xref: /petsc/src/dm/impls/plex/tests/ex53.c (revision 89f8043cca62f388f56c6240fe932552a256252a)
127e88307SToby Isaac const char help[] = "Test Berend's example";
227e88307SToby Isaac 
327e88307SToby Isaac #include <petscdmplex.h>
427e88307SToby Isaac #include <petscdmforest.h>
527e88307SToby Isaac 
627e88307SToby Isaac int main(int argc, char **argv)
727e88307SToby Isaac {
827e88307SToby Isaac   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
927e88307SToby Isaac   MPI_Comm       comm            = PETSC_COMM_WORLD;
1027e88307SToby Isaac   PetscInt       dim             = 3;
1127e88307SToby Isaac   PetscInt       cells_per_dir[] = {3, 3, 3};
1227e88307SToby Isaac   PetscReal      dir_min[]       = {0.0, 0.0, 0.0};
1327e88307SToby Isaac   PetscReal      dir_max[]       = {1.0, 1.0, 1.0};
1427e88307SToby Isaac   PetscInt       overlap         = 1;
1527e88307SToby Isaac   DMBoundaryType bcs[]           = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
1627e88307SToby Isaac   DM             forest, plex;
1727e88307SToby Isaac   Vec            CoordVec;
1827e88307SToby Isaac 
1927e88307SToby Isaac   PetscCall(DMCreate(comm, &forest));
2027e88307SToby Isaac   PetscCall(DMSetType(forest, DMP8EST));
2127e88307SToby Isaac   PetscCall(DMSetBasicAdjacency(forest, PETSC_TRUE, PETSC_TRUE));
2227e88307SToby Isaac   {
2327e88307SToby Isaac     DM dm_base, pdm, idm;
2427e88307SToby Isaac     PetscCall(DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE /* simplex */
2527e88307SToby Isaac                                   ,
2627e88307SToby Isaac                                   cells_per_dir, dir_min, dir_max, bcs, PETSC_TRUE /* interpolate */
2727e88307SToby Isaac                                   ,
2827e88307SToby Isaac                                   &dm_base));
2927e88307SToby Isaac     PetscCall(DMSetBasicAdjacency(dm_base, PETSC_TRUE, PETSC_TRUE));
3027e88307SToby Isaac     PetscCall(DMPlexDistribute(dm_base, overlap, NULL, &pdm));
3127e88307SToby Isaac     if (pdm) {
3227e88307SToby Isaac       PetscCall(DMDestroy(&dm_base));
3327e88307SToby Isaac       dm_base = pdm;
3427e88307SToby Isaac     }
3527e88307SToby Isaac     PetscCall(DMPlexInterpolate(dm_base, &idm));
3627e88307SToby Isaac     PetscCall(DMDestroy(&dm_base));
3727e88307SToby Isaac     dm_base = idm;
3827e88307SToby Isaac     PetscCall(DMLocalizeCoordinates(dm_base));
3927e88307SToby Isaac     PetscCall(DMPlexDistributeSetDefault(dm_base, PETSC_FALSE));
4027e88307SToby Isaac     PetscCall(DMSetFromOptions(dm_base));
4127e88307SToby Isaac     PetscCall(DMViewFromOptions(dm_base, NULL, "-dm_base_view"));
4227e88307SToby Isaac     PetscCall(DMCopyFields(dm_base, forest));
4327e88307SToby Isaac     PetscCall(DMForestSetBaseDM(forest, dm_base));
4427e88307SToby Isaac     PetscCall(DMDestroy(&dm_base));
4527e88307SToby Isaac   }
4627e88307SToby Isaac   PetscCall(DMForestSetPartitionOverlap(forest, 1));
4727e88307SToby Isaac   PetscCall(DMSetFromOptions(forest));
4827e88307SToby Isaac   PetscCall(DMSetUp(forest));
4927e88307SToby Isaac   PetscCall(DMViewFromOptions(forest, NULL, "-dm_forest_view"));
5027e88307SToby Isaac 
5127e88307SToby Isaac   PetscCall(DMConvert(forest, DMPLEX, &plex));
5227e88307SToby Isaac   PetscCall(DMGetCoordinates(plex, &CoordVec));
5327e88307SToby Isaac   PetscCall(DMViewFromOptions(forest, NULL, "-dm_plex_view"));
5427e88307SToby Isaac 
5527e88307SToby Isaac   PetscCall(DMDestroy(&plex));
5627e88307SToby Isaac   PetscCall(DMDestroy(&forest));
5727e88307SToby Isaac   PetscCall(PetscFinalize());
5827e88307SToby Isaac   return 0;
5927e88307SToby Isaac }
6027e88307SToby Isaac 
6127e88307SToby Isaac /*TEST
6227e88307SToby Isaac 
6327e88307SToby Isaac   test:
64*89f8043cSToby Isaac     requires: p4est
6527e88307SToby Isaac     suffix: 0
6627e88307SToby Isaac 
6727e88307SToby Isaac TEST*/
68