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