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