Actual source code: patch.c

petsc-3.4.5 2014-06-29
  1: #include <petsc-private/dmpatchimpl.h>   /*I      "petscdmpatch.h"   I*/
  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: */

 23: /*
 24:   DMPatchZoom - Create a version of the coarse patch (identified by rank) with halo on communicator commz

 26:   Collective on DM

 28:   Input Parameters:
 29:   + dm - the DM
 30:   . rank - the rank which holds the given patch
 31:   - commz - the new communicator for the patch

 33:   Output Parameters:
 34:   + dmz  - the patch DM
 35:   . sfz  - the PetscSF mapping the patch+halo to the zoomed version
 36:   . sfzr - the PetscSF mapping the patch to the restricted zoomed version

 38:   Level: intermediate

 40:   Note: All processes in commz should have the same rank (could autosplit comm)

 42: .seealso: DMPatchSolve()
 43: */
 44: PetscErrorCode DMPatchZoom(DM dm, Vec X, MatStencil lower, MatStencil upper, MPI_Comm commz, DM *dmz, PetscSF *sfz, PetscSF *sfzr)
 45: {
 46:   DMDAStencilType st;
 47:   MatStencil      blower, bupper, loclower, locupper;
 48:   IS              is;
 49:   const PetscInt  *ranges, *indices;
 50:   PetscInt        *localPoints  = NULL;
 51:   PetscSFNode     *remotePoints = NULL;
 52:   PetscInt        dim, dof;
 53:   PetscInt        M, N, P, rM, rN, rP, halo = 1, sxb, syb, szb, sxr, syr, szr, exr, eyr, ezr, mxb, myb, mzb, i, j, k, q;
 54:   PetscMPIInt     size;
 55:   PetscErrorCode  ierr;

 58:   MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);
 59:   /* Create patch DM */
 60:   DMDAGetInfo(dm, &dim, &M, &N, &P, 0,0,0, &dof, 0,0,0,0, &st);

 62:   /* Get piece for rank r, expanded by halo */
 63:   bupper.i = PetscMin(M, upper.i + halo); blower.i = PetscMax(lower.i - halo, 0);
 64:   bupper.j = PetscMin(N, upper.j + halo); blower.j = PetscMax(lower.j - halo, 0);
 65:   bupper.k = PetscMin(P, upper.k + halo); blower.k = PetscMax(lower.k - halo, 0);
 66:   rM       = bupper.i - blower.i;
 67:   rN       = bupper.j - blower.j;
 68:   rP       = bupper.k - blower.k;

 70:   if (commz != MPI_COMM_NULL) {
 71:     DMDACreate(commz, dmz);
 72:     DMDASetDim(*dmz, dim);
 73:     DMDASetSizes(*dmz, rM, rN, rP);
 74:     DMDASetNumProcs(*dmz, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE);
 75:     DMDASetBoundaryType(*dmz, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE);
 76:     DMDASetDof(*dmz, dof);
 77:     DMDASetStencilType(*dmz, st);
 78:     DMDASetStencilWidth(*dmz, 0);
 79:     DMDASetOwnershipRanges(*dmz, NULL, NULL, NULL);
 80:     DMSetFromOptions(*dmz);
 81:     DMSetUp(*dmz);
 82:     DMDAGetCorners(*dmz, &sxb, &syb, &szb, &mxb, &myb, &mzb);
 83:     sxr  = PetscMax(sxb,     lower.i - blower.i);
 84:     syr  = PetscMax(syb,     lower.j - blower.j);
 85:     szr  = PetscMax(szb,     lower.k - blower.k);
 86:     exr  = PetscMin(sxb+mxb, upper.i - blower.i);
 87:     eyr  = PetscMin(syb+myb, upper.j - blower.j);
 88:     ezr  = PetscMin(szb+mzb, upper.k - blower.k);
 89:     PetscMalloc2(rM*rN*rP,PetscInt,&localPoints,rM*rN*rP,PetscSFNode,&remotePoints);
 90:   } else {
 91:     sxr = syr = szr = exr = eyr = ezr = sxb = syb = szb = mxb = myb = mzb = 0;
 92:   }

 94:   /* Create SF for restricted map */
 95:   VecGetOwnershipRanges(X,&ranges);

 97:   loclower.i = blower.i + sxr; locupper.i = blower.i + exr;
 98:   loclower.j = blower.j + syr; locupper.j = blower.j + eyr;
 99:   loclower.k = blower.k + szr; locupper.k = blower.k + ezr;

101:   DMDACreatePatchIS(dm, &loclower, &locupper, &is);
102:   ISGetIndices(is, &indices);

104:   q = 0;
105:   for (k = szb; k < szb+mzb; ++k) {
106:     if ((k < szr) || (k >= ezr)) continue;
107:     for (j = syb; j < syb+myb; ++j) {
108:       if ((j < syr) || (j >= eyr)) continue;
109:       for (i = sxb; i < sxb+mxb; ++i) {
110:         const PetscInt lp = ((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:   ISRestoreIndices(is, &indices);
124:   ISDestroy(&is);
125:   PetscSFCreate(PetscObjectComm((PetscObject)dm), sfzr);
126:   PetscObjectSetName((PetscObject) *sfzr, "Restricted Map");
127:   PetscSFSetGraph(*sfzr, M*N*P, q, localPoints, PETSC_COPY_VALUES, remotePoints, PETSC_COPY_VALUES);

129:   /* Create SF for buffered map */
130:   loclower.i = blower.i + sxb; locupper.i = blower.i + sxb+mxb;
131:   loclower.j = blower.j + syb; locupper.j = blower.j + syb+myb;
132:   loclower.k = blower.k + szb; locupper.k = blower.k + szb+mzb;

134:   DMDACreatePatchIS(dm, &loclower, &locupper, &is);
135:   ISGetIndices(is, &indices);

137:   q = 0;
138:   for (k = szb; k < szb+mzb; ++k) {
139:     for (j = syb; j < syb+myb; ++j) {
140:       for (i = sxb; i < sxb+mxb; ++i, ++q) {
141:         PetscInt r;

143:         localPoints[q]        = q;
144:         PetscFindInt(indices[q], size+1, ranges, &r);
145:         remotePoints[q].rank  = r < 0 ? -(r+1) - 1 : r;
146:         remotePoints[q].index = indices[q] - ranges[remotePoints[q].rank];
147:       }
148:     }
149:   }
150:   ISRestoreIndices(is, &indices);
151:   ISDestroy(&is);
152:   PetscSFCreate(PetscObjectComm((PetscObject)dm), sfz);
153:   PetscObjectSetName((PetscObject) *sfz, "Buffered Map");
154:   PetscSFSetGraph(*sfz, M*N*P, q, localPoints, PETSC_COPY_VALUES, remotePoints, PETSC_COPY_VALUES);

156:   PetscFree2(localPoints, remotePoints);
157:   return(0);
158: }

160: typedef enum {PATCH_COMM_TYPE_WORLD = 0, PATCH_COMM_TYPE_SELF = 1} PatchCommType;

164: PetscErrorCode DMPatchSolve(DM dm)
165: {
166:   MPI_Comm       comm;
167:   MPI_Comm       commz;
168:   DM             dmc;
169:   PetscSF        sfz, sfzr;
170:   Vec            XC;
171:   MatStencil     patchSize, commSize, gridRank, lower, upper;
172:   PetscInt       M, N, P, i, j, k, l, m, n, p = 0;
173:   PetscMPIInt    rank, size;
174:   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, 0, &M, &N, &P, &l, &m, &n, 0,0,0,0,0,0);
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 %d\n", rank, newComm, newRank, *((PetscMPIInt*) &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, XC, 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:         /* Debug */
244:         PetscPrintf(comm, "Patch %d: (%d, %d, %d)--(%d, %d, %d)\n", p, lower.i, lower.j, lower.k, upper.i, upper.j, upper.k);
245:         if (dmz) {DMView(dmz, PETSC_VIEWER_STDOUT_(commz));}
246:         PetscSFView(sfz,  PETSC_VIEWER_STDOUT_(comm));
247:         PetscSFView(sfzr, PETSC_VIEWER_STDOUT_(comm));
248:         /* Scatter Xcoarse -> Xzoom */
249:         if (dmz) {DMGetGlobalVector(dmz, &XZ);}
250:         if (XZ)  {VecGetArray(XZ, &xzarray);}
251:         VecGetArray(XC, &xcarray);
252:         PetscSFBcastBegin(sfz, MPIU_SCALAR, xcarray, xzarray);
253:         PetscSFBcastEnd(sfz, MPIU_SCALAR, xcarray, xzarray);
254:         VecRestoreArray(XC, &xcarray);
255:         if (XZ)  {VecRestoreArray(XZ, &xzarray);}
256: #if 0
257:         /* Interpolate Xzoom -> Xfine, note that this may be on subcomms */
258:         DMRefine(dmz, MPI_COMM_NULL, &dmf);
259:         DMCreateInterpolation(dmz, dmf, &interpz, NULL);
260:         DMInterpolate(dmz, interpz, dmf);
261:         /* Smooth Xfine using two-step smoother, normal smoother plus Kaczmarz---moves back and forth from dmzoom to dmfine */
262:         /* Compute residual Rfine */
263:         /* Restrict Rfine to Rzoom_restricted */
264: #endif
265:         /* Scatter Rzoom_restricted -> Rcoarse_restricted */
266:         if (XZ)  {VecGetArray(XZ, &xzarray);}
267:         VecGetArray(XC, &xcarray);
268:         PetscSFReduceBegin(sfzr, MPIU_SCALAR, xzarray, xcarray, MPI_SUM);
269:         PetscSFReduceEnd(sfzr, MPIU_SCALAR, xzarray, xcarray, MPI_SUM);
270:         VecRestoreArray(XC, &xcarray);
271:         if (XZ)  {VecRestoreArray(XZ, &xzarray);}
272:         if (dmz) {DMRestoreGlobalVector(dmz, &XZ);}
273:         /* Compute global residual Rcoarse */
274:         /* TauCoarse = Rcoarse - Rcoarse_restricted */

276:         PetscSFDestroy(&sfz);
277:         PetscSFDestroy(&sfzr);
278:         DMDestroy(&dmz);
279:       }
280:     }
281:   }
282:   DMRestoreGlobalVector(dmc, &XC);
283:   return(0);
284: }

288: PetscErrorCode DMPatchView_Ascii(DM dm, PetscViewer viewer)
289: {
290:   DM_Patch          *mesh = (DM_Patch*) dm->data;
291:   PetscViewerFormat format;
292:   const char        *name;
293:   PetscErrorCode    ierr;

296:   PetscViewerGetFormat(viewer, &format);
297:   /* if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) */
298:   PetscObjectGetName((PetscObject) dm, &name);
299:   PetscViewerASCIIPrintf(viewer, "Patch DM %s\n", name);
300:   PetscViewerASCIIPushTab(viewer);
301:   PetscViewerASCIIPrintf(viewer, "Coarse DM\n");
302:   DMView(mesh->dmCoarse, viewer);
303:   PetscViewerASCIIPopTab(viewer);
304:   return(0);
305: }

309: PetscErrorCode DMView_Patch(DM dm, PetscViewer viewer)
310: {
311:   PetscBool      iascii, isbinary;

317:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
318:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERBINARY, &isbinary);
319:   if (iascii) {
320:     DMPatchView_Ascii(dm, viewer);
321: #if 0
322:   } else if (isbinary) {
323:     DMPatchView_Binary(dm, viewer);
324: #endif
325:   }
326:   return(0);
327: }

331: PetscErrorCode DMDestroy_Patch(DM dm)
332: {
333:   DM_Patch       *mesh = (DM_Patch*) dm->data;

337:   if (--mesh->refct > 0) return(0);
338:   DMDestroy(&mesh->dmCoarse);
339:   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
340:   PetscFree(mesh);
341:   return(0);
342: }

346: PetscErrorCode DMSetUp_Patch(DM dm)
347: {
348:   DM_Patch       *mesh = (DM_Patch*) dm->data;

353:   DMSetUp(mesh->dmCoarse);
354:   return(0);
355: }

359: PetscErrorCode DMCreateGlobalVector_Patch(DM dm, Vec *g)
360: {
361:   DM_Patch       *mesh = (DM_Patch*) dm->data;

366:   DMCreateGlobalVector(mesh->dmCoarse, g);
367:   return(0);
368: }

372: PetscErrorCode DMCreateLocalVector_Patch(DM dm, Vec *l)
373: {
374:   DM_Patch       *mesh = (DM_Patch*) dm->data;

379:   DMCreateLocalVector(mesh->dmCoarse, l);
380:   return(0);
381: }

385: PetscErrorCode DMCreateSubDM_Patch(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
386: {
389:   SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Tell me to code this");
390:   return(0);
391: }

395: PetscErrorCode DMPatchGetCoarse(DM dm, DM *dmCoarse)
396: {
397:   DM_Patch *mesh = (DM_Patch*) dm->data;

401:   *dmCoarse = mesh->dmCoarse;
402:   return(0);
403: }

407: PetscErrorCode DMPatchGetPatchSize(DM dm, MatStencil *patchSize)
408: {
409:   DM_Patch *mesh = (DM_Patch*) dm->data;

414:   *patchSize = mesh->patchSize;
415:   return(0);
416: }

420: PetscErrorCode DMPatchSetPatchSize(DM dm, MatStencil patchSize)
421: {
422:   DM_Patch *mesh = (DM_Patch*) dm->data;

426:   mesh->patchSize = patchSize;
427:   return(0);
428: }

432: PetscErrorCode DMPatchGetCommSize(DM dm, MatStencil *commSize)
433: {
434:   DM_Patch *mesh = (DM_Patch*) dm->data;

439:   *commSize = mesh->commSize;
440:   return(0);
441: }

445: PetscErrorCode DMPatchSetCommSize(DM dm, MatStencil commSize)
446: {
447:   DM_Patch *mesh = (DM_Patch*) dm->data;

451:   mesh->commSize = commSize;
452:   return(0);
453: }