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

This module contains the ZDG package methods. More...

Data Types

type  swfzdgtype
 

Functions/Subroutines

subroutine, public zdg_create (packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath, dis, cxs, unitconv)
 @ brief Create a new package object More...
 
subroutine zdg_allocate_scalars (this)
 @ brief Allocate scalars More...
 
subroutine zdg_allocate_arrays (this, nodelist, auxvar)
 @ brief Allocate arrays More...
 
subroutine zdg_da (this)
 @ brief Deallocate package memory More...
 
subroutine zdg_options (this)
 @ brief Source additional options for package More...
 
subroutine log_zdg_options (this, found)
 @ brief Log SWF specific package options More...
 
subroutine zdg_cf (this)
 @ brief Formulate the package hcof and rhs terms. More...
 
real(dp) function qcalc (this, i, depth, unitconv)
 
subroutine zdg_fc (this, rhs, ia, idxglo, matrix_sln)
 @ brief Copy hcof and rhs terms into solution. More...
 
subroutine define_listlabel (this)
 @ brief Define the list label for the package More...
 
logical function zdg_obs_supported (this)
 Determine if observations are supported. More...
 
subroutine zdg_df_obs (this)
 Define the observation types available in the package. More...
 
subroutine zdg_bd_obs (this)
 Save observations for the package. More...
 
subroutine zdg_rp_ts (this)
 Assign time series links for the package. More...
 
real(dp) function zdg_bound_value (this, col, row)
 @ brief Return a bound value More...
 

Variables

character(len=lenftype) ftype = 'ZDG'
 package ftype More...
 
character(len=16) text = ' ZDG'
 package flow text string More...
 

Detailed Description

This module can be used to represent outflow from streams using a zero-depth-gradient boundary.

Function/Subroutine Documentation

◆ define_listlabel()

subroutine swfzdgmodule::define_listlabel ( class(swfzdgtype), intent(inout)  this)
private

Method defined the list label for the ZDG package. The list label is the heading that is written to iout when PRINT_INPUT option is used.

Parameters
[in,out]thisSwfZdgType object

Definition at line 389 of file swf-zdg.f90.

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

◆ log_zdg_options()

subroutine swfzdgmodule::log_zdg_options ( class(swfzdgtype), intent(inout)  this,
type(swfzdgparamfoundtype), intent(in)  found 
)

Definition at line 232 of file swf-zdg.f90.

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'

◆ qcalc()

real(dp) function swfzdgmodule::qcalc ( class(swfzdgtype this,
integer(i4b), intent(in)  i,
real(dp), intent(in)  depth,
real(dp), intent(in)  unitconv 
)
Parameters
[in]iboundary number
[in]depthsimulated depth (stage - elevation) in reach n for this iteration
[in]unitconvconversion factor for roughness to length and time units of meters and seconds

Definition at line 321 of file swf-zdg.f90.

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 

◆ zdg_allocate_arrays()

subroutine swfzdgmodule::zdg_allocate_arrays ( class(swfzdgtype this,
integer(i4b), dimension(:), optional, pointer, contiguous  nodelist,
real(dp), dimension(:, :), optional, pointer, contiguous  auxvar 
)

Allocate and initialize arrays for the SWF package

Definition at line 151 of file swf-zdg.f90.

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)

◆ zdg_allocate_scalars()

subroutine swfzdgmodule::zdg_allocate_scalars ( class(swfzdgtype this)
private

Allocate and initialize scalars for the ZDG package. The base model allocate scalars method is also called.

Parameters
thisSwfZdgType object

Definition at line 130 of file swf-zdg.f90.

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

◆ zdg_bd_obs()

subroutine swfzdgmodule::zdg_bd_obs ( class(swfzdgtype this)
private

Method to save simulated values for the ZDG package.

Parameters
thisSwfZdgType object

Definition at line 455 of file swf-zdg.f90.

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
Here is the call graph for this function:

◆ zdg_bound_value()

real(dp) function swfzdgmodule::zdg_bound_value ( class(swfzdgtype), intent(inout)  this,
integer(i4b), intent(in)  col,
integer(i4b), intent(in)  row 
)
private

Return a bound value associated with an ncolbnd index and row.

Definition at line 530 of file swf-zdg.f90.

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
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
Here is the call graph for this function:

◆ zdg_cf()

subroutine swfzdgmodule::zdg_cf ( class(swfzdgtype this)

Formulate the hcof and rhs terms for the ZDG package that will be added to the coefficient matrix and right-hand side vector.

Parameters
thisSwfZdgType object

Definition at line 260 of file swf-zdg.f90.

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
real(dp) function, public get_perturbation(x)
Calculate a numerical perturbation given the value of x.
Definition: MathUtil.f90:372
subroutine squadratic(x, range, dydx, y)
@ brief sQuadratic
Here is the call graph for this function:

◆ zdg_create()

subroutine, public swfzdgmodule::zdg_create ( class(bndtype), pointer  packobj,
integer(i4b), intent(in)  id,
integer(i4b), intent(in)  ibcnum,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout,
character(len=*), intent(in)  namemodel,
character(len=*), intent(in)  pakname,
character(len=*), intent(in)  mempath,
class(disbasetype), intent(inout), pointer  dis,
type(swfcxstype), intent(in), pointer  cxs,
real(dp), intent(in)  unitconv 
)

Create a new ZDG Package object

Parameters
packobjpointer to default package type
[in]idpackage id
[in]ibcnumboundary condition number
[in]inunitunit number of ZDG package input file
[in]ioutunit number of model listing file
[in]namemodelmodel name
[in]paknamepackage name
[in]mempathinput mempath
[in,out]disthe pointer to the discretization
[in]cxsthe pointer to the cxs package
[in]unitconvunit conversion for roughness

Definition at line 72 of file swf-zdg.f90.

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
Here is the call graph for this function:
Here is the caller graph for this function:

◆ zdg_da()

subroutine swfzdgmodule::zdg_da ( class(swfzdgtype this)

Deallocate SWF package scalars and arrays.

Parameters
thisSwfZdgType object

Definition at line 185 of file swf-zdg.f90.

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)

◆ zdg_df_obs()

subroutine swfzdgmodule::zdg_df_obs ( class(swfzdgtype this)
private

Method to define the observation types available in the ZDG package.

Parameters
thisSwfZdgType object

Definition at line 434 of file swf-zdg.f90.

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
Here is the call graph for this function:

◆ zdg_fc()

subroutine swfzdgmodule::zdg_fc ( class(swfzdgtype this,
real(dp), dimension(:), intent(inout)  rhs,
integer(i4b), dimension(:), intent(in)  ia,
integer(i4b), dimension(:), intent(in)  idxglo,
class(matrixbasetype), pointer  matrix_sln 
)
private

Add the hcof and rhs terms for the ZDG package to the coefficient matrix and right-hand side vector.

Parameters
thisSwfZdgType object
[in,out]rhsright-hand side vector for model
[in]iasolution CRS row pointers
[in]idxglomapping vector for model (local) to solution (global)
matrix_slnsolution coefficient matrix

Definition at line 351 of file swf-zdg.f90.

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

◆ zdg_obs_supported()

logical function swfzdgmodule::zdg_obs_supported ( class(swfzdgtype this)
private

Function to determine if observations are supported by the ZDG package. Observations are supported by the ZDG package.

Returns
zdg_obs_supported boolean indicating if observations are supported
Parameters
thisSwfZdgType object

Definition at line 421 of file swf-zdg.f90.

422  ! -- dummy variables
423  class(SwfZdgType) :: this !< SwfZdgType object
424  !
425  ! -- set boolean
426  zdg_obs_supported = .true.

◆ zdg_options()

subroutine swfzdgmodule::zdg_options ( class(swfzdgtype), intent(inout)  this)

Source additional options for SWF package.

Parameters
[in,out]thisSwfZdgType object

Definition at line 209 of file swf-zdg.f90.

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)
subroutine, public urword(line, icol, istart, istop, ncode, n, r, iout, in)
Extract a word from a string.
Here is the call graph for this function:

◆ zdg_rp_ts()

subroutine swfzdgmodule::zdg_rp_ts ( class(swfzdgtype), intent(inout)  this)
private

Assign the time series links for the ZDG package. Only the Q variable can be defined with time series.

Parameters
[in,out]thisSwfZdgType object

Definition at line 505 of file swf-zdg.f90.

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
Here is the call graph for this function:

Variable Documentation

◆ ftype

character(len=lenftype) swfzdgmodule::ftype = 'ZDG'
private

Definition at line 30 of file swf-zdg.f90.

30  character(len=LENFTYPE) :: ftype = 'ZDG' !< package ftype

◆ text

character(len=16) swfzdgmodule::text = ' ZDG'
private

Definition at line 31 of file swf-zdg.f90.

31  character(len=16) :: text = ' ZDG' !< package flow text string