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

Data Types

type  ghbtype
 

Functions/Subroutines

subroutine, public ghb_create (packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
 Create a New Ghb Package and point bndobj to the new package. More...
 
subroutine ghb_da (this)
 Deallocate memory. More...
 
subroutine ghb_options (this)
 Set options specific to GhbType. More...
 
subroutine log_ghb_options (this, found_mover)
 Log options specific to GhbType. More...
 
subroutine ghb_allocate_arrays (this, nodelist, auxvar)
 Allocate arrays. More...
 
subroutine ghb_rp (this)
 Read and prepare. More...
 
subroutine ghb_ck (this)
 Check ghb boundary condition data. More...
 
subroutine ghb_cf (this)
 Formulate the HCOF and RHS terms. More...
 
subroutine ghb_fc (this, rhs, ia, idxglo, matrix_sln)
 Copy rhs and hcof into solution rhs and amat. More...
 
subroutine define_listlabel (this)
 Define the list heading that is written to iout when PRINT_INPUT option is used. More...
 
logical function ghb_obs_supported (this)
 Return true because GHB package supports observations. More...
 
subroutine ghb_df_obs (this)
 Store observation type supported by GHB package. More...
 
subroutine ghb_store_user_cond (this)
 Store user-specified conductance for GHB boundary type. More...
 
real(dp) function cond_mult (this, row)
 Apply multiplier to GHB conductance if option AUXMULTCOL is used. More...
 
real(dp) function ghb_bound_value (this, col, row)
 Return requested boundary value. More...
 

Variables

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

Function/Subroutine Documentation

◆ cond_mult()

real(dp) function ghbmodule::cond_mult ( class(ghbtype), intent(inout)  this,
integer(i4b), intent(in)  row 
)
private

Definition at line 376 of file gwf-ghb.f90.

377  ! -- modules
378  use constantsmodule, only: dzero
379  ! -- dummy variables
380  class(GhbType), intent(inout) :: this
381  integer(I4B), intent(in) :: row
382  ! -- result
383  real(DP) :: cond
384  !
385  if (this%iauxmultcol > 0) then
386  cond = this%cond(row) * this%auxvar(this%iauxmultcol, row)
387  else
388  cond = this%cond(row)
389  end if
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65

◆ define_listlabel()

subroutine ghbmodule::define_listlabel ( class(ghbtype), intent(inout)  this)
private

Definition at line 303 of file gwf-ghb.f90.

304  ! -- dummy
305  class(GhbType), intent(inout) :: this
306  !
307  ! -- create the header list label
308  this%listlabel = trim(this%filtyp)//' NO.'
309  if (this%dis%ndim == 3) then
310  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
311  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
312  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
313  elseif (this%dis%ndim == 2) then
314  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
315  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
316  else
317  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
318  end if
319  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'STAGE'
320  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'CONDUCTANCE'
321  if (this%inamedbound == 1) then
322  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
323  end if

◆ ghb_allocate_arrays()

subroutine ghbmodule::ghb_allocate_arrays ( class(ghbtype this,
integer(i4b), dimension(:), optional, pointer, contiguous  nodelist,
real(dp), dimension(:, :), optional, pointer, contiguous  auxvar 
)
private

Definition at line 142 of file gwf-ghb.f90.

143  ! -- modules
145  ! -- dummy
146  class(GhbType) :: this
147  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
148  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
149  !
150  ! -- call base type allocate arrays
151  call this%BndExtType%allocate_arrays(nodelist, auxvar)
152  !
153  ! -- set ghb input context pointers
154  call mem_setptr(this%bhead, 'BHEAD', this%input_mempath)
155  call mem_setptr(this%cond, 'COND', this%input_mempath)
156  !
157  ! --checkin ghb input context pointers
158  call mem_checkin(this%bhead, 'BHEAD', this%memoryPath, &
159  'BHEAD', this%input_mempath)
160  call mem_checkin(this%cond, 'COND', this%memoryPath, &
161  'COND', this%input_mempath)

◆ ghb_bound_value()

real(dp) function ghbmodule::ghb_bound_value ( class(ghbtype), intent(inout)  this,
integer(i4b), intent(in)  col,
integer(i4b), intent(in)  row 
)

Definition at line 394 of file gwf-ghb.f90.

395  ! -- modules
396  use constantsmodule, only: dzero
397  ! -- dummy
398  class(GhbType), intent(inout) :: this
399  integer(I4B), intent(in) :: col
400  integer(I4B), intent(in) :: row
401  ! -- result
402  real(DP) :: bndval
403  !
404  select case (col)
405  case (1)
406  bndval = this%bhead(row)
407  case (2)
408  bndval = this%cond_mult(row)
409  case default
410  errmsg = 'Programming error. GHB bound value requested column '&
411  &'outside range of ncolbnd (2).'
412  call store_error(errmsg)
413  call store_error_filename(this%input_fname)
414  end select
Here is the call graph for this function:

◆ ghb_cf()

subroutine ghbmodule::ghb_cf ( class(ghbtype this)

Skip if no GHBs

Definition at line 240 of file gwf-ghb.f90.

241  ! -- dummy
242  class(GhbType) :: this
243  ! -- local
244  integer(I4B) :: i, node
245  !
246  ! -- Return if no ghbs
247  if (this%nbound .eq. 0) return
248  !
249  ! -- Calculate hcof and rhs for each ghb entry
250  do i = 1, this%nbound
251  node = this%nodelist(i)
252  if (node == 0) cycle
253  if (this%ibound(node) .le. 0) then
254  this%hcof(i) = dzero
255  this%rhs(i) = dzero
256  cycle
257  end if
258  this%hcof(i) = -this%cond_mult(i)
259  this%rhs(i) = -this%cond_mult(i) * this%bhead(i)
260  end do

◆ ghb_ck()

subroutine ghbmodule::ghb_ck ( class(ghbtype), intent(inout)  this)

Definition at line 185 of file gwf-ghb.f90.

186  ! -- modules
187  use constantsmodule, only: linelength
189  ! -- dummy
190  class(GhbType), intent(inout) :: this
191  ! -- local
192  character(len=LINELENGTH) :: errmsg
193  integer(I4B) :: i
194  integer(I4B) :: node
195  real(DP) :: bt
196  ! -- formats
197  character(len=*), parameter :: fmtghberr = &
198  "('GHB BOUNDARY (',i0,') HEAD (',f10.3,') IS LESS THAN CELL &
199  &BOTTOM (',f10.3,')')"
200  character(len=*), parameter :: fmtcondmulterr = &
201  "('GHB BOUNDARY (',i0,') CONDUCTANCE MULTIPLIER (',g10.3,') IS &
202  &LESS THAN ZERO')"
203  character(len=*), parameter :: fmtconderr = &
204  "('GHB BOUNDARY (',i0,') CONDUCTANCE (',g10.3,') IS LESS THAN &
205  &ZERO')"
206  !
207  ! -- check stress period data
208  do i = 1, this%nbound
209  node = this%nodelist(i)
210  if (node == 0) cycle
211  bt = this%dis%bot(node)
212  ! -- accumulate errors
213  if (this%bhead(i) < bt .and. this%icelltype(node) /= 0) then
214  write (errmsg, fmt=fmtghberr) i, this%bhead(i), bt
215  call store_error(errmsg)
216  end if
217  if (this%iauxmultcol > 0) then
218  if (this%auxvar(this%iauxmultcol, i) < dzero) then
219  write (errmsg, fmt=fmtcondmulterr) &
220  i, this%auxvar(this%iauxmultcol, i)
221  call store_error(errmsg)
222  end if
223  end if
224  if (this%cond(i) < dzero) then
225  write (errmsg, fmt=fmtconderr) i, this%cond(i)
226  call store_error(errmsg)
227  end if
228  end do
229  !
230  !write summary of ghb package error messages
231  if (count_errors() > 0) then
232  call store_error_unit(this%inunit)
233  end if
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
subroutine, public store_error_unit(iunit, terminate)
Store the file unit number.
Definition: Sim.f90:168
Here is the call graph for this function:

◆ ghb_create()

subroutine, public ghbmodule::ghb_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 
)

Definition at line 46 of file gwf-ghb.f90.

48  ! -- dummy
49  class(BndType), pointer :: packobj
50  integer(I4B), intent(in) :: id
51  integer(I4B), intent(in) :: ibcnum
52  integer(I4B), intent(in) :: inunit
53  integer(I4B), intent(in) :: iout
54  character(len=*), intent(in) :: namemodel
55  character(len=*), intent(in) :: pakname
56  character(len=*), intent(in) :: mempath
57  ! -- local
58  type(GhbType), pointer :: ghbobj
59  !
60  ! -- allocate the object and assign values to object variables
61  allocate (ghbobj)
62  packobj => ghbobj
63  !
64  ! -- create name and memory path
65  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
66  packobj%text = text
67  !
68  ! -- allocate scalars
69  call packobj%allocate_scalars()
70  !
71  ! -- initialize package
72  call packobj%pack_initialize()
73  !
74  packobj%inunit = inunit
75  packobj%iout = iout
76  packobj%id = id
77  packobj%ibcnum = ibcnum
78  packobj%ncolbnd = 2
79  packobj%ictMemPath = create_mem_path(namemodel, 'NPF')
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ghb_da()

subroutine ghbmodule::ghb_da ( class(ghbtype this)
private

Definition at line 84 of file gwf-ghb.f90.

85  ! -- modules
87  ! -- dummy
88  class(GhbType) :: this
89  !
90  ! -- Deallocate parent package
91  call this%BndExtType%bnd_da()
92  !
93  ! -- arrays
94  call mem_deallocate(this%bhead, 'BHEAD', this%memoryPath)
95  call mem_deallocate(this%cond, 'COND', this%memoryPath)

◆ ghb_df_obs()

subroutine ghbmodule::ghb_df_obs ( class(ghbtype this)
private

Overrides BndTypebnd_df_obs

Definition at line 344 of file gwf-ghb.f90.

345  implicit none
346  ! -- dummy
347  class(GhbType) :: this
348  ! -- local
349  integer(I4B) :: indx
350  !
351  call this%obs%StoreObsType('ghb', .true., indx)
352  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
353  !
354  ! -- Store obs type and assign procedure pointer
355  ! for to-mvr observation type.
356  call this%obs%StoreObsType('to-mvr', .true., indx)
357  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
Here is the call graph for this function:

◆ ghb_fc()

subroutine ghbmodule::ghb_fc ( class(ghbtype 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

Definition at line 265 of file gwf-ghb.f90.

266  ! -- dummy
267  class(GhbType) :: this
268  real(DP), dimension(:), intent(inout) :: rhs
269  integer(I4B), dimension(:), intent(in) :: ia
270  integer(I4B), dimension(:), intent(in) :: idxglo
271  class(MatrixBaseType), pointer :: matrix_sln
272  ! -- local
273  integer(I4B) :: i, n, ipos
274  real(DP) :: cond, bhead, qghb
275  !
276  ! -- pakmvrobj fc
277  if (this%imover == 1) then
278  call this%pakmvrobj%fc()
279  end if
280  !
281  ! -- Copy package rhs and hcof into solution rhs and amat
282  do i = 1, this%nbound
283  n = this%nodelist(i)
284  if (n == 0) cycle
285  rhs(n) = rhs(n) + this%rhs(i)
286  ipos = ia(n)
287  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
288  !
289  ! -- If mover is active and this boundary is discharging,
290  ! store available water (as positive value).
291  bhead = this%bhead(i)
292  if (this%imover == 1 .and. this%xnew(n) > bhead) then
293  cond = this%cond_mult(i)
294  qghb = cond * (this%xnew(n) - bhead)
295  call this%pakmvrobj%accumulate_qformvr(i, qghb)
296  end if
297  end do

◆ ghb_obs_supported()

logical function ghbmodule::ghb_obs_supported ( class(ghbtype this)
private

Overrides BndTypebnd_obs_supported()

Definition at line 332 of file gwf-ghb.f90.

333  implicit none
334  ! -- dummy
335  class(GhbType) :: this
336  !
337  ghb_obs_supported = .true.

◆ ghb_options()

subroutine ghbmodule::ghb_options ( class(ghbtype), intent(inout)  this)

Definition at line 100 of file gwf-ghb.f90.

101  ! -- modules
104  ! -- dummy
105  class(GhbType), intent(inout) :: this
106  ! -- local
107  logical(LGP) :: found_mover
108  !
109  ! -- source base class options
110  call this%BndExtType%source_options()
111  !
112  ! -- source options from input context
113  call mem_set_value(this%imover, 'MOVER', this%input_mempath, found_mover)
114  !
115  ! -- log ghb specific options
116  call this%log_ghb_options(found_mover)
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23

◆ ghb_rp()

subroutine ghbmodule::ghb_rp ( class(ghbtype), intent(inout)  this)

Definition at line 166 of file gwf-ghb.f90.

167  ! -- modules
168  use tdismodule, only: kper
169  ! -- dummy
170  class(GhbType), intent(inout) :: this
171  !
172  if (this%iper /= kper) return
173  !
174  ! -- Call the parent class read and prepare
175  call this%BndExtType%bnd_rp()
176  !
177  ! -- store user cond
178  if (this%ivsc == 1) then
179  call this%ghb_store_user_cond()
180  end if
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23

◆ ghb_store_user_cond()

subroutine ghbmodule::ghb_store_user_cond ( class(ghbtype), intent(inout)  this)
private

Definition at line 362 of file gwf-ghb.f90.

363  ! -- dummy
364  class(GhbType), intent(inout) :: this
365  ! -- local
366  integer(I4B) :: n
367  !
368  ! -- store backup copy of conductance values
369  do n = 1, this%nbound
370  this%condinput(n) = this%cond_mult(n)
371  end do

◆ log_ghb_options()

subroutine ghbmodule::log_ghb_options ( class(ghbtype), intent(inout)  this,
logical(lgp), intent(in)  found_mover 
)

Definition at line 121 of file gwf-ghb.f90.

122  ! -- modules
123  ! -- dummy
124  class(GhbType), intent(inout) :: this
125  logical(LGP), intent(in) :: found_mover
126  !
127  ! -- log found options
128  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
129  //' OPTIONS'
130  !
131  if (found_mover) then
132  write (this%iout, '(4x,A)') 'MOVER OPTION ENABLED'
133  end if
134  !
135  ! -- close logging block
136  write (this%iout, '(1x,a)') &
137  'END OF '//trim(adjustl(this%text))//' OPTIONS'

Variable Documentation

◆ ftype

character(len=lenftype) ghbmodule::ftype = 'GHB'
private

Definition at line 18 of file gwf-ghb.f90.

18  character(len=LENFTYPE) :: ftype = 'GHB'

◆ text

character(len=lenpackagename) ghbmodule::text = ' GHB'
private

Definition at line 19 of file gwf-ghb.f90.

19  character(len=LENPACKAGENAME) :: text = ' GHB'