1c4762a1bSJed Brown static char help[] = "Evaluate the shape quality of a mesh\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown typedef struct { 6c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */ 7c4762a1bSJed Brown PetscBool report; /* Print a quality report */ 8*f108dbd7SJacob Faibussowitsch PetscReal condLimit, tol; /* Condition number limit for cell output */ 9*f108dbd7SJacob Faibussowitsch PetscBool interpolate, distribute, simplex; 10*f108dbd7SJacob Faibussowitsch PetscInt dim, overlap; 11c4762a1bSJed Brown } AppCtx; 12c4762a1bSJed Brown 13*f108dbd7SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 14c4762a1bSJed Brown { 15c4762a1bSJed Brown PetscErrorCode ierr; 16c4762a1bSJed Brown 17*f108dbd7SJacob Faibussowitsch PetscFunctionBeginUser; 18c4762a1bSJed Brown options->filename[0] = '\0'; 19c4762a1bSJed Brown options->report = PETSC_FALSE; 20*f108dbd7SJacob Faibussowitsch options->interpolate = PETSC_FALSE; 21*f108dbd7SJacob Faibussowitsch options->distribute = PETSC_FALSE; 22*f108dbd7SJacob Faibussowitsch options->simplex = PETSC_FALSE; 23*f108dbd7SJacob Faibussowitsch options->overlap = 0; 24*f108dbd7SJacob Faibussowitsch options->dim = 2; 25*f108dbd7SJacob Faibussowitsch options->tol = 0.5; 26c4762a1bSJed Brown options->condLimit = PETSC_DETERMINE; 27c4762a1bSJed Brown 28c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Mesh Quality Evaluation Options", "DMPLEX");CHKERRQ(ierr); 29589a23caSBarry Smith ierr = PetscOptionsString("-filename", "The mesh file", "ex9.c", options->filename, options->filename, sizeof(options->filename), NULL);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = PetscOptionsBool("-report", "Output a mesh quality report", "ex9.c", options->report, &options->report, NULL);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = PetscOptionsReal("-cond_limit", "Condition number limit for cell output", "ex9.c", options->condLimit, &options->condLimit, NULL);CHKERRQ(ierr); 32*f108dbd7SJacob Faibussowitsch ierr = PetscOptionsReal("-orth_qual_atol", "Absolute tolerance for Orthogonal Quality", "ex9.c", options->tol, &options->tol, NULL);CHKERRQ(ierr); 33*f108dbd7SJacob Faibussowitsch ierr = PetscOptionsBool("-interpolate", "Interpolate mesh", "ex9.c", options->interpolate, &options->interpolate, NULL);CHKERRQ(ierr); 34*f108dbd7SJacob Faibussowitsch ierr = PetscOptionsBool("-distribute", "Distribute mesh", "ex9.c", options->distribute, &options->distribute, NULL);CHKERRQ(ierr); 35*f108dbd7SJacob Faibussowitsch ierr = PetscOptionsBool("-simplex", "Create simplex mesh elements", "ex9.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr); 36*f108dbd7SJacob Faibussowitsch ierr = PetscOptionsInt("-overlap", "Number of overlap levels for dm", "ex9.c", options->overlap, &options->overlap, NULL);CHKERRQ(ierr); 37*f108dbd7SJacob Faibussowitsch ierr = PetscOptionsInt("-dim", "Dimension of mesh if generated", "ex9.c", options->dim, &options->dim, NULL);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 39*f108dbd7SJacob Faibussowitsch PetscFunctionReturn(0); 40*f108dbd7SJacob Faibussowitsch } 41c4762a1bSJed Brown 42*f108dbd7SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 43*f108dbd7SJacob Faibussowitsch { 44*f108dbd7SJacob Faibussowitsch PetscErrorCode ierr; 45*f108dbd7SJacob Faibussowitsch 46*f108dbd7SJacob Faibussowitsch PetscFunctionBeginUser; 47*f108dbd7SJacob Faibussowitsch if (user->filename[0]) { 48*f108dbd7SJacob Faibussowitsch ierr = DMPlexCreateFromFile(comm, user->filename, PETSC_TRUE, dm);CHKERRQ(ierr); 49*f108dbd7SJacob Faibussowitsch } else { 50*f108dbd7SJacob Faibussowitsch ierr = DMPlexCreateBoxMesh(comm, user->dim, user->simplex, NULL, NULL, NULL, NULL, user->interpolate, dm);CHKERRQ(ierr); 51*f108dbd7SJacob Faibussowitsch ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr); 52*f108dbd7SJacob Faibussowitsch } 53*f108dbd7SJacob Faibussowitsch ierr = DMSetUp(*dm);CHKERRQ(ierr); 54*f108dbd7SJacob Faibussowitsch ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 55*f108dbd7SJacob Faibussowitsch if (user->distribute) { 56*f108dbd7SJacob Faibussowitsch DM dmDist = NULL; 57*f108dbd7SJacob Faibussowitsch ierr = DMPlexDistribute(*dm, user->overlap, NULL, &dmDist);CHKERRQ(ierr); 58*f108dbd7SJacob Faibussowitsch if (dmDist) { 59*f108dbd7SJacob Faibussowitsch ierr = DMDestroy(dm);CHKERRQ(ierr); 60*f108dbd7SJacob Faibussowitsch *dm = dmDist; 61*f108dbd7SJacob Faibussowitsch } 62*f108dbd7SJacob Faibussowitsch } 63*f108dbd7SJacob Faibussowitsch ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 64c4762a1bSJed Brown PetscFunctionReturn(0); 65c4762a1bSJed Brown } 66c4762a1bSJed Brown 67c4762a1bSJed Brown int main(int argc, char **argv) 68c4762a1bSJed Brown { 69c4762a1bSJed Brown DM dm; 70*f108dbd7SJacob Faibussowitsch DMLabel OQLabel; 71*f108dbd7SJacob Faibussowitsch Vec OQ; 72c4762a1bSJed Brown AppCtx ctx; 73c4762a1bSJed Brown PetscErrorCode ierr; 74c4762a1bSJed Brown 75c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 76c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr); 77*f108dbd7SJacob Faibussowitsch ierr = CreateMesh(PETSC_COMM_WORLD, &ctx, &dm);CHKERRQ(ierr); 78c4762a1bSJed Brown ierr = DMPlexCheckCellShape(dm, ctx.report, ctx.condLimit);CHKERRQ(ierr); 79*f108dbd7SJacob Faibussowitsch ierr = DMPlexComputeOrthogonalQuality(dm, NULL, ctx.tol, &OQ, &OQLabel);CHKERRQ(ierr); 80*f108dbd7SJacob Faibussowitsch ierr = VecDestroy(&OQ);CHKERRQ(ierr); 81c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 82c4762a1bSJed Brown ierr = PetscFinalize(); 83c4762a1bSJed Brown return ierr; 84c4762a1bSJed Brown } 85c4762a1bSJed Brown 86c4762a1bSJed Brown /*TEST 87c4762a1bSJed Brown 88c4762a1bSJed Brown test: 89c4762a1bSJed Brown suffix: 0 90c4762a1bSJed Brown requires: exodusii 91c4762a1bSJed Brown nsize: {{1 2}} 92c4762a1bSJed Brown args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -report 93c4762a1bSJed Brown 94c4762a1bSJed Brown test: 95c4762a1bSJed Brown suffix: 1 96c4762a1bSJed Brown args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.msh -report 97c4762a1bSJed Brown 98*f108dbd7SJacob Faibussowitsch testset: 99*f108dbd7SJacob Faibussowitsch requires: exodusii 100*f108dbd7SJacob Faibussowitsch args: -dm_plex_orthogonal_quality_label_view -dm_plex_orthogonal_quality_vec_view 101*f108dbd7SJacob Faibussowitsch 102*f108dbd7SJacob Faibussowitsch test: 103*f108dbd7SJacob Faibussowitsch suffix: box_1 104*f108dbd7SJacob Faibussowitsch nsize: 1 105*f108dbd7SJacob Faibussowitsch args: -interpolate -distribute -dim 2 -dm_plex_box_faces 2,2 -orth_qual_atol 1.0 106*f108dbd7SJacob Faibussowitsch 107*f108dbd7SJacob Faibussowitsch test: 108*f108dbd7SJacob Faibussowitsch suffix: box_2 109*f108dbd7SJacob Faibussowitsch nsize: 2 110*f108dbd7SJacob Faibussowitsch args: -interpolate -distribute -dim 2 -dm_plex_box_faces 2,2 -orth_qual_atol 1.0 111*f108dbd7SJacob Faibussowitsch 112*f108dbd7SJacob Faibussowitsch test: 113*f108dbd7SJacob Faibussowitsch suffix: mesh_1 114*f108dbd7SJacob Faibussowitsch nsize: 1 115*f108dbd7SJacob Faibussowitsch args: -interpolate -distribute -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -orth_qual_atol 0.95 116*f108dbd7SJacob Faibussowitsch 117*f108dbd7SJacob Faibussowitsch test: 118*f108dbd7SJacob Faibussowitsch suffix: mesh_2 119*f108dbd7SJacob Faibussowitsch nsize: 2 120*f108dbd7SJacob Faibussowitsch args: -interpolate -distribute -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -orth_qual_atol 0.95 121c4762a1bSJed Brown TEST*/ 122