xref: /petsc/src/dm/impls/plex/tests/ex101.c (revision 92a26154bc948bc3af330cd082baf461ef1aa548)
153b735f8SJames Wright static char help[] = "Verify isoperiodic cone corrections";
253b735f8SJames Wright 
353b735f8SJames Wright #include <petscdmplex.h>
4*92a26154SJames Wright #include <petscsf.h>
5*92a26154SJames Wright #define EX "ex101.c"
653b735f8SJames Wright 
753b735f8SJames Wright // Creates periodic solution on a [0,1] x D domain for D dimension
853b735f8SJames Wright static PetscErrorCode project_function(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
953b735f8SJames Wright {
1053b735f8SJames Wright   PetscReal x_tot = 0;
1153b735f8SJames Wright 
1253b735f8SJames Wright   PetscFunctionBeginUser;
1353b735f8SJames Wright   for (PetscInt d = 0; d < dim; d++) x_tot += x[d];
1453b735f8SJames Wright   for (PetscInt c = 0; c < Nc; c++) {
15*92a26154SJames Wright     PetscScalar value = c % 2 ? PetscSinReal(2 * M_PI * x_tot) : PetscCosReal(2 * M_PI * x_tot);
1653b735f8SJames Wright     if (PetscAbsScalar(value) < 1e-7) value = 0.;
1753b735f8SJames Wright     u[c] = value;
1853b735f8SJames Wright   }
1953b735f8SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2053b735f8SJames Wright }
2153b735f8SJames Wright 
22*92a26154SJames Wright PetscErrorCode CreateFEField(DM dm, PetscBool use_natural_fe, PetscInt num_comps)
2353b735f8SJames Wright {
2453b735f8SJames Wright   PetscInt degree;
2553b735f8SJames Wright 
2653b735f8SJames Wright   PetscFunctionBeginUser;
2753b735f8SJames Wright   { // Get degree of the coords section
28*92a26154SJames Wright     PetscFE    fe;
29*92a26154SJames Wright     PetscSpace basis_space;
3053b735f8SJames Wright 
31*92a26154SJames Wright     if (use_natural_fe) {
32*92a26154SJames Wright       PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fe));
33*92a26154SJames Wright     } else {
34*92a26154SJames Wright       DM cdm;
3553b735f8SJames Wright       PetscCall(DMGetCoordinateDM(dm, &cdm));
36*92a26154SJames Wright       PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe));
37*92a26154SJames Wright     }
38*92a26154SJames Wright     PetscCall(PetscFEGetBasisSpace(fe, &basis_space));
39*92a26154SJames Wright     PetscCall(PetscSpaceGetDegree(basis_space, &degree, NULL));
4053b735f8SJames Wright   }
4153b735f8SJames Wright 
4253b735f8SJames Wright   PetscCall(DMClearFields(dm));
4353b735f8SJames Wright   PetscCall(DMSetLocalSection(dm, NULL)); // See https://gitlab.com/petsc/petsc/-/issues/1669
4453b735f8SJames Wright 
4553b735f8SJames Wright   { // Setup fe to load in the initial condition data
4653b735f8SJames Wright     PetscFE  fe;
4753b735f8SJames Wright     PetscInt dim;
4853b735f8SJames Wright 
4953b735f8SJames Wright     PetscCall(DMGetDimension(dm, &dim));
50*92a26154SJames Wright     PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, num_comps, PETSC_FALSE, degree, PETSC_DETERMINE, &fe));
5153b735f8SJames Wright     PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
5253b735f8SJames Wright     PetscCall(DMCreateDS(dm));
5353b735f8SJames Wright     PetscCall(PetscFEDestroy(&fe));
5453b735f8SJames Wright   }
5553b735f8SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
5653b735f8SJames Wright }
5753b735f8SJames Wright 
5853b735f8SJames Wright int main(int argc, char **argv)
5953b735f8SJames Wright {
6053b735f8SJames Wright   MPI_Comm  comm;
6153b735f8SJames Wright   DM        dm = NULL;
6253b735f8SJames Wright   Vec       V, V_G2L, V_local;
6353b735f8SJames Wright   PetscReal norm;
64*92a26154SJames Wright   PetscBool test_cgns_load = PETSC_FALSE;
65*92a26154SJames Wright   PetscInt  num_comps      = 1;
66*92a26154SJames Wright 
6753b735f8SJames Wright   PetscErrorCode (*funcs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) = {project_function};
6853b735f8SJames Wright 
6953b735f8SJames Wright   PetscFunctionBeginUser;
7053b735f8SJames Wright   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
7153b735f8SJames Wright   comm = PETSC_COMM_WORLD;
7253b735f8SJames Wright 
73*92a26154SJames Wright   PetscOptionsBegin(comm, "", "ex101.c Options", "DMPLEX");
74*92a26154SJames Wright   PetscCall(PetscOptionsBool("-test_cgns_load", "Test VecLoad using CGNS file", EX, test_cgns_load, &test_cgns_load, NULL));
75*92a26154SJames Wright   PetscCall(PetscOptionsInt("-num_comps", "Number of components in FE field", EX, num_comps, &num_comps, NULL));
76*92a26154SJames Wright   PetscOptionsEnd();
77*92a26154SJames Wright 
7853b735f8SJames Wright   PetscCall(DMCreate(comm, &dm));
7953b735f8SJames Wright   PetscCall(DMSetType(dm, DMPLEX));
8053b735f8SJames Wright   PetscCall(PetscObjectSetName((PetscObject)dm, "ex101_dm"));
8153b735f8SJames Wright   PetscCall(DMSetFromOptions(dm));
8253b735f8SJames Wright   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
83*92a26154SJames Wright   PetscCall(CreateFEField(dm, PETSC_FALSE, num_comps));
8453b735f8SJames Wright 
8553b735f8SJames Wright   // Verify that projected function on global vector (then projected onto local vector) is equal to projected function onto a local vector
8653b735f8SJames Wright   PetscCall(DMGetLocalVector(dm, &V_G2L));
8753b735f8SJames Wright   PetscCall(DMGetGlobalVector(dm, &V));
8853b735f8SJames Wright   PetscCall(DMProjectFunction(dm, 0, &funcs, NULL, INSERT_VALUES, V));
8953b735f8SJames Wright   PetscCall(DMGlobalToLocal(dm, V, INSERT_VALUES, V_G2L));
9053b735f8SJames Wright 
9153b735f8SJames Wright   PetscCall(DMGetLocalVector(dm, &V_local));
9253b735f8SJames Wright   PetscCall(DMProjectFunctionLocal(dm, 0, &funcs, NULL, INSERT_VALUES, V_local));
93*92a26154SJames Wright   PetscCall(VecViewFromOptions(V_local, NULL, "-local_view"));
9453b735f8SJames Wright 
95*92a26154SJames Wright   PetscCall(VecAXPY(V_G2L, -1, V_local));
96*92a26154SJames Wright   PetscCall(VecNorm(V_G2L, NORM_MAX, &norm));
9753b735f8SJames Wright   PetscReal tol = PetscDefined(USE_REAL___FLOAT128) ? 1e-12 : 1e4 * PETSC_MACHINE_EPSILON;
98*92a26154SJames Wright   if (norm > tol) PetscCall(PetscPrintf(comm, "Error! GlobalToLocal result does not match Local projection by norm %g\n", (double)norm));
99*92a26154SJames Wright 
100*92a26154SJames Wright   if (test_cgns_load) {
101*92a26154SJames Wright #ifndef PETSC_HAVE_CGNS
102*92a26154SJames Wright     SETERRQ(comm, PETSC_ERR_SUP, "PETSc not compiled with CGNS support");
103*92a26154SJames Wright #else
104*92a26154SJames Wright     PetscViewer viewer;
105*92a26154SJames Wright     DM          dm_read, dm_read_output;
106*92a26154SJames Wright     Vec         V_read, V_read_project2local, V_read_output2local;
107*92a26154SJames Wright     const char *filename = "test_file.cgns";
108*92a26154SJames Wright 
109*92a26154SJames Wright     // Write previous solution to filename
110*92a26154SJames Wright     PetscCall(PetscViewerCGNSOpen(comm, filename, FILE_MODE_WRITE, &viewer));
111*92a26154SJames Wright     PetscCall(VecView(V_local, viewer));
112*92a26154SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
113*92a26154SJames Wright 
114*92a26154SJames Wright     // Write DM from written CGNS file
115*92a26154SJames Wright     PetscCall(DMPlexCreateFromFile(comm, filename, "ex101_dm_read", PETSC_TRUE, &dm_read));
116*92a26154SJames Wright     PetscCall(DMSetFromOptions(dm_read));
117*92a26154SJames Wright     PetscCall(DMViewFromOptions(dm_read, NULL, "-dm_view"));
118*92a26154SJames Wright 
119*92a26154SJames Wright     { // Force isoperiodic point SF to be created to update sfNatural.
120*92a26154SJames Wright       // Needs to be done before removing the field corresponding to sfNatural
121*92a26154SJames Wright       PetscSection dummy_section;
122*92a26154SJames Wright       PetscCall(DMGetGlobalSection(dm_read, &dummy_section));
123*92a26154SJames Wright     }
124*92a26154SJames Wright     PetscCall(CreateFEField(dm_read, PETSC_TRUE, num_comps));
125*92a26154SJames Wright 
126*92a26154SJames Wright     // Use OutputDM as this doesn't use the isoperiodic pointSF and is compatible with sfNatural
127*92a26154SJames Wright     PetscCall(DMGetOutputDM(dm_read, &dm_read_output));
128*92a26154SJames Wright     PetscCall(DMGetLocalVector(dm_read, &V_read_project2local));
129*92a26154SJames Wright     PetscCall(DMGetLocalVector(dm_read, &V_read_output2local));
130*92a26154SJames Wright     PetscCall(PetscObjectSetName((PetscObject)V_read_output2local, "V_read_output2local"));
131*92a26154SJames Wright     PetscCall(PetscObjectSetName((PetscObject)V_read_project2local, "V_read_project2local"));
132*92a26154SJames Wright 
133*92a26154SJames Wright     PetscCall(DMProjectFunctionLocal(dm_read_output, 0, &funcs, NULL, INSERT_VALUES, V_read_project2local));
134*92a26154SJames Wright     PetscCall(VecViewFromOptions(V_read_project2local, NULL, "-project2local_view"));
135*92a26154SJames Wright 
136*92a26154SJames Wright     { // Read solution from file and communicate to local Vec
137*92a26154SJames Wright       PetscCall(DMGetGlobalVector(dm_read_output, &V_read));
138*92a26154SJames Wright       PetscCall(PetscObjectSetName((PetscObject)V_read, "V_read"));
139*92a26154SJames Wright 
140*92a26154SJames Wright       PetscCall(PetscViewerCGNSOpen(comm, filename, FILE_MODE_READ, &viewer));
141*92a26154SJames Wright       PetscCall(VecLoad(V_read, viewer));
142*92a26154SJames Wright       PetscCall(DMGlobalToLocal(dm_read_output, V_read, INSERT_VALUES, V_read_output2local));
143*92a26154SJames Wright 
144*92a26154SJames Wright       PetscCall(PetscViewerDestroy(&viewer));
145*92a26154SJames Wright       PetscCall(DMRestoreGlobalVector(dm_read_output, &V_read));
146*92a26154SJames Wright     }
147*92a26154SJames Wright     PetscCall(VecViewFromOptions(V_read_output2local, NULL, "-output2local_view"));
148*92a26154SJames Wright 
149*92a26154SJames Wright     PetscCall(VecAXPY(V_read_output2local, -1, V_read_project2local));
150*92a26154SJames Wright     PetscCall(VecNorm(V_read_output2local, NORM_MAX, &norm));
151*92a26154SJames Wright     if (norm > tol) PetscCall(PetscPrintf(comm, "Error! CGNS VecLoad result does not match Local projection by norm %g\n", (double)norm));
152*92a26154SJames Wright 
153*92a26154SJames Wright     PetscCall(DMRestoreLocalVector(dm_read, &V_read_project2local));
154*92a26154SJames Wright     PetscCall(DMRestoreLocalVector(dm_read, &V_read_output2local));
155*92a26154SJames Wright     PetscCall(DMDestroy(&dm_read));
156*92a26154SJames Wright #endif
157*92a26154SJames Wright   }
15853b735f8SJames Wright 
15953b735f8SJames Wright   PetscCall(DMRestoreGlobalVector(dm, &V));
16053b735f8SJames Wright   PetscCall(DMRestoreLocalVector(dm, &V_G2L));
16153b735f8SJames Wright   PetscCall(DMRestoreLocalVector(dm, &V_local));
16253b735f8SJames Wright   PetscCall(DMDestroy(&dm));
16353b735f8SJames Wright   PetscCall(PetscFinalize());
16453b735f8SJames Wright   return 0;
16553b735f8SJames Wright }
16653b735f8SJames Wright 
16753b735f8SJames Wright /*TEST
16853b735f8SJames Wright 
16953b735f8SJames Wright   test:
17053b735f8SJames Wright     nsize: 2
17153b735f8SJames Wright     args: -dm_plex_shape zbox -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 1,2,3 -dm_coord_space -dm_coord_petscspace_degree 3
17253b735f8SJames Wright     args: -dm_plex_box_bd periodic,periodic,periodic -dm_view ::ascii_info_detail -petscpartitioner_type simple
17353b735f8SJames Wright 
17453b735f8SJames Wright   test:
17553b735f8SJames Wright     requires: cgns
17653b735f8SJames Wright     suffix: cgns
17753b735f8SJames Wright     nsize: 3
178*92a26154SJames Wright     args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/2x2x2_Q3_wave.cgns -dm_plex_cgns_parallel -dm_view ::ascii_info_detail -dm_plex_box_label true -dm_plex_box_label_bd periodic,periodic,periodic -petscpartitioner_type simple -test_cgns_load -num_comps 2
179*92a26154SJames Wright 
180*92a26154SJames Wright   test:
181*92a26154SJames Wright     requires: cgns parmetis
182*92a26154SJames Wright     suffix: cgns_parmetis
183*92a26154SJames Wright     nsize: 3
184*92a26154SJames Wright     args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/2x2x2_Q3_wave.cgns -dm_plex_cgns_parallel -dm_plex_box_label true -dm_plex_box_label_bd periodic,periodic,periodic -petscpartitioner_type parmetis -test_cgns_load -num_comps 2
18553b735f8SJames Wright 
18653b735f8SJames Wright TEST*/
187