19b5c7ca5SMatthew G. Knepley static char help[] = "Tests adaptive refinement using DMForest, and uses HDF5.\n\n"; 29b5c7ca5SMatthew G. Knepley 39b5c7ca5SMatthew G. Knepley #include <petscdmforest.h> 49b5c7ca5SMatthew G. Knepley #include <petscdmplex.h> 59b5c7ca5SMatthew G. Knepley #include <petscviewerhdf5.h> 69b5c7ca5SMatthew G. Knepley 79b5c7ca5SMatthew G. Knepley int main (int argc, char **argv) 89b5c7ca5SMatthew G. Knepley { 99b5c7ca5SMatthew G. Knepley DM base, forest, plex; 109b5c7ca5SMatthew G. Knepley PetscSection s; 119b5c7ca5SMatthew G. Knepley PetscViewer viewer; 129b5c7ca5SMatthew G. Knepley Vec g = NULL, g2 = NULL; 139b5c7ca5SMatthew G. Knepley PetscReal nrm; 149b5c7ca5SMatthew G. Knepley PetscBool adapt = PETSC_FALSE, userSection = PETSC_FALSE; 159b5c7ca5SMatthew G. Knepley PetscInt vStart, vEnd, v, i; 169b5c7ca5SMatthew G. Knepley PetscErrorCode ierr; 179b5c7ca5SMatthew G. Knepley 189b5c7ca5SMatthew G. Knepley ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 19*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL, NULL, "-adapt", &adapt, NULL)); 20*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL, NULL, "-user_section", &userSection, NULL)); 219b5c7ca5SMatthew G. Knepley 229b5c7ca5SMatthew G. Knepley /* Create a base DMPlex mesh */ 23*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(PETSC_COMM_WORLD, &base)); 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(base, DMPLEX)); 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(base)); 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(base, NULL, "-dm_view")); 279b5c7ca5SMatthew G. Knepley 289b5c7ca5SMatthew G. Knepley /* Covert Plex mesh to Forest and destroy base */ 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(PETSC_COMM_WORLD, &forest)); 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(forest, DMP4EST)); 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMForestSetBaseDM(forest, base)); 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(forest)); 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&base)); 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(forest, NULL, "-dm_view")); 359b5c7ca5SMatthew G. Knepley 369b5c7ca5SMatthew G. Knepley if (adapt) { 379b5c7ca5SMatthew G. Knepley /* Adaptively refine the cell 0 of the mesh */ 389b5c7ca5SMatthew G. Knepley for (i = 0; i < 3; ++i) { 399b5c7ca5SMatthew G. Knepley DM postforest; 409b5c7ca5SMatthew G. Knepley DMLabel adaptLabel = NULL; 419b5c7ca5SMatthew G. Knepley 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel)); 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(adaptLabel, 0, DM_ADAPT_REFINE)); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMForestTemplate(forest, PETSC_COMM_WORLD, &postforest)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMForestSetAdaptivityLabel(postforest, adaptLabel)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelDestroy(&adaptLabel)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(postforest)); 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&forest)); 499b5c7ca5SMatthew G. Knepley forest = postforest; 509b5c7ca5SMatthew G. Knepley } 519b5c7ca5SMatthew G. Knepley } else { 529b5c7ca5SMatthew G. Knepley /* Adaptively refine all cells of the mesh */ 539b5c7ca5SMatthew G. Knepley PetscInt cStart, cEnd, c; 549b5c7ca5SMatthew G. Knepley 559b5c7ca5SMatthew G. Knepley for (i = 0; i < 3; ++i) { 569b5c7ca5SMatthew G. Knepley DM postforest; 579b5c7ca5SMatthew G. Knepley DMLabel adaptLabel = NULL; 589b5c7ca5SMatthew G. Knepley 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel)); 609b5c7ca5SMatthew G. Knepley 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMForestGetCellChart(forest, &cStart, &cEnd)); 629b5c7ca5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE)); 649b5c7ca5SMatthew G. Knepley } 659b5c7ca5SMatthew G. Knepley 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMForestTemplate(forest, PETSC_COMM_WORLD, &postforest)); 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMForestSetAdaptivityLabel(postforest, adaptLabel)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelDestroy(&adaptLabel)); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(postforest)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&forest)); 719b5c7ca5SMatthew G. Knepley forest = postforest; 729b5c7ca5SMatthew G. Knepley } 739b5c7ca5SMatthew G. Knepley } 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(forest, NULL, "-dm_view")); 759b5c7ca5SMatthew G. Knepley 769b5c7ca5SMatthew G. Knepley /* Setup the section*/ 779b5c7ca5SMatthew G. Knepley if (userSection) { 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMConvert(forest, DMPLEX, &plex)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(plex, NULL, "-dm_view")); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd)); 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "Vertices [%D, %D)\n", vStart, vEnd)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionCreate(PetscObjectComm((PetscObject) forest), &s)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetNumFields(s, 1)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetChart(s, vStart, vEnd)); 869b5c7ca5SMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetDof(s, v, 1)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldDof(s, v, 0, 1)); 899b5c7ca5SMatthew G. Knepley } 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetUp(s)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetLocalSection(forest, s)); 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectViewFromOptions((PetscObject) s, NULL, "-my_section_view")); 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionDestroy(&s)); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&plex)); 959b5c7ca5SMatthew G. Knepley } else { 969b5c7ca5SMatthew G. Knepley PetscFE fe; 979b5c7ca5SMatthew G. Knepley PetscInt dim; 989b5c7ca5SMatthew G. Knepley 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(forest, &dim)); 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, 1, PETSC_DETERMINE, &fe)); 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddField(forest, NULL, (PetscObject) fe)); 102*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(forest)); 1049b5c7ca5SMatthew G. Knepley } 1059b5c7ca5SMatthew G. Knepley 1069b5c7ca5SMatthew G. Knepley /* Create the global vector*/ 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(forest, &g)); 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) g, "g")); 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(g, 1.0)); 1109b5c7ca5SMatthew G. Knepley 1119b5c7ca5SMatthew G. Knepley /* Test global to local*/ 1129b5c7ca5SMatthew G. Knepley Vec l; 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(forest, &l)); 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(l)); 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocal(forest, g, INSERT_VALUES, l)); 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(g)); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobal(forest, l, INSERT_VALUES, g)); 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&l)); 1199b5c7ca5SMatthew G. Knepley 1209b5c7ca5SMatthew G. Knepley /* Save a vector*/ 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forestHDF.h5", FILE_MODE_WRITE, &viewer)); 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(g, viewer)); 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 1249b5c7ca5SMatthew G. Knepley 1259b5c7ca5SMatthew G. Knepley /* Load another vector to load into*/ 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(forest, &g2)); 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) g2, "g")); 128*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(g2)); 1299b5c7ca5SMatthew G. Knepley 1309b5c7ca5SMatthew G. Knepley /* Load a vector*/ 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forestHDF.h5", FILE_MODE_READ, &viewer)); 132*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecLoad(g2, viewer)); 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 1349b5c7ca5SMatthew G. Knepley 1359b5c7ca5SMatthew G. Knepley /* Check if the data is the same*/ 136*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(g2, -1.0, g)); 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(g2, NORM_INFINITY, &nrm)); 1382c71b3e2SJacob Faibussowitsch PetscCheckFalse(PetscAbsReal(nrm) > PETSC_SMALL,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Invalid difference norm %g", (double) nrm); 1399b5c7ca5SMatthew G. Knepley 140*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&g)); 141*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&g2)); 142*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&forest)); 1439b5c7ca5SMatthew G. Knepley ierr = PetscFinalize(); 1449b5c7ca5SMatthew G. Knepley return ierr; 1459b5c7ca5SMatthew G. Knepley } 1469b5c7ca5SMatthew G. Knepley 1479b5c7ca5SMatthew G. Knepley /*TEST 1489b5c7ca5SMatthew G. Knepley 1499b5c7ca5SMatthew G. Knepley build: 1509b5c7ca5SMatthew G. Knepley requires: hdf5 p4est 1519b5c7ca5SMatthew G. Knepley 1529b5c7ca5SMatthew G. Knepley test: 1539b5c7ca5SMatthew G. Knepley suffix: 0 1549b5c7ca5SMatthew G. Knepley nsize: {{1 2 5}} 15530602db0SMatthew G. Knepley args: -adapt -dm_plex_simplex 0 -dm_plex_box_faces 2,2 1569b5c7ca5SMatthew G. Knepley 1579b5c7ca5SMatthew G. Knepley TEST*/ 158