24 integer(I4B),
public :: iout = 0
26 character(len=LENMODELNAME) :: modelname
28 type(
listtype),
pointer,
private :: boundtaslinks => null()
29 character(len=LINELENGTH),
allocatable,
dimension(:) :: tasfiles
31 character(len=LENTIMESERIESNAME),
allocatable,
dimension(:) :: tasnames
61 character(len=*),
intent(in) :: modelname
62 integer(I4B),
intent(in) :: iout
64 if (
present(dis))
then
68 this%modelname = modelname
70 allocate (this%boundTasLinks)
71 allocate (this%tasfiles(0))
81 integer(I4B) :: nfiles
86 nfiles =
size(this%tasfiles)
87 allocate (this%taslist(nfiles))
88 allocate (this%tasnames(nfiles))
92 tasptr => this%taslist(i)
93 call tasptr%tas_init(this%tasfiles(i), this%modelname, &
94 this%iout, this%tasnames(i))
108 integer(I4B) :: i, j, nlinks, nvals, isize1, isize2, inunit
109 real(DP) :: begintime, endtime
111 character(len=*),
parameter :: fmt5 = &
112 "(/,'Time-array-series controlled arrays in stress period ', &
113 &i0, ', time step ', i0, ':')"
114 10
format(
'"', a,
'" package: ', a,
' array obtained from time-array series "', &
119 endtime = begintime +
delt
124 if (
associated(this%boundTasLinks))
then
125 nlinks = this%boundTasLinks%Count()
128 if (taslink%Iprpak == 1 .and. i == 1)
then
131 if (taslink%UseDefaultProc)
then
132 timearrayseries => taslink%timeArraySeries
133 nvals =
size(taslink%BndArray)
136 call timearrayseries%GetAverageValues(nvals, taslink%BndArray, &
140 if (taslink%ConvertFlux)
then
141 call this%tasmgr_convert_flux(taslink)
146 if (taslink%Iprpak == 1)
then
147 write (this%iout, 10) trim(taslink%PackageName), &
148 trim(taslink%Text), &
149 trim(taslink%timeArraySeries%Name)
152 if (i == nlinks)
then
153 write (this%iout,
'()')
161 if (taslink%UseDefaultProc)
then
162 if (
associated(taslink%RMultArray))
then
163 isize1 =
size(taslink%BndArray)
164 isize2 =
size(taslink%RMultArray)
165 if (isize1 == isize2 .and. isize1 == nvals)
then
167 taslink%BndArray(j) = taslink%BndArray(j) * taslink%RMultArray(j)
170 errmsg =
'Size mismatch between boundary and multiplier arrays'// &
171 ' using time-array series: '// &
172 trim(taslink%TimeArraySeries%Name)
174 inunit = taslink%TimeArraySeries%GetInunit()
194 n = this%boundTasLinks%Count()
201 do i = 1,
size(this%taslist)
202 call this%taslist(i)%da()
206 call this%boundTasLinks%Clear(.true.)
207 deallocate (this%boundTasLinks)
208 deallocate (this%tasfiles)
211 deallocate (this%taslist)
212 deallocate (this%tasnames)
216 this%boundTasLinks => null()
226 character(len=*),
intent(in) :: fname
231 indx =
size(this%tasfiles)
232 this%tasfiles(indx) = fname
243 character(len=*),
intent(in) :: pkgName
245 integer(I4B) :: i, j, nlinks
249 nlinks = this%boundTasLinks%Count()
252 if (
associated(taslink))
then
253 do j = 1,
size(taslink%BndArray)
254 taslink%BndArray(j) =
dzero
260 if (
associated(this%boundTasLinks))
then
262 nlinks = this%boundTasLinks%Count()
265 if (
associated(taslink))
then
267 call this%boundTasLinks%RemoveNode(i, .true.)
276 tasName, text, convertFlux, nodelist, inunit)
279 character(len=*),
intent(in) :: pkgName
280 real(DP),
dimension(:),
pointer :: bndArray
281 integer(I4B),
intent(in) :: iprpak
282 character(len=*),
intent(in) :: tasName
283 character(len=*),
intent(in) :: text
284 logical,
intent(in) :: convertFlux
285 integer(I4B),
dimension(:),
pointer,
contiguous,
intent(in) :: nodelist
286 integer(I4B),
intent(in) :: inunit
288 integer(I4B) :: i, nfiles, iloc
289 character(LINELENGTH) :: ermsg
294 nfiles =
size(this%tasnames)
297 if (this%tasnames(i) == tasname)
then
303 ermsg =
'Error: Time-array series "'//trim(tasname)//
'" not found.'
307 tasptr => this%taslist(iloc)
312 pkgname, bndarray, iprpak, &
314 newtaslink%ConvertFlux = convertflux
315 newtaslink%nodelist => nodelist
318 call this%tasmgr_add_link(newtaslink)
326 integer(I4B),
intent(in) :: indx
332 if (
associated(this%boundTasLinks))
then
345 if (
associated(this%boundtaslinks))
then
362 integer(I4B) :: i, n, noder
365 if (.not. (
associated(this%dis) .and. &
366 associated(taslink%nodelist)))
then
367 errmsg =
'Programming error. Cannot convert flux. Verify that '&
368 &
'a valid DIS instance and nodelist were provided.'
373 n =
size(taslink%BndArray)
375 noder = taslink%nodelist(i)
377 area = this%dis%get_area(noder)
378 taslink%BndArray(i) = taslink%BndArray(i) * area
This module contains simulation constants.
integer(i4b), parameter linelength
maximum length of a standard line
integer(i4b), parameter lenmodelname
maximum length of the model name
integer(i4b), parameter lenhugeline
maximum length of a huge line
integer(i4b), parameter lentimeseriesname
maximum length of a time series name
real(dp), parameter dzero
real constant zero
This module defines variable data types.
This module contains simulation methods.
subroutine, public store_error(msg, terminate)
Store an error message.
subroutine, public store_error_unit(iunit, terminate)
Store the file unit number.
This module contains simulation variables.
character(len=maxcharlen) errmsg
error message string
real(dp), pointer, public totimc
simulation time at start of time step
integer(i4b), pointer, public kstp
current time step number
integer(i4b), pointer, public kper
current stress period number
real(dp), pointer, public delt
length of the current time step
type(timearrayserieslinktype) function, pointer, public gettimearrayserieslinkfromlist(list, idx)
Get time-array series from a list and return.
subroutine, public addtimearrayserieslinktolist(list, tasLink)
Add time-array series to list.
subroutine, public constructtimearrayserieslink(newTasLink, timeArraySeries, pkgName, bndArray, iprpak, text)
Construct a time series of arrays that are linked.
subroutine tasmanager_df(this)
Define the time-array series manager.
subroutine tasmgr_add_link(this, tasLink)
Add a time arrays series link.
subroutine tasmgr_ad(this)
Time step (or subtime step) advance.
subroutine reset(this, pkgName)
Zero out arrays that are represented with time series.
subroutine maketaslink(this, pkgName, bndArray, iprpak, tasName, text, convertFlux, nodelist, inunit)
Make link from time-array series to package array.
subroutine tasmgr_convert_flux(this, tasLink)
Convert the array from a flux to a flow rate by multiplying by the cell area.
subroutine, public tasmanager_cr(this, dis, modelname, iout)
Create the time-array series manager.
integer(i4b) function countlinks(this)
Count number of links in the boundtaslinks list.
subroutine add_tasfile(this, fname)
Add a time-array series file.
type(timearrayserieslinktype) function, pointer getlink(this, indx)
Get link from the boundtaslinks list.
subroutine tasmgr_da(this)
Deallocate.
A generic heterogeneous doubly-linked list.