MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
gwf-rch.f90
Go to the documentation of this file.
1 module rchmodule
2  !
3  use kindmodule, only: dp, i4b, lgp
7  use bndmodule, only: bndtype
8  use bndextmodule, only: bndexttype
10  use simvariablesmodule, only: errmsg
13  !
14  implicit none
15  !
16  private
17  public :: rch_create
18  !
19  character(len=LENFTYPE) :: ftype = 'RCH'
20  character(len=LENPACKAGENAME) :: text = ' RCH'
21  character(len=LENPACKAGENAME) :: texta = ' RCHA'
22  !
23  type, extends(bndexttype) :: rchtype
24  real(dp), dimension(:), pointer, contiguous :: recharge => null() !< boundary recharge array
25  integer(I4B), dimension(:), pointer, contiguous :: nodesontop => null() ! User provided cell numbers; nodelist is cells where recharge is applied)
26  logical, pointer, private :: fixed_cell
27 
28  contains
29 
30  procedure :: rch_allocate_scalars
31  procedure :: allocate_arrays => rch_allocate_arrays
32  procedure :: source_options => rch_source_options
33  procedure :: log_rch_options
34  procedure :: read_initial_attr => rch_read_initial_attr
35  procedure :: bnd_rp => rch_rp
36  procedure :: bnd_cf => rch_cf
37  procedure :: bnd_fc => rch_fc
38  procedure :: bnd_da => rch_da
39  procedure :: set_nodesontop
40  procedure :: define_listlabel => rch_define_listlabel
41  procedure :: bound_value => rch_bound_value
42  ! -- for observations
43  procedure, public :: bnd_obs_supported => rch_obs_supported
44  procedure, public :: bnd_df_obs => rch_df_obs
45 
46  end type rchtype
47 
48 contains
49 
50  !> @brief Create a New Recharge Package
51  !!
52  !! Create new RCH package and point packobj to the new package
53  !<
54  subroutine rch_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
55  mempath)
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(rchtype), pointer :: rchobj
67  !
68  ! -- allocate recharge object and scalar variables
69  allocate (rchobj)
70  packobj => rchobj
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 rchobj%rch_allocate_scalars()
78  !
79  ! -- initialize package
80  call packobj%pack_initialize()
81  !
82  packobj%inunit = inunit
83  packobj%iout = iout
84  packobj%id = id
85  packobj%ibcnum = ibcnum
86  packobj%ictMemPath = create_mem_path(namemodel, 'NPF')
87  end subroutine rch_create
88 
89  !> @brief Allocate scalar members
90  !<
91  subroutine rch_allocate_scalars(this)
92  ! -- dummy
93  class(rchtype), intent(inout) :: this
94  !
95  ! -- allocate base scalars
96  call this%BndExtType%allocate_scalars()
97  !
98  ! -- allocate internal members
99  allocate (this%fixed_cell)
100  !
101  ! -- Set values
102  this%fixed_cell = .false.
103  end subroutine rch_allocate_scalars
104 
105  !> @brief Allocate package arrays
106  !<
107  subroutine rch_allocate_arrays(this, nodelist, auxvar)
108  ! -- modules
110  ! -- dummy
111  class(rchtype) :: this
112  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
113  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
114  !
115  ! -- allocate base arrays
116  call this%BndExtType%allocate_arrays(nodelist, auxvar)
117  !
118  ! -- set recharge input context pointer
119  call mem_setptr(this%recharge, 'RECHARGE', this%input_mempath)
120  !
121  ! -- checkin recharge input context pointer
122  call mem_checkin(this%recharge, 'RECHARGE', this%memoryPath, &
123  'RECHARGE', this%input_mempath)
124  end subroutine rch_allocate_arrays
125 
126  !> @brief Source options specific to RchType
127  !<
128  subroutine rch_source_options(this)
129  ! -- modules
131  implicit none
132  ! -- dummy
133  class(rchtype), intent(inout) :: this
134  ! -- local
135  logical(LGP) :: found_fixed_cell = .false.
136  !
137  ! -- source common bound options
138  call this%BndExtType%source_options()
139  !
140  ! -- update defaults with idm sourced values
141  call mem_set_value(this%fixed_cell, 'FIXED_CELL', this%input_mempath, &
142  found_fixed_cell)
143  !
144  if (this%readasarrays) then
145  this%text = texta
146  end if
147  !
148  ! -- log rch params
149  call this%log_rch_options(found_fixed_cell)
150  end subroutine rch_source_options
151 
152  !> @brief Log options specific to RchType
153  !<
154  subroutine log_rch_options(this, found_fixed_cell)
155  implicit none
156  ! -- dummy
157  class(rchtype), intent(inout) :: this
158  logical(LGP), intent(in) :: found_fixed_cell
159  ! -- formats
160  character(len=*), parameter :: fmtfixedcell = &
161  &"(4x, 'RECHARGE WILL BE APPLIED TO SPECIFIED CELL.')"
162  !
163  ! -- log found options
164  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
165  //' OPTIONS'
166  !
167  if (found_fixed_cell) then
168  write (this%iout, fmtfixedcell)
169  end if
170  !
171  ! -- close logging block
172  write (this%iout, '(1x,a)') &
173  'END OF '//trim(adjustl(this%text))//' OPTIONS'
174  end subroutine log_rch_options
175 
176  !> @brief Part of allocate and read
177  !<
178  subroutine rch_read_initial_attr(this)
179  ! -- dummy
180  class(rchtype), intent(inout) :: this
181  !
182  if (this%readasarrays) then
183  call this%default_nodelist()
184  end if
185  !
186  ! -- if fixed_cell option not set, then need to store nodelist
187  ! in the nodesontop array
188  if (.not. this%fixed_cell) call this%set_nodesontop()
189  end subroutine rch_read_initial_attr
190 
191  !> @brief Read and Prepare
192  !!
193  !! Read itmp and read new boundaries if itmp > 0
194  !<
195  subroutine rch_rp(this)
196  ! -- modules
197  use tdismodule, only: kper
198  implicit none
199  ! -- dummy
200  class(rchtype), intent(inout) :: this
201  !
202  if (this%iper /= kper) return
203  !
204  call this%BndExtType%bnd_rp()
205  !
206  ! -- copy nodelist to nodesontop if not fixed cell
207  if (.not. this%fixed_cell) call this%set_nodesontop()
208  !
209  if (this%iprpak /= 0) then
210  if (this%readasarrays) then
211  ! no-op
212  else
213  ! -- Write the list to iout
214  call this%write_list()
215  end if
216  end if
217  end subroutine rch_rp
218 
219  !> @brief Store nodelist in nodesontop
220  !<
221  subroutine set_nodesontop(this)
222  implicit none
223  ! -- dummy
224  class(rchtype), intent(inout) :: this
225  ! -- local
226  integer(I4B) :: n
227  !
228  ! -- allocate if necessary
229  if (.not. associated(this%nodesontop)) then
230  allocate (this%nodesontop(this%maxbound))
231  end if
232  !
233  ! -- copy nodelist into nodesontop
234  do n = 1, this%nbound
235  this%nodesontop(n) = this%nodelist(n)
236  end do
237  end subroutine set_nodesontop
238 
239  !> @brief Formulate the HCOF and RHS terms
240  !!
241  !! Skip if no recharge. Otherwise, calculate hcof and rhs
242  !<
243  subroutine rch_cf(this)
244  ! -- dummy
245  class(rchtype) :: this
246  ! -- local
247  integer(I4B) :: i, node
248  !
249  ! -- Return if no recharge
250  if (this%nbound == 0) return
251  !
252  ! -- Calculate hcof and rhs for each recharge entry
253  do i = 1, this%nbound
254  !
255  ! -- Find the node number
256  if (this%fixed_cell) then
257  node = this%nodelist(i)
258  else
259  node = this%nodesontop(i)
260  end if
261  !
262  ! -- cycle if nonexistent bound
263  if (node <= 0) then
264  this%hcof(i) = dzero
265  this%rhs(i) = dzero
266  cycle
267  end if
268  !
269  ! -- reset nodelist to highest active
270  if (.not. this%fixed_cell) then
271  if (this%ibound(node) == 0) &
272  call this%dis%highest_active(node, this%ibound)
273  this%nodelist(i) = node
274  end if
275  !
276  ! -- Set rhs and hcof
277  this%hcof(i) = dzero
278  if (this%iauxmultcol > 0) then
279  this%rhs(i) = -this%recharge(i) * this%dis%get_area(node) * &
280  this%auxvar(this%iauxmultcol, i)
281  else
282  this%rhs(i) = -this%recharge(i) * this%dis%get_area(node)
283  end if
284  if (this%ibound(node) <= 0) then
285  this%rhs(i) = dzero
286  cycle
287  end if
288  if (this%ibound(node) == iwetlake) then
289  this%rhs(i) = dzero
290  cycle
291  end if
292  end do
293  end subroutine rch_cf
294 
295  !> @brief Copy rhs and hcof into solution rhs and amat
296  !<
297  subroutine rch_fc(this, rhs, ia, idxglo, matrix_sln)
298  ! -- dummy
299  class(rchtype) :: this
300  real(DP), dimension(:), intent(inout) :: rhs
301  integer(I4B), dimension(:), intent(in) :: ia
302  integer(I4B), dimension(:), intent(in) :: idxglo
303  class(matrixbasetype), pointer :: matrix_sln
304  ! -- local
305  integer(I4B) :: i, n, ipos
306  !
307  ! -- Copy package rhs and hcof into solution rhs and amat
308  do i = 1, this%nbound
309  n = this%nodelist(i)
310  if (n <= 0) cycle
311  ! -- reset hcof and rhs for excluded cells
312  if (this%ibound(n) == iwetlake) then
313  this%hcof(i) = dzero
314  this%rhs(i) = dzero
315  cycle
316  end if
317  rhs(n) = rhs(n) + this%rhs(i)
318  ipos = ia(n)
319  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
320  end do
321  end subroutine rch_fc
322 
323  !> @brief Deallocate memory
324  !<
325  subroutine rch_da(this)
326  ! -- modules
328  ! -- dummy
329  class(rchtype) :: this
330  !
331  ! -- Deallocate parent package
332  call this%BndExtType%bnd_da()
333  !
334  ! -- scalars
335  deallocate (this%fixed_cell)
336  !
337  ! -- arrays
338  if (associated(this%nodesontop)) deallocate (this%nodesontop)
339  call mem_deallocate(this%recharge, 'RECHARGE', this%memoryPath)
340  end subroutine rch_da
341 
342  !> @brief Define the list heading that is written to iout when PRINT_INPUT
343  !! option is used.
344  !<
345  subroutine rch_define_listlabel(this)
346  ! -- dummy
347  class(rchtype), intent(inout) :: this
348  !
349  ! -- create the header list label
350  this%listlabel = trim(this%filtyp)//' NO.'
351  if (this%dis%ndim == 3) then
352  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
353  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
354  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
355  elseif (this%dis%ndim == 2) then
356  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
357  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
358  else
359  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
360  end if
361  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'RECHARGE'
362 ! if(this%multindex > 0) &
363 ! write(this%listlabel, '(a, a16)') trim(this%listlabel), 'MULTIPLIER'
364  if (this%inamedbound == 1) then
365  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
366  end if
367  end subroutine rch_define_listlabel
368 
369  ! -- Procedures related to observations
370 
371  !> @brief
372  !!
373  !! Overrides BndType%bnd_obs_supported()
374  !<
375  logical function rch_obs_supported(this)
376  implicit none
377  ! -- dummy
378  class(rchtype) :: this
379  !
380  rch_obs_supported = .true.
381  end function rch_obs_supported
382 
383  !> @brief Implements bnd_df_obs
384  !!
385  !! Store observation type supported by RCH package. Overrides
386  !! BndType%bnd_df_obs
387  !<
388  subroutine rch_df_obs(this)
389  implicit none
390  ! -- dummy
391  class(rchtype) :: this
392  ! -- local
393  integer(I4B) :: indx
394  !
395  call this%obs%StoreObsType('rch', .true., indx)
396  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
397  end subroutine rch_df_obs
398 
399  !> @brief Return requested boundary value
400  !<
401  function rch_bound_value(this, col, row) result(bndval)
402  ! -- modules
403  use constantsmodule, only: dzero
404  ! -- dummy
405  class(rchtype), intent(inout) :: this !< BndExtType object
406  integer(I4B), intent(in) :: col
407  integer(I4B), intent(in) :: row
408  ! -- result
409  real(dp) :: bndval
410  !
411  select case (col)
412  case (1)
413  if (this%iauxmultcol > 0) then
414  bndval = this%recharge(row) * this%auxvar(this%iauxmultcol, row)
415  else
416  bndval = this%recharge(row)
417  end if
418  case default
419  errmsg = 'Programming error. RCH bound value requested column '&
420  &'outside range of ncolbnd (1).'
421  call store_error(errmsg)
422  call store_error_filename(this%input_fname)
423  end select
424  end function rch_bound_value
425 
426 end module rchmodule
427 
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 iwetlake
integer constant for a dry lake
Definition: Constants.f90:52
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
integer(i4b), parameter maxcharlen
maximum length of char string
Definition: Constants.f90:47
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
subroutine rch_allocate_arrays(this, nodelist, auxvar)
Allocate package arrays.
Definition: gwf-rch.f90:108
logical function rch_obs_supported(this)
Overrides BndTypebnd_obs_supported()
Definition: gwf-rch.f90:376
subroutine rch_df_obs(this)
Implements bnd_df_obs.
Definition: gwf-rch.f90:389
subroutine, public rch_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
Create a New Recharge Package.
Definition: gwf-rch.f90:56
subroutine rch_fc(this, rhs, ia, idxglo, matrix_sln)
Copy rhs and hcof into solution rhs and amat.
Definition: gwf-rch.f90:298
real(dp) function rch_bound_value(this, col, row)
Return requested boundary value.
Definition: gwf-rch.f90:402
subroutine rch_allocate_scalars(this)
Allocate scalar members.
Definition: gwf-rch.f90:92
character(len=lenpackagename) texta
Definition: gwf-rch.f90:21
subroutine rch_read_initial_attr(this)
Part of allocate and read.
Definition: gwf-rch.f90:179
subroutine rch_rp(this)
Read and Prepare.
Definition: gwf-rch.f90:196
subroutine rch_define_listlabel(this)
Define the list heading that is written to iout when PRINT_INPUT option is used.
Definition: gwf-rch.f90:346
subroutine rch_cf(this)
Formulate the HCOF and RHS terms.
Definition: gwf-rch.f90:244
subroutine rch_source_options(this)
Source options specific to RchType.
Definition: gwf-rch.f90:129
subroutine rch_da(this)
Deallocate memory.
Definition: gwf-rch.f90:326
character(len=lenpackagename) text
Definition: gwf-rch.f90:20
subroutine log_rch_options(this, found_fixed_cell)
Log options specific to RchType.
Definition: gwf-rch.f90:155
subroutine set_nodesontop(this)
Store nodelist in nodesontop.
Definition: gwf-rch.f90:222
character(len=lenftype) ftype
Definition: gwf-rch.f90:19
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
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