1*c4762a1bSJed Brown static char help[] = "Tests for coarsening\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscdmplex.h> 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown typedef struct { 6*c4762a1bSJed Brown PetscInt debug; /* The debugging level */ 7*c4762a1bSJed Brown PetscInt dim; /* The topological mesh dimension */ 8*c4762a1bSJed Brown PetscInt testNum; /* The particular mesh to test */ 9*c4762a1bSJed Brown PetscBool uninterpolate; /* Uninterpolate the mesh at the end */ 10*c4762a1bSJed Brown char filename[2048]; /* The optional mesh file */ 11*c4762a1bSJed Brown } AppCtx; 12*c4762a1bSJed Brown 13*c4762a1bSJed Brown PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 14*c4762a1bSJed Brown { 15*c4762a1bSJed Brown PetscErrorCode ierr; 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown PetscFunctionBegin; 18*c4762a1bSJed Brown options->debug = 0; 19*c4762a1bSJed Brown options->dim = 2; 20*c4762a1bSJed Brown options->testNum = 0; 21*c4762a1bSJed Brown options->uninterpolate = PETSC_FALSE; 22*c4762a1bSJed Brown options->filename[0] = '\0'; 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr); 25*c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-debug", "The debugging level", "ex14.c", options->debug, &options->debug, NULL,0);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex14.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-test_num", "The particular mesh to test", "ex14.c", options->testNum, &options->testNum, NULL,0);CHKERRQ(ierr); 28*c4762a1bSJed Brown ierr = PetscOptionsBool("-uninterpolate", "Uninterpolate the mesh at the end", "ex14.c", options->uninterpolate, &options->uninterpolate, NULL);CHKERRQ(ierr); 29*c4762a1bSJed Brown ierr = PetscOptionsString("-f", "Filename to read", "ex14.c", options->filename, options->filename, sizeof(options->filename), NULL);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = PetscOptionsEnd(); 31*c4762a1bSJed Brown PetscFunctionReturn(0); 32*c4762a1bSJed Brown } 33*c4762a1bSJed Brown 34*c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 35*c4762a1bSJed Brown { 36*c4762a1bSJed Brown const char *filename = user->filename; 37*c4762a1bSJed Brown PetscInt dim = user->dim; 38*c4762a1bSJed Brown size_t len; 39*c4762a1bSJed Brown PetscErrorCode ierr; 40*c4762a1bSJed Brown 41*c4762a1bSJed Brown PetscFunctionBegin; 42*c4762a1bSJed Brown ierr = PetscStrlen(filename, &len);CHKERRQ(ierr); 43*c4762a1bSJed Brown if (!len) {ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);} 44*c4762a1bSJed Brown else {ierr = DMPlexCreateFromFile(comm, filename, PETSC_TRUE, dm);CHKERRQ(ierr);} 45*c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-orig_dm_view");CHKERRQ(ierr); 46*c4762a1bSJed Brown { 47*c4762a1bSJed Brown DM distributedMesh = NULL; 48*c4762a1bSJed Brown 49*c4762a1bSJed Brown ierr = DMPlexDistribute(*dm, 0, NULL, &distributedMesh);CHKERRQ(ierr); 50*c4762a1bSJed Brown if (distributedMesh) { 51*c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 52*c4762a1bSJed Brown *dm = distributedMesh; 53*c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dist_dm_view");CHKERRQ(ierr); 54*c4762a1bSJed Brown } 55*c4762a1bSJed Brown } 56*c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 57*c4762a1bSJed Brown if (user->uninterpolate) { 58*c4762a1bSJed Brown DM udm = NULL; 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown ierr = DMPlexUninterpolate(*dm, &udm);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 62*c4762a1bSJed Brown *dm = udm; 63*c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-un_dm_view");CHKERRQ(ierr); 64*c4762a1bSJed Brown } 65*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) *dm, "Test Mesh");CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 67*c4762a1bSJed Brown PetscFunctionReturn(0); 68*c4762a1bSJed Brown } 69*c4762a1bSJed Brown 70*c4762a1bSJed Brown int main(int argc, char **argv) 71*c4762a1bSJed Brown { 72*c4762a1bSJed Brown DM dm; 73*c4762a1bSJed Brown AppCtx user; 74*c4762a1bSJed Brown PetscErrorCode ierr; 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 77*c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = PetscFinalize(); 81*c4762a1bSJed Brown return ierr; 82*c4762a1bSJed Brown } 83*c4762a1bSJed Brown 84*c4762a1bSJed Brown /*TEST 85*c4762a1bSJed Brown test: 86*c4762a1bSJed Brown suffix: 0 87*c4762a1bSJed Brown requires: triangle 88*c4762a1bSJed Brown args: -dm_view -dm_refine 1 -dm_coarsen -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_faces 89*c4762a1bSJed Brown 90*c4762a1bSJed Brown TEST*/ 91