xref: /petsc/src/dm/impls/plex/tests/ex2.c (revision 12b8a6da77924c4c4db30ad685c4614c4146ffba)
1*12b8a6daSStefano Zampini static char help[] = "Benchmark DMPlexInterpolate.\n\n";
2*12b8a6daSStefano Zampini 
3*12b8a6daSStefano Zampini #include <petscdmplex.h>
4*12b8a6daSStefano Zampini 
5*12b8a6daSStefano Zampini int main(int argc, char **argv)
6*12b8a6daSStefano Zampini {
7*12b8a6daSStefano Zampini   DM dm, dm2;
8*12b8a6daSStefano Zampini #if defined(PETSC_USE_LOG)
9*12b8a6daSStefano Zampini   PetscLogStage stage;
10*12b8a6daSStefano Zampini #endif
11*12b8a6daSStefano Zampini 
12*12b8a6daSStefano Zampini   PetscFunctionBeginUser;
13*12b8a6daSStefano Zampini   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
14*12b8a6daSStefano Zampini   PetscCall(PetscLogStageRegister("Interpolate", &stage));
15*12b8a6daSStefano Zampini 
16*12b8a6daSStefano Zampini   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
17*12b8a6daSStefano Zampini   PetscCall(DMSetType(dm, DMPLEX));
18*12b8a6daSStefano Zampini   PetscCall(DMSetFromOptions(dm));
19*12b8a6daSStefano Zampini   PetscCall(DMViewFromOptions(dm, NULL, "-init_dm_view"));
20*12b8a6daSStefano Zampini 
21*12b8a6daSStefano Zampini   PetscCall(DMPlexUninterpolate(dm, &dm2));
22*12b8a6daSStefano Zampini   PetscCall(DMDestroy(&dm));
23*12b8a6daSStefano Zampini   dm = dm2;
24*12b8a6daSStefano Zampini   PetscCall(DMViewFromOptions(dm, NULL, "-unint_dm_view"));
25*12b8a6daSStefano Zampini 
26*12b8a6daSStefano Zampini   PetscCall(PetscLogStagePush(stage));
27*12b8a6daSStefano Zampini   PetscCall(DMPlexInterpolate(dm, &dm2));
28*12b8a6daSStefano Zampini   PetscCall(PetscLogStagePop());
29*12b8a6daSStefano Zampini 
30*12b8a6daSStefano Zampini   PetscCall(DMViewFromOptions(dm2, NULL, "-interp_dm_view"));
31*12b8a6daSStefano Zampini 
32*12b8a6daSStefano Zampini   PetscCall(DMDestroy(&dm2));
33*12b8a6daSStefano Zampini   PetscCall(DMDestroy(&dm));
34*12b8a6daSStefano Zampini   PetscCall(PetscFinalize());
35*12b8a6daSStefano Zampini   return 0;
36*12b8a6daSStefano Zampini }
37*12b8a6daSStefano Zampini 
38*12b8a6daSStefano Zampini /*TEST
39*12b8a6daSStefano Zampini 
40*12b8a6daSStefano Zampini   test:
41*12b8a6daSStefano Zampini     suffix: 0
42*12b8a6daSStefano Zampini     args: -dm_plex_simplex 0
43*12b8a6daSStefano Zampini 
44*12b8a6daSStefano Zampini TEST*/
45