Actual source code: gasm.c

petsc-master 2020-07-09
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,sizeof(fname),&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 = NULL;
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 (PetscDefined(USE_DEBUG)) {
943:     PetscInt        j,rstart,rend,*covered,lsize;
944:     const PetscInt  *indices;
945:     /* check if the inner indices cover and only cover the local portion of the preconditioning matrix */
946:     MatGetOwnershipRange(pc->pmat,&rstart,&rend);
947:     PetscCalloc1(rend-rstart,&covered);
948:     /* check if the current processor owns indices from others */
949:     for (i=0; i<n; i++) {
950:       ISGetIndices(osm->iis[i],&indices);
951:       ISGetLocalSize(osm->iis[i],&lsize);
952:       for (j=0; j<lsize; j++) {
953:         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]);
954:         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]);
955:         else covered[indices[j]-rstart] = 1;
956:       }
957:     ISRestoreIndices(osm->iis[i],&indices);
958:     }
959:     /* check if we miss any indices */
960:     for (i=rstart; i<rend; i++) {
961:       if (!covered[i-rstart]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"local entity %d was not covered by inner subdomains",i);
962:     }
963:     PetscFree(covered);
964:   }
965:   if (iis)  osm->user_subdomains = PETSC_TRUE;
966:   return(0);
967: }


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

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

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

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

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

996:   osm->sort_indices = doSort;
997:   return(0);
998: }

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

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

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

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

1030:     Collective on pc

1032:     Input Parameters:
1033: +   pc  - the preconditioner object
1034: .   n   - the number of subdomains for this processor
1035: .   iis - the index sets that define the inner subdomains (or NULL for PETSc to determine subdomains)
1036: -   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)

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

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

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


1054:     Level: advanced

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

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


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

1077:     Logically Collective on pc

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

1083:     Options Database Key:
1084: .   -pc_gasm_overlap <overlap> - Sets overlap

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

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

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

1102:     Level: intermediate

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

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

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

1124:     Logically Collective on PC

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

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

1139:     Level: intermediate

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

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

1155: /*@
1156:     PCGASMSetSortIndices - Determines whether subdomain indices are sorted.

1158:     Logically Collective on PC

1160:     Input Parameters:
1161: +   pc  - the preconditioner context
1162: -   doSort - sort the subdomain indices

1164:     Level: intermediate

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

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

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

1184:    Collective on PC iff first_local is requested

1186:    Input Parameter:
1187: .  pc - the preconditioner context

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

1195:    Note:
1196:    After PCGASMGetSubKSP() the array of KSPes is not to be freed

1198:    Currently for some matrix implementations only 1 block per processor
1199:    is supported.

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

1203:    Level: advanced

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

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

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

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

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

1234:    Notes:
1235:     Blocks can be shared by multiple processes.

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

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


1245:    Level: beginner

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

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

1257: M*/

1259: PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc)
1260: {
1262:   PC_GASM        *osm;

1265:   PetscNewLog(pc,&osm);

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

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

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


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

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

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

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

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

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

1413:   } else {

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

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

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

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

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

1465:   MatSubdomainsCreateCoalesce(A,N,n,iis);
1466:   return(0);
1467: }



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

1475:    Collective

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

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

1485:    Level: advanced

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


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


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

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

1521:    Collective

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

1528:    Level: intermediate

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

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

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


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



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

1602:    Collective

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

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


1616:    Level: advanced

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

1638:   PetscObjectGetComm((PetscObject)pc, &comm);
1639:   MPI_Comm_size(comm, &size);
1640:   MPI_Comm_rank(comm, &rank);
1641:   MatGetOwnershipRange(pc->pmat, &first, &last);
1642:   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) "
1643:                                       "does not respect the number of degrees of freedom per grid point %D", first, last, dof);

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

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

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

1777:     Not Collective

1779:     Input Parameter:
1780: .   pc - the preconditioner context

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


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

1792:     Level: advanced

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

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

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

1829:     Not Collective

1831:     Input Parameter:
1832: .   pc - the preconditioner context

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

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

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

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

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

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

1873:     Options Database Key:
1874: .   -pc_gasm_dm_subdomains -pc_gasm_overlap -pc_gasm_total_subdomains

1876:     Level: intermediate

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

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

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

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

1909:     Input Parameter:
1910: .   pc  - the preconditioner

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

1915:     Level: intermediate

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

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