xref: /petsc/src/dm/impls/patch/patch.c (revision 5d83a8b16d06840f96948f1a43aa9c83c769a60a)
1af0996ceSBarry Smith #include <petsc/private/dmpatchimpl.h> /*I      "petscdmpatch.h"   I*/
23a19ef87SMatthew G Knepley #include <petscdmda.h>
30c312b8eSJed Brown #include <petscsf.h>
43a19ef87SMatthew G Knepley 
55ee0e871SMatthew G Knepley /*
65ee0e871SMatthew G Knepley Solver loop to update \tau:
75ee0e871SMatthew G Knepley 
85ee0e871SMatthew G Knepley   DMZoom(dmc, &dmz)
95ee0e871SMatthew G Knepley   DMRefine(dmz, &dmf),
105ee0e871SMatthew G Knepley   Scatter Xcoarse -> Xzoom,
115ee0e871SMatthew G Knepley   Interpolate Xzoom -> Xfine (note that this may be on subcomms),
125ee0e871SMatthew G Knepley   Smooth Xfine using two-step smoother
135ee0e871SMatthew G Knepley     normal smoother plus Kaczmarz---moves back and forth from dmzoom to dmfine
145ee0e871SMatthew G Knepley   Compute residual Rfine
155ee0e871SMatthew G Knepley   Restrict Rfine to Rzoom_restricted
165ee0e871SMatthew G Knepley   Scatter Rzoom_restricted -> Rcoarse_restricted
175ee0e871SMatthew G Knepley   Compute global residual Rcoarse
185ee0e871SMatthew G Knepley   TauCoarse = Rcoarse - Rcoarse_restricted
195ee0e871SMatthew G Knepley */
205ee0e871SMatthew G Knepley 
21*5d83a8b1SBarry Smith /*@
2220f4b53cSBarry Smith   DMPatchZoom - Create patches of a `DMDA` on subsets of processes, indicated by `commz`
2359583ba0SMatthew G Knepley 
2420f4b53cSBarry Smith   Collective
25a487c981SJed Brown 
2659583ba0SMatthew G Knepley   Input Parameters:
2720f4b53cSBarry Smith + dm    - the `DM`
28c969029dSPierre Jolivet . lower - the lower left corner of the requested patch
29c969029dSPierre Jolivet . upper - the upper right corner of the requested patch
3020f4b53cSBarry Smith - commz - the new communicator for the patch, `MPI_COMM_NULL` indicates that the given rank will not own a patch
3159583ba0SMatthew G Knepley 
3259583ba0SMatthew G Knepley   Output Parameters:
3320f4b53cSBarry Smith + dmz  - the patch `DM`
3420f4b53cSBarry Smith . sfz  - the `PetscSF` mapping the patch+halo to the zoomed version (optional)
3520f4b53cSBarry Smith - sfzr - the `PetscSF` mapping the patch to the restricted zoomed version
3659583ba0SMatthew G Knepley 
3759583ba0SMatthew G Knepley   Level: intermediate
3859583ba0SMatthew G Knepley 
39db781477SPatrick Sanan .seealso: `DMPatchSolve()`, `DMDACreatePatchIS()`
4060c22052SBarry Smith @*/
41d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPatchZoom(DM dm, MatStencil lower, MatStencil upper, MPI_Comm commz, DM *dmz, PetscSF *sfz, PetscSF *sfzr)
42d71ae5a4SJacob Faibussowitsch {
4359583ba0SMatthew G Knepley   DMDAStencilType st;
44bd1050a3SMatthew G Knepley   MatStencil      blower, bupper, loclower, locupper;
45bd1050a3SMatthew G Knepley   IS              is;
46bd1050a3SMatthew G Knepley   const PetscInt *ranges, *indices;
470298fd71SBarry Smith   PetscInt       *localPoints  = NULL;
480298fd71SBarry Smith   PetscSFNode    *remotePoints = NULL;
49392fa6c6SMatthew G Knepley   PetscInt        dim, dof;
5060c22052SBarry Smith   PetscInt        M, N, P, rM, rN, rP, halo = 1, sxb, syb, szb, sxr, syr, szr, exr, eyr, ezr, mxb, myb, mzb, i, j, k, l, q;
51bd1050a3SMatthew G Knepley   PetscMPIInt     size;
523b5e53cdSSajid Ali   PetscBool       patchis_offproc = PETSC_TRUE;
5360c22052SBarry Smith   Vec             X;
5459583ba0SMatthew G Knepley 
5559583ba0SMatthew G Knepley   PetscFunctionBegin;
5660c22052SBarry Smith   if (!sfz) halo = 0;
579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
5859583ba0SMatthew G Knepley   /* Create patch DM */
599566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, &dim, &M, &N, &P, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, &st));
60a07761baSMatthew G Knepley 
61a07761baSMatthew G Knepley   /* Get piece for rank r, expanded by halo */
629371c9d4SSatish Balay   bupper.i = PetscMin(M, upper.i + halo);
639371c9d4SSatish Balay   blower.i = PetscMax(lower.i - halo, 0);
649371c9d4SSatish Balay   bupper.j = PetscMin(N, upper.j + halo);
659371c9d4SSatish Balay   blower.j = PetscMax(lower.j - halo, 0);
669371c9d4SSatish Balay   bupper.k = PetscMin(P, upper.k + halo);
679371c9d4SSatish Balay   blower.k = PetscMax(lower.k - halo, 0);
68392fa6c6SMatthew G Knepley   rM       = bupper.i - blower.i;
69392fa6c6SMatthew G Knepley   rN       = bupper.j - blower.j;
70392fa6c6SMatthew G Knepley   rP       = bupper.k - blower.k;
7159583ba0SMatthew G Knepley 
72bb71ef15SMatthew G Knepley   if (commz != MPI_COMM_NULL) {
739566063dSJacob Faibussowitsch     PetscCall(DMDACreate(commz, dmz));
749566063dSJacob Faibussowitsch     PetscCall(DMSetDimension(*dmz, dim));
759566063dSJacob Faibussowitsch     PetscCall(DMDASetSizes(*dmz, rM, rN, rP));
769566063dSJacob Faibussowitsch     PetscCall(DMDASetNumProcs(*dmz, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE));
779566063dSJacob Faibussowitsch     PetscCall(DMDASetBoundaryType(*dmz, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE));
789566063dSJacob Faibussowitsch     PetscCall(DMDASetDof(*dmz, dof));
799566063dSJacob Faibussowitsch     PetscCall(DMDASetStencilType(*dmz, st));
809566063dSJacob Faibussowitsch     PetscCall(DMDASetStencilWidth(*dmz, 0));
819566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(*dmz, NULL, NULL, NULL));
829566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(*dmz));
839566063dSJacob Faibussowitsch     PetscCall(DMSetUp(*dmz));
849566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(*dmz, &sxb, &syb, &szb, &mxb, &myb, &mzb));
85bd1050a3SMatthew G Knepley     sxr = PetscMax(sxb, lower.i - blower.i);
86bd1050a3SMatthew G Knepley     syr = PetscMax(syb, lower.j - blower.j);
87bd1050a3SMatthew G Knepley     szr = PetscMax(szb, lower.k - blower.k);
88bd1050a3SMatthew G Knepley     exr = PetscMin(sxb + mxb, upper.i - blower.i);
89bd1050a3SMatthew G Knepley     eyr = PetscMin(syb + myb, upper.j - blower.j);
90bd1050a3SMatthew G Knepley     ezr = PetscMin(szb + mzb, upper.k - blower.k);
919566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(dof * rM * rN * PetscMax(rP, 1), &localPoints, dof * rM * rN * PetscMax(rP, 1), &remotePoints));
92bb71ef15SMatthew G Knepley   } else {
93bb71ef15SMatthew G Knepley     sxr = syr = szr = exr = eyr = ezr = sxb = syb = szb = mxb = myb = mzb = 0;
94bb71ef15SMatthew G Knepley   }
95392fa6c6SMatthew G Knepley 
9659583ba0SMatthew G Knepley   /* Create SF for restricted map */
979566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &X));
989566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRanges(X, &ranges));
998865f1eaSKarl Rupp 
1009371c9d4SSatish Balay   loclower.i = blower.i + sxr;
1019371c9d4SSatish Balay   locupper.i = blower.i + exr;
1029371c9d4SSatish Balay   loclower.j = blower.j + syr;
1039371c9d4SSatish Balay   locupper.j = blower.j + eyr;
1049371c9d4SSatish Balay   loclower.k = blower.k + szr;
1059371c9d4SSatish Balay   locupper.k = blower.k + ezr;
1068865f1eaSKarl Rupp 
1079566063dSJacob Faibussowitsch   PetscCall(DMDACreatePatchIS(dm, &loclower, &locupper, &is, patchis_offproc));
1089566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is, &indices));
1098865f1eaSKarl Rupp 
1109371c9d4SSatish Balay   if (dim < 3) {
1119371c9d4SSatish Balay     mzb = 1;
1129371c9d4SSatish Balay     ezr = 1;
1139371c9d4SSatish Balay   }
114392fa6c6SMatthew G Knepley   q = 0;
115bd1050a3SMatthew G Knepley   for (k = szb; k < szb + mzb; ++k) {
116bd1050a3SMatthew G Knepley     if ((k < szr) || (k >= ezr)) continue;
117bd1050a3SMatthew G Knepley     for (j = syb; j < syb + myb; ++j) {
118bd1050a3SMatthew G Knepley       if ((j < syr) || (j >= eyr)) continue;
119bd1050a3SMatthew G Knepley       for (i = sxb; i < sxb + mxb; ++i) {
12060c22052SBarry Smith         for (l = 0; l < dof; l++) {
12160c22052SBarry Smith           const PetscInt lp = l + dof * (((k - szb) * rN + (j - syb)) * rM + i - sxb);
122bd1050a3SMatthew G Knepley           PetscInt       r;
123392fa6c6SMatthew G Knepley 
124bd1050a3SMatthew G Knepley           if ((i < sxr) || (i >= exr)) continue;
125392fa6c6SMatthew G Knepley           localPoints[q] = lp;
1269566063dSJacob Faibussowitsch           PetscCall(PetscFindInt(indices[q], size + 1, ranges, &r));
1278865f1eaSKarl Rupp 
128bd1050a3SMatthew G Knepley           remotePoints[q].rank  = r < 0 ? -(r + 1) - 1 : r;
129bd1050a3SMatthew G Knepley           remotePoints[q].index = indices[q] - ranges[remotePoints[q].rank];
130bd1050a3SMatthew G Knepley           ++q;
131392fa6c6SMatthew G Knepley         }
132392fa6c6SMatthew G Knepley       }
133392fa6c6SMatthew G Knepley     }
13460c22052SBarry Smith   }
1359566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is, &indices));
1369566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
1379566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), sfzr));
1389566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*sfzr, "Restricted Map"));
1399566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*sfzr, dof * M * N * P, q, localPoints, PETSC_COPY_VALUES, remotePoints, PETSC_COPY_VALUES));
140392fa6c6SMatthew G Knepley 
14160c22052SBarry Smith   if (sfz) {
142392fa6c6SMatthew G Knepley     /* Create SF for buffered map */
1439371c9d4SSatish Balay     loclower.i = blower.i + sxb;
1449371c9d4SSatish Balay     locupper.i = blower.i + sxb + mxb;
1459371c9d4SSatish Balay     loclower.j = blower.j + syb;
1469371c9d4SSatish Balay     locupper.j = blower.j + syb + myb;
1479371c9d4SSatish Balay     loclower.k = blower.k + szb;
1489371c9d4SSatish Balay     locupper.k = blower.k + szb + mzb;
1498865f1eaSKarl Rupp 
1509566063dSJacob Faibussowitsch     PetscCall(DMDACreatePatchIS(dm, &loclower, &locupper, &is, patchis_offproc));
1519566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is, &indices));
1528865f1eaSKarl Rupp 
153392fa6c6SMatthew G Knepley     q = 0;
154392fa6c6SMatthew G Knepley     for (k = szb; k < szb + mzb; ++k) {
155392fa6c6SMatthew G Knepley       for (j = syb; j < syb + myb; ++j) {
156392fa6c6SMatthew G Knepley         for (i = sxb; i < sxb + mxb; ++i, ++q) {
157bd1050a3SMatthew G Knepley           PetscInt r;
158392fa6c6SMatthew G Knepley 
159bd1050a3SMatthew G Knepley           localPoints[q] = q;
1609566063dSJacob Faibussowitsch           PetscCall(PetscFindInt(indices[q], size + 1, ranges, &r));
161bd1050a3SMatthew G Knepley           remotePoints[q].rank  = r < 0 ? -(r + 1) - 1 : r;
162bd1050a3SMatthew G Knepley           remotePoints[q].index = indices[q] - ranges[remotePoints[q].rank];
163392fa6c6SMatthew G Knepley         }
164392fa6c6SMatthew G Knepley       }
165392fa6c6SMatthew G Knepley     }
1669566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is, &indices));
1679566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is));
1689566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), sfz));
1699566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)*sfz, "Buffered Map"));
1709566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(*sfz, M * N * P, q, localPoints, PETSC_COPY_VALUES, remotePoints, PETSC_COPY_VALUES));
17160c22052SBarry Smith   }
172392fa6c6SMatthew G Knepley 
1739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
1749566063dSJacob Faibussowitsch   PetscCall(PetscFree2(localPoints, remotePoints));
1753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17659583ba0SMatthew G Knepley }
17759583ba0SMatthew G Knepley 
1789371c9d4SSatish Balay typedef enum {
1799371c9d4SSatish Balay   PATCH_COMM_TYPE_WORLD = 0,
1809371c9d4SSatish Balay   PATCH_COMM_TYPE_SELF  = 1
1819371c9d4SSatish Balay } PatchCommType;
182bd1050a3SMatthew G Knepley 
183d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPatchSolve(DM dm)
184d71ae5a4SJacob Faibussowitsch {
18582f516ccSBarry Smith   MPI_Comm    comm;
186bb71ef15SMatthew G Knepley   MPI_Comm    commz;
187bb71ef15SMatthew G Knepley   DM          dmc;
188392fa6c6SMatthew G Knepley   PetscSF     sfz, sfzr;
189bb71ef15SMatthew G Knepley   Vec         XC;
190bb71ef15SMatthew G Knepley   MatStencil  patchSize, commSize, gridRank, lower, upper;
191bb71ef15SMatthew G Knepley   PetscInt    M, N, P, i, j, k, l, m, n, p = 0;
192bb71ef15SMatthew G Knepley   PetscMPIInt rank, size;
193bb71ef15SMatthew G Knepley   PetscInt    debug = 0;
19459583ba0SMatthew G Knepley 
19559583ba0SMatthew G Knepley   PetscFunctionBegin;
1969566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1979566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1989566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
1999566063dSJacob Faibussowitsch   PetscCall(DMPatchGetCoarse(dm, &dmc));
2009566063dSJacob Faibussowitsch   PetscCall(DMPatchGetPatchSize(dm, &patchSize));
2019566063dSJacob Faibussowitsch   PetscCall(DMPatchGetCommSize(dm, &commSize));
2029566063dSJacob Faibussowitsch   PetscCall(DMPatchGetCommSize(dm, &commSize));
2039566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmc, &XC));
2049566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dmc, NULL, &M, &N, &P, &l, &m, &n, NULL, NULL, NULL, NULL, NULL, NULL));
2059371c9d4SSatish Balay   M = PetscMax(M, 1);
2069371c9d4SSatish Balay   l = PetscMax(l, 1);
2079371c9d4SSatish Balay   N = PetscMax(N, 1);
2089371c9d4SSatish Balay   m = PetscMax(m, 1);
2099371c9d4SSatish Balay   P = PetscMax(P, 1);
2109371c9d4SSatish Balay   n = PetscMax(n, 1);
2118865f1eaSKarl Rupp 
212bb71ef15SMatthew G Knepley   gridRank.i = rank % l;
213bb71ef15SMatthew G Knepley   gridRank.j = rank / l % m;
214bb71ef15SMatthew G Knepley   gridRank.k = rank / (l * m) % n;
2158865f1eaSKarl Rupp 
216bb71ef15SMatthew G Knepley   if (commSize.i * commSize.j * commSize.k == size || commSize.i * commSize.j * commSize.k == 0) {
2179371c9d4SSatish Balay     commSize.i = l;
2189371c9d4SSatish Balay     commSize.j = m;
2199371c9d4SSatish Balay     commSize.k = n;
220bb71ef15SMatthew G Knepley     commz      = comm;
221bb71ef15SMatthew G Knepley   } else if (commSize.i * commSize.j * commSize.k == 1) {
222bb71ef15SMatthew G Knepley     commz = PETSC_COMM_SELF;
223bb71ef15SMatthew G Knepley   } else {
224bb71ef15SMatthew G Knepley     const PetscMPIInt newComm = ((gridRank.k / commSize.k) * (m / commSize.j) + gridRank.j / commSize.j) * (l / commSize.i) + (gridRank.i / commSize.i);
225bb71ef15SMatthew G Knepley     const PetscMPIInt newRank = ((gridRank.k % commSize.k) * commSize.j + gridRank.j % commSize.j) * commSize.i + (gridRank.i % commSize.i);
226bb71ef15SMatthew G Knepley 
2279566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_split(comm, newComm, newRank, &commz));
2289566063dSJacob Faibussowitsch     if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Rank %d color %d key %d commz %p\n", rank, newComm, newRank, (void *)(MPI_Aint)commz));
229bb71ef15SMatthew G Knepley   }
230bb71ef15SMatthew G Knepley   /*
231bb71ef15SMatthew G Knepley    Assumptions:
232bb71ef15SMatthew G Knepley      - patchSize divides gridSize
233bb71ef15SMatthew G Knepley      - commSize divides gridSize
234bb71ef15SMatthew G Knepley      - commSize divides l,m,n
235bb71ef15SMatthew G Knepley    Ignore multiple patches per rank for now
236bb71ef15SMatthew G Knepley 
237bb71ef15SMatthew G Knepley    Multiple ranks per patch:
238bb71ef15SMatthew G Knepley      - l,m,n divides patchSize
239bb71ef15SMatthew G Knepley      - commSize divides patchSize
240bb71ef15SMatthew G Knepley    */
241392fa6c6SMatthew G Knepley   for (k = 0; k < P; k += PetscMax(patchSize.k, 1)) {
242392fa6c6SMatthew G Knepley     for (j = 0; j < N; j += PetscMax(patchSize.j, 1)) {
243392fa6c6SMatthew G Knepley       for (i = 0; i < M; i += PetscMax(patchSize.i, 1), ++p) {
244bb71ef15SMatthew G Knepley         MPI_Comm commp = MPI_COMM_NULL;
2450298fd71SBarry Smith         DM       dmz   = NULL;
2462f857c34SHong Zhang #if 0
2470298fd71SBarry Smith         DM          dmf     = NULL;
2480298fd71SBarry Smith         Mat         interpz = NULL;
2492f857c34SHong Zhang #endif
2500298fd71SBarry Smith         Vec          XZ      = NULL;
2510298fd71SBarry Smith         PetscScalar *xcarray = NULL;
2520298fd71SBarry Smith         PetscScalar *xzarray = NULL;
253bb87e006SMatthew G Knepley 
2549371c9d4SSatish Balay         if ((gridRank.k / commSize.k == p / (l / commSize.i * m / commSize.j) % n / commSize.k) && (gridRank.j / commSize.j == p / (l / commSize.i) % m / commSize.j) && (gridRank.i / commSize.i == p % l / commSize.i)) {
25563a3b9bcSJacob Faibussowitsch           if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Rank %d is accepting Patch %" PetscInt_FMT "\n", rank, p));
256bb71ef15SMatthew G Knepley           commp = commz;
257bb71ef15SMatthew G Knepley         }
258bb87e006SMatthew G Knepley         /* Zoom to coarse patch */
2599371c9d4SSatish Balay         lower.i = i;
2609371c9d4SSatish Balay         lower.j = j;
2619371c9d4SSatish Balay         lower.k = k;
2629371c9d4SSatish Balay         upper.i = i + patchSize.i;
2639371c9d4SSatish Balay         upper.j = j + patchSize.j;
2649371c9d4SSatish Balay         upper.k = k + patchSize.k;
2659566063dSJacob Faibussowitsch         PetscCall(DMPatchZoom(dmc, lower, upper, commp, &dmz, &sfz, &sfzr));
266da80777bSKarl Rupp         lower.c = 0; /* initialize member, otherwise compiler issues warnings */
267da80777bSKarl Rupp         upper.c = 0; /* initialize member, otherwise compiler issues warnings */
2689371c9d4SSatish Balay         if (debug)
2699371c9d4SSatish Balay           PetscCall(PetscPrintf(comm, "Patch %" PetscInt_FMT ": (%" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ")--(%" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ")\n", p, lower.i, lower.j, lower.k, upper.i, upper.j, upper.k));
2709566063dSJacob Faibussowitsch         if (dmz) PetscCall(DMView(dmz, PETSC_VIEWER_STDOUT_(commz)));
2719566063dSJacob Faibussowitsch         PetscCall(PetscSFView(sfz, PETSC_VIEWER_STDOUT_(comm)));
2729566063dSJacob Faibussowitsch         PetscCall(PetscSFView(sfzr, PETSC_VIEWER_STDOUT_(comm)));
27359583ba0SMatthew G Knepley         /* Scatter Xcoarse -> Xzoom */
2749566063dSJacob Faibussowitsch         if (dmz) PetscCall(DMGetGlobalVector(dmz, &XZ));
2759566063dSJacob Faibussowitsch         if (XZ) PetscCall(VecGetArray(XZ, &xzarray));
2769566063dSJacob Faibussowitsch         PetscCall(VecGetArray(XC, &xcarray));
2779566063dSJacob Faibussowitsch         PetscCall(PetscSFBcastBegin(sfz, MPIU_SCALAR, xcarray, xzarray, MPI_REPLACE));
2789566063dSJacob Faibussowitsch         PetscCall(PetscSFBcastEnd(sfz, MPIU_SCALAR, xcarray, xzarray, MPI_REPLACE));
2799566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(XC, &xcarray));
2809566063dSJacob Faibussowitsch         if (XZ) PetscCall(VecRestoreArray(XZ, &xzarray));
281bd1050a3SMatthew G Knepley #if 0
28259583ba0SMatthew G Knepley         /* Interpolate Xzoom -> Xfine, note that this may be on subcomms */
2839566063dSJacob Faibussowitsch         PetscCall(DMRefine(dmz, MPI_COMM_NULL, &dmf));
2849566063dSJacob Faibussowitsch         PetscCall(DMCreateInterpolation(dmz, dmf, &interpz, NULL));
2859566063dSJacob Faibussowitsch         PetscCall(DMInterpolate(dmz, interpz, dmf));
28659583ba0SMatthew G Knepley         /* Smooth Xfine using two-step smoother, normal smoother plus Kaczmarz---moves back and forth from dmzoom to dmfine */
28759583ba0SMatthew G Knepley         /* Compute residual Rfine */
28859583ba0SMatthew G Knepley         /* Restrict Rfine to Rzoom_restricted */
289bd1050a3SMatthew G Knepley #endif
29059583ba0SMatthew G Knepley         /* Scatter Rzoom_restricted -> Rcoarse_restricted */
2919566063dSJacob Faibussowitsch         if (XZ) PetscCall(VecGetArray(XZ, &xzarray));
2929566063dSJacob Faibussowitsch         PetscCall(VecGetArray(XC, &xcarray));
2939566063dSJacob Faibussowitsch         PetscCall(PetscSFReduceBegin(sfzr, MPIU_SCALAR, xzarray, xcarray, MPIU_SUM));
2949566063dSJacob Faibussowitsch         PetscCall(PetscSFReduceEnd(sfzr, MPIU_SCALAR, xzarray, xcarray, MPIU_SUM));
2959566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(XC, &xcarray));
2969566063dSJacob Faibussowitsch         if (XZ) PetscCall(VecRestoreArray(XZ, &xzarray));
2979566063dSJacob Faibussowitsch         if (dmz) PetscCall(DMRestoreGlobalVector(dmz, &XZ));
29859583ba0SMatthew G Knepley         /* Compute global residual Rcoarse */
29959583ba0SMatthew G Knepley         /* TauCoarse = Rcoarse - Rcoarse_restricted */
300bb87e006SMatthew G Knepley 
3019566063dSJacob Faibussowitsch         PetscCall(PetscSFDestroy(&sfz));
3029566063dSJacob Faibussowitsch         PetscCall(PetscSFDestroy(&sfzr));
3039566063dSJacob Faibussowitsch         PetscCall(DMDestroy(&dmz));
304a07761baSMatthew G Knepley       }
305392fa6c6SMatthew G Knepley     }
306392fa6c6SMatthew G Knepley   }
3079566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmc, &XC));
3083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30959583ba0SMatthew G Knepley }
31059583ba0SMatthew G Knepley 
31166976f2fSJacob Faibussowitsch static PetscErrorCode DMPatchView_ASCII(DM dm, PetscViewer viewer)
312d71ae5a4SJacob Faibussowitsch {
3133a19ef87SMatthew G Knepley   DM_Patch         *mesh = (DM_Patch *)dm->data;
3143a19ef87SMatthew G Knepley   PetscViewerFormat format;
315392fa6c6SMatthew G Knepley   const char       *name;
3163a19ef87SMatthew G Knepley 
3173a19ef87SMatthew G Knepley   PetscFunctionBegin;
3189566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
3193a19ef87SMatthew G Knepley   /* if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) */
3209566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
3219566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Patch DM %s\n", name));
3229566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
3239566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Coarse DM\n"));
3249566063dSJacob Faibussowitsch   PetscCall(DMView(mesh->dmCoarse, viewer));
3259566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
3263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3273a19ef87SMatthew G Knepley }
3283a19ef87SMatthew G Knepley 
329d71ae5a4SJacob Faibussowitsch PetscErrorCode DMView_Patch(DM dm, PetscViewer viewer)
330d71ae5a4SJacob Faibussowitsch {
3313a19ef87SMatthew G Knepley   PetscBool iascii, isbinary;
3323a19ef87SMatthew G Knepley 
3333a19ef87SMatthew G Knepley   PetscFunctionBegin;
3343a19ef87SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3353a19ef87SMatthew G Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
3369566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
3379566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
3381baa6e33SBarry Smith   if (iascii) PetscCall(DMPatchView_ASCII(dm, viewer));
3393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3403a19ef87SMatthew G Knepley }
3413a19ef87SMatthew G Knepley 
342d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDestroy_Patch(DM dm)
343d71ae5a4SJacob Faibussowitsch {
3443a19ef87SMatthew G Knepley   DM_Patch *mesh = (DM_Patch *)dm->data;
3453a19ef87SMatthew G Knepley 
3463a19ef87SMatthew G Knepley   PetscFunctionBegin;
3473ba16761SJacob Faibussowitsch   if (--mesh->refct > 0) PetscFunctionReturn(PETSC_SUCCESS);
3489566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&mesh->dmCoarse));
3493a19ef87SMatthew G Knepley   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
3509566063dSJacob Faibussowitsch   PetscCall(PetscFree(mesh));
3513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3523a19ef87SMatthew G Knepley }
3533a19ef87SMatthew G Knepley 
354d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetUp_Patch(DM dm)
355d71ae5a4SJacob Faibussowitsch {
3563a19ef87SMatthew G Knepley   DM_Patch *mesh = (DM_Patch *)dm->data;
3573a19ef87SMatthew G Knepley 
3583a19ef87SMatthew G Knepley   PetscFunctionBegin;
3593a19ef87SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3609566063dSJacob Faibussowitsch   PetscCall(DMSetUp(mesh->dmCoarse));
3613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3623a19ef87SMatthew G Knepley }
3633a19ef87SMatthew G Knepley 
364d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateGlobalVector_Patch(DM dm, Vec *g)
365d71ae5a4SJacob Faibussowitsch {
3663a19ef87SMatthew G Knepley   DM_Patch *mesh = (DM_Patch *)dm->data;
3673a19ef87SMatthew G Knepley 
3683a19ef87SMatthew G Knepley   PetscFunctionBegin;
3693a19ef87SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3709566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(mesh->dmCoarse, g));
3713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3723a19ef87SMatthew G Knepley }
3733a19ef87SMatthew G Knepley 
374d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateLocalVector_Patch(DM dm, Vec *l)
375d71ae5a4SJacob Faibussowitsch {
3763a19ef87SMatthew G Knepley   DM_Patch *mesh = (DM_Patch *)dm->data;
3773a19ef87SMatthew G Knepley 
3783a19ef87SMatthew G Knepley   PetscFunctionBegin;
3793a19ef87SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3809566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(mesh->dmCoarse, l));
3813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3823a19ef87SMatthew G Knepley }
3833a19ef87SMatthew G Knepley 
384d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateSubDM_Patch(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
385d71ae5a4SJacob Faibussowitsch {
38682f516ccSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Tell me to code this");
3873a19ef87SMatthew G Knepley }
388392fa6c6SMatthew G Knepley 
389d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPatchGetCoarse(DM dm, DM *dmCoarse)
390d71ae5a4SJacob Faibussowitsch {
391392fa6c6SMatthew G Knepley   DM_Patch *mesh = (DM_Patch *)dm->data;
392392fa6c6SMatthew G Knepley 
393392fa6c6SMatthew G Knepley   PetscFunctionBegin;
394392fa6c6SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
395392fa6c6SMatthew G Knepley   *dmCoarse = mesh->dmCoarse;
3963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
397392fa6c6SMatthew G Knepley }
398392fa6c6SMatthew G Knepley 
399d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPatchGetPatchSize(DM dm, MatStencil *patchSize)
400d71ae5a4SJacob Faibussowitsch {
401392fa6c6SMatthew G Knepley   DM_Patch *mesh = (DM_Patch *)dm->data;
402392fa6c6SMatthew G Knepley 
403392fa6c6SMatthew G Knepley   PetscFunctionBegin;
404392fa6c6SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4054f572ea9SToby Isaac   PetscAssertPointer(patchSize, 2);
406392fa6c6SMatthew G Knepley   *patchSize = mesh->patchSize;
4073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
408392fa6c6SMatthew G Knepley }
409392fa6c6SMatthew G Knepley 
410d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPatchSetPatchSize(DM dm, MatStencil patchSize)
411d71ae5a4SJacob Faibussowitsch {
412392fa6c6SMatthew G Knepley   DM_Patch *mesh = (DM_Patch *)dm->data;
413392fa6c6SMatthew G Knepley 
414392fa6c6SMatthew G Knepley   PetscFunctionBegin;
415392fa6c6SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
416392fa6c6SMatthew G Knepley   mesh->patchSize = patchSize;
4173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
418392fa6c6SMatthew G Knepley }
419bb71ef15SMatthew G Knepley 
420d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPatchGetCommSize(DM dm, MatStencil *commSize)
421d71ae5a4SJacob Faibussowitsch {
422bb71ef15SMatthew G Knepley   DM_Patch *mesh = (DM_Patch *)dm->data;
423bb71ef15SMatthew G Knepley 
424bb71ef15SMatthew G Knepley   PetscFunctionBegin;
425bb71ef15SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4264f572ea9SToby Isaac   PetscAssertPointer(commSize, 2);
427bb71ef15SMatthew G Knepley   *commSize = mesh->commSize;
4283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
429bb71ef15SMatthew G Knepley }
430bb71ef15SMatthew G Knepley 
431d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPatchSetCommSize(DM dm, MatStencil commSize)
432d71ae5a4SJacob Faibussowitsch {
433bb71ef15SMatthew G Knepley   DM_Patch *mesh = (DM_Patch *)dm->data;
434bb71ef15SMatthew G Knepley 
435bb71ef15SMatthew G Knepley   PetscFunctionBegin;
436bb71ef15SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
437bb71ef15SMatthew G Knepley   mesh->commSize = commSize;
4383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
439bb71ef15SMatthew G Knepley }
440