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 = NULL; 50 PetscSFNode *remotePoints = 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, 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,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 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);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(((PetscObject) dm)->comm, 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);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(((PetscObject) dm)->comm, 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 #undef __FUNCT__ 162 #define __FUNCT__ "DMPatchSolve" 163 PetscErrorCode DMPatchSolve(DM dm) 164 { 165 MPI_Comm comm = ((PetscObject) dm)->comm; 166 MPI_Comm commz; 167 DM dmc; 168 PetscSF sfz, sfzr; 169 Vec XC; 170 MatStencil patchSize, commSize, gridRank, lower, upper; 171 PetscInt M, N, P, i, j, k, l, m, n, p = 0; 172 PetscMPIInt rank, size; 173 PetscInt debug = 0; 174 PetscErrorCode ierr; 175 176 PetscFunctionBegin; 177 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 178 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 179 ierr = DMPatchGetCoarse(dm, &dmc);CHKERRQ(ierr); 180 ierr = DMPatchGetPatchSize(dm, &patchSize);CHKERRQ(ierr); 181 ierr = DMPatchGetCommSize(dm, &commSize);CHKERRQ(ierr); 182 ierr = DMPatchGetCommSize(dm, &commSize);CHKERRQ(ierr); 183 ierr = DMGetGlobalVector(dmc, &XC);CHKERRQ(ierr); 184 ierr = DMDAGetInfo(dmc, 0, &M, &N, &P, &l, &m, &n, 0,0,0,0,0,0);CHKERRQ(ierr); 185 M = PetscMax(M, 1); l = PetscMax(l, 1); 186 N = PetscMax(N, 1); m = PetscMax(m, 1); 187 P = PetscMax(P, 1); n = PetscMax(n, 1); 188 189 gridRank.i = rank % l; 190 gridRank.j = rank/l % m; 191 gridRank.k = rank/(l*m) % n; 192 193 if (commSize.i*commSize.j*commSize.k == size || commSize.i*commSize.j*commSize.k == 0) { 194 commSize.i = l; commSize.j = m; commSize.k = n; 195 commz = comm; 196 } else if (commSize.i*commSize.j*commSize.k == 1) { 197 commz = PETSC_COMM_SELF; 198 } else { 199 const PetscMPIInt newComm = ((gridRank.k/commSize.k)*(m/commSize.j) + gridRank.j/commSize.j)*(l/commSize.i) + (gridRank.i/commSize.i); 200 const PetscMPIInt newRank = ((gridRank.k%commSize.k)*commSize.j + gridRank.j%commSize.j)*commSize.i + (gridRank.i%commSize.i); 201 202 ierr = MPI_Comm_split(comm, newComm, newRank, &commz);CHKERRQ(ierr); 203 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "Rank %d color %d key %d commz %d\n", rank, newComm, newRank, *((PetscMPIInt*) &commz));CHKERRQ(ierr);} 204 } 205 /* 206 Assumptions: 207 - patchSize divides gridSize 208 - commSize divides gridSize 209 - commSize divides l,m,n 210 Ignore multiple patches per rank for now 211 212 Multiple ranks per patch: 213 - l,m,n divides patchSize 214 - commSize divides patchSize 215 */ 216 for (k = 0; k < P; k += PetscMax(patchSize.k, 1)) { 217 for (j = 0; j < N; j += PetscMax(patchSize.j, 1)) { 218 for (i = 0; i < M; i += PetscMax(patchSize.i, 1), ++p) { 219 MPI_Comm commp = MPI_COMM_NULL; 220 DM dmz = NULL; 221 #if 0 222 DM dmf = NULL; 223 Mat interpz = NULL; 224 #endif 225 Vec XZ = NULL; 226 PetscScalar *xcarray = NULL; 227 PetscScalar *xzarray = NULL; 228 229 if ((gridRank.k/commSize.k == p/(l/commSize.i * m/commSize.j) % n/commSize.k) && 230 (gridRank.j/commSize.j == p/(l/commSize.i) % m/commSize.j) && 231 (gridRank.i/commSize.i == p % l/commSize.i)) { 232 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "Rank %d is accepting Patch %d\n", rank, p);CHKERRQ(ierr);} 233 commp = commz; 234 } 235 /* Zoom to coarse patch */ 236 lower.i = i; lower.j = j; lower.k = k; 237 upper.i = i + patchSize.i; upper.j = j + patchSize.j; upper.k = k + patchSize.k; 238 ierr = DMPatchZoom(dmc, XC, lower, upper, commp, &dmz, &sfz, &sfzr);CHKERRQ(ierr); 239 lower.c = 0; /* initialize member, otherwise compiler issues warnings */ 240 upper.c = 0; /* initialize member, otherwise compiler issues warnings */ 241 /* Debug */ 242 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); 243 if (dmz) {ierr = DMView(dmz, PETSC_VIEWER_STDOUT_(commz));CHKERRQ(ierr);} 244 ierr = PetscSFView(sfz, PETSC_VIEWER_STDOUT_(comm));CHKERRQ(ierr); 245 ierr = PetscSFView(sfzr, PETSC_VIEWER_STDOUT_(comm));CHKERRQ(ierr); 246 /* Scatter Xcoarse -> Xzoom */ 247 if (dmz) {ierr = DMGetGlobalVector(dmz, &XZ);CHKERRQ(ierr);} 248 if (XZ) {ierr = VecGetArray(XZ, &xzarray);CHKERRQ(ierr);} 249 ierr = VecGetArray(XC, &xcarray);CHKERRQ(ierr); 250 ierr = PetscSFBcastBegin(sfz, MPIU_SCALAR, xcarray, xzarray);CHKERRQ(ierr); 251 ierr = PetscSFBcastEnd(sfz, MPIU_SCALAR, xcarray, xzarray);CHKERRQ(ierr); 252 ierr = VecRestoreArray(XC, &xcarray);CHKERRQ(ierr); 253 if (XZ) {ierr = VecRestoreArray(XZ, &xzarray);CHKERRQ(ierr);} 254 #if 0 255 /* Interpolate Xzoom -> Xfine, note that this may be on subcomms */ 256 ierr = DMRefine(dmz, MPI_COMM_NULL, &dmf);CHKERRQ(ierr); 257 ierr = DMCreateInterpolation(dmz, dmf, &interpz, NULL);CHKERRQ(ierr); 258 ierr = DMInterpolate(dmz, interpz, dmf);CHKERRQ(ierr); 259 /* Smooth Xfine using two-step smoother, normal smoother plus Kaczmarz---moves back and forth from dmzoom to dmfine */ 260 /* Compute residual Rfine */ 261 /* Restrict Rfine to Rzoom_restricted */ 262 #endif 263 /* Scatter Rzoom_restricted -> Rcoarse_restricted */ 264 if (XZ) {ierr = VecGetArray(XZ, &xzarray);CHKERRQ(ierr);} 265 ierr = VecGetArray(XC, &xcarray);CHKERRQ(ierr); 266 ierr = PetscSFReduceBegin(sfzr, MPIU_SCALAR, xzarray, xcarray, MPI_SUM);CHKERRQ(ierr); 267 ierr = PetscSFReduceEnd(sfzr, MPIU_SCALAR, xzarray, xcarray, MPI_SUM);CHKERRQ(ierr); 268 ierr = VecRestoreArray(XC, &xcarray);CHKERRQ(ierr); 269 if (XZ) {ierr = VecRestoreArray(XZ, &xzarray);CHKERRQ(ierr);} 270 if (dmz) {ierr = DMRestoreGlobalVector(dmz, &XZ);CHKERRQ(ierr);} 271 /* Compute global residual Rcoarse */ 272 /* TauCoarse = Rcoarse - Rcoarse_restricted */ 273 274 ierr = PetscSFDestroy(&sfz);CHKERRQ(ierr); 275 ierr = PetscSFDestroy(&sfzr);CHKERRQ(ierr); 276 ierr = DMDestroy(&dmz);CHKERRQ(ierr); 277 } 278 } 279 } 280 ierr = DMRestoreGlobalVector(dmc, &XC);CHKERRQ(ierr); 281 PetscFunctionReturn(0); 282 } 283 284 #undef __FUNCT__ 285 #define __FUNCT__ "DMPatchView_Ascii" 286 PetscErrorCode DMPatchView_Ascii(DM dm, PetscViewer viewer) 287 { 288 DM_Patch *mesh = (DM_Patch*) dm->data; 289 PetscViewerFormat format; 290 const char *name; 291 PetscErrorCode ierr; 292 293 PetscFunctionBegin; 294 ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 295 /* if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) */ 296 ierr = PetscObjectGetName((PetscObject) dm, &name);CHKERRQ(ierr); 297 ierr = PetscViewerASCIIPrintf(viewer, "Patch DM %s\n", name);CHKERRQ(ierr); 298 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 299 ierr = PetscViewerASCIIPrintf(viewer, "Coarse DM\n");CHKERRQ(ierr); 300 ierr = DMView(mesh->dmCoarse, viewer);CHKERRQ(ierr); 301 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 302 PetscFunctionReturn(0); 303 } 304 305 #undef __FUNCT__ 306 #define __FUNCT__ "DMView_Patch" 307 PetscErrorCode DMView_Patch(DM dm, PetscViewer viewer) 308 { 309 PetscBool iascii, isbinary; 310 PetscErrorCode ierr; 311 312 PetscFunctionBegin; 313 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 314 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 315 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 316 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERBINARY, &isbinary);CHKERRQ(ierr); 317 if (iascii) { 318 ierr = DMPatchView_Ascii(dm, viewer);CHKERRQ(ierr); 319 #if 0 320 } else if (isbinary) { 321 ierr = DMPatchView_Binary(dm, viewer);CHKERRQ(ierr); 322 #endif 323 } 324 PetscFunctionReturn(0); 325 } 326 327 #undef __FUNCT__ 328 #define __FUNCT__ "DMDestroy_Patch" 329 PetscErrorCode DMDestroy_Patch(DM dm) 330 { 331 DM_Patch *mesh = (DM_Patch*) dm->data; 332 PetscErrorCode ierr; 333 334 PetscFunctionBegin; 335 if (--mesh->refct > 0) PetscFunctionReturn(0); 336 ierr = DMDestroy(&mesh->dmCoarse);CHKERRQ(ierr); 337 /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 338 ierr = PetscFree(mesh);CHKERRQ(ierr); 339 PetscFunctionReturn(0); 340 } 341 342 #undef __FUNCT__ 343 #define __FUNCT__ "DMSetUp_Patch" 344 PetscErrorCode DMSetUp_Patch(DM dm) 345 { 346 DM_Patch *mesh = (DM_Patch*) dm->data; 347 PetscErrorCode ierr; 348 349 PetscFunctionBegin; 350 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 351 ierr = DMSetUp(mesh->dmCoarse);CHKERRQ(ierr); 352 PetscFunctionReturn(0); 353 } 354 355 #undef __FUNCT__ 356 #define __FUNCT__ "DMCreateGlobalVector_Patch" 357 PetscErrorCode DMCreateGlobalVector_Patch(DM dm, Vec *g) 358 { 359 DM_Patch *mesh = (DM_Patch*) dm->data; 360 PetscErrorCode ierr; 361 362 PetscFunctionBegin; 363 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 364 ierr = DMCreateGlobalVector(mesh->dmCoarse, g);CHKERRQ(ierr); 365 PetscFunctionReturn(0); 366 } 367 368 #undef __FUNCT__ 369 #define __FUNCT__ "DMCreateLocalVector_Patch" 370 PetscErrorCode DMCreateLocalVector_Patch(DM dm, Vec *l) 371 { 372 DM_Patch *mesh = (DM_Patch*) dm->data; 373 PetscErrorCode ierr; 374 375 PetscFunctionBegin; 376 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 377 ierr = DMCreateLocalVector(mesh->dmCoarse, l);CHKERRQ(ierr); 378 PetscFunctionReturn(0); 379 } 380 381 #undef __FUNCT__ 382 #define __FUNCT__ "DMCreateSubDM_Patch" 383 PetscErrorCode DMCreateSubDM_Patch(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm) 384 { 385 PetscFunctionBegin; 386 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 387 SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Tell me to code this"); 388 PetscFunctionReturn(0); 389 } 390 391 #undef __FUNCT__ 392 #define __FUNCT__ "DMPatchGetCoarse" 393 PetscErrorCode DMPatchGetCoarse(DM dm, DM *dmCoarse) 394 { 395 DM_Patch *mesh = (DM_Patch*) dm->data; 396 397 PetscFunctionBegin; 398 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 399 *dmCoarse = mesh->dmCoarse; 400 PetscFunctionReturn(0); 401 } 402 403 #undef __FUNCT__ 404 #define __FUNCT__ "DMPatchGetPatchSize" 405 PetscErrorCode DMPatchGetPatchSize(DM dm, MatStencil *patchSize) 406 { 407 DM_Patch *mesh = (DM_Patch*) dm->data; 408 409 PetscFunctionBegin; 410 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 411 PetscValidPointer(patchSize, 2); 412 *patchSize = mesh->patchSize; 413 PetscFunctionReturn(0); 414 } 415 416 #undef __FUNCT__ 417 #define __FUNCT__ "DMPatchSetPatchSize" 418 PetscErrorCode DMPatchSetPatchSize(DM dm, MatStencil patchSize) 419 { 420 DM_Patch *mesh = (DM_Patch*) dm->data; 421 422 PetscFunctionBegin; 423 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 424 mesh->patchSize = patchSize; 425 PetscFunctionReturn(0); 426 } 427 428 #undef __FUNCT__ 429 #define __FUNCT__ "DMPatchGetCommSize" 430 PetscErrorCode DMPatchGetCommSize(DM dm, MatStencil *commSize) 431 { 432 DM_Patch *mesh = (DM_Patch*) dm->data; 433 434 PetscFunctionBegin; 435 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 436 PetscValidPointer(commSize, 2); 437 *commSize = mesh->commSize; 438 PetscFunctionReturn(0); 439 } 440 441 #undef __FUNCT__ 442 #define __FUNCT__ "DMPatchSetCommSize" 443 PetscErrorCode DMPatchSetCommSize(DM dm, MatStencil commSize) 444 { 445 DM_Patch *mesh = (DM_Patch*) dm->data; 446 447 PetscFunctionBegin; 448 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 449 mesh->commSize = commSize; 450 PetscFunctionReturn(0); 451 } 452