MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
swf-pcp.f90
Go to the documentation of this file.
1 !> @brief This module contains the precipitation (PCP) package methods
2 !!
3 !! This module can be used to represent precipitation onto streams and
4 !! overland flow cells.
5 !<
7 
8  use kindmodule, only: dp, i4b, lgp
12  use bndmodule, only: bndtype
13  use bndextmodule, only: bndexttype
15  use simvariablesmodule, only: errmsg
21  use geomutilmodule, only: get_node
22  use basedismodule, only: disbasetype
23  use disv1dmodule, only: disv1dtype
24  use swfdfwmodule, only: swfdfwtype
25  use swfcxsmodule, only: swfcxstype
26 
27  implicit none
28 
29  private
30  public :: pcp_create
31 
32  character(len=LENFTYPE) :: ftype = 'PCP'
33  character(len=LENPACKAGENAME) :: text = ' PCP'
34  ! character(len=LENPACKAGENAME) :: texta = ' PCPA'
35 
36  type, extends(bndexttype) :: swfpcptype
37  real(dp), dimension(:), pointer, contiguous :: precipitation => null() !< boundary precipitation array
38  logical, pointer, private :: read_as_arrays
39 
40  ! pointers to other objects
41  type(swfdfwtype), pointer :: dfw
42  type(swfcxstype), pointer :: cxs
43 
44  contains
45 
46  procedure :: pcp_allocate_scalars
47  procedure :: allocate_arrays => pcp_allocate_arrays
48  procedure :: source_options => pcp_source_options
49  procedure :: source_dimensions => pcp_source_dimensions
50  procedure :: log_pcp_options
51  procedure :: read_initial_attr => pcp_read_initial_attr
52  procedure :: bnd_rp => pcp_rp
53  procedure :: bnd_ck => pcp_ck
54  procedure :: bnd_cf => pcp_cf
55  procedure :: bnd_fc => pcp_fc
56  procedure :: bnd_da => pcp_da
57  procedure :: define_listlabel => pcp_define_listlabel
58  procedure :: bound_value => pcp_bound_value
59  procedure, private :: default_nodelist
60  procedure, private :: reach_length_pointer
61  ! for observations
62  procedure, public :: bnd_obs_supported => pcp_obs_supported
63  procedure, public :: bnd_df_obs => pcp_df_obs
64 
65  end type swfpcptype
66 
67 contains
68 
69  !> @brief Create a Precipitation Package
70  !<
71  subroutine pcp_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
72  mempath, dis, dfw, cxs)
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
116  end subroutine pcp_create
117 
118  !> @brief Allocate scalar members
119  !<
120  subroutine pcp_allocate_scalars(this)
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.
132  end subroutine pcp_allocate_scalars
133 
134  !> @brief Allocate package arrays
135  !<
136  subroutine pcp_allocate_arrays(this, nodelist, auxvar)
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)
153  end subroutine pcp_allocate_arrays
154 
155  !> @brief Source options specific to PCPType
156  !<
157  subroutine pcp_source_options(this)
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)
175  end subroutine pcp_source_options
176 
177  !> @brief Log options specific to SwfPcpType
178  !<
179  subroutine log_pcp_options(this, found_readasarrays)
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'
199  end subroutine log_pcp_options
200 
201  !> @brief Source the dimensions for this package
202  !<
203  subroutine pcp_source_dimensions(this)
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()
231  end subroutine pcp_source_dimensions
232 
233  !> @brief Part of allocate and read
234  !<
235  subroutine pcp_read_initial_attr(this)
236  ! dummy
237  class(swfpcptype), intent(inout) :: this
238 
239  if (this%read_as_arrays) then
240  call this%default_nodelist()
241  end if
242  end subroutine pcp_read_initial_attr
243 
244  !> @brief Read and Prepare
245  !!
246  !! Read itmp and read new boundaries if itmp > 0
247  !<
248  subroutine pcp_rp(this)
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
268  end subroutine pcp_rp
269 
270  !> @brief Ensure precipitation is positive
271  !<
272  subroutine pcp_ck(this)
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
297  end subroutine pcp_ck
298 
299  !> @brief Formulate the HCOF and RHS terms
300  !!
301  !! Skip if no precipitation. Otherwise, calculate hcof and rhs
302  !<
303  subroutine pcp_cf(this)
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
370  end subroutine pcp_cf
371 
372  !> @brief Copy rhs and hcof into solution rhs and amat
373  !<
374  subroutine pcp_fc(this, rhs, ia, idxglo, matrix_sln)
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
392  end subroutine pcp_fc
393 
394  !> @brief Deallocate memory
395  !<
396  subroutine pcp_da(this)
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)
415  end subroutine pcp_da
416 
417  !> @brief Define the list heading that is written to iout when PRINT_INPUT
418  !! option is used.
419  !<
420  subroutine pcp_define_listlabel(this)
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
442  end subroutine pcp_define_listlabel
443 
444  !> @brief Assign default nodelist when READASARRAYS is specified.
445  !<
446  subroutine default_nodelist(this)
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 
463  end subroutine default_nodelist
464 
465  ! Procedures related to observations
466 
467  !> @brief
468  !!
469  !! Overrides BndType%bnd_obs_supported()
470  !<
471  logical function pcp_obs_supported(this)
472  implicit none
473  ! dummy
474  class(swfpcptype) :: this
475  pcp_obs_supported = .true.
476  end function pcp_obs_supported
477 
478  !> @brief Implements bnd_df_obs
479  !!
480  !! Store observation type supported by PCP package. Overrides
481  !! BndType%bnd_df_obs
482  !<
483  subroutine pcp_df_obs(this)
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
492  end subroutine pcp_df_obs
493 
494  !> @brief Return requested boundary value
495  !<
496  function pcp_bound_value(this, col, row) result(bndval)
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
519  end function pcp_bound_value
520 
521  function reach_length_pointer(this) result(ptr)
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 
536  end function reach_length_pointer
537 
538 end module swfpcpmodule
539 
This module contains block parser methods.
Definition: BlockParser.f90:7
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
integer(i4b), parameter maxcharlen
maximum length of char string
Definition: Constants.f90:47
integer(i4b) function, public get_node(ilay, irow, icol, nlay, nrow, ncol)
Get node number, given layer, row, and column indices for a structured grid. If any argument is inval...
Definition: GeomUtil.f90:83
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
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
This module contains the precipitation (PCP) package methods.
Definition: swf-pcp.f90:6
subroutine pcp_source_options(this)
Source options specific to PCPType.
Definition: swf-pcp.f90:158
real(dp) function, dimension(:), pointer reach_length_pointer(this)
Definition: swf-pcp.f90:522
subroutine pcp_allocate_arrays(this, nodelist, auxvar)
Allocate package arrays.
Definition: swf-pcp.f90:137
subroutine pcp_rp(this)
Read and Prepare.
Definition: swf-pcp.f90:249
subroutine log_pcp_options(this, found_readasarrays)
Log options specific to SwfPcpType.
Definition: swf-pcp.f90:180
subroutine pcp_fc(this, rhs, ia, idxglo, matrix_sln)
Copy rhs and hcof into solution rhs and amat.
Definition: swf-pcp.f90:375
character(len=lenftype) ftype
Definition: swf-pcp.f90:32
real(dp) function pcp_bound_value(this, col, row)
Return requested boundary value.
Definition: swf-pcp.f90:497
subroutine pcp_da(this)
Deallocate memory.
Definition: swf-pcp.f90:397
logical function pcp_obs_supported(this)
Overrides BndTypebnd_obs_supported()
Definition: swf-pcp.f90:472
subroutine default_nodelist(this)
Assign default nodelist when READASARRAYS is specified.
Definition: swf-pcp.f90:447
subroutine pcp_allocate_scalars(this)
Allocate scalar members.
Definition: swf-pcp.f90:121
subroutine pcp_read_initial_attr(this)
Part of allocate and read.
Definition: swf-pcp.f90:236
subroutine pcp_define_listlabel(this)
Define the list heading that is written to iout when PRINT_INPUT option is used.
Definition: swf-pcp.f90:421
subroutine pcp_df_obs(this)
Implements bnd_df_obs.
Definition: swf-pcp.f90:484
subroutine pcp_source_dimensions(this)
Source the dimensions for this package.
Definition: swf-pcp.f90:204
subroutine, public pcp_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath, dis, dfw, cxs)
Create a Precipitation Package.
Definition: swf-pcp.f90:73
subroutine pcp_ck(this)
Ensure precipitation is positive.
Definition: swf-pcp.f90:273
character(len=lenpackagename) text
Definition: swf-pcp.f90:33
subroutine pcp_cf(this)
Formulate the HCOF and RHS terms.
Definition: swf-pcp.f90:304
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