MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
gwf-ghb.f90
Go to the documentation of this file.
1 module ghbmodule
2  use kindmodule, only: dp, i4b, lgp
4  use simvariablesmodule, only: errmsg
7  use bndmodule, only: bndtype
8  use bndextmodule, only: bndexttype
11  !
12  implicit none
13  !
14  private
15  public :: ghb_create
16  public :: ghbtype
17  !
18  character(len=LENFTYPE) :: ftype = 'GHB'
19  character(len=LENPACKAGENAME) :: text = ' GHB'
20  !
21  type, extends(bndexttype) :: ghbtype
22  real(dp), dimension(:), pointer, contiguous :: bhead => null() !< GHB boundary head
23  real(dp), dimension(:), pointer, contiguous :: cond => null() !< GHB hydraulic conductance
24  contains
25  procedure :: allocate_arrays => ghb_allocate_arrays
26  procedure :: source_options => ghb_options
27  procedure :: log_ghb_options
28  procedure :: bnd_rp => ghb_rp
29  procedure :: bnd_ck => ghb_ck
30  procedure :: bnd_cf => ghb_cf
31  procedure :: bnd_fc => ghb_fc
32  procedure :: bnd_da => ghb_da
33  procedure :: define_listlabel
34  procedure :: bound_value => ghb_bound_value
35  procedure :: cond_mult
36  ! -- methods for observations
37  procedure, public :: bnd_obs_supported => ghb_obs_supported
38  procedure, public :: bnd_df_obs => ghb_df_obs
39  procedure, public :: ghb_store_user_cond
40  end type ghbtype
41 
42 contains
43 
44  !> @brief Create a New Ghb Package and point bndobj to the new package
45  !<
46  subroutine ghb_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
47  mempath)
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')
80  end subroutine ghb_create
81 
82  !> @brief Deallocate memory
83  !<
84  subroutine ghb_da(this)
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)
96  end subroutine ghb_da
97 
98  !> @brief Set options specific to GhbType
99  !<
100  subroutine ghb_options(this)
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)
117  end subroutine ghb_options
118 
119  !> @brief Log options specific to GhbType
120  !<
121  subroutine log_ghb_options(this, found_mover)
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'
138  end subroutine log_ghb_options
139 
140  !> @brief Allocate arrays
141  !<
142  subroutine ghb_allocate_arrays(this, nodelist, auxvar)
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)
162  end subroutine ghb_allocate_arrays
163 
164  !> @brief Read and prepare
165  !<
166  subroutine ghb_rp(this)
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
181  end subroutine ghb_rp
182 
183  !> @brief Check ghb boundary condition data
184  !<
185  subroutine ghb_ck(this)
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
234  end subroutine ghb_ck
235 
236  !> @brief Formulate the HCOF and RHS terms
237  !!
238  !! Skip if no GHBs
239  !<
240  subroutine ghb_cf(this)
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
261  end subroutine ghb_cf
262 
263  !> @brief Copy rhs and hcof into solution rhs and amat
264  !<
265  subroutine ghb_fc(this, rhs, ia, idxglo, matrix_sln)
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
298  end subroutine ghb_fc
299 
300  !> @brief Define the list heading that is written to iout when PRINT_INPUT
301  !! option is used
302  !<
303  subroutine define_listlabel(this)
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
324  end subroutine define_listlabel
325 
326  ! -- Procedures related to observations
327 
328  !> @brief Return true because GHB package supports observations
329  !!
330  !! Overrides BndType%bnd_obs_supported()
331  !<
332  logical function ghb_obs_supported(this)
333  implicit none
334  ! -- dummy
335  class(ghbtype) :: this
336  !
337  ghb_obs_supported = .true.
338  end function ghb_obs_supported
339 
340  !> @brief Store observation type supported by GHB package
341  !!
342  !! Overrides BndType%bnd_df_obs
343  !<
344  subroutine ghb_df_obs(this)
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
358  end subroutine ghb_df_obs
359 
360  !> @brief Store user-specified conductance for GHB boundary type
361  !<
362  subroutine ghb_store_user_cond(this)
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
372  end subroutine ghb_store_user_cond
373 
374  !> @brief Apply multiplier to GHB conductance if option AUXMULTCOL is used
375  !<
376  function cond_mult(this, row) result(cond)
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
390  end function cond_mult
391 
392  !> @brief Return requested boundary value
393  !<
394  function ghb_bound_value(this, col, row) result(bndval)
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
415  end function ghb_bound_value
416 
417 end module ghbmodule
This module contains the extended boundary package.
This module contains the base boundary package.
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:23
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
Definition: Constants.f90:39
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
subroutine ghb_cf(this)
Formulate the HCOF and RHS terms.
Definition: gwf-ghb.f90:241
subroutine ghb_allocate_arrays(this, nodelist, auxvar)
Allocate arrays.
Definition: gwf-ghb.f90:143
subroutine ghb_ck(this)
Check ghb boundary condition data.
Definition: gwf-ghb.f90:186
subroutine ghb_da(this)
Deallocate memory.
Definition: gwf-ghb.f90:85
subroutine ghb_fc(this, rhs, ia, idxglo, matrix_sln)
Copy rhs and hcof into solution rhs and amat.
Definition: gwf-ghb.f90:266
subroutine ghb_df_obs(this)
Store observation type supported by GHB package.
Definition: gwf-ghb.f90:345
character(len=lenpackagename) text
Definition: gwf-ghb.f90:19
real(dp) function ghb_bound_value(this, col, row)
Return requested boundary value.
Definition: gwf-ghb.f90:395
real(dp) function cond_mult(this, row)
Apply multiplier to GHB conductance if option AUXMULTCOL is used.
Definition: gwf-ghb.f90:377
character(len=lenftype) ftype
Definition: gwf-ghb.f90:18
subroutine ghb_store_user_cond(this)
Store user-specified conductance for GHB boundary type.
Definition: gwf-ghb.f90:363
logical function ghb_obs_supported(this)
Return true because GHB package supports observations.
Definition: gwf-ghb.f90:333
subroutine ghb_rp(this)
Read and prepare.
Definition: gwf-ghb.f90:167
subroutine, public ghb_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
Create a New Ghb Package and point bndobj to the new package.
Definition: gwf-ghb.f90:48
subroutine log_ghb_options(this, found_mover)
Log options specific to GhbType.
Definition: gwf-ghb.f90:122
subroutine ghb_options(this)
Set options specific to GhbType.
Definition: gwf-ghb.f90:101
subroutine define_listlabel(this)
Define the list heading that is written to iout when PRINT_INPUT option is used.
Definition: gwf-ghb.f90:304
This module defines variable data types.
Definition: kind.f90:8
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains the derived type ObsType.
Definition: Obs.f90:127
subroutine, public defaultobsidprocessor(obsrv, dis, inunitobs, iout)
@ brief Process IDstring provided for each observation
Definition: Obs.f90:246
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_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
subroutine, public store_error_unit(iunit, terminate)
Store the file unit number.
Definition: Sim.f90:168
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
@ brief BndType
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23