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