Actual source code: patch.c
1: #include <petsc/private/dmpatchimpl.h>
2: #include <petscdmda.h>
3: #include <petscsf.h>
5: /*
6: Solver loop to update \tau:
8: DMZoom(dmc, &dmz)
9: DMRefine(dmz, &dmf),
10: Scatter Xcoarse -> Xzoom,
11: Interpolate Xzoom -> Xfine (note that this may be on subcomms),
12: Smooth Xfine using two-step smoother
13: normal smoother plus Kaczmarz---moves back and forth from dmzoom to dmfine
14: Compute residual Rfine
15: Restrict Rfine to Rzoom_restricted
16: Scatter Rzoom_restricted -> Rcoarse_restricted
17: Compute global residual Rcoarse
18: TauCoarse = Rcoarse - Rcoarse_restricted
19: */
21: /*@C
22: DMPatchZoom - Create patches of a DMDA on subsets of processes, indicated by commz
24: Collective on dm
26: Input Parameters:
27: + dm - the DM
28: . lower,upper - the upper right corner and the lower left corner of the requested patch
29: - commz - the new communicator for the patch, MPI_COMM_NULL indicates that the given rank will not own a patch
31: Output Parameters:
32: + dmz - the patch DM
33: . sfz - the PetscSF mapping the patch+halo to the zoomed version (optional)
34: - sfzr - the PetscSF mapping the patch to the restricted zoomed version
36: Level: intermediate
38: .seealso: DMPatchSolve(), DMDACreatePatchIS()
39: @*/
40: PetscErrorCode DMPatchZoom(DM dm, MatStencil lower, MatStencil upper, MPI_Comm commz, DM *dmz, PetscSF *sfz, PetscSF *sfzr)
41: {
42: DMDAStencilType st;
43: MatStencil blower, bupper, loclower, locupper;
44: IS is;
45: const PetscInt *ranges, *indices;
46: PetscInt *localPoints = NULL;
47: PetscSFNode *remotePoints = NULL;
48: PetscInt dim, dof;
49: 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;
50: PetscMPIInt size;
51: PetscBool patchis_offproc = PETSC_TRUE;
52: Vec X;
54: if (!sfz) halo = 0;
55: MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);
56: /* Create patch DM */
57: DMDAGetInfo(dm, &dim, &M, &N, &P, NULL,NULL,NULL, &dof, NULL,NULL,NULL,NULL, &st);
59: /* Get piece for rank r, expanded by halo */
60: bupper.i = PetscMin(M, upper.i + halo); blower.i = PetscMax(lower.i - halo, 0);
61: bupper.j = PetscMin(N, upper.j + halo); blower.j = PetscMax(lower.j - halo, 0);
62: bupper.k = PetscMin(P, upper.k + halo); blower.k = PetscMax(lower.k - halo, 0);
63: rM = bupper.i - blower.i;
64: rN = bupper.j - blower.j;
65: rP = bupper.k - blower.k;
67: if (commz != MPI_COMM_NULL) {
68: DMDACreate(commz, dmz);
69: DMSetDimension(*dmz, dim);
70: DMDASetSizes(*dmz, rM, rN, rP);
71: DMDASetNumProcs(*dmz, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE);
72: DMDASetBoundaryType(*dmz, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE);
73: DMDASetDof(*dmz, dof);
74: DMDASetStencilType(*dmz, st);
75: DMDASetStencilWidth(*dmz, 0);
76: DMDASetOwnershipRanges(*dmz, NULL, NULL, NULL);
77: DMSetFromOptions(*dmz);
78: DMSetUp(*dmz);
79: DMDAGetCorners(*dmz, &sxb, &syb, &szb, &mxb, &myb, &mzb);
80: sxr = PetscMax(sxb, lower.i - blower.i);
81: syr = PetscMax(syb, lower.j - blower.j);
82: szr = PetscMax(szb, lower.k - blower.k);
83: exr = PetscMin(sxb+mxb, upper.i - blower.i);
84: eyr = PetscMin(syb+myb, upper.j - blower.j);
85: ezr = PetscMin(szb+mzb, upper.k - blower.k);
86: PetscMalloc2(dof*rM*rN*PetscMax(rP,1),&localPoints,dof*rM*rN*PetscMax(rP,1),&remotePoints);
87: } else {
88: sxr = syr = szr = exr = eyr = ezr = sxb = syb = szb = mxb = myb = mzb = 0;
89: }
91: /* Create SF for restricted map */
92: DMCreateGlobalVector(dm,&X);
93: VecGetOwnershipRanges(X,&ranges);
95: loclower.i = blower.i + sxr; locupper.i = blower.i + exr;
96: loclower.j = blower.j + syr; locupper.j = blower.j + eyr;
97: loclower.k = blower.k + szr; locupper.k = blower.k + ezr;
99: DMDACreatePatchIS(dm, &loclower, &locupper, &is, patchis_offproc);
100: ISGetIndices(is, &indices);
102: if (dim < 3) {mzb = 1; ezr = 1;}
103: q = 0;
104: for (k = szb; k < szb+mzb; ++k) {
105: if ((k < szr) || (k >= ezr)) continue;
106: for (j = syb; j < syb+myb; ++j) {
107: if ((j < syr) || (j >= eyr)) continue;
108: for (i = sxb; i < sxb+mxb; ++i) {
109: for (l=0; l<dof; l++) {
110: const PetscInt lp = l + dof*(((k-szb)*rN + (j-syb))*rM + i-sxb);
111: PetscInt r;
113: if ((i < sxr) || (i >= exr)) continue;
114: localPoints[q] = lp;
115: PetscFindInt(indices[q], size+1, ranges, &r);
117: remotePoints[q].rank = r < 0 ? -(r+1) - 1 : r;
118: remotePoints[q].index = indices[q] - ranges[remotePoints[q].rank];
119: ++q;
120: }
121: }
122: }
123: }
124: ISRestoreIndices(is, &indices);
125: ISDestroy(&is);
126: PetscSFCreate(PetscObjectComm((PetscObject)dm), sfzr);
127: PetscObjectSetName((PetscObject) *sfzr, "Restricted Map");
128: PetscSFSetGraph(*sfzr, dof*M*N*P, q, localPoints, PETSC_COPY_VALUES, remotePoints, PETSC_COPY_VALUES);
130: if (sfz) {
131: /* Create SF for buffered map */
132: loclower.i = blower.i + sxb; locupper.i = blower.i + sxb+mxb;
133: loclower.j = blower.j + syb; locupper.j = blower.j + syb+myb;
134: loclower.k = blower.k + szb; locupper.k = blower.k + szb+mzb;
136: DMDACreatePatchIS(dm, &loclower, &locupper, &is, patchis_offproc);
137: ISGetIndices(is, &indices);
139: q = 0;
140: for (k = szb; k < szb+mzb; ++k) {
141: for (j = syb; j < syb+myb; ++j) {
142: for (i = sxb; i < sxb+mxb; ++i, ++q) {
143: PetscInt r;
145: localPoints[q] = q;
146: PetscFindInt(indices[q], size+1, ranges, &r);
147: remotePoints[q].rank = r < 0 ? -(r+1) - 1 : r;
148: remotePoints[q].index = indices[q] - ranges[remotePoints[q].rank];
149: }
150: }
151: }
152: ISRestoreIndices(is, &indices);
153: ISDestroy(&is);
154: PetscSFCreate(PetscObjectComm((PetscObject)dm), sfz);
155: PetscObjectSetName((PetscObject) *sfz, "Buffered Map");
156: PetscSFSetGraph(*sfz, M*N*P, q, localPoints, PETSC_COPY_VALUES, remotePoints, PETSC_COPY_VALUES);
157: }
159: VecDestroy(&X);
160: PetscFree2(localPoints, remotePoints);
161: return 0;
162: }
164: typedef enum {PATCH_COMM_TYPE_WORLD = 0, PATCH_COMM_TYPE_SELF = 1} PatchCommType;
166: PetscErrorCode DMPatchSolve(DM dm)
167: {
168: MPI_Comm comm;
169: MPI_Comm commz;
170: DM dmc;
171: PetscSF sfz, sfzr;
172: Vec XC;
173: MatStencil patchSize, commSize, gridRank, lower, upper;
174: PetscInt M, N, P, i, j, k, l, m, n, p = 0;
175: PetscMPIInt rank, size;
176: PetscInt debug = 0;
178: PetscObjectGetComm((PetscObject)dm,&comm);
179: MPI_Comm_rank(comm, &rank);
180: MPI_Comm_size(comm, &size);
181: DMPatchGetCoarse(dm, &dmc);
182: DMPatchGetPatchSize(dm, &patchSize);
183: DMPatchGetCommSize(dm, &commSize);
184: DMPatchGetCommSize(dm, &commSize);
185: DMGetGlobalVector(dmc, &XC);
186: DMDAGetInfo(dmc, NULL, &M, &N, &P, &l, &m, &n, NULL,NULL,NULL,NULL,NULL,NULL);
187: M = PetscMax(M, 1); l = PetscMax(l, 1);
188: N = PetscMax(N, 1); m = PetscMax(m, 1);
189: P = PetscMax(P, 1); n = PetscMax(n, 1);
191: gridRank.i = rank % l;
192: gridRank.j = rank/l % m;
193: gridRank.k = rank/(l*m) % n;
195: if (commSize.i*commSize.j*commSize.k == size || commSize.i*commSize.j*commSize.k == 0) {
196: commSize.i = l; commSize.j = m; commSize.k = n;
197: commz = comm;
198: } else if (commSize.i*commSize.j*commSize.k == 1) {
199: commz = PETSC_COMM_SELF;
200: } else {
201: const PetscMPIInt newComm = ((gridRank.k/commSize.k)*(m/commSize.j) + gridRank.j/commSize.j)*(l/commSize.i) + (gridRank.i/commSize.i);
202: const PetscMPIInt newRank = ((gridRank.k%commSize.k)*commSize.j + gridRank.j%commSize.j)*commSize.i + (gridRank.i%commSize.i);
204: MPI_Comm_split(comm, newComm, newRank, &commz);
205: if (debug) PetscPrintf(PETSC_COMM_SELF, "Rank %d color %d key %d commz %p\n", rank, newComm, newRank, (void*)(MPI_Aint)commz);
206: }
207: /*
208: Assumptions:
209: - patchSize divides gridSize
210: - commSize divides gridSize
211: - commSize divides l,m,n
212: Ignore multiple patches per rank for now
214: Multiple ranks per patch:
215: - l,m,n divides patchSize
216: - commSize divides patchSize
217: */
218: for (k = 0; k < P; k += PetscMax(patchSize.k, 1)) {
219: for (j = 0; j < N; j += PetscMax(patchSize.j, 1)) {
220: for (i = 0; i < M; i += PetscMax(patchSize.i, 1), ++p) {
221: MPI_Comm commp = MPI_COMM_NULL;
222: DM dmz = NULL;
223: #if 0
224: DM dmf = NULL;
225: Mat interpz = NULL;
226: #endif
227: Vec XZ = NULL;
228: PetscScalar *xcarray = NULL;
229: PetscScalar *xzarray = NULL;
231: if ((gridRank.k/commSize.k == p/(l/commSize.i * m/commSize.j) % n/commSize.k) &&
232: (gridRank.j/commSize.j == p/(l/commSize.i) % m/commSize.j) &&
233: (gridRank.i/commSize.i == p % l/commSize.i)) {
234: if (debug) PetscPrintf(PETSC_COMM_SELF, "Rank %d is accepting Patch %d\n", rank, p);
235: commp = commz;
236: }
237: /* Zoom to coarse patch */
238: lower.i = i; lower.j = j; lower.k = k;
239: upper.i = i + patchSize.i; upper.j = j + patchSize.j; upper.k = k + patchSize.k;
240: DMPatchZoom(dmc, lower, upper, commp, &dmz, &sfz, &sfzr);
241: lower.c = 0; /* initialize member, otherwise compiler issues warnings */
242: upper.c = 0; /* initialize member, otherwise compiler issues warnings */
243: if (debug) PetscPrintf(comm, "Patch %d: (%d, %d, %d)--(%d, %d, %d)\n", p, lower.i, lower.j, lower.k, upper.i, upper.j, upper.k);
244: if (dmz) DMView(dmz, PETSC_VIEWER_STDOUT_(commz));
245: PetscSFView(sfz, PETSC_VIEWER_STDOUT_(comm));
246: PetscSFView(sfzr, PETSC_VIEWER_STDOUT_(comm));
247: /* Scatter Xcoarse -> Xzoom */
248: if (dmz) DMGetGlobalVector(dmz, &XZ);
249: if (XZ) VecGetArray(XZ, &xzarray);
250: VecGetArray(XC, &xcarray);
251: PetscSFBcastBegin(sfz, MPIU_SCALAR, xcarray, xzarray,MPI_REPLACE);
252: PetscSFBcastEnd(sfz, MPIU_SCALAR, xcarray, xzarray,MPI_REPLACE);
253: VecRestoreArray(XC, &xcarray);
254: if (XZ) VecRestoreArray(XZ, &xzarray);
255: #if 0
256: /* Interpolate Xzoom -> Xfine, note that this may be on subcomms */
257: DMRefine(dmz, MPI_COMM_NULL, &dmf);
258: DMCreateInterpolation(dmz, dmf, &interpz, NULL);
259: DMInterpolate(dmz, interpz, dmf);
260: /* Smooth Xfine using two-step smoother, normal smoother plus Kaczmarz---moves back and forth from dmzoom to dmfine */
261: /* Compute residual Rfine */
262: /* Restrict Rfine to Rzoom_restricted */
263: #endif
264: /* Scatter Rzoom_restricted -> Rcoarse_restricted */
265: if (XZ) VecGetArray(XZ, &xzarray);
266: VecGetArray(XC, &xcarray);
267: PetscSFReduceBegin(sfzr, MPIU_SCALAR, xzarray, xcarray, MPIU_SUM);
268: PetscSFReduceEnd(sfzr, MPIU_SCALAR, xzarray, xcarray, MPIU_SUM);
269: VecRestoreArray(XC, &xcarray);
270: if (XZ) VecRestoreArray(XZ, &xzarray);
271: if (dmz) DMRestoreGlobalVector(dmz, &XZ);
272: /* Compute global residual Rcoarse */
273: /* TauCoarse = Rcoarse - Rcoarse_restricted */
275: PetscSFDestroy(&sfz);
276: PetscSFDestroy(&sfzr);
277: DMDestroy(&dmz);
278: }
279: }
280: }
281: DMRestoreGlobalVector(dmc, &XC);
282: return 0;
283: }
285: PetscErrorCode DMPatchView_ASCII(DM dm, PetscViewer viewer)
286: {
287: DM_Patch *mesh = (DM_Patch*) dm->data;
288: PetscViewerFormat format;
289: const char *name;
291: PetscViewerGetFormat(viewer, &format);
292: /* if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) */
293: PetscObjectGetName((PetscObject) dm, &name);
294: PetscViewerASCIIPrintf(viewer, "Patch DM %s\n", name);
295: PetscViewerASCIIPushTab(viewer);
296: PetscViewerASCIIPrintf(viewer, "Coarse DM\n");
297: DMView(mesh->dmCoarse, viewer);
298: PetscViewerASCIIPopTab(viewer);
299: return 0;
300: }
302: PetscErrorCode DMView_Patch(DM dm, PetscViewer viewer)
303: {
304: PetscBool iascii, isbinary;
308: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
309: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERBINARY, &isbinary);
310: if (iascii) {
311: DMPatchView_ASCII(dm, viewer);
312: }
313: return 0;
314: }
316: PetscErrorCode DMDestroy_Patch(DM dm)
317: {
318: DM_Patch *mesh = (DM_Patch*) dm->data;
320: if (--mesh->refct > 0) return 0;
321: DMDestroy(&mesh->dmCoarse);
322: /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
323: PetscFree(mesh);
324: return 0;
325: }
327: PetscErrorCode DMSetUp_Patch(DM dm)
328: {
329: DM_Patch *mesh = (DM_Patch*) dm->data;
332: DMSetUp(mesh->dmCoarse);
333: return 0;
334: }
336: PetscErrorCode DMCreateGlobalVector_Patch(DM dm, Vec *g)
337: {
338: DM_Patch *mesh = (DM_Patch*) dm->data;
341: DMCreateGlobalVector(mesh->dmCoarse, g);
342: return 0;
343: }
345: PetscErrorCode DMCreateLocalVector_Patch(DM dm, Vec *l)
346: {
347: DM_Patch *mesh = (DM_Patch*) dm->data;
350: DMCreateLocalVector(mesh->dmCoarse, l);
351: return 0;
352: }
354: PetscErrorCode DMCreateSubDM_Patch(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
355: {
356: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Tell me to code this");
357: }
359: PetscErrorCode DMPatchGetCoarse(DM dm, DM *dmCoarse)
360: {
361: DM_Patch *mesh = (DM_Patch*) dm->data;
364: *dmCoarse = mesh->dmCoarse;
365: return 0;
366: }
368: PetscErrorCode DMPatchGetPatchSize(DM dm, MatStencil *patchSize)
369: {
370: DM_Patch *mesh = (DM_Patch*) dm->data;
374: *patchSize = mesh->patchSize;
375: return 0;
376: }
378: PetscErrorCode DMPatchSetPatchSize(DM dm, MatStencil patchSize)
379: {
380: DM_Patch *mesh = (DM_Patch*) dm->data;
383: mesh->patchSize = patchSize;
384: return 0;
385: }
387: PetscErrorCode DMPatchGetCommSize(DM dm, MatStencil *commSize)
388: {
389: DM_Patch *mesh = (DM_Patch*) dm->data;
393: *commSize = mesh->commSize;
394: return 0;
395: }
397: PetscErrorCode DMPatchSetCommSize(DM dm, MatStencil commSize)
398: {
399: DM_Patch *mesh = (DM_Patch*) dm->data;
402: mesh->commSize = commSize;
403: return 0;
404: }