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

419  ! -- modules
420  use constantsmodule, only: dzero
421  ! -- dummy variables
422  class(GhbType), intent(inout) :: this !< BndExtType object
423  integer(I4B), intent(in) :: row
424  ! -- result
425  real(DP) :: cond
426  !
427  if (this%iauxmultcol > 0) then
428  cond = this%cond(row) * this%auxvar(this%iauxmultcol, row)
429  else
430  cond = this%cond(row)
431  end if
432  !
433  ! -- Return
434  return
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64

◆ define_listlabel()

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

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

334  ! -- dummy
335  class(GhbType), intent(inout) :: this
336  !
337  ! -- create the header list label
338  this%listlabel = trim(this%filtyp)//' NO.'
339  if (this%dis%ndim == 3) then
340  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
341  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
342  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
343  elseif (this%dis%ndim == 2) then
344  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
345  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
346  else
347  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
348  end if
349  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'STAGE'
350  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'CONDUCTANCE'
351  if (this%inamedbound == 1) then
352  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
353  end if
354  !
355  ! -- Return
356  return

◆ 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 
)

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

156  ! -- modules
158  ! -- dummy
159  class(GhbType) :: this
160  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
161  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
162  !
163  ! -- call base type allocate arrays
164  call this%BndExtType%allocate_arrays(nodelist, auxvar)
165  !
166  ! -- set ghb input context pointers
167  call mem_setptr(this%bhead, 'BHEAD', this%input_mempath)
168  call mem_setptr(this%cond, 'COND', this%input_mempath)
169  !
170  ! --checkin ghb input context pointers
171  call mem_checkin(this%bhead, 'BHEAD', this%memoryPath, &
172  'BHEAD', this%input_mempath)
173  call mem_checkin(this%cond, 'COND', this%memoryPath, &
174  'COND', this%input_mempath)
175  !
176  ! -- Return
177  return

◆ 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 439 of file gwf-ghb.f90.

440  ! -- modules
441  use constantsmodule, only: dzero
442  ! -- dummy
443  class(GhbType), intent(inout) :: this !< BndExtType object
444  integer(I4B), intent(in) :: col
445  integer(I4B), intent(in) :: row
446  ! -- result
447  real(DP) :: bndval
448  !
449  select case (col)
450  case (1)
451  bndval = this%bhead(row)
452  case (2)
453  bndval = this%cond_mult(row)
454  case default
455  errmsg = 'Programming error. GHB bound value requested column '&
456  &'outside range of ncolbnd (2).'
457  call store_error(errmsg)
458  call store_error_filename(this%input_fname)
459  end select
460  !
461  ! -- Return
462  return
Here is the call graph for this function:

◆ ghb_cf()

subroutine ghbmodule::ghb_cf ( class(ghbtype this)

Skip if no GHBs

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

267  ! -- dummy
268  class(GhbType) :: this
269  ! -- local
270  integer(I4B) :: i, node
271  !
272  ! -- Return if no ghbs
273  if (this%nbound .eq. 0) return
274  !
275  ! -- Calculate hcof and rhs for each ghb entry
276  do i = 1, this%nbound
277  node = this%nodelist(i)
278  if (this%ibound(node) .le. 0) then
279  this%hcof(i) = dzero
280  this%rhs(i) = dzero
281  cycle
282  end if
283  this%hcof(i) = -this%cond_mult(i)
284  this%rhs(i) = -this%cond_mult(i) * this%bhead(i)
285  end do
286  !
287  ! -- Return
288  return

◆ ghb_ck()

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

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

210  ! -- modules
211  use constantsmodule, only: linelength
213  ! -- dummy
214  class(GhbType), intent(inout) :: this
215  ! -- local
216  character(len=LINELENGTH) :: errmsg
217  integer(I4B) :: i
218  integer(I4B) :: node
219  real(DP) :: bt
220  ! -- formats
221  character(len=*), parameter :: fmtghberr = &
222  "('GHB BOUNDARY (',i0,') HEAD (',f10.3,') IS LESS THAN CELL &
223  &BOTTOM (',f10.3,')')"
224  character(len=*), parameter :: fmtcondmulterr = &
225  "('GHB BOUNDARY (',i0,') CONDUCTANCE MULTIPLIER (',g10.3,') IS &
226  &LESS THAN ZERO')"
227  character(len=*), parameter :: fmtconderr = &
228  "('GHB BOUNDARY (',i0,') CONDUCTANCE (',g10.3,') IS LESS THAN &
229  &ZERO')"
230  !
231  ! -- check stress period data
232  do i = 1, this%nbound
233  node = this%nodelist(i)
234  bt = this%dis%bot(node)
235  ! -- accumulate errors
236  if (this%bhead(i) < bt .and. this%icelltype(node) /= 0) then
237  write (errmsg, fmt=fmtghberr) i, this%bhead(i), bt
238  call store_error(errmsg)
239  end if
240  if (this%iauxmultcol > 0) then
241  if (this%auxvar(this%iauxmultcol, i) < dzero) then
242  write (errmsg, fmt=fmtcondmulterr) &
243  i, this%auxvar(this%iauxmultcol, i)
244  call store_error(errmsg)
245  end if
246  end if
247  if (this%cond(i) < dzero) then
248  write (errmsg, fmt=fmtconderr) i, this%cond(i)
249  call store_error(errmsg)
250  end if
251  end do
252  !
253  !write summary of ghb package error messages
254  if (count_errors() > 0) then
255  call store_error_unit(this%inunit)
256  end if
257  !
258  ! -- Return
259  return
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
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')
79  !
80  ! -- Return
81  return
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 86 of file gwf-ghb.f90.

87  ! -- modules
89  ! -- dummy
90  class(GhbType) :: this
91  !
92  ! -- Deallocate parent package
93  call this%BndExtType%bnd_da()
94  !
95  ! -- arrays
96  call mem_deallocate(this%bhead, 'BHEAD', this%memoryPath)
97  call mem_deallocate(this%cond, 'COND', this%memoryPath)
98  !
99  ! -- Return
100  return

◆ ghb_df_obs()

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

Overrides BndTypebnd_df_obs

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

381  implicit none
382  ! -- dummy
383  class(GhbType) :: this
384  ! -- local
385  integer(I4B) :: indx
386  !
387  call this%obs%StoreObsType('ghb', .true., indx)
388  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
389  !
390  ! -- Store obs type and assign procedure pointer
391  ! for to-mvr observation type.
392  call this%obs%StoreObsType('to-mvr', .true., indx)
393  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
394  !
395  ! -- Return
396  return
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 293 of file gwf-ghb.f90.

294  ! -- dummy
295  class(GhbType) :: this
296  real(DP), dimension(:), intent(inout) :: rhs
297  integer(I4B), dimension(:), intent(in) :: ia
298  integer(I4B), dimension(:), intent(in) :: idxglo
299  class(MatrixBaseType), pointer :: matrix_sln
300  ! -- local
301  integer(I4B) :: i, n, ipos
302  real(DP) :: cond, bhead, qghb
303  !
304  ! -- pakmvrobj fc
305  if (this%imover == 1) then
306  call this%pakmvrobj%fc()
307  end if
308  !
309  ! -- Copy package rhs and hcof into solution rhs and amat
310  do i = 1, this%nbound
311  n = this%nodelist(i)
312  rhs(n) = rhs(n) + this%rhs(i)
313  ipos = ia(n)
314  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
315  !
316  ! -- If mover is active and this boundary is discharging,
317  ! store available water (as positive value).
318  bhead = this%bhead(i)
319  if (this%imover == 1 .and. this%xnew(n) > bhead) then
320  cond = this%cond_mult(i)
321  qghb = cond * (this%xnew(n) - bhead)
322  call this%pakmvrobj%accumulate_qformvr(i, qghb)
323  end if
324  end do
325  !
326  ! -- Return
327  return

◆ ghb_obs_supported()

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

Overrides BndTypebnd_obs_supported()

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

366  implicit none
367  ! -- dummy
368  class(GhbType) :: this
369  !
370  ghb_obs_supported = .true.
371  !
372  ! -- Return
373  return

◆ ghb_options()

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

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

106  ! -- modules
110  ! -- dummy
111  class(GhbType), intent(inout) :: this
112  ! -- local
113  type(GwfGhbParamFoundType) :: found
114  !
115  ! -- source base class options
116  call this%BndExtType%source_options()
117  !
118  ! -- source options from input context
119  call mem_set_value(this%imover, 'MOVER', this%input_mempath, found%mover)
120  !
121  ! -- log ghb specific options
122  call this%log_ghb_options(found)
123  !
124  ! -- Return
125  return
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 182 of file gwf-ghb.f90.

183  ! -- modules
184  use tdismodule, only: kper
185  ! -- dummy
186  class(GhbType), intent(inout) :: this
187  !
188  if (this%iper /= kper) return
189  !
190  ! -- Call the parent class read and prepare
191  call this%BndExtType%bnd_rp()
192  !
193  ! -- store user cond
194  if (this%ivsc == 1) then
195  call this%ghb_store_user_cond()
196  end if
197  !
198  ! -- Write the list to iout if requested
199  if (this%iprpak /= 0) then
200  call this%write_list()
201  end if
202  !
203  ! -- Return
204  return
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 401 of file gwf-ghb.f90.

402  ! -- dummy
403  class(GhbType), intent(inout) :: this !< BndExtType object
404  ! -- local
405  integer(I4B) :: n
406  !
407  ! -- store backup copy of conductance values
408  do n = 1, this%nbound
409  this%condinput(n) = this%cond_mult(n)
410  end do
411  !
412  ! -- Return
413  return

◆ log_ghb_options()

subroutine ghbmodule::log_ghb_options ( class(ghbtype), intent(inout)  this,
type(gwfghbparamfoundtype), intent(in)  found 
)
Parameters
[in,out]thisBndExtType object

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

131  ! -- modules
133  ! -- dummy
134  class(GhbType), intent(inout) :: this !< BndExtType object
135  type(GwfGhbParamFoundType), intent(in) :: found
136  !
137  ! -- log found options
138  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
139  //' OPTIONS'
140  !
141  if (found%mover) then
142  write (this%iout, '(4x,A)') 'MOVER OPTION ENABLED'
143  end if
144  !
145  ! -- close logging block
146  write (this%iout, '(1x,a)') &
147  'END OF '//trim(adjustl(this%text))//' OPTIONS'
148  !
149  ! -- Return
150  return

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'