xref: /petsc/src/dm/impls/forest/tests/ex3.c (revision 9b5c7ca5e0ec6bfcaa81b66c7955351a760cc116)
1*9b5c7ca5SMatthew G. Knepley static char help[] = "Tests adaptive refinement using DMForest, and uses HDF5.\n\n";
2*9b5c7ca5SMatthew G. Knepley 
3*9b5c7ca5SMatthew G. Knepley #include <petscdmforest.h>
4*9b5c7ca5SMatthew G. Knepley #include <petscdmplex.h>
5*9b5c7ca5SMatthew G. Knepley #include <petscviewerhdf5.h>
6*9b5c7ca5SMatthew G. Knepley 
7*9b5c7ca5SMatthew G. Knepley int main (int argc, char **argv)
8*9b5c7ca5SMatthew G. Knepley {
9*9b5c7ca5SMatthew G. Knepley   DM             base, forest, plex;
10*9b5c7ca5SMatthew G. Knepley   PetscSection   s;
11*9b5c7ca5SMatthew G. Knepley   PetscViewer    viewer;
12*9b5c7ca5SMatthew G. Knepley   Vec            g = NULL, g2 = NULL;
13*9b5c7ca5SMatthew G. Knepley   PetscReal      nrm;
14*9b5c7ca5SMatthew G. Knepley   PetscBool      adapt = PETSC_FALSE, userSection = PETSC_FALSE;
15*9b5c7ca5SMatthew G. Knepley   PetscInt       vStart, vEnd, v, i;
16*9b5c7ca5SMatthew G. Knepley   PetscErrorCode ierr;
17*9b5c7ca5SMatthew G. Knepley 
18*9b5c7ca5SMatthew G. Knepley   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
19*9b5c7ca5SMatthew G. Knepley   ierr = PetscOptionsGetBool(NULL, NULL, "-adapt", &adapt, NULL);CHKERRQ(ierr);
20*9b5c7ca5SMatthew G. Knepley   ierr = PetscOptionsGetBool(NULL, NULL, "-user_section", &userSection, NULL);CHKERRQ(ierr);
21*9b5c7ca5SMatthew G. Knepley 
22*9b5c7ca5SMatthew G. Knepley   /* Create a base DMPlex mesh */
23*9b5c7ca5SMatthew G. Knepley   ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, 2, PETSC_FALSE, NULL, NULL, NULL, NULL, PETSC_TRUE, &base);CHKERRQ(ierr);
24*9b5c7ca5SMatthew G. Knepley   ierr = DMSetFromOptions(base);CHKERRQ(ierr);
25*9b5c7ca5SMatthew G. Knepley   ierr = DMViewFromOptions(base, NULL, "-dm_view");CHKERRQ(ierr);
26*9b5c7ca5SMatthew G. Knepley 
27*9b5c7ca5SMatthew G. Knepley   /* Covert Plex mesh to Forest and destroy base */
28*9b5c7ca5SMatthew G. Knepley   ierr = DMCreate(PETSC_COMM_WORLD, &forest);CHKERRQ(ierr);
29*9b5c7ca5SMatthew G. Knepley   ierr = DMSetType(forest, DMP4EST);CHKERRQ(ierr);
30*9b5c7ca5SMatthew G. Knepley   ierr = DMForestSetBaseDM(forest, base);CHKERRQ(ierr);
31*9b5c7ca5SMatthew G. Knepley   ierr = DMSetUp(forest);CHKERRQ(ierr);
32*9b5c7ca5SMatthew G. Knepley   ierr = DMDestroy(&base);CHKERRQ(ierr);
33*9b5c7ca5SMatthew G. Knepley   ierr = DMViewFromOptions(forest, NULL, "-dm_view");CHKERRQ(ierr);
34*9b5c7ca5SMatthew G. Knepley 
35*9b5c7ca5SMatthew G. Knepley   if (adapt) {
36*9b5c7ca5SMatthew G. Knepley     /* Adaptively refine the cell 0 of the mesh */
37*9b5c7ca5SMatthew G. Knepley     for (i = 0; i < 3; ++i) {
38*9b5c7ca5SMatthew G. Knepley       DM      postforest;
39*9b5c7ca5SMatthew G. Knepley       DMLabel adaptLabel = NULL;
40*9b5c7ca5SMatthew G. Knepley 
41*9b5c7ca5SMatthew G. Knepley       ierr = DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel);CHKERRQ(ierr);
42*9b5c7ca5SMatthew G. Knepley       ierr = DMLabelSetValue(adaptLabel, 0, DM_ADAPT_REFINE);CHKERRQ(ierr);
43*9b5c7ca5SMatthew G. Knepley       ierr = DMForestTemplate(forest, PETSC_COMM_WORLD, &postforest);CHKERRQ(ierr);
44*9b5c7ca5SMatthew G. Knepley       ierr = DMForestSetAdaptivityLabel(postforest, adaptLabel);CHKERRQ(ierr);
45*9b5c7ca5SMatthew G. Knepley       ierr = DMLabelDestroy(&adaptLabel);CHKERRQ(ierr);
46*9b5c7ca5SMatthew G. Knepley       ierr = DMSetUp(postforest);CHKERRQ(ierr);
47*9b5c7ca5SMatthew G. Knepley       ierr = DMDestroy(&forest);CHKERRQ(ierr);
48*9b5c7ca5SMatthew G. Knepley       forest = postforest;
49*9b5c7ca5SMatthew G. Knepley     }
50*9b5c7ca5SMatthew G. Knepley   } else {
51*9b5c7ca5SMatthew G. Knepley     /* Adaptively refine all cells of the mesh */
52*9b5c7ca5SMatthew G. Knepley     PetscInt cStart, cEnd, c;
53*9b5c7ca5SMatthew G. Knepley 
54*9b5c7ca5SMatthew G. Knepley     for (i = 0; i < 3; ++i) {
55*9b5c7ca5SMatthew G. Knepley       DM      postforest;
56*9b5c7ca5SMatthew G. Knepley       DMLabel adaptLabel = NULL;
57*9b5c7ca5SMatthew G. Knepley 
58*9b5c7ca5SMatthew G. Knepley       ierr = DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel);CHKERRQ(ierr);
59*9b5c7ca5SMatthew G. Knepley 
60*9b5c7ca5SMatthew G. Knepley       ierr = DMForestGetCellChart(forest, &cStart, &cEnd);CHKERRQ(ierr);
61*9b5c7ca5SMatthew G. Knepley       for (c = cStart; c < cEnd; ++c) {
62*9b5c7ca5SMatthew G. Knepley         ierr = DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE);CHKERRQ(ierr);
63*9b5c7ca5SMatthew G. Knepley       }
64*9b5c7ca5SMatthew G. Knepley 
65*9b5c7ca5SMatthew G. Knepley       ierr = DMForestTemplate(forest, PETSC_COMM_WORLD, &postforest);CHKERRQ(ierr);
66*9b5c7ca5SMatthew G. Knepley       ierr = DMForestSetAdaptivityLabel(postforest, adaptLabel);CHKERRQ(ierr);
67*9b5c7ca5SMatthew G. Knepley       ierr = DMLabelDestroy(&adaptLabel);CHKERRQ(ierr);
68*9b5c7ca5SMatthew G. Knepley       ierr = DMSetUp(postforest);CHKERRQ(ierr);
69*9b5c7ca5SMatthew G. Knepley       ierr = DMDestroy(&forest);CHKERRQ(ierr);
70*9b5c7ca5SMatthew G. Knepley       forest = postforest;
71*9b5c7ca5SMatthew G. Knepley     }
72*9b5c7ca5SMatthew G. Knepley   }
73*9b5c7ca5SMatthew G. Knepley   ierr = DMViewFromOptions(forest, NULL, "-dm_view");CHKERRQ(ierr);
74*9b5c7ca5SMatthew G. Knepley 
75*9b5c7ca5SMatthew G. Knepley   /*  Setup the section*/
76*9b5c7ca5SMatthew G. Knepley   if (userSection) {
77*9b5c7ca5SMatthew G. Knepley     ierr = DMConvert(forest, DMPLEX, &plex);CHKERRQ(ierr);
78*9b5c7ca5SMatthew G. Knepley     ierr = DMViewFromOptions(plex, NULL, "-dm_view");CHKERRQ(ierr);
79*9b5c7ca5SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd);CHKERRQ(ierr);
80*9b5c7ca5SMatthew G. Knepley     ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "Vertices [%D, %D)\n", vStart, vEnd);CHKERRQ(ierr);
81*9b5c7ca5SMatthew G. Knepley     ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL);CHKERRQ(ierr);
82*9b5c7ca5SMatthew G. Knepley     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) forest), &s);CHKERRQ(ierr);
83*9b5c7ca5SMatthew G. Knepley     ierr = PetscSectionSetNumFields(s, 1);CHKERRQ(ierr);
84*9b5c7ca5SMatthew G. Knepley     ierr = PetscSectionSetChart(s, vStart, vEnd);CHKERRQ(ierr);
85*9b5c7ca5SMatthew G. Knepley     for (v = vStart; v < vEnd; ++v) {
86*9b5c7ca5SMatthew G. Knepley       ierr = PetscSectionSetDof(s, v, 1);CHKERRQ(ierr);
87*9b5c7ca5SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(s, v, 0, 1);CHKERRQ(ierr);
88*9b5c7ca5SMatthew G. Knepley     }
89*9b5c7ca5SMatthew G. Knepley     ierr = PetscSectionSetUp(s);CHKERRQ(ierr);
90*9b5c7ca5SMatthew G. Knepley     ierr = DMSetLocalSection(forest, s);CHKERRQ(ierr);
91*9b5c7ca5SMatthew G. Knepley     ierr = PetscObjectViewFromOptions((PetscObject) s, NULL, "-my_section_view");CHKERRQ(ierr);
92*9b5c7ca5SMatthew G. Knepley     ierr = PetscSectionDestroy(&s);CHKERRQ(ierr);
93*9b5c7ca5SMatthew G. Knepley     ierr = DMDestroy(&plex);CHKERRQ(ierr);
94*9b5c7ca5SMatthew G. Knepley   } else {
95*9b5c7ca5SMatthew G. Knepley     PetscFE  fe;
96*9b5c7ca5SMatthew G. Knepley     PetscInt dim;
97*9b5c7ca5SMatthew G. Knepley 
98*9b5c7ca5SMatthew G. Knepley     ierr = DMGetDimension(forest, &dim);CHKERRQ(ierr);
99*9b5c7ca5SMatthew G. Knepley     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
100*9b5c7ca5SMatthew G. Knepley     ierr = DMAddField(forest, NULL, (PetscObject) fe);CHKERRQ(ierr);
101*9b5c7ca5SMatthew G. Knepley     ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
102*9b5c7ca5SMatthew G. Knepley     ierr = DMCreateDS(forest);CHKERRQ(ierr);
103*9b5c7ca5SMatthew G. Knepley   }
104*9b5c7ca5SMatthew G. Knepley 
105*9b5c7ca5SMatthew G. Knepley   /* Create the global vector*/
106*9b5c7ca5SMatthew G. Knepley   ierr = DMCreateGlobalVector(forest, &g);CHKERRQ(ierr);
107*9b5c7ca5SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) g, "g");CHKERRQ(ierr);
108*9b5c7ca5SMatthew G. Knepley   ierr = VecSet(g, 1.0);CHKERRQ(ierr);
109*9b5c7ca5SMatthew G. Knepley 
110*9b5c7ca5SMatthew G. Knepley   /* Test global to local*/
111*9b5c7ca5SMatthew G. Knepley   Vec l;
112*9b5c7ca5SMatthew G. Knepley   ierr = DMCreateLocalVector(forest, &l);CHKERRQ(ierr);
113*9b5c7ca5SMatthew G. Knepley   ierr = VecZeroEntries(l);CHKERRQ(ierr);
114*9b5c7ca5SMatthew G. Knepley   ierr = DMGlobalToLocal(forest, g, INSERT_VALUES, l);CHKERRQ(ierr);
115*9b5c7ca5SMatthew G. Knepley   ierr = VecZeroEntries(g);CHKERRQ(ierr);
116*9b5c7ca5SMatthew G. Knepley   ierr = DMLocalToGlobal(forest, l, INSERT_VALUES, g);CHKERRQ(ierr);
117*9b5c7ca5SMatthew G. Knepley   ierr = VecDestroy(&l);CHKERRQ(ierr);
118*9b5c7ca5SMatthew G. Knepley 
119*9b5c7ca5SMatthew G. Knepley   /*  Save a vector*/
120*9b5c7ca5SMatthew G. Knepley   ierr = PetscViewerHDF5Open(PETSC_COMM_WORLD, "forestHDF.h5", FILE_MODE_WRITE, &viewer);
121*9b5c7ca5SMatthew G. Knepley   ierr = VecView(g, viewer);CHKERRQ(ierr);
122*9b5c7ca5SMatthew G. Knepley   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
123*9b5c7ca5SMatthew G. Knepley 
124*9b5c7ca5SMatthew G. Knepley   /* Load another vector to load into*/
125*9b5c7ca5SMatthew G. Knepley   ierr = DMCreateGlobalVector(forest, &g2);CHKERRQ(ierr);
126*9b5c7ca5SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) g2, "g");CHKERRQ(ierr);
127*9b5c7ca5SMatthew G. Knepley   ierr = VecZeroEntries(g2);CHKERRQ(ierr);
128*9b5c7ca5SMatthew G. Knepley 
129*9b5c7ca5SMatthew G. Knepley   /*  Load a vector*/
130*9b5c7ca5SMatthew G. Knepley   ierr = PetscViewerHDF5Open(PETSC_COMM_WORLD, "forestHDF.h5", FILE_MODE_READ, &viewer);
131*9b5c7ca5SMatthew G. Knepley   ierr = VecLoad(g2, viewer);CHKERRQ(ierr);
132*9b5c7ca5SMatthew G. Knepley   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
133*9b5c7ca5SMatthew G. Knepley 
134*9b5c7ca5SMatthew G. Knepley   /*  Check if the data is the same*/
135*9b5c7ca5SMatthew G. Knepley   ierr = VecAXPY(g2, -1.0, g);CHKERRQ(ierr);
136*9b5c7ca5SMatthew G. Knepley   ierr = VecNorm(g2, NORM_INFINITY, &nrm);CHKERRQ(ierr);
137*9b5c7ca5SMatthew G. Knepley   if (PetscAbsReal(nrm) > PETSC_SMALL) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Invalid difference norm %g", (double) nrm);
138*9b5c7ca5SMatthew G. Knepley 
139*9b5c7ca5SMatthew G. Knepley   ierr = VecDestroy(&g);CHKERRQ(ierr);
140*9b5c7ca5SMatthew G. Knepley   ierr = VecDestroy(&g2);CHKERRQ(ierr);
141*9b5c7ca5SMatthew G. Knepley   ierr = DMDestroy(&forest);CHKERRQ(ierr);
142*9b5c7ca5SMatthew G. Knepley   ierr = PetscFinalize();
143*9b5c7ca5SMatthew G. Knepley   return ierr;
144*9b5c7ca5SMatthew G. Knepley }
145*9b5c7ca5SMatthew G. Knepley 
146*9b5c7ca5SMatthew G. Knepley /*TEST
147*9b5c7ca5SMatthew G. Knepley 
148*9b5c7ca5SMatthew G. Knepley   build:
149*9b5c7ca5SMatthew G. Knepley     requires: hdf5 p4est
150*9b5c7ca5SMatthew G. Knepley 
151*9b5c7ca5SMatthew G. Knepley   test:
152*9b5c7ca5SMatthew G. Knepley     suffix: 0
153*9b5c7ca5SMatthew G. Knepley     nsize: {{1 2 5}}
154*9b5c7ca5SMatthew G. Knepley     args: -adapt -dm_plex_box_faces 2,2
155*9b5c7ca5SMatthew G. Knepley 
156*9b5c7ca5SMatthew G. Knepley TEST*/
157