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 80bb0ac0dSMatthew G. Knepley int main (int argc, char **argv) 90bb0ac0dSMatthew G. Knepley { 100bb0ac0dSMatthew G. Knepley 110bb0ac0dSMatthew G. Knepley DM base, forest, plex; 120bb0ac0dSMatthew G. Knepley Vec g, g2; 130bb0ac0dSMatthew G. Knepley PetscSection s; 140bb0ac0dSMatthew G. Knepley PetscViewer viewer; 150bb0ac0dSMatthew G. Knepley PetscReal diff; 160bb0ac0dSMatthew G. Knepley PetscInt min_refine = 2, overlap = 0; 170bb0ac0dSMatthew G. Knepley PetscInt vStart, vEnd, v; 180bb0ac0dSMatthew G. Knepley 19*327415f7SBarry Smith PetscFunctionBeginUser; 209566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 210bb0ac0dSMatthew G. Knepley 229566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &base)); 239566063dSJacob Faibussowitsch PetscCall(DMSetType(base, DMPLEX)); 249566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(base)); 250bb0ac0dSMatthew G. Knepley 269566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &forest)); 279566063dSJacob Faibussowitsch PetscCall(DMSetType(forest, DMP4EST)); 289566063dSJacob Faibussowitsch PetscCall(DMForestSetBaseDM(forest, base)); 299566063dSJacob Faibussowitsch PetscCall(DMForestSetInitialRefinement(forest, min_refine)); 309566063dSJacob Faibussowitsch PetscCall(DMForestSetPartitionOverlap(forest, overlap)); 319566063dSJacob Faibussowitsch PetscCall(DMSetUp(forest)); 329566063dSJacob Faibussowitsch PetscCall(DMDestroy(&base)); 339566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(forest, NULL, "-dm_view")); 340bb0ac0dSMatthew G. Knepley 359566063dSJacob Faibussowitsch PetscCall(DMConvert(forest, DMPLEX, &plex)); 369566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd)); 379566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex)); 389566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) forest), &s)); 399566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(s, vStart, vEnd)); 409566063dSJacob Faibussowitsch for (v = vStart; v < vEnd; ++v) PetscCall(PetscSectionSetDof(s, v, 1)); 419566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(s)); 429566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(forest, s)); 439566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&s)); 440bb0ac0dSMatthew G. Knepley 459566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(forest, &g)); 469566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) g, "g")); 479566063dSJacob Faibussowitsch PetscCall(VecSet(g, 1.0)); 489566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forest.h5", FILE_MODE_WRITE, &viewer)); 499566063dSJacob Faibussowitsch PetscCall(VecView(g, viewer)); 509566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 510bb0ac0dSMatthew G. Knepley 529566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(forest, &g2)); 539566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) g2, "g")); 549566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forest.h5", FILE_MODE_READ, &viewer)); 559566063dSJacob Faibussowitsch PetscCall(VecLoad(g2, viewer)); 569566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 570bb0ac0dSMatthew G. Knepley 580bb0ac0dSMatthew G. Knepley /* Check if the data is the same*/ 599566063dSJacob Faibussowitsch PetscCall(VecAXPY(g2, -1.0, g)); 609566063dSJacob Faibussowitsch PetscCall(VecNorm(g2, NORM_INFINITY, &diff)); 619566063dSJacob Faibussowitsch if (diff > PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Check failed: %g\n", (double) diff)); 620bb0ac0dSMatthew G. Knepley 639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&g)); 649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&g2)); 659566063dSJacob Faibussowitsch PetscCall(DMDestroy(&forest)); 669566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 67b122ec5aSJacob Faibussowitsch return 0; 680bb0ac0dSMatthew G. Knepley } 690bb0ac0dSMatthew G. Knepley 700bb0ac0dSMatthew G. Knepley /*TEST 710bb0ac0dSMatthew G. Knepley 720bb0ac0dSMatthew G. Knepley build: 730bb0ac0dSMatthew G. Knepley requires: hdf5 p4est 740bb0ac0dSMatthew G. Knepley 750bb0ac0dSMatthew G. Knepley test: 7630602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3 770bb0ac0dSMatthew G. Knepley 780bb0ac0dSMatthew G. Knepley TEST*/ 79