1472b9844SJames Wright #include "petscsf.h" 2472b9844SJames Wright static char help[] = "Demonstrate CGNS parallel load-save-reload cycle, including data\n\n"; 3472b9844SJames Wright 4472b9844SJames Wright #include <petscdmplex.h> 5472b9844SJames Wright #include <petscviewerhdf5.h> 6472b9844SJames Wright #define EX "ex15.c" 7472b9844SJames Wright 8472b9844SJames Wright typedef struct { 9472b9844SJames Wright char infile[PETSC_MAX_PATH_LEN]; /* Input mesh filename */ 10472b9844SJames Wright char outfile[PETSC_MAX_PATH_LEN]; /* Dump/reload mesh filename */ 11472b9844SJames Wright PetscBool heterogeneous; /* Test save on N / load on M */ 12472b9844SJames Wright PetscInt ntimes; /* How many times do the cycle */ 13472b9844SJames Wright } AppCtx; 14472b9844SJames Wright 15472b9844SJames Wright static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 16472b9844SJames Wright { 17472b9844SJames Wright PetscBool flg; 18472b9844SJames Wright 19472b9844SJames Wright PetscFunctionBeginUser; 20472b9844SJames Wright options->infile[0] = '\0'; 21472b9844SJames Wright options->outfile[0] = '\0'; 22472b9844SJames Wright options->heterogeneous = PETSC_FALSE; 23472b9844SJames Wright options->ntimes = 2; 24472b9844SJames Wright PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX"); 25472b9844SJames Wright PetscCall(PetscOptionsString("-infile", "The input CGNS file", EX, options->infile, options->infile, sizeof(options->infile), &flg)); 26472b9844SJames Wright PetscCall(PetscOptionsString("-outfile", "The output CGNS file", EX, options->outfile, options->outfile, sizeof(options->outfile), &flg)); 27472b9844SJames Wright PetscCall(PetscOptionsBool("-heterogeneous", "Test save on N / load on M", EX, options->heterogeneous, &options->heterogeneous, NULL)); 28472b9844SJames Wright PetscCall(PetscOptionsInt("-ntimes", "How many times do the cycle", EX, options->ntimes, &options->ntimes, NULL)); 29472b9844SJames Wright PetscOptionsEnd(); 30472b9844SJames Wright PetscCheck(flg, comm, PETSC_ERR_USER_INPUT, "-infile needs to be specified"); 31472b9844SJames Wright PetscCheck(flg, comm, PETSC_ERR_USER_INPUT, "-outfile needs to be specified"); 32472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 33472b9844SJames Wright } 34472b9844SJames Wright 35472b9844SJames Wright // @brief Create DM from CGNS file and setup PetscFE to VecLoad solution from that file 36472b9844SJames Wright PetscErrorCode ReadCGNSDM(MPI_Comm comm, const char filename[], DM *dm) 37472b9844SJames Wright { 38472b9844SJames Wright PetscInt degree; 39472b9844SJames Wright 40472b9844SJames Wright PetscFunctionBeginUser; 41472b9844SJames Wright PetscCall(DMPlexCreateFromFile(comm, filename, "ex15_plex", PETSC_TRUE, dm)); 42472b9844SJames Wright PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 43472b9844SJames Wright PetscCall(DMSetFromOptions(*dm)); 44472b9844SJames Wright PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 45472b9844SJames Wright 46472b9844SJames Wright /* Redistribute */ 47472b9844SJames Wright PetscCall(DMSetOptionsPrefix(*dm, "redistributed_")); 48472b9844SJames Wright PetscCall(DMSetFromOptions(*dm)); 49472b9844SJames Wright PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 50472b9844SJames Wright 51472b9844SJames Wright { // Get degree of the natural section 52472b9844SJames Wright PetscFE fe_natural; 53472b9844SJames Wright PetscDualSpace dual_space_natural; 54472b9844SJames Wright 55472b9844SJames Wright PetscCall(DMGetField(*dm, 0, NULL, (PetscObject *)&fe_natural)); 56472b9844SJames Wright PetscCall(PetscFEGetDualSpace(fe_natural, &dual_space_natural)); 57472b9844SJames Wright PetscCall(PetscDualSpaceGetOrder(dual_space_natural, °ree)); 58472b9844SJames Wright PetscCall(DMClearFields(*dm)); 59472b9844SJames Wright PetscCall(DMSetLocalSection(*dm, NULL)); 60472b9844SJames Wright } 61472b9844SJames Wright 62472b9844SJames Wright { // Setup fe to load in the initial condition data 63472b9844SJames Wright PetscFE fe; 64472b9844SJames Wright PetscInt dim; 65472b9844SJames Wright 66472b9844SJames Wright PetscCall(DMGetDimension(*dm, &dim)); 67472b9844SJames Wright PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 5, PETSC_FALSE, degree, PETSC_DETERMINE, &fe)); 68472b9844SJames Wright PetscCall(PetscObjectSetName((PetscObject)fe, "FE for VecLoad")); 69472b9844SJames Wright PetscCall(DMAddField(*dm, NULL, (PetscObject)fe)); 70472b9844SJames Wright PetscCall(DMCreateDS(*dm)); 71472b9844SJames Wright PetscCall(PetscFEDestroy(&fe)); 72472b9844SJames Wright } 73472b9844SJames Wright 74472b9844SJames Wright // Set section component names, used when writing out CGNS files 75472b9844SJames Wright PetscSection section; 76472b9844SJames Wright PetscCall(DMGetLocalSection(*dm, §ion)); 77472b9844SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 78472b9844SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure")); 79472b9844SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX")); 80472b9844SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY")); 81472b9844SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ")); 82472b9844SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature")); 83472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 84472b9844SJames Wright } 85472b9844SJames Wright 86472b9844SJames Wright // Verify that V_load is equivalent to V_serial, even if V_load is distributed 87472b9844SJames Wright PetscErrorCode VerifyLoadedSolution(DM dm_serial, Vec V_serial, DM dm_load, Vec V_load, PetscScalar tol) 88472b9844SJames Wright { 89472b9844SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm_load); 90472b9844SJames Wright PetscSF load_to_serial_sf; 91472b9844SJames Wright PetscScalar *array_load_bcast = NULL; 92472b9844SJames Wright PetscInt num_comps = 5; 93472b9844SJames Wright 94472b9844SJames Wright PetscFunctionBeginUser; 95472b9844SJames Wright { // Create SF to broadcast loaded vec nodes to serial vec nodes 96472b9844SJames Wright PetscInt dim, num_local_serial = 0, num_local_load; 97472b9844SJames Wright Vec coord_Vec_serial, coord_Vec_load; 98472b9844SJames Wright const PetscScalar *coord_serial = NULL, *coord_load; 99472b9844SJames Wright 100472b9844SJames Wright PetscCall(DMGetCoordinateDim(dm_load, &dim)); 101472b9844SJames Wright PetscCall(DMGetCoordinates(dm_load, &coord_Vec_load)); 102472b9844SJames Wright PetscCall(VecGetLocalSize(coord_Vec_load, &num_local_load)); 103472b9844SJames Wright num_local_load /= dim; 104472b9844SJames Wright 105472b9844SJames Wright PetscCall(VecGetArrayRead(coord_Vec_load, &coord_load)); 106472b9844SJames Wright 107472b9844SJames Wright if (dm_serial) { 108472b9844SJames Wright PetscCall(DMGetCoordinates(dm_serial, &coord_Vec_serial)); 109472b9844SJames Wright PetscCall(VecGetLocalSize(coord_Vec_serial, &num_local_serial)); 110472b9844SJames Wright num_local_serial /= dim; 111472b9844SJames Wright PetscCall(VecGetArrayRead(coord_Vec_serial, &coord_serial)); 112472b9844SJames Wright } 113472b9844SJames Wright 114472b9844SJames Wright PetscCall(PetscSFCreate(comm, &load_to_serial_sf)); 115472b9844SJames Wright PetscCall(PetscSFSetGraphFromCoordinates(load_to_serial_sf, num_local_load, num_local_serial, dim, 100 * PETSC_MACHINE_EPSILON, coord_load, coord_serial)); 116472b9844SJames Wright PetscCall(PetscSFViewFromOptions(load_to_serial_sf, NULL, "-verify_sf_view")); 117472b9844SJames Wright 118472b9844SJames Wright PetscCall(VecRestoreArrayRead(coord_Vec_load, &coord_load)); 119472b9844SJames Wright if (dm_serial) PetscCall(VecRestoreArrayRead(coord_Vec_serial, &coord_serial)); 120472b9844SJames Wright } 121472b9844SJames Wright 122472b9844SJames Wright { // Broadcast the loaded vector data into a serial-sized array 123472b9844SJames Wright PetscInt size_local_serial = 0; 124472b9844SJames Wright const PetscScalar *array_load; 125472b9844SJames Wright MPI_Datatype unit; 126472b9844SJames Wright 127472b9844SJames Wright PetscCall(VecGetArrayRead(V_load, &array_load)); 128472b9844SJames Wright if (V_serial) { 129472b9844SJames Wright PetscCall(VecGetLocalSize(V_serial, &size_local_serial)); 130472b9844SJames Wright PetscCall(PetscMalloc1(size_local_serial, &array_load_bcast)); 131472b9844SJames Wright } 132472b9844SJames Wright 133472b9844SJames Wright PetscCallMPI(MPI_Type_contiguous(num_comps, MPIU_REAL, &unit)); 134472b9844SJames Wright PetscCallMPI(MPI_Type_commit(&unit)); 135472b9844SJames Wright PetscCall(PetscSFBcastBegin(load_to_serial_sf, unit, array_load, array_load_bcast, MPI_REPLACE)); 136472b9844SJames Wright PetscCall(PetscSFBcastEnd(load_to_serial_sf, unit, array_load, array_load_bcast, MPI_REPLACE)); 137472b9844SJames Wright PetscCallMPI(MPI_Type_free(&unit)); 138472b9844SJames Wright PetscCall(VecRestoreArrayRead(V_load, &array_load)); 139472b9844SJames Wright } 140472b9844SJames Wright 141472b9844SJames Wright if (V_serial) { 142472b9844SJames Wright const PetscScalar *array_serial; 143472b9844SJames Wright PetscInt size_local_serial; 144472b9844SJames Wright 145472b9844SJames Wright PetscCall(VecGetArrayRead(V_serial, &array_serial)); 146472b9844SJames Wright PetscCall(VecGetLocalSize(V_serial, &size_local_serial)); 147472b9844SJames Wright for (PetscInt i = 0; i < size_local_serial; i++) { 148472b9844SJames 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])); 149472b9844SJames Wright } 150472b9844SJames Wright PetscCall(VecRestoreArrayRead(V_serial, &array_serial)); 151472b9844SJames Wright } 152472b9844SJames Wright 153472b9844SJames Wright PetscCall(PetscFree(array_load_bcast)); 154472b9844SJames Wright PetscCall(PetscSFDestroy(&load_to_serial_sf)); 155472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 156472b9844SJames Wright } 157472b9844SJames Wright 158472b9844SJames Wright int main(int argc, char **argv) 159472b9844SJames Wright { 160472b9844SJames Wright AppCtx user; 161472b9844SJames Wright MPI_Comm comm; 162472b9844SJames Wright PetscMPIInt gsize, grank, mycolor; 163472b9844SJames Wright PetscBool flg; 164472b9844SJames Wright const char *infilename; 165472b9844SJames Wright DM dm_serial = NULL; 166472b9844SJames Wright Vec V_serial = NULL; 167472b9844SJames Wright 168472b9844SJames Wright PetscFunctionBeginUser; 169472b9844SJames Wright PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 170472b9844SJames Wright PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 171472b9844SJames Wright PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &gsize)); 172472b9844SJames Wright PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &grank)); 173472b9844SJames Wright 174472b9844SJames Wright { // Read infile in serial 175472b9844SJames Wright PetscViewer viewer; 176472b9844SJames Wright PetscMPIInt gsize_serial; 177472b9844SJames Wright 178472b9844SJames Wright mycolor = grank == 0 ? 0 : 1; 179472b9844SJames Wright PetscCallMPI(MPI_Comm_split(PETSC_COMM_WORLD, mycolor, grank, &comm)); 180472b9844SJames Wright 181472b9844SJames Wright if (grank == 0) { 182472b9844SJames Wright PetscCallMPI(MPI_Comm_size(comm, &gsize_serial)); 183472b9844SJames Wright 184472b9844SJames Wright PetscCall(ReadCGNSDM(comm, user.infile, &dm_serial)); 185472b9844SJames Wright PetscCall(DMSetOptionsPrefix(dm_serial, "serial_")); 186472b9844SJames Wright 187472b9844SJames Wright /* We just test/demonstrate DM is indeed distributed - unneeded in the application code */ 188472b9844SJames Wright PetscCall(DMPlexIsDistributed(dm_serial, &flg)); 189472b9844SJames Wright PetscCall(PetscPrintf(comm, "Loaded mesh distributed? %s\n", PetscBools[flg])); 190472b9844SJames Wright 191472b9844SJames Wright PetscCall(DMViewFromOptions(dm_serial, NULL, "-dm_view")); 192472b9844SJames Wright PetscCall(PetscViewerCGNSOpen(comm, user.infile, FILE_MODE_READ, &viewer)); 193472b9844SJames Wright PetscCall(DMGetGlobalVector(dm_serial, &V_serial)); 194472b9844SJames Wright PetscCall(VecLoad(V_serial, viewer)); 195472b9844SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 196472b9844SJames Wright PetscCallMPI(MPI_Comm_free(&comm)); 197472b9844SJames Wright } 198472b9844SJames Wright PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 199472b9844SJames Wright } 200472b9844SJames Wright 201472b9844SJames Wright for (PetscInt i = 0; i < user.ntimes; i++) { 202472b9844SJames Wright if (i == 0) { 203472b9844SJames Wright /* Use infile for the initial load */ 204472b9844SJames Wright infilename = user.infile; 205472b9844SJames Wright } else { 206472b9844SJames Wright /* Use outfile for all I/O except the very initial load */ 207472b9844SJames Wright infilename = user.outfile; 208472b9844SJames Wright } 209472b9844SJames Wright 210472b9844SJames Wright if (user.heterogeneous) { 211472b9844SJames Wright mycolor = (PetscMPIInt)(grank > user.ntimes - i); 212472b9844SJames Wright } else { 213472b9844SJames Wright mycolor = (PetscMPIInt)0; 214472b9844SJames Wright } 215472b9844SJames Wright PetscCallMPI(MPI_Comm_split(PETSC_COMM_WORLD, mycolor, grank, &comm)); 216472b9844SJames Wright 217472b9844SJames Wright if (mycolor == 0) { 218472b9844SJames Wright /* Load/Save only on processes with mycolor == 0 */ 219472b9844SJames Wright DM dm; 220472b9844SJames Wright Vec V; 221472b9844SJames Wright PetscViewer viewer; 222472b9844SJames Wright PetscMPIInt comm_size; 223472b9844SJames Wright const char *name; 224472b9844SJames Wright PetscReal time; 225472b9844SJames Wright PetscBool set; 226472b9844SJames Wright 227472b9844SJames Wright PetscCallMPI(MPI_Comm_size(comm, &comm_size)); 228*67f7d742SJames Wright PetscCall(PetscPrintf(comm, "Begin cycle %" PetscInt_FMT ", comm size %d\n", i, comm_size)); 229472b9844SJames Wright 230472b9844SJames Wright // Load DM from CGNS file 231472b9844SJames Wright PetscCall(ReadCGNSDM(comm, infilename, &dm)); 232472b9844SJames Wright PetscCall(DMSetOptionsPrefix(dm, "loaded_")); 233472b9844SJames Wright PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 234472b9844SJames Wright 235472b9844SJames Wright // Load solution from CGNS file 236472b9844SJames Wright PetscCall(PetscViewerCGNSOpen(comm, infilename, FILE_MODE_READ, &viewer)); 237472b9844SJames Wright PetscCall(DMGetGlobalVector(dm, &V)); 238472b9844SJames Wright PetscCall(PetscViewerCGNSSetSolutionIndex(viewer, 1)); 239472b9844SJames Wright { // Test GetSolutionIndex, not needed in application code 240472b9844SJames Wright PetscInt solution_index; 241472b9844SJames Wright PetscCall(PetscViewerCGNSGetSolutionIndex(viewer, &solution_index)); 242472b9844SJames Wright PetscCheck(solution_index == 1, comm, PETSC_ERR_ARG_INCOMP, "Returned solution index wrong."); 243472b9844SJames Wright } 244472b9844SJames Wright PetscCall(PetscViewerCGNSGetSolutionName(viewer, &name)); 245472b9844SJames Wright PetscCall(PetscViewerCGNSGetSolutionTime(viewer, &time, &set)); 246472b9844SJames Wright PetscCheck(set, comm, PETSC_ERR_RETURN, "Time data wasn't set!"); 247472b9844SJames Wright PetscCall(PetscPrintf(comm, "Solution Name: %s, and time %g\n", name, time)); 248472b9844SJames Wright PetscCall(VecLoad(V, viewer)); 249472b9844SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 250472b9844SJames Wright 251472b9844SJames Wright // Verify loaded solution against serial solution 252472b9844SJames Wright PetscCall(VerifyLoadedSolution(dm_serial, V_serial, dm, V, 100 * PETSC_MACHINE_EPSILON)); 253472b9844SJames Wright 254472b9844SJames Wright // Write loaded solution to CGNS file 255472b9844SJames Wright PetscCall(PetscViewerCGNSOpen(comm, user.outfile, FILE_MODE_WRITE, &viewer)); 256472b9844SJames Wright PetscCall(VecView(V, viewer)); 257472b9844SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 258472b9844SJames Wright 259472b9844SJames Wright PetscCall(DMRestoreGlobalVector(dm, &V)); 260472b9844SJames Wright PetscCall(DMDestroy(&dm)); 261472b9844SJames Wright PetscCall(PetscPrintf(comm, "End cycle %" PetscInt_FMT "\n--------\n", i)); 262472b9844SJames Wright } 263472b9844SJames Wright PetscCallMPI(MPI_Comm_free(&comm)); 264472b9844SJames Wright PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 265472b9844SJames Wright } 266472b9844SJames Wright 267472b9844SJames Wright if (V_serial && dm_serial) PetscCall(DMRestoreGlobalVector(dm_serial, &V_serial)); 268472b9844SJames Wright if (dm_serial) PetscCall(DMDestroy(&dm_serial)); 269472b9844SJames Wright PetscCall(PetscFinalize()); 270472b9844SJames Wright return 0; 271472b9844SJames Wright } 272472b9844SJames Wright 273472b9844SJames Wright /*TEST 274472b9844SJames Wright build: 275472b9844SJames Wright requires: cgns 276472b9844SJames Wright testset: 277472b9844SJames Wright suffix: cgns 278472b9844SJames Wright requires: !complex 279472b9844SJames Wright nsize: 4 280472b9844SJames Wright args: -infile ${wPETSC_DIR}/share/petsc/datafiles/meshes/2x2x2_Q3_wave.cgns -outfile 2x2x2_Q3_wave_output.cgns 281f819935bSJames Wright args: -dm_plex_cgns_parallel -ntimes 3 -heterogeneous -serial_dm_view -loaded_dm_view -redistributed_dm_distribute 282472b9844SJames Wright test: 283472b9844SJames Wright # this partitioner should not shuffle anything, it should yield the same partitioning as the XDMF reader - added just for testing 284472b9844SJames Wright suffix: simple 285472b9844SJames Wright args: -petscpartitioner_type simple 286472b9844SJames Wright test: 287472b9844SJames Wright suffix: parmetis 288472b9844SJames Wright requires: parmetis 289472b9844SJames Wright args: -petscpartitioner_type parmetis 290472b9844SJames Wright test: 291472b9844SJames Wright suffix: ptscotch 292472b9844SJames Wright requires: ptscotch 293472b9844SJames Wright args: -petscpartitioner_type ptscotch 294472b9844SJames Wright 295472b9844SJames Wright TEST*/ 296