153b735f8SJames Wright static char help[] = "Verify isoperiodic cone corrections"; 253b735f8SJames Wright 353b735f8SJames Wright #include <petscdmplex.h> 492a26154SJames Wright #include <petscsf.h> 592a26154SJames 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++) { 1592a26154SJames 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 2292a26154SJames 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 2892a26154SJames Wright PetscFE fe; 2992a26154SJames Wright PetscSpace basis_space; 3053b735f8SJames Wright 3192a26154SJames Wright if (use_natural_fe) { 3292a26154SJames Wright PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fe)); 3392a26154SJames Wright } else { 3492a26154SJames Wright DM cdm; 3553b735f8SJames Wright PetscCall(DMGetCoordinateDM(dm, &cdm)); 3692a26154SJames Wright PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe)); 3792a26154SJames Wright } 3892a26154SJames Wright PetscCall(PetscFEGetBasisSpace(fe, &basis_space)); 3992a26154SJames Wright PetscCall(PetscSpaceGetDegree(basis_space, °ree, 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)); 5092a26154SJames 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; 6492a26154SJames Wright PetscBool test_cgns_load = PETSC_FALSE; 6592a26154SJames Wright PetscInt num_comps = 1; 6692a26154SJames 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 7392a26154SJames Wright PetscOptionsBegin(comm, "", "ex101.c Options", "DMPLEX"); 7492a26154SJames Wright PetscCall(PetscOptionsBool("-test_cgns_load", "Test VecLoad using CGNS file", EX, test_cgns_load, &test_cgns_load, NULL)); 7592a26154SJames Wright PetscCall(PetscOptionsInt("-num_comps", "Number of components in FE field", EX, num_comps, &num_comps, NULL)); 7692a26154SJames Wright PetscOptionsEnd(); 7792a26154SJames 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")); 8392a26154SJames 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)); 9392a26154SJames Wright PetscCall(VecViewFromOptions(V_local, NULL, "-local_view")); 9453b735f8SJames Wright 9592a26154SJames Wright PetscCall(VecAXPY(V_G2L, -1, V_local)); 9692a26154SJames Wright PetscCall(VecNorm(V_G2L, NORM_MAX, &norm)); 9753b735f8SJames Wright PetscReal tol = PetscDefined(USE_REAL___FLOAT128) ? 1e-12 : 1e4 * PETSC_MACHINE_EPSILON; 9892a26154SJames Wright if (norm > tol) PetscCall(PetscPrintf(comm, "Error! GlobalToLocal result does not match Local projection by norm %g\n", (double)norm)); 9992a26154SJames Wright 10092a26154SJames Wright if (test_cgns_load) { 10192a26154SJames Wright #ifndef PETSC_HAVE_CGNS 10292a26154SJames Wright SETERRQ(comm, PETSC_ERR_SUP, "PETSc not compiled with CGNS support"); 10392a26154SJames Wright #else 10492a26154SJames Wright PetscViewer viewer; 10592a26154SJames Wright DM dm_read, dm_read_output; 10692a26154SJames Wright Vec V_read, V_read_project2local, V_read_output2local; 10792a26154SJames Wright const char *filename = "test_file.cgns"; 10892a26154SJames Wright 10992a26154SJames Wright // Write previous solution to filename 11092a26154SJames Wright PetscCall(PetscViewerCGNSOpen(comm, filename, FILE_MODE_WRITE, &viewer)); 11192a26154SJames Wright PetscCall(VecView(V_local, viewer)); 11292a26154SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 11392a26154SJames Wright 11492a26154SJames Wright // Write DM from written CGNS file 11592a26154SJames Wright PetscCall(DMPlexCreateFromFile(comm, filename, "ex101_dm_read", PETSC_TRUE, &dm_read)); 11692a26154SJames Wright PetscCall(DMSetFromOptions(dm_read)); 11792a26154SJames Wright PetscCall(DMViewFromOptions(dm_read, NULL, "-dm_view")); 11892a26154SJames Wright 11992a26154SJames Wright { // Force isoperiodic point SF to be created to update sfNatural. 12092a26154SJames Wright // Needs to be done before removing the field corresponding to sfNatural 12192a26154SJames Wright PetscSection dummy_section; 12292a26154SJames Wright PetscCall(DMGetGlobalSection(dm_read, &dummy_section)); 12392a26154SJames Wright } 12492a26154SJames Wright PetscCall(CreateFEField(dm_read, PETSC_TRUE, num_comps)); 12592a26154SJames Wright 12692a26154SJames Wright // Use OutputDM as this doesn't use the isoperiodic pointSF and is compatible with sfNatural 12792a26154SJames Wright PetscCall(DMGetOutputDM(dm_read, &dm_read_output)); 12892a26154SJames Wright PetscCall(DMGetLocalVector(dm_read, &V_read_project2local)); 12992a26154SJames Wright PetscCall(DMGetLocalVector(dm_read, &V_read_output2local)); 13092a26154SJames Wright PetscCall(PetscObjectSetName((PetscObject)V_read_output2local, "V_read_output2local")); 13192a26154SJames Wright PetscCall(PetscObjectSetName((PetscObject)V_read_project2local, "V_read_project2local")); 13292a26154SJames Wright 13392a26154SJames Wright PetscCall(DMProjectFunctionLocal(dm_read_output, 0, &funcs, NULL, INSERT_VALUES, V_read_project2local)); 13492a26154SJames Wright PetscCall(VecViewFromOptions(V_read_project2local, NULL, "-project2local_view")); 13592a26154SJames Wright 13692a26154SJames Wright { // Read solution from file and communicate to local Vec 13792a26154SJames Wright PetscCall(DMGetGlobalVector(dm_read_output, &V_read)); 13892a26154SJames Wright PetscCall(PetscObjectSetName((PetscObject)V_read, "V_read")); 13992a26154SJames Wright 14092a26154SJames Wright PetscCall(PetscViewerCGNSOpen(comm, filename, FILE_MODE_READ, &viewer)); 14192a26154SJames Wright PetscCall(VecLoad(V_read, viewer)); 14292a26154SJames Wright PetscCall(DMGlobalToLocal(dm_read_output, V_read, INSERT_VALUES, V_read_output2local)); 14392a26154SJames Wright 14492a26154SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 14592a26154SJames Wright PetscCall(DMRestoreGlobalVector(dm_read_output, &V_read)); 14692a26154SJames Wright } 14792a26154SJames Wright PetscCall(VecViewFromOptions(V_read_output2local, NULL, "-output2local_view")); 14892a26154SJames Wright 14992a26154SJames Wright PetscCall(VecAXPY(V_read_output2local, -1, V_read_project2local)); 15092a26154SJames Wright PetscCall(VecNorm(V_read_output2local, NORM_MAX, &norm)); 15192a26154SJames Wright if (norm > tol) PetscCall(PetscPrintf(comm, "Error! CGNS VecLoad result does not match Local projection by norm %g\n", (double)norm)); 15292a26154SJames Wright 15392a26154SJames Wright PetscCall(DMRestoreLocalVector(dm_read, &V_read_project2local)); 15492a26154SJames Wright PetscCall(DMRestoreLocalVector(dm_read, &V_read_output2local)); 15592a26154SJames Wright PetscCall(DMDestroy(&dm_read)); 15692a26154SJames Wright #endif 15792a26154SJames 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 17892a26154SJames 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 17992a26154SJames Wright 18092a26154SJames Wright test: 18192a26154SJames Wright requires: cgns parmetis 18292a26154SJames Wright suffix: cgns_parmetis 18392a26154SJames Wright nsize: 3 18492a26154SJames 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 185*e0008caeSPierre Jolivet output_file: output/empty.out 18653b735f8SJames Wright 18753b735f8SJames Wright TEST*/ 188