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 179566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 189566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-adapt", &adapt, NULL)); 199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-user_section", &userSection, NULL)); 209b5c7ca5SMatthew G. Knepley 219b5c7ca5SMatthew G. Knepley /* Create a base DMPlex mesh */ 229566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &base)); 239566063dSJacob Faibussowitsch PetscCall(DMSetType(base, DMPLEX)); 249566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(base)); 259566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(base, NULL, "-dm_view")); 269b5c7ca5SMatthew G. Knepley 279b5c7ca5SMatthew G. Knepley /* Covert Plex mesh to Forest and destroy base */ 289566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &forest)); 299566063dSJacob Faibussowitsch PetscCall(DMSetType(forest, DMP4EST)); 309566063dSJacob Faibussowitsch PetscCall(DMForestSetBaseDM(forest, base)); 319566063dSJacob Faibussowitsch PetscCall(DMSetUp(forest)); 329566063dSJacob Faibussowitsch PetscCall(DMDestroy(&base)); 339566063dSJacob Faibussowitsch PetscCall(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 419566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel)); 429566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(adaptLabel, 0, DM_ADAPT_REFINE)); 439566063dSJacob Faibussowitsch PetscCall(DMForestTemplate(forest, PETSC_COMM_WORLD, &postforest)); 449566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityLabel(postforest, adaptLabel)); 459566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&adaptLabel)); 469566063dSJacob Faibussowitsch PetscCall(DMSetUp(postforest)); 479566063dSJacob Faibussowitsch PetscCall(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 589566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel)); 599b5c7ca5SMatthew G. Knepley 609566063dSJacob Faibussowitsch PetscCall(DMForestGetCellChart(forest, &cStart, &cEnd)); 619b5c7ca5SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 629566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE)); 639b5c7ca5SMatthew G. Knepley } 649b5c7ca5SMatthew G. Knepley 659566063dSJacob Faibussowitsch PetscCall(DMForestTemplate(forest, PETSC_COMM_WORLD, &postforest)); 669566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityLabel(postforest, adaptLabel)); 679566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&adaptLabel)); 689566063dSJacob Faibussowitsch PetscCall(DMSetUp(postforest)); 699566063dSJacob Faibussowitsch PetscCall(DMDestroy(&forest)); 709b5c7ca5SMatthew G. Knepley forest = postforest; 719b5c7ca5SMatthew G. Knepley } 729b5c7ca5SMatthew G. Knepley } 739566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(forest, NULL, "-dm_view")); 749b5c7ca5SMatthew G. Knepley 759b5c7ca5SMatthew G. Knepley /* Setup the section*/ 769b5c7ca5SMatthew G. Knepley if (userSection) { 779566063dSJacob Faibussowitsch PetscCall(DMConvert(forest, DMPLEX, &plex)); 789566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(plex, NULL, "-dm_view")); 799566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd)); 809566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "Vertices [%D, %D)\n", vStart, vEnd)); 819566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 829566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) forest), &s)); 839566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(s, 1)); 849566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(s, vStart, vEnd)); 859b5c7ca5SMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 869566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(s, v, 1)); 879566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(s, v, 0, 1)); 889b5c7ca5SMatthew G. Knepley } 899566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(s)); 909566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(forest, s)); 919566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject) s, NULL, "-my_section_view")); 929566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&s)); 939566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex)); 949b5c7ca5SMatthew G. Knepley } else { 959b5c7ca5SMatthew G. Knepley PetscFE fe; 969b5c7ca5SMatthew G. Knepley PetscInt dim; 979b5c7ca5SMatthew G. Knepley 989566063dSJacob Faibussowitsch PetscCall(DMGetDimension(forest, &dim)); 999566063dSJacob Faibussowitsch PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, 1, PETSC_DETERMINE, &fe)); 1009566063dSJacob Faibussowitsch PetscCall(DMAddField(forest, NULL, (PetscObject) fe)); 1019566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 1029566063dSJacob Faibussowitsch PetscCall(DMCreateDS(forest)); 1039b5c7ca5SMatthew G. Knepley } 1049b5c7ca5SMatthew G. Knepley 1059b5c7ca5SMatthew G. Knepley /* Create the global vector*/ 1069566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(forest, &g)); 1079566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) g, "g")); 1089566063dSJacob Faibussowitsch PetscCall(VecSet(g, 1.0)); 1099b5c7ca5SMatthew G. Knepley 1109b5c7ca5SMatthew G. Knepley /* Test global to local*/ 1119b5c7ca5SMatthew G. Knepley Vec l; 1129566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(forest, &l)); 1139566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(l)); 1149566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocal(forest, g, INSERT_VALUES, l)); 1159566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(g)); 1169566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobal(forest, l, INSERT_VALUES, g)); 1179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&l)); 1189b5c7ca5SMatthew G. Knepley 1199b5c7ca5SMatthew G. Knepley /* Save a vector*/ 1209566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forestHDF.h5", FILE_MODE_WRITE, &viewer)); 1219566063dSJacob Faibussowitsch PetscCall(VecView(g, viewer)); 1229566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1239b5c7ca5SMatthew G. Knepley 1249b5c7ca5SMatthew G. Knepley /* Load another vector to load into*/ 1259566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(forest, &g2)); 1269566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) g2, "g")); 1279566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(g2)); 1289b5c7ca5SMatthew G. Knepley 1299b5c7ca5SMatthew G. Knepley /* Load a vector*/ 1309566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forestHDF.h5", FILE_MODE_READ, &viewer)); 1319566063dSJacob Faibussowitsch PetscCall(VecLoad(g2, viewer)); 1329566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1339b5c7ca5SMatthew G. Knepley 1349b5c7ca5SMatthew G. Knepley /* Check if the data is the same*/ 1359566063dSJacob Faibussowitsch PetscCall(VecAXPY(g2, -1.0, g)); 1369566063dSJacob Faibussowitsch PetscCall(VecNorm(g2, NORM_INFINITY, &nrm)); 137*08401ef6SPierre Jolivet PetscCheck(PetscAbsReal(nrm) <= PETSC_SMALL,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Invalid difference norm %g", (double) nrm); 1389b5c7ca5SMatthew G. Knepley 1399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&g)); 1409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&g2)); 1419566063dSJacob Faibussowitsch PetscCall(DMDestroy(&forest)); 1429566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 143b122ec5aSJacob 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