xref: /petsc/src/dm/impls/patch/patch.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a) !
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 
2160c22052SBarry Smith /*@C
2260c22052SBarry Smith   DMPatchZoom - Create patches of a DMDA on subsets of processes, indicated by commz
2359583ba0SMatthew G Knepley 
24d083f849SBarry Smith   Collective on dm
25a487c981SJed Brown 
2659583ba0SMatthew G Knepley   Input Parameters:
2759583ba0SMatthew G Knepley + 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
3060c22052SBarry 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:
3359583ba0SMatthew G Knepley + dmz  - the patch DM
3460c22052SBarry Smith . sfz  - the PetscSF mapping the patch+halo to the zoomed version (optional)
3560c22052SBarry 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 @*/
41*9371c9d4SSatish Balay PetscErrorCode DMPatchZoom(DM dm, MatStencil lower, MatStencil upper, MPI_Comm commz, DM *dmz, PetscSF *sfz, PetscSF *sfzr) {
4259583ba0SMatthew G Knepley   DMDAStencilType st;
43bd1050a3SMatthew G Knepley   MatStencil      blower, bupper, loclower, locupper;
44bd1050a3SMatthew G Knepley   IS              is;
45bd1050a3SMatthew G Knepley   const PetscInt *ranges, *indices;
460298fd71SBarry Smith   PetscInt       *localPoints  = NULL;
470298fd71SBarry Smith   PetscSFNode    *remotePoints = NULL;
48392fa6c6SMatthew G Knepley   PetscInt        dim, dof;
4960c22052SBarry 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;
50bd1050a3SMatthew G Knepley   PetscMPIInt     size;
513b5e53cdSSajid Ali   PetscBool       patchis_offproc = PETSC_TRUE;
5260c22052SBarry Smith   Vec             X;
5359583ba0SMatthew G Knepley 
5459583ba0SMatthew G Knepley   PetscFunctionBegin;
5560c22052SBarry Smith   if (!sfz) halo = 0;
569566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
5759583ba0SMatthew G Knepley   /* Create patch DM */
589566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, &dim, &M, &N, &P, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, &st));
59a07761baSMatthew G Knepley 
60a07761baSMatthew G Knepley   /* Get piece for rank r, expanded by halo */
61*9371c9d4SSatish Balay   bupper.i = PetscMin(M, upper.i + halo);
62*9371c9d4SSatish Balay   blower.i = PetscMax(lower.i - halo, 0);
63*9371c9d4SSatish Balay   bupper.j = PetscMin(N, upper.j + halo);
64*9371c9d4SSatish Balay   blower.j = PetscMax(lower.j - halo, 0);
65*9371c9d4SSatish Balay   bupper.k = PetscMin(P, upper.k + halo);
66*9371c9d4SSatish Balay   blower.k = PetscMax(lower.k - halo, 0);
67392fa6c6SMatthew G Knepley   rM       = bupper.i - blower.i;
68392fa6c6SMatthew G Knepley   rN       = bupper.j - blower.j;
69392fa6c6SMatthew G Knepley   rP       = bupper.k - blower.k;
7059583ba0SMatthew G Knepley 
71bb71ef15SMatthew G Knepley   if (commz != MPI_COMM_NULL) {
729566063dSJacob Faibussowitsch     PetscCall(DMDACreate(commz, dmz));
739566063dSJacob Faibussowitsch     PetscCall(DMSetDimension(*dmz, dim));
749566063dSJacob Faibussowitsch     PetscCall(DMDASetSizes(*dmz, rM, rN, rP));
759566063dSJacob Faibussowitsch     PetscCall(DMDASetNumProcs(*dmz, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE));
769566063dSJacob Faibussowitsch     PetscCall(DMDASetBoundaryType(*dmz, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE));
779566063dSJacob Faibussowitsch     PetscCall(DMDASetDof(*dmz, dof));
789566063dSJacob Faibussowitsch     PetscCall(DMDASetStencilType(*dmz, st));
799566063dSJacob Faibussowitsch     PetscCall(DMDASetStencilWidth(*dmz, 0));
809566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(*dmz, NULL, NULL, NULL));
819566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(*dmz));
829566063dSJacob Faibussowitsch     PetscCall(DMSetUp(*dmz));
839566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(*dmz, &sxb, &syb, &szb, &mxb, &myb, &mzb));
84bd1050a3SMatthew G Knepley     sxr = PetscMax(sxb, lower.i - blower.i);
85bd1050a3SMatthew G Knepley     syr = PetscMax(syb, lower.j - blower.j);
86bd1050a3SMatthew G Knepley     szr = PetscMax(szb, lower.k - blower.k);
87bd1050a3SMatthew G Knepley     exr = PetscMin(sxb + mxb, upper.i - blower.i);
88bd1050a3SMatthew G Knepley     eyr = PetscMin(syb + myb, upper.j - blower.j);
89bd1050a3SMatthew G Knepley     ezr = PetscMin(szb + mzb, upper.k - blower.k);
909566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(dof * rM * rN * PetscMax(rP, 1), &localPoints, dof * rM * rN * PetscMax(rP, 1), &remotePoints));
91bb71ef15SMatthew G Knepley   } else {
92bb71ef15SMatthew G Knepley     sxr = syr = szr = exr = eyr = ezr = sxb = syb = szb = mxb = myb = mzb = 0;
93bb71ef15SMatthew G Knepley   }
94392fa6c6SMatthew G Knepley 
9559583ba0SMatthew G Knepley   /* Create SF for restricted map */
969566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &X));
979566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRanges(X, &ranges));
988865f1eaSKarl Rupp 
99*9371c9d4SSatish Balay   loclower.i = blower.i + sxr;
100*9371c9d4SSatish Balay   locupper.i = blower.i + exr;
101*9371c9d4SSatish Balay   loclower.j = blower.j + syr;
102*9371c9d4SSatish Balay   locupper.j = blower.j + eyr;
103*9371c9d4SSatish Balay   loclower.k = blower.k + szr;
104*9371c9d4SSatish Balay   locupper.k = blower.k + ezr;
1058865f1eaSKarl Rupp 
1069566063dSJacob Faibussowitsch   PetscCall(DMDACreatePatchIS(dm, &loclower, &locupper, &is, patchis_offproc));
1079566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is, &indices));
1088865f1eaSKarl Rupp 
109*9371c9d4SSatish Balay   if (dim < 3) {
110*9371c9d4SSatish Balay     mzb = 1;
111*9371c9d4SSatish Balay     ezr = 1;
112*9371c9d4SSatish Balay   }
113392fa6c6SMatthew G Knepley   q = 0;
114bd1050a3SMatthew G Knepley   for (k = szb; k < szb + mzb; ++k) {
115bd1050a3SMatthew G Knepley     if ((k < szr) || (k >= ezr)) continue;
116bd1050a3SMatthew G Knepley     for (j = syb; j < syb + myb; ++j) {
117bd1050a3SMatthew G Knepley       if ((j < syr) || (j >= eyr)) continue;
118bd1050a3SMatthew G Knepley       for (i = sxb; i < sxb + mxb; ++i) {
11960c22052SBarry Smith         for (l = 0; l < dof; l++) {
12060c22052SBarry Smith           const PetscInt lp = l + dof * (((k - szb) * rN + (j - syb)) * rM + i - sxb);
121bd1050a3SMatthew G Knepley           PetscInt       r;
122392fa6c6SMatthew G Knepley 
123bd1050a3SMatthew G Knepley           if ((i < sxr) || (i >= exr)) continue;
124392fa6c6SMatthew G Knepley           localPoints[q] = lp;
1259566063dSJacob Faibussowitsch           PetscCall(PetscFindInt(indices[q], size + 1, ranges, &r));
1268865f1eaSKarl Rupp 
127bd1050a3SMatthew G Knepley           remotePoints[q].rank  = r < 0 ? -(r + 1) - 1 : r;
128bd1050a3SMatthew G Knepley           remotePoints[q].index = indices[q] - ranges[remotePoints[q].rank];
129bd1050a3SMatthew G Knepley           ++q;
130392fa6c6SMatthew G Knepley         }
131392fa6c6SMatthew G Knepley       }
132392fa6c6SMatthew G Knepley     }
13360c22052SBarry Smith   }
1349566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is, &indices));
1359566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
1369566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), sfzr));
1379566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*sfzr, "Restricted Map"));
1389566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*sfzr, dof * M * N * P, q, localPoints, PETSC_COPY_VALUES, remotePoints, PETSC_COPY_VALUES));
139392fa6c6SMatthew G Knepley 
14060c22052SBarry Smith   if (sfz) {
141392fa6c6SMatthew G Knepley     /* Create SF for buffered map */
142*9371c9d4SSatish Balay     loclower.i = blower.i + sxb;
143*9371c9d4SSatish Balay     locupper.i = blower.i + sxb + mxb;
144*9371c9d4SSatish Balay     loclower.j = blower.j + syb;
145*9371c9d4SSatish Balay     locupper.j = blower.j + syb + myb;
146*9371c9d4SSatish Balay     loclower.k = blower.k + szb;
147*9371c9d4SSatish Balay     locupper.k = blower.k + szb + mzb;
1488865f1eaSKarl Rupp 
1499566063dSJacob Faibussowitsch     PetscCall(DMDACreatePatchIS(dm, &loclower, &locupper, &is, patchis_offproc));
1509566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is, &indices));
1518865f1eaSKarl Rupp 
152392fa6c6SMatthew G Knepley     q = 0;
153392fa6c6SMatthew G Knepley     for (k = szb; k < szb + mzb; ++k) {
154392fa6c6SMatthew G Knepley       for (j = syb; j < syb + myb; ++j) {
155392fa6c6SMatthew G Knepley         for (i = sxb; i < sxb + mxb; ++i, ++q) {
156bd1050a3SMatthew G Knepley           PetscInt r;
157392fa6c6SMatthew G Knepley 
158bd1050a3SMatthew G Knepley           localPoints[q] = q;
1599566063dSJacob Faibussowitsch           PetscCall(PetscFindInt(indices[q], size + 1, ranges, &r));
160bd1050a3SMatthew G Knepley           remotePoints[q].rank  = r < 0 ? -(r + 1) - 1 : r;
161bd1050a3SMatthew G Knepley           remotePoints[q].index = indices[q] - ranges[remotePoints[q].rank];
162392fa6c6SMatthew G Knepley         }
163392fa6c6SMatthew G Knepley       }
164392fa6c6SMatthew G Knepley     }
1659566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is, &indices));
1669566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is));
1679566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), sfz));
1689566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)*sfz, "Buffered Map"));
1699566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(*sfz, M * N * P, q, localPoints, PETSC_COPY_VALUES, remotePoints, PETSC_COPY_VALUES));
17060c22052SBarry Smith   }
171392fa6c6SMatthew G Knepley 
1729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
1739566063dSJacob Faibussowitsch   PetscCall(PetscFree2(localPoints, remotePoints));
17459583ba0SMatthew G Knepley   PetscFunctionReturn(0);
17559583ba0SMatthew G Knepley }
17659583ba0SMatthew G Knepley 
177*9371c9d4SSatish Balay typedef enum {
178*9371c9d4SSatish Balay   PATCH_COMM_TYPE_WORLD = 0,
179*9371c9d4SSatish Balay   PATCH_COMM_TYPE_SELF  = 1
180*9371c9d4SSatish Balay } PatchCommType;
181bd1050a3SMatthew G Knepley 
182*9371c9d4SSatish Balay PetscErrorCode DMPatchSolve(DM dm) {
18382f516ccSBarry Smith   MPI_Comm    comm;
184bb71ef15SMatthew G Knepley   MPI_Comm    commz;
185bb71ef15SMatthew G Knepley   DM          dmc;
186392fa6c6SMatthew G Knepley   PetscSF     sfz, sfzr;
187bb71ef15SMatthew G Knepley   Vec         XC;
188bb71ef15SMatthew G Knepley   MatStencil  patchSize, commSize, gridRank, lower, upper;
189bb71ef15SMatthew G Knepley   PetscInt    M, N, P, i, j, k, l, m, n, p = 0;
190bb71ef15SMatthew G Knepley   PetscMPIInt rank, size;
191bb71ef15SMatthew G Knepley   PetscInt    debug = 0;
19259583ba0SMatthew G Knepley 
19359583ba0SMatthew G Knepley   PetscFunctionBegin;
1949566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1959566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1969566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
1979566063dSJacob Faibussowitsch   PetscCall(DMPatchGetCoarse(dm, &dmc));
1989566063dSJacob Faibussowitsch   PetscCall(DMPatchGetPatchSize(dm, &patchSize));
1999566063dSJacob Faibussowitsch   PetscCall(DMPatchGetCommSize(dm, &commSize));
2009566063dSJacob Faibussowitsch   PetscCall(DMPatchGetCommSize(dm, &commSize));
2019566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmc, &XC));
2029566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dmc, NULL, &M, &N, &P, &l, &m, &n, NULL, NULL, NULL, NULL, NULL, NULL));
203*9371c9d4SSatish Balay   M = PetscMax(M, 1);
204*9371c9d4SSatish Balay   l = PetscMax(l, 1);
205*9371c9d4SSatish Balay   N = PetscMax(N, 1);
206*9371c9d4SSatish Balay   m = PetscMax(m, 1);
207*9371c9d4SSatish Balay   P = PetscMax(P, 1);
208*9371c9d4SSatish Balay   n = PetscMax(n, 1);
2098865f1eaSKarl Rupp 
210bb71ef15SMatthew G Knepley   gridRank.i = rank % l;
211bb71ef15SMatthew G Knepley   gridRank.j = rank / l % m;
212bb71ef15SMatthew G Knepley   gridRank.k = rank / (l * m) % n;
2138865f1eaSKarl Rupp 
214bb71ef15SMatthew G Knepley   if (commSize.i * commSize.j * commSize.k == size || commSize.i * commSize.j * commSize.k == 0) {
215*9371c9d4SSatish Balay     commSize.i = l;
216*9371c9d4SSatish Balay     commSize.j = m;
217*9371c9d4SSatish Balay     commSize.k = n;
218bb71ef15SMatthew G Knepley     commz      = comm;
219bb71ef15SMatthew G Knepley   } else if (commSize.i * commSize.j * commSize.k == 1) {
220bb71ef15SMatthew G Knepley     commz = PETSC_COMM_SELF;
221bb71ef15SMatthew G Knepley   } else {
222bb71ef15SMatthew G Knepley     const PetscMPIInt newComm = ((gridRank.k / commSize.k) * (m / commSize.j) + gridRank.j / commSize.j) * (l / commSize.i) + (gridRank.i / commSize.i);
223bb71ef15SMatthew G Knepley     const PetscMPIInt newRank = ((gridRank.k % commSize.k) * commSize.j + gridRank.j % commSize.j) * commSize.i + (gridRank.i % commSize.i);
224bb71ef15SMatthew G Knepley 
2259566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_split(comm, newComm, newRank, &commz));
2269566063dSJacob Faibussowitsch     if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Rank %d color %d key %d commz %p\n", rank, newComm, newRank, (void *)(MPI_Aint)commz));
227bb71ef15SMatthew G Knepley   }
228bb71ef15SMatthew G Knepley   /*
229bb71ef15SMatthew G Knepley    Assumptions:
230bb71ef15SMatthew G Knepley      - patchSize divides gridSize
231bb71ef15SMatthew G Knepley      - commSize divides gridSize
232bb71ef15SMatthew G Knepley      - commSize divides l,m,n
233bb71ef15SMatthew G Knepley    Ignore multiple patches per rank for now
234bb71ef15SMatthew G Knepley 
235bb71ef15SMatthew G Knepley    Multiple ranks per patch:
236bb71ef15SMatthew G Knepley      - l,m,n divides patchSize
237bb71ef15SMatthew G Knepley      - commSize divides patchSize
238bb71ef15SMatthew G Knepley    */
239392fa6c6SMatthew G Knepley   for (k = 0; k < P; k += PetscMax(patchSize.k, 1)) {
240392fa6c6SMatthew G Knepley     for (j = 0; j < N; j += PetscMax(patchSize.j, 1)) {
241392fa6c6SMatthew G Knepley       for (i = 0; i < M; i += PetscMax(patchSize.i, 1), ++p) {
242bb71ef15SMatthew G Knepley         MPI_Comm commp = MPI_COMM_NULL;
2430298fd71SBarry Smith         DM       dmz   = NULL;
2442f857c34SHong Zhang #if 0
2450298fd71SBarry Smith         DM          dmf     = NULL;
2460298fd71SBarry Smith         Mat         interpz = NULL;
2472f857c34SHong Zhang #endif
2480298fd71SBarry Smith         Vec          XZ      = NULL;
2490298fd71SBarry Smith         PetscScalar *xcarray = NULL;
2500298fd71SBarry Smith         PetscScalar *xzarray = NULL;
251bb87e006SMatthew G Knepley 
252*9371c9d4SSatish 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)) {
25363a3b9bcSJacob Faibussowitsch           if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Rank %d is accepting Patch %" PetscInt_FMT "\n", rank, p));
254bb71ef15SMatthew G Knepley           commp = commz;
255bb71ef15SMatthew G Knepley         }
256bb87e006SMatthew G Knepley         /* Zoom to coarse patch */
257*9371c9d4SSatish Balay         lower.i = i;
258*9371c9d4SSatish Balay         lower.j = j;
259*9371c9d4SSatish Balay         lower.k = k;
260*9371c9d4SSatish Balay         upper.i = i + patchSize.i;
261*9371c9d4SSatish Balay         upper.j = j + patchSize.j;
262*9371c9d4SSatish Balay         upper.k = k + patchSize.k;
2639566063dSJacob Faibussowitsch         PetscCall(DMPatchZoom(dmc, lower, upper, commp, &dmz, &sfz, &sfzr));
264da80777bSKarl Rupp         lower.c = 0; /* initialize member, otherwise compiler issues warnings */
265da80777bSKarl Rupp         upper.c = 0; /* initialize member, otherwise compiler issues warnings */
266*9371c9d4SSatish Balay         if (debug)
267*9371c9d4SSatish 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));
2689566063dSJacob Faibussowitsch         if (dmz) PetscCall(DMView(dmz, PETSC_VIEWER_STDOUT_(commz)));
2699566063dSJacob Faibussowitsch         PetscCall(PetscSFView(sfz, PETSC_VIEWER_STDOUT_(comm)));
2709566063dSJacob Faibussowitsch         PetscCall(PetscSFView(sfzr, PETSC_VIEWER_STDOUT_(comm)));
27159583ba0SMatthew G Knepley         /* Scatter Xcoarse -> Xzoom */
2729566063dSJacob Faibussowitsch         if (dmz) PetscCall(DMGetGlobalVector(dmz, &XZ));
2739566063dSJacob Faibussowitsch         if (XZ) PetscCall(VecGetArray(XZ, &xzarray));
2749566063dSJacob Faibussowitsch         PetscCall(VecGetArray(XC, &xcarray));
2759566063dSJacob Faibussowitsch         PetscCall(PetscSFBcastBegin(sfz, MPIU_SCALAR, xcarray, xzarray, MPI_REPLACE));
2769566063dSJacob Faibussowitsch         PetscCall(PetscSFBcastEnd(sfz, MPIU_SCALAR, xcarray, xzarray, MPI_REPLACE));
2779566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(XC, &xcarray));
2789566063dSJacob Faibussowitsch         if (XZ) PetscCall(VecRestoreArray(XZ, &xzarray));
279bd1050a3SMatthew G Knepley #if 0
28059583ba0SMatthew G Knepley         /* Interpolate Xzoom -> Xfine, note that this may be on subcomms */
2819566063dSJacob Faibussowitsch         PetscCall(DMRefine(dmz, MPI_COMM_NULL, &dmf));
2829566063dSJacob Faibussowitsch         PetscCall(DMCreateInterpolation(dmz, dmf, &interpz, NULL));
2839566063dSJacob Faibussowitsch         PetscCall(DMInterpolate(dmz, interpz, dmf));
28459583ba0SMatthew G Knepley         /* Smooth Xfine using two-step smoother, normal smoother plus Kaczmarz---moves back and forth from dmzoom to dmfine */
28559583ba0SMatthew G Knepley         /* Compute residual Rfine */
28659583ba0SMatthew G Knepley         /* Restrict Rfine to Rzoom_restricted */
287bd1050a3SMatthew G Knepley #endif
28859583ba0SMatthew G Knepley         /* Scatter Rzoom_restricted -> Rcoarse_restricted */
2899566063dSJacob Faibussowitsch         if (XZ) PetscCall(VecGetArray(XZ, &xzarray));
2909566063dSJacob Faibussowitsch         PetscCall(VecGetArray(XC, &xcarray));
2919566063dSJacob Faibussowitsch         PetscCall(PetscSFReduceBegin(sfzr, MPIU_SCALAR, xzarray, xcarray, MPIU_SUM));
2929566063dSJacob Faibussowitsch         PetscCall(PetscSFReduceEnd(sfzr, MPIU_SCALAR, xzarray, xcarray, MPIU_SUM));
2939566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(XC, &xcarray));
2949566063dSJacob Faibussowitsch         if (XZ) PetscCall(VecRestoreArray(XZ, &xzarray));
2959566063dSJacob Faibussowitsch         if (dmz) PetscCall(DMRestoreGlobalVector(dmz, &XZ));
29659583ba0SMatthew G Knepley         /* Compute global residual Rcoarse */
29759583ba0SMatthew G Knepley         /* TauCoarse = Rcoarse - Rcoarse_restricted */
298bb87e006SMatthew G Knepley 
2999566063dSJacob Faibussowitsch         PetscCall(PetscSFDestroy(&sfz));
3009566063dSJacob Faibussowitsch         PetscCall(PetscSFDestroy(&sfzr));
3019566063dSJacob Faibussowitsch         PetscCall(DMDestroy(&dmz));
302a07761baSMatthew G Knepley       }
303392fa6c6SMatthew G Knepley     }
304392fa6c6SMatthew G Knepley   }
3059566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmc, &XC));
30659583ba0SMatthew G Knepley   PetscFunctionReturn(0);
30759583ba0SMatthew G Knepley }
30859583ba0SMatthew G Knepley 
309*9371c9d4SSatish Balay PetscErrorCode DMPatchView_ASCII(DM dm, PetscViewer viewer) {
3103a19ef87SMatthew G Knepley   DM_Patch         *mesh = (DM_Patch *)dm->data;
3113a19ef87SMatthew G Knepley   PetscViewerFormat format;
312392fa6c6SMatthew G Knepley   const char       *name;
3133a19ef87SMatthew G Knepley 
3143a19ef87SMatthew G Knepley   PetscFunctionBegin;
3159566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
3163a19ef87SMatthew G Knepley   /* if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) */
3179566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
3189566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Patch DM %s\n", name));
3199566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
3209566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Coarse DM\n"));
3219566063dSJacob Faibussowitsch   PetscCall(DMView(mesh->dmCoarse, viewer));
3229566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
3233a19ef87SMatthew G Knepley   PetscFunctionReturn(0);
3243a19ef87SMatthew G Knepley }
3253a19ef87SMatthew G Knepley 
326*9371c9d4SSatish Balay PetscErrorCode DMView_Patch(DM dm, PetscViewer viewer) {
3273a19ef87SMatthew G Knepley   PetscBool iascii, isbinary;
3283a19ef87SMatthew G Knepley 
3293a19ef87SMatthew G Knepley   PetscFunctionBegin;
3303a19ef87SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3313a19ef87SMatthew G Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
3329566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
3339566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
3341baa6e33SBarry Smith   if (iascii) PetscCall(DMPatchView_ASCII(dm, viewer));
3353a19ef87SMatthew G Knepley   PetscFunctionReturn(0);
3363a19ef87SMatthew G Knepley }
3373a19ef87SMatthew G Knepley 
338*9371c9d4SSatish Balay PetscErrorCode DMDestroy_Patch(DM dm) {
3393a19ef87SMatthew G Knepley   DM_Patch *mesh = (DM_Patch *)dm->data;
3403a19ef87SMatthew G Knepley 
3413a19ef87SMatthew G Knepley   PetscFunctionBegin;
3428865f1eaSKarl Rupp   if (--mesh->refct > 0) PetscFunctionReturn(0);
3439566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&mesh->dmCoarse));
3443a19ef87SMatthew G Knepley   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
3459566063dSJacob Faibussowitsch   PetscCall(PetscFree(mesh));
3463a19ef87SMatthew G Knepley   PetscFunctionReturn(0);
3473a19ef87SMatthew G Knepley }
3483a19ef87SMatthew G Knepley 
349*9371c9d4SSatish Balay PetscErrorCode DMSetUp_Patch(DM dm) {
3503a19ef87SMatthew G Knepley   DM_Patch *mesh = (DM_Patch *)dm->data;
3513a19ef87SMatthew G Knepley 
3523a19ef87SMatthew G Knepley   PetscFunctionBegin;
3533a19ef87SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3549566063dSJacob Faibussowitsch   PetscCall(DMSetUp(mesh->dmCoarse));
3553a19ef87SMatthew G Knepley   PetscFunctionReturn(0);
3563a19ef87SMatthew G Knepley }
3573a19ef87SMatthew G Knepley 
358*9371c9d4SSatish Balay PetscErrorCode DMCreateGlobalVector_Patch(DM dm, Vec *g) {
3593a19ef87SMatthew G Knepley   DM_Patch *mesh = (DM_Patch *)dm->data;
3603a19ef87SMatthew G Knepley 
3613a19ef87SMatthew G Knepley   PetscFunctionBegin;
3623a19ef87SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3639566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(mesh->dmCoarse, g));
3643a19ef87SMatthew G Knepley   PetscFunctionReturn(0);
3653a19ef87SMatthew G Knepley }
3663a19ef87SMatthew G Knepley 
367*9371c9d4SSatish Balay PetscErrorCode DMCreateLocalVector_Patch(DM dm, Vec *l) {
3683a19ef87SMatthew G Knepley   DM_Patch *mesh = (DM_Patch *)dm->data;
3693a19ef87SMatthew G Knepley 
3703a19ef87SMatthew G Knepley   PetscFunctionBegin;
3713a19ef87SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3729566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(mesh->dmCoarse, l));
3733a19ef87SMatthew G Knepley   PetscFunctionReturn(0);
3743a19ef87SMatthew G Knepley }
3753a19ef87SMatthew G Knepley 
376*9371c9d4SSatish Balay PetscErrorCode DMCreateSubDM_Patch(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) {
37782f516ccSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Tell me to code this");
3783a19ef87SMatthew G Knepley }
379392fa6c6SMatthew G Knepley 
380*9371c9d4SSatish Balay PetscErrorCode DMPatchGetCoarse(DM dm, DM *dmCoarse) {
381392fa6c6SMatthew G Knepley   DM_Patch *mesh = (DM_Patch *)dm->data;
382392fa6c6SMatthew G Knepley 
383392fa6c6SMatthew G Knepley   PetscFunctionBegin;
384392fa6c6SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
385392fa6c6SMatthew G Knepley   *dmCoarse = mesh->dmCoarse;
386392fa6c6SMatthew G Knepley   PetscFunctionReturn(0);
387392fa6c6SMatthew G Knepley }
388392fa6c6SMatthew G Knepley 
389*9371c9d4SSatish Balay PetscErrorCode DMPatchGetPatchSize(DM dm, MatStencil *patchSize) {
390392fa6c6SMatthew G Knepley   DM_Patch *mesh = (DM_Patch *)dm->data;
391392fa6c6SMatthew G Knepley 
392392fa6c6SMatthew G Knepley   PetscFunctionBegin;
393392fa6c6SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
394392fa6c6SMatthew G Knepley   PetscValidPointer(patchSize, 2);
395392fa6c6SMatthew G Knepley   *patchSize = mesh->patchSize;
396392fa6c6SMatthew G Knepley   PetscFunctionReturn(0);
397392fa6c6SMatthew G Knepley }
398392fa6c6SMatthew G Knepley 
399*9371c9d4SSatish Balay PetscErrorCode DMPatchSetPatchSize(DM dm, MatStencil patchSize) {
400392fa6c6SMatthew G Knepley   DM_Patch *mesh = (DM_Patch *)dm->data;
401392fa6c6SMatthew G Knepley 
402392fa6c6SMatthew G Knepley   PetscFunctionBegin;
403392fa6c6SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
404392fa6c6SMatthew G Knepley   mesh->patchSize = patchSize;
405392fa6c6SMatthew G Knepley   PetscFunctionReturn(0);
406392fa6c6SMatthew G Knepley }
407bb71ef15SMatthew G Knepley 
408*9371c9d4SSatish Balay PetscErrorCode DMPatchGetCommSize(DM dm, MatStencil *commSize) {
409bb71ef15SMatthew G Knepley   DM_Patch *mesh = (DM_Patch *)dm->data;
410bb71ef15SMatthew G Knepley 
411bb71ef15SMatthew G Knepley   PetscFunctionBegin;
412bb71ef15SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
413bb71ef15SMatthew G Knepley   PetscValidPointer(commSize, 2);
414bb71ef15SMatthew G Knepley   *commSize = mesh->commSize;
415bb71ef15SMatthew G Knepley   PetscFunctionReturn(0);
416bb71ef15SMatthew G Knepley }
417bb71ef15SMatthew G Knepley 
418*9371c9d4SSatish Balay PetscErrorCode DMPatchSetCommSize(DM dm, MatStencil commSize) {
419bb71ef15SMatthew G Knepley   DM_Patch *mesh = (DM_Patch *)dm->data;
420bb71ef15SMatthew G Knepley 
421bb71ef15SMatthew G Knepley   PetscFunctionBegin;
422bb71ef15SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
423bb71ef15SMatthew G Knepley   mesh->commSize = commSize;
424bb71ef15SMatthew G Knepley   PetscFunctionReturn(0);
425bb71ef15SMatthew G Knepley }
426