spml/wt_module レファレンスマニュアル

spml/wt_module モジュールは球面上での流体運動を スペクトル法によって数値計算するための Fortran90 関数を提供するものである. 水平方向に球面調和函数変換および 上下の境界壁を扱うためのチェビシェフ変換を用いる場合の スペクトル計算のためのさまざまな関数を提供する. 内部で wa_module, at_module を用いている. 最下部では球面調和変換およびチェビシェフ変換のエンジンとして ISPACK の Fortran77 サブルーチンを用いている.

目次


サブルーチン・関数一覧(機能別)

, wt_BoundariesGrid, wt_Boundaries, wt_TorBoundariesGrid, wt_TorBoundaries, wt_TorMagBoundariesGrid, wt_TorMagBoundaries
初期化 機能
wt_Initial スペクトル変換の格子点数, 波数, 領域の大きさの設定
座標変数 機能
x_Lon, y_Lat, z_Rad 格子点座標(経度, 緯度, 動径座標)を格納した 1 次元配列.
x_Lon_weight, y_Lat_Weight x_Rad_Weight 重み座標を格納した 1 次元配列.
xyz_Lon, xyz_Lat, xyz_Rad 格子点データの経度・緯度・動径座標(X,Y,Z)(格子点データ型 3 次元配列)
wz_Rad 格子点データの動径座標(水平スペクトル・動径格子点データ型 2 次元配列)
wt_VMiss 欠損値として使われる値(初期値は -999.0)
基本変換 機能
xyz_wt, wt_xyz スペクトルデータと 3 次元格子データの間の変換 (球面調和函数, チェビシェフ変換)
xyz_wz, wz_xyz 3 次元格子データと水平スペクトル・動径格子データとの間の変換 (球面調和変換)
wz_wt, wt_wz 水平スペクトル・動径格子データとスペクトルデータとの間の変換 (チェビシェフ変換)
xy_w, w_xy スペクトルデータと 2 次元水平格子データの間の変換 (球面調和函数変換)
az_at, at_az 格子データとチェビシェフデータの間の変換を同時に複数個行う (チェビシェフ変換)
l_nm, nm_l スペクトルデータの格納位置と全波数・帯状波数の変換
微分 機能
wt_DRad_wt スペクトルデータに動径微分 ∂/∂r を作用させる
wt_DivRad_wt スペクトルデータに発散型動径微分 1/r^2 ∂/∂r r^2 = ∂/∂r + 2/r を作用させる
wt_RotRad_wt スペクトルデータに回転型動径微分 1/r ∂/∂r r = ∂/∂r + 1/r を作用させる
wt_Lapla_wt スペクトルデータにラプラシアンを作用させる
xyz_GradLon_wt スペクトルデータに勾配型経度微分 1/r cosφ・∂/∂λ を作用させる
xyz_GradLat_wt スペクトルデータに勾配型緯度微分 1/r ∂/∂φ を作用させる
wt_DivLon_xyz 格子データに発散型経度微分 1/r cosφ・∂/∂λ を作用させる
wt_DivLat_xyz 格子データに発散型緯度微分 1/r cosφ・∂(g cosφ)/∂φ を作用させる
wt_Div_xyt_xyt_xyz, xyz_Div_xyt_xyt_xyz ベクトル成分である 3 つの格子データに発散を作用させる
xyz_RotLon_wt_wt, xyz_RotLat_wt_wt, wt_RotRad_xyz_xyz, ベクトル場の回転の各成分を計算する.
トロイダルポロイダル計算用微分 機能
wt_KxRGrad_wt スペクトルデータに経度微分 k×r・▽ = ∂/∂λ を作用させる
xyz_KGrad_wt スペクトルデータに軸方向微分 k・▽ = cosφ/r ∂/∂φ + sinφ∂/∂r を作用させる
wt_L2_wt スペクトルデータに L2 演算子 = -水平ラプラシアン を作用させる.
wt_L2Inv_wt スペクトルデータに L2 演算子の逆 = -逆水平ラプラシアン を作用させる
wt_QOperator_wt スペクトルデータに演算子 Q=(k・▽-1/2(L2 k・▽+ k・▽L2)) を 作用させる
wt_RadRot_xyz_xyz ベクトル v の渦度と動径ベクトルrの内積 r・(▽×v) を計算する
wt_RadRotRot_xyz_xyz_xyz ベクトルの vr・(▽×▽×v) を計算する
wt_Potential2Vector トロイダルポロイダルポテンシャルからベクトル場を計算する
wt_Potential2Rotation トロイダルポロイダルポテンシャルで表される非発散ベクトル場の回転の各成分を計算する
非線形計算 機能
wt_VGradV ベクトル v から v▽v を計算する
ポロイダル/トロイダルモデル用スペクトル解析 機能
nmz_ToroidalEnergySpectrum_wt nz_ToroidalEnergySpectrum_wt トロイダルポテンシャルからエネルギーの球面調和函数各成分を計算する
nmz_PoloidalEnergySpectrum_wt nz_PoloidalEnergySpectrum_wt ポロイダルポテンシャルからエネルギーの球面調和函数各成分を計算する
境界値問題 機能
wt_BoundariesTau ディリクレ, ノイマン境界条件を適用する (タウ法, 選点法)
wt_TorBoundariesTau 速度トロイダルポテンシャルの境界条件を適用する(タウ法,選点法)
wz_LaplaPol2Pol_wz, wt_LaplaPol2Pol_wt 速度ポロイダルポテンシャル Φ を ▽2Φ から求める (入出力がそれぞれチェビシェフ格子点, チェビシェフ係数)
wt_TorMagBoundariesTau 磁場トロイダルポテンシャルの境界条件を適用する(タウ法, 選点法)
wt_PolMagBoundariesTau, wt_PolMagBoundariesGrid, wt_PolMagBoundaries 磁場トロイダルポテンシャル境界の境界条件を適用する(タウ法, 選点法)
積分・平均(3 次元データ) 機能
IntLonLatRad_xyz, AvrLonLatRad_xyz 3 次元格子点データの全領域積分および平均
z_IntLonLat_xyz, z_AvrLonLat_xyz 3 次元格子点データの緯度経度(水平・球面)積分および平均
y_IntLonRad_xyz, y_AvrLonRad_xyz 3 次元格子点データの緯度動径積分および平均
z_IntLatRad_xyz, z_AvrLatRad_xyz 3 次元格子点データの経度動径(子午面)積分および平均
yz_IntLon_xyz, yz_AvrLon_xyz 3 次元格子点データの経度方向積分および平均
xz_IntLat_xyz, xz_AvrLat_xyz 3 次元格子点データの緯度方向積分および平均
xz_IntRad_xyz, xz_AvrRad_xyz 3 次元格子点データの動径方向積分および平均
積分・平均(2 次元データ) 機能
IntLonLat_xy, AvrLonLat_xy 2 次元格子点データの水平(球面)積分および平均
IntLonRad_xz, AvrLonRad_xz 2 次元(XZ)格子点データの経度動径積分および平均
IntLatRad_yz, AvrLatRad_yz 2 次元(YZ)格子点データの緯度動径(子午面)積分および平均
y_IntLon_xy, y_AvrLon_xy 水平 2 次元(球面)格子点データの経度方向積分および平均.
x_IntLat_xy, x_AvrLat_xy 水平2 次元(球面)格子点データの緯度方向積分および平均.
z_IntLon_xz, z_AvrLon_xz 2 次元(XZ)格子点データの経度方向積分および平均
x_IntRad_xz, x_AvrRad_xz 2 次元(XZ)格子点データの動径方向積分および平均
z_IntLat_yz, z_AvrLat_yz 2 次元(YZ)格子点データの緯度方向積分および平均
y_IntRad_yz, y_AvrRad_yz 2 次元(YZ)格子点データの動径方向積分および平均
積分・平均(1 次元データ) 機能
IntLon_x, AvrLon_x 1 次元(X)格子点データの経度方向積分および平均
IntLat_y, AvrLat_y 1 次元(Y)格子点データの緯度方向積分および平均.
IntRad_z, AvrRad_z 1 次元(Z)格子点データの動径方向積分および平均.


関数・変数の名前と型について

命名法

各データの種類の説明


サブルーチンの説明

subroutine wt_Initial(im,jm,km,nm,lm,rmin,rmax)

  1. 機能 : スペクトル変換の格子点数, 波数, 動径座標の範囲を設定する.
  2. 引数の説明
          integer,intent(in)      :: im, jm, km  ! 格子点数(経度λ, 緯度φ, 動径 r)  
          integer,intent(in)      :: nm, lm      ! 切断波数(水平, 動径)
          real(8)                 :: rmin, rmax  ! 動径座標の範囲(球殻内外半径)
        
  3. 備考
    他の関数を呼ぶ前に, 最初にこのサブルーチンを呼んで初期設定を しなければならない.

subroutine wt_Potential2Vector (xyz_VLon,xyz_VLat,xyz_VRad,wt_TorPot,wt_PolPot)

  1. 機能 : トロイダル・ポロイダルポテンシャルからベクトル場を計算する.
  2. 引数の説明
          real(8), dimension(im,jm,0:km)     :: xyz_VLon   ! ベクトル場(経度成分)
          real(8), dimension(im,jm,0:km)     :: xyz_VLat   ! ベクトル場(緯度成分)
          real(8), dimension(im,jm,0:km)     :: xyz_VRad   ! ベクトル場(動径成分)
    
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)  :: wt_TorPot ! トロイダルポテンシャル
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)  :: wt_PorPot ! ポロイダルポテンシャル
        
  3. 備考
    トロイダル・ポロイダルポテンシャル Ψ, Φ からベクトル場の各成分 vλ, vφ, vr は次のように計算される.
    vλ = 1/r cosφ・∂/∂r (∂Φ/∂λ)
    vφ = - 1/r cosφ・∂Ψ/∂λ + 1/r ∂/∂r (r∂Ψ/∂φ)
    vr = L2Φ/r

subroutine wt_Potential2Rotation (xyz_RotVLON,xyz_RotVLAT,xyz_RotVRAD,wt_TORPOT,wt_POLPOT)

  1. 機能 : トロイダル・ポロイダルポテンシャルからベクトル場を計算する.
  2. 引数の説明
          ! ベクトル場の回転
          real(8), dimension(im,jm,0:km), intent(OUT) :: xyz_RotVLON   ! 経度成分
          real(8), dimension(im,jm,0:km), intent(OUT) :: xyz_RotVLAT   ! 緯度成分
          real(8), dimension(im,jm,0:km), intent(OUT) :: xyz_RotVRAD   ! 動径成分
    
          ! 入力ベクトル場を表すポテンシャル
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) &
                                             :: wt_TORPOT ! トロイダルポテンシャル
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) &
                                             :: wt_POLPOT ! ポロイダルポテンシャル
        
  3. 備考
    トロイダルポロイダルポテンシャルΨ,Φで表される非発散ベクトル場 v = ▽x(Ψr) + ▽x▽x(Φr) に対して, その回転は ▽xv = ▽x▽x(Ψr) + ▽x▽x▽x(Φr) = ▽x▽x(Ψr) - ▽x((▽^2Φ)r) となるので, 回転に対するトロイダルポロイダルポテンシャルが それぞれ-▽^2Φ, Ψとなる.

subroutine wt_VGradV(xyz_VGradV_Lon,xyz_VGradV_Lat,xyz_VGradV_Rad,xyz_VLon,xyz_VLat,xyz_VRad )

  1. 機能 : ベクトル場の v・▽v を計算する.
  2. 引数の説明
          real(8), dimension(im,jm,0:km),intent(out)   :: xyz_VGradV_Lon    ! v・▽v の経度成分(出力)
          real(8), dimension(im,jm,0:km),intent(out)   :: xyz_VGradV_Lat    ! v・▽v の緯度成分(出力)
          real(8), dimension(im,jm,0:km),intent(out)   :: xyz_VGradV_Rad    ! v・▽v の動径成分(出力)
          real(8), dimension(im,jm,0:km),intent(in)    :: xyz_VLon          ! ベクトル場経度成分(入力)
          real(8), dimension(im,jm,0:km),intent(in)    :: xyz_VLat          ! ベクトル場緯度成分(入力)
          real(8), dimension(im,jm,0:km),intent(in)    :: xyz_VRad          ! ベクトル場動径成分(入力)
        
  3. 備考
    ベクトル場 v=(vλ,vφ,vr) に対するv・▽vの各成分は次のように計算される.
    (v・▽v)λ = ▽・(vλv) + vλvr/r - vλvφtan(φ)/r
    (v・▽v)φ = ▽・(vφv) + vφvr/r - vλ2tan(φ)/r
    (v・▽v)r = ▽・(vrv) + (vλ2+vφ2)/r
    非発散速度場に対してはポテンシャルから wt_Potential2Rotation を用いて回転を計算し, 恒等式 v・▽v = ▽(v2/2) - vx▽xv を用いる方がよいだろう.

function nmz_ToroidalEnergySpectrum_wt(wt_TORPOT)

  1. 機能 : トロイダルポテンシャルから, トロイダルエネルギーの球面調和函数全波数 n, 帯状波数 m の 各成分を計算する
  2. 引数の説明
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) &
               :: wt_TORPOT                       ! トロイダルポテンシャル
          real(8), dimension(0:nm,-nm:nm,0:km) &
               :: nmz_ToroidalEnergySpectrum_wz   ! スペクトルトロイダル成分
        
  3. 備考

function nz_ToroidalEnergySpectrum_wt(wt_TORPOT)

  1. 機能 : トロイダルポテンシャルから, トロイダルエネルギーの球面調和函数全波数の各成分を計算する
  2. 引数の説明
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) &
               :: wt_TORPOT                       ! トロイダルポテンシャル
          real(8), dimension(0:nm,0:km) &
               :: nz_ToroidalEnergySpectrum_wt   ! スペクトルトロイダル成分
        
  3. 備考

function nmz_PoloidamEnergySpectrum_wt(wt_TORPOT)

  1. 機能 : ポロイダルポテンシャルから, ポロイダルエネルギーの球面調和函数全波数 n, 帯状波数 m の 各成分を計算する
  2. 引数の説明
          real(8), dimension((nm+1)*(nm+1),0:km), intent(in) &
               :: wt_POLPOT                      ! ポロイダルポテンシャル
          real(8), dimension(0:nm,-nm:nm,0:km) &
               :: nmz_ToroidalEnergySpectrum_wz   ! スペクトルトロイダル成分
        
  3. 備考

function nz_PoloidamEnergySpectrum_wt(wt_TORPOT)

  1. 機能 : ポロイダルポテンシャルから, ポロイダルエネルギーの球面調和函数全波数の各成分を計算する
  2. 引数の説明
          real(8), dimension((nm+1)*(nm+1),0:km), intent(in) &
               :: wt_POLPOT                      ! ポロイダルポテンシャル
          real(8), dimension(0:nm,0:km) &
               :: nz_PoloidalEnergySpectrum_wt   ! スペクトルトロイダル成分
        
  3. 備考

subroutine wt_Boundaries(wt_TorPot,cond)

  1. 備考 : このルーチンは wt_BoundariesTau に同じ. 以前の版との整合性のために残しており, 将来廃止予定である.

subroutine wt_BoundariesTau(wt,values,cond)

  1. 機能 : スペクトルデータにディリクレ・ノイマン境界条件を適用する
  2. 引数の説明
          real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)      :: wt
                  ! 境界条件を適用するデータ. 修正された値を返す. 
    
          real(8), dimension((nm+1)*(nm+1),2), intent(in), optional :: values
                  ! 境界での 値/勾配 分布を水平スペクトル変換したものを与える. 
                  ! 省略時は値/勾配 0 となる. 
    
          character(len=2), intent(in), optional             :: cond
                  ! 境界条件. 省略時は 'DD'
                  !   DD : 両端ディリクレ
                  !   DN,ND : ディリクレ/ノイマン条件
                  !   NN : 両端ノイマン
        
  3. 備考
    チェビシェフ空間において境界条件を満たすべく高次の係数を定める方法をとっている(タウ法).

subroutine wt_BoundariesGrid(wt,values,cond)

  1. 機能 : スペクトルデータにディリクレ・ノイマン境界条件を適用する
  2. 引数の説明
          real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)      :: wt
                  ! 境界条件を適用するデータ. 修正された値を返す. 
    
          real(8), dimension((nm+1)*(nm+1),2), intent(in), optional :: values
                  ! 境界での 値/勾配 分布を水平スペクトル変換したものを与える. 
                  ! 省略時は値/勾配 0 となる. 
    
          character(len=2), intent(in), optional             :: cond
                  ! 境界条件. 省略時は 'DD'
                  !   DD : 両端ディリクレ
                  !   DN,ND : ディリクレ/ノイマン条件
                  !   NN : 両端ノイマン
        
  3. 備考
    鉛直実格子点空間において内部領域の値と境界条件を満たすように条件を課している(選点法). このルーチンを用いるためには wt_Initial にて設定するチェビシェフ切断波数(lm)と鉛直格子点数(km)を等しくしておく必要がある.

subroutine wt_TorBoundaries(wt_TorPot,cond)

  1. 備考 : このルーチンは wt_TorBoundariesTau に同じ. 以前の版との整合性のために残しており, 将来廃止予定である.

subroutine wt_TorBoundariesTau(wt_TorPot,values,cond)

  1. 機能 : 速度トロイダルポテンシャルに対して境界条件を適用する.
  2. 引数の説明
          real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)      :: wt_TorPot
                  ! 境界条件を適用するデータ. 修正された値を返す. 
    
          real(8), dimension((nm+1)*(nm+1),2), intent(in), optional :: values
                  ! 両端境界でのトロイダルポテンシャル
                  ! 粘着条件の時のみ有効
    
          character(len=2), intent(in), optional  :: cond
                  ! 境界条件スイッチ. 省略時は 'RR'
                  !   RR : 両端粘着条件
                  !   RF,FR : 粘着/応力なし条件
                  !   FF : 両端応力なし条件
        
  3. 備考
    チェビシェフ空間において境界条件を満たすべく高次の係数を定める方法をとっている(タウ法). トロイダルポテンシャルΨに対して与えられる境界条件は
    粘着条件 : Ψ = Ψb(lon,lat). Ψb は境界球面での速度分布. default は 0 (静止状態).
    応力なし条件 : ∂(Ψ/r)/∂r = 0.
    最初に呼ばれるときの境界条件で以後計算される(要仕様変更).

subroutine wt_TorBoundariesGrid(wt_TorPot,values,cond)

  1. 機能 : 速度トロイダルポテンシャルに対して境界条件を適用する.
  2. 引数の説明
          real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)      :: wt_TorPot
                  ! 境界条件を適用するデータ. 修正された値を返す. 
    
          real(8), dimension((nm+1)*(nm+1),2), intent(in), optional :: values
                  ! 両端境界でのトロイダルポテンシャル
                  ! 粘着条件の時のみ有効
    
          character(len=2), intent(in), optional  :: cond
                  ! 境界条件スイッチ. 省略時は 'RR'
                  !   RR : 両端粘着条件
                  !   RF,FR : 粘着/応力なし条件
                  !   FF : 両端応力なし条件
        
  3. 備考
    鉛直実格子点空間において内部領域の値と境界条件を満たすように条件を課している(選点法). このルーチンを用いるためには wt_Initial にて設定するチェビシェフ切断波数(lm)と鉛直格子点数(km)を等しくしておく必要がある. トロイダルポテンシャルΨに対して与えられる境界条件は
    粘着条件 : Ψ = Ψb(lon,lat). Ψb は境界球面での速度分布. default は 0 (静止状態).
    応力なし条件 : ∂(Ψ/r)/∂r = 0.
    最初に呼ばれるときの境界条件で以後計算される(要仕様変更).

subroutine wt_TorMagBoundaries(wt_TorMag)

  1. 備考 : このルーチンは wt_TorMagBoundariesTau に同じ. 以前の版との整合性のために残しており, 将来廃止予定である.

subroutine wt_TorMagBoundariesTau(wt_TorMag)

  1. 機能 : 磁場トロイダルポテンシャルに対して境界条件を適用する.
  2. 引数の説明
          real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)   :: wt_TorMag
                  ! 境界条件を適用するデータ. 修正された値を返す. 
        
  3. 備考
    チェビシェフ空間において境界条件を満たすべく高次の係数を定める方法をとっている(タウ法). 現在のところ境界物質が非電気伝導体の場合のみ対応している. その場合, 磁場トロイダルポテンシャル g の境界条件は g=0 であるので wt_BoundariesTau で対応可能だが, 将来のため作成してある. 最初に呼ばれるときの境界条件で以後計算される(要仕様変更).

subroutine wt_TorMagBoundariesGrid(wt_TorMag)

  1. 機能 : 磁場トロイダルポテンシャルに対して境界条件を適用する.
  2. 引数の説明
          real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)   :: wt_TorMag
                  ! 境界条件を適用するデータ. 修正された値を返す. 
        
  3. 備考
    鉛直実格子点空間において内部領域の値と境界条件を満たすように条件を課している(選点法). このルーチンを用いるためには wt_Initial にて設定するチェビシェフ切断波数(lm)と鉛直格子点数(km)を等しくしておく必要がある. 現在のところ境界物質が非電気伝導体の場合のみ対応している. その場合, 磁場トロイダルポテンシャル g の境界条件は g=0 であるので wt_BoundariesTau で対応可能だが, 将来のため作成してある. 最初に呼ばれるときの境界条件で以後計算される(要仕様変更).

subroutine wt_PolMagBoundaries(wt_TorMag)

  1. 備考 : このルーチンは wt_PolMagBoundariesTau に同じ. 以前の版との整合性のために残しており, 将来廃止予定である.

subroutine wt_PolMagBoundariesTau(wt_PolMag)

  1. 機能 : 磁場ポロイダルポテンシャルに対して境界条件を適用する.
  2. 引数の説明
          real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)   :: wt_PolMag
                  ! 境界条件を適用するデータ. 修正された値を返す. 
        
  3. 備考
    チェビシェフ空間において境界条件を満たすべく高次の係数を定める方法をとっている(タウ法). 現在のところ境界物質が非電気伝導体の場合のみ対応している. その場合, 磁場ポロイダルポテンシャルの各水平スペクトル成分 h に たいして境界条件が与えられ,
    外側境界 : dh/dr + (n+1)h/r = 0
    内側境界 : dh/dr - nh/r = 0
    である. ここで n は h の水平全波数である. 最初に呼ばれるときの境界条件で以後計算される(要仕様変更).

subroutine wt_PolMagBoundariesGrid(wt_PolMag)

  1. 機能 : 磁場ポロイダルポテンシャルに対して境界条件を適用する.
  2. 引数の説明
          real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)   :: wt_PolMag
                  ! 境界条件を適用するデータ. 修正された値を返す. 
        
  3. 備考
    鉛直実格子点空間において内部領域の値と境界条件を満たすように条件を課している(選点法). このルーチンを用いるためには wt_Initial にて設定するチェビシェフ切断波数(lm)と鉛直格子点数(km)を等しくしておく必要がある. 現在のところ境界物質が非電気伝導体の場合のみ対応している. その場合, 磁場ポロイダルポテンシャルの各水平スペクトル成分 h に たいして境界条件が与えられ,
    外側境界 : dh/dr + (n+1)h/r = 0
    内側境界 : dh/dr - nh/r = 0
    である. ここで n は h の水平全波数である. 最初に呼ばれるときの境界条件で以後計算される(要仕様変更).


変数の説明

wa_module と共通な変数は省略する. wa_module の変数の説明の項を参照されたい.

z_Rad

  1. 説明 : 動径格子点座標を格納した 1 次元配列.
  2. 変数の型
          real(8), dimension(0:km) :: z_Rad
        
  3. 備考
    チェビシェフ変換用の格子点座標. at_module を参照のこと.

z_Rad_Weigtht

  1. 説明 : 動径重み座標を格納した 1 次元配列.
  2. 変数の型
          real(8), dimension(0:km) :: z_Rad_Weigtht
        
  3. 備考
    z_Rad_Weight には動径方向の積分の各格子点での重み(r2dr)を格納してある.

xyz_Lon, xyz_Lat, xyz_Rad

  1. 説明 : 各格子点(i,j,k)の位置の経度・緯度・動径座標を格納した格子データ.
  2. 変数の型
          real(8), dimension(im,jm,0:km) :: xyz_Lon, xyz_Lat, xyz_Rad
        
  3. 備考

wz_Rad

  1. 説明 : 水平スペクトル・動径格子点空間の各点(n,k)での動径座標を格納した格子データ.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:km) :: wz_Rad
        
  3. 備考

各関数の説明

以下では w_module.htm と共通でない関数のみ説明してある. w_module と共通な関数は w_module の関数の説明の項を参照されたい.

function xyz_wt(wt_data)

  1. 機能 : スペクトルデータから 3 次元格子データへ(逆)変換する.
  2. 変数の型
          real(8), dimension(im,jm,0:km)                     :: xyz_wt
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt
        
  3. 備考

function wt_xyz(xyz_data)

  1. 機能 : 3 次元格子データからスペクトルデータへ(正)変換する.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:lm)             :: wt_xyz
          real(8), dimension(im,jm,0:km), intent(in)         :: xyz
        
  3. 備考

function xyz_wz(wz_data)

  1. 機能 : 水平スペクトル・動径格子点データから 3 次元格子データへ(逆)変換する.
  2. 変数の型
          real(8), dimension(im,jm,0:km)                     :: xyz_wz
          real(8), dimension((nm+1)*(nm+1),0:km), intent(in) :: wz_data
        
  3. 備考

function wz_xyz(xyz_data)

  1. 機能 : 3 次元格子データから水平スペクトル・動径格子点データへ(正)変換する.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:km)             :: wz_xyz
          real(8), dimension(im,jm,0:km), intent(in)         :: xyz
        
  3. 備考

function wt_wz(wz_data)

  1. 機能 : 水平スペクトル・動径格子点データからスペクトルデータへ(正)変換する.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:km)             :: wt_wz
          real(8), dimension((nm+1)*(nm+1),0:km), intent(in) :: wz
        
  3. 備考

function wz_wt(wt_data)

  1. 機能 : スペクトルデータから水平スペクトル・動径格子点データへ(正)変換する.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:km)             :: wz_wt
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_data
        
  3. 備考

function wt_DRad_wt(wt_data)

  1. 機能 : 入力スペクトルデータに動径微分を作用する.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_data
          real(8), dimension((nm+1)*(nm+1),0:lm)             :: wt_DRad_wt
        
  3. 備考
    スペクトルデータの動径微分とは, 対応する格子点データに動径微分 ∂/∂r を作用させたデータのスペクトル変換のことである.

function wt_DivRad_wt(wt_data)

  1. 機能 : 入力スペクトルデータに発散型動径微分を作用する.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_data
          real(8), dimension((nm+1)*(nm+1),0:lm)             :: wt_DivRad_wt
        
  3. 備考
    スペクトルデータの発散型動径微分とは, 対応する格子点データに発散型動径微分を作用させたデータのスペクトル変換のことである. 発散型動径微分は 1/r^2 ∂/∂r r^2 = ∂/∂r + 2/r と計算される.

function wt_RotRad_wt(wt_data)

  1. 機能 : 入力スペクトルデータに回転型動径微分を作用する.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_data
          real(8), dimension((nm+1)*(nm+1),0:lm)             :: wt_RotRad_wt
        
  3. 備考
    スペクトルデータの回転型動径微分とは, 対応する格子点データに回転型動径微分を作用させたデータのスペクトル変換のことである. 回転型動径微分は 1/r ∂/∂r r = ∂/∂r + 1/r と計算される.

function wt_Lapla_wt(wt_data)

  1. 機能 : 入力スペクトルデータにラプラシアンを作用する.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt
          real(8), dimension((nm+1)*(nm+1),0:lm)             :: wt_Lapla_wt
        
  3. 備考
    スペクトルデータのラプラシアンとは, 対応する格子点データにラプラシアンを作用させたデータのスペクトル変換のことである. ラプラシアンは ▽^2 = 1/r^2 cos^2φ・∂^2/∂λ^2 + 1/r^2 cosφ・∂/∂φ(cosφ∂/∂φ) + 1/r^2 ∂/∂r (r^2 ∂/∂r) と計算される.

function xyz_GradLon_wt(wt_data)

  1. 機能 : スペクトルデータに勾配型経度微分を作用させる.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_data
          real(8), dimension(im,jm,0:km)                     :: xyz_GradLon_wt
        
  3. 備考
    スペクトルデータに対応する格子点データに 勾配型経度微分 1/r cosφ・∂/∂λ を作用させた格子データが返される

function xyz_GradLat_wt(wt_data)

  1. 機能 : スペクトルデータに勾配型経度微分を作用させる.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_data
          real(8), dimension(im,jm,0:km)                     :: xyz_GradLat_wt
        
  3. 備考
    入力スペクトルデータに対応する格子データに 勾配型緯度微分 1/r ∂/∂φ を作用させた格子データが返される

function wt_DivLon_xyz(xyz_data)

  1. 機能 : 格子データに発散型経度微分を作用させる.
  2. 変数の型
          real(8), dimension(im,jm,0:km), intent(in)   :: xyz
          real(8), dimension((nm+1)*(nm+1),0:lm)       :: wt_DivLon_xyz
        
  3. 備考
    入力格子データに発散型経度微分 1/r cosφ・∂/∂λ を作用させた スペクトルデータが返される

function wt_DivLat_xyz(xyz_data)

  1. 機能 : 格子データに発散型緯度微分を作用させる.
  2. 変数の型
          real(8), dimension(im,jm,0:km), intent(in)   :: xyz_data
          real(8), dimension((nm+1)*(nm+1),0:lm)       :: wt_DivLat_xyz
        
  3. 備考
    入力格子データに発散型経度微分 1/rcosφ・∂(cosφ)/∂φ を作用させた スペクトルデータが返される

function wt_Div_xyz_xyz(xyz_u,xyz_v,xyz_w)

  1. 機能 : ベクトル成分である 3 つの格子データに発散を作用させる.
  2. 変数の型
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_Vlon   ! 経度成分
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_Vlat   ! 緯度成分
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_Vrad   ! 動径成分
          real(8), dimension((nm+1)*(nm+1),0:lm)     :: wt_Div_xyz_xyz_xyz
        
  3. 備考
    第 1, 2 ,3 引数(u,v,w)がそれぞれベクトルの経度成分, 緯度成分, 動径成分を表し, 発散は 1/rcosφ・∂u/∂λ + 1/rcosφ・∂(v cosφ)/∂φ + 1/r^2 ∂/∂r (r^2 w)と計算される.

function xyz_Div_xyz_xyz(xyz_u,xyz_v,xyz_w)

  1. 機能 : ベクトル成分である 3 つの格子データに発散を作用させる.
  2. 変数の型
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_Vlon   ! 経度成分
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_Vlat   ! 緯度成分
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_Vrad   ! 動径成分
          real(8), dimension(im,jm,0:km)             :: xyz_Div_xyz_xyz_xyz
        
  3. 備考
    第 1, 2 ,3 引数(u,v,w)がそれぞれベクトルの経度成分, 緯度成分, 動径成分を表す. 極の特異性を回避するためにベクトル場に cosφ/r の重みをかけて, div V = (r/cosφ)・div (Vcosφ/r) + V_φtanφ/r + V_r/r なる恒等式を用いて計算している.

function xyz_RotLon_wt_wt(wt_Vrad,wt_Vlat)

  1. 機能 : ベクトル場の回転の経度成分を計算する.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_Vrad ! 動径成分
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_Vlat ! 緯度成分
          real(8), dimension(im,jm,0:km)                     :: xyz_RotLon_wt_wt
        
  3. 備考
    ベクトル場のの動径成分, 緯度成分である第 1, 2 引数 Vrad, Vlat から 回転の経度成分 1/rcosφ・∂Vrad/∂λ-1/r ∂(r Vlat)/∂r を計算する.

function xyz_RotLat_wt_wt(wt_Vlon,wt_Vrad)

  1. 機能 : ベクトル場の回転の緯度成分を計算する.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_Vlon ! 緯度成分
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_Vrad ! 動径成分
          real(8), dimension(im,jm,0:km)                     :: xyz_RotLat_wt_wt
        
  3. 備考
    ベクトル場の経度成分, 動径成分である第 1, 2 引数 Vlon, Vrad から 回転の緯度成分 1/r ∂(r Vlon)/∂r - 1/r ∂Vrad/∂φ を計算する.

function wt_RotRad_xyz_xyz(xyz_Vlat,xyz_Vlon)

  1. 機能 : ベクトル場の回転の動径成分を計算する.
  2. 変数の型
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_Vlat ! 緯度成分
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_Vlon ! 経度成分
          real(8), dimension((nm+1)*(nm+1),0:lm)     :: wt_RotRad_xyz_xyz
        
  3. 備考
    ベクトルの緯度成分, 経度成分である第 1, 2 引数 Vlat, Vlon に対して ベクトル場の回転の動径成分 1/rcosφ・∂Vlat/∂λ - 1/rcosφ・∂(Vlon cosφ)/∂φ を計算する.

function wt_KxRGrad_wt(wt_data)

  1. 機能 : 入力スペクトルデータに経度微分 k×r・▽ = ∂/∂λ を作用する.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_data
          real(8), dimension((nm+1)*(nm+1),0:lm)             :: wt_KxRGrad_wt
        
  3. 備考
    スペクトルデータの経度微分とは, 対応する格子点データに経度微分 k×r・▽ = ∂/∂λ を作用させた データのスペクトル変換のことである. ここでベクトル k は球の中心から北極向きの単位ベクトルである.

function xyz_KGrad_wt(wt_data)

  1. 機能 : 入力スペクトルデータに軸方向微分 k・▽ = cosφ/r ∂/∂φ + sinφ∂/∂r を作用する.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_data
          real(8), dimension(im,jm,0:km)                     :: xyz_KGrad_wt
        
  3. 備考
    入力スペクトルデータに対応する格子データに 軸方向微分 k・▽ = cosφ/r ∂/∂φ + sinφ∂/∂r を作用させた 格子データが返される ここでベクトル k は球の中心から北極向きの単位ベクトルである.

function wt_L2_wt(wt_data)

  1. 機能 : 入力スペクトルデータに L2 演算子を作用する.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_data
          real(8), dimension((nm+1)*(nm+1),0:lm)             :: wt_L2_wt
        
  3. 備考
    L2 演算子は単位球面上の水平ラプラシアンの逆符号にあたる. 入力スペクトルデータに対応する格子点データに 演算子 L2 = -1/cos^2φ・∂^2/∂λ^2 - 1/cosφ・∂/∂φ(cosφ∂/∂φ) を作用させたデータのスペクトル変換が返される.

function wt_L2Inv_wt(wt_data)

  1. 機能 : 入力スペクトルデータに L2 演算子の逆演算を作用する.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_data
          real(8), dimension((nm+1)*(nm+1),0:lm)             :: wt_L2Inv_wt
        
  3. 備考
    スペクトルデータに L2 演算子を作用させる関数 wt_L2_wt の逆計算を行う関数である.

function wt_QOperator_wt(wt_data)

  1. 機能 : 入力スペクトルデータに Q 演算子を作用する.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_data
          real(8), dimension((nm+1)*(nm+1),0:lm)             :: wt_QOperator_wt
        
  3. 備考
    入力スペクトルデータに対応する格子点データに演算子 Q=(k・▽-1/2(L2 k・▽+ k・▽L2)) を作用させたデータの スペクトル変換が返される.

function wt_RadRot_xyz_xyz(xyz_Vlon,xyz_Vlat)

  1. 機能 : ベクトルの渦度と動径ベクトルの内積 r・(▽×v) を計算する.
  2. 変数の型
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_VLon
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_VLat
          real(8), dimension((nm+1)*(nm+1),0:lm)     :: wt_RadRot_xyz_xyz
        
  3. 備考
    第 1, 2 引数(vλ, vφ)が それぞれベクトルの経度成分, 緯度成分を表す. r・(▽×v) = 1/cosφ・∂vφ/∂λ - 1/cosφ・∂(vλ cosφ)/∂φ のスペクトルデータが返される.

function wt_RadRotRot_xyz_xyz_xyz(xyz_VLon,xyz_VLat,xyz_VRad)

  1. 機能 : ベクトル v に対して r・(▽×▽×v) を計算する.
  2. 変数の型
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_VLon
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_VLat
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_VRad
          real(8), dimension((nm+1)*(nm+1),0:lm)     :: wt_RadRotRot_xyz_xyz_xyz
        
  3. 備考
    第 1, 2, 3 引数(vλ, vφ, vr)が それぞれベクトルの 経度成分, 緯度成分, 動径成分を表す. r・(▽×▽×v) = 1/r ∂/∂r (r・( 1/cosφ・∂vλ/∂λ + 1/cosφ・∂(vφ cosφ)/∂φ ) ) + L2 vr/r のスペクトルデータが返される.

function wz_LaplaPol2Pol_wz(wz_data,cond)

  1. 機能 : 速度ポロイダルポテンシャル Φ を ▽2Φ から計算する.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:km),intent(in)  :: wz_data
                  ! 入力▽^2Φ分布
    
          real(8), dimension((nm+1)*(nm+1),0:km)             :: wz_LaplaPol2Pol_wz
                  ! 出力ポロイダルポテンシャル分布
    
          character(len=2), intent(in), optional  :: cond
                  ! 境界条件スイッチ. 省略時は 'RR'
                  !   RR : 両端粘着条件
                  !   RF,FR : 粘着/応力なし条件
                  !   FF : 両端応力なし条件
        
  3. 備考
    チェビシェフ格子点空間で境界条件を適用している. この関数を用いるためには "wt_Initial にて設定するチェビシェフ切断波数(lm)と鉛直格子点数(km)を等しくしておく必要がある. 速度ポロイダルポテンシャル Φ を f = ▽2Φ から定める式は
    2Φ = f
    Φ = const. at boundaries.
    ∂Φ/∂r = 0 at boundaries (粘着条件) or ∂2Φ/∂r2 = 0 at boundaries (応力なし条件) .
    最初に呼ばれるときの境界条件で以後計算される(要仕様変更).

function wt_LaplaPol2Pol_wt(wt_data,cond)

  1. 機能 : 速度ポロイダルポテンシャル Φ を ▽2Φ から計算する.
  2. 変数の型
          real(8), dimension((nm+1)*(nm+1),0:lm),intent(in)  :: wt_data
                  ! 入力▽^2Φ分布
    
          real(8), dimension((nm+1)*(nm+1),0:lm)             :: wt_LaplaPol2Pol_wt
                  ! 出力ポロイダルポテンシャル分布
    
          character(len=2), intent(in), optional  :: cond
                  ! 境界条件スイッチ. 省略時は 'RR'
                  !   RR : 両端粘着条件
                  !   RF,FR : 粘着/応力なし条件
                  !   FF : 両端応力なし条件
        
  3. 備考
    チェビシェフ格子点空間で境界条件を適用している. この関数を用いるためには "wt_Initial にて設定するチェビシェフ切断波数(lm)と鉛直格子点数(km)を等しくしておく必要がある. 速度ポロイダルポテンシャル Φ を f = ▽2Φ から定める式は
    2Φ = f
    Φ = const. at boundaries.
    ∂Φ/∂r = 0 at boundaries (粘着条件) or ∂2Φ/∂r2 = 0 at boundaries (応力なし条件) .
    最初に呼ばれるときの境界条件で以後計算される(要仕様変更). 最終的にチェビシェフ係数の解が欲しい場合には, wz_LaplaPol2Pol_wz に比べて チェビシェフ -- 格子点変換が 1 回分少なくて済む.

function IntLonLatRad_xyz(xyz_data), function AvrLonLatRad_xyz(xyz_data)

  1. 機能 : 3 次元格子点データの緯度経度動径(全球)積分および平均
  2. 変数の型
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_data  ! 3 次元格子点データ
          real(8)                     :: IntLonLatRad_xyz         ! 全球積分値
    
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_data  ! 3 次元格子点データ
          real(8)                     :: AvrLonLatRad_xyz         ! 全球平均値
        
  3. 備考
    3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r) r2cosφ dλdφdr および ∫f(λ,φ,r) r2cosφ dλdφdr /(4π(ro3-ri3)/3) を計算する.

function z_IntLonLat_xyz(xyz_data), function z_AvrLonLat_xyz(xyz_data)

  1. 機能 : 3 次元格子点データの緯度経度(水平, 球面)積分および平均
  2. 変数の型
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_data  ! 3 次元格子点データ
          real(8), dimension(0:km)     :: z_IntLonLat_xyz         ! 動径格子点データ
    
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_data  ! 3 次元格子点データ
          real(8), dimension(0:km)     :: z_AvrLonLat_xyz         ! 動径格子点データ
        
  3. 備考
    3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r) cosφ dλdφ および ∫f(λ,φ,r) cosφ dλdφdr /4πを計算する.

function y_IntLonRad_xyz(xyz_data), function y_AvrLonRad_xyz(xyz_data)

  1. 機能 : 3 次元格子点データの経度動径積分および平均
  2. 変数の型
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_data ! 3 次元格子点データ
          real(8), dimension(jm)       :: y_IntLonRad_xyz        ! 緯度格子点データ
    
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_daa  ! 3 次元格子点データ
          real(8), dimension(jm)       :: y_AvrLonRad_xyz        ! 緯度格子点データ
        
  3. 備考
    3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r) r2dλdr および ∫f(λ,φ,r) r2dλdr /(2π(ro3-ri3)/3) を計算する.

function x_IntLatRad_xyz(xyz_data), function x_AvrLatRad_xyz(xyz_data)

  1. 機能 : 3 次元格子点データの緯度動径(子午面)積分および平均
  2. 変数の型
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_data ! 3 次元格子点データ
          real(8), dimension(im)     :: x_IntLatRad_xyz          ! 経度格子点データ
    
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_data ! 3 次元格子点データ
          real(8), dimension(im)     :: x_AvrLatRad_xyz          ! 経度格子点データ
        
  3. 備考
    3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r) r2cosφ dφdr および ∫f(λ,φ,r) r2cosφ dφdr /(2(ro3-ri3)/3) を計算する.

function yz_IntLon_xyz(xya_data), function yz_AvrLon_xyz(xya_data)

  1. 機能 : 3 次元格子点データの経度方向積分および平均.
  2. 変数の型
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_data  ! 3 次元格子点データ
          real(8), dimension(jm,0:km)  :: yz_IntLon_xyz           ! 子午面格子点データ
    
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_data  ! 3 次元格子点データ
          real(8), dimension(jm,0:km)  :: yz_AvrLon_xyz           ! 子午面格子点データ
        
  3. 備考
    3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r)dλ および ∫f(λ,φ,r)dλ/2π を計算する.

function xz_IntLat_xyz(xyz_data), function xz_AvrLat_xyz(xyz_data)

  1. 機能 : 3 次元格子点データの緯度方向域積分および平均.
  2. 変数の型
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_data  ! 3 次元格子点データ
          real(8), dimension(im,0:km)  :: xz_IntLat_xyz           ! 緯度円格子点データ
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_data  ! 3 次元格子点データ
          real(8), dimension(im,0:km)  :: xz_AvrLat_xyz           ! 緯度円格子点データ
        
  3. 備考
    3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r) cosφ dφ および ∫f(λ,φ,r) cosφ dφ/2 を計算する.

function xy_IntRad_xyz(xyz_data), function xy_AvrRad_xyz(xyz_data)

  1. 機能 : 3 次元格子点データの動径方向域積分および平均.
  2. 変数の型
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_data ! 3 次元格子点データ
          real(8), dimension(im,jm)  :: xy_IntRad_xyz            ! 水平格子点データ
    
          real(8), dimension(im,jm,0:km), intent(in) :: xyz_data ! 3 次元格子点データ
          real(8), dimension(im,jm)  :: xy_AvrRad_xyz            ! 水平格子点データ
        
  3. 備考
    3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r) r2dr および ∫f(λ,φ,r) r2dr /((ro3-ri3)/3) を計算する.

function IntLonRad_xz(xz_data), function AvrLonRad_xz(xz_data)

  1. 機能 : 2 次元(XZ)格子点データの経度動径積分および平均
  2. 変数の型
          real(8), dimension(im,0:km), intent(in) :: xz_data ! 2 次元格子点データ
          real(8) :: y_IntLonRad_xyz                         ! 積分値
    
          real(8), dimension(im,0:km), intent(in) :: xz_data ! 2 次元格子点データ
          real(8) :: y_AvrLonRad_xyz                         ! 平均値
        
  3. 備考
    2 次元データ f(λ,r) に対して ∫f(λ,r) r2dλdr および ∫f(λ,r) r2dλdr /(2π(ro3-ri3)/3) を計算する.

function z_IntLon_xz(xz_data), function z_AvrLon_xz(xz_data)

  1. 機能 : 2 次元(XZ)格子点データの経度方向積分および平均.
  2. 変数の型
          real(8), dimension(im,0:km), intent(in) :: xz_data  ! 3 次元格子点データ
          real(8), dimension(0:km)  :: yz_IntLon_xyz           ! 子午面格子点データ
    
          real(8), dimension(im,0:km), intent(in) :: xz_data  ! 3 次元格子点データ
          real(8), dimension(0:km)  :: z_AvrLon_xyz           ! 子午面格子点データ
        
  3. 備考
    2 次元データ f(λ,r) に対して ∫f(λ,r)dλ および ∫f(λ,r)dλ/2π を計算する.

function x_IntRad_xz(xz_data), function x_AvrRad_xz(xz_data)

  1. 機能 : 2 次元(XZ)格子点データの動径方向域積分および平均.
  2. 変数の型
          real(8), dimension(im,0:km), intent(in) :: xz_data ! 2 次元格子点データ
          real(8), dimension(im)  :: x_IntRad_xz             ! 水平格子点データ
    
          real(8), dimension(im,0:km), intent(in) :: xz_data ! 2 次元格子点データ
          real(8), dimension(im)  :: x_AvrRad_xz             ! 水平格子点データ
        
  3. 備考
    2 次元データ f(λ,r) に対して ∫f(λ,r) r2dr および ∫f(λ,r) r2dr /((ro3-ri3)/3) を計算する.

function IntLatRad_yz(yz_data), function AvrLatRad_yz(yz_data)

  1. 機能 : 2 次元(YZ)格子点データの緯度動径積分および平均
  2. 変数の型
          real(8), dimension(jm,0:km), intent(in) :: yz_data ! 2 次元格子点データ
          real(8)  :: IntLatRad_yz                           ! 積分値
    
          real(8), dimension(jm,0:km), intent(in) :: yz_daa  ! 3 次元格子点データ
          real(8) :: AvrLatRad_yz                            ! 平均値
        
  3. 備考
    2 次元データ f(φ,r) に対して ∫f(φ,r) r2cosφ dφdr および ∫f(φ,r) r2cosφ dφdr /(2(ro3-ri3)/3) を計算する.

function z_IntLat_yz(yz_data), function z_AvrLat_yz(yz_data)

  1. 機能 : 2 次元(YZ)格子点データの緯度方向域積分および平均.
  2. 変数の型
          real(8), dimension(jm,0:km), intent(in) :: yz_data  ! 2 次元格子点データ
          real(8), dimension(0:km)  :: z_IntLat_xz            ! 動径格子点データ
          real(8), dimension(im,0:km), intent(in) :: xz_data  ! 2 次元格子点データ
          real(8), dimension(0:km)  :: z_AvrLat_xz            ! 動径格子点データ
        
  3. 備考
    3 次元データ f(φ,r) に対して ∫f(φ,r) cosφ dφ および ∫f(φ,r) cosφ dφ/2 を計算する.

function y_IntRad_yz(yz_data), function y_AvrRad_yz(yz_data)

  1. 機能 : 2 次元(YZ)格子点データの動径方向域積分および平均.
  2. 変数の型
          real(8), dimension(jm,0:km), intent(in) :: yz_data ! 2 次元格子点データ
          real(8), dimension(jm)  :: xy_IntRad_xyz           ! 水平格子点データ
    
          real(8), dimension(jm,0:km), intent(in) :: yz_data ! 2 次元格子点データ
          real(8), dimension(jm)  :: y_AvrRad_xyz            ! 水平格子点データ
        
  3. 備考
    2 次元データ f(φ,r) に対して ∫f(φ,r) r2dr および ∫f(φ,r) r2dr /((ro3-ri3)/3) を計算する.

function IntRad_z(z_data), function AvrRad_z(z_data)

  1. 機能 : 1 次元(Z)格子点データの動径方向域積分および平均.
  2. 変数の型
          real(8), dimension(0:km), intent(in) :: z_data ! 1 次元格子点データ
          real(8)  :: IntRad_z                           ! 積分値
    
          real(8), dimension(0:km), intent(in) :: z_data ! 1 次元格子点データ
          real(8) :: AvrRad_z                            ! 平均値
        
  3. 備考
    1 次元データ f(r) に対して ∫f(r) r2dr および ∫f(r) r2dr /((ro3-ri3)/3) を計算する.


地球流体電脳倶楽部 SPMODEL プロジェクト
spmodel@gfd-dennou.org

2002/08/24 作成 (竹広真一)
2005/04/24 更新 (竹広真一)