MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
TVDScheme.f90
Go to the documentation of this file.
2  use kindmodule, only: dp, i4b
6  use basedismodule, only: disbasetype
7  use tspfmimodule, only: tspfmitype
8 
9  implicit none
10  private
11 
12  public :: tvdschemetype
13 
15  private
16  class(disbasetype), pointer :: dis
17  type(tspfmitype), pointer :: fmi
18  real(dp), dimension(:), pointer :: phi
19  integer(I4B), dimension(:), pointer, contiguous :: ibound => null() !< pointer to model ibound
20  contains
21  procedure :: compute
22  procedure :: set_field
23  end type tvdschemetype
24 
25  interface tvdschemetype
26  module procedure constructor
27  end interface tvdschemetype
28 
29 contains
30  function constructor(dis, fmi, ibound) result(interpolation_scheme)
31  ! -- return
32  type(tvdschemetype) :: interpolation_scheme
33  ! --dummy
34  class(disbasetype), pointer, intent(in) :: dis
35  type(tspfmitype), pointer, intent(in) :: fmi
36  integer(I4B), dimension(:), pointer, contiguous, intent(in) :: ibound
37 
38  interpolation_scheme%dis => dis
39  interpolation_scheme%fmi => fmi
40  interpolation_scheme%ibound => ibound
41 
42  end function constructor
43 
44  !> @brief Set the scalar field for which interpolation will be computed
45  !!
46  !! This method establishes a pointer to the scalar field data for
47  !! subsequent TVD interpolation computations.
48  !<
49  subroutine set_field(this, phi)
50  ! -- dummy
51  class(tvdschemetype), target :: this
52  real(DP), intent(in), dimension(:), pointer :: phi
53 
54  this%phi => phi
55  end subroutine set_field
56 
57  function compute(this, n, m, iposnm) result(phi_face)
58  !-- return
59  type(coefficientstype), target :: phi_face
60  ! -- dummy
61  class(tvdschemetype), target :: this
62  integer(I4B), intent(in) :: n
63  integer(I4B), intent(in) :: m
64  integer(I4B), intent(in) :: iposnm
65  ! -- local
66  integer(I4B) :: ipos, isympos, iup, idn, i2up, j
67  real(dp) :: qnm, qmax, qupj, elupdn, elup2up
68  real(dp) :: smooth, cdiff, alimiter
69  real(dp), pointer :: coef_up, coef_dn
70  !
71  ! -- Find upstream node
72  isympos = this%dis%con%jas(iposnm)
73  qnm = this%fmi%gwfflowja(iposnm)
74  if (qnm > dzero) then
75  ! -- positive flow into n means m is upstream
76  iup = m
77  idn = n
78 
79  coef_up => phi_face%c_m
80  coef_dn => phi_face%c_n
81  else
82  iup = n
83  idn = m
84 
85  coef_up => phi_face%c_n
86  coef_dn => phi_face%c_m
87  end if
88  elupdn = this%dis%con%cl1(isympos) + this%dis%con%cl2(isympos)
89  !
90  ! -- Add low order terms
91  coef_up = done
92  !
93  ! -- Add high order terms
94  !
95  ! -- Find second node upstream to iup
96  i2up = 0
97  qmax = dzero
98  do ipos = this%dis%con%ia(iup) + 1, this%dis%con%ia(iup + 1) - 1
99  j = this%dis%con%ja(ipos)
100  if (this%ibound(j) == 0) cycle
101  qupj = this%fmi%gwfflowja(ipos)
102  isympos = this%dis%con%jas(ipos)
103  if (qupj > qmax) then
104  qmax = qupj
105  i2up = j
106  elup2up = this%dis%con%cl1(isympos) + this%dis%con%cl2(isympos)
107  end if
108  end do
109  !
110  ! -- Calculate flux limiting term
111  if (i2up > 0) then
112  smooth = dzero
113  cdiff = abs(this%phi(idn) - this%phi(iup))
114  if (cdiff > dprec) then
115  smooth = (this%phi(iup) - this%phi(i2up)) / elup2up * &
116  elupdn / (this%phi(idn) - this%phi(iup))
117  end if
118  if (smooth > dzero) then
119  alimiter = dtwo * smooth / (done + smooth)
120  phi_face%rhs = -dhalf * alimiter * (this%phi(idn) - this%phi(iup))
121  end if
122  end if
123 
124  end function compute
125 
126 end module tvdschememodule
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:68
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter dprec
real constant machine precision
Definition: Constants.f90:120
real(dp), parameter dtwo
real constant 2
Definition: Constants.f90:79
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
This module defines variable data types.
Definition: kind.f90:8
type(coefficientstype) function, target compute(this, n, m, iposnm)
Definition: TVDScheme.f90:58
type(tvdschemetype) function constructor(dis, fmi, ibound)
Definition: TVDScheme.f90:31
subroutine set_field(this, phi)
Set the scalar field for which interpolation will be computed.
Definition: TVDScheme.f90:50