MODFLOW 6  version 6.8.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 442 of file swf-pcp.f90.

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

◆ log_pcp_options()

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

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

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

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

◆ pcp_allocate_scalars()

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

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

122  ! dummy
123  class(SwfPcpType), intent(inout) :: this
124 
125  ! allocate base scalars
126  call this%BndExtType%allocate_scalars()
127 
128  ! allocate internal members
129  allocate (this%read_as_arrays)
130 
131  ! Set values
132  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

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

493  ! modules
494  use constantsmodule, only: dzero
495  ! dummy
496  class(SwfPcpType), intent(inout) :: this
497  integer(I4B), intent(in) :: col
498  integer(I4B), intent(in) :: row
499  ! result
500  real(DP) :: bndval
501 
502  select case (col)
503  case (1)
504  if (this%iauxmultcol > 0) then
505  bndval = this%precipitation(row) * this%auxvar(this%iauxmultcol, row)
506  else
507  bndval = this%precipitation(row)
508  end if
509  case default
510  errmsg = 'Programming error. PCP bound value requested column '&
511  &'outside range of ncolbnd (1).'
512  call store_error(errmsg)
513  call store_error_filename(this%input_fname)
514  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 299 of file swf-pcp.f90.

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

◆ pcp_ck()

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

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

269  ! dummy
270  class(SwfPcpType), intent(inout) :: this
271  ! local
272  character(len=30) :: nodestr
273  integer(I4B) :: i, nr
274  character(len=*), parameter :: fmterr = &
275  &"('Specified stress ',i0, &
276  &' precipitation (',g0,') is less than zero for cell', a)"
277 
278  ! Ensure precipitation rates are positive
279  do i = 1, this%nbound
280  nr = this%nodelist(i)
281  if (nr <= 0) cycle
282  if (this%precipitation(i) < dzero) then
283  call this%dis%noder_to_string(nr, nodestr)
284  write (errmsg, fmt=fmterr) i, this%precipitation(i), trim(nodestr)
285  call store_error(errmsg)
286  end if
287  end do
288 
289  ! write summary of package error messages
290  if (count_errors() > 0) then
291  call store_error_filename(this%input_fname)
292  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%ncolbnd = 1
107  packobj%ictMemPath = create_mem_path(namemodel, 'DFW')
108 
109  ! store pointer to dis
110  pcpobj%dis => dis
111 
112  ! store pointer to dfw
113  pcpobj%dfw => dfw
114 
115  ! store pointer to cxs
116  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 392 of file swf-pcp.f90.

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

◆ pcp_define_listlabel()

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

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

417  ! dummy
418  class(SwfPcpType), intent(inout) :: this
419  !
420  ! create the header list label
421  this%listlabel = trim(this%filtyp)//' NO.'
422  if (this%dis%ndim == 3) then
423  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
424  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
425  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
426  elseif (this%dis%ndim == 2) then
427  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
428  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
429  else
430  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
431  end if
432  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'PRECIPITATION'
433 ! if(this%multindex > 0) &
434 ! write(this%listlabel, '(a, a16)') trim(this%listlabel), 'MULTIPLIER'
435  if (this%inamedbound == 1) then
436  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
437  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 479 of file swf-pcp.f90.

480  implicit none
481  ! dummy
482  class(SwfPcpType) :: this
483  ! local
484  integer(I4B) :: indx
485 
486  call this%obs%StoreObsType('pcp', .true., indx)
487  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 370 of file swf-pcp.f90.

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

◆ pcp_obs_supported()

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

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

468  implicit none
469  ! dummy
470  class(SwfPcpType) :: this
471  pcp_obs_supported = .true.

◆ pcp_read_initial_attr()

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

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

237  ! dummy
238  class(SwfPcpType), intent(inout) :: this
239 
240  if (this%read_as_arrays) then
241  call this%default_nodelist()
242  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 249 of file swf-pcp.f90.

250  ! modules
251  use tdismodule, only: kper
252  implicit none
253  ! dummy
254  class(SwfPcpType), intent(inout) :: this
255 
256  if (this%iper /= kper) return
257 
258  if (this%read_as_arrays) then
259  ! no need to do anything because this%precipitation points directly to
260  ! the input context precipitation, which is automatically updated by idm
261  else
262  call this%BndExtType%bnd_rp()
263  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 204 of file swf-pcp.f90.

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

159  ! modules
161  implicit none
162  ! dummy
163  class(SwfPcpType), intent(inout) :: this
164  ! local
165  logical(LGP) :: found_readasarrays = .false.
166 
167  ! source common bound options
168  call this%BndExtType%source_options()
169 
170  ! update defaults with idm sourced values
171  call mem_set_value(this%read_as_arrays, 'READASARRAYS', this%input_mempath, &
172  found_readasarrays)
173 
174  ! log pcp params
175  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 517 of file swf-pcp.f90.

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

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'