xref: /petsc/src/dm/impls/plex/tests/ex97f90.F90 (revision 9b88ac225e01f016352a5f4cd90e158abe5f5675)
13ead50e0SBlaise Bourdin#include "petsc/finclude/petsc.h"
2c5e229c2SMartin Diehlprogram ex97f90
33ead50e0SBlaise Bourdin  use petsc
43ead50e0SBlaise Bourdin  implicit none
53ead50e0SBlaise Bourdin
63ead50e0SBlaise Bourdin  ! Get the fortran kind associated with PetscInt and PetscReal so that we can use literal constants.
73ead50e0SBlaise Bourdin  PetscInt                           :: dummyPetscInt
83ead50e0SBlaise Bourdin  PetscReal                          :: dummyPetscreal
93ead50e0SBlaise Bourdin  integer, parameter                  :: kPI = kind(dummyPetscInt)
103ead50e0SBlaise Bourdin  integer, parameter                  :: kPR = kind(dummyPetscReal)
113ead50e0SBlaise Bourdin
123ead50e0SBlaise Bourdin  type(tDM)                          :: dm
133ead50e0SBlaise Bourdin  type(tDMLabel)                     :: label
143ead50e0SBlaise Bourdin  character(len=PETSC_MAX_PATH_LEN)  :: ifilename, iobuffer
153ead50e0SBlaise Bourdin  DMPolytopeType                     :: cellType
163ead50e0SBlaise Bourdin  PetscInt                           :: pStart, pEnd, p
173ead50e0SBlaise Bourdin  PetscErrorCode                     :: ierr
183ead50e0SBlaise Bourdin  PetscBool                          :: flg
193ead50e0SBlaise Bourdin
203ead50e0SBlaise Bourdin  PetscCallA(PetscInitialize(ierr))
213ead50e0SBlaise Bourdin
22dcb3e689SBarry Smith  PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-i', ifilename, flg, ierr))
23dcb3e689SBarry Smith  PetscCheckA(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, 'missing input file name -i <input file name>')
243ead50e0SBlaise Bourdin
253ead50e0SBlaise Bourdin  PetscCallA(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename, PETSC_NULL_CHARACTER, PETSC_TRUE, dm, ierr))
263ead50e0SBlaise Bourdin  PetscCallA(DMPlexDistributeSetDefault(dm, PETSC_FALSE, ierr))
27dcb3e689SBarry Smith  PetscCallA(PetscObjectSetName(dm, 'ex97f90', ierr))
283ead50e0SBlaise Bourdin  PetscCallA(DMSetFromOptions(dm, ierr))
29ce78bad3SBarry Smith  PetscCallA(DMViewFromOptions(dm, PETSC_NULL_OBJECT, '-dm_view', ierr))
303ead50e0SBlaise Bourdin
313ead50e0SBlaise Bourdin  PetscCallA(DMGetLabel(dm, 'celltype', label, ierr))
323ead50e0SBlaise Bourdin  PetscCallA(DMLabelView(label, PETSC_VIEWER_STDOUT_WORLD, ierr))
333ead50e0SBlaise Bourdin  PetscCallA(DMPlexGetHeightStratum(dm, 0_kPI, pStart, pEnd, ierr))
34*02c639afSMartin Diehl  do p = pStart, pEnd - 1
353ead50e0SBlaise Bourdin    PetscCallA(DMPlexGetCellType(dm, p, cellType, ierr))
36dcb3e689SBarry Smith    write (IOBuffer, '("cell: ",i3," type: ",i3,"\n")') p, cellType
373ead50e0SBlaise Bourdin    PetscCallA(PetscPrintf(PETSC_COMM_SELF, IOBuffer, ierr))
38*02c639afSMartin Diehl  end do
393ead50e0SBlaise Bourdin  PetscCallA(DMDestroy(dm, ierr))
403ead50e0SBlaise Bourdin
413ead50e0SBlaise Bourdin  PetscCallA(PetscFinalize(ierr))
423ead50e0SBlaise Bourdinend program ex97f90
433ead50e0SBlaise Bourdin
443ead50e0SBlaise Bourdin! /*TEST
453ead50e0SBlaise Bourdin!   build:
463ead50e0SBlaise Bourdin!     requires: !complex
473ead50e0SBlaise Bourdin!   testset:
483ead50e0SBlaise Bourdin!     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh -dm_view
493ead50e0SBlaise Bourdin!     nsize: 1
503ead50e0SBlaise Bourdin!     test:
513ead50e0SBlaise Bourdin!       suffix: 0
523ead50e0SBlaise Bourdin!       args:
533ead50e0SBlaise Bourdin! TEST*/
54