xref: /petsc/src/dm/impls/plex/tests/ex55.c (revision bd412c90fba8895e9763ccee76c7dd2e09a982d7)
15cca0b87SVaclav Hapla static char help[] = "Load and save the mesh and fields to HDF5 and ExodusII\n\n";
25cca0b87SVaclav Hapla 
3*bd412c90SVaclav 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 */
115cca0b87SVaclav Hapla   PetscBool         interpolate;                  /* Generate intermediate mesh elements */
12*bd412c90SVaclav Hapla   char              fname[PETSC_MAX_PATH_LEN];    /* Mesh filename */
13*bd412c90SVaclav Hapla   char              ofname[PETSC_MAX_PATH_LEN];   /* Output mesh filename */
14ed5e4e85SVaclav Hapla   char              meshname[PETSC_MAX_PATH_LEN]; /* Mesh name */
155cca0b87SVaclav Hapla   PetscViewerFormat format;                       /* Format to write and read */
165cca0b87SVaclav Hapla   PetscBool         second_write_read;            /* Write and read for the 2nd time */
17806c51deSksagiyam   PetscBool         use_low_level_functions;      /* Use low level functions for viewing and loading */
185cca0b87SVaclav Hapla } AppCtx;
195cca0b87SVaclav Hapla 
20d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
21d71ae5a4SJacob Faibussowitsch {
225cca0b87SVaclav Hapla   PetscFunctionBeginUser;
235cca0b87SVaclav Hapla   options->compare                 = PETSC_FALSE;
2438883f1bSVaclav Hapla   options->compare_labels          = PETSC_FALSE;
255cca0b87SVaclav Hapla   options->distribute              = PETSC_TRUE;
265cca0b87SVaclav Hapla   options->interpolate             = PETSC_FALSE;
27*bd412c90SVaclav Hapla   options->fname[0]                = '\0';
28ed5e4e85SVaclav Hapla   options->meshname[0]             = '\0';
295cca0b87SVaclav Hapla   options->format                  = PETSC_VIEWER_DEFAULT;
305cca0b87SVaclav Hapla   options->second_write_read       = PETSC_FALSE;
31806c51deSksagiyam   options->use_low_level_functions = PETSC_FALSE;
32*bd412c90SVaclav Hapla   PetscCall(PetscStrcpy(options->ofname, "ex55.h5"));
335cca0b87SVaclav Hapla 
34d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-compare", "Compare the meshes using DMPlexEqual()", "ex55.c", options->compare, &options->compare, NULL));
369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-compare_labels", "Compare labels in the meshes using DMCompareLabels()", "ex55.c", options->compare_labels, &options->compare_labels, NULL));
379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-distribute", "Distribute the mesh", "ex55.c", options->distribute, &options->distribute, NULL));
389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex55.c", options->interpolate, &options->interpolate, NULL));
39*bd412c90SVaclav Hapla   PetscCall(PetscOptionsString("-filename", "The mesh file", "ex55.c", options->fname, options->fname, sizeof(options->fname), NULL));
40*bd412c90SVaclav Hapla   PetscCall(PetscOptionsString("-ofilename", "The output mesh file", "ex55.c", options->ofname, options->ofname, sizeof(options->ofname), NULL));
419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-meshname", "The mesh file", "ex55.c", options->meshname, options->meshname, sizeof(options->meshname), NULL));
429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-format", "Format to write and read", "ex55.c", PetscViewerFormats, (PetscEnum)options->format, (PetscEnum *)&options->format, NULL));
439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-second_write_read", "Write and read for the 2nd time", "ex55.c", options->second_write_read, &options->second_write_read, NULL));
449566063dSJacob 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));
45d0609cedSBarry Smith   PetscOptionsEnd();
465cca0b87SVaclav Hapla   PetscFunctionReturn(0);
47*bd412c90SVaclav Hapla }
485cca0b87SVaclav Hapla 
49*bd412c90SVaclav Hapla static PetscErrorCode CheckDistributed(DM dm, PetscBool expectedDistributed)
50*bd412c90SVaclav Hapla {
51*bd412c90SVaclav Hapla   PetscMPIInt size;
52*bd412c90SVaclav Hapla   PetscBool   distributed;
53*bd412c90SVaclav Hapla   const char  YES[] = "DISTRIBUTED";
54*bd412c90SVaclav Hapla   const char  NO[]  = "NOT DISTRIBUTED";
55*bd412c90SVaclav Hapla 
56*bd412c90SVaclav Hapla   PetscFunctionBeginUser;
57*bd412c90SVaclav Hapla   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
58*bd412c90SVaclav Hapla   if (size < 2) PetscFunctionReturn(0);
59*bd412c90SVaclav Hapla   PetscCall(DMPlexIsDistributed(dm, &distributed));
60*bd412c90SVaclav Hapla   PetscCheck(distributed == expectedDistributed, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Expected DM being %s but actually is %s", expectedDistributed ? YES : NO, distributed ? YES : NO);
61*bd412c90SVaclav Hapla   PetscFunctionReturn(0);
62*bd412c90SVaclav Hapla }
63*bd412c90SVaclav Hapla 
64*bd412c90SVaclav Hapla static PetscErrorCode CheckInterpolated(DM dm, PetscBool expectedInterpolated)
65*bd412c90SVaclav Hapla {
66*bd412c90SVaclav Hapla   DMPlexInterpolatedFlag iflg;
67*bd412c90SVaclav Hapla   PetscBool              interpolated;
68*bd412c90SVaclav Hapla   const char             YES[] = "INTERPOLATED";
69*bd412c90SVaclav Hapla   const char             NO[]  = "NOT INTERPOLATED";
70*bd412c90SVaclav Hapla 
71*bd412c90SVaclav Hapla   PetscFunctionBeginUser;
72*bd412c90SVaclav Hapla   PetscCall(DMPlexIsInterpolatedCollective(dm, &iflg));
73*bd412c90SVaclav Hapla   interpolated = (PetscBool)(iflg == DMPLEX_INTERPOLATED_FULL);
74*bd412c90SVaclav Hapla   PetscCheck(interpolated == expectedInterpolated, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Expected DM being %s but actually is %s", expectedInterpolated ? YES : NO, interpolated ? YES : NO);
75*bd412c90SVaclav Hapla   PetscFunctionReturn(0);
76*bd412c90SVaclav Hapla }
77*bd412c90SVaclav Hapla 
78*bd412c90SVaclav Hapla static PetscErrorCode CheckDistributedInterpolated(DM dm, PetscViewer v, AppCtx *user)
79*bd412c90SVaclav Hapla {
80*bd412c90SVaclav Hapla   PetscViewerFormat format;
81*bd412c90SVaclav Hapla   PetscBool         distributed, interpolated;
82*bd412c90SVaclav Hapla 
83*bd412c90SVaclav Hapla   PetscFunctionBeginUser;
84*bd412c90SVaclav Hapla   PetscCall(PetscViewerGetFormat(v, &format));
85*bd412c90SVaclav Hapla   switch (format) {
86*bd412c90SVaclav Hapla   case PETSC_VIEWER_HDF5_XDMF:
87*bd412c90SVaclav Hapla   case PETSC_VIEWER_HDF5_VIZ: {
88*bd412c90SVaclav Hapla     distributed  = PETSC_TRUE;
89*bd412c90SVaclav Hapla     interpolated = PETSC_FALSE;
90*bd412c90SVaclav Hapla   }; break;
91*bd412c90SVaclav Hapla   case PETSC_VIEWER_HDF5_PETSC:
92*bd412c90SVaclav Hapla   case PETSC_VIEWER_DEFAULT:
93*bd412c90SVaclav Hapla   case PETSC_VIEWER_NATIVE: {
94*bd412c90SVaclav Hapla     DMPlexStorageVersion version;
95*bd412c90SVaclav Hapla 
96*bd412c90SVaclav Hapla     PetscCall(PetscViewerHDF5GetDMPlexStorageVersionReading(v, &version));
97*bd412c90SVaclav Hapla     distributed  = (PetscBool)(version->major >= 3);
98*bd412c90SVaclav Hapla     interpolated = user->interpolate;
99*bd412c90SVaclav Hapla   }; break;
100*bd412c90SVaclav Hapla   default: {
101*bd412c90SVaclav Hapla     distributed  = PETSC_FALSE;
102*bd412c90SVaclav Hapla     interpolated = user->interpolate;
103*bd412c90SVaclav Hapla   }
104*bd412c90SVaclav Hapla   }
105*bd412c90SVaclav Hapla   PetscCall(CheckDistributed(dm, distributed));
106*bd412c90SVaclav Hapla   PetscCall(CheckInterpolated(dm, interpolated));
107*bd412c90SVaclav Hapla   PetscFunctionReturn(0);
108*bd412c90SVaclav Hapla }
109*bd412c90SVaclav Hapla 
110*bd412c90SVaclav Hapla static PetscErrorCode DMPlexWriteAndReadHDF5(DM dm, const char filename[], const char prefix[], AppCtx *user, DM *dm_new)
111d71ae5a4SJacob Faibussowitsch {
1125cca0b87SVaclav Hapla   DM          dmnew;
113ed5e4e85SVaclav Hapla   const char  savedName[]  = "Mesh";
114ed5e4e85SVaclav Hapla   const char  loadedName[] = "Mesh_new";
1155cca0b87SVaclav Hapla   PetscViewer v;
1165cca0b87SVaclav Hapla 
1175cca0b87SVaclav Hapla   PetscFunctionBeginUser;
1189566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Open(PetscObjectComm((PetscObject)dm), filename, FILE_MODE_WRITE, &v));
119*bd412c90SVaclav Hapla   PetscCall(PetscViewerPushFormat(v, user->format));
1209566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)dm, savedName));
121*bd412c90SVaclav Hapla   if (user->use_low_level_functions) {
1229566063dSJacob Faibussowitsch     PetscCall(DMPlexTopologyView(dm, v));
1239566063dSJacob Faibussowitsch     PetscCall(DMPlexCoordinatesView(dm, v));
1249566063dSJacob Faibussowitsch     PetscCall(DMPlexLabelsView(dm, v));
125806c51deSksagiyam   } else {
1269566063dSJacob Faibussowitsch     PetscCall(DMView(dm, v));
127806c51deSksagiyam   }
1289566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(v, FILE_MODE_READ));
1299566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dmnew));
1309566063dSJacob Faibussowitsch   PetscCall(DMSetType(dmnew, DMPLEX));
131*bd412c90SVaclav Hapla   PetscCall(DMPlexDistributeSetDefault(dmnew, PETSC_FALSE));
1329566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)dmnew, savedName));
1339566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(dmnew, prefix));
134*bd412c90SVaclav Hapla   if (user->use_low_level_functions) {
135c9ad657eSksagiyam     PetscSF sfXC;
136c9ad657eSksagiyam 
1379566063dSJacob Faibussowitsch     PetscCall(DMPlexTopologyLoad(dmnew, v, &sfXC));
1389566063dSJacob Faibussowitsch     PetscCall(DMPlexCoordinatesLoad(dmnew, v, sfXC));
1399566063dSJacob Faibussowitsch     PetscCall(DMPlexLabelsLoad(dmnew, v, sfXC));
1409566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sfXC));
141806c51deSksagiyam   } else {
1429566063dSJacob Faibussowitsch     PetscCall(DMLoad(dmnew, v));
143806c51deSksagiyam   }
144*bd412c90SVaclav Hapla   PetscCall(CheckDistributedInterpolated(dmnew, v, user));
1459566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)dmnew, loadedName));
1469566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(v));
1479566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&v));
1485cca0b87SVaclav Hapla   *dm_new = dmnew;
1495cca0b87SVaclav Hapla   PetscFunctionReturn(0);
1505cca0b87SVaclav Hapla }
1515cca0b87SVaclav Hapla 
152d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
153d71ae5a4SJacob Faibussowitsch {
1545cca0b87SVaclav Hapla   DM               dm, dmnew;
1555cca0b87SVaclav Hapla   PetscPartitioner part;
1565cca0b87SVaclav Hapla   AppCtx           user;
1575cca0b87SVaclav Hapla   PetscBool        flg;
1585cca0b87SVaclav Hapla 
159327415f7SBarry Smith   PetscFunctionBeginUser;
1609566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1619566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
162*bd412c90SVaclav Hapla   PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, user.fname, user.meshname, user.interpolate, &dm));
1639566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(dm, "orig_"));
1649566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
165*bd412c90SVaclav Hapla   PetscCall(CheckInterpolated(dm, user.interpolate));
1665cca0b87SVaclav Hapla 
1675cca0b87SVaclav Hapla   if (user.distribute) {
1685cca0b87SVaclav Hapla     DM dmdist;
1695cca0b87SVaclav Hapla 
1709566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(dm, &part));
171*bd412c90SVaclav Hapla     PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSIMPLE));
1729566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetFromOptions(part));
1739566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(dm, 0, NULL, &dmdist));
1745cca0b87SVaclav Hapla     if (dmdist) {
1759566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dm));
1765cca0b87SVaclav Hapla       dm = dmdist;
177*bd412c90SVaclav Hapla       PetscCall(CheckDistributed(dm, PETSC_TRUE));
178*bd412c90SVaclav Hapla       PetscCall(CheckInterpolated(dm, user.interpolate));
1795cca0b87SVaclav Hapla     }
1805cca0b87SVaclav Hapla   }
1815cca0b87SVaclav Hapla 
1829566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(dm, NULL));
1839566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE));
1849566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
1859566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
1865cca0b87SVaclav Hapla 
187*bd412c90SVaclav Hapla   PetscCall(DMPlexWriteAndReadHDF5(dm, user.ofname, "new_", &user, &dmnew));
1885cca0b87SVaclav Hapla 
1895cca0b87SVaclav Hapla   if (user.second_write_read) {
1909566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dm));
1915cca0b87SVaclav Hapla     dm = dmnew;
192*bd412c90SVaclav Hapla     PetscCall(DMPlexWriteAndReadHDF5(dm, user.ofname, "new_", &user, &dmnew));
1935cca0b87SVaclav Hapla   }
1945cca0b87SVaclav Hapla 
1959566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dmnew, NULL, "-dm_view"));
1965cca0b87SVaclav Hapla 
1975cca0b87SVaclav Hapla   /* This currently makes sense only for sequential meshes. */
1985cca0b87SVaclav Hapla   if (user.compare) {
1999566063dSJacob Faibussowitsch     PetscCall(DMPlexEqual(dmnew, dm, &flg));
20028b400f6SJacob Faibussowitsch     PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "DMs are not equal");
2019566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "DMs equal\n"));
20238883f1bSVaclav Hapla   }
20338883f1bSVaclav Hapla   if (user.compare_labels) {
2049566063dSJacob Faibussowitsch     PetscCall(DMCompareLabels(dmnew, dm, NULL, NULL));
2059566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "DMLabels equal\n"));
2065cca0b87SVaclav Hapla   }
2075cca0b87SVaclav Hapla 
2089566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
2099566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmnew));
2109566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
211b122ec5aSJacob Faibussowitsch   return 0;
2125cca0b87SVaclav Hapla }
2135cca0b87SVaclav Hapla 
2145cca0b87SVaclav Hapla /*TEST
2155cca0b87SVaclav Hapla   build:
2165cca0b87SVaclav Hapla     requires: hdf5
2175cca0b87SVaclav Hapla   # Idempotence of saving/loading
2185cca0b87SVaclav Hapla   #   Have to replace Exodus file, which is creating uninterpolated edges
2195cca0b87SVaclav Hapla   test:
2205cca0b87SVaclav Hapla     suffix: 0
22138883f1bSVaclav Hapla     TODO: broken
22238883f1bSVaclav Hapla     requires: exodusii
2235cca0b87SVaclav Hapla     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail
2245cca0b87SVaclav Hapla     args: -format hdf5_petsc -compare
2255cca0b87SVaclav Hapla   test:
2265cca0b87SVaclav Hapla     suffix: 1
22738883f1bSVaclav Hapla     TODO: broken
22838883f1bSVaclav Hapla     requires: exodusii parmetis !defined(PETSC_USE_64BIT_INDICES)
2295cca0b87SVaclav Hapla     nsize: 2
2305cca0b87SVaclav Hapla     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail
2315cca0b87SVaclav Hapla     args: -petscpartitioner_type parmetis
2325cca0b87SVaclav Hapla     args: -format hdf5_petsc -new_dm_view ascii::ascii_info_detail
233*bd412c90SVaclav Hapla 
2345cca0b87SVaclav Hapla   testset:
2355cca0b87SVaclav Hapla     requires: exodusii
2365cca0b87SVaclav Hapla     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo
2375cca0b87SVaclav Hapla     test:
2385cca0b87SVaclav Hapla       suffix: 2
2395cca0b87SVaclav Hapla       nsize: {{1 2 4 8}separate output}
2405cca0b87SVaclav Hapla       args: -format {{default hdf5_petsc}separate output}
241*bd412c90SVaclav Hapla       args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0}}
2425cca0b87SVaclav Hapla       args: -interpolate {{0 1}separate output}
2435cca0b87SVaclav Hapla     test:
2445cca0b87SVaclav Hapla       suffix: 2a
2455cca0b87SVaclav Hapla       nsize: {{1 2 4 8}separate output}
2465cca0b87SVaclav Hapla       args: -format {{hdf5_xdmf hdf5_viz}separate output}
247*bd412c90SVaclav Hapla 
2485cca0b87SVaclav Hapla   test:
2495cca0b87SVaclav Hapla     suffix: 3
2505cca0b87SVaclav Hapla     requires: exodusii
25138883f1bSVaclav Hapla     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -compare -compare_labels
252*bd412c90SVaclav Hapla     args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0}}
2535cca0b87SVaclav Hapla 
2545cca0b87SVaclav Hapla   # Load HDF5 file in XDMF format in parallel, write, read dm1, write, read dm2, and compare dm1 and dm2
25538883f1bSVaclav Hapla   testset:
2565cca0b87SVaclav Hapla     suffix: 4
2575cca0b87SVaclav Hapla     requires: !complex
25838883f1bSVaclav Hapla     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf
25938883f1bSVaclav Hapla     args: -distribute 0 -second_write_read -compare
26038883f1bSVaclav Hapla     test:
26138883f1bSVaclav Hapla       suffix: hdf5_petsc
26238883f1bSVaclav Hapla       nsize: {{1 2}}
26338883f1bSVaclav Hapla       args: -format hdf5_petsc -compare_labels
264*bd412c90SVaclav Hapla       args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0}}
26538883f1bSVaclav Hapla     test:
26638883f1bSVaclav Hapla       suffix: hdf5_xdmf
26738883f1bSVaclav Hapla       nsize: {{1 3 8}}
26838883f1bSVaclav Hapla       args: -format hdf5_xdmf
2695cca0b87SVaclav Hapla 
270806c51deSksagiyam   # Use low level functions, DMPlexTopologyView()/Load(), DMPlexCoordinatesView()/Load(), and DMPlexLabelsView()/Load()
271*bd412c90SVaclav Hapla   # TODO: The output is very long so keeping just 1.0.0 version. This test should be redesigned or removed.
272806c51deSksagiyam   test:
273806c51deSksagiyam     suffix: 5
274806c51deSksagiyam     requires: exodusii
275806c51deSksagiyam     nsize: 2
276806c51deSksagiyam     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo
277806c51deSksagiyam     args: -dm_view ascii::ascii_info_detail
278806c51deSksagiyam     args: -new_dm_view ascii::ascii_info_detail
279*bd412c90SVaclav Hapla     args: -format hdf5_petsc -use_low_level_functions {{0 1}}
280*bd412c90SVaclav Hapla     args: -dm_plex_view_hdf5_storage_version 1.0.0
281806c51deSksagiyam 
2825cca0b87SVaclav Hapla   testset:
28338883f1bSVaclav Hapla     suffix: 6
28438883f1bSVaclav Hapla     requires: hdf5 !complex datafilespath
28538883f1bSVaclav Hapla     nsize: {{1 3}}
28638883f1bSVaclav Hapla     args: -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_geometry
28738883f1bSVaclav Hapla     args: -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
28838883f1bSVaclav Hapla     args: -format hdf5_petsc -second_write_read -compare -compare_labels
289*bd412c90SVaclav Hapla     args: -interpolate {{0 1}} -distribute {{0 1}}
290*bd412c90SVaclav Hapla     args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0}}
29138883f1bSVaclav Hapla 
29238883f1bSVaclav Hapla   testset:
2935cca0b87SVaclav Hapla     # the same data and settings as dm_impls_plex_tests-ex18_9%
29438883f1bSVaclav Hapla     suffix: 9
2955cca0b87SVaclav Hapla     requires: hdf5 !complex datafilespath
2964d056bb6SVaclav Hapla     nsize: {{1 2 4}}
2975cca0b87SVaclav Hapla     args: -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_geometry
2985cca0b87SVaclav Hapla     args: -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
299*bd412c90SVaclav Hapla     args: -format {{hdf5_petsc hdf5_xdmf}} -second_write_read -compare
3005cca0b87SVaclav Hapla     test:
30138883f1bSVaclav Hapla       suffix: hdf5_seqload
302*bd412c90SVaclav Hapla       args: -distribute
3035cca0b87SVaclav Hapla       args: -interpolate {{0 1}}
3045cca0b87SVaclav Hapla       args: -dm_plex_hdf5_force_sequential
3055cca0b87SVaclav Hapla     test:
30638883f1bSVaclav Hapla       suffix: hdf5_seqload_metis
3075cca0b87SVaclav Hapla       requires: parmetis
3085cca0b87SVaclav Hapla       args: -distribute -petscpartitioner_type parmetis
3095cca0b87SVaclav Hapla       args: -interpolate 1
3105cca0b87SVaclav Hapla       args: -dm_plex_hdf5_force_sequential
3115cca0b87SVaclav Hapla     test:
31238883f1bSVaclav Hapla       suffix: hdf5
313*bd412c90SVaclav Hapla       args: -interpolate 1
3145cca0b87SVaclav Hapla     test:
31538883f1bSVaclav Hapla       suffix: hdf5_repart
3165cca0b87SVaclav Hapla       requires: parmetis
3175cca0b87SVaclav Hapla       args: -distribute -petscpartitioner_type parmetis
3185cca0b87SVaclav Hapla       args: -interpolate 1
3195cca0b87SVaclav Hapla     test:
3205cca0b87SVaclav Hapla       TODO: Parallel partitioning of uninterpolated meshes not supported
32138883f1bSVaclav Hapla       suffix: hdf5_repart_ppu
3225cca0b87SVaclav Hapla       requires: parmetis
3235cca0b87SVaclav Hapla       args: -distribute -petscpartitioner_type parmetis
3245cca0b87SVaclav Hapla       args: -interpolate 0
3255cca0b87SVaclav Hapla 
3265cca0b87SVaclav Hapla   # reproduce PetscSFView() crash - fixed, left as regression test
3275cca0b87SVaclav Hapla   test:
3285cca0b87SVaclav Hapla     suffix: new_dm_view
3295cca0b87SVaclav Hapla     requires: exodusii
3305cca0b87SVaclav Hapla     nsize: 2
3315cca0b87SVaclav Hapla     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/TwoQuads.exo -new_dm_view ascii:ex5_new_dm_view.log:ascii_info_detail
332*bd412c90SVaclav Hapla     args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0}}
33396b7d01aSVaclav Hapla 
334*bd412c90SVaclav Hapla   # test backward compatibility with petsc_hdf5 format version 1.0.0, serial idempotence
33596b7d01aSVaclav Hapla   testset:
33696b7d01aSVaclav Hapla     suffix: 10-v3.16.0-v1.0.0
33796b7d01aSVaclav Hapla     requires: hdf5 !complex datafilespath
338e6ef0f20SVaclav Hapla     args: -dm_plex_check_all -compare -compare_labels
33996b7d01aSVaclav Hapla     args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0}} -use_low_level_functions {{0 1}}
34096b7d01aSVaclav Hapla     test:
34196b7d01aSVaclav Hapla       suffix: a
34296b7d01aSVaclav Hapla       args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
34396b7d01aSVaclav Hapla     test:
34496b7d01aSVaclav Hapla       suffix: b
34596b7d01aSVaclav Hapla       TODO: broken
34696b7d01aSVaclav Hapla       args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
34796b7d01aSVaclav Hapla     test:
34896b7d01aSVaclav Hapla       suffix: c
34996b7d01aSVaclav Hapla       args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
35096b7d01aSVaclav Hapla     test:
35196b7d01aSVaclav Hapla       suffix: d
35296b7d01aSVaclav Hapla       args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5
35396b7d01aSVaclav Hapla     test:
35496b7d01aSVaclav Hapla       suffix: e
35596b7d01aSVaclav Hapla       args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
35696b7d01aSVaclav Hapla     test:
35796b7d01aSVaclav Hapla       suffix: f
35896b7d01aSVaclav Hapla       args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
35996b7d01aSVaclav Hapla 
3605cca0b87SVaclav Hapla TEST*/
361