Actual source code: ex97f90.F90

  1: #include "petsc/finclude/petsc.h"
  2: program ex97f90
  3:   use petsc
  4:   implicit none

  6:   type(tDM)                          :: dm
  7:   type(tDMLabel)                     :: label
  8:   character(len=PETSC_MAX_PATH_LEN)  :: ifilename, iobuffer
  9:   DMPolytopeType                     :: cellType
 10:   PetscInt                           :: pStart, pEnd, p
 11:   PetscErrorCode                     :: ierr
 12:   PetscBool                          :: flg

 14:   PetscCallA(PetscInitialize(ierr))

 16:   PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-i', ifilename, flg, ierr))
 17:   PetscCheckA(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, 'missing input file name -i ')

 19:   PetscCallA(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename, PETSC_NULL_CHARACTER, PETSC_TRUE, dm, ierr))
 20:   PetscCallA(DMPlexDistributeSetDefault(dm, PETSC_FALSE, ierr))
 21:   PetscCallA(PetscObjectSetName(dm, 'ex97f90', ierr))
 22:   PetscCallA(DMSetFromOptions(dm, ierr))
 23:   PetscCallA(DMViewFromOptions(dm, PETSC_NULL_OBJECT, '-dm_view', ierr))

 25:   PetscCallA(DMGetLabel(dm, 'celltype', label, ierr))
 26:   PetscCallA(DMLabelView(label, PETSC_VIEWER_STDOUT_WORLD, ierr))
 27:   PetscCallA(DMPlexGetHeightStratum(dm, 0_PETSC_INT_KIND, pStart, pEnd, ierr))
 28:   do p = pStart, pEnd - 1
 29:     PetscCallA(DMPlexGetCellType(dm, p, cellType, ierr))
 30:     write (IOBuffer, '("cell: ",i3," type: ",i3,"\n")') p, cellType
 31:     PetscCallA(PetscPrintf(PETSC_COMM_SELF, IOBuffer, ierr))
 32:   end do
 33:   PetscCallA(DMDestroy(dm, ierr))

 35:   PetscCallA(PetscFinalize(ierr))
 36: end program ex97f90

 38: ! /*TEST
 39: !   build:
 40: !     requires: !complex
 41: !   testset:
 42: !     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh -dm_view
 43: !     nsize: 1
 44: !     test:
 45: !       suffix: 0
 46: !       args:
 47: ! TEST*/