Actual source code: stag1d.c

  1: /*
  2:    Functions specific to the 1-dimensional implementation of DMStag
  3: */
  4: #include <petsc/private/dmstagimpl.h>

  6: /*@
  7:   DMStagCreate1d - Create an object to manage data living on the elements and vertices of a parallelized regular 1D grid.

  9:   Collective

 11:   Input Parameters:
 12: + comm         - MPI communicator
 13: . bndx         - boundary type: `DM_BOUNDARY_NONE`, `DM_BOUNDARY_PERIODIC`, or `DM_BOUNDARY_GHOSTED`
 14: . M            - global number of elements
 15: . dof0         - number of degrees of freedom per vertex/0-cell
 16: . dof1         - number of degrees of freedom per element/1-cell
 17: . stencilType  - ghost/halo region type: `DMSTAG_STENCIL_BOX` or `DMSTAG_STENCIL_NONE`
 18: . stencilWidth - width, in elements, of halo/ghost region
 19: - lx           - array of local sizes, of length equal to the comm size, summing to `M` or `NULL`

 21:   Output Parameter:
 22: . dm - the new `DMSTAG` object

 24:   Options Database Keys:
 25: + -dm_view                                      - calls `DMViewFromOptions()` at the conclusion of `DMSetUp()`
 26: . -stag_grid_x nx                               - number of elements in the x direction
 27: . -stag_ghost_stencil_width                     - width of ghost region, in elements
 28: - -stag_boundary_type_x (none|ghosted|periodic) - `DMBoundaryType` value

 30:   Level: beginner

 32:   Notes:
 33:   You must call `DMSetUp()` after this call before using the `DM`.
 34:   If you wish to use the options database (see the keys above) to change values in the `DMSTAG`, you must call
 35:   `DMSetFromOptions()` after this function but before `DMSetUp()`.

 37: .seealso: [](ch_stag), `DMSTAG`, `DMStagCreate2d()`, `DMStagCreate3d()`, `DMDestroy()`, `DMView()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()`, `DMLocalToGlobalBegin()`, `DMDACreate1d()`
 38: @*/
 39: PetscErrorCode DMStagCreate1d(MPI_Comm comm, DMBoundaryType bndx, PetscInt M, PetscInt dof0, PetscInt dof1, DMStagStencilType stencilType, PetscInt stencilWidth, const PetscInt lx[], DM *dm)
 40: {
 41:   PetscMPIInt size;

 43:   PetscFunctionBegin;
 44:   PetscCallMPI(MPI_Comm_size(comm, &size));
 45:   PetscCall(DMCreate(comm, dm));
 46:   PetscCall(DMSetDimension(*dm, 1));
 47:   PetscCall(DMStagInitialize(bndx, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, M, 0, 0, size, 0, 0, dof0, dof1, 0, 0, stencilType, stencilWidth, lx, NULL, NULL, *dm));
 48:   PetscFunctionReturn(PETSC_SUCCESS);
 49: }

 51: PETSC_INTERN PetscErrorCode DMStagRestrictSimple_1d(DM dmf, Vec xf_local, DM dmc, Vec xc_local)
 52: {
 53:   PetscInt            Mf, Mc, factorx, dof[2];
 54:   PetscInt            xc, mc, nExtraxc, i, d;
 55:   PetscInt            ileftf, ielemf, ileftc, ielemc;
 56:   const PetscScalar **arrf;
 57:   PetscScalar       **arrc;

 59:   PetscFunctionBegin;
 60:   PetscCall(DMStagGetGlobalSizes(dmf, &Mf, NULL, NULL));
 61:   PetscCall(DMStagGetGlobalSizes(dmc, &Mc, NULL, NULL));
 62:   factorx = Mf / Mc;
 63:   PetscCall(DMStagGetDOF(dmc, &dof[0], &dof[1], NULL, NULL));

 65:   PetscCall(DMStagGetCorners(dmc, &xc, NULL, NULL, &mc, NULL, NULL, &nExtraxc, NULL, NULL));
 66:   PetscCall(VecZeroEntries(xc_local));
 67:   PetscCall(DMStagVecGetArray(dmf, xf_local, &arrf));
 68:   PetscCall(DMStagVecGetArray(dmc, xc_local, &arrc));
 69:   PetscCall(DMStagGetLocationSlot(dmf, DMSTAG_LEFT, 0, &ileftf));
 70:   PetscCall(DMStagGetLocationSlot(dmf, DMSTAG_ELEMENT, 0, &ielemf));
 71:   PetscCall(DMStagGetLocationSlot(dmc, DMSTAG_LEFT, 0, &ileftc));
 72:   PetscCall(DMStagGetLocationSlot(dmc, DMSTAG_ELEMENT, 0, &ielemc));

 74:   for (d = 0; d < dof[0]; ++d)
 75:     for (i = xc; i < xc + mc + nExtraxc; ++i) arrc[i][ileftc + d] = arrf[factorx * i][ileftf + d];

 77:   for (d = 0; d < dof[1]; ++d)
 78:     for (i = xc; i < xc + mc; ++i) {
 79:       if (factorx % 2 == 0) arrc[i][ielemc + d] = 0.5 * (arrf[factorx * i + factorx / 2 - 1][ielemf + d] + arrf[factorx * i + factorx / 2][ielemf + d]);
 80:       else arrc[i][ielemc + d] = arrf[factorx * i + factorx / 2][ielemf + d];
 81:     }

 83:   PetscCall(DMStagVecRestoreArray(dmf, xf_local, &arrf));
 84:   PetscCall(DMStagVecRestoreArray(dmc, xc_local, &arrc));
 85:   PetscFunctionReturn(PETSC_SUCCESS);
 86: }

 88: PETSC_INTERN PetscErrorCode DMStagSetUniformCoordinatesExplicit_1d(DM dm, PetscReal xmin, PetscReal xmax)
 89: {
 90:   DM_Stag      *stagCoord;
 91:   DM            dmCoord;
 92:   Vec           coordLocal;
 93:   PetscReal     h, min;
 94:   PetscScalar **arr;
 95:   PetscInt      start_ghost, n_ghost, s;
 96:   PetscInt      ileft, ielement;

 98:   PetscFunctionBegin;
 99:   PetscCall(DMGetCoordinateDM(dm, &dmCoord));
100:   stagCoord = (DM_Stag *)dmCoord->data;
101:   for (s = 0; s < 2; ++s) {
102:     PetscCheck(stagCoord->dof[s] == 0 || stagCoord->dof[s] == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Coordinate DM in 1 dimensions must have 0 or 1 dof on each stratum, but stratum %" PetscInt_FMT " has %" PetscInt_FMT " dof", s,
103:                stagCoord->dof[s]);
104:   }
105:   PetscCall(DMCreateLocalVector(dmCoord, &coordLocal));

107:   PetscCall(DMStagVecGetArray(dmCoord, coordLocal, &arr));
108:   if (stagCoord->dof[0]) PetscCall(DMStagGetLocationSlot(dmCoord, DMSTAG_LEFT, 0, &ileft));
109:   if (stagCoord->dof[1]) PetscCall(DMStagGetLocationSlot(dmCoord, DMSTAG_ELEMENT, 0, &ielement));
110:   PetscCall(DMStagGetGhostCorners(dmCoord, &start_ghost, NULL, NULL, &n_ghost, NULL, NULL));

112:   min = xmin;
113:   h   = (xmax - xmin) / stagCoord->N[0];

115:   for (PetscInt ind = start_ghost; ind < start_ghost + n_ghost; ++ind) {
116:     if (stagCoord->dof[0]) {
117:       const PetscReal off = 0.0;
118:       arr[ind][ileft]     = min + ((PetscReal)ind + off) * h;
119:     }
120:     if (stagCoord->dof[1]) {
121:       const PetscReal off = 0.5;
122:       arr[ind][ielement]  = min + ((PetscReal)ind + off) * h;
123:     }
124:   }
125:   PetscCall(DMStagVecRestoreArray(dmCoord, coordLocal, &arr));
126:   PetscCall(DMSetCoordinatesLocal(dm, coordLocal));
127:   PetscCall(VecDestroy(&coordLocal));
128:   PetscFunctionReturn(PETSC_SUCCESS);
129: }

131: /* Helper functions used in DMSetUp_Stag() */
132: static PetscErrorCode DMStagComputeLocationOffsets_1d(DM);

134: PETSC_INTERN PetscErrorCode DMSetUp_Stag_1d(DM dm)
135: {
136:   DM_Stag *const stag = (DM_Stag *)dm->data;
137:   PetscMPIInt    size, rank;
138:   MPI_Comm       comm;

140:   PetscFunctionBegin;
141:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
142:   PetscCallMPI(MPI_Comm_size(comm, &size));
143:   PetscCallMPI(MPI_Comm_rank(comm, &rank));

145:   /* Check Global size */
146:   PetscCheck(stag->N[0] >= 1, comm, PETSC_ERR_ARG_OUTOFRANGE, "Global grid size of %" PetscInt_FMT " < 1 specified", stag->N[0]);

148:   /* Local sizes */
149:   PetscCheck(stag->N[0] >= size, comm, PETSC_ERR_ARG_OUTOFRANGE, "More ranks (%d) than elements (%" PetscInt_FMT ") specified", size, stag->N[0]);
150:   if (!stag->l[0]) {
151:     /* Divide equally, giving an extra elements to higher ranks */
152:     PetscCall(PetscMalloc1(stag->nRanks[0], &stag->l[0]));
153:     for (PetscInt j = 0; j < stag->nRanks[0]; ++j) stag->l[0][j] = stag->N[0] / stag->nRanks[0] + (stag->N[0] % stag->nRanks[0] > j ? 1 : 0);
154:   }
155:   {
156:     PetscInt Nchk = 0;
157:     for (PetscInt j = 0; j < size; ++j) Nchk += stag->l[0][j];
158:     PetscCheck(Nchk == stag->N[0], comm, PETSC_ERR_ARG_OUTOFRANGE, "Sum of specified local sizes (%" PetscInt_FMT ") is not equal to global size (%" PetscInt_FMT ")", Nchk, stag->N[0]);
159:   }
160:   stag->n[0] = stag->l[0][rank];

162:   /* Rank (trivial in 1d) */
163:   stag->rank[0]      = rank;
164:   stag->firstRank[0] = (PetscBool)(rank == 0);
165:   stag->lastRank[0]  = (PetscBool)(rank == size - 1);

167:   /* Local (unghosted) numbers of entries */
168:   stag->entriesPerElement = stag->dof[0] + stag->dof[1];
169:   switch (stag->boundaryType[0]) {
170:   case DM_BOUNDARY_NONE:
171:   case DM_BOUNDARY_GHOSTED:
172:     stag->entries = stag->n[0] * stag->entriesPerElement + (stag->lastRank[0] ? stag->dof[0] : 0);
173:     break;
174:   case DM_BOUNDARY_PERIODIC:
175:     stag->entries = stag->n[0] * stag->entriesPerElement;
176:     break;
177:   default:
178:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported x boundary type %s", DMBoundaryTypes[stag->boundaryType[0]]);
179:   }

181:   /* Starting element */
182:   stag->start[0] = 0;
183:   for (PetscInt j = 0; j < stag->rank[0]; ++j) stag->start[0] += stag->l[0][j];

185:   /* Local/ghosted size and starting element */
186:   switch (stag->boundaryType[0]) {
187:   case DM_BOUNDARY_NONE:
188:     switch (stag->stencilType) {
189:     case DMSTAG_STENCIL_NONE: /* Only dummy cells on the right */
190:       stag->startGhost[0] = stag->start[0];
191:       stag->nGhost[0]     = stag->n[0] + (stag->lastRank[0] ? 1 : 0);
192:       break;
193:     case DMSTAG_STENCIL_STAR:
194:     case DMSTAG_STENCIL_BOX:
195:       stag->startGhost[0] = stag->firstRank[0] ? stag->start[0] : stag->start[0] - stag->stencilWidth;
196:       stag->nGhost[0]     = stag->n[0];
197:       stag->nGhost[0] += stag->firstRank[0] ? 0 : stag->stencilWidth;
198:       stag->nGhost[0] += stag->lastRank[0] ? 1 : stag->stencilWidth;
199:       break;
200:     default:
201:       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unrecognized ghost stencil type %d", stag->stencilType);
202:     }
203:     break;
204:   case DM_BOUNDARY_GHOSTED:
205:     switch (stag->stencilType) {
206:     case DMSTAG_STENCIL_NONE:
207:       stag->startGhost[0] = stag->start[0];
208:       stag->nGhost[0]     = stag->n[0] + (stag->lastRank[0] ? 1 : 0);
209:       break;
210:     case DMSTAG_STENCIL_STAR:
211:     case DMSTAG_STENCIL_BOX:
212:       stag->startGhost[0] = stag->start[0] - stag->stencilWidth; /* This value may be negative */
213:       stag->nGhost[0]     = stag->n[0] + 2 * stag->stencilWidth + (stag->lastRank[0] && stag->stencilWidth == 0 ? 1 : 0);
214:       break;
215:     default:
216:       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unrecognized ghost stencil type %d", stag->stencilType);
217:     }
218:     break;
219:   case DM_BOUNDARY_PERIODIC:
220:     switch (stag->stencilType) {
221:     case DMSTAG_STENCIL_NONE:
222:       stag->startGhost[0] = stag->start[0];
223:       stag->nGhost[0]     = stag->n[0];
224:       break;
225:     case DMSTAG_STENCIL_STAR:
226:     case DMSTAG_STENCIL_BOX:
227:       stag->startGhost[0] = stag->start[0] - stag->stencilWidth; /* This value may be negative */
228:       stag->nGhost[0]     = stag->n[0] + 2 * stag->stencilWidth;
229:       break;
230:     default:
231:       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unrecognized ghost stencil type %d", stag->stencilType);
232:     }
233:     break;
234:   default:
235:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported x boundary type %s", DMBoundaryTypes[stag->boundaryType[0]]);
236:   }

238:   /* Total size of ghosted/local representation */
239:   stag->entriesGhost = stag->nGhost[0] * stag->entriesPerElement;

241:   /* Define neighbors */
242:   PetscCall(PetscMalloc1(3, &stag->neighbors));
243:   if (stag->firstRank[0]) {
244:     switch (stag->boundaryType[0]) {
245:     case DM_BOUNDARY_GHOSTED:
246:     case DM_BOUNDARY_NONE:
247:       stag->neighbors[0] = -1;
248:       break;
249:     case DM_BOUNDARY_PERIODIC:
250:       stag->neighbors[0] = stag->nRanks[0] - 1;
251:       break;
252:     default:
253:       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported x boundary type %s", DMBoundaryTypes[stag->boundaryType[0]]);
254:     }
255:   } else {
256:     stag->neighbors[0] = stag->rank[0] - 1;
257:   }
258:   stag->neighbors[1] = stag->rank[0];
259:   if (stag->lastRank[0]) {
260:     switch (stag->boundaryType[0]) {
261:     case DM_BOUNDARY_GHOSTED:
262:     case DM_BOUNDARY_NONE:
263:       stag->neighbors[2] = -1;
264:       break;
265:     case DM_BOUNDARY_PERIODIC:
266:       stag->neighbors[2] = 0;
267:       break;
268:     default:
269:       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported x boundary type %s", DMBoundaryTypes[stag->boundaryType[0]]);
270:     }
271:   } else {
272:     stag->neighbors[2] = stag->rank[0] + 1;
273:   }

275:   PetscCheck(stag->n[0] >= stag->stencilWidth, PETSC_COMM_SELF, PETSC_ERR_SUP, "DMStag 1d setup does not support local sizes (%" PetscInt_FMT ") smaller than the elementwise stencil width (%" PetscInt_FMT ")", stag->n[0], stag->stencilWidth);

277:   /* Create global->local VecScatter and ISLocalToGlobalMapping */
278:   {
279:     PetscInt *idxLocal, *idxGlobal, *idxGlobalAll;
280:     PetscInt  i, iLocal, d, entriesToTransferTotal, ghostOffsetStart, ghostOffsetEnd, nNonDummyGhost;
281:     IS        isLocal, isGlobal;

283:     /* The offset on the right (may not be equal to the stencil width, as we
284:        always have at least one ghost element, to account for the boundary
285:        point, and may with ghosted boundaries), and the number of non-dummy ghost elements */
286:     ghostOffsetStart = stag->start[0] - stag->startGhost[0];
287:     ghostOffsetEnd   = stag->startGhost[0] + stag->nGhost[0] - (stag->start[0] + stag->n[0]);
288:     nNonDummyGhost   = stag->nGhost[0] - (stag->lastRank[0] ? ghostOffsetEnd : 0) - (stag->firstRank[0] ? ghostOffsetStart : 0);

290:     /* Compute the number of non-dummy entries in the local representation
291:        This is equal to the number of non-dummy elements in the local (ghosted) representation,
292:        plus some extra entries on the right boundary on the last rank*/
293:     switch (stag->boundaryType[0]) {
294:     case DM_BOUNDARY_GHOSTED:
295:     case DM_BOUNDARY_NONE:
296:       entriesToTransferTotal = nNonDummyGhost * stag->entriesPerElement + (stag->lastRank[0] ? stag->dof[0] : 0);
297:       break;
298:     case DM_BOUNDARY_PERIODIC:
299:       entriesToTransferTotal = stag->entriesGhost; /* No dummy points */
300:       break;
301:     default:
302:       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported x boundary type %s", DMBoundaryTypes[stag->boundaryType[0]]);
303:     }

305:     PetscCall(PetscMalloc1(entriesToTransferTotal, &idxLocal));
306:     PetscCall(PetscMalloc1(entriesToTransferTotal, &idxGlobal));
307:     PetscCall(PetscMalloc1(stag->entriesGhost, &idxGlobalAll));
308:     if (stag->boundaryType[0] == DM_BOUNDARY_NONE) {
309:       PetscInt count = 0, countAll = 0;
310:       /* Left ghost points and native points */
311:       for (i = stag->startGhost[0], iLocal = 0; iLocal < nNonDummyGhost; ++i, ++iLocal) {
312:         for (d = 0; d < stag->entriesPerElement; ++d, ++count, ++countAll) {
313:           idxLocal[count]        = iLocal * stag->entriesPerElement + d;
314:           idxGlobal[count]       = i * stag->entriesPerElement + d;
315:           idxGlobalAll[countAll] = i * stag->entriesPerElement + d;
316:         }
317:       }
318:       /* Ghost points on the right
319:          Special case for last (partial dummy) element on the last rank */
320:       if (stag->lastRank[0]) {
321:         i      = stag->N[0];
322:         iLocal = (stag->nGhost[0] - ghostOffsetEnd);
323:         /* Only vertex (0-cell) dofs in global representation */
324:         for (d = 0; d < stag->dof[0]; ++d, ++count, ++countAll) {
325:           idxGlobal[count]       = i * stag->entriesPerElement + d;
326:           idxLocal[count]        = iLocal * stag->entriesPerElement + d;
327:           idxGlobalAll[countAll] = i * stag->entriesPerElement + d;
328:         }
329:         for (d = stag->dof[0]; d < stag->entriesPerElement; ++d, ++countAll) { /* Additional dummy entries */
330:           idxGlobalAll[countAll] = -1;
331:         }
332:       }
333:     } else if (stag->boundaryType[0] == DM_BOUNDARY_PERIODIC) {
334:       PetscInt       count = 0, iLocal = 0; /* No dummy points, so idxGlobal and idxGlobalAll are identical */
335:       const PetscInt iMin = stag->firstRank[0] ? stag->start[0] : stag->startGhost[0];
336:       const PetscInt iMax = stag->lastRank[0] ? stag->startGhost[0] + stag->nGhost[0] - stag->stencilWidth : stag->startGhost[0] + stag->nGhost[0];
337:       /* Ghost points on the left */
338:       if (stag->firstRank[0]) {
339:         for (i = stag->N[0] - stag->stencilWidth; iLocal < stag->stencilWidth; ++i, ++iLocal) {
340:           for (d = 0; d < stag->entriesPerElement; ++d, ++count) {
341:             idxGlobal[count]    = i * stag->entriesPerElement + d;
342:             idxLocal[count]     = iLocal * stag->entriesPerElement + d;
343:             idxGlobalAll[count] = idxGlobal[count];
344:           }
345:         }
346:       }
347:       /* Native points */
348:       for (i = iMin; i < iMax; ++i, ++iLocal) {
349:         for (d = 0; d < stag->entriesPerElement; ++d, ++count) {
350:           idxGlobal[count]    = i * stag->entriesPerElement + d;
351:           idxLocal[count]     = iLocal * stag->entriesPerElement + d;
352:           idxGlobalAll[count] = idxGlobal[count];
353:         }
354:       }
355:       /* Ghost points on the right */
356:       if (stag->lastRank[0]) {
357:         for (i = 0; iLocal < stag->nGhost[0]; ++i, ++iLocal) {
358:           for (d = 0; d < stag->entriesPerElement; ++d, ++count) {
359:             idxGlobal[count]    = i * stag->entriesPerElement + d;
360:             idxLocal[count]     = iLocal * stag->entriesPerElement + d;
361:             idxGlobalAll[count] = idxGlobal[count];
362:           }
363:         }
364:       }
365:     } else if (stag->boundaryType[0] == DM_BOUNDARY_GHOSTED) {
366:       PetscInt count = 0, countAll = 0;
367:       /* Dummy elements on the left, on the first rank */
368:       if (stag->firstRank[0]) {
369:         for (iLocal = 0; iLocal < ghostOffsetStart; ++iLocal) {
370:           /* Complete elements full of dummy entries */
371:           for (d = 0; d < stag->entriesPerElement; ++d, ++countAll) idxGlobalAll[countAll] = -1;
372:         }
373:         i = 0; /* nonDummy entries start with global entry 0 */
374:       } else {
375:         /* nonDummy entries start as usual */
376:         i      = stag->startGhost[0];
377:         iLocal = 0;
378:       }

380:       /* non-Dummy entries */
381:       {
382:         PetscInt iLocalNonDummyMax = stag->firstRank[0] ? nNonDummyGhost + ghostOffsetStart : nNonDummyGhost;
383:         for (; iLocal < iLocalNonDummyMax; ++i, ++iLocal) {
384:           for (d = 0; d < stag->entriesPerElement; ++d, ++count, ++countAll) {
385:             idxLocal[count]        = iLocal * stag->entriesPerElement + d;
386:             idxGlobal[count]       = i * stag->entriesPerElement + d;
387:             idxGlobalAll[countAll] = i * stag->entriesPerElement + d;
388:           }
389:         }
390:       }

392:       /* (partial) dummy elements on the right, on the last rank */
393:       if (stag->lastRank[0]) {
394:         /* First one is partial dummy */
395:         i      = stag->N[0];
396:         iLocal = (stag->nGhost[0] - ghostOffsetEnd);
397:         for (d = 0; d < stag->dof[0]; ++d, ++count, ++countAll) { /* Only vertex (0-cell) dofs in global representation */
398:           idxLocal[count]        = iLocal * stag->entriesPerElement + d;
399:           idxGlobal[count]       = i * stag->entriesPerElement + d;
400:           idxGlobalAll[countAll] = i * stag->entriesPerElement + d;
401:         }
402:         for (d = stag->dof[0]; d < stag->entriesPerElement; ++d, ++countAll) { /* Additional dummy entries */
403:           idxGlobalAll[countAll] = -1;
404:         }
405:         for (iLocal = stag->nGhost[0] - ghostOffsetEnd + 1; iLocal < stag->nGhost[0]; ++iLocal) {
406:           /* Additional dummy elements */
407:           for (d = 0; d < stag->entriesPerElement; ++d, ++countAll) idxGlobalAll[countAll] = -1;
408:         }
409:       }
410:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported x boundary type %s", DMBoundaryTypes[stag->boundaryType[0]]);

412:     /* Create Local IS (transferring pointer ownership) */
413:     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), entriesToTransferTotal, idxLocal, PETSC_OWN_POINTER, &isLocal));

415:     /* Create Global IS (transferring pointer ownership) */
416:     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), entriesToTransferTotal, idxGlobal, PETSC_OWN_POINTER, &isGlobal));

418:     /* Create stag->gtol, which doesn't include dummy entries */
419:     {
420:       Vec local, global;
421:       PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)dm), 1, stag->entries, PETSC_DECIDE, NULL, &global));
422:       PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, PetscMax(stag->entriesPerElement, 1), stag->entriesGhost, NULL, &local));
423:       PetscCall(VecScatterCreate(global, isGlobal, local, isLocal, &stag->gtol));
424:       PetscCall(VecDestroy(&global));
425:       PetscCall(VecDestroy(&local));
426:     }

428:     /* In special cases, create a dedicated injective local-to-global map */
429:     if (stag->boundaryType[0] == DM_BOUNDARY_PERIODIC && stag->nRanks[0] == 1) PetscCall(DMStagPopulateLocalToGlobalInjective(dm));

431:     /* Destroy ISs */
432:     PetscCall(ISDestroy(&isLocal));
433:     PetscCall(ISDestroy(&isGlobal));

435:     /* Create local-to-global map (transferring pointer ownership) */
436:     PetscCall(ISLocalToGlobalMappingCreate(comm, 1, stag->entriesGhost, idxGlobalAll, PETSC_OWN_POINTER, &dm->ltogmap));
437:   }

439:   /* Precompute location offsets */
440:   PetscCall(DMStagComputeLocationOffsets_1d(dm));

442:   /* View from Options */
443:   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
444:   PetscFunctionReturn(PETSC_SUCCESS);
445: }

447: static PetscErrorCode DMStagComputeLocationOffsets_1d(DM dm)
448: {
449:   DM_Stag *const stag = (DM_Stag *)dm->data;
450:   const PetscInt epe  = stag->entriesPerElement;

452:   PetscFunctionBegin;
453:   PetscCall(PetscMalloc1(DMSTAG_NUMBER_LOCATIONS, &stag->locationOffsets));
454:   stag->locationOffsets[DMSTAG_LEFT]    = 0;
455:   stag->locationOffsets[DMSTAG_ELEMENT] = stag->locationOffsets[DMSTAG_LEFT] + stag->dof[0];
456:   stag->locationOffsets[DMSTAG_RIGHT]   = stag->locationOffsets[DMSTAG_LEFT] + epe;
457:   PetscFunctionReturn(PETSC_SUCCESS);
458: }

460: PETSC_INTERN PetscErrorCode DMStagPopulateLocalToGlobalInjective_1d(DM dm)
461: {
462:   DM_Stag *const stag = (DM_Stag *)dm->data;
463:   PetscInt      *idxLocal, *idxGlobal;
464:   PetscInt       i, iLocal, d, count;
465:   IS             isLocal, isGlobal;

467:   PetscFunctionBegin;
468:   PetscCall(PetscMalloc1(stag->entries, &idxLocal));
469:   PetscCall(PetscMalloc1(stag->entries, &idxGlobal));
470:   count  = 0;
471:   iLocal = stag->start[0] - stag->startGhost[0];
472:   for (i = stag->start[0]; i < stag->start[0] + stag->n[0]; ++i, ++iLocal) {
473:     for (d = 0; d < stag->entriesPerElement; ++d, ++count) {
474:       idxGlobal[count] = i * stag->entriesPerElement + d;
475:       idxLocal[count]  = iLocal * stag->entriesPerElement + d;
476:     }
477:   }
478:   if (stag->lastRank[0] && stag->boundaryType[0] != DM_BOUNDARY_PERIODIC) {
479:     i      = stag->start[0] + stag->n[0];
480:     iLocal = stag->start[0] - stag->startGhost[0] + stag->n[0];
481:     for (d = 0; d < stag->dof[0]; ++d, ++count) {
482:       idxGlobal[count] = i * stag->entriesPerElement + d;
483:       idxLocal[count]  = iLocal * stag->entriesPerElement + d;
484:     }
485:   }
486:   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), stag->entries, idxLocal, PETSC_OWN_POINTER, &isLocal));
487:   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), stag->entries, idxGlobal, PETSC_OWN_POINTER, &isGlobal));
488:   {
489:     Vec local, global;
490:     PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)dm), 1, stag->entries, PETSC_DECIDE, NULL, &global));
491:     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, PetscMax(stag->entriesPerElement, 1), stag->entriesGhost, NULL, &local));
492:     PetscCall(VecScatterCreate(local, isLocal, global, isGlobal, &stag->ltog_injective));
493:     PetscCall(VecDestroy(&global));
494:     PetscCall(VecDestroy(&local));
495:   }
496:   PetscCall(ISDestroy(&isLocal));
497:   PetscCall(ISDestroy(&isGlobal));
498:   PetscFunctionReturn(PETSC_SUCCESS);
499: }

501: PETSC_INTERN PetscErrorCode DMStagPopulateLocalToLocal1d_Internal(DM dm)
502: {
503:   DM_Stag *const stag = (DM_Stag *)dm->data;
504:   PetscInt      *idxRemap;
505:   PetscInt       leftGhostEntries;

507:   PetscFunctionBegin;
508:   PetscCall(VecScatterCopy(stag->gtol, &stag->ltol));
509:   PetscCall(PetscMalloc1(stag->entries, &idxRemap));

511:   leftGhostEntries = (stag->start[0] - stag->startGhost[0]) * stag->entriesPerElement;
512:   for (PetscInt i = 0; i < stag->entries; ++i) idxRemap[i] = leftGhostEntries + i;

514:   PetscCall(VecScatterRemap(stag->ltol, idxRemap, NULL));
515:   PetscCall(PetscFree(idxRemap));
516:   PetscFunctionReturn(PETSC_SUCCESS);
517: }

519: PETSC_INTERN PetscErrorCode DMCreateMatrix_Stag_1D_AIJ_Assemble(DM dm, Mat A)
520: {
521:   DMStagStencilType stencil_type;
522:   PetscInt          dof[2], start, n, n_extra, stencil_width, N, epe;
523:   DMBoundaryType    boundary_type_x;

525:   PetscFunctionBegin;
526:   PetscCall(DMStagGetDOF(dm, &dof[0], &dof[1], NULL, NULL));
527:   PetscCall(DMStagGetStencilType(dm, &stencil_type));
528:   PetscCall(DMStagGetStencilWidth(dm, &stencil_width));
529:   PetscCall(DMStagGetCorners(dm, &start, NULL, NULL, &n, NULL, NULL, &n_extra, NULL, NULL));
530:   PetscCall(DMStagGetGlobalSizes(dm, &N, NULL, NULL));
531:   PetscCall(DMStagGetEntriesPerElement(dm, &epe));
532:   PetscCall(DMStagGetBoundaryTypes(dm, &boundary_type_x, NULL, NULL));
533:   if (stencil_type == DMSTAG_STENCIL_NONE) {
534:     /* Couple all DOF at each location to each other */
535:     DMStagStencil *row_vertex, *row_element;

537:     PetscCall(PetscMalloc1(dof[0], &row_vertex));
538:     for (PetscInt c = 0; c < dof[0]; ++c) {
539:       row_vertex[c].loc = DMSTAG_LEFT;
540:       row_vertex[c].c   = c;
541:     }

543:     PetscCall(PetscMalloc1(dof[1], &row_element));
544:     for (PetscInt c = 0; c < dof[1]; ++c) {
545:       row_element[c].loc = DMSTAG_ELEMENT;
546:       row_element[c].c   = c;
547:     }

549:     for (PetscInt e = start; e < start + n + n_extra; ++e) {
550:       {
551:         for (PetscInt c = 0; c < dof[0]; ++c) row_vertex[c].i = e;
552:         PetscCall(DMStagMatSetValuesStencil(dm, A, dof[0], row_vertex, dof[0], row_vertex, NULL, INSERT_VALUES));
553:       }
554:       if (e < N) {
555:         for (PetscInt c = 0; c < dof[1]; ++c) row_element[c].i = e;
556:         PetscCall(DMStagMatSetValuesStencil(dm, A, dof[1], row_element, dof[1], row_element, NULL, INSERT_VALUES));
557:       }
558:     }
559:     PetscCall(PetscFree(row_vertex));
560:     PetscCall(PetscFree(row_element));
561:   } else if (stencil_type == DMSTAG_STENCIL_STAR || stencil_type == DMSTAG_STENCIL_BOX) {
562:     DMStagStencil *col, *row;

564:     PetscCall(PetscMalloc1(epe, &row));
565:     {
566:       PetscInt nrows = 0;
567:       for (PetscInt c = 0; c < dof[0]; ++c) {
568:         row[nrows].c   = c;
569:         row[nrows].loc = DMSTAG_LEFT;
570:         ++nrows;
571:       }
572:       for (PetscInt c = 0; c < dof[1]; ++c) {
573:         row[nrows].c   = c;
574:         row[nrows].loc = DMSTAG_ELEMENT;
575:         ++nrows;
576:       }
577:     }
578:     PetscCall(PetscMalloc1(epe, &col));
579:     {
580:       PetscInt ncols = 0;
581:       for (PetscInt c = 0; c < dof[0]; ++c) {
582:         col[ncols].c   = c;
583:         col[ncols].loc = DMSTAG_LEFT;
584:         ++ncols;
585:       }
586:       for (PetscInt c = 0; c < dof[1]; ++c) {
587:         col[ncols].c   = c;
588:         col[ncols].loc = DMSTAG_ELEMENT;
589:         ++ncols;
590:       }
591:     }
592:     for (PetscInt e = start; e < start + n + n_extra; ++e) {
593:       for (PetscInt i = 0; i < epe; ++i) row[i].i = e;
594:       for (PetscInt offset = -stencil_width; offset <= stencil_width; ++offset) {
595:         const PetscInt e_offset = e + offset;

597:         /* Only set values corresponding to elements which can have non-dummy entries,
598:            meaning those that map to unknowns in the global representation. In the periodic
599:            case, this is the entire stencil, but in all other cases, only includes a single
600:            "extra" element which is partially outside the physical domain (those points in the
601:            global representation */
602:         if (boundary_type_x == DM_BOUNDARY_PERIODIC || (e_offset < N + 1 && e_offset >= 0)) {
603:           for (PetscInt i = 0; i < epe; ++i) col[i].i = e_offset;
604:           PetscCall(DMStagMatSetValuesStencil(dm, A, epe, row, epe, col, NULL, INSERT_VALUES));
605:         }
606:       }
607:     }
608:     PetscCall(PetscFree(row));
609:     PetscCall(PetscFree(col));
610:   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported stencil type %s", DMStagStencilTypes[stencil_type]);
611:   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
612:   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
613:   PetscFunctionReturn(PETSC_SUCCESS);
614: }