190df3356SJames Wright #include "petscsf.h"
290df3356SJames Wright static char help[] = "Test CGNS writing output with isoperiodic boundaries\n\n";
390df3356SJames Wright // Also tests DMSetCoordinateDisc() for isoperiodic boundaries and projection = true
490df3356SJames Wright
590df3356SJames Wright #include <petscdmplex.h>
690df3356SJames Wright #include <petscviewerhdf5.h>
790df3356SJames Wright #define EX "ex100.c"
890df3356SJames Wright
project_function(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)9*2a8381b2SBarry Smith static PetscErrorCode project_function(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
1090df3356SJames Wright {
1190df3356SJames Wright PetscReal x_tot = 0;
1290df3356SJames Wright
1390df3356SJames Wright PetscFunctionBeginUser;
1490df3356SJames Wright for (PetscInt d = 0; d < dim; d++) x_tot += x[d];
15ac530a7eSPierre Jolivet for (PetscInt c = 0; c < Nc; c++) u[c] = sin(2 * M_PI * x_tot);
1690df3356SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
1790df3356SJames Wright }
1890df3356SJames Wright
main(int argc,char ** argv)1990df3356SJames Wright int main(int argc, char **argv)
2090df3356SJames Wright {
2190df3356SJames Wright DM dm_create, dm_read;
2290df3356SJames Wright Vec V;
2390df3356SJames Wright char file[PETSC_MAX_PATH_LEN];
2490df3356SJames Wright MPI_Comm comm;
2590df3356SJames Wright PetscInt solution_degree = 2;
2690df3356SJames Wright
2790df3356SJames Wright PetscFunctionBeginUser;
2890df3356SJames Wright PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2990df3356SJames Wright comm = PETSC_COMM_WORLD;
3090df3356SJames Wright
3190df3356SJames Wright PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
3290df3356SJames Wright PetscCall(PetscOptionsInt("-solution_degree", "The input CGNS file", EX, solution_degree, &solution_degree, NULL));
3390df3356SJames Wright PetscCall(PetscOptionsString("-file", "The input CGNS file", EX, file, file, sizeof(file), NULL));
3490df3356SJames Wright PetscOptionsEnd();
3590df3356SJames Wright
3690df3356SJames Wright { // Create DM
3790df3356SJames Wright PetscCall(DMCreate(comm, &dm_create));
3890df3356SJames Wright PetscCall(DMSetType(dm_create, DMPLEX));
3990df3356SJames Wright PetscCall(DMSetOptionsPrefix(dm_create, "create_"));
4090df3356SJames Wright PetscCall(DMSetFromOptions(dm_create));
4190df3356SJames Wright PetscCall(DMViewFromOptions(dm_create, NULL, "-dm_create_view"));
4290df3356SJames Wright
4390df3356SJames Wright { // Setup fe to load in the initial condition data
4490df3356SJames Wright PetscFE fe;
4590df3356SJames Wright PetscInt dim;
4690df3356SJames Wright
4790df3356SJames Wright PetscCall(DMGetDimension(dm_create, &dim));
4890df3356SJames Wright PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, solution_degree, PETSC_DETERMINE, &fe));
4990df3356SJames Wright PetscCall(DMAddField(dm_create, NULL, (PetscObject)fe));
5090df3356SJames Wright PetscCall(DMCreateDS(dm_create));
5190df3356SJames Wright PetscCall(PetscFEDestroy(&fe));
5290df3356SJames Wright }
5390df3356SJames Wright
54*2a8381b2SBarry Smith PetscErrorCode (*funcs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx) = {project_function};
5590df3356SJames Wright PetscCall(DMGetGlobalVector(dm_create, &V));
5690df3356SJames Wright PetscCall(DMProjectFunction(dm_create, 0, &funcs, NULL, INSERT_VALUES, V));
5790df3356SJames Wright PetscViewer viewer;
5890df3356SJames Wright PetscCall(PetscViewerCGNSOpen(comm, file, FILE_MODE_WRITE, &viewer));
5990df3356SJames Wright PetscCall(VecView(V, viewer));
6090df3356SJames Wright PetscCall(PetscViewerDestroy(&viewer));
6190df3356SJames Wright PetscCall(DMRestoreGlobalVector(dm_create, &V));
6290df3356SJames Wright
6390df3356SJames Wright PetscCall(DMDestroy(&dm_create));
6490df3356SJames Wright }
6590df3356SJames Wright
6690df3356SJames Wright {
6790df3356SJames Wright PetscSection coord_section;
6890df3356SJames Wright PetscInt cStart, cEnd;
6990df3356SJames Wright Vec coords;
7090df3356SJames Wright PetscInt cdim;
7190df3356SJames Wright PetscReal ref_cell_bounding_box_size[3];
7290df3356SJames Wright
7390df3356SJames Wright PetscCall(DMPlexCreateFromFile(comm, file, "ex100_plex", PETSC_TRUE, &dm_read));
7490df3356SJames Wright PetscCall(DMViewFromOptions(dm_read, NULL, "-dm_read_view"));
7590df3356SJames Wright PetscCall(DMGetCoordinateDim(dm_read, &cdim));
7690df3356SJames Wright PetscCall(DMGetCoordinateSection(dm_read, &coord_section));
7790df3356SJames Wright PetscCall(DMGetCoordinates(dm_read, &coords));
7890df3356SJames Wright PetscCall(DMPlexGetHeightStratum(dm_read, 0, &cStart, &cEnd));
7990df3356SJames Wright for (PetscInt cell = cStart; cell < cEnd; cell++) {
8090df3356SJames Wright PetscInt num_closure = 0;
8190df3356SJames Wright PetscScalar *cell_coords = NULL;
8290df3356SJames Wright PetscReal cell_bounding_box_size[3], difference;
8390df3356SJames Wright PetscReal min[3] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL};
8490df3356SJames Wright PetscReal max[3] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
8590df3356SJames Wright
8690df3356SJames Wright PetscCall(DMPlexVecGetClosure(dm_read, coord_section, coords, cell, &num_closure, &cell_coords));
8790df3356SJames Wright for (PetscInt n = 0; n < num_closure / cdim; n++) {
8890df3356SJames Wright for (PetscInt d = 0; d < cdim; ++d) {
8990df3356SJames Wright min[d] = PetscMin(min[d], PetscRealPart(cell_coords[n * cdim + d]));
9090df3356SJames Wright max[d] = PetscMax(max[d], PetscRealPart(cell_coords[n * cdim + d]));
9190df3356SJames Wright }
9290df3356SJames Wright }
9390df3356SJames Wright
9490df3356SJames Wright for (PetscInt d = 0; d < cdim; d++) cell_bounding_box_size[d] = max[d] - min[d];
9551af5445SJames Wright if (cell == cStart) PetscCall(PetscArraycpy(ref_cell_bounding_box_size, cell_bounding_box_size, cdim));
9690df3356SJames Wright
9790df3356SJames Wright for (PetscInt d = 0; d < cdim; d++) {
9890df3356SJames Wright difference = PetscAbsReal((ref_cell_bounding_box_size[d] - cell_bounding_box_size[d]) / ref_cell_bounding_box_size[d]);
9990df3356SJames Wright if (difference > PETSC_MACHINE_EPSILON * 100) {
10090df3356SJames 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);
10190df3356SJames Wright }
10290df3356SJames Wright }
10390df3356SJames Wright PetscCall(DMPlexVecRestoreClosure(dm_read, coord_section, coords, cell, &num_closure, &cell_coords));
10490df3356SJames Wright }
10590df3356SJames Wright PetscCall(DMDestroy(&dm_read));
10690df3356SJames Wright }
10790df3356SJames Wright
10890df3356SJames Wright PetscCall(PetscFinalize());
10990df3356SJames Wright return 0;
11090df3356SJames Wright }
11190df3356SJames Wright
11290df3356SJames Wright /*TEST
11390df3356SJames Wright build:
11490df3356SJames Wright requires: cgns
11590df3356SJames Wright test:
11690df3356SJames Wright # Coordinates of dm_create are linear, but the CGNS writer projects them to quadratic to match the solution.
11790df3356SJames Wright # The verification checks that all the cells of the resulting file are the same size
11890df3356SJames 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
1193886731fSPierre Jolivet output_file: output/empty.out
12090df3356SJames Wright TEST*/
121