Actual source code: gasm.c

petsc-3.13.6 2020-09-29
Report Typos and Errors
  1: /*
  2:   This file defines an "generalized" additive Schwarz preconditioner for any Mat implementation.
  3:   In this version each processor may intersect multiple subdomains and any subdomain may
  4:   intersect multiple processors.  Intersections of subdomains with processors are called *local
  5:   subdomains*.

  7:        N    - total number of distinct global subdomains          (set explicitly in PCGASMSetTotalSubdomains() or implicitly PCGASMSetSubdomains() and then calculated in PCSetUp_GASM())
  8:        n    - actual number of local subdomains on this processor (set in PCGASMSetSubdomains() or calculated in PCGASMSetTotalSubdomains())
  9:        nmax - maximum number of local subdomains per processor    (calculated in PCSetUp_GASM())
 10: */
 11:  #include <petsc/private/pcimpl.h>
 12:  #include <petscdm.h>

 14: typedef struct {
 15:   PetscInt    N,n,nmax;
 16:   PetscInt    overlap;                  /* overlap requested by user */
 17:   PCGASMType  type;                     /* use reduced interpolation, restriction or both */
 18:   PetscBool   type_set;                 /* if user set this value (so won't change it for symmetric problems) */
 19:   PetscBool   same_subdomain_solvers;   /* flag indicating whether all local solvers are same */
 20:   PetscBool   sort_indices;             /* flag to sort subdomain indices */
 21:   PetscBool   user_subdomains;          /* whether the user set explicit subdomain index sets -- keep them on PCReset() */
 22:   PetscBool   dm_subdomains;            /* whether DM is allowed to define subdomains */
 23:   PetscBool   hierarchicalpartitioning;
 24:   IS          *ois;                     /* index sets that define the outer (conceptually, overlapping) subdomains */
 25:   IS          *iis;                     /* index sets that define the inner (conceptually, nonoverlapping) subdomains */
 26:   KSP         *ksp;                     /* linear solvers for each subdomain */
 27:   Mat         *pmat;                    /* subdomain block matrices */
 28:   Vec         gx,gy;                    /* Merged work vectors */
 29:   Vec         *x,*y;                    /* Split work vectors; storage aliases pieces of storage of the above merged vectors. */
 30:   VecScatter  gorestriction;            /* merged restriction to disjoint union of outer subdomains */
 31:   VecScatter  girestriction;            /* merged restriction to disjoint union of inner subdomains */
 32:   VecScatter  pctoouter;
 33:   IS          permutationIS;
 34:   Mat         permutationP;
 35:   Mat         pcmat;
 36:   Vec         pcx,pcy;
 37: } PC_GASM;

 39: static PetscErrorCode  PCGASMComputeGlobalSubdomainNumbering_Private(PC pc,PetscInt **numbering,PetscInt **permutation)
 40: {
 41:   PC_GASM        *osm = (PC_GASM*)pc->data;
 42:   PetscInt       i;

 46:   /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */
 47:   PetscMalloc2(osm->n,numbering,osm->n,permutation);
 48:   PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->iis,NULL,*numbering);
 49:   for (i = 0; i < osm->n; ++i) (*permutation)[i] = i;
 50:   PetscSortIntWithPermutation(osm->n,*numbering,*permutation);
 51:   return(0);
 52: }

 54: static PetscErrorCode  PCGASMSubdomainView_Private(PC pc, PetscInt i, PetscViewer viewer)
 55: {
 56:   PC_GASM        *osm = (PC_GASM*)pc->data;
 57:   PetscInt       j,nidx;
 58:   const PetscInt *idx;
 59:   PetscViewer    sviewer;
 60:   char           *cidx;

 64:   if (i < 0 || i > osm->n) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Invalid subdomain %D: must nonnegative and less than %D", i, osm->n);
 65:   /* Inner subdomains. */
 66:   ISGetLocalSize(osm->iis[i], &nidx);
 67:   /*
 68:    No more than 15 characters per index plus a space.
 69:    PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
 70:    in case nidx == 0. That will take care of the space for the trailing '\0' as well.
 71:    For nidx == 0, the whole string 16 '\0'.
 72:    */
 73: #define len  16*(nidx+1)+1
 74:   PetscMalloc1(len, &cidx);
 75:   PetscViewerStringOpen(PETSC_COMM_SELF, cidx, len, &sviewer);
 76: #undef len
 77:   ISGetIndices(osm->iis[i], &idx);
 78:   for (j = 0; j < nidx; ++j) {
 79:     PetscViewerStringSPrintf(sviewer, "%D ", idx[j]);
 80:   }
 81:   ISRestoreIndices(osm->iis[i],&idx);
 82:   PetscViewerDestroy(&sviewer);
 83:   PetscViewerASCIIPrintf(viewer, "Inner subdomain:\n");
 84:   PetscViewerFlush(viewer);
 85:   PetscViewerASCIIPushSynchronized(viewer);
 86:   PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);
 87:   PetscViewerFlush(viewer);
 88:   PetscViewerASCIIPopSynchronized(viewer);
 89:   PetscViewerASCIIPrintf(viewer, "\n");
 90:   PetscViewerFlush(viewer);
 91:   PetscFree(cidx);
 92:   /* Outer subdomains. */
 93:   ISGetLocalSize(osm->ois[i], &nidx);
 94:   /*
 95:    No more than 15 characters per index plus a space.
 96:    PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
 97:    in case nidx == 0. That will take care of the space for the trailing '\0' as well.
 98:    For nidx == 0, the whole string 16 '\0'.
 99:    */
100: #define len  16*(nidx+1)+1
101:   PetscMalloc1(len, &cidx);
102:   PetscViewerStringOpen(PETSC_COMM_SELF, cidx, len, &sviewer);
103: #undef len
104:   ISGetIndices(osm->ois[i], &idx);
105:   for (j = 0; j < nidx; ++j) {
106:     PetscViewerStringSPrintf(sviewer,"%D ", idx[j]);
107:   }
108:   PetscViewerDestroy(&sviewer);
109:   ISRestoreIndices(osm->ois[i],&idx);
110:   PetscViewerASCIIPrintf(viewer, "Outer subdomain:\n");
111:   PetscViewerFlush(viewer);
112:   PetscViewerASCIIPushSynchronized(viewer);
113:   PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);
114:   PetscViewerFlush(viewer);
115:   PetscViewerASCIIPopSynchronized(viewer);
116:   PetscViewerASCIIPrintf(viewer, "\n");
117:   PetscViewerFlush(viewer);
118:   PetscFree(cidx);
119:   return(0);
120: }

122: static PetscErrorCode  PCGASMPrintSubdomains(PC pc)
123: {
124:   PC_GASM        *osm = (PC_GASM*)pc->data;
125:   const char     *prefix;
126:   char           fname[PETSC_MAX_PATH_LEN+1];
127:   PetscInt       l, d, count;
128:   PetscBool      doprint,found;
129:   PetscViewer    viewer, sviewer = NULL;
130:   PetscInt       *numbering,*permutation;/* global numbering of locally-supported subdomains and the permutation from the local ordering */

134:   PCGetOptionsPrefix(pc,&prefix);
135:   doprint  = PETSC_FALSE;
136:   PetscOptionsGetBool(NULL,prefix,"-pc_gasm_print_subdomains",&doprint,NULL);
137:   if (!doprint) return(0);
138:   PetscOptionsGetString(NULL,prefix,"-pc_gasm_print_subdomains",fname,PETSC_MAX_PATH_LEN,&found);
139:   if (!found) { PetscStrcpy(fname,"stdout"); };
140:   PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);
141:   /*
142:    Make sure the viewer has a name. Otherwise this may cause a deadlock or other weird errors when creating a subcomm viewer:
143:    the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed collectively on the comm.
144:   */
145:   PetscObjectName((PetscObject)viewer);
146:   l    = 0;
147:   PCGASMComputeGlobalSubdomainNumbering_Private(pc,&numbering,&permutation);
148:   for (count = 0; count < osm->N; ++count) {
149:     /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */
150:     if (l<osm->n) {
151:       d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */
152:       if (numbering[d] == count) {
153:         PetscViewerGetSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);
154:         PCGASMSubdomainView_Private(pc,d,sviewer);
155:         PetscViewerRestoreSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);
156:         ++l;
157:       }
158:     }
159:     MPI_Barrier(PetscObjectComm((PetscObject)pc));
160:   }
161:   PetscFree2(numbering,permutation);
162:   PetscViewerDestroy(&viewer);
163:   return(0);
164: }


167: static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer)
168: {
169:   PC_GASM        *osm = (PC_GASM*)pc->data;
170:   const char     *prefix;
172:   PetscMPIInt    rank, size;
173:   PetscInt       bsz;
174:   PetscBool      iascii,view_subdomains=PETSC_FALSE;
175:   PetscViewer    sviewer;
176:   PetscInt       count, l;
177:   char           overlap[256]     = "user-defined overlap";
178:   char           gsubdomains[256] = "unknown total number of subdomains";
179:   char           msubdomains[256] = "unknown max number of local subdomains";
180:   PetscInt       *numbering,*permutation;/* global numbering of locally-supported subdomains and the permutation from the local ordering */

183:   MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);
184:   MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);

186:   if (osm->overlap >= 0) {
187:     PetscSNPrintf(overlap,sizeof(overlap),"requested amount of overlap = %D",osm->overlap);
188:   }
189:   if (osm->N != PETSC_DETERMINE) {
190:     PetscSNPrintf(gsubdomains, sizeof(gsubdomains), "total number of subdomains = %D",osm->N);
191:   }
192:   if (osm->nmax != PETSC_DETERMINE) {
193:     PetscSNPrintf(msubdomains,sizeof(msubdomains),"max number of local subdomains = %D",osm->nmax);
194:   }

196:   PCGetOptionsPrefix(pc,&prefix);
197:   PetscOptionsGetBool(NULL,prefix,"-pc_gasm_view_subdomains",&view_subdomains,NULL);

199:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
200:   if (iascii) {
201:     /*
202:      Make sure the viewer has a name. Otherwise this may cause a deadlock when creating a subcomm viewer:
203:      the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed
204:      collectively on the comm.
205:      */
206:     PetscObjectName((PetscObject)viewer);
207:     PetscViewerASCIIPrintf(viewer,"  Restriction/interpolation type: %s\n",PCGASMTypes[osm->type]);
208:     PetscViewerASCIIPrintf(viewer,"  %s\n",overlap);
209:     PetscViewerASCIIPrintf(viewer,"  %s\n",gsubdomains);
210:     PetscViewerASCIIPrintf(viewer,"  %s\n",msubdomains);
211:     PetscViewerASCIIPushSynchronized(viewer);
212:     PetscViewerASCIISynchronizedPrintf(viewer,"  [%d|%d] number of locally-supported subdomains = %D\n",rank,size,osm->n);
213:     PetscViewerFlush(viewer);
214:     PetscViewerASCIIPopSynchronized(viewer);
215:     /* Cannot take advantage of osm->same_subdomain_solvers without a global numbering of subdomains. */
216:     PetscViewerASCIIPrintf(viewer,"  Subdomain solver info is as follows:\n");
217:     PetscViewerASCIIPushTab(viewer);
218:     PetscViewerASCIIPrintf(viewer,"  - - - - - - - - - - - - - - - - - -\n");
219:     /* Make sure that everybody waits for the banner to be printed. */
220:     MPI_Barrier(PetscObjectComm((PetscObject)viewer));
221:     /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */
222:     PCGASMComputeGlobalSubdomainNumbering_Private(pc,&numbering,&permutation);
223:     l = 0;
224:     for (count = 0; count < osm->N; ++count) {
225:       PetscMPIInt srank, ssize;
226:       if (l<osm->n) {
227:         PetscInt d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */
228:         if (numbering[d] == count) {
229:           MPI_Comm_size(((PetscObject)osm->ois[d])->comm, &ssize);
230:           MPI_Comm_rank(((PetscObject)osm->ois[d])->comm, &srank);
231:           PetscViewerGetSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);
232:           ISGetLocalSize(osm->ois[d],&bsz);
233:           PetscViewerASCIISynchronizedPrintf(sviewer,"  [%d|%d] (subcomm [%d|%d]) local subdomain number %D, local size = %D\n",rank,size,srank,ssize,d,bsz);
234:           PetscViewerFlush(sviewer);
235:           if (view_subdomains) {
236:             PCGASMSubdomainView_Private(pc,d,sviewer);
237:           }
238:           if (!pc->setupcalled) {
239:             PetscViewerASCIIPrintf(sviewer, "  Solver not set up yet: PCSetUp() not yet called\n");
240:           } else {
241:             KSPView(osm->ksp[d],sviewer);
242:           }
243:           PetscViewerASCIIPrintf(sviewer,"  - - - - - - - - - - - - - - - - - -\n");
244:           PetscViewerFlush(sviewer);
245:           PetscViewerRestoreSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);
246:           ++l;
247:         }
248:       }
249:       MPI_Barrier(PetscObjectComm((PetscObject)pc));
250:     }
251:     PetscFree2(numbering,permutation);
252:     PetscViewerASCIIPopTab(viewer);
253:   }
254:   PetscViewerFlush(viewer);
255:   PetscViewerASCIIPopSynchronized(viewer);
256:   return(0);
257: }

259: PETSC_INTERN PetscErrorCode  PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]);



263: PetscErrorCode PCGASMSetHierarchicalPartitioning(PC pc)
264: {
265:    PC_GASM              *osm = (PC_GASM*)pc->data;
266:    MatPartitioning       part;
267:    MPI_Comm              comm;
268:    PetscMPIInt           size;
269:    PetscInt              nlocalsubdomains,fromrows_localsize;
270:    IS                    partitioning,fromrows,isn;
271:    Vec                   outervec;
272:    PetscErrorCode        ierr;

275:    PetscObjectGetComm((PetscObject)pc,&comm);
276:    MPI_Comm_size(comm,&size);
277:    /* we do not need a hierarchical partitioning when
278:     * the total number of subdomains is consistent with
279:     * the number of MPI tasks.
280:     * For the following cases, we do not need to use HP
281:     * */
282:    if(osm->N==PETSC_DETERMINE || osm->N>=size || osm->N==1) return(0);
283:    if(size%osm->N != 0) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"have to specify the total number of subdomains %D to be a factor of the number of processors %d \n",osm->N,size);
284:    nlocalsubdomains = size/osm->N;
285:    osm->n           = 1;
286:    MatPartitioningCreate(comm,&part);
287:    MatPartitioningSetAdjacency(part,pc->pmat);
288:    MatPartitioningSetType(part,MATPARTITIONINGHIERARCH);
289:    MatPartitioningHierarchicalSetNcoarseparts(part,osm->N);
290:    MatPartitioningHierarchicalSetNfineparts(part,nlocalsubdomains);
291:    MatPartitioningSetFromOptions(part);
292:    /* get new processor owner number of each vertex */
293:    MatPartitioningApply(part,&partitioning);
294:    ISBuildTwoSided(partitioning,NULL,&fromrows);
295:    ISPartitioningToNumbering(partitioning,&isn);
296:    ISDestroy(&isn);
297:    ISGetLocalSize(fromrows,&fromrows_localsize);
298:    MatPartitioningDestroy(&part);
299:    MatCreateVecs(pc->pmat,&outervec,NULL);
300:    VecCreateMPI(comm,fromrows_localsize,PETSC_DETERMINE,&(osm->pcx));
301:    VecDuplicate(osm->pcx,&(osm->pcy));
302:    VecScatterCreate(osm->pcx,NULL,outervec,fromrows,&(osm->pctoouter));
303:    MatCreateSubMatrix(pc->pmat,fromrows,fromrows,MAT_INITIAL_MATRIX,&(osm->permutationP));
304:    PetscObjectReference((PetscObject)fromrows);
305:    osm->permutationIS = fromrows;
306:    osm->pcmat =  pc->pmat;
307:    PetscObjectReference((PetscObject)osm->permutationP);
308:    pc->pmat = osm->permutationP;
309:    VecDestroy(&outervec);
310:    ISDestroy(&fromrows);
311:    ISDestroy(&partitioning);
312:    osm->n           = PETSC_DETERMINE;
313:    return(0);
314: }



318: static PetscErrorCode PCSetUp_GASM(PC pc)
319: {
320:   PC_GASM        *osm = (PC_GASM*)pc->data;
322:   PetscInt       i,nInnerIndices,nTotalInnerIndices;
323:   PetscMPIInt    rank, size;
324:   MatReuse       scall = MAT_REUSE_MATRIX;
325:   KSP            ksp;
326:   PC             subpc;
327:   const char     *prefix,*pprefix;
328:   Vec            x,y;
329:   PetscInt       oni;       /* Number of indices in the i-th local outer subdomain.               */
330:   const PetscInt *oidxi;    /* Indices from the i-th subdomain local outer subdomain.             */
331:   PetscInt       on;        /* Number of indices in the disjoint union of local outer subdomains. */
332:   PetscInt       *oidx;     /* Indices in the disjoint union of local outer subdomains. */
333:   IS             gois;      /* Disjoint union the global indices of outer subdomains.             */
334:   IS             goid;      /* Identity IS of the size of the disjoint union of outer subdomains. */
335:   PetscScalar    *gxarray, *gyarray;
336:   PetscInt       gostart;   /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- over the disjoint union of outer subdomains. */
337:   PetscInt       num_subdomains    = 0;
338:   DM             *subdomain_dm     = NULL;
339:   char           **subdomain_names = NULL;
340:   PetscInt       *numbering;


344:   MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
345:   MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
346:   if (!pc->setupcalled) {
347:         /* use a hierarchical partitioning */
348:     if(osm->hierarchicalpartitioning){
349:       PCGASMSetHierarchicalPartitioning(pc);
350:     }
351:     if (osm->n == PETSC_DETERMINE) {
352:       if (osm->N != PETSC_DETERMINE) {
353:            /* No local subdomains given, but the desired number of total subdomains is known, so construct them accordingly. */
354:            PCGASMCreateSubdomains(pc->pmat,osm->N,&osm->n,&osm->iis);
355:       } else if (osm->dm_subdomains && pc->dm) {
356:         /* try pc->dm next, if allowed */
357:         PetscInt  d;
358:         IS       *inner_subdomain_is, *outer_subdomain_is;
359:         DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm);
360:         if (num_subdomains) {
361:           PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is);
362:         }
363:         for (d = 0; d < num_subdomains; ++d) {
364:           if (inner_subdomain_is) {ISDestroy(&inner_subdomain_is[d]);}
365:           if (outer_subdomain_is) {ISDestroy(&outer_subdomain_is[d]);}
366:         }
367:         PetscFree(inner_subdomain_is);
368:         PetscFree(outer_subdomain_is);
369:       } else {
370:         /* still no subdomains; use one per processor */
371:         osm->nmax = osm->n = 1;
372:         MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
373:         osm->N    = size;
374:         PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);
375:       }
376:     }
377:     if (!osm->iis) {
378:       /*
379:        osm->n was set in PCGASMSetSubdomains(), but the actual subdomains have not been supplied.
380:        We create the requisite number of local inner subdomains and then expand them into
381:        out subdomains, if necessary.
382:        */
383:       PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);
384:     }
385:     if (!osm->ois) {
386:       /*
387:             Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap
388:             has been requested, copy the inner subdomains over so they can be modified.
389:       */
390:       PetscMalloc1(osm->n,&osm->ois);
391:       for (i=0; i<osm->n; ++i) {
392:         if (osm->overlap > 0 && osm->N>1) { /* With positive overlap, osm->iis[i] will be modified */
393:           ISDuplicate(osm->iis[i],(osm->ois)+i);
394:           ISCopy(osm->iis[i],osm->ois[i]);
395:         } else {
396:           PetscObjectReference((PetscObject)((osm->iis)[i]));
397:           osm->ois[i] = osm->iis[i];
398:         }
399:       }
400:       if (osm->overlap>0 && osm->N>1) {
401:            /* Extend the "overlapping" regions by a number of steps */
402:            MatIncreaseOverlapSplit(pc->pmat,osm->n,osm->ois,osm->overlap);
403:       }
404:     }

406:     /* Now the subdomains are defined.  Determine their global and max local numbers, if necessary. */
407:     if (osm->nmax == PETSC_DETERMINE) {
408:       PetscMPIInt inwork,outwork;
409:       /* determine global number of subdomains and the max number of local subdomains */
410:       inwork = osm->n;
411:       MPIU_Allreduce(&inwork,&outwork,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)pc));
412:       osm->nmax  = outwork;
413:     }
414:     if (osm->N == PETSC_DETERMINE) {
415:       /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */
416:       PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,&osm->N,NULL);
417:     }


420:     if (osm->sort_indices) {
421:       for (i=0; i<osm->n; i++) {
422:         ISSort(osm->ois[i]);
423:         ISSort(osm->iis[i]);
424:       }
425:     }
426:     PCGetOptionsPrefix(pc,&prefix);
427:     PCGASMPrintSubdomains(pc);

429:     /*
430:        Merge the ISs, create merged vectors and restrictions.
431:      */
432:     /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */
433:     on = 0;
434:     for (i=0; i<osm->n; i++) {
435:       ISGetLocalSize(osm->ois[i],&oni);
436:       on  += oni;
437:     }
438:     PetscMalloc1(on, &oidx);
439:     on   = 0;
440:     /* Merge local indices together */
441:     for (i=0; i<osm->n; i++) {
442:       ISGetLocalSize(osm->ois[i],&oni);
443:       ISGetIndices(osm->ois[i],&oidxi);
444:       PetscArraycpy(oidx+on,oidxi,oni);
445:       ISRestoreIndices(osm->ois[i],&oidxi);
446:       on  += oni;
447:     }
448:     ISCreateGeneral(((PetscObject)(pc))->comm,on,oidx,PETSC_OWN_POINTER,&gois);
449:     nTotalInnerIndices = 0;
450:     for(i=0; i<osm->n; i++){
451:       ISGetLocalSize(osm->iis[i],&nInnerIndices);
452:       nTotalInnerIndices += nInnerIndices;
453:     }
454:     VecCreateMPI(((PetscObject)(pc))->comm,nTotalInnerIndices,PETSC_DETERMINE,&x);
455:     VecDuplicate(x,&y);

457:     VecCreateMPI(PetscObjectComm((PetscObject)pc),on,PETSC_DECIDE,&osm->gx);
458:     VecDuplicate(osm->gx,&osm->gy);
459:     VecGetOwnershipRange(osm->gx, &gostart, NULL);
460:     ISCreateStride(PetscObjectComm((PetscObject)pc),on,gostart,1, &goid);
461:     /* gois might indices not on local */
462:     VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction));
463:     PetscMalloc1(osm->n,&numbering);
464:     PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,NULL,numbering);
465:     VecDestroy(&x);
466:     ISDestroy(&gois);

468:     /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */
469:     {
470:       PetscInt        ini;           /* Number of indices the i-th a local inner subdomain. */
471:       PetscInt        in;            /* Number of indices in the disjoint union of local inner subdomains. */
472:       PetscInt       *iidx;          /* Global indices in the merged local inner subdomain. */
473:       PetscInt       *ioidx;         /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
474:       IS              giis;          /* IS for the disjoint union of inner subdomains. */
475:       IS              giois;         /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
476:       PetscScalar    *array;
477:       const PetscInt *indices;
478:       PetscInt        k;
479:       on = 0;
480:       for (i=0; i<osm->n; i++) {
481:         ISGetLocalSize(osm->ois[i],&oni);
482:         on  += oni;
483:       }
484:       PetscMalloc1(on, &iidx);
485:       PetscMalloc1(on, &ioidx);
486:       VecGetArray(y,&array);
487:       /* set communicator id to determine where overlap is */
488:       in   = 0;
489:       for (i=0; i<osm->n; i++) {
490:         ISGetLocalSize(osm->iis[i],&ini);
491:         for (k = 0; k < ini; ++k){
492:           array[in+k] = numbering[i];
493:         }
494:         in += ini;
495:       }
496:       VecRestoreArray(y,&array);
497:       VecScatterBegin(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD);
498:       VecScatterEnd(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD);
499:       VecGetOwnershipRange(osm->gy,&gostart, NULL);
500:       VecGetArray(osm->gy,&array);
501:       on  = 0;
502:       in  = 0;
503:       for (i=0; i<osm->n; i++) {
504:             ISGetLocalSize(osm->ois[i],&oni);
505:             ISGetIndices(osm->ois[i],&indices);
506:             for (k=0; k<oni; k++) {
507:           /*  skip overlapping indices to get inner domain */
508:           if(PetscRealPart(array[on+k]) != numbering[i]) continue;
509:           iidx[in]    = indices[k];
510:           ioidx[in++] = gostart+on+k;
511:             }
512:             ISRestoreIndices(osm->ois[i], &indices);
513:             on += oni;
514:       }
515:       VecRestoreArray(osm->gy,&array);
516:       ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,iidx,PETSC_OWN_POINTER,&giis);
517:       ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,ioidx,PETSC_OWN_POINTER,&giois);
518:       VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction);
519:       VecDestroy(&y);
520:       ISDestroy(&giis);
521:       ISDestroy(&giois);
522:     }
523:     ISDestroy(&goid);
524:     PetscFree(numbering);

526:     /* Create the subdomain work vectors. */
527:     PetscMalloc1(osm->n,&osm->x);
528:     PetscMalloc1(osm->n,&osm->y);
529:     VecGetArray(osm->gx, &gxarray);
530:     VecGetArray(osm->gy, &gyarray);
531:     for (i=0, on=0; i<osm->n; ++i, on += oni) {
532:       PetscInt oNi;
533:       ISGetLocalSize(osm->ois[i],&oni);
534:       /* on a sub communicator */
535:       ISGetSize(osm->ois[i],&oNi);
536:       VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i]);
537:       VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i]);
538:     }
539:     VecRestoreArray(osm->gx, &gxarray);
540:     VecRestoreArray(osm->gy, &gyarray);
541:     /* Create the subdomain solvers */
542:     PetscMalloc1(osm->n,&osm->ksp);
543:     for (i=0; i<osm->n; i++) {
544:       char subprefix[PETSC_MAX_PATH_LEN+1];
545:       KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp);
546:       KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);
547:       PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);
548:       PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);
549:       KSPSetType(ksp,KSPPREONLY);
550:       KSPGetPC(ksp,&subpc); /* Why do we need this here? */
551:       if (subdomain_dm) {
552:             KSPSetDM(ksp,subdomain_dm[i]);
553:             DMDestroy(subdomain_dm+i);
554:       }
555:       PCGetOptionsPrefix(pc,&prefix);
556:       KSPSetOptionsPrefix(ksp,prefix);
557:       if (subdomain_names && subdomain_names[i]) {
558:              PetscSNPrintf(subprefix,PETSC_MAX_PATH_LEN,"sub_%s_",subdomain_names[i]);
559:              KSPAppendOptionsPrefix(ksp,subprefix);
560:              PetscFree(subdomain_names[i]);
561:       }
562:       KSPAppendOptionsPrefix(ksp,"sub_");
563:       osm->ksp[i] = ksp;
564:     }
565:     PetscFree(subdomain_dm);
566:     PetscFree(subdomain_names);
567:     scall = MAT_INITIAL_MATRIX;

569:   } else { /* if (pc->setupcalled) */
570:     /*
571:        Destroy the submatrices from the previous iteration
572:     */
573:     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
574:       MatDestroyMatrices(osm->n,&osm->pmat);
575:       scall = MAT_INITIAL_MATRIX;
576:     }
577:     if(osm->permutationIS){
578:       MatCreateSubMatrix(pc->pmat,osm->permutationIS,osm->permutationIS,scall,&osm->permutationP);
579:       PetscObjectReference((PetscObject)osm->permutationP);
580:       osm->pcmat = pc->pmat;
581:       pc->pmat   = osm->permutationP;
582:     }

584:   }


587:   /*
588:      Extract out the submatrices.
589:   */
590:   if (size > 1) {
591:     MatCreateSubMatricesMPI(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);
592:   } else {
593:     MatCreateSubMatrices(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);
594:   }
595:   if (scall == MAT_INITIAL_MATRIX) {
596:     PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);
597:     for (i=0; i<osm->n; i++) {
598:       PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);
599:       PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);
600:     }
601:   }

603:   /* Return control to the user so that the submatrices can be modified (e.g., to apply
604:      different boundary conditions for the submatrices than for the global problem) */
605:   PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP);

607:   /*
608:      Loop over submatrices putting them into local ksps
609:   */
610:   for (i=0; i<osm->n; i++) {
611:     KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);
612:     if (!pc->setupcalled) {
613:       KSPSetFromOptions(osm->ksp[i]);
614:     }
615:   }
616:   if(osm->pcmat){
617:     MatDestroy(&pc->pmat);
618:     pc->pmat   = osm->pcmat;
619:     osm->pcmat = 0;
620:   }
621:   return(0);
622: }

624: static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc)
625: {
626:   PC_GASM        *osm = (PC_GASM*)pc->data;
628:   PetscInt       i;

631:   for (i=0; i<osm->n; i++) {
632:     KSPSetUp(osm->ksp[i]);
633:   }
634:   return(0);
635: }

637: static PetscErrorCode PCApply_GASM(PC pc,Vec xin,Vec yout)
638: {
639:   PC_GASM        *osm = (PC_GASM*)pc->data;
641:   PetscInt       i;
642:   Vec            x,y;
643:   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;

646:   if(osm->pctoouter){
647:     VecScatterBegin(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);
648:     VecScatterEnd(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);
649:     x = osm->pcx;
650:     y = osm->pcy;
651:   }else{
652:         x = xin;
653:         y = yout;
654:   }
655:   /*
656:      Support for limiting the restriction or interpolation only to the inner
657:      subdomain values (leaving the other values 0).
658:   */
659:   if (!(osm->type & PC_GASM_RESTRICT)) {
660:     /* have to zero the work RHS since scatter may leave some slots empty */
661:     VecZeroEntries(osm->gx);
662:     VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);
663:   } else {
664:     VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);
665:   }
666:   VecZeroEntries(osm->gy);
667:   if (!(osm->type & PC_GASM_RESTRICT)) {
668:     VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);
669:   } else {
670:     VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);
671:   }
672:   /* do the subdomain solves */
673:   for (i=0; i<osm->n; ++i) {
674:     KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);
675:     KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);
676:   }
677:   /* Do we need to zero y ?? */
678:   VecZeroEntries(y);
679:   if (!(osm->type & PC_GASM_INTERPOLATE)) {
680:     VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);
681:     VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);
682:   } else {
683:     VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);
684:     VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);
685:   }
686:   if(osm->pctoouter){
687:     VecScatterBegin(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);
688:     VecScatterEnd(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);
689:   }
690:   return(0);
691: }

693: static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec xin,Vec yout)
694: {
695:   PC_GASM        *osm = (PC_GASM*)pc->data;
697:   PetscInt       i;
698:   Vec            x,y;
699:   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;

702:   if(osm->pctoouter){
703:    VecScatterBegin(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);
704:    VecScatterEnd(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);
705:    x = osm->pcx;
706:    y = osm->pcy;
707:   }else{
708:         x = xin;
709:         y = yout;
710:   }
711:   /*
712:      Support for limiting the restriction or interpolation to only local
713:      subdomain values (leaving the other values 0).

715:      Note: these are reversed from the PCApply_GASM() because we are applying the
716:      transpose of the three terms
717:   */
718:   if (!(osm->type & PC_GASM_INTERPOLATE)) {
719:     /* have to zero the work RHS since scatter may leave some slots empty */
720:     VecZeroEntries(osm->gx);
721:     VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);
722:   } else {
723:     VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);
724:   }
725:   VecZeroEntries(osm->gy);
726:   if (!(osm->type & PC_GASM_INTERPOLATE)) {
727:     VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);
728:   } else {
729:     VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);
730:   }
731:   /* do the local solves */
732:   for (i=0; i<osm->n; ++i) { /* Note that the solves are local, so we can go to osm->n, rather than osm->nmax. */
733:     KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);
734:     KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);
735:   }
736:   VecZeroEntries(y);
737:   if (!(osm->type & PC_GASM_RESTRICT)) {
738:     VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);
739:     VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);
740:   } else {
741:     VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);
742:     VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);
743:   }
744:   if(osm->pctoouter){
745:    VecScatterBegin(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);
746:    VecScatterEnd(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);
747:   }
748:   return(0);
749: }

751: static PetscErrorCode PCReset_GASM(PC pc)
752: {
753:   PC_GASM        *osm = (PC_GASM*)pc->data;
755:   PetscInt       i;

758:   if (osm->ksp) {
759:     for (i=0; i<osm->n; i++) {
760:       KSPReset(osm->ksp[i]);
761:     }
762:   }
763:   if (osm->pmat) {
764:     if (osm->n > 0) {
765:       PetscMPIInt size;
766:       MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
767:       if (size > 1) {
768:         /* osm->pmat is created by MatCreateSubMatricesMPI(), cannot use MatDestroySubMatrices() */
769:         MatDestroyMatrices(osm->n,&osm->pmat);
770:       } else {
771:         MatDestroySubMatrices(osm->n,&osm->pmat);
772:       }
773:     }
774:   }
775:   if (osm->x) {
776:     for (i=0; i<osm->n; i++) {
777:       VecDestroy(&osm->x[i]);
778:       VecDestroy(&osm->y[i]);
779:     }
780:   }
781:   VecDestroy(&osm->gx);
782:   VecDestroy(&osm->gy);

784:   VecScatterDestroy(&osm->gorestriction);
785:   VecScatterDestroy(&osm->girestriction);
786:   if (!osm->user_subdomains) {
787:     PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);
788:     osm->N    = PETSC_DETERMINE;
789:     osm->nmax = PETSC_DETERMINE;
790:   }
791:   if(osm->pctoouter){
792:         VecScatterDestroy(&(osm->pctoouter));
793:   }
794:   if(osm->permutationIS){
795:         ISDestroy(&(osm->permutationIS));
796:   }
797:   if(osm->pcx){
798:         VecDestroy(&(osm->pcx));
799:   }
800:   if(osm->pcy){
801:         VecDestroy(&(osm->pcy));
802:   }
803:   if(osm->permutationP){
804:     MatDestroy(&(osm->permutationP));
805:   }
806:   if(osm->pcmat){
807:         MatDestroy(&osm->pcmat);
808:   }
809:   return(0);
810: }

812: static PetscErrorCode PCDestroy_GASM(PC pc)
813: {
814:   PC_GASM        *osm = (PC_GASM*)pc->data;
816:   PetscInt       i;

819:   PCReset_GASM(pc);
820:   /* PCReset will not destroy subdomains, if user_subdomains is true. */
821:   PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);
822:   if (osm->ksp) {
823:     for (i=0; i<osm->n; i++) {
824:       KSPDestroy(&osm->ksp[i]);
825:     }
826:     PetscFree(osm->ksp);
827:   }
828:   PetscFree(osm->x);
829:   PetscFree(osm->y);
830:   PetscFree(pc->data);
831:   return(0);
832: }

834: static PetscErrorCode PCSetFromOptions_GASM(PetscOptionItems *PetscOptionsObject,PC pc)
835: {
836:   PC_GASM        *osm = (PC_GASM*)pc->data;
838:   PetscInt       blocks,ovl;
839:   PetscBool      flg;
840:   PCGASMType     gasmtype;

843:   PetscOptionsHead(PetscOptionsObject,"Generalized additive Schwarz options");
844:   PetscOptionsBool("-pc_gasm_use_dm_subdomains","If subdomains aren't set, use DMCreateDomainDecomposition() to define subdomains.","PCGASMSetUseDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);
845:   PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->N,&blocks,&flg);
846:   if (flg) {
847:     PCGASMSetTotalSubdomains(pc,blocks);
848:   }
849:   PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);
850:   if (flg) {
851:     PCGASMSetOverlap(pc,ovl);
852:     osm->dm_subdomains = PETSC_FALSE;
853:   }
854:   flg  = PETSC_FALSE;
855:   PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);
856:   if (flg) {PCGASMSetType(pc,gasmtype);}
857:   PetscOptionsBool("-pc_gasm_use_hierachical_partitioning","use hierarchical partitioning",NULL,osm->hierarchicalpartitioning,&osm->hierarchicalpartitioning,&flg);
858:   PetscOptionsTail();
859:   return(0);
860: }

862: /*------------------------------------------------------------------------------------*/

864: /*@
865:     PCGASMSetTotalSubdomains - sets the total number of subdomains to use across the
866:                                communicator.
867:     Logically collective on pc

869:     Input Parameters:
870: +   pc  - the preconditioner
871: -   N   - total number of subdomains


874:     Level: beginner

876: .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap()
877:           PCGASMCreateSubdomains2D()
878: @*/
879: PetscErrorCode  PCGASMSetTotalSubdomains(PC pc,PetscInt N)
880: {
881:   PC_GASM        *osm = (PC_GASM*)pc->data;
882:   PetscMPIInt    size,rank;

886:   if (N < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Total number of subdomains must be 1 or more, got N = %D",N);
887:   if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before calling PCSetUp().");

889:   PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);
890:   osm->ois = osm->iis = NULL;

892:   MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
893:   MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
894:   osm->N    = N;
895:   osm->n    = PETSC_DETERMINE;
896:   osm->nmax = PETSC_DETERMINE;
897:   osm->dm_subdomains = PETSC_FALSE;
898:   return(0);
899: }


902: static PetscErrorCode  PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[])
903: {
904:   PC_GASM         *osm = (PC_GASM*)pc->data;
905:   PetscErrorCode  ierr;
906:   PetscInt        i;

909:   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, got n = %D",n);
910:   if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp().");

912:   PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);
913:   osm->iis  = osm->ois = NULL;
914:   osm->n    = n;
915:   osm->N    = PETSC_DETERMINE;
916:   osm->nmax = PETSC_DETERMINE;
917:   if (ois) {
918:     PetscMalloc1(n,&osm->ois);
919:     for (i=0; i<n; i++) {
920:       PetscObjectReference((PetscObject)ois[i]);
921:       osm->ois[i] = ois[i];
922:     }
923:     /*
924:        Since the user set the outer subdomains, even if nontrivial overlap was requested via PCGASMSetOverlap(),
925:        it will be ignored.  To avoid confusion later on (e.g., when viewing the PC), the overlap size is set to -1.
926:     */
927:     osm->overlap = -1;
928:     /* inner subdomains must be provided  */
929:     if (!iis) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"inner indices have to be provided \n");
930:   }/* end if */
931:   if (iis) {
932:     PetscMalloc1(n,&osm->iis);
933:     for (i=0; i<n; i++) {
934:       PetscObjectReference((PetscObject)iis[i]);
935:       osm->iis[i] = iis[i];
936:     }
937:     if (!ois) {
938:       osm->ois = NULL;
939:       /* if user does not provide outer indices, we will create the corresponding outer indices using  osm->overlap =1 in PCSetUp_GASM */
940:     }
941:   }
942: #if defined(PETSC_USE_DEBUG)
943:   {
944:     PetscInt        j,rstart,rend,*covered,lsize;
945:     const PetscInt  *indices;
946:     /* check if the inner indices cover and only cover the local portion of the preconditioning matrix */
947:     MatGetOwnershipRange(pc->pmat,&rstart,&rend);
948:     PetscCalloc1(rend-rstart,&covered);
949:     /* check if the current processor owns indices from others */
950:     for (i=0; i<n; i++) {
951:       ISGetIndices(osm->iis[i],&indices);
952:       ISGetLocalSize(osm->iis[i],&lsize);
953:       for (j=0; j<lsize; j++) {
954:         if (indices[j]<rstart || indices[j]>=rend) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"inner subdomains can not own an index %d from other processors", indices[j]);
955:         else if (covered[indices[j]-rstart]==1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"inner subdomains can not have an overlapping index %d ",indices[j]);
956:         else covered[indices[j]-rstart] = 1;
957:       }
958:     ISRestoreIndices(osm->iis[i],&indices);
959:     }
960:     /* check if we miss any indices */
961:     for (i=rstart; i<rend; i++) {
962:       if (!covered[i-rstart]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"local entity %d was not covered by inner subdomains",i);
963:     }
964:     PetscFree(covered);
965:   }
966: #endif
967:   if (iis)  osm->user_subdomains = PETSC_TRUE;
968:   return(0);
969: }


972: static PetscErrorCode  PCGASMSetOverlap_GASM(PC pc,PetscInt ovl)
973: {
974:   PC_GASM *osm = (PC_GASM*)pc->data;

977:   if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
978:   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp().");
979:   if (!pc->setupcalled) osm->overlap = ovl;
980:   return(0);
981: }

983: static PetscErrorCode  PCGASMSetType_GASM(PC pc,PCGASMType type)
984: {
985:   PC_GASM *osm = (PC_GASM*)pc->data;

988:   osm->type     = type;
989:   osm->type_set = PETSC_TRUE;
990:   return(0);
991: }

993: static PetscErrorCode  PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort)
994: {
995:   PC_GASM *osm = (PC_GASM*)pc->data;

998:   osm->sort_indices = doSort;
999:   return(0);
1000: }

1002: /*
1003:    FIXME: This routine might need to be modified now that multiple ranks per subdomain are allowed.
1004:         In particular, it would upset the global subdomain number calculation.
1005: */
1006: static PetscErrorCode  PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp)
1007: {
1008:   PC_GASM        *osm = (PC_GASM*)pc->data;

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

1014:   if (n) *n = osm->n;
1015:   if (first) {
1016:     MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));
1017:     *first -= osm->n;
1018:   }
1019:   if (ksp) {
1020:     /* Assume that local solves are now different; not necessarily
1021:        true, though!  This flag is used only for PCView_GASM() */
1022:     *ksp                        = osm->ksp;
1023:     osm->same_subdomain_solvers = PETSC_FALSE;
1024:   }
1025:   return(0);
1026: } /* PCGASMGetSubKSP_GASM() */

1028: /*@C
1029:     PCGASMSetSubdomains - Sets the subdomains for this processor
1030:     for the additive Schwarz preconditioner.

1032:     Collective on pc

1034:     Input Parameters:
1035: +   pc  - the preconditioner object
1036: .   n   - the number of subdomains for this processor
1037: .   iis - the index sets that define the inner subdomains (or NULL for PETSc to determine subdomains)
1038: -   ois - the index sets that define the outer subdomains (or NULL to use the same as iis, or to construct by expanding iis by the requested overlap)

1040:     Notes:
1041:     The IS indices use the parallel, global numbering of the vector entries.
1042:     Inner subdomains are those where the correction is applied.
1043:     Outer subdomains are those where the residual necessary to obtain the
1044:     corrections is obtained (see PCGASMType for the use of inner/outer subdomains).
1045:     Both inner and outer subdomains can extend over several processors.
1046:     This processor's portion of a subdomain is known as a local subdomain.

1048:     Inner subdomains can not overlap with each other, do not have any entities from remote processors,
1049:     and  have to cover the entire local subdomain owned by the current processor. The index sets on each
1050:     process should be ordered such that the ith local subdomain is connected to the ith remote subdomain
1051:     on another MPI process.

1053:     By default the GASM preconditioner uses 1 (local) subdomain per processor.


1056:     Level: advanced

1058: .seealso: PCGASMSetNumSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
1059:           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
1060: @*/
1061: PetscErrorCode  PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[])
1062: {
1063:   PC_GASM *osm = (PC_GASM*)pc->data;

1068:   PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));
1069:   osm->dm_subdomains = PETSC_FALSE;
1070:   return(0);
1071: }


1074: /*@
1075:     PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the
1076:     additive Schwarz preconditioner.  Either all or no processors in the
1077:     pc communicator must call this routine.

1079:     Logically Collective on pc

1081:     Input Parameters:
1082: +   pc  - the preconditioner context
1083: -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 0)

1085:     Options Database Key:
1086: .   -pc_gasm_overlap <overlap> - Sets overlap

1088:     Notes:
1089:     By default the GASM preconditioner uses 1 subdomain per processor.  To use
1090:     multiple subdomain per perocessor or "straddling" subdomains that intersect
1091:     multiple processors use PCGASMSetSubdomains() (or option -pc_gasm_total_subdomains <n>).

1093:     The overlap defaults to 0, so if one desires that no additional
1094:     overlap be computed beyond what may have been set with a call to
1095:     PCGASMSetSubdomains(), then ovl must be set to be 0.  In particular, if one does
1096:     not explicitly set the subdomains in Section 1.5 Writing Application Codes with PETSc code, then all overlap would be computed
1097:     internally by PETSc, and using an overlap of 0 would result in an GASM
1098:     variant that is equivalent to the block Jacobi preconditioner.

1100:     Note that one can define initial index sets with any overlap via
1101:     PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows
1102:     PETSc to extend that overlap further, if desired.

1104:     Level: intermediate

1106: .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1107:           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
1108: @*/
1109: PetscErrorCode  PCGASMSetOverlap(PC pc,PetscInt ovl)
1110: {
1112:   PC_GASM *osm = (PC_GASM*)pc->data;

1117:   PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));
1118:   osm->dm_subdomains = PETSC_FALSE;
1119:   return(0);
1120: }

1122: /*@
1123:     PCGASMSetType - Sets the type of restriction and interpolation used
1124:     for local problems in the additive Schwarz method.

1126:     Logically Collective on PC

1128:     Input Parameters:
1129: +   pc  - the preconditioner context
1130: -   type - variant of GASM, one of
1131: .vb
1132:       PC_GASM_BASIC       - full interpolation and restriction
1133:       PC_GASM_RESTRICT    - full restriction, local processor interpolation
1134:       PC_GASM_INTERPOLATE - full interpolation, local processor restriction
1135:       PC_GASM_NONE        - local processor restriction and interpolation
1136: .ve

1138:     Options Database Key:
1139: .   -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type

1141:     Level: intermediate

1143: .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1144:           PCGASMCreateSubdomains2D()
1145: @*/
1146: PetscErrorCode  PCGASMSetType(PC pc,PCGASMType type)
1147: {

1153:   PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));
1154:   return(0);
1155: }

1157: /*@
1158:     PCGASMSetSortIndices - Determines whether subdomain indices are sorted.

1160:     Logically Collective on PC

1162:     Input Parameters:
1163: +   pc  - the preconditioner context
1164: -   doSort - sort the subdomain indices

1166:     Level: intermediate

1168: .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1169:           PCGASMCreateSubdomains2D()
1170: @*/
1171: PetscErrorCode  PCGASMSetSortIndices(PC pc,PetscBool doSort)
1172: {

1178:   PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));
1179:   return(0);
1180: }

1182: /*@C
1183:    PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on
1184:    this processor.

1186:    Collective on PC iff first_local is requested

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

1191:    Output Parameters:
1192: +  n_local - the number of blocks on this processor or NULL
1193: .  first_local - the global number of the first block on this processor or NULL,
1194:                  all processors must request or all must pass NULL
1195: -  ksp - the array of KSP contexts

1197:    Note:
1198:    After PCGASMGetSubKSP() the array of KSPes is not to be freed

1200:    Currently for some matrix implementations only 1 block per processor
1201:    is supported.

1203:    You must call KSPSetUp() before calling PCGASMGetSubKSP().

1205:    Level: advanced

1207: .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(),
1208:           PCGASMCreateSubdomains2D(),
1209: @*/
1210: PetscErrorCode  PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
1211: {

1216:   PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));
1217:   return(0);
1218: }

1220: /* -------------------------------------------------------------------------------------*/
1221: /*MC
1222:    PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
1223:            its own KSP object.

1225:    Options Database Keys:
1226: +  -pc_gasm_total_subdomains <n>  - Sets total number of local subdomains to be distributed among processors
1227: .  -pc_gasm_view_subdomains       - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view
1228: .  -pc_gasm_print_subdomains      - activates the printing of subdomain indices in PCSetUp()
1229: .  -pc_gasm_overlap <ovl>         - Sets overlap by which to (automatically) extend local subdomains
1230: -  -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type

1232:      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
1233:       will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use
1234:       -pc_gasm_type basic to use the standard GASM.

1236:    Notes:
1237:     Blocks can be shared by multiple processes.

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

1242:      To set the options on the solvers separate for each block call PCGASMGetSubKSP()
1243:          and set the options directly on the resulting KSP object (you can access its PC
1244:          with KSPGetPC())


1247:    Level: beginner

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

1255: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1256:            PCBJACOBI,  PCGASMGetSubKSP(), PCGASMSetSubdomains(),
1257:            PCSetModifySubMatrices(), PCGASMSetOverlap(), PCGASMSetType()

1259: M*/

1261: PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc)
1262: {
1264:   PC_GASM        *osm;

1267:   PetscNewLog(pc,&osm);

1269:   osm->N                        = PETSC_DETERMINE;
1270:   osm->n                        = PETSC_DECIDE;
1271:   osm->nmax                     = PETSC_DETERMINE;
1272:   osm->overlap                  = 0;
1273:   osm->ksp                      = 0;
1274:   osm->gorestriction            = 0;
1275:   osm->girestriction            = 0;
1276:   osm->pctoouter                = 0;
1277:   osm->gx                       = 0;
1278:   osm->gy                       = 0;
1279:   osm->x                        = 0;
1280:   osm->y                        = 0;
1281:   osm->pcx                      = 0;
1282:   osm->pcy                      = 0;
1283:   osm->permutationIS            = 0;
1284:   osm->permutationP             = 0;
1285:   osm->pcmat                    = 0;
1286:   osm->ois                      = 0;
1287:   osm->iis                      = 0;
1288:   osm->pmat                     = 0;
1289:   osm->type                     = PC_GASM_RESTRICT;
1290:   osm->same_subdomain_solvers   = PETSC_TRUE;
1291:   osm->sort_indices             = PETSC_TRUE;
1292:   osm->dm_subdomains            = PETSC_FALSE;
1293:   osm->hierarchicalpartitioning = PETSC_FALSE;

1295:   pc->data                 = (void*)osm;
1296:   pc->ops->apply           = PCApply_GASM;
1297:   pc->ops->applytranspose  = PCApplyTranspose_GASM;
1298:   pc->ops->setup           = PCSetUp_GASM;
1299:   pc->ops->reset           = PCReset_GASM;
1300:   pc->ops->destroy         = PCDestroy_GASM;
1301:   pc->ops->setfromoptions  = PCSetFromOptions_GASM;
1302:   pc->ops->setuponblocks   = PCSetUpOnBlocks_GASM;
1303:   pc->ops->view            = PCView_GASM;
1304:   pc->ops->applyrichardson = 0;

1306:   PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSubdomains_C",PCGASMSetSubdomains_GASM);
1307:   PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetOverlap_C",PCGASMSetOverlap_GASM);
1308:   PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetType_C",PCGASMSetType_GASM);
1309:   PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSortIndices_C",PCGASMSetSortIndices_GASM);
1310:   PetscObjectComposeFunction((PetscObject)pc,"PCGASMGetSubKSP_C",PCGASMGetSubKSP_GASM);
1311:   return(0);
1312: }


1315: PetscErrorCode  PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[])
1316: {
1317:   MatPartitioning mpart;
1318:   const char      *prefix;
1319:   PetscInt        i,j,rstart,rend,bs;
1320:   PetscBool       hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
1321:   Mat             Ad     = NULL, adj;
1322:   IS              ispart,isnumb,*is;
1323:   PetscErrorCode  ierr;

1326:   if (nloc < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local subdomains must > 0, got nloc = %D",nloc);

1328:   /* Get prefix, row distribution, and block size */
1329:   MatGetOptionsPrefix(A,&prefix);
1330:   MatGetOwnershipRange(A,&rstart,&rend);
1331:   MatGetBlockSize(A,&bs);
1332:   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);

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

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

1407:     PetscInt mbs   = (rend-rstart)/bs;
1408:     PetscInt start = rstart;
1409:     for (i=0; i<nloc; i++) {
1410:       PetscInt count = (mbs/nloc + ((mbs % nloc) > i)) * bs;
1411:       ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);
1412:       start += count;
1413:     }

1415:   } else {

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

1419:     const PetscInt *numbering;
1420:     PetscInt       *count,nidx,*indices,*newidx,start=0;
1421:     /* Get node count in each partition */
1422:     PetscMalloc1(nloc,&count);
1423:     ISPartitioningCount(ispart,nloc,count);
1424:     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1425:       for (i=0; i<nloc; i++) count[i] *= bs;
1426:     }
1427:     /* Build indices from node numbering */
1428:     ISGetLocalSize(isnumb,&nidx);
1429:     PetscMalloc1(nidx,&indices);
1430:     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1431:     ISGetIndices(isnumb,&numbering);
1432:     PetscSortIntWithPermutation(nidx,numbering,indices);
1433:     ISRestoreIndices(isnumb,&numbering);
1434:     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1435:       PetscMalloc1(nidx*bs,&newidx);
1436:       for (i=0; i<nidx; i++) {
1437:         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
1438:       }
1439:       PetscFree(indices);
1440:       nidx   *= bs;
1441:       indices = newidx;
1442:     }
1443:     /* Shift to get global indices */
1444:     for (i=0; i<nidx; i++) indices[i] += rstart;

1446:     /* Build the index sets for each block */
1447:     for (i=0; i<nloc; i++) {
1448:       ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);
1449:       ISSort(is[i]);
1450:       start += count[i];
1451:     }

1453:     PetscFree(count);
1454:     PetscFree(indices);
1455:     ISDestroy(&isnumb);
1456:     ISDestroy(&ispart);
1457:   }
1458:   *iis = is;
1459:   return(0);
1460: }

1462: PETSC_INTERN PetscErrorCode  PCGASMCreateStraddlingSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[])
1463: {
1464:   PetscErrorCode  ierr;

1467:   MatSubdomainsCreateCoalesce(A,N,n,iis);
1468:   return(0);
1469: }



1473: /*@C
1474:    PCGASMCreateSubdomains - Creates n index sets defining n nonoverlapping subdomains for the additive
1475:    Schwarz preconditioner for a any problem based on its matrix.

1477:    Collective

1479:    Input Parameters:
1480: +  A       - The global matrix operator
1481: -  N       - the number of global subdomains requested

1483:    Output Parameters:
1484: +  n   - the number of subdomains created on this processor
1485: -  iis - the array of index sets defining the local inner subdomains (on which the correction is applied)

1487:    Level: advanced

1489:    Note: When N >= A's communicator size, each subdomain is local -- contained within a single processor.
1490:          When N < size, the subdomains are 'straddling' (processor boundaries) and are no longer local.
1491:          The resulting subdomains can be use in PCGASMSetSubdomains(pc,n,iss,NULL).  The overlapping
1492:          outer subdomains will be automatically generated from these according to the requested amount of
1493:          overlap; this is currently supported only with local subdomains.


1496: .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains()
1497: @*/
1498: PetscErrorCode  PCGASMCreateSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[])
1499: {
1500:   PetscMPIInt     size;
1501:   PetscErrorCode  ierr;


1507:   if (N < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of subdomains must be > 0, N = %D",N);
1508:   MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
1509:   if (N >= size) {
1510:     *n = N/size + (N%size);
1511:     PCGASMCreateLocalSubdomains(A,*n,iis);
1512:   } else {
1513:     PCGASMCreateStraddlingSubdomains(A,N,n,iis);
1514:   }
1515:   return(0);
1516: }

1518: /*@C
1519:    PCGASMDestroySubdomains - Destroys the index sets created with
1520:    PCGASMCreateSubdomains() or PCGASMCreateSubdomains2D. Should be
1521:    called after setting subdomains with PCGASMSetSubdomains().

1523:    Collective

1525:    Input Parameters:
1526: +  n   - the number of index sets
1527: .  iis - the array of inner subdomains,
1528: -  ois - the array of outer subdomains, can be NULL

1530:    Level: intermediate

1532:    Notes:
1533:     this is merely a convenience subroutine that walks each list,
1534:    destroys each IS on the list, and then frees the list. At the end the
1535:    list pointers are set to NULL.

1537: .seealso: PCGASMCreateSubdomains(), PCGASMSetSubdomains()
1538: @*/
1539: PetscErrorCode  PCGASMDestroySubdomains(PetscInt n,IS **iis,IS **ois)
1540: {
1541:   PetscInt       i;

1545:   if (n <= 0) return(0);
1546:   if (ois) {
1548:     if (*ois) {
1550:       for (i=0; i<n; i++) {
1551:         ISDestroy(&(*ois)[i]);
1552:       }
1553:       PetscFree((*ois));
1554:     }
1555:   }
1556:   if (iis) {
1558:     if (*iis) {
1560:       for (i=0; i<n; i++) {
1561:         ISDestroy(&(*iis)[i]);
1562:       }
1563:       PetscFree((*iis));
1564:     }
1565:   }
1566:   return(0);
1567: }


1570: #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \
1571:   {                                                                                                       \
1572:     PetscInt first_row = first/M, last_row = last/M+1;                                                     \
1573:     /*                                                                                                    \
1574:      Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners       \
1575:      of the bounding box of the intersection of the subdomain with the local ownership range (local       \
1576:      subdomain).                                                                                          \
1577:      Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows       \
1578:      of the intersection.                                                                                 \
1579:     */                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                             \
1580:     /* ylow_loc is the grid row containing the first element of the local sumbdomain */                   \
1581:     *ylow_loc = PetscMax(first_row,ylow);                                                                    \
1582:     /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \
1583:     *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft;                                                                            \
1584:     /* yhigh_loc is the grid row above the last local subdomain element */                                                                    \
1585:     *yhigh_loc = PetscMin(last_row,yhigh);                                                                                                     \
1586:     /* xright is the offset of the end of the  local subdomain within its grid row (might actually be outside the local subdomain) */         \
1587:     *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright;                                                                          \
1588:     /* Now compute the size of the local subdomain n. */ \
1589:     *n = 0;                                               \
1590:     if (*ylow_loc < *yhigh_loc) {                           \
1591:       PetscInt width = xright-xleft;                     \
1592:       *n += width*(*yhigh_loc-*ylow_loc-1);                 \
1593:       *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \
1594:       *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \
1595:     } \
1596:   }



1600: /*@
1601:    PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1602:    preconditioner for a two-dimensional problem on a regular grid.

1604:    Collective

1606:    Input Parameters:
1607: +  M, N               - the global number of grid points in the x and y directions
1608: .  Mdomains, Ndomains - the global number of subdomains in the x and y directions
1609: .  dof                - degrees of freedom per node
1610: -  overlap            - overlap in mesh lines

1612:    Output Parameters:
1613: +  Nsub - the number of local subdomains created
1614: .  iis  - array of index sets defining inner (nonoverlapping) subdomains
1615: -  ois  - array of index sets defining outer (overlapping, if overlap > 0) subdomains


1618:    Level: advanced

1620: .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), PCGASMSetOverlap()
1621: @*/
1622: PetscErrorCode  PCGASMCreateSubdomains2D(PC pc,PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap,PetscInt *nsub,IS **iis,IS **ois)
1623: {
1625:   PetscMPIInt    size, rank;
1626:   PetscInt       i, j;
1627:   PetscInt       maxheight, maxwidth;
1628:   PetscInt       xstart, xleft, xright, xleft_loc, xright_loc;
1629:   PetscInt       ystart, ylow,  yhigh,  ylow_loc,  yhigh_loc;
1630:   PetscInt       x[2][2], y[2][2], n[2];
1631:   PetscInt       first, last;
1632:   PetscInt       nidx, *idx;
1633:   PetscInt       ii,jj,s,q,d;
1634:   PetscInt       k,kk;
1635:   PetscMPIInt    color;
1636:   MPI_Comm       comm, subcomm;
1637:   IS             **xis = 0, **is = ois, **is_local = iis;

1640:   PetscObjectGetComm((PetscObject)pc, &comm);
1641:   MPI_Comm_size(comm, &size);
1642:   MPI_Comm_rank(comm, &rank);
1643:   MatGetOwnershipRange(pc->pmat, &first, &last);
1644:   if (first%dof || last%dof) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Matrix row partitioning unsuitable for domain decomposition: local row range (%D,%D) "
1645:                                       "does not respect the number of degrees of freedom per grid point %D", first, last, dof);

1647:   /* Determine the number of domains with nonzero intersections with the local ownership range. */
1648:   s      = 0;
1649:   ystart = 0;
1650:   for (j=0; j<Ndomains; ++j) {
1651:     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1652:     if (maxheight < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the vertical directon for mesh height %D", Ndomains, N);
1653:     /* Vertical domain limits with an overlap. */
1654:     ylow   = PetscMax(ystart - overlap,0);
1655:     yhigh  = PetscMin(ystart + maxheight + overlap,N);
1656:     xstart = 0;
1657:     for (i=0; i<Mdomains; ++i) {
1658:       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1659:       if (maxwidth < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the horizontal direction for mesh width %D", Mdomains, M);
1660:       /* Horizontal domain limits with an overlap. */
1661:       xleft  = PetscMax(xstart - overlap,0);
1662:       xright = PetscMin(xstart + maxwidth + overlap,M);
1663:       /*
1664:          Determine whether this subdomain intersects this processor's ownership range of pc->pmat.
1665:       */
1666:       PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1667:       if (nidx) ++s;
1668:       xstart += maxwidth;
1669:     } /* for (i = 0; i < Mdomains; ++i) */
1670:     ystart += maxheight;
1671:   } /* for (j = 0; j < Ndomains; ++j) */

1673:   /* Now we can allocate the necessary number of ISs. */
1674:   *nsub  = s;
1675:   PetscMalloc1(*nsub,is);
1676:   PetscMalloc1(*nsub,is_local);
1677:   s      = 0;
1678:   ystart = 0;
1679:   for (j=0; j<Ndomains; ++j) {
1680:     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1681:     if (maxheight < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the vertical directon for mesh height %D", Ndomains, N);
1682:     /* Vertical domain limits with an overlap. */
1683:     y[0][0] = PetscMax(ystart - overlap,0);
1684:     y[0][1] = PetscMin(ystart + maxheight + overlap,N);
1685:     /* Vertical domain limits without an overlap. */
1686:     y[1][0] = ystart;
1687:     y[1][1] = PetscMin(ystart + maxheight,N);
1688:     xstart  = 0;
1689:     for (i=0; i<Mdomains; ++i) {
1690:       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1691:       if (maxwidth < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the horizontal direction for mesh width %D", Mdomains, M);
1692:       /* Horizontal domain limits with an overlap. */
1693:       x[0][0] = PetscMax(xstart - overlap,0);
1694:       x[0][1] = PetscMin(xstart + maxwidth + overlap,M);
1695:       /* Horizontal domain limits without an overlap. */
1696:       x[1][0] = xstart;
1697:       x[1][1] = PetscMin(xstart+maxwidth,M);
1698:       /*
1699:          Determine whether this domain intersects this processor's ownership range of pc->pmat.
1700:          Do this twice: first for the domains with overlaps, and once without.
1701:          During the first pass create the subcommunicators, and use them on the second pass as well.
1702:       */
1703:       for (q = 0; q < 2; ++q) {
1704:         PetscBool split = PETSC_FALSE;
1705:         /*
1706:           domain limits, (xleft, xright) and (ylow, yheigh) are adjusted
1707:           according to whether the domain with an overlap or without is considered.
1708:         */
1709:         xleft = x[q][0]; xright = x[q][1];
1710:         ylow  = y[q][0]; yhigh  = y[q][1];
1711:         PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1712:         nidx *= dof;
1713:         n[q]  = nidx;
1714:         /*
1715:          Based on the counted number of indices in the local domain *with an overlap*,
1716:          construct a subcommunicator of all the processors supporting this domain.
1717:          Observe that a domain with an overlap might have nontrivial local support,
1718:          while the domain without an overlap might not.  Hence, the decision to participate
1719:          in the subcommunicator must be based on the domain with an overlap.
1720:          */
1721:         if (q == 0) {
1722:           if (nidx) color = 1;
1723:           else color = MPI_UNDEFINED;
1724:           MPI_Comm_split(comm, color, rank, &subcomm);
1725:           split = PETSC_TRUE;
1726:         }
1727:         /*
1728:          Proceed only if the number of local indices *with an overlap* is nonzero.
1729:          */
1730:         if (n[0]) {
1731:           if (q == 0) xis = is;
1732:           if (q == 1) {
1733:             /*
1734:              The IS for the no-overlap subdomain shares a communicator with the overlapping domain.
1735:              Moreover, if the overlap is zero, the two ISs are identical.
1736:              */
1737:             if (overlap == 0) {
1738:               (*is_local)[s] = (*is)[s];
1739:               PetscObjectReference((PetscObject)(*is)[s]);
1740:               continue;
1741:             } else {
1742:               xis     = is_local;
1743:               subcomm = ((PetscObject)(*is)[s])->comm;
1744:             }
1745:           } /* if (q == 1) */
1746:           idx  = NULL;
1747:           PetscMalloc1(nidx,&idx);
1748:           if (nidx) {
1749:             k = 0;
1750:             for (jj=ylow_loc; jj<yhigh_loc; ++jj) {
1751:               PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft;
1752:               PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright;
1753:               kk = dof*(M*jj + x0);
1754:               for (ii=x0; ii<x1; ++ii) {
1755:                 for (d = 0; d < dof; ++d) {
1756:                   idx[k++] = kk++;
1757:                 }
1758:               }
1759:             }
1760:           }
1761:           ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);
1762:           if (split) {
1763:             MPI_Comm_free(&subcomm);
1764:           }
1765:         }/* if (n[0]) */
1766:       }/* for (q = 0; q < 2; ++q) */
1767:       if (n[0]) ++s;
1768:       xstart += maxwidth;
1769:     } /* for (i = 0; i < Mdomains; ++i) */
1770:     ystart += maxheight;
1771:   } /* for (j = 0; j < Ndomains; ++j) */
1772:   return(0);
1773: }

1775: /*@C
1776:     PCGASMGetSubdomains - Gets the subdomains supported on this processor
1777:     for the additive Schwarz preconditioner.

1779:     Not Collective

1781:     Input Parameter:
1782: .   pc - the preconditioner context

1784:     Output Parameters:
1785: +   n   - the number of subdomains for this processor (default value = 1)
1786: .   iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be NULL)
1787: -   ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be NULL)


1790:     Notes:
1791:     The user is responsible for destroying the ISs and freeing the returned arrays.
1792:     The IS numbering is in the parallel, global numbering of the vector.

1794:     Level: advanced

1796: .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), PCGASMCreateSubdomains2D(),
1797:           PCGASMSetSubdomains(), PCGASMGetSubmatrices()
1798: @*/
1799: PetscErrorCode  PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[])
1800: {
1801:   PC_GASM        *osm;
1803:   PetscBool      match;
1804:   PetscInt       i;

1808:   PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);
1809:   if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1810:   osm = (PC_GASM*)pc->data;
1811:   if (n) *n = osm->n;
1812:   if (iis) {
1813:     PetscMalloc1(osm->n, iis);
1814:   }
1815:   if (ois) {
1816:     PetscMalloc1(osm->n, ois);
1817:   }
1818:   if (iis || ois) {
1819:     for (i = 0; i < osm->n; ++i) {
1820:       if (iis) (*iis)[i] = osm->iis[i];
1821:       if (ois) (*ois)[i] = osm->ois[i];
1822:     }
1823:   }
1824:   return(0);
1825: }

1827: /*@C
1828:     PCGASMGetSubmatrices - Gets the local submatrices (for this processor
1829:     only) for the additive Schwarz preconditioner.

1831:     Not Collective

1833:     Input Parameter:
1834: .   pc - the preconditioner context

1836:     Output Parameters:
1837: +   n   - the number of matrices for this processor (default value = 1)
1838: -   mat - the matrices

1840:     Notes:
1841:     matrices returned by this routine have the same communicators as the index sets (IS)
1842:            used to define subdomains in PCGASMSetSubdomains()
1843:     Level: advanced

1845: .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(),
1846:           PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains()
1847: @*/
1848: PetscErrorCode  PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1849: {
1850:   PC_GASM        *osm;
1852:   PetscBool      match;

1858:   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUp() or PCSetUp().");
1859:   PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);
1860:   if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1861:   osm = (PC_GASM*)pc->data;
1862:   if (n) *n = osm->n;
1863:   if (mat) *mat = osm->pmat;
1864:   return(0);
1865: }

1867: /*@
1868:     PCGASMSetUseDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1869:     Logically Collective

1871:     Input Parameter:
1872: +   pc  - the preconditioner
1873: -   flg - boolean indicating whether to use subdomains defined by the DM

1875:     Options Database Key:
1876: .   -pc_gasm_dm_subdomains -pc_gasm_overlap -pc_gasm_total_subdomains

1878:     Level: intermediate

1880:     Notes:
1881:     PCGASMSetSubdomains(), PCGASMSetTotalSubdomains() or PCGASMSetOverlap() take precedence over PCGASMSetUseDMSubdomains(),
1882:     so setting PCGASMSetSubdomains() with nontrivial subdomain ISs or any of PCGASMSetTotalSubdomains() and PCGASMSetOverlap()
1883:     automatically turns the latter off.

1885: .seealso: PCGASMGetUseDMSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap()
1886:           PCGASMCreateSubdomains2D()
1887: @*/
1888: PetscErrorCode  PCGASMSetUseDMSubdomains(PC pc,PetscBool flg)
1889: {
1890:   PC_GASM        *osm = (PC_GASM*)pc->data;
1892:   PetscBool      match;

1897:   if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
1898:   PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);
1899:   if (match) {
1900:     if (!osm->user_subdomains && osm->N == PETSC_DETERMINE && osm->overlap < 0) {
1901:       osm->dm_subdomains = flg;
1902:     }
1903:   }
1904:   return(0);
1905: }

1907: /*@
1908:     PCGASMGetUseDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1909:     Not Collective

1911:     Input Parameter:
1912: .   pc  - the preconditioner

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

1917:     Level: intermediate

1919: .seealso: PCGASMSetUseDMSubdomains(), PCGASMSetOverlap()
1920:           PCGASMCreateSubdomains2D()
1921: @*/
1922: PetscErrorCode  PCGASMGetUseDMSubdomains(PC pc,PetscBool* flg)
1923: {
1924:   PC_GASM        *osm = (PC_GASM*)pc->data;
1926:   PetscBool      match;

1931:   PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);
1932:   if (match) {
1933:     if (flg) *flg = osm->dm_subdomains;
1934:   }
1935:   return(0);
1936: }