Actual source code: asm.c

petsc-3.13.6 2020-09-29
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 solve is same for all blocks, in the following KSP and PC objects:\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,PETSC_MAX_PATH_LEN,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,0);

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:       }
481:       else { /* interpolate the overlapping i-block solution to the local solution */
482:         VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);
483:         VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);
484:       }

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

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

507: static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y)
508: {
509:   PC_ASM         *osm = (PC_ASM*)pc->data;
511:   PetscInt       i,n_local_true = osm->n_local_true;
512:   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;

515:   /*
516:      Support for limiting the restriction or interpolation to only local
517:      subdomain values (leaving the other values 0).

519:      Note: these are reversed from the PCApply_ASM() because we are applying the
520:      transpose of the three terms
521:   */

523:   if (!(osm->type & PC_ASM_INTERPOLATE)) {
524:     forward = SCATTER_FORWARD_LOCAL;
525:     /* have to zero the work RHS since scatter may leave some slots empty */
526:     VecSet(osm->lx, 0.0);
527:   }
528:   if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL;

530:   /* zero the global and the local solutions */
531:   VecSet(y, 0.0);
532:   VecSet(osm->ly, 0.0);

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

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

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

545:     /* solve the overlapping i-block */
546:     PetscLogEventBegin(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);
547:     KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i]);
548:     KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);
549:     PetscLogEventEnd(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);

551:     if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution */
552:       VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);
553:       VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);
554:     } else { /* interpolate the overlapping i-block solution to the local solution */
555:       VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);
556:       VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);
557:     }

559:     if (i < n_local_true-1) {
560:       /* Restrict local RHS to the overlapping (i+1)-block RHS */
561:       VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);
562:       VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);
563:     }
564:   }
565:   /* Add the local solution to the global solution including the ghost nodes */
566:   VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse);
567:   VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse);

569:   return(0);

571: }

573: static PetscErrorCode PCReset_ASM(PC pc)
574: {
575:   PC_ASM         *osm = (PC_ASM*)pc->data;
577:   PetscInt       i;

580:   if (osm->ksp) {
581:     for (i=0; i<osm->n_local_true; i++) {
582:       KSPReset(osm->ksp[i]);
583:     }
584:   }
585:   if (osm->pmat) {
586:     if (osm->n_local_true > 0) {
587:       MatDestroySubMatrices(osm->n_local_true,&osm->pmat);
588:     }
589:   }
590:   if (osm->lrestriction) {
591:     VecScatterDestroy(&osm->restriction);
592:     for (i=0; i<osm->n_local_true; i++) {
593:       VecScatterDestroy(&osm->lrestriction[i]);
594:       if (osm->lprolongation) {VecScatterDestroy(&osm->lprolongation[i]);}
595:       VecDestroy(&osm->x[i]);
596:       VecDestroy(&osm->y[i]);
597:     }
598:     PetscFree(osm->lrestriction);
599:     if (osm->lprolongation) {PetscFree(osm->lprolongation);}
600:     PetscFree(osm->x);
601:     PetscFree(osm->y);

603:   }
604:   PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);
605:   ISDestroy(&osm->lis);
606:   VecDestroy(&osm->lx);
607:   VecDestroy(&osm->ly);
608:   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
609:     MatDestroyMatrices(osm->n_local_true, &osm->lmats);
610:   }

612:   PetscFree(osm->sub_mat_type);

614:   osm->is       = 0;
615:   osm->is_local = 0;
616:   return(0);
617: }

619: static PetscErrorCode PCDestroy_ASM(PC pc)
620: {
621:   PC_ASM         *osm = (PC_ASM*)pc->data;
623:   PetscInt       i;

626:   PCReset_ASM(pc);
627:   if (osm->ksp) {
628:     for (i=0; i<osm->n_local_true; i++) {
629:       KSPDestroy(&osm->ksp[i]);
630:     }
631:     PetscFree(osm->ksp);
632:   }
633:   PetscFree(pc->data);
634:   return(0);
635: }

637: static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc)
638: {
639:   PC_ASM         *osm = (PC_ASM*)pc->data;
641:   PetscInt       blocks,ovl;
642:   PetscBool      flg;
643:   PCASMType      asmtype;
644:   PCCompositeType loctype;
645:   char           sub_mat_type[256];

648:   PetscOptionsHead(PetscOptionsObject,"Additive Schwarz options");
649:   PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);
650:   PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);
651:   if (flg) {
652:     PCASMSetTotalSubdomains(pc,blocks,NULL,NULL);
653:     osm->dm_subdomains = PETSC_FALSE;
654:   }
655:   PetscOptionsInt("-pc_asm_local_blocks","Number of local subdomains","PCASMSetLocalSubdomains",osm->n_local_true,&blocks,&flg);
656:   if (flg) {
657:     PCASMSetLocalSubdomains(pc,blocks,NULL,NULL);
658:     osm->dm_subdomains = PETSC_FALSE;
659:   }
660:   PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);
661:   if (flg) {
662:     PCASMSetOverlap(pc,ovl);
663:     osm->dm_subdomains = PETSC_FALSE;
664:   }
665:   flg  = PETSC_FALSE;
666:   PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);
667:   if (flg) {PCASMSetType(pc,asmtype); }
668:   flg  = PETSC_FALSE;
669:   PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg);
670:   if (flg) {PCASMSetLocalType(pc,loctype); }
671:   PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg);
672:   if (flg) {
673:     PCASMSetSubMatType(pc,sub_mat_type);
674:   }
675:   PetscOptionsTail();
676:   return(0);
677: }

679: /*------------------------------------------------------------------------------------*/

681: static PetscErrorCode  PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[])
682: {
683:   PC_ASM         *osm = (PC_ASM*)pc->data;
685:   PetscInt       i;

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

691:   if (!pc->setupcalled) {
692:     if (is) {
693:       for (i=0; i<n; i++) {PetscObjectReference((PetscObject)is[i]);}
694:     }
695:     if (is_local) {
696:       for (i=0; i<n; i++) {PetscObjectReference((PetscObject)is_local[i]);}
697:     }
698:     PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);

700:     osm->n_local_true = n;
701:     osm->is           = 0;
702:     osm->is_local     = 0;
703:     if (is) {
704:       PetscMalloc1(n,&osm->is);
705:       for (i=0; i<n; i++) osm->is[i] = is[i];
706:       /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */
707:       osm->overlap = -1;
708:     }
709:     if (is_local) {
710:       PetscMalloc1(n,&osm->is_local);
711:       for (i=0; i<n; i++) osm->is_local[i] = is_local[i];
712:       if (!is) {
713:         PetscMalloc1(osm->n_local_true,&osm->is);
714:         for (i=0; i<osm->n_local_true; i++) {
715:           if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
716:             ISDuplicate(osm->is_local[i],&osm->is[i]);
717:             ISCopy(osm->is_local[i],osm->is[i]);
718:           } else {
719:             PetscObjectReference((PetscObject)osm->is_local[i]);
720:             osm->is[i] = osm->is_local[i];
721:           }
722:         }
723:       }
724:     }
725:   }
726:   return(0);
727: }

729: static PetscErrorCode  PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local)
730: {
731:   PC_ASM         *osm = (PC_ASM*)pc->data;
733:   PetscMPIInt    rank,size;
734:   PetscInt       n;

737:   if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N);
738:   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.");

740:   /*
741:      Split the subdomains equally among all processors
742:   */
743:   MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
744:   MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
745:   n    = N/size + ((N % size) > rank);
746:   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);
747:   if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp().");
748:   if (!pc->setupcalled) {
749:     PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);

751:     osm->n_local_true = n;
752:     osm->is           = 0;
753:     osm->is_local     = 0;
754:   }
755:   return(0);
756: }

758: static PetscErrorCode  PCASMSetOverlap_ASM(PC pc,PetscInt ovl)
759: {
760:   PC_ASM *osm = (PC_ASM*)pc->data;

763:   if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
764:   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp().");
765:   if (!pc->setupcalled) osm->overlap = ovl;
766:   return(0);
767: }

769: static PetscErrorCode  PCASMSetType_ASM(PC pc,PCASMType type)
770: {
771:   PC_ASM *osm = (PC_ASM*)pc->data;

774:   osm->type     = type;
775:   osm->type_set = PETSC_TRUE;
776:   return(0);
777: }

779: static PetscErrorCode  PCASMGetType_ASM(PC pc,PCASMType *type)
780: {
781:   PC_ASM *osm = (PC_ASM*)pc->data;

784:   *type = osm->type;
785:   return(0);
786: }

788: static PetscErrorCode  PCASMSetLocalType_ASM(PC pc, PCCompositeType type)
789: {
790:   PC_ASM *osm = (PC_ASM *) pc->data;

793:   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");
794:   osm->loctype = type;
795:   return(0);
796: }

798: static PetscErrorCode  PCASMGetLocalType_ASM(PC pc, PCCompositeType *type)
799: {
800:   PC_ASM *osm = (PC_ASM *) pc->data;

803:   *type = osm->loctype;
804:   return(0);
805: }

807: static PetscErrorCode  PCASMSetSortIndices_ASM(PC pc,PetscBool  doSort)
808: {
809:   PC_ASM *osm = (PC_ASM*)pc->data;

812:   osm->sort_indices = doSort;
813:   return(0);
814: }

816: static PetscErrorCode  PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
817: {
818:   PC_ASM         *osm = (PC_ASM*)pc->data;

822:   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");

824:   if (n_local) *n_local = osm->n_local_true;
825:   if (first_local) {
826:     MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));
827:     *first_local -= osm->n_local_true;
828:   }
829:   if (ksp) {
830:     /* Assume that local solves are now different; not necessarily
831:        true though!  This flag is used only for PCView_ASM() */
832:     *ksp                   = osm->ksp;
833:     osm->same_local_solves = PETSC_FALSE;
834:   }
835:   return(0);
836: }

838: static PetscErrorCode  PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type)
839: {
840:   PC_ASM         *osm = (PC_ASM*)pc->data;

845:   *sub_mat_type = osm->sub_mat_type;
846:   return(0);
847: }

849: static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type)
850: {
851:   PetscErrorCode    ierr;
852:   PC_ASM            *osm = (PC_ASM*)pc->data;

856:   PetscFree(osm->sub_mat_type);
857:   PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type);
858:   return(0);
859: }

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

864:     Collective on pc

866:     Input Parameters:
867: +   pc - the preconditioner context
868: .   n - the number of subdomains for this processor (default value = 1)
869: .   is - the index set that defines the subdomains for this processor
870:          (or NULL for PETSc to determine subdomains)
871: -   is_local - the index sets that define the local part of the subdomains for this processor, not used unless PCASMType is PC_ASM_RESTRICT
872:          (or NULL to not provide these)

874:     Options Database Key:
875:     To set the total number of subdomain blocks rather than specify the
876:     index sets, use the option
877: .    -pc_asm_local_blocks <blks> - Sets local blocks

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

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

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

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

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

892:     Level: advanced

894: .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
895:           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), PCASMType, PCASMSetType()
896: @*/
897: PetscErrorCode  PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
898: {

903:   PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));
904:   return(0);
905: }

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

912:     Collective on pc

914:     Input Parameters:
915: +   pc - the preconditioner context
916: .   N  - the number of subdomains for all processors
917: .   is - the index sets that define the subdomains for all processors
918:          (or NULL to ask PETSc to determine the subdomains)
919: -   is_local - the index sets that define the local part of the subdomains for this processor
920:          (or NULL to not provide this information)

922:     Options Database Key:
923:     To set the total number of subdomain blocks rather than specify the
924:     index sets, use the option
925: .    -pc_asm_blocks <blks> - Sets total blocks

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

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

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

935:     Use PCASMSetLocalSubdomains() to set local subdomains.

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

939:     Level: advanced

941: .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
942:           PCASMCreateSubdomains2D()
943: @*/
944: PetscErrorCode  PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[])
945: {

950:   PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local));
951:   return(0);
952: }

954: /*@
955:     PCASMSetOverlap - Sets the overlap between a pair of subdomains for the
956:     additive Schwarz preconditioner.  Either all or no processors in the
957:     PC communicator must call this routine.

959:     Logically Collective on pc

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

965:     Options Database Key:
966: .   -pc_asm_overlap <ovl> - Sets overlap

968:     Notes:
969:     By default the ASM preconditioner uses 1 block per processor.  To use
970:     multiple blocks per perocessor, see PCASMSetTotalSubdomains() and
971:     PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>).

973:     The overlap defaults to 1, so if one desires that no additional
974:     overlap be computed beyond what may have been set with a call to
975:     PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl
976:     must be set to be 0.  In particular, if one does not explicitly set
977:     the subdomains an Section 1.5 Writing Application Codes with PETSc code, then all overlap would be computed
978:     internally by PETSc, and using an overlap of 0 would result in an ASM
979:     variant that is equivalent to the block Jacobi preconditioner.

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

984:     Note that one can define initial index sets with any overlap via
985:     PCASMSetLocalSubdomains(); the routine
986:     PCASMSetOverlap() merely allows PETSc to extend that overlap further
987:     if desired.

989:     Level: intermediate

991: .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
992:           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), MatIncreaseOverlap()
993: @*/
994: PetscErrorCode  PCASMSetOverlap(PC pc,PetscInt ovl)
995: {

1001:   PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl));
1002:   return(0);
1003: }

1005: /*@
1006:     PCASMSetType - Sets the type of restriction and interpolation used
1007:     for local problems in the additive Schwarz method.

1009:     Logically Collective on pc

1011:     Input Parameters:
1012: +   pc  - the preconditioner context
1013: -   type - variant of ASM, one of
1014: .vb
1015:       PC_ASM_BASIC       - full interpolation and restriction
1016:       PC_ASM_RESTRICT    - full restriction, local processor interpolation
1017:       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
1018:       PC_ASM_NONE        - local processor restriction and interpolation
1019: .ve

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

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

1028:     Level: intermediate

1030: .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1031:           PCASMCreateSubdomains2D(), PCASMType, PCASMSetLocalType(), PCASMGetLocalType()
1032: @*/
1033: PetscErrorCode  PCASMSetType(PC pc,PCASMType type)
1034: {

1040:   PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type));
1041:   return(0);
1042: }

1044: /*@
1045:     PCASMGetType - Gets the type of restriction and interpolation used
1046:     for local problems in the additive Schwarz method.

1048:     Logically Collective on pc

1050:     Input Parameter:
1051: .   pc  - the preconditioner context

1053:     Output Parameter:
1054: .   type - variant of ASM, one of

1056: .vb
1057:       PC_ASM_BASIC       - full interpolation and restriction
1058:       PC_ASM_RESTRICT    - full restriction, local processor interpolation
1059:       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
1060:       PC_ASM_NONE        - local processor restriction and interpolation
1061: .ve

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

1066:     Level: intermediate

1068: .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1069:           PCASMCreateSubdomains2D(), PCASMType, PCASMSetType(), PCASMSetLocalType(), PCASMGetLocalType()
1070: @*/
1071: PetscErrorCode  PCASMGetType(PC pc,PCASMType *type)
1072: {

1077:   PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type));
1078:   return(0);
1079: }

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

1084:   Logically Collective on pc

1086:   Input Parameters:
1087: + pc  - the preconditioner context
1088: - type - type of composition, one of
1089: .vb
1090:   PC_COMPOSITE_ADDITIVE       - local additive combination
1091:   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
1092: .ve

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

1097:   Level: intermediate

1099: .seealso: PCASMSetType(), PCASMGetType(), PCASMGetLocalType(), PCASM, PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
1100: @*/
1101: PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type)
1102: {

1108:   PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));
1109:   return(0);
1110: }

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

1115:   Logically Collective on pc

1117:   Input Parameter:
1118: . pc  - the preconditioner context

1120:   Output Parameter:
1121: . type - type of composition, one of
1122: .vb
1123:   PC_COMPOSITE_ADDITIVE       - local additive combination
1124:   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
1125: .ve

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

1130:   Level: intermediate

1132: .seealso: PCASMSetType(), PCASMGetType(), PCASMSetLocalType(), PCASMCreate(), PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
1133: @*/
1134: PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type)
1135: {

1141:   PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));
1142:   return(0);
1143: }

1145: /*@
1146:     PCASMSetSortIndices - Determines whether subdomain indices are sorted.

1148:     Logically Collective on pc

1150:     Input Parameters:
1151: +   pc  - the preconditioner context
1152: -   doSort - sort the subdomain indices

1154:     Level: intermediate

1156: .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1157:           PCASMCreateSubdomains2D()
1158: @*/
1159: PetscErrorCode  PCASMSetSortIndices(PC pc,PetscBool doSort)
1160: {

1166:   PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));
1167:   return(0);
1168: }

1170: /*@C
1171:    PCASMGetSubKSP - Gets the local KSP contexts for all blocks on
1172:    this processor.

1174:    Collective on pc iff first_local is requested

1176:    Input Parameter:
1177: .  pc - the preconditioner context

1179:    Output Parameters:
1180: +  n_local - the number of blocks on this processor or NULL
1181: .  first_local - the global number of the first block on this processor or NULL,
1182:                  all processors must request or all must pass NULL
1183: -  ksp - the array of KSP contexts

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

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

1190:    Fortran note:
1191:    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.

1193:    Level: advanced

1195: .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(),
1196:           PCASMCreateSubdomains2D(),
1197: @*/
1198: PetscErrorCode  PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
1199: {

1204:   PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));
1205:   return(0);
1206: }

1208: /* -------------------------------------------------------------------------------------*/
1209: /*MC
1210:    PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
1211:            its own KSP object.

1213:    Options Database Keys:
1214: +  -pc_asm_blocks <blks> - Sets total blocks
1215: .  -pc_asm_overlap <ovl> - Sets overlap
1216: .  -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict
1217: -  -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive

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

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

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

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

1234:    Level: beginner

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

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

1246: M*/

1248: PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc)
1249: {
1251:   PC_ASM         *osm;

1254:   PetscNewLog(pc,&osm);

1256:   osm->n                 = PETSC_DECIDE;
1257:   osm->n_local           = 0;
1258:   osm->n_local_true      = PETSC_DECIDE;
1259:   osm->overlap           = 1;
1260:   osm->ksp               = 0;
1261:   osm->restriction       = 0;
1262:   osm->lprolongation     = 0;
1263:   osm->lrestriction      = 0;
1264:   osm->x                 = 0;
1265:   osm->y                 = 0;
1266:   osm->is                = 0;
1267:   osm->is_local          = 0;
1268:   osm->mat               = 0;
1269:   osm->pmat              = 0;
1270:   osm->type              = PC_ASM_RESTRICT;
1271:   osm->loctype           = PC_COMPOSITE_ADDITIVE;
1272:   osm->same_local_solves = PETSC_TRUE;
1273:   osm->sort_indices      = PETSC_TRUE;
1274:   osm->dm_subdomains     = PETSC_FALSE;
1275:   osm->sub_mat_type      = NULL;

1277:   pc->data                 = (void*)osm;
1278:   pc->ops->apply           = PCApply_ASM;
1279:   pc->ops->applytranspose  = PCApplyTranspose_ASM;
1280:   pc->ops->setup           = PCSetUp_ASM;
1281:   pc->ops->reset           = PCReset_ASM;
1282:   pc->ops->destroy         = PCDestroy_ASM;
1283:   pc->ops->setfromoptions  = PCSetFromOptions_ASM;
1284:   pc->ops->setuponblocks   = PCSetUpOnBlocks_ASM;
1285:   pc->ops->view            = PCView_ASM;
1286:   pc->ops->applyrichardson = 0;

1288:   PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM);
1289:   PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM);
1290:   PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM);
1291:   PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM);
1292:   PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM);
1293:   PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM);
1294:   PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM);
1295:   PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM);
1296:   PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM);
1297:   PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM);
1298:   PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM);
1299:   return(0);
1300: }

1302: /*@C
1303:    PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
1304:    preconditioner for a any problem on a general grid.

1306:    Collective

1308:    Input Parameters:
1309: +  A - The global matrix operator
1310: -  n - the number of local blocks

1312:    Output Parameters:
1313: .  outis - the array of index sets defining the subdomains

1315:    Level: advanced

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

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

1322: .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains()
1323: @*/
1324: PetscErrorCode  PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[])
1325: {
1326:   MatPartitioning mpart;
1327:   const char      *prefix;
1328:   PetscInt        i,j,rstart,rend,bs;
1329:   PetscBool       hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
1330:   Mat             Ad     = NULL, adj;
1331:   IS              ispart,isnumb,*is;
1332:   PetscErrorCode  ierr;

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

1339:   /* Get prefix, row distribution, and block size */
1340:   MatGetOptionsPrefix(A,&prefix);
1341:   MatGetOwnershipRange(A,&rstart,&rend);
1342:   MatGetBlockSize(A,&bs);
1343:   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);

1345:   /* Get diagonal block from matrix if possible */
1346:   MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop);
1347:   if (hasop) {
1348:     MatGetDiagonalBlock(A,&Ad);
1349:   }
1350:   if (Ad) {
1351:     PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);
1352:     if (!isbaij) {PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);}
1353:   }
1354:   if (Ad && n > 1) {
1355:     PetscBool match,done;
1356:     /* Try to setup a good matrix partitioning if available */
1357:     MatPartitioningCreate(PETSC_COMM_SELF,&mpart);
1358:     PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);
1359:     MatPartitioningSetFromOptions(mpart);
1360:     PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);
1361:     if (!match) {
1362:       PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);
1363:     }
1364:     if (!match) { /* assume a "good" partitioner is available */
1365:       PetscInt       na;
1366:       const PetscInt *ia,*ja;
1367:       MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);
1368:       if (done) {
1369:         /* Build adjacency matrix by hand. Unfortunately a call to
1370:            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1371:            remove the block-aij structure and we cannot expect
1372:            MatPartitioning to split vertices as we need */
1373:         PetscInt       i,j,len,nnz,cnt,*iia=0,*jja=0;
1374:         const PetscInt *row;
1375:         nnz = 0;
1376:         for (i=0; i<na; i++) { /* count number of nonzeros */
1377:           len = ia[i+1] - ia[i];
1378:           row = ja + ia[i];
1379:           for (j=0; j<len; j++) {
1380:             if (row[j] == i) { /* don't count diagonal */
1381:               len--; break;
1382:             }
1383:           }
1384:           nnz += len;
1385:         }
1386:         PetscMalloc1(na+1,&iia);
1387:         PetscMalloc1(nnz,&jja);
1388:         nnz    = 0;
1389:         iia[0] = 0;
1390:         for (i=0; i<na; i++) { /* fill adjacency */
1391:           cnt = 0;
1392:           len = ia[i+1] - ia[i];
1393:           row = ja + ia[i];
1394:           for (j=0; j<len; j++) {
1395:             if (row[j] != i) { /* if not diagonal */
1396:               jja[nnz+cnt++] = row[j];
1397:             }
1398:           }
1399:           nnz     += cnt;
1400:           iia[i+1] = nnz;
1401:         }
1402:         /* Partitioning of the adjacency matrix */
1403:         MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);
1404:         MatPartitioningSetAdjacency(mpart,adj);
1405:         MatPartitioningSetNParts(mpart,n);
1406:         MatPartitioningApply(mpart,&ispart);
1407:         ISPartitioningToNumbering(ispart,&isnumb);
1408:         MatDestroy(&adj);
1409:         foundpart = PETSC_TRUE;
1410:       }
1411:       MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);
1412:     }
1413:     MatPartitioningDestroy(&mpart);
1414:   }

1416:   PetscMalloc1(n,&is);
1417:   *outis = is;

1419:   if (!foundpart) {

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

1423:     PetscInt mbs   = (rend-rstart)/bs;
1424:     PetscInt start = rstart;
1425:     for (i=0; i<n; i++) {
1426:       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
1427:       ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);
1428:       start += count;
1429:     }

1431:   } else {

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

1435:     const PetscInt *numbering;
1436:     PetscInt       *count,nidx,*indices,*newidx,start=0;
1437:     /* Get node count in each partition */
1438:     PetscMalloc1(n,&count);
1439:     ISPartitioningCount(ispart,n,count);
1440:     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1441:       for (i=0; i<n; i++) count[i] *= bs;
1442:     }
1443:     /* Build indices from node numbering */
1444:     ISGetLocalSize(isnumb,&nidx);
1445:     PetscMalloc1(nidx,&indices);
1446:     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1447:     ISGetIndices(isnumb,&numbering);
1448:     PetscSortIntWithPermutation(nidx,numbering,indices);
1449:     ISRestoreIndices(isnumb,&numbering);
1450:     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1451:       PetscMalloc1(nidx*bs,&newidx);
1452:       for (i=0; i<nidx; i++) {
1453:         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
1454:       }
1455:       PetscFree(indices);
1456:       nidx   *= bs;
1457:       indices = newidx;
1458:     }
1459:     /* Shift to get global indices */
1460:     for (i=0; i<nidx; i++) indices[i] += rstart;

1462:     /* Build the index sets for each block */
1463:     for (i=0; i<n; i++) {
1464:       ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);
1465:       ISSort(is[i]);
1466:       start += count[i];
1467:     }

1469:     PetscFree(count);
1470:     PetscFree(indices);
1471:     ISDestroy(&isnumb);
1472:     ISDestroy(&ispart);

1474:   }
1475:   return(0);
1476: }

1478: /*@C
1479:    PCASMDestroySubdomains - Destroys the index sets created with
1480:    PCASMCreateSubdomains(). Should be called after setting subdomains
1481:    with PCASMSetLocalSubdomains().

1483:    Collective

1485:    Input Parameters:
1486: +  n - the number of index sets
1487: .  is - the array of index sets
1488: -  is_local - the array of local index sets, can be NULL

1490:    Level: advanced

1492: .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains()
1493: @*/
1494: PetscErrorCode  PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[])
1495: {
1496:   PetscInt       i;

1500:   if (n <= 0) return(0);
1501:   if (is) {
1503:     for (i=0; i<n; i++) { ISDestroy(&is[i]); }
1504:     PetscFree(is);
1505:   }
1506:   if (is_local) {
1508:     for (i=0; i<n; i++) { ISDestroy(&is_local[i]); }
1509:     PetscFree(is_local);
1510:   }
1511:   return(0);
1512: }

1514: /*@
1515:    PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1516:    preconditioner for a two-dimensional problem on a regular grid.

1518:    Not Collective

1520:    Input Parameters:
1521: +  m, n - the number of mesh points in the x and y directions
1522: .  M, N - the number of subdomains in the x and y directions
1523: .  dof - degrees of freedom per node
1524: -  overlap - overlap in mesh lines

1526:    Output Parameters:
1527: +  Nsub - the number of subdomains created
1528: .  is - array of index sets defining overlapping (if overlap > 0) subdomains
1529: -  is_local - array of index sets defining non-overlapping subdomains

1531:    Note:
1532:    Presently PCAMSCreateSubdomains2d() is valid only for sequential
1533:    preconditioners.  More general related routines are
1534:    PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains().

1536:    Level: advanced

1538: .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
1539:           PCASMSetOverlap()
1540: @*/
1541: PetscErrorCode  PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local)
1542: {
1543:   PetscInt       i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer;
1545:   PetscInt       nidx,*idx,loc,ii,jj,count;

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

1550:   *Nsub     = N*M;
1551:   PetscMalloc1(*Nsub,is);
1552:   PetscMalloc1(*Nsub,is_local);
1553:   ystart    = 0;
1554:   loc_outer = 0;
1555:   for (i=0; i<N; i++) {
1556:     height = n/N + ((n % N) > i); /* height of subdomain */
1557:     if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n");
1558:     yleft  = ystart - overlap; if (yleft < 0) yleft = 0;
1559:     yright = ystart + height + overlap; if (yright > n) yright = n;
1560:     xstart = 0;
1561:     for (j=0; j<M; j++) {
1562:       width = m/M + ((m % M) > j); /* width of subdomain */
1563:       if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m");
1564:       xleft  = xstart - overlap; if (xleft < 0) xleft = 0;
1565:       xright = xstart + width + overlap; if (xright > m) xright = m;
1566:       nidx   = (xright - xleft)*(yright - yleft);
1567:       PetscMalloc1(nidx,&idx);
1568:       loc    = 0;
1569:       for (ii=yleft; ii<yright; ii++) {
1570:         count = m*ii + xleft;
1571:         for (jj=xleft; jj<xright; jj++) idx[loc++] = count++;
1572:       }
1573:       ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);
1574:       if (overlap == 0) {
1575:         PetscObjectReference((PetscObject)(*is)[loc_outer]);

1577:         (*is_local)[loc_outer] = (*is)[loc_outer];
1578:       } else {
1579:         for (loc=0,ii=ystart; ii<ystart+height; ii++) {
1580:           for (jj=xstart; jj<xstart+width; jj++) {
1581:             idx[loc++] = m*ii + jj;
1582:           }
1583:         }
1584:         ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);
1585:       }
1586:       PetscFree(idx);
1587:       xstart += width;
1588:       loc_outer++;
1589:     }
1590:     ystart += height;
1591:   }
1592:   for (i=0; i<*Nsub; i++) { ISSort((*is)[i]); }
1593:   return(0);
1594: }

1596: /*@C
1597:     PCASMGetLocalSubdomains - Gets the local subdomains (for this processor
1598:     only) for the additive Schwarz preconditioner.

1600:     Not Collective

1602:     Input Parameter:
1603: .   pc - the preconditioner context

1605:     Output Parameters:
1606: +   n - the number of subdomains for this processor (default value = 1)
1607: .   is - the index sets that define the subdomains for this processor
1608: -   is_local - the index sets that define the local part of the subdomains for this processor (can be NULL)


1611:     Notes:
1612:     The IS numbering is in the parallel, global numbering of the vector.

1614:     Level: advanced

1616: .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
1617:           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices()
1618: @*/
1619: PetscErrorCode  PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[])
1620: {
1621:   PC_ASM         *osm = (PC_ASM*)pc->data;
1623:   PetscBool      match;

1629:   PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);
1630:   if (!match) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM");
1631:   if (n) *n = osm->n_local_true;
1632:   if (is) *is = osm->is;
1633:   if (is_local) *is_local = osm->is_local;
1634:   return(0);
1635: }

1637: /*@C
1638:     PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor
1639:     only) for the additive Schwarz preconditioner.

1641:     Not Collective

1643:     Input Parameter:
1644: .   pc - the preconditioner context

1646:     Output Parameters:
1647: +   n - the number of matrices for this processor (default value = 1)
1648: -   mat - the matrices

1650:     Level: advanced

1652:     Notes:
1653:     Call after PCSetUp() (or KSPSetUp()) but before PCApply() and before PCSetUpOnBlocks())

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

1657: .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
1658:           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubMatrices()
1659: @*/
1660: PetscErrorCode  PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1661: {
1662:   PC_ASM         *osm;
1664:   PetscBool      match;

1670:   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUp() or PCSetUp().");
1671:   PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);
1672:   if (!match) {
1673:     if (n) *n = 0;
1674:     if (mat) *mat = NULL;
1675:   } else {
1676:     osm = (PC_ASM*)pc->data;
1677:     if (n) *n = osm->n_local_true;
1678:     if (mat) *mat = osm->pmat;
1679:   }
1680:   return(0);
1681: }

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

1686:     Logically Collective

1688:     Input Parameter:
1689: +   pc  - the preconditioner
1690: -   flg - boolean indicating whether to use subdomains defined by the DM

1692:     Options Database Key:
1693: .   -pc_asm_dm_subdomains

1695:     Level: intermediate

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

1701: .seealso: PCASMGetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1702:           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1703: @*/
1704: PetscErrorCode  PCASMSetDMSubdomains(PC pc,PetscBool flg)
1705: {
1706:   PC_ASM         *osm = (PC_ASM*)pc->data;
1708:   PetscBool      match;

1713:   if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
1714:   PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);
1715:   if (match) {
1716:     osm->dm_subdomains = flg;
1717:   }
1718:   return(0);
1719: }

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

1725:     Input Parameter:
1726: .   pc  - the preconditioner

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

1731:     Level: intermediate

1733: .seealso: PCASMSetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1734:           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1735: @*/
1736: PetscErrorCode  PCASMGetDMSubdomains(PC pc,PetscBool* flg)
1737: {
1738:   PC_ASM         *osm = (PC_ASM*)pc->data;
1740:   PetscBool      match;

1745:   PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);
1746:   if (match) {
1747:     if (flg) *flg = osm->dm_subdomains;
1748:   }
1749:   return(0);
1750: }

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

1755:    Not Collective

1757:    Input Parameter:
1758: .  pc - the PC

1760:    Output Parameter:
1761: .  -pc_asm_sub_mat_type - name of matrix type

1763:    Level: advanced

1765: .seealso: PCASMSetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
1766: @*/
1767: PetscErrorCode  PCASMGetSubMatType(PC pc,MatType *sub_mat_type){

1770:   PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type));
1771:   return(0);
1772: }

1774: /*@
1775:      PCASMSetSubMatType - Set the type of matrix used for ASM subsolves

1777:    Collective on Mat

1779:    Input Parameters:
1780: +  pc             - the PC object
1781: -  sub_mat_type   - matrix type

1783:    Options Database Key:
1784: .  -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.

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

1789:   Level: advanced

1791: .seealso: PCASMGetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
1792: @*/
1793: PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type)
1794: {

1797:   PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type));
1798:   return(0);
1799: }