MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
evtmodule Module Reference

Data Types

type  evttype
 

Functions/Subroutines

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. More...
 
subroutine evt_allocate_scalars (this)
 Allocate package scalar members. More...
 
subroutine evt_allocate_arrays (this, nodelist, auxvar)
 Allocate package arrays. More...
 
subroutine evt_source_options (this)
 Source options specific to EvtType. More...
 
subroutine evt_log_options (this, found_fixed_cell, found_readasarrays, found_surfratespec)
 Source options specific to EvtType. More...
 
subroutine evt_source_dimensions (this)
 Source the dimensions for this package. More...
 
subroutine evt_read_initial_attr (this)
 Part of allocate and read. More...
 
subroutine evt_rp (this)
 Read and Prepare. More...
 
subroutine check_pxdp (this)
 Subroutine to check pxdp. More...
 
subroutine set_nodesontop (this)
 Store nodelist in nodesontop. More...
 
subroutine evt_cf (this)
 Formulate the HCOF and RHS terms. More...
 
subroutine evt_fc (this, rhs, ia, idxglo, matrix_sln)
 Copy rhs and hcof into solution rhs and amat. More...
 
subroutine evt_da (this)
 Deallocate. More...
 
subroutine evt_define_listlabel (this)
 Define the list heading that is written to iout when PRINT_INPUT option is used. More...
 
subroutine default_nodelist (this)
 Assign default nodelist when READASARRAYS is specified. More...
 
logical function evt_obs_supported (this)
 Return true because EVT package supports observations. More...
 
subroutine evt_df_obs (this)
 Store observation type supported by EVT package. More...
 
real(dp) function evt_bound_value (this, col, row)
 Return requested boundary value. More...
 
subroutine nodelist_update (nodelist, nbound, maxbound, dis, input_mempath)
 Update the nodelist based on IEVT input. More...
 

Variables

character(len=lenftype) ftype = 'EVT'
 
character(len=lenpackagename) text = ' EVT'
 
character(len=lenpackagename) texta = ' EVTA'
 

Function/Subroutine Documentation

◆ check_pxdp()

subroutine evtmodule::check_pxdp ( class(evttype), intent(inout)  this)

If the number of EVT segments (nseg) is greater than one, then pxdp must be monotically increasing from zero to one. Check to make sure this is the case.

Parameters
[in,out]thisEvtType

Definition at line 406 of file gwf-evt.f90.

407  ! -- dummy
408  class(EvtType), intent(inout) :: this !< EvtType
409  ! -- local
410  integer(I4B) :: n
411  integer(I4B) :: node
412  integer(I4B) :: i
413  integer(I4B) :: ierrmono
414  real(DP) :: pxdp1
415  real(DP) :: pxdp2
416  character(len=15) :: nodestr
417  ! -- formats
418  character(len=*), parameter :: fmtpxdp0 = &
419  &"('PXDP must be between 0 and 1. Found ', G0, ' for cell ', A)"
420  character(len=*), parameter :: fmtpxdp = &
421  &"('PXDP is not monotonically increasing for cell ', A)"
422  !
423  ! -- check and make sure that pxdp is monotonically increasing and
424  ! that pxdp values are between 0 and 1
425  do n = 1, this%nbound
426  node = this%nodelist(n)
427  pxdp1 = dzero
428  ierrmono = 0
429  segloop: do i = 1, this%nseg
430  !
431  ! -- set and check pxdp2
432  if (i < this%nseg) then
433  pxdp2 = this%pxdp(i, n)
434  if (pxdp2 <= dzero .or. pxdp2 >= done) then
435  call this%dis%noder_to_string(node, nodestr)
436  write (errmsg, fmtpxdp0) pxdp2, trim(nodestr)
437  call store_error(errmsg)
438  end if
439  else
440  pxdp2 = done
441  end if
442  !
443  ! -- check for monotonically increasing condition
444  if (pxdp2 - pxdp1 < dzero) then
445  if (ierrmono == 0) then
446  ! -- only store mono error once for each node
447  call this%dis%noder_to_string(node, nodestr)
448  write (errmsg, fmtpxdp) trim(nodestr)
449  call store_error(errmsg)
450  end if
451  ierrmono = 1
452  end if
453  pxdp1 = pxdp2
454  end do segloop
455  end do
456  !
457  ! -- terminate if errors encountered
458  if (count_errors() > 0) then
459  call store_error_filename(this%input_fname)
460  end if
Here is the call graph for this function:

◆ default_nodelist()

subroutine evtmodule::default_nodelist ( class(evttype this)
private

Equivalent to reading IEVT as CONSTANT 1

Definition at line 740 of file gwf-evt.f90.

741  ! -- modules
742  use simmodule, only: store_error
743  use constantsmodule, only: linelength
744  ! -- dummy
745  class(EvtType) :: this
746  ! -- local
747  integer(I4B) :: il, ir, ic, ncol, nrow, nlay, nodeu, noder, ipos
748  !
749  ! -- set variables
750  if (this%dis%ndim == 3) then
751  nlay = this%dis%mshape(1)
752  nrow = this%dis%mshape(2)
753  ncol = this%dis%mshape(3)
754  elseif (this%dis%ndim == 2) then
755  nlay = this%dis%mshape(1)
756  nrow = 1
757  ncol = this%dis%mshape(2)
758  end if
759  !
760  ! -- Populate nodelist
761  ipos = 1
762  il = 1
763  do ir = 1, nrow
764  do ic = 1, ncol
765  nodeu = get_node(il, ir, ic, nlay, nrow, ncol)
766  noder = this%dis%get_nodenumber(nodeu, 0)
767  this%nodelist(ipos) = noder
768  ipos = ipos + 1
769  end do
770  end do
771  !
772  ! -- assign nbound.
773  this%nbound = ipos - 1
774  !
775  ! -- if fixed_cell option not set, then need to store nodelist
776  ! in the nodesontop array
777  if (.not. this%fixed_cell) call this%set_nodesontop()
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
Here is the call graph for this function:

◆ evt_allocate_arrays()

subroutine evtmodule::evt_allocate_arrays ( class(evttype this,
integer(i4b), dimension(:), optional, pointer, contiguous  nodelist,
real(dp), dimension(:, :), optional, pointer, contiguous  auxvar 
)

Definition at line 132 of file gwf-evt.f90.

133  ! -- modules
135  ! -- dummy
136  class(EvtType) :: this
137  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
138  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
139  !
140  ! -- call standard BndType allocate scalars
141  call this%BndExtType%allocate_arrays(nodelist, auxvar)
142  !
143  ! -- set EVT input context pointers
144  call mem_setptr(this%surface, 'SURFACE', this%input_mempath)
145  call mem_setptr(this%rate, 'RATE', this%input_mempath)
146  call mem_setptr(this%depth, 'DEPTH', this%input_mempath)
147  !
148  ! -- checkin EVT input context pointers
149  call mem_checkin(this%surface, 'SURFACE', this%memoryPath, &
150  'SURFACE', this%input_mempath)
151  call mem_checkin(this%rate, 'RATE', this%memoryPath, &
152  'RATE', this%input_mempath)
153  call mem_checkin(this%depth, 'DEPTH', this%memoryPath, &
154  'DEPTH', this%input_mempath)
155  !
156  ! -- set list input segment descriptors
157  if (.not. this%read_as_arrays) then
158  if (this%nseg > 1) then
159  !
160  ! -- set pxdp and petm input context pointers
161  call mem_setptr(this%pxdp, 'PXDP', this%input_mempath)
162  call mem_setptr(this%petm, 'PETM', this%input_mempath)
163  !
164  ! -- checkin pxdp and petm input context pointers
165  call mem_checkin(this%pxdp, 'PXDP', this%memoryPath, &
166  'PXDP', this%input_mempath)
167  call mem_checkin(this%petm, 'PETM', this%memoryPath, &
168  'PETM', this%input_mempath)
169  end if
170  !
171  if (this%surfratespecified) then
172  !
173  ! -- set petm0 input context pointer
174  call mem_setptr(this%petm0, 'PETM0', this%input_mempath)
175  !
176  ! -- cehckin petm0 input context pointer
177  call mem_checkin(this%petm0, 'PETM0', this%memoryPath, &
178  'PETM0', this%input_mempath)
179  end if
180  end if

◆ evt_allocate_scalars()

subroutine evtmodule::evt_allocate_scalars ( class(evttype), intent(inout)  this)
private

Definition at line 104 of file gwf-evt.f90.

105  ! -- modules
107  ! -- dummy
108  class(EvtType), intent(inout) :: this
109  !
110  ! -- call standard BndType allocate scalars
111  call this%BndExtType%allocate_scalars()
112  !
113  ! -- allocate the object and assign values to object variables
114  call mem_allocate(this%nseg, 'NSEG', this%memoryPath)
115  !
116  ! -- allocate internal members
117  allocate (this%segsdefined)
118  allocate (this%fixed_cell)
119  allocate (this%read_as_arrays)
120  allocate (this%surfratespecified)
121  !
122  ! -- Set values
123  this%nseg = 1
124  this%segsdefined = .true.
125  this%fixed_cell = .false.
126  this%read_as_arrays = .false.
127  this%surfratespecified = .false.

◆ evt_bound_value()

real(dp) function evtmodule::evt_bound_value ( class(evttype), intent(inout)  this,
integer(i4b), intent(in)  col,
integer(i4b), intent(in)  row 
)
private
Parameters
[in,out]thisBndExtType object

Definition at line 809 of file gwf-evt.f90.

810  ! -- modules
811  use constantsmodule, only: dzero
812  ! -- dummy variables
813  class(EvtType), intent(inout) :: this !< BndExtType object
814  integer(I4B), intent(in) :: col
815  integer(I4B), intent(in) :: row
816  ! -- result
817  real(DP) :: bndval
818  ! -- local
819  integer(I4B) :: idx
820  !
821  ! -- initialize
822  idx = 0
823  !
824  select case (col)
825  case (1)
826  bndval = this%surface(row)
827  case (2)
828  if (this%iauxmultcol > 0) then
829  bndval = this%rate(row) * this%auxvar(this%iauxmultcol, row)
830  else
831  bndval = this%rate(row)
832  end if
833  case (3)
834  bndval = this%depth(row)
835  case default
836  if (col > 0) then
837  if (this%nseg > 1) then
838  if (col < (3 + this%nseg)) then
839  idx = col - 3
840  bndval = this%pxdp(idx, row)
841  else if (col < (3 + (2 * this%nseg) - 1)) then
842  idx = col - (3 + this%nseg - 1)
843  bndval = this%petm(idx, row)
844  else if (col == (3 + (2 * this%nseg) - 1)) then
845  if (this%surfratespecified) then
846  idx = 1
847  bndval = this%petm0(row)
848  end if
849  end if
850  else if (this%surfratespecified) then
851  if (col == 4) then
852  idx = 1
853  bndval = this%petm0(row)
854  end if
855  end if
856  end if
857  !
858  ! -- set error if idx not found
859  if (idx == 0) then
860  write (errmsg, '(a,i0,a)') &
861  'Programming error. EVT bound value requested column '&
862  &'outside range of ncolbnd (', this%ncolbnd, ').'
863  call store_error(errmsg)
864  call store_error_filename(this%input_fname)
865  end if
866  !
867  end select
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
Here is the call graph for this function:

◆ evt_cf()

subroutine evtmodule::evt_cf ( class(evttype this)
private

Definition at line 484 of file gwf-evt.f90.

485  ! -- dummy
486  class(EvtType) :: this
487  ! -- local
488  integer(I4B) :: i, iseg, node
489  integer(I4B) :: idxdepth, idxrate
490  real(DP) :: c, d, h, s, x
491  real(DP) :: petm0
492  real(DP) :: petm1, petm2, pxdp1, pxdp2, thcof, trhs
493  !
494  ! -- Return if no ET nodes
495  if (this%nbound == 0) return
496  !
497  ! -- Calculate hcof and rhs for each ET node
498  do i = 1, this%nbound
499  !
500  ! -- Find the node number
501  if (this%fixed_cell) then
502  node = this%nodelist(i)
503  else
504  node = this%nodesontop(i)
505  end if
506  !
507  ! -- cycle if nonexistent bound
508  if (node <= 0) then
509  this%hcof(i) = dzero
510  this%rhs(i) = dzero
511  cycle
512  end if
513  !
514  ! -- reset nodelist to highest active
515  if (.not. this%fixed_cell) then
516  if (this%ibound(node) == 0) &
517  call this%dis%highest_active(node, this%ibound)
518  this%nodelist(i) = node
519  end if
520  !
521  ! -- set rhs and hcof to zero
522  this%rhs(i) = dzero
523  this%hcof(i) = dzero
524  !
525  ! -- if ibound is positive and not overlain by a lake, then add terms
526  if (this%ibound(node) > 0 .and. this%ibound(node) /= iwetlake) then
527  !
528  if (this%iauxmultcol > 0) then
529  c = this%rate(i) * this%dis%get_area(node) * &
530  this%auxvar(this%iauxmultcol, i)
531  else
532  c = this%rate(i) * this%dis%get_area(node)
533  end if
534  s = this%surface(i)
535  h = this%xnew(node)
536  if (this%surfratespecified) then
537  petm0 = this%petm0(i)
538  end if
539  !
540  ! -- If head in cell is greater than or equal to SURFACE, ET is constant
541  if (h >= s) then
542  if (this%surfratespecified) then
543  ! -- Subtract -PETM0 * max rate from RHS
544  this%rhs(i) = this%rhs(i) + petm0 * c
545  else
546  ! -- Subtract -RATE from RHS
547  this%rhs(i) = this%rhs(i) + c
548  end if
549  else
550  ! -- If depth to water >= extinction depth, then ET is 0
551  d = s - h
552  x = this%depth(i)
553  if (d < x) then
554  ! -- Variable range. add ET terms to both RHS and HCOF.
555  if (this%nseg > 1) then
556  ! -- Determine which segment applies based on head, and
557  ! calculate terms to add to RHS and HCOF
558  !
559  ! -- Set proportions corresponding to surface elevation
560  ! and initial indices for depth and rate.
561  ! -- Idxdepth will point to the elements of bound containing
562  ! proportion of depth at the bottom of each segment.
563  ! Idxrate will point to the elements of bound containing
564  ! proportion of ET rate at the bottom of each segment.
565  ! If surfratespecified is true, skip over the elements
566  ! containing pxdp0 (=0.0) and petm0.
567  pxdp1 = dzero
568  if (this%surfratespecified) then
569  petm1 = petm0
570  else
571  petm1 = done
572  end if
573  ! -- Initialize indices to point to elements preceding
574  ! pxdp1 and petm1 (values for lower end of segment 1).
575  idxdepth = 0
576  idxrate = 0
577  ! -- Iterate through segments to find segment that contains
578  ! current depth of head below ET surface.
579  segloop: do iseg = 1, this%nseg
580  ! -- Set proportions corresponding to lower end of
581  ! segment
582  if (iseg < this%nseg) then
583  ! -- Increment the indices for depth and rate
584  idxdepth = idxdepth + 1
585  idxrate = idxrate + 1
586  ! -- Get proportions for lower end of segment
587  pxdp2 = this%pxdp(idxdepth, i)
588  petm2 = this%petm(idxrate, i)
589  else
590  pxdp2 = done
591  petm2 = dzero
592  end if
593  if (d <= pxdp2 * x) then
594  ! -- head is in domain of this segment
595  exit segloop
596  end if
597  ! -- Proportions at lower end of segment will be for
598  ! upper end of segment next time through loop
599  pxdp1 = pxdp2
600  petm1 = petm2
601  end do segloop
602  ! -- Calculate terms to add to RHS and HCOF based on
603  ! segment that applies at head elevation
604  thcof = -(petm1 - petm2) * c / ((pxdp2 - pxdp1) * x)
605  trhs = thcof * (s - pxdp1 * x) + petm1 * c
606  else
607  ! -- Calculate terms to add to RHS and HCOF based on simple
608  ! linear relation of ET vs. head for single segment
609  trhs = c - c * s / x
610  thcof = -c / x
611  end if
612  this%rhs(i) = this%rhs(i) + trhs
613  this%hcof(i) = this%hcof(i) + thcof
614  end if
615  end if
616  end if
617  !
618  end do

◆ evt_create()

subroutine, public evtmodule::evt_create ( class(bndtype), pointer  packobj,
integer(i4b), intent(in)  id,
integer(i4b), intent(in)  ibcnum,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout,
character(len=*), intent(in)  namemodel,
character(len=*), intent(in)  pakname,
character(len=*), intent(in)  mempath 
)

Definition at line 67 of file gwf-evt.f90.

69  ! -- dummy
70  class(BndType), pointer :: packobj
71  integer(I4B), intent(in) :: id
72  integer(I4B), intent(in) :: ibcnum
73  integer(I4B), intent(in) :: inunit
74  integer(I4B), intent(in) :: iout
75  character(len=*), intent(in) :: namemodel
76  character(len=*), intent(in) :: pakname
77  character(len=*), intent(in) :: mempath
78  ! -- local
79  type(EvtType), pointer :: evtobj
80  !
81  ! -- allocate evt object and scalar variables
82  allocate (evtobj)
83  packobj => evtobj
84  !
85  ! -- create name and memory path
86  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
87  packobj%text = text
88  !
89  ! -- allocate scalars
90  call evtobj%evt_allocate_scalars()
91  !
92  ! -- initialize package
93  call packobj%pack_initialize()
94 
95  packobj%inunit = inunit
96  packobj%iout = iout
97  packobj%id = id
98  packobj%ibcnum = ibcnum
99  packobj%ictMemPath = create_mem_path(namemodel, 'NPF')
Here is the call graph for this function:
Here is the caller graph for this function:

◆ evt_da()

subroutine evtmodule::evt_da ( class(evttype this)
private

Definition at line 651 of file gwf-evt.f90.

652  ! -- modules
654  ! -- dummy
655  class(EvtType) :: this
656  !
657  ! -- arrays
658  if (associated(this%nodesontop)) deallocate (this%nodesontop)
659  call mem_deallocate(this%surface, 'SURFACE', this%memoryPath)
660  call mem_deallocate(this%rate, 'RATE', this%memoryPath)
661  call mem_deallocate(this%depth, 'DEPTH', this%memoryPath)
662  !
663  if (.not. this%read_as_arrays) then
664  if (this%nseg > 1) then
665  call mem_deallocate(this%pxdp, 'PXDP', this%memoryPath)
666  call mem_deallocate(this%petm, 'PETM', this%memoryPath)
667  end if
668  !
669  if (this%surfratespecified) then
670  call mem_deallocate(this%petm0, 'PETM0', this%memoryPath)
671  end if
672  end if
673  !
674  ! -- scalars
675  call mem_deallocate(this%nseg)
676  deallocate (this%segsdefined)
677  deallocate (this%fixed_cell)
678  deallocate (this%read_as_arrays)
679  deallocate (this%surfratespecified)
680  !
681  ! -- Deallocate parent package
682  call this%BndExtType%bnd_da()

◆ evt_define_listlabel()

subroutine evtmodule::evt_define_listlabel ( class(evttype), intent(inout)  this)

Definition at line 688 of file gwf-evt.f90.

689  ! -- dummy
690  class(EvtType), intent(inout) :: this
691  ! -- local
692  integer(I4B) :: nsegm1, i
693  !
694  ! -- create the header list label
695  this%listlabel = trim(this%filtyp)//' NO.'
696  if (this%dis%ndim == 3) then
697  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
698  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
699  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
700  elseif (this%dis%ndim == 2) then
701  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
702  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
703  else
704  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
705  end if
706  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'SURFACE'
707  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'MAX. RATE'
708  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'EXT. DEPTH'
709  !
710  ! -- add headings for as many PXDP and PETM columns as needed
711  nsegm1 = this%nseg - 1
712  if (nsegm1 > 0) then
713  do i = 1, nsegm1
714  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'PXDP'
715  end do
716  do i = 1, nsegm1
717  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'PETM'
718  end do
719  end if
720  !
721  ! -- PETM0, if SURF_RATE_SPECIFIED is used
722  if (this%surfratespecified) then
723  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'PETM0'
724  end if
725  !
726 ! ! -- multiplier
727 ! if(this%multindex > 0) &
728 ! write(this%listlabel, '(a, a16)') trim(this%listlabel), 'MULTIPLIER'
729  !
730  ! -- boundary name
731  if (this%inamedbound == 1) then
732  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
733  end if

◆ evt_df_obs()

subroutine evtmodule::evt_df_obs ( class(evttype this)
private

Overrides BndTypebnd_df_obs

Definition at line 797 of file gwf-evt.f90.

798  ! -- dummy
799  class(EvtType) :: this
800  ! -- local
801  integer(I4B) :: indx
802  !
803  call this%obs%StoreObsType('evt', .true., indx)
804  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
Here is the call graph for this function:

◆ evt_fc()

subroutine evtmodule::evt_fc ( class(evttype this,
real(dp), dimension(:), intent(inout)  rhs,
integer(i4b), dimension(:), intent(in)  ia,
integer(i4b), dimension(:), intent(in)  idxglo,
class(matrixbasetype), pointer  matrix_sln 
)
private

Definition at line 623 of file gwf-evt.f90.

624  ! -- dummy
625  class(EvtType) :: this
626  real(DP), dimension(:), intent(inout) :: rhs
627  integer(I4B), dimension(:), intent(in) :: ia
628  integer(I4B), dimension(:), intent(in) :: idxglo
629  class(MatrixBaseType), pointer :: matrix_sln
630  ! -- local
631  integer(I4B) :: i, n, ipos
632  !
633  ! -- Copy package rhs and hcof into solution rhs and amat
634  do i = 1, this%nbound
635  n = this%nodelist(i)
636  if (n <= 0) cycle
637  ! -- reset hcof and rhs for excluded cells
638  if (this%ibound(n) == iwetlake) then
639  this%hcof(i) = dzero
640  this%rhs(i) = dzero
641  cycle
642  end if
643  rhs(n) = rhs(n) + this%rhs(i)
644  ipos = ia(n)
645  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
646  end do

◆ evt_log_options()

subroutine evtmodule::evt_log_options ( class(evttype), intent(inout)  this,
logical(lgp), intent(in)  found_fixed_cell,
logical(lgp), intent(in)  found_readasarrays,
logical(lgp), intent(in)  found_surfratespec 
)

Definition at line 233 of file gwf-evt.f90.

235  ! -- modules
239  ! -- dummy
240  class(EvtType), intent(inout) :: this
241  logical(LGP), intent(in) :: found_fixed_cell
242  logical(LGP), intent(in) :: found_readasarrays
243  logical(LGP), intent(in) :: found_surfratespec
244  ! -- formats
245  character(len=*), parameter :: fmtihact = &
246  &"(4x, 'EVAPOTRANSPIRATION WILL BE APPLIED TO HIGHEST ACTIVE CELL.')"
247  character(len=*), parameter :: fmtfixedcell = &
248  &"(4x, 'EVAPOTRANSPIRATION WILL BE APPLIED TO SPECIFIED CELL.')"
249  character(len=*), parameter :: fmtreadasarrays = &
250  &"(4x, 'EVAPOTRANSPIRATION INPUT WILL BE READ AS ARRAYS.')"
251  character(len=*), parameter :: fmtsrz = &
252  &"(4x, 'ET RATE AT SURFACE WILL BE ZERO.')"
253  character(len=*), parameter :: fmtsrs = &
254  &"(4x, 'ET RATE AT SURFACE WILL BE AS SPECIFIED BY PETM0.')"
255  !
256  ! -- log found options
257  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
258  //' OPTIONS'
259  !
260  if (found_fixed_cell) then
261  write (this%iout, fmtfixedcell)
262  end if
263  !
264  if (found_readasarrays) then
265  write (this%iout, fmtreadasarrays)
266  end if
267  !
268  if (found_surfratespec) then
269  write (this%iout, fmtsrs)
270  end if
271  !
272  ! -- close logging block
273  write (this%iout, '(1x,a)') &
274  'END OF '//trim(adjustl(this%text))//' OPTIONS'
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23

◆ evt_obs_supported()

logical function evtmodule::evt_obs_supported ( class(evttype this)

Overrides BndTypebnd_obs_supported()

Definition at line 786 of file gwf-evt.f90.

787  ! -- dummy
788  class(EvtType) :: this
789  !
790  evt_obs_supported = .true.

◆ evt_read_initial_attr()

subroutine evtmodule::evt_read_initial_attr ( class(evttype), intent(inout)  this)

If READASARRAYS has been specified, assign default IEVT = 1

Definition at line 352 of file gwf-evt.f90.

353  ! -- dummy
354  class(EvtType), intent(inout) :: this
355  !
356  if (this%read_as_arrays) then
357  call this%default_nodelist()
358  end if

◆ evt_rp()

subroutine evtmodule::evt_rp ( class(evttype), intent(inout)  this)
private

Read itmp and new boundaries if itmp > 0

Definition at line 365 of file gwf-evt.f90.

366  use tdismodule, only: kper
367  implicit none
368  ! -- dummy
369  class(EvtType), intent(inout) :: this
370  !
371  if (this%iper /= kper) return
372  !
373  if (this%read_as_arrays) then
374  !
375  ! -- update nodelist based on IRCH input
376  call nodelist_update(this%nodelist, this%nbound, this%maxbound, &
377  this%dis, this%input_mempath)
378  !
379  else
380  !
381  ! -- process the input list arrays
382  call this%BndExtType%bnd_rp()
383  !
384  ! -- ensure pxdp is monotonically increasing
385  if (this%nseg > 1) then
386  call this%check_pxdp()
387  end if
388  !
389  ! -- Write the list to iout if requested
390  if (this%iprpak /= 0) then
391  call this%write_list()
392  end if
393  !
394  end if
395  !
396  ! -- copy nodelist to nodesontop if not fixed cell
397  if (.not. this%fixed_cell) call this%set_nodesontop()
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
Here is the call graph for this function:

◆ evt_source_dimensions()

subroutine evtmodule::evt_source_dimensions ( class(evttype), intent(inout)  this)

Definition at line 279 of file gwf-evt.f90.

280  ! -- modules
282  ! -- dummy
283  class(EvtType), intent(inout) :: this
284  ! -- local
285  logical(LGP) :: found_nseg = .false.
286  ! -- format
287  character(len=*), parameter :: fmtnsegerr = &
288  &"('Error: In EVT, NSEG must be > 0 but is specified as ',i0)"
289  !
290  ! Dimensions block is not required if:
291  ! (1) discretization is DIS or DISV, and
292  ! (2) READASARRAYS option has been specified.
293  if (this%read_as_arrays) then
294  this%maxbound = this%dis%get_ncpl()
295  !
296  ! -- verify dimensions were set
297  if (this%maxbound <= 0) then
298  write (errmsg, '(a)') &
299  'MAXBOUND must be an integer greater than zero.'
300  call store_error(errmsg)
301  call store_error_filename(this%input_fname)
302  end if
303  !
304  else
305  !
306  ! -- source maxbound
307  call this%BndExtType%source_dimensions()
308  !
309  ! -- log found options
310  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
311  //' DIMENSIONS'
312  !
313  ! -- update defaults with idm sourced values
314  call mem_set_value(this%nseg, 'NSEG', this%input_mempath, found_nseg)
315  !
316  if (found_nseg) then
317  !
318  write (this%iout, '(4x,a,i0)') 'NSEG = ', this%nseg
319  !
320  if (this%nseg < 1) then
321  write (errmsg, fmtnsegerr) this%nseg
322  call store_error(errmsg)
323  call store_error_filename(this%input_fname)
324  !
325  elseif (this%nseg > 1) then
326  ! NSEG>1 is supported only if readasarrays is false
327  if (this%read_as_arrays) then
328  errmsg = 'In the EVT package, NSEG cannot be greater than 1'// &
329  ' when READASARRAYS is used.'
330  call store_error(errmsg)
331  call store_error_filename(this%input_fname)
332  end if
333  !
334  end if
335  end if
336  !
337  ! -- close logging block
338  write (this%iout, '(1x,a)') &
339  'END OF '//trim(adjustl(this%text))//' DIMENSIONS'
340  !
341  end if
342  !
343  ! -- Call define_listlabel to construct the list label that is written
344  ! when PRINT_INPUT option is used.
345  call this%define_listlabel()
Here is the call graph for this function:

◆ evt_source_options()

subroutine evtmodule::evt_source_options ( class(evttype), intent(inout)  this)

Definition at line 185 of file gwf-evt.f90.

186  ! -- modules
188  ! -- dummy
189  class(EvtType), intent(inout) :: this
190  ! -- local
191  logical(LGP) :: found_fixed_cell = .false.
192  logical(LGP) :: found_readasarrays = .false.
193  logical(LGP) :: found_surfratespec = .false.
194  !
195  ! -- source common bound options
196  call this%BndExtType%source_options()
197  !
198  ! -- update defaults with idm sourced values
199  call mem_set_value(this%fixed_cell, 'FIXED_CELL', &
200  this%input_mempath, found_fixed_cell)
201  call mem_set_value(this%read_as_arrays, 'READASARRAYS', &
202  this%input_mempath, found_readasarrays)
203  call mem_set_value(this%surfratespecified, 'SURFRATESPEC', &
204  this%input_mempath, found_surfratespec)
205  !
206  if (found_readasarrays) then
207  if (this%dis%supports_layers()) then
208  this%text = texta
209  else
210  errmsg = 'READASARRAYS option is not compatible with selected'// &
211  ' discretization type.'
212  call store_error(errmsg)
213  call store_error_filename(this%input_fname)
214  end if
215  end if
216  !
217  if (found_readasarrays .and. found_surfratespec) then
218  if (this%read_as_arrays) then
219  errmsg = 'READASARRAYS option is not compatible with the'// &
220  ' SURF_RATE_SPECIFIED option.'
221  call store_error(errmsg)
222  call store_error_filename(this%input_fname)
223  end if
224  end if
225  !
226  ! -- log evt specific options
227  call this%evt_log_options(found_fixed_cell, found_readasarrays, &
228  found_surfratespec)
Here is the call graph for this function:

◆ nodelist_update()

subroutine evtmodule::nodelist_update ( integer(i4b), dimension(:), intent(inout), pointer, contiguous  nodelist,
integer(i4b), intent(inout)  nbound,
integer(i4b), intent(in)  maxbound,
class(disbasetype), intent(in), pointer  dis,
character(len=*), intent(in)  input_mempath 
)

This is a module scoped routine to check for IEVT input. If array input was provided, INIEVT and IEVT will be allocated in the input context. If the read state variable INIEVT is set to 1 during this period update, IEVT input was read and is used here to update the nodelist.

Definition at line 877 of file gwf-evt.f90.

879  ! -- modules
881  use basedismodule, only: disbasetype
882  ! -- dummy
883  integer(I4B), dimension(:), contiguous, &
884  pointer, intent(inout) :: nodelist
885  class(DisBaseType), pointer, intent(in) :: dis
886  character(len=*), intent(in) :: input_mempath
887  integer(I4B), intent(inout) :: nbound
888  integer(I4B), intent(in) :: maxbound
889  ! -- format
890  character(len=24) :: aname = ' LAYER OR NODE INDEX'
891  ! -- local
892  integer(I4B), dimension(:), contiguous, pointer :: ievt => null()
893  integer(I4B), pointer :: inievt => null()
894  !
895  ! -- set pointer to input context INIEVT
896  call mem_setptr(inievt, 'INIEVT', input_mempath)
897  !
898  ! -- check INIEVT read state
899  if (inievt == 1) then
900  ! -- ievt was read this period
901  !
902  ! -- set pointer to input context IEVT
903  call mem_setptr(ievt, 'IEVT', input_mempath)
904  !
905  ! -- update nodelist
906  call dis%nlarray_to_nodelist(ievt, nodelist, &
907  maxbound, nbound, aname)
908  end if
Here is the caller graph for this function:

◆ set_nodesontop()

subroutine evtmodule::set_nodesontop ( class(evttype), intent(inout)  this)
private

Definition at line 465 of file gwf-evt.f90.

466  ! -- dummy
467  class(EvtType), intent(inout) :: this
468  ! -- local
469  integer(I4B) :: n
470  !
471  ! -- allocate if necessary
472  if (.not. associated(this%nodesontop)) then
473  allocate (this%nodesontop(this%maxbound))
474  end if
475  !
476  ! -- copy nodelist into nodesontop
477  do n = 1, this%nbound
478  this%nodesontop(n) = this%nodelist(n)
479  end do

Variable Documentation

◆ ftype

character(len=lenftype) evtmodule::ftype = 'EVT'
private

Definition at line 21 of file gwf-evt.f90.

21  character(len=LENFTYPE) :: ftype = 'EVT'

◆ text

character(len=lenpackagename) evtmodule::text = ' EVT'
private

Definition at line 22 of file gwf-evt.f90.

22  character(len=LENPACKAGENAME) :: text = ' EVT'

◆ texta

character(len=lenpackagename) evtmodule::texta = ' EVTA'
private

Definition at line 23 of file gwf-evt.f90.

23  character(len=LENPACKAGENAME) :: texta = ' EVTA'