24 use swfdfwmodule,
only: swfdfwtype
32 character(len=LENFTYPE) ::
ftype =
'PCP'
33 character(len=LENPACKAGENAME) ::
text =
' PCP'
37 real(dp),
dimension(:),
pointer,
contiguous :: precipitation => null()
38 logical,
pointer,
private :: read_as_arrays
41 type(swfdfwtype),
pointer :: dfw
71 subroutine pcp_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
72 mempath, dis, dfw, cxs)
74 class(
bndtype),
pointer :: packobj
75 integer(I4B),
intent(in) :: id
76 integer(I4B),
intent(in) :: ibcnum
77 integer(I4B),
intent(in) :: inunit
78 integer(I4B),
intent(in) :: iout
79 character(len=*),
intent(in) :: namemodel
80 character(len=*),
intent(in) :: pakname
81 character(len=*),
intent(in) :: mempath
83 type(swfdfwtype),
pointer,
intent(in) :: dfw
93 call packobj%set_names(ibcnum, namemodel, pakname,
ftype, mempath)
97 call pcpobj%pcp_allocate_scalars()
100 call packobj%pack_initialize()
102 packobj%inunit = inunit
105 packobj%ibcnum = ibcnum
125 call this%BndExtType%allocate_scalars()
128 allocate (this%read_as_arrays)
131 this%read_as_arrays = .false.
141 integer(I4B),
dimension(:),
pointer,
contiguous,
optional :: nodelist
142 real(DP),
dimension(:, :),
pointer,
contiguous,
optional :: auxvar
145 call this%BndExtType%allocate_arrays(nodelist, auxvar)
148 call mem_setptr(this%precipitation,
'PRECIPITATION', this%input_mempath)
151 call mem_checkin(this%precipitation,
'PRECIPITATION', this%memoryPath, &
152 'PRECIPITATION', this%input_mempath)
164 logical(LGP) :: found_readasarrays = .false.
167 call this%BndExtType%source_options()
170 call mem_set_value(this%read_as_arrays,
'READASARRAYS', this%input_mempath, &
174 call this%log_pcp_options(found_readasarrays)
183 logical(LGP),
intent(in) :: found_readasarrays
185 character(len=*),
parameter :: fmtreadasarrays = &
186 &
"(4x, 'PRECIPITATION INPUT WILL BE READ AS ARRAY(S).')"
189 write (this%iout,
'(/1x,a)')
'PROCESSING '//trim(adjustl(this%text)) &
192 if (found_readasarrays)
then
193 write (this%iout, fmtreadasarrays)
197 write (this%iout,
'(1x,a)') &
198 'END OF '//trim(adjustl(this%text))//
' OPTIONS'
207 if (this%read_as_arrays)
then
211 this%maxbound = this%dis%get_ncpl()
214 if (this%maxbound <= 0)
then
216 'MAXBOUND must be an integer greater than zero.'
224 call this%BndExtType%source_dimensions()
230 call this%define_listlabel()
239 if (this%read_as_arrays)
then
240 call this%default_nodelist()
255 if (this%iper /=
kper)
return
257 if (this%read_as_arrays)
then
261 call this%BndExtType%bnd_rp()
265 if (this%iprpak /= 0)
then
266 call this%write_list()
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)"
283 do i = 1, this%nbound
284 nr = this%nodelist(i)
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)
309 integer(I4B) :: idcxs
312 real(DP) :: width_channel
313 real(DP) :: top_width
315 real(DP),
dimension(:),
pointer :: reach_length
318 if (this%nbound == 0)
return
321 reach_length => this%reach_length_pointer()
324 do i = 1, this%nbound
327 node = this%nodelist(i)
340 if (this%dis%is_2d())
then
342 area = this%dis%get_area(node)
343 else if (this%dis%is_1d())
then
345 idcxs = this%dfw%idcxs(node)
346 call this%dis%get_flow_width(node, node, 0, width_channel, &
348 top_width = this%cxs%get_maximum_top_width(idcxs, width_channel)
349 area = reach_length(node) * top_width
353 qpcp = this%precipitation(i) * area
356 if (this%iauxmultcol > 0)
then
357 qpcp = qpcp * this%auxvar(this%iauxmultcol, i)
364 if (this%ibound(node) <= 0)
then
374 subroutine pcp_fc(this, rhs, ia, idxglo, matrix_sln)
377 real(DP),
dimension(:),
intent(inout) :: rhs
378 integer(I4B),
dimension(:),
intent(in) :: ia
379 integer(I4B),
dimension(:),
intent(in) :: idxglo
382 integer(I4B) :: i, n, ipos
385 do i = 1, this%nbound
388 rhs(n) = rhs(n) + this%rhs(i)
390 call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
403 call this%BndExtType%bnd_da()
406 deallocate (this%read_as_arrays)
409 call mem_deallocate(this%precipitation,
'PRECIPITATION', this%memoryPath)
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'
434 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'NODE'
436 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'PRECIPITATION'
439 if (this%inamedbound == 1)
then
440 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'BOUNDARY NAME'
450 integer(I4B) :: nodeu, noder
455 do nodeu = 1, this%maxbound
456 noder = this%dis%get_nodenumber(nodeu, 0)
457 this%nodelist(nodeu) = noder
461 this%nbound = this%maxbound
490 call this%obs%StoreObsType(
'pcp', .true., indx)
501 integer(I4B),
intent(in) :: col
502 integer(I4B),
intent(in) :: row
508 if (this%iauxmultcol > 0)
then
509 bndval = this%precipitation(row) * this%auxvar(this%iauxmultcol, row)
511 bndval = this%precipitation(row)
514 errmsg =
'Programming error. PCP bound value requested column '&
515 &
'outside range of ncolbnd (1).'
525 real(dp),
dimension(:),
pointer :: ptr
This module contains block parser methods.
This module contains the extended boundary package.
This module contains the base boundary package.
This module contains simulation constants.
integer(i4b), parameter linelength
maximum length of a standard line
integer(i4b), parameter lenpackagename
maximum length of the package name
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
real(dp), parameter dzero
real constant zero
integer(i4b), parameter maxcharlen
maximum length of char string
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...
This module defines variable data types.
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains the derived type ObsType.
subroutine, public defaultobsidprocessor(obsrv, dis, inunitobs, iout)
@ brief Process IDstring provided for each observation
This module contains simulation methods.
subroutine, public store_error(msg, terminate)
Store an error message.
integer(i4b) function, public count_errors()
Return number of errors.
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
This module contains simulation variables.
character(len=maxcharlen) errmsg
error message string
This module contains the precipitation (PCP) package methods.
subroutine pcp_source_options(this)
Source options specific to PCPType.
real(dp) function, dimension(:), pointer reach_length_pointer(this)
subroutine pcp_allocate_arrays(this, nodelist, auxvar)
Allocate package arrays.
subroutine pcp_rp(this)
Read and Prepare.
subroutine log_pcp_options(this, found_readasarrays)
Log options specific to SwfPcpType.
subroutine pcp_fc(this, rhs, ia, idxglo, matrix_sln)
Copy rhs and hcof into solution rhs and amat.
character(len=lenftype) ftype
real(dp) function pcp_bound_value(this, col, row)
Return requested boundary value.
subroutine pcp_da(this)
Deallocate memory.
logical function pcp_obs_supported(this)
Overrides BndTypebnd_obs_supported()
subroutine default_nodelist(this)
Assign default nodelist when READASARRAYS is specified.
subroutine pcp_allocate_scalars(this)
Allocate scalar members.
subroutine pcp_read_initial_attr(this)
Part of allocate and read.
subroutine pcp_define_listlabel(this)
Define the list heading that is written to iout when PRINT_INPUT option is used.
subroutine pcp_df_obs(this)
Implements bnd_df_obs.
subroutine pcp_source_dimensions(this)
Source the dimensions for this package.
subroutine, public pcp_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath, dis, dfw, cxs)
Create a Precipitation Package.
subroutine pcp_ck(this)
Ensure precipitation is positive.
character(len=lenpackagename) text
subroutine pcp_cf(this)
Formulate the HCOF and RHS terms.
integer(i4b), pointer, public kper
current stress period number
This class is used to store a single deferred-length character string. It was designed to work in an ...