xref: /petsc/src/dm/impls/plex/tests/ex53.c (revision 4bbe91384854fec5c16b677958fdcd4ec4bdaf6b)
121399927SMatthew G. Knepley const char help[] = "Test distribution with overlap using DMForest";
227e88307SToby Isaac 
327e88307SToby Isaac #include <petscdmplex.h>
427e88307SToby Isaac #include <petscdmforest.h>
527e88307SToby Isaac 
main(int argc,char ** argv)627e88307SToby Isaac int main(int argc, char **argv)
727e88307SToby Isaac {
827e88307SToby Isaac   DM       forest, plex;
927e88307SToby Isaac   Vec      CoordVec;
1021399927SMatthew G. Knepley   MPI_Comm comm;
1127e88307SToby Isaac 
1221399927SMatthew G. Knepley   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1321399927SMatthew G. Knepley   comm = PETSC_COMM_WORLD;
1427e88307SToby Isaac   PetscCall(DMCreate(comm, &forest));
1521399927SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)forest, "forest"));
1627e88307SToby Isaac   PetscCall(DMSetType(forest, DMP8EST));
1727e88307SToby Isaac   PetscCall(DMSetBasicAdjacency(forest, PETSC_TRUE, PETSC_TRUE));
1827e88307SToby Isaac   {
1921399927SMatthew G. Knepley     DM dm_base;
2021399927SMatthew G. Knepley 
2121399927SMatthew G. Knepley     PetscCall(DMCreate(comm, &dm_base));
2221399927SMatthew G. Knepley     PetscCall(DMSetType(dm_base, DMPLEX));
2321399927SMatthew G. Knepley     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm_base, "base_"));
2427e88307SToby Isaac     PetscCall(DMSetFromOptions(dm_base));
2521399927SMatthew G. Knepley     PetscCall(DMViewFromOptions(dm_base, NULL, "-dm_view"));
26*bb4b53efSMatthew G. Knepley     PetscCall(DMCopyFields(dm_base, PETSC_DETERMINE, PETSC_DETERMINE, forest));
2727e88307SToby Isaac     PetscCall(DMForestSetBaseDM(forest, dm_base));
2827e88307SToby Isaac     PetscCall(DMDestroy(&dm_base));
2927e88307SToby Isaac   }
3027e88307SToby Isaac   PetscCall(DMForestSetPartitionOverlap(forest, 1));
3127e88307SToby Isaac   PetscCall(DMSetFromOptions(forest));
3227e88307SToby Isaac   PetscCall(DMSetUp(forest));
3327e88307SToby Isaac   PetscCall(DMViewFromOptions(forest, NULL, "-dm_forest_view"));
3427e88307SToby Isaac 
3527e88307SToby Isaac   PetscCall(DMConvert(forest, DMPLEX, &plex));
3627e88307SToby Isaac   PetscCall(DMGetCoordinates(plex, &CoordVec));
3727e88307SToby Isaac   PetscCall(DMViewFromOptions(forest, NULL, "-dm_plex_view"));
3827e88307SToby Isaac 
3927e88307SToby Isaac   PetscCall(DMDestroy(&plex));
4027e88307SToby Isaac   PetscCall(DMDestroy(&forest));
4127e88307SToby Isaac   PetscCall(PetscFinalize());
4227e88307SToby Isaac   return 0;
4327e88307SToby Isaac }
4427e88307SToby Isaac 
4527e88307SToby Isaac /*TEST
4627e88307SToby Isaac 
4721399927SMatthew G. Knepley   testset:
4889f8043cSToby Isaac     requires: p4est
4921399927SMatthew G. Knepley     args: -base_dm_plex_dim 3 -base_dm_plex_simplex 0 -base_dm_plex_box_faces 3,3,3 -base_dm_distribute_overlap 1 \
5021399927SMatthew G. Knepley           -base_dm_plex_adj_cone true -base_dm_plex_adj_closure true \
5121399927SMatthew G. Knepley           -base_dm_view -dm_forest_view -dm_plex_view
5221399927SMatthew G. Knepley 
5321399927SMatthew G. Knepley     test:
5427e88307SToby Isaac       suffix: 0
5527e88307SToby Isaac 
5621399927SMatthew G. Knepley     test:
5721399927SMatthew G. Knepley       suffix: 1
5821399927SMatthew G. Knepley       nsize: 3
5921399927SMatthew G. Knepley       args: -petscpartitioner_type simple
6021399927SMatthew G. Knepley 
6127e88307SToby Isaac TEST*/
62