xref: /petsc/src/dm/impls/plex/tutorials/ex7.c (revision e600fa544e2bb197ca2af9b6e65ea465976dec56)
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   PetscErrorCode ierr;
10c4762a1bSJed Brown 
11c4762a1bSJed Brown   PetscFunctionBeginUser;
12c4762a1bSJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
13c4762a1bSJed Brown   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &s);CHKERRQ(ierr);
14c4762a1bSJed Brown   ierr = PetscSectionSetNumFields(s, 1);CHKERRQ(ierr);
15c4762a1bSJed Brown   ierr = PetscSectionSetFieldComponents(s, 0, 1);CHKERRQ(ierr);
16c4762a1bSJed Brown   ierr = PetscSectionSetChart(s, vStart, vEnd);CHKERRQ(ierr);
17c4762a1bSJed Brown   for (v = vStart; v < vEnd; ++v) {
18c4762a1bSJed Brown     ierr = PetscSectionSetDof(s, v, 1);CHKERRQ(ierr);
19c4762a1bSJed Brown     ierr = PetscSectionSetFieldDof(s, v, 0, 1);CHKERRQ(ierr);
20c4762a1bSJed Brown   }
21c4762a1bSJed Brown   ierr = PetscSectionSetUp(s);CHKERRQ(ierr);
22c4762a1bSJed Brown   ierr = DMSetLocalSection(dm, s);CHKERRQ(ierr);
23c4762a1bSJed Brown   ierr = PetscSectionDestroy(&s);CHKERRQ(ierr);
24c4762a1bSJed Brown   PetscFunctionReturn(0);
25c4762a1bSJed Brown }
26c4762a1bSJed Brown 
27c4762a1bSJed Brown int main(int argc, char **argv)
28c4762a1bSJed Brown {
29c4762a1bSJed Brown   DM             dm;
30c4762a1bSJed Brown   Vec            u;
31c4762a1bSJed Brown   PetscErrorCode ierr;
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
3430602db0SMatthew G. Knepley   ierr = DMCreate(PETSC_COMM_WORLD, &dm);CHKERRQ(ierr);
3530602db0SMatthew G. Knepley   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
36c4762a1bSJed Brown   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
3730602db0SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) dm, "Sphere");CHKERRQ(ierr);
38c4762a1bSJed Brown   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
3930602db0SMatthew G. Knepley 
40c4762a1bSJed Brown   ierr = SetupSection(dm);CHKERRQ(ierr);
41c4762a1bSJed Brown   ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr);
42c4762a1bSJed Brown   ierr = VecSet(u, 2);CHKERRQ(ierr);
43c4762a1bSJed Brown   ierr = VecViewFromOptions(u, NULL, "-vec_view");CHKERRQ(ierr);
44c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr);
45c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
46c4762a1bSJed Brown   ierr = PetscFinalize();
47c4762a1bSJed Brown   return ierr;
48c4762a1bSJed Brown }
49c4762a1bSJed Brown 
50c4762a1bSJed Brown /*TEST
51c4762a1bSJed Brown 
5230602db0SMatthew G. Knepley   testset:
53c4762a1bSJed Brown     requires: !__float128
5430602db0SMatthew G. Knepley     args: -dm_plex_shape sphere -dm_view
55c4762a1bSJed Brown 
56c4762a1bSJed Brown     test:
5730602db0SMatthew G. Knepley       suffix: 2d_quad
5830602db0SMatthew G. Knepley       args: -dm_plex_simplex 0
59c4762a1bSJed Brown 
60c4762a1bSJed Brown     test:
61c4762a1bSJed Brown       suffix: 2d_tri
6230602db0SMatthew G. Knepley       args:
63c4762a1bSJed Brown 
64c4762a1bSJed Brown     test:
65c4762a1bSJed Brown       suffix: 3d_tri
6630602db0SMatthew G. Knepley       args: -dm_plex_dim 3
6730602db0SMatthew G. Knepley 
6830602db0SMatthew G. Knepley   testset:
69c4762a1bSJed Brown     requires: !__float128
70*e600fa54SMatthew G. Knepley     args: -dm_plex_shape sphere -petscpartitioner_type simple -dm_view
7130602db0SMatthew G. Knepley 
7230602db0SMatthew G. Knepley     test:
7330602db0SMatthew G. Knepley       suffix: 2d_quad_parallel
7430602db0SMatthew G. Knepley       nsize: 2
7530602db0SMatthew G. Knepley       args: -dm_plex_simplex 0
7630602db0SMatthew G. Knepley 
7730602db0SMatthew G. Knepley     test:
7830602db0SMatthew G. Knepley       suffix: 2d_tri_parallel
7930602db0SMatthew G. Knepley       nsize: 2
80c4762a1bSJed Brown 
81c4762a1bSJed Brown     test:
82c4762a1bSJed Brown       suffix: 3d_tri_parallel
83c4762a1bSJed Brown       nsize: 2
8430602db0SMatthew G. Knepley       args: -dm_plex_dim 3
85c4762a1bSJed Brown 
86c4762a1bSJed Brown TEST*/
87