MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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 431 of file gwf-evt.f90.

432  ! -- dummy
433  class(EvtType), intent(inout) :: this !< EvtType
434  ! -- local
435  integer(I4B) :: n
436  integer(I4B) :: node
437  integer(I4B) :: i
438  integer(I4B) :: ierrmono
439  real(DP) :: pxdp1
440  real(DP) :: pxdp2
441  character(len=15) :: nodestr
442  ! -- formats
443  character(len=*), parameter :: fmtpxdp0 = &
444  &"('PXDP must be between 0 and 1. Found ', G0, ' for cell ', A)"
445  character(len=*), parameter :: fmtpxdp = &
446  &"('PXDP is not monotonically increasing for cell ', A)"
447  !
448  ! -- check and make sure that pxdp is monotonically increasing and
449  ! that pxdp values are between 0 and 1
450  do n = 1, this%nbound
451  node = this%nodelist(n)
452  pxdp1 = dzero
453  ierrmono = 0
454  segloop: do i = 1, this%nseg
455  !
456  ! -- set and check pxdp2
457  if (i < this%nseg) then
458  pxdp2 = this%pxdp(i, n)
459  if (pxdp2 <= dzero .or. pxdp2 >= done) then
460  call this%dis%noder_to_string(node, nodestr)
461  write (errmsg, fmtpxdp0) pxdp2, trim(nodestr)
462  call store_error(errmsg)
463  end if
464  else
465  pxdp2 = done
466  end if
467  !
468  ! -- check for monotonically increasing condition
469  if (pxdp2 - pxdp1 < dzero) then
470  if (ierrmono == 0) then
471  ! -- only store mono error once for each node
472  call this%dis%noder_to_string(node, nodestr)
473  write (errmsg, fmtpxdp) trim(nodestr)
474  call store_error(errmsg)
475  end if
476  ierrmono = 1
477  end if
478  pxdp1 = pxdp2
479  end do segloop
480  end do
481  !
482  ! -- terminate if errors encountered
483  if (count_errors() > 0) then
484  call this%parser%StoreErrorUnit()
485  end if
486  !
487  ! -- Return
488  return
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 783 of file gwf-evt.f90.

784  ! -- modules
785  use simmodule, only: store_error
786  use constantsmodule, only: linelength
787  ! -- dummy
788  class(EvtType) :: this
789  ! -- local
790  integer(I4B) :: il, ir, ic, ncol, nrow, nlay, nodeu, noder, ipos
791  !
792  ! -- set variables
793  if (this%dis%ndim == 3) then
794  nlay = this%dis%mshape(1)
795  nrow = this%dis%mshape(2)
796  ncol = this%dis%mshape(3)
797  elseif (this%dis%ndim == 2) then
798  nlay = this%dis%mshape(1)
799  nrow = 1
800  ncol = this%dis%mshape(2)
801  end if
802  !
803  ! -- Populate nodelist
804  ipos = 1
805  il = 1
806  do ir = 1, nrow
807  do ic = 1, ncol
808  nodeu = get_node(il, ir, ic, nlay, nrow, ncol)
809  noder = this%dis%get_nodenumber(nodeu, 0)
810  this%nodelist(ipos) = noder
811  ipos = ipos + 1
812  end do
813  end do
814  !
815  ! -- assign nbound.
816  this%nbound = ipos - 1
817  !
818  ! -- if fixed_cell option not set, then need to store nodelist
819  ! in the nodesontop array
820  if (.not. this%fixed_cell) call this%set_nodesontop()
821  !
822  ! -- Return
823  return
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
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 139 of file gwf-evt.f90.

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

◆ evt_allocate_scalars()

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

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

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

◆ 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 861 of file gwf-evt.f90.

862  ! -- modules
863  use constantsmodule, only: dzero
864  ! -- dummy variables
865  class(EvtType), intent(inout) :: this !< BndExtType object
866  integer(I4B), intent(in) :: col
867  integer(I4B), intent(in) :: row
868  ! -- result
869  real(DP) :: bndval
870  ! -- local
871  integer(I4B) :: idx
872  !
873  ! -- initialize
874  idx = 0
875  !
876  select case (col)
877  case (1)
878  bndval = this%surface(row)
879  case (2)
880  if (this%iauxmultcol > 0) then
881  bndval = this%rate(row) * this%auxvar(this%iauxmultcol, row)
882  else
883  bndval = this%rate(row)
884  end if
885  case (3)
886  bndval = this%depth(row)
887  case default
888  if (col > 0) then
889  if (this%nseg > 1) then
890  if (col < (3 + this%nseg)) then
891  idx = col - 3
892  bndval = this%pxdp(idx, row)
893  else if (col < (3 + (2 * this%nseg) - 1)) then
894  idx = col - (3 + this%nseg - 1)
895  bndval = this%petm(idx, row)
896  else if (col == (3 + (2 * this%nseg) - 1)) then
897  if (this%surfratespecified) then
898  idx = 1
899  bndval = this%petm0(row)
900  end if
901  end if
902  else if (this%surfratespecified) then
903  if (col == 4) then
904  idx = 1
905  bndval = this%petm0(row)
906  end if
907  end if
908  end if
909  !
910  ! -- set error if idx not found
911  if (idx == 0) then
912  write (errmsg, '(a,i0,a)') &
913  'Programming error. EVT bound value requested column '&
914  &'outside range of ncolbnd (', this%ncolbnd, ').'
915  call store_error(errmsg)
916  call store_error_filename(this%input_fname)
917  end if
918  !
919  end select
920  !
921  ! -- Return
922  return
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
Here is the call graph for this function:

◆ evt_cf()

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

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

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

◆ 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 68 of file gwf-evt.f90.

70  ! -- dummy
71  class(BndType), pointer :: packobj
72  integer(I4B), intent(in) :: id
73  integer(I4B), intent(in) :: ibcnum
74  integer(I4B), intent(in) :: inunit
75  integer(I4B), intent(in) :: iout
76  character(len=*), intent(in) :: namemodel
77  character(len=*), intent(in) :: pakname
78  character(len=*), intent(in) :: mempath
79  ! -- local
80  type(EvtType), pointer :: evtobj
81  !
82  ! -- allocate evt object and scalar variables
83  allocate (evtobj)
84  packobj => evtobj
85  !
86  ! -- create name and memory path
87  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
88  packobj%text = text
89  !
90  ! -- allocate scalars
91  call evtobj%evt_allocate_scalars()
92  !
93  ! -- initialize package
94  call packobj%pack_initialize()
95 
96  packobj%inunit = inunit
97  packobj%iout = iout
98  packobj%id = id
99  packobj%ibcnum = ibcnum
100  packobj%ictMemPath = create_mem_path(namemodel, 'NPF')
101  !
102  ! -- Return
103  return
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 688 of file gwf-evt.f90.

689  ! -- modules
691  ! -- dummy
692  class(EvtType) :: this
693  !
694  ! -- arrays
695  if (associated(this%nodesontop)) deallocate (this%nodesontop)
696  call mem_deallocate(this%surface, 'SURFACE', this%memoryPath)
697  call mem_deallocate(this%rate, 'RATE', this%memoryPath)
698  call mem_deallocate(this%depth, 'DEPTH', this%memoryPath)
699  !
700  if (.not. this%read_as_arrays) then
701  if (this%nseg > 1) then
702  call mem_deallocate(this%pxdp, 'PXDP', this%memoryPath)
703  call mem_deallocate(this%petm, 'PETM', this%memoryPath)
704  end if
705  !
706  if (this%surfratespecified) then
707  call mem_deallocate(this%petm0, 'PETM0', this%memoryPath)
708  end if
709  end if
710  !
711  ! -- scalars
712  call mem_deallocate(this%nseg)
713  deallocate (this%segsdefined)
714  deallocate (this%fixed_cell)
715  deallocate (this%read_as_arrays)
716  deallocate (this%surfratespecified)
717  !
718  ! -- Deallocate parent package
719  call this%BndExtType%bnd_da()
720  !
721  ! -- Return
722  return

◆ evt_define_listlabel()

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

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

729  ! -- dummy
730  class(EvtType), intent(inout) :: this
731  ! -- local
732  integer(I4B) :: nsegm1, i
733  !
734  ! -- create the header list label
735  this%listlabel = trim(this%filtyp)//' NO.'
736  if (this%dis%ndim == 3) then
737  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
738  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
739  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
740  elseif (this%dis%ndim == 2) then
741  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
742  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
743  else
744  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
745  end if
746  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'SURFACE'
747  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'MAX. RATE'
748  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'EXT. DEPTH'
749  !
750  ! -- add headings for as many PXDP and PETM columns as needed
751  nsegm1 = this%nseg - 1
752  if (nsegm1 > 0) then
753  do i = 1, nsegm1
754  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'PXDP'
755  end do
756  do i = 1, nsegm1
757  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'PETM'
758  end do
759  end if
760  !
761  ! -- PETM0, if SURF_RATE_SPECIFIED is used
762  if (this%surfratespecified) then
763  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'PETM0'
764  end if
765  !
766 ! ! -- multiplier
767 ! if(this%multindex > 0) &
768 ! write(this%listlabel, '(a, a16)') trim(this%listlabel), 'MULTIPLIER'
769  !
770  ! -- boundary name
771  if (this%inamedbound == 1) then
772  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
773  end if
774  !
775  ! -- Return
776  return

◆ evt_df_obs()

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

Overrides BndTypebnd_df_obs

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

847  ! -- dummy
848  class(EvtType) :: this
849  ! -- local
850  integer(I4B) :: indx
851  !
852  call this%obs%StoreObsType('evt', .true., indx)
853  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
854  !
855  ! -- Return
856  return
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 657 of file gwf-evt.f90.

658  ! -- dummy
659  class(EvtType) :: this
660  real(DP), dimension(:), intent(inout) :: rhs
661  integer(I4B), dimension(:), intent(in) :: ia
662  integer(I4B), dimension(:), intent(in) :: idxglo
663  class(MatrixBaseType), pointer :: matrix_sln
664  ! -- local
665  integer(I4B) :: i, n, ipos
666  !
667  ! -- Copy package rhs and hcof into solution rhs and amat
668  do i = 1, this%nbound
669  n = this%nodelist(i)
670  if (n <= 0) cycle
671  ! -- reset hcof and rhs for excluded cells
672  if (this%ibound(n) == iwetlake) then
673  this%hcof(i) = dzero
674  this%rhs(i) = dzero
675  cycle
676  end if
677  rhs(n) = rhs(n) + this%rhs(i)
678  ipos = ia(n)
679  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
680  end do
681  !
682  ! -- Return
683  return

◆ 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 246 of file gwf-evt.f90.

248  ! -- modules
252  ! -- dummy
253  class(EvtType), intent(inout) :: this
254  logical(LGP), intent(in) :: found_fixed_cell
255  logical(LGP), intent(in) :: found_readasarrays
256  logical(LGP), intent(in) :: found_surfratespec
257  ! -- formats
258  character(len=*), parameter :: fmtihact = &
259  &"(4x, 'EVAPOTRANSPIRATION WILL BE APPLIED TO HIGHEST ACTIVE CELL.')"
260  character(len=*), parameter :: fmtfixedcell = &
261  &"(4x, 'EVAPOTRANSPIRATION WILL BE APPLIED TO SPECIFIED CELL.')"
262  character(len=*), parameter :: fmtreadasarrays = &
263  &"(4x, 'EVAPOTRANSPIRATION INPUT WILL BE READ AS ARRAYS.')"
264  character(len=*), parameter :: fmtsrz = &
265  &"(4x, 'ET RATE AT SURFACE WILL BE ZERO.')"
266  character(len=*), parameter :: fmtsrs = &
267  &"(4x, 'ET RATE AT SURFACE WILL BE AS SPECIFIED BY PETM0.')"
268  !
269  ! -- log found options
270  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
271  //' OPTIONS'
272  !
273  if (found_fixed_cell) then
274  write (this%iout, fmtfixedcell)
275  end if
276  !
277  if (found_readasarrays) then
278  write (this%iout, fmtreadasarrays)
279  end if
280  !
281  if (found_surfratespec) then
282  write (this%iout, fmtsrs)
283  end if
284  !
285  ! -- close logging block
286  write (this%iout, '(1x,a)') &
287  'END OF '//trim(adjustl(this%text))//' OPTIONS'
288  !
289  ! -- Return
290  return
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 832 of file gwf-evt.f90.

833  ! -- dummy
834  class(EvtType) :: this
835  !
836  evt_obs_supported = .true.
837  !
838  ! -- Return
839  return

◆ 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 371 of file gwf-evt.f90.

372  ! -- dummy
373  class(EvtType), intent(inout) :: this
374  !
375  if (this%read_as_arrays) then
376  call this%default_nodelist()
377  end if
378  !
379  ! -- Return
380  return

◆ evt_rp()

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

Read itmp and new boundaries if itmp > 0

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

388  use tdismodule, only: kper
389  implicit none
390  ! -- dummy
391  class(EvtType), intent(inout) :: this
392  !
393  if (this%iper /= kper) return
394  !
395  if (this%read_as_arrays) then
396  !
397  ! -- update nodelist based on IRCH input
398  call nodelist_update(this%nodelist, this%nbound, this%maxbound, &
399  this%dis, this%input_mempath)
400  !
401  else
402  !
403  ! -- process the input list arrays
404  call this%BndExtType%bnd_rp()
405  !
406  ! -- ensure pxdp is monotonically increasing
407  if (this%nseg > 1) then
408  call this%check_pxdp()
409  end if
410  !
411  ! -- Write the list to iout if requested
412  if (this%iprpak /= 0) then
413  call this%write_list()
414  end if
415  !
416  end if
417  !
418  ! -- copy nodelist to nodesontop if not fixed cell
419  if (.not. this%fixed_cell) call this%set_nodesontop()
420  !
421  ! -- Return
422  return
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 295 of file gwf-evt.f90.

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

◆ evt_source_options()

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

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

196  ! -- modules
198  ! -- dummy
199  class(EvtType), intent(inout) :: this
200  ! -- local
201  logical(LGP) :: found_fixed_cell = .false.
202  logical(LGP) :: found_readasarrays = .false.
203  logical(LGP) :: found_surfratespec = .false.
204  !
205  ! -- source common bound options
206  call this%BndExtType%source_options()
207  !
208  ! -- update defaults with idm sourced values
209  call mem_set_value(this%fixed_cell, 'FIXED_CELL', &
210  this%input_mempath, found_fixed_cell)
211  call mem_set_value(this%read_as_arrays, 'READASARRAYS', &
212  this%input_mempath, found_readasarrays)
213  call mem_set_value(this%surfratespecified, 'SURFRATESPEC', &
214  this%input_mempath, found_surfratespec)
215  !
216  if (found_readasarrays) then
217  if (this%dis%supports_layers()) then
218  this%text = texta
219  else
220  errmsg = 'READASARRAYS option is not compatible with selected'// &
221  ' discretization type.'
222  call store_error(errmsg)
223  call store_error_filename(this%input_fname)
224  end if
225  end if
226  !
227  if (found_readasarrays .and. found_surfratespec) then
228  if (this%read_as_arrays) then
229  errmsg = 'READASARRAYS option is not compatible with the'// &
230  ' SURF_RATE_SPECIFIED option.'
231  call store_error(errmsg)
232  call store_error_filename(this%input_fname)
233  end if
234  end if
235  !
236  ! -- log evt specific options
237  call this%evt_log_options(found_fixed_cell, found_readasarrays, &
238  found_surfratespec)
239  !
240  ! -- Return
241  return
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 932 of file gwf-evt.f90.

934  ! -- modules
936  use basedismodule, only: disbasetype
937  ! -- dummy
938  integer(I4B), dimension(:), contiguous, &
939  pointer, intent(inout) :: nodelist
940  class(DisBaseType), pointer, intent(in) :: dis
941  character(len=*), intent(in) :: input_mempath
942  integer(I4B), intent(inout) :: nbound
943  integer(I4B), intent(in) :: maxbound
944  ! -- format
945  character(len=24) :: aname = ' LAYER OR NODE INDEX'
946  ! -- local
947  integer(I4B), dimension(:), contiguous, pointer :: ievt => null()
948  integer(I4B), pointer :: inievt => null()
949  !
950  ! -- set pointer to input context INIEVT
951  call mem_setptr(inievt, 'INIEVT', input_mempath)
952  !
953  ! -- check INIEVT read state
954  if (inievt == 1) then
955  ! -- ievt was read this period
956  !
957  ! -- set pointer to input context IEVT
958  call mem_setptr(ievt, 'IEVT', input_mempath)
959  !
960  ! -- update nodelist
961  call dis%nlarray_to_nodelist(ievt, nodelist, &
962  maxbound, nbound, aname)
963  end if
964  !
965  ! -- Return
966  return
Here is the caller graph for this function:

◆ set_nodesontop()

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

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

494  ! -- dummy
495  class(EvtType), intent(inout) :: this
496  ! -- local
497  integer(I4B) :: n
498  !
499  ! -- allocate if necessary
500  if (.not. associated(this%nodesontop)) then
501  allocate (this%nodesontop(this%maxbound))
502  end if
503  !
504  ! -- copy nodelist into nodesontop
505  do n = 1, this%nbound
506  this%nodesontop(n) = this%nodelist(n)
507  end do
508  !
509  ! -- Return
510  return

Variable Documentation

◆ ftype

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

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

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

◆ text

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

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

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

◆ texta

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

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

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