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