Class | dc_date |
In: |
dc_date.f90
|
日付および時刻を扱うための手続きを提供します
以下の手続きは構造型 dc_date_types#DC_DATETIME または dc_date_types#DC_DIFFTIME 変数 (日時, 時刻に関する情報を格納) を対象とします.
Create : | 初期化 |
assignment(=) : | 〃 |
Eval : | 日時, 時刻情報を個別に取得 |
toChar : | 日時, 時刻情報を文字型変数へ変換 |
Put_Line : | 日時, 時刻情報の印字を行います. |
EvalDay : | 日数 (実数型) に換算して取得 |
EvalHour : | 時間 (実数型) に換算して取得 |
EvalMin : | 分 (実数型) に換算して取得 |
EvalSec : | 秒 (実数型) に換算して取得 |
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 型同士) |
SetZone : | タイムゾーンを変更 |
以下の手続きは dc_date_types 内部の変数を変更します.
SetCaltype : | 暦法のデフォルトを変更 |
SetSecOfDay : | 1 日の秒数のデフォルトを変更 |
その他の手続き
ValidCaltype : | 暦法が有効なものかをチェック |
ValidZone : | タイムゾーンとして有効化をチェック |
ZoneToDiff : | タイムゾーンを dc_date_types#DC_DIFFTIME 変数へと変換 |
DC_DATETIME 型の変数にサブルーチン Create を用いると, 時刻が設定されます. 下記のように特に年月日を指定しないと現在時刻が設定されます. 設定された時刻は 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: Create, toChar implicit none type(DC_DATETIME) :: time call Create(time) 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_string, only: Printf use dc_date_types, only: DC_DATETIME, DC_DIFFTIME use dc_date, only: Create, toChar, operator(+) implicit none type(DC_DATETIME) :: time_before, time_after type(DC_DIFFTIME) :: diff call Create(time_before, year=2006, mon=6, day=10, hour=14, min=15, sec=0.0d0) call Create(diff, day=25, hour=12, min=50) 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: Create, 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 Create(DelTimef, 5.0d0, 'sec') ! Δt = 5.0 (秒) call Create(intervalf, 1.0d0, 'min') ! データの出力間隔 = 1.0 (分) call Create(calctimef, 12.0d0, 'hour') ! 計算時間 = 12.0 (時間) 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
Subroutine : | |||
time : | type(DC_DATETIME), intent(out) | ||
year : | integer, intent(in), optional
| ||
mon : | integer, intent(in), optional
| ||
day : | integer, intent(in), optional
| ||
hour : | integer, intent(in), optional
| ||
min : | integer, intent(in), optional
| ||
sec : | real(DP),intent(in), optional
| ||
zone : | character(*), intent(in), optional
| ||
caltype : | integer, intent(in), optional
| ||
day_seconds : | real(DP),intent(in), optional
| ||
err : | logical, intent(out), optional |
dc_date_types#DC_DATETIME 型変数の生成を行います. 引数 year, mon, day, hour, min, sec の全てを与えない場合, このサブルーチンが呼ばれた際の時刻が使用されます.
引数 caltype には暦法を設定します. dc_date_types#CAL_CYCLIC, dc_date_types#CAL_NOLEAP, dc_date_types#CAL_JULIAN, dc_date_types#CAL_GREGORIAN のいづれかを与えてください. 引数 caltype を指定しない場合, 暦法は dc_date_types#CAL_GREGORIAN に設定されます.
引数 zone には UTC からの時差を設定します. ’+09:00’ や ’-13:00’ のように時差を 6 文字で指定してください. 引数 zone を指定しない場合, date_and_time 組み込みサブルーチン によって得られる時差を設定します.
引数 day_seconds には 1 日何秒かを設定します. この引数を 指定しない場合, dc_date_types#day_seconds の値が用いられます. dc_date_types#day_seconds は SetSecOfDay で変更可能です.
引数 caltype および, zone に不適切な値が与えられた場合, エラーを発生させます. 引数 err を与えている場合には err に .true. が返り, プログラムは続行します.
Original external subprogram is dcdatetimecreate.f90#DCDateTimeCreate1
Subroutine : | |||
diff : | type(DC_DIFFTIME), intent(out) | ||
year : | integer, intent(in), optional
| ||
mon : | integer, intent(in), optional
| ||
day : | integer, intent(in), optional
| ||
hour : | integer, intent(in), optional
| ||
min : | integer, intent(in), optional
| ||
sec : | real(DP),intent(in), optional
| ||
day_seconds : | real(DP),intent(in), optional
|
dc_date_types#DC_DIFFTIME 型変数の生成を行います. 引数 year, mon, day, hour, min, sec を与えない場合, 0 が与えられたことになります.
引数 day_seconds には 1 日何秒かを設定します. この引数を 指定しない場合, dc_date_types#day_seconds の値が用いられます. dc_date_types#day_seconds は SetSecOfDay で変更可能です.
Original external subprogram is dcdatetimecreate.f90#DCDiffTimeCreate1
Subroutine : | |
diff : | type(DC_DIFFTIME), intent(out) |
value : | real(DP), intent(in) |
unit : | character(*), intent(in) |
err : | logical, intent(out), optional |
dc_date_types#DC_DIFFTIME 型変数の生成を行います. 引数 value に数値を, unit に単位を表す文字列を与えてください. unit に指定できるのは以下の文字列です. (大文字小文字は区別しません).
年 : | 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 |
これらに該当しない文字列を unit に与えた場合, エラーを発生させます. 引数 err を与えている場合には err に .true. が返り, プログラムは続行します.
Original external subprogram is dcdatetimecreate.f90#DCDiffTimeCreate2
Subroutine : | |||
time : | type(DC_DATETIME), intent(in) | ||
year : | integer, intent(out), optional
| ||
mon : | integer, intent(out), optional
| ||
day : | integer, intent(out), optional
| ||
hour : | integer, intent(out), optional
| ||
min : | integer, intent(out), optional
| ||
sec : | real(DP),intent(out), optional
| ||
caltype : | integer, intent(out), optional
| ||
zone : | character(*), intent(out), optional
|
dc_date_types#DC_DATETIME 型変数 time を 年 year, 月 mon, 日 day, 時間 hour, 分 min, 秒 sec, 暦法 caltype, タイムゾーン zone に変換して返します.
Original external subprogram is dcdatetimeeval.f90#DCDateTimeEval1
Subroutine : | |||
diff : | type(DC_DIFFTIME), intent(in) | ||
year : | integer, intent(out), optional
| ||
mon : | integer, intent(out), optional
| ||
day : | integer, intent(out), optional
| ||
hour : | integer, intent(out), optional
| ||
min : | integer, intent(out), optional
| ||
sec : | real(DP),intent(out), optional
|
Original external subprogram is dcdatetimeeval.f90#DCDiffTimeEval1
Function : | |
result : | real(DP) |
diff : | type(DC_DIFFTIME), intent(in) |
unit : | character(*), intent(in) |
dc_date_types#DC_DIFFTIME 型変数の日時を unit の単位 に換算して倍精度実数型変数で返します. unit には 日 dc_date_types#UNIT_DAY, 時 dc_date_types#UNIT_HOUR, 分 dc_date_types#UNIT_MIN, 秒 dc_date_types#UNIT_SEC を与えることが可能です. これらに該当しない文字列を unit に与えた場合 0.0d0 が返ります.
Original external subprogram is dcdatetimeeval.f90#DCDiffTimeEvalByUnit
Function : | |
result : | real(DP) |
time : | type(DC_DATETIME), intent(in) |
unit : | character(*), intent(in) |
dc_date_types#DC_DATETIME 型変数の日時を unit の単位 に換算して倍精度実数型変数で返します. unit には 日 dc_date_types#UNIT_DAY, 時 dc_date_types#UNIT_HOUR, 分 dc_date_types#UNIT_MIN, 秒 dc_date_types#UNIT_SEC を与えることが可能です. これらに該当しない文字列を unit に与えた場合 0.0d0 が返ります.
Original external subprogram is dcdatetimeeval.f90#DCDateTimeEvalByUnit
Function : | |
result : | real(DP) |
time : | type(DC_DATETIME), intent(in) |
dc_date_types#DC_DATETIME 型変数の日時を日数に換算して 倍精度実数型変数で返します. (例えば 12 時間は 0.5 日と換算されます).
Original external subprogram is dcdatetimeeval.f90#DCDateTimeEvalDay
Function : | |
result : | real(DP) |
diff : | type(DC_DIFFTIME), intent(in) |
dc_date_types#DC_DIFFTIME 型変数の日時を日数に換算して 倍精度実数型変数で返します. (例えば 12 時間は 0.5 日と換算されます).
1 ヶ月は dc_date_types#CYCLIC_MDAYS と換算します.
Original external subprogram is dcdatetimeeval.f90#DCDiffTimeEvalDay
Function : | |
result : | real(DP) |
time : | type(DC_DATETIME), intent(in) |
dc_date_types#DC_DATETIME 型変数の日時を時間に換算して 倍精度実数型変数で返します. (例えば 2 日は 48 時間に, 30 分 は 0.5 時間と換算されます).
Original external subprogram is dcdatetimeeval.f90#DCDateTimeEvalHour
Function : | |
result : | real(DP) |
diff : | type(DC_DIFFTIME), intent(in) |
dc_date_types#DC_DIFFTIME 型変数の日時を時間に換算して 倍精度実数型変数で返します. (例えば 2 日は 48 時間に, 30 分 は 0.5 時間と換算されます).
1 ヶ月は dc_date_types#CYCLIC_MDAYS と換算します.
Original external subprogram is dcdatetimeeval.f90#DCDiffTimeEvalHour
Function : | |
result : | real(DP) |
time : | type(DC_DATETIME), intent(in) |
dc_date_types#DC_DATETIME 型変数の日時を分に換算して 倍精度実数型変数で返します. (例えば 1 日は 3600 分に, 30 秒 は 0.5 分と換算されます).
Original external subprogram is dcdatetimeeval.f90#DCDateTimeEvalMin
Function : | |
result : | real(DP) |
diff : | type(DC_DIFFTIME), intent(in) |
dc_date_types#DC_DIFFTIME 型変数の日時を分に換算して 倍精度実数型変数で返します. (例えば 1 日は 3600 分に, 30 秒 は 0.5 分と換算されます).
1 ヶ月は dc_date_types#CYCLIC_MDAYS と換算します.
Original external subprogram is dcdatetimeeval.f90#DCDiffTimeEvalMin
Function : | |
result : | real(DP) |
time : | type(DC_DATETIME), intent(in) |
dc_date_types#DC_DATETIME 型変数の日時を秒に換算して 倍精度実数型変数で返します.
Original external subprogram is dcdatetimeeval.f90#DCDateTimeEvalSec
Function : | |
result : | real(DP) |
diff : | type(DC_DIFFTIME), intent(in) |
dc_date_types#DC_DIFFTIME 型変数の日時を秒に換算して 倍精度実数型変数で返します.
1 ヶ月は dc_date_types#CYCLIC_MDAYS と換算します.
Original external subprogram is dcdatetimeeval.f90#DCDiffTimeEvalSec
Subroutine : | |
diff : | type(DC_DIFFTIME), intent(in) |
unit : | integer, intent(in), optional |
dc_date_types#DC_DIFFTIME 型変数の印字を行います. unit には出力先の装置番号を 与えてください. unit を与えない場合, 標準出力へ表示されます.
Original external subprogram is dcdatetimeputline.f90#DCDiffTimePutLine
Subroutine : | |
time : | type(DC_DATETIME), intent(in) |
unit : | integer, intent(in), optional |
dc_date_types#DC_DATETIME 型変数の印字を行います. unit には出力先の装置番号を 与えてください. unit を与えない場合, 標準出力へ表示されます.
Original external subprogram is dcdatetimeputline.f90#DCDateTimePutLine
Subroutine : | |
caltype : | integer, intent(in) |
暦法のデフォルトを設定します. dc_date_types#CAL_CYCLIC, dc_date_types#CAL_NOLEAP, dc_date_types#CAL_JULIAN, dc_date_types#CAL_GREGORIAN のいづれかを引数 caltype に与えてください.
なお, この手続きを呼ばない場合, デフォルトの暦法は dc_date_types#CAL_GREGORIAN に設定されています.
Original external subprogram is dcdatetimesetcaltype.f90#DCDateTimeSetCaltype
Subroutine : | |
sec : | real(DP), intent(in) |
1 日の秒数のデフォルトを設定します.
なお, この手続きを呼ばない場合, デフォルトの 1 日の秒数は dc_date_types#DAY_SECONDS_EARTH に設定されています.
Original external subprogram is dcdatetimesetsecofday.f90#DCDateTimeSetSecOfDay
Subroutine : | |
time : | type(DC_DATETIME), intent(inout) |
zone : | character(*), intent(in) |
err : | logical, intent(out), optional |
引数 time のタイムゾーンを zone へと変更します. 実質的な日時は変更しません.
引数 zone に不適切な値が与えられた場合, エラーを発生させます. 引数 err を与えている場合には err に .true. が返り, プログラムは続行します.
Original external subprogram is dcdatetimezone.f90#DCDateTimeSetZone
Function : | |
result : | logical |
caltype : | integer, intent(in) |
与えられる暦法が dc_date_types 内で有効であれば .true. を, それ以外の場合は .false. を返します.
Original external subprogram is dcdatetimevalidcaltype.f90#DCDateTimeValidCaltype
Function : | |
result : | logical |
zone : | character(*), intent(in) |
与えられるタイムゾーンの表記が有効であれば .true. を, それ以外の場合は .false. を返します.
タイムゾーンの表記は ’+09:00’ のように, 1 文字目が ’+’ または ’-’, 2〜3, 5〜6 文字目が数値で, 4 文字目が ’:’ となります.
Original external subprogram is dcdatetimezone.f90#DCDateTimeValidZone
Function : | |
diff : | type(DC_DIFFTIME) |
zone : | character(*), intent(in) |
与えられるタイムゾーンを dc_date_types#DC_DIFFTIME 変数へと 変換して返します. タイムゾーンの表記が無効な場合は ’+00:00’ が与えられたと解釈します.
Original external subprogram is dcdatetimezone.f90#DCDateTimeZoneToDiff
Subroutine : | |
time : | type(DC_DATETIME), intent(out) |
sec : | real, intent(in) |
dc_date_types#DC_DATETIME 型変数の生成を行います. 引数 sec には秒数を与えてください. 年月日, 時分を使って 指定を行いたい場合は Create を利用してください.
Original external subprogram is dcdatetimecreate.f90#DCDateTimeCreateR
Subroutine : | |
time : | type(DC_DATETIME), intent(out) |
sec : | real(DP), intent(in) |
Original external subprogram is dcdatetimecreate.f90#DCDateTimeCreateD
Subroutine : | |
diff : | type(DC_DIFFTIME), intent(out) |
sec : | real, intent(in) |
dc_date_types#DC_DIFFTIME 型変数の生成を行います. 引数 sec には秒数を与えてください. 年月日, 時分を使って 指定を行いたい場合は Create を利用してください.
Original external subprogram is dcdatetimecreate.f90#DCDiffTimeCreateR
Subroutine : | |
diff : | type(DC_DIFFTIME), intent(out) |
sec : | real(DP), intent(in) |
Original external subprogram is dcdatetimecreate.f90#DCDiffTimeCreateD
Subroutine : | |
day : | integer, intent(inout) |
sec : | real(DP), intent(inout) |
day_seconds : | real(DP), intent(in), optional |
このサブルーチンは内部向けなので dc_date モジュール外では 極力使用しないでください.
日付 day と秒数 sec の正規化を行います. sec が day_seconds (省略される場合は dc_date_types#day_seconds) を超える場合, day に繰上げを行います. また, sec と day の符号が逆の場合, 同符号になるよう 設定します.
subroutine dcdate_normalize(day, sec, day_seconds) ! !=== 日と秒の正規化 ! ! このサブルーチンは内部向けなので dc_date モジュール外では ! 極力使用しないでください. ! ! 日付 *day* と秒数 *sec* の正規化を行います. *sec* が *day_seconds* ! (省略される場合は dc_date_types#day_seconds) を超える場合, *day* ! に繰上げを行います. ! また, *sec* と *day* の符号が逆の場合, 同符号になるよう ! 設定します. ! use dc_date_types, only: day_seconds_default => day_seconds implicit none integer, intent(inout):: day real(DP), intent(inout):: sec real(DP), intent(in), optional:: day_seconds integer:: sgn real(DP):: day_sec continue if (present(day_seconds)) then day_sec = day_seconds else day_sec = day_seconds_default end if if (abs(sec) > day_sec) then day = day + int(sec / day_sec) sec = modulo(sec, day_sec) end if if ((sec > 0.0 .and. day < 0) .or. (sec < 0.0 .and. day > 0)) then sgn = sign(day, 1) day = day - sgn sec = sec + sgn * day_sec endif end subroutine dcdate_normalize
Function : | |
unit : | character(TOKEN) |
str : | character(*), intent(in) |
このサブルーチンは内部向けなので 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 |
これらに該当しない文字列を str に与えた場合, 空文字が返ります.
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 ! ! これらに該当しない文字列を *str* に与えた場合, 空文字が返ります. ! use dc_types, only: TOKEN use dc_date_types, only: UNIT_YEAR, UNIT_MONTH, UNIT_DAY, UNIT_HOUR, UNIT_MIN, UNIT_SEC 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_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
Function : | |
result : | type(DC_DIFFTIME) |
diff1 : | type(DC_DIFFTIME), intent(in) |
diff2 : | type(DC_DIFFTIME), intent(in) |
引数 diff1 を diff2 で除算した際の余りを返します.
※ 注意: 月差と日時の混在する除算は近似的結果になるおそれがあります
type(DC_DIFFTIME) function dcdate_mod_ff(diff1, diff2) result(result) ! ! 引数 <b>diff1</b> を <b>diff2</b> で除算した際の余りを返します. ! ! ※ 注意: 月差と日時の混在する除算は近似的結果になるおそれがあります ! use dc_date_types, only: CYCLIC_MDAYS implicit none type(DC_DIFFTIME), intent(in):: diff1, diff2 real(DP):: sec1, sec2 continue result % day_seconds = diff1 % day_seconds if (diff1 % day == 0 .and. diff2 % day == 0 .and. diff1 % sec == 0.0 .and. diff2 % sec == 0.0) then result % mon = mod(diff1 % mon, diff2 % mon) result % day = 0 result % sec = 0.0 else if (diff1 % sec == 0.0 .and. diff2 % sec == 0.0) then result % mon = 0 result % day = mod((CYCLIC_MDAYS * diff1 % mon + diff1 % day), (CYCLIC_MDAYS * diff2 % mon + diff2 % day)) result % sec = 0.0 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 = mod(sec1, sec2) result % day = 0.0 result % mon = 0.0 call dcdate_normalize(result % day, result % sec, result % day_seconds) endif end function dcdate_mod_ff
Function : | |
result : | type(DC_DIFFTIME) |
factor : | integer, intent(in) |
diff : | type(DC_DIFFTIME), intent(in) |
日時差 diff と facter とを乗算した結果を返します.
type(DC_DIFFTIME) function dcdate_mul_if(factor, diff) result(result) ! ! 日時差 *diff* と *facter* とを乗算した結果を返します. ! 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 call dcdate_normalize(result % day, result % sec, result % day_seconds) end function dcdate_mul_if
Function : | |
result : | type(DC_DIFFTIME) |
diff : | type(DC_DIFFTIME), intent(in) |
factor : | integer, intent(in) |
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
Function : | |
result : | type(DC_DIFFTIME) |
factor : | real, intent(in) |
diff : | type(DC_DIFFTIME), intent(in) |
※ 注意 : 月差を非整数倍すると近似的結果になるおそれがあります
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
Function : | |
result : | type(DC_DIFFTIME) |
diff : | type(DC_DIFFTIME), intent(in) |
factor : | real, intent(in) |
※ 注意 : 月差を非整数倍すると近似的結果になるおそれがあります
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
Function : | |
result : | type(DC_DIFFTIME) |
factor : | real(DP), intent(in) |
diff : | type(DC_DIFFTIME), intent(in) |
※ 注意 : 月差を非整数倍すると近似的結果になるおそれがあります
type(DC_DIFFTIME) function dcdate_mul_df(factor, diff) result(result) ! ! ※ 注意 : 月差を非整数倍すると近似的結果になるおそれがあります use dc_types, only: DP use dc_date_types, only: CYCLIC_MDAYS implicit none real(DP), intent(in):: factor type(DC_DIFFTIME), intent(in):: diff real(DP):: 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 call dcdate_normalize(result % day, result % sec, result % day_seconds) end function dcdate_mul_df
Function : | |
result : | type(DC_DIFFTIME) |
diff : | type(DC_DIFFTIME), intent(in) |
factor : | real(DP), intent(in) |
※ 注意 : 月差を非整数倍すると近似的結果になるおそれがあります
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
Function : | |
result : | type(DC_DATETIME) |
diff : | type(DC_DIFFTIME), intent(in) |
time : | type(DC_DATETIME), intent(in) |
2 つの日時 (DC_DATETIME 型) もしくは 日時差 (DC_DIFFTIME 型)の加算を行います.
type(DC_DATETIME) function dcdate_add_ft(diff, time) result(result) ! ! 2 つの日時 (DC_DATETIME 型) もしくは ! 日時差 (DC_DIFFTIME 型)の加算を行います. ! implicit none type(DC_DIFFTIME), intent(in):: diff type(DC_DATETIME), intent(in):: time integer:: time_year, time_mon, time_day, time_caltype real(DP):: time_sec character(6):: time_zone continue call Eval(time, year = time_year, mon = time_mon, day = time_day, sec = time_sec, caltype = time_caltype, zone = time_zone) call Create(result, year=time_year, mon = time_mon + diff % mon, day = time_day + diff % day, sec = time_sec + diff % sec, caltype = time_caltype, zone = time_zone) end function dcdate_add_ft
Function : | |
result : | type(DC_DIFFTIME) |
diff1 : | type(DC_DIFFTIME), intent(in) |
diff2 : | type(DC_DIFFTIME), intent(in) |
type(DC_DIFFTIME) function dcdate_add_ff(diff1, diff2) result(result) 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_normalize(result % day, result % sec, result % day_seconds) end function dcdate_add_ff
Function : | |
result : | type(DC_DATETIME) |
time : | type(DC_DATETIME), intent(in) |
diff : | type(DC_DIFFTIME), intent(in) |
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
Function : | |
result : | type(DC_DATETIME) |
time : | type(DC_DATETIME), intent(in) |
diff : | type(DC_DIFFTIME), intent(in) |
2 つの日時 (DC_DATETIME 型) もしくは 日時差 (DC_DIFFTIME 型)の減算を行います.
type(DC_DATETIME) function dcdate_sub_tf(time, diff) result(result) ! ! 2 つの日時 (DC_DATETIME 型) もしくは ! 日時差 (DC_DIFFTIME 型)の減算を行います. ! implicit none type(DC_DATETIME), intent(in):: time type(DC_DIFFTIME), intent(in):: diff integer:: time_year, time_mon, time_day, time_caltype real(DP):: time_sec character(6):: time_zone continue call Eval(time, year = time_year, mon = time_mon, day = time_day, sec = time_sec, caltype = time_caltype, zone = time_zone) call Create(result, year=time_year, mon = time_mon - diff % mon, day = time_day - diff % day, sec = time_sec - diff % sec, caltype = time_caltype, zone = time_zone) end function dcdate_sub_tf
Function : | |
result : | type(DC_DIFFTIME) |
diff1 : | type(DC_DIFFTIME), intent(in) |
diff2 : | type(DC_DIFFTIME), intent(in) |
type(DC_DIFFTIME) function dcdate_sub_ff(diff1, diff2) result(result) 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_normalize(result % day, result % sec, result % day_seconds) end function dcdate_sub_ff
Function : | |
result : | type(DC_DIFFTIME) |
time1 : | type(DC_DATETIME), intent(in) |
time2 : | type(DC_DATETIME), intent(in) |
type(DC_DIFFTIME) function dcdate_sub_tt(time1, time2) result(result) implicit none type(DC_DATETIME), intent(in):: time1, time2 continue result % day = time1 % day - time2 % day result % sec = time1 % sec - time2 % sec + EvalSec(ZoneToDiff(time1 % zone) - ZoneToDiff(time2 % zone)) result % day_seconds = time1 % day_seconds call dcdate_normalize(result % day, result % sec, result % day_seconds) end function dcdate_sub_tt
Function : | |
result : | type(DC_DIFFTIME) |
diff : | type(DC_DIFFTIME), intent(in) |
denominator : | integer, intent(in) |
日時差 diff を denominator で除算した結果を返します.
※ 注意 : 月差を除算すると近似的結果になるおそれがあります
type(DC_DIFFTIME) function dcdate_div_fi(diff, denominator) result(result) ! ! 日時差 *diff* を *denominator* で除算した結果を返します. ! ! ※ 注意 : 月差を除算すると近似的結果になるおそれがあります use dc_date_types, only: CYCLIC_MDAYS implicit none type(DC_DIFFTIME), intent(in):: diff integer, intent(in):: denominator continue result % mon = diff % mon / denominator ! 月からの近似的繰り下がりは日単位でしか行わない result % day = diff % day / denominator + int((CYCLIC_MDAYS * mod(diff % mon, denominator)) / denominator) result % sec = diff % sec / denominator + (diff % day_seconds * mod(diff % day, denominator)) / denominator end function dcdate_div_fi
Function : | |
result : | type(DC_DIFFTIME) |
diff : | type(DC_DIFFTIME), intent(in) |
denominator : | real(DP), intent(in) |
※ 注意 : 月差を除算すると近似的結果になるおそれがあります
type(DC_DIFFTIME) function dcdate_div_fd(diff, denominator) result(result) ! ! ※ 注意 : 月差を除算すると近似的結果になるおそれがあります use dc_date_types, only: CYCLIC_MDAYS implicit none type(DC_DIFFTIME), intent(in):: diff real(DP), intent(in):: denominator real(DP):: month, day continue month = diff % mon / denominator result % mon = int(month) day = 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 call dcdate_normalize(result % day, result % sec, result % day_seconds) end function dcdate_div_fd
Function : | |
result : | type(DC_DIFFTIME) |
diff : | type(DC_DIFFTIME), intent(in) |
denominator : | real, intent(in) |
※ 注意 : 月差を除算すると近似的結果になるおそれがあります
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
Function : | |
result : | real(DP) |
diff1 : | type(DC_DIFFTIME), intent(in) |
diff2 : | type(DC_DIFFTIME), intent(in) |
※ 注意 : 月差と日時の混在する除算は近似的結果になるおそれがあります
real(DP) function dcdate_div_ff(diff1, diff2) result(result) ! ! ※ 注意 : 月差と日時の混在する除算は近似的結果になるおそれがあります use dc_date_types, only: CYCLIC_MDAYS 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
Function : | |
result : | logical |
time1 : | type(DC_DATETIME), intent(in) |
time2 : | type(DC_DATETIME), intent(in) |
2 つの引数の日時を比較します. 2 つ目の引数に格納される日時が 1 つ目の引数に格納される日時 よりも進んでいる場合, .true. が返ります.
logical function dcdate_lt_tt(time1, time2) result(result) ! ! 2 つの引数の日時を比較します. ! 2 つ目の引数に格納される日時が 1 つ目の引数に格納される日時 ! よりも進んでいる場合, .true. が返ります. ! implicit none type(DC_DATETIME), intent(in):: time1, time2 continue result = .not. dcdate_gt_tt(time1, time2) end function dcdate_lt_tt
Function : | |
result : | logical |
diff1 : | type(DC_DIFFTIME), intent(in) |
diff2 : | type(DC_DIFFTIME), intent(in) |
2 つの引数の日時差を比較します. 2 つ目の引数に格納される日時差が 1 つ目の引数に格納される日時差 よりも大きい場合, .true. が返ります.
logical function dcdate_lt_ff(diff1, diff2) result(result) ! ! 2 つの引数の日時差を比較します. ! 2 つ目の引数に格納される日時差が 1 つ目の引数に格納される日時差 ! よりも大きい場合, .true. が返ります. ! implicit none type(DC_DIFFTIME), intent(in):: diff1, diff2 continue result = .not. dcdate_gt_ff(diff1, diff2) end function dcdate_lt_ff
Function : | |
result : | logical |
time1 : | type(DC_DATETIME), intent(in) |
time2 : | type(DC_DATETIME), intent(in) |
2 つの引数の日時を比較します. 1 つ目の引数に格納される日時が 2 つ目の引数に格納される日時 と同じ場合, .true. が返ります.
logical function dcdate_eq_tt(time1, time2) result(result) ! ! 2 つの引数の日時を比較します. ! 1 つ目の引数に格納される日時が 2 つ目の引数に格納される日時 ! と同じ場合, .true. が返ります. ! implicit none type(DC_DATETIME), intent(in):: time1, time2 real(DP):: time1_sec, time2_sec continue time1_sec = EvalSec(time1) + EvalSec(ZoneToDiff(time1 % zone)) time2_sec = EvalSec(time2) + EvalSec(ZoneToDiff(time2 % zone)) if (time1_sec == time2_sec) then result = .true. else result = .false. end if end function dcdate_eq_tt
Function : | |
result : | logical |
diff1 : | type(DC_DIFFTIME), intent(in) |
diff2 : | type(DC_DIFFTIME), intent(in) |
2 つの引数の日時差を比較します. 1 つ目の引数に格納される日時差が 2 つ目の引数に格納される日時差 と同じ場合, .true. が返ります.
logical function dcdate_eq_ff(diff1, diff2) result(result) ! ! 2 つの引数の日時差を比較します. ! 1 つ目の引数に格納される日時差が 2 つ目の引数に格納される日時差 ! と同じ場合, .true. が返ります. ! implicit none type(DC_DIFFTIME), intent(in):: diff1, diff2 continue if (EvalSec(diff1) == EvalSec(diff2)) then result = .true. else result = .false. end if end function dcdate_eq_ff
Function : | |
result : | logical |
i : | integer, intent(in) |
diff : | type(DC_DIFFTIME), intent(in) |
引数 diff の日時差が i と等しいかどうかを比較します. diff を秒数に換算した値と i とが等しい場合, .true. が返ります.
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
Function : | |
result : | logical |
diff : | type(DC_DIFFTIME), intent(in) |
i : | integer, intent(in) |
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
Function : | |
result : | logical |
r : | real, intent(in) |
diff : | type(DC_DIFFTIME), intent(in) |
引数 diff の日時差が r と等しいかどうかを比較します. diff を秒数に換算した値と r とが等しい場合, .true. が返ります.
logical function dcdate_eq_rf(r, diff) result(result) ! ! 引数 *diff* の日時差が *r* と等しいかどうかを比較します. *diff* ! を秒数に換算した値と *r* とが等しい場合, .true. が返ります. ! implicit none type(DC_DIFFTIME), intent(in):: diff real, intent(in):: r continue if (real(EvalSec(diff)) == r) then result = .true. else result = .false. end if end function dcdate_eq_rf
Function : | |
result : | logical |
diff : | type(DC_DIFFTIME), intent(in) |
r : | real, intent(in) |
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
Function : | |
result : | logical |
d : | real(DP), intent(in) |
diff : | type(DC_DIFFTIME), intent(in) |
引数 diff の日時差が d と等しいかどうかを比較します. diff を秒数に換算した値と d とが等しい場合, .true. が返ります.
logical function dcdate_eq_df(d, diff) result(result) ! ! 引数 *diff* の日時差が *d* と等しいかどうかを比較します. *diff* ! を秒数に換算した値と *d* とが等しい場合, .true. が返ります. ! use dc_types, only: DP implicit none type(DC_DIFFTIME), intent(in):: diff real(DP), intent(in):: d continue if (real(EvalSec(diff)) == d) then result = .true. else result = .false. end if end function dcdate_eq_df
Function : | |
result : | logical |
diff : | type(DC_DIFFTIME), intent(in) |
d : | real(DP), intent(in) |
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
Function : | |
result : | logical |
time1 : | type(DC_DATETIME), intent(in) |
time2 : | type(DC_DATETIME), intent(in) |
2 つの引数の日時を比較します. 1 つ目の引数に格納される日時が 2 つ目の引数に格納される日時 よりも進んでいる場合, .true. が返ります.
logical function dcdate_gt_tt(time1, time2) result(result) ! ! 2 つの引数の日時を比較します. ! 1 つ目の引数に格納される日時が 2 つ目の引数に格納される日時 ! よりも進んでいる場合, .true. が返ります. ! implicit none type(DC_DATETIME), intent(in):: time1, time2 real(DP):: time1_sec, time2_sec continue time1_sec = EvalSec(time1) + EvalSec(ZoneToDiff(time1 % zone)) time2_sec = EvalSec(time2) + EvalSec(ZoneToDiff(time2 % zone)) if (time1_sec > time2_sec) then result = .true. else result = .false. end if end function dcdate_gt_tt
Function : | |
result : | logical |
diff1 : | type(DC_DIFFTIME), intent(in) |
diff2 : | type(DC_DIFFTIME), intent(in) |
2 つの引数の日時差を比較します. 1 つ目の引数に格納される日時差が 2 つ目の引数に格納される日時差 よりも大きい場合, .true. が返ります.
logical function dcdate_gt_ff(diff1, diff2) result(result) ! ! 2 つの引数の日時差を比較します. ! 1 つ目の引数に格納される日時差が 2 つ目の引数に格納される日時差 ! よりも大きい場合, .true. が返ります. ! use dc_date_types, only: CYCLIC_MDAYS implicit none type(DC_DIFFTIME), intent(in):: diff1, diff2 continue if (EvalSec(diff1) > EvalSec(diff2)) then result = .true. else result = .false. end if end function dcdate_gt_ff
Function : | |
result : | character(STRING) |
time : | type(DC_DATETIME), intent(in) |
dc_date_types#DC_DATETIME 型変数を文字型変数へ変換して返します. 書式は下記のように JIS X 0301 の完全表記です.
YYYY-MM-DDThh:mm:ss.sTZD
YYYY は年, MM は月, DD は日, hh は時, mm は分, ss.s は秒, TZD はタイムゾーンを表します.
Original external subprogram is dcdatetimetochar.f90#DCDateTimeToChar
Function : | |
result : | character(STRING) |
diff : | type(DC_DIFFTIME), intent(in) |
dc_date_types#DC_DIFFTIME 型変数を文字型変数へ変換して返します. 書式は以下のようになります.
YYYY-MM-DDThh:mm:ss.s
YYYY は年, MM は月, DD は日, hh は時, mm は分, ss.s は秒を表します. ただし, DD は 2 桁を超える場合があります. (dc_date_types#DC_DIFFTIME は X ヶ月後, X 日前, などを表現するため のデータ型なので, 日を月に繰り上げたり, 月を日に繰り下げることを しません. また「年」の情報も持ちません. 1 年の日数や 1 月の日数は dc_date_types#DC_DATETIME 側で決まります).
Original external subprogram is dcdatetimetochar.f90#DCDiffTimeToChar