MODFLOW 6  version 6.8.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 757 of file TimeArraySeries.f90.

758  ! -- dummy
759  class(*), pointer, intent(inout) :: obj
760  ! -- return
761  type(TimeArraySeriesType), pointer :: res
762  !
763  res => null()
764  if (.not. associated(obj)) return
765  !
766  select type (obj)
767  type is (timearrayseriestype)
768  res => obj
769  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 629 of file TimeArraySeries.f90.

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

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

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

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

775  ! -- dummy
776  type(ListType), intent(inout) :: list
777  integer, intent(in) :: indx
778  ! -- return
779  type(TimeArraySeriesType), pointer :: res
780  ! -- local
781  class(*), pointer :: obj
782  !
783  obj => list%GetItem(indx)
784  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
342  ! -- dummy
343  class(TimeArraySeriesType), intent(inout) :: this
344  ! -- local
345  integer(I4B) :: i, ierr, istart, istat, istop, lloc, nrow, ncol, &
346  nodesperlayer
347  logical :: lopen, isFound
348  type(TimeArrayType), pointer :: ta => null()
349  character(len=LENMEMPATH) :: mempath
350  integer(I4B), dimension(:), contiguous, pointer :: mshape
351  !
352  ! -- initialize
353  istart = 1
354  istat = 0
355  istop = 1
356  lloc = 1
357  nullify (mshape)
358  !
359  ! -- create mempath to model input context where MODEL_SHAPE is
360  ! -- stored during static DIS loading
361  mempath = create_mem_path(component=this%modelname, context=idm_context)
362  !
363  ! -- set mshape pointer
364  call mem_setptr(mshape, 'MODEL_SHAPE', mempath)
365  !
366  ! Get dimensions for supported discretization type
367  if (size(mshape) == 2) then
368  nodesperlayer = mshape(2)
369  nrow = 1
370  ncol = mshape(2)
371  else if (size(mshape) == 3) then
372  nodesperlayer = mshape(2) * mshape(3)
373  nrow = mshape(2)
374  ncol = mshape(3)
375  else
376  errmsg = 'Time array series is not supported for selected &
377  &discretization type.'
378  call store_error(errmsg)
379  call this%parser%StoreErrorUnit()
380  end if
381  !
382  read_next_array = .false.
383  inquire (unit=this%inunit, opened=lopen)
384  if (lopen) then
385  call constructtimearray(ta, this%modelname)
386  ! -- read a time and an array from the input file
387  ! -- Get a TIME block and read the time
388  call this%parser%GetBlock('TIME', isfound, ierr, &
389  supportopenclose=.false.)
390  if (isfound) then
391  ta%taTime = this%parser%GetDouble()
392  ! -- Read the array
393  call readarray(this%parser%iuactive, ta%taArray, this%Name, &
394  size(mshape), ncol, nrow, 1, nodesperlayer, &
395  this%iout, 0, 0)
396  !
397  ! -- multiply values by sfac
398  do i = 1, nodesperlayer
399  ta%taArray(i) = ta%taArray(i) * this%sfac
400  end do
401  !
402  ! -- append the new time array to the list
403  call addtimearraytolist(this%list, ta)
404  read_next_array = .true.
405  !
406  ! -- make sure block is closed
407  call this%parser%terminateblock()
408  end if
409  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
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=linelength) idm_context
Here is the call graph for this function:

◆ tas_da()

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

Definition at line 734 of file TimeArraySeries.f90.

735  ! -- dummy
736  class(TimeArraySeriesType), intent(inout) :: this
737  ! -- local
738  integer :: i, n
739  type(TimeArrayType), pointer :: ta => null()
740  !
741  ! -- Deallocate contents of each time array in list
742  n = this%list%Count()
743  do i = 1, n
744  ta => gettimearrayfromlist(this%list, i)
745  call ta%da()
746  end do
747  !
748  ! -- Deallocate the list of time arrays
749  call this%list%Clear(.true.)
750  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: