MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
utvdschememodule Module Reference

Data Types

interface  utvdschemetype
 Total Variation Diminishing (TVD) interpolation scheme. More...
 

Functions/Subroutines

type(utvdschemetype) function constructor (dis, fmi, gradient)
 
subroutine destructor (this)
 
subroutine set_field (this, phi)
 Set the scalar field for which interpolation will be computed. More...
 
type(coefficientstype) function, target compute (this, n, m, iposnm)
 
real(dp) function limiter (this, r)
 
subroutine compute_node_distance (this)
 

Function/Subroutine Documentation

◆ compute()

type(coefficientstype) function, target utvdschememodule::compute ( class(utvdschemetype), target  this,
integer(i4b), intent(in)  n,
integer(i4b), intent(in)  m,
integer(i4b), intent(in)  iposnm 
)
private

Definition at line 106 of file UTVDScheme.f90.

107  !-- return
108  type(CoefficientsType), target :: phi_face ! Output: coefficients for the face between cells n and m
109  ! -- dummy
110  class(UTVDSchemeType), target :: this
111  integer(I4B), intent(in) :: n ! Index of the first cell
112  integer(I4B), intent(in) :: m ! Index of the second cell
113  integer(I4B), intent(in) :: iposnm ! Position in the connectivity array for the n-m connection
114  ! -- local
115  integer(I4B) :: iup, idn, isympos ! Indices for upwind, downwind, and symmetric position in connectivity
116  real(DP) :: qnm ! Flow rate across the face between n and m
117  real(DP), pointer :: coef_up, coef_dn ! Pointers to upwind and downwind coefficients in phi_face
118  real(DP), dimension(3) :: grad_c ! gradient at upwind cell
119  real(DP), dimension(3) :: dnm ! vector from upwind to downwind cell
120  real(DP) :: smooth ! Smoothness indicator for limiter i.e. ratio of gradients
121  real(DP) :: alimiter ! Value of the TVD limiter
122  real(DP) :: cl1, cl2 ! Connection length of cell n to the face and of cell m to the face
123  real(DP) :: cl_up, cl_dn ! Connection length of cell iup to the face and of cell idn to the face
124  real(DP) :: relative_distance ! Relative distance factor for high-order term
125  real(DP) :: c_virtual ! Virtual node concentration (Darwish method)
126  real(DP), pointer :: min_phi, max_phi ! Local minimum and maximum among cell and neighbors
127 
128  isympos = this%dis%con%jas(iposnm)
129  qnm = this%fmi%gwfflowja(iposnm)
130  !
131  ! -- Get the distance from node n to the face (cl1) and from node m to the face (cl2).
132  ! The distances are dependent on the the node numbering convention and the direction of the connection.
133  if (n < m) then
134  cl1 = this%dis%con%cl1(isympos)
135  cl2 = this%dis%con%cl2(isympos)
136  else
137  cl1 = this%dis%con%cl2(isympos)
138  cl2 = this%dis%con%cl1(isympos)
139  end if
140  !
141  ! -- Find upstream node
142  if (qnm > dzero) then
143  ! -- positive flow into n means m is upstream
144  iup = m
145  idn = n
146 
147  cl_up = cl2
148  cl_dn = cl1
149 
150  coef_up => phi_face%c_m
151  coef_dn => phi_face%c_n
152  else
153  iup = n
154  idn = m
155 
156  cl_up = cl1
157  cl_dn = cl2
158 
159  coef_up => phi_face%c_n
160  coef_dn => phi_face%c_m
161  end if
162  !
163  ! Determine direction of distance vector from upwind to downwind cell
164  ! The cached_node_distance always stores vector from lower-numbered node to higher-numbered node.
165  ! Since we need dnm to point from upwind (iup) to downwind (idn), we must adjust the sign:
166  ! - If iup > idn: the cached vector points from idn to iup, so we negate it to get iup to idn
167  ! - If iup < idn: the cached vector already points from iup to idn, so use it as-is
168  if (iup > idn) then
169  dnm = -this%cached_node_distance(isympos, :)
170  else
171  dnm = this%cached_node_distance(isympos, :)
172  end if
173  !
174  ! -- Add low order terms
175  coef_up = done
176  !
177  ! -- Add high order terms
178  !
179  ! -- Return if straddled cells have same value
180  if (abs(this%phi(idn) - this%phi(iup)) < dsame) return
181  !
182  ! -- Compute cell concentration gradient
183  grad_c = this%gradient%get(iup)
184  !
185  ! Darwish's method to compute virtual node concentration
186  c_virtual = this%phi(idn) - 2.0_dp * (dot_product(grad_c, dnm))
187  !
188  ! Enforce local TVD condition.
189  ! This is done by limiting the virtual concentration to the range of
190  ! the max and min concentration of the neighbouring cells.
191  min_phi => this%min_max_phi%get_min(iup)
192  max_phi => this%min_max_phi%get_max(iup)
193 
194  if (c_virtual > max_phi) then
195  c_virtual = max_phi
196  end if
197 
198  if (c_virtual < max(min_phi, dzero)) then
199  c_virtual = max(min_phi, dzero)
200  end if
201  !
202  ! -- Compute smoothness factor
203  smooth = (this%phi(iup) - c_virtual) / (this%phi(idn) - this%phi(iup))
204  !
205  ! -- Compute limiter
206  alimiter = this%limiter(smooth)
207 
208  ! High order term is:
209  relative_distance = cl_up / (cl_up + cl_dn)
210  phi_face%rhs = -relative_distance * alimiter * (this%phi(idn) - this%phi(iup))
211 
212  ! Alternative way of writing the high order term by adding it to the
213  ! coefficients matrix. The equation to be added is:
214  ! high_order = cl_up / (cl_up + cl_dn) * alimiter * qnm * (phi(idn) - phi(iup))
215  ! This is split into two parts:
216  ! coef_up = coef_up - relative_distance * alimiter
217  ! coef_dn = coef_dn + relative_distance * alimiter
218 

◆ compute_node_distance()

subroutine utvdschememodule::compute_node_distance ( class(utvdschemetype), target  this)
private

Definition at line 246 of file UTVDScheme.f90.

247  ! -- dummy
248  class(UTVDSchemeType), target :: this
249  ! -- local
250  integer(I4B) :: n, m, ipos, isympos
251 
252  this%cached_node_distance = 0.0_dp
253  do n = 1, this%dis%nodes
254  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
255  m = this%dis%con%ja(ipos)
256  if (m <= n) cycle
257 
258  isympos = this%dis%con%jas(ipos)
259  this%cached_node_distance(isympos, :) = node_distance(this%dis, n, m)
260  end do
261  end do
262 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ constructor()

type(utvdschemetype) function utvdschememodule::constructor ( class(disbasetype), intent(in), pointer  dis,
type(tspfmitype), intent(in), pointer  fmi,
class(igradienttype), intent(in), allocatable, target  gradient 
)
private

Definition at line 62 of file UTVDScheme.f90.

64  ! -- return
65  type(UTVDSchemeType) :: interpolation_scheme
66  ! -- dummy
67  class(DisBaseType), pointer, intent(in) :: dis
68  type(TspFmiType), pointer, intent(in) :: fmi
69  class(IGradientType), allocatable, intent(in), target :: gradient
70 
71  interpolation_scheme%dis => dis
72  interpolation_scheme%fmi => fmi
73 
74  interpolation_scheme%gradient => gradient
75  interpolation_scheme%min_max_phi = localcellextrematype(dis)
76 
77  allocate (interpolation_scheme%cached_node_distance(dis%njas, 3))
78  call compute_node_distance(interpolation_scheme)
79 
Here is the caller graph for this function:

◆ destructor()

subroutine utvdschememodule::destructor ( type(utvdschemetype), intent(inout)  this)
private

Definition at line 82 of file UTVDScheme.f90.

83  ! -- dummy
84  type(UTVDSchemeType), intent(inout) :: this
85 
86  deallocate (this%cached_node_distance)
Here is the caller graph for this function:

◆ limiter()

real(dp) function utvdschememodule::limiter ( class(utvdschemetype this,
real(dp)  r 
)
private

Definition at line 221 of file UTVDScheme.f90.

222  ! -- return
223  real(DP) :: theta ! limited slope
224  ! -- dummy
225  class(UTVDSchemeType) :: this
226  real(DP) :: r ! ratio of successive gradients
227 
228  select case (this%limiter_id)
229  case (2) ! van Leer
230  theta = max(0.0_dp, min((r + dabs(r)) / (1.0_dp + dabs(r)), 2.0_dp))
231  case (3) ! Koren
232  theta = max(0.0_dp, min(2.0_dp * r, &
233  1.0_dp / 3.0_dp + 2.0_dp / 3.0_dp * r, 2.0_dp))
234  case (4) ! Superbee
235  theta = max(0.0_dp, min(2.0_dp * r, 1.0_dp), min(r, 2.0_dp))
236  case (5) ! van Albada
237  theta = max(0.0_dp, (r * r + r) / (r * r + 1.0_dp))
238  case (6) ! Koren modified
239  theta = max(0.0_dp, min(4.0_dp * r * r + r, &
240  1.0_dp / 3.0_dp + 2.0_dp / 3.0_dp * r, 2.0_dp))
241  case default
242  theta = dzero
243  end select

◆ set_field()

subroutine utvdschememodule::set_field ( class(utvdschemetype), target  this,
real(dp), dimension(:), intent(in), pointer  phi 
)
private

This method establishes a pointer to the scalar field data and updates any dependent cached data (gradients and local extrema) to ensure subsequent interpolation computations use the current field values.

Definition at line 96 of file UTVDScheme.f90.

97  ! -- dummy
98  class(UTVDSchemeType), target :: this
99  real(DP), intent(in), dimension(:), pointer :: phi
100 
101  this%phi => phi
102  call this%gradient%set_field(phi)
103  call this%min_max_phi%set_field(phi)