xref: /petsc/src/dm/impls/forest/tests/ex1.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
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 
main(int argc,char ** argv)8d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
9d71ae5a4SJacob Faibussowitsch {
100bb0ac0dSMatthew G. Knepley   DM           base, forest, plex;
110bb0ac0dSMatthew G. Knepley   Vec          g, g2;
120bb0ac0dSMatthew G. Knepley   PetscSection s;
130bb0ac0dSMatthew G. Knepley   PetscViewer  viewer;
140bb0ac0dSMatthew G. Knepley   PetscReal    diff;
150bb0ac0dSMatthew G. Knepley   PetscInt     min_refine = 2, overlap = 0;
160bb0ac0dSMatthew G. Knepley   PetscInt     vStart, vEnd, v;
170bb0ac0dSMatthew G. Knepley 
18327415f7SBarry Smith   PetscFunctionBeginUser;
199566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
200bb0ac0dSMatthew G. Knepley 
219566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &base));
229566063dSJacob Faibussowitsch   PetscCall(DMSetType(base, DMPLEX));
239566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(base));
240bb0ac0dSMatthew G. Knepley 
259566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &forest));
269566063dSJacob Faibussowitsch   PetscCall(DMSetType(forest, DMP4EST));
279566063dSJacob Faibussowitsch   PetscCall(DMForestSetBaseDM(forest, base));
289566063dSJacob Faibussowitsch   PetscCall(DMForestSetInitialRefinement(forest, min_refine));
299566063dSJacob Faibussowitsch   PetscCall(DMForestSetPartitionOverlap(forest, overlap));
309566063dSJacob Faibussowitsch   PetscCall(DMSetUp(forest));
319566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&base));
329566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(forest, NULL, "-dm_view"));
330bb0ac0dSMatthew G. Knepley 
349566063dSJacob Faibussowitsch   PetscCall(DMConvert(forest, DMPLEX, &plex));
359566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd));
369566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
379566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)forest), &s));
389566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(s, vStart, vEnd));
399566063dSJacob Faibussowitsch   for (v = vStart; v < vEnd; ++v) PetscCall(PetscSectionSetDof(s, v, 1));
409566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(s));
419566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(forest, s));
429566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&s));
430bb0ac0dSMatthew G. Knepley 
449566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(forest, &g));
459566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)g, "g"));
469566063dSJacob Faibussowitsch   PetscCall(VecSet(g, 1.0));
479566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forest.h5", FILE_MODE_WRITE, &viewer));
489566063dSJacob Faibussowitsch   PetscCall(VecView(g, viewer));
499566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
500bb0ac0dSMatthew G. Knepley 
519566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(forest, &g2));
529566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)g2, "g"));
539566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forest.h5", FILE_MODE_READ, &viewer));
549566063dSJacob Faibussowitsch   PetscCall(VecLoad(g2, viewer));
559566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
560bb0ac0dSMatthew G. Knepley 
570bb0ac0dSMatthew G. Knepley   /*  Check if the data is the same*/
589566063dSJacob Faibussowitsch   PetscCall(VecAXPY(g2, -1.0, g));
599566063dSJacob Faibussowitsch   PetscCall(VecNorm(g2, NORM_INFINITY, &diff));
609566063dSJacob Faibussowitsch   if (diff > PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Check failed: %g\n", (double)diff));
610bb0ac0dSMatthew G. Knepley 
629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&g));
639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&g2));
649566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&forest));
659566063dSJacob 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
76*3886731fSPierre Jolivet     output_file: output/empty.out
770bb0ac0dSMatthew G. Knepley 
780bb0ac0dSMatthew G. Knepley TEST*/
79