19 character(len=LENFTYPE) ::
ftype =
'EVT'
20 character(len=LENPACKAGENAME) ::
text =
' EVT'
21 character(len=LENPACKAGENAME) ::
texta =
' EVTA'
25 logical,
pointer,
private :: segsdefined
26 logical,
pointer,
private :: fixed_cell
27 logical,
pointer,
private :: surfratespecified
29 integer(I4B),
pointer,
private :: nseg => null()
31 real(dp),
dimension(:),
pointer,
contiguous :: surface => null()
32 real(dp),
dimension(:),
pointer,
contiguous :: rate => null()
33 real(dp),
dimension(:),
pointer,
contiguous :: depth => null()
34 real(dp),
dimension(:, :),
pointer,
contiguous :: pxdp => null()
35 real(dp),
dimension(:, :),
pointer,
contiguous :: petm => null()
36 real(dp),
dimension(:),
pointer,
contiguous :: petm0 => null()
37 integer(I4B),
dimension(:),
pointer,
contiguous :: nodesontop => null()
63 subroutine evt_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
66 class(
bndtype),
pointer :: packobj
67 integer(I4B),
intent(in) :: id
68 integer(I4B),
intent(in) :: ibcnum
69 integer(I4B),
intent(in) :: inunit
70 integer(I4B),
intent(in) :: iout
71 character(len=*),
intent(in) :: namemodel
72 character(len=*),
intent(in) :: pakname
73 character(len=*),
intent(in) :: mempath
75 type(
evttype),
pointer :: evtobj
82 call packobj%set_names(ibcnum, namemodel, pakname,
ftype, mempath)
86 call evtobj%evt_allocate_scalars()
89 call packobj%pack_initialize()
91 packobj%inunit = inunit
94 packobj%ibcnum = ibcnum
105 class(
evttype),
intent(inout) :: this
108 call this%BndExtType%allocate_scalars()
114 allocate (this%segsdefined)
115 allocate (this%fixed_cell)
116 allocate (this%surfratespecified)
120 this%segsdefined = .true.
121 this%fixed_cell = .false.
122 this%surfratespecified = .false.
132 integer(I4B),
dimension(:),
pointer,
contiguous,
optional :: nodelist
133 real(DP),
dimension(:, :),
pointer,
contiguous,
optional :: auxvar
136 call this%BndExtType%allocate_arrays(nodelist, auxvar)
139 call mem_setptr(this%surface,
'SURFACE', this%input_mempath)
140 call mem_setptr(this%rate,
'RATE', this%input_mempath)
141 call mem_setptr(this%depth,
'DEPTH', this%input_mempath)
144 call mem_checkin(this%surface,
'SURFACE', this%memoryPath, &
145 'SURFACE', this%input_mempath)
146 call mem_checkin(this%rate,
'RATE', this%memoryPath, &
147 'RATE', this%input_mempath)
148 call mem_checkin(this%depth,
'DEPTH', this%memoryPath, &
149 'DEPTH', this%input_mempath)
152 if (.not. this%readasarrays)
then
153 if (this%nseg > 1)
then
156 call mem_setptr(this%pxdp,
'PXDP', this%input_mempath)
157 call mem_setptr(this%petm,
'PETM', this%input_mempath)
160 call mem_checkin(this%pxdp,
'PXDP', this%memoryPath, &
161 'PXDP', this%input_mempath)
162 call mem_checkin(this%petm,
'PETM', this%memoryPath, &
163 'PETM', this%input_mempath)
166 if (this%surfratespecified)
then
169 call mem_setptr(this%petm0,
'PETM0', this%input_mempath)
172 call mem_checkin(this%petm0,
'PETM0', this%memoryPath, &
173 'PETM0', this%input_mempath)
184 class(
evttype),
intent(inout) :: this
186 logical(LGP) :: found_fixed_cell = .false.
187 logical(LGP) :: found_surfratespec = .false.
190 call this%BndExtType%source_options()
194 this%input_mempath, found_fixed_cell)
196 this%input_mempath, found_surfratespec)
198 if (this%readasarrays)
then
199 if (found_surfratespec)
then
200 errmsg =
'READASARRAYS option is not compatible with the'// &
201 ' SURF_RATE_SPECIFIED option.'
210 call this%evt_log_options(found_fixed_cell, found_surfratespec)
220 class(
evttype),
intent(inout) :: this
221 logical(LGP),
intent(in) :: found_fixed_cell
222 logical(LGP),
intent(in) :: found_surfratespec
224 character(len=*),
parameter :: fmtihact = &
225 &
"(4x, 'EVAPOTRANSPIRATION WILL BE APPLIED TO HIGHEST ACTIVE CELL.')"
226 character(len=*),
parameter :: fmtfixedcell = &
227 &
"(4x, 'EVAPOTRANSPIRATION WILL BE APPLIED TO SPECIFIED CELL.')"
228 character(len=*),
parameter :: fmtsrz = &
229 &
"(4x, 'ET RATE AT SURFACE WILL BE ZERO.')"
230 character(len=*),
parameter :: fmtsrs = &
231 &
"(4x, 'ET RATE AT SURFACE WILL BE AS SPECIFIED BY PETM0.')"
234 write (this%iout,
'(/1x,a)')
'PROCESSING '//trim(adjustl(this%text)) &
237 if (found_fixed_cell)
then
238 write (this%iout, fmtfixedcell)
241 if (found_surfratespec)
then
242 write (this%iout, fmtsrs)
246 write (this%iout,
'(1x,a)') &
247 'END OF '//trim(adjustl(this%text))//
' OPTIONS'
256 class(
evttype),
intent(inout) :: this
258 logical(LGP) :: found_nseg = .false.
260 character(len=*),
parameter :: fmtnsegerr = &
261 &
"('Error: In EVT, NSEG must be > 0 but is specified as ',i0)"
264 call this%BndExtType%source_dimensions()
266 if (this%readasarrays)
then
271 write (this%iout,
'(/1x,a)')
'PROCESSING '//trim(adjustl(this%text)) &
275 call mem_set_value(this%nseg,
'NSEG', this%input_mempath, found_nseg)
279 write (this%iout,
'(4x,a,i0)')
'NSEG = ', this%nseg
281 if (this%nseg < 1)
then
282 write (
errmsg, fmtnsegerr) this%nseg
286 elseif (this%nseg > 1)
then
288 if (this%readasarrays)
then
289 errmsg =
'In the EVT package, NSEG cannot be greater than 1'// &
290 ' when READASARRAYS is used.'
299 write (this%iout,
'(1x,a)') &
300 'END OF '//trim(adjustl(this%text))//
' DIMENSIONS'
311 class(
evttype),
intent(inout) :: this
313 if (this%readasarrays)
then
314 call this%default_nodelist()
319 if (.not. this%fixed_cell)
call this%set_nodesontop()
330 class(
evttype),
intent(inout) :: this
332 if (this%iper /=
kper)
return
335 call this%BndExtType%bnd_rp()
338 if (this%nseg > 1)
then
339 call this%check_pxdp()
343 if (.not. this%fixed_cell)
call this%set_nodesontop()
354 class(
evttype),
intent(inout) :: this
359 integer(I4B) :: ierrmono
362 character(len=15) :: nodestr
364 character(len=*),
parameter :: fmtpxdp0 = &
365 &
"('PXDP must be between 0 and 1. Found ', G0, ' for cell ', A)"
366 character(len=*),
parameter :: fmtpxdp = &
367 &
"('PXDP is not monotonically increasing for cell ', A)"
371 do n = 1, this%nbound
372 node = this%nodelist(n)
375 segloop:
do i = 1, this%nseg
378 if (i < this%nseg)
then
379 pxdp2 = this%pxdp(i, n)
380 if (pxdp2 <=
dzero .or. pxdp2 >=
done)
then
381 call this%dis%noder_to_string(node, nodestr)
382 write (
errmsg, fmtpxdp0) pxdp2, trim(nodestr)
390 if (pxdp2 - pxdp1 <
dzero)
then
391 if (ierrmono == 0)
then
393 call this%dis%noder_to_string(node, nodestr)
394 write (
errmsg, fmtpxdp) trim(nodestr)
413 class(
evttype),
intent(inout) :: this
418 if (.not.
associated(this%nodesontop))
then
419 allocate (this%nodesontop(this%maxbound))
423 do n = 1, this%nbound
424 this%nodesontop(n) = this%nodelist(n)
434 integer(I4B) :: i, iseg, node
435 integer(I4B) :: idxdepth, idxrate
436 real(DP) :: c, d, h, s, x
438 real(DP) :: petm1, petm2, pxdp1, pxdp2, thcof, trhs
441 if (this%nbound == 0)
return
444 do i = 1, this%nbound
447 if (this%fixed_cell)
then
448 node = this%nodelist(i)
450 node = this%nodesontop(i)
461 if (.not. this%fixed_cell)
then
462 if (this%ibound(node) == 0) &
463 call this%dis%highest_active(node, this%ibound)
464 this%nodelist(i) = node
472 if (this%ibound(node) > 0 .and. this%ibound(node) /=
iwetlake)
then
474 if (this%iauxmultcol > 0)
then
475 c = this%rate(i) * this%dis%get_area(node) * &
476 this%auxvar(this%iauxmultcol, i)
478 c = this%rate(i) * this%dis%get_area(node)
482 if (this%surfratespecified)
then
483 petm0 = this%petm0(i)
488 if (this%surfratespecified)
then
490 this%rhs(i) = this%rhs(i) + petm0 * c
493 this%rhs(i) = this%rhs(i) + c
501 if (this%nseg > 1)
then
514 if (this%surfratespecified)
then
525 segloop:
do iseg = 1, this%nseg
528 if (iseg < this%nseg)
then
530 idxdepth = idxdepth + 1
531 idxrate = idxrate + 1
533 pxdp2 = this%pxdp(idxdepth, i)
534 petm2 = this%petm(idxrate, i)
539 if (d <= pxdp2 * x)
then
550 thcof = -(petm1 - petm2) * c / ((pxdp2 - pxdp1) * x)
551 trhs = thcof * (s - pxdp1 * x) + petm1 * c
558 this%rhs(i) = this%rhs(i) + trhs
559 this%hcof(i) = this%hcof(i) + thcof
569 subroutine evt_fc(this, rhs, ia, idxglo, matrix_sln)
572 real(DP),
dimension(:),
intent(inout) :: rhs
573 integer(I4B),
dimension(:),
intent(in) :: ia
574 integer(I4B),
dimension(:),
intent(in) :: idxglo
577 integer(I4B) :: i, n, ipos
580 do i = 1, this%nbound
584 if (this%ibound(n) ==
iwetlake)
then
589 rhs(n) = rhs(n) + this%rhs(i)
591 call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
604 if (
associated(this%nodesontop))
deallocate (this%nodesontop)
609 if (.not. this%readasarrays)
then
610 if (this%nseg > 1)
then
615 if (this%surfratespecified)
then
622 deallocate (this%segsdefined)
623 deallocate (this%fixed_cell)
624 deallocate (this%surfratespecified)
627 call this%BndExtType%bnd_da()
635 class(
evttype),
intent(inout) :: this
637 integer(I4B) :: nsegm1, i
640 this%listlabel = trim(this%filtyp)//
' NO.'
641 if (this%dis%ndim == 3)
then
642 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'LAYER'
643 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'ROW'
644 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'COL'
645 elseif (this%dis%ndim == 2)
then
646 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'LAYER'
647 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'CELL2D'
649 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'NODE'
651 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'SURFACE'
652 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'MAX. RATE'
653 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'EXT. DEPTH'
656 nsegm1 = this%nseg - 1
659 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'PXDP'
662 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'PETM'
667 if (this%surfratespecified)
then
668 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'PETM0'
676 if (this%inamedbound == 1)
then
677 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'BOUNDARY NAME'
704 call this%obs%StoreObsType(
'evt', .true., indx)
714 class(
evttype),
intent(inout) :: this
715 integer(I4B),
intent(in) :: col
716 integer(I4B),
intent(in) :: row
727 bndval = this%surface(row)
729 if (this%iauxmultcol > 0)
then
730 bndval = this%rate(row) * this%auxvar(this%iauxmultcol, row)
732 bndval = this%rate(row)
735 bndval = this%depth(row)
738 if (this%nseg > 1)
then
739 if (col < (3 + this%nseg))
then
741 bndval = this%pxdp(idx, row)
742 else if (col < (3 + (2 * this%nseg) - 1))
then
743 idx = col - (3 + this%nseg - 1)
744 bndval = this%petm(idx, row)
745 else if (col == (3 + (2 * this%nseg) - 1))
then
746 if (this%surfratespecified)
then
748 bndval = this%petm0(row)
751 else if (this%surfratespecified)
then
754 bndval = this%petm0(row)
761 write (
errmsg,
'(a,i0,a)') &
762 'Programming error. EVT bound value requested column '&
763 &
'outside range of ncolbnd (', this%ncolbnd,
').'
This module contains the extended boundary package.
This module contains the base boundary package.
This module contains simulation constants.
integer(i4b), parameter lenpackagename
maximum length of the package name
integer(i4b), parameter iwetlake
integer constant for a dry lake
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
real(dp), parameter done
real constant 1
subroutine evt_source_dimensions(this)
Source the dimensions for this package.
subroutine evt_rp(this)
Read and Prepare.
subroutine set_nodesontop(this)
Store nodelist in nodesontop.
subroutine evt_cf(this)
Formulate the HCOF and RHS terms.
character(len=lenpackagename) text
subroutine evt_da(this)
Deallocate.
subroutine evt_df_obs(this)
Store observation type supported by EVT package.
subroutine evt_allocate_arrays(this, nodelist, auxvar)
Allocate package arrays.
real(dp) function evt_bound_value(this, col, row)
Return requested boundary value.
subroutine evt_log_options(this, found_fixed_cell, found_surfratespec)
Source options specific to EvtType.
subroutine evt_read_initial_attr(this)
Part of allocate and read.
subroutine evt_source_options(this)
Source options specific to EvtType.
subroutine evt_allocate_scalars(this)
Allocate package scalar members.
subroutine, public evt_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
Create a new Evapotranspiration Segments Package and point pakobj to the new package.
character(len=lenftype) ftype
subroutine evt_fc(this, rhs, ia, idxglo, matrix_sln)
Copy rhs and hcof into solution rhs and amat.
logical function evt_obs_supported(this)
Return true because EVT package supports observations.
character(len=lenpackagename) texta
subroutine evt_define_listlabel(this)
Define the list heading that is written to iout when PRINT_INPUT option is used.
subroutine check_pxdp(this)
Subroutine to check pxdp.
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
integer(i4b), pointer, public kper
current stress period number