MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
chdmodule Module Reference

Data Types

type  chdtype
 

Functions/Subroutines

subroutine, public chd_create (packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
 Create a new constant head package. More...
 
subroutine chd_allocate_arrays (this, nodelist, auxvar)
 Allocate arrays specific to the constant head package. More...
 
subroutine chd_rp (this)
 Constant concentration/temperature read and prepare (rp) routine. More...
 
subroutine chd_ad (this)
 Constant head package advance routine. More...
 
subroutine chd_ck (this)
 Check constant concentration/temperature boundary condition data. More...
 
subroutine chd_fc (this, rhs, ia, idxglo, matrix_sln)
 Override bnd_fc and do nothing. More...
 
subroutine chd_cq (this, x, flowja, iadv)
 Calculate flow associated with constant head boundary. More...
 
subroutine calc_chd_rate (this)
 Calculate the CHD cell rates, to be called. More...
 
subroutine chd_bd (this, model_budget)
 Add package ratin/ratout to model budget. More...
 
subroutine chd_da (this)
 Deallocate memory. More...
 
subroutine define_listlabel (this)
 Define the list heading that is written to iout when PRINT_INPUT option is used. More...
 
logical function chd_obs_supported (this)
 Overrides bnd_obs_supported from bndType class. More...
 
subroutine chd_df_obs (this)
 Overrides bnd_df_obs from bndType class. More...
 
real(dp) function head_mult (this, row)
 Apply auxiliary multiplier to specified head if appropriate. More...
 
real(dp) function chd_bound_value (this, col, row)
 @ brief Return a bound value More...
 

Variables

character(len=lenftype) ftype = 'CHD'
 
character(len=lenpackagename) text = ' CHD'
 

Function/Subroutine Documentation

◆ calc_chd_rate()

subroutine chdmodule::calc_chd_rate ( class(chdtype), intent(inout)  this)
private

Definition at line 260 of file gwf-chd.f90.

261  ! -- modules
262  ! -- dummy
263  class(ChdType), intent(inout) :: this
264  ! -- local
265  integer(I4B) :: i
266  integer(I4B) :: ipos
267  integer(I4B) :: node
268  integer(I4B) :: n2
269  integer(I4B) :: idiag
270  real(DP) :: rate
271  real(DP) :: ratein, rateout
272  real(DP) :: q
273  !
274  ! -- If no boundaries, skip flow calculations.
275  if (this%nbound > 0) then
276  !
277  ! -- Loop through each boundary calculating flow.
278  do i = 1, this%nbound
279  node = this%nodelist(i)
280  idiag = this%dis%con%ia(node)
281  rate = dzero
282  ratein = dzero
283  rateout = dzero
284  !
285  ! -- Calculate the flow rate into the cell.
286  do ipos = this%dis%con%ia(node) + 1, &
287  this%dis%con%ia(node + 1) - 1
288  q = this%flowja(ipos)
289  rate = rate - q
290  ! -- only accumulate chin and chout for active
291  ! connected cells
292  n2 = this%dis%con%ja(ipos)
293  if (this%ibound(n2) > 0) then
294  if (q < dzero) then
295  ratein = ratein - q
296  else
297  rateout = rateout + q
298  end if
299  end if
300  end do
301  !
302  ! -- For chd, store total flow in rhs so it is available for other
303  ! calculations
304  this%rhs(i) = -rate
305  this%hcof(i) = dzero
306  !
307  ! -- Save simulated value to simvals array.
308  this%simvals(i) = rate
309  this%ratechdin(i) = ratein
310  this%ratechdout(i) = rateout
311  this%flowja(idiag) = this%flowja(idiag) + rate
312  !
313  end do
314  !
315  end if

◆ chd_ad()

subroutine chdmodule::chd_ad ( class(chdtype this)

Add package connections to matrix

Definition at line 171 of file gwf-chd.f90.

172  ! -- modules
173  ! -- dummy
174  class(ChdType) :: this
175  ! -- local
176  integer(I4B) :: i, node
177  real(DP) :: hb
178  ! -- formats
179  !
180  ! -- Process each entry in the specified-head cell list
181  do i = 1, this%nbound
182  node = this%nodelist(i)
183  hb = this%head_mult(i)
184  !
185  this%xnew(node) = hb
186  this%xold(node) = this%xnew(node)
187  end do
188  !
189  ! -- For each observation, push simulated value and corresponding
190  ! simulation time from "current" to "preceding" and reset
191  ! "current" value.
192  call this%obs%obs_ad()

◆ chd_allocate_arrays()

subroutine chdmodule::chd_allocate_arrays ( class(chdtype this,
integer(i4b), dimension(:), optional, pointer, contiguous  nodelist,
real(dp), dimension(:, :), optional, pointer, contiguous  auxvar 
)
private

Definition at line 93 of file gwf-chd.f90.

94  ! -- modules
96  ! -- dummy
97  class(ChdType) :: this
98  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
99  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
100  ! -- local
101  integer(I4B) :: i
102  !
103  ! -- call standard BndType allocate scalars
104  call this%BndExtType%allocate_arrays(nodelist, auxvar)
105  !
106  ! -- allocate ratechdex
107  call mem_allocate(this%ratechdin, this%maxbound, 'RATECHDIN', this%memoryPath)
108  call mem_allocate(this%ratechdout, this%maxbound, 'RATECHDOUT', &
109  this%memoryPath)
110  do i = 1, this%maxbound
111  this%ratechdin(i) = dzero
112  this%ratechdout(i) = dzero
113  end do
114  !
115  ! -- set constant head array input context pointer
116  call mem_setptr(this%head, 'HEAD', this%input_mempath)
117  !
118  ! -- checkin constant head array input context pointer
119  call mem_checkin(this%head, 'HEAD', this%memoryPath, &
120  'HEAD', this%input_mempath)

◆ chd_bd()

subroutine chdmodule::chd_bd ( class(chdtype this,
type(budgettype), intent(inout)  model_budget 
)
private

Definition at line 320 of file gwf-chd.f90.

321  ! -- add package ratin/ratout to model budget
322  use tdismodule, only: delt
324  class(ChdType) :: this
325  type(BudgetType), intent(inout) :: model_budget
326  real(DP) :: ratin
327  real(DP) :: ratout
328  real(DP) :: dum
329  integer(I4B) :: isuppress_output
330 
331  ! For CHDs at an exchange, under some conditions
332  ! (XT3D), the model flowja into the cell is not
333  ! finalized until after exg_cq. So we calculate
334  ! the CHD rate here
335  call this%calc_chd_rate()
336 
337  isuppress_output = 0
338  call rate_accumulator(this%ratechdin(1:this%nbound), ratin, dum)
339  call rate_accumulator(this%ratechdout(1:this%nbound), ratout, dum)
340  call model_budget%addentry(ratin, ratout, delt, this%text, &
341  isuppress_output, this%packName)
This module contains the BudgetModule.
Definition: Budget.f90:20
subroutine, public rate_accumulator(flow, rin, rout)
@ brief Rate accumulator subroutine
Definition: Budget.f90:632
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29
Derived type for the Budget object.
Definition: Budget.f90:39
Here is the call graph for this function:

◆ chd_bound_value()

real(dp) function chdmodule::chd_bound_value ( class(chdtype), intent(inout)  this,
integer(i4b), intent(in)  col,
integer(i4b), intent(in)  row 
)

Return a bound value associated with an ncolbnd index and row.

Definition at line 438 of file gwf-chd.f90.

439  ! -- modules
440  use constantsmodule, only: dzero
441  ! -- dummy variables
442  class(ChdType), intent(inout) :: this
443  integer(I4B), intent(in) :: col
444  integer(I4B), intent(in) :: row
445  ! -- result
446  real(DP) :: bndval
447  !
448  select case (col)
449  case (1)
450  bndval = this%head_mult(row)
451  case default
452  errmsg = 'Programming error. CHD bound value requested column '&
453  &'outside range of ncolbnd (1).'
454  call store_error(errmsg)
455  call store_error_filename(this%input_fname)
456  end select
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
Here is the call graph for this function:

◆ chd_ck()

subroutine chdmodule::chd_ck ( class(chdtype), intent(inout)  this)
private

Definition at line 197 of file gwf-chd.f90.

198  ! -- modules
199  ! -- dummy
200  class(ChdType), intent(inout) :: this
201  ! -- local
202  character(len=30) :: nodestr
203  integer(I4B) :: i
204  integer(I4B) :: node
205  real(DP) :: bt
206  ! -- formats
207  character(len=*), parameter :: fmtchderr = &
208  "('CHD BOUNDARY ',i0,' HEAD (',g0,') IS LESS THAN CELL &
209  &BOTTOM (',g0,')',' FOR CELL ',a)"
210  !
211  ! -- check stress period data
212  do i = 1, this%nbound
213  node = this%nodelist(i)
214  bt = this%dis%bot(node)
215  ! -- accumulate errors
216  if (this%head_mult(i) < bt .and. this%icelltype(node) /= 0) then
217  call this%dis%noder_to_string(node, nodestr)
218  write (errmsg, fmt=fmtchderr) i, this%head_mult(i), bt, trim(nodestr)
219  call store_error(errmsg)
220  end if
221  end do
222  !
223  ! write summary of chd package error messages
224  if (count_errors() > 0) then
225  call store_error_filename(this%input_fname)
226  end if
Here is the call graph for this function:

◆ chd_cq()

subroutine chdmodule::chd_cq ( class(chdtype), intent(inout)  this,
real(dp), dimension(:), intent(in)  x,
real(dp), dimension(:), intent(inout), contiguous  flowja,
integer(i4b), intent(in), optional  iadv 
)
private

This method overrides bnd_cq()

Definition at line 248 of file gwf-chd.f90.

249  class(ChdType), intent(inout) :: this
250  real(DP), dimension(:), intent(in) :: x
251  real(DP), dimension(:), contiguous, intent(inout) :: flowja
252  integer(I4B), optional, intent(in) :: iadv
253 
254  ! NB: the rate calculation cannot be done until chd_bd below
255 

◆ chd_create()

subroutine, public chdmodule::chd_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 
)

Routine points packobj to the newly created package

Definition at line 54 of file gwf-chd.f90.

56  ! -- dummy
57  class(BndType), pointer :: packobj
58  integer(I4B), intent(in) :: id
59  integer(I4B), intent(in) :: ibcnum
60  integer(I4B), intent(in) :: inunit
61  integer(I4B), intent(in) :: iout
62  character(len=*), intent(in) :: namemodel
63  character(len=*), intent(in) :: pakname
64  character(len=*), intent(in) :: mempath
65  ! -- local
66  type(ChdType), pointer :: chdobj
67  !
68  ! -- allocate the object and assign values to object variables
69  allocate (chdobj)
70  packobj => chdobj
71  !
72  ! -- create name and memory path
73  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
74  packobj%text = text
75  !
76  ! -- allocate scalars
77  call chdobj%allocate_scalars()
78  !
79  ! -- initialize package
80  call packobj%pack_initialize()
81  !
82  ! -- store values
83  packobj%inunit = inunit
84  packobj%iout = iout
85  packobj%id = id
86  packobj%ibcnum = ibcnum
87  packobj%ncolbnd = 1
88  packobj%ictMemPath = create_mem_path(namemodel, 'NPF')
Here is the call graph for this function:
Here is the caller graph for this function:

◆ chd_da()

subroutine chdmodule::chd_da ( class(chdtype this)

Definition at line 346 of file gwf-chd.f90.

347  ! -- modules
349  ! -- dummy
350  class(ChdType) :: this
351  !
352  ! -- Deallocate parent package
353  call this%BndExtType%bnd_da()
354  !
355  ! -- arrays
356  call mem_deallocate(this%ratechdin)
357  call mem_deallocate(this%ratechdout)
358  call mem_deallocate(this%head, 'HEAD', this%memoryPath)

◆ chd_df_obs()

subroutine chdmodule::chd_df_obs ( class(chdtype this)
private

(1) Store observation type supported by CHD package and (2) override BndTypebnd_df_obs

Definition at line 404 of file gwf-chd.f90.

405  implicit none
406  ! -- dummy
407  class(ChdType) :: this
408  ! -- local
409  integer(I4B) :: indx
410  !
411  call this%obs%StoreObsType('chd', .true., indx)
412  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
Here is the call graph for this function:

◆ chd_fc()

subroutine chdmodule::chd_fc ( class(chdtype 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

For constant head boundary type, the call to bnd_fc needs to be overwritten to do nothing

Definition at line 234 of file gwf-chd.f90.

235  ! -- dummy
236  class(ChdType) :: this
237  real(DP), dimension(:), intent(inout) :: rhs
238  integer(I4B), dimension(:), intent(in) :: ia
239  integer(I4B), dimension(:), intent(in) :: idxglo
240  class(MatrixBaseType), pointer :: matrix_sln
241  ! -- local

◆ chd_obs_supported()

logical function chdmodule::chd_obs_supported ( class(chdtype this)
private

Return true since CHD package supports observations

Definition at line 391 of file gwf-chd.f90.

392  implicit none
393  !
394  class(ChdType) :: this
395  !
396  chd_obs_supported = .true.

◆ chd_rp()

subroutine chdmodule::chd_rp ( class(chdtype), intent(inout)  this)

Definition at line 125 of file gwf-chd.f90.

126  !
127  use tdismodule, only: kper
128  ! -- dummy
129  class(ChdType), intent(inout) :: this
130  ! -- local
131  character(len=30) :: nodestr
132  integer(I4B) :: i, node, ibd, ierr
133  !
134  if (this%iper /= kper) return
135  !
136  ! -- Reset previous CHDs to active cell
137  do i = 1, this%nbound
138  node = this%nodelist(i)
139  this%ibound(node) = this%ibcnum
140  end do
141  !
142  ! -- Call the parent class read and prepare
143  call this%BndExtType%bnd_rp()
144  !
145  ! -- Set ibound to -(ibcnum + 1) for constant head cells
146  ierr = 0
147  do i = 1, this%nbound
148  node = this%nodelist(i)
149  ibd = this%ibound(node)
150  if (ibd < 0) then
151  call this%dis%noder_to_string(node, nodestr)
152  write (errmsg, '(3a)') &
153  'Cell is already a constant head (', trim(adjustl(nodestr)), ').'
154  call store_error(errmsg)
155  ierr = ierr + 1
156  else
157  this%ibound(node) = -this%ibcnum
158  end if
159  end do
160  !
161  ! -- Stop if errors detected
162  if (ierr > 0) then
163  call store_error_filename(this%input_fname)
164  end if
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
Here is the call graph for this function:

◆ define_listlabel()

subroutine chdmodule::define_listlabel ( class(chdtype), intent(inout)  this)

Definition at line 364 of file gwf-chd.f90.

365  class(ChdType), intent(inout) :: this
366  !
367  ! -- create the header list label
368  this%listlabel = trim(this%filtyp)//' NO.'
369  if (this%dis%ndim == 3) then
370  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
371  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
372  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
373  elseif (this%dis%ndim == 2) then
374  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
375  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
376  else
377  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
378  end if
379  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'HEAD'
380  if (this%inamedbound == 1) then
381  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
382  end if

◆ head_mult()

real(dp) function chdmodule::head_mult ( class(chdtype), intent(inout)  this,
integer(i4b), intent(in)  row 
)
private

Definition at line 417 of file gwf-chd.f90.

418  ! -- modules
419  use constantsmodule, only: dzero
420  ! -- dummy variables
421  class(ChdType), intent(inout) :: this
422  integer(I4B), intent(in) :: row
423  ! -- result
424  real(DP) :: head
425  !
426  if (this%iauxmultcol > 0) then
427  head = this%head(row) * this%auxvar(this%iauxmultcol, row)
428  else
429  head = this%head(row)
430  end if

Variable Documentation

◆ ftype

character(len=lenftype) chdmodule::ftype = 'CHD'
private

Definition at line 22 of file gwf-chd.f90.

22  character(len=LENFTYPE) :: ftype = 'CHD'

◆ text

character(len=lenpackagename) chdmodule::text = ' CHD'
private

Definition at line 23 of file gwf-chd.f90.

23  character(len=LENPACKAGENAME) :: text = ' CHD'