Actual source code: ex1f90.F90

  1: program ex1f90
  2: #include <petsc/finclude/petscdmlabel.h>
  3:   use petscdm
  4:   implicit NONE

  6:   type(tDM)                         :: dm, dmDist
  7:   character(len=PETSC_MAX_PATH_LEN) :: filename
  8:   PetscBool                         :: interpolate = PETSC_FALSE
  9:   PetscBool                         :: flg
 10:   PetscErrorCode                    :: ierr
 11:   PetscInt                          :: izero
 12:   izero = 0

 14:   PetscCallA(PetscInitialize(ierr))
 15:   PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-i', filename, flg, ierr))
 16:   PetscCallA(PetscOptionsGetBool(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-interpolate', interpolate, flg, ierr))

 18:   PetscCallA(DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, 'ex1f90_plex', interpolate, dm, ierr))
 19:   PetscCallA(DMPlexDistribute(dm, izero, PETSC_NULL_SF, dmDist, ierr))
 20:   if (.not. PetscObjectIsNull(dmDist)) then
 21:     PetscCallA(DMDestroy(dm, ierr))
 22:     dm = dmDist
 23:   end if

 25:   PetscCallA(ViewLabels(dm, PETSC_VIEWER_STDOUT_WORLD, ierr))
 26:   PetscCallA(DMDestroy(dm, ierr))
 27:   PetscCallA(PetscFinalize(ierr))

 29: contains
 30:   subroutine ViewLabels(dm, viewer, ierr)
 31:     type(tDM)                        :: dm
 32:     type(tPetscViewer)               :: viewer
 33:     PetscErrorCode                   :: ierr

 35:     DMLabel                          :: label
 36:     type(tIS)                        :: labelIS
 37:     character(len=PETSC_MAX_PATH_LEN):: labelName, IObuffer
 38:     PetscInt                         :: numLabels, l

 40:     PetscCall(DMGetNumLabels(dm, numLabels, ierr))
 41:     write (IObuffer, *) 'Number of labels: ', numLabels, '\n'
 42:     PetscCall(PetscViewerASCIIPrintf(viewer, IObuffer, ierr))
 43:     do l = 0, numLabels - 1
 44:       PetscCall(DMGetLabelName(dm, l, labelName, ierr))
 45:       write (IObuffer, *) 'label ', l, ' name: ', trim(labelName), '\n'
 46:       PetscCall(PetscViewerASCIIPrintf(viewer, IObuffer, ierr))

 48:       PetscCall(PetscViewerASCIIPrintf(viewer, 'IS of values\n', ierr))
 49:       PetscCall(DMGetLabel(dm, labelName, label, ierr))
 50:       PetscCall(DMLabelGetValueIS(label, labelIS, ierr))
 51: !      PetscCall(PetscViewerASCIIPushTab(viewer,ierr))
 52:       PetscCall(ISView(labelIS, viewer, ierr))
 53: !      PetscCall(PetscViewerASCIIPopTab(viewer,ierr))
 54:       PetscCall(ISDestroy(labelIS, ierr))
 55:       PetscCall(PetscViewerASCIIPrintf(viewer, '\n', ierr))
 56:     end do

 58:     PetscCall(PetscViewerASCIIPrintf(viewer, '\n\nCell Set label IS\n', ierr))
 59:     PetscCall(DMGetLabel(dm, 'Cell Sets', label, ierr))
 60:     PetscCall(DMLabelGetValueIS(label, labelIS, ierr))
 61:     PetscCall(ISView(labelIS, viewer, ierr))
 62:     PetscCall(ISDestroy(labelIS, ierr))
 63:   end subroutine viewLabels
 64: end program ex1F90

 66: !/*TEST
 67: !
 68: !  test:
 69: !    suffix: 0
 70: !    args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -interpolate
 71: !    requires: exodusii
 72: !
 73: !TEST*/