Class | dyn_spectral_as83 |
In: |
dynamics/dyn_spectral_as83.f90
|
Note that Japanese and English are described in parallel.
力学過程を演算するモジュールです. 水平離散化にスペクトル法を, 鉛直離散化には Arakawa and Suarez (1983) を用いています.
時間積分法にはリープフロッグスキームを用いています. デフォルトでは, $ Delta t $ を大きくとるために, 重力波項に セミインプリシット法を適用しています. 重力波項をエクスプリシット法によって解くことも可能です. 詳しくは Create および NAMELIST#dyn_as83_nml を参照してください.
This is a dynamical core module. Spectral method (for horizontal) and Arakawa and Suarez (1983) method (for vertical) are used.
Leap-frog scheme is used as time integration. By default, semi-implicit scheme is applied to gravitational terms for extension of $ Delta t $ . Explicit scheme can be applied to gravitational terms. For details, see "Create" or "NAMELIST#dyn_as83_nml".
DynSpAsCreate : | DYNSPAS83 型変数の初期化 |
Dynamics : | 力学過程の演算 |
VorDiv2UV : | 渦度発散を風速に変換 |
UV2VorDiv : | 風速を渦度発散に変換 |
DynSpAsEqualAxes : | DYNSPAS83 型変数に格納されている座標値の確認 |
DynSpAsGetAxes : | 座標値の取得 |
DynSpAsClose : | DYNSPAS83 型変数の終了処理 |
DynSpAsPutLine : | DYNSPAS83 型変数に格納されている情報の印字 |
DynSpAsInitialized : | DYNSPAS83 型変数が初期設定されているか否か |
DynSpAsSetTime : | 時刻の設定 |
———— : | ———— |
DynSpAsCreate : | Constructor of "DYNSPAS83" |
Dynamics : | Calculation of dynamical core |
VorDiv2UV : | Convert vorticity and divergence into wind |
UV2VorDiv : | Convert wind into vorticity and divergence |
DynSpAsEqualAxes : | Confirm data of axes in "DYNSPAS83" |
DynSpAsGetAxes : | Get data of axes |
DynSpAsClose : | Deconstructor of "DYNSPAS83" |
DynSpAsPutLine : | Print information of "DYNSPAS83" |
DynSpAsInitialized : | Check initialization of "DYNSPAS83" |
DynSpAsSetTime : | Configure time |
主プログラムの始めに, DYNSPAS83 型の変数を Create で初期化します. 次に, 時間積分ループ内で渦度, 発散, 温度, 比湿, 地表面気圧 ( $ t $ および $ t-\Delta t $ ) を Dynamics に与えてください. それらの値を用いて力学過程の演算が行われ, $ t+Delta t $ の東西風速, 南北風速, 渦度, 発散, 温度, 比湿, 地表面気圧が得られます. 最後に, DYNSPAS83 型の変数の終了処理を Close にて行います. 主プログラムで座標値が必要となる場合には, GetAxes を用いて 座標値を取得してください.
In the begging of program, initialize "DYNSPAS83" by "Create". Next, in the time integration loop, provide vorticity, divergence, temperature, specific humidity, surface pressure on $ t $ and $ t-\Delta t $ to Dynamics. Then the physical values at $ t+Delta t $ is returned. Finally, terminate "DYNSPAS83" by "Close". If you need data of axes, get that by GetAxes.
Derived Type : | |||
initialized = .false. : | logical
| ||
nmax : | integer
| ||
imax : | integer
| ||
jmax : | integer
| ||
kmax : | integer
| ||
dyn_sp : | type(DYNSP) | ||
dyn_as : | type(DYNAS83) | ||
current_time : | type(DC_DIFFTIME)
| ||
delta_time : | type(DC_DIFFTIME)
| ||
gthstnml : | type(GTHST_NMLINFO)
|
まず, Create で "DYNSPAS83" 型の変数を初期設定して下さい. 初期化された "DYNSPAS83" 型の変数を再度利用する際には, Close によって終了処理を行ってください.
Initialize "DYNSPAS83" variable by "Create" before usage. If you reuse "DYNSPAS83" variable again for another application, terminate by "Close".
Subroutine : | |||
dyn_sp_as : | type(DYNSPAS83), intent(inout) | ||
err : | logical, intent(out), optional
|
DYNSPAS83 型の変数の終了処理を行います. なお, 与えられた dyn_sp_as が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Deconstructor of "DYNSPAS83". Note that if dyn_sp_as is not initialized by "Create" yet, error is occurred.
subroutine DynSpAsClose( dyn_sp_as, err ) ! ! DYNSPAS83 型の変数の終了処理を行います. ! なお, 与えられた *dyn_sp_as* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Deconstructor of "DYNSPAS83". ! Note that if *dyn_sp_as* is not initialized by "Create" yet, ! error is occurred. ! use dc_trace, only: BeginSub, EndSub use dc_types, only: STRING, STDOUT use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT use dyn_spectral, only: Close use dyn_as83, only: Close use gt4_history_nmlinfo, only: HstNmlInfoClose, HstNmlInfoNames, HstNmlInfoAssocGtHist, HstNmlInfoPutLine use gt4_history, only: HistoryClose, HistoryInitialized implicit none type(DYNSPAS83), intent(inout):: dyn_sp_as 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. !----------------------------------- ! ヒストリファイルへのデータ出力設定 ! Configure the settings for history data output character(STRING):: name = '' ! 変数名. Variable identifier character(STRING):: varnames ! 変数名リスト. ! List of variables character(TOKEN), pointer:: varnames_array(:) =>null() ! 変数名リスト配列. ! List of variables (array) integer:: i, vnmax type(GT_HISTORY), pointer:: gthist =>null() ! gt4_history モジュール用構造体. ! Derived type for "gt4_history" module !----------------------------------- ! 作業変数 ! Work variables integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpAsClose' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp_as % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSPAS83' goto 999 end if !----------------------------------------------------------------- ! スペクトル法と Arakawa and Suarez (1983) の設定の消去 ! Clear the settings for spectral method and Arakawa and Suarez (1983) !----------------------------------------------------------------- call Close(dyn_sp = dyn_sp_as % dyn_sp, err = err) ! (out) call Close(dyn_as = dyn_sp_as % dyn_as, err = err) ! (out) !----------------------------------------------------------------- ! ヒストリファイルへのデータ出力の終了処理 ! Terminate the settings for history data output !----------------------------------------------------------------- varnames = HstNmlInfoNames( dyn_sp_as % gthstnml ) call Split( str = varnames, sep = ',', carray = varnames_array ) ! (out) vnmax = size( varnames_array ) do i = 1, vnmax name = varnames_array(i) if ( trim( name ) == '' ) exit nullify( gthist ) call HstNmlInfoAssocGtHist( gthstnml = dyn_sp_as % gthstnml, name = name, history = gthist, err = err ) ! (out) if ( HistoryInitialized( gthist ) ) then call HistoryClose( history = gthist, err = err ) ! (out) end if end do !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- dyn_sp_as % initialized = .false. 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpAsClose
Subroutine : | |||||
dyn_sp_as : | type(DYNSPAS83), intent(inout) | ||||
nmax : | integer, intent(in)
| ||||
imax : | integer, intent(in)
| ||||
jmax : | integer, intent(in)
| ||||
kmax : | integer, intent(in)
| ||||
PI : | real(DP), intent(in)
| ||||
PI = 3.1415926535897930_DP : | real(DP), parameter
| ||||
RPlanet : | real(DP), intent(in)
| ||||
Omega : | real(DP), intent(in)
| ||||
Cp : | real(DP), intent(in)
| ||||
RAir : | real(DP), intent(in)
| ||||
EpsV : | real(DP), intent(in)
| ||||
VisOrder : | integer, intent(in)
| ||||
EFoldTime : | real(DP), intent(in)
| ||||
DelTime : | real(DP), intent(in)
| ||||
xy_Phis(0:imax-1, 0:jmax-1) : | real(DP), intent(in), optional
| ||||
openmp_threads : | integer, intent(in), optional
| ||||
r_SigmaSet(:) : | real(DP), intent(in), optional
| ||||
wa_module_initialized : | logical, intent(in), optional
| ||||
time_integration_scheme : | character(*), intent(in), optional
| ||||
current_time_value : | real, intent(in), optional
| ||||
current_time_unit : | character(*), intent(in), optional
| ||||
history_varlist : | character(*), intent(in), optional
| ||||
history_interval_value : | real, intent(in), optional
| ||||
history_interval_unit : | character(*), intent(in), optional
| ||||
history_precision : | character(*), intent(in), optional
| ||||
history_fileprefix : | character(*), intent(in), optional
| ||||
nmlfile : | character(*), intent(in), optional
| ||||
err : | logical, intent(out), optional
|
引数 dyn_sp_as に力学過程の設定を行う. 他のサブルーチンを使用する前に必ずこのサブルーチンによって DYNSPAS83 型の変数を初期設定してください.
NAMELIST を利用する場合には引数 nmlfile に NAMELIST ファイル名 を与えてください. NAMELIST 変数群の詳細に関しては NAMELIST#dyn_spectral_as83_nml を参照してください.
Configure the settings for dynamical core to dyn_sp_as. Before other subroutines are used, initialize "DYNSPAS83" variable.
In order to use NAMELIST, specify a NAMELIST filename to argument nmlfile. See "NAMELIST#dyn_spectral_as83_nml" for details about a NAMELIST group.
subroutine DynSpAsCreate( dyn_sp_as, nmax, imax, jmax, kmax, PI, RPlanet, Omega, Cp, RAir, EpsV, VisOrder, EFoldTime, DelTime, xy_Phis, openmp_threads, r_SigmaSet, wa_module_initialized, time_integration_scheme, current_time_value, current_time_unit, history_varlist, history_interval_value, history_interval_unit, history_precision, history_fileprefix, nmlfile, err) ! ! 引数 *dyn_sp_as* に力学過程の設定を行う. ! 他のサブルーチンを使用する前に必ずこのサブルーチンによって ! DYNSPAS83 型の変数を初期設定してください. ! ! NAMELIST を利用する場合には引数 *nmlfile* に NAMELIST ファイル名 ! を与えてください. NAMELIST 変数群の詳細に関しては ! NAMELIST#dyn_spectral_as83_nml を参照してください. ! ! Configure the settings for dynamical core to *dyn_sp_as*. ! Before other subroutines are used, initialize "DYNSPAS83" ! variable. ! ! In order to use NAMELIST, specify a NAMELIST filename to ! argument *nmlfile*. See "NAMELIST#dyn_spectral_as83_nml" ! for details about a NAMELIST group. ! use dc_trace, only: BeginSub, EndSub use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar use dc_types, only: STRING, DP, STDOUT use dc_present, only: present_and_not_empty, present_and_true use dc_message, only: MessageNotify use dc_date, only: DCDiffTimeCreate, EvalSec, EvalByUnit use dc_error, only: StoreError, DC_NOERR, DC_EALREADYINIT, DC_EARGLACK, DC_ENEGATIVE, DC_ENOFILEREAD, USR_ERRNO use dyn_spectral, only: Create, GetCoriolis, GetSpectralCoeff, GetDiffCoeff, GetPhis, GetAxes use dyn_as83, only: Create, GetAxes use dc_hash, only: HASH, DCHashPut, DCHashRewind, DCHashNext, DCHashDelete use gt4_history_nmlinfo, only: HstNmlInfoCreate, HstNmlInfoAdd, HstNmlInfoEndDefine, HstNmlInfoPutLine use gt4_history, only: GT_HISTORY, HistoryCreate, HistoryAddVariable, HistoryPut, HistoryAddAttr implicit none type(DYNSPAS83), intent(inout):: dyn_sp_as integer, intent(in):: nmax ! 最大全波数. ! Maximum truncated wavenumber integer, intent(in):: imax ! 経度格子点数. ! Number of grid points in longitude integer, intent(in):: jmax ! 緯度格子点数. ! Number of grid points in latitude integer, intent(in):: kmax ! 鉛直層数. ! Number of vertical level real(DP), intent(in):: PI ! $ \pi $ . 円周率. Circular constant real(DP), intent(in):: RPlanet ! $ a $ . 惑星半径. Radius of planet real(DP), intent(in):: Omega ! $ \Omega $ . 回転角速度. Angular velocity !!$ real(DP), intent(in):: Grav ! $ g $ . 重力加速度. Gravitational acceleration real(DP), intent(in):: Cp ! $ C_p $ . 大気定圧比熱. Specific heat of air at constant pressure real(DP), intent(in):: RAir ! $ R $ . 大気気体定数. Gas constant of air real(DP), intent(in):: EpsV ! $ \epsilon_v $ . 水蒸気分子量比. Molecular weight ratio of water vapor integer, intent(in):: VisOrder ! 超粘性の次数. Order of hyper-viscosity real(DP), intent(in):: EFoldTime ! 最大波数に対する e-folding time. E-folding time for maximum wavenumber real(DP), intent(in):: DelTime ! $ \Delta t $ . タイムステップ. Time step real(DP), intent(in), optional:: xy_Phis (0:imax-1, 0:jmax-1) ! $ \Phi_s $ . 地表ジオポテンシャル. ! Surface geo-potential integer, intent(in), optional:: openmp_threads ! OPENMP での最大スレッド数. ! openmp_threads に 1 より大きな値を指定すれば ! ISPACK[http://www.gfd-dennou.org/library/ispack/] ! の球面調和函数変換 OPENMP 並列計算 ! ルーチンが用いられる. 並列計算を実行するには, ! 実行時に環境変数 OMP_NUM_THREADS ! を openmp_threads 以下の数字に設定する ! 等のシステムに応じた準備が必要となる. ! ! openmp_threads に 1 より大きな値を ! 指定しなければ並列計算ルーチンは呼ばれない. ! ! Maximum number of threads in OPENMP. ! If integer more than 1 is given, ! spherical harmonics transration subroutines ! with OPENMP parallel computation in ! ISPACK[http://www.gfd-dennou.org/library/ispack/] ! is used. In practice, specify environment ! variable (OMP_NUM_THREADS etc.) for parallel ! computation. ! ! If number less than 2 is given, ! subroutines of parallel computation ! is not called. real(DP), intent(in), optional:: r_SigmaSet (:) ! デフォルトでは, kmax の値に応じ, ! 自動的に $ \sigma $ レベル (半整数) ! は設定される (kmax がある特定の値のみ). ! $ \sigma $ レベル (半整数). ! を明示的に設定する必要がある場合, ! この引数を与える. ! ! By default, half $ \sigma $ level is ! specified automatically according to ! the value of kmax (only the value is ! same as particular values). ! If the half $ \sigma $ level is specified ! manually, give this argument. logical, intent(in), optional :: wa_module_initialized ! wa_module (SPMODEL ライブラリ) 初期化フラグ. ! SPMODEL ライブラリの wa_module が ! ! 既にプログラム上で ! wa_Initial によって初期化されている場合, ! 再度 wa_Initial を呼ぶとエラーが生じます. ! この引数に .true. を与えることで, ! wa_Initial を呼ばないようにします. ! ! "wa_module" (SPMODEL library) ! initialization flag. ! ! When "wa_module" (SPMODEL library) ! is initialized by ! "wa_Initial" already, second "wa_Initial" ! causes an error. ! If .true. is specified to this argument, ! "wa_Initial" is not called. ! character(*), intent(in), optional:: time_integration_scheme ! 時間積分法. ! 以下の方法を選択可能. ! ! Time integration scheme. ! Available schemes are as follows. ! ! * "Semi-implicit" ! * "Explicit" ! real, intent(in), optional:: current_time_value ! 現在時刻の数値. Numerical value of current time character(*), intent(in), optional:: current_time_unit ! 現在時刻の単位. Unit of current time character(*), intent(in), optional:: history_varlist ! ヒストリデータの出力変数リスト. ! カンマで区切って並べる. ! (例: "Data1,Data2" ). ! ! List of variables output to history data. ! Delimiter is comma. ! (exp. "Data1,Data2" ). ! real, intent(in), optional:: history_interval_value ! ヒストリデータの出力間隔の数値. ! Numerical value for interval of history data output character(*), intent(in), optional:: history_interval_unit ! ヒストリデータの出力間隔の単位. ! Unit for interval of history data output character(*), intent(in), optional:: history_precision ! ヒストリデータの精度. ! Precision of history data character(*), intent(in), optional:: history_fileprefix ! ヒストリデータのファイル名の接頭詞. ! Prefix of history data filenames character(*), intent(in), optional :: nmlfile ! NAMELIST ファイルの名称. ! この引数に空文字以外を与えた場合, ! 指定されたファイルから ! NAMELIST 変数群を読み込みます. ! ファイルを読み込めない場合にはエラーを ! 生じます. ! ! NAMELIST 変数群の詳細に関しては ! NAMELIST#dyn_spectral_as83_nml ! を参照してください. ! ! NAMELIST file name. ! If nonnull character is specified to ! this argument, ! NAMELIST group name is loaded from the ! file. ! If the file can not be read, ! an error occurs. ! ! See "NAMELIST#dyn_spectral_as83_nml" ! for details about a NAMELIST group. ! 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. !----------------------------------- ! ヒストリファイルへのデータ出力設定 ! Configure the settings for history data output character(STRING):: name = '' ! 変数名. Variable identifier character(STRING):: longname = '' ! 変数の記述的名称. Descriptive name of variables character(STRING), allocatable:: dims(:) ! 座標軸の名称. Name of axes character(STRING):: units = '' ! 単位. Units type(GT_HISTORY), pointer:: gthist =>null() ! gt4_history モジュール用構造体. ! Derived type for "gt4_history" module character(TOKEN):: precision ! ヒストリデータの精度. ! Precision of history data logical:: average ! 出力データの平均化フラグ. ! Flag for average of output data type(HASH):: registered_varnames ! このモジュールから出力される変数名のリスト. ! ! List of names of variables output ! from this module. logical:: end !----------------------------------- ! 座標変数 ! Coordinate variables real(DP):: x_Lon (0:imax-1) ! 経度. Longitude real(DP):: y_Lat (0:jmax-1) ! 緯度. Latitude real(DP):: x_Lon_Weight(0:imax-1) ! 経度座標重み. ! Weight of longitude real(DP):: y_Lat_Weight(0:jmax-1) ! 緯度座標重み. ! Weight of latitude real(DP):: z_Sigma (0:kmax-1) ! $ \sigma $ レベル (整数). ! Full $ \sigma $ level real(DP):: r_Sigma (0:kmax) ! $ \sigma $ レベル (半整数). ! Half $ \sigma $ level !----------------------------------- ! 作業変数 ! Work variables real(DP):: xy_Cori (0:imax-1, 0:jmax-1) ! $ f\equiv 2\Omega\sin\varphi $ . ! コリオリパラメータ. Coriolis parameter integer:: nmo (1:2, 0:nmax, 0:nmax) ! スペクトルの添字順番. ! Spectral subscript expression real(DP):: wz_rn ((nmax+1)**2, 0:kmax-1) ! $ -n \times (n+1) $ . ラプラシアンの係数. ! Laplacian coefficient real(DP):: wz_DiffVorDiv ((nmax+1)**2, 0:kmax-1) ! $ -K_{HD} [(-1)^{N_D/2}\nabla^{N_D}- (2/a^2)^{N_D/2}] $ . ! 運動量水平拡散係数. ! Coefficient of horizontal momentum diffusion real(DP):: wz_DiffTherm ((nmax+1)**2, 0:kmax-1) ! $ -(-1)^{N_D/2}K_{HD}\nabla^{N_D} $ . ! 熱, 水水平拡散係数. ! Coefficient of horizontal thermal and water diffusion real(DP):: w_Phis ((nmax+1)**2) ! $ \Phi_s $ . 地表ジオポテンシャル. ! Surface geo-potential integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpAsCreate' continue call BeginSub(subname, version) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (dyn_sp_as % initialized) then stat = DC_EALREADYINIT cause_c = 'DYNSPAS83' goto 999 end if !----------------------------------------------------------------- ! 引数の正当性のチェック ! Validation of arguments !----------------------------------------------------------------- if (nmax < 1) then stat = DC_ENEGATIVE cause_c = 'nmax' goto 999 end if if (imax < 1) then stat = DC_ENEGATIVE cause_c = 'imax' goto 999 end if if (jmax < 1) then stat = DC_ENEGATIVE cause_c = 'jmax' goto 999 end if if (kmax < 1) then stat = DC_ENEGATIVE cause_c = 'kmax' goto 999 end if !----------------------------------------------------------------- ! 波数・格子点の設定 ! Configure wave number and grid point !----------------------------------------------------------------- dyn_sp_as % nmax = nmax dyn_sp_as % imax = imax dyn_sp_as % jmax = jmax dyn_sp_as % kmax = kmax !----------------------------------------------------------------- ! 時刻管理 ! Time control !----------------------------------------------------------------- if ( present(current_time_value) .and. present(current_time_unit) ) then call DCDiffTimeCreate( diff = dyn_sp_as % current_time, value = real(current_time_value, DP), unit = current_time_unit ) ! (in) else call DCDiffTimeCreate( diff = dyn_sp_as % current_time, value = 0.0_DP, unit = 'sec' ) ! (in) end if call DCDiffTimeCreate( diff = dyn_sp_as % delta_time, value = DelTime, unit = 'sec' ) ! (in) !----------------------------------------------------------------- ! スペクトル法の設定 ! Configure the settings for spectral method !----------------------------------------------------------------- call Create( dyn_sp = dyn_sp_as % dyn_sp, nmax = nmax, imax = imax, jmax = jmax, kmax = kmax, PI = PI, RPlanet = RPlanet, Omega = Omega, Cp = Cp, VisOrder = VisOrder, EFoldTime = EFoldTime, DelTime = DelTime, xy_Phis = xy_Phis, openmp_threads = openmp_threads, wa_module_initialized = wa_module_initialized, nmlfile = nmlfile ) ! (in) !----------------------------------------------------------------- ! Arakawa and Suarez (1983) の設定 ! Configure the settings for Arakawa and Suarez (1983) !----------------------------------------------------------------- !--------------------------------------------- ! コリオリパラメータの取得 ! Get Coriolis parameter call GetCoriolis( dyn_sp = dyn_sp_as % dyn_sp, xy_Cori = xy_Cori ) ! (out) !--------------------------------------------- ! スペクトル係数の取得 ! Get spectral coefficient call GetSpectralCoeff( dyn_sp = dyn_sp_as % dyn_sp, nmo = nmo, wz_rn = wz_rn ) ! (out) !--------------------------------------------- ! 拡散係数の取得 ! Get diffusion coefficient call GetDiffCoeff( dyn_sp = dyn_sp_as % dyn_sp, wz_DiffVorDiv = wz_DiffVorDiv, wz_DiffTherm = wz_DiffTherm ) ! (out) !--------------------------------------------- ! 地形データ (地表 $ \Phi $ )の取得 ! Get geography data (surface $ \Phi $ ) call GetPhis( dyn_sp = dyn_sp_as % dyn_sp, w_Phis = w_Phis ) ! (out) call Create( dyn_as = dyn_sp_as % dyn_as, nmax = nmax, imax = imax, jmax = jmax, kmax = kmax, PI = PI, RPlanet = RPlanet, Cp = Cp, RAir = RAir, EpsV = EpsV, DelTime = DelTime, xy_Cori = xy_Cori, nmo = nmo, wz_rn = wz_rn, wz_DiffVorDiv = wz_DiffVorDiv, wz_DiffTherm = wz_DiffTherm, r_SigmaSet = r_SigmaSet, w_Phis = w_Phis, time_integration_scheme = time_integration_scheme, nmlfile = nmlfile ) ! (in) !----------------------------------------------------------------- ! 緯度経度座標値の取得 ! Get data of latitude and longitude !----------------------------------------------------------------- call GetAxes( dyn_sp = dyn_sp_as % dyn_sp, x_Lon = x_Lon, y_Lat = y_Lat, x_Lon_Weight = x_Lon_Weight, y_Lat_Weight = y_Lat_Weight ) ! (out) !----------------------------------------------------------------- ! 鉛直レベル座標値の取得 ! Get data of vertical level !----------------------------------------------------------------- call GetAxes( dyn_as = dyn_sp_as % dyn_as, z_Sigma = z_Sigma, r_Sigma = r_Sigma ) ! (out) !----------------------------------------------------------------- ! ヒストリーデータに関する時刻管理 ! Time control for history data !----------------------------------------------------------------- if ( present(current_time_value) .and. present(current_time_unit) ) then call DCDiffTimeCreate( diff = dyn_sp_as % current_time, value = real( current_time_value, DP ), unit = current_time_unit ) ! (in) else call DCDiffTimeCreate( diff = dyn_sp_as % current_time, value = 0.0_DP, unit = 'sec' ) ! (in) end if call DCDiffTimeCreate( diff = dyn_sp_as % delta_time, value = DelTime, unit = 'sec' ) ! (in) !----------------------------------------------------------------- ! ヒストリファイルへのデータ出力設定 ! Configure the settings for history data output !----------------------------------------------------------------- !------------------------- ! デフォルト値 ! Default values call HstNmlInfoCreate( gthstnml = dyn_sp_as % gthstnml ) ! (inout) !------------------------- ! オプショナル引数からの値 ! Values from optional arguments call HstNmlInfoAdd( gthstnml = dyn_sp_as % gthstnml, name = '', interval_value = history_interval_value, interval_unit = history_interval_unit, precision = history_precision, average = .false., fileprefix = history_fileprefix ) ! (in) if ( present(history_varlist) ) then call HstNmlInfoAdd( gthstnml = dyn_sp_as % gthstnml, name = history_varlist ) ! (in) end if !----------------------------------------------------------------- ! NAMELIST からの値の読み込み ! Load values from NAMELIST !----------------------------------------------------------------- if ( present_and_not_empty(nmlfile) ) then call MessageNotify( 'M', subname, 'Loading NAMELIST file "%c" ...', c1 = trim(nmlfile) ) call NmlRead ( nmlfile = nmlfile, gthstnml = dyn_sp_as % gthstnml, err = err ) ! (out) if ( present_and_true(err) ) then call MessageNotify( 'W', subname, '"%c" can not be read.', c1 = trim(nmlfile) ) stat = DC_ENOFILEREAD cause_c = nmlfile goto 999 end if end if call HstNmlInfoEndDefine( gthstnml = dyn_sp_as % gthstnml, err = err ) ! (out) if ( present_and_true( err ) ) then stat = USR_ERRNO goto 999 end if !----------------------------------------------------------------- ! データ出力の初期設定 ! Initialize data output !----------------------------------------------------------------- !------------------------- ! xyr_SigmaDot の出力設定 ! Configure the settings for "xyr_SigmaDot" output name = 'SigmaDot' longname = 'sigma-vertical velocity' units = '1 s-1' allocate( dims(4) ) dims = StoA( 'lon', 'lat', 'sigm', 'time' ) ! 出力ファイルの初期設定. ! * gthist (gt4_history#GT_HISTORY) が設定される. ! Initialize output file. ! * "gthist" (gt4_history#GT_HISTORY) is configured. call output_init ! これは内部サブルーチン. This is internal subroutine ! 属性の付加などを行う場合には以下のようにする. ! Describe codes as follows in order to add attributes etc. !!$ if ( associated( gthist ) ) then !!$ call HistoryAddAttr( & !!$ & history = gthist, & ! (inout) !!$ & varname = name, attrname = 'missing_value', & ! (in) !!$ & value = 9999.0 ) ! (in) !!$ end if deallocate( dims ) !------------------------- ! xy_DPiDt の出力設定 ! Configure the settings for "xy_DPiDt" output name = 'DPiDt' longname = 'Pi (log Ps) tendency)' units = 'Pa s-1' allocate( dims(3) ) dims = StoA( 'lon', 'lat', 'time' ) ! 出力ファイルの初期設定. ! * gthist (gt4_history#GT_HISTORY) が設定される. ! Initialize output file. ! * "gthist" (gt4_history#GT_HISTORY) is configured. call output_init ! これは内部サブルーチン. This is internal subroutine deallocate( dims ) !----------------------------------------------------------------- ! このモジュールから出力される変数名のリストを表示 ! Print list of names of variables output from this module !----------------------------------------------------------------- call Printf( STDOUT, ' *** MESSAGE *** +---- "%c" module output varnames list -----', c1 = 'dcmodel_sample_code' ) call DCHashRewind( hashv = registered_varnames ) ! (inout) do call DCHashNext( hashv = registered_varnames, key = name, value = longname, end = end ) ! (out) if ( end ) exit call Printf( STDOUT, ' *** MESSAGE *** | "%c" (%c)', c1 = trim(name), c2 = trim(longname) ) enddo call DCHashDelete( hashv = registered_varnames ) ! (inout) call Printf( STDOUT, ' *** MESSAGE *** `----------------------------------------' ) !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- dyn_sp_as % initialized = .true. 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) contains subroutine output_init ! ! 変数 *name* に関して出力ファイルの初期設定を行います. ! 出力ファイル名や出力間隔などの情報は dyn_sp_as % gthstnml ! から取り出されます. ! ! 変数 *name* に関して出力が行われる場合には, ! *gthist* に出力先ファイルの gt4_history#GT_HISTORY ! 型変数を結合させます. そうでない場合は, *gthist* を空状態にします. ! ! また, 出力データの精度を precision に, ! 出力データ平均化の可否を average に設定します. ! ! 標準出力に表示される変数リスト *registered_varnames* に ! *name*, *longname*, *dims*, *units* が登録されます. ! ! An output file is initialized for a variable *name*. ! Information such as the output filename and output intervals ! is taken out of "dyn_sp_as % gthstnml". ! ! When output is done for the variable *name*, *gthist* is ! associated with the "gt4_history#GT_HISTORY" variable of ! the output file. Otherwise, *gthist* is nullified. ! ! Moreover, the accuracy of output data is set to *precision*, and ! right or wrong of averaging the output data is set to *average*. ! ! *name*, *longname*, *dims*, *units* are registered to ! a list of variables *registered_varnames* that is printed to ! standard output. ! use dc_date, only: DCDiffTimeCreate, EvalSec, EvalByUnit use gt4_history_nmlinfo, only: GTHST_NMLINFO, HstNmlInfoInitialized, HstNmlInfoInquire, HstNmlInfoOutputValid, HstNmlInfoAssocGtHist, HstNmlInfoPutLine use gt4_history, only: GT_HISTORY, HistoryCreate, HistoryAddVariable, HistoryPut, HistoryAddAttr, HistoryInitialized !----------------------------------- ! 作業変数 ! Work variables character(STRING):: file ! ヒストリデータのファイル名. ! History data filenames character(STRING):: dims_str ! 座標軸のリスト. ! List of axes real:: interval_value ! ヒストリデータの出力間隔の数値. ! Numerical value for interval of history data output character(TOKEN):: interval_unit ! ヒストリデータの出力間隔の単位. ! Unit for interval of history data output real(DP), parameter:: PI = 3.1415926535897930_DP ! $ \pi $ . 円周率. Circular constant continue !----------------------------------------------------------------- ! 標準出力に表示される変数の登録 ! Register a variable name for print to standard output !----------------------------------------------------------------- if ( allocated(dims) ) then dims_str = JoinChar( dims, ',' ) else dims_str = '' end if call DCHashPut( hashv = registered_varnames, key = name, value = trim( longname ) // ' [' // trim( units ) // '] {' // trim( dims_str ) // '}' ) ! (in) !----------------------------------------------------------------- ! 変数の初期化 ! Initialize variable !----------------------------------------------------------------- nullify( gthist ) precision = 'float' average = .false. !----------------------------------------------------------------- ! 出力が有効かどうかを確認する ! Confirm whether the output is effective !----------------------------------------------------------------- if ( .not. HstNmlInfoOutputValid( dyn_sp_as % gthstnml, name ) ) then return end if !----------------------------------------------------------------- ! GT_HISTORY 変数の取得 ! Get "GT_HISTORY" variable !----------------------------------------------------------------- call HstNmlInfoAssocGtHist( gthstnml = dyn_sp_as % gthstnml, name = name, history = gthist, err = err ) ! (out) if ( present_and_true( err ) ) return call HstNmlInfoInquire( gthstnml = dyn_sp_as % gthstnml, name = name, precision = precision, average = average, err = err ) ! (out) if ( present_and_true( err ) ) return !----------------------------------------------------------------- ! GT_HISTORY 変数の初期設定の確認 ! Check initialization of "GT_HISTORY" variable !----------------------------------------------------------------- if ( HistoryInitialized( gthist ) ) then !--------------------------------------------------------------- ! HistoryAddVariable による変数作成 ! A variable is created by "HistoryAddVariable" !--------------------------------------------------------------- call HistoryAddVariable( history = gthist, varname = name, dims = dims, longname = longname, units = units, xtype = precision, average = average ) ! (in) return end if !----------------------------------------------------------------- ! HistoryCreate のための設定値の取得 ! Get the settings for "HistoryCreate" !----------------------------------------------------------------- call HstNmlInfoInquire( gthstnml = dyn_sp_as % gthstnml, name = name, file = file, interval_unit = interval_unit, interval_value = interval_value, err = err ) ! (out) if ( present_and_true( err ) ) return !----------------------------------------------------------------- ! HistoryCreate によるファイル作成 ! Files are created by "HistoryCreate" !----------------------------------------------------------------- call HistoryCreate( history = gthist, file = file, title = 'Dynamical core calculation', source = 'dcpam4 : ' // trim(version), institution = 'GFD Dennou Club', dims = StoA('lon', 'lat', 'sig', 'sigm', 'time'), dimsizes = (/ dyn_sp_as % imax, dyn_sp_as % jmax, dyn_sp_as % kmax, dyn_sp_as % kmax+1, 0 /), longnames = StoA('longitude', 'latitude', 'sigma at layer midpoints', 'sigma at layer end-points (half level)', 'time'), units = StoA( 'degree_east', 'degree_north', '1', '1', interval_unit ), origin = real( EvalbyUnit( dyn_sp_as % current_time, interval_unit) ), interval = interval_value, err = err ) ! (out) if ( present_and_true( err ) ) then nullify( gthist ) return end if call HistoryAddAttr( history = gthist, varname = 'lon', attrname = 'standard_name', value = 'longitude' ) ! (in) call HistoryAddAttr( history = gthist, varname = 'lat', attrname = 'standard_name', value = 'latitude' ) ! (in) call HistoryAddAttr( history = gthist, varname = 'sig', attrname = 'standard_name', value = 'atmosphere_sigma_coordinate' ) ! (in) call HistoryAddAttr( history = gthist, varname = 'sigm', attrname = 'standard_name', value = 'atmosphere_sigma_coordinate' ) ! (in) call HistoryAddAttr( history = gthist, varname = 'time', attrname = 'standard_name', value = 'time' ) ! (in) call HistoryAddAttr( history = gthist, varname = 'sig', attrname = 'positive', value = 'down' ) ! (in) call HistoryAddAttr( history = gthist, varname = 'sigm', attrname = 'positive', value = 'down' ) ! (in) call HistoryPut( history = gthist, varname = 'lon', array = x_Lon / PI * 180.0_DP ) ! (in) call HistoryPut( history = gthist, varname = 'lat', array = y_Lat / PI * 180.0_DP ) ! (in) call HistoryPut( history = gthist, varname = 'sig', array = z_Sigma ) ! (in) call HistoryPut( history = gthist, varname = 'sigm', array = r_Sigma ) ! (in) call HistoryAddVariable( history = gthist, varname = 'lon_weight', dims = StoA('lon'), longname = 'weight for integration in longitude', units = 'radian', xtype = 'double' ) ! (in) call HistoryAddAttr( history = gthist, varname = 'lon', attrname = 'gt_calc_weight', value = 'lon_weight' ) ! (in) call HistoryPut( history = gthist, varname = 'lon_weight', array = x_Lon_Weight ) ! (in) call HistoryAddVariable( history = gthist, varname = 'lat_weight', dims = StoA('lat'), longname = 'weight for integration in latitude', units = 'radian', xtype = 'double' ) ! (in) call HistoryAddAttr( history = gthist, varname = 'lat', attrname = 'gt_calc_weight', value = 'lat_weight' ) ! (in) call HistoryPut( history = gthist, varname = 'lat_weight', array = y_Lat_Weight ) ! (in) !----------------------------------------------------------------- ! HistoryAddVariable による変数作成 ! A variable is created by "HistoryAddVariable" !----------------------------------------------------------------- if ( HistoryInitialized( gthist ) ) then call HistoryAddVariable( history = gthist, varname = name, dims = dims, longname = longname, units = units, xtype = precision, average = average ) ! (in) else nullify( gthist ) end if end subroutine output_init end subroutine DynSpAsCreate
Subroutine : | |||
dyn_sp_as : | type(DYNSPAS83), intent(inout) | ||
x_Lon(0:dyn_sp_as%imax-1) : | real(DP), intent(in), optional
| ||
y_Lat(0:dyn_sp_as%jmax-1) : | real(DP), intent(in), optional
| ||
z_Sigma(0:dyn_sp_as%kmax-1) : | real(DP), intent(in), optional
| ||
r_Sigma(0:dyn_sp_as%kmax) : | real(DP), intent(in), optional
| ||
err : | logical, intent(out), optional
|
与えられた座標軸データ と dyn_sp_as 内に 保持される座標データが等しいことを確認します. 異なる場合にはエラーを発生させます.
なお, 与えられた dyn_sp_as が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Confirm equality between given axes data and axes data stored in dyn_sp_as. If the equality is not confirmed, error is occurred.
If dyn_sp_as is not initialized by "Create" yet, error is occurred.
subroutine DynSpAsEqualAxes( dyn_sp_as, x_Lon, y_Lat, z_Sigma, r_Sigma, err ) ! ! 与えられた座標軸データ と *dyn_sp_as* 内に ! 保持される座標データが等しいことを確認します. ! 異なる場合にはエラーを発生させます. ! ! なお, 与えられた *dyn_sp_as* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Confirm equality between given axes data and ! axes data stored in *dyn_sp_as*. ! If the equality is not confirmed, error is occurred. ! ! If *dyn_sp_as* is not initialized by "Create" yet, ! error is occurred. ! use dc_trace, only: BeginSub, EndSub use dc_types, only: STRING, STDOUT, DP use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT, DCPAM_EAXISMISMATCH use dyn_spectral, only: EqualAxes use dyn_as83, only: EqualAxes implicit none type(DYNSPAS83), intent(inout):: dyn_sp_as real(DP), intent(in), optional:: x_Lon (0:dyn_sp_as%imax-1) ! 経度. Longitude real(DP), intent(in), optional:: y_Lat (0:dyn_sp_as%jmax-1) ! 緯度. Latitude real(DP), intent(in), optional:: z_Sigma (0:dyn_sp_as%kmax-1) ! $ \sigma $ レベル (整数). ! Full $ \sigma $ level real(DP), intent(in), optional:: r_Sigma (0:dyn_sp_as%kmax) ! $ \sigma $ レベル (半整数). ! Half $ \sigma $ level 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. integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpAsEqualAxes' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp_as % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSPAS83' goto 999 end if !----------------------------------------------------------------- ! 軸データの比較 ! Compare axes data !----------------------------------------------------------------- call EqualAxes( dyn_sp = dyn_sp_as % dyn_sp, x_Lon = x_Lon, y_Lat = y_Lat, err = err ) ! (out) call EqualAxes( dyn_as = dyn_sp_as % dyn_as, z_Sigma = z_Sigma, r_Sigma = r_Sigma, err = err ) ! (out) !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpAsEqualAxes
Subroutine : | |||
dyn_sp_as : | type(DYNSPAS83), intent(inout) | ||
x_Lon(0:dyn_sp_as%imax-1) : | real(DP), intent(out)
| ||
y_Lat(0:dyn_sp_as%jmax-1) : | real(DP), intent(out)
| ||
z_Sigma(0:dyn_sp_as%kmax-1) : | real(DP), intent(out)
| ||
r_Sigma(0:dyn_sp_as%kmax) : | real(DP), intent(out)
| ||
z_DelSigma(0:dyn_sp_as%kmax-1) : | real(DP), intent(out), optional
| ||
err : | logical, intent(out), optional
|
dyn_sp_as に格納されている座標値を返します.
なお, 与えられた dyn_sp_as が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Return data of axes stored in dyn_sp_as.
If dyn_sp_as is not initialized by "Create" yet, error is occurred.
subroutine DynSpAsGetAxes( dyn_sp_as, x_Lon, y_Lat, z_Sigma, r_Sigma, z_DelSigma, err ) ! ! *dyn_sp_as* に格納されている座標値を返します. ! ! なお, 与えられた *dyn_sp_as* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Return data of axes stored in *dyn_sp_as*. ! ! If *dyn_sp_as* is not initialized by "Create" yet, ! error is occurred. ! use dc_trace, only: BeginSub, EndSub use dc_types, only: STRING, STDOUT, DP use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT use dyn_spectral, only: GetAxes use dyn_as83, only: GetAxes implicit none type(DYNSPAS83), intent(inout):: dyn_sp_as real(DP), intent(out):: x_Lon (0:dyn_sp_as%imax-1) ! 経度. Longitude real(DP), intent(out):: y_Lat (0:dyn_sp_as%jmax-1) ! 緯度. Latitude real(DP), intent(out):: z_Sigma (0:dyn_sp_as%kmax-1) ! $ \sigma $ レベル (整数). ! Full $ \sigma $ level real(DP), intent(out):: r_Sigma (0:dyn_sp_as%kmax) ! $ \sigma $ レベル (半整数). ! Half $ \sigma $ level real(DP), intent(out), optional:: z_DelSigma (0:dyn_sp_as%kmax-1) ! $ \Delta \sigma $ (整数). ! $ \Delta \sigma $ (Full) 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. integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpAsSample' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp_as % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSPAS83' goto 999 end if !----------------------------------------------------------------- ! 緯度経度座標値の取得 ! Get data of latitude and longitude !----------------------------------------------------------------- call GetAxes( dyn_sp = dyn_sp_as % dyn_sp, x_Lon = x_Lon, y_Lat = y_Lat ) ! (out) !----------------------------------------------------------------- ! 鉛直レベル座標値の取得 ! Get data of vertical level !----------------------------------------------------------------- call GetAxes( dyn_as = dyn_sp_as % dyn_as, z_Sigma = z_Sigma, r_Sigma = r_Sigma, z_DelSigma = z_DelSigma ) ! (out) !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpAsGetAxes
Function : | |
result : | logical |
dyn_sp_as : | type(DYNSPAS83), intent(in) |
dyn_sp_as が初期設定されている場合には .true. が, 初期設定されていない場合には .false. が返ります.
If dyn_sp_as is initialized, .true. is returned. If dyn_sp_as is not initialized, .false. is returned.
logical function DynSpAsInitialized( dyn_sp_as ) result(result) ! ! *dyn_sp_as* が初期設定されている場合には .true. が, ! 初期設定されていない場合には .false. が返ります. ! ! If *dyn_sp_as* is initialized, .true. is returned. ! If *dyn_sp_as* is not initialized, .false. is returned. ! implicit none type(DYNSPAS83), intent(in):: dyn_sp_as continue result = dyn_sp_as % initialized end function DynSpAsInitialized
Subroutine : | |||
dyn_sp_as : | type(DYNSPAS83), intent(inout) | ||
unit : | integer, intent(in), optional
| ||
indent : | character(*), intent(in), optional
| ||
err : | logical, intent(out), optional
|
引数 dyn_sp_as に設定されている情報を印字します. デフォルトではメッセージは標準出力に出力されます. unit に装置番号を指定することで, 出力先を変更することが可能です.
Print information of dyn_sp_as. By default messages are output to standard output. Unit number for output can be changed by unit argument.
subroutine DynSpAsPutLine( dyn_sp_as, unit, indent, err ) ! ! 引数 *dyn_sp_as* に設定されている情報を印字します. ! デフォルトではメッセージは標準出力に出力されます. ! *unit* に装置番号を指定することで, 出力先を変更することが可能です. ! ! Print information of *dyn_sp_as*. ! 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_types, only: STRING, STDOUT use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar use dyn_spectral, only: PutLine use dyn_as83, only: PutLine use gt4_history_nmlinfo, only: GTHST_NMLINFO, HstNmlInfoPutLine implicit none type(DYNSPAS83), intent(inout):: dyn_sp_as 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 integer:: stat character(STRING):: cause_c integer:: out_unit integer:: indent_len character(STRING):: indent_str character(*), parameter:: subname = 'DynSpAsPutLine' 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 !----------------------------------------------------------------- ! "DYNSPAS83" の設定の印字 ! Print the settings for "DYNSPAS83" !----------------------------------------------------------------- if (dyn_sp_as % initialized) then call Printf(out_unit, indent_str(1:indent_len) // '#<DYNSPAS83:: @initialized=%y', l=(/dyn_sp_as % initialized/)) call PutLine(dyn_sp_as % dyn_sp, out_unit, indent_str(1:indent_len) // ' ', err) call PutLine(dyn_sp_as % dyn_as, out_unit, indent_str(1:indent_len) // ' ', err) call HstNmlInfoPutLine( gthstnml = dyn_sp_as % gthstnml, unit = out_unit, indent = indent_str(1:indent_len) // ' ' ) ! (in) call Printf( out_unit, indent_str(1:indent_len) // '>' ) else call Printf(out_unit, indent_str(1:indent_len) // '#<DYNSPAS83:: @initialized=%y>', l=(/dyn_sp_as % initialized/)) end if !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpAsPutLine
Subroutine : | |||
dyn_sp_as : | type(DYNSPAS83), intent(inout) | ||
xyz_VorB(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xyz_DivB(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xyz_TempB(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xyz_QVapB(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xy_PsB(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1) : | real(DP), intent(in)
| ||
xyz_VorN(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xyz_DivN(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xyz_TempN(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xyz_QVapN(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xy_PsN(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1) : | real(DP), intent(in)
| ||
xyz_VorA(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(out)
| ||
xyz_DivA(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(out)
| ||
xyz_TempA(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(out)
| ||
xyz_QVapA(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(out)
| ||
xy_PsA(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1) : | real(DP), intent(out)
| ||
xyz_DVorDt(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in), optional
| ||
xyz_DDivDt(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in), optional
| ||
xyz_DTempDt(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in), optional
| ||
xyz_DQVapDt(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in), optional
| ||
historyput_flag : | logical, intent(in), optional
| ||
err : | logical, intent(out), optional
|
力学過程の演算を行い, 与えられた $ t-\Delta t $ および $ t $ の 渦度 (xyz_VorB, xyz_VorN), 発散 (xyz_DivB, xyz_DivN), 温度 (xyz_TempB, xyz_TempN), 比湿 (xyz_QVapB, xyz_QVapN), 地表面気圧 (xyz_PsB, xyz_PsN) から, $ t+Delta t $ の 渦度 (xyz_VorA), 発散 (xyz_DivA), 温度 (xyz_TempA), 比湿 (xyz_QVapA), 地表面気圧 (xyz_PsA) を返します.
別の物理プロセスによる渦度, 発散, 温度, 比湿の変化を, 力学過程の変化に足して次のステップを計算する場合には, それらの変化を xyz_DVorDt, xyz_DDivDt, xyz_DTempDt, xyz_DQVapDt に与えてください.
時間積分法にはリープフロッグスキームを用いています. デフォルトでは, $ Delta t $ を大きくとるために, 重力波項に セミインプリシット法を適用しています. ただし, Create の time_integration_scheme に "Explicit" を指定した場合には重力波項もエクスプリシット法によって 解きます.
なお, 与えられた dyn_sp_as が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculating dynamical core, from vorticity (xyz_VorB, xyz_VorN), divergence (xyz_DivB, xyz_DivN), temperature (xyz_TempB, xyz_TempN), specific humidity (xyz_QVapB, xyz_QVapN), surface pressure (xyz_PsB, xyz_PsN) at $ t-\Delta t $ and $ t $, vorticity (xyz_VorA), divergence (xyz_DivA), temperature (xyz_TempA), specific humidity (xyz_QVapA), surface pressure (xyz_PsA) are returned.
In order to add tendencies of vorticity, divergence, temperature and specific humidity by another physical process to tendencies of this dynamical process and to calculate the values at next time step, give these tendencies to "xyz_DVorDt", "xyz_DDivDt", "xyz_DTempDt", "xyz_DQVapDt"
Leap-frog scheme is used as time integration. By default, semi-implicit scheme is applied to gravitational terms for extension of $ Delta t $ . If "Explicit" is specified to time_integration_scheme in "Create", explicit scheme is applied to gravitational terms as well as other terms.
If dyn_sp_as is not initialized by "Create" yet, error is occurred.
Alias for DynSpAsDynamics
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
xyz_U(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_V(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_Vor(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
xyz_Div(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
与えられた東西風速と南北風速から渦度と発散を計算します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculate vorticity and divergence from given zonal and meridional wind.
If dyn_sp is not initialized by "Create" yet, error is occurred.
Original external subprogram is dyn_spectral#UV2VorDiv
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
xyz_U(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_V(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_Vor(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
xyz_Div(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
与えられた東西風速と南北風速から渦度と発散を計算します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculate vorticity and divergence from given zonal and meridional wind.
If dyn_sp is not initialized by "Create" yet, error is occurred.
Original external subprogram is dyn_spectral#UV2VorDiv
Subroutine : | |||
dyn_sp_as : | type(DYNSPAS83), intent(inout) | ||
xyz_Vor(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(out)
| ||
xyz_Div(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(out)
| ||
xyz_U(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xyz_V(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
err : | logical, intent(out), optional
|
与えられた東西風速と南北風速から渦度と発散を計算します.
なお, 与えられた dyn_sp_as が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculate vorticity and divergence from given zonal and meridional wind.
If dyn_sp_as is not initialized by "Create" yet, error is occurred.
Alias for DynSpAsUV2VorDiv
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
xyz_Vor(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_Div(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_U(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
xyz_V(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
与えられた渦度と発散から東西風速と南北風速を計算します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculate zonal and meridional wind from given vorticity and divergence.
If dyn_sp is not initialized by "Create" yet, error is occurred.
Original external subprogram is dyn_spectral#VorDiv2UV
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
xyz_Vor(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_Div(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_U(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
xyz_V(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
与えられた渦度と発散から東西風速と南北風速を計算します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculate zonal and meridional wind from given vorticity and divergence.
If dyn_sp is not initialized by "Create" yet, error is occurred.
Original external subprogram is dyn_spectral#VorDiv2UV
Subroutine : | |||
dyn_sp_as : | type(DYNSPAS83), intent(inout) | ||
xyz_Vor(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xyz_Div(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xyz_U(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(out)
| ||
xyz_V(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
与えられた渦度と発散から東西風速と南北風速を計算します.
なお, 与えられた dyn_sp_as が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculate zonal and meridional wind from given vorticity and divergence.
If dyn_sp_as is not initialized by "Create" yet, error is occurred.
Alias for DynSpAsVorDiv2UV
Subroutine : | |||
dyn_sp_as : | type(DYNSPAS83), intent(inout) | ||
xyz_VorB(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xyz_DivB(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xyz_TempB(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xyz_QVapB(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xy_PsB(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1) : | real(DP), intent(in)
| ||
xyz_VorN(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xyz_DivN(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xyz_TempN(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xyz_QVapN(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xy_PsN(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1) : | real(DP), intent(in)
| ||
xyz_VorA(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(out)
| ||
xyz_DivA(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(out)
| ||
xyz_TempA(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(out)
| ||
xyz_QVapA(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(out)
| ||
xy_PsA(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1) : | real(DP), intent(out)
| ||
xyz_DVorDt(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in), optional
| ||
xyz_DDivDt(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in), optional
| ||
xyz_DTempDt(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in), optional
| ||
xyz_DQVapDt(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in), optional
| ||
historyput_flag : | logical, intent(in), optional
| ||
err : | logical, intent(out), optional
|
力学過程の演算を行い, 与えられた $ t-\Delta t $ および $ t $ の 渦度 (xyz_VorB, xyz_VorN), 発散 (xyz_DivB, xyz_DivN), 温度 (xyz_TempB, xyz_TempN), 比湿 (xyz_QVapB, xyz_QVapN), 地表面気圧 (xyz_PsB, xyz_PsN) から, $ t+Delta t $ の 渦度 (xyz_VorA), 発散 (xyz_DivA), 温度 (xyz_TempA), 比湿 (xyz_QVapA), 地表面気圧 (xyz_PsA) を返します.
別の物理プロセスによる渦度, 発散, 温度, 比湿の変化を, 力学過程の変化に足して次のステップを計算する場合には, それらの変化を xyz_DVorDt, xyz_DDivDt, xyz_DTempDt, xyz_DQVapDt に与えてください.
時間積分法にはリープフロッグスキームを用いています. デフォルトでは, $ Delta t $ を大きくとるために, 重力波項に セミインプリシット法を適用しています. ただし, Create の time_integration_scheme に "Explicit" を指定した場合には重力波項もエクスプリシット法によって 解きます.
なお, 与えられた dyn_sp_as が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculating dynamical core, from vorticity (xyz_VorB, xyz_VorN), divergence (xyz_DivB, xyz_DivN), temperature (xyz_TempB, xyz_TempN), specific humidity (xyz_QVapB, xyz_QVapN), surface pressure (xyz_PsB, xyz_PsN) at $ t-\Delta t $ and $ t $, vorticity (xyz_VorA), divergence (xyz_DivA), temperature (xyz_TempA), specific humidity (xyz_QVapA), surface pressure (xyz_PsA) are returned.
In order to add tendencies of vorticity, divergence, temperature and specific humidity by another physical process to tendencies of this dynamical process and to calculate the values at next time step, give these tendencies to "xyz_DVorDt", "xyz_DDivDt", "xyz_DTempDt", "xyz_DQVapDt"
Leap-frog scheme is used as time integration. By default, semi-implicit scheme is applied to gravitational terms for extension of $ Delta t $ . If "Explicit" is specified to time_integration_scheme in "Create", explicit scheme is applied to gravitational terms as well as other terms.
If dyn_sp_as is not initialized by "Create" yet, error is occurred.
subroutine DynSpAsDynamics( dyn_sp_as, xyz_VorB, xyz_DivB, xyz_TempB, xyz_QVapB, xy_PsB, xyz_VorN, xyz_DivN, xyz_TempN, xyz_QVapN, xy_PsN, xyz_VorA, xyz_DivA, xyz_TempA, xyz_QVapA, xy_PsA, xyz_DVorDt, xyz_DDivDt, xyz_DTempDt, xyz_DQVapDt, historyput_flag, err ) ! ! 力学過程の演算を行い, 与えられた $ t-\Delta t $ および $ t $ の ! 渦度 (xyz_VorB, xyz_VorN), 発散 (xyz_DivB, xyz_DivN), ! 温度 (xyz_TempB, xyz_TempN), 比湿 (xyz_QVapB, xyz_QVapN), ! 地表面気圧 (xyz_PsB, xyz_PsN) から, ! $ t+\Delta t $ の ! 渦度 (xyz_VorA), 発散 (xyz_DivA), 温度 (xyz_TempA), ! 比湿 (xyz_QVapA), 地表面気圧 (xyz_PsA) を返します. ! ! 別の物理プロセスによる渦度, 発散, 温度, 比湿の変化を, ! 力学過程の変化に足して次のステップを計算する場合には, ! それらの変化を xyz_DVorDt, xyz_DDivDt, xyz_DTempDt, xyz_DQVapDt ! に与えてください. ! ! 時間積分法にはリープフロッグスキームを用いています. ! デフォルトでは, $ \Delta t $ を大きくとるために, 重力波項に ! セミインプリシット法を適用しています. ! ただし, Create の *time_integration_scheme* に ! "Explicit" を指定した場合には重力波項もエクスプリシット法によって ! 解きます. ! ! なお, 与えられた *dyn_sp_as* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Calculating dynamical core, ! from vorticity (xyz_VorB, xyz_VorN), divergence (xyz_DivB, xyz_DivN), ! temperature (xyz_TempB, xyz_TempN), ! specific humidity (xyz_QVapB, xyz_QVapN), ! surface pressure (xyz_PsB, xyz_PsN) at $ t-\Delta t $ and $ t $, ! vorticity (xyz_VorA), divergence (xyz_DivA), temperature (xyz_TempA), ! specific humidity (xyz_QVapA), surface pressure (xyz_PsA) ! are returned. ! ! In order to add tendencies of vorticity, divergence, ! temperature and specific humidity by another physical process to ! tendencies of this dynamical process ! and to calculate the values at next time step, ! give these tendencies to ! "xyz_DVorDt", "xyz_DDivDt", "xyz_DTempDt", "xyz_DQVapDt" ! ! Leap-frog scheme is used as time integration. ! By default, semi-implicit scheme is applied to gravitational terms ! for extension of $ \Delta t $ . ! If "Explicit" is specified to *time_integration_scheme* in "Create", ! explicit scheme is applied to gravitational terms as well as ! other terms. ! ! If *dyn_sp_as* is not initialized by "Create" yet, ! error is occurred. ! use dc_trace, only: BeginSub, EndSub use dc_types, only: STRING, STDOUT, DP, TOKEN use dc_string, only: LChar, PutLine, Printf, Split, StrInclude, StoA, JoinChar use dc_present, only: present_and_true use dc_date, only: mod, operator(+), operator(==), EvalbyUnit use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT use dyn_spectral, only: GradPi, VorDiv2UV, Tendency, TendencyExplicit, Spectral2Grid, Grid2Spectral, DiffusionCorrectTemp, DiffusionVorDiv use dyn_as83, only: NonLinearOnGrid, TimeIntegration, Inquire, WTplusGPiOnGrid, HDivOnGrid use gt4_history_nmlinfo, only: HstNmlInfoInquire, HstNmlInfoOutputValid, HstNmlInfoAssocGtHist, HstNmlInfoPutLine use gt4_history, only: HistoryPut, HistoryInitialized implicit none type(DYNSPAS83), intent(inout):: dyn_sp_as real(DP), intent(in):: xyz_VorB (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ \zeta (t-\Delta t) $ . 渦度. Vorticity real(DP), intent(in):: xyz_DivB (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ D (t-\Delta t) $ . 発散. Divergence real(DP), intent(in):: xyz_TempB (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ T (t-\Delta t) $ . 温度. Temperature real(DP), intent(in):: xy_PsB (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1) ! $ p_s (t-\Delta t) $ . 地表面気圧. Surface pressure real(DP), intent(in):: xyz_QVapB (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ q (t-\Delta t) $ . 比湿. Specific humidity real(DP), intent(in):: xyz_VorN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ \zeta (t) $ . 渦度. Vorticity real(DP), intent(in):: xyz_DivN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ D (t) $ . 発散. Divergence real(DP), intent(in):: xyz_TempN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ T (t) $ . 温度. Temperature real(DP), intent(in):: xy_PsN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1) ! $ p_s (t) $ . 地表面気圧. Surface pressure real(DP), intent(in):: xyz_QVapN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ q (t) $ . 比湿. Specific humidity real(DP), intent(out):: xyz_VorA (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ \zeta (t+\Delta t) $ . 渦度. Vorticity real(DP), intent(out):: xyz_DivA (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ D (t+\Delta t) $ . 発散. Divergence real(DP), intent(out):: xyz_TempA (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ T (t+\Delta t) $ . 温度. Temperature real(DP), intent(out):: xy_PsA (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1) ! $ p_s (t+\Delta t) $ . 地表面気圧. Surface pressure real(DP), intent(out):: xyz_QVapA (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ q (t+\Delta t) $ . 比湿. Specific humidity real(DP), intent(in), optional:: xyz_DVorDt (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ \DP{\zeta}{t} $ . ! 渦度変化. ! Vorticity tendency real(DP), intent(in), optional:: xyz_DDivDt (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ \DP{D}{t} $ . ! 発散変化. ! Divergence tendency real(DP), intent(in), optional:: xyz_DTempDt (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ \DP{T}{t} $ . ! 温度変化. ! Temperature tendency real(DP), intent(in), optional:: xyz_DQVapDt (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ \DP{q}{t} $ . ! 比湿変化. ! Specific humidity tendency logical, intent(in), optional:: historyput_flag ! データ出力のフラグ. ! SetTime によって時刻を明示的に ! 指定した場合には, この引数に ! .true. または .false. を指定する ! ことでデータ出力のオンオフを ! 明示的に指定する必要があります. ! デフォルトは .false. です. ! ! Data output flag. ! When time is specified by "SetTime", ! explicit specification of data output ! on/off by specifying ".true." or ".false." ! to this argument. ! Default value is ".false.". ! 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 real(DP):: xy_GradLonPiN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1) ! $ \frac{1}{\cos \varphi} \DP{\pi}{\lambda} $ real(DP):: xy_GradLatPiN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1) ! $ \DP{\pi}{\varphi} $ real(DP):: xyz_UN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ u (t) $ . 東西風速. Zonal wind. real(DP):: xyz_VN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ v (t) $ . 南北風速. Meridional wind. real(DP):: xyz_UA (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ u (t+\Delta t) $ . 東西風速. Zonal wind real(DP):: xyz_VA (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ v (t+\Delta t) $ . 南北風速. Meridional wind real(DP):: xyz_UAdvN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ u_A (t) $ . 東西運動量移流項. ! Zonal advection of momentum real(DP):: xyz_VAdvN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ v_A (t) $ . 南北運動量移流項. ! Meridional advection of momentum real(DP):: xyz_DTempDtN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ H (t) $ . 温度時間変化項. ! Temperature tendency real(DP):: xyz_DQVapDtN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ R (t) $ . 比湿時間変化項. ! Specific humidity tendency real(DP):: xyz_KEN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ KE (t) $ . 運動エネルギー項. ! Kinematic energy real(DP):: xyz_TempUAdvN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ uT (t) $ . 温度東西移流項. ! Zonal advection of temperature real(DP):: xyz_TempVAdvN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ vT (t) $ . 温度南北移流項. ! Meridional advection of temperature real(DP):: xyr_SigmaDotN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax) ! $ \dot{\sigma} $ . ! 鉛直流. Vertical flow real(DP):: xy_DPiDtN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1) ! $ Z $ . 地表面気圧時間変化項. ! Surface pressure tendency real(DP):: xyz_QVapUAdvN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ uq (t) $ . 比湿東西移流項. ! Zonal advection of specific humidity real(DP):: xyz_QVapVAdvN (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ vq (t) $ . 比湿南北移流項. ! Meridional advection of specific humidity real(DP):: wz_DVorDtN ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1) ! $ \DD{\zeta}{t} (t) $ . 渦度変化 (スペクトル). ! Vorticity tendency (spectral) real(DP):: wz_DDivDtN ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1) ! $ \DD{D}{t} (t) $ . 発散変化 (スペクトル). ! Divergence tendency (spectral) real(DP):: wz_DTempDtN ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1) ! $ \DD{T}{t} (t) $ . 温度変化 (スペクトル). ! Temperature tendency (spectral) real(DP):: wz_DQVapDtN ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1) ! $ \DD{q}{t} (t) $ . 比湿変化 (スペクトル). ! Specific humidity tendency (spectral) real(DP):: w_DPiDtN ((dyn_sp_as%nmax+1)**2) ! $ \DD{p_s}{t} (t) $ . 地表面気圧変化 (スペクトル). ! Surface pressure tendency (spectral) real(DP):: wz_VorB ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1) ! $ \zeta (t-\Delta t) $ . 渦度 (スペクトル). ! Vorticity (spectral) real(DP):: wz_DivB ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1) ! $ D (t-\Delta t) $ . 発散 (スペクトル). ! Divergence (spectral) real(DP):: wz_TempB ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1) ! $ T (t-\Delta t) $ . 温度 (スペクトル). ! Temperature (spectral) real(DP):: w_PiB ((dyn_sp_as%nmax+1)**2) ! $ \pi = \ln p_s (t-\Delta t) $ . 地表面気圧 (スペクトル). ! Surface pressure (spectral) real(DP):: wz_QVapB ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1) ! $ q (t-\Delta t) $ . 比湿 (スペクトル). ! Specific humidity (spectral) real(DP):: wz_VorA ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1) ! $ \zeta (t+\Delta t) $ . 渦度 (スペクトル). ! Vorticity (spectral) real(DP):: wz_DivA ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1) ! $ D (t+\Delta t) $ . 発散 (スペクトル). ! Divergence (spectral) real(DP):: wz_TempA ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1) ! $ T (t+\Delta t) $ . 温度 (スペクトル). ! Temperature (spectral) real(DP):: w_PiA ((dyn_sp_as%nmax+1)**2) ! $ \pi = \ln p_s (t+\Delta t) $ . 地表面気圧 (スペクトル). ! Surface pressure (spectral) real(DP):: wz_QVapA ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1) ! $ q (t+\Delta t) $ . 比湿 (スペクトル). ! Specific humidity (spectral) real(DP):: wz_VorDiffA ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1) ! $ \mathscr{D}(\zeta) $ . ! 運動量水平拡散による渦度変化 (スペクトル). ! Vorticity tendency by ! horizontal momentum diffusion (spectral) real(DP):: wz_DivDiffA ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1) ! $ \mathscr{D}(D) $ . ! 運動量水平拡散による発散変化 (スペクトル). ! Divergence tendency by ! horizontal momentum diffusion (spectral) !----------------------------------- ! 外部のプロセスの時間変化項のための作業変数 ! Work variables for tendency of external processes real(DP):: xyz_DVorDtWork (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ \DP{\zeta}{t} $ . ! 渦度変化. ! Vorticity tendency real(DP):: xyz_DDivDtWork (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ \DP{D}{t} $ . ! 発散変化. ! Divergence tendency real(DP):: xyz_DTempDtWork (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ \DP{T}{t} $ . ! 温度変化. ! Temperature tendency real(DP):: xyz_DQVapDtWork (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ \DP{q}{t} $ . ! 比湿変化. ! Specific humidity tendency real(DP):: wz_DVorDtWork ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1) ! $ \DD{\zeta}{t} (t) $ . 渦度変化 (スペクトル). ! Vorticity tendency (spectral) real(DP):: wz_DDivDtWork ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1) ! $ \DD{D}{t} (t) $ . 発散変化 (スペクトル). ! Divergence tendency (spectral) real(DP):: wz_DTempDtWork ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1) ! $ \DD{T}{t} (t) $ . 温度変化 (スペクトル). ! Temperature tendency (spectral) real(DP):: wz_DQVapDtWork ((dyn_sp_as%nmax+1)**2, 0:dyn_sp_as%kmax-1) ! $ \DD{q}{t} (t) $ . 比湿変化 (スペクトル). ! Specific humidity tendency (spectral) !----------------------------------- ! エクスプリシット法のための作業変数 ! Work variables for explicit scheme real(DP):: xyz_exWTGPi (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ \underline{W} \Dvect{T} + \Dvect{G} \pi $ . real(DP):: xyz_exHDiv (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! $ \underline{h} \Dvect{D} $ . character(TOKEN):: time_integration_scheme ! 時間積分法. ! 以下の方法を選択可能. ! ! Time integration scheme. ! Available schemes are as follows. ! ! * "Semi-implicit" ! * "Explicit" ! !----------------------------------- ! ヒストリファイルへのデータ出力設定 ! Configure the settings for history data output character(STRING):: name = '' ! 変数名. Variable identifier real:: time ! 時刻. Time type(GT_HISTORY), pointer:: gthist =>null() ! gt4_history モジュール用構造体. ! Derived type for "gt4_history" module !----------------------------------- ! 作業変数 ! Work variables integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpAsDynamics' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp_as % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSPAS83' goto 999 end if !----------------------------------------------------------------- ! 時間積分法の問い合わせ ! Inquire time integration scheme !----------------------------------------------------------------- call Inquire( dyn_as = dyn_sp_as % dyn_as, time_integration_scheme = time_integration_scheme ) ! (out) !----------------------------------------------------------------- ! 時間変化項の計算 ! Calculate tendency terms !----------------------------------------------------------------- !--------------------------------------------- ! 地表面気圧の空間変化項の計算 ! Calculate spatial surface pressure tendency ! call GradPi( dyn_sp = dyn_sp_as % dyn_sp, xy_Ps = xy_PsN, xy_GradLonPi = xy_GradLonPiN, xy_GradLatPi = xy_GradLatPiN ) ! (out) !--------------------------------------------- ! 渦度発散から風速の計算 ! Calculate wind velocity from vorticity and divergence ! call VorDiv2UV( dyn_sp = dyn_sp_as % dyn_sp, xyz_Vor = xyz_VorN, xyz_Div = xyz_DivN, xyz_U = xyz_UN, xyz_V = xyz_VN ) ! (out) !--------------------------------------------- ! 格子点上での非線形力学項の計算 ! Calculate non-linear dynamical terms on grid points ! call NonLinearOnGrid( dyn_as = dyn_sp_as % dyn_as, xyz_U = xyz_UN, xyz_V = xyz_VN, xyz_Vor = xyz_VorN, xyz_Div = xyz_DivN, xyz_Temp = xyz_TempN, xyz_QVap = xyz_QVapN, xy_GradLonPi = xy_GradLonPiN, xy_GradLatPi = xy_GradLatPiN, xyz_UAdv = xyz_UAdvN, xyz_VAdv = xyz_VAdvN, xyz_DTempDt = xyz_DTempDtN, xyz_DQVapDt = xyz_DQVapDtN, xyz_KE = xyz_KEN, xyz_TempUAdv = xyz_TempUAdvN, xyz_TempVAdv = xyz_TempVAdvN, xyr_SigmaDot = xyr_SigmaDotN, xy_DPiDt = xy_DPiDtN, xyz_QVapUAdv = xyz_QVapUAdvN, xyz_QVapVAdv = xyz_QVapVAdvN ) ! (out) !--------------------------------------------- ! スペクトル時間変化項の計算 ! Calculate spectral tendency terms ! call Tendency( dyn_sp = dyn_sp_as % dyn_sp, xyz_UAdv = xyz_UAdvN, xyz_VAdv = xyz_VAdvN, xyz_KE = xyz_KEN, xyz_TempUAdv = xyz_TempUAdvN, xyz_TempVAdv = xyz_TempVAdvN, xyz_DTempDt = xyz_DTempDtN, xy_DPiDt = xy_DPiDtN, xyz_QVapUAdv = xyz_QVapUAdvN, xyz_QVapVAdv = xyz_QVapVAdvN, xyz_DQVapDt = xyz_DQVapDtN, wz_DVorDt = wz_DVorDtN, wz_DDivDt = wz_DDivDtN, wz_DTempDt = wz_DTempDtN, w_DPiDt = w_DPiDtN, wz_DQVapDt = wz_DQVapDtN ) ! (out) !--------------------------------------------- ! エクスプリシット法を用いる際の計算 ! Calculate for explicit scheme ! select case ( LChar( trim( time_integration_scheme ) ) ) case ('explicit') !--------------------------------------------- ! $ \underline{W} \Dvect{T} + \Dvect{G} \pi $ の格子点値の計算 ! Calculate $ \underline{W} \Dvect{T} + \Dvect{G} \pi $ on grid ! call WTplusGPiOnGrid( dyn_as = dyn_sp_as % dyn_as, xyz_Temp = xyz_TempN, xy_Ps = xy_PsN, xyz_exWTGPi = xyz_exWTGPi ) ! (out) !--------------------------------------------- ! $ \underline{h} \Dvect{D} $ の格子点値の計算 ! Calculate $ \underline{h} \Dvect{D} $ on grid ! call HDivOnGrid( dyn_as = dyn_sp_as % dyn_as, xyz_Div = xyz_DivN, xyz_exHDiv = xyz_exHDiv ) ! (out) !--------------------------------------------- ! 発散, 温度 (スペクトル) の時間変化項の修正 ! Modify divergence and temperature tencency (spectral) ! call TendencyExplicit( dyn_sp = dyn_sp_as % dyn_sp, xyz_exWTGPi = xyz_exWTGPi, xyz_exHDiv = xyz_exHDiv, wz_DDivDt = wz_DDivDtN, wz_DTempDt = wz_DTempDtN ) ! (inout) end select !--------------------------------------------- ! 格子点値をスペクトル値へ ( $ t-\Delta t$ ) ! Exchange grid values to spectral values ( $ t-\Delta t$ ) ! call Grid2Spectral( dyn_sp = dyn_sp_as % dyn_sp, xyz_Data = xyz_VorB, wz_Data = wz_VorB ) ! (out) call Grid2Spectral( dyn_sp = dyn_sp_as % dyn_sp, xyz_Data = xyz_DivB, wz_Data = wz_DivB ) ! (out) call Grid2Spectral( dyn_sp = dyn_sp_as % dyn_sp, xyz_Data = xyz_TempB, wz_Data = wz_TempB ) ! (out) call Grid2Spectral( dyn_sp = dyn_sp_as % dyn_sp, xy_Data = xy_PsB, math_func = 'log', w_Data = w_PiB ) ! (out) call Grid2Spectral( dyn_sp = dyn_sp_as % dyn_sp, xyz_Data = xyz_QVapB, wz_Data = wz_QVapB ) ! (out) !----------------------------------------------------------------- ! 外部プロセスによる時間変化格子データをスペクトルデータへ変換 ! Convert tendency data on grid by external processes into spectral data !----------------------------------------------------------------- xyz_DVorDtWork = 0.0_DP xyz_DDivDtWork = 0.0_DP xyz_DTempDtWork = 0.0_DP xyz_DQVapDtWork = 0.0_DP if ( present(xyz_DVorDt ) ) xyz_DVorDtWork = xyz_DVorDt if ( present(xyz_DDivDt ) ) xyz_DDivDtWork = xyz_DDivDt if ( present(xyz_DTempDt) ) xyz_DTempDtWork = xyz_DTempDt if ( present(xyz_DQVapDt) ) xyz_DQVapDtWork = xyz_DQVapDt !--------------------------------------------- ! 格子点値をスペクトル値へ ( $ t-\Delta t$ ) ! Exchange grid values to spectral values ( $ t-\Delta t$ ) ! call Grid2Spectral( dyn_sp = dyn_sp_as % dyn_sp, xyz_Data = xyz_DVorDtWork, wz_Data = wz_DVorDtWork ) ! (out) call Grid2Spectral( dyn_sp = dyn_sp_as % dyn_sp, xyz_Data = xyz_DDivDtWork, wz_Data = wz_DDivDtWork ) ! (out) call Grid2Spectral( dyn_sp = dyn_sp_as % dyn_sp, xyz_Data = xyz_DTempDtWork, wz_Data = wz_DTempDtWork ) ! (out) call Grid2Spectral( dyn_sp = dyn_sp_as % dyn_sp, xyz_Data = xyz_DQVapDtWork, wz_Data = wz_DQVapDtWork ) ! (out) !----------------------------------------------------------------- ! 時間積分 ! Time integration !----------------------------------------------------------------- call TimeIntegration( dyn_as = dyn_sp_as % dyn_as, wz_DVorDtN = wz_DVorDtN + wz_DVorDtWork, wz_DDivDtN = wz_DDivDtN + wz_DDivDtWork, wz_DTempDtN = wz_DTempDtN + wz_DTempDtWork, wz_DQVapDtN = wz_DQVapDtN + wz_DQVapDtWork, w_DPiDtN = w_DPiDtN, wz_VorB = wz_VorB, wz_DivB = wz_DivB, wz_TempB = wz_TempB, wz_QVapB = wz_QVapB, w_PiB = w_PiB, wz_VorA = wz_VorA, wz_DivA = wz_DivA, wz_TempA = wz_TempA, wz_QVapA = wz_QVapA, w_PiA = w_PiA ) ! (out) !--------------------------------------------- ! スペクトル値を格子点値へ ( $ t+\Delta t$ ) ! Exchange spectral values to grid values ( $ t+\Delta t$ ) ! call Spectral2Grid( dyn_sp = dyn_sp_as % dyn_sp, wz_Data = wz_VorA, xyz_Data = xyz_VorA ) ! (out) call Spectral2Grid( dyn_sp = dyn_sp_as % dyn_sp, wz_Data = wz_DivA, xyz_Data = xyz_DivA ) ! (out) call Spectral2Grid( dyn_sp = dyn_sp_as % dyn_sp, wz_Data = wz_TempA, xyz_Data = xyz_TempA ) ! (out) call Spectral2Grid( dyn_sp = dyn_sp_as % dyn_sp, w_Data = w_PiA, math_func = 'exp', xy_Data = xy_PsA ) ! (out) call Spectral2Grid( dyn_sp = dyn_sp_as % dyn_sp, wz_Data = wz_QVapA, xyz_Data = xyz_QVapA ) ! (out) call VorDiv2UV( dyn_sp = dyn_sp_as % dyn_sp, xyz_Vor = xyz_VorA, xyz_Div = xyz_DivA, xyz_U = xyz_UA, xyz_V = xyz_VA ) ! (out) !----------------------------------------------------------------- ! 拡散による補正 ! Correction by diffusion !----------------------------------------------------------------- !--------------------------------------------- ! 運動量水平拡散による渦度発散の時間変化 ! Vorticity and divergence tendency by ! horizontal diffusion of momentum ! call DiffusionVorDiv( dyn_sp = dyn_sp_as % dyn_sp, wz_Vor = wz_VorA, wz_Div = wz_DivA, wz_VorDiff = wz_VorDiffA, wz_DivDiff = wz_DivDiffA ) ! (out) !--------------------------------------------- ! 運動量水平拡散による摩擦熱補正 ! Frictional thermal correction by horizontal momentum diffusion ! call DiffusionCorrectTemp( dyn_sp = dyn_sp_as % dyn_sp, wz_VorDiff = wz_VorDiffA, wz_DivDiff = wz_DivDiffA, xyz_U = xyz_UA, xyz_V = xyz_VA, xyz_Temp = xyz_TempA ) ! (inout) !---------------------------------------------------------------- ! ヒストリファイルへのデータ出力 ! History data output !---------------------------------------------------------------- !------------------------- ! xyr_SigmaDot の出力 ! Output "xyr_SigmaDot" name = 'SigmaDot' ! 出力のチェック. ! * gthist (gt4_history#GT_HISTORY), time (単精度実数型) が設定される. ! Check for output. ! * "gthist" (gt4_history#GT_HISTORY), time (real) are configured. call output_check ! これは内部サブルーチン. This is internal subroutine ! 出力データを引数 array に渡す. ! Give output data to argument "array" if ( associated( gthist ) ) then call HistoryPut( history = gthist, varname = name, array = xyr_SigmaDotN, time = time, quiet = .false., err = err ) ! (out) end if !------------------------- ! xy_DPiDt の出力 ! Output "xy_DPiDt" name = 'DPiDt' ! 出力のチェック. ! * gthist (gt4_history#GT_HISTORY), time (単精度実数型) が設定される. ! Check for output. ! * "gthist" (gt4_history#GT_HISTORY), time (real) are configured. call output_check ! これは内部サブルーチン. This is internal subroutine ! 出力データを引数 array に渡す. ! Give output data to argument "array" if ( associated( gthist ) ) then call HistoryPut( history = gthist, varname = name, array = xy_DPiDtN, time = time, quiet = .false., err = err ) ! (out) end if !----------------------------------------------------------------- ! 時刻の更新 ! Update time !----------------------------------------------------------------- dyn_sp_as % current_time = dyn_sp_as % current_time + dyn_sp_as % delta_time !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) contains subroutine output_check ! ! 変数 *name* を出力するかどうかをチェックします. ! 出力に関する情報は dyn_sp_as % gthstnml から取り出されます. ! ! 変数 *name* に関して出力するよう設定されている場合には, ! *gthist* に出力先ファイルの gt4_history#GT_HISTORY ! 型変数を結合させます. そうでない場合は, *gthist* を空状態にします. ! ! また, 現在時刻を *time* に設定します. ! ! Check whether to output variable *name*. ! Information about output is taken out of "dyn_sp_as % gthstnml". ! ! When output is done for the variable *name*, *gthist* is ! associated with "gt4_history#GT_HISTORY" variable of ! the output file. Otherwise, *gthist* is nullified. ! ! Moreover, current time is set to *time*. ! character(TOKEN):: interval_unit ! ヒストリデータの出力間隔の単位. ! Unit for interval of history data output continue nullify( gthist ) time = 0.0 if ( HstNmlInfoOutputValid( dyn_sp_as % gthstnml, name ) ) then call HstNmlInfoInquire( gthstnml = dyn_sp_as % gthstnml, name = name, interval_unit = interval_unit ) ! (out) time = real( EvalbyUnit( dyn_sp_as % current_time, interval_unit ) ) if ( present_and_true( historyput_flag ) ) time = 0.0 call HstNmlInfoAssocGtHist( gthstnml = dyn_sp_as % gthstnml, name = name, history = gthist, err = err ) ! (out) if ( present_and_true( err ) ) then nullify( gthist ) return end if if ( .not. HistoryInitialized( gthist ) ) nullify( gthist ) end if end subroutine output_check end subroutine DynSpAsDynamics
Subroutine : | |||
nmlfile : | character(*), intent(in)
| ||
gthstnml : | type(GTHST_NMLINFO), intent(inout)
| ||
err : | logical, intent(out), optional
|
NAMELIST ファイル nmlfile から値を入力するための 内部サブルーチンです. DynSpAsCreate 内で呼び出されることを 想定しています.
値が NAMELIST ファイル内で指定されていない場合には, 入力された値がそのまま返ります.
なお, nmlfile に空文字が与えられた場合, または 与えられた nmlfile を読み込むことができない場合, プログラムはエラーを発生させます.
This is an internal subroutine to input values from NAMELIST file nmlfile. This subroutine is expected to be called by "DynSpAsCreate".
A value not specified in NAMELIST file is returned without change.
If nmlfile is empty, or nmlfile can not be read, error is occurred.
This procedure input/output NAMELIST#dyn_spectral_as83_history_nml .
subroutine DynSpAsNmlRead( nmlfile, gthstnml, err ) ! ! NAMELIST ファイル *nmlfile* から値を入力するための ! 内部サブルーチンです. DynSpAsCreate 内で呼び出されることを ! 想定しています. ! ! 値が NAMELIST ファイル内で指定されていない場合には, ! 入力された値がそのまま返ります. ! ! なお, *nmlfile* に空文字が与えられた場合, または ! 与えられた *nmlfile* を読み込むことができない場合, ! プログラムはエラーを発生させます. ! ! This is an internal subroutine to input values from ! NAMELIST file *nmlfile*. This subroutine is expected to be ! called by "DynSpAsCreate". ! ! A value not specified in NAMELIST file is returned ! without change. ! ! If *nmlfile* is empty, or *nmlfile* can not be read, ! 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_iounit, only: FileOpen use dc_message, only: MessageNotify use dc_present, only: present_and_true use dc_date, only: DCDiffTimeCreate use dc_error, only: StoreError, DC_NOERR, DC_ENOFILEREAD, DC_ENOASSOC use gt4_history_nmlinfo, only: GTHST_NMLINFO, HstNmlInfoAdd, HstNmlInfoInquire, HstNmlInfoInitialized, HstNmlInfoPutLine implicit none character(*), intent(in):: nmlfile ! NAMELIST ファイルの名称. ! NAMELIST file name !!$ real(DP), intent(inout):: CoefAlpha !!$ ! $ \alpha $ . 係数. Coefficient !!$ !!$ character(*), intent(inout):: key00_ !!$ character(TOKEN):: key00 !!$ ! キーワード. Keyword !!$ type(GTHST_NMLINFO), intent(inout):: gthstnml ! NAMELIST#dcmodel_sample_code_history_nml ! から入手される個別のデータ出力情報. ! ! 初期設定やデフォルト値の設定などを ! 行った後に与えること. ! ! Individual data output information from ! "NAMELIST#dcmodel_sample_code_history_nml". ! ! Before this argument is given to ! this procedure, initialize and ! configure the defaut settings. 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. !!$ namelist /dyn_spectral_as83_nml/ & !!$ & CoefAlpha, key00 ! dyn_spectral_as83 モジュール用 ! NAMELIST 変数群名. ! ! dyn_spectral_as83#DynSpAsCreate を使用する際に, ! オプショナル引数 *nmlfile* へ NAMELIST ! ファイル名を指定することで, そのファイルから ! この NAMELIST 変数群を読み込みます. ! ! NAMELIST group name for ! "dyn_spectral_as83" module. ! ! If a NAMELIST filename is specified to ! an optional argument *nmlfile* ! when "dyn_spectral_as83#DynSpAsCreate" is used, ! this NAMELIST group is loaded from ! the file. character(STRING):: name ! 変数名. ! 空白の場合には, この他の設定値は ! dcmodel_sample_code モジュールにおいて ! 出力されるデータ全ての ! デフォルト値となります. ! ! "Data1,Data2" のようにカンマで区切って複数 ! の変数を指定することも可能です. ! ! Variable identifier. ! If blank is given, other values are ! used as default values of output data ! in "dcmodel_sample_code". ! ! Multiple variables can be specified ! as "Data1,Data2" too. Delimiter is comma. character(STRING):: file ! 出力ファイル名. ! これはデフォルト値としては使用されません. ! *name* に値が設定されている時のみ有効です. ! ! Output file name. ! This is not used as default value. ! This value is valid only when *name* is ! specified. real:: interval_value ! ヒストリデータの出力間隔の数値. ! 負の値を与えると, 出力を抑止します. ! Numerical value for interval of history data output ! Negative values suppresses output. character(TOKEN):: interval_unit ! ヒストリデータの出力間隔の単位. ! Unit for interval of history data output character(TOKEN):: precision ! ヒストリデータの精度. ! Precision of history data logical:: average ! 出力データの平均化フラグ. ! Flag for average of output data character(STRING):: fileprefix ! ヒストリデータのファイル名の接頭詞. ! Prefixes of history data filenames namelist /dyn_spectral_as83_history_nml/ name, file, interval_value, interval_unit, precision, fileprefix, average ! dyn_spectral_as83 モジュールのヒストリデータ用 ! NAMELIST 変数群名. ! ! dyn_spectral_as83#DynSpAsCreate を使用する際に, ! オプショナル引数 *nmlfile* へ NAMELIST ! ファイル名を指定することで, そのファイルから ! この NAMELIST 変数群を読み込みます. ! ! NAMELIST group name for ! history data of "dyn_spectral_as83" module. ! ! If a NAMELIST filename is specified to ! an optional argument *nmlfile* ! when "dyn_spectral_as83#DynSpAsCreate" is used, ! this NAMELIST group is loaded from ! the file. !----------------------------------- ! 作業変数 ! Work variables integer:: stat character(STRING):: cause_c integer:: unit_nml ! NAMELIST ファイルオープン用装置番号. ! Unit number for NAMELIST file open integer:: iostat_nml ! NAMELIST 読み込み時の IOSTAT. ! IOSTAT of NAMELIST read character(TOKEN):: pos_nml ! NAMELIST 読み込み時のファイル位置. ! File position of NAMELIST read character(*), parameter:: subname = 'DynSpAsNmlRead' continue call BeginSub( subname ) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- !----------------------------------------------------------------- ! 文字型引数を NAMELIST 変数群へ代入 ! Substitute character arguments to NAMELIST group !----------------------------------------------------------------- !!$ key00 = key00_ !---------------------------------------------------------------- ! NAMELIST ファイルのオープン ! Open NAMELIST file !---------------------------------------------------------------- call FileOpen( unit = unit_nml, file = nmlfile, mode = 'r', err = err ) ! (out) if ( present_and_true(err) ) then stat = DC_ENOFILEREAD cause_c = nmlfile goto 999 end if !----------------------------------------------------------------- ! NAMELIST 変数群の取得 ! Get NAMELIST group !----------------------------------------------------------------- !------------------------- ! 係数などの取得 ! Get coefficients etc. !!$ rewind( unit_nml ) !!$ read( unit = unit_nml, & ! (in) !!$ & nml = dyn_spectral_as83_nml, iostat = iostat_nml ) ! (out) !!$ if ( iostat_nml == 0 ) then !!$ call MessageNotify( 'M', subname, & !!$ & 'NAMELIST group "%c" is loaded from "%c".', & !!$ & c1 = 'dyn_spectral_as83_nml', c2 = trim(nmlfile) ) !!$ write(STDOUT, nml = dyn_spectral_as83_nml) !!$ else !!$ call MessageNotify( 'W', subname, & !!$ & 'NAMELIST group "%c" is not found in "%c" (iostat=%d).', & !!$ & c1 = 'dyn_spectral_as83_nml', c2 = trim(nmlfile), & !!$ & i = (/iostat_nml/) ) !!$ end if !------------------------- ! 出力データの個別情報の取得 ! Get individual information of output data rewind( unit_nml ) iostat_nml = 0 pos_nml = '' do while ( trim(pos_nml) /= 'APPEND' .and. iostat_nml == 0 ) name = '' file = '' call HstNmlInfoInquire( gthstnml = gthstnml, interval_value = interval_value, interval_unit = interval_unit, precision = precision, average = average, fileprefix = fileprefix ) ! (out) read( unit = unit_nml, nml = dyn_spectral_as83_history_nml, iostat = iostat_nml ) ! (out) inquire( unit = unit_nml, position = pos_nml ) ! (out) if ( iostat_nml == 0 ) then call MessageNotify( 'M', subname, 'NAMELIST group "%c" is loaded from "%c".', c1='dcmodel_sample_code_history_nml', c2=trim(nmlfile) ) write(STDOUT, nml = dyn_spectral_as83_history_nml) call HstNmlInfoAdd( gthstnml = gthstnml, name = name, file = file, interval_value = interval_value, interval_unit = interval_unit, precision = precision, average = average, fileprefix = fileprefix ) ! (in) else call MessageNotify( 'W', subname, 'NAMELIST group "%c" is not found in "%c" any more (iostat=%d).', c1='dcmodel_sample_code_history_nml', c2=trim(nmlfile), i = (/iostat_nml/) ) end if end do close( unit_nml ) !----------------------------------------------------------------- ! NAMELIST 変数群を文字型引数へ代入 ! Substitute NAMELIST group to character arguments !----------------------------------------------------------------- !!$ key00_ = key00 !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError( stat, subname, err, cause_c ) call EndSub( subname ) end subroutine DynSpAsNmlRead
Subroutine : | |||
dyn_sp_as : | type(DYNSPAS83), intent(inout) | ||
current_time_value : | real(DP), intent(in)
| ||
current_time_unit : | character(*), intent(in)
| ||
err : | logical, intent(out), optional
|
dyn_sp_as に対して時刻の設定を行います.
ヒストリデータを出力している場合には, ヒストリデータの 出力時刻も設定します. 一度でもこのサブルーチンを呼んだ場合には, それ以後のヒストリデータ出力前にこのサブルーチンを呼び出し, 時刻の設定を行ってください. また, データを出力するサブルーチンに対しても オプショナル引数 historyput_flag に .true. を与えてください.
なお, 与えられた dyn_sp_as が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Set time to dyn_sp_as.
When history data are output, the output time of history data are specified. Once this subroutine is called, the time of history data must be specified by this routine before history data output In additional, give ".true." to an optional argument "historyput_flag" of a data output subroutine.
If dyn_sp_as is not initialized by "Create" yet, error is occurred.
subroutine DynSpAsSetTime( dyn_sp_as, current_time_value, current_time_unit, err ) ! ! *dyn_sp_as* に対して時刻の設定を行います. ! ! ヒストリデータを出力している場合には, ヒストリデータの ! 出力時刻も設定します. 一度でもこのサブルーチンを呼んだ場合には, ! それ以後のヒストリデータ出力前にこのサブルーチンを呼び出し, ! 時刻の設定を行ってください. ! また, データを出力するサブルーチンに対しても ! オプショナル引数 historyput_flag に .true. を与えてください. ! ! なお, 与えられた *dyn_sp_as* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Set time to *dyn_sp_as*. ! ! When history data are output, the output time of history data ! are specified. ! Once this subroutine is called, the time of history data must be ! specified by this routine before history data output ! In additional, give ".true." to an optional argument ! "historyput_flag" of a data output subroutine. ! ! If *dyn_sp_as* 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_date, only: DCDiffTimeCreate, EvalbyUnit use dc_error, only: StoreError, DC_NOERR, DC_ENOTINIT use gt4_history_nmlinfo, only: HstNmlInfoAdd, HstNmlInfoInquire, HstNmlInfoNames, HstNmlInfoAssocGtHist, HstNmlInfoOutputStepDisable, HstNmlInfoPutLine use gt4_history, only: HistorySetTime, HistoryInitialized implicit none type(DYNSPAS83), intent(inout):: dyn_sp_as real(DP), intent(in):: current_time_value ! 現在時刻の数値. Numerical value of current time character(*), intent(in):: current_time_unit ! 現在時刻の単位. Unit of current time 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. !----------------------------------- ! ヒストリファイルへのデータ出力設定 ! Configure the settings for history data output character(STRING):: name = '' ! 変数名. Variable identifier character(TOKEN):: interval_unit ! ヒストリデータの出力間隔の単位. ! Unit for interval of history data output character(STRING):: varnames ! 変数名リスト. ! List of variables character(TOKEN), pointer:: varnames_array(:) =>null() ! 変数名リスト配列. ! List of variables (array) integer:: i, vnmax type(GT_HISTORY), pointer:: gthist =>null() ! gt4_history モジュール用構造体. ! Derived type for "gt4_history" module !----------------------------------- ! 作業変数 ! Work variables integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpAsSetTime' continue call BeginSub( subname ) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if ( .not. dyn_sp_as % initialized ) then stat = DC_ENOTINIT cause_c = 'DYNSPAS83' goto 999 end if !----------------------------------------------------------------- ! 時刻設定 ! Configure time !----------------------------------------------------------------- call DCDiffTimeCreate( diff = dyn_sp_as % current_time, value = current_time_value, unit = current_time_unit ) ! (in) !----------------------------------------------------------------- ! ヒストリファイルへのデータの時刻設定 ! Configure the time of history data !----------------------------------------------------------------- varnames = HstNmlInfoNames( dyn_sp_as % gthstnml ) call Split( str = varnames, sep = ',', carray = varnames_array ) ! (out) vnmax = size( varnames_array ) do i = 1, vnmax name = varnames_array(i) if ( trim( name ) == '' ) exit call HstNmlInfoOutputStepDisable( gthstnml = dyn_sp_as % gthstnml, name = name, err = err ) ! (out) nullify( gthist ) call HstNmlInfoAssocGtHist( gthstnml = dyn_sp_as % gthstnml, name = name, history = gthist, err = err ) ! (out) if ( HistoryInitialized( gthist ) ) then call HstNmlInfoInquire( gthstnml = dyn_sp_as % gthstnml, name = name, interval_unit = interval_unit ) ! (out) call HistorySetTime( history = gthist, time = real( EvalbyUnit( dyn_sp_as % current_time, interval_unit) ) ) ! (in) end if end do !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError( stat, subname, err, cause_c ) call EndSub( subname ) end subroutine DynSpAsSetTime
Subroutine : | |||
dyn_sp_as : | type(DYNSPAS83), intent(inout) | ||
xyz_Vor(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(out)
| ||
xyz_Div(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(out)
| ||
xyz_U(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xyz_V(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
err : | logical, intent(out), optional
|
与えられた東西風速と南北風速から渦度と発散を計算します.
なお, 与えられた dyn_sp_as が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculate vorticity and divergence from given zonal and meridional wind.
If dyn_sp_as is not initialized by "Create" yet, error is occurred.
subroutine DynSpAsUV2VorDiv( dyn_sp_as, xyz_Vor, xyz_Div, xyz_U, xyz_V, err ) ! ! 与えられた東西風速と南北風速から渦度と発散を計算します. ! ! なお, 与えられた *dyn_sp_as* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Calculate vorticity and divergence ! from given zonal and meridional wind. ! ! If *dyn_sp_as* is not initialized by "Create" yet, ! error is occurred. ! use dc_trace, only: BeginSub, EndSub use dc_types, only: STRING, STDOUT, DP use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT use dyn_spectral, only: UV2VorDiv implicit none type(DYNSPAS83), intent(inout):: dyn_sp_as real(DP), intent(in):: xyz_U (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! 東西風速 $ u $ real(DP), intent(in):: xyz_V (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! 南北風速 $ v $ real(DP), intent(out):: xyz_Vor (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! 渦度 $ \zeta $ real(DP), intent(out):: xyz_Div (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! 発散 $ D $ 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. integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpAsUV2VorDiv' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp_as % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSPAS83' goto 999 end if !----------------------------------------------------------------- ! 風速から渦度発散の計算 ! Calculate vorticity and divergence from wind velocity !----------------------------------------------------------------- call UV2VorDiv( dyn_sp = dyn_sp_as % dyn_sp, xyz_U = xyz_U, xyz_V = xyz_V, xyz_Vor = xyz_Vor, xyz_Div = xyz_Div ) ! (out) !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpAsUV2VorDiv
Subroutine : | |||
dyn_sp_as : | type(DYNSPAS83), intent(inout) | ||
xyz_Vor(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xyz_Div(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(in)
| ||
xyz_U(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(out)
| ||
xyz_V(0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
与えられた渦度と発散から東西風速と南北風速を計算します.
なお, 与えられた dyn_sp_as が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculate zonal and meridional wind from given vorticity and divergence.
If dyn_sp_as is not initialized by "Create" yet, error is occurred.
subroutine DynSpAsVorDiv2UV( dyn_sp_as, xyz_Vor, xyz_Div, xyz_U, xyz_V, err ) ! ! 与えられた渦度と発散から東西風速と南北風速を計算します. ! ! なお, 与えられた *dyn_sp_as* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Calculate zonal and meridional wind from given ! vorticity and divergence. ! ! If *dyn_sp_as* is not initialized by "Create" yet, ! error is occurred. ! use dc_trace, only: BeginSub, EndSub use dc_types, only: STRING, STDOUT, DP use dc_string, only: PutLine, Printf, Split, StrInclude, StoA, JoinChar use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT use dyn_spectral, only: VorDiv2UV implicit none type(DYNSPAS83), intent(inout):: dyn_sp_as real(DP), intent(in):: xyz_Vor (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! 渦度 $ \zeta $ real(DP), intent(in):: xyz_Div (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! 発散 $ D $ real(DP), intent(out):: xyz_U (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! 東西風速 $ u $ real(DP), intent(out):: xyz_V (0:dyn_sp_as%imax-1, 0:dyn_sp_as%jmax-1, 0:dyn_sp_as%kmax-1) ! 南北風速 $ v $ 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. integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpAsVorDiv2UV' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp_as % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSPAS83' goto 999 end if !----------------------------------------------------------------- ! 渦度発散から風速の計算 ! Calculate wind velocity from vorticity and divergence !----------------------------------------------------------------- call VorDiv2UV( dyn_sp = dyn_sp_as % dyn_sp, xyz_Vor = xyz_Vor, xyz_Div = xyz_Div, xyz_U = xyz_U, xyz_V = xyz_V ) ! (out) !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpAsVorDiv2UV
Subroutine : | |||
nmlfile : | character(*), intent(in)
| ||
gthstnml : | type(GTHST_NMLINFO), intent(inout)
| ||
err : | logical, intent(out), optional
|
NAMELIST ファイル nmlfile から値を入力するための 内部サブルーチンです. DynSpAsCreate 内で呼び出されることを 想定しています.
値が NAMELIST ファイル内で指定されていない場合には, 入力された値がそのまま返ります.
なお, nmlfile に空文字が与えられた場合, または 与えられた nmlfile を読み込むことができない場合, プログラムはエラーを発生させます.
This is an internal subroutine to input values from NAMELIST file nmlfile. This subroutine is expected to be called by "DynSpAsCreate".
A value not specified in NAMELIST file is returned without change.
If nmlfile is empty, or nmlfile can not be read, error is occurred.
This procedure input/output NAMELIST#dyn_spectral_as83_history_nml .
Alias for DynSpAsNmlRead
Subroutine : |
変数 name を出力するかどうかをチェックします. 出力に関する情報は dyn_sp_as % gthstnml から取り出されます.
変数 name に関して出力するよう設定されている場合には, gthist に出力先ファイルの gt4_history#GT_HISTORY 型変数を結合させます. そうでない場合は, gthist を空状態にします.
また, 現在時刻を time に設定します.
Check whether to output variable name. Information about output is taken out of "dyn_sp_as % gthstnml".
When output is done for the variable name, gthist is associated with "gt4_history#GT_HISTORY" variable of the output file. Otherwise, gthist is nullified.
Moreover, current time is set to time.
subroutine output_check ! ! 変数 *name* を出力するかどうかをチェックします. ! 出力に関する情報は dyn_sp_as % gthstnml から取り出されます. ! ! 変数 *name* に関して出力するよう設定されている場合には, ! *gthist* に出力先ファイルの gt4_history#GT_HISTORY ! 型変数を結合させます. そうでない場合は, *gthist* を空状態にします. ! ! また, 現在時刻を *time* に設定します. ! ! Check whether to output variable *name*. ! Information about output is taken out of "dyn_sp_as % gthstnml". ! ! When output is done for the variable *name*, *gthist* is ! associated with "gt4_history#GT_HISTORY" variable of ! the output file. Otherwise, *gthist* is nullified. ! ! Moreover, current time is set to *time*. ! character(TOKEN):: interval_unit ! ヒストリデータの出力間隔の単位. ! Unit for interval of history data output continue nullify( gthist ) time = 0.0 if ( HstNmlInfoOutputValid( dyn_sp_as % gthstnml, name ) ) then call HstNmlInfoInquire( gthstnml = dyn_sp_as % gthstnml, name = name, interval_unit = interval_unit ) ! (out) time = real( EvalbyUnit( dyn_sp_as % current_time, interval_unit ) ) if ( present_and_true( historyput_flag ) ) time = 0.0 call HstNmlInfoAssocGtHist( gthstnml = dyn_sp_as % gthstnml, name = name, history = gthist, err = err ) ! (out) if ( present_and_true( err ) ) then nullify( gthist ) return end if if ( .not. HistoryInitialized( gthist ) ) nullify( gthist ) end if end subroutine output_check
Subroutine : |
変数 name に関して出力ファイルの初期設定を行います. 出力ファイル名や出力間隔などの情報は dyn_sp_as % gthstnml から取り出されます.
変数 name に関して出力が行われる場合には, gthist に出力先ファイルの gt4_history#GT_HISTORY 型変数を結合させます. そうでない場合は, gthist を空状態にします.
また, 出力データの精度を precision に, 出力データ平均化の可否を average に設定します.
標準出力に表示される変数リスト registered_varnames に name, longname, dims, units が登録されます.
An output file is initialized for a variable name. Information such as the output filename and output intervals is taken out of "dyn_sp_as % gthstnml".
When output is done for the variable name, gthist is associated with the "gt4_history#GT_HISTORY" variable of the output file. Otherwise, gthist is nullified.
Moreover, the accuracy of output data is set to precision, and right or wrong of averaging the output data is set to average.
name, longname, dims, units are registered to a list of variables registered_varnames that is printed to standard output.
subroutine output_init ! ! 変数 *name* に関して出力ファイルの初期設定を行います. ! 出力ファイル名や出力間隔などの情報は dyn_sp_as % gthstnml ! から取り出されます. ! ! 変数 *name* に関して出力が行われる場合には, ! *gthist* に出力先ファイルの gt4_history#GT_HISTORY ! 型変数を結合させます. そうでない場合は, *gthist* を空状態にします. ! ! また, 出力データの精度を precision に, ! 出力データ平均化の可否を average に設定します. ! ! 標準出力に表示される変数リスト *registered_varnames* に ! *name*, *longname*, *dims*, *units* が登録されます. ! ! An output file is initialized for a variable *name*. ! Information such as the output filename and output intervals ! is taken out of "dyn_sp_as % gthstnml". ! ! When output is done for the variable *name*, *gthist* is ! associated with the "gt4_history#GT_HISTORY" variable of ! the output file. Otherwise, *gthist* is nullified. ! ! Moreover, the accuracy of output data is set to *precision*, and ! right or wrong of averaging the output data is set to *average*. ! ! *name*, *longname*, *dims*, *units* are registered to ! a list of variables *registered_varnames* that is printed to ! standard output. ! use dc_date, only: DCDiffTimeCreate, EvalSec, EvalByUnit use gt4_history_nmlinfo, only: GTHST_NMLINFO, HstNmlInfoInitialized, HstNmlInfoInquire, HstNmlInfoOutputValid, HstNmlInfoAssocGtHist, HstNmlInfoPutLine use gt4_history, only: GT_HISTORY, HistoryCreate, HistoryAddVariable, HistoryPut, HistoryAddAttr, HistoryInitialized !----------------------------------- ! 作業変数 ! Work variables character(STRING):: file ! ヒストリデータのファイル名. ! History data filenames character(STRING):: dims_str ! 座標軸のリスト. ! List of axes real:: interval_value ! ヒストリデータの出力間隔の数値. ! Numerical value for interval of history data output character(TOKEN):: interval_unit ! ヒストリデータの出力間隔の単位. ! Unit for interval of history data output real(DP), parameter:: PI = 3.1415926535897930_DP ! $ \pi $ . 円周率. Circular constant continue !----------------------------------------------------------------- ! 標準出力に表示される変数の登録 ! Register a variable name for print to standard output !----------------------------------------------------------------- if ( allocated(dims) ) then dims_str = JoinChar( dims, ',' ) else dims_str = '' end if call DCHashPut( hashv = registered_varnames, key = name, value = trim( longname ) // ' [' // trim( units ) // '] {' // trim( dims_str ) // '}' ) ! (in) !----------------------------------------------------------------- ! 変数の初期化 ! Initialize variable !----------------------------------------------------------------- nullify( gthist ) precision = 'float' average = .false. !----------------------------------------------------------------- ! 出力が有効かどうかを確認する ! Confirm whether the output is effective !----------------------------------------------------------------- if ( .not. HstNmlInfoOutputValid( dyn_sp_as % gthstnml, name ) ) then return end if !----------------------------------------------------------------- ! GT_HISTORY 変数の取得 ! Get "GT_HISTORY" variable !----------------------------------------------------------------- call HstNmlInfoAssocGtHist( gthstnml = dyn_sp_as % gthstnml, name = name, history = gthist, err = err ) ! (out) if ( present_and_true( err ) ) return call HstNmlInfoInquire( gthstnml = dyn_sp_as % gthstnml, name = name, precision = precision, average = average, err = err ) ! (out) if ( present_and_true( err ) ) return !----------------------------------------------------------------- ! GT_HISTORY 変数の初期設定の確認 ! Check initialization of "GT_HISTORY" variable !----------------------------------------------------------------- if ( HistoryInitialized( gthist ) ) then !--------------------------------------------------------------- ! HistoryAddVariable による変数作成 ! A variable is created by "HistoryAddVariable" !--------------------------------------------------------------- call HistoryAddVariable( history = gthist, varname = name, dims = dims, longname = longname, units = units, xtype = precision, average = average ) ! (in) return end if !----------------------------------------------------------------- ! HistoryCreate のための設定値の取得 ! Get the settings for "HistoryCreate" !----------------------------------------------------------------- call HstNmlInfoInquire( gthstnml = dyn_sp_as % gthstnml, name = name, file = file, interval_unit = interval_unit, interval_value = interval_value, err = err ) ! (out) if ( present_and_true( err ) ) return !----------------------------------------------------------------- ! HistoryCreate によるファイル作成 ! Files are created by "HistoryCreate" !----------------------------------------------------------------- call HistoryCreate( history = gthist, file = file, title = 'Dynamical core calculation', source = 'dcpam4 : ' // trim(version), institution = 'GFD Dennou Club', dims = StoA('lon', 'lat', 'sig', 'sigm', 'time'), dimsizes = (/ dyn_sp_as % imax, dyn_sp_as % jmax, dyn_sp_as % kmax, dyn_sp_as % kmax+1, 0 /), longnames = StoA('longitude', 'latitude', 'sigma at layer midpoints', 'sigma at layer end-points (half level)', 'time'), units = StoA( 'degree_east', 'degree_north', '1', '1', interval_unit ), origin = real( EvalbyUnit( dyn_sp_as % current_time, interval_unit) ), interval = interval_value, err = err ) ! (out) if ( present_and_true( err ) ) then nullify( gthist ) return end if call HistoryAddAttr( history = gthist, varname = 'lon', attrname = 'standard_name', value = 'longitude' ) ! (in) call HistoryAddAttr( history = gthist, varname = 'lat', attrname = 'standard_name', value = 'latitude' ) ! (in) call HistoryAddAttr( history = gthist, varname = 'sig', attrname = 'standard_name', value = 'atmosphere_sigma_coordinate' ) ! (in) call HistoryAddAttr( history = gthist, varname = 'sigm', attrname = 'standard_name', value = 'atmosphere_sigma_coordinate' ) ! (in) call HistoryAddAttr( history = gthist, varname = 'time', attrname = 'standard_name', value = 'time' ) ! (in) call HistoryAddAttr( history = gthist, varname = 'sig', attrname = 'positive', value = 'down' ) ! (in) call HistoryAddAttr( history = gthist, varname = 'sigm', attrname = 'positive', value = 'down' ) ! (in) call HistoryPut( history = gthist, varname = 'lon', array = x_Lon / PI * 180.0_DP ) ! (in) call HistoryPut( history = gthist, varname = 'lat', array = y_Lat / PI * 180.0_DP ) ! (in) call HistoryPut( history = gthist, varname = 'sig', array = z_Sigma ) ! (in) call HistoryPut( history = gthist, varname = 'sigm', array = r_Sigma ) ! (in) call HistoryAddVariable( history = gthist, varname = 'lon_weight', dims = StoA('lon'), longname = 'weight for integration in longitude', units = 'radian', xtype = 'double' ) ! (in) call HistoryAddAttr( history = gthist, varname = 'lon', attrname = 'gt_calc_weight', value = 'lon_weight' ) ! (in) call HistoryPut( history = gthist, varname = 'lon_weight', array = x_Lon_Weight ) ! (in) call HistoryAddVariable( history = gthist, varname = 'lat_weight', dims = StoA('lat'), longname = 'weight for integration in latitude', units = 'radian', xtype = 'double' ) ! (in) call HistoryAddAttr( history = gthist, varname = 'lat', attrname = 'gt_calc_weight', value = 'lat_weight' ) ! (in) call HistoryPut( history = gthist, varname = 'lat_weight', array = y_Lat_Weight ) ! (in) !----------------------------------------------------------------- ! HistoryAddVariable による変数作成 ! A variable is created by "HistoryAddVariable" !----------------------------------------------------------------- if ( HistoryInitialized( gthist ) ) then call HistoryAddVariable( history = gthist, varname = name, dims = dims, longname = longname, units = units, xtype = precision, average = average ) ! (in) else nullify( gthist ) end if end subroutine output_init
Constant : | |
version = ’$Name: dcpam4-20080626 $’ // ’$Id: dyn_spectral_as83.f90,v 1.34 2008-06-14 13:32:29 morikawa Exp $’ : | character(*), parameter |