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