xref: /petsc/src/dm/impls/plex/tutorials/ex1.c (revision a69119a591a03a9d906b29c0a4e9802e4d7c9795)
1 static char help[] = "Define a simple field over the mesh\n\n";
2 
3 #include <petscdmplex.h>
4 
5 int main(int argc, char **argv) {
6   DM           dm;
7   Vec          u;
8   PetscSection section;
9   PetscViewer  viewer;
10   PetscInt     dim, numFields, numBC, i;
11   PetscInt     numComp[3];
12   PetscInt     numDof[12];
13   PetscInt     bcField[1];
14   IS           bcPointIS[1];
15 
16   PetscFunctionBeginUser;
17   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
18   /* Create a mesh */
19   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
20   PetscCall(DMSetType(dm, DMPLEX));
21   PetscCall(DMSetFromOptions(dm));
22   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
23   PetscCall(DMGetDimension(dm, &dim));
24   /* Create a scalar field u, a vector field v, and a surface vector field w */
25   numFields  = 3;
26   numComp[0] = 1;
27   numComp[1] = dim;
28   numComp[2] = dim - 1;
29   for (i = 0; i < numFields * (dim + 1); ++i) numDof[i] = 0;
30   /* Let u be defined on vertices */
31   numDof[0 * (dim + 1) + 0]       = 1;
32   /* Let v be defined on cells */
33   numDof[1 * (dim + 1) + dim]     = dim;
34   /* Let w be defined on faces */
35   numDof[2 * (dim + 1) + dim - 1] = dim - 1;
36   /* Setup boundary conditions */
37   numBC                           = 1;
38   /* Prescribe a Dirichlet condition on u on the boundary
39        Label "marker" is made by the mesh creation routine */
40   bcField[0]                      = 0;
41   PetscCall(DMGetStratumIS(dm, "marker", 1, &bcPointIS[0]));
42   /* Create a PetscSection with this data layout */
43   PetscCall(DMSetNumFields(dm, numFields));
44   PetscCall(DMPlexCreateSection(dm, NULL, numComp, numDof, numBC, bcField, NULL, bcPointIS, NULL, &section));
45   PetscCall(ISDestroy(&bcPointIS[0]));
46   /* Name the Field variables */
47   PetscCall(PetscSectionSetFieldName(section, 0, "u"));
48   PetscCall(PetscSectionSetFieldName(section, 1, "v"));
49   PetscCall(PetscSectionSetFieldName(section, 2, "w"));
50   PetscCall(PetscSectionView(section, PETSC_VIEWER_STDOUT_WORLD));
51   /* Tell the DM to use this data layout */
52   PetscCall(DMSetLocalSection(dm, section));
53   /* Create a Vec with this layout and view it */
54   PetscCall(DMGetGlobalVector(dm, &u));
55   PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer));
56   PetscCall(PetscViewerSetType(viewer, PETSCVIEWERVTK));
57   PetscCall(PetscViewerFileSetName(viewer, "sol.vtu"));
58   PetscCall(VecView(u, viewer));
59   PetscCall(PetscViewerDestroy(&viewer));
60   PetscCall(DMRestoreGlobalVector(dm, &u));
61   /* Cleanup */
62   PetscCall(PetscSectionDestroy(&section));
63   PetscCall(DMDestroy(&dm));
64   PetscCall(PetscFinalize());
65   return 0;
66 }
67 
68 /*TEST
69 
70   test:
71     suffix: 0
72     requires: triangle
73     args: -info :~sys,mat
74   test:
75     suffix: 1
76     requires: ctetgen
77     args: -dm_plex_dim 3 -info :~sys,mat
78 
79 TEST*/
80