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