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