Actual source code: gasm.c
1: /*
2: This file defines an "generalized" additive Schwarz preconditioner for any Mat implementation.
3: In this version each MPI rank may intersect multiple subdomains and any subdomain may
4: intersect multiple MPI ranks. Intersections of subdomains with MPI ranks 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 rank (set in PCGASMSetSubdomains() or calculated in PCGASMSetTotalSubdomains())
9: nmax - maximum number of local subdomains per rank (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;
44: PetscFunctionBegin;
45: /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */
46: PetscCall(PetscMalloc2(osm->n, numbering, osm->n, permutation));
47: PetscCall(PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc), osm->n, (PetscObject *)osm->iis, NULL, *numbering));
48: for (i = 0; i < osm->n; ++i) (*permutation)[i] = i;
49: PetscCall(PetscSortIntWithPermutation(osm->n, *numbering, *permutation));
50: PetscFunctionReturn(PETSC_SUCCESS);
51: }
53: static PetscErrorCode PCGASMSubdomainView_Private(PC pc, PetscInt i, PetscViewer viewer)
54: {
55: PC_GASM *osm = (PC_GASM *)pc->data;
56: PetscInt j, nidx;
57: const PetscInt *idx;
58: PetscViewer sviewer;
59: char *cidx;
61: PetscFunctionBegin;
62: PetscCheck(i >= 0 && i < osm->n, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Invalid subdomain %" PetscInt_FMT ": must nonnegative and less than %" PetscInt_FMT, i, osm->n);
63: /* Inner subdomains. */
64: PetscCall(ISGetLocalSize(osm->iis[i], &nidx));
65: /*
66: No more than 15 characters per index plus a space.
67: PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
68: in case nidx == 0. That will take care of the space for the trailing '\0' as well.
69: For nidx == 0, the whole string 16 '\0'.
70: */
71: #define len 16 * (nidx + 1) + 1
72: PetscCall(PetscMalloc1(len, &cidx));
73: PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, cidx, len, &sviewer));
74: #undef len
75: PetscCall(ISGetIndices(osm->iis[i], &idx));
76: for (j = 0; j < nidx; ++j) PetscCall(PetscViewerStringSPrintf(sviewer, "%" PetscInt_FMT " ", idx[j]));
77: PetscCall(ISRestoreIndices(osm->iis[i], &idx));
78: PetscCall(PetscViewerDestroy(&sviewer));
79: PetscCall(PetscViewerASCIIPrintf(viewer, "Inner subdomain:\n"));
80: PetscCall(PetscViewerFlush(viewer));
81: PetscCall(PetscViewerASCIIPushSynchronized(viewer));
82: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx));
83: PetscCall(PetscViewerFlush(viewer));
84: PetscCall(PetscViewerASCIIPopSynchronized(viewer));
85: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
86: PetscCall(PetscViewerFlush(viewer));
87: PetscCall(PetscFree(cidx));
88: /* Outer subdomains. */
89: PetscCall(ISGetLocalSize(osm->ois[i], &nidx));
90: /*
91: No more than 15 characters per index plus a space.
92: PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
93: in case nidx == 0. That will take care of the space for the trailing '\0' as well.
94: For nidx == 0, the whole string 16 '\0'.
95: */
96: #define len 16 * (nidx + 1) + 1
97: PetscCall(PetscMalloc1(len, &cidx));
98: PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, cidx, len, &sviewer));
99: #undef len
100: PetscCall(ISGetIndices(osm->ois[i], &idx));
101: for (j = 0; j < nidx; ++j) PetscCall(PetscViewerStringSPrintf(sviewer, "%" PetscInt_FMT " ", idx[j]));
102: PetscCall(PetscViewerDestroy(&sviewer));
103: PetscCall(ISRestoreIndices(osm->ois[i], &idx));
104: PetscCall(PetscViewerASCIIPrintf(viewer, "Outer subdomain:\n"));
105: PetscCall(PetscViewerFlush(viewer));
106: PetscCall(PetscViewerASCIIPushSynchronized(viewer));
107: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx));
108: PetscCall(PetscViewerFlush(viewer));
109: PetscCall(PetscViewerASCIIPopSynchronized(viewer));
110: PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
111: PetscCall(PetscViewerFlush(viewer));
112: PetscCall(PetscFree(cidx));
113: PetscFunctionReturn(PETSC_SUCCESS);
114: }
116: static PetscErrorCode PCGASMPrintSubdomains(PC pc)
117: {
118: PC_GASM *osm = (PC_GASM *)pc->data;
119: const char *prefix;
120: char fname[PETSC_MAX_PATH_LEN + 1];
121: PetscInt l, d, count;
122: PetscBool found;
123: PetscViewer viewer, sviewer = NULL;
124: PetscInt *numbering, *permutation; /* global numbering of locally-supported subdomains and the permutation from the local ordering */
126: PetscFunctionBegin;
127: PetscCall(PCGetOptionsPrefix(pc, &prefix));
128: PetscCall(PetscOptionsHasName(NULL, prefix, "-pc_gasm_print_subdomains", &found));
129: if (!found) PetscFunctionReturn(PETSC_SUCCESS);
130: PetscCall(PetscOptionsGetString(NULL, prefix, "-pc_gasm_print_subdomains", fname, sizeof(fname), &found));
131: if (!found) PetscCall(PetscStrncpy(fname, "stdout", sizeof(fname)));
132: PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc), fname, &viewer));
133: /*
134: Make sure the viewer has a name. Otherwise this may cause a deadlock or other weird errors when creating a subcomm viewer:
135: the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed collectively on the comm.
136: */
137: PetscCall(PetscObjectName((PetscObject)viewer));
138: l = 0;
139: PetscCall(PCGASMComputeGlobalSubdomainNumbering_Private(pc, &numbering, &permutation));
140: for (count = 0; count < osm->N; ++count) {
141: /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */
142: if (l < osm->n) {
143: d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */
144: if (numbering[d] == count) {
145: PetscCall(PetscViewerGetSubViewer(viewer, ((PetscObject)osm->ois[d])->comm, &sviewer));
146: PetscCall(PCGASMSubdomainView_Private(pc, d, sviewer));
147: PetscCall(PetscViewerRestoreSubViewer(viewer, ((PetscObject)osm->ois[d])->comm, &sviewer));
148: ++l;
149: }
150: }
151: PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)pc)));
152: }
153: PetscCall(PetscFree2(numbering, permutation));
154: PetscCall(PetscViewerDestroy(&viewer));
155: PetscFunctionReturn(PETSC_SUCCESS);
156: }
158: static PetscErrorCode PCView_GASM(PC pc, PetscViewer viewer)
159: {
160: PC_GASM *osm = (PC_GASM *)pc->data;
161: const char *prefix;
162: PetscMPIInt rank, size;
163: PetscInt bsz;
164: PetscBool iascii, view_subdomains = PETSC_FALSE;
165: PetscViewer sviewer;
166: PetscInt count, l;
167: char overlap[256] = "user-defined overlap";
168: char gsubdomains[256] = "unknown total number of subdomains";
169: char msubdomains[256] = "unknown max number of local subdomains";
170: PetscInt *numbering, *permutation; /* global numbering of locally-supported subdomains and the permutation from the local ordering */
172: PetscFunctionBegin;
173: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
174: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
176: if (osm->overlap >= 0) PetscCall(PetscSNPrintf(overlap, sizeof(overlap), "requested amount of overlap = %" PetscInt_FMT, osm->overlap));
177: if (osm->N != PETSC_DETERMINE) PetscCall(PetscSNPrintf(gsubdomains, sizeof(gsubdomains), "total number of subdomains = %" PetscInt_FMT, osm->N));
178: if (osm->nmax != PETSC_DETERMINE) PetscCall(PetscSNPrintf(msubdomains, sizeof(msubdomains), "max number of local subdomains = %" PetscInt_FMT, osm->nmax));
180: PetscCall(PCGetOptionsPrefix(pc, &prefix));
181: PetscCall(PetscOptionsGetBool(NULL, prefix, "-pc_gasm_view_subdomains", &view_subdomains, NULL));
183: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
184: if (iascii) {
185: /*
186: Make sure the viewer has a name. Otherwise this may cause a deadlock when creating a subcomm viewer:
187: the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed
188: collectively on the comm.
189: */
190: PetscCall(PetscObjectName((PetscObject)viewer));
191: PetscCall(PetscViewerASCIIPrintf(viewer, " Restriction/interpolation type: %s\n", PCGASMTypes[osm->type]));
192: PetscCall(PetscViewerASCIIPrintf(viewer, " %s\n", overlap));
193: PetscCall(PetscViewerASCIIPrintf(viewer, " %s\n", gsubdomains));
194: PetscCall(PetscViewerASCIIPrintf(viewer, " %s\n", msubdomains));
195: PetscCall(PetscViewerASCIIPushSynchronized(viewer));
196: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " [%d|%d] number of locally-supported subdomains = %" PetscInt_FMT "\n", rank, size, osm->n));
197: PetscCall(PetscViewerFlush(viewer));
198: PetscCall(PetscViewerASCIIPopSynchronized(viewer));
199: /* Cannot take advantage of osm->same_subdomain_solvers without a global numbering of subdomains. */
200: PetscCall(PetscViewerASCIIPrintf(viewer, " Subdomain solver info is as follows:\n"));
201: PetscCall(PetscViewerASCIIPushTab(viewer));
202: PetscCall(PetscViewerASCIIPrintf(viewer, " - - - - - - - - - - - - - - - - - -\n"));
203: /* Make sure that everybody waits for the banner to be printed. */
204: PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer)));
205: /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */
206: PetscCall(PCGASMComputeGlobalSubdomainNumbering_Private(pc, &numbering, &permutation));
207: l = 0;
208: for (count = 0; count < osm->N; ++count) {
209: PetscMPIInt srank, ssize;
210: if (l < osm->n) {
211: PetscInt d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */
212: if (numbering[d] == count) {
213: PetscCallMPI(MPI_Comm_size(((PetscObject)osm->ois[d])->comm, &ssize));
214: PetscCallMPI(MPI_Comm_rank(((PetscObject)osm->ois[d])->comm, &srank));
215: PetscCall(PetscViewerGetSubViewer(viewer, ((PetscObject)osm->ois[d])->comm, &sviewer));
216: PetscCall(ISGetLocalSize(osm->ois[d], &bsz));
217: PetscCall(PetscViewerASCIISynchronizedPrintf(sviewer, " [%d|%d] (subcomm [%d|%d]) local subdomain number %" PetscInt_FMT ", local size = %" PetscInt_FMT "\n", rank, size, srank, ssize, d, bsz));
218: PetscCall(PetscViewerFlush(sviewer));
219: if (view_subdomains) PetscCall(PCGASMSubdomainView_Private(pc, d, sviewer));
220: if (!pc->setupcalled) {
221: PetscCall(PetscViewerASCIIPrintf(sviewer, " Solver not set up yet: PCSetUp() not yet called\n"));
222: } else {
223: PetscCall(KSPView(osm->ksp[d], sviewer));
224: }
225: PetscCall(PetscViewerASCIIPrintf(sviewer, " - - - - - - - - - - - - - - - - - -\n"));
226: PetscCall(PetscViewerFlush(sviewer));
227: PetscCall(PetscViewerRestoreSubViewer(viewer, ((PetscObject)osm->ois[d])->comm, &sviewer));
228: ++l;
229: }
230: }
231: PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)pc)));
232: }
233: PetscCall(PetscFree2(numbering, permutation));
234: PetscCall(PetscViewerASCIIPopTab(viewer));
235: PetscCall(PetscViewerFlush(viewer));
236: /* this line is needed to match the extra PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
237: PetscCall(PetscViewerASCIIPopSynchronized(viewer));
238: }
239: PetscFunctionReturn(PETSC_SUCCESS);
240: }
242: PETSC_INTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]);
244: PetscErrorCode PCGASMSetHierarchicalPartitioning(PC pc)
245: {
246: PC_GASM *osm = (PC_GASM *)pc->data;
247: MatPartitioning part;
248: MPI_Comm comm;
249: PetscMPIInt size;
250: PetscInt nlocalsubdomains, fromrows_localsize;
251: IS partitioning, fromrows, isn;
252: Vec outervec;
254: PetscFunctionBegin;
255: PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
256: PetscCallMPI(MPI_Comm_size(comm, &size));
257: /* we do not need a hierarchical partitioning when
258: * the total number of subdomains is consistent with
259: * the number of MPI tasks.
260: * For the following cases, we do not need to use HP
261: * */
262: if (osm->N == PETSC_DETERMINE || osm->N >= size || osm->N == 1) PetscFunctionReturn(PETSC_SUCCESS);
263: PetscCheck(size % osm->N == 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "have to specify the total number of subdomains %" PetscInt_FMT " to be a factor of the number of ranks %d ", osm->N, size);
264: nlocalsubdomains = size / osm->N;
265: osm->n = 1;
266: PetscCall(MatPartitioningCreate(comm, &part));
267: PetscCall(MatPartitioningSetAdjacency(part, pc->pmat));
268: PetscCall(MatPartitioningSetType(part, MATPARTITIONINGHIERARCH));
269: PetscCall(MatPartitioningHierarchicalSetNcoarseparts(part, osm->N));
270: PetscCall(MatPartitioningHierarchicalSetNfineparts(part, nlocalsubdomains));
271: PetscCall(MatPartitioningSetFromOptions(part));
272: /* get new rank owner number of each vertex */
273: PetscCall(MatPartitioningApply(part, &partitioning));
274: PetscCall(ISBuildTwoSided(partitioning, NULL, &fromrows));
275: PetscCall(ISPartitioningToNumbering(partitioning, &isn));
276: PetscCall(ISDestroy(&isn));
277: PetscCall(ISGetLocalSize(fromrows, &fromrows_localsize));
278: PetscCall(MatPartitioningDestroy(&part));
279: PetscCall(MatCreateVecs(pc->pmat, &outervec, NULL));
280: PetscCall(VecCreateMPI(comm, fromrows_localsize, PETSC_DETERMINE, &(osm->pcx)));
281: PetscCall(VecDuplicate(osm->pcx, &(osm->pcy)));
282: PetscCall(VecScatterCreate(osm->pcx, NULL, outervec, fromrows, &(osm->pctoouter)));
283: PetscCall(MatCreateSubMatrix(pc->pmat, fromrows, fromrows, MAT_INITIAL_MATRIX, &(osm->permutationP)));
284: PetscCall(PetscObjectReference((PetscObject)fromrows));
285: osm->permutationIS = fromrows;
286: osm->pcmat = pc->pmat;
287: PetscCall(PetscObjectReference((PetscObject)osm->permutationP));
288: pc->pmat = osm->permutationP;
289: PetscCall(VecDestroy(&outervec));
290: PetscCall(ISDestroy(&fromrows));
291: PetscCall(ISDestroy(&partitioning));
292: osm->n = PETSC_DETERMINE;
293: PetscFunctionReturn(PETSC_SUCCESS);
294: }
296: static PetscErrorCode PCSetUp_GASM(PC pc)
297: {
298: PC_GASM *osm = (PC_GASM *)pc->data;
299: PetscInt i, nInnerIndices, nTotalInnerIndices;
300: PetscMPIInt rank, size;
301: MatReuse scall = MAT_REUSE_MATRIX;
302: KSP ksp;
303: PC subpc;
304: const char *prefix, *pprefix;
305: Vec x, y;
306: PetscInt oni; /* Number of indices in the i-th local outer subdomain. */
307: const PetscInt *oidxi; /* Indices from the i-th subdomain local outer subdomain. */
308: PetscInt on; /* Number of indices in the disjoint union of local outer subdomains. */
309: PetscInt *oidx; /* Indices in the disjoint union of local outer subdomains. */
310: IS gois; /* Disjoint union the global indices of outer subdomains. */
311: IS goid; /* Identity IS of the size of the disjoint union of outer subdomains. */
312: PetscScalar *gxarray, *gyarray;
313: PetscInt gostart; /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- over the disjoint union of outer subdomains. */
314: PetscInt num_subdomains = 0;
315: DM *subdomain_dm = NULL;
316: char **subdomain_names = NULL;
317: PetscInt *numbering;
319: PetscFunctionBegin;
320: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
321: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
322: if (!pc->setupcalled) {
323: /* use a hierarchical partitioning */
324: if (osm->hierarchicalpartitioning) PetscCall(PCGASMSetHierarchicalPartitioning(pc));
325: if (osm->n == PETSC_DETERMINE) {
326: if (osm->N != PETSC_DETERMINE) {
327: /* No local subdomains given, but the desired number of total subdomains is known, so construct them accordingly. */
328: PetscCall(PCGASMCreateSubdomains(pc->pmat, osm->N, &osm->n, &osm->iis));
329: } else if (osm->dm_subdomains && pc->dm) {
330: /* try pc->dm next, if allowed */
331: PetscInt d;
332: IS *inner_subdomain_is, *outer_subdomain_is;
333: PetscCall(DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm));
334: if (num_subdomains) PetscCall(PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is));
335: for (d = 0; d < num_subdomains; ++d) {
336: if (inner_subdomain_is) PetscCall(ISDestroy(&inner_subdomain_is[d]));
337: if (outer_subdomain_is) PetscCall(ISDestroy(&outer_subdomain_is[d]));
338: }
339: PetscCall(PetscFree(inner_subdomain_is));
340: PetscCall(PetscFree(outer_subdomain_is));
341: } else {
342: /* still no subdomains; use one per rank */
343: osm->nmax = osm->n = 1;
344: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
345: osm->N = size;
346: PetscCall(PCGASMCreateLocalSubdomains(pc->pmat, osm->n, &osm->iis));
347: }
348: }
349: if (!osm->iis) {
350: /*
351: osm->n was set in PCGASMSetSubdomains(), but the actual subdomains have not been supplied.
352: We create the requisite number of local inner subdomains and then expand them into
353: out subdomains, if necessary.
354: */
355: PetscCall(PCGASMCreateLocalSubdomains(pc->pmat, osm->n, &osm->iis));
356: }
357: if (!osm->ois) {
358: /*
359: Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap
360: has been requested, copy the inner subdomains over so they can be modified.
361: */
362: PetscCall(PetscMalloc1(osm->n, &osm->ois));
363: for (i = 0; i < osm->n; ++i) {
364: if (osm->overlap > 0 && osm->N > 1) { /* With positive overlap, osm->iis[i] will be modified */
365: PetscCall(ISDuplicate(osm->iis[i], (osm->ois) + i));
366: PetscCall(ISCopy(osm->iis[i], osm->ois[i]));
367: } else {
368: PetscCall(PetscObjectReference((PetscObject)((osm->iis)[i])));
369: osm->ois[i] = osm->iis[i];
370: }
371: }
372: if (osm->overlap > 0 && osm->N > 1) {
373: /* Extend the "overlapping" regions by a number of steps */
374: PetscCall(MatIncreaseOverlapSplit(pc->pmat, osm->n, osm->ois, osm->overlap));
375: }
376: }
378: /* Now the subdomains are defined. Determine their global and max local numbers, if necessary. */
379: if (osm->nmax == PETSC_DETERMINE) {
380: PetscMPIInt inwork, outwork;
381: /* determine global number of subdomains and the max number of local subdomains */
382: inwork = osm->n;
383: PetscCall(MPIU_Allreduce(&inwork, &outwork, 1, MPI_INT, MPI_MAX, PetscObjectComm((PetscObject)pc)));
384: osm->nmax = outwork;
385: }
386: if (osm->N == PETSC_DETERMINE) {
387: /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */
388: PetscCall(PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc), osm->n, (PetscObject *)osm->ois, &osm->N, NULL));
389: }
391: if (osm->sort_indices) {
392: for (i = 0; i < osm->n; i++) {
393: PetscCall(ISSort(osm->ois[i]));
394: PetscCall(ISSort(osm->iis[i]));
395: }
396: }
397: PetscCall(PCGetOptionsPrefix(pc, &prefix));
398: PetscCall(PCGASMPrintSubdomains(pc));
400: /*
401: Merge the ISs, create merged vectors and restrictions.
402: */
403: /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */
404: on = 0;
405: for (i = 0; i < osm->n; i++) {
406: PetscCall(ISGetLocalSize(osm->ois[i], &oni));
407: on += oni;
408: }
409: PetscCall(PetscMalloc1(on, &oidx));
410: on = 0;
411: /* Merge local indices together */
412: for (i = 0; i < osm->n; i++) {
413: PetscCall(ISGetLocalSize(osm->ois[i], &oni));
414: PetscCall(ISGetIndices(osm->ois[i], &oidxi));
415: PetscCall(PetscArraycpy(oidx + on, oidxi, oni));
416: PetscCall(ISRestoreIndices(osm->ois[i], &oidxi));
417: on += oni;
418: }
419: PetscCall(ISCreateGeneral(((PetscObject)(pc))->comm, on, oidx, PETSC_OWN_POINTER, &gois));
420: nTotalInnerIndices = 0;
421: for (i = 0; i < osm->n; i++) {
422: PetscCall(ISGetLocalSize(osm->iis[i], &nInnerIndices));
423: nTotalInnerIndices += nInnerIndices;
424: }
425: PetscCall(VecCreateMPI(((PetscObject)(pc))->comm, nTotalInnerIndices, PETSC_DETERMINE, &x));
426: PetscCall(VecDuplicate(x, &y));
428: PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)pc), on, PETSC_DECIDE, &osm->gx));
429: PetscCall(VecDuplicate(osm->gx, &osm->gy));
430: PetscCall(VecGetOwnershipRange(osm->gx, &gostart, NULL));
431: PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pc), on, gostart, 1, &goid));
432: /* gois might indices not on local */
433: PetscCall(VecScatterCreate(x, gois, osm->gx, goid, &(osm->gorestriction)));
434: PetscCall(PetscMalloc1(osm->n, &numbering));
435: PetscCall(PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc), osm->n, (PetscObject *)osm->ois, NULL, numbering));
436: PetscCall(VecDestroy(&x));
437: PetscCall(ISDestroy(&gois));
439: /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */
440: {
441: PetscInt ini; /* Number of indices the i-th a local inner subdomain. */
442: PetscInt in; /* Number of indices in the disjoint union of local inner subdomains. */
443: PetscInt *iidx; /* Global indices in the merged local inner subdomain. */
444: PetscInt *ioidx; /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
445: IS giis; /* IS for the disjoint union of inner subdomains. */
446: IS giois; /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
447: PetscScalar *array;
448: const PetscInt *indices;
449: PetscInt k;
450: on = 0;
451: for (i = 0; i < osm->n; i++) {
452: PetscCall(ISGetLocalSize(osm->ois[i], &oni));
453: on += oni;
454: }
455: PetscCall(PetscMalloc1(on, &iidx));
456: PetscCall(PetscMalloc1(on, &ioidx));
457: PetscCall(VecGetArray(y, &array));
458: /* set communicator id to determine where overlap is */
459: in = 0;
460: for (i = 0; i < osm->n; i++) {
461: PetscCall(ISGetLocalSize(osm->iis[i], &ini));
462: for (k = 0; k < ini; ++k) array[in + k] = numbering[i];
463: in += ini;
464: }
465: PetscCall(VecRestoreArray(y, &array));
466: PetscCall(VecScatterBegin(osm->gorestriction, y, osm->gy, INSERT_VALUES, SCATTER_FORWARD));
467: PetscCall(VecScatterEnd(osm->gorestriction, y, osm->gy, INSERT_VALUES, SCATTER_FORWARD));
468: PetscCall(VecGetOwnershipRange(osm->gy, &gostart, NULL));
469: PetscCall(VecGetArray(osm->gy, &array));
470: on = 0;
471: in = 0;
472: for (i = 0; i < osm->n; i++) {
473: PetscCall(ISGetLocalSize(osm->ois[i], &oni));
474: PetscCall(ISGetIndices(osm->ois[i], &indices));
475: for (k = 0; k < oni; k++) {
476: /* skip overlapping indices to get inner domain */
477: if (PetscRealPart(array[on + k]) != numbering[i]) continue;
478: iidx[in] = indices[k];
479: ioidx[in++] = gostart + on + k;
480: }
481: PetscCall(ISRestoreIndices(osm->ois[i], &indices));
482: on += oni;
483: }
484: PetscCall(VecRestoreArray(osm->gy, &array));
485: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc), in, iidx, PETSC_OWN_POINTER, &giis));
486: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc), in, ioidx, PETSC_OWN_POINTER, &giois));
487: PetscCall(VecScatterCreate(y, giis, osm->gy, giois, &osm->girestriction));
488: PetscCall(VecDestroy(&y));
489: PetscCall(ISDestroy(&giis));
490: PetscCall(ISDestroy(&giois));
491: }
492: PetscCall(ISDestroy(&goid));
493: PetscCall(PetscFree(numbering));
495: /* Create the subdomain work vectors. */
496: PetscCall(PetscMalloc1(osm->n, &osm->x));
497: PetscCall(PetscMalloc1(osm->n, &osm->y));
498: PetscCall(VecGetArray(osm->gx, &gxarray));
499: PetscCall(VecGetArray(osm->gy, &gyarray));
500: for (i = 0, on = 0; i < osm->n; ++i, on += oni) {
501: PetscInt oNi;
502: PetscCall(ISGetLocalSize(osm->ois[i], &oni));
503: /* on a sub communicator */
504: PetscCall(ISGetSize(osm->ois[i], &oNi));
505: PetscCall(VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm, 1, oni, oNi, gxarray + on, &osm->x[i]));
506: PetscCall(VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm, 1, oni, oNi, gyarray + on, &osm->y[i]));
507: }
508: PetscCall(VecRestoreArray(osm->gx, &gxarray));
509: PetscCall(VecRestoreArray(osm->gy, &gyarray));
510: /* Create the subdomain solvers */
511: PetscCall(PetscMalloc1(osm->n, &osm->ksp));
512: for (i = 0; i < osm->n; i++) {
513: char subprefix[PETSC_MAX_PATH_LEN + 1];
514: PetscCall(KSPCreate(((PetscObject)(osm->ois[i]))->comm, &ksp));
515: PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
516: PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
517: PetscCall(KSPSetType(ksp, KSPPREONLY));
518: PetscCall(KSPGetPC(ksp, &subpc)); /* Why do we need this here? */
519: if (subdomain_dm) {
520: PetscCall(KSPSetDM(ksp, subdomain_dm[i]));
521: PetscCall(DMDestroy(subdomain_dm + i));
522: }
523: PetscCall(PCGetOptionsPrefix(pc, &prefix));
524: PetscCall(KSPSetOptionsPrefix(ksp, prefix));
525: if (subdomain_names && subdomain_names[i]) {
526: PetscCall(PetscSNPrintf(subprefix, PETSC_MAX_PATH_LEN, "sub_%s_", subdomain_names[i]));
527: PetscCall(KSPAppendOptionsPrefix(ksp, subprefix));
528: PetscCall(PetscFree(subdomain_names[i]));
529: }
530: PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
531: osm->ksp[i] = ksp;
532: }
533: PetscCall(PetscFree(subdomain_dm));
534: PetscCall(PetscFree(subdomain_names));
535: scall = MAT_INITIAL_MATRIX;
536: } else { /* if (pc->setupcalled) */
537: /*
538: Destroy the submatrices from the previous iteration
539: */
540: if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
541: PetscCall(MatDestroyMatrices(osm->n, &osm->pmat));
542: scall = MAT_INITIAL_MATRIX;
543: }
544: if (osm->permutationIS) {
545: PetscCall(MatCreateSubMatrix(pc->pmat, osm->permutationIS, osm->permutationIS, scall, &osm->permutationP));
546: PetscCall(PetscObjectReference((PetscObject)osm->permutationP));
547: osm->pcmat = pc->pmat;
548: pc->pmat = osm->permutationP;
549: }
550: }
552: /*
553: Extract the submatrices.
554: */
555: if (size > 1) {
556: PetscCall(MatCreateSubMatricesMPI(pc->pmat, osm->n, osm->ois, osm->ois, scall, &osm->pmat));
557: } else {
558: PetscCall(MatCreateSubMatrices(pc->pmat, osm->n, osm->ois, osm->ois, scall, &osm->pmat));
559: }
560: if (scall == MAT_INITIAL_MATRIX) {
561: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc->pmat, &pprefix));
562: for (i = 0; i < osm->n; i++) PetscCall(PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i], pprefix));
563: }
565: /* Return control to the user so that the submatrices can be modified (e.g., to apply
566: different boundary conditions for the submatrices than for the global problem) */
567: PetscCall(PCModifySubMatrices(pc, osm->n, osm->ois, osm->ois, osm->pmat, pc->modifysubmatricesP));
569: /*
570: Loop over submatrices putting them into local ksps
571: */
572: for (i = 0; i < osm->n; i++) {
573: PetscCall(KSPSetOperators(osm->ksp[i], osm->pmat[i], osm->pmat[i]));
574: PetscCall(KSPGetOptionsPrefix(osm->ksp[i], &prefix));
575: PetscCall(MatSetOptionsPrefix(osm->pmat[i], prefix));
576: if (!pc->setupcalled) PetscCall(KSPSetFromOptions(osm->ksp[i]));
577: }
578: if (osm->pcmat) {
579: PetscCall(MatDestroy(&pc->pmat));
580: pc->pmat = osm->pcmat;
581: osm->pcmat = NULL;
582: }
583: PetscFunctionReturn(PETSC_SUCCESS);
584: }
586: static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc)
587: {
588: PC_GASM *osm = (PC_GASM *)pc->data;
589: PetscInt i;
591: PetscFunctionBegin;
592: for (i = 0; i < osm->n; i++) PetscCall(KSPSetUp(osm->ksp[i]));
593: PetscFunctionReturn(PETSC_SUCCESS);
594: }
596: static PetscErrorCode PCApply_GASM(PC pc, Vec xin, Vec yout)
597: {
598: PC_GASM *osm = (PC_GASM *)pc->data;
599: PetscInt i;
600: Vec x, y;
601: ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE;
603: PetscFunctionBegin;
604: if (osm->pctoouter) {
605: PetscCall(VecScatterBegin(osm->pctoouter, xin, osm->pcx, INSERT_VALUES, SCATTER_REVERSE));
606: PetscCall(VecScatterEnd(osm->pctoouter, xin, osm->pcx, INSERT_VALUES, SCATTER_REVERSE));
607: x = osm->pcx;
608: y = osm->pcy;
609: } else {
610: x = xin;
611: y = yout;
612: }
613: /*
614: support for limiting the restriction or interpolation only to the inner
615: subdomain values (leaving the other values 0).
616: */
617: if (!(osm->type & PC_GASM_RESTRICT)) {
618: /* have to zero the work RHS since scatter may leave some slots empty */
619: PetscCall(VecZeroEntries(osm->gx));
620: PetscCall(VecScatterBegin(osm->girestriction, x, osm->gx, INSERT_VALUES, forward));
621: } else {
622: PetscCall(VecScatterBegin(osm->gorestriction, x, osm->gx, INSERT_VALUES, forward));
623: }
624: PetscCall(VecZeroEntries(osm->gy));
625: if (!(osm->type & PC_GASM_RESTRICT)) {
626: PetscCall(VecScatterEnd(osm->girestriction, x, osm->gx, INSERT_VALUES, forward));
627: } else {
628: PetscCall(VecScatterEnd(osm->gorestriction, x, osm->gx, INSERT_VALUES, forward));
629: }
630: /* do the subdomain solves */
631: for (i = 0; i < osm->n; ++i) {
632: PetscCall(KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]));
633: PetscCall(KSPCheckSolve(osm->ksp[i], pc, osm->y[i]));
634: }
635: /* do we need to zero y? */
636: PetscCall(VecZeroEntries(y));
637: if (!(osm->type & PC_GASM_INTERPOLATE)) {
638: PetscCall(VecScatterBegin(osm->girestriction, osm->gy, y, ADD_VALUES, reverse));
639: PetscCall(VecScatterEnd(osm->girestriction, osm->gy, y, ADD_VALUES, reverse));
640: } else {
641: PetscCall(VecScatterBegin(osm->gorestriction, osm->gy, y, ADD_VALUES, reverse));
642: PetscCall(VecScatterEnd(osm->gorestriction, osm->gy, y, ADD_VALUES, reverse));
643: }
644: if (osm->pctoouter) {
645: PetscCall(VecScatterBegin(osm->pctoouter, y, yout, INSERT_VALUES, SCATTER_FORWARD));
646: PetscCall(VecScatterEnd(osm->pctoouter, y, yout, INSERT_VALUES, SCATTER_FORWARD));
647: }
648: PetscFunctionReturn(PETSC_SUCCESS);
649: }
651: static PetscErrorCode PCMatApply_GASM(PC pc, Mat Xin, Mat Yout)
652: {
653: PC_GASM *osm = (PC_GASM *)pc->data;
654: Mat X, Y, O = NULL, Z, W;
655: Vec x, y;
656: PetscInt i, m, M, N;
657: ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE;
659: PetscFunctionBegin;
660: PetscCheck(osm->n == 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not yet implemented");
661: PetscCall(MatGetSize(Xin, NULL, &N));
662: if (osm->pctoouter) {
663: PetscCall(VecGetLocalSize(osm->pcx, &m));
664: PetscCall(VecGetSize(osm->pcx, &M));
665: PetscCall(MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]), m, PETSC_DECIDE, M, N, NULL, &O));
666: for (i = 0; i < N; ++i) {
667: PetscCall(MatDenseGetColumnVecRead(Xin, i, &x));
668: PetscCall(MatDenseGetColumnVecWrite(O, i, &y));
669: PetscCall(VecScatterBegin(osm->pctoouter, x, y, INSERT_VALUES, SCATTER_REVERSE));
670: PetscCall(VecScatterEnd(osm->pctoouter, x, y, INSERT_VALUES, SCATTER_REVERSE));
671: PetscCall(MatDenseRestoreColumnVecWrite(O, i, &y));
672: PetscCall(MatDenseRestoreColumnVecRead(Xin, i, &x));
673: }
674: X = Y = O;
675: } else {
676: X = Xin;
677: Y = Yout;
678: }
679: /*
680: support for limiting the restriction or interpolation only to the inner
681: subdomain values (leaving the other values 0).
682: */
683: PetscCall(VecGetLocalSize(osm->x[0], &m));
684: PetscCall(VecGetSize(osm->x[0], &M));
685: PetscCall(MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]), m, PETSC_DECIDE, M, N, NULL, &Z));
686: for (i = 0; i < N; ++i) {
687: PetscCall(MatDenseGetColumnVecRead(X, i, &x));
688: PetscCall(MatDenseGetColumnVecWrite(Z, i, &y));
689: if (!(osm->type & PC_GASM_RESTRICT)) {
690: /* have to zero the work RHS since scatter may leave some slots empty */
691: PetscCall(VecZeroEntries(y));
692: PetscCall(VecScatterBegin(osm->girestriction, x, y, INSERT_VALUES, forward));
693: PetscCall(VecScatterEnd(osm->girestriction, x, y, INSERT_VALUES, forward));
694: } else {
695: PetscCall(VecScatterBegin(osm->gorestriction, x, y, INSERT_VALUES, forward));
696: PetscCall(VecScatterEnd(osm->gorestriction, x, y, INSERT_VALUES, forward));
697: }
698: PetscCall(MatDenseRestoreColumnVecWrite(Z, i, &y));
699: PetscCall(MatDenseRestoreColumnVecRead(X, i, &x));
700: }
701: PetscCall(MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]), m, PETSC_DECIDE, M, N, NULL, &W));
702: PetscCall(MatSetOption(Z, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
703: PetscCall(MatAssemblyBegin(Z, MAT_FINAL_ASSEMBLY));
704: PetscCall(MatAssemblyEnd(Z, MAT_FINAL_ASSEMBLY));
705: /* do the subdomain solve */
706: PetscCall(KSPMatSolve(osm->ksp[0], Z, W));
707: PetscCall(KSPCheckSolve(osm->ksp[0], pc, NULL));
708: PetscCall(MatDestroy(&Z));
709: /* do we need to zero y? */
710: PetscCall(MatZeroEntries(Y));
711: for (i = 0; i < N; ++i) {
712: PetscCall(MatDenseGetColumnVecWrite(Y, i, &y));
713: PetscCall(MatDenseGetColumnVecRead(W, i, &x));
714: if (!(osm->type & PC_GASM_INTERPOLATE)) {
715: PetscCall(VecScatterBegin(osm->girestriction, x, y, ADD_VALUES, reverse));
716: PetscCall(VecScatterEnd(osm->girestriction, x, y, ADD_VALUES, reverse));
717: } else {
718: PetscCall(VecScatterBegin(osm->gorestriction, x, y, ADD_VALUES, reverse));
719: PetscCall(VecScatterEnd(osm->gorestriction, x, y, ADD_VALUES, reverse));
720: }
721: PetscCall(MatDenseRestoreColumnVecRead(W, i, &x));
722: if (osm->pctoouter) {
723: PetscCall(MatDenseGetColumnVecWrite(Yout, i, &x));
724: PetscCall(VecScatterBegin(osm->pctoouter, y, x, INSERT_VALUES, SCATTER_FORWARD));
725: PetscCall(VecScatterEnd(osm->pctoouter, y, x, INSERT_VALUES, SCATTER_FORWARD));
726: PetscCall(MatDenseRestoreColumnVecRead(Yout, i, &x));
727: }
728: PetscCall(MatDenseRestoreColumnVecWrite(Y, i, &y));
729: }
730: PetscCall(MatDestroy(&W));
731: PetscCall(MatDestroy(&O));
732: PetscFunctionReturn(PETSC_SUCCESS);
733: }
735: static PetscErrorCode PCApplyTranspose_GASM(PC pc, Vec xin, Vec yout)
736: {
737: PC_GASM *osm = (PC_GASM *)pc->data;
738: PetscInt i;
739: Vec x, y;
740: ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE;
742: PetscFunctionBegin;
743: if (osm->pctoouter) {
744: PetscCall(VecScatterBegin(osm->pctoouter, xin, osm->pcx, INSERT_VALUES, SCATTER_REVERSE));
745: PetscCall(VecScatterEnd(osm->pctoouter, xin, osm->pcx, INSERT_VALUES, SCATTER_REVERSE));
746: x = osm->pcx;
747: y = osm->pcy;
748: } else {
749: x = xin;
750: y = yout;
751: }
752: /*
753: Support for limiting the restriction or interpolation to only local
754: subdomain values (leaving the other values 0).
756: Note: these are reversed from the PCApply_GASM() because we are applying the
757: transpose of the three terms
758: */
759: if (!(osm->type & PC_GASM_INTERPOLATE)) {
760: /* have to zero the work RHS since scatter may leave some slots empty */
761: PetscCall(VecZeroEntries(osm->gx));
762: PetscCall(VecScatterBegin(osm->girestriction, x, osm->gx, INSERT_VALUES, forward));
763: } else {
764: PetscCall(VecScatterBegin(osm->gorestriction, x, osm->gx, INSERT_VALUES, forward));
765: }
766: PetscCall(VecZeroEntries(osm->gy));
767: if (!(osm->type & PC_GASM_INTERPOLATE)) {
768: PetscCall(VecScatterEnd(osm->girestriction, x, osm->gx, INSERT_VALUES, forward));
769: } else {
770: PetscCall(VecScatterEnd(osm->gorestriction, x, osm->gx, INSERT_VALUES, forward));
771: }
772: /* do the local solves */
773: for (i = 0; i < osm->n; ++i) { /* Note that the solves are local, so we can go to osm->n, rather than osm->nmax. */
774: PetscCall(KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i]));
775: PetscCall(KSPCheckSolve(osm->ksp[i], pc, osm->y[i]));
776: }
777: PetscCall(VecZeroEntries(y));
778: if (!(osm->type & PC_GASM_RESTRICT)) {
779: PetscCall(VecScatterBegin(osm->girestriction, osm->gy, y, ADD_VALUES, reverse));
780: PetscCall(VecScatterEnd(osm->girestriction, osm->gy, y, ADD_VALUES, reverse));
781: } else {
782: PetscCall(VecScatterBegin(osm->gorestriction, osm->gy, y, ADD_VALUES, reverse));
783: PetscCall(VecScatterEnd(osm->gorestriction, osm->gy, y, ADD_VALUES, reverse));
784: }
785: if (osm->pctoouter) {
786: PetscCall(VecScatterBegin(osm->pctoouter, y, yout, INSERT_VALUES, SCATTER_FORWARD));
787: PetscCall(VecScatterEnd(osm->pctoouter, y, yout, INSERT_VALUES, SCATTER_FORWARD));
788: }
789: PetscFunctionReturn(PETSC_SUCCESS);
790: }
792: static PetscErrorCode PCReset_GASM(PC pc)
793: {
794: PC_GASM *osm = (PC_GASM *)pc->data;
795: PetscInt i;
797: PetscFunctionBegin;
798: if (osm->ksp) {
799: for (i = 0; i < osm->n; i++) PetscCall(KSPReset(osm->ksp[i]));
800: }
801: if (osm->pmat) {
802: if (osm->n > 0) {
803: PetscMPIInt size;
804: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
805: if (size > 1) {
806: /* osm->pmat is created by MatCreateSubMatricesMPI(), cannot use MatDestroySubMatrices() */
807: PetscCall(MatDestroyMatrices(osm->n, &osm->pmat));
808: } else {
809: PetscCall(MatDestroySubMatrices(osm->n, &osm->pmat));
810: }
811: }
812: }
813: if (osm->x) {
814: for (i = 0; i < osm->n; i++) {
815: PetscCall(VecDestroy(&osm->x[i]));
816: PetscCall(VecDestroy(&osm->y[i]));
817: }
818: }
819: PetscCall(VecDestroy(&osm->gx));
820: PetscCall(VecDestroy(&osm->gy));
822: PetscCall(VecScatterDestroy(&osm->gorestriction));
823: PetscCall(VecScatterDestroy(&osm->girestriction));
824: if (!osm->user_subdomains) {
825: PetscCall(PCGASMDestroySubdomains(osm->n, &osm->ois, &osm->iis));
826: osm->N = PETSC_DETERMINE;
827: osm->nmax = PETSC_DETERMINE;
828: }
829: if (osm->pctoouter) PetscCall(VecScatterDestroy(&(osm->pctoouter)));
830: if (osm->permutationIS) PetscCall(ISDestroy(&(osm->permutationIS)));
831: if (osm->pcx) PetscCall(VecDestroy(&(osm->pcx)));
832: if (osm->pcy) PetscCall(VecDestroy(&(osm->pcy)));
833: if (osm->permutationP) PetscCall(MatDestroy(&(osm->permutationP)));
834: if (osm->pcmat) PetscCall(MatDestroy(&osm->pcmat));
835: PetscFunctionReturn(PETSC_SUCCESS);
836: }
838: static PetscErrorCode PCDestroy_GASM(PC pc)
839: {
840: PC_GASM *osm = (PC_GASM *)pc->data;
841: PetscInt i;
843: PetscFunctionBegin;
844: PetscCall(PCReset_GASM(pc));
845: /* PCReset will not destroy subdomains, if user_subdomains is true. */
846: PetscCall(PCGASMDestroySubdomains(osm->n, &osm->ois, &osm->iis));
847: if (osm->ksp) {
848: for (i = 0; i < osm->n; i++) PetscCall(KSPDestroy(&osm->ksp[i]));
849: PetscCall(PetscFree(osm->ksp));
850: }
851: PetscCall(PetscFree(osm->x));
852: PetscCall(PetscFree(osm->y));
853: PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMSetSubdomains_C", NULL));
854: PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMSetOverlap_C", NULL));
855: PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMSetType_C", NULL));
856: PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMSetSortIndices_C", NULL));
857: PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMGetSubKSP_C", NULL));
858: PetscCall(PetscFree(pc->data));
859: PetscFunctionReturn(PETSC_SUCCESS);
860: }
862: static PetscErrorCode PCSetFromOptions_GASM(PC pc, PetscOptionItems *PetscOptionsObject)
863: {
864: PC_GASM *osm = (PC_GASM *)pc->data;
865: PetscInt blocks, ovl;
866: PetscBool flg;
867: PCGASMType gasmtype;
869: PetscFunctionBegin;
870: PetscOptionsHeadBegin(PetscOptionsObject, "Generalized additive Schwarz options");
871: PetscCall(PetscOptionsBool("-pc_gasm_use_dm_subdomains", "If subdomains aren't set, use DMCreateDomainDecomposition() to define subdomains.", "PCGASMSetUseDMSubdomains", osm->dm_subdomains, &osm->dm_subdomains, &flg));
872: PetscCall(PetscOptionsInt("-pc_gasm_total_subdomains", "Total number of subdomains across communicator", "PCGASMSetTotalSubdomains", osm->N, &blocks, &flg));
873: if (flg) PetscCall(PCGASMSetTotalSubdomains(pc, blocks));
874: PetscCall(PetscOptionsInt("-pc_gasm_overlap", "Number of overlapping degrees of freedom", "PCGASMSetOverlap", osm->overlap, &ovl, &flg));
875: if (flg) {
876: PetscCall(PCGASMSetOverlap(pc, ovl));
877: osm->dm_subdomains = PETSC_FALSE;
878: }
879: flg = PETSC_FALSE;
880: PetscCall(PetscOptionsEnum("-pc_gasm_type", "Type of restriction/extension", "PCGASMSetType", PCGASMTypes, (PetscEnum)osm->type, (PetscEnum *)&gasmtype, &flg));
881: if (flg) PetscCall(PCGASMSetType(pc, gasmtype));
882: PetscCall(PetscOptionsBool("-pc_gasm_use_hierachical_partitioning", "use hierarchical partitioning", NULL, osm->hierarchicalpartitioning, &osm->hierarchicalpartitioning, &flg));
883: PetscOptionsHeadEnd();
884: PetscFunctionReturn(PETSC_SUCCESS);
885: }
887: /*@
888: PCGASMSetTotalSubdomains - sets the total number of subdomains to use across the communicator for `PCGASM`
890: Logically Collective
892: Input Parameters:
893: + pc - the preconditioner
894: - N - total number of subdomains
896: Level: beginner
898: .seealso: `PCGASM`, `PCGASMSetSubdomains()`, `PCGASMSetOverlap()`
899: `PCGASMCreateSubdomains2D()`
900: @*/
901: PetscErrorCode PCGASMSetTotalSubdomains(PC pc, PetscInt N)
902: {
903: PC_GASM *osm = (PC_GASM *)pc->data;
904: PetscMPIInt size, rank;
906: PetscFunctionBegin;
907: PetscCheck(N >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Total number of subdomains must be 1 or more, got N = %" PetscInt_FMT, N);
908: PetscCheck(!pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCGASMSetTotalSubdomains() should be called before calling PCSetUp().");
910: PetscCall(PCGASMDestroySubdomains(osm->n, &osm->iis, &osm->ois));
911: osm->ois = osm->iis = NULL;
913: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
914: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
915: osm->N = N;
916: osm->n = PETSC_DETERMINE;
917: osm->nmax = PETSC_DETERMINE;
918: osm->dm_subdomains = PETSC_FALSE;
919: PetscFunctionReturn(PETSC_SUCCESS);
920: }
922: static PetscErrorCode PCGASMSetSubdomains_GASM(PC pc, PetscInt n, IS iis[], IS ois[])
923: {
924: PC_GASM *osm = (PC_GASM *)pc->data;
925: PetscInt i;
927: PetscFunctionBegin;
928: PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Each MPI rank must have 1 or more subdomains, got n = %" PetscInt_FMT, n);
929: PetscCheck(!pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCGASMSetSubdomains() should be called before calling PCSetUp().");
931: PetscCall(PCGASMDestroySubdomains(osm->n, &osm->iis, &osm->ois));
932: osm->iis = osm->ois = NULL;
933: osm->n = n;
934: osm->N = PETSC_DETERMINE;
935: osm->nmax = PETSC_DETERMINE;
936: if (ois) {
937: PetscCall(PetscMalloc1(n, &osm->ois));
938: for (i = 0; i < n; i++) {
939: PetscCall(PetscObjectReference((PetscObject)ois[i]));
940: osm->ois[i] = ois[i];
941: }
942: /*
943: Since the user set the outer subdomains, even if nontrivial overlap was requested via PCGASMSetOverlap(),
944: it will be ignored. To avoid confusion later on (e.g., when viewing the PC), the overlap size is set to -1.
945: */
946: osm->overlap = -1;
947: /* inner subdomains must be provided */
948: PetscCheck(iis, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "inner indices have to be provided ");
949: } /* end if */
950: if (iis) {
951: PetscCall(PetscMalloc1(n, &osm->iis));
952: for (i = 0; i < n; i++) {
953: PetscCall(PetscObjectReference((PetscObject)iis[i]));
954: osm->iis[i] = iis[i];
955: }
956: if (!ois) {
957: osm->ois = NULL;
958: /* if user does not provide outer indices, we will create the corresponding outer indices using osm->overlap =1 in PCSetUp_GASM */
959: }
960: }
961: if (PetscDefined(USE_DEBUG)) {
962: PetscInt j, rstart, rend, *covered, lsize;
963: const PetscInt *indices;
964: /* check if the inner indices cover and only cover the local portion of the preconditioning matrix */
965: PetscCall(MatGetOwnershipRange(pc->pmat, &rstart, &rend));
966: PetscCall(PetscCalloc1(rend - rstart, &covered));
967: /* check if the current MPI rank owns indices from others */
968: for (i = 0; i < n; i++) {
969: PetscCall(ISGetIndices(osm->iis[i], &indices));
970: PetscCall(ISGetLocalSize(osm->iis[i], &lsize));
971: for (j = 0; j < lsize; j++) {
972: PetscCheck(indices[j] >= rstart && indices[j] < rend, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "inner subdomains can not own an index %" PetscInt_FMT " from other ranks", indices[j]);
973: PetscCheck(covered[indices[j] - rstart] != 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "inner subdomains can not have an overlapping index %" PetscInt_FMT " ", indices[j]);
974: covered[indices[j] - rstart] = 1;
975: }
976: PetscCall(ISRestoreIndices(osm->iis[i], &indices));
977: }
978: /* check if we miss any indices */
979: for (i = rstart; i < rend; i++) PetscCheck(covered[i - rstart], PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "local entity %" PetscInt_FMT " was not covered by inner subdomains", i);
980: PetscCall(PetscFree(covered));
981: }
982: if (iis) osm->user_subdomains = PETSC_TRUE;
983: PetscFunctionReturn(PETSC_SUCCESS);
984: }
986: static PetscErrorCode PCGASMSetOverlap_GASM(PC pc, PetscInt ovl)
987: {
988: PC_GASM *osm = (PC_GASM *)pc->data;
990: PetscFunctionBegin;
991: PetscCheck(ovl >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap value requested");
992: PetscCheck(!pc->setupcalled || ovl == osm->overlap, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCGASMSetOverlap() should be called before PCSetUp().");
993: if (!pc->setupcalled) osm->overlap = ovl;
994: PetscFunctionReturn(PETSC_SUCCESS);
995: }
997: static PetscErrorCode PCGASMSetType_GASM(PC pc, PCGASMType type)
998: {
999: PC_GASM *osm = (PC_GASM *)pc->data;
1001: PetscFunctionBegin;
1002: osm->type = type;
1003: osm->type_set = PETSC_TRUE;
1004: PetscFunctionReturn(PETSC_SUCCESS);
1005: }
1007: static PetscErrorCode PCGASMSetSortIndices_GASM(PC pc, PetscBool doSort)
1008: {
1009: PC_GASM *osm = (PC_GASM *)pc->data;
1011: PetscFunctionBegin;
1012: osm->sort_indices = doSort;
1013: PetscFunctionReturn(PETSC_SUCCESS);
1014: }
1016: /*
1017: FIXME: This routine might need to be modified now that multiple ranks per subdomain are allowed.
1018: In particular, it would upset the global subdomain number calculation.
1019: */
1020: static PetscErrorCode PCGASMGetSubKSP_GASM(PC pc, PetscInt *n, PetscInt *first, KSP **ksp)
1021: {
1022: PC_GASM *osm = (PC_GASM *)pc->data;
1024: PetscFunctionBegin;
1025: PetscCheck(osm->n >= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Need to call PCSetUp() on PC (or KSPSetUp() on the outer KSP object) before calling here");
1027: if (n) *n = osm->n;
1028: if (first) {
1029: PetscCallMPI(MPI_Scan(&osm->n, first, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)pc)));
1030: *first -= osm->n;
1031: }
1032: if (ksp) {
1033: /* Assume that local solves are now different; not necessarily
1034: true, though! This flag is used only for PCView_GASM() */
1035: *ksp = osm->ksp;
1036: osm->same_subdomain_solvers = PETSC_FALSE;
1037: }
1038: PetscFunctionReturn(PETSC_SUCCESS);
1039: } /* PCGASMGetSubKSP_GASM() */
1041: /*@C
1042: PCGASMSetSubdomains - Sets the subdomains for this MPI rank
1043: for the additive Schwarz preconditioner with multiple MPI ranks per subdomain, `PCGASM`
1045: Collective
1047: Input Parameters:
1048: + pc - the preconditioner object
1049: . n - the number of subdomains for this MPI rank
1050: . iis - the index sets that define the inner subdomains (or `NULL` for PETSc to determine subdomains)
1051: - ois - the index sets that define the outer subdomains (or `NULL` to use the same as `iis`, or to construct by expanding `iis` by
1052: the requested overlap)
1054: Level: advanced
1056: Notes:
1057: The `IS` indices use the parallel, global numbering of the vector entries.
1059: Inner subdomains are those where the correction is applied.
1061: Outer subdomains are those where the residual necessary to obtain the
1062: corrections is obtained (see `PCGASMType` for the use of inner/outer subdomains).
1064: Both inner and outer subdomains can extend over several MPI ranks.
1065: This rank's portion of a subdomain is known as a local subdomain.
1067: Inner subdomains can not overlap with each other, do not have any entities from remote ranks,
1068: and have to cover the entire local subdomain owned by the current rank. The index sets on each
1069: rank should be ordered such that the ith local subdomain is connected to the ith remote subdomain
1070: on another MPI rank.
1072: By default the `PGASM` preconditioner uses 1 (local) subdomain per MPI rank.
1074: The `iis` and `ois` arrays may be freed after this call using `PCGASMDestroySubdomains()`
1076: .seealso: `PCGASM`, `PCGASMSetOverlap()`, `PCGASMGetSubKSP()`, `PCGASMDestroySubdomains()`,
1077: `PCGASMCreateSubdomains2D()`, `PCGASMGetSubdomains()`
1078: @*/
1079: PetscErrorCode PCGASMSetSubdomains(PC pc, PetscInt n, IS iis[], IS ois[])
1080: {
1081: PC_GASM *osm = (PC_GASM *)pc->data;
1083: PetscFunctionBegin;
1085: PetscTryMethod(pc, "PCGASMSetSubdomains_C", (PC, PetscInt, IS[], IS[]), (pc, n, iis, ois));
1086: osm->dm_subdomains = PETSC_FALSE;
1087: PetscFunctionReturn(PETSC_SUCCESS);
1088: }
1090: /*@
1091: PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the
1092: additive Schwarz preconditioner `PCGASM`. Either all or no MPI ranks in the
1093: pc communicator must call this routine.
1095: Logically Collective
1097: Input Parameters:
1098: + pc - the preconditioner context
1099: - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 0)
1101: Options Database Key:
1102: . -pc_gasm_overlap <overlap> - Sets overlap
1104: Level: intermediate
1106: Notes:
1107: By default the `PCGASM` preconditioner uses 1 subdomain per rank. To use
1108: multiple subdomain per perocessor or "straddling" subdomains that intersect
1109: multiple ranks use `PCGASMSetSubdomains()` (or option `-pc_gasm_total_subdomains` <n>).
1111: The overlap defaults to 0, so if one desires that no additional
1112: overlap be computed beyond what may have been set with a call to
1113: `PCGASMSetSubdomains()`, then `ovl` must be set to be 0. In particular, if one does
1114: not explicitly set the subdomains in application code, then all overlap would be computed
1115: internally by PETSc, and using an overlap of 0 would result in an `PCGASM`
1116: variant that is equivalent to the block Jacobi preconditioner.
1118: One can define initial index sets with any overlap via
1119: `PCGASMSetSubdomains()`; the routine `PCGASMSetOverlap()` merely allows
1120: PETSc to extend that overlap further, if desired.
1122: .seealso: `PCGASM`, `PCGASMSetSubdomains()`, `PCGASMGetSubKSP()`,
1123: `PCGASMCreateSubdomains2D()`, `PCGASMGetSubdomains()`
1124: @*/
1125: PetscErrorCode PCGASMSetOverlap(PC pc, PetscInt ovl)
1126: {
1127: PC_GASM *osm = (PC_GASM *)pc->data;
1129: PetscFunctionBegin;
1132: PetscTryMethod(pc, "PCGASMSetOverlap_C", (PC, PetscInt), (pc, ovl));
1133: osm->dm_subdomains = PETSC_FALSE;
1134: PetscFunctionReturn(PETSC_SUCCESS);
1135: }
1137: /*@
1138: PCGASMSetType - Sets the type of restriction and interpolation used
1139: for local problems in the `PCGASM` additive Schwarz method.
1141: Logically Collective
1143: Input Parameters:
1144: + pc - the preconditioner context
1145: - type - variant of `PCGASM`, one of
1146: .vb
1147: `PC_GASM_BASIC` - full interpolation and restriction
1148: `PC_GASM_RESTRICT` - full restriction, local MPI rank interpolation
1149: `PC_GASM_INTERPOLATE` - full interpolation, local MPI rank restriction
1150: `PC_GASM_NONE` - local MPI rank restriction and interpolation
1151: .ve
1153: Options Database Key:
1154: . -pc_gasm_type [basic,restrict,interpolate,none] - Sets `PCGASM` type
1156: Level: intermediate
1158: .seealso: `PCGASM`, `PCGASMSetSubdomains()`, `PCGASMGetSubKSP()`,
1159: `PCGASMCreateSubdomains2D()`, `PCASM`, `PCASMSetType()`
1160: @*/
1161: PetscErrorCode PCGASMSetType(PC pc, PCGASMType type)
1162: {
1163: PetscFunctionBegin;
1166: PetscTryMethod(pc, "PCGASMSetType_C", (PC, PCGASMType), (pc, type));
1167: PetscFunctionReturn(PETSC_SUCCESS);
1168: }
1170: /*@
1171: PCGASMSetSortIndices - Determines whether subdomain indices are sorted.
1173: Logically Collective
1175: Input Parameters:
1176: + pc - the preconditioner context
1177: - doSort - sort the subdomain indices
1179: Level: intermediate
1181: .seealso: `PCGASM`, `PCGASMSetSubdomains()`, `PCGASMGetSubKSP()`,
1182: `PCGASMCreateSubdomains2D()`
1183: @*/
1184: PetscErrorCode PCGASMSetSortIndices(PC pc, PetscBool doSort)
1185: {
1186: PetscFunctionBegin;
1189: PetscTryMethod(pc, "PCGASMSetSortIndices_C", (PC, PetscBool), (pc, doSort));
1190: PetscFunctionReturn(PETSC_SUCCESS);
1191: }
1193: /*@C
1194: PCGASMGetSubKSP - Gets the local `KSP` contexts for all subdomains on
1195: this MPI rank.
1197: Collective iff first_local is requested
1199: Input Parameter:
1200: . pc - the preconditioner context
1202: Output Parameters:
1203: + n_local - the number of blocks on this MPI rank or `NULL`
1204: . first_local - the global number of the first block on this rank or `NULL`,
1205: all ranks must request or all must pass `NULL`
1206: - ksp - the array of `KSP` contexts
1208: Level: advanced
1210: Note:
1211: After `PCGASMGetSubKSP()` the array of `KSP`es is not to be freed
1213: Currently for some matrix implementations only 1 block per MPI process
1214: is supported.
1216: You must call `KSPSetUp()` before calling `PCGASMGetSubKSP()`.
1218: .seealso: `PCGASM`, `PCGASMSetSubdomains()`, `PCGASMSetOverlap()`,
1219: `PCGASMCreateSubdomains2D()`,
1220: @*/
1221: PetscErrorCode PCGASMGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[])
1222: {
1223: PetscFunctionBegin;
1225: PetscUseMethod(pc, "PCGASMGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp));
1226: PetscFunctionReturn(PETSC_SUCCESS);
1227: }
1229: /*MC
1230: PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
1231: its own `KSP` object on a subset of MPI ranks
1233: Options Database Keys:
1234: + -pc_gasm_total_subdomains <n> - Sets total number of local subdomains to be distributed among MPI ranks
1235: . -pc_gasm_view_subdomains - activates the printing of subdomain indices in `PCView()`, -ksp_view or -snes_view
1236: . -pc_gasm_print_subdomains - activates the printing of subdomain indices in `PCSetUp()`
1237: . -pc_gasm_overlap <ovl> - Sets overlap by which to (automatically) extend local subdomains
1238: - -pc_gasm_type [basic,restrict,interpolate,none] - Sets `PCGASMType`
1240: Level: beginner
1242: Notes:
1243: To set options on the solvers for each block append `-sub_` to all the `KSP`, and `PC`
1244: options database keys. For example, `-sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly`
1246: To set the options on the solvers separate for each block call `PCGASMGetSubKSP()`
1247: and set the options directly on the resulting `KSP` object (you can access its `PC`
1248: with `KSPGetPC()`)
1250: References:
1251: + * - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions
1252: Courant Institute, New York University Technical report
1253: - * - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
1254: Cambridge University Press.
1256: .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCASM`, `PCGASMType`, `PCGASMSetType()`,
1257: `PCBJACOBI`, `PCGASMGetSubKSP()`, `PCGASMSetSubdomains()`,
1258: `PCSetModifySubMatrices()`, `PCGASMSetOverlap()`, `PCGASMSetType()`
1259: M*/
1261: PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc)
1262: {
1263: PC_GASM *osm;
1265: PetscFunctionBegin;
1266: PetscCall(PetscNew(&osm));
1268: osm->N = PETSC_DETERMINE;
1269: osm->n = PETSC_DECIDE;
1270: osm->nmax = PETSC_DETERMINE;
1271: osm->overlap = 0;
1272: osm->ksp = NULL;
1273: osm->gorestriction = NULL;
1274: osm->girestriction = NULL;
1275: osm->pctoouter = NULL;
1276: osm->gx = NULL;
1277: osm->gy = NULL;
1278: osm->x = NULL;
1279: osm->y = NULL;
1280: osm->pcx = NULL;
1281: osm->pcy = NULL;
1282: osm->permutationIS = NULL;
1283: osm->permutationP = NULL;
1284: osm->pcmat = NULL;
1285: osm->ois = NULL;
1286: osm->iis = NULL;
1287: osm->pmat = NULL;
1288: osm->type = PC_GASM_RESTRICT;
1289: osm->same_subdomain_solvers = PETSC_TRUE;
1290: osm->sort_indices = PETSC_TRUE;
1291: osm->dm_subdomains = PETSC_FALSE;
1292: osm->hierarchicalpartitioning = PETSC_FALSE;
1294: pc->data = (void *)osm;
1295: pc->ops->apply = PCApply_GASM;
1296: pc->ops->matapply = PCMatApply_GASM;
1297: pc->ops->applytranspose = PCApplyTranspose_GASM;
1298: pc->ops->setup = PCSetUp_GASM;
1299: pc->ops->reset = PCReset_GASM;
1300: pc->ops->destroy = PCDestroy_GASM;
1301: pc->ops->setfromoptions = PCSetFromOptions_GASM;
1302: pc->ops->setuponblocks = PCSetUpOnBlocks_GASM;
1303: pc->ops->view = PCView_GASM;
1304: pc->ops->applyrichardson = NULL;
1306: PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMSetSubdomains_C", PCGASMSetSubdomains_GASM));
1307: PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMSetOverlap_C", PCGASMSetOverlap_GASM));
1308: PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMSetType_C", PCGASMSetType_GASM));
1309: PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMSetSortIndices_C", PCGASMSetSortIndices_GASM));
1310: PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMGetSubKSP_C", PCGASMGetSubKSP_GASM));
1311: PetscFunctionReturn(PETSC_SUCCESS);
1312: }
1314: PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[])
1315: {
1316: MatPartitioning mpart;
1317: const char *prefix;
1318: PetscInt i, j, rstart, rend, bs;
1319: PetscBool hasop, isbaij = PETSC_FALSE, foundpart = PETSC_FALSE;
1320: Mat Ad = NULL, adj;
1321: IS ispart, isnumb, *is;
1323: PetscFunctionBegin;
1324: PetscCheck(nloc >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "number of local subdomains must > 0, got nloc = %" PetscInt_FMT, nloc);
1326: /* Get prefix, row distribution, and block size */
1327: PetscCall(MatGetOptionsPrefix(A, &prefix));
1328: PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
1329: PetscCall(MatGetBlockSize(A, &bs));
1330: PetscCheck(rstart / bs * bs == rstart && rend / bs * bs == rend, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "bad row distribution [%" PetscInt_FMT ",%" PetscInt_FMT ") for matrix block size %" PetscInt_FMT, rstart, rend, bs);
1332: /* Get diagonal block from matrix if possible */
1333: PetscCall(MatHasOperation(A, MATOP_GET_DIAGONAL_BLOCK, &hasop));
1334: if (hasop) PetscCall(MatGetDiagonalBlock(A, &Ad));
1335: if (Ad) {
1336: PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad, MATSEQBAIJ, &isbaij));
1337: if (!isbaij) PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad, MATSEQSBAIJ, &isbaij));
1338: }
1339: if (Ad && nloc > 1) {
1340: PetscBool match, done;
1341: /* Try to setup a good matrix partitioning if available */
1342: PetscCall(MatPartitioningCreate(PETSC_COMM_SELF, &mpart));
1343: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart, prefix));
1344: PetscCall(MatPartitioningSetFromOptions(mpart));
1345: PetscCall(PetscObjectTypeCompare((PetscObject)mpart, MATPARTITIONINGCURRENT, &match));
1346: if (!match) PetscCall(PetscObjectTypeCompare((PetscObject)mpart, MATPARTITIONINGSQUARE, &match));
1347: if (!match) { /* assume a "good" partitioner is available */
1348: PetscInt na;
1349: const PetscInt *ia, *ja;
1350: PetscCall(MatGetRowIJ(Ad, 0, PETSC_TRUE, isbaij, &na, &ia, &ja, &done));
1351: if (done) {
1352: /* Build adjacency matrix by hand. Unfortunately a call to
1353: MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1354: remove the block-aij structure and we cannot expect
1355: MatPartitioning to split vertices as we need */
1356: PetscInt i, j, len, nnz, cnt, *iia = NULL, *jja = NULL;
1357: const PetscInt *row;
1358: nnz = 0;
1359: for (i = 0; i < na; i++) { /* count number of nonzeros */
1360: len = ia[i + 1] - ia[i];
1361: row = ja + ia[i];
1362: for (j = 0; j < len; j++) {
1363: if (row[j] == i) { /* don't count diagonal */
1364: len--;
1365: break;
1366: }
1367: }
1368: nnz += len;
1369: }
1370: PetscCall(PetscMalloc1(na + 1, &iia));
1371: PetscCall(PetscMalloc1(nnz, &jja));
1372: nnz = 0;
1373: iia[0] = 0;
1374: for (i = 0; i < na; i++) { /* fill adjacency */
1375: cnt = 0;
1376: len = ia[i + 1] - ia[i];
1377: row = ja + ia[i];
1378: for (j = 0; j < len; j++) {
1379: if (row[j] != i) jja[nnz + cnt++] = row[j]; /* if not diagonal */
1380: }
1381: nnz += cnt;
1382: iia[i + 1] = nnz;
1383: }
1384: /* Partitioning of the adjacency matrix */
1385: PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, na, na, iia, jja, NULL, &adj));
1386: PetscCall(MatPartitioningSetAdjacency(mpart, adj));
1387: PetscCall(MatPartitioningSetNParts(mpart, nloc));
1388: PetscCall(MatPartitioningApply(mpart, &ispart));
1389: PetscCall(ISPartitioningToNumbering(ispart, &isnumb));
1390: PetscCall(MatDestroy(&adj));
1391: foundpart = PETSC_TRUE;
1392: }
1393: PetscCall(MatRestoreRowIJ(Ad, 0, PETSC_TRUE, isbaij, &na, &ia, &ja, &done));
1394: }
1395: PetscCall(MatPartitioningDestroy(&mpart));
1396: }
1397: PetscCall(PetscMalloc1(nloc, &is));
1398: if (!foundpart) {
1399: /* Partitioning by contiguous chunks of rows */
1401: PetscInt mbs = (rend - rstart) / bs;
1402: PetscInt start = rstart;
1403: for (i = 0; i < nloc; i++) {
1404: PetscInt count = (mbs / nloc + ((mbs % nloc) > i)) * bs;
1405: PetscCall(ISCreateStride(PETSC_COMM_SELF, count, start, 1, &is[i]));
1406: start += count;
1407: }
1409: } else {
1410: /* Partitioning by adjacency of diagonal block */
1412: const PetscInt *numbering;
1413: PetscInt *count, nidx, *indices, *newidx, start = 0;
1414: /* Get node count in each partition */
1415: PetscCall(PetscMalloc1(nloc, &count));
1416: PetscCall(ISPartitioningCount(ispart, nloc, count));
1417: if (isbaij && bs > 1) { /* adjust for the block-aij case */
1418: for (i = 0; i < nloc; i++) count[i] *= bs;
1419: }
1420: /* Build indices from node numbering */
1421: PetscCall(ISGetLocalSize(isnumb, &nidx));
1422: PetscCall(PetscMalloc1(nidx, &indices));
1423: for (i = 0; i < nidx; i++) indices[i] = i; /* needs to be initialized */
1424: PetscCall(ISGetIndices(isnumb, &numbering));
1425: PetscCall(PetscSortIntWithPermutation(nidx, numbering, indices));
1426: PetscCall(ISRestoreIndices(isnumb, &numbering));
1427: if (isbaij && bs > 1) { /* adjust for the block-aij case */
1428: PetscCall(PetscMalloc1(nidx * bs, &newidx));
1429: for (i = 0; i < nidx; i++) {
1430: for (j = 0; j < bs; j++) newidx[i * bs + j] = indices[i] * bs + j;
1431: }
1432: PetscCall(PetscFree(indices));
1433: nidx *= bs;
1434: indices = newidx;
1435: }
1436: /* Shift to get global indices */
1437: for (i = 0; i < nidx; i++) indices[i] += rstart;
1439: /* Build the index sets for each block */
1440: for (i = 0; i < nloc; i++) {
1441: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, count[i], &indices[start], PETSC_COPY_VALUES, &is[i]));
1442: PetscCall(ISSort(is[i]));
1443: start += count[i];
1444: }
1446: PetscCall(PetscFree(count));
1447: PetscCall(PetscFree(indices));
1448: PetscCall(ISDestroy(&isnumb));
1449: PetscCall(ISDestroy(&ispart));
1450: }
1451: *iis = is;
1452: PetscFunctionReturn(PETSC_SUCCESS);
1453: }
1455: PETSC_INTERN PetscErrorCode PCGASMCreateStraddlingSubdomains(Mat A, PetscInt N, PetscInt *n, IS *iis[])
1456: {
1457: PetscFunctionBegin;
1458: PetscCall(MatSubdomainsCreateCoalesce(A, N, n, iis));
1459: PetscFunctionReturn(PETSC_SUCCESS);
1460: }
1462: /*@C
1463: PCGASMCreateSubdomains - Creates `n` index sets defining `n` nonoverlapping subdomains on this MPI process for the `PCGASM` additive
1464: Schwarz preconditioner for a any problem based on its matrix.
1466: Collective
1468: Input Parameters:
1469: + A - The global matrix operator
1470: - N - the number of global subdomains requested
1472: Output Parameters:
1473: + n - the number of subdomains created on this MPI rank
1474: - iis - the array of index sets defining the local inner subdomains (on which the correction is applied)
1476: Level: advanced
1478: Notes:
1479: When `N` >= A's communicator size, each subdomain is local -- contained within a single MPI process.
1480: When `N` < size, the subdomains are 'straddling' (rank boundaries) and are no longer local.
1481: The resulting subdomains can be use in `PCGASMSetSubdomains`(pc,n,iss,`NULL`). The overlapping
1482: outer subdomains will be automatically generated from these according to the requested amount of
1483: overlap; this is currently supported only with local subdomains.
1485: Use `PCGASMDestroySubdomains()` to free the array and the list of index sets.
1487: .seealso: `PCGASM`, `PCGASMSetSubdomains()`, `PCGASMDestroySubdomains()`
1488: @*/
1489: PetscErrorCode PCGASMCreateSubdomains(Mat A, PetscInt N, PetscInt *n, IS *iis[])
1490: {
1491: PetscMPIInt size;
1493: PetscFunctionBegin;
1497: PetscCheck(N >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of subdomains must be > 0, N = %" PetscInt_FMT, N);
1498: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1499: if (N >= size) {
1500: *n = N / size + (N % size);
1501: PetscCall(PCGASMCreateLocalSubdomains(A, *n, iis));
1502: } else {
1503: PetscCall(PCGASMCreateStraddlingSubdomains(A, N, n, iis));
1504: }
1505: PetscFunctionReturn(PETSC_SUCCESS);
1506: }
1508: /*@C
1509: PCGASMDestroySubdomains - Destroys the index sets created with
1510: `PCGASMCreateSubdomains()` or `PCGASMCreateSubdomains2D()`. Should be
1511: called after setting subdomains with `PCGASMSetSubdomains()`.
1513: Collective
1515: Input Parameters:
1516: + n - the number of index sets
1517: . iis - the array of inner subdomains
1518: - ois - the array of outer subdomains, can be `NULL`
1520: Level: intermediate
1522: Note:
1523: This is a convenience subroutine that walks each list,
1524: destroys each `IS` on the list, and then frees the list. At the end the
1525: list pointers are set to `NULL`.
1527: .seealso: `PCGASM`, `PCGASMCreateSubdomains()`, `PCGASMSetSubdomains()`, `PCGASMSetSubdomains()`
1528: @*/
1529: PetscErrorCode PCGASMDestroySubdomains(PetscInt n, IS **iis, IS **ois)
1530: {
1531: PetscInt i;
1533: PetscFunctionBegin;
1534: if (n <= 0) PetscFunctionReturn(PETSC_SUCCESS);
1535: if (ois) {
1537: if (*ois) {
1539: for (i = 0; i < n; i++) PetscCall(ISDestroy(&(*ois)[i]));
1540: PetscCall(PetscFree((*ois)));
1541: }
1542: }
1543: if (iis) {
1545: if (*iis) {
1547: for (i = 0; i < n; i++) PetscCall(ISDestroy(&(*iis)[i]));
1548: PetscCall(PetscFree((*iis)));
1549: }
1550: }
1551: PetscFunctionReturn(PETSC_SUCCESS);
1552: }
1554: #define PCGASMLocalSubdomainBounds2D(M, N, xleft, ylow, xright, yhigh, first, last, xleft_loc, ylow_loc, xright_loc, yhigh_loc, n) \
1555: { \
1556: PetscInt first_row = first / M, last_row = last / M + 1; \
1557: /* \
1558: Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners \
1559: of the bounding box of the intersection of the subdomain with the local ownership range (local \
1560: subdomain). \
1561: Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows \
1562: of the intersection. \
1563: */ \
1564: /* ylow_loc is the grid row containing the first element of the local sumbdomain */ \
1565: *ylow_loc = PetscMax(first_row, ylow); \
1566: /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \
1567: *xleft_loc = *ylow_loc == first_row ? PetscMax(first % M, xleft) : xleft; \
1568: /* yhigh_loc is the grid row above the last local subdomain element */ \
1569: *yhigh_loc = PetscMin(last_row, yhigh); \
1570: /* xright is the offset of the end of the local subdomain within its grid row (might actually be outside the local subdomain) */ \
1571: *xright_loc = *yhigh_loc == last_row ? PetscMin(xright, last % M) : xright; \
1572: /* Now compute the size of the local subdomain n. */ \
1573: *n = 0; \
1574: if (*ylow_loc < *yhigh_loc) { \
1575: PetscInt width = xright - xleft; \
1576: *n += width * (*yhigh_loc - *ylow_loc - 1); \
1577: *n += PetscMin(PetscMax(*xright_loc - xleft, 0), width); \
1578: *n -= PetscMin(PetscMax(*xleft_loc - xleft, 0), width); \
1579: } \
1580: }
1582: /*@C
1583: PCGASMCreateSubdomains2D - Creates the index sets for the `PCGASM` overlapping Schwarz
1584: preconditioner for a two-dimensional problem on a regular grid.
1586: Collective
1588: Input Parameters:
1589: + pc - the preconditioner context
1590: . M - the global number of grid points in the x direction
1591: . N - the global number of grid points in the y direction
1592: . Mdomains - the global number of subdomains in the x direction
1593: . Ndomains - the global number of subdomains in the y direction
1594: . dof - degrees of freedom per node
1595: - overlap - overlap in mesh lines
1597: Output Parameters:
1598: + Nsub - the number of local subdomains created
1599: . iis - array of index sets defining inner (nonoverlapping) subdomains
1600: - ois - array of index sets defining outer (overlapping, if overlap > 0) subdomains
1602: Level: advanced
1604: Note:
1605: Use `PCGASMDestroySubdomains()` to free the index sets and the arrays
1607: Fortran Note:
1608: The `IS` must be declared as an array of length long enough to hold `Nsub` entries
1610: .seealso: `PCGASM`, `PCGASMSetSubdomains()`, `PCGASMGetSubKSP()`, `PCGASMSetOverlap()`, `PCASMCreateSubdomains2D()`,
1611: `PCGASMDestroySubdomains()`
1612: @*/
1613: PetscErrorCode PCGASMCreateSubdomains2D(PC pc, PetscInt M, PetscInt N, PetscInt Mdomains, PetscInt Ndomains, PetscInt dof, PetscInt overlap, PetscInt *nsub, IS **iis, IS **ois)
1614: {
1615: PetscMPIInt size, rank;
1616: PetscInt i, j;
1617: PetscInt maxheight, maxwidth;
1618: PetscInt xstart, xleft, xright, xleft_loc, xright_loc;
1619: PetscInt ystart, ylow, yhigh, ylow_loc, yhigh_loc;
1620: PetscInt x[2][2], y[2][2], n[2];
1621: PetscInt first, last;
1622: PetscInt nidx, *idx;
1623: PetscInt ii, jj, s, q, d;
1624: PetscInt k, kk;
1625: PetscMPIInt color;
1626: MPI_Comm comm, subcomm;
1627: IS **xis = NULL, **is = ois, **is_local = iis;
1629: PetscFunctionBegin;
1630: PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
1631: PetscCallMPI(MPI_Comm_size(comm, &size));
1632: PetscCallMPI(MPI_Comm_rank(comm, &rank));
1633: PetscCall(MatGetOwnershipRange(pc->pmat, &first, &last));
1634: PetscCheck((first % dof) == 0 && (last % dof) == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE,
1635: "Matrix row partitioning unsuitable for domain decomposition: local row range (%" PetscInt_FMT ",%" PetscInt_FMT ") "
1636: "does not respect the number of degrees of freedom per grid point %" PetscInt_FMT,
1637: first, last, dof);
1639: /* Determine the number of domains with nonzero intersections with the local ownership range. */
1640: s = 0;
1641: ystart = 0;
1642: for (j = 0; j < Ndomains; ++j) {
1643: maxheight = N / Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1644: PetscCheck(maxheight >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many %" PetscInt_FMT " subdomains in the vertical direction for mesh height %" PetscInt_FMT, Ndomains, N);
1645: /* Vertical domain limits with an overlap. */
1646: ylow = PetscMax(ystart - overlap, 0);
1647: yhigh = PetscMin(ystart + maxheight + overlap, N);
1648: xstart = 0;
1649: for (i = 0; i < Mdomains; ++i) {
1650: maxwidth = M / Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1651: PetscCheck(maxwidth >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many %" PetscInt_FMT " subdomains in the horizontal direction for mesh width %" PetscInt_FMT, Mdomains, M);
1652: /* Horizontal domain limits with an overlap. */
1653: xleft = PetscMax(xstart - overlap, 0);
1654: xright = PetscMin(xstart + maxwidth + overlap, M);
1655: /*
1656: Determine whether this subdomain intersects this rank's ownership range of pc->pmat.
1657: */
1658: PCGASMLocalSubdomainBounds2D(M, N, xleft, ylow, xright, yhigh, first, last, (&xleft_loc), (&ylow_loc), (&xright_loc), (&yhigh_loc), (&nidx));
1659: if (nidx) ++s;
1660: xstart += maxwidth;
1661: } /* for (i = 0; i < Mdomains; ++i) */
1662: ystart += maxheight;
1663: } /* for (j = 0; j < Ndomains; ++j) */
1665: /* Now we can allocate the necessary number of ISs. */
1666: *nsub = s;
1667: PetscCall(PetscMalloc1(*nsub, is));
1668: PetscCall(PetscMalloc1(*nsub, is_local));
1669: s = 0;
1670: ystart = 0;
1671: for (j = 0; j < Ndomains; ++j) {
1672: maxheight = N / Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1673: PetscCheck(maxheight >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many %" PetscInt_FMT " subdomains in the vertical direction for mesh height %" PetscInt_FMT, Ndomains, N);
1674: /* Vertical domain limits with an overlap. */
1675: y[0][0] = PetscMax(ystart - overlap, 0);
1676: y[0][1] = PetscMin(ystart + maxheight + overlap, N);
1677: /* Vertical domain limits without an overlap. */
1678: y[1][0] = ystart;
1679: y[1][1] = PetscMin(ystart + maxheight, N);
1680: xstart = 0;
1681: for (i = 0; i < Mdomains; ++i) {
1682: maxwidth = M / Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1683: PetscCheck(maxwidth >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many %" PetscInt_FMT " subdomains in the horizontal direction for mesh width %" PetscInt_FMT, Mdomains, M);
1684: /* Horizontal domain limits with an overlap. */
1685: x[0][0] = PetscMax(xstart - overlap, 0);
1686: x[0][1] = PetscMin(xstart + maxwidth + overlap, M);
1687: /* Horizontal domain limits without an overlap. */
1688: x[1][0] = xstart;
1689: x[1][1] = PetscMin(xstart + maxwidth, M);
1690: /*
1691: Determine whether this domain intersects this rank's ownership range of pc->pmat.
1692: Do this twice: first for the domains with overlaps, and once without.
1693: During the first pass create the subcommunicators, and use them on the second pass as well.
1694: */
1695: for (q = 0; q < 2; ++q) {
1696: PetscBool split = PETSC_FALSE;
1697: /*
1698: domain limits, (xleft, xright) and (ylow, yheigh) are adjusted
1699: according to whether the domain with an overlap or without is considered.
1700: */
1701: xleft = x[q][0];
1702: xright = x[q][1];
1703: ylow = y[q][0];
1704: yhigh = y[q][1];
1705: PCGASMLocalSubdomainBounds2D(M, N, xleft, ylow, xright, yhigh, first, last, (&xleft_loc), (&ylow_loc), (&xright_loc), (&yhigh_loc), (&nidx));
1706: nidx *= dof;
1707: n[q] = nidx;
1708: /*
1709: Based on the counted number of indices in the local domain *with an overlap*,
1710: construct a subcommunicator of all the MPI ranks supporting this domain.
1711: Observe that a domain with an overlap might have nontrivial local support,
1712: while the domain without an overlap might not. Hence, the decision to participate
1713: in the subcommunicator must be based on the domain with an overlap.
1714: */
1715: if (q == 0) {
1716: if (nidx) color = 1;
1717: else color = MPI_UNDEFINED;
1718: PetscCallMPI(MPI_Comm_split(comm, color, rank, &subcomm));
1719: split = PETSC_TRUE;
1720: }
1721: /*
1722: Proceed only if the number of local indices *with an overlap* is nonzero.
1723: */
1724: if (n[0]) {
1725: if (q == 0) xis = is;
1726: if (q == 1) {
1727: /*
1728: The IS for the no-overlap subdomain shares a communicator with the overlapping domain.
1729: Moreover, if the overlap is zero, the two ISs are identical.
1730: */
1731: if (overlap == 0) {
1732: (*is_local)[s] = (*is)[s];
1733: PetscCall(PetscObjectReference((PetscObject)(*is)[s]));
1734: continue;
1735: } else {
1736: xis = is_local;
1737: subcomm = ((PetscObject)(*is)[s])->comm;
1738: }
1739: } /* if (q == 1) */
1740: idx = NULL;
1741: PetscCall(PetscMalloc1(nidx, &idx));
1742: if (nidx) {
1743: k = 0;
1744: for (jj = ylow_loc; jj < yhigh_loc; ++jj) {
1745: PetscInt x0 = (jj == ylow_loc) ? xleft_loc : xleft;
1746: PetscInt x1 = (jj == yhigh_loc - 1) ? xright_loc : xright;
1747: kk = dof * (M * jj + x0);
1748: for (ii = x0; ii < x1; ++ii) {
1749: for (d = 0; d < dof; ++d) idx[k++] = kk++;
1750: }
1751: }
1752: }
1753: PetscCall(ISCreateGeneral(subcomm, nidx, idx, PETSC_OWN_POINTER, (*xis) + s));
1754: if (split) PetscCallMPI(MPI_Comm_free(&subcomm));
1755: } /* if (n[0]) */
1756: } /* for (q = 0; q < 2; ++q) */
1757: if (n[0]) ++s;
1758: xstart += maxwidth;
1759: } /* for (i = 0; i < Mdomains; ++i) */
1760: ystart += maxheight;
1761: } /* for (j = 0; j < Ndomains; ++j) */
1762: PetscFunctionReturn(PETSC_SUCCESS);
1763: }
1765: /*@C
1766: PCGASMGetSubdomains - Gets the subdomains supported on this MPI rank
1767: for the `PCGASM` additive Schwarz preconditioner.
1769: Not Collective
1771: Input Parameter:
1772: . pc - the preconditioner context
1774: Output Parameters:
1775: + n - the number of subdomains for this MPI rank (default value = 1)
1776: . iis - the index sets that define the inner subdomains (without overlap) supported on this rank (can be `NULL`)
1777: - ois - the index sets that define the outer subdomains (with overlap) supported on this rank (can be `NULL`)
1779: Level: advanced
1781: Notes:
1782: The user is responsible for destroying the `IS`s and freeing the returned arrays, this can be done with
1783: `PCGASMDestroySubdomains()`
1785: The `IS` numbering is in the parallel, global numbering of the vector.
1787: .seealso: `PCGASM`, `PCGASMSetOverlap()`, `PCGASMGetSubKSP()`, `PCGASMCreateSubdomains2D()`,
1788: `PCGASMSetSubdomains()`, `PCGASMGetSubmatrices()`, `PCGASMDestroySubdomains()`
1789: @*/
1790: PetscErrorCode PCGASMGetSubdomains(PC pc, PetscInt *n, IS *iis[], IS *ois[])
1791: {
1792: PC_GASM *osm;
1793: PetscBool match;
1794: PetscInt i;
1796: PetscFunctionBegin;
1798: PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCGASM, &match));
1799: PetscCheck(match, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1800: osm = (PC_GASM *)pc->data;
1801: if (n) *n = osm->n;
1802: if (iis) PetscCall(PetscMalloc1(osm->n, iis));
1803: if (ois) PetscCall(PetscMalloc1(osm->n, ois));
1804: if (iis || ois) {
1805: for (i = 0; i < osm->n; ++i) {
1806: if (iis) (*iis)[i] = osm->iis[i];
1807: if (ois) (*ois)[i] = osm->ois[i];
1808: }
1809: }
1810: PetscFunctionReturn(PETSC_SUCCESS);
1811: }
1813: /*@C
1814: PCGASMGetSubmatrices - Gets the local submatrices (for this MPI rank
1815: only) for the `PCGASM` additive Schwarz preconditioner.
1817: Not Collective
1819: Input Parameter:
1820: . pc - the preconditioner context
1822: Output Parameters:
1823: + n - the number of matrices for this MPI rank (default value = 1)
1824: - mat - the matrices
1826: Level: advanced
1828: Note:
1829: Matrices returned by this routine have the same communicators as the index sets (`IS`)
1830: used to define subdomains in `PCGASMSetSubdomains()`
1832: .seealso: `PCGASM`, `PCGASMSetOverlap()`, `PCGASMGetSubKSP()`,
1833: `PCGASMCreateSubdomains2D()`, `PCGASMSetSubdomains()`, `PCGASMGetSubdomains()`
1834: @*/
1835: PetscErrorCode PCGASMGetSubmatrices(PC pc, PetscInt *n, Mat *mat[])
1836: {
1837: PC_GASM *osm;
1838: PetscBool match;
1840: PetscFunctionBegin;
1844: PetscCheck(pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call after KSPSetUp() or PCSetUp().");
1845: PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCGASM, &match));
1846: PetscCheck(match, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1847: osm = (PC_GASM *)pc->data;
1848: if (n) *n = osm->n;
1849: if (mat) *mat = osm->pmat;
1850: PetscFunctionReturn(PETSC_SUCCESS);
1851: }
1853: /*@
1854: PCGASMSetUseDMSubdomains - Indicates whether to use `DMCreateDomainDecomposition()` to define the subdomains, whenever possible for `PCGASM`
1856: Logically Collective
1858: Input Parameters:
1859: + pc - the preconditioner
1860: - flg - boolean indicating whether to use subdomains defined by the `DM`
1862: Options Database Key:
1863: . -pc_gasm_dm_subdomains -pc_gasm_overlap -pc_gasm_total_subdomains
1865: Level: intermediate
1867: Note:
1868: `PCGASMSetSubdomains()`, `PCGASMSetTotalSubdomains()` or `PCGASMSetOverlap()` take precedence over `PCGASMSetUseDMSubdomains()`,
1869: so setting `PCGASMSetSubdomains()` with nontrivial subdomain ISs or any of `PCGASMSetTotalSubdomains()` and `PCGASMSetOverlap()`
1870: automatically turns the latter off.
1872: .seealso: `PCGASM`, `PCGASMGetUseDMSubdomains()`, `PCGASMSetSubdomains()`, `PCGASMSetOverlap()`
1873: `PCGASMCreateSubdomains2D()`
1874: @*/
1875: PetscErrorCode PCGASMSetUseDMSubdomains(PC pc, PetscBool flg)
1876: {
1877: PC_GASM *osm = (PC_GASM *)pc->data;
1878: PetscBool match;
1880: PetscFunctionBegin;
1883: PetscCheck(!pc->setupcalled, ((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Not for a setup PC.");
1884: PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCGASM, &match));
1885: if (match) {
1886: if (!osm->user_subdomains && osm->N == PETSC_DETERMINE && osm->overlap < 0) osm->dm_subdomains = flg;
1887: }
1888: PetscFunctionReturn(PETSC_SUCCESS);
1889: }
1891: /*@
1892: PCGASMGetUseDMSubdomains - Returns flag indicating whether to use `DMCreateDomainDecomposition()` to define the subdomains, whenever possible with `PCGASM`
1894: Not Collective
1896: Input Parameter:
1897: . pc - the preconditioner
1899: Output Parameter:
1900: . flg - boolean indicating whether to use subdomains defined by the `DM`
1902: Level: intermediate
1904: .seealso: `PCGASM`, `PCGASMSetUseDMSubdomains()`, `PCGASMSetOverlap()`
1905: `PCGASMCreateSubdomains2D()`
1906: @*/
1907: PetscErrorCode PCGASMGetUseDMSubdomains(PC pc, PetscBool *flg)
1908: {
1909: PC_GASM *osm = (PC_GASM *)pc->data;
1910: PetscBool match;
1912: PetscFunctionBegin;
1915: PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCGASM, &match));
1916: if (match) {
1917: if (flg) *flg = osm->dm_subdomains;
1918: }
1919: PetscFunctionReturn(PETSC_SUCCESS);
1920: }