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