MODFLOW 6  version 6.6.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_rp (this)
 @ brief SWF read and prepare 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 414 of file swf-zdg.f90.

415  ! -- dummy variables
416  class(SwfZdgType), intent(inout) :: this !< SwfZdgType object
417  !
418  ! -- create the header list label
419  this%listlabel = trim(this%filtyp)//' NO.'
420  if (this%dis%ndim == 3) then
421  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
422  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
423  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
424  elseif (this%dis%ndim == 2) then
425  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
426  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
427  else
428  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
429  end if
430  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'FLOW RATE'
431  if (this%inamedbound == 1) then
432  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
433  end if

◆ log_zdg_options()

subroutine swfzdgmodule::log_zdg_options ( class(swfzdgtype), intent(inout)  this,
type(swfzdgparamfoundtype), intent(in)  found 
)
Parameters
[in,out]thisBndExtType object

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

238  ! -- modules
240  ! -- dummy variables
241  class(SwfZdgType), intent(inout) :: this !< BndExtType object
242  type(SwfZdgParamFoundType), intent(in) :: found
243  ! -- local variables
244  ! -- format
245  !
246  ! -- log found options
247  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
248  //' OPTIONS'
249  !
250  ! if (found%mover) then
251  ! write (this%iout, '(4x,A)') 'MOVER OPTION ENABLED'
252  ! end if
253  !
254  ! -- close logging block
255  write (this%iout, '(1x,a)') &
256  '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 346 of file swf-zdg.f90.

347  ! dummy
348  class(SwfZdgType) :: this
349  integer(I4B), intent(in) :: i !< boundary number
350  real(DP), intent(in) :: depth !< simulated depth (stage - elevation) in reach n for this iteration
351  real(DP), intent(in) :: unitconv !< conversion factor for roughness to length and time units of meters and seconds
352  ! return
353  real(DP) :: q
354  ! local
355  integer(I4B) :: idcxs
356  real(DP) :: width
357  real(DP) :: rough
358  real(DP) :: slope
359  real(DP) :: conveyance
360 
361  idcxs = this%idcxs(i)
362  width = this%width(i)
363  rough = this%rough(i)
364  slope = this%slope(i)
365  conveyance = this%cxs%get_conveyance(idcxs, width, depth, rough)
366  q = conveyance * slope**dhalf * unitconv
367 

◆ 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 156 of file swf-zdg.f90.

157  ! -- modules
159  ! -- dummy
160  class(SwfZdgType) :: this
161  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
162  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
163  ! -- local
164  !
165  ! -- call BndExtType allocate scalars
166  call this%BndExtType%allocate_arrays(nodelist, auxvar)
167  !
168  ! -- set array input context pointer
169  call mem_setptr(this%idcxs, 'IDCXS', this%input_mempath)
170  call mem_setptr(this%width, 'WIDTH', this%input_mempath)
171  call mem_setptr(this%slope, 'SLOPE', this%input_mempath)
172  call mem_setptr(this%rough, 'ROUGH', this%input_mempath)
173  !
174  ! -- checkin array input context pointer
175  call mem_checkin(this%idcxs, 'IDCXS', this%memoryPath, &
176  'IDCXS', this%input_mempath)
177  call mem_checkin(this%width, 'WIDTH', this%memoryPath, &
178  'WIDTH', this%input_mempath)
179  call mem_checkin(this%slope, 'SLOPE', this%memoryPath, &
180  'SLOPE', this%input_mempath)
181  call mem_checkin(this%rough, 'ROUGH', this%memoryPath, &
182  '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 135 of file swf-zdg.f90.

136  ! -- modules
138  ! -- dummy variables
139  class(SwfZdgType) :: this !< SwfZdgType object
140  !
141  ! -- call base type allocate scalars
142  call this%BndExtType%allocate_scalars()
143  !
144  ! -- allocate the object and assign values to object variables
145  call mem_allocate(this%unitconv, 'UNITCONV', this%memoryPath)
146  !
147  ! -- Set values
148  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 480 of file swf-zdg.f90.

481  ! -- dummy variables
482  class(SwfZdgType) :: this !< SwfZdgType object
483  ! -- local variables
484  integer(I4B) :: i
485  integer(I4B) :: n
486  integer(I4B) :: jj
487  real(DP) :: v
488  type(ObserveType), pointer :: obsrv => null()
489  !
490  ! -- clear the observations
491  call this%obs%obs_bd_clear()
492  !
493  ! -- Save simulated values for all of package's observations.
494  do i = 1, this%obs%npakobs
495  obsrv => this%obs%pakobs(i)%obsrv
496  if (obsrv%BndFound) then
497  do n = 1, obsrv%indxbnds_count
498  v = dnodata
499  jj = obsrv%indxbnds(n)
500  select case (obsrv%ObsTypeId)
501  case ('TO-MVR')
502  if (this%imover == 1) then
503  v = this%pakmvrobj%get_qtomvr(jj)
504  if (v > dzero) then
505  v = -v
506  end if
507  end if
508  case ('ZDG')
509  v = this%simvals(jj)
510  case default
511  errmsg = 'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
512  call store_error(errmsg)
513  end select
514  call this%obs%SaveOneSimval(obsrv, v)
515  end do
516  else
517  call this%obs%SaveOneSimval(obsrv, dnodata)
518  end if
519  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.

Parameters
[in,out]thisBndExtType object

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

556  ! -- modules
557  use constantsmodule, only: dzero
558  ! -- dummy variables
559  class(SwfZdgType), intent(inout) :: this !< BndExtType object
560  integer(I4B), intent(in) :: col
561  integer(I4B), intent(in) :: row
562  ! -- result
563  real(DP) :: bndval
564  !
565  select case (col)
566  case (1)
567  bndval = this%idcxs(row)
568  case (2)
569  bndval = this%width(row)
570  case (3)
571  bndval = this%slope(row)
572  case (4)
573  bndval = this%rough(row)
574  case default
575  errmsg = 'Programming error. ZDG bound value requested column '&
576  &'outside range of ncolbnd (1).'
577  call store_error(errmsg)
578  call store_error_filename(this%input_fname)
579  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 285 of file swf-zdg.f90.

286  ! modules
288  use smoothingmodule, only: squadratic
289  ! dummy variables
290  class(SwfZdgType) :: this !< SwfZdgType object
291  ! local variables
292  integer(I4B) :: i, node
293  real(DP) :: q
294  real(DP) :: qeps
295  real(DP) :: absdhdxsq
296  real(DP) :: depth
297  real(DP) :: derv
298  real(DP) :: eps
299  real(DP) :: range = 1.d-6
300  real(DP) :: dydx
301  real(DP) :: smooth_factor
302  !
303  ! -- Return if no inflows
304  if (this%nbound == 0) return
305  !
306  ! -- Calculate hcof and rhs for each zdg entry
307  do i = 1, this%nbound
308 
309  node = this%nodelist(i)
310  if (this%ibound(node) <= 0) then
311  this%hcof(i) = dzero
312  this%rhs(i) = dzero
313  cycle
314  end if
315 
316  ! -- calculate terms and add to hcof and rhs
317  absdhdxsq = this%slope(i)**dhalf
318  depth = this%xnew(node) - this%dis%bot(node)
319 
320  ! smooth the depth
321  call squadratic(depth, range, dydx, smooth_factor)
322  depth = depth * smooth_factor
323 
324  ! -- calculate unperturbed q
325  q = -this%qcalc(i, depth, this%unitconv)
326 
327  ! -- calculate perturbed q
328  eps = get_perturbation(depth)
329  qeps = -this%qcalc(i, depth + eps, this%unitconv)
330 
331  ! -- calculate derivative
332  derv = (qeps - q) / eps
333 
334  ! -- add terms to hcof and rhs
335  this%hcof(i) = derv
336  this%rhs(i) = -q + derv * this%xnew(node)
337 
338  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 77 of file swf-zdg.f90.

79  ! -- dummy variables
80  class(BndType), pointer :: packobj !< pointer to default package type
81  integer(I4B), intent(in) :: id !< package id
82  integer(I4B), intent(in) :: ibcnum !< boundary condition number
83  integer(I4B), intent(in) :: inunit !< unit number of ZDG package input file
84  integer(I4B), intent(in) :: iout !< unit number of model listing file
85  character(len=*), intent(in) :: namemodel !< model name
86  character(len=*), intent(in) :: pakname !< package name
87  character(len=*), intent(in) :: mempath !< input mempath
88  class(DisBaseType), pointer, intent(inout) :: dis !< the pointer to the discretization
89  type(SwfCxsType), pointer, intent(in) :: cxs !< the pointer to the cxs package
90  real(DP), intent(in) :: unitconv !< unit conversion for roughness
91  ! -- local variables
92  type(SwfZdgType), pointer :: zdgobj
93  !
94  ! -- allocate the object and assign values to object variables
95  allocate (zdgobj)
96  packobj => zdgobj
97  !
98  ! -- create name and memory path
99  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
100  packobj%text = text
101  !
102  ! -- allocate scalars
103  call zdgobj%allocate_scalars()
104  !
105  ! -- initialize package
106  call packobj%pack_initialize()
107 
108  packobj%inunit = inunit
109  packobj%iout = iout
110  packobj%id = id
111  packobj%ibcnum = ibcnum
112  packobj%ncolbnd = 1
113  packobj%iscloc = 1
114  packobj%ictMemPath = create_mem_path(namemodel, 'DFW')
115  !
116  ! -- store pointer to disv1d
117  select type (dis)
118  type is (disv1dtype)
119  zdgobj%disv1d => dis
120  end select
121  !
122  ! -- store pointer to cxs
123  zdgobj%cxs => cxs
124  !
125  ! -- store unit conversion
126  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 190 of file swf-zdg.f90.

191  ! -- modules
193  ! -- dummy variables
194  class(SwfZdgType) :: this !< SwfZdgType object
195  !
196  ! -- Deallocate parent package
197  call this%BndExtType%bnd_da()
198  !
199  ! -- arrays
200  call mem_deallocate(this%idcxs, 'IDCXS', this%memoryPath)
201  call mem_deallocate(this%width, 'WIDTH', this%memoryPath)
202  call mem_deallocate(this%slope, 'SLOPE', this%memoryPath)
203  call mem_deallocate(this%rough, 'ROUGH', this%memoryPath)
204  !
205  ! -- scalars
206  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 459 of file swf-zdg.f90.

460  ! -- dummy variables
461  class(SwfZdgType) :: this !< SwfZdgType object
462  ! -- local variables
463  integer(I4B) :: indx
464  !
465  ! -- initialize observations
466  call this%obs%StoreObsType('zdg', .true., indx)
467  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
468  !
469  ! -- Store obs type and assign procedure pointer
470  ! for to-mvr observation type.
471  call this%obs%StoreObsType('to-mvr', .true., indx)
472  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 376 of file swf-zdg.f90.

377  ! -- dummy variables
378  class(SwfZdgType) :: this !< SwfZdgType object
379  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
380  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
381  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
382  class(MatrixBaseType), pointer :: matrix_sln !< solution coefficient matrix
383  ! -- local variables
384  integer(I4B) :: i
385  integer(I4B) :: n
386  integer(I4B) :: ipos
387  !
388  ! -- pakmvrobj fc
389  if (this%imover == 1) then
390  call this%pakmvrobj%fc()
391  end if
392  !
393  ! -- Copy package rhs and hcof into solution rhs and amat
394  do i = 1, this%nbound
395  n = this%nodelist(i)
396  rhs(n) = rhs(n) + this%rhs(i)
397  ipos = ia(n)
398  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
399  !
400  ! -- If mover is active and this zdg item is discharging,
401  ! store available water (as positive value).
402  if (this%imover == 1 .and. this%rhs(i) > dzero) then
403  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
404  end if
405  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 446 of file swf-zdg.f90.

447  ! -- dummy variables
448  class(SwfZdgType) :: this !< SwfZdgType object
449  !
450  ! -- set boolean
451  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 214 of file swf-zdg.f90.

215  ! -- modules
216  use inputoutputmodule, only: urword
219  ! -- dummy variables
220  class(SwfZdgType), intent(inout) :: this !< SwfZdgType object
221  ! -- local variables
222  type(SwfZdgParamFoundType) :: found
223  ! -- formats
224  !
225  ! -- source base BndExtType options
226  call this%BndExtType%source_options()
227  !
228  ! -- source options from input context
229  ! none
230  !
231  ! -- log SWF specific options
232  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()

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

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

263  use tdismodule, only: kper
264  ! -- dummy
265  class(SwfZdgType), intent(inout) :: this
266  ! -- local
267  !
268  if (this%iper /= kper) return
269  !
270  ! -- Call the parent class read and prepare
271  call this%BndExtType%bnd_rp()
272  !
273  ! -- Write the list to iout if requested
274  if (this%iprpak /= 0) then
275  call this%write_list()
276  end if
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23

◆ 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 530 of file swf-zdg.f90.

531  ! -- dummy variables
532  class(SwfZdgType), intent(inout) :: this !< SwfZdgType object
533  ! -- local variables
534  integer(I4B) :: i, nlinks
535  type(TimeSeriesLinkType), pointer :: tslink => null()
536  !
537  ! -- set up the time series links
538  nlinks = this%TsManager%boundtslinks%Count()
539  do i = 1, nlinks
540  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
541  if (associated(tslink)) then
542  if (tslink%JCol == 1) then
543  tslink%Text = 'Q'
544  end if
545  end if
546  end do
Here is the call graph for this function:

Variable Documentation

◆ ftype

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

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

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

◆ text

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

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

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