xref: /petsc/src/ksp/pc/impls/gasm/gasm.c (revision fe8fb07485f33b84ebb214f2984459f87e22ea81)
1b1a0a8a3SJed Brown /*
2f746d493SDmitry Karpeev   This file defines an "generalized" additive Schwarz preconditioner for any Mat implementation.
3f1580f4eSBarry Smith   In this version each MPI rank may intersect multiple subdomains and any subdomain may
4f1580f4eSBarry Smith   intersect multiple MPI ranks.  Intersections of subdomains with MPI ranks are called *local
56a4f0f73SDmitry Karpeev   subdomains*.
6b1a0a8a3SJed Brown 
78f3b4b4dSDmitry Karpeev        N    - total number of distinct global subdomains  (set explicitly in PCGASMSetTotalSubdomains() or implicitly PCGASMSetSubdomains() and then calculated in PCSetUp_GASM())
8f1580f4eSBarry Smith        n    - actual number of local subdomains on this rank (set in PCGASMSetSubdomains() or calculated in PCGASMSetTotalSubdomains())
9f1580f4eSBarry Smith        nmax - maximum number of local subdomains per rank    (calculated in PCSetUp_GASM())
10b1a0a8a3SJed Brown */
11af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
121e25c274SJed Brown #include <petscdm.h>
13b1a0a8a3SJed Brown 
14b1a0a8a3SJed Brown typedef struct {
15f746d493SDmitry Karpeev   PetscInt   N, n, nmax;
16b1a0a8a3SJed Brown   PetscInt   overlap;                /* overlap requested by user */
178f3b4b4dSDmitry Karpeev   PCGASMType type;                   /* use reduced interpolation, restriction or both */
188f3b4b4dSDmitry Karpeev   PetscBool  type_set;               /* if user set this value (so won't change it for symmetric problems) */
198f3b4b4dSDmitry Karpeev   PetscBool  same_subdomain_solvers; /* flag indicating whether all local solvers are same */
208f3b4b4dSDmitry Karpeev   PetscBool  sort_indices;           /* flag to sort subdomain indices */
218f3b4b4dSDmitry Karpeev   PetscBool  user_subdomains;        /* whether the user set explicit subdomain index sets -- keep them on PCReset() */
228f3b4b4dSDmitry Karpeev   PetscBool  dm_subdomains;          /* whether DM is allowed to define subdomains */
23ea91fabdSFande Kong   PetscBool  hierarchicalpartitioning;
248f3b4b4dSDmitry Karpeev   IS        *ois;           /* index sets that define the outer (conceptually, overlapping) subdomains */
258f3b4b4dSDmitry Karpeev   IS        *iis;           /* index sets that define the inner (conceptually, nonoverlapping) subdomains */
268f3b4b4dSDmitry Karpeev   KSP       *ksp;           /* linear solvers for each subdomain */
278f3b4b4dSDmitry Karpeev   Mat       *pmat;          /* subdomain block matrices */
28f746d493SDmitry Karpeev   Vec        gx, gy;        /* Merged work vectors */
29f746d493SDmitry Karpeev   Vec       *x, *y;         /* Split work vectors; storage aliases pieces of storage of the above merged vectors. */
306a4f0f73SDmitry Karpeev   VecScatter gorestriction; /* merged restriction to disjoint union of outer subdomains */
316a4f0f73SDmitry Karpeev   VecScatter girestriction; /* merged restriction to disjoint union of inner subdomains */
32ea91fabdSFande Kong   VecScatter pctoouter;
33ea91fabdSFande Kong   IS         permutationIS;
34ea91fabdSFande Kong   Mat        permutationP;
35ea91fabdSFande Kong   Mat        pcmat;
36ea91fabdSFande Kong   Vec        pcx, pcy;
37f746d493SDmitry Karpeev } PC_GASM;
38b1a0a8a3SJed Brown 
39d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGASMComputeGlobalSubdomainNumbering_Private(PC pc, PetscInt **numbering, PetscInt **permutation)
40d71ae5a4SJacob Faibussowitsch {
418f3b4b4dSDmitry Karpeev   PC_GASM *osm = (PC_GASM *)pc->data;
428f3b4b4dSDmitry Karpeev   PetscInt i;
438f3b4b4dSDmitry Karpeev 
448f3b4b4dSDmitry Karpeev   PetscFunctionBegin;
458f3b4b4dSDmitry Karpeev   /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */
469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(osm->n, numbering, osm->n, permutation));
479566063dSJacob Faibussowitsch   PetscCall(PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc), osm->n, (PetscObject *)osm->iis, NULL, *numbering));
488f3b4b4dSDmitry Karpeev   for (i = 0; i < osm->n; ++i) (*permutation)[i] = i;
499566063dSJacob Faibussowitsch   PetscCall(PetscSortIntWithPermutation(osm->n, *numbering, *permutation));
503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
518f3b4b4dSDmitry Karpeev }
528f3b4b4dSDmitry Karpeev 
53d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGASMSubdomainView_Private(PC pc, PetscInt i, PetscViewer viewer)
54d71ae5a4SJacob Faibussowitsch {
55af538404SDmitry Karpeev   PC_GASM        *osm = (PC_GASM *)pc->data;
5606b43e7eSDmitry Karpeev   PetscInt        j, nidx;
57af538404SDmitry Karpeev   const PetscInt *idx;
5806b43e7eSDmitry Karpeev   PetscViewer     sviewer;
59af538404SDmitry Karpeev   char           *cidx;
60af538404SDmitry Karpeev 
61af538404SDmitry Karpeev   PetscFunctionBegin;
62*fe8fb074SBarry Smith   PetscCheck(i >= -1 && 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*fe8fb074SBarry Smith 
644bde246dSDmitry Karpeev   /* Inner subdomains. */
654bde246dSDmitry Karpeev   /*
664bde246dSDmitry Karpeev    No more than 15 characters per index plus a space.
674bde246dSDmitry Karpeev    PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
684bde246dSDmitry Karpeev    in case nidx == 0. That will take care of the space for the trailing '\0' as well.
694bde246dSDmitry Karpeev    For nidx == 0, the whole string 16 '\0'.
704bde246dSDmitry Karpeev    */
71*fe8fb074SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "Inner subdomain:\n"));
72*fe8fb074SBarry Smith   PetscCall(PetscViewerFlush(viewer));
73*fe8fb074SBarry Smith   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
74*fe8fb074SBarry Smith   if (i > -1) {
75*fe8fb074SBarry Smith     PetscCall(ISGetLocalSize(osm->iis[i], &nidx));
76*fe8fb074SBarry Smith     PetscCall(PetscMalloc1(16 * (nidx + 1) + 1, &cidx));
77*fe8fb074SBarry Smith     PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16 * (nidx + 1) + 1, &sviewer));
789566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(osm->iis[i], &idx));
7948a46eb9SPierre Jolivet     for (j = 0; j < nidx; ++j) PetscCall(PetscViewerStringSPrintf(sviewer, "%" PetscInt_FMT " ", idx[j]));
809566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(osm->iis[i], &idx));
819566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&sviewer));
829566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx));
83*fe8fb074SBarry Smith     PetscCall(PetscFree(cidx));
84*fe8fb074SBarry Smith   }
859566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
869566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
879566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
889566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
89*fe8fb074SBarry Smith 
904bde246dSDmitry Karpeev   /* Outer subdomains. */
91eec7e3faSDmitry Karpeev   /*
92eec7e3faSDmitry Karpeev    No more than 15 characters per index plus a space.
93eec7e3faSDmitry Karpeev    PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
94eec7e3faSDmitry Karpeev    in case nidx == 0. That will take care of the space for the trailing '\0' as well.
95eec7e3faSDmitry Karpeev    For nidx == 0, the whole string 16 '\0'.
96eec7e3faSDmitry Karpeev    */
97*fe8fb074SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "Outer subdomain:\n"));
98*fe8fb074SBarry Smith   PetscCall(PetscViewerFlush(viewer));
99*fe8fb074SBarry Smith   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
100*fe8fb074SBarry Smith   if (i > -1) {
101*fe8fb074SBarry Smith     PetscCall(ISGetLocalSize(osm->ois[i], &nidx));
102*fe8fb074SBarry Smith     PetscCall(PetscMalloc1(16 * (nidx + 1) + 1, &cidx));
103*fe8fb074SBarry Smith     PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16 * (nidx + 1) + 1, &sviewer));
1049566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(osm->ois[i], &idx));
10548a46eb9SPierre Jolivet     for (j = 0; j < nidx; ++j) PetscCall(PetscViewerStringSPrintf(sviewer, "%" PetscInt_FMT " ", idx[j]));
1069566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&sviewer));
1079566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(osm->ois[i], &idx));
1089566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx));
109*fe8fb074SBarry Smith     PetscCall(PetscFree(cidx));
110*fe8fb074SBarry Smith   }
1119566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
1129566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1139566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1149566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
1153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11606b43e7eSDmitry Karpeev }
11706b43e7eSDmitry Karpeev 
118d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGASMPrintSubdomains(PC pc)
119d71ae5a4SJacob Faibussowitsch {
12006b43e7eSDmitry Karpeev   PC_GASM    *osm = (PC_GASM *)pc->data;
12106b43e7eSDmitry Karpeev   const char *prefix;
12206b43e7eSDmitry Karpeev   char        fname[PETSC_MAX_PATH_LEN + 1];
1238f3b4b4dSDmitry Karpeev   PetscInt    l, d, count;
1249140fbc5SPierre Jolivet   PetscBool   found;
125*fe8fb074SBarry Smith   PetscViewer viewer;
1268f3b4b4dSDmitry Karpeev   PetscInt   *numbering, *permutation; /* global numbering of locally-supported subdomains and the permutation from the local ordering */
12706b43e7eSDmitry Karpeev 
12806b43e7eSDmitry Karpeev   PetscFunctionBegin;
1299566063dSJacob Faibussowitsch   PetscCall(PCGetOptionsPrefix(pc, &prefix));
1309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, prefix, "-pc_gasm_print_subdomains", &found));
1313ba16761SJacob Faibussowitsch   if (!found) PetscFunctionReturn(PETSC_SUCCESS);
1329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, prefix, "-pc_gasm_print_subdomains", fname, sizeof(fname), &found));
133c6a7a370SJeremy L Thompson   if (!found) PetscCall(PetscStrncpy(fname, "stdout", sizeof(fname)));
1349566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc), fname, &viewer));
1354bde246dSDmitry Karpeev   /*
1364bde246dSDmitry Karpeev    Make sure the viewer has a name. Otherwise this may cause a deadlock or other weird errors when creating a subcomm viewer:
1374bde246dSDmitry Karpeev    the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed collectively on the comm.
1384bde246dSDmitry Karpeev   */
1399566063dSJacob Faibussowitsch   PetscCall(PetscObjectName((PetscObject)viewer));
1404bde246dSDmitry Karpeev   l = 0;
1419566063dSJacob Faibussowitsch   PetscCall(PCGASMComputeGlobalSubdomainNumbering_Private(pc, &numbering, &permutation));
1428f3b4b4dSDmitry Karpeev   for (count = 0; count < osm->N; ++count) {
1434bde246dSDmitry Karpeev     /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */
1444bde246dSDmitry Karpeev     if (l < osm->n) {
1454bde246dSDmitry Karpeev       d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */
146*fe8fb074SBarry Smith       if (numbering[d] == count) l++;
147*fe8fb074SBarry Smith       else d = -1;
148*fe8fb074SBarry Smith     } else d = -1;
149*fe8fb074SBarry Smith     PetscCall(PCGASMSubdomainView_Private(pc, d, viewer));
1504bde246dSDmitry Karpeev   }
1519566063dSJacob Faibussowitsch   PetscCall(PetscFree2(numbering, permutation));
1529566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
1533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
154af538404SDmitry Karpeev }
155af538404SDmitry Karpeev 
156d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_GASM(PC pc, PetscViewer viewer)
157d71ae5a4SJacob Faibussowitsch {
158f746d493SDmitry Karpeev   PC_GASM    *osm = (PC_GASM *)pc->data;
159af538404SDmitry Karpeev   const char *prefix;
160af538404SDmitry Karpeev   PetscMPIInt rank, size;
1618f3b4b4dSDmitry Karpeev   PetscInt    bsz;
16206b43e7eSDmitry Karpeev   PetscBool   iascii, view_subdomains = PETSC_FALSE;
163b1a0a8a3SJed Brown   PetscViewer sviewer;
1648f3b4b4dSDmitry Karpeev   PetscInt    count, l;
1656a4f0f73SDmitry Karpeev   char        overlap[256]     = "user-defined overlap";
1666a4f0f73SDmitry Karpeev   char        gsubdomains[256] = "unknown total number of subdomains";
16706b43e7eSDmitry Karpeev   char        msubdomains[256] = "unknown max number of local subdomains";
1688f3b4b4dSDmitry Karpeev   PetscInt   *numbering, *permutation; /* global numbering of locally-supported subdomains and the permutation from the local ordering */
16911aeaf0aSBarry Smith 
170b1a0a8a3SJed Brown   PetscFunctionBegin;
1719566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
1729566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
17306b43e7eSDmitry Karpeev 
17448a46eb9SPierre Jolivet   if (osm->overlap >= 0) PetscCall(PetscSNPrintf(overlap, sizeof(overlap), "requested amount of overlap = %" PetscInt_FMT, osm->overlap));
17548a46eb9SPierre Jolivet   if (osm->N != PETSC_DETERMINE) PetscCall(PetscSNPrintf(gsubdomains, sizeof(gsubdomains), "total number of subdomains = %" PetscInt_FMT, osm->N));
17648a46eb9SPierre Jolivet   if (osm->nmax != PETSC_DETERMINE) PetscCall(PetscSNPrintf(msubdomains, sizeof(msubdomains), "max number of local subdomains = %" PetscInt_FMT, osm->nmax));
17706b43e7eSDmitry Karpeev 
1789566063dSJacob Faibussowitsch   PetscCall(PCGetOptionsPrefix(pc, &prefix));
1799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, prefix, "-pc_gasm_view_subdomains", &view_subdomains, NULL));
18006b43e7eSDmitry Karpeev 
1819566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
182b1a0a8a3SJed Brown   if (iascii) {
18306b43e7eSDmitry Karpeev     /*
18406b43e7eSDmitry Karpeev      Make sure the viewer has a name. Otherwise this may cause a deadlock when creating a subcomm viewer:
18506b43e7eSDmitry Karpeev      the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed
18606b43e7eSDmitry Karpeev      collectively on the comm.
18706b43e7eSDmitry Karpeev      */
1889566063dSJacob Faibussowitsch     PetscCall(PetscObjectName((PetscObject)viewer));
1899566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Restriction/interpolation type: %s\n", PCGASMTypes[osm->type]));
1909566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  %s\n", overlap));
1919566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  %s\n", gsubdomains));
1929566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  %s\n", msubdomains));
1939566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
19463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  [%d|%d] number of locally-supported subdomains = %" PetscInt_FMT "\n", rank, size, osm->n));
1959566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
1969566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1976a4f0f73SDmitry Karpeev     /* Cannot take advantage of osm->same_subdomain_solvers without a global numbering of subdomains. */
1989566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Subdomain solver info is as follows:\n"));
1999566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
2009566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  - - - - - - - - - - - - - - - - - -\n"));
20106b43e7eSDmitry Karpeev     /* Make sure that everybody waits for the banner to be printed. */
2029566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer)));
2034bde246dSDmitry Karpeev     /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */
2049566063dSJacob Faibussowitsch     PetscCall(PCGASMComputeGlobalSubdomainNumbering_Private(pc, &numbering, &permutation));
20506b43e7eSDmitry Karpeev     l = 0;
2068f3b4b4dSDmitry Karpeev     for (count = 0; count < osm->N; ++count) {
20706b43e7eSDmitry Karpeev       PetscMPIInt srank, ssize;
20806b43e7eSDmitry Karpeev       if (l < osm->n) {
20906b43e7eSDmitry Karpeev         PetscInt d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */
21006b43e7eSDmitry Karpeev         if (numbering[d] == count) {
2119566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Comm_size(((PetscObject)osm->ois[d])->comm, &ssize));
2129566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Comm_rank(((PetscObject)osm->ois[d])->comm, &srank));
2139566063dSJacob Faibussowitsch           PetscCall(PetscViewerGetSubViewer(viewer, ((PetscObject)osm->ois[d])->comm, &sviewer));
2149566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(osm->ois[d], &bsz));
21563a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(sviewer, "  [%d|%d] (subcomm [%d|%d]) local subdomain number %" PetscInt_FMT ", local size = %" PetscInt_FMT "\n", rank, size, srank, ssize, d, bsz));
2169566063dSJacob Faibussowitsch           PetscCall(PetscViewerFlush(sviewer));
2171baa6e33SBarry Smith           if (view_subdomains) PetscCall(PCGASMSubdomainView_Private(pc, d, sviewer));
2186a4f0f73SDmitry Karpeev           if (!pc->setupcalled) {
2199566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(sviewer, "  Solver not set up yet: PCSetUp() not yet called\n"));
2202fa5cd67SKarl Rupp           } else {
2219566063dSJacob Faibussowitsch             PetscCall(KSPView(osm->ksp[d], sviewer));
2226a4f0f73SDmitry Karpeev           }
2239566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(sviewer, "  - - - - - - - - - - - - - - - - - -\n"));
2249566063dSJacob Faibussowitsch           PetscCall(PetscViewerFlush(sviewer));
2259566063dSJacob Faibussowitsch           PetscCall(PetscViewerRestoreSubViewer(viewer, ((PetscObject)osm->ois[d])->comm, &sviewer));
22606b43e7eSDmitry Karpeev           ++l;
227b1a0a8a3SJed Brown         }
22806b43e7eSDmitry Karpeev       }
2299566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)pc)));
230b1a0a8a3SJed Brown     }
2319566063dSJacob Faibussowitsch     PetscCall(PetscFree2(numbering, permutation));
2329566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
2339566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
2349530cbd7SBarry Smith     /* this line is needed to match the extra PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
2359566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2369530cbd7SBarry Smith   }
2373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
238b1a0a8a3SJed Brown }
239b1a0a8a3SJed Brown 
2408f3b4b4dSDmitry Karpeev PETSC_INTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]);
241af538404SDmitry Karpeev 
24266976f2fSJacob Faibussowitsch static PetscErrorCode PCGASMSetHierarchicalPartitioning(PC pc)
243d71ae5a4SJacob Faibussowitsch {
244ea91fabdSFande Kong   PC_GASM        *osm = (PC_GASM *)pc->data;
245ea91fabdSFande Kong   MatPartitioning part;
246ea91fabdSFande Kong   MPI_Comm        comm;
247ea91fabdSFande Kong   PetscMPIInt     size;
248ea91fabdSFande Kong   PetscInt        nlocalsubdomains, fromrows_localsize;
249ea91fabdSFande Kong   IS              partitioning, fromrows, isn;
250ea91fabdSFande Kong   Vec             outervec;
251ea91fabdSFande Kong 
252ea91fabdSFande Kong   PetscFunctionBegin;
2539566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
2549566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
255ea91fabdSFande Kong   /* we do not need a hierarchical partitioning when
256ea91fabdSFande Kong     * the total number of subdomains is consistent with
257ea91fabdSFande Kong     * the number of MPI tasks.
258ea91fabdSFande Kong     * For the following cases, we do not need to use HP
259ea91fabdSFande Kong     * */
2603ba16761SJacob Faibussowitsch   if (osm->N == PETSC_DETERMINE || osm->N >= size || osm->N == 1) PetscFunctionReturn(PETSC_SUCCESS);
261f1580f4eSBarry Smith   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);
262ea91fabdSFande Kong   nlocalsubdomains = size / osm->N;
263ea91fabdSFande Kong   osm->n           = 1;
2649566063dSJacob Faibussowitsch   PetscCall(MatPartitioningCreate(comm, &part));
2659566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetAdjacency(part, pc->pmat));
2669566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetType(part, MATPARTITIONINGHIERARCH));
2679566063dSJacob Faibussowitsch   PetscCall(MatPartitioningHierarchicalSetNcoarseparts(part, osm->N));
2689566063dSJacob Faibussowitsch   PetscCall(MatPartitioningHierarchicalSetNfineparts(part, nlocalsubdomains));
2699566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetFromOptions(part));
270f1580f4eSBarry Smith   /* get new rank owner number of each vertex */
2719566063dSJacob Faibussowitsch   PetscCall(MatPartitioningApply(part, &partitioning));
2729566063dSJacob Faibussowitsch   PetscCall(ISBuildTwoSided(partitioning, NULL, &fromrows));
2739566063dSJacob Faibussowitsch   PetscCall(ISPartitioningToNumbering(partitioning, &isn));
2749566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isn));
2759566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(fromrows, &fromrows_localsize));
2769566063dSJacob Faibussowitsch   PetscCall(MatPartitioningDestroy(&part));
2779566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pc->pmat, &outervec, NULL));
2789566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(comm, fromrows_localsize, PETSC_DETERMINE, &(osm->pcx)));
2799566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(osm->pcx, &(osm->pcy)));
2809566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(osm->pcx, NULL, outervec, fromrows, &(osm->pctoouter)));
2819566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix(pc->pmat, fromrows, fromrows, MAT_INITIAL_MATRIX, &(osm->permutationP)));
2829566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)fromrows));
283ea91fabdSFande Kong   osm->permutationIS = fromrows;
284ea91fabdSFande Kong   osm->pcmat         = pc->pmat;
2859566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)osm->permutationP));
286ea91fabdSFande Kong   pc->pmat = osm->permutationP;
2879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&outervec));
2889566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&fromrows));
2899566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&partitioning));
290ea91fabdSFande Kong   osm->n = PETSC_DETERMINE;
2913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
292ea91fabdSFande Kong }
293ea91fabdSFande Kong 
294d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_GASM(PC pc)
295d71ae5a4SJacob Faibussowitsch {
296f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM *)pc->data;
297ea91fabdSFande Kong   PetscInt        i, nInnerIndices, nTotalInnerIndices;
298c10bc1a1SDmitry Karpeev   PetscMPIInt     rank, size;
299b1a0a8a3SJed Brown   MatReuse        scall = MAT_REUSE_MATRIX;
300b1a0a8a3SJed Brown   KSP             ksp;
301b1a0a8a3SJed Brown   PC              subpc;
302b1a0a8a3SJed Brown   const char     *prefix, *pprefix;
303f746d493SDmitry Karpeev   Vec             x, y;
3046a4f0f73SDmitry Karpeev   PetscInt        oni;   /* Number of indices in the i-th local outer subdomain.               */
3056a4f0f73SDmitry Karpeev   const PetscInt *oidxi; /* Indices from the i-th subdomain local outer subdomain.             */
3066a4f0f73SDmitry Karpeev   PetscInt        on;    /* Number of indices in the disjoint union of local outer subdomains. */
3076a4f0f73SDmitry Karpeev   PetscInt       *oidx;  /* Indices in the disjoint union of local outer subdomains. */
3086a4f0f73SDmitry Karpeev   IS              gois;  /* Disjoint union the global indices of outer subdomains.             */
3096a4f0f73SDmitry Karpeev   IS              goid;  /* Identity IS of the size of the disjoint union of outer subdomains. */
310f746d493SDmitry Karpeev   PetscScalar    *gxarray, *gyarray;
311930d09c1SFande Kong   PetscInt        gostart; /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- over the disjoint union of outer subdomains. */
3128f3b4b4dSDmitry Karpeev   PetscInt        num_subdomains  = 0;
3130298fd71SBarry Smith   DM             *subdomain_dm    = NULL;
3148f3b4b4dSDmitry Karpeev   char          **subdomain_names = NULL;
315f771a274SFande Kong   PetscInt       *numbering;
3168f3b4b4dSDmitry Karpeev 
317b1a0a8a3SJed Brown   PetscFunctionBegin;
3189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
3199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
320b1a0a8a3SJed Brown   if (!pc->setupcalled) {
321ea91fabdSFande Kong     /* use a hierarchical partitioning */
3221baa6e33SBarry Smith     if (osm->hierarchicalpartitioning) PetscCall(PCGASMSetHierarchicalPartitioning(pc));
3238f3b4b4dSDmitry Karpeev     if (osm->n == PETSC_DETERMINE) {
3248f3b4b4dSDmitry Karpeev       if (osm->N != PETSC_DETERMINE) {
3258f3b4b4dSDmitry Karpeev         /* No local subdomains given, but the desired number of total subdomains is known, so construct them accordingly. */
3269566063dSJacob Faibussowitsch         PetscCall(PCGASMCreateSubdomains(pc->pmat, osm->N, &osm->n, &osm->iis));
3278f3b4b4dSDmitry Karpeev       } else if (osm->dm_subdomains && pc->dm) {
3288f3b4b4dSDmitry Karpeev         /* try pc->dm next, if allowed */
3298f3b4b4dSDmitry Karpeev         PetscInt d;
330a35b7b57SDmitry Karpeev         IS      *inner_subdomain_is, *outer_subdomain_is;
3319566063dSJacob Faibussowitsch         PetscCall(DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm));
3321baa6e33SBarry Smith         if (num_subdomains) PetscCall(PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is));
333a35b7b57SDmitry Karpeev         for (d = 0; d < num_subdomains; ++d) {
3349566063dSJacob Faibussowitsch           if (inner_subdomain_is) PetscCall(ISDestroy(&inner_subdomain_is[d]));
3359566063dSJacob Faibussowitsch           if (outer_subdomain_is) PetscCall(ISDestroy(&outer_subdomain_is[d]));
336fdc48646SDmitry Karpeev         }
3379566063dSJacob Faibussowitsch         PetscCall(PetscFree(inner_subdomain_is));
3389566063dSJacob Faibussowitsch         PetscCall(PetscFree(outer_subdomain_is));
3398f3b4b4dSDmitry Karpeev       } else {
340f1580f4eSBarry Smith         /* still no subdomains; use one per rank */
34144fe18e5SDmitry Karpeev         osm->nmax = osm->n = 1;
3429566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
343f746d493SDmitry Karpeev         osm->N = size;
3449566063dSJacob Faibussowitsch         PetscCall(PCGASMCreateLocalSubdomains(pc->pmat, osm->n, &osm->iis));
345fdc48646SDmitry Karpeev       }
34606b43e7eSDmitry Karpeev     }
347a35b7b57SDmitry Karpeev     if (!osm->iis) {
348a35b7b57SDmitry Karpeev       /*
3498f3b4b4dSDmitry Karpeev        osm->n was set in PCGASMSetSubdomains(), but the actual subdomains have not been supplied.
3508f3b4b4dSDmitry Karpeev        We create the requisite number of local inner subdomains and then expand them into
3518f3b4b4dSDmitry Karpeev        out subdomains, if necessary.
352a35b7b57SDmitry Karpeev        */
3539566063dSJacob Faibussowitsch       PetscCall(PCGASMCreateLocalSubdomains(pc->pmat, osm->n, &osm->iis));
354a35b7b57SDmitry Karpeev     }
3558f3b4b4dSDmitry Karpeev     if (!osm->ois) {
3568f3b4b4dSDmitry Karpeev       /*
3578f3b4b4dSDmitry Karpeev             Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap
3588f3b4b4dSDmitry Karpeev             has been requested, copy the inner subdomains over so they can be modified.
3598f3b4b4dSDmitry Karpeev       */
3609566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(osm->n, &osm->ois));
3618f3b4b4dSDmitry Karpeev       for (i = 0; i < osm->n; ++i) {
362ea91fabdSFande Kong         if (osm->overlap > 0 && osm->N > 1) { /* With positive overlap, osm->iis[i] will be modified */
3639566063dSJacob Faibussowitsch           PetscCall(ISDuplicate(osm->iis[i], (osm->ois) + i));
3649566063dSJacob Faibussowitsch           PetscCall(ISCopy(osm->iis[i], osm->ois[i]));
3658f3b4b4dSDmitry Karpeev         } else {
3669566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)((osm->iis)[i])));
3678f3b4b4dSDmitry Karpeev           osm->ois[i] = osm->iis[i];
3688f3b4b4dSDmitry Karpeev         }
3698f3b4b4dSDmitry Karpeev       }
370ea91fabdSFande Kong       if (osm->overlap > 0 && osm->N > 1) {
3718f3b4b4dSDmitry Karpeev         /* Extend the "overlapping" regions by a number of steps */
3729566063dSJacob Faibussowitsch         PetscCall(MatIncreaseOverlapSplit(pc->pmat, osm->n, osm->ois, osm->overlap));
3738f3b4b4dSDmitry Karpeev       }
374b1a0a8a3SJed Brown     }
3756a4f0f73SDmitry Karpeev 
3768f3b4b4dSDmitry Karpeev     /* Now the subdomains are defined.  Determine their global and max local numbers, if necessary. */
3778f3b4b4dSDmitry Karpeev     if (osm->nmax == PETSC_DETERMINE) {
3788f3b4b4dSDmitry Karpeev       PetscMPIInt inwork, outwork;
3798f3b4b4dSDmitry Karpeev       /* determine global number of subdomains and the max number of local subdomains */
3808f3b4b4dSDmitry Karpeev       inwork = osm->n;
3811c2dc1cbSBarry Smith       PetscCall(MPIU_Allreduce(&inwork, &outwork, 1, MPI_INT, MPI_MAX, PetscObjectComm((PetscObject)pc)));
3828f3b4b4dSDmitry Karpeev       osm->nmax = outwork;
3838f3b4b4dSDmitry Karpeev     }
3848f3b4b4dSDmitry Karpeev     if (osm->N == PETSC_DETERMINE) {
3858f3b4b4dSDmitry Karpeev       /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */
3869566063dSJacob Faibussowitsch       PetscCall(PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc), osm->n, (PetscObject *)osm->ois, &osm->N, NULL));
3878f3b4b4dSDmitry Karpeev     }
3888f3b4b4dSDmitry Karpeev 
389b1a0a8a3SJed Brown     if (osm->sort_indices) {
390f746d493SDmitry Karpeev       for (i = 0; i < osm->n; i++) {
3919566063dSJacob Faibussowitsch         PetscCall(ISSort(osm->ois[i]));
3929566063dSJacob Faibussowitsch         PetscCall(ISSort(osm->iis[i]));
393b1a0a8a3SJed Brown       }
394b1a0a8a3SJed Brown     }
3959566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc, &prefix));
3969566063dSJacob Faibussowitsch     PetscCall(PCGASMPrintSubdomains(pc));
3978f3b4b4dSDmitry Karpeev 
3986a4f0f73SDmitry Karpeev     /*
3996a4f0f73SDmitry Karpeev        Merge the ISs, create merged vectors and restrictions.
4006a4f0f73SDmitry Karpeev      */
4016a4f0f73SDmitry Karpeev     /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */
4026a4f0f73SDmitry Karpeev     on = 0;
403f746d493SDmitry Karpeev     for (i = 0; i < osm->n; i++) {
4049566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(osm->ois[i], &oni));
4056a4f0f73SDmitry Karpeev       on += oni;
406f746d493SDmitry Karpeev     }
4079566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(on, &oidx));
4086a4f0f73SDmitry Karpeev     on = 0;
409930d09c1SFande Kong     /* Merge local indices together */
410f746d493SDmitry Karpeev     for (i = 0; i < osm->n; i++) {
4119566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(osm->ois[i], &oni));
4129566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(osm->ois[i], &oidxi));
4139566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(oidx + on, oidxi, oni));
4149566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(osm->ois[i], &oidxi));
4156a4f0f73SDmitry Karpeev       on += oni;
416f746d493SDmitry Karpeev     }
4179566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(((PetscObject)(pc))->comm, on, oidx, PETSC_OWN_POINTER, &gois));
418ea91fabdSFande Kong     nTotalInnerIndices = 0;
419ea91fabdSFande Kong     for (i = 0; i < osm->n; i++) {
4209566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(osm->iis[i], &nInnerIndices));
421ea91fabdSFande Kong       nTotalInnerIndices += nInnerIndices;
422ea91fabdSFande Kong     }
4239566063dSJacob Faibussowitsch     PetscCall(VecCreateMPI(((PetscObject)(pc))->comm, nTotalInnerIndices, PETSC_DETERMINE, &x));
4249566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &y));
425ea91fabdSFande Kong 
4269566063dSJacob Faibussowitsch     PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)pc), on, PETSC_DECIDE, &osm->gx));
4279566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(osm->gx, &osm->gy));
4289566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(osm->gx, &gostart, NULL));
4299566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pc), on, gostart, 1, &goid));
430930d09c1SFande Kong     /* gois might indices not on local */
4319566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(x, gois, osm->gx, goid, &(osm->gorestriction)));
4329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(osm->n, &numbering));
4339566063dSJacob Faibussowitsch     PetscCall(PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc), osm->n, (PetscObject *)osm->ois, NULL, numbering));
4349566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&x));
4359566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&gois));
4362fa5cd67SKarl Rupp 
4376a4f0f73SDmitry Karpeev     /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */
4382fa5cd67SKarl Rupp     {
4392fa5cd67SKarl Rupp       PetscInt        ini;   /* Number of indices the i-th a local inner subdomain. */
4408966356dSPierre Jolivet       PetscInt        in;    /* Number of indices in the disjoint union of local inner subdomains. */
4416a4f0f73SDmitry Karpeev       PetscInt       *iidx;  /* Global indices in the merged local inner subdomain. */
4426a4f0f73SDmitry Karpeev       PetscInt       *ioidx; /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
4436a4f0f73SDmitry Karpeev       IS              giis;  /* IS for the disjoint union of inner subdomains. */
4446a4f0f73SDmitry Karpeev       IS              giois; /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
445f771a274SFande Kong       PetscScalar    *array;
446f771a274SFande Kong       const PetscInt *indices;
447f771a274SFande Kong       PetscInt        k;
4486a4f0f73SDmitry Karpeev       on = 0;
449f746d493SDmitry Karpeev       for (i = 0; i < osm->n; i++) {
4509566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(osm->ois[i], &oni));
4516a4f0f73SDmitry Karpeev         on += oni;
452f746d493SDmitry Karpeev       }
4539566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(on, &iidx));
4549566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(on, &ioidx));
4559566063dSJacob Faibussowitsch       PetscCall(VecGetArray(y, &array));
456e12b4729SFande Kong       /* set communicator id to determine where overlap is */
457f771a274SFande Kong       in = 0;
458f771a274SFande Kong       for (i = 0; i < osm->n; i++) {
4599566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(osm->iis[i], &ini));
460ad540459SPierre Jolivet         for (k = 0; k < ini; ++k) array[in + k] = numbering[i];
461f771a274SFande Kong         in += ini;
462f771a274SFande Kong       }
4639566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(y, &array));
4649566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->gorestriction, y, osm->gy, INSERT_VALUES, SCATTER_FORWARD));
4659566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->gorestriction, y, osm->gy, INSERT_VALUES, SCATTER_FORWARD));
4669566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(osm->gy, &gostart, NULL));
4679566063dSJacob Faibussowitsch       PetscCall(VecGetArray(osm->gy, &array));
468f771a274SFande Kong       on = 0;
469f771a274SFande Kong       in = 0;
470f771a274SFande Kong       for (i = 0; i < osm->n; i++) {
4719566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(osm->ois[i], &oni));
4729566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(osm->ois[i], &indices));
473f771a274SFande Kong         for (k = 0; k < oni; k++) {
474e12b4729SFande Kong           /*  skip overlapping indices to get inner domain */
47543081b6cSDmitry Karpeev           if (PetscRealPart(array[on + k]) != numbering[i]) continue;
476f771a274SFande Kong           iidx[in]    = indices[k];
477f771a274SFande Kong           ioidx[in++] = gostart + on + k;
478f771a274SFande Kong         }
4799566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(osm->ois[i], &indices));
480f771a274SFande Kong         on += oni;
481f771a274SFande Kong       }
4829566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(osm->gy, &array));
4839566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc), in, iidx, PETSC_OWN_POINTER, &giis));
4849566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc), in, ioidx, PETSC_OWN_POINTER, &giois));
4859566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(y, giis, osm->gy, giois, &osm->girestriction));
4869566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&y));
4879566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&giis));
4889566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&giois));
489b1a0a8a3SJed Brown     }
4909566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&goid));
4919566063dSJacob Faibussowitsch     PetscCall(PetscFree(numbering));
4922fa5cd67SKarl Rupp 
493f746d493SDmitry Karpeev     /* Create the subdomain work vectors. */
4949566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(osm->n, &osm->x));
4959566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(osm->n, &osm->y));
4969566063dSJacob Faibussowitsch     PetscCall(VecGetArray(osm->gx, &gxarray));
4979566063dSJacob Faibussowitsch     PetscCall(VecGetArray(osm->gy, &gyarray));
4986a4f0f73SDmitry Karpeev     for (i = 0, on = 0; i < osm->n; ++i, on += oni) {
4996a4f0f73SDmitry Karpeev       PetscInt oNi;
5009566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(osm->ois[i], &oni));
501930d09c1SFande Kong       /* on a sub communicator */
5029566063dSJacob Faibussowitsch       PetscCall(ISGetSize(osm->ois[i], &oNi));
5039566063dSJacob Faibussowitsch       PetscCall(VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm, 1, oni, oNi, gxarray + on, &osm->x[i]));
5049566063dSJacob Faibussowitsch       PetscCall(VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm, 1, oni, oNi, gyarray + on, &osm->y[i]));
505b1a0a8a3SJed Brown     }
5069566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(osm->gx, &gxarray));
5079566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(osm->gy, &gyarray));
5088f3b4b4dSDmitry Karpeev     /* Create the subdomain solvers */
5099566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(osm->n, &osm->ksp));
5108f3b4b4dSDmitry Karpeev     for (i = 0; i < osm->n; i++) {
5118f3b4b4dSDmitry Karpeev       char subprefix[PETSC_MAX_PATH_LEN + 1];
5129566063dSJacob Faibussowitsch       PetscCall(KSPCreate(((PetscObject)(osm->ois[i]))->comm, &ksp));
5133821be0aSBarry Smith       PetscCall(KSPSetNestLevel(ksp, pc->kspnestlevel));
5149566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
5159566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
5169566063dSJacob Faibussowitsch       PetscCall(KSPSetType(ksp, KSPPREONLY));
5179566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(ksp, &subpc)); /* Why do we need this here? */
5188f3b4b4dSDmitry Karpeev       if (subdomain_dm) {
5199566063dSJacob Faibussowitsch         PetscCall(KSPSetDM(ksp, subdomain_dm[i]));
5209566063dSJacob Faibussowitsch         PetscCall(DMDestroy(subdomain_dm + i));
5218f3b4b4dSDmitry Karpeev       }
5229566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
5239566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(ksp, prefix));
5248f3b4b4dSDmitry Karpeev       if (subdomain_names && subdomain_names[i]) {
5259566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(subprefix, PETSC_MAX_PATH_LEN, "sub_%s_", subdomain_names[i]));
5269566063dSJacob Faibussowitsch         PetscCall(KSPAppendOptionsPrefix(ksp, subprefix));
5279566063dSJacob Faibussowitsch         PetscCall(PetscFree(subdomain_names[i]));
5288f3b4b4dSDmitry Karpeev       }
5299566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
530b1a0a8a3SJed Brown       osm->ksp[i] = ksp;
531b1a0a8a3SJed Brown     }
5329566063dSJacob Faibussowitsch     PetscCall(PetscFree(subdomain_dm));
5339566063dSJacob Faibussowitsch     PetscCall(PetscFree(subdomain_names));
534b1a0a8a3SJed Brown     scall = MAT_INITIAL_MATRIX;
5358f3b4b4dSDmitry Karpeev   } else { /* if (pc->setupcalled) */
536b1a0a8a3SJed Brown     /*
5378f3b4b4dSDmitry Karpeev        Destroy the submatrices from the previous iteration
538b1a0a8a3SJed Brown     */
539b1a0a8a3SJed Brown     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
5409566063dSJacob Faibussowitsch       PetscCall(MatDestroyMatrices(osm->n, &osm->pmat));
541b1a0a8a3SJed Brown       scall = MAT_INITIAL_MATRIX;
542b1a0a8a3SJed Brown     }
543ea91fabdSFande Kong     if (osm->permutationIS) {
5449566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(pc->pmat, osm->permutationIS, osm->permutationIS, scall, &osm->permutationP));
5459566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)osm->permutationP));
546ea91fabdSFande Kong       osm->pcmat = pc->pmat;
547ea91fabdSFande Kong       pc->pmat   = osm->permutationP;
548b1a0a8a3SJed Brown     }
549ea91fabdSFande Kong   }
550ea91fabdSFande Kong 
551b1a0a8a3SJed Brown   /*
5522da392ccSBarry Smith      Extract the submatrices.
553b1a0a8a3SJed Brown   */
55481a5c4aaSDmitry Karpeev   if (size > 1) {
5559566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatricesMPI(pc->pmat, osm->n, osm->ois, osm->ois, scall, &osm->pmat));
5562fa5cd67SKarl Rupp   } else {
5579566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrices(pc->pmat, osm->n, osm->ois, osm->ois, scall, &osm->pmat));
55881a5c4aaSDmitry Karpeev   }
559b1a0a8a3SJed Brown   if (scall == MAT_INITIAL_MATRIX) {
5609566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc->pmat, &pprefix));
561aa624791SPierre Jolivet     for (i = 0; i < osm->n; i++) PetscCall(PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i], pprefix));
562b1a0a8a3SJed Brown   }
563b1a0a8a3SJed Brown 
564b1a0a8a3SJed Brown   /* Return control to the user so that the submatrices can be modified (e.g., to apply
565b1a0a8a3SJed Brown      different boundary conditions for the submatrices than for the global problem) */
5669566063dSJacob Faibussowitsch   PetscCall(PCModifySubMatrices(pc, osm->n, osm->ois, osm->ois, osm->pmat, pc->modifysubmatricesP));
567b1a0a8a3SJed Brown 
568b1a0a8a3SJed Brown   /*
5696a4f0f73SDmitry Karpeev      Loop over submatrices putting them into local ksps
570b1a0a8a3SJed Brown   */
571f746d493SDmitry Karpeev   for (i = 0; i < osm->n; i++) {
5729566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(osm->ksp[i], osm->pmat[i], osm->pmat[i]));
5739566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(osm->ksp[i], &prefix));
5749566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(osm->pmat[i], prefix));
57548a46eb9SPierre Jolivet     if (!pc->setupcalled) PetscCall(KSPSetFromOptions(osm->ksp[i]));
576b1a0a8a3SJed Brown   }
577ea91fabdSFande Kong   if (osm->pcmat) {
5789566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&pc->pmat));
579ea91fabdSFande Kong     pc->pmat   = osm->pcmat;
5800a545947SLisandro Dalcin     osm->pcmat = NULL;
581ea91fabdSFande Kong   }
5823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
583b1a0a8a3SJed Brown }
584b1a0a8a3SJed Brown 
585d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc)
586d71ae5a4SJacob Faibussowitsch {
587f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM *)pc->data;
588b1a0a8a3SJed Brown   PetscInt i;
589b1a0a8a3SJed Brown 
590b1a0a8a3SJed Brown   PetscFunctionBegin;
59148a46eb9SPierre Jolivet   for (i = 0; i < osm->n; i++) PetscCall(KSPSetUp(osm->ksp[i]));
5923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
593b1a0a8a3SJed Brown }
594b1a0a8a3SJed Brown 
595d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_GASM(PC pc, Vec xin, Vec yout)
596d71ae5a4SJacob Faibussowitsch {
597f746d493SDmitry Karpeev   PC_GASM    *osm = (PC_GASM *)pc->data;
598f746d493SDmitry Karpeev   PetscInt    i;
599ea91fabdSFande Kong   Vec         x, y;
600b1a0a8a3SJed Brown   ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE;
601b1a0a8a3SJed Brown 
602b1a0a8a3SJed Brown   PetscFunctionBegin;
603ea91fabdSFande Kong   if (osm->pctoouter) {
6049566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(osm->pctoouter, xin, osm->pcx, INSERT_VALUES, SCATTER_REVERSE));
6059566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(osm->pctoouter, xin, osm->pcx, INSERT_VALUES, SCATTER_REVERSE));
606ea91fabdSFande Kong     x = osm->pcx;
607ea91fabdSFande Kong     y = osm->pcy;
608ea91fabdSFande Kong   } else {
609ea91fabdSFande Kong     x = xin;
610ea91fabdSFande Kong     y = yout;
611ea91fabdSFande Kong   }
612b1a0a8a3SJed Brown   /*
61348e38a8aSPierre Jolivet      support for limiting the restriction or interpolation only to the inner
614b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
615b1a0a8a3SJed Brown   */
616f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
617b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
6189566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(osm->gx));
6199566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(osm->girestriction, x, osm->gx, INSERT_VALUES, forward));
6202fa5cd67SKarl Rupp   } else {
6219566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(osm->gorestriction, x, osm->gx, INSERT_VALUES, forward));
622b1a0a8a3SJed Brown   }
6239566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(osm->gy));
6246a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
6259566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(osm->girestriction, x, osm->gx, INSERT_VALUES, forward));
6262fa5cd67SKarl Rupp   } else {
6279566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(osm->gorestriction, x, osm->gx, INSERT_VALUES, forward));
6286a4f0f73SDmitry Karpeev   }
62986b96d36SDmitry Karpeev   /* do the subdomain solves */
63086b96d36SDmitry Karpeev   for (i = 0; i < osm->n; ++i) {
6319566063dSJacob Faibussowitsch     PetscCall(KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]));
6329566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(osm->ksp[i], pc, osm->y[i]));
633b1a0a8a3SJed Brown   }
63448e38a8aSPierre Jolivet   /* do we need to zero y? */
6359566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(y));
6366a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
6379566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(osm->girestriction, osm->gy, y, ADD_VALUES, reverse));
6389566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(osm->girestriction, osm->gy, y, ADD_VALUES, reverse));
6392fa5cd67SKarl Rupp   } else {
6409566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(osm->gorestriction, osm->gy, y, ADD_VALUES, reverse));
6419566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(osm->gorestriction, osm->gy, y, ADD_VALUES, reverse));
6426a4f0f73SDmitry Karpeev   }
643ea91fabdSFande Kong   if (osm->pctoouter) {
6449566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(osm->pctoouter, y, yout, INSERT_VALUES, SCATTER_FORWARD));
6459566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(osm->pctoouter, y, yout, INSERT_VALUES, SCATTER_FORWARD));
646ea91fabdSFande Kong   }
6473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
648b1a0a8a3SJed Brown }
649b1a0a8a3SJed Brown 
650d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_GASM(PC pc, Mat Xin, Mat Yout)
651d71ae5a4SJacob Faibussowitsch {
65248e38a8aSPierre Jolivet   PC_GASM    *osm = (PC_GASM *)pc->data;
65348e38a8aSPierre Jolivet   Mat         X, Y, O = NULL, Z, W;
65448e38a8aSPierre Jolivet   Vec         x, y;
65548e38a8aSPierre Jolivet   PetscInt    i, m, M, N;
65648e38a8aSPierre Jolivet   ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE;
65748e38a8aSPierre Jolivet 
65848e38a8aSPierre Jolivet   PetscFunctionBegin;
65908401ef6SPierre Jolivet   PetscCheck(osm->n == 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not yet implemented");
6609566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Xin, NULL, &N));
66148e38a8aSPierre Jolivet   if (osm->pctoouter) {
6629566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(osm->pcx, &m));
6639566063dSJacob Faibussowitsch     PetscCall(VecGetSize(osm->pcx, &M));
6649566063dSJacob Faibussowitsch     PetscCall(MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]), m, PETSC_DECIDE, M, N, NULL, &O));
66548e38a8aSPierre Jolivet     for (i = 0; i < N; ++i) {
6669566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecRead(Xin, i, &x));
6679566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecWrite(O, i, &y));
6689566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->pctoouter, x, y, INSERT_VALUES, SCATTER_REVERSE));
6699566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->pctoouter, x, y, INSERT_VALUES, SCATTER_REVERSE));
6709566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecWrite(O, i, &y));
6719566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecRead(Xin, i, &x));
67248e38a8aSPierre Jolivet     }
67348e38a8aSPierre Jolivet     X = Y = O;
67448e38a8aSPierre Jolivet   } else {
67548e38a8aSPierre Jolivet     X = Xin;
67648e38a8aSPierre Jolivet     Y = Yout;
67748e38a8aSPierre Jolivet   }
67848e38a8aSPierre Jolivet   /*
67948e38a8aSPierre Jolivet      support for limiting the restriction or interpolation only to the inner
68048e38a8aSPierre Jolivet      subdomain values (leaving the other values 0).
68148e38a8aSPierre Jolivet   */
6829566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(osm->x[0], &m));
6839566063dSJacob Faibussowitsch   PetscCall(VecGetSize(osm->x[0], &M));
6849566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]), m, PETSC_DECIDE, M, N, NULL, &Z));
68548e38a8aSPierre Jolivet   for (i = 0; i < N; ++i) {
6869566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumnVecRead(X, i, &x));
6879566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumnVecWrite(Z, i, &y));
68848e38a8aSPierre Jolivet     if (!(osm->type & PC_GASM_RESTRICT)) {
68948e38a8aSPierre Jolivet       /* have to zero the work RHS since scatter may leave some slots empty */
6909566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(y));
6919566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->girestriction, x, y, INSERT_VALUES, forward));
6929566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->girestriction, x, y, INSERT_VALUES, forward));
69348e38a8aSPierre Jolivet     } else {
6949566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->gorestriction, x, y, INSERT_VALUES, forward));
6959566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->gorestriction, x, y, INSERT_VALUES, forward));
69648e38a8aSPierre Jolivet     }
6979566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumnVecWrite(Z, i, &y));
6989566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumnVecRead(X, i, &x));
69948e38a8aSPierre Jolivet   }
7009566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]), m, PETSC_DECIDE, M, N, NULL, &W));
7019566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Z, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
7029566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Z, MAT_FINAL_ASSEMBLY));
7039566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Z, MAT_FINAL_ASSEMBLY));
70448e38a8aSPierre Jolivet   /* do the subdomain solve */
7059566063dSJacob Faibussowitsch   PetscCall(KSPMatSolve(osm->ksp[0], Z, W));
7069566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(osm->ksp[0], pc, NULL));
7079566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Z));
70848e38a8aSPierre Jolivet   /* do we need to zero y? */
7099566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(Y));
71048e38a8aSPierre Jolivet   for (i = 0; i < N; ++i) {
7119566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumnVecWrite(Y, i, &y));
7129566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumnVecRead(W, i, &x));
71348e38a8aSPierre Jolivet     if (!(osm->type & PC_GASM_INTERPOLATE)) {
7149566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->girestriction, x, y, ADD_VALUES, reverse));
7159566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->girestriction, x, y, ADD_VALUES, reverse));
71648e38a8aSPierre Jolivet     } else {
7179566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->gorestriction, x, y, ADD_VALUES, reverse));
7189566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->gorestriction, x, y, ADD_VALUES, reverse));
71948e38a8aSPierre Jolivet     }
7209566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumnVecRead(W, i, &x));
72148e38a8aSPierre Jolivet     if (osm->pctoouter) {
7229566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecWrite(Yout, i, &x));
7239566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->pctoouter, y, x, INSERT_VALUES, SCATTER_FORWARD));
7249566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->pctoouter, y, x, INSERT_VALUES, SCATTER_FORWARD));
7259566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecRead(Yout, i, &x));
72648e38a8aSPierre Jolivet     }
7279566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumnVecWrite(Y, i, &y));
72848e38a8aSPierre Jolivet   }
7299566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&W));
7309566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&O));
7313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
73248e38a8aSPierre Jolivet }
73348e38a8aSPierre Jolivet 
734d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_GASM(PC pc, Vec xin, Vec yout)
735d71ae5a4SJacob Faibussowitsch {
736f746d493SDmitry Karpeev   PC_GASM    *osm = (PC_GASM *)pc->data;
737f746d493SDmitry Karpeev   PetscInt    i;
738ea91fabdSFande Kong   Vec         x, y;
739b1a0a8a3SJed Brown   ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE;
740b1a0a8a3SJed Brown 
741b1a0a8a3SJed Brown   PetscFunctionBegin;
742ea91fabdSFande Kong   if (osm->pctoouter) {
7439566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(osm->pctoouter, xin, osm->pcx, INSERT_VALUES, SCATTER_REVERSE));
7449566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(osm->pctoouter, xin, osm->pcx, INSERT_VALUES, SCATTER_REVERSE));
745ea91fabdSFande Kong     x = osm->pcx;
746ea91fabdSFande Kong     y = osm->pcy;
747ea91fabdSFande Kong   } else {
748ea91fabdSFande Kong     x = xin;
749ea91fabdSFande Kong     y = yout;
750ea91fabdSFande Kong   }
751b1a0a8a3SJed Brown   /*
752b1a0a8a3SJed Brown      Support for limiting the restriction or interpolation to only local
753b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
754b1a0a8a3SJed Brown 
755f746d493SDmitry Karpeev      Note: these are reversed from the PCApply_GASM() because we are applying the
756b1a0a8a3SJed Brown      transpose of the three terms
757b1a0a8a3SJed Brown   */
758f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
759b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
7609566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(osm->gx));
7619566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(osm->girestriction, x, osm->gx, INSERT_VALUES, forward));
7622fa5cd67SKarl Rupp   } else {
7639566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(osm->gorestriction, x, osm->gx, INSERT_VALUES, forward));
764b1a0a8a3SJed Brown   }
7659566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(osm->gy));
7666a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
7679566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(osm->girestriction, x, osm->gx, INSERT_VALUES, forward));
7682fa5cd67SKarl Rupp   } else {
7699566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(osm->gorestriction, x, osm->gx, INSERT_VALUES, forward));
7706a4f0f73SDmitry Karpeev   }
771b1a0a8a3SJed Brown   /* do the local solves */
772ab3e923fSDmitry Karpeev   for (i = 0; i < osm->n; ++i) { /* Note that the solves are local, so we can go to osm->n, rather than osm->nmax. */
7739566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i]));
7749566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(osm->ksp[i], pc, osm->y[i]));
775b1a0a8a3SJed Brown   }
7769566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(y));
7776a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
7789566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(osm->girestriction, osm->gy, y, ADD_VALUES, reverse));
7799566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(osm->girestriction, osm->gy, y, ADD_VALUES, reverse));
7802fa5cd67SKarl Rupp   } else {
7819566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(osm->gorestriction, osm->gy, y, ADD_VALUES, reverse));
7829566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(osm->gorestriction, osm->gy, y, ADD_VALUES, reverse));
7836a4f0f73SDmitry Karpeev   }
784ea91fabdSFande Kong   if (osm->pctoouter) {
7859566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(osm->pctoouter, y, yout, INSERT_VALUES, SCATTER_FORWARD));
7869566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(osm->pctoouter, y, yout, INSERT_VALUES, SCATTER_FORWARD));
787ea91fabdSFande Kong   }
7883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
789b1a0a8a3SJed Brown }
790b1a0a8a3SJed Brown 
791d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_GASM(PC pc)
792d71ae5a4SJacob Faibussowitsch {
793f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM *)pc->data;
794b1a0a8a3SJed Brown   PetscInt i;
795b1a0a8a3SJed Brown 
796b1a0a8a3SJed Brown   PetscFunctionBegin;
797b1a0a8a3SJed Brown   if (osm->ksp) {
79848a46eb9SPierre Jolivet     for (i = 0; i < osm->n; i++) PetscCall(KSPReset(osm->ksp[i]));
799b1a0a8a3SJed Brown   }
800b1a0a8a3SJed Brown   if (osm->pmat) {
801f746d493SDmitry Karpeev     if (osm->n > 0) {
802df750dc8SHong Zhang       PetscMPIInt size;
8039566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
804df750dc8SHong Zhang       if (size > 1) {
8057dae84e0SHong Zhang         /* osm->pmat is created by MatCreateSubMatricesMPI(), cannot use MatDestroySubMatrices() */
8069566063dSJacob Faibussowitsch         PetscCall(MatDestroyMatrices(osm->n, &osm->pmat));
807df750dc8SHong Zhang       } else {
8089566063dSJacob Faibussowitsch         PetscCall(MatDestroySubMatrices(osm->n, &osm->pmat));
809df750dc8SHong Zhang       }
810b1a0a8a3SJed Brown     }
811b1a0a8a3SJed Brown   }
812d34fcf5fSBarry Smith   if (osm->x) {
813f746d493SDmitry Karpeev     for (i = 0; i < osm->n; i++) {
8149566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&osm->x[i]));
8159566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&osm->y[i]));
816b1a0a8a3SJed Brown     }
817d34fcf5fSBarry Smith   }
8189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&osm->gx));
8199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&osm->gy));
820ab3e923fSDmitry Karpeev 
8219566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&osm->gorestriction));
8229566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&osm->girestriction));
8238f3b4b4dSDmitry Karpeev   if (!osm->user_subdomains) {
8249566063dSJacob Faibussowitsch     PetscCall(PCGASMDestroySubdomains(osm->n, &osm->ois, &osm->iis));
8258f3b4b4dSDmitry Karpeev     osm->N    = PETSC_DETERMINE;
8268f3b4b4dSDmitry Karpeev     osm->nmax = PETSC_DETERMINE;
8278f3b4b4dSDmitry Karpeev   }
82848a46eb9SPierre Jolivet   if (osm->pctoouter) PetscCall(VecScatterDestroy(&(osm->pctoouter)));
82948a46eb9SPierre Jolivet   if (osm->permutationIS) PetscCall(ISDestroy(&(osm->permutationIS)));
83048a46eb9SPierre Jolivet   if (osm->pcx) PetscCall(VecDestroy(&(osm->pcx)));
83148a46eb9SPierre Jolivet   if (osm->pcy) PetscCall(VecDestroy(&(osm->pcy)));
83248a46eb9SPierre Jolivet   if (osm->permutationP) PetscCall(MatDestroy(&(osm->permutationP)));
83348a46eb9SPierre Jolivet   if (osm->pcmat) PetscCall(MatDestroy(&osm->pcmat));
8343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
835a06653b4SBarry Smith }
836a06653b4SBarry Smith 
837d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_GASM(PC pc)
838d71ae5a4SJacob Faibussowitsch {
839a06653b4SBarry Smith   PC_GASM *osm = (PC_GASM *)pc->data;
840a06653b4SBarry Smith   PetscInt i;
841a06653b4SBarry Smith 
842a06653b4SBarry Smith   PetscFunctionBegin;
8439566063dSJacob Faibussowitsch   PetscCall(PCReset_GASM(pc));
844135757f6SDmitry Karpeev   /* PCReset will not destroy subdomains, if user_subdomains is true. */
8459566063dSJacob Faibussowitsch   PetscCall(PCGASMDestroySubdomains(osm->n, &osm->ois, &osm->iis));
846a06653b4SBarry Smith   if (osm->ksp) {
84748a46eb9SPierre Jolivet     for (i = 0; i < osm->n; i++) PetscCall(KSPDestroy(&osm->ksp[i]));
8489566063dSJacob Faibussowitsch     PetscCall(PetscFree(osm->ksp));
849a06653b4SBarry Smith   }
8509566063dSJacob Faibussowitsch   PetscCall(PetscFree(osm->x));
8519566063dSJacob Faibussowitsch   PetscCall(PetscFree(osm->y));
8522e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMSetSubdomains_C", NULL));
8532e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMSetOverlap_C", NULL));
8542e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMSetType_C", NULL));
8552e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMSetSortIndices_C", NULL));
8562e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMGetSubKSP_C", NULL));
8579566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
8583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
859b1a0a8a3SJed Brown }
860b1a0a8a3SJed Brown 
861d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_GASM(PC pc, PetscOptionItems *PetscOptionsObject)
862d71ae5a4SJacob Faibussowitsch {
863f746d493SDmitry Karpeev   PC_GASM   *osm = (PC_GASM *)pc->data;
864b1a0a8a3SJed Brown   PetscInt   blocks, ovl;
86557501b6eSBarry Smith   PetscBool  flg;
866f746d493SDmitry Karpeev   PCGASMType gasmtype;
867b1a0a8a3SJed Brown 
868b1a0a8a3SJed Brown   PetscFunctionBegin;
869d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Generalized additive Schwarz options");
8709566063dSJacob Faibussowitsch   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));
8719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_gasm_total_subdomains", "Total number of subdomains across communicator", "PCGASMSetTotalSubdomains", osm->N, &blocks, &flg));
8721baa6e33SBarry Smith   if (flg) PetscCall(PCGASMSetTotalSubdomains(pc, blocks));
8739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_gasm_overlap", "Number of overlapping degrees of freedom", "PCGASMSetOverlap", osm->overlap, &ovl, &flg));
87465db9045SDmitry Karpeev   if (flg) {
8759566063dSJacob Faibussowitsch     PetscCall(PCGASMSetOverlap(pc, ovl));
876d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
87765db9045SDmitry Karpeev   }
878b1a0a8a3SJed Brown   flg = PETSC_FALSE;
8799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_gasm_type", "Type of restriction/extension", "PCGASMSetType", PCGASMTypes, (PetscEnum)osm->type, (PetscEnum *)&gasmtype, &flg));
8809566063dSJacob Faibussowitsch   if (flg) PetscCall(PCGASMSetType(pc, gasmtype));
8819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_gasm_use_hierachical_partitioning", "use hierarchical partitioning", NULL, osm->hierarchicalpartitioning, &osm->hierarchicalpartitioning, &flg));
882d0609cedSBarry Smith   PetscOptionsHeadEnd();
8833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
884b1a0a8a3SJed Brown }
885b1a0a8a3SJed Brown 
8868f3b4b4dSDmitry Karpeev /*@
887f1580f4eSBarry Smith   PCGASMSetTotalSubdomains - sets the total number of subdomains to use across the communicator for `PCGASM`
888c3339decSBarry Smith 
889c3339decSBarry Smith   Logically Collective
8908f3b4b4dSDmitry Karpeev 
8918f3b4b4dSDmitry Karpeev   Input Parameters:
8928f3b4b4dSDmitry Karpeev + pc - the preconditioner
8938f3b4b4dSDmitry Karpeev - N  - total number of subdomains
8948f3b4b4dSDmitry Karpeev 
8958f3b4b4dSDmitry Karpeev   Level: beginner
8968f3b4b4dSDmitry Karpeev 
897f1580f4eSBarry Smith .seealso: `PCGASM`, `PCGASMSetSubdomains()`, `PCGASMSetOverlap()`
898db781477SPatrick Sanan           `PCGASMCreateSubdomains2D()`
8998f3b4b4dSDmitry Karpeev @*/
900d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGASMSetTotalSubdomains(PC pc, PetscInt N)
901d71ae5a4SJacob Faibussowitsch {
9028f3b4b4dSDmitry Karpeev   PC_GASM    *osm = (PC_GASM *)pc->data;
9038f3b4b4dSDmitry Karpeev   PetscMPIInt size, rank;
9048f3b4b4dSDmitry Karpeev 
9058f3b4b4dSDmitry Karpeev   PetscFunctionBegin;
90663a3b9bcSJacob Faibussowitsch   PetscCheck(N >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Total number of subdomains must be 1 or more, got N = %" PetscInt_FMT, N);
90728b400f6SJacob Faibussowitsch   PetscCheck(!pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCGASMSetTotalSubdomains() should be called before calling PCSetUp().");
9088f3b4b4dSDmitry Karpeev 
9099566063dSJacob Faibussowitsch   PetscCall(PCGASMDestroySubdomains(osm->n, &osm->iis, &osm->ois));
9108f3b4b4dSDmitry Karpeev   osm->ois = osm->iis = NULL;
9118f3b4b4dSDmitry Karpeev 
9129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
9139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
9148f3b4b4dSDmitry Karpeev   osm->N             = N;
9158f3b4b4dSDmitry Karpeev   osm->n             = PETSC_DETERMINE;
9168f3b4b4dSDmitry Karpeev   osm->nmax          = PETSC_DETERMINE;
9178f3b4b4dSDmitry Karpeev   osm->dm_subdomains = PETSC_FALSE;
9183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9198f3b4b4dSDmitry Karpeev }
9208f3b4b4dSDmitry Karpeev 
921d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGASMSetSubdomains_GASM(PC pc, PetscInt n, IS iis[], IS ois[])
922d71ae5a4SJacob Faibussowitsch {
923f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM *)pc->data;
924b1a0a8a3SJed Brown   PetscInt i;
925b1a0a8a3SJed Brown 
926b1a0a8a3SJed Brown   PetscFunctionBegin;
927f1580f4eSBarry Smith   PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Each MPI rank must have 1 or more subdomains, got n = %" PetscInt_FMT, n);
92828b400f6SJacob Faibussowitsch   PetscCheck(!pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCGASMSetSubdomains() should be called before calling PCSetUp().");
929b1a0a8a3SJed Brown 
9309566063dSJacob Faibussowitsch   PetscCall(PCGASMDestroySubdomains(osm->n, &osm->iis, &osm->ois));
9318f3b4b4dSDmitry Karpeev   osm->iis = osm->ois = NULL;
9328f3b4b4dSDmitry Karpeev   osm->n              = n;
9338f3b4b4dSDmitry Karpeev   osm->N              = PETSC_DETERMINE;
9348f3b4b4dSDmitry Karpeev   osm->nmax           = PETSC_DETERMINE;
935a35b7b57SDmitry Karpeev   if (ois) {
9369566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &osm->ois));
9378f3b4b4dSDmitry Karpeev     for (i = 0; i < n; i++) {
9389566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)ois[i]));
9398f3b4b4dSDmitry Karpeev       osm->ois[i] = ois[i];
9408f3b4b4dSDmitry Karpeev     }
9418f3b4b4dSDmitry Karpeev     /*
9428f3b4b4dSDmitry Karpeev        Since the user set the outer subdomains, even if nontrivial overlap was requested via PCGASMSetOverlap(),
9438f3b4b4dSDmitry Karpeev        it will be ignored.  To avoid confusion later on (e.g., when viewing the PC), the overlap size is set to -1.
9448f3b4b4dSDmitry Karpeev     */
945b1a0a8a3SJed Brown     osm->overlap = -1;
946670417b2SFande Kong     /* inner subdomains must be provided  */
94728b400f6SJacob Faibussowitsch     PetscCheck(iis, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "inner indices have to be provided ");
948670417b2SFande Kong   } /* end if */
949a35b7b57SDmitry Karpeev   if (iis) {
9509566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &osm->iis));
9518f3b4b4dSDmitry Karpeev     for (i = 0; i < n; i++) {
9529566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)iis[i]));
9538f3b4b4dSDmitry Karpeev       osm->iis[i] = iis[i];
9548f3b4b4dSDmitry Karpeev     }
955a35b7b57SDmitry Karpeev     if (!ois) {
956390e1bf2SBarry Smith       osm->ois = NULL;
957670417b2SFande Kong       /* if user does not provide outer indices, we will create the corresponding outer indices using  osm->overlap =1 in PCSetUp_GASM */
958670417b2SFande Kong     }
959670417b2SFande Kong   }
96076bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
961670417b2SFande Kong     PetscInt        j, rstart, rend, *covered, lsize;
962670417b2SFande Kong     const PetscInt *indices;
963670417b2SFande Kong     /* check if the inner indices cover and only cover the local portion of the preconditioning matrix */
9649566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(pc->pmat, &rstart, &rend));
9659566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(rend - rstart, &covered));
966f1580f4eSBarry Smith     /* check if the current MPI rank owns indices from others */
967a35b7b57SDmitry Karpeev     for (i = 0; i < n; i++) {
9689566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(osm->iis[i], &indices));
9699566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(osm->iis[i], &lsize));
970670417b2SFande Kong       for (j = 0; j < lsize; j++) {
971f1580f4eSBarry Smith         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]);
9722472a847SBarry Smith         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]);
97363a3b9bcSJacob Faibussowitsch         covered[indices[j] - rstart] = 1;
974a35b7b57SDmitry Karpeev       }
9759566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(osm->iis[i], &indices));
9768f3b4b4dSDmitry Karpeev     }
977670417b2SFande Kong     /* check if we miss any indices */
978ad540459SPierre Jolivet     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);
9799566063dSJacob Faibussowitsch     PetscCall(PetscFree(covered));
980a35b7b57SDmitry Karpeev   }
9818f3b4b4dSDmitry Karpeev   if (iis) osm->user_subdomains = PETSC_TRUE;
9823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
983b1a0a8a3SJed Brown }
984b1a0a8a3SJed Brown 
985d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGASMSetOverlap_GASM(PC pc, PetscInt ovl)
986d71ae5a4SJacob Faibussowitsch {
987f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM *)pc->data;
988b1a0a8a3SJed Brown 
989b1a0a8a3SJed Brown   PetscFunctionBegin;
99008401ef6SPierre Jolivet   PetscCheck(ovl >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap value requested");
9912472a847SBarry Smith   PetscCheck(!pc->setupcalled || ovl == osm->overlap, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCGASMSetOverlap() should be called before PCSetUp().");
9922fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
9933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
994b1a0a8a3SJed Brown }
995b1a0a8a3SJed Brown 
996d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGASMSetType_GASM(PC pc, PCGASMType type)
997d71ae5a4SJacob Faibussowitsch {
998f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM *)pc->data;
999b1a0a8a3SJed Brown 
1000b1a0a8a3SJed Brown   PetscFunctionBegin;
1001b1a0a8a3SJed Brown   osm->type     = type;
1002b1a0a8a3SJed Brown   osm->type_set = PETSC_TRUE;
10033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1004b1a0a8a3SJed Brown }
1005b1a0a8a3SJed Brown 
1006d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGASMSetSortIndices_GASM(PC pc, PetscBool doSort)
1007d71ae5a4SJacob Faibussowitsch {
1008f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM *)pc->data;
1009b1a0a8a3SJed Brown 
1010b1a0a8a3SJed Brown   PetscFunctionBegin;
1011b1a0a8a3SJed Brown   osm->sort_indices = doSort;
10123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1013b1a0a8a3SJed Brown }
1014b1a0a8a3SJed Brown 
101544fe18e5SDmitry Karpeev /*
10168f3b4b4dSDmitry Karpeev    FIXME: This routine might need to be modified now that multiple ranks per subdomain are allowed.
101744fe18e5SDmitry Karpeev         In particular, it would upset the global subdomain number calculation.
101844fe18e5SDmitry Karpeev */
1019d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGASMGetSubKSP_GASM(PC pc, PetscInt *n, PetscInt *first, KSP **ksp)
1020d71ae5a4SJacob Faibussowitsch {
1021f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM *)pc->data;
1022b1a0a8a3SJed Brown 
1023b1a0a8a3SJed Brown   PetscFunctionBegin;
102408401ef6SPierre Jolivet   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");
1025b1a0a8a3SJed Brown 
10262fa5cd67SKarl Rupp   if (n) *n = osm->n;
1027ab3e923fSDmitry Karpeev   if (first) {
10289566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Scan(&osm->n, first, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)pc)));
1029ab3e923fSDmitry Karpeev     *first -= osm->n;
1030b1a0a8a3SJed Brown   }
1031b1a0a8a3SJed Brown   if (ksp) {
1032b1a0a8a3SJed Brown     /* Assume that local solves are now different; not necessarily
103306b43e7eSDmitry Karpeev        true, though!  This flag is used only for PCView_GASM() */
1034b1a0a8a3SJed Brown     *ksp                        = osm->ksp;
10356a4f0f73SDmitry Karpeev     osm->same_subdomain_solvers = PETSC_FALSE;
1036b1a0a8a3SJed Brown   }
10373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1038ab3e923fSDmitry Karpeev } /* PCGASMGetSubKSP_GASM() */
1039b1a0a8a3SJed Brown 
1040b1a0a8a3SJed Brown /*@C
1041f1580f4eSBarry Smith   PCGASMSetSubdomains - Sets the subdomains for this MPI rank
1042f1580f4eSBarry Smith   for the additive Schwarz preconditioner with multiple MPI ranks per subdomain, `PCGASM`
1043b1a0a8a3SJed Brown 
1044c3339decSBarry Smith   Collective
1045b1a0a8a3SJed Brown 
1046b1a0a8a3SJed Brown   Input Parameters:
10478f3b4b4dSDmitry Karpeev + pc  - the preconditioner object
1048f1580f4eSBarry Smith . n   - the number of subdomains for this MPI rank
104920f4b53cSBarry Smith . iis - the index sets that define the inner subdomains (or `NULL` for PETSc to determine subdomains)
1050f36f9100SBarry Smith - ois - the index sets that define the outer subdomains (or `NULL` to use the same as `iis`, or to construct by expanding `iis` by
1051f36f9100SBarry Smith           the requested overlap)
105220f4b53cSBarry Smith 
105320f4b53cSBarry Smith   Level: advanced
1054b1a0a8a3SJed Brown 
1055b1a0a8a3SJed Brown   Notes:
1056f1580f4eSBarry Smith   The `IS` indices use the parallel, global numbering of the vector entries.
1057f36f9100SBarry Smith 
10586a4f0f73SDmitry Karpeev   Inner subdomains are those where the correction is applied.
1059f36f9100SBarry Smith 
10606a4f0f73SDmitry Karpeev   Outer subdomains are those where the residual necessary to obtain the
1061f1580f4eSBarry Smith   corrections is obtained (see `PCGASMType` for the use of inner/outer subdomains).
1062f36f9100SBarry Smith 
1063f1580f4eSBarry Smith   Both inner and outer subdomains can extend over several MPI ranks.
1064f1580f4eSBarry Smith   This rank's portion of a subdomain is known as a local subdomain.
10656a4f0f73SDmitry Karpeev 
1066f1580f4eSBarry Smith   Inner subdomains can not overlap with each other, do not have any entities from remote ranks,
1067f1580f4eSBarry Smith   and  have to cover the entire local subdomain owned by the current rank. The index sets on each
1068f1580f4eSBarry Smith   rank should be ordered such that the ith local subdomain is connected to the ith remote subdomain
1069f1580f4eSBarry Smith   on another MPI rank.
1070670417b2SFande Kong 
1071f1580f4eSBarry Smith   By default the `PGASM` preconditioner uses 1 (local) subdomain per MPI rank.
10726a4f0f73SDmitry Karpeev 
1073f36f9100SBarry Smith   The `iis` and `ois` arrays may be freed after this call using `PCGASMDestroySubdomains()`
1074f36f9100SBarry Smith 
1075f36f9100SBarry Smith .seealso: `PCGASM`, `PCGASMSetOverlap()`, `PCGASMGetSubKSP()`, `PCGASMDestroySubdomains()`,
1076db781477SPatrick Sanan           `PCGASMCreateSubdomains2D()`, `PCGASMGetSubdomains()`
1077b1a0a8a3SJed Brown @*/
1078d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGASMSetSubdomains(PC pc, PetscInt n, IS iis[], IS ois[])
1079d71ae5a4SJacob Faibussowitsch {
10808f3b4b4dSDmitry Karpeev   PC_GASM *osm = (PC_GASM *)pc->data;
1081b1a0a8a3SJed Brown 
1082b1a0a8a3SJed Brown   PetscFunctionBegin;
1083b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1084cac4c232SBarry Smith   PetscTryMethod(pc, "PCGASMSetSubdomains_C", (PC, PetscInt, IS[], IS[]), (pc, n, iis, ois));
10858f3b4b4dSDmitry Karpeev   osm->dm_subdomains = PETSC_FALSE;
10863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1087b1a0a8a3SJed Brown }
1088b1a0a8a3SJed Brown 
1089b1a0a8a3SJed Brown /*@
1090f746d493SDmitry Karpeev   PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the
1091f1580f4eSBarry Smith   additive Schwarz preconditioner `PCGASM`.  Either all or no MPI ranks in the
10928f3b4b4dSDmitry Karpeev   pc communicator must call this routine.
1093b1a0a8a3SJed Brown 
1094c3339decSBarry Smith   Logically Collective
1095b1a0a8a3SJed Brown 
1096b1a0a8a3SJed Brown   Input Parameters:
1097b1a0a8a3SJed Brown + pc  - the preconditioner context
10988f3b4b4dSDmitry Karpeev - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 0)
1099b1a0a8a3SJed Brown 
1100b1a0a8a3SJed Brown   Options Database Key:
110106b43e7eSDmitry Karpeev . -pc_gasm_overlap <overlap> - Sets overlap
1102b1a0a8a3SJed Brown 
110320f4b53cSBarry Smith   Level: intermediate
110420f4b53cSBarry Smith 
1105b1a0a8a3SJed Brown   Notes:
1106f1580f4eSBarry Smith   By default the `PCGASM` preconditioner uses 1 subdomain per rank.  To use
11078f3b4b4dSDmitry Karpeev   multiple subdomain per perocessor or "straddling" subdomains that intersect
1108f36f9100SBarry Smith   multiple ranks use `PCGASMSetSubdomains()` (or option `-pc_gasm_total_subdomains` <n>).
1109b1a0a8a3SJed Brown 
11108f3b4b4dSDmitry Karpeev   The overlap defaults to 0, so if one desires that no additional
1111b1a0a8a3SJed Brown   overlap be computed beyond what may have been set with a call to
1112f36f9100SBarry Smith   `PCGASMSetSubdomains()`, then `ovl` must be set to be 0.  In particular, if one does
11138f3b4b4dSDmitry Karpeev   not explicitly set the subdomains in application code, then all overlap would be computed
1114f1580f4eSBarry Smith   internally by PETSc, and using an overlap of 0 would result in an `PCGASM`
1115b1a0a8a3SJed Brown   variant that is equivalent to the block Jacobi preconditioner.
1116b1a0a8a3SJed Brown 
1117f1580f4eSBarry Smith   One can define initial index sets with any overlap via
1118f1580f4eSBarry Smith   `PCGASMSetSubdomains()`; the routine `PCGASMSetOverlap()` merely allows
111906b43e7eSDmitry Karpeev   PETSc to extend that overlap further, if desired.
1120b1a0a8a3SJed Brown 
1121f1580f4eSBarry Smith .seealso: `PCGASM`, `PCGASMSetSubdomains()`, `PCGASMGetSubKSP()`,
1122db781477SPatrick Sanan           `PCGASMCreateSubdomains2D()`, `PCGASMGetSubdomains()`
1123b1a0a8a3SJed Brown @*/
1124d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGASMSetOverlap(PC pc, PetscInt ovl)
1125d71ae5a4SJacob Faibussowitsch {
11268f3b4b4dSDmitry Karpeev   PC_GASM *osm = (PC_GASM *)pc->data;
1127b1a0a8a3SJed Brown 
1128b1a0a8a3SJed Brown   PetscFunctionBegin;
1129b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1130b1a0a8a3SJed Brown   PetscValidLogicalCollectiveInt(pc, ovl, 2);
1131cac4c232SBarry Smith   PetscTryMethod(pc, "PCGASMSetOverlap_C", (PC, PetscInt), (pc, ovl));
11328f3b4b4dSDmitry Karpeev   osm->dm_subdomains = PETSC_FALSE;
11333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1134b1a0a8a3SJed Brown }
1135b1a0a8a3SJed Brown 
1136b1a0a8a3SJed Brown /*@
1137f746d493SDmitry Karpeev   PCGASMSetType - Sets the type of restriction and interpolation used
1138f1580f4eSBarry Smith   for local problems in the `PCGASM` additive Schwarz method.
1139b1a0a8a3SJed Brown 
1140c3339decSBarry Smith   Logically Collective
1141b1a0a8a3SJed Brown 
1142b1a0a8a3SJed Brown   Input Parameters:
1143b1a0a8a3SJed Brown + pc   - the preconditioner context
1144f1580f4eSBarry Smith - type - variant of `PCGASM`, one of
1145b1a0a8a3SJed Brown .vb
1146f1580f4eSBarry Smith       `PC_GASM_BASIC`       - full interpolation and restriction
1147f1580f4eSBarry Smith       `PC_GASM_RESTRICT`    - full restriction, local MPI rank interpolation
1148f1580f4eSBarry Smith       `PC_GASM_INTERPOLATE` - full interpolation, local MPI rank restriction
1149f1580f4eSBarry Smith       `PC_GASM_NONE`        - local MPI rank restriction and interpolation
1150b1a0a8a3SJed Brown .ve
1151b1a0a8a3SJed Brown 
1152b1a0a8a3SJed Brown   Options Database Key:
1153f36f9100SBarry Smith . -pc_gasm_type [basic,restrict,interpolate,none] - Sets `PCGASM` type
1154b1a0a8a3SJed Brown 
1155b1a0a8a3SJed Brown   Level: intermediate
1156b1a0a8a3SJed Brown 
1157f1580f4eSBarry Smith .seealso: `PCGASM`, `PCGASMSetSubdomains()`, `PCGASMGetSubKSP()`,
1158f1580f4eSBarry Smith           `PCGASMCreateSubdomains2D()`, `PCASM`, `PCASMSetType()`
1159b1a0a8a3SJed Brown @*/
1160d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGASMSetType(PC pc, PCGASMType type)
1161d71ae5a4SJacob Faibussowitsch {
1162b1a0a8a3SJed Brown   PetscFunctionBegin;
1163b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1164b1a0a8a3SJed Brown   PetscValidLogicalCollectiveEnum(pc, type, 2);
1165cac4c232SBarry Smith   PetscTryMethod(pc, "PCGASMSetType_C", (PC, PCGASMType), (pc, type));
11663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1167b1a0a8a3SJed Brown }
1168b1a0a8a3SJed Brown 
1169b1a0a8a3SJed Brown /*@
1170f746d493SDmitry Karpeev   PCGASMSetSortIndices - Determines whether subdomain indices are sorted.
1171b1a0a8a3SJed Brown 
1172c3339decSBarry Smith   Logically Collective
1173b1a0a8a3SJed Brown 
1174b1a0a8a3SJed Brown   Input Parameters:
1175b1a0a8a3SJed Brown + pc     - the preconditioner context
1176b1a0a8a3SJed Brown - doSort - sort the subdomain indices
1177b1a0a8a3SJed Brown 
1178b1a0a8a3SJed Brown   Level: intermediate
1179b1a0a8a3SJed Brown 
1180f1580f4eSBarry Smith .seealso: `PCGASM`, `PCGASMSetSubdomains()`, `PCGASMGetSubKSP()`,
1181db781477SPatrick Sanan           `PCGASMCreateSubdomains2D()`
1182b1a0a8a3SJed Brown @*/
1183d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGASMSetSortIndices(PC pc, PetscBool doSort)
1184d71ae5a4SJacob Faibussowitsch {
1185b1a0a8a3SJed Brown   PetscFunctionBegin;
1186b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1187b1a0a8a3SJed Brown   PetscValidLogicalCollectiveBool(pc, doSort, 2);
1188cac4c232SBarry Smith   PetscTryMethod(pc, "PCGASMSetSortIndices_C", (PC, PetscBool), (pc, doSort));
11893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1190b1a0a8a3SJed Brown }
1191b1a0a8a3SJed Brown 
1192b1a0a8a3SJed Brown /*@C
1193f1580f4eSBarry Smith   PCGASMGetSubKSP - Gets the local `KSP` contexts for all subdomains on
1194f1580f4eSBarry Smith   this MPI rank.
1195b1a0a8a3SJed Brown 
1196c3339decSBarry Smith   Collective iff first_local is requested
1197b1a0a8a3SJed Brown 
1198b1a0a8a3SJed Brown   Input Parameter:
1199b1a0a8a3SJed Brown . pc - the preconditioner context
1200b1a0a8a3SJed Brown 
1201b1a0a8a3SJed Brown   Output Parameters:
1202f36f9100SBarry Smith + n_local     - the number of blocks on this MPI rank or `NULL`
1203f36f9100SBarry Smith . first_local - the global number of the first block on this rank or `NULL`,
1204f36f9100SBarry Smith                  all ranks must request or all must pass `NULL`
1205f1580f4eSBarry Smith - ksp         - the array of `KSP` contexts
1206b1a0a8a3SJed Brown 
1207f36f9100SBarry Smith   Level: advanced
1208f36f9100SBarry Smith 
1209b1a0a8a3SJed Brown   Note:
1210f1580f4eSBarry Smith   After `PCGASMGetSubKSP()` the array of `KSP`es is not to be freed
1211b1a0a8a3SJed Brown 
1212f36f9100SBarry Smith   Currently for some matrix implementations only 1 block per MPI process
1213b1a0a8a3SJed Brown   is supported.
1214b1a0a8a3SJed Brown 
1215f1580f4eSBarry Smith   You must call `KSPSetUp()` before calling `PCGASMGetSubKSP()`.
1216b1a0a8a3SJed Brown 
1217f1580f4eSBarry Smith .seealso: `PCGASM`, `PCGASMSetSubdomains()`, `PCGASMSetOverlap()`,
1218db781477SPatrick Sanan           `PCGASMCreateSubdomains2D()`,
1219b1a0a8a3SJed Brown @*/
1220d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGASMGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[])
1221d71ae5a4SJacob Faibussowitsch {
1222b1a0a8a3SJed Brown   PetscFunctionBegin;
1223b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1224cac4c232SBarry Smith   PetscUseMethod(pc, "PCGASMGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp));
12253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1226b1a0a8a3SJed Brown }
1227b1a0a8a3SJed Brown 
1228b1a0a8a3SJed Brown /*MC
1229f746d493SDmitry Karpeev    PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
1230f1580f4eSBarry Smith            its own `KSP` object on a subset of MPI ranks
1231b1a0a8a3SJed Brown 
1232b1a0a8a3SJed Brown    Options Database Keys:
1233f1580f4eSBarry Smith +  -pc_gasm_total_subdomains <n>  - Sets total number of local subdomains to be distributed among MPI ranks
1234f1580f4eSBarry Smith .  -pc_gasm_view_subdomains       - activates the printing of subdomain indices in `PCView()`, -ksp_view or -snes_view
1235f1580f4eSBarry Smith .  -pc_gasm_print_subdomains      - activates the printing of subdomain indices in `PCSetUp()`
123606b43e7eSDmitry Karpeev .  -pc_gasm_overlap <ovl>         - Sets overlap by which to (automatically) extend local subdomains
1237f1580f4eSBarry Smith -  -pc_gasm_type [basic,restrict,interpolate,none] - Sets `PCGASMType`
1238b1a0a8a3SJed Brown 
1239f36f9100SBarry Smith    Level: beginner
1240f36f9100SBarry Smith 
124195452b02SPatrick Sanan    Notes:
1242f36f9100SBarry Smith    To set options on the solvers for each block append `-sub_` to all the `KSP`, and `PC`
1243f36f9100SBarry Smith    options database keys. For example, `-sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly`
1244b1a0a8a3SJed Brown 
1245f1580f4eSBarry Smith    To set the options on the solvers separate for each block call `PCGASMGetSubKSP()`
1246f1580f4eSBarry Smith    and set the options directly on the resulting `KSP` object (you can access its `PC`
12470462cc06SPierre Jolivet    with `KSPGetPC()`)
1248b1a0a8a3SJed Brown 
1249b1a0a8a3SJed Brown     References:
1250606c0280SSatish Balay +   * - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions
125196a0c994SBarry Smith      Courant Institute, New York University Technical report
1252606c0280SSatish Balay -   * - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
125396a0c994SBarry Smith     Cambridge University Press.
1254b1a0a8a3SJed Brown 
1255f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCASM`, `PCGASMType`, `PCGASMSetType()`,
1256db781477SPatrick Sanan           `PCBJACOBI`, `PCGASMGetSubKSP()`, `PCGASMSetSubdomains()`,
1257db781477SPatrick Sanan           `PCSetModifySubMatrices()`, `PCGASMSetOverlap()`, `PCGASMSetType()`
1258b1a0a8a3SJed Brown M*/
1259b1a0a8a3SJed Brown 
1260d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc)
1261d71ae5a4SJacob Faibussowitsch {
1262f746d493SDmitry Karpeev   PC_GASM *osm;
1263b1a0a8a3SJed Brown 
1264b1a0a8a3SJed Brown   PetscFunctionBegin;
12654dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&osm));
12662fa5cd67SKarl Rupp 
12678f3b4b4dSDmitry Karpeev   osm->N                        = PETSC_DETERMINE;
126806b43e7eSDmitry Karpeev   osm->n                        = PETSC_DECIDE;
12698f3b4b4dSDmitry Karpeev   osm->nmax                     = PETSC_DETERMINE;
12708f3b4b4dSDmitry Karpeev   osm->overlap                  = 0;
12710a545947SLisandro Dalcin   osm->ksp                      = NULL;
12720a545947SLisandro Dalcin   osm->gorestriction            = NULL;
12730a545947SLisandro Dalcin   osm->girestriction            = NULL;
12740a545947SLisandro Dalcin   osm->pctoouter                = NULL;
12750a545947SLisandro Dalcin   osm->gx                       = NULL;
12760a545947SLisandro Dalcin   osm->gy                       = NULL;
12770a545947SLisandro Dalcin   osm->x                        = NULL;
12780a545947SLisandro Dalcin   osm->y                        = NULL;
12790a545947SLisandro Dalcin   osm->pcx                      = NULL;
12800a545947SLisandro Dalcin   osm->pcy                      = NULL;
12810a545947SLisandro Dalcin   osm->permutationIS            = NULL;
12820a545947SLisandro Dalcin   osm->permutationP             = NULL;
12830a545947SLisandro Dalcin   osm->pcmat                    = NULL;
12840a545947SLisandro Dalcin   osm->ois                      = NULL;
12850a545947SLisandro Dalcin   osm->iis                      = NULL;
12860a545947SLisandro Dalcin   osm->pmat                     = NULL;
1287f746d493SDmitry Karpeev   osm->type                     = PC_GASM_RESTRICT;
12886a4f0f73SDmitry Karpeev   osm->same_subdomain_solvers   = PETSC_TRUE;
1289b1a0a8a3SJed Brown   osm->sort_indices             = PETSC_TRUE;
1290d709ea83SDmitry Karpeev   osm->dm_subdomains            = PETSC_FALSE;
1291ea91fabdSFande Kong   osm->hierarchicalpartitioning = PETSC_FALSE;
1292b1a0a8a3SJed Brown 
1293b1a0a8a3SJed Brown   pc->data                 = (void *)osm;
1294f746d493SDmitry Karpeev   pc->ops->apply           = PCApply_GASM;
129548e38a8aSPierre Jolivet   pc->ops->matapply        = PCMatApply_GASM;
1296f746d493SDmitry Karpeev   pc->ops->applytranspose  = PCApplyTranspose_GASM;
1297f746d493SDmitry Karpeev   pc->ops->setup           = PCSetUp_GASM;
1298a06653b4SBarry Smith   pc->ops->reset           = PCReset_GASM;
1299f746d493SDmitry Karpeev   pc->ops->destroy         = PCDestroy_GASM;
1300f746d493SDmitry Karpeev   pc->ops->setfromoptions  = PCSetFromOptions_GASM;
1301f746d493SDmitry Karpeev   pc->ops->setuponblocks   = PCSetUpOnBlocks_GASM;
1302f746d493SDmitry Karpeev   pc->ops->view            = PCView_GASM;
13030a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
1304b1a0a8a3SJed Brown 
13059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMSetSubdomains_C", PCGASMSetSubdomains_GASM));
13069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMSetOverlap_C", PCGASMSetOverlap_GASM));
13079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMSetType_C", PCGASMSetType_GASM));
13089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMSetSortIndices_C", PCGASMSetSortIndices_GASM));
13099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGASMGetSubKSP_C", PCGASMGetSubKSP_GASM));
13103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1311b1a0a8a3SJed Brown }
1312b1a0a8a3SJed Brown 
1313d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[])
1314d71ae5a4SJacob Faibussowitsch {
1315b1a0a8a3SJed Brown   MatPartitioning mpart;
1316b1a0a8a3SJed Brown   const char     *prefix;
1317b1a0a8a3SJed Brown   PetscInt        i, j, rstart, rend, bs;
1318976e8c5aSLisandro Dalcin   PetscBool       hasop, isbaij = PETSC_FALSE, foundpart = PETSC_FALSE;
13190298fd71SBarry Smith   Mat             Ad = NULL, adj;
1320b1a0a8a3SJed Brown   IS              ispart, isnumb, *is;
1321b1a0a8a3SJed Brown 
1322b1a0a8a3SJed Brown   PetscFunctionBegin;
132363a3b9bcSJacob Faibussowitsch   PetscCheck(nloc >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "number of local subdomains must > 0, got nloc = %" PetscInt_FMT, nloc);
1324b1a0a8a3SJed Brown 
1325b1a0a8a3SJed Brown   /* Get prefix, row distribution, and block size */
13269566063dSJacob Faibussowitsch   PetscCall(MatGetOptionsPrefix(A, &prefix));
13279566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
13289566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
13292472a847SBarry Smith   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);
1330b1a0a8a3SJed Brown 
1331b1a0a8a3SJed Brown   /* Get diagonal block from matrix if possible */
13329566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(A, MATOP_GET_DIAGONAL_BLOCK, &hasop));
133348a46eb9SPierre Jolivet   if (hasop) PetscCall(MatGetDiagonalBlock(A, &Ad));
1334b1a0a8a3SJed Brown   if (Ad) {
13359566063dSJacob Faibussowitsch     PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad, MATSEQBAIJ, &isbaij));
13369566063dSJacob Faibussowitsch     if (!isbaij) PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad, MATSEQSBAIJ, &isbaij));
1337b1a0a8a3SJed Brown   }
13388f3b4b4dSDmitry Karpeev   if (Ad && nloc > 1) {
1339b1a0a8a3SJed Brown     PetscBool match, done;
1340b1a0a8a3SJed Brown     /* Try to setup a good matrix partitioning if available */
13419566063dSJacob Faibussowitsch     PetscCall(MatPartitioningCreate(PETSC_COMM_SELF, &mpart));
13429566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart, prefix));
13439566063dSJacob Faibussowitsch     PetscCall(MatPartitioningSetFromOptions(mpart));
13449566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)mpart, MATPARTITIONINGCURRENT, &match));
134548a46eb9SPierre Jolivet     if (!match) PetscCall(PetscObjectTypeCompare((PetscObject)mpart, MATPARTITIONINGSQUARE, &match));
1346b1a0a8a3SJed Brown     if (!match) { /* assume a "good" partitioner is available */
13471a83f524SJed Brown       PetscInt        na;
13481a83f524SJed Brown       const PetscInt *ia, *ja;
13499566063dSJacob Faibussowitsch       PetscCall(MatGetRowIJ(Ad, 0, PETSC_TRUE, isbaij, &na, &ia, &ja, &done));
1350b1a0a8a3SJed Brown       if (done) {
1351b1a0a8a3SJed Brown         /* Build adjacency matrix by hand. Unfortunately a call to
1352b1a0a8a3SJed Brown            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1353b1a0a8a3SJed Brown            remove the block-aij structure and we cannot expect
1354b1a0a8a3SJed Brown            MatPartitioning to split vertices as we need */
13550a545947SLisandro Dalcin         PetscInt        i, j, len, nnz, cnt, *iia = NULL, *jja = NULL;
13561a83f524SJed Brown         const PetscInt *row;
1357b1a0a8a3SJed Brown         nnz = 0;
1358b1a0a8a3SJed Brown         for (i = 0; i < na; i++) { /* count number of nonzeros */
1359b1a0a8a3SJed Brown           len = ia[i + 1] - ia[i];
1360b1a0a8a3SJed Brown           row = ja + ia[i];
1361b1a0a8a3SJed Brown           for (j = 0; j < len; j++) {
1362b1a0a8a3SJed Brown             if (row[j] == i) { /* don't count diagonal */
13639371c9d4SSatish Balay               len--;
13649371c9d4SSatish Balay               break;
1365b1a0a8a3SJed Brown             }
1366b1a0a8a3SJed Brown           }
1367b1a0a8a3SJed Brown           nnz += len;
1368b1a0a8a3SJed Brown         }
13699566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(na + 1, &iia));
13709566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(nnz, &jja));
1371b1a0a8a3SJed Brown         nnz    = 0;
1372b1a0a8a3SJed Brown         iia[0] = 0;
1373b1a0a8a3SJed Brown         for (i = 0; i < na; i++) { /* fill adjacency */
1374b1a0a8a3SJed Brown           cnt = 0;
1375b1a0a8a3SJed Brown           len = ia[i + 1] - ia[i];
1376b1a0a8a3SJed Brown           row = ja + ia[i];
1377b1a0a8a3SJed Brown           for (j = 0; j < len; j++) {
13782fa5cd67SKarl Rupp             if (row[j] != i) jja[nnz + cnt++] = row[j]; /* if not diagonal */
1379b1a0a8a3SJed Brown           }
1380b1a0a8a3SJed Brown           nnz += cnt;
1381b1a0a8a3SJed Brown           iia[i + 1] = nnz;
1382b1a0a8a3SJed Brown         }
1383b1a0a8a3SJed Brown         /* Partitioning of the adjacency matrix */
13849566063dSJacob Faibussowitsch         PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, na, na, iia, jja, NULL, &adj));
13859566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetAdjacency(mpart, adj));
13869566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetNParts(mpart, nloc));
13879566063dSJacob Faibussowitsch         PetscCall(MatPartitioningApply(mpart, &ispart));
13889566063dSJacob Faibussowitsch         PetscCall(ISPartitioningToNumbering(ispart, &isnumb));
13899566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&adj));
1390b1a0a8a3SJed Brown         foundpart = PETSC_TRUE;
1391b1a0a8a3SJed Brown       }
13929566063dSJacob Faibussowitsch       PetscCall(MatRestoreRowIJ(Ad, 0, PETSC_TRUE, isbaij, &na, &ia, &ja, &done));
1393b1a0a8a3SJed Brown     }
13949566063dSJacob Faibussowitsch     PetscCall(MatPartitioningDestroy(&mpart));
1395b1a0a8a3SJed Brown   }
13969566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nloc, &is));
1397b1a0a8a3SJed Brown   if (!foundpart) {
1398b1a0a8a3SJed Brown     /* Partitioning by contiguous chunks of rows */
1399b1a0a8a3SJed Brown 
1400b1a0a8a3SJed Brown     PetscInt mbs   = (rend - rstart) / bs;
1401b1a0a8a3SJed Brown     PetscInt start = rstart;
14028f3b4b4dSDmitry Karpeev     for (i = 0; i < nloc; i++) {
14038f3b4b4dSDmitry Karpeev       PetscInt count = (mbs / nloc + ((mbs % nloc) > i)) * bs;
14049566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, count, start, 1, &is[i]));
1405b1a0a8a3SJed Brown       start += count;
1406b1a0a8a3SJed Brown     }
1407b1a0a8a3SJed Brown 
1408b1a0a8a3SJed Brown   } else {
1409b1a0a8a3SJed Brown     /* Partitioning by adjacency of diagonal block  */
1410b1a0a8a3SJed Brown 
1411b1a0a8a3SJed Brown     const PetscInt *numbering;
1412b1a0a8a3SJed Brown     PetscInt       *count, nidx, *indices, *newidx, start = 0;
1413b1a0a8a3SJed Brown     /* Get node count in each partition */
14149566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &count));
14159566063dSJacob Faibussowitsch     PetscCall(ISPartitioningCount(ispart, nloc, count));
1416b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
14178f3b4b4dSDmitry Karpeev       for (i = 0; i < nloc; i++) count[i] *= bs;
1418b1a0a8a3SJed Brown     }
1419b1a0a8a3SJed Brown     /* Build indices from node numbering */
14209566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(isnumb, &nidx));
14219566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nidx, &indices));
1422b1a0a8a3SJed Brown     for (i = 0; i < nidx; i++) indices[i] = i; /* needs to be initialized */
14239566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(isnumb, &numbering));
14249566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithPermutation(nidx, numbering, indices));
14259566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(isnumb, &numbering));
1426b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
14279566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nidx * bs, &newidx));
14282fa5cd67SKarl Rupp       for (i = 0; i < nidx; i++) {
14292fa5cd67SKarl Rupp         for (j = 0; j < bs; j++) newidx[i * bs + j] = indices[i] * bs + j;
14302fa5cd67SKarl Rupp       }
14319566063dSJacob Faibussowitsch       PetscCall(PetscFree(indices));
1432b1a0a8a3SJed Brown       nidx *= bs;
1433b1a0a8a3SJed Brown       indices = newidx;
1434b1a0a8a3SJed Brown     }
1435b1a0a8a3SJed Brown     /* Shift to get global indices */
1436b1a0a8a3SJed Brown     for (i = 0; i < nidx; i++) indices[i] += rstart;
1437b1a0a8a3SJed Brown 
1438b1a0a8a3SJed Brown     /* Build the index sets for each block */
14398f3b4b4dSDmitry Karpeev     for (i = 0; i < nloc; i++) {
14409566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, count[i], &indices[start], PETSC_COPY_VALUES, &is[i]));
14419566063dSJacob Faibussowitsch       PetscCall(ISSort(is[i]));
1442b1a0a8a3SJed Brown       start += count[i];
1443b1a0a8a3SJed Brown     }
1444b1a0a8a3SJed Brown 
14459566063dSJacob Faibussowitsch     PetscCall(PetscFree(count));
14469566063dSJacob Faibussowitsch     PetscCall(PetscFree(indices));
14479566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isnumb));
14489566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ispart));
1449b1a0a8a3SJed Brown   }
14506a4f0f73SDmitry Karpeev   *iis = is;
14513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14528f3b4b4dSDmitry Karpeev }
14538f3b4b4dSDmitry Karpeev 
1454d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PCGASMCreateStraddlingSubdomains(Mat A, PetscInt N, PetscInt *n, IS *iis[])
1455d71ae5a4SJacob Faibussowitsch {
14568f3b4b4dSDmitry Karpeev   PetscFunctionBegin;
14579566063dSJacob Faibussowitsch   PetscCall(MatSubdomainsCreateCoalesce(A, N, n, iis));
14583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14598f3b4b4dSDmitry Karpeev }
14608f3b4b4dSDmitry Karpeev 
14618f3b4b4dSDmitry Karpeev /*@C
1462f36f9100SBarry Smith   PCGASMCreateSubdomains - Creates `n` index sets defining `n` nonoverlapping subdomains on this MPI process for the `PCGASM` additive
14638f3b4b4dSDmitry Karpeev   Schwarz preconditioner for a any problem based on its matrix.
14648f3b4b4dSDmitry Karpeev 
14658f3b4b4dSDmitry Karpeev   Collective
14668f3b4b4dSDmitry Karpeev 
14678f3b4b4dSDmitry Karpeev   Input Parameters:
14688f3b4b4dSDmitry Karpeev + A - The global matrix operator
14698f3b4b4dSDmitry Karpeev - N - the number of global subdomains requested
14708f3b4b4dSDmitry Karpeev 
14718f3b4b4dSDmitry Karpeev   Output Parameters:
1472f1580f4eSBarry Smith + n   - the number of subdomains created on this MPI rank
14738f3b4b4dSDmitry Karpeev - iis - the array of index sets defining the local inner subdomains (on which the correction is applied)
14748f3b4b4dSDmitry Karpeev 
14758f3b4b4dSDmitry Karpeev   Level: advanced
14768f3b4b4dSDmitry Karpeev 
1477f36f9100SBarry Smith   Notes:
1478f36f9100SBarry Smith   When `N` >= A's communicator size, each subdomain is local -- contained within a single MPI process.
147920f4b53cSBarry Smith   When `N` < size, the subdomains are 'straddling' (rank boundaries) and are no longer local.
148020f4b53cSBarry Smith   The resulting subdomains can be use in `PCGASMSetSubdomains`(pc,n,iss,`NULL`).  The overlapping
14818f3b4b4dSDmitry Karpeev   outer subdomains will be automatically generated from these according to the requested amount of
14828f3b4b4dSDmitry Karpeev   overlap; this is currently supported only with local subdomains.
14838f3b4b4dSDmitry Karpeev 
1484f36f9100SBarry Smith   Use `PCGASMDestroySubdomains()` to free the array and the list of index sets.
1485f36f9100SBarry Smith 
1486f1580f4eSBarry Smith .seealso: `PCGASM`, `PCGASMSetSubdomains()`, `PCGASMDestroySubdomains()`
14878f3b4b4dSDmitry Karpeev @*/
1488d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGASMCreateSubdomains(Mat A, PetscInt N, PetscInt *n, IS *iis[])
1489d71ae5a4SJacob Faibussowitsch {
14908f3b4b4dSDmitry Karpeev   PetscMPIInt size;
14918f3b4b4dSDmitry Karpeev 
14928f3b4b4dSDmitry Karpeev   PetscFunctionBegin;
14938f3b4b4dSDmitry Karpeev   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
14944f572ea9SToby Isaac   PetscAssertPointer(iis, 4);
14958f3b4b4dSDmitry Karpeev 
149663a3b9bcSJacob Faibussowitsch   PetscCheck(N >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of subdomains must be > 0, N = %" PetscInt_FMT, N);
14979566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
14988f3b4b4dSDmitry Karpeev   if (N >= size) {
14998f3b4b4dSDmitry Karpeev     *n = N / size + (N % size);
15009566063dSJacob Faibussowitsch     PetscCall(PCGASMCreateLocalSubdomains(A, *n, iis));
15016a4f0f73SDmitry Karpeev   } else {
15029566063dSJacob Faibussowitsch     PetscCall(PCGASMCreateStraddlingSubdomains(A, N, n, iis));
15036a4f0f73SDmitry Karpeev   }
15043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1505b1a0a8a3SJed Brown }
1506b1a0a8a3SJed Brown 
1507b1a0a8a3SJed Brown /*@C
1508f746d493SDmitry Karpeev   PCGASMDestroySubdomains - Destroys the index sets created with
1509f1580f4eSBarry Smith   `PCGASMCreateSubdomains()` or `PCGASMCreateSubdomains2D()`. Should be
1510f1580f4eSBarry Smith   called after setting subdomains with `PCGASMSetSubdomains()`.
1511b1a0a8a3SJed Brown 
1512b1a0a8a3SJed Brown   Collective
1513b1a0a8a3SJed Brown 
1514b1a0a8a3SJed Brown   Input Parameters:
1515b1a0a8a3SJed Brown + n   - the number of index sets
1516f36f9100SBarry Smith . iis - the array of inner subdomains
1517f36f9100SBarry Smith - ois - the array of outer subdomains, can be `NULL`
1518b1a0a8a3SJed Brown 
15196a4f0f73SDmitry Karpeev   Level: intermediate
15206a4f0f73SDmitry Karpeev 
1521f1580f4eSBarry Smith   Note:
1522f1580f4eSBarry Smith   This is a convenience subroutine that walks each list,
1523f1580f4eSBarry Smith   destroys each `IS` on the list, and then frees the list. At the end the
152420f4b53cSBarry Smith   list pointers are set to `NULL`.
1525b1a0a8a3SJed Brown 
1526feefa0e1SJacob Faibussowitsch .seealso: `PCGASM`, `PCGASMCreateSubdomains()`, `PCGASMSetSubdomains()`
1527b1a0a8a3SJed Brown @*/
1528d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGASMDestroySubdomains(PetscInt n, IS **iis, IS **ois)
1529d71ae5a4SJacob Faibussowitsch {
1530b1a0a8a3SJed Brown   PetscInt i;
15315fd66863SKarl Rupp 
1532b1a0a8a3SJed Brown   PetscFunctionBegin;
15333ba16761SJacob Faibussowitsch   if (n <= 0) PetscFunctionReturn(PETSC_SUCCESS);
15346a4f0f73SDmitry Karpeev   if (ois) {
15354f572ea9SToby Isaac     PetscAssertPointer(ois, 3);
15362c112581SDmitry Karpeev     if (*ois) {
15374f572ea9SToby Isaac       PetscAssertPointer(*ois, 3);
153848a46eb9SPierre Jolivet       for (i = 0; i < n; i++) PetscCall(ISDestroy(&(*ois)[i]));
15399566063dSJacob Faibussowitsch       PetscCall(PetscFree((*ois)));
15402c112581SDmitry Karpeev     }
1541b1a0a8a3SJed Brown   }
1542b9d0fdaaSFande Kong   if (iis) {
15434f572ea9SToby Isaac     PetscAssertPointer(iis, 2);
1544b9d0fdaaSFande Kong     if (*iis) {
15454f572ea9SToby Isaac       PetscAssertPointer(*iis, 2);
154648a46eb9SPierre Jolivet       for (i = 0; i < n; i++) PetscCall(ISDestroy(&(*iis)[i]));
15479566063dSJacob Faibussowitsch       PetscCall(PetscFree((*iis)));
1548b9d0fdaaSFande Kong     }
1549b9d0fdaaSFande Kong   }
15503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1551b1a0a8a3SJed Brown }
1552b1a0a8a3SJed Brown 
1553af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M, N, xleft, ylow, xright, yhigh, first, last, xleft_loc, ylow_loc, xright_loc, yhigh_loc, n) \
1554a8f51744SPierre Jolivet   do { \
1555af538404SDmitry Karpeev     PetscInt first_row = first / M, last_row = last / M + 1; \
1556af538404SDmitry Karpeev     /*                                                                                                    \
1557af538404SDmitry Karpeev      Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners       \
1558af538404SDmitry Karpeev      of the bounding box of the intersection of the subdomain with the local ownership range (local       \
1559af538404SDmitry Karpeev      subdomain).                                                                                          \
1560af538404SDmitry Karpeev      Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows       \
1561af538404SDmitry Karpeev      of the intersection.                                                                                 \
1562af538404SDmitry Karpeev     */ \
1563af538404SDmitry Karpeev     /* ylow_loc is the grid row containing the first element of the local sumbdomain */ \
1564eec7e3faSDmitry Karpeev     *ylow_loc = PetscMax(first_row, ylow); \
1565af538404SDmitry Karpeev     /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \
1566eec7e3faSDmitry Karpeev     *xleft_loc = *ylow_loc == first_row ? PetscMax(first % M, xleft) : xleft; \
1567af538404SDmitry Karpeev     /* yhigh_loc is the grid row above the last local subdomain element */ \
1568eec7e3faSDmitry Karpeev     *yhigh_loc = PetscMin(last_row, yhigh); \
1569af538404SDmitry Karpeev     /* xright is the offset of the end of the  local subdomain within its grid row (might actually be outside the local subdomain) */ \
1570eec7e3faSDmitry Karpeev     *xright_loc = *yhigh_loc == last_row ? PetscMin(xright, last % M) : xright; \
1571af538404SDmitry Karpeev     /* Now compute the size of the local subdomain n. */ \
1572c3518bceSDmitry Karpeev     *n = 0; \
1573eec7e3faSDmitry Karpeev     if (*ylow_loc < *yhigh_loc) { \
1574af538404SDmitry Karpeev       PetscInt width = xright - xleft; \
1575eec7e3faSDmitry Karpeev       *n += width * (*yhigh_loc - *ylow_loc - 1); \
1576eec7e3faSDmitry Karpeev       *n += PetscMin(PetscMax(*xright_loc - xleft, 0), width); \
1577eec7e3faSDmitry Karpeev       *n -= PetscMin(PetscMax(*xleft_loc - xleft, 0), width); \
1578af538404SDmitry Karpeev     } \
1579a8f51744SPierre Jolivet   } while (0)
1580af538404SDmitry Karpeev 
1581f36f9100SBarry Smith /*@C
1582f1580f4eSBarry Smith   PCGASMCreateSubdomains2D - Creates the index sets for the `PCGASM` overlapping Schwarz
1583b1a0a8a3SJed Brown   preconditioner for a two-dimensional problem on a regular grid.
1584b1a0a8a3SJed Brown 
1585af538404SDmitry Karpeev   Collective
1586b1a0a8a3SJed Brown 
1587b1a0a8a3SJed Brown   Input Parameters:
15886b867d5aSJose E. Roman + pc       - the preconditioner context
15896b867d5aSJose E. Roman . M        - the global number of grid points in the x direction
15906b867d5aSJose E. Roman . N        - the global number of grid points in the y direction
15916b867d5aSJose E. Roman . Mdomains - the global number of subdomains in the x direction
15926b867d5aSJose E. Roman . Ndomains - the global number of subdomains in the y direction
1593b1a0a8a3SJed Brown . dof      - degrees of freedom per node
1594b1a0a8a3SJed Brown - overlap  - overlap in mesh lines
1595b1a0a8a3SJed Brown 
1596b1a0a8a3SJed Brown   Output Parameters:
1597feefa0e1SJacob Faibussowitsch + nsub - the number of local subdomains created
15986a4f0f73SDmitry Karpeev . iis  - array of index sets defining inner (nonoverlapping) subdomains
15996a4f0f73SDmitry Karpeev - ois  - array of index sets defining outer (overlapping, if overlap > 0) subdomains
1600b1a0a8a3SJed Brown 
1601b1a0a8a3SJed Brown   Level: advanced
1602b1a0a8a3SJed Brown 
1603f36f9100SBarry Smith   Note:
1604f36f9100SBarry Smith   Use `PCGASMDestroySubdomains()` to free the index sets and the arrays
1605f36f9100SBarry Smith 
1606feefa0e1SJacob Faibussowitsch   Fortran Notes:
1607f36f9100SBarry Smith   The `IS` must be declared as an array of length long enough to hold `Nsub` entries
1608f36f9100SBarry Smith 
1609f36f9100SBarry Smith .seealso: `PCGASM`, `PCGASMSetSubdomains()`, `PCGASMGetSubKSP()`, `PCGASMSetOverlap()`, `PCASMCreateSubdomains2D()`,
1610f36f9100SBarry Smith           `PCGASMDestroySubdomains()`
1611b1a0a8a3SJed Brown @*/
1612d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGASMCreateSubdomains2D(PC pc, PetscInt M, PetscInt N, PetscInt Mdomains, PetscInt Ndomains, PetscInt dof, PetscInt overlap, PetscInt *nsub, IS **iis, IS **ois)
1613d71ae5a4SJacob Faibussowitsch {
16142bbb417fSDmitry Karpeev   PetscMPIInt size, rank;
16152bbb417fSDmitry Karpeev   PetscInt    i, j;
16162bbb417fSDmitry Karpeev   PetscInt    maxheight, maxwidth;
16172bbb417fSDmitry Karpeev   PetscInt    xstart, xleft, xright, xleft_loc, xright_loc;
16182bbb417fSDmitry Karpeev   PetscInt    ystart, ylow, yhigh, ylow_loc, yhigh_loc;
1619eec7e3faSDmitry Karpeev   PetscInt    x[2][2], y[2][2], n[2];
16202bbb417fSDmitry Karpeev   PetscInt    first, last;
16212bbb417fSDmitry Karpeev   PetscInt    nidx, *idx;
16222bbb417fSDmitry Karpeev   PetscInt    ii, jj, s, q, d;
1623af538404SDmitry Karpeev   PetscInt    k, kk;
16242bbb417fSDmitry Karpeev   PetscMPIInt color;
16252bbb417fSDmitry Karpeev   MPI_Comm    comm, subcomm;
16260a545947SLisandro Dalcin   IS        **xis = NULL, **is = ois, **is_local = iis;
1627d34fcf5fSBarry Smith 
1628b1a0a8a3SJed Brown   PetscFunctionBegin;
16299566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
16309566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
16319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
16329566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(pc->pmat, &first, &last));
16339371c9d4SSatish Balay   PetscCheck((first % dof) == 0 && (last % dof) == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE,
16349371c9d4SSatish Balay              "Matrix row partitioning unsuitable for domain decomposition: local row range (%" PetscInt_FMT ",%" PetscInt_FMT ") "
16359371c9d4SSatish Balay              "does not respect the number of degrees of freedom per grid point %" PetscInt_FMT,
16369371c9d4SSatish Balay              first, last, dof);
1637d34fcf5fSBarry Smith 
1638af538404SDmitry Karpeev   /* Determine the number of domains with nonzero intersections with the local ownership range. */
16392bbb417fSDmitry Karpeev   s      = 0;
1640af538404SDmitry Karpeev   ystart = 0;
1641af538404SDmitry Karpeev   for (j = 0; j < Ndomains; ++j) {
1642af538404SDmitry Karpeev     maxheight = N / Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1643da81f932SPierre Jolivet     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);
1644eec7e3faSDmitry Karpeev     /* Vertical domain limits with an overlap. */
1645eec7e3faSDmitry Karpeev     ylow   = PetscMax(ystart - overlap, 0);
1646eec7e3faSDmitry Karpeev     yhigh  = PetscMin(ystart + maxheight + overlap, N);
1647b1a0a8a3SJed Brown     xstart = 0;
1648af538404SDmitry Karpeev     for (i = 0; i < Mdomains; ++i) {
1649af538404SDmitry Karpeev       maxwidth = M / Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
165063a3b9bcSJacob Faibussowitsch       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);
1651eec7e3faSDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1652eec7e3faSDmitry Karpeev       xleft  = PetscMax(xstart - overlap, 0);
1653eec7e3faSDmitry Karpeev       xright = PetscMin(xstart + maxwidth + overlap, M);
1654af538404SDmitry Karpeev       /*
1655f1580f4eSBarry Smith          Determine whether this subdomain intersects this rank's ownership range of pc->pmat.
1656af538404SDmitry Karpeev       */
1657c3518bceSDmitry Karpeev       PCGASMLocalSubdomainBounds2D(M, N, xleft, ylow, xright, yhigh, first, last, (&xleft_loc), (&ylow_loc), (&xright_loc), (&yhigh_loc), (&nidx));
16582fa5cd67SKarl Rupp       if (nidx) ++s;
1659af538404SDmitry Karpeev       xstart += maxwidth;
1660af538404SDmitry Karpeev     } /* for (i = 0; i < Mdomains; ++i) */
1661af538404SDmitry Karpeev     ystart += maxheight;
1662af538404SDmitry Karpeev   } /* for (j = 0; j < Ndomains; ++j) */
16632fa5cd67SKarl Rupp 
1664af538404SDmitry Karpeev   /* Now we can allocate the necessary number of ISs. */
1665af538404SDmitry Karpeev   *nsub = s;
16669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(*nsub, is));
16679566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(*nsub, is_local));
1668af538404SDmitry Karpeev   s      = 0;
1669af538404SDmitry Karpeev   ystart = 0;
1670af538404SDmitry Karpeev   for (j = 0; j < Ndomains; ++j) {
1671af538404SDmitry Karpeev     maxheight = N / Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1672da81f932SPierre Jolivet     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);
1673af538404SDmitry Karpeev     /* Vertical domain limits with an overlap. */
1674eec7e3faSDmitry Karpeev     y[0][0] = PetscMax(ystart - overlap, 0);
1675eec7e3faSDmitry Karpeev     y[0][1] = PetscMin(ystart + maxheight + overlap, N);
1676eec7e3faSDmitry Karpeev     /* Vertical domain limits without an overlap. */
1677eec7e3faSDmitry Karpeev     y[1][0] = ystart;
1678eec7e3faSDmitry Karpeev     y[1][1] = PetscMin(ystart + maxheight, N);
1679eec7e3faSDmitry Karpeev     xstart  = 0;
1680af538404SDmitry Karpeev     for (i = 0; i < Mdomains; ++i) {
1681af538404SDmitry Karpeev       maxwidth = M / Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
168263a3b9bcSJacob Faibussowitsch       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);
1683af538404SDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1684eec7e3faSDmitry Karpeev       x[0][0] = PetscMax(xstart - overlap, 0);
1685eec7e3faSDmitry Karpeev       x[0][1] = PetscMin(xstart + maxwidth + overlap, M);
1686eec7e3faSDmitry Karpeev       /* Horizontal domain limits without an overlap. */
1687eec7e3faSDmitry Karpeev       x[1][0] = xstart;
1688eec7e3faSDmitry Karpeev       x[1][1] = PetscMin(xstart + maxwidth, M);
16892bbb417fSDmitry Karpeev       /*
1690f1580f4eSBarry Smith          Determine whether this domain intersects this rank's ownership range of pc->pmat.
16912bbb417fSDmitry Karpeev          Do this twice: first for the domains with overlaps, and once without.
16922bbb417fSDmitry Karpeev          During the first pass create the subcommunicators, and use them on the second pass as well.
16932bbb417fSDmitry Karpeev       */
16942bbb417fSDmitry Karpeev       for (q = 0; q < 2; ++q) {
1695cc96b2e5SBarry Smith         PetscBool split = PETSC_FALSE;
16962bbb417fSDmitry Karpeev         /*
16972bbb417fSDmitry Karpeev           domain limits, (xleft, xright) and (ylow, yheigh) are adjusted
16982bbb417fSDmitry Karpeev           according to whether the domain with an overlap or without is considered.
16992bbb417fSDmitry Karpeev         */
17009371c9d4SSatish Balay         xleft  = x[q][0];
17019371c9d4SSatish Balay         xright = x[q][1];
17029371c9d4SSatish Balay         ylow   = y[q][0];
17039371c9d4SSatish Balay         yhigh  = y[q][1];
1704c3518bceSDmitry Karpeev         PCGASMLocalSubdomainBounds2D(M, N, xleft, ylow, xright, yhigh, first, last, (&xleft_loc), (&ylow_loc), (&xright_loc), (&yhigh_loc), (&nidx));
1705af538404SDmitry Karpeev         nidx *= dof;
1706eec7e3faSDmitry Karpeev         n[q] = nidx;
17072bbb417fSDmitry Karpeev         /*
1708eec7e3faSDmitry Karpeev          Based on the counted number of indices in the local domain *with an overlap*,
1709f1580f4eSBarry Smith          construct a subcommunicator of all the MPI ranks supporting this domain.
1710eec7e3faSDmitry Karpeev          Observe that a domain with an overlap might have nontrivial local support,
1711eec7e3faSDmitry Karpeev          while the domain without an overlap might not.  Hence, the decision to participate
1712eec7e3faSDmitry Karpeev          in the subcommunicator must be based on the domain with an overlap.
17132bbb417fSDmitry Karpeev          */
17142bbb417fSDmitry Karpeev         if (q == 0) {
17152fa5cd67SKarl Rupp           if (nidx) color = 1;
17162fa5cd67SKarl Rupp           else color = MPI_UNDEFINED;
17179566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Comm_split(comm, color, rank, &subcomm));
1718cc96b2e5SBarry Smith           split = PETSC_TRUE;
17192bbb417fSDmitry Karpeev         }
1720af538404SDmitry Karpeev         /*
1721eec7e3faSDmitry Karpeev          Proceed only if the number of local indices *with an overlap* is nonzero.
1722af538404SDmitry Karpeev          */
1723eec7e3faSDmitry Karpeev         if (n[0]) {
17242fa5cd67SKarl Rupp           if (q == 0) xis = is;
1725af538404SDmitry Karpeev           if (q == 1) {
1726af538404SDmitry Karpeev             /*
1727eec7e3faSDmitry Karpeev              The IS for the no-overlap subdomain shares a communicator with the overlapping domain.
1728af538404SDmitry Karpeev              Moreover, if the overlap is zero, the two ISs are identical.
1729af538404SDmitry Karpeev              */
1730b1a0a8a3SJed Brown             if (overlap == 0) {
1731eec7e3faSDmitry Karpeev               (*is_local)[s] = (*is)[s];
17329566063dSJacob Faibussowitsch               PetscCall(PetscObjectReference((PetscObject)(*is)[s]));
17332bbb417fSDmitry Karpeev               continue;
1734d34fcf5fSBarry Smith             } else {
17356a4f0f73SDmitry Karpeev               xis     = is_local;
1736eec7e3faSDmitry Karpeev               subcomm = ((PetscObject)(*is)[s])->comm;
17372bbb417fSDmitry Karpeev             }
1738af538404SDmitry Karpeev           } /* if (q == 1) */
17390298fd71SBarry Smith           idx = NULL;
17409566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(nidx, &idx));
1741eec7e3faSDmitry Karpeev           if (nidx) {
17422bbb417fSDmitry Karpeev             k = 0;
17432bbb417fSDmitry Karpeev             for (jj = ylow_loc; jj < yhigh_loc; ++jj) {
1744af538404SDmitry Karpeev               PetscInt x0 = (jj == ylow_loc) ? xleft_loc : xleft;
1745af538404SDmitry Karpeev               PetscInt x1 = (jj == yhigh_loc - 1) ? xright_loc : xright;
1746af538404SDmitry Karpeev               kk          = dof * (M * jj + x0);
17472bbb417fSDmitry Karpeev               for (ii = x0; ii < x1; ++ii) {
1748ad540459SPierre Jolivet                 for (d = 0; d < dof; ++d) idx[k++] = kk++;
1749b1a0a8a3SJed Brown               }
1750b1a0a8a3SJed Brown             }
1751eec7e3faSDmitry Karpeev           }
17529566063dSJacob Faibussowitsch           PetscCall(ISCreateGeneral(subcomm, nidx, idx, PETSC_OWN_POINTER, (*xis) + s));
175348a46eb9SPierre Jolivet           if (split) PetscCallMPI(MPI_Comm_free(&subcomm));
1754eec7e3faSDmitry Karpeev         } /* if (n[0]) */
17552bbb417fSDmitry Karpeev       }   /* for (q = 0; q < 2; ++q) */
17562fa5cd67SKarl Rupp       if (n[0]) ++s;
1757af538404SDmitry Karpeev       xstart += maxwidth;
1758af538404SDmitry Karpeev     } /* for (i = 0; i < Mdomains; ++i) */
17592bbb417fSDmitry Karpeev     ystart += maxheight;
1760af538404SDmitry Karpeev   } /* for (j = 0; j < Ndomains; ++j) */
17613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1762b1a0a8a3SJed Brown }
1763b1a0a8a3SJed Brown 
1764b1a0a8a3SJed Brown /*@C
1765f1580f4eSBarry Smith   PCGASMGetSubdomains - Gets the subdomains supported on this MPI rank
1766f1580f4eSBarry Smith   for the `PCGASM` additive Schwarz preconditioner.
1767b1a0a8a3SJed Brown 
1768b1a0a8a3SJed Brown   Not Collective
1769b1a0a8a3SJed Brown 
1770b1a0a8a3SJed Brown   Input Parameter:
1771b1a0a8a3SJed Brown . pc - the preconditioner context
1772b1a0a8a3SJed Brown 
1773b1a0a8a3SJed Brown   Output Parameters:
1774f1580f4eSBarry Smith + n   - the number of subdomains for this MPI rank (default value = 1)
177520f4b53cSBarry Smith . iis - the index sets that define the inner subdomains (without overlap) supported on this rank (can be `NULL`)
177620f4b53cSBarry Smith - ois - the index sets that define the outer subdomains (with overlap) supported on this rank (can be `NULL`)
177720f4b53cSBarry Smith 
177820f4b53cSBarry Smith   Level: advanced
1779b1a0a8a3SJed Brown 
1780f36f9100SBarry Smith   Notes:
1781f36f9100SBarry Smith   The user is responsible for destroying the `IS`s and freeing the returned arrays, this can be done with
1782f36f9100SBarry Smith   `PCGASMDestroySubdomains()`
1783f36f9100SBarry Smith 
1784f1580f4eSBarry Smith   The `IS` numbering is in the parallel, global numbering of the vector.
1785b1a0a8a3SJed Brown 
1786f1580f4eSBarry Smith .seealso: `PCGASM`, `PCGASMSetOverlap()`, `PCGASMGetSubKSP()`, `PCGASMCreateSubdomains2D()`,
1787f36f9100SBarry Smith           `PCGASMSetSubdomains()`, `PCGASMGetSubmatrices()`, `PCGASMDestroySubdomains()`
1788b1a0a8a3SJed Brown @*/
1789d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGASMGetSubdomains(PC pc, PetscInt *n, IS *iis[], IS *ois[])
1790d71ae5a4SJacob Faibussowitsch {
1791f746d493SDmitry Karpeev   PC_GASM  *osm;
1792b1a0a8a3SJed Brown   PetscBool match;
17936a4f0f73SDmitry Karpeev   PetscInt  i;
17945fd66863SKarl Rupp 
1795b1a0a8a3SJed Brown   PetscFunctionBegin;
1796b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
17979566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCGASM, &match));
179828b400f6SJacob Faibussowitsch   PetscCheck(match, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1799f746d493SDmitry Karpeev   osm = (PC_GASM *)pc->data;
1800ab3e923fSDmitry Karpeev   if (n) *n = osm->n;
18011baa6e33SBarry Smith   if (iis) PetscCall(PetscMalloc1(osm->n, iis));
18021baa6e33SBarry Smith   if (ois) PetscCall(PetscMalloc1(osm->n, ois));
18036a4f0f73SDmitry Karpeev   if (iis || ois) {
18046a4f0f73SDmitry Karpeev     for (i = 0; i < osm->n; ++i) {
18056a4f0f73SDmitry Karpeev       if (iis) (*iis)[i] = osm->iis[i];
18066a4f0f73SDmitry Karpeev       if (ois) (*ois)[i] = osm->ois[i];
18076a4f0f73SDmitry Karpeev     }
1808b1a0a8a3SJed Brown   }
18093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1810b1a0a8a3SJed Brown }
1811b1a0a8a3SJed Brown 
1812b1a0a8a3SJed Brown /*@C
1813f1580f4eSBarry Smith   PCGASMGetSubmatrices - Gets the local submatrices (for this MPI rank
1814f1580f4eSBarry Smith   only) for the `PCGASM` additive Schwarz preconditioner.
1815b1a0a8a3SJed Brown 
1816b1a0a8a3SJed Brown   Not Collective
1817b1a0a8a3SJed Brown 
1818b1a0a8a3SJed Brown   Input Parameter:
1819b1a0a8a3SJed Brown . pc - the preconditioner context
1820b1a0a8a3SJed Brown 
1821b1a0a8a3SJed Brown   Output Parameters:
1822f1580f4eSBarry Smith + n   - the number of matrices for this MPI rank (default value = 1)
1823b1a0a8a3SJed Brown - mat - the matrices
1824b1a0a8a3SJed Brown 
1825b1a0a8a3SJed Brown   Level: advanced
1826b1a0a8a3SJed Brown 
1827f36f9100SBarry Smith   Note:
1828f36f9100SBarry Smith   Matrices returned by this routine have the same communicators as the index sets (`IS`)
1829f36f9100SBarry Smith   used to define subdomains in `PCGASMSetSubdomains()`
1830f36f9100SBarry Smith 
1831f1580f4eSBarry Smith .seealso: `PCGASM`, `PCGASMSetOverlap()`, `PCGASMGetSubKSP()`,
1832db781477SPatrick Sanan           `PCGASMCreateSubdomains2D()`, `PCGASMSetSubdomains()`, `PCGASMGetSubdomains()`
1833b1a0a8a3SJed Brown @*/
1834d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGASMGetSubmatrices(PC pc, PetscInt *n, Mat *mat[])
1835d71ae5a4SJacob Faibussowitsch {
1836f746d493SDmitry Karpeev   PC_GASM  *osm;
1837b1a0a8a3SJed Brown   PetscBool match;
1838b1a0a8a3SJed Brown 
1839b1a0a8a3SJed Brown   PetscFunctionBegin;
1840b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
18414f572ea9SToby Isaac   PetscAssertPointer(n, 2);
18424f572ea9SToby Isaac   if (mat) PetscAssertPointer(mat, 3);
184328b400f6SJacob Faibussowitsch   PetscCheck(pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call after KSPSetUp() or PCSetUp().");
18449566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCGASM, &match));
184528b400f6SJacob Faibussowitsch   PetscCheck(match, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1846f746d493SDmitry Karpeev   osm = (PC_GASM *)pc->data;
1847ab3e923fSDmitry Karpeev   if (n) *n = osm->n;
1848b1a0a8a3SJed Brown   if (mat) *mat = osm->pmat;
18493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1850b1a0a8a3SJed Brown }
1851d709ea83SDmitry Karpeev 
1852d709ea83SDmitry Karpeev /*@
1853f1580f4eSBarry Smith   PCGASMSetUseDMSubdomains - Indicates whether to use `DMCreateDomainDecomposition()` to define the subdomains, whenever possible for `PCGASM`
1854f1580f4eSBarry Smith 
1855d709ea83SDmitry Karpeev   Logically Collective
1856d709ea83SDmitry Karpeev 
1857d8d19677SJose E. Roman   Input Parameters:
1858d709ea83SDmitry Karpeev + pc  - the preconditioner
1859f1580f4eSBarry Smith - flg - boolean indicating whether to use subdomains defined by the `DM`
1860d709ea83SDmitry Karpeev 
1861d709ea83SDmitry Karpeev   Options Database Key:
1862feefa0e1SJacob Faibussowitsch + -pc_gasm_dm_subdomains    - configure subdomains
1863feefa0e1SJacob Faibussowitsch . -pc_gasm_overlap          - set overlap
1864feefa0e1SJacob Faibussowitsch - -pc_gasm_total_subdomains - set number of subdomains
1865d709ea83SDmitry Karpeev 
1866d709ea83SDmitry Karpeev   Level: intermediate
1867d709ea83SDmitry Karpeev 
1868f1580f4eSBarry Smith   Note:
1869f1580f4eSBarry Smith   `PCGASMSetSubdomains()`, `PCGASMSetTotalSubdomains()` or `PCGASMSetOverlap()` take precedence over `PCGASMSetUseDMSubdomains()`,
1870f1580f4eSBarry Smith   so setting `PCGASMSetSubdomains()` with nontrivial subdomain ISs or any of `PCGASMSetTotalSubdomains()` and `PCGASMSetOverlap()`
18718f3b4b4dSDmitry Karpeev   automatically turns the latter off.
1872d709ea83SDmitry Karpeev 
1873f1580f4eSBarry Smith .seealso: `PCGASM`, `PCGASMGetUseDMSubdomains()`, `PCGASMSetSubdomains()`, `PCGASMSetOverlap()`
1874db781477SPatrick Sanan           `PCGASMCreateSubdomains2D()`
1875d709ea83SDmitry Karpeev @*/
1876d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGASMSetUseDMSubdomains(PC pc, PetscBool flg)
1877d71ae5a4SJacob Faibussowitsch {
1878d709ea83SDmitry Karpeev   PC_GASM  *osm = (PC_GASM *)pc->data;
1879d709ea83SDmitry Karpeev   PetscBool match;
1880d709ea83SDmitry Karpeev 
1881d709ea83SDmitry Karpeev   PetscFunctionBegin;
1882d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1883d709ea83SDmitry Karpeev   PetscValidLogicalCollectiveBool(pc, flg, 2);
188428b400f6SJacob Faibussowitsch   PetscCheck(!pc->setupcalled, ((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Not for a setup PC.");
18859566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCGASM, &match));
1886d709ea83SDmitry Karpeev   if (match) {
1887ad540459SPierre Jolivet     if (!osm->user_subdomains && osm->N == PETSC_DETERMINE && osm->overlap < 0) osm->dm_subdomains = flg;
18888f3b4b4dSDmitry Karpeev   }
18893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1890d709ea83SDmitry Karpeev }
1891d709ea83SDmitry Karpeev 
1892d709ea83SDmitry Karpeev /*@
1893f1580f4eSBarry Smith   PCGASMGetUseDMSubdomains - Returns flag indicating whether to use `DMCreateDomainDecomposition()` to define the subdomains, whenever possible with `PCGASM`
1894f1580f4eSBarry Smith 
1895d709ea83SDmitry Karpeev   Not Collective
1896d709ea83SDmitry Karpeev 
1897d709ea83SDmitry Karpeev   Input Parameter:
1898d709ea83SDmitry Karpeev . pc - the preconditioner
1899d709ea83SDmitry Karpeev 
1900d709ea83SDmitry Karpeev   Output Parameter:
1901f1580f4eSBarry Smith . flg - boolean indicating whether to use subdomains defined by the `DM`
1902d709ea83SDmitry Karpeev 
1903d709ea83SDmitry Karpeev   Level: intermediate
1904d709ea83SDmitry Karpeev 
1905f1580f4eSBarry Smith .seealso: `PCGASM`, `PCGASMSetUseDMSubdomains()`, `PCGASMSetOverlap()`
1906db781477SPatrick Sanan           `PCGASMCreateSubdomains2D()`
1907d709ea83SDmitry Karpeev @*/
1908d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGASMGetUseDMSubdomains(PC pc, PetscBool *flg)
1909d71ae5a4SJacob Faibussowitsch {
1910d709ea83SDmitry Karpeev   PC_GASM  *osm = (PC_GASM *)pc->data;
1911d709ea83SDmitry Karpeev   PetscBool match;
1912d709ea83SDmitry Karpeev 
1913d709ea83SDmitry Karpeev   PetscFunctionBegin;
1914d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
19154f572ea9SToby Isaac   PetscAssertPointer(flg, 2);
19169566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCGASM, &match));
1917d709ea83SDmitry Karpeev   if (match) {
1918d709ea83SDmitry Karpeev     if (flg) *flg = osm->dm_subdomains;
1919d709ea83SDmitry Karpeev   }
19203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1921d709ea83SDmitry Karpeev }
1922