xref: /petsc/src/dm/impls/plex/tutorials/ex7.c (revision 9566063d113dddea24716c546802770db7481bc0)
1c4762a1bSJed Brown static char help[] = "Create a Plex sphere from quads and create a P1 section\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown static PetscErrorCode SetupSection(DM dm)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   PetscSection   s;
8c4762a1bSJed Brown   PetscInt       vStart, vEnd, v;
9c4762a1bSJed Brown 
10c4762a1bSJed Brown   PetscFunctionBeginUser;
11*9566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
12*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) dm), &s));
13*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetNumFields(s, 1));
14*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldComponents(s, 0, 1));
15*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(s, vStart, vEnd));
16c4762a1bSJed Brown   for (v = vStart; v < vEnd; ++v) {
17*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(s, v, 1));
18*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldDof(s, v, 0, 1));
19c4762a1bSJed Brown   }
20*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(s));
21*9566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dm, s));
22*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&s));
23c4762a1bSJed Brown   PetscFunctionReturn(0);
24c4762a1bSJed Brown }
25c4762a1bSJed Brown 
26c4762a1bSJed Brown int main(int argc, char **argv)
27c4762a1bSJed Brown {
28c4762a1bSJed Brown   DM             dm;
29c4762a1bSJed Brown   Vec            u;
30c4762a1bSJed Brown 
31*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
32*9566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
33*9566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
34*9566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
35*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) dm, "Sphere"));
36*9566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
3730602db0SMatthew G. Knepley 
38*9566063dSJacob Faibussowitsch   PetscCall(SetupSection(dm));
39*9566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &u));
40*9566063dSJacob Faibussowitsch   PetscCall(VecSet(u, 2));
41*9566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-vec_view"));
42*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &u));
43*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
44*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
45b122ec5aSJacob Faibussowitsch   return 0;
46c4762a1bSJed Brown }
47c4762a1bSJed Brown 
48c4762a1bSJed Brown /*TEST
49c4762a1bSJed Brown 
5030602db0SMatthew G. Knepley   testset:
51c4762a1bSJed Brown     requires: !__float128
5230602db0SMatthew G. Knepley     args: -dm_plex_shape sphere -dm_view
53c4762a1bSJed Brown 
54c4762a1bSJed Brown     test:
5530602db0SMatthew G. Knepley       suffix: 2d_quad
5630602db0SMatthew G. Knepley       args: -dm_plex_simplex 0
57c4762a1bSJed Brown 
58c4762a1bSJed Brown     test:
59c4762a1bSJed Brown       suffix: 2d_tri
6030602db0SMatthew G. Knepley       args:
61c4762a1bSJed Brown 
62c4762a1bSJed Brown     test:
63c4762a1bSJed Brown       suffix: 3d_tri
6430602db0SMatthew G. Knepley       args: -dm_plex_dim 3
6530602db0SMatthew G. Knepley 
6630602db0SMatthew G. Knepley   testset:
67c4762a1bSJed Brown     requires: !__float128
68e600fa54SMatthew G. Knepley     args: -dm_plex_shape sphere -petscpartitioner_type simple -dm_view
6930602db0SMatthew G. Knepley 
7030602db0SMatthew G. Knepley     test:
7130602db0SMatthew G. Knepley       suffix: 2d_quad_parallel
7230602db0SMatthew G. Knepley       nsize: 2
7330602db0SMatthew G. Knepley       args: -dm_plex_simplex 0
7430602db0SMatthew G. Knepley 
7530602db0SMatthew G. Knepley     test:
7630602db0SMatthew G. Knepley       suffix: 2d_tri_parallel
7730602db0SMatthew G. Knepley       nsize: 2
78c4762a1bSJed Brown 
79c4762a1bSJed Brown     test:
80c4762a1bSJed Brown       suffix: 3d_tri_parallel
81c4762a1bSJed Brown       nsize: 2
8230602db0SMatthew G. Knepley       args: -dm_plex_dim 3
83c4762a1bSJed Brown 
84c4762a1bSJed Brown TEST*/
85