MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
swfevpmodule Module Reference

This module contains the evaporation (EVP) package methods. More...

Data Types

type  swfevptype
 

Functions/Subroutines

subroutine, public evp_create (packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath, dis, dfw, cxs)
 Create a Evaporation Package. More...
 
subroutine evp_allocate_scalars (this)
 Allocate scalar members. More...
 
subroutine evp_allocate_arrays (this, nodelist, auxvar)
 Allocate package arrays. More...
 
subroutine evp_source_options (this)
 Source options specific to EVPType. More...
 
subroutine log_evp_options (this, found_readasarrays)
 Log options specific to SwfEvpType. More...
 
subroutine evp_source_dimensions (this)
 Source the dimensions for this package. More...
 
subroutine evp_read_initial_attr (this)
 Part of allocate and read. More...
 
subroutine evp_rp (this)
 Read and Prepare. More...
 
subroutine evp_ck (this)
 Ensure evaporation is positive. More...
 
subroutine evp_cf (this)
 Formulate the HCOF and RHS terms. More...
 
real(dp) function get_qevp (this, node, rlen, snew, sold, evaporation)
 Calculate qevp. More...
 
real(dp) function get_evap_reduce_mult (this, stage, bottom)
 Calculate multiplier to reduce evap as depth goes to zero. More...
 
subroutine evp_fc (this, rhs, ia, idxglo, matrix_sln)
 Copy rhs and hcof into solution rhs and amat. More...
 
subroutine evp_da (this)
 Deallocate memory. More...
 
subroutine evp_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 evp_obs_supported (this)
 Overrides BndTypebnd_obs_supported() More...
 
subroutine evp_df_obs (this)
 Implements bnd_df_obs. More...
 
real(dp) function evp_bound_value (this, col, row)
 Return requested boundary value. More...
 
real(dp) function, dimension(:), pointer reach_length_pointer (this)
 

Variables

character(len=lenftype) ftype = 'EVP'
 
character(len=lenpackagename) text = ' EVP'
 

Detailed Description

This module can be used to represent evaporation onto streams and overland flow cells.

Function/Subroutine Documentation

◆ default_nodelist()

subroutine swfevpmodule::default_nodelist ( class(swfevptype this)
private

Definition at line 550 of file swf-evp.f90.

551  ! dummy
552  class(SwfEvpType) :: this
553  ! local
554  integer(I4B) :: nodeu, noder
555 
556  ! This is only called for readasarrays, so nodelist will be the size of
557  ! the user grid, and will have a value of 0 for any entries where idomain
558  ! is not 1
559  do nodeu = 1, this%maxbound
560  noder = this%dis%get_nodenumber(nodeu, 0)
561  this%nodelist(nodeu) = noder
562  end do
563 
564  ! Assign nbound
565  this%nbound = this%maxbound
566 

◆ evp_allocate_arrays()

subroutine swfevpmodule::evp_allocate_arrays ( class(swfevptype this,
integer(i4b), dimension(:), optional, pointer, contiguous  nodelist,
real(dp), dimension(:, :), optional, pointer, contiguous  auxvar 
)
private

Definition at line 148 of file swf-evp.f90.

149  ! modules
151  ! dummy
152  class(SwfEvpType) :: this
153  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
154  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
155 
156  ! allocate base arrays
157  call this%BndExtType%allocate_arrays(nodelist, auxvar)
158 
159  ! set input context pointers
160  call mem_setptr(this%evaporation, 'EVAPORATION', this%input_mempath)
161 
162  ! checkin input context pointers
163  call mem_checkin(this%evaporation, 'EVAPORATION', this%memoryPath, &
164  'EVAPORATION', this%input_mempath)

◆ evp_allocate_scalars()

subroutine swfevpmodule::evp_allocate_scalars ( class(swfevptype), intent(inout)  this)
private

Definition at line 128 of file swf-evp.f90.

129  ! dummy
130  class(SwfEvpType), intent(inout) :: this
131 
132  ! allocate base scalars
133  call this%BndExtType%allocate_scalars()
134 
135  ! allocate internal members
136  call mem_allocate(this%iflowred, 'IFLOWRED', this%memoryPath)
137  call mem_allocate(this%reduction_depth, 'REDUCTION_DEPTH', this%memoryPath)
138  allocate (this%read_as_arrays)
139 
140  ! Set values
141  this%iflowred = 1
142  this%reduction_depth = dem6
143  this%read_as_arrays = .false.

◆ evp_bound_value()

real(dp) function swfevpmodule::evp_bound_value ( class(swfevptype), intent(inout)  this,
integer(i4b), intent(in)  col,
integer(i4b), intent(in)  row 
)
private

Definition at line 600 of file swf-evp.f90.

601  ! modules
602  use constantsmodule, only: dzero
603  ! dummy
604  class(SwfEvpType), intent(inout) :: this
605  integer(I4B), intent(in) :: col
606  integer(I4B), intent(in) :: row
607  ! result
608  real(DP) :: bndval
609 
610  select case (col)
611  case (1)
612  if (this%iauxmultcol > 0) then
613  bndval = this%evaporation(row) * this%auxvar(this%iauxmultcol, row)
614  else
615  bndval = this%evaporation(row)
616  end if
617  case default
618  errmsg = 'Programming error. EVP bound value requested column '&
619  &'outside range of ncolbnd (1).'
620  call store_error(errmsg)
621  call store_error_filename(this%input_fname)
622  end select
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
Here is the call graph for this function:

◆ evp_cf()

subroutine swfevpmodule::evp_cf ( class(swfevptype this)
private

Skip if no evaporation. Otherwise, calculate hcof and rhs

Definition at line 310 of file swf-evp.f90.

311  ! dummy
312  class(SwfEvpType) :: this
313  ! local
314  integer(I4B) :: i
315  integer(I4B) :: node
316  integer(I4B) :: inwt
317  real(DP) :: q
318  real(DP) :: qeps
319  real(DP) :: eps
320  real(DP) :: derv
321  real(DP) :: evap
322  real(DP) :: rlen
323  real(DP), dimension(:), pointer :: reach_length
324 
325  ! Return if no evaporation
326  if (this%nbound == 0) return
327 
328  ! Set pointer to reach_length for 1d
329  reach_length => this%reach_length_pointer()
330  rlen = dzero
331 
332  ! Calculate hcof and rhs for each evaporation entry
333  do i = 1, this%nbound
334 
335  ! Find the node number
336  node = this%nodelist(i)
337 
338  ! cycle if nonexistent bound
339  if (node <= 0) then
340  this%hcof(i) = dzero
341  this%rhs(i) = dzero
342  cycle
343  end if
344 
345  ! cycle if dry or constant head
346  if (this%ibound(node) <= 0) then
347  this%hcof(i) = dzero
348  this%rhs(i) = dzero
349  cycle
350  end if
351 
352  ! Initialize hcof
353  this%hcof(i) = dzero
354 
355  ! assign evap rate in length per time and multiply by auxvar
356  evap = this%evaporation(i)
357  if (this%iauxmultcol > 0) then
358  evap = evap * this%auxvar(this%iauxmultcol, i)
359  end if
360 
361  ! get reach length for 1d channel
362  if (this%dis%is_1d()) then
363  rlen = reach_length(node)
364  end if
365 
366  ! Calculate volumetric evaporation flow in L^3/Td and add to rhs
367  q = -this%get_qevp(node, rlen, this%xnew(node), this%xold(node), evap)
368  this%rhs(i) = -q
369 
370  ! Code for adding newton terms
371  inwt = 1
372  if (inwt == 1) then
373 
374  ! calculate perturbed q
375  eps = get_perturbation(this%xnew(node))
376  qeps = -this%get_qevp(node, rlen, this%xnew(node) + eps, &
377  this%xold(node), evap)
378 
379  ! calculate derivative
380  derv = (qeps - q) / eps
381 
382  ! add derivative to hcof and update rhs with derivate contribution
383  this%hcof(i) = derv
384  this%rhs(i) = this%rhs(i) + derv * this%xnew(node)
385  end if
386 
387  end do
Here is the call graph for this function:

◆ evp_ck()

subroutine swfevpmodule::evp_ck ( class(swfevptype), intent(inout)  this)

Definition at line 279 of file swf-evp.f90.

280  ! dummy
281  class(SwfEvpType), intent(inout) :: this
282  ! local
283  character(len=30) :: nodestr
284  integer(I4B) :: i, nr
285  character(len=*), parameter :: fmterr = &
286  &"('Specified stress ',i0, &
287  &' evaporation (',g0,') is less than zero for cell', a)"
288 
289  ! Ensure evaporation rates are positive
290  do i = 1, this%nbound
291  nr = this%nodelist(i)
292  if (nr <= 0) cycle
293  if (this%evaporation(i) < dzero) then
294  call this%dis%noder_to_string(nr, nodestr)
295  write (errmsg, fmt=fmterr) i, this%evaporation(i), trim(nodestr)
296  call store_error(errmsg)
297  end if
298  end do
299 
300  ! write summary of package error messages
301  if (count_errors() > 0) then
302  call store_error_filename(this%input_fname)
303  end if
Here is the call graph for this function:

◆ evp_create()

subroutine, public swfevpmodule::evp_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,
class(disbasetype), intent(inout), pointer  dis,
type(swfdfwtype), intent(in), pointer  dfw,
type(swfcxstype), intent(in), pointer  cxs 
)
Parameters
packobjpointer to default package type
[in]idpackage id
[in]ibcnumboundary condition number
[in]inunitunit number of CDB package input file
[in]ioutunit number of model listing file
[in]namemodelmodel name
[in]paknamepackage name
[in]mempathinput mempath
[in,out]disthe pointer to the discretization
[in]dfwthe pointer to the dfw package
[in]cxsthe pointer to the cxs package

Definition at line 78 of file swf-evp.f90.

80  ! dummy
81  class(BndType), pointer :: packobj !< pointer to default package type
82  integer(I4B), intent(in) :: id !< package id
83  integer(I4B), intent(in) :: ibcnum !< boundary condition number
84  integer(I4B), intent(in) :: inunit !< unit number of CDB package input file
85  integer(I4B), intent(in) :: iout !< unit number of model listing file
86  character(len=*), intent(in) :: namemodel !< model name
87  character(len=*), intent(in) :: pakname !< package name
88  character(len=*), intent(in) :: mempath !< input mempath
89  class(DisBaseType), pointer, intent(inout) :: dis !< the pointer to the discretization
90  type(SwfDfwType), pointer, intent(in) :: dfw !< the pointer to the dfw package
91  type(SwfCxsType), pointer, intent(in) :: cxs !< the pointer to the cxs package
92  ! local
93  type(SwfEvpType), pointer :: evpobj
94 
95  ! allocate evaporation object and scalar variables
96  allocate (evpobj)
97  packobj => evpobj
98 
99  ! create name and memory path
100  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
101  packobj%text = text
102 
103  ! allocate scalars
104  call evpobj%evp_allocate_scalars()
105 
106  ! initialize package
107  call packobj%pack_initialize()
108 
109  packobj%inunit = inunit
110  packobj%iout = iout
111  packobj%id = id
112  packobj%ibcnum = ibcnum
113  packobj%ncolbnd = 1
114  packobj%ictMemPath = create_mem_path(namemodel, 'DFW')
115 
116  ! store pointer to dis
117  evpobj%dis => dis
118 
119  ! store pointer to dfw
120  evpobj%dfw => dfw
121 
122  ! store pointer to cxs
123  evpobj%cxs => cxs
Here is the call graph for this function:
Here is the caller graph for this function:

◆ evp_da()

subroutine swfevpmodule::evp_da ( class(swfevptype this)
private

Definition at line 499 of file swf-evp.f90.

500  ! modules
501  ! dummy
502  class(SwfEvpType) :: this
503 
504  ! Deallocate parent package
505  call this%BndExtType%bnd_da()
506 
507  ! scalars
508  call mem_deallocate(this%iflowred)
509  call mem_deallocate(this%reduction_depth)
510  deallocate (this%read_as_arrays)
511 
512  ! arrays
513  call mem_deallocate(this%evaporation, 'EVAPORATION', this%memoryPath)
514 
515  ! pointers
516  nullify (this%dis)
517  nullify (this%dfw)
518  nullify (this%cxs)

◆ evp_define_listlabel()

subroutine swfevpmodule::evp_define_listlabel ( class(swfevptype), intent(inout)  this)
private

Definition at line 524 of file swf-evp.f90.

525  ! dummy
526  class(SwfEvpType), intent(inout) :: this
527  !
528  ! create the header list label
529  this%listlabel = trim(this%filtyp)//' NO.'
530  if (this%dis%ndim == 3) then
531  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
532  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
533  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
534  elseif (this%dis%ndim == 2) then
535  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
536  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
537  else
538  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
539  end if
540  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'EVAPORATION'
541 ! if(this%multindex > 0) &
542 ! write(this%listlabel, '(a, a16)') trim(this%listlabel), 'MULTIPLIER'
543  if (this%inamedbound == 1) then
544  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
545  end if

◆ evp_df_obs()

subroutine swfevpmodule::evp_df_obs ( class(swfevptype this)
private

Store observation type supported by EVP package. Overrides BndTypebnd_df_obs

Definition at line 587 of file swf-evp.f90.

588  implicit none
589  ! dummy
590  class(SwfEvpType) :: this
591  ! local
592  integer(I4B) :: indx
593 
594  call this%obs%StoreObsType('evp', .true., indx)
595  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
Here is the call graph for this function:

◆ evp_fc()

subroutine swfevpmodule::evp_fc ( class(swfevptype 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 477 of file swf-evp.f90.

478  ! dummy
479  class(SwfEvpType) :: this
480  real(DP), dimension(:), intent(inout) :: rhs
481  integer(I4B), dimension(:), intent(in) :: ia
482  integer(I4B), dimension(:), intent(in) :: idxglo
483  class(MatrixBaseType), pointer :: matrix_sln
484  ! local
485  integer(I4B) :: i, n, ipos
486 
487  ! Copy package rhs and hcof into solution rhs and amat
488  do i = 1, this%nbound
489  n = this%nodelist(i)
490  if (n <= 0) cycle
491  rhs(n) = rhs(n) + this%rhs(i)
492  ipos = ia(n)
493  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
494  end do

◆ evp_obs_supported()

logical function swfevpmodule::evp_obs_supported ( class(swfevptype this)
private

Definition at line 575 of file swf-evp.f90.

576  implicit none
577  ! dummy
578  class(SwfEvpType) :: this
579  evp_obs_supported = .true.

◆ evp_read_initial_attr()

subroutine swfevpmodule::evp_read_initial_attr ( class(swfevptype), intent(inout)  this)
private

Definition at line 247 of file swf-evp.f90.

248  ! dummy
249  class(SwfEvpType), intent(inout) :: this
250 
251  if (this%read_as_arrays) then
252  call this%default_nodelist()
253  end if

◆ evp_rp()

subroutine swfevpmodule::evp_rp ( class(swfevptype), intent(inout)  this)
private

Read itmp and read new boundaries if itmp > 0

Definition at line 260 of file swf-evp.f90.

261  ! modules
262  use tdismodule, only: kper
263  implicit none
264  ! dummy
265  class(SwfEvpType), intent(inout) :: this
266 
267  if (this%iper /= kper) return
268 
269  if (this%read_as_arrays) then
270  ! no need to do anything because this%evaporation points directly to
271  ! the input context evaporation, which is automatically updated by idm
272  else
273  call this%BndExtType%bnd_rp()
274  end if
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23

◆ evp_source_dimensions()

subroutine swfevpmodule::evp_source_dimensions ( class(swfevptype), intent(inout)  this)
private

Definition at line 215 of file swf-evp.f90.

216  ! dummy
217  class(SwfEvpType), intent(inout) :: this
218 
219  if (this%read_as_arrays) then
220 
221  ! Set maxbound to the number of cells per layer, which is simply
222  ! nrow * ncol for a dis2d grid, and nodesuser for disv2d and disv1d
223  this%maxbound = this%dis%get_ncpl()
224 
225  ! verify dimensions were set
226  if (this%maxbound <= 0) then
227  write (errmsg, '(a)') &
228  'MAXBOUND must be an integer greater than zero.'
229  call store_error(errmsg)
230  call store_error_filename(this%input_fname)
231  end if
232 
233  else
234 
235  ! source maxbound
236  call this%BndExtType%source_dimensions()
237 
238  end if
239 
240  ! Call define_listlabel to construct the list label that is written
241  ! when PRINT_INPUT option is used.
242  call this%define_listlabel()
Here is the call graph for this function:

◆ evp_source_options()

subroutine swfevpmodule::evp_source_options ( class(swfevptype), intent(inout)  this)

Definition at line 169 of file swf-evp.f90.

170  ! modules
172  implicit none
173  ! dummy
174  class(SwfEvpType), intent(inout) :: this
175  ! local
176  logical(LGP) :: found_readasarrays = .false.
177 
178  ! source common bound options
179  call this%BndExtType%source_options()
180 
181  ! update defaults with idm sourced values
182  call mem_set_value(this%read_as_arrays, 'READASARRAYS', this%input_mempath, &
183  found_readasarrays)
184 
185  ! log evp params
186  call this%log_evp_options(found_readasarrays)

◆ get_evap_reduce_mult()

real(dp) function swfevpmodule::get_evap_reduce_mult ( class(swfevptype this,
real(dp), intent(in)  stage,
real(dp), intent(in)  bottom 
)
private

Definition at line 457 of file swf-evp.f90.

458  ! dummy
459  class(SwfEvpType) :: this
460  real(DP), intent(in) :: stage
461  real(DP), intent(in) :: bottom
462  ! return
463  real(DP) :: qmult
464  ! local
465  real(DP) :: tp
466 
467  qmult = done
468  if (this%iflowred == 1) then
469  tp = bottom + this%reduction_depth
470  qmult = sqsaturation(tp, bottom, stage)
471  end if
472 
Here is the call graph for this function:

◆ get_qevp()

real(dp) function swfevpmodule::get_qevp ( class(swfevptype this,
integer(i4b), intent(in)  node,
real(dp), intent(in)  rlen,
real(dp), intent(in)  snew,
real(dp), intent(in)  sold,
real(dp), intent(in)  evaporation 
)
private

Calculate qevp for both channel and overland flow grids. Approximate the average water surface width of the channel as wavg = delta A over delta h, and then multiply wavg by reach length to come up with surface water area for the channel. Reduce evaporation when depths are small and shut it off when there is no water in the cell.

Parameters
thisthis instance
[in]nodereduced node number
[in]rlenlength of reach
[in]snewcurrent stage in reach
[in]soldprevious stage in reach
[in]evaporationevaporation rate in length per time

Definition at line 399 of file swf-evp.f90.

400  ! dummy
401  class(SwfEvpType) :: this !< this instance
402  integer(I4B), intent(in) :: node !< reduced node number
403  real(DP), intent(in) :: rlen !< length of reach
404  real(DP), intent(in) :: snew !< current stage in reach
405  real(DP), intent(in) :: sold !< previous stage in reach
406  real(DP), intent(in) :: evaporation !< evaporation rate in length per time
407  ! return
408  real(DP) :: qevp
409  ! local
410  integer(I4B) :: idcxs
411  real(DP) :: depth
412  real(DP) :: bt
413  real(DP) :: area
414  real(DP) :: anew
415  real(DP) :: aold
416  real(DP) :: denom
417  real(DP) :: wavg
418  real(DP) :: width_channel
419  real(DP) :: dummy
420  real(DP) :: qmult
421 
422  ! calculate depth of water
423  bt = this%dis%bot(node)
424 
425  ! Determine the water surface area
426  if (this%dis%is_2d()) then
427  ! overland flow case
428  area = this%dis%get_area(node)
429  else if (this%dis%is_1d()) then
430  ! channel case
431  idcxs = this%dfw%idcxs(node)
432  call this%dis%get_flow_width(node, node, 0, width_channel, dummy)
433 
434  depth = snew - bt
435  anew = this%cxs%get_area(idcxs, width_channel, depth)
436  depth = sold - bt
437  aold = this%cxs%get_area(idcxs, width_channel, depth)
438  wavg = this%cxs%get_wetted_top_width(idcxs, width_channel, depth)
439  denom = snew - sold
440  if (abs(denom) > dprec) then
441  wavg = (anew - aold) / (snew - sold)
442  end if
443  area = rlen * wavg
444 
445  end if
446 
447  ! Reduce the evap rate as cell goes dry
448  qmult = this%get_evap_reduce_mult(snew, bt)
449 
450  ! calculate volumetric evaporation flow in L^3/T
451  qevp = evaporation * area * qmult
452 

◆ log_evp_options()

subroutine swfevpmodule::log_evp_options ( class(swfevptype), intent(inout)  this,
logical(lgp), intent(in)  found_readasarrays 
)

Definition at line 191 of file swf-evp.f90.

192  implicit none
193  ! dummy
194  class(SwfEvpType), intent(inout) :: this
195  logical(LGP), intent(in) :: found_readasarrays
196  ! formats
197  character(len=*), parameter :: fmtreadasarrays = &
198  &"(4x, 'EVAPORATION INPUT WILL BE READ AS ARRAY(S).')"
199 
200  ! log found options
201  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
202  //' OPTIONS'
203 
204  if (found_readasarrays) then
205  write (this%iout, fmtreadasarrays)
206  end if
207 
208  ! close logging block
209  write (this%iout, '(1x,a)') &
210  'END OF '//trim(adjustl(this%text))//' OPTIONS'

◆ reach_length_pointer()

real(dp) function, dimension(:), pointer swfevpmodule::reach_length_pointer ( class(swfevptype this)
Parameters
thisthis instance

Definition at line 625 of file swf-evp.f90.

626  ! dummy
627  class(SwfEvpType) :: this !< this instance
628  ! return
629  real(DP), dimension(:), pointer :: ptr
630  ! local
631  class(DisBaseType), pointer :: dis
632 
633  ptr => null()
634  dis => this%dis
635  select type (dis)
636  type is (disv1dtype)
637  ptr => dis%length
638  end select
639 

Variable Documentation

◆ ftype

character(len=lenftype) swfevpmodule::ftype = 'EVP'
private

Definition at line 35 of file swf-evp.f90.

35  character(len=LENFTYPE) :: ftype = 'EVP'

◆ text

character(len=lenpackagename) swfevpmodule::text = ' EVP'
private

Definition at line 36 of file swf-evp.f90.

36  character(len=LENPACKAGENAME) :: text = ' EVP'