MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
swf-cdb.f90
Go to the documentation of this file.
1 !> @brief This module contains the CDB package methods
2 !!
3 !! This module can be used to represent outflow from streams using
4 !! a critical depth boundary.
5 !!
6 !<
8  ! -- modules
9  use kindmodule, only: dp, i4b
10  use constantsmodule, only: dzero, 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 swfcxsmodule, only: swfcxstype
27  !
28  implicit none
29  !
30  private
31  public :: cdb_create
32  !
33  character(len=LENFTYPE) :: ftype = 'CDB' !< package ftype
34  character(len=16) :: text = ' CDB' !< package flow text string
35  !
36  type, extends(bndexttype) :: swfcdbtype
37 
38  integer(I4B), dimension(:), pointer, contiguous :: idcxs => null() !< cross section id
39  real(dp), dimension(:), pointer, contiguous :: width => null() !< channel width
40  real(dp), pointer :: gravconv => null() !< conversion factor gravity in m/s^2 to model units
41 
42  ! -- pointers other objects
43  type(swfcxstype), pointer :: cxs
44 
45  contains
46  procedure :: allocate_scalars => cdb_allocate_scalars
47  procedure :: allocate_arrays => cdb_allocate_arrays
48  procedure :: source_options => cdb_options
49  procedure :: log_cdb_options
50  procedure :: bnd_rp => cdb_rp
51  procedure :: bnd_cf => cdb_cf
52  procedure :: bnd_fc => cdb_fc
53  procedure :: bnd_da => cdb_da
54  procedure :: define_listlabel
55  procedure :: bound_value => cdb_bound_value
56  ! -- methods for observations
57  procedure, public :: bnd_obs_supported => cdb_obs_supported
58  procedure, public :: bnd_df_obs => cdb_df_obs
59  procedure, public :: bnd_bd_obs => cdb_bd_obs
60  ! -- methods for time series
61  procedure, public :: bnd_rp_ts => cdb_rp_ts
62  ! -- private
63  procedure, private :: qcalc
64  end type swfcdbtype
65 
66 contains
67 
68  !> @ brief Create a new package object
69  !!
70  !! Create a new CDB Package object
71  !!
72  !<
73  subroutine cdb_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
74  mempath, dis, cxs, lengthconv, timeconv)
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 CDB 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) :: lengthconv !< conversion factor from model length to meters
87  real(dp), intent(in) :: timeconv !< conversion factor from model time units to seconds
88  ! -- local variables
89  type(swfcdbtype), pointer :: cdbobj
90  !
91  ! -- allocate the object and assign values to object variables
92  allocate (cdbobj)
93  packobj => cdbobj
94  !
95  ! -- create name and memory path
96  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
97  packobj%text = text
98  !
99  ! -- allocate scalars
100  call cdbobj%allocate_scalars()
101  !
102  ! -- initialize package
103  call packobj%pack_initialize()
104 
105  packobj%inunit = inunit
106  packobj%iout = iout
107  packobj%id = id
108  packobj%ibcnum = ibcnum
109  packobj%ncolbnd = 1
110  packobj%iscloc = 1
111  packobj%ictMemPath = create_mem_path(namemodel, 'DFW')
112 
113  ! -- store pointer to dis
114  cdbobj%dis => dis
115 
116  ! -- store pointer to cxs
117  cdbobj%cxs => cxs
118  !
119  ! -- store unit conversion
120  cdbobj%gravconv = dgravity * lengthconv * timeconv**2
121  end subroutine cdb_create
122 
123  !> @ brief Allocate scalars
124  !!
125  !! Allocate and initialize scalars for the CDB package. The base model
126  !! allocate scalars method is also called.
127  !!
128  !<
129  subroutine cdb_allocate_scalars(this)
130  ! -- modules
132  ! -- dummy variables
133  class(swfcdbtype) :: this !< SwfcdbType object
134  !
135  ! -- call base type allocate scalars
136  call this%BndExtType%allocate_scalars()
137  !
138  ! -- allocate the object and assign values to object variables
139  call mem_allocate(this%gravconv, 'GRAVCONV', this%memoryPath)
140  !
141  ! -- Set values
142  this%gravconv = dzero
143  end subroutine cdb_allocate_scalars
144 
145  !> @ brief Allocate arrays
146  !!
147  !! Allocate and initialize arrays for the SWF package
148  !!
149  !<
150  subroutine cdb_allocate_arrays(this, nodelist, auxvar)
151  ! -- modules
153  ! -- dummy
154  class(swfcdbtype) :: this
155  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
156  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
157  ! -- local
158  !
159  ! -- call BndExtType allocate scalars
160  call this%BndExtType%allocate_arrays(nodelist, auxvar)
161  !
162  ! -- set array input context pointer
163  call mem_setptr(this%idcxs, 'IDCXS', this%input_mempath)
164  call mem_setptr(this%width, 'WIDTH', this%input_mempath)
165  !
166  ! -- checkin array input context pointer
167  call mem_checkin(this%idcxs, 'IDCXS', this%memoryPath, &
168  'IDCXS', this%input_mempath)
169  call mem_checkin(this%width, 'WIDTH', this%memoryPath, &
170  'WIDTH', this%input_mempath)
171  end subroutine cdb_allocate_arrays
172 
173  !> @ brief Deallocate package memory
174  !!
175  !! Deallocate SWF package scalars and arrays.
176  !!
177  !<
178  subroutine cdb_da(this)
179  ! -- modules
181  ! -- dummy variables
182  class(swfcdbtype) :: this !< SwfcdbType object
183  !
184  ! -- Deallocate parent package
185  call this%BndExtType%bnd_da()
186  !
187  ! -- arrays
188  call mem_deallocate(this%idcxs, 'IDCXS', this%memoryPath)
189  call mem_deallocate(this%width, 'WIDTH', this%memoryPath)
190  !
191  ! -- scalars
192  call mem_deallocate(this%gravconv)
193  end subroutine cdb_da
194 
195  !> @ brief Source additional options for package
196  !!
197  !! Source additional options for SWF package.
198  !!
199  !<
200  subroutine cdb_options(this)
201  ! -- modules
202  use inputoutputmodule, only: urword
205  ! -- dummy variables
206  class(swfcdbtype), intent(inout) :: this !< SwfCdbType object
207  ! -- local variables
208  type(swfcdbparamfoundtype) :: found
209  ! -- formats
210  !
211  ! -- source base BndExtType options
212  call this%BndExtType%source_options()
213  !
214  ! -- source options from input context
215  ! none
216  !
217  ! -- log SWF specific options
218  call this%log_cdb_options(found)
219  end subroutine cdb_options
220 
221  !> @ brief Log SWF specific package options
222  !<
223  subroutine log_cdb_options(this, found)
224  ! -- modules
226  ! -- dummy variables
227  class(swfcdbtype), intent(inout) :: this !< BndExtType object
228  type(swfcdbparamfoundtype), intent(in) :: found
229  ! -- local variables
230  ! -- format
231  !
232  ! -- log found options
233  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
234  //' OPTIONS'
235  !
236  ! if (found%mover) then
237  ! write (this%iout, '(4x,A)') 'MOVER OPTION ENABLED'
238  ! end if
239  !
240  ! -- close logging block
241  write (this%iout, '(1x,a)') &
242  'END OF '//trim(adjustl(this%text))//' OPTIONS'
243  end subroutine log_cdb_options
244 
245  !> @ brief SWF read and prepare
246  !!
247  !<
248  subroutine cdb_rp(this)
249  use tdismodule, only: kper
250  ! -- dummy
251  class(swfcdbtype), intent(inout) :: this
252  ! -- local
253  !
254  if (this%iper /= kper) return
255  !
256  ! -- Call the parent class read and prepare
257  call this%BndExtType%bnd_rp()
258  !
259  ! -- Write the list to iout if requested
260  if (this%iprpak /= 0) then
261  call this%write_list()
262  end if
263  end subroutine cdb_rp
264 
265  !> @ brief Formulate the package hcof and rhs terms.
266  !!
267  !! Formulate the hcof and rhs terms for the CDB package that will be
268  !! added to the coefficient matrix and right-hand side vector.
269  !!
270  !<
271  subroutine cdb_cf(this)
272  ! modules
274  ! -- dummy variables
275  class(swfcdbtype) :: this !< SwfCdbType object
276  ! -- local variables
277  integer(I4B) :: i, node
278  real(DP) :: q
279  real(DP) :: qeps
280  real(DP) :: depth
281  real(DP) :: derv
282  real(DP) :: eps
283  !
284  ! -- Return if no inflows
285  if (this%nbound == 0) return
286  !
287  ! -- Calculate hcof and rhs for each cdb entry
288  do i = 1, this%nbound
289 
290  node = this%nodelist(i)
291  if (this%ibound(node) <= 0) then
292  this%hcof(i) = dzero
293  this%rhs(i) = dzero
294  cycle
295  end if
296 
297  ! -- calculate terms and add to hcof and rhs
298  depth = this%xnew(node) - this%dis%bot(node)
299 
300  ! -- calculate unperturbed q
301  q = this%qcalc(i, depth)
302 
303  ! -- calculate perturbed q
304  eps = get_perturbation(depth)
305  qeps = this%qcalc(i, depth + eps)
306 
307  ! -- calculate derivative
308  derv = (qeps - q) / eps
309 
310  ! -- add terms to hcof and rhs
311  this%hcof(i) = derv
312  this%rhs(i) = -q + derv * this%xnew(node)
313 
314  end do
315  end subroutine cdb_cf
316 
317  !> @brief Calculate critical depth boundary flow
318  !<
319  function qcalc(this, i, depth) result(q)
320  ! modules
321  ! dummy
322  class(swfcdbtype) :: this
323  integer(I4B), intent(in) :: i !< boundary number
324  real(dp), intent(in) :: depth !< simulated depth (stage - elevation) in reach n for this iteration
325  ! return
326  real(dp) :: q
327  ! local
328  integer(I4B) :: idcxs
329  real(dp) :: width
330  real(dp) :: a
331  real(dp) :: r
332 
333  idcxs = this%idcxs(i)
334  width = this%width(i)
335  a = this%cxs%get_area(idcxs, width, depth)
336  r = this%cxs%get_hydraulic_radius(idcxs, width, depth, area=a)
337 
338  q = this%gravconv * a**dtwo * r
339  if (q > dprec) then
340  q = q**dhalf
341  else
342  q = dzero
343  end if
344  q = -q
345 
346  end function qcalc
347 
348  !> @ brief Copy hcof and rhs terms into solution.
349  !!
350  !! Add the hcof and rhs terms for the CDB package to the
351  !! coefficient matrix and right-hand side vector.
352  !!
353  !<
354  subroutine cdb_fc(this, rhs, ia, idxglo, matrix_sln)
355  ! -- dummy variables
356  class(swfcdbtype) :: this !< SwfCdbType object
357  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
358  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
359  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
360  class(matrixbasetype), pointer :: matrix_sln !< solution coefficient matrix
361  ! -- local variables
362  integer(I4B) :: i
363  integer(I4B) :: n
364  integer(I4B) :: ipos
365  !
366  ! -- pakmvrobj fc
367  if (this%imover == 1) then
368  call this%pakmvrobj%fc()
369  end if
370  !
371  ! -- Copy package rhs and hcof into solution rhs and amat
372  do i = 1, this%nbound
373  n = this%nodelist(i)
374  rhs(n) = rhs(n) + this%rhs(i)
375  ipos = ia(n)
376  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
377  !
378  ! -- If mover is active and this cdb item is discharging,
379  ! store available water (as positive value).
380  if (this%imover == 1 .and. this%rhs(i) > dzero) then
381  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
382  end if
383  end do
384  end subroutine cdb_fc
385 
386  !> @ brief Define the list label for the package
387  !!
388  !! Method defined the list label for the CDB package. The list label is
389  !! the heading that is written to iout when PRINT_INPUT option is used.
390  !!
391  !<
392  subroutine define_listlabel(this)
393  ! -- dummy variables
394  class(swfcdbtype), intent(inout) :: this !< SwfCdbType object
395  !
396  ! -- create the header list label
397  this%listlabel = trim(this%filtyp)//' NO.'
398  if (this%dis%ndim == 3) then
399  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
400  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
401  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
402  elseif (this%dis%ndim == 2) then
403  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
404  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
405  else
406  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
407  end if
408  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'FLOW RATE'
409  if (this%inamedbound == 1) then
410  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
411  end if
412  end subroutine define_listlabel
413 
414  ! -- Procedures related to observations
415 
416  !> @brief Determine if observations are supported.
417  !!
418  !! Function to determine if observations are supported by the CDB package.
419  !! Observations are supported by the CDB package.
420  !!
421  !! @return cdb_obs_supported boolean indicating if observations are supported
422  !!
423  !<
424  logical function cdb_obs_supported(this)
425  ! -- dummy variables
426  class(swfcdbtype) :: this !< SwfCdbType object
427  !
428  ! -- set boolean
429  cdb_obs_supported = .true.
430  end function cdb_obs_supported
431 
432  !> @brief Define the observation types available in the package
433  !!
434  !! Method to define the observation types available in the CDB package.
435  !!
436  !<
437  subroutine cdb_df_obs(this)
438  ! -- dummy variables
439  class(swfcdbtype) :: this !< SwfCdbType object
440  ! -- local variables
441  integer(I4B) :: indx
442  !
443  ! -- initialize observations
444  call this%obs%StoreObsType('cdb', .true., indx)
445  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
446  !
447  ! -- Store obs type and assign procedure pointer
448  ! for to-mvr observation type.
449  call this%obs%StoreObsType('to-mvr', .true., indx)
450  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
451  end subroutine cdb_df_obs
452 
453  !> @brief Save observations for the package
454  !!
455  !! Method to save simulated values for the CDB package.
456  !!
457  !<
458  subroutine cdb_bd_obs(this)
459  ! -- dummy variables
460  class(swfcdbtype) :: this !< SwfCdbType object
461  ! -- local variables
462  integer(I4B) :: i
463  integer(I4B) :: n
464  integer(I4B) :: jj
465  real(DP) :: v
466  type(observetype), pointer :: obsrv => null()
467  !
468  ! -- clear the observations
469  call this%obs%obs_bd_clear()
470  !
471  ! -- Save simulated values for all of package's observations.
472  do i = 1, this%obs%npakobs
473  obsrv => this%obs%pakobs(i)%obsrv
474  if (obsrv%BndFound) then
475  do n = 1, obsrv%indxbnds_count
476  v = dnodata
477  jj = obsrv%indxbnds(n)
478  select case (obsrv%ObsTypeId)
479  case ('TO-MVR')
480  if (this%imover == 1) then
481  v = this%pakmvrobj%get_qtomvr(jj)
482  if (v > dzero) then
483  v = -v
484  end if
485  end if
486  case ('CDB')
487  v = this%simvals(jj)
488  case default
489  errmsg = 'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
490  call store_error(errmsg)
491  end select
492  call this%obs%SaveOneSimval(obsrv, v)
493  end do
494  else
495  call this%obs%SaveOneSimval(obsrv, dnodata)
496  end if
497  end do
498  end subroutine cdb_bd_obs
499 
500  ! -- Procedure related to time series
501 
502  !> @brief Assign time series links for the package
503  !!
504  !! Assign the time series links for the CDB package. Only
505  !! the Q variable can be defined with time series.
506  !!
507  !<
508  subroutine cdb_rp_ts(this)
509  ! -- dummy variables
510  class(swfcdbtype), intent(inout) :: this !< SwfCdbType object
511  ! -- local variables
512  integer(I4B) :: i, nlinks
513  type(timeserieslinktype), pointer :: tslink => null()
514  !
515  ! -- set up the time series links
516  nlinks = this%TsManager%boundtslinks%Count()
517  do i = 1, nlinks
518  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
519  if (associated(tslink)) then
520  if (tslink%JCol == 1) then
521  tslink%Text = 'Q'
522  end if
523  end if
524  end do
525  end subroutine cdb_rp_ts
526 
527  !> @ brief Return a bound value
528  !!
529  !! Return a bound value associated with an ncolbnd index
530  !! and row.
531  !!
532  !<
533  function cdb_bound_value(this, col, row) result(bndval)
534  ! -- modules
535  use constantsmodule, only: dzero
536  ! -- dummy variables
537  class(swfcdbtype), intent(inout) :: this !< BndExtType object
538  integer(I4B), intent(in) :: col
539  integer(I4B), intent(in) :: row
540  ! -- result
541  real(dp) :: bndval
542  !
543  select case (col)
544  case (1)
545  bndval = this%idcxs(row)
546  case (2)
547  bndval = this%width(row)
548  case default
549  errmsg = 'Programming error. CDB bound value requested column '&
550  &'outside range of ncolbnd (1).'
551  call store_error(errmsg)
552  call store_error_filename(this%input_fname)
553  end select
554  end function cdb_bound_value
555 
556 end module swfcdbmodule
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
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:95
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:68
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
Definition: Constants.f90:39
real(dp), parameter dgravity
real constant gravitational acceleration (m/(s s))
Definition: Constants.f90:132
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter dprec
real constant machine precision
Definition: Constants.f90:120
real(dp), parameter dtwo
real constant 2
Definition: Constants.f90:79
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:372
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:246
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 CDB package methods.
Definition: swf-cdb.f90:7
subroutine cdb_allocate_arrays(this, nodelist, auxvar)
@ brief Allocate arrays
Definition: swf-cdb.f90:151
subroutine cdb_rp(this)
@ brief SWF read and prepare
Definition: swf-cdb.f90:249
subroutine, public cdb_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath, dis, cxs, lengthconv, timeconv)
@ brief Create a new package object
Definition: swf-cdb.f90:75
real(dp) function qcalc(this, i, depth)
Calculate critical depth boundary flow.
Definition: swf-cdb.f90:320
logical function cdb_obs_supported(this)
Determine if observations are supported.
Definition: swf-cdb.f90:425
subroutine log_cdb_options(this, found)
@ brief Log SWF specific package options
Definition: swf-cdb.f90:224
subroutine cdb_bd_obs(this)
Save observations for the package.
Definition: swf-cdb.f90:459
character(len=16) text
package flow text string
Definition: swf-cdb.f90:34
character(len=lenftype) ftype
package ftype
Definition: swf-cdb.f90:33
subroutine cdb_rp_ts(this)
Assign time series links for the package.
Definition: swf-cdb.f90:509
subroutine cdb_da(this)
@ brief Deallocate package memory
Definition: swf-cdb.f90:179
subroutine cdb_allocate_scalars(this)
@ brief Allocate scalars
Definition: swf-cdb.f90:130
subroutine cdb_options(this)
@ brief Source additional options for package
Definition: swf-cdb.f90:201
subroutine cdb_cf(this)
@ brief Formulate the package hcof and rhs terms.
Definition: swf-cdb.f90:272
subroutine cdb_df_obs(this)
Define the observation types available in the package.
Definition: swf-cdb.f90:438
subroutine cdb_fc(this, rhs, ia, idxglo, matrix_sln)
@ brief Copy hcof and rhs terms into solution.
Definition: swf-cdb.f90:355
subroutine define_listlabel(this)
@ brief Define the list label for the package
Definition: swf-cdb.f90:393
real(dp) function cdb_bound_value(this, col, row)
@ brief Return a bound value
Definition: swf-cdb.f90:534
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