xref: /petsc/src/dm/impls/plex/tutorials/ex7.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
12*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionCreate(PetscObjectComm((PetscObject) dm), &s));
13*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionSetNumFields(s, 1));
14*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionSetFieldComponents(s, 0, 1));
15*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionSetChart(s, vStart, vEnd));
16c4762a1bSJed Brown   for (v = vStart; v < vEnd; ++v) {
17*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetDof(s, v, 1));
18*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetFieldDof(s, v, 0, 1));
19c4762a1bSJed Brown   }
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionSetUp(s));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetLocalSection(dm, s));
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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   PetscErrorCode ierr;
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(PETSC_COMM_WORLD, &dm));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(dm, DMPLEX));
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(dm));
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) dm, "Sphere"));
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view"));
3830602db0SMatthew G. Knepley 
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SetupSection(dm));
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(dm, &u));
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(u, 2));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(u, NULL, "-vec_view"));
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(dm, &u));
44*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
45c4762a1bSJed Brown   ierr = PetscFinalize();
46c4762a1bSJed Brown   return ierr;
47c4762a1bSJed Brown }
48c4762a1bSJed Brown 
49c4762a1bSJed Brown /*TEST
50c4762a1bSJed Brown 
5130602db0SMatthew G. Knepley   testset:
52c4762a1bSJed Brown     requires: !__float128
5330602db0SMatthew G. Knepley     args: -dm_plex_shape sphere -dm_view
54c4762a1bSJed Brown 
55c4762a1bSJed Brown     test:
5630602db0SMatthew G. Knepley       suffix: 2d_quad
5730602db0SMatthew G. Knepley       args: -dm_plex_simplex 0
58c4762a1bSJed Brown 
59c4762a1bSJed Brown     test:
60c4762a1bSJed Brown       suffix: 2d_tri
6130602db0SMatthew G. Knepley       args:
62c4762a1bSJed Brown 
63c4762a1bSJed Brown     test:
64c4762a1bSJed Brown       suffix: 3d_tri
6530602db0SMatthew G. Knepley       args: -dm_plex_dim 3
6630602db0SMatthew G. Knepley 
6730602db0SMatthew G. Knepley   testset:
68c4762a1bSJed Brown     requires: !__float128
69e600fa54SMatthew G. Knepley     args: -dm_plex_shape sphere -petscpartitioner_type simple -dm_view
7030602db0SMatthew G. Knepley 
7130602db0SMatthew G. Knepley     test:
7230602db0SMatthew G. Knepley       suffix: 2d_quad_parallel
7330602db0SMatthew G. Knepley       nsize: 2
7430602db0SMatthew G. Knepley       args: -dm_plex_simplex 0
7530602db0SMatthew G. Knepley 
7630602db0SMatthew G. Knepley     test:
7730602db0SMatthew G. Knepley       suffix: 2d_tri_parallel
7830602db0SMatthew G. Knepley       nsize: 2
79c4762a1bSJed Brown 
80c4762a1bSJed Brown     test:
81c4762a1bSJed Brown       suffix: 3d_tri_parallel
82c4762a1bSJed Brown       nsize: 2
8330602db0SMatthew G. Knepley       args: -dm_plex_dim 3
84c4762a1bSJed Brown 
85c4762a1bSJed Brown TEST*/
86