Actual source code: ex9.c

petsc-3.14.6 2021-03-30
Report Typos and Errors
  1: static char help[] = "Evaluate the shape quality of a mesh\n\n";

  3: #include <petscdmplex.h>

  5: typedef struct {
  6:   char      filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */
  7:   PetscBool report;                       /* Print a quality report */
  8:   PetscReal condLimit, tol;               /* Condition number limit for cell output */
  9: } AppCtx;

 11: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
 12: {

 16:   options->filename[0] = '\0';
 17:   options->report      = PETSC_FALSE;
 18:   options->tol         = 0.5;
 19:   options->condLimit   = PETSC_DETERMINE;

 21:   PetscOptionsBegin(comm, "", "Mesh Quality Evaluation Options", "DMPLEX");
 22:   PetscOptionsString("-filename", "The mesh file", "ex9.c", options->filename, options->filename, sizeof(options->filename), NULL);
 23:   PetscOptionsBool("-report", "Output a mesh quality report", "ex9.c", options->report, &options->report, NULL);
 24:   PetscOptionsReal("-cond_limit", "Condition number limit for cell output", "ex9.c", options->condLimit, &options->condLimit, NULL);
 25:   PetscOptionsReal("-orth_qual_atol", "Absolute tolerance for Orthogonal Quality", "ex9.c", options->tol, &options->tol, NULL);
 26:   PetscOptionsEnd();
 27:   return(0);
 28: }

 30: static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
 31: {
 32:   PetscErrorCode  ierr;

 35:   if (user->filename[0]) {DMPlexCreateFromFile(comm, user->filename, PETSC_TRUE, dm);}
 36:   else                   {DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);}
 37:   PetscObjectSetName((PetscObject) *dm, "Mesh");
 38:   DMSetFromOptions(*dm);
 39:   DMViewFromOptions(*dm, NULL, "-dm_view");
 40:   return(0);
 41: }

 43: int main(int argc, char **argv)
 44: {
 45:   DM             dm;
 46:   DMLabel        OQLabel;
 47:   Vec            OQ;
 48:   AppCtx         ctx;

 51:   PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
 52:   ProcessOptions(PETSC_COMM_WORLD, &ctx);
 53:   CreateMesh(PETSC_COMM_WORLD, &ctx, &dm);
 54:   DMPlexCheckCellShape(dm, ctx.report, ctx.condLimit);
 55:   DMPlexComputeOrthogonalQuality(dm, NULL, ctx.tol, &OQ, &OQLabel);
 56:   VecDestroy(&OQ);
 57:   DMDestroy(&dm);
 58:   PetscFinalize();
 59:   return ierr;
 60: }

 62: /*TEST

 64:   test:
 65:     suffix: 0
 66:     requires: exodusii
 67:     nsize: {{1 2}}
 68:     args: -dm_distribute -petscpartitioner_type simple -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -report

 70:   test:
 71:     suffix: 1
 72:     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.msh -report

 74:   testset:
 75:     requires: exodusii
 76:     args: -dm_plex_orthogonal_quality_label_view -dm_plex_orthogonal_quality_vec_view

 78:     test:
 79:       suffix: box_1
 80:       nsize: 1
 81:       args: -dm_plex_box_simplex 0 -dm_plex_box_faces 2,2 -orth_qual_atol 1.0

 83:     test:
 84:       suffix: box_2
 85:       nsize: 2
 86:       args: -dm_distribute -petscpartitioner_type simple -dm_plex_box_simplex 0 -dm_plex_box_faces 2,2 -orth_qual_atol 1.0

 88:     test:
 89:       suffix: mesh_1
 90:       nsize: 1
 91:       args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -orth_qual_atol 0.95

 93:     test:
 94:       suffix: mesh_2
 95:       nsize: 2
 96:       args: -dm_distribute -petscpartitioner_type simple -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -orth_qual_atol 0.95
 97: TEST*/