1c77c71ffSToby Isaac const char help[] = "Test clearing stale AMR data (example contributed by Berend van Wachem)"; 2c77c71ffSToby Isaac 3c77c71ffSToby Isaac #include <petscdmplex.h> 4c77c71ffSToby Isaac #include <petscdmforest.h> 5c77c71ffSToby Isaac 6c77c71ffSToby Isaac PetscErrorCode CloneDMWithNewSection(DM OriginalDM, DM *NewDM, PetscInt NFields) 7c77c71ffSToby Isaac { 8c77c71ffSToby Isaac PetscSection section; 9c77c71ffSToby Isaac PetscInt *NumComp, *NumDof; 104d86920dSPierre Jolivet 114d86920dSPierre Jolivet PetscFunctionBegin; 12c77c71ffSToby Isaac PetscCall(DMClone(OriginalDM, NewDM)); 13c77c71ffSToby Isaac PetscCall(DMPlexDistributeSetDefault(*NewDM, PETSC_FALSE)); 14c77c71ffSToby Isaac PetscCall(DMClearDS(*NewDM)); 15c77c71ffSToby Isaac PetscCall(PetscCalloc2(1, &NumComp, 3, &NumDof)); 16c77c71ffSToby Isaac NumComp[0] = 1; 17c77c71ffSToby Isaac NumDof[2] = NFields; 18c77c71ffSToby Isaac PetscCall(DMSetNumFields(*NewDM, 1)); 19c77c71ffSToby Isaac PetscCall(DMSetFromOptions(*NewDM)); 20c77c71ffSToby Isaac PetscCall(DMPlexCreateSection(*NewDM, NULL, NumComp, NumDof, 0, NULL, NULL, NULL, NULL, §ion)); 21c77c71ffSToby Isaac PetscCall(DMSetLocalSection(*NewDM, section)); 22c77c71ffSToby Isaac PetscCall(PetscFree2(NumComp, NumDof)); 23c77c71ffSToby Isaac PetscCall(PetscSectionDestroy(§ion)); 24c77c71ffSToby Isaac PetscFE fe; 25c77c71ffSToby Isaac PetscCall(PetscFECreateDefault(PETSC_COMM_WORLD, 2, 1, PETSC_FALSE, NULL, PETSC_DEFAULT, &fe)); 26c77c71ffSToby Isaac PetscCall(DMSetField(*NewDM, 0, NULL, (PetscObject)fe)); 27c77c71ffSToby Isaac PetscCall(PetscFEDestroy(&fe)); 28c77c71ffSToby Isaac PetscCall(DMCreateDS(*NewDM)); 29c77c71ffSToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 30c77c71ffSToby Isaac } 31c77c71ffSToby Isaac 32c77c71ffSToby Isaac int main(int argc, char **argv) 33c77c71ffSToby Isaac { 34c77c71ffSToby Isaac PetscInt dim = 2; 35c77c71ffSToby Isaac PetscInt cells_per_dir[] = {1, 1}; 36c77c71ffSToby Isaac PetscReal dir_min[] = {0.0, 0.0}; 37c77c71ffSToby Isaac PetscReal dir_max[] = {1.0, 1.0}; 38c77c71ffSToby Isaac DMBoundaryType bcs[] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE}; 394d86920dSPierre Jolivet DM forest, plex, NewDM; 40c77c71ffSToby Isaac Vec NewDMVecGlobal, NewDMVecLocal; 41c77c71ffSToby Isaac 424d86920dSPierre Jolivet PetscFunctionBeginUser; 434d86920dSPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 444d86920dSPierre Jolivet PetscCall(DMCreate(PETSC_COMM_WORLD, &forest)); 45c77c71ffSToby Isaac PetscCall(DMSetType(forest, DMP4EST)); 46c77c71ffSToby Isaac { 47c77c71ffSToby Isaac DM dm_base; 484d86920dSPierre Jolivet PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, /* simplex */ PETSC_FALSE, cells_per_dir, dir_min, dir_max, bcs, /* interpolate */ PETSC_TRUE, &dm_base)); 49c77c71ffSToby Isaac PetscCall(DMSetFromOptions(dm_base)); 50c77c71ffSToby Isaac PetscCall(DMViewFromOptions(dm_base, NULL, "-dm_base_view")); 51*bb4b53efSMatthew G. Knepley PetscCall(DMCopyFields(dm_base, PETSC_DETERMINE, PETSC_DETERMINE, forest)); 52c77c71ffSToby Isaac PetscCall(DMForestSetBaseDM(forest, dm_base)); 53c77c71ffSToby Isaac PetscCall(DMDestroy(&dm_base)); 54c77c71ffSToby Isaac } 55c77c71ffSToby Isaac PetscCall(DMSetFromOptions(forest)); 56c77c71ffSToby Isaac PetscCall(DMSetUp(forest)); 57c77c71ffSToby Isaac 58c77c71ffSToby Isaac PetscCall(DMViewFromOptions(forest, NULL, "-dm_forest_view")); 59c77c71ffSToby Isaac 60c77c71ffSToby Isaac PetscCall(DMConvert(forest, DMPLEX, &plex)); 61c77c71ffSToby Isaac 62c77c71ffSToby Isaac PetscInt numFields = 2; 63c77c71ffSToby Isaac PetscInt numComp[2] = {1, 1}; 64c77c71ffSToby Isaac PetscInt numDof[6] = {0}; 65c77c71ffSToby Isaac for (PetscInt i = 0; i < numFields; i++) numDof[i * (dim + 1) + dim] = 1; 66c77c71ffSToby Isaac 67c77c71ffSToby Isaac PetscCall(DMSetNumFields(plex, numFields)); 68c77c71ffSToby Isaac PetscCall(DMSetNumFields(forest, numFields)); 69c77c71ffSToby Isaac 70c77c71ffSToby Isaac PetscSection section; 71c77c71ffSToby Isaac PetscCall(DMPlexCreateSection(plex, NULL, numComp, numDof, 0, NULL, NULL, NULL, NULL, §ion)); 72c77c71ffSToby Isaac 73c77c71ffSToby Isaac const char *names[] = {"field 0", "field 1"}; 74c77c71ffSToby Isaac for (PetscInt i = 0; i < numFields; i++) PetscCall(PetscSectionSetFieldName(section, i, names[i])); 75c77c71ffSToby Isaac PetscCall(DMSetLocalSection(forest, section)); 76c77c71ffSToby Isaac PetscCall(DMSetLocalSection(plex, section)); 77c77c71ffSToby Isaac PetscCall(PetscSectionDestroy(§ion)); 78c77c71ffSToby Isaac 79c77c71ffSToby Isaac PetscFE fe; 804d86920dSPierre Jolivet PetscCall(PetscFECreateDefault(PETSC_COMM_WORLD, dim, 1, PETSC_FALSE, NULL, PETSC_DEFAULT, &fe)); 81c77c71ffSToby Isaac for (PetscInt i = 0; i < numFields; i++) { 82c77c71ffSToby Isaac PetscCall(DMSetField(plex, i, NULL, (PetscObject)fe)); 83c77c71ffSToby Isaac PetscCall(DMSetField(forest, i, NULL, (PetscObject)fe)); 84c77c71ffSToby Isaac } 85c77c71ffSToby Isaac PetscCall(PetscFEDestroy(&fe)); 86c77c71ffSToby Isaac 87c77c71ffSToby Isaac PetscCall(DMCreateDS(plex)); 88c77c71ffSToby Isaac PetscCall(DMCreateDS(forest)); 89c77c71ffSToby Isaac 90c77c71ffSToby Isaac /* Make another DM, based on the layout of the previous DM, but with a different number of fields */ 91c77c71ffSToby Isaac PetscCall(CloneDMWithNewSection(plex, &NewDM, 1)); 92c77c71ffSToby Isaac PetscCall(DMCreateDS(NewDM)); 93c77c71ffSToby Isaac PetscCall(DMCreateGlobalVector(NewDM, &NewDMVecGlobal)); 94c77c71ffSToby Isaac PetscCall(DMCreateLocalVector(NewDM, &NewDMVecLocal)); 95c77c71ffSToby Isaac PetscCall(VecSet(NewDMVecGlobal, 3.141592)); 96c77c71ffSToby Isaac PetscCall(DMGlobalToLocal(NewDM, NewDMVecGlobal, INSERT_VALUES, NewDMVecLocal)); 97c77c71ffSToby Isaac PetscCall(VecDestroy(&NewDMVecLocal)); 98c77c71ffSToby Isaac PetscCall(VecDestroy(&NewDMVecGlobal)); 99c77c71ffSToby Isaac PetscCall(DMClearDS(NewDM)); 100c77c71ffSToby Isaac PetscCall(DMDestroy(&NewDM)); 101c77c71ffSToby Isaac 102c77c71ffSToby Isaac Vec g_vec, l_vec; 103c77c71ffSToby Isaac PetscCall(DMCreateGlobalVector(plex, &g_vec)); 104c77c71ffSToby Isaac PetscCall(VecSet(g_vec, 1.0)); 105c77c71ffSToby Isaac PetscCall(DMCreateLocalVector(plex, &l_vec)); 106c77c71ffSToby Isaac PetscCall(DMGlobalToLocal(plex, g_vec, INSERT_VALUES, l_vec)); 107c77c71ffSToby Isaac PetscCall(VecViewFromOptions(l_vec, NULL, "-local_vec_view")); 108c77c71ffSToby Isaac PetscCall(VecDestroy(&l_vec)); 109c77c71ffSToby Isaac PetscCall(VecDestroy(&g_vec)); 110c77c71ffSToby Isaac 111c77c71ffSToby Isaac PetscCall(DMViewFromOptions(forest, NULL, "-dm_plex_view")); 112c77c71ffSToby Isaac PetscCall(DMDestroy(&plex)); 113c77c71ffSToby Isaac PetscCall(DMDestroy(&forest)); 114c77c71ffSToby Isaac PetscCall(PetscFinalize()); 115c77c71ffSToby Isaac return 0; 116c77c71ffSToby Isaac } 117c77c71ffSToby Isaac 118c77c71ffSToby Isaac /*TEST 119c77c71ffSToby Isaac 120c77c71ffSToby Isaac test: 121c77c71ffSToby Isaac suffix: 0 122c77c71ffSToby Isaac requires: p4est 123c77c71ffSToby Isaac args: -dm_forest_initial_refinement 1 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash 124c77c71ffSToby Isaac 125c77c71ffSToby Isaac TEST*/ 126