MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
UTVDScheme.f90
Go to the documentation of this file.
2  use kindmodule, only: dp, i4b
3  use constantsmodule, only: done, dzero, dsame, dhalf
6  use basedismodule, only: disbasetype
7  use tspfmimodule, only: tspfmitype
8  use igradient, only: igradienttype
9 
10  use disutilsmodule, only: node_distance
12 
13  implicit none
14  private
15 
16  public :: utvdschemetype
17 
18  !> @brief Total Variation Diminishing (TVD) interpolation scheme.
19  !!
20  !! This class implements a high-resolution, TVD interpolation scheme for use in transport modeling.
21  !! It extends a generic interpolation scheme interface and supports multiple TVD limiters (van Leer,
22  !! Koren, Superbee, van Albada, etc.) for controlling numerical diffusion and oscillations.
23  !! The default limiter is van Leer, but others can be selected by changing the `limiter_id` member.
24  !!
25  !! The scheme uses a combination of low-order upwind and high-order limited terms to compute
26  !! face concentrations. The high-order term is constructed using a gradient-based virtual node
27  !! value, following the approach described by Darwish et al. An additional TVD clamp is applied
28  !! to the virtual node value to enforce TVD compliance, especially on grids where the original
29  !! method may not guarantee monotonicity.
30  !!
31  !! - Supports both structured and unstructured grids via polymorphic discretization and gradient objects.
32  !! - The limiter can be selected via the `limiter_id` member (default is van Leer).
33  !! - The method `find_local_extrema` finds the minimum and maximum values among the current cell and its neighbors,
34  !! which is used to enforce the TVD condition.
35  !! - The `compute` method calculates the face coefficients for the transport equation.
36  !<
38  private
39  class(disbasetype), pointer :: dis
40  type(tspfmitype), pointer :: fmi
41  class(igradienttype), pointer :: gradient
42 
43  real(dp), dimension(:), pointer :: phi
44  type(localcellextrematype), allocatable :: min_max_phi ! local minimum values at nodes
45  integer(I4B) :: limiter_id = 2 ! default to van Leer limiter
46  real(dp), dimension(:, :), allocatable :: cached_node_distance ! distance vectors
47  contains
48  procedure :: compute
49  procedure :: set_field
50  final :: destructor
51 
52  procedure, private :: limiter
53  procedure, private :: compute_node_distance
54  end type utvdschemetype
55 
56  interface utvdschemetype
57  module procedure constructor
58  end interface utvdschemetype
59 
60 contains
61 
62  function constructor(dis, fmi, gradient) &
63  result(interpolation_scheme)
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 
80  end function constructor
81 
82  subroutine destructor(this)
83  ! -- dummy
84  type(utvdschemetype), intent(inout) :: this
85 
86  deallocate (this%cached_node_distance)
87  end subroutine destructor
88 
89  !> @brief Set the scalar field for which interpolation will be computed
90  !!
91  !! This method establishes a pointer to the scalar field data and updates
92  !! any dependent cached data (gradients and local extrema) to ensure
93  !! subsequent interpolation computations use the current field values.
94  !!
95  !<
96  subroutine set_field(this, phi)
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)
104  end subroutine set_field
105 
106  function compute(this, n, m, iposnm) result(phi_face)
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 
219  end function compute
220 
221  function limiter(this, r) result(theta)
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
244  end function
245 
246  subroutine compute_node_distance(this)
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 
263  end subroutine compute_node_distance
264 
265 end module utvdschememodule
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dsame
real constant for values that are considered the same based on machine precision
Definition: Constants.f90:122
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 done
real constant 1
Definition: Constants.f90:76
real(dp) function, dimension(3), public node_distance(dis, n, m)
Returns the vector distance from cell n to cell m.
Definition: DisUtils.f90:39
This module defines variable data types.
Definition: kind.f90:8
type(coefficientstype) function, target compute(this, n, m, iposnm)
Definition: UTVDScheme.f90:107
subroutine set_field(this, phi)
Set the scalar field for which interpolation will be computed.
Definition: UTVDScheme.f90:97
subroutine destructor(this)
Definition: UTVDScheme.f90:83
type(utvdschemetype) function constructor(dis, fmi, gradient)
Definition: UTVDScheme.f90:64
subroutine compute_node_distance(this)
Definition: UTVDScheme.f90:247
real(dp) function limiter(this, r)
Definition: UTVDScheme.f90:222
Abstract interface for cell-based gradient computation.
Definition: IGradient.f90:15
Computes and caches local extrema for each cell and its connected neighbors.
Total Variation Diminishing (TVD) interpolation scheme.
Definition: UTVDScheme.f90:37