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