xref: /petsc/src/dm/impls/plex/tests/ex53.c (revision 27e88307fb5906948719cedfd0799d3aeca9c5c7)
1*27e88307SToby Isaac const char help[] = "Test Berend's example";
2*27e88307SToby Isaac 
3*27e88307SToby Isaac #include <petscdmplex.h>
4*27e88307SToby Isaac #include <petscdmforest.h>
5*27e88307SToby Isaac 
6*27e88307SToby Isaac int main(int argc, char **argv)
7*27e88307SToby Isaac {
8*27e88307SToby Isaac   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
9*27e88307SToby Isaac   MPI_Comm       comm            = PETSC_COMM_WORLD;
10*27e88307SToby Isaac   PetscInt       dim             = 3;
11*27e88307SToby Isaac   PetscInt       cells_per_dir[] = {3, 3, 3};
12*27e88307SToby Isaac   PetscReal      dir_min[]       = {0.0, 0.0, 0.0};
13*27e88307SToby Isaac   PetscReal      dir_max[]       = {1.0, 1.0, 1.0};
14*27e88307SToby Isaac   PetscInt       overlap         = 1;
15*27e88307SToby Isaac   DMBoundaryType bcs[]           = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
16*27e88307SToby Isaac   DM             forest, plex;
17*27e88307SToby Isaac   Vec            CoordVec;
18*27e88307SToby Isaac 
19*27e88307SToby Isaac   PetscCall(DMCreate(comm, &forest));
20*27e88307SToby Isaac   PetscCall(DMSetType(forest, DMP8EST));
21*27e88307SToby Isaac   PetscCall(DMSetBasicAdjacency(forest, PETSC_TRUE, PETSC_TRUE));
22*27e88307SToby Isaac   {
23*27e88307SToby Isaac     DM dm_base, pdm, idm;
24*27e88307SToby Isaac     PetscCall(DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE /* simplex */
25*27e88307SToby Isaac                                   ,
26*27e88307SToby Isaac                                   cells_per_dir, dir_min, dir_max, bcs, PETSC_TRUE /* interpolate */
27*27e88307SToby Isaac                                   ,
28*27e88307SToby Isaac                                   &dm_base));
29*27e88307SToby Isaac     PetscCall(DMSetBasicAdjacency(dm_base, PETSC_TRUE, PETSC_TRUE));
30*27e88307SToby Isaac     PetscCall(DMPlexDistribute(dm_base, overlap, NULL, &pdm));
31*27e88307SToby Isaac     if (pdm) {
32*27e88307SToby Isaac       PetscCall(DMDestroy(&dm_base));
33*27e88307SToby Isaac       dm_base = pdm;
34*27e88307SToby Isaac     }
35*27e88307SToby Isaac     PetscCall(DMPlexInterpolate(dm_base, &idm));
36*27e88307SToby Isaac     PetscCall(DMDestroy(&dm_base));
37*27e88307SToby Isaac     dm_base = idm;
38*27e88307SToby Isaac     PetscCall(DMLocalizeCoordinates(dm_base));
39*27e88307SToby Isaac     PetscCall(DMPlexDistributeSetDefault(dm_base, PETSC_FALSE));
40*27e88307SToby Isaac     PetscCall(DMSetFromOptions(dm_base));
41*27e88307SToby Isaac     PetscCall(DMViewFromOptions(dm_base, NULL, "-dm_base_view"));
42*27e88307SToby Isaac     PetscCall(DMCopyFields(dm_base, forest));
43*27e88307SToby Isaac     PetscCall(DMForestSetBaseDM(forest, dm_base));
44*27e88307SToby Isaac     PetscCall(DMDestroy(&dm_base));
45*27e88307SToby Isaac   }
46*27e88307SToby Isaac   PetscCall(DMForestSetPartitionOverlap(forest, 1));
47*27e88307SToby Isaac   PetscCall(DMSetFromOptions(forest));
48*27e88307SToby Isaac   PetscCall(DMSetUp(forest));
49*27e88307SToby Isaac   PetscCall(DMViewFromOptions(forest, NULL, "-dm_forest_view"));
50*27e88307SToby Isaac 
51*27e88307SToby Isaac   PetscCall(DMConvert(forest, DMPLEX, &plex));
52*27e88307SToby Isaac   PetscCall(DMGetCoordinates(plex, &CoordVec));
53*27e88307SToby Isaac   PetscCall(DMViewFromOptions(forest, NULL, "-dm_plex_view"));
54*27e88307SToby Isaac 
55*27e88307SToby Isaac   PetscCall(DMDestroy(&plex));
56*27e88307SToby Isaac   PetscCall(DMDestroy(&forest));
57*27e88307SToby Isaac   PetscCall(PetscFinalize());
58*27e88307SToby Isaac   return 0;
59*27e88307SToby Isaac }
60*27e88307SToby Isaac 
61*27e88307SToby Isaac /*TEST
62*27e88307SToby Isaac 
63*27e88307SToby Isaac   test:
64*27e88307SToby Isaac     suffix: 0
65*27e88307SToby Isaac 
66*27e88307SToby Isaac TEST*/
67