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