program ex62f90
#include "petsc/finclude/petscdmplex.h"
use petscdmplex
implicit none
#include "exodusII.inc"
! Get the Fortran kind associated with PetscInt and PetscReal so that we can use literal constants.
PetscInt :: dummyPetscInt
PetscReal :: dummyPetscreal
integer, parameter :: kPI = kind(dummyPetscInt)
integer, parameter :: kPR = kind(dummyPetscReal)
type(tDM) :: dm, dmU, dmA, dmS, dmUA, dmUA2, pDM
type(tDM), dimension(:), pointer :: dmList
type(tVec) :: X, U, A, S, UA, UA2
type(tIS) :: isU, isA, isS, isUA
type(tPetscSection) :: section
PetscInt :: fieldU = 0
PetscInt :: fieldA = 2
PetscInt :: fieldS = 1
PetscInt, dimension(2) :: fieldUA = [0, 2]
character(len=PETSC_MAX_PATH_LEN) :: ifilename, ofilename, IOBuffer
integer :: exoid = -1
type(tIS) :: csIS
PetscInt, dimension(:), pointer :: csID
PetscInt, dimension(:), pointer :: pStartDepth, pEndDepth
PetscInt :: order = 1
Integer :: i
PetscInt :: sdim, d, pStart, pEnd, p, numCS, set, j
PetscMPIInt :: rank, numProc
PetscBool :: flg
PetscErrorCode :: ierr
MPI_Comm :: comm
type(tPetscViewer) :: viewer
Character(len=MXSTLN) :: sJunk
PetscInt :: numstep = 3, step
Integer :: numNodalVar, numZonalVar
character(len=MXNAME), dimension(4) :: nodalVarName2D = ["U_x ", &
"U_y ", &
"Alpha", &
"Beta "]
character(len=MXNAME), dimension(3) :: zonalVarName2D = ["Sigma_11", &
"Sigma_12", &
"Sigma_22"]
character(len=MXNAME), dimension(5) :: nodalVarName3D = ["U_x ", &
"U_y ", &
"U_z ", &
"Alpha", &
"Beta "]
character(len=MXNAME), dimension(6) :: zonalVarName3D = ["Sigma_11", &
"Sigma_22", &
"Sigma_33", &
"Sigma_23", &
"Sigma_13", &
"Sigma_12"]
logical, dimension(:, :), pointer :: truthtable
type(tIS) :: cellIS
PetscInt, dimension(:), pointer :: cellID
PetscInt :: numCells, cell, closureSize
PetscInt, dimension(:), pointer :: closureA, closure
type(tPetscSection) :: sectionUA, coordSection
type(tVec) :: UALoc, coord
PetscScalar, dimension(:), pointer :: cval, xyz
PetscInt :: dofUA, offUA, c
! dof layout ordered by increasing height in the DAG: cell, face, edge, vertex
PetscInt, dimension(3), target :: dofS2D = [0, 0, 3]
PetscInt, dimension(3), target :: dofUP1Tri = [2, 0, 0]
PetscInt, dimension(3), target :: dofAP1Tri = [1, 0, 0]
PetscInt, dimension(3), target :: dofUP2Tri = [2, 2, 0]
PetscInt, dimension(3), target :: dofAP2Tri = [1, 1, 0]
PetscInt, dimension(3), target :: dofUP1Quad = [2, 0, 0]
PetscInt, dimension(3), target :: dofAP1Quad = [1, 0, 0]
PetscInt, dimension(3), target :: dofUP2Quad = [2, 2, 2]
PetscInt, dimension(3), target :: dofAP2Quad = [1, 1, 1]
PetscInt, dimension(4), target :: dofS3D = [0, 0, 0, 6]
PetscInt, dimension(4), target :: dofUP1Tet = [3, 0, 0, 0]
PetscInt, dimension(4), target :: dofAP1Tet = [1, 0, 0, 0]
PetscInt, dimension(4), target :: dofUP2Tet = [3, 3, 0, 0]
PetscInt, dimension(4), target :: dofAP2Tet = [1, 1, 0, 0]
PetscInt, dimension(4), target :: dofUP1Hex = [3, 0, 0, 0]
PetscInt, dimension(4), target :: dofAP1Hex = [1, 0, 0, 0]
PetscInt, dimension(4), target :: dofUP2Hex = [3, 3, 3, 3]
PetscInt, dimension(4), target :: dofAP2Hex = [1, 1, 1, 1]
PetscInt, dimension(:), pointer :: dofU, dofA, dofS
type(tPetscSF) :: migrationSF
PetscPartitioner :: part
PetscLayout :: layout
type(tVec) :: tmpVec
PetscReal :: norm
PetscReal :: time = 1.234_kPR
PetscCallA(PetscInitialize(PETSC_NULL_CHARACTER, ierr))
if (ierr /= 0) then
print *, 'Unable to initialize PETSc'
stop
end if
PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr))
PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD, numProc, ierr))
PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-i', ifilename, flg, ierr))
PetscCheckA(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, 'missing input file name -i ')
PetscCallA(PetscOptionsGetString(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-o', ofilename, flg, ierr))
PetscCheckA(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, 'missing output file name -o