!== dc_date.f90 - 日付および時刻に関する手続きを提供するモジュール
!
! Authors:: Yasuhiro MORIKAWA, Eizi TOYODA
! Version:: $Id: dc_date.f90,v 1.11 2009-01-14 10:30:59 morikawa Exp $
! Tag Name:: $Name: gtool5-20090217 $
! Copyright:: Copyright (C) GFD Dennou Club, 2000-2005. All rights reserved.
! License:: See COPYRIGHT[link:../../COPYRIGHT]
module dc_date
!== Overview
!
! 日付および時刻を扱うための手続きを提供します
!
!== Procedures List
!
! 以下の手続きは構造型 dc_date_types#DC_DATETIME または
! dc_date_types#DC_DIFFTIME 変数 (日時, 時刻に関する情報を格納)
! を対象とします.
!
! DCDateTimeCreate :: dc_date_types#DC_DATETIME 型変数の初期設定
! DCDiffTimeCreate :: dc_date_types#DC_DIFFTIME 型変数の初期設定
!
! #assignment(=) :: dc_date_types#DC_DATETIME 型変数および
! dc_date_types#DC_DIFFTIME 型変数の初期設定
!
! Eval :: 日時, 時刻情報を個別に取得
!
! toChar :: 日時, 時刻情報を文字型変数へ変換
!
! EvalDay :: 日数 (実数型) に換算して取得
! EvalHour :: 時間 (実数型) に換算して取得
! EvalMin :: 分 (実数型) に換算して取得
! EvalSec :: 秒 (実数型) に換算して取得
! EvalNondim :: 無次元時間 (実数型) に換算して取得
! EvalByUnit :: 単位を指定し, 日, 時, 分, 秒のいづれか (実数型)
! に換算して取得
!
! #operator(+) :: 加算 (dc_date_types#DC_DATETIME 型 および
! dc_date_types#DC_DIFFTIME 型 同士)
! #operator(-) :: 減算 (dc_date_types#DC_DATETIME 型 および
! dc_date_types#DC_DIFFTIME 型 同士)
! #operator(*) :: 乗算 (dc_date_types#DC_DIFFTIME 型と数値型)
! #operator(/) :: 除算 (dc_date_types#DC_DIFFTIME 型と数値型)
! mod :: 余り (dc_date_types#DC_DIFFTIME 型同士)
! #operator(==) :: 比較 (dc_date_types#DC_DATETIME 型同士)
! #operator(>) :: 比較 (dc_date_types#DC_DATETIME 型同士)
! #operator(<) :: 比較 (dc_date_types#DC_DATETIME 型同士)
! max :: 大きい値を返す
! min :: 小さい値を返す
!
! SetZone :: タイムゾーンを変更
!
! DCDateTimePutLine :: dc_date_types#DC_DATETIME 型変数に格納されている
! 日時, 時刻情報の印字
! DCDiffTimePutLine :: dc_date_types#DC_DIFFTIME 型変数に格納されている
! 日時, 時刻情報の印字
!
!
! 以下の手続きは dc_date_types 内部の変数を変更します.
!
! SetCaltype :: 暦法のデフォルトを変更
! SetSecOfDay :: 1 日の秒数のデフォルトを変更
!
! その他の手続き
!
! ValidCaltype :: 暦法が有効なものかをチェック
! ValidZone :: タイムゾーンとして有効化をチェック
! ZoneToDiff :: タイムゾーンを dc_date_types#DC_DIFFTIME 変数へと変換
! ParseTimeUnits :: 時間の単位を解析し, 単位のシンボルを返します.
!
!== Usage
!
!=== 現在時刻の表示
!
! DC_DATETIME 型の変数にサブルーチン DCDateTimeCreate
! を用いると, 時刻が設定されます.
! 下記のように特に年月日を指定しないと現在時刻が設定されます.
! 設定された時刻は toChar によって文字型変数へと変換できます.
! サブルーチン Printf に関しては dc_string#Printf を参照ください.
!
! program dc_date_sapmle1
! use dc_string, only: Printf
! use dc_date_types, only: DC_DATETIME
! use dc_date, only: DCDateTimeCreate, toChar
! implicit none
! type(DC_DATETIME) :: time
!
! call DCDateTimeCreate( time = time ) ! (out)
! call Printf( fmt = 'current date and time is %c', c1 = trim( toChar(time) ) )
! end program dc_date_sapmle1
!
!=== 日時, 時刻情報の加算
!
! DC_DIFFTIME 型の変数は日時差を表現します. 下記の例では,
! 日時差を表現するための変数として *diff*
! を用意し, サブルーチン Create によって 25 日 + 12 時間 + 50 分の日時差
! を設定しています. DC_DATETIME 型の変数 *time_before* と *diff* とを
! #operator(+) によって加算することで *time_before* から
! 25 日 + 12 時間 + 50 分を進めた日時 *time_after* を取得しています.
!
! program dc_date_sapmle2
! use dc_types, only: DP
! use dc_string, only: Printf
! use dc_date_types, only: DC_DATETIME, DC_DIFFTIME
! use dc_date, only: DCDateTimeCreate, DCDiffTimeCreate, toChar, operator(+)
! implicit none
! type(DC_DATETIME) :: time_before, time_after
! type(DC_DIFFTIME) :: diff
!
! call DCDateTimeCreate( time = time_before, & ! (out)
! & year = 2006, mon = 6, day = 10, & ! (in)
! & hour = 14, min = 15, sec = 0.0_DP ) ! (in)
! call DCDiffTimeCreate( diff = diff, & ! (out)
! & day = 25, hour = 12, min = 50) ! (in)
!
! time_after = time_before + diff
!
! call Printf( fmt = '%c + %c = %c', &
! & c1 = trim( toChar(time_before) ), c2 = trim( toChar(diff) ), &
! & c3 = trim( toChar(time_after) ) )
! end program dc_date_sapmle2
!
!
!=== 時間積分のループへの応用
!
! 以下は dA/dt = - αA (初期値 1, α=0.0001) を t = 12 (時間)
! まで解くプログラムの例です. 時間積分には前進差分を用いています.
! Δt, データの出力間隔, 計算時間に DC_DIFFTIME を用いることで,
! ループの終了処理や
! データ出力の際の時刻の比較が容易となります.
!
! program dc_date_sapmle3
! use dc_types, only: DP
! use dc_date, only: DCDiffTimeCreate, EvalSec, EvalByUnit, mod, &
! & operator(*), operator(==), operator(>)
! use dc_date_types, only: DC_DIFFTIME
! implicit none
! real(DP) :: func_a = 1.0d0 ! 関数 A の初期値
! real(DP), parameter :: alph = 0.0001d0 ! 係数 α
! character(*), parameter :: out_unit = 'hour' ! 出力される時刻の単位
! type(DC_DIFFTIME):: DelTimef, intervalf, calctimef
! integer :: i
! continue
! call DCDiffTimeCreate( & ! Δt = 5.0 (秒)
! & diff = DelTimef, & ! (out)
! & value = 5.0_DP, unit = 'sec') ! (in)
! call DCDiffTimeCreate( & ! データ出力間隔 = 1.0 (分)
! & diff = intervalf, & ! (out)
! & value = 1.0_DP, unit = 'min') ! (in)
! call DCDiffTimeCreate( & ! 計算時間 = 12.0 (時間)
! & diff = calctimef, & ! (out)
! & value = 12.0_DP, unit = 'hour') ! (in)
!
! open( 10, file='dc_date_sample.dat' )
! write(10,'(A,A,A)') '# ', out_unit, ' value'
!
! i = 1
! do
! if (DelTimef * i > calctimef) exit ! 計算時間を過ぎたら終了
!
! !---------------------------------------------
! ! A_(n+1) = (1 - αΔt) * A_(n)
! !---------------------------------------------
! func_a = (1.0 - alph * EvalSec(DelTimef)) * func_a
!
! !---------------------------------------------
! ! intervalf (1 分) 毎にデータを出力
! !---------------------------------------------
! if (mod(DelTimef * i, intervalf) == 0) then
! write(10,*) ' ', EvalByUnit( DelTimef * i, out_unit ), func_a
! end if
! i = i + 1
! end do
! end program dc_date_sapmle3
!
!
! 日付時刻と時間間隔を区別する。
! 型宣言によって自明に定まるサブルーチンは dc_date_types に置く。
use dc_date_types, only: DC_DATETIME, DC_DIFFTIME
use dc_types, only: DP, STRING, TOKEN
use dc_present, only: present_and_not_empty
implicit none
private
public:: DCDateTimeCreate, DCDiffTimeCreate
public:: DCDateTimePutLine, DCDiffTimePutLine
public:: Eval
public:: SetCaltype, SetZone, SetSecOfDay
public:: ValidCaltype, ValidZone, ZoneToDiff, ParseTimeUnits
public:: assignment(=)
public:: mod, operator(/), operator(-), operator(+), operator(*)
public:: operator(<), operator(>), operator(>=), operator(<=)
public:: operator(==), max, min
public:: toChar, toCharCal
public:: EvalDay, EvalHour, EvalMin, EvalSec, EvalNondim, EvalByUnit
public:: dcdate_normalize, dcdate_parse_unit, dcdate_set_day_seconds_scl
!-----------------------------------------------
! 後方互換用
! For backward compatibility
public:: Create, PutLine
interface DCDateTimeCreate
subroutine DCDateTimeCreate1(time, &
& year, mon, day, hour, min, sec, &
& zone, zone_hour, zone_min, caltype, caltype_str, day_seconds, &
& sclyear, sclmon, sclday, sclsec, err) !:doc-priority 40:
use dc_types, only: DP
use dc_date_types, only: DC_DATETIME
use dc_scaledsec, only: DC_SCALED_SEC
type(DC_DATETIME), intent(out):: time
integer, intent(in), optional:: year, mon, day, hour, min
real(DP),intent(in), optional:: sec, day_seconds
character(*), intent(in), optional :: zone
integer, intent(in), optional :: zone_hour
integer, intent(in), optional :: zone_min
integer, intent(in), optional:: caltype
character(*), intent(in), optional:: caltype_str
type(DC_SCALED_SEC), intent(in), optional:: sclyear, sclmon, sclday, sclsec
logical, intent(out), optional:: err
end subroutine DCDateTimeCreate1
end interface
interface DCDiffTimeCreate
subroutine DCDiffTimeCreate1(diff, &
& year, mon, day, hour, min, sec, day_seconds, nondim, &
& sclyear, sclmon, sclday, sclsec ) !:doc-priority 60:
use dc_types, only: DP
use dc_date_types, only: DC_DIFFTIME
use dc_scaledsec, only: DC_SCALED_SEC
type(DC_DIFFTIME), intent(out) :: diff
integer, intent(in), optional:: year, mon, day, hour, min
real(DP),intent(in), optional:: sec, day_seconds, nondim
type(DC_SCALED_SEC), intent(in), optional:: sclyear, sclmon, sclday, sclsec
end subroutine DCDiffTimeCreate1
subroutine DCDiffTimeCreate2D(diff, value, unit, unit_symbol, err) !:doc-priority 70:
use dc_types, only: DP
use dc_date_types, only: DC_DIFFTIME
type(DC_DIFFTIME), intent(out) :: diff
real(DP), intent(in) :: value
character(*), intent(in) :: unit
integer, intent(in), optional :: unit_symbol
logical, intent(out), optional :: err
end subroutine DCDiffTimeCreate2D
subroutine DCDiffTimeCreate2R(diff, value, unit, unit_symbol, err) !:doc-priority 80:
use dc_date_types, only: DC_DIFFTIME
type(DC_DIFFTIME), intent(out) :: diff
real, intent(in) :: value
character(*), intent(in) :: unit
integer, intent(in), optional :: unit_symbol
logical, intent(out), optional :: err
end subroutine DCDiffTimeCreate2R
subroutine DCDiffTimeCreate2I(diff, value, unit, unit_symbol, err) !:doc-priority 90:
use dc_date_types, only: DC_DIFFTIME
type(DC_DIFFTIME), intent(out) :: diff
integer, intent(in) :: value
character(*), intent(in) :: unit
integer, intent(in), optional :: unit_symbol
logical, intent(out), optional :: err
end subroutine DCDiffTimeCreate2I
end interface
interface DCDateTimePutLine
subroutine DCDateTimePutLine(time, unit, indent)
use dc_date_types, only: DC_DATETIME
type(DC_DATETIME), intent(in) :: time
integer, intent(in), optional :: unit
character(*), intent(in), optional:: indent
end subroutine DCDateTimePutLine
end interface
interface DCDiffTimePutLine
subroutine DCDiffTimePutLine(diff, unit, indent)
use dc_date_types, only: DC_DIFFTIME
type(DC_DIFFTIME), intent(in) :: diff
integer, intent(in), optional :: unit
character(*), intent(in), optional:: indent
end subroutine DCDiffTimePutLine
end interface
interface assignment(=)
subroutine DCDateTimeCreateI(time, sec) !:doc-priority 20:
use dc_date_types, only: DC_DATETIME
type(DC_DATETIME), intent(out):: time
integer, intent(in):: sec
end subroutine DCDateTimeCreateI
subroutine DCDateTimeCreateR(time, sec) !:doc-priority 30:
use dc_date_types, only: DC_DATETIME
type(DC_DATETIME), intent(out):: time
real, intent(in):: sec
end subroutine DCDateTimeCreateR
subroutine DCDateTimeCreateD(time, sec) !:doc-priority 40:
use dc_types, only: DP
use dc_date_types, only: DC_DATETIME
type(DC_DATETIME), intent(out):: time
real(DP), intent(in):: sec
end subroutine DCDateTimeCreateD
subroutine DCDiffTimeCreateI(diff, sec) !:doc-priority 60:
use dc_date_types, only: DC_DIFFTIME
type(DC_DIFFTIME), intent(out):: diff
integer, intent(in):: sec
end subroutine DCDiffTimeCreateI
subroutine DCDiffTimeCreateR(diff, sec) !:doc-priority 70:
use dc_date_types, only: DC_DIFFTIME
type(DC_DIFFTIME), intent(out):: diff
real, intent(in):: sec
end subroutine DCDiffTimeCreateR
subroutine DCDiffTimeCreateD(diff, sec) !:doc-priority 80:
use dc_types, only: DP
use dc_date_types, only: DC_DIFFTIME
type(DC_DIFFTIME), intent(out):: diff
real(DP), intent(in):: sec
end subroutine DCDiffTimeCreateD
!!$ subroutine DCDateLetFC(diff, string)
!!$ use dc_date_types, only: DC_DIFFTIME
!!$ type(DC_DIFFTIME), intent(out):: diff
!!$ character(len = *), intent(in):: string
!!$ end subroutine DCDateLetFC
!!$
!!$ subroutine DCDateLetTC(time, string)
!!$ use dc_date_types, only: DC_DATETIME
!!$ type(DC_DATETIME), intent(out):: time
!!$ character(len = *), intent(in):: string
!!$ end subroutine DCDateLetTC
end interface
interface SetCaltype
subroutine DCDateTimeSetCaltype(caltype)
integer, intent(in):: caltype
end subroutine DCDateTimeSetCaltype
end interface
interface SetSecOfDay
subroutine DCDateTimeSetSecOfDay(sec)
use dc_types, only: DP
real(DP), intent(in):: sec
end subroutine DCDateTimeSetSecOfDay
end interface
interface ValidCaltype
function DCDateTimeValidCaltype(caltype) result(result)
integer, intent(in):: caltype
logical:: result
end function DCDateTimeValidCaltype
end interface
interface ValidZone
function DCDateTimeValidZone(zone) result(result)
character(*), intent(in):: zone
logical:: result
end function DCDateTimeValidZone
end interface
interface ZoneToDiff
function DCDateTimeZoneToDiff(zone) result(diff)
use dc_date_types, only: DC_DIFFTIME
type(DC_DIFFTIME):: diff
character(*), intent(in):: zone
end function DCDateTimeZoneToDiff
end interface
interface SetZone
subroutine DCDateTimeSetZone(time, zone, err)
use dc_date_types, only: DC_DATETIME
type(DC_DATETIME), intent(inout):: time
character(*), intent(in):: zone
logical, intent(out), optional:: err
end subroutine DCDateTimeSetZone
end interface
interface Eval
subroutine DCDateTimeEval1(time, year, mon, day, hour, min, &
& sec, caltype, zone, sclyear, sclmon, sclday, sclsec) !:doc-priority 40:
use dc_types, only: DP
use dc_date_types, only: DC_DATETIME
use dc_scaledsec, only: DC_SCALED_SEC
type(DC_DATETIME), intent(in):: time
integer, intent(out), optional:: year, mon, day, hour, min, caltype
real(DP), intent(out), optional:: sec
character(*), intent(out), optional:: zone
type(DC_SCALED_SEC), intent(out), optional:: sclyear, sclmon, sclday, sclsec
end subroutine DCDateTimeEval1
!!$ subroutine DCDateTimeEval0(time, mon, day, sec)
!!$ use dc_date_types, only: DC_DATETIME
!!$ use dc_types, only: DP
!!$ type(DC_DATETIME), intent(in):: time
!!$ integer, intent(out):: mon, day
!!$ real(DP), intent(out):: sec
!!$ end subroutine DCDateTimeEval0
subroutine DCDiffTimeEval1(diff, &
& year, mon, day, hour, min, sec, nondim, &
& sclyear, sclmon, sclday, sclsec, sclnondim, err) !:doc-priority 60:
use dc_types, only: DP
use dc_date_types, only: DC_DIFFTIME
use dc_scaledsec, only: DC_SCALED_SEC
type(DC_DIFFTIME), intent(in):: diff
integer, intent(out), optional:: year, mon, day, hour, min
real(DP), intent(out), optional:: sec, nondim
type(DC_SCALED_SEC), intent(out), optional:: sclyear, sclmon, sclday, sclsec, sclnondim
logical, intent(out), optional :: err
end subroutine DCDiffTimeEval1
end interface
interface EvalDay
function DCDateTimeEvalDay(time) result(result) !:doc-priority 40:
use dc_types, only: DP
use dc_date_types, only: DC_DATETIME
real(DP):: result
type(DC_DATETIME), intent(in):: time
end function DCDateTimeEvalDay
function DCDiffTimeEvalDay(diff) result(result) !:doc-priority 60:
use dc_types, only: DP
use dc_date_types, only: DC_DIFFTIME
real(DP):: result
type(DC_DIFFTIME), intent(in):: diff
end function DCDiffTimeEvalDay
end interface
interface EvalHour
function DCDateTimeEvalHour(time) result(result) !:doc-priority 40:
use dc_types, only: DP
use dc_date_types, only: DC_DATETIME
real(DP):: result
type(DC_DATETIME), intent(in):: time
end function DCDateTimeEvalHour
function DCDifftimeEvalHour(diff) result(result) !:doc-priority 60:
use dc_types, only: DP
use dc_date_types, only: DC_DIFFTIME
real(DP):: result
type(DC_DIFFTIME), intent(in):: diff
end function DCDifftimeEvalHour
end interface
interface EvalMin
function DCDateTimeEvalMin(time) result(result) !:doc-priority 40:
use dc_types, only: DP
use dc_date_types, only: DC_DATETIME
real(DP):: result
type(DC_DATETIME), intent(in):: time
end function DCDateTimeEvalMin
function DCDifftimeEvalMin(diff) result(result) !:doc-priority 60:
use dc_types, only: DP
use dc_date_types, only: DC_DIFFTIME
real(DP):: result
type(DC_DIFFTIME), intent(in):: diff
end function DCDifftimeEvalMin
end interface
interface EvalSec
function DCDateTimeEvalSec(time) result(result) !:doc-priority 40:
use dc_types, only: DP
use dc_date_types, only: DC_DATETIME
real(DP):: result
type(DC_DATETIME), intent(in):: time
end function DCDateTimeEvalSec
function DCDifftimeEvalSec(diff) result(result) !:doc-priority 60:
use dc_types, only: DP
use dc_date_types, only: DC_DIFFTIME
real(DP):: result
type(DC_DIFFTIME), intent(in):: diff
end function DCDifftimeEvalSec
end interface
interface EvalNondim
function DCDiffTimeEvalNondim(diff) result(result)
use dc_types, only: DP
use dc_date_types, only: DC_DIFFTIME
real(DP):: result
type(DC_DIFFTIME), intent(in):: diff
end function DCDiffTimeEvalNondim
end interface
interface EvalSclSec
function DCDateTimeEvalSclSec(time) result(result) !:doc-priority 40:
use dc_date_types, only: DC_DATETIME
use dc_scaledsec, only: DC_SCALED_SEC
type(DC_SCALED_SEC):: result
type(DC_DATETIME), intent(in):: time
end function DCDateTimeEvalSclSec
function DCDifftimeEvalSclSec(diff) result(result) !:doc-priority 60:
use dc_date_types, only: DC_DIFFTIME
use dc_scaledsec, only: DC_SCALED_SEC
type(DC_SCALED_SEC):: result
type(DC_DIFFTIME), intent(in):: diff
end function DCDifftimeEvalSclSec
end interface
interface EvalByUnit
function DCDateTimeEvalByUnit(time, unit, unit_symbol) result(result)
use dc_types, only: DP, TOKEN
use dc_date_types, only: DC_DATETIME
real(DP):: result
type(DC_DATETIME), intent(in):: time
character(*), intent(in), optional:: unit
integer, intent(in), optional:: unit_symbol
end function DCDateTimeEvalByUnit
function DCDiffTimeEvalByUnit(diff, unit, unit_symbol) result(result)
use dc_types, only: DP, TOKEN
use dc_date_types, only: DC_DIFFTIME
real(DP):: result
type(DC_DIFFTIME), intent(in):: diff
character(*), intent(in), optional:: unit
integer, intent(in), optional:: unit_symbol
end function DCDiffTimeEvalByUnit
end interface
interface toChar
function DCDateTimeToChar(time) result(result) !:doc-priority 40:
use dc_types, only: STRING
use dc_date_types, only: DC_DATETIME
character(STRING) :: result
type(DC_DATETIME), intent(in):: time
end function DCDateTimeToChar
function DCDiffTimeToChar(diff) result(result) !:doc-priority 60:
use dc_types, only: STRING
use dc_date_types, only: DC_DIFFTIME
character(STRING) :: result
type(DC_DIFFTIME), intent(in):: diff
end function DCDiffTimeToChar
end interface
interface toCharCal
function DCDateTimeToCharCal(time, upcase) result(result)
use dc_types, only: TOKEN
use dc_date_types, only: DC_DATETIME
character(TOKEN) :: result
type(DC_DATETIME), intent(in):: time
logical, intent(in), optional:: upcase
end function DCDateTimeToCharCal
end interface
interface operator(+)
module procedure dcdate_add_ft
module procedure dcdate_add_tf
module procedure dcdate_add_ff
module procedure dcdate_add_fr
module procedure dcdate_add_fd
module procedure dcdate_add_fi
end interface
interface operator(-)
module procedure dcdate_sub_tf !:doc-priority 40:
module procedure dcdate_sub_tt
module procedure dcdate_sub_ff
module procedure dcdate_sub_fr
module procedure dcdate_sub_fd
module procedure dcdate_sub_fi
end interface
interface operator(*)
module procedure dcdate_mul_if !:doc-priority 51:
module procedure dcdate_mul_fi !:doc-priority 52:
module procedure dcdate_mul_rf !:doc-priority 61:
module procedure dcdate_mul_fr !:doc-priority 62:
module procedure dcdate_mul_df !:doc-priority 71:
module procedure dcdate_mul_fd !:doc-priority 72:
end interface
interface operator(/)
module procedure dcdate_div_fi
module procedure dcdate_div_fr
module procedure dcdate_div_fd
module procedure dcdate_div_ff
end interface
interface mod
module procedure dcdate_mod_ff
end interface
interface operator(==)
module procedure dcdate_eq_tt !:doc-priority 30:
module procedure dcdate_eq_ff !:doc-priority 40:
module procedure dcdate_eq_if !:doc-priority 51:
module procedure dcdate_eq_fi !:doc-priority 52:
module procedure dcdate_eq_rf !:doc-priority 61:
module procedure dcdate_eq_fr !:doc-priority 62:
module procedure dcdate_eq_df !:doc-priority 71:
module procedure dcdate_eq_fd !:doc-priority 72:
end interface
interface operator(>)
module procedure dcdate_gt_tt !:doc-priority 30:
module procedure dcdate_gt_ff !:doc-priority 40:
module procedure dcdate_gt_fi !:doc-priority 42:
module procedure dcdate_gt_if !:doc-priority 44:
end interface
interface operator(<)
module procedure dcdate_lt_tt !:doc-priority 30:
module procedure dcdate_lt_ff !:doc-priority 40:
module procedure dcdate_lt_fi !:doc-priority 42:
module procedure dcdate_lt_if !:doc-priority 44:
end interface
interface operator(>=)
module procedure dcdate_ge_tt !:doc-priority 30:
module procedure dcdate_ge_ff !:doc-priority 40:
module procedure dcdate_ge_fi !:doc-priority 42:
module procedure dcdate_ge_if !:doc-priority 44:
end interface
interface operator(<=)
module procedure dcdate_le_tt !:doc-priority 30:
module procedure dcdate_le_ff !:doc-priority 40:
module procedure dcdate_le_fi !:doc-priority 42:
module procedure dcdate_le_if !:doc-priority 44:
end interface
interface max
module procedure dcdate_max_tt !:doc-priority 30:
module procedure dcdate_max_ff !:doc-priority 40:
end interface
interface min
module procedure dcdate_min_tt !:doc-priority 30:
module procedure dcdate_min_ff !:doc-priority 40:
end interface
!-----------------------------------------------
! 後方互換用
! For backward compatibility
interface Create
subroutine DCDateTimeCreate1_bc(time, &
& year, mon, day, hour, min, sec, &
& zone, caltype, day_seconds, err) !:doc-priority 40:
use dc_types, only: DP
use dc_date_types, only: DC_DATETIME
use dc_scaledsec, only: DC_SCALED_SEC
type(DC_DATETIME), intent(out):: time
integer, intent(in), optional:: year, mon, day, hour, min
real(DP),intent(in), optional:: sec, day_seconds
character(*), intent(in), optional :: zone
integer, intent(in), optional:: caltype
logical, intent(out), optional:: err
end subroutine DCDateTimeCreate1_bc
subroutine DCDiffTimeCreate1_bc(diff, &
& year, mon, day, hour, min, sec, day_seconds ) !:doc-priority 60:
use dc_types, only: DP
use dc_date_types, only: DC_DIFFTIME
use dc_scaledsec, only: DC_SCALED_SEC
type(DC_DIFFTIME), intent(out) :: diff
integer, intent(in), optional:: year, mon, day, hour, min
real(DP),intent(in), optional:: sec, day_seconds
end subroutine DCDiffTimeCreate1_bc
subroutine DCDiffTimeCreate2_bc(diff, value, unit, unit_symbol, err) !:doc-priority 70:
use dc_types, only: DP
use dc_date_types, only: DC_DIFFTIME
type(DC_DIFFTIME), intent(out) :: diff
real(DP), intent(in) :: value
character(*), intent(in) :: unit
integer, intent(in), optional :: unit_symbol
logical, intent(out), optional :: err
end subroutine DCDiffTimeCreate2_bc
end interface
interface PutLine
subroutine DCDateTimePutLine_bc(time, unit)
use dc_date_types, only: DC_DATETIME
type(DC_DATETIME), intent(in) :: time
integer, intent(in), optional :: unit
end subroutine DCDateTimePutLine_bc
subroutine DCDiffTimePutLine_bc(diff, unit)
use dc_date_types, only: DC_DIFFTIME
type(DC_DIFFTIME), intent(in) :: diff
integer, intent(in), optional :: unit
end subroutine DCDiffTimePutLine_bc
end interface
contains
subroutine dcdate_normalize(day, sec, day_seconds, nondim_flag)
!
!=== 日と秒の正規化
!
! このサブルーチンは内部向けなので dc_date モジュール外では
! 極力使用しないでください.
!
! 日付 *day* と秒数 *sec* の正規化を行います. *sec* が *day_seconds*
! (省略される場合は dc_date_types#day_seconds) を超える場合, *day*
! に繰上げを行います.
! また, *sec* と *day* の符号が逆の場合, 同符号になるよう
! 設定します.
!
use dc_date_types, only: &
& flag_set_day_seconds_scl, day_seconds_scl
use dc_scaledsec, only: DC_SCALED_SEC, &
& operator(<), operator(>), operator(<=), operator(>=), &
& operator(+), operator(-), operator(*), operator(/), &
& modulo, int, abs, sign
implicit none
type(DC_SCALED_SEC), intent(inout):: day
type(DC_SCALED_SEC), intent(inout):: sec
type(DC_SCALED_SEC), intent(in), optional:: day_seconds
logical, intent(in):: nondim_flag
type(DC_SCALED_SEC):: sgn, day_sec, zero_sec
continue
if ( nondim_flag ) return
if (present(day_seconds)) then
day_sec = day_seconds
else
if ( .not. flag_set_day_seconds_scl ) call dcdate_set_day_seconds_scl
day_sec = day_seconds_scl
end if
if (abs(sec) >= day_sec) then
day = day + int(sec / day_sec)
sec = modulo(sec, day_sec)
end if
!! zero_sec = 0 (デフォルト値 = 0 を使用する).
if ( ( sec > zero_sec .and. day < zero_sec ) &
& .or. ( sec < zero_sec .and. day > zero_sec ) ) then
sgn = sign(day, 1)
day = day - sgn
sec = sec + sgn * day_sec
endif
end subroutine dcdate_normalize
subroutine dcdate_set_day_seconds_scl
use dc_scaledsec, only: DC_SCALED_SEC, assignment(=)
use dc_date_types, only: day_seconds, &
& flag_set_day_seconds_scl, day_seconds_scl
continue
if ( .not. flag_set_day_seconds_scl ) then
flag_set_day_seconds_scl = .true.
day_seconds_scl = day_seconds
end if
end subroutine dcdate_set_day_seconds_scl
subroutine dcdate_nondimcheck(opr, diff1, diff2, rslt)
!
! このサブルーチンは内部向けなので dc_date モジュール外では
! 極力使用しないでください.
!
! diff1 と diff2 が両方とも有次元もしくは無次元かをチェックし,
! 両方が同じであれば, その結果を rslt に適用します.
! 2つの引数で片方が有次元, もう片方が無次元の場合には
! エラーを発生させます.
!
use dc_error, only: StoreError, DC_EDIMTIME
implicit none
character(*), intent(in):: opr ! 演算子の名称
type(DC_DIFFTIME), intent(in):: diff1, diff2
type(DC_DIFFTIME), intent(inout):: rslt
continue
if ( ( diff1 % nondim_flag .and. .not. diff2 % nondim_flag ) &
& .or. ( .not. diff1 % nondim_flag .and. diff2 % nondim_flag ) ) then
call StoreError(DC_EDIMTIME, opr)
end if
rslt % nondim_flag = diff1 % nondim_flag
end subroutine dcdate_nondimcheck
function ParseTimeUnits(str) result(symbol)
!
! 引数 *str* に与えられた文字列を解釈し, 日時の単位を示す
! シンボルを返します. それぞれ以下の文字列が日時の単位として解釈されます.
! 大文字と小文字は区別されません.
!
! 年 :: dc_date_types#UNIT_YEAR
! 月 :: dc_date_types#UNIT_MONTH
! 日 :: dc_date_types#UNIT_DAY
! 時 :: dc_date_types#UNIT_HOUR
! 分 :: dc_date_types#UNIT_MIN
! 秒 :: dc_date_types#UNIT_SEC
! 無次元時間 :: dc_date_types#UNIT_NONDIM
!
! 返るシンボル (整数型) は以下の通りです.
!
! 年 :: dc_date_types#UNIT_SYMBOL_YEAR
! 月 :: dc_date_types#UNIT_SYMBOL_MONTH
! 日 :: dc_date_types#UNIT_SYMBOL_DAY
! 時 :: dc_date_types#UNIT_SYMBOL_HOUR
! 分 :: dc_date_types#UNIT_SYMBOL_MIN
! 秒 :: dc_date_types#UNIT_SYMBOL_SEC
! 無次元時間 :: dc_date_types#UNIT_SYMBOL_NONDIM
!
! これらに該当しない文字列を *str* に与えた場合,
! dc_date_types#UNIT_SYMBOL_ERR が返ります.
!
use dc_types, only: TOKEN
use dc_date_types, only: UNIT_YEAR, UNIT_MONTH, UNIT_DAY, &
& UNIT_HOUR, UNIT_MIN, UNIT_SEC, UNIT_NONDIM, &
& UNIT_SYMBOL_YEAR, UNIT_SYMBOL_MONTH, UNIT_SYMBOL_DAY, &
& UNIT_SYMBOL_HOUR, UNIT_SYMBOL_MIN, UNIT_SYMBOL_SEC, &
& UNIT_SYMBOL_NONDIM, UNIT_SYMBOL_ERR
use dc_string, only: StriEq
implicit none
character(*), intent(in):: str
integer:: symbol
integer:: unit_str_size, i
character(TOKEN):: unit
continue
unit = adjustl(str)
unit_str_size = size(UNIT_NONDIM)
do i = 1, unit_str_size
if (StriEq(trim(unit), trim(UNIT_NONDIM(i)))) then
symbol = UNIT_SYMBOL_NONDIM
return
end if
end do
unit_str_size = size(UNIT_SEC)
do i = 1, unit_str_size
if (StriEq(trim(unit), trim(UNIT_SEC(i)))) then
symbol = UNIT_SYMBOL_SEC
return
end if
end do
unit_str_size = size(UNIT_MIN)
do i = 1, unit_str_size
if (StriEq(trim(unit), trim(UNIT_MIN(i)))) then
symbol = UNIT_SYMBOL_MIN
return
end if
end do
unit_str_size = size(UNIT_HOUR)
do i = 1, unit_str_size
if (StriEq(trim(unit), trim(UNIT_HOUR(i)))) then
symbol = UNIT_SYMBOL_HOUR
return
end if
end do
unit_str_size = size(UNIT_DAY)
do i = 1, unit_str_size
if (StriEq(trim(unit), trim(UNIT_DAY(i)))) then
symbol = UNIT_SYMBOL_DAY
return
end if
end do
unit_str_size = size(UNIT_MONTH)
do i = 1, unit_str_size
if (StriEq(trim(unit), trim(UNIT_MONTH(i)))) then
symbol = UNIT_SYMBOL_MONTH
return
end if
end do
unit_str_size = size(UNIT_YEAR)
do i = 1, unit_str_size
if (StriEq(trim(unit), trim(UNIT_YEAR(i)))) then
symbol = UNIT_SYMBOL_YEAR
return
end if
end do
symbol = UNIT_SYMBOL_ERR
end function ParseTimeUnits
character(TOKEN) function dcdate_parse_unit(str) result(unit)
!
! このサブルーチンは内部向けなので dc_date モジュール外では
! 極力使用しないでください.
!
! 引数 *str* に与えられた文字列を解釈し, 日時の単位を
! 返します. それぞれ以下の文字列が日時の単位として解釈されます.
! 大文字と小文字は区別されません.
! 返る文字列は以下の文字型の配列の先頭の文字列です.
! (例: *str* に 'hrs.' が与えられる場合, dc_date_types#UNIT_HOUR
! 配列の先頭の文字列 UNIT_HOUR(1) が返ります.)
!
! 年 :: dc_date_types#UNIT_YEAR
! 月 :: dc_date_types#UNIT_MONTH
! 日 :: dc_date_types#UNIT_DAY
! 時 :: dc_date_types#UNIT_HOUR
! 分 :: dc_date_types#UNIT_MIN
! 秒 :: dc_date_types#UNIT_SEC
! 無次元時間 :: dc_date_types#UNIT_NONDIM
!
! これらに該当しない文字列を *str* に与えた場合, 空文字が返ります.
!
use dc_types, only: TOKEN
use dc_date_types, only: UNIT_YEAR, UNIT_MONTH, UNIT_DAY, &
& UNIT_HOUR, UNIT_MIN, UNIT_SEC, UNIT_NONDIM
use dc_string, only: StriEq
implicit none
character(*), intent(in):: str
integer :: unit_str_size, i
continue
unit = adjustl(str)
unit_str_size = size(UNIT_NONDIM)
do i = 1, unit_str_size
if (StriEq(trim(unit), trim(UNIT_NONDIM(i)))) then
unit = UNIT_NONDIM(1)
return
end if
end do
unit_str_size = size(UNIT_SEC)
do i = 1, unit_str_size
if (StriEq(trim(unit), trim(UNIT_SEC(i)))) then
unit = UNIT_SEC(1)
return
end if
end do
unit_str_size = size(UNIT_MIN)
do i = 1, unit_str_size
if (StriEq(trim(unit), trim(UNIT_MIN(i)))) then
unit = UNIT_MIN(1)
return
end if
end do
unit_str_size = size(UNIT_HOUR)
do i = 1, unit_str_size
if (StriEq(trim(unit), trim(UNIT_HOUR(i)))) then
unit = UNIT_HOUR(1)
return
end if
end do
unit_str_size = size(UNIT_DAY)
do i = 1, unit_str_size
if (StriEq(trim(unit), trim(UNIT_DAY(i)))) then
unit = UNIT_DAY(1)
return
end if
end do
unit_str_size = size(UNIT_MONTH)
do i = 1, unit_str_size
if (StriEq(trim(unit), trim(UNIT_MONTH(i)))) then
unit = UNIT_MONTH(1)
return
end if
end do
unit_str_size = size(UNIT_YEAR)
do i = 1, unit_str_size
if (StriEq(trim(unit), trim(UNIT_YEAR(i)))) then
unit = UNIT_YEAR(1)
return
end if
end do
unit = ''
end function dcdate_parse_unit
type(DC_DATETIME) function dcdate_add_ft(diff, time) result(result)
!
! 2 つの日時 (DC_DATETIME 型) もしくは
! 日時差 (DC_DIFFTIME 型)の加算を行います.
!
use dc_scaledsec, only: DC_SCALED_SEC, &
& operator(<), operator(>), operator(<=), operator(>=), &
& operator(+), operator(-), operator(*), operator(/), &
& modulo, int, abs, sign
implicit none
type(DC_DIFFTIME), intent(in):: diff
type(DC_DATETIME), intent(in):: time
type(DC_SCALED_SEC):: time_year, time_mon, time_day, time_sec
integer:: time_caltype
character(6):: time_zone
continue
call Eval(time, sclyear = time_year, sclmon = time_mon, &
& sclday = time_day, sclsec = time_sec, &
& caltype = time_caltype, zone = time_zone)
call DCDateTimeCreate(result, sclyear = time_year, &
& sclmon = time_mon + diff % mon, &
& sclday = time_day + diff % day, &
& sclsec = time_sec + diff % sec, &
& caltype = time_caltype, zone = time_zone)
end function dcdate_add_ft
type(DC_DATETIME) function dcdate_add_tf(time, diff) result(result)
use dc_types, only: DP
implicit none
type(DC_DATETIME), intent(in):: time
type(DC_DIFFTIME), intent(in):: diff
continue
result = dcdate_add_ft(diff, time)
end function dcdate_add_tf
type(DC_DIFFTIME) function dcdate_add_ff(diff1, diff2) result(result)
use dc_scaledsec, only: operator(+)
implicit none
type(DC_DIFFTIME), intent(in):: diff1, diff2
continue
result % mon = diff1 % mon + diff2 % mon
result % day = diff1 % day + diff2 % day
result % sec = diff1 % sec + diff2 % sec
result % day_seconds = diff1 % day_seconds
call dcdate_nondimcheck('dc_date#operator(+)', diff1, diff2, result)
call dcdate_normalize(result % day, result % sec, result % day_seconds, result % nondim_flag)
end function dcdate_add_ff
type(DC_DIFFTIME) function dcdate_add_fd(diff, sec) result(result)
use dc_scaledsec, only: operator(+)
implicit none
type(DC_DIFFTIME), intent(in):: diff
real(DP), intent(in):: sec
continue
result % mon = diff % mon
result % day = diff % day
result % sec = diff % sec + sec
result % day_seconds = diff % day_seconds
call dcdate_normalize(result % day, result % sec, result % day_seconds, result % nondim_flag)
end function dcdate_add_fd
type(DC_DIFFTIME) function dcdate_add_fr(diff, sec) result(result)
use dc_scaledsec, only: operator(+)
implicit none
type(DC_DIFFTIME), intent(in):: diff
real, intent(in):: sec
continue
result = diff + real( sec, DP )
end function dcdate_add_fr
type(DC_DIFFTIME) function dcdate_add_fi(diff, sec) result(result)
use dc_scaledsec, only: operator(+)
implicit none
type(DC_DIFFTIME), intent(in):: diff
integer, intent(in):: sec
continue
result = diff + real( sec, DP )
end function dcdate_add_fi
type(DC_DATETIME) function dcdate_sub_tf(time, diff) result(result)
!
! 2 つの日時 (DC_DATETIME 型) もしくは
! 日時差 (DC_DIFFTIME 型)の減算を行います.
!
use dc_scaledsec, only: DC_SCALED_SEC, &
& operator(<), operator(>), operator(<=), operator(>=), &
& operator(+), operator(-), operator(*), operator(/), &
& modulo, int, abs, sign
implicit none
type(DC_DATETIME), intent(in):: time
type(DC_DIFFTIME), intent(in):: diff
type(DC_SCALED_SEC):: time_year, time_mon, time_day, time_sec
integer:: time_caltype
character(6):: time_zone
continue
call Eval(time, &
& sclyear = time_year, sclmon = time_mon, sclday = time_day, &
& sclsec = time_sec, caltype = time_caltype, zone = time_zone )
call DCDateTimeCreate(result, &
& sclyear = time_year, &
& sclmon = time_mon - diff % mon, &
& sclday = time_day - diff % day, &
& sclsec = time_sec - diff % sec, &
& caltype = time_caltype, zone = time_zone)
end function dcdate_sub_tf
type(DC_DIFFTIME) function dcdate_sub_tt(time1, time2) result(result)
use dc_scaledsec, only: DC_SCALED_SEC, &
& operator(<), operator(>), operator(<=), operator(>=), &
& operator(+), operator(-), operator(*), operator(/), &
& modulo, int, abs, sign
implicit none
type(DC_DATETIME), intent(in):: time1, time2
continue
result % day = time1 % day - time2 % day
result % sec = time1 % sec - time2 % sec &
& + EvalSclSec(ZoneToDiff(time1 % zone) - ZoneToDiff(time2 % zone))
result % day_seconds = time1 % day_seconds
call dcdate_normalize(result % day, result % sec, result % day_seconds, result % nondim_flag)
end function dcdate_sub_tt
type(DC_DIFFTIME) function dcdate_sub_ff(diff1, diff2) result(result)
use dc_scaledsec, only: operator(-)
implicit none
type(DC_DIFFTIME), intent(in):: diff1, diff2
continue
result % mon = diff1 % mon - diff2 % mon
result % day = diff1 % day - diff2 % day
result % sec = diff1 % sec - diff2 % sec
result % day_seconds = diff1 % day_seconds
call dcdate_nondimcheck('dc_date#operator(-)', diff1, diff2, result)
call dcdate_normalize(result % day, result % sec, result % day_seconds, result % nondim_flag)
end function dcdate_sub_ff
type(DC_DIFFTIME) function dcdate_sub_fd(diff, sec) result(result)
use dc_scaledsec, only: operator(-)
implicit none
type(DC_DIFFTIME), intent(in):: diff
real(DP), intent(in):: sec
continue
result % mon = diff % mon
result % day = diff % day
result % sec = diff % sec - sec
result % day_seconds = diff % day_seconds
call dcdate_normalize(result % day, result % sec, result % day_seconds, result % nondim_flag)
end function dcdate_sub_fd
type(DC_DIFFTIME) function dcdate_sub_fr(diff, sec) result(result)
use dc_scaledsec, only: operator(-)
implicit none
type(DC_DIFFTIME), intent(in):: diff
real, intent(in):: sec
continue
result = diff - real( sec, DP )
end function dcdate_sub_fr
type(DC_DIFFTIME) function dcdate_sub_fi(diff, sec) result(result)
use dc_scaledsec, only: operator(-)
implicit none
type(DC_DIFFTIME), intent(in):: diff
integer, intent(in):: sec
continue
result = diff - real( sec, DP )
end function dcdate_sub_fi
type(DC_DIFFTIME) function dcdate_mul_if(factor, diff) result(result)
!
! 日時差 *diff* と *facter* とを乗算した結果を返します.
!
use dc_scaledsec, only: operator(*)
implicit none
integer, intent(in):: factor
type(DC_DIFFTIME), intent(in):: diff
continue
result % mon = factor * diff % mon
result % day = factor * diff % day
result % sec = factor * diff % sec
result % day_seconds = diff % day_seconds
result % nondim_flag = diff % nondim_flag
call dcdate_normalize(result % day, result % sec, result % day_seconds, result % nondim_flag)
end function dcdate_mul_if
type(DC_DIFFTIME) function dcdate_mul_fi(diff, factor) result(result)
implicit none
type(DC_DIFFTIME), intent(in):: diff
integer, intent(in):: factor
continue
result = dcdate_mul_if(factor, diff)
end function dcdate_mul_fi
type(DC_DIFFTIME) function dcdate_mul_rf(factor, diff) result(result)
!
! ※ 注意 : 月差を非整数倍すると近似的結果になるおそれがあります
use dc_types, only: DP
implicit none
real, intent(in):: factor
type(DC_DIFFTIME), intent(in):: diff
continue
result = dcdate_mul_df(real(factor, DP), diff)
end function dcdate_mul_rf
type(DC_DIFFTIME) function dcdate_mul_fr(diff, factor) result(result)
!
! ※ 注意 : 月差を非整数倍すると近似的結果になるおそれがあります
implicit none
type(DC_DIFFTIME), intent(in):: diff
real, intent(in):: factor
continue
result = dcdate_mul_rf(factor, diff)
end function dcdate_mul_fr
type(DC_DIFFTIME) function dcdate_mul_df(factor, diff) result(result)
!
! ※ 注意 : 月差を非整数倍すると近似的結果になるおそれがあります
use dc_types, only: DP
use dc_date_types, only: CYCLIC_MDAYS
use dc_scaledsec, only: DC_SCALED_SEC, &
& operator(<), operator(>), operator(<=), operator(>=), &
& operator(+), operator(-), operator(*), operator(/), &
& modulo, int, abs, sign
implicit none
real(DP), intent(in):: factor
type(DC_DIFFTIME), intent(in):: diff
type(DC_SCALED_SEC):: month, day
continue
month = factor * diff % mon
result % mon = int(month)
day = factor * diff % day + int(CYCLIC_MDAYS * (month - result % mon))
result % day = int(day)
result % sec = &
& factor * diff % sec + (day - result % day) * diff % day_seconds
result % day_seconds = diff % day_seconds
result % nondim_flag = diff % nondim_flag
call dcdate_normalize(result % day, result % sec, result % day_seconds, result % nondim_flag)
end function dcdate_mul_df
type(DC_DIFFTIME) function dcdate_mul_fd(diff, factor) result(result)
!
! ※ 注意 : 月差を非整数倍すると近似的結果になるおそれがあります
use dc_types, only: DP
implicit none
type(DC_DIFFTIME), intent(in):: diff
real(DP), intent(in):: factor
continue
result = dcdate_mul_df(factor, diff)
end function dcdate_mul_fd
type(DC_DIFFTIME) function dcdate_div_fi(diff, denominator) result(result)
!
! 日時差 *diff* を *denominator* で除算した結果を返します.
!
! ※ 注意 : 月差を除算すると近似的結果になるおそれがあります
use dc_date_types, only: CYCLIC_MDAYS
use dc_scaledsec, only: DC_SCALED_SEC, &
& operator(<), operator(>), operator(<=), operator(>=), &
& operator(+), operator(-), operator(*), operator(/), &
& modscl => mod, modulo, int, abs, sign
implicit none
type(DC_DIFFTIME), intent(in):: diff
integer, intent(in):: denominator
continue
result % mon = int( diff % mon / denominator )
! 月からの近似的繰り下がりは日単位でしか行わない
result % day = &
& int( diff % day / denominator ) &
& + int( (CYCLIC_MDAYS * modscl(diff % mon, denominator)) / denominator )
result % sec = diff % sec / denominator + &
& (diff % day_seconds * modscl(diff % day, denominator)) / &
& denominator
result % nondim_flag = diff % nondim_flag
end function dcdate_div_fi
type(DC_DIFFTIME) function dcdate_div_fr(diff, denominator) result(result)
!
! ※ 注意 : 月差を除算すると近似的結果になるおそれがあります
implicit none
type(DC_DIFFTIME), intent(in):: diff
real, intent(in):: denominator
continue
result = dcdate_div_fd(diff, real(denominator, DP))
end function dcdate_div_fr
type(DC_DIFFTIME) function dcdate_div_fd(diff, denominator) result(result)
!
! ※ 注意 : 月差を除算すると近似的結果になるおそれがあります
use dc_date_types, only: CYCLIC_MDAYS
use dc_scaledsec, only: DC_SCALED_SEC, &
& operator(<), operator(>), operator(<=), operator(>=), &
& operator(+), operator(-), operator(*), operator(/), &
& modulo, int, abs, sign
implicit none
type(DC_DIFFTIME), intent(in):: diff
real(DP), intent(in):: denominator
type(DC_SCALED_SEC):: month, day
continue
month = int( diff % mon / denominator )
result % mon = int(month)
day = int( diff % day / denominator ) &
& + int(CYCLIC_MDAYS * (month - result % mon))
result % day = int(day)
result % sec = &
& diff % sec / denominator + (day - result % day) * diff % day_seconds
result % day_seconds = diff % day_seconds
result % nondim_flag = diff % nondim_flag
call dcdate_normalize(result % day, result % sec, result % day_seconds, result % nondim_flag)
end function dcdate_div_fd
real(DP) function dcdate_div_ff(diff1, diff2) result(result)
!
! ※ 注意 : 月差と日時の混在する除算は近似的結果になるおそれがあります
use dc_date_types, only: CYCLIC_MDAYS
use dc_scaledsec, only: DC_SCALED_SEC, assignment(=), &
& operator(<), operator(>), operator(<=), operator(>=), &
& operator(+), operator(-), operator(*), operator(/), &
& modulo, int, abs, sign
implicit none
type(DC_DIFFTIME), intent(in):: diff1, diff2
continue
! ゼロ割対応コードが必要か?
result = &
& (diff1 % day_seconds * (CYCLIC_MDAYS * diff1 % mon + diff1 % day) &
& + diff1 % sec) / &
& (diff2 % day_seconds * (CYCLIC_MDAYS * diff2 % mon + diff2 % day) &
& + diff2 % sec)
end function dcdate_div_ff
type(DC_DIFFTIME) function dcdate_mod_ff(diff1, diff2) result(result)
!
! 引数 diff1 を diff2 で除算した際の余りを返します.
!
! ※ 注意: 月差と日時の混在する除算は近似的結果になるおそれがあります
!
use dc_date_types, only: CYCLIC_MDAYS
use dc_scaledsec, only: DC_SCALED_SEC, &
& operator(==), operator(<), operator(>), operator(<=), operator(>=), &
& operator(+), operator(-), operator(*), operator(/), &
& modscl => mod, modulo, int, abs, sign
implicit none
type(DC_DIFFTIME), intent(in):: diff1, diff2
type(DC_SCALED_SEC):: sec1, sec2
type(DC_SCALED_SEC):: zero_sec
continue
result % day_seconds = diff1 % day_seconds
if (diff1 % day == zero_sec .and. diff2 % day == zero_sec .and. &
& diff1 % sec == zero_sec .and. diff2 % sec == zero_sec) then
result % mon = modscl(diff1 % mon, diff2 % mon)
result % day = zero_sec
result % sec = zero_sec
else if (diff1 % sec == zero_sec .and. diff2 % sec == zero_sec) then
result % mon = zero_sec
result % day = modscl((CYCLIC_MDAYS * diff1 % mon + diff1 % day), &
& (CYCLIC_MDAYS * diff2 % mon + diff2 % day))
result % sec = zero_sec
else
sec1 = diff1 % day_seconds * (CYCLIC_MDAYS * diff1 % mon + diff1 % day) &
& + diff1 % sec
sec2 = diff2 % day_seconds * (CYCLIC_MDAYS * diff2 % mon + diff2 % day) &
& + diff2 % sec
result % sec = modscl(sec1, sec2)
result % day = zero_sec
result % mon = zero_sec
call dcdate_normalize(result % day, result % sec, result % day_seconds, result % nondim_flag)
endif
call dcdate_nondimcheck('dc_date#mod', diff1, diff2, result)
end function dcdate_mod_ff
logical function dcdate_gt_tt(time1, time2) result(result)
!
! 2 つの引数の日時を比較します.
! 1 つ目の引数に格納される日時が 2 つ目の引数に格納される日時
! よりも進んでいる場合, .true. が返ります.
!
use dc_scaledsec, only: DC_SCALED_SEC, &
& operator(==), operator(<), operator(>), operator(<=), operator(>=), &
& operator(+), operator(-), operator(*), operator(/), &
& modulo, int, abs, sign
implicit none
type(DC_DATETIME), intent(in):: time1, time2
type(DC_SCALED_SEC):: year1, year2, time1_sec, time2_sec
continue
call Eval(time1, sclyear=year1)
call Eval(time2, sclyear=year2)
if (year1 > year2) then
result = .true.
elseif (year1 < year2) then
result = .false.
else
time1_sec = EvalSclSec(time1) + EvalSclSec(ZoneToDiff(time1 % zone))
time2_sec = EvalSclSec(time2) + EvalSclSec(ZoneToDiff(time2 % zone))
if (time1_sec > time2_sec) then
result = .true.
else
result = .false.
end if
end if
end function dcdate_gt_tt
logical function dcdate_gt_ff(diff1, diff2) result(result)
!
! 2 つの引数の日時差を比較します.
! 1 つ目の引数に格納される日時差が 2 つ目の引数に格納される日時差
! よりも大きい場合, .true. が返ります.
!
use dc_date_types, only: CYCLIC_MDAYS
use dc_scaledsec, only: &
& operator(<), operator(>), operator(<=), operator(>=), operator(==)
implicit none
type(DC_DIFFTIME), intent(in):: diff1, diff2
continue
if ( diff1 % day_seconds == diff2 % day_seconds ) then
if ( diff1 % mon > diff2 % mon ) then
result = .true. ; return
elseif ( diff1 % mon < diff2 % mon ) then
result = .false. ; return
end if
if ( diff1 % day > diff2 % day ) then
result = .true. ; return
elseif ( diff1 % day < diff2 % day ) then
result = .false. ; return
end if
if ( diff1 % sec > diff2 % sec ) then
result = .true. ; return
elseif ( diff1 % sec < diff2 % sec ) then
result = .false. ; return
end if
result = .false.
else
if (EvalSec(diff1) > EvalSec(diff2)) then
result = .true.
else
result = .false.
end if
end if
end function dcdate_gt_ff
logical function dcdate_gt_fi(diff, factor) result(result)
!
! 2 つの引数の日時差を比較します.
! 1 つ目の引数に格納される日時差が 2 つ目の引数に格納される日時差
! よりも大きい場合, .true. が返ります.
!
use dc_date_types, only: CYCLIC_MDAYS
implicit none
type(DC_DIFFTIME), intent(in):: diff
integer, intent(in):: factor
continue
result = EvalSec(diff) > factor
end function dcdate_gt_fi
logical function dcdate_gt_if(factor, diff) result(result)
!
! 2 つの引数の日時差を比較します.
! 1 つ目の引数に格納される日時差が 2 つ目の引数に格納される日時差
! よりも大きい場合, .true. が返ります.
!
use dc_date_types, only: CYCLIC_MDAYS
implicit none
integer, intent(in):: factor
type(DC_DIFFTIME), intent(in):: diff
continue
result = factor > EvalSec(diff)
end function dcdate_gt_if
logical function dcdate_lt_tt(time1, time2) result(result)
!
! 2 つの引数の日時を比較します.
! 2 つ目の引数に格納される日時が 1 つ目の引数に格納される日時
! よりも進んでいる場合, .true. が返ります.
!
use dc_scaledsec, only: DC_SCALED_SEC, &
& operator(==), operator(<), operator(>), operator(<=), operator(>=), &
& operator(+)
implicit none
type(DC_DATETIME), intent(in):: time1, time2
type(DC_SCALED_SEC):: year1, year2, time1_sec, time2_sec
continue
call Eval(time1, sclyear=year1)
call Eval(time2, sclyear=year2)
if (year1 < year2) then
result = .true.
elseif (year1 > year2) then
result = .false.
else
time1_sec = EvalSclSec(time1) + EvalSclSec(ZoneToDiff(time1 % zone))
time2_sec = EvalSclSec(time2) + EvalSclSec(ZoneToDiff(time2 % zone))
if (time1_sec < time2_sec) then
result = .true.
else
result = .false.
end if
end if
end function dcdate_lt_tt
logical function dcdate_lt_ff(diff1, diff2) result(result)
!
! 2 つの引数の日時差を比較します.
! 2 つ目の引数に格納される日時差が 1 つ目の引数に格納される日時差
! よりも大きい場合, .true. が返ります.
!
use dc_scaledsec, only: &
& operator(<), operator(>), operator(<=), operator(>=), operator(==)
implicit none
type(DC_DIFFTIME), intent(in):: diff1, diff2
continue
if ( diff1 % day_seconds == diff2 % day_seconds ) then
if ( diff1 % mon < diff2 % mon ) then
result = .true. ; return
elseif ( diff1 % mon > diff2 % mon ) then
result = .false. ; return
end if
if ( diff1 % day < diff2 % day ) then
result = .true. ; return
elseif ( diff1 % day > diff2 % day ) then
result = .false. ; return
end if
if ( diff1 % sec < diff2 % sec ) then
result = .true. ; return
elseif ( diff1 % sec > diff2 % sec ) then
result = .false. ; return
end if
result = .false.
else
if (EvalSec(diff1) < EvalSec(diff2)) then
result = .true.
else
result = .false.
end if
end if
end function dcdate_lt_ff
logical function dcdate_lt_fi(diff, factor) result(result)
!
! 2 つの引数の日時差を比較します.
! 2 つ目の引数に格納される日時差が 1 つ目の引数に格納される日時差
! よりも大きい場合, .true. が返ります.
!
implicit none
type(DC_DIFFTIME), intent(in):: diff
integer, intent(in):: factor
continue
result = EvalSec(diff) < factor
end function dcdate_lt_fi
logical function dcdate_lt_if(factor, diff) result(result)
!
! 2 つの引数の日時差を比較します.
! 2 つ目の引数に格納される日時差が 1 つ目の引数に格納される日時差
! よりも大きい場合, .true. が返ります.
!
implicit none
integer, intent(in):: factor
type(DC_DIFFTIME), intent(in):: diff
continue
result = factor < EvalSec(diff)
end function dcdate_lt_if
logical function dcdate_ge_tt(time1, time2) result(result)
!
! 2 つの引数の日時を比較します.
! 1 つ目の引数に格納される日時が 2 つ目の引数に格納される日時
! よりも進んでいる場合かもしくは等しい場合, .true. が返ります.
!
implicit none
type(DC_DATETIME), intent(in):: time1, time2
continue
result = .not. time1 < time2
end function dcdate_ge_tt
logical function dcdate_ge_ff(diff1, diff2) result(result)
!
! 2 つの引数の日時差を比較します.
! 1 つ目の引数に格納される日時差が 2 つ目の引数に格納される日時差
! よりも大きい場合かもしくは等しい場合, .true. が返ります.
!
implicit none
type(DC_DIFFTIME), intent(in):: diff1, diff2
continue
result = .not. diff1 < diff2
end function dcdate_ge_ff
logical function dcdate_ge_fi(diff, factor) result(result)
!
! 2 つの引数の日時差を比較します.
! 1 つ目の引数に格納される日時差が 2 つ目の引数に格納される日時差
! よりも大きい場合かもしくは等しい場合, .true. が返ります.
!
implicit none
type(DC_DIFFTIME), intent(in):: diff
integer, intent(in):: factor
continue
result = .not. diff < factor
end function dcdate_ge_fi
logical function dcdate_ge_if(factor, diff) result(result)
!
! 2 つの引数の日時差を比較します.
! 1 つ目の引数に格納される日時差が 2 つ目の引数に格納される日時差
! よりも大きい場合かもしくは等しい場合, .true. が返ります.
!
implicit none
integer, intent(in):: factor
type(DC_DIFFTIME), intent(in):: diff
continue
result = .not. factor < diff
end function dcdate_ge_if
logical function dcdate_le_tt(time1, time2) result(result)
!
! 2 つの引数の日時を比較します.
! 2 つ目の引数に格納される日時が 1 つ目の引数に格納される日時
! よりも進んでいる場合かもしくは等しい場合, .true. が返ります.
!
implicit none
type(DC_DATETIME), intent(in):: time1, time2
continue
result = .not. time1 > time2
end function dcdate_le_tt
logical function dcdate_le_ff(diff1, diff2) result(result)
!
! 2 つの引数の日時差を比較します.
! 2 つ目の引数に格納される日時差が 1 つ目の引数に格納される日時差
! よりも大きい場合かもしくは等しい場合, .true. が返ります.
!
implicit none
type(DC_DIFFTIME), intent(in):: diff1, diff2
continue
result = .not. diff1 > diff2
end function dcdate_le_ff
logical function dcdate_le_fi(diff, factor) result(result)
!
! 2 つの引数の日時差を比較します.
! 2 つ目の引数に格納される日時差が 1 つ目の引数に格納される日時差
! よりも大きい場合かもしくは等しい場合, .true. が返ります.
!
implicit none
type(DC_DIFFTIME), intent(in):: diff
integer, intent(in):: factor
continue
result = .not. diff > factor
end function dcdate_le_fi
logical function dcdate_le_if(factor, diff) result(result)
!
! 2 つの引数の日時差を比較します.
! 2 つ目の引数に格納される日時差が 1 つ目の引数に格納される日時差
! よりも大きい場合かもしくは等しい場合, .true. が返ります.
!
implicit none
integer, intent(in):: factor
type(DC_DIFFTIME), intent(in):: diff
continue
result = .not. factor > diff
end function dcdate_le_if
type(DC_DATETIME) function dcdate_max_tt(time1, time2) result(result)
!
! 2 つの引数の日時を比較し, より日時が進んでいる方を返します.
!
implicit none
type(DC_DATETIME), intent(in):: time1, time2
continue
if ( time1 > time2 ) then
result = time1
else
result = time2
end if
end function dcdate_max_tt
type(DC_DIFFTIME) function dcdate_max_ff(diff1, diff2) result(result)
!
! 2 つの引数の日時差を比較し, より大きい方を返します.
!
implicit none
type(DC_DIFFTIME), intent(in):: diff1, diff2
continue
if ( diff1 > diff2 ) then
result = diff1
else
result = diff2
end if
call dcdate_nondimcheck('dc_date#max', diff1, diff2, result)
end function dcdate_max_ff
type(DC_DATETIME) function dcdate_min_tt(time1, time2) result(result)
!
! 2 つの引数の日時を比較し, より日時が遅れている方を返します.
!
implicit none
type(DC_DATETIME), intent(in):: time1, time2
continue
if ( time1 < time2 ) then
result = time1
else
result = time2
end if
end function dcdate_min_tt
type(DC_DIFFTIME) function dcdate_min_ff(diff1, diff2) result(result)
!
! 2 つの引数の日時差を比較し, より小さい方を返します.
!
implicit none
type(DC_DIFFTIME), intent(in):: diff1, diff2
continue
if ( diff1 < diff2 ) then
result = diff1
else
result = diff2
end if
call dcdate_nondimcheck('dc_date#min', diff1, diff2, result)
end function dcdate_min_ff
logical function dcdate_eq_tt(time1, time2) result(result)
!
! 2 つの引数の日時を比較します.
! 1 つ目の引数に格納される日時が 2 つ目の引数に格納される日時
! と同じ場合, .true. が返ります.
!
use dc_scaledsec, only: DC_SCALED_SEC, &
& operator(==), operator(<), operator(>), operator(<=), operator(>=), &
& operator(+), assignment(=)
implicit none
type(DC_DATETIME), intent(in):: time1, time2
type(DC_SCALED_SEC):: year1, year2, time1_sec, time2_sec
continue
call Eval(time1, sclyear=year1)
call Eval(time2, sclyear=year2)
time1_sec = EvalSclSec(time1) + EvalSclSec(ZoneToDiff(time1 % zone))
time2_sec = EvalSclSec(time2) + EvalSclSec(ZoneToDiff(time2 % zone))
if (year1 == year2 .and. time1_sec == time2_sec) then
result = .true.
else
result = .false.
end if
end function dcdate_eq_tt
logical function dcdate_eq_ff(diff1, diff2) result(result)
!
! 2 つの引数の日時差を比較します.
! 1 つ目の引数に格納される日時差が 2 つ目の引数に格納される日時差
! と同じ場合, .true. が返ります.
!
use dc_scaledsec, only: operator(==)
implicit none
type(DC_DIFFTIME), intent(in):: diff1, diff2
continue
if ( diff1 % mon == diff2 % mon &
& .and. diff1 % day == diff2 % day &
& .and. diff1 % sec == diff2 % sec ) then
result = .true.
else
result = .false.
end if
end function dcdate_eq_ff
logical function dcdate_eq_if(i, diff) result(result)
!
! 引数 *diff* の日時差が *i* と等しいかどうかを比較します. *diff*
! を秒数に換算した値と *i* とが等しい場合, .true. が返ります.
!
implicit none
type(DC_DIFFTIME), intent(in):: diff
integer, intent(in):: i
continue
result = dcdate_eq_rf(real(i), diff)
end function dcdate_eq_if
logical function dcdate_eq_fi(diff, i) result(result)
implicit none
type(DC_DIFFTIME), intent(in):: diff
integer, intent(in):: i
continue
result = dcdate_eq_if(i, diff)
end function dcdate_eq_fi
logical function dcdate_eq_rf(r, diff) result(result)
!
! 引数 *diff* の日時差が *r* と等しいかどうかを比較します. *diff*
! を秒数に換算した値と *r* とが等しい場合, .true. が返ります.
!
use dc_scaledsec, only: operator(==)
implicit none
type(DC_DIFFTIME), intent(in):: diff
real, intent(in):: r
continue
if (EvalSclSec(diff) == r) then
result = .true.
else
result = .false.
end if
end function dcdate_eq_rf
logical function dcdate_eq_fr(diff, r) result(result)
implicit none
type(DC_DIFFTIME), intent(in):: diff
real, intent(in):: r
continue
result = dcdate_eq_rf(r, diff)
end function dcdate_eq_fr
logical function dcdate_eq_df(d, diff) result(result)
!
! 引数 *diff* の日時差が *d* と等しいかどうかを比較します. *diff*
! を秒数に換算した値と *d* とが等しい場合, .true. が返ります.
!
use dc_types, only: DP
use dc_scaledsec, only: operator(==)
implicit none
type(DC_DIFFTIME), intent(in):: diff
real(DP), intent(in):: d
continue
if (EvalSclSec(diff) == d) then
result = .true.
else
result = .false.
end if
end function dcdate_eq_df
logical function dcdate_eq_fd(diff, d) result(result)
use dc_types, only: DP
implicit none
type(DC_DIFFTIME), intent(in):: diff
real(DP), intent(in):: d
continue
result = dcdate_eq_df(d, diff)
end function dcdate_eq_fd
end module dc_date