1*472b9844SJames Wright #include "petscsf.h" 2*472b9844SJames Wright static char help[] = "Demonstrate CGNS parallel load-save-reload cycle, including data\n\n"; 3*472b9844SJames Wright 4*472b9844SJames Wright #include <petscdmplex.h> 5*472b9844SJames Wright #include <petscviewerhdf5.h> 6*472b9844SJames Wright #define EX "ex15.c" 7*472b9844SJames Wright 8*472b9844SJames Wright typedef struct { 9*472b9844SJames Wright char infile[PETSC_MAX_PATH_LEN]; /* Input mesh filename */ 10*472b9844SJames Wright char outfile[PETSC_MAX_PATH_LEN]; /* Dump/reload mesh filename */ 11*472b9844SJames Wright PetscBool heterogeneous; /* Test save on N / load on M */ 12*472b9844SJames Wright PetscInt ntimes; /* How many times do the cycle */ 13*472b9844SJames Wright } AppCtx; 14*472b9844SJames Wright 15*472b9844SJames Wright static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 16*472b9844SJames Wright { 17*472b9844SJames Wright PetscBool flg; 18*472b9844SJames Wright 19*472b9844SJames Wright PetscFunctionBeginUser; 20*472b9844SJames Wright options->infile[0] = '\0'; 21*472b9844SJames Wright options->outfile[0] = '\0'; 22*472b9844SJames Wright options->heterogeneous = PETSC_FALSE; 23*472b9844SJames Wright options->ntimes = 2; 24*472b9844SJames Wright PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX"); 25*472b9844SJames Wright PetscCall(PetscOptionsString("-infile", "The input CGNS file", EX, options->infile, options->infile, sizeof(options->infile), &flg)); 26*472b9844SJames Wright PetscCall(PetscOptionsString("-outfile", "The output CGNS file", EX, options->outfile, options->outfile, sizeof(options->outfile), &flg)); 27*472b9844SJames Wright PetscCall(PetscOptionsBool("-heterogeneous", "Test save on N / load on M", EX, options->heterogeneous, &options->heterogeneous, NULL)); 28*472b9844SJames Wright PetscCall(PetscOptionsInt("-ntimes", "How many times do the cycle", EX, options->ntimes, &options->ntimes, NULL)); 29*472b9844SJames Wright PetscOptionsEnd(); 30*472b9844SJames Wright PetscCheck(flg, comm, PETSC_ERR_USER_INPUT, "-infile needs to be specified"); 31*472b9844SJames Wright PetscCheck(flg, comm, PETSC_ERR_USER_INPUT, "-outfile needs to be specified"); 32*472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 33*472b9844SJames Wright } 34*472b9844SJames Wright 35*472b9844SJames Wright // @brief Create DM from CGNS file and setup PetscFE to VecLoad solution from that file 36*472b9844SJames Wright PetscErrorCode ReadCGNSDM(MPI_Comm comm, const char filename[], DM *dm) 37*472b9844SJames Wright { 38*472b9844SJames Wright PetscInt degree; 39*472b9844SJames Wright 40*472b9844SJames Wright PetscFunctionBeginUser; 41*472b9844SJames Wright PetscCall(DMPlexCreateFromFile(comm, filename, "ex15_plex", PETSC_TRUE, dm)); 42*472b9844SJames Wright PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 43*472b9844SJames Wright PetscCall(DMSetFromOptions(*dm)); 44*472b9844SJames Wright PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 45*472b9844SJames Wright 46*472b9844SJames Wright /* Redistribute */ 47*472b9844SJames Wright PetscCall(DMSetOptionsPrefix(*dm, "redistributed_")); 48*472b9844SJames Wright PetscCall(DMSetFromOptions(*dm)); 49*472b9844SJames Wright PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 50*472b9844SJames Wright 51*472b9844SJames Wright { // Get degree of the natural section 52*472b9844SJames Wright PetscFE fe_natural; 53*472b9844SJames Wright PetscDualSpace dual_space_natural; 54*472b9844SJames Wright 55*472b9844SJames Wright PetscCall(DMGetField(*dm, 0, NULL, (PetscObject *)&fe_natural)); 56*472b9844SJames Wright PetscCall(PetscFEGetDualSpace(fe_natural, &dual_space_natural)); 57*472b9844SJames Wright PetscCall(PetscDualSpaceGetOrder(dual_space_natural, °ree)); 58*472b9844SJames Wright PetscCall(DMClearFields(*dm)); 59*472b9844SJames Wright PetscCall(DMSetLocalSection(*dm, NULL)); 60*472b9844SJames Wright } 61*472b9844SJames Wright 62*472b9844SJames Wright { // Setup fe to load in the initial condition data 63*472b9844SJames Wright PetscFE fe; 64*472b9844SJames Wright PetscInt dim; 65*472b9844SJames Wright 66*472b9844SJames Wright PetscCall(DMGetDimension(*dm, &dim)); 67*472b9844SJames Wright PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 5, PETSC_FALSE, degree, PETSC_DETERMINE, &fe)); 68*472b9844SJames Wright PetscCall(PetscObjectSetName((PetscObject)fe, "FE for VecLoad")); 69*472b9844SJames Wright PetscCall(DMAddField(*dm, NULL, (PetscObject)fe)); 70*472b9844SJames Wright PetscCall(DMCreateDS(*dm)); 71*472b9844SJames Wright PetscCall(PetscFEDestroy(&fe)); 72*472b9844SJames Wright } 73*472b9844SJames Wright 74*472b9844SJames Wright // Set section component names, used when writing out CGNS files 75*472b9844SJames Wright PetscSection section; 76*472b9844SJames Wright PetscCall(DMGetLocalSection(*dm, §ion)); 77*472b9844SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 78*472b9844SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure")); 79*472b9844SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX")); 80*472b9844SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY")); 81*472b9844SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ")); 82*472b9844SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature")); 83*472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 84*472b9844SJames Wright } 85*472b9844SJames Wright 86*472b9844SJames Wright // Verify that V_load is equivalent to V_serial, even if V_load is distributed 87*472b9844SJames Wright PetscErrorCode VerifyLoadedSolution(DM dm_serial, Vec V_serial, DM dm_load, Vec V_load, PetscScalar tol) 88*472b9844SJames Wright { 89*472b9844SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm_load); 90*472b9844SJames Wright PetscSF load_to_serial_sf; 91*472b9844SJames Wright PetscScalar *array_load_bcast = NULL; 92*472b9844SJames Wright PetscInt num_comps = 5; 93*472b9844SJames Wright 94*472b9844SJames Wright PetscFunctionBeginUser; 95*472b9844SJames Wright { // Create SF to broadcast loaded vec nodes to serial vec nodes 96*472b9844SJames Wright PetscInt dim, num_local_serial = 0, num_local_load; 97*472b9844SJames Wright Vec coord_Vec_serial, coord_Vec_load; 98*472b9844SJames Wright const PetscScalar *coord_serial = NULL, *coord_load; 99*472b9844SJames Wright 100*472b9844SJames Wright PetscCall(DMGetCoordinateDim(dm_load, &dim)); 101*472b9844SJames Wright PetscCall(DMGetCoordinates(dm_load, &coord_Vec_load)); 102*472b9844SJames Wright PetscCall(VecGetLocalSize(coord_Vec_load, &num_local_load)); 103*472b9844SJames Wright num_local_load /= dim; 104*472b9844SJames Wright 105*472b9844SJames Wright PetscCall(VecGetArrayRead(coord_Vec_load, &coord_load)); 106*472b9844SJames Wright 107*472b9844SJames Wright if (dm_serial) { 108*472b9844SJames Wright PetscCall(DMGetCoordinates(dm_serial, &coord_Vec_serial)); 109*472b9844SJames Wright PetscCall(VecGetLocalSize(coord_Vec_serial, &num_local_serial)); 110*472b9844SJames Wright num_local_serial /= dim; 111*472b9844SJames Wright PetscCall(VecGetArrayRead(coord_Vec_serial, &coord_serial)); 112*472b9844SJames Wright } 113*472b9844SJames Wright 114*472b9844SJames Wright PetscCall(PetscSFCreate(comm, &load_to_serial_sf)); 115*472b9844SJames Wright PetscCall(PetscSFSetGraphFromCoordinates(load_to_serial_sf, num_local_load, num_local_serial, dim, 100 * PETSC_MACHINE_EPSILON, coord_load, coord_serial)); 116*472b9844SJames Wright PetscCall(PetscSFViewFromOptions(load_to_serial_sf, NULL, "-verify_sf_view")); 117*472b9844SJames Wright 118*472b9844SJames Wright PetscCall(VecRestoreArrayRead(coord_Vec_load, &coord_load)); 119*472b9844SJames Wright if (dm_serial) PetscCall(VecRestoreArrayRead(coord_Vec_serial, &coord_serial)); 120*472b9844SJames Wright } 121*472b9844SJames Wright 122*472b9844SJames Wright { // Broadcast the loaded vector data into a serial-sized array 123*472b9844SJames Wright PetscInt size_local_serial = 0; 124*472b9844SJames Wright const PetscScalar *array_load; 125*472b9844SJames Wright MPI_Datatype unit; 126*472b9844SJames Wright 127*472b9844SJames Wright PetscCall(VecGetArrayRead(V_load, &array_load)); 128*472b9844SJames Wright if (V_serial) { 129*472b9844SJames Wright PetscCall(VecGetLocalSize(V_serial, &size_local_serial)); 130*472b9844SJames Wright PetscCall(PetscMalloc1(size_local_serial, &array_load_bcast)); 131*472b9844SJames Wright } 132*472b9844SJames Wright 133*472b9844SJames Wright PetscCallMPI(MPI_Type_contiguous(num_comps, MPIU_REAL, &unit)); 134*472b9844SJames Wright PetscCallMPI(MPI_Type_commit(&unit)); 135*472b9844SJames Wright PetscCall(PetscSFBcastBegin(load_to_serial_sf, unit, array_load, array_load_bcast, MPI_REPLACE)); 136*472b9844SJames Wright PetscCall(PetscSFBcastEnd(load_to_serial_sf, unit, array_load, array_load_bcast, MPI_REPLACE)); 137*472b9844SJames Wright PetscCallMPI(MPI_Type_free(&unit)); 138*472b9844SJames Wright PetscCall(VecRestoreArrayRead(V_load, &array_load)); 139*472b9844SJames Wright } 140*472b9844SJames Wright 141*472b9844SJames Wright if (V_serial) { 142*472b9844SJames Wright const PetscScalar *array_serial; 143*472b9844SJames Wright PetscInt size_local_serial; 144*472b9844SJames Wright 145*472b9844SJames Wright PetscCall(VecGetArrayRead(V_serial, &array_serial)); 146*472b9844SJames Wright PetscCall(VecGetLocalSize(V_serial, &size_local_serial)); 147*472b9844SJames Wright for (PetscInt i = 0; i < size_local_serial; i++) { 148*472b9844SJames Wright if (PetscAbs(array_serial[i] - array_load_bcast[i]) > tol) PetscCall(PetscPrintf(comm, "DoF %" PetscInt_FMT " is inconsistent. Found %g, expected %g\n", i, array_load_bcast[i], array_serial[i])); 149*472b9844SJames Wright } 150*472b9844SJames Wright PetscCall(VecRestoreArrayRead(V_serial, &array_serial)); 151*472b9844SJames Wright } 152*472b9844SJames Wright 153*472b9844SJames Wright PetscCall(PetscFree(array_load_bcast)); 154*472b9844SJames Wright PetscCall(PetscSFDestroy(&load_to_serial_sf)); 155*472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 156*472b9844SJames Wright } 157*472b9844SJames Wright 158*472b9844SJames Wright int main(int argc, char **argv) 159*472b9844SJames Wright { 160*472b9844SJames Wright AppCtx user; 161*472b9844SJames Wright MPI_Comm comm; 162*472b9844SJames Wright PetscMPIInt gsize, grank, mycolor; 163*472b9844SJames Wright PetscBool flg; 164*472b9844SJames Wright const char *infilename; 165*472b9844SJames Wright DM dm_serial = NULL; 166*472b9844SJames Wright Vec V_serial = NULL; 167*472b9844SJames Wright 168*472b9844SJames Wright PetscFunctionBeginUser; 169*472b9844SJames Wright PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 170*472b9844SJames Wright PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 171*472b9844SJames Wright PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &gsize)); 172*472b9844SJames Wright PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &grank)); 173*472b9844SJames Wright 174*472b9844SJames Wright { // Read infile in serial 175*472b9844SJames Wright PetscViewer viewer; 176*472b9844SJames Wright PetscMPIInt gsize_serial; 177*472b9844SJames Wright 178*472b9844SJames Wright mycolor = grank == 0 ? 0 : 1; 179*472b9844SJames Wright PetscCallMPI(MPI_Comm_split(PETSC_COMM_WORLD, mycolor, grank, &comm)); 180*472b9844SJames Wright 181*472b9844SJames Wright if (grank == 0) { 182*472b9844SJames Wright PetscCallMPI(MPI_Comm_size(comm, &gsize_serial)); 183*472b9844SJames Wright 184*472b9844SJames Wright PetscCall(ReadCGNSDM(comm, user.infile, &dm_serial)); 185*472b9844SJames Wright PetscCall(DMSetOptionsPrefix(dm_serial, "serial_")); 186*472b9844SJames Wright 187*472b9844SJames Wright /* We just test/demonstrate DM is indeed distributed - unneeded in the application code */ 188*472b9844SJames Wright PetscCall(DMPlexIsDistributed(dm_serial, &flg)); 189*472b9844SJames Wright PetscCall(PetscPrintf(comm, "Loaded mesh distributed? %s\n", PetscBools[flg])); 190*472b9844SJames Wright 191*472b9844SJames Wright PetscCall(DMViewFromOptions(dm_serial, NULL, "-dm_view")); 192*472b9844SJames Wright PetscCall(PetscViewerCGNSOpen(comm, user.infile, FILE_MODE_READ, &viewer)); 193*472b9844SJames Wright PetscCall(DMGetGlobalVector(dm_serial, &V_serial)); 194*472b9844SJames Wright PetscCall(VecLoad(V_serial, viewer)); 195*472b9844SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 196*472b9844SJames Wright PetscCallMPI(MPI_Comm_free(&comm)); 197*472b9844SJames Wright } 198*472b9844SJames Wright PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 199*472b9844SJames Wright } 200*472b9844SJames Wright 201*472b9844SJames Wright for (PetscInt i = 0; i < user.ntimes; i++) { 202*472b9844SJames Wright if (i == 0) { 203*472b9844SJames Wright /* Use infile for the initial load */ 204*472b9844SJames Wright infilename = user.infile; 205*472b9844SJames Wright } else { 206*472b9844SJames Wright /* Use outfile for all I/O except the very initial load */ 207*472b9844SJames Wright infilename = user.outfile; 208*472b9844SJames Wright } 209*472b9844SJames Wright 210*472b9844SJames Wright if (user.heterogeneous) { 211*472b9844SJames Wright mycolor = (PetscMPIInt)(grank > user.ntimes - i); 212*472b9844SJames Wright } else { 213*472b9844SJames Wright mycolor = (PetscMPIInt)0; 214*472b9844SJames Wright } 215*472b9844SJames Wright PetscCallMPI(MPI_Comm_split(PETSC_COMM_WORLD, mycolor, grank, &comm)); 216*472b9844SJames Wright 217*472b9844SJames Wright if (mycolor == 0) { 218*472b9844SJames Wright /* Load/Save only on processes with mycolor == 0 */ 219*472b9844SJames Wright DM dm; 220*472b9844SJames Wright Vec V; 221*472b9844SJames Wright PetscViewer viewer; 222*472b9844SJames Wright PetscMPIInt comm_size; 223*472b9844SJames Wright const char *name; 224*472b9844SJames Wright PetscReal time; 225*472b9844SJames Wright PetscBool set; 226*472b9844SJames Wright 227*472b9844SJames Wright PetscCallMPI(MPI_Comm_size(comm, &comm_size)); 228*472b9844SJames Wright PetscCall(PetscPrintf(comm, "Begin cycle %" PetscInt_FMT ", comm size %" PetscInt_FMT "\n", i, comm_size)); 229*472b9844SJames Wright 230*472b9844SJames Wright // Load DM from CGNS file 231*472b9844SJames Wright PetscCall(ReadCGNSDM(comm, infilename, &dm)); 232*472b9844SJames Wright PetscCall(DMSetOptionsPrefix(dm, "loaded_")); 233*472b9844SJames Wright PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 234*472b9844SJames Wright 235*472b9844SJames Wright // Load solution from CGNS file 236*472b9844SJames Wright PetscCall(PetscViewerCGNSOpen(comm, infilename, FILE_MODE_READ, &viewer)); 237*472b9844SJames Wright PetscCall(DMGetGlobalVector(dm, &V)); 238*472b9844SJames Wright PetscCall(PetscViewerCGNSSetSolutionIndex(viewer, 1)); 239*472b9844SJames Wright { // Test GetSolutionIndex, not needed in application code 240*472b9844SJames Wright PetscInt solution_index; 241*472b9844SJames Wright PetscCall(PetscViewerCGNSGetSolutionIndex(viewer, &solution_index)); 242*472b9844SJames Wright PetscCheck(solution_index == 1, comm, PETSC_ERR_ARG_INCOMP, "Returned solution index wrong."); 243*472b9844SJames Wright } 244*472b9844SJames Wright PetscCall(PetscViewerCGNSGetSolutionName(viewer, &name)); 245*472b9844SJames Wright PetscCall(PetscViewerCGNSGetSolutionTime(viewer, &time, &set)); 246*472b9844SJames Wright PetscCheck(set, comm, PETSC_ERR_RETURN, "Time data wasn't set!"); 247*472b9844SJames Wright PetscCall(PetscPrintf(comm, "Solution Name: %s, and time %g\n", name, time)); 248*472b9844SJames Wright PetscCall(VecLoad(V, viewer)); 249*472b9844SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 250*472b9844SJames Wright 251*472b9844SJames Wright // Verify loaded solution against serial solution 252*472b9844SJames Wright PetscCall(VerifyLoadedSolution(dm_serial, V_serial, dm, V, 100 * PETSC_MACHINE_EPSILON)); 253*472b9844SJames Wright 254*472b9844SJames Wright // Write loaded solution to CGNS file 255*472b9844SJames Wright PetscCall(PetscViewerCGNSOpen(comm, user.outfile, FILE_MODE_WRITE, &viewer)); 256*472b9844SJames Wright PetscCall(VecView(V, viewer)); 257*472b9844SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 258*472b9844SJames Wright 259*472b9844SJames Wright PetscCall(DMRestoreGlobalVector(dm, &V)); 260*472b9844SJames Wright PetscCall(DMDestroy(&dm)); 261*472b9844SJames Wright PetscCall(PetscPrintf(comm, "End cycle %" PetscInt_FMT "\n--------\n", i)); 262*472b9844SJames Wright } 263*472b9844SJames Wright PetscCallMPI(MPI_Comm_free(&comm)); 264*472b9844SJames Wright PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 265*472b9844SJames Wright } 266*472b9844SJames Wright 267*472b9844SJames Wright if (V_serial && dm_serial) PetscCall(DMRestoreGlobalVector(dm_serial, &V_serial)); 268*472b9844SJames Wright if (dm_serial) PetscCall(DMDestroy(&dm_serial)); 269*472b9844SJames Wright PetscCall(PetscFinalize()); 270*472b9844SJames Wright return 0; 271*472b9844SJames Wright } 272*472b9844SJames Wright 273*472b9844SJames Wright /*TEST 274*472b9844SJames Wright build: 275*472b9844SJames Wright requires: cgns 276*472b9844SJames Wright testset: 277*472b9844SJames Wright suffix: cgns 278*472b9844SJames Wright requires: !complex 279*472b9844SJames Wright nsize: 4 280*472b9844SJames Wright args: -infile ${wPETSC_DIR}/share/petsc/datafiles/meshes/2x2x2_Q3_wave.cgns -outfile 2x2x2_Q3_wave_output.cgns 281*472b9844SJames Wright args: -dm_plex_cgns_parallel -ntimes 3 -heterogeneous -serial_dm_view -loaded_dm_view 282*472b9844SJames Wright test: 283*472b9844SJames Wright # this partitioner should not shuffle anything, it should yield the same partitioning as the XDMF reader - added just for testing 284*472b9844SJames Wright suffix: simple 285*472b9844SJames Wright args: -petscpartitioner_type simple 286*472b9844SJames Wright test: 287*472b9844SJames Wright suffix: parmetis 288*472b9844SJames Wright requires: parmetis 289*472b9844SJames Wright args: -petscpartitioner_type parmetis 290*472b9844SJames Wright test: 291*472b9844SJames Wright suffix: ptscotch 292*472b9844SJames Wright requires: ptscotch 293*472b9844SJames Wright args: -petscpartitioner_type ptscotch 294*472b9844SJames Wright 295*472b9844SJames Wright TEST*/ 296