xref: /petsc/src/dm/impls/plex/tests/ex38.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
13c3287b0SJed Brown const char help[] = "Test DMPlexInsertBoundaryValues with DMPlexSetClosurePermutationTensor.\n";
23c3287b0SJed Brown 
33c3287b0SJed Brown #include <petscdmplex.h>
43c3287b0SJed Brown 
53c3287b0SJed Brown static PetscErrorCode bc_func(PetscInt dim, PetscReal time,
63c3287b0SJed Brown                      const PetscReal coords[], PetscInt num_comp_u,
73c3287b0SJed Brown                      PetscScalar *u, void *ctx)
83c3287b0SJed Brown {
93c3287b0SJed Brown   PetscFunctionBeginUser;
103c3287b0SJed Brown   for (PetscInt i=0; i<num_comp_u; i++) u[i] = coords[i];
113c3287b0SJed Brown   PetscFunctionReturn(0);
123c3287b0SJed Brown }
133c3287b0SJed Brown 
143c3287b0SJed Brown int main(int argc,char **argv)
153c3287b0SJed Brown {
163c3287b0SJed Brown   DM        dm;
173c3287b0SJed Brown   PetscFE   fe;
183c3287b0SJed Brown   Vec       U_loc;
193c3287b0SJed Brown   PetscInt  dim, order = 1;
203c3287b0SJed Brown   PetscBool tensorCoords = PETSC_TRUE;
213c3287b0SJed Brown 
223c3287b0SJed Brown   /* Initialize PETSc */
23*327415f7SBarry Smith   PetscFunctionBeginUser;
243c3287b0SJed Brown   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
253c3287b0SJed Brown   PetscCall(PetscOptionsGetBool(NULL, NULL, "-tensor_coords", &tensorCoords, NULL));
263c3287b0SJed Brown   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
273c3287b0SJed Brown   PetscCall(DMSetType(dm, DMPLEX));
283c3287b0SJed Brown   PetscCall(DMSetFromOptions(dm));
293c3287b0SJed Brown 
303c3287b0SJed Brown   PetscCall(DMGetDimension(dm, &dim));
313c3287b0SJed Brown   PetscCall(PetscFECreateLagrange(PETSC_COMM_WORLD, dim, dim, PETSC_FALSE, order, order, &fe));
323c3287b0SJed Brown   PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
333c3287b0SJed Brown   PetscCall(PetscFEDestroy(&fe));
343c3287b0SJed Brown   PetscCall(DMCreateDS(dm));
353c3287b0SJed Brown 
363c3287b0SJed Brown   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
373c3287b0SJed Brown 
383c3287b0SJed Brown   DMLabel label;
393c3287b0SJed Brown   PetscInt marker_ids[] = {1};
403c3287b0SJed Brown   PetscCall(DMGetLabel(dm, "marker", &label));
413c3287b0SJed Brown   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "mms", label, 1, marker_ids, 0, 0, NULL, (void(*)(void))bc_func, NULL, NULL, NULL));
423c3287b0SJed Brown   PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
433c3287b0SJed Brown   {
443c3287b0SJed Brown     DM cdm;
453c3287b0SJed Brown     PetscCall(DMGetCoordinateDM(dm, &cdm));
463c3287b0SJed Brown     if (tensorCoords) PetscCall(DMPlexSetClosurePermutationTensor(cdm, PETSC_DETERMINE, NULL));
473c3287b0SJed Brown   }
483c3287b0SJed Brown 
493c3287b0SJed Brown   PetscCall(DMCreateLocalVector(dm, &U_loc));
503c3287b0SJed Brown   PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, U_loc, 1., NULL, NULL, NULL));
513c3287b0SJed Brown   PetscCall(VecViewFromOptions(U_loc, NULL, "-u_loc_vec_view"));
523c3287b0SJed Brown   PetscCall(VecDestroy(&U_loc));
533c3287b0SJed Brown   PetscCall(DMDestroy(&dm));
543c3287b0SJed Brown   PetscCall(PetscFinalize());
553c3287b0SJed Brown   return 0;
563c3287b0SJed Brown }
573c3287b0SJed Brown 
583c3287b0SJed Brown /*TEST
593c3287b0SJed Brown   test:
603c3287b0SJed Brown     suffix: 2d
613c3287b0SJed Brown     args: -dm_plex_simplex 0 -dm_plex_dim 2 -dm_plex_box_faces 3,3 -u_loc_vec_view
623c3287b0SJed Brown   test:
633c3287b0SJed Brown     suffix: 3d
643c3287b0SJed Brown     args: -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -u_loc_vec_view
653c3287b0SJed Brown TEST*/
66