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

This module contains the WEL package methods. More...

Data Types

type  weltype
 

Functions/Subroutines

subroutine, public wel_create (packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
 @ brief Create a new package object More...
 
subroutine wel_da (this)
 @ brief Deallocate package memory More...
 
subroutine wel_allocate_scalars (this)
 @ brief Allocate scalars More...
 
subroutine wel_allocate_arrays (this, nodelist, auxvar)
 @ brief Allocate arrays More...
 
subroutine wel_options (this)
 @ brief Source additional options for package More...
 
subroutine log_wel_options (this, found)
 @ brief Log WEL specific package options More...
 
subroutine wel_cf (this)
 @ brief Formulate the package hcof and rhs terms. More...
 
subroutine wel_fc (this, rhs, ia, idxglo, matrix_sln)
 @ brief Copy hcof and rhs terms into solution. More...
 
subroutine wel_fn (this, rhs, ia, idxglo, matrix_sln)
 @ brief Add Newton-Raphson terms for package into solution. More...
 
subroutine wel_afr_csv_init (this, fname)
 Initialize the auto flow reduce csv output file. More...
 
subroutine wel_afr_csv_write (this)
 Write out auto flow reductions only when & where they occur. More...
 
subroutine define_listlabel (this)
 @ brief Define the list label for the package More...
 
logical function wel_obs_supported (this)
 Determine if observations are supported. More...
 
subroutine wel_df_obs (this)
 Define the observation types available in the package. More...
 
subroutine wel_bd_obs (this)
 Save observations for the package. More...
 
real(dp) function q_mult (this, row)
 
real(dp) function wel_bound_value (this, col, row)
 @ brief Return a bound value More...
 

Variables

character(len=lenftype) ftype = 'WEL'
 package ftype More...
 
character(len=16) text = ' WEL'
 package flow text string More...
 

Detailed Description

This module contains the overridden methods for the standard WEL package. Several methods need to be overridden because of the AUTO_FLOW_REDUCE option. Overridden methods include:

  • bnd_cf (AUTO_FLOW_REDUCE)
  • bnd_fc (AUTO_FLOW_REDUCE)
  • bnd_fn (AUTO_FLOW_REDUCE Newton-Raphson terms)
  • bnd_ot_package_flows (write AUTO_FLOW_REDUCE terms to csv file)
  • bnd_da (deallocate AUTO_FLOW_REDUCE variables)
  • bnd_bd_obs (wel-reduction observation added)

Function/Subroutine Documentation

◆ define_listlabel()

subroutine welmodule::define_listlabel ( class(weltype), intent(inout)  this)

Method defined the list label for the WEL package. The list label is the heading that is written to iout when PRINT_INPUT option is used.

Parameters
[in,out]thisWelType object

Definition at line 495 of file gwf-wel.f90.

496  ! -- dummy variables
497  class(WelType), intent(inout) :: this !< WelType object
498  !
499  ! -- create the header list label
500  this%listlabel = trim(this%filtyp)//' NO.'
501  if (this%dis%ndim == 3) then
502  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
503  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
504  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
505  elseif (this%dis%ndim == 2) then
506  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
507  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
508  else
509  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
510  end if
511  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'STRESS RATE'
512  if (this%inamedbound == 1) then
513  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
514  end if

◆ log_wel_options()

subroutine welmodule::log_wel_options ( class(weltype), intent(inout)  this,
type(gwfwelparamfoundtype), intent(in)  found 
)

Definition at line 256 of file gwf-wel.f90.

257  ! -- modules
259  ! -- dummy variables
260  class(WelType), intent(inout) :: this
261  type(GwfWelParamFoundType), intent(in) :: found
262  ! -- local variables
263  ! -- format
264  character(len=*), parameter :: fmtflowred = &
265  &"(4x, 'AUTOMATIC FLOW REDUCTION OF WELLS IMPLEMENTED.')"
266  character(len=*), parameter :: fmtflowredv = &
267  &"(4x, 'AUTOMATIC FLOW REDUCTION FRACTION (',g15.7,').')"
268  character(len=*), parameter :: fmtflowredl = &
269  &"(4x, 'AUTOMATIC FLOW REDUCTION LENGTH (',g15.7,').')"
270  !
271  ! -- log found options
272  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
273  //' OPTIONS'
274 
275  if (found%iflowredlen) then
276  write (this%iout, fmtflowred)
277  write (this%iout, '(4x,A)') &
278  'AUTOMATIC FLOW REDUCTION FRACTION INTERPRETED AS A LENGTH'
279  end if
280 
281  if (found%flowred) then
282  if (this%iflowredlen == 0) then
283  write (this%iout, fmtflowredv) this%flowred
284  else
285  write (this%iout, fmtflowredl) this%flowred
286  end if
287  end if
288  !
289  if (found%afrcsvfile) then
290  ! -- currently no-op
291  end if
292 
293  if (found%mover) then
294  write (this%iout, '(4x,A)') 'MOVER OPTION ENABLED'
295  end if
296  !
297  ! -- close logging block
298  write (this%iout, '(1x,a)') &
299  'END OF '//trim(adjustl(this%text))//' OPTIONS'

◆ q_mult()

real(dp) function welmodule::q_mult ( class(weltype), intent(inout)  this,
integer(i4b), intent(in)  row 
)
private

Definition at line 617 of file gwf-wel.f90.

618  ! -- modules
619  use constantsmodule, only: dzero
620  ! -- dummy variables
621  class(WelType), intent(inout) :: this
622  integer(I4B), intent(in) :: row
623  ! -- result
624  real(DP) :: q
625  !
626  if (this%iauxmultcol > 0) then
627  q = this%q(row) * this%auxvar(this%iauxmultcol, row)
628  else
629  q = this%q(row)
630  end if
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65

◆ wel_afr_csv_init()

subroutine welmodule::wel_afr_csv_init ( class(weltype), intent(inout)  this,
character(len=*), intent(in)  fname 
)
private
Parameters
[in,out]thisWelType object

Definition at line 442 of file gwf-wel.f90.

443  ! -- dummy variables
444  class(WelType), intent(inout) :: this !< WelType object
445  character(len=*), intent(in) :: fname
446  ! -- format
447  character(len=*), parameter :: fmtafrcsv = &
448  "(4x, 'AUTO FLOW REDUCE INFORMATION WILL BE SAVED TO FILE: ', a, /4x, &
449  &'OPENED ON UNIT: ', I0)"
450 
451  this%ioutafrcsv = getunit()
452  call openfile(this%ioutafrcsv, this%iout, fname, 'CSV', &
453  filstat_opt='REPLACE')
454  write (this%iout, fmtafrcsv) trim(adjustl(fname)), &
455  this%ioutafrcsv
456  write (this%ioutafrcsv, '(a)') &
457  'time,period,step,boundnumber,cellnumber,rate-requested,&
458  &rate-actual,wel-reduction'
Here is the call graph for this function:

◆ wel_afr_csv_write()

subroutine welmodule::wel_afr_csv_write ( class(weltype), intent(inout)  this)
private
Parameters
[in,out]thisWelType object

Definition at line 462 of file gwf-wel.f90.

463  ! -- modules
464  use tdismodule, only: totim, kstp, kper
465  ! -- dummy variables
466  class(WelType), intent(inout) :: this !< WelType object
467  ! -- local
468  integer(I4B) :: i
469  integer(I4B) :: nodereduced
470  integer(I4B) :: nodeuser
471  real(DP) :: v
472  ! -- format
473  do i = 1, this%nbound
474  nodereduced = this%nodelist(i)
475  !
476  ! -- test if node is constant or inactive
477  if (this%ibound(nodereduced) <= 0) then
478  cycle
479  end if
480  v = this%q_mult(i) + this%rhs(i)
481  if (v < dzero) then
482  nodeuser = this%dis%get_nodeuser(nodereduced)
483  write (this%ioutafrcsv, '(*(G0,:,","))') &
484  totim, kper, kstp, i, nodeuser, this%q_mult(i), this%simvals(i), v
485  end if
486  end do
real(dp), pointer, public totim
time relative to start of simulation
Definition: tdis.f90:32
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23

◆ wel_allocate_arrays()

subroutine welmodule::wel_allocate_arrays ( class(weltype this,
integer(i4b), dimension(:), optional, pointer, contiguous  nodelist,
real(dp), dimension(:, :), optional, pointer, contiguous  auxvar 
)

Allocate and initialize arrays for the WEL package

Definition at line 164 of file gwf-wel.f90.

165  ! -- modules
167  ! -- dummy
168  class(WelType) :: this
169  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
170  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
171  ! -- local
172  !
173  ! -- call BndExtType allocate scalars
174  call this%BndExtType%allocate_arrays(nodelist, auxvar)
175  !
176  ! -- set constant head array input context pointer
177  call mem_setptr(this%q, 'Q', this%input_mempath)
178  !
179  ! -- checkin constant head array input context pointer
180  call mem_checkin(this%q, 'Q', this%memoryPath, &
181  'Q', this%input_mempath)

◆ wel_allocate_scalars()

subroutine welmodule::wel_allocate_scalars ( class(weltype this)

Allocate and initialize scalars for the WEL package. The base model allocate scalars method is also called.

Parameters
thisWelType object

Definition at line 137 of file gwf-wel.f90.

138  ! -- modules
140  ! -- dummy variables
141  class(WelType) :: this !< WelType object
142  !
143  ! -- call base type allocate scalars
144  call this%BndExtType%allocate_scalars()
145  !
146  ! -- allocate the object and assign values to object variables
147  call mem_allocate(this%iflowred, 'IFLOWRED', this%memoryPath)
148  call mem_allocate(this%flowred, 'FLOWRED', this%memoryPath)
149  call mem_allocate(this%ioutafrcsv, 'IOUTAFRCSV', this%memoryPath)
150  call mem_allocate(this%iflowredlen, 'IFLOWREDLEN', this%memoryPath)
151  !
152  ! -- Set values
153  this%iflowred = 0
154  this%ioutafrcsv = 0
155  this%flowred = dzero
156  this%iflowredlen = 0

◆ wel_bd_obs()

subroutine welmodule::wel_bd_obs ( class(weltype this)
private

Method to save simulated values for the WEL package.

Parameters
thisWelType object

Definition at line 566 of file gwf-wel.f90.

567  ! -- dummy variables
568  class(WelType) :: this !< WelType object
569  ! -- local variables
570  integer(I4B) :: i
571  integer(I4B) :: n
572  integer(I4B) :: jj
573  real(DP) :: v
574  type(ObserveType), pointer :: obsrv => null()
575  !
576  ! -- clear the observations
577  call this%obs%obs_bd_clear()
578  !
579  ! -- Save simulated values for all of package's observations.
580  do i = 1, this%obs%npakobs
581  obsrv => this%obs%pakobs(i)%obsrv
582  if (obsrv%BndFound) then
583  do n = 1, obsrv%indxbnds_count
584  v = dnodata
585  jj = obsrv%indxbnds(n)
586  select case (obsrv%ObsTypeId)
587  case ('TO-MVR')
588  if (this%imover == 1) then
589  v = this%pakmvrobj%get_qtomvr(jj)
590  if (v > dzero) then
591  v = -v
592  end if
593  end if
594  case ('WEL')
595  v = this%simvals(jj)
596  case ('WEL-REDUCTION')
597  if (this%iflowred > 0) then
598  v = this%q_mult(jj) + this%rhs(jj)
599  end if
600  case default
601  errmsg = 'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
602  call store_error(errmsg)
603  end select
604  call this%obs%SaveOneSimval(obsrv, v)
605  end do
606  else
607  call this%obs%SaveOneSimval(obsrv, dnodata)
608  end if
609  end do
610  !
611  ! -- Write the auto flow reduce csv file entries for this step
612  if (this%ioutafrcsv > 0) then
613  call this%wel_afr_csv_write()
614  end if
Here is the call graph for this function:

◆ wel_bound_value()

real(dp) function welmodule::wel_bound_value ( class(weltype), intent(inout)  this,
integer(i4b), intent(in)  col,
integer(i4b), intent(in)  row 
)

Return a bound value associated with an ncolbnd index and row.

Definition at line 639 of file gwf-wel.f90.

640  ! -- modules
641  use constantsmodule, only: dzero
642  ! -- dummy variables
643  class(WelType), intent(inout) :: this
644  integer(I4B), intent(in) :: col
645  integer(I4B), intent(in) :: row
646  ! -- result
647  real(DP) :: bndval
648  !
649  select case (col)
650  case (1)
651  bndval = this%q_mult(row)
652  case default
653  errmsg = 'Programming error. WEL bound value requested column '&
654  &'outside range of ncolbnd (1).'
655  call store_error(errmsg)
656  call store_error_filename(this%input_fname)
657  end select
Here is the call graph for this function:

◆ wel_cf()

subroutine welmodule::wel_cf ( class(weltype this)

Formulate the hcof and rhs terms for the WEL package that will be added to the coefficient matrix and right-hand side vector.

Parameters
thisWelType object

Definition at line 308 of file gwf-wel.f90.

309  ! -- dummy variables
310  class(WelType) :: this !< WelType object
311  ! -- local variables
312  integer(I4B) :: i, node, ict
313  real(DP) :: qmult
314  real(DP) :: q
315  real(DP) :: tp
316  real(DP) :: bt
317  real(DP) :: thick
318  !
319  ! -- Return if no wells
320  if (this%nbound == 0) return
321  !
322  ! -- Calculate hcof and rhs for each well entry
323  do i = 1, this%nbound
324  node = this%nodelist(i)
325  this%hcof(i) = dzero
326  if (this%ibound(node) <= 0) then
327  this%rhs(i) = dzero
328  cycle
329  end if
330  q = this%q_mult(i)
331  if (this%iflowred /= 0 .and. q < dzero) then
332  ict = this%icelltype(node)
333  if (ict /= 0) then
334  bt = this%dis%bot(node)
335  if (this%iflowredlen == 0) then
336  thick = this%dis%top(node) - bt
337  else
338  thick = done
339  end if
340  tp = bt + this%flowred * thick
341  qmult = sqsaturation(tp, bt, this%xnew(node))
342  q = q * qmult
343  end if
344  end if
345  this%rhs(i) = -q
346  end do
Here is the call graph for this function:

◆ wel_create()

subroutine, public welmodule::wel_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 
)

Create a new WEL Package object

Parameters
packobjpointer to default package type
[in]idpackage id
[in]ibcnumboundary condition number
[in]inunitunit number of WEL package input file
[in]ioutunit number of model listing file
[in]namemodelmodel name
[in]paknamepackage name
[in]mempathinput mempath

Definition at line 73 of file gwf-wel.f90.

75  ! -- dummy variables
76  class(BndType), pointer :: packobj !< pointer to default package type
77  integer(I4B), intent(in) :: id !< package id
78  integer(I4B), intent(in) :: ibcnum !< boundary condition number
79  integer(I4B), intent(in) :: inunit !< unit number of WEL package input file
80  integer(I4B), intent(in) :: iout !< unit number of model listing file
81  character(len=*), intent(in) :: namemodel !< model name
82  character(len=*), intent(in) :: pakname !< package name
83  character(len=*), intent(in) :: mempath !< input mempath
84  ! -- local variables
85  type(WelType), pointer :: welobj
86  !
87  ! -- allocate the object and assign values to object variables
88  allocate (welobj)
89  packobj => welobj
90  !
91  ! -- create name and memory path
92  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
93  packobj%text = text
94  !
95  ! -- allocate scalars
96  call welobj%allocate_scalars()
97  !
98  ! -- initialize package
99  call packobj%pack_initialize()
100 
101  packobj%inunit = inunit
102  packobj%iout = iout
103  packobj%id = id
104  packobj%ibcnum = ibcnum
105  packobj%ncolbnd = 1
106  packobj%ictMemPath = create_mem_path(namemodel, 'NPF')
Here is the call graph for this function:
Here is the caller graph for this function:

◆ wel_da()

subroutine welmodule::wel_da ( class(weltype this)
private

Deallocate WEL package scalars and arrays.

Parameters
thisWelType object

Definition at line 114 of file gwf-wel.f90.

115  ! -- modules
117  ! -- dummy variables
118  class(WelType) :: this !< WelType object
119  !
120  ! -- Deallocate parent package
121  call this%BndExtType%bnd_da()
122  !
123  ! -- scalars
124  call mem_deallocate(this%iflowred)
125  call mem_deallocate(this%flowred)
126  call mem_deallocate(this%ioutafrcsv)
127  call mem_deallocate(this%iflowredlen)
128  call mem_deallocate(this%q, 'Q', this%memoryPath)

◆ wel_df_obs()

subroutine welmodule::wel_df_obs ( class(weltype this)
private

Method to define the observation types available in the WEL package.

Parameters
thisWelType object

Definition at line 540 of file gwf-wel.f90.

541  ! -- dummy variables
542  class(WelType) :: this !< WelType object
543  ! -- local variables
544  integer(I4B) :: indx
545  !
546  ! -- initialize observations
547  call this%obs%StoreObsType('wel', .true., indx)
548  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
549  !
550  ! -- Store obs type and assign procedure pointer
551  ! for to-mvr observation type.
552  call this%obs%StoreObsType('to-mvr', .true., indx)
553  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
554  !
555  ! -- Store obs type and assign procedure pointer
556  ! for wel-reduction observation type.
557  call this%obs%StoreObsType('wel-reduction', .true., indx)
558  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
Here is the call graph for this function:

◆ wel_fc()

subroutine welmodule::wel_fc ( class(weltype 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

Add the hcof and rhs terms for the WEL package to the coefficient matrix and right-hand side vector.

Parameters
thisWelType object
[in,out]rhsright-hand side vector for model
[in]iasolution CRS row pointers
[in]idxglomapping vector for model (local) to solution (global)
matrix_slnsolution coefficient matrix

Definition at line 355 of file gwf-wel.f90.

356  ! -- dummy variables
357  class(WelType) :: this !< WelType object
358  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
359  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
360  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
361  class(MatrixBaseType), pointer :: matrix_sln !< solution coefficient matrix
362  ! -- local variables
363  integer(I4B) :: i
364  integer(I4B) :: n
365  integer(I4B) :: ipos
366  !
367  ! -- pakmvrobj fc
368  if (this%imover == 1) then
369  call this%pakmvrobj%fc()
370  end if
371  !
372  ! -- Copy package rhs and hcof into solution rhs and amat
373  do i = 1, this%nbound
374  n = this%nodelist(i)
375  rhs(n) = rhs(n) + this%rhs(i)
376  ipos = ia(n)
377  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
378  !
379  ! -- If mover is active and this well is discharging,
380  ! store available water (as positive value).
381  if (this%imover == 1 .and. this%rhs(i) > dzero) then
382  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
383  end if
384  end do

◆ wel_fn()

subroutine welmodule::wel_fn ( class(weltype 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

Calculate and add the Newton-Raphson terms for the WEL package to the coefficient matrix and right-hand side vector.

Parameters
thisWelType object
[in,out]rhsright-hand side vector for model
[in]iasolution CRS row pointers
[in]idxglomapping vector for model (local) to solution (global)
matrix_slnsolution coefficient matrix

Definition at line 393 of file gwf-wel.f90.

394  ! -- dummy variables
395  class(WelType) :: this !< WelType object
396  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
397  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
398  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
399  class(MatrixBaseType), pointer :: matrix_sln !< solution coefficient matrix
400  ! -- local variables
401  integer(I4B) :: i
402  integer(I4B) :: node
403  integer(I4B) :: ipos
404  integer(I4B) :: ict
405  real(DP) :: drterm
406  real(DP) :: q
407  real(DP) :: tp
408  real(DP) :: bt
409  real(DP) :: thick
410  !
411  ! -- Copy package rhs and hcof into solution rhs and amat
412  do i = 1, this%nbound
413  node = this%nodelist(i)
414  !
415  ! -- test if node is constant or inactive
416  if (this%ibound(node) <= 0) then
417  cycle
418  end if
419  !
420  ! -- well rate is possibly head dependent
421  ict = this%icelltype(node)
422  if (this%iflowred /= 0 .and. ict /= 0) then
423  ipos = ia(node)
424  q = -this%rhs(i)
425  if (q < dzero) then
426  ! -- calculate derivative for well
427  tp = this%dis%top(node)
428  bt = this%dis%bot(node)
429  thick = tp - bt
430  tp = bt + this%flowred * thick
431  drterm = sqsaturationderivative(tp, bt, this%xnew(node))
432  drterm = drterm * this%q_mult(i)
433  !--fill amat and rhs with newton-raphson terms
434  call matrix_sln%add_value_pos(idxglo(ipos), drterm)
435  rhs(node) = rhs(node) + drterm * this%xnew(node)
436  end if
437  end if
438  end do
Here is the call graph for this function:

◆ wel_obs_supported()

logical function welmodule::wel_obs_supported ( class(weltype this)
private

Function to determine if observations are supported by the WEL package. Observations are supported by the WEL package.

Returns
wel_obs_supported boolean indicating if observations are supported
Parameters
thisWelType object

Definition at line 527 of file gwf-wel.f90.

528  ! -- dummy variables
529  class(WelType) :: this !< WelType object
530  !
531  ! -- set boolean
532  wel_obs_supported = .true.

◆ wel_options()

subroutine welmodule::wel_options ( class(weltype), intent(inout)  this)

Source additional options for WEL package.

Parameters
[in,out]thisWelType object

Definition at line 189 of file gwf-wel.f90.

190  ! -- modules
191  use inputoutputmodule, only: urword
194  ! -- dummy variables
195  class(WelType), intent(inout) :: this !< WelType object
196  ! -- local variables
197  character(len=LINELENGTH) :: fname
198  type(GwfWelParamFoundType) :: found
199  ! -- formats
200  character(len=*), parameter :: fmtflowred = &
201  &"(4x, 'AUTOMATIC FLOW REDUCTION OF WELLS IMPLEMENTED.')"
202  character(len=*), parameter :: fmtflowredv = &
203  &"(4x, 'AUTOMATIC FLOW REDUCTION FRACTION (',g15.7,').')"
204  !
205  ! -- source base BndExtType options
206  call this%BndExtType%source_options()
207  !
208  ! -- source well options from input context
209  call mem_set_value(this%flowred, 'FLOWRED', this%input_mempath, found%flowred)
210  call mem_set_value(fname, 'AFRCSVFILE', this%input_mempath, found%afrcsvfile)
211  call mem_set_value(this%imover, 'MOVER', this%input_mempath, found%mover)
212  call mem_set_value(this%iflowredlen, 'IFLOWREDLEN', this%input_mempath, &
213  found%iflowredlen)
214 
215  if (found%iflowredlen) then
216  if (found%flowred .eqv. .false.) then
217  write (warnmsg, '(a)') &
218  'FLOW_REDUCTION_LENGTH option specified but a AUTO_FLOW_REDUCTION value &
219  &is not specified. The FLOW_REDUCTION_LENGTH option will be ignored.'
220  call store_warning(warnmsg)
221  else
222  this%iflowredlen = 1
223  end if
224  end if
225 
226  if (found%flowred) then
227  this%iflowred = 1
228  if (this%flowred <= dzero) then
229  if (found%iflowredlen) then
230  write (errmsg, '(a)') &
231  'An AUTO_FLOW_REDUCTION value less than or equal to zero cannot be &
232  &specified if the FLOW_REDUCTION_LENGTH option is specified.'
233  call store_error(errmsg)
234  else
235  this%flowred = dem1
236  end if
237  else if (this%flowred > done .and. this%iflowredlen == 0) then
238  this%flowred = done
239  end if
240  end if
241 
242  if (found%afrcsvfile) then
243  call this%wel_afr_csv_init(fname)
244  end if
245 
246  if (found%mover) then
247  this%imover = 1
248  end if
249 
250  ! -- log WEL specific options
251  call this%log_wel_options(found)
subroutine, public urword(line, icol, istart, istop, ncode, n, r, iout, in)
Extract a word from a string.
Here is the call graph for this function:

Variable Documentation

◆ ftype

character(len=lenftype) welmodule::ftype = 'WEL'
private

Definition at line 36 of file gwf-wel.f90.

36  character(len=LENFTYPE) :: ftype = 'WEL' !< package ftype

◆ text

character(len=16) welmodule::text = ' WEL'
private

Definition at line 37 of file gwf-wel.f90.

37  character(len=16) :: text = ' WEL' !< package flow text string