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()
73 subroutine wel_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
76 class(
bndtype),
pointer :: packobj
77 integer(I4B),
intent(in) :: id
78 integer(I4B),
intent(in) :: ibcnum
79 integer(I4B),
intent(in) :: inunit
80 integer(I4B),
intent(in) :: iout
81 character(len=*),
intent(in) :: namemodel
82 character(len=*),
intent(in) :: pakname
83 character(len=*),
intent(in) :: mempath
85 type(
weltype),
pointer :: welobj
92 call packobj%set_names(ibcnum, namemodel, pakname,
ftype, mempath)
96 call welobj%allocate_scalars()
99 call packobj%pack_initialize()
101 packobj%inunit = inunit
104 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'
312 integer(I4B) :: i, node, ict
320 if (this%nbound == 0)
return
323 do i = 1, this%nbound
324 node = this%nodelist(i)
326 if (this%ibound(node) <= 0)
then
331 if (this%iflowred /= 0 .and. q <
dzero)
then
332 ict = this%icelltype(node)
334 bt = this%dis%bot(node)
335 if (this%iflowredlen == 0)
then
336 thick = this%dis%top(node) - bt
340 tp = bt + this%flowred * thick
355 subroutine wel_fc(this, rhs, ia, idxglo, matrix_sln)
358 real(DP),
dimension(:),
intent(inout) :: rhs
359 integer(I4B),
dimension(:),
intent(in) :: ia
360 integer(I4B),
dimension(:),
intent(in) :: idxglo
368 if (this%imover == 1)
then
369 call this%pakmvrobj%fc()
373 do i = 1, this%nbound
375 rhs(n) = rhs(n) + this%rhs(i)
377 call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
381 if (this%imover == 1 .and. this%rhs(i) >
dzero)
then
382 call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
393 subroutine wel_fn(this, rhs, ia, idxglo, matrix_sln)
396 real(DP),
dimension(:),
intent(inout) :: rhs
397 integer(I4B),
dimension(:),
intent(in) :: ia
398 integer(I4B),
dimension(:),
intent(in) :: idxglo
412 do i = 1, this%nbound
413 node = this%nodelist(i)
416 if (this%ibound(node) <= 0)
then
421 ict = this%icelltype(node)
422 if (this%iflowred /= 0 .and. ict /= 0)
then
427 tp = this%dis%top(node)
428 bt = this%dis%bot(node)
429 if (this%iflowredlen == 0)
then
434 tp = bt + this%flowred * thick
436 drterm = drterm * this%q_mult(i)
438 call matrix_sln%add_value_pos(idxglo(ipos), drterm)
439 rhs(node) = rhs(node) + drterm * this%xnew(node)
448 class(
weltype),
intent(inout) :: this
449 character(len=*),
intent(in) :: fname
451 character(len=*),
parameter :: fmtafrcsv = &
452 "(4x, 'AUTO FLOW REDUCE INFORMATION WILL BE SAVED TO FILE: ', a, /4x, &
453 &'OPENED ON UNIT: ', I0)"
456 call openfile(this%ioutafrcsv, this%iout, fname,
'CSV', &
457 filstat_opt=
'REPLACE')
458 write (this%iout, fmtafrcsv) trim(adjustl(fname)), &
460 write (this%ioutafrcsv,
'(a)') &
461 'time,period,step,boundnumber,cellnumber,rate-requested,&
462 &rate-actual,wel-reduction'
470 class(
weltype),
intent(inout) :: this
473 integer(I4B) :: nodereduced
474 integer(I4B) :: nodeuser
477 do i = 1, this%nbound
478 nodereduced = this%nodelist(i)
481 if (this%ibound(nodereduced) <= 0)
then
484 v = this%q_mult(i) + this%rhs(i)
486 nodeuser = this%dis%get_nodeuser(nodereduced)
487 write (this%ioutafrcsv,
'(*(G0,:,","))') &
488 totim,
kper,
kstp, i, nodeuser, this%q_mult(i), this%simvals(i), v
501 class(
weltype),
intent(inout) :: this
504 this%listlabel = trim(this%filtyp)//
' NO.'
505 if (this%dis%ndim == 3)
then
506 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'LAYER'
507 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'ROW'
508 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'COL'
509 elseif (this%dis%ndim == 2)
then
510 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'LAYER'
511 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'CELL2D'
513 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'NODE'
515 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'STRESS RATE'
516 if (this%inamedbound == 1)
then
517 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'BOUNDARY NAME'
551 call this%obs%StoreObsType(
'wel', .true., indx)
556 call this%obs%StoreObsType(
'to-mvr', .true., indx)
561 call this%obs%StoreObsType(
'wel-reduction', .true., indx)
581 call this%obs%obs_bd_clear()
584 do i = 1, this%obs%npakobs
585 obsrv => this%obs%pakobs(i)%obsrv
586 if (obsrv%BndFound)
then
587 do n = 1, obsrv%indxbnds_count
589 jj = obsrv%indxbnds(n)
590 select case (obsrv%ObsTypeId)
592 if (this%imover == 1)
then
593 v = this%pakmvrobj%get_qtomvr(jj)
600 case (
'WEL-REDUCTION')
601 if (this%iflowred > 0)
then
602 v = this%q_mult(jj) + this%rhs(jj)
605 errmsg =
'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
608 call this%obs%SaveOneSimval(obsrv, v)
611 call this%obs%SaveOneSimval(obsrv,
dnodata)
616 if (this%ioutafrcsv > 0)
then
617 call this%wel_afr_csv_write()
625 class(
weltype),
intent(inout) :: this
626 integer(I4B),
intent(in) :: row
630 if (this%iauxmultcol > 0)
then
631 q = this%q(row) * this%auxvar(this%iauxmultcol, row)
647 class(
weltype),
intent(inout) :: this
648 integer(I4B),
intent(in) :: col
649 integer(I4B),
intent(in) :: row
655 bndval = this%q_mult(row)
657 errmsg =
'Programming error. WEL bound value requested column '&
658 &
'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_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.