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