Actual source code: gasm.c
petsc-3.14.0 2020-09-29
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: }
253: PetscViewerFlush(viewer);
254: PetscViewerASCIIPopSynchronized(viewer);
255: return(0);
256: }
258: PETSC_INTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]);
260: PetscErrorCode PCGASMSetHierarchicalPartitioning(PC pc)
261: {
262: PC_GASM *osm = (PC_GASM*)pc->data;
263: MatPartitioning part;
264: MPI_Comm comm;
265: PetscMPIInt size;
266: PetscInt nlocalsubdomains,fromrows_localsize;
267: IS partitioning,fromrows,isn;
268: Vec outervec;
269: PetscErrorCode ierr;
272: PetscObjectGetComm((PetscObject)pc,&comm);
273: MPI_Comm_size(comm,&size);
274: /* we do not need a hierarchical partitioning when
275: * the total number of subdomains is consistent with
276: * the number of MPI tasks.
277: * For the following cases, we do not need to use HP
278: * */
279: if (osm->N==PETSC_DETERMINE || osm->N>=size || osm->N==1) return(0);
280: 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);
281: nlocalsubdomains = size/osm->N;
282: osm->n = 1;
283: MatPartitioningCreate(comm,&part);
284: MatPartitioningSetAdjacency(part,pc->pmat);
285: MatPartitioningSetType(part,MATPARTITIONINGHIERARCH);
286: MatPartitioningHierarchicalSetNcoarseparts(part,osm->N);
287: MatPartitioningHierarchicalSetNfineparts(part,nlocalsubdomains);
288: MatPartitioningSetFromOptions(part);
289: /* get new processor owner number of each vertex */
290: MatPartitioningApply(part,&partitioning);
291: ISBuildTwoSided(partitioning,NULL,&fromrows);
292: ISPartitioningToNumbering(partitioning,&isn);
293: ISDestroy(&isn);
294: ISGetLocalSize(fromrows,&fromrows_localsize);
295: MatPartitioningDestroy(&part);
296: MatCreateVecs(pc->pmat,&outervec,NULL);
297: VecCreateMPI(comm,fromrows_localsize,PETSC_DETERMINE,&(osm->pcx));
298: VecDuplicate(osm->pcx,&(osm->pcy));
299: VecScatterCreate(osm->pcx,NULL,outervec,fromrows,&(osm->pctoouter));
300: MatCreateSubMatrix(pc->pmat,fromrows,fromrows,MAT_INITIAL_MATRIX,&(osm->permutationP));
301: PetscObjectReference((PetscObject)fromrows);
302: osm->permutationIS = fromrows;
303: osm->pcmat = pc->pmat;
304: PetscObjectReference((PetscObject)osm->permutationP);
305: pc->pmat = osm->permutationP;
306: VecDestroy(&outervec);
307: ISDestroy(&fromrows);
308: ISDestroy(&partitioning);
309: osm->n = PETSC_DETERMINE;
310: return(0);
311: }
313: static PetscErrorCode PCSetUp_GASM(PC pc)
314: {
315: PC_GASM *osm = (PC_GASM*)pc->data;
317: PetscInt i,nInnerIndices,nTotalInnerIndices;
318: PetscMPIInt rank, size;
319: MatReuse scall = MAT_REUSE_MATRIX;
320: KSP ksp;
321: PC subpc;
322: const char *prefix,*pprefix;
323: Vec x,y;
324: PetscInt oni; /* Number of indices in the i-th local outer subdomain. */
325: const PetscInt *oidxi; /* Indices from the i-th subdomain local outer subdomain. */
326: PetscInt on; /* Number of indices in the disjoint union of local outer subdomains. */
327: PetscInt *oidx; /* Indices in the disjoint union of local outer subdomains. */
328: IS gois; /* Disjoint union the global indices of outer subdomains. */
329: IS goid; /* Identity IS of the size of the disjoint union of outer subdomains. */
330: PetscScalar *gxarray, *gyarray;
331: PetscInt gostart; /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- over the disjoint union of outer subdomains. */
332: PetscInt num_subdomains = 0;
333: DM *subdomain_dm = NULL;
334: char **subdomain_names = NULL;
335: PetscInt *numbering;
339: MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
340: MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
341: if (!pc->setupcalled) {
342: /* use a hierarchical partitioning */
343: if (osm->hierarchicalpartitioning){
344: PCGASMSetHierarchicalPartitioning(pc);
345: }
346: if (osm->n == PETSC_DETERMINE) {
347: if (osm->N != PETSC_DETERMINE) {
348: /* No local subdomains given, but the desired number of total subdomains is known, so construct them accordingly. */
349: PCGASMCreateSubdomains(pc->pmat,osm->N,&osm->n,&osm->iis);
350: } else if (osm->dm_subdomains && pc->dm) {
351: /* try pc->dm next, if allowed */
352: PetscInt d;
353: IS *inner_subdomain_is, *outer_subdomain_is;
354: DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm);
355: if (num_subdomains) {
356: PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is);
357: }
358: for (d = 0; d < num_subdomains; ++d) {
359: if (inner_subdomain_is) {ISDestroy(&inner_subdomain_is[d]);}
360: if (outer_subdomain_is) {ISDestroy(&outer_subdomain_is[d]);}
361: }
362: PetscFree(inner_subdomain_is);
363: PetscFree(outer_subdomain_is);
364: } else {
365: /* still no subdomains; use one per processor */
366: osm->nmax = osm->n = 1;
367: MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
368: osm->N = size;
369: PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);
370: }
371: }
372: if (!osm->iis) {
373: /*
374: osm->n was set in PCGASMSetSubdomains(), but the actual subdomains have not been supplied.
375: We create the requisite number of local inner subdomains and then expand them into
376: out subdomains, if necessary.
377: */
378: PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);
379: }
380: if (!osm->ois) {
381: /*
382: Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap
383: has been requested, copy the inner subdomains over so they can be modified.
384: */
385: PetscMalloc1(osm->n,&osm->ois);
386: for (i=0; i<osm->n; ++i) {
387: if (osm->overlap > 0 && osm->N>1) { /* With positive overlap, osm->iis[i] will be modified */
388: ISDuplicate(osm->iis[i],(osm->ois)+i);
389: ISCopy(osm->iis[i],osm->ois[i]);
390: } else {
391: PetscObjectReference((PetscObject)((osm->iis)[i]));
392: osm->ois[i] = osm->iis[i];
393: }
394: }
395: if (osm->overlap>0 && osm->N>1) {
396: /* Extend the "overlapping" regions by a number of steps */
397: MatIncreaseOverlapSplit(pc->pmat,osm->n,osm->ois,osm->overlap);
398: }
399: }
401: /* Now the subdomains are defined. Determine their global and max local numbers, if necessary. */
402: if (osm->nmax == PETSC_DETERMINE) {
403: PetscMPIInt inwork,outwork;
404: /* determine global number of subdomains and the max number of local subdomains */
405: inwork = osm->n;
406: MPIU_Allreduce(&inwork,&outwork,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)pc));
407: osm->nmax = outwork;
408: }
409: if (osm->N == PETSC_DETERMINE) {
410: /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */
411: PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,&osm->N,NULL);
412: }
414: if (osm->sort_indices) {
415: for (i=0; i<osm->n; i++) {
416: ISSort(osm->ois[i]);
417: ISSort(osm->iis[i]);
418: }
419: }
420: PCGetOptionsPrefix(pc,&prefix);
421: PCGASMPrintSubdomains(pc);
423: /*
424: Merge the ISs, create merged vectors and restrictions.
425: */
426: /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */
427: on = 0;
428: for (i=0; i<osm->n; i++) {
429: ISGetLocalSize(osm->ois[i],&oni);
430: on += oni;
431: }
432: PetscMalloc1(on, &oidx);
433: on = 0;
434: /* Merge local indices together */
435: for (i=0; i<osm->n; i++) {
436: ISGetLocalSize(osm->ois[i],&oni);
437: ISGetIndices(osm->ois[i],&oidxi);
438: PetscArraycpy(oidx+on,oidxi,oni);
439: ISRestoreIndices(osm->ois[i],&oidxi);
440: on += oni;
441: }
442: ISCreateGeneral(((PetscObject)(pc))->comm,on,oidx,PETSC_OWN_POINTER,&gois);
443: nTotalInnerIndices = 0;
444: for (i=0; i<osm->n; i++){
445: ISGetLocalSize(osm->iis[i],&nInnerIndices);
446: nTotalInnerIndices += nInnerIndices;
447: }
448: VecCreateMPI(((PetscObject)(pc))->comm,nTotalInnerIndices,PETSC_DETERMINE,&x);
449: VecDuplicate(x,&y);
451: VecCreateMPI(PetscObjectComm((PetscObject)pc),on,PETSC_DECIDE,&osm->gx);
452: VecDuplicate(osm->gx,&osm->gy);
453: VecGetOwnershipRange(osm->gx, &gostart, NULL);
454: ISCreateStride(PetscObjectComm((PetscObject)pc),on,gostart,1, &goid);
455: /* gois might indices not on local */
456: VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction));
457: PetscMalloc1(osm->n,&numbering);
458: PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,NULL,numbering);
459: VecDestroy(&x);
460: ISDestroy(&gois);
462: /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */
463: {
464: PetscInt ini; /* Number of indices the i-th a local inner subdomain. */
465: PetscInt in; /* Number of indices in the disjoint union of local inner subdomains. */
466: PetscInt *iidx; /* Global indices in the merged local inner subdomain. */
467: PetscInt *ioidx; /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
468: IS giis; /* IS for the disjoint union of inner subdomains. */
469: IS giois; /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
470: PetscScalar *array;
471: const PetscInt *indices;
472: PetscInt k;
473: on = 0;
474: for (i=0; i<osm->n; i++) {
475: ISGetLocalSize(osm->ois[i],&oni);
476: on += oni;
477: }
478: PetscMalloc1(on, &iidx);
479: PetscMalloc1(on, &ioidx);
480: VecGetArray(y,&array);
481: /* set communicator id to determine where overlap is */
482: in = 0;
483: for (i=0; i<osm->n; i++) {
484: ISGetLocalSize(osm->iis[i],&ini);
485: for (k = 0; k < ini; ++k){
486: array[in+k] = numbering[i];
487: }
488: in += ini;
489: }
490: VecRestoreArray(y,&array);
491: VecScatterBegin(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD);
492: VecScatterEnd(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD);
493: VecGetOwnershipRange(osm->gy,&gostart, NULL);
494: VecGetArray(osm->gy,&array);
495: on = 0;
496: in = 0;
497: for (i=0; i<osm->n; i++) {
498: ISGetLocalSize(osm->ois[i],&oni);
499: ISGetIndices(osm->ois[i],&indices);
500: for (k=0; k<oni; k++) {
501: /* skip overlapping indices to get inner domain */
502: if (PetscRealPart(array[on+k]) != numbering[i]) continue;
503: iidx[in] = indices[k];
504: ioidx[in++] = gostart+on+k;
505: }
506: ISRestoreIndices(osm->ois[i], &indices);
507: on += oni;
508: }
509: VecRestoreArray(osm->gy,&array);
510: ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,iidx,PETSC_OWN_POINTER,&giis);
511: ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,ioidx,PETSC_OWN_POINTER,&giois);
512: VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction);
513: VecDestroy(&y);
514: ISDestroy(&giis);
515: ISDestroy(&giois);
516: }
517: ISDestroy(&goid);
518: PetscFree(numbering);
520: /* Create the subdomain work vectors. */
521: PetscMalloc1(osm->n,&osm->x);
522: PetscMalloc1(osm->n,&osm->y);
523: VecGetArray(osm->gx, &gxarray);
524: VecGetArray(osm->gy, &gyarray);
525: for (i=0, on=0; i<osm->n; ++i, on += oni) {
526: PetscInt oNi;
527: ISGetLocalSize(osm->ois[i],&oni);
528: /* on a sub communicator */
529: ISGetSize(osm->ois[i],&oNi);
530: VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i]);
531: VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i]);
532: }
533: VecRestoreArray(osm->gx, &gxarray);
534: VecRestoreArray(osm->gy, &gyarray);
535: /* Create the subdomain solvers */
536: PetscMalloc1(osm->n,&osm->ksp);
537: for (i=0; i<osm->n; i++) {
538: char subprefix[PETSC_MAX_PATH_LEN+1];
539: KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp);
540: KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);
541: PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);
542: PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);
543: KSPSetType(ksp,KSPPREONLY);
544: KSPGetPC(ksp,&subpc); /* Why do we need this here? */
545: if (subdomain_dm) {
546: KSPSetDM(ksp,subdomain_dm[i]);
547: DMDestroy(subdomain_dm+i);
548: }
549: PCGetOptionsPrefix(pc,&prefix);
550: KSPSetOptionsPrefix(ksp,prefix);
551: if (subdomain_names && subdomain_names[i]) {
552: PetscSNPrintf(subprefix,PETSC_MAX_PATH_LEN,"sub_%s_",subdomain_names[i]);
553: KSPAppendOptionsPrefix(ksp,subprefix);
554: PetscFree(subdomain_names[i]);
555: }
556: KSPAppendOptionsPrefix(ksp,"sub_");
557: osm->ksp[i] = ksp;
558: }
559: PetscFree(subdomain_dm);
560: PetscFree(subdomain_names);
561: scall = MAT_INITIAL_MATRIX;
562: } else { /* if (pc->setupcalled) */
563: /*
564: Destroy the submatrices from the previous iteration
565: */
566: if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
567: MatDestroyMatrices(osm->n,&osm->pmat);
568: scall = MAT_INITIAL_MATRIX;
569: }
570: if (osm->permutationIS){
571: MatCreateSubMatrix(pc->pmat,osm->permutationIS,osm->permutationIS,scall,&osm->permutationP);
572: PetscObjectReference((PetscObject)osm->permutationP);
573: osm->pcmat = pc->pmat;
574: pc->pmat = osm->permutationP;
575: }
576: }
578: /*
579: Extract the submatrices.
580: */
581: if (size > 1) {
582: MatCreateSubMatricesMPI(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);
583: } else {
584: MatCreateSubMatrices(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);
585: }
586: if (scall == MAT_INITIAL_MATRIX) {
587: PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);
588: for (i=0; i<osm->n; i++) {
589: PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);
590: PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);
591: }
592: }
594: /* Return control to the user so that the submatrices can be modified (e.g., to apply
595: different boundary conditions for the submatrices than for the global problem) */
596: PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP);
598: /*
599: Loop over submatrices putting them into local ksps
600: */
601: for (i=0; i<osm->n; i++) {
602: KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);
603: if (!pc->setupcalled) {
604: KSPSetFromOptions(osm->ksp[i]);
605: }
606: }
607: if (osm->pcmat){
608: MatDestroy(&pc->pmat);
609: pc->pmat = osm->pcmat;
610: osm->pcmat = NULL;
611: }
612: return(0);
613: }
615: static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc)
616: {
617: PC_GASM *osm = (PC_GASM*)pc->data;
619: PetscInt i;
622: for (i=0; i<osm->n; i++) {
623: KSPSetUp(osm->ksp[i]);
624: }
625: return(0);
626: }
628: static PetscErrorCode PCApply_GASM(PC pc,Vec xin,Vec yout)
629: {
630: PC_GASM *osm = (PC_GASM*)pc->data;
632: PetscInt i;
633: Vec x,y;
634: ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
637: if (osm->pctoouter) {
638: VecScatterBegin(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);
639: VecScatterEnd(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);
640: x = osm->pcx;
641: y = osm->pcy;
642: } else {
643: x = xin;
644: y = yout;
645: }
646: /*
647: support for limiting the restriction or interpolation only to the inner
648: subdomain values (leaving the other values 0).
649: */
650: if (!(osm->type & PC_GASM_RESTRICT)) {
651: /* have to zero the work RHS since scatter may leave some slots empty */
652: VecZeroEntries(osm->gx);
653: VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);
654: } else {
655: VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);
656: }
657: VecZeroEntries(osm->gy);
658: if (!(osm->type & PC_GASM_RESTRICT)) {
659: VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);
660: } else {
661: VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);
662: }
663: /* do the subdomain solves */
664: for (i=0; i<osm->n; ++i) {
665: KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);
666: KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);
667: }
668: /* do we need to zero y? */
669: VecZeroEntries(y);
670: if (!(osm->type & PC_GASM_INTERPOLATE)) {
671: VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);
672: VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);
673: } else {
674: VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);
675: VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);
676: }
677: if (osm->pctoouter) {
678: VecScatterBegin(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);
679: VecScatterEnd(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);
680: }
681: return(0);
682: }
684: static PetscErrorCode PCMatApply_GASM(PC pc,Mat Xin,Mat Yout)
685: {
686: PC_GASM *osm = (PC_GASM*)pc->data;
687: Mat X,Y,O=NULL,Z,W;
688: Vec x,y;
689: PetscInt i,m,M,N;
690: ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
694: if (osm->n != 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not yet implemented");
695: MatGetSize(Xin,NULL,&N);
696: if (osm->pctoouter) {
697: VecGetLocalSize(osm->pcx,&m);
698: VecGetSize(osm->pcx,&M);
699: MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]),m,PETSC_DECIDE,M,N,NULL,&O);
700: for (i = 0; i < N; ++i) {
701: MatDenseGetColumnVecRead(Xin,i,&x);
702: MatDenseGetColumnVecWrite(O,i,&y);
703: VecScatterBegin(osm->pctoouter,x,y,INSERT_VALUES,SCATTER_REVERSE);
704: VecScatterEnd(osm->pctoouter,x,y,INSERT_VALUES,SCATTER_REVERSE);
705: MatDenseRestoreColumnVecWrite(O,i,&y);
706: MatDenseRestoreColumnVecRead(Xin,i,&x);
707: }
708: X = Y = O;
709: } else {
710: X = Xin;
711: Y = Yout;
712: }
713: /*
714: support for limiting the restriction or interpolation only to the inner
715: subdomain values (leaving the other values 0).
716: */
717: VecGetLocalSize(osm->x[0],&m);
718: VecGetSize(osm->x[0],&M);
719: MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]),m,PETSC_DECIDE,M,N,NULL,&Z);
720: for (i = 0; i < N; ++i) {
721: MatDenseGetColumnVecRead(X,i,&x);
722: MatDenseGetColumnVecWrite(Z,i,&y);
723: if (!(osm->type & PC_GASM_RESTRICT)) {
724: /* have to zero the work RHS since scatter may leave some slots empty */
725: VecZeroEntries(y);
726: VecScatterBegin(osm->girestriction,x,y,INSERT_VALUES,forward);
727: VecScatterEnd(osm->girestriction,x,y,INSERT_VALUES,forward);
728: } else {
729: VecScatterBegin(osm->gorestriction,x,y,INSERT_VALUES,forward);
730: VecScatterEnd(osm->gorestriction,x,y,INSERT_VALUES,forward);
731: }
732: MatDenseRestoreColumnVecWrite(Z,i,&y);
733: MatDenseRestoreColumnVecRead(X,i,&x);
734: }
735: MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]),m,PETSC_DECIDE,M,N,NULL,&W);
736: MatSetOption(Z,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);
737: MatAssemblyBegin(Z,MAT_FINAL_ASSEMBLY);
738: MatAssemblyEnd(Z,MAT_FINAL_ASSEMBLY);
739: /* do the subdomain solve */
740: KSPMatSolve(osm->ksp[0],Z,W);
741: KSPCheckSolve(osm->ksp[0],pc,NULL);
742: MatDestroy(&Z);
743: /* do we need to zero y? */
744: MatZeroEntries(Y);
745: for (i = 0; i < N; ++i) {
746: MatDenseGetColumnVecWrite(Y,i,&y);
747: MatDenseGetColumnVecRead(W,i,&x);
748: if (!(osm->type & PC_GASM_INTERPOLATE)) {
749: VecScatterBegin(osm->girestriction,x,y,ADD_VALUES,reverse);
750: VecScatterEnd(osm->girestriction,x,y,ADD_VALUES,reverse);
751: } else {
752: VecScatterBegin(osm->gorestriction,x,y,ADD_VALUES,reverse);
753: VecScatterEnd(osm->gorestriction,x,y,ADD_VALUES,reverse);
754: }
755: MatDenseRestoreColumnVecRead(W,i,&x);
756: if (osm->pctoouter) {
757: MatDenseGetColumnVecWrite(Yout,i,&x);
758: VecScatterBegin(osm->pctoouter,y,x,INSERT_VALUES,SCATTER_FORWARD);
759: VecScatterEnd(osm->pctoouter,y,x,INSERT_VALUES,SCATTER_FORWARD);
760: MatDenseRestoreColumnVecRead(Yout,i,&x);
761: }
762: MatDenseRestoreColumnVecWrite(Y,i,&y);
763: }
764: MatDestroy(&W);
765: MatDestroy(&O);
766: return(0);
767: }
769: static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec xin,Vec yout)
770: {
771: PC_GASM *osm = (PC_GASM*)pc->data;
773: PetscInt i;
774: Vec x,y;
775: ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
778: if (osm->pctoouter){
779: VecScatterBegin(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);
780: VecScatterEnd(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);
781: x = osm->pcx;
782: y = osm->pcy;
783: }else{
784: x = xin;
785: y = yout;
786: }
787: /*
788: Support for limiting the restriction or interpolation to only local
789: subdomain values (leaving the other values 0).
791: Note: these are reversed from the PCApply_GASM() because we are applying the
792: transpose of the three terms
793: */
794: if (!(osm->type & PC_GASM_INTERPOLATE)) {
795: /* have to zero the work RHS since scatter may leave some slots empty */
796: VecZeroEntries(osm->gx);
797: VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);
798: } else {
799: VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);
800: }
801: VecZeroEntries(osm->gy);
802: if (!(osm->type & PC_GASM_INTERPOLATE)) {
803: VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);
804: } else {
805: VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);
806: }
807: /* do the local solves */
808: for (i=0; i<osm->n; ++i) { /* Note that the solves are local, so we can go to osm->n, rather than osm->nmax. */
809: KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);
810: KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);
811: }
812: VecZeroEntries(y);
813: if (!(osm->type & PC_GASM_RESTRICT)) {
814: VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);
815: VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);
816: } else {
817: VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);
818: VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);
819: }
820: if (osm->pctoouter){
821: VecScatterBegin(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);
822: VecScatterEnd(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);
823: }
824: return(0);
825: }
827: static PetscErrorCode PCReset_GASM(PC pc)
828: {
829: PC_GASM *osm = (PC_GASM*)pc->data;
831: PetscInt i;
834: if (osm->ksp) {
835: for (i=0; i<osm->n; i++) {
836: KSPReset(osm->ksp[i]);
837: }
838: }
839: if (osm->pmat) {
840: if (osm->n > 0) {
841: PetscMPIInt size;
842: MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
843: if (size > 1) {
844: /* osm->pmat is created by MatCreateSubMatricesMPI(), cannot use MatDestroySubMatrices() */
845: MatDestroyMatrices(osm->n,&osm->pmat);
846: } else {
847: MatDestroySubMatrices(osm->n,&osm->pmat);
848: }
849: }
850: }
851: if (osm->x) {
852: for (i=0; i<osm->n; i++) {
853: VecDestroy(&osm->x[i]);
854: VecDestroy(&osm->y[i]);
855: }
856: }
857: VecDestroy(&osm->gx);
858: VecDestroy(&osm->gy);
860: VecScatterDestroy(&osm->gorestriction);
861: VecScatterDestroy(&osm->girestriction);
862: if (!osm->user_subdomains) {
863: PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);
864: osm->N = PETSC_DETERMINE;
865: osm->nmax = PETSC_DETERMINE;
866: }
867: if (osm->pctoouter){
868: VecScatterDestroy(&(osm->pctoouter));
869: }
870: if (osm->permutationIS){
871: ISDestroy(&(osm->permutationIS));
872: }
873: if (osm->pcx){
874: VecDestroy(&(osm->pcx));
875: }
876: if (osm->pcy){
877: VecDestroy(&(osm->pcy));
878: }
879: if (osm->permutationP){
880: MatDestroy(&(osm->permutationP));
881: }
882: if (osm->pcmat){
883: MatDestroy(&osm->pcmat);
884: }
885: return(0);
886: }
888: static PetscErrorCode PCDestroy_GASM(PC pc)
889: {
890: PC_GASM *osm = (PC_GASM*)pc->data;
892: PetscInt i;
895: PCReset_GASM(pc);
896: /* PCReset will not destroy subdomains, if user_subdomains is true. */
897: PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);
898: if (osm->ksp) {
899: for (i=0; i<osm->n; i++) {
900: KSPDestroy(&osm->ksp[i]);
901: }
902: PetscFree(osm->ksp);
903: }
904: PetscFree(osm->x);
905: PetscFree(osm->y);
906: PetscFree(pc->data);
907: return(0);
908: }
910: static PetscErrorCode PCSetFromOptions_GASM(PetscOptionItems *PetscOptionsObject,PC pc)
911: {
912: PC_GASM *osm = (PC_GASM*)pc->data;
914: PetscInt blocks,ovl;
915: PetscBool flg;
916: PCGASMType gasmtype;
919: PetscOptionsHead(PetscOptionsObject,"Generalized additive Schwarz options");
920: PetscOptionsBool("-pc_gasm_use_dm_subdomains","If subdomains aren't set, use DMCreateDomainDecomposition() to define subdomains.","PCGASMSetUseDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);
921: PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->N,&blocks,&flg);
922: if (flg) {
923: PCGASMSetTotalSubdomains(pc,blocks);
924: }
925: PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);
926: if (flg) {
927: PCGASMSetOverlap(pc,ovl);
928: osm->dm_subdomains = PETSC_FALSE;
929: }
930: flg = PETSC_FALSE;
931: PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);
932: if (flg) {PCGASMSetType(pc,gasmtype);}
933: PetscOptionsBool("-pc_gasm_use_hierachical_partitioning","use hierarchical partitioning",NULL,osm->hierarchicalpartitioning,&osm->hierarchicalpartitioning,&flg);
934: PetscOptionsTail();
935: return(0);
936: }
938: /*------------------------------------------------------------------------------------*/
940: /*@
941: PCGASMSetTotalSubdomains - sets the total number of subdomains to use across the
942: communicator.
943: Logically collective on pc
945: Input Parameters:
946: + pc - the preconditioner
947: - N - total number of subdomains
950: Level: beginner
952: .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap()
953: PCGASMCreateSubdomains2D()
954: @*/
955: PetscErrorCode PCGASMSetTotalSubdomains(PC pc,PetscInt N)
956: {
957: PC_GASM *osm = (PC_GASM*)pc->data;
958: PetscMPIInt size,rank;
962: if (N < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Total number of subdomains must be 1 or more, got N = %D",N);
963: if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before calling PCSetUp().");
965: PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);
966: osm->ois = osm->iis = NULL;
968: MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
969: MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
970: osm->N = N;
971: osm->n = PETSC_DETERMINE;
972: osm->nmax = PETSC_DETERMINE;
973: osm->dm_subdomains = PETSC_FALSE;
974: return(0);
975: }
977: static PetscErrorCode PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[])
978: {
979: PC_GASM *osm = (PC_GASM*)pc->data;
980: PetscErrorCode ierr;
981: PetscInt i;
984: if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, got n = %D",n);
985: if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp().");
987: PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);
988: osm->iis = osm->ois = NULL;
989: osm->n = n;
990: osm->N = PETSC_DETERMINE;
991: osm->nmax = PETSC_DETERMINE;
992: if (ois) {
993: PetscMalloc1(n,&osm->ois);
994: for (i=0; i<n; i++) {
995: PetscObjectReference((PetscObject)ois[i]);
996: osm->ois[i] = ois[i];
997: }
998: /*
999: Since the user set the outer subdomains, even if nontrivial overlap was requested via PCGASMSetOverlap(),
1000: it will be ignored. To avoid confusion later on (e.g., when viewing the PC), the overlap size is set to -1.
1001: */
1002: osm->overlap = -1;
1003: /* inner subdomains must be provided */
1004: if (!iis) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"inner indices have to be provided \n");
1005: }/* end if */
1006: if (iis) {
1007: PetscMalloc1(n,&osm->iis);
1008: for (i=0; i<n; i++) {
1009: PetscObjectReference((PetscObject)iis[i]);
1010: osm->iis[i] = iis[i];
1011: }
1012: if (!ois) {
1013: osm->ois = NULL;
1014: /* if user does not provide outer indices, we will create the corresponding outer indices using osm->overlap =1 in PCSetUp_GASM */
1015: }
1016: }
1017: if (PetscDefined(USE_DEBUG)) {
1018: PetscInt j,rstart,rend,*covered,lsize;
1019: const PetscInt *indices;
1020: /* check if the inner indices cover and only cover the local portion of the preconditioning matrix */
1021: MatGetOwnershipRange(pc->pmat,&rstart,&rend);
1022: PetscCalloc1(rend-rstart,&covered);
1023: /* check if the current processor owns indices from others */
1024: for (i=0; i<n; i++) {
1025: ISGetIndices(osm->iis[i],&indices);
1026: ISGetLocalSize(osm->iis[i],&lsize);
1027: for (j=0; j<lsize; j++) {
1028: 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]);
1029: 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]);
1030: else covered[indices[j]-rstart] = 1;
1031: }
1032: ISRestoreIndices(osm->iis[i],&indices);
1033: }
1034: /* check if we miss any indices */
1035: for (i=rstart; i<rend; i++) {
1036: if (!covered[i-rstart]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"local entity %d was not covered by inner subdomains",i);
1037: }
1038: PetscFree(covered);
1039: }
1040: if (iis) osm->user_subdomains = PETSC_TRUE;
1041: return(0);
1042: }
1044: static PetscErrorCode PCGASMSetOverlap_GASM(PC pc,PetscInt ovl)
1045: {
1046: PC_GASM *osm = (PC_GASM*)pc->data;
1049: if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
1050: if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp().");
1051: if (!pc->setupcalled) osm->overlap = ovl;
1052: return(0);
1053: }
1055: static PetscErrorCode PCGASMSetType_GASM(PC pc,PCGASMType type)
1056: {
1057: PC_GASM *osm = (PC_GASM*)pc->data;
1060: osm->type = type;
1061: osm->type_set = PETSC_TRUE;
1062: return(0);
1063: }
1065: static PetscErrorCode PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort)
1066: {
1067: PC_GASM *osm = (PC_GASM*)pc->data;
1070: osm->sort_indices = doSort;
1071: return(0);
1072: }
1074: /*
1075: FIXME: This routine might need to be modified now that multiple ranks per subdomain are allowed.
1076: In particular, it would upset the global subdomain number calculation.
1077: */
1078: static PetscErrorCode PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp)
1079: {
1080: PC_GASM *osm = (PC_GASM*)pc->data;
1084: 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");
1086: if (n) *n = osm->n;
1087: if (first) {
1088: MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));
1089: *first -= osm->n;
1090: }
1091: if (ksp) {
1092: /* Assume that local solves are now different; not necessarily
1093: true, though! This flag is used only for PCView_GASM() */
1094: *ksp = osm->ksp;
1095: osm->same_subdomain_solvers = PETSC_FALSE;
1096: }
1097: return(0);
1098: } /* PCGASMGetSubKSP_GASM() */
1100: /*@C
1101: PCGASMSetSubdomains - Sets the subdomains for this processor
1102: for the additive Schwarz preconditioner.
1104: Collective on pc
1106: Input Parameters:
1107: + pc - the preconditioner object
1108: . n - the number of subdomains for this processor
1109: . iis - the index sets that define the inner subdomains (or NULL for PETSc to determine subdomains)
1110: - 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)
1112: Notes:
1113: The IS indices use the parallel, global numbering of the vector entries.
1114: Inner subdomains are those where the correction is applied.
1115: Outer subdomains are those where the residual necessary to obtain the
1116: corrections is obtained (see PCGASMType for the use of inner/outer subdomains).
1117: Both inner and outer subdomains can extend over several processors.
1118: This processor's portion of a subdomain is known as a local subdomain.
1120: Inner subdomains can not overlap with each other, do not have any entities from remote processors,
1121: and have to cover the entire local subdomain owned by the current processor. The index sets on each
1122: process should be ordered such that the ith local subdomain is connected to the ith remote subdomain
1123: on another MPI process.
1125: By default the GASM preconditioner uses 1 (local) subdomain per processor.
1128: Level: advanced
1130: .seealso: PCGASMSetNumSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
1131: PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
1132: @*/
1133: PetscErrorCode PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[])
1134: {
1135: PC_GASM *osm = (PC_GASM*)pc->data;
1140: PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));
1141: osm->dm_subdomains = PETSC_FALSE;
1142: return(0);
1143: }
1145: /*@
1146: PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the
1147: additive Schwarz preconditioner. Either all or no processors in the
1148: pc communicator must call this routine.
1150: Logically Collective on pc
1152: Input Parameters:
1153: + pc - the preconditioner context
1154: - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 0)
1156: Options Database Key:
1157: . -pc_gasm_overlap <overlap> - Sets overlap
1159: Notes:
1160: By default the GASM preconditioner uses 1 subdomain per processor. To use
1161: multiple subdomain per perocessor or "straddling" subdomains that intersect
1162: multiple processors use PCGASMSetSubdomains() (or option -pc_gasm_total_subdomains <n>).
1164: The overlap defaults to 0, so if one desires that no additional
1165: overlap be computed beyond what may have been set with a call to
1166: PCGASMSetSubdomains(), then ovl must be set to be 0. In particular, if one does
1167: not explicitly set the subdomains in application code, then all overlap would be computed
1168: internally by PETSc, and using an overlap of 0 would result in an GASM
1169: variant that is equivalent to the block Jacobi preconditioner.
1171: Note that one can define initial index sets with any overlap via
1172: PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows
1173: PETSc to extend that overlap further, if desired.
1175: Level: intermediate
1177: .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1178: PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
1179: @*/
1180: PetscErrorCode PCGASMSetOverlap(PC pc,PetscInt ovl)
1181: {
1183: PC_GASM *osm = (PC_GASM*)pc->data;
1188: PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));
1189: osm->dm_subdomains = PETSC_FALSE;
1190: return(0);
1191: }
1193: /*@
1194: PCGASMSetType - Sets the type of restriction and interpolation used
1195: for local problems in the additive Schwarz method.
1197: Logically Collective on PC
1199: Input Parameters:
1200: + pc - the preconditioner context
1201: - type - variant of GASM, one of
1202: .vb
1203: PC_GASM_BASIC - full interpolation and restriction
1204: PC_GASM_RESTRICT - full restriction, local processor interpolation
1205: PC_GASM_INTERPOLATE - full interpolation, local processor restriction
1206: PC_GASM_NONE - local processor restriction and interpolation
1207: .ve
1209: Options Database Key:
1210: . -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1212: Level: intermediate
1214: .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1215: PCGASMCreateSubdomains2D()
1216: @*/
1217: PetscErrorCode PCGASMSetType(PC pc,PCGASMType type)
1218: {
1224: PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));
1225: return(0);
1226: }
1228: /*@
1229: PCGASMSetSortIndices - Determines whether subdomain indices are sorted.
1231: Logically Collective on PC
1233: Input Parameters:
1234: + pc - the preconditioner context
1235: - doSort - sort the subdomain indices
1237: Level: intermediate
1239: .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1240: PCGASMCreateSubdomains2D()
1241: @*/
1242: PetscErrorCode PCGASMSetSortIndices(PC pc,PetscBool doSort)
1243: {
1249: PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));
1250: return(0);
1251: }
1253: /*@C
1254: PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on
1255: this processor.
1257: Collective on PC iff first_local is requested
1259: Input Parameter:
1260: . pc - the preconditioner context
1262: Output Parameters:
1263: + n_local - the number of blocks on this processor or NULL
1264: . first_local - the global number of the first block on this processor or NULL,
1265: all processors must request or all must pass NULL
1266: - ksp - the array of KSP contexts
1268: Note:
1269: After PCGASMGetSubKSP() the array of KSPes is not to be freed
1271: Currently for some matrix implementations only 1 block per processor
1272: is supported.
1274: You must call KSPSetUp() before calling PCGASMGetSubKSP().
1276: Level: advanced
1278: .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(),
1279: PCGASMCreateSubdomains2D(),
1280: @*/
1281: PetscErrorCode PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
1282: {
1287: PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));
1288: return(0);
1289: }
1291: /* -------------------------------------------------------------------------------------*/
1292: /*MC
1293: PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
1294: its own KSP object.
1296: Options Database Keys:
1297: + -pc_gasm_total_subdomains <n> - Sets total number of local subdomains to be distributed among processors
1298: . -pc_gasm_view_subdomains - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view
1299: . -pc_gasm_print_subdomains - activates the printing of subdomain indices in PCSetUp()
1300: . -pc_gasm_overlap <ovl> - Sets overlap by which to (automatically) extend local subdomains
1301: - -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1303: IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
1304: will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use
1305: -pc_gasm_type basic to use the standard GASM.
1307: Notes:
1308: Blocks can be shared by multiple processes.
1310: To set options on the solvers for each block append -sub_ to all the KSP, and PC
1311: options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
1313: To set the options on the solvers separate for each block call PCGASMGetSubKSP()
1314: and set the options directly on the resulting KSP object (you can access its PC
1315: with KSPGetPC())
1318: Level: beginner
1320: References:
1321: + 1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions
1322: Courant Institute, New York University Technical report
1323: - 2. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
1324: Cambridge University Press.
1326: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
1327: PCBJACOBI, PCGASMGetSubKSP(), PCGASMSetSubdomains(),
1328: PCSetModifySubMatrices(), PCGASMSetOverlap(), PCGASMSetType()
1330: M*/
1332: PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc)
1333: {
1335: PC_GASM *osm;
1338: PetscNewLog(pc,&osm);
1340: osm->N = PETSC_DETERMINE;
1341: osm->n = PETSC_DECIDE;
1342: osm->nmax = PETSC_DETERMINE;
1343: osm->overlap = 0;
1344: osm->ksp = NULL;
1345: osm->gorestriction = NULL;
1346: osm->girestriction = NULL;
1347: osm->pctoouter = NULL;
1348: osm->gx = NULL;
1349: osm->gy = NULL;
1350: osm->x = NULL;
1351: osm->y = NULL;
1352: osm->pcx = NULL;
1353: osm->pcy = NULL;
1354: osm->permutationIS = NULL;
1355: osm->permutationP = NULL;
1356: osm->pcmat = NULL;
1357: osm->ois = NULL;
1358: osm->iis = NULL;
1359: osm->pmat = NULL;
1360: osm->type = PC_GASM_RESTRICT;
1361: osm->same_subdomain_solvers = PETSC_TRUE;
1362: osm->sort_indices = PETSC_TRUE;
1363: osm->dm_subdomains = PETSC_FALSE;
1364: osm->hierarchicalpartitioning = PETSC_FALSE;
1366: pc->data = (void*)osm;
1367: pc->ops->apply = PCApply_GASM;
1368: pc->ops->matapply = PCMatApply_GASM;
1369: pc->ops->applytranspose = PCApplyTranspose_GASM;
1370: pc->ops->setup = PCSetUp_GASM;
1371: pc->ops->reset = PCReset_GASM;
1372: pc->ops->destroy = PCDestroy_GASM;
1373: pc->ops->setfromoptions = PCSetFromOptions_GASM;
1374: pc->ops->setuponblocks = PCSetUpOnBlocks_GASM;
1375: pc->ops->view = PCView_GASM;
1376: pc->ops->applyrichardson = NULL;
1378: PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSubdomains_C",PCGASMSetSubdomains_GASM);
1379: PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetOverlap_C",PCGASMSetOverlap_GASM);
1380: PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetType_C",PCGASMSetType_GASM);
1381: PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSortIndices_C",PCGASMSetSortIndices_GASM);
1382: PetscObjectComposeFunction((PetscObject)pc,"PCGASMGetSubKSP_C",PCGASMGetSubKSP_GASM);
1383: return(0);
1384: }
1386: PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[])
1387: {
1388: MatPartitioning mpart;
1389: const char *prefix;
1390: PetscInt i,j,rstart,rend,bs;
1391: PetscBool hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
1392: Mat Ad = NULL, adj;
1393: IS ispart,isnumb,*is;
1394: PetscErrorCode ierr;
1397: if (nloc < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local subdomains must > 0, got nloc = %D",nloc);
1399: /* Get prefix, row distribution, and block size */
1400: MatGetOptionsPrefix(A,&prefix);
1401: MatGetOwnershipRange(A,&rstart,&rend);
1402: MatGetBlockSize(A,&bs);
1403: 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);
1405: /* Get diagonal block from matrix if possible */
1406: MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop);
1407: if (hasop) {
1408: MatGetDiagonalBlock(A,&Ad);
1409: }
1410: if (Ad) {
1411: PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);
1412: if (!isbaij) {PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);}
1413: }
1414: if (Ad && nloc > 1) {
1415: PetscBool match,done;
1416: /* Try to setup a good matrix partitioning if available */
1417: MatPartitioningCreate(PETSC_COMM_SELF,&mpart);
1418: PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);
1419: MatPartitioningSetFromOptions(mpart);
1420: PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);
1421: if (!match) {
1422: PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);
1423: }
1424: if (!match) { /* assume a "good" partitioner is available */
1425: PetscInt na;
1426: const PetscInt *ia,*ja;
1427: MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);
1428: if (done) {
1429: /* Build adjacency matrix by hand. Unfortunately a call to
1430: MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1431: remove the block-aij structure and we cannot expect
1432: MatPartitioning to split vertices as we need */
1433: PetscInt i,j,len,nnz,cnt,*iia=NULL,*jja=NULL;
1434: const PetscInt *row;
1435: nnz = 0;
1436: for (i=0; i<na; i++) { /* count number of nonzeros */
1437: len = ia[i+1] - ia[i];
1438: row = ja + ia[i];
1439: for (j=0; j<len; j++) {
1440: if (row[j] == i) { /* don't count diagonal */
1441: len--; break;
1442: }
1443: }
1444: nnz += len;
1445: }
1446: PetscMalloc1(na+1,&iia);
1447: PetscMalloc1(nnz,&jja);
1448: nnz = 0;
1449: iia[0] = 0;
1450: for (i=0; i<na; i++) { /* fill adjacency */
1451: cnt = 0;
1452: len = ia[i+1] - ia[i];
1453: row = ja + ia[i];
1454: for (j=0; j<len; j++) {
1455: if (row[j] != i) jja[nnz+cnt++] = row[j]; /* if not diagonal */
1456: }
1457: nnz += cnt;
1458: iia[i+1] = nnz;
1459: }
1460: /* Partitioning of the adjacency matrix */
1461: MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);
1462: MatPartitioningSetAdjacency(mpart,adj);
1463: MatPartitioningSetNParts(mpart,nloc);
1464: MatPartitioningApply(mpart,&ispart);
1465: ISPartitioningToNumbering(ispart,&isnumb);
1466: MatDestroy(&adj);
1467: foundpart = PETSC_TRUE;
1468: }
1469: MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);
1470: }
1471: MatPartitioningDestroy(&mpart);
1472: }
1473: PetscMalloc1(nloc,&is);
1474: if (!foundpart) {
1476: /* Partitioning by contiguous chunks of rows */
1478: PetscInt mbs = (rend-rstart)/bs;
1479: PetscInt start = rstart;
1480: for (i=0; i<nloc; i++) {
1481: PetscInt count = (mbs/nloc + ((mbs % nloc) > i)) * bs;
1482: ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);
1483: start += count;
1484: }
1486: } else {
1488: /* Partitioning by adjacency of diagonal block */
1490: const PetscInt *numbering;
1491: PetscInt *count,nidx,*indices,*newidx,start=0;
1492: /* Get node count in each partition */
1493: PetscMalloc1(nloc,&count);
1494: ISPartitioningCount(ispart,nloc,count);
1495: if (isbaij && bs > 1) { /* adjust for the block-aij case */
1496: for (i=0; i<nloc; i++) count[i] *= bs;
1497: }
1498: /* Build indices from node numbering */
1499: ISGetLocalSize(isnumb,&nidx);
1500: PetscMalloc1(nidx,&indices);
1501: for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1502: ISGetIndices(isnumb,&numbering);
1503: PetscSortIntWithPermutation(nidx,numbering,indices);
1504: ISRestoreIndices(isnumb,&numbering);
1505: if (isbaij && bs > 1) { /* adjust for the block-aij case */
1506: PetscMalloc1(nidx*bs,&newidx);
1507: for (i=0; i<nidx; i++) {
1508: for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
1509: }
1510: PetscFree(indices);
1511: nidx *= bs;
1512: indices = newidx;
1513: }
1514: /* Shift to get global indices */
1515: for (i=0; i<nidx; i++) indices[i] += rstart;
1517: /* Build the index sets for each block */
1518: for (i=0; i<nloc; i++) {
1519: ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);
1520: ISSort(is[i]);
1521: start += count[i];
1522: }
1524: PetscFree(count);
1525: PetscFree(indices);
1526: ISDestroy(&isnumb);
1527: ISDestroy(&ispart);
1528: }
1529: *iis = is;
1530: return(0);
1531: }
1533: PETSC_INTERN PetscErrorCode PCGASMCreateStraddlingSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[])
1534: {
1535: PetscErrorCode ierr;
1538: MatSubdomainsCreateCoalesce(A,N,n,iis);
1539: return(0);
1540: }
1542: /*@C
1543: PCGASMCreateSubdomains - Creates n index sets defining n nonoverlapping subdomains for the additive
1544: Schwarz preconditioner for a any problem based on its matrix.
1546: Collective
1548: Input Parameters:
1549: + A - The global matrix operator
1550: - N - the number of global subdomains requested
1552: Output Parameters:
1553: + n - the number of subdomains created on this processor
1554: - iis - the array of index sets defining the local inner subdomains (on which the correction is applied)
1556: Level: advanced
1558: Note: When N >= A's communicator size, each subdomain is local -- contained within a single processor.
1559: When N < size, the subdomains are 'straddling' (processor boundaries) and are no longer local.
1560: The resulting subdomains can be use in PCGASMSetSubdomains(pc,n,iss,NULL). The overlapping
1561: outer subdomains will be automatically generated from these according to the requested amount of
1562: overlap; this is currently supported only with local subdomains.
1565: .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains()
1566: @*/
1567: PetscErrorCode PCGASMCreateSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[])
1568: {
1569: PetscMPIInt size;
1570: PetscErrorCode ierr;
1576: if (N < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of subdomains must be > 0, N = %D",N);
1577: MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
1578: if (N >= size) {
1579: *n = N/size + (N%size);
1580: PCGASMCreateLocalSubdomains(A,*n,iis);
1581: } else {
1582: PCGASMCreateStraddlingSubdomains(A,N,n,iis);
1583: }
1584: return(0);
1585: }
1587: /*@C
1588: PCGASMDestroySubdomains - Destroys the index sets created with
1589: PCGASMCreateSubdomains() or PCGASMCreateSubdomains2D. Should be
1590: called after setting subdomains with PCGASMSetSubdomains().
1592: Collective
1594: Input Parameters:
1595: + n - the number of index sets
1596: . iis - the array of inner subdomains,
1597: - ois - the array of outer subdomains, can be NULL
1599: Level: intermediate
1601: Notes:
1602: this is merely a convenience subroutine that walks each list,
1603: destroys each IS on the list, and then frees the list. At the end the
1604: list pointers are set to NULL.
1606: .seealso: PCGASMCreateSubdomains(), PCGASMSetSubdomains()
1607: @*/
1608: PetscErrorCode PCGASMDestroySubdomains(PetscInt n,IS **iis,IS **ois)
1609: {
1610: PetscInt i;
1614: if (n <= 0) return(0);
1615: if (ois) {
1617: if (*ois) {
1619: for (i=0; i<n; i++) {
1620: ISDestroy(&(*ois)[i]);
1621: }
1622: PetscFree((*ois));
1623: }
1624: }
1625: if (iis) {
1627: if (*iis) {
1629: for (i=0; i<n; i++) {
1630: ISDestroy(&(*iis)[i]);
1631: }
1632: PetscFree((*iis));
1633: }
1634: }
1635: return(0);
1636: }
1638: #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \
1639: { \
1640: PetscInt first_row = first/M, last_row = last/M+1; \
1641: /* \
1642: Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners \
1643: of the bounding box of the intersection of the subdomain with the local ownership range (local \
1644: subdomain). \
1645: Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows \
1646: of the intersection. \
1647: */ \
1648: /* ylow_loc is the grid row containing the first element of the local sumbdomain */ \
1649: *ylow_loc = PetscMax(first_row,ylow); \
1650: /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \
1651: *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft; \
1652: /* yhigh_loc is the grid row above the last local subdomain element */ \
1653: *yhigh_loc = PetscMin(last_row,yhigh); \
1654: /* xright is the offset of the end of the local subdomain within its grid row (might actually be outside the local subdomain) */ \
1655: *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright; \
1656: /* Now compute the size of the local subdomain n. */ \
1657: *n = 0; \
1658: if (*ylow_loc < *yhigh_loc) { \
1659: PetscInt width = xright-xleft; \
1660: *n += width*(*yhigh_loc-*ylow_loc-1); \
1661: *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \
1662: *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \
1663: } \
1664: }
1666: /*@
1667: PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1668: preconditioner for a two-dimensional problem on a regular grid.
1670: Collective
1672: Input Parameters:
1673: + M, N - the global number of grid points in the x and y directions
1674: . Mdomains, Ndomains - the global number of subdomains in the x and y directions
1675: . dof - degrees of freedom per node
1676: - overlap - overlap in mesh lines
1678: Output Parameters:
1679: + Nsub - the number of local subdomains created
1680: . iis - array of index sets defining inner (nonoverlapping) subdomains
1681: - ois - array of index sets defining outer (overlapping, if overlap > 0) subdomains
1684: Level: advanced
1686: .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), PCGASMSetOverlap()
1687: @*/
1688: PetscErrorCode PCGASMCreateSubdomains2D(PC pc,PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap,PetscInt *nsub,IS **iis,IS **ois)
1689: {
1691: PetscMPIInt size, rank;
1692: PetscInt i, j;
1693: PetscInt maxheight, maxwidth;
1694: PetscInt xstart, xleft, xright, xleft_loc, xright_loc;
1695: PetscInt ystart, ylow, yhigh, ylow_loc, yhigh_loc;
1696: PetscInt x[2][2], y[2][2], n[2];
1697: PetscInt first, last;
1698: PetscInt nidx, *idx;
1699: PetscInt ii,jj,s,q,d;
1700: PetscInt k,kk;
1701: PetscMPIInt color;
1702: MPI_Comm comm, subcomm;
1703: IS **xis = NULL, **is = ois, **is_local = iis;
1706: PetscObjectGetComm((PetscObject)pc, &comm);
1707: MPI_Comm_size(comm, &size);
1708: MPI_Comm_rank(comm, &rank);
1709: MatGetOwnershipRange(pc->pmat, &first, &last);
1710: 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) "
1711: "does not respect the number of degrees of freedom per grid point %D", first, last, dof);
1713: /* Determine the number of domains with nonzero intersections with the local ownership range. */
1714: s = 0;
1715: ystart = 0;
1716: for (j=0; j<Ndomains; ++j) {
1717: maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1718: 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);
1719: /* Vertical domain limits with an overlap. */
1720: ylow = PetscMax(ystart - overlap,0);
1721: yhigh = PetscMin(ystart + maxheight + overlap,N);
1722: xstart = 0;
1723: for (i=0; i<Mdomains; ++i) {
1724: maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1725: 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);
1726: /* Horizontal domain limits with an overlap. */
1727: xleft = PetscMax(xstart - overlap,0);
1728: xright = PetscMin(xstart + maxwidth + overlap,M);
1729: /*
1730: Determine whether this subdomain intersects this processor's ownership range of pc->pmat.
1731: */
1732: PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1733: if (nidx) ++s;
1734: xstart += maxwidth;
1735: } /* for (i = 0; i < Mdomains; ++i) */
1736: ystart += maxheight;
1737: } /* for (j = 0; j < Ndomains; ++j) */
1739: /* Now we can allocate the necessary number of ISs. */
1740: *nsub = s;
1741: PetscMalloc1(*nsub,is);
1742: PetscMalloc1(*nsub,is_local);
1743: s = 0;
1744: ystart = 0;
1745: for (j=0; j<Ndomains; ++j) {
1746: maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1747: 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);
1748: /* Vertical domain limits with an overlap. */
1749: y[0][0] = PetscMax(ystart - overlap,0);
1750: y[0][1] = PetscMin(ystart + maxheight + overlap,N);
1751: /* Vertical domain limits without an overlap. */
1752: y[1][0] = ystart;
1753: y[1][1] = PetscMin(ystart + maxheight,N);
1754: xstart = 0;
1755: for (i=0; i<Mdomains; ++i) {
1756: maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1757: 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);
1758: /* Horizontal domain limits with an overlap. */
1759: x[0][0] = PetscMax(xstart - overlap,0);
1760: x[0][1] = PetscMin(xstart + maxwidth + overlap,M);
1761: /* Horizontal domain limits without an overlap. */
1762: x[1][0] = xstart;
1763: x[1][1] = PetscMin(xstart+maxwidth,M);
1764: /*
1765: Determine whether this domain intersects this processor's ownership range of pc->pmat.
1766: Do this twice: first for the domains with overlaps, and once without.
1767: During the first pass create the subcommunicators, and use them on the second pass as well.
1768: */
1769: for (q = 0; q < 2; ++q) {
1770: PetscBool split = PETSC_FALSE;
1771: /*
1772: domain limits, (xleft, xright) and (ylow, yheigh) are adjusted
1773: according to whether the domain with an overlap or without is considered.
1774: */
1775: xleft = x[q][0]; xright = x[q][1];
1776: ylow = y[q][0]; yhigh = y[q][1];
1777: PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1778: nidx *= dof;
1779: n[q] = nidx;
1780: /*
1781: Based on the counted number of indices in the local domain *with an overlap*,
1782: construct a subcommunicator of all the processors supporting this domain.
1783: Observe that a domain with an overlap might have nontrivial local support,
1784: while the domain without an overlap might not. Hence, the decision to participate
1785: in the subcommunicator must be based on the domain with an overlap.
1786: */
1787: if (q == 0) {
1788: if (nidx) color = 1;
1789: else color = MPI_UNDEFINED;
1790: MPI_Comm_split(comm, color, rank, &subcomm);
1791: split = PETSC_TRUE;
1792: }
1793: /*
1794: Proceed only if the number of local indices *with an overlap* is nonzero.
1795: */
1796: if (n[0]) {
1797: if (q == 0) xis = is;
1798: if (q == 1) {
1799: /*
1800: The IS for the no-overlap subdomain shares a communicator with the overlapping domain.
1801: Moreover, if the overlap is zero, the two ISs are identical.
1802: */
1803: if (overlap == 0) {
1804: (*is_local)[s] = (*is)[s];
1805: PetscObjectReference((PetscObject)(*is)[s]);
1806: continue;
1807: } else {
1808: xis = is_local;
1809: subcomm = ((PetscObject)(*is)[s])->comm;
1810: }
1811: } /* if (q == 1) */
1812: idx = NULL;
1813: PetscMalloc1(nidx,&idx);
1814: if (nidx) {
1815: k = 0;
1816: for (jj=ylow_loc; jj<yhigh_loc; ++jj) {
1817: PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft;
1818: PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright;
1819: kk = dof*(M*jj + x0);
1820: for (ii=x0; ii<x1; ++ii) {
1821: for (d = 0; d < dof; ++d) {
1822: idx[k++] = kk++;
1823: }
1824: }
1825: }
1826: }
1827: ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);
1828: if (split) {
1829: MPI_Comm_free(&subcomm);
1830: }
1831: }/* if (n[0]) */
1832: }/* for (q = 0; q < 2; ++q) */
1833: if (n[0]) ++s;
1834: xstart += maxwidth;
1835: } /* for (i = 0; i < Mdomains; ++i) */
1836: ystart += maxheight;
1837: } /* for (j = 0; j < Ndomains; ++j) */
1838: return(0);
1839: }
1841: /*@C
1842: PCGASMGetSubdomains - Gets the subdomains supported on this processor
1843: for the additive Schwarz preconditioner.
1845: Not Collective
1847: Input Parameter:
1848: . pc - the preconditioner context
1850: Output Parameters:
1851: + n - the number of subdomains for this processor (default value = 1)
1852: . iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be NULL)
1853: - ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be NULL)
1856: Notes:
1857: The user is responsible for destroying the ISs and freeing the returned arrays.
1858: The IS numbering is in the parallel, global numbering of the vector.
1860: Level: advanced
1862: .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), PCGASMCreateSubdomains2D(),
1863: PCGASMSetSubdomains(), PCGASMGetSubmatrices()
1864: @*/
1865: PetscErrorCode PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[])
1866: {
1867: PC_GASM *osm;
1869: PetscBool match;
1870: PetscInt i;
1874: PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);
1875: if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1876: osm = (PC_GASM*)pc->data;
1877: if (n) *n = osm->n;
1878: if (iis) {
1879: PetscMalloc1(osm->n, iis);
1880: }
1881: if (ois) {
1882: PetscMalloc1(osm->n, ois);
1883: }
1884: if (iis || ois) {
1885: for (i = 0; i < osm->n; ++i) {
1886: if (iis) (*iis)[i] = osm->iis[i];
1887: if (ois) (*ois)[i] = osm->ois[i];
1888: }
1889: }
1890: return(0);
1891: }
1893: /*@C
1894: PCGASMGetSubmatrices - Gets the local submatrices (for this processor
1895: only) for the additive Schwarz preconditioner.
1897: Not Collective
1899: Input Parameter:
1900: . pc - the preconditioner context
1902: Output Parameters:
1903: + n - the number of matrices for this processor (default value = 1)
1904: - mat - the matrices
1906: Notes:
1907: matrices returned by this routine have the same communicators as the index sets (IS)
1908: used to define subdomains in PCGASMSetSubdomains()
1909: Level: advanced
1911: .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(),
1912: PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains()
1913: @*/
1914: PetscErrorCode PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1915: {
1916: PC_GASM *osm;
1918: PetscBool match;
1924: if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUp() or PCSetUp().");
1925: PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);
1926: if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1927: osm = (PC_GASM*)pc->data;
1928: if (n) *n = osm->n;
1929: if (mat) *mat = osm->pmat;
1930: return(0);
1931: }
1933: /*@
1934: PCGASMSetUseDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1935: Logically Collective
1937: Input Parameter:
1938: + pc - the preconditioner
1939: - flg - boolean indicating whether to use subdomains defined by the DM
1941: Options Database Key:
1942: . -pc_gasm_dm_subdomains -pc_gasm_overlap -pc_gasm_total_subdomains
1944: Level: intermediate
1946: Notes:
1947: PCGASMSetSubdomains(), PCGASMSetTotalSubdomains() or PCGASMSetOverlap() take precedence over PCGASMSetUseDMSubdomains(),
1948: so setting PCGASMSetSubdomains() with nontrivial subdomain ISs or any of PCGASMSetTotalSubdomains() and PCGASMSetOverlap()
1949: automatically turns the latter off.
1951: .seealso: PCGASMGetUseDMSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap()
1952: PCGASMCreateSubdomains2D()
1953: @*/
1954: PetscErrorCode PCGASMSetUseDMSubdomains(PC pc,PetscBool flg)
1955: {
1956: PC_GASM *osm = (PC_GASM*)pc->data;
1958: PetscBool match;
1963: if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
1964: PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);
1965: if (match) {
1966: if (!osm->user_subdomains && osm->N == PETSC_DETERMINE && osm->overlap < 0) {
1967: osm->dm_subdomains = flg;
1968: }
1969: }
1970: return(0);
1971: }
1973: /*@
1974: PCGASMGetUseDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1975: Not Collective
1977: Input Parameter:
1978: . pc - the preconditioner
1980: Output Parameter:
1981: . flg - boolean indicating whether to use subdomains defined by the DM
1983: Level: intermediate
1985: .seealso: PCGASMSetUseDMSubdomains(), PCGASMSetOverlap()
1986: PCGASMCreateSubdomains2D()
1987: @*/
1988: PetscErrorCode PCGASMGetUseDMSubdomains(PC pc,PetscBool* flg)
1989: {
1990: PC_GASM *osm = (PC_GASM*)pc->data;
1992: PetscBool match;
1997: PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);
1998: if (match) {
1999: if (flg) *flg = osm->dm_subdomains;
2000: }
2001: return(0);
2002: }