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
104 class(
evttype),
intent(inout) :: this
107 call this%BndExtType%allocate_scalars()
113 allocate (this%segsdefined)
114 allocate (this%fixed_cell)
115 allocate (this%surfratespecified)
119 this%segsdefined = .true.
120 this%fixed_cell = .false.
121 this%surfratespecified = .false.
131 integer(I4B),
dimension(:),
pointer,
contiguous,
optional :: nodelist
132 real(DP),
dimension(:, :),
pointer,
contiguous,
optional :: auxvar
135 call this%BndExtType%allocate_arrays(nodelist, auxvar)
138 call mem_setptr(this%surface,
'SURFACE', this%input_mempath)
139 call mem_setptr(this%rate,
'RATE', this%input_mempath)
140 call mem_setptr(this%depth,
'DEPTH', this%input_mempath)
143 call mem_checkin(this%surface,
'SURFACE', this%memoryPath, &
144 'SURFACE', this%input_mempath)
145 call mem_checkin(this%rate,
'RATE', this%memoryPath, &
146 'RATE', this%input_mempath)
147 call mem_checkin(this%depth,
'DEPTH', this%memoryPath, &
148 'DEPTH', this%input_mempath)
151 if (.not. this%readasarrays)
then
152 if (this%nseg > 1)
then
155 call mem_setptr(this%pxdp,
'PXDP', this%input_mempath)
156 call mem_setptr(this%petm,
'PETM', this%input_mempath)
159 call mem_checkin(this%pxdp,
'PXDP', this%memoryPath, &
160 'PXDP', this%input_mempath)
161 call mem_checkin(this%petm,
'PETM', this%memoryPath, &
162 'PETM', this%input_mempath)
165 if (this%surfratespecified)
then
168 call mem_setptr(this%petm0,
'PETM0', this%input_mempath)
171 call mem_checkin(this%petm0,
'PETM0', this%memoryPath, &
172 'PETM0', this%input_mempath)
183 class(
evttype),
intent(inout) :: this
185 logical(LGP) :: found_fixed_cell = .false.
186 logical(LGP) :: found_surfratespec = .false.
189 call this%BndExtType%source_options()
193 this%input_mempath, found_fixed_cell)
195 this%input_mempath, found_surfratespec)
197 if (this%readasarrays)
then
198 if (found_surfratespec)
then
199 errmsg =
'READASARRAYS option is not compatible with the'// &
200 ' SURF_RATE_SPECIFIED option.'
209 call this%evt_log_options(found_fixed_cell, found_surfratespec)
219 class(
evttype),
intent(inout) :: this
220 logical(LGP),
intent(in) :: found_fixed_cell
221 logical(LGP),
intent(in) :: found_surfratespec
223 character(len=*),
parameter :: fmtihact = &
224 &
"(4x, 'EVAPOTRANSPIRATION WILL BE APPLIED TO HIGHEST ACTIVE CELL.')"
225 character(len=*),
parameter :: fmtfixedcell = &
226 &
"(4x, 'EVAPOTRANSPIRATION WILL BE APPLIED TO SPECIFIED CELL.')"
227 character(len=*),
parameter :: fmtsrz = &
228 &
"(4x, 'ET RATE AT SURFACE WILL BE ZERO.')"
229 character(len=*),
parameter :: fmtsrs = &
230 &
"(4x, 'ET RATE AT SURFACE WILL BE AS SPECIFIED BY PETM0.')"
233 write (this%iout,
'(/1x,a)')
'PROCESSING '//trim(adjustl(this%text)) &
236 if (found_fixed_cell)
then
237 write (this%iout, fmtfixedcell)
240 if (found_surfratespec)
then
241 write (this%iout, fmtsrs)
245 write (this%iout,
'(1x,a)') &
246 'END OF '//trim(adjustl(this%text))//
' OPTIONS'
255 class(
evttype),
intent(inout) :: this
257 logical(LGP) :: found_nseg = .false.
259 character(len=*),
parameter :: fmtnsegerr = &
260 &
"('Error: In EVT, NSEG must be > 0 but is specified as ',i0)"
263 call this%BndExtType%source_dimensions()
265 if (this%readasarrays)
then
270 write (this%iout,
'(/1x,a)')
'PROCESSING '//trim(adjustl(this%text)) &
274 call mem_set_value(this%nseg,
'NSEG', this%input_mempath, found_nseg)
278 write (this%iout,
'(4x,a,i0)')
'NSEG = ', this%nseg
280 if (this%nseg < 1)
then
281 write (
errmsg, fmtnsegerr) this%nseg
285 elseif (this%nseg > 1)
then
287 if (this%readasarrays)
then
288 errmsg =
'In the EVT package, NSEG cannot be greater than 1'// &
289 ' when READASARRAYS is used.'
298 write (this%iout,
'(1x,a)') &
299 'END OF '//trim(adjustl(this%text))//
' DIMENSIONS'
310 class(
evttype),
intent(inout) :: this
312 if (this%readasarrays)
then
313 call this%default_nodelist()
318 if (.not. this%fixed_cell)
call this%set_nodesontop()
329 class(
evttype),
intent(inout) :: this
331 if (this%iper /=
kper)
return
334 call this%BndExtType%bnd_rp()
337 if (this%nseg > 1)
then
338 call this%check_pxdp()
341 if (this%iprpak /= 0)
then
342 if (this%readasarrays)
then
346 call this%write_list()
351 if (.not. this%fixed_cell)
call this%set_nodesontop()
362 class(
evttype),
intent(inout) :: this
367 integer(I4B) :: ierrmono
370 character(len=15) :: nodestr
372 character(len=*),
parameter :: fmtpxdp0 = &
373 &
"('PXDP must be between 0 and 1. Found ', G0, ' for cell ', A)"
374 character(len=*),
parameter :: fmtpxdp = &
375 &
"('PXDP is not monotonically increasing for cell ', A)"
379 do n = 1, this%nbound
380 node = this%nodelist(n)
383 segloop:
do i = 1, this%nseg
386 if (i < this%nseg)
then
387 pxdp2 = this%pxdp(i, n)
388 if (pxdp2 <=
dzero .or. pxdp2 >=
done)
then
389 call this%dis%noder_to_string(node, nodestr)
390 write (
errmsg, fmtpxdp0) pxdp2, trim(nodestr)
398 if (pxdp2 - pxdp1 <
dzero)
then
399 if (ierrmono == 0)
then
401 call this%dis%noder_to_string(node, nodestr)
402 write (
errmsg, fmtpxdp) trim(nodestr)
421 class(
evttype),
intent(inout) :: this
426 if (.not.
associated(this%nodesontop))
then
427 allocate (this%nodesontop(this%maxbound))
431 do n = 1, this%nbound
432 this%nodesontop(n) = this%nodelist(n)
442 integer(I4B) :: i, iseg, node
443 integer(I4B) :: idxdepth, idxrate
444 real(DP) :: c, d, h, s, x
446 real(DP) :: petm1, petm2, pxdp1, pxdp2, thcof, trhs
449 if (this%nbound == 0)
return
452 do i = 1, this%nbound
455 if (this%fixed_cell)
then
456 node = this%nodelist(i)
458 node = this%nodesontop(i)
469 if (.not. this%fixed_cell)
then
470 if (this%ibound(node) == 0) &
471 call this%dis%highest_active(node, this%ibound)
472 this%nodelist(i) = node
480 if (this%ibound(node) > 0 .and. this%ibound(node) /=
iwetlake)
then
482 if (this%iauxmultcol > 0)
then
483 c = this%rate(i) * this%dis%get_area(node) * &
484 this%auxvar(this%iauxmultcol, i)
486 c = this%rate(i) * this%dis%get_area(node)
490 if (this%surfratespecified)
then
491 petm0 = this%petm0(i)
496 if (this%surfratespecified)
then
498 this%rhs(i) = this%rhs(i) + petm0 * c
501 this%rhs(i) = this%rhs(i) + c
509 if (this%nseg > 1)
then
522 if (this%surfratespecified)
then
533 segloop:
do iseg = 1, this%nseg
536 if (iseg < this%nseg)
then
538 idxdepth = idxdepth + 1
539 idxrate = idxrate + 1
541 pxdp2 = this%pxdp(idxdepth, i)
542 petm2 = this%petm(idxrate, i)
547 if (d <= pxdp2 * x)
then
558 thcof = -(petm1 - petm2) * c / ((pxdp2 - pxdp1) * x)
559 trhs = thcof * (s - pxdp1 * x) + petm1 * c
566 this%rhs(i) = this%rhs(i) + trhs
567 this%hcof(i) = this%hcof(i) + thcof
577 subroutine evt_fc(this, rhs, ia, idxglo, matrix_sln)
580 real(DP),
dimension(:),
intent(inout) :: rhs
581 integer(I4B),
dimension(:),
intent(in) :: ia
582 integer(I4B),
dimension(:),
intent(in) :: idxglo
585 integer(I4B) :: i, n, ipos
588 do i = 1, this%nbound
592 if (this%ibound(n) ==
iwetlake)
then
597 rhs(n) = rhs(n) + this%rhs(i)
599 call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
612 if (
associated(this%nodesontop))
deallocate (this%nodesontop)
617 if (.not. this%readasarrays)
then
618 if (this%nseg > 1)
then
623 if (this%surfratespecified)
then
630 deallocate (this%segsdefined)
631 deallocate (this%fixed_cell)
632 deallocate (this%surfratespecified)
635 call this%BndExtType%bnd_da()
643 class(
evttype),
intent(inout) :: this
645 integer(I4B) :: nsegm1, i
648 this%listlabel = trim(this%filtyp)//
' NO.'
649 if (this%dis%ndim == 3)
then
650 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'LAYER'
651 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'ROW'
652 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'COL'
653 elseif (this%dis%ndim == 2)
then
654 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'LAYER'
655 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'CELL2D'
657 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'NODE'
659 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'SURFACE'
660 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'MAX. RATE'
661 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'EXT. DEPTH'
664 nsegm1 = this%nseg - 1
667 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'PXDP'
670 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'PETM'
675 if (this%surfratespecified)
then
676 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'PETM0'
684 if (this%inamedbound == 1)
then
685 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'BOUNDARY NAME'
712 call this%obs%StoreObsType(
'evt', .true., indx)
722 class(
evttype),
intent(inout) :: this
723 integer(I4B),
intent(in) :: col
724 integer(I4B),
intent(in) :: row
735 bndval = this%surface(row)
737 if (this%iauxmultcol > 0)
then
738 bndval = this%rate(row) * this%auxvar(this%iauxmultcol, row)
740 bndval = this%rate(row)
743 bndval = this%depth(row)
746 if (this%nseg > 1)
then
747 if (col < (3 + this%nseg))
then
749 bndval = this%pxdp(idx, row)
750 else if (col < (3 + (2 * this%nseg) - 1))
then
751 idx = col - (3 + this%nseg - 1)
752 bndval = this%petm(idx, row)
753 else if (col == (3 + (2 * this%nseg) - 1))
then
754 if (this%surfratespecified)
then
756 bndval = this%petm0(row)
759 else if (this%surfratespecified)
then
762 bndval = this%petm0(row)
769 write (
errmsg,
'(a,i0,a)') &
770 'Programming error. EVT bound value requested column '&
771 &
'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