MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
TimeArraySeries.f90
Go to the documentation of this file.
2 
7  use mathutilmodule, only: is_close
9  use kindmodule, only: dp, i4b
10  use listmodule, only: listtype, listnodetype
11  use simvariablesmodule, only: errmsg
16  use, intrinsic :: iso_fortran_env, only: iostat_end
17 
18  implicit none
19  private
22 
24  ! -- Public members
25  character(len=LENTIMESERIESNAME), public :: name = ''
26  ! -- Private members
27  integer(I4B), private :: inunit = 0
28  integer(I4B), private :: iout = 0
29  integer(I4B), private :: imethod = undefined
30  real(dp), private :: sfac = done
31  character(len=LINELENGTH), private :: datafile = ''
32  logical, private :: autodeallocate = .true.
33  type(listtype), pointer, private :: list => null()
34  character(len=LENMODELNAME) :: modelname
35  type(blockparsertype), private :: parser
36 
37  contains
38 
39  ! -- Public procedures
40  procedure, public :: tas_init
41  procedure, public :: getaveragevalues
42  procedure, public :: getinunit
43  procedure, public :: da => tas_da
44  ! -- Private procedures
45  procedure, private :: get_integrated_values
46  procedure, private :: get_latest_preceding_node
47  procedure, private :: get_values_at_time
48  procedure, private :: get_surrounding_records
49  procedure, private :: read_next_array
50  procedure, private :: deallocatebackward
51  end type timearrayseriestype
52 
53 contains
54 
55  ! -- Constructor for TimeArraySeriesType
56 
57  !> @brief Allocate a new TimeArraySeriesType object.
58  !<
59  subroutine constructtimearrayseries(newTas, filename)
60  ! -- dummy
61  type(timearrayseriestype), pointer, intent(out) :: newtas
62  character(len=*), intent(in) :: filename
63  ! -- local
64  logical :: lex
65  ! -- formats
66 10 format('Error: Time-array-series file "', a, '" does not exist.')
67  !
68  ! -- Allocate a new object of type TimeArraySeriesType
69  allocate (newtas)
70  allocate (newtas%list)
71  !
72  ! -- Ensure that input file exists
73  inquire (file=filename, exist=lex)
74  if (.not. lex) then
75  write (errmsg, 10) trim(filename)
76  call store_error(errmsg, terminate=.true.)
77  end if
78  newtas%datafile = filename
79  !
80  ! -- Return
81  return
82  end subroutine constructtimearrayseries
83 
84  ! -- Public procedures
85 
86  !> @brief Initialize the time array series
87  !<
88  subroutine tas_init(this, fname, modelname, iout, tasname, autoDeallocate)
89  ! -- dummy
90  class(timearrayseriestype), intent(inout) :: this
91  character(len=*), intent(in) :: fname
92  character(len=*), intent(in) :: modelname
93  integer(I4B), intent(in) :: iout
94  character(len=*), intent(inout) :: tasname
95  logical, optional, intent(in) :: autoDeallocate
96  ! -- local
97  integer(I4B) :: istatus
98  integer(I4B) :: ierr
99  integer(I4B) :: inunit
100  character(len=40) :: keyword, keyvalue
101  logical :: found, continueread, endOfBlock
102  !
103  ! -- initialize some variables
104  if (present(autodeallocate)) this%autoDeallocate = autodeallocate
105  this%dataFile = fname
106  allocate (this%list)
107  !
108  ! -- assign members
109  this%modelname = modelname
110  this%iout = iout
111  !
112  ! -- open time-array series input file
113  inunit = getunit()
114  this%inunit = inunit
115  call openfile(inunit, 0, fname, 'TAS6')
116  !
117  ! -- initialize block parser
118  call this%parser%Initialize(this%inunit, this%iout)
119  !
120  ! -- read ATTRIBUTES block
121  continueread = .false.
122  ierr = 0
123  !
124  ! -- get BEGIN line of ATTRIBUTES block
125  call this%parser%GetBlock('ATTRIBUTES', found, ierr, &
126  supportopenclose=.true.)
127  if (.not. found) then
128  errmsg = 'Error: Attributes block not found in file: '// &
129  trim(fname)
130  call store_error(errmsg)
131  call this%parser%StoreErrorUnit()
132  end if
133  !
134  ! -- parse ATTRIBUTES entries
135  do
136  ! -- read line from input
137  call this%parser%GetNextLine(endofblock)
138  if (endofblock) exit
139  !
140  ! -- get the keyword
141  call this%parser%GetStringCaps(keyword)
142  !
143  ! -- get the word following the keyword (the key value)
144  call this%parser%GetStringCaps(keyvalue)
145  select case (keyword)
146  case ('NAME')
147  this%Name = keyvalue
148  tasname = keyvalue
149  case ('METHOD')
150  select case (keyvalue)
151  case ('STEPWISE')
152  this%iMethod = stepwise
153  case ('LINEAR')
154  this%iMethod = linear
155  case default
156  errmsg = 'Unknown interpolation method: "'//trim(keyvalue)//'"'
157  call store_error(errmsg)
158  call this%parser%StoreErrorUnit()
159  end select
160  case ('AUTODEALLOCATE')
161  this%autoDeallocate = (keyvalue == 'TRUE')
162  case ('SFAC')
163  read (keyvalue, *, iostat=istatus) this%sfac
164  if (istatus /= 0) then
165  errmsg = 'Error reading numeric SFAC value from "'//trim(keyvalue) &
166  //'"'
167  call store_error(errmsg)
168  call this%parser%StoreErrorUnit()
169  end if
170  case default
171  errmsg = 'Unknown option found in ATTRIBUTES block: "'// &
172  trim(keyword)//'"'
173  call store_error(errmsg)
174  call this%parser%StoreErrorUnit()
175  end select
176  end do
177  !
178  ! -- ensure that NAME and METHOD have been specified
179  if (this%Name == '') then
180  errmsg = 'Name not specified for time array series in file: '// &
181  trim(this%dataFile)
182  call store_error(errmsg)
183  call this%parser%StoreErrorUnit()
184  end if
185  if (this%iMethod == undefined) then
186  errmsg = 'Interpolation method not specified for time'// &
187  ' array series in file: '//trim(this%dataFile)
188  call store_error(errmsg)
189  call this%parser%StoreErrorUnit()
190  end if
191  !
192  ! -- handle any errors encountered so far
193  if (count_errors() > 0) then
194  errmsg = 'Error(s) encountered initializing time array series from file: ' &
195  //trim(this%dataFile)
196  call store_error(errmsg)
197  call this%parser%StoreErrorUnit()
198  end if
199  !
200  ! -- try to read first time array into linked list
201  if (.not. this%read_next_array()) then
202  errmsg = 'Error encountered reading time-array data from file: '// &
203  trim(this%dataFile)
204  call store_error(errmsg)
205  call this%parser%StoreErrorUnit()
206  end if
207  !
208  ! -- Return
209  return
210  end subroutine tas_init
211 
212  !> @brief Populate an array time-weighted average value for a specified time
213  !! span
214  !<
215  subroutine getaveragevalues(this, nvals, values, time0, time1)
216  ! -- dummy
217  class(timearrayseriestype), intent(inout) :: this
218  integer(I4B), intent(in) :: nvals
219  real(DP), dimension(nvals), intent(inout) :: values
220  real(DP), intent(in) :: time0
221  real(DP), intent(in) :: time1
222  ! -- local
223  integer(I4B) :: i
224  real(DP) :: timediff
225  !
226  timediff = time1 - time0
227  if (timediff > 0) then
228  call this%get_integrated_values(nvals, values, time0, time1)
229  do i = 1, nvals
230  values(i) = values(i) / timediff
231  end do
232  else
233  ! -- time0 and time1 are the same, so skip the integration step.
234  call this%get_values_at_time(nvals, values, time0)
235  end if
236  !
237  ! -- Return
238  return
239  end subroutine getaveragevalues
240 
241  !> @brief Return unit number
242  !<
243  function getinunit(this)
244  ! -- return
245  integer(I4B) :: getinunit
246  ! -- dummy
247  class(timearrayseriestype) :: this
248  !
249  getinunit = this%inunit
250  !
251  ! -- Return
252  return
253  end function getinunit
254 
255  ! -- Private procedures
256 
257  !> @brief Get surrounding records
258  !<
259  subroutine get_surrounding_records(this, time, taEarlier, taLater)
260  ! -- dummy
261  class(timearrayseriestype), intent(inout) :: this
262  real(DP), intent(in) :: time
263  type(timearraytype), pointer, intent(inout) :: taEarlier
264  type(timearraytype), pointer, intent(inout) :: taLater
265  ! -- local
266  real(DP) :: time0, time1
267  type(listnodetype), pointer :: currNode => null()
268  type(listnodetype), pointer :: node0 => null()
269  type(listnodetype), pointer :: node1 => null()
270  type(timearraytype), pointer :: ta => null(), ta0 => null(), ta1 => null()
271  class(*), pointer :: obj
272  !
273  taearlier => null()
274  talater => null()
275  !
276  if (associated(this%list%firstNode)) then
277  currnode => this%list%firstNode
278  end if
279  !
280  ! -- If the next node is earlier than time of interest, advance along
281  ! linked list until the next node is later than time of interest.
282  do
283  if (associated(currnode)) then
284  if (associated(currnode%nextNode)) then
285  obj => currnode%nextNode%GetItem()
286  ta => castastimearraytype(obj)
287  if (ta%taTime <= time) then
288  currnode => currnode%nextNode
289  else
290  exit
291  end if
292  else
293  ! -- read another array
294  if (.not. this%read_next_array()) exit
295  end if
296  else
297  exit
298  end if
299  end do
300  !
301  if (associated(currnode)) then
302  !
303  ! -- find earlier record
304  node0 => currnode
305  obj => node0%GetItem()
306  ta0 => castastimearraytype(obj)
307  time0 = ta0%taTime
308  do while (time0 > time)
309  if (associated(node0%prevNode)) then
310  node0 => node0%prevNode
311  obj => node0%GetItem()
312  ta0 => castastimearraytype(obj)
313  time0 = ta0%taTime
314  else
315  exit
316  end if
317  end do
318  !
319  ! -- find later record
320  node1 => currnode
321  obj => node1%GetItem()
322  ta1 => castastimearraytype(obj)
323  time1 = ta1%taTime
324  do while (time1 < time)
325  if (associated(node1%nextNode)) then
326  node1 => node1%nextNode
327  obj => node1%GetItem()
328  ta1 => castastimearraytype(obj)
329  time1 = ta1%taTime
330  else
331  ! -- get next array
332  if (.not. this%read_next_array()) then
333  ! -- end of file reached, so exit loop
334  exit
335  end if
336  end if
337  end do
338  !
339  end if
340  !
341  if (time0 <= time) taearlier => ta0
342  if (time1 >= time) talater => ta1
343  !
344  ! -- Return
345  return
346  end subroutine get_surrounding_records
347 
348  !> @brief Read next time array from input file and append to list
349  !<
350  logical function read_next_array(this)
351  ! -- modules
352  use constantsmodule, only: lenmempath
355  ! -- dummy
356  class(timearrayseriestype), intent(inout) :: this
357  ! -- local
358  integer(I4B) :: i, ierr, istart, istat, istop, lloc, nrow, ncol, &
359  nodesperlayer
360  logical :: lopen, isfound
361  type(timearraytype), pointer :: ta => null()
362  character(len=LENMEMPATH) :: mempath
363  integer(I4B), dimension(:), contiguous, pointer :: mshape
364  !
365  ! -- initialize
366  istart = 1
367  istat = 0
368  istop = 1
369  lloc = 1
370  nullify (mshape)
371  !
372  ! -- create mempath
373  mempath = create_mem_path(component=this%modelname, subcomponent='DIS')
374  !
375  ! -- set mshape pointer
376  call mem_setptr(mshape, 'MSHAPE', mempath)
377  !
378  ! Get dimensions for supported discretization type
379  if (size(mshape) == 2) then
380  nodesperlayer = mshape(2)
381  nrow = 1
382  ncol = mshape(2)
383  else if (size(mshape) == 3) then
384  nodesperlayer = mshape(2) * mshape(3)
385  nrow = mshape(2)
386  ncol = mshape(3)
387  else
388  errmsg = 'Time array series is not supported for selected &
389  &discretization type.'
390  call store_error(errmsg)
391  call this%parser%StoreErrorUnit()
392  end if
393  !
394  read_next_array = .false.
395  inquire (unit=this%inunit, opened=lopen)
396  if (lopen) then
397  call constructtimearray(ta, this%modelname)
398  ! -- read a time and an array from the input file
399  ! -- Get a TIME block and read the time
400  call this%parser%GetBlock('TIME', isfound, ierr, &
401  supportopenclose=.false.)
402  if (isfound) then
403  ta%taTime = this%parser%GetDouble()
404  ! -- Read the array
405  call readarray(this%parser%iuactive, ta%taArray, this%Name, &
406  size(mshape), ncol, nrow, 1, nodesperlayer, &
407  this%iout, 0, 0)
408  !
409  ! -- multiply values by sfac
410  do i = 1, nodesperlayer
411  ta%taArray(i) = ta%taArray(i) * this%sfac
412  end do
413  !
414  ! -- append the new time array to the list
415  call addtimearraytolist(this%list, ta)
416  read_next_array = .true.
417  !
418  ! -- make sure block is closed
419  call this%parser%terminateblock()
420  end if
421  end if
422  !
423  ! -- Return
424  return
425  end function read_next_array
426 
427  !> @brief Return an array of values for a specified time, same units as
428  !! time-series values
429  !<
430  subroutine get_values_at_time(this, nvals, values, time)
431  ! -- dummy
432  class(timearrayseriestype), intent(inout) :: this
433  integer(I4B), intent(in) :: nvals
434  real(DP), dimension(nvals), intent(inout) :: values
435  real(DP), intent(in) :: time ! time of interest
436  ! -- local
437  integer(I4B) :: i, ierr
438  real(DP) :: ratio, time0, time1, timediff, timediffi, val0, val1, &
439  valdiff
440  type(timearraytype), pointer :: taEarlier => null()
441  type(timearraytype), pointer :: taLater => null()
442  ! -- formats
443 10 format('Error getting array at time ', g10.3, &
444  ' for time-array series "', a, '"')
445  !
446  ierr = 0
447  call this%get_surrounding_records(time, taearlier, talater)
448  if (associated(taearlier)) then
449  if (associated(talater)) then
450  ! -- values are available for both earlier and later times
451  if (this%iMethod == stepwise) then
452  ! -- Just populate values from elements of earlier time array
453  values = taearlier%taArray
454  elseif (this%iMethod == linear) then
455  ! -- perform linear interpolation
456  time0 = taearlier%taTime
457  time1 = talater%tatime
458  timediff = time1 - time0
459  timediffi = time - time0
460  if (timediff > 0) then
461  ratio = timediffi / timediff
462  else
463  ! -- should not happen if TS does not contain duplicate times
464  ratio = 0.5d0
465  end if
466  ! -- Iterate through all elements and perform interpolation.
467  do i = 1, nvals
468  val0 = taearlier%taArray(i)
469  val1 = talater%taArray(i)
470  valdiff = val1 - val0
471  values(i) = val0 + (ratio * valdiff)
472  end do
473  else
474  ierr = 1
475  end if
476  else
477  if (is_close(taearlier%taTime, time)) then
478  values = taearlier%taArray
479  else
480  ! -- Only earlier time is available, and it is not time of interest;
481  ! however, if method is STEPWISE, use value for earlier time.
482  if (this%iMethod == stepwise) then
483  values = taearlier%taArray
484  else
485  ierr = 1
486  end if
487  end if
488  end if
489  else
490  if (associated(talater)) then
491  if (is_close(talater%taTime, time)) then
492  values = talater%taArray
493  else
494  ! -- only later time is available, and it is not time of interest
495  ierr = 1
496  end if
497  else
498  ! -- Neither earlier nor later time is available.
499  ! This should never happen!
500  ierr = 1
501  end if
502  end if
503  !
504  if (ierr > 0) then
505  write (errmsg, 10) time, trim(this%Name)
506  call store_error(errmsg)
507  call store_error_unit(this%inunit)
508  end if
509  !
510  ! -- Return
511  return
512  end subroutine get_values_at_time
513 
514  !> @brief Populates an array with integrated values for a specified time span
515  !!
516  !! Units: (ts-value-unit)*time
517  !<
518  subroutine get_integrated_values(this, nvals, values, time0, time1)
519  ! -- dummy
520  class(timearrayseriestype), intent(inout) :: this
521  integer(I4B), intent(in) :: nvals
522  real(DP), dimension(nvals), intent(inout) :: values
523  real(DP), intent(in) :: time0
524  real(DP), intent(in) :: time1
525  ! -- local
526  integer(I4B) :: i
527  real(DP) :: area, currTime, nextTime, ratio0, ratio1, t0, &
528  t01, t1, timediff, value, value0, value1, valuediff
529  logical :: ldone
530  type(listnodetype), pointer :: precNode => null()
531  type(listnodetype), pointer :: currNode => null(), nextnode => null()
532  type(timearraytype), pointer :: currRecord => null(), nextrecord => null()
533  class(*), pointer :: currObj => null(), nextobj => null()
534  ! -- formats
535 10 format('Error encountered while performing integration', &
536  ' for time-array series "', a, '" for time interval: ', &
537  g12.5, ' to ', g12.5)
538  !
539  values = dzero
540  value = dzero
541  ldone = .false.
542  t1 = -done
543  call this%get_latest_preceding_node(time0, precnode)
544  if (associated(precnode)) then
545  currnode => precnode
546  do while (.not. ldone)
547  currobj => currnode%GetItem()
548  currrecord => castastimearraytype(currobj)
549  currtime = currrecord%taTime
550  if (currtime < time1) then
551  if (.not. associated(currnode%nextNode)) then
552  ! -- try to read the next array
553  if (.not. this%read_next_array()) then
554  write (errmsg, 10) trim(this%Name), time0, time1
555  call store_error(errmsg)
556  call store_error_unit(this%inunit)
557  end if
558  end if
559  if (associated(currnode%nextNode)) then
560  nextnode => currnode%nextNode
561  nextobj => nextnode%GetItem()
562  nextrecord => castastimearraytype(nextobj)
563  nexttime = nextrecord%taTime
564  ! -- determine lower and upper limits of time span of interest
565  ! within current interval
566  if (currtime >= time0) then
567  t0 = currtime
568  else
569  t0 = time0
570  end if
571  if (nexttime <= time1) then
572  t1 = nexttime
573  else
574  t1 = time1
575  end if
576  ! -- For each element, find area of rectangle
577  ! or trapezoid delimited by t0 and t1.
578  t01 = t1 - t0
579  select case (this%iMethod)
580  case (stepwise)
581  do i = 1, nvals
582  ! -- compute area of a rectangle
583  value0 = currrecord%taArray(i)
584  area = value0 * t01
585  ! -- add area to integrated value
586  values(i) = values(i) + area
587  end do
588  case (linear)
589  do i = 1, nvals
590  ! -- compute area of a trapezoid
591  timediff = nexttime - currtime
592  ratio0 = (t0 - currtime) / timediff
593  ratio1 = (t1 - currtime) / timediff
594  valuediff = nextrecord%taArray(i) - currrecord%taArray(i)
595  value0 = currrecord%taArray(i) + ratio0 * valuediff
596  value1 = currrecord%taArray(i) + ratio1 * valuediff
597  area = 0.5d0 * t01 * (value0 + value1)
598  ! -- add area to integrated value
599  values(i) = values(i) + area
600  end do
601  end select
602  else
603  write (errmsg, 10) trim(this%Name), time0, time1
604  call store_error(errmsg)
605  call store_error('(Probable programming error)', terminate=.true.)
606  end if
607  else
608  ! Current node time = time1 so should be done
609  ldone = .true.
610  end if
611  !
612  ! -- Are we done yet?
613  if (t1 >= time1) then
614  ldone = .true.
615  else
616  if (.not. associated(currnode%nextNode)) then
617  ! -- try to read the next array
618  if (.not. this%read_next_array()) then
619  write (errmsg, 10) trim(this%Name), time0, time1
620  call store_error(errmsg)
621  call this%parser%StoreErrorUnit()
622  end if
623  end if
624  if (associated(currnode%nextNode)) then
625  currnode => currnode%nextNode
626  else
627  write (errmsg, 10) trim(this%Name), time0, time1
628  call store_error(errmsg)
629  call store_error('(Probable programming error)', terminate=.true.)
630  end if
631  end if
632  end do
633  end if
634  !
635  if (this%autoDeallocate) then
636  if (associated(precnode)) then
637  if (associated(precnode%prevNode)) then
638  call this%DeallocateBackward(precnode%prevNode)
639  end if
640  end if
641  end if
642  !
643  ! -- Return
644  return
645  end subroutine get_integrated_values
646 
647  !> @brief Deallocate fromNode and all previous nodes in list;
648  !! reassign firstNode
649  !<
650  subroutine deallocatebackward(this, fromNode)
651  ! -- dummy
652  class(timearrayseriestype), intent(inout) :: this
653  type(listnodetype), pointer, intent(inout) :: fromNode
654  !
655  ! -- local
656  type(listnodetype), pointer :: current => null()
657  type(listnodetype), pointer :: prev => null()
658  type(timearraytype), pointer :: ta => null()
659  class(*), pointer :: obj => null()
660  !
661  if (associated(fromnode)) then
662  ! -- reassign firstNode
663  if (associated(fromnode%nextNode)) then
664  this%list%firstNode => fromnode%nextNode
665  else
666  this%list%firstNode => null()
667  end if
668  ! -- deallocate fromNode and all previous nodes
669  current => fromnode
670  do while (associated(current))
671  prev => current%prevNode
672  obj => current%GetItem()
673  ta => castastimearraytype(obj)
674  ! -- Deallocate the contents of this time array,
675  ! then remove it from the list
676  call ta%da()
677  call this%list%RemoveNode(current, .true.)
678  current => prev
679  end do
680  fromnode => null()
681  end if
682  !
683  ! -- Return
684  return
685  end subroutine deallocatebackward
686 
687  !> @brief Return pointer to ListNodeType object for the node representing
688  !! the latest preceding time in the time series
689  !<
690  subroutine get_latest_preceding_node(this, time, tslNode)
691  ! -- dummy
692  class(timearrayseriestype), intent(inout) :: this
693  real(DP), intent(in) :: time
694  type(listnodetype), pointer, intent(inout) :: tslNode
695  ! -- local
696  real(DP) :: time0
697  type(listnodetype), pointer :: currNode => null()
698  type(listnodetype), pointer :: node0 => null()
699  type(timearraytype), pointer :: ta => null()
700  type(timearraytype), pointer :: ta0 => null()
701  class(*), pointer :: obj => null()
702  !
703  tslnode => null()
704  if (associated(this%list%firstNode)) then
705  currnode => this%list%firstNode
706  else
707  call store_error('probable programming error in &
708  &get_latest_preceding_node', &
709  terminate=.true.)
710  end if
711  !
712  continue
713  ! -- If the next node is earlier than time of interest, advance along
714  ! linked list until the next node is later than time of interest.
715  do
716  if (associated(currnode)) then
717  if (associated(currnode%nextNode)) then
718  obj => currnode%nextNode%GetItem()
719  ta => castastimearraytype(obj)
720  if (ta%taTime < time .or. is_close(ta%taTime, time)) then
721  currnode => currnode%nextNode
722  else
723  exit
724  end if
725  else
726  ! -- read another record
727  if (.not. this%read_next_array()) exit
728  end if
729  else
730  exit
731  end if
732  end do
733  !
734  if (associated(currnode)) then
735  !
736  ! -- find earlier record
737  node0 => currnode
738  obj => node0%GetItem()
739  ta0 => castastimearraytype(obj)
740  time0 = ta0%taTime
741  do while (time0 > time)
742  if (associated(node0%prevNode)) then
743  node0 => node0%prevNode
744  obj => node0%GetItem()
745  ta0 => castastimearraytype(obj)
746  time0 = ta0%taTime
747  else
748  exit
749  end if
750  end do
751  end if
752  !
753  if (time0 <= time) tslnode => node0
754  !
755  ! -- Return
756  return
757  end subroutine get_latest_preceding_node
758 
759  !> @brief Deallocate memory
760  !<
761  subroutine tas_da(this)
762  ! -- dummy
763  class(timearrayseriestype), intent(inout) :: this
764  ! -- local
765  integer :: i, n
766  type(timearraytype), pointer :: ta => null()
767  !
768  ! -- Deallocate contents of each time array in list
769  n = this%list%Count()
770  do i = 1, n
771  ta => gettimearrayfromlist(this%list, i)
772  call ta%da()
773  end do
774  !
775  ! -- Deallocate the list of time arrays
776  call this%list%Clear(.true.)
777  deallocate (this%list)
778  !
779  ! -- Return
780  return
781  end subroutine tas_da
782 
783  ! -- Procedures not type-bound
784 
785  !> @brief Cast an unlimited polymorphic object as class(TimeArraySeriesType)
786  !<
787  function castastimearrayseriestype(obj) result(res)
788  ! -- dummy
789  class(*), pointer, intent(inout) :: obj
790  ! -- return
791  type(timearrayseriestype), pointer :: res
792  !
793  res => null()
794  if (.not. associated(obj)) return
795  !
796  select type (obj)
797  type is (timearrayseriestype)
798  res => obj
799  end select
800  !
801  ! -- Return
802  return
803  end function castastimearrayseriestype
804 
805  !> @brief Get time array from list
806  !<
807  function gettimearrayseriesfromlist(list, indx) result(res)
808  ! -- dummy
809  type(listtype), intent(inout) :: list
810  integer, intent(in) :: indx
811  ! -- return
812  type(timearrayseriestype), pointer :: res
813  ! -- local
814  class(*), pointer :: obj
815  !
816  obj => list%GetItem(indx)
817  res => castastimearrayseriestype(obj)
818  !
819  ! -- Return
820  return
821  end function gettimearrayseriesfromlist
822 
823 end module timearrayseriesmodule
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
integer(i4b), parameter lenmodelname
maximum length of the model name
Definition: Constants.f90:21
integer(i4b), parameter lentimeseriesname
maximum length of a time series name
Definition: Constants.f90:41
@ 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
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:26
real(dp), parameter done
real constant 1
Definition: Constants.f90:75
integer(i4b) function, public getunit()
Get a free unit number.
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
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
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(timearraytype) function, pointer, public castastimearraytype(obj)
Cast an unlimited polymorphic object as TimeArrayType.
Definition: TimeArray.f90:78
type(timearraytype) function, pointer, public gettimearrayfromlist(list, indx)
Retrieve a time array from a list.
Definition: TimeArray.f90:114
subroutine, public constructtimearray(newTa, modelname)
Construct time array.
Definition: TimeArray.f90:36
subroutine, public addtimearraytolist(list, timearray)
Add a time array to a to list.
Definition: TimeArray.f90:98
subroutine, public constructtimearrayseries(newTas, filename)
Allocate a new TimeArraySeriesType object.
subroutine tas_init(this, fname, modelname, iout, tasname, autoDeallocate)
Initialize the time array series.
subroutine get_latest_preceding_node(this, time, tslNode)
Return pointer to ListNodeType object for the node representing the latest preceding time in the time...
subroutine getaveragevalues(this, nvals, values, time0, time1)
Populate an array time-weighted average value for a specified time span.
subroutine tas_da(this)
Deallocate memory.
subroutine get_values_at_time(this, nvals, values, time)
Return an array of values for a specified time, same units as time-series values.
subroutine deallocatebackward(this, fromNode)
Deallocate fromNode and all previous nodes in list; reassign firstNode.
subroutine get_surrounding_records(this, time, taEarlier, taLater)
Get surrounding records.
type(timearrayseriestype) function, pointer, public gettimearrayseriesfromlist(list, indx)
Get time array from list.
subroutine get_integrated_values(this, nvals, values, time0, time1)
Populates an array with integrated values for a specified time span.
type(timearrayseriestype) function, pointer, public castastimearrayseriestype(obj)
Cast an unlimited polymorphic object as class(TimeArraySeriesType)
integer(i4b) function getinunit(this)
Return unit number.
logical function read_next_array(this)
Read next time array from input file and append to list.
A generic heterogeneous doubly-linked list.
Definition: List.f90:10