Actual source code: asm.c

  1: #define PETSCKSP_DLL

  3: /*
  4:   This file defines an additive Schwarz preconditioner for any Mat implementation.

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

 10:        n - total number of true subdomains on all processors
 11:        n_local_true - actual number of subdomains on this processor
 12:        n_local = maximum over all processors of n_local_true
 13: */
 14:  #include private/pcimpl.h

 16: typedef struct {
 17:   PetscInt   n,n_local,n_local_true;
 18:   PetscInt   overlap;             /* overlap requested by user */
 19:   KSP        *ksp;                /* linear solvers for each block */
 20:   VecScatter *scat;               /* mapping to subregion */
 21:   Vec        *x,*y;               /* work vectors */
 22:   IS         *is;                 /* index set that defines each subdomain */
 23:   Mat        *mat,*pmat;          /* mat is not currently used */
 24:   PCASMType  type;                /* use reduced interpolation, restriction or both */
 25:   PetscTruth type_set;            /* if user set this value (so won't change it for symmetric problems) */
 26:   PetscTruth same_local_solves;   /* flag indicating whether all local solvers are same */
 27:   PetscTruth inplace;             /* indicates that the sub-matrices are deleted after 
 28:                                      PCSetUpOnBlocks() is done. Similar to inplace 
 29:                                      factorization in the case of LU and ILU */
 30: } PC_ASM;

 34: static PetscErrorCode PCView_ASM(PC pc,PetscViewer viewer)
 35: {
 36:   PC_ASM         *osm = (PC_ASM*)pc->data;
 38:   PetscMPIInt    rank;
 39:   PetscInt       i,bsz;
 40:   PetscTruth     iascii,isstring;
 41:   PetscViewer    sviewer;


 45:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
 46:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
 47:   if (iascii) {
 48:     if (osm->n > 0) {
 49:       PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: total subdomain blocks = %D, amount of overlap = %D\n",osm->n,osm->overlap);
 50:     } else {
 51:       PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: total subdomain blocks not yet set, amount of overlap = %D\n",osm->overlap);
 52:     }
 53:     PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: restriction/interpolation type - %s\n",PCASMTypes[osm->type]);
 54:     MPI_Comm_rank(((PetscObject)pc)->comm,&rank);
 55:     if (osm->same_local_solves) {
 56:       if (osm->ksp) {
 57:         PetscViewerASCIIPrintf(viewer,"  Local solve is same for all blocks, in the following KSP and PC objects:\n");
 58:         PetscViewerGetSingleton(viewer,&sviewer);
 59:         if (!rank) {
 60:           PetscViewerASCIIPushTab(viewer);
 61:           KSPView(osm->ksp[0],sviewer);
 62:           PetscViewerASCIIPopTab(viewer);
 63:         }
 64:         PetscViewerRestoreSingleton(viewer,&sviewer);
 65:       }
 66:     } else {
 67:       PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,osm->n_local_true);
 68:       PetscViewerFlush(viewer);
 69:       PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following KSP and PC objects:\n");
 70:       PetscViewerASCIIPushTab(viewer);
 71:       PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");
 72:       for (i=0; i<osm->n_local; i++) {
 73:         PetscViewerGetSingleton(viewer,&sviewer);
 74:         if (i < osm->n_local_true) {
 75:           ISGetLocalSize(osm->is[i],&bsz);
 76:           PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);
 77:           KSPView(osm->ksp[i],sviewer);
 78:           PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");
 79:         }
 80:         PetscViewerRestoreSingleton(viewer,&sviewer);
 81:       }
 82:       PetscViewerASCIIPopTab(viewer);
 83:       PetscViewerFlush(viewer);
 84:     }
 85:   } else if (isstring) {
 86:     PetscViewerStringSPrintf(viewer," blocks=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCASMTypes[osm->type]);
 87:     PetscViewerGetSingleton(viewer,&sviewer);
 88:     if (osm->ksp) {KSPView(osm->ksp[0],sviewer);}
 89:     PetscViewerRestoreSingleton(viewer,&sviewer);
 90:   } else {
 91:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCASM",((PetscObject)viewer)->type_name);
 92:   }
 93:   return(0);
 94: }

 98: static PetscErrorCode PCASMPrintSubdomains(PC pc)
 99: {
100:   PC_ASM         *osm  = (PC_ASM*)pc->data;
101:   const char     *prefix;
102:   char           fname[PETSC_MAX_PATH_LEN+1];
103:   PetscViewer    viewer;
104:   PetscInt       i,j,nidx;
105:   const PetscInt *idx;
108:   PCGetOptionsPrefix(pc,&prefix);
109:   PetscOptionsGetString(prefix,"-pc_asm_print_subdomains",
110:                                fname,PETSC_MAX_PATH_LEN,PETSC_NULL);
111:   if (fname[0] == 0) { PetscStrcpy(fname,"stdout"); };
112:   PetscViewerASCIIOpen(((PetscObject)pc)->comm,fname,&viewer);
113:   for (i=0;i<osm->n_local_true;i++) {
114:     ISGetLocalSize(osm->is[i],&nidx);
115:     ISGetIndices(osm->is[i],&idx);
116:     for (j=0; j<nidx; j++) {
117:       PetscViewerASCIISynchronizedPrintf(viewer,"%D ",idx[j]);
118:     }
119:     ISRestoreIndices(osm->is[i],&idx);
120:     PetscViewerASCIISynchronizedPrintf(viewer,"\n");
121:   }
122:   PetscViewerFlush(viewer);
123:   PetscViewerDestroy(viewer);
124:   return(0);
125: }

129: static PetscErrorCode PCSetUp_ASM(PC pc)
130: {
131:   PC_ASM         *osm  = (PC_ASM*)pc->data;
133:   PetscTruth     symset,flg;
134:   PetscInt       i,m;
135:   PetscMPIInt    size;
136:   MatReuse       scall = MAT_REUSE_MATRIX;
137:   IS             isl;
138:   KSP            ksp;
139:   PC             subpc;
140:   const char     *prefix,*pprefix;
141:   Vec            vec;

144:   if (!pc->setupcalled) {

146:     if (!osm->type_set) {
147:       MatIsSymmetricKnown(pc->pmat,&symset,&flg);
148:       if (symset && flg) { osm->type = PC_ASM_BASIC; }
149:     }

151:     if (osm->n == PETSC_DECIDE && osm->n_local_true < 1) {
152:       /* no subdomains given, use one per processor */
153:       osm->n_local = osm->n_local_true = 1;
154:       MPI_Comm_size(((PetscObject)pc)->comm,&size);
155:       osm->n = size;
156:     } else if (osm->n == PETSC_DECIDE) {
157:       /* determine global number of subdomains */
158:       PetscInt inwork[2],outwork[2];
159:       inwork[0] = inwork[1] = osm->n_local_true;
160:       MPI_Allreduce(inwork,outwork,1,MPIU_2INT,PetscMaxSum_Op,((PetscObject)pc)->comm);
161:       osm->n_local = outwork[0];
162:       osm->n       = outwork[1];
163:     }

165:     if (!osm->is){ /* create the index sets */
166:       PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is);
167:     }
168:     PCGetOptionsPrefix(pc,&prefix);
169:     PetscOptionsHasName(prefix,"-pc_asm_print_subdomains",&flg);
170:     if (flg) { PCASMPrintSubdomains(pc); }

172:     /*  Extend the "overlapping" regions by a number of steps  */
173:     MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap);
174:     for (i=0; i<osm->n_local_true; i++) {
175:       ISSort(osm->is[i]);
176:     }

178:     /* Create the local work vectors and scatter contexts */
179:     MatGetVecs(pc->pmat,&vec,0);
180:     PetscMalloc(osm->n_local*sizeof(VecScatter *),&osm->scat);
181:     PetscMalloc(osm->n_local*sizeof(Vec *),&osm->x);
182:     PetscMalloc(osm->n_local*sizeof(Vec *),&osm->y);
183:     for (i=0; i<osm->n_local_true; i++) {
184:       ISGetLocalSize(osm->is[i],&m);
185:       VecCreateSeq(PETSC_COMM_SELF,m,&osm->x[i]);
186:       VecDuplicate(osm->x[i],&osm->y[i]);
187:       ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);
188:       VecScatterCreate(vec,osm->is[i],osm->x[i],isl,&osm->scat[i]);
189:       ISDestroy(isl);
190:     }
191:     for (i=osm->n_local_true; i<osm->n_local; i++) {
192:       VecCreateSeq(PETSC_COMM_SELF,0,&osm->x[i]);
193:       VecDuplicate(osm->x[i],&osm->y[i]);
194:       ISCreateStride(PETSC_COMM_SELF,0,0,1,&isl);
195:       VecScatterCreate(vec,isl,osm->x[i],isl,&osm->scat[i]);
196:       ISDestroy(isl);
197:     }
198:     VecDestroy(vec);

200:     /* Create the local solvers */
201:     PetscMalloc(osm->n_local_true*sizeof(KSP *),&osm->ksp);
202:     for (i=0; i<osm->n_local_true; i++) {
203:       KSPCreate(PETSC_COMM_SELF,&ksp);
204:       PetscLogObjectParent(pc,ksp);
205:       PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);
206:       KSPSetType(ksp,KSPPREONLY);
207:       KSPGetPC(ksp,&subpc);
208:       PCGetOptionsPrefix(pc,&prefix);
209:       KSPSetOptionsPrefix(ksp,prefix);
210:       KSPAppendOptionsPrefix(ksp,"sub_");
211:       osm->ksp[i] = ksp;
212:     }
213:     scall = MAT_INITIAL_MATRIX;

215:   } else {
216:     /* 
217:        Destroy the blocks from the previous iteration
218:     */
219:     if (pc->flag == DIFFERENT_NONZERO_PATTERN || osm->inplace) {
220:       if (!osm->inplace) {
221:         MatDestroyMatrices(osm->n_local_true,&osm->pmat);
222:       }
223:       scall = MAT_INITIAL_MATRIX;
224:     }
225:   }

227:   /* 
228:      Extract out the submatrices
229:   */
230:   MatGetSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);
231:   if (scall == MAT_INITIAL_MATRIX) {
232:     PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);
233:     for (i=0; i<osm->n_local_true; i++) {
234:       PetscLogObjectParent(pc,osm->pmat[i]);
235:       PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);
236:     }
237:   }

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

243:   /* 
244:      Loop over subdomains putting them into local ksp
245:   */
246:   for (i=0; i<osm->n_local_true; i++) {
247:     KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i],pc->flag);
248:     if (!pc->setupcalled) {
249:       KSPSetFromOptions(osm->ksp[i]);
250:     }
251:   }

253:   return(0);
254: }

258: static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc)
259: {
260:   PC_ASM         *osm = (PC_ASM*)pc->data;
262:   PetscInt       i;

265:   for (i=0; i<osm->n_local_true; i++) {
266:     KSPSetUp(osm->ksp[i]);
267:   }
268:   /* 
269:      If inplace flag is set, then destroy the matrix after the setup
270:      on blocks is done.
271:   */
272:   if (osm->inplace && osm->n_local_true > 0) {
273:     MatDestroyMatrices(osm->n_local_true,&osm->pmat);
274:   }
275:   return(0);
276: }

280: static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y)
281: {
282:   PC_ASM         *osm = (PC_ASM*)pc->data;
284:   PetscInt       i,n_local = osm->n_local,n_local_true = osm->n_local_true;
285:   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;

288:   /*
289:      Support for limiting the restriction or interpolation to only local 
290:      subdomain values (leaving the other values 0). 
291:   */
292:   if (!(osm->type & PC_ASM_RESTRICT)) {
293:     forward = SCATTER_FORWARD_LOCAL;
294:     /* have to zero the work RHS since scatter may leave some slots empty */
295:     for (i=0; i<n_local_true; i++) {
296:       VecSet(osm->x[i],0.0);
297:     }
298:   }
299:   if (!(osm->type & PC_ASM_INTERPOLATE)) {
300:     reverse = SCATTER_REVERSE_LOCAL;
301:   }

303:   for (i=0; i<n_local; i++) {
304:     VecScatterBegin(osm->scat[i],x,osm->x[i],INSERT_VALUES,forward);
305:   }
306:   VecSet(y,0.0);
307:   /* do the local solves */
308:   for (i=0; i<n_local_true; i++) {
309:     VecScatterEnd(osm->scat[i],x,osm->x[i],INSERT_VALUES,forward);
310:     KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);
311:     VecScatterBegin(osm->scat[i],osm->y[i],y,ADD_VALUES,reverse);
312:   }
313:   /* handle the rest of the scatters that do not have local solves */
314:   for (i=n_local_true; i<n_local; i++) {
315:     VecScatterEnd(osm->scat[i],x,osm->x[i],INSERT_VALUES,forward);
316:     VecScatterBegin(osm->scat[i],osm->y[i],y,ADD_VALUES,reverse);
317:   }
318:   for (i=0; i<n_local; i++) {
319:     VecScatterEnd(osm->scat[i],osm->y[i],y,ADD_VALUES,reverse);
320:   }
321:   return(0);
322: }

326: static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y)
327: {
328:   PC_ASM         *osm = (PC_ASM*)pc->data;
330:   PetscInt       i,n_local = osm->n_local,n_local_true = osm->n_local_true;
331:   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;

334:   /*
335:      Support for limiting the restriction or interpolation to only local 
336:      subdomain values (leaving the other values 0).

338:      Note: these are reversed from the PCApply_ASM() because we are applying the 
339:      transpose of the three terms 
340:   */
341:   if (!(osm->type & PC_ASM_INTERPOLATE)) {
342:     forward = SCATTER_FORWARD_LOCAL;
343:     /* have to zero the work RHS since scatter may leave some slots empty */
344:     for (i=0; i<n_local_true; i++) {
345:       VecSet(osm->x[i],0.0);
346:     }
347:   }
348:   if (!(osm->type & PC_ASM_RESTRICT)) {
349:     reverse = SCATTER_REVERSE_LOCAL;
350:   }

352:   for (i=0; i<n_local; i++) {
353:     VecScatterBegin(osm->scat[i],x,osm->x[i],INSERT_VALUES,forward);
354:   }
355:   VecSet(y,0.0);
356:   /* do the local solves */
357:   for (i=0; i<n_local_true; i++) {
358:     VecScatterEnd(osm->scat[i],x,osm->x[i],INSERT_VALUES,forward);
359:     KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);
360:     VecScatterBegin(osm->scat[i],osm->y[i],y,ADD_VALUES,reverse);
361:   }
362:   /* handle the rest of the scatters that do not have local solves */
363:   for (i=n_local_true; i<n_local; i++) {
364:     VecScatterEnd(osm->scat[i],x,osm->x[i],INSERT_VALUES,forward);
365:     VecScatterBegin(osm->scat[i],osm->y[i],y,ADD_VALUES,reverse);
366:   }
367:   for (i=0; i<n_local; i++) {
368:     VecScatterEnd(osm->scat[i],osm->y[i],y,ADD_VALUES,reverse);
369:   }
370:   return(0);
371: }

375: static PetscErrorCode PCDestroy_ASM(PC pc)
376: {
377:   PC_ASM         *osm = (PC_ASM*)pc->data;
379:   PetscInt       i;

382:   if (osm->ksp) {
383:     for (i=0; i<osm->n_local_true; i++) {
384:       KSPDestroy(osm->ksp[i]);
385:     }
386:     PetscFree(osm->ksp);
387:   }
388:   if (osm->pmat && !osm->inplace) {
389:     if (osm->n_local_true > 0) {
390:       MatDestroyMatrices(osm->n_local_true,&osm->pmat);
391:     }
392:   }
393:   if (osm->scat) {
394:     for (i=0; i<osm->n_local; i++) {
395:       VecScatterDestroy(osm->scat[i]);
396:       VecDestroy(osm->x[i]);
397:       VecDestroy(osm->y[i]);
398:     }
399:     PetscFree(osm->scat);
400:     PetscFree(osm->x);
401:     PetscFree(osm->y);
402:   }
403:   if (osm->is) {
404:     PCASMDestroySubdomains(osm->n_local_true,osm->is);
405:   }
406:   PetscFree(osm);
407:   return(0);
408: }

412: static PetscErrorCode PCSetFromOptions_ASM(PC pc)
413: {
414:   PC_ASM         *osm = (PC_ASM*)pc->data;
416:   PetscInt       blocks,ovl;
417:   PetscTruth     symset,flg;
418:   PCASMType      asmtype;

421:   /* set the type to symmetric if matrix is symmetric */
422:   if (!osm->type_set && pc->pmat) {
423:     MatIsSymmetricKnown(pc->pmat,&symset,&flg);
424:     if (symset && flg) { osm->type = PC_ASM_BASIC; }
425:   }
426:   PetscOptionsHead("Additive Schwarz options");
427:     PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",
428:                            osm->n,&blocks,&flg);
429:     if (flg) {PCASMSetTotalSubdomains(pc,blocks,PETSC_NULL); }
430:     PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",
431:                            osm->overlap,&ovl,&flg);
432:     if (flg) {PCASMSetOverlap(pc,ovl); }
433:     PetscOptionsName("-pc_asm_in_place","Perform matrix factorization inplace","PCASMSetUseInPlace",&flg);
434:     if (flg) {PCASMSetUseInPlace(pc); }
435:     PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",
436:                             PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);
437:     if (flg) {PCASMSetType(pc,asmtype); }
438:   PetscOptionsTail();
439:   return(0);
440: }

442: /*------------------------------------------------------------------------------------*/

447: PetscErrorCode  PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[])
448: {
449:   PC_ASM         *osm = (PC_ASM*)pc->data;
451:   PetscInt       i;

454:   if (n < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n);
455:   if (pc->setupcalled && (n != osm->n_local_true || is)) {
456:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"PCASMSetLocalSubdomains() should be called before calling PCSetup().");
457:   }
458:   if (!pc->setupcalled) {
459:     if (is) {
460:       for (i=0; i<n; i++) {PetscObjectReference((PetscObject)is[i]);}
461:     }
462:     if (osm->is) {
463:       PCASMDestroySubdomains(osm->n_local_true,osm->is);
464:     }
465:     osm->n_local_true = n;
466:     osm->is           = 0;
467:     if (is) {
468:       PetscMalloc(n*sizeof(IS *),&osm->is);
469:       for (i=0; i<n; i++) { osm->is[i] = is[i]; }
470:     }
471:   }
472:   return(0);
473: }

479: PetscErrorCode  PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is)
480: {
481:   PC_ASM         *osm = (PC_ASM*)pc->data;
483:   PetscMPIInt    rank,size;
484:   PetscInt       n;

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

490:   /*
491:      Split the subdomains equally amoung all processors 
492:   */
493:   MPI_Comm_rank(((PetscObject)pc)->comm,&rank);
494:   MPI_Comm_size(((PetscObject)pc)->comm,&size);
495:   n = N/size + ((N % size) > rank);
496:   if (!n) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"Process %d must have at least one block: total processors %d total blocks %D",(int)rank,(int)size,N);
497:   if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetup().");
498:   if (!pc->setupcalled) {
499:     if (osm->is) {
500:       PCASMDestroySubdomains(osm->n_local_true,osm->is);
501:     }
502:     osm->n_local_true = n;
503:     osm->is           = 0;
504:   }
505:   return(0);
506: }

512: PetscErrorCode  PCASMSetOverlap_ASM(PC pc,PetscInt ovl)
513: {
514:   PC_ASM *osm = (PC_ASM*)pc->data;

517:   if (ovl < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
518:   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetup().");
519:   if (!pc->setupcalled) {
520:     osm->overlap = ovl;
521:   }
522:   return(0);
523: }

529: PetscErrorCode  PCASMSetType_ASM(PC pc,PCASMType type)
530: {
531:   PC_ASM *osm = (PC_ASM*)pc->data;

534:   osm->type     = type;
535:   osm->type_set = PETSC_TRUE;
536:   return(0);
537: }

543: PetscErrorCode  PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
544: {
545:   PC_ASM         *osm = (PC_ASM*)pc->data;

549:   if (osm->n_local_true < 1) {
550:     SETERRQ(PETSC_ERR_ORDER,"Need to call PCSetUP() on PC (or KSPSetUp() on the outer KSP object) before calling here");
551:   }

553:   if (n_local) {
554:     *n_local = osm->n_local_true;
555:   }
556:   if (first_local) {
557:     MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);
558:     *first_local -= osm->n_local_true;
559:   }
560:   if (ksp) {
561:     /* Assume that local solves are now different; not necessarily
562:        true though!  This flag is used only for PCView_ASM() */
563:     *ksp                   = osm->ksp;
564:     osm->same_local_solves = PETSC_FALSE;
565:   }
566:   return(0);
567: }

573: PetscErrorCode  PCASMSetUseInPlace_ASM(PC pc)
574: {
575:   PC_ASM *osm = (PC_ASM*)pc->data;

578:   osm->inplace = PETSC_TRUE;
579:   return(0);
580: }

583: /*----------------------------------------------------------------------------*/
586: /*@
587:    PCASMSetUseInPlace - Tells the system to destroy the matrix after setup is done.

589:    Collective on PC

591:    Input Parameters:
592: .  pc - the preconditioner context

594:    Options Database Key:
595: .  -pc_asm_in_place - Activates in-place factorization

597:    Note:
598:    PCASMSetUseInplace() can only be used with the KSP method KSPPREONLY, and
599:    when the original matrix is not required during the Solve process.
600:    This destroys the matrix, early thus, saving on memory usage.

602:    Level: intermediate

604: .keywords: PC, set, factorization, direct, inplace, in-place, ASM

606: .seealso: PCFactorSetUseInPlace()
607: @*/
608: PetscErrorCode  PCASMSetUseInPlace(PC pc)
609: {
610:   PetscErrorCode ierr,(*f)(PC);

614:   PetscObjectQueryFunction((PetscObject)pc,"PCASMSetUseInPlace_C",(void (**)(void))&f);
615:   if (f) {
616:     (*f)(pc);
617:   }
618:   return(0);
619: }
620: /*----------------------------------------------------------------------------*/

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

628:     Collective on PC 

630:     Input Parameters:
631: +   pc - the preconditioner context
632: .   n - the number of subdomains for this processor (default value = 1)
633: -   is - the index sets that define the subdomains for this processor
634:          (or PETSC_NULL for PETSc to determine subdomains)

636:     Notes:
637:     The IS numbering is in the parallel, global numbering of the vector.

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

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

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

646:     Level: advanced

648: .keywords: PC, ASM, set, local, subdomains, additive Schwarz

650: .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
651:           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains()
652: @*/
653: PetscErrorCode  PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[])
654: {
655:   PetscErrorCode ierr,(*f)(PC,PetscInt,IS[]);

659:   PetscObjectQueryFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",(void (**)(void))&f);
660:   if (f) {
661:     (*f)(pc,n,is);
662:   }
663:   return(0);
664: }

668: /*@C
669:     PCASMSetTotalSubdomains - Sets the subdomains for all processor for the 
670:     additive Schwarz preconditioner.  Either all or no processors in the
671:     PC communicator must call this routine, with the same index sets.

673:     Collective on PC

675:     Input Parameters:
676: +   pc - the preconditioner context
677: .   n - the number of subdomains for all processors
678: -   is - the index sets that define the subdomains for all processor
679:          (or PETSC_NULL for PETSc to determine subdomains)

681:     Options Database Key:
682:     To set the total number of subdomain blocks rather than specify the
683:     index sets, use the option
684: .    -pc_asm_blocks <blks> - Sets total blocks

686:     Notes:
687:     Currently you cannot use this to set the actual subdomains with the argument is.

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

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

694:     Use PCASMSetLocalSubdomains() to set local subdomains.

696:     Level: advanced

698: .keywords: PC, ASM, set, total, global, subdomains, additive Schwarz

700: .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
701:           PCASMCreateSubdomains2D()
702: @*/
703: PetscErrorCode  PCASMSetTotalSubdomains(PC pc,PetscInt N,IS *is)
704: {
705:   PetscErrorCode ierr,(*f)(PC,PetscInt,IS *);

709:   PetscObjectQueryFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",(void (**)(void))&f);
710:   if (f) {
711:     (*f)(pc,N,is);
712:   }
713:   return(0);
714: }

718: /*@
719:     PCASMSetOverlap - Sets the overlap between a pair of subdomains for the
720:     additive Schwarz preconditioner.  Either all or no processors in the
721:     PC communicator must call this routine. 

723:     Collective on PC

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

729:     Options Database Key:
730: .   -pc_asm_overlap <ovl> - Sets overlap

732:     Notes:
733:     By default the ASM preconditioner uses 1 block per processor.  To use
734:     multiple blocks per perocessor, see PCASMSetTotalSubdomains() and
735:     PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>).

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

745:     Note that one can define initial index sets with any overlap via
746:     PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(); the routine
747:     PCASMSetOverlap() merely allows PETSc to extend that overlap further
748:     if desired.

750:     Level: intermediate

752: .keywords: PC, ASM, set, overlap

754: .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
755:           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains()
756: @*/
757: PetscErrorCode  PCASMSetOverlap(PC pc,PetscInt ovl)
758: {
759:   PetscErrorCode ierr,(*f)(PC,PetscInt);

763:   PetscObjectQueryFunction((PetscObject)pc,"PCASMSetOverlap_C",(void (**)(void))&f);
764:   if (f) {
765:     (*f)(pc,ovl);
766:   }
767:   return(0);
768: }

772: /*@
773:     PCASMSetType - Sets the type of restriction and interpolation used
774:     for local problems in the additive Schwarz method.

776:     Collective on PC

778:     Input Parameters:
779: +   pc  - the preconditioner context
780: -   type - variant of ASM, one of
781: .vb
782:       PC_ASM_BASIC       - full interpolation and restriction
783:       PC_ASM_RESTRICT    - full restriction, local processor interpolation
784:       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
785:       PC_ASM_NONE        - local processor restriction and interpolation
786: .ve

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

791:     Level: intermediate

793: .keywords: PC, ASM, set, type

795: .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
796:           PCASMCreateSubdomains2D()
797: @*/
798: PetscErrorCode  PCASMSetType(PC pc,PCASMType type)
799: {
800:   PetscErrorCode ierr,(*f)(PC,PCASMType);

804:   PetscObjectQueryFunction((PetscObject)pc,"PCASMSetType_C",(void (**)(void))&f);
805:   if (f) {
806:     (*f)(pc,type);
807:   }
808:   return(0);
809: }

813: /*@C
814:    PCASMGetSubKSP - Gets the local KSP contexts for all blocks on
815:    this processor.
816:    
817:    Collective on PC iff first_local is requested

819:    Input Parameter:
820: .  pc - the preconditioner context

822:    Output Parameters:
823: +  n_local - the number of blocks on this processor or PETSC_NULL
824: .  first_local - the global number of the first block on this processor or PETSC_NULL,
825:                  all processors must request or all must pass PETSC_NULL
826: -  ksp - the array of KSP contexts

828:    Note:  
829:    After PCASMGetSubKSP() the array of KSPes is not to be freed

831:    Currently for some matrix implementations only 1 block per processor 
832:    is supported.
833:    
834:    You must call KSPSetUp() before calling PCASMGetSubKSP().

836:    Level: advanced

838: .keywords: PC, ASM, additive Schwarz, get, sub, KSP, context

840: .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(),
841:           PCASMCreateSubdomains2D(),
842: @*/
843: PetscErrorCode  PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
844: {
845:   PetscErrorCode ierr,(*f)(PC,PetscInt*,PetscInt*,KSP **);

849:   PetscObjectQueryFunction((PetscObject)pc,"PCASMGetSubKSP_C",(void (**)(void))&f);
850:   if (f) {
851:     (*f)(pc,n_local,first_local,ksp);
852:   } else {
853:     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
854:   }

856:  return(0);
857: }

859: /* -------------------------------------------------------------------------------------*/
860: /*MC
861:    PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 
862:            its own KSP object.

864:    Options Database Keys:
865: +  -pc_asm_truelocal - Activates PCASMSetUseTrueLocal()
866: .  -pc_asm_in_place - Activates in-place factorization
867: .  -pc_asm_blocks <blks> - Sets total blocks
868: .  -pc_asm_overlap <ovl> - Sets overlap
869: -  -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type

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

875:    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
876:      than one processor. Defaults to one block per processor.

878:      To set options on the solvers for each block append -sub_ to all the KSP, and PC
879:         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
880:         
881:      To set the options on the solvers separate for each block call PCASMGetSubKSP()
882:          and set the options directly on the resulting KSP object (you can access its PC
883:          with KSPGetPC())


886:    Level: beginner

888:    Concepts: additive Schwarz method

890:     References:
891:     An additive variant of the Schwarz alternating method for the case of many subregions
892:     M Dryja, OB Widlund - Courant Institute, New York University Technical report

894:     Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 
895:     Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X.

897: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
898:            PCBJACOBI, PCASMSetUseTrueLocal(), PCASMGetSubKSP(), PCASMSetLocalSubdomains(),
899:            PCASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCASMSetOverlap(), PCASMSetType(),
900:            PCASMSetUseInPlace()
901: M*/

906: PetscErrorCode  PCCreate_ASM(PC pc)
907: {
909:   PC_ASM         *osm;

912:   PetscNewLog(pc,PC_ASM,&osm);
913:   osm->n                 = PETSC_DECIDE;
914:   osm->n_local           = 0;
915:   osm->n_local_true      = 0;
916:   osm->overlap           = 1;
917:   osm->ksp               = 0;
918:   osm->scat              = 0;
919:   osm->x                 = 0;
920:   osm->y                 = 0;
921:   osm->is                = 0;
922:   osm->mat               = 0;
923:   osm->pmat              = 0;
924:   osm->type              = PC_ASM_RESTRICT;
925:   osm->same_local_solves = PETSC_TRUE;
926:   osm->inplace           = PETSC_FALSE;

928:   pc->data                   = (void*)osm;
929:   pc->ops->apply             = PCApply_ASM;
930:   pc->ops->applytranspose    = PCApplyTranspose_ASM;
931:   pc->ops->setup             = PCSetUp_ASM;
932:   pc->ops->destroy           = PCDestroy_ASM;
933:   pc->ops->setfromoptions    = PCSetFromOptions_ASM;
934:   pc->ops->setuponblocks     = PCSetUpOnBlocks_ASM;
935:   pc->ops->view              = PCView_ASM;
936:   pc->ops->applyrichardson   = 0;

938:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetLocalSubdomains_C","PCASMSetLocalSubdomains_ASM",
939:                     PCASMSetLocalSubdomains_ASM);
940:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetTotalSubdomains_C","PCASMSetTotalSubdomains_ASM",
941:                     PCASMSetTotalSubdomains_ASM);
942:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetOverlap_C","PCASMSetOverlap_ASM",
943:                     PCASMSetOverlap_ASM);
944:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetType_C","PCASMSetType_ASM",
945:                     PCASMSetType_ASM);
946:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMGetSubKSP_C","PCASMGetSubKSP_ASM",
947:                     PCASMGetSubKSP_ASM);
948:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetUseInPlace_C","PCASMSetUseInPlace_ASM",
949:                     PCASMSetUseInPlace_ASM);
950:   return(0);
951: }


957: /*@C
958:    PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz 
959:    preconditioner for a any problem on a general grid.

961:    Collective

963:    Input Parameters:
964: +  A - The global matrix operator
965: -  n - the number of local blocks

967:    Output Parameters:
968: .  outis - the array of index sets defining the subdomains

970:    Level: advanced

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

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

977: .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid

979: .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains()
980: @*/
981: PetscErrorCode  PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[])
982: {
983:   MatPartitioning           mpart;
984:   const char                *prefix;
985:   PetscErrorCode            (*f)(Mat,PetscTruth*,MatReuse,Mat*);
986:   PetscMPIInt               size;
987:   PetscInt                  i,j,rstart,rend,bs;
988:   PetscTruth                iscopy = PETSC_FALSE,isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
989:   Mat                       Ad = PETSC_NULL, adj;
990:   IS                        ispart,isnumb,*is;
991:   PetscErrorCode            ierr;

996:   if (n < 1) {SETERRQ1(PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);}
997: 
998:   /* Get prefix, row distribution, and block size */
999:   MatGetOptionsPrefix(A,&prefix);
1000:   MatGetOwnershipRange(A,&rstart,&rend);
1001:   MatGetBlockSize(A,&bs);
1002:   if (rstart/bs*bs != rstart || rend/bs*bs != rend) {
1003:     SETERRQ3(PETSC_ERR_ARG_WRONG,"bad row distribution [%D,%D) for matrix block size %D",rstart,rend,bs);
1004:   }
1005:   /* Get diagonal block from matrix if possible */
1006:   MPI_Comm_size(((PetscObject)A)->comm,&size);
1007:   PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",(void (**)(void))&f);
1008:   if (f) {
1009:     (*f)(A,&iscopy,MAT_INITIAL_MATRIX,&Ad);
1010:   } else if (size == 1) {
1011:     iscopy = PETSC_FALSE; Ad = A;
1012:   } else {
1013:     iscopy = PETSC_FALSE; Ad = PETSC_NULL;
1014:   }
1015:   if (Ad) {
1016:     PetscTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);
1017:     if (!isbaij) {PetscTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);}
1018:   }
1019:   if (Ad && n > 1) {
1020:     PetscTruth match,done;
1021:     /* Try to setup a good matrix partitioning if available */
1022:     MatPartitioningCreate(PETSC_COMM_SELF,&mpart);
1023:     PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);
1024:     MatPartitioningSetFromOptions(mpart);
1025:     PetscTypeCompare((PetscObject)mpart,MAT_PARTITIONING_CURRENT,&match);
1026:     if (!match) {
1027:       PetscTypeCompare((PetscObject)mpart,MAT_PARTITIONING_SQUARE,&match);
1028:     }
1029:     if (!match) { /* assume a "good" partitioner is available */
1030:       PetscInt na,*ia,*ja;
1031:       MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);
1032:       if (done) {
1033:         /* Build adjacency matrix by hand. Unfortunately a call to
1034:            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1035:            remove the block-aij structure and we cannot expect
1036:            MatPartitioning to split vertices as we need */
1037:         PetscInt i,j,*row,len,nnz,cnt,*iia=0,*jja=0;
1038:         nnz = 0;
1039:         for (i=0; i<na; i++) { /* count number of nonzeros */
1040:           len = ia[i+1] - ia[i];
1041:           row = ja + ia[i];
1042:           for (j=0; j<len; j++) {
1043:             if (row[j] == i) { /* don't count diagonal */
1044:               len--; break;
1045:             }
1046:           }
1047:           nnz += len;
1048:         }
1049:         PetscMalloc((na+1)*sizeof(PetscInt),&iia);
1050:         PetscMalloc((nnz)*sizeof(PetscInt),&jja);
1051:         nnz    = 0;
1052:         iia[0] = 0;
1053:         for (i=0; i<na; i++) { /* fill adjacency */
1054:           cnt = 0;
1055:           len = ia[i+1] - ia[i];
1056:           row = ja + ia[i];
1057:           for (j=0; j<len; j++) {
1058:             if (row[j] != i) { /* if not diagonal */
1059:               jja[nnz+cnt++] = row[j];
1060:             }
1061:           }
1062:           nnz += cnt;
1063:           iia[i+1] = nnz;
1064:         }
1065:         /* Partitioning of the adjacency matrix */
1066:         MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,PETSC_NULL,&adj);
1067:         MatPartitioningSetAdjacency(mpart,adj);
1068:         MatPartitioningSetNParts(mpart,n);
1069:         MatPartitioningApply(mpart,&ispart);
1070:         ISPartitioningToNumbering(ispart,&isnumb);
1071:         MatDestroy(adj);
1072:         foundpart = PETSC_TRUE;
1073:       }
1074:       MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);
1075:     }
1076:     MatPartitioningDestroy(mpart);
1077:   }
1078:   if (iscopy) {MatDestroy(Ad);}
1079: 
1080:   PetscMalloc(n*sizeof(IS),&is);
1081:   *outis = is;

1083:   if (!foundpart) {

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

1087:     PetscInt mbs   = (rend-rstart)/bs;
1088:     PetscInt start = rstart;
1089:     for (i=0; i<n; i++) {
1090:       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
1091:       ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);
1092:       start += count;
1093:     }
1094: 
1095:   } else {

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

1099:     const PetscInt *numbering;
1100:     PetscInt       *count,nidx,*indices,*newidx,start=0;
1101:     /* Get node count in each partition */
1102:     PetscMalloc(n*sizeof(PetscInt),&count);
1103:     ISPartitioningCount(ispart,n,count);
1104:     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1105:       for (i=0; i<n; i++) count[i] *= bs;
1106:     }
1107:     /* Build indices from node numbering */
1108:     ISGetLocalSize(isnumb,&nidx);
1109:     PetscMalloc(nidx*sizeof(PetscInt),&indices);
1110:     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1111:     ISGetIndices(isnumb,&numbering);
1112:     PetscSortIntWithPermutation(nidx,numbering,indices);
1113:     ISRestoreIndices(isnumb,&numbering);
1114:     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1115:       PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);
1116:       for (i=0; i<nidx; i++)
1117:         for (j=0; j<bs; j++)
1118:           newidx[i*bs+j] = indices[i]*bs + j;
1119:       PetscFree(indices);
1120:       nidx   *= bs;
1121:       indices = newidx;
1122:     }
1123:     /* Shift to get global indices */
1124:     for (i=0; i<nidx; i++) indices[i] += rstart;
1125: 
1126:     /* Build the index sets for each block */
1127:     for (i=0; i<n; i++) {
1128:       ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],&is[i]);
1129:       ISSort(is[i]);
1130:       start += count[i];
1131:     }

1133:     PetscFree(count);
1134:     PetscFree(indices);
1135:     ISDestroy(isnumb);
1136:     ISDestroy(ispart);

1138:   }
1139: 
1140:   return(0);
1141: }

1145: /*@C
1146:    PCASMDestroySubdomains - Destroys the index sets created with
1147:    PCASMCreateSubdomains(). Should be called after setting subdomains
1148:    with PCASMSetLocalSubdomains().

1150:    Collective

1152:    Input Parameters:
1153: +  n - the number of index sets
1154: -  is - the array of index sets

1156:    Level: advanced

1158: .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid

1160: .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains()
1161: @*/
1162: PetscErrorCode  PCASMDestroySubdomains(PetscInt n, IS is[])
1163: {
1164:   PetscInt       i;
1167:   if (n <= 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"n must be > 0: n = %D",n);
1169:   for (i=0; i<n; i++) { ISDestroy(is[i]); }
1170:   PetscFree(is);
1171:   return(0);
1172: }

1176: /*@
1177:    PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 
1178:    preconditioner for a two-dimensional problem on a regular grid.

1180:    Not Collective

1182:    Input Parameters:
1183: +  m, n - the number of mesh points in the x and y directions
1184: .  M, N - the number of subdomains in the x and y directions
1185: .  dof - degrees of freedom per node
1186: -  overlap - overlap in mesh lines

1188:    Output Parameters:
1189: +  Nsub - the number of subdomains created
1190: -  is - the array of index sets defining the subdomains

1192:    Note:
1193:    Presently PCAMSCreateSubdomains2d() is valid only for sequential
1194:    preconditioners.  More general related routines are
1195:    PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains().

1197:    Level: advanced

1199: .keywords: PC, ASM, additive Schwarz, create, subdomains, 2D, regular grid

1201: .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
1202:           PCASMSetOverlap()
1203: @*/
1204: PetscErrorCode  PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is)
1205: {
1206:   PetscInt       i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outter;
1208:   PetscInt       nidx,*idx,loc,ii,jj,count;

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

1213:   *Nsub = N*M;
1214:   PetscMalloc((*Nsub)*sizeof(IS *),is);
1215:   ystart = 0;
1216:   loc_outter = 0;
1217:   for (i=0; i<N; i++) {
1218:     height = n/N + ((n % N) > i); /* height of subdomain */
1219:     if (height < 2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n");
1220:     yleft  = ystart - overlap; if (yleft < 0) yleft = 0;
1221:     yright = ystart + height + overlap; if (yright > n) yright = n;
1222:     xstart = 0;
1223:     for (j=0; j<M; j++) {
1224:       width = m/M + ((m % M) > j); /* width of subdomain */
1225:       if (width < 2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m");
1226:       xleft  = xstart - overlap; if (xleft < 0) xleft = 0;
1227:       xright = xstart + width + overlap; if (xright > m) xright = m;
1228:       nidx   = (xright - xleft)*(yright - yleft);
1229:       PetscMalloc(nidx*sizeof(PetscInt),&idx);
1230:       loc    = 0;
1231:       for (ii=yleft; ii<yright; ii++) {
1232:         count = m*ii + xleft;
1233:         for (jj=xleft; jj<xright; jj++) {
1234:           idx[loc++] = count++;
1235:         }
1236:       }
1237:       ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,(*is)+loc_outter++);
1238:       PetscFree(idx);
1239:       xstart += width;
1240:     }
1241:     ystart += height;
1242:   }
1243:   for (i=0; i<*Nsub; i++) { ISSort((*is)[i]); }
1244:   return(0);
1245: }

1249: /*@C
1250:     PCASMGetLocalSubdomains - Gets the local subdomains (for this processor
1251:     only) for the additive Schwarz preconditioner. 

1253:     Collective on PC 

1255:     Input Parameter:
1256: .   pc - the preconditioner context

1258:     Output Parameters:
1259: +   n - the number of subdomains for this processor (default value = 1)
1260: -   is - the index sets that define the subdomains for this processor
1261:          

1263:     Notes:
1264:     The IS numbering is in the parallel, global numbering of the vector.

1266:     Level: advanced

1268: .keywords: PC, ASM, set, local, subdomains, additive Schwarz

1270: .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
1271:           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices()
1272: @*/
1273: PetscErrorCode  PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[])
1274: {
1275:   PC_ASM         *osm;
1277:   PetscTruth     match;

1283:   PetscTypeCompare((PetscObject)pc,PCASM,&match);
1284:   if (!match) {
1285:     if (n)  *n  = 0;
1286:     if (is) *is = PETSC_NULL;
1287:   } else {
1288:     osm = (PC_ASM*)pc->data;
1289:     if (n)  *n  = osm->n_local_true;
1290:     if (is) *is = osm->is;
1291:   }
1292:   return(0);
1293: }

1297: /*@C
1298:     PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor
1299:     only) for the additive Schwarz preconditioner. 

1301:     Collective on PC 

1303:     Input Parameter:
1304: .   pc - the preconditioner context

1306:     Output Parameters:
1307: +   n - the number of matrices for this processor (default value = 1)
1308: -   mat - the matrices
1309:          

1311:     Level: advanced

1313: .keywords: PC, ASM, set, local, subdomains, additive Schwarz, block Jacobi

1315: .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
1316:           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1317: @*/
1318: PetscErrorCode  PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1319: {
1320:   PC_ASM         *osm;
1322:   PetscTruth     match;

1328:   if (!pc->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1329:   PetscTypeCompare((PetscObject)pc,PCASM,&match);
1330:   if (!match) {
1331:     if (n)   *n   = 0;
1332:     if (mat) *mat = PETSC_NULL;
1333:   } else {
1334:     osm = (PC_ASM*)pc->data;
1335:     if (n)   *n   = osm->n_local_true;
1336:     if (mat) *mat = osm->pmat;
1337:   }
1338:   return(0);
1339: }