MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
timearrayseriesmodule Module Reference

Data Types

type  timearrayseriestype
 

Functions/Subroutines

subroutine, public constructtimearrayseries (newTas, filename)
 Allocate a new TimeArraySeriesType object. More...
 
subroutine tas_init (this, fname, modelname, iout, tasname, autoDeallocate)
 Initialize the time array series. More...
 
subroutine getaveragevalues (this, nvals, values, time0, time1)
 Populate an array time-weighted average value for a specified time span. More...
 
integer(i4b) function getinunit (this)
 Return unit number. More...
 
subroutine get_surrounding_records (this, time, taEarlier, taLater)
 Get surrounding records. More...
 
logical function read_next_array (this)
 Read next time array from input file and append to list. More...
 
subroutine get_values_at_time (this, nvals, values, time)
 Return an array of values for a specified time, same units as time-series values. More...
 
subroutine get_integrated_values (this, nvals, values, time0, time1)
 Populates an array with integrated values for a specified time span. More...
 
subroutine deallocatebackward (this, fromNode)
 Deallocate fromNode and all previous nodes in list; reassign firstNode. More...
 
subroutine get_latest_preceding_node (this, time, tslNode)
 Return pointer to ListNodeType object for the node representing the latest preceding time in the time series. More...
 
subroutine tas_da (this)
 Deallocate memory. More...
 
type(timearrayseriestype) function, pointer, public castastimearrayseriestype (obj)
 Cast an unlimited polymorphic object as class(TimeArraySeriesType) More...
 
type(timearrayseriestype) function, pointer, public gettimearrayseriesfromlist (list, indx)
 Get time array from list. More...
 

Function/Subroutine Documentation

◆ castastimearrayseriestype()

type(timearrayseriestype) function, pointer, public timearrayseriesmodule::castastimearrayseriestype ( class(*), intent(inout), pointer  obj)

Definition at line 755 of file TimeArraySeries.f90.

756  ! -- dummy
757  class(*), pointer, intent(inout) :: obj
758  ! -- return
759  type(TimeArraySeriesType), pointer :: res
760  !
761  res => null()
762  if (.not. associated(obj)) return
763  !
764  select type (obj)
765  type is (timearrayseriestype)
766  res => obj
767  end select
Here is the caller graph for this function:

◆ constructtimearrayseries()

subroutine, public timearrayseriesmodule::constructtimearrayseries ( type(timearrayseriestype), intent(out), pointer  newTas,
character(len=*), intent(in)  filename 
)

Definition at line 60 of file TimeArraySeries.f90.

61  ! -- dummy
62  type(TimeArraySeriesType), pointer, intent(out) :: newTas
63  character(len=*), intent(in) :: filename
64  ! -- local
65  logical :: lex
66  ! -- formats
67 10 format('Error: Time-array-series file "', a, '" does not exist.')
68  !
69  ! -- Allocate a new object of type TimeArraySeriesType
70  allocate (newtas)
71  allocate (newtas%list)
72  !
73  ! -- Ensure that input file exists
74  inquire (file=filename, exist=lex)
75  if (.not. lex) then
76  write (errmsg, 10) trim(filename)
77  call store_error(errmsg, terminate=.true.)
78  end if
79  newtas%datafile = filename
Here is the call graph for this function:

◆ deallocatebackward()

subroutine timearrayseriesmodule::deallocatebackward ( class(timearrayseriestype), intent(inout)  this,
type(listnodetype), intent(inout), pointer  fromNode 
)
private

Definition at line 627 of file TimeArraySeries.f90.

628  ! -- dummy
629  class(TimeArraySeriesType), intent(inout) :: this
630  type(ListNodeType), pointer, intent(inout) :: fromNode
631  !
632  ! -- local
633  type(ListNodeType), pointer :: current => null()
634  type(ListNodeType), pointer :: prev => null()
635  type(TimeArrayType), pointer :: ta => null()
636  class(*), pointer :: obj => null()
637  !
638  if (associated(fromnode)) then
639  ! -- reassign firstNode
640  if (associated(fromnode%nextNode)) then
641  this%list%firstNode => fromnode%nextNode
642  else
643  this%list%firstNode => null()
644  end if
645  ! -- deallocate fromNode and all previous nodes
646  current => fromnode
647  do while (associated(current))
648  prev => current%prevNode
649  obj => current%GetItem()
650  ta => castastimearraytype(obj)
651  ! -- Deallocate the contents of this time array,
652  ! then remove it from the list
653  call ta%da()
654  call this%list%RemoveNode(current, .true.)
655  current => prev
656  end do
657  fromnode => null()
658  end if
Here is the call graph for this function:

◆ get_integrated_values()

subroutine timearrayseriesmodule::get_integrated_values ( class(timearrayseriestype), intent(inout)  this,
integer(i4b), intent(in)  nvals,
real(dp), dimension(nvals), intent(inout)  values,
real(dp), intent(in)  time0,
real(dp), intent(in)  time1 
)
private

Units: (ts-value-unit)*time

Definition at line 498 of file TimeArraySeries.f90.

499  ! -- dummy
500  class(TimeArraySeriesType), intent(inout) :: this
501  integer(I4B), intent(in) :: nvals
502  real(DP), dimension(nvals), intent(inout) :: values
503  real(DP), intent(in) :: time0
504  real(DP), intent(in) :: time1
505  ! -- local
506  integer(I4B) :: i
507  real(DP) :: area, currTime, nextTime, ratio0, ratio1, t0, &
508  t01, t1, timediff, value, value0, value1, valuediff
509  logical :: ldone
510  type(ListNodeType), pointer :: precNode => null()
511  type(ListNodeType), pointer :: currNode => null(), nextnode => null()
512  type(TimeArrayType), pointer :: currRecord => null(), nextrecord => null()
513  class(*), pointer :: currObj => null(), nextobj => null()
514  ! -- formats
515 10 format('Error encountered while performing integration', &
516  ' for time-array series "', a, '" for time interval: ', &
517  g12.5, ' to ', g12.5)
518  !
519  values = dzero
520  value = dzero
521  ldone = .false.
522  t1 = -done
523  call this%get_latest_preceding_node(time0, precnode)
524  if (associated(precnode)) then
525  currnode => precnode
526  do while (.not. ldone)
527  currobj => currnode%GetItem()
528  currrecord => castastimearraytype(currobj)
529  currtime = currrecord%taTime
530  if (currtime < time1) then
531  if (.not. associated(currnode%nextNode)) then
532  ! -- try to read the next array
533  if (.not. this%read_next_array()) then
534  write (errmsg, 10) trim(this%Name), time0, time1
535  call store_error(errmsg)
536  call store_error_unit(this%inunit)
537  end if
538  end if
539  if (associated(currnode%nextNode)) then
540  nextnode => currnode%nextNode
541  nextobj => nextnode%GetItem()
542  nextrecord => castastimearraytype(nextobj)
543  nexttime = nextrecord%taTime
544  ! -- determine lower and upper limits of time span of interest
545  ! within current interval
546  if (currtime >= time0) then
547  t0 = currtime
548  else
549  t0 = time0
550  end if
551  if (nexttime <= time1) then
552  t1 = nexttime
553  else
554  t1 = time1
555  end if
556  ! -- For each element, find area of rectangle
557  ! or trapezoid delimited by t0 and t1.
558  t01 = t1 - t0
559  select case (this%iMethod)
560  case (stepwise)
561  do i = 1, nvals
562  ! -- compute area of a rectangle
563  value0 = currrecord%taArray(i)
564  area = value0 * t01
565  ! -- add area to integrated value
566  values(i) = values(i) + area
567  end do
568  case (linear)
569  do i = 1, nvals
570  ! -- compute area of a trapezoid
571  timediff = nexttime - currtime
572  ratio0 = (t0 - currtime) / timediff
573  ratio1 = (t1 - currtime) / timediff
574  valuediff = nextrecord%taArray(i) - currrecord%taArray(i)
575  value0 = currrecord%taArray(i) + ratio0 * valuediff
576  value1 = currrecord%taArray(i) + ratio1 * valuediff
577  area = 0.5d0 * t01 * (value0 + value1)
578  ! -- add area to integrated value
579  values(i) = values(i) + area
580  end do
581  end select
582  else
583  write (errmsg, 10) trim(this%Name), time0, time1
584  call store_error(errmsg)
585  call store_error('(Probable programming error)', terminate=.true.)
586  end if
587  else
588  ! Current node time = time1 so should be done
589  ldone = .true.
590  end if
591  !
592  ! -- Are we done yet?
593  if (t1 >= time1) then
594  ldone = .true.
595  else
596  if (.not. associated(currnode%nextNode)) then
597  ! -- try to read the next array
598  if (.not. this%read_next_array()) then
599  write (errmsg, 10) trim(this%Name), time0, time1
600  call store_error(errmsg)
601  call this%parser%StoreErrorUnit()
602  end if
603  end if
604  if (associated(currnode%nextNode)) then
605  currnode => currnode%nextNode
606  else
607  write (errmsg, 10) trim(this%Name), time0, time1
608  call store_error(errmsg)
609  call store_error('(Probable programming error)', terminate=.true.)
610  end if
611  end if
612  end do
613  end if
614  !
615  if (this%autoDeallocate) then
616  if (associated(precnode)) then
617  if (associated(precnode%prevNode)) then
618  call this%DeallocateBackward(precnode%prevNode)
619  end if
620  end if
621  end if
Here is the call graph for this function:

◆ get_latest_preceding_node()

subroutine timearrayseriesmodule::get_latest_preceding_node ( class(timearrayseriestype), intent(inout)  this,
real(dp), intent(in)  time,
type(listnodetype), intent(inout), pointer  tslNode 
)
private

Definition at line 664 of file TimeArraySeries.f90.

665  ! -- dummy
666  class(TimeArraySeriesType), intent(inout) :: this
667  real(DP), intent(in) :: time
668  type(ListNodeType), pointer, intent(inout) :: tslNode
669  ! -- local
670  real(DP) :: time0
671  type(ListNodeType), pointer :: currNode => null()
672  type(ListNodeType), pointer :: node0 => null()
673  type(TimeArrayType), pointer :: ta => null()
674  type(TimeArrayType), pointer :: ta0 => null()
675  class(*), pointer :: obj => null()
676  !
677  tslnode => null()
678  if (associated(this%list%firstNode)) then
679  currnode => this%list%firstNode
680  else
681  call store_error('probable programming error in &
682  &get_latest_preceding_node', &
683  terminate=.true.)
684  end if
685  !
686  continue
687  ! -- If the next node is earlier than time of interest, advance along
688  ! linked list until the next node is later than time of interest.
689  do
690  if (associated(currnode)) then
691  if (associated(currnode%nextNode)) then
692  obj => currnode%nextNode%GetItem()
693  ta => castastimearraytype(obj)
694  if (ta%taTime < time .or. is_close(ta%taTime, time)) then
695  currnode => currnode%nextNode
696  else
697  exit
698  end if
699  else
700  ! -- read another record
701  if (.not. this%read_next_array()) exit
702  end if
703  else
704  exit
705  end if
706  end do
707  !
708  if (associated(currnode)) then
709  !
710  ! -- find earlier record
711  node0 => currnode
712  obj => node0%GetItem()
713  ta0 => castastimearraytype(obj)
714  time0 = ta0%taTime
715  do while (time0 > time)
716  if (associated(node0%prevNode)) then
717  node0 => node0%prevNode
718  obj => node0%GetItem()
719  ta0 => castastimearraytype(obj)
720  time0 = ta0%taTime
721  else
722  exit
723  end if
724  end do
725  end if
726  !
727  if (time0 <= time) tslnode => node0
Here is the call graph for this function:

◆ get_surrounding_records()

subroutine timearrayseriesmodule::get_surrounding_records ( class(timearrayseriestype), intent(inout)  this,
real(dp), intent(in)  time,
type(timearraytype), intent(inout), pointer  taEarlier,
type(timearraytype), intent(inout), pointer  taLater 
)
private

Definition at line 248 of file TimeArraySeries.f90.

249  ! -- dummy
250  class(TimeArraySeriesType), intent(inout) :: this
251  real(DP), intent(in) :: time
252  type(TimeArrayType), pointer, intent(inout) :: taEarlier
253  type(TimeArrayType), pointer, intent(inout) :: taLater
254  ! -- local
255  real(DP) :: time0, time1
256  type(ListNodeType), pointer :: currNode => null()
257  type(ListNodeType), pointer :: node0 => null()
258  type(ListNodeType), pointer :: node1 => null()
259  type(TimeArrayType), pointer :: ta => null(), ta0 => null(), ta1 => null()
260  class(*), pointer :: obj
261  !
262  taearlier => null()
263  talater => null()
264  !
265  if (associated(this%list%firstNode)) then
266  currnode => this%list%firstNode
267  end if
268  !
269  ! -- If the next node is earlier than time of interest, advance along
270  ! linked list until the next node is later than time of interest.
271  do
272  if (associated(currnode)) then
273  if (associated(currnode%nextNode)) then
274  obj => currnode%nextNode%GetItem()
275  ta => castastimearraytype(obj)
276  if (ta%taTime <= time) then
277  currnode => currnode%nextNode
278  else
279  exit
280  end if
281  else
282  ! -- read another array
283  if (.not. this%read_next_array()) exit
284  end if
285  else
286  exit
287  end if
288  end do
289  !
290  if (associated(currnode)) then
291  !
292  ! -- find earlier record
293  node0 => currnode
294  obj => node0%GetItem()
295  ta0 => castastimearraytype(obj)
296  time0 = ta0%taTime
297  do while (time0 > time)
298  if (associated(node0%prevNode)) then
299  node0 => node0%prevNode
300  obj => node0%GetItem()
301  ta0 => castastimearraytype(obj)
302  time0 = ta0%taTime
303  else
304  exit
305  end if
306  end do
307  !
308  ! -- find later record
309  node1 => currnode
310  obj => node1%GetItem()
311  ta1 => castastimearraytype(obj)
312  time1 = ta1%taTime
313  do while (time1 < time)
314  if (associated(node1%nextNode)) then
315  node1 => node1%nextNode
316  obj => node1%GetItem()
317  ta1 => castastimearraytype(obj)
318  time1 = ta1%taTime
319  else
320  ! -- get next array
321  if (.not. this%read_next_array()) then
322  ! -- end of file reached, so exit loop
323  exit
324  end if
325  end if
326  end do
327  !
328  end if
329  !
330  if (time0 <= time) taearlier => ta0
331  if (time1 >= time) talater => ta1
Here is the call graph for this function:

◆ get_values_at_time()

subroutine timearrayseriesmodule::get_values_at_time ( class(timearrayseriestype), intent(inout)  this,
integer(i4b), intent(in)  nvals,
real(dp), dimension(nvals), intent(inout)  values,
real(dp), intent(in)  time 
)

Definition at line 413 of file TimeArraySeries.f90.

414  ! -- dummy
415  class(TimeArraySeriesType), intent(inout) :: this
416  integer(I4B), intent(in) :: nvals
417  real(DP), dimension(nvals), intent(inout) :: values
418  real(DP), intent(in) :: time ! time of interest
419  ! -- local
420  integer(I4B) :: i, ierr
421  real(DP) :: ratio, time0, time1, timediff, timediffi, val0, val1, &
422  valdiff
423  type(TimeArrayType), pointer :: taEarlier => null()
424  type(TimeArrayType), pointer :: taLater => null()
425  ! -- formats
426 10 format('Error getting array at time ', g10.3, &
427  ' for time-array series "', a, '"')
428  !
429  ierr = 0
430  call this%get_surrounding_records(time, taearlier, talater)
431  if (associated(taearlier)) then
432  if (associated(talater)) then
433  ! -- values are available for both earlier and later times
434  if (this%iMethod == stepwise) then
435  ! -- Just populate values from elements of earlier time array
436  values = taearlier%taArray
437  elseif (this%iMethod == linear) then
438  ! -- perform linear interpolation
439  time0 = taearlier%taTime
440  time1 = talater%tatime
441  timediff = time1 - time0
442  timediffi = time - time0
443  if (timediff > 0) then
444  ratio = timediffi / timediff
445  else
446  ! -- should not happen if TS does not contain duplicate times
447  ratio = 0.5d0
448  end if
449  ! -- Iterate through all elements and perform interpolation.
450  do i = 1, nvals
451  val0 = taearlier%taArray(i)
452  val1 = talater%taArray(i)
453  valdiff = val1 - val0
454  values(i) = val0 + (ratio * valdiff)
455  end do
456  else
457  ierr = 1
458  end if
459  else
460  if (is_close(taearlier%taTime, time)) then
461  values = taearlier%taArray
462  else
463  ! -- Only earlier time is available, and it is not time of interest;
464  ! however, if method is STEPWISE, use value for earlier time.
465  if (this%iMethod == stepwise) then
466  values = taearlier%taArray
467  else
468  ierr = 1
469  end if
470  end if
471  end if
472  else
473  if (associated(talater)) then
474  if (is_close(talater%taTime, time)) then
475  values = talater%taArray
476  else
477  ! -- only later time is available, and it is not time of interest
478  ierr = 1
479  end if
480  else
481  ! -- Neither earlier nor later time is available.
482  ! This should never happen!
483  ierr = 1
484  end if
485  end if
486  !
487  if (ierr > 0) then
488  write (errmsg, 10) time, trim(this%Name)
489  call store_error(errmsg)
490  call store_error_unit(this%inunit)
491  end if
Here is the call graph for this function:

◆ getaveragevalues()

subroutine timearrayseriesmodule::getaveragevalues ( class(timearrayseriestype), intent(inout)  this,
integer(i4b), intent(in)  nvals,
real(dp), dimension(nvals), intent(inout)  values,
real(dp), intent(in)  time0,
real(dp), intent(in)  time1 
)
private

Definition at line 210 of file TimeArraySeries.f90.

211  ! -- dummy
212  class(TimeArraySeriesType), intent(inout) :: this
213  integer(I4B), intent(in) :: nvals
214  real(DP), dimension(nvals), intent(inout) :: values
215  real(DP), intent(in) :: time0
216  real(DP), intent(in) :: time1
217  ! -- local
218  integer(I4B) :: i
219  real(DP) :: timediff
220  !
221  timediff = time1 - time0
222  if (timediff > 0) then
223  call this%get_integrated_values(nvals, values, time0, time1)
224  do i = 1, nvals
225  values(i) = values(i) / timediff
226  end do
227  else
228  ! -- time0 and time1 are the same, so skip the integration step.
229  call this%get_values_at_time(nvals, values, time0)
230  end if

◆ getinunit()

integer(i4b) function timearrayseriesmodule::getinunit ( class(timearrayseriestype this)
private

Definition at line 235 of file TimeArraySeries.f90.

236  ! -- return
237  integer(I4B) :: GetInunit
238  ! -- dummy
239  class(TimeArraySeriesType) :: this
240  !
241  getinunit = this%inunit

◆ gettimearrayseriesfromlist()

type(timearrayseriestype) function, pointer, public timearrayseriesmodule::gettimearrayseriesfromlist ( type(listtype), intent(inout)  list,
integer, intent(in)  indx 
)

Definition at line 772 of file TimeArraySeries.f90.

773  ! -- dummy
774  type(ListType), intent(inout) :: list
775  integer, intent(in) :: indx
776  ! -- return
777  type(TimeArraySeriesType), pointer :: res
778  ! -- local
779  class(*), pointer :: obj
780  !
781  obj => list%GetItem(indx)
782  res => castastimearrayseriestype(obj)
Here is the call graph for this function:

◆ read_next_array()

logical function timearrayseriesmodule::read_next_array ( class(timearrayseriestype), intent(inout)  this)
private

Definition at line 336 of file TimeArraySeries.f90.

337  ! -- modules
338  use constantsmodule, only: lenmempath
341  ! -- dummy
342  class(TimeArraySeriesType), intent(inout) :: this
343  ! -- local
344  integer(I4B) :: i, ierr, istart, istat, istop, lloc, nrow, ncol, &
345  nodesperlayer
346  logical :: lopen, isFound
347  type(TimeArrayType), pointer :: ta => null()
348  character(len=LENMEMPATH) :: mempath
349  integer(I4B), dimension(:), contiguous, pointer :: mshape
350  !
351  ! -- initialize
352  istart = 1
353  istat = 0
354  istop = 1
355  lloc = 1
356  nullify (mshape)
357  !
358  ! -- create mempath
359  mempath = create_mem_path(component=this%modelname, subcomponent='DIS')
360  !
361  ! -- set mshape pointer
362  call mem_setptr(mshape, 'MSHAPE', mempath)
363  !
364  ! Get dimensions for supported discretization type
365  if (size(mshape) == 2) then
366  nodesperlayer = mshape(2)
367  nrow = 1
368  ncol = mshape(2)
369  else if (size(mshape) == 3) then
370  nodesperlayer = mshape(2) * mshape(3)
371  nrow = mshape(2)
372  ncol = mshape(3)
373  else
374  errmsg = 'Time array series is not supported for selected &
375  &discretization type.'
376  call store_error(errmsg)
377  call this%parser%StoreErrorUnit()
378  end if
379  !
380  read_next_array = .false.
381  inquire (unit=this%inunit, opened=lopen)
382  if (lopen) then
383  call constructtimearray(ta, this%modelname)
384  ! -- read a time and an array from the input file
385  ! -- Get a TIME block and read the time
386  call this%parser%GetBlock('TIME', isfound, ierr, &
387  supportopenclose=.false.)
388  if (isfound) then
389  ta%taTime = this%parser%GetDouble()
390  ! -- Read the array
391  call readarray(this%parser%iuactive, ta%taArray, this%Name, &
392  size(mshape), ncol, nrow, 1, nodesperlayer, &
393  this%iout, 0, 0)
394  !
395  ! -- multiply values by sfac
396  do i = 1, nodesperlayer
397  ta%taArray(i) = ta%taArray(i) * this%sfac
398  end do
399  !
400  ! -- append the new time array to the list
401  call addtimearraytolist(this%list, ta)
402  read_next_array = .true.
403  !
404  ! -- make sure block is closed
405  call this%parser%terminateblock()
406  end if
407  end if
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
Here is the call graph for this function:

◆ tas_da()

subroutine timearrayseriesmodule::tas_da ( class(timearrayseriestype), intent(inout)  this)
private

Definition at line 732 of file TimeArraySeries.f90.

733  ! -- dummy
734  class(TimeArraySeriesType), intent(inout) :: this
735  ! -- local
736  integer :: i, n
737  type(TimeArrayType), pointer :: ta => null()
738  !
739  ! -- Deallocate contents of each time array in list
740  n = this%list%Count()
741  do i = 1, n
742  ta => gettimearrayfromlist(this%list, i)
743  call ta%da()
744  end do
745  !
746  ! -- Deallocate the list of time arrays
747  call this%list%Clear(.true.)
748  deallocate (this%list)
Here is the call graph for this function:

◆ tas_init()

subroutine timearrayseriesmodule::tas_init ( class(timearrayseriestype), intent(inout)  this,
character(len=*), intent(in)  fname,
character(len=*), intent(in)  modelname,
integer(i4b), intent(in)  iout,
character(len=*), intent(inout)  tasname,
logical, intent(in), optional  autoDeallocate 
)
private

Definition at line 86 of file TimeArraySeries.f90.

87  ! -- dummy
88  class(TimeArraySeriesType), intent(inout) :: this
89  character(len=*), intent(in) :: fname
90  character(len=*), intent(in) :: modelname
91  integer(I4B), intent(in) :: iout
92  character(len=*), intent(inout) :: tasname
93  logical, optional, intent(in) :: autoDeallocate
94  ! -- local
95  integer(I4B) :: istatus
96  integer(I4B) :: ierr
97  integer(I4B) :: inunit
98  character(len=40) :: keyword, keyvalue
99  logical :: found, continueread, endOfBlock
100  !
101  ! -- initialize some variables
102  if (present(autodeallocate)) this%autoDeallocate = autodeallocate
103  this%dataFile = fname
104  allocate (this%list)
105  !
106  ! -- assign members
107  this%modelname = modelname
108  this%iout = iout
109  !
110  ! -- open time-array series input file
111  inunit = getunit()
112  this%inunit = inunit
113  call openfile(inunit, 0, fname, 'TAS6')
114  !
115  ! -- initialize block parser
116  call this%parser%Initialize(this%inunit, this%iout)
117  !
118  ! -- read ATTRIBUTES block
119  continueread = .false.
120  ierr = 0
121  !
122  ! -- get BEGIN line of ATTRIBUTES block
123  call this%parser%GetBlock('ATTRIBUTES', found, ierr, &
124  supportopenclose=.true.)
125  if (.not. found) then
126  errmsg = 'Error: Attributes block not found in file: '// &
127  trim(fname)
128  call store_error(errmsg)
129  call this%parser%StoreErrorUnit()
130  end if
131  !
132  ! -- parse ATTRIBUTES entries
133  do
134  ! -- read line from input
135  call this%parser%GetNextLine(endofblock)
136  if (endofblock) exit
137  !
138  ! -- get the keyword
139  call this%parser%GetStringCaps(keyword)
140  !
141  ! -- get the word following the keyword (the key value)
142  call this%parser%GetStringCaps(keyvalue)
143  select case (keyword)
144  case ('NAME')
145  this%Name = keyvalue
146  tasname = keyvalue
147  case ('METHOD')
148  select case (keyvalue)
149  case ('STEPWISE')
150  this%iMethod = stepwise
151  case ('LINEAR')
152  this%iMethod = linear
153  case default
154  errmsg = 'Unknown interpolation method: "'//trim(keyvalue)//'"'
155  call store_error(errmsg)
156  call this%parser%StoreErrorUnit()
157  end select
158  case ('AUTODEALLOCATE')
159  this%autoDeallocate = (keyvalue == 'TRUE')
160  case ('SFAC')
161  read (keyvalue, *, iostat=istatus) this%sfac
162  if (istatus /= 0) then
163  errmsg = 'Error reading numeric SFAC value from "'//trim(keyvalue) &
164  //'"'
165  call store_error(errmsg)
166  call this%parser%StoreErrorUnit()
167  end if
168  case default
169  errmsg = 'Unknown option found in ATTRIBUTES block: "'// &
170  trim(keyword)//'"'
171  call store_error(errmsg)
172  call this%parser%StoreErrorUnit()
173  end select
174  end do
175  !
176  ! -- ensure that NAME and METHOD have been specified
177  if (this%Name == '') then
178  errmsg = 'Name not specified for time array series in file: '// &
179  trim(this%dataFile)
180  call store_error(errmsg)
181  call this%parser%StoreErrorUnit()
182  end if
183  if (this%iMethod == undefined) then
184  errmsg = 'Interpolation method not specified for time'// &
185  ' array series in file: '//trim(this%dataFile)
186  call store_error(errmsg)
187  call this%parser%StoreErrorUnit()
188  end if
189  !
190  ! -- handle any errors encountered so far
191  if (count_errors() > 0) then
192  errmsg = 'Error(s) encountered initializing time array series from file: ' &
193  //trim(this%dataFile)
194  call store_error(errmsg)
195  call this%parser%StoreErrorUnit()
196  end if
197  !
198  ! -- try to read first time array into linked list
199  if (.not. this%read_next_array()) then
200  errmsg = 'Error encountered reading time-array data from file: '// &
201  trim(this%dataFile)
202  call store_error(errmsg)
203  call this%parser%StoreErrorUnit()
204  end if
Here is the call graph for this function: