| Class | wq_module_sjpack |
| In: |
libsrc/wq_module_sjpack/wq_module_sjpack.f90
|
| Authors: | Shin-ichi Takehiro, Youhei SASAKI |
| Version: | $Id: wq_module_sjpack.f90 619 2013-10-02 13:27:43Z takepiro $ |
| Copyright&License: | See COPYRIGHT |
spml/wq_module_sjpack モジュールは球内での流体運動をスペクトル法に よって数値計算するための Fortran90 関数を提供するものである.
水平方向に球面調和函数変換および動径方向に Matsushima and Marcus (1994) で提唱された多項式を用いた スペクトル計算のためのさまざまな関数を提供する.
内部で wa_module, aq_module を用いている. 最下部では球面調和変換のエンジンとして ISPACK の Fortran77 サブルーチンを用いている.
| wq_ : | スペクトルデータ(球面調和函数・チェビシェフ変換) |
| nmr_ : | 水平スペクトルデータ(全波数 n, 帯状波数各成分, 動径) |
| nr_ : | 水平スペクトルデータ(全波数 n, 動径) |
| xyr_ : | 3 次元格子点データ(経度・緯度・動径) |
| wr_ : | 水平スペクトル, 動径格子点データ |
| _wq : | スペクトルデータ |
| _xyr : | 3 次元格子点データ |
| _xyr_xyr : | 2 つの3 次元格子点データ, … |
| wq_Initial : | スペクトル変換の格子点数, 波数, 領域の大きさの設定 |
| x_Lon, y_Lat, r_Rad : | 格子点座標(緯度, 経度, 動径座標)を 格納した1 次元配列 |
| x_Lon_Weight, y_Lat_Weight, r_Rad_Weight : | 重み座標を格納した 1 次元配列 |
| xyr_Lon, xyr_Lat, xyr_Rad : | 格子点データの経度・緯度・動径座標(X,Y,Z) (格子点データ型 3 次元配列) |
| xyr_wq, wq_xyr : | スペクトルデータと 3 次元格子データの間の変換 (球面調和函数, チェビシェフ変換) |
| xyr_wr, wr_xyr : | 3 次元格子データと水平スペクトル・動径格子データとの 間の変換 (球面調和函数) |
| wr_wq, wq_wr : | スペクトルデータと水平スペクトル・動径格子データとの 間の変換 (チェビシェフ変換) |
| w_xy, xy_w : | スペクトルデータと 2 次元水平格子データの 間の変換(球面調和函数変換) |
| l_nm, nm_l : | スペクトルデータの格納位置と全波数・帯状波数の変換 |
| wq_RadDRad_wq : | スペクトルデータに動径微分 r∂/∂r を作用させる |
| wr_DivRad_wq : | スペクトルデータに発散型動径微分 1/r^2 ∂/∂r r^2 = ∂/∂r + 2/r を作用させる |
| wr_DivRad_xyr : | 格子点データに発散型動径微分 1/r^2 ∂/∂r r^2 = ∂/∂r + 2/r を作用させる |
| wr_RotDRad_wq : | スペクトルデータに回転型動径微分 1/r ∂/∂rr = ∂/∂r + 1/r を作用させる |
| wr_RotDRad_wr : | スペクトルデータに回転型動径微分 1/r ∂/∂rr = ∂/∂r + 1/r を作用させる |
| wq_RotDRad_wr : | スペクトルデータに回転型動径微分 1/r ∂/∂rr = ∂/∂r + 1/r を作用させる |
| wq_Lapla_wq : | スペクトルデータにラプラシアンを作用させる |
| xyr_GradLon_wq : | スペクトルデータに勾配型経度微分 1/rcosφ・∂/∂λを作用させる |
| xyr_GradLat_wq : | スペクトルデータに勾配型緯度微分 1/r・∂/∂φを作用させる |
| wr_DivLon_xyr : | 格子データに発散型経度微分 1/rcosφ・∂/∂λを作用させる |
| wr_DivLat_xyr : | 格子データに発散型緯度微分 1/rcosφ・∂(g cosφ)/∂φを作用させる |
| wr_Div_xyr_xyr_xyr : | ベクトル成分である 3 つの格子データに 発散を作用させる |
| xyr_Div_xyr_xyr_xyr : | ベクトル成分である 3 つの格子データに 発散を作用させる |
| xyr_RotLon_wq_wq : | ベクトル場の回転の経度成分を計算する |
| xyr_RotLat_wq_wq : | ベクトル場の回転の緯度成分を計算する |
| wr_RotRad_xyr_xyr : | ベクトル場の回転の動径成分を計算する |
| wq_KxRGrad_wq : | スペクトルデータに経度微分 k×r・▽ = ∂/∂λを作用させる |
| xyr_KGrad_wq : | スペクトルデータに軸方向微分 k・▽ = cosφ/r ∂/∂φ + sinφ∂/∂r を作用させる |
| wq_L2_wq : | スペクトルデータに L2 演算子 = -水平ラプラシアンを 作用させる |
| wq_L2Inv_wq : | スペクトルデータに L2 演算子の逆 = -逆水平ラプラシアン を作用させる |
| wq_QOperator_wq : | スペクトルデータに演算子 Q=(k・▽-1/2(L2 k・▽+ k・▽L2)) を作用させる |
| wr_RadRot_xyr_xyr : | ベクトル v の渦度と動径ベクトル r の内積 r・(▽×v) を 計算する |
| wr_RadRotRot_xyr_xyr_xyr : | ベクトルの v の r・(▽×▽×v) を計算する |
| wq_RadRotRot_xyr_xyr_xyr : | ベクトルの v の r・(▽×▽×v) を計算する |
| wq_Potential2Vector : | トロイダルポロイダルポテンシャルから ベクトル場を計算する |
| wq_Potential2Rotation : | トロイダルポロイダルポテンシャルで表される 非発散ベクトル場の回転の各成分を計算する |
| nmr_ToroidalEnergySpectrum_wq, nr_ToroidalEnergySpectrum_wq : | トロイダルポテンシャルからエネルギーの球面調和函数各成分を計算する |
| nmr_PoloidalEnergySpectrum_wq, nr_PoloidalEnergySpectrum_wq : | ポロイダルポテンシャルからエネルギーの球面調和函数各成分を計算する |
| wq_BoundaryTau, wr_BoundaryGrid, wq_Boundary : | ディリクレ, ノイマン境界条件を適用する(タウ法, 選点法) |
| wq_TorBoundaryTau, wr_TorBoundaryGrid, wq_TorBoundary : | 速度トロイダルポテンシャルの境界条件を適用する(タウ法,選点法) |
| wq_LaplaPol2Pol_wq, wq_LaplaPol2PolTau_wq : | 速度ポロイダルポテンシャルΦを▽^2Φから求める (入出力がそれぞれ格子点およびスペクトル係数) |
| wq_TorMagBoundaryTau, wr_TorMagBoundaryGrid, wq_TorMagBoundary : | 磁場トロイダルポテンシャルの境界条件を適用する(タウ法, 選点法) |
| wq_PolMagBoundaryTau, wr_PolMagBoundaryGrid, wq_PolMagBoundary : | 磁場トロイダルポテンシャル境界の境界条件を適用する(タウ法, 選点法) |
| IntLonLatRad_xyr, AvrLonLatRad_xyr : | 3 次元格子点データの 全領域積分および平均 |
| r_IntLonLat_xyr, r_AvrLonLat_xyr : | 3 次元格子点データの 緯度経度(水平・球面)積分および平均 |
| y_IntLonRad_xyr, y_AvrLonRad_xyr : | 3 次元格子点データの 緯度動径積分および平均 |
| r_IntLatRad_xyr, r_AvrLatRad_xyr : | 3 次元格子点データの 経度動径(子午面)積分および平均 |
| yr_IntLon_xyr, yr_AvrLon_xyr : | 3 次元格子点データの 経度方向積分および平均 |
| xr_IntLat_xyr, xr_AvrLat_xyr : | 3 次元格子点データの 緯度方向積分および平均 |
| xr_IntRad_xyr, xr_AvrRad_xyr : | 3 次元格子点データの 動径方向積分および平均 |
| IntLonLat_xy, AvrLonLat_xy : | 2 次元格子点データの水平(球面)積分および平均 |
| IntLonRad_xr, AvrLonRad_xr : | 2 次元(XZ)格子点データの経度動径積分 および平均 |
| IntLatRad_yr, AvrLatRad_yr : | 2 次元(YZ)格子点データの緯度動径(子午面) 積分および平均 |
| y_IntLon_xy, y_AvrLon_xy : | 水平 2 次元(球面)格子点データの経度方向 積分および平均 |
| x_IntLat_xy, x_AvrLat_xy : | 水平2 次元(球面)格子点データの緯度方向積分 および平均 |
| r_IntLon_xr, r_AvrLon_xr : | 2 次元(XZ)格子点データの経度方向積分および 平均 |
| x_IntRad_xr, x_AvrRad_xr : | 2 次元(XZ)格子点データの動径方向積分および 平均 |
| r_IntLat_yr, r_AvrLat_yr : | 2 次元(YZ)格子点データの緯度方向積分および 平均 |
| y_IntRad_yr, y_AvrRad_yr : | 2 次元(YZ)格子点データの動径方向積分および 平均 |
| IntLon_x, AvrLon_x : | 1 次元(X)格子点データの経度方向積分および平均 |
| IntLat_y, AvrLat_y : | 1 次元(Y)格子点データの緯度方向積分および平均 |
| IntRad_r, AvrRad_r : | 1 次元(Z)格子点データの動径方向積分および平均 |
| Function : | |||
| AvrLatRad_yr : | real(8)
| ||
| yr : | real(8), dimension(1:jm,km), intent(in)
|
緯度動径(子午面)積分
2 次元(YR)格子点データの緯度動径(子午面)平均
2 次元データ f(φ,r) に対して
∫f(φ,r) r^2cosφ dφdr /(2(r[o]^3-r[i]^3)/3)
を計算する.
function AvrLatRad_yr(yr) ! 緯度動径(子午面)積分
!
! 2 次元(YR)格子点データの緯度動径(子午面)平均
!
! 2 次元データ f(φ,r) に対して
!
! ∫f(φ,r) r^2cosφ dφdr /(2(r[o]^3-r[i]^3)/3)
!
! を計算する.
!
real(8), dimension(1:jm,km), intent(in) :: yr
!(in) 2 次元緯度動径(子午面)格子点データ
real(8) :: AvrLatRad_yr
!(out) 平均値
AvrLatRad_yr = IntLatRad_yr(yr)/(sum(y_Lat_Weight)*sum(r_Rad_Weight))
end function AvrLatRad_yr
| Function : | |||
| AvrLat_y : | real(8)
| ||
| y_data(1:jm) : | real(8), intent(in)
|
1 次元(Y)格子点データの緯度(Y)方向平均(1 層用).
実際には格子点データ各点毎に y_Y_Weight をかけた総和を計算し, y_Y_Weight の総和で割ることで平均している.
Original external subprogram is wa_module_sjpack#AvrLat_y
| Function : | |||
| AvrLonLatRad_xyr : | real(8)
| ||
| xyr : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
緯度経度動径(全球)積分
3 次元格子点データの緯度経度動径(全球)積分
3 次元データ f(λ,φ,r) に対して
∫f(λ,φ,r) r^2cosφ dλdφdr /(4π(r[o]^3-r[i]^3)/3)
を計算する.
function AvrLonLatRad_xyr(xyr) ! 緯度経度動径(全球)積分
!
! 3 次元格子点データの緯度経度動径(全球)積分
!
! 3 次元データ f(λ,φ,r) に対して
!
! ∫f(λ,φ,r) r^2cosφ dλdφdr /(4π(r[o]^3-r[i]^3)/3)
!
! を計算する.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr
!(in) 3 次元経度緯度動径格子点データ
real(8) :: AvrLonLatRad_xyr
!(out) 全球平均値
AvrLonLatRad_xyr = IntLonLatRad_xyr(xyr) /(sum(x_Lon_Weight)*sum(y_Lat_Weight) * sum(r_Rad_Weight))
end function AvrLonLatRad_xyr
| Function : | |||
| AvrLonLat_xy : | real(8)
| ||
| xy_data(0:im-1,1:jm) : | real(8), intent(in)
|
2 次元緯度経度格子点データの全領域平均(1 層用).
実際には格子点データ各点毎に x_X_Weight, y_Y_Weight をかけた 総和を計算し, x_X_Weight*y_Y_Weight の総和で割ることで平均している.
Original external subprogram is wa_module_sjpack#AvrLonLat_xy
| Function : | |||
| AvrLonRad_xr : | real(8)
| ||
| xr : | real(8), dimension(0:im-1,km), intent(in)
|
経度動径(緯度円)積分
2 次元(XR)格子点データの経度動径平均
2 次元データ f(λ,r) に対して
∫f(λ,r) r^2dλdr /(2π(r[o]^3-r[i]^3)/3)
を計算する.
function AvrLonRad_xr(xr) ! 経度動径(緯度円)積分
!
! 2 次元(XR)格子点データの経度動径平均
!
! 2 次元データ f(λ,r) に対して
!
! ∫f(λ,r) r^2dλdr /(2π(r[o]^3-r[i]^3)/3)
!
! を計算する.
!
real(8), dimension(0:im-1,km), intent(in) :: xr
! (in) 2 次元格子点データ
real(8) :: AvrLonRad_xr
! 積分値
AvrLonRad_xr = IntLonRad_xr(xr)/(sum(x_Lon_Weight)*sum(r_Rad_Weight))
end function AvrLonRad_xr
| Function : | |||
| AvrLon_x : | real(8)
| ||
| x_data(0:im-1) : | real(8), intent(in)
|
1 次元(X)格子点データの経度(X)方向平均(1 層用).
実際には格子点データ各点毎に x_X_Weight をかけた総和を計算し, x_X_Weight の総和で割ることで平均している.
Original external subprogram is wa_module_sjpack#AvrLon_x
| Function : | |||
| AvrRad_r : | real(8)
| ||
| z : | real(8), dimension(km), intent(in)
|
1 次元(Z)格子点データの動径方向域平均.
1 次元データ f(r) に対して ∫f(r) r^2dr /((r[o]^3-r[i]^3)/3) を 計算する.
function AvrRad_r(z)
!
! 1 次元(Z)格子点データの動径方向域平均.
!
! 1 次元データ f(r) に対して ∫f(r) r^2dr /((r[o]^3-r[i]^3)/3) を
! 計算する.
!
real(8), dimension(km), intent(in) :: z
!(in) 1 次元動径格子点データ
real(8) :: AvrRad_r
!(out) 平均値
AvrRad_r = IntRad_r(z)/sum(r_Rad_Weight)
end function AvrRad_r
| Function : | |||
| IntLatRad_yr : | real(8)
| ||
| yr : | real(8), dimension(1:jm,km), intent(in)
|
2 次元(YR)格子点データの緯度動径積分(子午面)および平均
2 次元データ f(φ,r) に対して ∫f(φ,r) r^2cosφ dφdr を計算する.
function IntLatRad_yr(yr)
!
! 2 次元(YR)格子点データの緯度動径積分(子午面)および平均
!
! 2 次元データ f(φ,r) に対して ∫f(φ,r) r^2cosφ dφdr を計算する.
!
real(8), dimension(1:jm,km), intent(in) :: yr
!(in) 2 次元緯度動径(子午面)格子点データ
real(8) :: IntLatRad_yr
!(out) 積分値
integer :: j, k
IntLatRad_yr = 0.0D0
do k=1,km
do j=1,jm
IntLatRad_yr = IntLatRad_yr + yr(j,k) * y_Lat_Weight(j) * r_Rad_Weight(k)
enddo
enddo
end function IntLatRad_yr
| Function : | |||
| IntLat_y : | real(8)
| ||
| y_data(1:jm) : | real(8), intent(in)
|
1 次元緯度(Y)格子点データの Y 方向積分(1 層用).
実際には格子点データ各点毎に y_Y_Weight をかけた総和を計算している.
Original external subprogram is wa_module_sjpack#IntLat_y
| Function : | |||
| IntLonLatRad_xyr : | real(8)
| ||
| xyr : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
緯度経度動径(全球)積分
3 次元格子点データの緯度経度動径(全球)積分
3 次元データ f(λ,φ,r) に対して
∫f(λ,φ,r) r^2cosφ dλdφdr
を計算する.
function IntLonLatRad_xyr(xyr) ! 緯度経度動径(全球)積分
!
! 3 次元格子点データの緯度経度動径(全球)積分
!
! 3 次元データ f(λ,φ,r) に対して
!
! ∫f(λ,φ,r) r^2cosφ dλdφdr
!
! を計算する.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr
!(in) 3 次元経度緯度動径格子点データ
real(8) :: IntLonLatRad_xyr
!(out) 全球積分値
integer :: i, j, k
IntLonLatRad_xyr = 0.0D0
do k=1,km
do j=1,jm
do i=0,im-1
IntLonLatRad_xyr = IntLonLatRad_xyr + xyr(i,j,k) * x_Lon_Weight(i) * y_Lat_Weight(j) * r_Rad_Weight(k)
enddo
enddo
enddo
end function IntLonLatRad_xyr
| Function : | |||
| IntLonLat_xy : | real(8)
| ||
| xy_data(0:im-1,1:jm) : | real(8), intent(in)
|
2 次元緯度経度格子点データの全領域積分(1 層用).
実際には格子点データ各点毎に x_X_Weight, y_Y_Weight をかけた 総和を計算している.
Original external subprogram is wa_module_sjpack#IntLonLat_xy
| Function : | |||
| IntLonRad_xr : | real(8)
| ||
| xr : | real(8), dimension(0:im-1,km), intent(in)
|
経度動径(緯度円)積分
2 次元(XR)格子点データの経度動径積分
2 次元データ f(λ,r) に対して∫f(λ,r) r^2dλdr を計算する.
function IntLonRad_xr(xr) ! 経度動径(緯度円)積分
!
! 2 次元(XR)格子点データの経度動径積分
!
! 2 次元データ f(λ,r) に対して∫f(λ,r) r^2dλdr を計算する.
!
real(8), dimension(0:im-1,km), intent(in) :: xr
!(in) 2 次元緯度動径格子点データ
real(8) :: IntLonRad_xr
!(out) 積分値
integer :: i, k
IntLonRad_xr = 0.0D0
do k=1,km
do i=0,im-1
IntLonRad_xr = IntLonRad_xr + xr(i,k) * x_Lon_Weight(i) * r_Rad_Weight(k)
enddo
enddo
end function IntLonRad_xr
| Function : | |||
| IntLon_x : | real(8)
| ||
| x_data(0:im-1) : | real(8), intent(in)
|
1 次元経度(X)格子点データの X 方向積分(1 層用).
実際には格子点データ各点毎に x_X_Weight をかけた総和を計算している.
Original external subprogram is wa_module_sjpack#IntLon_x
| Function : | |||
| IntRad_r : | real(8)
| ||
| z : | real(8), dimension(km), intent(in)
|
動径積分
1 次元(Z)格子点データの動径方向域積分.
1 次元データ f(r) に対して ∫f(r) r^2dr を計算する.
function IntRad_r(z) ! 動径積分
!
! 1 次元(Z)格子点データの動径方向域積分.
!
! 1 次元データ f(r) に対して ∫f(r) r^2dr を計算する.
!
real(8), dimension(km), intent(in) :: z
!(in) 1 次元動径格子点データ
real(8) :: IntRad_r
!(out) 積分値
integer :: k
IntRad_r = 0.0d0
do k=1,km
IntRad_r = IntRad_r + z(k) * r_Rad_Weight(k)
enddo
end function IntRad_r
| Function : | |||
| l_nm_array00 : | integer
| ||
| n : | integer, intent(in)
| ||
| m : | integer, intent(in)
|
全波数(n)と東西波数(m)からそのスペクトルデータの格納位置を返す.
引数 n,m がともに整数値の場合, 整数値を返す.
Original external subprogram is wa_module_sjpack#l_nm
| Function : | |||
| l_nm_array01(size(marray)) : | integer
| ||
| n : | integer, intent(in)
| ||
| marray(:) : | integer, intent(in)
|
スペクトルデータの格納位置
全波数(n)と東西波数(m)からそのスペクトルデータの格納位置を返す.
第 1 引数 n が整数, 第 2 引数 marray が整数 1 次元配列の場合, marray と同じ大きさの 1 次元整数配列を返す.
Original external subprogram is wa_module_sjpack#l_nm
| Function : | |||
| l_nm_array10(size(narray)) : | integer
| ||
| narray(:) : | integer, intent(in)
| ||
| m : | integer, intent(in)
|
全波数(n)と東西波数(m)からそのスペクトルデータの格納位置を返す.
第 1 引数 narray が整数 1 次元配列, 第 2 引数 m が整数の場合, narray と同じ大きさの 1 次元整数配列を返す.
Original external subprogram is wa_module_sjpack#l_nm
| Function : | |||
| l_nm_array11(size(narray)) : | integer
| ||
| narray(:) : | integer, intent(in)
| ||
| marray(:) : | integer, intent(in)
|
全波数(n)と東西波数(m)からそのスペクトルデータの格納位置を返す.
第 1,2 引数 narray, marray がともに整数 1 次元配列の場合, narray, marray と同じ大きさの 1 次元整数配列を返す. narray, marray は同じ大きさでなければならない.
Original external subprogram is wa_module_sjpack#l_nm
| Function : | |||
| nm_l_int(2) : | integer
| ||
| l : | integer, intent(in)
|
スペクトルデータの格納位置(l)から全波数(n)と東西波数(m)を返す.
引数 l が整数値の場合, 対応する全波数と帯状波数を 長さ 2 の 1 次元整数値を返す. nm_l(1) が全波数, nm_l(2) が帯状波数である.
Original external subprogram is wa_module_sjpack#nm_l
| Function : | |||
| nm_l_array(size(larray),2) : | integer
| ||
| larray(:) : | integer, intent(in)
|
スペクトルデータの格納位置(l)から全波数(n)と東西波数(m)を返す.
引数 larray が整数 1 次元配列の場合, larray に対応する n, m を格納した 2 次元整数配列を返す. nm_l_array(:,1) が全波数, nm_l_array(:,2) が帯状波数である.
Original external subprogram is wa_module_sjpack#nm_l
| Function : | |||
| nmr_PoloidalEnergySpectrum_wq : | real(8), dimension(0:nm,-nm:nm,km)
| ||
| wq_POLPOT : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
ポロイダルポテンシャルから, ポロイダルエネルギーの 球面調和函数全波数 n, 帯状波数 m の各成分を計算する.
* 全波数 n, 帯状波数 m のポロイダルポテンシャルのスペクトル成分
φ(n,m,r)から全波数 n, 帯状波数 m 成分のポロイダルエネルギー
スペクトルは
(1/2)n(n+1)4π{|d(rφ(n,m,r))/dr|^2 + n(n+1)|φ(n,m,r)|^2}
と計算される.
* 全てのエネルギースペクトル成分の和を動径積分したもの(r^2の重み無し)
が球殻内での全エネルギーに等しい.
* データの存在しない全波数 n, 帯状波数 m の配列には欠損値が格納される.
欠損値の値はモジュール変数 wq_VMiss によって設定できる
(初期値は -999.0)
function nmr_PoloidalEnergySpectrum_wq(wq_POLPOT)
!
! ポロイダルポテンシャルから, ポロイダルエネルギーの
! 球面調和函数全波数 n, 帯状波数 m の各成分を計算する.
!
! * 全波数 n, 帯状波数 m のポロイダルポテンシャルのスペクトル成分
! φ(n,m,r)から全波数 n, 帯状波数 m 成分のポロイダルエネルギー
! スペクトルは
!
! (1/2)n(n+1)4π{|d(rφ(n,m,r))/dr|^2 + n(n+1)|φ(n,m,r)|^2}
!
! と計算される.
!
! * 全てのエネルギースペクトル成分の和を動径積分したもの(r^2の重み無し)
! が球殻内での全エネルギーに等しい.
!
! * データの存在しない全波数 n, 帯状波数 m の配列には欠損値が格納される.
! 欠損値の値はモジュール変数 wq_VMiss によって設定できる
! (初期値は -999.0)
!
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq_POLPOT
!(in) ポロイダルポテンシャル
real(8), dimension(0:nm,-nm:nm,km) :: nmr_PoloidalEnergySpectrum_wq
!(out) エネルギースペクトルポロイダル成分
real(8), dimension((nm+1)*(nm+1),km) ::wr_DATA1 ! 作業領域
real(8), dimension((nm+1)*(nm+1),km) ::wr_DATA2 ! 作業領域
integer :: n, m
nmr_PoloidalEnergySpectrum_wq = wq_VMiss
wr_Data1 = wr_wq(wq_POLPOT)
wr_Data2 = wr_wq(wq_RadDRad_wq(wq_POLPOT)+wq_POLPOT) ! d(rφ)/dr
! = rdφ/dr+φ
do n=0,nm
nmr_PoloidalEnergySpectrum_wq(n,0,:) = + 0.5* n*(n+1)* (4*pi) *( wr_Data2(l_nm(n,0),:)**2 + n*(n+1)*wr_Data1(l_nm(n,0),:)**2 )
do m=1,n
nmr_PoloidalEnergySpectrum_wq(n,m,:) = + 0.5* n*(n+1)* (4*pi) *( wr_Data2(l_nm(n,m),:)**2 + wr_Data2(l_nm(n,-m),:)**2 + n*(n+1)* ( wr_Data1(l_nm(n,m),:)**2 + wr_Data1(l_nm(n,-m),:)**2))
nmr_PoloidalEnergySpectrum_wq(n,-m,:) = nmr_PoloidalEnergySpectrum_wq(n,m,:)
enddo
enddo
end function nmr_PoloidalEnergySpectrum_wq
| Function : | |||
| nmr_ToroidalEnergySpectrum_wq : | real(8), dimension(0:nm,-nm:nm,km)
| ||
| wq_TORPOT : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
トロイダルポテンシャルから, トロイダルエネルギーの 球面調和函数全波数 n, 帯状波数 m の各成分を計算する
* 全波数 n, 帯状波数 m のトロイダルポテンシャルのスペクトル成分 ψ(n,m,r)から全波数 n, 帯状波数 m 成分のトロイダルエネルギー スペクトルは (1/2)n(n+1)4πr^2|ψ(n,m,r)|^2 と計算される. * 全てのエネルギースペクトル成分の和を動径積分したもの(r^2の重み無し) が球殻内での全エネルギーに等しい. * データの存在しない全波数 n, 帯状波数 m の配列には欠損値が格納される. wq_VMiss によって設定できる (初期値は -999.0)
function nmr_ToroidalEnergySpectrum_wq(wq_TORPOT)
!
! トロイダルポテンシャルから, トロイダルエネルギーの
! 球面調和函数全波数 n, 帯状波数 m の各成分を計算する
!
! * 全波数 n, 帯状波数 m のトロイダルポテンシャルのスペクトル成分
! ψ(n,m,r)から全波数 n, 帯状波数 m 成分のトロイダルエネルギー
! スペクトルは (1/2)n(n+1)4πr^2|ψ(n,m,r)|^2 と計算される.
!
! * 全てのエネルギースペクトル成分の和を動径積分したもの(r^2の重み無し)
! が球殻内での全エネルギーに等しい.
!
! * データの存在しない全波数 n, 帯状波数 m の配列には欠損値が格納される.
! wq_VMiss によって設定できる (初期値は -999.0)
!
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq_TORPOT
!(in) トロイダルポテンシャル
real(8), dimension(0:nm,-nm:nm,km) :: nmr_ToroidalEnergySpectrum_wq
!(out) エネルギースペクトルトロイダル成分
real(8), dimension((nm+1)*(nm+1),km) ::wr_DATA ! 作業領域
integer :: n, m
nmr_ToroidalEnergySpectrum_wq = wq_VMiss
wr_DATA = wr_wq(wq_TORPOT)
do n=0,nm
nmr_ToroidalEnergySpectrum_wq(n,0,:) = 0.5 * n*(n+1)* (4*pi) * r_Rad**2 * wr_DATA(l_nm(n,0),:)**2
do m=1,n
nmr_ToroidalEnergySpectrum_wq(n,m,:) = 0.5 * n*(n+1)* (4*pi) * r_Rad**2 * (wr_DATA(l_nm(n,m),:)**2+wr_DATA(l_nm(n,-m),:)**2)
nmr_ToroidalEnergySpectrum_wq(n,-m,:) = nmr_ToroidalEnergySpectrum_wq(n,m,:)
enddo
enddo
end function nmr_ToroidalEnergySpectrum_wq
| Function : | |||
| nr_PoloidalEnergySpectrum_wq : | real(8), dimension(0:nm,km)
| ||
| wq_POLPOT : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
ポロイダルポテンシャルから, ポロイダルエネルギーの 球面調和函数全波数の各成分を計算する
* 全波数 n, 帯状波数 m のポロイダルポテンシャルのスペクトル成分
φ(n,m,r)から全波数 n 成分のポロイダルエネルギースペクトルは
Σ[m=-n]^n ((1/2)n(n+1)4π{|d(rφ(n,m,r))/dr|^2
+ n(n+1)|φ(n,m,r)|^2}
と計算される.
* 全ての全波数に対してのエネルギースペクトル成分の和を動径積分したもの
(r^2の重み無し)が球殻内での全エネルギーに等しい.
function nr_PoloidalEnergySpectrum_wq(wq_POLPOT)
!
! ポロイダルポテンシャルから, ポロイダルエネルギーの
! 球面調和函数全波数の各成分を計算する
!
! * 全波数 n, 帯状波数 m のポロイダルポテンシャルのスペクトル成分
! φ(n,m,r)から全波数 n 成分のポロイダルエネルギースペクトルは
!
! Σ[m=-n]^n ((1/2)n(n+1)4π{|d(rφ(n,m,r))/dr|^2
! + n(n+1)|φ(n,m,r)|^2}
!
! と計算される.
!
! * 全ての全波数に対してのエネルギースペクトル成分の和を動径積分したもの
! (r^2の重み無し)が球殻内での全エネルギーに等しい.
!
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq_POLPOT
!(in) ポロイダルポテンシャル
real(8), dimension(0:nm,km) :: nr_PoloidalEnergySpectrum_wq
!(out) エネルギースペクトルポロイダル成分
real(8), dimension((nm+1)*(nm+1),km) ::wr_DATA1 ! 作業領域
real(8), dimension((nm+1)*(nm+1),km) ::wr_DATA2 ! 作業領域
integer :: n, m
wr_Data1 = wr_wq(wq_POLPOT)
wr_Data2 = wr_wq(wq_RadDRad_wq(wq_POLPOT)+wq_POLPOT) ! d(rφ)/dr
! = rdφ/dr+φ
do n=0,nm
nr_PoloidalEnergySpectrum_wq(n,:) = 0.5* n*(n+1)* (4*pi) *( wr_Data2(l_nm(n,0),:)**2 + n*(n+1)*wr_Data1(l_nm(n,0),:)**2 )
do m=1,n
nr_PoloidalEnergySpectrum_wq(n,:) = nr_PoloidalEnergySpectrum_wq(n,:) + 2 * 0.5* n*(n+1)* (4*pi) *( wr_Data2(l_nm(n,m),:)**2 + wr_Data2(l_nm(n,-m),:)**2 + n*(n+1)*(wr_Data1(l_nm(n,m),:)**2 +wr_Data1(l_nm(n,-m),:)**2))
enddo
enddo
end function nr_PoloidalEnergySpectrum_wq
| Function : | |||
| nr_ToroidalEnergySpectrum_wq : | real(8), dimension(0:nm,km)
| ||
| wq_TORPOT : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
トロイダルポテンシャルから, トロイダルエネルギーの 球面調和函数全波数の各成分を計算する.
* 全波数 n, 帯状波数 m のトロイダルポテンシャルのスペクトル成分 ψ(n,m,r)から全波数 n 成分のトロイダルエネルギースペクトルは Σ[m=-n]^n(1/2)n(n+1)4πr^2|ψ(n,m,r)|^2 と計算される.
が球殻内での全エネルギーに等しい.
function nr_ToroidalEnergySpectrum_wq(wq_TORPOT)
!
! トロイダルポテンシャルから, トロイダルエネルギーの
! 球面調和函数全波数の各成分を計算する.
!
! * 全波数 n, 帯状波数 m のトロイダルポテンシャルのスペクトル成分
! ψ(n,m,r)から全波数 n 成分のトロイダルエネルギースペクトルは
! Σ[m=-n]^n(1/2)n(n+1)4πr^2|ψ(n,m,r)|^2 と計算される.
!
! * 全てのエネルギースペクトル成分の和を動径積分したもの(r^2の重み無し)
! が球殻内での全エネルギーに等しい.
!
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq_TORPOT
!(in) トロイダルポテンシャル
real(8), dimension(0:nm,km) :: nr_ToroidalEnergySpectrum_wq
!(out) エネルギースペクトルトロイダル成分
real(8), dimension((nm+1)*(nm+1),km) ::wr_DATA ! 作業領域
integer :: n, m
wr_DATA = wr_wq(wq_TORPOT)
do n=0,nm
nr_ToroidalEnergySpectrum_wq(n,:) = 0.5 * n*(n+1)* (4*pi) * r_Rad**2 * wr_Data(l_nm(n,0),:)**2
do m=1,n
nr_ToroidalEnergySpectrum_wq(n,:) = nr_ToroidalEnergySpectrum_wq(n,:) + 0.5 * n*(n+1)* (4*pi) * r_Rad**2 * 2* (wr_Data(l_nm(n,m),:)**2 + wr_Data(l_nm(n,-m),:)**2)
enddo
enddo
end function nr_ToroidalEnergySpectrum_wq
| Function : | |||
| q_r2Inv_q : | real(8), dimension(0:km)
| ||
| q_data : | real(8), dimension(:), intent(in)
|
入力スペクトルデータに対して積 r^2 のスペクトル係数 を計算する(1 次元配列用).
Original external subprogram is aq_module#q_r2Inv_q
| Function : | |||
| q_r2_q : | real(8), dimension(0:km)
| ||
| q_data : | real(8), dimension(:), intent(in)
|
入力スペクトルデータに対して積 r^2 のスペクトル係数 を計算する(1 次元配列用).
Original external subprogram is aq_module#q_r2_q
| Function : | |||
| q_rDr_q : | real(8), dimension(0:km)
| ||
| q_data : | real(8), dimension(:), intent(in)
|
入力スペクトルデータに r(d/dR) 微分を作用する(1 次元配列用).
スペクトルデータの r(d/dR) 微分とは, 対応する格子点データに R 微分を 作用させたデータのスペクトル変換のことである.
Original external subprogram is aq_module#q_rDr_q
| Function : | |||
| r_AvrLat_yr : | real(8), dimension(km)
| ||
| yr : | real(8), dimension(1:jm,km), intent(in)
|
2 次元(YR)格子点データの緯度方向域平均.
2 次元データ f(φ,r) に対して ∫f(φ,r) cosφ dφ/2 を計算する.
function r_AvrLat_yr(yr)
!
! 2 次元(YR)格子点データの緯度方向域平均.
!
! 2 次元データ f(φ,r) に対して ∫f(φ,r) cosφ dφ/2 を計算する.
!
real(8), dimension(1:jm,km), intent(in) :: yr
!(in) 2 次元緯度動径(子午面)格子点データ
real(8), dimension(km) :: r_AvrLat_yr
!(out) 緯度平均された 1 次元動径格子点データ
r_AvrLat_yr = r_IntLat_yr(yr)/sum(y_Lat_Weight)
end function r_AvrLat_yr
| Function : | |||
| r_AvrLonLat_xyr : | real(8), dimension(km)
| ||
| xyr : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
緯度経度(水平)積分
3 次元格子点データの緯度経度(水平, 球面)積分
3 次元データ f(λ,φ,r) に対して
∫f(λ,φ,r) cosφ dλdφ /4π
を計算する.
function r_AvrLonLat_xyr(xyr) ! 緯度経度(水平)積分
!
! 3 次元格子点データの緯度経度(水平, 球面)積分
!
! 3 次元データ f(λ,φ,r) に対して
!
! ∫f(λ,φ,r) cosφ dλdφ /4π
!
! を計算する.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr
!(in) 3 次元経度緯度動径格子点データ
real(8), dimension(km) :: r_AvrLonLat_xyr
!(out) 緯度経度(水平, 球面)平均された 1 次元動径格子点データ
r_AvrLonLat_xyr = r_IntLonLat_xyr(xyr) /(sum(x_Lon_Weight)*sum(y_Lat_Weight))
end function r_AvrLonLat_xyr
| Function : | |||
| r_AvrLon_xr : | real(8), dimension(km)
| ||
| xr : | real(8), dimension(0:im-1,km), intent(in)
|
経度(帯状)積分
2 次元(XR)格子点データの経度方向平均.
2 次元データ f(λ,r) に対して ∫f(λ,r)dλ/2π を計算する.
function r_AvrLon_xr(xr) ! 経度(帯状)積分
!
! 2 次元(XR)格子点データの経度方向平均.
!
! 2 次元データ f(λ,r) に対して ∫f(λ,r)dλ/2π を計算する.
!
real(8), dimension(0:im-1,km), intent(in) :: xr
!(in) 2 次元緯度動径格子点データ
real(8), dimension(km) :: r_AvrLon_xr
!(out) 経度平均された 1 次元動径格子点データ
r_AvrLon_xr = r_IntLon_xr(xr)/sum(x_Lon_Weight)
end function r_AvrLon_xr
| Function : | |||
| r_IntLat_yr : | real(8), dimension(km)
| ||
| yr : | real(8), dimension(1:jm,km), intent(in)
|
緯度積分
2 次元(YR)格子点データの緯度方向域積分.
2 次元データ f(φ,r) に対して∫f(φ,r) cosφ dφ を計算する.
function r_IntLat_yr(yr) ! 緯度積分
!
! 2 次元(YR)格子点データの緯度方向域積分.
!
! 2 次元データ f(φ,r) に対して∫f(φ,r) cosφ dφ を計算する.
!
real(8), dimension(1:jm,km), intent(in) :: yr
!(in) 2 次元緯度動径(子午面)格子点データ
real(8), dimension(km) :: r_IntLat_yr
!(out) 緯度積分された 1 次元動径格子点データ
integer :: j
r_IntLat_yr = 0.0d0
do j=1,jm
r_IntLat_yr(:) = r_IntLat_yr(:) + yr(j,:) * y_Lat_Weight(j)
enddo
end function r_IntLat_yr
| Function : | |||
| r_IntLonLat_xyr : | real(8), dimension(km)
| ||
| xyr : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
緯度経度(水平)積分
3 次元格子点データの緯度経度(水平, 球面)積分
3 次元データ f(λ,φ,r) に対して
∫f(λ,φ,r) cosφ dλdφ
を計算する.
function r_IntLonLat_xyr(xyr) ! 緯度経度(水平)積分
!
! 3 次元格子点データの緯度経度(水平, 球面)積分
!
! 3 次元データ f(λ,φ,r) に対して
!
! ∫f(λ,φ,r) cosφ dλdφ
!
! を計算する.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr
!(in) 3 次元経度緯度動径格子点データ
real(8), dimension(km) :: r_IntLonLat_xyr
!(out) 緯度経度(水平, 球面)積分された 1 次元動径格子点データ
integer :: i, j
r_IntLonLat_xyr = 0.0D0
do j=1,jm
do i=0,im-1
r_IntLonLat_xyr = r_IntLonLat_xyr + xyr(i,j,:) * x_Lon_Weight(i) * y_Lat_Weight(j)
enddo
enddo
end function r_IntLonLat_xyr
| Function : | |||
| r_IntLon_xr : | real(8), dimension(km)
| ||
| xr : | real(8), dimension(0:im-1,km), intent(in)
|
2 次元(XR)格子点データの経度方向積分.
2 次元データ f(λ,r) に対して ∫f(λ,r)dλ を計算する.
function r_IntLon_xr(xr)
!
! 2 次元(XR)格子点データの経度方向積分.
!
! 2 次元データ f(λ,r) に対して ∫f(λ,r)dλ を計算する.
!
real(8), dimension(0:im-1,km), intent(in) :: xr
!(in) 2 次元緯度動径格子点データ
real(8), dimension(km) :: r_IntLon_xr
!(out) 経度積分された 1 次元動径格子点データ
integer :: i
r_IntLon_xr = 0.0d0
do i=0,im-1
r_IntLon_xr(:) = r_IntLon_xr(:) + xr(i,:) * x_Lon_Weight(i)
enddo
end function r_IntLon_xr
| Variable : | |||
| g_R_Weight(:) : | real(8), allocatable
|
Original external subprogram is aq_module#g_R_Weight
| Variable : | |||
| g_R(:) : | real(8), allocatable
|
Original external subprogram is aq_module#g_R
| Function : | |||
| w_xy((mm+1)*(mm+1)) : | real(8)
| ||
| xy_data(0:im-1,1:jm) : | real(8), intent(in)
| ||
| ipow : | integer, intent(in), optional
| ||
| iflag : | integer, intent(in), optional
|
格子データからスペクトルデータへ(正)変換する(1 層用).
Original external subprogram is wa_module_sjpack#w_xy
| Subroutine : | |||
| wq : | real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
| ||
| value : | real(8), dimension((nm+1)*(nm+1)), intent(in), optional
| ||
| cond : | character(len=1), intent(in), optional
|
スペクトルデータにディリクレ・ノイマン境界条件を適用する Chebyshev 空間での境界条件適用(タウ法)
チェビシェフ空間において境界条件を満たすべく高次の係数を 定める方法をとっている(タウ法).
Alias for wq_BoundaryTau
| Subroutine : | |||
| wq : | real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
| ||
| value : | real(8), dimension((nm+1)*(nm+1)), intent(in), optional
| ||
| cond : | character(len=1), intent(in), optional
|
スペクトルデータにディリクレ・ノイマン境界条件を適用する Chebyshev 空間での境界条件適用(タウ法)
チェビシェフ空間において境界条件を満たすべく高次の係数を 定める方法をとっている(タウ法).
subroutine wq_BoundaryTau(wq,value,cond)
!
! スペクトルデータにディリクレ・ノイマン境界条件を適用する
! Chebyshev 空間での境界条件適用(タウ法)
!
! チェビシェフ空間において境界条件を満たすべく高次の係数を
! 定める方法をとっている(タウ法).
!
real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout) :: wq
!(inout) 境界条件を適用するデータ. 修正された値を返す.
real(8), dimension((nm+1)*(nm+1)), intent(in), optional :: value
!(in) 境界での 値/勾配 分布を水平スペクトル変換したものを与える.
! 省略時は値/勾配 0 となる.
character(len=1), intent(in), optional :: cond
!(in) 境界条件. 省略時は 'D'
! D : 外側ディリクレ条件
! N : 外側ノイマン条件
if (.not. present(cond)) then
if (present(value)) then
call aq_BoundaryTau_D(wq,value)
else
call aq_BoundaryTau_D(wq)
endif
return
endif
select case(cond)
case ('N')
if (present(value)) then
call aq_BoundaryTau_N(wq,value)
else
call aq_BoundaryTau_N(wq)
endif
case ('D')
if (present(value)) then
call aq_BoundaryTau_D(wq,value)
else
call aq_BoundaryTau_D(wq)
endif
case default
call MessageNotify('E','wq_BoundaryTau','B.C. not supported')
end select
end subroutine wq_BoundaryTau
| Subroutine : | |||
| i : | integer,intent(in)
| ||
| j : | integer,intent(in)
| ||
| k : | integer,intent(in)
| ||
| n : | integer,intent(in)
| ||
| l : | integer,intent(in)
| ||
| r : | real(8),intent(in)
| ||
| np : | integer,intent(in), optional
| ||
| wa_init : | logical,intent(in), optional
|
スペクトル変換の格子点数, 波数, 動径座標の範囲を設定する.
他の関数を呼ぶ前に, 最初にこのサブルーチンを呼んで初期設定を しなければならない.
np に 1 より大きな値を指定すれば ISPACK の球面調和函数変換 OPENMP 並列計算ルーチンが用いられる. 並列計算を実行するには, 実行時に環境変数 OMP_NUM_THREADS を np 以下の数字に設定する等の システムに応じた準備が必要となる.
np に 1 より大きな値を指定しなければ並列計算ルーチンは呼ばれない.
subroutine wq_Initial(i,j,k,n,l,r,np,wa_init)
!
! スペクトル変換の格子点数, 波数, 動径座標の範囲を設定する.
!
! 他の関数を呼ぶ前に, 最初にこのサブルーチンを呼んで初期設定を
! しなければならない.
!
! np に 1 より大きな値を指定すれば ISPACK の球面調和函数変換
! OPENMP 並列計算ルーチンが用いられる. 並列計算を実行するには,
! 実行時に環境変数 OMP_NUM_THREADS を np 以下の数字に設定する等の
! システムに応じた準備が必要となる.
!
! np に 1 より大きな値を指定しなければ並列計算ルーチンは呼ばれない.
!
!
integer,intent(in) :: i ! 格子点数(経度λ)
integer,intent(in) :: j ! 格子点数(緯度φ)
integer,intent(in) :: k ! 格子点数(動径 r)
integer,intent(in) :: n ! 切断波数(水平全波数)
integer,intent(in) :: l ! 切断波数(動径波数)
real(8),intent(in) :: r ! 球半径
integer,intent(in), optional :: np ! OPENMP での最大スレッド数
logical,intent(in), optional :: wa_init ! wa_initial スイッチ
logical :: wa_initialize=.true. ! wa_initial スイッチ
integer :: nn, mm
im = i
jm = j
km = k
nm = n
lm = l
ra = r
if ( present(wa_init) ) then
wa_initialize = wa_init
else
wa_initialize = .true.
endif
if ( wa_initialize ) then
if ( present(np) ) then
call wa_Initial(nm,im,jm,km,np)
else
call wa_Initial(nm,im,jm,km)
endif
endif
allocate(nd((nm+1)*(nm+1)))
do nn=0,nm
do mm=-nn,nn
nd(l_nm(nn,mm)) = nn ! 水平 n 次は r^n から
enddo
enddo
call aq_Initial(km,lm,ra,alpha,beta,nd)
allocate(xyr_Lon(0:im-1,1:jm,km))
allocate(xyr_Lat(0:im-1,1:jm,km))
allocate(xyr_Rad(0:im-1,1:jm,km))
allocate(wr_Rad((nm+1)*(nm+1),km))
xyr_Lon = spread(xy_Lon,3,km)
xyr_Lat = spread(xy_Lat,3,km)
xyr_Rad = spread(spread(r_Rad,1,jm),1,im)
wr_Rad = spread(r_Rad,1,(nm+1)*(nm+1))
call MessageNotify('M','wq_initial','wq_module_sjpack (2009/12/06) is initialized')
end subroutine wq_initial
| Function : | |||
| wq_KxRGrad_wq : | real(8), dimension((nm+1)*(nm+1),0:lm)
| ||
| wq : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
入力スペクトルデータに経度微分 k×r・▽ = ∂/∂λを作用する.
function wq_KxRGrad_wq(wq)
!
! 入力スペクトルデータに経度微分 k×r・▽ = ∂/∂λを作用する.
!
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq
!(in) 2 次元球面調和函数チェビシェフスペクトルデータ
real(8), dimension((nm+1)*(nm+1),0:lm) :: wq_KxRGrad_wq
!(out) 経度微分を作用された 2 次元スペクトルデータ
wq_KxRGrad_wq = wa_Dlon_wa(wq)
end function wq_KxRGrad_wq
| Function : | |||
| wq_L2Inv_wq : | real(8), dimension((nm+1)*(nm+1),0:lm)
| ||
| wq : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
入力スペクトルデータに L^2 演算子の逆演算(-逆水平ラプラシアン)を 作用する.
スペクトルデータに L^2 演算子を作用させる関数 wq_L2_wq の逆計算を 行う関数である.
function wq_L2Inv_wq(wq)
!
! 入力スペクトルデータに L^2 演算子の逆演算(-逆水平ラプラシアン)を
! 作用する.
!
! スペクトルデータに L^2 演算子を作用させる関数 wq_L2_wq の逆計算を
! 行う関数である.
!
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq
!(in) 2 次元球面調和函数チェビシェフスペクトルデータ
real(8), dimension((nm+1)*(nm+1),0:lm) :: wq_L2Inv_wq
!(out) L^2 演算子の逆演算を作用された 2 次元スペクトルデータ
wq_L2Inv_wq = -wa_LaplaInv_wa(wq)
end function wq_L2Inv_wq
| Function : | |||
| wq_L2_wq : | real(8), dimension((nm+1)*(nm+1),0:lm)
| ||
| wq : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
入力スペクトルデータに L^2 演算子(=-水平ラプラシアン)を作用する.
L^2 演算子は単位球面上の水平ラプラシアンの逆符号にあたる.
入力スペクトルデ ータに対応する格子点データに演算子
L^2 = -1/cos^2φ・∂^2/∂λ^2 - 1/cosφ・∂/∂φ(cosφ∂/∂φ)
を作用させたデータのスペクトル変換が返される.
function wq_L2_wq(wq)
!
! 入力スペクトルデータに L^2 演算子(=-水平ラプラシアン)を作用する.
!
! L^2 演算子は単位球面上の水平ラプラシアンの逆符号にあたる.
! 入力スペクトルデ ータに対応する格子点データに演算子
!
! L^2 = -1/cos^2φ・∂^2/∂λ^2 - 1/cosφ・∂/∂φ(cosφ∂/∂φ)
!
! を作用させたデータのスペクトル変換が返される.
!
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq
!(in) 2 次元球面調和函数チェビシェフスペクトルデータ
real(8), dimension((nm+1)*(nm+1),0:lm) :: wq_L2_wq
!(out) L^2 演算子を作用された 2 次元スペクトルデータ
wq_L2_wq = -wa_Lapla_wa(wq)
end function wq_L2_wq
| Function : | |||
| wq_LaplaPol2PolTau_wq : | real(8), dimension((nm+1)*(nm+1),0:lm)
| ||
| wq : | real(8), dimension((nm+1)*(nm+1),0:lm),intent(in)
| ||
| cond : | character(len=1), intent(in), optional
| ||
| new : | logical, intent(IN), optional
|
速度ポロイダルポテンシャルΦを▽^2Φから計算する.
スペクトル空間で境界条件を適用している(タウ法).
速度ポロイダルポテンシャルΦを f = ▽^2Φから定める式は
▽^2Φ = f
Φ = const. at Boundary.
∂Φ/∂r = 0 at Boundary (粘着条件)
or ∂^2Φ/∂r^2 = 0 at Boundary (応力なし条件)
最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.
function wq_LaplaPol2PolTau_wq(wq,cond,new)
!
! 速度ポロイダルポテンシャルΦを▽^2Φから計算する.
!
! スペクトル空間で境界条件を適用している(タウ法).
!
! 速度ポロイダルポテンシャルΦを f = ▽^2Φから定める式は
!
! ▽^2Φ = f
! Φ = const. at Boundary.
! ∂Φ/∂r = 0 at Boundary (粘着条件)
! or ∂^2Φ/∂r^2 = 0 at Boundary (応力なし条件)
!
! 最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.
!
real(8), dimension((nm+1)*(nm+1),0:lm),intent(in) :: wq
!(in) 入力▽^2φ分布
real(8), dimension((nm+1)*(nm+1),0:lm) :: wq_LaplaPol2PolTau_wq
!(out) 出力ポロイダルポテンシャル分布
character(len=1), intent(in), optional :: cond
!(in) 境界条件スイッチ. 省略時は 'R'
! R : 上側粘着条件
! F : 上側応力なし条件
logical, intent(IN), optional :: new
!(in) true だと境界条件計算用行列を強制的に新たに作る.
! default は false.
real(8), dimension(:,:,:), allocatable :: alu ! 内部領域計算用
integer, dimension(:,:), allocatable :: kp ! 内部領域計算用
real(8), dimension(:,:,:), allocatable :: alub ! 境界条件計算用
integer, dimension(:,:), allocatable :: kpb ! 境界条件計算用
real(8), dimension((nm+1)*(nm+1),km) :: wr_work
real(8), dimension((nm+1)*(nm+1),0:lm) :: wq_work
logical :: rigid ! 境界条件
logical :: first = .true.
logical :: new_matrix = .false.
integer :: l,n, lend
save :: alu, kp, first
save :: alub, kpb
if (.not. present(cond)) then
rigid=.TRUE.
else
select case (cond)
case ('R')
rigid = .TRUE.
case ('F')
rigid = .FALSE.
case default
call MessageNotify('E','wq_laplapol2pol_wq','B.C. not supported')
end select
endif
if (.not. present(new)) then
new_matrix=.false.
else
new_matrix=new
endif
if ( first .OR. new_matrix ) then
first = .false.
if ( allocated(alu) ) deallocate(alu)
if ( allocated(kp) ) deallocate(kp)
if ( allocated(alub) ) deallocate(alub)
if ( allocated(kpb) ) deallocate(kpb)
allocate(alu((nm+1)*(nm+1),0:lm,0:lm),kp((nm+1)*(nm+1),0:lm))
allocate(alub((nm+1)*(nm+1),0:lm,0:lm),kpb((nm+1)*(nm+1),0:lm))
!---- 内部領域計算用行列 -----
do l=0,lm
wq_work = 0.0
wq_work(:,l) = 1.0D0
alu(:,:,l) = wq_Lapla_wq(wq_work)
enddo
! 0 成分のところを 1 で埋める.
do n=1,(nm+1)**2
do l=0,nd(n)-1
alu(n,l,l) = 1.0D0
enddo
do l=nd(n)+1,lm,2
alu(n,l,l) = 1.0D0
enddo
enddo
! alu(:,:,nd(n)) 列は 0 なので 1 をいれておく.
! l=nd(n) 成分は境界条件で決める.
do n=1,(nm+1)**2
if ( mod(nd(n),2) .eq. mod(lm,2) ) then
alu(n,lm,nd(n)) = 1.0D0
else
alu(n,lm-1,nd(n)) = 1.0D0
endif
enddo
call ludecomp(alu,kp)
!---- 境界条件計算用行列 -----
alub = 0.0D0
do l=0,lm
alub(:,l,l) = 1.0D0
enddo
do l=0,lm
wq_work=0.0D0
wq_work(:,l)=1.0D0
wr_work = wr_wq(wq_work)
! 運動学的条件. 流線は境界で一定
! l=nd(n) 成分を境界条件で決める.
do n=1,(nm+1)**2
alub(n,nd(n),l) = wr_work(n,km)
enddo
! 力学的条件粘着条件
! l=lend 成分を境界条件で決める.
if ( rigid ) then
wr_work=wr_wq(wq_RadDRad_wq(wq_work))/wr_Rad
else
wr_work=wr_wq(wq_RadDRad_wq(wq_RadDRad_wq(wq_work)) -wq_RadDRad_wq(wq_work) ) /wr_Rad**2
endif
do n=1,(nm+1)**2
if ( mod(nd(n),2) .eq. mod(lm,2) ) then
lend = lm
else
lend = lm-1
endif
alub(n,lend,l) = wr_work(n,km)
enddo
enddo
! 関係ないところを 0 で埋める.
do n=1,(nm+1)**2
if ( mod(nd(n),2) .eq. mod(lm,2) ) then
lend = lm
else
lend = lm-1
endif
do l=0,nd(n)-1
alub(n,nd(n),l) = 0.0D0
alub(n,lend,l) = 0.0D0
enddo
do l=nd(n)+1,lm,2
alub(n,nd(n),l) = 0.0D0
alub(n,lend,l) = 0.0D0
enddo
enddo
call ludecomp(alub,kpb)
if ( rigid ) then
call MessageNotify('M','wq_LaplaPol2PolTau_wq', 'Matrix to apply rigid b.c. newly produced.')
else
call MessageNotify('M','wq_LaplaPol2PolTau_wq', 'Matrix to apply stress-free b.c. newly produced.')
endif
endif
! 内部領域計算
wq_work = wq
wq_work = lusolve(alu,kp,wq_work)
! 境界条件計算
do n=1,(nm+1)**2
wq_work(n,nd(n)) = 0
if ( mod(nd(n),2) .eq. mod(lm,2) ) then
wq_work(n,lm) = 0
else
wq_work(n,lm-1) = 0
endif
enddo
wq_laplapol2polTau_wq = lusolve(alub,kpb,wq_work)
end function wq_LaplaPol2PolTau_wq
| Function : | |||
| wq_LaplaPol2PolTau_wq : | real(8), dimension((nm+1)*(nm+1),0:lm)
| ||
| wq : | real(8), dimension((nm+1)*(nm+1),0:lm),intent(in)
| ||
| cond : | character(len=1), intent(in), optional
| ||
| new : | logical, intent(IN), optional
|
速度ポロイダルポテンシャルΦを▽^2Φから計算する.
スペクトル空間で境界条件を適用している(タウ法).
速度ポロイダルポテンシャルΦを f = ▽^2Φから定める式は
▽^2Φ = f
Φ = const. at Boundary.
∂Φ/∂r = 0 at Boundary (粘着条件)
or ∂^2Φ/∂r^2 = 0 at Boundary (応力なし条件)
最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.
Alias for wq_LaplaPol2PolTau_wq
| Function : | |||
| wq_Lapla_wq : | real(8), dimension((nm+1)*(nm+1),0:lm)
| ||
| wq : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
入力スペクトルデータにラプラシアン
▽^2 = 1/r^2 cos^2φ・∂^2/∂λ^2
+ 1/r^2 cosφ・∂/∂φ(cosφ∂/∂φ)
+ 1/r^2 ∂/∂r (r^2 ∂/∂r)
[1/r^2 (r∂/∂r)(r∂/∂r) + 1/r^2(r∂/∂r)]
を作用する.
スペクトルデータのラプラシアンとは, 対応する格子点データに ラプラシアンを作用させたデータのスペクトル変換のことである.
function wq_Lapla_wq(wq)
! 入力スペクトルデータにラプラシアン
!
! ▽^2 = 1/r^2 cos^2φ・∂^2/∂λ^2
! + 1/r^2 cosφ・∂/∂φ(cosφ∂/∂φ)
! + 1/r^2 ∂/∂r (r^2 ∂/∂r)
! [1/r^2 (r∂/∂r)(r∂/∂r) + 1/r^2(r∂/∂r)]
!
! を作用する.
!
! スペクトルデータのラプラシアンとは, 対応する格子点データに
! ラプラシアンを作用させたデータのスペクトル変換のことである.
!
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq
!(in) 2 次元球面調和函数チェビシェフスペクトルデータ
real(8), dimension((nm+1)*(nm+1),0:lm) :: wq_Lapla_wq
!(out) ラプラシアンを作用された水平スペクトルデータ
wq_Lapla_wq = wq_Rad2Inv_wq( wq_RadDRad_wq(wq_RadDRad_wq(wq)) + wq_RadDRad_wq(wq)+ wa_Lapla_wa(wq) )
end function wq_Lapla_wq
| Subroutine : | |||
| wq_POL : | real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
| ||
| new : | logical, intent(IN), optional
|
磁場ポロイダルポテンシャルに対して境界条件を適用する. Chebyshev 空間での境界条件適用
チェビシェフ空間において境界条件を満たすべく高次の係数を定める方法を とっている(タウ法). 現在のところ境界物質が非電気伝導体の場合のみ 対応している. その場合, 磁場ポロイダルポテンシャルの各水平スペクトル 成分 h にたいして境界条件が与えられ,
* 外側境界 : dh/dr + (n+1)h/r = 0
である. ここで n は h の水平全波数である.
最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.
Alias for wq_PolMagBoundaryTau
| Subroutine : | |||
| wq_POL : | real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
| ||
| new : | logical, intent(IN), optional
|
磁場ポロイダルポテンシャルに対して境界条件を適用する. Chebyshev 空間での境界条件適用
チェビシェフ空間において境界条件を満たすべく高次の係数を定める方法を とっている(タウ法). 現在のところ境界物質が非電気伝導体の場合のみ 対応している. その場合, 磁場ポロイダルポテンシャルの各水平スペクトル 成分 h にたいして境界条件が与えられ,
* 外側境界 : dh/dr + (n+1)h/r = 0
である. ここで n は h の水平全波数である.
最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.
subroutine wq_PolmagBoundaryTau(wq_POL,new)
!
! 磁場ポロイダルポテンシャルに対して境界条件を適用する.
! Chebyshev 空間での境界条件適用
!
! チェビシェフ空間において境界条件を満たすべく高次の係数を定める方法を
! とっている(タウ法). 現在のところ境界物質が非電気伝導体の場合のみ
! 対応している. その場合, 磁場ポロイダルポテンシャルの各水平スペクトル
! 成分 h にたいして境界条件が与えられ,
!
! * 外側境界 : dh/dr + (n+1)h/r = 0
!
! である. ここで n は h の水平全波数である.
!
! 最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.
!
real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout) :: wq_POL
!(inout) 境界条件を適用するデータ. 修正された値を返す.
logical, intent(IN), optional :: new
!(in) true だと境界条件計算用行列を強制的に新たに作る.
! default は false.
real(8), dimension(:,:,:), allocatable :: alu
integer, dimension(:,:), allocatable :: kp
real(8), dimension(:,:), allocatable :: wq_I
real(8), dimension(:,:), allocatable :: wr_PSI
real(8), dimension(:,:), allocatable :: wr_DPSIDR
logical :: first = .true.
logical :: new_matrix = .false.
integer :: l, n, nn(2), lend
save :: alu, kp, first
if (.not. present(new)) then
new_matrix=.false.
else
new_matrix=new
endif
if ( first .OR. new_matrix ) then
first = .false.
if ( allocated(alu) ) deallocate(alu)
if ( allocated(kp) ) deallocate(kp)
if ( allocated(wq_I) ) deallocate(wq_I)
if ( allocated(wr_PSI) ) deallocate(wr_PSI)
if ( allocated(wr_DPSIDR) ) deallocate(wr_DPSIDR)
allocate(alu((nm+1)*(nm+1),0:lm,0:lm),kp((nm+1)*(nm+1),0:lm))
allocate(wq_I((nm+1)*(nm+1),0:lm))
allocate(wr_PSI((nm+1)*(nm+1),km),wr_DPSIDR((nm+1)*(nm+1),km))
alu = 0.0D0
do l=0,lm
alu(:,l,l) = 1.0D0
enddo
! 非電気伝導体
do l=0,lm
wq_I = 0.0D0
wq_I(:,l) = 1.0D0
wr_PSI = wr_wq(wq_I)
wr_DPSIDR = wr_wq(wq_RadDRad_wq(wq_I))/wr_Rad
do n=1,(nm+1)**2
if ( mod(nd(n),2) .eq. mod(lm,2) ) then
lend = lm
else
lend = lm-1
endif
nn=nm_l(n)
alu(n,lend,l) = wr_DPSIDR(n,km) + (nn(1)+1) * wr_PSI(n,km)/r_RAD(km)
enddo
enddo
! 関係ないところを 0 で埋める.
do n=1,(nm+1)**2
if ( mod(nd(n),2) .eq. mod(lm,2) ) then
lend = lm
else
lend = lm-1
endif
do l=0,nd(n)-1
alu(n,lend,l) = 0.0D0
enddo
do l=nd(n)+1,lm,2
alu(n,lend,l) = 0.0D0
enddo
enddo
call ludecomp(alu,kp)
deallocate(wq_I,wr_PSI,wr_DPSIDR)
call MessageNotify('M','PolmagBoundaryTau', 'Matrix to apply b.c. newly produced.')
endif
do n=1,(nm+1)**2
if ( mod(nd(n),2) .eq. mod(lm,2) ) then
wq_POL(n,lm) = 0.0
else
wq_POL(n,lm-1) = 0.0
endif
enddo
wq_POL = lusolve(alu,kp,wq_POL)
end subroutine wq_PolmagBoundaryTau
| Subroutine : | |||
| xyr_RotVLON : | real(8), dimension(0:im-1,1:jm,km), intent(OUT)
| ||
| xyr_RotVLAT : | real(8), dimension(0:im-1,1:jm,km), intent(OUT)
| ||
| xyr_RotVRAD : | real(8), dimension(0:im-1,1:jm,km), intent(OUT)
| ||
| wq_TORPOT : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
| ||
| wq_POLPOT : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
トロイダルポロイダルポテンシャルΨ,Φで表される非発散ベクトル場
v = ▽x(Ψr) + ▽x▽x(Φr)
に対して, その回転
▽xv = ▽x▽x(Ψr) + ▽x▽x▽x(Φr) = ▽x▽x(Ψr) - ▽x((▽^2Φ)r)
を計算する.
subroutine wq_Potential2Rotation( xyr_RotVLON,xyr_RotVLAT,xyr_RotVRAD,wq_TORPOT,wq_POLPOT)
!
! トロイダルポロイダルポテンシャルΨ,Φで表される非発散ベクトル場
!
! v = ▽x(Ψr) + ▽x▽x(Φr)
!
! に対して, その回転
!
! ▽xv = ▽x▽x(Ψr) + ▽x▽x▽x(Φr) = ▽x▽x(Ψr) - ▽x((▽^2Φ)r)
!
! を計算する.
! ベクトル場の回転
real(8), dimension(0:im-1,1:jm,km), intent(OUT) :: xyr_RotVLON
!(out) 回転の経度成分
real(8), dimension(0:im-1,1:jm,km), intent(OUT) :: xyr_RotVLAT
!(out) 回転の緯度成分
real(8), dimension(0:im-1,1:jm,km), intent(OUT) :: xyr_RotVRAD
!(out) 回転の動径成分
! 入力ベクトル場を表すポテンシャル
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq_TORPOT
!(in) トロイダルポテンシャル
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq_POLPOT
!(in) ポロイダルポテンシャル
call wq_Potential2Vector( xyr_RotVLON,xyr_RotVLAT,xyr_RotVRAD, -wq_Lapla_wq(wq_POLPOT), wq_TORPOT)
end subroutine wq_Potential2Rotation
| Subroutine : | |||
| xyr_VLON : | real(8), dimension(0:im-1,1:jm,km)
| ||
| xyr_VLAT : | real(8), dimension(0:im-1,1:jm,km)
| ||
| xyr_VRAD : | real(8), dimension(0:im-1,1:jm,km)
| ||
| wq_TORPOT : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
| ||
| wq_POLPOT : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
トロイダルポロイダルポテンシャルΨ,Φで表される非発散ベクトル場
v = ▽x(Ψr) + ▽x▽x(Φr)
の各成分を計算する
subroutine wq_Potential2Vector( xyr_VLON,xyr_VLAT,xyr_VRAD,wq_TORPOT,wq_POLPOT)
!
! トロイダルポロイダルポテンシャルΨ,Φで表される非発散ベクトル場
!
! v = ▽x(Ψr) + ▽x▽x(Φr)
!
! の各成分を計算する
!
real(8), dimension(0:im-1,1:jm,km) :: xyr_VLON
!(out) ベクトル場の経度成分
real(8), dimension(0:im-1,1:jm,km) :: xyr_VLAT
!(out) ベクトル場の緯度成分
real(8), dimension(0:im-1,1:jm,km) :: xyr_VRAD
!(out) ベクトル場の動径成分
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq_TORPOT
!(in) トロイダルポテンシャル
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq_POLPOT
!(in) ポロイダルポテンシャル
xyr_VLON = xyr_RAD * xyr_GradLat_wq(wq_TORPOT) + xya_GradLon_wa(wr_RotDRad_wq(wq_POLPOT))
xyr_VLAT = - xyr_RAD * xyr_GradLon_wq(wq_TORPOT) + xya_GradLat_wa(wr_RotDRad_wq(wq_POLPOT))
xyr_VRAD = xyr_wq(wq_L2_wq(wq_POLPOT))/xyr_RAD
end subroutine wq_Potential2Vector
| Function : | |||
| wq_QOperator_wq : | real(8), dimension((nm+1)*(nm+1),0:lm)
| ||
| wq : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
入力スペクトルデータに対応する格子点データに演算子
Q=(k・▽-1/2(L2 k・▽+ k・▽L2))
を作用させたデータのスペクトル変換が返される.
function wq_QOperator_wq(wq)
!
! 入力スペクトルデータに対応する格子点データに演算子
!
! Q=(k・▽-1/2(L2 k・▽+ k・▽L2))
!
! を作用させたデータのスペクトル変換が返される.
!
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq
!(in) 2 次元球面調和函数チェビシェフスペクトルデータ
real(8), dimension((nm+1)*(nm+1),0:lm) :: wq_QOperator_wq
!(out) Q 演算子を作用された 2 次元スペクトルデータ
wq_QOperator_wq = wq_xyr(xyr_KGrad_wq(wq) - xyr_KGrad_wq(wq_L2_wq(wq))/2) - wq_L2_wq(wq_xyr(xyr_KGrad_wq(wq)))/2
end function wq_QOperator_wq
| Function : | |||
| aq_r2Inv_aq : | real(8), dimension(size(aq_data,1),0:size(aq_data,2)-1)
| ||
| aq_data : | real(8), dimension(:,0:), intent(in)
|
入力スペクトルデータに対して積 r^-2 のスペクトル係数 を計算する(2 次元配列用).
a_n^m = aq_r2Inv_aq/sqrt(Inm), b_n^m = aq_data/sqrt(Inm)
b_n^m = (n-|m|)(n+|m|+beta-1)/((2n+gamma-5)(2n+gamma-3)) a_n-2^m
+ (2n(n+gamma-1) + 2|m|(|m|+beta-1)+(gamma-3)(beta+1)
/((2n+gamma+1)(2n+gamma-3)) a_n^m
+ (n-|m|+gamma-beta)(n+|m|+gamma-1)
/((2n+gamma+3)(2n+gamma+1)) a_n+2^m
Original external subprogram is aq_module#aq_r2Inv_aq
| Function : | |||
| aq_r2_aq : | real(8), dimension(size(aq_data,1),0:size(aq_data,2)-1)
| ||
| aq_data : | real(8), dimension(:,0:), intent(in)
|
入力スペクトルデータに対して積 r^2 のスペクトル係数 を計算する(2 次元配列用).
a_n^m = aq_data/sqrt(Inm), b_n^m = aq_rDr_aq/sqrt(Inm)
b_n^m = (n-|m|)(n+|m|+beta-1)/((2n+gamma-5)(2n+gamma-3)) a_n-2^m
+ (2n(n+gamma-1) + 2|m|(|m|+beta-1)+(gamma-3)(beta+1)
/((2n+gamma+1)(2n+gamma-3)) a_n^m
+ (n-|m|+gamma-beta)(n+|m|+gamma-1)
/((2n+gamma+3)(2n+gamma+1)) a_n+2^m
Original external subprogram is aq_module#aq_r2_aq
| Function : | |||
| aq_rDr_aq : | real(8), dimension(size(aq_data,1),0:size(aq_data,2)-1)
| ||
| aq_data : | real(8), dimension(:,0:), intent(in)
|
入力スペクトルデータに対して微分 r(d/dr) のスペクトル係数 を計算する(2 次元配列用).
a_n = aq_data/sqrt(Inm), b_n = aq_rDr_aq/sqrt(Inm)
b_n = (2n+gamma-1)/(2n+gamma+3)b_n+2
+ (2n+gamma-1)(n+gamma+1)/(2n+gamma+3)a_n+2 + n a_n
Original external subprogram is aq_module#aq_rDr_aq
| Function : | |||
| wq_RadRotRot_xyr_xyr_xyr : | real(8), dimension((nm+1)*(nm+1),0:lm)
| ||
| xyr_VLON : | real(8), dimension(0:im-1,1:jm,km), intent(in)
| ||
| xyr_VLAT : | real(8), dimension(0:im-1,1:jm,km), intent(in)
| ||
| xyr_VRAD : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
ベクトル v に対して r・(▽×▽×v) を計算する.
第 1, 2, 3 引数(v[λ], v[φ], v[r])がそれぞれベクトルの経度成分, 緯度成分, 動径成分を表す.
r・(▽×▽×v) = 1/r ∂/∂r (r・( 1/cosφ・∂v[λ]/∂λ
+ 1/cosφ・∂(v[φ] cosφ)/∂φ ) )
+ L^2 v[r]/r
のスペクトルデータが返される.
function wq_RadRotRot_xyr_xyr_xyr(xyr_VLON,xyr_VLAT,xyr_VRAD)
!
! ベクトル v に対して r・(▽×▽×v) を計算する.
!
! 第 1, 2, 3 引数(v[λ], v[φ], v[r])がそれぞれベクトルの経度成分,
! 緯度成分, 動径成分を表す.
!
! r・(▽×▽×v) = 1/r ∂/∂r (r・( 1/cosφ・∂v[λ]/∂λ
! + 1/cosφ・∂(v[φ] cosφ)/∂φ ) )
! + L^2 v[r]/r
!
! のスペクトルデータが返される.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr_VLON
!(in) ベクトルの経度成分
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr_VLAT
!(in) ベクトルの緯度成分
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr_VRAD
!(in) ベクトルの動径成分
real(8), dimension((nm+1)*(nm+1),0:lm) :: wq_RadRotRot_xyr_xyr_xyr
!(out) ベクトル v の r・(▽×▽×v)
wq_RadRotRot_xyr_xyr_xyr = wq_RotDRad_wr( wa_DivLon_xya(xyr_VLON)+ wa_DivLat_xya(xyr_VLAT)) - wa_Lapla_wa(wq_xyr(xyr_VRAD/xyr_RAD))
end function wq_RadRotRot_xyr_xyr_xyr
| Function : | |||
| wq_RotDRad_wr : | real(8), dimension((nm+1)*(nm+1),0:lm)
| ||
| wr : | real(8), dimension((nm+1)*(nm+1),km), intent(in)
|
入力スペクトルデータに回転型動径微分
1/r ∂(r.)/∂r = ∂(.)/∂r + (.)/r
を作用する.
スペクトルデータの回転型動径微分とは, 対応する格子点データに 回転型動径微分を作用させたデータのスペクトル変換のことである.
function wq_RotDRad_wr(wr)
!
! 入力スペクトルデータに回転型動径微分
!
! 1/r ∂(r.)/∂r = ∂(.)/∂r + (.)/r
!
! を作用する.
!
! スペクトルデータの回転型動径微分とは, 対応する格子点データに
! 回転型動径微分を作用させたデータのスペクトル変換のことである.
!
real(8), dimension((nm+1)*(nm+1),km), intent(in) :: wr
!(in) 2 次元球面調和函数チェビシェフスペクトルデータ
real(8), dimension((nm+1)*(nm+1),0:lm) :: wq_RotDRad_wr
!(out) 回転型動径微分を作用された水平スペクトル動径格子点データ
wq_RotDRad_wr = wq_Rad2Inv_wq(wq_RadDRad_wq(wq_wr(wr*wr_Rad)))
end function wq_RotDRad_wr
| Subroutine : | |||
| wq_TORPOT : | real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
| ||
| value : | real(8), dimension((nm+1)*(nm+1)), intent(in), optional
| ||
| cond : | character(len=1), intent(in), optional
| ||
| new : | logical, intent(IN), optional
|
速度トロイダルポテンシャルに対して境界条件を適用する. Chebyshev 空間での境界条件適用.
速度トロイダルポテンシャルΨに対して与えられる境界条件は
* 粘着条件 : Ψ = Ψb(lon,lat). Ψb は境界球面での速度分布.
default は 0(静止状態).
* 応力なし条件 : ∂(Ψ/r)/∂r = 0.
最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.
Alias for wq_TorBoundaryTau
| Subroutine : | |||
| wq_TORPOT : | real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
| ||
| value : | real(8), dimension((nm+1)*(nm+1)), intent(in), optional
| ||
| cond : | character(len=1), intent(in), optional
| ||
| new : | logical, intent(IN), optional
|
速度トロイダルポテンシャルに対して境界条件を適用する. Chebyshev 空間での境界条件適用.
速度トロイダルポテンシャルΨに対して与えられる境界条件は
* 粘着条件 : Ψ = Ψb(lon,lat). Ψb は境界球面での速度分布.
default は 0(静止状態).
* 応力なし条件 : ∂(Ψ/r)/∂r = 0.
最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.
subroutine wq_TorBoundaryTau(wq_TORPOT,value,cond,new)
!
! 速度トロイダルポテンシャルに対して境界条件を適用する.
! Chebyshev 空間での境界条件適用.
!
! 速度トロイダルポテンシャルΨに対して与えられる境界条件は
!
! * 粘着条件 : Ψ = Ψb(lon,lat). Ψb は境界球面での速度分布.
! default は 0(静止状態).
!
! * 応力なし条件 : ∂(Ψ/r)/∂r = 0.
!
! 最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.
!
real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout) :: wq_TORPOT
!(inout) 境界条件を適用するデータ. 修正された値を返す.
real(8), dimension((nm+1)*(nm+1)), intent(in), optional :: value
!(in) 両端境界でのトロイダルポテンシャル
! 粘着条件の時のみ有効
character(len=1), intent(in), optional :: cond
!(in) 境界条件スイッチ. 省略時は 'R'
! R : 上側粘着条件
! F : 上側応力なし条件
logical, intent(IN), optional :: new
!(in) true だと境界条件計算用行列を強制的に新たに作る.
! default は false.
real(8), dimension(:,:,:), allocatable :: alu
integer, dimension(:,:), allocatable :: kp
real(8), dimension((nm+1)*(nm+1),0:lm) :: wq_data
real(8), dimension((nm+1)*(nm+1),km) :: wr_data
logical :: rigid ! 境界条件
logical :: first = .true.
logical :: new_matrix = .false.
integer :: n, l, lend
save :: alu, kp, first
if (.not. present(cond)) then
rigid=.TRUE.
else
select case (cond)
case ('R')
rigid = .TRUE.
case ('F')
rigid = .FALSE.
case default
call MessageNotify('E','wq_TorBoundaryTau','B.C. not supported')
end select
endif
if (.not. present(new)) then
new_matrix=.false.
else
new_matrix=new
endif
if ( first .OR. new_matrix ) then
first = .false.
if ( allocated(alu) ) deallocate(alu)
if ( allocated(kp) ) deallocate(kp)
allocate(alu((nm+1)*(nm+1),0:lm,0:lm),kp((nm+1)*(nm+1),0:lm))
alu = 0.0D0
do l=0,lm
alu(:,l,l)=1.0D0
enddo
! 力学的条件
do l=0,lm
wq_data = 0.0
wq_data(:,l) = 1.0D0
if ( rigid ) then ! 力学的条件粘着
wr_data = wr_wq(wq_data)
else ! 力学的条件自由すべり
wr_data = wr_wq(wq_RadDRad_wq(wq_data)- wq_data)/wr_Rad
endif
do n=1,(nm+1)**2
if ( mod(nd(n),2) .eq. mod(lm,2) ) then
alu(n,lm,l) = wr_data(n,km)
else
alu(n,lm-1,l) = wr_data(n,km)
endif
end do
enddo
! 関係ないところを 0 で埋める.
do n=1,(nm+1)**2
if ( mod(nd(n),2) .eq. mod(lm,2) ) then
lend = lm
else
lend = lm-1
endif
do l=0,nd(n)-1
alu(n,lend,l) = 0.0D0
enddo
do l=nd(n)+1,lm,2
alu(n,lend,l) = 0.0D0
enddo
enddo
call ludecomp(alu,kp)
if ( rigid .AND. present(value) ) then
call MessageNotify('M','wq_TorBoundaryTau', 'Toroidal potential at k=km was given by the optional variable.')
else if ( rigid .AND. (.NOT.present(value)) ) then
call MessageNotify('M','wq_TorBoundaryTau', 'Toroidal potential at k=km was set to zero.')
else if ( (.NOT. rigid) .AND. present(value) ) then
call MessageNotify('W','wq_TorBoundaryTau', 'Boundary value k=km cannot be set under stress-free condition.')
endif
call MessageNotify('M','wq_TorBoundaryTau', 'Matrix to apply b.c. newly produced.')
endif
do n=1,(nm+1)**2
if ( mod(nd(n),2) .eq. mod(lm,2) ) then
lend = lm
else
lend = lm-1
endif
if ( rigid .AND. present(value) ) then
wq_torpot(n,lend) = value(n)
else
wq_torpot(n,lend) = 0.0D0
endif
enddo
wq_torpot = lusolve(alu,kp,wq_TORPOT)
end subroutine wq_TorBoundaryTau
| Subroutine : | |||
| wq_TOR : | real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
| ||
| new : | logical, intent(IN), optional
|
Alias for wq_TorMagBoundaryTau
| Subroutine : | |||
| wq_TOR : | real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
| ||
| new : | logical, intent(IN), optional
|
subroutine wq_TormagBoundaryTau(wq_TOR,new)
! 磁場トロイダルポテンシャルに対して境界条件を適用する.
! Chebyshev 空間での境界条件適用
!
! チェビシェフ空間において境界条件を満たすべく高次の係数を定める方法を
! とっている(タウ法). 現在のところ境界物質が非電気伝導体の場合のみ
! 対応している. その場合, 磁場トロイダルポテンシャルの境界条件は
!
! 外側
! wq_psi = 0 at the outer boundary
!
! であるから wq_Boundary で対応可能だが, 将来のため別途作成しておく.
!
! 最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.
!
real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout) :: wq_TOR
!(inout) 境界条件を適用するデータ. 修正された値を返す.
logical, intent(IN), optional :: new
!(in) true だと境界条件計算用行列を強制的に新たに作る.
! default は false.
real(8), dimension(:,:,:), allocatable :: alu
integer, dimension(:,:), allocatable :: kp
real(8), dimension(:,:), allocatable :: wq_I
real(8), dimension(:,:), allocatable :: wr_PSI
logical :: first = .true.
logical :: new_matrix = .false.
integer :: n, l, lend
save :: alu, kp, first
if (.not. present(new)) then
new_matrix=.false.
else
new_matrix=new
endif
if ( first .OR. new_matrix ) then
first = .false.
if ( allocated(alu) ) deallocate(alu)
if ( allocated(kp) ) deallocate(kp)
if ( allocated(wq_I) ) deallocate(wq_I)
if ( allocated(wr_PSI) ) deallocate(wr_PSI)
allocate(alu((nm+1)*(nm+1),0:lm,0:lm),kp((nm+1)*(nm+1),0:lm))
allocate(wq_I((nm+1)*(nm+1),0:lm),wr_PSI((nm+1)*(nm+1),km))
alu = 0.0D0
do l=0,lm
alu(:,l,l) = 1.0D0
enddo
do l=0,lm
wq_I = 0.0
wq_I(:,l) = 1.0
! 非電気伝導体
wr_PSI = wr_wq(wq_I)
do n=1,(nm+1)**2
if ( mod(nd(n),2) .eq. mod(lm,2) ) then
alu(n,lm,l) = wr_Psi(n,km)
lend = lm
else
alu(n,lm-1,l) = wr_Psi(n,km)
endif
enddo
enddo
! 関係ないところを 0 で埋める.
do n=1,(nm+1)**2
if ( mod(nd(n),2) .eq. mod(lm,2) ) then
lend = lm
else
lend = lm-1
endif
do l=0,nd(n)-1
alu(n,lend,l) = 0.0D0
enddo
do l=nd(n)+1,lm,2
alu(n,lend,l) = 0.0D0
enddo
enddo
call ludecomp(alu,kp)
deallocate(wq_I,wr_PSI)
call MessageNotify('M','TormagBoundaryTau', 'Matrix to apply b.c. newly produced.')
endif
do n=1,(nm+1)**2
if ( mod(nd(n),2) .eq. mod(lm,2) ) then
wq_TOR(n,lm) = 0.0
else
wq_TOR(n,lm-1) = 0.0
endif
enddo
wq_TOR = lusolve(alu,kp,wq_TOR)
end subroutine wq_TormagBoundaryTau
| Function : | |||
| wq_wr : | real(8), dimension((nm+1)*(nm+1),0:lm)
| ||
| wr : | real(8), dimension((nm+1)*(nm+1),km), intent(in)
|
水平スペクトル・動径格子点データからスペクトルデータへ(正)変換する.
function wq_wr(wr)
!
! 水平スペクトル・動径格子点データからスペクトルデータへ(正)変換する.
!
real(8), dimension((nm+1)*(nm+1),0:lm) :: wq_wr
!(out) 2 次元球面調和函数チェビシェフスペクトルデータ
real(8), dimension((nm+1)*(nm+1),km), intent(in) :: wr
!(in) 2 次元球面調和函数スペクトル・動径格子点データ
wq_wr = aq_ar(wr)
end function wq_wr
| Function : | |||
| wq_xyr : | real(8), dimension((nm+1)*(nm+1),0:lm)
| ||
| xyr : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
3 次元格子点データからスペクトルデータへ(正)変換する.
function wq_xyr(xyr)
!
! 3 次元格子点データからスペクトルデータへ(正)変換する.
!
real(8), dimension((nm+1)*(nm+1),0:lm) :: wq_xyr
!(out) 2 次元球面調和函数チェビシェフスペクトルデータ
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr
!(in) 3 次元経度緯度動径格子点データ
wq_xyr = wq_wr(wa_xya(xyr))
end function wq_xyr
| Subroutine : | |||
| wr : | real(8), dimension((nm+1)*(nm+1),km),intent(inout)
| ||
| value : | real(8), dimension((nm+1)*(nm+1)), intent(in), optional
| ||
| cond : | character(len=1), intent(in), optional
|
スペクトルデータにディリクレ・ノイマン境界条件を適用する 実空間での境界条件適用
鉛直実格子点空間において内部領域の値と境界条件を満たすように 条件を課している(選点法).
subroutine wr_BoundaryGrid(wr,value,cond)
!
! スペクトルデータにディリクレ・ノイマン境界条件を適用する
! 実空間での境界条件適用
!
! 鉛直実格子点空間において内部領域の値と境界条件を満たすように
! 条件を課している(選点法).
!
real(8), dimension((nm+1)*(nm+1),km),intent(inout) :: wr
!(inout) 境界条件を適用するデータ. 修正された値を返す.
real(8), dimension((nm+1)*(nm+1)), intent(in), optional :: value
!(in) 境界での 値/勾配 分布を水平スペクトル変換したものを与える.
! 省略時は値/勾配 0 となる.
character(len=1), intent(in), optional :: cond
!(in) 境界条件. 省略時は 'D'
! D : 外側ディリクレ条件
! N : 外側ノイマン条件
if (.not. present(cond)) then
if (present(value)) then
call ag_BoundaryGrid_D(wr,value)
else
call ag_BoundaryGrid_D(wr)
endif
return
endif
select case(cond)
case ('N')
if (present(value)) then
call ag_BoundaryGrid_N(wr,value)
else
call ag_BoundaryGrid_N(wr)
endif
case ('D')
if (present(value)) then
call ag_BoundaryGrid_D(wr,value)
else
call ag_BoundaryGrid_D(wr)
endif
case default
call MessageNotify('E','wr_BoundaryGrid','B.C. not supported')
end select
end subroutine wr_BoundaryGrid
| Function : | |||
| wr_DivLat_xyr : | real(8), dimension((nm+1)*(nm+1),km)
| ||
| xyr : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
格子データに発散型緯度微分 1/rcosφ・∂(f cosφ)/∂φ を 作用させたスペクトルデータを返す.
function wr_DivLat_xyr(xyr)
!
! 格子データに発散型緯度微分 1/rcosφ・∂(f cosφ)/∂φ を
! 作用させたスペクトルデータを返す.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr
!(in) 3 次元経度緯度動径格子点データ
real(8), dimension((nm+1)*(nm+1),km) :: wr_DivLat_xyr
!(out) 発散型緯度微分を作用された水平スペクトル動径格子点データ
wr_DivLat_xyr = wa_DivLat_xya(xyr/xyr_Rad)
end function wr_DivLat_xyr
| Function : | |||
| wr_DivLon_xyr : | real(8), dimension((nm+1)*(nm+1),km)
| ||
| xyr : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
格子点データに発散型経度微分 1/rcosφ・∂/∂λ を作用させた スペクトルデータを返す.
function wr_DivLon_xyr(xyr)
!
! 格子点データに発散型経度微分 1/rcosφ・∂/∂λ を作用させた
! スペクトルデータを返す.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr
!(in) 3 次元経度緯度動径格子点データ
real(8), dimension((nm+1)*(nm+1),km) :: wr_DivLon_xyr
!(out) 発散型経度微分を作用された水平スペクトル動径格子点データ
wr_DivLon_xyr = wa_DivLon_xya(xyr/xyr_Rad)
end function wr_DivLon_xyr
| Function : | |||
| wr_DivRad_wq : | real(8), dimension((nm+1)*(nm+1),km)
| ||
| wq : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
入力スペクトルデータに発散型動径微分
1/r^2 ∂/∂r (r^2 .)= ∂/∂r + 2/r
を作用する.
スペクトルデータの発散型動径微分とは, 対応する格子点データに 発散型動径微分を作用させたデータのスペクトル変換のことである.
function wr_DivRad_wq(wq)
!
! 入力スペクトルデータに発散型動径微分
!
! 1/r^2 ∂/∂r (r^2 .)= ∂/∂r + 2/r
!
! を作用する.
!
! スペクトルデータの発散型動径微分とは, 対応する格子点データに
! 発散型動径微分を作用させたデータのスペクトル変換のことである.
!
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq
!(in) 2 次元球面調和函数チェビシェフスペクトルデータ
real(8), dimension((nm+1)*(nm+1),km) :: wr_DivRad_wq
!(out) 発散型動径微分を作用された水平スペクトル動径格子点データ
wr_DivRad_wq = wr_wq(wq_RadDRad_wq(wq))/wr_Rad + 2/wr_Rad * wr_wq(wq)
end function wr_DivRad_wq
| Function : | |||
| wr_DivRad_xyr : | real(8), dimension((nm+1)*(nm+1),km)
| ||
| xyr : | real(8), dimension(0:im-1,jm,km), intent(in)
|
格子点データに発散型動径微分
1/r^2 ∂/∂r (r^2 = ∂/∂r + 2/= 1/r∂/∂(r.) + 1/r
を作用する.
この関数の入力は動径次数と水平全波数の偶奇性が異なっていることを 仮定している. 偶奇性が一致している場合には wr_DigRad_wq を使用すること.
function wr_DivRad_xyr(xyr)
!
! 格子点データに発散型動径微分
!
! 1/r^2 ∂/∂r (r^2 = ∂/∂r + 2/= 1/r∂/∂(r.) + 1/r
!
! を作用する.
!
! この関数の入力は動径次数と水平全波数の偶奇性が異なっていることを
! 仮定している. 偶奇性が一致している場合には wr_DigRad_wq を使用すること.
!
real(8), dimension(0:im-1,jm,km), intent(in) :: xyr
!(in) 2 次元球面調和函数チェビシェフスペクトルデータ
real(8), dimension((nm+1)*(nm+1),km) :: wr_DivRad_xyr
!(out) 発散型動径微分を作用された水平スペクトル動径格子点データ
wr_DivRad_xyr = wr_wq(wq_RadDRad_wq(wq_xyr(xyr_Rad*xyr)))/wr_Rad**2 + 1/wr_Rad * wr_xyr(xyr)
end function wr_DivRad_xyr
| Function : | |||
| wr_Div_xyr_xyr_xyr : | real(8), dimension((nm+1)*(nm+1),km)
| ||
| xyr_Vlon : | real(8), dimension(0:im-1,1:jm,km), intent(in)
| ||
| xyr_Vlat : | real(8), dimension(0:im-1,1:jm,km), intent(in)
| ||
| xyr_Vrad : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
べクトル成分である 3 つの格子データに発散を作用させた スペクトルデータを返す.
第 1, 2 ,3 引数(u,v,w)がそれぞれベクトルの経度成分, 緯度成分, 動径成分を表し, 発散は
1/rcosφ・∂u/∂λ + 1/rcosφ・∂(v cosφ)/∂φ
+ 1/r^2 ∂/∂r (r^2 w)
と計算される.
function wr_Div_xyr_xyr_xyr(xyr_Vlon,xyr_Vlat,xyr_Vrad)
!
! べクトル成分である 3 つの格子データに発散を作用させた
! スペクトルデータを返す.
!
! 第 1, 2 ,3 引数(u,v,w)がそれぞれベクトルの経度成分, 緯度成分,
! 動径成分を表し, 発散は
!
! 1/rcosφ・∂u/∂λ + 1/rcosφ・∂(v cosφ)/∂φ
! + 1/r^2 ∂/∂r (r^2 w)
!
! と計算される.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr_Vlon
!(in) ベクトル場の経度成分
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr_Vlat
!(in) ベクトル場の緯度成分
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr_Vrad
!(in) ベクトル場の動径成分
real(8), dimension((nm+1)*(nm+1),km) :: wr_Div_xyr_xyr_xyr
!(out) ベクトル場の発散
wr_Div_xyr_xyr_xyr = wr_DivLon_xyr(xyr_Vlon) + wr_DivLat_xyr(xyr_Vlat) + wr_DivRad_xyr(xyr_Vrad)
end function wr_Div_xyr_xyr_xyr
| Subroutine : | |||
| wr_POL : | real(8), dimension((nm+1)*(nm+1),km),intent(inout)
| ||
| new : | logical, intent(IN), optional
|
磁場ポロイダルポテンシャルに対して境界条件を適用する. 鉛直実空間での境界条件適用.
鉛直実格子点空間において内部領域の値と境界条件を満たすように 条件を課している(選点法).
現在のところ境界物質が非電気伝導体の場合のみ対応している. その場合, 磁場ポロイダルポテンシャルの各水平スペクトル成分 h に たいして境界条件が与えられ,
* 外側境界 : dh/dr + (n+1)h/r = 0
である. ここで n は h の水平全波数である.
最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.
subroutine wr_PolmagBoundaryGrid(wr_POL,new)
!
! 磁場ポロイダルポテンシャルに対して境界条件を適用する.
! 鉛直実空間での境界条件適用.
!
! 鉛直実格子点空間において内部領域の値と境界条件を満たすように
! 条件を課している(選点法).
!
! 現在のところ境界物質が非電気伝導体の場合のみ対応している.
! その場合, 磁場ポロイダルポテンシャルの各水平スペクトル成分 h に
! たいして境界条件が与えられ,
!
! * 外側境界 : dh/dr + (n+1)h/r = 0
!
! である. ここで n は h の水平全波数である.
!
! 最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.
!
real(8), dimension((nm+1)*(nm+1),km),intent(inout) :: wr_POL
!(inout) 境界条件を適用するデータ. 修正された値を返す.
logical, intent(IN), optional :: new
!(in) true だと境界条件計算用行列を強制的に新たに作る.
! default は false.
real(8), dimension(:,:,:), allocatable :: alu
integer, dimension(:,:), allocatable :: kp
real(8), dimension(:,:), allocatable :: wr_I
real(8), dimension(:,:), allocatable :: wr_PSI
real(8), dimension(:,:), allocatable :: wr_DPSIDR
logical :: first = .true.
logical :: new_matrix = .false.
integer :: k, n, nn(2)
save :: alu, kp, first
if (.not. present(new)) then
new_matrix=.false.
else
new_matrix=new
endif
if ( first .OR. new_matrix ) then
first = .false.
if ( allocated(alu) ) deallocate(alu)
if ( allocated(kp) ) deallocate(kp)
if ( allocated(wr_I) ) deallocate(wr_I)
if ( allocated(wr_PSI) ) deallocate(wr_PSI)
if ( allocated(wr_DPSIDR) ) deallocate(wr_DPSIDR)
allocate(alu((nm+1)*(nm+1),km,km),kp((nm+1)*(nm+1),km))
allocate(wr_I((nm+1)*(nm+1),km))
allocate(wr_PSI((nm+1)*(nm+1),km),wr_DPSIDR((nm+1)*(nm+1),km))
do k=1,km
wr_I = 0.0D0
wr_I(:,k)=1.0D0
alu(:,:,k) = wr_I ! 内部領域は値そのまま.
enddo
! 非電気伝導体
do k=1,km
wr_I = 0.0D0
wr_I(:,k) = 1.0D0
wr_PSI = wr_I
wr_DPSIDR = wr_wq(wq_RadDRad_wq(wq_wr(wr_I)))/wr_Rad
do n=1,(nm+1)*(nm+1)
nn=nm_l(n)
alu(n,km,k) = wr_DPSIDR(n,km) + (nn(1)+1) * wr_PSI(n,km)/r_RAD(km)
enddo
end do
call ludecomp(alu,kp)
deallocate(wr_I,wr_PSI,wr_DPSIDR)
call MessageNotify('M','PolmagBoundaryGrid', 'Matrix to apply b.c. newly produced.')
endif
wr_POL(:,km) = 0.0D0
wr_POL = lusolve(alu,kp,wr_POL)
end subroutine wr_PolmagBoundaryGrid
| Function : | |||
| wr_RadRotRot_xyr_xyr_xyr : | real(8), dimension((nm+1)*(nm+1),km)
| ||
| xyr_VLON : | real(8), dimension(0:im-1,1:jm,km), intent(in)
| ||
| xyr_VLAT : | real(8), dimension(0:im-1,1:jm,km), intent(in)
| ||
| xyr_VRAD : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
ベクトル v に対して r・(▽×▽×v) を計算する.
第 1, 2, 3 引数(v[λ], v[φ], v[r])がそれぞれベクトルの経度成分, 緯度成分, 動径成分を表す.
r・(▽×▽×v) = 1/r ∂/∂r (r・( 1/cosφ・∂v[λ]/∂λ
+ 1/cosφ・∂(v[φ] cosφ)/∂φ ) )
+ L^2 v[r]/r
のスペクトルデータが返される.
function wr_RadRotRot_xyr_xyr_xyr(xyr_VLON,xyr_VLAT,xyr_VRAD)
!
! ベクトル v に対して r・(▽×▽×v) を計算する.
!
! 第 1, 2, 3 引数(v[λ], v[φ], v[r])がそれぞれベクトルの経度成分,
! 緯度成分, 動径成分を表す.
!
! r・(▽×▽×v) = 1/r ∂/∂r (r・( 1/cosφ・∂v[λ]/∂λ
! + 1/cosφ・∂(v[φ] cosφ)/∂φ ) )
! + L^2 v[r]/r
!
! のスペクトルデータが返される.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr_VLON
!(in) ベクトルの経度成分
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr_VLAT
!(in) ベクトルの緯度成分
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr_VRAD
!(in) ベクトルの動径成分
real(8), dimension((nm+1)*(nm+1),km) :: wr_RadRotRot_xyr_xyr_xyr
!(out) ベクトル v の r・(▽×▽×v)
wr_RadRotRot_xyr_xyr_xyr = wr_RotDRad_wr( wa_DivLon_xya(xyr_VLON)+ wa_DivLat_xya(xyr_VLAT)) - wa_Lapla_wa(wr_xyr(xyr_VRAD/xyr_RAD))
end function wr_RadRotRot_xyr_xyr_xyr
| Function : | |||
| wr_RadRot_xyr_xyr : | real(8), dimension((nm+1)*(nm+1),km)
| ||
| xyr_VLON : | real(8), dimension(0:im-1,1:jm,km), intent(in)
| ||
| xyr_VLAT : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
r・(▽×v)
ベクトルの渦度と動径ベクトルの内積 r・(▽×v) を計算する.
第 1, 2 引数(v[λ], v[φ])がそれぞれベクトルの経度成分, 緯度成分を表す.
r・(▽×v) = 1/cosφ・∂v[φ]/∂λ - 1/cosφ・∂(v[λ] cosφ)/∂φ
のスペクトル データが返される.
function wr_RadRot_xyr_xyr(xyr_VLON,xyr_VLAT) ! r・(▽×v)
!
! ベクトルの渦度と動径ベクトルの内積 r・(▽×v) を計算する.
!
! 第 1, 2 引数(v[λ], v[φ])がそれぞれベクトルの経度成分, 緯度成分を表す.
!
! r・(▽×v) = 1/cosφ・∂v[φ]/∂λ - 1/cosφ・∂(v[λ] cosφ)/∂φ
!
! のスペクトル データが返される.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr_VLON
!(in) ベクトルの経度成分
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr_VLAT
!(in) ベクトルの緯度成分
real(8), dimension((nm+1)*(nm+1),km) :: wr_RadRot_xyr_xyr
!(out) ベクトルの渦度と動径ベクトルの内積
wr_RadRot_xyr_xyr = wa_DivLon_xya(xyr_VLAT) - wa_DivLat_xya(xyr_VLON)
end function wr_RadRot_xyr_xyr
| Function : | |||
| wr_RotDRad_wq : | real(8), dimension((nm+1)*(nm+1),km)
| ||
| wq : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
入力スペクトルデータに回転型動径微分
1/r ∂(r.)/∂r = ∂(.)/∂r + (.)/r
を作用する.
スペクトルデータの回転型動径微分とは, 対応する格子点データに 回転型動径微分を作用させたデータのスペクトル変換のことである.
function wr_RotDRad_wq(wq)
!
! 入力スペクトルデータに回転型動径微分
!
! 1/r ∂(r.)/∂r = ∂(.)/∂r + (.)/r
!
! を作用する.
!
! スペクトルデータの回転型動径微分とは, 対応する格子点データに
! 回転型動径微分を作用させたデータのスペクトル変換のことである.
!
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq
!(in) 2 次元球面調和函数チェビシェフスペクトルデータ
real(8), dimension((nm+1)*(nm+1),km) :: wr_RotDRad_wq
!(out) 回転型動径微分を作用された水平スペクトル動径格子点データ
wr_RotDRad_wq = wr_wq(wq_RadDrad_wq(wq))/wr_Rad + wr_wq(wq)/wr_Rad
end function wr_RotDRad_wq
| Function : | |||
| wr_RotDRad_wr : | real(8), dimension((nm+1)*(nm+1),km)
| ||
| wr : | real(8), dimension((nm+1)*(nm+1),km), intent(in)
|
入力スペクトルデータに回転型動径微分
1/r ∂(r.)/∂r = ∂(.)/∂r + (.)/r
を作用する.
スペクトルデータの回転型動径微分とは, 対応する格子点データに 回転型動径微分を作用させたデータのスペクトル変換のことである.
function wr_RotDRad_wr(wr)
!
! 入力スペクトルデータに回転型動径微分
!
! 1/r ∂(r.)/∂r = ∂(.)/∂r + (.)/r
!
! を作用する.
!
! スペクトルデータの回転型動径微分とは, 対応する格子点データに
! 回転型動径微分を作用させたデータのスペクトル変換のことである.
!
real(8), dimension((nm+1)*(nm+1),km), intent(in) :: wr
!(in) 2 次元球面調和函数チェビシェフスペクトルデータ
real(8), dimension((nm+1)*(nm+1),km) :: wr_RotDRad_wr
!(out) 回転型動径微分を作用された水平スペクトル動径格子点データ
wr_RotDRad_wr = wr_wq(wq_RadDRad_wq(wq_wr(wr*wr_Rad)))/wr_Rad**2
end function wr_RotDRad_wr
| Function : | |||
| wr_RotRad_xyr_xyr : | real(8), dimension((nm+1)*(nm+1),km)
| ||
| xyr_Vlat : | real(8), dimension(0:im-1,1:jm,km), intent(in)
| ||
| xyr_Vlon : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
ベクトルの緯度成分, 経度成分である第 1, 2 引数 Vlat, Vlon に対して ベクトル場の回転の動径成分
1/rcosφ・∂Vlat/∂λ - 1/rcosφ・∂(Vlon cosφ)/∂φ
を計算する.
function wr_RotRad_xyr_xyr(xyr_Vlat,xyr_Vlon)
!
! ベクトルの緯度成分, 経度成分である第 1, 2 引数 Vlat, Vlon に対して
! ベクトル場の回転の動径成分
!
! 1/rcosφ・∂Vlat/∂λ - 1/rcosφ・∂(Vlon cosφ)/∂φ
!
! を計算する.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr_Vlat
!(in) ベクトル場の緯度成分
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr_Vlon
!(in) ベクトル場の経度成分
real(8), dimension((nm+1)*(nm+1),km) :: wr_RotRad_xyr_xyr
!(out) ベクトル場の回転の動径成分
wr_RotRad_xyr_xyr = wr_DivLon_xyr(xyr_Vlat) - wr_DivLat_xyr(xyr_Vlon)
end function wr_RotRad_xyr_xyr
| Subroutine : | |||
| wr_TORPOT : | real(8), dimension((nm+1)*(nm+1),km),intent(inout)
| ||
| value : | real(8), dimension((nm+1)*(nm+1)), intent(in), optional
| ||
| cond : | character(len=1), intent(in), optional
| ||
| new : | logical, intent(IN), optional
|
速度トロイダルポテンシャルに対して境界条件を適用する. 実空間での境界条件適用
鉛直実格子点空間において内部領域の値と境界条件を満たすように 条件を課している(選点法).
速度トロイダルポテンシャルΨに対して与えられる境界条件は
* 粘着条件 : Ψ = Ψb(lon,lat). Ψb は境界球面での速度分布.
default は 0 (静止状態).
* 応力なし条件 : ∂(Ψ/r)/∂r = 0.
最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.
subroutine wr_TorBoundaryGrid(wr_TORPOT,value,cond,new)
!
! 速度トロイダルポテンシャルに対して境界条件を適用する.
! 実空間での境界条件適用
!
! 鉛直実格子点空間において内部領域の値と境界条件を満たすように
! 条件を課している(選点法).
!
! 速度トロイダルポテンシャルΨに対して与えられる境界条件は
!
! * 粘着条件 : Ψ = Ψb(lon,lat). Ψb は境界球面での速度分布.
! default は 0 (静止状態).
!
! * 応力なし条件 : ∂(Ψ/r)/∂r = 0.
!
! 最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.
!
real(8), dimension((nm+1)*(nm+1),km),intent(inout) :: wr_TORPOT
!(inout) 境界条件を適用するデータ. 修正された値を返す.
real(8), dimension((nm+1)*(nm+1)), intent(in), optional :: value
!(in) 両端境界でのトロイダルポテンシャル
! 粘着条件の時のみ有効
character(len=1), intent(in), optional :: cond
!(in) 境界条件スイッチ. 省略時は 'R'
! R : 上側粘着条件
! F : 上側応力なし条件
logical, intent(IN), optional :: new
!(in) true だと境界条件計算用行列を強制的に新たに作る.
! default は false.
real(8), dimension(:,:,:), allocatable :: alu
integer, dimension(:,:), allocatable :: kp
real(8), dimension((nm+1)*(nm+1),0:lm) :: wq_data
real(8), dimension((nm+1)*(nm+1),km) :: wr_data
logical :: rigid ! 境界条件
logical :: first = .true.
logical :: new_matrix = .false.
integer :: k
save :: alu, kp, first
if (.not. present(cond)) then
rigid=.TRUE.
else
select case (cond)
case ('R')
rigid = .TRUE.
case ('F')
rigid = .FALSE.
case default
call MessageNotify('E','wr_TorBoundaryGrid','B.C. not supported')
end select
endif
if (.not. present(new)) then
new_matrix=.false.
else
new_matrix=new
endif
if ( first .OR. new_matrix ) then
first = .false.
if ( allocated(alu) ) deallocate(alu)
if ( allocated(kp) ) deallocate(kp)
allocate(alu((nm+1)*(nm+1),km,km),kp((nm+1)*(nm+1),km))
alu = 0.0D0
do k=1,km
wr_data = 0.0D0
wr_data(:,k)=1.0D0
alu(:,:,k) = wr_data
enddo
if ( .not. rigid ) then
do k=1,km
wr_data = 0.0D0
wr_data(:,k)=1.0D0
wq_data = wq_wr(wr_data)
wr_data = wr_wq(wq_RadDRad_wq(wq_data) - wq_data)/wr_Rad
alu(:,km,k) = wr_data(:,km)
enddo
endif
call ludecomp(alu,kp)
if ( rigid .AND. present(value) ) then
call MessageNotify('M','wr_TorBoundaryGrid', 'Toroidal potential at k=km was given by the optional variable.')
else if ( rigid .AND. (.NOT.present(value)) ) then
call MessageNotify('M','wr_TorBoundaryGrid', 'Toroidal potential at k=km was set to zero.')
else if ( (.NOT. rigid) .AND. present(value) ) then
call MessageNotify('W','wr_TorBoundaryGrid', 'Boundary value at k=km cannot be set under stress-free condition.')
endif
call MessageNotify('M','wr_TorBoundaryGrid', 'Matrix to apply b.c. newly produced.')
endif
if ( rigid .AND. present(value) ) then
wr_TorPot(:,km) = value
else
wr_TorPot(:,km) = 0.0D0
endif
wr_TorPot = lusolve(alu,kp,wr_TorPot)
end subroutine wr_TorBoundaryGrid
| Subroutine : | |||
| wr_TOR : | real(8), dimension((nm+1)*(nm+1),km),intent(inout)
| ||
| new : | logical, intent(IN), optional
|
磁場トロイダルポテンシャルに対して境界条件を適用する. 鉛直実空間での境界条件適用.
鉛直実格子点空間において内部領域の値と境界条件を満たすように 条件を課している(選点法).
現在のところ境界物質が非電気伝導体の場合のみ対応している. その場合, 磁場トロイダルポテンシャルの境界条件は
外側
wq_psi = 0 at the outer boundary
であるので wq_Boundary で対応可能だが, 将来のため別途作成しておく
最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.
subroutine wr_TormagBoundaryGrid(wr_TOR,new)
!
! 磁場トロイダルポテンシャルに対して境界条件を適用する.
! 鉛直実空間での境界条件適用.
!
! 鉛直実格子点空間において内部領域の値と境界条件を満たすように
! 条件を課している(選点法).
!
! 現在のところ境界物質が非電気伝導体の場合のみ対応している.
! その場合, 磁場トロイダルポテンシャルの境界条件は
!
! 外側
! wq_psi = 0 at the outer boundary
!
! であるので wq_Boundary で対応可能だが, 将来のため別途作成しておく
!
! 最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.
!
real(8), dimension((nm+1)*(nm+1),km),intent(inout) :: wr_TOR
!(inout) 境界条件を適用するデータ. 修正された値を返す.
logical, intent(IN), optional :: new
!(in) (ダミー) true だと境界条件計算用行列を強制的に新たに作る.
! default は false.
wr_TOR(:,km) = 0.0D0
end subroutine wr_TormagBoundaryGrid
| Function : | |||
| wr_wq : | real(8), dimension((nm+1)*(nm+1),km)
| ||
| wq : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
スペクトルデータから水平スペクトル・動径格子点データへ(正)変換する.
function wr_wq(wq)
!
! スペクトルデータから水平スペクトル・動径格子点データへ(正)変換する.
!
real(8), dimension((nm+1)*(nm+1),km) :: wr_wq
!(out) 2 次元球面調和函数スペクトル・動径格子点データ
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq
!(in) 2 次元球面調和函数チェビシェフスペクトルデータ
wr_wq = ar_aq(wq)
end function wr_wq
| Function : | |||
| wr_xyr : | real(8), dimension((nm+1)*(nm+1),km)
| ||
| xyr : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
3 次元格子データから水平スペクトル・動径格子点データへ(正)変換する.
function wr_xyr(xyr)
!
! 3 次元格子データから水平スペクトル・動径格子点データへ(正)変換する.
!
real(8), dimension((nm+1)*(nm+1),km) :: wr_xyr
!(out) 2 次元球面調和函数スペクトル・動径格子点データ
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr
!(in) 3 次元経度緯度動径格子点データ
wr_xyr = wa_xya(xyr)
end function wr_xyr
| Function : | |||
| x_AvrLatRad_xyr : | real(8), dimension(0:im-1)
| ||
| xyr : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
緯度動径(子午面)積分
3 次元格子点データの緯度動径(子午面)平均
3 次元データ f(λ,φ,r) に対して
∫f(λ,,r) r^2cosφ dφdr /(2(r[o]^3-r[i]^3)/3)
を計算する.
function x_AvrLatRad_xyr(xyr) ! 緯度動径(子午面)積分
!
! 3 次元格子点データの緯度動径(子午面)平均
!
! 3 次元データ f(λ,φ,r) に対して
!
! ∫f(λ,,r) r^2cosφ dφdr /(2(r[o]^3-r[i]^3)/3)
!
! を計算する.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr
!(in) 3 次元経度緯度動径格子点データ
real(8), dimension(0:im-1) :: x_AvrLatRad_xyr
!(out) 緯度動径(子午面)平均された 1 次元経度格子点データ
x_AvrLatRad_xyr = x_IntLatRad_xyr(xyr) /( sum(y_Lat_Weight)*sum(r_Rad_Weight) )
end function x_AvrLatRad_xyr
| Function : | |||
| x_AvrLat_xy(im) : | real(8)
| ||
| xy_data(0:im-1,1:jm) : | real(8), intent(in)
|
2 次元緯度経度格子点データの緯度(Y)方向平均(1 層用).
実際には格子点データ各点毎に y_Y_Weight をかけた総和を計算し, y_Y_Weight の総和で割ることで平均している.
Original external subprogram is wa_module_sjpack#x_AvrLat_xy
| Function : | |||
| x_AvrRad_xr : | real(8), dimension(0:im-1)
| ||
| xr : | real(8), dimension(0:im-1,km), intent(in)
|
動径積分
2 次元(XR)格子点データの動径方向域平均.
2 次元データ f(λ,r) に対して
∫f(λ,r) r^2dr /((r[o]^3-r[i]^3)/3)
を計算する.
function x_AvrRad_xr(xr) ! 動径積分
!
! 2 次元(XR)格子点データの動径方向域平均.
!
! 2 次元データ f(λ,r) に対して
!
! ∫f(λ,r) r^2dr /((r[o]^3-r[i]^3)/3)
!
! を計算する.
!
real(8), dimension(0:im-1,km), intent(in) :: xr
!(in) 2 次元緯度動径格子点データ
real(8), dimension(0:im-1) :: x_AvrRad_xr
!(out) 動径平均された 1 次元経度格子点データ
x_AvrRad_xr = x_IntRad_xr(xr)/sum(r_Rad_Weight)
end function x_AvrRad_xr
| Function : | |||
| x_IntLatRad_xyr : | real(8), dimension(0:im-1)
| ||
| xyr : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
3 次元格子点データの緯度動径(子午面)積分
3 次元データ f(λ,φ,r) に対して
∫f(λ,φ,r) r^2cosφ dφdr
を計算する.
function x_IntLatRad_xyr(xyr)
!
! 3 次元格子点データの緯度動径(子午面)積分
!
! 3 次元データ f(λ,φ,r) に対して
!
! ∫f(λ,φ,r) r^2cosφ dφdr
!
! を計算する.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr
!(in) 3 次元経度緯度動径格子点データ
real(8), dimension(0:im-1) :: x_IntLatRad_xyr
!(out) 緯度動径(子午面)積分された 1 次元経度格子点データ
integer :: j, k
x_IntLatRad_xyr = 0.0D0
do k=1,km
do j=1,jm
x_IntLatRad_xyr = x_IntLatRad_xyr + xyr(:,j,k) * y_Lat_Weight(j) * r_Rad_Weight(k)
enddo
enddo
end function x_IntLatRad_xyr
| Function : | |||
| x_IntLat_xy(0:im-1) : | real(8)
| ||
| xy_data(0:im-1,1:jm) : | real(8), intent(in)
|
2 次元緯度経度格子点データの緯度(Y)方向積分(1 層用).
実際には格子点データ各点毎に y_Y_Weight をかけた総和を計算している.
Original external subprogram is wa_module_sjpack#x_IntLat_xy
| Function : | |||
| x_IntRad_xr : | real(8), dimension(0:im-1)
| ||
| xr : | real(8), dimension(0:im-1,km), intent(in)
|
2 次元(XR)格子点データの動径方向域積分.
2 次元データ f(λ,r) に対して ∫f(λ,r) r^2dr を計算する.
function x_IntRad_xr(xr)
!
! 2 次元(XR)格子点データの動径方向域積分.
!
! 2 次元データ f(λ,r) に対して ∫f(λ,r) r^2dr を計算する.
!
real(8), dimension(0:im-1,km), intent(in) :: xr
!(in) 2 次元緯度動径格子点データ
real(8), dimension(0:im-1) :: x_IntRad_xr
!(out) 動径積分された 1 次元経度格子点データ
integer :: k
x_IntRad_xr = 0.0d0
do k=1,km
x_IntRad_xr(:) = x_IntRad_xr(:) + xr(:,k) * r_Rad_Weight(k)
enddo
end function x_IntRad_xr
| Variable : | |||
| x_Lon(:) : | real(8), allocatable
|
Original external subprogram is wa_module_sjpack#x_Lon
| Variable : | |||
| x_Lon_Weight(:) : | real(8), allocatable
|
Original external subprogram is wa_module_sjpack#x_Lon_Weight
| Function : | |||
| xr_AvrLat_xyr : | real(8), dimension(0:im-1,km)
| ||
| xyr : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
緯度積分
3 次元格子点データの緯度方向域平均.
3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r)cosφ dφ/2 を計算する.
function xr_AvrLat_xyr(xyr) ! 緯度積分
!
! 3 次元格子点データの緯度方向域平均.
!
! 3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r)cosφ dφ/2 を計算する.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr
!(in) 3 次元経度緯度動径格子点データ
real(8), dimension(0:im-1,km) :: xr_AvrLat_xyr
!(out) 緯度平均された 2 次元緯度動径格子点データ
xr_AvrLat_xyr = xr_IntLat_xyr(xyr)/sum(y_Lat_Weight)
end function xr_AvrLat_xyr
| Function : | |||
| xr_IntLat_xyr : | real(8), dimension(0:im-1,km)
| ||
| xyr : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
3 次元格子点データの緯度方向域積分.
3 次元データ f(λ,φ,r) に対して∫f(λ,φ,r) cosφ dφ を計算する.
function xr_IntLat_xyr(xyr)
!
! 3 次元格子点データの緯度方向域積分.
!
! 3 次元データ f(λ,φ,r) に対して∫f(λ,φ,r) cosφ dφ を計算する.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr
!(in) 3 次元経度緯度動径格子点データ
real(8), dimension(0:im-1,km) :: xr_IntLat_xyr
!(out) 緯度積分された 2 次元緯度動径格子点データ
! 緯度円格子点データ
integer :: j
xr_IntLat_xyr = 0.0d0
do j=1,jm
xr_IntLat_xyr(:,:) = xr_IntLat_xyr(:,:) + xyr(:,j,:) * y_Lat_Weight(j)
enddo
end function xr_IntLat_xyr
| Function : | |||
| xy_AvrRad_xyr : | real(8), dimension(0:im-1,1:jm)
| ||
| xyr : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
3 次元格子点データの動径方向域平均.
3 次元データ f(λ,φ,r) に対して
∫f(λ,φ,r) r^2dr/((r[o]^3-r[i]^3)/3)
を計算する.
function xy_AvrRad_xyr(xyr)
!
! 3 次元格子点データの動径方向域平均.
!
! 3 次元データ f(λ,φ,r) に対して
!
! ∫f(λ,φ,r) r^2dr/((r[o]^3-r[i]^3)/3)
!
! を計算する.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr
!(in) 3 次元経度緯度動径格子点データ
real(8), dimension(0:im-1,1:jm) :: xy_AvrRad_xyr
!(out) 動径平均された 2 次元経度緯度(水平, 球面)格子点データ
! 水平格子点データ
xy_AvrRad_xyr = xy_IntRad_xyr(xyr)/sum(r_Rad_Weight)
end function xy_AvrRad_xyr
| Function : | |||
| xy_IntRad_xyr : | real(8), dimension(0:im-1,1:jm)
| ||
| xyr : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
動径積分
3 次元格子点データの動径方向域積分.
3 次元データ f(λ,φ,r) に対して∫f(λ,φ,r) r^2dr を計算する.
function xy_IntRad_xyr(xyr) ! 動径積分
!
! 3 次元格子点データの動径方向域積分.
!
! 3 次元データ f(λ,φ,r) に対して∫f(λ,φ,r) r^2dr を計算する.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr
!(in) 3 次元経度緯度動径格子点データ
real(8), dimension(0:im-1,1:jm) :: xy_IntRad_xyr
!(out) 動径積分された 2 次元経度緯度(水平, 球面)格子点データ
integer :: k
xy_IntRad_xyr = 0.0d0
do k=1,km
xy_IntRad_xyr(:,:) = xy_IntRad_xyr(:,:) + xyr(:,:,k) * r_Rad_Weight(k)
enddo
end function xy_IntRad_xyr
| Variable : | |
| xy_Lat(:,:) : | real(8), allocatable |
Original external subprogram is wa_module_sjpack#xy_Lat
| Variable : | |
| xy_Lon(:,:) : | real(8), allocatable |
Original external subprogram is wa_module_sjpack#xy_Lon
| Function : | |||
| xy_w(0:im-1,1:jm) : | real(8)
| ||
| w_data((mm+1)*(mm+1)) : | real(8), intent(in)
| ||
| ipow : | integer, intent(in), optional
| ||
| iflag : | integer, intent(in), optional
|
スペクトルデータから格子データへ変換する(1 層用).
Original external subprogram is wa_module_sjpack#xy_w
| Function : | |||
| xyr_Div_xyr_xyr_xyr : | real(8), dimension(0:im-1,1:jm,km)
| ||
| xyr_Vlon : | real(8), dimension(0:im-1,1:jm,km), intent(in)
| ||
| xyr_Vlat : | real(8), dimension(0:im-1,1:jm,km), intent(in)
| ||
| xyr_Vrad : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
ベクトル成分である 3 つの格子データに発散を作用させる.
第 1, 2 ,3 引数(u,v,w)がそれぞれベクトルの経度成分, 緯度成分, 動径成分を表す.
極の特異性を回避するためにベクトル場に cosφ/r の重みをかけて 計算している.
div V = (r/cosφ)・div (Vcosφ/r) + V_φtanφ/r + V_r/r
function xyr_Div_xyr_xyr_xyr(xyr_Vlon,xyr_Vlat,xyr_Vrad)
!
! ベクトル成分である 3 つの格子データに発散を作用させる.
!
! 第 1, 2 ,3 引数(u,v,w)がそれぞれベクトルの経度成分, 緯度成分,
! 動径成分を表す.
!
! 極の特異性を回避するためにベクトル場に cosφ/r の重みをかけて
! 計算している.
!
! div V = (r/cosφ)・div (Vcosφ/r) + V_φtanφ/r + V_r/r
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr_Vlon
!(in) ベクトル場の経度成分
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr_Vlat
!(in) ベクトル場の緯度成分
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr_Vrad
!(in) ベクトル場の動径成分
real(8), dimension(0:im-1,1:jm,km) :: xyr_Div_xyr_xyr_xyr
!(out) ベクトル場の発散
xyr_Div_xyr_xyr_xyr = xyr_Rad/cos(xyr_Lat) * xyr_wr(wr_Div_xyr_xyr_xyr(xyr_VLon*cos(xyr_Lat)/xyr_Rad, xyr_VLat*cos(xyr_Lat)/xyr_Rad, xyr_VRad*cos(xyr_Lat)/xyr_Rad )) + xyr_VLat*tan(xyr_Lat)/xyr_Rad + xyr_VRad/xyr_Rad
end function xyr_Div_xyr_xyr_xyr
| Function : | |||
| xyr_GradLat_wq : | real(8), dimension(0:im-1,1:jm,km)
| ||
| wq : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
スペクトルデータに勾配型経度微分 1/r ∂/∂φ を作用させる.
function xyr_GradLat_wq(wq)
!
! スペクトルデータに勾配型経度微分 1/r ∂/∂φ を作用させる.
!
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq
!(in) 2 次元球面調和函数チェビシェフスペクトルデータ
real(8), dimension(0:im-1,1:jm,km) :: xyr_GradLat_wq
!(out) 勾配型緯度微分を作用された 2 次元スペクトルデータ
xyr_GradLat_wq = xya_GradLat_wa(wr_wq(wq))/xyr_Rad
end function xyr_GradLat_wq
| Function : | |||
| xyr_GradLon_wq : | real(8), dimension(0:im-1,1:jm,km)
| ||
| wq : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
スペクトルデータに勾配型経度微分 1/rcosφ・∂/∂λ を作用させる.
function xyr_GradLon_wq(wq)
!
! スペクトルデータに勾配型経度微分 1/rcosφ・∂/∂λ
! を作用させる.
!
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq
!(in) 2 次元球面調和函数チェビシェフスペクトルデータ
real(8), dimension(0:im-1,1:jm,km) :: xyr_GradLon_wq
!(out) 勾配型経度微分を作用された 2 次元スペクトルデータ
xyr_GradLon_wq = xya_GradLon_wa(wr_wq(wq))/xyr_Rad
end function xyr_GradLon_wq
| Function : | |||
| xyr_KGrad_wq : | real(8), dimension(0:im-1,1:jm,km)
| ||
| wq : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
k・▽ = cosφ/r ∂/∂φ + sinφ∂/∂r
入力スペクトルデータに対応する格子データに軸方向微分
k・▽ = cosφ/r ∂/∂φ + sinφ∂/∂r
を作用させた格子データが返される. ここでベクトル k は球の中心から北極向きの単位ベクトルである.
function xyr_KGrad_wq(wq) ! k・▽ = cosφ/r ∂/∂φ + sinφ∂/∂r
!
! 入力スペクトルデータに対応する格子データに軸方向微分
!
! k・▽ = cosφ/r ∂/∂φ + sinφ∂/∂r
!
! を作用させた格子データが返される.
! ここでベクトル k は球の中心から北極向きの単位ベクトルである.
!
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq
!(in) 2 次元球面調和函数チェビシェフスペクトルデータ
real(8), dimension(0:im-1,1:jm,km) :: xyr_KGrad_wq
!(out) 軸方向微分を作用された 2 次元スペクトルデータ
xyr_KGrad_wq = cos(xyr_Lat)*xyr_GradLat_wq(wq) + sin(xyr_Lat)*xyr_wq(wq_RadDRad_wq(wq))/xyr_Rad
end function xyr_KGrad_wq
| Function : | |||
| xyr_RotLat_wq_wq : | real(8), dimension(0:im-1,1:jm,km)
| ||
| wq_Vlon : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
| ||
| wq_Vrad : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
ベクトル場の経度成分, 動径成分である第 1, 2 引数 Vlon, Vrad から 回転の緯度成分
1/r ∂(r Vlon)/∂r - 1/rcosφ・∂Vrad/∂λ
を計算する.
function xyr_RotLat_wq_wq(wq_Vlon,wq_Vrad)
!
! ベクトル場の経度成分, 動径成分である第 1, 2 引数 Vlon, Vrad から
! 回転の緯度成分
!
! 1/r ∂(r Vlon)/∂r - 1/rcosφ・∂Vrad/∂λ
!
! を計算する.
!
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq_Vlon
!(in) ベクトル場の経度成分
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq_Vrad
!(in) ベクトル場の動径成分
real(8), dimension(0:im-1,1:jm,km) :: xyr_RotLat_wq_wq
!(out) ベクトル場の回転の緯度成分
xyr_RotLat_wq_wq = xyr_wr(wr_RotDRad_wq(wq_Vlon)) - xyr_GradLon_wq(wq_Vrad)
end function xyr_RotLat_wq_wq
| Function : | |||
| xyr_RotLon_wq_wq : | real(8), dimension(0:im-1,1:jm,km)
| ||
| wq_Vrad : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
| ||
| wq_Vlat : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
ベクトル場の動径成分, 緯度成分である第 1, 2 引数 Vrad, Vlat から 回転の経度成分
1/r ∂Vrad/∂φ-1/r ∂(r Vlat)/∂r を計算する.
を計算する
function xyr_RotLon_wq_wq(wq_Vrad,wq_Vlat)
!
! ベクトル場の動径成分, 緯度成分である第 1, 2 引数 Vrad, Vlat から
! 回転の経度成分
!
! 1/r ∂Vrad/∂φ-1/r ∂(r Vlat)/∂r を計算する.
!
! を計算する
!
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq_Vrad
!(in) ベクトル場の動径成分
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq_Vlat
!(in) ベクトル場の緯度成分
real(8), dimension(0:im-1,1:jm,km) :: xyr_RotLon_wq_wq
!(out) ベクトル場の回転の経度成分
xyr_RotLon_wq_wq = xyr_GradLat_wq(wq_Vrad) - xyr_wr(wr_RotDRad_wq(wq_Vlat))
end function xyr_RotLon_wq_wq
| Function : | |||
| xyr_wq : | real(8), dimension(0:im-1,1:jm,km)
| ||
| wq : | real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
|
スペクトルデータから 3 次元格子点データへ(逆)変換する.
function xyr_wq(wq)
!
! スペクトルデータから 3 次元格子点データへ(逆)変換する.
!
real(8), dimension(0:im-1,1:jm,km) :: xyr_wq
!(out) 3 次元経度緯度動径格子点データ
real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wq
!(in) 2 次元球面調和函数チェビシェフスペクトルデータ
xyr_wq = xya_wa(wr_wq(wq))
end function xyr_wq
| Function : | |||
| xyr_wr : | real(8), dimension(0:im-1,1:jm,km)
| ||
| wr : | real(8), dimension((nm+1)*(nm+1),km), intent(in)
|
水平スペクトル・動径格子点データから 3 次元格子点データへ(逆)変換する.
function xyr_wr(wr)
!
! 水平スペクトル・動径格子点データから 3 次元格子点データへ(逆)変換する.
!
real(8), dimension(0:im-1,1:jm,km) :: xyr_wr
!(out) 3 次元経度緯度動径格子点データ
real(8), dimension((nm+1)*(nm+1),km), intent(in) :: wr
!(in) 2 次元球面調和函数スペクトル・動径格子点データ
xyr_wr = xya_wa(wr)
end function xyr_wr
| Function : | |||
| y_AvrLonRad_xyr : | real(8), dimension(1:jm)
| ||
| xyr : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
経度動径(緯度円)積分
3 次元格子点データの経度動径(緯度円)平均.
3 次元データ f(λ,φ,r) に対して
∫f(λ,φ,r) r^2dλdr /(2π(r[o]^3-r[i]^3)/3)
を計算する.
function y_AvrLonRad_xyr(xyr) ! 経度動径(緯度円)積分
!
! 3 次元格子点データの経度動径(緯度円)平均.
!
! 3 次元データ f(λ,φ,r) に対して
!
! ∫f(λ,φ,r) r^2dλdr /(2π(r[o]^3-r[i]^3)/3)
!
! を計算する.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr
!(in) 3 次元経度緯度動径格子点データ
real(8), dimension(1:jm) :: y_AvrLonRad_xyr
!(out) 経度動径(緯度円)平均された 1 次元緯度格子点データ
y_AvrLonRad_xyr = y_IntLonRad_xyr(xyr) /(sum(x_Lon_Weight)*sum(r_Rad_Weight))
end function y_AvrLonRad_xyr
| Function : | |||
| y_AvrLon_xy(1:jm) : | real(8)
| ||
| xy_data(0:im-1,1:jm) : | real(8), intent(in)
|
2 次元緯度経度格子点データの経度(X)方向平均(1 層用).
実際には格子点データ各点毎に x_X_Weight をかけた総和を計算し, x_X_Weight の総和で割ることで平均している.
Original external subprogram is wa_module_sjpack#y_AvrLon_xy
| Function : | |||
| y_AvrRad_yr : | real(8), dimension(1:jm)
| ||
| yr : | real(8), dimension(1:jm,km), intent(in)
|
2 次元(YR)格子点データの動径方向域平均.
2 次元データ f(φ,r) に対して ∫f(φ,r) r^2dr /((r[o]^3-r[i]^3)/3) を計算する.
function y_AvrRad_yr(yr)
!
! 2 次元(YR)格子点データの動径方向域平均.
!
! 2 次元データ f(φ,r) に対して ∫f(φ,r) r^2dr /((r[o]^3-r[i]^3)/3)
! を計算する.
!
real(8), dimension(1:jm,km), intent(in) :: yr
!(in) 2 次元緯度動径(子午面)格子点データ
real(8), dimension(1:jm) :: y_AvrRad_yr
!(out) 動径平均された 1 次元緯度格子点データ
y_AvrRad_yr = y_IntRad_yr(yr)/sum(r_Rad_Weight)
end function y_AvrRad_yr
| Function : | |||
| y_IntLonRad_xyr : | real(8), dimension(1:jm)
| ||
| xyr : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
3 次元格子点データの経度動径(緯度円)積分.
3 次元データ f(λ,φ,r) に対して∫f(λ,φ,r) r^2dλdr を計算する.
function y_IntLonRad_xyr(xyr)
!
! 3 次元格子点データの経度動径(緯度円)積分.
!
! 3 次元データ f(λ,φ,r) に対して∫f(λ,φ,r) r^2dλdr を計算する.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr
!(in) 3 次元経度緯度動径格子点データ
real(8), dimension(1:jm) :: y_IntLonRad_xyr
!(out) 経度動径(緯度円)積分された 1 次元緯度格子点データ
integer :: i, k
y_IntLonRad_xyr = 0.0D0
do k=1,km
do i=0,im-1
y_IntLonRad_xyr = y_IntLonRad_xyr + xyr(i,:,k) * x_Lon_Weight(i) * r_Rad_Weight(k)
enddo
enddo
end function y_IntLonRad_xyr
| Function : | |||
| y_IntLon_xy(1:jm) : | real(8)
| ||
| xy_data(0:im-1,1:jm) : | real(8), intent(in)
|
2 次元緯度経度格子点データの経度(X)方向積分(1 層用).
実際には格子点データ各点毎に x_X_Weight をかけた総和を計算している.
Original external subprogram is wa_module_sjpack#y_IntLon_xy
| Function : | |||
| y_IntRad_yr : | real(8), dimension(1:jm)
| ||
| yr : | real(8), dimension(1:jm,km), intent(in)
|
動径積分
2 次元(YR)格子点データの動径方向域積分.
2 次元データ f(φ,r) に対して∫f(φ,r) r^2dr を計算する.
function y_IntRad_yr(yr) ! 動径積分
!
! 2 次元(YR)格子点データの動径方向域積分.
!
! 2 次元データ f(φ,r) に対して∫f(φ,r) r^2dr を計算する.
!
real(8), dimension(1:jm,km), intent(in) :: yr
!(in) 2 次元緯度動径(子午面)格子点データ
real(8), dimension(1:jm) :: y_IntRad_yr
!(out) 動径積分された 1 次元緯度格子点データ
integer :: k
y_IntRad_yr = 0.0d0
do k=1,km
y_IntRad_yr(:) = y_IntRad_yr(:) + yr(:,k) * r_Rad_Weight(k)
enddo
end function y_IntRad_yr
| Variable : | |||
| y_Lat(:) : | real(8), allocatable
|
Original external subprogram is wa_module_sjpack#y_Lat
| Variable : | |||
| y_Lat_Weight(:) : | real(8), allocatable
|
Original external subprogram is wa_module_sjpack#y_Lat_Weight
| Function : | |||
| yr_AvrLon_xyr : | real(8), dimension(1:jm,km)
| ||
| xyr : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
経度(帯状)積分
3 次元格子点データの経度方向(帯状)平均.
3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r)dλ/2π を計算する.
function yr_AvrLon_xyr(xyr) ! 経度(帯状)積分
!
! 3 次元格子点データの経度方向(帯状)平均.
!
! 3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r)dλ/2π を計算する.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr
!(in) 3 次元経度緯度動径格子点データ
real(8), dimension(1:jm,km) :: yr_AvrLon_xyr
!(out) 経度方向(帯状)平均された 2 次元子午面格子点データ
yr_AvrLon_xyr = yr_IntLon_xyr(xyr)/sum(x_Lon_Weight)
end function yr_AvrLon_xyr
| Function : | |||
| yr_IntLon_xyr : | real(8), dimension(1:jm,km)
| ||
| xyr : | real(8), dimension(0:im-1,1:jm,km), intent(in)
|
経度(帯状)積分
3 次元格子点データの経度方向(帯状)積分.
3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r)dλ を計算する.
function yr_IntLon_xyr(xyr) ! 経度(帯状)積分
!
! 3 次元格子点データの経度方向(帯状)積分.
!
! 3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r)dλ を計算する.
!
real(8), dimension(0:im-1,1:jm,km), intent(in) :: xyr
!(in) 3 次元経度緯度動径格子点データ
real(8), dimension(1:jm,km) :: yr_IntLon_xyr
!(out) 経度方向(帯状)積分された 2 次元子午面格子点データ
integer :: i
yr_IntLon_xyr = 0.0d0
do i=0,im-1
yr_IntLon_xyr(:,:) = yr_IntLon_xyr(:,:) + xyr(i,:,:) * x_Lon_Weight(i)
enddo
end function yr_IntLon_xyr