xref: /petsc/src/dm/impls/plex/tests/ex30.c (revision dfd57a172ac9fa6c7b5fe6de6ab5df85cefc2996)
1c4762a1bSJed Brown const char help[] = "Test memory allocation in DMPlex refinement.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petsc.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown int main(int argc, char **argv)
6c4762a1bSJed Brown {
7b5a892a1SMatthew G. Knepley   DM             dm;
8c4762a1bSJed Brown   PetscErrorCode ierr;
9c4762a1bSJed Brown 
102da392ccSBarry Smith   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
11b5a892a1SMatthew G. Knepley   ierr = DMCreate(PETSC_COMM_WORLD, &dm);CHKERRQ(ierr);
12c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) dm, "BaryDM");CHKERRQ(ierr);
13b5a892a1SMatthew G. Knepley   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
14b5a892a1SMatthew G. Knepley   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
15c4762a1bSJed Brown   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
16b5a892a1SMatthew G. Knepley   //ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr);
17b5a892a1SMatthew G. Knepley   //ierr = DMRefine(dm, comm, &rdm);CHKERRQ(ierr);
18b5a892a1SMatthew G. Knepley   //ierr = DMPlexConvertOldOrientations_Internal(dm);CHKERRQ(ierr);
19b5a892a1SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) dm, "RefinedDM");CHKERRQ(ierr);
20b5a892a1SMatthew G. Knepley   ierr = PetscObjectSetOptionsPrefix((PetscObject) dm, "ref_");CHKERRQ(ierr);
21b5a892a1SMatthew G. Knepley   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
22b5a892a1SMatthew G. Knepley   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
23c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
24c4762a1bSJed Brown   ierr = PetscFinalize();
25c4762a1bSJed Brown   return ierr;
26c4762a1bSJed Brown }
27c4762a1bSJed Brown 
28c4762a1bSJed Brown /*TEST
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   test:
31*dfd57a17SPierre Jolivet     requires: hdf5 double !complex !defined(PETSC_USE_64BIT_INDICES)
32b5a892a1SMatthew G. Knepley     args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/barycentricallyrefinedcube.h5 -dm_view ascii::ASCII_INFO_DETAIL -ref_dm_refine 1 -ref_dm_view ascii::ASCII_INFO_DETAIL
33c4762a1bSJed Brown 
34c4762a1bSJed Brown TEST*/
35