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