Actual source code: asm.c

petsc-master 2020-09-23
Report Typos and Errors
  1: /*
  2:   This file defines an additive Schwarz preconditioner for any Mat implementation.

  4:   Note that each processor may have any number of subdomains. But in order to
  5:   deal easily with the VecScatter(), we treat each processor as if it has the
  6:   same number of subdomains.

  8:        n - total number of true subdomains on all processors
  9:        n_local_true - actual number of subdomains on this processor
 10:        n_local = maximum over all processors of n_local_true
 11: */

 13: #include <../src/ksp/pc/impls/asm/asm.h>

 15: static PetscErrorCode PCView_ASM(PC pc,PetscViewer viewer)
 16: {
 17:   PC_ASM         *osm = (PC_ASM*)pc->data;
 19:   PetscMPIInt    rank;
 20:   PetscInt       i,bsz;
 21:   PetscBool      iascii,isstring;
 22:   PetscViewer    sviewer;

 25:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
 26:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
 27:   if (iascii) {
 28:     char overlaps[256] = "user-defined overlap",blocks[256] = "total subdomain blocks not yet set";
 29:     if (osm->overlap >= 0) {PetscSNPrintf(overlaps,sizeof(overlaps),"amount of overlap = %D",osm->overlap);}
 30:     if (osm->n > 0) {PetscSNPrintf(blocks,sizeof(blocks),"total subdomain blocks = %D",osm->n);}
 31:     PetscViewerASCIIPrintf(viewer,"  %s, %s\n",blocks,overlaps);
 32:     PetscViewerASCIIPrintf(viewer,"  restriction/interpolation type - %s\n",PCASMTypes[osm->type]);
 33:     if (osm->dm_subdomains) {PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: using DM to define subdomains\n");}
 34:     if (osm->loctype != PC_COMPOSITE_ADDITIVE) {PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: local solve composition type - %s\n",PCCompositeTypes[osm->loctype]);}
 35:     MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
 36:     if (osm->same_local_solves) {
 37:       if (osm->ksp) {
 38:         PetscViewerASCIIPrintf(viewer,"  Local solver is the same for all blocks, as in the following KSP and PC objects on rank 0:\n");
 39:         PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
 40:         if (!rank) {
 41:           PetscViewerASCIIPushTab(viewer);
 42:           KSPView(osm->ksp[0],sviewer);
 43:           PetscViewerASCIIPopTab(viewer);
 44:         }
 45:         PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
 46:       }
 47:     } else {
 48:       PetscViewerASCIIPushSynchronized(viewer);
 49:       PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,osm->n_local_true);
 50:       PetscViewerFlush(viewer);
 51:       PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following KSP and PC objects:\n");
 52:       PetscViewerASCIIPushTab(viewer);
 53:       PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");
 54:       PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
 55:       for (i=0; i<osm->n_local_true; i++) {
 56:         ISGetLocalSize(osm->is[i],&bsz);
 57:         PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);
 58:         KSPView(osm->ksp[i],sviewer);
 59:         PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");
 60:       }
 61:       PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
 62:       PetscViewerASCIIPopTab(viewer);
 63:       PetscViewerFlush(viewer);
 64:       PetscViewerASCIIPopSynchronized(viewer);
 65:     }
 66:   } else if (isstring) {
 67:     PetscViewerStringSPrintf(viewer," blocks=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCASMTypes[osm->type]);
 68:     PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
 69:     if (osm->ksp) {KSPView(osm->ksp[0],sviewer);}
 70:     PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
 71:   }
 72:   return(0);
 73: }

 75: static PetscErrorCode PCASMPrintSubdomains(PC pc)
 76: {
 77:   PC_ASM         *osm = (PC_ASM*)pc->data;
 78:   const char     *prefix;
 79:   char           fname[PETSC_MAX_PATH_LEN+1];
 80:   PetscViewer    viewer, sviewer;
 81:   char           *s;
 82:   PetscInt       i,j,nidx;
 83:   const PetscInt *idx;
 84:   PetscMPIInt    rank, size;

 88:   MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);
 89:   MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);
 90:   PCGetOptionsPrefix(pc,&prefix);
 91:   PetscOptionsGetString(NULL,prefix,"-pc_asm_print_subdomains",fname,sizeof(fname),NULL);
 92:   if (fname[0] == 0) { PetscStrcpy(fname,"stdout"); };
 93:   PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);
 94:   for (i=0; i<osm->n_local; i++) {
 95:     if (i < osm->n_local_true) {
 96:       ISGetLocalSize(osm->is[i],&nidx);
 97:       ISGetIndices(osm->is[i],&idx);
 98:       /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
 99: #define len  16*(nidx+1)+512
100:       PetscMalloc1(len,&s);
101:       PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer);
102: #undef len
103:       PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D with overlap:\n", rank, size, i);
104:       for (j=0; j<nidx; j++) {
105:         PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);
106:       }
107:       ISRestoreIndices(osm->is[i],&idx);
108:       PetscViewerStringSPrintf(sviewer,"\n");
109:       PetscViewerDestroy(&sviewer);
110:       PetscViewerASCIIPushSynchronized(viewer);
111:       PetscViewerASCIISynchronizedPrintf(viewer, s);
112:       PetscViewerFlush(viewer);
113:       PetscViewerASCIIPopSynchronized(viewer);
114:       PetscFree(s);
115:       if (osm->is_local) {
116:         /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
117: #define len  16*(nidx+1)+512
118:         PetscMalloc1(len, &s);
119:         PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer);
120: #undef len
121:         PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D without overlap:\n", rank, size, i);
122:         ISGetLocalSize(osm->is_local[i],&nidx);
123:         ISGetIndices(osm->is_local[i],&idx);
124:         for (j=0; j<nidx; j++) {
125:           PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);
126:         }
127:         ISRestoreIndices(osm->is_local[i],&idx);
128:         PetscViewerStringSPrintf(sviewer,"\n");
129:         PetscViewerDestroy(&sviewer);
130:         PetscViewerASCIIPushSynchronized(viewer);
131:         PetscViewerASCIISynchronizedPrintf(viewer, s);
132:         PetscViewerFlush(viewer);
133:         PetscViewerASCIIPopSynchronized(viewer);
134:         PetscFree(s);
135:       }
136:     } else {
137:       /* Participate in collective viewer calls. */
138:       PetscViewerASCIIPushSynchronized(viewer);
139:       PetscViewerFlush(viewer);
140:       PetscViewerASCIIPopSynchronized(viewer);
141:       /* Assume either all ranks have is_local or none do. */
142:       if (osm->is_local) {
143:         PetscViewerASCIIPushSynchronized(viewer);
144:         PetscViewerFlush(viewer);
145:         PetscViewerASCIIPopSynchronized(viewer);
146:       }
147:     }
148:   }
149:   PetscViewerFlush(viewer);
150:   PetscViewerDestroy(&viewer);
151:   return(0);
152: }

154: static PetscErrorCode PCSetUp_ASM(PC pc)
155: {
156:   PC_ASM         *osm = (PC_ASM*)pc->data;
158:   PetscBool      flg;
159:   PetscInt       i,m,m_local;
160:   MatReuse       scall = MAT_REUSE_MATRIX;
161:   IS             isl;
162:   KSP            ksp;
163:   PC             subpc;
164:   const char     *prefix,*pprefix;
165:   Vec            vec;
166:   DM             *domain_dm = NULL;

169:   if (!pc->setupcalled) {
170:     PetscInt m;

172:     /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */
173:     if (osm->n_local_true == PETSC_DECIDE) {
174:       /* no subdomains given */
175:       /* try pc->dm first, if allowed */
176:       if (osm->dm_subdomains && pc->dm) {
177:         PetscInt  num_domains, d;
178:         char      **domain_names;
179:         IS        *inner_domain_is, *outer_domain_is;
180:         DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm);
181:         osm->overlap = -1; /* We do not want to increase the overlap of the IS.
182:                               A future improvement of this code might allow one to use
183:                               DM-defined subdomains and also increase the overlap,
184:                               but that is not currently supported */
185:         if (num_domains) {
186:           PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is);
187:         }
188:         for (d = 0; d < num_domains; ++d) {
189:           if (domain_names)    {PetscFree(domain_names[d]);}
190:           if (inner_domain_is) {ISDestroy(&inner_domain_is[d]);}
191:           if (outer_domain_is) {ISDestroy(&outer_domain_is[d]);}
192:         }
193:         PetscFree(domain_names);
194:         PetscFree(inner_domain_is);
195:         PetscFree(outer_domain_is);
196:       }
197:       if (osm->n_local_true == PETSC_DECIDE) {
198:         /* still no subdomains; use one subdomain per processor */
199:         osm->n_local_true = 1;
200:       }
201:     }
202:     { /* determine the global and max number of subdomains */
203:       struct {PetscInt max,sum;} inwork,outwork;
204:       PetscMPIInt size;

206:       inwork.max   = osm->n_local_true;
207:       inwork.sum   = osm->n_local_true;
208:       MPIU_Allreduce(&inwork,&outwork,1,MPIU_2INT,MPIU_MAXSUM_OP,PetscObjectComm((PetscObject)pc));
209:       osm->n_local = outwork.max;
210:       osm->n       = outwork.sum;

212:       MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
213:       if (outwork.max == 1 && outwork.sum == size) {
214:         /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */
215:         MatSetOption(pc->pmat,MAT_SUBMAT_SINGLEIS,PETSC_TRUE);
216:       }
217:     }
218:     if (!osm->is) { /* create the index sets */
219:       PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is);
220:     }
221:     if (osm->n_local_true > 1 && !osm->is_local) {
222:       PetscMalloc1(osm->n_local_true,&osm->is_local);
223:       for (i=0; i<osm->n_local_true; i++) {
224:         if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
225:           ISDuplicate(osm->is[i],&osm->is_local[i]);
226:           ISCopy(osm->is[i],osm->is_local[i]);
227:         } else {
228:           PetscObjectReference((PetscObject)osm->is[i]);
229:           osm->is_local[i] = osm->is[i];
230:         }
231:       }
232:     }
233:     PCGetOptionsPrefix(pc,&prefix);
234:     flg  = PETSC_FALSE;
235:     PetscOptionsGetBool(NULL,prefix,"-pc_asm_print_subdomains",&flg,NULL);
236:     if (flg) { PCASMPrintSubdomains(pc); }

238:     if (osm->overlap > 0) {
239:       /* Extend the "overlapping" regions by a number of steps */
240:       MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap);
241:     }
242:     if (osm->sort_indices) {
243:       for (i=0; i<osm->n_local_true; i++) {
244:         ISSort(osm->is[i]);
245:         if (osm->is_local) {
246:           ISSort(osm->is_local[i]);
247:         }
248:       }
249:     }

251:     if (!osm->ksp) {
252:       /* Create the local solvers */
253:       PetscMalloc1(osm->n_local_true,&osm->ksp);
254:       if (domain_dm) {
255:         PetscInfo(pc,"Setting up ASM subproblems using the embedded DM\n");
256:       }
257:       for (i=0; i<osm->n_local_true; i++) {
258:         KSPCreate(PETSC_COMM_SELF,&ksp);
259:         KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);
260:         PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);
261:         PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);
262:         KSPSetType(ksp,KSPPREONLY);
263:         KSPGetPC(ksp,&subpc);
264:         PCGetOptionsPrefix(pc,&prefix);
265:         KSPSetOptionsPrefix(ksp,prefix);
266:         KSPAppendOptionsPrefix(ksp,"sub_");
267:         if (domain_dm) {
268:           KSPSetDM(ksp, domain_dm[i]);
269:           KSPSetDMActive(ksp, PETSC_FALSE);
270:           DMDestroy(&domain_dm[i]);
271:         }
272:         osm->ksp[i] = ksp;
273:       }
274:       if (domain_dm) {
275:         PetscFree(domain_dm);
276:       }
277:     }

279:     ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis);
280:     ISSortRemoveDups(osm->lis);
281:     ISGetLocalSize(osm->lis, &m);
282:     VecCreateSeq(PETSC_COMM_SELF, m, &osm->lx);
283:     VecDuplicate(osm->lx, &osm->ly);

285:     scall = MAT_INITIAL_MATRIX;
286:   } else {
287:     /*
288:        Destroy the blocks from the previous iteration
289:     */
290:     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
291:       MatDestroyMatrices(osm->n_local_true,&osm->pmat);
292:       scall = MAT_INITIAL_MATRIX;
293:     }
294:   }

296:   /*
297:      Extract out the submatrices
298:   */
299:   MatCreateSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);
300:   if (scall == MAT_INITIAL_MATRIX) {
301:     PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);
302:     for (i=0; i<osm->n_local_true; i++) {
303:       PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);
304:       PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);
305:     }
306:   }

308:   /* Convert the types of the submatrices (if needbe) */
309:   if (osm->sub_mat_type) {
310:     for (i=0; i<osm->n_local_true; i++) {
311:       MatConvert(osm->pmat[i],osm->sub_mat_type,MAT_INPLACE_MATRIX,&(osm->pmat[i]));
312:     }
313:   }

315:   if (!pc->setupcalled) {
316:     /* Create the local work vectors (from the local matrices) and scatter contexts */
317:     MatCreateVecs(pc->pmat,&vec,NULL);

319:     if (osm->is_local && (osm->type == PC_ASM_INTERPOLATE || osm->type == PC_ASM_NONE)) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot use interpolate or none PCASMType if is_local was provided to PCASMSetLocalSubdomains()");
320:     if (osm->is_local && osm->type == PC_ASM_RESTRICT && osm->loctype == PC_COMPOSITE_ADDITIVE) {
321:       PetscMalloc1(osm->n_local_true,&osm->lprolongation);
322:     }
323:     PetscMalloc1(osm->n_local_true,&osm->lrestriction);
324:     PetscMalloc1(osm->n_local_true,&osm->x);
325:     PetscMalloc1(osm->n_local_true,&osm->y);

327:     ISGetLocalSize(osm->lis,&m);
328:     ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);
329:     VecScatterCreate(vec,osm->lis,osm->lx,isl,&osm->restriction);
330:     ISDestroy(&isl);


333:     for (i=0; i<osm->n_local_true; ++i) {
334:       ISLocalToGlobalMapping ltog;
335:       IS                     isll;
336:       const PetscInt         *idx_is;
337:       PetscInt               *idx_lis,nout;

339:       ISGetLocalSize(osm->is[i],&m);
340:       MatCreateVecs(osm->pmat[i],&osm->x[i],NULL);
341:       VecDuplicate(osm->x[i],&osm->y[i]);

343:       /* generate a scatter from ly to y[i] picking all the overlapping is[i] entries */
344:       ISLocalToGlobalMappingCreateIS(osm->lis,&ltog);
345:       ISGetLocalSize(osm->is[i],&m);
346:       ISGetIndices(osm->is[i], &idx_is);
347:       PetscMalloc1(m,&idx_lis);
348:       ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m,idx_is,&nout,idx_lis);
349:       if (nout != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is not a subset of lis");
350:       ISRestoreIndices(osm->is[i], &idx_is);
351:       ISCreateGeneral(PETSC_COMM_SELF,m,idx_lis,PETSC_OWN_POINTER,&isll);
352:       ISLocalToGlobalMappingDestroy(&ltog);
353:       ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);
354:       VecScatterCreate(osm->ly,isll,osm->y[i],isl,&osm->lrestriction[i]);
355:       ISDestroy(&isll);
356:       ISDestroy(&isl);
357:       if (osm->lprolongation) { /* generate a scatter from y[i] to ly picking only the the non-overlapping is_local[i] entries */
358:         ISLocalToGlobalMapping ltog;
359:         IS                     isll,isll_local;
360:         const PetscInt         *idx_local;
361:         PetscInt               *idx1, *idx2, nout;

363:         ISGetLocalSize(osm->is_local[i],&m_local);
364:         ISGetIndices(osm->is_local[i], &idx_local);

366:         ISLocalToGlobalMappingCreateIS(osm->is[i],&ltog);
367:         PetscMalloc1(m_local,&idx1);
368:         ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx1);
369:         ISLocalToGlobalMappingDestroy(&ltog);
370:         if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is");
371:         ISCreateGeneral(PETSC_COMM_SELF,m_local,idx1,PETSC_OWN_POINTER,&isll);

373:         ISLocalToGlobalMappingCreateIS(osm->lis,&ltog);
374:         PetscMalloc1(m_local,&idx2);
375:         ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx2);
376:         ISLocalToGlobalMappingDestroy(&ltog);
377:         if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of lis");
378:         ISCreateGeneral(PETSC_COMM_SELF,m_local,idx2,PETSC_OWN_POINTER,&isll_local);

380:         ISRestoreIndices(osm->is_local[i], &idx_local);
381:         VecScatterCreate(osm->y[i],isll,osm->ly,isll_local,&osm->lprolongation[i]);

383:         ISDestroy(&isll);
384:         ISDestroy(&isll_local);
385:       }
386:     }
387:     VecDestroy(&vec);
388:   }

390:   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
391:     IS      *cis;
392:     PetscInt c;

394:     PetscMalloc1(osm->n_local_true, &cis);
395:     for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis;
396:     MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats);
397:     PetscFree(cis);
398:   }

400:   /* Return control to the user so that the submatrices can be modified (e.g., to apply
401:      different boundary conditions for the submatrices than for the global problem) */
402:   PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);

404:   /*
405:      Loop over subdomains putting them into local ksp
406:   */
407:   for (i=0; i<osm->n_local_true; i++) {
408:     KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);
409:     if (!pc->setupcalled) {
410:       KSPSetFromOptions(osm->ksp[i]);
411:     }
412:   }
413:   return(0);
414: }

416: static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc)
417: {
418:   PC_ASM             *osm = (PC_ASM*)pc->data;
419:   PetscErrorCode     ierr;
420:   PetscInt           i;
421:   KSPConvergedReason reason;

424:   for (i=0; i<osm->n_local_true; i++) {
425:     KSPSetUp(osm->ksp[i]);
426:     KSPGetConvergedReason(osm->ksp[i],&reason);
427:     if (reason == KSP_DIVERGED_PC_FAILED) {
428:       pc->failedreason = PC_SUBPC_ERROR;
429:     }
430:   }
431:   return(0);
432: }

434: static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y)
435: {
436:   PC_ASM         *osm = (PC_ASM*)pc->data;
438:   PetscInt       i,n_local_true = osm->n_local_true;
439:   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;

442:   /*
443:      support for limiting the restriction or interpolation to only local
444:      subdomain values (leaving the other values 0).
445:   */
446:   if (!(osm->type & PC_ASM_RESTRICT)) {
447:     forward = SCATTER_FORWARD_LOCAL;
448:     /* have to zero the work RHS since scatter may leave some slots empty */
449:     VecSet(osm->lx, 0.0);
450:   }
451:   if (!(osm->type & PC_ASM_INTERPOLATE)) {
452:     reverse = SCATTER_REVERSE_LOCAL;
453:   }

455:   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) {
456:     /* zero the global and the local solutions */
457:     VecSet(y, 0.0);
458:     VecSet(osm->ly, 0.0);

460:     /* copy the global RHS to local RHS including the ghost nodes */
461:     VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);
462:     VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);

464:     /* restrict local RHS to the overlapping 0-block RHS */
465:     VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);
466:     VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);

468:     /* do the local solves */
469:     for (i = 0; i < n_local_true; ++i) {

471:       /* solve the overlapping i-block */
472:       PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i],0);
473:       KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);
474:       KSPCheckSolve(osm->ksp[i], pc, osm->y[i]);
475:       PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0);

477:       if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution (only for restrictive additive) */
478:         VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);
479:         VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);
480:       } else { /* interpolate the overlapping i-block solution to the local solution */
481:         VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);
482:         VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);
483:       }

485:       if (i < n_local_true-1) {
486:         /* restrict local RHS to the overlapping (i+1)-block RHS */
487:         VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);
488:         VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);

490:         if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
491:           /* update the overlapping (i+1)-block RHS using the current local solution */
492:           MatMult(osm->lmats[i+1], osm->ly, osm->y[i+1]);
493:           VecAXPBY(osm->x[i+1],-1.,1., osm->y[i+1]);
494:         }
495:       }
496:     }
497:     /* add the local solution to the global solution including the ghost nodes */
498:     VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse);
499:     VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse);
500:   } else SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
501:   return(0);
502: }

504: static PetscErrorCode PCMatApply_ASM(PC pc,Mat X,Mat Y)
505: {
506:   PC_ASM         *osm = (PC_ASM*)pc->data;
507:   Mat            Z,W;
508:   Vec            x;
509:   PetscInt       i,m,N;
510:   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;

514:   if (osm->n_local_true > 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not yet implemented");
515:   /*
516:      support for limiting the restriction or interpolation to only local
517:      subdomain values (leaving the other values 0).
518:   */
519:   if (!(osm->type & PC_ASM_RESTRICT)) {
520:     forward = SCATTER_FORWARD_LOCAL;
521:     /* have to zero the work RHS since scatter may leave some slots empty */
522:     VecSet(osm->lx, 0.0);
523:   }
524:   if (!(osm->type & PC_ASM_INTERPOLATE)) {
525:     reverse = SCATTER_REVERSE_LOCAL;
526:   }
527:   VecGetLocalSize(osm->x[0], &m);
528:   MatGetSize(X, NULL, &N);
529:   MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &Z);
530:   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) {
531:     /* zero the global and the local solutions */
532:     MatZeroEntries(Y);
533:     VecSet(osm->ly, 0.0);

535:     for (i = 0; i < N; ++i) {
536:       MatDenseGetColumnVecRead(X, i, &x);
537:       /* copy the global RHS to local RHS including the ghost nodes */
538:       VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);
539:       VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);
540:       MatDenseRestoreColumnVecRead(X, i, &x);

542:       MatDenseGetColumnVecWrite(Z, i, &x);
543:       /* restrict local RHS to the overlapping 0-block RHS */
544:       VecScatterBegin(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward);
545:       VecScatterEnd(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward);
546:       MatDenseRestoreColumnVecWrite(Z, i, &x);
547:     }
548:     MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &W);
549:     /* solve the overlapping 0-block */
550:     PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[0], Z, W, 0);
551:     KSPMatSolve(osm->ksp[0], Z, W);
552:     KSPCheckSolve(osm->ksp[0], pc, NULL);
553:     PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[0], Z, W,0);
554:     MatDestroy(&Z);

556:     for (i = 0; i < N; ++i) {
557:       VecSet(osm->ly, 0.0);
558:       MatDenseGetColumnVecRead(W, i, &x);
559:       if (osm->lprolongation) { /* interpolate the non-overlapping 0-block solution to the local solution (only for restrictive additive) */
560:         VecScatterBegin(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward);
561:         VecScatterEnd(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward);
562:       } else { /* interpolate the overlapping 0-block solution to the local solution */
563:         VecScatterBegin(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse);
564:         VecScatterEnd(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse);
565:       }
566:       MatDenseRestoreColumnVecRead(W, i, &x);

568:       MatDenseGetColumnVecWrite(Y, i, &x);
569:       /* add the local solution to the global solution including the ghost nodes */
570:       VecScatterBegin(osm->restriction, osm->ly, x, ADD_VALUES, reverse);
571:       VecScatterEnd(osm->restriction, osm->ly, x, ADD_VALUES, reverse);
572:       MatDenseRestoreColumnVecWrite(Y, i, &x);
573:     }
574:     MatDestroy(&W);
575:   } else SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
576:   return(0);
577: }

579: static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y)
580: {
581:   PC_ASM         *osm = (PC_ASM*)pc->data;
583:   PetscInt       i,n_local_true = osm->n_local_true;
584:   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;

587:   /*
588:      Support for limiting the restriction or interpolation to only local
589:      subdomain values (leaving the other values 0).

591:      Note: these are reversed from the PCApply_ASM() because we are applying the
592:      transpose of the three terms
593:   */

595:   if (!(osm->type & PC_ASM_INTERPOLATE)) {
596:     forward = SCATTER_FORWARD_LOCAL;
597:     /* have to zero the work RHS since scatter may leave some slots empty */
598:     VecSet(osm->lx, 0.0);
599:   }
600:   if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL;

602:   /* zero the global and the local solutions */
603:   VecSet(y, 0.0);
604:   VecSet(osm->ly, 0.0);

606:   /* Copy the global RHS to local RHS including the ghost nodes */
607:   VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);
608:   VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);

610:   /* Restrict local RHS to the overlapping 0-block RHS */
611:   VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);
612:   VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);

614:   /* do the local solves */
615:   for (i = 0; i < n_local_true; ++i) {

617:     /* solve the overlapping i-block */
618:     PetscLogEventBegin(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);
619:     KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i]);
620:     KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);
621:     PetscLogEventEnd(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);

623:     if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution */
624:       VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);
625:       VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);
626:     } else { /* interpolate the overlapping i-block solution to the local solution */
627:       VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);
628:       VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);
629:     }

631:     if (i < n_local_true-1) {
632:       /* Restrict local RHS to the overlapping (i+1)-block RHS */
633:       VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);
634:       VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);
635:     }
636:   }
637:   /* Add the local solution to the global solution including the ghost nodes */
638:   VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse);
639:   VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse);

641:   return(0);

643: }

645: static PetscErrorCode PCReset_ASM(PC pc)
646: {
647:   PC_ASM         *osm = (PC_ASM*)pc->data;
649:   PetscInt       i;

652:   if (osm->ksp) {
653:     for (i=0; i<osm->n_local_true; i++) {
654:       KSPReset(osm->ksp[i]);
655:     }
656:   }
657:   if (osm->pmat) {
658:     if (osm->n_local_true > 0) {
659:       MatDestroySubMatrices(osm->n_local_true,&osm->pmat);
660:     }
661:   }
662:   if (osm->lrestriction) {
663:     VecScatterDestroy(&osm->restriction);
664:     for (i=0; i<osm->n_local_true; i++) {
665:       VecScatterDestroy(&osm->lrestriction[i]);
666:       if (osm->lprolongation) {VecScatterDestroy(&osm->lprolongation[i]);}
667:       VecDestroy(&osm->x[i]);
668:       VecDestroy(&osm->y[i]);
669:     }
670:     PetscFree(osm->lrestriction);
671:     if (osm->lprolongation) {PetscFree(osm->lprolongation);}
672:     PetscFree(osm->x);
673:     PetscFree(osm->y);

675:   }
676:   PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);
677:   ISDestroy(&osm->lis);
678:   VecDestroy(&osm->lx);
679:   VecDestroy(&osm->ly);
680:   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
681:     MatDestroyMatrices(osm->n_local_true, &osm->lmats);
682:   }

684:   PetscFree(osm->sub_mat_type);

686:   osm->is       = NULL;
687:   osm->is_local = NULL;
688:   return(0);
689: }

691: static PetscErrorCode PCDestroy_ASM(PC pc)
692: {
693:   PC_ASM         *osm = (PC_ASM*)pc->data;
695:   PetscInt       i;

698:   PCReset_ASM(pc);
699:   if (osm->ksp) {
700:     for (i=0; i<osm->n_local_true; i++) {
701:       KSPDestroy(&osm->ksp[i]);
702:     }
703:     PetscFree(osm->ksp);
704:   }
705:   PetscFree(pc->data);
706:   return(0);
707: }

709: static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc)
710: {
711:   PC_ASM         *osm = (PC_ASM*)pc->data;
713:   PetscInt       blocks,ovl;
714:   PetscBool      flg;
715:   PCASMType      asmtype;
716:   PCCompositeType loctype;
717:   char           sub_mat_type[256];

720:   PetscOptionsHead(PetscOptionsObject,"Additive Schwarz options");
721:   PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);
722:   PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);
723:   if (flg) {
724:     PCASMSetTotalSubdomains(pc,blocks,NULL,NULL);
725:     osm->dm_subdomains = PETSC_FALSE;
726:   }
727:   PetscOptionsInt("-pc_asm_local_blocks","Number of local subdomains","PCASMSetLocalSubdomains",osm->n_local_true,&blocks,&flg);
728:   if (flg) {
729:     PCASMSetLocalSubdomains(pc,blocks,NULL,NULL);
730:     osm->dm_subdomains = PETSC_FALSE;
731:   }
732:   PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);
733:   if (flg) {
734:     PCASMSetOverlap(pc,ovl);
735:     osm->dm_subdomains = PETSC_FALSE;
736:   }
737:   flg  = PETSC_FALSE;
738:   PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);
739:   if (flg) {PCASMSetType(pc,asmtype); }
740:   flg  = PETSC_FALSE;
741:   PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg);
742:   if (flg) {PCASMSetLocalType(pc,loctype); }
743:   PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg);
744:   if (flg) {
745:     PCASMSetSubMatType(pc,sub_mat_type);
746:   }
747:   PetscOptionsTail();
748:   return(0);
749: }

751: /*------------------------------------------------------------------------------------*/

753: static PetscErrorCode  PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[])
754: {
755:   PC_ASM         *osm = (PC_ASM*)pc->data;
757:   PetscInt       i;

760:   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n);
761:   if (pc->setupcalled && (n != osm->n_local_true || is)) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetLocalSubdomains() should be called before calling PCSetUp().");

763:   if (!pc->setupcalled) {
764:     if (is) {
765:       for (i=0; i<n; i++) {PetscObjectReference((PetscObject)is[i]);}
766:     }
767:     if (is_local) {
768:       for (i=0; i<n; i++) {PetscObjectReference((PetscObject)is_local[i]);}
769:     }
770:     PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);

772:     osm->n_local_true = n;
773:     osm->is           = NULL;
774:     osm->is_local     = NULL;
775:     if (is) {
776:       PetscMalloc1(n,&osm->is);
777:       for (i=0; i<n; i++) osm->is[i] = is[i];
778:       /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */
779:       osm->overlap = -1;
780:     }
781:     if (is_local) {
782:       PetscMalloc1(n,&osm->is_local);
783:       for (i=0; i<n; i++) osm->is_local[i] = is_local[i];
784:       if (!is) {
785:         PetscMalloc1(osm->n_local_true,&osm->is);
786:         for (i=0; i<osm->n_local_true; i++) {
787:           if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
788:             ISDuplicate(osm->is_local[i],&osm->is[i]);
789:             ISCopy(osm->is_local[i],osm->is[i]);
790:           } else {
791:             PetscObjectReference((PetscObject)osm->is_local[i]);
792:             osm->is[i] = osm->is_local[i];
793:           }
794:         }
795:       }
796:     }
797:   }
798:   return(0);
799: }

801: static PetscErrorCode  PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local)
802: {
803:   PC_ASM         *osm = (PC_ASM*)pc->data;
805:   PetscMPIInt    rank,size;
806:   PetscInt       n;

809:   if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N);
810:   if (is || is_local) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Use PCASMSetLocalSubdomains() to set specific index sets\n\they cannot be set globally yet.");

812:   /*
813:      Split the subdomains equally among all processors
814:   */
815:   MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
816:   MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
817:   n    = N/size + ((N % size) > rank);
818:   if (!n) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Process %d must have at least one block: total processors %d total blocks %D",(int)rank,(int)size,N);
819:   if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp().");
820:   if (!pc->setupcalled) {
821:     PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);

823:     osm->n_local_true = n;
824:     osm->is           = NULL;
825:     osm->is_local     = NULL;
826:   }
827:   return(0);
828: }

830: static PetscErrorCode  PCASMSetOverlap_ASM(PC pc,PetscInt ovl)
831: {
832:   PC_ASM *osm = (PC_ASM*)pc->data;

835:   if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
836:   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp().");
837:   if (!pc->setupcalled) osm->overlap = ovl;
838:   return(0);
839: }

841: static PetscErrorCode  PCASMSetType_ASM(PC pc,PCASMType type)
842: {
843:   PC_ASM *osm = (PC_ASM*)pc->data;

846:   osm->type     = type;
847:   osm->type_set = PETSC_TRUE;
848:   return(0);
849: }

851: static PetscErrorCode  PCASMGetType_ASM(PC pc,PCASMType *type)
852: {
853:   PC_ASM *osm = (PC_ASM*)pc->data;

856:   *type = osm->type;
857:   return(0);
858: }

860: static PetscErrorCode  PCASMSetLocalType_ASM(PC pc, PCCompositeType type)
861: {
862:   PC_ASM *osm = (PC_ASM *) pc->data;

865:   if (type != PC_COMPOSITE_ADDITIVE && type != PC_COMPOSITE_MULTIPLICATIVE) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Only supports additive or multiplicative as the local type");
866:   osm->loctype = type;
867:   return(0);
868: }

870: static PetscErrorCode  PCASMGetLocalType_ASM(PC pc, PCCompositeType *type)
871: {
872:   PC_ASM *osm = (PC_ASM *) pc->data;

875:   *type = osm->loctype;
876:   return(0);
877: }

879: static PetscErrorCode  PCASMSetSortIndices_ASM(PC pc,PetscBool  doSort)
880: {
881:   PC_ASM *osm = (PC_ASM*)pc->data;

884:   osm->sort_indices = doSort;
885:   return(0);
886: }

888: static PetscErrorCode  PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
889: {
890:   PC_ASM         *osm = (PC_ASM*)pc->data;

894:   if (osm->n_local_true < 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Need to call PCSetUp() on PC (or KSPSetUp() on the outer KSP object) before calling here");

896:   if (n_local) *n_local = osm->n_local_true;
897:   if (first_local) {
898:     MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));
899:     *first_local -= osm->n_local_true;
900:   }
901:   if (ksp) {
902:     /* Assume that local solves are now different; not necessarily
903:        true though!  This flag is used only for PCView_ASM() */
904:     *ksp                   = osm->ksp;
905:     osm->same_local_solves = PETSC_FALSE;
906:   }
907:   return(0);
908: }

910: static PetscErrorCode  PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type)
911: {
912:   PC_ASM         *osm = (PC_ASM*)pc->data;

917:   *sub_mat_type = osm->sub_mat_type;
918:   return(0);
919: }

921: static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type)
922: {
923:   PetscErrorCode    ierr;
924:   PC_ASM            *osm = (PC_ASM*)pc->data;

928:   PetscFree(osm->sub_mat_type);
929:   PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type);
930:   return(0);
931: }

933: /*@C
934:     PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner.

936:     Collective on pc

938:     Input Parameters:
939: +   pc - the preconditioner context
940: .   n - the number of subdomains for this processor (default value = 1)
941: .   is - the index set that defines the subdomains for this processor
942:          (or NULL for PETSc to determine subdomains)
943: -   is_local - the index sets that define the local part of the subdomains for this processor, not used unless PCASMType is PC_ASM_RESTRICT
944:          (or NULL to not provide these)

946:     Options Database Key:
947:     To set the total number of subdomain blocks rather than specify the
948:     index sets, use the option
949: .    -pc_asm_local_blocks <blks> - Sets local blocks

951:     Notes:
952:     The IS numbering is in the parallel, global numbering of the vector for both is and is_local

954:     By default the ASM preconditioner uses 1 block per processor.

956:     Use PCASMSetTotalSubdomains() to set the subdomains for all processors.

958:     If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated
959:     back to form the global solution (this is the standard restricted additive Schwarz method)

961:     If the is_local is provided and PCASMType is PC_ASM_INTERPOLATE or PC_ASM_NONE then an error is generated since there is
962:     no code to handle that case.

964:     Level: advanced

966: .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
967:           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), PCASMType, PCASMSetType()
968: @*/
969: PetscErrorCode  PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
970: {

975:   PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));
976:   return(0);
977: }

979: /*@C
980:     PCASMSetTotalSubdomains - Sets the subdomains for all processors for the
981:     additive Schwarz preconditioner.  Either all or no processors in the
982:     PC communicator must call this routine, with the same index sets.

984:     Collective on pc

986:     Input Parameters:
987: +   pc - the preconditioner context
988: .   N  - the number of subdomains for all processors
989: .   is - the index sets that define the subdomains for all processors
990:          (or NULL to ask PETSc to determine the subdomains)
991: -   is_local - the index sets that define the local part of the subdomains for this processor
992:          (or NULL to not provide this information)

994:     Options Database Key:
995:     To set the total number of subdomain blocks rather than specify the
996:     index sets, use the option
997: .    -pc_asm_blocks <blks> - Sets total blocks

999:     Notes:
1000:     Currently you cannot use this to set the actual subdomains with the argument is or is_local.

1002:     By default the ASM preconditioner uses 1 block per processor.

1004:     These index sets cannot be destroyed until after completion of the
1005:     linear solves for which the ASM preconditioner is being used.

1007:     Use PCASMSetLocalSubdomains() to set local subdomains.

1009:     The IS numbering is in the parallel, global numbering of the vector for both is and is_local

1011:     Level: advanced

1013: .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
1014:           PCASMCreateSubdomains2D()
1015: @*/
1016: PetscErrorCode  PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[])
1017: {

1022:   PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local));
1023:   return(0);
1024: }

1026: /*@
1027:     PCASMSetOverlap - Sets the overlap between a pair of subdomains for the
1028:     additive Schwarz preconditioner.  Either all or no processors in the
1029:     PC communicator must call this routine.

1031:     Logically Collective on pc

1033:     Input Parameters:
1034: +   pc  - the preconditioner context
1035: -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)

1037:     Options Database Key:
1038: .   -pc_asm_overlap <ovl> - Sets overlap

1040:     Notes:
1041:     By default the ASM preconditioner uses 1 block per processor.  To use
1042:     multiple blocks per perocessor, see PCASMSetTotalSubdomains() and
1043:     PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>).

1045:     The overlap defaults to 1, so if one desires that no additional
1046:     overlap be computed beyond what may have been set with a call to
1047:     PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl
1048:     must be set to be 0.  In particular, if one does not explicitly set
1049:     the subdomains an application code, then all overlap would be computed
1050:     internally by PETSc, and using an overlap of 0 would result in an ASM
1051:     variant that is equivalent to the block Jacobi preconditioner.

1053:     The default algorithm used by PETSc to increase overlap is fast, but not scalable,
1054:     use the option -mat_increase_overlap_scalable when the problem and number of processes is large.

1056:     Note that one can define initial index sets with any overlap via
1057:     PCASMSetLocalSubdomains(); the routine
1058:     PCASMSetOverlap() merely allows PETSc to extend that overlap further
1059:     if desired.

1061:     Level: intermediate

1063: .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
1064:           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), MatIncreaseOverlap()
1065: @*/
1066: PetscErrorCode  PCASMSetOverlap(PC pc,PetscInt ovl)
1067: {

1073:   PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl));
1074:   return(0);
1075: }

1077: /*@
1078:     PCASMSetType - Sets the type of restriction and interpolation used
1079:     for local problems in the additive Schwarz method.

1081:     Logically Collective on pc

1083:     Input Parameters:
1084: +   pc  - the preconditioner context
1085: -   type - variant of ASM, one of
1086: .vb
1087:       PC_ASM_BASIC       - full interpolation and restriction
1088:       PC_ASM_RESTRICT    - full restriction, local processor interpolation
1089:       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
1090:       PC_ASM_NONE        - local processor restriction and interpolation
1091: .ve

1093:     Options Database Key:
1094: .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type

1096:     Notes:
1097:     if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected
1098:     to limit the local processor interpolation

1100:     Level: intermediate

1102: .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1103:           PCASMCreateSubdomains2D(), PCASMType, PCASMSetLocalType(), PCASMGetLocalType()
1104: @*/
1105: PetscErrorCode  PCASMSetType(PC pc,PCASMType type)
1106: {

1112:   PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type));
1113:   return(0);
1114: }

1116: /*@
1117:     PCASMGetType - Gets the type of restriction and interpolation used
1118:     for local problems in the additive Schwarz method.

1120:     Logically Collective on pc

1122:     Input Parameter:
1123: .   pc  - the preconditioner context

1125:     Output Parameter:
1126: .   type - variant of ASM, one of

1128: .vb
1129:       PC_ASM_BASIC       - full interpolation and restriction
1130:       PC_ASM_RESTRICT    - full restriction, local processor interpolation
1131:       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
1132:       PC_ASM_NONE        - local processor restriction and interpolation
1133: .ve

1135:     Options Database Key:
1136: .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type

1138:     Level: intermediate

1140: .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1141:           PCASMCreateSubdomains2D(), PCASMType, PCASMSetType(), PCASMSetLocalType(), PCASMGetLocalType()
1142: @*/
1143: PetscErrorCode  PCASMGetType(PC pc,PCASMType *type)
1144: {

1149:   PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type));
1150:   return(0);
1151: }

1153: /*@
1154:   PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method.

1156:   Logically Collective on pc

1158:   Input Parameters:
1159: + pc  - the preconditioner context
1160: - type - type of composition, one of
1161: .vb
1162:   PC_COMPOSITE_ADDITIVE       - local additive combination
1163:   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
1164: .ve

1166:   Options Database Key:
1167: . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type

1169:   Level: intermediate

1171: .seealso: PCASMSetType(), PCASMGetType(), PCASMGetLocalType(), PCASM, PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
1172: @*/
1173: PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type)
1174: {

1180:   PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));
1181:   return(0);
1182: }

1184: /*@
1185:   PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method.

1187:   Logically Collective on pc

1189:   Input Parameter:
1190: . pc  - the preconditioner context

1192:   Output Parameter:
1193: . type - type of composition, one of
1194: .vb
1195:   PC_COMPOSITE_ADDITIVE       - local additive combination
1196:   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
1197: .ve

1199:   Options Database Key:
1200: . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type

1202:   Level: intermediate

1204: .seealso: PCASMSetType(), PCASMGetType(), PCASMSetLocalType(), PCASMCreate(), PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
1205: @*/
1206: PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type)
1207: {

1213:   PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));
1214:   return(0);
1215: }

1217: /*@
1218:     PCASMSetSortIndices - Determines whether subdomain indices are sorted.

1220:     Logically Collective on pc

1222:     Input Parameters:
1223: +   pc  - the preconditioner context
1224: -   doSort - sort the subdomain indices

1226:     Level: intermediate

1228: .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1229:           PCASMCreateSubdomains2D()
1230: @*/
1231: PetscErrorCode  PCASMSetSortIndices(PC pc,PetscBool doSort)
1232: {

1238:   PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));
1239:   return(0);
1240: }

1242: /*@C
1243:    PCASMGetSubKSP - Gets the local KSP contexts for all blocks on
1244:    this processor.

1246:    Collective on pc iff first_local is requested

1248:    Input Parameter:
1249: .  pc - the preconditioner context

1251:    Output Parameters:
1252: +  n_local - the number of blocks on this processor or NULL
1253: .  first_local - the global number of the first block on this processor or NULL,
1254:                  all processors must request or all must pass NULL
1255: -  ksp - the array of KSP contexts

1257:    Note:
1258:    After PCASMGetSubKSP() the array of KSPes is not to be freed.

1260:    You must call KSPSetUp() before calling PCASMGetSubKSP().

1262:    Fortran note:
1263:    The output argument 'ksp' must be an array of sufficient length or PETSC_NULL_KSP. The latter can be used to learn the necessary length.

1265:    Level: advanced

1267: .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(),
1268:           PCASMCreateSubdomains2D(),
1269: @*/
1270: PetscErrorCode  PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
1271: {

1276:   PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));
1277:   return(0);
1278: }

1280: /* -------------------------------------------------------------------------------------*/
1281: /*MC
1282:    PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
1283:            its own KSP object.

1285:    Options Database Keys:
1286: +  -pc_asm_blocks <blks> - Sets total blocks
1287: .  -pc_asm_overlap <ovl> - Sets overlap
1288: .  -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict
1289: -  -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive

1291:      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
1292:       will get a different convergence rate due to the default option of -pc_asm_type restrict. Use
1293:       -pc_asm_type basic to use the standard ASM.

1295:    Notes:
1296:     Each processor can have one or more blocks, but a block cannot be shared by more
1297:      than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor.

1299:      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1300:         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly

1302:      To set the options on the solvers separate for each block call PCASMGetSubKSP()
1303:          and set the options directly on the resulting KSP object (you can access its PC
1304:          with KSPGetPC())

1306:    Level: beginner

1308:     References:
1309: +   1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions
1310:      Courant Institute, New York University Technical report
1311: -   2. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
1312:     Cambridge University Press.

1314: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1315:            PCBJACOBI, PCASMGetSubKSP(), PCASMSetLocalSubdomains(), PCASMType, PCASMGetType(), PCASMSetLocalType(), PCASMGetLocalType()
1316:            PCASMSetTotalSubdomains(), PCSetModifySubMatrices(), PCASMSetOverlap(), PCASMSetType(), PCCompositeType

1318: M*/

1320: PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc)
1321: {
1323:   PC_ASM         *osm;

1326:   PetscNewLog(pc,&osm);

1328:   osm->n                 = PETSC_DECIDE;
1329:   osm->n_local           = 0;
1330:   osm->n_local_true      = PETSC_DECIDE;
1331:   osm->overlap           = 1;
1332:   osm->ksp               = NULL;
1333:   osm->restriction       = NULL;
1334:   osm->lprolongation     = NULL;
1335:   osm->lrestriction      = NULL;
1336:   osm->x                 = NULL;
1337:   osm->y                 = NULL;
1338:   osm->is                = NULL;
1339:   osm->is_local          = NULL;
1340:   osm->mat               = NULL;
1341:   osm->pmat              = NULL;
1342:   osm->type              = PC_ASM_RESTRICT;
1343:   osm->loctype           = PC_COMPOSITE_ADDITIVE;
1344:   osm->same_local_solves = PETSC_TRUE;
1345:   osm->sort_indices      = PETSC_TRUE;
1346:   osm->dm_subdomains     = PETSC_FALSE;
1347:   osm->sub_mat_type      = NULL;

1349:   pc->data                 = (void*)osm;
1350:   pc->ops->apply           = PCApply_ASM;
1351:   pc->ops->matapply        = PCMatApply_ASM;
1352:   pc->ops->applytranspose  = PCApplyTranspose_ASM;
1353:   pc->ops->setup           = PCSetUp_ASM;
1354:   pc->ops->reset           = PCReset_ASM;
1355:   pc->ops->destroy         = PCDestroy_ASM;
1356:   pc->ops->setfromoptions  = PCSetFromOptions_ASM;
1357:   pc->ops->setuponblocks   = PCSetUpOnBlocks_ASM;
1358:   pc->ops->view            = PCView_ASM;
1359:   pc->ops->applyrichardson = NULL;

1361:   PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM);
1362:   PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM);
1363:   PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM);
1364:   PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM);
1365:   PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM);
1366:   PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM);
1367:   PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM);
1368:   PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM);
1369:   PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM);
1370:   PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM);
1371:   PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM);
1372:   return(0);
1373: }

1375: /*@C
1376:    PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
1377:    preconditioner for a any problem on a general grid.

1379:    Collective

1381:    Input Parameters:
1382: +  A - The global matrix operator
1383: -  n - the number of local blocks

1385:    Output Parameters:
1386: .  outis - the array of index sets defining the subdomains

1388:    Level: advanced

1390:    Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap
1391:     from these if you use PCASMSetLocalSubdomains()

1393:     In the Fortran version you must provide the array outis[] already allocated of length n.

1395: .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains()
1396: @*/
1397: PetscErrorCode  PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[])
1398: {
1399:   MatPartitioning mpart;
1400:   const char      *prefix;
1401:   PetscInt        i,j,rstart,rend,bs;
1402:   PetscBool       hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
1403:   Mat             Ad     = NULL, adj;
1404:   IS              ispart,isnumb,*is;
1405:   PetscErrorCode  ierr;

1410:   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);

1412:   /* Get prefix, row distribution, and block size */
1413:   MatGetOptionsPrefix(A,&prefix);
1414:   MatGetOwnershipRange(A,&rstart,&rend);
1415:   MatGetBlockSize(A,&bs);
1416:   if (rstart/bs*bs != rstart || rend/bs*bs != rend) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"bad row distribution [%D,%D) for matrix block size %D",rstart,rend,bs);

1418:   /* Get diagonal block from matrix if possible */
1419:   MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop);
1420:   if (hasop) {
1421:     MatGetDiagonalBlock(A,&Ad);
1422:   }
1423:   if (Ad) {
1424:     PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);
1425:     if (!isbaij) {PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);}
1426:   }
1427:   if (Ad && n > 1) {
1428:     PetscBool match,done;
1429:     /* Try to setup a good matrix partitioning if available */
1430:     MatPartitioningCreate(PETSC_COMM_SELF,&mpart);
1431:     PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);
1432:     MatPartitioningSetFromOptions(mpart);
1433:     PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);
1434:     if (!match) {
1435:       PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);
1436:     }
1437:     if (!match) { /* assume a "good" partitioner is available */
1438:       PetscInt       na;
1439:       const PetscInt *ia,*ja;
1440:       MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);
1441:       if (done) {
1442:         /* Build adjacency matrix by hand. Unfortunately a call to
1443:            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1444:            remove the block-aij structure and we cannot expect
1445:            MatPartitioning to split vertices as we need */
1446:         PetscInt       i,j,len,nnz,cnt,*iia=NULL,*jja=NULL;
1447:         const PetscInt *row;
1448:         nnz = 0;
1449:         for (i=0; i<na; i++) { /* count number of nonzeros */
1450:           len = ia[i+1] - ia[i];
1451:           row = ja + ia[i];
1452:           for (j=0; j<len; j++) {
1453:             if (row[j] == i) { /* don't count diagonal */
1454:               len--; break;
1455:             }
1456:           }
1457:           nnz += len;
1458:         }
1459:         PetscMalloc1(na+1,&iia);
1460:         PetscMalloc1(nnz,&jja);
1461:         nnz    = 0;
1462:         iia[0] = 0;
1463:         for (i=0; i<na; i++) { /* fill adjacency */
1464:           cnt = 0;
1465:           len = ia[i+1] - ia[i];
1466:           row = ja + ia[i];
1467:           for (j=0; j<len; j++) {
1468:             if (row[j] != i) { /* if not diagonal */
1469:               jja[nnz+cnt++] = row[j];
1470:             }
1471:           }
1472:           nnz     += cnt;
1473:           iia[i+1] = nnz;
1474:         }
1475:         /* Partitioning of the adjacency matrix */
1476:         MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);
1477:         MatPartitioningSetAdjacency(mpart,adj);
1478:         MatPartitioningSetNParts(mpart,n);
1479:         MatPartitioningApply(mpart,&ispart);
1480:         ISPartitioningToNumbering(ispart,&isnumb);
1481:         MatDestroy(&adj);
1482:         foundpart = PETSC_TRUE;
1483:       }
1484:       MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);
1485:     }
1486:     MatPartitioningDestroy(&mpart);
1487:   }

1489:   PetscMalloc1(n,&is);
1490:   *outis = is;

1492:   if (!foundpart) {

1494:     /* Partitioning by contiguous chunks of rows */

1496:     PetscInt mbs   = (rend-rstart)/bs;
1497:     PetscInt start = rstart;
1498:     for (i=0; i<n; i++) {
1499:       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
1500:       ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);
1501:       start += count;
1502:     }

1504:   } else {

1506:     /* Partitioning by adjacency of diagonal block  */

1508:     const PetscInt *numbering;
1509:     PetscInt       *count,nidx,*indices,*newidx,start=0;
1510:     /* Get node count in each partition */
1511:     PetscMalloc1(n,&count);
1512:     ISPartitioningCount(ispart,n,count);
1513:     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1514:       for (i=0; i<n; i++) count[i] *= bs;
1515:     }
1516:     /* Build indices from node numbering */
1517:     ISGetLocalSize(isnumb,&nidx);
1518:     PetscMalloc1(nidx,&indices);
1519:     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1520:     ISGetIndices(isnumb,&numbering);
1521:     PetscSortIntWithPermutation(nidx,numbering,indices);
1522:     ISRestoreIndices(isnumb,&numbering);
1523:     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1524:       PetscMalloc1(nidx*bs,&newidx);
1525:       for (i=0; i<nidx; i++) {
1526:         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
1527:       }
1528:       PetscFree(indices);
1529:       nidx   *= bs;
1530:       indices = newidx;
1531:     }
1532:     /* Shift to get global indices */
1533:     for (i=0; i<nidx; i++) indices[i] += rstart;

1535:     /* Build the index sets for each block */
1536:     for (i=0; i<n; i++) {
1537:       ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);
1538:       ISSort(is[i]);
1539:       start += count[i];
1540:     }

1542:     PetscFree(count);
1543:     PetscFree(indices);
1544:     ISDestroy(&isnumb);
1545:     ISDestroy(&ispart);

1547:   }
1548:   return(0);
1549: }

1551: /*@C
1552:    PCASMDestroySubdomains - Destroys the index sets created with
1553:    PCASMCreateSubdomains(). Should be called after setting subdomains
1554:    with PCASMSetLocalSubdomains().

1556:    Collective

1558:    Input Parameters:
1559: +  n - the number of index sets
1560: .  is - the array of index sets
1561: -  is_local - the array of local index sets, can be NULL

1563:    Level: advanced

1565: .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains()
1566: @*/
1567: PetscErrorCode  PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[])
1568: {
1569:   PetscInt       i;

1573:   if (n <= 0) return(0);
1574:   if (is) {
1576:     for (i=0; i<n; i++) { ISDestroy(&is[i]); }
1577:     PetscFree(is);
1578:   }
1579:   if (is_local) {
1581:     for (i=0; i<n; i++) { ISDestroy(&is_local[i]); }
1582:     PetscFree(is_local);
1583:   }
1584:   return(0);
1585: }

1587: /*@
1588:    PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1589:    preconditioner for a two-dimensional problem on a regular grid.

1591:    Not Collective

1593:    Input Parameters:
1594: +  m, n - the number of mesh points in the x and y directions
1595: .  M, N - the number of subdomains in the x and y directions
1596: .  dof - degrees of freedom per node
1597: -  overlap - overlap in mesh lines

1599:    Output Parameters:
1600: +  Nsub - the number of subdomains created
1601: .  is - array of index sets defining overlapping (if overlap > 0) subdomains
1602: -  is_local - array of index sets defining non-overlapping subdomains

1604:    Note:
1605:    Presently PCAMSCreateSubdomains2d() is valid only for sequential
1606:    preconditioners.  More general related routines are
1607:    PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains().

1609:    Level: advanced

1611: .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
1612:           PCASMSetOverlap()
1613: @*/
1614: PetscErrorCode  PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local)
1615: {
1616:   PetscInt       i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer;
1618:   PetscInt       nidx,*idx,loc,ii,jj,count;

1621:   if (dof != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," ");

1623:   *Nsub     = N*M;
1624:   PetscMalloc1(*Nsub,is);
1625:   PetscMalloc1(*Nsub,is_local);
1626:   ystart    = 0;
1627:   loc_outer = 0;
1628:   for (i=0; i<N; i++) {
1629:     height = n/N + ((n % N) > i); /* height of subdomain */
1630:     if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n");
1631:     yleft  = ystart - overlap; if (yleft < 0) yleft = 0;
1632:     yright = ystart + height + overlap; if (yright > n) yright = n;
1633:     xstart = 0;
1634:     for (j=0; j<M; j++) {
1635:       width = m/M + ((m % M) > j); /* width of subdomain */
1636:       if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m");
1637:       xleft  = xstart - overlap; if (xleft < 0) xleft = 0;
1638:       xright = xstart + width + overlap; if (xright > m) xright = m;
1639:       nidx   = (xright - xleft)*(yright - yleft);
1640:       PetscMalloc1(nidx,&idx);
1641:       loc    = 0;
1642:       for (ii=yleft; ii<yright; ii++) {
1643:         count = m*ii + xleft;
1644:         for (jj=xleft; jj<xright; jj++) idx[loc++] = count++;
1645:       }
1646:       ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);
1647:       if (overlap == 0) {
1648:         PetscObjectReference((PetscObject)(*is)[loc_outer]);

1650:         (*is_local)[loc_outer] = (*is)[loc_outer];
1651:       } else {
1652:         for (loc=0,ii=ystart; ii<ystart+height; ii++) {
1653:           for (jj=xstart; jj<xstart+width; jj++) {
1654:             idx[loc++] = m*ii + jj;
1655:           }
1656:         }
1657:         ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);
1658:       }
1659:       PetscFree(idx);
1660:       xstart += width;
1661:       loc_outer++;
1662:     }
1663:     ystart += height;
1664:   }
1665:   for (i=0; i<*Nsub; i++) { ISSort((*is)[i]); }
1666:   return(0);
1667: }

1669: /*@C
1670:     PCASMGetLocalSubdomains - Gets the local subdomains (for this processor
1671:     only) for the additive Schwarz preconditioner.

1673:     Not Collective

1675:     Input Parameter:
1676: .   pc - the preconditioner context

1678:     Output Parameters:
1679: +   n - the number of subdomains for this processor (default value = 1)
1680: .   is - the index sets that define the subdomains for this processor
1681: -   is_local - the index sets that define the local part of the subdomains for this processor (can be NULL)


1684:     Notes:
1685:     The IS numbering is in the parallel, global numbering of the vector.

1687:     Level: advanced

1689: .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
1690:           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices()
1691: @*/
1692: PetscErrorCode  PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[])
1693: {
1694:   PC_ASM         *osm = (PC_ASM*)pc->data;
1696:   PetscBool      match;

1702:   PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);
1703:   if (!match) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM");
1704:   if (n) *n = osm->n_local_true;
1705:   if (is) *is = osm->is;
1706:   if (is_local) *is_local = osm->is_local;
1707:   return(0);
1708: }

1710: /*@C
1711:     PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor
1712:     only) for the additive Schwarz preconditioner.

1714:     Not Collective

1716:     Input Parameter:
1717: .   pc - the preconditioner context

1719:     Output Parameters:
1720: +   n - the number of matrices for this processor (default value = 1)
1721: -   mat - the matrices

1723:     Level: advanced

1725:     Notes:
1726:     Call after PCSetUp() (or KSPSetUp()) but before PCApply() and before PCSetUpOnBlocks())

1728:            Usually one would use PCSetModifySubMatrices() to change the submatrices in building the preconditioner.

1730: .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
1731:           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubMatrices()
1732: @*/
1733: PetscErrorCode  PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1734: {
1735:   PC_ASM         *osm;
1737:   PetscBool      match;

1743:   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUp() or PCSetUp().");
1744:   PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);
1745:   if (!match) {
1746:     if (n) *n = 0;
1747:     if (mat) *mat = NULL;
1748:   } else {
1749:     osm = (PC_ASM*)pc->data;
1750:     if (n) *n = osm->n_local_true;
1751:     if (mat) *mat = osm->pmat;
1752:   }
1753:   return(0);
1754: }

1756: /*@
1757:     PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.

1759:     Logically Collective

1761:     Input Parameter:
1762: +   pc  - the preconditioner
1763: -   flg - boolean indicating whether to use subdomains defined by the DM

1765:     Options Database Key:
1766: .   -pc_asm_dm_subdomains

1768:     Level: intermediate

1770:     Notes:
1771:     PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(),
1772:     so setting either of the first two effectively turns the latter off.

1774: .seealso: PCASMGetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1775:           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1776: @*/
1777: PetscErrorCode  PCASMSetDMSubdomains(PC pc,PetscBool flg)
1778: {
1779:   PC_ASM         *osm = (PC_ASM*)pc->data;
1781:   PetscBool      match;

1786:   if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
1787:   PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);
1788:   if (match) {
1789:     osm->dm_subdomains = flg;
1790:   }
1791:   return(0);
1792: }

1794: /*@
1795:     PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1796:     Not Collective

1798:     Input Parameter:
1799: .   pc  - the preconditioner

1801:     Output Parameter:
1802: .   flg - boolean indicating whether to use subdomains defined by the DM

1804:     Level: intermediate

1806: .seealso: PCASMSetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1807:           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1808: @*/
1809: PetscErrorCode  PCASMGetDMSubdomains(PC pc,PetscBool* flg)
1810: {
1811:   PC_ASM         *osm = (PC_ASM*)pc->data;
1813:   PetscBool      match;

1818:   PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);
1819:   if (match) {
1820:     if (flg) *flg = osm->dm_subdomains;
1821:   }
1822:   return(0);
1823: }

1825: /*@
1826:      PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string.

1828:    Not Collective

1830:    Input Parameter:
1831: .  pc - the PC

1833:    Output Parameter:
1834: .  -pc_asm_sub_mat_type - name of matrix type

1836:    Level: advanced

1838: .seealso: PCASMSetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
1839: @*/
1840: PetscErrorCode  PCASMGetSubMatType(PC pc,MatType *sub_mat_type) {

1843:   PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type));
1844:   return(0);
1845: }

1847: /*@
1848:      PCASMSetSubMatType - Set the type of matrix used for ASM subsolves

1850:    Collective on Mat

1852:    Input Parameters:
1853: +  pc             - the PC object
1854: -  sub_mat_type   - matrix type

1856:    Options Database Key:
1857: .  -pc_asm_sub_mat_type  <sub_mat_type> - Sets the matrix type used for subsolves, for example, seqaijviennacl. If you specify a base name like aijviennacl, the corresponding sequential type is assumed.

1859:    Notes:
1860:    See "${PETSC_DIR}/include/petscmat.h" for available types

1862:   Level: advanced

1864: .seealso: PCASMGetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
1865: @*/
1866: PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type)
1867: {

1870:   PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type));
1871:   return(0);
1872: }