xref: /petsc/src/dm/impls/forest/tests/ex1.c (revision 9566063d113dddea24716c546802770db7481bc0)
10bb0ac0dSMatthew G. Knepley static char help[] = "Test HDF5 input and output.\n\
20bb0ac0dSMatthew G. Knepley This exposed a bug with sharing discretizations.\n\n\n";
30bb0ac0dSMatthew G. Knepley 
40bb0ac0dSMatthew G. Knepley #include <petscdmforest.h>
50bb0ac0dSMatthew G. Knepley #include <petscdmplex.h>
60bb0ac0dSMatthew G. Knepley #include <petscviewerhdf5.h>
70bb0ac0dSMatthew G. Knepley 
80bb0ac0dSMatthew G. Knepley int main (int argc, char **argv)
90bb0ac0dSMatthew G. Knepley {
100bb0ac0dSMatthew G. Knepley 
110bb0ac0dSMatthew G. Knepley   DM             base, forest, plex;
120bb0ac0dSMatthew G. Knepley   Vec            g, g2;
130bb0ac0dSMatthew G. Knepley   PetscSection   s;
140bb0ac0dSMatthew G. Knepley   PetscViewer    viewer;
150bb0ac0dSMatthew G. Knepley   PetscReal      diff;
160bb0ac0dSMatthew G. Knepley   PetscInt       min_refine = 2, overlap = 0;
170bb0ac0dSMatthew G. Knepley   PetscInt       vStart, vEnd, v;
180bb0ac0dSMatthew G. Knepley 
19*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
200bb0ac0dSMatthew G. Knepley 
21*9566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &base));
22*9566063dSJacob Faibussowitsch   PetscCall(DMSetType(base, DMPLEX));
23*9566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(base));
240bb0ac0dSMatthew G. Knepley 
25*9566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &forest));
26*9566063dSJacob Faibussowitsch   PetscCall(DMSetType(forest, DMP4EST));
27*9566063dSJacob Faibussowitsch   PetscCall(DMForestSetBaseDM(forest, base));
28*9566063dSJacob Faibussowitsch   PetscCall(DMForestSetInitialRefinement(forest, min_refine));
29*9566063dSJacob Faibussowitsch   PetscCall(DMForestSetPartitionOverlap(forest, overlap));
30*9566063dSJacob Faibussowitsch   PetscCall(DMSetUp(forest));
31*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&base));
32*9566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(forest, NULL, "-dm_view"));
330bb0ac0dSMatthew G. Knepley 
34*9566063dSJacob Faibussowitsch   PetscCall(DMConvert(forest, DMPLEX, &plex));
35*9566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd));
36*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
37*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) forest), &s));
38*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(s, vStart, vEnd));
39*9566063dSJacob Faibussowitsch   for (v = vStart; v < vEnd; ++v) PetscCall(PetscSectionSetDof(s, v, 1));
40*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(s));
41*9566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(forest, s));
42*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&s));
430bb0ac0dSMatthew G. Knepley 
44*9566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(forest, &g));
45*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) g, "g"));
46*9566063dSJacob Faibussowitsch   PetscCall(VecSet(g, 1.0));
47*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forest.h5", FILE_MODE_WRITE, &viewer));
48*9566063dSJacob Faibussowitsch   PetscCall(VecView(g, viewer));
49*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
500bb0ac0dSMatthew G. Knepley 
51*9566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(forest, &g2));
52*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) g2, "g"));
53*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forest.h5", FILE_MODE_READ, &viewer));
54*9566063dSJacob Faibussowitsch   PetscCall(VecLoad(g2, viewer));
55*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
560bb0ac0dSMatthew G. Knepley 
570bb0ac0dSMatthew G. Knepley   /*  Check if the data is the same*/
58*9566063dSJacob Faibussowitsch   PetscCall(VecAXPY(g2, -1.0, g));
59*9566063dSJacob Faibussowitsch   PetscCall(VecNorm(g2, NORM_INFINITY, &diff));
60*9566063dSJacob Faibussowitsch   if (diff > PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Check failed: %g\n", (double) diff));
610bb0ac0dSMatthew G. Knepley 
62*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&g));
63*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&g2));
64*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&forest));
65*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
66b122ec5aSJacob Faibussowitsch   return 0;
670bb0ac0dSMatthew G. Knepley }
680bb0ac0dSMatthew G. Knepley 
690bb0ac0dSMatthew G. Knepley /*TEST
700bb0ac0dSMatthew G. Knepley 
710bb0ac0dSMatthew G. Knepley   build:
720bb0ac0dSMatthew G. Knepley     requires: hdf5 p4est
730bb0ac0dSMatthew G. Knepley 
740bb0ac0dSMatthew G. Knepley   test:
7530602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3
760bb0ac0dSMatthew G. Knepley 
770bb0ac0dSMatthew G. Knepley TEST*/
78