Actual source code: plexhdf5.c
petsc-master 2020-08-25
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, §ion);
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, §ionGlobal);
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: PetscObject isZero;
240: const char *name;
241: PetscReal time;
245: VecGetDM(v, &dm);
246: DMGetLocalVector(dm, &locv);
247: PetscObjectGetName((PetscObject) v, &name);
248: PetscObjectSetName((PetscObject) locv, name);
249: PetscObjectQuery((PetscObject) v, "__Vec_bc_zero__", &isZero);
250: PetscObjectCompose((PetscObject) locv, "__Vec_bc_zero__", isZero);
251: DMGlobalToLocalBegin(dm, v, INSERT_VALUES, locv);
252: DMGlobalToLocalEnd(dm, v, INSERT_VALUES, locv);
253: DMGetOutputSequenceNumber(dm, NULL, &time);
254: DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locv, time, NULL, NULL, NULL);
255: PetscViewerHDF5PushGroup(viewer, "/fields");
256: PetscViewerPushFormat(viewer, PETSC_VIEWER_HDF5_VIZ);
257: VecView_Plex_Local_HDF5_Internal(locv, viewer);
258: PetscViewerPopFormat(viewer);
259: PetscViewerHDF5PopGroup(viewer);
260: PetscObjectCompose((PetscObject) locv, "__Vec_bc_zero__", NULL);
261: DMRestoreLocalVector(dm, &locv);
262: return(0);
263: }
265: PetscErrorCode VecView_Plex_HDF5_Native_Internal(Vec v, PetscViewer viewer)
266: {
267: PetscBool isseq;
271: PetscObjectTypeCompare((PetscObject) v, VECSEQ, &isseq);
272: PetscViewerHDF5PushGroup(viewer, "/fields");
273: if (isseq) {VecView_Seq(v, viewer);}
274: else {VecView_MPI(v, viewer);}
275: PetscViewerHDF5PopGroup(viewer);
276: return(0);
277: }
279: PetscErrorCode VecLoad_Plex_HDF5_Internal(Vec v, PetscViewer viewer)
280: {
281: DM dm;
282: Vec locv;
283: const char *name;
284: PetscInt seqnum;
288: VecGetDM(v, &dm);
289: DMGetLocalVector(dm, &locv);
290: PetscObjectGetName((PetscObject) v, &name);
291: PetscObjectSetName((PetscObject) locv, name);
292: DMGetOutputSequenceNumber(dm, &seqnum, NULL);
293: PetscViewerHDF5SetTimestep(viewer, seqnum);
294: PetscViewerHDF5PushGroup(viewer, "/fields");
295: VecLoad_Plex_Local(locv, viewer);
296: PetscViewerHDF5PopGroup(viewer);
297: DMLocalToGlobalBegin(dm, locv, INSERT_VALUES, v);
298: DMLocalToGlobalEnd(dm, locv, INSERT_VALUES, v);
299: DMRestoreLocalVector(dm, &locv);
300: return(0);
301: }
303: PetscErrorCode VecLoad_Plex_HDF5_Native_Internal(Vec v, PetscViewer viewer)
304: {
305: DM dm;
306: PetscInt seqnum;
310: VecGetDM(v, &dm);
311: DMGetOutputSequenceNumber(dm, &seqnum, NULL);
312: PetscViewerHDF5SetTimestep(viewer, seqnum);
313: PetscViewerHDF5PushGroup(viewer, "/fields");
314: VecLoad_Default(v, viewer);
315: PetscViewerHDF5PopGroup(viewer);
316: return(0);
317: }
319: static PetscErrorCode DMPlexWriteTopology_HDF5_Static(DM dm, IS globalPointNumbers, PetscViewer viewer)
320: {
321: IS orderIS, conesIS, cellsIS, orntsIS;
322: const PetscInt *gpoint;
323: PetscInt *order, *sizes, *cones, *ornts;
324: PetscInt dim, pStart, pEnd, p, conesSize = 0, cellsSize = 0, c = 0, s = 0;
325: PetscErrorCode ierr;
328: ISGetIndices(globalPointNumbers, &gpoint);
329: DMGetDimension(dm, &dim);
330: DMPlexGetChart(dm, &pStart, &pEnd);
331: for (p = pStart; p < pEnd; ++p) {
332: if (gpoint[p] >= 0) {
333: PetscInt coneSize;
335: DMPlexGetConeSize(dm, p, &coneSize);
336: conesSize += 1;
337: cellsSize += coneSize;
338: }
339: }
340: PetscMalloc1(conesSize, &order);
341: PetscMalloc1(conesSize, &sizes);
342: PetscMalloc1(cellsSize, &cones);
343: PetscMalloc1(cellsSize, &ornts);
344: for (p = pStart; p < pEnd; ++p) {
345: if (gpoint[p] >= 0) {
346: const PetscInt *cone, *ornt;
347: PetscInt coneSize, cp;
349: DMPlexGetConeSize(dm, p, &coneSize);
350: DMPlexGetCone(dm, p, &cone);
351: DMPlexGetConeOrientation(dm, p, &ornt);
352: order[s] = gpoint[p];
353: sizes[s++] = coneSize;
354: for (cp = 0; cp < coneSize; ++cp, ++c) {cones[c] = gpoint[cone[cp]] < 0 ? -(gpoint[cone[cp]]+1) : gpoint[cone[cp]]; ornts[c] = ornt[cp];}
355: }
356: }
357: if (s != conesSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_LIB, "Total number of points %d != %d", s, conesSize);
358: if (c != cellsSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_LIB, "Total number of cone points %d != %d", c, cellsSize);
359: ISCreateGeneral(PetscObjectComm((PetscObject) dm), conesSize, order, PETSC_OWN_POINTER, &orderIS);
360: PetscObjectSetName((PetscObject) orderIS, "order");
361: ISCreateGeneral(PetscObjectComm((PetscObject) dm), conesSize, sizes, PETSC_OWN_POINTER, &conesIS);
362: PetscObjectSetName((PetscObject) conesIS, "cones");
363: ISCreateGeneral(PetscObjectComm((PetscObject) dm), cellsSize, cones, PETSC_OWN_POINTER, &cellsIS);
364: PetscObjectSetName((PetscObject) cellsIS, "cells");
365: ISCreateGeneral(PetscObjectComm((PetscObject) dm), cellsSize, ornts, PETSC_OWN_POINTER, &orntsIS);
366: PetscObjectSetName((PetscObject) orntsIS, "orientation");
367: PetscViewerHDF5PushGroup(viewer, "/topology");
368: ISView(orderIS, viewer);
369: ISView(conesIS, viewer);
370: ISView(cellsIS, viewer);
371: PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) cellsIS, "cell_dim", PETSC_INT, (void *) &dim);
372: ISView(orntsIS, viewer);
373: PetscViewerHDF5PopGroup(viewer);
374: ISDestroy(&orderIS);
375: ISDestroy(&conesIS);
376: ISDestroy(&cellsIS);
377: ISDestroy(&orntsIS);
378: ISRestoreIndices(globalPointNumbers, &gpoint);
379: return(0);
380: }
382: static PetscErrorCode CreateConesIS_Private(DM dm, PetscInt cStart, PetscInt cEnd, PetscInt *numCorners, IS *cellIS)
383: {
384: PetscSF sfPoint;
385: DMLabel cutLabel, cutVertexLabel = NULL;
386: IS globalVertexNumbers, cutvertices = NULL;
387: const PetscInt *gvertex, *cutverts = NULL;
388: PetscInt *vertices;
389: PetscInt conesSize = 0;
390: PetscInt dim, numCornersLocal = 0, cell, vStart, vEnd, vExtra = 0, v;
391: PetscErrorCode ierr;
394: *numCorners = 0;
395: DMGetDimension(dm, &dim);
396: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
398: for (cell = cStart; cell < cEnd; ++cell) {
399: PetscInt *closure = NULL;
400: PetscInt closureSize, v, Nc = 0;
402: DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);
403: for (v = 0; v < closureSize*2; v += 2) {
404: if ((closure[v] >= vStart) && (closure[v] < vEnd)) ++Nc;
405: }
406: DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);
407: conesSize += Nc;
408: if (!numCornersLocal) numCornersLocal = Nc;
409: else if (numCornersLocal != Nc) numCornersLocal = 1;
410: }
411: MPIU_Allreduce(&numCornersLocal, numCorners, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject) dm));
412: if (numCornersLocal && (numCornersLocal != *numCorners || *numCorners == 1)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Visualization topology currently only supports identical cell shapes");
413: /* Handle periodic cuts by identifying vertices which should be duplicated */
414: DMGetLabel(dm, "periodic_cut", &cutLabel);
415: DMPlexCreateCutVertexLabel_Private(dm, cutLabel, &cutVertexLabel);
416: if (cutVertexLabel) {DMLabelGetStratumIS(cutVertexLabel, 1, &cutvertices);}
417: if (cutvertices) {
418: ISGetIndices(cutvertices, &cutverts);
419: ISGetLocalSize(cutvertices, &vExtra);
420: }
421: DMGetPointSF(dm, &sfPoint);
422: if (cutLabel) {
423: const PetscInt *ilocal;
424: const PetscSFNode *iremote;
425: PetscInt nroots, nleaves;
427: PetscSFGetGraph(sfPoint, &nroots, &nleaves, &ilocal, &iremote);
428: if (nleaves < 0) {
429: PetscObjectReference((PetscObject) sfPoint);
430: } else {
431: PetscSFCreate(PetscObjectComm((PetscObject) sfPoint), &sfPoint);
432: PetscSFSetGraph(sfPoint, nroots+vExtra, nleaves, ilocal, PETSC_USE_POINTER, iremote, PETSC_USE_POINTER);
433: }
434: } else {
435: PetscObjectReference((PetscObject) sfPoint);
436: }
437: /* Number all vertices */
438: DMPlexCreateNumbering_Plex(dm, vStart, vEnd+vExtra, 0, NULL, sfPoint, &globalVertexNumbers);
439: PetscSFDestroy(&sfPoint);
440: /* Create cones */
441: ISGetIndices(globalVertexNumbers, &gvertex);
442: PetscMalloc1(conesSize, &vertices);
443: for (cell = cStart, v = 0; cell < cEnd; ++cell) {
444: PetscInt *closure = NULL;
445: PetscInt closureSize, Nc = 0, p, value = -1;
446: PetscBool replace;
448: if (cutLabel) {DMLabelGetValue(cutLabel, cell, &value);}
449: replace = (value == 2) ? PETSC_TRUE : PETSC_FALSE;
450: DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);
451: for (p = 0; p < closureSize*2; p += 2) {
452: if ((closure[p] >= vStart) && (closure[p] < vEnd)) {
453: closure[Nc++] = closure[p];
454: }
455: }
456: DMPlexReorderCell(dm, cell, closure);
457: for (p = 0; p < Nc; ++p) {
458: PetscInt nv, gv = gvertex[closure[p] - vStart];
460: if (replace) {
461: PetscFindInt(closure[p], vExtra, cutverts, &nv);
462: if (nv >= 0) gv = gvertex[vEnd - vStart + nv];
463: }
464: vertices[v++] = gv < 0 ? -(gv+1) : gv;
465: }
466: DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);
467: }
468: ISRestoreIndices(globalVertexNumbers, &gvertex);
469: ISDestroy(&globalVertexNumbers);
470: if (cutvertices) {ISRestoreIndices(cutvertices, &cutverts);}
471: ISDestroy(&cutvertices);
472: DMLabelDestroy(&cutVertexLabel);
473: if (v != conesSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_LIB, "Total number of cell vertices %d != %d", v, conesSize);
474: ISCreateGeneral(PetscObjectComm((PetscObject) dm), conesSize, vertices, PETSC_OWN_POINTER, cellIS);
475: PetscLayoutSetBlockSize((*cellIS)->map, *numCorners);
476: PetscObjectSetName((PetscObject) *cellIS, "cells");
477: return(0);
478: }
480: static PetscErrorCode DMPlexWriteTopology_Vertices_HDF5_Static(DM dm, PetscViewer viewer)
481: {
482: DM cdm;
483: DMLabel depthLabel, ctLabel;
484: IS cellIS;
485: PetscInt dim, depth, cellHeight, c;
486: hid_t fileId, groupId;
487: PetscErrorCode ierr;
490: PetscViewerHDF5PushGroup(viewer, "/viz");
491: PetscViewerHDF5OpenGroup(viewer, &fileId, &groupId);
492: PetscStackCallHDF5(H5Gclose,(groupId));
494: PetscViewerHDF5PopGroup(viewer);
495: DMGetDimension(dm, &dim);
496: DMPlexGetDepth(dm, &depth);
497: DMGetCoordinateDM(dm, &cdm);
498: DMPlexGetVTKCellHeight(dm, &cellHeight);
499: DMPlexGetDepthLabel(dm, &depthLabel);
500: DMPlexGetCellTypeLabel(dm, &ctLabel);
501: for (c = 0; c < DM_NUM_POLYTOPES; ++c) {
502: const DMPolytopeType ict = (DMPolytopeType) c;
503: PetscInt pStart, pEnd, dep, numCorners, n = 0;
504: PetscBool output = PETSC_FALSE, doOutput;
506: DMLabelGetStratumBounds(ctLabel, ict, &pStart, &pEnd);
507: if (pStart >= 0) {
508: DMLabelGetValue(depthLabel, pStart, &dep);
509: if (dep == depth - cellHeight) output = PETSC_TRUE;
510: }
511: MPI_Allreduce(&output, &doOutput, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject) dm));
512: if (!doOutput) continue;
513: CreateConesIS_Private(dm, pStart, pEnd, &numCorners, &cellIS);
514: if (!n) {
515: PetscViewerHDF5PushGroup(viewer, "/viz/topology");
516: } else {
517: char group[PETSC_MAX_PATH_LEN];
519: PetscSNPrintf(group, PETSC_MAX_PATH_LEN, "/viz/topology_%D", n);
520: PetscViewerHDF5PushGroup(viewer, group);
521: }
522: ISView(cellIS, viewer);
523: PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) cellIS, "cell_corners", PETSC_INT, (void *) &numCorners);
524: PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) cellIS, "cell_dim", PETSC_INT, (void *) &dim);
525: ISDestroy(&cellIS);
526: PetscViewerHDF5PopGroup(viewer);
527: ++n;
528: }
529: return(0);
530: }
532: static PetscErrorCode DMPlexWriteCoordinates_HDF5_Static(DM dm, PetscViewer viewer)
533: {
534: DM cdm;
535: Vec coordinates, newcoords;
536: PetscReal lengthScale;
537: PetscInt m, M, bs;
541: DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);
542: DMGetCoordinateDM(dm, &cdm);
543: DMGetCoordinates(dm, &coordinates);
544: VecCreate(PetscObjectComm((PetscObject) coordinates), &newcoords);
545: PetscObjectSetName((PetscObject) newcoords, "vertices");
546: VecGetSize(coordinates, &M);
547: VecGetLocalSize(coordinates, &m);
548: VecSetSizes(newcoords, m, M);
549: VecGetBlockSize(coordinates, &bs);
550: VecSetBlockSize(newcoords, bs);
551: VecSetType(newcoords,VECSTANDARD);
552: VecCopy(coordinates, newcoords);
553: VecScale(newcoords, lengthScale);
554: /* Did not use DMGetGlobalVector() in order to bypass default group assignment */
555: PetscViewerHDF5PushGroup(viewer, "/geometry");
556: PetscViewerPushFormat(viewer, PETSC_VIEWER_NATIVE);
557: VecView(newcoords, viewer);
558: PetscViewerPopFormat(viewer);
559: PetscViewerHDF5PopGroup(viewer);
560: VecDestroy(&newcoords);
561: return(0);
562: }
564: static PetscErrorCode DMPlexWriteCoordinates_Vertices_HDF5_Static(DM dm, PetscViewer viewer)
565: {
566: DM cdm;
567: Vec coordinatesLocal, newcoords;
568: PetscSection cSection, cGlobalSection;
569: PetscScalar *coords, *ncoords;
570: DMLabel cutLabel, cutVertexLabel = NULL;
571: const PetscReal *L;
572: const DMBoundaryType *bd;
573: PetscReal lengthScale;
574: PetscInt vStart, vEnd, v, bs, N, coordSize, dof, off, d;
575: PetscBool localized, embedded;
576: hid_t fileId, groupId;
577: PetscErrorCode ierr;
580: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
581: DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);
582: DMGetCoordinatesLocal(dm, &coordinatesLocal);
583: VecGetBlockSize(coordinatesLocal, &bs);
584: DMGetCoordinatesLocalized(dm, &localized);
585: if (localized == PETSC_FALSE) return(0);
586: DMGetPeriodicity(dm, NULL, NULL, &L, &bd);
587: DMGetCoordinateDM(dm, &cdm);
588: DMGetLocalSection(cdm, &cSection);
589: DMGetGlobalSection(cdm, &cGlobalSection);
590: DMGetLabel(dm, "periodic_cut", &cutLabel);
591: N = 0;
593: DMPlexCreateCutVertexLabel_Private(dm, cutLabel, &cutVertexLabel);
594: VecCreate(PetscObjectComm((PetscObject) dm), &newcoords);
595: PetscSectionGetDof(cSection, vStart, &dof);
596: PetscPrintf(PETSC_COMM_SELF, "DOF: %D\n", dof);
597: embedded = (PetscBool) (L && dof == 2 && !cutLabel);
598: if (cutVertexLabel) {
599: DMLabelGetStratumSize(cutVertexLabel, 1, &v);
600: N += dof*v;
601: }
602: for (v = vStart; v < vEnd; ++v) {
603: PetscSectionGetDof(cGlobalSection, v, &dof);
604: if (dof < 0) continue;
605: if (embedded) N += dof+1;
606: else N += dof;
607: }
608: if (embedded) {VecSetBlockSize(newcoords, bs+1);}
609: else {VecSetBlockSize(newcoords, bs);}
610: VecSetSizes(newcoords, N, PETSC_DETERMINE);
611: VecSetType(newcoords, VECSTANDARD);
612: VecGetArray(coordinatesLocal, &coords);
613: VecGetArray(newcoords, &ncoords);
614: coordSize = 0;
615: for (v = vStart; v < vEnd; ++v) {
616: PetscSectionGetDof(cGlobalSection, v, &dof);
617: PetscSectionGetOffset(cSection, v, &off);
618: if (dof < 0) continue;
619: if (embedded) {
620: if ((bd[0] == DM_BOUNDARY_PERIODIC) && (bd[1] == DM_BOUNDARY_PERIODIC)) {
621: PetscReal theta, phi, r, R;
622: /* XY-periodic */
623: /* Suppose its an y-z circle, then
624: \hat r = (0, cos(th), sin(th)) \hat x = (1, 0, 0)
625: and the circle in that plane is
626: \hat r cos(phi) + \hat x sin(phi) */
627: theta = 2.0*PETSC_PI*PetscRealPart(coords[off+1])/L[1];
628: phi = 2.0*PETSC_PI*PetscRealPart(coords[off+0])/L[0];
629: r = L[0]/(2.0*PETSC_PI * 2.0*L[1]);
630: R = L[1]/(2.0*PETSC_PI);
631: ncoords[coordSize++] = PetscSinReal(phi) * r;
632: ncoords[coordSize++] = -PetscCosReal(theta) * (R + r * PetscCosReal(phi));
633: ncoords[coordSize++] = PetscSinReal(theta) * (R + r * PetscCosReal(phi));
634: } else if ((bd[0] == DM_BOUNDARY_PERIODIC)) {
635: /* X-periodic */
636: ncoords[coordSize++] = -PetscCosReal(2.0*PETSC_PI*PetscRealPart(coords[off+0])/L[0])*(L[0]/(2.0*PETSC_PI));
637: ncoords[coordSize++] = coords[off+1];
638: ncoords[coordSize++] = PetscSinReal(2.0*PETSC_PI*PetscRealPart(coords[off+0])/L[0])*(L[0]/(2.0*PETSC_PI));
639: } else if ((bd[1] == DM_BOUNDARY_PERIODIC)) {
640: /* Y-periodic */
641: ncoords[coordSize++] = coords[off+0];
642: ncoords[coordSize++] = PetscSinReal(2.0*PETSC_PI*PetscRealPart(coords[off+1])/L[1])*(L[1]/(2.0*PETSC_PI));
643: ncoords[coordSize++] = -PetscCosReal(2.0*PETSC_PI*PetscRealPart(coords[off+1])/L[1])*(L[1]/(2.0*PETSC_PI));
644: } else if ((bd[0] == DM_BOUNDARY_TWIST)) {
645: PetscReal phi, r, R;
646: /* Mobius strip */
647: /* Suppose its an x-z circle, then
648: \hat r = (-cos(phi), 0, sin(phi)) \hat y = (0, 1, 0)
649: and in that plane we rotate by pi as we go around the circle
650: \hat r cos(phi/2) + \hat y sin(phi/2) */
651: phi = 2.0*PETSC_PI*PetscRealPart(coords[off+0])/L[0];
652: R = L[0];
653: r = PetscRealPart(coords[off+1]) - L[1]/2.0;
654: ncoords[coordSize++] = -PetscCosReal(phi) * (R + r * PetscCosReal(phi/2.0));
655: ncoords[coordSize++] = PetscSinReal(phi/2.0) * r;
656: ncoords[coordSize++] = PetscSinReal(phi) * (R + r * PetscCosReal(phi/2.0));
657: } else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Cannot handle periodicity in this domain");
658: } else {
659: for (d = 0; d < dof; ++d, ++coordSize) ncoords[coordSize] = coords[off+d];
660: }
661: }
662: if (cutVertexLabel) {
663: IS vertices;
664: const PetscInt *verts;
665: PetscInt n;
667: DMLabelGetStratumIS(cutVertexLabel, 1, &vertices);
668: if (vertices) {
669: ISGetIndices(vertices, &verts);
670: ISGetLocalSize(vertices, &n);
671: for (v = 0; v < n; ++v) {
672: PetscSectionGetDof(cSection, verts[v], &dof);
673: PetscSectionGetOffset(cSection, verts[v], &off);
674: for (d = 0; d < dof; ++d) ncoords[coordSize++] = coords[off+d] + ((bd[d] == DM_BOUNDARY_PERIODIC) ? L[d] : 0.0);
675: }
676: ISRestoreIndices(vertices, &verts);
677: ISDestroy(&vertices);
678: }
679: }
680: if (coordSize != N) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Mismatched sizes: %D != %D", coordSize, N);
681: DMLabelDestroy(&cutVertexLabel);
682: VecRestoreArray(coordinatesLocal, &coords);
683: VecRestoreArray(newcoords, &ncoords);
684: PetscObjectSetName((PetscObject) newcoords, "vertices");
685: VecScale(newcoords, lengthScale);
686: PetscViewerHDF5PushGroup(viewer, "/viz");
687: PetscViewerHDF5OpenGroup(viewer, &fileId, &groupId);
688: PetscStackCallHDF5(H5Gclose,(groupId));
689: PetscViewerHDF5PopGroup(viewer);
690: PetscViewerHDF5PushGroup(viewer, "/viz/geometry");
691: VecView(newcoords, viewer);
692: PetscViewerHDF5PopGroup(viewer);
693: VecDestroy(&newcoords);
694: return(0);
695: }
698: static PetscErrorCode DMPlexWriteLabels_HDF5_Static(DM dm, IS globalPointNumbers, PetscViewer viewer)
699: {
700: const PetscInt *gpoint;
701: PetscInt numLabels, l;
702: hid_t fileId, groupId;
703: PetscErrorCode ierr;
706: ISGetIndices(globalPointNumbers, &gpoint);
707: PetscViewerHDF5PushGroup(viewer, "/labels");
708: PetscViewerHDF5OpenGroup(viewer, &fileId, &groupId);
709: if (groupId != fileId) PetscStackCallHDF5(H5Gclose,(groupId));
710: PetscViewerHDF5PopGroup(viewer);
711: DMGetNumLabels(dm, &numLabels);
712: for (l = 0; l < numLabels; ++l) {
713: DMLabel label;
714: const char *name;
715: IS valueIS, pvalueIS, globalValueIS;
716: const PetscInt *values;
717: PetscInt numValues, v;
718: PetscBool isDepth, output;
719: char group[PETSC_MAX_PATH_LEN];
721: DMGetLabelName(dm, l, &name);
722: DMGetLabelOutput(dm, name, &output);
723: PetscStrncmp(name, "depth", 10, &isDepth);
724: if (isDepth || !output) continue;
725: DMGetLabel(dm, name, &label);
726: PetscSNPrintf(group, PETSC_MAX_PATH_LEN, "/labels/%s", name);
727: PetscViewerHDF5PushGroup(viewer, group);
728: PetscViewerHDF5OpenGroup(viewer, &fileId, &groupId);
729: if (groupId != fileId) PetscStackCallHDF5(H5Gclose,(groupId));
730: PetscViewerHDF5PopGroup(viewer);
731: DMLabelGetValueIS(label, &valueIS);
732: /* Must copy to a new IS on the global comm */
733: ISGetLocalSize(valueIS, &numValues);
734: ISGetIndices(valueIS, &values);
735: ISCreateGeneral(PetscObjectComm((PetscObject) dm), numValues, values, PETSC_COPY_VALUES, &pvalueIS);
736: ISRestoreIndices(valueIS, &values);
737: ISAllGather(pvalueIS, &globalValueIS);
738: ISDestroy(&pvalueIS);
739: ISSortRemoveDups(globalValueIS);
740: ISGetLocalSize(globalValueIS, &numValues);
741: ISGetIndices(globalValueIS, &values);
742: for (v = 0; v < numValues; ++v) {
743: IS stratumIS, globalStratumIS;
744: const PetscInt *spoints = NULL;
745: PetscInt *gspoints, n = 0, gn, p;
746: const char *iname = "indices";
748: PetscSNPrintf(group, PETSC_MAX_PATH_LEN, "/labels/%s/%d", name, values[v]);
749: DMLabelGetStratumIS(label, values[v], &stratumIS);
751: if (stratumIS) {ISGetLocalSize(stratumIS, &n);}
752: if (stratumIS) {ISGetIndices(stratumIS, &spoints);}
753: for (gn = 0, p = 0; p < n; ++p) if (gpoint[spoints[p]] >= 0) ++gn;
754: PetscMalloc1(gn,&gspoints);
755: for (gn = 0, p = 0; p < n; ++p) if (gpoint[spoints[p]] >= 0) gspoints[gn++] = gpoint[spoints[p]];
756: if (stratumIS) {ISRestoreIndices(stratumIS, &spoints);}
757: ISCreateGeneral(PetscObjectComm((PetscObject) dm), gn, gspoints, PETSC_OWN_POINTER, &globalStratumIS);
758: if (stratumIS) {PetscObjectGetName((PetscObject) stratumIS, &iname);}
759: PetscObjectSetName((PetscObject) globalStratumIS, iname);
761: PetscViewerHDF5PushGroup(viewer, group);
762: ISView(globalStratumIS, viewer);
763: PetscViewerHDF5PopGroup(viewer);
764: ISDestroy(&globalStratumIS);
765: ISDestroy(&stratumIS);
766: }
767: ISRestoreIndices(globalValueIS, &values);
768: ISDestroy(&globalValueIS);
769: ISDestroy(&valueIS);
770: }
771: ISRestoreIndices(globalPointNumbers, &gpoint);
772: return(0);
773: }
775: /* We only write cells and vertices. Does this screw up parallel reading? */
776: PetscErrorCode DMPlexView_HDF5_Internal(DM dm, PetscViewer viewer)
777: {
778: IS globalPointNumbers;
779: PetscViewerFormat format;
780: PetscBool viz_geom=PETSC_FALSE, xdmf_topo=PETSC_FALSE, petsc_topo=PETSC_FALSE;
781: PetscErrorCode ierr;
784: DMPlexCreatePointNumbering(dm, &globalPointNumbers);
785: DMPlexWriteCoordinates_HDF5_Static(dm, viewer);
786: DMPlexWriteLabels_HDF5_Static(dm, globalPointNumbers, viewer);
788: PetscViewerGetFormat(viewer, &format);
789: switch (format) {
790: case PETSC_VIEWER_HDF5_VIZ:
791: viz_geom = PETSC_TRUE;
792: xdmf_topo = PETSC_TRUE;
793: break;
794: case PETSC_VIEWER_HDF5_XDMF:
795: xdmf_topo = PETSC_TRUE;
796: break;
797: case PETSC_VIEWER_HDF5_PETSC:
798: petsc_topo = PETSC_TRUE;
799: break;
800: case PETSC_VIEWER_DEFAULT:
801: case PETSC_VIEWER_NATIVE:
802: viz_geom = PETSC_TRUE;
803: xdmf_topo = PETSC_TRUE;
804: petsc_topo = PETSC_TRUE;
805: break;
806: default:
807: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 output.", PetscViewerFormats[format]);
808: }
810: if (viz_geom) {DMPlexWriteCoordinates_Vertices_HDF5_Static(dm, viewer);}
811: if (xdmf_topo) {DMPlexWriteTopology_Vertices_HDF5_Static(dm, viewer);}
812: if (petsc_topo) {DMPlexWriteTopology_HDF5_Static(dm, globalPointNumbers, viewer);}
814: ISDestroy(&globalPointNumbers);
815: return(0);
816: }
818: typedef struct {
819: PetscMPIInt rank;
820: DM dm;
821: PetscViewer viewer;
822: DMLabel label;
823: } LabelCtx;
825: static herr_t ReadLabelStratumHDF5_Static(hid_t g_id, const char *name, const H5L_info_t *info, void *op_data)
826: {
827: PetscViewer viewer = ((LabelCtx *) op_data)->viewer;
828: DMLabel label = ((LabelCtx *) op_data)->label;
829: IS stratumIS;
830: const PetscInt *ind;
831: PetscInt value, N, i;
832: const char *lname;
833: char group[PETSC_MAX_PATH_LEN];
834: PetscErrorCode ierr;
836: PetscOptionsStringToInt(name, &value);
837: ISCreate(PetscObjectComm((PetscObject) viewer), &stratumIS);
838: PetscObjectSetName((PetscObject) stratumIS, "indices");
839: PetscObjectGetName((PetscObject) label, &lname);
840: PetscSNPrintf(group, PETSC_MAX_PATH_LEN, "/labels/%s/%s", lname, name);
841: PetscViewerHDF5PushGroup(viewer, group);
842: {
843: /* Force serial load */
844: PetscViewerHDF5ReadSizes(viewer, "indices", NULL, &N);
845: PetscLayoutSetLocalSize(stratumIS->map, !((LabelCtx *) op_data)->rank ? N : 0);
846: PetscLayoutSetSize(stratumIS->map, N);
847: }
848: ISLoad(stratumIS, viewer);
849: PetscViewerHDF5PopGroup(viewer);
850: ISGetLocalSize(stratumIS, &N);
851: ISGetIndices(stratumIS, &ind);
852: for (i = 0; i < N; ++i) {DMLabelSetValue(label, ind[i], value);}
853: ISRestoreIndices(stratumIS, &ind);
854: ISDestroy(&stratumIS);
855: return 0;
856: }
858: static herr_t ReadLabelHDF5_Static(hid_t g_id, const char *name, const H5L_info_t *info, void *op_data)
859: {
860: DM dm = ((LabelCtx *) op_data)->dm;
861: hsize_t idx = 0;
863: herr_t err;
865: DMCreateLabel(dm, name); if (ierr) return (herr_t) ierr;
866: DMGetLabel(dm, name, &((LabelCtx *) op_data)->label); if (ierr) return (herr_t) ierr;
867: PetscStackCall("H5Literate_by_name",err = H5Literate_by_name(g_id, name, H5_INDEX_NAME, H5_ITER_NATIVE, &idx, ReadLabelStratumHDF5_Static, op_data, 0));
868: return err;
869: }
871: PetscErrorCode DMPlexLoadLabels_HDF5_Internal(DM dm, PetscViewer viewer)
872: {
873: LabelCtx ctx;
874: hid_t fileId, groupId;
875: hsize_t idx = 0;
876: PetscErrorCode ierr;
879: MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &ctx.rank);
880: ctx.dm = dm;
881: ctx.viewer = viewer;
882: PetscViewerHDF5PushGroup(viewer, "/labels");
883: PetscViewerHDF5OpenGroup(viewer, &fileId, &groupId);
884: PetscStackCallHDF5(H5Literate,(groupId, H5_INDEX_NAME, H5_ITER_NATIVE, &idx, ReadLabelHDF5_Static, &ctx));
885: PetscStackCallHDF5(H5Gclose,(groupId));
886: PetscViewerHDF5PopGroup(viewer);
887: return(0);
888: }
890: /* The first version will read everything onto proc 0, letting the user distribute
891: The next will create a naive partition, and then rebalance after reading
892: */
893: PetscErrorCode DMPlexLoad_HDF5_Internal(DM dm, PetscViewer viewer)
894: {
895: PetscSection coordSection;
896: Vec coordinates;
897: IS orderIS, conesIS, cellsIS, orntsIS;
898: const PetscInt *order, *cones, *cells, *ornts;
899: PetscReal lengthScale;
900: PetscInt *cone, *ornt;
901: PetscInt dim, spatialDim, N, numVertices, vStart, vEnd, v, pEnd, p, q, maxConeSize = 0, c;
902: PetscMPIInt rank;
903: PetscErrorCode ierr;
906: MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
907: /* Read toplogy */
908: PetscViewerHDF5PushGroup(viewer, "/topology");
909: ISCreate(PetscObjectComm((PetscObject) dm), &orderIS);
910: PetscObjectSetName((PetscObject) orderIS, "order");
911: ISCreate(PetscObjectComm((PetscObject) dm), &conesIS);
912: PetscObjectSetName((PetscObject) conesIS, "cones");
913: ISCreate(PetscObjectComm((PetscObject) dm), &cellsIS);
914: PetscObjectSetName((PetscObject) cellsIS, "cells");
915: ISCreate(PetscObjectComm((PetscObject) dm), &orntsIS);
916: PetscObjectSetName((PetscObject) orntsIS, "orientation");
917: PetscViewerHDF5ReadObjectAttribute(viewer, (PetscObject) cellsIS, "cell_dim", PETSC_INT, (void *) &dim);
918: DMSetDimension(dm, dim);
919: {
920: /* Force serial load */
921: PetscViewerHDF5ReadSizes(viewer, "order", NULL, &pEnd);
922: PetscLayoutSetLocalSize(orderIS->map, !rank ? pEnd : 0);
923: PetscLayoutSetSize(orderIS->map, pEnd);
924: PetscViewerHDF5ReadSizes(viewer, "cones", NULL, &pEnd);
925: PetscLayoutSetLocalSize(conesIS->map, !rank ? pEnd : 0);
926: PetscLayoutSetSize(conesIS->map, pEnd);
927: pEnd = !rank ? pEnd : 0;
928: PetscViewerHDF5ReadSizes(viewer, "cells", NULL, &N);
929: PetscLayoutSetLocalSize(cellsIS->map, !rank ? N : 0);
930: PetscLayoutSetSize(cellsIS->map, N);
931: PetscViewerHDF5ReadSizes(viewer, "orientation", NULL, &N);
932: PetscLayoutSetLocalSize(orntsIS->map, !rank ? N : 0);
933: PetscLayoutSetSize(orntsIS->map, N);
934: }
935: ISLoad(orderIS, viewer);
936: ISLoad(conesIS, viewer);
937: ISLoad(cellsIS, viewer);
938: ISLoad(orntsIS, viewer);
939: PetscViewerHDF5PopGroup(viewer);
940: /* Read geometry */
941: PetscViewerHDF5PushGroup(viewer, "/geometry");
942: VecCreate(PetscObjectComm((PetscObject) dm), &coordinates);
943: PetscObjectSetName((PetscObject) coordinates, "vertices");
944: {
945: /* Force serial load */
946: PetscViewerHDF5ReadSizes(viewer, "vertices", &spatialDim, &N);
947: VecSetSizes(coordinates, !rank ? N : 0, N);
948: VecSetBlockSize(coordinates, spatialDim);
949: }
950: VecLoad(coordinates, viewer);
951: PetscViewerHDF5PopGroup(viewer);
952: DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);
953: VecScale(coordinates, 1.0/lengthScale);
954: VecGetLocalSize(coordinates, &numVertices);
955: VecGetBlockSize(coordinates, &spatialDim);
956: numVertices /= spatialDim;
957: /* Create Plex */
958: DMPlexSetChart(dm, 0, pEnd);
959: ISGetIndices(orderIS, &order);
960: ISGetIndices(conesIS, &cones);
961: for (p = 0; p < pEnd; ++p) {
962: DMPlexSetConeSize(dm, order[p], cones[p]);
963: maxConeSize = PetscMax(maxConeSize, cones[p]);
964: }
965: DMSetUp(dm);
966: ISGetIndices(cellsIS, &cells);
967: ISGetIndices(orntsIS, &ornts);
968: PetscMalloc2(maxConeSize,&cone,maxConeSize,&ornt);
969: for (p = 0, q = 0; p < pEnd; ++p) {
970: for (c = 0; c < cones[p]; ++c, ++q) {cone[c] = cells[q]; ornt[c] = ornts[q];}
971: DMPlexSetCone(dm, order[p], cone);
972: DMPlexSetConeOrientation(dm, order[p], ornt);
973: }
974: PetscFree2(cone,ornt);
975: ISRestoreIndices(orderIS, &order);
976: ISRestoreIndices(conesIS, &cones);
977: ISRestoreIndices(cellsIS, &cells);
978: ISRestoreIndices(orntsIS, &ornts);
979: ISDestroy(&orderIS);
980: ISDestroy(&conesIS);
981: ISDestroy(&cellsIS);
982: ISDestroy(&orntsIS);
983: DMPlexSymmetrize(dm);
984: DMPlexStratify(dm);
985: /* Create coordinates */
986: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
987: 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);
988: DMGetCoordinateSection(dm, &coordSection);
989: PetscSectionSetNumFields(coordSection, 1);
990: PetscSectionSetFieldComponents(coordSection, 0, spatialDim);
991: PetscSectionSetChart(coordSection, vStart, vEnd);
992: for (v = vStart; v < vEnd; ++v) {
993: PetscSectionSetDof(coordSection, v, spatialDim);
994: PetscSectionSetFieldDof(coordSection, v, 0, spatialDim);
995: }
996: PetscSectionSetUp(coordSection);
997: DMSetCoordinates(dm, coordinates);
998: VecDestroy(&coordinates);
999: /* Read Labels */
1000: DMPlexLoadLabels_HDF5_Internal(dm, viewer);
1001: return(0);
1002: }
1003: #endif