Actual source code: plexhdf5.c

  1: #include <petsc/private/dmpleximpl.h>
  2: #include <petsc/private/isimpl.h>
  3: #include <petsc/private/vecimpl.h>
  4: #include <petsc/private/viewerhdf5impl.h>
  5: #include <petsclayouthdf5.h>

  7: /* Logging support */
  8: PetscLogEvent DMPLEX_DistributionView, DMPLEX_DistributionLoad;

 10: #if defined(PETSC_HAVE_HDF5)
 11: static PetscErrorCode PetscViewerParseVersion_Private(PetscViewer, const char[], DMPlexStorageVersion *);
 12: static PetscErrorCode PetscViewerCheckVersion_Private(PetscViewer, DMPlexStorageVersion);
 13: static PetscErrorCode PetscViewerAttachVersion_Private(PetscViewer, const char[], DMPlexStorageVersion);
 14: static PetscErrorCode PetscViewerGetAttachedVersion_Private(PetscViewer, const char[], DMPlexStorageVersion *);

 16: PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer);

 18: static PetscErrorCode PetscViewerPrintVersion_Private(PetscViewer viewer, DMPlexStorageVersion version, char str[], size_t len)
 19: {
 20:   PetscFunctionBegin;
 21:   PetscCall(PetscViewerCheckVersion_Private(viewer, version));
 22:   PetscCall(PetscSNPrintf(str, len, "%d.%d.%d", version->major, version->minor, version->subminor));
 23:   PetscFunctionReturn(PETSC_SUCCESS);
 24: }

 26: static PetscErrorCode PetscViewerParseVersion_Private(PetscViewer viewer, const char str[], DMPlexStorageVersion *version)
 27: {
 28:   PetscToken           t;
 29:   char                *ts;
 30:   PetscInt             i;
 31:   PetscInt             ti[3];
 32:   DMPlexStorageVersion v;

 34:   PetscFunctionBegin;
 35:   PetscCall(PetscTokenCreate(str, '.', &t));
 36:   for (i = 0; i < 3; i++) {
 37:     PetscCall(PetscTokenFind(t, &ts));
 38:     PetscCheck(ts, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Malformed version string %s", str);
 39:     PetscCall(PetscOptionsStringToInt(ts, &ti[i]));
 40:   }
 41:   PetscCall(PetscTokenFind(t, &ts));
 42:   PetscCheck(!ts, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Malformed version string %s", str);
 43:   PetscCall(PetscTokenDestroy(&t));
 44:   PetscCall(PetscNew(&v));
 45:   v->major    = (int)ti[0];
 46:   v->minor    = (int)ti[1];
 47:   v->subminor = (int)ti[2];
 48:   PetscCall(PetscViewerCheckVersion_Private(viewer, v));
 49:   *version = v;
 50:   PetscFunctionReturn(PETSC_SUCCESS);
 51: }

 53: static PetscErrorCode PetscViewerAttachVersion_Private(PetscViewer viewer, const char key[], DMPlexStorageVersion v)
 54: {
 55:   PetscFunctionBegin;
 56:   PetscCall(PetscObjectContainerCompose((PetscObject)viewer, key, v, PetscContainerUserDestroyDefault));
 57:   PetscFunctionReturn(PETSC_SUCCESS);
 58: }

 60: static PetscErrorCode PetscViewerGetAttachedVersion_Private(PetscViewer viewer, const char key[], DMPlexStorageVersion *v)
 61: {
 62:   PetscContainer cont;

 64:   PetscFunctionBegin;
 65:   PetscCall(PetscObjectQuery((PetscObject)viewer, key, (PetscObject *)&cont));
 66:   *v = NULL;
 67:   if (cont) PetscCall(PetscContainerGetPointer(cont, (void **)v));
 68:   PetscFunctionReturn(PETSC_SUCCESS);
 69: }

 71: /*
 72:   Version log:
 73:   1.0.0 legacy version (default if no "dmplex_storage_version" attribute found in file)
 74:   2.0.0 introduce versioning and multiple topologies
 75:   2.1.0 introduce distributions
 76: */
 77: static PetscErrorCode PetscViewerCheckVersion_Private(PetscViewer viewer, DMPlexStorageVersion version)
 78: {
 79:   PetscBool valid = PETSC_FALSE;

 81:   PetscFunctionBegin;
 82:   switch (version->major) {
 83:   case 1:
 84:     switch (version->minor) {
 85:     case 0:
 86:       switch (version->subminor) {
 87:       case 0:
 88:         valid = PETSC_TRUE;
 89:         break;
 90:       }
 91:       break;
 92:     }
 93:     break;
 94:   case 2:
 95:     switch (version->minor) {
 96:     case 0:
 97:       switch (version->subminor) {
 98:       case 0:
 99:         valid = PETSC_TRUE;
100:         break;
101:       }
102:       break;
103:     case 1:
104:       switch (version->subminor) {
105:       case 0:
106:         valid = PETSC_TRUE;
107:         break;
108:       }
109:       break;
110:     }
111:     break;
112:   case 3:
113:     switch (version->minor) {
114:     case 0:
115:       switch (version->subminor) {
116:       case 0:
117:         valid = PETSC_TRUE;
118:         break;
119:       }
120:       break;
121:     }
122:     break;
123:   }
124:   PetscCheck(valid, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "DMPlexStorageVersion %d.%d.%d not supported", version->major, version->minor, version->subminor);
125:   PetscFunctionReturn(PETSC_SUCCESS);
126: }

128: static inline PetscBool DMPlexStorageVersionGE(DMPlexStorageVersion version, int major, int minor, int subminor)
129: {
130:   return (PetscBool)((version->major == major && version->minor == minor && version->subminor >= subminor) || (version->major == major && version->minor > minor) || (version->major > major));
131: }

133: /*@C
134:   PetscViewerHDF5SetDMPlexStorageVersionWriting - Set the storage version for writing

136:   Logically collective

138:   Input Parameters:
139: + viewer  - The `PetscViewer`
140: - version - The storage format version

142:   Level: advanced

144:   Note:
145:   The version has major, minor, and subminor integers. Parallel operations are only available for version 3.0.0.

147: .seealso: [](ch_dmbase), `DM`, `PetscViewerHDF5GetDMPlexStorageVersionWriting()`, `PetscViewerHDF5GetDMPlexStorageVersionReading()`, `PetscViewerHDF5SetDMPlexStorageVersionReading()`
148: @*/
149: PetscErrorCode PetscViewerHDF5SetDMPlexStorageVersionWriting(PetscViewer viewer, DMPlexStorageVersion version)
150: {
151:   const char           ATTR_NAME[] = "dmplex_storage_version";
152:   DMPlexStorageVersion viewerVersion;
153:   PetscBool            fileHasVersion;
154:   char                 fileVersion[16], versionStr[16], viewerVersionStr[16];

156:   PetscFunctionBegin;
158:   PetscAssertPointer(version, 2);
159:   PetscCall(PetscViewerPrintVersion_Private(viewer, version, versionStr, 16));
160:   PetscCall(PetscViewerGetAttachedVersion_Private(viewer, DMPLEX_STORAGE_VERSION_WRITING_KEY, &viewerVersion));
161:   if (viewerVersion) {
162:     PetscBool flg;

164:     PetscCall(PetscViewerPrintVersion_Private(viewer, viewerVersion, viewerVersionStr, 16));
165:     PetscCall(PetscStrcmp(versionStr, viewerVersionStr, &flg));
166:     PetscCheck(flg, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "User requested DMPlex storage version %s but viewer already has version %s - cannot mix versions", versionStr, viewerVersionStr);
167:   }

169:   PetscCall(PetscViewerHDF5HasAttribute(viewer, NULL, ATTR_NAME, &fileHasVersion));
170:   if (fileHasVersion) {
171:     PetscBool flg;
172:     char     *tmp;

174:     PetscCall(PetscViewerHDF5ReadAttribute(viewer, "/", ATTR_NAME, PETSC_STRING, NULL, &tmp));
175:     PetscCall(PetscStrncpy(fileVersion, tmp, sizeof(fileVersion)));
176:     PetscCall(PetscFree(tmp));
177:     PetscCall(PetscStrcmp(fileVersion, versionStr, &flg));
178:     PetscCheck(flg, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "User requested DMPlex storage version %s but file already has version %s - cannot mix versions", versionStr, fileVersion);
179:   } else {
180:     PetscCall(PetscViewerHDF5WriteAttribute(viewer, "/", ATTR_NAME, PETSC_STRING, versionStr));
181:   }
182:   PetscCall(PetscNew(&viewerVersion));
183:   viewerVersion->major    = version->major;
184:   viewerVersion->minor    = version->minor;
185:   viewerVersion->subminor = version->subminor;
186:   PetscCall(PetscViewerAttachVersion_Private(viewer, DMPLEX_STORAGE_VERSION_WRITING_KEY, viewerVersion));
187:   PetscFunctionReturn(PETSC_SUCCESS);
188: }

190: /*@C
191:   PetscViewerHDF5GetDMPlexStorageVersionWriting - Get the storage version for writing

193:   Logically collective

195:   Input Parameter:
196: . viewer - The `PetscViewer`

198:   Output Parameter:
199: . version - The storage format version

201:   Options Database Keys:
202: . -dm_plex_view_hdf5_storage_version <num> - Overrides the storage format version

204:   Level: advanced

206:   Note:
207:   The version has major, minor, and subminor integers. Parallel operations are only available for version 3.0.0.

209: .seealso: [](ch_dmbase), `DM`, `PetscViewerHDF5SetDMPlexStorageVersionWriting()`, `PetscViewerHDF5GetDMPlexStorageVersionReading()`, `PetscViewerHDF5SetDMPlexStorageVersionReading()`
210: @*/
211: PetscErrorCode PetscViewerHDF5GetDMPlexStorageVersionWriting(PetscViewer viewer, DMPlexStorageVersion *version)
212: {
213:   const char ATTR_NAME[] = "dmplex_storage_version";
214:   PetscBool  fileHasVersion;
215:   char       optVersion[16], fileVersion[16];

217:   PetscFunctionBegin;
219:   PetscAssertPointer(version, 2);
220:   PetscCall(PetscViewerGetAttachedVersion_Private(viewer, DMPLEX_STORAGE_VERSION_WRITING_KEY, version));
221:   if (*version) PetscFunctionReturn(PETSC_SUCCESS);

223:   PetscCall(PetscStrncpy(fileVersion, DMPLEX_STORAGE_VERSION_STABLE, sizeof(fileVersion)));
224:   PetscCall(PetscViewerHDF5HasAttribute(viewer, NULL, ATTR_NAME, &fileHasVersion));
225:   if (fileHasVersion) {
226:     char *tmp;

228:     PetscCall(PetscViewerHDF5ReadAttribute(viewer, "/", ATTR_NAME, PETSC_STRING, NULL, &tmp));
229:     PetscCall(PetscStrncpy(fileVersion, tmp, sizeof(fileVersion)));
230:     PetscCall(PetscFree(tmp));
231:   }
232:   PetscCall(PetscStrncpy(optVersion, fileVersion, sizeof(optVersion)));
233:   PetscObjectOptionsBegin((PetscObject)viewer);
234:   PetscCall(PetscOptionsString("-dm_plex_view_hdf5_storage_version", "DMPlex HDF5 viewer storage version", NULL, optVersion, optVersion, sizeof(optVersion), NULL));
235:   PetscOptionsEnd();
236:   if (!fileHasVersion) {
237:     PetscCall(PetscViewerHDF5WriteAttribute(viewer, "/", ATTR_NAME, PETSC_STRING, optVersion));
238:   } else {
239:     PetscBool flg;

241:     PetscCall(PetscStrcmp(fileVersion, optVersion, &flg));
242:     PetscCheck(flg, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "User requested DMPlex storage version %s but file already has version %s - cannot mix versions", optVersion, fileVersion);
243:   }
244:   PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "petsc_version_git", PETSC_STRING, PETSC_VERSION_GIT));
245:   PetscCall(PetscViewerParseVersion_Private(viewer, optVersion, version));
246:   PetscCall(PetscViewerAttachVersion_Private(viewer, DMPLEX_STORAGE_VERSION_WRITING_KEY, *version));
247:   PetscFunctionReturn(PETSC_SUCCESS);
248: }

250: /*@C
251:   PetscViewerHDF5SetDMPlexStorageVersionReading - Set the storage version for reading

253:   Logically collective

255:   Input Parameters:
256: + viewer  - The `PetscViewer`
257: - version - The storage format version

259:   Level: advanced

261:   Note:
262:   The version has major, minor, and subminor integers. Parallel operations are only available for version 3.0.0.

264: .seealso: [](ch_dmbase), `DM`, `PetscViewerHDF5GetDMPlexStorageVersionReading()`, `PetscViewerHDF5GetDMPlexStorageVersionWriting()`, `PetscViewerHDF5SetDMPlexStorageVersionWriting()`
265: @*/
266: PetscErrorCode PetscViewerHDF5SetDMPlexStorageVersionReading(PetscViewer viewer, DMPlexStorageVersion version)
267: {
268:   const char           ATTR_NAME[] = "dmplex_storage_version";
269:   DMPlexStorageVersion viewerVersion;
270:   PetscBool            fileHasVersion;
271:   char                 versionStr[16], viewerVersionStr[16];

273:   PetscFunctionBegin;
275:   PetscAssertPointer(version, 2);
276:   PetscCall(PetscViewerPrintVersion_Private(viewer, version, versionStr, 16));
277:   PetscCall(PetscViewerGetAttachedVersion_Private(viewer, DMPLEX_STORAGE_VERSION_READING_KEY, &viewerVersion));
278:   if (viewerVersion) {
279:     PetscBool flg;

281:     PetscCall(PetscViewerPrintVersion_Private(viewer, viewerVersion, viewerVersionStr, 16));
282:     PetscCall(PetscStrcmp(versionStr, viewerVersionStr, &flg));
283:     PetscCheck(flg, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "User requested DMPlex storage version %s but viewer already has version %s - cannot mix versions", versionStr, viewerVersionStr);
284:   }

286:   PetscCall(PetscViewerHDF5HasAttribute(viewer, NULL, ATTR_NAME, &fileHasVersion));
287:   if (fileHasVersion) {
288:     char     *fileVersion;
289:     PetscBool flg;

291:     PetscCall(PetscViewerHDF5ReadAttribute(viewer, "/", ATTR_NAME, PETSC_STRING, NULL, &fileVersion));
292:     PetscCall(PetscStrcmp(fileVersion, versionStr, &flg));
293:     PetscCheck(flg, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "User requested DMPlex storage version %s but file already has version %s - cannot mix versions", versionStr, fileVersion);
294:     PetscCall(PetscFree(fileVersion));
295:   }
296:   PetscCall(PetscNew(&viewerVersion));
297:   viewerVersion->major    = version->major;
298:   viewerVersion->minor    = version->minor;
299:   viewerVersion->subminor = version->subminor;
300:   PetscCall(PetscViewerAttachVersion_Private(viewer, DMPLEX_STORAGE_VERSION_READING_KEY, viewerVersion));
301:   PetscFunctionReturn(PETSC_SUCCESS);
302: }

304: /*@C
305:   PetscViewerHDF5GetDMPlexStorageVersionReading - Get the storage version for reading

307:   Logically collective

309:   Input Parameter:
310: . viewer - The `PetscViewer`

312:   Output Parameter:
313: . version - The storage format version

315:   Options Database Keys:
316: . -dm_plex_view_hdf5_storage_version <num> - Overrides the storage format version

318:   Level: advanced

320:   Note:
321:   The version has major, minor, and subminor integers. Parallel operations are only available for version 3.0.0.

323: .seealso: [](ch_dmbase), `DM`, `PetscViewerHDF5SetDMPlexStorageVersionReading()`, `PetscViewerHDF5GetDMPlexStorageVersionWriting()`, `PetscViewerHDF5SetDMPlexStorageVersionWriting()`
324: @*/
325: PetscErrorCode PetscViewerHDF5GetDMPlexStorageVersionReading(PetscViewer viewer, DMPlexStorageVersion *version)
326: {
327:   const char ATTR_NAME[] = "dmplex_storage_version";
328:   char      *defaultVersion;
329:   char      *versionString;

331:   PetscFunctionBegin;
333:   PetscAssertPointer(version, 2);
334:   PetscCall(PetscViewerGetAttachedVersion_Private(viewer, DMPLEX_STORAGE_VERSION_READING_KEY, version));
335:   if (*version) PetscFunctionReturn(PETSC_SUCCESS);

337:   //TODO string HDF5 attribute handling is terrible and should be redesigned
338:   PetscCall(PetscStrallocpy(DMPLEX_STORAGE_VERSION_FIRST, &defaultVersion));
339:   PetscCall(PetscViewerHDF5ReadAttribute(viewer, "/", ATTR_NAME, PETSC_STRING, &defaultVersion, &versionString));
340:   PetscCall(PetscViewerParseVersion_Private(viewer, versionString, version));
341:   PetscCall(PetscViewerAttachVersion_Private(viewer, DMPLEX_STORAGE_VERSION_READING_KEY, *version));
342:   PetscCall(PetscFree(versionString));
343:   PetscCall(PetscFree(defaultVersion));
344:   PetscFunctionReturn(PETSC_SUCCESS);
345: }

347: static PetscErrorCode DMPlexGetHDF5Name_Private(DM dm, const char *name[])
348: {
349:   PetscFunctionBegin;
350:   if (((PetscObject)dm)->name) {
351:     PetscCall(PetscObjectGetName((PetscObject)dm, name));
352:   } else {
353:     *name = "plex";
354:   }
355:   PetscFunctionReturn(PETSC_SUCCESS);
356: }

358: PetscErrorCode DMSequenceGetLength_HDF5_Internal(DM dm, const char seqname[], PetscInt *seqlen, PetscViewer viewer)
359: {
360:   hid_t     file, group, dset, dspace;
361:   hsize_t   rdim, *dims;
362:   char     *groupname;
363:   PetscBool has;

365:   PetscFunctionBegin;
366:   PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &groupname));
367:   PetscCall(PetscViewerHDF5HasDataset(viewer, seqname, &has));
368:   PetscCheck(has, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", seqname, groupname);

370:   PetscCall(PetscViewerHDF5OpenGroup(viewer, NULL, &file, &group));
371:   PetscCallHDF5Return(dset, H5Dopen2, (group, seqname, H5P_DEFAULT));
372:   PetscCallHDF5Return(dspace, H5Dget_space, (dset));
373:   PetscCallHDF5ReturnNoCheck(rdim, H5Sget_simple_extent_dims, (dspace, NULL, NULL));
374:   PetscCall(PetscMalloc1(rdim, &dims));
375:   PetscCallHDF5ReturnNoCheck(rdim, H5Sget_simple_extent_dims, (dspace, dims, NULL));
376:   *seqlen = (PetscInt)dims[0];
377:   PetscCall(PetscFree(dims));
378:   PetscCallHDF5(H5Dclose, (dset));
379:   PetscCallHDF5(H5Gclose, (group));
380:   PetscFunctionReturn(PETSC_SUCCESS);
381: }

383: static PetscErrorCode DMSequenceView_HDF5(DM dm, const char seqname[], PetscInt seqnum, PetscScalar value, PetscViewer viewer)
384: {
385:   Vec         stamp;
386:   PetscMPIInt rank;

388:   PetscFunctionBegin;
389:   if (seqnum < 0) PetscFunctionReturn(PETSC_SUCCESS);
390:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
391:   PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)viewer), rank ? 0 : 1, 1, &stamp));
392:   PetscCall(VecSetBlockSize(stamp, 1));
393:   PetscCall(PetscObjectSetName((PetscObject)stamp, seqname));
394:   if (rank == 0) {
395:     PetscReal timeScale;
396:     PetscBool istime;

398:     PetscCall(PetscStrncmp(seqname, "time", 5, &istime));
399:     if (istime) {
400:       PetscCall(DMPlexGetScale(dm, PETSC_UNIT_TIME, &timeScale));
401:       value *= timeScale;
402:     }
403:     PetscCall(VecSetValue(stamp, 0, value, INSERT_VALUES));
404:   }
405:   PetscCall(VecAssemblyBegin(stamp));
406:   PetscCall(VecAssemblyEnd(stamp));
407:   PetscCall(PetscViewerHDF5PushGroup(viewer, "/"));
408:   PetscCall(PetscViewerHDF5PushTimestepping(viewer));
409:   PetscCall(PetscViewerHDF5SetTimestep(viewer, seqnum)); /* seqnum < 0 jumps out above */
410:   PetscCall(VecView(stamp, viewer));
411:   PetscCall(PetscViewerHDF5PopTimestepping(viewer));
412:   PetscCall(PetscViewerHDF5PopGroup(viewer));
413:   PetscCall(VecDestroy(&stamp));
414:   PetscFunctionReturn(PETSC_SUCCESS);
415: }

417: PetscErrorCode DMSequenceLoad_HDF5_Internal(DM dm, const char seqname[], PetscInt seqnum, PetscScalar *value, PetscViewer viewer)
418: {
419:   Vec         stamp;
420:   PetscMPIInt rank;

422:   PetscFunctionBegin;
423:   if (seqnum < 0) PetscFunctionReturn(PETSC_SUCCESS);
424:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
425:   PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)viewer), rank ? 0 : 1, 1, &stamp));
426:   PetscCall(VecSetBlockSize(stamp, 1));
427:   PetscCall(PetscObjectSetName((PetscObject)stamp, seqname));
428:   PetscCall(PetscViewerHDF5PushGroup(viewer, "/"));
429:   PetscCall(PetscViewerHDF5PushTimestepping(viewer));
430:   PetscCall(PetscViewerHDF5SetTimestep(viewer, seqnum)); /* seqnum < 0 jumps out above */
431:   PetscCall(VecLoad(stamp, viewer));
432:   PetscCall(PetscViewerHDF5PopTimestepping(viewer));
433:   PetscCall(PetscViewerHDF5PopGroup(viewer));
434:   if (rank == 0) {
435:     const PetscScalar *a;
436:     PetscReal          timeScale;
437:     PetscBool          istime;

439:     PetscCall(VecGetArrayRead(stamp, &a));
440:     *value = a[0];
441:     PetscCall(VecRestoreArrayRead(stamp, &a));
442:     PetscCall(PetscStrncmp(seqname, "time", 5, &istime));
443:     if (istime) {
444:       PetscCall(DMPlexGetScale(dm, PETSC_UNIT_TIME, &timeScale));
445:       *value /= timeScale;
446:     }
447:   }
448:   PetscCall(VecDestroy(&stamp));
449:   PetscFunctionReturn(PETSC_SUCCESS);
450: }

452: static PetscErrorCode DMPlexCreateCutVertexLabel_Private(DM dm, DMLabel cutLabel, DMLabel *cutVertexLabel)
453: {
454:   IS              cutcells = NULL;
455:   const PetscInt *cutc;
456:   PetscInt        cellHeight, vStart, vEnd, cStart, cEnd, c;

458:   PetscFunctionBegin;
459:   if (!cutLabel) PetscFunctionReturn(PETSC_SUCCESS);
460:   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
461:   PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd));
462:   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
463:   /* Label vertices that should be duplicated */
464:   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Cut Vertices", cutVertexLabel));
465:   PetscCall(DMLabelGetStratumIS(cutLabel, 2, &cutcells));
466:   if (cutcells) {
467:     PetscInt n;

469:     PetscCall(ISGetIndices(cutcells, &cutc));
470:     PetscCall(ISGetLocalSize(cutcells, &n));
471:     for (c = 0; c < n; ++c) {
472:       if ((cutc[c] >= cStart) && (cutc[c] < cEnd)) {
473:         PetscInt *closure = NULL;
474:         PetscInt  closureSize, cl, value;

476:         PetscCall(DMPlexGetTransitiveClosure(dm, cutc[c], PETSC_TRUE, &closureSize, &closure));
477:         for (cl = 0; cl < closureSize * 2; cl += 2) {
478:           if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) {
479:             PetscCall(DMLabelGetValue(cutLabel, closure[cl], &value));
480:             if (value == 1) PetscCall(DMLabelSetValue(*cutVertexLabel, closure[cl], 1));
481:           }
482:         }
483:         PetscCall(DMPlexRestoreTransitiveClosure(dm, cutc[c], PETSC_TRUE, &closureSize, &closure));
484:       }
485:     }
486:     PetscCall(ISRestoreIndices(cutcells, &cutc));
487:   }
488:   PetscCall(ISDestroy(&cutcells));
489:   PetscFunctionReturn(PETSC_SUCCESS);
490: }

492: PetscErrorCode VecView_Plex_Local_HDF5_Internal(Vec v, PetscViewer viewer)
493: {
494:   DM                dm;
495:   DM                dmBC;
496:   PetscSection      section, sectionGlobal;
497:   Vec               gv;
498:   const char       *name;
499:   PetscViewerFormat format;
500:   PetscInt          seqnum;
501:   PetscReal         seqval;
502:   PetscBool         isseq;

504:   PetscFunctionBegin;
505:   PetscCall(PetscObjectTypeCompare((PetscObject)v, VECSEQ, &isseq));
506:   PetscCall(VecGetDM(v, &dm));
507:   PetscCall(DMGetLocalSection(dm, &section));
508:   PetscCall(DMGetOutputSequenceNumber(dm, &seqnum, &seqval));
509:   PetscCall(DMSequenceView_HDF5(dm, "time", seqnum, (PetscScalar)seqval, viewer));
510:   if (seqnum >= 0) {
511:     PetscCall(PetscViewerHDF5PushTimestepping(viewer));
512:     PetscCall(PetscViewerHDF5SetTimestep(viewer, seqnum));
513:   }
514:   PetscCall(PetscViewerGetFormat(viewer, &format));
515:   PetscCall(DMGetOutputDM(dm, &dmBC));
516:   PetscCall(DMGetGlobalSection(dmBC, &sectionGlobal));
517:   PetscCall(DMGetGlobalVector(dmBC, &gv));
518:   PetscCall(PetscObjectGetName((PetscObject)v, &name));
519:   PetscCall(PetscObjectSetName((PetscObject)gv, name));
520:   PetscCall(DMLocalToGlobalBegin(dmBC, v, INSERT_VALUES, gv));
521:   PetscCall(DMLocalToGlobalEnd(dmBC, v, INSERT_VALUES, gv));
522:   PetscCall(PetscObjectTypeCompare((PetscObject)gv, VECSEQ, &isseq));
523:   if (format == PETSC_VIEWER_HDF5_VIZ) {
524:     /* Output visualization representation */
525:     PetscInt numFields, f;
526:     DMLabel  cutLabel, cutVertexLabel = NULL;

528:     PetscCall(PetscSectionGetNumFields(section, &numFields));
529:     PetscCall(DMGetLabel(dm, "periodic_cut", &cutLabel));
530:     for (f = 0; f < numFields; ++f) {
531:       Vec                     subv;
532:       IS                      is;
533:       const char             *fname, *fgroup, *componentName, *fname_def = "unnamed";
534:       char                    subname[PETSC_MAX_PATH_LEN];
535:       PetscInt                Nc, c;
536:       PetscInt                pStart, pEnd;
537:       PetscViewerVTKFieldType ft;

539:       PetscCall(DMPlexGetFieldType_Internal(dm, section, f, &pStart, &pEnd, &ft));
540:       if (ft == PETSC_VTK_INVALID) continue;
541:       fgroup = (ft == PETSC_VTK_POINT_VECTOR_FIELD) || (ft == PETSC_VTK_POINT_FIELD) ? "/vertex_fields" : "/cell_fields";
542:       PetscCall(PetscSectionGetFieldName(section, f, &fname));
543:       if (!fname) fname = fname_def;

545:       PetscCall(PetscViewerHDF5PushGroup(viewer, fgroup));

547:       if (cutLabel) {
548:         const PetscScalar *ga;
549:         PetscScalar       *suba;
550:         PetscInt           gstart, subSize = 0, extSize = 0, subOff = 0, newOff = 0, p;

552:         PetscCall(DMPlexCreateCutVertexLabel_Private(dm, cutLabel, &cutVertexLabel));
553:         PetscCall(PetscSectionGetFieldComponents(section, f, &Nc));
554:         for (p = pStart; p < pEnd; ++p) {
555:           PetscInt gdof, fdof = 0, val;

557:           PetscCall(PetscSectionGetDof(sectionGlobal, p, &gdof));
558:           if (gdof > 0) PetscCall(PetscSectionGetFieldDof(section, p, f, &fdof));
559:           subSize += fdof;
560:           PetscCall(DMLabelGetValue(cutVertexLabel, p, &val));
561:           if (val == 1) extSize += fdof;
562:         }
563:         PetscCall(VecCreate(PetscObjectComm((PetscObject)gv), &subv));
564:         PetscCall(VecSetSizes(subv, subSize + extSize, PETSC_DETERMINE));
565:         PetscCall(VecSetBlockSize(subv, Nc));
566:         PetscCall(VecSetType(subv, VECSTANDARD));
567:         PetscCall(VecGetOwnershipRange(gv, &gstart, NULL));
568:         PetscCall(VecGetArrayRead(gv, &ga));
569:         PetscCall(VecGetArray(subv, &suba));
570:         for (p = pStart; p < pEnd; ++p) {
571:           PetscInt gdof, goff, val;

573:           PetscCall(PetscSectionGetDof(sectionGlobal, p, &gdof));
574:           if (gdof > 0) {
575:             PetscInt fdof, fc, f2, poff = 0;

577:             PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
578:             /* Can get rid of this loop by storing field information in the global section */
579:             for (f2 = 0; f2 < f; ++f2) {
580:               PetscCall(PetscSectionGetFieldDof(section, p, f2, &fdof));
581:               poff += fdof;
582:             }
583:             PetscCall(PetscSectionGetFieldDof(section, p, f, &fdof));
584:             for (fc = 0; fc < fdof; ++fc, ++subOff) suba[subOff] = ga[goff + poff + fc - gstart];
585:             PetscCall(DMLabelGetValue(cutVertexLabel, p, &val));
586:             if (val == 1) {
587:               for (fc = 0; fc < fdof; ++fc, ++newOff) suba[subSize + newOff] = ga[goff + poff + fc - gstart];
588:             }
589:           }
590:         }
591:         PetscCall(VecRestoreArrayRead(gv, &ga));
592:         PetscCall(VecRestoreArray(subv, &suba));
593:         PetscCall(DMLabelDestroy(&cutVertexLabel));
594:       } else {
595:         PetscCall(PetscSectionGetField_Internal(section, sectionGlobal, gv, f, pStart, pEnd, &is, &subv));
596:       }
597:       PetscCall(PetscStrncpy(subname, name, sizeof(subname)));
598:       PetscCall(PetscStrlcat(subname, "_", sizeof(subname)));
599:       PetscCall(PetscStrlcat(subname, fname, sizeof(subname)));
600:       PetscCall(PetscObjectSetName((PetscObject)subv, subname));
601:       if (isseq) PetscCall(VecView_Seq(subv, viewer));
602:       else PetscCall(VecView_MPI(subv, viewer));
603:       if (ft == PETSC_VTK_POINT_VECTOR_FIELD || ft == PETSC_VTK_CELL_VECTOR_FIELD) {
604:         PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)subv, "vector_field_type", PETSC_STRING, "vector"));
605:       } else {
606:         PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)subv, "vector_field_type", PETSC_STRING, "scalar"));
607:       }

609:       /* Output the component names in the field if available */
610:       PetscCall(PetscSectionGetFieldComponents(section, f, &Nc));
611:       for (c = 0; c < Nc; ++c) {
612:         char componentNameLabel[PETSC_MAX_PATH_LEN];
613:         PetscCall(PetscSectionGetComponentName(section, f, c, &componentName));
614:         PetscCall(PetscSNPrintf(componentNameLabel, sizeof(componentNameLabel), "componentName%" PetscInt_FMT, c));
615:         PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)subv, componentNameLabel, PETSC_STRING, componentName));
616:       }

618:       if (cutLabel) PetscCall(VecDestroy(&subv));
619:       else PetscCall(PetscSectionRestoreField_Internal(section, sectionGlobal, gv, f, pStart, pEnd, &is, &subv));
620:       PetscCall(PetscViewerHDF5PopGroup(viewer));
621:     }
622:   } else {
623:     /* Output full vector */
624:     PetscCall(PetscViewerHDF5PushGroup(viewer, "/fields"));
625:     if (isseq) PetscCall(VecView_Seq(gv, viewer));
626:     else PetscCall(VecView_MPI(gv, viewer));
627:     PetscCall(PetscViewerHDF5PopGroup(viewer));
628:   }
629:   if (seqnum >= 0) PetscCall(PetscViewerHDF5PopTimestepping(viewer));
630:   PetscCall(DMRestoreGlobalVector(dmBC, &gv));
631:   PetscFunctionReturn(PETSC_SUCCESS);
632: }

634: PetscErrorCode VecView_Plex_HDF5_Internal(Vec v, PetscViewer viewer)
635: {
636:   DM          dm;
637:   Vec         locv;
638:   PetscObject isZero;
639:   const char *name;
640:   PetscReal   time;

642:   PetscFunctionBegin;
643:   PetscCall(VecGetDM(v, &dm));
644:   PetscCall(DMGetLocalVector(dm, &locv));
645:   PetscCall(PetscObjectGetName((PetscObject)v, &name));
646:   PetscCall(PetscObjectSetName((PetscObject)locv, name));
647:   PetscCall(PetscObjectQuery((PetscObject)v, "__Vec_bc_zero__", &isZero));
648:   PetscCall(PetscObjectCompose((PetscObject)locv, "__Vec_bc_zero__", isZero));
649:   PetscCall(DMGlobalToLocalBegin(dm, v, INSERT_VALUES, locv));
650:   PetscCall(DMGlobalToLocalEnd(dm, v, INSERT_VALUES, locv));
651:   PetscCall(DMGetOutputSequenceNumber(dm, NULL, &time));
652:   PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locv, time, NULL, NULL, NULL));
653:   PetscCall(VecView_Plex_Local_HDF5_Internal(locv, viewer));
654:   PetscCall(PetscObjectCompose((PetscObject)locv, "__Vec_bc_zero__", NULL));
655:   PetscCall(DMRestoreLocalVector(dm, &locv));
656:   PetscFunctionReturn(PETSC_SUCCESS);
657: }

659: PetscErrorCode VecView_Plex_HDF5_Native_Internal(Vec v, PetscViewer viewer)
660: {
661:   PetscBool isseq;

663:   PetscFunctionBegin;
664:   PetscCall(PetscObjectTypeCompare((PetscObject)v, VECSEQ, &isseq));
665:   PetscCall(PetscViewerHDF5PushGroup(viewer, "/fields"));
666:   if (isseq) PetscCall(VecView_Seq(v, viewer));
667:   else PetscCall(VecView_MPI(v, viewer));
668:   PetscCall(PetscViewerHDF5PopGroup(viewer));
669:   PetscFunctionReturn(PETSC_SUCCESS);
670: }

672: PetscErrorCode VecLoad_Plex_HDF5_Internal(Vec v, PetscViewer viewer)
673: {
674:   DM          dm;
675:   Vec         locv;
676:   const char *name;
677:   PetscInt    seqnum;

679:   PetscFunctionBegin;
680:   PetscCall(VecGetDM(v, &dm));
681:   PetscCall(DMGetLocalVector(dm, &locv));
682:   PetscCall(PetscObjectGetName((PetscObject)v, &name));
683:   PetscCall(PetscObjectSetName((PetscObject)locv, name));
684:   PetscCall(DMGetOutputSequenceNumber(dm, &seqnum, NULL));
685:   PetscCall(PetscViewerHDF5PushGroup(viewer, "/fields"));
686:   if (seqnum >= 0) {
687:     PetscCall(PetscViewerHDF5PushTimestepping(viewer));
688:     PetscCall(PetscViewerHDF5SetTimestep(viewer, seqnum));
689:   }
690:   PetscCall(VecLoad_Plex_Local(locv, viewer));
691:   if (seqnum >= 0) PetscCall(PetscViewerHDF5PopTimestepping(viewer));
692:   PetscCall(PetscViewerHDF5PopGroup(viewer));
693:   PetscCall(DMLocalToGlobalBegin(dm, locv, INSERT_VALUES, v));
694:   PetscCall(DMLocalToGlobalEnd(dm, locv, INSERT_VALUES, v));
695:   PetscCall(DMRestoreLocalVector(dm, &locv));
696:   PetscFunctionReturn(PETSC_SUCCESS);
697: }

699: PetscErrorCode VecLoad_Plex_HDF5_Native_Internal(Vec v, PetscViewer viewer)
700: {
701:   DM       dm;
702:   PetscInt seqnum;

704:   PetscFunctionBegin;
705:   PetscCall(VecGetDM(v, &dm));
706:   PetscCall(DMGetOutputSequenceNumber(dm, &seqnum, NULL));
707:   PetscCall(PetscViewerHDF5PushGroup(viewer, "/fields"));
708:   if (seqnum >= 0) {
709:     PetscCall(PetscViewerHDF5PushTimestepping(viewer));
710:     PetscCall(PetscViewerHDF5SetTimestep(viewer, seqnum));
711:   }
712:   PetscCall(VecLoad_Default(v, viewer));
713:   if (seqnum >= 0) PetscCall(PetscViewerHDF5PopTimestepping(viewer));
714:   PetscCall(PetscViewerHDF5PopGroup(viewer));
715:   PetscFunctionReturn(PETSC_SUCCESS);
716: }

718: static PetscErrorCode DMPlexDistributionView_HDF5_Private(DM dm, IS globalPointNumbers, PetscViewer viewer)
719: {
720:   MPI_Comm           comm;
721:   PetscMPIInt        size, rank;
722:   PetscInt           size_petsc_int;
723:   const char        *topologydm_name, *distribution_name;
724:   const PetscInt    *gpoint;
725:   PetscInt           pStart, pEnd, p;
726:   PetscSF            pointSF;
727:   PetscInt           nroots, nleaves;
728:   const PetscInt    *ilocal;
729:   const PetscSFNode *iremote;
730:   IS                 chartSizesIS, ownersIS, gpointsIS;
731:   PetscInt          *chartSize, *owners, *gpoints;

733:   PetscFunctionBegin;
734:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
735:   PetscCallMPI(MPI_Comm_size(comm, &size));
736:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
737:   PetscCall(DMPlexGetHDF5Name_Private(dm, &topologydm_name));
738:   PetscCall(DMPlexDistributionGetName(dm, &distribution_name));
739:   if (!distribution_name) PetscFunctionReturn(PETSC_SUCCESS);
740:   PetscCall(PetscLogEventBegin(DMPLEX_DistributionView, viewer, 0, 0, 0));
741:   size_petsc_int = (PetscInt)size;
742:   PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "comm_size", PETSC_INT, (void *)&size_petsc_int));
743:   PetscCall(ISGetIndices(globalPointNumbers, &gpoint));
744:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
745:   PetscCall(PetscMalloc1(1, &chartSize));
746:   *chartSize = pEnd - pStart;
747:   PetscCall(PetscMalloc1(*chartSize, &owners));
748:   PetscCall(PetscMalloc1(*chartSize, &gpoints));
749:   PetscCall(DMGetPointSF(dm, &pointSF));
750:   PetscCall(PetscSFGetGraph(pointSF, &nroots, &nleaves, &ilocal, &iremote));
751:   for (p = pStart; p < pEnd; ++p) {
752:     PetscInt gp = gpoint[p - pStart];

754:     owners[p - pStart]  = rank;
755:     gpoints[p - pStart] = (gp < 0 ? -(gp + 1) : gp);
756:   }
757:   for (p = 0; p < nleaves; ++p) {
758:     PetscInt ilocalp = (ilocal ? ilocal[p] : p);

760:     owners[ilocalp] = iremote[p].rank;
761:   }
762:   PetscCall(ISCreateGeneral(comm, 1, chartSize, PETSC_OWN_POINTER, &chartSizesIS));
763:   PetscCall(ISCreateGeneral(comm, *chartSize, owners, PETSC_OWN_POINTER, &ownersIS));
764:   PetscCall(ISCreateGeneral(comm, *chartSize, gpoints, PETSC_OWN_POINTER, &gpointsIS));
765:   PetscCall(PetscObjectSetName((PetscObject)chartSizesIS, "chart_sizes"));
766:   PetscCall(PetscObjectSetName((PetscObject)ownersIS, "owners"));
767:   PetscCall(PetscObjectSetName((PetscObject)gpointsIS, "global_point_numbers"));
768:   PetscCall(ISView(chartSizesIS, viewer));
769:   PetscCall(ISView(ownersIS, viewer));
770:   PetscCall(ISView(gpointsIS, viewer));
771:   PetscCall(ISDestroy(&chartSizesIS));
772:   PetscCall(ISDestroy(&ownersIS));
773:   PetscCall(ISDestroy(&gpointsIS));
774:   PetscCall(ISRestoreIndices(globalPointNumbers, &gpoint));
775:   PetscCall(PetscLogEventEnd(DMPLEX_DistributionView, viewer, 0, 0, 0));
776:   PetscFunctionReturn(PETSC_SUCCESS);
777: }

779: static PetscErrorCode DMPlexTopologyView_HDF5_Inner_Private(DM dm, IS globalPointNumbers, PetscViewer viewer, PetscInt pStart, PetscInt pEnd, const char pointsName[], const char coneSizesName[], const char conesName[], const char orientationsName[])
780: {
781:   IS              coneSizesIS, conesIS, orientationsIS;
782:   PetscInt       *coneSizes, *cones, *orientations;
783:   const PetscInt *gpoint;
784:   PetscInt        nPoints = 0, conesSize = 0;
785:   PetscInt        p, c, s;
786:   MPI_Comm        comm;

788:   PetscFunctionBegin;
789:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
790:   PetscCall(ISGetIndices(globalPointNumbers, &gpoint));
791:   for (p = pStart; p < pEnd; ++p) {
792:     if (gpoint[p] >= 0) {
793:       PetscInt coneSize;

795:       PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
796:       nPoints += 1;
797:       conesSize += coneSize;
798:     }
799:   }
800:   PetscCall(PetscMalloc1(nPoints, &coneSizes));
801:   PetscCall(PetscMalloc1(conesSize, &cones));
802:   PetscCall(PetscMalloc1(conesSize, &orientations));
803:   for (p = pStart, c = 0, s = 0; p < pEnd; ++p) {
804:     if (gpoint[p] >= 0) {
805:       const PetscInt *cone, *ornt;
806:       PetscInt        coneSize, cp;

808:       PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
809:       PetscCall(DMPlexGetCone(dm, p, &cone));
810:       PetscCall(DMPlexGetConeOrientation(dm, p, &ornt));
811:       coneSizes[s] = coneSize;
812:       for (cp = 0; cp < coneSize; ++cp, ++c) {
813:         cones[c]        = gpoint[cone[cp]] < 0 ? -(gpoint[cone[cp]] + 1) : gpoint[cone[cp]];
814:         orientations[c] = ornt[cp];
815:       }
816:       ++s;
817:     }
818:   }
819:   PetscCheck(s == nPoints, PETSC_COMM_SELF, PETSC_ERR_LIB, "Total number of points %" PetscInt_FMT " != %" PetscInt_FMT, s, nPoints);
820:   PetscCheck(c == conesSize, PETSC_COMM_SELF, PETSC_ERR_LIB, "Total number of cone points %" PetscInt_FMT " != %" PetscInt_FMT, c, conesSize);
821:   PetscCall(ISCreateGeneral(comm, nPoints, coneSizes, PETSC_OWN_POINTER, &coneSizesIS));
822:   PetscCall(ISCreateGeneral(comm, conesSize, cones, PETSC_OWN_POINTER, &conesIS));
823:   PetscCall(ISCreateGeneral(comm, conesSize, orientations, PETSC_OWN_POINTER, &orientationsIS));
824:   PetscCall(PetscObjectSetName((PetscObject)coneSizesIS, coneSizesName));
825:   PetscCall(PetscObjectSetName((PetscObject)conesIS, conesName));
826:   PetscCall(PetscObjectSetName((PetscObject)orientationsIS, orientationsName));
827:   PetscCall(ISView(coneSizesIS, viewer));
828:   PetscCall(ISView(conesIS, viewer));
829:   PetscCall(ISView(orientationsIS, viewer));
830:   PetscCall(ISDestroy(&coneSizesIS));
831:   PetscCall(ISDestroy(&conesIS));
832:   PetscCall(ISDestroy(&orientationsIS));
833:   if (pointsName) {
834:     IS        pointsIS;
835:     PetscInt *points;

837:     PetscCall(PetscMalloc1(nPoints, &points));
838:     for (p = pStart, c = 0, s = 0; p < pEnd; ++p) {
839:       if (gpoint[p] >= 0) {
840:         points[s] = gpoint[p];
841:         ++s;
842:       }
843:     }
844:     PetscCall(ISCreateGeneral(comm, nPoints, points, PETSC_OWN_POINTER, &pointsIS));
845:     PetscCall(PetscObjectSetName((PetscObject)pointsIS, pointsName));
846:     PetscCall(ISView(pointsIS, viewer));
847:     PetscCall(ISDestroy(&pointsIS));
848:   }
849:   PetscCall(ISRestoreIndices(globalPointNumbers, &gpoint));
850:   PetscFunctionReturn(PETSC_SUCCESS);
851: }

853: static PetscErrorCode DMPlexTopologyView_HDF5_Legacy_Private(DM dm, IS globalPointNumbers, PetscViewer viewer)
854: {
855:   const char *pointsName, *coneSizesName, *conesName, *orientationsName;
856:   PetscInt    pStart, pEnd, dim;

858:   PetscFunctionBegin;
859:   pointsName       = "order";
860:   coneSizesName    = "cones";
861:   conesName        = "cells";
862:   orientationsName = "orientation";
863:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
864:   PetscCall(DMPlexTopologyView_HDF5_Inner_Private(dm, globalPointNumbers, viewer, pStart, pEnd, pointsName, coneSizesName, conesName, orientationsName));
865:   PetscCall(DMGetDimension(dm, &dim));
866:   PetscCall(PetscViewerHDF5WriteAttribute(viewer, conesName, "cell_dim", PETSC_INT, (void *)&dim));
867:   PetscFunctionReturn(PETSC_SUCCESS);
868: }

870: //TODO get this numbering right away without needing this function
871: /* Renumber global point numbers so that they are 0-based per stratum */
872: static PetscErrorCode RenumberGlobalPointNumbersPerStratum_Private(DM dm, IS globalPointNumbers, IS *newGlobalPointNumbers, IS *strataPermutation)
873: {
874:   PetscInt        d, depth, p, n;
875:   PetscInt       *offsets;
876:   const PetscInt *gpn;
877:   PetscInt       *ngpn;
878:   MPI_Comm        comm;

880:   PetscFunctionBegin;
881:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
882:   PetscCall(ISGetLocalSize(globalPointNumbers, &n));
883:   PetscCall(ISGetIndices(globalPointNumbers, &gpn));
884:   PetscCall(PetscMalloc1(n, &ngpn));
885:   PetscCall(DMPlexGetDepth(dm, &depth));
886:   PetscCall(PetscMalloc1(depth + 1, &offsets));
887:   for (d = 0; d <= depth; d++) {
888:     PetscInt pStart, pEnd;

890:     PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd));
891:     offsets[d] = PETSC_INT_MAX;
892:     for (p = pStart; p < pEnd; p++) {
893:       if (gpn[p] >= 0 && gpn[p] < offsets[d]) offsets[d] = gpn[p];
894:     }
895:   }
896:   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, offsets, depth + 1, MPIU_INT, MPI_MIN, comm));
897:   for (d = 0; d <= depth; d++) {
898:     PetscInt pStart, pEnd;

900:     PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd));
901:     for (p = pStart; p < pEnd; p++) ngpn[p] = gpn[p] - PetscSign(gpn[p]) * offsets[d];
902:   }
903:   PetscCall(ISRestoreIndices(globalPointNumbers, &gpn));
904:   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)globalPointNumbers), n, ngpn, PETSC_OWN_POINTER, newGlobalPointNumbers));
905:   {
906:     PetscInt *perm;

908:     PetscCall(PetscMalloc1(depth + 1, &perm));
909:     for (d = 0; d <= depth; d++) perm[d] = d;
910:     PetscCall(PetscSortIntWithPermutation(depth + 1, offsets, perm));
911:     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, depth + 1, perm, PETSC_OWN_POINTER, strataPermutation));
912:   }
913:   PetscCall(PetscFree(offsets));
914:   PetscFunctionReturn(PETSC_SUCCESS);
915: }

917: static PetscErrorCode DMPlexTopologyView_HDF5_Private(DM dm, IS globalPointNumbers, PetscViewer viewer)
918: {
919:   IS          globalPointNumbers0, strataPermutation;
920:   const char *coneSizesName, *conesName, *orientationsName;
921:   PetscInt    depth, d;
922:   MPI_Comm    comm;

924:   PetscFunctionBegin;
925:   coneSizesName    = "cone_sizes";
926:   conesName        = "cones";
927:   orientationsName = "orientations";
928:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
929:   PetscCall(DMPlexGetDepth(dm, &depth));
930:   {
931:     PetscInt dim;

933:     PetscCall(DMGetDimension(dm, &dim));
934:     PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "cell_dim", PETSC_INT, &dim));
935:     PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "depth", PETSC_INT, &depth));
936:   }

938:   PetscCall(RenumberGlobalPointNumbersPerStratum_Private(dm, globalPointNumbers, &globalPointNumbers0, &strataPermutation));
939:   /* TODO dirty trick to save serial IS using the same parallel viewer */
940:   {
941:     IS              spOnComm;
942:     PetscInt        n   = 0, N;
943:     const PetscInt *idx = NULL;
944:     const PetscInt *old;
945:     PetscMPIInt     rank;

947:     PetscCallMPI(MPI_Comm_rank(comm, &rank));
948:     PetscCall(ISGetLocalSize(strataPermutation, &N));
949:     PetscCall(ISGetIndices(strataPermutation, &old));
950:     if (!rank) {
951:       n   = N;
952:       idx = old;
953:     }
954:     PetscCall(ISCreateGeneral(comm, n, idx, PETSC_COPY_VALUES, &spOnComm));
955:     PetscCall(ISRestoreIndices(strataPermutation, &old));
956:     PetscCall(ISDestroy(&strataPermutation));
957:     strataPermutation = spOnComm;
958:   }
959:   PetscCall(PetscObjectSetName((PetscObject)strataPermutation, "permutation"));
960:   PetscCall(ISView(strataPermutation, viewer));
961:   PetscCall(PetscViewerHDF5PushGroup(viewer, "strata"));
962:   for (d = 0; d <= depth; d++) {
963:     PetscInt pStart, pEnd;
964:     char     group[128];

966:     PetscCall(PetscSNPrintf(group, sizeof(group), "%" PetscInt_FMT, d));
967:     PetscCall(PetscViewerHDF5PushGroup(viewer, group));
968:     PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd));
969:     PetscCall(DMPlexTopologyView_HDF5_Inner_Private(dm, globalPointNumbers0, viewer, pStart, pEnd, NULL, coneSizesName, conesName, orientationsName));
970:     PetscCall(PetscViewerHDF5PopGroup(viewer));
971:   }
972:   PetscCall(PetscViewerHDF5PopGroup(viewer)); /* strata */
973:   PetscCall(ISDestroy(&globalPointNumbers0));
974:   PetscCall(ISDestroy(&strataPermutation));
975:   PetscFunctionReturn(PETSC_SUCCESS);
976: }

978: PetscErrorCode DMPlexTopologyView_HDF5_Internal(DM dm, IS globalPointNumbers, PetscViewer viewer)
979: {
980:   DMPlexStorageVersion version;
981:   const char          *topologydm_name;
982:   char                 group[PETSC_MAX_PATH_LEN];

984:   PetscFunctionBegin;
985:   PetscCall(PetscViewerHDF5GetDMPlexStorageVersionWriting(viewer, &version));
986:   PetscCall(PetscInfo(dm, "Writing DM %s storage version %d.%d.%d\n", dm->hdr.name, version->major, version->minor, version->subminor));
987:   PetscCall(DMPlexGetHDF5Name_Private(dm, &topologydm_name));
988:   if (DMPlexStorageVersionGE(version, 2, 0, 0)) {
989:     PetscCall(PetscSNPrintf(group, sizeof(group), "topologies/%s", topologydm_name));
990:   } else {
991:     PetscCall(PetscStrncpy(group, "/", sizeof(group)));
992:   }
993:   PetscCall(PetscViewerHDF5PushGroup(viewer, group));

995:   PetscCall(PetscViewerHDF5PushGroup(viewer, "topology"));
996:   if (version->major < 3) {
997:     PetscCall(DMPlexTopologyView_HDF5_Legacy_Private(dm, globalPointNumbers, viewer));
998:   } else {
999:     /* since DMPlexStorageVersion 3.0.0 */
1000:     PetscCall(DMPlexTopologyView_HDF5_Private(dm, globalPointNumbers, viewer));
1001:   }
1002:   PetscCall(PetscViewerHDF5PopGroup(viewer)); /* "topology" */

1004:   if (DMPlexStorageVersionGE(version, 2, 1, 0)) {
1005:     const char *distribution_name;

1007:     PetscCall(DMPlexDistributionGetName(dm, &distribution_name));
1008:     PetscCall(PetscViewerHDF5PushGroup(viewer, "distributions"));
1009:     PetscCall(PetscViewerHDF5WriteGroup(viewer, NULL));
1010:     PetscCall(PetscViewerHDF5PushGroup(viewer, distribution_name));
1011:     PetscCall(DMPlexDistributionView_HDF5_Private(dm, globalPointNumbers, viewer));
1012:     PetscCall(PetscViewerHDF5PopGroup(viewer)); /* distribution_name */
1013:     PetscCall(PetscViewerHDF5PopGroup(viewer)); /* "distributions" */
1014:   }

1016:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1017:   PetscFunctionReturn(PETSC_SUCCESS);
1018: }

1020: static PetscErrorCode CreateConesIS_Private(DM dm, PetscInt cStart, PetscInt cEnd, IS globalCellNumbers, PetscInt *numCorners, IS *cellIS)
1021: {
1022:   PetscSF         sfPoint;
1023:   DMLabel         cutLabel, cutVertexLabel         = NULL;
1024:   IS              globalVertexNumbers, cutvertices = NULL;
1025:   const PetscInt *gcell, *gvertex, *cutverts = NULL;
1026:   PetscInt       *vertices;
1027:   PetscInt        conesSize = 0;
1028:   PetscInt        dim, numCornersLocal = 0, cell, vStart, vEnd, vExtra = 0, v;

1030:   PetscFunctionBegin;
1031:   *numCorners = 0;
1032:   PetscCall(DMGetDimension(dm, &dim));
1033:   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1034:   PetscCall(ISGetIndices(globalCellNumbers, &gcell));

1036:   for (cell = cStart; cell < cEnd; ++cell) {
1037:     PetscInt *closure = NULL;
1038:     PetscInt  closureSize, v, Nc = 0;

1040:     if (gcell[cell] < 0) continue;
1041:     PetscCall(DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure));
1042:     for (v = 0; v < closureSize * 2; v += 2) {
1043:       if ((closure[v] >= vStart) && (closure[v] < vEnd)) ++Nc;
1044:     }
1045:     PetscCall(DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure));
1046:     conesSize += Nc;
1047:     if (!numCornersLocal) numCornersLocal = Nc;
1048:     else if (numCornersLocal != Nc) numCornersLocal = 1;
1049:   }
1050:   PetscCallMPI(MPIU_Allreduce(&numCornersLocal, numCorners, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
1051:   PetscCheck(!numCornersLocal || !(numCornersLocal != *numCorners || *numCorners == 0), PETSC_COMM_SELF, PETSC_ERR_SUP, "Visualization topology currently only supports identical cell shapes");
1052:   /* Handle periodic cuts by identifying vertices which should be duplicated */
1053:   PetscCall(DMGetLabel(dm, "periodic_cut", &cutLabel));
1054:   PetscCall(DMPlexCreateCutVertexLabel_Private(dm, cutLabel, &cutVertexLabel));
1055:   if (cutVertexLabel) PetscCall(DMLabelGetStratumIS(cutVertexLabel, 1, &cutvertices));
1056:   if (cutvertices) {
1057:     PetscCall(ISGetIndices(cutvertices, &cutverts));
1058:     PetscCall(ISGetLocalSize(cutvertices, &vExtra));
1059:   }
1060:   PetscCall(DMGetPointSF(dm, &sfPoint));
1061:   if (cutLabel) {
1062:     const PetscInt    *ilocal;
1063:     const PetscSFNode *iremote;
1064:     PetscInt           nroots, nleaves;

1066:     PetscCall(PetscSFGetGraph(sfPoint, &nroots, &nleaves, &ilocal, &iremote));
1067:     if (nleaves < 0) {
1068:       PetscCall(PetscObjectReference((PetscObject)sfPoint));
1069:     } else {
1070:       PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sfPoint), &sfPoint));
1071:       PetscCall(PetscSFSetGraph(sfPoint, nroots + vExtra, nleaves, (PetscInt *)ilocal, PETSC_COPY_VALUES, (PetscSFNode *)iremote, PETSC_COPY_VALUES));
1072:     }
1073:   } else {
1074:     PetscCall(PetscObjectReference((PetscObject)sfPoint));
1075:   }
1076:   /* Number all vertices */
1077:   PetscCall(DMPlexCreateNumbering_Plex(dm, vStart, vEnd + vExtra, 0, NULL, sfPoint, &globalVertexNumbers));
1078:   PetscCall(PetscSFDestroy(&sfPoint));
1079:   /* Create cones */
1080:   PetscCall(ISGetIndices(globalVertexNumbers, &gvertex));
1081:   PetscCall(PetscMalloc1(conesSize, &vertices));
1082:   for (cell = cStart, v = 0; cell < cEnd; ++cell) {
1083:     PetscInt *closure = NULL;
1084:     PetscInt  closureSize, Nc = 0, p, value = -1;
1085:     PetscBool replace;

1087:     if (gcell[cell] < 0) continue;
1088:     if (cutLabel) PetscCall(DMLabelGetValue(cutLabel, cell, &value));
1089:     replace = (value == 2) ? PETSC_TRUE : PETSC_FALSE;
1090:     PetscCall(DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure));
1091:     for (p = 0; p < closureSize * 2; p += 2) {
1092:       if ((closure[p] >= vStart) && (closure[p] < vEnd)) closure[Nc++] = closure[p];
1093:     }
1094:     PetscCall(DMPlexReorderCell(dm, cell, closure));
1095:     for (p = 0; p < Nc; ++p) {
1096:       PetscInt nv, gv = gvertex[closure[p] - vStart];

1098:       if (replace) {
1099:         PetscCall(PetscFindInt(closure[p], vExtra, cutverts, &nv));
1100:         if (nv >= 0) gv = gvertex[vEnd - vStart + nv];
1101:       }
1102:       vertices[v++] = gv < 0 ? -(gv + 1) : gv;
1103:     }
1104:     PetscCall(DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure));
1105:   }
1106:   PetscCall(ISRestoreIndices(globalVertexNumbers, &gvertex));
1107:   PetscCall(ISDestroy(&globalVertexNumbers));
1108:   PetscCall(ISRestoreIndices(globalCellNumbers, &gcell));
1109:   if (cutvertices) PetscCall(ISRestoreIndices(cutvertices, &cutverts));
1110:   PetscCall(ISDestroy(&cutvertices));
1111:   PetscCall(DMLabelDestroy(&cutVertexLabel));
1112:   PetscCheck(v == conesSize, PETSC_COMM_SELF, PETSC_ERR_LIB, "Total number of cell vertices %" PetscInt_FMT " != %" PetscInt_FMT, v, conesSize);
1113:   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), conesSize, vertices, PETSC_OWN_POINTER, cellIS));
1114:   PetscCall(PetscLayoutSetBlockSize((*cellIS)->map, *numCorners));
1115:   PetscCall(PetscObjectSetName((PetscObject)*cellIS, "cells"));
1116:   PetscFunctionReturn(PETSC_SUCCESS);
1117: }

1119: static PetscErrorCode DMPlexTopologyView_HDF5_XDMF_Private(DM dm, IS globalCellNumbers, PetscViewer viewer)
1120: {
1121:   DM       cdm;
1122:   DMLabel  depthLabel, ctLabel;
1123:   IS       cellIS;
1124:   PetscInt dim, depth, cellHeight, c, n = 0;

1126:   PetscFunctionBegin;
1127:   PetscCall(PetscViewerHDF5PushGroup(viewer, "/viz"));
1128:   PetscCall(PetscViewerHDF5WriteGroup(viewer, NULL));

1130:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1131:   PetscCall(DMGetDimension(dm, &dim));
1132:   PetscCall(DMPlexGetDepth(dm, &depth));
1133:   PetscCall(DMGetCoordinateDM(dm, &cdm));
1134:   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
1135:   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
1136:   PetscCall(DMPlexGetCellTypeLabel(dm, &ctLabel));
1137:   for (c = 0; c < DM_NUM_POLYTOPES; ++c) {
1138:     const DMPolytopeType ict = (DMPolytopeType)c;
1139:     PetscInt             pStart, pEnd, dep, numCorners;
1140:     PetscBool            output = PETSC_FALSE, doOutput;

1142:     if (ict == DM_POLYTOPE_FV_GHOST) continue;
1143:     PetscCall(DMLabelGetStratumBounds(ctLabel, ict, &pStart, &pEnd));
1144:     if (pStart >= 0) {
1145:       PetscCall(DMLabelGetValue(depthLabel, pStart, &dep));
1146:       if (dep == depth - cellHeight) output = PETSC_TRUE;
1147:     }
1148:     PetscCallMPI(MPIU_Allreduce(&output, &doOutput, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)dm)));
1149:     if (!doOutput) continue;
1150:     PetscCall(CreateConesIS_Private(dm, pStart, pEnd, globalCellNumbers, &numCorners, &cellIS));
1151:     if (!n) {
1152:       PetscCall(PetscViewerHDF5PushGroup(viewer, "/viz/topology"));
1153:     } else {
1154:       char group[PETSC_MAX_PATH_LEN];

1156:       PetscCall(PetscSNPrintf(group, PETSC_MAX_PATH_LEN, "/viz/topology_%" PetscInt_FMT, n));
1157:       PetscCall(PetscViewerHDF5PushGroup(viewer, group));
1158:     }
1159:     PetscCall(ISView(cellIS, viewer));
1160:     PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)cellIS, "cell_corners", PETSC_INT, (void *)&numCorners));
1161:     PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)cellIS, "cell_dim", PETSC_INT, (void *)&dim));
1162:     PetscCall(ISDestroy(&cellIS));
1163:     PetscCall(PetscViewerHDF5PopGroup(viewer));
1164:     ++n;
1165:   }
1166:   PetscFunctionReturn(PETSC_SUCCESS);
1167: }

1169: static PetscErrorCode DMPlexCoordinatesView_HDF5_Legacy_Private(DM dm, PetscViewer viewer)
1170: {
1171:   DM        cdm;
1172:   Vec       coordinates, newcoords;
1173:   PetscReal lengthScale;
1174:   PetscInt  m, M, bs;

1176:   PetscFunctionBegin;
1177:   PetscCall(DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale));
1178:   PetscCall(DMGetCoordinateDM(dm, &cdm));
1179:   PetscCall(DMGetCoordinates(dm, &coordinates));
1180:   PetscCall(VecCreate(PetscObjectComm((PetscObject)coordinates), &newcoords));
1181:   PetscCall(PetscObjectSetName((PetscObject)newcoords, "vertices"));
1182:   PetscCall(VecGetSize(coordinates, &M));
1183:   PetscCall(VecGetLocalSize(coordinates, &m));
1184:   PetscCall(VecSetSizes(newcoords, m, M));
1185:   PetscCall(VecGetBlockSize(coordinates, &bs));
1186:   PetscCall(VecSetBlockSize(newcoords, bs));
1187:   PetscCall(VecSetType(newcoords, VECSTANDARD));
1188:   PetscCall(VecCopy(coordinates, newcoords));
1189:   PetscCall(VecScale(newcoords, lengthScale));
1190:   /* Did not use DMGetGlobalVector() in order to bypass default group assignment */
1191:   PetscCall(PetscViewerHDF5PushGroup(viewer, "/geometry"));
1192:   PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_NATIVE));
1193:   PetscCall(VecView(newcoords, viewer));
1194:   PetscCall(PetscViewerPopFormat(viewer));
1195:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1196:   PetscCall(VecDestroy(&newcoords));
1197:   PetscFunctionReturn(PETSC_SUCCESS);
1198: }

1200: PetscErrorCode DMPlexCoordinatesView_HDF5_Internal(DM dm, PetscViewer viewer)
1201: {
1202:   DM          cdm;
1203:   Vec         coords, newcoords;
1204:   PetscInt    m, M, bs;
1205:   PetscReal   lengthScale;
1206:   const char *topologydm_name, *coordinatedm_name, *coordinates_name;

1208:   PetscFunctionBegin;
1209:   {
1210:     PetscViewerFormat    format;
1211:     DMPlexStorageVersion version;

1213:     PetscCall(PetscViewerGetFormat(viewer, &format));
1214:     PetscCall(PetscViewerHDF5GetDMPlexStorageVersionWriting(viewer, &version));
1215:     if (!DMPlexStorageVersionGE(version, 2, 0, 0) || format == PETSC_VIEWER_HDF5_XDMF || format == PETSC_VIEWER_HDF5_VIZ) {
1216:       PetscCall(DMPlexCoordinatesView_HDF5_Legacy_Private(dm, viewer));
1217:       PetscFunctionReturn(PETSC_SUCCESS);
1218:     }
1219:   }
1220:   /* since 2.0.0 */
1221:   PetscCall(DMGetCoordinateDM(dm, &cdm));
1222:   PetscCall(DMGetCoordinates(dm, &coords));
1223:   PetscCall(PetscObjectGetName((PetscObject)cdm, &coordinatedm_name));
1224:   PetscCall(PetscObjectGetName((PetscObject)coords, &coordinates_name));
1225:   PetscCall(DMPlexGetHDF5Name_Private(dm, &topologydm_name));
1226:   PetscCall(PetscViewerHDF5PushGroup(viewer, "topologies"));
1227:   PetscCall(PetscViewerHDF5PushGroup(viewer, topologydm_name));
1228:   PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "coordinateDMName", PETSC_STRING, coordinatedm_name));
1229:   PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "coordinatesName", PETSC_STRING, coordinates_name));
1230:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1231:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1232:   PetscCall(DMPlexSectionView(dm, viewer, cdm));
1233:   PetscCall(VecCreate(PetscObjectComm((PetscObject)coords), &newcoords));
1234:   PetscCall(PetscObjectSetName((PetscObject)newcoords, coordinates_name));
1235:   PetscCall(VecGetSize(coords, &M));
1236:   PetscCall(VecGetLocalSize(coords, &m));
1237:   PetscCall(VecSetSizes(newcoords, m, M));
1238:   PetscCall(VecGetBlockSize(coords, &bs));
1239:   PetscCall(VecSetBlockSize(newcoords, bs));
1240:   PetscCall(VecSetType(newcoords, VECSTANDARD));
1241:   PetscCall(VecCopy(coords, newcoords));
1242:   PetscCall(DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale));
1243:   PetscCall(VecScale(newcoords, lengthScale));
1244:   PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_NATIVE));
1245:   PetscCall(DMPlexGlobalVectorView(dm, viewer, cdm, newcoords));
1246:   PetscCall(PetscViewerPopFormat(viewer));
1247:   PetscCall(VecDestroy(&newcoords));
1248:   PetscFunctionReturn(PETSC_SUCCESS);
1249: }

1251: static PetscErrorCode DMPlexCoordinatesView_HDF5_XDMF_Private(DM dm, PetscViewer viewer)
1252: {
1253:   DM               cdm;
1254:   Vec              coordinatesLocal, newcoords;
1255:   PetscSection     cSection, cGlobalSection;
1256:   PetscScalar     *coords, *ncoords;
1257:   DMLabel          cutLabel, cutVertexLabel = NULL;
1258:   const PetscReal *L;
1259:   PetscReal        lengthScale;
1260:   PetscInt         vStart, vEnd, v, bs, N, coordSize, dof, off, d;
1261:   PetscBool        localized, embedded;

1263:   PetscFunctionBegin;
1264:   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1265:   PetscCall(DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale));
1266:   PetscCall(DMGetCoordinatesLocal(dm, &coordinatesLocal));
1267:   PetscCall(VecGetBlockSize(coordinatesLocal, &bs));
1268:   PetscCall(DMGetCoordinatesLocalized(dm, &localized));
1269:   if (localized == PETSC_FALSE) PetscFunctionReturn(PETSC_SUCCESS);
1270:   PetscCall(DMGetPeriodicity(dm, NULL, NULL, &L));
1271:   PetscCall(DMGetCoordinateDM(dm, &cdm));
1272:   PetscCall(DMGetLocalSection(cdm, &cSection));
1273:   PetscCall(DMGetGlobalSection(cdm, &cGlobalSection));
1274:   PetscCall(DMGetLabel(dm, "periodic_cut", &cutLabel));
1275:   N = 0;

1277:   PetscCall(DMPlexCreateCutVertexLabel_Private(dm, cutLabel, &cutVertexLabel));
1278:   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), &newcoords));
1279:   PetscCall(PetscSectionGetDof(cSection, vStart, &dof));
1280:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "DOF: %" PetscInt_FMT "\n", dof));
1281:   embedded = (PetscBool)(L && dof == 2 && !cutLabel);
1282:   if (cutVertexLabel) {
1283:     PetscCall(DMLabelGetStratumSize(cutVertexLabel, 1, &v));
1284:     N += dof * v;
1285:   }
1286:   for (v = vStart; v < vEnd; ++v) {
1287:     PetscCall(PetscSectionGetDof(cGlobalSection, v, &dof));
1288:     if (dof < 0) continue;
1289:     if (embedded) N += dof + 1;
1290:     else N += dof;
1291:   }
1292:   if (embedded) PetscCall(VecSetBlockSize(newcoords, bs + 1));
1293:   else PetscCall(VecSetBlockSize(newcoords, bs));
1294:   PetscCall(VecSetSizes(newcoords, N, PETSC_DETERMINE));
1295:   PetscCall(VecSetType(newcoords, VECSTANDARD));
1296:   PetscCall(VecGetArray(coordinatesLocal, &coords));
1297:   PetscCall(VecGetArray(newcoords, &ncoords));
1298:   coordSize = 0;
1299:   for (v = vStart; v < vEnd; ++v) {
1300:     PetscCall(PetscSectionGetDof(cGlobalSection, v, &dof));
1301:     PetscCall(PetscSectionGetOffset(cSection, v, &off));
1302:     if (dof < 0) continue;
1303:     if (embedded) {
1304:       if (L && (L[0] > 0.0) && (L[1] > 0.0)) {
1305:         PetscReal theta, phi, r, R;
1306:         /* XY-periodic */
1307:         /* Suppose its an y-z circle, then
1308:              \hat r = (0, cos(th), sin(th)) \hat x = (1, 0, 0)
1309:            and the circle in that plane is
1310:              \hat r cos(phi) + \hat x sin(phi) */
1311:         theta                = 2.0 * PETSC_PI * PetscRealPart(coords[off + 1]) / L[1];
1312:         phi                  = 2.0 * PETSC_PI * PetscRealPart(coords[off + 0]) / L[0];
1313:         r                    = L[0] / (2.0 * PETSC_PI * 2.0 * L[1]);
1314:         R                    = L[1] / (2.0 * PETSC_PI);
1315:         ncoords[coordSize++] = PetscSinReal(phi) * r;
1316:         ncoords[coordSize++] = -PetscCosReal(theta) * (R + r * PetscCosReal(phi));
1317:         ncoords[coordSize++] = PetscSinReal(theta) * (R + r * PetscCosReal(phi));
1318:       } else if (L && (L[0] > 0.0)) {
1319:         /* X-periodic */
1320:         ncoords[coordSize++] = -PetscCosReal(2.0 * PETSC_PI * PetscRealPart(coords[off + 0]) / L[0]) * (L[0] / (2.0 * PETSC_PI));
1321:         ncoords[coordSize++] = coords[off + 1];
1322:         ncoords[coordSize++] = PetscSinReal(2.0 * PETSC_PI * PetscRealPart(coords[off + 0]) / L[0]) * (L[0] / (2.0 * PETSC_PI));
1323:       } else if (L && (L[1] > 0.0)) {
1324:         /* Y-periodic */
1325:         ncoords[coordSize++] = coords[off + 0];
1326:         ncoords[coordSize++] = PetscSinReal(2.0 * PETSC_PI * PetscRealPart(coords[off + 1]) / L[1]) * (L[1] / (2.0 * PETSC_PI));
1327:         ncoords[coordSize++] = -PetscCosReal(2.0 * PETSC_PI * PetscRealPart(coords[off + 1]) / L[1]) * (L[1] / (2.0 * PETSC_PI));
1328:   #if 0
1329:       } else if ((bd[0] == DM_BOUNDARY_TWIST)) {
1330:         PetscReal phi, r, R;
1331:         /* Mobius strip */
1332:         /* Suppose its an x-z circle, then
1333:              \hat r = (-cos(phi), 0, sin(phi)) \hat y = (0, 1, 0)
1334:            and in that plane we rotate by pi as we go around the circle
1335:              \hat r cos(phi/2) + \hat y sin(phi/2) */
1336:         phi   = 2.0*PETSC_PI*PetscRealPart(coords[off+0])/L[0];
1337:         R     = L[0];
1338:         r     = PetscRealPart(coords[off+1]) - L[1]/2.0;
1339:         ncoords[coordSize++] = -PetscCosReal(phi) * (R + r * PetscCosReal(phi/2.0));
1340:         ncoords[coordSize++] =  PetscSinReal(phi/2.0) * r;
1341:         ncoords[coordSize++] =  PetscSinReal(phi) * (R + r * PetscCosReal(phi/2.0));
1342:   #endif
1343:       } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot handle periodicity in this domain");
1344:     } else {
1345:       for (d = 0; d < dof; ++d, ++coordSize) ncoords[coordSize] = coords[off + d];
1346:     }
1347:   }
1348:   if (cutVertexLabel) {
1349:     IS              vertices;
1350:     const PetscInt *verts;
1351:     PetscInt        n;

1353:     PetscCall(DMLabelGetStratumIS(cutVertexLabel, 1, &vertices));
1354:     if (vertices) {
1355:       PetscCall(ISGetIndices(vertices, &verts));
1356:       PetscCall(ISGetLocalSize(vertices, &n));
1357:       for (v = 0; v < n; ++v) {
1358:         PetscCall(PetscSectionGetDof(cSection, verts[v], &dof));
1359:         PetscCall(PetscSectionGetOffset(cSection, verts[v], &off));
1360:         for (d = 0; d < dof; ++d) ncoords[coordSize++] = coords[off + d] + ((L[d] > 0.) ? L[d] : 0.0);
1361:       }
1362:       PetscCall(ISRestoreIndices(vertices, &verts));
1363:       PetscCall(ISDestroy(&vertices));
1364:     }
1365:   }
1366:   PetscCheck(coordSize == N, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Mismatched sizes: %" PetscInt_FMT " != %" PetscInt_FMT, coordSize, N);
1367:   PetscCall(DMLabelDestroy(&cutVertexLabel));
1368:   PetscCall(VecRestoreArray(coordinatesLocal, &coords));
1369:   PetscCall(VecRestoreArray(newcoords, &ncoords));
1370:   PetscCall(PetscObjectSetName((PetscObject)newcoords, "vertices"));
1371:   PetscCall(VecScale(newcoords, lengthScale));
1372:   PetscCall(PetscViewerHDF5PushGroup(viewer, "/viz"));
1373:   PetscCall(PetscViewerHDF5WriteGroup(viewer, NULL));
1374:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1375:   PetscCall(PetscViewerHDF5PushGroup(viewer, "/viz/geometry"));
1376:   PetscCall(VecView(newcoords, viewer));
1377:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1378:   PetscCall(VecDestroy(&newcoords));
1379:   PetscFunctionReturn(PETSC_SUCCESS);
1380: }

1382: PetscErrorCode DMPlexLabelsView_HDF5_Internal(DM dm, IS globalPointNumbers, PetscViewer viewer)
1383: {
1384:   const char          *topologydm_name;
1385:   const PetscInt      *gpoint;
1386:   PetscInt             numLabels, l;
1387:   DMPlexStorageVersion version;
1388:   char                 group[PETSC_MAX_PATH_LEN];

1390:   PetscFunctionBegin;
1391:   PetscCall(PetscViewerHDF5GetDMPlexStorageVersionWriting(viewer, &version));
1392:   PetscCall(ISGetIndices(globalPointNumbers, &gpoint));
1393:   PetscCall(DMPlexGetHDF5Name_Private(dm, &topologydm_name));
1394:   if (DMPlexStorageVersionGE(version, 2, 0, 0)) {
1395:     PetscCall(PetscSNPrintf(group, sizeof(group), "topologies/%s/labels", topologydm_name));
1396:   } else {
1397:     PetscCall(PetscStrncpy(group, "/labels", sizeof(group)));
1398:   }
1399:   PetscCall(PetscViewerHDF5PushGroup(viewer, group));
1400:   PetscCall(DMGetNumLabels(dm, &numLabels));
1401:   for (l = 0; l < numLabels; ++l) {
1402:     DMLabel         label;
1403:     const char     *name;
1404:     IS              valueIS, pvalueIS, globalValueIS;
1405:     const PetscInt *values;
1406:     PetscInt        numValues, v;
1407:     PetscBool       isDepth, output;

1409:     PetscCall(DMGetLabelByNum(dm, l, &label));
1410:     PetscCall(PetscObjectGetName((PetscObject)label, &name));
1411:     PetscCall(DMGetLabelOutput(dm, name, &output));
1412:     PetscCall(PetscStrncmp(name, "depth", 10, &isDepth));
1413:     if (isDepth || !output) continue;
1414:     PetscCall(PetscViewerHDF5PushGroup(viewer, name));
1415:     PetscCall(DMLabelGetValueIS(label, &valueIS));
1416:     /* Must copy to a new IS on the global comm */
1417:     PetscCall(ISGetLocalSize(valueIS, &numValues));
1418:     PetscCall(ISGetIndices(valueIS, &values));
1419:     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), numValues, values, PETSC_COPY_VALUES, &pvalueIS));
1420:     PetscCall(ISRestoreIndices(valueIS, &values));
1421:     PetscCall(ISAllGather(pvalueIS, &globalValueIS));
1422:     PetscCall(ISDestroy(&pvalueIS));
1423:     PetscCall(ISSortRemoveDups(globalValueIS));
1424:     PetscCall(ISGetLocalSize(globalValueIS, &numValues));
1425:     PetscCall(ISGetIndices(globalValueIS, &values));
1426:     for (v = 0; v < numValues; ++v) {
1427:       IS              stratumIS, globalStratumIS;
1428:       const PetscInt *spoints = NULL;
1429:       PetscInt       *gspoints, n = 0, gn, p;
1430:       const char     *iname = "indices";
1431:       char            group[PETSC_MAX_PATH_LEN];

1433:       PetscCall(PetscSNPrintf(group, sizeof(group), "%" PetscInt_FMT, values[v]));
1434:       PetscCall(PetscViewerHDF5PushGroup(viewer, group));
1435:       PetscCall(DMLabelGetStratumIS(label, values[v], &stratumIS));

1437:       if (stratumIS) PetscCall(ISGetLocalSize(stratumIS, &n));
1438:       if (stratumIS) PetscCall(ISGetIndices(stratumIS, &spoints));
1439:       for (gn = 0, p = 0; p < n; ++p)
1440:         if (gpoint[spoints[p]] >= 0) ++gn;
1441:       PetscCall(PetscMalloc1(gn, &gspoints));
1442:       for (gn = 0, p = 0; p < n; ++p)
1443:         if (gpoint[spoints[p]] >= 0) gspoints[gn++] = gpoint[spoints[p]];
1444:       if (stratumIS) PetscCall(ISRestoreIndices(stratumIS, &spoints));
1445:       PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), gn, gspoints, PETSC_OWN_POINTER, &globalStratumIS));
1446:       PetscCall(PetscObjectSetName((PetscObject)globalStratumIS, iname));

1448:       PetscCall(ISView(globalStratumIS, viewer));
1449:       PetscCall(ISDestroy(&globalStratumIS));
1450:       PetscCall(ISDestroy(&stratumIS));
1451:       PetscCall(PetscViewerHDF5PopGroup(viewer));
1452:     }
1453:     PetscCall(ISRestoreIndices(globalValueIS, &values));
1454:     PetscCall(ISDestroy(&globalValueIS));
1455:     PetscCall(ISDestroy(&valueIS));
1456:     PetscCall(PetscViewerHDF5PopGroup(viewer));
1457:   }
1458:   PetscCall(ISRestoreIndices(globalPointNumbers, &gpoint));
1459:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1460:   PetscFunctionReturn(PETSC_SUCCESS);
1461: }

1463: /* We only write cells and vertices. Does this screw up parallel reading? */
1464: PetscErrorCode DMPlexView_HDF5_Internal(DM dm, PetscViewer viewer)
1465: {
1466:   IS                globalPointNumbers;
1467:   PetscViewerFormat format;
1468:   PetscBool         viz_geom = PETSC_FALSE, xdmf_topo = PETSC_FALSE, petsc_topo = PETSC_FALSE;

1470:   PetscFunctionBegin;
1471:   PetscCall(DMPlexCreatePointNumbering(dm, &globalPointNumbers));
1472:   PetscCall(DMPlexCoordinatesView_HDF5_Internal(dm, viewer));

1474:   PetscCall(PetscViewerGetFormat(viewer, &format));
1475:   switch (format) {
1476:   case PETSC_VIEWER_HDF5_VIZ:
1477:     viz_geom  = PETSC_TRUE;
1478:     xdmf_topo = PETSC_TRUE;
1479:     break;
1480:   case PETSC_VIEWER_HDF5_XDMF:
1481:     xdmf_topo = PETSC_TRUE;
1482:     break;
1483:   case PETSC_VIEWER_HDF5_PETSC:
1484:     petsc_topo = PETSC_TRUE;
1485:     break;
1486:   case PETSC_VIEWER_DEFAULT:
1487:   case PETSC_VIEWER_NATIVE:
1488:     viz_geom   = PETSC_TRUE;
1489:     xdmf_topo  = PETSC_TRUE;
1490:     petsc_topo = PETSC_TRUE;
1491:     break;
1492:   default:
1493:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 output.", PetscViewerFormats[format]);
1494:   }

1496:   if (viz_geom) PetscCall(DMPlexCoordinatesView_HDF5_XDMF_Private(dm, viewer));
1497:   if (xdmf_topo) PetscCall(DMPlexTopologyView_HDF5_XDMF_Private(dm, globalPointNumbers, viewer));
1498:   if (petsc_topo) {
1499:     PetscCall(DMPlexTopologyView_HDF5_Internal(dm, globalPointNumbers, viewer));
1500:     PetscCall(DMPlexLabelsView_HDF5_Internal(dm, globalPointNumbers, viewer));
1501:   }

1503:   PetscCall(ISDestroy(&globalPointNumbers));
1504:   PetscFunctionReturn(PETSC_SUCCESS);
1505: }

1507: PetscErrorCode DMPlexSectionView_HDF5_Internal(DM dm, PetscViewer viewer, DM sectiondm)
1508: {
1509:   MPI_Comm     comm;
1510:   const char  *topologydm_name;
1511:   const char  *sectiondm_name;
1512:   PetscSection gsection;

1514:   PetscFunctionBegin;
1515:   PetscCall(PetscObjectGetComm((PetscObject)sectiondm, &comm));
1516:   PetscCall(DMPlexGetHDF5Name_Private(dm, &topologydm_name));
1517:   PetscCall(PetscObjectGetName((PetscObject)sectiondm, &sectiondm_name));
1518:   PetscCall(PetscViewerHDF5PushGroup(viewer, "topologies"));
1519:   PetscCall(PetscViewerHDF5PushGroup(viewer, topologydm_name));
1520:   PetscCall(PetscViewerHDF5PushGroup(viewer, "dms"));
1521:   PetscCall(PetscViewerHDF5PushGroup(viewer, sectiondm_name));
1522:   PetscCall(DMGetGlobalSection(sectiondm, &gsection));
1523:   /* Save raw section */
1524:   PetscCall(PetscSectionView(gsection, viewer));
1525:   /* Save plex wrapper */
1526:   {
1527:     PetscInt        pStart, pEnd, p, n;
1528:     IS              globalPointNumbers;
1529:     const PetscInt *gpoints;
1530:     IS              orderIS;
1531:     PetscInt       *order;

1533:     PetscCall(PetscSectionGetChart(gsection, &pStart, &pEnd));
1534:     PetscCall(DMPlexCreatePointNumbering(dm, &globalPointNumbers));
1535:     PetscCall(ISGetIndices(globalPointNumbers, &gpoints));
1536:     for (p = pStart, n = 0; p < pEnd; ++p)
1537:       if (gpoints[p] >= 0) n++;
1538:     /* "order" is an array of global point numbers.
1539:        When loading, it is used with topology/order array
1540:        to match section points with plex topology points. */
1541:     PetscCall(PetscMalloc1(n, &order));
1542:     for (p = pStart, n = 0; p < pEnd; ++p)
1543:       if (gpoints[p] >= 0) order[n++] = gpoints[p];
1544:     PetscCall(ISRestoreIndices(globalPointNumbers, &gpoints));
1545:     PetscCall(ISDestroy(&globalPointNumbers));
1546:     PetscCall(ISCreateGeneral(comm, n, order, PETSC_OWN_POINTER, &orderIS));
1547:     PetscCall(PetscObjectSetName((PetscObject)orderIS, "order"));
1548:     PetscCall(ISView(orderIS, viewer));
1549:     PetscCall(ISDestroy(&orderIS));
1550:   }
1551:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1552:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1553:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1554:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1555:   PetscFunctionReturn(PETSC_SUCCESS);
1556: }

1558: PetscErrorCode DMPlexGlobalVectorView_HDF5_Internal(DM dm, PetscViewer viewer, DM sectiondm, Vec vec)
1559: {
1560:   const char *topologydm_name;
1561:   const char *sectiondm_name;
1562:   const char *vec_name;
1563:   PetscInt    bs;

1565:   PetscFunctionBegin;
1566:   /* Check consistency */
1567:   {
1568:     PetscSF pointsf, pointsf1;

1570:     PetscCall(DMGetPointSF(dm, &pointsf));
1571:     PetscCall(DMGetPointSF(sectiondm, &pointsf1));
1572:     PetscCheck(pointsf1 == pointsf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Mismatching point SFs for dm and sectiondm");
1573:   }
1574:   PetscCall(DMPlexGetHDF5Name_Private(dm, &topologydm_name));
1575:   PetscCall(PetscObjectGetName((PetscObject)sectiondm, &sectiondm_name));
1576:   PetscCall(PetscObjectGetName((PetscObject)vec, &vec_name));
1577:   PetscCall(PetscViewerHDF5PushGroup(viewer, "topologies"));
1578:   PetscCall(PetscViewerHDF5PushGroup(viewer, topologydm_name));
1579:   PetscCall(PetscViewerHDF5PushGroup(viewer, "dms"));
1580:   PetscCall(PetscViewerHDF5PushGroup(viewer, sectiondm_name));
1581:   PetscCall(PetscViewerHDF5PushGroup(viewer, "vecs"));
1582:   PetscCall(PetscViewerHDF5PushGroup(viewer, vec_name));
1583:   PetscCall(VecGetBlockSize(vec, &bs));
1584:   PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "blockSize", PETSC_INT, (void *)&bs));
1585:   PetscCall(VecSetBlockSize(vec, 1));
1586:   /* VecView(vec, viewer) would call (*vec->opt->view)(vec, viewer), but,    */
1587:   /* if vec was created with DMGet{Global, Local}Vector(), vec->opt->view    */
1588:   /* is set to VecView_Plex, which would save vec in a predefined location.  */
1589:   /* To save vec in where we want, we create a new Vec (temp) with           */
1590:   /* VecCreate(), wrap the vec data in temp, and call VecView(temp, viewer). */
1591:   {
1592:     Vec                temp;
1593:     const PetscScalar *array;
1594:     PetscLayout        map;

1596:     PetscCall(VecCreate(PetscObjectComm((PetscObject)vec), &temp));
1597:     PetscCall(PetscObjectSetName((PetscObject)temp, vec_name));
1598:     PetscCall(VecGetLayout(vec, &map));
1599:     PetscCall(VecSetLayout(temp, map));
1600:     PetscCall(VecSetUp(temp));
1601:     PetscCall(VecGetArrayRead(vec, &array));
1602:     PetscCall(VecPlaceArray(temp, array));
1603:     PetscCall(VecView(temp, viewer));
1604:     PetscCall(VecResetArray(temp));
1605:     PetscCall(VecRestoreArrayRead(vec, &array));
1606:     PetscCall(VecDestroy(&temp));
1607:   }
1608:   PetscCall(VecSetBlockSize(vec, bs));
1609:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1610:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1611:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1612:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1613:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1614:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1615:   PetscFunctionReturn(PETSC_SUCCESS);
1616: }

1618: PetscErrorCode DMPlexLocalVectorView_HDF5_Internal(DM dm, PetscViewer viewer, DM sectiondm, Vec vec)
1619: {
1620:   MPI_Comm     comm;
1621:   const char  *topologydm_name;
1622:   const char  *sectiondm_name;
1623:   const char  *vec_name;
1624:   PetscSection section;
1625:   PetscBool    includesConstraints;
1626:   Vec          gvec;
1627:   PetscInt     m, bs;

1629:   PetscFunctionBegin;
1630:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1631:   /* Check consistency */
1632:   {
1633:     PetscSF pointsf, pointsf1;

1635:     PetscCall(DMGetPointSF(dm, &pointsf));
1636:     PetscCall(DMGetPointSF(sectiondm, &pointsf1));
1637:     PetscCheck(pointsf1 == pointsf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Mismatching point SFs for dm and sectiondm");
1638:   }
1639:   PetscCall(DMPlexGetHDF5Name_Private(dm, &topologydm_name));
1640:   PetscCall(PetscObjectGetName((PetscObject)sectiondm, &sectiondm_name));
1641:   PetscCall(PetscObjectGetName((PetscObject)vec, &vec_name));
1642:   PetscCall(PetscViewerHDF5PushGroup(viewer, "topologies"));
1643:   PetscCall(PetscViewerHDF5PushGroup(viewer, topologydm_name));
1644:   PetscCall(PetscViewerHDF5PushGroup(viewer, "dms"));
1645:   PetscCall(PetscViewerHDF5PushGroup(viewer, sectiondm_name));
1646:   PetscCall(PetscViewerHDF5PushGroup(viewer, "vecs"));
1647:   PetscCall(PetscViewerHDF5PushGroup(viewer, vec_name));
1648:   PetscCall(VecGetBlockSize(vec, &bs));
1649:   PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "blockSize", PETSC_INT, (void *)&bs));
1650:   PetscCall(VecCreate(comm, &gvec));
1651:   PetscCall(PetscObjectSetName((PetscObject)gvec, vec_name));
1652:   PetscCall(DMGetGlobalSection(sectiondm, &section));
1653:   PetscCall(PetscSectionGetIncludesConstraints(section, &includesConstraints));
1654:   if (includesConstraints) PetscCall(PetscSectionGetStorageSize(section, &m));
1655:   else PetscCall(PetscSectionGetConstrainedStorageSize(section, &m));
1656:   PetscCall(VecSetSizes(gvec, m, PETSC_DECIDE));
1657:   PetscCall(VecSetUp(gvec));
1658:   PetscCall(DMLocalToGlobalBegin(sectiondm, vec, INSERT_VALUES, gvec));
1659:   PetscCall(DMLocalToGlobalEnd(sectiondm, vec, INSERT_VALUES, gvec));
1660:   PetscCall(VecView(gvec, viewer));
1661:   PetscCall(VecDestroy(&gvec));
1662:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1663:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1664:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1665:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1666:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1667:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1668:   PetscFunctionReturn(PETSC_SUCCESS);
1669: }

1671: struct _n_LoadLabelsCtx {
1672:   MPI_Comm    comm;
1673:   PetscMPIInt rank;
1674:   DM          dm;
1675:   PetscViewer viewer;
1676:   DMLabel     label;
1677:   PetscSF     sfXC;
1678:   PetscLayout layoutX;
1679: };
1680: typedef struct _n_LoadLabelsCtx *LoadLabelsCtx;

1682: static PetscErrorCode LoadLabelsCtxCreate(DM dm, PetscViewer viewer, PetscSF sfXC, LoadLabelsCtx *ctx)
1683: {
1684:   PetscFunctionBegin;
1685:   PetscCall(PetscNew(ctx));
1686:   PetscCall(PetscObjectReference((PetscObject)((*ctx)->dm = dm)));
1687:   PetscCall(PetscObjectReference((PetscObject)((*ctx)->viewer = viewer)));
1688:   PetscCall(PetscObjectGetComm((PetscObject)dm, &(*ctx)->comm));
1689:   PetscCallMPI(MPI_Comm_rank((*ctx)->comm, &(*ctx)->rank));
1690:   (*ctx)->sfXC = sfXC;
1691:   if (sfXC) {
1692:     PetscInt nX;

1694:     PetscCall(PetscObjectReference((PetscObject)sfXC));
1695:     PetscCall(PetscSFGetGraph(sfXC, &nX, NULL, NULL, NULL));
1696:     PetscCall(PetscLayoutCreateFromSizes((*ctx)->comm, nX, PETSC_DECIDE, 1, &(*ctx)->layoutX));
1697:   }
1698:   PetscFunctionReturn(PETSC_SUCCESS);
1699: }

1701: static PetscErrorCode LoadLabelsCtxDestroy(LoadLabelsCtx *ctx)
1702: {
1703:   PetscFunctionBegin;
1704:   if (!*ctx) PetscFunctionReturn(PETSC_SUCCESS);
1705:   PetscCall(DMDestroy(&(*ctx)->dm));
1706:   PetscCall(PetscViewerDestroy(&(*ctx)->viewer));
1707:   PetscCall(PetscSFDestroy(&(*ctx)->sfXC));
1708:   PetscCall(PetscLayoutDestroy(&(*ctx)->layoutX));
1709:   PetscCall(PetscFree(*ctx));
1710:   PetscFunctionReturn(PETSC_SUCCESS);
1711: }

1713: /*
1714:     A: on-disk points
1715:     X: global points [0, NX)
1716:     C: distributed plex points
1717: */
1718: static herr_t ReadLabelStratumHDF5_Distribute_Private(IS stratumIS, LoadLabelsCtx ctx, IS *newStratumIS)
1719: {
1720:   MPI_Comm        comm    = ctx->comm;
1721:   PetscSF         sfXC    = ctx->sfXC;
1722:   PetscLayout     layoutX = ctx->layoutX;
1723:   PetscSF         sfXA;
1724:   const PetscInt *A_points;
1725:   PetscInt        nX, nC;
1726:   PetscInt        n;

1728:   PetscFunctionBegin;
1729:   PetscCall(PetscSFGetGraph(sfXC, &nX, &nC, NULL, NULL));
1730:   PetscCall(ISGetLocalSize(stratumIS, &n));
1731:   PetscCall(ISGetIndices(stratumIS, &A_points));
1732:   PetscCall(PetscSFCreate(comm, &sfXA));
1733:   PetscCall(PetscSFSetGraphLayout(sfXA, layoutX, n, NULL, PETSC_USE_POINTER, A_points));
1734:   PetscCall(ISCreate(comm, newStratumIS));
1735:   PetscCall(ISSetType(*newStratumIS, ISGENERAL));
1736:   {
1737:     PetscInt   i;
1738:     PetscBool *A_mask, *X_mask, *C_mask;

1740:     PetscCall(PetscCalloc3(n, &A_mask, nX, &X_mask, nC, &C_mask));
1741:     for (i = 0; i < n; i++) A_mask[i] = PETSC_TRUE;
1742:     PetscCall(PetscSFReduceBegin(sfXA, MPIU_BOOL, A_mask, X_mask, MPI_REPLACE));
1743:     PetscCall(PetscSFReduceEnd(sfXA, MPIU_BOOL, A_mask, X_mask, MPI_REPLACE));
1744:     PetscCall(PetscSFBcastBegin(sfXC, MPIU_BOOL, X_mask, C_mask, MPI_LOR));
1745:     PetscCall(PetscSFBcastEnd(sfXC, MPIU_BOOL, X_mask, C_mask, MPI_LOR));
1746:     PetscCall(ISGeneralSetIndicesFromMask(*newStratumIS, 0, nC, C_mask));
1747:     PetscCall(PetscFree3(A_mask, X_mask, C_mask));
1748:   }
1749:   PetscCall(PetscSFDestroy(&sfXA));
1750:   PetscCall(ISRestoreIndices(stratumIS, &A_points));
1751:   PetscFunctionReturn(PETSC_SUCCESS);
1752: }

1754: static herr_t ReadLabelStratumHDF5_Static(hid_t g_id, const char *vname, const H5L_info_t *info, void *op_data)
1755: {
1756:   LoadLabelsCtx   ctx    = (LoadLabelsCtx)op_data;
1757:   PetscViewer     viewer = ctx->viewer;
1758:   DMLabel         label  = ctx->label;
1759:   MPI_Comm        comm   = ctx->comm;
1760:   IS              stratumIS;
1761:   const PetscInt *ind;
1762:   PetscInt        value, N, i;

1764:   PetscCall(PetscOptionsStringToInt(vname, &value));
1765:   PetscCall(ISCreate(comm, &stratumIS));
1766:   PetscCall(PetscObjectSetName((PetscObject)stratumIS, "indices"));
1767:   PetscCall(PetscViewerHDF5PushGroup(viewer, vname)); /* labels/<lname>/<vname> */

1769:   if (!ctx->sfXC) {
1770:     /* Force serial load */
1771:     PetscCall(PetscViewerHDF5ReadSizes(viewer, "indices", NULL, &N));
1772:     PetscCall(PetscLayoutSetLocalSize(stratumIS->map, !ctx->rank ? N : 0));
1773:     PetscCall(PetscLayoutSetSize(stratumIS->map, N));
1774:   }
1775:   PetscCall(ISLoad(stratumIS, viewer));

1777:   if (ctx->sfXC) {
1778:     IS newStratumIS;

1780:     PetscCallHDF5(ReadLabelStratumHDF5_Distribute_Private, (stratumIS, ctx, &newStratumIS));
1781:     PetscCall(ISDestroy(&stratumIS));
1782:     stratumIS = newStratumIS;
1783:   }

1785:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1786:   PetscCall(ISGetLocalSize(stratumIS, &N));
1787:   PetscCall(ISGetIndices(stratumIS, &ind));
1788:   for (i = 0; i < N; ++i) PetscCall(DMLabelSetValue(label, ind[i], value));
1789:   PetscCall(ISRestoreIndices(stratumIS, &ind));
1790:   PetscCall(ISDestroy(&stratumIS));
1791:   return 0;
1792: }

1794: /* TODO: Fix this code, it is returning PETSc error codes when it should be translating them to herr_t codes */
1795: static herr_t ReadLabelHDF5_Static(hid_t g_id, const char *lname, const H5L_info_t *info, void *op_data)
1796: {
1797:   LoadLabelsCtx  ctx = (LoadLabelsCtx)op_data;
1798:   DM             dm  = ctx->dm;
1799:   hsize_t        idx = 0;
1800:   PetscErrorCode ierr;
1801:   PetscBool      flg;
1802:   herr_t         err;

1804:   PetscCall(DMHasLabel(dm, lname, &flg));
1805:   if (flg) PetscCall(DMRemoveLabel(dm, lname, NULL));
1806:   ierr = DMCreateLabel(dm, lname);
1807:   if (ierr) return (herr_t)ierr;
1808:   ierr = DMGetLabel(dm, lname, &ctx->label);
1809:   if (ierr) return (herr_t)ierr;
1810:   ierr = PetscViewerHDF5PushGroup(ctx->viewer, lname);
1811:   if (ierr) return (herr_t)ierr;
1812:   /* Iterate over the label's strata */
1813:   PetscCallHDF5Return(err, H5Literate_by_name, (g_id, lname, H5_INDEX_NAME, H5_ITER_NATIVE, &idx, ReadLabelStratumHDF5_Static, op_data, 0));
1814:   ierr = PetscViewerHDF5PopGroup(ctx->viewer);
1815:   if (ierr) return (herr_t)ierr;
1816:   return err;
1817: }

1819: PetscErrorCode DMPlexLabelsLoad_HDF5_Internal(DM dm, PetscViewer viewer, PetscSF sfXC)
1820: {
1821:   const char          *topologydm_name;
1822:   LoadLabelsCtx        ctx;
1823:   hsize_t              idx = 0;
1824:   char                 group[PETSC_MAX_PATH_LEN];
1825:   DMPlexStorageVersion version;
1826:   PetscBool            distributed, hasGroup;

1828:   PetscFunctionBegin;
1829:   PetscCall(DMPlexIsDistributed(dm, &distributed));
1830:   if (distributed) PetscCheck(sfXC, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_NULL, "PetscSF must be given for parallel load");
1831:   PetscCall(LoadLabelsCtxCreate(dm, viewer, sfXC, &ctx));
1832:   PetscCall(DMPlexGetHDF5Name_Private(dm, &topologydm_name));
1833:   PetscCall(PetscViewerHDF5GetDMPlexStorageVersionReading(viewer, &version));
1834:   if (DMPlexStorageVersionGE(version, 2, 0, 0)) {
1835:     PetscCall(PetscSNPrintf(group, sizeof(group), "topologies/%s/labels", topologydm_name));
1836:   } else {
1837:     PetscCall(PetscStrncpy(group, "labels", sizeof(group)));
1838:   }
1839:   PetscCall(PetscViewerHDF5PushGroup(viewer, group));
1840:   PetscCall(PetscViewerHDF5HasGroup(viewer, NULL, &hasGroup));
1841:   if (hasGroup) {
1842:     hid_t fileId, groupId;

1844:     PetscCall(PetscViewerHDF5OpenGroup(viewer, NULL, &fileId, &groupId));
1845:     /* Iterate over labels */
1846:     PetscCallHDF5(H5Literate, (groupId, H5_INDEX_NAME, H5_ITER_NATIVE, &idx, ReadLabelHDF5_Static, ctx));
1847:     PetscCallHDF5(H5Gclose, (groupId));
1848:   }
1849:   PetscCall(PetscViewerHDF5PopGroup(viewer));
1850:   PetscCall(LoadLabelsCtxDestroy(&ctx));
1851:   PetscFunctionReturn(PETSC_SUCCESS);
1852: }

1854: static PetscErrorCode DMPlexDistributionLoad_HDF5_Private(DM dm, PetscViewer viewer, PetscSF sf, PetscSF *distsf, DM *distdm)
1855: {
1856:   MPI_Comm        comm;
1857:   PetscMPIInt     size, rank;
1858:   PetscInt        dist_size;
1859:   const char     *distribution_name;
1860:   PetscInt        p, lsize;
1861:   IS              chartSizesIS, ownersIS, gpointsIS;
1862:   const PetscInt *chartSize, *owners, *gpoints;
1863:   PetscLayout     layout;
1864:   PetscBool       has;

1866:   PetscFunctionBegin;
1867:   *distsf = NULL;
1868:   *distdm = NULL;
1869:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1870:   PetscCallMPI(MPI_Comm_size(comm, &size));
1871:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1872:   PetscCall(DMPlexDistributionGetName(dm, &distribution_name));
1873:   if (!distribution_name) PetscFunctionReturn(PETSC_SUCCESS);
1874:   PetscCall(PetscLogEventBegin(DMPLEX_DistributionLoad, viewer, 0, 0, 0));
1875:   PetscCall(PetscViewerHDF5HasGroup(viewer, NULL, &has));
1876:   if (!has) {
1877:     char *full_group;

1879:     PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &full_group));
1880:     PetscCheck(has, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Distribution %s cannot be found: HDF5 group %s not found in file", distribution_name, full_group);
1881:   }
1882:   PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "comm_size", PETSC_INT, NULL, (void *)&dist_size));
1883:   PetscCheck(dist_size == (PetscInt)size, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mismatching comm sizes: comm size of this session (%d) != comm size used for %s (%" PetscInt_FMT ")", size, distribution_name, dist_size);
1884:   PetscCall(ISCreate(comm, &chartSizesIS));
1885:   PetscCall(PetscObjectSetName((PetscObject)chartSizesIS, "chart_sizes"));
1886:   PetscCall(ISCreate(comm, &ownersIS));
1887:   PetscCall(PetscObjectSetName((PetscObject)ownersIS, "owners"));
1888:   PetscCall(ISCreate(comm, &gpointsIS));
1889:   PetscCall(PetscObjectSetName((PetscObject)gpointsIS, "global_point_numbers"));
1890:   PetscCall(PetscLayoutSetLocalSize(chartSizesIS->map, 1));
1891:   PetscCall(ISLoad(chartSizesIS, viewer));
1892:   PetscCall(ISGetIndices(chartSizesIS, &chartSize));
1893:   PetscCall(PetscLayoutSetLocalSize(ownersIS->map, *chartSize));
1894:   PetscCall(PetscLayoutSetLocalSize(gpointsIS->map, *chartSize));
1895:   PetscCall(ISLoad(ownersIS, viewer));
1896:   PetscCall(ISLoad(gpointsIS, viewer));
1897:   PetscCall(ISGetIndices(ownersIS, &owners));
1898:   PetscCall(ISGetIndices(gpointsIS, &gpoints));
1899:   PetscCall(PetscSFCreate(comm, distsf));
1900:   PetscCall(PetscSFSetFromOptions(*distsf));
1901:   PetscCall(PetscLayoutCreate(comm, &layout));
1902:   PetscCall(PetscSFGetGraph(sf, &lsize, NULL, NULL, NULL));
1903:   PetscCall(PetscLayoutSetLocalSize(layout, lsize));
1904:   PetscCall(PetscLayoutSetBlockSize(layout, 1));
1905:   PetscCall(PetscLayoutSetUp(layout));
1906:   PetscCall(PetscSFSetGraphLayout(*distsf, layout, *chartSize, NULL, PETSC_OWN_POINTER, gpoints));
1907:   PetscCall(PetscLayoutDestroy(&layout));
1908:   /* Migrate DM */
1909:   {
1910:     PetscInt     pStart, pEnd;
1911:     PetscSFNode *buffer0, *buffer1, *buffer2;

1913:     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1914:     PetscCall(PetscMalloc2(pEnd - pStart, &buffer0, lsize, &buffer1));
1915:     PetscCall(PetscMalloc1(*chartSize, &buffer2));
1916:     {
1917:       PetscSF            workPointSF;
1918:       PetscInt           workNroots, workNleaves;
1919:       const PetscInt    *workIlocal;
1920:       const PetscSFNode *workIremote;

1922:       for (p = pStart; p < pEnd; ++p) {
1923:         buffer0[p - pStart].rank  = rank;
1924:         buffer0[p - pStart].index = p - pStart;
1925:       }
1926:       PetscCall(DMGetPointSF(dm, &workPointSF));
1927:       PetscCall(PetscSFGetGraph(workPointSF, &workNroots, &workNleaves, &workIlocal, &workIremote));
1928:       for (p = 0; p < workNleaves; ++p) {
1929:         PetscInt workIlocalp = (workIlocal ? workIlocal[p] : p);

1931:         buffer0[workIlocalp].rank = -1;
1932:       }
1933:     }
1934:     for (p = 0; p < lsize; ++p) buffer1[p].rank = -1;
1935:     for (p = 0; p < *chartSize; ++p) buffer2[p].rank = -1;
1936:     PetscCall(PetscSFReduceBegin(sf, MPIU_SF_NODE, buffer0, buffer1, MPI_MAXLOC));
1937:     PetscCall(PetscSFReduceEnd(sf, MPIU_SF_NODE, buffer0, buffer1, MPI_MAXLOC));
1938:     PetscCall(PetscSFBcastBegin(*distsf, MPIU_SF_NODE, buffer1, buffer2, MPI_REPLACE));
1939:     PetscCall(PetscSFBcastEnd(*distsf, MPIU_SF_NODE, buffer1, buffer2, MPI_REPLACE));
1940:     if (PetscDefined(USE_DEBUG)) {
1941:       for (p = 0; p < *chartSize; ++p) {
1942:         PetscCheck(buffer2[p].rank >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Found negative root rank %" PetscInt_FMT " at local point %" PetscInt_FMT " on rank %d when making migrationSF", buffer2[p].rank, p, rank);
1943:       }
1944:     }
1945:     PetscCall(PetscFree2(buffer0, buffer1));
1946:     PetscCall(DMCreate(comm, distdm));
1947:     PetscCall(DMSetType(*distdm, DMPLEX));
1948:     PetscCall(PetscObjectSetName((PetscObject)*distdm, ((PetscObject)dm)->name));
1949:     PetscCall(DMPlexDistributionSetName(*distdm, distribution_name));
1950:     {
1951:       PetscSF migrationSF;

1953:       PetscCall(PetscSFCreate(comm, &migrationSF));
1954:       PetscCall(PetscSFSetFromOptions(migrationSF));
1955:       PetscCall(PetscSFSetGraph(migrationSF, pEnd - pStart, *chartSize, NULL, PETSC_OWN_POINTER, buffer2, PETSC_OWN_POINTER));
1956:       PetscCall(PetscSFSetUp(migrationSF));
1957:       PetscCall(DMPlexMigrate(dm, migrationSF, *distdm));
1958:       PetscCall(PetscSFDestroy(&migrationSF));
1959:     }
1960:   }
1961:   /* Set pointSF */
1962:   {
1963:     PetscSF      pointSF;
1964:     PetscInt    *ilocal, nleaves, q;
1965:     PetscSFNode *iremote, *buffer0, *buffer1;

1967:     PetscCall(PetscMalloc2(*chartSize, &buffer0, lsize, &buffer1));
1968:     for (p = 0, nleaves = 0; p < *chartSize; ++p) {
1969:       if (owners[p] == rank) {
1970:         buffer0[p].rank = rank;
1971:       } else {
1972:         buffer0[p].rank = -1;
1973:         nleaves++;
1974:       }
1975:       buffer0[p].index = p;
1976:     }
1977:     for (p = 0; p < lsize; ++p) buffer1[p].rank = -1;
1978:     PetscCall(PetscSFReduceBegin(*distsf, MPIU_SF_NODE, buffer0, buffer1, MPI_MAXLOC));
1979:     PetscCall(PetscSFReduceEnd(*distsf, MPIU_SF_NODE, buffer0, buffer1, MPI_MAXLOC));
1980:     for (p = 0; p < *chartSize; ++p) buffer0[p].rank = -1;
1981:     PetscCall(PetscSFBcastBegin(*distsf, MPIU_SF_NODE, buffer1, buffer0, MPI_REPLACE));
1982:     PetscCall(PetscSFBcastEnd(*distsf, MPIU_SF_NODE, buffer1, buffer0, MPI_REPLACE));
1983:     if (PetscDefined(USE_DEBUG)) {
1984:       for (p = 0; p < *chartSize; ++p) PetscCheck(buffer0[p].rank >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Found negative root rank %" PetscInt_FMT " at local point %" PetscInt_FMT " on rank %d when making pointSF", buffer0[p].rank, p, rank);
1985:     }
1986:     PetscCall(PetscMalloc1(nleaves, &ilocal));
1987:     PetscCall(PetscMalloc1(nleaves, &iremote));
1988:     for (p = 0, q = 0; p < *chartSize; ++p) {
1989:       if (buffer0[p].rank != rank) {
1990:         ilocal[q]        = p;
1991:         iremote[q].rank  = buffer0[p].rank;
1992:         iremote[q].index = buffer0[p].index;
1993:         q++;
1994:       }
1995:     }
1996:     PetscCheck(q == nleaves, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Mismatching leaf sizes: %" PetscInt_FMT " != %" PetscInt_FMT, q, nleaves);
1997:     PetscCall(PetscFree2(buffer0, buffer1));
1998:     PetscCall(PetscSFCreate(comm, &pointSF));
1999:     PetscCall(PetscSFSetFromOptions(pointSF));
2000:     PetscCall(PetscSFSetGraph(pointSF, *chartSize, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
2001:     PetscCall(DMSetPointSF(*distdm, pointSF));
2002:     {
2003:       DM cdm;

2005:       PetscCall(DMGetCoordinateDM(*distdm, &cdm));
2006:       PetscCall(DMSetPointSF(cdm, pointSF));
2007:     }
2008:     PetscCall(PetscSFDestroy(&pointSF));
2009:   }
2010:   PetscCall(ISRestoreIndices(chartSizesIS, &chartSize));
2011:   PetscCall(ISRestoreIndices(ownersIS, &owners));
2012:   PetscCall(ISRestoreIndices(gpointsIS, &gpoints));
2013:   PetscCall(ISDestroy(&chartSizesIS));
2014:   PetscCall(ISDestroy(&ownersIS));
2015:   PetscCall(ISDestroy(&gpointsIS));
2016:   /* Record that overlap has been manually created.               */
2017:   /* This is to pass `DMPlexCheckPointSF()`, which checks that    */
2018:   /* pointSF does not contain cells in the leaves if overlap = 0. */
2019:   PetscCall(DMPlexSetOverlap_Plex(*distdm, NULL, DMPLEX_OVERLAP_MANUAL));
2020:   PetscCall(DMPlexDistributeSetDefault(*distdm, PETSC_FALSE));
2021:   PetscCall(DMPlexReorderSetDefault(*distdm, DM_REORDER_DEFAULT_FALSE));
2022:   PetscCall(PetscLogEventEnd(DMPLEX_DistributionLoad, viewer, 0, 0, 0));
2023:   PetscFunctionReturn(PETSC_SUCCESS);
2024: }

2026: // Serial load of topology
2027: static PetscErrorCode DMPlexTopologyLoad_HDF5_Legacy_Private(DM dm, PetscViewer viewer, PetscSF *sf)
2028: {
2029:   MPI_Comm        comm;
2030:   const char     *pointsName, *coneSizesName, *conesName, *orientationsName;
2031:   IS              pointsIS, coneSizesIS, conesIS, orientationsIS;
2032:   const PetscInt *points, *coneSizes, *cones, *orientations;
2033:   PetscInt       *cone, *ornt;
2034:   PetscInt        dim, N, Np, pEnd, p, q, maxConeSize = 0, c;
2035:   PetscMPIInt     size, rank;

2037:   PetscFunctionBegin;
2038:   pointsName       = "order";
2039:   coneSizesName    = "cones";
2040:   conesName        = "cells";
2041:   orientationsName = "orientation";
2042:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2043:   PetscCallMPI(MPI_Comm_size(comm, &size));
2044:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2045:   PetscCall(ISCreate(comm, &pointsIS));
2046:   PetscCall(PetscObjectSetName((PetscObject)pointsIS, pointsName));
2047:   PetscCall(ISCreate(comm, &coneSizesIS));
2048:   PetscCall(PetscObjectSetName((PetscObject)coneSizesIS, coneSizesName));
2049:   PetscCall(ISCreate(comm, &conesIS));
2050:   PetscCall(PetscObjectSetName((PetscObject)conesIS, conesName));
2051:   PetscCall(ISCreate(comm, &orientationsIS));
2052:   PetscCall(PetscObjectSetName((PetscObject)orientationsIS, orientationsName));
2053:   PetscCall(PetscViewerHDF5ReadObjectAttribute(viewer, (PetscObject)conesIS, "cell_dim", PETSC_INT, NULL, &dim));
2054:   PetscCall(DMSetDimension(dm, dim));
2055:   {
2056:     /* Force serial load */
2057:     PetscCall(PetscInfo(dm, "Loading DM %s in serial\n", dm->hdr.name));
2058:     PetscCall(PetscViewerHDF5ReadSizes(viewer, pointsName, NULL, &Np));
2059:     PetscCall(PetscLayoutSetLocalSize(pointsIS->map, rank == 0 ? Np : 0));
2060:     PetscCall(PetscLayoutSetSize(pointsIS->map, Np));
2061:     pEnd = rank == 0 ? Np : 0;
2062:     PetscCall(PetscViewerHDF5ReadSizes(viewer, coneSizesName, NULL, &Np));
2063:     PetscCall(PetscLayoutSetLocalSize(coneSizesIS->map, rank == 0 ? Np : 0));
2064:     PetscCall(PetscLayoutSetSize(coneSizesIS->map, Np));
2065:     PetscCall(PetscViewerHDF5ReadSizes(viewer, conesName, NULL, &N));
2066:     PetscCall(PetscLayoutSetLocalSize(conesIS->map, rank == 0 ? N : 0));
2067:     PetscCall(PetscLayoutSetSize(conesIS->map, N));
2068:     PetscCall(PetscViewerHDF5ReadSizes(viewer, orientationsName, NULL, &N));
2069:     PetscCall(PetscLayoutSetLocalSize(orientationsIS->map, rank == 0 ? N : 0));
2070:     PetscCall(PetscLayoutSetSize(orientationsIS->map, N));
2071:   }
2072:   PetscCall(ISLoad(pointsIS, viewer));
2073:   PetscCall(ISLoad(coneSizesIS, viewer));
2074:   PetscCall(ISLoad(conesIS, viewer));
2075:   PetscCall(ISLoad(orientationsIS, viewer));
2076:   /* Create Plex */
2077:   PetscCall(DMPlexSetChart(dm, 0, pEnd));
2078:   PetscCall(ISGetIndices(pointsIS, &points));
2079:   PetscCall(ISGetIndices(coneSizesIS, &coneSizes));
2080:   for (p = 0; p < pEnd; ++p) {
2081:     PetscCall(DMPlexSetConeSize(dm, points[p], coneSizes[p]));
2082:     maxConeSize = PetscMax(maxConeSize, coneSizes[p]);
2083:   }
2084:   PetscCall(DMSetUp(dm));
2085:   PetscCall(ISGetIndices(conesIS, &cones));
2086:   PetscCall(ISGetIndices(orientationsIS, &orientations));
2087:   PetscCall(PetscMalloc2(maxConeSize, &cone, maxConeSize, &ornt));
2088:   for (p = 0, q = 0; p < pEnd; ++p) {
2089:     for (c = 0; c < coneSizes[p]; ++c, ++q) {
2090:       cone[c] = cones[q];
2091:       ornt[c] = orientations[q];
2092:     }
2093:     PetscCall(DMPlexSetCone(dm, points[p], cone));
2094:     PetscCall(DMPlexSetConeOrientation(dm, points[p], ornt));
2095:   }
2096:   PetscCall(PetscFree2(cone, ornt));
2097:   /* Create global section migration SF */
2098:   if (sf) {
2099:     PetscLayout layout;
2100:     PetscInt   *globalIndices;

2102:     PetscCall(PetscMalloc1(pEnd, &globalIndices));
2103:     /* plex point == globalPointNumber in this case */
2104:     for (p = 0; p < pEnd; ++p) globalIndices[p] = p;
2105:     PetscCall(PetscLayoutCreate(comm, &layout));
2106:     PetscCall(PetscLayoutSetSize(layout, Np));
2107:     PetscCall(PetscLayoutSetBlockSize(layout, 1));
2108:     PetscCall(PetscLayoutSetUp(layout));
2109:     PetscCall(PetscSFCreate(comm, sf));
2110:     PetscCall(PetscSFSetFromOptions(*sf));
2111:     PetscCall(PetscSFSetGraphLayout(*sf, layout, pEnd, NULL, PETSC_OWN_POINTER, globalIndices));
2112:     PetscCall(PetscLayoutDestroy(&layout));
2113:     PetscCall(PetscFree(globalIndices));
2114:   }
2115:   /* Clean-up */
2116:   PetscCall(ISRestoreIndices(pointsIS, &points));
2117:   PetscCall(ISRestoreIndices(coneSizesIS, &coneSizes));
2118:   PetscCall(ISRestoreIndices(conesIS, &cones));
2119:   PetscCall(ISRestoreIndices(orientationsIS, &orientations));
2120:   PetscCall(ISDestroy(&pointsIS));
2121:   PetscCall(ISDestroy(&coneSizesIS));
2122:   PetscCall(ISDestroy(&conesIS));
2123:   PetscCall(ISDestroy(&orientationsIS));
2124:   /* Fill in the rest of the topology structure */
2125:   PetscCall(DMPlexSymmetrize(dm));
2126:   PetscCall(DMPlexStratify(dm));
2127:   PetscFunctionReturn(PETSC_SUCCESS);
2128: }

2130: /* Representation of two DMPlex strata in 0-based global numbering */
2131: struct _n_PlexLayer {
2132:   PetscInt     d;
2133:   IS           conesIS, orientationsIS;
2134:   PetscSection coneSizesSection;
2135:   PetscLayout  vertexLayout;
2136:   PetscSF      overlapSF, l2gSF; //TODO maybe confusing names (in DMPlex in general)
2137:   PetscInt     offset, conesOffset, leafOffset;
2138: };
2139: typedef struct _n_PlexLayer *PlexLayer;

2141: static PetscErrorCode PlexLayerDestroy(PlexLayer *layer)
2142: {
2143:   PetscFunctionBegin;
2144:   if (!*layer) PetscFunctionReturn(PETSC_SUCCESS);
2145:   PetscCall(PetscSectionDestroy(&(*layer)->coneSizesSection));
2146:   PetscCall(ISDestroy(&(*layer)->conesIS));
2147:   PetscCall(ISDestroy(&(*layer)->orientationsIS));
2148:   PetscCall(PetscSFDestroy(&(*layer)->overlapSF));
2149:   PetscCall(PetscSFDestroy(&(*layer)->l2gSF));
2150:   PetscCall(PetscLayoutDestroy(&(*layer)->vertexLayout));
2151:   PetscCall(PetscFree(*layer));
2152:   PetscFunctionReturn(PETSC_SUCCESS);
2153: }

2155: static PetscErrorCode PlexLayerCreate_Private(PlexLayer *layer)
2156: {
2157:   PetscFunctionBegin;
2158:   PetscCall(PetscNew(layer));
2159:   (*layer)->d           = -1;
2160:   (*layer)->offset      = -1;
2161:   (*layer)->conesOffset = -1;
2162:   (*layer)->leafOffset  = -1;
2163:   PetscFunctionReturn(PETSC_SUCCESS);
2164: }

2166: // Parallel load of a depth stratum
2167: static PetscErrorCode PlexLayerLoad_Private(PlexLayer layer, PetscViewer viewer, PetscInt d, PetscLayout pointsLayout)
2168: {
2169:   char         path[128];
2170:   MPI_Comm     comm;
2171:   const char  *coneSizesName, *conesName, *orientationsName;
2172:   IS           coneSizesIS, conesIS, orientationsIS;
2173:   PetscSection coneSizesSection;
2174:   PetscLayout  vertexLayout = NULL;
2175:   PetscInt     s;

2177:   PetscFunctionBegin;
2178:   coneSizesName    = "cone_sizes";
2179:   conesName        = "cones";
2180:   orientationsName = "orientations";
2181:   PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));

2183:   /* query size of next lower depth stratum (next lower dimension) */
2184:   if (d > 0) {
2185:     PetscInt NVertices;

2187:     PetscCall(PetscSNPrintf(path, sizeof(path), "%" PetscInt_FMT "/%s", d - 1, coneSizesName));
2188:     PetscCall(PetscViewerHDF5ReadSizes(viewer, path, NULL, &NVertices));
2189:     PetscCall(PetscLayoutCreate(comm, &vertexLayout));
2190:     PetscCall(PetscLayoutSetSize(vertexLayout, NVertices));
2191:     PetscCall(PetscLayoutSetUp(vertexLayout));
2192:   }

2194:   PetscCall(PetscSNPrintf(path, sizeof(path), "%" PetscInt_FMT, d));
2195:   PetscCall(PetscViewerHDF5PushGroup(viewer, path));

2197:   /* create coneSizesSection from stored IS coneSizes */
2198:   {
2199:     const PetscInt *coneSizes;

2201:     PetscCall(ISCreate(comm, &coneSizesIS));
2202:     PetscCall(PetscObjectSetName((PetscObject)coneSizesIS, coneSizesName));
2203:     if (pointsLayout) PetscCall(ISSetLayout(coneSizesIS, pointsLayout));
2204:     PetscCall(ISLoad(coneSizesIS, viewer));
2205:     if (!pointsLayout) PetscCall(ISGetLayout(coneSizesIS, &pointsLayout));
2206:     PetscCall(ISGetIndices(coneSizesIS, &coneSizes));
2207:     PetscCall(PetscSectionCreate(comm, &coneSizesSection));
2208:     //TODO different start ?
2209:     PetscCall(PetscSectionSetChart(coneSizesSection, 0, pointsLayout->n));
2210:     for (s = 0; s < pointsLayout->n; ++s) PetscCall(PetscSectionSetDof(coneSizesSection, s, coneSizes[s]));
2211:     PetscCall(PetscSectionSetUp(coneSizesSection));
2212:     PetscCall(ISRestoreIndices(coneSizesIS, &coneSizes));
2213:     {
2214:       PetscLayout tmp = NULL;

2216:       /* We need to keep the layout until the end of function */
2217:       PetscCall(PetscLayoutReference((PetscLayout)pointsLayout, &tmp));
2218:     }
2219:     PetscCall(ISDestroy(&coneSizesIS));
2220:   }

2222:   /* use value layout of coneSizesSection as layout of cones and orientations */
2223:   {
2224:     PetscLayout conesLayout;

2226:     PetscCall(PetscSectionGetValueLayout(comm, coneSizesSection, &conesLayout));
2227:     PetscCall(ISCreate(comm, &conesIS));
2228:     PetscCall(ISCreate(comm, &orientationsIS));
2229:     PetscCall(PetscObjectSetName((PetscObject)conesIS, conesName));
2230:     PetscCall(PetscObjectSetName((PetscObject)orientationsIS, orientationsName));
2231:     PetscCall(PetscLayoutDuplicate(conesLayout, &conesIS->map));
2232:     PetscCall(PetscLayoutDuplicate(conesLayout, &orientationsIS->map));
2233:     PetscCall(ISLoad(conesIS, viewer));
2234:     PetscCall(ISLoad(orientationsIS, viewer));
2235:     PetscCall(PetscLayoutDestroy(&conesLayout));
2236:   }

2238:   /* check assertion that layout of points is the same as point layout of coneSizesSection */
2239:   {
2240:     PetscLayout pointsLayout0;
2241:     PetscBool   flg;

2243:     PetscCall(PetscSectionGetPointLayout(comm, coneSizesSection, &pointsLayout0));
2244:     PetscCall(PetscLayoutCompare(pointsLayout, pointsLayout0, &flg));
2245:     PetscCheck(flg, comm, PETSC_ERR_PLIB, "points layout != coneSizesSection point layout");
2246:     PetscCall(PetscLayoutDestroy(&pointsLayout0));
2247:   }
2248:   PetscCall(PetscViewerHDF5PopGroup(viewer));
2249:   PetscCall(PetscLayoutDestroy(&pointsLayout));

2251:   layer->d                = d;
2252:   layer->conesIS          = conesIS;
2253:   layer->coneSizesSection = coneSizesSection;
2254:   layer->orientationsIS   = orientationsIS;
2255:   layer->vertexLayout     = vertexLayout;
2256:   PetscFunctionReturn(PETSC_SUCCESS);
2257: }

2259: static PetscErrorCode PlexLayerDistribute_Private(PlexLayer layer, PetscSF cellLocalToGlobalSF)
2260: {
2261:   IS           newConesIS, newOrientationsIS;
2262:   PetscSection newConeSizesSection;
2263:   MPI_Comm     comm;

2265:   PetscFunctionBegin;
2266:   PetscCall(PetscObjectGetComm((PetscObject)cellLocalToGlobalSF, &comm));
2267:   PetscCall(PetscSectionCreate(comm, &newConeSizesSection));
2268:   //TODO rename to something like ISDistribute() with optional PetscSection argument
2269:   PetscCall(DMPlexDistributeFieldIS(NULL, cellLocalToGlobalSF, layer->coneSizesSection, layer->conesIS, newConeSizesSection, &newConesIS));
2270:   PetscCall(DMPlexDistributeFieldIS(NULL, cellLocalToGlobalSF, layer->coneSizesSection, layer->orientationsIS, newConeSizesSection, &newOrientationsIS));

2272:   PetscCall(PetscObjectSetName((PetscObject)newConeSizesSection, ((PetscObject)layer->coneSizesSection)->name));
2273:   PetscCall(PetscObjectSetName((PetscObject)newConesIS, ((PetscObject)layer->conesIS)->name));
2274:   PetscCall(PetscObjectSetName((PetscObject)newOrientationsIS, ((PetscObject)layer->orientationsIS)->name));
2275:   PetscCall(PetscSectionDestroy(&layer->coneSizesSection));
2276:   PetscCall(ISDestroy(&layer->conesIS));
2277:   PetscCall(ISDestroy(&layer->orientationsIS));
2278:   layer->coneSizesSection = newConeSizesSection;
2279:   layer->conesIS          = newConesIS;
2280:   layer->orientationsIS   = newOrientationsIS;
2281:   PetscFunctionReturn(PETSC_SUCCESS);
2282: }

2284:   //TODO share code with DMPlexBuildFromCellListParallel()
2285: #include <petsc/private/hashseti.h>
2286: static PetscErrorCode PlexLayerCreateSFs_Private(PlexLayer layer, PetscSF *vertexOverlapSF, PetscSF *sfXC)
2287: {
2288:   PetscLayout     vertexLayout     = layer->vertexLayout;
2289:   PetscSection    coneSection      = layer->coneSizesSection;
2290:   IS              cellVertexData   = layer->conesIS;
2291:   IS              coneOrientations = layer->orientationsIS;
2292:   PetscSF         vl2gSF, vOverlapSF;
2293:   PetscInt       *verticesAdj;
2294:   PetscInt        i, n, numVerticesAdj;
2295:   const PetscInt *cvd, *co = NULL;
2296:   MPI_Comm        comm;

2298:   PetscFunctionBegin;
2299:   PetscCall(PetscObjectGetComm((PetscObject)coneSection, &comm));
2300:   PetscCall(PetscSectionGetStorageSize(coneSection, &n));
2301:   {
2302:     PetscInt n0;

2304:     PetscCall(ISGetLocalSize(cellVertexData, &n0));
2305:     PetscCheck(n == n0, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local size of IS cellVertexData = %" PetscInt_FMT " != %" PetscInt_FMT " = storage size of PetscSection coneSection", n0, n);
2306:     PetscCall(ISGetIndices(cellVertexData, &cvd));
2307:   }
2308:   if (coneOrientations) {
2309:     PetscInt n0;

2311:     PetscCall(ISGetLocalSize(coneOrientations, &n0));
2312:     PetscCheck(n == n0, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local size of IS coneOrientations = %" PetscInt_FMT " != %" PetscInt_FMT " = storage size of PetscSection coneSection", n0, n);
2313:     PetscCall(ISGetIndices(coneOrientations, &co));
2314:   }
2315:   /* Get/check global number of vertices */
2316:   {
2317:     PetscInt NVerticesInCells = PETSC_INT_MIN;

2319:     /* NVerticesInCells = max(cellVertexData) + 1 */
2320:     for (i = 0; i < n; i++)
2321:       if (cvd[i] > NVerticesInCells) NVerticesInCells = cvd[i];
2322:     ++NVerticesInCells;
2323:     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &NVerticesInCells, 1, MPIU_INT, MPI_MAX, comm));

2325:     if (vertexLayout->n == PETSC_DECIDE && vertexLayout->N == PETSC_DECIDE) vertexLayout->N = NVerticesInCells;
2326:     else
2327:       PetscCheck(vertexLayout->N == PETSC_DECIDE || vertexLayout->N >= NVerticesInCells, comm, PETSC_ERR_ARG_SIZ, "Specified global number of vertices %" PetscInt_FMT " must be greater than or equal to the number of unique vertices in the cell-vertex dataset %" PetscInt_FMT,
2328:                  vertexLayout->N, NVerticesInCells);
2329:     PetscCall(PetscLayoutSetUp(vertexLayout));
2330:   }
2331:   /* Find locally unique vertices in cellVertexData */
2332:   {
2333:     PetscHSetI vhash;
2334:     PetscInt   off = 0;

2336:     PetscCall(PetscHSetICreate(&vhash));
2337:     for (i = 0; i < n; ++i) PetscCall(PetscHSetIAdd(vhash, cvd[i]));
2338:     PetscCall(PetscHSetIGetSize(vhash, &numVerticesAdj));
2339:     PetscCall(PetscMalloc1(numVerticesAdj, &verticesAdj));
2340:     PetscCall(PetscHSetIGetElems(vhash, &off, verticesAdj));
2341:     PetscAssert(off == numVerticesAdj, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Assertion failed: off == numVerticesAdj (%" PetscInt_FMT " != %" PetscInt_FMT ")", off, numVerticesAdj);
2342:     PetscCall(PetscHSetIDestroy(&vhash));
2343:   }
2344:   /* We must sort vertices to preserve numbering */
2345:   PetscCall(PetscSortInt(numVerticesAdj, verticesAdj));
2346:   /* Connect locally unique vertices with star forest */
2347:   PetscCall(PetscSFCreateByMatchingIndices(vertexLayout, numVerticesAdj, verticesAdj, NULL, 0, numVerticesAdj, verticesAdj, NULL, 0, &vl2gSF, &vOverlapSF));
2348:   PetscCall(PetscObjectSetName((PetscObject)vOverlapSF, "overlapSF"));
2349:   PetscCall(PetscObjectSetName((PetscObject)vl2gSF, "localToGlobalSF"));

2351:   PetscCall(PetscFree(verticesAdj));
2352:   *vertexOverlapSF = vOverlapSF;
2353:   *sfXC            = vl2gSF;
2354:   PetscFunctionReturn(PETSC_SUCCESS);
2355: }

2357: static PetscErrorCode PlexLayerCreateCellSFs_Private(PlexLayer layer, PetscSF *cellOverlapSF, PetscSF *cellLocalToGlobalSF)
2358: {
2359:   PetscSection coneSection = layer->coneSizesSection;
2360:   PetscInt     nCells;
2361:   MPI_Comm     comm;

2363:   PetscFunctionBegin;
2364:   PetscCall(PetscObjectGetComm((PetscObject)coneSection, &comm));
2365:   {
2366:     PetscInt cStart;

2368:     PetscCall(PetscSectionGetChart(coneSection, &cStart, &nCells));
2369:     PetscCheck(cStart == 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "coneSection must start at 0");
2370:   }
2371:   /* Create overlapSF as empty SF with the right number of roots */
2372:   PetscCall(PetscSFCreate(comm, cellOverlapSF));
2373:   PetscCall(PetscSFSetGraph(*cellOverlapSF, nCells, 0, NULL, PETSC_USE_POINTER, NULL, PETSC_USE_POINTER));
2374:   PetscCall(PetscSFSetUp(*cellOverlapSF));
2375:   /* Create localToGlobalSF as identity mapping */
2376:   {
2377:     PetscLayout map;

2379:     PetscCall(PetscLayoutCreateFromSizes(comm, nCells, PETSC_DECIDE, 1, &map));
2380:     PetscCall(PetscSFCreateFromLayouts(map, map, cellLocalToGlobalSF));
2381:     PetscCall(PetscSFSetUp(*cellLocalToGlobalSF));
2382:     PetscCall(PetscLayoutDestroy(&map));
2383:   }
2384:   PetscFunctionReturn(PETSC_SUCCESS);
2385: }

2387: static PetscErrorCode DMPlexTopologyBuildFromLayers_Private(DM dm, PetscInt depth, PlexLayer *layers, IS strataPermutation)
2388: {
2389:   const PetscInt *permArr;
2390:   PetscInt        d, nPoints;
2391:   MPI_Comm        comm;

2393:   PetscFunctionBegin;
2394:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2395:   PetscCall(ISGetIndices(strataPermutation, &permArr));

2397:   /* Count points, strata offsets and cones offsets (taking strataPermutation into account) */
2398:   {
2399:     PetscInt stratumOffset = 0;
2400:     PetscInt conesOffset   = 0;

2402:     for (d = 0; d <= depth; d++) {
2403:       const PetscInt  e = permArr[d];
2404:       const PlexLayer l = layers[e];
2405:       PetscInt        lo, n, size;

2407:       PetscCall(PetscSectionGetChart(l->coneSizesSection, &lo, &n));
2408:       PetscCall(PetscSectionGetStorageSize(l->coneSizesSection, &size));
2409:       PetscCheck(lo == 0, comm, PETSC_ERR_PLIB, "starting point should be 0 in coneSizesSection %" PetscInt_FMT, d);
2410:       l->offset      = stratumOffset;
2411:       l->conesOffset = conesOffset;
2412:       stratumOffset += n;
2413:       conesOffset += size;
2414:     }
2415:     nPoints = stratumOffset;
2416:   }

2418:   /* Set interval for all plex points */
2419:   //TODO we should store starting point of plex
2420:   PetscCall(DMPlexSetChart(dm, 0, nPoints));

2422:   /* Set up plex coneSection from layer coneSections */
2423:   {
2424:     PetscSection coneSection;

2426:     PetscCall(DMPlexGetConeSection(dm, &coneSection));
2427:     for (d = 0; d <= depth; d++) {
2428:       const PlexLayer l = layers[d];
2429:       PetscInt        n, q;

2431:       PetscCall(PetscSectionGetChart(l->coneSizesSection, NULL, &n));
2432:       for (q = 0; q < n; q++) {
2433:         const PetscInt p = l->offset + q;
2434:         PetscInt       coneSize;

2436:         PetscCall(PetscSectionGetDof(l->coneSizesSection, q, &coneSize));
2437:         PetscCall(PetscSectionSetDof(coneSection, p, coneSize));
2438:       }
2439:     }
2440:   }
2441:   //TODO this is terrible, DMSetUp_Plex() should be DMPlexSetUpSections() or so
2442:   PetscCall(DMSetUp(dm));

2444:   /* Renumber cones points from layer-global numbering to plex-local numbering */
2445:   {
2446:     PetscInt *cones, *ornts;

2448:     PetscCall(DMPlexGetCones(dm, &cones));
2449:     PetscCall(DMPlexGetConeOrientations(dm, &ornts));
2450:     for (d = 1; d <= depth; d++) {
2451:       const PlexLayer l = layers[d];
2452:       PetscInt        i, lConesSize;
2453:       PetscInt       *lCones;
2454:       const PetscInt *lOrnts;
2455:       PetscInt       *pCones = &cones[l->conesOffset];
2456:       PetscInt       *pOrnts = &ornts[l->conesOffset];

2458:       PetscCall(PetscSectionGetStorageSize(l->coneSizesSection, &lConesSize));
2459:       /* Get cones in local plex numbering */
2460:       {
2461:         ISLocalToGlobalMapping l2g;
2462:         PetscLayout            vertexLayout = l->vertexLayout;
2463:         PetscSF                vertexSF     = layers[d - 1]->l2gSF; /* vertices of this layer are cells of previous layer */
2464:         const PetscInt        *gCones;
2465:         PetscInt               lConesSize0;

2467:         PetscCall(ISGetLocalSize(l->conesIS, &lConesSize0));
2468:         PetscCheck(lConesSize0 == lConesSize, comm, PETSC_ERR_PLIB, "layer %" PetscInt_FMT " size(conesIS) = %" PetscInt_FMT " != %" PetscInt_FMT " = storageSize(coneSizesSection)", d, lConesSize0, lConesSize);
2469:         PetscCall(ISGetLocalSize(l->orientationsIS, &lConesSize0));
2470:         PetscCheck(lConesSize0 == lConesSize, comm, PETSC_ERR_PLIB, "layer %" PetscInt_FMT " size(orientationsIS) = %" PetscInt_FMT " != %" PetscInt_FMT " = storageSize(coneSizesSection)", d, lConesSize0, lConesSize);

2472:         PetscCall(PetscMalloc1(lConesSize, &lCones));
2473:         PetscCall(ISGetIndices(l->conesIS, &gCones));
2474:         PetscCall(ISLocalToGlobalMappingCreateSF(vertexSF, vertexLayout->rstart, &l2g));
2475:         PetscCall(ISGlobalToLocalMappingApply(l2g, IS_GTOLM_MASK, lConesSize, gCones, &lConesSize0, lCones));
2476:         PetscCheck(lConesSize0 == lConesSize, comm, PETSC_ERR_PLIB, "global to local does not cover all indices (%" PetscInt_FMT " of %" PetscInt_FMT ")", lConesSize0, lConesSize);
2477:         PetscCall(ISLocalToGlobalMappingDestroy(&l2g));
2478:         PetscCall(ISRestoreIndices(l->conesIS, &gCones));
2479:       }
2480:       PetscCall(ISGetIndices(l->orientationsIS, &lOrnts));
2481:       /* Set cones, need to add stratum offset */
2482:       for (i = 0; i < lConesSize; i++) {
2483:         pCones[i] = lCones[i] + layers[d - 1]->offset; /* cone points of current layer are points of previous layer */
2484:         pOrnts[i] = lOrnts[i];
2485:       }
2486:       PetscCall(PetscFree(lCones));
2487:       PetscCall(ISRestoreIndices(l->orientationsIS, &lOrnts));
2488:     }
2489:   }
2490:   PetscCall(DMPlexSymmetrize(dm));
2491:   PetscCall(DMPlexStratify(dm));
2492:   PetscCall(ISRestoreIndices(strataPermutation, &permArr));
2493:   PetscFunctionReturn(PETSC_SUCCESS);
2494: }

2496: static PetscErrorCode PlexLayerConcatenateSFs_Private(MPI_Comm comm, PetscInt depth, PlexLayer layers[], IS strataPermutation, PetscSF *overlapSF, PetscSF *l2gSF)
2497: {
2498:   PetscInt        d;
2499:   PetscSF        *osfs, *lsfs;
2500:   PetscInt       *leafOffsets;
2501:   const PetscInt *permArr;

2503:   PetscFunctionBegin;
2504:   PetscCall(ISGetIndices(strataPermutation, &permArr));
2505:   PetscCall(PetscCalloc3(depth + 1, &osfs, depth + 1, &lsfs, depth + 1, &leafOffsets));
2506:   for (d = 0; d <= depth; d++) {
2507:     const PetscInt e = permArr[d];

2509:     PetscAssert(e == layers[e]->d, PETSC_COMM_SELF, PETSC_ERR_PLIB, "assertion: e == layers[e]->d");
2510:     osfs[d]        = layers[e]->overlapSF;
2511:     lsfs[d]        = layers[e]->l2gSF;
2512:     leafOffsets[d] = layers[e]->offset;
2513:   }
2514:   PetscCall(PetscSFConcatenate(comm, depth + 1, osfs, PETSCSF_CONCATENATE_ROOTMODE_LOCAL, leafOffsets, overlapSF));
2515:   PetscCall(PetscSFConcatenate(comm, depth + 1, lsfs, PETSCSF_CONCATENATE_ROOTMODE_GLOBAL, leafOffsets, l2gSF));
2516:   PetscCall(PetscFree3(osfs, lsfs, leafOffsets));
2517:   PetscCall(ISRestoreIndices(strataPermutation, &permArr));
2518:   PetscFunctionReturn(PETSC_SUCCESS);
2519: }

2521: // Parallel load of topology
2522: static PetscErrorCode DMPlexTopologyLoad_HDF5_Private(DM dm, PetscViewer viewer, PetscSF *sfXC)
2523: {
2524:   PlexLayer  *layers;
2525:   IS          strataPermutation;
2526:   PetscLayout pointsLayout;
2527:   PetscInt    depth;
2528:   PetscInt    d;
2529:   MPI_Comm    comm;

2531:   PetscFunctionBegin;
2532:   {
2533:     PetscInt dim;

2535:     PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "depth", PETSC_INT, NULL, &depth));
2536:     PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "cell_dim", PETSC_INT, NULL, &dim));
2537:     PetscCall(DMSetDimension(dm, dim));
2538:   }
2539:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));

2541:   PetscCall(PetscInfo(dm, "Loading DM %s in parallel\n", dm->hdr.name));
2542:   {
2543:     IS spOnComm;

2545:     PetscCall(ISCreate(comm, &spOnComm));
2546:     PetscCall(PetscObjectSetName((PetscObject)spOnComm, "permutation"));
2547:     PetscCall(ISLoad(spOnComm, viewer));
2548:     /* have the same serial IS on every rank */
2549:     PetscCall(ISAllGather(spOnComm, &strataPermutation));
2550:     PetscCall(PetscObjectSetName((PetscObject)strataPermutation, ((PetscObject)spOnComm)->name));
2551:     PetscCall(ISDestroy(&spOnComm));
2552:   }

2554:   /* Create layers, load raw data for each layer */
2555:   PetscCall(PetscViewerHDF5PushGroup(viewer, "strata"));
2556:   PetscCall(PetscMalloc1(depth + 1, &layers));
2557:   for (d = depth, pointsLayout = NULL; d >= 0; pointsLayout = layers[d]->vertexLayout, d--) {
2558:     PetscCall(PlexLayerCreate_Private(&layers[d]));
2559:     PetscCall(PlexLayerLoad_Private(layers[d], viewer, d, pointsLayout));
2560:   }
2561:   PetscCall(PetscViewerHDF5PopGroup(viewer)); /* strata */

2563:   for (d = depth; d >= 0; d--) {
2564:     /* Redistribute cells and vertices for each applicable layer */
2565:     if (d < depth) PetscCall(PlexLayerDistribute_Private(layers[d], layers[d]->l2gSF));
2566:     /* Create vertex overlap SF and vertex localToGlobal SF for each applicable layer */
2567:     if (d > 0) PetscCall(PlexLayerCreateSFs_Private(layers[d], &layers[d - 1]->overlapSF, &layers[d - 1]->l2gSF));
2568:   }
2569:   /* Build trivial SFs for the cell layer as well */
2570:   PetscCall(PlexLayerCreateCellSFs_Private(layers[depth], &layers[depth]->overlapSF, &layers[depth]->l2gSF));

2572:   /* Build DMPlex topology from the layers */
2573:   PetscCall(DMPlexTopologyBuildFromLayers_Private(dm, depth, layers, strataPermutation));

2575:   /* Build overall point SF alias overlap SF */
2576:   {
2577:     PetscSF overlapSF;

2579:     PetscCall(PlexLayerConcatenateSFs_Private(comm, depth, layers, strataPermutation, &overlapSF, sfXC));
2580:     PetscCall(DMSetPointSF(dm, overlapSF));
2581:     PetscCall(PetscSFDestroy(&overlapSF));
2582:   }

2584:   for (d = depth; d >= 0; d--) PetscCall(PlexLayerDestroy(&layers[d]));
2585:   PetscCall(PetscFree(layers));
2586:   PetscCall(ISDestroy(&strataPermutation));
2587:   PetscFunctionReturn(PETSC_SUCCESS);
2588: }

2590: PetscErrorCode DMPlexTopologyLoad_HDF5_Internal(DM dm, PetscViewer viewer, PetscSF *sfXC)
2591: {
2592:   DMPlexStorageVersion version;
2593:   const char          *topologydm_name;
2594:   char                 group[PETSC_MAX_PATH_LEN];
2595:   PetscSF              sfwork = NULL;

2597:   PetscFunctionBegin;
2598:   PetscCall(DMPlexGetHDF5Name_Private(dm, &topologydm_name));
2599:   PetscCall(PetscViewerHDF5GetDMPlexStorageVersionReading(viewer, &version));
2600:   if (DMPlexStorageVersionGE(version, 2, 0, 0)) {
2601:     PetscCall(PetscSNPrintf(group, sizeof(group), "topologies/%s", topologydm_name));
2602:   } else {
2603:     PetscCall(PetscStrncpy(group, "/", sizeof(group)));
2604:   }
2605:   PetscCall(PetscViewerHDF5PushGroup(viewer, group));

2607:   PetscCall(PetscViewerHDF5PushGroup(viewer, "topology"));
2608:   if (version->major < 3) {
2609:     PetscCall(DMPlexTopologyLoad_HDF5_Legacy_Private(dm, viewer, &sfwork));
2610:   } else {
2611:     /* since DMPlexStorageVersion 3.0.0 */
2612:     PetscCall(DMPlexTopologyLoad_HDF5_Private(dm, viewer, &sfwork));
2613:   }
2614:   PetscCall(PetscViewerHDF5PopGroup(viewer)); /* "topology" */

2616:   if (DMPlexStorageVersionGE(version, 2, 0, 0)) {
2617:     DM          distdm;
2618:     PetscSF     distsf;
2619:     const char *distribution_name;
2620:     PetscBool   exists;

2622:     PetscCall(DMPlexDistributionGetName(dm, &distribution_name));
2623:     /* this group is guaranteed to be present since DMPlexStorageVersion 2.1.0 */
2624:     PetscCall(PetscViewerHDF5PushGroup(viewer, "distributions"));
2625:     PetscCall(PetscViewerHDF5HasGroup(viewer, NULL, &exists));
2626:     if (exists) {
2627:       PetscCall(PetscViewerHDF5PushGroup(viewer, distribution_name));
2628:       PetscCall(DMPlexDistributionLoad_HDF5_Private(dm, viewer, sfwork, &distsf, &distdm));
2629:       if (distdm) {
2630:         PetscCall(DMPlexReplace_Internal(dm, &distdm));
2631:         PetscCall(PetscSFDestroy(&sfwork));
2632:         sfwork = distsf;
2633:       }
2634:       PetscCall(PetscViewerHDF5PopGroup(viewer)); /* distribution_name */
2635:     }
2636:     PetscCall(PetscViewerHDF5PopGroup(viewer)); /* "distributions" */
2637:   }
2638:   if (sfXC) {
2639:     *sfXC = sfwork;
2640:   } else {
2641:     PetscCall(PetscSFDestroy(&sfwork));
2642:   }

2644:   PetscCall(PetscViewerHDF5PopGroup(viewer));
2645:   PetscFunctionReturn(PETSC_SUCCESS);
2646: }

2648: /* If the file is old, it not only has different path to the coordinates, but   */
2649: /* does not contain coordinateDMs, so must fall back to the old implementation. */
2650: static PetscErrorCode DMPlexCoordinatesLoad_HDF5_Legacy_Private(DM dm, PetscViewer viewer)
2651: {
2652:   PetscSection coordSection;
2653:   Vec          coordinates;
2654:   PetscReal    lengthScale;
2655:   PetscInt     spatialDim, N, numVertices, vStart, vEnd, v;
2656:   PetscMPIInt  rank;

2658:   PetscFunctionBegin;
2659:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
2660:   /* Read geometry */
2661:   PetscCall(PetscViewerHDF5PushGroup(viewer, "/geometry"));
2662:   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), &coordinates));
2663:   PetscCall(PetscObjectSetName((PetscObject)coordinates, "vertices"));
2664:   {
2665:     /* Force serial load */
2666:     PetscCall(PetscViewerHDF5ReadSizes(viewer, "vertices", &spatialDim, &N));
2667:     PetscCall(VecSetSizes(coordinates, rank == 0 ? N : 0, N));
2668:     PetscCall(VecSetBlockSize(coordinates, spatialDim));
2669:   }
2670:   PetscCall(VecLoad(coordinates, viewer));
2671:   PetscCall(PetscViewerHDF5PopGroup(viewer));
2672:   PetscCall(DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale));
2673:   PetscCall(VecScale(coordinates, 1.0 / lengthScale));
2674:   PetscCall(VecGetLocalSize(coordinates, &numVertices));
2675:   PetscCall(VecGetBlockSize(coordinates, &spatialDim));
2676:   numVertices /= spatialDim;
2677:   /* Create coordinates */
2678:   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2679:   PetscCheck(numVertices == vEnd - vStart, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of coordinates loaded %" PetscInt_FMT " does not match number of vertices %" PetscInt_FMT, numVertices, vEnd - vStart);
2680:   PetscCall(DMGetCoordinateSection(dm, &coordSection));
2681:   PetscCall(PetscSectionSetNumFields(coordSection, 1));
2682:   PetscCall(PetscSectionSetFieldComponents(coordSection, 0, spatialDim));
2683:   PetscCall(PetscSectionSetChart(coordSection, vStart, vEnd));
2684:   for (v = vStart; v < vEnd; ++v) {
2685:     PetscCall(PetscSectionSetDof(coordSection, v, spatialDim));
2686:     PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, spatialDim));
2687:   }
2688:   PetscCall(PetscSectionSetUp(coordSection));
2689:   PetscCall(DMSetCoordinates(dm, coordinates));
2690:   PetscCall(VecDestroy(&coordinates));
2691:   PetscFunctionReturn(PETSC_SUCCESS);
2692: }

2694: PetscErrorCode DMPlexCoordinatesLoad_HDF5_Internal(DM dm, PetscViewer viewer, PetscSF sfXC)
2695: {
2696:   DMPlexStorageVersion version;
2697:   DM                   cdm;
2698:   Vec                  coords;
2699:   PetscInt             blockSize;
2700:   PetscReal            lengthScale;
2701:   PetscSF              lsf;
2702:   const char          *topologydm_name;
2703:   char                *coordinatedm_name, *coordinates_name;

2705:   PetscFunctionBegin;
2706:   PetscCall(PetscViewerHDF5GetDMPlexStorageVersionReading(viewer, &version));
2707:   if (!DMPlexStorageVersionGE(version, 2, 0, 0)) {
2708:     PetscCall(DMPlexCoordinatesLoad_HDF5_Legacy_Private(dm, viewer));
2709:     PetscFunctionReturn(PETSC_SUCCESS);
2710:   }
2711:   /* else: since DMPlexStorageVersion 2.0.0 */
2712:   PetscCheck(sfXC, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_NULL, "PetscSF must be given for parallel load");
2713:   PetscCall(DMPlexGetHDF5Name_Private(dm, &topologydm_name));
2714:   PetscCall(PetscViewerHDF5PushGroup(viewer, "topologies"));
2715:   PetscCall(PetscViewerHDF5PushGroup(viewer, topologydm_name));
2716:   PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "coordinateDMName", PETSC_STRING, NULL, &coordinatedm_name));
2717:   PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "coordinatesName", PETSC_STRING, NULL, &coordinates_name));
2718:   PetscCall(PetscViewerHDF5PopGroup(viewer));
2719:   PetscCall(PetscViewerHDF5PopGroup(viewer));
2720:   PetscCall(DMGetCoordinateDM(dm, &cdm));
2721:   PetscCall(PetscObjectSetName((PetscObject)cdm, coordinatedm_name));
2722:   PetscCall(PetscFree(coordinatedm_name));
2723:   /* lsf: on-disk data -> in-memory local vector associated with cdm's local section */
2724:   PetscCall(DMPlexSectionLoad(dm, viewer, cdm, sfXC, NULL, &lsf));
2725:   PetscCall(DMCreateLocalVector(cdm, &coords));
2726:   PetscCall(PetscObjectSetName((PetscObject)coords, coordinates_name));
2727:   PetscCall(PetscFree(coordinates_name));
2728:   PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_NATIVE));
2729:   PetscCall(DMPlexLocalVectorLoad(dm, viewer, cdm, lsf, coords));
2730:   PetscCall(PetscViewerPopFormat(viewer));
2731:   PetscCall(DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale));
2732:   PetscCall(VecScale(coords, 1.0 / lengthScale));
2733:   PetscCall(DMSetCoordinatesLocal(dm, coords));
2734:   PetscCall(VecGetBlockSize(coords, &blockSize));
2735:   PetscCall(DMSetCoordinateDim(dm, blockSize));
2736:   PetscCall(VecDestroy(&coords));
2737:   PetscCall(PetscSFDestroy(&lsf));
2738:   PetscFunctionReturn(PETSC_SUCCESS);
2739: }

2741: PetscErrorCode DMPlexLoad_HDF5_Internal(DM dm, PetscViewer viewer)
2742: {
2743:   DMPlexStorageVersion version;

2745:   PetscFunctionBegin;
2746:   PetscCall(PetscViewerHDF5GetDMPlexStorageVersionReading(viewer, &version));
2747:   PetscCall(PetscInfo(dm, "Loading DM %s storage version %d.%d.%d\n", dm->hdr.name, version->major, version->minor, version->subminor));
2748:   if (!DMPlexStorageVersionGE(version, 2, 0, 0)) {
2749:     PetscCall(DMPlexTopologyLoad_HDF5_Internal(dm, viewer, NULL));
2750:     PetscCall(DMPlexLabelsLoad_HDF5_Internal(dm, viewer, NULL));
2751:     PetscCall(DMPlexCoordinatesLoad_HDF5_Legacy_Private(dm, viewer));
2752:   } else {
2753:     PetscSF sfXC;

2755:     /* since DMPlexStorageVersion 2.0.0 */
2756:     PetscCall(DMPlexTopologyLoad_HDF5_Internal(dm, viewer, &sfXC));
2757:     PetscCall(DMPlexLabelsLoad_HDF5_Internal(dm, viewer, sfXC));
2758:     PetscCall(DMPlexCoordinatesLoad_HDF5_Internal(dm, viewer, sfXC));
2759:     PetscCall(PetscSFDestroy(&sfXC));
2760:   }
2761:   PetscFunctionReturn(PETSC_SUCCESS);
2762: }

2764: static PetscErrorCode DMPlexSectionLoad_HDF5_Internal_CreateDataSF(PetscSection rootSection, PetscLayout layout, PetscInt globalOffsets[], PetscSection leafSection, PetscSF *sectionSF)
2765: {
2766:   MPI_Comm  comm;
2767:   PetscInt  pStart, pEnd, p, m;
2768:   PetscInt *goffs, *ilocal;
2769:   PetscBool rootIncludeConstraints, leafIncludeConstraints;

2771:   PetscFunctionBegin;
2772:   PetscCall(PetscObjectGetComm((PetscObject)leafSection, &comm));
2773:   PetscCall(PetscSectionGetChart(leafSection, &pStart, &pEnd));
2774:   PetscCall(PetscSectionGetIncludesConstraints(rootSection, &rootIncludeConstraints));
2775:   PetscCall(PetscSectionGetIncludesConstraints(leafSection, &leafIncludeConstraints));
2776:   if (rootIncludeConstraints && leafIncludeConstraints) PetscCall(PetscSectionGetStorageSize(leafSection, &m));
2777:   else PetscCall(PetscSectionGetConstrainedStorageSize(leafSection, &m));
2778:   PetscCall(PetscMalloc1(m, &ilocal));
2779:   PetscCall(PetscMalloc1(m, &goffs));
2780:   /* Currently, PetscSFDistributeSection() returns globalOffsets[] only */
2781:   /* for the top-level section (not for each field), so one must have   */
2782:   /* rootSection->pointMajor == PETSC_TRUE.                             */
2783:   PetscCheck(rootSection->pointMajor, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for field major ordering");
2784:   /* Currently, we also assume that leafSection->pointMajor == PETSC_TRUE. */
2785:   PetscCheck(leafSection->pointMajor, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for field major ordering");
2786:   for (p = pStart, m = 0; p < pEnd; ++p) {
2787:     PetscInt        dof, cdof, i, j, off, goff;
2788:     const PetscInt *cinds;

2790:     PetscCall(PetscSectionGetDof(leafSection, p, &dof));
2791:     if (dof < 0) continue;
2792:     goff = globalOffsets[p - pStart];
2793:     PetscCall(PetscSectionGetOffset(leafSection, p, &off));
2794:     PetscCall(PetscSectionGetConstraintDof(leafSection, p, &cdof));
2795:     PetscCall(PetscSectionGetConstraintIndices(leafSection, p, &cinds));
2796:     for (i = 0, j = 0; i < dof; ++i) {
2797:       PetscBool constrained = (PetscBool)(j < cdof && i == cinds[j]);

2799:       if (!constrained || (leafIncludeConstraints && rootIncludeConstraints)) {
2800:         ilocal[m]  = off++;
2801:         goffs[m++] = goff++;
2802:       } else if (leafIncludeConstraints && !rootIncludeConstraints) ++off;
2803:       else if (!leafIncludeConstraints && rootIncludeConstraints) ++goff;
2804:       if (constrained) ++j;
2805:     }
2806:   }
2807:   PetscCall(PetscSFCreate(comm, sectionSF));
2808:   PetscCall(PetscSFSetFromOptions(*sectionSF));
2809:   PetscCall(PetscSFSetGraphLayout(*sectionSF, layout, m, ilocal, PETSC_OWN_POINTER, goffs));
2810:   PetscCall(PetscFree(goffs));
2811:   PetscFunctionReturn(PETSC_SUCCESS);
2812: }

2814: PetscErrorCode DMPlexSectionLoad_HDF5_Internal(DM dm, PetscViewer viewer, DM sectiondm, PetscSF sfXB, PetscSF *gsf, PetscSF *lsf)
2815: {
2816:   MPI_Comm     comm;
2817:   PetscMPIInt  size, rank;
2818:   const char  *topologydm_name;
2819:   const char  *sectiondm_name;
2820:   PetscSection sectionA, sectionB;
2821:   PetscInt     nX, n, i;
2822:   PetscSF      sfAB;

2824:   PetscFunctionBegin;
2825:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2826:   PetscCallMPI(MPI_Comm_size(comm, &size));
2827:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2828:   PetscCall(DMPlexGetHDF5Name_Private(dm, &topologydm_name));
2829:   PetscCall(PetscObjectGetName((PetscObject)sectiondm, &sectiondm_name));
2830:   PetscCall(PetscViewerHDF5PushGroup(viewer, "topologies"));
2831:   PetscCall(PetscViewerHDF5PushGroup(viewer, topologydm_name));
2832:   PetscCall(PetscViewerHDF5PushGroup(viewer, "dms"));
2833:   PetscCall(PetscViewerHDF5PushGroup(viewer, sectiondm_name));
2834:   /* A: on-disk points                        */
2835:   /* X: list of global point numbers, [0, NX) */
2836:   /* B: plex points                           */
2837:   /* Load raw section (sectionA)              */
2838:   PetscCall(PetscSectionCreate(comm, &sectionA));
2839:   PetscCall(PetscSectionLoad(sectionA, viewer));
2840:   PetscCall(PetscSectionGetChart(sectionA, NULL, &n));
2841:   /* Create sfAB: A -> B */
2842:   #if defined(PETSC_USE_DEBUG)
2843:   {
2844:     PetscInt N, N1;

2846:     PetscCall(PetscViewerHDF5ReadSizes(viewer, "order", NULL, &N1));
2847:     PetscCallMPI(MPIU_Allreduce(&n, &N, 1, MPIU_INT, MPI_SUM, comm));
2848:     PetscCheck(N1 == N, comm, PETSC_ERR_ARG_SIZ, "Mismatching sizes: on-disk order array size (%" PetscInt_FMT ") != number of loaded section points (%" PetscInt_FMT ")", N1, N);
2849:   }
2850:   #endif
2851:   {
2852:     IS              orderIS;
2853:     const PetscInt *gpoints;
2854:     PetscSF         sfXA, sfAX;
2855:     PetscLayout     layout;
2856:     PetscSFNode    *owners, *buffer;
2857:     PetscInt        nleaves;
2858:     PetscInt       *ilocal;
2859:     PetscSFNode    *iremote;

2861:     /* Create sfAX: A -> X */
2862:     PetscCall(ISCreate(comm, &orderIS));
2863:     PetscCall(PetscObjectSetName((PetscObject)orderIS, "order"));
2864:     PetscCall(PetscLayoutSetLocalSize(orderIS->map, n));
2865:     PetscCall(ISLoad(orderIS, viewer));
2866:     PetscCall(PetscLayoutCreate(comm, &layout));
2867:     PetscCall(PetscSFGetGraph(sfXB, &nX, NULL, NULL, NULL));
2868:     PetscCall(PetscLayoutSetLocalSize(layout, nX));
2869:     PetscCall(PetscLayoutSetBlockSize(layout, 1));
2870:     PetscCall(PetscLayoutSetUp(layout));
2871:     PetscCall(PetscSFCreate(comm, &sfXA));
2872:     PetscCall(ISGetIndices(orderIS, &gpoints));
2873:     PetscCall(PetscSFSetGraphLayout(sfXA, layout, n, NULL, PETSC_OWN_POINTER, gpoints));
2874:     PetscCall(ISRestoreIndices(orderIS, &gpoints));
2875:     PetscCall(ISDestroy(&orderIS));
2876:     PetscCall(PetscLayoutDestroy(&layout));
2877:     PetscCall(PetscMalloc1(n, &owners));
2878:     PetscCall(PetscMalloc1(nX, &buffer));
2879:     for (i = 0; i < n; ++i) {
2880:       owners[i].rank  = rank;
2881:       owners[i].index = i;
2882:     }
2883:     for (i = 0; i < nX; ++i) {
2884:       buffer[i].rank  = -1;
2885:       buffer[i].index = -1;
2886:     }
2887:     PetscCall(PetscSFReduceBegin(sfXA, MPIU_SF_NODE, owners, buffer, MPI_MAXLOC));
2888:     PetscCall(PetscSFReduceEnd(sfXA, MPIU_SF_NODE, owners, buffer, MPI_MAXLOC));
2889:     PetscCall(PetscSFDestroy(&sfXA));
2890:     PetscCall(PetscFree(owners));
2891:     for (i = 0, nleaves = 0; i < nX; ++i)
2892:       if (buffer[i].rank >= 0) nleaves++;
2893:     PetscCall(PetscMalloc1(nleaves, &ilocal));
2894:     PetscCall(PetscMalloc1(nleaves, &iremote));
2895:     for (i = 0, nleaves = 0; i < nX; ++i) {
2896:       if (buffer[i].rank >= 0) {
2897:         ilocal[nleaves]        = i;
2898:         iremote[nleaves].rank  = buffer[i].rank;
2899:         iremote[nleaves].index = buffer[i].index;
2900:         nleaves++;
2901:       }
2902:     }
2903:     PetscCall(PetscFree(buffer));
2904:     PetscCall(PetscSFCreate(comm, &sfAX));
2905:     PetscCall(PetscSFSetFromOptions(sfAX));
2906:     PetscCall(PetscSFSetGraph(sfAX, n, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
2907:     PetscCall(PetscSFCompose(sfAX, sfXB, &sfAB));
2908:     PetscCall(PetscSFDestroy(&sfAX));
2909:   }
2910:   PetscCall(PetscViewerHDF5PopGroup(viewer));
2911:   PetscCall(PetscViewerHDF5PopGroup(viewer));
2912:   PetscCall(PetscViewerHDF5PopGroup(viewer));
2913:   PetscCall(PetscViewerHDF5PopGroup(viewer));
2914:   /* Create plex section (sectionB) */
2915:   PetscCall(DMGetLocalSection(sectiondm, &sectionB));
2916:   if (lsf || gsf) {
2917:     PetscLayout layout;
2918:     PetscInt    M, m;
2919:     PetscInt   *offsetsA;
2920:     PetscBool   includesConstraintsA;

2922:     PetscCall(PetscSFDistributeSection(sfAB, sectionA, &offsetsA, sectionB));
2923:     PetscCall(PetscSectionGetIncludesConstraints(sectionA, &includesConstraintsA));
2924:     if (includesConstraintsA) PetscCall(PetscSectionGetStorageSize(sectionA, &m));
2925:     else PetscCall(PetscSectionGetConstrainedStorageSize(sectionA, &m));
2926:     PetscCallMPI(MPIU_Allreduce(&m, &M, 1, MPIU_INT, MPI_SUM, comm));
2927:     PetscCall(PetscLayoutCreate(comm, &layout));
2928:     PetscCall(PetscLayoutSetSize(layout, M));
2929:     PetscCall(PetscLayoutSetUp(layout));
2930:     if (lsf) {
2931:       PetscSF lsfABdata;

2933:       PetscCall(DMPlexSectionLoad_HDF5_Internal_CreateDataSF(sectionA, layout, offsetsA, sectionB, &lsfABdata));
2934:       *lsf = lsfABdata;
2935:     }
2936:     if (gsf) {
2937:       PetscSection gsectionB, gsectionB1;
2938:       PetscBool    includesConstraintsB;
2939:       PetscSF      gsfABdata, pointsf;

2941:       PetscCall(DMGetGlobalSection(sectiondm, &gsectionB1));
2942:       PetscCall(PetscSectionGetIncludesConstraints(gsectionB1, &includesConstraintsB));
2943:       PetscCall(DMGetPointSF(sectiondm, &pointsf));
2944:       PetscCall(PetscSectionCreateGlobalSection(sectionB, pointsf, PETSC_TRUE, includesConstraintsB, PETSC_TRUE, &gsectionB));
2945:       PetscCall(DMPlexSectionLoad_HDF5_Internal_CreateDataSF(sectionA, layout, offsetsA, gsectionB, &gsfABdata));
2946:       PetscCall(PetscSectionDestroy(&gsectionB));
2947:       *gsf = gsfABdata;
2948:     }
2949:     PetscCall(PetscLayoutDestroy(&layout));
2950:     PetscCall(PetscFree(offsetsA));
2951:   } else {
2952:     PetscCall(PetscSFDistributeSection(sfAB, sectionA, NULL, sectionB));
2953:   }
2954:   PetscCall(PetscSFDestroy(&sfAB));
2955:   PetscCall(PetscSectionDestroy(&sectionA));
2956:   PetscFunctionReturn(PETSC_SUCCESS);
2957: }

2959: PetscErrorCode DMPlexVecLoad_HDF5_Internal(DM dm, PetscViewer viewer, DM sectiondm, PetscSF sf, Vec vec)
2960: {
2961:   MPI_Comm           comm;
2962:   const char        *topologydm_name;
2963:   const char        *sectiondm_name;
2964:   const char        *vec_name;
2965:   Vec                vecA;
2966:   PetscInt           mA, m, bs;
2967:   const PetscInt    *ilocal;
2968:   const PetscScalar *src;
2969:   PetscScalar       *dest;

2971:   PetscFunctionBegin;
2972:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2973:   PetscCall(DMPlexGetHDF5Name_Private(dm, &topologydm_name));
2974:   PetscCall(PetscObjectGetName((PetscObject)sectiondm, &sectiondm_name));
2975:   PetscCall(PetscObjectGetName((PetscObject)vec, &vec_name));
2976:   PetscCall(PetscViewerHDF5PushGroup(viewer, "topologies"));
2977:   PetscCall(PetscViewerHDF5PushGroup(viewer, topologydm_name));
2978:   PetscCall(PetscViewerHDF5PushGroup(viewer, "dms"));
2979:   PetscCall(PetscViewerHDF5PushGroup(viewer, sectiondm_name));
2980:   PetscCall(PetscViewerHDF5PushGroup(viewer, "vecs"));
2981:   PetscCall(PetscViewerHDF5PushGroup(viewer, vec_name));
2982:   PetscCall(VecCreate(comm, &vecA));
2983:   PetscCall(PetscObjectSetName((PetscObject)vecA, vec_name));
2984:   PetscCall(PetscSFGetGraph(sf, &mA, &m, &ilocal, NULL));
2985:   /* Check consistency */
2986:   {
2987:     PetscSF  pointsf, pointsf1;
2988:     PetscInt m1, i, j;

2990:     PetscCall(DMGetPointSF(dm, &pointsf));
2991:     PetscCall(DMGetPointSF(sectiondm, &pointsf1));
2992:     PetscCheck(pointsf1 == pointsf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Mismatching point SFs for dm and sectiondm");
2993:   #if defined(PETSC_USE_DEBUG)
2994:     {
2995:       PetscInt MA, MA1;

2997:       PetscCallMPI(MPIU_Allreduce(&mA, &MA, 1, MPIU_INT, MPI_SUM, comm));
2998:       PetscCall(PetscViewerHDF5ReadSizes(viewer, vec_name, NULL, &MA1));
2999:       PetscCheck(MA1 == MA, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total SF root size (%" PetscInt_FMT ") != On-disk vector data size (%" PetscInt_FMT ")", MA, MA1);
3000:     }
3001:   #endif
3002:     PetscCall(VecGetLocalSize(vec, &m1));
3003:     PetscCheck(m1 >= m, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Target vector size (%" PetscInt_FMT ") < SF leaf size (%" PetscInt_FMT ")", m1, m);
3004:     for (i = 0; i < m; ++i) {
3005:       j = ilocal ? ilocal[i] : i;
3006:       PetscCheck(j >= 0 && j < m1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Leaf's %" PetscInt_FMT "-th index, %" PetscInt_FMT ", not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", i, j, (PetscInt)0, m1);
3007:     }
3008:   }
3009:   PetscCall(VecSetSizes(vecA, mA, PETSC_DECIDE));
3010:   PetscCall(VecLoad(vecA, viewer));
3011:   PetscCall(VecGetArrayRead(vecA, &src));
3012:   PetscCall(VecGetArray(vec, &dest));
3013:   PetscCall(PetscSFBcastBegin(sf, MPIU_SCALAR, src, dest, MPI_REPLACE));
3014:   PetscCall(PetscSFBcastEnd(sf, MPIU_SCALAR, src, dest, MPI_REPLACE));
3015:   PetscCall(VecRestoreArray(vec, &dest));
3016:   PetscCall(VecRestoreArrayRead(vecA, &src));
3017:   PetscCall(VecDestroy(&vecA));
3018:   PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "blockSize", PETSC_INT, NULL, (void *)&bs));
3019:   PetscCall(VecSetBlockSize(vec, bs));
3020:   PetscCall(PetscViewerHDF5PopGroup(viewer));
3021:   PetscCall(PetscViewerHDF5PopGroup(viewer));
3022:   PetscCall(PetscViewerHDF5PopGroup(viewer));
3023:   PetscCall(PetscViewerHDF5PopGroup(viewer));
3024:   PetscCall(PetscViewerHDF5PopGroup(viewer));
3025:   PetscCall(PetscViewerHDF5PopGroup(viewer));
3026:   PetscFunctionReturn(PETSC_SUCCESS);
3027: }
3028: #endif