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