xref: /petsc/src/dm/impls/plex/tests/ex30.c (revision 05552d4bd0342c207cc7c527d55aae516dc232d5)
1 const char help[] = "Test memory allocation in DMPlex refinement.\n\n";
2 
3 #include <petsc.h>
4 
5 int main(int argc, char **argv)
6 {
7   PetscErrorCode ierr;
8   DM             dm, rdm;
9   PetscViewer    vwr;
10   PetscBool      flg;
11   char           datafile[PETSC_MAX_PATH_LEN];
12   MPI_Comm       comm;
13 
14   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
15   comm = PETSC_COMM_WORLD;
16   ierr = PetscViewerCreate(comm, &vwr);CHKERRQ(ierr);
17   ierr = PetscViewerSetType(vwr, PETSCVIEWERHDF5);CHKERRQ(ierr);
18   ierr = PetscViewerFileSetMode(vwr, FILE_MODE_READ);CHKERRQ(ierr);
19   ierr = PetscOptionsGetString(NULL, NULL, "-f", datafile, sizeof(datafile), &flg);CHKERRQ(ierr);
20   if (!flg) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Must provide meshfile");
21   ierr = PetscViewerFileSetName(vwr, datafile);CHKERRQ(ierr);
22   ierr = DMCreate(comm, &dm);CHKERRQ(ierr);
23   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
24   ierr = DMLoad(dm, vwr);CHKERRQ(ierr);
25   ierr = PetscViewerDestroy(&vwr);CHKERRQ(ierr);
26   ierr = PetscObjectSetName((PetscObject)dm, "BaryDM");CHKERRQ(ierr);
27   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
28   ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr);
29   ierr = DMRefine(dm, comm, &rdm);CHKERRQ(ierr);
30   ierr = DMDestroy(&dm);CHKERRQ(ierr);
31   ierr = PetscObjectSetName((PetscObject)rdm, "RefinedDM");CHKERRQ(ierr);
32   ierr = DMViewFromOptions(rdm, NULL, "-refined_dm_view");CHKERRQ(ierr);
33   ierr = DMDestroy(&rdm);CHKERRQ(ierr);
34   ierr = PetscFinalize();
35   return ierr;
36 }
37 
38 /*TEST
39 
40   test:
41     requires: hdf5 double !complex !define(PETSC_USE_64BIT_INDICES)
42     args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/barycentricallyrefinedcube.h5 -dm_view ascii::ASCII_INFO_DETAIL -refined_dm_view ascii::ASCII_INFO_DETAIL
43 
44 TEST*/
45