Actual source code: ex1.c
petsc-dev 2014-02-02
1: static char help[] = "Define a simple field over the mesh\n\n";
3: #include <petscdmplex.h>
7: int main(int argc, char **argv)
8: {
9: DM dm;
10: Vec u;
11: PetscSection section;
12: PetscViewer viewer;
13: PetscInt dim = 2, numFields, numBC, i;
14: PetscInt numComp[3];
15: PetscInt numDof[12];
16: PetscInt bcField[1];
17: IS bcPointIS[1];
18: PetscBool interpolate = PETSC_TRUE;
21: PetscInitialize(&argc, &argv, NULL, help);
22: PetscOptionsGetInt(NULL, "-dim", &dim, NULL);
23: /* Create a mesh */
24: DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, interpolate, &dm);
25: /* Create a scalar field u, a vector field v, and a surface vector field w */
26: numFields = 3;
27: numComp[0] = 1;
28: numComp[1] = dim;
29: numComp[2] = dim-1;
30: for (i = 0; i < numFields*(dim+1); ++i) numDof[i] = 0;
31: /* Let u be defined on vertices */
32: numDof[0*(dim+1)+0] = 1;
33: /* Let v be defined on cells */
34: numDof[1*(dim+1)+dim] = dim;
35: /* Let v be defined on faces */
36: numDof[2*(dim+1)+dim-1] = dim-1;
37: /* Setup boundary conditions */
38: numBC = 1;
39: /* Prescribe a Dirichlet condition on u on the boundary
40: Label "marker" is made by the mesh creation routine */
41: bcField[0] = 0;
42: DMPlexGetStratumIS(dm, "marker", 1, &bcPointIS[0]);
43: /* Create a PetscSection with this data layout */
44: DMPlexCreateSection(dm, dim, numFields, numComp, numDof, numBC, bcField, bcPointIS, §ion);
45: ISDestroy(&bcPointIS[0]);
46: /* Name the Field variables */
47: PetscSectionSetFieldName(section, 0, "u");
48: PetscSectionSetFieldName(section, 1, "v");
49: PetscSectionSetFieldName(section, 2, "w");
50: PetscSectionView(section, PETSC_VIEWER_STDOUT_WORLD);
51: /* Tell the DM to use this data layout */
52: DMSetDefaultSection(dm, section);
53: /* Create a Vec with this layout and view it */
54: DMGetGlobalVector(dm, &u);
55: PetscViewerCreate(PETSC_COMM_WORLD, &viewer);
56: PetscViewerSetType(viewer, PETSCVIEWERVTK);
57: PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_VTK);
58: PetscViewerFileSetName(viewer, "sol.vtk");
59: VecView(u, viewer);
60: PetscViewerDestroy(&viewer);
61: DMRestoreGlobalVector(dm, &u);
62: /* Cleanup */
63: PetscSectionDestroy(§ion);
64: DMDestroy(&dm);
65: PetscFinalize();
66: return 0;
67: }