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)); 3228b400f6SJacob 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)); 3428b400f6SJacob 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; 52c9ad657eSksagiyam const char exampleDMPlexName[] = "DMPlex Object"; 53efa12513Sksagiyam const char *infilename; 54efa12513Sksagiyam PetscViewerFormat informat; 557d26aeb3SVaclav Hapla 56*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL,help)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user)); 585f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&gsize)); 595f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&grank)); 607d26aeb3SVaclav Hapla 617d26aeb3SVaclav Hapla for (i=0; i<user.ntimes; i++) { 62efa12513Sksagiyam if (i==0) { 63efa12513Sksagiyam /* Use infile/informat for the initial load */ 64efa12513Sksagiyam infilename = user.infile; 65efa12513Sksagiyam informat = user.informat; 66efa12513Sksagiyam } else { 67efa12513Sksagiyam /* Use outfile/outformat for all I/O except the very initial load */ 68efa12513Sksagiyam infilename = user.outfile; 69efa12513Sksagiyam informat = user.outformat; 70efa12513Sksagiyam } 71efa12513Sksagiyam 72efa12513Sksagiyam if (user.heterogeneous) { 73efa12513Sksagiyam mycolor = (PetscMPIInt)(grank > user.ntimes-i); 74efa12513Sksagiyam } else { 75efa12513Sksagiyam mycolor = (PetscMPIInt)0; 76efa12513Sksagiyam /* comm = PETSC_COMM_WORLD; */ 77efa12513Sksagiyam } 785f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_split(PETSC_COMM_WORLD,mycolor,grank,&comm)); 79efa12513Sksagiyam 80efa12513Sksagiyam if (mycolor == 0) { 81efa12513Sksagiyam /* Load/Save only on processes with mycolor == 0 */ 827d26aeb3SVaclav Hapla DM dm; 837d26aeb3SVaclav Hapla PetscPartitioner part; 847d26aeb3SVaclav Hapla PetscViewer v; 857d26aeb3SVaclav Hapla 865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(comm, "Begin cycle %D\n",i)); 877d26aeb3SVaclav Hapla 887d26aeb3SVaclav Hapla /* Load data from XDMF into dm in parallel */ 897d26aeb3SVaclav Hapla /* We could also use 905f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromFile(PETSC_COMM_WORLD, user.filename, "ex5_plex", PETSC_TRUE, &dm)); 917d26aeb3SVaclav Hapla This currently support a few more formats than DMLoad(). 927d26aeb3SVaclav Hapla */ 935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerHDF5Open(comm, infilename, FILE_MODE_READ, &v)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(v, informat)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, &dm)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(dm, DMPLEX)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) dm, exampleDMPlexName)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOptionsPrefix(dm,"loaded_")); 995f80ce2aSJacob Faibussowitsch CHKERRQ(DMLoad(dm, v)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeSetDefault(dm, PETSC_FALSE)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(dm)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view")); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(v)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&v)); 1057d26aeb3SVaclav Hapla 1067d26aeb3SVaclav Hapla /* We just test/demonstrate DM is indeed distributed - unneeded in the application code */ 1075f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsDistributed(dm, &flg)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(comm, "Loaded mesh distributed? %s\n", PetscBools[flg])); 1097d26aeb3SVaclav Hapla 1107d26aeb3SVaclav Hapla /* Interpolate */ 1117d26aeb3SVaclav Hapla //TODO we want to be able to do this from options in DMSetFromOptions() probably 1127d26aeb3SVaclav Hapla //TODO we want to be able to do this in-place 1137d26aeb3SVaclav Hapla { 1147d26aeb3SVaclav Hapla DM idm; 1157d26aeb3SVaclav Hapla 1165f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(dm, &idm)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 1187d26aeb3SVaclav Hapla dm = idm; 1195f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOptionsPrefix(dm,"interpolated_")); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(dm)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view")); 1227d26aeb3SVaclav Hapla } 1237d26aeb3SVaclav Hapla 1247d26aeb3SVaclav Hapla /* Redistribute */ 1257d26aeb3SVaclav Hapla //TODO we want to be able to do this from options in DMSetFromOptions() probably 1267d26aeb3SVaclav Hapla if (user.redistribute) { 1277d26aeb3SVaclav Hapla DM dmdist; 1287d26aeb3SVaclav Hapla 1295f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetPartitioner(dm, &part)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPartitionerSetFromOptions(part)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistribute(dm, 0, NULL, &dmdist)); 1327d26aeb3SVaclav Hapla //TODO we want to be able to do this in-place 1337d26aeb3SVaclav Hapla if (dmdist) { 1345f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 1357d26aeb3SVaclav Hapla dm = dmdist; 1365f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOptionsPrefix(dm,"redistributed_")); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(dm)); 1385f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view")); 1397d26aeb3SVaclav Hapla } 1407d26aeb3SVaclav Hapla } 1417d26aeb3SVaclav Hapla 1427d26aeb3SVaclav Hapla /* Save redistributed dm to XDMF in parallel and destroy it */ 1435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerHDF5Open(comm, user.outfile, FILE_MODE_WRITE, &v)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(v, user.outformat)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) dm, exampleDMPlexName)); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(DMView(dm, v)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(v)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&v)); 1495f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 1507d26aeb3SVaclav Hapla 1515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(comm, "End cycle %D\n--------\n",i)); 1527d26aeb3SVaclav Hapla } 1535f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_free(&comm)); 1545f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD)); 155efa12513Sksagiyam } 1567d26aeb3SVaclav Hapla 1577d26aeb3SVaclav Hapla /* Final clean-up */ 158*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 159*b122ec5aSJacob Faibussowitsch return 0; 1607d26aeb3SVaclav Hapla } 1617d26aeb3SVaclav Hapla 1627d26aeb3SVaclav Hapla /*TEST 1637d26aeb3SVaclav Hapla build: 1647d26aeb3SVaclav Hapla requires: hdf5 1657d26aeb3SVaclav Hapla testset: 1667d26aeb3SVaclav Hapla suffix: 0 1677d26aeb3SVaclav Hapla requires: !complex 168efa12513Sksagiyam nsize: 4 1697d26aeb3SVaclav Hapla args: -infile ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -informat hdf5_xdmf 1707d26aeb3SVaclav Hapla args: -outfile ex5_dump.h5 -outformat {{hdf5_xdmf hdf5_petsc}separate output} 1717d26aeb3SVaclav Hapla args: -ntimes 3 172efa12513Sksagiyam args: -loaded_dm_view -interpolated_dm_view -redistributed_dm_view 1737d26aeb3SVaclav Hapla test: 1747d26aeb3SVaclav Hapla # this partitioner should not shuffle anything, it should yield the same partititioning as the XDMF reader - added just for testing 1757d26aeb3SVaclav Hapla suffix: simple 1767d26aeb3SVaclav Hapla args: -petscpartitioner_type simple 1777d26aeb3SVaclav Hapla test: 1787d26aeb3SVaclav Hapla suffix: parmetis 1797d26aeb3SVaclav Hapla requires: parmetis 1807d26aeb3SVaclav Hapla args: -petscpartitioner_type parmetis 1817d26aeb3SVaclav Hapla test: 1827d26aeb3SVaclav Hapla suffix: ptscotch 1837d26aeb3SVaclav Hapla requires: ptscotch 1847d26aeb3SVaclav Hapla args: -petscpartitioner_type ptscotch 185efa12513Sksagiyam 186efa12513Sksagiyam testset: 187efa12513Sksagiyam suffix: 1 188efa12513Sksagiyam requires: !complex 189efa12513Sksagiyam nsize: 4 190efa12513Sksagiyam args: -infile ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -informat hdf5_xdmf 191efa12513Sksagiyam args: -outfile ex5_dump.h5 -outformat {{hdf5_xdmf hdf5_petsc}separate output} 192efa12513Sksagiyam args: -ntimes 3 193efa12513Sksagiyam args: -heterogeneous True 194efa12513Sksagiyam args: -loaded_dm_view -interpolated_dm_view -redistributed_dm_view 195efa12513Sksagiyam test: 196efa12513Sksagiyam suffix: simple 197efa12513Sksagiyam args: -petscpartitioner_type simple 198efa12513Sksagiyam test: 199efa12513Sksagiyam suffix: parmetis 200efa12513Sksagiyam requires: parmetis 201efa12513Sksagiyam args: -petscpartitioner_type parmetis 202efa12513Sksagiyam test: 203efa12513Sksagiyam suffix: ptscotch 204efa12513Sksagiyam requires: ptscotch 205efa12513Sksagiyam args: -petscpartitioner_type ptscotch 206efa12513Sksagiyam 2077d26aeb3SVaclav Hapla TEST*/ 208