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 17*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL, help)); 185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL, NULL, "-adapt", &adapt, NULL)); 195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL, NULL, "-user_section", &userSection, NULL)); 209b5c7ca5SMatthew G. Knepley 219b5c7ca5SMatthew G. Knepley /* Create a base DMPlex mesh */ 225f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(PETSC_COMM_WORLD, &base)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(base, DMPLEX)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(base)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(base, NULL, "-dm_view")); 269b5c7ca5SMatthew G. Knepley 279b5c7ca5SMatthew G. Knepley /* Covert Plex mesh to Forest and destroy base */ 285f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(PETSC_COMM_WORLD, &forest)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(forest, DMP4EST)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(DMForestSetBaseDM(forest, base)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(forest)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&base)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(forest, NULL, "-dm_view")); 349b5c7ca5SMatthew G. Knepley 359b5c7ca5SMatthew G. Knepley if (adapt) { 369b5c7ca5SMatthew G. Knepley /* Adaptively refine the cell 0 of the mesh */ 379b5c7ca5SMatthew G. Knepley for (i = 0; i < 3; ++i) { 389b5c7ca5SMatthew G. Knepley DM postforest; 399b5c7ca5SMatthew G. Knepley DMLabel adaptLabel = NULL; 409b5c7ca5SMatthew G. Knepley 415f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(adaptLabel, 0, DM_ADAPT_REFINE)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(DMForestTemplate(forest, PETSC_COMM_WORLD, &postforest)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(DMForestSetAdaptivityLabel(postforest, adaptLabel)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelDestroy(&adaptLabel)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(postforest)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&forest)); 489b5c7ca5SMatthew G. Knepley forest = postforest; 499b5c7ca5SMatthew G. Knepley } 509b5c7ca5SMatthew G. Knepley } else { 519b5c7ca5SMatthew G. Knepley /* Adaptively refine all cells of the mesh */ 529b5c7ca5SMatthew G. Knepley PetscInt cStart, cEnd, c; 539b5c7ca5SMatthew G. Knepley 549b5c7ca5SMatthew G. Knepley for (i = 0; i < 3; ++i) { 559b5c7ca5SMatthew G. Knepley DM postforest; 569b5c7ca5SMatthew G. Knepley DMLabel adaptLabel = NULL; 579b5c7ca5SMatthew G. Knepley 585f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel)); 599b5c7ca5SMatthew G. Knepley 605f80ce2aSJacob Faibussowitsch CHKERRQ(DMForestGetCellChart(forest, &cStart, &cEnd)); 619b5c7ca5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 625f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE)); 639b5c7ca5SMatthew G. Knepley } 649b5c7ca5SMatthew G. Knepley 655f80ce2aSJacob Faibussowitsch CHKERRQ(DMForestTemplate(forest, PETSC_COMM_WORLD, &postforest)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(DMForestSetAdaptivityLabel(postforest, adaptLabel)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelDestroy(&adaptLabel)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(postforest)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&forest)); 709b5c7ca5SMatthew G. Knepley forest = postforest; 719b5c7ca5SMatthew G. Knepley } 729b5c7ca5SMatthew G. Knepley } 735f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(forest, NULL, "-dm_view")); 749b5c7ca5SMatthew G. Knepley 759b5c7ca5SMatthew G. Knepley /* Setup the section*/ 769b5c7ca5SMatthew G. Knepley if (userSection) { 775f80ce2aSJacob Faibussowitsch CHKERRQ(DMConvert(forest, DMPLEX, &plex)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(plex, NULL, "-dm_view")); 795f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "Vertices [%D, %D)\n", vStart, vEnd)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionCreate(PetscObjectComm((PetscObject) forest), &s)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetNumFields(s, 1)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetChart(s, vStart, vEnd)); 859b5c7ca5SMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetDof(s, v, 1)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldDof(s, v, 0, 1)); 889b5c7ca5SMatthew G. Knepley } 895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetUp(s)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetLocalSection(forest, s)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectViewFromOptions((PetscObject) s, NULL, "-my_section_view")); 925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionDestroy(&s)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&plex)); 949b5c7ca5SMatthew G. Knepley } else { 959b5c7ca5SMatthew G. Knepley PetscFE fe; 969b5c7ca5SMatthew G. Knepley PetscInt dim; 979b5c7ca5SMatthew G. Knepley 985f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(forest, &dim)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, 1, PETSC_DETERMINE, &fe)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddField(forest, NULL, (PetscObject) fe)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(forest)); 1039b5c7ca5SMatthew G. Knepley } 1049b5c7ca5SMatthew G. Knepley 1059b5c7ca5SMatthew G. Knepley /* Create the global vector*/ 1065f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(forest, &g)); 1075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) g, "g")); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(g, 1.0)); 1099b5c7ca5SMatthew G. Knepley 1109b5c7ca5SMatthew G. Knepley /* Test global to local*/ 1119b5c7ca5SMatthew G. Knepley Vec l; 1125f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(forest, &l)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(l)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocal(forest, g, INSERT_VALUES, l)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(g)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobal(forest, l, INSERT_VALUES, g)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&l)); 1189b5c7ca5SMatthew G. Knepley 1199b5c7ca5SMatthew G. Knepley /* Save a vector*/ 1205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forestHDF.h5", FILE_MODE_WRITE, &viewer)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(g, viewer)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 1239b5c7ca5SMatthew G. Knepley 1249b5c7ca5SMatthew G. Knepley /* Load another vector to load into*/ 1255f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(forest, &g2)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) g2, "g")); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(g2)); 1289b5c7ca5SMatthew G. Knepley 1299b5c7ca5SMatthew G. Knepley /* Load a vector*/ 1305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forestHDF.h5", FILE_MODE_READ, &viewer)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(VecLoad(g2, viewer)); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 1339b5c7ca5SMatthew G. Knepley 1349b5c7ca5SMatthew G. Knepley /* Check if the data is the same*/ 1355f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(g2, -1.0, g)); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(g2, NORM_INFINITY, &nrm)); 1372c71b3e2SJacob Faibussowitsch PetscCheckFalse(PetscAbsReal(nrm) > PETSC_SMALL,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Invalid difference norm %g", (double) nrm); 1389b5c7ca5SMatthew G. Knepley 1395f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&g)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&g2)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&forest)); 142*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 143*b122ec5aSJacob Faibussowitsch return 0; 1449b5c7ca5SMatthew G. Knepley } 1459b5c7ca5SMatthew G. Knepley 1469b5c7ca5SMatthew G. Knepley /*TEST 1479b5c7ca5SMatthew G. Knepley 1489b5c7ca5SMatthew G. Knepley build: 1499b5c7ca5SMatthew G. Knepley requires: hdf5 p4est 1509b5c7ca5SMatthew G. Knepley 1519b5c7ca5SMatthew G. Knepley test: 1529b5c7ca5SMatthew G. Knepley suffix: 0 1539b5c7ca5SMatthew G. Knepley nsize: {{1 2 5}} 15430602db0SMatthew G. Knepley args: -adapt -dm_plex_simplex 0 -dm_plex_box_faces 2,2 1559b5c7ca5SMatthew G. Knepley 1569b5c7ca5SMatthew G. Knepley TEST*/ 157