MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
gwf-evt.f90
Go to the documentation of this file.
1 module evtmodule
2  !
3  use kindmodule, only: dp, i4b, lgp
7  use bndmodule, only: bndtype
8  use bndextmodule, only: bndexttype
10  use simvariablesmodule, only: errmsg
13  !
14  implicit none
15  !
16  private
17  public :: evt_create
18  !
19  character(len=LENFTYPE) :: ftype = 'EVT'
20  character(len=LENPACKAGENAME) :: text = ' EVT'
21  character(len=LENPACKAGENAME) :: texta = ' EVTA'
22  !
23  type, extends(bndexttype) :: evttype
24  ! -- logicals
25  logical, pointer, private :: segsdefined
26  logical, pointer, private :: fixed_cell
27  logical, pointer, private :: surfratespecified
28  ! -- integers
29  integer(I4B), pointer, private :: nseg => null() !< number of ET segments
30  ! -- arrays
31  real(dp), dimension(:), pointer, contiguous :: surface => null() !< elevation of the ET surface
32  real(dp), dimension(:), pointer, contiguous :: rate => null() !< maximum ET flux rate
33  real(dp), dimension(:), pointer, contiguous :: depth => null() !< ET extinction depth
34  real(dp), dimension(:, :), pointer, contiguous :: pxdp => null() !< proportion of ET extinction depth at bottom of segment
35  real(dp), dimension(:, :), pointer, contiguous :: petm => null() !< proportion of max ET flux rate at bottom of segment
36  real(dp), dimension(:), pointer, contiguous :: petm0 => null() !< proportion of max ET flux rate that will apply when head is at or above ET surface
37  integer(I4B), dimension(:), pointer, contiguous :: nodesontop => null()
38  contains
39  procedure :: evt_allocate_scalars
40  procedure :: allocate_arrays => evt_allocate_arrays
41  procedure :: source_options => evt_source_options
42  procedure :: source_dimensions => evt_source_dimensions
43  procedure :: evt_log_options
44  procedure :: read_initial_attr => evt_read_initial_attr
45  procedure :: bnd_rp => evt_rp
46  procedure :: set_nodesontop
47  procedure :: bnd_cf => evt_cf
48  procedure :: bnd_fc => evt_fc
49  procedure :: bnd_da => evt_da
50  procedure :: define_listlabel => evt_define_listlabel
51  procedure :: bound_value => evt_bound_value
52  procedure, private :: check_pxdp
53  ! -- for observations
54  procedure, public :: bnd_obs_supported => evt_obs_supported
55  procedure, public :: bnd_df_obs => evt_df_obs
56  end type evttype
57 
58 contains
59 
60  !> @brief Create a new Evapotranspiration Segments Package and point pakobj
61  !! to the new package
62  !<
63  subroutine evt_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
64  mempath)
65  ! -- dummy
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
74  ! -- local
75  type(evttype), pointer :: evtobj
76  !
77  ! -- allocate evt object and scalar variables
78  allocate (evtobj)
79  packobj => evtobj
80  !
81  ! -- create name and memory path
82  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
83  packobj%text = text
84  !
85  ! -- allocate scalars
86  call evtobj%evt_allocate_scalars()
87  !
88  ! -- initialize package
89  call packobj%pack_initialize()
90 
91  packobj%inunit = inunit
92  packobj%iout = iout
93  packobj%id = id
94  packobj%ibcnum = ibcnum
95  packobj%ncolbnd = 3
96  packobj%ictMemPath = create_mem_path(namemodel, 'NPF')
97  end subroutine evt_create
98 
99  !> @brief Allocate package scalar members
100  !<
101  subroutine evt_allocate_scalars(this)
102  ! -- modules
104  ! -- dummy
105  class(evttype), intent(inout) :: this
106  !
107  ! -- call standard BndType allocate scalars
108  call this%BndExtType%allocate_scalars()
109  !
110  ! -- allocate the object and assign values to object variables
111  call mem_allocate(this%nseg, 'NSEG', this%memoryPath)
112  !
113  ! -- allocate internal members
114  allocate (this%segsdefined)
115  allocate (this%fixed_cell)
116  allocate (this%surfratespecified)
117  !
118  ! -- Set values
119  this%nseg = 1
120  this%segsdefined = .true.
121  this%fixed_cell = .false.
122  this%surfratespecified = .false.
123  end subroutine evt_allocate_scalars
124 
125  !> @brief Allocate package arrays
126  !<
127  subroutine evt_allocate_arrays(this, nodelist, auxvar)
128  ! -- modules
130  ! -- dummy
131  class(evttype) :: this
132  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
133  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
134  !
135  ! -- call standard BndType allocate scalars
136  call this%BndExtType%allocate_arrays(nodelist, auxvar)
137  !
138  ! -- set EVT input context pointers
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)
142  !
143  ! -- checkin EVT input context pointers
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)
150  !
151  ! -- set list input segment descriptors
152  if (.not. this%readasarrays) then
153  if (this%nseg > 1) then
154  !
155  ! -- set pxdp and petm input context pointers
156  call mem_setptr(this%pxdp, 'PXDP', this%input_mempath)
157  call mem_setptr(this%petm, 'PETM', this%input_mempath)
158  !
159  ! -- checkin pxdp and petm input context pointers
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)
164  end if
165  !
166  if (this%surfratespecified) then
167  !
168  ! -- set petm0 input context pointer
169  call mem_setptr(this%petm0, 'PETM0', this%input_mempath)
170  !
171  ! -- cehckin petm0 input context pointer
172  call mem_checkin(this%petm0, 'PETM0', this%memoryPath, &
173  'PETM0', this%input_mempath)
174  end if
175  end if
176  end subroutine evt_allocate_arrays
177 
178  !> @brief Source options specific to EvtType
179  !<
180  subroutine evt_source_options(this)
181  ! -- modules
183  ! -- dummy
184  class(evttype), intent(inout) :: this
185  ! -- local
186  logical(LGP) :: found_fixed_cell = .false.
187  logical(LGP) :: found_surfratespec = .false.
188  !
189  ! -- source common bound options
190  call this%BndExtType%source_options()
191  !
192  ! -- update defaults with idm sourced values
193  call mem_set_value(this%fixed_cell, 'FIXED_CELL', &
194  this%input_mempath, found_fixed_cell)
195  call mem_set_value(this%surfratespecified, 'SURFRATESPEC', &
196  this%input_mempath, found_surfratespec)
197  !
198  if (this%readasarrays) then
199  if (found_surfratespec) then
200  errmsg = 'READASARRAYS option is not compatible with the'// &
201  ' SURF_RATE_SPECIFIED option.'
202  call store_error(errmsg)
203  call store_error_filename(this%input_fname)
204  end if
205  !
206  this%text = texta
207  end if
208  !
209  ! -- log evt specific options
210  call this%evt_log_options(found_fixed_cell, found_surfratespec)
211  end subroutine evt_source_options
212 
213  !> @brief Source options specific to EvtType
214  !<
215  subroutine evt_log_options(this, found_fixed_cell, found_surfratespec)
216  ! -- modules
219  ! -- dummy
220  class(evttype), intent(inout) :: this
221  logical(LGP), intent(in) :: found_fixed_cell
222  logical(LGP), intent(in) :: found_surfratespec
223  ! -- formats
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.')"
232  !
233  ! -- log found options
234  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
235  //' OPTIONS'
236  !
237  if (found_fixed_cell) then
238  write (this%iout, fmtfixedcell)
239  end if
240  !
241  if (found_surfratespec) then
242  write (this%iout, fmtsrs)
243  end if
244  !
245  ! -- close logging block
246  write (this%iout, '(1x,a)') &
247  'END OF '//trim(adjustl(this%text))//' OPTIONS'
248  end subroutine evt_log_options
249 
250  !> @brief Source the dimensions for this package
251  !<
252  subroutine evt_source_dimensions(this)
253  ! -- modules
255  ! -- dummy
256  class(evttype), intent(inout) :: this
257  ! -- local
258  logical(LGP) :: found_nseg = .false.
259  ! -- format
260  character(len=*), parameter :: fmtnsegerr = &
261  &"('Error: In EVT, NSEG must be > 0 but is specified as ',i0)"
262 
263  ! -- set maxbound
264  call this%BndExtType%source_dimensions()
265 
266  if (this%readasarrays) then
267  ! -- base class call is sufficient
268  else
269  !
270  ! -- log found options
271  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
272  //' DIMENSIONS'
273  !
274  ! -- update defaults with idm sourced values
275  call mem_set_value(this%nseg, 'NSEG', this%input_mempath, found_nseg)
276  !
277  if (found_nseg) then
278  !
279  write (this%iout, '(4x,a,i0)') 'NSEG = ', this%nseg
280  !
281  if (this%nseg < 1) then
282  write (errmsg, fmtnsegerr) this%nseg
283  call store_error(errmsg)
284  call store_error_filename(this%input_fname)
285  !
286  elseif (this%nseg > 1) then
287  ! NSEG>1 is supported only if readasarrays is false
288  if (this%readasarrays) then
289  errmsg = 'In the EVT package, NSEG cannot be greater than 1'// &
290  ' when READASARRAYS is used.'
291  call store_error(errmsg)
292  call store_error_filename(this%input_fname)
293  end if
294  !
295  end if
296  end if
297  !
298  ! -- close logging block
299  write (this%iout, '(1x,a)') &
300  'END OF '//trim(adjustl(this%text))//' DIMENSIONS'
301  !
302  end if
303  end subroutine evt_source_dimensions
304 
305  !> @brief Part of allocate and read
306  !!
307  !! If READASARRAYS has been specified, assign default IEVT = 1
308  !<
309  subroutine evt_read_initial_attr(this)
310  ! -- dummy
311  class(evttype), intent(inout) :: this
312  !
313  if (this%readasarrays) then
314  call this%default_nodelist()
315  end if
316  !
317  ! -- if fixed_cell option not set, then need to store nodelist
318  ! in the nodesontop array
319  if (.not. this%fixed_cell) call this%set_nodesontop()
320  end subroutine evt_read_initial_attr
321 
322  !> @brief Read and Prepare
323  !!
324  !! Read itmp and new boundaries if itmp > 0
325  !<
326  subroutine evt_rp(this)
327  use tdismodule, only: kper
328  implicit none
329  ! -- dummy
330  class(evttype), intent(inout) :: this
331  !
332  if (this%iper /= kper) return
333  !
334  ! -- process the input list arrays
335  call this%BndExtType%bnd_rp()
336  !
337  ! -- ensure pxdp is monotonically increasing
338  if (this%nseg > 1) then
339  call this%check_pxdp()
340  end if
341  !
342  ! -- copy nodelist to nodesontop if not fixed cell
343  if (.not. this%fixed_cell) call this%set_nodesontop()
344  end subroutine evt_rp
345 
346  !> @brief Subroutine to check pxdp
347  !!
348  !! If the number of EVT segments (nseg) is greater than one, then
349  !! pxdp must be monotically increasing from zero to one. Check
350  !! to make sure this is the case.
351  !<
352  subroutine check_pxdp(this)
353  ! -- dummy
354  class(evttype), intent(inout) :: this !< EvtType
355  ! -- local
356  integer(I4B) :: n
357  integer(I4B) :: node
358  integer(I4B) :: i
359  integer(I4B) :: ierrmono
360  real(DP) :: pxdp1
361  real(DP) :: pxdp2
362  character(len=15) :: nodestr
363  ! -- formats
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)"
368  !
369  ! -- check and make sure that pxdp is monotonically increasing and
370  ! that pxdp values are between 0 and 1
371  do n = 1, this%nbound
372  node = this%nodelist(n)
373  pxdp1 = dzero
374  ierrmono = 0
375  segloop: do i = 1, this%nseg
376  !
377  ! -- set and check pxdp2
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)
383  call store_error(errmsg)
384  end if
385  else
386  pxdp2 = done
387  end if
388  !
389  ! -- check for monotonically increasing condition
390  if (pxdp2 - pxdp1 < dzero) then
391  if (ierrmono == 0) then
392  ! -- only store mono error once for each node
393  call this%dis%noder_to_string(node, nodestr)
394  write (errmsg, fmtpxdp) trim(nodestr)
395  call store_error(errmsg)
396  end if
397  ierrmono = 1
398  end if
399  pxdp1 = pxdp2
400  end do segloop
401  end do
402  !
403  ! -- terminate if errors encountered
404  if (count_errors() > 0) then
405  call store_error_filename(this%input_fname)
406  end if
407  end subroutine check_pxdp
408 
409  !> @brief Store nodelist in nodesontop
410  !<
411  subroutine set_nodesontop(this)
412  ! -- dummy
413  class(evttype), intent(inout) :: this
414  ! -- local
415  integer(I4B) :: n
416  !
417  ! -- allocate if necessary
418  if (.not. associated(this%nodesontop)) then
419  allocate (this%nodesontop(this%maxbound))
420  end if
421  !
422  ! -- copy nodelist into nodesontop
423  do n = 1, this%nbound
424  this%nodesontop(n) = this%nodelist(n)
425  end do
426  end subroutine set_nodesontop
427 
428  !> @brief Formulate the HCOF and RHS terms
429  !<
430  subroutine evt_cf(this)
431  ! -- dummy
432  class(evttype) :: this
433  ! -- local
434  integer(I4B) :: i, iseg, node
435  integer(I4B) :: idxdepth, idxrate
436  real(DP) :: c, d, h, s, x
437  real(DP) :: petm0
438  real(DP) :: petm1, petm2, pxdp1, pxdp2, thcof, trhs
439  !
440  ! -- Return if no ET nodes
441  if (this%nbound == 0) return
442  !
443  ! -- Calculate hcof and rhs for each ET node
444  do i = 1, this%nbound
445  !
446  ! -- Find the node number
447  if (this%fixed_cell) then
448  node = this%nodelist(i)
449  else
450  node = this%nodesontop(i)
451  end if
452  !
453  ! -- cycle if nonexistent bound
454  if (node <= 0) then
455  this%hcof(i) = dzero
456  this%rhs(i) = dzero
457  cycle
458  end if
459  !
460  ! -- reset nodelist to highest active
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
465  end if
466  !
467  ! -- set rhs and hcof to zero
468  this%rhs(i) = dzero
469  this%hcof(i) = dzero
470  !
471  ! -- if ibound is positive and not overlain by a lake, then add terms
472  if (this%ibound(node) > 0 .and. this%ibound(node) /= iwetlake) then
473  !
474  if (this%iauxmultcol > 0) then
475  c = this%rate(i) * this%dis%get_area(node) * &
476  this%auxvar(this%iauxmultcol, i)
477  else
478  c = this%rate(i) * this%dis%get_area(node)
479  end if
480  s = this%surface(i)
481  h = this%xnew(node)
482  if (this%surfratespecified) then
483  petm0 = this%petm0(i)
484  end if
485  !
486  ! -- If head in cell is greater than or equal to SURFACE, ET is constant
487  if (h >= s) then
488  if (this%surfratespecified) then
489  ! -- Subtract -PETM0 * max rate from RHS
490  this%rhs(i) = this%rhs(i) + petm0 * c
491  else
492  ! -- Subtract -RATE from RHS
493  this%rhs(i) = this%rhs(i) + c
494  end if
495  else
496  ! -- If depth to water >= extinction depth, then ET is 0
497  d = s - h
498  x = this%depth(i)
499  if (d < x) then
500  ! -- Variable range. add ET terms to both RHS and HCOF.
501  if (this%nseg > 1) then
502  ! -- Determine which segment applies based on head, and
503  ! calculate terms to add to RHS and HCOF
504  !
505  ! -- Set proportions corresponding to surface elevation
506  ! and initial indices for depth and rate.
507  ! -- Idxdepth will point to the elements of bound containing
508  ! proportion of depth at the bottom of each segment.
509  ! Idxrate will point to the elements of bound containing
510  ! proportion of ET rate at the bottom of each segment.
511  ! If surfratespecified is true, skip over the elements
512  ! containing pxdp0 (=0.0) and petm0.
513  pxdp1 = dzero
514  if (this%surfratespecified) then
515  petm1 = petm0
516  else
517  petm1 = done
518  end if
519  ! -- Initialize indices to point to elements preceding
520  ! pxdp1 and petm1 (values for lower end of segment 1).
521  idxdepth = 0
522  idxrate = 0
523  ! -- Iterate through segments to find segment that contains
524  ! current depth of head below ET surface.
525  segloop: do iseg = 1, this%nseg
526  ! -- Set proportions corresponding to lower end of
527  ! segment
528  if (iseg < this%nseg) then
529  ! -- Increment the indices for depth and rate
530  idxdepth = idxdepth + 1
531  idxrate = idxrate + 1
532  ! -- Get proportions for lower end of segment
533  pxdp2 = this%pxdp(idxdepth, i)
534  petm2 = this%petm(idxrate, i)
535  else
536  pxdp2 = done
537  petm2 = dzero
538  end if
539  if (d <= pxdp2 * x) then
540  ! -- head is in domain of this segment
541  exit segloop
542  end if
543  ! -- Proportions at lower end of segment will be for
544  ! upper end of segment next time through loop
545  pxdp1 = pxdp2
546  petm1 = petm2
547  end do segloop
548  ! -- Calculate terms to add to RHS and HCOF based on
549  ! segment that applies at head elevation
550  thcof = -(petm1 - petm2) * c / ((pxdp2 - pxdp1) * x)
551  trhs = thcof * (s - pxdp1 * x) + petm1 * c
552  else
553  ! -- Calculate terms to add to RHS and HCOF based on simple
554  ! linear relation of ET vs. head for single segment
555  trhs = c - c * s / x
556  thcof = -c / x
557  end if
558  this%rhs(i) = this%rhs(i) + trhs
559  this%hcof(i) = this%hcof(i) + thcof
560  end if
561  end if
562  end if
563  !
564  end do
565  end subroutine evt_cf
566 
567  !> @brief Copy rhs and hcof into solution rhs and amat
568  !<
569  subroutine evt_fc(this, rhs, ia, idxglo, matrix_sln)
570  ! -- dummy
571  class(evttype) :: this
572  real(DP), dimension(:), intent(inout) :: rhs
573  integer(I4B), dimension(:), intent(in) :: ia
574  integer(I4B), dimension(:), intent(in) :: idxglo
575  class(matrixbasetype), pointer :: matrix_sln
576  ! -- local
577  integer(I4B) :: i, n, ipos
578  !
579  ! -- Copy package rhs and hcof into solution rhs and amat
580  do i = 1, this%nbound
581  n = this%nodelist(i)
582  if (n <= 0) cycle
583  ! -- reset hcof and rhs for excluded cells
584  if (this%ibound(n) == iwetlake) then
585  this%hcof(i) = dzero
586  this%rhs(i) = dzero
587  cycle
588  end if
589  rhs(n) = rhs(n) + this%rhs(i)
590  ipos = ia(n)
591  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
592  end do
593  end subroutine evt_fc
594 
595  !> @brief Deallocate
596  !<
597  subroutine evt_da(this)
598  ! -- modules
600  ! -- dummy
601  class(evttype) :: this
602  !
603  ! -- arrays
604  if (associated(this%nodesontop)) deallocate (this%nodesontop)
605  call mem_deallocate(this%surface, 'SURFACE', this%memoryPath)
606  call mem_deallocate(this%rate, 'RATE', this%memoryPath)
607  call mem_deallocate(this%depth, 'DEPTH', this%memoryPath)
608  !
609  if (.not. this%readasarrays) then
610  if (this%nseg > 1) then
611  call mem_deallocate(this%pxdp, 'PXDP', this%memoryPath)
612  call mem_deallocate(this%petm, 'PETM', this%memoryPath)
613  end if
614  !
615  if (this%surfratespecified) then
616  call mem_deallocate(this%petm0, 'PETM0', this%memoryPath)
617  end if
618  end if
619  !
620  ! -- scalars
621  call mem_deallocate(this%nseg)
622  deallocate (this%segsdefined)
623  deallocate (this%fixed_cell)
624  deallocate (this%surfratespecified)
625  !
626  ! -- Deallocate parent package
627  call this%BndExtType%bnd_da()
628  end subroutine evt_da
629 
630  !> @brief Define the list heading that is written to iout when PRINT_INPUT
631  !! option is used
632  !<
633  subroutine evt_define_listlabel(this)
634  ! -- dummy
635  class(evttype), intent(inout) :: this
636  ! -- local
637  integer(I4B) :: nsegm1, i
638  !
639  ! -- create the header list label
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'
648  else
649  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
650  end if
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'
654  !
655  ! -- add headings for as many PXDP and PETM columns as needed
656  nsegm1 = this%nseg - 1
657  if (nsegm1 > 0) then
658  do i = 1, nsegm1
659  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'PXDP'
660  end do
661  do i = 1, nsegm1
662  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'PETM'
663  end do
664  end if
665  !
666  ! -- PETM0, if SURF_RATE_SPECIFIED is used
667  if (this%surfratespecified) then
668  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'PETM0'
669  end if
670  !
671 ! ! -- multiplier
672 ! if(this%multindex > 0) &
673 ! write(this%listlabel, '(a, a16)') trim(this%listlabel), 'MULTIPLIER'
674  !
675  ! -- boundary name
676  if (this%inamedbound == 1) then
677  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
678  end if
679  end subroutine evt_define_listlabel
680 
681  ! -- Procedures related to observations
682 
683  !> @brief Return true because EVT package supports observations
684  !!
685  !! Overrides BndType%bnd_obs_supported()
686  !<
687  logical function evt_obs_supported(this)
688  ! -- dummy
689  class(evttype) :: this
690  !
691  evt_obs_supported = .true.
692  end function evt_obs_supported
693 
694  !> @brief Store observation type supported by EVT package
695  !!
696  !! Overrides BndType%bnd_df_obs
697  !<
698  subroutine evt_df_obs(this)
699  ! -- dummy
700  class(evttype) :: this
701  ! -- local
702  integer(I4B) :: indx
703  !
704  call this%obs%StoreObsType('evt', .true., indx)
705  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
706  end subroutine evt_df_obs
707 
708  !> @brief Return requested boundary value
709  !<
710  function evt_bound_value(this, col, row) result(bndval)
711  ! -- modules
712  use constantsmodule, only: dzero
713  ! -- dummy variables
714  class(evttype), intent(inout) :: this
715  integer(I4B), intent(in) :: col
716  integer(I4B), intent(in) :: row
717  ! -- result
718  real(dp) :: bndval
719  ! -- local
720  integer(I4B) :: idx
721  !
722  ! -- initialize
723  idx = 0
724  !
725  select case (col)
726  case (1)
727  bndval = this%surface(row)
728  case (2)
729  if (this%iauxmultcol > 0) then
730  bndval = this%rate(row) * this%auxvar(this%iauxmultcol, row)
731  else
732  bndval = this%rate(row)
733  end if
734  case (3)
735  bndval = this%depth(row)
736  case default
737  if (col > 0) then
738  if (this%nseg > 1) then
739  if (col < (3 + this%nseg)) then
740  idx = col - 3
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
747  idx = 1
748  bndval = this%petm0(row)
749  end if
750  end if
751  else if (this%surfratespecified) then
752  if (col == 4) then
753  idx = 1
754  bndval = this%petm0(row)
755  end if
756  end if
757  end if
758  !
759  ! -- set error if idx not found
760  if (idx == 0) then
761  write (errmsg, '(a,i0,a)') &
762  'Programming error. EVT bound value requested column '&
763  &'outside range of ncolbnd (', this%ncolbnd, ').'
764  call store_error(errmsg)
765  call store_error_filename(this%input_fname)
766  end if
767  !
768  end select
769  end function evt_bound_value
770 
771 end module evtmodule
772 
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 lenpackagename
maximum length of the package name
Definition: Constants.f90:23
integer(i4b), parameter iwetlake
integer constant for a dry lake
Definition: Constants.f90:52
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
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
subroutine evt_source_dimensions(this)
Source the dimensions for this package.
Definition: gwf-evt.f90:253
subroutine evt_rp(this)
Read and Prepare.
Definition: gwf-evt.f90:327
subroutine set_nodesontop(this)
Store nodelist in nodesontop.
Definition: gwf-evt.f90:412
subroutine evt_cf(this)
Formulate the HCOF and RHS terms.
Definition: gwf-evt.f90:431
character(len=lenpackagename) text
Definition: gwf-evt.f90:20
subroutine evt_da(this)
Deallocate.
Definition: gwf-evt.f90:598
subroutine evt_df_obs(this)
Store observation type supported by EVT package.
Definition: gwf-evt.f90:699
subroutine evt_allocate_arrays(this, nodelist, auxvar)
Allocate package arrays.
Definition: gwf-evt.f90:128
real(dp) function evt_bound_value(this, col, row)
Return requested boundary value.
Definition: gwf-evt.f90:711
subroutine evt_log_options(this, found_fixed_cell, found_surfratespec)
Source options specific to EvtType.
Definition: gwf-evt.f90:216
subroutine evt_read_initial_attr(this)
Part of allocate and read.
Definition: gwf-evt.f90:310
subroutine evt_source_options(this)
Source options specific to EvtType.
Definition: gwf-evt.f90:181
subroutine evt_allocate_scalars(this)
Allocate package scalar members.
Definition: gwf-evt.f90:102
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.
Definition: gwf-evt.f90:65
character(len=lenftype) ftype
Definition: gwf-evt.f90:19
subroutine evt_fc(this, rhs, ia, idxglo, matrix_sln)
Copy rhs and hcof into solution rhs and amat.
Definition: gwf-evt.f90:570
logical function evt_obs_supported(this)
Return true because EVT package supports observations.
Definition: gwf-evt.f90:688
character(len=lenpackagename) texta
Definition: gwf-evt.f90:21
subroutine evt_define_listlabel(this)
Define the list heading that is written to iout when PRINT_INPUT option is used.
Definition: gwf-evt.f90:634
subroutine check_pxdp(this)
Subroutine to check pxdp.
Definition: gwf-evt.f90:353
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
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
@ brief BndType