1*d2522c19Sksagiyam static char help[] = "Tests save/load plex with distribution in HDF5.\n\n"; 2*d2522c19Sksagiyam 3*d2522c19Sksagiyam #include <petscdmplex.h> 4*d2522c19Sksagiyam #include <petscsf.h> 5*d2522c19Sksagiyam #include <petsclayouthdf5.h> 6*d2522c19Sksagiyam 7*d2522c19Sksagiyam typedef struct { 8*d2522c19Sksagiyam char fname[PETSC_MAX_PATH_LEN]; /* Output mesh filename */ 9*d2522c19Sksagiyam } AppCtx; 10*d2522c19Sksagiyam 11*d2522c19Sksagiyam PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 12*d2522c19Sksagiyam PetscBool flg; 13*d2522c19Sksagiyam 14*d2522c19Sksagiyam PetscFunctionBegin; 15*d2522c19Sksagiyam options->fname[0] = '\0'; 16*d2522c19Sksagiyam PetscOptionsBegin(comm, "", "DMPlex View/Load Test Options", "DMPLEX"); 17*d2522c19Sksagiyam PetscCall(PetscOptionsString("-fname", "The output mesh file", "ex51.c", options->fname, options->fname, sizeof(options->fname), &flg)); 18*d2522c19Sksagiyam PetscOptionsEnd(); 19*d2522c19Sksagiyam PetscFunctionReturn(0); 20*d2522c19Sksagiyam } 21*d2522c19Sksagiyam 22*d2522c19Sksagiyam int main(int argc, char **argv) { 23*d2522c19Sksagiyam const char exampleDMPlexName[] = "exampleDMPlex"; 24*d2522c19Sksagiyam const char exampleDistributionName[] = "exampleDistribution"; 25*d2522c19Sksagiyam PetscViewerFormat format = PETSC_VIEWER_HDF5_PETSC; 26*d2522c19Sksagiyam AppCtx user; 27*d2522c19Sksagiyam 28*d2522c19Sksagiyam PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 29*d2522c19Sksagiyam PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 30*d2522c19Sksagiyam /* Save */ 31*d2522c19Sksagiyam { 32*d2522c19Sksagiyam DM dm; 33*d2522c19Sksagiyam PetscViewer viewer; 34*d2522c19Sksagiyam 35*d2522c19Sksagiyam PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, user.fname, FILE_MODE_WRITE, &viewer)); 36*d2522c19Sksagiyam /* Save exampleDMPlex */ 37*d2522c19Sksagiyam { 38*d2522c19Sksagiyam DM pdm; 39*d2522c19Sksagiyam PetscPartitioner part; 40*d2522c19Sksagiyam const PetscInt faces[2] = {6, 1}; 41*d2522c19Sksagiyam PetscSF sf; 42*d2522c19Sksagiyam PetscInt overlap = 1; 43*d2522c19Sksagiyam 44*d2522c19Sksagiyam PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, 2, PETSC_FALSE, faces, NULL, NULL, NULL, PETSC_TRUE, &dm)); 45*d2522c19Sksagiyam PetscCall(DMPlexGetPartitioner(dm, &part)); 46*d2522c19Sksagiyam PetscCall(PetscPartitionerSetFromOptions(part)); 47*d2522c19Sksagiyam PetscCall(DMPlexDistribute(dm, overlap, &sf, &pdm)); 48*d2522c19Sksagiyam if (pdm) { 49*d2522c19Sksagiyam PetscCall(DMDestroy(&dm)); 50*d2522c19Sksagiyam dm = pdm; 51*d2522c19Sksagiyam } 52*d2522c19Sksagiyam PetscCall(PetscSFDestroy(&sf)); 53*d2522c19Sksagiyam PetscCall(PetscObjectSetName((PetscObject)dm, exampleDMPlexName)); 54*d2522c19Sksagiyam PetscCall(DMPlexDistributionSetName(dm, exampleDistributionName)); 55*d2522c19Sksagiyam PetscCall(PetscViewerPushFormat(viewer, format)); 56*d2522c19Sksagiyam PetscCall(DMPlexTopologyView(dm, viewer)); 57*d2522c19Sksagiyam PetscCall(DMPlexLabelsView(dm, viewer)); 58*d2522c19Sksagiyam PetscCall(PetscViewerPopFormat(viewer)); 59*d2522c19Sksagiyam } 60*d2522c19Sksagiyam /* Save coordinates */ 61*d2522c19Sksagiyam PetscCall(PetscViewerPushFormat(viewer, format)); 62*d2522c19Sksagiyam PetscCall(DMPlexCoordinatesView(dm, viewer)); 63*d2522c19Sksagiyam PetscCall(PetscViewerPopFormat(viewer)); 64*d2522c19Sksagiyam PetscCall(PetscObjectSetName((PetscObject)dm, "Save: DM")); 65*d2522c19Sksagiyam PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 66*d2522c19Sksagiyam PetscCall(DMDestroy(&dm)); 67*d2522c19Sksagiyam PetscCall(PetscViewerDestroy(&viewer)); 68*d2522c19Sksagiyam } 69*d2522c19Sksagiyam /* Load */ 70*d2522c19Sksagiyam { 71*d2522c19Sksagiyam DM dm; 72*d2522c19Sksagiyam PetscSF sfXC; 73*d2522c19Sksagiyam PetscViewer viewer; 74*d2522c19Sksagiyam 75*d2522c19Sksagiyam PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, user.fname, FILE_MODE_READ, &viewer)); 76*d2522c19Sksagiyam /* Load exampleDMPlex */ 77*d2522c19Sksagiyam PetscCall(DMCreate(PETSC_COMM_WORLD, &dm)); 78*d2522c19Sksagiyam PetscCall(DMSetType(dm, DMPLEX)); 79*d2522c19Sksagiyam PetscCall(PetscObjectSetName((PetscObject)dm, exampleDMPlexName)); 80*d2522c19Sksagiyam PetscCall(DMPlexDistributionSetName(dm, exampleDistributionName)); 81*d2522c19Sksagiyam /* sfXC: X -> C */ 82*d2522c19Sksagiyam /* X: set of globalPointNumbers, [0, N) */ 83*d2522c19Sksagiyam /* C: loaded in-memory plex */ 84*d2522c19Sksagiyam PetscCall(PetscViewerPushFormat(viewer, format)); 85*d2522c19Sksagiyam PetscCall(DMPlexTopologyLoad(dm, viewer, &sfXC)); 86*d2522c19Sksagiyam PetscCall(PetscViewerPopFormat(viewer)); 87*d2522c19Sksagiyam /* Do not distribute (Already distributed just like the saved plex) */ 88*d2522c19Sksagiyam /* Load labels */ 89*d2522c19Sksagiyam PetscCall(DMPlexLabelsLoad(dm, viewer, sfXC)); 90*d2522c19Sksagiyam /* Load coordinates */ 91*d2522c19Sksagiyam PetscCall(PetscViewerPushFormat(viewer, format)); 92*d2522c19Sksagiyam PetscCall(DMPlexCoordinatesLoad(dm, viewer, sfXC)); 93*d2522c19Sksagiyam PetscCall(PetscViewerPopFormat(viewer)); 94*d2522c19Sksagiyam PetscCall(DMSetFromOptions(dm)); 95*d2522c19Sksagiyam /* Print the exact same plex as the saved one */ 96*d2522c19Sksagiyam PetscCall(PetscObjectSetName((PetscObject)dm, "Load: DM")); 97*d2522c19Sksagiyam PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 98*d2522c19Sksagiyam PetscCall(PetscSFDestroy(&sfXC)); 99*d2522c19Sksagiyam PetscCall(DMDestroy(&dm)); 100*d2522c19Sksagiyam PetscCall(PetscViewerDestroy(&viewer)); 101*d2522c19Sksagiyam } 102*d2522c19Sksagiyam /* Finalize */ 103*d2522c19Sksagiyam PetscCall(PetscFinalize()); 104*d2522c19Sksagiyam return 0; 105*d2522c19Sksagiyam } 106*d2522c19Sksagiyam 107*d2522c19Sksagiyam /*TEST 108*d2522c19Sksagiyam 109*d2522c19Sksagiyam build: 110*d2522c19Sksagiyam requires: hdf5 111*d2522c19Sksagiyam requires: !complex 112*d2522c19Sksagiyam test: 113*d2522c19Sksagiyam suffix: 0 114*d2522c19Sksagiyam requires: parmetis 115*d2522c19Sksagiyam nsize: {{1 2 4}separate output} 116*d2522c19Sksagiyam args: -fname ex51_dump.h5 -dm_view ascii::ascii_info_detail 117*d2522c19Sksagiyam args: -petscpartitioner_type parmetis 118*d2522c19Sksagiyam args: -dm_plex_view_hdf5_storage_version 2.0.0 119*d2522c19Sksagiyam 120*d2522c19Sksagiyam TEST*/ 121