36 character(len=LENFTYPE) ::
ftype =
'WEL'
37 character(len=16) ::
text =
' WEL'
40 real(dp),
dimension(:),
pointer,
contiguous :: q => null()
41 integer(I4B),
pointer :: iflowred => null()
42 real(dp),
pointer :: flowred => null()
43 integer(I4B),
pointer :: ioutafrcsv => null()
44 integer(I4B),
pointer :: iflowredlen => null()
74 subroutine wel_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
77 class(
bndtype),
pointer :: packobj
78 integer(I4B),
intent(in) :: id
79 integer(I4B),
intent(in) :: ibcnum
80 integer(I4B),
intent(in) :: inunit
81 integer(I4B),
intent(in) :: iout
82 character(len=*),
intent(in) :: namemodel
83 character(len=*),
intent(in) :: pakname
84 character(len=*),
intent(in) :: mempath
86 type(
weltype),
pointer :: welobj
93 call packobj%set_names(ibcnum, namemodel, pakname,
ftype, mempath)
97 call welobj%allocate_scalars()
100 call packobj%pack_initialize()
102 packobj%inunit = inunit
105 packobj%ibcnum = ibcnum
121 call this%BndExtType%bnd_da()
144 call this%BndExtType%allocate_scalars()
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)
169 integer(I4B),
dimension(:),
pointer,
contiguous,
optional :: nodelist
170 real(DP),
dimension(:, :),
pointer,
contiguous,
optional :: auxvar
174 call this%BndExtType%allocate_arrays(nodelist, auxvar)
177 call mem_setptr(this%q,
'Q', this%input_mempath)
181 'Q', this%input_mempath)
195 class(
weltype),
intent(inout) :: this
197 character(len=LINELENGTH) :: fname
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,').')"
206 call this%BndExtType%source_options()
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, &
215 if (found%iflowredlen)
then
216 if (found%flowred .eqv. .false.)
then
218 'FLOW_REDUCTION_LENGTH option specified but a AUTO_FLOW_REDUCTION value &
219 &is not specified. The FLOW_REDUCTION_LENGTH option will be ignored.'
226 if (found%flowred)
then
228 if (this%flowred <=
dzero)
then
229 if (found%iflowredlen)
then
231 'An AUTO_FLOW_REDUCTION value less than or equal to zero cannot be &
232 &specified if the FLOW_REDUCTION_LENGTH option is specified.'
237 else if (this%flowred >
done .and. this%iflowredlen == 0)
then
242 if (found%afrcsvfile)
then
243 call this%wel_afr_csv_init(fname)
246 if (found%mover)
then
251 call this%log_wel_options(found)
260 class(
weltype),
intent(inout) :: this
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,').')"
272 write (this%iout,
'(/1x,a)')
'PROCESSING '//trim(adjustl(this%text)) &
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'
281 if (found%flowred)
then
282 if (this%iflowredlen == 0)
then
283 write (this%iout, fmtflowredv) this%flowred
285 write (this%iout, fmtflowredl) this%flowred
289 if (found%afrcsvfile)
then
293 if (found%mover)
then
294 write (this%iout,
'(4x,A)')
'MOVER OPTION ENABLED'
298 write (this%iout,
'(1x,a)') &
299 'END OF '//trim(adjustl(this%text))//
' OPTIONS'
308 class(
weltype),
intent(inout) :: this
311 if (this%iper /=
kper)
return
314 call this%BndExtType%bnd_rp()
317 if (this%iprpak /= 0)
then
318 call this%write_list()
332 integer(I4B) :: i, node, ict
340 if (this%nbound == 0)
return
343 do i = 1, this%nbound
344 node = this%nodelist(i)
346 if (this%ibound(node) <= 0)
then
351 if (this%iflowred /= 0 .and. q <
dzero)
then
352 ict = this%icelltype(node)
354 bt = this%dis%bot(node)
355 if (this%iflowredlen == 0)
then
356 thick = this%dis%top(node) - bt
360 tp = bt + this%flowred * thick
375 subroutine wel_fc(this, rhs, ia, idxglo, matrix_sln)
378 real(DP),
dimension(:),
intent(inout) :: rhs
379 integer(I4B),
dimension(:),
intent(in) :: ia
380 integer(I4B),
dimension(:),
intent(in) :: idxglo
388 if (this%imover == 1)
then
389 call this%pakmvrobj%fc()
393 do i = 1, this%nbound
395 rhs(n) = rhs(n) + this%rhs(i)
397 call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
401 if (this%imover == 1 .and. this%rhs(i) >
dzero)
then
402 call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
413 subroutine wel_fn(this, rhs, ia, idxglo, matrix_sln)
416 real(DP),
dimension(:),
intent(inout) :: rhs
417 integer(I4B),
dimension(:),
intent(in) :: ia
418 integer(I4B),
dimension(:),
intent(in) :: idxglo
432 do i = 1, this%nbound
433 node = this%nodelist(i)
436 if (this%ibound(node) <= 0)
then
441 ict = this%icelltype(node)
442 if (this%iflowred /= 0 .and. ict /= 0)
then
447 tp = this%dis%top(node)
448 bt = this%dis%bot(node)
450 tp = bt + this%flowred * thick
452 drterm = drterm * this%q_mult(i)
454 call matrix_sln%add_value_pos(idxglo(ipos), drterm)
455 rhs(node) = rhs(node) + drterm * this%xnew(node)
464 class(
weltype),
intent(inout) :: this
465 character(len=*),
intent(in) :: fname
467 character(len=*),
parameter :: fmtafrcsv = &
468 "(4x, 'AUTO FLOW REDUCE INFORMATION WILL BE SAVED TO FILE: ', a, /4x, &
469 &'OPENED ON UNIT: ', I0)"
472 call openfile(this%ioutafrcsv, this%iout, fname,
'CSV', &
473 filstat_opt=
'REPLACE')
474 write (this%iout, fmtafrcsv) trim(adjustl(fname)), &
476 write (this%ioutafrcsv,
'(a)') &
477 'time,period,step,boundnumber,cellnumber,rate-requested,&
478 &rate-actual,wel-reduction'
486 class(
weltype),
intent(inout) :: this
489 integer(I4B) :: nodereduced
490 integer(I4B) :: nodeuser
493 do i = 1, this%nbound
494 nodereduced = this%nodelist(i)
497 if (this%ibound(nodereduced) <= 0)
then
500 v = this%q_mult(i) + this%rhs(i)
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
517 class(
weltype),
intent(inout) :: this
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'
529 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'NODE'
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'
567 call this%obs%StoreObsType(
'wel', .true., indx)
572 call this%obs%StoreObsType(
'to-mvr', .true., indx)
577 call this%obs%StoreObsType(
'wel-reduction', .true., indx)
597 call this%obs%obs_bd_clear()
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
605 jj = obsrv%indxbnds(n)
606 select case (obsrv%ObsTypeId)
608 if (this%imover == 1)
then
609 v = this%pakmvrobj%get_qtomvr(jj)
616 case (
'WEL-REDUCTION')
617 if (this%iflowred > 0)
then
618 v = this%q_mult(jj) + this%rhs(jj)
621 errmsg =
'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
624 call this%obs%SaveOneSimval(obsrv, v)
627 call this%obs%SaveOneSimval(obsrv,
dnodata)
632 if (this%ioutafrcsv > 0)
then
633 call this%wel_afr_csv_write()
641 class(
weltype),
intent(inout) :: this
642 integer(I4B),
intent(in) :: row
646 if (this%iauxmultcol > 0)
then
647 q = this%q(row) * this%auxvar(this%iauxmultcol, row)
663 class(
weltype),
intent(inout) :: this
664 integer(I4B),
intent(in) :: col
665 integer(I4B),
intent(in) :: row
671 bndval = this%q_mult(row)
673 errmsg =
'Programming error. WEL bound value requested column '&
674 &
'outside range of ncolbnd (1).'
This module contains block parser methods.
This module contains the extended boundary package.
This module contains the base boundary package.
This module contains simulation constants.
integer(i4b), parameter linelength
maximum length of a standard line
real(dp), parameter dnodata
real no data constant
real(dp), parameter dem1
real constant 1e-1
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
real(dp), parameter dzero
real constant zero
real(dp), parameter done
real constant 1
This module defines variable data types.
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.
This module contains the derived type ObsType.
subroutine, public defaultobsidprocessor(obsrv, dis, inunitobs, iout)
@ brief Process IDstring provided for each observation
This module contains simulation methods.
subroutine, public store_warning(msg, substring)
Store warning message.
subroutine, public store_error(msg, terminate)
Store an error message.
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
This module contains simulation variables.
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
integer(i4b), pointer, public kstp
current time step number
integer(i4b), pointer, public kper
current stress period number
This module contains the WEL package methods.
subroutine, public wel_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
@ brief Create a new package object
subroutine wel_allocate_scalars(this)
@ brief Allocate scalars
subroutine wel_fc(this, rhs, ia, idxglo, matrix_sln)
@ brief Copy hcof and rhs terms into solution.
subroutine define_listlabel(this)
@ brief Define the list label for the package
subroutine wel_rp(this)
@ brief WEL read and prepare
subroutine wel_allocate_arrays(this, nodelist, auxvar)
@ brief Allocate arrays
subroutine wel_options(this)
@ brief Source additional options for package
subroutine wel_afr_csv_write(this)
Write out auto flow reductions only when & where they occur.
real(dp) function q_mult(this, row)
real(dp) function wel_bound_value(this, col, row)
@ brief Return a bound value
subroutine wel_da(this)
@ brief Deallocate package memory
subroutine wel_fn(this, rhs, ia, idxglo, matrix_sln)
@ brief Add Newton-Raphson terms for package into solution.
subroutine wel_afr_csv_init(this, fname)
Initialize the auto flow reduce csv output file.
character(len=lenftype) ftype
package ftype
subroutine wel_cf(this)
@ brief Formulate the package hcof and rhs terms.
subroutine wel_bd_obs(this)
Save observations for the package.
logical function wel_obs_supported(this)
Determine if observations are supported.
subroutine log_wel_options(this, found)
@ brief Log WEL specific package options
character(len=16) text
package flow text string
subroutine wel_df_obs(this)
Define the observation types available in the package.