Class gt4_history_nmlinfo
In: gt4_history_nmlinfo.f90

NAMELIST の使用を想定したヒストリデータ出力情報管理用ユーティリティ

Utilities for history data output information management assuming use of NAMELIST

Note that Japanese and English are described in parallel.

比較的大規模な数値モデルにおいて, データ出力の情報管理を 支援するためのモジュールです. 個別のモジュールがそれぞれ独立にデータ出力を行うことと, NAMELIST を用いて出力ファイルや出力間隔などを変更すること 想定して設計されています. ただし, このモジュール自体の主目的は情報の管理であり, 実際のデータ出力には gt4_history モジュールを 使用することに注意してください.

This module supports information management of data output in a comparatively large-scale numerical model. This module is designed expecting that individual modules perform data output independently, and output filename or output interval is changed from NAMELIST. Note that the purpose of this module is information management, therefore actual data output is performed by "gt4_history" module.

Procedures List

Create :GT_HISTORY_NMLINFO 型変数の初期設定
Close :GT_HISTORY_NMLINFO 型変数の終了処理
PutLine :GT_HISTORY_NMLINFO 型変数に格納されている情報の印字
initialized :GT_HISTORY_NMLINFO 型変数が初期設定されているか否か
define_mode :定義モードの場合に真を返す
EndDefine :変数情報定義モードから出力モードへ変更
Add :変数情報の追加
Delete :変数情報の削除
ResetDefault :デフォルト設定のみに戻す
Inquire :変数情報の問い合わせ
AssocGtHist :変数に応じた gt4_history#GT_HISTORY 型変数を返す
OutputStepDisable :output_step が常に .false. を返すよう設定する
output_valid :出力の設定が有効である場合に真を返す
output_step :現在の時刻が出力のタイミングの場合に真を返す
names :登録されている変数名リストを返す
———— :————
Create :Constructor of "GT_HISTORY_NMLINFO"
Close :Deconstructor of "GT_HISTORY_NMLINFO"
PutLine :Print information of "GT_HISTORY_NMLINFO"
initialized :Check initialization of "GT_HISTORY_NMLINFO"
define_mode :True is returned if current state is define mode
EndDefine :Change define mode about information of variables to output mode
Add :Add information of variables
Delete :Delete information of variables
ResetDefault :Reset to default settings
Inquire :Inquire information of variables
AssocGtHist :"gt4_history#GT_HISTORY" correspond to variable is returned
OutputStepDisable :Configure that "output_step" returns .false. already
output_valid :True is returned when a configuration of output is valid
output_step :True is returned when current time is output timing
names :Return list of registerd variable identifiers

Usage

このモジュールは以下のような手順で用いてください.

このモジュールを使用したサンプル Fortran プログラム 作成スクリプトが www.gfd-dennou.org/library/dcpam/dcpam4/dcpam4_current/script/f90/dcmodel_f90sample_maker.rb から入手できます. Ruby で記述されており, 実行することで サンプルとなる Fortran プログラムが作成されます. 下記の解説のみでは実際の利用法やご利益が分かりにくいため, サンプル Fortran プログラムを実際に見てみることをオススメします.

  1. モジュール内で, GT_HISTORY_NMLINFO 型の変数を定義しておきます.
  2. Create を用いて, GT_HISTORY_NMLINFO 型の変数の初期設定を行います. この際にデフォルトの出力間隔 interval_value, interval_unit, 精度 precision, 平均化 average, 出力ファイル名接頭詞 fileprefix を設定します.
  3. プログラムがデフォルトで出力する変数がある場合, Add を 使用して登録してください. name には変数名を与えます. name は変数を識別するためのキーと して利用するため, 異なる変数に対して同じ name を指定しないで ください. file には出力ファイル名を与えます. 与えない場合, 適当なファイル名が設定されます. その他の情報は上記と同様です.
  4. NAMELIST から得られた出力変数の情報を, Add を使用して登録してください. Add で既に登録済みの name を再度登録することで, 設定が上書きされます.
  5. 登録が完了したら, EndDefine を用いて, 定義モードから出力モードへ 移行してください.
  6. gt4_history#HistoryCreate, gt4_history#HistoryAddVariable gt4_history#HistoryPut 等で出力設定およびデータ出力を行う際には, AssocGtHist に対し, 変数名 namegt4_history#GT_HISTORY 型のポインタ history を渡してください. GT_HISTORY_NMLINFO 型の変数に登録されている name に関する gt4_history#GT_HISTORY 型変数に結合された history が返ります. この history を上記 gt4_history のサブルーチン群の引数 history に渡して 出力設定およびデータ出力を行ってください. gt4_history#HistoryCreate に必要な 出力間隔や精度は Inquire を用いて得ることができます. 使用後は, NULLIFY によって history を空状態にしてください. (DEALLOCATE を使用すると出力に関する情報が失われるため, 使用しないでください).

    それぞれの変数に関して, 出力設定が有効かどうかについては, output_valid で知ることが可能です.

    また, 時間積分中に gt4_history#HistoryPut を使用する際 に, 現在時刻が出力タイミングかどうかについては, output_step で知ることが可能です.

  7. ファイルの出力が終了したら, 上記の手順と同様に gt4_history#GT_HISTORY 型の history を取得し, gt4_history#HistoryClose によって終了処理を行ってください.
  8. 最後に, Close によって, GT_HISTORY_NMLINFO 型の変数の 終了処理を行います.

Use this module as follows.

Sample Fortran programs generator (Ruby script) is available from www.gfd-dennou.org/library/dcpam/dcpam4/dcpam4_current/script/f90/dcmodel_f90sample_maker.rb . Sample Fortran programs are created by executing this script. Because neither actual usage nor the profit are understood easily only from the following explanations, It is recommended to see sample Fortran programs actually.

  1. Declare "GT_HISTORY_NMLINFO" variable in the module.
  2. Initialize "GT_HISTORY_NMLINFO" variable by "Create". On this occasion, configure default interval_value, interval_unit, precision, average, fileprefix (prefix of output file).
  3. Register by using "Add" when there are variables that the program outputs by default. variable identifier is given to name. Do not specify same name for different variables because name is used as a key to identify the variable. The output file name is given to file. A suitable file name is set when not giving it. The extra information is similar to the above-mentioned.
  4. Register information of output variables obtained from NAMELIST by using "Add". When registered name is registered again, the setting concerning the name has been overwritten.
  5. Shift from the define mode to output mode by using "EndDefine" when registration is completed.
  6. Pass "AssocGtHist" variable identifier name and history of "gt4_history#GT_HISTORY" pointer when setting output and data output is performed with "gt4_history#HistoryCreate", "gt4_history#HistoryAddVariable" "gt4_history#HistoryPut" etc. history is associated to "gt4_history#GT_HISTORY" correspond to name stored in "GT_HISTORY_NMLINFO" variable. Pass the history to subroutines in "gt4_history" above-mentioned, and configure output setting and output data. Necessary output interval and precision for "gt4_history#HistoryCreate" can be obtained by using "Inquire". Please put history into a null state by NULLIFY after use. (Information of output is lost when DEALLOCATE is used, so do not use it).

    It can know whether the output setting is effective for each variable with "output_valid".

    Moreover, it can know time now to be whether output timing when "gt4_history#HistoryPut" is used while integrating time with "output_step".

  7. Acquire history of "gt4_history#GT_HISTORY" type as well as the above-mentioned procedure, and terminate that by "gt4_history#HistoryClose" when the output of the file ends.
  8. Finally, the termination of the variable of "GT_HISTORY_NMLINFO" type is done by "Close".

Methods

Included Modules

dc_types gt4_history dc_trace dc_string dc_present dc_message dc_error dc_date_types dc_date dc_hash

Attributes

Public Instance methods

Add( gthstnml, [name], [file], [interval_value], [interval_unit], [precision], [average], [fileprefix], [err] )
Subroutine :recursive
gthstnml :type(GT_HISTORY_NMLINFO), intent(inout)
name :character(*), intent(in), optional
: 変数名.

先頭の空白は無視されます.

"Data1,Data2" のようにカンマで区切って複数 の変数を指定することも可能です.

Variable identifier.

Blanks at the head of the name are ignored.

Multiple variables can be specified as "Data1,Data2" too. Delimiter is comma.

file :character(*), intent(in), optional
: ヒストリデータのファイル名. History data filenames
interval_value :real, intent(in), optional
: ヒストリデータの出力間隔の数値. 負の値を与えると, 出力を抑止します.

Numerical value for interval of history data output. Negative values suppresses output.

interval_unit :character(*), intent(in), optional
: ヒストリデータの出力間隔の単位. Unit for interval of history data output
precision :character(*), intent(in), optional
: ヒストリデータの精度. Precision of history data
average :logical, intent(in), optional
: 出力データの平均化フラグ. Flag for average of output data.
fileprefix :character(*), intent(in), optional
: ヒストリデータのファイル名の接頭詞. Prefixes of history data filenames
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

変数の出力情報を加えます.

デフォルト値を設定するには, name を与えないか, または name に空白を与えてください. デフォルト値を与える場合, file に与えられる情報は無視されます. fileprefix はデフォルト値に与える場合のみ有効です.

name に変数名が指定され, その際に file が与えられない, または空白が与えられる場合, file には "<name に与えられた文字>.nc" が指定されます.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合, プログラムはエラーを発生させます.

Add output information of a variable.

In order to set default values, specify blank to name or do not specify name. When default values are specified, file is ignored. fileprefix is valid only when default values are specified.

When a variable identifier is specified to name and file is not specified or blanks are specified to file, "<string given to name>.nc" is specified to file.

If gthstnml is not initialized by "Create" yet, error is occurred.

Alias for HistoryNmlInfoAdd

AssocGtHist( gthstnml, name, history, [err] )
Subroutine :
gthstnml :type(GT_HISTORY_NMLINFO), intent(in)
name :character(*), intent(in)
: 変数名. 先頭の空白は無視されます.

Variable identifier. Blanks at the head of the name are ignored.

history :type(GT_HISTORY), pointer
: (out)

gt4_history モジュール用構造体. Derived type for "gt4_history" module

err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

与えられた gt4_history#GT_HISTORY 型のポインタ history に対し, gthstnml 内の name に関する gt4_history#GT_HISTORY 型変数を 結合します. 空状態の history を与えてください.

EndDefine で定義モードから出力モードに 移行した後に呼び出してください. EndDefine を呼ぶ前にこのサブルーチンを使用すると, プログラムはエラーを発生させます.

name に関する情報が見当たらない場合, プログラムはエラーを発生させます. name が空文字の場合にも, プログラムはエラーを発生させます.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合にも, プログラムはエラーを発生させます.

This subroutine associates given "gt4_history#GT_HISTORY" pointer history to "gt4_history#GT_HISTORY" correspond to name in gthstnml. Give null history.

Use after state is changed from define mode to output mode by "EndDefile". If this subroutine is used before "EndDefine" is used, error is occurred.

When data correspond to name is not found, error is occurred. When name is blank, error is occurred too.

If gthstnml is not initialized by "Create" yet, error is occurred.

Alias for HistoryNmlInfoAssocGtHist

Close( gthstnml, [err] )
Subroutine :
gthstnml :type(GT_HISTORY_NMLINFO), intent(inout)
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

GT_HISTORY_NMLINFO 型の変数の終了処理を行います.

このサブルーチンを使用する前に, gthstnml に格納されている gt4_history#GT_HISTORY 型の全ての変数に対して, gt4_history#HistoryClose を用いて終了処理を行ってください. 終了処理されていないものがある場合, プログラムはエラーを発生させます.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合, プログラムはエラーを発生させます.

Deconstructor of "GT_HISTORY_NMLINFO".

Terminate all "gt4_history#GT_HISTORY" variables in gthstnml by "gt4_history#HistoryClose" before this subroutine is used. If unterminated variables remain, error is occurred.

Note that if gthstnml is not initialized by "Create" yet, error is occurred.

Alias for HistoryNmlInfoClose

Create( gthstnml, [interval_value], [interval_unit], [precision], [average], [fileprefix], [err] )
Subroutine :
gthstnml :type(GT_HISTORY_NMLINFO), intent(inout)
interval_value :real, intent(in), optional
: ヒストリデータの出力間隔の数値. 負の値を与えると, 出力を抑止します.

Numerical value for interval of history data output. Negative values suppresses output.

interval_unit :character(*), intent(in), optional
: ヒストリデータの出力間隔の単位. Unit for interval of history data output
precision :character(*), intent(in), optional
: ヒストリデータの精度. Precision of history data
average :logical, intent(in), optional
: 出力データの平均化フラグ. Flag for average of output data.
fileprefix :character(*), intent(in), optional
: ヒストリデータのファイル名の接頭詞. Prefixes of history data filenames
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

GT_HISTORY_NMLINFO 型の変数の初期設定を行います. 他のサブルーチンを使用する前に必ずこのサブルーチンによって GT_HISTORY_NMLINFO 型の変数を初期設定してください.

interval_value, interval_unit, precision, average はデフォルト値として設定されます. fileprefix は各変数の出力ファイル名の接頭詞として 使用されます.

なお, 与えられた gthstnml が既に初期設定されている場合, プログラムはエラーを発生させます.

Constructor of "GT_HISTORY_NMLINFO". Initialize gthstnml by this subroutine, before other procedures are used,

interval_value, interval_unit, precision, average are set as default values. fileprefix is used as prefixes of output filenames of each variable.

Note that if gthstnml is already initialized by this procedure, error is occurred.

Alias for HistoryNmlInfoCreate

Delete( gthstnml, name, [err] )
Subroutine :recursive
gthstnml :type(GT_HISTORY_NMLINFO), intent(inout)
name :character(*), intent(in)
: 変数名.

先頭の空白は無視されます.

"Data1,Data2" のようにカンマで区切って複数 の変数を指定することが可能です.

Variable identifier.

Blanks at the head of the name are ignored.

Multiple variables can be specified as "Data1,Data2". Delimiter is comma.

err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

変数の出力情報を削除します.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合, プログラムはエラーを発生させます.

Delete output information of a variable.

If gthstnml is not initialized by "Create" yet, error is occurred.

Alias for HistoryNmlInfoDelete

EndDefine( gthstnml, [err] )
Subroutine :
gthstnml :type(GT_HISTORY_NMLINFO), intent(inout)
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

定義モードから出力モードに移行し, gthstnml に設定した情報を確定します. AssocGTHist サブルーチンを呼び出す前に, 必ずこのサブルーチンを呼び出してください. このサブルーチンを呼んだ後に Add, Delete, ResetDefault を呼ぶとプログラムはエラーを発生させます.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合にも, プログラムはエラーを発生させます.

Change from define mode to output mode, and determine information configured in gthstnml. Use this subroutine before "AssocGTHist" is used. If "Add", "Delete", "ResetDefault" are used after this subroutine is used, error is occurred.

If gthstnml is not initialized by "Create" yet, error is occurred.

Alias for HistoryNmlInfoEndDefine

GT_HISTORY_NMLINFO
Derived Type :
initialized = .false. :logical
: 初期設定フラグ. Initialization flag
define_mode = .true. :logical
: 定義状態を表すフラグ. Flag that represents define mode
gthstnml_list =>null() :type(GT_HISTORY_NMLINFO_ENTRY), pointer
: 変数ごとの情報リスト. 格納される情報については GT_HISTORY_NMLINFO_ENTRY を参照のこと.

Information list about individual variable See "GT_HISTORY_NMLINFO_ENTRY" about stored information.

NAMELIST から取得したヒストリデータの出力情報 を格納するための構造データ型です. まず, Create で "GT_HISTORY_NMLINFO" 型の変数を初期設定して下さい. 初期設定された "GT_HISTORY_NMLINFO" 型の変数を再度利用する際には, Close によって終了処理を行ってください.

This derived type is worked in order to store information about data output from NAMELIST. Initialize "GT_HISTORY_NMLINFO" variable by "Create" before usage. If you reuse "GT_HISTORY_NMLINFO" variable again for another application, terminate by "Close".

Inquire( gthstnml, [name], [file], [interval_value], [interval_unit], [precision], [average], [fileprefix], [err] )
Subroutine :
gthstnml :type(GT_HISTORY_NMLINFO), intent(in)
name :character(*), intent(in), optional
: 変数名. 先頭の空白は無視されます.

Variable identifier. Blanks at the head of the name are ignored.

file :character(*), intent(out), optional
: ヒストリデータのファイル名. History data filenames
interval_value :real, intent(out), optional
: ヒストリデータの出力間隔の数値. 負の値を与えると, 出力を抑止します.

Numerical value for interval of history data output. Negative values suppresses output.

interval_unit :character(*), intent(out), optional
: ヒストリデータの出力間隔の単位. Unit for interval of history data output
precision :character(*), intent(out), optional
: ヒストリデータの精度. Precision of history data
average :logical, intent(out), optional
: 出力データの平均化フラグ. Flag for average of output data.
fileprefix :character(*), intent(out), optional
: ヒストリデータのファイル名の接頭詞. Prefixes of history data filenames
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

変数の出力情報を取得します.

デフォルト値を取得するには, name を与えないか, または name に空白を与えてください.

name に関するデータが存在しない場合, エラーを発生させます.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合, プログラムはエラーを発生させます.

Inquire output information of a variable.

If data correspond to name is not found, error is occurred.

If gthstnml is not initialized by "Create" yet, error is occurred.

Alias for HistoryNmlInfoInquire

OutputStepDisable( gthstnml, name, [err] )
Subroutine :
gthstnml :type(GT_HISTORY_NMLINFO), intent(in)
name :character(*), intent(in)
: 変数名. 先頭の空白は無視されます.

Variable identifier. Blanks at the head of the name are ignored.

err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

このサブルーチンを使用すると, name に関して, 以降は output_step が常に .false. を返すようになります.

データ出力間隔を出力の初期設定から変更し, データを出力するたびに時刻を明示的に指定する場合に利用することを 想定しています.

EndDefine で定義モードから出力モードに 移行した後に呼び出してください. EndDefine を呼ぶ前にこのサブルーチンを使用すると, プログラムはエラーを発生させます.

name に関する情報が見当たらない場合, プログラムはエラーを発生させます. name が空文字の場合にも, プログラムはエラーを発生させます.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合にも, プログラムはエラーを発生させます.

After this subroutine is used, "output_step" returns .false. already corresponding to the name.

This subroutine expected to use when interval of data output is changed from initialization of output, and time is specified explicitly whenever data is output.

Use after state is changed from define mode to output mode by "EndDefile". If this subroutine is used before "EndDefine" is used, error is occurred.

When data correspond to name is not found, error is occurred. When name is blank, error is occurred too.

If gthstnml is not initialized by "Create" yet, error is occurred.

Alias for HistoryNmlInfoOutputStepDisable

PutLine( gthstnml, [unit], [indent], [err] )
Subroutine :
gthstnml :type(GT_HISTORY_NMLINFO), intent(in)
unit :integer, intent(in), optional
: 出力先の装置番号. デフォルトの出力先は標準出力.

Unit number for output. Default value is standard output.

indent :character(*), intent(in), optional
: 表示されるメッセージの字下げ.

Indent of displayed messages.

err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

引数 gthstnml に設定されている情報を印字します. デフォルトではメッセージは標準出力に出力されます. unit に装置番号を指定することで, 出力先を変更することが可能です.

Print information of gthstnml. By default messages are output to standard output. Unit number for output can be changed by unit argument.

Alias for HistoryNmlInfoPutLine

ResetDefault( gthstnml, [err] )
Subroutine :
gthstnml :type(GT_HISTORY_NMLINFO), intent(inout)
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

デフォルト値を残し, 登録したデータを削除します.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合, プログラムはエラーを発生させます.

Stored data is deleted without default settings.

If gthstnml is not initialized by "Create" yet, error is occurred.

Alias for HistoryNmlInfoResetDefault

define_mode( gthstnml ) result(result)
Function :
result :logical
gthstnml :type(GT_HISTORY_NMLINFO), intent(in)

gthstnml が定義モードであれば .true. が, 定義モードでなければ .false. が返ります.

If gthstnml is define mode, .true. is returned. If gthstnml is not define mode, .false. is returned.

Alias for HistoryNmlInfoDefineMode

initialized( gthstnml ) result(result)
Function :
result :logical
gthstnml :type(GT_HISTORY_NMLINFO), intent(in)

gthstnml が初期設定されている場合には .true. が, 初期設定されていない場合には .false. が返ります.

If gthstnml is initialized, .true. is returned. If gthstnml is not initialized, .false. is returned.

Alias for HistoryNmlInfoInitialized

names( gthstnml ) result(result)
Function :
result :character(STRING)
gthstnml :type(GT_HISTORY_NMLINFO), intent(in)

gthstnml が設定されている変数リストをカンマでつなげて 返します.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合, 空文字が返ります.

List of variables registered in gthstnml is join with ’,’, and returned.

If gthstnml is not initialized by "Create" yet, blank is returned.

Alias for HistoryNmlInfoNames

output_step( gthstnml, name, time ) result(result)
Function :
result :logical
gthstnml :type(GT_HISTORY_NMLINFO), intent(in)
name :character(*), intent(in)
: 変数名. 先頭の空白は無視されます.

Variable identifier. Blanks at the head of the name are ignored.

time :type(DC_DIFFTIME), intent(in)
: 現在時刻. Current time

time が変数 name の出力されるタイミングであれば .true. を, そうでなければ .false. を返します. gthstnml が初期設定されていない場合にも .false. が返ります. name に関するデータが存在しない場合にも .false. が返ります.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合, プログラムはエラーを発生させます.

If time is the time that a variable name is output, .true. is returned, otherwise .false. is returned When gthstnml is not initialized, .false. is returned too. When data correspond to name is not found, .false. is returned too.

If gthstnml is not initialized by "Create" yet, error is occurred.

Alias for HistoryNmlInfoOutputStep

output_valid( gthstnml, name ) result(result)
Function :
result :logical
gthstnml :type(GT_HISTORY_NMLINFO), intent(in)
name :character(*), intent(in)
: 変数名. 先頭の空白は無視されます.

Variable identifier. Blanks at the head of the name are ignored.

変数 name の出力が有効であれば, .true. を, そうでなければ .false. を返します. gthstnml が初期設定されていない場合にも .false. が返ります. name に関するデータが存在しない場合にも .false. が返ります.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合, プログラムはエラーを発生させます.

If output of a variable name is valid, .true. is returned, otherwise .false. is returned When gthstnml is not initialized, .false. is returned too. When data correspond to name is not found, .false. is returned too.

If gthstnml is not initialized by "Create" yet, error is occurred.

Alias for HistoryNmlInfoOutputValid

Private Instance methods

GT_HISTORY_NMLINFO_ENTRY
Derived Type :
name :character(TOKEN)
: 変数名. Variable identifier
file :character(STRING)
: ヒストリデータのファイル名. History data filenames
interval_value =>null() :real, pointer
: ヒストリデータの出力間隔の数値. 負の値を与えると, 出力を抑止します.

Numerical value for interval of history data output. Negative values suppresses output.

interval_unit =>null() :character(TOKEN), pointer
: ヒストリデータの出力間隔の単位. Unit for interval of history data output
precision =>null() :character(TOKEN), pointer
: ヒストリデータの精度. Precision of history data
average =>null() :logical, pointer
: 出力データの平均化フラグ. Flag for average of output data.
fileprefix =>null() :character(STRING), pointer
: ヒストリデータのファイル名の接頭詞. Prefixes of history data filenames
output_step_disable = .false. :logical
: output_step 無効化フラグ. "output_step" disable flag
dummy0 :logical
: 8 ビット境界用のダミー変数. Dummy variable for 8 bit boundary
history =>null() :type(GT_HISTORY), pointer
: gt4_history モジュール用構造体. Derived type for "gt4_history" module
next =>null() :type(GT_HISTORY_NMLINFO_ENTRY), pointer
: リスト構造のための変数. A variable for a list structure

出力変数ごとの情報を格納するための構造体です. この構造体はモジュール内で使用されることを想定しているため, モジュール外部からは使用しないでください.

Information about individual output variable is stored in this derived type. It is expected that this derived type is used internally, so do not refer from the outside.

Subroutine :recursive
gthstnml :type(GT_HISTORY_NMLINFO), intent(inout)
name :character(*), intent(in), optional
: 変数名.

先頭の空白は無視されます.

"Data1,Data2" のようにカンマで区切って複数 の変数を指定することも可能です.

Variable identifier.

Blanks at the head of the name are ignored.

Multiple variables can be specified as "Data1,Data2" too. Delimiter is comma.

file :character(*), intent(in), optional
: ヒストリデータのファイル名. History data filenames
interval_value :real, intent(in), optional
: ヒストリデータの出力間隔の数値. 負の値を与えると, 出力を抑止します.

Numerical value for interval of history data output. Negative values suppresses output.

interval_unit :character(*), intent(in), optional
: ヒストリデータの出力間隔の単位. Unit for interval of history data output
precision :character(*), intent(in), optional
: ヒストリデータの精度. Precision of history data
average :logical, intent(in), optional
: 出力データの平均化フラグ. Flag for average of output data.
fileprefix :character(*), intent(in), optional
: ヒストリデータのファイル名の接頭詞. Prefixes of history data filenames
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

変数の出力情報を加えます.

デフォルト値を設定するには, name を与えないか, または name に空白を与えてください. デフォルト値を与える場合, file に与えられる情報は無視されます. fileprefix はデフォルト値に与える場合のみ有効です.

name に変数名が指定され, その際に file が与えられない, または空白が与えられる場合, file には "<name に与えられた文字>.nc" が指定されます.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合, プログラムはエラーを発生させます.

Add output information of a variable.

In order to set default values, specify blank to name or do not specify name. When default values are specified, file is ignored. fileprefix is valid only when default values are specified.

When a variable identifier is specified to name and file is not specified or blanks are specified to file, "<string given to name>.nc" is specified to file.

If gthstnml is not initialized by "Create" yet, error is occurred.

[Source]

  recursive subroutine HistoryNmlInfoAdd( gthstnml, name, file, interval_value, interval_unit, precision, average, fileprefix, err )
    !
    ! 変数の出力情報を加えます. 
    !
    ! デフォルト値を設定するには, *name* を与えないか, または
    ! *name* に空白を与えてください. 
    ! デフォルト値を与える場合, *file* に与えられる情報は無視されます. 
    ! *fileprefix* はデフォルト値に与える場合のみ有効です. 
    !
    ! *name* に変数名が指定され, その際に *file* が与えられない, 
    ! または空白が与えられる場合, *file* には 
    ! "<i><*name* に与えられた文字></i>.nc" が指定されます. 
    !
    ! なお, 与えられた *gthstnml* が Create によって初期設定
    ! されていない場合, プログラムはエラーを発生させます. 
    !
    ! Add output information of a variable.
    ! 
    ! In order to set default values, specify blank to *name* or
    ! do not specify *name*.
    ! When default values are specified, *file* is ignored. 
    ! *fileprefix* is valid only when default values are specified. 
    !
    ! When a variable identifier is specified to *name* and 
    ! *file* is not specified or blanks are specified to *file*,
    ! "<i><string given to *name*></i>.nc" is specified to *file*.
    !
    ! If *gthstnml* is not initialized by "Create" yet, 
    ! error is occurred.
    !
    use dc_trace, only: BeginSub, EndSub, DbgMessage
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_present, only: present_and_not_empty, present_and_true, present_select
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    use dc_date_types, only: DC_DIFFTIME
    use dc_date, only: Create
    use dc_message, only: MessageNotify
    use dc_error, only: StoreError, DC_NOERR, DC_ENOTINIT, DC_EARGLACK, USR_ERRNO, HST_ENOTINDEFINE
    implicit none
    type(GT_HISTORY_NMLINFO), intent(inout):: gthstnml
    character(*), intent(in), optional:: name
                              ! 変数名. 
                              ! 
                              ! 先頭の空白は無視されます. 
                              ! 
                              ! "Data1,Data2" のようにカンマで区切って複数
                              ! の変数を指定することも可能です. 
                              !
                              ! 
                              ! Variable identifier. 
                              ! 
                              ! Blanks at the head of the name are ignored. 
                              ! 
                              ! Multiple variables can be specified 
                              ! as "Data1,Data2" too. Delimiter is comma. 
                              !
                              ! 
    character(*), intent(in), optional:: file
                              ! ヒストリデータのファイル名. 
                              ! History data filenames
    real, intent(in), optional:: interval_value
                              ! ヒストリデータの出力間隔の数値. 
                              ! 負の値を与えると, 出力を抑止します. 
                              ! 
                              ! Numerical value for interval of history data output. 
                              ! Negative values suppresses output.
    character(*), intent(in), optional:: interval_unit
                              ! ヒストリデータの出力間隔の単位. 
                              ! Unit for interval of history data output
    character(*), intent(in), optional:: precision
                              ! ヒストリデータの精度. 
                              ! Precision of history data
    logical, intent(in), optional:: average
                              ! 出力データの平均化フラグ. 
                              ! Flag for average of output data.
    character(*), intent(in), optional:: fileprefix
                              ! ヒストリデータのファイル名の接頭詞. 
                              ! Prefixes of history data filenames
    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    !-----------------------------------
    !  作業変数
    !  Work variables
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: hptr =>null()
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: hptr_last =>null()
    type(DC_DIFFTIME):: interval_time
    character(TOKEN), pointer:: varnames_array(:) =>null()
    integer:: i, vnmax
    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = 'HistoryNmlInfoAdd'
  continue
    call BeginSub( subname, fmt = '@name=%a @file=%a @interval_value=%r @interval_unit=%a @precision=%a @average=%y @fileprefix=%a', r  = (/ present_select(.true., -1.0, interval_value) /), l  = (/ present_and_true(average) /), ca = StoA( present_select(.true., '<no>', name), present_select(.true., '<no>', file), present_select(.true., '<no>', interval_unit), present_select(.true., '<no>', precision), present_select(.true., '<no>', fileprefix) ) )

    stat = DC_NOERR
    cause_c = ''

    !-----------------------------------------------------------------
    !  初期設定のチェック
    !  Check initialization
    !-----------------------------------------------------------------
    if ( .not. gthstnml % initialized ) then
      stat = DC_ENOTINIT
      cause_c = 'GT_HISTORY_NMLINFO'
      goto 999
    end if

    if ( .not. gthstnml % define_mode ) then
      stat = HST_ENOTINDEFINE
      cause_c = 'Add'
      goto 999
    end if

    !-----------------------------------------------------------------
    !  複数の変数を設定する場合
    !  Configure multiple variables
    !-----------------------------------------------------------------
    if ( present_and_not_empty(name) ) then
      if ( index(name, ',') > 0 ) then
        call DbgMessage( 'multiple entries (%c) will be created', c1 = trim(name) )
!!$        if ( present(file) ) call DbgMessage( 'argument @file=%c is ignored', c1 = trim(file) )

        call Split( str = name, sep = ',', carray = varnames_array )        ! (out)
        vnmax = size( varnames_array )

        do i = 1, vnmax
          call Add( gthstnml = gthstnml, name = varnames_array(i), file = file, interval_value = interval_value, interval_unit = interval_unit, precision = precision, average = average, err = err )                        ! (out)
          if ( present_and_true( err ) ) then
            deallocate( varnames_array )
            stat = USR_ERRNO
            goto 999
          end if
        end do
        deallocate( varnames_array )
        goto 999
      end if
    end if

    !-----------------------------------------------------------------
    !  *gthstnml* へ情報を追加.
    !  Add information to *gthstnml*
    !-----------------------------------------------------------------
    if ( .not. present_and_not_empty(name) ) then
      if ( present(interval_value) ) gthstnml % gthstnml_list % interval_value = interval_value
      if ( present(interval_unit)  ) gthstnml % gthstnml_list % interval_unit  = interval_unit 
      if ( present(precision)      ) gthstnml % gthstnml_list % precision      = precision     
      if ( present(average)        ) gthstnml % gthstnml_list % average        = average       
      if ( present(fileprefix)     ) gthstnml % gthstnml_list % fileprefix     = fileprefix    

      hptr => gthstnml % gthstnml_list

    else
      hptr => gthstnml % gthstnml_list
      call ListSearch( gthstnml_list = hptr, name = name )           ! (in)
      if ( .not. associated(hptr) ) then
        call DbgMessage( 'new entry (%c) is created', c1 = trim( adjustl( name ) ) )

        hptr_last => gthstnml % gthstnml_list
        call ListLast( gthstnml_list = hptr_last ) ! (inout)
        allocate( hptr )

        nullify( hptr % next )

        hptr % interval_value => gthstnml % gthstnml_list % interval_value 
        hptr % interval_unit  => gthstnml % gthstnml_list % interval_unit  
        hptr % precision      => gthstnml % gthstnml_list % precision      
        hptr % average        => gthstnml % gthstnml_list % average        
        hptr % fileprefix     => gthstnml % gthstnml_list % fileprefix     

        hptr_last % next => hptr
      else
        call DbgMessage( 'entry (%c) is overwritten', c1 = trim( adjustl( name ) ) )
      end if

      hptr % name  = adjustl( name )
      if ( present_and_not_empty(file) ) then
        hptr % file = file
        nullify(  hptr % fileprefix )
        allocate( hptr % fileprefix )
        hptr % fileprefix = ''
      else
        hptr % file = trim(name) // '.nc'
      end if

      if ( present(interval_value) ) then
        nullify(  hptr % interval_value )
        allocate( hptr % interval_value )
        hptr % interval_value = interval_value 
      end if
      if ( present(interval_unit)  ) then
        nullify(  hptr % interval_unit  )
        allocate( hptr % interval_unit  )
        hptr % interval_unit  = interval_unit  
      end if
      if ( present(precision)      ) then
        nullify(  hptr % precision      )
        allocate( hptr % precision      )
        hptr % precision      = precision      
      end if
      if ( present(average)        ) then
        nullify(  hptr % average        )
        allocate( hptr % average        )
        hptr % average        = average        
      end if

    end if

    !---------------------------------------------------------------
    !  時間の単位のチェック
    !  Check unit of time
    !---------------------------------------------------------------
    call Create( diff = interval_time, value = real( hptr % interval_value, DP ), unit = hptr % interval_unit, err = err )                                  ! (out)
    if ( present_and_true( err ) ) then
      call Delete( gthstnml = gthstnml, name = name )                   ! (in)
      stat = USR_ERRNO
      goto 999
    end if

    nullify( hptr )

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
999 continue
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )
  end subroutine HistoryNmlInfoAdd
Subroutine :
gthstnml :type(GT_HISTORY_NMLINFO), intent(in)
name :character(*), intent(in)
: 変数名. 先頭の空白は無視されます.

Variable identifier. Blanks at the head of the name are ignored.

history :type(GT_HISTORY), pointer
: (out)

gt4_history モジュール用構造体. Derived type for "gt4_history" module

err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

与えられた gt4_history#GT_HISTORY 型のポインタ history に対し, gthstnml 内の name に関する gt4_history#GT_HISTORY 型変数を 結合します. 空状態の history を与えてください.

EndDefine で定義モードから出力モードに 移行した後に呼び出してください. EndDefine を呼ぶ前にこのサブルーチンを使用すると, プログラムはエラーを発生させます.

name に関する情報が見当たらない場合, プログラムはエラーを発生させます. name が空文字の場合にも, プログラムはエラーを発生させます.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合にも, プログラムはエラーを発生させます.

This subroutine associates given "gt4_history#GT_HISTORY" pointer history to "gt4_history#GT_HISTORY" correspond to name in gthstnml. Give null history.

Use after state is changed from define mode to output mode by "EndDefile". If this subroutine is used before "EndDefine" is used, error is occurred.

When data correspond to name is not found, error is occurred. When name is blank, error is occurred too.

If gthstnml is not initialized by "Create" yet, error is occurred.

[Source]

  subroutine HistoryNmlInfoAssocGtHist( gthstnml, name, history, err )
    !
    ! 与えられた gt4_history#GT_HISTORY 型のポインタ *history* に対し, 
    ! *gthstnml* 内の *name* に関する gt4_history#GT_HISTORY 型変数を
    ! 結合します. 
    ! 空状態の *history* を与えてください. 
    !
    ! EndDefine で定義モードから出力モードに
    ! 移行した後に呼び出してください. 
    ! EndDefine を呼ぶ前にこのサブルーチンを使用すると, 
    ! プログラムはエラーを発生させます. 
    !
    ! *name* に関する情報が見当たらない場合, 
    ! プログラムはエラーを発生させます. 
    ! *name* が空文字の場合にも, 
    ! プログラムはエラーを発生させます. 
    !
    ! なお, 与えられた *gthstnml* が Create によって初期設定
    ! されていない場合にも, プログラムはエラーを発生させます. 
    !
    ! This subroutine associates given "gt4_history#GT_HISTORY" 
    ! pointer *history* to 
    ! "gt4_history#GT_HISTORY" correspond to *name* in *gthstnml*. 
    ! Give null *history*. 
    !
    ! Use after state is changed from define mode to
    ! output mode by "EndDefile". 
    ! If this subroutine is used before 
    ! "EndDefine" is used, error is occurred. 
    !
    ! When data correspond to *name* is not found, 
    ! error is occurred.
    ! When *name* is blank, 
    ! error is occurred too.
    ! 
    ! If *gthstnml* is not initialized by "Create" yet, 
    ! error is occurred.
    !
    use dc_trace, only: BeginSub, EndSub
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    use dc_error, only: StoreError, DC_NOERR, DC_ENOTINIT, DC_ENOENTRY, HST_EBADNAME, HST_EINDEFINE
    implicit none
    type(GT_HISTORY_NMLINFO), intent(in):: gthstnml
    character(*), intent(in):: name
                              ! 変数名. 
                              ! 先頭の空白は無視されます. 
                              ! 
                              ! Variable identifier. 
                              ! Blanks at the head of the name are ignored. 
    type(GT_HISTORY), pointer:: history
                              ! (out)
                              ! 
                              ! gt4_history モジュール用構造体. 
                              ! Derived type for "gt4_history" module
    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    !-----------------------------------
    !  作業変数
    !  Work variables
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: hptr =>null()
    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = 'HistoryNmlInfoAssocGtHist'
  continue
    call BeginSub( subname )
    stat = DC_NOERR
    cause_c = ''

    !-----------------------------------------------------------------
    !  初期設定のチェック
    !  Check initialization
    !-----------------------------------------------------------------
    if ( .not. gthstnml % initialized ) then
      stat = DC_ENOTINIT
      cause_c = 'GT_HISTORY_NMLINFO'
      goto 999
    end if

    if ( trim( name ) == '' ) then
      stat = HST_EBADNAME
      cause_c = ''
      goto 999
    end if

    if ( gthstnml % define_mode ) then
      stat = HST_EINDEFINE
      cause_c = 'AssocGtHist'
      goto 999
    end if

    !-----------------------------------------------------------------
    !  *gthstnml* 内から, *name* に関する history を探査.
    !  Search "history" correspond to *name* in *gthstnml*
    !-----------------------------------------------------------------
    hptr => gthstnml % gthstnml_list
    call ListSearch( gthstnml_list = hptr, name = name )           ! (in)

    if ( .not. associated( hptr ) ) then
      stat = DC_ENOENTRY
      cause_c = adjustl( name )
      goto 999
    end if

    nullify( history )
    history => hptr % history

    nullify( hptr )

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
999 continue
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )
  end subroutine HistoryNmlInfoAssocGtHist
Subroutine :
gthstnml :type(GT_HISTORY_NMLINFO), intent(inout)
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

GT_HISTORY_NMLINFO 型の変数の終了処理を行います.

このサブルーチンを使用する前に, gthstnml に格納されている gt4_history#GT_HISTORY 型の全ての変数に対して, gt4_history#HistoryClose を用いて終了処理を行ってください. 終了処理されていないものがある場合, プログラムはエラーを発生させます.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合, プログラムはエラーを発生させます.

Deconstructor of "GT_HISTORY_NMLINFO".

Terminate all "gt4_history#GT_HISTORY" variables in gthstnml by "gt4_history#HistoryClose" before this subroutine is used. If unterminated variables remain, error is occurred.

Note that if gthstnml is not initialized by "Create" yet, error is occurred.

[Source]

  subroutine HistoryNmlInfoClose( gthstnml, err )
    !
    ! GT_HISTORY_NMLINFO 型の変数の終了処理を行います. 
    !
    ! このサブルーチンを使用する前に, *gthstnml* に格納されている
    ! gt4_history#GT_HISTORY 型の全ての変数に対して, 
    ! gt4_history#HistoryClose を用いて終了処理を行ってください. 
    ! 終了処理されていないものがある場合, 
    ! プログラムはエラーを発生させます. 
    !
    ! なお, 与えられた *gthstnml* が Create によって初期設定
    ! されていない場合, プログラムはエラーを発生させます. 
    !
    ! Deconstructor of "GT_HISTORY_NMLINFO". 
    !
    ! Terminate all "gt4_history#GT_HISTORY" variables in *gthstnml*
    ! by "gt4_history#HistoryClose" before this subroutine is used. 
    ! If unterminated variables remain, 
    ! error is occurred. 
    !
    ! Note that if *gthstnml* is not initialized by "Create" yet, 
    ! error is occurred. 
    !
    use dc_trace, only: BeginSub, EndSub, DbgMessage
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    use dc_error, only: StoreError, DC_NOERR, DC_ENOTINIT, HST_ENOTTERMGTHIST
    use gt4_history, only: initialized
    implicit none
    type(GT_HISTORY_NMLINFO), intent(inout):: gthstnml
    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    !-----------------------------------
    !  作業変数
    !  Work variables
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: hptr =>null()
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: hptr_prev =>null()
    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = 'HistoryNmlInfoClose'
  continue
    call BeginSub( subname )
    stat = DC_NOERR
    cause_c = ''

    !-----------------------------------------------------------------
    !  初期設定のチェック
    !  Check initialization
    !-----------------------------------------------------------------
    if ( .not. gthstnml % initialized ) then
      stat = DC_ENOTINIT
      cause_c = 'GT_HISTORY_NMLINFO'
      goto 999
    end if

    !-----------------------------------------------------------------
    !  "GT_HISTORY_NMLINFO" の設定の消去
    !  Clear the settings for "GT_HISTORY_NMLINFO"
    !-----------------------------------------------------------------
    do 
      hptr => gthstnml % gthstnml_list
      call ListLast( gthstnml_list = hptr, previous = hptr_prev )             ! (out)
      call DbgMessage( 'remove entry (%c)', c1 = trim(hptr % name) )
      if ( trim( hptr % name ) == '' ) exit
      if ( .not. gthstnml % define_mode ) then
        if ( initialized( hptr % history ) ) then
          stat = HST_ENOTTERMGTHIST
          cause_c = hptr % name
          goto 999
        end if
      end if
      deallocate( hptr )
      nullify( hptr_prev % next )
    end do
    deallocate( gthstnml % gthstnml_list )

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
    gthstnml % initialized = .false.
    gthstnml % define_mode = .true.
999 continue
    nullify( hptr )
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )
  end subroutine HistoryNmlInfoClose
Subroutine :
gthstnml :type(GT_HISTORY_NMLINFO), intent(inout)
interval_value :real, intent(in), optional
: ヒストリデータの出力間隔の数値. 負の値を与えると, 出力を抑止します.

Numerical value for interval of history data output. Negative values suppresses output.

interval_unit :character(*), intent(in), optional
: ヒストリデータの出力間隔の単位. Unit for interval of history data output
precision :character(*), intent(in), optional
: ヒストリデータの精度. Precision of history data
average :logical, intent(in), optional
: 出力データの平均化フラグ. Flag for average of output data.
fileprefix :character(*), intent(in), optional
: ヒストリデータのファイル名の接頭詞. Prefixes of history data filenames
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

GT_HISTORY_NMLINFO 型の変数の初期設定を行います. 他のサブルーチンを使用する前に必ずこのサブルーチンによって GT_HISTORY_NMLINFO 型の変数を初期設定してください.

interval_value, interval_unit, precision, average はデフォルト値として設定されます. fileprefix は各変数の出力ファイル名の接頭詞として 使用されます.

なお, 与えられた gthstnml が既に初期設定されている場合, プログラムはエラーを発生させます.

Constructor of "GT_HISTORY_NMLINFO". Initialize gthstnml by this subroutine, before other procedures are used,

interval_value, interval_unit, precision, average are set as default values. fileprefix is used as prefixes of output filenames of each variable.

Note that if gthstnml is already initialized by this procedure, error is occurred.

[Source]

  subroutine HistoryNmlInfoCreate( gthstnml, interval_value, interval_unit, precision, average, fileprefix, err )
    !
    ! GT_HISTORY_NMLINFO 型の変数の初期設定を行います. 
    ! 他のサブルーチンを使用する前に必ずこのサブルーチンによって 
    ! GT_HISTORY_NMLINFO 型の変数を初期設定してください. 
    !
    ! *interval_value*, 
    ! *interval_unit*, 
    ! *precision*, 
    ! *average*
    ! はデフォルト値として設定されます. 
    ! *fileprefix* は各変数の出力ファイル名の接頭詞として
    ! 使用されます. 
    !
    ! なお, 与えられた *gthstnml* が既に初期設定されている場合, 
    ! プログラムはエラーを発生させます. 
    !
    ! Constructor of "GT_HISTORY_NMLINFO". 
    ! Initialize *gthstnml* by this subroutine, 
    ! before other procedures are used, 
    !
    ! *interval_value*, 
    ! *interval_unit*, 
    ! *precision*, 
    ! *average* 
    ! are set as default values. 
    ! *fileprefix* is used as prefixes of output filenames of 
    ! each variable. 
    !
    ! Note that if *gthstnml* is already initialized 
    ! by this procedure, error is occurred. 
    !
    use dc_trace, only: BeginSub, EndSub
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    use dc_present, only: present_and_not_empty, present_and_true, present_select
    use dc_message, only: MessageNotify
    use dc_error, only: StoreError, DC_NOERR, DC_EALREADYINIT, DC_EARGLACK, DC_ENEGATIVE, DC_ENOFILEREAD, USR_ERRNO
    use dc_date_types, only: DC_DIFFTIME
    use dc_date, only: Create
    implicit none
    type(GT_HISTORY_NMLINFO), intent(inout):: gthstnml
    real, intent(in), optional:: interval_value
                              ! ヒストリデータの出力間隔の数値. 
                              ! 負の値を与えると, 出力を抑止します. 
                              ! 
                              ! Numerical value for interval of history data output. 
                              ! Negative values suppresses output.
    character(*), intent(in), optional:: interval_unit
                              ! ヒストリデータの出力間隔の単位. 
                              ! Unit for interval of history data output
    character(*), intent(in), optional:: precision
                              ! ヒストリデータの精度. 
                              ! Precision of history data
    logical, intent(in), optional:: average
                              ! 出力データの平均化フラグ. 
                              ! Flag for average of output data.
    character(*), intent(in), optional:: fileprefix
                              ! ヒストリデータのファイル名の接頭詞. 
                              ! Prefixes of history data filenames
    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    !-----------------------------------
    !  作業変数
    !  Work variables
    type(DC_DIFFTIME):: interval_time
    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = 'HistoryNmlInfoCreate'
  continue
    call BeginSub( subname, fmt = '@interval_value=%r @interval_unit=%c @precision=%c @average=%y @fileprefix=%c', r  = (/ present_select(.true., -1.0, interval_value) /), c1 = trim( present_select(.true., '<no>', interval_unit) ), c2 = trim( present_select(.true., '<no>', precision) ), l  = (/ present_and_true(average) /), c3 = trim( present_select(.true., '<no>', fileprefix) ), version = version )
    stat = DC_NOERR
    cause_c = ''

    !-----------------------------------------------------------------
    !  初期設定のチェック
    !  Check initialization
    !-----------------------------------------------------------------
    if ( gthstnml % initialized ) then
      stat = DC_EALREADYINIT
      cause_c = 'GT_HISTORY_NMLINFO'
      goto 999
    end if

    !-----------------------------------------------------------------
    !  割付
    !  Allocate
    !-----------------------------------------------------------------
    allocate( gthstnml % gthstnml_list )
    nullify( gthstnml % gthstnml_list % next )

    !-----------------------------------------------------------------
    !  デフォルト値の設定
    !  Configure default values
    !-----------------------------------------------------------------
    gthstnml % gthstnml_list % name = ''
    gthstnml % gthstnml_list % file = ''

    allocate( gthstnml % gthstnml_list % interval_value )
    allocate( gthstnml % gthstnml_list % interval_unit  )
    allocate( gthstnml % gthstnml_list % precision      )
    allocate( gthstnml % gthstnml_list % average        )
    allocate( gthstnml % gthstnml_list % fileprefix    )
    gthstnml % gthstnml_list % interval_value = -1.0
    gthstnml % gthstnml_list % interval_unit  = 'sec'
    gthstnml % gthstnml_list % precision      = 'float'
    gthstnml % gthstnml_list % average        = .false.
    gthstnml % gthstnml_list % fileprefix     = ''

    if ( present(interval_value) ) gthstnml % gthstnml_list % interval_value = interval_value
    if ( present(interval_unit)  ) gthstnml % gthstnml_list % interval_unit  = interval_unit 
    if ( present(precision)      ) gthstnml % gthstnml_list % precision      = precision     

    if ( present(average)        ) gthstnml % gthstnml_list % average        = average       
    if ( present(fileprefix)     ) gthstnml % gthstnml_list % fileprefix     = fileprefix    

    !-----------------------------------------------------------------
    !  時間の単位のチェック
    !  Check unit of time
    !-----------------------------------------------------------------
    call Create( diff = interval_time, value = real( gthstnml % gthstnml_list % interval_value, DP ), unit = gthstnml % gthstnml_list % interval_unit, err = err )                                                ! (out)
    if ( present_and_true( err ) ) then
      stat = USR_ERRNO
      goto 999
    end if

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
    gthstnml % initialized = .true.
    gthstnml % define_mode = .true.
999 continue
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )
  end subroutine HistoryNmlInfoCreate
Function :
result :logical
gthstnml :type(GT_HISTORY_NMLINFO), intent(in)

gthstnml が定義モードであれば .true. が, 定義モードでなければ .false. が返ります.

If gthstnml is define mode, .true. is returned. If gthstnml is not define mode, .false. is returned.

[Source]

  logical function HistoryNmlInfoDefineMode( gthstnml ) result(result)
    !
    ! *gthstnml* が定義モードであれば .true. が, 
    ! 定義モードでなければ .false. が返ります. 
    !
    ! If *gthstnml* is define mode, .true. is returned. 
    ! If *gthstnml* is not define mode, .false. is returned. 
    !
    implicit none
    type(GT_HISTORY_NMLINFO), intent(in):: gthstnml
  continue
    result = gthstnml % define_mode
  end function HistoryNmlInfoDefineMode
Subroutine :recursive
gthstnml :type(GT_HISTORY_NMLINFO), intent(inout)
name :character(*), intent(in)
: 変数名.

先頭の空白は無視されます.

"Data1,Data2" のようにカンマで区切って複数 の変数を指定することが可能です.

Variable identifier.

Blanks at the head of the name are ignored.

Multiple variables can be specified as "Data1,Data2". Delimiter is comma.

err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

変数の出力情報を削除します.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合, プログラムはエラーを発生させます.

Delete output information of a variable.

If gthstnml is not initialized by "Create" yet, error is occurred.

[Source]

  recursive subroutine HistoryNmlInfoDelete( gthstnml, name, err )
    !
    ! 変数の出力情報を削除します. 
    !
    ! なお, 与えられた *gthstnml* が Create によって初期設定
    ! されていない場合, プログラムはエラーを発生させます. 
    !
    ! Delete output information of a variable.
    ! 
    ! If *gthstnml* is not initialized by "Create" yet, 
    ! error is occurred.
    !
    use dc_trace, only: BeginSub, EndSub, DbgMessage
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_present, only: present_and_not_empty, present_and_true
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    use dc_error, only: StoreError, DC_NOERR, DC_ENOTINIT, DC_EARGLACK, USR_ERRNO, HST_ENOTINDEFINE
    implicit none
    type(GT_HISTORY_NMLINFO), intent(inout):: gthstnml
    character(*), intent(in):: name
                              ! 変数名. 
                              ! 
                              ! 先頭の空白は無視されます. 
                              ! 
                              ! "Data1,Data2" のようにカンマで区切って複数
                              ! の変数を指定することが可能です. 
                              ! 
                              ! Variable identifier. 
                              ! 
                              ! Blanks at the head of the name are ignored. 
                              ! 
                              ! Multiple variables can be specified 
                              ! as "Data1,Data2". Delimiter is comma. 
                              ! 
    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    !-----------------------------------
    !  作業変数
    !  Work variables
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: hptr =>null()
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: hptr_prev =>null()
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: hptr_next =>null()
    character(TOKEN), pointer:: varnames_array(:) =>null()
    integer:: i, vnmax
    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = 'HistoryNmlInfoDelete'
  continue
    call BeginSub( subname, fmt = '@name=%c', c1 = trim( name ) )
    stat = DC_NOERR
    cause_c = ''

    !-----------------------------------------------------------------
    !  初期設定のチェック
    !  Check initialization
    !-----------------------------------------------------------------
    if ( .not. gthstnml % initialized ) then
      stat = DC_ENOTINIT
      cause_c = 'GT_HISTORY_NMLINFO'
      goto 999
    end if

    if ( .not. gthstnml % define_mode ) then
      stat = HST_ENOTINDEFINE
      cause_c = 'Delete'
      goto 999
    end if

    !-----------------------------------------------------------------
    !  複数の変数を削除する場合
    !  Delete multiple variables
    !-----------------------------------------------------------------
    if ( present_and_not_empty(name) ) then
      if ( index(name, ',') > 0 ) then
        call DbgMessage( 'multiple entries (%c) will be deleted', c1 = trim(name) )
        call Split( str = name, sep = ',', carray = varnames_array )        ! (out)
        vnmax = size( varnames_array )

        do i = 1, vnmax
          call Delete( gthstnml = gthstnml, name = varnames_array(i), err = err )                        ! (out)
          if ( present_and_true( err ) ) then
            deallocate( varnames_array )
            stat = USR_ERRNO
            goto 999
          end if
        end do
        deallocate( varnames_array )
        goto 999
      end if
    end if

    !-----------------------------------------------------------------
    !  *gthstnml* の情報を削除.
    !  Delete information in *gthstnml*
    !-----------------------------------------------------------------
    hptr => gthstnml % gthstnml_list
    call ListSearch( gthstnml_list = hptr, name = name, previous = hptr_prev, next = hptr_next )      ! (out)

    if ( .not. associated( hptr ) ) goto 999
    if ( ( trim(hptr % name) /= '' ) .and. associated( hptr_prev ) ) then
      call DbgMessage( 'entry (%c) is deleted', c1 = trim( adjustl( name ) ) )
      hptr_prev % next => hptr_next
      deallocate( hptr )
    end if

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
999 continue
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )
  end subroutine HistoryNmlInfoDelete
Subroutine :
gthstnml :type(GT_HISTORY_NMLINFO), intent(inout)
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

定義モードから出力モードに移行し, gthstnml に設定した情報を確定します. AssocGTHist サブルーチンを呼び出す前に, 必ずこのサブルーチンを呼び出してください. このサブルーチンを呼んだ後に Add, Delete, ResetDefault を呼ぶとプログラムはエラーを発生させます.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合にも, プログラムはエラーを発生させます.

Change from define mode to output mode, and determine information configured in gthstnml. Use this subroutine before "AssocGTHist" is used. If "Add", "Delete", "ResetDefault" are used after this subroutine is used, error is occurred.

If gthstnml is not initialized by "Create" yet, error is occurred.

[Source]

  subroutine HistoryNmlInfoEndDefine( gthstnml, err )
    !
    ! 定義モードから出力モードに移行し, 
    ! *gthstnml* に設定した情報を確定します. 
    ! AssocGTHist サブルーチンを呼び出す前に, 
    ! 必ずこのサブルーチンを呼び出してください. 
    ! このサブルーチンを呼んだ後に Add, Delete, ResetDefault 
    ! を呼ぶとプログラムはエラーを発生させます. 
    !
    ! なお, 与えられた *gthstnml* が Create によって初期設定
    ! されていない場合にも, プログラムはエラーを発生させます. 
    !
    ! Change from define mode to output mode, 
    ! and determine information configured in *gthstnml*. 
    ! Use this subroutine before "AssocGTHist" is used. 
    ! If "Add", "Delete", "ResetDefault" are used after 
    ! this subroutine is used, error is occurred. 
    !
    ! If *gthstnml* is not initialized by "Create" yet, 
    ! error is occurred.
    !
    use dc_trace, only: BeginSub, EndSub
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    use dc_hash, only: HASH, Put, Get, Rewind, Next, number
    use dc_error, only: StoreError, DC_NOERR, DC_ENOTINIT, HST_ENOTINDEFINE, HST_EINTFILE
    use dc_message, only: MessageNotify
    implicit none
    type(GT_HISTORY_NMLINFO), intent(inout):: gthstnml
    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    !-----------------------------------
    !  ファイルの重複をチェックするための変数
    !  Variables for checking duplication of files
    type(HASH):: opened_files
    character(STRING):: opname, opfile
    logical:: end

    !-----------------------------------
    !  作業変数
    !  Work variables
    character(STRING):: fullfilename
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: hptr =>null()
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: hptr_prev =>null()
    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = 'HistoryNmlInfoEndDefine'
  continue
    call BeginSub( subname )
    stat = DC_NOERR
    cause_c = ''

    !-----------------------------------------------------------------
    !  初期設定のチェック
    !  Check initialization
    !-----------------------------------------------------------------
    if ( .not. gthstnml % initialized ) then
      stat = DC_ENOTINIT
      cause_c = 'GT_HISTORY_NMLINFO'
      goto 999
    end if

    if ( .not. gthstnml % define_mode ) then
      stat = HST_ENOTINDEFINE
      cause_c = 'EndDefine'
      goto 999
    end if

    !-----------------------------------------------------------------
    !  gt4_history#GT_HISTORY 変数の割付
    !  Allocate "gt4_history#GT_HISTORY" variables
    !-----------------------------------------------------------------
    hptr => gthstnml % gthstnml_list
    allocate( hptr % history )
    WholeLoop : do while ( associated( hptr % next ) )
      call ListNext( gthstnml_list = hptr ) ! (inout)
      if ( trim(hptr % name) == '' .or. trim(hptr % file) == '' ) cycle WholeLoop

      fullfilename = trim( hptr % fileprefix ) // hptr % file

      !---------------------------------------------------------------
      !  以前に同一ファイル名の gt4_history#GT_HISTORY 変数がある場合, そちらに結合
      !  If "gt4_history#GT_HISTORY" that has same filename exist already, associate to it
      !---------------------------------------------------------------
      nullify( hptr_prev )
      call Rewind(opened_files) ! (inout)
      SearchLoop : do
        call Next( opened_files, opname, opfile, end )  ! (out)
        if ( end ) exit SearchLoop
        if ( trim(opfile) /= trim(fullfilename) ) cycle SearchLoop
        hptr_prev => gthstnml % gthstnml_list
        call ListSearch( gthstnml_list = hptr_prev, name = opname )              ! (in)
        if ( .not. associated( hptr_prev ) ) cycle SearchLoop
        !
        ! interval_value, interval_unit の同一性をチェック
        ! Check consistency of "interval_value", "interval_unit"
        !
        if ( hptr % interval_value /= hptr_prev % interval_value ) then
          call MessageNotify( 'W', subname, '@interval_value=%r (var=%a) and @interval_value=%r (var=%a) are applied to a file "%a"', r = (/hptr % interval_value, hptr_prev % interval_value/), ca = StoA(hptr % name, hptr_prev % name, fullfilename) )
          stat = HST_EINTFILE
          cause_c = fullfilename
          goto 999
        elseif ( hptr % interval_unit /= hptr_prev % interval_unit ) then
          call MessageNotify( 'W', subname, '@interval_unit=%a (var=%a) and @interval_unit=%a (var=%a) are applied to a file "%a"', ca = StoA(hptr % interval_unit, hptr % name, hptr_prev % interval_unit, hptr_prev % name, fullfilename) )
          stat = HST_EINTFILE
          cause_c = fullfilename
          goto 999
        end if
        hptr % history => hptr_prev % history
        exit SearchLoop
      end do SearchLoop

      !---------------------------------------------------------------
      !  新規に割付
      !  Allocate newly
      !---------------------------------------------------------------
      if ( .not. associated( hptr % history ) ) then
        allocate( hptr % history )
      end if

      !---------------------------------------------------------------
      !  割り付けられた名前とファイル名を登録
      !  Regist allocated name and filename
      !---------------------------------------------------------------
      call Put( opened_files, hptr % name, fullfilename ) ! (in)

    end do WholeLoop

    nullify( hptr )
    nullify( hptr_prev )

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
    gthstnml % define_mode = .false.
999 continue
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )
  end subroutine HistoryNmlInfoEndDefine
Function :
result :logical
gthstnml :type(GT_HISTORY_NMLINFO), intent(in)

gthstnml が初期設定されている場合には .true. が, 初期設定されていない場合には .false. が返ります.

If gthstnml is initialized, .true. is returned. If gthstnml is not initialized, .false. is returned.

[Source]

  logical function HistoryNmlInfoInitialized( gthstnml ) result(result)
    !
    ! *gthstnml* が初期設定されている場合には .true. が, 
    ! 初期設定されていない場合には .false. が返ります. 
    !
    ! If *gthstnml* is initialized, .true. is returned. 
    ! If *gthstnml* is not initialized, .false. is returned. 
    !
    implicit none
    type(GT_HISTORY_NMLINFO), intent(in):: gthstnml
  continue
    result = gthstnml % initialized
  end function HistoryNmlInfoInitialized
Subroutine :
gthstnml :type(GT_HISTORY_NMLINFO), intent(in)
name :character(*), intent(in), optional
: 変数名. 先頭の空白は無視されます.

Variable identifier. Blanks at the head of the name are ignored.

file :character(*), intent(out), optional
: ヒストリデータのファイル名. History data filenames
interval_value :real, intent(out), optional
: ヒストリデータの出力間隔の数値. 負の値を与えると, 出力を抑止します.

Numerical value for interval of history data output. Negative values suppresses output.

interval_unit :character(*), intent(out), optional
: ヒストリデータの出力間隔の単位. Unit for interval of history data output
precision :character(*), intent(out), optional
: ヒストリデータの精度. Precision of history data
average :logical, intent(out), optional
: 出力データの平均化フラグ. Flag for average of output data.
fileprefix :character(*), intent(out), optional
: ヒストリデータのファイル名の接頭詞. Prefixes of history data filenames
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

変数の出力情報を取得します.

デフォルト値を取得するには, name を与えないか, または name に空白を与えてください.

name に関するデータが存在しない場合, エラーを発生させます.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合, プログラムはエラーを発生させます.

Inquire output information of a variable.

If data correspond to name is not found, error is occurred.

If gthstnml is not initialized by "Create" yet, error is occurred.

[Source]

  subroutine HistoryNmlInfoInquire( gthstnml, name, file, interval_value, interval_unit, precision, average, fileprefix, err )
    !
    ! 変数の出力情報を取得します. 
    !
    ! デフォルト値を取得するには, *name* を与えないか, または
    ! *name* に空白を与えてください. 
    !
    ! *name* に関するデータが存在しない場合, エラーを発生させます. 
    !
    ! なお, 与えられた *gthstnml* が Create によって初期設定
    ! されていない場合, プログラムはエラーを発生させます. 
    !
    ! Inquire output information of a variable.
    ! 
    ! If data correspond to *name* is not found, 
    ! error is occurred. 
    !
    ! If *gthstnml* is not initialized by "Create" yet, 
    ! error is occurred.
    !
    use dc_trace, only: BeginSub, EndSub
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_present, only: present_and_not_empty, present_and_true
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    use dc_error, only: StoreError, DC_NOERR, DC_ENOTINIT, DC_EARGLACK, DC_ENOENTRY
    implicit none
    type(GT_HISTORY_NMLINFO), intent(in):: gthstnml
    character(*), intent(in), optional:: name
                              ! 変数名. 
                              ! 先頭の空白は無視されます. 
                              ! 
                              ! Variable identifier. 
                              ! Blanks at the head of the name are ignored. 
    character(*), intent(out), optional:: file
                              ! ヒストリデータのファイル名. 
                              ! History data filenames
    real, intent(out), optional:: interval_value
                              ! ヒストリデータの出力間隔の数値. 
                              ! 負の値を与えると, 出力を抑止します. 
                              ! 
                              ! Numerical value for interval of history data output. 
                              ! Negative values suppresses output.
    character(*), intent(out), optional:: interval_unit
                              ! ヒストリデータの出力間隔の単位. 
                              ! Unit for interval of history data output
    character(*), intent(out), optional:: precision
                              ! ヒストリデータの精度. 
                              ! Precision of history data
    logical, intent(out), optional:: average
                              ! 出力データの平均化フラグ. 
                              ! Flag for average of output data.
    character(*), intent(out), optional:: fileprefix
                              ! ヒストリデータのファイル名の接頭詞. 
                              ! Prefixes of history data filenames

    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    !-----------------------------------
    !  作業変数
    !  Work variables
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: hptr =>null()
    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = 'HistoryNmlInfoInquire'
  continue
    call BeginSub( subname )
    stat = DC_NOERR
    cause_c = ''

    !-----------------------------------------------------------------
    !  初期設定のチェック
    !  Check initialization
    !-----------------------------------------------------------------
    if ( .not. gthstnml % initialized ) then
      stat = DC_ENOTINIT
      cause_c = 'GT_HISTORY_NMLINFO'
      goto 999
    end if

    !-----------------------------------------------------------------
    !  *gthstnml* 内から, *name* に関する情報を探査.
    !  Search information correspond to *name* in *gthstnml*
    !-----------------------------------------------------------------
    hptr => gthstnml % gthstnml_list
    call ListSearch( gthstnml_list = hptr, name = name )           ! (in)

    if ( .not. associated( hptr ) ) then
      stat = DC_ENOENTRY
      cause_c = adjustl( name )
      goto 999
    end if

    if ( hptr % name == '' ) then
      if ( present(file)           ) file           = ''
    else
      if ( present(file)           ) file           = trim( hptr % fileprefix ) // hptr % file
    end if
    if ( present(interval_value) ) interval_value = hptr % interval_value
    if ( present(interval_unit)  ) interval_unit  = hptr % interval_unit 
    if ( present(precision)      ) precision      = hptr % precision     
    if ( present(average)        ) average        = hptr % average       
    if ( present(fileprefix)     ) fileprefix     = hptr % fileprefix    

    nullify( hptr )

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
999 continue
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )
  end subroutine HistoryNmlInfoInquire
Subroutine :
gthstnml_list :type(GT_HISTORY_NMLINFO_ENTRY), pointer
: (inout)
previous :type(GT_HISTORY_NMLINFO_ENTRY), pointer, optional
: (out)
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

リスト構造である gthstnml_list (GT_HISTORY_NMLINFO_ENTRY 型) を受け取り, 最後のエントリに再結合して返します. gthstnml_list が始めから空の場合には空状態を返します.

previous が与えられる場合, 当該エントリの一つ前の エントリに結合します.

gthstnml_list (type "GT_HISTORY_NMLINFO_ENTRY") that is a list structure is recieved, and gthstnml_list is reassociated to last entry, and returned. If gthstnml_list is null from the beginning, null is returned.

If previous is given, an entry previous to the above entry is associated.

[Source]

  subroutine HistoryNmlInfoListLast( gthstnml_list, previous, err )
    !
    ! リスト構造である *gthstnml_list* (GT_HISTORY_NMLINFO_ENTRY 型) を受け取り, 
    ! 最後のエントリに再結合して返します. 
    ! *gthstnml_list* が始めから空の場合には空状態を返します. 
    !
    ! *previous* が与えられる場合, 当該エントリの一つ前の
    ! エントリに結合します. 
    !
    ! *gthstnml_list* (type "GT_HISTORY_NMLINFO_ENTRY") that is a list structure 
    ! is recieved, and *gthstnml_list* is reassociated to 
    ! last entry, and returned. 
    ! If *gthstnml_list* is null from the beginning, null is returned.
    !
    ! If *previous* is given, an entry previous to the above entry
    ! is associated. 
    !
    use dc_trace, only: BeginSub, EndSub
    use dc_error, only: StoreError, DC_NOERR
    implicit none
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: gthstnml_list
                              ! (inout)
    type(GT_HISTORY_NMLINFO_ENTRY), pointer, optional:: previous
                              ! (out)
    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    !-----------------------------------
    !  作業変数
    !  Work variables
    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = 'HistoryNmlInfoListLast'
  continue
    call BeginSub( subname )
    stat = DC_NOERR
    cause_c = ''

    if ( present( previous ) ) nullify( previous )

    !-----------------------------------------------------------------
    !  空状態の場合は何もしないで返す
    !  If null, return without change
    !-----------------------------------------------------------------
    if ( .not. associated( gthstnml_list ) ) goto 999

    !-----------------------------------------------------------------
    !  最後のエントリの *next* に結合して返す
    !  "*next*" in last entry is associated, and returned
    !-----------------------------------------------------------------
    do while ( associated( gthstnml_list % next ) )
      if ( present( previous ) ) previous => gthstnml_list
      call ListNext( gthstnml_list = gthstnml_list ) ! (inout)
    end do

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
999 continue
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )
  end subroutine HistoryNmlInfoListLast
Subroutine :
gthstnml_list :type(GT_HISTORY_NMLINFO_ENTRY), pointer
: (inout)
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

リスト構造である gthstnml_list (GT_HISTORY_NMLINFO_ENTRY 型) を受け取り, 次のエントリを gthstnml_list に再結合して返します. 次のエントリが無い場合, gthstnml_list の最後のエントリの next (空状態) に接続して返します. gthstnml_list が始めから空の場合には空状態を返します.

gthstnml_list (type "GT_HISTORY_NMLINFO_ENTRY") that is a list structure is recieved, and gthstnml_list is reassociated to next entry, and is returned. If next entry is not found, gthstnml_list is associated to next in last entry (null), and returned. If gthstnml_list is null from the beginning, null is returned.

[Source]

  subroutine HistoryNmlInfoListNext( gthstnml_list, err )
    !
    ! リスト構造である *gthstnml_list* (GT_HISTORY_NMLINFO_ENTRY 型) を受け取り, 
    ! 次のエントリを *gthstnml_list* に再結合して返します. 
    ! 次のエントリが無い場合, *gthstnml_list* の最後のエントリの 
    ! *next* (空状態) に接続して返します. 
    ! *gthstnml_list* が始めから空の場合には空状態を返します. 
    !
    ! *gthstnml_list* (type "GT_HISTORY_NMLINFO_ENTRY") that is a list structure 
    ! is recieved, and *gthstnml_list* is reassociated to next entry, and
    ! is returned. 
    ! If next entry is not found, *gthstnml_list* is associated to 
    ! *next* in last entry (null), and returned. 
    ! If *gthstnml_list* is null from the beginning, null is returned.
    !
    use dc_trace, only: BeginSub, EndSub
    use dc_error, only: StoreError, DC_NOERR
    implicit none
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: gthstnml_list
                              ! (inout)
    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    !-----------------------------------
    !  作業変数
    !  Work variables
    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = 'HistoryNmlInfoListNext'
  continue
    call BeginSub( subname )
    stat = DC_NOERR
    cause_c = ''

    !-----------------------------------------------------------------
    !  空状態の場合は何もしないで返す
    !  If null, return without change
    !-----------------------------------------------------------------
    if ( .not. associated( gthstnml_list ) ) goto 999

    !-----------------------------------------------------------------
    !  次のエントリに結合して返す
    !  Next entry is associated, and returned
    !-----------------------------------------------------------------
    gthstnml_list => gthstnml_list % next

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
999 continue
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )
  end subroutine HistoryNmlInfoListNext
Subroutine :
gthstnml_list :type(GT_HISTORY_NMLINFO_ENTRY), pointer
: (inout)
name :character(*), intent(in)
: 変数名. 先頭の空白は無視されます.

Variable identifier. Blanks at the head of the name are ignored.

previous :type(GT_HISTORY_NMLINFO_ENTRY), pointer, optional
: (out)
next :type(GT_HISTORY_NMLINFO_ENTRY), pointer, optional
: (out)
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

リスト構造である gthstnml_list (GT_HISTORY_NMLINFO_ENTRY 型) を受け取り, 引数 name と同じ値を持つエントリに再結合して返します. 見つからない場合は空状態を返します. gthstnml_list が始めから空の場合には空状態を返します.

previous が与えられる場合, 当該エントリの一つ前の エントリに結合します. 前のエントリが無い場合には 空状態を返します.

next が与えられる場合, 当該エントリの一つ後ろの エントリに結合します. 後ろのエントリが無い場合には 空状態を返します.

gthstnml_list (type "GT_HISTORY_NMLINFO_ENTRY") that is a list structure is recieved, and gthstnml_list is reassociated to the entry that has a value that is same as argument name, and returned. If the entry is not found, null is returned. If gthstnml_list is null from the beginning, null is returned.

If previous is given, an entry previous to the above entry is associated. If previous entries are not found, null is returned.

If next is given, an entry next to the above entry is associated. If next entries are not found, null is returned.

[Source]

  subroutine HistoryNmlInfoListSearch( gthstnml_list, name, previous, next, err )
    !
    ! リスト構造である *gthstnml_list* (GT_HISTORY_NMLINFO_ENTRY 型) を受け取り, 
    ! 引数 *name* と同じ値を持つエントリに再結合して返します. 
    ! 見つからない場合は空状態を返します.
    ! *gthstnml_list* が始めから空の場合には空状態を返します. 
    ! 
    ! *previous* が与えられる場合, 当該エントリの一つ前の
    ! エントリに結合します. 前のエントリが無い場合には
    ! 空状態を返します. 
    !
    ! *next* が与えられる場合, 当該エントリの一つ後ろの
    ! エントリに結合します. 後ろのエントリが無い場合には
    ! 空状態を返します. 
    ! 
    ! *gthstnml_list* (type "GT_HISTORY_NMLINFO_ENTRY") that is a list structure 
    ! is recieved, and *gthstnml_list* is reassociated to 
    ! the entry that has a value that is same as argument *name*, 
    ! and returned. 
    ! If the entry is not found, null is returned. 
    ! If *gthstnml_list* is null from the beginning, null is returned. 
    !
    ! If *previous* is given, an entry previous to the above entry
    ! is associated. If previous entries are not found, 
    ! null is returned.
    !
    ! If *next* is given, an entry next to the above entry
    ! is associated. If next entries are not found, 
    ! null is returned.
    !
    use dc_trace, only: BeginSub, EndSub
    use dc_error, only: StoreError, DC_NOERR
    implicit none
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: gthstnml_list
                              ! (inout)
    character(*), intent(in):: name
                              ! 変数名. 
                              ! 先頭の空白は無視されます. 
                              ! 
                              ! Variable identifier. 
                              ! Blanks at the head of the name are ignored. 
    type(GT_HISTORY_NMLINFO_ENTRY), pointer, optional:: previous
                              ! (out)
    type(GT_HISTORY_NMLINFO_ENTRY), pointer, optional:: next
                              ! (out)
    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    !-----------------------------------
    !  作業変数
    !  Work variables
    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = 'HistoryNmlInfoListSearch'
  continue
    call BeginSub( subname )
    stat = DC_NOERR
    cause_c = ''

    !-----------------------------------------------------------------
    !  空状態の場合は何もしないで返す
    !  If null, return without change
    !-----------------------------------------------------------------
    if ( .not. associated( gthstnml_list ) ) goto 999

    !-----------------------------------------------------------------
    !  引数 *name* と同じ *name* を持つエントリを探査
    !  The entry that has *name* that is same as argument *name* is searched
    !-----------------------------------------------------------------
    if ( present( previous ) ) nullify( previous )
    if ( present( next ) ) nullify( next )
    if ( trim( adjustl( gthstnml_list % name ) ) == trim( adjustl( name ) ) ) then
      if ( present( next ) ) then
        next => gthstnml_list % next
      end if
      goto 999
    end if

    do while ( associated( gthstnml_list ) )
      if ( present( previous ) ) previous => gthstnml_list
      call ListNext( gthstnml_list = gthstnml_list ) ! (inout)
      if ( .not. associated( gthstnml_list ) ) goto 999
      if ( trim( adjustl( gthstnml_list % name ) ) == trim( adjustl( name ) ) ) then
        if ( present( next ) ) then
          next => gthstnml_list % next
        end if
        goto 999
      end if
    end do

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
999 continue
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )
  end subroutine HistoryNmlInfoListSearch
Function :
result :character(STRING)
gthstnml :type(GT_HISTORY_NMLINFO), intent(in)

gthstnml が設定されている変数リストをカンマでつなげて 返します.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合, 空文字が返ります.

List of variables registered in gthstnml is join with ’,’, and returned.

If gthstnml is not initialized by "Create" yet, blank is returned.

[Source]

  character(STRING) function HistoryNmlInfoNames( gthstnml ) result(result)
    !
    ! *gthstnml* が設定されている変数リストをカンマでつなげて
    ! 返します. 
    !
    ! なお, 与えられた *gthstnml* が Create によって初期設定
    ! されていない場合, 空文字が返ります. 
    !
    ! List of variables registered in *gthstnml* is join with ',', 
    ! and returned. 
    ! 
    ! If *gthstnml* is not initialized by "Create" yet, 
    ! blank is returned. 
    !
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    implicit none
    type(GT_HISTORY_NMLINFO), intent(in):: gthstnml

    !-----------------------------------
    !  作業変数
    !  Work variables
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: hptr =>null()
    logical:: first
!!$    character(*), parameter:: subname = 'HistoryNmlInfoNames'
  continue

    result = ''
    first = .true.

    !-----------------------------------------------------------------
    !  初期設定のチェック
    !  Check initialization
    !-----------------------------------------------------------------
    if ( .not. gthstnml % initialized ) goto 999

    !-----------------------------------------------------------------
    !  情報の取り出し
    !  Fetch information
    !-----------------------------------------------------------------
    hptr => gthstnml % gthstnml_list
    do while ( associated( hptr % next ) )
      call ListNext( gthstnml_list = hptr ) ! (inout)
      if ( first ) then
        result = adjustl( hptr % name )
        first = .false.
      else
        result = trim( result ) // ',' // adjustl( hptr % name )
      end if
    end do

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
999 continue
    nullify( hptr )
  end function HistoryNmlInfoNames
Function :
result :logical
gthstnml :type(GT_HISTORY_NMLINFO), intent(in)
name :character(*), intent(in)
: 変数名. 先頭の空白は無視されます.

Variable identifier. Blanks at the head of the name are ignored.

time :type(DC_DIFFTIME), intent(in)
: 現在時刻. Current time

time が変数 name の出力されるタイミングであれば .true. を, そうでなければ .false. を返します. gthstnml が初期設定されていない場合にも .false. が返ります. name に関するデータが存在しない場合にも .false. が返ります.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合, プログラムはエラーを発生させます.

If time is the time that a variable name is output, .true. is returned, otherwise .false. is returned When gthstnml is not initialized, .false. is returned too. When data correspond to name is not found, .false. is returned too.

If gthstnml is not initialized by "Create" yet, error is occurred.

[Source]

  logical function HistoryNmlInfoOutputStep( gthstnml, name, time ) result(result)
    !
    ! *time* が変数 *name* の出力されるタイミングであれば
    ! .true. を, そうでなければ .false. を返します. 
    ! *gthstnml* が初期設定されていない場合にも .false. が返ります. 
    ! *name* に関するデータが存在しない場合にも .false. が返ります. 
    !
    ! なお, 与えられた *gthstnml* が Create によって初期設定
    ! されていない場合, プログラムはエラーを発生させます. 
    !
    ! If *time* is the time that a variable *name* is output,
    ! .true. is returned, otherwise .false. is returned
    ! When *gthstnml* is not initialized, .false. is returned too.
    ! When data correspond to *name* is not found, .false. is returned too.
    ! 
    ! If *gthstnml* is not initialized by "Create" yet, 
    ! error is occurred.
    !
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_present, only: present_and_not_empty, present_and_true
    use dc_date_types, only: DC_DIFFTIME
    use dc_date, only: Create, mod, operator(==), toChar
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    implicit none
    type(GT_HISTORY_NMLINFO), intent(in):: gthstnml
    character(*), intent(in):: name
                              ! 変数名. 
                              ! 先頭の空白は無視されます. 
                              ! 
                              ! Variable identifier. 
                              ! Blanks at the head of the name are ignored. 
    type(DC_DIFFTIME), intent(in):: time
                              ! 現在時刻. Current time

    !-----------------------------------
    !  作業変数
    !  Work variables
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: hptr =>null()
    type(DC_DIFFTIME):: interval_time
!!$    character(*), parameter:: subname = 'HistoryNmlInfoOutputStep'
  continue

    result = .false.

    !-----------------------------------------------------------------
    !  初期設定のチェック
    !  Check initialization
    !-----------------------------------------------------------------
    if ( .not. gthstnml % initialized ) goto 999

    !-----------------------------------------------------------------
    !  情報格納変数への結合
    !  Associate a variable storing information
    !-----------------------------------------------------------------
    hptr => gthstnml % gthstnml_list
    call ListSearch( gthstnml_list = hptr, name = name )           ! (in)

    if ( .not. associated( hptr ) ) goto 999
    if ( hptr % output_step_disable ) goto 999

    !-----------------------------------------------------------------
    !  時刻のチェック
    !  Check time
    !-----------------------------------------------------------------
    if ( .not. hptr % interval_value > 0.0 ) goto 999

    call Create( diff = interval_time, value = real( hptr % interval_value, DP ), unit = hptr % interval_unit )                ! (in)

    if ( mod( time, interval_time ) == 0 ) then
      result = .true.
    end if

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
999 continue
    nullify( hptr )
  end function HistoryNmlInfoOutputStep
Subroutine :
gthstnml :type(GT_HISTORY_NMLINFO), intent(in)
name :character(*), intent(in)
: 変数名. 先頭の空白は無視されます.

Variable identifier. Blanks at the head of the name are ignored.

err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

このサブルーチンを使用すると, name に関して, 以降は output_step が常に .false. を返すようになります.

データ出力間隔を出力の初期設定から変更し, データを出力するたびに時刻を明示的に指定する場合に利用することを 想定しています.

EndDefine で定義モードから出力モードに 移行した後に呼び出してください. EndDefine を呼ぶ前にこのサブルーチンを使用すると, プログラムはエラーを発生させます.

name に関する情報が見当たらない場合, プログラムはエラーを発生させます. name が空文字の場合にも, プログラムはエラーを発生させます.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合にも, プログラムはエラーを発生させます.

After this subroutine is used, "output_step" returns .false. already corresponding to the name.

This subroutine expected to use when interval of data output is changed from initialization of output, and time is specified explicitly whenever data is output.

Use after state is changed from define mode to output mode by "EndDefile". If this subroutine is used before "EndDefine" is used, error is occurred.

When data correspond to name is not found, error is occurred. When name is blank, error is occurred too.

If gthstnml is not initialized by "Create" yet, error is occurred.

[Source]

  subroutine HistoryNmlInfoOutputStepDisable( gthstnml, name, err )
    !
    ! このサブルーチンを使用すると, *name* に関して, 
    ! 以降は output_step が常に .false. を返すようになります. 
    !
    ! データ出力間隔を出力の初期設定から変更し, 
    ! データを出力するたびに時刻を明示的に指定する場合に利用することを
    ! 想定しています. 
    !
    ! EndDefine で定義モードから出力モードに
    ! 移行した後に呼び出してください. 
    ! EndDefine を呼ぶ前にこのサブルーチンを使用すると, 
    ! プログラムはエラーを発生させます. 
    !
    ! *name* に関する情報が見当たらない場合, 
    ! プログラムはエラーを発生させます. 
    ! *name* が空文字の場合にも, 
    ! プログラムはエラーを発生させます. 
    !
    ! なお, 与えられた *gthstnml* が Create によって初期設定
    ! されていない場合にも, プログラムはエラーを発生させます. 
    !
    ! After this subroutine is used, 
    ! "output_step" returns .false. already 
    ! corresponding to the *name*. 
    !
    ! This subroutine expected to use when 
    ! interval of data output is changed from initialization of output, 
    ! and time is specified explicitly whenever data is output. 
    !
    ! Use after state is changed from define mode to
    ! output mode by "EndDefile". 
    ! If this subroutine is used before 
    ! "EndDefine" is used, error is occurred. 
    !
    ! When data correspond to *name* is not found, error is occurred.
    ! When *name* is blank, error is occurred too.
    !
    ! If *gthstnml* is not initialized by "Create" yet, 
    ! error is occurred.
    !
    use dc_trace, only: BeginSub, EndSub
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    use dc_error, only: StoreError, DC_NOERR, DC_ENOTINIT, DC_ENOENTRY, HST_EBADNAME, HST_EINDEFINE
    implicit none
    type(GT_HISTORY_NMLINFO), intent(in):: gthstnml
    character(*), intent(in):: name
                              ! 変数名. 
                              ! 先頭の空白は無視されます. 
                              ! 
                              ! Variable identifier. 
                              ! Blanks at the head of the name are ignored. 
    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    !-----------------------------------
    !  作業変数
    !  Work variables
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: hptr =>null()
    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = 'HistoryNmlInfoOutputStepDisable'
  continue
    call BeginSub( subname )
    stat = DC_NOERR
    cause_c = ''

    !-----------------------------------------------------------------
    !  初期設定のチェック
    !  Check initialization
    !-----------------------------------------------------------------
    if ( .not. gthstnml % initialized ) then
      stat = DC_ENOTINIT
      cause_c = 'GT_HISTORY_NMLINFO'
      goto 999
    end if

    if ( trim( name ) == '' ) then
      stat = HST_EBADNAME
      cause_c = ''
      goto 999
    end if

    if ( gthstnml % define_mode ) then
      stat = HST_EINDEFINE
      cause_c = 'OutputStepDisable'
      goto 999
    end if

    !-----------------------------------------------------------------
    !  *gthstnml* 内から, *name* に関する history を探査.
    !  Search "history" correspond to *name* in *gthstnml*
    !-----------------------------------------------------------------
    hptr => gthstnml % gthstnml_list
    call ListSearch( gthstnml_list = hptr, name = name )           ! (in)

    if ( .not. associated( hptr ) ) then
      stat = DC_ENOENTRY
      cause_c = adjustl( name )
      goto 999
    end if

    hptr % output_step_disable = .true.

    nullify( hptr )

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
999 continue
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )
  end subroutine HistoryNmlInfoOutputStepDisable
Function :
result :logical
gthstnml :type(GT_HISTORY_NMLINFO), intent(in)
name :character(*), intent(in)
: 変数名. 先頭の空白は無視されます.

Variable identifier. Blanks at the head of the name are ignored.

変数 name の出力が有効であれば, .true. を, そうでなければ .false. を返します. gthstnml が初期設定されていない場合にも .false. が返ります. name に関するデータが存在しない場合にも .false. が返ります.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合, プログラムはエラーを発生させます.

If output of a variable name is valid, .true. is returned, otherwise .false. is returned When gthstnml is not initialized, .false. is returned too. When data correspond to name is not found, .false. is returned too.

If gthstnml is not initialized by "Create" yet, error is occurred.

[Source]

  logical function HistoryNmlInfoOutputValid( gthstnml, name ) result(result)
    !
    ! 変数 *name* の出力が有効であれば, 
    ! .true. を, そうでなければ .false. を返します. 
    ! *gthstnml* が初期設定されていない場合にも .false. が返ります. 
    ! *name* に関するデータが存在しない場合にも .false. が返ります. 
    !
    ! なお, 与えられた *gthstnml* が Create によって初期設定
    ! されていない場合, プログラムはエラーを発生させます. 
    !
    ! If output of a variable *name* is valid, 
    ! .true. is returned, otherwise .false. is returned
    ! When *gthstnml* is not initialized, .false. is returned too.
    ! When data correspond to *name* is not found, .false. is returned too.
    ! 
    ! If *gthstnml* is not initialized by "Create" yet, 
    ! error is occurred.
    !
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_present, only: present_and_not_empty, present_and_true
    use dc_date_types, only: DC_DIFFTIME
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    implicit none
    type(GT_HISTORY_NMLINFO), intent(in):: gthstnml
    character(*), intent(in):: name
                              ! 変数名. 
                              ! 先頭の空白は無視されます. 
                              ! 
                              ! Variable identifier. 
                              ! Blanks at the head of the name are ignored. 

    !-----------------------------------
    !  作業変数
    !  Work variables
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: hptr =>null()
!!$    character(*), parameter:: subname = 'HistoryNmlInfoOutputValid'
  continue

    result = .false.

    !-----------------------------------------------------------------
    !  初期設定のチェック
    !  Check initialization
    !-----------------------------------------------------------------
    if ( .not. gthstnml % initialized ) goto 999

    !-----------------------------------------------------------------
    !  情報格納変数への結合
    !  Associate a variable storing information
    !-----------------------------------------------------------------
    hptr => gthstnml % gthstnml_list
    call ListSearch( gthstnml_list = hptr, name = name )           ! (in)

    if ( .not. associated( hptr ) ) goto 999

    !-----------------------------------------------------------------
    !  出力の有効性のチェック
    !  Check validity of output
    !-----------------------------------------------------------------
    if ( hptr % interval_value > 0.0 ) then
      result = .true.
      goto 999
    end if

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
999 continue
    nullify( hptr )
  end function HistoryNmlInfoOutputValid
Subroutine :
gthstnml :type(GT_HISTORY_NMLINFO), intent(in)
unit :integer, intent(in), optional
: 出力先の装置番号. デフォルトの出力先は標準出力.

Unit number for output. Default value is standard output.

indent :character(*), intent(in), optional
: 表示されるメッセージの字下げ.

Indent of displayed messages.

err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

引数 gthstnml に設定されている情報を印字します. デフォルトではメッセージは標準出力に出力されます. unit に装置番号を指定することで, 出力先を変更することが可能です.

Print information of gthstnml. By default messages are output to standard output. Unit number for output can be changed by unit argument.

[Source]

  subroutine HistoryNmlInfoPutLine( gthstnml, unit, indent, err )
    !
    ! 引数 *gthstnml* に設定されている情報を印字します. 
    ! デフォルトではメッセージは標準出力に出力されます. 
    ! *unit* に装置番号を指定することで, 出力先を変更することが可能です. 
    !
    ! Print information of *gthstnml*. 
    ! By default messages are output to standard output. 
    ! Unit number for output can be changed by *unit* argument. 
    !
    use dc_trace, only: BeginSub, EndSub
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    use dc_error, only: StoreError, DC_NOERR, DC_ENOTINIT
    use gt4_history, only: HistoryPutLine
    implicit none
    type(GT_HISTORY_NMLINFO), intent(in):: gthstnml
    integer, intent(in), optional:: unit
                              ! 出力先の装置番号. 
                              ! デフォルトの出力先は標準出力. 
                              !
                              ! Unit number for output. 
                              ! Default value is standard output. 
    character(*), intent(in), optional:: indent
                              ! 表示されるメッセージの字下げ. 
                              !
                              ! Indent of displayed messages. 
    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    !-----------------------------------
    !  作業変数
    !  Work variables
    type(GT_HISTORY_NMLINFO_ENTRY), pointer:: hptr =>null()
    integer:: stat
    character(STRING):: cause_c
    integer:: out_unit
    integer:: indent_len
    character(STRING):: indent_str
    character(*), parameter:: subname = 'HistoryNmlInfoPutLine'
  continue
    call BeginSub( subname )
    stat = DC_NOERR
    cause_c = ''

    !-----------------------------------------------------------------
    !  初期設定のチェック
    !  Check initialization
    !-----------------------------------------------------------------
    if ( present(unit) ) then
      out_unit = unit
    else
      out_unit = STDOUT
    end if

    indent_len = 0
    indent_str = ''
    if ( present(indent) ) then
      if ( len(indent) /= 0 ) then
        indent_len = len(indent)
        indent_str(1:indent_len) = indent
      end if
    end if

    !-----------------------------------------------------------------
    !  "GT_HISTORY_NMLINFO" の設定の印字
    !  Print the settings for "GT_HISTORY_NMLINFO"
    !-----------------------------------------------------------------
    if ( gthstnml % initialized ) then
      call Printf( out_unit, indent_str(1:indent_len) // '#<GT_HISTORY_NMLINFO:: @initialized=%y define_mode=%y', l = (/gthstnml % initialized, gthstnml % define_mode/) )

      hptr => gthstnml % gthstnml_list

      do while ( associated( hptr ) )

        call Printf( out_unit, indent_str(1:indent_len) // ' #<GT_HISTORY_NMLINFO_ENTRY:: @name=%c @file=%c', c1 = trim(hptr % name), c2 = trim(hptr % file) )

        call Printf( out_unit, indent_str(1:indent_len) // '  @interval_value=%r @interval_unit=%c', r = (/hptr % interval_value/), c1 = trim(hptr % interval_unit) )

        call Printf( out_unit, indent_str(1:indent_len) // '  @output_step_disable=%y', l = (/hptr % output_step_disable/) )

        call Printf( out_unit, indent_str(1:indent_len) // '  @precision=%c @average=%y', c1 = trim(hptr % precision), l = (/ hptr % average /) )

        call Printf( out_unit, indent_str(1:indent_len) // '  @fileprefix=%c', c1 = trim(hptr % fileprefix) )

        if ( .not. gthstnml % define_mode ) then
          call Printf( out_unit, indent_str(1:indent_len) // '  @history=' )

          call HistoryPutLine( hptr % history, unit = out_unit, indent = indent_str(1:indent_len) // '   ' )
        end if

        call ListNext( gthstnml_list = hptr ) ! (inout)
      end do

      call Printf( out_unit, indent_str(1:indent_len) // ' >' )

      call Printf( out_unit, indent_str(1:indent_len) // '>' )
    else
      call Printf( out_unit, indent_str(1:indent_len) // '#<GT_HISTORY_NMLINFO:: @initialized=%y>', l = (/gthstnml % initialized/) )
    end if

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
999 continue
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )
  end subroutine HistoryNmlInfoPutLine
Subroutine :
gthstnml :type(GT_HISTORY_NMLINFO), intent(inout)
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

デフォルト値を残し, 登録したデータを削除します.

なお, 与えられた gthstnmlCreate によって初期設定 されていない場合, プログラムはエラーを発生させます.

Stored data is deleted without default settings.

If gthstnml is not initialized by "Create" yet, error is occurred.

[Source]

  subroutine HistoryNmlInfoResetDefault( gthstnml, err )
    !
    ! デフォルト値を残し, 登録したデータを削除します.
    !
    ! なお, 与えられた *gthstnml* が Create によって初期設定
    ! されていない場合, プログラムはエラーを発生させます. 
    !
    ! Stored data is deleted without default settings. 
    !
    ! If *gthstnml* is not initialized by "Create" yet, 
    ! error is occurred.
    !
    use dc_trace, only: BeginSub, EndSub
    use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar
    use dc_types, only: DP, STRING, TOKEN, STDOUT
    use dc_error, only: StoreError, DC_NOERR, DC_ENOTINIT, HST_ENOTINDEFINE
    implicit none
    type(GT_HISTORY_NMLINFO), intent(inout):: gthstnml
    logical, intent(out), optional:: err
                              ! 例外処理用フラグ. 
                              ! デフォルトでは, この手続き内でエラーが
                              ! 生じた場合, プログラムは強制終了します. 
                              ! 引数 *err* が与えられる場合, 
                              ! プログラムは強制終了せず, 代わりに
                              ! *err* に .true. が代入されます. 
                              !
                              ! Exception handling flag. 
                              ! By default, when error occur in 
                              ! this procedure, the program aborts. 
                              ! If this *err* argument is given, 
                              ! .true. is substituted to *err* and 
                              ! the program does not abort. 

    !-----------------------------------
    !  作業変数
    !  Work variables
    character(STRING):: varnames
    character(TOKEN), pointer:: varnames_array(:) =>null()
    integer:: i, vnmax
    integer:: stat
    character(STRING):: cause_c
    character(*), parameter:: subname = 'HistoryNmlInfoResetDefault'
  continue
    call BeginSub( subname )
    stat = DC_NOERR
    cause_c = ''

    !-----------------------------------------------------------------
    !  初期設定のチェック
    !  Check initialization
    !-----------------------------------------------------------------
    if ( .not. gthstnml % initialized ) then
      stat = DC_ENOTINIT
      cause_c = 'GT_HISTORY_NMLINFO'
      goto 999
    end if

    if ( .not. gthstnml % define_mode ) then
      stat = HST_ENOTINDEFINE
      cause_c = 'ResetDefault'
      goto 999
    end if

    !-----------------------------------------------------------------
    !  変数名リストの取得
    !  Get varnames list
    !-----------------------------------------------------------------
    varnames = names( gthstnml )
    call Split( str = varnames, sep = ',', carray = varnames_array )            ! (out)
    vnmax = size( varnames_array )

    do i = 1, vnmax
      call Delete( gthstnml = gthstnml, name = varnames_array(i) )      ! (in)
    end do

    !-----------------------------------------------------------------
    !  終了処理, 例外処理
    !  Termination and Exception handling
    !-----------------------------------------------------------------
999 continue
    call StoreError( stat, subname, err, cause_c )
    call EndSub( subname )
  end subroutine HistoryNmlInfoResetDefault
ListLast( gthstnml_list, [previous], [err] )
Subroutine :
gthstnml_list :type(GT_HISTORY_NMLINFO_ENTRY), pointer
: (inout)
previous :type(GT_HISTORY_NMLINFO_ENTRY), pointer, optional
: (out)
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

リスト構造である gthstnml_list (GT_HISTORY_NMLINFO_ENTRY 型) を受け取り, 最後のエントリに再結合して返します. gthstnml_list が始めから空の場合には空状態を返します.

previous が与えられる場合, 当該エントリの一つ前の エントリに結合します.

gthstnml_list (type "GT_HISTORY_NMLINFO_ENTRY") that is a list structure is recieved, and gthstnml_list is reassociated to last entry, and returned. If gthstnml_list is null from the beginning, null is returned.

If previous is given, an entry previous to the above entry is associated.

Alias for HistoryNmlInfoListLast

ListNext( gthstnml_list, [err] )
Subroutine :
gthstnml_list :type(GT_HISTORY_NMLINFO_ENTRY), pointer
: (inout)
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

リスト構造である gthstnml_list (GT_HISTORY_NMLINFO_ENTRY 型) を受け取り, 次のエントリを gthstnml_list に再結合して返します. 次のエントリが無い場合, gthstnml_list の最後のエントリの next (空状態) に接続して返します. gthstnml_list が始めから空の場合には空状態を返します.

gthstnml_list (type "GT_HISTORY_NMLINFO_ENTRY") that is a list structure is recieved, and gthstnml_list is reassociated to next entry, and is returned. If next entry is not found, gthstnml_list is associated to next in last entry (null), and returned. If gthstnml_list is null from the beginning, null is returned.

Alias for HistoryNmlInfoListNext

ListSearch( gthstnml_list, name, [previous], [next], [err] )
Subroutine :
gthstnml_list :type(GT_HISTORY_NMLINFO_ENTRY), pointer
: (inout)
name :character(*), intent(in)
: 変数名. 先頭の空白は無視されます.

Variable identifier. Blanks at the head of the name are ignored.

previous :type(GT_HISTORY_NMLINFO_ENTRY), pointer, optional
: (out)
next :type(GT_HISTORY_NMLINFO_ENTRY), pointer, optional
: (out)
err :logical, intent(out), optional
: 例外処理用フラグ. デフォルトでは, この手続き内でエラーが 生じた場合, プログラムは強制終了します. 引数 err が与えられる場合, プログラムは強制終了せず, 代わりに err に .true. が代入されます.

Exception handling flag. By default, when error occur in this procedure, the program aborts. If this err argument is given, .true. is substituted to err and the program does not abort.

リスト構造である gthstnml_list (GT_HISTORY_NMLINFO_ENTRY 型) を受け取り, 引数 name と同じ値を持つエントリに再結合して返します. 見つからない場合は空状態を返します. gthstnml_list が始めから空の場合には空状態を返します.

previous が与えられる場合, 当該エントリの一つ前の エントリに結合します. 前のエントリが無い場合には 空状態を返します.

next が与えられる場合, 当該エントリの一つ後ろの エントリに結合します. 後ろのエントリが無い場合には 空状態を返します.

gthstnml_list (type "GT_HISTORY_NMLINFO_ENTRY") that is a list structure is recieved, and gthstnml_list is reassociated to the entry that has a value that is same as argument name, and returned. If the entry is not found, null is returned. If gthstnml_list is null from the beginning, null is returned.

If previous is given, an entry previous to the above entry is associated. If previous entries are not found, null is returned.

If next is given, an entry next to the above entry is associated. If next entries are not found, null is returned.

Alias for HistoryNmlInfoListSearch

version
Constant :
version = ’$Name: gt4f90io-20080210 $’ // ’$Id: gt4_history_nmlinfo.f90,v 1.5 2008/02/10 16:51:15 morikawa Exp $’ :character(*), parameter

[Validate]