xref: /petsc/src/dm/impls/plex/tests/ex61.c (revision f1580f4e3ce5d5b2393648fd039d0d41b440385d)
1 const char help[] = "Test boundary condition insertion";
2 
3 #include <petscdmplex.h>
4 
5 static PetscErrorCode set_one(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[], void *ctx) {
6   PetscFunctionBegin;
7   bcval[0] = 1.;
8   PetscFunctionReturn(0);
9 }
10 
11 static PetscErrorCode set_two(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[], void *ctx) {
12   PetscFunctionBegin;
13   bcval[0] = 2.;
14   PetscFunctionReturn(0);
15 }
16 
17 int main(int argc, char **argv) {
18   DM       dm;
19   DMLabel  label;
20   PetscInt in_value  = 1;
21   PetscInt out_value = 3;
22   PetscInt comps[]   = {0};
23   PetscFE  fe;
24   Vec      localVec;
25 
26   PetscFunctionBeginUser;
27   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
28   PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, 2, PETSC_FALSE, NULL, NULL, NULL, NULL, PETSC_TRUE, &dm));
29   PetscCall(DMGetLabel(dm, "Face Sets", &label));
30   PetscCall(PetscFECreateLagrange(PETSC_COMM_WORLD, 2, 1, PETSC_FALSE, 1, PETSC_DETERMINE, &fe));
31   PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
32   PetscCall(PetscFEDestroy(&fe));
33   PetscCall(DMCreateDS(dm));
34   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "inflow condition", label, 1, &in_value, 0, 1, comps, (void (*)(void))set_one, NULL, NULL, NULL));
35   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "outflow condition", label, 1, &out_value, 0, 1, comps, (void (*)(void))set_two, NULL, NULL, NULL));
36   PetscCall(DMCreateLocalVector(dm, &localVec));
37   PetscCall(VecSet(localVec, 0.));
38   PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localVec, 0.0, NULL, NULL, NULL));
39   PetscCall(VecView(localVec, NULL));
40   PetscCall(VecDestroy(&localVec));
41   PetscCall(DMDestroy(&dm));
42   PetscCall(PetscFinalize());
43   return 0;
44 }
45 
46 /*TEST
47 
48   test:
49     suffix: 0
50 
51 TEST*/
52