xref: /petsc/src/dm/impls/plex/tests/ex97f90.F90 (revision 3ead50e02836b3e0f420ae7d56092ac5c1313c0c)
1*3ead50e0SBlaise Bourdinprogram ex97f90
2*3ead50e0SBlaise Bourdin#include "petsc/finclude/petsc.h"
3*3ead50e0SBlaise Bourdin    use petsc
4*3ead50e0SBlaise Bourdin    implicit none
5*3ead50e0SBlaise Bourdin
6*3ead50e0SBlaise Bourdin    ! Get the fortran kind associated with PetscInt and PetscReal so that we can use literal constants.
7*3ead50e0SBlaise Bourdin    PetscInt                           :: dummyPetscInt
8*3ead50e0SBlaise Bourdin    PetscReal                          :: dummyPetscreal
9*3ead50e0SBlaise Bourdin    integer,parameter                  :: kPI = kind(dummyPetscInt)
10*3ead50e0SBlaise Bourdin    integer,parameter                  :: kPR = kind(dummyPetscReal)
11*3ead50e0SBlaise Bourdin
12*3ead50e0SBlaise Bourdin    type(tDM)                          :: dm
13*3ead50e0SBlaise Bourdin    type(tDMLabel)                     :: label
14*3ead50e0SBlaise Bourdin    character(len=PETSC_MAX_PATH_LEN)  :: ifilename,iobuffer
15*3ead50e0SBlaise Bourdin    DMPolytopeType                     :: cellType
16*3ead50e0SBlaise Bourdin    PetscInt                           :: pStart,pEnd,p
17*3ead50e0SBlaise Bourdin    PetscErrorCode                     :: ierr
18*3ead50e0SBlaise Bourdin    PetscBool                          :: flg
19*3ead50e0SBlaise Bourdin
20*3ead50e0SBlaise Bourdin    PetscCallA(PetscInitialize(ierr))
21*3ead50e0SBlaise Bourdin
22*3ead50e0SBlaise Bourdin    PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,"-i",ifilename,flg,ierr))
23*3ead50e0SBlaise Bourdin    if (.not. flg) then
24*3ead50e0SBlaise Bourdin        SETERRA(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"missing input file name -i <input file name>")
25*3ead50e0SBlaise Bourdin    end if
26*3ead50e0SBlaise Bourdin
27*3ead50e0SBlaise Bourdin    PetscCallA(DMPlexCreateFromFile(PETSC_COMM_WORLD,ifilename,PETSC_NULL_CHARACTER,PETSC_TRUE,dm,ierr))
28*3ead50e0SBlaise Bourdin    PetscCallA(DMPlexDistributeSetDefault(dm,PETSC_FALSE,ierr))
29*3ead50e0SBlaise Bourdin    PetscCallA(PetscObjectSetName(dm,"ex97f90",ierr))
30*3ead50e0SBlaise Bourdin    PetscCallA(DMSetFromOptions(dm,ierr))
31*3ead50e0SBlaise Bourdin    PetscCallA(DMViewFromOptions(dm,PETSC_NULL_OPTIONS,"-dm_view",ierr))
32*3ead50e0SBlaise Bourdin
33*3ead50e0SBlaise Bourdin    PetscCallA(DMGetLabel(dm,'celltype',label,ierr))
34*3ead50e0SBlaise Bourdin    PetscCallA(DMLabelView(label,PETSC_VIEWER_STDOUT_WORLD,ierr))
35*3ead50e0SBlaise Bourdin    PetscCallA(DMPlexGetHeightStratum(dm,0_kPI,pStart,pEnd,ierr))
36*3ead50e0SBlaise Bourdin    Do p = pStart,pEnd-1
37*3ead50e0SBlaise Bourdin        PetscCallA(DMPlexGetCellType(dm,p,cellType,ierr))
38*3ead50e0SBlaise Bourdin        Write(IOBuffer,'("cell: ",i3," type: ",i3,"\n")' ) p,cellType
39*3ead50e0SBlaise Bourdin        PetscCallA(PetscPrintf(PETSC_COMM_SELF,IOBuffer,ierr))
40*3ead50e0SBlaise Bourdin    End Do
41*3ead50e0SBlaise Bourdin    PetscCallA(DMDestroy(dm,ierr))
42*3ead50e0SBlaise Bourdin
43*3ead50e0SBlaise Bourdin    PetscCallA(PetscFinalize(ierr))
44*3ead50e0SBlaise Bourdinend program ex97f90
45*3ead50e0SBlaise Bourdin
46*3ead50e0SBlaise Bourdin! /*TEST
47*3ead50e0SBlaise Bourdin!   build:
48*3ead50e0SBlaise Bourdin!     requires: !complex
49*3ead50e0SBlaise Bourdin!   testset:
50*3ead50e0SBlaise Bourdin!     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh -dm_view
51*3ead50e0SBlaise Bourdin!     nsize: 1
52*3ead50e0SBlaise Bourdin!     test:
53*3ead50e0SBlaise Bourdin!       suffix: 0
54*3ead50e0SBlaise Bourdin!       args:
55*3ead50e0SBlaise Bourdin! TEST*/
56