Class | dyn_spectral |
In: |
dynamics/dyn_spectral.f90
|
Note that Japanese and English are described in parallel.
力学コアモジュールです. 水平離散化には以下の手法を用いています.
This is a dynamical core module. Used horizontal discretization method is as follows.
Create : | DYNSP 型変数の初期設定 |
Close : | DYNSP 型変数の終了処理 |
PutLine : | DYNSP 型変数に格納されている情報の印字 |
initialized : | DYNSP 型変数が初期設定されているか否か |
EqualAxes : | DYNSP 型変数に格納されている座標値の確認 |
GetAxes : | DYNSP 型変数に格納されている座標値の取得 |
GetCoriolis : | コリオリパラメータの取得 |
GetSpectralCoeff : | スペクトルの添字順序とラプラシアンの係数の取得 |
GetDiffCoeff : | 運動量水平拡散係数と熱, 水水平拡散係数の取得 |
GetPhis : | 地形データ (地表 $ Phi $ ) のスペクトル値を取得 |
GradPi : | 地表面気圧の空間変化の計算 |
VorDiv2UV : | 渦度と発散から東西風速と南北風速を計算 |
UV2VorDiv : | 東西風速と南北風速から渦度と発散を計算 |
Tendency : | 渦度や発散などの予報変数のスペクトル時間変化項の計算 |
TendencyExplicit : | 発散, 温度のスペクトル時間変化項へ地形効果などの追加 (エクスプリシット法で使用) |
Grid2Spectral : | 格子点値からスペクトル値を計算 |
Spectral2Grid : | スペクトル値から格子点値を計算 |
DiffusionVorDiv : | 運動量水平拡散による渦度発散の時間変化の計算 |
DiffusionCorrectTemp : | 運動量水平拡散による摩擦熱補正の計算 |
———— : | ———— |
Create : | Constructor of "DYNSP" |
Close : | Deconstructor of "DYNSP" |
PutLine : | Print information of "DYNSP" |
initialized : | Check initialization of "DYNSP" |
EqualAxes : | Confirm data of axes in "DYNSP" |
GetAxes : | Get data of axes in "DYNSP" |
GetCoriolis : | Get Coriolis parameter |
GetSpectralCoeff : | Get spectral subscript expression, and
Laplacian coefficient |
GetDiffCoeff : | Get coefficient of horizontal momentum diffusion, and coefficient of horizontal thermal and water diffusion. |
GetPhis : | Get spectral geography data (surface $ Phi $ ) |
GradPi : | Calculate spatial variations of surface pressure |
VorDiv2UV : | Calculate zonal and meridional wind from given vorticity and divergence |
UV2VorDiv : | Calculate vorticity and divergence from given zonal and meridional wind. |
Tendency : | Calculate spectral tendency of predictive variables like vorticity and divergence |
TendencyExplicit : | Add surface effect etc. to divergence and temperature tendency (For explicit scheme). |
Grid2Spectral : | Calculate spectral values from given grid values |
Spectral2Grid : | Calculate grid values from given spectral values |
DiffusionVorDiv : | Caluculate vorticity and divergence tendency by horizontal diffusion of momentum |
DiffusionCorrectTemp : | Caluculate frictional thermal correction by
horizontal momentum diffusion |
始めに, DYNSP 型の変数を定義し, Create で初期化します. このモジュールではスペクトル変換による計算を行うための サブルーチン群が用意されています. 上記から選んで使用してください. DYNSP 型の変数の終了処理には Close を用います.
First, initialize "DYNSP" by "Create". This module provides subroutines that calculate by spectral transform. Select from above list. In order to terminate "DYNSP", use "Close".
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
err : | logical, intent(out), optional
|
引数 dyn_sp に設定された内容を空にします. DYNSP 型の 変数に関して再利用する際に使用します. この手続きによって 空となった変数は再度 Create によって初期化して使用します. なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Deconstructor of "DYNSP". Note that if dyn_sp is not initialized by "Create" yet, error is occurred.
Alias for DynSpectralClose
Subroutine : | |||||
dyn_sp : | type(DYNSP), intent(inout) | ||||
nmax : | integer, intent(in)
| ||||
imax : | integer, intent(in)
| ||||
jmax : | integer, intent(in)
| ||||
kmax : | integer, intent(in)
| ||||
PI : | real(DP), intent(in)
| ||||
RPlanet : | real(DP), intent(in)
| ||||
Omega : | real(DP), intent(in)
| ||||
Cp : | real(DP), intent(in)
| ||||
VisOrder : | integer, intent(in)
| ||||
EFoldTime : | real(DP), intent(in)
| ||||
DelTime : | real(DP), intent(in)
| ||||
viscous_effect : | logical, intent(in), optional
| ||||
xy_Phis(0:imax-1, 0:jmax-1) : | real(DP), intent(in), optional
| ||||
openmp_threads : | integer, intent(in), optional
| ||||
wa_module_initialized : | logical, intent(in), optional
| ||||
nmlfile : | character(*), intent(in), optional
| ||||
err : | logical, intent(out), optional
|
引数 dyn_sp にスペクトル法を用いた演算の設定を行います. 他のサブルーチンを使用する前に必ずこのサブルーチンによって DYNSP 型の変数を初期設定してください. なお, 与えられた dyn_sp が既に初期設定されている場合, プログラムはエラーを発生させます.
NAMELIST を利用する場合には引数 nmlfile に NAMELIST ファイル名 を与えてください. NAMELIST 変数群の詳細に関しては NAMELIST#dyn_spectral_nml を参照してください.
Configure the settings for spectral method to dyn_sp. Initialize dyn_sp by this subroutine, before other procedures are used, Note that if dyn_sp is already initialized by this procedure, error is occurred.
In order to use NAMELIST, specify a NAMELIST filename to argument nmlfile. See "NAMELIST#dyn_spectral_nml" for details about a NAMELIST group.
Alias for DynSpectralCreate
Derived Type : | |||||
initialized = .false. : | logical
| ||||
nmax : | integer
| ||||
imax : | integer
| ||||
jmax : | integer
| ||||
kmax : | integer
| ||||
PI : | real(DP)
| ||||
RPlanet : | real(DP)
| ||||
Omega : | real(DP)
| ||||
Cp : | real(DP)
| ||||
VisOrder : | integer
| ||||
EFoldTime : | real(DP)
| ||||
viscous_effect : | logical
| ||||
DelTime : | real(DP)
| ||||
xy_Cori(:,:) =>null() : | real(DP), pointer
| ||||
xy_UVFact(:,:) =>null() : | real(DP), pointer
| ||||
wz_DiffVorDiv(:,:) =>null() : | real(DP), pointer
| ||||
wz_DiffTherm(:,:) =>null() : | real(DP), pointer
| ||||
xy_Phis(:,:) =>null() : | real(DP), pointer
| ||||
openmp_threads = 0 : | integer
|
DYNSP 型の変数を使用する際には必ず Create によって初期設定を 行ってください. 一度 Create で初期化された変数を別途再利用 する際は, Close によって終了処理を行ってください.
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
wz_VorDiff((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
wz_DivDiff((dyn_sp%nmax+1)**2, 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(in)
| ||
xyz_V(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_Temp(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(inout)
| ||
err : | logical, intent(out), optional
|
運動量水平拡散による摩擦熱補正を計算します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Caluculate frictional thermal correction by horizontal momentum diffusion
If dyn_sp is not initialized by "Create" yet, error is occurred.
Alias for DynSpectralDiffusionCorrectTemp
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
wz_Vor((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
wz_Div((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
wz_VorDiff((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
wz_DivDiff((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
運動量水平拡散による渦度発散の時間変化を計算します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Caluculate vorticity and divergence tendency by horizontal diffusion of momentum
If dyn_sp is not initialized by "Create" yet, error is occurred.
Alias for DynSpectralDiffusionVorDiv
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
x_Lon(0:dyn_sp%imax-1) : | real(DP), intent(in), optional
| ||
y_Lat(0:dyn_sp%jmax-1) : | real(DP), intent(in), optional
| ||
err : | logical, intent(out), optional
|
与えられた x_Lon および y_Lat と dyn_sp 内に 保持される緯度経度データが等しいことを確認します. 異なる場合にはエラーを発生させます.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Confirm equality between given x_Lon, y_lat and longitude and latitude data stored in dyn_sp. If the equality is not confirmed, error is occurred.
If dyn_sp is not initialized by "Create" yet, error is occurred.
Alias for DynSpectralEqualAxes
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
x_Lon(0:dyn_sp%imax-1) : | real(DP), intent(out)
| ||
y_Lat(0:dyn_sp%jmax-1) : | real(DP), intent(out)
| ||
x_Lon_Weight(0:dyn_sp%imax-1) : | real(DP), intent(out), optional
| ||
y_Lat_Weight(0:dyn_sp%jmax-1) : | real(DP), intent(out), optional
| ||
err : | logical, intent(out), optional
|
dyn_sp に格納されている座標値を返します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Return data of axes stored in dyn_sp.
If dyn_sp is not initialized by "Create" yet, error is occurred.
Alias for DynSpectralGetAxes
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
xy_Cori(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
コリオリパラメータ $ f\equiv 2\Omega\sin\varphi $ を返します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Return Coriolis parameter $ f\equiv 2\Omega\sin\varphi $ .
If dyn_sp is not initialized by "Create" yet, error is occurred.
Alias for DynSpectralGetCoriolis
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
wz_DiffVorDiv((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
wz_DiffTherm((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
運動量水平拡散係数と熱, 水水平拡散係数を返します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Return coefficient of horizontal momentum diffusion, and coefficient of horizontal thermal and water diffusion.
If dyn_sp is not initialized by "Create" yet, error is occurred.
Alias for DynSpectralGetDiffCoeff
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
w_Phis((dyn_sp%nmax+1)**2) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
地形データ (地表 $ Phi $ ) のスペクトル値を返します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Return spectral geography data (surface $ Phi $ ).
If dyn_sp is not initialized by "Create" yet, error is occurred.
Alias for DynSpectralGetPhis
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
nmo(1:2, 0:dyn_sp%nmax, 0:dyn_sp%nmax) : | integer, intent(out)
| ||
wz_rn((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
スペクトルの添字順序とラプラシアンの係数を返します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Return spectral subscript expression, and Laplacian coefficient.
If dyn_sp is not initialized by "Create" yet, error is occurred.
Alias for DynSpectralGetSpectralCoeff
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
xy_Ps(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) : | real(DP), intent(in)
| ||
xy_GradLonPi(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) : | real(DP), intent(out)
| ||
xy_GradLatPi(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
与えられた地表面気圧 $ P_s $ から, 地表面気圧の空間変化 $ DD{pi}{x} $ および $ DD{pi}{y} $ を返す. なお, $ pi = ln P_s $ である.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
From given surface presure $ P_s $, spatial variations of surface pressure $ DD{pi}{x} $ and $ DD{pi}{y} $ are returned, where $ pi = ln P_s $.
If dyn_sp is not initialized by "Create" yet, error is occurred.
Alias for DynSpectralGradPi
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
xy_Data(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) : | real(DP), intent(in)
| ||
w_Data((dyn_sp%nmax+1)**2) : | real(DP), intent(out)
| ||
math_func : | character(*), intent(in), optional
| ||
err : | logical, intent(out), optional
|
与えられた格子点値からスペクトル値を計算します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculate spectral values from given grid values.
If dyn_sp is not initialized by "Create" yet, error is occurred.
Alias for DynSpectralGrid2Spectral2
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
xyz_Data(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
wz_Data((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
math_func : | character(*), intent(in), optional
| ||
err : | logical, intent(out), optional
|
与えられた格子点値からスペクトル値を計算します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculate spectral values from given grid values.
If dyn_sp is not initialized by "Create" yet, error is occurred.
Alias for DynSpectralGrid2Spectral3
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(in) | ||
unit : | integer, intent(in), optional
| ||
indent : | character(*), intent(in), optional
| ||
err : | logical, intent(out), optional
|
引数 dyn_sp に設定されている情報を印字します. デフォルトではメッセージは標準出力に出力されます. unit に装置番号を指定することで, 出力先を変更することが可能です.
Print information of dyn_sp. By default messages are output to standard output. Unit number for output can be changed by unit argument.
Alias for DynSpectralPutLine
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
w_Data((dyn_sp%nmax+1)**2) : | real(DP), intent(in)
| ||
xy_Data(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) : | real(DP), intent(out)
| ||
math_func : | character(*), intent(in), optional
| ||
err : | logical, intent(out), optional
|
与えられたスペクトル値から格子点値を計算します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculate grid values from given spectral values.
If dyn_sp is not initialized by "Create" yet, error is occurred.
Alias for DynSpectralSpectral2Grid2
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
wz_Data((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_Data(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
math_func : | character(*), intent(in), optional
| ||
err : | logical, intent(out), optional
|
与えられたスペクトル値から格子点値を計算します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculate grid values from given spectral values.
If dyn_sp is not initialized by "Create" yet, error is occurred.
Alias for DynSpectralSpectral2Grid3
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
xyz_UAdv(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_VAdv(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_KE(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_TempUAdv(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_TempVAdv(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_DTempDt(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xy_DPiDt(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) : | real(DP), intent(in)
| ||
xyz_QVapUAdv(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_QVapVAdv(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_DQVapDt(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
wz_DVorDt((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
wz_DDivDt((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
wz_DTempDt((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
w_DPiDt((dyn_sp%nmax+1)**2) : | real(DP), intent(out)
| ||
wz_DQVapDt((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
与えられた格子点値 (移流項, 運動エネルギー項, 温度変化項 $ H $ , 比湿変化項 $ R $ ) から, スペクトルの時間変化項 (渦度, 発散, 温度, 地表面気圧, 比湿) を計算します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculate spectral tendency terms (vorticity, divergence, temperature, surface pressure, specific humidity) from grid points values (advection, kinetic energy, and temperature tendency $ H $ , specific humidity tendency $ R $ ).
If dyn_sp is not initialized by "Create" yet, error is occurred.
Alias for DynSpectralTendency
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
xyz_exWTGPi(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_exHDiv(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
wz_DDivDt((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(inout)
| ||
wz_DTempDt((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(inout)
| ||
err : | logical, intent(out), optional
|
発散 (スペクトル) の時間変化項に $ Phi_s underline{W} Dvect{T} + Dvect{G} pi $ の効果を, 発散 (スペクトル) の時間変化項に $ underline{h} Dvect{D} $ の効果を追加します. 時間積分法にエクスプリシット法を使用している場合に呼び出します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Add effect of $ Phi_s + underline{W} Dvect{T} + Dvect{G} pi $ to divergence (spectral) tendency, and add effect of $ underline{h} Dvect{D} $ to temperature (spectral) tendency. This routine called then explicit scheme is used as time integration scheme.
If dyn_sp is not initialized by "Create" yet, error is occurred.
Alias for DynSpectralTendencyExplicit
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.
Alias for DynSpectralUV2VorDiv
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.
Alias for DynSpectralVorDiv2UV
Function : | |
result : | logical |
dyn_sp : | type(DYNSP), intent(in) |
dyn_sp が初期設定されている場合には .true. が, 初期設定されていない場合には .false. が返ります.
If dyn_sp is initialized, .true. is returned. If dyn_sp is not initialized, .false. is returned.
Alias for DynSpectralInitialized
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
err : | logical, intent(out), optional
|
引数 dyn_sp に設定された内容を空にします. DYNSP 型の 変数に関して再利用する際に使用します. この手続きによって 空となった変数は再度 Create によって初期化して使用します. なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Deconstructor of "DYNSP". Note that if dyn_sp is not initialized by "Create" yet, error is occurred.
subroutine DynSpectralClose( dyn_sp, err ) ! ! 引数 *dyn_sp* に設定された内容を空にします. DYNSP 型の ! 変数に関して再利用する際に使用します. この手続きによって ! 空となった変数は再度 Create によって初期化して使用します. ! なお, 与えられた *dyn_sp* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Deconstructor of "DYNSP". ! Note that if *dyn_sp* 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 use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT implicit none type(DYNSP), intent(inout):: dyn_sp 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 = 'DynSpectralClose' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSP' goto 999 end if !----------------------------------------------------------------- ! "DYNSP" の設定の消去 ! Clear the settings for "DYNSP" !----------------------------------------------------------------- deallocate( dyn_sp % xy_Cori ) deallocate( dyn_sp % xy_UVFact ) deallocate( dyn_sp % wz_DiffVorDiv ) deallocate( dyn_sp % wz_DiffTherm ) deallocate( dyn_sp % xy_Phis ) !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- dyn_sp % initialized = .false. 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralClose
Subroutine : | |||||
dyn_sp : | type(DYNSP), intent(inout) | ||||
nmax : | integer, intent(in)
| ||||
imax : | integer, intent(in)
| ||||
jmax : | integer, intent(in)
| ||||
kmax : | integer, intent(in)
| ||||
PI : | real(DP), intent(in)
| ||||
RPlanet : | real(DP), intent(in)
| ||||
Omega : | real(DP), intent(in)
| ||||
Cp : | real(DP), intent(in)
| ||||
VisOrder : | integer, intent(in)
| ||||
EFoldTime : | real(DP), intent(in)
| ||||
DelTime : | real(DP), intent(in)
| ||||
viscous_effect : | logical, intent(in), optional
| ||||
xy_Phis(0:imax-1, 0:jmax-1) : | real(DP), intent(in), optional
| ||||
openmp_threads : | integer, intent(in), optional
| ||||
wa_module_initialized : | logical, intent(in), optional
| ||||
nmlfile : | character(*), intent(in), optional
| ||||
err : | logical, intent(out), optional
|
引数 dyn_sp にスペクトル法を用いた演算の設定を行います. 他のサブルーチンを使用する前に必ずこのサブルーチンによって DYNSP 型の変数を初期設定してください. なお, 与えられた dyn_sp が既に初期設定されている場合, プログラムはエラーを発生させます.
NAMELIST を利用する場合には引数 nmlfile に NAMELIST ファイル名 を与えてください. NAMELIST 変数群の詳細に関しては NAMELIST#dyn_spectral_nml を参照してください.
Configure the settings for spectral method to dyn_sp. Initialize dyn_sp by this subroutine, before other procedures are used, Note that if dyn_sp is already initialized by this procedure, error is occurred.
In order to use NAMELIST, specify a NAMELIST filename to argument nmlfile. See "NAMELIST#dyn_spectral_nml" for details about a NAMELIST group.
subroutine DynSpectralCreate( dyn_sp, nmax, imax, jmax, kmax, PI, RPlanet, Omega, Cp, VisOrder, EFoldTime, DelTime, viscous_effect, xy_Phis, openmp_threads, wa_module_initialized, nmlfile, err ) ! ! 引数 *dyn_sp* にスペクトル法を用いた演算の設定を行います. ! 他のサブルーチンを使用する前に必ずこのサブルーチンによって ! DYNSP 型の変数を初期設定してください. ! なお, 与えられた *dyn_sp* が既に初期設定されている場合, ! プログラムはエラーを発生させます. ! ! NAMELIST を利用する場合には引数 *nmlfile* に NAMELIST ファイル名 ! を与えてください. NAMELIST 変数群の詳細に関しては ! NAMELIST#dyn_spectral_nml を参照してください. ! ! Configure the settings for spectral method to *dyn_sp*. ! Initialize *dyn_sp* by this subroutine, ! before other procedures are used, ! Note that if *dyn_sp* is already initialized ! by this procedure, error is occurred. ! ! In order to use NAMELIST, specify a NAMELIST filename to ! argument *nmlfile*. See "NAMELIST#dyn_spectral_nml" ! for details about a NAMELIST group. ! use wa_module, only: wa_Initial, xy_Lat use w_module, only: rn use dc_types, only: DP, STRING use dc_present, only: present_and_not_empty, present_and_true use dc_message, only: MessageNotify use dc_string, only: PutLine use dc_trace, only: BeginSub, EndSub, DbgMessage use dc_error, only: DC_ENOFILEREAD, DC_ENOTINIT use dcpam_error, only: StoreError, DC_NOERR, DCPAM_EALREADYINIT, DCPAM_ENEGATIVE implicit none type(DYNSP), intent(inout):: dyn_sp 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 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 logical, intent(in), optional:: viscous_effect ! 粘性効果のスイッチ. ! デフォルトでは *VisOrder* と *EFoldTime* ! とで設定される水平粘性の効果が ! 運動量, 熱, 水に対して働きます. ! この引数に .false. を与えることで ! 水平粘性を無効にします. ! ! Switch of viscous effect. ! By default, horizontal diffusion set ! with *VisOrder* and *EFoldTime* make ! an effect on momentum, heating, and water. ! If .false. is specified to this argument, ! the horizontal diffusion becomes invalid. ! 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 より大きな値を ! 指定しなければ並列計算ルーチンは呼ばれない. 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 :: nmlfile ! NAMELIST ファイルの名称. ! この引数に空文字以外を与えた場合, ! 指定されたファイルから ! NAMELIST 変数群を読み込みます. ! ファイルを読み込めない場合にはエラーを ! 生じます. ! ! NAMELIST 変数群の詳細に関しては ! NAMELIST#dyn_spectral_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_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. integer:: k ! DO ループ用作業変数 ! Work variables for DO loop integer:: stat character(STRING):: cause_c real(DP):: wz_rn ((nmax+1)**2, 0:kmax-1) ! $ -n \times (n+1) $ . ラプラシアンの係数. ! Laplacian coefficient real(DP):: VisCoef ! 超粘性係数. Hyper-viscosity coefficient character(*), parameter:: subname = 'DynSpectralCreate' continue call BeginSub(subname, version) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (dyn_sp % initialized) then stat = DCPAM_EALREADYINIT cause_c = 'DYNSP' goto 999 end if !----------------------------------------------------------------- ! 引数の正当性のチェック ! Validation of arguments !----------------------------------------------------------------- if (nmax < 1) then stat = DCPAM_ENEGATIVE cause_c = 'nmax' goto 999 end if if (imax < 1) then stat = DCPAM_ENEGATIVE cause_c = 'imax' goto 999 end if if (jmax < 1) then stat = DCPAM_ENEGATIVE cause_c = 'jmax' goto 999 end if if (kmax < 1) then stat = DCPAM_ENEGATIVE cause_c = 'kmax' goto 999 end if !----------------------------------------------------------------- ! 波数・格子点の設定 ! Configure wave number and grid point !----------------------------------------------------------------- dyn_sp % nmax = nmax dyn_sp % imax = imax dyn_sp % jmax = jmax dyn_sp % kmax = kmax !----------------------------------------------------------------- ! 物理定数の設定 ! Configure physical constants !----------------------------------------------------------------- dyn_sp % PI = PI dyn_sp % RPlanet = RPlanet dyn_sp % Omega = Omega !!$ dyn_sp % Grav = Grav dyn_sp % Cp = Cp !!$ dyn_sp % RAir = RAir dyn_sp % VisOrder = VisOrder dyn_sp % EFoldTime = EFoldTime dyn_sp % DelTime = DelTime if ( present(viscous_effect) ) then dyn_sp % viscous_effect = viscous_effect else dyn_sp % viscous_effect = .true. end if !!$ dyn_sp % Kappa = RAir / Cp !----------------------------------------------------------------- ! OpenMP の設定 ! Configure OpenMP !----------------------------------------------------------------- if (present(openmp_threads)) then if (openmp_threads > 1) then dyn_sp % openmp_threads = openmp_threads end if end if !----------------------------------------------------------------- ! NAMELIST の読み込み ! Load NAMELIST !----------------------------------------------------------------- if ( present_and_not_empty(nmlfile) ) then call MessageNotify( 'M', subname, 'Loading NAMELIST file "%c" ...', c1=trim(nmlfile) ) call NmlRead ( nmlfile = nmlfile, VisOrder = dyn_sp % VisOrder, EFoldTime = dyn_sp % EFoldTime, viscous_effect = dyn_sp % viscous_effect, openmp_threads = dyn_sp % openmp_threads, 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 !----------------------------------------------------------------- ! spml ライブラリ wa_module の初期設定 ! Initialization of "wa_module" of "spml" library !----------------------------------------------------------------- if ( .not. present_and_true(wa_module_initialized) ) then if (wa_module_initialized_stored) then call MessageNotify('W', subname, 'spml library wa_module can not be initialized more than once.' // ' Object (DYNSP) was not initialized.') stat = DC_ENOTINIT cause_c = 'DYNSP' goto 999 else if (dyn_sp % openmp_threads < 2) then call wa_Initial(nmax, imax, jmax, kmax) else call wa_Initial(nmax, imax, jmax, kmax, openmp_threads) end if wa_module_initialized_stored = .true. end if elseif ( present_and_true(wa_module_initialized) ) then wa_module_initialized_stored = .true. end if !----------------------------------------------------------------- ! コリオリ力, $ u \rightarrow U, v \rightarrow V $ 係数の計算 ! Calculate Coriolis force, coefficients for ! $ u \rightarrow U, v \rightarrow V $ !----------------------------------------------------------------- allocate( dyn_sp % xy_Cori (0:imax-1, 0:jmax-1) ) dyn_sp % xy_Cori = 2.0_DP * Omega * sin(xy_Lat) allocate( dyn_sp % xy_UVFact (0:imax-1, 0:jmax-1) ) dyn_sp % xy_UVFact = cos( xy_Lat ) !----------------------------------------------------------------- ! スペクトル演算用係数の設定 ! Calculate coefficients for spectral calculation !----------------------------------------------------------------- allocate( dyn_sp % wz_DiffVorDiv ((nmax+1)**2, 0:kmax-1) ) allocate( dyn_sp % wz_DiffTherm ((nmax+1)**2, 0:kmax-1) ) do k = 0, kmax-1 wz_rn(:,k) = rn(:,1) enddo ! 粘性係数の計算 (最大波数の e-folding time が EFoldTime となるように) ! Calculate viscosity coefficient if ( dyn_sp % viscous_effect ) then VisCoef = (real((nmax*(nmax+1)), DP) / RPlanet**2 )**(-VisOrder / 2) / EFoldTime else VisCoef = 0.0_DP end if call DbgMessage('VisCoef=<%f>', d=(/VisCoef/)) dyn_sp % wz_DiffTherm = VisCoef * ( (-wz_rn / RPlanet**2)**(VisOrder / 2) ) dyn_sp % wz_DiffVorDiv = dyn_sp % wz_DiffTherm + VisCoef * ( - (2.0_DP / RPlanet**2)**(VisOrder / 2)) !------------------------------------------------------------------- ! 地形データ (地表 $ \Phi $ ) の設定 ! Configure geography data (surface $ \Phi $ ) !------------------------------------------------------------------- allocate( dyn_sp % xy_Phis (0:imax-1, 0:jmax-1) ) if (present(xy_Phis)) then dyn_sp % xy_Phis = xy_Phis else dyn_sp % xy_Phis = 0.0_DP end if !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- dyn_sp % initialized = .true. 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralCreate
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
wz_VorDiff((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
wz_DivDiff((dyn_sp%nmax+1)**2, 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(in)
| ||
xyz_V(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_Temp(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(inout)
| ||
err : | logical, intent(out), optional
|
運動量水平拡散による摩擦熱補正を計算します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Caluculate frictional thermal correction by horizontal momentum diffusion
If dyn_sp is not initialized by "Create" yet, error is occurred.
subroutine DynSpectralDiffusionCorrectTemp( dyn_sp, wz_VorDiff, wz_DivDiff, xyz_U, xyz_V, xyz_Temp, err ) ! ! 運動量水平拡散による摩擦熱補正を計算します. ! ! なお, 与えられた *dyn_sp* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Caluculate frictional thermal correction by horizontal ! momentum diffusion ! ! If *dyn_sp* 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 use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT use wa_module, only: xya_wa implicit none type(DYNSP), intent(inout):: dyn_sp real(DP), intent(in):: wz_VorDiff ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! $ \mathscr{D}(\zeta) $ . ! 運動量水平拡散による渦度変化 (スペクトル). ! Vorticity tendency by ! horizontal momentum diffusion (spectral) real(DP), intent(in):: wz_DivDiff ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! $ \mathscr{D}(D) $ . ! 運動量水平拡散による発散変化 (スペクトル). ! Divergence tendency by ! horizontal momentum diffusion (spectral) real(DP), intent(in):: xyz_U (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ U $ . 東西風速. Zonal wind. real(DP), intent(in):: xyz_V (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ V $ . 南北風速. Meridional wind. real(DP), intent(inout):: xyz_Temp (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ T $ . 温度. Temperature 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):: xyz_UDiff (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! 運動量水平拡散による東西風変化. ! Zonal wind tendency by ! horizontal momentum diffusion real(DP):: xyz_VDiff (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! 運動量水平拡散による南北風変化. ! Meridional wind tendency by ! horizontal momentum diffusion real(DP):: Cp ! $ C_p $ . 大気定圧比熱. Specific heat of air at constant pressure real(DP):: DelTime ! $ \Delta t $ . タイムステップ. Time step integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpectralDiffusionCorrectTemp' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSP' goto 999 end if !----------------------------------------------------------------- ! *dyn_sp* に格納される物理定数, $ \Delta t $ の取り出し ! Fetch physical constant, $ \Delta t $ from *dyn_sp* !----------------------------------------------------------------- Cp = dyn_sp % Cp DelTime = dyn_sp % DelTime !----------------------------------------------------------------- ! 運動量水平拡散による東西風速, 南北風速の時間変化の計算 ! Caluculate zonal and meridional wind tendency by ! horizontal diffusion of momentum !----------------------------------------------------------------- call VorDiv2UV( dyn_sp, xya_wa( wz_VorDiff ), xya_wa( wz_DivDiff ), xyz_UDiff, xyz_VDiff ) ! (out) !----------------------------------------------------------------- ! 運動量水平拡散による摩擦熱補正の計算 ! Caluculate frictional thermal correction by ! horizontal diffusion of momentum !----------------------------------------------------------------- xyz_Temp = xyz_Temp - ( xyz_U * xyz_UDiff + xyz_V * xyz_VDiff ) / Cp * 2.0_DP * DelTime !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralDiffusionCorrectTemp
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
wz_Vor((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
wz_Div((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
wz_VorDiff((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
wz_DivDiff((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
運動量水平拡散による渦度発散の時間変化を計算します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Caluculate vorticity and divergence tendency by horizontal diffusion of momentum
If dyn_sp is not initialized by "Create" yet, error is occurred.
subroutine DynSpectralDiffusionVorDiv( dyn_sp, wz_Vor, wz_Div, wz_VorDiff, wz_DivDiff, err ) ! ! 運動量水平拡散による渦度発散の時間変化を計算します. ! ! なお, 与えられた *dyn_sp* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Caluculate vorticity and divergence tendency by ! horizontal diffusion of momentum ! ! If *dyn_sp* 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 use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT implicit none type(DYNSP), intent(inout):: dyn_sp real(DP), intent(in):: wz_Vor ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! $ \zeta $ . 渦度 (スペクトル). ! Vorticity (spectral) real(DP), intent(in):: wz_Div ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! $ D $ . 発散 (スペクトル). ! Divergence (spectral) real(DP), intent(out):: wz_VorDiff ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! $ \mathscr{D}(\zeta) $ . ! 運動量水平拡散による渦度変化 (スペクトル). ! Vorticity tendency by ! horizontal momentum diffusion (spectral) real(DP), intent(out):: wz_DivDiff ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! $ \mathscr{D}(D) $ . ! 運動量水平拡散による発散変化 (スペクトル). ! Divergence tendency by ! horizontal momentum diffusion (spectral) 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):: wz_DiffVorDiv((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! $ -K_{HD} [(-1)^{N_D/2}\nabla^{N_D}- (2/a^2)^{N_D/2}] $ . ! 運動量水平拡散係数. ! Coefficient of horizontal momentum diffusion integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpectralDiffusionVorDiv' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSP' goto 999 end if !----------------------------------------------------------------- ! *dyn_sp* に格納される拡散係数の取り出し ! Fetch diffusion coefficients from *dyn_sp* !----------------------------------------------------------------- wz_DiffVorDiv = dyn_sp % wz_DiffVorDiv !----------------------------------------------------------------- ! 運動量水平拡散による渦度発散の時間変化の計算 ! Caluculate vorticity and divergence tendency by ! horizontal diffusion of momentum !----------------------------------------------------------------- wz_VorDiff = - wz_Vor * wz_DiffVorDiv wz_DivDiff = - wz_Div * wz_DiffVorDiv !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralDiffusionVorDiv
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
x_Lon(0:dyn_sp%imax-1) : | real(DP), intent(in), optional
| ||
y_Lat(0:dyn_sp%jmax-1) : | real(DP), intent(in), optional
| ||
err : | logical, intent(out), optional
|
与えられた x_Lon および y_Lat と dyn_sp 内に 保持される緯度経度データが等しいことを確認します. 異なる場合にはエラーを発生させます.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Confirm equality between given x_Lon, y_lat and longitude and latitude data stored in dyn_sp. If the equality is not confirmed, error is occurred.
If dyn_sp is not initialized by "Create" yet, error is occurred.
subroutine DynSpectralEqualAxes( dyn_sp, x_Lon, y_Lat, err ) ! ! 与えられた *x_Lon* および *y_Lat* と *dyn_sp* 内に ! 保持される緯度経度データが等しいことを確認します. ! 異なる場合にはエラーを発生させます. ! ! なお, 与えられた *dyn_sp* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Confirm equality between given *x_Lon*, *y_lat* and ! longitude and latitude data stored in *dyn_sp*. ! If the equality is not confirmed, error is occurred. ! ! If *dyn_sp* is not initialized by "Create" yet, ! error is occurred. ! use dc_trace, only: BeginSub, EndSub use dc_types, only: STRING, STDOUT use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT, DCPAM_EAXISMISMATCH use dc_message, only: MessageNotify use dc_string, only: JoinChar, toChar, PutLine implicit none type(DYNSP), intent(inout):: dyn_sp real(DP), intent(in), optional:: x_Lon (0:dyn_sp%imax-1) ! 経度. Longitude real(DP), intent(in), optional:: y_Lat (0:dyn_sp%jmax-1) ! 緯度. Latitude 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):: x_LonStored (0:dyn_sp%imax-1) ! 経度. Longitude real(DP):: y_LatStored (0:dyn_sp%jmax-1) ! 緯度. Latitude logical:: lon_gt_check (0:dyn_sp%imax-1), lon_lt_check (0:dyn_sp%imax-1) logical:: lat_gt_check (0:dyn_sp%jmax-1), lat_lt_check (0:dyn_sp%jmax-1) real(DP), parameter:: check_accuracy = 1.0e-6_DP integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpectralEqualAxes' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSP' goto 999 end if !----------------------------------------------------------------- ! 軸データの比較 ! Compare axes data !----------------------------------------------------------------- call GetAxes ( dyn_sp = dyn_sp, x_Lon = x_LonStored, y_Lat = y_LatStored ) ! (out) if ( present(x_Lon) ) then lon_gt_check = x_Lon > ( x_LonStored - check_accuracy ) lon_lt_check = x_Lon < ( x_LonStored + check_accuracy ) if ( any(.not. lon_gt_check) .or. any(.not. lon_lt_check) ) then !!$ call MessageNotify ( 'W', subname, & !!$ & '%c (argument)=%*f', & !!$ & c1 = 'x_Lon', d = x_Lon, n = (/dyn_sp % imax/) ) !!$ call MessageNotify ( 'W', subname, & !!$ & 'dyn_sp %% %c =%*f', & !!$ & c1 = 'x_Lon', d = x_LonStored, n = (/dyn_sp % imax/) ) stat = DCPAM_EAXISMISMATCH cause_c = 'x_Lon' goto 999 end if end if if ( present(y_Lat) ) then lat_gt_check = y_Lat > ( y_LatStored - check_accuracy ) lat_lt_check = y_Lat < ( y_LatStored + check_accuracy ) if ( any(.not. lat_gt_check) .or. any(.not. lat_lt_check) ) then !!$ call MessageNotify ( 'W', subname, & !!$ & '%c (argument)=%*f', & !!$ & c1 = 'y_Lat', d = y_Lat, n = (/dyn_sp % jmax/) ) !!$ call MessageNotify ( 'W', subname, & !!$ & 'dyn_sp %% %c =%*f', & !!$ & c1 = 'y_Lat', d = y_LatStored, n = (/dyn_sp % jmax /) ) stat = DCPAM_EAXISMISMATCH cause_c = 'y_Lat' goto 999 end if end if !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralEqualAxes
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
x_Lon(0:dyn_sp%imax-1) : | real(DP), intent(out)
| ||
y_Lat(0:dyn_sp%jmax-1) : | real(DP), intent(out)
| ||
x_Lon_Weight(0:dyn_sp%imax-1) : | real(DP), intent(out), optional
| ||
y_Lat_Weight(0:dyn_sp%jmax-1) : | real(DP), intent(out), optional
| ||
err : | logical, intent(out), optional
|
dyn_sp に格納されている座標値を返します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Return data of axes stored in dyn_sp.
If dyn_sp is not initialized by "Create" yet, error is occurred.
subroutine DynSpectralGetAxes( dyn_sp, x_Lon, y_Lat, x_Lon_Weight, y_Lat_Weight, err ) ! ! *dyn_sp* に格納されている座標値を返します. ! ! なお, 与えられた *dyn_sp* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Return data of axes stored in *dyn_sp*. ! ! If *dyn_sp* 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 use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT use wa_module, only: x_Lon_spml => x_Lon, y_Lat_spml => y_Lat, x_Lon_Weight_spml => x_Lon_Weight, y_Lat_Weight_spml => y_Lat_Weight implicit none type(DYNSP), intent(inout):: dyn_sp real(DP), intent(out):: x_Lon (0:dyn_sp%imax-1) ! 経度. Longitude real(DP), intent(out):: y_Lat (0:dyn_sp%jmax-1) ! 緯度. Latitude real(DP), intent(out), optional:: x_Lon_Weight(0:dyn_sp%imax-1) ! 経度座標重み. ! Weight of longitude real(DP), intent(out), optional:: y_Lat_Weight(0:dyn_sp%jmax-1) ! 緯度座標重み. ! Weight of latitude 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 = 'DynSpectralGetAxes' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSP' goto 999 end if !----------------------------------------------------------------- ! 緯度経度座標値の設定 ! Configure data of latitude and longitude !----------------------------------------------------------------- x_Lon = x_Lon_spml y_Lat = y_Lat_spml if ( present( x_Lon_Weight ) ) x_Lon_Weight = x_Lon_Weight_spml if ( present( y_Lat_Weight ) ) y_Lat_Weight = y_Lat_Weight_spml !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralGetAxes
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
xy_Cori(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
コリオリパラメータ $ f\equiv 2\Omega\sin\varphi $ を返します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Return Coriolis parameter $ f\equiv 2\Omega\sin\varphi $ .
If dyn_sp is not initialized by "Create" yet, error is occurred.
subroutine DynSpectralGetCoriolis( dyn_sp, xy_Cori, err ) ! ! コリオリパラメータ $ f\equiv 2\Omega\sin\varphi $ を返します. ! ! なお, 与えられた *dyn_sp* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Return Coriolis parameter $ f\equiv 2\Omega\sin\varphi $ . ! ! If *dyn_sp* 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 use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT implicit none type(DYNSP), intent(inout):: dyn_sp real(DP), intent(out):: xy_Cori (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) ! $ f\equiv 2\Omega\sin\varphi $ . ! コリオリパラメータ. Coriolis parameter 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 = 'DynSpectralGetCoriolis' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSP' goto 999 end if !----------------------------------------------------------------- ! *dyn_sp* に格納されるコリオリパラメータの取り出し ! Fetch Coriolis parameter from *dyn_sp* !----------------------------------------------------------------- xy_Cori = dyn_sp % xy_Cori !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralGetCoriolis
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
wz_DiffVorDiv((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
wz_DiffTherm((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
運動量水平拡散係数と熱, 水水平拡散係数を返します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Return coefficient of horizontal momentum diffusion, and coefficient of horizontal thermal and water diffusion.
If dyn_sp is not initialized by "Create" yet, error is occurred.
subroutine DynSpectralGetDiffCoeff( dyn_sp, wz_DiffVorDiv, wz_DiffTherm, err ) ! ! 運動量水平拡散係数と熱, 水水平拡散係数を返します. ! ! なお, 与えられた *dyn_sp* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Return coefficient of horizontal momentum diffusion, ! and coefficient of horizontal thermal and water diffusion. ! ! If *dyn_sp* 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 use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT implicit none type(DYNSP), intent(inout):: dyn_sp real(DP), intent(out):: wz_DiffVorDiv ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! $ -K_{HD} [(-1)^{N_D/2}\nabla^{N_D}- (2/a^2)^{N_D/2}] $ . ! 運動量水平拡散係数. ! Coefficient of horizontal momentum diffusion real(DP), intent(out):: wz_DiffTherm ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! $ -(-1)^{N_D/2}K_{HD}\nabla^{N_D} $ . ! 熱, 水水平拡散係数. ! Coefficient of horizontal thermal and water diffusion 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 = 'DynSpectralGetDiffCoeff' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSP' goto 999 end if !----------------------------------------------------------------- ! *dyn_sp* に格納される拡散係数の取り出し ! Fetch diffusion coefficients from *dyn_sp* !----------------------------------------------------------------- wz_DiffVorDiv = dyn_sp % wz_DiffVorDiv wz_DiffTherm = dyn_sp % wz_DiffTherm !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralGetDiffCoeff
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
w_Phis((dyn_sp%nmax+1)**2) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
地形データ (地表 $ Phi $ ) のスペクトル値を返します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Return spectral geography data (surface $ Phi $ ).
If dyn_sp is not initialized by "Create" yet, error is occurred.
subroutine DynSpectralGetPhis( dyn_sp, w_Phis, err ) ! ! 地形データ (地表 $ \Phi $ ) のスペクトル値を返します. ! ! なお, 与えられた *dyn_sp* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Return spectral geography data (surface $ \Phi $ ). ! ! If *dyn_sp* 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 use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT use wa_module, only: w_xy implicit none type(DYNSP), intent(inout):: dyn_sp real(DP), intent(out):: w_Phis ((dyn_sp%nmax+1)**2) ! $ \Phi_s $ . 地表ジオポテンシャル. ! Surface geo-potential 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 = 'DynSpectralGetPhis' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSP' goto 999 end if !------------------------------------------------------------------- ! 地形データ (地表 $ \Phi $ ) のスペクトル値の計算 ! Calculate spectral geography data (surface $ \Phi $ ) !------------------------------------------------------------------- w_Phis = w_xy( dyn_sp % xy_Phis ) !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralGetPhis
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
nmo(1:2, 0:dyn_sp%nmax, 0:dyn_sp%nmax) : | integer, intent(out)
| ||
wz_rn((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
スペクトルの添字順序とラプラシアンの係数を返します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Return spectral subscript expression, and Laplacian coefficient.
If dyn_sp is not initialized by "Create" yet, error is occurred.
subroutine DynSpectralGetSpectralCoeff( dyn_sp, nmo, wz_rn, err ) ! ! スペクトルの添字順序とラプラシアンの係数を返します. ! ! なお, 与えられた *dyn_sp* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Return spectral subscript expression, and ! Laplacian coefficient. ! ! If *dyn_sp* 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 use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT use wa_module, only: l_nm use w_module, only: rn implicit none type(DYNSP), intent(inout):: dyn_sp integer, intent(out):: nmo (1:2, 0:dyn_sp%nmax, 0:dyn_sp%nmax) ! スペクトルの添字順番. ! Spectral subscript expression real(DP), intent(out):: wz_rn ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! $ -n \times (n+1) $ . ラプラシアンの係数. ! Laplacian coefficient 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:: kmax ! 鉛直層数. ! Number of vertical level integer:: nmax ! 最大全波数. ! Maximum truncated wavenumber integer:: mmax ! 最大東西波数. ! Maximum truncated zonal wavenumber integer:: lmax ! 最大南北波数. ! Maximum truncated meridional wavenumber integer:: k, l, m ! DO ループ用作業変数 ! Work variables for DO loop integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpectralGetSpectralCoeff' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSP' goto 999 end if !----------------------------------------------------------------- ! スペクトルの添字順番の取り出し ! Fetch spectral subscript expression !----------------------------------------------------------------- nmax = dyn_sp % nmax mmax = dyn_sp % nmax lmax = dyn_sp % nmax nmo = 0 do l = 0, lmax do m = 0, min(mmax, nmax-l) nmo(1,m,l) = l_nm(m+l, m) nmo(2,m,l) = l_nm(m+l, -m) enddo enddo !----------------------------------------------------------------- ! ラプラシアンの係数 $ -n \times (n+1) $ の取り出し ! Fetch Laplacian coefficient $ -n \times (n+1) $ !----------------------------------------------------------------- kmax = dyn_sp % kmax do k = 0, kmax-1 wz_rn(:,k) = rn(:,1) enddo !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralGetSpectralCoeff
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
xy_Ps(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) : | real(DP), intent(in)
| ||
xy_GradLonPi(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) : | real(DP), intent(out)
| ||
xy_GradLatPi(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
与えられた地表面気圧 $ P_s $ から, 地表面気圧の空間変化 $ DD{pi}{x} $ および $ DD{pi}{y} $ を返す. なお, $ pi = ln P_s $ である.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
From given surface presure $ P_s $, spatial variations of surface pressure $ DD{pi}{x} $ and $ DD{pi}{y} $ are returned, where $ pi = ln P_s $.
If dyn_sp is not initialized by "Create" yet, error is occurred.
subroutine DynSpectralGradPi( dyn_sp, xy_Ps, xy_GradLonPi, xy_GradLatPi, err ) ! ! 与えられた地表面気圧 $ P_s $ から, ! 地表面気圧の空間変化 $ \DD{\pi}{x} $ および $ \DD{\pi}{y} $ を返す. ! なお, $ \pi = \ln P_s $ である. ! ! なお, 与えられた *dyn_sp* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! From given surface presure $ P_s $, ! spatial variations of surface pressure ! $ \DD{\pi}{x} $ and $ \DD{\pi}{y} $ are returned, where ! $ \pi = \ln P_s $. ! ! If *dyn_sp* 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 use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT use wa_module, only: xy_GradLon_w, xy_GradLat_w, w_xy implicit none type(DYNSP), intent(inout):: dyn_sp real(DP), intent(in):: xy_Ps (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) ! $ P_s $ 地表面気圧. ! Surface pressure real(DP), intent(out):: xy_GradLonPi (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) ! $ \DD{\pi}{x} $ real(DP), intent(out):: xy_GradLatPi (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) ! $ \DD{\pi}{y} $ 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. real(DP):: xy_Pi (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) ! $ \pi = \ln P_s $ real(DP):: w_Pi ((dyn_sp%nmax+1)**2) ! $ \pi $ スペクトル integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpectralGradPi' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSP' goto 999 end if !----------------------------------------------------------------- ! 地表面気圧の空間変化の計算 ! Calculate spatial surface pressure tendency !----------------------------------------------------------------- xy_Pi = log( xy_Ps ) w_Pi = w_xy( xy_Pi ) xy_GradLonPi = xy_GradLon_w( w_Pi ) xy_GradLatPi = xy_GradLat_w( w_Pi ) !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralGradPi
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
xy_Data(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) : | real(DP), intent(in)
| ||
w_Data((dyn_sp%nmax+1)**2) : | real(DP), intent(out)
| ||
math_func : | character(*), intent(in), optional
| ||
err : | logical, intent(out), optional
|
与えられた格子点値からスペクトル値を計算します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculate spectral values from given grid values.
If dyn_sp is not initialized by "Create" yet, error is occurred.
subroutine DynSpectralGrid2Spectral2( dyn_sp, xy_Data, w_Data, math_func, err ) ! ! 与えられた格子点値からスペクトル値を計算します. ! ! なお, 与えられた *dyn_sp* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Calculate spectral values from given grid values. ! ! If *dyn_sp* 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 use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT, DCPAM_EBADMATHFUNC use wa_module, only: w_xy implicit none type(DYNSP), intent(inout):: dyn_sp real(DP), intent(in):: xy_Data (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) ! 格子点データ. Grid data real(DP), intent(out):: w_Data ((dyn_sp%nmax+1)**2) ! スペクトルデータ. Spectral Data character(*), intent(in), optional:: math_func ! 数学関数. ! 以下の数学関数を選択可能. ! ! Mathematical function. ! Available functions are as follows. ! ! * "log" ! * "exp" ! 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 = 'DynSpectralGrid2Spectral2' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSP' goto 999 end if !------------------------------------------------------------------- ! スペクトル正変換 ! Forward spectral transformation !------------------------------------------------------------------- if ( present(math_func) ) then select case ( trim(math_func) ) case ('log') w_Data = w_xy( log( xy_Data ) ) case ('exp') w_Data = exp( w_xy( xy_Data ) ) case default stat = DCPAM_EBADMATHFUNC cause_c = trim(math_func) goto 999 end select else w_Data = w_xy( xy_Data ) end if !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralGrid2Spectral2
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
xyz_Data(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
wz_Data((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
math_func : | character(*), intent(in), optional
| ||
err : | logical, intent(out), optional
|
与えられた格子点値からスペクトル値を計算します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculate spectral values from given grid values.
If dyn_sp is not initialized by "Create" yet, error is occurred.
subroutine DynSpectralGrid2Spectral3( dyn_sp, xyz_Data, wz_Data, math_func, err ) ! ! 与えられた格子点値からスペクトル値を計算します. ! ! なお, 与えられた *dyn_sp* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Calculate spectral values from given grid values. ! ! If *dyn_sp* 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 use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT, DCPAM_EBADMATHFUNC use wa_module, only: wa_xya implicit none type(DYNSP), intent(inout):: dyn_sp real(DP), intent(in):: xyz_Data (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! 格子点データ. Grid data real(DP), intent(out):: wz_Data ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! スペクトルデータ. Spectral Data character(*), intent(in), optional:: math_func ! 数学関数. ! 以下の数学関数を選択可能. ! ! Mathematical function. ! Available functions are as follows. ! ! * "log" ! * "exp" ! 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 = 'DynSpectralGrid2Spectral3' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSP' goto 999 end if !------------------------------------------------------------------- ! スペクトル正変換 ! Forward spectral transformation !------------------------------------------------------------------- if ( present(math_func) ) then select case ( trim(math_func) ) case ('log') wz_Data = wa_xya( log( xyz_Data ) ) case ('exp') wz_Data = exp( wa_xya( xyz_Data ) ) case default stat = DCPAM_EBADMATHFUNC cause_c = trim(math_func) goto 999 end select else wz_Data = wa_xya( xyz_Data ) end if !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralGrid2Spectral3
Function : | |
result : | logical |
dyn_sp : | type(DYNSP), intent(in) |
dyn_sp が初期設定されている場合には .true. が, 初期設定されていない場合には .false. が返ります.
If dyn_sp is initialized, .true. is returned. If dyn_sp is not initialized, .false. is returned.
logical function DynSpectralInitialized( dyn_sp ) result(result) ! ! *dyn_sp* が初期設定されている場合には .true. が, ! 初期設定されていない場合には .false. が返ります. ! ! If *dyn_sp* is initialized, .true. is returned. ! If *dyn_sp* is not initialized, .false. is returned. ! implicit none type(DYNSP), intent(in):: dyn_sp continue result = dyn_sp % initialized end function DynSpectralInitialized
Subroutine : | |||
nmlfile : | character(*), intent(in)
| ||
VisOrder : | integer, intent(inout)
| ||
EFoldTime : | real(DP), intent(inout)
| ||
viscous_effect : | logical, intent(inout)
| ||
openmp_threads : | integer, intent(inout)
| ||
err : | logical, intent(out), optional
|
NAMELIST ファイル nmlfile から値を入力するための 内部サブルーチンです. Create 内で呼び出されることを 想定しています.
値が NAMELIST ファイル内で指定されていない場合には, 入力された値がそのまま返ります.
なお, nmlfile に空文字が与えられた場合, または 与えられた nmlfile を読み込むことができない場合, プログラムはエラーを発生させます.
This is an internal subroutine to input values from NAMELIST file nmlfile. This subroutine is expected to be called by "Create".
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_nml .
subroutine DynSpectralNmlRead( nmlfile, VisOrder, EFoldTime, viscous_effect, openmp_threads, err ) ! ! NAMELIST ファイル *nmlfile* から値を入力するための ! 内部サブルーチンです. Create 内で呼び出されることを ! 想定しています. ! ! 値が NAMELIST ファイル内で指定されていない場合には, ! 入力された値がそのまま返ります. ! ! なお, *nmlfile* に空文字が与えられた場合, または ! 与えられた *nmlfile* を読み込むことができない場合, ! プログラムはエラーを発生させます. ! ! This is an internal subroutine to input values from ! NAMELIST file *nmlfile*. This subroutine is expected to be ! called by "Create". ! ! 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_types, only: DP, STRING, TOKEN, STDOUT use dc_string, only: PutLine use dc_iounit, only: FileOpen use dc_message, only: MessageNotify use dc_present, only: present_and_true use dc_error, only: StoreError, DC_NOERR, DC_ENOFILEREAD implicit none character(*), intent(in):: nmlfile ! NAMELIST ファイルの名称. ! NAMELIST file name integer, intent(inout):: VisOrder ! 超粘性の次数. Order of hyper-viscosity real(DP), intent(inout):: EFoldTime ! 最大波数に対する e-folding time. E-folding time for maximum wavenumber logical, intent(inout):: viscous_effect ! 粘性効果のスイッチ. ! デフォルトでは *VisOrder* と *EFoldTime* ! とで設定される水平粘性の効果が ! 運動量, 熱, 水に対して働きます. ! この引数に .false. を与えることで ! 水平粘性を無効にします. ! ! Switch of viscous effect. ! By default, horizontal diffusion set ! with *VisOrder* and *EFoldTime* make ! an effect on momentum, heating, and water. ! If .false. is specified to this argument, ! the horizontal diffusion becomes invalid. ! integer, intent(inout):: openmp_threads ! OPENMP での最大スレッド数. ! openmp_threads に 1 より大きな値を指定すれば ! ISPACK[http://www.gfd-dennou.org/library/ispack/] ! の球面調和函数変換 OPENMP 並列計算 ! ルーチンが用いられる. 並列計算を実行するには, ! 実行時に環境変数 OMP_NUM_THREADS ! を openmp_threads 以下の数字に設定する ! 等のシステムに応じた準備が必要となる. ! ! openmp_threads に 1 より大きな値を ! 指定しなければ並列計算ルーチンは呼ばれない. 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_nml/ VisOrder, EFoldTime, viscous_effect, openmp_threads ! dyn_spectral モジュール用 ! NAMELIST 変数群名. ! ! dyn_spectral#Create を使用する際に, ! オプショナル引数 *nmlfile* へ NAMELIST ! ファイル名を指定することで, そのファイルから ! この NAMELIST 変数群を読み込みます. ! ! NAMELIST group name for ! "dyn_spectral" module. ! ! If a NAMELIST filename is specified to ! an optional argument *nmlfile* ! when "dyn_spectral#Create" 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(*), parameter:: subname = 'DynSpectralNmlRead' continue call BeginSub( subname ) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 文字型引数を NAMELIST 変数群へ代入 ! Substitute character arguments to NAMELIST group !----------------------------------------------------------------- !!$ param_c = param_c_ !---------------------------------------------------------------- ! 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 !----------------------------------------------------------------- read( unit = unit_nml, nml = dyn_spectral_nml, iostat = iostat_nml ) ! (out) if ( iostat_nml == 0 ) then call MessageNotify( 'M', subname, 'NAMELIST group "%c" is loaded from "%c".', c1='dyn_spectral_nml', c2=trim(nmlfile) ) write(STDOUT, nml = dyn_spectral_nml) else call MessageNotify( 'W', subname, 'NAMELIST group "%c" is not found in "%c" (iostat=%d).', c1='dyn_spectral_nml', c2=trim(nmlfile), i=(/iostat_nml/) ) end if close( unit_nml ) !----------------------------------------------------------------- ! NAMELIST 変数群を文字型引数へ代入 ! Substitute NAMELIST group to character arguments !----------------------------------------------------------------- !!$ param_c_ = param_c !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError( stat, subname, err, cause_c ) call EndSub( subname ) end subroutine DynSpectralNmlRead
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(in) | ||
unit : | integer, intent(in), optional
| ||
indent : | character(*), intent(in), optional
| ||
err : | logical, intent(out), optional
|
引数 dyn_sp に設定されている情報を印字します. デフォルトではメッセージは標準出力に出力されます. unit に装置番号を指定することで, 出力先を変更することが可能です.
Print information of dyn_sp. By default messages are output to standard output. Unit number for output can be changed by unit argument.
subroutine DynSpectralPutLine( dyn_sp, unit, indent, err ) ! ! 引数 *dyn_sp* に設定されている情報を印字します. ! デフォルトではメッセージは標準出力に出力されます. ! *unit* に装置番号を指定することで, 出力先を変更することが可能です. ! ! Print information of *dyn_sp*. ! By default messages are output to standard output. ! Unit number for output can be changed by *unit* argument. ! use dc_types, only: STRING, STDOUT use dc_date, only: toChar use dc_trace, only: BeginSub, EndSub use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT use dc_string, only: Printf, PutLine implicit none type(DYNSP), intent(in):: dyn_sp 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. integer:: stat character(STRING):: cause_c integer:: out_unit integer:: indent_len character(STRING):: indent_str integer:: xy_minpos(1:2), xy_maxpos(1:2) integer:: xy_lbound(1:2), xy_ubound(1:2) integer:: wz_minpos(1:2), wz_maxpos(1:2) integer:: wz_lbound(1:2), wz_ubound(1:2) character(len = *), parameter:: subname = 'DynSpectralPutLine' 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 !----------------------------------------------------------------- ! "DYNSP" の設定の印字 ! Print the settings for "DYNSP" !----------------------------------------------------------------- if (dyn_sp % initialized) then call Printf(out_unit, indent_str(1:indent_len) // '#<DYNSP:: @initialized=%y @nmax=%d @imax=%d @jmax=%d @kmax=%d', i=(/dyn_sp % nmax, dyn_sp % imax, dyn_sp % jmax, dyn_sp % kmax/), l=(/dyn_sp % initialized/)) call Printf(out_unit, indent_str(1:indent_len) // ' @PI=%f @RPlanet=%f @Omega=%f @Cp=%f', d=(/dyn_sp % PI, dyn_sp % RPlanet, dyn_sp % Omega, dyn_sp % Cp/)) !!$ call Printf(out_unit, & !!$ & indent_str(1:indent_len) // & !!$ & ' @Grav=%f @RAir=%f @Kappa=%f ', & !!$ & d=(/dyn_sp % Grav, dyn_sp % RAir, dyn_sp % Kappa/)) call Printf(out_unit, indent_str(1:indent_len) // ' @VisOrder=%d @EFoldTime=%f @DelTime=%f', i=(/dyn_sp % VisOrder/), d=(/dyn_sp % EFoldTime, dyn_sp % DelTime/)) call Printf(out_unit, indent_str(1:indent_len) // ' @openmp_threads=%d', i=(/dyn_sp % openmp_threads/)) xy_lbound = lbound(dyn_sp % xy_Cori) xy_ubound = ubound(dyn_sp % xy_Cori) xy_minpos = minloc(dyn_sp % xy_Cori) xy_maxpos = maxloc(dyn_sp % xy_Cori) xy_minpos(1) = xy_minpos(1) + xy_lbound(1) - 1 xy_minpos(2) = xy_minpos(2) + xy_lbound(2) - 1 xy_maxpos(1) = xy_maxpos(1) + xy_lbound(1) - 1 xy_maxpos(2) = xy_maxpos(2) + xy_lbound(2) - 1 call Printf(out_unit, indent_str(1:indent_len) // ' @xy_Cori(%d:%d,%d:%d)=' // '[min:%f, max:%f]', i=(/xy_lbound(1), xy_ubound(1), xy_lbound(2), xy_ubound(2)/), d=(/dyn_sp % xy_Cori(xy_minpos(1), xy_minpos(2)), dyn_sp % xy_Cori(xy_maxpos(1), xy_maxpos(2))/)) xy_lbound = lbound(dyn_sp % xy_UVFact) xy_ubound = ubound(dyn_sp % xy_UVFact) xy_minpos = minloc(dyn_sp % xy_UVFact) xy_maxpos = maxloc(dyn_sp % xy_UVFact) xy_minpos(1) = xy_minpos(1) + xy_lbound(1) - 1 xy_minpos(2) = xy_minpos(2) + xy_lbound(2) - 1 xy_maxpos(1) = xy_maxpos(1) + xy_lbound(1) - 1 xy_maxpos(2) = xy_maxpos(2) + xy_lbound(2) - 1 call Printf(out_unit, indent_str(1:indent_len) // ' @xy_UVFact(%d:%d,%d:%d)=' // '[min:%f, max:%f]', i=(/xy_lbound(1), xy_ubound(1), xy_lbound(2), xy_ubound(2)/), d=(/dyn_sp % xy_UVFact(xy_minpos(1), xy_minpos(2)), dyn_sp % xy_UVFact(xy_maxpos(1), xy_maxpos(2))/)) wz_lbound = lbound(dyn_sp % wz_DiffVorDiv) wz_ubound = ubound(dyn_sp % wz_DiffVorDiv) wz_minpos = minloc(dyn_sp % wz_DiffVorDiv) wz_maxpos = maxloc(dyn_sp % wz_DiffVorDiv) wz_minpos(1) = wz_minpos(1) + wz_lbound(1) - 1 wz_minpos(2) = wz_minpos(2) + wz_lbound(2) - 1 wz_maxpos(1) = wz_maxpos(1) + wz_lbound(1) - 1 wz_maxpos(2) = wz_maxpos(2) + wz_lbound(2) - 1 call Printf(out_unit, indent_str(1:indent_len) // ' @wz_DiffVorDiv(%d:%d,%d:%d)=' // '[min:%f, max:%f]', i=(/wz_lbound(1), wz_ubound(1), wz_lbound(2), wz_ubound(2)/), d=(/dyn_sp % wz_DiffVorDiv(wz_lbound(1), wz_lbound(2)), dyn_sp % wz_DiffVorDiv(wz_ubound(1), wz_ubound(2))/)) wz_lbound = lbound(dyn_sp % wz_DiffTherm) wz_ubound = ubound(dyn_sp % wz_DiffTherm) wz_minpos = minloc(dyn_sp % wz_DiffTherm) wz_maxpos = maxloc(dyn_sp % wz_DiffTherm) wz_minpos(1) = wz_minpos(1) + wz_lbound(1) - 1 wz_minpos(2) = wz_minpos(2) + wz_lbound(2) - 1 wz_maxpos(1) = wz_maxpos(1) + wz_lbound(1) - 1 wz_maxpos(2) = wz_maxpos(2) + wz_lbound(2) - 1 call Printf(out_unit, indent_str(1:indent_len) // ' @wz_DiffTherm(%d:%d,%d:%d)=' // '[min:%f, max:%f]', i=(/wz_lbound(1), wz_ubound(1), wz_lbound(2), wz_ubound(2)/), d=(/dyn_sp % wz_DiffTherm(wz_lbound(1), wz_lbound(2)), dyn_sp % wz_DiffTherm(wz_ubound(1), wz_ubound(2))/)) call Printf(out_unit, indent_str(1:indent_len) // '>') else call Printf(out_unit, indent_str(1:indent_len) // '#<DYNSP:: @initialized=%y>', l=(/dyn_sp % initialized/)) end if !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralPutLine
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
w_Data((dyn_sp%nmax+1)**2) : | real(DP), intent(in)
| ||
xy_Data(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) : | real(DP), intent(out)
| ||
math_func : | character(*), intent(in), optional
| ||
err : | logical, intent(out), optional
|
与えられたスペクトル値から格子点値を計算します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculate grid values from given spectral values.
If dyn_sp is not initialized by "Create" yet, error is occurred.
subroutine DynSpectralSpectral2Grid2( dyn_sp, w_Data, xy_Data, math_func, err ) ! ! 与えられたスペクトル値から格子点値を計算します. ! ! なお, 与えられた *dyn_sp* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Calculate grid values from given spectral values. ! ! If *dyn_sp* 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 use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT, DCPAM_EBADMATHFUNC use wa_module, only: xy_w implicit none type(DYNSP), intent(inout):: dyn_sp real(DP), intent(in):: w_Data ((dyn_sp%nmax+1)**2) ! スペクトルデータ. Spectral Data real(DP), intent(out):: xy_Data (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) ! 格子点データ. Grid data character(*), intent(in), optional:: math_func ! 数学関数. ! 以下の数学関数を選択可能. ! ! Mathematical function. ! Available functions are as follows. ! ! * "log" ! * "exp" ! logical, intent(out), optional:: err ! 例外処理用フラグ. ! デフォルトでは, この手続き内でエラーが ! 生じた場合, プログラムは強制終了します. ! 引数 *err* が与えられる場合, ! プログラムは強制終了せず, 代わりに ! *err* に .true. が代入されます. ! ! Exception handling flag. ! By default, when error occur in ! this procedure, the program aborts. ! If this *err* argument is given, ! .true. is substituted to *err* and ! the program does not abort. !----------------------------------- ! 作業変数 ! Work variables integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpectralSpectral2Grid2' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSP' goto 999 end if !------------------------------------------------------------------- ! スペクトル逆変換 ! Backward spectral transformation !------------------------------------------------------------------- if ( present(math_func) ) then select case ( trim(math_func) ) case ('log') xy_Data = xy_w( log( w_Data ) ) case ('exp') xy_Data = exp( xy_w( w_Data ) ) case default stat = DCPAM_EBADMATHFUNC cause_c = trim(math_func) goto 999 end select else xy_Data = xy_w( w_Data ) end if !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralSpectral2Grid2
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
wz_Data((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_Data(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
math_func : | character(*), intent(in), optional
| ||
err : | logical, intent(out), optional
|
与えられたスペクトル値から格子点値を計算します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculate grid values from given spectral values.
If dyn_sp is not initialized by "Create" yet, error is occurred.
subroutine DynSpectralSpectral2Grid3( dyn_sp, wz_Data, xyz_Data, math_func, err ) ! ! 与えられたスペクトル値から格子点値を計算します. ! ! なお, 与えられた *dyn_sp* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Calculate grid values from given spectral values. ! ! If *dyn_sp* 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 use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT, DCPAM_EBADMATHFUNC use wa_module, only: xya_wa implicit none type(DYNSP), intent(inout):: dyn_sp real(DP), intent(in):: wz_Data ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! スペクトルデータ. Spectral Data real(DP), intent(out):: xyz_Data (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! 格子点データ. Grid data character(*), intent(in), optional:: math_func ! 数学関数. ! 以下の数学関数を選択可能. ! ! Mathematical function. ! Available functions are as follows. ! ! * "log" ! * "exp" ! logical, intent(out), optional:: err ! 例外処理用フラグ. ! デフォルトでは, この手続き内でエラーが ! 生じた場合, プログラムは強制終了します. ! 引数 *err* が与えられる場合, ! プログラムは強制終了せず, 代わりに ! *err* に .true. が代入されます. ! ! Exception handling flag. ! By default, when error occur in ! this procedure, the program aborts. ! If this *err* argument is given, ! .true. is substituted to *err* and ! the program does not abort. !----------------------------------- ! 作業変数 ! Work variables integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpectralSpectral2Grid3' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSP' goto 999 end if !------------------------------------------------------------------- ! スペクトル逆変換 ! Backward spectral transformation !------------------------------------------------------------------- if ( present(math_func) ) then select case ( trim(math_func) ) case ('log') xyz_Data = xya_wa( log( wz_Data ) ) case ('exp') xyz_Data = exp( xya_wa( wz_Data ) ) case default stat = DCPAM_EBADMATHFUNC cause_c = trim(math_func) goto 999 end select else xyz_Data = xya_wa( wz_Data ) end if !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralSpectral2Grid3
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
xyz_UAdv(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_VAdv(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_KE(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_TempUAdv(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_TempVAdv(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_DTempDt(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xy_DPiDt(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) : | real(DP), intent(in)
| ||
xyz_QVapUAdv(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_QVapVAdv(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_DQVapDt(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
wz_DVorDt((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
wz_DDivDt((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
wz_DTempDt((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
w_DPiDt((dyn_sp%nmax+1)**2) : | real(DP), intent(out)
| ||
wz_DQVapDt((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(out)
| ||
err : | logical, intent(out), optional
|
与えられた格子点値 (移流項, 運動エネルギー項, 温度変化項 $ H $ , 比湿変化項 $ R $ ) から, スペクトルの時間変化項 (渦度, 発散, 温度, 地表面気圧, 比湿) を計算します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Calculate spectral tendency terms (vorticity, divergence, temperature, surface pressure, specific humidity) from grid points values (advection, kinetic energy, and temperature tendency $ H $ , specific humidity tendency $ R $ ).
If dyn_sp is not initialized by "Create" yet, error is occurred.
subroutine DynSpectralTendency( dyn_sp, xyz_UAdv, xyz_VAdv, xyz_KE, xyz_TempUAdv, xyz_TempVAdv, xyz_DTempDt, xy_DPiDt, xyz_QVapUAdv, xyz_QVapVAdv, xyz_DQVapDt, wz_DVorDt, wz_DDivDt, wz_DTempDt, w_DPiDt, wz_DQVapDt, err ) ! ! 与えられた格子点値 (移流項, 運動エネルギー項, ! 温度変化項 $ H $ , 比湿変化項 $ R $ ) から, ! スペクトルの時間変化項 (渦度, 発散, 温度, 地表面気圧, 比湿) ! を計算します. ! ! なお, 与えられた *dyn_sp* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Calculate spectral tendency terms (vorticity, divergence, ! temperature, surface pressure, specific humidity) from ! grid points values (advection, kinetic energy, and ! temperature tendency $ H $ , specific humidity tendency $ R $ ). ! ! If *dyn_sp* 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 use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT use wa_module, only: wa_Div_xya_xya, wa_Lapla_wa, wa_xya, w_xy implicit none type(DYNSP), intent(inout):: dyn_sp 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. real(DP), intent(in):: xyz_UAdv (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ UA $ . 東西運動量移流項. ! Zonal advection of momentum real(DP), intent(in):: xyz_VAdv (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ VA $ . 南北運動量移流項. ! Meridional advection of momentum real(DP), intent(in):: xyz_KE (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ E $ . 運動エネルギー項. ! Kinematic energy real(DP), intent(in):: xyz_TempUAdv (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ UT $ . 温度東西移流項. ! Zonal advection of temperature real(DP), intent(in):: xyz_TempVAdv (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ VT $ . 温度南北移流項. ! Meridional advection of temperature real(DP), intent(in):: xyz_DTempDt (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ H $ . 温度時間変化項. ! Temperature tendency real(DP), intent(in):: xy_DPiDt (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) ! $ Z $ . 地表面気圧時間変化項. ! Surface pressure tendency real(DP), intent(in):: xyz_QVapUAdv (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ Uq $ . 比湿東西移流項. ! Zonal advection of specific humidity real(DP), intent(in):: xyz_QVapVAdv (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ Vq $ . 比湿南北移流項. ! Meridional advection of specific humidity real(DP), intent(in):: xyz_DQVapDt (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ R $ . 比湿時間変化項. ! Specific humidity tendency real(DP), intent(out):: wz_DVorDt ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! $ \DD{\zeta}{t} $ . 渦度変化 (スペクトル). ! Vorticity tendency (spectral) real(DP), intent(out):: wz_DDivDt ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! $ \DD{D}{t} $ . 発散変化 (スペクトル). ! Divergence tendency (spectral) real(DP), intent(out):: wz_DTempDt ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! $ \DD{T}{t} $ . 温度変化 (スペクトル). ! Temperature tendency (spectral) real(DP), intent(out):: w_DPiDt ((dyn_sp%nmax+1)**2) ! $ \DD{Ps}{t} $ . 地表面気圧変化 (スペクトル). ! Surface pressure tendency (spectral) real(DP), intent(out):: wz_DQVapDt ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! $ \DD{q}{t} $ . 比湿変化 (スペクトル). ! Specific humidity tendency (spectral) real(DP):: RPlanet ! $ a $ . 惑星半径. Radius of planet integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpectralTendency' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSP' goto 999 end if !----------------------------------------------------------------- ! *dyn_sp* に格納される物理定数の取り出し ! Fetch physical constants from *dyn_sp* !----------------------------------------------------------------- RPlanet = dyn_sp % RPlanet !------------------------------------------------------------------- ! 渦度の時間変化 (スペクトル) の計算 ! Calculate vorticity tendency (spectral) !------------------------------------------------------------------- wz_DVorDt = wa_Div_xya_xya( xyz_VAdv, - xyz_UAdv ) / RPlanet !------------------------------------------------------------------- ! 発散の時間変化 (スペクトル) の計算 ! Calculate divergence tendency (spectral) !------------------------------------------------------------------- wz_DDivDt = wa_Div_xya_xya( xyz_UAdv, xyz_VAdv ) / RPlanet - wa_Lapla_wa( wa_xya(xyz_KE) ) / RPlanet**2 !------------------------------------------------------------------- ! 温度の時間変化 (スペクトル) の計算 ! Calculate temperature tendency (spectral) !------------------------------------------------------------------- wz_DTempDt = - wa_Div_xya_xya( xyz_TempUAdv, xyz_TempVAdv ) / RPlanet + wa_xya( xyz_DTempDt ) !------------------------------------------------------------------- ! 地表面気圧の時間変化 (スペクトル) の計算 ! Calculate surface pressure tendency (spectral) !------------------------------------------------------------------- w_DPiDt = w_xy( xy_DPiDt ) !------------------------------------------------------------------- ! 比湿の時間変化 (スペクトル) の計算 ! Calculate specific humidity tendency (spectral) !------------------------------------------------------------------- wz_DQVapDt = - wa_Div_xya_xya( xyz_QVapUAdv, xyz_QVapVAdv ) / RPlanet + wa_xya( xyz_DQVapDt ) !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralTendency
Subroutine : | |||
dyn_sp : | type(DYNSP), intent(inout) | ||
xyz_exWTGPi(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
xyz_exHDiv(0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) : | real(DP), intent(in)
| ||
wz_DDivDt((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(inout)
| ||
wz_DTempDt((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) : | real(DP), intent(inout)
| ||
err : | logical, intent(out), optional
|
発散 (スペクトル) の時間変化項に $ Phi_s underline{W} Dvect{T} + Dvect{G} pi $ の効果を, 発散 (スペクトル) の時間変化項に $ underline{h} Dvect{D} $ の効果を追加します. 時間積分法にエクスプリシット法を使用している場合に呼び出します.
なお, 与えられた dyn_sp が Create によって初期設定 されていない場合, プログラムはエラーを発生させます.
Add effect of $ Phi_s + underline{W} Dvect{T} + Dvect{G} pi $ to divergence (spectral) tendency, and add effect of $ underline{h} Dvect{D} $ to temperature (spectral) tendency. This routine called then explicit scheme is used as time integration scheme.
If dyn_sp is not initialized by "Create" yet, error is occurred.
subroutine DynSpectralTendencyExplicit( dyn_sp, xyz_exWTGPi, xyz_exHDiv, wz_DDivDt, wz_DTempDt, err ) ! ! 発散 (スペクトル) の時間変化項に ! $ \Phi_s \underline{W} \Dvect{T} + \Dvect{G} \pi $ の効果を, ! 発散 (スペクトル) の時間変化項に $ \underline{h} \Dvect{D} $ ! の効果を追加します. ! 時間積分法にエクスプリシット法を使用している場合に呼び出します. ! ! なお, 与えられた *dyn_sp* が Create によって初期設定 ! されていない場合, プログラムはエラーを発生させます. ! ! Add effect of $ \Phi_s + \underline{W} \Dvect{T} + \Dvect{G} \pi $ to ! divergence (spectral) tendency, and ! add effect of $ \underline{h} \Dvect{D} $ to ! temperature (spectral) tendency. ! This routine called then explicit scheme is used as ! time integration scheme. ! ! If *dyn_sp* 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 use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT use wa_module, only: wa_Lapla_wa, wa_xya, w_Lapla_w, w_xy implicit none type(DYNSP), intent(inout):: dyn_sp real(DP), intent(in):: xyz_exWTGPi (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ \underline{W} \Dvect{T} + \Dvect{G} \pi $ . real(DP), intent(in):: xyz_exHDiv (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ \underline{h} \Dvect{D} $ . real(DP), intent(inout):: wz_DDivDt ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! $ \DD{D}{t} $ . 発散変化 (スペクトル). ! Divergence tendency (spectral) real(DP), intent(inout):: wz_DTempDt ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! $ \DD{T}{t} $ . 温度変化 (スペクトル). ! Temperature tendency (spectral) 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:: k ! DO ループ用作業変数 ! Work variables for DO loop integer:: kmax ! 鉛直層数. ! Number of vertical level real(DP):: RPlanet ! $ a $ . 惑星半径. Radius of planet real(DP):: xy_Phis (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1) ! $ \Phi_s $ . 地表ジオポテンシャル. ! Surface geo-potential integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpectralTendencyExplicit' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSP' goto 999 end if !----------------------------------------------------------------- ! *dyn_sp* に格納される格子点の取り出し ! Fetch grid point from *dyn_sp* !----------------------------------------------------------------- kmax = dyn_sp % kmax !----------------------------------------------------------------- ! *dyn_sp* に格納される物理定数の取り出し ! Fetch physical constants from *dyn_sp* !----------------------------------------------------------------- RPlanet = dyn_sp % RPlanet !----------------------------------------------------------------- ! *dyn_sp* に格納される地表ジオポテンシャルの取り出し ! Fetch surface geo-potential from *dyn_sp* !----------------------------------------------------------------- xy_Phis = dyn_sp % xy_Phis !------------------------------------------------------------------- ! 発散の時間変化 (スペクトル) の計算 ! Calculate divergence tendency (spectral) !------------------------------------------------------------------- wz_DDivDt = wz_DDivDt - wa_Lapla_wa( wa_xya( xyz_exWTGPi ) ) / RPlanet**2 do k = 0, kmax - 1 wz_DDivDt(:,k) = wz_DDivDt(:,k) - w_Lapla_w( w_xy( xy_Phis ) ) / RPlanet**2 end do !------------------------------------------------------------------- ! 温度の時間変化 (スペクトル) の計算 ! Calculate temperature tendency (spectral) !------------------------------------------------------------------- wz_DTempDt = wz_DTempDt - wa_xya( xyz_exHDiv ) !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralTendencyExplicit
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.
subroutine DynSpectralUV2VorDiv( dyn_sp, xyz_U, xyz_V, xyz_Vor, xyz_Div, err ) ! ! 与えられた東西風速と南北風速から渦度と発散を計算します. ! ! なお, 与えられた *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. ! use dc_trace, only: BeginSub, EndSub use dc_types, only: STRING, STDOUT use dc_string, only: PutLine use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT use wa_module, only: xya_wa, wa_Div_xya_xya implicit none type(DYNSP), intent(inout):: dyn_sp real(DP), intent(in):: xyz_U (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ U $ . 東西風速. ! Zonal wind real(DP), intent(in):: xyz_V (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ V $ . 南北風速. ! Meridional wind real(DP), intent(out):: xyz_Vor (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ \zeta $ . 渦度. ! Vorticity real(DP), intent(out):: xyz_Div (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ D $ . 発散. ! Divergence 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. real(DP):: RPlanet ! 惑星半径 integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpectralUV2VorDiv' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSP' goto 999 end if !----------------------------------------------------------------- ! *dyn_sp* に格納される物理定数の取り出し ! Fetch physical constants from *dyn_sp* !----------------------------------------------------------------- RPlanet = dyn_sp % RPlanet !----------------------------------------------------------------- ! 渦度 $ \zeta $ と発散 $ D $ の計算 ! Caluculate vorticity $ \zeta $ and divergence $ D $ !----------------------------------------------------------------- xyz_Vor = xya_wa( wa_Div_xya_xya( xyz_V , - xyz_U ) / RPlanet ) xyz_Div = xya_wa( wa_Div_xya_xya( xyz_U , xyz_V ) / RPlanet ) !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralUV2VorDiv
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.
subroutine DynSpectralVorDiv2UV( dyn_sp, xyz_Vor, xyz_Div, xyz_U, xyz_V, err ) ! ! 与えられた渦度と発散から東西風速と南北風速を計算します. ! ! なお, 与えられた *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. ! use dc_trace, only: BeginSub, EndSub use dc_types, only: STRING, STDOUT use dc_string, only: PutLine use dcpam_error, only: StoreError, DC_NOERR, DCPAM_ENOTINIT use wa_module, only: wa_xya, wa_LaplaInv_wa, xya_GradLon_wa, xya_GradLat_wa implicit none type(DYNSP), intent(inout):: dyn_sp real(DP), intent(in):: xyz_Vor (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ \zeta $ . 渦度. ! Vorticity real(DP), intent(in):: xyz_Div (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ D $ . 発散. ! Divergence real(DP), intent(out):: xyz_U (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ U $ . 東西風速. ! Zonal wind real(DP), intent(out):: xyz_V (0:dyn_sp%imax-1, 0:dyn_sp%jmax-1, 0:dyn_sp%kmax-1) ! $ V $ . 南北風速. ! Meridional wind 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. real(DP):: RPlanet ! $ a $ . 惑星半径. Radius of planet real(DP):: wz_Psi ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! 流線関数 $ \psi $ スペクトル real(DP):: wz_Chi ((dyn_sp%nmax+1)**2, 0:dyn_sp%kmax-1) ! ポテンシャル $ \chi $ スペクトル integer:: stat character(STRING):: cause_c character(*), parameter:: subname = 'DynSpectralVorDiv2UV' continue call BeginSub(subname) stat = DC_NOERR cause_c = '' !----------------------------------------------------------------- ! 初期設定のチェック ! Check initialization !----------------------------------------------------------------- if (.not. dyn_sp % initialized) then stat = DCPAM_ENOTINIT cause_c = 'DYNSP' goto 999 end if !----------------------------------------------------------------- ! *dyn_sp* に格納される物理定数の取り出し ! Fetch physical constants from *dyn_sp* !----------------------------------------------------------------- RPlanet = dyn_sp % RPlanet !----------------------------------------------------------------- ! 東西風速 $ U $ と南北風速 $ V $ の計算 ! Caluculate zonal wind $ U $ and meridional wind $ V $ !----------------------------------------------------------------- wz_Psi = wa_LaplaInv_wa( wa_xya( xyz_Vor ) ) * RPlanet**2 wz_Chi = wa_LaplaInv_wa( wa_xya( xyz_Div ) ) * RPlanet**2 xyz_U = ( xya_GradLon_wa( wz_Chi ) - xya_GradLat_wa( wz_Psi ) ) / RPlanet xyz_V = ( xya_GradLon_wa( wz_Psi ) + xya_GradLat_wa( wz_Chi ) ) / RPlanet !----------------------------------------------------------------- ! 終了処理, 例外処理 ! Termination and Exception handling !----------------------------------------------------------------- 999 continue call StoreError(stat, subname, err, cause_c) call EndSub(subname) end subroutine DynSpectralVorDiv2UV
Subroutine : | |||
nmlfile : | character(*), intent(in)
| ||
VisOrder : | integer, intent(inout)
| ||
EFoldTime : | real(DP), intent(inout)
| ||
viscous_effect : | logical, intent(inout)
| ||
openmp_threads : | integer, intent(inout)
| ||
err : | logical, intent(out), optional
|
NAMELIST ファイル nmlfile から値を入力するための 内部サブルーチンです. Create 内で呼び出されることを 想定しています.
値が NAMELIST ファイル内で指定されていない場合には, 入力された値がそのまま返ります.
なお, nmlfile に空文字が与えられた場合, または 与えられた nmlfile を読み込むことができない場合, プログラムはエラーを発生させます.
This is an internal subroutine to input values from NAMELIST file nmlfile. This subroutine is expected to be called by "Create".
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_nml .
Alias for DynSpectralNmlRead
Constant : | |
version = ’$Name: dcpam4-20080609-1 $’ // ’$Id: dyn_spectral.f90,v 1.19 2007-10-10 19:43:20 morikawa Exp $’ : | character(*), parameter |