MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
tvdschememodule Module Reference

Data Types

interface  tvdschemetype
 

Functions/Subroutines

type(tvdschemetype) function constructor (dis, fmi, ibound)
 
type(coefficientstype) function, target compute (this, n, m, iposnm, phi)
 

Function/Subroutine Documentation

◆ compute()

type(coefficientstype) function, target tvdschememodule::compute ( class(tvdschemetype), target  this,
integer(i4b), intent(in)  n,
integer(i4b), intent(in)  m,
integer(i4b), intent(in)  iposnm,
real(dp), dimension(:), intent(in)  phi 
)
private

Definition at line 42 of file TVDScheme.f90.

43  !-- return
44  type(CoefficientsType), target :: phi_face
45  ! -- dummy
46  class(TVDSchemeType), target :: this
47  integer(I4B), intent(in) :: n
48  integer(I4B), intent(in) :: m
49  integer(I4B), intent(in) :: iposnm
50  real(DP), intent(in), dimension(:) :: phi
51  ! -- local
52  integer(I4B) :: ipos, isympos, iup, idn, i2up, j
53  real(DP) :: qnm, qmax, qupj, elupdn, elup2up
54  real(DP) :: smooth, cdiff, alimiter
55  real(DP), pointer :: coef_up, coef_dn
56  !
57  ! -- Find upstream node
58  isympos = this%dis%con%jas(iposnm)
59  qnm = this%fmi%gwfflowja(iposnm)
60  if (qnm > dzero) then
61  ! -- positive flow into n means m is upstream
62  iup = m
63  idn = n
64 
65  coef_up => phi_face%c_m
66  coef_dn => phi_face%c_n
67  else
68  iup = n
69  idn = m
70 
71  coef_up => phi_face%c_n
72  coef_dn => phi_face%c_m
73  end if
74  elupdn = this%dis%con%cl1(isympos) + this%dis%con%cl2(isympos)
75  !
76  ! -- Add low order terms
77  coef_up = done
78  !
79  ! -- Add high order terms
80  !
81  ! -- Find second node upstream to iup
82  i2up = 0
83  qmax = dzero
84  do ipos = this%dis%con%ia(iup) + 1, this%dis%con%ia(iup + 1) - 1
85  j = this%dis%con%ja(ipos)
86  if (this%ibound(j) == 0) cycle
87  qupj = this%fmi%gwfflowja(ipos)
88  isympos = this%dis%con%jas(ipos)
89  if (qupj > qmax) then
90  qmax = qupj
91  i2up = j
92  elup2up = this%dis%con%cl1(isympos) + this%dis%con%cl2(isympos)
93  end if
94  end do
95  !
96  ! -- Calculate flux limiting term
97  if (i2up > 0) then
98  smooth = dzero
99  cdiff = abs(phi(idn) - phi(iup))
100  if (cdiff > dprec) then
101  smooth = (phi(iup) - phi(i2up)) / elup2up * &
102  elupdn / (phi(idn) - phi(iup))
103  end if
104  if (smooth > dzero) then
105  alimiter = dtwo * smooth / (done + smooth)
106  phi_face%rhs = -dhalf * alimiter * (phi(idn) - phi(iup))
107  end if
108  end if
109 

◆ constructor()

type(tvdschemetype) function tvdschememodule::constructor ( class(disbasetype), intent(in), pointer  dis,
type(tspfmitype), intent(in), pointer  fmi,
integer(i4b), dimension(:), intent(in), pointer, contiguous  ibound 
)
private

Definition at line 28 of file TVDScheme.f90.

29  ! -- return
30  type(TVDSchemeType) :: interpolation_scheme
31  ! --dummy
32  class(DisBaseType), pointer, intent(in) :: dis
33  type(TspFmiType), pointer, intent(in) :: fmi
34  integer(I4B), dimension(:), pointer, contiguous, intent(in) :: ibound
35 
36  interpolation_scheme%dis => dis
37  interpolation_scheme%fmi => fmi
38  interpolation_scheme%ibound => ibound
39 
Here is the caller graph for this function: