MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
tsp-adv.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b
4  use constantsmodule, only: done, dzero, dhalf, dtwo
6  use basedismodule, only: disbasetype
7  use tspfmimodule, only: tspfmitype
10 
11  implicit none
12  private
13  public :: tspadvtype
14  public :: adv_cr
15 
16  type, extends(numericalpackagetype) :: tspadvtype
17 
18  integer(I4B), pointer :: iadvwt => null() !< advection scheme (0 up, 1 central, 2 tvd)
19  integer(I4B), dimension(:), pointer, contiguous :: ibound => null() !< pointer to model ibound
20  type(tspfmitype), pointer :: fmi => null() !< pointer to fmi object
21  real(dp), pointer :: eqnsclfac => null() !< governing equation scale factor; =1. for solute; =rhow*cpw for energy
22 
23  contains
24 
25  procedure :: adv_df
26  procedure :: adv_ar
27  procedure :: adv_fc
28  procedure :: adv_cq
29  procedure :: adv_da
30 
31  procedure :: allocate_scalars
32  procedure, private :: read_options
33  procedure, private :: advqtvd
34  procedure, private :: advtvd_bd
35  procedure :: adv_weight
36  procedure :: advtvd
37 
38  end type tspadvtype
39 
40 contains
41 
42  !> @ brief Create a new ADV object
43  !!
44  !! Create a new ADV package
45  !<
46  subroutine adv_cr(advobj, name_model, inunit, iout, fmi, eqnsclfac)
47  ! -- dummy
48  type(tspadvtype), pointer :: advobj
49  character(len=*), intent(in) :: name_model
50  integer(I4B), intent(in) :: inunit
51  integer(I4B), intent(in) :: iout
52  type(tspfmitype), intent(in), target :: fmi
53  real(dp), intent(in), pointer :: eqnsclfac !< governing equation scale factor
54  !
55  ! -- Create the object
56  allocate (advobj)
57  !
58  ! -- create name and memory path
59  call advobj%set_names(1, name_model, 'ADV', 'ADV')
60  !
61  ! -- Allocate scalars
62  call advobj%allocate_scalars()
63  !
64  ! -- Set variables
65  advobj%inunit = inunit
66  advobj%iout = iout
67  advobj%fmi => fmi
68  advobj%eqnsclfac => eqnsclfac
69  !
70  ! -- Return
71  return
72  end subroutine adv_cr
73 
74  !> @brief Define ADV object
75  !!
76  !! Define the ADV package
77  !<
78  subroutine adv_df(this, adv_options)
79  ! -- dummy
80  class(tspadvtype) :: this
81  type(tspadvoptionstype), optional, intent(in) :: adv_options !< the optional options, for when not constructing from file
82  ! -- local
83  character(len=*), parameter :: fmtadv = &
84  "(1x,/1x,'ADV-- ADVECTION PACKAGE, VERSION 1, 8/25/2017', &
85  &' INPUT READ FROM UNIT ', i0, //)"
86  !
87  ! -- Read or set advection options
88  if (.not. present(adv_options)) then
89  !
90  ! -- Initialize block parser (adv has no define, so it's
91  ! not done until here)
92  call this%parser%Initialize(this%inunit, this%iout)
93  !
94  ! --print a message identifying the advection package.
95  write (this%iout, fmtadv) this%inunit
96  !
97  ! --read options from file
98  call this%read_options()
99  else
100  !
101  ! --set options from input arg
102  this%iadvwt = adv_options%iAdvScheme
103  end if
104  !
105  ! -- Return
106  return
107  end subroutine adv_df
108 
109  !> @brief Allocate and read method for package
110  !!
111  !! Method to allocate and read static data for the ADV package.
112  !<
113  subroutine adv_ar(this, dis, ibound)
114  ! -- modules
115  ! -- dummy
116  class(tspadvtype) :: this
117  class(disbasetype), pointer, intent(in) :: dis
118  integer(I4B), dimension(:), pointer, contiguous, intent(in) :: ibound
119  ! -- local
120  ! -- formats
121  !
122  ! -- adv pointers to arguments that were passed in
123  this%dis => dis
124  this%ibound => ibound
125  !
126  ! -- Return
127  return
128  end subroutine adv_ar
129 
130  !> @brief Fill coefficient method for ADV package
131  !!
132  !! Method to calculate coefficients and fill amat and rhs.
133  !<
134  subroutine adv_fc(this, nodes, matrix_sln, idxglo, cnew, rhs)
135  ! -- modules
136  ! -- dummy
137  class(tspadvtype) :: this
138  integer, intent(in) :: nodes
139  class(matrixbasetype), pointer :: matrix_sln
140  integer(I4B), intent(in), dimension(:) :: idxglo
141  real(DP), intent(in), dimension(:) :: cnew
142  real(DP), dimension(:), intent(inout) :: rhs
143  ! -- local
144  integer(I4B) :: n, m, idiag, ipos
145  real(DP) :: omega, qnm
146  !
147  ! -- Calculate advection terms and add to solution rhs and hcof. qnm
148  ! is the volumetric flow rate and has dimensions of L^/T.
149  do n = 1, nodes
150  if (this%ibound(n) == 0) cycle
151  idiag = this%dis%con%ia(n)
152  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
153  if (this%dis%con%mask(ipos) == 0) cycle
154  m = this%dis%con%ja(ipos)
155  if (this%ibound(m) == 0) cycle
156  qnm = this%fmi%gwfflowja(ipos) * this%eqnsclfac
157  omega = this%adv_weight(this%iadvwt, ipos, n, m, qnm)
158  call matrix_sln%add_value_pos(idxglo(ipos), qnm * (done - omega))
159  call matrix_sln%add_value_pos(idxglo(idiag), qnm * omega)
160  end do
161  end do
162  !
163  ! -- TVD
164  if (this%iadvwt == 2) then
165  do n = 1, nodes
166  if (this%ibound(n) == 0) cycle
167  call this%advtvd(n, cnew, rhs)
168  end do
169  end if
170  !
171  ! -- Return
172  return
173  end subroutine adv_fc
174 
175  !> @brief Calculate TVD
176  !!
177  !! Use explicit scheme to calculate the advective component of transport.
178  !! TVD is an acronym for Total-Variation Diminishing
179  !<
180  subroutine advtvd(this, n, cnew, rhs)
181  ! -- modules
182  ! -- dummy
183  class(tspadvtype) :: this
184  integer(I4B), intent(in) :: n
185  real(DP), dimension(:), intent(in) :: cnew
186  real(DP), dimension(:), intent(inout) :: rhs
187  ! -- local
188  real(DP) :: qtvd
189  integer(I4B) :: m, ipos
190  !
191  ! -- Loop through each n connection. This will
192  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
193  if (this%dis%con%mask(ipos) == 0) cycle
194  m = this%dis%con%ja(ipos)
195  if (m > n .and. this%ibound(m) /= 0) then
196  qtvd = this%advqtvd(n, m, ipos, cnew)
197  rhs(n) = rhs(n) - qtvd
198  rhs(m) = rhs(m) + qtvd
199  end if
200  end do
201  !
202  ! -- Return
203  return
204  end subroutine advtvd
205 
206  !> @brief Calculate TVD
207  !!
208  !! Use explicit scheme to calculate the advective component of transport.
209  !! TVD is an acronym for Total-Variation Diminishing
210  !<
211  function advqtvd(this, n, m, iposnm, cnew) result(qtvd)
212  ! -- modules
213  use constantsmodule, only: dprec
214  ! -- return
215  real(dp) :: qtvd
216  ! -- dummy
217  class(tspadvtype) :: this
218  integer(I4B), intent(in) :: n
219  integer(I4B), intent(in) :: m
220  integer(I4B), intent(in) :: iposnm
221  real(dp), dimension(:), intent(in) :: cnew
222  ! -- local
223  integer(I4B) :: ipos, isympos, iup, idn, i2up, j
224  real(dp) :: qnm, qmax, qupj, elupdn, elup2up
225  real(dp) :: smooth, cdiff, alimiter
226  !
227  ! -- initialize
228  qtvd = dzero
229  !
230  ! -- Find upstream node
231  isympos = this%dis%con%jas(iposnm)
232  qnm = this%fmi%gwfflowja(iposnm)
233  if (qnm > dzero) then
234  ! -- positive flow into n means m is upstream
235  iup = m
236  idn = n
237  else
238  iup = n
239  idn = m
240  end if
241  elupdn = this%dis%con%cl1(isympos) + this%dis%con%cl2(isympos)
242  !
243  ! -- Find second node upstream to iup
244  i2up = 0
245  qmax = dzero
246  do ipos = this%dis%con%ia(iup) + 1, this%dis%con%ia(iup + 1) - 1
247  j = this%dis%con%ja(ipos)
248  if (this%ibound(j) == 0) cycle
249  qupj = this%fmi%gwfflowja(ipos)
250  isympos = this%dis%con%jas(ipos)
251  if (qupj > qmax) then
252  qmax = qupj
253  i2up = j
254  elup2up = this%dis%con%cl1(isympos) + this%dis%con%cl2(isympos)
255  end if
256  end do
257  !
258  ! -- Calculate flux limiting term
259  if (i2up > 0) then
260  smooth = dzero
261  cdiff = abs(cnew(idn) - cnew(iup))
262  if (cdiff > dprec) then
263  smooth = (cnew(iup) - cnew(i2up)) / elup2up * &
264  elupdn / (cnew(idn) - cnew(iup))
265  end if
266  if (smooth > dzero) then
267  alimiter = dtwo * smooth / (done + smooth)
268  qtvd = dhalf * alimiter * qnm * (cnew(idn) - cnew(iup))
269  qtvd = qtvd * this%eqnsclfac
270  end if
271  end if
272  !
273  ! -- Return
274  return
275  end function advqtvd
276 
277  !> @brief Calculate advection contribution to flowja
278  !<
279  subroutine adv_cq(this, cnew, flowja)
280  ! -- modules
281  ! -- dummy
282  class(tspadvtype) :: this
283  real(DP), intent(in), dimension(:) :: cnew
284  real(DP), intent(inout), dimension(:) :: flowja
285  ! -- local
286  integer(I4B) :: nodes
287  integer(I4B) :: n, m, idiag, ipos
288  real(DP) :: omega, qnm
289  !
290  ! -- Calculate advection and add to flowja. qnm is the volumetric flow
291  ! rate and has dimensions of L^/T.
292  nodes = this%dis%nodes
293  do n = 1, nodes
294  if (this%ibound(n) == 0) cycle
295  idiag = this%dis%con%ia(n)
296  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
297  m = this%dis%con%ja(ipos)
298  if (this%ibound(m) == 0) cycle
299  qnm = this%fmi%gwfflowja(ipos) * this%eqnsclfac
300  omega = this%adv_weight(this%iadvwt, ipos, n, m, qnm)
301  flowja(ipos) = flowja(ipos) + qnm * omega * cnew(n) + &
302  qnm * (done - omega) * cnew(m)
303  end do
304  end do
305  !
306  ! -- TVD
307  if (this%iadvwt == 2) call this%advtvd_bd(cnew, flowja)
308  !
309  ! -- Return
310  return
311  end subroutine adv_cq
312 
313  !> @brief Add TVD contribution to flowja
314  subroutine advtvd_bd(this, cnew, flowja)
315  ! -- modules
316  ! -- dummy
317  class(tspadvtype) :: this
318  real(DP), dimension(:), intent(in) :: cnew
319  real(DP), dimension(:), intent(inout) :: flowja
320  ! -- local
321  real(DP) :: qtvd, qnm
322  integer(I4B) :: nodes, n, m, ipos
323  !
324  nodes = this%dis%nodes
325  do n = 1, nodes
326  if (this%ibound(n) == 0) cycle
327  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
328  m = this%dis%con%ja(ipos)
329  if (this%ibound(m) /= 0) then
330  qnm = this%fmi%gwfflowja(ipos)
331  qtvd = this%advqtvd(n, m, ipos, cnew)
332  flowja(ipos) = flowja(ipos) + qtvd
333  end if
334  end do
335  end do
336  !
337  ! -- Return
338  return
339  end subroutine advtvd_bd
340 
341  !> @brief Deallocate memory
342  !<
343  subroutine adv_da(this)
344  ! -- modules
346  ! -- dummy
347  class(tspadvtype) :: this
348  !
349  ! -- Deallocate arrays if package was active
350  if (this%inunit > 0) then
351  end if
352  !
353  ! -- nullify pointers
354  this%ibound => null()
355  !
356  ! -- Scalars
357  call mem_deallocate(this%iadvwt)
358  !
359  ! -- deallocate parent
360  call this%NumericalPackageType%da()
361  !
362  ! -- Return
363  return
364  end subroutine adv_da
365 
366  !> @brief Allocate scalars specific to the streamflow energy transport (SFE)
367  !! package.
368  !<
369  subroutine allocate_scalars(this)
370  ! -- modules
372  ! -- dummy
373  class(tspadvtype) :: this
374  ! -- local
375 ! ------------------------------------------------------------------------------
376  !
377  ! -- allocate scalars in NumericalPackageType
378  call this%NumericalPackageType%allocate_scalars()
379  !
380  ! -- Allocate
381  call mem_allocate(this%iadvwt, 'IADVWT', this%memoryPath)
382  !
383  ! -- Initialize
384  this%iadvwt = 0
385  !
386  ! -- Advection creates an asymmetric coefficient matrix
387  this%iasym = 1
388  !
389  ! -- Return
390  return
391  end subroutine allocate_scalars
392 
393  !> @brief Read options
394  !!
395  !! Read the options block
396  !<
397  subroutine read_options(this)
398  ! -- modules
399  use constantsmodule, only: linelength
400  use simmodule, only: store_error
401  ! -- dummy
402  class(tspadvtype) :: this
403  ! -- local
404  character(len=LINELENGTH) :: errmsg, keyword
405  integer(I4B) :: ierr
406  logical :: isfound, endOfBlock
407  ! -- formats
408  character(len=*), parameter :: fmtiadvwt = &
409  &"(4x,'ADVECTION WEIGHTING SCHEME HAS BEEN SET TO: ', a)"
410  !
411  ! -- get options block
412  call this%parser%GetBlock('OPTIONS', isfound, ierr, blockrequired=.false., &
413  supportopenclose=.true.)
414  !
415  ! -- parse options block if detected
416  if (isfound) then
417  write (this%iout, '(1x,a)') 'PROCESSING ADVECTION OPTIONS'
418  do
419  call this%parser%GetNextLine(endofblock)
420  if (endofblock) exit
421  call this%parser%GetStringCaps(keyword)
422  select case (keyword)
423  case ('SCHEME')
424  call this%parser%GetStringCaps(keyword)
425  select case (keyword)
426  case ('UPSTREAM')
427  this%iadvwt = 0
428  write (this%iout, fmtiadvwt) 'UPSTREAM'
429  case ('CENTRAL')
430  this%iadvwt = 1
431  write (this%iout, fmtiadvwt) 'CENTRAL'
432  case ('TVD')
433  this%iadvwt = 2
434  write (this%iout, fmtiadvwt) 'TVD'
435  case default
436  write (errmsg, '(a, a)') &
437  'Unknown scheme: ', trim(keyword)
438  call store_error(errmsg)
439  write (errmsg, '(a, a)') &
440  'Scheme must be "UPSTREAM", "CENTRAL" or "TVD"'
441  call store_error(errmsg)
442  call this%parser%StoreErrorUnit()
443  end select
444  case default
445  write (errmsg, '(a,a)') 'Unknown ADVECTION option: ', &
446  trim(keyword)
447  call store_error(errmsg, terminate=.true.)
448  end select
449  end do
450  write (this%iout, '(1x,a)') 'END OF ADVECTION OPTIONS'
451  end if
452  !
453  ! -- Return
454  return
455  end subroutine read_options
456 
457  !> @ brief Advection weight
458  !!
459  !! Calculate the advection weight
460  !<
461  function adv_weight(this, iadvwt, ipos, n, m, qnm) result(omega)
462  ! -- return
463  real(dp) :: omega
464  ! -- dummy
465  class(tspadvtype) :: this
466  integer, intent(in) :: iadvwt
467  integer, intent(in) :: ipos
468  integer, intent(in) :: n
469  integer, intent(in) :: m
470  real(dp), intent(in) :: qnm
471  ! -- local
472  real(dp) :: lnm, lmn
473 ! ------------------------------------------------------------------------------
474  select case (iadvwt)
475  case (1)
476  ! -- calculate weight based on distances between nodes and the shared
477  ! face of the connection
478  if (this%dis%con%ihc(this%dis%con%jas(ipos)) == 0) then
479  ! -- vertical connection; assume cell is fully saturated
480  lnm = dhalf * (this%dis%top(n) - this%dis%bot(n))
481  lmn = dhalf * (this%dis%top(m) - this%dis%bot(m))
482  else
483  ! -- horizontal connection
484  lnm = this%dis%con%cl1(this%dis%con%jas(ipos))
485  lmn = this%dis%con%cl2(this%dis%con%jas(ipos))
486  end if
487  omega = lmn / (lnm + lmn)
488  case (0, 2)
489  ! -- use upstream weighting for upstream and tvd schemes
490  if (qnm > dzero) then
491  omega = dzero
492  else
493  omega = done
494  end if
495  end select
496  !
497  ! -- return
498  return
499  end function adv_weight
500 
501 end module tspadvmodule
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:67
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
real(dp), parameter dprec
real constant machine precision
Definition: Constants.f90:119
real(dp), parameter dtwo
real constant 2
Definition: Constants.f90:78
real(dp), parameter done
real constant 1
Definition: Constants.f90:75
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 advtvd_bd(this, cnew, flowja)
Add TVD contribution to flowja.
Definition: tsp-adv.f90:315
real(dp) function advqtvd(this, n, m, iposnm, cnew)
Calculate TVD.
Definition: tsp-adv.f90:212
subroutine adv_df(this, adv_options)
Define ADV object.
Definition: tsp-adv.f90:79
subroutine adv_fc(this, nodes, matrix_sln, idxglo, cnew, rhs)
Fill coefficient method for ADV package.
Definition: tsp-adv.f90:135
subroutine adv_cq(this, cnew, flowja)
Calculate advection contribution to flowja.
Definition: tsp-adv.f90:280
subroutine adv_da(this)
Deallocate memory.
Definition: tsp-adv.f90:344
subroutine advtvd(this, n, cnew, rhs)
Calculate TVD.
Definition: tsp-adv.f90:181
subroutine adv_ar(this, dis, ibound)
Allocate and read method for package.
Definition: tsp-adv.f90:114
real(dp) function adv_weight(this, iadvwt, ipos, n, m, qnm)
@ brief Advection weight
Definition: tsp-adv.f90:462
subroutine read_options(this)
Read options.
Definition: tsp-adv.f90:398
subroutine, public adv_cr(advobj, name_model, inunit, iout, fmi, eqnsclfac)
@ brief Create a new ADV object
Definition: tsp-adv.f90:47
subroutine allocate_scalars(this)
Allocate scalars specific to the streamflow energy transport (SFE) package.
Definition: tsp-adv.f90:370