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 7d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 8d71ae5a4SJacob Faibussowitsch { 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 17327415f7SBarry Smith PetscFunctionBeginUser; 189566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-adapt", &adapt, NULL)); 209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-user_section", &userSection, NULL)); 219b5c7ca5SMatthew G. Knepley 229b5c7ca5SMatthew G. Knepley /* Create a base DMPlex mesh */ 239566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &base)); 249566063dSJacob Faibussowitsch PetscCall(DMSetType(base, DMPLEX)); 259566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(base)); 269566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(base, NULL, "-dm_view")); 279b5c7ca5SMatthew G. Knepley 28*d5b43468SJose E. Roman /* Convert Plex mesh to Forest and destroy base */ 299566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &forest)); 309566063dSJacob Faibussowitsch PetscCall(DMSetType(forest, DMP4EST)); 319566063dSJacob Faibussowitsch PetscCall(DMForestSetBaseDM(forest, base)); 329566063dSJacob Faibussowitsch PetscCall(DMSetUp(forest)); 339566063dSJacob Faibussowitsch PetscCall(DMDestroy(&base)); 349566063dSJacob Faibussowitsch PetscCall(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 429566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel)); 439566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(adaptLabel, 0, DM_ADAPT_REFINE)); 449566063dSJacob Faibussowitsch PetscCall(DMForestTemplate(forest, PETSC_COMM_WORLD, &postforest)); 459566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityLabel(postforest, adaptLabel)); 469566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&adaptLabel)); 479566063dSJacob Faibussowitsch PetscCall(DMSetUp(postforest)); 489566063dSJacob Faibussowitsch PetscCall(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 599566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel)); 609b5c7ca5SMatthew G. Knepley 619566063dSJacob Faibussowitsch PetscCall(DMForestGetCellChart(forest, &cStart, &cEnd)); 6248a46eb9SPierre Jolivet for (c = cStart; c < cEnd; ++c) PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE)); 639b5c7ca5SMatthew G. Knepley 649566063dSJacob Faibussowitsch PetscCall(DMForestTemplate(forest, PETSC_COMM_WORLD, &postforest)); 659566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityLabel(postforest, adaptLabel)); 669566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&adaptLabel)); 679566063dSJacob Faibussowitsch PetscCall(DMSetUp(postforest)); 689566063dSJacob Faibussowitsch PetscCall(DMDestroy(&forest)); 699b5c7ca5SMatthew G. Knepley forest = postforest; 709b5c7ca5SMatthew G. Knepley } 719b5c7ca5SMatthew G. Knepley } 729566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(forest, NULL, "-dm_view")); 739b5c7ca5SMatthew G. Knepley 749b5c7ca5SMatthew G. Knepley /* Setup the section*/ 759b5c7ca5SMatthew G. Knepley if (userSection) { 769566063dSJacob Faibussowitsch PetscCall(DMConvert(forest, DMPLEX, &plex)); 779566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(plex, NULL, "-dm_view")); 789566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd)); 7963a3b9bcSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "Vertices [%" PetscInt_FMT ", %" PetscInt_FMT ")\n", vStart, vEnd)); 809566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 819566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)forest), &s)); 829566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(s, 1)); 839566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(s, vStart, vEnd)); 849b5c7ca5SMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 859566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(s, v, 1)); 869566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(s, v, 0, 1)); 879b5c7ca5SMatthew G. Knepley } 889566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(s)); 899566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(forest, s)); 909566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)s, NULL, "-my_section_view")); 919566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&s)); 929566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex)); 939b5c7ca5SMatthew G. Knepley } else { 949b5c7ca5SMatthew G. Knepley PetscFE fe; 959b5c7ca5SMatthew G. Knepley PetscInt dim; 969b5c7ca5SMatthew G. Knepley 979566063dSJacob Faibussowitsch PetscCall(DMGetDimension(forest, &dim)); 989566063dSJacob Faibussowitsch PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, 1, PETSC_DETERMINE, &fe)); 999566063dSJacob Faibussowitsch PetscCall(DMAddField(forest, NULL, (PetscObject)fe)); 1009566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 1019566063dSJacob Faibussowitsch PetscCall(DMCreateDS(forest)); 1029b5c7ca5SMatthew G. Knepley } 1039b5c7ca5SMatthew G. Knepley 1049b5c7ca5SMatthew G. Knepley /* Create the global vector*/ 1059566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(forest, &g)); 1069566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)g, "g")); 1079566063dSJacob Faibussowitsch PetscCall(VecSet(g, 1.0)); 1089b5c7ca5SMatthew G. Knepley 1099b5c7ca5SMatthew G. Knepley /* Test global to local*/ 1109b5c7ca5SMatthew G. Knepley Vec l; 1119566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(forest, &l)); 1129566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(l)); 1139566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocal(forest, g, INSERT_VALUES, l)); 1149566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(g)); 1159566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobal(forest, l, INSERT_VALUES, g)); 1169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&l)); 1179b5c7ca5SMatthew G. Knepley 1189b5c7ca5SMatthew G. Knepley /* Save a vector*/ 1199566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forestHDF.h5", FILE_MODE_WRITE, &viewer)); 1209566063dSJacob Faibussowitsch PetscCall(VecView(g, viewer)); 1219566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1229b5c7ca5SMatthew G. Knepley 1239b5c7ca5SMatthew G. Knepley /* Load another vector to load into*/ 1249566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(forest, &g2)); 1259566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)g2, "g")); 1269566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(g2)); 1279b5c7ca5SMatthew G. Knepley 1289b5c7ca5SMatthew G. Knepley /* Load a vector*/ 1299566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forestHDF.h5", FILE_MODE_READ, &viewer)); 1309566063dSJacob Faibussowitsch PetscCall(VecLoad(g2, viewer)); 1319566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1329b5c7ca5SMatthew G. Knepley 1339b5c7ca5SMatthew G. Knepley /* Check if the data is the same*/ 1349566063dSJacob Faibussowitsch PetscCall(VecAXPY(g2, -1.0, g)); 1359566063dSJacob Faibussowitsch PetscCall(VecNorm(g2, NORM_INFINITY, &nrm)); 13608401ef6SPierre Jolivet PetscCheck(PetscAbsReal(nrm) <= PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Invalid difference norm %g", (double)nrm); 1379b5c7ca5SMatthew G. Knepley 1389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&g)); 1399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&g2)); 1409566063dSJacob Faibussowitsch PetscCall(DMDestroy(&forest)); 1419566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 142b122ec5aSJacob Faibussowitsch return 0; 1439b5c7ca5SMatthew G. Knepley } 1449b5c7ca5SMatthew G. Knepley 1459b5c7ca5SMatthew G. Knepley /*TEST 1469b5c7ca5SMatthew G. Knepley 1479b5c7ca5SMatthew G. Knepley build: 1489b5c7ca5SMatthew G. Knepley requires: hdf5 p4est 1499b5c7ca5SMatthew G. Knepley 1509b5c7ca5SMatthew G. Knepley test: 1519b5c7ca5SMatthew G. Knepley suffix: 0 1529b5c7ca5SMatthew G. Knepley nsize: {{1 2 5}} 15330602db0SMatthew G. Knepley args: -adapt -dm_plex_simplex 0 -dm_plex_box_faces 2,2 1549b5c7ca5SMatthew G. Knepley 1559b5c7ca5SMatthew G. Knepley TEST*/ 156