MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
gwf-wel.f90
Go to the documentation of this file.
1 !> @brief This module contains the WEL package methods
2 !!
3 !! This module contains the overridden methods for the standard WEL package.
4 !! Several methods need to be overridden because of the AUTO_FLOW_REDUCE
5 !! option. Overridden methods include:
6 !! - bnd_cf (AUTO_FLOW_REDUCE)
7 !! - bnd_fc (AUTO_FLOW_REDUCE)
8 !! - bnd_fn (AUTO_FLOW_REDUCE Newton-Raphson terms)
9 !! - bnd_ot_package_flows (write AUTO_FLOW_REDUCE terms to csv file)
10 !! - bnd_da (deallocate AUTO_FLOW_REDUCE variables)
11 !! - bnd_bd_obs (wel-reduction observation added)
12 !!
13 !<
14 
15 module welmodule
16  ! -- modules used by WelModule methods
17  use kindmodule, only: dp, i4b
22  use bndmodule, only: bndtype
23  use bndextmodule, only: bndexttype
26  use observemodule, only: observetype
30  !
31  implicit none
32  !
33  private
34  public :: wel_create
35  !
36  character(len=LENFTYPE) :: ftype = 'WEL' !< package ftype
37  character(len=16) :: text = ' WEL' !< package flow text string
38  !
39  type, extends(bndexttype) :: weltype
40  real(dp), dimension(:), pointer, contiguous :: q => null() !< volumetric well rate
41  integer(I4B), pointer :: iflowred => null() !< flag indicating if the AUTO_FLOW_REDUCE option is active
42  real(dp), pointer :: flowred => null() !< AUTO_FLOW_REDUCE variable
43  integer(I4B), pointer :: ioutafrcsv => null() !< unit number for CSV output file containing wells with reduced puping rates
44  integer(I4B), pointer :: iflowredlen => null() !< flag indicating flowred variable is a length value
45  contains
46  procedure :: allocate_scalars => wel_allocate_scalars
47  procedure :: allocate_arrays => wel_allocate_arrays
48  procedure :: source_options => wel_options
49  procedure :: log_wel_options
50  procedure :: bnd_rp => wel_rp
51  procedure :: bnd_cf => wel_cf
52  procedure :: bnd_fc => wel_fc
53  procedure :: bnd_fn => wel_fn
54  procedure :: bnd_da => wel_da
55  procedure :: define_listlabel
56  procedure :: bound_value => wel_bound_value
57  procedure :: q_mult
58  ! -- methods for observations
59  procedure, public :: bnd_obs_supported => wel_obs_supported
60  procedure, public :: bnd_df_obs => wel_df_obs
61  procedure, public :: bnd_bd_obs => wel_bd_obs
62  ! -- afr
63  procedure, private :: wel_afr_csv_init
64  procedure, private :: wel_afr_csv_write
65  end type weltype
66 
67 contains
68 
69  !> @ brief Create a new package object
70  !!
71  !! Create a new WEL Package object
72  !!
73  !<
74  subroutine wel_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
75  mempath)
76  ! -- dummy variables
77  class(bndtype), pointer :: packobj !< pointer to default package type
78  integer(I4B), intent(in) :: id !< package id
79  integer(I4B), intent(in) :: ibcnum !< boundary condition number
80  integer(I4B), intent(in) :: inunit !< unit number of WEL package input file
81  integer(I4B), intent(in) :: iout !< unit number of model listing file
82  character(len=*), intent(in) :: namemodel !< model name
83  character(len=*), intent(in) :: pakname !< package name
84  character(len=*), intent(in) :: mempath !< input mempath
85  ! -- local variables
86  type(weltype), pointer :: welobj
87  !
88  ! -- allocate the object and assign values to object variables
89  allocate (welobj)
90  packobj => welobj
91  !
92  ! -- create name and memory path
93  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
94  packobj%text = text
95  !
96  ! -- allocate scalars
97  call welobj%allocate_scalars()
98  !
99  ! -- initialize package
100  call packobj%pack_initialize()
101 
102  packobj%inunit = inunit
103  packobj%iout = iout
104  packobj%id = id
105  packobj%ibcnum = ibcnum
106  packobj%ictMemPath = create_mem_path(namemodel, 'NPF')
107  end subroutine wel_create
108 
109  !> @ brief Deallocate package memory
110  !!
111  !! Deallocate WEL package scalars and arrays.
112  !!
113  !<
114  subroutine wel_da(this)
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)
129  end subroutine wel_da
130 
131  !> @ brief Allocate scalars
132  !!
133  !! Allocate and initialize scalars for the WEL package. The base model
134  !! allocate scalars method is also called.
135  !!
136  !<
137  subroutine wel_allocate_scalars(this)
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
157  end subroutine wel_allocate_scalars
158 
159  !> @ brief Allocate arrays
160  !!
161  !! Allocate and initialize arrays for the WEL package
162  !!
163  !<
164  subroutine wel_allocate_arrays(this, nodelist, auxvar)
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)
182  end subroutine wel_allocate_arrays
183 
184  !> @ brief Source additional options for package
185  !!
186  !! Source additional options for WEL package.
187  !!
188  !<
189  subroutine wel_options(this)
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)
252  end subroutine wel_options
253 
254  !> @ brief Log WEL specific package options
255  !<
256  subroutine log_wel_options(this, found)
257  ! -- modules
259  ! -- dummy variables
260  class(weltype), intent(inout) :: this !< BndExtType object
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'
300  end subroutine log_wel_options
301 
302  !> @ brief WEL read and prepare
303  !!
304  !<
305  subroutine wel_rp(this)
306  use tdismodule, only: kper
307  ! -- dummy
308  class(weltype), intent(inout) :: this
309  ! -- local
310  !
311  if (this%iper /= kper) return
312  !
313  ! -- Call the parent class read and prepare
314  call this%BndExtType%bnd_rp()
315  !
316  ! -- Write the list to iout if requested
317  if (this%iprpak /= 0) then
318  call this%write_list()
319  end if
320  end subroutine wel_rp
321 
322  !> @ brief Formulate the package hcof and rhs terms.
323  !!
324  !! Formulate the hcof and rhs terms for the WEL package that will be
325  !! added to the coefficient matrix and right-hand side vector.
326  !!
327  !<
328  subroutine wel_cf(this)
329  ! -- dummy variables
330  class(weltype) :: this !< WelType object
331  ! -- local variables
332  integer(I4B) :: i, node, ict
333  real(DP) :: qmult
334  real(DP) :: q
335  real(DP) :: tp
336  real(DP) :: bt
337  real(DP) :: thick
338  !
339  ! -- Return if no wells
340  if (this%nbound == 0) return
341  !
342  ! -- Calculate hcof and rhs for each well entry
343  do i = 1, this%nbound
344  node = this%nodelist(i)
345  this%hcof(i) = dzero
346  if (this%ibound(node) <= 0) then
347  this%rhs(i) = dzero
348  cycle
349  end if
350  q = this%q_mult(i)
351  if (this%iflowred /= 0 .and. q < dzero) then
352  ict = this%icelltype(node)
353  if (ict /= 0) then
354  bt = this%dis%bot(node)
355  if (this%iflowredlen == 0) then
356  thick = this%dis%top(node) - bt
357  else
358  thick = done
359  end if
360  tp = bt + this%flowred * thick
361  qmult = sqsaturation(tp, bt, this%xnew(node))
362  q = q * qmult
363  end if
364  end if
365  this%rhs(i) = -q
366  end do
367  end subroutine wel_cf
368 
369  !> @ brief Copy hcof and rhs terms into solution.
370  !!
371  !! Add the hcof and rhs terms for the WEL package to the
372  !! coefficient matrix and right-hand side vector.
373  !!
374  !<
375  subroutine wel_fc(this, rhs, ia, idxglo, matrix_sln)
376  ! -- dummy variables
377  class(weltype) :: this !< WelType object
378  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
379  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
380  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
381  class(matrixbasetype), pointer :: matrix_sln !< solution coefficient matrix
382  ! -- local variables
383  integer(I4B) :: i
384  integer(I4B) :: n
385  integer(I4B) :: ipos
386  !
387  ! -- pakmvrobj fc
388  if (this%imover == 1) then
389  call this%pakmvrobj%fc()
390  end if
391  !
392  ! -- Copy package rhs and hcof into solution rhs and amat
393  do i = 1, this%nbound
394  n = this%nodelist(i)
395  rhs(n) = rhs(n) + this%rhs(i)
396  ipos = ia(n)
397  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
398  !
399  ! -- If mover is active and this well is discharging,
400  ! store available water (as positive value).
401  if (this%imover == 1 .and. this%rhs(i) > dzero) then
402  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
403  end if
404  end do
405  end subroutine wel_fc
406 
407  !> @ brief Add Newton-Raphson terms for package into solution.
408  !!
409  !! Calculate and add the Newton-Raphson terms for the WEL package to the
410  !! coefficient matrix and right-hand side vector.
411  !!
412  !<
413  subroutine wel_fn(this, rhs, ia, idxglo, matrix_sln)
414  ! -- dummy variables
415  class(weltype) :: this !< WelType object
416  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
417  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
418  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
419  class(matrixbasetype), pointer :: matrix_sln !< solution coefficient matrix
420  ! -- local variables
421  integer(I4B) :: i
422  integer(I4B) :: node
423  integer(I4B) :: ipos
424  integer(I4B) :: ict
425  real(DP) :: drterm
426  real(DP) :: q
427  real(DP) :: tp
428  real(DP) :: bt
429  real(DP) :: thick
430  !
431  ! -- Copy package rhs and hcof into solution rhs and amat
432  do i = 1, this%nbound
433  node = this%nodelist(i)
434  !
435  ! -- test if node is constant or inactive
436  if (this%ibound(node) <= 0) then
437  cycle
438  end if
439  !
440  ! -- well rate is possibly head dependent
441  ict = this%icelltype(node)
442  if (this%iflowred /= 0 .and. ict /= 0) then
443  ipos = ia(node)
444  q = -this%rhs(i)
445  if (q < dzero) then
446  ! -- calculate derivative for well
447  tp = this%dis%top(node)
448  bt = this%dis%bot(node)
449  thick = tp - bt
450  tp = bt + this%flowred * thick
451  drterm = sqsaturationderivative(tp, bt, this%xnew(node))
452  drterm = drterm * this%q_mult(i)
453  !--fill amat and rhs with newton-raphson terms
454  call matrix_sln%add_value_pos(idxglo(ipos), drterm)
455  rhs(node) = rhs(node) + drterm * this%xnew(node)
456  end if
457  end if
458  end do
459  end subroutine wel_fn
460 
461  !> @brief Initialize the auto flow reduce csv output file
462  subroutine wel_afr_csv_init(this, fname)
463  ! -- dummy variables
464  class(weltype), intent(inout) :: this !< WelType object
465  character(len=*), intent(in) :: fname
466  ! -- format
467  character(len=*), parameter :: fmtafrcsv = &
468  "(4x, 'AUTO FLOW REDUCE INFORMATION WILL BE SAVED TO FILE: ', a, /4x, &
469  &'OPENED ON UNIT: ', I0)"
470 
471  this%ioutafrcsv = getunit()
472  call openfile(this%ioutafrcsv, this%iout, fname, 'CSV', &
473  filstat_opt='REPLACE')
474  write (this%iout, fmtafrcsv) trim(adjustl(fname)), &
475  this%ioutafrcsv
476  write (this%ioutafrcsv, '(a)') &
477  'time,period,step,boundnumber,cellnumber,rate-requested,&
478  &rate-actual,wel-reduction'
479  end subroutine wel_afr_csv_init
480 
481  !> @brief Write out auto flow reductions only when & where they occur
482  subroutine wel_afr_csv_write(this)
483  ! -- modules
484  use tdismodule, only: totim, kstp, kper
485  ! -- dummy variables
486  class(weltype), intent(inout) :: this !< WelType object
487  ! -- local
488  integer(I4B) :: i
489  integer(I4B) :: nodereduced
490  integer(I4B) :: nodeuser
491  real(DP) :: v
492  ! -- format
493  do i = 1, this%nbound
494  nodereduced = this%nodelist(i)
495  !
496  ! -- test if node is constant or inactive
497  if (this%ibound(nodereduced) <= 0) then
498  cycle
499  end if
500  v = this%q_mult(i) + this%rhs(i)
501  if (v < dzero) then
502  nodeuser = this%dis%get_nodeuser(nodereduced)
503  write (this%ioutafrcsv, '(*(G0,:,","))') &
504  totim, kper, kstp, i, nodeuser, this%q_mult(i), this%simvals(i), v
505  end if
506  end do
507  end subroutine wel_afr_csv_write
508 
509  !> @ brief Define the list label for the package
510  !!
511  !! Method defined the list label for the WEL package. The list label is
512  !! the heading that is written to iout when PRINT_INPUT option is used.
513  !!
514  !<
515  subroutine define_listlabel(this)
516  ! -- dummy variables
517  class(weltype), intent(inout) :: this !< WelType object
518  !
519  ! -- create the header list label
520  this%listlabel = trim(this%filtyp)//' NO.'
521  if (this%dis%ndim == 3) then
522  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
523  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
524  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
525  elseif (this%dis%ndim == 2) then
526  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
527  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
528  else
529  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
530  end if
531  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'STRESS RATE'
532  if (this%inamedbound == 1) then
533  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
534  end if
535  end subroutine define_listlabel
536 
537  ! -- Procedures related to observations
538 
539  !> @brief Determine if observations are supported.
540  !!
541  !! Function to determine if observations are supported by the WEL package.
542  !! Observations are supported by the WEL package.
543  !!
544  !! @return wel_obs_supported boolean indicating if observations are supported
545  !!
546  !<
547  logical function wel_obs_supported(this)
548  ! -- dummy variables
549  class(weltype) :: this !< WelType object
550  !
551  ! -- set boolean
552  wel_obs_supported = .true.
553  end function wel_obs_supported
554 
555  !> @brief Define the observation types available in the package
556  !!
557  !! Method to define the observation types available in the WEL package.
558  !!
559  !<
560  subroutine wel_df_obs(this)
561  ! -- dummy variables
562  class(weltype) :: this !< WelType object
563  ! -- local variables
564  integer(I4B) :: indx
565  !
566  ! -- initialize observations
567  call this%obs%StoreObsType('wel', .true., indx)
568  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
569  !
570  ! -- Store obs type and assign procedure pointer
571  ! for to-mvr observation type.
572  call this%obs%StoreObsType('to-mvr', .true., indx)
573  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
574  !
575  ! -- Store obs type and assign procedure pointer
576  ! for wel-reduction observation type.
577  call this%obs%StoreObsType('wel-reduction', .true., indx)
578  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
579  end subroutine wel_df_obs
580 
581  !> @brief Save observations for the package
582  !!
583  !! Method to save simulated values for the WEL package.
584  !!
585  !<
586  subroutine wel_bd_obs(this)
587  ! -- dummy variables
588  class(weltype) :: this !< WelType object
589  ! -- local variables
590  integer(I4B) :: i
591  integer(I4B) :: n
592  integer(I4B) :: jj
593  real(DP) :: v
594  type(observetype), pointer :: obsrv => null()
595  !
596  ! -- clear the observations
597  call this%obs%obs_bd_clear()
598  !
599  ! -- Save simulated values for all of package's observations.
600  do i = 1, this%obs%npakobs
601  obsrv => this%obs%pakobs(i)%obsrv
602  if (obsrv%BndFound) then
603  do n = 1, obsrv%indxbnds_count
604  v = dnodata
605  jj = obsrv%indxbnds(n)
606  select case (obsrv%ObsTypeId)
607  case ('TO-MVR')
608  if (this%imover == 1) then
609  v = this%pakmvrobj%get_qtomvr(jj)
610  if (v > dzero) then
611  v = -v
612  end if
613  end if
614  case ('WEL')
615  v = this%simvals(jj)
616  case ('WEL-REDUCTION')
617  if (this%iflowred > 0) then
618  v = this%q_mult(jj) + this%rhs(jj)
619  end if
620  case default
621  errmsg = 'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
622  call store_error(errmsg)
623  end select
624  call this%obs%SaveOneSimval(obsrv, v)
625  end do
626  else
627  call this%obs%SaveOneSimval(obsrv, dnodata)
628  end if
629  end do
630  !
631  ! -- Write the auto flow reduce csv file entries for this step
632  if (this%ioutafrcsv > 0) then
633  call this%wel_afr_csv_write()
634  end if
635  end subroutine wel_bd_obs
636 
637  function q_mult(this, row) result(q)
638  ! -- modules
639  use constantsmodule, only: dzero
640  ! -- dummy variables
641  class(weltype), intent(inout) :: this !< BndExtType object
642  integer(I4B), intent(in) :: row
643  ! -- result
644  real(dp) :: q
645  !
646  if (this%iauxmultcol > 0) then
647  q = this%q(row) * this%auxvar(this%iauxmultcol, row)
648  else
649  q = this%q(row)
650  end if
651  end function q_mult
652 
653  !> @ brief Return a bound value
654  !!
655  !! Return a bound value associated with an ncolbnd index
656  !! and row.
657  !!
658  !<
659  function wel_bound_value(this, col, row) result(bndval)
660  ! -- modules
661  use constantsmodule, only: dzero
662  ! -- dummy variables
663  class(weltype), intent(inout) :: this !< BndExtType object
664  integer(I4B), intent(in) :: col
665  integer(I4B), intent(in) :: row
666  ! -- result
667  real(dp) :: bndval
668  !
669  select case (col)
670  case (1)
671  bndval = this%q_mult(row)
672  case default
673  errmsg = 'Programming error. WEL bound value requested column '&
674  &'outside range of ncolbnd (1).'
675  call store_error(errmsg)
676  call store_error_filename(this%input_fname)
677  end select
678  end function wel_bound_value
679 
680 end module welmodule
This module contains block parser methods.
Definition: BlockParser.f90:7
This module contains the extended boundary package.
This module contains the base boundary package.
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:95
real(dp), parameter dem1
real constant 1e-1
Definition: Constants.f90:103
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
Definition: Constants.f90:39
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
integer(i4b) function, public getunit()
Get a free unit number.
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
subroutine, public urword(line, icol, istart, istop, ncode, n, r, iout, in)
Extract a word from a string.
This module defines variable data types.
Definition: kind.f90:8
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains the derived types ObserveType and ObsDataType.
Definition: Observe.f90:15
This module contains the derived type ObsType.
Definition: Obs.f90:127
subroutine, public defaultobsidprocessor(obsrv, dis, inunitobs, iout)
@ brief Process IDstring provided for each observation
Definition: Obs.f90:246
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_warning(msg, substring)
Store warning message.
Definition: Sim.f90:236
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
character(len=maxcharlen) warnmsg
warning message string
real(dp) function sqsaturationderivative(top, bot, x, c1, c2)
@ brief sQSaturationDerivative
real(dp) function sqsaturation(top, bot, x, c1, c2)
@ brief sQSaturation
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
This module contains the WEL package methods.
Definition: gwf-wel.f90:15
subroutine, public wel_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
@ brief Create a new package object
Definition: gwf-wel.f90:76
subroutine wel_allocate_scalars(this)
@ brief Allocate scalars
Definition: gwf-wel.f90:138
subroutine wel_fc(this, rhs, ia, idxglo, matrix_sln)
@ brief Copy hcof and rhs terms into solution.
Definition: gwf-wel.f90:376
subroutine define_listlabel(this)
@ brief Define the list label for the package
Definition: gwf-wel.f90:516
subroutine wel_rp(this)
@ brief WEL read and prepare
Definition: gwf-wel.f90:306
subroutine wel_allocate_arrays(this, nodelist, auxvar)
@ brief Allocate arrays
Definition: gwf-wel.f90:165
subroutine wel_options(this)
@ brief Source additional options for package
Definition: gwf-wel.f90:190
subroutine wel_afr_csv_write(this)
Write out auto flow reductions only when & where they occur.
Definition: gwf-wel.f90:483
real(dp) function q_mult(this, row)
Definition: gwf-wel.f90:638
real(dp) function wel_bound_value(this, col, row)
@ brief Return a bound value
Definition: gwf-wel.f90:660
subroutine wel_da(this)
@ brief Deallocate package memory
Definition: gwf-wel.f90:115
subroutine wel_fn(this, rhs, ia, idxglo, matrix_sln)
@ brief Add Newton-Raphson terms for package into solution.
Definition: gwf-wel.f90:414
subroutine wel_afr_csv_init(this, fname)
Initialize the auto flow reduce csv output file.
Definition: gwf-wel.f90:463
character(len=lenftype) ftype
package ftype
Definition: gwf-wel.f90:36
subroutine wel_cf(this)
@ brief Formulate the package hcof and rhs terms.
Definition: gwf-wel.f90:329
subroutine wel_bd_obs(this)
Save observations for the package.
Definition: gwf-wel.f90:587
logical function wel_obs_supported(this)
Determine if observations are supported.
Definition: gwf-wel.f90:548
subroutine log_wel_options(this, found)
@ brief Log WEL specific package options
Definition: gwf-wel.f90:257
character(len=16) text
package flow text string
Definition: gwf-wel.f90:37
subroutine wel_df_obs(this)
Define the observation types available in the package.
Definition: gwf-wel.f90:561
@ brief BndType