xref: /petsc/src/dm/impls/plex/tests/ex38.c (revision a69119a591a03a9d906b29c0a4e9802e4d7c9795)
1 const char help[] = "Test DMPlexInsertBoundaryValues with DMPlexSetClosurePermutationTensor.\n";
2 
3 #include <petscdmplex.h>
4 
5 static PetscErrorCode bc_func(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt num_comp_u, PetscScalar *u, void *ctx) {
6   PetscFunctionBeginUser;
7   for (PetscInt i = 0; i < num_comp_u; i++) u[i] = coords[i];
8   PetscFunctionReturn(0);
9 }
10 
11 int main(int argc, char **argv) {
12   DM        dm;
13   PetscFE   fe;
14   Vec       U_loc;
15   PetscInt  dim, order = 1;
16   PetscBool tensorCoords = PETSC_TRUE;
17 
18   /* Initialize PETSc */
19   PetscFunctionBeginUser;
20   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
21   PetscCall(PetscOptionsGetBool(NULL, NULL, "-tensor_coords", &tensorCoords, NULL));
22   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
23   PetscCall(DMSetType(dm, DMPLEX));
24   PetscCall(DMSetFromOptions(dm));
25 
26   PetscCall(DMGetDimension(dm, &dim));
27   PetscCall(PetscFECreateLagrange(PETSC_COMM_WORLD, dim, dim, PETSC_FALSE, order, order, &fe));
28   PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
29   PetscCall(PetscFEDestroy(&fe));
30   PetscCall(DMCreateDS(dm));
31 
32   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
33 
34   DMLabel  label;
35   PetscInt marker_ids[] = {1};
36   PetscCall(DMGetLabel(dm, "marker", &label));
37   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "mms", label, 1, marker_ids, 0, 0, NULL, (void (*)(void))bc_func, NULL, NULL, NULL));
38   PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
39   {
40     DM cdm;
41     PetscCall(DMGetCoordinateDM(dm, &cdm));
42     if (tensorCoords) PetscCall(DMPlexSetClosurePermutationTensor(cdm, PETSC_DETERMINE, NULL));
43   }
44 
45   PetscCall(DMCreateLocalVector(dm, &U_loc));
46   PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, U_loc, 1., NULL, NULL, NULL));
47   PetscCall(VecViewFromOptions(U_loc, NULL, "-u_loc_vec_view"));
48   PetscCall(VecDestroy(&U_loc));
49   PetscCall(DMDestroy(&dm));
50   PetscCall(PetscFinalize());
51   return 0;
52 }
53 
54 /*TEST
55   test:
56     suffix: 2d
57     args: -dm_plex_simplex 0 -dm_plex_dim 2 -dm_plex_box_faces 3,3 -u_loc_vec_view
58   test:
59     suffix: 3d
60     args: -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -u_loc_vec_view
61 TEST*/
62