xref: /petsc/src/dm/impls/forest/tests/ex1.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 static char help[] = "Test HDF5 input and output.\n\
2 This exposed a bug with sharing discretizations.\n\n\n";
3 
4 #include <petscdmforest.h>
5 #include <petscdmplex.h>
6 #include <petscviewerhdf5.h>
7 
8 int main (int argc, char **argv)
9 {
10 
11   DM             base, forest, plex;
12   Vec            g, g2;
13   PetscSection   s;
14   PetscViewer    viewer;
15   PetscReal      diff;
16   PetscInt       min_refine = 2, overlap = 0;
17   PetscInt       vStart, vEnd, v;
18 
19   PetscFunctionBeginUser;
20   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
21 
22   PetscCall(DMCreate(PETSC_COMM_WORLD, &base));
23   PetscCall(DMSetType(base, DMPLEX));
24   PetscCall(DMSetFromOptions(base));
25 
26   PetscCall(DMCreate(PETSC_COMM_WORLD, &forest));
27   PetscCall(DMSetType(forest, DMP4EST));
28   PetscCall(DMForestSetBaseDM(forest, base));
29   PetscCall(DMForestSetInitialRefinement(forest, min_refine));
30   PetscCall(DMForestSetPartitionOverlap(forest, overlap));
31   PetscCall(DMSetUp(forest));
32   PetscCall(DMDestroy(&base));
33   PetscCall(DMViewFromOptions(forest, NULL, "-dm_view"));
34 
35   PetscCall(DMConvert(forest, DMPLEX, &plex));
36   PetscCall(DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd));
37   PetscCall(DMDestroy(&plex));
38   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) forest), &s));
39   PetscCall(PetscSectionSetChart(s, vStart, vEnd));
40   for (v = vStart; v < vEnd; ++v) PetscCall(PetscSectionSetDof(s, v, 1));
41   PetscCall(PetscSectionSetUp(s));
42   PetscCall(DMSetLocalSection(forest, s));
43   PetscCall(PetscSectionDestroy(&s));
44 
45   PetscCall(DMCreateGlobalVector(forest, &g));
46   PetscCall(PetscObjectSetName((PetscObject) g, "g"));
47   PetscCall(VecSet(g, 1.0));
48   PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forest.h5", FILE_MODE_WRITE, &viewer));
49   PetscCall(VecView(g, viewer));
50   PetscCall(PetscViewerDestroy(&viewer));
51 
52   PetscCall(DMCreateGlobalVector(forest, &g2));
53   PetscCall(PetscObjectSetName((PetscObject) g2, "g"));
54   PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forest.h5", FILE_MODE_READ, &viewer));
55   PetscCall(VecLoad(g2, viewer));
56   PetscCall(PetscViewerDestroy(&viewer));
57 
58   /*  Check if the data is the same*/
59   PetscCall(VecAXPY(g2, -1.0, g));
60   PetscCall(VecNorm(g2, NORM_INFINITY, &diff));
61   if (diff > PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Check failed: %g\n", (double) diff));
62 
63   PetscCall(VecDestroy(&g));
64   PetscCall(VecDestroy(&g2));
65   PetscCall(DMDestroy(&forest));
66   PetscCall(PetscFinalize());
67   return 0;
68 }
69 
70 /*TEST
71 
72   build:
73     requires: hdf5 p4est
74 
75   test:
76     args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3
77 
78 TEST*/
79