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