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