MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
tsp-adv.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b
6  use basedismodule, only: disbasetype
7  use tspfmimodule, only: tspfmitype
10  ! -- Gradient schemes
11  use igradient, only: igradienttype
14  ! -- Interpolation schemes
22 
23  implicit none
24  private
25  public :: tspadvtype
26  public :: adv_cr
27 
28  type, extends(numericalpackagetype) :: tspadvtype
29  integer(I4B), pointer :: iadvwt => null() !< advection scheme. See ADV_SCHEME_* constants
30  real(dp), pointer :: ats_percel => null() !< user-specified fractional number of cells advection can move a particle during one time step
31  integer(I4B), dimension(:), pointer, contiguous :: ibound => null() !< pointer to model ibound
32  type(tspfmitype), pointer :: fmi => null() !< pointer to fmi object
33  real(dp), pointer :: eqnsclfac => null() !< governing equation scale factor; =1. for solute; =rhow*cpw for energy
34 
35  class(interpolationschemeinterface), allocatable :: face_interpolation !< interpolation scheme for face values
36  class(igradienttype), allocatable :: gradient !< cell centered gradient
37  contains
38 
39  procedure :: adv_df
40  procedure :: adv_ar
41  procedure :: adv_dt
42  procedure :: adv_fc
43  procedure :: adv_cq
44  procedure :: adv_da
45 
46  procedure :: allocate_scalars
47  procedure, private :: source_options
48 
49  end type tspadvtype
50 
51 contains
52 
53  !> @ brief Create a new ADV object
54  !!
55  !! Create a new ADV package
56  !<
57  subroutine adv_cr(advobj, name_model, input_mempath, inunit, iout, fmi, &
58  eqnsclfac)
59  ! -- dummy
60  type(tspadvtype), pointer :: advobj
61  character(len=*), intent(in) :: name_model
62  character(len=*), intent(in) :: input_mempath
63  integer(I4B), intent(in) :: inunit
64  integer(I4B), intent(in) :: iout
65  type(tspfmitype), intent(in), target :: fmi
66  real(dp), intent(in), pointer :: eqnsclfac !< governing equation scale factor
67  !
68  ! -- Create the object
69  allocate (advobj)
70  !
71  ! -- create name and memory path
72  call advobj%set_names(1, name_model, 'ADV', 'ADV', input_mempath)
73  !
74  ! -- Allocate scalars
75  call advobj%allocate_scalars()
76  !
77  ! -- Set variables
78  advobj%inunit = inunit
79  advobj%iout = iout
80  advobj%fmi => fmi
81  advobj%eqnsclfac => eqnsclfac
82  end subroutine adv_cr
83 
84  !> @brief Define ADV object
85  !!
86  !! Define the ADV package
87  !<
88  subroutine adv_df(this, adv_options)
89  ! -- dummy
90  class(tspadvtype) :: this
91  type(tspadvoptionstype), optional, intent(in) :: adv_options !< the optional options, for when not constructing from file
92  ! -- local
93  character(len=*), parameter :: fmtadv = &
94  "(1x,/1x,'ADV -- ADVECTION PACKAGE, VERSION 1, 8/25/2017', &
95  &' INPUT READ FROM MEMPATH: ', A, //)"
96  !
97  ! -- Read or set advection options
98  if (.not. present(adv_options)) then
99  !
100  ! --print a message identifying the advection package.
101  write (this%iout, fmtadv) this%input_mempath
102  !
103  ! --read options from file
104  call this%source_options()
105  else
106  !
107  ! --set options from input arg
108  this%iadvwt = adv_options%iAdvScheme
109  end if
110  end subroutine adv_df
111 
112  !> @brief Allocate and read method for package
113  !!
114  !! Method to allocate and read static data for the ADV package.
115  !<
116  subroutine adv_ar(this, dis, ibound)
117  ! -- modules
118  use simmodule, only: store_error
119  ! -- dummy
120  class(tspadvtype) :: this
121  class(disbasetype), pointer, intent(in) :: dis
122  integer(I4B), dimension(:), pointer, contiguous, intent(in) :: ibound
123  ! -- local
124  integer(I4B) :: iadvwt_value
125  class(igradienttype), allocatable :: gradient
126  ! -- adv pointers to arguments that were passed in
127  this%dis => dis
128  this%ibound => ibound
129  !
130  ! -- Create interpolation scheme
131  iadvwt_value = this%iadvwt ! Dereference iadvwt to work with case statement
132  select case (iadvwt_value)
133  case (adv_scheme_upstream)
134  this%face_interpolation = &
135  upstreamschemetype(this%dis, this%fmi)
136  case (adv_scheme_central)
137  this%face_interpolation = &
138  centraldifferenceschemetype(this%dis, this%fmi)
139  case (adv_scheme_tvd)
140  this%face_interpolation = &
141  tvdschemetype(this%dis, this%fmi, this%ibound)
142  case (adv_scheme_utvd)
143  gradient = leastsquaresgradienttype(this%dis)
144  this%gradient = cachedgradienttype(gradient, this%dis)
145  this%face_interpolation = &
146  utvdschemetype(this%dis, this%fmi, this%gradient)
147  case default
148  call store_error("Unknown advection scheme", terminate=.true.)
149  end select
150  end subroutine adv_ar
151 
152  !> @brief Calculate maximum time step length
153  !!
154  !! Return the largest time step that meets stability constraints
155  !<
156  subroutine adv_dt(this, dtmax, msg, thetam)
157  ! dummy
158  class(tspadvtype) :: this !< this instance
159  real(DP), intent(out) :: dtmax !< maximum allowable dt subject to stability constraint
160  character(len=*), intent(inout) :: msg !< package/cell dt constraint message
161  real(DP), dimension(:), intent(in) :: thetam !< porosity
162  ! local
163  integer(I4B) :: n
164  integer(I4B) :: m
165  integer(I4B) :: ipos
166  integer(I4B) :: nrmax
167  character(len=LINELENGTH) :: cellstr
168  real(DP) :: dt
169  real(DP) :: flowmax
170  real(DP) :: flowsumpos
171  real(DP) :: flowsumneg
172  real(DP) :: flownm
173  real(DP) :: cell_volume
174  dtmax = dnodata
175  nrmax = 0
176  msg = ''
177 
178  ! If ats_percel not specified by user, then return without making
179  ! the courant time step calculation
180  if (this%ats_percel == dnodata) then
181  return
182  end if
183 
184  ! Calculate time step lengths based on stability constraint for each cell
185  ! and store the smallest one
186  do n = 1, this%dis%nodes
187  if (this%ibound(n) == 0) cycle
188  flowsumneg = dzero
189  flowsumpos = dzero
190  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
191  if (this%dis%con%mask(ipos) == 0) cycle
192  m = this%dis%con%ja(ipos)
193  if (this%ibound(m) == 0) cycle
194  flownm = this%fmi%gwfflowja(ipos)
195  if (flownm < dzero) then
196  flowsumneg = flowsumneg - flownm
197  else
198  flowsumpos = flowsumpos + flownm
199  end if
200  end do
201  flowmax = max(flowsumneg, flowsumpos)
202  if (flowmax < dprec) cycle
203  cell_volume = this%dis%get_cell_volume(n, this%dis%top(n))
204  dt = cell_volume * this%fmi%gwfsat(n) * thetam(n) / flowmax
205  dt = dt * this%ats_percel
206  if (dt < dtmax) then
207  dtmax = dt
208  nrmax = n
209  end if
210  end do
211  if (nrmax > 0) then
212  call this%dis%noder_to_string(nrmax, cellstr)
213  write (msg, *) adjustl(trim(this%memoryPath))//'-'//trim(cellstr)
214  end if
215  end subroutine adv_dt
216 
217  !> @brief Fill coefficient method for ADV package
218  !!
219  !! Method to calculate coefficients and fill amat and rhs.
220  !<
221  subroutine adv_fc(this, nodes, matrix_sln, idxglo, cnew, rhs)
222  ! -- modules
223  ! -- dummy
224  class(tspadvtype) :: this !< this instance
225  integer(I4B), intent(in) :: nodes !< number of nodes
226  class(matrixbasetype), pointer :: matrix_sln !< pointer to solution matrix
227  integer(I4B), intent(in), dimension(:) :: idxglo !< global indices for matrix
228  real(DP), intent(in), dimension(:), target :: cnew !< new concentration/temperature values
229  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector
230  ! -- local
231  integer(I4B) :: n, m, idiag, ipos
232  real(DP) :: qnm !< volumetric flow rate
233  type(coefficientstype) :: coefficients
234 
235  ! Calculate internal domain fluxes and add to matrix_sln and rhs.
236  call this%face_interpolation%set_field(cnew)
237  do n = 1, nodes
238  if (this%ibound(n) == 0) cycle ! skip inactive nodes
239  idiag = this%dis%con%ia(n)
240  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
241  m = this%dis%con%ja(ipos)
242  if (this%ibound(m) == 0) cycle ! skip inactive nodes
243  if (this%dis%con%mask(ipos) == 0) cycle ! skip masked connections
244 
245  qnm = this%fmi%gwfflowja(ipos) * this%eqnsclfac
246  coefficients = this%face_interpolation%compute(n, m, ipos)
247 
248  call matrix_sln%add_value_pos(idxglo(idiag), qnm * coefficients%c_n)
249  call matrix_sln%add_value_pos(idxglo(ipos), qnm * coefficients%c_m)
250  rhs(n) = rhs(n) + qnm * coefficients%rhs
251  end do
252  end do
253  end subroutine adv_fc
254 
255  !> @brief Calculate advection contribution to flowja
256  !<
257  subroutine adv_cq(this, cnew, flowja)
258  ! -- modules
259  ! -- dummy
260  class(tspadvtype) :: this
261  real(DP), intent(in), dimension(:), target :: cnew
262  real(DP), intent(inout), dimension(:) :: flowja
263  ! -- local
264  integer(I4B) :: nodes
265  integer(I4B) :: n, m, ipos
266  real(DP) :: qnm
267  type(coefficientstype) :: coefficients
268  !
269  ! -- Calculate advection and add to flowja. qnm is the volumetric flow
270  ! rate and has dimensions of L^/T.
271  nodes = this%dis%nodes
272 
273  call this%face_interpolation%set_field(cnew)
274  do n = 1, nodes
275  if (this%ibound(n) == 0) cycle
276  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
277  m = this%dis%con%ja(ipos)
278  if (this%ibound(m) == 0) cycle
279  qnm = this%fmi%gwfflowja(ipos) * this%eqnsclfac
280 
281  coefficients = this%face_interpolation%compute(n, m, ipos)
282  flowja(ipos) = flowja(ipos) &
283  + qnm * coefficients%c_n * cnew(n) &
284  + qnm * coefficients%c_m * cnew(m) &
285  - qnm * coefficients%rhs
286  end do
287  end do
288 
289  end subroutine adv_cq
290 
291  !> @brief Deallocate memory
292  !<
293  subroutine adv_da(this)
294  ! -- modules
296  ! -- dummy
297  class(tspadvtype) :: this
298  !
299  ! -- Deallocate arrays if package was active
300  if (this%inunit > 0) then
301  end if
302  !
303  ! -- nullify pointers
304  this%ibound => null()
305  !
306  ! -- Scalars
307  call mem_deallocate(this%iadvwt)
308  call mem_deallocate(this%ats_percel)
309  !
310  ! -- deallocate parent
311  call this%NumericalPackageType%da()
312  end subroutine adv_da
313 
314  !> @brief Allocate scalars specific to the streamflow energy transport (SFE)
315  !! package.
316  !<
317  subroutine allocate_scalars(this)
318  ! -- modules
320  ! -- dummy
321  class(tspadvtype) :: this
322  ! -- local
323  !
324  ! -- allocate scalars in NumericalPackageType
325  call this%NumericalPackageType%allocate_scalars()
326  !
327  ! -- Allocate
328  call mem_allocate(this%iadvwt, 'IADVWT', this%memoryPath)
329  call mem_allocate(this%ats_percel, 'ATS_PERCEL', this%memoryPath)
330  !
331  ! -- Initialize
332  this%iadvwt = adv_scheme_upstream
333  this%ats_percel = dnodata
334  !
335  ! -- Advection creates an asymmetric coefficient matrix
336  this%iasym = 1
337  end subroutine allocate_scalars
338 
339  !> @brief Source input options
340  !<
341  subroutine source_options(this)
342  ! -- modules
343  use kindmodule, only: lgp
344  use constantsmodule, only: lenvarname
345  use simvariablesmodule, only: errmsg
348  ! -- dummy
349  class(tspadvtype) :: this
350  ! -- locals
351  character(len=LENVARNAME), dimension(4) :: scheme = &
352  &[character(len=LENVARNAME) :: 'UPSTREAM', 'CENTRAL', 'TVD', 'UTVD']
353  logical(lgp) :: found_scheme, found_atspercel
354  ! -- formats
355  character(len=*), parameter :: fmtiadvwt = &
356  &"(4x,'ADVECTION WEIGHTING SCHEME HAS BEEN SET TO: ', a)"
357 
358  ! update defaults with input sourced values
359  call mem_set_value(this%iadvwt, 'SCHEME', this%input_mempath, &
360  scheme, found_scheme)
361  call mem_set_value(this%ats_percel, 'ATS_PERCEL', this%input_mempath, &
362  found_atspercel)
363 
364  if (found_scheme) then
365  ! should currently be set to index of scheme names
366  if (this%iadvwt == 0) then
367  write (errmsg, '(a, a)') &
368  'Unknown scheme, must be "UPSTREAM", "CENTRAL", "TVD" or "UTVD"'
369  call store_error(errmsg)
370  call store_error_filename(this%input_fname)
371  else
372  ! scheme parameters are 0 based
373  this%iadvwt = this%iadvwt - 1
374  end if
375  end if
376 
377  if (found_atspercel) then
378  if (this%ats_percel == dzero) this%ats_percel = dnodata
379  end if
380 
381  ! log options
382  write (this%iout, '(1x,a)') 'PROCESSING ADVECTION OPTIONS'
383  if (found_scheme) then
384  write (this%iout, fmtiadvwt) trim(scheme(this%iadvwt + 1))
385  end if
386  if (found_atspercel) then
387  write (this%iout, '(4x,a,1pg15.6)') &
388  'User-specified fractional cell distance for adaptive time &
389  &steps: ', this%ats_percel
390  end if
391  write (this%iout, '(1x,a)') 'END OF ADVECTION OPTIONS'
392  end subroutine source_options
393 
394 end module tspadvmodule
integer(i4b), parameter adv_scheme_tvd
integer(i4b), parameter adv_scheme_upstream
integer(i4b), parameter adv_scheme_utvd
integer(i4b), parameter adv_scheme_central
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:95
integer(i4b), parameter lenvarname
maximum length of a variable name
Definition: Constants.f90:17
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 done
real constant 1
Definition: Constants.f90:76
This module defines variable data types.
Definition: kind.f90:8
This module contains the base numerical package type.
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
subroutine, public adv_cr(advobj, name_model, input_mempath, inunit, iout, fmi, eqnsclfac)
@ brief Create a new ADV object
Definition: tsp-adv.f90:59
subroutine source_options(this)
Source input options.
Definition: tsp-adv.f90:342
subroutine adv_df(this, adv_options)
Define ADV object.
Definition: tsp-adv.f90:89
subroutine adv_fc(this, nodes, matrix_sln, idxglo, cnew, rhs)
Fill coefficient method for ADV package.
Definition: tsp-adv.f90:222
subroutine adv_cq(this, cnew, flowja)
Calculate advection contribution to flowja.
Definition: tsp-adv.f90:258
subroutine adv_da(this)
Deallocate memory.
Definition: tsp-adv.f90:294
subroutine adv_ar(this, dis, ibound)
Allocate and read method for package.
Definition: tsp-adv.f90:117
subroutine adv_dt(this, dtmax, msg, thetam)
Calculate maximum time step length.
Definition: tsp-adv.f90:157
subroutine allocate_scalars(this)
Allocate scalars specific to the streamflow energy transport (SFE) package.
Definition: tsp-adv.f90:318
Decorator that adds caching to any gradient computation implementation.
Abstract interface for cell-based gradient computation.
Definition: IGradient.f90:15
Weighted least-squares gradient method for structured and unstructured grids.
Total Variation Diminishing (TVD) interpolation scheme.
Definition: UTVDScheme.f90:37