xref: /petsc/src/dm/impls/forest/tests/ex1.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   CHKERRQ(PetscInitialize(&argc, &argv, NULL, help));
20 
21   CHKERRQ(DMCreate(PETSC_COMM_WORLD, &base));
22   CHKERRQ(DMSetType(base, DMPLEX));
23   CHKERRQ(DMSetFromOptions(base));
24 
25   CHKERRQ(DMCreate(PETSC_COMM_WORLD, &forest));
26   CHKERRQ(DMSetType(forest, DMP4EST));
27   CHKERRQ(DMForestSetBaseDM(forest, base));
28   CHKERRQ(DMForestSetInitialRefinement(forest, min_refine));
29   CHKERRQ(DMForestSetPartitionOverlap(forest, overlap));
30   CHKERRQ(DMSetUp(forest));
31   CHKERRQ(DMDestroy(&base));
32   CHKERRQ(DMViewFromOptions(forest, NULL, "-dm_view"));
33 
34   CHKERRQ(DMConvert(forest, DMPLEX, &plex));
35   CHKERRQ(DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd));
36   CHKERRQ(DMDestroy(&plex));
37   CHKERRQ(PetscSectionCreate(PetscObjectComm((PetscObject) forest), &s));
38   CHKERRQ(PetscSectionSetChart(s, vStart, vEnd));
39   for (v = vStart; v < vEnd; ++v) CHKERRQ(PetscSectionSetDof(s, v, 1));
40   CHKERRQ(PetscSectionSetUp(s));
41   CHKERRQ(DMSetLocalSection(forest, s));
42   CHKERRQ(PetscSectionDestroy(&s));
43 
44   CHKERRQ(DMCreateGlobalVector(forest, &g));
45   CHKERRQ(PetscObjectSetName((PetscObject) g, "g"));
46   CHKERRQ(VecSet(g, 1.0));
47   CHKERRQ(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forest.h5", FILE_MODE_WRITE, &viewer));
48   CHKERRQ(VecView(g, viewer));
49   CHKERRQ(PetscViewerDestroy(&viewer));
50 
51   CHKERRQ(DMCreateGlobalVector(forest, &g2));
52   CHKERRQ(PetscObjectSetName((PetscObject) g2, "g"));
53   CHKERRQ(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forest.h5", FILE_MODE_READ, &viewer));
54   CHKERRQ(VecLoad(g2, viewer));
55   CHKERRQ(PetscViewerDestroy(&viewer));
56 
57   /*  Check if the data is the same*/
58   CHKERRQ(VecAXPY(g2, -1.0, g));
59   CHKERRQ(VecNorm(g2, NORM_INFINITY, &diff));
60   if (diff > PETSC_MACHINE_EPSILON) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Check failed: %g\n", (double) diff));
61 
62   CHKERRQ(VecDestroy(&g));
63   CHKERRQ(VecDestroy(&g2));
64   CHKERRQ(DMDestroy(&forest));
65   CHKERRQ(PetscFinalize());
66   return 0;
67 }
68 
69 /*TEST
70 
71   build:
72     requires: hdf5 p4est
73 
74   test:
75     args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3
76 
77 TEST*/
78