Actual source code: plexhdf5.c

petsc-3.12.5 2020-03-29
Report Typos and Errors
  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: PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer);

  9: #if defined(PETSC_HAVE_HDF5)
 10: static PetscErrorCode DMSequenceView_HDF5(DM dm, const char *seqname, PetscInt seqnum, PetscScalar value, PetscViewer viewer)
 11: {
 12:   Vec            stamp;
 13:   PetscMPIInt    rank;

 17:   if (seqnum < 0) return(0);
 18:   MPI_Comm_rank(PetscObjectComm((PetscObject) viewer), &rank);
 19:   VecCreateMPI(PetscObjectComm((PetscObject) viewer), rank ? 0 : 1, 1, &stamp);
 20:   VecSetBlockSize(stamp, 1);
 21:   PetscObjectSetName((PetscObject) stamp, seqname);
 22:   if (!rank) {
 23:     PetscReal timeScale;
 24:     PetscBool istime;

 26:     PetscStrncmp(seqname, "time", 5, &istime);
 27:     if (istime) {DMPlexGetScale(dm, PETSC_UNIT_TIME, &timeScale); value *= timeScale;}
 28:     VecSetValue(stamp, 0, value, INSERT_VALUES);
 29:   }
 30:   VecAssemblyBegin(stamp);
 31:   VecAssemblyEnd(stamp);
 32:   PetscViewerHDF5PushGroup(viewer, "/");
 33:   PetscViewerHDF5SetTimestep(viewer, seqnum);
 34:   VecView(stamp, viewer);
 35:   PetscViewerHDF5PopGroup(viewer);
 36:   VecDestroy(&stamp);
 37:   return(0);
 38: }

 40: PetscErrorCode DMSequenceLoad_HDF5_Internal(DM dm, const char *seqname, PetscInt seqnum, PetscScalar *value, PetscViewer viewer)
 41: {
 42:   Vec            stamp;
 43:   PetscMPIInt    rank;

 47:   if (seqnum < 0) return(0);
 48:   MPI_Comm_rank(PetscObjectComm((PetscObject) viewer), &rank);
 49:   VecCreateMPI(PetscObjectComm((PetscObject) viewer), rank ? 0 : 1, 1, &stamp);
 50:   VecSetBlockSize(stamp, 1);
 51:   PetscObjectSetName((PetscObject) stamp, seqname);
 52:   PetscViewerHDF5PushGroup(viewer, "/");
 53:   PetscViewerHDF5SetTimestep(viewer, seqnum);
 54:   VecLoad(stamp, viewer);
 55:   PetscViewerHDF5PopGroup(viewer);
 56:   if (!rank) {
 57:     const PetscScalar *a;
 58:     PetscReal timeScale;
 59:     PetscBool istime;

 61:     VecGetArrayRead(stamp, &a);
 62:     *value = a[0];
 63:     VecRestoreArrayRead(stamp, &a);
 64:     PetscStrncmp(seqname, "time", 5, &istime);
 65:     if (istime) {DMPlexGetScale(dm, PETSC_UNIT_TIME, &timeScale); *value /= timeScale;}
 66:   }
 67:   VecDestroy(&stamp);
 68:   return(0);
 69: }

 71: static PetscErrorCode DMPlexCreateCutVertexLabel_Private(DM dm, DMLabel cutLabel, DMLabel *cutVertexLabel)
 72: {
 73:   IS              cutcells = NULL;
 74:   const PetscInt *cutc;
 75:   PetscInt        cellHeight, vStart, vEnd, cStart, cEnd, c;
 76:   PetscErrorCode  ierr;

 79:   if (!cutLabel) return(0);
 80:   DMPlexGetVTKCellHeight(dm, &cellHeight);
 81:   DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);
 82:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
 83:   /* Label vertices that should be duplicated */
 84:   DMLabelCreate(PETSC_COMM_SELF, "Cut Vertices", cutVertexLabel);
 85:   DMLabelGetStratumIS(cutLabel, 2, &cutcells);
 86:   if (cutcells) {
 87:     PetscInt n;

 89:     ISGetIndices(cutcells, &cutc);
 90:     ISGetLocalSize(cutcells, &n);
 91:     for (c = 0; c < n; ++c) {
 92:       if ((cutc[c] >= cStart) && (cutc[c] < cEnd)) {
 93:         PetscInt *closure = NULL;
 94:         PetscInt  closureSize, cl, value;

 96:         DMPlexGetTransitiveClosure(dm, cutc[c], PETSC_TRUE, &closureSize, &closure);
 97:         for (cl = 0; cl < closureSize*2; cl += 2) {
 98:           if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) {
 99:             DMLabelGetValue(cutLabel, closure[cl], &value);
100:             if (value == 1) {
101:               DMLabelSetValue(*cutVertexLabel, closure[cl], 1);
102:             }
103:           }
104:         }
105:         DMPlexRestoreTransitiveClosure(dm, cutc[c], PETSC_TRUE, &closureSize, &closure);
106:       }
107:     }
108:     ISRestoreIndices(cutcells, &cutc);
109:   }
110:   ISDestroy(&cutcells);
111:   return(0);
112: }

114: PetscErrorCode VecView_Plex_Local_HDF5_Internal(Vec v, PetscViewer viewer)
115: {
116:   DM                      dm;
117:   DM                      dmBC;
118:   PetscSection            section, sectionGlobal;
119:   Vec                     gv;
120:   const char             *name;
121:   PetscViewerVTKFieldType ft;
122:   PetscViewerFormat       format;
123:   PetscInt                seqnum;
124:   PetscReal               seqval;
125:   PetscBool               isseq;
126:   PetscErrorCode          ierr;

129:   PetscObjectTypeCompare((PetscObject) v, VECSEQ, &isseq);
130:   VecGetDM(v, &dm);
131:   DMGetLocalSection(dm, &section);
132:   DMGetOutputSequenceNumber(dm, &seqnum, &seqval);
133:   PetscViewerHDF5SetTimestep(viewer, seqnum);
134:   DMSequenceView_HDF5(dm, "time", seqnum, (PetscScalar) seqval, viewer);
135:   PetscViewerGetFormat(viewer, &format);
136:   DMGetOutputDM(dm, &dmBC);
137:   DMGetGlobalSection(dmBC, &sectionGlobal);
138:   DMGetGlobalVector(dmBC, &gv);
139:   PetscObjectGetName((PetscObject) v, &name);
140:   PetscObjectSetName((PetscObject) gv, name);
141:   DMLocalToGlobalBegin(dmBC, v, INSERT_VALUES, gv);
142:   DMLocalToGlobalEnd(dmBC, v, INSERT_VALUES, gv);
143:   PetscObjectTypeCompare((PetscObject) gv, VECSEQ, &isseq);
144:   if (isseq) {VecView_Seq(gv, viewer);}
145:   else       {VecView_MPI(gv, viewer);}
146:   if (format == PETSC_VIEWER_HDF5_VIZ) {
147:     /* Output visualization representation */
148:     PetscInt numFields, f;
149:     DMLabel  cutLabel, cutVertexLabel = NULL;

151:     PetscSectionGetNumFields(section, &numFields);
152:     DMGetLabel(dm, "periodic_cut", &cutLabel);
153:     for (f = 0; f < numFields; ++f) {
154:       Vec         subv;
155:       IS          is;
156:       const char *fname, *fgroup;
157:       char        subname[PETSC_MAX_PATH_LEN];
158:       PetscInt    pStart, pEnd;

160:       DMPlexGetFieldType_Internal(dm, section, f, &pStart, &pEnd, &ft);
161:       fgroup = (ft == PETSC_VTK_POINT_VECTOR_FIELD) || (ft == PETSC_VTK_POINT_FIELD) ? "/vertex_fields" : "/cell_fields";
162:       PetscSectionGetFieldName(section, f, &fname);
163:       if (!fname) continue;
164:       PetscViewerHDF5PushGroup(viewer, fgroup);
165:       if (cutLabel) {
166:         const PetscScalar *ga;
167:         PetscScalar       *suba;
168:         PetscInt           Nc, gstart, subSize = 0, extSize = 0, subOff = 0, newOff = 0, p;

170:         DMPlexCreateCutVertexLabel_Private(dm, cutLabel, &cutVertexLabel);
171:         PetscSectionGetFieldComponents(section, f, &Nc);
172:         for (p = pStart; p < pEnd; ++p) {
173:           PetscInt gdof, fdof = 0, val;

175:           PetscSectionGetDof(sectionGlobal, p, &gdof);
176:           if (gdof > 0) {PetscSectionGetFieldDof(section, p, f, &fdof);}
177:           subSize += fdof;
178:           DMLabelGetValue(cutVertexLabel, p, &val);
179:           if (val == 1) extSize += fdof;
180:         }
181:         VecCreate(PetscObjectComm((PetscObject) gv), &subv);
182:         VecSetSizes(subv, subSize+extSize, PETSC_DETERMINE);
183:         VecSetBlockSize(subv, Nc);
184:         VecSetType(subv, VECSTANDARD);
185:         VecGetOwnershipRange(gv, &gstart, NULL);
186:         VecGetArrayRead(gv, &ga);
187:         VecGetArray(subv, &suba);
188:         for (p = pStart; p < pEnd; ++p) {
189:           PetscInt gdof, goff, val;

191:           PetscSectionGetDof(sectionGlobal, p, &gdof);
192:           if (gdof > 0) {
193:             PetscInt fdof, fc, f2, poff = 0;

195:             PetscSectionGetOffset(sectionGlobal, p, &goff);
196:             /* Can get rid of this loop by storing field information in the global section */
197:             for (f2 = 0; f2 < f; ++f2) {
198:               PetscSectionGetFieldDof(section, p, f2, &fdof);
199:               poff += fdof;
200:             }
201:             PetscSectionGetFieldDof(section, p, f, &fdof);
202:             for (fc = 0; fc < fdof; ++fc, ++subOff) suba[subOff] = ga[goff+poff+fc - gstart];
203:             DMLabelGetValue(cutVertexLabel, p, &val);
204:             if (val == 1) {
205:               for (fc = 0; fc < fdof; ++fc, ++newOff) suba[subSize+newOff] = ga[goff+poff+fc - gstart];
206:             }
207:           }
208:         }
209:         VecRestoreArrayRead(gv, &ga);
210:         VecRestoreArray(subv, &suba);
211:         DMLabelDestroy(&cutVertexLabel);
212:       } else {
213:         PetscSectionGetField_Internal(section, sectionGlobal, gv, f, pStart, pEnd, &is, &subv);
214:       }
215:       PetscStrncpy(subname, name,sizeof(subname));
216:       PetscStrlcat(subname, "_",sizeof(subname));
217:       PetscStrlcat(subname, fname,sizeof(subname));
218:       PetscObjectSetName((PetscObject) subv, subname);
219:       if (isseq) {VecView_Seq(subv, viewer);}
220:       else       {VecView_MPI(subv, viewer);}
221:       if ((ft == PETSC_VTK_POINT_VECTOR_FIELD) || (ft == PETSC_VTK_CELL_VECTOR_FIELD)) {
222:         PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) subv, "vector_field_type", PETSC_STRING, "vector");
223:       } else {
224:         PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) subv, "vector_field_type", PETSC_STRING, "scalar");
225:       }
226:       if (cutLabel) {VecDestroy(&subv);}
227:       else          {PetscSectionRestoreField_Internal(section, sectionGlobal, gv, f, pStart, pEnd, &is, &subv);}
228:       PetscViewerHDF5PopGroup(viewer);
229:     }
230:   }
231:   DMRestoreGlobalVector(dmBC, &gv);
232:   return(0);
233: }

235: PetscErrorCode VecView_Plex_HDF5_Internal(Vec v, PetscViewer viewer)
236: {
237:   DM             dm;
238:   Vec            locv;
239:   const char    *name;
240:   PetscReal      time;

244:   VecGetDM(v, &dm);
245:   DMGetLocalVector(dm, &locv);
246:   PetscObjectGetName((PetscObject) v, &name);
247:   PetscObjectSetName((PetscObject) locv, name);
248:   DMGlobalToLocalBegin(dm, v, INSERT_VALUES, locv);
249:   DMGlobalToLocalEnd(dm, v, INSERT_VALUES, locv);
250:   DMGetOutputSequenceNumber(dm, NULL, &time);
251:   DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locv, time, NULL, NULL, NULL);
252:   PetscViewerHDF5PushGroup(viewer, "/fields");
253:   PetscViewerPushFormat(viewer, PETSC_VIEWER_HDF5_VIZ);
254:   VecView_Plex_Local_HDF5_Internal(locv, viewer);
255:   PetscViewerPopFormat(viewer);
256:   PetscViewerHDF5PopGroup(viewer);
257:   DMRestoreLocalVector(dm, &locv);
258:   return(0);
259: }

261: PetscErrorCode VecView_Plex_HDF5_Native_Internal(Vec v, PetscViewer viewer)
262: {
263:   PetscBool      isseq;

267:   PetscObjectTypeCompare((PetscObject) v, VECSEQ, &isseq);
268:   PetscViewerHDF5PushGroup(viewer, "/fields");
269:   if (isseq) {VecView_Seq(v, viewer);}
270:   else       {VecView_MPI(v, viewer);}
271:   PetscViewerHDF5PopGroup(viewer);
272:   return(0);
273: }

275: PetscErrorCode VecLoad_Plex_HDF5_Internal(Vec v, PetscViewer viewer)
276: {
277:   DM             dm;
278:   Vec            locv;
279:   const char    *name;
280:   PetscInt       seqnum;

284:   VecGetDM(v, &dm);
285:   DMGetLocalVector(dm, &locv);
286:   PetscObjectGetName((PetscObject) v, &name);
287:   PetscObjectSetName((PetscObject) locv, name);
288:   DMGetOutputSequenceNumber(dm, &seqnum, NULL);
289:   PetscViewerHDF5SetTimestep(viewer, seqnum);
290:   PetscViewerHDF5PushGroup(viewer, "/fields");
291:   VecLoad_Plex_Local(locv, viewer);
292:   PetscViewerHDF5PopGroup(viewer);
293:   DMLocalToGlobalBegin(dm, locv, INSERT_VALUES, v);
294:   DMLocalToGlobalEnd(dm, locv, INSERT_VALUES, v);
295:   DMRestoreLocalVector(dm, &locv);
296:   return(0);
297: }

299: PetscErrorCode VecLoad_Plex_HDF5_Native_Internal(Vec v, PetscViewer viewer)
300: {
301:   DM             dm;
302:   PetscInt       seqnum;

306:   VecGetDM(v, &dm);
307:   DMGetOutputSequenceNumber(dm, &seqnum, NULL);
308:   PetscViewerHDF5SetTimestep(viewer, seqnum);
309:   PetscViewerHDF5PushGroup(viewer, "/fields");
310:   VecLoad_Default(v, viewer);
311:   PetscViewerHDF5PopGroup(viewer);
312:   return(0);
313: }

315: static PetscErrorCode DMPlexWriteTopology_HDF5_Static(DM dm, IS globalPointNumbers, PetscViewer viewer)
316: {
317:   IS              orderIS, conesIS, cellsIS, orntsIS;
318:   const PetscInt *gpoint;
319:   PetscInt       *order, *sizes, *cones, *ornts;
320:   PetscInt        dim, pStart, pEnd, p, conesSize = 0, cellsSize = 0, c = 0, s = 0;
321:   PetscErrorCode  ierr;

324:   ISGetIndices(globalPointNumbers, &gpoint);
325:   DMGetDimension(dm, &dim);
326:   DMPlexGetChart(dm, &pStart, &pEnd);
327:   for (p = pStart; p < pEnd; ++p) {
328:     if (gpoint[p] >= 0) {
329:       PetscInt coneSize;

331:       DMPlexGetConeSize(dm, p, &coneSize);
332:       conesSize += 1;
333:       cellsSize += coneSize;
334:     }
335:   }
336:   PetscMalloc1(conesSize, &order);
337:   PetscMalloc1(conesSize, &sizes);
338:   PetscMalloc1(cellsSize, &cones);
339:   PetscMalloc1(cellsSize, &ornts);
340:   for (p = pStart; p < pEnd; ++p) {
341:     if (gpoint[p] >= 0) {
342:       const PetscInt *cone, *ornt;
343:       PetscInt        coneSize, cp;

345:       DMPlexGetConeSize(dm, p, &coneSize);
346:       DMPlexGetCone(dm, p, &cone);
347:       DMPlexGetConeOrientation(dm, p, &ornt);
348:       order[s]   = gpoint[p];
349:       sizes[s++] = coneSize;
350:       for (cp = 0; cp < coneSize; ++cp, ++c) {cones[c] = gpoint[cone[cp]] < 0 ? -(gpoint[cone[cp]]+1) : gpoint[cone[cp]]; ornts[c] = ornt[cp];}
351:     }
352:   }
353:   if (s != conesSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_LIB, "Total number of points %d != %d", s, conesSize);
354:   if (c != cellsSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_LIB, "Total number of cone points %d != %d", c, cellsSize);
355:   ISCreateGeneral(PetscObjectComm((PetscObject) dm), conesSize, order, PETSC_OWN_POINTER, &orderIS);
356:   PetscObjectSetName((PetscObject) orderIS, "order");
357:   ISCreateGeneral(PetscObjectComm((PetscObject) dm), conesSize, sizes, PETSC_OWN_POINTER, &conesIS);
358:   PetscObjectSetName((PetscObject) conesIS, "cones");
359:   ISCreateGeneral(PetscObjectComm((PetscObject) dm), cellsSize, cones, PETSC_OWN_POINTER, &cellsIS);
360:   PetscObjectSetName((PetscObject) cellsIS, "cells");
361:   ISCreateGeneral(PetscObjectComm((PetscObject) dm), cellsSize, ornts, PETSC_OWN_POINTER, &orntsIS);
362:   PetscObjectSetName((PetscObject) orntsIS, "orientation");
363:   PetscViewerHDF5PushGroup(viewer, "/topology");
364:   ISView(orderIS, viewer);
365:   ISView(conesIS, viewer);
366:   ISView(cellsIS, viewer);
367:   PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) cellsIS, "cell_dim", PETSC_INT, (void *) &dim);
368:   ISView(orntsIS, viewer);
369:   PetscViewerHDF5PopGroup(viewer);
370:   ISDestroy(&orderIS);
371:   ISDestroy(&conesIS);
372:   ISDestroy(&cellsIS);
373:   ISDestroy(&orntsIS);
374:   ISRestoreIndices(globalPointNumbers, &gpoint);
375:   return(0);
376: }

378: static PetscErrorCode CreateConesIS_Private(DM dm, PetscInt cStart, PetscInt cEnd, PetscInt *numCorners, IS *cellIS)
379: {
380:   PetscSF         sfPoint;
381:   DMLabel         cutLabel, cutVertexLabel = NULL;
382:   IS              globalVertexNumbers, cutvertices = NULL;
383:   const PetscInt *gvertex, *cutverts = NULL;
384:   PetscInt       *vertices;
385:   PetscInt        conesSize = 0;
386:   PetscInt        dim, numCornersLocal = 0, cell, vStart, vEnd, vExtra = 0, v;
387:   PetscErrorCode  ierr;

390:   *numCorners = 0;
391:   DMGetDimension(dm, &dim);
392:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);

394:   for (cell = cStart; cell < cEnd; ++cell) {
395:     PetscInt *closure = NULL;
396:     PetscInt  closureSize, v, Nc = 0;

398:     DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);
399:     for (v = 0; v < closureSize*2; v += 2) {
400:       if ((closure[v] >= vStart) && (closure[v] < vEnd)) ++Nc;
401:     }
402:     DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);
403:     conesSize += Nc;
404:     if (!numCornersLocal)           numCornersLocal = Nc;
405:     else if (numCornersLocal != Nc) numCornersLocal = 1;
406:   }
407:   MPIU_Allreduce(&numCornersLocal, numCorners, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject) dm));
408:   if (numCornersLocal && (numCornersLocal != *numCorners || *numCorners == 1)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Visualization topology currently only supports identical cell shapes");
409:   /* Handle periodic cuts by identifying vertices which should be duplicated */
410:   DMGetLabel(dm, "periodic_cut", &cutLabel);
411:   DMPlexCreateCutVertexLabel_Private(dm, cutLabel, &cutVertexLabel);
412:   if (cutVertexLabel) {DMLabelGetStratumIS(cutVertexLabel, 1, &cutvertices);}
413:   if (cutvertices) {
414:     ISGetIndices(cutvertices, &cutverts);
415:     ISGetLocalSize(cutvertices, &vExtra);
416:   }
417:   DMGetPointSF(dm, &sfPoint);
418:   if (cutLabel) {
419:     const PetscInt    *ilocal;
420:     const PetscSFNode *iremote;
421:     PetscInt           nroots, nleaves;

423:     PetscSFGetGraph(sfPoint, &nroots, &nleaves, &ilocal, &iremote);
424:     if (nleaves < 0) {
425:       PetscObjectReference((PetscObject) sfPoint);
426:     } else {
427:       PetscSFCreate(PetscObjectComm((PetscObject) sfPoint), &sfPoint);
428:       PetscSFSetGraph(sfPoint, nroots+vExtra, nleaves, ilocal, PETSC_USE_POINTER, iremote, PETSC_USE_POINTER);
429:     }
430:   } else {
431:     PetscObjectReference((PetscObject) sfPoint);
432:   }
433:   /* Number all vertices */
434:   DMPlexCreateNumbering_Internal(dm, vStart, vEnd+vExtra, 0, NULL, sfPoint, &globalVertexNumbers);
435:   PetscSFDestroy(&sfPoint);
436:   /* Create cones */
437:   ISGetIndices(globalVertexNumbers, &gvertex);
438:   PetscMalloc1(conesSize, &vertices);
439:   for (cell = cStart, v = 0; cell < cEnd; ++cell) {
440:     PetscInt *closure = NULL;
441:     PetscInt  closureSize, Nc = 0, p, value = -1;
442:     PetscBool replace;

444:     if (cutLabel) {DMLabelGetValue(cutLabel, cell, &value);}
445:     replace = (value == 2) ? PETSC_TRUE : PETSC_FALSE;
446:     DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);
447:     for (p = 0; p < closureSize*2; p += 2) {
448:       if ((closure[p] >= vStart) && (closure[p] < vEnd)) {
449:         closure[Nc++] = closure[p];
450:       }
451:     }
452:     DMPlexInvertCell_Internal(dim, Nc, closure);
453:     for (p = 0; p < Nc; ++p) {
454:       PetscInt nv, gv = gvertex[closure[p] - vStart];

456:       if (replace) {
457:         PetscFindInt(closure[p], vExtra, cutverts, &nv);
458:         if (nv >= 0) gv = gvertex[vEnd - vStart + nv];
459:       }
460:       vertices[v++] = gv < 0 ? -(gv+1) : gv;
461:     }
462:     DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);
463:   }
464:   ISRestoreIndices(globalVertexNumbers, &gvertex);
465:   ISDestroy(&globalVertexNumbers);
466:   if (cutvertices) {ISRestoreIndices(cutvertices, &cutverts);}
467:   ISDestroy(&cutvertices);
468:   DMLabelDestroy(&cutVertexLabel);
469:   if (v != conesSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_LIB, "Total number of cell vertices %d != %d", v, conesSize);
470:   ISCreateGeneral(PetscObjectComm((PetscObject) dm), conesSize, vertices, PETSC_OWN_POINTER, cellIS);
471:   PetscLayoutSetBlockSize((*cellIS)->map, *numCorners);
472:   PetscObjectSetName((PetscObject) *cellIS, "cells");
473:   return(0);
474: }

476: static PetscErrorCode DMPlexWriteTopology_Vertices_HDF5_Static(DM dm, PetscViewer viewer)
477: {
478:   DM              cdm;
479:   IS              cellIS = NULL, hcellIS = NULL;
480:   PetscInt        dim, cellHeight, cStart, cMax, cEnd, numCorners, numHCorners;
481:   hid_t           fileId, groupId;
482:   PetscErrorCode  ierr;

485:   DMGetDimension(dm, &dim);
486:   DMGetCoordinateDM(dm, &cdm);
487:   DMPlexGetVTKCellHeight(dm, &cellHeight);
488:   DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);
489:   DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);
490:   cMax = cMax >= 0 ? cMax : cEnd;
491:   CreateConesIS_Private(dm, cStart, cMax, &numCorners,  &cellIS);
492:   if (cMax < cEnd) {CreateConesIS_Private(dm, cMax,   cEnd, &numHCorners, &hcellIS);}

494:   PetscViewerHDF5PushGroup(viewer, "/viz");
495:   PetscViewerHDF5OpenGroup(viewer, &fileId, &groupId);
496:   PetscStackCallHDF5(H5Gclose,(groupId));
497:   PetscViewerHDF5PopGroup(viewer);
498:   if (cellIS) {
499:     PetscViewerHDF5PushGroup(viewer, "/viz/topology");
500:     ISView(cellIS, viewer);
501:     PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) cellIS, "cell_corners", PETSC_INT, (void *) &numCorners);
502:     PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) cellIS, "cell_dim",     PETSC_INT, (void *) &dim);
503:     ISDestroy(&cellIS);
504:     PetscViewerHDF5PopGroup(viewer);
505:   }
506:   if (hcellIS) {
507:     PetscViewerHDF5PushGroup(viewer, "/viz/hybrid_topology");
508:     ISView(hcellIS, viewer);
509:     PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) hcellIS, "cell_corners", PETSC_INT, (void *) &numHCorners);
510:     PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) hcellIS, "cell_dim",     PETSC_INT, (void *) &dim);
511:     ISDestroy(&hcellIS);
512:     PetscViewerHDF5PopGroup(viewer);
513:   }
514:   return(0);
515: }

517: static PetscErrorCode DMPlexWriteCoordinates_HDF5_Static(DM dm, PetscViewer viewer)
518: {
519:   DM             cdm;
520:   Vec            coordinates, newcoords;
521:   PetscReal      lengthScale;
522:   PetscInt       m, M, bs;

526:   DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);
527:   DMGetCoordinateDM(dm, &cdm);
528:   DMGetCoordinates(dm, &coordinates);
529:   VecCreate(PetscObjectComm((PetscObject) coordinates), &newcoords);
530:   PetscObjectSetName((PetscObject) newcoords, "vertices");
531:   VecGetSize(coordinates, &M);
532:   VecGetLocalSize(coordinates, &m);
533:   VecSetSizes(newcoords, m, M);
534:   VecGetBlockSize(coordinates, &bs);
535:   VecSetBlockSize(newcoords, bs);
536:   VecSetType(newcoords,VECSTANDARD);
537:   VecCopy(coordinates, newcoords);
538:   VecScale(newcoords, lengthScale);
539:   /* Did not use DMGetGlobalVector() in order to bypass default group assignment */
540:   PetscViewerHDF5PushGroup(viewer, "/geometry");
541:   PetscViewerPushFormat(viewer, PETSC_VIEWER_NATIVE);
542:   VecView(newcoords, viewer);
543:   PetscViewerPopFormat(viewer);
544:   PetscViewerHDF5PopGroup(viewer);
545:   VecDestroy(&newcoords);
546:   return(0);
547: }

549: static PetscErrorCode DMPlexWriteCoordinates_Vertices_HDF5_Static(DM dm, PetscViewer viewer)
550: {
551:   DM               cdm;
552:   Vec              coordinatesLocal, newcoords;
553:   PetscSection     cSection, cGlobalSection;
554:   PetscScalar     *coords, *ncoords;
555:   DMLabel          cutLabel, cutVertexLabel = NULL;
556:   const PetscReal *L;
557:   const DMBoundaryType *bd;
558:   PetscReal        lengthScale;
559:   PetscInt         vStart, vEnd, v, bs, N, coordSize, dof, off, d;
560:   PetscBool        localized, embedded;
561:   hid_t            fileId, groupId;
562:   PetscErrorCode   ierr;

565:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
566:   DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);
567:   DMGetCoordinatesLocal(dm, &coordinatesLocal);
568:   VecGetBlockSize(coordinatesLocal, &bs);
569:   DMGetCoordinatesLocalized(dm, &localized);
570:   if (localized == PETSC_FALSE) return(0);
571:   DMGetPeriodicity(dm, NULL, NULL, &L, &bd);
572:   DMGetCoordinateDM(dm, &cdm);
573:   DMGetLocalSection(cdm, &cSection);
574:   DMGetGlobalSection(cdm, &cGlobalSection);
575:   DMGetLabel(dm, "periodic_cut", &cutLabel);
576:   N    = 0;

578:   DMPlexCreateCutVertexLabel_Private(dm, cutLabel, &cutVertexLabel);
579:   VecCreate(PetscObjectComm((PetscObject) dm), &newcoords);
580:   PetscSectionGetDof(cSection, vStart, &dof);
581:   PetscPrintf(PETSC_COMM_SELF, "DOF: %D\n", dof);
582:   embedded  = (PetscBool) (L && dof == 2 && !cutLabel);
583:   if (cutVertexLabel) {
584:     DMLabelGetStratumSize(cutVertexLabel, 1, &v);
585:     N   += dof*v;
586:   }
587:   for (v = vStart; v < vEnd; ++v) {
588:     PetscSectionGetDof(cGlobalSection, v, &dof);
589:     if (dof < 0) continue;
590:     if (embedded) N += dof+1;
591:     else          N += dof;
592:   }
593:   if (embedded) {VecSetBlockSize(newcoords, bs+1);}
594:   else          {VecSetBlockSize(newcoords, bs);}
595:   VecSetSizes(newcoords, N, PETSC_DETERMINE);
596:   VecSetType(newcoords, VECSTANDARD);
597:   VecGetArray(coordinatesLocal, &coords);
598:   VecGetArray(newcoords,        &ncoords);
599:   coordSize = 0;
600:   for (v = vStart; v < vEnd; ++v) {
601:     PetscSectionGetDof(cGlobalSection, v, &dof);
602:     PetscSectionGetOffset(cSection, v, &off);
603:     if (dof < 0) continue;
604:     if (embedded) {
605:       if ((bd[0] == DM_BOUNDARY_PERIODIC) && (bd[1] == DM_BOUNDARY_PERIODIC)) {
606:         PetscReal theta, phi, r, R;
607:         /* XY-periodic */
608:         /* Suppose its an y-z circle, then
609:              \hat r = (0, cos(th), sin(th)) \hat x = (1, 0, 0)
610:            and the circle in that plane is
611:              \hat r cos(phi) + \hat x sin(phi) */
612:         theta = 2.0*PETSC_PI*PetscRealPart(coords[off+1])/L[1];
613:         phi   = 2.0*PETSC_PI*PetscRealPart(coords[off+0])/L[0];
614:         r     = L[0]/(2.0*PETSC_PI * 2.0*L[1]);
615:         R     = L[1]/(2.0*PETSC_PI);
616:         ncoords[coordSize++] =  PetscSinReal(phi) * r;
617:         ncoords[coordSize++] = -PetscCosReal(theta) * (R + r * PetscCosReal(phi));
618:         ncoords[coordSize++] =  PetscSinReal(theta) * (R + r * PetscCosReal(phi));
619:       } else if ((bd[0] == DM_BOUNDARY_PERIODIC)) {
620:         /* X-periodic */
621:         ncoords[coordSize++] = -PetscCosReal(2.0*PETSC_PI*PetscRealPart(coords[off+0])/L[0])*(L[0]/(2.0*PETSC_PI));
622:         ncoords[coordSize++] = coords[off+1];
623:         ncoords[coordSize++] = PetscSinReal(2.0*PETSC_PI*PetscRealPart(coords[off+0])/L[0])*(L[0]/(2.0*PETSC_PI));
624:       } else if ((bd[1] == DM_BOUNDARY_PERIODIC)) {
625:         /* Y-periodic */
626:         ncoords[coordSize++] = coords[off+0];
627:         ncoords[coordSize++] = PetscSinReal(2.0*PETSC_PI*PetscRealPart(coords[off+1])/L[1])*(L[1]/(2.0*PETSC_PI));
628:         ncoords[coordSize++] = -PetscCosReal(2.0*PETSC_PI*PetscRealPart(coords[off+1])/L[1])*(L[1]/(2.0*PETSC_PI));
629:       } else if ((bd[0] == DM_BOUNDARY_TWIST)) {
630:         PetscReal phi, r, R;
631:         /* Mobius strip */
632:         /* Suppose its an x-z circle, then
633:              \hat r = (-cos(phi), 0, sin(phi)) \hat y = (0, 1, 0)
634:            and in that plane we rotate by pi as we go around the circle
635:              \hat r cos(phi/2) + \hat y sin(phi/2) */
636:         phi   = 2.0*PETSC_PI*PetscRealPart(coords[off+0])/L[0];
637:         R     = L[0];
638:         r     = PetscRealPart(coords[off+1]) - L[1]/2.0;
639:         ncoords[coordSize++] = -PetscCosReal(phi) * (R + r * PetscCosReal(phi/2.0));
640:         ncoords[coordSize++] =  PetscSinReal(phi/2.0) * r;
641:         ncoords[coordSize++] =  PetscSinReal(phi) * (R + r * PetscCosReal(phi/2.0));
642:       } else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Cannot handle periodicity in this domain");
643:     } else {
644:       for (d = 0; d < dof; ++d, ++coordSize) ncoords[coordSize] = coords[off+d];
645:     }
646:   }
647:   if (cutVertexLabel) {
648:     IS              vertices;
649:     const PetscInt *verts;
650:     PetscInt        n;

652:     DMLabelGetStratumIS(cutVertexLabel, 1, &vertices);
653:     if (vertices) {
654:       ISGetIndices(vertices, &verts);
655:       ISGetLocalSize(vertices, &n);
656:       for (v = 0; v < n; ++v) {
657:         PetscSectionGetDof(cSection, verts[v], &dof);
658:         PetscSectionGetOffset(cSection, verts[v], &off);
659:         for (d = 0; d < dof; ++d) ncoords[coordSize++] = coords[off+d] + ((bd[d] == DM_BOUNDARY_PERIODIC) ? L[d] : 0.0);
660:       }
661:       ISRestoreIndices(vertices, &verts);
662:       ISDestroy(&vertices);
663:     }
664:   }
665:   if (coordSize != N) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Mismatched sizes: %D != %D", coordSize, N);
666:   DMLabelDestroy(&cutVertexLabel);
667:   VecRestoreArray(coordinatesLocal, &coords);
668:   VecRestoreArray(newcoords,        &ncoords);
669:   PetscObjectSetName((PetscObject) newcoords, "vertices");
670:   VecScale(newcoords, lengthScale);
671:   PetscViewerHDF5PushGroup(viewer, "/viz");
672:   PetscViewerHDF5OpenGroup(viewer, &fileId, &groupId);
673:   PetscStackCallHDF5(H5Gclose,(groupId));
674:   PetscViewerHDF5PopGroup(viewer);
675:   PetscViewerHDF5PushGroup(viewer, "/viz/geometry");
676:   VecView(newcoords, viewer);
677:   PetscViewerHDF5PopGroup(viewer);
678:   VecDestroy(&newcoords);
679:   return(0);
680: }


683: static PetscErrorCode DMPlexWriteLabels_HDF5_Static(DM dm, IS globalPointNumbers, PetscViewer viewer)
684: {
685:   const PetscInt   *gpoint;
686:   PetscInt          numLabels, l;
687:   hid_t             fileId, groupId;
688:   PetscErrorCode    ierr;

691:   ISGetIndices(globalPointNumbers, &gpoint);
692:   PetscViewerHDF5PushGroup(viewer, "/labels");
693:   PetscViewerHDF5OpenGroup(viewer, &fileId, &groupId);
694:   if (groupId != fileId) PetscStackCallHDF5(H5Gclose,(groupId));
695:   PetscViewerHDF5PopGroup(viewer);
696:   DMGetNumLabels(dm, &numLabels);
697:   for (l = 0; l < numLabels; ++l) {
698:     DMLabel         label;
699:     const char     *name;
700:     IS              valueIS, pvalueIS, globalValueIS;
701:     const PetscInt *values;
702:     PetscInt        numValues, v;
703:     PetscBool       isDepth, output;
704:     char            group[PETSC_MAX_PATH_LEN];

706:     DMGetLabelName(dm, l, &name);
707:     DMGetLabelOutput(dm, name, &output);
708:     PetscStrncmp(name, "depth", 10, &isDepth);
709:     if (isDepth || !output) continue;
710:     DMGetLabel(dm, name, &label);
711:     PetscSNPrintf(group, PETSC_MAX_PATH_LEN, "/labels/%s", name);
712:     PetscViewerHDF5PushGroup(viewer, group);
713:     PetscViewerHDF5OpenGroup(viewer, &fileId, &groupId);
714:     if (groupId != fileId) PetscStackCallHDF5(H5Gclose,(groupId));
715:     PetscViewerHDF5PopGroup(viewer);
716:     DMLabelGetValueIS(label, &valueIS);
717:     /* Must copy to a new IS on the global comm */
718:     ISGetLocalSize(valueIS, &numValues);
719:     ISGetIndices(valueIS, &values);
720:     ISCreateGeneral(PetscObjectComm((PetscObject) dm), numValues, values, PETSC_COPY_VALUES, &pvalueIS);
721:     ISRestoreIndices(valueIS, &values);
722:     ISAllGather(pvalueIS, &globalValueIS);
723:     ISDestroy(&pvalueIS);
724:     ISSortRemoveDups(globalValueIS);
725:     ISGetLocalSize(globalValueIS, &numValues);
726:     ISGetIndices(globalValueIS, &values);
727:     for (v = 0; v < numValues; ++v) {
728:       IS              stratumIS, globalStratumIS;
729:       const PetscInt *spoints = NULL;
730:       PetscInt       *gspoints, n = 0, gn, p;
731:       const char     *iname = "indices";

733:       PetscSNPrintf(group, PETSC_MAX_PATH_LEN, "/labels/%s/%d", name, values[v]);
734:       DMLabelGetStratumIS(label, values[v], &stratumIS);

736:       if (stratumIS) {ISGetLocalSize(stratumIS, &n);}
737:       if (stratumIS) {ISGetIndices(stratumIS, &spoints);}
738:       for (gn = 0, p = 0; p < n; ++p) if (gpoint[spoints[p]] >= 0) ++gn;
739:       PetscMalloc1(gn,&gspoints);
740:       for (gn = 0, p = 0; p < n; ++p) if (gpoint[spoints[p]] >= 0) gspoints[gn++] = gpoint[spoints[p]];
741:       if (stratumIS) {ISRestoreIndices(stratumIS, &spoints);}
742:       ISCreateGeneral(PetscObjectComm((PetscObject) dm), gn, gspoints, PETSC_OWN_POINTER, &globalStratumIS);
743:       if (stratumIS) {PetscObjectGetName((PetscObject) stratumIS, &iname);}
744:       PetscObjectSetName((PetscObject) globalStratumIS, iname);

746:       PetscViewerHDF5PushGroup(viewer, group);
747:       ISView(globalStratumIS, viewer);
748:       PetscViewerHDF5PopGroup(viewer);
749:       ISDestroy(&globalStratumIS);
750:       ISDestroy(&stratumIS);
751:     }
752:     ISRestoreIndices(globalValueIS, &values);
753:     ISDestroy(&globalValueIS);
754:     ISDestroy(&valueIS);
755:   }
756:   ISRestoreIndices(globalPointNumbers, &gpoint);
757:   return(0);
758: }

760: /* We only write cells and vertices. Does this screw up parallel reading? */
761: PetscErrorCode DMPlexView_HDF5_Internal(DM dm, PetscViewer viewer)
762: {
763:   IS                globalPointNumbers;
764:   PetscViewerFormat format;
765:   PetscBool         viz_geom=PETSC_FALSE, xdmf_topo=PETSC_FALSE, petsc_topo=PETSC_FALSE;
766:   PetscErrorCode    ierr;

769:   DMPlexCreatePointNumbering(dm, &globalPointNumbers);
770:   DMPlexWriteCoordinates_HDF5_Static(dm, viewer);
771:   DMPlexWriteLabels_HDF5_Static(dm, globalPointNumbers, viewer);

773:   PetscViewerGetFormat(viewer, &format);
774:   switch (format) {
775:     case PETSC_VIEWER_HDF5_VIZ:
776:       viz_geom    = PETSC_TRUE;
777:       xdmf_topo   = PETSC_TRUE;
778:       break;
779:     case PETSC_VIEWER_HDF5_XDMF:
780:       xdmf_topo   = PETSC_TRUE;
781:       break;
782:     case PETSC_VIEWER_HDF5_PETSC:
783:       petsc_topo  = PETSC_TRUE;
784:       break;
785:     case PETSC_VIEWER_DEFAULT:
786:     case PETSC_VIEWER_NATIVE:
787:       viz_geom    = PETSC_TRUE;
788:       xdmf_topo   = PETSC_TRUE;
789:       petsc_topo  = PETSC_TRUE;
790:       break;
791:     default:
792:       SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 output.", PetscViewerFormats[format]);
793:   }

795:   if (viz_geom)   {DMPlexWriteCoordinates_Vertices_HDF5_Static(dm, viewer);}
796:   if (xdmf_topo)  {DMPlexWriteTopology_Vertices_HDF5_Static(dm, viewer);}
797:   if (petsc_topo) {DMPlexWriteTopology_HDF5_Static(dm, globalPointNumbers, viewer);}

799:   ISDestroy(&globalPointNumbers);
800:   return(0);
801: }

803: typedef struct {
804:   PetscMPIInt rank;
805:   DM          dm;
806:   PetscViewer viewer;
807:   DMLabel     label;
808: } LabelCtx;

810: static herr_t ReadLabelStratumHDF5_Static(hid_t g_id, const char *name, const H5L_info_t *info, void *op_data)
811: {
812:   PetscViewer     viewer = ((LabelCtx *) op_data)->viewer;
813:   DMLabel         label  = ((LabelCtx *) op_data)->label;
814:   IS              stratumIS;
815:   const PetscInt *ind;
816:   PetscInt        value, N, i;
817:   const char     *lname;
818:   char            group[PETSC_MAX_PATH_LEN];
819:   PetscErrorCode  ierr;

821:   PetscOptionsStringToInt(name, &value);
822:   ISCreate(PetscObjectComm((PetscObject) viewer), &stratumIS);
823:   PetscObjectSetName((PetscObject) stratumIS, "indices");
824:   PetscObjectGetName((PetscObject) label, &lname);
825:   PetscSNPrintf(group, PETSC_MAX_PATH_LEN, "/labels/%s/%s", lname, name);
826:   PetscViewerHDF5PushGroup(viewer, group);
827:   {
828:     /* Force serial load */
829:     PetscViewerHDF5ReadSizes(viewer, "indices", NULL, &N);
830:     PetscLayoutSetLocalSize(stratumIS->map, !((LabelCtx *) op_data)->rank ? N : 0);
831:     PetscLayoutSetSize(stratumIS->map, N);
832:   }
833:   ISLoad(stratumIS, viewer);
834:   PetscViewerHDF5PopGroup(viewer);
835:   ISGetLocalSize(stratumIS, &N);
836:   ISGetIndices(stratumIS, &ind);
837:   for (i = 0; i < N; ++i) {DMLabelSetValue(label, ind[i], value);}
838:   ISRestoreIndices(stratumIS, &ind);
839:   ISDestroy(&stratumIS);
840:   return 0;
841: }

843: static herr_t ReadLabelHDF5_Static(hid_t g_id, const char *name, const H5L_info_t *info, void *op_data)
844: {
845:   DM             dm  = ((LabelCtx *) op_data)->dm;
846:   hsize_t        idx = 0;
848:   herr_t         err;

850:   DMCreateLabel(dm, name); if (ierr) return (herr_t) ierr;
851:   DMGetLabel(dm, name, &((LabelCtx *) op_data)->label); if (ierr) return (herr_t) ierr;
852:   PetscStackCall("H5Literate_by_name",err = H5Literate_by_name(g_id, name, H5_INDEX_NAME, H5_ITER_NATIVE, &idx, ReadLabelStratumHDF5_Static, op_data, 0));
853:   return err;
854: }

856: PetscErrorCode DMPlexLoadLabels_HDF5_Internal(DM dm, PetscViewer viewer)
857: {
858:   LabelCtx        ctx;
859:   hid_t           fileId, groupId;
860:   hsize_t         idx = 0;
861:   PetscErrorCode  ierr;

864:   MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &ctx.rank);
865:   ctx.dm     = dm;
866:   ctx.viewer = viewer;
867:   PetscViewerHDF5PushGroup(viewer, "/labels");
868:   PetscViewerHDF5OpenGroup(viewer, &fileId, &groupId);
869:   PetscStackCallHDF5(H5Literate,(groupId, H5_INDEX_NAME, H5_ITER_NATIVE, &idx, ReadLabelHDF5_Static, &ctx));
870:   PetscStackCallHDF5(H5Gclose,(groupId));
871:   PetscViewerHDF5PopGroup(viewer);
872:   return(0);
873: }

875: /* The first version will read everything onto proc 0, letting the user distribute
876:    The next will create a naive partition, and then rebalance after reading
877: */
878: PetscErrorCode DMPlexLoad_HDF5_Internal(DM dm, PetscViewer viewer)
879: {
880:   PetscSection    coordSection;
881:   Vec             coordinates;
882:   IS              orderIS, conesIS, cellsIS, orntsIS;
883:   const PetscInt *order, *cones, *cells, *ornts;
884:   PetscReal       lengthScale;
885:   PetscInt       *cone, *ornt;
886:   PetscInt        dim, spatialDim, N, numVertices, vStart, vEnd, v, pEnd, p, q, maxConeSize = 0, c;
887:   PetscMPIInt     rank;
888:   PetscErrorCode  ierr;

891:   MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
892:   /* Read toplogy */
893:   PetscViewerHDF5PushGroup(viewer, "/topology");
894:   ISCreate(PetscObjectComm((PetscObject) dm), &orderIS);
895:   PetscObjectSetName((PetscObject) orderIS, "order");
896:   ISCreate(PetscObjectComm((PetscObject) dm), &conesIS);
897:   PetscObjectSetName((PetscObject) conesIS, "cones");
898:   ISCreate(PetscObjectComm((PetscObject) dm), &cellsIS);
899:   PetscObjectSetName((PetscObject) cellsIS, "cells");
900:   ISCreate(PetscObjectComm((PetscObject) dm), &orntsIS);
901:   PetscObjectSetName((PetscObject) orntsIS, "orientation");
902:   PetscViewerHDF5ReadObjectAttribute(viewer, (PetscObject) cellsIS, "cell_dim", PETSC_INT, (void *) &dim);
903:   DMSetDimension(dm, dim);
904:   {
905:     /* Force serial load */
906:     PetscViewerHDF5ReadSizes(viewer, "order", NULL, &pEnd);
907:     PetscLayoutSetLocalSize(orderIS->map, !rank ? pEnd : 0);
908:     PetscLayoutSetSize(orderIS->map, pEnd);
909:     PetscViewerHDF5ReadSizes(viewer, "cones", NULL, &pEnd);
910:     PetscLayoutSetLocalSize(conesIS->map, !rank ? pEnd : 0);
911:     PetscLayoutSetSize(conesIS->map, pEnd);
912:     pEnd = !rank ? pEnd : 0;
913:     PetscViewerHDF5ReadSizes(viewer, "cells", NULL, &N);
914:     PetscLayoutSetLocalSize(cellsIS->map, !rank ? N : 0);
915:     PetscLayoutSetSize(cellsIS->map, N);
916:     PetscViewerHDF5ReadSizes(viewer, "orientation", NULL, &N);
917:     PetscLayoutSetLocalSize(orntsIS->map, !rank ? N : 0);
918:     PetscLayoutSetSize(orntsIS->map, N);
919:   }
920:   ISLoad(orderIS, viewer);
921:   ISLoad(conesIS, viewer);
922:   ISLoad(cellsIS, viewer);
923:   ISLoad(orntsIS, viewer);
924:   PetscViewerHDF5PopGroup(viewer);
925:   /* Read geometry */
926:   PetscViewerHDF5PushGroup(viewer, "/geometry");
927:   VecCreate(PetscObjectComm((PetscObject) dm), &coordinates);
928:   PetscObjectSetName((PetscObject) coordinates, "vertices");
929:   {
930:     /* Force serial load */
931:     PetscViewerHDF5ReadSizes(viewer, "vertices", &spatialDim, &N);
932:     VecSetSizes(coordinates, !rank ? N : 0, N);
933:     VecSetBlockSize(coordinates, spatialDim);
934:   }
935:   VecLoad(coordinates, viewer);
936:   PetscViewerHDF5PopGroup(viewer);
937:   DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);
938:   VecScale(coordinates, 1.0/lengthScale);
939:   VecGetLocalSize(coordinates, &numVertices);
940:   VecGetBlockSize(coordinates, &spatialDim);
941:   numVertices /= spatialDim;
942:   /* Create Plex */
943:   DMPlexSetChart(dm, 0, pEnd);
944:   ISGetIndices(orderIS, &order);
945:   ISGetIndices(conesIS, &cones);
946:   for (p = 0; p < pEnd; ++p) {
947:     DMPlexSetConeSize(dm, order[p], cones[p]);
948:     maxConeSize = PetscMax(maxConeSize, cones[p]);
949:   }
950:   DMSetUp(dm);
951:   ISGetIndices(cellsIS, &cells);
952:   ISGetIndices(orntsIS, &ornts);
953:   PetscMalloc2(maxConeSize,&cone,maxConeSize,&ornt);
954:   for (p = 0, q = 0; p < pEnd; ++p) {
955:     for (c = 0; c < cones[p]; ++c, ++q) {cone[c] = cells[q]; ornt[c] = ornts[q];}
956:     DMPlexSetCone(dm, order[p], cone);
957:     DMPlexSetConeOrientation(dm, order[p], ornt);
958:   }
959:   PetscFree2(cone,ornt);
960:   ISRestoreIndices(orderIS, &order);
961:   ISRestoreIndices(conesIS, &cones);
962:   ISRestoreIndices(cellsIS, &cells);
963:   ISRestoreIndices(orntsIS, &ornts);
964:   ISDestroy(&orderIS);
965:   ISDestroy(&conesIS);
966:   ISDestroy(&cellsIS);
967:   ISDestroy(&orntsIS);
968:   DMPlexSymmetrize(dm);
969:   DMPlexStratify(dm);
970:   /* Create coordinates */
971:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
972:   if (numVertices != vEnd - vStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of coordinates loaded %d does not match number of vertices %d", numVertices, vEnd - vStart);
973:   DMGetCoordinateSection(dm, &coordSection);
974:   PetscSectionSetNumFields(coordSection, 1);
975:   PetscSectionSetFieldComponents(coordSection, 0, spatialDim);
976:   PetscSectionSetChart(coordSection, vStart, vEnd);
977:   for (v = vStart; v < vEnd; ++v) {
978:     PetscSectionSetDof(coordSection, v, spatialDim);
979:     PetscSectionSetFieldDof(coordSection, v, 0, spatialDim);
980:   }
981:   PetscSectionSetUp(coordSection);
982:   DMSetCoordinates(dm, coordinates);
983:   VecDestroy(&coordinates);
984:   /* Read Labels */
985:   DMPlexLoadLabels_HDF5_Internal(dm, viewer);
986:   return(0);
987: }
988: #endif