MODFLOW 6  version 6.7.0.dev2
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
Parameters
[in,out]thisBndExtType object

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

381  ! -- modules
382  use constantsmodule, only: dzero
383  ! -- dummy variables
384  class(GhbType), intent(inout) :: this !< BndExtType object
385  integer(I4B), intent(in) :: row
386  ! -- result
387  real(DP) :: cond
388  !
389  if (this%iauxmultcol > 0) then
390  cond = this%cond(row) * this%auxvar(this%iauxmultcol, row)
391  else
392  cond = this%cond(row)
393  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 307 of file gwf-ghb.f90.

308  ! -- dummy
309  class(GhbType), intent(inout) :: this
310  !
311  ! -- create the header list label
312  this%listlabel = trim(this%filtyp)//' NO.'
313  if (this%dis%ndim == 3) then
314  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
315  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
316  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
317  elseif (this%dis%ndim == 2) then
318  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
319  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
320  else
321  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
322  end if
323  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'STAGE'
324  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'CONDUCTANCE'
325  if (this%inamedbound == 1) then
326  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
327  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 141 of file gwf-ghb.f90.

142  ! -- modules
144  ! -- dummy
145  class(GhbType) :: this
146  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
147  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
148  !
149  ! -- call base type allocate arrays
150  call this%BndExtType%allocate_arrays(nodelist, auxvar)
151  !
152  ! -- set ghb input context pointers
153  call mem_setptr(this%bhead, 'BHEAD', this%input_mempath)
154  call mem_setptr(this%cond, 'COND', this%input_mempath)
155  !
156  ! --checkin ghb input context pointers
157  call mem_checkin(this%bhead, 'BHEAD', this%memoryPath, &
158  'BHEAD', this%input_mempath)
159  call mem_checkin(this%cond, 'COND', this%memoryPath, &
160  '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 
)
Parameters
[in,out]thisBndExtType object

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

399  ! -- modules
400  use constantsmodule, only: dzero
401  ! -- dummy
402  class(GhbType), intent(inout) :: this !< BndExtType object
403  integer(I4B), intent(in) :: col
404  integer(I4B), intent(in) :: row
405  ! -- result
406  real(DP) :: bndval
407  !
408  select case (col)
409  case (1)
410  bndval = this%bhead(row)
411  case (2)
412  bndval = this%cond_mult(row)
413  case default
414  errmsg = 'Programming error. GHB bound value requested column '&
415  &'outside range of ncolbnd (2).'
416  call store_error(errmsg)
417  call store_error_filename(this%input_fname)
418  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 244 of file gwf-ghb.f90.

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

◆ ghb_ck()

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

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

190  ! -- modules
191  use constantsmodule, only: linelength
193  ! -- dummy
194  class(GhbType), intent(inout) :: this
195  ! -- local
196  character(len=LINELENGTH) :: errmsg
197  integer(I4B) :: i
198  integer(I4B) :: node
199  real(DP) :: bt
200  ! -- formats
201  character(len=*), parameter :: fmtghberr = &
202  "('GHB BOUNDARY (',i0,') HEAD (',f10.3,') IS LESS THAN CELL &
203  &BOTTOM (',f10.3,')')"
204  character(len=*), parameter :: fmtcondmulterr = &
205  "('GHB BOUNDARY (',i0,') CONDUCTANCE MULTIPLIER (',g10.3,') IS &
206  &LESS THAN ZERO')"
207  character(len=*), parameter :: fmtconderr = &
208  "('GHB BOUNDARY (',i0,') CONDUCTANCE (',g10.3,') IS LESS THAN &
209  &ZERO')"
210  !
211  ! -- check stress period data
212  do i = 1, this%nbound
213  node = this%nodelist(i)
214  if (node == 0) cycle
215  bt = this%dis%bot(node)
216  ! -- accumulate errors
217  if (this%bhead(i) < bt .and. this%icelltype(node) /= 0) then
218  write (errmsg, fmt=fmtghberr) i, this%bhead(i), bt
219  call store_error(errmsg)
220  end if
221  if (this%iauxmultcol > 0) then
222  if (this%auxvar(this%iauxmultcol, i) < dzero) then
223  write (errmsg, fmt=fmtcondmulterr) &
224  i, this%auxvar(this%iauxmultcol, i)
225  call store_error(errmsg)
226  end if
227  end if
228  if (this%cond(i) < dzero) then
229  write (errmsg, fmt=fmtconderr) i, this%cond(i)
230  call store_error(errmsg)
231  end if
232  end do
233  !
234  !write summary of ghb package error messages
235  if (count_errors() > 0) then
236  call store_error_unit(this%inunit)
237  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%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 83 of file gwf-ghb.f90.

84  ! -- modules
86  ! -- dummy
87  class(GhbType) :: this
88  !
89  ! -- Deallocate parent package
90  call this%BndExtType%bnd_da()
91  !
92  ! -- arrays
93  call mem_deallocate(this%bhead, 'BHEAD', this%memoryPath)
94  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 348 of file gwf-ghb.f90.

349  implicit none
350  ! -- dummy
351  class(GhbType) :: this
352  ! -- local
353  integer(I4B) :: indx
354  !
355  call this%obs%StoreObsType('ghb', .true., indx)
356  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
357  !
358  ! -- Store obs type and assign procedure pointer
359  ! for to-mvr observation type.
360  call this%obs%StoreObsType('to-mvr', .true., indx)
361  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 269 of file gwf-ghb.f90.

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

◆ ghb_obs_supported()

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

Overrides BndTypebnd_obs_supported()

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

337  implicit none
338  ! -- dummy
339  class(GhbType) :: this
340  !
341  ghb_obs_supported = .true.

◆ ghb_options()

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

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

100  ! -- modules
103  ! -- dummy
104  class(GhbType), intent(inout) :: this
105  ! -- local
106  logical(LGP) :: found_mover
107  !
108  ! -- source base class options
109  call this%BndExtType%source_options()
110  !
111  ! -- source options from input context
112  call mem_set_value(this%imover, 'MOVER', this%input_mempath, found_mover)
113  !
114  ! -- log ghb specific options
115  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 165 of file gwf-ghb.f90.

166  ! -- modules
167  use tdismodule, only: kper
168  ! -- dummy
169  class(GhbType), intent(inout) :: this
170  !
171  if (this%iper /= kper) return
172  !
173  ! -- Call the parent class read and prepare
174  call this%BndExtType%bnd_rp()
175  !
176  ! -- store user cond
177  if (this%ivsc == 1) then
178  call this%ghb_store_user_cond()
179  end if
180  !
181  ! -- Write the list to iout if requested
182  if (this%iprpak /= 0) then
183  call this%write_list()
184  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
Parameters
[in,out]thisBndExtType object

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

367  ! -- dummy
368  class(GhbType), intent(inout) :: this !< BndExtType object
369  ! -- local
370  integer(I4B) :: n
371  !
372  ! -- store backup copy of conductance values
373  do n = 1, this%nbound
374  this%condinput(n) = this%cond_mult(n)
375  end do

◆ log_ghb_options()

subroutine ghbmodule::log_ghb_options ( class(ghbtype), intent(inout)  this,
logical(lgp), intent(in)  found_mover 
)
Parameters
[in,out]thisBndExtType object

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

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