15cca0b87SVaclav Hapla static char help[] = "Load and save the mesh and fields to HDF5 and ExodusII\n\n";
25cca0b87SVaclav Hapla
3bd412c90SVaclav Hapla #include <petsc/private/dmpleximpl.h>
45cca0b87SVaclav Hapla #include <petscviewerhdf5.h>
55cca0b87SVaclav Hapla #include <petscsf.h>
65cca0b87SVaclav Hapla
75cca0b87SVaclav Hapla typedef struct {
85cca0b87SVaclav Hapla PetscBool compare; /* Compare the meshes using DMPlexEqual() */
938883f1bSVaclav Hapla PetscBool compare_labels; /* Compare labels in the meshes using DMCompareLabels() */
105cca0b87SVaclav Hapla PetscBool distribute; /* Distribute the mesh */
11bcde9acaSMatthew G. Knepley PetscBool field; /* Layout a field over the mesh */
12bcde9acaSMatthew G. Knepley PetscBool reorder; /* Reorder the points in the section */
13bd412c90SVaclav Hapla char ofname[PETSC_MAX_PATH_LEN]; /* Output mesh filename */
145cca0b87SVaclav Hapla PetscViewerFormat format; /* Format to write and read */
155cca0b87SVaclav Hapla PetscBool second_write_read; /* Write and read for the 2nd time */
16806c51deSksagiyam PetscBool use_low_level_functions; /* Use low level functions for viewing and loading */
175cca0b87SVaclav Hapla } AppCtx;
185cca0b87SVaclav Hapla
ProcessOptions(MPI_Comm comm,AppCtx * options)19d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
20d71ae5a4SJacob Faibussowitsch {
215cca0b87SVaclav Hapla PetscFunctionBeginUser;
225cca0b87SVaclav Hapla options->compare = PETSC_FALSE;
2338883f1bSVaclav Hapla options->compare_labels = PETSC_FALSE;
245cca0b87SVaclav Hapla options->distribute = PETSC_TRUE;
25bcde9acaSMatthew G. Knepley options->field = PETSC_FALSE;
26bcde9acaSMatthew G. Knepley options->reorder = PETSC_FALSE;
275cca0b87SVaclav Hapla options->format = PETSC_VIEWER_DEFAULT;
285cca0b87SVaclav Hapla options->second_write_read = PETSC_FALSE;
29806c51deSksagiyam options->use_low_level_functions = PETSC_FALSE;
30c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(options->ofname, "ex55.h5", sizeof(options->ofname)));
315cca0b87SVaclav Hapla
32d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
339566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-compare", "Compare the meshes using DMPlexEqual()", "ex55.c", options->compare, &options->compare, NULL));
349566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-compare_labels", "Compare labels in the meshes using DMCompareLabels()", "ex55.c", options->compare_labels, &options->compare_labels, NULL));
359566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-distribute", "Distribute the mesh", "ex55.c", options->distribute, &options->distribute, NULL));
36bcde9acaSMatthew G. Knepley PetscCall(PetscOptionsBool("-field", "Layout a field over the mesh", "ex55.c", options->field, &options->field, NULL));
37bcde9acaSMatthew G. Knepley PetscCall(PetscOptionsBool("-reorder", "Reorder the points in the section", "ex55.c", options->reorder, &options->reorder, NULL));
38bd412c90SVaclav Hapla PetscCall(PetscOptionsString("-ofilename", "The output mesh file", "ex55.c", options->ofname, options->ofname, sizeof(options->ofname), NULL));
399566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-format", "Format to write and read", "ex55.c", PetscViewerFormats, (PetscEnum)options->format, (PetscEnum *)&options->format, NULL));
409566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-second_write_read", "Write and read for the 2nd time", "ex55.c", options->second_write_read, &options->second_write_read, NULL));
419566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_low_level_functions", "Use low level functions for viewing and loading", "ex55.c", options->use_low_level_functions, &options->use_low_level_functions, NULL));
42d0609cedSBarry Smith PetscOptionsEnd();
433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
44bd412c90SVaclav Hapla }
455cca0b87SVaclav Hapla
CreateMesh(MPI_Comm comm,DM * dm)46bcde9acaSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm)
47bcde9acaSMatthew G. Knepley {
48bcde9acaSMatthew G. Knepley PetscFunctionBeginUser;
49bcde9acaSMatthew G. Knepley PetscCall(DMCreate(comm, dm));
50bcde9acaSMatthew G. Knepley PetscCall(DMSetType(*dm, DMPLEX));
51bcde9acaSMatthew G. Knepley PetscCall(DMSetOptionsPrefix(*dm, "orig_"));
52bcde9acaSMatthew G. Knepley PetscCall(DMSetFromOptions(*dm));
53bcde9acaSMatthew G. Knepley PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
54bcde9acaSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
55bcde9acaSMatthew G. Knepley }
56bcde9acaSMatthew G. Knepley
CreateDiscretization(DM dm)57bcde9acaSMatthew G. Knepley static PetscErrorCode CreateDiscretization(DM dm)
58bcde9acaSMatthew G. Knepley {
59bcde9acaSMatthew G. Knepley PetscFE fe;
60bcde9acaSMatthew G. Knepley DMPolytopeType ct;
61bcde9acaSMatthew G. Knepley PetscInt dim, cStart;
62bcde9acaSMatthew G. Knepley
63bcde9acaSMatthew G. Knepley PetscFunctionBeginUser;
64bcde9acaSMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim));
65bcde9acaSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
66bcde9acaSMatthew G. Knepley PetscCall(DMPlexGetCellType(dm, cStart, &ct));
67bcde9acaSMatthew G. Knepley PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, 1, ct, 1, PETSC_DETERMINE, &fe));
68bcde9acaSMatthew G. Knepley PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
69bcde9acaSMatthew G. Knepley PetscCall(PetscFEDestroy(&fe));
70bcde9acaSMatthew G. Knepley PetscCall(DMCreateDS(dm));
71bcde9acaSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
72bcde9acaSMatthew G. Knepley }
73bcde9acaSMatthew G. Knepley
CheckDistributed(DM dm,PetscBool expectedDistributed)74bd412c90SVaclav Hapla static PetscErrorCode CheckDistributed(DM dm, PetscBool expectedDistributed)
75bd412c90SVaclav Hapla {
76bd412c90SVaclav Hapla PetscMPIInt size;
77bd412c90SVaclav Hapla PetscBool distributed;
78bd412c90SVaclav Hapla const char YES[] = "DISTRIBUTED";
79bd412c90SVaclav Hapla const char NO[] = "NOT DISTRIBUTED";
80bd412c90SVaclav Hapla
81bd412c90SVaclav Hapla PetscFunctionBeginUser;
82bd412c90SVaclav Hapla PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
833ba16761SJacob Faibussowitsch if (size < 2) PetscFunctionReturn(PETSC_SUCCESS);
84bd412c90SVaclav Hapla PetscCall(DMPlexIsDistributed(dm, &distributed));
85bd412c90SVaclav Hapla PetscCheck(distributed == expectedDistributed, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Expected DM being %s but actually is %s", expectedDistributed ? YES : NO, distributed ? YES : NO);
863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
87bd412c90SVaclav Hapla }
88bd412c90SVaclav Hapla
CheckInterpolated(DM dm,PetscBool expectedInterpolated)89bd412c90SVaclav Hapla static PetscErrorCode CheckInterpolated(DM dm, PetscBool expectedInterpolated)
90bd412c90SVaclav Hapla {
91bd412c90SVaclav Hapla DMPlexInterpolatedFlag iflg;
92bd412c90SVaclav Hapla PetscBool interpolated;
93bd412c90SVaclav Hapla const char YES[] = "INTERPOLATED";
94bd412c90SVaclav Hapla const char NO[] = "NOT INTERPOLATED";
95bd412c90SVaclav Hapla
96bd412c90SVaclav Hapla PetscFunctionBeginUser;
97bd412c90SVaclav Hapla PetscCall(DMPlexIsInterpolatedCollective(dm, &iflg));
98bd412c90SVaclav Hapla interpolated = (PetscBool)(iflg == DMPLEX_INTERPOLATED_FULL);
99bd412c90SVaclav Hapla PetscCheck(interpolated == expectedInterpolated, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Expected DM being %s but actually is %s", expectedInterpolated ? YES : NO, interpolated ? YES : NO);
1003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
101bd412c90SVaclav Hapla }
102bd412c90SVaclav Hapla
CheckDistributedInterpolated(DM dm,PetscBool expectedInterpolated,PetscViewer v,AppCtx * user)103bcde9acaSMatthew G. Knepley static PetscErrorCode CheckDistributedInterpolated(DM dm, PetscBool expectedInterpolated, PetscViewer v, AppCtx *user)
104bd412c90SVaclav Hapla {
105bd412c90SVaclav Hapla PetscViewerFormat format;
106bcde9acaSMatthew G. Knepley PetscBool distributed, interpolated = expectedInterpolated;
107bd412c90SVaclav Hapla
108bd412c90SVaclav Hapla PetscFunctionBeginUser;
109bd412c90SVaclav Hapla PetscCall(PetscViewerGetFormat(v, &format));
110bd412c90SVaclav Hapla switch (format) {
111bd412c90SVaclav Hapla case PETSC_VIEWER_HDF5_XDMF:
112bd412c90SVaclav Hapla case PETSC_VIEWER_HDF5_VIZ: {
113bd412c90SVaclav Hapla distributed = PETSC_TRUE;
114bd412c90SVaclav Hapla interpolated = PETSC_FALSE;
115bd412c90SVaclav Hapla }; break;
116bd412c90SVaclav Hapla case PETSC_VIEWER_HDF5_PETSC:
117bd412c90SVaclav Hapla case PETSC_VIEWER_DEFAULT:
118bd412c90SVaclav Hapla case PETSC_VIEWER_NATIVE: {
119bd412c90SVaclav Hapla DMPlexStorageVersion version;
120bd412c90SVaclav Hapla
121bd412c90SVaclav Hapla PetscCall(PetscViewerHDF5GetDMPlexStorageVersionReading(v, &version));
122bd412c90SVaclav Hapla distributed = (PetscBool)(version->major >= 3);
123bd412c90SVaclav Hapla }; break;
124bcde9acaSMatthew G. Knepley default:
125bd412c90SVaclav Hapla distributed = PETSC_FALSE;
126bd412c90SVaclav Hapla }
127bd412c90SVaclav Hapla PetscCall(CheckDistributed(dm, distributed));
128bd412c90SVaclav Hapla PetscCall(CheckInterpolated(dm, interpolated));
1293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
130bd412c90SVaclav Hapla }
131bd412c90SVaclav Hapla
DMPlexWriteAndReadHDF5(DM dm,Vec vec,const char filename[],const char prefix[],PetscBool expectedInterpolated,AppCtx * user,DM * dm_new,Vec * v_new)132bcde9acaSMatthew G. Knepley static PetscErrorCode DMPlexWriteAndReadHDF5(DM dm, Vec vec, const char filename[], const char prefix[], PetscBool expectedInterpolated, AppCtx *user, DM *dm_new, Vec *v_new)
133d71ae5a4SJacob Faibussowitsch {
1345cca0b87SVaclav Hapla DM dmnew;
135bcde9acaSMatthew G. Knepley Vec vnew = NULL;
136ed5e4e85SVaclav Hapla const char savedName[] = "Mesh";
137ed5e4e85SVaclav Hapla const char loadedName[] = "Mesh_new";
1385cca0b87SVaclav Hapla PetscViewer v;
1395cca0b87SVaclav Hapla
1405cca0b87SVaclav Hapla PetscFunctionBeginUser;
1419566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Open(PetscObjectComm((PetscObject)dm), filename, FILE_MODE_WRITE, &v));
142bd412c90SVaclav Hapla PetscCall(PetscViewerPushFormat(v, user->format));
1439566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dm, savedName));
144bd412c90SVaclav Hapla if (user->use_low_level_functions) {
1459566063dSJacob Faibussowitsch PetscCall(DMPlexTopologyView(dm, v));
1469566063dSJacob Faibussowitsch PetscCall(DMPlexCoordinatesView(dm, v));
1479566063dSJacob Faibussowitsch PetscCall(DMPlexLabelsView(dm, v));
148806c51deSksagiyam } else {
1499566063dSJacob Faibussowitsch PetscCall(DMView(dm, v));
150bcde9acaSMatthew G. Knepley if (vec) PetscCall(VecView(vec, v));
151806c51deSksagiyam }
1529566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(v, FILE_MODE_READ));
1539566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dmnew));
1549566063dSJacob Faibussowitsch PetscCall(DMSetType(dmnew, DMPLEX));
155bd412c90SVaclav Hapla PetscCall(DMPlexDistributeSetDefault(dmnew, PETSC_FALSE));
1569566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dmnew, savedName));
1579566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(dmnew, prefix));
158bd412c90SVaclav Hapla if (user->use_low_level_functions) {
159c9ad657eSksagiyam PetscSF sfXC;
160c9ad657eSksagiyam
1619566063dSJacob Faibussowitsch PetscCall(DMPlexTopologyLoad(dmnew, v, &sfXC));
1629566063dSJacob Faibussowitsch PetscCall(DMPlexCoordinatesLoad(dmnew, v, sfXC));
1639566063dSJacob Faibussowitsch PetscCall(DMPlexLabelsLoad(dmnew, v, sfXC));
1649566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfXC));
165806c51deSksagiyam } else {
1669566063dSJacob Faibussowitsch PetscCall(DMLoad(dmnew, v));
167bcde9acaSMatthew G. Knepley if (vec) {
168bcde9acaSMatthew G. Knepley PetscCall(CreateDiscretization(dmnew));
169bcde9acaSMatthew G. Knepley PetscCall(DMCreateGlobalVector(dmnew, &vnew));
170bcde9acaSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)vnew, "solution"));
171bcde9acaSMatthew G. Knepley PetscCall(VecLoad(vnew, v));
172806c51deSksagiyam }
173bcde9acaSMatthew G. Knepley }
17407dc484cSMatthew G. Knepley DMLabel celltypes;
17507dc484cSMatthew G. Knepley PetscCall(DMPlexGetCellTypeLabel(dmnew, &celltypes));
176bcde9acaSMatthew G. Knepley PetscCall(CheckDistributedInterpolated(dmnew, expectedInterpolated, v, user));
1779566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dmnew, loadedName));
1789566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(v));
1799566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&v));
1805cca0b87SVaclav Hapla *dm_new = dmnew;
181bcde9acaSMatthew G. Knepley *v_new = vnew;
1823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1835cca0b87SVaclav Hapla }
1845cca0b87SVaclav Hapla
main(int argc,char ** argv)185d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
186d71ae5a4SJacob Faibussowitsch {
1875cca0b87SVaclav Hapla DM dm, dmnew;
188bcde9acaSMatthew G. Knepley Vec v = NULL, vnew = NULL;
1895cca0b87SVaclav Hapla PetscPartitioner part;
1905cca0b87SVaclav Hapla AppCtx user;
191bcde9acaSMatthew G. Knepley PetscBool interpolated = PETSC_TRUE, flg;
1925cca0b87SVaclav Hapla
193327415f7SBarry Smith PetscFunctionBeginUser;
1949566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1959566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
196bcde9acaSMatthew G. Knepley PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm));
197bcde9acaSMatthew G. Knepley PetscCall(PetscOptionsGetBool(NULL, "orig_", "-dm_plex_interpolate", &interpolated, NULL));
198bcde9acaSMatthew G. Knepley PetscCall(CheckInterpolated(dm, interpolated));
1995cca0b87SVaclav Hapla
2005cca0b87SVaclav Hapla if (user.distribute) {
2015cca0b87SVaclav Hapla DM dmdist;
2025cca0b87SVaclav Hapla
2039566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(dm, &part));
204bd412c90SVaclav Hapla PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSIMPLE));
2059566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part));
2069566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(dm, 0, NULL, &dmdist));
2075cca0b87SVaclav Hapla if (dmdist) {
2089566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
2095cca0b87SVaclav Hapla dm = dmdist;
210bd412c90SVaclav Hapla PetscCall(CheckDistributed(dm, PETSC_TRUE));
211bcde9acaSMatthew G. Knepley PetscCall(CheckInterpolated(dm, interpolated));
2125cca0b87SVaclav Hapla }
2135cca0b87SVaclav Hapla }
214bcde9acaSMatthew G. Knepley if (user.field) {
215bcde9acaSMatthew G. Knepley PetscSection gs;
216bcde9acaSMatthew G. Knepley PetscScalar *a;
217bcde9acaSMatthew G. Knepley PetscInt pStart, pEnd, rStart;
218bcde9acaSMatthew G. Knepley
219bcde9acaSMatthew G. Knepley PetscCall(CreateDiscretization(dm));
220bcde9acaSMatthew G. Knepley
221bcde9acaSMatthew G. Knepley PetscCall(DMCreateGlobalVector(dm, &v));
222bcde9acaSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)v, "solution"));
223bcde9acaSMatthew G. Knepley PetscCall(DMGetGlobalSection(dm, &gs));
224bcde9acaSMatthew G. Knepley PetscCall(PetscSectionGetChart(gs, &pStart, &pEnd));
225bcde9acaSMatthew G. Knepley PetscCall(VecGetOwnershipRange(v, &rStart, NULL));
226bcde9acaSMatthew G. Knepley PetscCall(VecGetArrayWrite(v, &a));
227bcde9acaSMatthew G. Knepley for (PetscInt p = pStart; p < pEnd; ++p) {
228bcde9acaSMatthew G. Knepley PetscInt dof, off;
229bcde9acaSMatthew G. Knepley
230bcde9acaSMatthew G. Knepley PetscCall(PetscSectionGetDof(gs, p, &dof));
231bcde9acaSMatthew G. Knepley PetscCall(PetscSectionGetOffset(gs, p, &off));
232bcde9acaSMatthew G. Knepley if (off < 0) continue;
233bcde9acaSMatthew G. Knepley for (PetscInt d = 0; d < dof; ++d) a[off + d] = p;
234bcde9acaSMatthew G. Knepley }
235bcde9acaSMatthew G. Knepley PetscCall(VecRestoreArrayWrite(v, &a));
236bcde9acaSMatthew G. Knepley }
2375cca0b87SVaclav Hapla
2389566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(dm, NULL));
2399566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE));
2409566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm));
2419566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
2425cca0b87SVaclav Hapla
243bcde9acaSMatthew G. Knepley PetscCall(DMPlexWriteAndReadHDF5(dm, v, user.ofname, "new_", interpolated, &user, &dmnew, &vnew));
2445cca0b87SVaclav Hapla
2455cca0b87SVaclav Hapla if (user.second_write_read) {
2469566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
2475cca0b87SVaclav Hapla dm = dmnew;
248bcde9acaSMatthew G. Knepley PetscCall(DMPlexWriteAndReadHDF5(dm, v, user.ofname, "new_", interpolated, &user, &dmnew, &vnew));
2495cca0b87SVaclav Hapla }
2505cca0b87SVaclav Hapla
2519566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dmnew, NULL, "-dm_view"));
2525cca0b87SVaclav Hapla
2535cca0b87SVaclav Hapla /* This currently makes sense only for sequential meshes. */
2545cca0b87SVaclav Hapla if (user.compare) {
2559566063dSJacob Faibussowitsch PetscCall(DMPlexEqual(dmnew, dm, &flg));
25628b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "DMs are not equal");
2579566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "DMs equal\n"));
258bcde9acaSMatthew G. Knepley if (v) {
259bcde9acaSMatthew G. Knepley PetscCall(VecEqual(vnew, v, &flg));
260bcde9acaSMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Vecs %s\n", flg ? "equal" : "are not equal"));
261bcde9acaSMatthew G. Knepley PetscCall(VecViewFromOptions(v, NULL, "-old_vec_view"));
262bcde9acaSMatthew G. Knepley PetscCall(VecViewFromOptions(vnew, NULL, "-new_vec_view"));
263bcde9acaSMatthew G. Knepley }
26438883f1bSVaclav Hapla }
26538883f1bSVaclav Hapla if (user.compare_labels) {
2669566063dSJacob Faibussowitsch PetscCall(DMCompareLabels(dmnew, dm, NULL, NULL));
2679566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "DMLabels equal\n"));
2685cca0b87SVaclav Hapla }
2695cca0b87SVaclav Hapla
270bcde9acaSMatthew G. Knepley PetscCall(VecDestroy(&v));
2719566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
272bcde9acaSMatthew G. Knepley PetscCall(VecDestroy(&vnew));
2739566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmnew));
2749566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
275b122ec5aSJacob Faibussowitsch return 0;
2765cca0b87SVaclav Hapla }
2775cca0b87SVaclav Hapla
2785cca0b87SVaclav Hapla /*TEST
2795cca0b87SVaclav Hapla build:
2805cca0b87SVaclav Hapla requires: hdf5
2815cca0b87SVaclav Hapla # Idempotence of saving/loading
2825cca0b87SVaclav Hapla # Have to replace Exodus file, which is creating uninterpolated edges
2835cca0b87SVaclav Hapla test:
2845cca0b87SVaclav Hapla suffix: 0
28538883f1bSVaclav Hapla TODO: broken
28638883f1bSVaclav Hapla requires: exodusii
2875cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail
2885cca0b87SVaclav Hapla args: -format hdf5_petsc -compare
2895cca0b87SVaclav Hapla test:
2905cca0b87SVaclav Hapla suffix: 1
29138883f1bSVaclav Hapla TODO: broken
29238883f1bSVaclav Hapla requires: exodusii parmetis !defined(PETSC_USE_64BIT_INDICES)
2935cca0b87SVaclav Hapla nsize: 2
2945cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail
2955cca0b87SVaclav Hapla args: -petscpartitioner_type parmetis
2965cca0b87SVaclav Hapla args: -format hdf5_petsc -new_dm_view ascii::ascii_info_detail
297bd412c90SVaclav Hapla
2985cca0b87SVaclav Hapla testset:
2995cca0b87SVaclav Hapla requires: exodusii
300bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo
3013886731fSPierre Jolivet output_file: output/empty.out
3025cca0b87SVaclav Hapla test:
3035cca0b87SVaclav Hapla suffix: 2
3045cca0b87SVaclav Hapla nsize: {{1 2 4 8}separate output}
3055cca0b87SVaclav Hapla args: -format {{default hdf5_petsc}separate output}
3068886822aSVaclav Hapla args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0 3.0.0}}
307bcde9acaSMatthew G. Knepley args: -orig_dm_plex_interpolate {{0 1}separate output}
3085cca0b87SVaclav Hapla test:
3095cca0b87SVaclav Hapla suffix: 2a
3105cca0b87SVaclav Hapla nsize: {{1 2 4 8}separate output}
3115cca0b87SVaclav Hapla args: -format {{hdf5_xdmf hdf5_viz}separate output}
312bd412c90SVaclav Hapla
3135cca0b87SVaclav Hapla test:
3145cca0b87SVaclav Hapla suffix: 3
3155cca0b87SVaclav Hapla requires: exodusii
316bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -compare -compare_labels
3178886822aSVaclav Hapla args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0 3.0.0}}
3185cca0b87SVaclav Hapla
3195cca0b87SVaclav Hapla # Load HDF5 file in XDMF format in parallel, write, read dm1, write, read dm2, and compare dm1 and dm2
32038883f1bSVaclav Hapla testset:
3215cca0b87SVaclav Hapla suffix: 4
3225cca0b87SVaclav Hapla requires: !complex
323bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf
32438883f1bSVaclav Hapla args: -distribute 0 -second_write_read -compare
32538883f1bSVaclav Hapla test:
32638883f1bSVaclav Hapla suffix: hdf5_petsc
32738883f1bSVaclav Hapla nsize: {{1 2}}
32838883f1bSVaclav Hapla args: -format hdf5_petsc -compare_labels
3298886822aSVaclav Hapla args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0 3.0.0}}
33038883f1bSVaclav Hapla test:
33138883f1bSVaclav Hapla suffix: hdf5_xdmf
33238883f1bSVaclav Hapla nsize: {{1 3 8}}
33338883f1bSVaclav Hapla args: -format hdf5_xdmf
3345cca0b87SVaclav Hapla
335806c51deSksagiyam # Use low level functions, DMPlexTopologyView()/Load(), DMPlexCoordinatesView()/Load(), and DMPlexLabelsView()/Load()
336bd412c90SVaclav Hapla # TODO: The output is very long so keeping just 1.0.0 version. This test should be redesigned or removed.
337806c51deSksagiyam test:
338806c51deSksagiyam suffix: 5
339806c51deSksagiyam requires: exodusii
340806c51deSksagiyam nsize: 2
341bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo
342bcde9acaSMatthew G. Knepley args: -orig_dm_plex_interpolate 0 -orig_dm_distribute 0
343806c51deSksagiyam args: -dm_view ascii::ascii_info_detail
344806c51deSksagiyam args: -new_dm_view ascii::ascii_info_detail
345bd412c90SVaclav Hapla args: -format hdf5_petsc -use_low_level_functions {{0 1}}
346bd412c90SVaclav Hapla args: -dm_plex_view_hdf5_storage_version 1.0.0
347806c51deSksagiyam
3485cca0b87SVaclav Hapla testset:
34938883f1bSVaclav Hapla suffix: 6
35038883f1bSVaclav Hapla requires: hdf5 !complex datafilespath
35138883f1bSVaclav Hapla nsize: {{1 3}}
35238883f1bSVaclav Hapla args: -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_geometry
353bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${DATAFILESPATH}/meshes/cube-hexahedra-refined.h5 -dm_plex_create_from_hdf5_xdmf -dm_plex_hdf5_topology_path /cells -dm_plex_hdf5_geometry_path /coordinates
354bcde9acaSMatthew G. Knepley args: -orig_dm_distribute 0
35538883f1bSVaclav Hapla args: -format hdf5_petsc -second_write_read -compare -compare_labels
356bcde9acaSMatthew G. Knepley args: -orig_dm_plex_interpolate {{0 1}} -distribute {{0 1}}
3578886822aSVaclav Hapla args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0 3.0.0}}
35838883f1bSVaclav Hapla
35938883f1bSVaclav Hapla testset:
3605cca0b87SVaclav Hapla # the same data and settings as dm_impls_plex_tests-ex18_9%
36138883f1bSVaclav Hapla suffix: 9
3625cca0b87SVaclav Hapla requires: hdf5 !complex datafilespath
3634d056bb6SVaclav Hapla nsize: {{1 2 4}}
3645cca0b87SVaclav Hapla args: -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_geometry
365bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${DATAFILESPATH}/meshes/cube-hexahedra-refined.h5 -dm_plex_create_from_hdf5_xdmf -dm_plex_hdf5_topology_path /cells -dm_plex_hdf5_geometry_path /coordinates
366bcde9acaSMatthew G. Knepley args: -orig_dm_distribute 0
367bd412c90SVaclav Hapla args: -format {{hdf5_petsc hdf5_xdmf}} -second_write_read -compare
3688886822aSVaclav Hapla args: -dm_plex_view_hdf5_storage_version 3.0.0
3695cca0b87SVaclav Hapla test:
37038883f1bSVaclav Hapla suffix: hdf5_seqload
371bd412c90SVaclav Hapla args: -distribute
372bcde9acaSMatthew G. Knepley args: -orig_dm_plex_interpolate {{0 1}}
3735cca0b87SVaclav Hapla args: -dm_plex_hdf5_force_sequential
3745cca0b87SVaclav Hapla test:
37538883f1bSVaclav Hapla suffix: hdf5_seqload_metis
3765cca0b87SVaclav Hapla requires: parmetis
3775cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type parmetis
378bcde9acaSMatthew G. Knepley args: -orig_dm_plex_interpolate 1
3795cca0b87SVaclav Hapla args: -dm_plex_hdf5_force_sequential
3805cca0b87SVaclav Hapla test:
38138883f1bSVaclav Hapla suffix: hdf5
382bcde9acaSMatthew G. Knepley args: -orig_dm_plex_interpolate 1
3835cca0b87SVaclav Hapla test:
38438883f1bSVaclav Hapla suffix: hdf5_repart
3855cca0b87SVaclav Hapla requires: parmetis
3865cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type parmetis
387bcde9acaSMatthew G. Knepley args: -orig_dm_plex_interpolate 1
3885cca0b87SVaclav Hapla test:
3895cca0b87SVaclav Hapla TODO: Parallel partitioning of uninterpolated meshes not supported
39038883f1bSVaclav Hapla suffix: hdf5_repart_ppu
3915cca0b87SVaclav Hapla requires: parmetis
3925cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type parmetis
393bcde9acaSMatthew G. Knepley args: -orig_dm_plex_interpolate 0
3945cca0b87SVaclav Hapla
3955cca0b87SVaclav Hapla # reproduce PetscSFView() crash - fixed, left as regression test
3965cca0b87SVaclav Hapla test:
3975cca0b87SVaclav Hapla suffix: new_dm_view
3985cca0b87SVaclav Hapla requires: exodusii
3995cca0b87SVaclav Hapla nsize: 2
400bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/TwoQuads.exo -new_dm_view ascii:ex5_new_dm_view.log:ascii_info_detail
4018886822aSVaclav Hapla args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0 3.0.0}}
402*e0008caeSPierre Jolivet output_file: output/empty.out
40396b7d01aSVaclav Hapla
404bd412c90SVaclav Hapla # test backward compatibility with petsc_hdf5 format version 1.0.0, serial idempotence
40596b7d01aSVaclav Hapla testset:
40696b7d01aSVaclav Hapla suffix: 10-v3.16.0-v1.0.0
40796b7d01aSVaclav Hapla requires: hdf5 !complex datafilespath
408e6ef0f20SVaclav Hapla args: -dm_plex_check_all -compare -compare_labels
4099ed1248dSVaclav Hapla args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0 3.0.0}} -use_low_level_functions {{0 1}}
41096b7d01aSVaclav Hapla test:
41196b7d01aSVaclav Hapla suffix: a
412bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
41396b7d01aSVaclav Hapla test:
41496b7d01aSVaclav Hapla suffix: b
41596b7d01aSVaclav Hapla TODO: broken
416bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
41796b7d01aSVaclav Hapla test:
41896b7d01aSVaclav Hapla suffix: c
419bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
42096b7d01aSVaclav Hapla test:
42196b7d01aSVaclav Hapla suffix: d
422bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5
42396b7d01aSVaclav Hapla test:
42496b7d01aSVaclav Hapla suffix: e
425bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
42696b7d01aSVaclav Hapla test:
42796b7d01aSVaclav Hapla suffix: f
428bcde9acaSMatthew G. Knepley args: -orig_dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
429bcde9acaSMatthew G. Knepley
430bcde9acaSMatthew G. Knepley # test permuted sections with petsc_hdf5 format version 1.0.0
431bcde9acaSMatthew G. Knepley testset:
432bcde9acaSMatthew G. Knepley suffix: 11
433bcde9acaSMatthew G. Knepley requires: hdf5 triangle
434bcde9acaSMatthew G. Knepley args: -field
435bcde9acaSMatthew G. Knepley args: -dm_plex_check_all -compare -compare_labels
436bcde9acaSMatthew G. Knepley args: -orig_dm_plex_box_faces 3,3 -dm_plex_view_hdf5_storage_version 1.0.0
437bcde9acaSMatthew G. Knepley args: -orig_dm_reorder_section -orig_dm_reorder_section_type reverse
438bcde9acaSMatthew G. Knepley
439bcde9acaSMatthew G. Knepley test:
440bcde9acaSMatthew G. Knepley suffix: serial
441bcde9acaSMatthew G. Knepley test:
442bcde9acaSMatthew G. Knepley suffix: serial_no_perm
443bcde9acaSMatthew G. Knepley args: -orig_dm_ignore_perm_output
44496b7d01aSVaclav Hapla
4455cca0b87SVaclav Hapla TEST*/
446