MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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
10  use constantsmodule, only: dzero, dem1, done, lenftype, dnodata, &
12  use simvariablesmodule, only: errmsg
15  use bndmodule, only: bndtype
16  use bndextmodule, only: bndexttype
19  use observemodule, only: observetype
25  use basedismodule, only: disbasetype
26  use disv1dmodule, only: disv1dtype
27  use swfcxsmodule, only: swfcxstype
28  !
29  implicit none
30  !
31  private
32  public :: zdg_create
33  !
34  character(len=LENFTYPE) :: ftype = 'ZDG' !< package ftype
35  character(len=16) :: text = ' ZDG' !< package flow text string
36  !
37  type, extends(bndexttype) :: swfzdgtype
38 
39  integer(I4B), dimension(:), pointer, contiguous :: idcxs => null() !< cross section id
40  real(dp), dimension(:), pointer, contiguous :: width => null() !< channel width
41  real(dp), dimension(:), pointer, contiguous :: slope => null() !< channel slope
42  real(dp), dimension(:), pointer, contiguous :: rough => null() !< channel roughness
43  real(dp), pointer :: unitconv => null() !< conversion factor for roughness to length and time units of meters and seconds
44 
45  ! -- pointers other objects
46  type(disv1dtype), pointer :: disv1d
47  type(swfcxstype), pointer :: cxs
48 
49  contains
50  procedure :: allocate_scalars => zdg_allocate_scalars
51  procedure :: allocate_arrays => zdg_allocate_arrays
52  procedure :: source_options => zdg_options
53  procedure :: log_zdg_options
54  procedure :: bnd_rp => zdg_rp
55  procedure :: bnd_cf => zdg_cf
56  procedure :: bnd_fc => zdg_fc
57  procedure :: bnd_da => zdg_da
58  procedure :: define_listlabel
59  procedure :: bound_value => zdg_bound_value
60  ! -- methods for observations
61  procedure, public :: bnd_obs_supported => zdg_obs_supported
62  procedure, public :: bnd_df_obs => zdg_df_obs
63  procedure, public :: bnd_bd_obs => zdg_bd_obs
64  ! -- methods for time series
65  procedure, public :: bnd_rp_ts => zdg_rp_ts
66  ! -- private
67  procedure, private :: get_cond
68  end type swfzdgtype
69 
70 contains
71 
72  !> @ brief Create a new package object
73  !!
74  !! Create a new ZDG Package object
75  !!
76  !<
77  subroutine zdg_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
78  mempath, dis, cxs, unitconv)
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
130  end subroutine zdg_create
131 
132  !> @ brief Allocate scalars
133  !!
134  !! Allocate and initialize scalars for the ZDG package. The base model
135  !! allocate scalars method is also called.
136  !!
137  !<
138  subroutine zdg_allocate_scalars(this)
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
155  end subroutine zdg_allocate_scalars
156 
157  !> @ brief Allocate arrays
158  !!
159  !! Allocate and initialize arrays for the SWF package
160  !!
161  !<
162  subroutine zdg_allocate_arrays(this, nodelist, auxvar)
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
192  end subroutine zdg_allocate_arrays
193 
194  !> @ brief Deallocate package memory
195  !!
196  !! Deallocate SWF package scalars and arrays.
197  !!
198  !<
199  subroutine zdg_da(this)
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
219  end subroutine zdg_da
220 
221  !> @ brief Source additional options for package
222  !!
223  !! Source additional options for SWF package.
224  !!
225  !<
226  subroutine zdg_options(this)
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
248  end subroutine zdg_options
249 
250  !> @ brief Log SWF specific package options
251  !<
252  subroutine log_zdg_options(this, found)
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
275  end subroutine log_zdg_options
276 
277  !> @ brief SWF read and prepare
278  !!
279  !<
280  subroutine zdg_rp(this)
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
298  end subroutine zdg_rp
299 
300  !> @ brief Formulate the package hcof and rhs terms.
301  !!
302  !! Formulate the hcof and rhs terms for the ZDG package that will be
303  !! added to the coefficient matrix and right-hand side vector.
304  !!
305  !<
306  subroutine zdg_cf(this)
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
358  end subroutine zdg_cf
359 
360  !> @brief Calculate conductance-like term
361  !!
362  !! Conductance normally has a dx term in the denominator
363  !! but that is not included here, as the flow is calculated
364  !! using Q = C * slope. The returned c value from this
365  !! function has dimensions of L3/T.
366  !!
367  !<
368  function get_cond(this, i, depth, absdhdxsq, unitconv) result(c)
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 
400  end function get_cond
401 
402  !> @ brief Copy hcof and rhs terms into solution.
403  !!
404  !! Add the hcof and rhs terms for the ZDG package to the
405  !! coefficient matrix and right-hand side vector.
406  !!
407  !<
408  subroutine zdg_fc(this, rhs, ia, idxglo, matrix_sln)
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
441  end subroutine zdg_fc
442 
443  !> @ brief Define the list label for the package
444  !!
445  !! Method defined the list label for the ZDG package. The list label is
446  !! the heading that is written to iout when PRINT_INPUT option is used.
447  !!
448  !<
449  subroutine define_listlabel(this)
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
472  end subroutine define_listlabel
473 
474  ! -- Procedures related to observations
475 
476  !> @brief Determine if observations are supported.
477  !!
478  !! Function to determine if observations are supported by the ZDG package.
479  !! Observations are supported by the ZDG package.
480  !!
481  !! @return zdg_obs_supported boolean indicating if observations are supported
482  !!
483  !<
484  logical function zdg_obs_supported(this)
485  ! -- dummy variables
486  class(swfzdgtype) :: this !< SwfZdgType object
487  !
488  ! -- set boolean
489  zdg_obs_supported = .true.
490  !
491  ! -- return
492  return
493  end function zdg_obs_supported
494 
495  !> @brief Define the observation types available in the package
496  !!
497  !! Method to define the observation types available in the ZDG package.
498  !!
499  !<
500  subroutine zdg_df_obs(this)
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
517  end subroutine zdg_df_obs
518 
519  !> @brief Save observations for the package
520  !!
521  !! Method to save simulated values for the ZDG package.
522  !!
523  !<
524  subroutine zdg_bd_obs(this)
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
567  end subroutine zdg_bd_obs
568 
569  ! -- Procedure related to time series
570 
571  !> @brief Assign time series links for the package
572  !!
573  !! Assign the time series links for the ZDG package. Only
574  !! the Q variable can be defined with time series.
575  !!
576  !<
577  subroutine zdg_rp_ts(this)
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
597  end subroutine zdg_rp_ts
598 
599  !> @ brief Return a bound value
600  !!
601  !! Return a bound value associated with an ncolbnd index
602  !! and row.
603  !!
604  !<
605  function zdg_bound_value(this, col, row) result(bndval)
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
633  end function zdg_bound_value
634 
635 end module swfzdgmodule
This module contains block parser methods.
Definition: BlockParser.f90:7
This module contains the extended boundary package.
This module contains the base boundary package.
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
real(dp), parameter dtwothirds
real constant 2/3
Definition: Constants.f90:69
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:94
real(dp), parameter dem1
real constant 1e-1
Definition: Constants.f90:102
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:67
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
Definition: Constants.f90:38
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
real(dp), parameter done
real constant 1
Definition: Constants.f90:75
integer(i4b) function, public getunit()
Get a free unit number.
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
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:368
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:249
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
real(dp) function sqsaturationderivative(top, bot, x, c1, c2)
@ brief sQSaturationDerivative
real(dp) function sqsaturation(top, bot, x, c1, c2)
@ brief sQSaturation
This module contains the ZDG package methods.
Definition: swf-zdg.f90:7
subroutine zdg_allocate_scalars(this)
@ brief Allocate scalars
Definition: swf-zdg.f90:139
subroutine define_listlabel(this)
@ brief Define the list label for the package
Definition: swf-zdg.f90:450
logical function zdg_obs_supported(this)
Determine if observations are supported.
Definition: swf-zdg.f90:485
subroutine zdg_cf(this)
@ brief Formulate the package hcof and rhs terms.
Definition: swf-zdg.f90:307
subroutine zdg_fc(this, rhs, ia, idxglo, matrix_sln)
@ brief Copy hcof and rhs terms into solution.
Definition: swf-zdg.f90:409
real(dp) function zdg_bound_value(this, col, row)
@ brief Return a bound value
Definition: swf-zdg.f90:606
subroutine log_zdg_options(this, found)
@ brief Log SWF specific package options
Definition: swf-zdg.f90:253
real(dp) function get_cond(this, i, depth, absdhdxsq, unitconv)
Calculate conductance-like term.
Definition: swf-zdg.f90:369
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:79
subroutine zdg_rp_ts(this)
Assign time series links for the package.
Definition: swf-zdg.f90:578
character(len=16) text
package flow text string
Definition: swf-zdg.f90:35
subroutine zdg_da(this)
@ brief Deallocate package memory
Definition: swf-zdg.f90:200
subroutine zdg_options(this)
@ brief Source additional options for package
Definition: swf-zdg.f90:227
subroutine zdg_rp(this)
@ brief SWF read and prepare
Definition: swf-zdg.f90:281
character(len=lenftype) ftype
package ftype
Definition: swf-zdg.f90:34
subroutine zdg_allocate_arrays(this, nodelist, auxvar)
@ brief Allocate arrays
Definition: swf-zdg.f90:163
subroutine zdg_bd_obs(this)
Save observations for the package.
Definition: swf-zdg.f90:525
subroutine zdg_df_obs(this)
Define the observation types available in the package.
Definition: swf-zdg.f90:501
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
type(timeserieslinktype) function, pointer, public gettimeserieslinkfromlist(list, indx)
Get time series link from a list.
@ brief BndType