MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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 get_cond (this, i, depth, absdhdxsq, unitconv)
 Calculate conductance-like term. More...
 
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 449 of file swf-zdg.f90.

450  ! -- dummy variables
451  class(SwfZdgType), intent(inout) :: this !< SwfZdgType object
452  !
453  ! -- create the header list label
454  this%listlabel = trim(this%filtyp)//' NO.'
455  if (this%dis%ndim == 3) then
456  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
457  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
458  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
459  elseif (this%dis%ndim == 2) then
460  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
461  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
462  else
463  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
464  end if
465  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'FLOW RATE'
466  if (this%inamedbound == 1) then
467  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
468  end if
469  !
470  ! -- return
471  return

◆ get_cond()

real(dp) function swfzdgmodule::get_cond ( class(swfzdgtype this,
integer(i4b), intent(in)  i,
real(dp), intent(in)  depth,
real(dp), intent(in)  absdhdxsq,
real(dp), intent(in)  unitconv 
)

Conductance normally has a dx term in the denominator but that is not included here, as the flow is calculated using Q = C * slope. The returned c value from this function has dimensions of L3/T.

Parameters
[in]iboundary number
[in]depthsimulated depth (stage - elevation) in reach n for this iteration
[in]absdhdxsqabsolute value of simulated hydraulic gradient
[in]unitconvconversion factor for roughness to length and time units of meters and seconds

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

369  ! -- modules
370  ! -- dummy
371  class(SwfZdgType) :: this
372  integer(I4B), intent(in) :: i !< boundary number
373  real(DP), intent(in) :: depth !< simulated depth (stage - elevation) in reach n for this iteration
374  real(DP), intent(in) :: absdhdxsq !< absolute value of simulated hydraulic gradient
375  real(DP), intent(in) :: unitconv !< conversion factor for roughness to length and time units of meters and seconds
376  ! -- local
377  integer(I4B) :: idcxs
378  real(DP) :: c
379  real(DP) :: width
380  real(DP) :: rough
381  real(DP) :: slope
382  real(DP) :: roughc
383  real(DP) :: a
384  real(DP) :: r
385  real(DP) :: conveyance
386  !
387  idcxs = this%idcxs(i)
388  width = this%width(i)
389  rough = this%rough(i)
390  slope = this%slope(i)
391  roughc = this%cxs%get_roughness(idcxs, width, depth, rough, &
392  slope)
393  a = this%cxs%get_area(idcxs, width, depth)
394  r = this%cxs%get_hydraulic_radius(idcxs, width, depth, area=a)
395 
396  !conveyance = a * r**DTWOTHIRDS / roughc
397  conveyance = this%cxs%get_conveyance(idcxs, width, depth, rough)
398  c = conveyance / absdhdxsq
399 

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

253  ! -- modules
255  ! -- dummy variables
256  class(SwfZdgType), intent(inout) :: this !< BndExtType object
257  type(SwfZdgParamFoundType), intent(in) :: found
258  ! -- local variables
259  ! -- format
260  !
261  ! -- log found options
262  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
263  //' OPTIONS'
264  !
265  ! if (found%mover) then
266  ! write (this%iout, '(4x,A)') 'MOVER OPTION ENABLED'
267  ! end if
268  !
269  ! -- close logging block
270  write (this%iout, '(1x,a)') &
271  'END OF '//trim(adjustl(this%text))//' OPTIONS'
272  !
273  ! -- return
274  return

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

163  ! -- modules
165  ! -- dummy
166  class(SwfZdgType) :: this
167  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
168  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
169  ! -- local
170  !
171  ! -- call BndExtType allocate scalars
172  call this%BndExtType%allocate_arrays(nodelist, auxvar)
173  !
174  ! -- set array input context pointer
175  call mem_setptr(this%idcxs, 'IDCXS', this%input_mempath)
176  call mem_setptr(this%width, 'WIDTH', this%input_mempath)
177  call mem_setptr(this%slope, 'SLOPE', this%input_mempath)
178  call mem_setptr(this%rough, 'ROUGH', this%input_mempath)
179  !
180  ! -- checkin array input context pointer
181  call mem_checkin(this%idcxs, 'IDCXS', this%memoryPath, &
182  'IDCXS', this%input_mempath)
183  call mem_checkin(this%width, 'WIDTH', this%memoryPath, &
184  'WIDTH', this%input_mempath)
185  call mem_checkin(this%slope, 'SLOPE', this%memoryPath, &
186  'SLOPE', this%input_mempath)
187  call mem_checkin(this%rough, 'ROUGH', this%memoryPath, &
188  'ROUGH', this%input_mempath)
189  !
190  ! -- return
191  return

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

139  ! -- modules
141  ! -- dummy variables
142  class(SwfZdgType) :: this !< SwfZdgType object
143  !
144  ! -- call base type allocate scalars
145  call this%BndExtType%allocate_scalars()
146  !
147  ! -- allocate the object and assign values to object variables
148  call mem_allocate(this%unitconv, 'UNITCONV', this%memoryPath)
149  !
150  ! -- Set values
151  this%unitconv = dzero
152  !
153  ! -- return
154  return

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

525  ! -- dummy variables
526  class(SwfZdgType) :: this !< SwfZdgType object
527  ! -- local variables
528  integer(I4B) :: i
529  integer(I4B) :: n
530  integer(I4B) :: jj
531  real(DP) :: v
532  type(ObserveType), pointer :: obsrv => null()
533  !
534  ! -- clear the observations
535  call this%obs%obs_bd_clear()
536  !
537  ! -- Save simulated values for all of package's observations.
538  do i = 1, this%obs%npakobs
539  obsrv => this%obs%pakobs(i)%obsrv
540  if (obsrv%BndFound) then
541  do n = 1, obsrv%indxbnds_count
542  v = dnodata
543  jj = obsrv%indxbnds(n)
544  select case (obsrv%ObsTypeId)
545  case ('TO-MVR')
546  if (this%imover == 1) then
547  v = this%pakmvrobj%get_qtomvr(jj)
548  if (v > dzero) then
549  v = -v
550  end if
551  end if
552  case ('ZDG')
553  v = this%simvals(jj)
554  case default
555  errmsg = 'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
556  call store_error(errmsg)
557  end select
558  call this%obs%SaveOneSimval(obsrv, v)
559  end do
560  else
561  call this%obs%SaveOneSimval(obsrv, dnodata)
562  end if
563  end do
564  !
565  ! -- return
566  return
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 605 of file swf-zdg.f90.

606  ! -- modules
607  use constantsmodule, only: dzero
608  ! -- dummy variables
609  class(SwfZdgType), intent(inout) :: this !< BndExtType object
610  integer(I4B), intent(in) :: col
611  integer(I4B), intent(in) :: row
612  ! -- result
613  real(DP) :: bndval
614  !
615  select case (col)
616  case (1)
617  bndval = this%idcxs(row)
618  case (2)
619  bndval = this%width(row)
620  case (3)
621  bndval = this%slope(row)
622  case (4)
623  bndval = this%rough(row)
624  case default
625  errmsg = 'Programming error. ZDG bound value requested column '&
626  &'outside range of ncolbnd (1).'
627  call store_error(errmsg)
628  call store_error_filename(this%input_fname)
629  end select
630  !
631  ! -- return
632  return
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
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 306 of file swf-zdg.f90.

307  ! modules
309  ! -- dummy variables
310  class(SwfZdgType) :: this !< SwfZdgType object
311  ! -- local variables
312  integer(I4B) :: i, node
313  real(DP) :: q
314  real(DP) :: qeps
315  real(DP) :: absdhdxsq
316  real(DP) :: depth
317  real(DP) :: cond
318  real(DP) :: derv
319  real(DP) :: eps
320  !
321  ! -- Return if no inflows
322  if (this%nbound == 0) return
323  !
324  ! -- Calculate hcof and rhs for each zdg entry
325  do i = 1, this%nbound
326 
327  node = this%nodelist(i)
328  if (this%ibound(node) <= 0) then
329  this%hcof(i) = dzero
330  this%rhs(i) = dzero
331  cycle
332  end if
333 
334  ! -- calculate terms and add to hcof and rhs
335  absdhdxsq = this%slope(i)**dhalf
336  depth = this%xnew(node) - this%dis%bot(node)
337 
338  ! -- calculate unperturbed q
339  ! TODO: UNITCONV?!
340  cond = this%get_cond(i, depth, absdhdxsq, this%unitconv)
341  q = -cond * this%slope(i)
342 
343  ! -- calculate perturbed q
344  eps = get_perturbation(depth)
345  cond = this%get_cond(i, depth + eps, absdhdxsq, this%unitconv)
346  qeps = -cond * this%slope(i)
347 
348  ! -- calculate derivative
349  derv = (qeps - q) / eps
350 
351  ! -- add terms to hcof and rhs
352  this%hcof(i) = derv
353  this%rhs(i) = -q + derv * this%xnew(node)
354 
355  end do
356  !
357  return
real(dp) function, public get_perturbation(x)
Calculate a numerical perturbation given the value of x.
Definition: MathUtil.f90:446
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
127  !
128  ! -- return
129  return
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 199 of file swf-zdg.f90.

200  ! -- modules
202  ! -- dummy variables
203  class(SwfZdgType) :: this !< SwfZdgType object
204  !
205  ! -- Deallocate parent package
206  call this%BndExtType%bnd_da()
207  !
208  ! -- arrays
209  call mem_deallocate(this%idcxs, 'IDCXS', this%memoryPath)
210  call mem_deallocate(this%width, 'WIDTH', this%memoryPath)
211  call mem_deallocate(this%slope, 'SLOPE', this%memoryPath)
212  call mem_deallocate(this%rough, 'ROUGH', this%memoryPath)
213  !
214  ! -- scalars
215  call mem_deallocate(this%unitconv)
216  !
217  ! -- return
218  return

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

501  ! -- dummy variables
502  class(SwfZdgType) :: this !< SwfZdgType object
503  ! -- local variables
504  integer(I4B) :: indx
505  !
506  ! -- initialize observations
507  call this%obs%StoreObsType('zdg', .true., indx)
508  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
509  !
510  ! -- Store obs type and assign procedure pointer
511  ! for to-mvr observation type.
512  call this%obs%StoreObsType('to-mvr', .true., indx)
513  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
514  !
515  ! -- return
516  return
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 408 of file swf-zdg.f90.

409  ! -- dummy variables
410  class(SwfZdgType) :: this !< SwfZdgType object
411  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
412  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
413  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
414  class(MatrixBaseType), pointer :: matrix_sln !< solution coefficient matrix
415  ! -- local variables
416  integer(I4B) :: i
417  integer(I4B) :: n
418  integer(I4B) :: ipos
419  !
420  ! -- pakmvrobj fc
421  if (this%imover == 1) then
422  call this%pakmvrobj%fc()
423  end if
424  !
425  ! -- Copy package rhs and hcof into solution rhs and amat
426  do i = 1, this%nbound
427  n = this%nodelist(i)
428  rhs(n) = rhs(n) + this%rhs(i)
429  ipos = ia(n)
430  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
431  !
432  ! -- If mover is active and this zdg item is discharging,
433  ! store available water (as positive value).
434  if (this%imover == 1 .and. this%rhs(i) > dzero) then
435  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
436  end if
437  end do
438  !
439  ! -- return
440  return

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

485  ! -- dummy variables
486  class(SwfZdgType) :: this !< SwfZdgType object
487  !
488  ! -- set boolean
489  zdg_obs_supported = .true.
490  !
491  ! -- return
492  return

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

227  ! -- modules
228  use inputoutputmodule, only: urword
231  ! -- dummy variables
232  class(SwfZdgType), intent(inout) :: this !< SwfZdgType object
233  ! -- local variables
234  type(SwfZdgParamFoundType) :: found
235  ! -- formats
236  !
237  ! -- source base BndExtType options
238  call this%BndExtType%source_options()
239  !
240  ! -- source options from input context
241  ! none
242  !
243  ! -- log SWF specific options
244  call this%log_zdg_options(found)
245  !
246  ! -- return
247  return
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 280 of file swf-zdg.f90.

281  use tdismodule, only: kper
282  ! -- dummy
283  class(SwfZdgType), intent(inout) :: this
284  ! -- local
285  !
286  if (this%iper /= kper) return
287  !
288  ! -- Call the parent class read and prepare
289  call this%BndExtType%bnd_rp()
290  !
291  ! -- Write the list to iout if requested
292  if (this%iprpak /= 0) then
293  call this%write_list()
294  end if
295  !
296  ! -- return
297  return
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 577 of file swf-zdg.f90.

578  ! -- dummy variables
579  class(SwfZdgType), intent(inout) :: this !< SwfZdgType object
580  ! -- local variables
581  integer(I4B) :: i, nlinks
582  type(TimeSeriesLinkType), pointer :: tslink => null()
583  !
584  ! -- set up the time series links
585  nlinks = this%TsManager%boundtslinks%Count()
586  do i = 1, nlinks
587  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
588  if (associated(tslink)) then
589  if (tslink%JCol == 1) then
590  tslink%Text = 'Q'
591  end if
592  end if
593  end do
594  !
595  ! -- return
596  return
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