xref: /petsc/src/dm/impls/plex/tutorials/ex5.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
17d26aeb3SVaclav Hapla static char help[] = "Demonstrate HDF5/XDMF load-save-reload cycle\n\n";
27d26aeb3SVaclav Hapla 
37d26aeb3SVaclav Hapla #include <petscdmplex.h>
47d26aeb3SVaclav Hapla #include <petscviewerhdf5.h>
57d26aeb3SVaclav Hapla #define EX "ex5.c"
67d26aeb3SVaclav Hapla 
77d26aeb3SVaclav Hapla typedef struct {
87d26aeb3SVaclav Hapla   char              infile[PETSC_MAX_PATH_LEN];  /* Input mesh filename */
97d26aeb3SVaclav Hapla   char              outfile[PETSC_MAX_PATH_LEN]; /* Dump/reload mesh filename */
107d26aeb3SVaclav Hapla   PetscViewerFormat informat;                    /* Input mesh format */
117d26aeb3SVaclav Hapla   PetscViewerFormat outformat;                   /* Dump/reload mesh format */
127d26aeb3SVaclav Hapla   PetscBool         redistribute;                /* Redistribute the mesh */
13efa12513Sksagiyam   PetscBool         heterogeneous;               /* Test save on N / load on M */
147d26aeb3SVaclav Hapla   PetscInt          ntimes;                      /* How many times do the cycle */
157d26aeb3SVaclav Hapla } AppCtx;
167d26aeb3SVaclav Hapla 
177d26aeb3SVaclav Hapla static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
187d26aeb3SVaclav Hapla {
197d26aeb3SVaclav Hapla   PetscBool      flg;
207d26aeb3SVaclav Hapla   PetscErrorCode ierr;
217d26aeb3SVaclav Hapla 
227d26aeb3SVaclav Hapla   PetscFunctionBeginUser;
237d26aeb3SVaclav Hapla   options->infile[0]     = '\0';
247d26aeb3SVaclav Hapla   options->outfile[0]    = '\0';
257d26aeb3SVaclav Hapla   options->informat      = PETSC_VIEWER_HDF5_XDMF;
267d26aeb3SVaclav Hapla   options->outformat     = PETSC_VIEWER_HDF5_XDMF;
277d26aeb3SVaclav Hapla   options->redistribute  = PETSC_TRUE;
28efa12513Sksagiyam   options->heterogeneous = PETSC_FALSE;
297d26aeb3SVaclav Hapla   options->ntimes        = 2;
307d26aeb3SVaclav Hapla   ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsString("-infile", "The input mesh file", EX, options->infile, options->infile, sizeof(options->infile), &flg));
32*28b400f6SJacob Faibussowitsch   PetscCheck(flg,comm, PETSC_ERR_USER_INPUT, "-infile needs to be specified");
335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsString("-outfile", "The output mesh file (by default it's the same as infile)", EX, options->outfile, options->outfile, sizeof(options->outfile), &flg));
34*28b400f6SJacob Faibussowitsch   PetscCheck(flg,comm, PETSC_ERR_USER_INPUT, "-outfile needs to be specified");
355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsEnum("-informat", "Input mesh format", EX, PetscViewerFormats, (PetscEnum)options->informat, (PetscEnum*)&options->informat, NULL));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsEnum("-outformat", "Dump/reload mesh format", EX, PetscViewerFormats, (PetscEnum)options->outformat, (PetscEnum*)&options->outformat, NULL));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-redistribute", "Redistribute the mesh", EX, options->redistribute, &options->redistribute, NULL));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-heterogeneous", "Test save on N / load on M", EX, options->heterogeneous, &options->heterogeneous, NULL));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-ntimes", "How many times do the cycle", EX, options->ntimes, &options->ntimes, NULL));
401e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
417d26aeb3SVaclav Hapla   PetscFunctionReturn(0);
427d26aeb3SVaclav Hapla };
437d26aeb3SVaclav Hapla 
447d26aeb3SVaclav Hapla //TODO test DMLabel I/O (not yet working for PETSC_VIEWER_HDF5_XDMF)
457d26aeb3SVaclav Hapla int main(int argc, char **argv)
467d26aeb3SVaclav Hapla {
477d26aeb3SVaclav Hapla   AppCtx            user;
48efa12513Sksagiyam   MPI_Comm          comm;
49efa12513Sksagiyam   PetscMPIInt       gsize, grank, mycolor;
507d26aeb3SVaclav Hapla   PetscInt          i;
517d26aeb3SVaclav Hapla   PetscBool         flg;
527d26aeb3SVaclav Hapla   PetscErrorCode    ierr;
53c9ad657eSksagiyam   const char        exampleDMPlexName[] = "DMPlex Object";
54efa12513Sksagiyam   const char        *infilename;
55efa12513Sksagiyam   PetscViewerFormat informat;
567d26aeb3SVaclav Hapla 
577d26aeb3SVaclav Hapla   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
585f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user));
595f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&gsize));
605f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&grank));
617d26aeb3SVaclav Hapla 
627d26aeb3SVaclav Hapla   for (i=0; i<user.ntimes; i++) {
63efa12513Sksagiyam     if (i==0) {
64efa12513Sksagiyam       /* Use infile/informat for the initial load */
65efa12513Sksagiyam       infilename = user.infile;
66efa12513Sksagiyam       informat   = user.informat;
67efa12513Sksagiyam     } else {
68efa12513Sksagiyam       /* Use outfile/outformat for all I/O except the very initial load */
69efa12513Sksagiyam       infilename = user.outfile;
70efa12513Sksagiyam       informat   = user.outformat;
71efa12513Sksagiyam     }
72efa12513Sksagiyam 
73efa12513Sksagiyam     if (user.heterogeneous) {
74efa12513Sksagiyam       mycolor = (PetscMPIInt)(grank > user.ntimes-i);
75efa12513Sksagiyam     } else {
76efa12513Sksagiyam       mycolor = (PetscMPIInt)0;
77efa12513Sksagiyam       /* comm = PETSC_COMM_WORLD; */
78efa12513Sksagiyam     }
795f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_split(PETSC_COMM_WORLD,mycolor,grank,&comm));
80efa12513Sksagiyam 
81efa12513Sksagiyam     if (mycolor == 0) {
82efa12513Sksagiyam       /* Load/Save only on processes with mycolor == 0 */
837d26aeb3SVaclav Hapla       DM                dm;
847d26aeb3SVaclav Hapla       PetscPartitioner  part;
857d26aeb3SVaclav Hapla       PetscViewer       v;
867d26aeb3SVaclav Hapla 
875f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm, "Begin cycle %D\n",i));
887d26aeb3SVaclav Hapla 
897d26aeb3SVaclav Hapla       /* Load data from XDMF into dm in parallel */
907d26aeb3SVaclav Hapla       /* We could also use
915f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexCreateFromFile(PETSC_COMM_WORLD, user.filename, "ex5_plex", PETSC_TRUE, &dm));
927d26aeb3SVaclav Hapla         This currently support a few more formats than DMLoad().
937d26aeb3SVaclav Hapla       */
945f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerHDF5Open(comm, infilename, FILE_MODE_READ, &v));
955f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerPushFormat(v, informat));
965f80ce2aSJacob Faibussowitsch       CHKERRQ(DMCreate(comm, &dm));
975f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSetType(dm, DMPLEX));
985f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectSetName((PetscObject) dm, exampleDMPlexName));
995f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSetOptionsPrefix(dm,"loaded_"));
1005f80ce2aSJacob Faibussowitsch       CHKERRQ(DMLoad(dm, v));
1015f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexDistributeSetDefault(dm, PETSC_FALSE));
1025f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSetFromOptions(dm));
1035f80ce2aSJacob Faibussowitsch       CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view"));
1045f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerPopFormat(v));
1055f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerDestroy(&v));
1067d26aeb3SVaclav Hapla 
1077d26aeb3SVaclav Hapla       /* We just test/demonstrate DM is indeed distributed - unneeded in the application code */
1085f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexIsDistributed(dm, &flg));
1095f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm, "Loaded mesh distributed? %s\n", PetscBools[flg]));
1107d26aeb3SVaclav Hapla 
1117d26aeb3SVaclav Hapla       /* Interpolate */
1127d26aeb3SVaclav Hapla       //TODO we want to be able to do this from options in DMSetFromOptions() probably
1137d26aeb3SVaclav Hapla       //TODO we want to be able to do this in-place
1147d26aeb3SVaclav Hapla       {
1157d26aeb3SVaclav Hapla         DM idm;
1167d26aeb3SVaclav Hapla 
1175f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexInterpolate(dm, &idm));
1185f80ce2aSJacob Faibussowitsch         CHKERRQ(DMDestroy(&dm));
1197d26aeb3SVaclav Hapla         dm   = idm;
1205f80ce2aSJacob Faibussowitsch           CHKERRQ(DMSetOptionsPrefix(dm,"interpolated_"));
1215f80ce2aSJacob Faibussowitsch           CHKERRQ(DMSetFromOptions(dm));
1225f80ce2aSJacob Faibussowitsch           CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view"));
1237d26aeb3SVaclav Hapla       }
1247d26aeb3SVaclav Hapla 
1257d26aeb3SVaclav Hapla       /* Redistribute */
1267d26aeb3SVaclav Hapla       //TODO we want to be able to do this from options in DMSetFromOptions() probably
1277d26aeb3SVaclav Hapla       if (user.redistribute) {
1287d26aeb3SVaclav Hapla         DM dmdist;
1297d26aeb3SVaclav Hapla 
1305f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexGetPartitioner(dm, &part));
1315f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPartitionerSetFromOptions(part));
1325f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexDistribute(dm, 0, NULL, &dmdist));
1337d26aeb3SVaclav Hapla         //TODO we want to be able to do this in-place
1347d26aeb3SVaclav Hapla         if (dmdist) {
1355f80ce2aSJacob Faibussowitsch           CHKERRQ(DMDestroy(&dm));
1367d26aeb3SVaclav Hapla           dm   = dmdist;
1375f80ce2aSJacob Faibussowitsch           CHKERRQ(DMSetOptionsPrefix(dm,"redistributed_"));
1385f80ce2aSJacob Faibussowitsch           CHKERRQ(DMSetFromOptions(dm));
1395f80ce2aSJacob Faibussowitsch           CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view"));
1407d26aeb3SVaclav Hapla         }
1417d26aeb3SVaclav Hapla       }
1427d26aeb3SVaclav Hapla 
1437d26aeb3SVaclav Hapla       /* Save redistributed dm to XDMF in parallel and destroy it */
1445f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerHDF5Open(comm, user.outfile, FILE_MODE_WRITE, &v));
1455f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerPushFormat(v, user.outformat));
1465f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectSetName((PetscObject) dm, exampleDMPlexName));
1475f80ce2aSJacob Faibussowitsch       CHKERRQ(DMView(dm, v));
1485f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerPopFormat(v));
1495f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerDestroy(&v));
1505f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDestroy(&dm));
1517d26aeb3SVaclav Hapla 
1525f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm, "End   cycle %D\n--------\n",i));
1537d26aeb3SVaclav Hapla     }
1545f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_free(&comm));
1555f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD));
156efa12513Sksagiyam   }
1577d26aeb3SVaclav Hapla 
1587d26aeb3SVaclav Hapla   /* Final clean-up */
1597d26aeb3SVaclav Hapla   ierr = PetscFinalize();
1607d26aeb3SVaclav Hapla   return ierr;
1617d26aeb3SVaclav Hapla }
1627d26aeb3SVaclav Hapla 
1637d26aeb3SVaclav Hapla /*TEST
1647d26aeb3SVaclav Hapla   build:
1657d26aeb3SVaclav Hapla     requires: hdf5
1667d26aeb3SVaclav Hapla   testset:
1677d26aeb3SVaclav Hapla     suffix: 0
1687d26aeb3SVaclav Hapla     requires: !complex
169efa12513Sksagiyam     nsize: 4
1707d26aeb3SVaclav Hapla     args: -infile ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -informat hdf5_xdmf
1717d26aeb3SVaclav Hapla     args: -outfile ex5_dump.h5 -outformat {{hdf5_xdmf hdf5_petsc}separate output}
1727d26aeb3SVaclav Hapla     args: -ntimes 3
173efa12513Sksagiyam     args: -loaded_dm_view -interpolated_dm_view -redistributed_dm_view
1747d26aeb3SVaclav Hapla     test:
1757d26aeb3SVaclav Hapla       # this partitioner should not shuffle anything, it should yield the same partititioning as the XDMF reader - added just for testing
1767d26aeb3SVaclav Hapla       suffix: simple
1777d26aeb3SVaclav Hapla       args: -petscpartitioner_type simple
1787d26aeb3SVaclav Hapla     test:
1797d26aeb3SVaclav Hapla       suffix: parmetis
1807d26aeb3SVaclav Hapla       requires: parmetis
1817d26aeb3SVaclav Hapla       args: -petscpartitioner_type parmetis
1827d26aeb3SVaclav Hapla     test:
1837d26aeb3SVaclav Hapla       suffix: ptscotch
1847d26aeb3SVaclav Hapla       requires: ptscotch
1857d26aeb3SVaclav Hapla       args: -petscpartitioner_type ptscotch
186efa12513Sksagiyam 
187efa12513Sksagiyam   testset:
188efa12513Sksagiyam     suffix: 1
189efa12513Sksagiyam     requires: !complex
190efa12513Sksagiyam     nsize: 4
191efa12513Sksagiyam     args: -infile ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -informat hdf5_xdmf
192efa12513Sksagiyam     args: -outfile ex5_dump.h5 -outformat {{hdf5_xdmf hdf5_petsc}separate output}
193efa12513Sksagiyam     args: -ntimes 3
194efa12513Sksagiyam     args: -heterogeneous True
195efa12513Sksagiyam     args: -loaded_dm_view -interpolated_dm_view -redistributed_dm_view
196efa12513Sksagiyam     test:
197efa12513Sksagiyam       suffix: simple
198efa12513Sksagiyam       args: -petscpartitioner_type simple
199efa12513Sksagiyam     test:
200efa12513Sksagiyam       suffix: parmetis
201efa12513Sksagiyam       requires: parmetis
202efa12513Sksagiyam       args: -petscpartitioner_type parmetis
203efa12513Sksagiyam     test:
204efa12513Sksagiyam       suffix: ptscotch
205efa12513Sksagiyam       requires: ptscotch
206efa12513Sksagiyam       args: -petscpartitioner_type ptscotch
207efa12513Sksagiyam 
2087d26aeb3SVaclav Hapla TEST*/
209