MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
TimeArraySeriesManager.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b
4  use simvariablesmodule, only: errmsg
7  use listmodule, only: listtype
9  use tdismodule, only: delt, totimc, kper, kstp
15  use basedismodule, only: disbasetype
16 
17  implicit none
18 
19  private
21 
23  ! -- Public members
24  integer(I4B), public :: iout = 0 ! output unit num
25  class(disbasetype), pointer :: dis => null() ! pointer to dis
26  character(len=LENMODELNAME) :: modelname
27  ! -- Private members
28  type(listtype), pointer, private :: boundtaslinks => null() ! list of TAS links
29  character(len=LINELENGTH), allocatable, dimension(:) :: tasfiles ! list of TA file names
30  type(timearrayseriestype), dimension(:), pointer, contiguous :: taslist ! array of TA pointers
31  character(len=LENTIMESERIESNAME), allocatable, dimension(:) :: tasnames ! array of TA names
32 
33  contains
34 
35  ! -- Public procedures
36  procedure, public :: tasmanager_df
37  procedure, public :: ad => tasmgr_ad
38  procedure, public :: da => tasmgr_da
39  procedure, public :: add_tasfile
40  procedure, public :: countlinks
41  procedure, public :: getlink
42  procedure, public :: maketaslink
43  procedure, public :: reset
44  ! -- Private procedures
45  procedure, private :: tasmgr_add_link
46  procedure, private :: tasmgr_convert_flux
48 
49 contains
50 
51 ! -- Type-bound procedures of TimeArraySeriesManagerType
52 
53  ! -- Public procedures
54 
55  !> @brief Create the time-array series manager
56  !<
57  subroutine tasmanager_cr(this, dis, modelname, iout)
58  ! -- dummy
59  type(timearrayseriesmanagertype) :: this
60  class(disbasetype), pointer, optional :: dis
61  character(len=*), intent(in) :: modelname
62  integer(I4B), intent(in) :: iout
63  !
64  if (present(dis)) then
65  this%dis => dis
66  end if
67  !
68  this%modelname = modelname
69  this%iout = iout
70  allocate (this%boundTasLinks)
71  allocate (this%tasfiles(0))
72  end subroutine tasmanager_cr
73 
74  !> @brief Define the time-array series manager
75  !<
76  subroutine tasmanager_df(this)
77  ! -- dummy
78  class(timearrayseriesmanagertype) :: this
79  ! -- local
80  type(timearrayseriestype), pointer :: tasptr => null()
81  integer(I4B) :: nfiles
82  integer(I4B) :: i
83  !
84  ! -- determine how many tasfiles. This is the number of time array series
85  ! so allocate arrays to store them
86  nfiles = size(this%tasfiles)
87  allocate (this%taslist(nfiles))
88  allocate (this%tasnames(nfiles))
89  !
90  ! -- Setup a time array series for each file specified
91  do i = 1, nfiles
92  tasptr => this%taslist(i)
93  call tasptr%tas_init(this%tasfiles(i), this%modelname, &
94  this%iout, this%tasnames(i))
95  end do
96  end subroutine tasmanager_df
97 
98  !> @brief Time step (or subtime step) advance.
99  !!
100  !! Call this each time step or subtime step.
101  !<
102  subroutine tasmgr_ad(this)
103  ! -- dummy
104  class(timearrayseriesmanagertype) :: this
105  ! -- local
106  type(timearrayserieslinktype), pointer :: tasLink => null()
107  type(timearrayseriestype), pointer :: timearrayseries => null()
108  integer(I4B) :: i, j, nlinks, nvals, isize1, isize2, inunit
109  real(DP) :: begintime, endtime
110  ! -- formats
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 "', &
115  a, '"')
116  !
117  ! -- Initialize time variables
118  begintime = totimc
119  endtime = begintime + delt
120  !
121  ! -- Iterate through boundtaslinks and update specified
122  ! array with array of average values obtained from
123  ! appropriate time series.
124  if (associated(this%boundTasLinks)) then
125  nlinks = this%boundTasLinks%Count()
126  do i = 1, nlinks
127  taslink => gettimearrayserieslinkfromlist(this%boundTasLinks, i)
128  if (taslink%Iprpak == 1 .and. i == 1) then
129  write (this%iout, fmt5) kper, kstp
130  end if
131  if (taslink%UseDefaultProc) then
132  timearrayseries => taslink%timeArraySeries
133  nvals = size(taslink%BndArray)
134  !
135  ! -- Fill the package array with integrated values
136  call timearrayseries%GetAverageValues(nvals, taslink%BndArray, &
137  begintime, endtime)
138  !
139  ! -- If conversion from flux to flow is required, multiply by cell area
140  if (taslink%ConvertFlux) then
141  call this%tasmgr_convert_flux(taslink)
142  end if
143  !
144  ! -- If PRINT_INPUT is specified, write information
145  ! regarding source of time-array series data
146  if (taslink%Iprpak == 1) then
147  write (this%iout, 10) trim(taslink%PackageName), &
148  trim(taslink%Text), &
149  trim(taslink%timeArraySeries%Name)
150  end if
151  end if
152  if (i == nlinks) then
153  write (this%iout, '()')
154  end if
155  end do
156  !
157  ! -- Now that all array values have been substituted, can now multiply
158  ! an array by a multiplier array
159  do i = 1, nlinks
160  taslink => gettimearrayserieslinkfromlist(this%boundTasLinks, i)
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
166  do j = 1, nvals
167  taslink%BndArray(j) = taslink%BndArray(j) * taslink%RMultArray(j)
168  end do
169  else
170  errmsg = 'Size mismatch between boundary and multiplier arrays'// &
171  ' using time-array series: '// &
172  trim(taslink%TimeArraySeries%Name)
173  call store_error(errmsg)
174  inunit = taslink%TimeArraySeries%GetInunit()
175  call store_error_unit(inunit)
176  end if
177  end if
178  end if
179  end do
180  end if
181  end subroutine tasmgr_ad
182 
183  !> @brief Deallocate
184  !<
185  subroutine tasmgr_da(this)
186  ! -- dummy
187  class(timearrayseriesmanagertype) :: this
188  ! -- local
189  integer :: i, n
190  type(timearrayserieslinktype), pointer :: tasLink => null()
191  !
192  ! -- Deallocate contents of each TimeArraySeriesType object in list
193  ! of time-array series links.
194  n = this%boundTasLinks%Count()
195  do i = 1, n
196  taslink => gettimearrayserieslinkfromlist(this%boundTasLinks, i)
197  call taslink%da()
198  end do
199  !
200  ! -- Go through and deallocate individual time array series
201  do i = 1, size(this%taslist)
202  call this%taslist(i)%da()
203  end do
204  !
205  ! -- Deallocate the list of time-array series links.
206  call this%boundTasLinks%Clear(.true.)
207  deallocate (this%boundTasLinks)
208  deallocate (this%tasfiles)
209  !
210  ! -- Deallocate the time array series
211  deallocate (this%taslist)
212  deallocate (this%tasnames)
213  !
214  ! -- nullify pointers
215  this%dis => null()
216  this%boundTasLinks => null()
217  end subroutine tasmgr_da
218 
219  !> @brief Add a time-array series file
220  !<
221  subroutine add_tasfile(this, fname)
222  ! -- modules
224  ! -- dummy
225  class(timearrayseriesmanagertype) :: this
226  character(len=*), intent(in) :: fname
227  ! -- local
228  integer(I4B) :: indx
229  !
230  call expandarray(this%tasfiles, 1)
231  indx = size(this%tasfiles)
232  this%tasfiles(indx) = fname
233  end subroutine add_tasfile
234 
235  !> @brief Zero out arrays that are represented with time series
236  !!
237  !! Delete all existing links from time array series to package arrays as they
238  !! will need to be created with a new BEGIN PERIOD block
239  !<
240  subroutine reset(this, pkgName)
241  ! -- dummy
242  class(timearrayseriesmanagertype) :: this
243  character(len=*), intent(in) :: pkgName
244  ! -- local
245  integer(I4B) :: i, j, nlinks
246  type(timearrayserieslinktype), pointer :: taslink
247  !
248  ! -- Reassign all linked elements to zero
249  nlinks = this%boundTasLinks%Count()
250  do i = 1, nlinks
251  taslink => gettimearrayserieslinkfromlist(this%boundTasLinks, i)
252  if (associated(taslink)) then
253  do j = 1, size(taslink%BndArray)
254  taslink%BndArray(j) = dzero
255  end do
256  end if
257  end do
258  !
259  ! -- Delete all existing time array links
260  if (associated(this%boundTasLinks)) then
261  ! Deallocate and remove all links belonging to package
262  nlinks = this%boundTasLinks%Count()
263  do i = nlinks, 1, -1
264  taslink => gettimearrayserieslinkfromlist(this%boundTasLinks, i)
265  if (associated(taslink)) then
266  call taslink%da()
267  call this%boundTasLinks%RemoveNode(i, .true.)
268  end if
269  end do
270  end if
271  end subroutine reset
272 
273  !> @brief Make link from time-array series to package array
274  !<
275  subroutine maketaslink(this, pkgName, bndArray, iprpak, &
276  tasName, text, convertFlux, nodelist, inunit)
277  ! -- dummy
278  class(timearrayseriesmanagertype) :: this
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
287  ! -- local
288  integer(I4B) :: i, nfiles, iloc
289  character(LINELENGTH) :: ermsg
290  type(timearrayserieslinktype), pointer :: newTasLink
291  type(timearrayseriestype), pointer :: tasptr => null()
292  !
293  ! -- Find the time array series
294  nfiles = size(this%tasnames)
295  iloc = 0
296  do i = 1, nfiles
297  if (this%tasnames(i) == tasname) then
298  iloc = i
299  exit
300  end if
301  end do
302  if (iloc == 0) then
303  ermsg = 'Error: Time-array series "'//trim(tasname)//'" not found.'
304  call store_error(ermsg)
305  call store_error_unit(inunit)
306  end if
307  tasptr => this%taslist(iloc)
308  !
309  ! -- Construct a time-array series link
310  newtaslink => null()
311  call constructtimearrayserieslink(newtaslink, tasptr, &
312  pkgname, bndarray, iprpak, &
313  text)
314  newtaslink%ConvertFlux = convertflux
315  newtaslink%nodelist => nodelist
316  !
317  ! -- Add link to list of links
318  call this%tasmgr_add_link(newtaslink)
319  end subroutine maketaslink
320 
321  !> @brief Get link from the boundtaslinks list
322  !<
323  function getlink(this, indx) result(tasLink)
324  ! -- dummy
325  class(timearrayseriesmanagertype) :: this
326  integer(I4B), intent(in) :: indx
327  ! -- return
328  type(timearrayserieslinktype), pointer :: taslink
329  !
330  taslink => null()
331  !
332  if (associated(this%boundTasLinks)) then
333  taslink => gettimearrayserieslinkfromlist(this%boundTasLinks, indx)
334  end if
335  end function getlink
336 
337  !> @brief Count number of links in the boundtaslinks list
338  !<
339  function countlinks(this)
340  ! -- return
341  integer(I4B) :: countlinks
342  ! -- dummy
343  class(timearrayseriesmanagertype) :: this
344  !
345  if (associated(this%boundtaslinks)) then
346  countlinks = this%boundTasLinks%Count()
347  else
348  countlinks = 0
349  end if
350  end function countlinks
351 
352  ! -- Private procedures
353 
354  !> @brief Convert the array from a flux to a flow rate by multiplying by the
355  !! cell area
356  !<
357  subroutine tasmgr_convert_flux(this, tasLink)
358  ! -- dummy
359  class(timearrayseriesmanagertype) :: this
360  type(timearrayserieslinktype), pointer, intent(inout) :: tasLink
361  ! -- local
362  integer(I4B) :: i, n, noder
363  real(DP) :: area
364  !
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.'
369  call store_error(errmsg)
370  call store_error_unit(taslink%TimeArraySeries%GetInunit())
371  end if
372  !
373  n = size(taslink%BndArray)
374  do i = 1, n
375  noder = taslink%nodelist(i)
376  if (noder > 0) then
377  area = this%dis%get_area(noder)
378  taslink%BndArray(i) = taslink%BndArray(i) * area
379  end if
380  end do
381  end subroutine tasmgr_convert_flux
382 
383  !> @brief Add a time arrays series link
384  !<
385  subroutine tasmgr_add_link(this, tasLink)
386  ! -- dummy
387  class(timearrayseriesmanagertype) :: this
388  type(timearrayserieslinktype), pointer :: tasLink
389  ! -- local
390  !
391  call addtimearrayserieslinktolist(this%boundTasLinks, taslink)
392  end subroutine tasmgr_add_link
393 
395 
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
integer(i4b), parameter lenmodelname
maximum length of the model name
Definition: Constants.f90:22
integer(i4b), parameter lenhugeline
maximum length of a huge line
Definition: Constants.f90:16
integer(i4b), parameter lentimeseriesname
maximum length of a time series name
Definition: Constants.f90:42
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
This module defines variable data types.
Definition: kind.f90:8
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
subroutine, public store_error_unit(iunit, terminate)
Store the file unit number.
Definition: Sim.f90:168
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
real(dp), pointer, public totimc
simulation time at start of time step
Definition: tdis.f90:33
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
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29
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.
Definition: List.f90:14