Actual source code: ex7.c
1: static char help[] = "Create a Plex sphere from quads and create a P1 section\n\n";
3: #include <petscdmplex.h>
5: typedef struct {
6: PetscInt dim; /* Topological problem dimension */
7: PetscBool simplex; /* Mesh with simplices */
8: } AppCtx;
10: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
11: {
15: options->dim = 2;
16: options->simplex = PETSC_FALSE;
18: PetscOptionsBegin(comm, "", "Sphere Mesh Options", "DMPLEX");
19: PetscOptionsRangeInt("-dim", "Problem dimension", "ex7.c", options->dim, &options->dim, NULL,1,3);
20: PetscOptionsBool("-simplex", "Use simplices, or tensor product cells", "ex7.c", options->simplex, &options->simplex, NULL);
21: PetscOptionsEnd();
22: return(0);
23: }
25: static PetscErrorCode SetupSection(DM dm)
26: {
27: PetscSection s;
28: PetscInt vStart, vEnd, v;
32: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
33: PetscSectionCreate(PetscObjectComm((PetscObject) dm), &s);
34: PetscSectionSetNumFields(s, 1);
35: PetscSectionSetFieldComponents(s, 0, 1);
36: PetscSectionSetChart(s, vStart, vEnd);
37: for (v = vStart; v < vEnd; ++v) {
38: PetscSectionSetDof(s, v, 1);
39: PetscSectionSetFieldDof(s, v, 0, 1);
40: }
41: PetscSectionSetUp(s);
42: DMSetLocalSection(dm, s);
43: PetscSectionDestroy(&s);
44: return(0);
45: }
47: int main(int argc, char **argv)
48: {
49: DM dm;
50: Vec u;
51: AppCtx ctx;
54: PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
55: ProcessOptions(PETSC_COMM_WORLD, &ctx);
56: DMPlexCreateSphereMesh(PETSC_COMM_WORLD, ctx.dim, ctx.simplex, 1.0, &dm);
57: PetscObjectSetName((PetscObject) dm, "Sphere");
58: /* Distribute mesh over processes */
59: {
60: DM dmDist = NULL;
61: PetscPartitioner part;
62: DMPlexGetPartitioner(dm, &part);
63: PetscPartitionerSetFromOptions(part);
64: DMPlexDistribute(dm, 0, NULL, &dmDist);
65: if (dmDist) {
66: DMDestroy(&dm);
67: dm = dmDist;
68: }
69: }
70: DMSetFromOptions(dm);
71: DMViewFromOptions(dm, NULL, "-dm_view");
72: SetupSection(dm);
73: DMGetGlobalVector(dm, &u);
74: VecSet(u, 2);
75: VecViewFromOptions(u, NULL, "-vec_view");
76: DMRestoreGlobalVector(dm, &u);
77: DMDestroy(&dm);
78: PetscFinalize();
79: return ierr;
80: }
82: /*TEST
84: test:
85: suffix: 2d_quad
86: requires: !__float128
87: args: -dm_view
89: test:
90: suffix: 2d_quad_parallel
91: requires: !__float128
92: args: -dm_view -petscpartitioner_type simple
93: nsize: 2
95: test:
96: suffix: 2d_tri
97: requires: !__float128
98: args: -simplex -dm_view
100: test:
101: suffix: 2d_tri_parallel
102: requires: !__float128
103: args: -simplex -dm_view -petscpartitioner_type simple
104: nsize: 2
106: test:
107: suffix: 3d_tri
108: requires: !__float128
109: args: -dim 3 -simplex -dm_view
111: test:
112: suffix: 3d_tri_parallel
113: requires: !__float128
114: args: -dim 3 -simplex -dm_view -petscpartitioner_type simple
115: nsize: 2
117: TEST*/