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); 15 if (ierr) return ierr; 16 comm = PETSC_COMM_WORLD; 17 ierr = PetscViewerCreate(comm, &vwr);CHKERRQ(ierr); 18 ierr = PetscViewerSetType(vwr, PETSCVIEWERHDF5);CHKERRQ(ierr); 19 ierr = PetscViewerFileSetMode(vwr, FILE_MODE_READ);CHKERRQ(ierr); 20 ierr = PetscOptionsGetString(NULL, NULL, "-f", datafile, PETSC_MAX_PATH_LEN, &flg);CHKERRQ(ierr); 21 if (!flg) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Must provide meshfile"); 22 ierr = PetscViewerFileSetName(vwr, datafile);CHKERRQ(ierr); 23 ierr = DMCreate(comm, &dm);CHKERRQ(ierr); 24 ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr); 25 ierr = DMLoad(dm, vwr);CHKERRQ(ierr); 26 ierr = PetscViewerDestroy(&vwr);CHKERRQ(ierr); 27 ierr = PetscObjectSetName((PetscObject)dm, "BaryDM");CHKERRQ(ierr); 28 ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr); 29 ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr); 30 ierr = DMRefine(dm, comm, &rdm);CHKERRQ(ierr); 31 ierr = DMDestroy(&dm);CHKERRQ(ierr); 32 ierr = PetscObjectSetName((PetscObject)rdm, "RefinedDM");CHKERRQ(ierr); 33 ierr = DMViewFromOptions(rdm, NULL, "-refined_dm_view");CHKERRQ(ierr); 34 ierr = DMDestroy(&rdm);CHKERRQ(ierr); 35 ierr = PetscFinalize(); 36 return ierr; 37 } 38 39 /*TEST 40 41 test: 42 requires: hdf5 double !complex !define(PETSC_USE_64BIT_INDICES) 43 args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/barycentricallyrefinedcube.h5 -dm_view ascii::ASCII_INFO_DETAIL -refined_dm_view ascii::ASCII_INFO_DETAIL 44 45 TEST*/ 46