MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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  !
73  ! -- Return
74  return
75  end subroutine tasmanager_cr
76 
77  !> @brief Define the time-array series manager
78  !<
79  subroutine tasmanager_df(this)
80  ! -- dummy
81  class(timearrayseriesmanagertype) :: this
82  ! -- local
83  type(timearrayseriestype), pointer :: tasptr => null()
84  integer(I4B) :: nfiles
85  integer(I4B) :: i
86  !
87  ! -- determine how many tasfiles. This is the number of time array series
88  ! so allocate arrays to store them
89  nfiles = size(this%tasfiles)
90  allocate (this%taslist(nfiles))
91  allocate (this%tasnames(nfiles))
92  !
93  ! -- Setup a time array series for each file specified
94  do i = 1, nfiles
95  tasptr => this%taslist(i)
96  call tasptr%tas_init(this%tasfiles(i), this%modelname, &
97  this%iout, this%tasnames(i))
98  end do
99  !
100  ! -- Return
101  return
102  end subroutine tasmanager_df
103 
104  !> @brief Time step (or subtime step) advance.
105  !!
106  !! Call this each time step or subtime step.
107  !<
108  subroutine tasmgr_ad(this)
109  ! -- dummy
110  class(timearrayseriesmanagertype) :: this
111  ! -- local
112  type(timearrayserieslinktype), pointer :: tasLink => null()
113  type(timearrayseriestype), pointer :: timearrayseries => null()
114  integer(I4B) :: i, j, nlinks, nvals, isize1, isize2, inunit
115  real(DP) :: begintime, endtime
116  ! -- formats
117  character(len=*), parameter :: fmt5 = &
118  "(/,'Time-array-series controlled arrays in stress period ', &
119  &i0, ', time step ', i0, ':')"
120 10 format('"', a, '" package: ', a, ' array obtained from time-array series "', &
121  a, '"')
122  !
123  ! -- Initialize time variables
124  begintime = totimc
125  endtime = begintime + delt
126  !
127  ! -- Iterate through boundtaslinks and update specified
128  ! array with array of average values obtained from
129  ! appropriate time series.
130  if (associated(this%boundTasLinks)) then
131  nlinks = this%boundTasLinks%Count()
132  do i = 1, nlinks
133  taslink => gettimearrayserieslinkfromlist(this%boundTasLinks, i)
134  if (taslink%Iprpak == 1 .and. i == 1) then
135  write (this%iout, fmt5) kper, kstp
136  end if
137  if (taslink%UseDefaultProc) then
138  timearrayseries => taslink%timeArraySeries
139  nvals = size(taslink%BndArray)
140  !
141  ! -- Fill the package array with integrated values
142  call timearrayseries%GetAverageValues(nvals, taslink%BndArray, &
143  begintime, endtime)
144  !
145  ! -- If conversion from flux to flow is required, multiply by cell area
146  if (taslink%ConvertFlux) then
147  call this%tasmgr_convert_flux(taslink)
148  end if
149  !
150  ! -- If PRINT_INPUT is specified, write information
151  ! regarding source of time-array series data
152  if (taslink%Iprpak == 1) then
153  write (this%iout, 10) trim(taslink%PackageName), &
154  trim(taslink%Text), &
155  trim(taslink%timeArraySeries%Name)
156  end if
157  end if
158  if (i == nlinks) then
159  write (this%iout, '()')
160  end if
161  end do
162  !
163  ! -- Now that all array values have been substituted, can now multiply
164  ! an array by a multiplier array
165  do i = 1, nlinks
166  taslink => gettimearrayserieslinkfromlist(this%boundTasLinks, i)
167  if (taslink%UseDefaultProc) then
168  if (associated(taslink%RMultArray)) then
169  isize1 = size(taslink%BndArray)
170  isize2 = size(taslink%RMultArray)
171  if (isize1 == isize2 .and. isize1 == nvals) then
172  do j = 1, nvals
173  taslink%BndArray(j) = taslink%BndArray(j) * taslink%RMultArray(j)
174  end do
175  else
176  errmsg = 'Size mismatch between boundary and multiplier arrays'// &
177  ' using time-array series: '// &
178  trim(taslink%TimeArraySeries%Name)
179  call store_error(errmsg)
180  inunit = taslink%TimeArraySeries%GetInunit()
181  call store_error_unit(inunit)
182  end if
183  end if
184  end if
185  end do
186  end if
187  !
188  ! -- Return
189  return
190  end subroutine tasmgr_ad
191 
192  !> @brief Deallocate
193  !<
194  subroutine tasmgr_da(this)
195  ! -- dummy
196  class(timearrayseriesmanagertype) :: this
197  ! -- local
198  integer :: i, n
199  type(timearrayserieslinktype), pointer :: tasLink => null()
200  !
201  ! -- Deallocate contents of each TimeArraySeriesType object in list
202  ! of time-array series links.
203  n = this%boundTasLinks%Count()
204  do i = 1, n
205  taslink => gettimearrayserieslinkfromlist(this%boundTasLinks, i)
206  call taslink%da()
207  end do
208  !
209  ! -- Go through and deallocate individual time array series
210  do i = 1, size(this%taslist)
211  call this%taslist(i)%da()
212  end do
213  !
214  ! -- Deallocate the list of time-array series links.
215  call this%boundTasLinks%Clear(.true.)
216  deallocate (this%boundTasLinks)
217  deallocate (this%tasfiles)
218  !
219  ! -- Deallocate the time array series
220  deallocate (this%taslist)
221  deallocate (this%tasnames)
222  !
223  ! -- nullify pointers
224  this%dis => null()
225  this%boundTasLinks => null()
226  !
227  ! -- Return
228  return
229  end subroutine tasmgr_da
230 
231  !> @brief Add a time-array series file
232  !<
233  subroutine add_tasfile(this, fname)
234  ! -- modules
236  ! -- dummy
237  class(timearrayseriesmanagertype) :: this
238  character(len=*), intent(in) :: fname
239  ! -- local
240  integer(I4B) :: indx
241  !
242  call expandarray(this%tasfiles, 1)
243  indx = size(this%tasfiles)
244  this%tasfiles(indx) = fname
245  !
246  ! -- Return
247  return
248  end subroutine add_tasfile
249 
250  !> @brief Zero out arrays that are represented with time series
251  !!
252  !! Delete all existing links from time array series to package arrays as they
253  !! will need to be created with a new BEGIN PERIOD block
254  !<
255  subroutine reset(this, pkgName)
256  ! -- dummy
257  class(timearrayseriesmanagertype) :: this
258  character(len=*), intent(in) :: pkgName
259  ! -- local
260  integer(I4B) :: i, j, nlinks
261  type(timearrayserieslinktype), pointer :: taslink
262  !
263  ! -- Reassign all linked elements to zero
264  nlinks = this%boundTasLinks%Count()
265  do i = 1, nlinks
266  taslink => gettimearrayserieslinkfromlist(this%boundTasLinks, i)
267  if (associated(taslink)) then
268  do j = 1, size(taslink%BndArray)
269  taslink%BndArray(j) = dzero
270  end do
271  end if
272  end do
273  !
274  ! -- Delete all existing time array links
275  if (associated(this%boundTasLinks)) then
276  ! Deallocate and remove all links belonging to package
277  nlinks = this%boundTasLinks%Count()
278  do i = nlinks, 1, -1
279  taslink => gettimearrayserieslinkfromlist(this%boundTasLinks, i)
280  if (associated(taslink)) then
281  call taslink%da()
282  call this%boundTasLinks%RemoveNode(i, .true.)
283  end if
284  end do
285  end if
286  !
287  ! -- Return
288  return
289  end subroutine reset
290 
291  !> @brief Make link from time-array series to package array
292  !<
293  subroutine maketaslink(this, pkgName, bndArray, iprpak, &
294  tasName, text, convertFlux, nodelist, inunit)
295  ! -- dummy
296  class(timearrayseriesmanagertype) :: this
297  character(len=*), intent(in) :: pkgName
298  real(DP), dimension(:), pointer :: bndArray
299  integer(I4B), intent(in) :: iprpak
300  character(len=*), intent(in) :: tasName
301  character(len=*), intent(in) :: text
302  logical, intent(in) :: convertFlux
303  integer(I4B), dimension(:), pointer, contiguous, intent(in) :: nodelist
304  integer(I4B), intent(in) :: inunit
305  ! -- local
306  integer(I4B) :: i, nfiles, iloc
307  character(LINELENGTH) :: ermsg
308  type(timearrayserieslinktype), pointer :: newTasLink
309  type(timearrayseriestype), pointer :: tasptr => null()
310  !
311  ! -- Find the time array series
312  nfiles = size(this%tasnames)
313  iloc = 0
314  do i = 1, nfiles
315  if (this%tasnames(i) == tasname) then
316  iloc = i
317  exit
318  end if
319  end do
320  if (iloc == 0) then
321  ermsg = 'Error: Time-array series "'//trim(tasname)//'" not found.'
322  call store_error(ermsg)
323  call store_error_unit(inunit)
324  end if
325  tasptr => this%taslist(iloc)
326  !
327  ! -- Construct a time-array series link
328  newtaslink => null()
329  call constructtimearrayserieslink(newtaslink, tasptr, &
330  pkgname, bndarray, iprpak, &
331  text)
332  newtaslink%ConvertFlux = convertflux
333  newtaslink%nodelist => nodelist
334  !
335  ! -- Add link to list of links
336  call this%tasmgr_add_link(newtaslink)
337  !
338  ! -- Return
339  return
340  end subroutine maketaslink
341 
342  !> @brief Get link from the boundtaslinks list
343  !<
344  function getlink(this, indx) result(tasLink)
345  ! -- dummy
346  class(timearrayseriesmanagertype) :: this
347  integer(I4B), intent(in) :: indx
348  ! -- return
349  type(timearrayserieslinktype), pointer :: taslink
350  !
351  taslink => null()
352  !
353  if (associated(this%boundTasLinks)) then
354  taslink => gettimearrayserieslinkfromlist(this%boundTasLinks, indx)
355  end if
356  !
357  ! -- Return
358  return
359  end function getlink
360 
361  !> @brief Count number of links in the boundtaslinks list
362  !<
363  function countlinks(this)
364  ! -- return
365  integer(I4B) :: countlinks
366  ! -- dummy
367  class(timearrayseriesmanagertype) :: this
368  !
369  if (associated(this%boundtaslinks)) then
370  countlinks = this%boundTasLinks%Count()
371  else
372  countlinks = 0
373  end if
374  !
375  ! -- Return
376  return
377  end function countlinks
378 
379  ! -- Private procedures
380 
381  !> @brief Convert the array from a flux to a flow rate by multiplying by the
382  !! cell area
383  !<
384  subroutine tasmgr_convert_flux(this, tasLink)
385  ! -- dummy
386  class(timearrayseriesmanagertype) :: this
387  type(timearrayserieslinktype), pointer, intent(inout) :: tasLink
388  ! -- local
389  integer(I4B) :: i, n, noder
390  real(DP) :: area
391  !
392  if (.not. (associated(this%dis) .and. &
393  associated(taslink%nodelist))) then
394  errmsg = 'Programming error. Cannot convert flux. Verify that '&
395  &'a valid DIS instance and nodelist were provided.'
396  call store_error(errmsg)
397  call store_error_unit(taslink%TimeArraySeries%GetInunit())
398  end if
399  !
400  n = size(taslink%BndArray)
401  do i = 1, n
402  noder = taslink%nodelist(i)
403  if (noder > 0) then
404  area = this%dis%get_area(noder)
405  taslink%BndArray(i) = taslink%BndArray(i) * area
406  end if
407  end do
408  !
409  ! -- Return
410  return
411  end subroutine tasmgr_convert_flux
412 
413  !> @brief Add a time arrays series link
414  !<
415  subroutine tasmgr_add_link(this, tasLink)
416  ! -- dummy
417  class(timearrayseriesmanagertype) :: this
418  type(timearrayserieslinktype), pointer :: tasLink
419  ! -- local
420  !
421  call addtimearrayserieslinktolist(this%boundTasLinks, taslink)
422  !
423  ! -- Return
424  return
425  end subroutine tasmgr_add_link
426 
428 
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
integer(i4b), parameter lenmodelname
maximum length of the model name
Definition: Constants.f90:21
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:41
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
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:10