Actual source code: ex1.c

petsc-3.11.4 2019-09-28
Report Typos and Errors
  1: static char help[] = "Tests DMLabel operations.\n\n";

  3:  #include <petscdm.h>
  4:  #include <petscdmplex.h>

  6: PetscErrorCode ViewLabels(DM dm, PetscViewer viewer)
  7: {
  8:   DMLabel        label;
  9:   IS             labelIS;
 10:   const char    *labelName;
 11:   PetscInt       numLabels, l;

 15:   /* query the number and name of labels*/
 16:   DMGetNumLabels(dm, &numLabels);
 17:   PetscViewerASCIIPrintf(viewer, "Number of labels: %d\n", numLabels);
 18:   for (l = 0; l < numLabels; ++l) {
 19:     DMGetLabelName(dm, l, &labelName);
 20:     PetscViewerASCIIPrintf(viewer, "Label %d: name: %s\n", l, labelName);
 21:     PetscViewerASCIIPrintf(viewer, "IS of values\n");
 22:     DMGetLabel(dm, labelName, &label);
 23:     DMLabelGetValueIS(label, &labelIS);
 24:     PetscViewerASCIIPushTab(viewer);
 25:     ISView(labelIS, viewer);
 26:     PetscViewerASCIIPopTab(viewer);
 27:     ISDestroy(&labelIS);
 28:     PetscViewerASCIIPrintf(viewer, "\n");
 29:   }
 30:   /* Making sure that string literals work */
 31:   PetscViewerASCIIPrintf(viewer,"\n\nCell Set label IS\n");
 32:   DMGetLabel(dm, "Cell Sets", &label);
 33:   if (label) {
 34:     DMLabelGetValueIS(label, &labelIS);
 35:     ISView(labelIS, viewer);
 36:     ISDestroy(&labelIS);
 37:   }
 38:   return(0);
 39: }

 41: int main(int argc, char **argv)
 42: {
 43:   DM             dm, dmDist;
 44:   char           filename[PETSC_MAX_PATH_LEN];
 45:   PetscBool      interpolate = PETSC_FALSE;

 48:   /* initialize and get options */
 49:   PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
 50:   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "DMLabel ex1 Options", "DMLabel");
 51:   PetscOptionsString("-i", "filename to read", "ex1.c", filename, filename, sizeof(filename), NULL);
 52:   PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex1.c", interpolate, &interpolate, NULL);
 53:   PetscOptionsEnd();

 55:   /* create and distribute DM */
 56:   DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, interpolate, &dm);
 57:   DMPlexDistribute(dm, 0, NULL, &dmDist);
 58:   if (dmDist) {
 59:     DMDestroy(&dm);
 60:     dm   = dmDist;
 61:   }
 62:   DMSetFromOptions(dm);
 63:   ViewLabels(dm, PETSC_VIEWER_STDOUT_WORLD);
 64:   DMDestroy(&dm);
 65:   PetscFinalize();
 66:   return ierr;
 67: }

 69: /*TEST

 71:   test:
 72:     suffix: 0
 73:     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -interpolate
 74:     requires: exodusii

 76: TEST*/