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