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*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL, help)); 200bb0ac0dSMatthew G. Knepley 215f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(PETSC_COMM_WORLD, &base)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(base, DMPLEX)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(base)); 240bb0ac0dSMatthew G. Knepley 255f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(PETSC_COMM_WORLD, &forest)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(forest, DMP4EST)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(DMForestSetBaseDM(forest, base)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(DMForestSetInitialRefinement(forest, min_refine)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(DMForestSetPartitionOverlap(forest, overlap)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(forest)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&base)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(forest, NULL, "-dm_view")); 330bb0ac0dSMatthew G. Knepley 345f80ce2aSJacob Faibussowitsch CHKERRQ(DMConvert(forest, DMPLEX, &plex)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&plex)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionCreate(PetscObjectComm((PetscObject) forest), &s)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetChart(s, vStart, vEnd)); 395f80ce2aSJacob Faibussowitsch for (v = vStart; v < vEnd; ++v) CHKERRQ(PetscSectionSetDof(s, v, 1)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetUp(s)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetLocalSection(forest, s)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionDestroy(&s)); 430bb0ac0dSMatthew G. Knepley 445f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(forest, &g)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) g, "g")); 465f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(g, 1.0)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forest.h5", FILE_MODE_WRITE, &viewer)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(g, viewer)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 500bb0ac0dSMatthew G. Knepley 515f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(forest, &g2)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) g2, "g")); 535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forest.h5", FILE_MODE_READ, &viewer)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(VecLoad(g2, viewer)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 560bb0ac0dSMatthew G. Knepley 570bb0ac0dSMatthew G. Knepley /* Check if the data is the same*/ 585f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(g2, -1.0, g)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(g2, NORM_INFINITY, &diff)); 605f80ce2aSJacob Faibussowitsch if (diff > PETSC_MACHINE_EPSILON) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Check failed: %g\n", (double) diff)); 610bb0ac0dSMatthew G. Knepley 625f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&g)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&g2)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&forest)); 65*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 66*b122ec5aSJacob 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 760bb0ac0dSMatthew G. Knepley 770bb0ac0dSMatthew G. Knepley TEST*/ 78