MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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 282 of file gwf-chd.f90.

283  ! -- modules
284  ! -- dummy
285  class(ChdType), intent(inout) :: this
286  ! -- local
287  integer(I4B) :: i
288  integer(I4B) :: ipos
289  integer(I4B) :: node
290  integer(I4B) :: n2
291  integer(I4B) :: idiag
292  real(DP) :: rate
293  real(DP) :: ratein, rateout
294  real(DP) :: q
295  !
296  ! -- If no boundaries, skip flow calculations.
297  if (this%nbound > 0) then
298  !
299  ! -- Loop through each boundary calculating flow.
300  do i = 1, this%nbound
301  node = this%nodelist(i)
302  idiag = this%dis%con%ia(node)
303  rate = dzero
304  ratein = dzero
305  rateout = dzero
306  !
307  ! -- Calculate the flow rate into the cell.
308  do ipos = this%dis%con%ia(node) + 1, &
309  this%dis%con%ia(node + 1) - 1
310  q = this%flowja(ipos)
311  rate = rate - q
312  ! -- only accumulate chin and chout for active
313  ! connected cells
314  n2 = this%dis%con%ja(ipos)
315  if (this%ibound(n2) > 0) then
316  if (q < dzero) then
317  ratein = ratein - q
318  else
319  rateout = rateout + q
320  end if
321  end if
322  end do
323  !
324  ! -- For chd, store total flow in rhs so it is available for other
325  ! calculations
326  this%rhs(i) = -rate
327  this%hcof(i) = dzero
328  !
329  ! -- Save simulated value to simvals array.
330  this%simvals(i) = rate
331  this%ratechdin(i) = ratein
332  this%ratechdout(i) = rateout
333  this%flowja(idiag) = this%flowja(idiag) + rate
334  !
335  end do
336  !
337  end if
338  !
339  ! -- Return
340  return

◆ chd_ad()

subroutine chdmodule::chd_ad ( class(chdtype this)

Add package connections to matrix

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

185  ! -- modules
186  ! -- dummy
187  class(ChdType) :: this
188  ! -- local
189  integer(I4B) :: i, node
190  real(DP) :: hb
191  ! -- formats
192  !
193  ! -- Process each entry in the specified-head cell list
194  do i = 1, this%nbound
195  node = this%nodelist(i)
196  hb = this%head_mult(i)
197  !
198  this%xnew(node) = hb
199  this%xold(node) = this%xnew(node)
200  end do
201  !
202  ! -- For each observation, push simulated value and corresponding
203  ! simulation time from "current" to "preceding" and reset
204  ! "current" value.
205  call this%obs%obs_ad()
206  !
207  ! -- Return
208  return

◆ 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 95 of file gwf-chd.f90.

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

◆ chd_bd()

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

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

346  ! -- add package ratin/ratout to model budget
347  use tdismodule, only: delt
349  class(ChdType) :: this
350  type(BudgetType), intent(inout) :: model_budget
351  real(DP) :: ratin
352  real(DP) :: ratout
353  real(DP) :: dum
354  integer(I4B) :: isuppress_output
355 
356  ! For CHDs at an exchange, under some conditions
357  ! (XT3D), the model flowja into the cell is not
358  ! finalized until after exg_cq. So we calculate
359  ! the CHD rate here
360  call this%calc_chd_rate()
361 
362  isuppress_output = 0
363  call rate_accumulator(this%ratechdin(1:this%nbound), ratin, dum)
364  call rate_accumulator(this%ratechdout(1:this%nbound), ratout, dum)
365  call model_budget%addentry(ratin, ratout, delt, this%text, &
366  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:664
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.

Parameters
[in,out]thisBndType object

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

479  ! -- modules
480  use constantsmodule, only: dzero
481  ! -- dummy variables
482  class(ChdType), intent(inout) :: this !< BndType object
483  integer(I4B), intent(in) :: col
484  integer(I4B), intent(in) :: row
485  ! -- result
486  real(DP) :: bndval
487  !
488  select case (col)
489  case (1)
490  bndval = this%head_mult(row)
491  case default
492  errmsg = 'Programming error. CHD bound value requested column '&
493  &'outside range of ncolbnd (1).'
494  call store_error(errmsg)
495  call store_error_filename(this%input_fname)
496  end select
497  !
498  ! -- Return
499  return
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
Here is the call graph for this function:

◆ chd_ck()

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

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

214  ! -- modules
215  ! -- dummy
216  class(ChdType), intent(inout) :: this
217  ! -- local
218  character(len=30) :: nodestr
219  integer(I4B) :: i
220  integer(I4B) :: node
221  real(DP) :: bt
222  ! -- formats
223  character(len=*), parameter :: fmtchderr = &
224  "('CHD BOUNDARY ',i0,' HEAD (',g0,') IS LESS THAN CELL &
225  &BOTTOM (',g0,')',' FOR CELL ',a)"
226  !
227  ! -- check stress period data
228  do i = 1, this%nbound
229  node = this%nodelist(i)
230  bt = this%dis%bot(node)
231  ! -- accumulate errors
232  if (this%head_mult(i) < bt .and. this%icelltype(node) /= 0) then
233  call this%dis%noder_to_string(node, nodestr)
234  write (errmsg, fmt=fmtchderr) i, this%head_mult(i), bt, trim(nodestr)
235  call store_error(errmsg)
236  end if
237  end do
238  !
239  ! write summary of chd package error messages
240  if (count_errors() > 0) then
241  call store_error_filename(this%input_fname)
242  end if
243  !
244  ! -- Return
245  return
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 270 of file gwf-chd.f90.

271  class(ChdType), intent(inout) :: this
272  real(DP), dimension(:), intent(in) :: x
273  real(DP), dimension(:), contiguous, intent(inout) :: flowja
274  integer(I4B), optional, intent(in) :: iadv
275 
276  ! NB: the rate calculation cannot be done until chd_bd below
277 

◆ 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%ictMemPath = create_mem_path(namemodel, 'NPF')
88  !
89  ! -- Return
90  return
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 371 of file gwf-chd.f90.

372  ! -- modules
374  ! -- dummy
375  class(ChdType) :: this
376  !
377  ! -- Deallocate parent package
378  call this%BndExtType%bnd_da()
379  !
380  ! -- arrays
381  call mem_deallocate(this%ratechdin)
382  call mem_deallocate(this%ratechdout)
383  call mem_deallocate(this%head, 'HEAD', this%memoryPath)
384  !
385  ! -- Return
386  return

◆ 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 438 of file gwf-chd.f90.

439  implicit none
440  ! -- dummy
441  class(ChdType) :: this
442  ! -- local
443  integer(I4B) :: indx
444  !
445  call this%obs%StoreObsType('chd', .true., indx)
446  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
447  !
448  ! -- Return
449  return
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 253 of file gwf-chd.f90.

254  ! -- dummy
255  class(ChdType) :: this
256  real(DP), dimension(:), intent(inout) :: rhs
257  integer(I4B), dimension(:), intent(in) :: ia
258  integer(I4B), dimension(:), intent(in) :: idxglo
259  class(MatrixBaseType), pointer :: matrix_sln
260  ! -- local
261  !
262  ! -- Return
263  return

◆ chd_obs_supported()

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

Return true since CHD package supports observations

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

423  implicit none
424  !
425  class(ChdType) :: this
426  !
427  chd_obs_supported = .true.
428  !
429  ! -- Return
430  return

◆ chd_rp()

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

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

131  !
132  use tdismodule, only: kper
133  ! -- dummy
134  class(ChdType), intent(inout) :: this
135  ! -- local
136  character(len=30) :: nodestr
137  integer(I4B) :: i, node, ibd, ierr
138  !
139  if (this%iper /= kper) return
140  !
141  ! -- Reset previous CHDs to active cell
142  do i = 1, this%nbound
143  node = this%nodelist(i)
144  this%ibound(node) = this%ibcnum
145  end do
146  !
147  ! -- Call the parent class read and prepare
148  call this%BndExtType%bnd_rp()
149  !
150  ! -- Set ibound to -(ibcnum + 1) for constant head cells
151  ierr = 0
152  do i = 1, this%nbound
153  node = this%nodelist(i)
154  ibd = this%ibound(node)
155  if (ibd < 0) then
156  call this%dis%noder_to_string(node, nodestr)
157  write (errmsg, '(3a)') &
158  'Cell is already a constant head (', trim(adjustl(nodestr)), ').'
159  call store_error(errmsg)
160  ierr = ierr + 1
161  else
162  this%ibound(node) = -this%ibcnum
163  end if
164  end do
165  !
166  ! -- Stop if errors detected
167  if (ierr > 0) then
168  call store_error_filename(this%input_fname)
169  end if
170  !
171  ! -- Write the list to iout if requested
172  if (this%iprpak /= 0) then
173  call this%write_list()
174  end if
175  !
176  ! -- Return
177  return
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 392 of file gwf-chd.f90.

393  class(ChdType), intent(inout) :: this
394  !
395  ! -- create the header list label
396  this%listlabel = trim(this%filtyp)//' NO.'
397  if (this%dis%ndim == 3) then
398  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
399  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
400  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
401  elseif (this%dis%ndim == 2) then
402  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
403  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
404  else
405  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
406  end if
407  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'HEAD'
408  if (this%inamedbound == 1) then
409  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
410  end if
411  !
412  ! -- Return
413  return

◆ head_mult()

real(dp) function chdmodule::head_mult ( class(chdtype), intent(inout)  this,
integer(i4b), intent(in)  row 
)
private
Parameters
[in,out]thisBndExtType object

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

455  ! -- modules
456  use constantsmodule, only: dzero
457  ! -- dummy variables
458  class(ChdType), intent(inout) :: this !< BndExtType object
459  integer(I4B), intent(in) :: row
460  ! -- result
461  real(DP) :: head
462  !
463  if (this%iauxmultcol > 0) then
464  head = this%head(row) * this%auxvar(this%iauxmultcol, row)
465  else
466  head = this%head(row)
467  end if
468  !
469  ! -- Return
470  return

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'