1 !***********************************************************************
2 !* GNU Lesser General Public License
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 !> @defgroup time_manager_mod time_manager_mod
20 !> @ingroup time_manager
21 !> @link http://www.gfdl.noaa.gov/fms-cgi-bin/cvsweb.cgi/FMS/
22 !> @brief A software package that provides a set of simple interfaces for
23 !! modelers to perform computations related to time and dates.
24 !> The changes between the lima revision and this revision are more
25 !! extensive that all those between antwerp and lima.
26 !! A brief description of these changes follows.
28 !! <LI> Added option to set the smallest time increment to something less than one second.
29 !! This is controlled by calling the pubic subroutine set_ticks_per_second. </LI>
31 !! <LI>Gregorian calendar fixed.</LI>
33 !! <LI>Optional error flag added to calling arguments of public routines.
34 !! This allows the using routine to terminate the program. It is likely that more
35 !! diagnostic information is available from the user than from time_manager alone.
36 !! If the error flag is present then it is the responsibility of the using
37 !! routine to test it and add additional information to the error message.</LI>
39 !! <LI>Removed the restriction that time increments be positive in routines that increment or decrement
40 !! time and date. The option to prohibit negative increments can be turned on via optional argument.</LI>
42 !! <LI>subroutine set_date_c modified to handle strings that include only hours or only hours and minutes.
43 !! This complies with CF convensions.</LI>
45 !! <LI>Made calendar specific routines private.
46 !! They are not used, and should not be used, by any using code.</LI>
48 !! <LI>Error messages made more informative.</LI>
50 !! The module defines a type that can be used to represent discrete
51 !! times (accurate to one second) and to map these times into dates
52 !! using a variety of calendars. A time is mapped to a date by
53 !! representing the time with respect to an arbitrary base date (refer
54 !! to <B>NOTES</B> section for the base date setting).
56 !! The time_manager provides a single defined type, time_type, which is
57 !! used to store time and date quantities. A time_type is a positive
58 !! definite quantity that represents an interval of time. It can be
59 !! most easily thought of as representing the number of seconds in some
60 !! time interval. A time interval can be mapped to a date under a given
61 !! calendar definition by using it to represent the time that has passed
62 !! since some base date. A number of interfaces are provided to operate
63 !! on time_type variables and their associated calendars. Time intervals
64 !! can be as large as n days where n is the largest number represented by
65 !! the default integer type on a compiler. This is typically considerably
66 !! greater than 10 million years (assuming 32 bit integer representation)
67 !! which is likely to be adequate for most applications. The description
68 !! of the interfaces is separated into two sections. The first deals with
69 !! operations on time intervals while the second deals with operations
70 !! that convert time intervals to dates for a given calendar.
72 !! The smallest increment of time is referred to as a tick.
73 !! A tick cannot be larger than 1 second, which also is the default.
74 !! The number of ticks per second is set via pubic subroutine set_ticks_per_second.
75 !! For example, ticks_per_second = 1000 will set the tick to one millisecond.
77 !! <DATA NAME="time_type" TYPE="derived type">
78 !! Derived-type data variable used to store time and date quantities. It
79 !! contains three PRIVATE variables: days, seconds and ticks.
83 !> @brief File for @ref time_manager_mod
85 !> @addtogroup time_manager_mod
87 module time_manager_mod
90 use platform_mod, only: r8_kind
91 use constants_mod, only: rseconds_per_day=>seconds_per_day
92 use fms_mod, only: error_mesg, FATAL, WARNING, write_version_number, stdout
97 ! Module defines a single type
100 ! Operators defined on time_type
101 public operator(+), operator(-), operator(*), operator(/), &
102 operator(>), operator(>=), operator(==), operator(/=), &
103 operator(<), operator(<=), operator(//), assignment(=)
105 ! Subroutines and functions operating on time_type
106 public set_time, increment_time, decrement_time, get_time, interval_alarm
107 public repeat_alarm, time_type_to_real, real_to_time_type
108 public time_list_error
110 ! List of available calendar types
111 public THIRTY_DAY_MONTHS, JULIAN, GREGORIAN, NOLEAP, NO_CALENDAR, INVALID_CALENDAR
113 ! Subroutines and functions involving relations between time and calendar
114 public set_calendar_type
115 public get_calendar_type
116 public set_ticks_per_second
117 public get_ticks_per_second
120 public increment_date
121 public decrement_date
124 public length_of_year
129 public valid_calendar_types
131 ! Subroutines for printing version number and time type
132 public :: time_manager_init, print_time, print_date
134 ! The following exist only for interpolator.F90
135 ! interpolator.F90 uses them to do a calendar conversion,
136 ! which is also done by get_cal_time. interpolator.F90
137 ! should be modified to use get_cal_time instead.
138 ! After interpolator.F90 is fixed, these can be removed
139 ! and the corresponding private routines can be renamed.
140 ! (e.g., rename set_date_julian_private to be just set_date_julian)
141 public :: set_date_julian, set_date_no_leap, get_date_julian, get_date_no_leap
143 public :: date_to_string
145 !====================================================================
147 ! Global data to define calendar type
148 integer, parameter :: THIRTY_DAY_MONTHS = 1, JULIAN = 2, &
149 GREGORIAN = 3, NOLEAP = 4, &
150 NO_CALENDAR = 0, INVALID_CALENDAR =-1
151 integer, private :: calendar_type = NO_CALENDAR
152 integer, parameter :: max_type = 4
154 ! Define number of days per month
155 integer, private :: days_per_month(12) = (/31,28,31,30,31,30,31,31,30,31,30,31/)
156 integer, parameter :: seconds_per_day = rseconds_per_day ! This should automatically cast real to integer
157 integer, parameter :: days_in_400_year_period = 146097 !< Used only for gregorian
158 integer,parameter :: do_floor = 0
159 integer,parameter :: do_nearest = 1
163 !> @brief Type to represent amounts of time.
164 !> Implemented as seconds and days to allow for larger intervals.
165 !> @ingroup time_manager_mod
173 !> Operator override interface for use with @ref time_type
174 !> @ingroup time_manager_mod
175 interface operator (+); module procedure time_plus; end interface
176 !> Operator override interface for use with @ref time_type
177 !> @ingroup time_manager_mod
178 interface operator (-); module procedure time_minus; end interface
179 !> Operator override interface for use with @ref time_type
180 !> @ingroup time_manager_mod
181 interface operator (*); module procedure time_scalar_mult
182 module procedure scalar_time_mult; end interface
183 !> Operator override interface for use with @ref time_type
184 !> @ingroup time_manager_mod
185 interface operator (/); module procedure time_scalar_divide
186 module procedure time_divide; end interface
187 !> Operator override interface for use with @ref time_type
188 !> @ingroup time_manager_mod
189 interface operator (>); module procedure time_gt; end interface
190 !> Operator override interface for use with @ref time_type
191 !> @ingroup time_manager_mod
192 interface operator (>=); module procedure time_ge; end interface
193 !> Operator override interface for use with @ref time_type
194 !> @ingroup time_manager_mod
195 interface operator (<); module procedure time_lt; end interface
196 !> Operator override interface for use with @ref time_type
197 !> @ingroup time_manager_mod
198 interface operator (<=); module procedure time_le; end interface
199 !> Operator override interface for use with @ref time_type
200 !> @ingroup time_manager_mod
201 interface operator (==); module procedure time_eq; end interface
202 !> Operator override interface for use with @ref time_type
203 !> @ingroup time_manager_mod
204 interface operator (/=); module procedure time_ne; end interface
205 !> Operator override interface for use with @ref time_type
206 !> @ingroup time_manager_mod
207 interface operator (//); module procedure time_real_divide; end interface
208 !> Operator override interface for use with @ref time_type
209 !> @ingroup time_manager_mod
210 interface assignment(=); module procedure time_assignment; end interface
212 !======================================================================
214 !> @brief Given some number of seconds and days, returns the
215 !! corresponding time_type.
217 !> Given some number of seconds and days, returns the
218 !! corresponding time_type. set_time has two forms;
219 !! one accepts integer input, the other a character string with the day and second counts.
220 !! For the first form, there are no restrictions on the range of the inputs,
221 !! except that the result must be positive time.
222 !! e.g. days=-1, seconds=86401 is acceptable.
223 !! For the second form, days and seconds must both be positive.
225 !! <br>Example usage:
227 !! type(time_type) :: time1, time2
228 !! time1 = set_time(seconds, days, ticks, err_msg)
229 !! time2 = set_time("100 43200", err_msg, allow_rounding)
231 !> @ingroup time_manager_mod
233 module procedure set_time_i, set_time_c
236 !> @brief Given an input date in year, month, days, etc., creates a
237 !! time_type that represents this time interval from the
238 !! internally defined base date.
240 !> Given a date, computes the corresponding time given the selected
241 !! date time mapping algorithm. Note that it is possible to specify
242 !! any number of illegal dates; these should be checked for and generate
243 !! errors as appropriate.
245 !! <br>Example usage:
246 !! <br> Integer input
247 !! @code{.F90} time = set_date(year, month, day, hours, minute, second, tick, err_msg) @endcode
249 !! @code{.F90} time = set_date_c(time_string, zero_year_warning, err_msg, allow_rounding) @endcode
251 !! @param time_string A character string containing a date formatted
252 !! according to CF conventions. e.g. '1980-12-31 23:59:59.9'
254 !! @param zero_year_warning If the year number is zero, it will be silently changed to one,
255 !! unless zero_year_warning=.true., in which case a WARNING message
256 !! will also be issued.
258 !! @param allow_rounding When .true., any fractions of a second will be rounded off to the nearest
259 !! tick. When .false., it is a fatal error if the second fraction cannot be exactly
260 !! represented by a number of ticks.
262 !! @param err_msg When present, and when non-blank, a fatal error condition as been detected.
263 !! The string itself is an error message.
264 !! It is recommended that, when err_msg is present in the call
265 !! to this routine, the next line of code should be something
268 !! if(err_msg /= '') call error_mesg('my_routine','additional info: '//trim(err_msg) ,FATAL)
271 !> @ingroup time_manager_mod
273 module procedure set_date_i, set_date_c
276 !> @addtogroup time_manager_mod
279 !======================================================================
281 ! Include variable "version" to be written to log file.
282 #include<file_version.h>
283 logical :: module_is_initialized = .false.
285 !======================================================================
287 ! A tick is the smallest increment of time.
288 ! That is, smallest increment of time = (1/ticks_per_second) seconds
290 integer :: ticks_per_second = 1
292 !======================================================================
295 ! First define all operations on time intervals independent of calendar
297 !> Returns a time interval corresponding to this number of days, seconds, and ticks.
298 !! days, seconds and ticks may be negative, but resulting time must be positive.
299 function set_time_private(seconds, days, ticks, Time_out, err_msg)
302 ! To understand why inputs may be negative,
303 ! one needs to understand the intrinsic function "modulo".
304 ! The expanation below is copied from a web page on fortran 90
306 ! In addition, CEILING, FLOOR and MODULO have been added to Fortran 90.
307 ! Only the last one is difficult to explain, which is most easily done with the examples from ISO (1991)
309 ! MOD (8,5) gives 3 MODULO (8,5) gives 3
310 ! MOD (-8,5) gives -3 MODULO (-8,5) gives 2
311 ! MOD (8,-5) gives 3 MODULO (8,-5) gives -2
312 ! MOD (-8,-5) gives -3 MODULO (-8,-5) gives -3
314 ! I don't think it is difficult to explain.
315 ! I think that is it sufficient to say this:
316 ! "The result of modulo(n,m) has the sign of m"
319 logical :: set_time_private
320 integer, intent(in) :: seconds, days, ticks
321 type(time_type), intent(out) :: Time_out
322 character(len=*), intent(out) :: err_msg
323 integer :: seconds_new, days_new, ticks_new
325 seconds_new = seconds + floor(ticks/real(ticks_per_second))
326 ticks_new = modulo(ticks,ticks_per_second)
327 days_new = days + floor(seconds_new/real(seconds_per_day))
328 seconds_new = modulo(seconds_new,seconds_per_day)
330 if ( seconds_new < 0 .or. ticks_new < 0) then
331 call error_mesg('function set_time_i','Bad result for time. Contact those responsible for maintaining time_manager'&
335 if(days_new < 0) then
336 write(err_msg,'(a,i6,a,i6,a,i6)') 'time is negative. days=',days_new,' seconds=',seconds_new,' ticks=',ticks_new
337 set_time_private = .false.
339 Time_out%days = days_new
340 Time_out%seconds = seconds_new
341 Time_out%ticks = ticks_new
343 set_time_private = .true.
346 end function set_time_private
347 !---------------------------------------------------------------------------
348 !> @brief Returns a time_type set to the given amount of time via integer amounts.
349 function set_time_i(seconds, days, ticks, err_msg)
350 type(time_type) :: set_time_i
351 integer, intent(in) :: seconds !< A number of seconds
352 integer, intent(in), optional :: days !< A number of days
353 integer, intent(in), optional :: ticks !< A number of ticks
354 character(len=*), intent(out), optional :: err_msg !< When present, and when non-blank, a fatal
355 !! error condition as been detected. The string itself is an error message.
356 !! It is recommended that, when err_msg is present in the call
357 !! to this routine, the next line of code should be something
359 !! if(err_msg /= '') call error_mesg('my_routine','additional info: '//trim(err_msg),FATAL)
360 character(len=128) :: err_msg_local
361 integer :: odays, oticks
363 if(.not.module_is_initialized) call time_manager_init
365 odays = 0; if(present(days)) odays = days
366 oticks = 0; if(present(ticks)) oticks = ticks
367 if(present(err_msg)) err_msg = ''
369 if(.not.set_time_private(seconds, odays, oticks, set_time_i, err_msg_local)) then
370 if(error_handler('function set_time_i', trim(err_msg_local), err_msg)) return
373 end function set_time_i
374 !---------------------------------------------------------------------------
376 !> @brief Returns a time_type set to the given amount of time via a string
377 function set_time_c(string, err_msg, allow_rounding)
379 type(time_type) :: set_time_c
380 character(len=*), intent(in) :: string !< Contains days and seconds separated by a single blank.
381 !! days must be integer, seconds may be integer or real.
382 !! Examples: '100 43200' '100 43200.50'
383 character(len=*), intent(out), optional :: err_msg !< When present, and when non-blank, a fatal
384 !! error condition as been detected. The string itself is an error message.
385 !! It is recommended that, when err_msg is present in the call
386 !! to this routine, the next line of code should be something
388 !! if(err_msg /= '') call error_mesg('my_routine','additional info: '//trim(err_msg),FATAL)
389 logical, intent(in), optional :: allow_rounding !< When .true., any fractions of a second will be
390 !! rounded off to the nearest tick. When .false., it is a fatal error
391 !! if the second fraction cannot be exactly represented by a number of ticks.
393 character(len=4) :: formt='(i )'
394 integer :: i1, i2, i3, day, second, tick, nsps
395 character(len=32) :: string_sifted_left
396 character(len=128) :: err_msg_local
397 logical :: allow_rounding_local
399 if(.not.module_is_initialized) call time_manager_init
400 if(present(err_msg)) err_msg = ''
401 allow_rounding_local=.true.; if(present(allow_rounding)) allow_rounding_local=allow_rounding
403 err_msg_local = 'Form of character time stamp is incorrect. The character time stamp is: '//trim(string)
405 string_sifted_left = adjustl(string)
406 i1 = index(trim(string_sifted_left),' ')
408 if(error_handler('function set_time_c', err_msg_local, err_msg)) return
410 if(index(string,'-') /= 0 .or. index(string,':') /= 0) then
411 if(error_handler('function set_time_c', err_msg_local, err_msg)) return
414 i2 = index(trim(string_sifted_left),'.')
415 i3 = len_trim(cut0(string_sifted_left))
417 if(i2 /= 0) then ! There is no decimal point
418 ! Check that decimal is on seconds (not days)
420 if(error_handler('function set_time_c', err_msg_local, err_msg)) return
423 write(formt(3:3),'(i1)') i1-1
424 read(string_sifted_left(1:i1-1),formt) day
426 if(i2 == 0) then ! There is no decimal point
427 write(formt(3:3),'(i1)') i3-i1
428 read(string_sifted_left(i1+1:i3),formt) second
430 else ! There is a decimal point
431 ! nsps = spaces occupied by whole number of seconds
436 write(formt(3:3),'(i1)') nsps
437 read(string_sifted_left(i1+1:i2-1),formt) second
440 if(.not.get_tick_from_string(string_sifted_left(i2:i3), err_msg_local, allow_rounding_local, tick)) then
441 if(error_handler('function set_time_c', err_msg_local, err_msg)) return
443 ! If tick has been rounded up to ticks_per_second, then bump up second.
444 if(tick == ticks_per_second) then
450 if(.not.set_time_private(second, day, tick, set_time_c, err_msg_local)) then
451 if(error_handler('function set_time_c', err_msg_local, err_msg)) return
454 end function set_time_c
455 !---------------------------------------------------------------------------
457 function get_tick_from_string(string, err_msg, allow_rounding, tick)
459 logical :: get_tick_from_string
460 character(len=*), intent(in) :: string
461 character(len=*), intent(out) :: err_msg
462 logical, intent(in) :: allow_rounding
463 integer, intent(out) :: tick
465 character(len=4) :: formt='(i )'
466 integer :: i3, nspf, fraction, magnitude, tpsfrac
469 get_tick_from_string = .true.
470 i3 = len_trim(string)
471 nspf = i3 - 1 ! nspf = spaces occupied by fractional seconds, excluding decimal point
473 tick = 0 ! Nothing to the right of the decimal point
475 write(formt(3:3),'(i1)') nspf
476 read(string(2:i3),formt) fraction
477 if(fraction == 0) then
478 tick = 0 ! All zeros to the right of the decimal point
481 tpsfrac = ticks_per_second*fraction
482 if(allow_rounding) then
483 tick = nint((real(tpsfrac)/magnitude))
485 if(modulo(tpsfrac,magnitude) == 0) then
486 tick = tpsfrac/magnitude
488 write(err_msg,'(a,i6)') 'Second fraction cannot be exactly represented with ticks. '// &
489 'fraction='//trim(string)//' ticks_per_second=',ticks_per_second
490 get_tick_from_string = .false.
496 end function get_tick_from_string
498 !> @brief Returns days and seconds ( < 86400 ) corresponding to a time.
499 !! <TT>err_msg</TT> should be checked for any errors.
501 !> @param time A @ref time_type interval to get days and seconds from
502 !! @param [out] seconds The number of seconds
503 !! @param [out] days The number of seconds
504 !! @param [out] ticks The number of ticks
505 !! @param [out] err_msg Contains an error message on failure
506 !! <br>Example usage:
507 !! @code{.F90} get_time(time, seconds, days, ticks, err_msg) @endcode
508 subroutine get_time(Time, seconds, days, ticks, err_msg)
510 ! Returns days and seconds ( < 86400 ) corresponding to a time.
512 type(time_type), intent(in) :: Time
513 integer, intent(out) :: seconds
514 integer, intent(out), optional :: days, ticks
515 character(len=*), intent(out), optional :: err_msg
516 character(len=128) :: err_msg_local
518 if(.not.module_is_initialized) call time_manager_init
519 if(present(err_msg)) err_msg = ''
521 seconds = Time%seconds
523 if(present(ticks)) then
526 if(Time%ticks /= 0) then
527 err_msg_local = 'subroutine get_time: ticks must be present when time has a second fraction'
528 if(error_handler('subroutine get_time', err_msg_local, err_msg)) return
532 if (present(days)) then
535 if (Time%days > (huge(seconds) - seconds)/seconds_per_day) then
536 err_msg_local = 'Integer overflow in seconds. Optional argument days must be present.'
537 if(error_handler('subroutine get_time', err_msg_local, err_msg)) return
539 seconds = seconds + Time%days * seconds_per_day
542 end subroutine get_time
544 !> @brief Increments a time by seconds and days.
546 !> Given a time and an increment of days and seconds, returns
547 !! a new time_type that represents the given time after the given increment.
548 !! @returns incremented @ref time_type
549 function increment_time(Time, seconds, days, ticks, err_msg, allow_neg_inc)
551 type(time_type) :: increment_time
552 type(time_type), intent(in) :: Time !< A time interval
553 integer, intent(in) :: seconds !< Increment of seconds
554 integer, intent(in), optional :: days, ticks !< Increment of days and ticks
555 character(len=*), intent(out), optional :: err_msg !< When present and non-blank, a fatal error
556 !! condition has been detected, with the string itself
557 !! as the error message.
558 logical, intent(in), optional :: allow_neg_inc !< When false, negative increments give fatal errors
560 integer :: odays, oticks
561 character(len=128) :: err_msg_local
562 logical :: allow_neg_inc_local
564 odays = 0; if(present(days)) odays = days
565 oticks = 0; if(present(ticks)) oticks = ticks
566 allow_neg_inc_local=.true.; if(present(allow_neg_inc)) allow_neg_inc_local=allow_neg_inc
568 if(.not.allow_neg_inc_local) then
569 if(seconds < 0 .or. odays < 0 .or. oticks < 0) then
570 write(err_msg_local,10) seconds, odays, oticks
571 10 format('One or more time increments are negative: seconds=',i6,' days=',i6,' ticks=',i6)
572 if(error_handler('function increment_time', err_msg_local, err_msg)) return
576 if(.not.increment_time_private(Time, seconds, odays, oticks, increment_time, err_msg_local)) then
577 if(error_handler('function increment_time', err_msg_local, err_msg)) return
580 end function increment_time
581 !--------------------------------------------------------------------------
583 !> Increments a time by seconds, days and ticks.
584 function increment_time_private(Time_in, seconds, days, ticks, Time_out, err_msg)
587 logical :: increment_time_private
588 type(time_type), intent(in) :: Time_in
589 integer, intent(in) :: seconds, days, ticks
590 type(time_type), intent(out) :: Time_out
591 character(len=*), intent(out) :: err_msg
593 ! Watch for immediate overflow on days or seconds
594 if(days >= huge(days) - Time_in%days) then
595 err_msg = 'Integer overflow in days in increment_time'
596 increment_time_private = .false.
599 if(seconds >= huge(seconds) - Time_in%seconds) then
600 err_msg = 'Integer overflow in seconds in increment_time'
601 increment_time_private = .false.
605 increment_time_private = set_time_private(Time_in%seconds+seconds, Time_in%days+days, Time_in%ticks+ticks, &
608 end function increment_time_private
610 !--------------------------------------------------------------------------
612 !> @brief Decrements a time by seconds and days.
614 !> Given a time and a decrement of days and seconds, returns
615 !! a time that subtracts this decrement from an input time.
616 !> @returns A time that suvtracts this decrement from an input time. A negative result is a fatal error.
617 function decrement_time(Time, seconds, days, ticks, err_msg, allow_neg_inc)
619 ! Decrements a time by seconds, days and ticks.
621 type(time_type) :: decrement_time
622 type(time_type), intent(in) :: Time !< A time interval
623 integer, intent(in) :: seconds !< Decrement of seconds
624 integer, intent(in), optional :: days, ticks !< Decrement of days and ticks
625 character(len=*), intent(out), optional :: err_msg !< Present and non-blank when a fatal error has
626 !! occured, holds the error message.
627 logical, intent(in), optional :: allow_neg_inc !< Throws fatal warning when set to false if
628 !! negative values are used to decrement. Default
631 integer :: odays, oticks
632 character(len=128) :: err_msg_local
633 logical :: allow_neg_inc_local
635 odays = 0; if (present(days)) odays = days
636 oticks = 0; if (present(ticks)) oticks = ticks
637 allow_neg_inc_local=.true.; if(present(allow_neg_inc)) allow_neg_inc_local=allow_neg_inc
639 if(.not.allow_neg_inc_local) then
640 if(seconds < 0 .or. odays < 0 .or. oticks < 0) then
641 write(err_msg_local,10) seconds,odays,oticks
642 10 format('One or more time increments are negative: seconds=',i6,' days=',i6,' ticks=',i6)
643 if(error_handler('function decrement_time', err_msg_local, err_msg)) return
647 if(.not.increment_time_private(Time, -seconds, -odays, -oticks, decrement_time, err_msg_local)) then
648 if(error_handler('function decrement_time', err_msg_local, err_msg)) return
651 end function decrement_time
653 !--------------------------------------------------------------------------
654 ! <FUNCTION NAME="time_gt operator(>)">
657 ! Returns true if time1 > time2.
660 ! Returns true if time1 > time2.
662 ! <IN NAME="time1" UNITS="" TYPE="time_type" DIM="">
665 ! <IN NAME="time2" UNITS="" TYPE="time_type" DIM="">
668 ! <OUT NAME="" UNITS="" TYPE="logical" DIM="" DEFAULT="">
669 ! Returns true if time1 > time2
672 ! time_gt(time1, time2)
675 !> @brief Returns true if time1 > time2
676 function time_gt(time1, time2)
679 type(time_type), intent(in) :: time1, time2
681 time_gt = (time1%days > time2%days)
682 if(time1%days == time2%days) then
683 if(time1%seconds == time2%seconds) then
684 time_gt = (time1%ticks > time2%ticks)
686 time_gt = (time1%seconds > time2%seconds)
693 !--------------------------------------------------------------------------
694 ! <FUNCTION NAME="time_ge; operator(>=)">
697 ! Returns true if time1 >= time2.
700 ! Returns true if time1 >= time2.
703 ! time_ge(time1, time2)
706 ! <IN NAME="time1" UNITS="" TYPE="time_type" DIM="">
709 ! <IN NAME="time2" UNITS="" TYPE="time_type" DIM="">
712 ! <OUT NAME="" UNITS="" TYPE="logical" DIM="" DEFAULT="">
713 ! Returns true if time1 >= time2
716 function time_ge(time1, time2)
718 ! Returns true if time1 >= time2
721 type(time_type), intent(in) :: time1, time2
723 time_ge = (time_gt(time1, time2) .or. time_eq(time1, time2))
728 !--------------------------------------------------------------------------
729 ! <FUNCTION NAME="time_lt; operator(<)">
732 ! Returns true if time1 < time2.
735 ! Returns true if time1 < time2.
738 ! time_lt(time1, time2)
741 ! <IN NAME="time1" UNITS="" TYPE="time_type" DIM="">
744 ! <IN NAME="time2" UNITS="" TYPE="time_type" DIM="">
747 ! <OUT NAME="" UNITS="" TYPE="logical" DIM="" DEFAULT="">
748 ! Returns true if time1 < time2
751 function time_lt(time1, time2)
753 ! Returns true if time1 < time2
756 type(time_type), intent(in) :: time1, time2
758 time_lt = (time1%days < time2%days)
759 if(time1%days == time2%days)then
760 if(time1%seconds == time2%seconds) then
761 time_lt = (time1%ticks < time2%ticks)
763 time_lt = (time1%seconds < time2%seconds)
769 !--------------------------------------------------------------------------
770 ! <FUNCTION NAME="time_le; operator(<=)">
773 ! Returns true if time1 <= time2.
776 ! Returns true if time1 <= time2.
779 ! time_le(time1, time2)
782 ! <IN NAME="time1" UNITS="" TYPE="time_type" DIM="">
785 ! <IN NAME="time2" UNITS="" TYPE="time_type" DIM="">
788 ! <OUT NAME="" UNITS="" TYPE="logical" DIM="" DEFAULT="">
789 ! Returns true if time1 <= time2
792 !> Returns true if time1 <= time2
793 function time_le(time1, time2)
797 type(time_type), intent(in) :: time1, time2
799 time_le = (time_lt(time1, time2) .or. time_eq(time1, time2))
804 !--------------------------------------------------------------------------
805 ! <FUNCTION NAME="time_eq; operator(==)">
808 ! Returns true if time1 == time2.
811 ! Returns true if time1 == time2.
814 ! time_eq(time1, time2)
817 ! <IN NAME="time1" UNITS="" TYPE="time_type" DIM="">
820 ! <IN NAME="time2" UNITS="" TYPE="time_type" DIM="">
823 ! <OUT NAME="" UNITS="" TYPE="logical" DIM="" DEFAULT="">
824 ! Returns true if time1 == time2
827 function time_eq(time1, time2)
829 ! Returns true if time1 == time2
832 type(time_type), intent(in) :: time1, time2
834 if(.not.module_is_initialized) call time_manager_init
836 time_eq = (time1%seconds == time2%seconds .and. time1%days == time2%days &
837 .and. time1%ticks == time2%ticks)
842 !--------------------------------------------------------------------------
843 ! <FUNCTION NAME="time_ne; operator(/=)">
846 ! Returns true if time1 /= time2.
849 ! Returns true if time1 /= time2.
852 ! time_ne(time1, time2)
855 ! <IN NAME="time1" UNITS="" TYPE="time_type" DIM="">
858 ! <IN NAME="time2" UNITS="" TYPE="time_type" DIM="">
861 ! <OUT NAME="" UNITS="" TYPE="logical" DIM="" DEFAULT="">
862 ! Returns true if time1 /= time2
865 !> Returns true if time1 /= time2
866 function time_ne(time1, time2)
869 type(time_type), intent(in) :: time1, time2
871 time_ne = (.not. time_eq(time1, time2))
876 !-------------------------------------------------------------------------
877 ! <FUNCTION NAME="time_plus; operator(+)">
880 ! Returns sum of two time_types.
886 ! Returns sum of two time_types.
889 ! <IN NAME="time1" UNITS="" TYPE="time_type" DIM="">
892 ! <IN NAME="time2" UNITS="" TYPE="time_type" DIM="">
895 ! <OUT NAME="" UNITS="" TYPE="time_type" DIM="" DEFAULT="">
896 ! Returns sum of two time_types.
899 !> Returns sum of two time_types
900 function time_plus(time1, time2)
903 type(time_type) :: time_plus
904 type(time_type), intent(in) :: time1, time2
906 if(.not.module_is_initialized) call time_manager_init
908 time_plus = increment_time(time1, time2%seconds, time2%days, time2%ticks)
910 end function time_plus
913 !-------------------------------------------------------------------------
914 ! <FUNCTION NAME="time_minus; operator(-)">
917 ! Returns difference of two time_types.
920 ! Returns difference of two time_types. WARNING: a time type is positive
921 ! so by definition time1 - time2 is the same as time2 - time1.
924 ! time_minus(time1, time2)
930 ! <IN NAME="time1" UNITS="" TYPE="time_type" DIM="">
933 ! <IN NAME="time2" UNITS="" TYPE="time_type" DIM="">
936 ! <OUT NAME="" UNITS="" TYPE="time_type" DIM="" DEFAULT="">
937 ! Returns difference of two time_types.
940 !> Returns difference of two time_types. WARNING: a time type is positive
941 !! so by definition time1 - time2 is the same as time2 - time1.
942 function time_minus(time1, time2)
945 type(time_type) :: time_minus
946 type(time_type), intent(in) :: time1, time2
948 if(.not.module_is_initialized) call time_manager_init
950 if(time1 > time2) then
951 time_minus = decrement_time(time1, time2%seconds, time2%days, time2%ticks)
953 time_minus = decrement_time(time2, time1%seconds, time1%days, time1%ticks)
956 end function time_minus
959 !--------------------------------------------------------------------------
960 ! <FUNCTION NAME="time_scalar_mult; operator(*)">
963 ! Returns time multiplied by integer factor n.
966 ! Returns time multiplied by integer factor n.
969 ! time_scalar_mult(time, n)
972 ! <IN NAME="time" UNITS="" TYPE="time_type" DIM="">
975 ! <IN NAME="n" UNITS="" TYPE="integer" DIM="">
978 ! <OUT NAME="" UNITS="" TYPE="time_type" DIM="" DEFAULT="">
979 ! Returns time multiplied by integer factor n.
982 !> Returns time multiplied by integer factor n
983 function time_scalar_mult(time, n)
986 type(time_type) :: time_scalar_mult
987 type(time_type), intent(in) :: time
988 integer, intent(in) :: n
989 integer :: days, seconds, ticks, num_sec
990 double precision :: sec_prod, tick_prod
992 if(.not.module_is_initialized) call time_manager_init
994 ! Multiplying here in a reasonable fashion to avoid overflow is tricky
995 ! Could multiply by some large factor n, and seconds could be up to 86399
996 ! Need to avoid overflowing integers and wrapping around to negatives
997 ! ticks could be up to ticks_per_second-1
999 tick_prod = dble(time%ticks) * dble(n)
1000 num_sec = int(tick_prod/dble(ticks_per_second))
1001 sec_prod = dble(time%seconds) * dble(n) + num_sec
1002 ticks = int(tick_prod - num_sec * ticks_per_second)
1004 ! If sec_prod is large compared to precision of double precision, things
1005 ! can go bad. Need to warn and abort on this.
1006 ! The same is true of tick_prod but is is more likely to happen to sec_prod,
1007 ! so let's just test sec_prod. (A test of tick_prod would be necessary only
1008 ! if ticks_per_second were greater than seconds_per_day)
1009 if(sec_prod /= 0.0) then
1010 if(log10(sec_prod) > precision(sec_prod) - 3) call error_mesg('time_scalar_mult', &
1011 'Insufficient precision to handle scalar product in time_scalar_mult; contact developer',FATAL)
1014 days = int(sec_prod / dble(seconds_per_day))
1015 seconds = int(sec_prod - dble(days) * dble(seconds_per_day))
1017 time_scalar_mult = set_time(seconds, time%days * n + days, ticks)
1019 end function time_scalar_mult
1022 !-------------------------------------------------------------------------
1023 ! <FUNCTION NAME="scalar_time_mult; operator(*)">
1026 ! Returns time multiplied by integer factor n.
1029 ! Returns time multiplied by integer factor n.
1033 ! scalar_time_mult(n, time)
1036 ! <IN NAME="time" UNITS="" TYPE="time_type" DIM="">A time interval.</IN>
1037 ! <IN NAME="n" UNITS="" TYPE="integer" DIM=""> An integer. </IN>
1038 ! <OUT NAME="" UNITS="" TYPE="time_type" DIM="" DEFAULT="">
1039 ! Returns time multiplied by integer factor n.
1042 !> Returns time multipled by integer factor n
1043 function scalar_time_mult(n, time)
1045 type(time_type) :: scalar_time_mult
1046 type(time_type), intent(in) :: time
1047 integer, intent(in) :: n
1049 scalar_time_mult = time_scalar_mult(time, n)
1051 end function scalar_time_mult
1054 !-------------------------------------------------------------------------
1055 ! <FUNCTION NAME="time_divide; operator(/)">
1058 ! Returns the largest integer, n, for which time1 >= time2 * n.
1061 ! Returns the largest integer, n, for which time1 >= time2 * n.
1065 ! time_divide(time1, time2)
1068 ! <IN NAME="time1" UNITS="" TYPE="time_type" DIM="">
1071 ! <IN NAME="time2" UNITS="" TYPE="time_type" DIM="">
1074 ! <OUT NAME="" UNITS="" TYPE="integer" DIM="" DEFAULT="">
1075 ! Returns the largest integer, n, for which time1 >= time2 * n.
1078 !> Returns the largest integer, n, for which time1 >= time2 * n.
1079 function time_divide(time1, time2)
1081 integer :: time_divide
1082 type(time_type), intent(in) :: time1, time2
1083 double precision :: d1, d2
1085 if(.not.module_is_initialized) call time_manager_init
1087 ! Convert time intervals to floating point days; risky for general performance?
1088 d1 = time1%days * dble(seconds_per_day) + dble(time1%seconds) + time1%ticks/dble(ticks_per_second)
1089 d2 = time2%days * dble(seconds_per_day) + dble(time2%seconds) + time2%ticks/dble(ticks_per_second)
1091 ! Get integer quotient of this, check carefully to avoid round-off problems.
1092 time_divide = int(d1 / d2)
1094 ! Verify time_divide*time2 is <= time1 and (time_divide + 1)*time2 is > time1
1095 if(time_divide * time2 > time1 .or. (time_divide + 1) * time2 <= time1) &
1096 call error_mesg('time_divide',' quotient error :: notify developer',FATAL)
1098 end function time_divide
1101 !-------------------------------------------------------------------------
1102 ! <FUNCTION NAME="time_real_divide; operator(//)">
1105 ! Returns the double precision quotient of two times.
1108 ! Returns the double precision quotient of two times.
1112 ! time_real_divide(time1, time2)
1115 ! <IN NAME="time1" UNITS="" TYPE="time_type" DIM="">
1118 ! <IN NAME="time2" UNITS="" TYPE="time_type" DIM="">
1121 ! <OUT NAME="" UNITS="" TYPE="integer" DIM="double precision" DEFAULT="">
1122 ! Returns the double precision quotient of two times
1125 !> Returns the double precision quotient of two times
1126 function time_real_divide(time1, time2)
1128 double precision :: time_real_divide
1129 type(time_type), intent(in) :: time1, time2
1130 double precision :: d1, d2
1132 if(.not.module_is_initialized) call time_manager_init
1134 ! Convert time intervals to floating point seconds; risky for general performance?
1135 d1 = time1%days * dble(seconds_per_day) + dble(time1%seconds) + dble(time1%ticks)/dble(ticks_per_second)
1136 d2 = time2%days * dble(seconds_per_day) + dble(time2%seconds) + dble(time2%ticks)/dble(ticks_per_second)
1138 time_real_divide = d1 / d2
1140 end function time_real_divide
1143 !-------------------------------------------------------------------------
1144 ! <SUBROUTINE NAME="time_assignment; assignment(=)">
1147 ! Assigns all components of the time_type variable on
1148 ! RHS to same components of time_type variable on LHS.
1151 ! Assigns all components of the time_type variable on
1152 ! RHS to same components of time_type variable on LHS.
1158 ! <OUT NAME="time1" UNITS="" TYPE="time_type" DIM="">
1159 ! A time type variable.
1161 ! <IN NAME="time2" UNITS="" TYPE="time_type" DIM="">
1162 ! A time type variable.
1165 !> Assigns all components of the time_type variable on
1166 !! RHS to same components of time_type variable on LHS.
1167 subroutine time_assignment(time1, time2)
1168 type(time_type), intent(out) :: time1
1169 type(time_type), intent(in) :: time2
1170 time1%seconds = time2%seconds
1171 time1%days = time2%days
1172 time1%ticks = time2%ticks
1173 end subroutine time_assignment
1176 !-------------------------------------------------------------------------
1177 ! <FUNCTION NAME="time_type_to_real">
1179 ! Converts time to seconds and returns it as a real number
1182 ! Converts time to seconds and returns it as a real number
1185 ! time_type_to_real(time)
1187 ! <IN NAME="time" UNITS="" TYPE="time_type" DIM="">
1191 function time_type_to_real(time)
1193 real(kind=r8_kind) :: time_type_to_real
1194 type(time_type), intent(in) :: time
1196 if(.not.module_is_initialized) call time_manager_init
1198 time_type_to_real = real(dble(time%days) * 86400.d0 + dble(time%seconds) + &
1199 dble(time%ticks)/dble(ticks_per_second), kind=r8_kind)
1201 end function time_type_to_real
1203 !> @brief Convert a real number of seconds into a time_type variable.
1204 !! @return A filled time type variable, and an error message if an
1206 function real_to_time_type(x,err_msg) result(t)
1207 real,intent(in) :: x !< Number of seconds.
1208 character(len=*),intent(out),optional :: err_msg !< Error message.
1209 type(time_type) :: t
1213 character(len=128) :: err_msg_local
1214 real,parameter :: spd = real(86400)
1217 tps = real(ticks_per_second)
1219 days = safe_rtoi(a,do_floor)
1220 a = x - real(days)*spd
1221 seconds = safe_rtoi(a,do_floor)
1222 a = (a - real(seconds))*tps
1223 ticks = safe_rtoi(a,do_nearest)
1224 if (.not. set_time_private(seconds,days,ticks,t,err_msg_local)) then
1225 if (error_handler('function real_to_time_type',err_msg_local,err_msg)) then
1229 end function real_to_time_type
1231 !> @brief Convert a floating point value to an integer value.
1232 !! @return The integer value, using the input rounding mode.
1233 function safe_rtoi(rval,mode) result(ival)
1234 real,intent(in) :: rval !< A floating point value.
1235 integer,intent(in) :: mode !< A rouding mode (either "do_floor" or
1239 big = real(huge(ival))
1240 if (rval .le. big .and. -1.*rval .ge. -1.*big) then
1241 if (mode .eq. do_floor) then
1243 elseif (mode .eq. do_nearest) then
1246 call error_mesg("safe_rtoi","mode must be either do_floor" &
1247 //" or do_nearest.",FATAL)
1250 call error_mesg("safe_rtoi","input value cannot be safely" &
1251 //" converted to a 32-bit integer.",FATAL)
1253 end function safe_rtoi
1255 !-------------------------------------------------------------------------
1256 ! <FUNCTION NAME="time_scalar_divide; operator(/)">
1259 ! Returns the largest time, t, for which n * t <= time.
1262 ! Returns the largest time, t, for which n * t <= time.
1265 ! time_scalar_divide(time, n)
1268 ! <IN NAME="time" UNITS="" TYPE="time_type" DIM="">
1271 ! <IN NAME="n" UNITS="" TYPE="integer" DIM="">
1272 ! An integer factor.
1274 ! <OUT NAME="" UNITS="" TYPE="integer" DIM="double precision" DEFAULT="">
1275 ! Returns the largest time, t, for which n * t <= time.
1278 !> Returns the largest time, t, for which n * t <= time.
1279 function time_scalar_divide(time, n)
1281 ! Returns the largest time, t, for which n * t <= time
1283 type(time_type) :: time_scalar_divide
1284 type(time_type), intent(in) :: time
1285 integer, intent(in) :: n
1286 double precision :: d, div, dseconds_per_day, dticks_per_second
1287 integer :: days, seconds, ticks
1288 type(time_type) :: prod1, prod2
1289 character(len=128) tmp1,tmp2
1292 ! Convert time interval to floating point days; risky for general performance?
1293 dseconds_per_day = dble(seconds_per_day)
1294 dticks_per_second = dble(ticks_per_second)
1295 d = time%days*dseconds_per_day*dticks_per_second + dble(time%seconds)*dticks_per_second + dble(time%ticks)
1298 days = int(div/(dseconds_per_day*dticks_per_second))
1299 seconds = int(div/dticks_per_second - days*dseconds_per_day)
1300 ticks = int(div - (days*dseconds_per_day + dble(seconds))*dticks_per_second)
1301 time_scalar_divide = set_time(seconds, days, ticks)
1303 ! Need to make sure that roundoff isn't killing this
1304 prod1 = n * time_scalar_divide
1305 prod2 = n * (increment_time(time_scalar_divide, days=0, seconds=0, ticks=1))
1306 if(prod1 > time .or. prod2 <= time) then
1307 call get_time(time, seconds, days, ticks)
1308 write(tmp1,20) days,seconds,ticks
1309 call get_time(time_scalar_divide, seconds, days, ticks)
1310 write(tmp2,30) n,days,seconds,ticks
1311 ltmp = error_handler('time_scalar_divide',' quotient error:'//trim(tmp1)//trim(tmp2))
1312 20 format('time=',i7,' days, ',i6,' seconds, ',i6,' ticks')
1313 30 format(' time divided by',i6,'=',i7,' days, ',i6,' seconds, ',i6,' ticks')
1316 end function time_scalar_divide
1319 !-------------------------------------------------------------------------
1320 ! <FUNCTION NAME="interval_alarm">
1323 ! Given a time, and a time interval, this function returns true
1324 ! if this is the closest time step to the alarm time.
1327 ! This is a specialized operation that is frequently performed in models.
1328 ! Given a time, and a time interval, this function is true if this is the
1329 ! closest time step to the alarm time. The actual computation is:
1331 ! if((alarm_time - time) <= (time_interval / 2))
1333 ! If the function is true, the alarm time is incremented by the
1334 ! alarm_interval; WARNING, this is a featured side effect. Otherwise, the
1335 ! function is false and there are no other effects. CAUTION: if the
1336 ! alarm_interval is smaller than the time_interval, the alarm may fail to
1337 ! return true ever again. Watch
1338 ! for problems if the new alarm time is less than time + time_interval
1341 ! interval_alarm(time, time_interval, alarm, alarm_interval)
1344 ! <IN NAME="time" TYPE="time_type"> Current time. </IN>
1345 ! <IN NAME="time_interval" TYPE="time_type"> A time interval. </IN>
1346 ! <IN NAME="alarm_interval" TYPE="time_type"> A time interval. </IN>
1347 ! <OUT NAME="interval_alarm" TYPE="logical">
1348 ! Returns either True or false.
1350 ! <INOUT NAME="alarm" TYPE="time_type">
1351 ! An alarm time, which is incremented by the alarm_interval
1352 ! if the function is true.
1355 !> Supports a commonly used type of test on times for models. Given the
1356 !! current time, and a time for an alarm, determines if this is the closest
1357 !! time to the alarm time given a time step of time_interval. If this
1358 !! is the closest time (alarm - time <= time_interval/2), the function
1359 !! returns true and the alarm is incremented by the alarm_interval. Watch
1360 !! for problems if the new alarm time is less than time + time_interval
1361 function interval_alarm(time, time_interval, alarm, alarm_interval)
1363 logical :: interval_alarm
1364 type(time_type), intent(in) :: time, time_interval, alarm_interval
1365 type(time_type), intent(inout) :: alarm
1367 if((alarm - time) <= (time_interval / 2)) then
1368 interval_alarm = .TRUE.
1369 alarm = alarm + alarm_interval
1371 interval_alarm = .FALSE.
1374 end function interval_alarm
1377 !--------------------------------------------------------------------------
1378 ! <FUNCTION NAME="repeat_alarm">
1381 ! Repeat_alarm supports an alarm that goes off with
1382 ! alarm_frequency and lasts for alarm_length.
1385 ! Repeat_alarm supports an alarm that goes off with alarm_frequency and
1386 ! lasts for alarm_length. If the nearest occurence of an alarm time
1387 ! is less than half an alarm_length from the input time, repeat_alarm
1388 ! is true. For instance, if the alarm_frequency is 1 day, and the
1389 ! alarm_length is 2 hours, then repeat_alarm is true from time 2300 on
1390 ! day n to time 0100 on day n + 1 for all n.
1393 ! repeat_alarm(time, alarm_frequency, alarm_length)
1396 ! <IN NAME="time" TYPE="time_type"> Current time. </IN>
1397 ! <IN NAME="alarm_frequency" TYPE="time_type">
1398 ! A time interval for alarm_frequency.
1400 ! <IN NAME="alarm_length" TYPE="time_type">
1401 ! A time interval for alarm_length.
1403 ! <OUT NAME="repeat_alarm" TYPE="logical">
1404 ! Returns either True or false.
1407 !> Repeat_alarm supports an alarm that goes off with alarm_frequency and
1408 !! lasts for alarm_length. If the nearest occurence of an alarm time
1409 !! is less than half an alarm_length from the input time, repeat_alarm
1410 !! is true. For instance, if the alarm_frequency is 1 day, and the
1411 !! alarm_length is 2 hours, then repeat_alarm is true from time 2300 on
1412 !! day n to time 0100 on day n + 1 for all n.
1413 function repeat_alarm(time, alarm_frequency, alarm_length)
1415 logical :: repeat_alarm
1416 type(time_type), intent(in) :: time, alarm_frequency, alarm_length
1417 type(time_type) :: prev, next
1419 prev = (time / alarm_frequency) * alarm_frequency
1420 next = prev + alarm_frequency
1421 if(time - prev <= alarm_length / 2 .or. next - time <= alarm_length / 2) then
1422 repeat_alarm = .TRUE.
1424 repeat_alarm = .FALSE.
1427 end function repeat_alarm
1430 !--------------------------------------------------------------------------
1432 !=========================================================================
1433 ! CALENDAR OPERATIONS BEGIN HERE
1434 !=========================================================================
1436 ! <SUBROUTINE NAME="set_calendar_type">
1439 ! Sets the default calendar type for mapping time intervals to dates.
1442 ! A constant number for setting the calendar type.
1444 ! <TEMPLATE> set_calendar_type(type, err_msg) </TEMPLATE>
1446 ! <IN NAME="type" TYPE="integer" DIM="(scalar)" DEFAULT="NO_CALENDAR">
1447 ! A constant number for setting the calendar type.
1449 ! <OUT NAME="err_msg" TYPE="character, optional" DIM="(scalar)">
1450 ! When present, and when non-blank, a fatal error condition as been detected.
1451 ! The string itself is an error message.
1452 ! It is recommended that, when err_msg is present in the call
1453 ! to this routine, the next line of code should be something
1455 ! if(err_msg /= '') call error_mesg('my_routine','additional info: '//trim(err_msg),FATAL)
1458 !> @brief Sets calendar_type.
1459 !! For the Gregorian calendar, negative years and the proleptic calendar are not used;
1460 !! and the discontinuity of days in October 1582 (when the Gregorian calendar was adopted by select groups in Europe)
1461 !! is also not taken into account.
1462 subroutine set_calendar_type(type, err_msg)
1464 ! Selects calendar for default mapping from time to date.
1466 integer, intent(in) :: type
1467 character(len=*), intent(out), optional :: err_msg
1468 character(len=256) :: err_msg_local
1470 if(.not.module_is_initialized) call time_manager_init()
1472 if(present(err_msg)) err_msg = ''
1474 if(type < 0 .or. type > max_type) then
1475 err_msg_local = 'Illegal calendar type'
1476 if(error_handler('subroutine set_calendar_type', err_msg_local, err_msg)) return
1479 if(seconds_per_day /= 86400 .and. type /= NO_CALENDAR ) then
1480 err_msg_local = 'Only calendar type NO_CALENDAR is allowed when seconds_per_day is not 86400.'// &
1481 ' You are using '//trim(valid_calendar_types(type))//' and seconds_per_day='
1482 write(err_msg_local(len_trim(err_msg_local)+1:len_trim(err_msg_local)+8),'(i8)') seconds_per_day
1483 if(error_handler('subroutine set_calendar_type', err_msg_local, err_msg)) return
1486 calendar_type = type
1488 end subroutine set_calendar_type
1491 !------------------------------------------------------------------------
1492 ! <FUNCTION NAME="get_calendar_type">
1495 ! Returns the value of the default calendar type for mapping
1496 ! from time to date.
1499 ! There are no arguments in this function. It returns the value of
1500 ! the default calendar type for mapping from time to date.
1503 ! get_calendar_type()
1506 !> Returns default calendar type for mapping from time to date.
1507 function get_calendar_type()
1509 integer :: get_calendar_type
1511 get_calendar_type = calendar_type
1513 end function get_calendar_type
1516 !------------------------------------------------------------------------
1517 ! <SUBROUTINE NAME="set_ticks_per_second">
1520 ! Sets the number of ticks per second.
1523 ! Sets the number of ticks per second.
1525 ! <TEMPLATE> call set_ticks_per_second(ticks_per_second) </TEMPLATE>
1526 ! <IN NAME="type" TYPE="integer" DIM="(scalar)" DEFAULT="1"> </IN>
1528 !> Sets the number of ticks per second.
1529 subroutine set_ticks_per_second(tps)
1530 integer, intent(in) :: tps
1532 ticks_per_second = tps
1534 end subroutine set_ticks_per_second
1538 !------------------------------------------------------------------------
1539 ! <FUNCTION NAME="get_ticks_per_second">
1542 ! Returns the number of ticks per second.
1545 ! Returns the number of ticks per second.
1548 ! ticks_per_second = get_ticks_per_second()
1551 !> Returns the number of ticks per second.
1552 function get_ticks_per_second()
1553 integer :: get_ticks_per_second
1555 get_ticks_per_second = ticks_per_second
1557 end function get_ticks_per_second
1560 !------------------------------------------------------------------------
1562 !========================================================================
1563 ! START OF get_date BLOCK
1564 ! <SUBROUTINE NAME="get_date">
1567 ! Given a time_interval, returns the corresponding date under
1568 ! the selected calendar.
1571 ! Given a time_interval, returns the corresponding date under
1572 ! the selected calendar.
1575 ! get_date(time, year, month, day, hour, minute, second, tick, err_msg)
1577 ! <IN NAME="time" TYPE="time_type"> A time interval.</IN>
1578 ! <OUT NAME="year" TYPE="integer"></OUT>
1579 ! <OUT NAME="month" TYPE="integer"></OUT>
1580 ! <OUT NAME="day" TYPE="integer"></OUT>
1581 ! <OUT NAME="hour" TYPE="integer"></OUT>
1582 ! <OUT NAME="minute" TYPE="integer"></OUT>
1583 ! <OUT NAME="second" TYPE="integer"></OUT>
1584 ! <OUT NAME="tick" TYPE="integer, optional"></OUT>
1585 ! <OUT NAME="err_msg" TYPE="character, optional" DIM="(scalar)">
1586 ! When present, and when non-blank, a fatal error condition as been detected.
1587 ! The string itself is an error message.
1588 ! It is recommended that, when err_msg is present in the call
1589 ! to this routine, the next line of code should be something
1591 ! if(err_msg /= '') call error_mesg('my_routine','additional info: '//trim(err_msg),FATAL)
1594 !> @brief Gets the date for different calendar types.
1595 subroutine get_date(time, year, month, day, hour, minute, second, tick, err_msg)
1597 ! Given a time, computes the corresponding date given the selected calendar
1599 type(time_type), intent(in) :: time
1600 integer, intent(out) :: second, minute, hour, day, month, year
1601 integer, intent(out), optional :: tick
1602 character(len=*), intent(out), optional :: err_msg
1603 character(len=128) :: err_msg_local
1606 if(.not.module_is_initialized) call time_manager_init
1607 if(present(err_msg)) err_msg = ''
1609 select case(calendar_type)
1610 case(THIRTY_DAY_MONTHS)
1611 call get_date_thirty(time, year, month, day, hour, minute, second, tick1)
1613 call get_date_gregorian(time, year, month, day, hour, minute, second, tick1)
1615 call get_date_julian_private(time, year, month, day, hour, minute, second, tick1)
1617 call get_date_no_leap_private(time, year, month, day, hour, minute, second, tick1)
1619 err_msg_local = 'Cannot produce a date when the calendar type is NO_CALENDAR'
1620 if(error_handler('subroutine get_date', err_msg_local, err_msg)) return
1622 err_msg_local = 'Invalid calendar type'
1623 if(error_handler('subroutine get_date', err_msg_local, err_msg)) return
1626 if(present(tick)) then
1630 err_msg_local = 'tick must be present when time has a second fraction'
1631 if(error_handler('subroutine get_date', err_msg_local, err_msg)) return
1635 end subroutine get_date
1637 !------------------------------------------------------------------------
1639 !> @brief Gets the date on a Gregorian calendar.
1640 !! Computes the year, month, day on the fly from the quantity time%days
1641 subroutine get_date_gregorian(time, year, month, day, hour, minute, second, tick)
1643 type(time_type), intent(in) :: time
1644 integer, intent(out) :: year, month, day, hour, minute, second
1645 integer, intent(out) :: tick
1646 integer :: iday, isec
1648 integer :: l !< value of 1 if leap year; value of 0 otherwise
1649 integer :: ncenturies !< number of centuries passed in the current 400 year period
1650 integer :: nlpyrs !< number of leap years in the current century
1651 integer :: yearx, monthx, dayx, idayx !< temporary values for year, month, day
1652 integer :: i !< counter, dummy variable
1654 ! Computes date corresponding to time for gregorian calendar
1656 !Carried over from the old subroutine
1657 if(Time%seconds >= 86400) then ! This check appears to be unecessary.
1658 call error_mesg('get_date','Time%seconds .ge. 86400 in subroutine get_date_gregorian',FATAL)
1661 iday = mod(Time%days+1,days_in_400_year_period)
1665 if( iday.eq.0 ) then !year 400
1668 else if( iday.gt.365 ) then
1669 yearx = int(iday/365) - 1 !approximation off by -1 year at most
1670 ncenturies = int(yearx/100)
1671 nlpyrs = int((yearx-ncenturies*100)/4)
1672 idayx = ncenturies*36524 + (yearx-ncenturies*100)*365 + nlpyrs !36524 days in a century
1673 if( ncenturies.eq.4 ) idayx = idayx + 1 !year 400 is a leap year
1674 l = 0 ; if ( leap_year_gregorian_int(yearx+1) ) l = 1
1675 if ( (iday-idayx).gt.365+l ) then !off by -1 year
1677 idayx = idayx + 365 + l
1682 year = 400*int((Time%days+1)/days_in_400_year_period) + yearx
1684 l = 0 ; if( leap_year_gregorian_int(year) ) l = 1
1686 if( dayx.le.31 ) then
1690 monthx = int(dayx/30)
1693 dayx = dayx - days_per_month(i)
1694 if(i.eq.2) dayx = dayx - l
1698 dayx = dayx - days_per_month(i)
1703 if( dayx.le.0 ) then
1705 day = dayx + days_per_month(monthx)
1706 if(monthx.eq.2) day = day + l
1710 hour = Time%seconds / 3600
1711 isec = Time%seconds - 3600*hour
1713 second = isec - 60*minute
1716 end subroutine get_date_gregorian
1717 !------------------------------------------------------------------------
1719 function cut0(string)
1720 character(len=256) :: cut0
1721 character(len=*), intent(in) :: string
1727 if(ichar(string(i:i)) == 0 ) then
1734 !------------------------------------------------------------------------
1736 !> Base date for Julian calendar is year 1 with all multiples of 4
1737 !! years being leap years.
1738 subroutine get_date_julian_private(time, year, month, day, hour, minute, second, tick)
1741 type(time_type), intent(in) :: time
1742 integer, intent(out) :: second, minute, hour, day, month, year
1743 integer, intent(out) :: tick
1744 integer :: m, t, nfour, nex, days_this_month
1747 ! find number of four year periods; also get modulo number of days
1748 nfour = time%days / (4 * 365 + 1)
1749 day = modulo(time%days, (4 * 365 + 1))
1751 ! Find out what year in four year chunk
1757 day=modulo(day, 365) + 1
1760 ! Is this a leap year?
1763 year = 1 + 4 * nfour + nex
1765 ! find month and day
1768 days_this_month = days_per_month(m)
1769 if(leap .and. m == 2) days_this_month = 29
1770 if(day <= days_this_month) exit
1771 day = day - days_this_month
1774 ! find hour,minute and second
1776 hour = t / (60 * 60)
1777 t = t - hour * (60 * 60)
1779 second = t - 60 * minute
1781 end subroutine get_date_julian_private
1783 !------------------------------------------------------------------------
1784 subroutine get_date_julian(time, year, month, day, hour, minute, second)
1786 ! No need to include tick in argument list because this routine
1787 ! exists only for interpolator.F90, which does not need it.
1789 type(time_type), intent(in) :: time
1790 integer, intent(out) :: second, minute, hour, day, month, year
1793 call get_date_julian_private(time, year, month, day, hour, minute, second, tick)
1795 end subroutine get_date_julian
1797 !------------------------------------------------------------------------
1799 subroutine get_date_thirty(time, year, month, day, hour, minute, second, tick)
1801 ! Computes date corresponding to time interval for 30 day months, 12
1804 type(time_type), intent(in) :: time
1805 integer, intent(out) :: second, minute, hour, day, month, year
1806 integer, intent(out) :: tick
1807 integer :: t, dmonth, dyear
1810 dyear = t / (30 * 12)
1812 t = t - dyear * (30 * 12)
1815 day = t -dmonth * 30 + 1
1818 hour = t / (60 * 60)
1819 t = t - hour * (60 * 60)
1821 second = t - 60 * minute
1824 end subroutine get_date_thirty
1825 !------------------------------------------------------------------------
1827 subroutine get_date_no_leap_private(time, year, month, day, hour, minute, second, tick)
1829 ! Base date for NOLEAP calendar is year 1.
1831 type(time_type), intent(in) :: time
1832 integer, intent(out) :: second, minute, hour, day, month, year
1833 integer, intent(out) :: tick
1836 ! get modulo number of days
1837 year = time%days / 365 + 1
1838 day = modulo(time%days, 365) + 1
1840 ! find month and day
1843 if(day <= days_per_month(m)) exit
1844 day = day - days_per_month(m)
1847 ! find hour,minute and second
1849 hour = t / (60 * 60)
1850 t = t - hour * (60 * 60)
1852 second = t - 60 * minute
1855 end subroutine get_date_no_leap_private
1857 !------------------------------------------------------------------------
1858 subroutine get_date_no_leap(time, year, month, day, hour, minute, second)
1860 ! No need to include tick in argument list because this routine
1861 ! exists only for interpolator.F90, which does not need it.
1863 type(time_type), intent(in) :: time
1864 integer, intent(out) :: second, minute, hour, day, month, year
1867 call get_date_no_leap_private(time, year, month, day, hour, minute, second, tick)
1869 end subroutine get_date_no_leap
1870 !------------------------------------------------------------------------
1872 ! END OF get_date BLOCK
1873 !========================================================================
1874 ! START OF set_date BLOCK
1875 ! <FUNCTION NAME="set_date">
1878 ! Given an input date in year, month, days, etc., creates a
1879 ! time_type that represents this time interval from the
1880 ! internally defined base date.
1883 ! Given a date, computes the corresponding time given the selected
1884 ! date time mapping algorithm. Note that it is possible to specify
1885 ! any number of illegal dates; these should be checked for and generate
1886 ! errors as appropriate.
1889 ! 1. set_date(year, month, day, hours, minute, second, tick, err_msg)
1892 ! 2. set_date_c(time_string, zero_year_warning, err_msg, allow_rounding)
1893 ! time_string is a character string containing a date formatted
1894 ! according to CF conventions. e.g. '1980-12-31 23:59:59.9'
1896 ! <IN NAME="time" TYPE="time_type"> A time interval.</IN>
1897 ! <IN NAME="year" TYPE="integer"></IN>
1898 ! <IN NAME="month" TYPE="integer"></IN>
1899 ! <IN NAME="day" TYPE="integer"></IN>
1900 ! <IN NAME="hour" TYPE="integer"></IN>
1901 ! <IN NAME="minute" TYPE="integer"></IN>
1902 ! <IN NAME="second" TYPE="integer"></IN>
1903 ! <IN NAME="tick" TYPE="integer"></IN>
1904 ! <IN NAME="zero_year_warning" TYPE="logical">
1905 ! If the year number is zero, it will be silently changed to one,
1906 ! unless zero_year_warning=.true., in which case a WARNING message
1907 ! will also be issued.
1909 ! <IN NAME="allow_rounding" TYPE="logical, optional" DEFAULT=".true.">
1910 ! When .true., any fractions of a second will be rounded off to the nearest tick.
1911 ! When .false., it is a fatal error if the second fraction cannot be exactly
1912 ! represented by a number of ticks.
1914 ! <OUT NAME="err_msg" TYPE="character, optional" DIM="(scalar)">
1915 ! When present, and when non-blank, a fatal error condition as been detected.
1916 ! The string itself is an error message.
1917 ! It is recommended that, when err_msg is present in the call
1918 ! to this routine, the next line of code should be something
1920 ! if(err_msg /= '') call error_mesg('my_routine','additional info: '//trim(err_msg),FATAL)
1922 ! <OUT NAME="set_date" TYPE="time_type"> A time interval.</OUT>
1924 !> @brief Sets days for different calendar types.
1925 function set_date_private(year, month, day, hour, minute, second, tick, Time_out, err_msg)
1927 ! Given a date, computes the corresponding time given the selected
1928 ! date time mapping algorithm. Note that it is possible to specify
1929 ! any number of illegal dates; these are checked for and generate
1930 ! errors as appropriate.
1932 logical :: set_date_private
1933 integer, intent(in) :: year, month, day, hour, minute, second, tick
1934 type(time_type) :: Time_out
1935 character(len=*), intent(out) :: err_msg
1937 if(.not.module_is_initialized) call time_manager_init
1941 select case(calendar_type)
1942 case(THIRTY_DAY_MONTHS)
1943 set_date_private = set_date_thirty (year, month, day, hour, minute, second, tick, Time_out, err_msg)
1945 set_date_private = set_date_gregorian (year, month, day, hour, minute, second, tick, Time_out, err_msg)
1947 set_date_private = set_date_julian_private (year, month, day, hour, minute, second, tick, Time_out, err_msg)
1949 set_date_private = set_date_no_leap_private(year, month, day, hour, minute, second, tick, Time_out, err_msg)
1951 err_msg = 'Cannot produce a date when calendar type is NO_CALENDAR'
1952 set_date_private = .false.
1954 err_msg = 'Invalid calendar type'
1955 set_date_private = .false.
1958 end function set_date_private
1961 !------------------------------------------------------------------------
1963 !> @brief Calls set_date_private to set days for different calendar types.
1964 function set_date_i(year, month, day, hour, minute, second, tick, err_msg)
1965 type(time_type) :: set_date_i
1966 integer, intent(in) :: day, month, year
1967 integer, intent(in), optional :: second, minute, hour, tick
1968 character(len=*), intent(out), optional :: err_msg
1969 integer :: osecond, ominute, ohour, otick
1970 character(len=128) :: err_msg_local
1972 if(.not.module_is_initialized) call time_manager_init
1973 if(present(err_msg)) err_msg = ''
1975 ! Missing optionals are set to 0
1976 osecond = 0; if(present(second)) osecond = second
1977 ominute = 0; if(present(minute)) ominute = minute
1978 ohour = 0; if(present(hour)) ohour = hour
1979 otick = 0; if(present(tick)) otick = tick
1981 if(.not.set_date_private(year, month, day, ohour, ominute, osecond, otick, set_date_i, err_msg_local)) then
1982 if(error_handler('function set_date_i', err_msg_local, err_msg)) return
1985 end function set_date_i
1986 !------------------------------------------------------------------------
1988 !> @brief Calls set_date_private for different calendar types when given a string input.
1989 function set_date_c(string, zero_year_warning, err_msg, allow_rounding)
1991 ! Examples of acceptable forms of string:
1993 ! 1980-01-01 00:00:00
1994 ! 1980-01-01 00:00:00.50
1998 ! year number must occupy 4 spaces.
1999 ! months, days, hours, minutes, seconds may occupy 1 or 2 spaces
2000 ! year, month and day must be separated by a '-'
2001 ! hour, minute, second must be separated by a ':'
2002 ! hour, minute, second are optional. If not present then zero is assumed.
2003 ! second may be a real number.
2005 ! zero_year_warning:
2006 ! If the year number is zero, it will be silently changed to one,
2007 ! unless zero_year_warning=.true., in which case a WARNING message
2008 ! will also be issued
2010 type(time_type) :: set_date_c
2011 character(len=*), intent(in) :: string
2012 logical, intent(in), optional :: zero_year_warning
2013 character(len=*), intent(out), optional :: err_msg
2014 logical, intent(in), optional :: allow_rounding
2015 character(len=4) :: formt='(i )'
2016 logical :: correct_form, zero_year_warning_local, allow_rounding_local
2017 integer :: i1, i2, i3, i4, i5, i6, i7
2018 character(len=32) :: string_sifted_left
2019 integer :: year, month, day, hour, minute, second, tick
2020 character(len=128) :: err_msg_local
2022 if(.not.module_is_initialized) call time_manager_init()
2023 if(present(err_msg)) err_msg = ''
2024 if(present(zero_year_warning)) then
2025 zero_year_warning_local = zero_year_warning
2027 zero_year_warning_local = .true.
2029 if(present(allow_rounding)) then
2030 allow_rounding_local = allow_rounding
2032 allow_rounding_local = .true.
2035 string_sifted_left = adjustl(string)
2036 i1 = index(string_sifted_left,'-')
2037 i2 = index(string_sifted_left,'-',back=.true.)
2038 i3 = index(string_sifted_left,':')
2039 i4 = index(string_sifted_left,':',back=.true.)
2040 i5 = len_trim(cut0(string_sifted_left))
2041 i6 = index(string_sifted_left,'.',back=.true.)
2042 correct_form = (i1 > 1) ! year number must occupy at least 1 space
2043 correct_form = correct_form .and. (i2-i1 == 2 .or. i2-i1 == 3) ! month number must occupy 1 or 2 spaces
2044 if(.not.correct_form) then
2045 err_msg_local = 'Form of character time stamp is incorrect. The character time stamp is: '//trim(string)
2046 if(error_handler('function set_date_c', err_msg_local, err_msg)) return
2048 write(formt(3:3),'(i1)') i1-1
2049 read(string_sifted_left(1:i1-1),formt) year
2052 if(zero_year_warning_local) then
2053 call error_mesg('set_date_c','Year zero is invalid. Resetting year to 1', WARNING)
2056 write(formt(3:3),'(i1)') i2-i1-1
2057 read(string_sifted_left(i1+1:i2-1),formt) month
2059 read(string_sifted_left(i2+1:i7),'(i2)') day
2062 ! There are no minutes or seconds in the string
2067 ! There is no clocktime in the string at all
2070 ! The clocktime includes only hours
2071 read(string_sifted_left(i5-1:i5),'(i2)') hour
2073 else if(i3 == i4) then
2074 ! The string includes hours and minutes, but no seconds
2075 read(string_sifted_left(i3-2:i3-1),'(i2)') hour
2076 write(formt(3:3),'(i1)') i5-i3
2077 read(string_sifted_left(i3+1:i5),formt) minute
2081 ! The string includes hours, minutes, and seconds
2082 read(string_sifted_left(i3-2:i3-1),'(i2)') hour
2083 write(formt(3:3),'(i1)') i4-i3-1
2084 read(string_sifted_left(i3+1:i4-1),formt) minute
2085 write(formt(3:3),'(i1)') i5-i4
2087 ! There are no fractional seconds
2088 read(string_sifted_left(i4+1:i5),formt) second
2091 read(string_sifted_left(i4+1:i6-1),formt) second
2092 if(.not.get_tick_from_string(string_sifted_left(i6:i5), err_msg_local, allow_rounding_local, tick)) then
2093 if(error_handler('function set_date_c', err_msg_local, err_msg)) return
2095 ! If tick has been rounded up to ticks_per_second, then bump up second.
2096 if(tick == ticks_per_second) then
2103 if(.not.set_date_private(year, month, day, hour, minute, second, tick, set_date_c, err_msg_local)) then
2104 if(error_handler('function set_date_c', err_msg_local, err_msg)) return
2107 end function set_date_c
2109 !------------------------------------------------------------------------
2111 !> @brief Sets Time_out%days on a Gregorian calendar
2112 !! Computes the total number of days between 1/1/0001 to the current month/day/year
2113 function set_date_gregorian(year, month, day, hour, minute, second, tick, Time_out, err_msg)
2114 logical :: set_date_gregorian
2116 ! Computes time corresponding to date for gregorian calendar.
2118 integer, intent(in) :: year, month, day, hour, minute, second, tick
2119 type(time_type), intent(out) :: Time_out
2120 character(len=*), intent(out) :: err_msg
2121 integer :: yearx, monthx, dayx, hrx, minx, secx, tickx, ncenturies, nlpyrs, l
2123 if( .not.valid_increments(year,month,day,hour,minute,second,tick,err_msg) ) then
2124 set_date_gregorian = .false.
2128 l = 0 ; if( leap_year_gregorian_int(year) ) l = 1
2130 ! Check if date is invalid
2132 if(day.gt.days_per_month(month)+l .or. day.lt.1) then
2133 err_msg = 'Invalid_date. Date='//convert_integer_date_to_char(year,month,day,hour,minute,second)
2134 set_date_gregorian = .false.
2138 if(day.gt.days_per_month(month) .or. day.lt.1) then
2139 err_msg = 'Invalid_date. Date='//convert_integer_date_to_char(year,month,day,hour,minute,second)
2140 set_date_gregorian = .false.
2145 Time_out%seconds = second + 60*(minute + 60*hour)
2147 yearx = mod(year-1,400)
2150 ncenturies = int( yearx/100 )
2151 nlpyrs = int( (yearx-ncenturies*100)/4 )
2152 dayx = ncenturies*36524 + (yearx-ncenturies*100)*365 + nlpyrs ! 36524 days in 100 years, year 100 not
2154 if(ncenturies.eq.4) dayx = dayx + 1 ! year 400 is a leap year
2157 select case( month )
2158 case(1) ; dayx = dayx
2159 case(2) ; dayx = dayx + 31
2160 case(3) ; dayx = dayx + 59 + l
2161 case(4) ; dayx = dayx + 90 + l
2162 case(5) ; dayx = dayx + 120 + l
2163 case(6) ; dayx = dayx + 151 + l
2164 case(7) ; dayx = dayx + 181 + l
2165 case(8) ; dayx = dayx + 212 + l
2166 case(9) ; dayx = dayx + 243 + l
2167 case(10) ; dayx = dayx + 273 + l
2168 case(11) ; dayx = dayx + 304 + l
2169 case(12) ; dayx = dayx + 334 + l
2172 dayx = int((year-1)/400)*days_in_400_year_period + dayx + day - 1
2173 Time_out%days = dayx
2174 Time_out%ticks = tick
2177 set_date_gregorian = .true.
2180 yearx = year ; monthx = month ; dayx = day
2181 hrx = hour ; minx = minute ; secx = second ; tickx = tick
2182 call get_date_gregorian(Time_out, yearx, monthx, dayx, hrx, minx, secx, tickx)
2183 l = 0 ; if( leap_year_gregorian_int(yearx) ) l = 1
2184 if( monthx.lt.1 .or. monthx.gt.12 ) then
2185 err_msg = 'Invalid_date. Date='//convert_integer_date_to_char(yearx,monthx,dayx,hour,minute,second)
2186 set_date_gregorian = .false.
2187 else if( dayx.lt.1 .or. dayx.gt.days_per_month(monthx) ) then
2188 if( monthx.eq.2 .and. dayx.le.days_per_month(monthx)+l ) return
2189 err_msg = 'Invalid_date. Date='//convert_integer_date_to_char(yearx,monthx,dayx,hour,minute,second)
2190 set_date_gregorian = .false.
2193 end function set_date_gregorian
2195 !------------------------------------------------------------------------
2197 function set_date_julian_private(year, month, day, hour, minute, second, tick, Time_out, err_msg)
2198 logical :: set_date_julian_private
2200 ! Returns time corresponding to date for julian calendar.
2202 integer, intent(in) :: year, month, day, hour, minute, second, tick
2203 type(time_type), intent(out) :: Time_out
2204 character(len=*), intent(out) :: err_msg
2205 integer :: ndays, m, nleapyr
2208 if( .not.valid_increments(year,month,day,hour,minute,second,tick,err_msg) ) then
2209 set_date_julian_private = .false.
2213 if(month /= 2 .and. day > days_per_month(month)) then
2214 err_msg = 'Invalid date. Date='//convert_integer_date_to_char(year,month,day,hour,minute,second)
2215 set_date_julian_private = .false.
2219 ! Is this a leap year?
2220 leap = (modulo(year,4) == 0)
2221 ! compute number of complete leap years from year 1
2222 nleapyr = (year - 1) / 4
2224 ! Finish checking for day specication errors
2225 if(month == 2 .and. (day > 29 .or. ((.not. leap) .and. day > 28))) then
2226 err_msg = 'Invalid date. Date='//convert_integer_date_to_char(year,month,day,hour,minute,second)
2227 set_date_julian_private = .false.
2233 ndays = ndays + days_per_month(m)
2234 if(leap .and. m == 2) ndays = ndays + 1
2237 Time_out%seconds = second + 60 * (minute + 60 * hour)
2238 Time_out%days = day -1 + ndays + 365*(year - nleapyr - 1) + 366*(nleapyr)
2239 Time_out%ticks = tick
2241 set_date_julian_private = .true.
2243 end function set_date_julian_private
2245 !------------------------------------------------------------------------
2246 function set_date_julian(year, month, day, hour, minute, second)
2248 ! No need to include tick or err_msg in argument list because this
2249 ! routine exists only for interpolator.F90, which does not need them.
2251 type(time_type) :: set_date_julian
2252 integer, intent(in) :: year, month, day, hour, minute, second
2253 character(len=128) :: err_msg
2255 if(.not.set_date_julian_private(year, month, day, hour, minute, second, 0, set_date_julian, err_msg)) then
2256 call error_mesg('set_date_julian',trim(err_msg),FATAL)
2259 end function set_date_julian
2260 !------------------------------------------------------------------------
2262 function set_date_thirty(year, month, day, hour, minute, second, tick, Time_out, err_msg)
2263 logical :: set_date_thirty
2265 ! Computes time corresponding to date for thirty day months.
2267 integer, intent(in) :: year, month, day, hour, minute, second, tick
2268 type(time_type), intent(out) :: Time_out
2269 character(len=*), intent(out) :: err_msg
2271 if( .not.valid_increments(year,month,day,hour,minute,second,tick,err_msg) ) then
2272 set_date_thirty = .false.
2277 err_msg = 'Invalid date. Date='//convert_integer_date_to_char(year,month,day,hour,minute,second)
2278 set_date_thirty = .false.
2282 Time_out%days = (day - 1) + 30 * ((month - 1) + 12 * (year - 1))
2283 Time_out%seconds = second + 60 * (minute + 60 * hour)
2284 Time_out%ticks = tick
2286 set_date_thirty = .true.
2288 end function set_date_thirty
2290 !------------------------------------------------------------------------
2292 function set_date_no_leap_private(year, month, day, hour, minute, second, tick, Time_out, err_msg)
2293 logical :: set_date_no_leap_private
2295 ! Computes time corresponding to date for fixed 365 day year calendar.
2297 integer, intent(in) :: year, month, day, hour, minute, second, tick
2298 type(time_type), intent(out) :: Time_out
2299 character(len=*), intent(out) :: err_msg
2302 if( .not.valid_increments(year,month,day,hour,minute,second,tick,err_msg) ) then
2303 set_date_no_leap_private = .false.
2307 if(day > days_per_month(month)) then
2308 err_msg = 'Invalid date. Date='//convert_integer_date_to_char(year,month,day,hour,minute,second)
2309 set_date_no_leap_private = .false.
2315 ndays = ndays + days_per_month(m)
2318 ! No need for err_msg in call to set_time because previous checks ensure positive value of time.
2319 Time_out = set_time(second + 60 * (minute + 60 * hour), day -1 + ndays + 365 * (year - 1), tick)
2321 set_date_no_leap_private = .true.
2323 end function set_date_no_leap_private
2324 !------------------------------------------------------------------------
2326 function set_date_no_leap(year, month, day, hour, minute, second)
2328 ! No need to include tick or err_msg in argument list because this
2329 ! routine exists only for interpolator.F90, which does not need them.
2331 type(time_type) :: set_date_no_leap
2332 integer, intent(in) :: year, month, day, hour, minute, second
2333 character(len=128) :: err_msg
2335 if(.not.set_date_no_leap_private(year, month, day, hour, minute, second, 0, set_date_no_leap, err_msg)) then
2336 call error_mesg('set_date_no_leap',trim(err_msg),FATAL)
2339 end function set_date_no_leap
2341 !=========================================================================
2343 function valid_increments(year, month, day, hour, minute, second, tick, err_msg)
2344 logical :: valid_increments
2345 integer, intent(in) :: year, month, day, hour, minute, second, tick
2346 character(len=128), intent(out) :: err_msg
2348 ! Check for invalid values
2351 valid_increments = .true.
2352 if(second > 59 .or. second < 0 .or. minute > 59 .or. minute < 0 &
2353 .or. hour > 23 .or. hour < 0 .or. day > 31 .or. day < 1 &
2354 .or. month > 12 .or. month < 1 .or. year < 1) then
2355 err_msg = 'Invalid date. Date='//convert_integer_date_to_char(year,month,day,hour,minute,second)
2356 valid_increments = .false.
2359 if(tick < 0 .or. tick >= ticks_per_second) then
2360 write(err_msg,'(a,i6)') 'Invalid number of ticks. tick=',tick
2361 valid_increments = .false.
2364 end function valid_increments
2366 !=========================================================================
2368 function convert_integer_date_to_char(year, month, day, hour, minute, second)
2369 character(len=19) :: convert_integer_date_to_char
2370 integer, intent(in) :: year, month, day
2371 integer, intent(in) :: hour, minute, second
2373 write(convert_integer_date_to_char,10) year,month,day,hour,minute,second
2374 10 format(i4.4, '-', i2.2, '-', i2.2, ' ', i2.2, ':', i2.2, ':', i2.2)
2376 end function convert_integer_date_to_char
2378 !=========================================================================
2379 ! END OF set_date BLOCK
2380 !=========================================================================
2382 ! <FUNCTION NAME="increment_date">
2385 ! Increments the date represented by a time interval and the
2386 ! default calendar type by a number of seconds, etc.
2389 ! Given a time and some date increment, computes a new time. Depending
2390 ! on the mapping algorithm from date to time, it may be possible to specify
2391 ! undefined increments (i.e. if one increments by 68 days and 3 months in
2392 ! a Julian calendar, it matters which order these operations are done and
2393 ! we don't want to deal with stuff like that, make it an error).
2396 ! increment_date(time, years, months, days, hours, minutes, seconds, ticks, err_msg)
2398 ! <IN NAME="time" TYPE="time_type"> A time interval.</IN>
2399 ! <IN NAME="years" TYPE="integer">An increment of years.</IN>
2400 ! <IN NAME="months" TYPE="integer">An increment of months.</IN>
2401 ! <IN NAME="days" TYPE="integer">An increment of days.</IN>
2402 ! <IN NAME="hours" TYPE="integer">An increment of hours.</IN>
2403 ! <IN NAME="minutes" TYPE="integer">An increment of minutes.</IN>
2404 ! <IN NAME="seconds" TYPE="integer">An increment of seconds.</IN>
2405 ! <IN NAME="ticks" TYPE="integer">An increment of ticks.</IN>
2406 ! <OUT NAME="err_msg" TYPE="character, optional" DIM="(scalar)">
2407 ! When present, and when non-blank, a fatal error condition as been detected.
2408 ! The string itself is an error message.
2409 ! It is recommended that, when err_msg is present in the call
2410 ! to this routine, the next line of code should be something
2412 ! if(err_msg /= '') call error_mesg('my_routine','additional info: '//trim(err_msg),FATAL)
2414 ! <OUT NAME="increment_date" TYPE="time_type"> A new time based on the input
2415 ! time interval and the calendar type.
2417 ! <IN NAME="allow_neg_inc" TYPE="logical, optional" DIM="(scalar)" DEFAULT=".true.">
2418 ! When .false., it is a fatal error if any of the input time increments are negative.
2419 ! This mimics the behavior of lima and earlier revisions.
2422 ! For all but the thirty_day_months calendar, increments to months
2423 ! and years must be made separately from other units because of the
2424 ! non-associative nature of addition.
2425 ! If the result is a negative time (i.e. date before the base date)
2426 ! it is considered a fatal error.
2429 function increment_date(Time, years, months, days, hours, minutes, seconds, ticks, err_msg, allow_neg_inc)
2431 ! Given a time and some date increment, computes a new time. Depending
2432 ! on the mapping algorithm from date to time, it may be possible to specify
2433 ! undefined increments (i.e. if one increments by 68 days and 3 months in
2434 ! a Julian calendar, it matters which order these operations are done and
2435 ! we don't want to deal with stuff like that, make it an error).
2437 ! This routine operates in one of two modes.
2438 ! 1. days, hours, minutes, seconds, ticks are incremented, years and months must be zero or absent arguments.
2439 ! 2. years and/or months are incremented, other time increments must be zero or absent arguments.
2441 type(time_type) :: increment_date
2442 type(time_type), intent(in) :: Time
2443 integer, intent(in), optional :: years, months, days, hours, minutes, seconds, ticks
2444 character(len=*), intent(out), optional :: err_msg
2445 logical, intent(in), optional :: allow_neg_inc
2447 integer :: oyears, omonths, odays, ohours, ominutes, oseconds, oticks
2448 character(len=128) :: err_msg_local
2449 logical :: allow_neg_inc_local
2451 if(.not.module_is_initialized) call time_manager_init
2452 if(present(err_msg)) err_msg = ''
2454 ! Missing optionals are set to 0
2455 oseconds = 0; if(present(seconds)) oseconds = seconds
2456 ominutes = 0; if(present(minutes)) ominutes = minutes
2457 ohours = 0; if(present(hours)) ohours = hours
2458 odays = 0; if(present(days)) odays = days
2459 omonths = 0; if(present(months)) omonths = months
2460 oyears = 0; if(present(years)) oyears = years
2461 oticks = 0; if(present(ticks)) oticks = ticks
2462 allow_neg_inc_local=.true.; if(present(allow_neg_inc)) allow_neg_inc_local=allow_neg_inc
2464 if(.not.allow_neg_inc_local) then
2465 if(oyears < 0 .or. omonths < 0 .or. odays < 0 .or. ohours < 0 .or. ominutes < 0 .or. oseconds < 0 .or. &
2467 write(err_msg_local,10) oyears, omonths, odays, ohours, ominutes, oseconds, oticks
2468 if(error_handler('function increment_time', err_msg_local, err_msg)) return
2471 10 format('One or more time increments are negative: '// &
2472 'years=',i6,' months=',i6,' days=',i6,' hours=',i6,' minutes=',i6,' seconds=',i6,' ticks=',i6)
2474 if(.not.increment_date_private( &
2475 Time, oyears, omonths, odays, ohours, ominutes, oseconds, oticks, increment_date, err_msg_local)) then
2476 if(error_handler('function increment_date', err_msg_local, err_msg)) return
2479 end function increment_date
2483 !=======================================================================
2485 function increment_date_private(Time, years, months, days, hours, minutes, seconds, ticks, Time_out, err_msg)
2487 ! Given a time and some date increment, computes a new time. Depending
2488 ! on the mapping algorithm from date to time, it may be possible to specify
2489 ! undefined increments (i.e. if one increments by 68 days and 3 months in
2490 ! a Julian calendar, it matters which order these operations are done and
2491 ! we don't want to deal with stuff like that, make it an error).
2493 ! This routine operates in one of two modes.
2494 ! 1. days, hours, minutes, seconds, ticks are incremented, years and months must be zero or absent arguments.
2495 ! 2. years and/or months are incremented, other time increments must be zero or absent arguments.
2497 ! Negative increments are always allowed in the private version of this routine.
2499 logical :: increment_date_private
2500 type(time_type), intent(in) :: Time
2501 integer, intent(in) :: years, months, days, hours, minutes, seconds, ticks
2502 type(time_type), intent(out) :: Time_out
2503 character(len=*), intent(out) :: err_msg
2504 integer :: cyear , cmonth , cday , chour , cminute , csecond , ctick
2505 logical :: mode_1, mode_2
2508 increment_date_private = .true.
2510 mode_1 = days /= 0 .or. hours /= 0 .or. minutes /= 0 .or. seconds /= 0 .or. ticks /= 0
2511 mode_2 = years /= 0 .or. months /= 0
2513 if(.not.mode_1 .and. .not.mode_2) then
2514 ! All time increments are zero
2519 if(mode_1 .and. mode_2) then
2520 err_msg = 'years and/or months must not be incremented with other time units'
2521 increment_date_private = .false.
2526 csecond = seconds + 60 * (minutes + 60 * hours)
2527 increment_date_private = increment_time_private(Time, csecond, days, ticks, Time_out, err_msg)
2531 ! Convert Time to a date
2532 select case(calendar_type)
2533 case(THIRTY_DAY_MONTHS)
2534 call get_date_thirty (Time, cyear, cmonth, cday, chour, cminute, csecond, ctick)
2536 call get_date_no_leap_private (Time, cyear, cmonth, cday, chour, cminute, csecond, ctick)
2538 call get_date_julian_private (Time, cyear, cmonth, cday, chour, cminute, csecond, ctick)
2540 call get_date_gregorian(Time, cyear, cmonth, cday, chour, cminute, csecond, ctick)
2542 err_msg = 'Cannot increment a date when the calendar type is NO_CALENDAR'
2543 increment_date_private = .false.
2546 err_msg = 'Invalid calendar type'
2547 increment_date_private = .false.
2551 ! Add month increment
2552 cmonth = cmonth + months
2554 ! Adjust year and month number when cmonth falls outside the range 1 to 12
2555 cyear = cyear + floor((cmonth-1)/12.)
2556 cmonth = modulo((cmonth-1),12) + 1
2558 ! Add year increment
2559 cyear = cyear + years
2561 ! Convert this back into a time.
2562 select case(calendar_type)
2563 case(THIRTY_DAY_MONTHS)
2564 increment_date_private = set_date_thirty (cyear, cmonth, cday, chour, cminute, csecond, ctick, Time_out, err_msg)
2566 increment_date_private = set_date_no_leap_private (cyear, cmonth, cday, chour, cminute, csecond, ctick, &
2567 & Time_out, err_msg)
2569 increment_date_private = set_date_julian_private (cyear, cmonth, cday, chour, cminute, csecond, ctick, &
2570 & Time_out, err_msg)
2572 increment_date_private = set_date_gregorian(cyear, cmonth, cday, chour, cminute, csecond, ctick, Time_out, &
2577 end function increment_date_private
2579 !=========================================================================
2580 ! <FUNCTION NAME="decrement_date">
2583 ! Decrements the date represented by a time interval and the
2584 ! default calendar type by a number of seconds, etc.
2587 ! Given a time and some date decrement, computes a new time. Depending
2588 ! on the mapping algorithm from date to time, it may be possible to specify
2589 ! undefined decrements (i.e. if one decrements by 68 days and 3 months in
2590 ! a Julian calendar, it matters which order these operations are done and
2591 ! we don't want to deal with stuff like that, make it an error).
2594 ! decrement_date(time, years, months, days, hours, minutes, seconds, ticks, err_msg))
2596 ! <IN NAME="time" TYPE="time_type"> A time interval.</IN>
2597 ! <IN NAME="years" TYPE="integer">An decrement of years.</IN>
2598 ! <IN NAME="months" TYPE="integer">An decrement of months.</IN>
2599 ! <IN NAME="days" TYPE="integer">An decrement of days.</IN>
2600 ! <IN NAME="hours" TYPE="integer">An decrement of hours.</IN>
2601 ! <IN NAME="minutes" TYPE="integer">An decrement of minutes.</IN>
2602 ! <IN NAME="seconds" TYPE="integer">An decrement of seconds.</IN>
2603 ! <IN NAME="ticks" TYPE="integer">An decrement of ticks.</IN>
2604 ! <OUT NAME="err_msg" TYPE="character, optional" DIM="(scalar)">
2605 ! When present, and when non-blank, a fatal error condition as been detected.
2606 ! The string itself is an error message.
2607 ! It is recommended that, when err_msg is present in the call
2608 ! to this routine, the next line of code should be something
2610 ! if(err_msg /= '') call error_mesg('my_routine','additional info: '//trim(err_msg),FATAL)
2612 ! <OUT NAME="decrement_date" TYPE="time_type"> A new time based on the input
2613 ! time interval and the calendar type.
2615 ! <IN NAME="allow_neg_inc" TYPE="logical, optional" DIM="(scalar)" DEFAULT=".true.">
2616 ! When .false., it is a fatal error if any of the input time increments are negative.
2617 ! This mimics the behavior of lima and earlier revisions.
2620 ! For all but the thirty_day_months calendar, decrements to months
2621 ! and years must be made separately from other units because of the
2622 ! non-associative nature of addition.
2623 ! If the result is a negative time (i.e. date before the base date)
2624 ! it is considered a fatal error.
2627 function decrement_date(Time, years, months, days, hours, minutes, seconds, ticks, err_msg, allow_neg_inc)
2629 type(time_type) :: decrement_date
2630 type(time_type), intent(in) :: Time
2631 integer, intent(in), optional :: seconds, minutes, hours, days, months, years, ticks
2632 character(len=*), intent(out), optional :: err_msg
2633 logical, intent(in), optional :: allow_neg_inc
2635 integer :: oseconds, ominutes, ohours, odays, omonths, oyears, oticks
2636 character(len=128) :: err_msg_local
2637 logical :: allow_neg_inc_local
2639 if(present(err_msg)) err_msg = ''
2641 ! Missing optionals are set to 0
2642 oseconds = 0; if(present(seconds)) oseconds = seconds
2643 ominutes = 0; if(present(minutes)) ominutes = minutes
2644 ohours = 0; if(present(hours)) ohours = hours
2645 odays = 0; if(present(days)) odays = days
2646 omonths = 0; if(present(months)) omonths = months
2647 oyears = 0; if(present(years)) oyears = years
2648 oticks = 0; if(present(ticks)) oticks = ticks
2649 allow_neg_inc_local=.true.; if(present(allow_neg_inc)) allow_neg_inc_local=allow_neg_inc
2651 if(.not.allow_neg_inc_local) then
2652 if(oyears < 0 .or. omonths < 0 .or. odays < 0 .or. ohours < 0 .or. ominutes < 0 .or. oseconds < 0 .or. &
2654 write(err_msg_local,10) oyears, omonths, odays, ohours, ominutes, oseconds, oticks
2655 if(error_handler('function decrement_date', err_msg_local, err_msg)) return
2658 10 format('One or more time increments are negative: '// &
2659 'years=',i6,' months=',i6,' days=',i6,' hours=',i6,' minutes=',i6,' seconds=',i6,' ticks=',i6)
2661 if(.not.increment_date_private( &
2662 Time, -oyears, -omonths, -odays, -ohours, -ominutes, -oseconds, -oticks, decrement_date, err_msg_local)) then
2663 if(error_handler('function decrement_date', err_msg_local, err_msg)) return
2666 end function decrement_date
2669 !=========================================================================
2670 ! START days_in_month BLOCK
2671 ! <FUNCTION NAME="days_in_month">
2674 ! Given a time interval, gives the number of days in the
2675 ! month corresponding to the default calendar.
2678 ! Given a time, computes the corresponding date given the selected
2679 ! date time mapping algorithm.
2681 ! <TEMPLATE> days_in_month(time) </TEMPLATE>
2683 ! <IN NAME="time" UNITS="" TYPE="time_type" DIM="">A time interval.</IN>
2684 ! <OUT NAME="days_in_month" UNITS="" TYPE="integer" DIM="" DEFAULT="">
2685 ! The number of days in the month given the selected time
2686 ! mapping algorithm.
2689 function days_in_month(Time, err_msg)
2691 ! Given a time, computes the corresponding date given the selected
2692 ! date time mapping algorithm
2694 integer :: days_in_month
2695 type(time_type), intent(in) :: Time
2696 character(len=*), intent(out), optional :: err_msg
2698 if(.not.module_is_initialized) call time_manager_init
2699 if(present(err_msg)) err_msg = ''
2701 select case(calendar_type)
2702 case(THIRTY_DAY_MONTHS)
2703 days_in_month = days_in_month_thirty(Time)
2705 days_in_month = days_in_month_gregorian(Time)
2707 days_in_month = days_in_month_julian(Time)
2709 days_in_month = days_in_month_no_leap(Time)
2711 if(error_handler('function days_in_month', &
2712 'days_in_month makes no sense when the calendar type is NO_CALENDAR', err_msg)) return
2714 if(error_handler('function days_in_month', 'Invalid calendar type', err_msg)) return
2716 end function days_in_month
2719 !--------------------------------------------------------------------------
2721 function days_in_month_gregorian(Time)
2723 ! Returns the number of days in a gregorian month.
2725 integer :: days_in_month_gregorian
2726 type(time_type), intent(in) :: Time
2727 integer :: year, month, day, hour, minute, second, ticks
2729 call get_date_gregorian(Time, year, month, day, hour, minute, second, ticks)
2730 days_in_month_gregorian = days_per_month(month)
2731 if(leap_year_gregorian_int(year) .and. month == 2) days_in_month_gregorian = 29
2733 end function days_in_month_gregorian
2735 !--------------------------------------------------------------------------
2736 function days_in_month_julian(Time)
2738 ! Returns the number of days in a julian month.
2740 integer :: days_in_month_julian
2741 type(time_type), intent(in) :: Time
2742 integer :: year, month, day, hour, minute, second, ticks
2744 call get_date_julian_private(Time, year, month, day, hour, minute, second, ticks)
2745 days_in_month_julian = days_per_month(month)
2746 if(leap_year_julian(Time) .and. month == 2) days_in_month_julian = 29
2748 end function days_in_month_julian
2750 !--------------------------------------------------------------------------
2751 function days_in_month_thirty(Time)
2753 ! Returns the number of days in a thirty day month (needed for transparent
2754 ! changes to calendar type).
2756 integer :: days_in_month_thirty
2757 type(time_type), intent(in) :: Time
2759 days_in_month_thirty = 30
2761 end function days_in_month_thirty
2763 !--------------------------------------------------------------------------
2764 function days_in_month_no_leap(Time)
2766 ! Returns the number of days in a 365 day year month.
2768 integer :: days_in_month_no_leap
2769 type(time_type), intent(in) :: Time
2770 integer :: year, month, day, hour, minute, second, ticks
2772 call get_date_no_leap_private(Time, year, month, day, hour, minute, second, ticks)
2773 days_in_month_no_leap= days_per_month(month)
2775 end function days_in_month_no_leap
2777 ! END OF days_in_month BLOCK
2778 !==========================================================================
2779 ! START OF leap_year BLOCK
2780 ! <FUNCTION NAME="leap_year">
2783 ! Returns true if the year corresponding to the input time is
2784 ! a leap year. Always returns false for THIRTY_DAY_MONTHS and NOLEAP.
2787 ! Returns true if the year corresponding to the input time is
2788 ! a leap year. Always returns false for THIRTY_DAY_MONTHS and NOLEAP.
2790 ! <TEMPLATE> leap_year(time) </TEMPLATE>
2792 ! <IN NAME="time" UNITS="" TYPE="time_type" DIM="">A time interval.</IN>
2793 ! <OUT NAME="leap_year" UNITS="" TYPE="calendar_type" DIM="" DEFAULT="">
2794 ! true if the year corresponding to the input time is a leap year.
2797 function leap_year(Time, err_msg)
2799 ! Is this date in a leap year for default calendar?
2801 logical :: leap_year
2802 type(time_type), intent(in) :: Time
2803 character(len=*), intent(out), optional :: err_msg
2805 if(.not.module_is_initialized) call time_manager_init
2806 if(present(err_msg)) err_msg=''
2808 select case(calendar_type)
2809 case(THIRTY_DAY_MONTHS)
2810 leap_year = leap_year_thirty(Time)
2812 leap_year = leap_year_gregorian(Time)
2814 leap_year = leap_year_julian(Time)
2816 leap_year = leap_year_no_leap(Time)
2818 if(error_handler('function leap_year', 'Invalid calendar type in leap_year', err_msg)) return
2820 end function leap_year
2823 !--------------------------------------------------------------------------
2825 function leap_year_gregorian(Time)
2827 ! Is this a leap year for gregorian calendar?
2829 logical :: leap_year_gregorian
2830 type(time_type), intent(in) :: Time
2831 integer :: seconds, minutes, hours, day, month, year
2833 call get_date(Time, year, month, day, hours, minutes, seconds)
2834 leap_year_gregorian = leap_year_gregorian_int(year)
2836 end function leap_year_gregorian
2838 !--------------------------------------------------------------------------
2840 function leap_year_gregorian_int(year)
2841 logical :: leap_year_gregorian_int
2842 integer, intent(in) :: year
2844 leap_year_gregorian_int = mod(year,4) == 0
2845 leap_year_gregorian_int = leap_year_gregorian_int .and. .not.mod(year,100) == 0
2846 leap_year_gregorian_int = leap_year_gregorian_int .or. mod(year,400) == 0
2848 end function leap_year_gregorian_int
2850 !--------------------------------------------------------------------------
2852 function leap_year_julian(Time)
2854 ! Returns the number of days in a julian month.
2856 logical :: leap_year_julian
2857 type(time_type), intent(in) :: Time
2858 integer :: seconds, minutes, hours, day, month, year
2860 call get_date(Time, year, month, day, hours, minutes, seconds)
2861 leap_year_julian = ((year / 4 * 4) == year)
2863 end function leap_year_julian
2865 !--------------------------------------------------------------------------
2867 function leap_year_thirty(Time)
2869 ! No leap years in thirty day months, included for transparency.
2871 logical :: leap_year_thirty
2872 type(time_type), intent(in) :: Time
2874 leap_year_thirty = .FALSE.
2876 end function leap_year_thirty
2878 !--------------------------------------------------------------------------
2880 function leap_year_no_leap(Time)
2882 ! Another tough one; no leap year returns false for leap year inquiry.
2884 logical :: leap_year_no_leap
2885 type(time_type), intent(in) :: Time
2887 leap_year_no_leap = .FALSE.
2889 end function leap_year_no_leap
2891 !END OF leap_year BLOCK
2892 !==========================================================================
2894 !> @brief Returns the mean length of the year in the default calendar setting.
2896 !> There are no arguments in this function. It returns the mean
2897 !! length of the year for the default calendar.
2898 function length_of_year()
2900 ! What is the length of the year for the default calendar type
2902 type(time_type) :: length_of_year
2904 if(.not.module_is_initialized) call time_manager_init
2906 select case(calendar_type)
2907 case(THIRTY_DAY_MONTHS)
2908 length_of_year = length_of_year_thirty()
2910 length_of_year = length_of_year_gregorian()
2912 length_of_year = length_of_year_julian()
2914 length_of_year = length_of_year_no_leap()
2916 call error_mesg('length_of_year','Invalid calendar type in length_of_year',FATAL)
2918 end function length_of_year
2921 !--------------------------------------------------------------------------
2923 function length_of_year_thirty()
2925 type(time_type) :: length_of_year_thirty
2927 length_of_year_thirty = set_time(0, 360)
2929 end function length_of_year_thirty
2931 !---------------------------------------------------------------------------
2933 function length_of_year_gregorian()
2935 type(time_type) :: length_of_year_gregorian
2937 length_of_year_gregorian = set_time(20952, 365) !20952 = 86500 * (days_in_400_yrs/400. - (days_in_400_yrs/400))
2939 end function length_of_year_gregorian
2941 !--------------------------------------------------------------------------
2943 function length_of_year_julian()
2945 type(time_type) :: length_of_year_julian
2947 length_of_year_julian = set_time(21600, 365) !21600 = (24/4) * 60 * 60
2949 end function length_of_year_julian
2951 !--------------------------------------------------------------------------
2953 function length_of_year_no_leap()
2955 type(time_type) :: length_of_year_no_leap
2957 length_of_year_no_leap = set_time(0, 365)
2959 end function length_of_year_no_leap
2961 !--------------------------------------------------------------------------
2963 ! END OF length_of_year BLOCK
2964 !==========================================================================
2966 !==========================================================================
2967 !> Returns number of day in year for given time. Jan 1st is day 1, not zero!
2968 function day_of_year(time)
2969 integer :: day_of_year
2970 type(time_type), intent(in) :: Time
2972 integer :: second, minute, hour, day, month, year
2973 type(time_type) :: t
2975 call get_date(time,year,month,day,hour,minute,second)
2976 t = time-set_date(year,1,1,0,0,0)
2977 day_of_year = t%days + 1
2980 ! START OF days_in_year BLOCK
2981 ! <FUNCTION NAME="days_in_year">
2983 !> @brief Returns the number of days in the calendar year corresponding to the date represented by
2984 !! time for the default calendar.
2985 !> @returns The number of days in this year for the default calendar type.
2986 function days_in_year(Time)
2988 ! What is the number of days in this year for the default calendar type
2990 integer :: days_in_year
2991 type(time_type), intent(in) :: Time !< A time interval
2993 if(.not.module_is_initialized) call time_manager_init
2995 select case(calendar_type)
2996 case(THIRTY_DAY_MONTHS)
2997 days_in_year = days_in_year_thirty(Time)
2999 days_in_year = days_in_year_gregorian(Time)
3001 days_in_year = days_in_year_julian(Time)
3003 days_in_year = days_in_year_no_leap(Time)
3005 call error_mesg('days_in_year','Invalid calendar type in days_in_year',FATAL)
3007 end function days_in_year
3010 !--------------------------------------------------------------------------
3012 function days_in_year_thirty(Time)
3014 integer :: days_in_year_thirty
3015 type(time_type), intent(in) :: Time
3017 days_in_year_thirty = 360
3019 end function days_in_year_thirty
3021 !---------------------------------------------------------------------------
3023 function days_in_year_gregorian(Time)
3025 integer :: days_in_year_gregorian
3026 type(time_type), intent(in) :: Time
3028 if(leap_year_gregorian(Time)) then
3029 days_in_year_gregorian = 366
3031 days_in_year_gregorian = 365
3034 end function days_in_year_gregorian
3036 !--------------------------------------------------------------------------
3037 function days_in_year_julian(Time)
3039 integer :: days_in_year_julian
3040 type(time_type), intent(in) :: Time
3042 if(leap_year_julian(Time)) then
3043 days_in_year_julian = 366
3045 days_in_year_julian = 365
3048 end function days_in_year_julian
3050 !--------------------------------------------------------------------------
3052 function days_in_year_no_leap(Time)
3054 integer :: days_in_year_no_leap
3055 type(time_type), intent(in) :: Time
3057 days_in_year_no_leap = 365
3059 end function days_in_year_no_leap
3061 !--------------------------------------------------------------------------
3063 ! END OF days_in_year BLOCK
3065 !> @brief Returns a character string containing the name of the
3066 !! month corresponding to month number n.
3068 !> Definition is the same for all calendar types.
3069 !! @returns The character string associated with a month. All calendars have 12 months and return
3070 !! full month names, not abreviations.
3071 function month_name(n)
3073 ! Returns character string associated with a month, for now, all calendars
3074 ! have 12 months and will return standard names.
3076 character (len=9) :: month_name
3077 integer, intent(in) :: n !< Month number
3078 character (len = 9), dimension(12) :: months = (/'January ', 'February ', &
3079 'March ', 'April ', 'May ', 'June ', 'July ', &
3080 'August ', 'September', 'October ', 'November ', 'December '/)
3082 if(.not.module_is_initialized) call time_manager_init
3084 if(n < 1 .or. n > 12) call error_mesg('month_name','Illegal month index',FATAL)
3086 month_name = months(n)
3088 end function month_name
3090 !==========================================================================
3092 !> The purpose of this routine is to prevent the addition of an excessive amount of code in order to implement
3093 !! the error handling scheme involving an optional error flag of type character.
3094 !! It allows one line of code to accomplish what would otherwise require 6 lines.
3095 !! A value of .true. for this function is a flag to the caller that it should immediately return to it's caller.
3096 function error_handler(routine, err_msg_local, err_msg)
3099 logical :: error_handler
3100 character(len=*), intent(in) :: routine, err_msg_local
3101 character(len=*), intent(out), optional :: err_msg
3103 error_handler = .false.
3104 if(present(err_msg)) then
3105 err_msg = err_msg_local
3106 error_handler = .true.
3108 call error_mesg(trim(routine),trim(err_msg_local),FATAL)
3111 end function error_handler
3113 !------------------------------------------------------------------------
3115 !> Initialization routine. Writes the version information to the log file
3116 subroutine time_manager_init ( )
3118 if (module_is_initialized) return ! silent return if already called
3120 call write_version_number("TIME_MANAGER_MOD", version)
3121 module_is_initialized = .true.
3123 end subroutine time_manager_init
3125 !------------------------------------------------------------------------
3127 !> @brief Prints the given time_type argument as a time (using days, seconds and ticks)
3129 !> @note There is no check for PE number.
3130 subroutine print_time (Time,str,unit)
3131 type(time_type) , intent(in) :: Time !< Time that will be printed
3132 character (len=*), intent(in), optional :: str !< Character string that precedes the printed time
3133 integer , intent(in), optional :: unit !< Unit number for printed output, defaults to stdout
3134 integer :: s,d,ticks, ns,nd,nt, unit_in
3135 character(len=19) :: fmt
3137 ! prints the time to standard output (or optional unit) as days and seconds
3138 ! NOTE: there is no check for PE number
3141 if (present(unit)) unit_in = unit
3143 call get_time (Time,s,d,ticks)
3146 ! get number of digits for days and seconds strings
3147 nd = int(log10(real(max(1,d))))+1
3148 ns = int(log10(real(max(1,s))))+1
3149 nt = int(log10(real(max(1,ticks))))+1
3150 write (fmt,10) nd, ns, nt
3151 10 format ('(a,i',i2.2,',a,i',i2.2,',a,i',i2.2,')')
3153 if (present(str)) then
3154 write (unit_in,fmt) trim(str)//' day=', d, ', sec=', s, ', ticks=', ticks
3156 write (unit_in,fmt) 'TIME: day=', d, ', sec=', s, ', ticks=', ticks
3159 end subroutine print_time
3161 !> @brief Prints the time to standard output (or optional unit) as a date.
3163 !! Prints the given time_type argument as a date (using year, month, day,
3164 !! hour, minutes, seconds and ticks). NOTE: there is no check for PE number.
3165 subroutine print_date (Time,str,unit)
3166 type(time_type) , intent(in) :: Time !< Time that will be printed
3167 character (len=*), intent(in), optional :: str !< Character string that precedes the printed time
3168 integer , intent(in), optional :: unit !< Unit number for printed output, defaults to stdout
3169 integer :: y,mo,d,h,m,s, unit_in
3170 character(len=9) :: mon
3172 ! prints the time to standard output (or optional unit) as a date
3173 ! NOTE: there is no check for PE number
3176 if (present(unit)) unit_in = unit
3178 call get_date (Time,y,mo,d,h,m,s)
3179 mon = month_name(mo)
3180 if (present(str)) then
3181 write (unit_in,10) trim(str)//' ', y,mon(1:3),' ',d,' ',h,':',m,':',s
3183 write (unit_in,10) 'DATE: ', y,mon(1:3),' ',d,' ',h,':',m,':',s
3185 10 format (a,i4,1x,a3,4(a1,i2.2))
3187 end subroutine print_date
3189 !------------------------------------------------------------------------
3191 !> @brief Returns a character string that describes the
3192 !! calendar type corresponding to the input integer.
3194 !> @returns A character string describing the calendar type
3195 function valid_calendar_types(ncal, err_msg)
3196 integer, intent(in) :: ncal !< Integer corresponding to a valid calendar type
3197 character(len=*), intent(out), optional :: err_msg !< Holds an error message when present
3198 character(len=24) :: valid_calendar_types
3199 character(len=128) :: err_msg_local
3201 if(.not.module_is_initialized) call time_manager_init
3203 if(present(err_msg)) err_msg = ''
3205 if(ncal == NO_CALENDAR) then
3206 valid_calendar_types = 'NO_CALENDAR '
3207 else if(ncal == THIRTY_DAY_MONTHS) then
3208 valid_calendar_types = '360_DAY '
3209 else if(ncal == JULIAN) then
3210 valid_calendar_types = 'JULIAN '
3211 else if(ncal == GREGORIAN) then
3212 valid_calendar_types = 'GREGORIAN '
3213 else if(ncal == NOLEAP) then
3214 valid_calendar_types = 'NOLEAP '
3216 write(err_msg_local,'(a,i4,a)') 'calendar type=',ncal,' is invalid.'
3217 if(error_handler('function valid_calendar_types', err_msg_local, err_msg)) return
3219 end function valid_calendar_types
3221 !------------------------------------------------------------------------
3223 !--- get the a character string that represents the time. The format will be
3224 !--- yyyymmdd.hhmmss
3225 function date_to_string(time, err_msg)
3226 type(time_type), intent(in) :: time
3227 character(len=*), intent(out), optional :: err_msg
3228 character(len=128) :: err_msg_local
3229 character(len=15) :: date_to_string
3230 integer :: yr,mon,day,hr,min,sec
3232 if(present(err_msg)) err_msg = ''
3233 call get_date(time,yr,mon,day,hr,min,sec)
3234 if (yr <= 9999) then
3235 write(date_to_string,'(I4.4,I2.2,I2.2,".",I2.2,I2.2,I2.2)') yr, mon, day, hr, min, sec
3237 write(err_msg_local, '(a,i4.4,a)') 'year = ', yr, ' should be less than 10000'
3238 if(error_handler('function date_to_string', err_msg_local, err_msg)) return
3241 end function date_to_string
3243 !> \author Tom Robinson thomas.robinson@noaa.gov
3244 !! \brief This routine converts the integer t%days to a string
3245 subroutine time_list_error (T,Terr)
3246 type(time_type), intent(in) :: t !< time_type input
3247 character(len=:), allocatable :: terr !< String holding the t%days
3248 !> Allocate the string
3249 allocate (character(len=10) :: terr)
3250 !> Write the integer to the string
3251 write (terr,'(I0)') t%days
3252 end subroutine time_list_error
3255 end module time_manager_mod
3259 ! <TESTPROGRAM NAME="time_main2">
3261 ! use time_manager_mod
3263 ! type(time_type) :: dt, init_date, astro_base_date, time, final_date
3264 ! type(time_type) :: next_rad_time, mid_date
3265 ! type(time_type) :: repeat_alarm_freq, repeat_alarm_length
3266 ! integer :: num_steps, i, days, months, years, seconds, minutes, hours
3267 ! integer :: months2, length
3268 ! real :: astro_days
3270 ! !Set calendar type
3271 ! ! call set_calendar_type(THIRTY_DAY_MONTHS)
3272 ! call set_calendar_type(JULIAN)
3273 ! ! call set_calendar_type(NOLEAP)
3276 ! dt = set_time(1100, 0)
3278 ! ! Set initial date
3279 ! init_date = set_date(1992, 1, 1)
3281 ! ! Set date for astronomy delta calculation
3282 ! astro_base_date = set_date(1970, 1, 1, 12, 0, 0)
3284 ! ! Copy initial time to model current time
3287 ! ! Determine how many steps to do to run one year
3288 ! final_date = increment_date(init_date, years = 1)
3289 ! num_steps = (final_date - init_date) / dt
3290 ! write(*, *) 'Number of steps is' , num_steps
3292 ! ! Want to compute radiation at initial step, then every two hours
3293 ! next_rad_time = time + set_time(7200, 0)
3295 ! ! Test repeat alarm
3296 ! repeat_alarm_freq = set_time(0, 1)
3297 ! repeat_alarm_length = set_time(7200, 0)
3299 ! ! Loop through a year
3300 ! do i = 1, num_steps
3305 ! ! Test repeat alarm
3306 ! if(repeat_alarm(time, repeat_alarm_freq, repeat_alarm_length)) &
3307 ! write(*, *) 'REPEAT ALARM IS TRUE'
3309 ! ! Should radiation be computed? Three possible tests.
3310 ! ! First test assumes exact interval; just ask if times are equal
3311 ! ! if(time == next_rad_time) then
3312 ! ! Second test computes rad on last time step that is <= radiation time
3313 ! ! if((next_rad_time - time) < dt .and. time < next_rad) then
3314 ! ! Third test computes rad on time step closest to radiation time
3315 ! if(interval_alarm(time, dt, next_rad_time, set_time(7200, 0))) then
3316 ! call get_date(time, years, months, days, hours, minutes, seconds)
3317 ! write(*, *) days, month_name(months), years, hours, minutes, seconds
3319 ! ! Need to compute real number of days between current time and astro_base
3320 ! call get_time(time - astro_base_date, seconds, days)
3321 ! astro_days = days + seconds / 86400.
3322 ! ! write(*, *) 'astro offset ', astro_days
3325 ! ! Can compute daily, monthly, yearly, hourly, etc. diagnostics as for rad
3327 ! ! Example: do diagnostics on last time step of this month
3328 ! call get_date(time + dt, years, months2, days, hours, minutes, seconds)
3329 ! call get_date(time, years, months, days, hours, minutes, seconds)
3330 ! if(months /= months2) then
3331 ! write(*, *) 'last timestep of month'
3332 ! write(*, *) days, months, years, hours, minutes, seconds
3335 ! ! Example: mid-month diagnostics; inefficient to make things clear
3336 ! length = days_in_month(time)
3337 ! call get_date(time, years, months, days, hours, minutes, seconds)
3338 ! mid_date = set_date(years, months, 1) + set_time(0, length) / 2
3340 ! if(time < mid_date .and. (mid_date - time) < dt) then
3341 ! write(*, *) 'mid-month time'
3342 ! write(*, *) days, months, years, hours, minutes, seconds
3348 ! end program time_main2
3352 ! The <a name="base date">base date</a> is implicitly defined so users don't
3353 ! need to be concerned with it. For the curious, the base date is defined as
3354 ! 0 seconds, 0 minutes, 0 hours, day 1, month 1, year 1
3357 ! Please note that a time is a positive definite quantity.
3360 ! See the <LINK SRC="TEST PROGRAM">Test Program </LINK> for a simple program
3361 ! that shows some of the capabilities of the time manager.
3365 ! close documentation grouping