1*90df3356SJames Wright #include "petscsf.h" 2*90df3356SJames Wright static char help[] = "Test CGNS writing output with isoperiodic boundaries\n\n"; 3*90df3356SJames Wright // Also tests DMSetCoordinateDisc() for isoperiodic boundaries and projection = true 4*90df3356SJames Wright 5*90df3356SJames Wright #include <petscdmplex.h> 6*90df3356SJames Wright #include <petscviewerhdf5.h> 7*90df3356SJames Wright #define EX "ex100.c" 8*90df3356SJames Wright 9*90df3356SJames Wright static PetscErrorCode project_function(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 10*90df3356SJames Wright { 11*90df3356SJames Wright PetscReal x_tot = 0; 12*90df3356SJames Wright 13*90df3356SJames Wright PetscFunctionBeginUser; 14*90df3356SJames Wright for (PetscInt d = 0; d < dim; d++) x_tot += x[d]; 15*90df3356SJames Wright for (PetscInt c = 0; c < Nc; c++) { u[c] = sin(2 * M_PI * x_tot); } 16*90df3356SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 17*90df3356SJames Wright } 18*90df3356SJames Wright 19*90df3356SJames Wright int main(int argc, char **argv) 20*90df3356SJames Wright { 21*90df3356SJames Wright DM dm_create, dm_read; 22*90df3356SJames Wright Vec V; 23*90df3356SJames Wright char file[PETSC_MAX_PATH_LEN]; 24*90df3356SJames Wright MPI_Comm comm; 25*90df3356SJames Wright PetscInt solution_degree = 2; 26*90df3356SJames Wright 27*90df3356SJames Wright PetscFunctionBeginUser; 28*90df3356SJames Wright PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 29*90df3356SJames Wright comm = PETSC_COMM_WORLD; 30*90df3356SJames Wright 31*90df3356SJames Wright PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX"); 32*90df3356SJames Wright PetscCall(PetscOptionsInt("-solution_degree", "The input CGNS file", EX, solution_degree, &solution_degree, NULL)); 33*90df3356SJames Wright PetscCall(PetscOptionsString("-file", "The input CGNS file", EX, file, file, sizeof(file), NULL)); 34*90df3356SJames Wright PetscOptionsEnd(); 35*90df3356SJames Wright 36*90df3356SJames Wright { // Create DM 37*90df3356SJames Wright PetscCall(DMCreate(comm, &dm_create)); 38*90df3356SJames Wright PetscCall(DMSetType(dm_create, DMPLEX)); 39*90df3356SJames Wright PetscCall(DMSetOptionsPrefix(dm_create, "create_")); 40*90df3356SJames Wright PetscCall(DMSetFromOptions(dm_create)); 41*90df3356SJames Wright PetscCall(DMViewFromOptions(dm_create, NULL, "-dm_create_view")); 42*90df3356SJames Wright 43*90df3356SJames Wright { // Setup fe to load in the initial condition data 44*90df3356SJames Wright PetscFE fe; 45*90df3356SJames Wright PetscInt dim; 46*90df3356SJames Wright 47*90df3356SJames Wright PetscCall(DMGetDimension(dm_create, &dim)); 48*90df3356SJames Wright PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, solution_degree, PETSC_DETERMINE, &fe)); 49*90df3356SJames Wright PetscCall(DMAddField(dm_create, NULL, (PetscObject)fe)); 50*90df3356SJames Wright PetscCall(DMCreateDS(dm_create)); 51*90df3356SJames Wright PetscCall(PetscFEDestroy(&fe)); 52*90df3356SJames Wright } 53*90df3356SJames Wright 54*90df3356SJames Wright PetscErrorCode (*funcs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) = {project_function}; 55*90df3356SJames Wright PetscCall(DMGetGlobalVector(dm_create, &V)); 56*90df3356SJames Wright PetscCall(DMProjectFunction(dm_create, 0, &funcs, NULL, INSERT_VALUES, V)); 57*90df3356SJames Wright PetscViewer viewer; 58*90df3356SJames Wright PetscCall(PetscViewerCGNSOpen(comm, file, FILE_MODE_WRITE, &viewer)); 59*90df3356SJames Wright PetscCall(VecView(V, viewer)); 60*90df3356SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 61*90df3356SJames Wright PetscCall(DMRestoreGlobalVector(dm_create, &V)); 62*90df3356SJames Wright 63*90df3356SJames Wright PetscCall(DMDestroy(&dm_create)); 64*90df3356SJames Wright } 65*90df3356SJames Wright 66*90df3356SJames Wright { 67*90df3356SJames Wright PetscSection coord_section; 68*90df3356SJames Wright PetscInt cStart, cEnd; 69*90df3356SJames Wright Vec coords; 70*90df3356SJames Wright PetscInt cdim; 71*90df3356SJames Wright PetscReal ref_cell_bounding_box_size[3]; 72*90df3356SJames Wright 73*90df3356SJames Wright PetscCall(DMPlexCreateFromFile(comm, file, "ex100_plex", PETSC_TRUE, &dm_read)); 74*90df3356SJames Wright PetscCall(DMViewFromOptions(dm_read, NULL, "-dm_read_view")); 75*90df3356SJames Wright PetscCall(DMGetCoordinateDim(dm_read, &cdim)); 76*90df3356SJames Wright PetscCall(DMGetCoordinateSection(dm_read, &coord_section)); 77*90df3356SJames Wright PetscCall(DMGetCoordinates(dm_read, &coords)); 78*90df3356SJames Wright PetscCall(DMPlexGetHeightStratum(dm_read, 0, &cStart, &cEnd)); 79*90df3356SJames Wright for (PetscInt cell = cStart; cell < cEnd; cell++) { 80*90df3356SJames Wright PetscInt num_closure = 0; 81*90df3356SJames Wright PetscScalar *cell_coords = NULL; 82*90df3356SJames Wright PetscReal cell_bounding_box_size[3], difference; 83*90df3356SJames Wright PetscReal min[3] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL}; 84*90df3356SJames Wright PetscReal max[3] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 85*90df3356SJames Wright 86*90df3356SJames Wright PetscCall(DMPlexVecGetClosure(dm_read, coord_section, coords, cell, &num_closure, &cell_coords)); 87*90df3356SJames Wright for (PetscInt n = 0; n < num_closure / cdim; n++) { 88*90df3356SJames Wright for (PetscInt d = 0; d < cdim; ++d) { 89*90df3356SJames Wright min[d] = PetscMin(min[d], PetscRealPart(cell_coords[n * cdim + d])); 90*90df3356SJames Wright max[d] = PetscMax(max[d], PetscRealPart(cell_coords[n * cdim + d])); 91*90df3356SJames Wright } 92*90df3356SJames Wright } 93*90df3356SJames Wright 94*90df3356SJames Wright for (PetscInt d = 0; d < cdim; d++) cell_bounding_box_size[d] = max[d] - min[d]; 95*90df3356SJames Wright if (cell == cStart) PetscCall(PetscArraycpy(ref_cell_bounding_box_size, cell_bounding_box_size, 3)); 96*90df3356SJames Wright 97*90df3356SJames Wright for (PetscInt d = 0; d < cdim; d++) { 98*90df3356SJames Wright difference = PetscAbsReal((ref_cell_bounding_box_size[d] - cell_bounding_box_size[d]) / ref_cell_bounding_box_size[d]); 99*90df3356SJames Wright if (difference > PETSC_MACHINE_EPSILON * 100) { 100*90df3356SJames Wright PetscPrintf(comm, "Cell %" PetscInt_FMT " doesn't match bounding box size of Cell %" PetscInt_FMT " in dimension %" PetscInt_FMT ". Relative difference = %g\n", cell, cStart, d, (double)difference); 101*90df3356SJames Wright } 102*90df3356SJames Wright } 103*90df3356SJames Wright PetscCall(DMPlexVecRestoreClosure(dm_read, coord_section, coords, cell, &num_closure, &cell_coords)); 104*90df3356SJames Wright } 105*90df3356SJames Wright PetscCall(DMDestroy(&dm_read)); 106*90df3356SJames Wright } 107*90df3356SJames Wright 108*90df3356SJames Wright PetscCall(PetscFinalize()); 109*90df3356SJames Wright return 0; 110*90df3356SJames Wright } 111*90df3356SJames Wright 112*90df3356SJames Wright /*TEST 113*90df3356SJames Wright build: 114*90df3356SJames Wright requires: cgns 115*90df3356SJames Wright test: 116*90df3356SJames Wright # Coordinates of dm_create are linear, but the CGNS writer projects them to quadratic to match the solution. 117*90df3356SJames Wright # The verification checks that all the cells of the resulting file are the same size 118*90df3356SJames Wright args: -create_dm_plex_shape zbox -create_dm_plex_box_faces 2,2 -create_dm_plex_box_bd periodic,periodic -file test.cgns -dm_plex_cgns_parallel 119*90df3356SJames Wright TEST*/ 120