MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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_rp (this)
 @ brief WEL read and prepare 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 510 of file gwf-wel.f90.

511  ! -- dummy variables
512  class(WelType), intent(inout) :: this !< WelType object
513  !
514  ! -- create the header list label
515  this%listlabel = trim(this%filtyp)//' NO.'
516  if (this%dis%ndim == 3) then
517  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
518  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
519  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
520  elseif (this%dis%ndim == 2) then
521  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
522  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
523  else
524  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
525  end if
526  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'STRESS RATE'
527  if (this%inamedbound == 1) then
528  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
529  end if
530  !
531  ! -- return
532  return

◆ log_wel_options()

subroutine welmodule::log_wel_options ( class(weltype), intent(inout)  this,
type(gwfwelparamfoundtype), intent(in)  found 
)
Parameters
[in,out]thisBndExtType object

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

250  ! -- modules
252  ! -- dummy variables
253  class(WelType), intent(inout) :: this !< BndExtType object
254  type(GwfWelParamFoundType), intent(in) :: found
255  ! -- local variables
256  ! -- format
257  character(len=*), parameter :: fmtflowred = &
258  &"(4x, 'AUTOMATIC FLOW REDUCTION OF WELLS IMPLEMENTED.')"
259  character(len=*), parameter :: fmtflowredv = &
260  &"(4x, 'AUTOMATIC FLOW REDUCTION FRACTION (',g15.7,').')"
261  !
262  ! -- log found options
263  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
264  //' OPTIONS'
265  !
266  if (found%flowred) then
267  if (this%iflowred > 0) &
268  write (this%iout, fmtflowred)
269  write (this%iout, fmtflowredv) this%flowred
270  end if
271  !
272  if (found%afrcsvfile) then
273  ! -- currently no-op
274  end if
275  !
276  if (found%mover) then
277  write (this%iout, '(4x,A)') 'MOVER OPTION ENABLED'
278  end if
279  !
280  ! -- close logging block
281  write (this%iout, '(1x,a)') &
282  'END OF '//trim(adjustl(this%text))//' OPTIONS'
283  !
284  ! -- return
285  return

◆ q_mult()

real(dp) function welmodule::q_mult ( class(weltype), intent(inout)  this,
integer(i4b), intent(in)  row 
)
private
Parameters
[in,out]thisBndExtType object

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

645  ! -- modules
646  use constantsmodule, only: dzero
647  ! -- dummy variables
648  class(WelType), intent(inout) :: this !< BndExtType object
649  integer(I4B), intent(in) :: row
650  ! -- result
651  real(DP) :: q
652  !
653  if (this%iauxmultcol > 0) then
654  q = this%q(row) * this%auxvar(this%iauxmultcol, row)
655  else
656  q = this%q(row)
657  end if
658  !
659  ! -- return
660  return
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64

◆ 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 456 of file gwf-wel.f90.

457  ! -- dummy variables
458  class(WelType), intent(inout) :: this !< WelType object
459  character(len=*), intent(in) :: fname
460  ! -- format
461  character(len=*), parameter :: fmtafrcsv = &
462  "(4x, 'AUTO FLOW REDUCE INFORMATION WILL BE SAVED TO FILE: ', a, /4x, &
463  &'OPENED ON UNIT: ', I0)"
464 
465  this%ioutafrcsv = getunit()
466  call openfile(this%ioutafrcsv, this%iout, fname, 'CSV', &
467  filstat_opt='REPLACE')
468  write (this%iout, fmtafrcsv) trim(adjustl(fname)), &
469  this%ioutafrcsv
470  write (this%ioutafrcsv, '(a)') &
471  'time,period,step,boundnumber,cellnumber,rate-requested,&
472  &rate-actual,wel-reduction'
473  return
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 477 of file gwf-wel.f90.

478  ! -- modules
479  use tdismodule, only: totim, kstp, kper
480  ! -- dummy variables
481  class(WelType), intent(inout) :: this !< WelType object
482  ! -- local
483  integer(I4B) :: i
484  integer(I4B) :: nodereduced
485  integer(I4B) :: nodeuser
486  real(DP) :: v
487  ! -- format
488  do i = 1, this%nbound
489  nodereduced = this%nodelist(i)
490  !
491  ! -- test if node is constant or inactive
492  if (this%ibound(nodereduced) <= 0) then
493  cycle
494  end if
495  v = this%q_mult(i) + this%rhs(i)
496  if (v < dzero) then
497  nodeuser = this%dis%get_nodeuser(nodereduced)
498  write (this%ioutafrcsv, '(*(G0,:,","))') &
499  totim, kper, kstp, i, nodeuser, this%q_mult(i), this%simvals(i), v
500  end if
501  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 169 of file gwf-wel.f90.

170  ! -- modules
172  ! -- dummy
173  class(WelType) :: this
174  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
175  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
176  ! -- local
177  !
178  ! -- call BndExtType allocate scalars
179  call this%BndExtType%allocate_arrays(nodelist, auxvar)
180  !
181  ! -- set constant head array input context pointer
182  call mem_setptr(this%q, 'Q', this%input_mempath)
183  !
184  ! -- checkin constant head array input context pointer
185  call mem_checkin(this%q, 'Q', this%memoryPath, &
186  'Q', this%input_mempath)
187  !
188  ! -- return
189  return

◆ 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 141 of file gwf-wel.f90.

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

◆ 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 590 of file gwf-wel.f90.

591  ! -- dummy variables
592  class(WelType) :: this !< WelType object
593  ! -- local variables
594  integer(I4B) :: i
595  integer(I4B) :: n
596  integer(I4B) :: jj
597  real(DP) :: v
598  type(ObserveType), pointer :: obsrv => null()
599  !
600  ! -- clear the observations
601  call this%obs%obs_bd_clear()
602  !
603  ! -- Save simulated values for all of package's observations.
604  do i = 1, this%obs%npakobs
605  obsrv => this%obs%pakobs(i)%obsrv
606  if (obsrv%BndFound) then
607  do n = 1, obsrv%indxbnds_count
608  v = dnodata
609  jj = obsrv%indxbnds(n)
610  select case (obsrv%ObsTypeId)
611  case ('TO-MVR')
612  if (this%imover == 1) then
613  v = this%pakmvrobj%get_qtomvr(jj)
614  if (v > dzero) then
615  v = -v
616  end if
617  end if
618  case ('WEL')
619  v = this%simvals(jj)
620  case ('WEL-REDUCTION')
621  if (this%iflowred > 0) then
622  v = this%q_mult(jj) + this%rhs(jj)
623  end if
624  case default
625  errmsg = 'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
626  call store_error(errmsg)
627  end select
628  call this%obs%SaveOneSimval(obsrv, v)
629  end do
630  else
631  call this%obs%SaveOneSimval(obsrv, dnodata)
632  end if
633  end do
634  !
635  ! -- Write the auto flow reduce csv file entries for this step
636  if (this%ioutafrcsv > 0) then
637  call this%wel_afr_csv_write()
638  end if
639  !
640  ! -- return
641  return
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.

Parameters
[in,out]thisBndExtType object

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

670  ! -- modules
671  use constantsmodule, only: dzero
672  ! -- dummy variables
673  class(WelType), intent(inout) :: this !< BndExtType object
674  integer(I4B), intent(in) :: col
675  integer(I4B), intent(in) :: row
676  ! -- result
677  real(DP) :: bndval
678  !
679  select case (col)
680  case (1)
681  bndval = this%q_mult(row)
682  case default
683  errmsg = 'Programming error. WEL bound value requested column '&
684  &'outside range of ncolbnd (1).'
685  call store_error(errmsg)
686  call store_error_filename(this%input_fname)
687  end select
688  !
689  ! -- return
690  return
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 317 of file gwf-wel.f90.

318  ! -- dummy variables
319  class(WelType) :: this !< WelType object
320  ! -- local variables
321  integer(I4B) :: i, node, ict
322  real(DP) :: qmult
323  real(DP) :: q
324  real(DP) :: tp
325  real(DP) :: bt
326  real(DP) :: thick
327  !
328  ! -- Return if no wells
329  if (this%nbound == 0) return
330  !
331  ! -- Calculate hcof and rhs for each well entry
332  do i = 1, this%nbound
333  node = this%nodelist(i)
334  this%hcof(i) = dzero
335  if (this%ibound(node) <= 0) then
336  this%rhs(i) = dzero
337  cycle
338  end if
339  q = this%q_mult(i)
340  if (this%iflowred /= 0 .and. q < dzero) then
341  ict = this%icelltype(node)
342  if (ict /= 0) then
343  tp = this%dis%top(node)
344  bt = this%dis%bot(node)
345  thick = tp - bt
346  tp = bt + this%flowred * thick
347  qmult = sqsaturation(tp, bt, this%xnew(node))
348  q = q * qmult
349  end if
350  end if
351  this%rhs(i) = -q
352  end do
353  !
354  return
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%ictMemPath = create_mem_path(namemodel, 'NPF')
106  !
107  ! -- return
108  return
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 116 of file gwf-wel.f90.

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

◆ 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 561 of file gwf-wel.f90.

562  ! -- dummy variables
563  class(WelType) :: this !< WelType object
564  ! -- local variables
565  integer(I4B) :: indx
566  !
567  ! -- initialize observations
568  call this%obs%StoreObsType('wel', .true., indx)
569  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
570  !
571  ! -- Store obs type and assign procedure pointer
572  ! for to-mvr observation type.
573  call this%obs%StoreObsType('to-mvr', .true., indx)
574  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
575  !
576  ! -- Store obs type and assign procedure pointer
577  ! for wel-reduction observation type.
578  call this%obs%StoreObsType('wel-reduction', .true., indx)
579  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
580  !
581  ! -- return
582  return
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 363 of file gwf-wel.f90.

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

◆ 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 404 of file gwf-wel.f90.

405  ! -- dummy variables
406  class(WelType) :: this !< WelType object
407  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
408  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
409  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
410  class(MatrixBaseType), pointer :: matrix_sln !< solution coefficient matrix
411  ! -- local variables
412  integer(I4B) :: i
413  integer(I4B) :: node
414  integer(I4B) :: ipos
415  integer(I4B) :: ict
416  real(DP) :: drterm
417  real(DP) :: q
418  real(DP) :: tp
419  real(DP) :: bt
420  real(DP) :: thick
421  !
422  ! -- Copy package rhs and hcof into solution rhs and amat
423  do i = 1, this%nbound
424  node = this%nodelist(i)
425  !
426  ! -- test if node is constant or inactive
427  if (this%ibound(node) <= 0) then
428  cycle
429  end if
430  !
431  ! -- well rate is possibly head dependent
432  ict = this%icelltype(node)
433  if (this%iflowred /= 0 .and. ict /= 0) then
434  ipos = ia(node)
435  q = -this%rhs(i)
436  if (q < dzero) then
437  ! -- calculate derivative for well
438  tp = this%dis%top(node)
439  bt = this%dis%bot(node)
440  thick = tp - bt
441  tp = bt + this%flowred * thick
442  drterm = sqsaturationderivative(tp, bt, this%xnew(node))
443  drterm = drterm * this%q_mult(i)
444  !--fill amat and rhs with newton-raphson terms
445  call matrix_sln%add_value_pos(idxglo(ipos), drterm)
446  rhs(node) = rhs(node) + drterm * this%xnew(node)
447  end if
448  end if
449  end do
450  !
451  ! -- return
452  return
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 545 of file gwf-wel.f90.

546  ! -- dummy variables
547  class(WelType) :: this !< WelType object
548  !
549  ! -- set boolean
550  wel_obs_supported = .true.
551  !
552  ! -- return
553  return

◆ 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 197 of file gwf-wel.f90.

198  ! -- modules
199  use inputoutputmodule, only: urword
202  ! -- dummy variables
203  class(WelType), intent(inout) :: this !< WelType object
204  ! -- local variables
205  character(len=LINELENGTH) :: fname
206  type(GwfWelParamFoundType) :: found
207  ! -- formats
208  character(len=*), parameter :: fmtflowred = &
209  &"(4x, 'AUTOMATIC FLOW REDUCTION OF WELLS IMPLEMENTED.')"
210  character(len=*), parameter :: fmtflowredv = &
211  &"(4x, 'AUTOMATIC FLOW REDUCTION FRACTION (',g15.7,').')"
212  !
213  ! -- source base BndExtType options
214  call this%BndExtType%source_options()
215  !
216  ! -- source well options from input context
217  call mem_set_value(this%flowred, 'FLOWRED', this%input_mempath, found%flowred)
218  call mem_set_value(fname, 'AFRCSVFILE', this%input_mempath, found%afrcsvfile)
219  call mem_set_value(this%imover, 'MOVER', this%input_mempath, found%mover)
220  !
221  if (found%flowred) then
222  !
223  this%iflowred = 1
224  !
225  if (this%flowred <= dzero) then
226  this%flowred = dem1
227  else if (this%flowred > done) then
228  this%flowred = done
229  end if
230  end if
231  !
232  if (found%afrcsvfile) then
233  call this%wel_afr_csv_init(fname)
234  end if
235  !
236  if (found%mover) then
237  this%imover = 1
238  end if
239  !
240  ! -- log WEL specific options
241  call this%log_wel_options(found)
242  !
243  ! -- return
244  return
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:

◆ wel_rp()

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

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

292  use tdismodule, only: kper
293  ! -- dummy
294  class(WelType), intent(inout) :: this
295  ! -- local
296  !
297  if (this%iper /= kper) return
298  !
299  ! -- Call the parent class read and prepare
300  call this%BndExtType%bnd_rp()
301  !
302  ! -- Write the list to iout if requested
303  if (this%iprpak /= 0) then
304  call this%write_list()
305  end if
306  !
307  ! -- return
308  return

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