MODFLOW 6  version 6.7.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 410 of file swf-zdg.f90.

411  ! -- dummy variables
412  class(SwfZdgType), intent(inout) :: this !< SwfZdgType object
413  !
414  ! -- create the header list label
415  this%listlabel = trim(this%filtyp)//' NO.'
416  if (this%dis%ndim == 3) then
417  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
418  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
419  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
420  elseif (this%dis%ndim == 2) then
421  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
422  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
423  else
424  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
425  end if
426  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'FLOW RATE'
427  if (this%inamedbound == 1) then
428  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
429  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 233 of file swf-zdg.f90.

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

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

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

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

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

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

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

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

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

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

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

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

443  ! -- dummy variables
444  class(SwfZdgType) :: this !< SwfZdgType object
445  !
446  ! -- set boolean
447  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 210 of file swf-zdg.f90.

211  ! -- modules
212  use inputoutputmodule, only: urword
215  ! -- dummy variables
216  class(SwfZdgType), intent(inout) :: this !< SwfZdgType object
217  ! -- local variables
218  type(SwfZdgParamFoundType) :: found
219  ! -- formats
220  !
221  ! -- source base BndExtType options
222  call this%BndExtType%source_options()
223  !
224  ! -- source options from input context
225  ! none
226  !
227  ! -- log SWF specific options
228  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 258 of file swf-zdg.f90.

259  use tdismodule, only: kper
260  ! -- dummy
261  class(SwfZdgType), intent(inout) :: this
262  ! -- local
263  !
264  if (this%iper /= kper) return
265  !
266  ! -- Call the parent class read and prepare
267  call this%BndExtType%bnd_rp()
268  !
269  ! -- Write the list to iout if requested
270  if (this%iprpak /= 0) then
271  call this%write_list()
272  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 526 of file swf-zdg.f90.

527  ! -- dummy variables
528  class(SwfZdgType), intent(inout) :: this !< SwfZdgType object
529  ! -- local variables
530  integer(I4B) :: i, nlinks
531  type(TimeSeriesLinkType), pointer :: tslink => null()
532  !
533  ! -- set up the time series links
534  nlinks = this%TsManager%boundtslinks%Count()
535  do i = 1, nlinks
536  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
537  if (associated(tslink)) then
538  if (tslink%JCol == 1) then
539  tslink%Text = 'Q'
540  end if
541  end if
542  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