MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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  !
122  ! -- return
123  return
124  end subroutine cdb_create
125 
126  !> @ brief Allocate scalars
127  !!
128  !! Allocate and initialize scalars for the CDB package. The base model
129  !! allocate scalars method is also called.
130  !!
131  !<
132  subroutine cdb_allocate_scalars(this)
133  ! -- modules
135  ! -- dummy variables
136  class(swfcdbtype) :: this !< SwfcdbType object
137  !
138  ! -- call base type allocate scalars
139  call this%BndExtType%allocate_scalars()
140  !
141  ! -- allocate the object and assign values to object variables
142  call mem_allocate(this%gravconv, 'GRAVCONV', this%memoryPath)
143  !
144  ! -- Set values
145  this%gravconv = dzero
146  !
147  ! -- return
148  return
149  end subroutine cdb_allocate_scalars
150 
151  !> @ brief Allocate arrays
152  !!
153  !! Allocate and initialize arrays for the SWF package
154  !!
155  !<
156  subroutine cdb_allocate_arrays(this, nodelist, auxvar)
157  ! -- modules
159  ! -- dummy
160  class(swfcdbtype) :: this
161  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
162  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
163  ! -- local
164  !
165  ! -- call BndExtType allocate scalars
166  call this%BndExtType%allocate_arrays(nodelist, auxvar)
167  !
168  ! -- set array input context pointer
169  call mem_setptr(this%idcxs, 'IDCXS', this%input_mempath)
170  call mem_setptr(this%width, 'WIDTH', this%input_mempath)
171  !
172  ! -- checkin array input context pointer
173  call mem_checkin(this%idcxs, 'IDCXS', this%memoryPath, &
174  'IDCXS', this%input_mempath)
175  call mem_checkin(this%width, 'WIDTH', this%memoryPath, &
176  'WIDTH', this%input_mempath)
177  !
178  ! -- return
179  return
180  end subroutine cdb_allocate_arrays
181 
182  !> @ brief Deallocate package memory
183  !!
184  !! Deallocate SWF package scalars and arrays.
185  !!
186  !<
187  subroutine cdb_da(this)
188  ! -- modules
190  ! -- dummy variables
191  class(swfcdbtype) :: this !< SwfcdbType object
192  !
193  ! -- Deallocate parent package
194  call this%BndExtType%bnd_da()
195  !
196  ! -- arrays
197  call mem_deallocate(this%idcxs, 'IDCXS', this%memoryPath)
198  call mem_deallocate(this%width, 'WIDTH', this%memoryPath)
199  !
200  ! -- scalars
201  call mem_deallocate(this%gravconv)
202  !
203  ! -- return
204  return
205  end subroutine cdb_da
206 
207  !> @ brief Source additional options for package
208  !!
209  !! Source additional options for SWF package.
210  !!
211  !<
212  subroutine cdb_options(this)
213  ! -- modules
214  use inputoutputmodule, only: urword
217  ! -- dummy variables
218  class(swfcdbtype), intent(inout) :: this !< SwfCdbType object
219  ! -- local variables
220  type(swfcdbparamfoundtype) :: found
221  ! -- formats
222  !
223  ! -- source base BndExtType options
224  call this%BndExtType%source_options()
225  !
226  ! -- source options from input context
227  ! none
228  !
229  ! -- log SWF specific options
230  call this%log_cdb_options(found)
231  !
232  ! -- return
233  return
234  end subroutine cdb_options
235 
236  !> @ brief Log SWF specific package options
237  !<
238  subroutine log_cdb_options(this, found)
239  ! -- modules
241  ! -- dummy variables
242  class(swfcdbtype), intent(inout) :: this !< BndExtType object
243  type(swfcdbparamfoundtype), intent(in) :: found
244  ! -- local variables
245  ! -- format
246  !
247  ! -- log found options
248  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
249  //' OPTIONS'
250  !
251  ! if (found%mover) then
252  ! write (this%iout, '(4x,A)') 'MOVER OPTION ENABLED'
253  ! end if
254  !
255  ! -- close logging block
256  write (this%iout, '(1x,a)') &
257  'END OF '//trim(adjustl(this%text))//' OPTIONS'
258  !
259  ! -- return
260  return
261  end subroutine log_cdb_options
262 
263  !> @ brief SWF read and prepare
264  !!
265  !<
266  subroutine cdb_rp(this)
267  use tdismodule, only: kper
268  ! -- dummy
269  class(swfcdbtype), intent(inout) :: this
270  ! -- local
271  !
272  if (this%iper /= kper) return
273  !
274  ! -- Call the parent class read and prepare
275  call this%BndExtType%bnd_rp()
276  !
277  ! -- Write the list to iout if requested
278  if (this%iprpak /= 0) then
279  call this%write_list()
280  end if
281  !
282  ! -- return
283  return
284  end subroutine cdb_rp
285 
286  !> @ brief Formulate the package hcof and rhs terms.
287  !!
288  !! Formulate the hcof and rhs terms for the CDB package that will be
289  !! added to the coefficient matrix and right-hand side vector.
290  !!
291  !<
292  subroutine cdb_cf(this)
293  ! modules
295  ! -- dummy variables
296  class(swfcdbtype) :: this !< SwfCdbType object
297  ! -- local variables
298  integer(I4B) :: i, node
299  real(DP) :: q
300  real(DP) :: qeps
301  real(DP) :: depth
302  real(DP) :: derv
303  real(DP) :: eps
304  !
305  ! -- Return if no inflows
306  if (this%nbound == 0) return
307  !
308  ! -- Calculate hcof and rhs for each cdb entry
309  do i = 1, this%nbound
310 
311  node = this%nodelist(i)
312  if (this%ibound(node) <= 0) then
313  this%hcof(i) = dzero
314  this%rhs(i) = dzero
315  cycle
316  end if
317 
318  ! -- calculate terms and add to hcof and rhs
319  depth = this%xnew(node) - this%dis%bot(node)
320 
321  ! -- calculate unperturbed q
322  q = this%qcalc(i, depth)
323 
324  ! -- calculate perturbed q
325  eps = get_perturbation(depth)
326  qeps = this%qcalc(i, depth + eps)
327 
328  ! -- calculate derivative
329  derv = (qeps - q) / eps
330 
331  ! -- add terms to hcof and rhs
332  this%hcof(i) = derv
333  this%rhs(i) = -q + derv * this%xnew(node)
334 
335  end do
336  !
337  return
338  end subroutine cdb_cf
339 
340  !> @brief Calculate critical depth boundary flow
341  !<
342  function qcalc(this, i, depth) result(q)
343  ! modules
344  ! dummy
345  class(swfcdbtype) :: this
346  integer(I4B), intent(in) :: i !< boundary number
347  real(dp), intent(in) :: depth !< simulated depth (stage - elevation) in reach n for this iteration
348  ! return
349  real(dp) :: q
350  ! local
351  integer(I4B) :: idcxs
352  real(dp) :: width
353  real(dp) :: a
354  real(dp) :: r
355 
356  idcxs = this%idcxs(i)
357  width = this%width(i)
358  a = this%cxs%get_area(idcxs, width, depth)
359  r = this%cxs%get_hydraulic_radius(idcxs, width, depth, area=a)
360 
361  q = this%gravconv * a**dtwo * r
362  if (q > dprec) then
363  q = q**dhalf
364  else
365  q = dzero
366  end if
367  q = -q
368 
369  end function qcalc
370 
371  !> @ brief Copy hcof and rhs terms into solution.
372  !!
373  !! Add the hcof and rhs terms for the CDB package to the
374  !! coefficient matrix and right-hand side vector.
375  !!
376  !<
377  subroutine cdb_fc(this, rhs, ia, idxglo, matrix_sln)
378  ! -- dummy variables
379  class(swfcdbtype) :: this !< SwfCdbType object
380  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
381  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
382  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
383  class(matrixbasetype), pointer :: matrix_sln !< solution coefficient matrix
384  ! -- local variables
385  integer(I4B) :: i
386  integer(I4B) :: n
387  integer(I4B) :: ipos
388  !
389  ! -- pakmvrobj fc
390  if (this%imover == 1) then
391  call this%pakmvrobj%fc()
392  end if
393  !
394  ! -- Copy package rhs and hcof into solution rhs and amat
395  do i = 1, this%nbound
396  n = this%nodelist(i)
397  rhs(n) = rhs(n) + this%rhs(i)
398  ipos = ia(n)
399  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
400  !
401  ! -- If mover is active and this cdb item is discharging,
402  ! store available water (as positive value).
403  if (this%imover == 1 .and. this%rhs(i) > dzero) then
404  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
405  end if
406  end do
407  !
408  ! -- return
409  return
410  end subroutine cdb_fc
411 
412  !> @ brief Define the list label for the package
413  !!
414  !! Method defined the list label for the CDB package. The list label is
415  !! the heading that is written to iout when PRINT_INPUT option is used.
416  !!
417  !<
418  subroutine define_listlabel(this)
419  ! -- dummy variables
420  class(swfcdbtype), intent(inout) :: this !< SwfCdbType object
421  !
422  ! -- create the header list label
423  this%listlabel = trim(this%filtyp)//' NO.'
424  if (this%dis%ndim == 3) then
425  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
426  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
427  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
428  elseif (this%dis%ndim == 2) then
429  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
430  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
431  else
432  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
433  end if
434  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'FLOW RATE'
435  if (this%inamedbound == 1) then
436  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
437  end if
438  !
439  ! -- return
440  return
441  end subroutine define_listlabel
442 
443  ! -- Procedures related to observations
444 
445  !> @brief Determine if observations are supported.
446  !!
447  !! Function to determine if observations are supported by the CDB package.
448  !! Observations are supported by the CDB package.
449  !!
450  !! @return cdb_obs_supported boolean indicating if observations are supported
451  !!
452  !<
453  logical function cdb_obs_supported(this)
454  ! -- dummy variables
455  class(swfcdbtype) :: this !< SwfCdbType object
456  !
457  ! -- set boolean
458  cdb_obs_supported = .true.
459  !
460  ! -- return
461  return
462  end function cdb_obs_supported
463 
464  !> @brief Define the observation types available in the package
465  !!
466  !! Method to define the observation types available in the CDB package.
467  !!
468  !<
469  subroutine cdb_df_obs(this)
470  ! -- dummy variables
471  class(swfcdbtype) :: this !< SwfCdbType object
472  ! -- local variables
473  integer(I4B) :: indx
474  !
475  ! -- initialize observations
476  call this%obs%StoreObsType('cdb', .true., indx)
477  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
478  !
479  ! -- Store obs type and assign procedure pointer
480  ! for to-mvr observation type.
481  call this%obs%StoreObsType('to-mvr', .true., indx)
482  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
483  !
484  ! -- return
485  return
486  end subroutine cdb_df_obs
487 
488  !> @brief Save observations for the package
489  !!
490  !! Method to save simulated values for the CDB package.
491  !!
492  !<
493  subroutine cdb_bd_obs(this)
494  ! -- dummy variables
495  class(swfcdbtype) :: this !< SwfCdbType object
496  ! -- local variables
497  integer(I4B) :: i
498  integer(I4B) :: n
499  integer(I4B) :: jj
500  real(DP) :: v
501  type(observetype), pointer :: obsrv => null()
502  !
503  ! -- clear the observations
504  call this%obs%obs_bd_clear()
505  !
506  ! -- Save simulated values for all of package's observations.
507  do i = 1, this%obs%npakobs
508  obsrv => this%obs%pakobs(i)%obsrv
509  if (obsrv%BndFound) then
510  do n = 1, obsrv%indxbnds_count
511  v = dnodata
512  jj = obsrv%indxbnds(n)
513  select case (obsrv%ObsTypeId)
514  case ('TO-MVR')
515  if (this%imover == 1) then
516  v = this%pakmvrobj%get_qtomvr(jj)
517  if (v > dzero) then
518  v = -v
519  end if
520  end if
521  case ('CDB')
522  v = this%simvals(jj)
523  case default
524  errmsg = 'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
525  call store_error(errmsg)
526  end select
527  call this%obs%SaveOneSimval(obsrv, v)
528  end do
529  else
530  call this%obs%SaveOneSimval(obsrv, dnodata)
531  end if
532  end do
533  !
534  ! -- return
535  return
536  end subroutine cdb_bd_obs
537 
538  ! -- Procedure related to time series
539 
540  !> @brief Assign time series links for the package
541  !!
542  !! Assign the time series links for the CDB package. Only
543  !! the Q variable can be defined with time series.
544  !!
545  !<
546  subroutine cdb_rp_ts(this)
547  ! -- dummy variables
548  class(swfcdbtype), intent(inout) :: this !< SwfCdbType object
549  ! -- local variables
550  integer(I4B) :: i, nlinks
551  type(timeserieslinktype), pointer :: tslink => null()
552  !
553  ! -- set up the time series links
554  nlinks = this%TsManager%boundtslinks%Count()
555  do i = 1, nlinks
556  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
557  if (associated(tslink)) then
558  if (tslink%JCol == 1) then
559  tslink%Text = 'Q'
560  end if
561  end if
562  end do
563  !
564  ! -- return
565  return
566  end subroutine cdb_rp_ts
567 
568  !> @ brief Return a bound value
569  !!
570  !! Return a bound value associated with an ncolbnd index
571  !! and row.
572  !!
573  !<
574  function cdb_bound_value(this, col, row) result(bndval)
575  ! -- modules
576  use constantsmodule, only: dzero
577  ! -- dummy variables
578  class(swfcdbtype), intent(inout) :: this !< BndExtType object
579  integer(I4B), intent(in) :: col
580  integer(I4B), intent(in) :: row
581  ! -- result
582  real(dp) :: bndval
583  !
584  select case (col)
585  case (1)
586  bndval = this%idcxs(row)
587  case (2)
588  bndval = this%width(row)
589  case default
590  errmsg = 'Programming error. CDB bound value requested column '&
591  &'outside range of ncolbnd (1).'
592  call store_error(errmsg)
593  call store_error_filename(this%input_fname)
594  end select
595  !
596  ! -- return
597  return
598  end function cdb_bound_value
599 
600 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:94
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 dgravity
real constant gravitational acceleration (m/(s s))
Definition: Constants.f90:131
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
real(dp), parameter dprec
real constant machine precision
Definition: Constants.f90:119
real(dp), parameter dtwo
real constant 2
Definition: Constants.f90:78
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 CDB package methods.
Definition: swf-cdb.f90:7
subroutine cdb_allocate_arrays(this, nodelist, auxvar)
@ brief Allocate arrays
Definition: swf-cdb.f90:157
subroutine cdb_rp(this)
@ brief SWF read and prepare
Definition: swf-cdb.f90:267
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:343
logical function cdb_obs_supported(this)
Determine if observations are supported.
Definition: swf-cdb.f90:454
subroutine log_cdb_options(this, found)
@ brief Log SWF specific package options
Definition: swf-cdb.f90:239
subroutine cdb_bd_obs(this)
Save observations for the package.
Definition: swf-cdb.f90:494
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:547
subroutine cdb_da(this)
@ brief Deallocate package memory
Definition: swf-cdb.f90:188
subroutine cdb_allocate_scalars(this)
@ brief Allocate scalars
Definition: swf-cdb.f90:133
subroutine cdb_options(this)
@ brief Source additional options for package
Definition: swf-cdb.f90:213
subroutine cdb_cf(this)
@ brief Formulate the package hcof and rhs terms.
Definition: swf-cdb.f90:293
subroutine cdb_df_obs(this)
Define the observation types available in the package.
Definition: swf-cdb.f90:470
subroutine cdb_fc(this, rhs, ia, idxglo, matrix_sln)
@ brief Copy hcof and rhs terms into solution.
Definition: swf-cdb.f90:378
subroutine define_listlabel(this)
@ brief Define the list label for the package
Definition: swf-cdb.f90:419
real(dp) function cdb_bound_value(this, col, row)
@ brief Return a bound value
Definition: swf-cdb.f90:575
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