27 use swfdfwmodule,
only: swfdfwtype
35 character(len=LENFTYPE) ::
ftype =
'EVP'
36 character(len=LENPACKAGENAME) ::
text =
' EVP'
40 real(dp),
dimension(:),
pointer,
contiguous :: evaporation => null()
41 integer(I4B),
pointer :: iflowred => null()
42 real(dp),
pointer :: reduction_depth => null()
43 logical,
pointer,
private :: read_as_arrays
46 type(swfdfwtype),
pointer :: dfw
78 subroutine evp_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
79 mempath, dis, dfw, cxs)
81 class(
bndtype),
pointer :: packobj
82 integer(I4B),
intent(in) :: id
83 integer(I4B),
intent(in) :: ibcnum
84 integer(I4B),
intent(in) :: inunit
85 integer(I4B),
intent(in) :: iout
86 character(len=*),
intent(in) :: namemodel
87 character(len=*),
intent(in) :: pakname
88 character(len=*),
intent(in) :: mempath
90 type(swfdfwtype),
pointer,
intent(in) :: dfw
100 call packobj%set_names(ibcnum, namemodel, pakname,
ftype, mempath)
104 call evpobj%evp_allocate_scalars()
107 call packobj%pack_initialize()
109 packobj%inunit = inunit
112 packobj%ibcnum = ibcnum
133 call this%BndExtType%allocate_scalars()
136 call mem_allocate(this%iflowred,
'IFLOWRED', this%memoryPath)
137 call mem_allocate(this%reduction_depth,
'REDUCTION_DEPTH', this%memoryPath)
138 allocate (this%read_as_arrays)
142 this%reduction_depth =
dem6
143 this%read_as_arrays = .false.
153 integer(I4B),
dimension(:),
pointer,
contiguous,
optional :: nodelist
154 real(DP),
dimension(:, :),
pointer,
contiguous,
optional :: auxvar
157 call this%BndExtType%allocate_arrays(nodelist, auxvar)
160 call mem_setptr(this%evaporation,
'EVAPORATION', this%input_mempath)
163 call mem_checkin(this%evaporation,
'EVAPORATION', this%memoryPath, &
164 'EVAPORATION', this%input_mempath)
176 logical(LGP) :: found_readasarrays = .false.
179 call this%BndExtType%source_options()
182 call mem_set_value(this%read_as_arrays,
'READASARRAYS', this%input_mempath, &
186 call this%log_evp_options(found_readasarrays)
195 logical(LGP),
intent(in) :: found_readasarrays
197 character(len=*),
parameter :: fmtreadasarrays = &
198 &
"(4x, 'EVAPORATION INPUT WILL BE READ AS ARRAY(S).')"
201 write (this%iout,
'(/1x,a)')
'PROCESSING '//trim(adjustl(this%text)) &
204 if (found_readasarrays)
then
205 write (this%iout, fmtreadasarrays)
209 write (this%iout,
'(1x,a)') &
210 'END OF '//trim(adjustl(this%text))//
' OPTIONS'
219 if (this%read_as_arrays)
then
223 this%maxbound = this%dis%get_ncpl()
226 if (this%maxbound <= 0)
then
228 'MAXBOUND must be an integer greater than zero.'
236 call this%BndExtType%source_dimensions()
242 call this%define_listlabel()
251 if (this%read_as_arrays)
then
252 call this%default_nodelist()
267 if (this%iper /=
kper)
return
269 if (this%read_as_arrays)
then
273 call this%BndExtType%bnd_rp()
283 character(len=30) :: nodestr
284 integer(I4B) :: i, nr
285 character(len=*),
parameter :: fmterr = &
286 &
"('Specified stress ',i0, &
287 &' evaporation (',g0,') is less than zero for cell', a)"
290 do i = 1, this%nbound
291 nr = this%nodelist(i)
293 if (this%evaporation(i) <
dzero)
then
294 call this%dis%noder_to_string(nr, nodestr)
295 write (
errmsg, fmt=fmterr) i, this%evaporation(i), trim(nodestr)
323 real(DP),
dimension(:),
pointer :: reach_length
326 if (this%nbound == 0)
return
329 reach_length => this%reach_length_pointer()
333 do i = 1, this%nbound
336 node = this%nodelist(i)
346 if (this%ibound(node) <= 0)
then
356 evap = this%evaporation(i)
357 if (this%iauxmultcol > 0)
then
358 evap = evap * this%auxvar(this%iauxmultcol, i)
362 if (this%dis%is_1d())
then
363 rlen = reach_length(node)
367 q = -this%get_qevp(node, rlen, this%xnew(node), this%xold(node), evap)
376 qeps = -this%get_qevp(node, rlen, this%xnew(node) + eps, &
377 this%xold(node), evap)
380 derv = (qeps - q) / eps
384 this%rhs(i) = this%rhs(i) + derv * this%xnew(node)
399 function get_qevp(this, node, rlen, snew, sold, evaporation)
result(qevp)
402 integer(I4B),
intent(in) :: node
403 real(dp),
intent(in) :: rlen
404 real(dp),
intent(in) :: snew
405 real(dp),
intent(in) :: sold
406 real(dp),
intent(in) :: evaporation
410 integer(I4B) :: idcxs
418 real(dp) :: width_channel
423 bt = this%dis%bot(node)
426 if (this%dis%is_2d())
then
428 area = this%dis%get_area(node)
429 else if (this%dis%is_1d())
then
431 idcxs = this%dfw%idcxs(node)
432 call this%dis%get_flow_width(node, node, 0, width_channel, dummy)
435 anew = this%cxs%get_area(idcxs, width_channel, depth)
437 aold = this%cxs%get_area(idcxs, width_channel, depth)
438 wavg = this%cxs%get_wetted_top_width(idcxs, width_channel, depth)
440 if (abs(denom) >
dprec)
then
441 wavg = (anew - aold) / (snew - sold)
448 qmult = this%get_evap_reduce_mult(snew, bt)
451 qevp = evaporation * area * qmult
460 real(dp),
intent(in) :: stage
461 real(dp),
intent(in) :: bottom
468 if (this%iflowred == 1)
then
469 tp = bottom + this%reduction_depth
477 subroutine evp_fc(this, rhs, ia, idxglo, matrix_sln)
480 real(DP),
dimension(:),
intent(inout) :: rhs
481 integer(I4B),
dimension(:),
intent(in) :: ia
482 integer(I4B),
dimension(:),
intent(in) :: idxglo
485 integer(I4B) :: i, n, ipos
488 do i = 1, this%nbound
491 rhs(n) = rhs(n) + this%rhs(i)
493 call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
505 call this%BndExtType%bnd_da()
510 deallocate (this%read_as_arrays)
513 call mem_deallocate(this%evaporation,
'EVAPORATION', this%memoryPath)
529 this%listlabel = trim(this%filtyp)//
' NO.'
530 if (this%dis%ndim == 3)
then
531 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'LAYER'
532 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'ROW'
533 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'COL'
534 elseif (this%dis%ndim == 2)
then
535 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'LAYER'
536 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'CELL2D'
538 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'NODE'
540 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'EVAPORATION'
543 if (this%inamedbound == 1)
then
544 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'BOUNDARY NAME'
554 integer(I4B) :: nodeu, noder
559 do nodeu = 1, this%maxbound
560 noder = this%dis%get_nodenumber(nodeu, 0)
561 this%nodelist(nodeu) = noder
565 this%nbound = this%maxbound
594 call this%obs%StoreObsType(
'evp', .true., indx)
605 integer(I4B),
intent(in) :: col
606 integer(I4B),
intent(in) :: row
612 if (this%iauxmultcol > 0)
then
613 bndval = this%evaporation(row) * this%auxvar(this%iauxmultcol, row)
615 bndval = this%evaporation(row)
618 errmsg =
'Programming error. EVP bound value requested column '&
619 &
'outside range of ncolbnd (1).'
629 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
real(dp), parameter dhalf
real constant 1/2
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
real(dp), parameter dem6
real constant 1e-6
real(dp), parameter dzero
real constant zero
real(dp), parameter dprec
real constant machine precision
integer(i4b), parameter maxcharlen
maximum length of char string
real(dp), parameter done
real constant 1
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.
real(dp) function, public get_perturbation(x)
Calculate a numerical perturbation given the value of x.
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
real(dp) function sqsaturation(top, bot, x, c1, c2)
@ brief sQSaturation
This module contains the evaporation (EVP) package methods.
character(len=lenpackagename) text
real(dp) function get_qevp(this, node, rlen, snew, sold, evaporation)
Calculate qevp.
subroutine evp_da(this)
Deallocate memory.
subroutine evp_allocate_scalars(this)
Allocate scalar members.
subroutine evp_df_obs(this)
Implements bnd_df_obs.
subroutine evp_ck(this)
Ensure evaporation is positive.
subroutine evp_cf(this)
Formulate the HCOF and RHS terms.
logical function evp_obs_supported(this)
Overrides BndTypebnd_obs_supported()
subroutine evp_fc(this, rhs, ia, idxglo, matrix_sln)
Copy rhs and hcof into solution rhs and amat.
subroutine, public evp_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath, dis, dfw, cxs)
Create a Evaporation Package.
subroutine evp_rp(this)
Read and Prepare.
subroutine evp_define_listlabel(this)
Define the list heading that is written to iout when PRINT_INPUT option is used.
subroutine evp_allocate_arrays(this, nodelist, auxvar)
Allocate package arrays.
subroutine evp_source_dimensions(this)
Source the dimensions for this package.
real(dp) function evp_bound_value(this, col, row)
Return requested boundary value.
real(dp) function, dimension(:), pointer reach_length_pointer(this)
subroutine evp_source_options(this)
Source options specific to EVPType.
subroutine evp_read_initial_attr(this)
Part of allocate and read.
subroutine default_nodelist(this)
Assign default nodelist when READASARRAYS is specified.
subroutine log_evp_options(this, found_readasarrays)
Log options specific to SwfEvpType.
real(dp) function get_evap_reduce_mult(this, stage, bottom)
Calculate multiplier to reduce evap as depth goes to zero.
character(len=lenftype) ftype
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 ...