Actual source code: gasm.c

petsc-3.14.6 2021-03-30
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: }

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

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

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

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

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

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

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

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

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


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

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

415:     if (osm->sort_indices) {
416:       for (i=0; i<osm->n; i++) {
417:         ISSort(osm->ois[i]);
418:         ISSort(osm->iis[i]);
419:       }
420:     }
421:     PCGetOptionsPrefix(pc,&prefix);
422:     PCGASMPrintSubdomains(pc);

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

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

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

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

579:   /*
580:      Extract the submatrices.
581:   */
582:   if (size > 1) {
583:     MatCreateSubMatricesMPI(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);
584:   } else {
585:     MatCreateSubMatrices(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);
586:   }
587:   if (scall == MAT_INITIAL_MATRIX) {
588:     PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);
589:     for (i=0; i<osm->n; i++) {
590:       PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);
591:       PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);
592:     }
593:   }

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

599:   /*
600:      Loop over submatrices putting them into local ksps
601:   */
602:   for (i=0; i<osm->n; i++) {
603:     KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);
604:     if (!pc->setupcalled) {
605:       KSPSetFromOptions(osm->ksp[i]);
606:     }
607:   }
608:   if (osm->pcmat){
609:     MatDestroy(&pc->pmat);
610:     pc->pmat   = osm->pcmat;
611:     osm->pcmat = NULL;
612:   }
613:   return(0);
614: }

616: static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc)
617: {
618:   PC_GASM        *osm = (PC_GASM*)pc->data;
620:   PetscInt       i;

623:   for (i=0; i<osm->n; i++) {
624:     KSPSetUp(osm->ksp[i]);
625:   }
626:   return(0);
627: }

629: static PetscErrorCode PCApply_GASM(PC pc,Vec xin,Vec yout)
630: {
631:   PC_GASM        *osm = (PC_GASM*)pc->data;
633:   PetscInt       i;
634:   Vec            x,y;
635:   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;

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

685: static PetscErrorCode PCMatApply_GASM(PC pc,Mat Xin,Mat Yout)
686: {
687:   PC_GASM        *osm = (PC_GASM*)pc->data;
688:   Mat            X,Y,O=NULL,Z,W;
689:   Vec            x,y;
690:   PetscInt       i,m,M,N;
691:   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;

695:   if (osm->n != 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not yet implemented");
696:   MatGetSize(Xin,NULL,&N);
697:   if (osm->pctoouter) {
698:     VecGetLocalSize(osm->pcx,&m);
699:     VecGetSize(osm->pcx,&M);
700:     MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]),m,PETSC_DECIDE,M,N,NULL,&O);
701:     for (i = 0; i < N; ++i) {
702:       MatDenseGetColumnVecRead(Xin,i,&x);
703:       MatDenseGetColumnVecWrite(O,i,&y);
704:       VecScatterBegin(osm->pctoouter,x,y,INSERT_VALUES,SCATTER_REVERSE);
705:       VecScatterEnd(osm->pctoouter,x,y,INSERT_VALUES,SCATTER_REVERSE);
706:       MatDenseRestoreColumnVecWrite(O,i,&y);
707:       MatDenseRestoreColumnVecRead(Xin,i,&x);
708:     }
709:     X = Y = O;
710:   } else {
711:     X = Xin;
712:     Y = Yout;
713:   }
714:   /*
715:      support for limiting the restriction or interpolation only to the inner
716:      subdomain values (leaving the other values 0).
717:   */
718:   VecGetLocalSize(osm->x[0],&m);
719:   VecGetSize(osm->x[0],&M);
720:   MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]),m,PETSC_DECIDE,M,N,NULL,&Z);
721:   for (i = 0; i < N; ++i) {
722:     MatDenseGetColumnVecRead(X,i,&x);
723:     MatDenseGetColumnVecWrite(Z,i,&y);
724:     if (!(osm->type & PC_GASM_RESTRICT)) {
725:       /* have to zero the work RHS since scatter may leave some slots empty */
726:       VecZeroEntries(y);
727:       VecScatterBegin(osm->girestriction,x,y,INSERT_VALUES,forward);
728:       VecScatterEnd(osm->girestriction,x,y,INSERT_VALUES,forward);
729:     } else {
730:       VecScatterBegin(osm->gorestriction,x,y,INSERT_VALUES,forward);
731:       VecScatterEnd(osm->gorestriction,x,y,INSERT_VALUES,forward);
732:     }
733:     MatDenseRestoreColumnVecWrite(Z,i,&y);
734:     MatDenseRestoreColumnVecRead(X,i,&x);
735:   }
736:   MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]),m,PETSC_DECIDE,M,N,NULL,&W);
737:   MatSetOption(Z,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);
738:   MatAssemblyBegin(Z,MAT_FINAL_ASSEMBLY);
739:   MatAssemblyEnd(Z,MAT_FINAL_ASSEMBLY);
740:   /* do the subdomain solve */
741:   KSPMatSolve(osm->ksp[0],Z,W);
742:   KSPCheckSolve(osm->ksp[0],pc,NULL);
743:   MatDestroy(&Z);
744:   /* do we need to zero y? */
745:   MatZeroEntries(Y);
746:   for (i = 0; i < N; ++i) {
747:     MatDenseGetColumnVecWrite(Y,i,&y);
748:     MatDenseGetColumnVecRead(W,i,&x);
749:     if (!(osm->type & PC_GASM_INTERPOLATE)) {
750:       VecScatterBegin(osm->girestriction,x,y,ADD_VALUES,reverse);
751:       VecScatterEnd(osm->girestriction,x,y,ADD_VALUES,reverse);
752:     } else {
753:       VecScatterBegin(osm->gorestriction,x,y,ADD_VALUES,reverse);
754:       VecScatterEnd(osm->gorestriction,x,y,ADD_VALUES,reverse);
755:     }
756:     MatDenseRestoreColumnVecRead(W,i,&x);
757:     if (osm->pctoouter) {
758:       MatDenseGetColumnVecWrite(Yout,i,&x);
759:       VecScatterBegin(osm->pctoouter,y,x,INSERT_VALUES,SCATTER_FORWARD);
760:       VecScatterEnd(osm->pctoouter,y,x,INSERT_VALUES,SCATTER_FORWARD);
761:       MatDenseRestoreColumnVecRead(Yout,i,&x);
762:     }
763:     MatDenseRestoreColumnVecWrite(Y,i,&y);
764:   }
765:   MatDestroy(&W);
766:   MatDestroy(&O);
767:   return(0);
768: }

770: static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec xin,Vec yout)
771: {
772:   PC_GASM        *osm = (PC_GASM*)pc->data;
774:   PetscInt       i;
775:   Vec            x,y;
776:   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;

779:   if (osm->pctoouter){
780:    VecScatterBegin(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);
781:    VecScatterEnd(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);
782:    x = osm->pcx;
783:    y = osm->pcy;
784:   }else{
785:         x = xin;
786:         y = yout;
787:   }
788:   /*
789:      Support for limiting the restriction or interpolation to only local
790:      subdomain values (leaving the other values 0).

792:      Note: these are reversed from the PCApply_GASM() because we are applying the
793:      transpose of the three terms
794:   */
795:   if (!(osm->type & PC_GASM_INTERPOLATE)) {
796:     /* have to zero the work RHS since scatter may leave some slots empty */
797:     VecZeroEntries(osm->gx);
798:     VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);
799:   } else {
800:     VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);
801:   }
802:   VecZeroEntries(osm->gy);
803:   if (!(osm->type & PC_GASM_INTERPOLATE)) {
804:     VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);
805:   } else {
806:     VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);
807:   }
808:   /* do the local solves */
809:   for (i=0; i<osm->n; ++i) { /* Note that the solves are local, so we can go to osm->n, rather than osm->nmax. */
810:     KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);
811:     KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);
812:   }
813:   VecZeroEntries(y);
814:   if (!(osm->type & PC_GASM_RESTRICT)) {
815:     VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);
816:     VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);
817:   } else {
818:     VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);
819:     VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);
820:   }
821:   if (osm->pctoouter){
822:    VecScatterBegin(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);
823:    VecScatterEnd(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);
824:   }
825:   return(0);
826: }

828: static PetscErrorCode PCReset_GASM(PC pc)
829: {
830:   PC_GASM        *osm = (PC_GASM*)pc->data;
832:   PetscInt       i;

835:   if (osm->ksp) {
836:     for (i=0; i<osm->n; i++) {
837:       KSPReset(osm->ksp[i]);
838:     }
839:   }
840:   if (osm->pmat) {
841:     if (osm->n > 0) {
842:       PetscMPIInt size;
843:       MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
844:       if (size > 1) {
845:         /* osm->pmat is created by MatCreateSubMatricesMPI(), cannot use MatDestroySubMatrices() */
846:         MatDestroyMatrices(osm->n,&osm->pmat);
847:       } else {
848:         MatDestroySubMatrices(osm->n,&osm->pmat);
849:       }
850:     }
851:   }
852:   if (osm->x) {
853:     for (i=0; i<osm->n; i++) {
854:       VecDestroy(&osm->x[i]);
855:       VecDestroy(&osm->y[i]);
856:     }
857:   }
858:   VecDestroy(&osm->gx);
859:   VecDestroy(&osm->gy);

861:   VecScatterDestroy(&osm->gorestriction);
862:   VecScatterDestroy(&osm->girestriction);
863:   if (!osm->user_subdomains) {
864:     PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);
865:     osm->N    = PETSC_DETERMINE;
866:     osm->nmax = PETSC_DETERMINE;
867:   }
868:   if (osm->pctoouter){
869:         VecScatterDestroy(&(osm->pctoouter));
870:   }
871:   if (osm->permutationIS){
872:         ISDestroy(&(osm->permutationIS));
873:   }
874:   if (osm->pcx){
875:         VecDestroy(&(osm->pcx));
876:   }
877:   if (osm->pcy){
878:         VecDestroy(&(osm->pcy));
879:   }
880:   if (osm->permutationP){
881:     MatDestroy(&(osm->permutationP));
882:   }
883:   if (osm->pcmat){
884:         MatDestroy(&osm->pcmat);
885:   }
886:   return(0);
887: }

889: static PetscErrorCode PCDestroy_GASM(PC pc)
890: {
891:   PC_GASM        *osm = (PC_GASM*)pc->data;
893:   PetscInt       i;

896:   PCReset_GASM(pc);
897:   /* PCReset will not destroy subdomains, if user_subdomains is true. */
898:   PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);
899:   if (osm->ksp) {
900:     for (i=0; i<osm->n; i++) {
901:       KSPDestroy(&osm->ksp[i]);
902:     }
903:     PetscFree(osm->ksp);
904:   }
905:   PetscFree(osm->x);
906:   PetscFree(osm->y);
907:   PetscFree(pc->data);
908:   return(0);
909: }

911: static PetscErrorCode PCSetFromOptions_GASM(PetscOptionItems *PetscOptionsObject,PC pc)
912: {
913:   PC_GASM        *osm = (PC_GASM*)pc->data;
915:   PetscInt       blocks,ovl;
916:   PetscBool      flg;
917:   PCGASMType     gasmtype;

920:   PetscOptionsHead(PetscOptionsObject,"Generalized additive Schwarz options");
921:   PetscOptionsBool("-pc_gasm_use_dm_subdomains","If subdomains aren't set, use DMCreateDomainDecomposition() to define subdomains.","PCGASMSetUseDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);
922:   PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->N,&blocks,&flg);
923:   if (flg) {
924:     PCGASMSetTotalSubdomains(pc,blocks);
925:   }
926:   PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);
927:   if (flg) {
928:     PCGASMSetOverlap(pc,ovl);
929:     osm->dm_subdomains = PETSC_FALSE;
930:   }
931:   flg  = PETSC_FALSE;
932:   PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);
933:   if (flg) {PCGASMSetType(pc,gasmtype);}
934:   PetscOptionsBool("-pc_gasm_use_hierachical_partitioning","use hierarchical partitioning",NULL,osm->hierarchicalpartitioning,&osm->hierarchicalpartitioning,&flg);
935:   PetscOptionsTail();
936:   return(0);
937: }

939: /*------------------------------------------------------------------------------------*/

941: /*@
942:     PCGASMSetTotalSubdomains - sets the total number of subdomains to use across the
943:                                communicator.
944:     Logically collective on pc

946:     Input Parameters:
947: +   pc  - the preconditioner
948: -   N   - total number of subdomains


951:     Level: beginner

953: .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap()
954:           PCGASMCreateSubdomains2D()
955: @*/
956: PetscErrorCode  PCGASMSetTotalSubdomains(PC pc,PetscInt N)
957: {
958:   PC_GASM        *osm = (PC_GASM*)pc->data;
959:   PetscMPIInt    size,rank;

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

966:   PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);
967:   osm->ois = osm->iis = NULL;

969:   MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
970:   MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
971:   osm->N    = N;
972:   osm->n    = PETSC_DETERMINE;
973:   osm->nmax = PETSC_DETERMINE;
974:   osm->dm_subdomains = PETSC_FALSE;
975:   return(0);
976: }

978: static PetscErrorCode  PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[])
979: {
980:   PC_GASM         *osm = (PC_GASM*)pc->data;
981:   PetscErrorCode  ierr;
982:   PetscInt        i;

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

988:   PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);
989:   osm->iis  = osm->ois = NULL;
990:   osm->n    = n;
991:   osm->N    = PETSC_DETERMINE;
992:   osm->nmax = PETSC_DETERMINE;
993:   if (ois) {
994:     PetscMalloc1(n,&osm->ois);
995:     for (i=0; i<n; i++) {
996:       PetscObjectReference((PetscObject)ois[i]);
997:       osm->ois[i] = ois[i];
998:     }
999:     /*
1000:        Since the user set the outer subdomains, even if nontrivial overlap was requested via PCGASMSetOverlap(),
1001:        it will be ignored.  To avoid confusion later on (e.g., when viewing the PC), the overlap size is set to -1.
1002:     */
1003:     osm->overlap = -1;
1004:     /* inner subdomains must be provided  */
1005:     if (!iis) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"inner indices have to be provided \n");
1006:   }/* end if */
1007:   if (iis) {
1008:     PetscMalloc1(n,&osm->iis);
1009:     for (i=0; i<n; i++) {
1010:       PetscObjectReference((PetscObject)iis[i]);
1011:       osm->iis[i] = iis[i];
1012:     }
1013:     if (!ois) {
1014:       osm->ois = NULL;
1015:       /* if user does not provide outer indices, we will create the corresponding outer indices using  osm->overlap =1 in PCSetUp_GASM */
1016:     }
1017:   }
1018:   if (PetscDefined(USE_DEBUG)) {
1019:     PetscInt        j,rstart,rend,*covered,lsize;
1020:     const PetscInt  *indices;
1021:     /* check if the inner indices cover and only cover the local portion of the preconditioning matrix */
1022:     MatGetOwnershipRange(pc->pmat,&rstart,&rend);
1023:     PetscCalloc1(rend-rstart,&covered);
1024:     /* check if the current processor owns indices from others */
1025:     for (i=0; i<n; i++) {
1026:       ISGetIndices(osm->iis[i],&indices);
1027:       ISGetLocalSize(osm->iis[i],&lsize);
1028:       for (j=0; j<lsize; j++) {
1029:         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]);
1030:         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]);
1031:         else covered[indices[j]-rstart] = 1;
1032:       }
1033:     ISRestoreIndices(osm->iis[i],&indices);
1034:     }
1035:     /* check if we miss any indices */
1036:     for (i=rstart; i<rend; i++) {
1037:       if (!covered[i-rstart]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"local entity %d was not covered by inner subdomains",i);
1038:     }
1039:     PetscFree(covered);
1040:   }
1041:   if (iis)  osm->user_subdomains = PETSC_TRUE;
1042:   return(0);
1043: }

1045: static PetscErrorCode  PCGASMSetOverlap_GASM(PC pc,PetscInt ovl)
1046: {
1047:   PC_GASM *osm = (PC_GASM*)pc->data;

1050:   if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
1051:   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp().");
1052:   if (!pc->setupcalled) osm->overlap = ovl;
1053:   return(0);
1054: }

1056: static PetscErrorCode  PCGASMSetType_GASM(PC pc,PCGASMType type)
1057: {
1058:   PC_GASM *osm = (PC_GASM*)pc->data;

1061:   osm->type     = type;
1062:   osm->type_set = PETSC_TRUE;
1063:   return(0);
1064: }

1066: static PetscErrorCode  PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort)
1067: {
1068:   PC_GASM *osm = (PC_GASM*)pc->data;

1071:   osm->sort_indices = doSort;
1072:   return(0);
1073: }

1075: /*
1076:    FIXME: This routine might need to be modified now that multiple ranks per subdomain are allowed.
1077:         In particular, it would upset the global subdomain number calculation.
1078: */
1079: static PetscErrorCode  PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp)
1080: {
1081:   PC_GASM        *osm = (PC_GASM*)pc->data;

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

1087:   if (n) *n = osm->n;
1088:   if (first) {
1089:     MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));
1090:     *first -= osm->n;
1091:   }
1092:   if (ksp) {
1093:     /* Assume that local solves are now different; not necessarily
1094:        true, though!  This flag is used only for PCView_GASM() */
1095:     *ksp                        = osm->ksp;
1096:     osm->same_subdomain_solvers = PETSC_FALSE;
1097:   }
1098:   return(0);
1099: } /* PCGASMGetSubKSP_GASM() */

1101: /*@C
1102:     PCGASMSetSubdomains - Sets the subdomains for this processor
1103:     for the additive Schwarz preconditioner.

1105:     Collective on pc

1107:     Input Parameters:
1108: +   pc  - the preconditioner object
1109: .   n   - the number of subdomains for this processor
1110: .   iis - the index sets that define the inner subdomains (or NULL for PETSc to determine subdomains)
1111: -   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)

1113:     Notes:
1114:     The IS indices use the parallel, global numbering of the vector entries.
1115:     Inner subdomains are those where the correction is applied.
1116:     Outer subdomains are those where the residual necessary to obtain the
1117:     corrections is obtained (see PCGASMType for the use of inner/outer subdomains).
1118:     Both inner and outer subdomains can extend over several processors.
1119:     This processor's portion of a subdomain is known as a local subdomain.

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

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


1129:     Level: advanced

1131: .seealso: PCGASMSetNumSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
1132:           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
1133: @*/
1134: PetscErrorCode  PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[])
1135: {
1136:   PC_GASM *osm = (PC_GASM*)pc->data;

1141:   PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));
1142:   osm->dm_subdomains = PETSC_FALSE;
1143:   return(0);
1144: }

1146: /*@
1147:     PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the
1148:     additive Schwarz preconditioner.  Either all or no processors in the
1149:     pc communicator must call this routine.

1151:     Logically Collective on pc

1153:     Input Parameters:
1154: +   pc  - the preconditioner context
1155: -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 0)

1157:     Options Database Key:
1158: .   -pc_gasm_overlap <overlap> - Sets overlap

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

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

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

1176:     Level: intermediate

1178: .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1179:           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
1180: @*/
1181: PetscErrorCode  PCGASMSetOverlap(PC pc,PetscInt ovl)
1182: {
1184:   PC_GASM *osm = (PC_GASM*)pc->data;

1189:   PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));
1190:   osm->dm_subdomains = PETSC_FALSE;
1191:   return(0);
1192: }

1194: /*@
1195:     PCGASMSetType - Sets the type of restriction and interpolation used
1196:     for local problems in the additive Schwarz method.

1198:     Logically Collective on PC

1200:     Input Parameters:
1201: +   pc  - the preconditioner context
1202: -   type - variant of GASM, one of
1203: .vb
1204:       PC_GASM_BASIC       - full interpolation and restriction
1205:       PC_GASM_RESTRICT    - full restriction, local processor interpolation
1206:       PC_GASM_INTERPOLATE - full interpolation, local processor restriction
1207:       PC_GASM_NONE        - local processor restriction and interpolation
1208: .ve

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

1213:     Level: intermediate

1215: .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1216:           PCGASMCreateSubdomains2D()
1217: @*/
1218: PetscErrorCode  PCGASMSetType(PC pc,PCGASMType type)
1219: {

1225:   PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));
1226:   return(0);
1227: }

1229: /*@
1230:     PCGASMSetSortIndices - Determines whether subdomain indices are sorted.

1232:     Logically Collective on PC

1234:     Input Parameters:
1235: +   pc  - the preconditioner context
1236: -   doSort - sort the subdomain indices

1238:     Level: intermediate

1240: .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1241:           PCGASMCreateSubdomains2D()
1242: @*/
1243: PetscErrorCode  PCGASMSetSortIndices(PC pc,PetscBool doSort)
1244: {

1250:   PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));
1251:   return(0);
1252: }

1254: /*@C
1255:    PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on
1256:    this processor.

1258:    Collective on PC iff first_local is requested

1260:    Input Parameter:
1261: .  pc - the preconditioner context

1263:    Output Parameters:
1264: +  n_local - the number of blocks on this processor or NULL
1265: .  first_local - the global number of the first block on this processor or NULL,
1266:                  all processors must request or all must pass NULL
1267: -  ksp - the array of KSP contexts

1269:    Note:
1270:    After PCGASMGetSubKSP() the array of KSPes is not to be freed

1272:    Currently for some matrix implementations only 1 block per processor
1273:    is supported.

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

1277:    Level: advanced

1279: .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(),
1280:           PCGASMCreateSubdomains2D(),
1281: @*/
1282: PetscErrorCode  PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
1283: {

1288:   PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));
1289:   return(0);
1290: }

1292: /* -------------------------------------------------------------------------------------*/
1293: /*MC
1294:    PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
1295:            its own KSP object.

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

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

1308:    Notes:
1309:     Blocks can be shared by multiple processes.

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

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


1319:    Level: beginner

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

1327: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1328:            PCBJACOBI,  PCGASMGetSubKSP(), PCGASMSetSubdomains(),
1329:            PCSetModifySubMatrices(), PCGASMSetOverlap(), PCGASMSetType()

1331: M*/

1333: PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc)
1334: {
1336:   PC_GASM        *osm;

1339:   PetscNewLog(pc,&osm);

1341:   osm->N                        = PETSC_DETERMINE;
1342:   osm->n                        = PETSC_DECIDE;
1343:   osm->nmax                     = PETSC_DETERMINE;
1344:   osm->overlap                  = 0;
1345:   osm->ksp                      = NULL;
1346:   osm->gorestriction            = NULL;
1347:   osm->girestriction            = NULL;
1348:   osm->pctoouter                = NULL;
1349:   osm->gx                       = NULL;
1350:   osm->gy                       = NULL;
1351:   osm->x                        = NULL;
1352:   osm->y                        = NULL;
1353:   osm->pcx                      = NULL;
1354:   osm->pcy                      = NULL;
1355:   osm->permutationIS            = NULL;
1356:   osm->permutationP             = NULL;
1357:   osm->pcmat                    = NULL;
1358:   osm->ois                      = NULL;
1359:   osm->iis                      = NULL;
1360:   osm->pmat                     = NULL;
1361:   osm->type                     = PC_GASM_RESTRICT;
1362:   osm->same_subdomain_solvers   = PETSC_TRUE;
1363:   osm->sort_indices             = PETSC_TRUE;
1364:   osm->dm_subdomains            = PETSC_FALSE;
1365:   osm->hierarchicalpartitioning = PETSC_FALSE;

1367:   pc->data                 = (void*)osm;
1368:   pc->ops->apply           = PCApply_GASM;
1369:   pc->ops->matapply        = PCMatApply_GASM;
1370:   pc->ops->applytranspose  = PCApplyTranspose_GASM;
1371:   pc->ops->setup           = PCSetUp_GASM;
1372:   pc->ops->reset           = PCReset_GASM;
1373:   pc->ops->destroy         = PCDestroy_GASM;
1374:   pc->ops->setfromoptions  = PCSetFromOptions_GASM;
1375:   pc->ops->setuponblocks   = PCSetUpOnBlocks_GASM;
1376:   pc->ops->view            = PCView_GASM;
1377:   pc->ops->applyrichardson = NULL;

1379:   PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSubdomains_C",PCGASMSetSubdomains_GASM);
1380:   PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetOverlap_C",PCGASMSetOverlap_GASM);
1381:   PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetType_C",PCGASMSetType_GASM);
1382:   PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSortIndices_C",PCGASMSetSortIndices_GASM);
1383:   PetscObjectComposeFunction((PetscObject)pc,"PCGASMGetSubKSP_C",PCGASMGetSubKSP_GASM);
1384:   return(0);
1385: }

1387: PetscErrorCode  PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[])
1388: {
1389:   MatPartitioning mpart;
1390:   const char      *prefix;
1391:   PetscInt        i,j,rstart,rend,bs;
1392:   PetscBool       hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
1393:   Mat             Ad     = NULL, adj;
1394:   IS              ispart,isnumb,*is;
1395:   PetscErrorCode  ierr;

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

1400:   /* Get prefix, row distribution, and block size */
1401:   MatGetOptionsPrefix(A,&prefix);
1402:   MatGetOwnershipRange(A,&rstart,&rend);
1403:   MatGetBlockSize(A,&bs);
1404:   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);

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

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

1479:     PetscInt mbs   = (rend-rstart)/bs;
1480:     PetscInt start = rstart;
1481:     for (i=0; i<nloc; i++) {
1482:       PetscInt count = (mbs/nloc + ((mbs % nloc) > i)) * bs;
1483:       ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);
1484:       start += count;
1485:     }

1487:   } else {

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

1491:     const PetscInt *numbering;
1492:     PetscInt       *count,nidx,*indices,*newidx,start=0;
1493:     /* Get node count in each partition */
1494:     PetscMalloc1(nloc,&count);
1495:     ISPartitioningCount(ispart,nloc,count);
1496:     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1497:       for (i=0; i<nloc; i++) count[i] *= bs;
1498:     }
1499:     /* Build indices from node numbering */
1500:     ISGetLocalSize(isnumb,&nidx);
1501:     PetscMalloc1(nidx,&indices);
1502:     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1503:     ISGetIndices(isnumb,&numbering);
1504:     PetscSortIntWithPermutation(nidx,numbering,indices);
1505:     ISRestoreIndices(isnumb,&numbering);
1506:     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1507:       PetscMalloc1(nidx*bs,&newidx);
1508:       for (i=0; i<nidx; i++) {
1509:         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
1510:       }
1511:       PetscFree(indices);
1512:       nidx   *= bs;
1513:       indices = newidx;
1514:     }
1515:     /* Shift to get global indices */
1516:     for (i=0; i<nidx; i++) indices[i] += rstart;

1518:     /* Build the index sets for each block */
1519:     for (i=0; i<nloc; i++) {
1520:       ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);
1521:       ISSort(is[i]);
1522:       start += count[i];
1523:     }

1525:     PetscFree(count);
1526:     PetscFree(indices);
1527:     ISDestroy(&isnumb);
1528:     ISDestroy(&ispart);
1529:   }
1530:   *iis = is;
1531:   return(0);
1532: }

1534: PETSC_INTERN PetscErrorCode  PCGASMCreateStraddlingSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[])
1535: {
1536:   PetscErrorCode  ierr;

1539:   MatSubdomainsCreateCoalesce(A,N,n,iis);
1540:   return(0);
1541: }

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

1547:    Collective

1549:    Input Parameters:
1550: +  A       - The global matrix operator
1551: -  N       - the number of global subdomains requested

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

1557:    Level: advanced

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


1566: .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains()
1567: @*/
1568: PetscErrorCode  PCGASMCreateSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[])
1569: {
1570:   PetscMPIInt     size;
1571:   PetscErrorCode  ierr;


1577:   if (N < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of subdomains must be > 0, N = %D",N);
1578:   MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
1579:   if (N >= size) {
1580:     *n = N/size + (N%size);
1581:     PCGASMCreateLocalSubdomains(A,*n,iis);
1582:   } else {
1583:     PCGASMCreateStraddlingSubdomains(A,N,n,iis);
1584:   }
1585:   return(0);
1586: }

1588: /*@C
1589:    PCGASMDestroySubdomains - Destroys the index sets created with
1590:    PCGASMCreateSubdomains() or PCGASMCreateSubdomains2D. Should be
1591:    called after setting subdomains with PCGASMSetSubdomains().

1593:    Collective

1595:    Input Parameters:
1596: +  n   - the number of index sets
1597: .  iis - the array of inner subdomains,
1598: -  ois - the array of outer subdomains, can be NULL

1600:    Level: intermediate

1602:    Notes:
1603:     this is merely a convenience subroutine that walks each list,
1604:    destroys each IS on the list, and then frees the list. At the end the
1605:    list pointers are set to NULL.

1607: .seealso: PCGASMCreateSubdomains(), PCGASMSetSubdomains()
1608: @*/
1609: PetscErrorCode  PCGASMDestroySubdomains(PetscInt n,IS **iis,IS **ois)
1610: {
1611:   PetscInt       i;

1615:   if (n <= 0) return(0);
1616:   if (ois) {
1618:     if (*ois) {
1620:       for (i=0; i<n; i++) {
1621:         ISDestroy(&(*ois)[i]);
1622:       }
1623:       PetscFree((*ois));
1624:     }
1625:   }
1626:   if (iis) {
1628:     if (*iis) {
1630:       for (i=0; i<n; i++) {
1631:         ISDestroy(&(*iis)[i]);
1632:       }
1633:       PetscFree((*iis));
1634:     }
1635:   }
1636:   return(0);
1637: }

1639: #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \
1640:   {                                                                                                       \
1641:     PetscInt first_row = first/M, last_row = last/M+1;                                                     \
1642:     /*                                                                                                    \
1643:      Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners       \
1644:      of the bounding box of the intersection of the subdomain with the local ownership range (local       \
1645:      subdomain).                                                                                          \
1646:      Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows       \
1647:      of the intersection.                                                                                 \
1648:     */                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                             \
1649:     /* ylow_loc is the grid row containing the first element of the local sumbdomain */                   \
1650:     *ylow_loc = PetscMax(first_row,ylow);                                                                    \
1651:     /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \
1652:     *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft;                                                                            \
1653:     /* yhigh_loc is the grid row above the last local subdomain element */                                                                    \
1654:     *yhigh_loc = PetscMin(last_row,yhigh);                                                                                                     \
1655:     /* xright is the offset of the end of the  local subdomain within its grid row (might actually be outside the local subdomain) */         \
1656:     *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright;                                                                          \
1657:     /* Now compute the size of the local subdomain n. */ \
1658:     *n = 0;                                               \
1659:     if (*ylow_loc < *yhigh_loc) {                           \
1660:       PetscInt width = xright-xleft;                     \
1661:       *n += width*(*yhigh_loc-*ylow_loc-1);                 \
1662:       *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \
1663:       *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \
1664:     } \
1665:   }

1667: /*@
1668:    PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1669:    preconditioner for a two-dimensional problem on a regular grid.

1671:    Collective

1673:    Input Parameters:
1674: +  M, N               - the global number of grid points in the x and y directions
1675: .  Mdomains, Ndomains - the global number of subdomains in the x and y directions
1676: .  dof                - degrees of freedom per node
1677: -  overlap            - overlap in mesh lines

1679:    Output Parameters:
1680: +  Nsub - the number of local subdomains created
1681: .  iis  - array of index sets defining inner (nonoverlapping) subdomains
1682: -  ois  - array of index sets defining outer (overlapping, if overlap > 0) subdomains


1685:    Level: advanced

1687: .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), PCGASMSetOverlap()
1688: @*/
1689: PetscErrorCode  PCGASMCreateSubdomains2D(PC pc,PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap,PetscInt *nsub,IS **iis,IS **ois)
1690: {
1692:   PetscMPIInt    size, rank;
1693:   PetscInt       i, j;
1694:   PetscInt       maxheight, maxwidth;
1695:   PetscInt       xstart, xleft, xright, xleft_loc, xright_loc;
1696:   PetscInt       ystart, ylow,  yhigh,  ylow_loc,  yhigh_loc;
1697:   PetscInt       x[2][2], y[2][2], n[2];
1698:   PetscInt       first, last;
1699:   PetscInt       nidx, *idx;
1700:   PetscInt       ii,jj,s,q,d;
1701:   PetscInt       k,kk;
1702:   PetscMPIInt    color;
1703:   MPI_Comm       comm, subcomm;
1704:   IS             **xis = NULL, **is = ois, **is_local = iis;

1707:   PetscObjectGetComm((PetscObject)pc, &comm);
1708:   MPI_Comm_size(comm, &size);
1709:   MPI_Comm_rank(comm, &rank);
1710:   MatGetOwnershipRange(pc->pmat, &first, &last);
1711:   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) "
1712:                                       "does not respect the number of degrees of freedom per grid point %D", first, last, dof);

1714:   /* Determine the number of domains with nonzero intersections with the local ownership range. */
1715:   s      = 0;
1716:   ystart = 0;
1717:   for (j=0; j<Ndomains; ++j) {
1718:     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1719:     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);
1720:     /* Vertical domain limits with an overlap. */
1721:     ylow   = PetscMax(ystart - overlap,0);
1722:     yhigh  = PetscMin(ystart + maxheight + overlap,N);
1723:     xstart = 0;
1724:     for (i=0; i<Mdomains; ++i) {
1725:       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1726:       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);
1727:       /* Horizontal domain limits with an overlap. */
1728:       xleft  = PetscMax(xstart - overlap,0);
1729:       xright = PetscMin(xstart + maxwidth + overlap,M);
1730:       /*
1731:          Determine whether this subdomain intersects this processor's ownership range of pc->pmat.
1732:       */
1733:       PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1734:       if (nidx) ++s;
1735:       xstart += maxwidth;
1736:     } /* for (i = 0; i < Mdomains; ++i) */
1737:     ystart += maxheight;
1738:   } /* for (j = 0; j < Ndomains; ++j) */

1740:   /* Now we can allocate the necessary number of ISs. */
1741:   *nsub  = s;
1742:   PetscMalloc1(*nsub,is);
1743:   PetscMalloc1(*nsub,is_local);
1744:   s      = 0;
1745:   ystart = 0;
1746:   for (j=0; j<Ndomains; ++j) {
1747:     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1748:     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);
1749:     /* Vertical domain limits with an overlap. */
1750:     y[0][0] = PetscMax(ystart - overlap,0);
1751:     y[0][1] = PetscMin(ystart + maxheight + overlap,N);
1752:     /* Vertical domain limits without an overlap. */
1753:     y[1][0] = ystart;
1754:     y[1][1] = PetscMin(ystart + maxheight,N);
1755:     xstart  = 0;
1756:     for (i=0; i<Mdomains; ++i) {
1757:       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1758:       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);
1759:       /* Horizontal domain limits with an overlap. */
1760:       x[0][0] = PetscMax(xstart - overlap,0);
1761:       x[0][1] = PetscMin(xstart + maxwidth + overlap,M);
1762:       /* Horizontal domain limits without an overlap. */
1763:       x[1][0] = xstart;
1764:       x[1][1] = PetscMin(xstart+maxwidth,M);
1765:       /*
1766:          Determine whether this domain intersects this processor's ownership range of pc->pmat.
1767:          Do this twice: first for the domains with overlaps, and once without.
1768:          During the first pass create the subcommunicators, and use them on the second pass as well.
1769:       */
1770:       for (q = 0; q < 2; ++q) {
1771:         PetscBool split = PETSC_FALSE;
1772:         /*
1773:           domain limits, (xleft, xright) and (ylow, yheigh) are adjusted
1774:           according to whether the domain with an overlap or without is considered.
1775:         */
1776:         xleft = x[q][0]; xright = x[q][1];
1777:         ylow  = y[q][0]; yhigh  = y[q][1];
1778:         PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1779:         nidx *= dof;
1780:         n[q]  = nidx;
1781:         /*
1782:          Based on the counted number of indices in the local domain *with an overlap*,
1783:          construct a subcommunicator of all the processors supporting this domain.
1784:          Observe that a domain with an overlap might have nontrivial local support,
1785:          while the domain without an overlap might not.  Hence, the decision to participate
1786:          in the subcommunicator must be based on the domain with an overlap.
1787:          */
1788:         if (q == 0) {
1789:           if (nidx) color = 1;
1790:           else color = MPI_UNDEFINED;
1791:           MPI_Comm_split(comm, color, rank, &subcomm);
1792:           split = PETSC_TRUE;
1793:         }
1794:         /*
1795:          Proceed only if the number of local indices *with an overlap* is nonzero.
1796:          */
1797:         if (n[0]) {
1798:           if (q == 0) xis = is;
1799:           if (q == 1) {
1800:             /*
1801:              The IS for the no-overlap subdomain shares a communicator with the overlapping domain.
1802:              Moreover, if the overlap is zero, the two ISs are identical.
1803:              */
1804:             if (overlap == 0) {
1805:               (*is_local)[s] = (*is)[s];
1806:               PetscObjectReference((PetscObject)(*is)[s]);
1807:               continue;
1808:             } else {
1809:               xis     = is_local;
1810:               subcomm = ((PetscObject)(*is)[s])->comm;
1811:             }
1812:           } /* if (q == 1) */
1813:           idx  = NULL;
1814:           PetscMalloc1(nidx,&idx);
1815:           if (nidx) {
1816:             k = 0;
1817:             for (jj=ylow_loc; jj<yhigh_loc; ++jj) {
1818:               PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft;
1819:               PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright;
1820:               kk = dof*(M*jj + x0);
1821:               for (ii=x0; ii<x1; ++ii) {
1822:                 for (d = 0; d < dof; ++d) {
1823:                   idx[k++] = kk++;
1824:                 }
1825:               }
1826:             }
1827:           }
1828:           ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);
1829:           if (split) {
1830:             MPI_Comm_free(&subcomm);
1831:           }
1832:         }/* if (n[0]) */
1833:       }/* for (q = 0; q < 2; ++q) */
1834:       if (n[0]) ++s;
1835:       xstart += maxwidth;
1836:     } /* for (i = 0; i < Mdomains; ++i) */
1837:     ystart += maxheight;
1838:   } /* for (j = 0; j < Ndomains; ++j) */
1839:   return(0);
1840: }

1842: /*@C
1843:     PCGASMGetSubdomains - Gets the subdomains supported on this processor
1844:     for the additive Schwarz preconditioner.

1846:     Not Collective

1848:     Input Parameter:
1849: .   pc - the preconditioner context

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


1857:     Notes:
1858:     The user is responsible for destroying the ISs and freeing the returned arrays.
1859:     The IS numbering is in the parallel, global numbering of the vector.

1861:     Level: advanced

1863: .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), PCGASMCreateSubdomains2D(),
1864:           PCGASMSetSubdomains(), PCGASMGetSubmatrices()
1865: @*/
1866: PetscErrorCode  PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[])
1867: {
1868:   PC_GASM        *osm;
1870:   PetscBool      match;
1871:   PetscInt       i;

1875:   PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);
1876:   if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1877:   osm = (PC_GASM*)pc->data;
1878:   if (n) *n = osm->n;
1879:   if (iis) {
1880:     PetscMalloc1(osm->n, iis);
1881:   }
1882:   if (ois) {
1883:     PetscMalloc1(osm->n, ois);
1884:   }
1885:   if (iis || ois) {
1886:     for (i = 0; i < osm->n; ++i) {
1887:       if (iis) (*iis)[i] = osm->iis[i];
1888:       if (ois) (*ois)[i] = osm->ois[i];
1889:     }
1890:   }
1891:   return(0);
1892: }

1894: /*@C
1895:     PCGASMGetSubmatrices - Gets the local submatrices (for this processor
1896:     only) for the additive Schwarz preconditioner.

1898:     Not Collective

1900:     Input Parameter:
1901: .   pc - the preconditioner context

1903:     Output Parameters:
1904: +   n   - the number of matrices for this processor (default value = 1)
1905: -   mat - the matrices

1907:     Notes:
1908:     matrices returned by this routine have the same communicators as the index sets (IS)
1909:            used to define subdomains in PCGASMSetSubdomains()
1910:     Level: advanced

1912: .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(),
1913:           PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains()
1914: @*/
1915: PetscErrorCode  PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1916: {
1917:   PC_GASM        *osm;
1919:   PetscBool      match;

1925:   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUp() or PCSetUp().");
1926:   PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);
1927:   if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1928:   osm = (PC_GASM*)pc->data;
1929:   if (n) *n = osm->n;
1930:   if (mat) *mat = osm->pmat;
1931:   return(0);
1932: }

1934: /*@
1935:     PCGASMSetUseDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1936:     Logically Collective

1938:     Input Parameter:
1939: +   pc  - the preconditioner
1940: -   flg - boolean indicating whether to use subdomains defined by the DM

1942:     Options Database Key:
1943: .   -pc_gasm_dm_subdomains -pc_gasm_overlap -pc_gasm_total_subdomains

1945:     Level: intermediate

1947:     Notes:
1948:     PCGASMSetSubdomains(), PCGASMSetTotalSubdomains() or PCGASMSetOverlap() take precedence over PCGASMSetUseDMSubdomains(),
1949:     so setting PCGASMSetSubdomains() with nontrivial subdomain ISs or any of PCGASMSetTotalSubdomains() and PCGASMSetOverlap()
1950:     automatically turns the latter off.

1952: .seealso: PCGASMGetUseDMSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap()
1953:           PCGASMCreateSubdomains2D()
1954: @*/
1955: PetscErrorCode  PCGASMSetUseDMSubdomains(PC pc,PetscBool flg)
1956: {
1957:   PC_GASM        *osm = (PC_GASM*)pc->data;
1959:   PetscBool      match;

1964:   if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
1965:   PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);
1966:   if (match) {
1967:     if (!osm->user_subdomains && osm->N == PETSC_DETERMINE && osm->overlap < 0) {
1968:       osm->dm_subdomains = flg;
1969:     }
1970:   }
1971:   return(0);
1972: }

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

1978:     Input Parameter:
1979: .   pc  - the preconditioner

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

1984:     Level: intermediate

1986: .seealso: PCGASMSetUseDMSubdomains(), PCGASMSetOverlap()
1987:           PCGASMCreateSubdomains2D()
1988: @*/
1989: PetscErrorCode  PCGASMGetUseDMSubdomains(PC pc,PetscBool* flg)
1990: {
1991:   PC_GASM        *osm = (PC_GASM*)pc->data;
1993:   PetscBool      match;

1998:   PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);
1999:   if (match) {
2000:     if (flg) *flg = osm->dm_subdomains;
2001:   }
2002:   return(0);
2003: }