MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
TimeSeries.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b
8  use mathutilmodule, only: is_close
10  use listmodule, only: listtype, listnodetype
11  use simvariablesmodule, only: errmsg
12  use simmodule, only: count_errors, store_error, &
18 
19  private
24 
26  ! -- Public members
27  integer(I4B), public :: imethod = undefined
28  character(len=LENTIMESERIESNAME), public :: name = ''
29  ! -- Private members
30  real(dp), private :: sfac = done
31  logical, public :: autodeallocate = .true.
32  type(listtype), pointer, private :: list => null()
33  class(timeseriesfiletype), pointer, private :: tsfile => null()
34 
35  contains
36 
37  ! -- Public procedures
38  procedure, public :: addtimeseriesrecord
39  procedure, public :: clear
40  procedure, public :: findlatesttime
41  procedure, public :: get_surrounding_records
42  procedure, public :: get_surrounding_nodes
43  procedure, public :: getcurrenttimeseriesrecord
44  procedure, public :: getnexttimeseriesrecord
45  procedure, public :: getprevioustimeseriesrecord
46  procedure, public :: gettimeseriesrecord
47  procedure, public :: getvalue
48  procedure, public :: initializetimeseries => initialize_time_series
49  procedure, public :: inserttsr
50  procedure, public :: reset
51  ! -- Private procedures
52  procedure, private :: da => ts_da
53  procedure, private :: get_average_value
54  procedure, private :: get_integrated_value
55  procedure, private :: get_latest_preceding_node
56  procedure, private :: get_value_at_time
57  procedure, private :: initialize_time_series
58  procedure, private :: read_next_record
59  end type timeseriestype
60 
62  ! -- Private members
63  integer(I4B), public :: inunit = 0
64  integer(I4B), public :: iout = 0
65  integer(I4B), public :: ntimeseries = 0
66  logical, public :: finishedreading = .false.
67  character(len=LINELENGTH), public :: datafile = ''
68  type(timeseriestype), dimension(:), &
69  pointer, contiguous, public :: timeseries => null()
70  type(blockparsertype), pointer, public :: parser
71 
72  contains
73 
74  ! -- Public procedures
75  procedure, public :: count
76  procedure, public :: initializetsfile
77  procedure, public :: gettimeseries
78  procedure, public :: da => tsf_da
79  ! -- Private procedures
80  procedure, private :: read_tsfile_line
81  end type timeseriesfiletype
82 
84  ! -- Public members
85  type(timeseriestype), pointer, public :: timeseries => null()
87 
88 contains
89 
90  ! -- non-type-bound procedures
91 
92  !> @brief Construct time series file
93  !<
94  subroutine constructtimeseriesfile(newTimeSeriesFile)
95  ! -- dummy
96  type(timeseriesfiletype), pointer, intent(inout) :: newtimeseriesfile
97  !
98  allocate (newtimeseriesfile)
99  allocate (newtimeseriesfile%parser)
100  !
101  ! -- Return
102  return
103  end subroutine constructtimeseriesfile
104 
105  !> @brief Cast an unlimited polymorphic object as class(TimeSeriesFileType)
106  !<
107  function castastimeseriesfiletype(obj) result(res)
108  ! -- dummy
109  class(*), pointer, intent(inout) :: obj
110  ! -- return
111  type(timeseriesfiletype), pointer :: res
112  !
113  res => null()
114  if (.not. associated(obj)) return
115  !
116  select type (obj)
117  type is (timeseriesfiletype)
118  res => obj
119  end select
120  !
121  ! -- Return
122  return
123  end function castastimeseriesfiletype
124 
125  !> @brief Cast an unlimited polymorphic object as class(TimeSeriesFileType)
126  !<
127  function castastimeseriesfileclass(obj) result(res)
128  ! -- dummy
129  class(*), pointer, intent(inout) :: obj
130  ! -- return
131  type(timeseriesfiletype), pointer :: res
132  !
133  res => null()
134  if (.not. associated(obj)) return
135  !
136  select type (obj)
137  class is (timeseriesfiletype)
138  res => obj
139  end select
140  !
141  ! -- Return
142  return
143  end function castastimeseriesfileclass
144 
145  !> @brief Add time series file to list
146  !<
147  subroutine addtimeseriesfiletolist(list, tsfile)
148  ! -- dummy
149  type(listtype), intent(inout) :: list
150  class(timeseriesfiletype), pointer, intent(inout) :: tsfile
151  ! -- local
152  class(*), pointer :: obj => null()
153  !
154  obj => tsfile
155  call list%Add(obj)
156  !
157  ! -- Return
158  return
159  end subroutine addtimeseriesfiletolist
160 
161  !> @brief Get time series from list
162  !<
163  function gettimeseriesfilefromlist(list, idx) result(res)
164  ! -- dummy
165  type(listtype), intent(inout) :: list
166  integer(I4B), intent(in) :: idx
167  ! -- return
168  type(timeseriesfiletype), pointer :: res
169  ! -- local
170  class(*), pointer :: obj => null()
171  !
172  obj => list%GetItem(idx)
173  res => castastimeseriesfiletype(obj)
174  !
175  if (.not. associated(res)) then
176  res => castastimeseriesfileclass(obj)
177  end if
178  !
179  ! -- Return
180  return
181  end function gettimeseriesfilefromlist
182 
183  !> @brief Compare two time series; if they are identical, return true
184  !<
185  function sametimeseries(ts1, ts2) result(same)
186  ! -- dummy
187  type(timeseriestype), intent(in) :: ts1
188  type(timeseriestype), intent(in) :: ts2
189  ! -- return
190  logical :: same
191  ! -- local
192  integer :: i, n1, n2
193  type(timeseriesrecordtype), pointer :: tsr1, tsr2
194  !
195  same = .false.
196  n1 = ts1%list%Count()
197  n2 = ts2%list%Count()
198  if (n1 /= n2) return
199  !
200  call ts1%Reset()
201  call ts2%Reset()
202  !
203  do i = 1, n1
204  tsr1 => ts1%GetNextTimeSeriesRecord()
205  tsr2 => ts2%GetNextTimeSeriesRecord()
206  if (tsr1%tsrTime /= tsr2%tsrTime) return
207  if (tsr1%tsrValue /= tsr2%tsrValue) return
208  end do
209  !
210  same = .true.
211  !
212  ! -- Return
213  return
214  end function sametimeseries
215 
216  ! Type-bound procedures of TimeSeriesType
217 
218  !> @brief Get time series value
219  !!
220  !! If iMethod is STEPWISE or LINEAR:
221  !! Return a time-weighted average value for a specified time span.
222  !! If iMethod is LINEAREND:
223  !! Return value at time1. Time0 argument is ignored.
224  !! Units: (ts-value-unit)
225  !<
226  function getvalue(this, time0, time1, extendToEndOfSimulation)
227  ! -- return
228  real(dp) :: getvalue
229  ! -- dummy
230  class(timeseriestype), intent(inout) :: this
231  real(dp), intent(in) :: time0
232  real(dp), intent(in) :: time1
233  logical, intent(in), optional :: extendtoendofsimulation
234  ! -- local
235  logical :: extend
236  !
237  if (present(extendtoendofsimulation)) then
238  extend = extendtoendofsimulation
239  else
240  extend = .false.
241  end if
242  !
243  select case (this%iMethod)
244  case (stepwise, linear)
245  getvalue = this%get_average_value(time0, time1, extend)
246  case (linearend)
247  getvalue = this%get_value_at_time(time1, extend)
248  end select
249  !
250  ! -- Return
251  return
252  end function getvalue
253 
254  !> @brief Initialize time series
255  !!
256  !! Open time-series file and read options and first time-series record.
257  !<
258  subroutine initialize_time_series(this, tsfile, name, autoDeallocate)
259  ! -- dummy
260  class(timeseriestype), intent(inout) :: this
261  class(timeseriesfiletype), target :: tsfile
262  character(len=*), intent(in) :: name
263  logical, intent(in), optional :: autoDeallocate
264  ! -- local
265  character(len=LENTIMESERIESNAME) :: tsNameTemp
266  !
267  ! -- Assign the time-series tsfile, name, and autoDeallocate
268  this%tsfile => tsfile
269  ! Store time-series name as all caps
270  tsnametemp = name
271  call upcase(tsnametemp)
272  this%Name = tsnametemp
273  !
274  this%iMethod = undefined
275  !
276  if (present(autodeallocate)) this%autoDeallocate = autodeallocate
277  !
278  ! -- allocate the list
279  allocate (this%list)
280  !
281  ! -- ensure that NAME has been specified
282  if (this%Name == '') then
283  errmsg = 'Name not specified for time series.'
284  call store_error(errmsg, terminate=.true.)
285  end if
286  !
287  ! -- Return
288  return
289  end subroutine initialize_time_series
290 
291  !> @brief Get surrounding records
292  !<
293  subroutine get_surrounding_records(this, time, tsrecEarlier, tsrecLater)
294  ! -- dummy
295  class(timeseriestype), intent(inout) :: this
296  real(DP), intent(in) :: time
297  type(timeseriesrecordtype), pointer, intent(inout) :: tsrecEarlier
298  type(timeseriesrecordtype), pointer, intent(inout) :: tsrecLater
299  ! -- local
300  real(DP) :: time0, time1
301  type(listnodetype), pointer :: currNode => null()
302  type(listnodetype), pointer :: tsNode0 => null()
303  type(listnodetype), pointer :: tsNode1 => null()
304  type(timeseriesrecordtype), pointer :: tsr => null(), tsrec0 => null()
305  type(timeseriesrecordtype), pointer :: tsrec1 => null()
306  class(*), pointer :: obj => null()
307  !
308  tsrecearlier => null()
309  tsreclater => null()
310  !
311  if (associated(this%list%firstNode)) then
312  currnode => this%list%firstNode
313  end if
314  !
315  ! -- If the next node is earlier than time of interest, advance along
316  ! linked list until the next node is later than time of interest.
317  do
318  if (associated(currnode)) then
319  if (associated(currnode%nextNode)) then
320  obj => currnode%nextNode%GetItem()
321  tsr => castastimeseriesrecordtype(obj)
322  if (tsr%tsrTime < time .and. .not. is_close(tsr%tsrTime, time)) then
323  currnode => currnode%nextNode
324  else
325  exit
326  end if
327  else
328  ! -- read another record
329  if (.not. this%read_next_record()) exit
330  end if
331  else
332  exit
333  end if
334  end do
335  !
336  if (associated(currnode)) then
337  !
338  ! -- find earlier record
339  tsnode0 => currnode
340  obj => tsnode0%GetItem()
341  tsrec0 => castastimeseriesrecordtype(obj)
342  time0 = tsrec0%tsrTime
343  do while (time0 > time)
344  if (associated(tsnode0%prevNode)) then
345  tsnode0 => tsnode0%prevNode
346  obj => tsnode0%GetItem()
347  tsrec0 => castastimeseriesrecordtype(obj)
348  time0 = tsrec0%tsrTime
349  else
350  exit
351  end if
352  end do
353  !
354  ! -- find later record
355  tsnode1 => currnode
356  obj => tsnode1%GetItem()
357  tsrec1 => castastimeseriesrecordtype(obj)
358  time1 = tsrec1%tsrTime
359  do while (time1 < time .and. .not. is_close(time1, time))
360  if (associated(tsnode1%nextNode)) then
361  tsnode1 => tsnode1%nextNode
362  obj => tsnode1%GetItem()
363  tsrec1 => castastimeseriesrecordtype(obj)
364  time1 = tsrec1%tsrTime
365  else
366  ! -- get next record
367  if (.not. this%read_next_record()) then
368  ! -- end of file reached, so exit loop
369  exit
370  end if
371  end if
372  end do
373  !
374  end if
375  !
376  if (time0 < time .or. is_close(time0, time)) tsrecearlier => tsrec0
377  if (time1 > time .or. is_close(time1, time)) tsreclater => tsrec1
378  !
379  ! -- Return
380  return
381  end subroutine get_surrounding_records
382 
383  !> @brief Get surrounding nodes
384  !!
385  !! This subroutine is for working with time series already entirely stored
386  !! in memory -- it does not read data from a file.
387  !<
388  subroutine get_surrounding_nodes(this, time, nodeEarlier, nodeLater)
389  ! -- dummy
390  class(timeseriestype), intent(inout) :: this
391  real(DP), intent(in) :: time
392  type(listnodetype), pointer, intent(inout) :: nodeEarlier
393  type(listnodetype), pointer, intent(inout) :: nodeLater
394  ! -- local
395  real(DP) :: time0, time1
396  type(listnodetype), pointer :: currNode => null()
397  type(listnodetype), pointer :: tsNode0 => null()
398  type(listnodetype), pointer :: tsNode1 => null()
399  type(timeseriesrecordtype), pointer :: tsr => null(), tsrec0 => null()
400  type(timeseriesrecordtype), pointer :: tsrec1 => null()
401  type(timeseriesrecordtype), pointer :: tsrecEarlier
402  type(timeseriesrecordtype), pointer :: tsrecLater
403  class(*), pointer :: obj => null()
404  !
405  tsrecearlier => null()
406  tsreclater => null()
407  nodeearlier => null()
408  nodelater => null()
409  !
410  if (associated(this%list%firstNode)) then
411  currnode => this%list%firstNode
412  end if
413  !
414  ! -- If the next node is earlier than time of interest, advance along
415  ! linked list until the next node is later than time of interest.
416  do
417  if (associated(currnode)) then
418  if (associated(currnode%nextNode)) then
419  obj => currnode%nextNode%GetItem()
420  tsr => castastimeseriesrecordtype(obj)
421  if (tsr%tsrTime < time .and. .not. is_close(tsr%tsrTime, time)) then
422  currnode => currnode%nextNode
423  else
424  exit
425  end if
426  else
427  exit
428  end if
429  else
430  exit
431  end if
432  end do
433  !
434  if (associated(currnode)) then
435  !
436  ! -- find earlier record
437  tsnode0 => currnode
438  obj => tsnode0%GetItem()
439  tsrec0 => castastimeseriesrecordtype(obj)
440  time0 = tsrec0%tsrTime
441  do while (time0 > time)
442  if (associated(tsnode0%prevNode)) then
443  tsnode0 => tsnode0%prevNode
444  obj => tsnode0%GetItem()
445  tsrec0 => castastimeseriesrecordtype(obj)
446  time0 = tsrec0%tsrTime
447  else
448  exit
449  end if
450  end do
451  !
452  ! -- find later record
453  tsnode1 => currnode
454  obj => tsnode1%GetItem()
455  tsrec1 => castastimeseriesrecordtype(obj)
456  time1 = tsrec1%tsrTime
457  do while (time1 < time .and. .not. is_close(time1, time))
458  if (associated(tsnode1%nextNode)) then
459  tsnode1 => tsnode1%nextNode
460  obj => tsnode1%GetItem()
461  tsrec1 => castastimeseriesrecordtype(obj)
462  time1 = tsrec1%tsrTime
463  else
464  exit
465  end if
466  end do
467  !
468  end if
469  !
470  if (time0 < time .or. is_close(time0, time)) then
471  tsrecearlier => tsrec0
472  nodeearlier => tsnode0
473  end if
474  if (time1 > time .or. is_close(time1, time)) then
475  tsreclater => tsrec1
476  nodelater => tsnode1
477  end if
478  !
479  ! -- Return
480  return
481  end subroutine get_surrounding_nodes
482 
483  !> @brief Read next record
484  !!
485  !! Read next time-series record from input file
486  !<
487  logical function read_next_record(this)
488  ! -- dummy
489  class(timeseriestype), intent(inout) :: this
490  !
491  ! -- If we have already encountered the end of the TIMESERIES block, do not try to read any further
492  if (this%tsfile%finishedReading) then
493  read_next_record = .false.
494  return
495  end if
496  !
497  read_next_record = this%tsfile%read_tsfile_line()
498  if (.not. read_next_record) then
499  this%tsfile%finishedReading = .true.
500  end if
501  !
502  ! -- Return
503  return
504  end function read_next_record
505 
506  !> @brief Get value for a time
507  !!
508  !! Return a value for a specified time, same units as time-series values
509  !<
510  function get_value_at_time(this, time, extendToEndOfSimulation)
511  ! -- return
512  real(dp) :: get_value_at_time
513  ! -- dummy
514  class(timeseriestype), intent(inout) :: this
515  real(dp), intent(in) :: time ! time of interest
516  logical, intent(in) :: extendtoendofsimulation
517  ! -- local
518  integer(I4B) :: ierr
519  real(dp) :: ratio, time0, time1, timediff, timediffi, val0, val1, &
520  valdiff
521  type(timeseriesrecordtype), pointer :: tsrearlier => null()
522  type(timeseriesrecordtype), pointer :: tsrlater => null()
523  ! -- formats
524 10 format('Error getting value at time ', g10.3, ' for time series "', a, '"')
525  !
526  ierr = 0
527  call this%get_surrounding_records(time, tsrearlier, tsrlater)
528  if (associated(tsrearlier)) then
529  if (associated(tsrlater)) then
530  ! -- values are available for both earlier and later times
531  if (this%iMethod == stepwise) then
532  get_value_at_time = tsrearlier%tsrValue
533  elseif (this%iMethod == linear .or. this%iMethod == linearend) then
534  ! -- For get_value_at_time, result is the same for either
535  ! linear method.
536  ! -- Perform linear interpolation.
537  time0 = tsrearlier%tsrTime
538  time1 = tsrlater%tsrtime
539  timediff = time1 - time0
540  timediffi = time - time0
541  if (timediff > 0) then
542  ratio = timediffi / timediff
543  else
544  ! -- should not happen if TS does not contain duplicate times
545  ratio = 0.5d0
546  end if
547  val0 = tsrearlier%tsrValue
548  val1 = tsrlater%tsrValue
549  valdiff = val1 - val0
550  get_value_at_time = val0 + (ratio * valdiff)
551  else
552  ierr = 1
553  end if
554  else
555  if (extendtoendofsimulation .or. is_close(tsrearlier%tsrTime, time)) then
556  get_value_at_time = tsrearlier%tsrValue
557  else
558  ! -- Only earlier time is available, and it is not time of interest;
559  ! however, if method is STEPWISE, use value for earlier time.
560  if (this%iMethod == stepwise) then
561  get_value_at_time = tsrearlier%tsrValue
562  else
563  ierr = 1
564  end if
565  end if
566  end if
567  else
568  if (associated(tsrlater)) then
569  if (is_close(tsrlater%tsrTime, time)) then
570  get_value_at_time = tsrlater%tsrValue
571  else
572  ! -- only later time is available, and it is not time of interest
573  ierr = 1
574  end if
575  else
576  ! -- Neither earlier nor later time is available.
577  ! This should never happen!
578  ierr = 1
579  end if
580  end if
581  !
582  if (ierr > 0) then
583  write (errmsg, 10) time, trim(this%Name)
584  call store_error(errmsg, terminate=.true.)
585  end if
586  !
587  ! -- Return
588  return
589  end function get_value_at_time
590 
591  !> @brief Get integrated value
592  !!
593  !! Return an integrated value for a specified time span.
594  !! Units: (ts-value-unit)*time
595  !<
596  function get_integrated_value(this, time0, time1, extendToEndOfSimulation)
597  ! -- return
598  real(dp) :: get_integrated_value
599  ! -- dummy
600  class(timeseriestype), intent(inout) :: this
601  real(dp), intent(in) :: time0
602  real(dp), intent(in) :: time1
603  logical, intent(in) :: extendtoendofsimulation
604  ! -- local
605  real(dp) :: area, currtime, nexttime, ratio0, ratio1, t0, t01, t1, &
606  timediff, value, value0, value1, valuediff, currval, nextval
607  logical :: ldone, lprocess
608  type(listnodetype), pointer :: tslnodepreceding => null()
609  type(listnodetype), pointer :: currnode => null(), nextnode => null()
610  type(timeseriesrecordtype), pointer :: currrecord => null()
611  type(timeseriesrecordtype), pointer :: nextrecord => null()
612  class(*), pointer :: currobj => null(), nextobj => null()
613  ! -- formats
614 10 format('Error encountered while performing integration', &
615  ' for time series "', a, '" for time interval: ', g12.5, ' to ', g12.5)
616  !
617  value = dzero
618  ldone = .false.
619  t1 = -done
620  call this%get_latest_preceding_node(time0, tslnodepreceding)
621  if (associated(tslnodepreceding)) then
622  currnode => tslnodepreceding
623  do while (.not. ldone)
624  currobj => currnode%GetItem()
625  currrecord => castastimeseriesrecordtype(currobj)
626  currtime = currrecord%tsrTime
627  if (is_close(currtime, time1)) then
628  ! Current node time = time1 so should be ldone
629  ldone = .true.
630  elseif (currtime < time1) then
631  if (.not. associated(currnode%nextNode)) then
632  ! -- try to read the next record
633  if (.not. this%read_next_record()) then
634  if (.not. extendtoendofsimulation) then
635  write (errmsg, 10) trim(this%Name), time0, time1
636  call store_error(errmsg, terminate=.true.)
637  end if
638  end if
639  end if
640  !
641  currval = currrecord%tsrValue
642  lprocess = .false.
643  if (associated(currnode%nextNode)) then
644  nextnode => currnode%nextNode
645  nextobj => nextnode%GetItem()
646  nextrecord => castastimeseriesrecordtype(nextobj)
647  nexttime = nextrecord%tsrTime
648  nextval = nextrecord%tsrValue
649  lprocess = .true.
650  elseif (extendtoendofsimulation) then
651  ! -- Last time series value extends forever, so integrate the final value over all simulation time after the end of the series
652  nexttime = time1
653  nextval = currval
654  lprocess = .true.
655  end if
656  !
657  if (lprocess) then
658  ! -- determine lower and upper limits of time span of interest
659  ! within current interval
660  if (currtime > time0 .or. is_close(currtime, time0)) then
661  t0 = currtime
662  else
663  t0 = time0
664  end if
665  if (nexttime < time1 .or. is_close(nexttime, time1)) then
666  t1 = nexttime
667  else
668  t1 = time1
669  end if
670  ! -- find area of rectangle or trapezoid delimited by t0 and t1
671  t01 = t1 - t0
672  select case (this%iMethod)
673  case (stepwise)
674  ! -- compute area of a rectangle
675  value0 = currval
676  area = value0 * t01
677  case (linear, linearend)
678  ! -- compute area of a trapezoid
679  timediff = nexttime - currtime
680  ratio0 = (t0 - currtime) / timediff
681  ratio1 = (t1 - currtime) / timediff
682  valuediff = nextval - currval
683  value0 = currval + ratio0 * valuediff
684  value1 = currval + ratio1 * valuediff
685  if (this%iMethod == linear) then
686  area = 0.5d0 * t01 * (value0 + value1)
687  elseif (this%iMethod == linearend) then
688  area = dzero
689  value = value1
690  end if
691  end select
692  ! -- add area to integrated value
693  value = value + area
694  end if
695  end if
696  !
697  ! -- Are we done yet?
698  if (t1 > time1) then
699  ldone = .true.
700  elseif (is_close(t1, time1)) then
701  ldone = .true.
702  else
703  ! -- We are not done yet
704  if (.not. associated(currnode%nextNode)) then
705  ! -- Not done and no more data, so try to read the next record
706  if (.not. this%read_next_record()) then
707  write (errmsg, 10) trim(this%Name), time0, time1
708  call store_error(errmsg, terminate=.true.)
709  end if
710  elseif (associated(currnode%nextNode)) then
711  currnode => currnode%nextNode
712  end if
713  end if
714  end do
715  end if
716  !
717  get_integrated_value = value
718  if (this%autoDeallocate) then
719  if (associated(tslnodepreceding)) then
720  if (associated(tslnodepreceding%prevNode)) then
721  call this%list%DeallocateBackward(tslnodepreceding%prevNode)
722  end if
723  end if
724  end if
725  !
726  ! -- Return
727  return
728  end function get_integrated_value
729 
730  !> @brief Get average value
731  !!
732  !! Return a time-weighted average value for a specified time span.
733  !! Units: (ts-value-unit)
734  !<
735  function get_average_value(this, time0, time1, extendToEndOfSimulation)
736  ! -- return
737  real(dp) :: get_average_value
738  ! -- dummy
739  class(timeseriestype), intent(inout) :: this
740  real(dp), intent(in) :: time0
741  real(dp), intent(in) :: time1
742  logical, intent(in) :: extendtoendofsimulation
743  ! -- local
744  real(dp) :: timediff, value, valueintegrated
745  !
746  timediff = time1 - time0
747  if (timediff > 0) then
748  valueintegrated = this%get_integrated_value(time0, time1, &
749  extendtoendofsimulation)
750  if (this%iMethod == linearend) then
751  value = valueintegrated
752  else
753  value = valueintegrated / timediff
754  end if
755  else
756  ! -- time0 and time1 are the same
757  value = this%get_value_at_time(time0, extendtoendofsimulation)
758  end if
759  get_average_value = value
760  !
761  ! -- Return
762  return
763  end function get_average_value
764 
765  !> @brief Get latest preceding node
766  !!
767  !! Return pointer to ListNodeType object for the node representing the
768  !! latest preceding time in the time series
769  !<
770  subroutine get_latest_preceding_node(this, time, tslNode)
771  ! -- dummy
772  class(timeseriestype), intent(inout) :: this
773  real(DP), intent(in) :: time
774  type(listnodetype), pointer, intent(inout) :: tslNode
775  ! -- local
776  real(DP) :: time0
777  type(listnodetype), pointer :: currNode => null()
778  type(listnodetype), pointer :: tsNode0 => null()
779  type(timeseriesrecordtype), pointer :: tsr => null()
780  type(timeseriesrecordtype), pointer :: tsrec0 => null()
781  class(*), pointer :: obj => null()
782  !
783  tslnode => null()
784  if (associated(this%list%firstNode)) then
785  currnode => this%list%firstNode
786  else
787  call store_error('probable programming error in &
788  &get_latest_preceding_node', &
789  terminate=.true.)
790  end if
791  !
792  ! -- If the next node is earlier than time of interest, advance along
793  ! linked list until the next node is later than time of interest.
794  do
795  if (associated(currnode)) then
796  if (associated(currnode%nextNode)) then
797  obj => currnode%nextNode%GetItem()
798  tsr => castastimeseriesrecordtype(obj)
799  if (tsr%tsrTime < time .or. is_close(tsr%tsrTime, time)) then
800  currnode => currnode%nextNode
801  else
802  exit
803  end if
804  else
805  ! -- read another record
806  if (.not. this%read_next_record()) exit
807  end if
808  else
809  exit
810  end if
811  end do
812  !
813  if (associated(currnode)) then
814  !
815  ! -- find earlier record
816  tsnode0 => currnode
817  obj => tsnode0%GetItem()
818  tsrec0 => castastimeseriesrecordtype(obj)
819  time0 = tsrec0%tsrTime
820  do while (time0 > time)
821  if (associated(tsnode0%prevNode)) then
822  tsnode0 => tsnode0%prevNode
823  obj => tsnode0%GetItem()
824  tsrec0 => castastimeseriesrecordtype(obj)
825  time0 = tsrec0%tsrTime
826  else
827  exit
828  end if
829  end do
830  end if
831  !
832  if (time0 < time .or. is_close(time0, time)) tslnode => tsnode0
833  !
834  ! -- Return
835  return
836  end subroutine get_latest_preceding_node
837 
838  !> @brief Deallocate
839  !<
840  subroutine ts_da(this)
841  ! -- dummy
842  class(timeseriestype), intent(inout) :: this
843  !
844  if (associated(this%list)) then
845  call this%list%Clear(.true.)
846  deallocate (this%list)
847  end if
848  !
849  ! -- Return
850  return
851  end subroutine ts_da
852 
853  !> @brief Add ts record
854  !<
855  subroutine addtimeseriesrecord(this, tsr)
856  ! -- dummy
857  class(timeseriestype) :: this
858  type(timeseriesrecordtype), pointer, intent(inout) :: tsr
859  ! -- local
860  class(*), pointer :: obj => null()
861  !
862  obj => tsr
863  call this%list%Add(obj)
864  !
865  ! -- Return
866  return
867  end subroutine addtimeseriesrecord
868 
869  !> @brief Get current ts record
870  !<
871  function getcurrenttimeseriesrecord(this) result(res)
872  ! -- dummy
873  class(timeseriestype) :: this
874  ! -- result
875  type(timeseriesrecordtype), pointer :: res
876  ! -- local
877  class(*), pointer :: obj => null()
878  !
879  obj => null()
880  res => null()
881  obj => this%list%GetItem()
882  if (associated(obj)) then
883  res => castastimeseriesrecordtype(obj)
884  end if
885  !
886  ! -- Return
887  return
888  end function getcurrenttimeseriesrecord
889 
890  !> @brief Get previous ts record
891  !<
892  function getprevioustimeseriesrecord(this) result(res)
893  ! -- dummy
894  class(timeseriestype) :: this
895  ! -- result
896  type(timeseriesrecordtype), pointer :: res
897  ! -- local
898  class(*), pointer :: obj => null()
899  !
900  obj => null()
901  res => null()
902  obj => this%list%GetPreviousItem()
903  if (associated(obj)) then
904  res => castastimeseriesrecordtype(obj)
905  end if
906  !
907  ! -- Return
908  return
909  end function getprevioustimeseriesrecord
910 
911  !> @brief Get next ts record
912  !<
913  function getnexttimeseriesrecord(this) result(res)
914  ! -- dummy
915  class(timeseriestype) :: this
916  ! -- result
917  type(timeseriesrecordtype), pointer :: res
918  ! -- local
919  class(*), pointer :: obj => null()
920  !
921  obj => null()
922  res => null()
923  obj => this%list%GetNextItem()
924  if (associated(obj)) then
925  res => castastimeseriesrecordtype(obj)
926  end if
927  !
928  ! -- Return
929  return
930  end function getnexttimeseriesrecord
931 
932  !> @brief Get ts record
933  !<
934  function gettimeseriesrecord(this, time, epsi) result(res)
935  ! -- dummy
936  class(timeseriestype) :: this
937  double precision, intent(in) :: time
938  double precision, intent(in) :: epsi
939  ! -- result
940  type(timeseriesrecordtype), pointer :: res
941  ! -- local
942  type(timeseriesrecordtype), pointer :: tsr
943  !
944  call this%list%Reset()
945  res => null()
946  do
947  tsr => this%GetNextTimeSeriesRecord()
948  if (associated(tsr)) then
949  if (is_close(tsr%tsrTime, time)) then
950  res => tsr
951  exit
952  end if
953  if (tsr%tsrTime > time) exit
954  else
955  exit
956  end if
957  end do
958  !
959  ! -- Return
960  return
961  end function gettimeseriesrecord
962 
963  !> @brief Reset
964  !<
965  subroutine reset(this)
966  ! -- dummy
967  class(timeseriestype) :: this
968  !
969  call this%list%Reset()
970  !
971  ! -- Return
972  return
973  end subroutine reset
974 
975  !> @brief Insert a time series record
976  !<
977  subroutine inserttsr(this, tsr)
978  ! -- dummy
979  class(timeseriestype), intent(inout) :: this
980  type(timeseriesrecordtype), pointer, intent(inout) :: tsr
981  ! -- local
982  double precision :: badtime, time, time0, time1
983  type(timeseriesrecordtype), pointer :: tsrEarlier, tsrLater
984  type(listnodetype), pointer :: nodeEarlier, nodeLater
985  class(*), pointer :: obj => null()
986  !
987  badtime = -9.0d30
988  time0 = badtime
989  time1 = badtime
990  time = tsr%tsrTime
991  call this%get_surrounding_nodes(time, nodeearlier, nodelater)
992  !
993  if (associated(nodeearlier)) then
994  obj => nodeearlier%GetItem()
995  tsrearlier => castastimeseriesrecordtype(obj)
996  if (associated(tsrearlier)) then
997  time0 = tsrearlier%tsrTime
998  end if
999  end if
1000  !
1001  if (associated(nodelater)) then
1002  obj => nodelater%GetItem()
1003  tsrlater => castastimeseriesrecordtype(obj)
1004  if (associated(tsrlater)) then
1005  time1 = tsrlater%tsrTime
1006  end if
1007  end if
1008  !
1009  if (time0 > badtime) then
1010  ! Time0 is valid
1011  if (time1 > badtime) then
1012  ! Both time0 and time1 are valid
1013  if (time > time0 .and. time < time1) then
1014  ! Insert record between two list nodes
1015  obj => tsr
1016  call this%list%InsertBefore(obj, nodelater)
1017  else
1018  ! No need to insert a time series record, but if existing record
1019  ! for time of interest has NODATA as tsrValue, replace tsrValue
1020  if (time == time0 .and. tsrearlier%tsrValue == dnodata .and. &
1021  tsr%tsrValue /= dnodata) then
1022  tsrearlier%tsrValue = tsr%tsrValue
1023  elseif (time == time1 .and. tsrlater%tsrValue == dnodata .and. &
1024  tsr%tsrValue /= dnodata) then
1025  tsrlater%tsrValue = tsr%tsrValue
1026  end if
1027  end if
1028  else
1029  ! Time0 is valid and time1 is invalid. Just add tsr to the list.
1030  call this%AddTimeSeriesRecord(tsr)
1031  end if
1032  else
1033  ! Time0 is invalid, so time1 must be for first node in list
1034  if (time1 > badtime) then
1035  ! Time 1 is valid
1036  if (time < time1) then
1037  ! Insert tsr at beginning of list
1038  obj => tsr
1039  call this%list%InsertBefore(obj, nodelater)
1040  elseif (time == time1) then
1041  ! No need to insert a time series record, but if existing record
1042  ! for time of interest has NODATA as tsrValue, replace tsrValue
1043  if (tsrlater%tsrValue == dnodata .and. tsr%tsrValue /= dnodata) then
1044  tsrlater%tsrValue = tsr%tsrValue
1045  end if
1046  end if
1047  else
1048  ! Both time0 and time1 are invalid. Just add tsr to the list.
1049  call this%AddTimeSeriesRecord(tsr)
1050  end if
1051  end if
1052  !
1053  ! -- Return
1054  return
1055  end subroutine inserttsr
1056 
1057  !> @brief Find latest time
1058  !<
1059  function findlatesttime(this, readToEnd) result(endtime)
1060  ! -- dummy
1061  class(timeseriestype), intent(inout) :: this
1062  logical, intent(in), optional :: readtoend
1063  ! -- local
1064  integer :: nrecords
1065  type(timeseriesrecordtype), pointer :: tsr
1066  class(*), pointer :: obj => null()
1067  ! -- return
1068  double precision :: endtime
1069  !
1070  ! -- If the caller requested the very last time in the series (readToEnd is true), check that we have first read all records
1071  if (present(readtoend)) then
1072  if (readtoend) then
1073  do while (this%read_next_record())
1074  end do
1075  end if
1076  end if
1077  !
1078  nrecords = this%list%Count()
1079  obj => this%list%GetItem(nrecords)
1080  tsr => castastimeseriesrecordtype(obj)
1081  endtime = tsr%tsrTime
1082  !
1083  ! -- Return
1084  return
1085  end function findlatesttime
1086 
1087  !> @brief Clear the list of time series records
1088  !<
1089  subroutine clear(this, destroy)
1090  ! -- dummy
1091  class(timeseriestype), intent(inout) :: this
1092  logical, optional, intent(in) :: destroy
1093  !
1094  call this%list%Clear(destroy)
1095  !
1096  ! -- Return
1097  return
1098  end subroutine clear
1099 
1100 ! Type-bound procedures of TimeSeriesFileType
1101 
1102  !> @brief Count number of time series
1103  !<
1104  function count(this)
1105  ! -- return
1106  integer(I4B) :: count
1107  ! -- dummy
1108  class(timeseriesfiletype) :: this
1109  !
1110  if (associated(this%timeSeries)) then
1111  count = size(this%timeSeries)
1112  else
1113  count = 0
1114  end if
1115  !
1116  ! -- Return
1117  return
1118  end function count
1119 
1120  !> @brief Get time series
1121  !<
1122  function gettimeseries(this, indx) result(res)
1123  ! -- dummy
1124  class(timeseriesfiletype) :: this
1125  integer(I4B), intent(in) :: indx
1126  ! -- return
1127  type(timeseriestype), pointer :: res
1128  !
1129  res => null()
1130  if (indx > 0 .and. indx <= this%nTimeSeries) then
1131  res => this%timeSeries(indx)
1132  end if
1133  !
1134  ! -- Return
1135  return
1136  end function gettimeseries
1137 
1138  !> @brief Open time-series tsfile file and read options and first record,
1139  !! which may contain data to define multiple time series.
1140  !<
1141  subroutine initializetsfile(this, filename, iout, autoDeallocate)
1142  ! -- dummy
1143  class(timeseriesfiletype), target, intent(inout) :: this
1144  character(len=*), intent(in) :: filename
1145  integer(I4B), intent(in) :: iout
1146  logical, optional, intent(in) :: autoDeallocate
1147  ! -- local
1148  integer(I4B) :: iMethod, istatus, j, nwords
1149  integer(I4B) :: ierr, inunit
1150  logical :: autoDeallocateLocal = .true.
1151  logical :: continueread, found, endOfBlock
1152  logical :: methodWasSet
1153  real(DP) :: sfaclocal
1154  character(len=40) :: keyword, keyvalue
1155  character(len=:), allocatable :: line
1156  character(len=LENTIMESERIESNAME), allocatable, dimension(:) :: words
1157  !
1158  ! -- Initialize some variables
1159  if (present(autodeallocate)) autodeallocatelocal = autodeallocate
1160  imethod = undefined
1161  methodwasset = .false.
1162  !
1163  ! -- Assign members
1164  this%iout = iout
1165  this%datafile = filename
1166  !
1167  ! -- Open the time-series tsfile input file
1168  this%inunit = getunit()
1169  inunit = this%inunit
1170  call openfile(inunit, 0, filename, 'TS6')
1171  !
1172  ! -- Initialize block parser
1173  call this%parser%Initialize(this%inunit, this%iout)
1174  !
1175  ! -- Read the ATTRIBUTES block and count time series
1176  continueread = .false.
1177  ierr = 0
1178  !
1179  ! -- get BEGIN line of ATTRIBUTES block
1180  call this%parser%GetBlock('ATTRIBUTES', found, ierr, &
1181  supportopenclose=.true.)
1182  if (ierr /= 0) then
1183  ! end of file
1184  errmsg = 'End-of-file encountered while searching for'// &
1185  ' ATTRIBUTES in time-series '// &
1186  'input file "'//trim(this%datafile)//'"'
1187  call store_error(errmsg)
1188  call this%parser%StoreErrorUnit()
1189  elseif (.not. found) then
1190  errmsg = 'ATTRIBUTES block not found in time-series '// &
1191  'tsfile input file "'//trim(this%datafile)//'"'
1192  call store_error(errmsg)
1193  call this%parser%StoreErrorUnit()
1194  end if
1195  !
1196  ! -- parse ATTRIBUTES entries
1197  do
1198  ! -- read a line from input
1199  call this%parser%GetNextLine(endofblock)
1200  if (endofblock) exit
1201  !
1202  ! -- get the keyword
1203  call this%parser%GetStringCaps(keyword)
1204  !
1205  ! support either NAME or NAMES as equivalent keywords
1206  if (keyword == 'NAMES') keyword = 'NAME'
1207  !
1208  if (keyword /= 'NAME' .and. keyword /= 'METHODS' .and. &
1209  keyword /= 'SFACS') then
1210  ! -- get the word following the keyword (the key value)
1211  call this%parser%GetStringCaps(keyvalue)
1212  end if
1213  !
1214  select case (keyword)
1215  case ('NAME')
1216 ! line = line(istart:linelen)
1217  call this%parser%GetRemainingLine(line)
1218  call parseline(line, nwords, words, this%parser%iuactive)
1219  this%nTimeSeries = nwords
1220  ! -- Allocate the timeSeries array and initialize each
1221  ! time series.
1222  allocate (this%timeSeries(this%nTimeSeries))
1223  do j = 1, this%nTimeSeries
1224  call this%timeSeries(j)%initialize_time_series(this, words(j), &
1225  autodeallocatelocal)
1226  end do
1227  case ('METHOD')
1228  methodwasset = .true.
1229  if (this%nTimeSeries == 0) then
1230  errmsg = 'Error: NAME attribute not provided before METHOD in file: ' &
1231  //trim(filename)
1232  call store_error(errmsg)
1233  call this%parser%StoreErrorUnit()
1234  end if
1235  select case (keyvalue)
1236  case ('STEPWISE')
1237  imethod = stepwise
1238  case ('LINEAR')
1239  imethod = linear
1240  case ('LINEAREND')
1241  imethod = linearend
1242  case default
1243  errmsg = 'Unknown interpolation method: "'//trim(keyvalue)//'"'
1244  call store_error(errmsg)
1245  end select
1246  do j = 1, this%nTimeSeries
1247  this%timeSeries(j)%iMethod = imethod
1248  end do
1249  case ('METHODS')
1250  methodwasset = .true.
1251  if (this%nTimeSeries == 0) then
1252  errmsg = 'Error: NAME attribute not provided before METHODS in file: ' &
1253  //trim(filename)
1254  call store_error(errmsg)
1255  call this%parser%StoreErrorUnit()
1256  end if
1257  call this%parser%GetRemainingLine(line)
1258  call parseline(line, nwords, words, this%parser%iuactive)
1259  if (nwords < this%nTimeSeries) then
1260  errmsg = 'METHODS attribute does not list a method for'// &
1261  ' all time series.'
1262  call store_error(errmsg)
1263  call this%parser%StoreErrorUnit()
1264  end if
1265  do j = 1, this%nTimeSeries
1266  call upcase(words(j))
1267  select case (words(j))
1268  case ('STEPWISE')
1269  imethod = stepwise
1270  case ('LINEAR')
1271  imethod = linear
1272  case ('LINEAREND')
1273  imethod = linearend
1274  case default
1275  errmsg = 'Unknown interpolation method: "'//trim(words(j))//'"'
1276  call store_error(errmsg)
1277  end select
1278  this%timeSeries(j)%iMethod = imethod
1279  end do
1280  case ('SFAC')
1281  if (this%nTimeSeries == 0) then
1282  errmsg = 'NAME attribute not provided before SFAC in file: ' &
1283  //trim(filename)
1284  call store_error(errmsg)
1285  call this%parser%StoreErrorUnit()
1286  end if
1287  read (keyvalue, *, iostat=istatus) sfaclocal
1288  if (istatus /= 0) then
1289  errmsg = 'Error reading numeric value from: "'//trim(keyvalue)//'"'
1290  call store_error(errmsg)
1291  end if
1292  do j = 1, this%nTimeSeries
1293  this%timeSeries(j)%sfac = sfaclocal
1294  end do
1295  case ('SFACS')
1296  if (this%nTimeSeries == 0) then
1297  errmsg = 'NAME attribute not provided before SFACS in file: ' &
1298  //trim(filename)
1299  call store_error(errmsg)
1300  call this%parser%StoreErrorUnit()
1301  end if
1302  do j = 1, this%nTimeSeries
1303  sfaclocal = this%parser%GetDouble()
1304  this%timeSeries(j)%sfac = sfaclocal
1305  end do
1306  case ('AUTODEALLOCATE')
1307  do j = 1, this%nTimeSeries
1308  this%timeSeries(j)%autoDeallocate = (keyvalue == 'TRUE')
1309  end do
1310  case default
1311  errmsg = 'Unknown option found in ATTRIBUTES block: "'// &
1312  trim(keyword)//'"'
1313  call store_error(errmsg)
1314  call this%parser%StoreErrorUnit()
1315  end select
1316  end do
1317  !
1318  ! -- Get TIMESERIES block
1319  call this%parser%GetBlock('TIMESERIES', found, ierr, &
1320  supportopenclose=.true.)
1321  !
1322  ! -- Read the first line of time-series data
1323  if (.not. this%read_tsfile_line()) then
1324  errmsg = 'Error: No time-series data contained in file: '// &
1325  trim(this%datafile)
1326  call store_error(errmsg)
1327  end if
1328  !
1329  ! -- Ensure method was set
1330  if (.not. methodwasset) then
1331  errmsg = 'Interpolation method was not set. METHOD or METHODS &
1332  &must be specified in the ATTRIBUTES block for this time series file.'
1333  call store_error(errmsg)
1334  end if
1335  !
1336  ! -- Clean up and return
1337  if (allocated(words)) deallocate (words)
1338  !
1339  if (count_errors() > 0) then
1340  call this%parser%StoreErrorUnit()
1341  end if
1342  !
1343  ! -- Return
1344  return
1345  end subroutine initializetsfile
1346 
1347  !> @brief Read time series file line
1348  logical function read_tsfile_line(this)
1349  ! -- dummy
1350  class(timeseriesfiletype), intent(inout) :: this
1351  ! -- local
1352  real(dp) :: tsrtime, tsrvalue
1353  integer(I4B) :: i
1354  logical :: endofblock
1355  type(timeseriesrecordtype), pointer :: tsrecord => null()
1356  !
1357  read_tsfile_line = .false.
1358  !
1359  ! -- Get an arbitrary length, non-comment, non-blank line
1360  ! from the input file.
1361  call this%parser%GetNextLine(endofblock)
1362  !
1363  ! -- Check if we've reached the end of the TIMESERIES block
1364  if (endofblock) then
1365  return
1366  end if
1367  !
1368  ! -- Get the time
1369  tsrtime = this%parser%GetDouble()
1370  !
1371  ! -- Construct a new record and append a new node to each time series
1372  tsloop: do i = 1, this%nTimeSeries
1373  tsrvalue = this%parser%GetDouble()
1374  if (tsrvalue == dnodata) cycle tsloop
1375  ! -- multiply value by sfac
1376  tsrvalue = tsrvalue * this%timeSeries(i)%sfac
1377  call constructtimeseriesrecord(tsrecord, tsrtime, tsrvalue)
1378  call addtimeseriesrecordtolist(this%timeSeries(i)%list, tsrecord)
1379  end do tsloop
1380  read_tsfile_line = .true.
1381  !
1382  ! -- Return
1383  return
1384  end function read_tsfile_line
1385 
1386  !> @brief Deallocate memory
1387  !<
1388  subroutine tsf_da(this)
1389  ! -- dummy
1390  class(timeseriesfiletype), intent(inout) :: this
1391  ! -- local
1392  integer :: i, n
1393  type(timeseriestype), pointer :: ts => null()
1394  !
1395  n = this%Count()
1396  do i = 1, n
1397  ts => this%GetTimeSeries(i)
1398  if (associated(ts)) then
1399  call ts%da()
1400 ! deallocate(ts)
1401  end if
1402  end do
1403  !
1404  deallocate (this%timeSeries)
1405  deallocate (this%parser)
1406  !
1407  ! -- Return
1408  return
1409  end subroutine tsf_da
1410 
1411 end module timeseriesmodule
This module contains block parser methods.
Definition: BlockParser.f90:7
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:94
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
@ linearend
linear end interpolation
Definition: Constants.f90:146
@ undefined
undefined interpolation
Definition: Constants.f90:143
@ linear
linear interpolation
Definition: Constants.f90:145
@ stepwise
stepwise interpolation
Definition: Constants.f90:144
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
real(dp), parameter done
real constant 1
Definition: Constants.f90:75
integer(i4b) function, public getunit()
Get a free unit number.
subroutine, public parseline(line, nwords, words, inunit, filename)
Parse a line into words.
subroutine, public upcase(word)
Convert to upper case.
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
This module defines variable data types.
Definition: kind.f90:8
pure logical function, public is_close(a, b, rtol, atol, symmetric)
Check if a real value is approximately equal to another.
Definition: MathUtil.f90:46
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
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
type(timeseriestype) function, pointer gettimeseries(this, indx)
Get time series.
type(timeseriesrecordtype) function, pointer getcurrenttimeseriesrecord(this)
Get current ts record.
Definition: TimeSeries.f90:872
subroutine get_surrounding_nodes(this, time, nodeEarlier, nodeLater)
Get surrounding nodes.
Definition: TimeSeries.f90:389
subroutine inserttsr(this, tsr)
Insert a time series record.
Definition: TimeSeries.f90:978
integer(i4b) function count(this)
Count number of time series.
subroutine initialize_time_series(this, tsfile, name, autoDeallocate)
Initialize time series.
Definition: TimeSeries.f90:259
logical function read_next_record(this)
Read next record.
Definition: TimeSeries.f90:488
type(timeseriesfiletype) function, pointer, public gettimeseriesfilefromlist(list, idx)
Get time series from list.
Definition: TimeSeries.f90:164
subroutine, public constructtimeseriesfile(newTimeSeriesFile)
Construct time series file.
Definition: TimeSeries.f90:95
subroutine get_latest_preceding_node(this, time, tslNode)
Get latest preceding node.
Definition: TimeSeries.f90:771
subroutine reset(this)
Reset.
Definition: TimeSeries.f90:966
subroutine initializetsfile(this, filename, iout, autoDeallocate)
Open time-series tsfile file and read options and first record, which may contain data to define mult...
type(timeseriesrecordtype) function, pointer gettimeseriesrecord(this, time, epsi)
Get ts record.
Definition: TimeSeries.f90:935
subroutine, public addtimeseriesfiletolist(list, tsfile)
Add time series file to list.
Definition: TimeSeries.f90:148
real(dp) function getvalue(this, time0, time1, extendToEndOfSimulation)
Get time series value.
Definition: TimeSeries.f90:227
type(timeseriesfiletype) function, pointer castastimeseriesfiletype(obj)
Cast an unlimited polymorphic object as class(TimeSeriesFileType)
Definition: TimeSeries.f90:108
type(timeseriesfiletype) function, pointer, public castastimeseriesfileclass(obj)
Cast an unlimited polymorphic object as class(TimeSeriesFileType)
Definition: TimeSeries.f90:128
real(dp) function get_integrated_value(this, time0, time1, extendToEndOfSimulation)
Get integrated value.
Definition: TimeSeries.f90:597
real(dp) function get_value_at_time(this, time, extendToEndOfSimulation)
Get value for a time.
Definition: TimeSeries.f90:511
logical function, public sametimeseries(ts1, ts2)
Compare two time series; if they are identical, return true.
Definition: TimeSeries.f90:186
type(timeseriesrecordtype) function, pointer getnexttimeseriesrecord(this)
Get next ts record.
Definition: TimeSeries.f90:914
real(dp) function get_average_value(this, time0, time1, extendToEndOfSimulation)
Get average value.
Definition: TimeSeries.f90:736
subroutine clear(this, destroy)
Clear the list of time series records.
type(timeseriesrecordtype) function, pointer getprevioustimeseriesrecord(this)
Get previous ts record.
Definition: TimeSeries.f90:893
subroutine addtimeseriesrecord(this, tsr)
Add ts record.
Definition: TimeSeries.f90:856
subroutine ts_da(this)
Deallocate.
Definition: TimeSeries.f90:841
subroutine tsf_da(this)
Deallocate memory.
subroutine get_surrounding_records(this, time, tsrecEarlier, tsrecLater)
Get surrounding records.
Definition: TimeSeries.f90:294
logical function read_tsfile_line(this)
Read time series file line.
double precision function findlatesttime(this, readToEnd)
Find latest time.
subroutine, public addtimeseriesrecordtolist(list, tsrecord)
Add time series record to list.
subroutine, public constructtimeseriesrecord(newTsRecord, time, value)
Allocate and assign members of a new TimeSeriesRecordType object.
type(timeseriesrecordtype) function, pointer, public castastimeseriesrecordtype(obj)
Cast an unlimited polymorphic object as TimeSeriesRecordType.
A generic heterogeneous doubly-linked list.
Definition: List.f90:10