MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
swfpcpmodule Module Reference

This module contains the precipitation (PCP) package methods. More...

Data Types

type  swfpcptype
 

Functions/Subroutines

subroutine, public pcp_create (packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath, dis, dfw, cxs)
 Create a Precipitation Package. More...
 
subroutine pcp_allocate_scalars (this)
 Allocate scalar members. More...
 
subroutine pcp_allocate_arrays (this, nodelist, auxvar)
 Allocate package arrays. More...
 
subroutine pcp_source_options (this)
 Source options specific to PCPType. More...
 
subroutine log_pcp_options (this, found_readasarrays)
 Log options specific to SwfPcpType. More...
 
subroutine pcp_source_dimensions (this)
 Source the dimensions for this package. More...
 
subroutine pcp_read_initial_attr (this)
 Part of allocate and read. More...
 
subroutine pcp_rp (this)
 Read and Prepare. More...
 
subroutine pcp_ck (this)
 Ensure precipitation is positive. More...
 
subroutine pcp_cf (this)
 Formulate the HCOF and RHS terms. More...
 
subroutine pcp_fc (this, rhs, ia, idxglo, matrix_sln)
 Copy rhs and hcof into solution rhs and amat. More...
 
subroutine pcp_da (this)
 Deallocate memory. More...
 
subroutine pcp_define_listlabel (this)
 Define the list heading that is written to iout when PRINT_INPUT option is used. More...
 
subroutine default_nodelist (this)
 Assign default nodelist when READASARRAYS is specified. More...
 
logical function pcp_obs_supported (this)
 Overrides BndTypebnd_obs_supported() More...
 
subroutine pcp_df_obs (this)
 Implements bnd_df_obs. More...
 
real(dp) function pcp_bound_value (this, col, row)
 Return requested boundary value. More...
 
real(dp) function, dimension(:), pointer reach_length_pointer (this)
 

Variables

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

Detailed Description

This module can be used to represent precipitation onto streams and overland flow cells.

Function/Subroutine Documentation

◆ default_nodelist()

subroutine swfpcpmodule::default_nodelist ( class(swfpcptype this)
private

Definition at line 446 of file swf-pcp.f90.

447  ! dummy
448  class(SwfPcpType) :: this
449  ! local
450  integer(I4B) :: nodeu, noder
451 
452  ! This is only called for readasarrays, so nodelist will be the size of
453  ! the user grid, and will have a value of 0 for any entries where idomain
454  ! is not 1
455  do nodeu = 1, this%maxbound
456  noder = this%dis%get_nodenumber(nodeu, 0)
457  this%nodelist(nodeu) = noder
458  end do
459 
460  ! Assign nbound
461  this%nbound = this%maxbound
462 

◆ log_pcp_options()

subroutine swfpcpmodule::log_pcp_options ( class(swfpcptype), intent(inout)  this,
logical(lgp), intent(in)  found_readasarrays 
)

Definition at line 179 of file swf-pcp.f90.

180  implicit none
181  ! dummy
182  class(SwfPcpType), intent(inout) :: this
183  logical(LGP), intent(in) :: found_readasarrays
184  ! formats
185  character(len=*), parameter :: fmtreadasarrays = &
186  &"(4x, 'PRECIPITATION INPUT WILL BE READ AS ARRAY(S).')"
187 
188  ! log found options
189  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
190  //' OPTIONS'
191 
192  if (found_readasarrays) then
193  write (this%iout, fmtreadasarrays)
194  end if
195 
196  ! close logging block
197  write (this%iout, '(1x,a)') &
198  'END OF '//trim(adjustl(this%text))//' OPTIONS'

◆ pcp_allocate_arrays()

subroutine swfpcpmodule::pcp_allocate_arrays ( class(swfpcptype this,
integer(i4b), dimension(:), optional, pointer, contiguous  nodelist,
real(dp), dimension(:, :), optional, pointer, contiguous  auxvar 
)
private

Definition at line 136 of file swf-pcp.f90.

137  ! modules
139  ! dummy
140  class(SwfPcpType) :: this
141  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
142  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
143 
144  ! allocate base arrays
145  call this%BndExtType%allocate_arrays(nodelist, auxvar)
146 
147  ! set input context pointers
148  call mem_setptr(this%precipitation, 'PRECIPITATION', this%input_mempath)
149 
150  ! checkin input context pointers
151  call mem_checkin(this%precipitation, 'PRECIPITATION', this%memoryPath, &
152  'PRECIPITATION', this%input_mempath)

◆ pcp_allocate_scalars()

subroutine swfpcpmodule::pcp_allocate_scalars ( class(swfpcptype), intent(inout)  this)
private

Definition at line 120 of file swf-pcp.f90.

121  ! dummy
122  class(SwfPcpType), intent(inout) :: this
123 
124  ! allocate base scalars
125  call this%BndExtType%allocate_scalars()
126 
127  ! allocate internal members
128  allocate (this%read_as_arrays)
129 
130  ! Set values
131  this%read_as_arrays = .false.

◆ pcp_bound_value()

real(dp) function swfpcpmodule::pcp_bound_value ( class(swfpcptype), intent(inout)  this,
integer(i4b), intent(in)  col,
integer(i4b), intent(in)  row 
)
private
Parameters
[in,out]thisBndExtType object

Definition at line 496 of file swf-pcp.f90.

497  ! modules
498  use constantsmodule, only: dzero
499  ! dummy
500  class(SwfPcpType), intent(inout) :: this !< BndExtType object
501  integer(I4B), intent(in) :: col
502  integer(I4B), intent(in) :: row
503  ! result
504  real(DP) :: bndval
505 
506  select case (col)
507  case (1)
508  if (this%iauxmultcol > 0) then
509  bndval = this%precipitation(row) * this%auxvar(this%iauxmultcol, row)
510  else
511  bndval = this%precipitation(row)
512  end if
513  case default
514  errmsg = 'Programming error. PCP bound value requested column '&
515  &'outside range of ncolbnd (1).'
516  call store_error(errmsg)
517  call store_error_filename(this%input_fname)
518  end select
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
Here is the call graph for this function:

◆ pcp_cf()

subroutine swfpcpmodule::pcp_cf ( class(swfpcptype this)
private

Skip if no precipitation. Otherwise, calculate hcof and rhs

Definition at line 303 of file swf-pcp.f90.

304  ! dummy
305  class(SwfPcpType) :: this
306  ! local
307  integer(I4B) :: i
308  integer(I4B) :: node
309  integer(I4B) :: idcxs
310  real(DP) :: qpcp
311  real(DP) :: area
312  real(DP) :: width_channel
313  real(DP) :: top_width
314  real(DP) :: dummy
315  real(DP), dimension(:), pointer :: reach_length
316 
317  ! Return if no precipitation
318  if (this%nbound == 0) return
319 
320  ! Set pointer to reach_length for 1d
321  reach_length => this%reach_length_pointer()
322 
323  ! Calculate hcof and rhs for each precipitation entry
324  do i = 1, this%nbound
325 
326  ! Find the node number
327  node = this%nodelist(i)
328 
329  ! cycle if nonexistent bound
330  if (node <= 0) then
331  this%hcof(i) = dzero
332  this%rhs(i) = dzero
333  cycle
334  end if
335 
336  ! Initialize hcof
337  this%hcof(i) = dzero
338 
339  ! Determine the water surface area
340  if (this%dis%is_2d()) then
341  ! this is for overland flow case
342  area = this%dis%get_area(node)
343  else if (this%dis%is_1d()) then
344  ! this is for channel case
345  idcxs = this%dfw%idcxs(node)
346  call this%dis%get_flow_width(node, node, 0, width_channel, &
347  dummy)
348  top_width = this%cxs%get_maximum_top_width(idcxs, width_channel)
349  area = reach_length(node) * top_width
350  end if
351 
352  ! calculate volumetric precipitation flow in L^3/T
353  qpcp = this%precipitation(i) * area
354 
355  ! multiplier
356  if (this%iauxmultcol > 0) then
357  qpcp = qpcp * this%auxvar(this%iauxmultcol, i)
358  end if
359 
360  ! rhs contribution
361  this%rhs(i) = -qpcp
362 
363  ! zero out contribution if cell is inactive or constant head
364  if (this%ibound(node) <= 0) then
365  this%rhs(i) = dzero
366  cycle
367  end if
368 
369  end do

◆ pcp_ck()

subroutine swfpcpmodule::pcp_ck ( class(swfpcptype), intent(inout)  this)

Definition at line 272 of file swf-pcp.f90.

273  ! dummy
274  class(SwfPcpType), intent(inout) :: this
275  ! local
276  character(len=30) :: nodestr
277  integer(I4B) :: i, nr
278  character(len=*), parameter :: fmterr = &
279  &"('Specified stress ',i0, &
280  &' precipitation (',g0,') is less than zero for cell', a)"
281 
282  ! Ensure precipitation rates are positive
283  do i = 1, this%nbound
284  nr = this%nodelist(i)
285  if (nr <= 0) cycle
286  if (this%precipitation(i) < dzero) then
287  call this%dis%noder_to_string(nr, nodestr)
288  write (errmsg, fmt=fmterr) i, this%precipitation(i), trim(nodestr)
289  call store_error(errmsg)
290  end if
291  end do
292 
293  ! write summary of package error messages
294  if (count_errors() > 0) then
295  call store_error_filename(this%input_fname)
296  end if
Here is the call graph for this function:

◆ pcp_create()

subroutine, public swfpcpmodule::pcp_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,
class(disbasetype), intent(inout), pointer  dis,
type(swfdfwtype), intent(in), pointer  dfw,
type(swfcxstype), intent(in), pointer  cxs 
)
Parameters
packobjpointer to default package type
[in]idpackage id
[in]ibcnumboundary condition number
[in]inunitunit number of CDB package input file
[in]ioutunit number of model listing file
[in]namemodelmodel name
[in]paknamepackage name
[in]mempathinput mempath
[in,out]disthe pointer to the discretization
[in]dfwthe pointer to the dfw package
[in]cxsthe pointer to the cxs package

Definition at line 71 of file swf-pcp.f90.

73  ! dummy
74  class(BndType), pointer :: packobj !< pointer to default package type
75  integer(I4B), intent(in) :: id !< package id
76  integer(I4B), intent(in) :: ibcnum !< boundary condition number
77  integer(I4B), intent(in) :: inunit !< unit number of CDB package input file
78  integer(I4B), intent(in) :: iout !< unit number of model listing file
79  character(len=*), intent(in) :: namemodel !< model name
80  character(len=*), intent(in) :: pakname !< package name
81  character(len=*), intent(in) :: mempath !< input mempath
82  class(DisBaseType), pointer, intent(inout) :: dis !< the pointer to the discretization
83  type(SwfDfwType), pointer, intent(in) :: dfw !< the pointer to the dfw package
84  type(SwfCxsType), pointer, intent(in) :: cxs !< the pointer to the cxs package
85  ! local
86  type(SwfPcpType), pointer :: pcpobj
87 
88  ! allocate precipitation object and scalar variables
89  allocate (pcpobj)
90  packobj => pcpobj
91 
92  ! create name and memory path
93  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
94  packobj%text = text
95 
96  ! allocate scalars
97  call pcpobj%pcp_allocate_scalars()
98 
99  ! initialize package
100  call packobj%pack_initialize()
101 
102  packobj%inunit = inunit
103  packobj%iout = iout
104  packobj%id = id
105  packobj%ibcnum = ibcnum
106  packobj%ictMemPath = create_mem_path(namemodel, 'DFW')
107 
108  ! store pointer to dis
109  pcpobj%dis => dis
110 
111  ! store pointer to dfw
112  pcpobj%dfw => dfw
113 
114  ! store pointer to cxs
115  pcpobj%cxs => cxs
Here is the call graph for this function:
Here is the caller graph for this function:

◆ pcp_da()

subroutine swfpcpmodule::pcp_da ( class(swfpcptype this)
private

Definition at line 396 of file swf-pcp.f90.

397  ! modules
399  ! dummy
400  class(SwfPcpType) :: this
401 
402  ! Deallocate parent package
403  call this%BndExtType%bnd_da()
404 
405  ! scalars
406  deallocate (this%read_as_arrays)
407 
408  ! arrays
409  call mem_deallocate(this%precipitation, 'PRECIPITATION', this%memoryPath)
410 
411  ! pointers
412  nullify (this%dis)
413  nullify (this%dfw)
414  nullify (this%cxs)

◆ pcp_define_listlabel()

subroutine swfpcpmodule::pcp_define_listlabel ( class(swfpcptype), intent(inout)  this)

Definition at line 420 of file swf-pcp.f90.

421  ! dummy
422  class(SwfPcpType), intent(inout) :: this
423  !
424  ! create the header list label
425  this%listlabel = trim(this%filtyp)//' NO.'
426  if (this%dis%ndim == 3) then
427  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
428  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
429  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
430  elseif (this%dis%ndim == 2) then
431  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
432  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
433  else
434  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
435  end if
436  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'PRECIPITATION'
437 ! if(this%multindex > 0) &
438 ! write(this%listlabel, '(a, a16)') trim(this%listlabel), 'MULTIPLIER'
439  if (this%inamedbound == 1) then
440  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
441  end if

◆ pcp_df_obs()

subroutine swfpcpmodule::pcp_df_obs ( class(swfpcptype this)
private

Store observation type supported by PCP package. Overrides BndTypebnd_df_obs

Definition at line 483 of file swf-pcp.f90.

484  implicit none
485  ! dummy
486  class(SwfPcpType) :: this
487  ! local
488  integer(I4B) :: indx
489 
490  call this%obs%StoreObsType('pcp', .true., indx)
491  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
Here is the call graph for this function:

◆ pcp_fc()

subroutine swfpcpmodule::pcp_fc ( class(swfpcptype 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 374 of file swf-pcp.f90.

375  ! dummy
376  class(SwfPcpType) :: this
377  real(DP), dimension(:), intent(inout) :: rhs
378  integer(I4B), dimension(:), intent(in) :: ia
379  integer(I4B), dimension(:), intent(in) :: idxglo
380  class(MatrixBaseType), pointer :: matrix_sln
381  ! local
382  integer(I4B) :: i, n, ipos
383 
384  ! Copy package rhs and hcof into solution rhs and amat
385  do i = 1, this%nbound
386  n = this%nodelist(i)
387  if (n <= 0) cycle
388  rhs(n) = rhs(n) + this%rhs(i)
389  ipos = ia(n)
390  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
391  end do

◆ pcp_obs_supported()

logical function swfpcpmodule::pcp_obs_supported ( class(swfpcptype this)
private

Definition at line 471 of file swf-pcp.f90.

472  implicit none
473  ! dummy
474  class(SwfPcpType) :: this
475  pcp_obs_supported = .true.

◆ pcp_read_initial_attr()

subroutine swfpcpmodule::pcp_read_initial_attr ( class(swfpcptype), intent(inout)  this)
private

Definition at line 235 of file swf-pcp.f90.

236  ! dummy
237  class(SwfPcpType), intent(inout) :: this
238 
239  if (this%read_as_arrays) then
240  call this%default_nodelist()
241  end if

◆ pcp_rp()

subroutine swfpcpmodule::pcp_rp ( class(swfpcptype), intent(inout)  this)
private

Read itmp and read new boundaries if itmp > 0

Definition at line 248 of file swf-pcp.f90.

249  ! modules
250  use tdismodule, only: kper
251  implicit none
252  ! dummy
253  class(SwfPcpType), intent(inout) :: this
254 
255  if (this%iper /= kper) return
256 
257  if (this%read_as_arrays) then
258  ! no need to do anything because this%precipitation points directly to
259  ! the input context precipitation, which is automatically updated by idm
260  else
261  call this%BndExtType%bnd_rp()
262  end if
263 
264  ! Write the list to iout if requested
265  if (this%iprpak /= 0) then
266  call this%write_list()
267  end if
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23

◆ pcp_source_dimensions()

subroutine swfpcpmodule::pcp_source_dimensions ( class(swfpcptype), intent(inout)  this)
private

Definition at line 203 of file swf-pcp.f90.

204  ! dummy
205  class(SwfPcpType), intent(inout) :: this
206 
207  if (this%read_as_arrays) then
208 
209  ! Set maxbound to the number of cells per layer, which is simply
210  ! nrow * ncol for a dis2d grid, and nodesuser for disv2d and disv1d
211  this%maxbound = this%dis%get_ncpl()
212 
213  ! verify dimensions were set
214  if (this%maxbound <= 0) then
215  write (errmsg, '(a)') &
216  'MAXBOUND must be an integer greater than zero.'
217  call store_error(errmsg)
218  call store_error_filename(this%input_fname)
219  end if
220 
221  else
222 
223  ! source maxbound
224  call this%BndExtType%source_dimensions()
225 
226  end if
227 
228  ! Call define_listlabel to construct the list label that is written
229  ! when PRINT_INPUT option is used.
230  call this%define_listlabel()
Here is the call graph for this function:

◆ pcp_source_options()

subroutine swfpcpmodule::pcp_source_options ( class(swfpcptype), intent(inout)  this)

Definition at line 157 of file swf-pcp.f90.

158  ! modules
160  implicit none
161  ! dummy
162  class(SwfPcpType), intent(inout) :: this
163  ! local
164  logical(LGP) :: found_readasarrays = .false.
165 
166  ! source common bound options
167  call this%BndExtType%source_options()
168 
169  ! update defaults with idm sourced values
170  call mem_set_value(this%read_as_arrays, 'READASARRAYS', this%input_mempath, &
171  found_readasarrays)
172 
173  ! log pcp params
174  call this%log_pcp_options(found_readasarrays)

◆ reach_length_pointer()

real(dp) function, dimension(:), pointer swfpcpmodule::reach_length_pointer ( class(swfpcptype this)
Parameters
thisthis instance

Definition at line 521 of file swf-pcp.f90.

522  ! dummy
523  class(SwfPcpType) :: this !< this instance
524  ! return
525  real(DP), dimension(:), pointer :: ptr
526  ! local
527  class(DisBaseType), pointer :: dis
528 
529  ptr => null()
530  dis => this%dis
531  select type (dis)
532  type is (disv1dtype)
533  ptr => dis%length
534  end select
535 

Variable Documentation

◆ ftype

character(len=lenftype) swfpcpmodule::ftype = 'PCP'
private

Definition at line 32 of file swf-pcp.f90.

32  character(len=LENFTYPE) :: ftype = 'PCP'

◆ text

character(len=lenpackagename) swfpcpmodule::text = ' PCP'
private

Definition at line 33 of file swf-pcp.f90.

33  character(len=LENPACKAGENAME) :: text = ' PCP'