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