Actual source code: ex5.c
petsc-3.8.4 2018-03-24
1: static char help[] = "Load and save the mesh and fields to HDF5 and ExodusII\n\n";
3: #include <petscdmplex.h>
4: #include <petscviewerhdf5.h>
5: #include <petscsf.h>
7: typedef struct {
8: PetscBool interpolate; /* Generate intermediate mesh elements */
9: char filename[PETSC_MAX_PATH_LEN]; /* Mesh filename */
10: } AppCtx;
12: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
13: {
17: options->interpolate = PETSC_FALSE;
18: options->filename[0] = '\0';
20: PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
21: PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex5.c", options->interpolate, &options->interpolate, NULL);
22: PetscOptionsString("-filename", "The mesh file", "ex5.c", options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL);
23: PetscOptionsEnd();
24: return(0);
25: };
27: int main(int argc, char **argv)
28: {
29: DM dm, dmdist, dmnew;
30: AppCtx user;
31: PetscViewer v;
32: PetscSF pointSF;
33: PetscBool flg;
36: PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
37: ProcessOptions(PETSC_COMM_WORLD, &user);
38: DMPlexCreateFromFile(PETSC_COMM_WORLD, user.filename, user.interpolate, &dm);
39: DMViewFromOptions(dm, NULL, "-orig_dm_view");
40: DMPlexDistribute(dm, 0, &pointSF, &dmdist);
41: if (dmdist) {
42: DMDestroy(&dm);
43: dm = dmdist;
44: }
45: PetscSFDestroy(&pointSF);
46: DMSetFromOptions(dm);
47: DMViewFromOptions(dm, NULL, "-dm_view");
49: PetscViewerHDF5Open(PetscObjectComm((PetscObject) dm), "dmdist.h5", FILE_MODE_WRITE, &v);
50: DMView(dm, v);
51: PetscViewerDestroy(&v);
53: DMCreate(PetscObjectComm((PetscObject) dm), &dmnew);
54: DMSetType(dmnew, DMPLEX);
55: PetscViewerHDF5Open(PETSC_COMM_WORLD, "dmdist.h5", FILE_MODE_READ, &v);
56: DMLoad(dmnew, v);
57: PetscViewerDestroy(&v);
58: DMViewFromOptions(dmnew, NULL, "-new_dm_view");
60: /* The NATIVE format for coordiante viewing is killing parallel output, since we have a local vector. Map it to global, and it will work. */
61: DMPlexEqual(dmnew, dm, &flg);
62: if (flg) {PetscPrintf(PETSC_COMM_WORLD,"DMs equal\n");}
63: else {PetscPrintf(PETSC_COMM_WORLD,"DMs are not equal\n");}
65: DMDestroy(&dm);
66: DMDestroy(&dmnew);
67: PetscFinalize();
68: return ierr;
69: }
71: /*TEST
72: build:
73: requires: exodusii
74: # Idempotence of saving/loading
75: test:
76: suffix: 0
77: requires: exodusii broken
78: args: -filename ${PETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail
79: test:
80: suffix: 1
81: requires: exodusii broken
82: nsize: 2
83: args: -filename ${PETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail
85: TEST*/