xref: /petsc/src/dm/impls/plex/tutorials/ex7.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1 static char help[] = "Create a Plex sphere from quads and create a P1 section\n\n";
2 
3 #include <petscdmplex.h>
4 
5 static PetscErrorCode SetupSection(DM dm) {
6   PetscSection s;
7   PetscInt     vStart, vEnd, v;
8 
9   PetscFunctionBeginUser;
10   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
11   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &s));
12   PetscCall(PetscSectionSetNumFields(s, 1));
13   PetscCall(PetscSectionSetFieldComponents(s, 0, 1));
14   PetscCall(PetscSectionSetChart(s, vStart, vEnd));
15   for (v = vStart; v < vEnd; ++v) {
16     PetscCall(PetscSectionSetDof(s, v, 1));
17     PetscCall(PetscSectionSetFieldDof(s, v, 0, 1));
18   }
19   PetscCall(PetscSectionSetUp(s));
20   PetscCall(DMSetLocalSection(dm, s));
21   PetscCall(PetscSectionDestroy(&s));
22   PetscFunctionReturn(0);
23 }
24 
25 int main(int argc, char **argv) {
26   DM  dm;
27   Vec u;
28 
29   PetscFunctionBeginUser;
30   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
31   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
32   PetscCall(DMSetType(dm, DMPLEX));
33   PetscCall(DMSetFromOptions(dm));
34   PetscCall(PetscObjectSetName((PetscObject)dm, "Sphere"));
35   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
36 
37   PetscCall(SetupSection(dm));
38   PetscCall(DMGetGlobalVector(dm, &u));
39   PetscCall(VecSet(u, 2));
40   PetscCall(VecViewFromOptions(u, NULL, "-vec_view"));
41   PetscCall(DMRestoreGlobalVector(dm, &u));
42   PetscCall(DMDestroy(&dm));
43   PetscCall(PetscFinalize());
44   return 0;
45 }
46 
47 /*TEST
48 
49   testset:
50     requires: !__float128
51     args: -dm_plex_shape sphere -dm_view
52 
53     test:
54       suffix: 2d_quad
55       args: -dm_plex_simplex 0
56 
57     test:
58       suffix: 2d_tri
59       args:
60 
61     test:
62       suffix: 3d_tri
63       args: -dm_plex_dim 3
64 
65   testset:
66     requires: !__float128
67     args: -dm_plex_shape sphere -petscpartitioner_type simple -dm_view
68 
69     test:
70       suffix: 2d_quad_parallel
71       nsize: 2
72       args: -dm_plex_simplex 0
73 
74     test:
75       suffix: 2d_tri_parallel
76       nsize: 2
77 
78     test:
79       suffix: 3d_tri_parallel
80       nsize: 2
81       args: -dm_plex_dim 3
82 
83 TEST*/
84