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