xref: /petsc/src/dm/impls/plex/tests/ex48.c (revision d5b43468fb8780a8feea140ccd6fa3e6a50411cc)
1 static char help[] = "Tests for VecGetValuesSection / VecSetValuesSection \n\n";
2 
3 #include <petscdmplex.h>
4 
5 int main(int argc, char **argv)
6 {
7   DM           dm;
8   Vec          v;
9   PetscSection section;
10   PetscScalar  val[2];
11   PetscInt     pStart, pEnd, p;
12 
13   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
14 
15   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
16   PetscCall(DMSetType(dm, DMPLEX));
17   PetscCall(DMSetFromOptions(dm));
18   PetscCall(DMViewFromOptions(dm, NULL, "-d_view"));
19 
20   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section));
21   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
22   PetscCall(PetscSectionSetChart(section, pStart, pEnd));
23   PetscCall(DMPlexGetHeightStratum(dm, 0, &pStart, &pEnd));
24   for (p = pStart; p < pEnd; p++) PetscCall(PetscSectionSetDof(section, p, 1));
25   PetscCall(DMPlexGetDepthStratum(dm, 0, &pStart, &pEnd));
26   for (p = pStart; p < pEnd; p++) PetscCall(PetscSectionSetDof(section, p, 2));
27   PetscCall(PetscSectionSetUp(section));
28   PetscCall(DMSetLocalSection(dm, section));
29   PetscCall(PetscSectionViewFromOptions(section, NULL, "-s_view"));
30 
31   PetscCall(DMCreateGlobalVector(dm, &v));
32   PetscCall(VecViewFromOptions(v, NULL, "-v_view"));
33 
34   /* look through all cells and change "cell values" */
35   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
36   for (p = pStart; p < pEnd; ++p) {
37     PetscInt dof;
38 
39     PetscCall(PetscSectionGetDof(section, p, &dof));
40     for (PetscInt d = 0; d < dof; ++d) val[d] = 100 * p + d;
41     PetscCall(VecSetValuesSection(v, section, p, val, INSERT_VALUES));
42   }
43   PetscCall(VecView(v, PETSC_VIEWER_STDOUT_WORLD));
44 
45   for (p = pStart; p < pEnd; ++p) {
46     PetscScalar *x;
47     PetscInt     dof;
48 
49     PetscCall(PetscSectionGetDof(section, p, &dof));
50     PetscCall(VecGetValuesSection(v, section, p, &x));
51     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Point #%" PetscInt_FMT " %" PetscInt_FMT " dof\n", p, dof));
52   }
53 
54   PetscCall(VecDestroy(&v));
55   PetscCall(PetscSectionDestroy(&section));
56   PetscCall(DMDestroy(&dm));
57   PetscCall(PetscFinalize());
58   return 0;
59 }
60 
61 /*TEST
62 
63   test:
64     suffix: 0
65     args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/quads-q2.msh
66 
67 TEST*/
68