MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
swf-zdg.f90
Go to the documentation of this file.
1 !> @brief This module contains the ZDG package methods
2 !!
3 !! This module can be used to represent outflow from streams using
4 !! a zero-depth-gradient boundary.
5 !!
6 !<
8  ! -- modules
9  use kindmodule, only: dp, i4b
11  use simvariablesmodule, only: errmsg
14  use bndmodule, only: bndtype
15  use bndextmodule, only: bndexttype
17  use observemodule, only: observetype
21  use basedismodule, only: disbasetype
22  use disv1dmodule, only: disv1dtype
23  use swfcxsmodule, only: swfcxstype
24  !
25  implicit none
26  !
27  private
28  public :: zdg_create
29  !
30  character(len=LENFTYPE) :: ftype = 'ZDG' !< package ftype
31  character(len=16) :: text = ' ZDG' !< package flow text string
32  !
33  type, extends(bndexttype) :: swfzdgtype
34 
35  integer(I4B), dimension(:), pointer, contiguous :: idcxs => null() !< cross section id
36  real(dp), dimension(:), pointer, contiguous :: width => null() !< channel width
37  real(dp), dimension(:), pointer, contiguous :: slope => null() !< channel slope
38  real(dp), dimension(:), pointer, contiguous :: rough => null() !< channel roughness
39  real(dp), pointer :: unitconv => null() !< conversion factor for roughness to length and time units of meters and seconds
40 
41  ! -- pointers other objects
42  type(disv1dtype), pointer :: disv1d
43  type(swfcxstype), pointer :: cxs
44 
45  contains
46  procedure :: allocate_scalars => zdg_allocate_scalars
47  procedure :: allocate_arrays => zdg_allocate_arrays
48  procedure :: source_options => zdg_options
49  procedure :: log_zdg_options
50  procedure :: bnd_cf => zdg_cf
51  procedure :: bnd_fc => zdg_fc
52  procedure :: bnd_da => zdg_da
53  procedure :: define_listlabel
54  procedure :: bound_value => zdg_bound_value
55  ! -- methods for observations
56  procedure, public :: bnd_obs_supported => zdg_obs_supported
57  procedure, public :: bnd_df_obs => zdg_df_obs
58  procedure, public :: bnd_bd_obs => zdg_bd_obs
59  ! -- methods for time series
60  procedure, public :: bnd_rp_ts => zdg_rp_ts
61  ! -- private
62  procedure, private :: qcalc
63  end type swfzdgtype
64 
65 contains
66 
67  !> @ brief Create a new package object
68  !!
69  !! Create a new ZDG Package object
70  !!
71  !<
72  subroutine zdg_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
73  mempath, dis, cxs, unitconv)
74  ! -- dummy variables
75  class(bndtype), pointer :: packobj !< pointer to default package type
76  integer(I4B), intent(in) :: id !< package id
77  integer(I4B), intent(in) :: ibcnum !< boundary condition number
78  integer(I4B), intent(in) :: inunit !< unit number of ZDG package input file
79  integer(I4B), intent(in) :: iout !< unit number of model listing file
80  character(len=*), intent(in) :: namemodel !< model name
81  character(len=*), intent(in) :: pakname !< package name
82  character(len=*), intent(in) :: mempath !< input mempath
83  class(disbasetype), pointer, intent(inout) :: dis !< the pointer to the discretization
84  type(swfcxstype), pointer, intent(in) :: cxs !< the pointer to the cxs package
85  real(dp), intent(in) :: unitconv !< unit conversion for roughness
86  ! -- local variables
87  type(swfzdgtype), pointer :: zdgobj
88  !
89  ! -- allocate the object and assign values to object variables
90  allocate (zdgobj)
91  packobj => zdgobj
92  !
93  ! -- create name and memory path
94  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
95  packobj%text = text
96  !
97  ! -- allocate scalars
98  call zdgobj%allocate_scalars()
99  !
100  ! -- initialize package
101  call packobj%pack_initialize()
102 
103  packobj%inunit = inunit
104  packobj%iout = iout
105  packobj%id = id
106  packobj%ibcnum = ibcnum
107  packobj%ncolbnd = 1
108  packobj%iscloc = 1
109  packobj%ictMemPath = create_mem_path(namemodel, 'DFW')
110  !
111  ! -- store pointer to disv1d
112  select type (dis)
113  type is (disv1dtype)
114  zdgobj%disv1d => dis
115  end select
116  !
117  ! -- store pointer to cxs
118  zdgobj%cxs => cxs
119  !
120  ! -- store unit conversion
121  zdgobj%unitconv = unitconv
122  end subroutine zdg_create
123 
124  !> @ brief Allocate scalars
125  !!
126  !! Allocate and initialize scalars for the ZDG package. The base model
127  !! allocate scalars method is also called.
128  !!
129  !<
130  subroutine zdg_allocate_scalars(this)
131  ! -- modules
133  ! -- dummy variables
134  class(swfzdgtype) :: this !< SwfZdgType object
135  !
136  ! -- call base type allocate scalars
137  call this%BndExtType%allocate_scalars()
138  !
139  ! -- allocate the object and assign values to object variables
140  call mem_allocate(this%unitconv, 'UNITCONV', this%memoryPath)
141  !
142  ! -- Set values
143  this%unitconv = dzero
144  end subroutine zdg_allocate_scalars
145 
146  !> @ brief Allocate arrays
147  !!
148  !! Allocate and initialize arrays for the SWF package
149  !!
150  !<
151  subroutine zdg_allocate_arrays(this, nodelist, auxvar)
152  ! -- modules
154  ! -- dummy
155  class(swfzdgtype) :: this
156  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
157  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
158  ! -- local
159  !
160  ! -- call BndExtType allocate scalars
161  call this%BndExtType%allocate_arrays(nodelist, auxvar)
162  !
163  ! -- set array input context pointer
164  call mem_setptr(this%idcxs, 'IDCXS', this%input_mempath)
165  call mem_setptr(this%width, 'WIDTH', this%input_mempath)
166  call mem_setptr(this%slope, 'SLOPE', this%input_mempath)
167  call mem_setptr(this%rough, 'ROUGH', this%input_mempath)
168  !
169  ! -- checkin array input context pointer
170  call mem_checkin(this%idcxs, 'IDCXS', this%memoryPath, &
171  'IDCXS', this%input_mempath)
172  call mem_checkin(this%width, 'WIDTH', this%memoryPath, &
173  'WIDTH', this%input_mempath)
174  call mem_checkin(this%slope, 'SLOPE', this%memoryPath, &
175  'SLOPE', this%input_mempath)
176  call mem_checkin(this%rough, 'ROUGH', this%memoryPath, &
177  'ROUGH', this%input_mempath)
178  end subroutine zdg_allocate_arrays
179 
180  !> @ brief Deallocate package memory
181  !!
182  !! Deallocate SWF package scalars and arrays.
183  !!
184  !<
185  subroutine zdg_da(this)
186  ! -- modules
188  ! -- dummy variables
189  class(swfzdgtype) :: this !< SwfZdgType object
190  !
191  ! -- Deallocate parent package
192  call this%BndExtType%bnd_da()
193  !
194  ! -- arrays
195  call mem_deallocate(this%idcxs, 'IDCXS', this%memoryPath)
196  call mem_deallocate(this%width, 'WIDTH', this%memoryPath)
197  call mem_deallocate(this%slope, 'SLOPE', this%memoryPath)
198  call mem_deallocate(this%rough, 'ROUGH', this%memoryPath)
199  !
200  ! -- scalars
201  call mem_deallocate(this%unitconv)
202  end subroutine zdg_da
203 
204  !> @ brief Source additional options for package
205  !!
206  !! Source additional options for SWF package.
207  !!
208  !<
209  subroutine zdg_options(this)
210  ! -- modules
211  use inputoutputmodule, only: urword
214  ! -- dummy variables
215  class(swfzdgtype), intent(inout) :: this !< SwfZdgType object
216  ! -- local variables
217  type(swfzdgparamfoundtype) :: found
218  ! -- formats
219  !
220  ! -- source base BndExtType options
221  call this%BndExtType%source_options()
222  !
223  ! -- source options from input context
224  ! none
225  !
226  ! -- log SWF specific options
227  call this%log_zdg_options(found)
228  end subroutine zdg_options
229 
230  !> @ brief Log SWF specific package options
231  !<
232  subroutine log_zdg_options(this, found)
233  ! -- modules
235  ! -- dummy variables
236  class(swfzdgtype), intent(inout) :: this
237  type(swfzdgparamfoundtype), intent(in) :: found
238  ! -- local variables
239  ! -- format
240  !
241  ! -- log found options
242  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
243  //' OPTIONS'
244  !
245  ! if (found%mover) then
246  ! write (this%iout, '(4x,A)') 'MOVER OPTION ENABLED'
247  ! end if
248  !
249  ! -- close logging block
250  write (this%iout, '(1x,a)') &
251  'END OF '//trim(adjustl(this%text))//' OPTIONS'
252  end subroutine log_zdg_options
253 
254  !> @ brief Formulate the package hcof and rhs terms.
255  !!
256  !! Formulate the hcof and rhs terms for the ZDG package that will be
257  !! added to the coefficient matrix and right-hand side vector.
258  !!
259  !<
260  subroutine zdg_cf(this)
261  ! modules
263  use smoothingmodule, only: squadratic
264  ! dummy variables
265  class(swfzdgtype) :: this !< SwfZdgType object
266  ! local variables
267  integer(I4B) :: i, node
268  real(DP) :: q
269  real(DP) :: qeps
270  real(DP) :: absdhdxsq
271  real(DP) :: depth
272  real(DP) :: derv
273  real(DP) :: eps
274  real(DP) :: range = 1.d-6
275  real(DP) :: dydx
276  real(DP) :: smooth_factor
277  !
278  ! -- Return if no inflows
279  if (this%nbound == 0) return
280  !
281  ! -- Calculate hcof and rhs for each zdg entry
282  do i = 1, this%nbound
283 
284  node = this%nodelist(i)
285  if (this%ibound(node) <= 0) then
286  this%hcof(i) = dzero
287  this%rhs(i) = dzero
288  cycle
289  end if
290 
291  ! -- calculate terms and add to hcof and rhs
292  absdhdxsq = this%slope(i)**dhalf
293  depth = this%xnew(node) - this%dis%bot(node)
294 
295  ! smooth the depth
296  call squadratic(depth, range, dydx, smooth_factor)
297  depth = depth * smooth_factor
298 
299  ! -- calculate unperturbed q
300  q = -this%qcalc(i, depth, this%unitconv)
301 
302  ! -- calculate perturbed q
303  eps = get_perturbation(depth)
304  qeps = -this%qcalc(i, depth + eps, this%unitconv)
305 
306  ! -- calculate derivative
307  derv = (qeps - q) / eps
308 
309  ! -- add terms to hcof and rhs
310  this%hcof(i) = derv
311  this%rhs(i) = -q + derv * this%xnew(node)
312 
313  end do
314  end subroutine zdg_cf
315 
316  ! !> @brief Calculate flow
317  !!
318  !! Calculate volumetric flow rate for the zero-depth gradient
319  !! condition. Flow is positive.
320  ! !<
321  function qcalc(this, i, depth, unitconv) result(q)
322  ! dummy
323  class(swfzdgtype) :: this
324  integer(I4B), intent(in) :: i !< boundary number
325  real(dp), intent(in) :: depth !< simulated depth (stage - elevation) in reach n for this iteration
326  real(dp), intent(in) :: unitconv !< conversion factor for roughness to length and time units of meters and seconds
327  ! return
328  real(dp) :: q
329  ! local
330  integer(I4B) :: idcxs
331  real(dp) :: width
332  real(dp) :: rough
333  real(dp) :: slope
334  real(dp) :: conveyance
335 
336  idcxs = this%idcxs(i)
337  width = this%width(i)
338  rough = this%rough(i)
339  slope = this%slope(i)
340  conveyance = this%cxs%get_conveyance(idcxs, width, depth, rough)
341  q = conveyance * slope**dhalf * unitconv
342 
343  end function qcalc
344 
345  !> @ brief Copy hcof and rhs terms into solution.
346  !!
347  !! Add the hcof and rhs terms for the ZDG package to the
348  !! coefficient matrix and right-hand side vector.
349  !!
350  !<
351  subroutine zdg_fc(this, rhs, ia, idxglo, matrix_sln)
352  ! -- dummy variables
353  class(swfzdgtype) :: this !< SwfZdgType object
354  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
355  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
356  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
357  class(matrixbasetype), pointer :: matrix_sln !< solution coefficient matrix
358  ! -- local variables
359  integer(I4B) :: i
360  integer(I4B) :: n
361  integer(I4B) :: ipos
362  !
363  ! -- pakmvrobj fc
364  if (this%imover == 1) then
365  call this%pakmvrobj%fc()
366  end if
367  !
368  ! -- Copy package rhs and hcof into solution rhs and amat
369  do i = 1, this%nbound
370  n = this%nodelist(i)
371  rhs(n) = rhs(n) + this%rhs(i)
372  ipos = ia(n)
373  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
374  !
375  ! -- If mover is active and this zdg item is discharging,
376  ! store available water (as positive value).
377  if (this%imover == 1 .and. this%rhs(i) > dzero) then
378  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
379  end if
380  end do
381  end subroutine zdg_fc
382 
383  !> @ brief Define the list label for the package
384  !!
385  !! Method defined the list label for the ZDG package. The list label is
386  !! the heading that is written to iout when PRINT_INPUT option is used.
387  !!
388  !<
389  subroutine define_listlabel(this)
390  ! -- dummy variables
391  class(swfzdgtype), intent(inout) :: this !< SwfZdgType object
392  !
393  ! -- create the header list label
394  this%listlabel = trim(this%filtyp)//' NO.'
395  if (this%dis%ndim == 3) then
396  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
397  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
398  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
399  elseif (this%dis%ndim == 2) then
400  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
401  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
402  else
403  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
404  end if
405  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'FLOW RATE'
406  if (this%inamedbound == 1) then
407  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
408  end if
409  end subroutine define_listlabel
410 
411  ! -- Procedures related to observations
412 
413  !> @brief Determine if observations are supported.
414  !!
415  !! Function to determine if observations are supported by the ZDG package.
416  !! Observations are supported by the ZDG package.
417  !!
418  !! @return zdg_obs_supported boolean indicating if observations are supported
419  !!
420  !<
421  logical function zdg_obs_supported(this)
422  ! -- dummy variables
423  class(swfzdgtype) :: this !< SwfZdgType object
424  !
425  ! -- set boolean
426  zdg_obs_supported = .true.
427  end function zdg_obs_supported
428 
429  !> @brief Define the observation types available in the package
430  !!
431  !! Method to define the observation types available in the ZDG package.
432  !!
433  !<
434  subroutine zdg_df_obs(this)
435  ! -- dummy variables
436  class(swfzdgtype) :: this !< SwfZdgType object
437  ! -- local variables
438  integer(I4B) :: indx
439  !
440  ! -- initialize observations
441  call this%obs%StoreObsType('zdg', .true., indx)
442  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
443  !
444  ! -- Store obs type and assign procedure pointer
445  ! for to-mvr observation type.
446  call this%obs%StoreObsType('to-mvr', .true., indx)
447  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
448  end subroutine zdg_df_obs
449 
450  !> @brief Save observations for the package
451  !!
452  !! Method to save simulated values for the ZDG package.
453  !!
454  !<
455  subroutine zdg_bd_obs(this)
456  ! -- dummy variables
457  class(swfzdgtype) :: this !< SwfZdgType object
458  ! -- local variables
459  integer(I4B) :: i
460  integer(I4B) :: n
461  integer(I4B) :: jj
462  real(DP) :: v
463  type(observetype), pointer :: obsrv => null()
464  !
465  ! -- clear the observations
466  call this%obs%obs_bd_clear()
467  !
468  ! -- Save simulated values for all of package's observations.
469  do i = 1, this%obs%npakobs
470  obsrv => this%obs%pakobs(i)%obsrv
471  if (obsrv%BndFound) then
472  do n = 1, obsrv%indxbnds_count
473  v = dnodata
474  jj = obsrv%indxbnds(n)
475  select case (obsrv%ObsTypeId)
476  case ('TO-MVR')
477  if (this%imover == 1) then
478  v = this%pakmvrobj%get_qtomvr(jj)
479  if (v > dzero) then
480  v = -v
481  end if
482  end if
483  case ('ZDG')
484  v = this%simvals(jj)
485  case default
486  errmsg = 'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
487  call store_error(errmsg)
488  end select
489  call this%obs%SaveOneSimval(obsrv, v)
490  end do
491  else
492  call this%obs%SaveOneSimval(obsrv, dnodata)
493  end if
494  end do
495  end subroutine zdg_bd_obs
496 
497  ! -- Procedure related to time series
498 
499  !> @brief Assign time series links for the package
500  !!
501  !! Assign the time series links for the ZDG package. Only
502  !! the Q variable can be defined with time series.
503  !!
504  !<
505  subroutine zdg_rp_ts(this)
506  ! -- dummy variables
507  class(swfzdgtype), intent(inout) :: this !< SwfZdgType object
508  ! -- local variables
509  integer(I4B) :: i, nlinks
510  type(timeserieslinktype), pointer :: tslink => null()
511  !
512  ! -- set up the time series links
513  nlinks = this%TsManager%boundtslinks%Count()
514  do i = 1, nlinks
515  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
516  if (associated(tslink)) then
517  if (tslink%JCol == 1) then
518  tslink%Text = 'Q'
519  end if
520  end if
521  end do
522  end subroutine zdg_rp_ts
523 
524  !> @ brief Return a bound value
525  !!
526  !! Return a bound value associated with an ncolbnd index
527  !! and row.
528  !!
529  !<
530  function zdg_bound_value(this, col, row) result(bndval)
531  ! -- modules
532  use constantsmodule, only: dzero
533  ! -- dummy variables
534  class(swfzdgtype), intent(inout) :: this
535  integer(I4B), intent(in) :: col
536  integer(I4B), intent(in) :: row
537  ! -- result
538  real(dp) :: bndval
539  !
540  select case (col)
541  case (1)
542  bndval = this%idcxs(row)
543  case (2)
544  bndval = this%width(row)
545  case (3)
546  bndval = this%slope(row)
547  case (4)
548  bndval = this%rough(row)
549  case default
550  errmsg = 'Programming error. ZDG bound value requested column '&
551  &'outside range of ncolbnd (1).'
552  call store_error(errmsg)
553  call store_error_filename(this%input_fname)
554  end select
555  end function zdg_bound_value
556 
557 end module swfzdgmodule
This module contains the extended boundary package.
This module contains the base boundary package.
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:95
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:68
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
Definition: Constants.f90:39
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
subroutine, public urword(line, icol, istart, istop, ncode, n, r, iout, in)
Extract a word from a string.
This module defines variable data types.
Definition: kind.f90:8
real(dp) function, public get_perturbation(x)
Calculate a numerical perturbation given the value of x.
Definition: MathUtil.f90:372
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains the derived types ObserveType and ObsDataType.
Definition: Observe.f90:15
This module contains the derived type ObsType.
Definition: Obs.f90:127
subroutine, public defaultobsidprocessor(obsrv, dis, inunitobs, iout)
@ brief Process IDstring provided for each observation
Definition: Obs.f90:246
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 squadratic(x, range, dydx, y)
@ brief sQuadratic
This module contains the ZDG package methods.
Definition: swf-zdg.f90:7
subroutine zdg_allocate_scalars(this)
@ brief Allocate scalars
Definition: swf-zdg.f90:131
subroutine define_listlabel(this)
@ brief Define the list label for the package
Definition: swf-zdg.f90:390
logical function zdg_obs_supported(this)
Determine if observations are supported.
Definition: swf-zdg.f90:422
subroutine zdg_cf(this)
@ brief Formulate the package hcof and rhs terms.
Definition: swf-zdg.f90:261
subroutine zdg_fc(this, rhs, ia, idxglo, matrix_sln)
@ brief Copy hcof and rhs terms into solution.
Definition: swf-zdg.f90:352
real(dp) function zdg_bound_value(this, col, row)
@ brief Return a bound value
Definition: swf-zdg.f90:531
subroutine log_zdg_options(this, found)
@ brief Log SWF specific package options
Definition: swf-zdg.f90:233
subroutine, public zdg_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath, dis, cxs, unitconv)
@ brief Create a new package object
Definition: swf-zdg.f90:74
subroutine zdg_rp_ts(this)
Assign time series links for the package.
Definition: swf-zdg.f90:506
character(len=16) text
package flow text string
Definition: swf-zdg.f90:31
subroutine zdg_da(this)
@ brief Deallocate package memory
Definition: swf-zdg.f90:186
subroutine zdg_options(this)
@ brief Source additional options for package
Definition: swf-zdg.f90:210
character(len=lenftype) ftype
package ftype
Definition: swf-zdg.f90:30
subroutine zdg_allocate_arrays(this, nodelist, auxvar)
@ brief Allocate arrays
Definition: swf-zdg.f90:152
subroutine zdg_bd_obs(this)
Save observations for the package.
Definition: swf-zdg.f90:456
real(dp) function qcalc(this, i, depth, unitconv)
Definition: swf-zdg.f90:322
subroutine zdg_df_obs(this)
Define the observation types available in the package.
Definition: swf-zdg.f90:435
type(timeserieslinktype) function, pointer, public gettimeserieslinkfromlist(list, indx)
Get time series link from a list.
@ brief BndType