Class wt_mpi_module
In: libsrc/wt_mpi_module/wt_mpi_module.f90

wt_mpi_module

Authors:Shin-ichi Takehiro, Youhei SASAKI
Version:$Id: wt_mpi_module.f90 598 2013-08-20 03:23:44Z takepiro $
Copyright&License:See COPYRIGHT

概要

spml/wt_mpi_module モジュールは球面上および球殻内での流体運動を スペクトル法と MPI 並列化によって数値計算するための Fortran90 関数を提供するものである.

水平方向に球面調和函数変換および上下の境界壁を扱うための チェビシェフ変換を用いる場合のスペクトル計算のためのさまざまな 関数を提供する.

内部で wt_module, wa_mpi_module を用いている. 最下部では球面調和変換 およびチェビシェフ変換のエンジンとして ISPACK の Fortran77 サブルーチンを用いている.

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

命名法

  • 関数名の先頭(wt_, nmz_, nz_, xyz_, xvz_ wz_, w_, xy_, x_, y_, z_, a_)は, 返す値の形を示している.
       wt_ :: スペクトルデータ(球面調和函数・チェビシェフ変換)
       nmz_:: 水平スペクトルデータ(全波数 n, 帯状波数各成分, 動径)
       nz_ :: 水平スペクトルデータ(全波数 n, 動径)
      xyz_ :: 3 次元格子点データ(経度・緯度・動径)
      xvz_ :: 3 次元分散格子点データ(経度・緯度・動径)
       wz_ :: 水平スペクトル, 動径格子点データ
    
  • 関数名の間の文字列(DLon, GradLat, GradLat, DivLon, DivLat, Lapla,..) は, その関数の作用を表している.
  • 関数名の最後 (_wt, _xyz, _xvz, wz, _w, _xy, _xv, _x, _y, _z, _a) は, 入力変数の形がスペクトルデータおよび格子点データであることを示している.
          _wt :: スペクトルデータ
         _xyz :: 3 次元格子点データ
     _xyz_xyz :: 2 つの3 次元格子点データ, ...
         _xvz :: 3 次元分散格子点データ
     _xvz_xvz :: 2 つの3 次元分散格子点データ, ...
    

各データの種類の説明

  • xyz : 3 次元格子点データ(経度・緯度・動径)
    • 変数の種類と次元は real(8), dimension(0:im-1,jm,0:km).
    • im, jm, km はそれぞれ経度, 緯度, 動径座標の格子点数であり, サブルーチン wt_mpi_Initial にてあらかじめ設定しておく.
  • xvz : 3 次元分散格子点データ(経度・緯度・動径)
    • 変数の種類と次元は real(8), dimension(0:im-1,jc,0:km).
    • im, km はそれぞれ経度, 動径座標の格子点数であり,
    • jc はこのプロセスで保有する緯度格子点数である. サブルーチン wt_mpi_Initial を呼ぶと jc が設定される.
  • wt : スペクトルデータ
    • 変数の種類と次元は real(8), dimension((nm+1)*(nm+1),0:lm).
    • nm は球面調和函数の最大全波数, lm はチェビシェフ多項式の最大次数 であり, サブルーチン wt_mpi_Initial にてあらかじめ設定しておく.
    • 水平スペクトルデータの格納のされ方は関数 l_nm, nm_l によって調べる ことができる.
  • nmz : 水平スペクトルデータの並んだ 3 次元配列.
    • 変数の種類と次元は real(8), dimension(0:nm,-nm:nm,0:km).
    • 第 1 次元が水平全波数, 第 2 次元が帯状波数, 第 3 次元が動径座標を表す.
    • nm は球面調和函数の最大全波数であり, サブルーチン wt_Initial にて あらかじめ設定しておく.
  • nz : スペクトルデータの並んだ 2 次元配列.
    • 変数の種類と次元は real(8), dimension(0:nm,0:km). 第 1 次元が水平全波数を表す.
    • nm は球面調和函数の最大全波数であり, サブルーチン wt_Initial にてあらかじめ設定しておく.
  • wz : 水平スペクトル, 動径格子点データ.
    • 変数の種類と次元は real(8), dimension((nm+1)*(nm+1),0:km).
  • wt_ で始まる関数が返す値はスペクトルデータに同じ.
  • xyz_ で始まる関数が返す値は 3 次元格子点データに同じ.
  • xvz_ で始まる関数が返す値は 3 次元分散格子点データに同じ.
  • wz_ で始まる関数が返す値は水平スペクトル, 動径格子点データに同じ.
  • スペクトルデータに対する微分等の作用とは, 対応する格子点データに 微分などを作用させたデータをスペクトル変換したものことである.

変数・手続き群の要約

初期化

wt_mpi_Initial :スペクトル変換の格子点数, 波数, 領域の大きさの設定

座標変数

x_Lon, x_Lon_Weight, :格子点座標(経度)と重みを格納した 1 次元配列
y_Lat, y_Lat_Weight :格子点座標(緯度)と重みを格納した 1 次元配列
v_Lat, v_Lat_Weight :分散格子点座標(緯度)と重みを格納した 1 次元配列
z_Rad, z_Rad_Weight :格子点座標(同形)と重みを格納した 1 次元配列
xyz_Lon, xyz_Lat, xyz_Rad :格子点データの経度・緯度・動径座標(X,Y,Z) (格子点データ型 3 次元配列)
xvz_Lon, xvz_Lat, xvz_Rad :格子点データの経度・緯度・動径座標(X,Y,Z) (分散格子点データ型 3 次元配列)

基本変換

xyz_wt, wt_xyz :スペクトルデータと 3 次元格子データの間の変換 (球面調和函数, チェビシェフ変換)
xvz_wt, wt_xvz :スペクトルデータと 3 次元分散格子データの間の変換 (球面調和函数, チェビシェフ変換)
xyz_wz, wz_xyz :3 次元格子データと水平スペクトル・動径格子データとの間の変換 (球面調和函数)
xvz_wz, wz_xvz :3 次元分散格子データと水平スペクトル・動径格子データとの間の変換 (球面調和函数)
wz_wt, wt_wz :スペクトルデータと水平スペクトル・動径格子データとの間の変換 (チェビシェフ変換)
w_xy, xy_w :スペクトルデータと 2 次元水平格子データの間の変換(球面調和函数変換)
w_xv, xv_w :スペクトルデータと 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 ∂/∂rr = ∂/∂r + 1/r を作用させる
wt_Lapla_wt :スペクトルデータにラプラシアンを作用させる
xyz_GradLon_wt :スペクトルデータに勾配型経度微分 1/rcosφ・∂/∂λを作用させる
xvz_GradLon_wt :スペクトルデータに勾配型経度微分 1/rcosφ・∂/∂λを作用させる
xyz_GradLat_wt :スペクトルデータに勾配型緯度微分 1/r・∂/∂φを作用させる
xvz_GradLat_wt :スペクトルデータに勾配型緯度微分 1/r・∂/∂φを作用させる(分散格子)
wt_DivLon_xyz :格子データに発散型経度微分 1/rcosφ・∂/∂λを作用させる
wt_DivLon_xvz :分散格子データに発散型経度微分 1/rcosφ・∂/∂λを作用させるa

Methods

AvrLatRad_vz   AvrLatRad_yz   AvrLat_v   AvrLat_y   AvrLonLatRad_xvz   AvrLonLatRad_xyz   AvrLonLat_xv   AvrLonLat_xy   AvrLonRad_xz   AvrLon_x   AvrRad_z   IntLatRad_vz   IntLatRad_yz   IntLat_v   IntLat_y   IntLonLatRad_xvz   IntLonLatRad_xyz   IntLonLat_xv   IntLonLat_xy   IntLonRad_xz   IntLon_x   IntRad_z   Interpolate_wt   at_Dr_at   at_az   az_at   jc   l_nm   l_nm   l_nm   l_nm   nm_l   nm_l   nmz_PoloidalEnergySpectrum_wt   nmz_ToroidalEnergySpectrum_wt   nz_PoloidalEnergySpectrum_wt   nz_ToroidalEnergySpectrum_wt   t_Dr_t   v_AvrLonRad_xvz   v_AvrLon_xv   v_AvrRad_vz   v_IntLonRad_xvz   v_IntLon_xv   v_IntRad_vz   v_Lat   v_Lat_Weight   vz_AvrLon_xvz   vz_IntLon_xvz   w_xv   w_xy   wt_Boundaries   wt_BoundariesGrid   wt_BoundariesTau   wt_DRad_wt   wt_DivLat_xvz   wt_DivLat_xyz   wt_DivLon_xvz   wt_DivLon_xyz   wt_DivRad_wt   wt_Div_xvz_xvz_xvz   wt_Div_xyz_xyz_xyz   wt_KxRGrad_wt   wt_L2Inv_wt   wt_L2_wt   wt_LaplaPol2PolGrid_wt   wt_LaplaPol2PolTau_wt   wt_LaplaPol2Pol_wt   wt_Lapla_wt   wt_PolMagBoundaries   wt_PolmagBoundariesGrid   wt_PolmagBoundariesTau   wt_Potential2Rotation   wt_Potential2RotationMPI   wt_Potential2Vector   wt_Potential2VectorMPI   wt_QOperatorMPI_wt   wt_QOperator_wt   wt_RadRotRot_xvz_xvz_xvz   wt_RadRotRot_xyz_xyz_xyz   wt_RadRot_xvz_xvz   wt_RadRot_xyz_xyz   wt_RotRad_wt   wt_RotRad_xvz_xvz   wt_RotRad_xyz_xyz   wt_TorBoundaries   wt_TorBoundariesGrid   wt_TorBoundariesTau   wt_TorMagBoundaries   wt_TormagBoundariesGrid   wt_TormagBoundariesTau   wt_VGradV   wt_VGradVMPI   wt_VMiss   wt_mpi_Initial   wt_wz   wt_xvz   wt_xyz   wz_LaplaPol2Pol_wz   wz_RAD   wz_wt   wz_xvz   wz_xyz   x_AvrLatRad_xvz   x_AvrLatRad_xyz   x_AvrLat_xv   x_AvrLat_xy   x_AvrRad_xz   x_IntLatRad_xvz   x_IntLatRad_xyz   x_IntLat_xv   x_IntLat_xy   x_IntRad_xz   x_Lon   x_Lon_Weight   xv_AvrRad_xvz   xv_IntRad_xvz   xv_Lat   xv_Lon   xv_w   xvz_Div_xvz_xvz_xvz   xvz_GradLat_wt   xvz_GradLon_wt   xvz_KGrad_wt   xvz_LAT   xvz_LON   xvz_RAD   xvz_RotLat_wt_wt   xvz_RotLon_wt_wt   xvz_wt   xvz_wz   xy_AvrRad_xyz   xy_IntRad_xyz   xy_Lat   xy_Lon   xy_w   xyz_Div_xyz_xyz_xyz   xyz_GradLat_wt   xyz_GradLon_wt   xyz_KGrad_wt   xyz_LAT   xyz_LON   xyz_RAD   xyz_RotLat_wt_wt   xyz_RotLon_wt_wt   xyz_wt   xyz_wz   xz_AvrLat_xvz   xz_AvrLat_xyz   xz_IntLat_xvz   xz_IntLat_xyz   y_AvrLonRad_xyz   y_AvrLon_xy   y_AvrRad_yz   y_IntLonRad_xyz   y_IntLon_xy   y_IntRad_yz   y_Lat   y_Lat_Weight   yz_AvrLon_xyz   yz_IntLon_xyz   z_AvrLat_vz   z_AvrLat_yz   z_AvrLonLat_xvz   z_AvrLonLat_xyz   z_AvrLon_xz   z_IntLat_vz   z_IntLat_yz   z_IntLonLat_xvz   z_IntLonLat_xyz   z_IntLon_xz   z_RAD   z_RAD_WEIGHT  

Included Modules

dc_message lumatrix mpi w_base_mpi_module w_deriv_mpi_module w_integral_mpi_module wa_base_mpi_module wa_deriv_mpi_module wt_module

Public Instance methods

Function :
AvrLatRad_vz :real(8)
: (out) 平均値
vz :real(8), dimension(jc,0:km), intent(in)
: (in) 2 次元緯度動径(子午面)格子点データ

緯度動径(子午面)積分

2 次元(VZ)格子点データの緯度動径(子午面)平均

2 次元データ f(φ,r) に対して

   ∫f(φ,r) r^2cosφ dφdr /(2(r[o]^3-r[i]^3)/3)

を計算する.

[Source]

    function AvrLatRad_vz(vz)  ! 緯度動径(子午面)積分
      !
      ! 2 次元(VZ)格子点データの緯度動径(子午面)平均
      !
      ! 2 次元データ f(φ,r) に対して
      !
      !    ∫f(φ,r) r^2cosφ dφdr /(2(r[o]^3-r[i]^3)/3) 
      !
      ! を計算する.
      !
      real(8), dimension(jc,0:km), intent(in) :: vz
      !(in) 2 次元緯度動径(子午面)格子点データ

      real(8)                   :: AvrLatRad_vz
      !(out) 平均値

      AvrLatRad_vz = IntLatRad_vz(vz)/(sum(y_Lat_Weight)*sum(z_Rad_Weight))

    end function AvrLatRad_vz
AvrLatRad_yz( yz ) result(AvrLatRad_yz)
Function :
AvrLatRad_yz :real(8)
: (out) 平均値
yz :real(8), dimension(1:jm,0:km), intent(in)
: (in) 2 次元緯度動径(子午面)格子点データ

緯度動径(子午面)積分

2 次元(YZ)格子点データの緯度動径(子午面)平均

2 次元データ f(φ,r) に対して

   ∫f(φ,r) r^2cosφ dφdr /(2(r[o]^3-r[i]^3)/3)

を計算する.

Original external subprogram is wt_module#AvrLatRad_yz

AvrLat_v( v_data ) result(AvrLat_v)
Function :
AvrLat_v :real(8)
: (out) 平均値
v_data(jc) :real(8), intent(in)
: (in) 1 次元緯度格子点データ

1 次元(Y)格子点データの緯度(Y)方向平均(1 層用).

実際には格子点データ各点毎に v_Y_Weight をかけた総和を計算し, v_Y_Weight の総和で割ることで平均している.

Original external subprogram is w_integral_mpi_module#AvrLat_v

AvrLat_y( y_data ) result(AvrLat_y)
Function :
AvrLat_y :real(8)
: (out) 平均値
y_data(1:jm) :real(8), intent(in)
: (in) 1 次元緯度格子点データ

1 次元(Y)格子点データの緯度(Y)方向平均(1 層用).

実際には格子点データ各点毎に y_Y_Weight をかけた総和を計算し, y_Y_Weight の総和で割ることで平均している.

Original external subprogram is wt_module#AvrLat_y

Function :
AvrLonLatRad_xvz :real(8)
: (out) 全球平均値
xvz :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

緯度経度動径(全球)積分

3 次元格子点データの緯度経度動径(全球)積分

3 次元データ f(λ,φ,r) に対して

   ∫f(λ,φ,r) r^2cosφ dλdφdr /(4π(r[o]^3-r[i]^3)/3)

を計算する.

[Source]

    function AvrLonLatRad_xvz(xvz) ! 緯度経度動径(全球)積分
      !
      ! 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,jc,0:km), intent(in) :: xvz
      !(in) 3 次元経度緯度動径格子点データ

      real(8)                     :: AvrLonLatRad_xvz
      !(out) 全球平均値

      AvrLonLatRad_xvz = IntLonLatRad_xvz(xvz) /(sum(x_Lon_Weight)*sum(y_Lat_Weight) * sum(z_Rad_Weight))

    end function AvrLonLatRad_xvz
AvrLonLatRad_xyz( xyz ) result(AvrLonLatRad_xyz)
Function :
AvrLonLatRad_xyz :real(8)
: (out) 全球平均値
xyz :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

緯度経度動径(全球)積分

3 次元格子点データの緯度経度動径(全球)積分

3 次元データ f(λ,φ,r) に対して

   ∫f(λ,φ,r) r^2cosφ dλdφdr /(4π(r[o]^3-r[i]^3)/3)

を計算する.

Original external subprogram is wt_module#AvrLonLatRad_xyz

AvrLonLat_xv( xv_data ) result(AvrLonLat_xv)
Function :
AvrLonLat_xv :real(8)
: (out) 平均値
xv_data(0:im-1,jc) :real(8), intent(in)
: (in) 2 次元経度緯度格子点データ(0:im-1,jc)

2 次元緯度経度格子点データの全領域平均(1 層用).

実際には格子点データ各点毎に x_X_Weight, v_Y_Weight をかけた 総和を計算し, x_X_Weight*v_Y_Weight の総和で割ることで平均している.

Original external subprogram is w_integral_mpi_module#AvrLonLat_xv

AvrLonLat_xy( xy_data ) result(AvrLonLat_xy)
Function :
AvrLonLat_xy :real(8)
: (out) 平均値
xy_data(0:im-1,1:jm) :real(8), intent(in)
: (in) 2 次元経度緯度格子点データ

2 次元緯度経度格子点データの全領域平均(1 層用).

実際には格子点データ各点毎に x_X_Weight, y_Y_Weight をかけた 総和を計算し, x_X_Weight*y_Y_Weight の総和で割ることで平均している.

Original external subprogram is wt_module#AvrLonLat_xy

AvrLonRad_xz( xz ) result(AvrLonRad_xz)
Function :
AvrLonRad_xz :real(8)
: 積分値
xz :real(8), dimension(0:im-1,0:km), intent(in)
: (in) 2 次元格子点データ

経度動径(緯度円)積分

2 次元(XZ)格子点データの経度動径平均

2 次元データ f(λ,r) に対して

   ∫f(λ,r) r^2dλdr /(2π(r[o]^3-r[i]^3)/3)

を計算する.

Original external subprogram is wt_module#AvrLonRad_xz

AvrLon_x( x_data ) result(AvrLon_x)
Function :
AvrLon_x :real(8)
: (out) 平均値
x_data(0:im-1) :real(8), intent(in)
: (in) 1 次元経度(X)格子点データ

1 次元(X)格子点データの経度(X)方向平均(1 層用).

実際には格子点データ各点毎に x_X_Weight をかけた総和を計算し, x_X_Weight の総和で割ることで平均している.

Original external subprogram is wt_module#AvrLon_x

AvrRad_z( z ) result(AvrRad_z)
Function :
AvrRad_z :real(8)
: (out) 平均値
z :real(8), dimension(0:km), intent(in)
: (in) 1 次元動径格子点データ

1 次元(Z)格子点データの動径方向域平均.

1 次元データ f(r) に対して ∫f(r) r^2dr /((r[o]^3-r[i]^3)/3) を 計算する.

Original external subprogram is wt_module#AvrRad_z

Function :
IntLatRad_vz :real(8)
: (out) 積分値
vz :real(8), dimension(jc,0:km), intent(in)
: (in) 2 次元緯度動径(子午面)格子点データ

2 次元(VZ)格子点データの緯度動径積分(子午面)および平均

2 次元データ f(φ,r) に対して ∫f(φ,r) r^2cosφ dφdr を計算する.

[Source]

    function IntLatRad_vz(vz)
      !
      ! 2 次元(VZ)格子点データの緯度動径積分(子午面)および平均
      !
      ! 2 次元データ f(φ,r) に対して ∫f(φ,r) r^2cosφ dφdr を計算する.
      !
      real(8), dimension(jc,0:km), intent(in) :: vz
      !(in) 2 次元緯度動径(子午面)格子点データ

      real(8)                   :: IntLatRad_vz
      !(out) 積分値

      real(8)                   :: IntLatRadTMP
      integer :: j, k

      IntLatRad_vz = 0
      do k=0,km
         do j=1,jc
            IntLatRad_vz = IntLatRad_vz + vz(j,k) * v_Lat_Weight(j) * z_Rad_Weight(k)
         enddo
      enddo

      IntLatRadTmp=IntLatRad_vz
      CALL MPI_ALLREDUCE(IntLatRadTMP,IntLatRad_vz,1,MPI_REAL8, MPI_SUM,MPI_COMM_WORLD,IERR)

    end function IntLatRad_vz
IntLatRad_yz( yz ) result(IntLatRad_yz)
Function :
IntLatRad_yz :real(8)
: (out) 積分値
yz :real(8), dimension(1:jm,0:km), intent(in)
: (in) 2 次元緯度動径(子午面)格子点データ

2 次元(YZ)格子点データの緯度動径積分(子午面)および平均

2 次元データ f(φ,r) に対して ∫f(φ,r) r^2cosφ dφdr を計算する.

Original external subprogram is wt_module#IntLatRad_yz

IntLat_v( v_data ) result(IntLat_v)
Function :
IntLat_v :real(8)
: (out) 積分値
v_data(jc) :real(8), intent(in)
: (in) 1 次元緯度(Y)格子点データ

1 次元緯度(Y)格子点データの Y 方向積分(1 層用).

実際には格子点データ各点毎に v_Y_Weight をかけた総和を計算している.

Original external subprogram is w_integral_mpi_module#IntLat_v

IntLat_y( y_data ) result(IntLat_y)
Function :
IntLat_y :real(8)
: (out) 積分値
y_data(1:jm) :real(8), intent(in)
: (in) 1 次元緯度(Y)格子点データ

1 次元緯度(Y)格子点データの Y 方向積分(1 層用).

実際には格子点データ各点毎に y_Y_Weight をかけた総和を計算している.

Original external subprogram is wt_module#IntLat_y

Function :
IntLonLatRad_xvz :real(8)
: (out) 全球積分値
xvz :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

緯度経度動径(全球)積分

3 次元格子点データの緯度経度動径(全球)積分

3 次元データ f(λ,φ,r) に対して

    ∫f(λ,φ,r) r^2cosφ dλdφdr

を計算する.

[Source]

    function IntLonLatRad_xvz(xvz) ! 緯度経度動径(全球)積分
      !
      ! 3 次元格子点データの緯度経度動径(全球)積分
      !
      ! 3 次元データ f(λ,φ,r) に対して
      !
      !     ∫f(λ,φ,r) r^2cosφ dλdφdr 
      !
      ! を計算する.
      !
      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz 
      !(in) 3 次元経度緯度動径格子点データ

      real(8)                     :: IntLonLatRad_xvz 
      !(out) 全球積分値

      real(8)                     :: IntLonLatRadTMP
      integer :: i, j, k

      IntLonLatRad_xvz = 0
      do k=0,km
         do j=1,jc
            do i=0,im-1
               IntLonLatRad_xvz = IntLonLatRad_xvz + xvz(i,j,k) * x_Lon_Weight(i) * v_Lat_Weight(j) * z_Rad_Weight(k)
            enddo
         enddo
      enddo

      IntLonLatRadTmp=IntLonLatRad_xvz
      CALL MPI_ALLREDUCE(IntLonLatRadTMP,IntLonLatRad_xvz,1,MPI_REAL8, MPI_SUM,MPI_COMM_WORLD,IERR)

    end function IntLonLatRad_xvz
IntLonLatRad_xyz( xyz ) result(IntLonLatRad_xyz)
Function :
IntLonLatRad_xyz :real(8)
: (out) 全球積分値
xyz :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

緯度経度動径(全球)積分

3 次元格子点データの緯度経度動径(全球)積分

3 次元データ f(λ,φ,r) に対して

    ∫f(λ,φ,r) r^2cosφ dλdφdr

を計算する.

Original external subprogram is wt_module#IntLonLatRad_xyz

IntLonLat_xv( xv_data ) result(IntLonLat_xv)
Function :
IntLonLat_xv :real(8)
: (out) 積分値
xv_data(0:im-1,jc) :real(8), intent(in)
: (in) 2 次元経度緯度格子点データ(0:im-1,jc)

2 次元緯度経度格子点データの全領域積分(1 層用).

実際には格子点データ各点毎に x_X_Weight, v_V_Weight をかけた 総和を計算している.

Original external subprogram is w_integral_mpi_module#IntLonLat_xv

IntLonLat_xy( xy_data ) result(IntLonLat_xy)
Function :
IntLonLat_xy :real(8)
: (out) 積分値
xy_data(0:im-1,1:jm) :real(8), intent(in)
: (in) 2 次元経度緯度格子点データ(0:im-1,1:jm)

2 次元緯度経度格子点データの全領域積分(1 層用).

実際には格子点データ各点毎に x_X_Weight, y_Y_Weight をかけた 総和を計算している.

Original external subprogram is wt_module#IntLonLat_xy

IntLonRad_xz( xz ) result(IntLonRad_xz)
Function :
IntLonRad_xz :real(8)
: (out) 積分値
xz :real(8), dimension(0:im-1,0:km), intent(in)
: (in) 2 次元緯度動径格子点データ

経度動径(緯度円)積分

2 次元(XZ)格子点データの経度動径積分

2 次元データ f(λ,r) に対して∫f(λ,r) r^2dλdr を計算する.

Original external subprogram is wt_module#IntLonRad_xz

IntLon_x( x_data ) result(IntLon_x)
Function :
IntLon_x :real(8)
: (out) 積分値
x_data(0:im-1) :real(8), intent(in)
: (in) 1 次元経度(X)格子点データ

1 次元経度(X)格子点データの X 方向積分(1 層用).

実際には格子点データ各点毎に x_X_Weight をかけた総和を計算している.

Original external subprogram is wt_module#IntLon_x

IntRad_z( z ) result(IntRad_z)
Function :
IntRad_z :real(8)
: (out) 積分値
z :real(8), dimension(0:km), intent(in)
: (in) 1 次元動径格子点データ

動径積分

1 次元(Z)格子点データの動径方向域積分.

1 次元データ f(r) に対して ∫f(r) r^2dr を計算する.

Original external subprogram is wt_module#IntRad_z

Interpolate_wt( wt_data, alon, alat, arad ) result(Interpolate_wt)
Function :
Interpolate_wt :real(8)
: 補間した値
wt_data((nm+1)**2,0:km) :real(8), intent(IN)
: スペクトルデータ
alon :real(8), intent(IN)
: 補間する位置(経度)
alat :real(8), intent(IN)
: 補間する位置(緯度)
arad :real(8), intent(IN)
: 補間する位置(動径)

緯度 alon, 経度 alat 動径 arad における関数値を その球面調和変換係数 wa_data から補間計算する

Original external subprogram is wt_module#Interpolate_wt

at_Dr_at( at_data ) result(at_Dx_at)
Function :
at_Dx_at :real(8), dimension(size(at_data,1),0:size(at_data,2)-1)
: (out) チェビシェフデータの X 微分
at_data :real(8), dimension(:,0:), intent(in)
: (in) 入力チェビシェフデータ

入力チェビシェフデータに X 微分を作用する(2 次元配列用).

チェビシェフデータの X 微分とは, 対応する格子点データに X 微分を 作用させたデータのチェビシェフ変換のことである.

Original external subprogram is wt_module#at_Dr_at

at_az( ag_data ) result(at_ag)
Function :
at_ag :real(8), dimension(size(ag_data,1),0:km)
: (out) チェビシェフデータ
ag_data :real(8), dimension(:,:), intent(in)
: (in) 格子点データ

格子データからチェビシェフデータへ変換する(2 次元配列用).

Original external subprogram is wt_module#at_az

az_at( at_data ) result(ag_at)
Function :
ag_at :real(8), dimension(size(at_data,1),0:im)
: (out) 格子点データ
at_data :real(8), dimension(:,:), intent(in)
: (in) チェビシェフデータ

チェビシェフデータから格子データへ変換する(2 次元配列用).

Original external subprogram is wt_module#az_at

jc
Variable :
jc :integer
: 分散配置用変数

Original external subprogram is w_base_mpi_module#jc

l_nm( n, m ) result(l_nm_array00)
Function :
l_nm_array00 :integer
: (out) スペクトルデータの格納位置
n :integer, intent(in)
: (in) 全波数
m :integer, intent(in)
: (in) 帯状波数

全波数(n)と東西波数(m)からそのスペクトルデータの格納位置を返す.

引数 n,m がともに整数値の場合, 整数値を返す.

Original external subprogram is wt_module#l_nm

l_nm( n, marray ) result(l_nm_array01)
Function :
l_nm_array01(size(marray)) :integer
: (out) スペクトルデータ位置
n :integer, intent(in)
: (in) 全波数
marray(:) :integer, intent(in)
: (in) 帯状波数

スペクトルデータの格納位置

全波数(n)と東西波数(m)からそのスペクトルデータの格納位置を返す.

第 1 引数 n が整数, 第 2 引数 marray が整数 1 次元配列の場合, marray と同じ大きさの 1 次元整数配列を返す.

Original external subprogram is wt_module#l_nm

l_nm( narray, m ) result(l_nm_array10)
Function :
l_nm_array10(size(narray)) :integer
: (out) スペクトルデータ位置
narray(:) :integer, intent(in)
: (in) 全波数
m :integer, intent(in)
: (in) 帯状波数

全波数(n)と東西波数(m)からそのスペクトルデータの格納位置を返す.

第 1 引数 narray が整数 1 次元配列, 第 2 引数 m が整数の場合, narray と同じ大きさの 1 次元整数配列を返す.

Original external subprogram is wt_module#l_nm

l_nm( narray, marray ) result(l_nm_array11)
Function :
l_nm_array11(size(narray)) :integer
: (out) スペクトルデータ位置
narray(:) :integer, intent(in)
: (in) 全波数
marray(:) :integer, intent(in)
: (in) 帯状波数

全波数(n)と東西波数(m)からそのスペクトルデータの格納位置を返す.

第 1,2 引数 narray, marray がともに整数 1 次元配列の場合, narray, marray と同じ大きさの 1 次元整数配列を返す. narray, marray は同じ大きさでなければならない.

Original external subprogram is wt_module#l_nm

nm_l( l ) result(nm_l_int)
Function :
nm_l_int(2) :integer
: (out) 全波数, 帯状波数
l :integer, intent(in)
: (in) スペクトルデータの格納位置

スペクトルデータの格納位置(l)から全波数(n)と東西波数(m)を返す.

引数 l が整数値の場合, 対応する全波数と帯状波数を 長さ 2 の 1 次元整数値を返す. nm_l(1) が全波数, nm_l(2) が帯状波数である.

Original external subprogram is wt_module#nm_l

nm_l( larray ) result(nm_l_array)
Function :
nm_l_array(size(larray),2) :integer
: (in) スペクトルデータの格納位置
larray(:) :integer, intent(in)
: (out) 全波数, 帯状波数

スペクトルデータの格納位置(l)から全波数(n)と東西波数(m)を返す.

引数 larray が整数 1 次元配列の場合, larray に対応する n, m を格納した 2 次元整数配列を返す. nm_l_array(:,1) が全波数, nm_l_array(:,2) が帯状波数である.

Original external subprogram is wt_module#nm_l

nmz_PoloidalEnergySpectrum_wt( wt_POLPOT ) result(nmz_PoloidalEnergySpectrum_wt)
Function :
nmz_PoloidalEnergySpectrum_wt :real(8), dimension(0:nm,-nm:nm,0:km)
: (out) エネルギースペクトルポロイダル成分
wt_POLPOT :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) ポロイダルポテンシャル

ポロイダルポテンシャルから, ポロイダルエネルギーの 球面調和函数全波数 n, 帯状波数 m の各成分を計算する.

 * 全波数 n, 帯状波数 m のポロイダルポテンシャルのスペクトル成分
   φ(n,m,r)から全波数 n, 帯状波数 m 成分のポロイダルエネルギー
   スペクトルは

     (1/2)n(n+1)4πr^2{[d(rφ(n,m,r))/dr]^2 + n(n+1)φ(n,m,r)^2}

   と計算される.

 * 全てのエネルギースペクトル成分の和を動径積分したもの(r^2の重み無し)
   が球殻内での全エネルギーに等しい.

 * データの存在しない全波数 n, 帯状波数 m の配列には欠損値が格納される.
   欠損値の値はモジュール変数 wt_VMiss によって設定できる
   (初期値は -999.0)

Original external subprogram is wt_module#nmz_PoloidalEnergySpectrum_wt

nmz_ToroidalEnergySpectrum_wt( wt_TORPOT ) result(nmz_ToroidalEnergySpectrum_wt)
Function :
nmz_ToroidalEnergySpectrum_wt :real(8), dimension(0:nm,-nm:nm,0:km)
: (out) エネルギースペクトルトロイダル成分
wt_TORPOT :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (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 の配列には欠損値が格納される.
   wt_VMiss によって設定できる (初期値は -999.0)

Original external subprogram is wt_module#nmz_ToroidalEnergySpectrum_wt

nz_PoloidalEnergySpectrum_wt( wt_POLPOT ) result(nz_PoloidalEnergySpectrum_wt)
Function :
nz_PoloidalEnergySpectrum_wt :real(8), dimension(0:nm,0:km)
: (out) エネルギースペクトルポロイダル成分
wt_POLPOT :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) ポロイダルポテンシャル

ポロイダルポテンシャルから, ポロイダルエネルギーの 球面調和函数全波数の各成分を計算する

 * 全波数 n, 帯状波数 m のポロイダルポテンシャルのスペクトル成分
   φ(n,m,r)から全波数 n 成分のポロイダルエネルギースペクトルは

     Σ[m=-n]^n ((1/2)n(n+1)4πr^2{[d(rφ(n,m,r))/dr]^2
                + n(n+1)φ(n,m,r)^2}

   と計算される.

 * 全ての全波数に対してのエネルギースペクトル成分の和を動径積分したもの
   (r^2の重み無し)が球殻内での全エネルギーに等しい.

Original external subprogram is wt_module#nz_PoloidalEnergySpectrum_wt

nz_ToroidalEnergySpectrum_wt( wt_TORPOT ) result(nz_ToroidalEnergySpectrum_wt)
Function :
nz_ToroidalEnergySpectrum_wt :real(8), dimension(0:nm,0:km)
: (out) エネルギースペクトルトロイダル成分
wt_TORPOT :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) トロイダルポテンシャル

トロイダルポテンシャルから, トロイダルエネルギーの 球面調和函数全波数の各成分を計算する.

 * 全波数 n, 帯状波数 m のトロイダルポテンシャルのスペクトル成分
   ψ(n,m,r)から全波数 n 成分のトロイダルエネルギースペクトルは
   Σ[m=-n]^n(1/2)n(n+1)4πr^2ψ(n,m,r)^2 と計算される.
  • 全てのエネルギースペクトル成分の和を動径積分したもの(r^2の重み無し)
     が球殻内での全エネルギーに等しい.
    

Original external subprogram is wt_module#nz_ToroidalEnergySpectrum_wt

t_Dr_t( t_data ) result(t_Dx_t)
Function :
t_Dx_t :real(8), dimension(size(t_data))
: (out) チェビシェフデータの X 微分
t_data :real(8), dimension(:), intent(in)
: (in) 入力チェビシェフデータ

入力チェビシェフデータに X 微分を作用する(1 次元配列用).

チェビシェフデータの X 微分とは, 対応する格子点データに X 微分を 作用させたデータのチェビシェフ変換のことである.

Original external subprogram is wt_module#t_Dr_t

Function :
v_AvrLonRad_xvz :real(8), dimension(jc)
: (out) 経度動径(緯度円)平均された 1 次元緯度格子点データ
xvz :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

経度動径(緯度円)積分

3 次元格子点データの経度動径(緯度円)平均.

3 次元データ f(λ,φ,r) に対して

    ∫f(λ,φ,r) r^2dλdr /(2π(r[o]^3-r[i]^3)/3)

を計算する.

[Source]

    function v_AvrLonRad_xvz(xvz)  ! 経度動径(緯度円)積分
      !
      ! 3 次元格子点データの経度動径(緯度円)平均.
      !
      ! 3 次元データ f(λ,φ,r) に対して
      !
      !     ∫f(λ,φ,r) r^2dλdr /(2π(r[o]^3-r[i]^3)/3) 
      !
      ! を計算する.
      !
      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz
      !(in) 3 次元経度緯度動径格子点データ

      real(8), dimension(jc)       :: v_AvrLonRad_xvz
      !(out) 経度動径(緯度円)平均された 1 次元緯度格子点データ

      v_AvrLonRad_xvz = v_IntLonRad_xvz(xvz) /(sum(x_Lon_Weight)*sum(z_Rad_Weight))

    end function v_AvrLonRad_xvz
v_AvrLon_xv( xv_data ) result(v_AvrLon_xv)
Function :
v_AvrLon_xv(jc) :real(8)
: (out) 平均された 1 次元緯度(Y)格子点
xv_data(0:im-1,jc) :real(8), intent(in)
: (in) 2 次元経度緯度格子点データ(0:im-1,jc)

2 次元緯度経度格子点データの経度(X)方向平均(1 層用).

実際には格子点データ各点毎に x_X_Weight をかけた総和を計算し, x_X_Weight の総和で割ることで平均している.

Original external subprogram is w_integral_mpi_module#v_AvrLon_xv

Function :
v_AvrRad_vz :real(8), dimension(jc)
: (out) 動径平均された 1 次元緯度格子点データ
vz :real(8), dimension(jc,0:km), intent(in)
: (in) 2 次元緯度動径(子午面)格子点データ

2 次元(VZ)格子点データの動径方向域平均.

2 次元データ f(φ,r) に対して ∫f(φ,r) r^2dr /((r[o]^3-r[i]^3)/3) を計算する.

[Source]

    function v_AvrRad_vz(vz)
      !
      ! 2 次元(VZ)格子点データの動径方向域平均.
      !
      ! 2 次元データ f(φ,r) に対して ∫f(φ,r) r^2dr /((r[o]^3-r[i]^3)/3) 
      ! を計算する.
      !
      real(8), dimension(jc,0:km), intent(in) :: vz
      !(in) 2 次元緯度動径(子午面)格子点データ

      real(8), dimension(jc)  :: v_AvrRad_vz
      !(out) 動径平均された 1 次元緯度格子点データ

      v_AvrRad_vz = v_IntRad_vz(vz)/sum(z_Rad_Weight)

    end function v_AvrRad_vz
Function :
v_IntLonRad_xvz :real(8), dimension(jc)
: (out) 経度動径(緯度円)積分された 1 次元緯度格子点データ
xvz :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

3 次元格子点データの経度動径(緯度円)積分.

3 次元データ f(λ,φ,r) に対して∫f(λ,φ,r) r^2dλdr を計算する.

[Source]

    function v_IntLonRad_xvz(xvz)
      !
      ! 3 次元格子点データの経度動径(緯度円)積分.
      !
      ! 3 次元データ f(λ,φ,r) に対して∫f(λ,φ,r) r^2dλdr を計算する.
      !
      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz
      !(in) 3 次元経度緯度動径格子点データ

      real(8), dimension(jc)       :: v_IntLonRad_xvz
      !(out) 経度動径(緯度円)積分された 1 次元緯度格子点データ

      integer :: i, k

      v_IntLonRad_xvz = 0
      do k=0,km
         do i=0,im-1
            v_IntLonRad_xvz = v_IntLonRad_xvz + xvz(i,:,k) * x_Lon_Weight(i) * z_Rad_Weight(k)
         enddo
      enddo
    end function v_IntLonRad_xvz
v_IntLon_xv( xv_data ) result(v_IntLon_xv)
Function :
v_IntLon_xv(jc) :real(8)
: (out) 積分された 1 次元緯度(Y)格子点データ
xv_data(0:im-1,jc) :real(8), intent(in)
: (in) 2 次元経度緯度格子点データ(0:im-1,jc)

2 次元緯度経度格子点データの経度(X)方向積分(1 層用).

実際には格子点データ各点毎に x_X_Weight をかけた総和を計算している.

Original external subprogram is w_integral_mpi_module#v_IntLon_xv

Function :
v_IntRad_vz :real(8), dimension(jc)
: (out) 動径積分された 1 次元緯度格子点データ
vz :real(8), dimension(jc,0:km), intent(in)
: (in) 2 次元緯度動径(子午面)格子点データ

動径積分

2 次元(VZ)格子点データの動径方向域積分.

2 次元データ f(φ,r) に対して∫f(φ,r) r^2dr を計算する.

[Source]

    function v_IntRad_vz(vz)  ! 動径積分
      !
      ! 2 次元(VZ)格子点データの動径方向域積分.
      !
      ! 2 次元データ f(φ,r) に対して∫f(φ,r) r^2dr を計算する.
      !
      real(8), dimension(jc,0:km), intent(in) :: vz
      !(in) 2 次元緯度動径(子午面)格子点データ

      real(8), dimension(jc)  :: v_IntRad_vz
      !(out) 動径積分された 1 次元緯度格子点データ

      integer :: k

      v_IntRad_vz = 0.0d0
      do k=0,km
         v_IntRad_vz(:) = v_IntRad_vz(:) + vz(:,k) * z_Rad_Weight(k) 
      enddo
    end function v_IntRad_vz
v_Lat
Variable :
v_Lat(:) :real(8), allocatable
: 緯度経度

Original external subprogram is w_base_mpi_module#v_Lat

v_Lat_Weight
Variable :
v_Lat_Weight(:) :real(8), allocatable
: 緯度経度

Original external subprogram is w_base_mpi_module#v_Lat_Weight

Function :
vz_AvrLon_xvz :real(8), dimension(jc,0:km)
: (out) 経度方向(帯状)平均された 2 次元子午面格子点データ
xvz :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

経度(帯状)積分

3 次元格子点データの経度方向(帯状)平均.

3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r)dλ/2π を計算する.

[Source]

    function vz_AvrLon_xvz(xvz)  ! 経度(帯状)積分
      !
      ! 3 次元格子点データの経度方向(帯状)平均.
      !
      ! 3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r)dλ/2π を計算する.
      !
      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz
      !(in) 3 次元経度緯度動径格子点データ

      real(8), dimension(jc,0:km)  :: vz_AvrLon_xvz
      !(out) 経度方向(帯状)平均された 2 次元子午面格子点データ

      vz_AvrLon_xvz = vz_IntLon_xvz(xvz)/sum(x_Lon_Weight)

    end function vz_AvrLon_xvz
Function :
vz_IntLon_xvz :real(8), dimension(jc,0:km)
: (out) 経度方向(帯状)積分された 2 次元子午面格子点データ
xvz :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

経度(帯状)積分

3 次元格子点データの経度方向(帯状)積分.

3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r)dλ を計算する.

[Source]

    function vz_IntLon_xvz(xvz)  ! 経度(帯状)積分
      !
      ! 3 次元格子点データの経度方向(帯状)積分.
      !
      ! 3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r)dλ を計算する.
      !
      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz
      !(in) 3 次元経度緯度動径格子点データ

      real(8), dimension(jc,0:km)  :: vz_IntLon_xvz
      !(out) 経度方向(帯状)積分された 2 次元子午面格子点データ

      integer :: i

      vz_IntLon_xvz = 0.0d0
      do i=0,im-1
         vz_IntLon_xvz(:,:) = vz_IntLon_xvz(:,:) + xvz(i,:,:) * x_Lon_Weight(i)
      enddo
    end function vz_IntLon_xvz
w_xv( xv_data, [ipow], [iflag] ) result(w_xv)
Function :
w_xv((nm+1)*(nm+1)) :real(8)
: (out) スペクトルデータ
xv_data(0:im-1,jc) :real(8), intent(in)
: (in) 格子点データ
ipow :integer, intent(in), optional
: (in) 変換時に同時に作用させる 1/cosφ の次数. 省略時は 0.
iflag :integer, intent(in), optional
: 変換の種類
    0 : 通常の正変換
    1 : 経度微分を作用させた正変換
   -1 : 緯度微分を作用させた正変換
    2 : sinφを作用させた正変換
  省略時は 0.

格子データからスペクトルデータへ(正)変換する(1 層用).

Original external subprogram is w_base_mpi_module#w_xv

w_xy( xy_data, [ipow], [iflag] ) result(w_xy)
Function :
w_xy((nm+1)*(nm+1)) :real(8)
: (out) スペクトルデータ
xy_data(0:im-1,1:jm) :real(8), intent(in)
: (in) 格子点データ
ipow :integer, intent(in), optional
: (in) 変換時に同時に作用させる 1/cosφ の次数. 省略時は 0.
iflag :integer, intent(in), optional
: 変換の種類
   0 : 通常の正変換
  -1 : 経度微分を作用させた正変換
   1 : 緯度微分 1/cosφ・∂(f cos^2φ)/∂φ を作用させた正変換
   2 : sinφを作用させた正変換
 省略時は 0.

格子データからスペクトルデータへ(正)変換する(1 層用).

Original external subprogram is wt_module#w_xy

wt_Boundaries( wt, [values], [cond] )
Subroutine :
wt :real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
: (inout) 境界条件を適用するデータ. 修正された値を返す.
values :real(8), dimension((nm+1)*(nm+1),2), intent(in), optional
: (in) 境界での 値/勾配 分布を水平スペクトル変換したものを与える.
    省略時は値/勾配 0 となる.
cond :character(len=2), intent(in), optional
: (in) 境界条件. 省略時は ‘DD‘
       DD : 両端ディリクレ条件
       DN : 上端ディリクレ, 下端ノイマン条件
       ND : 上端ノイマン, 下端ディリクレ条件
       NN : 両端ノイマン条件

スペクトルデータにディリクレ・ノイマン境界条件を適用する Chebyshev 空間での境界条件適用(タウ法)

チェビシェフ空間において境界条件を満たすべく高次の係数を 定める方法をとっている(タウ法).

Original external subprogram is wt_module#wt_Boundaries

wt_BoundariesGrid( wt, [values], [cond] )
Subroutine :
wt :real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
: (inout) 境界条件を適用するデータ. 修正された値を返す.
values :real(8), dimension((nm+1)*(nm+1),2), intent(in), optional
: (in) 境界での 値/勾配 分布を水平スペクトル変換したものを与える.
   省略時は値/勾配 0 となる.
cond :character(len=2), intent(in), optional
: (in) 境界条件. 省略時は ‘DD‘
       DD : 両端ディリクレ条件
       DN : 上端ディリクレ, 下端ノイマン条件
       ND : 上端ノイマン, 下端ディリクレ条件
       NN : 両端ノイマン条件

スペクトルデータにディリクレ・ノイマン境界条件を適用する 実空間での境界条件適用

鉛直実格子点空間において内部領域の値と境界条件を満たすように 条件を課している(選点法). このルーチンを用いるためには wt_Initial にて設定するチェビシェフ切断波数(lm)と鉛直格子点数(km)を 等しくしておく必要がある.

Original external subprogram is wt_module#wt_BoundariesGrid

wt_BoundariesTau( wt, [values], [cond] )
Subroutine :
wt :real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
: (inout) 境界条件を適用するデータ. 修正された値を返す.
values :real(8), dimension((nm+1)*(nm+1),2), intent(in), optional
: (in) 境界での 値/勾配 分布を水平スペクトル変換したものを与える.
    省略時は値/勾配 0 となる.
cond :character(len=2), intent(in), optional
: (in) 境界条件. 省略時は ‘DD‘
       DD : 両端ディリクレ条件
       DN : 上端ディリクレ, 下端ノイマン条件
       ND : 上端ノイマン, 下端ディリクレ条件
       NN : 両端ノイマン条件

スペクトルデータにディリクレ・ノイマン境界条件を適用する Chebyshev 空間での境界条件適用(タウ法)

チェビシェフ空間において境界条件を満たすべく高次の係数を 定める方法をとっている(タウ法).

Original external subprogram is wt_module#wt_BoundariesTau

wt_DRad_wt( wt ) result(wt_DRad_wt)
Function :
wt_DRad_wt :real(8), dimension((nm+1)*(nm+1),0:lm)
: (in) 動径微分された2 次元球面調和函数チェビシェフスペクトルデータ
wt :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) 2 次元球面調和函数チェビシェフスペクトルデータ

入力スペクトルデータに動径微分 ∂/∂r を作用する.

スペクトルデータの動径微分とは, 対応する格子点データに動径微分を 作用させたデータのスペクトル変換のことである.

Original external subprogram is wt_module#wt_DRad_wt

Function :
wt_DivLat_xvz :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) 発散型緯度微分を作用された 2 次元スペクトルデータ
xvz :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

格子データに発散型緯度微分 1/rcosφ・∂(f cosφ)/∂φ を 作用させたスペクトルデータを返す.

[Source]

    function wt_DivLat_xvz(xvz)
      !
      ! 格子データに発散型緯度微分 1/rcosφ・∂(f cosφ)/∂φ を
      ! 作用させたスペクトルデータを返す.
      !
      real(8), dimension(0:im-1,jc,0:km), intent(in)   :: xvz
      !(in) 3 次元経度緯度動径格子点データ

      real(8), dimension((nm+1)*(nm+1),0:lm)       :: wt_DivLat_xvz
      !(out) 発散型緯度微分を作用された 2 次元スペクトルデータ

      wt_DivLat_xvz = wt_wz(wa_divlat_xva(xvz/xvz_Rad))
    end function wt_DivLat_xvz
wt_DivLat_xyz( xyz ) result(wt_DivLat_xyz)
Function :
wt_DivLat_xyz :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) 発散型緯度微分を作用された 2 次元スペクトルデータ
xyz :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

格子データに発散型緯度微分 1/rcosφ・∂(f cosφ)/∂φ を 作用させたスペクトルデータを返す.

Original external subprogram is wt_module#wt_DivLat_xyz

Function :
wt_DivLon_xvz :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) 発散型経度微分を作用された 2 次元スペクトルデータ
xvz :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

格子点データに発散型経度微分 1/rcosφ・∂/∂λ を作用させた スペクトルデータを返す.

[Source]

    function wt_DivLon_xvz(xvz)
      ! 
      ! 格子点データに発散型経度微分 1/rcosφ・∂/∂λ を作用させた
      ! スペクトルデータを返す.
      !
      real(8), dimension(0:im-1,jc,0:km), intent(in)   :: xvz
      !(in) 3 次元経度緯度動径格子点データ

      real(8), dimension((nm+1)*(nm+1),0:lm)       :: wt_DivLon_xvz
      !(out) 発散型経度微分を作用された 2 次元スペクトルデータ

      wt_DivLon_xvz = wt_wz(wa_DivLon_xva(xvz/xvz_Rad))
    end function wt_DivLon_xvz
wt_DivLon_xyz( xyz ) result(wt_DivLon_xyz)
Function :
wt_DivLon_xyz :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) 発散型経度微分を作用された 2 次元スペクトルデータ
xyz :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

格子点データに発散型経度微分 1/rcosφ・∂/∂λ を作用させた スペクトルデータを返す.

Original external subprogram is wt_module#wt_DivLon_xyz

wt_DivRad_wt( wt ) result(wt_DivRad_wt)
Function :
wt_DivRad_wt :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) 発散型動径微分を作用された 2 次元スペクトルデータ
wt :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) 2 次元球面調和函数チェビシェフスペクトルデータ

入力スペクトルデータに発散型動径微分

      1/r^2 ∂/∂r (r^2 .)= ∂/∂r + 2/r

を作用する.

スペクトルデータの発散型動径微分とは, 対応する格子点データに 発散型動径微分を作用させたデータのスペクトル変換のことである.

Original external subprogram is wt_module#wt_DivRad_wt

Function :
wt_Div_xvz_xvz_xvz :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) ベクトル場の発散
xvz_Vlon :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) ベクトル場の経度成分
xvz_Vlat :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) ベクトル場の緯度成分
xvz_Vrad :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) ベクトル場の動径成分

べクトル成分である 3 つの格子データに発散を作用させた スペクトルデータを返す.

第 1, 2 ,3 引数(u,v,w)がそれぞれベクトルの経度成分, 緯度成分, 動径成分を表し, 発散は

     1/rcosφ・∂u/∂λ + 1/rcosφ・∂(v cosφ)/∂φ
   + 1/r^2 ∂/∂r (r^2 w)

と計算される.

[Source]

    function wt_Div_xvz_xvz_xvz(xvz_Vlon,xvz_Vlat,xvz_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,jc,0:km), intent(in) :: xvz_Vlon
      !(in) ベクトル場の経度成分

      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz_Vlat
      !(in) ベクトル場の緯度成分

      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz_Vrad
      !(in) ベクトル場の動径成分

      real(8), dimension((nm+1)*(nm+1),0:lm)     :: wt_Div_xvz_xvz_xvz
      !(out) ベクトル場の発散

      wt_Div_xvz_xvz_xvz =   wt_DivLon_xvz(xvz_Vlon) + wt_DivLat_xvz(xvz_Vlat) + wt_DivRad_wt(wt_xvz(xvz_Vrad))

    end function wt_Div_xvz_xvz_xvz
wt_Div_xyz_xyz_xyz( xyz_Vlon, xyz_Vlat, xyz_Vrad ) result(wt_Div_xyz_xyz_xyz)
Function :
wt_Div_xyz_xyz_xyz :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) ベクトル場の発散
xyz_Vlon :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) ベクトル場の経度成分
xyz_Vlat :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) ベクトル場の緯度成分
xyz_Vrad :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) ベクトル場の動径成分

べクトル成分である 3 つの格子データに発散を作用させた スペクトルデータを返す.

第 1, 2 ,3 引数(u,v,w)がそれぞれベクトルの経度成分, 緯度成分, 動径成分を表し, 発散は

     1/rcosφ・∂u/∂λ + 1/rcosφ・∂(v cosφ)/∂φ
   + 1/r^2 ∂/∂r (r^2 w)

と計算される.

Original external subprogram is wt_module#wt_Div_xyz_xyz_xyz

wt_KxRGrad_wt( wt ) result(wt_KxRGrad_wt)
Function :
wt_KxRGrad_wt :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) 経度微分を作用された 2 次元スペクトルデータ
wt :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) 2 次元球面調和函数チェビシェフスペクトルデータ

入力スペクトルデータに経度微分 k×r・▽ = ∂/∂λを作用する.

Original external subprogram is wt_module#wt_KxRGrad_wt

wt_L2Inv_wt( wt ) result(wt_L2Inv_wt)
Function :
wt_L2Inv_wt :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) L^2 演算子の逆演算を作用された 2 次元スペクトルデータ
wt :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) 2 次元球面調和函数チェビシェフスペクトルデータ

入力スペクトルデータに L^2 演算子の逆演算(-逆水平ラプラシアン)を 作用する.

スペクトルデータに L^2 演算子を作用させる関数 wt_L2_wt の逆計算を 行う関数である.

Original external subprogram is wt_module#wt_L2Inv_wt

wt_L2_wt( wt ) result(wt_L2_wt)
Function :
wt_L2_wt :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) L^2 演算子を作用された 2 次元スペクトルデータ
wt :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) 2 次元球面調和函数チェビシェフスペクトルデータ

入力スペクトルデータに L^2 演算子(=-水平ラプラシアン)を作用する.

L^2 演算子は単位球面上の水平ラプラシアンの逆符号にあたる.

 入力スペクトルデ ータに対応する格子点データに演算子

    L^2 = -1/cos^2φ・∂^2/∂λ^2 - 1/cosφ・∂/∂φ(cosφ∂/∂φ)

を作用させたデータのスペクトル変換が返される.

Original external subprogram is wt_module#wt_L2_wt

wt_LaplaPol2PolGrid_wt( wt, [cond], [new] ) result(wt_LaplaPol2PolGrid_wt)
Function :
wt_LaplaPol2PolGrid_wt :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) 出力ポロイダルポテンシャル分布
wt :real(8), dimension((nm+1)*(nm+1),0:lm),intent(in)
: (in) 入力▽^2φ分布
cond :character(len=2), intent(in), optional
: (in) 境界条件スイッチ. 省略時は ‘RR‘
    RR    : 両端粘着条件
    RF    : 上端粘着, 下端応力なし条件
    FR    : 上端応力なし, 下端粘着条件
    FF    : 両端応力なし条件
new :logical, intent(IN), optional
: (in) true だと境界条件計算用行列を強制的に新たに作る.
    default は false.

速度ポロイダルポテンシャルΦを▽^2Φから計算する. チェビシェフ格子点空間で境界条件を適用している.

この関数を用いるためには wt_Initial にて設定する チェビシェフ切断波数(lm)と鉛直格子点数(km)を等しく しておく必要がある.

速度ポロイダルポテンシャルΦを f = ▽^2Φから定める式は

   ▽^2Φ = f
     Φ = const. at boundaries.
     ∂Φ/∂r = 0 at boundaries          (粘着条件)
     or ∂^2Φ/∂r^2 = 0 at boundaries   (応力なし条件)

最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.

最終的にチェビシェフ係数の解が欲しい場合には, wz_LaplaPol2Pol_wz に 比べてチェビシェフ — 格子点変換が 1 回分少なくて済む.

Original external subprogram is wt_module#wt_LaplaPol2PolGrid_wt

wt_LaplaPol2PolTau_wt( wt, [cond], [new] ) result(wt_LaplaPol2PolTau_wt)
Function :
wt_LaplaPol2PolTau_wt :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) 出力ポロイダルポテンシャル分布
wt :real(8), dimension((nm+1)*(nm+1),0:lm),intent(in)
: (in) 入力▽^2φ分布
cond :character(len=2), intent(in), optional
: (in) 境界条件スイッチ. 省略時は ‘RR‘
    RR    : 両端粘着条件
    RF    : 上端粘着, 下端応力なし条件
    FR    : 上端応力なし, 下端粘着条件
    FF    : 両端応力なし条件
new :logical, intent(IN), optional
: (in) true だと境界条件計算用行列を強制的に新たに作る.
    default は false.

速度ポロイダルポテンシャルΦを▽^2Φから計算する. チェビシェフタウ法で境界条件を適用している.

速度ポロイダルポテンシャルΦを f = ▽^2Φから定める式は

   ▽^2Φ = f
     Φ = const. at boundaries.
     ∂Φ/∂r = 0 at boundaries          (粘着条件)
     or ∂^2Φ/∂r^2 = 0 at boundaries   (応力なし条件)

最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.

Original external subprogram is wt_module#wt_LaplaPol2PolTau_wt

wt_LaplaPol2Pol_wt( wt, [cond], [new] ) result(wt_LaplaPol2PolTau_wt)
Function :
wt_LaplaPol2PolTau_wt :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) 出力ポロイダルポテンシャル分布
wt :real(8), dimension((nm+1)*(nm+1),0:lm),intent(in)
: (in) 入力▽^2φ分布
cond :character(len=2), intent(in), optional
: (in) 境界条件スイッチ. 省略時は ‘RR‘
    RR    : 両端粘着条件
    RF    : 上端粘着, 下端応力なし条件
    FR    : 上端応力なし, 下端粘着条件
    FF    : 両端応力なし条件
new :logical, intent(IN), optional
: (in) true だと境界条件計算用行列を強制的に新たに作る.
    default は false.

速度ポロイダルポテンシャルΦを▽^2Φから計算する. チェビシェフタウ法で境界条件を適用している.

速度ポロイダルポテンシャルΦを f = ▽^2Φから定める式は

   ▽^2Φ = f
     Φ = const. at boundaries.
     ∂Φ/∂r = 0 at boundaries          (粘着条件)
     or ∂^2Φ/∂r^2 = 0 at boundaries   (応力なし条件)

最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.

Original external subprogram is wt_module#wt_LaplaPol2Pol_wt

wt_Lapla_wt( wt ) result(wt_Lapla_wt)
Function :
wt_Lapla_wt :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) ラプラシアンを作用された 2 次元スペクトルデータ
wt :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) 2 次元球面調和函数チェビシェフスペクトルデータ

入力スペクトルデータにラプラシアン

    ▽^2 =   1/r^2 cos^2φ・∂^2/∂λ^2
           + 1/r^2 cosφ・∂/∂φ(cosφ∂/∂φ)
           + 1/r^2 ∂/∂r (r^2 ∂/∂r)

を作用する.

スペクトルデータのラプラシアンとは, 対応する格子点データに ラプラシアンを作用させたデータのスペクトル変換のことである.

Original external subprogram is wt_module#wt_Lapla_wt

wt_PolMagBoundaries( wt_POL, [new] )
Subroutine :
wt_POL :real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
: (inout) 境界条件を適用するデータ. 修正された値を返す.
new :logical, intent(IN), optional
: (in) true だと境界条件計算用行列を強制的に新たに作る.
    default は false.

磁場ポロイダルポテンシャルに対して境界条件を適用する. Chebyshev 空間での境界条件適用

チェビシェフ空間において境界条件を満たすべく高次の係数を定める方法を とっている(タウ法). 現在のところ境界物質が非電気伝導体の場合のみ 対応している. その場合, 磁場ポロイダルポテンシャルの各水平スペクトル 成分 h にたいして境界条件が与えられ,

 * 外側境界 : dh/dr + (n+1)h/r = 0
 * 内側境界 : dh/dr - nh/r = 0

である. ここで n は h の水平全波数である.

最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.

Original external subprogram is wt_module#wt_PolMagBoundaries

wt_PolmagBoundariesGrid( wt_POL, [new] )
Subroutine :
wt_POL :real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
: (inout) 境界条件を適用するデータ. 修正された値を返す.
new :logical, intent(IN), optional
: (in) true だと境界条件計算用行列を強制的に新たに作る.
    default は false.

磁場ポロイダルポテンシャルに対して境界条件を適用する. 鉛直実空間での境界条件適用.

鉛直実格子点空間において内部領域の値と境界条件を満たすように 条件を課している(選点法). このルーチンを用いるためには wt_Initial にて設定するチェビシェフ切断波数(lm)と鉛直格子点数(km)を 等しくしておく必要がある.

現在のところ境界物質が非電気伝導体の場合のみ対応している. その場合, 磁場ポロイダルポテンシャルの各水平スペクトル成分 h に たいして境界条件が与えられ,

 * 外側境界 : dh/dr + (n+1)h/r = 0
 * 内側境界 : dh/dr - nh/r = 0

である. ここで n は h の水平全波数である.

最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.

Original external subprogram is wt_module#wt_PolmagBoundariesGrid

wt_PolmagBoundariesTau( wt_POL, [new] )
Subroutine :
wt_POL :real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
: (inout) 境界条件を適用するデータ. 修正された値を返す.
new :logical, intent(IN), optional
: (in) true だと境界条件計算用行列を強制的に新たに作る.
    default は false.

磁場ポロイダルポテンシャルに対して境界条件を適用する. Chebyshev 空間での境界条件適用

チェビシェフ空間において境界条件を満たすべく高次の係数を定める方法を とっている(タウ法). 現在のところ境界物質が非電気伝導体の場合のみ 対応している. その場合, 磁場ポロイダルポテンシャルの各水平スペクトル 成分 h にたいして境界条件が与えられ,

 * 外側境界 : dh/dr + (n+1)h/r = 0
 * 内側境界 : dh/dr - nh/r = 0

である. ここで n は h の水平全波数である.

最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.

Original external subprogram is wt_module#wt_PolmagBoundariesTau

wt_Potential2Rotation( xyz_RotVLON, xyz_RotVLAT, xyz_RotVRAD, wt_TORPOT, wt_POLPOT )
Subroutine :
xyz_RotVLON :real(8), dimension(0:im-1,1:jm,0:km), intent(OUT)
: (out) 回転の経度成分
xyz_RotVLAT :real(8), dimension(0:im-1,1:jm,0:km), intent(OUT)
: (out) 回転の緯度成分
xyz_RotVRAD :real(8), dimension(0:im-1,1:jm,0:km), intent(OUT)
: (out) 回転の動径成分
wt_TORPOT :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) トロイダルポテンシャル
wt_POLPOT :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) ポロイダルポテンシャル

トロイダルポロイダルポテンシャルΨ,Φで表される非発散ベクトル場

    v = ▽x(Ψr) + ▽x▽x(Φr)

に対して, その回転

    ▽xv = ▽x▽x(Ψr) + ▽x▽x▽x(Φr) = ▽x▽x(Ψr) - ▽x((▽^2Φ)r)

を計算する.

Original external subprogram is wt_module#wt_Potential2Rotation

Subroutine :
xvz_RotVLON :real(8), dimension(0:im-1,jc,0:km), intent(OUT)
: (out) 回転の経度成分
xvz_RotVLAT :real(8), dimension(0:im-1,jc,0:km), intent(OUT)
: (out) 回転の緯度成分
xvz_RotVRAD :real(8), dimension(0:im-1,jc,0:km), intent(OUT)
: (out) 回転の動径成分
wt_TORPOT :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) トロイダルポテンシャル
wt_POLPOT :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) ポロイダルポテンシャル

トロイダルポロイダルポテンシャルΨ,Φで表される非発散ベクトル場

    v = ▽x(Ψr) + ▽x▽x(Φr)

に対して, その回転

    ▽xv = ▽x▽x(Ψr) + ▽x▽x▽x(Φr) = ▽x▽x(Ψr) - ▽x((▽^2Φ)r)

を計算する.

[Source]

    subroutine wt_Potential2RotationMPI( xvz_RotVLON,xvz_RotVLAT,xvz_RotVRAD,wt_TORPOT,wt_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,jc,0:km), intent(OUT) :: xvz_RotVLON
      !(out) 回転の経度成分

      real(8), dimension(0:im-1,jc,0:km), intent(OUT) :: xvz_RotVLAT
      !(out) 回転の緯度成分

      real(8), dimension(0:im-1,jc,0:km), intent(OUT) :: xvz_RotVRAD
      !(out) 回転の動径成分

      ! 入力ベクトル場を表すポテンシャル
      real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_TORPOT
      !(in) トロイダルポテンシャル

      real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_POLPOT
      !(in) ポロイダルポテンシャル

      call wt_Potential2VectorMPI( xvz_RotVLON,xvz_RotVLAT,xvz_RotVRAD, -wt_Lapla_wt(wt_POLPOT), wt_TORPOT)

    end subroutine wt_Potential2RotationMPI
wt_Potential2Vector( xyz_VLON, xyz_VLAT, xyz_VRAD, wt_TORPOT, wt_POLPOT )
Subroutine :
xyz_VLON :real(8), dimension(0:im-1,1:jm,0:km)
: (out) ベクトル場の経度成分
xyz_VLAT :real(8), dimension(0:im-1,1:jm,0:km)
: (out) ベクトル場の緯度成分
xyz_VRAD :real(8), dimension(0:im-1,1:jm,0:km)
: (out) ベクトル場の動径成分
wt_TORPOT :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) トロイダルポテンシャル
wt_POLPOT :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) ポロイダルポテンシャル

トロイダルポロイダルポテンシャルΨ,Φで表される非発散ベクトル場

    v = ▽x(Ψr) + ▽x▽x(Φr)

の各成分を計算する

Original external subprogram is wt_module#wt_Potential2Vector

Subroutine :
xvz_VLON :real(8), dimension(0:im-1,jc,0:km)
: (out) ベクトル場の経度成分
xvz_VLAT :real(8), dimension(0:im-1,jc,0:km)
: (out) ベクトル場の緯度成分
xvz_VRAD :real(8), dimension(0:im-1,jc,0:km)
: (out) ベクトル場の動径成分
wt_TORPOT :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) トロイダルポテンシャル
wt_POLPOT :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) ポロイダルポテンシャル

トロイダルポロイダルポテンシャルΨ,Φで表される非発散ベクトル場

    v = ▽x(Ψr) + ▽x▽x(Φr)

の各成分を計算する

[Source]

    subroutine wt_Potential2VectorMPI( xvz_VLON,xvz_VLAT,xvz_VRAD,wt_TORPOT,wt_POLPOT)
      !
      ! トロイダルポロイダルポテンシャルΨ,Φで表される非発散ベクトル場
      !
      !     v = ▽x(Ψr) + ▽x▽x(Φr) 
      !
      ! の各成分を計算する
      !
      real(8), dimension(0:im-1,jc,0:km)     :: xvz_VLON
      !(out) ベクトル場の経度成分

      real(8), dimension(0:im-1,jc,0:km)     :: xvz_VLAT
      !(out) ベクトル場の緯度成分

      real(8), dimension(0:im-1,jc,0:km)     :: xvz_VRAD
      !(out) ベクトル場の動径成分

      real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_TORPOT
      !(in) トロイダルポテンシャル

      real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_POLPOT
      !(in) ポロイダルポテンシャル

      xvz_VLON =   xvz_RAD * xvz_GradLat_wt(wt_TORPOT) + xva_GradLon_wa(wz_wt(wt_RotRad_wt(wt_POLPOT)))
      xvz_VLAT = - xvz_RAD * xvz_GradLon_wt(wt_TORPOT) + xva_GradLat_wa(wz_wt(wt_RotRad_wt(wt_POLPOT)))
      xvz_VRAD = xvz_wt(wt_L2_wt(wt_POLPOT))/xvz_RAD

    end subroutine wt_Potential2VectorMPI
Function :
wt_QOperatorMPI_wt :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) Q 演算子を作用された 2 次元スペクトルデータ
wt :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) 2 次元球面調和函数チェビシェフスペクトルデータ

入力スペクトルデータに対応する格子点データに演算子

   Q=(k・▽-1/2(L2 k・▽+ k・▽L2))

を作用させたデータのスペクトル変換が返される.

[Source]

    function wt_QOperatorMPI_wt(wt)
      !
      ! 入力スペクトルデータに対応する格子点データに演算子
      !
      !    Q=(k・▽-1/2(L2 k・▽+ k・▽L2))
      !
      ! を作用させたデータのスペクトル変換が返される.
      !
      real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt
      !(in) 2 次元球面調和函数チェビシェフスペクトルデータ

      real(8), dimension((nm+1)*(nm+1),0:lm)             :: wt_QOperatorMPI_wt
      !(out) Q 演算子を作用された 2 次元スペクトルデータ

      wt_QOperatorMPI_wt = wt_xvz(xvz_KGrad_wt(wt) - xvz_KGrad_wt(wt_L2_wt(wt))/2) - wt_L2_wt(wt_xvz(xvz_KGrad_wt(wt)))/2

    end function wt_QOperatorMPI_wt
wt_QOperator_wt( wt ) result(wt_QOperator_wt)
Function :
wt_QOperator_wt :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) Q 演算子を作用された 2 次元スペクトルデータ
wt :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) 2 次元球面調和函数チェビシェフスペクトルデータ

入力スペクトルデータに対応する格子点データに演算子

   Q=(k・▽-1/2(L2 k・▽+ k・▽L2))

を作用させたデータのスペクトル変換が返される.

Original external subprogram is wt_module#wt_QOperator_wt

Function :
wt_RadRotRot_xvz_xvz_xvz :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) ベクトル v の r・(▽×▽×v)
xvz_VLON :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) ベクトルの経度成分
xvz_VLAT :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) ベクトルの緯度成分
xvz_VRAD :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (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

のスペクトルデータが返される.

[Source]

    function wt_RadRotRot_xvz_xvz_xvz(xvz_VLON,xvz_VLAT,xvz_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,jc,0:km), intent(in) :: xvz_VLON
      !(in) ベクトルの経度成分

      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz_VLAT
      !(in) ベクトルの緯度成分

      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz_VRAD
      !(in) ベクトルの動径成分

      real(8), dimension((nm+1)*(nm+1),0:lm)     :: wt_RadRotRot_xvz_xvz_xvz
      !(out) ベクトル v の r・(▽×▽×v) 

      wt_RadRotRot_xvz_xvz_xvz = wt_RotRad_wt(wt_wz( (wa_DivLon_xva(xvz_VLON)+ wa_DivLat_xva(xvz_VLAT)))) + wt_L2_wt(wt_xvz(xvz_VRAD/xvz_RAD))

    end function wt_RadRotRot_xvz_xvz_xvz
wt_RadRotRot_xyz_xyz_xyz( xyz_VLON, xyz_VLAT, xyz_VRAD ) result(wt_RadRotRot_xyz_xyz_xyz)
Function :
wt_RadRotRot_xyz_xyz_xyz :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) ベクトル v の r・(▽×▽×v)
xyz_VLON :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) ベクトルの経度成分
xyz_VLAT :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) ベクトルの緯度成分
xyz_VRAD :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (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

のスペクトルデータが返される.

Original external subprogram is wt_module#wt_RadRotRot_xyz_xyz_xyz

Function :
wt_RadRot_xvz_xvz :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) ベクトルの渦度と動径ベクトルの内積
xvz_VLON :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) ベクトルの経度成分
xvz_VLAT :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) ベクトルの緯度成分

r・(▽×v)

ベクトルの渦度と動径ベクトルの内積 r・(▽×v) を計算する.

第 1, 2 引数(v[λ], v[φ])がそれぞれベクトルの経度成分, 緯度成分を表す.

   r・(▽×v) = 1/cosφ・∂v[φ]/∂λ - 1/cosφ・∂(v[λ] cosφ)/∂φ

のスペクトル データが返される.

[Source]

    function wt_RadRot_xvz_xvz(xvz_VLON,xvz_VLAT)  ! r・(▽×v)
      !
      ! ベクトルの渦度と動径ベクトルの内積 r・(▽×v) を計算する.
      !
      ! 第 1, 2 引数(v[λ], v[φ])がそれぞれベクトルの経度成分, 緯度成分を表す.
      !
      !    r・(▽×v) = 1/cosφ・∂v[φ]/∂λ - 1/cosφ・∂(v[λ] cosφ)/∂φ
      !
      ! のスペクトル データが返される.
      !
      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz_VLON
      !(in) ベクトルの経度成分

      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz_VLAT
      !(in) ベクトルの緯度成分

      real(8), dimension((nm+1)*(nm+1),0:lm)     :: wt_RadRot_xvz_xvz
      !(out) ベクトルの渦度と動径ベクトルの内積

      wt_RadRot_xvz_xvz = wt_wz(wa_DivLon_xva(xvz_VLAT) - wa_DivLat_xva(xvz_VLON))
      
    end function wt_RadRot_xvz_xvz
wt_RadRot_xyz_xyz( xyz_VLON, xyz_VLAT ) result(wt_RadRot_xyz_xyz)
Function :
wt_RadRot_xyz_xyz :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) ベクトルの渦度と動径ベクトルの内積
xyz_VLON :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) ベクトルの経度成分
xyz_VLAT :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) ベクトルの緯度成分

r・(▽×v)

ベクトルの渦度と動径ベクトルの内積 r・(▽×v) を計算する.

第 1, 2 引数(v[λ], v[φ])がそれぞれベクトルの経度成分, 緯度成分を表す.

   r・(▽×v) = 1/cosφ・∂v[φ]/∂λ - 1/cosφ・∂(v[λ] cosφ)/∂φ

のスペクトル データが返される.

Original external subprogram is wt_module#wt_RadRot_xyz_xyz

wt_RotRad_wt( wt ) result(wt_RotRad_wt)
Function :
wt_RotRad_wt :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) 回転型動径微分を作用された 2 次元スペクトルデータ
wt :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) 2 次元球面調和函数チェビシェフスペクトルデータ

入力スペクトルデータに回転型動径微分

     1/r ∂(r.)/∂r = ∂(.)/∂r + (.)/r

を作用する.

スペクトルデータの回転型動径微分とは, 対応する格子点データに 回転型動径微分を作用させたデータのスペクトル変換のことである.

Original external subprogram is wt_module#wt_RotRad_wt

Function :
wt_RotRad_xvz_xvz :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) ベクトル場の回転の動径成分
xvz_Vlat :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) ベクトル場の緯度成分
xvz_Vlon :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) ベクトル場の経度成分

ベクトルの緯度成分, 経度成分である第 1, 2 引数 Vlat, Vlon に対して ベクトル場の回転の動径成分

   1/rcosφ・∂Vlat/∂λ - 1/rcosφ・∂(Vlon cosφ)/∂φ

を計算する.

[Source]

    function wt_RotRad_xvz_xvz(xvz_Vlat,xvz_Vlon) 
      !
      ! ベクトルの緯度成分, 経度成分である第 1, 2 引数 Vlat, Vlon に対して
      ! ベクトル場の回転の動径成分 
      !
      !    1/rcosφ・∂Vlat/∂λ - 1/rcosφ・∂(Vlon cosφ)/∂φ
      !
      ! を計算する.
      !
      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz_Vlat
      !(in) ベクトル場の緯度成分

      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz_Vlon
      !(in) ベクトル場の経度成分

      real(8), dimension((nm+1)*(nm+1),0:lm)     :: wt_RotRad_xvz_xvz
      !(out) ベクトル場の回転の動径成分

        wt_RotRad_xvz_xvz =   wt_DivLon_xvz(xvz_Vlat) - wt_DivLat_xvz(xvz_Vlon)

    end function wt_RotRad_xvz_xvz
wt_RotRad_xyz_xyz( xyz_Vlat, xyz_Vlon ) result(wt_RotRad_xyz_xyz)
Function :
wt_RotRad_xyz_xyz :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) ベクトル場の回転の動径成分
xyz_Vlat :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) ベクトル場の緯度成分
xyz_Vlon :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) ベクトル場の経度成分

ベクトルの緯度成分, 経度成分である第 1, 2 引数 Vlat, Vlon に対して ベクトル場の回転の動径成分

   1/rcosφ・∂Vlat/∂λ - 1/rcosφ・∂(Vlon cosφ)/∂φ

を計算する.

Original external subprogram is wt_module#wt_RotRad_xyz_xyz

wt_TorBoundaries( wt_TORPOT, [values], [cond], [new] )
Subroutine :
wt_TORPOT :real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
: (inout) 境界条件を適用するデータ. 修正された値を返す.
values :real(8), dimension((nm+1)*(nm+1),2), intent(in), optional
: (in) 両端境界でのトロイダルポテンシャル
    粘着条件の時のみ有効
cond :character(len=2), intent(in), optional
: (in) 境界条件スイッチ. 省略時は ‘RR‘
    RR    : 両端粘着条件
    RF    : 上端粘着, 下端応力なし条件
    FR    : 上端応力なし, 下端粘着条件
    FF    : 両端応力なし条件
new :logical, intent(IN), optional
: (in) true だと境界条件計算用行列を強制的に新たに作る.
    default は false.

速度トロイダルポテンシャルに対して境界条件を適用する. Chebyshev 空間での境界条件適用.

速度トロイダルポテンシャルΨに対して与えられる境界条件は

  * 粘着条件 : Ψ = Ψb(lon,lat). Ψb は境界球面での速度分布.
                                  default は 0(静止状態).

  * 応力なし条件 : ∂(Ψ/r)/∂r = 0.

最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.

Original external subprogram is wt_module#wt_TorBoundaries

wt_TorBoundariesGrid( wt_TORPOT, [values], [cond], [new] )
Subroutine :
wt_TORPOT :real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
: (inout) 境界条件を適用するデータ. 修正された値を返す.
values :real(8), dimension((nm+1)*(nm+1),2), intent(in), optional
: (in) 両端境界でのトロイダルポテンシャル
    粘着条件の時のみ有効
cond :character(len=2), intent(in), optional
: (in) 境界条件スイッチ. 省略時は ‘RR‘
    RR    : 両端粘着条件
    RF    : 上端粘着, 下端応力なし条件
    FR    : 上端応力なし, 下端粘着条件
    FF    : 両端応力なし条件
new :logical, intent(IN), optional
: (in) true だと境界条件計算用行列を強制的に新たに作る.
    default は false.

速度トロイダルポテンシャルに対して境界条件を適用する. 実空間での境界条件適用

鉛直実格子点空間において内部領域の値と境界条件を満たすように 条件を課している(選点法). このルーチンを用いるためには wt_Initial にて設定するチェビシェフ切断波数(lm)と鉛直格子点数(km)を 等しくしておく必要がある.

速度トロイダルポテンシャルΨに対して与えられる境界条件は

  * 粘着条件 : Ψ = Ψb(lon,lat). Ψb は境界球面での速度分布.
                                  default は 0 (静止状態).

  * 応力なし条件 : ∂(Ψ/r)/∂r = 0.

最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.

Original external subprogram is wt_module#wt_TorBoundariesGrid

wt_TorBoundariesTau( wt_TORPOT, [values], [cond], [new] )
Subroutine :
wt_TORPOT :real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
: (inout) 境界条件を適用するデータ. 修正された値を返す.
values :real(8), dimension((nm+1)*(nm+1),2), intent(in), optional
: (in) 両端境界でのトロイダルポテンシャル
    粘着条件の時のみ有効
cond :character(len=2), intent(in), optional
: (in) 境界条件スイッチ. 省略時は ‘RR‘
    RR    : 両端粘着条件
    RF    : 上端粘着, 下端応力なし条件
    FR    : 上端応力なし, 下端粘着条件
    FF    : 両端応力なし条件
new :logical, intent(IN), optional
: (in) true だと境界条件計算用行列を強制的に新たに作る.
    default は false.

速度トロイダルポテンシャルに対して境界条件を適用する. Chebyshev 空間での境界条件適用.

速度トロイダルポテンシャルΨに対して与えられる境界条件は

  * 粘着条件 : Ψ = Ψb(lon,lat). Ψb は境界球面での速度分布.
                                  default は 0(静止状態).

  * 応力なし条件 : ∂(Ψ/r)/∂r = 0.

最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.

Original external subprogram is wt_module#wt_TorBoundariesTau

wt_TorMagBoundaries( wt_TOR, [new] )
Subroutine :
wt_TOR :real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
: (inout) 境界条件を適用するデータ. 修正された値を返す.
new :logical, intent(IN), optional
: (in) true だと境界条件計算用行列を強制的に新たに作る.
    default は false.

Original external subprogram is wt_module#wt_TorMagBoundaries

wt_TormagBoundariesGrid( wt_TOR, [new] )
Subroutine :
wt_TOR :real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
: (inout) 境界条件を適用するデータ. 修正された値を返す.
new :logical, intent(IN), optional
: (in) true だと境界条件計算用行列を強制的に新たに作る.
    default は false.

磁場トロイダルポテンシャルに対して境界条件を適用する. 鉛直実空間での境界条件適用.

鉛直実格子点空間において内部領域の値と境界条件を満たすように 条件を課している(選点法). このルーチンを用いるためには wt_Initial にて設定するチェビシェフ切断波数(lm)と鉛直格子点数(km)を 等しくしておく必要がある.

現在のところ境界物質が非電気伝導体の場合のみ対応している. その場合, 磁場トロイダルポテンシャルの境界条件は

外側

   wt_psi = 0   at the outer boundary

内側

   wt_psi = 0       at the inner boundary

であるので wt_Boundaries で対応可能だが, 将来のため別途作成しておく

最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.

Original external subprogram is wt_module#wt_TormagBoundariesGrid

wt_TormagBoundariesTau( wt_TOR, [new] )
Subroutine :
wt_TOR :real(8), dimension((nm+1)*(nm+1),0:lm),intent(inout)
: (inout) 境界条件を適用するデータ. 修正された値を返す.
new :logical, intent(IN), optional
: (in) true だと境界条件計算用行列を強制的に新たに作る.
    default は false.

Original external subprogram is wt_module#wt_TormagBoundariesTau

wt_VGradV( xyz_VGRADV_LON, xyz_VGRADV_LAT, xyz_VGRADV_RAD, xyz_VLON, xyz_VLAT, xyz_VRAD )
Subroutine :
xyz_VGRADV_LON :real(8), dimension(0:im-1,1:jm,0:km),intent(out)
: (out) (v・▽v) 経度成分
xyz_VGRADV_LAT :real(8), dimension(0:im-1,1:jm,0:km),intent(out)
: (out) (v・▽v) 緯度成分
xyz_VGRADV_RAD :real(8), dimension(0:im-1,1:jm,0:km),intent(out)
: (out) (v・▽v) 動径成分
xyz_VLON :real(8), dimension(0:im-1,1:jm,0:km),intent(in)
: (in) ベクトル場 v の経度成分
xyz_VLAT :real(8), dimension(0:im-1,1:jm,0:km),intent(in)
: (in) ベクトル場 v の緯度成分
xyz_VRAD :real(8), dimension(0:im-1,1:jm,0:km),intent(in)
: (in) ベクトル場 v の動径成分

ベクトル場の v・▽v を計算する.

ベクトル場 v=(v[λ],v[φ],v[r]) に対するv・▽vの各成分は 次のように計算される.

  (v・▽v)[λ] = ▽・(v[λ]v) + v[λ]v[r]/r - v[λ]v[φ]tan(φ)/r
  (v・▽v)[φ] = ▽・(v[φ]v) + v[φ]v[r]/r - v[λ]^2tan(φ)/r
  (v・▽v)[r] = ▽・(v[r]v) + (v[λ]^2+v[φ]^2)/r

非発散速度場に対してはポテンシャルから wt_Potential2Rotation を 用いて回転を計算し, 恒等式 v・▽v = ▽(v[2^/2) - vx▽xv を 用いる方がよいだろう.

Original external subprogram is wt_module#wt_VGradV

Subroutine :
xvz_VGRADV_LON :real(8), dimension(0:im-1,jc,0:km),intent(out)
: (out) (v・▽v) 経度成分
xvz_VGRADV_LAT :real(8), dimension(0:im-1,jc,0:km),intent(out)
: (out) (v・▽v) 緯度成分
xvz_VGRADV_RAD :real(8), dimension(0:im-1,jc,0:km),intent(out)
: (out) (v・▽v) 動径成分
xvz_VLON :real(8), dimension(0:im-1,jc,0:km),intent(in)
: (in) ベクトル場 v の経度成分
xvz_VLAT :real(8), dimension(0:im-1,jc,0:km),intent(in)
: (in) ベクトル場 v の緯度成分
xvz_VRAD :real(8), dimension(0:im-1,jc,0:km),intent(in)
: (in) ベクトル場 v の動径成分

ベクトル場の v・▽v を計算する.

ベクトル場 v=(v[λ],v[φ],v[r]) に対するv・▽vの各成分は 次のように計算される.

  (v・▽v)[λ] = ▽・(v[λ]v) + v[λ]v[r]/r - v[λ]v[φ]tan(φ)/r
  (v・▽v)[φ] = ▽・(v[φ]v) + v[φ]v[r]/r - v[λ]^2tan(φ)/r
  (v・▽v)[r] = ▽・(v[r]v) + (v[λ]^2+v[φ]^2)/r

非発散速度場に対してはポテンシャルから wt_Potential2Rotation を 用いて回転を計算し, 恒等式 v・▽v = ▽(v[2^/2) - vx▽xv を 用いる方がよいだろう.

[Source]

    subroutine wt_VGradVMPI(xvz_VGRADV_LON,xvz_VGRADV_LAT,xvz_VGRADV_RAD, xvz_VLON,xvz_VLAT,xvz_VRAD )
      !
      ! ベクトル場の v・▽v を計算する.
      !
      ! ベクトル場 v=(v[λ],v[φ],v[r]) に対するv・▽vの各成分は
      ! 次のように計算される.
      !
      !   (v・▽v)[λ] = ▽・(v[λ]v) + v[λ]v[r]/r - v[λ]v[φ]tan(φ)/r
      !   (v・▽v)[φ] = ▽・(v[φ]v) + v[φ]v[r]/r - v[λ]^2tan(φ)/r
      !   (v・▽v)[r] = ▽・(v[r]v) + (v[λ]^2+v[φ]^2)/r
      !
      ! 非発散速度場に対してはポテンシャルから wt_Potential2Rotation を
      ! 用いて回転を計算し, 恒等式 v・▽v = ▽(v[2^/2) - vx▽xv を
      ! 用いる方がよいだろう.
      !
      real(8), dimension(0:im-1,jc,0:km),intent(out)   :: xvz_VGRADV_LON
      !(out) (v・▽v) 経度成分

      real(8), dimension(0:im-1,jc,0:km),intent(out)   :: xvz_VGRADV_LAT
      !(out) (v・▽v) 緯度成分

      real(8), dimension(0:im-1,jc,0:km),intent(out)   :: xvz_VGRADV_RAD
      !(out) (v・▽v) 動径成分

      real(8), dimension(0:im-1,jc,0:km),intent(in)    :: xvz_VLON
      !(in) ベクトル場 v の経度成分

      real(8), dimension(0:im-1,jc,0:km),intent(in)    :: xvz_VLAT
      !(in) ベクトル場 v の緯度成分

      real(8), dimension(0:im-1,jc,0:km),intent(in)    :: xvz_VRAD
      !(in) ベクトル場 v の動径成分

      xvz_VGRADV_LON = xvz_Div_xvz_xvz_xvz( xvz_VLON * xvz_VLON, xvz_VLON*xvz_VLAT, xvz_VLON*xvz_VRAD ) + xvz_VLON*xvz_VRAD/xvz_RAD - xvz_VLON*xvz_VLAT*tan(xvz_LAT)/xvz_RAD 

      xvz_VGRADV_LAT = xvz_Div_xvz_xvz_xvz( xvz_VLAT*xvz_VLON, xvz_VLAT*xvz_VLAT, xvz_VLAT*xvz_VRAD ) + xvz_VLAT*xvz_VRAD/xvz_RAD + xvz_VLON**2*tan(xvz_LAT)/xvz_RAD 

      xvz_VGRADV_RAD = xvz_Div_xvz_xvz_xvz( xvz_VRAD*xvz_VLON, xvz_VRAD*xvz_VLAT, xvz_VRAD*xvz_VRAD ) - (xvz_VLON**2 + xvz_VLAT**2)/xvz_RAD 

    end subroutine wt_VGradVMPI
wt_VMiss
Variable :
wt_VMiss = -999.0 :real(8)
: 欠損値

Original external subprogram is wt_module#wt_VMiss

Subroutine :
i :integer,intent(in)
: 格子点数(経度λ)
j :integer,intent(in)
: 格子点数(緯度φ)
k :integer,intent(in)
: 格子点数(動径 r)
n :integer,intent(in)
: 切断波数(水平全波数)
l :integer,intent(in)
: 切断波数(動径波数)
r_in :real(8),intent(in)
: 球殻内半径
r_out :real(8),intent(in)
: 球殻外半径
np :integer,intent(in), optional
: OPENMP での最大スレッド数
wa_init :logical,intent(in), optional
: wa_initial スイッチ

スペクトル変換の格子点数, 波数, 動径座標の範囲を設定する.

他の関数を呼ぶ前に, 最初にこのサブルーチンを呼んで初期設定を しなければならない.

[Source]

   subroutine wt_mpi_Initial(i,j,k,n,l,r_in,r_out,np,wa_init)
     !
     ! スペクトル変換の格子点数, 波数, 動径座標の範囲を設定する.
     !
     ! 他の関数を呼ぶ前に, 最初にこのサブルーチンを呼んで初期設定を
     ! しなければならない. 
     !
     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_in           ! 球殻内半径
     real(8),intent(in) :: r_out          ! 球殻外半径

     integer,intent(in), optional :: np   ! OPENMP での最大スレッド数
     logical,intent(in), optional :: wa_init   ! wa_initial スイッチ

     logical    :: wa_initialize=.true.   ! wa_initial スイッチ

     im = i  
      jm = j 
      km = k
     nm = n  
      lm = l
     ri = r_in 
      ro = r_out

     if ( present(wa_init) ) then
        wa_initialize = wa_init
     else
        wa_initialize = .true.
     endif

     if ( present(np) ) then
        call wt_Initial(im,jm,km,nm,lm,ri,ro,np=np,wa_init=wa_initialize)
     else
        call wt_Initial(im,jm,km,nm,lm,ri,ro,wa_init=wa_initialize)
     endif
     if ( wa_initialize ) then
        call w_base_mpi_Initial
        call w_deriv_mpi_Initial
        call wa_base_mpi_Initial
     endif

     allocate(xvz_Lon(0:im-1,jc,0:km))
     allocate(xvz_Lat(0:im-1,jc,0:km))
     allocate(xvz_Rad(0:im-1,jc,0:km))

     xvz_Lon = spread(xv_Lon,3,km+1)
     xvz_Lat = spread(xv_Lat,3,km+1)
     xvz_Rad = spread(spread(z_Rad,1,jc),1,im)

     call MessageNotify('M','wt_mpi_initial', 'wt_mpi_module (2013/08/20) is initialized')

   end subroutine wt_mpi_Initial
wt_wz( wz ) result(wt_wz)
Function :
wt_wz :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) 2 次元球面調和函数チェビシェフスペクトルデータ
wz :real(8), dimension((nm+1)*(nm+1),0:km), intent(in)
: (in) 2 次元球面調和函数スペクトル・動径格子点データ

水平スペクトル・動径格子点データからスペクトルデータへ(正)変換する.

Original external subprogram is wt_module#wt_wz

Function :
wt_xvz :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) 2 次元球面調和函数チェビシェフスペクトルデータ
xvz :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

3 次元格子点データからスペクトルデータへ(正)変換する.

[Source]

    function wt_xvz(xvz)
      !
      ! 3 次元格子点データからスペクトルデータへ(正)変換する.
      !
      real(8), dimension((nm+1)*(nm+1),0:lm)             :: wt_xvz
      !(out) 2 次元球面調和函数チェビシェフスペクトルデータ

      real(8), dimension(0:im-1,jc,0:km), intent(in)         :: xvz
      !(in) 3 次元経度緯度動径格子点データ

      wt_xvz = wt_wz(wa_xva(xvz))

    end function wt_xvz
wt_xyz( xyz ) result(wt_xyz)
Function :
wt_xyz :real(8), dimension((nm+1)*(nm+1),0:lm)
: (out) 2 次元球面調和函数チェビシェフスペクトルデータ
xyz :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

3 次元格子点データからスペクトルデータへ(正)変換する.

Original external subprogram is wt_module#wt_xyz

wz_LaplaPol2Pol_wz( wz, [cond], [new] ) result(wz_LaplaPol2Pol_wz)
Function :
wz_LaplaPol2Pol_wz :real(8), dimension((nm+1)*(nm+1),0:km)
: (out) 出力ポロイダルポテンシャル分布
wz :real(8), dimension((nm+1)*(nm+1),0:km),intent(in)
: (in) 入力▽^2φ分布
cond :character(len=2), intent(in), optional
: (in) 境界条件スイッチ. 省略時は ‘RR‘
    RR    : 両端粘着条件
    RF    : 上端粘着, 下端応力なし条件
    FR    : 上端応力なし, 下端粘着条件
    FF    : 両端応力なし条件
new :logical, intent(IN), optional
: (in) true だと境界条件計算用行列を強制的に新たに作る.
    default は false.

速度ポロイダルポテンシャルΦを▽^2Φから計算する.

チェビシェフ格子点空間で境界条件を適用している. この関数を用いるためには wt_Initial にて設定する チェビシェフ切断波数(lm)と鉛直格子点数(km)を等しく しておく必要がある.

速度ポロイダルポテンシャルΦを f = ▽^2Φから定める式は

  ▽^2Φ = f
    Φ = const. at boundaries.
    ∂Φ/∂r = 0 at boundaries           (粘着条件)
    or ∂^2Φ/∂r^2 = 0 at boundaries    (応力なし条件)

最初に呼ばれるときはオプショナル引数 new に関係なく行列が設定される.

Original external subprogram is wt_module#wz_LaplaPol2Pol_wz

wz_RAD
Variable :
wz_RAD :real(8), dimension(:,:), allocatable
: 座標

Original external subprogram is wt_module#wz_RAD

wz_wt( wt ) result(wz_wt)
Function :
wz_wt :real(8), dimension((nm+1)*(nm+1),0:km)
: (out) 2 次元球面調和函数スペクトル・動径格子点データ
wt :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) 2 次元球面調和函数チェビシェフスペクトルデータ

スペクトルデータから水平スペクトル・動径格子点データへ(正)変換する.

Original external subprogram is wt_module#wz_wt

Function :
wz_xvz :real(8), dimension((nm+1)*(nm+1),0:km)
: (out) 2 次元球面調和函数スペクトル・動径格子点データ
xvz :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

3 次元格子データから水平スペクトル・動径格子点データへ(正)変換する.

[Source]

    function wz_xvz(xvz)
      !
      ! 3 次元格子データから水平スペクトル・動径格子点データへ(正)変換する.
      !
      real(8), dimension((nm+1)*(nm+1),0:km)             :: wz_xvz
      !(out) 2 次元球面調和函数スペクトル・動径格子点データ

      real(8), dimension(0:im-1,jc,0:km), intent(in)         :: xvz
      !(in) 3 次元経度緯度動径格子点データ

      wz_xvz = wa_xva(xvz)

    end function wz_xvz
wz_xyz( xyz ) result(wz_xyz)
Function :
wz_xyz :real(8), dimension((nm+1)*(nm+1),0:km)
: (out) 2 次元球面調和函数スペクトル・動径格子点データ
xyz :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

3 次元格子データから水平スペクトル・動径格子点データへ(正)変換する.

Original external subprogram is wt_module#wz_xyz

Function :
x_AvrLatRad_xvz :real(8), dimension(0:im-1)
: (out) 緯度動径(子午面)平均された 1 次元経度格子点データ
xvz :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

緯度動径(子午面)積分

3 次元格子点データの緯度動径(子午面)平均

3 次元データ f(λ,φ,r) に対して

   ∫f(λ,,r) r^2cosφ dφdr /(2(r[o]^3-r[i]^3)/3)

を計算する.

[Source]

    function x_AvrLatRad_xvz(xvz)  ! 緯度動径(子午面)積分
      !
      ! 3 次元格子点データの緯度動径(子午面)平均
      !
      ! 3 次元データ f(λ,φ,r) に対して
      !
      !    ∫f(λ,,r) r^2cosφ dφdr /(2(r[o]^3-r[i]^3)/3) 
      !
      ! を計算する.
      !
      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz
      !(in) 3 次元経度緯度動径格子点データ

      real(8), dimension(0:im-1)     :: x_AvrLatRad_xvz
      !(out) 緯度動径(子午面)平均された 1 次元経度格子点データ

      x_AvrLatRad_xvz = x_IntLatRad_xvz(xvz) /( sum(y_Lat_Weight)*sum(z_Rad_Weight) )

    end function x_AvrLatRad_xvz
x_AvrLatRad_xyz( xyz ) result(x_AvrLatRad_xyz)
Function :
x_AvrLatRad_xyz :real(8), dimension(0:im-1)
: (out) 緯度動径(子午面)平均された 1 次元経度格子点データ
xyz :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

緯度動径(子午面)積分

3 次元格子点データの緯度動径(子午面)平均

3 次元データ f(λ,φ,r) に対して

   ∫f(λ,,r) r^2cosφ dφdr /(2(r[o]^3-r[i]^3)/3)

を計算する.

Original external subprogram is wt_module#x_AvrLatRad_xyz

x_AvrLat_xv( xv_data ) result(x_AvrLat_xv)
Function :
x_AvrLat_xv(im) :real(8)
: (out) 平均された 1 次元経度(X)格子点データ
xv_data(0:im-1,jc) :real(8), intent(in)
: 格子点(0:im-1,jc) (in) 2 次元経度緯度格子点データ(0:im-1,jc)

2 次元緯度経度格子点データの緯度(Y)方向平均(1 層用).

実際には格子点データ各点毎に v_Lat_Weight をかけた総和を計算し, y_Lat_Weight の総和で割ることで平均している.

Original external subprogram is w_integral_mpi_module#x_AvrLat_xv

x_AvrLat_xy( xy_data ) result(x_AvrLat_xy)
Function :
x_AvrLat_xy(im) :real(8)
: (out) 平均された 1 次元経度(X)格子点データ
xy_data(0:im-1,1:jm) :real(8), intent(in)
: (in) 2 次元経度緯度格子点データ(0:im-1,1:jm)

2 次元緯度経度格子点データの緯度(Y)方向平均(1 層用).

実際には格子点データ各点毎に y_Y_Weight をかけた総和を計算し, y_Y_Weight の総和で割ることで平均している.

Original external subprogram is wt_module#x_AvrLat_xy

x_AvrRad_xz( xz ) result(x_AvrRad_xz)
Function :
x_AvrRad_xz :real(8), dimension(0:im-1)
: (out) 動径平均された 1 次元経度格子点データ
xz :real(8), dimension(0:im-1,0:km), intent(in)
: (in) 2 次元緯度動径格子点データ

動径積分

2 次元(XZ)格子点データの動径方向域平均.

2 次元データ f(λ,r) に対して

  ∫f(λ,r) r^2dr /((r[o]^3-r[i]^3)/3)

を計算する.

Original external subprogram is wt_module#x_AvrRad_xz

Function :
x_IntLatRad_xvz :real(8), dimension(0:im-1)
: (out) 緯度動径(子午面)積分された 1 次元経度格子点データ
xvz :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

3 次元格子点データの緯度動径(子午面)積分

3 次元データ f(λ,φ,r) に対して

   ∫f(λ,φ,r) r^2cosφ dφdr

を計算する.

[Source]

    function x_IntLatRad_xvz(xvz)
      !
      ! 3 次元格子点データの緯度動径(子午面)積分
      !
      ! 3 次元データ f(λ,φ,r) に対して
      !
      !    ∫f(λ,φ,r) r^2cosφ dφdr 
      !
      ! を計算する.
      !
      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz
      !(in) 3 次元経度緯度動径格子点データ

      real(8), dimension(0:im-1)     :: x_IntLatRad_xvz
      !(out) 緯度動径(子午面)積分された 1 次元経度格子点データ

      real(8), dimension(0:im-1)     :: x_IntLatRadTMP
      integer :: j, k

      x_IntLatRad_xvz = 0
      do k=0,km
         do j=1,jc
            x_IntLatRad_xvz = x_IntLatRad_xvz + xvz(:,j,k) * v_Lat_Weight(j) * z_Rad_Weight(k)
         enddo
      enddo

      x_IntLatRadTmp=x_IntLatRad_xvz
      CALL MPI_ALLREDUCE(x_IntLatRadTMP,x_IntLatRad_xvz,im,MPI_REAL8, MPI_SUM,MPI_COMM_WORLD,IERR)

    end function x_IntLatRad_xvz
x_IntLatRad_xyz( xyz ) result(x_IntLatRad_xyz)
Function :
x_IntLatRad_xyz :real(8), dimension(0:im-1)
: (out) 緯度動径(子午面)積分された 1 次元経度格子点データ
xyz :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

3 次元格子点データの緯度動径(子午面)積分

3 次元データ f(λ,φ,r) に対して

   ∫f(λ,φ,r) r^2cosφ dφdr

を計算する.

Original external subprogram is wt_module#x_IntLatRad_xyz

x_IntLat_xv( xv_data ) result(x_IntLat_xv)
Function :
x_IntLat_xv(0:im-1) :real(8)
: (out) 積分された 1 次元経度(X)格子点データ
xv_data(0:im-1,jc) :real(8), intent(in)
: (in) 2 次元経度緯度格子点データ(0:im-1,jc)

2 次元緯度経度格子点データの緯度(Y)方向積分(1 層用).

実際には格子点データ各点毎に v_Y_Weight をかけた総和を計算している.

Original external subprogram is w_integral_mpi_module#x_IntLat_xv

x_IntLat_xy( xy_data ) result(x_IntLat_xy)
Function :
x_IntLat_xy(0:im-1) :real(8)
: (out) 積分された 1 次元経度(X)格子点データ
xy_data(0:im-1,1:jm) :real(8), intent(in)
: (in) 2 次元経度緯度格子点データ(0:im-1,1:jm)

2 次元緯度経度格子点データの緯度(Y)方向積分(1 層用).

実際には格子点データ各点毎に y_Y_Weight をかけた総和を計算している.

Original external subprogram is wt_module#x_IntLat_xy

x_IntRad_xz( xz ) result(x_IntRad_xz)
Function :
x_IntRad_xz :real(8), dimension(0:im-1)
: (out) 動径積分された 1 次元経度格子点データ
xz :real(8), dimension(0:im-1,0:km), intent(in)
: (in) 2 次元緯度動径格子点データ

2 次元(XZ)格子点データの動径方向域積分.

2 次元データ f(λ,r) に対して ∫f(λ,r) r^2dr を計算する.

Original external subprogram is wt_module#x_IntRad_xz

x_Lon
Variable :
x_Lon(:) :real(8), allocatable
: 緯度経度

Original external subprogram is wt_module#x_Lon

x_Lon_Weight
Variable :
x_Lon_Weight(:) :real(8), allocatable
: 座標重み

Original external subprogram is wt_module#x_Lon_Weight

Function :
xv_AvrRad_xvz :real(8), dimension(0:im-1,jc)
: 水平格子点データ (out) 動径平均された 2 次元経度緯度(水平, 球面)格子点データ
xvz :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

3 次元格子点データの動径方向域平均.

3 次元データ f(λ,φ,r) に対して

   ∫f(λ,φ,r) r^2dr/((r[o]^3-r[i]^3)/3)

を計算する.

[Source]

    function xv_AvrRad_xvz(xvz)
      !
      ! 3 次元格子点データの動径方向域平均.
      !
      ! 3 次元データ f(λ,φ,r) に対して 
      !
      !    ∫f(λ,φ,r) r^2dr/((r[o]^3-r[i]^3)/3) 
      !
      ! を計算する.
      !
      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz
      !(in) 3 次元経度緯度動径格子点データ

      real(8), dimension(0:im-1,jc)  :: xv_AvrRad_xvz          ! 水平格子点データ
      !(out) 動径平均された 2 次元経度緯度(水平, 球面)格子点データ

      xv_AvrRad_xvz = xv_IntRad_xvz(xvz)/sum(z_Rad_Weight)

    end function xv_AvrRad_xvz
Function :
xv_IntRad_xvz :real(8), dimension(0:im-1,jc)
: (out) 動径積分された 2 次元経度緯度(水平, 球面)格子点データ
xvz :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

動径積分

3 次元格子点データの動径方向域積分.

3 次元データ f(λ,φ,r) に対して∫f(λ,φ,r) r^2dr を計算する.

[Source]

    function xv_IntRad_xvz(xvz)  ! 動径積分
      !
      ! 3 次元格子点データの動径方向域積分.
      !
      ! 3 次元データ f(λ,φ,r) に対して∫f(λ,φ,r) r^2dr を計算する.
      !
      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz
      !(in) 3 次元経度緯度動径格子点データ

      real(8), dimension(0:im-1,jc)  :: xv_IntRad_xvz
      !(out) 動径積分された 2 次元経度緯度(水平, 球面)格子点データ

      integer :: k

      xv_IntRad_xvz = 0.0d0
      do k=0,km
         xv_IntRad_xvz(:,:) = xv_IntRad_xvz(:,:) + xvz(:,:,k) * z_Rad_Weight(k) 
      enddo
    end function xv_IntRad_xvz
xv_Lat
Variable :
xv_Lat(:,:) :real(8), allocatable

Original external subprogram is w_base_mpi_module#xv_Lat

xv_Lon
Variable :
xv_Lon(:,:) :real(8), allocatable

Original external subprogram is w_base_mpi_module#xv_Lon

xv_w( w_data, [ipow], [iflag] ) result(xv_w)
Function :
xv_w(0:im-1,jc) :real(8)
: (out) 格子点データ
w_data((nm+1)*(nm+1)) :real(8), intent(in)
: (in) スペクトルデータ
ipow :integer, intent(in), optional
: (in) 作用させる 1/cosφ の次数. 省略時は 0.
iflag :integer, intent(in), optional
: (in) 変換の種類
    0 : 通常の正変換
    1 : 経度微分を作用させた正変換
   -1 : 緯度微分を作用させた正変換
    2 : sinφを作用させた正変換
    省略時は 0.

スペクトルデータから分散格子データへ変換する(1 層用).

Original external subprogram is w_base_mpi_module#xv_w

Function :
xvz_Div_xvz_xvz_xvz :real(8), dimension(0:im-1,jc,0:km)
: (out) ベクトル場の発散
xvz_Vlon :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) ベクトル場の経度成分
xvz_Vlat :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) ベクトル場の緯度成分
xvz_Vrad :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) ベクトル場の動径成分

ベクトル成分である 3 つの格子データに発散を作用させる.

第 1, 2 ,3 引数(u,v,w)がそれぞれベクトルの経度成分, 緯度成分, 動径成分を表す.

極の特異性を回避するためにベクトル場に cosφ/r の重みをかけて 計算している.

     div V = (r/cosφ)・div (Vcosφ/r) + V_φtanφ/r + V_r/r

[Source]

    function xvz_Div_xvz_xvz_xvz(xvz_Vlon,xvz_Vlat,xvz_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,jc,0:km), intent(in) :: xvz_Vlon
      !(in) ベクトル場の経度成分

      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz_Vlat
      !(in) ベクトル場の緯度成分

      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz_Vrad
      !(in) ベクトル場の動径成分

      real(8), dimension(0:im-1,jc,0:km)             :: xvz_Div_xvz_xvz_xvz
      !(out) ベクトル場の発散

      xvz_Div_xvz_xvz_xvz = xvz_Rad/cos(xvz_Lat) * xvz_wt(wt_Div_xvz_xvz_xvz(xvz_VLon*cos(xvz_Lat)/xvz_Rad, xvz_VLat*cos(xvz_Lat)/xvz_Rad, xvz_VRad*cos(xvz_Lat)/xvz_Rad )) + xvz_VLat*tan(xvz_Lat)/xvz_Rad + xvz_VRad/xvz_Rad

    end function xvz_Div_xvz_xvz_xvz
Function :
xvz_GradLat_wt :real(8), dimension(0:im-1,jc,0:km)
: (out) 勾配型緯度微分を作用された 2 次元スペクトルデータ
wt :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) 2 次元球面調和函数チェビシェフスペクトルデータ

スペクトルデータに勾配型経度微分 1/r ∂/∂φ を作用させる.

[Source]

    function xvz_GradLat_wt(wt) 
      !
      ! スペクトルデータに勾配型経度微分 1/r ∂/∂φ を作用させる.
      !
      real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt
      !(in) 2 次元球面調和函数チェビシェフスペクトルデータ

      real(8), dimension(0:im-1,jc,0:km)                     :: xvz_GradLat_wt
      !(out) 勾配型緯度微分を作用された 2 次元スペクトルデータ

      xvz_GradLat_wt = xva_GradLat_wa(wz_wt(wt))/xvz_Rad
    end function xvz_GradLat_wt
Function :
xvz_GradLon_wt :real(8), dimension(0:im-1,jc,0:km)
: (out) 勾配型経度微分を作用された 2 次元スペクトルデータ
wt :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) 2 次元球面調和函数チェビシェフスペクトルデータ

スペクトルデータに勾配型経度微分 1/rcosφ・∂/∂λ を作用させる.

[Source]

    function xvz_GradLon_wt(wt)
      !
      ! スペクトルデータに勾配型経度微分 1/rcosφ・∂/∂λ
      ! を作用させる.
      !
      real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt
      !(in) 2 次元球面調和函数チェビシェフスペクトルデータ

      real(8), dimension(0:im-1,jc,0:km)                     :: xvz_GradLon_wt
      !(out) 勾配型経度微分を作用された 2 次元スペクトルデータ

      xvz_GradLon_wt = xva_GradLon_wa(wz_wt(wt))/xvz_Rad

    end function xvz_GradLon_wt
Function :
xvz_KGrad_wt :real(8), dimension(0:im-1,jc,0:km)
: (out) 軸方向微分を作用された 2 次元スペクトルデータ
wt :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) 2 次元球面調和函数チェビシェフスペクトルデータ

k・▽ = cosφ/r ∂/∂φ + sinφ∂/∂r

入力スペクトルデータに対応する格子データに軸方向微分

   k・▽ = cosφ/r ∂/∂φ + sinφ∂/∂r

を作用させた格子データが返される. ここでベクトル k は球の中心から北極向きの単位ベクトルである.

[Source]

    function xvz_KGrad_wt(wt)    ! k・▽ = cosφ/r ∂/∂φ + sinφ∂/∂r
      !
      ! 入力スペクトルデータに対応する格子データに軸方向微分 
      !
      !    k・▽ = cosφ/r ∂/∂φ + sinφ∂/∂r 
      !
      ! を作用させた格子データが返される. 
      ! ここでベクトル k は球の中心から北極向きの単位ベクトルである.
      !
      real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt
      !(in) 2 次元球面調和函数チェビシェフスペクトルデータ

      real(8), dimension(0:im-1,jc,0:km)                     :: xvz_KGrad_wt
      !(out) 軸方向微分を作用された 2 次元スペクトルデータ

      xvz_KGrad_wt =  cos(xvz_Lat)*xvz_GradLat_wt(wt) + sin(xvz_Lat)*xvz_wt(wt_Drad_wt(wt))
    end function xvz_KGrad_wt
xvz_LAT
Variable :
xvz_LAT :real(8), dimension(:,:,:), allocatable
: 座標
xvz_LON
Variable :
xvz_LON :real(8), dimension(:,:,:), allocatable
: 座標
xvz_RAD
Variable :
xvz_RAD :real(8), dimension(:,:,:), allocatable
: 座標
Function :
xvz_RotLat_wt_wt :real(8), dimension(0:im-1,jc,0:km)
: (out) ベクトル場の回転の緯度成分
wt_Vlon :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) ベクトル場の経度成分
wt_Vrad :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) ベクトル場の動径成分

ベクトル場の経度成分, 動径成分である第 1, 2 引数 Vlon, Vrad から 回転の緯度成分

   1/r ∂(r Vlon)/∂r - 1/rcosφ・∂Vrad/∂λ

を計算する.

[Source]

    function xvz_RotLat_wt_wt(wt_Vlon,wt_Vrad) 
      !
      ! ベクトル場の経度成分, 動径成分である第 1, 2 引数 Vlon, Vrad から
      ! 回転の緯度成分 
      !
      !    1/r ∂(r Vlon)/∂r - 1/rcosφ・∂Vrad/∂λ
      !
      ! を計算する.
      !
      real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_Vlon
      !(in) ベクトル場の経度成分

      real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_Vrad
      !(in) ベクトル場の動径成分

      real(8), dimension(0:im-1,jc,0:km)                     :: xvz_RotLat_wt_wt
      !(out) ベクトル場の回転の緯度成分

        xvz_RotLat_wt_wt =   xvz_wt(wt_RotRad_wt(wt_Vlon)) - xvz_GradLon_wt(wt_Vrad) 

    end function xvz_RotLat_wt_wt
Function :
xvz_RotLon_wt_wt :real(8), dimension(0:im-1,jc,0:km)
: (out) ベクトル場の回転の経度成分
wt_Vrad :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) ベクトル場の動径成分
wt_Vlat :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) ベクトル場の緯度成分

ベクトル場の動径成分, 緯度成分である第 1, 2 引数 Vrad, Vlat から 回転の経度成分

   1/r ∂Vrad/∂φ-1/r ∂(r Vlat)/∂r を計算する.

を計算する

[Source]

    function xvz_RotLon_wt_wt(wt_Vrad,wt_Vlat) 
      !
      ! ベクトル場の動径成分, 緯度成分である第 1, 2 引数 Vrad, Vlat から
      ! 回転の経度成分 
      !
      !    1/r ∂Vrad/∂φ-1/r ∂(r Vlat)/∂r を計算する.
      !
      ! を計算する
      !
      real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_Vrad
      !(in) ベクトル場の動径成分

      real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt_Vlat
      !(in) ベクトル場の緯度成分

      real(8), dimension(0:im-1,jc,0:km)                     :: xvz_RotLon_wt_wt
      !(out) ベクトル場の回転の経度成分

        xvz_RotLon_wt_wt =   xvz_GradLat_wt(wt_Vrad) - xvz_wt(wt_RotRad_wt(wt_Vlat))

    end function xvz_RotLon_wt_wt
Function :
xvz_wt :real(8), dimension(0:im-1,jc,0:km)
: (out) 3 次元経度緯度動径格子点データ
wt :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) 2 次元球面調和函数チェビシェフスペクトルデータ

スペクトルデータから 3 次元格子点データへ(逆)変換する.

[Source]

    function xvz_wt(wt)
      !
      ! スペクトルデータから 3 次元格子点データへ(逆)変換する.
      !
      real(8), dimension(0:im-1,jc,0:km)                     :: xvz_wt
      !(out) 3 次元経度緯度動径格子点データ

      real(8), dimension((nm+1)*(nm+1),0:lm), intent(in) :: wt
      !(in) 2 次元球面調和函数チェビシェフスペクトルデータ

      xvz_wt = xva_wa(wz_wt(wt))

    end function xvz_wt
Function :
xvz_wz :real(8), dimension(0:im-1,jc,0:km)
: (out) 3 次元経度緯度動径格子点データ
wz :real(8), dimension((nm+1)*(nm+1),0:km), intent(in)
: (in) 2 次元球面調和函数スペクトル・動径格子点データ

水平スペクトル・動径格子点データから 3 次元格子点データへ(逆)変換する.

[Source]

    function xvz_wz(wz)
      !
      ! 水平スペクトル・動径格子点データから 3 次元格子点データへ(逆)変換する.
      !
      real(8), dimension(0:im-1,jc,0:km)                     :: xvz_wz
      !(out) 3 次元経度緯度動径格子点データ

      real(8), dimension((nm+1)*(nm+1),0:km), intent(in) :: wz
      !(in) 2 次元球面調和函数スペクトル・動径格子点データ

      xvz_wz = xva_wa(wz)

    end function xvz_wz
xy_AvrRad_xyz( xyz ) result(xy_AvrRad_xyz)
Function :
xy_AvrRad_xyz :real(8), dimension(0:im-1,1:jm)
: (out) 動径平均された 2 次元経度緯度(水平, 球面)格子点データ 水平格子点データ
xyz :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

3 次元格子点データの動径方向域平均.

3 次元データ f(λ,φ,r) に対して

   ∫f(λ,φ,r) r^2dr/((r[o]^3-r[i]^3)/3)

を計算する.

Original external subprogram is wt_module#xy_AvrRad_xyz

xy_IntRad_xyz( xyz ) result(xy_IntRad_xyz)
Function :
xy_IntRad_xyz :real(8), dimension(0:im-1,1:jm)
: (out) 動径積分された 2 次元経度緯度(水平, 球面)格子点データ
xyz :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

動径積分

3 次元格子点データの動径方向域積分.

3 次元データ f(λ,φ,r) に対して∫f(λ,φ,r) r^2dr を計算する.

Original external subprogram is wt_module#xy_IntRad_xyz

xy_Lat
Variable :
xy_Lat(:,:) :real(8), allocatable

Original external subprogram is wt_module#xy_Lat

xy_Lon
Variable :
xy_Lon(:,:) :real(8), allocatable

Original external subprogram is wt_module#xy_Lon

xy_w( w_data, [ipow], [iflag] ) result(xy_w)
Function :
xy_w(0:im-1,1:jm) :real(8)
: (out) 格子点データ
w_data((nm+1)*(nm+1)) :real(8), intent(in)
: (in) スペクトルデータ
ipow :integer, intent(in), optional
: (in) 作用させる 1/cosφ の次数. 省略時は 0.
iflag :integer, intent(in), optional
: (in) 変換の種類
   0 : 通常の正変換
  -1 : 経度微分を作用させた逆変換
   1 : 緯度微分 cosφ・∂/∂φ を作用させた逆変換
   2 : sinφを作用させた逆変換
   省略時は 0.

スペクトルデータから格子データへ変換する(1 層用).

Original external subprogram is wt_module#xy_w

xyz_Div_xyz_xyz_xyz( xyz_Vlon, xyz_Vlat, xyz_Vrad ) result(xyz_Div_xyz_xyz_xyz)
Function :
xyz_Div_xyz_xyz_xyz :real(8), dimension(0:im-1,1:jm,0:km)
: (out) ベクトル場の発散
xyz_Vlon :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) ベクトル場の経度成分
xyz_Vlat :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) ベクトル場の緯度成分
xyz_Vrad :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) ベクトル場の動径成分

ベクトル成分である 3 つの格子データに発散を作用させる.

第 1, 2 ,3 引数(u,v,w)がそれぞれベクトルの経度成分, 緯度成分, 動径成分を表す.

極の特異性を回避するためにベクトル場に cosφ/r の重みをかけて 計算している.

     div V = (r/cosφ)・div (Vcosφ/r) + V_φtanφ/r + V_r/r

Original external subprogram is wt_module#xyz_Div_xyz_xyz_xyz

xyz_GradLat_wt( wt ) result(xyz_GradLat_wt)
Function :
xyz_GradLat_wt :real(8), dimension(0:im-1,1:jm,0:km)
: (out) 勾配型緯度微分を作用された 2 次元スペクトルデータ
wt :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) 2 次元球面調和函数チェビシェフスペクトルデータ

スペクトルデータに勾配型経度微分 1/r ∂/∂φ を作用させる.

Original external subprogram is wt_module#xyz_GradLat_wt

xyz_GradLon_wt( wt ) result(xyz_GradLon_wt)
Function :
xyz_GradLon_wt :real(8), dimension(0:im-1,1:jm,0:km)
: (out) 勾配型経度微分を作用された 2 次元スペクトルデータ
wt :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) 2 次元球面調和函数チェビシェフスペクトルデータ

スペクトルデータに勾配型経度微分 1/rcosφ・∂/∂λ を作用させる.

Original external subprogram is wt_module#xyz_GradLon_wt

xyz_KGrad_wt( wt ) result(xyz_KGrad_wt)
Function :
xyz_KGrad_wt :real(8), dimension(0:im-1,1:jm,0:km)
: (out) 軸方向微分を作用された 2 次元スペクトルデータ
wt :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) 2 次元球面調和函数チェビシェフスペクトルデータ

k・▽ = cosφ/r ∂/∂φ + sinφ∂/∂r

入力スペクトルデータに対応する格子データに軸方向微分

   k・▽ = cosφ/r ∂/∂φ + sinφ∂/∂r

を作用させた格子データが返される. ここでベクトル k は球の中心から北極向きの単位ベクトルである.

Original external subprogram is wt_module#xyz_KGrad_wt

xyz_LAT
Variable :
xyz_LAT :real(8), dimension(:,:,:), allocatable
: 座標

Original external subprogram is wt_module#xyz_LAT

xyz_LON
Variable :
xyz_LON :real(8), dimension(:,:,:), allocatable
: 座標

Original external subprogram is wt_module#xyz_LON

xyz_RAD
Variable :
xyz_RAD :real(8), dimension(:,:,:), allocatable
: 座標

Original external subprogram is wt_module#xyz_RAD

xyz_RotLat_wt_wt( wt_Vlon, wt_Vrad ) result(xyz_RotLat_wt_wt)
Function :
xyz_RotLat_wt_wt :real(8), dimension(0:im-1,1:jm,0:km)
: (out) ベクトル場の回転の緯度成分
wt_Vlon :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) ベクトル場の経度成分
wt_Vrad :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) ベクトル場の動径成分

ベクトル場の経度成分, 動径成分である第 1, 2 引数 Vlon, Vrad から 回転の緯度成分

   1/r ∂(r Vlon)/∂r - 1/rcosφ・∂Vrad/∂λ

を計算する.

Original external subprogram is wt_module#xyz_RotLat_wt_wt

xyz_RotLon_wt_wt( wt_Vrad, wt_Vlat ) result(xyz_RotLon_wt_wt)
Function :
xyz_RotLon_wt_wt :real(8), dimension(0:im-1,1:jm,0:km)
: (out) ベクトル場の回転の経度成分
wt_Vrad :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) ベクトル場の動径成分
wt_Vlat :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) ベクトル場の緯度成分

ベクトル場の動径成分, 緯度成分である第 1, 2 引数 Vrad, Vlat から 回転の経度成分

   1/r ∂Vrad/∂φ-1/r ∂(r Vlat)/∂r を計算する.

を計算する

Original external subprogram is wt_module#xyz_RotLon_wt_wt

xyz_wt( wt ) result(xyz_wt)
Function :
xyz_wt :real(8), dimension(0:im-1,1:jm,0:km)
: (out) 3 次元経度緯度動径格子点データ
wt :real(8), dimension((nm+1)*(nm+1),0:lm), intent(in)
: (in) 2 次元球面調和函数チェビシェフスペクトルデータ

スペクトルデータから 3 次元格子点データへ(逆)変換する.

Original external subprogram is wt_module#xyz_wt

xyz_wz( wz ) result(xyz_wz)
Function :
xyz_wz :real(8), dimension(0:im-1,1:jm,0:km)
: (out) 3 次元経度緯度動径格子点データ
wz :real(8), dimension((nm+1)*(nm+1),0:km), intent(in)
: (in) 2 次元球面調和函数スペクトル・動径格子点データ

水平スペクトル・動径格子点データから 3 次元格子点データへ(逆)変換する.

Original external subprogram is wt_module#xyz_wz

Function :
xz_AvrLat_xvz :real(8), dimension(im,0:km)
: (out) 緯度平均された 2 次元緯度動径格子点データ
xvz :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

緯度積分

3 次元格子点データの緯度方向域平均.

3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r)cosφ dφ/2 を計算する.

[Source]

    function xz_AvrLat_xvz(xvz)  ! 緯度積分
      !
      ! 3 次元格子点データの緯度方向域平均.
      !
      ! 3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r)cosφ dφ/2 を計算する.
      !
      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz
      !(in) 3 次元経度緯度動径格子点データ

      real(8), dimension(im,0:km)  :: xz_AvrLat_xvz
      !(out) 緯度平均された 2 次元緯度動径格子点データ

      xz_AvrLat_xvz = xz_IntLat_xvz(xvz)/sum(y_Lat_Weight)

    end function xz_AvrLat_xvz
xz_AvrLat_xyz( xyz ) result(xz_AvrLat_xyz)
Function :
xz_AvrLat_xyz :real(8), dimension(0:im-1,0:km)
: (out) 緯度平均された 2 次元緯度動径格子点データ
xyz :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

緯度積分

3 次元格子点データの緯度方向域平均.

3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r)cosφ dφ/2 を計算する.

Original external subprogram is wt_module#xz_AvrLat_xyz

Function :
xz_IntLat_xvz :real(8), dimension(0:im-1,0:km)
: 緯度円格子点データ (out) 緯度積分された 2 次元緯度動径格子点データ
xvz :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

3 次元格子点データの緯度方向域積分.

3 次元データ f(λ,φ,r) に対して∫f(λ,φ,r) cosφ dφ を計算する.

[Source]

    function xz_IntLat_xvz(xvz)
      !
      ! 3 次元格子点データの緯度方向域積分.
      !
      ! 3 次元データ f(λ,φ,r) に対して∫f(λ,φ,r) cosφ dφ を計算する.
      !
      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz
      !(in) 3 次元経度緯度動径格子点データ

      real(8), dimension(0:im-1,0:km)  :: xz_IntLat_xvz        ! 緯度円格子点データ
      !(out) 緯度積分された 2 次元緯度動径格子点データ

      real(8), dimension(0:im-1,0:km)  :: xz_IntLatTMP
      integer :: j

      xz_IntLat_xvz = 0.0d0
      do j=1,jc
         xz_IntLat_xvz(:,:) = xz_IntLat_xvz(:,:) + xvz(:,j,:) * v_Lat_Weight(j)
      enddo
      xz_IntLatTmp=xz_IntLat_xvz
      CALL MPI_ALLREDUCE(xz_IntLatTMP,xz_IntLat_xvz,im*(km+1),MPI_REAL8, MPI_SUM,MPI_COMM_WORLD,IERR)

    end function xz_IntLat_xvz
xz_IntLat_xyz( xyz ) result(xz_IntLat_xyz)
Function :
xz_IntLat_xyz :real(8), dimension(0:im-1,0:km)
: (out) 緯度積分された 2 次元緯度動径格子点データ. 緯度円格子点データ
xyz :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

3 次元格子点データの緯度方向域積分.

3 次元データ f(λ,φ,r) に対して∫f(λ,φ,r) cosφ dφ を計算する.

Original external subprogram is wt_module#xz_IntLat_xyz

y_AvrLonRad_xyz( xyz ) result(y_AvrLonRad_xyz)
Function :
y_AvrLonRad_xyz :real(8), dimension(1:jm)
: (out) 経度動径(緯度円)平均された 1 次元緯度格子点データ
xyz :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

経度動径(緯度円)積分

3 次元格子点データの経度動径(緯度円)平均.

3 次元データ f(λ,φ,r) に対して

    ∫f(λ,φ,r) r^2dλdr /(2π(r[o]^3-r[i]^3)/3)

を計算する.

Original external subprogram is wt_module#y_AvrLonRad_xyz

y_AvrLon_xy( xy_data ) result(y_AvrLon_xy)
Function :
y_AvrLon_xy(1:jm) :real(8)
: (out) 平均された 1 次元緯度(Y)格子点
xy_data(0:im-1,1:jm) :real(8), intent(in)
: (in) 2 次元経度緯度格子点データ(0:im-1,1:jm)

2 次元緯度経度格子点データの経度(X)方向平均(1 層用).

実際には格子点データ各点毎に x_X_Weight をかけた総和を計算し, x_X_Weight の総和で割ることで平均している.

Original external subprogram is wt_module#y_AvrLon_xy

y_AvrRad_yz( yz ) result(y_AvrRad_yz)
Function :
y_AvrRad_yz :real(8), dimension(1:jm)
: (out) 動径平均された 1 次元緯度格子点データ
yz :real(8), dimension(1:jm,0:km), intent(in)
: (in) 2 次元緯度動径(子午面)格子点データ

2 次元(YZ)格子点データの動径方向域平均.

2 次元データ f(φ,r) に対して ∫f(φ,r) r^2dr /((r[o]^3-r[i]^3)/3) を計算する.

Original external subprogram is wt_module#y_AvrRad_yz

y_IntLonRad_xyz( xyz ) result(y_IntLonRad_xyz)
Function :
y_IntLonRad_xyz :real(8), dimension(1:jm)
: (out) 経度動径(緯度円)積分された 1 次元緯度格子点データ
xyz :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

3 次元格子点データの経度動径(緯度円)積分.

3 次元データ f(λ,φ,r) に対して∫f(λ,φ,r) r^2dλdr を計算する.

Original external subprogram is wt_module#y_IntLonRad_xyz

y_IntLon_xy( xy_data ) result(y_IntLon_xy)
Function :
y_IntLon_xy(1:jm) :real(8)
: (out) 積分された 1 次元緯度(Y)格子点データ
xy_data(0:im-1,1:jm) :real(8), intent(in)
: (in) 2 次元経度緯度格子点データ(0:im-1,1:jm)

2 次元緯度経度格子点データの経度(X)方向積分(1 層用).

実際には格子点データ各点毎に x_X_Weight をかけた総和を計算している.

Original external subprogram is wt_module#y_IntLon_xy

y_IntRad_yz( yz ) result(y_IntRad_yz)
Function :
y_IntRad_yz :real(8), dimension(1:jm)
: (out) 動径積分された 1 次元緯度格子点データ
yz :real(8), dimension(1:jm,0:km), intent(in)
: (in) 2 次元緯度動径(子午面)格子点データ

動径積分

2 次元(YZ)格子点データの動径方向域積分.

2 次元データ f(φ,r) に対して∫f(φ,r) r^2dr を計算する.

Original external subprogram is wt_module#y_IntRad_yz

y_Lat
Variable :
y_Lat(:) :real(8), allocatable
: 緯度経度

Original external subprogram is wt_module#y_Lat

y_Lat_Weight
Variable :
y_Lat_Weight(:) :real(8), allocatable
: 座標重み

Original external subprogram is wt_module#y_Lat_Weight

yz_AvrLon_xyz( xyz ) result(yz_AvrLon_xyz)
Function :
yz_AvrLon_xyz :real(8), dimension(1:jm,0:km)
: (out) 経度方向(帯状)平均された 2 次元子午面格子点データ
xyz :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

経度(帯状)積分

3 次元格子点データの経度方向(帯状)平均.

3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r)dλ/2π を計算する.

Original external subprogram is wt_module#yz_AvrLon_xyz

yz_IntLon_xyz( xyz ) result(yz_IntLon_xyz)
Function :
yz_IntLon_xyz :real(8), dimension(1:jm,0:km)
: (out) 経度方向(帯状)積分された 2 次元子午面格子点データ
xyz :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

経度(帯状)積分

3 次元格子点データの経度方向(帯状)積分.

3 次元データ f(λ,φ,r) に対して ∫f(λ,φ,r)dλ を計算する.

Original external subprogram is wt_module#yz_IntLon_xyz

Function :
z_AvrLat_vz :real(8), dimension(0:km)
: (out) 緯度平均された 1 次元動径格子点データ
vz :real(8), dimension(jc,0:km), intent(in)
: (in) 2 次元緯度動径(子午面)格子点データ

2 次元(VZ)格子点データの緯度方向域平均.

2 次元データ f(φ,r) に対して ∫f(φ,r) cosφ dφ/2 を計算する.

[Source]

    function z_AvrLat_vz(vz)
      !
      ! 2 次元(VZ)格子点データの緯度方向域平均.
      !
      ! 2 次元データ f(φ,r) に対して ∫f(φ,r) cosφ dφ/2 を計算する.
      !
      real(8), dimension(jc,0:km), intent(in) :: vz
      !(in) 2 次元緯度動径(子午面)格子点データ

      real(8), dimension(0:km)  :: z_AvrLat_vz
      !(out) 緯度平均された 1 次元動径格子点データ

      z_AvrLat_vz = z_IntLat_vz(vz)/sum(y_Lat_Weight)
    end function z_AvrLat_vz
z_AvrLat_yz( yz ) result(z_AvrLat_yz)
Function :
z_AvrLat_yz :real(8), dimension(0:km)
: (out) 緯度平均された 1 次元動径格子点データ
yz :real(8), dimension(1:jm,0:km), intent(in)
: (in) 2 次元緯度動径(子午面)格子点データ

2 次元(YZ)格子点データの緯度方向域平均.

2 次元データ f(φ,r) に対して ∫f(φ,r) cosφ dφ/2 を計算する.

Original external subprogram is wt_module#z_AvrLat_yz

Function :
z_AvrLonLat_xvz :real(8), dimension(0:km)
: (out) 緯度経度(水平, 球面)平均された 1 次元動径格子点データ
xvz :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

緯度経度(水平)積分

3 次元格子点データの緯度経度(水平, 球面)積分

3 次元データ f(λ,φ,r) に対して

   ∫f(λ,φ,r) cosφ dλdφ /4π

を計算する.

[Source]

    function z_AvrLonLat_xvz(xvz)  ! 緯度経度(水平)積分
      !
      ! 3 次元格子点データの緯度経度(水平, 球面)積分
      ! 
      ! 3 次元データ f(λ,φ,r) に対して
      !
      !    ∫f(λ,φ,r) cosφ dλdφ /4π 
      !
      ! を計算する.
      !
      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz
      !(in) 3 次元経度緯度動径格子点データ

      real(8), dimension(0:km)     :: z_AvrLonLat_xvz
      !(out) 緯度経度(水平, 球面)平均された 1 次元動径格子点データ

      z_AvrLonLat_xvz = z_IntLonLat_xvz(xvz) /(sum(x_Lon_Weight)*sum(y_Lat_Weight))

    end function z_AvrLonLat_xvz
z_AvrLonLat_xyz( xyz ) result(z_AvrLonLat_xyz)
Function :
z_AvrLonLat_xyz :real(8), dimension(0:km)
: (out) 緯度経度(水平, 球面)平均された 1 次元動径格子点データ
xyz :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

緯度経度(水平)積分

3 次元格子点データの緯度経度(水平, 球面)積分

3 次元データ f(λ,φ,r) に対して

   ∫f(λ,φ,r) cosφ dλdφ /4π

を計算する.

Original external subprogram is wt_module#z_AvrLonLat_xyz

z_AvrLon_xz( xz ) result(z_AvrLon_xz)
Function :
z_AvrLon_xz :real(8), dimension(0:km)
: (out) 経度平均された 1 次元動径格子点データ
xz :real(8), dimension(0:im-1,0:km), intent(in)
: (in) 2 次元緯度動径格子点データ

経度(帯状)積分

2 次元(XZ)格子点データの経度方向平均.

2 次元データ f(λ,r) に対して ∫f(λ,r)dλ/2π を計算する.

Original external subprogram is wt_module#z_AvrLon_xz

Function :
z_IntLat_vz :real(8), dimension(0:km)
: (out) 緯度積分された 1 次元動径格子点データ
vz :real(8), dimension(jc,0:km), intent(in)
: (in) 2 次元緯度動径(子午面)格子点データ

緯度積分

2 次元(VZ)格子点データの緯度方向域積分.

2 次元データ f(φ,r) に対して∫f(φ,r) cosφ dφ を計算する.

[Source]

    function z_IntLat_vz(vz)  ! 緯度積分
      !
      ! 2 次元(VZ)格子点データの緯度方向域積分.
      !
      ! 2 次元データ f(φ,r) に対して∫f(φ,r) cosφ dφ を計算する.
      !
      real(8), dimension(jc,0:km), intent(in) :: vz
      !(in) 2 次元緯度動径(子午面)格子点データ

      real(8), dimension(0:km)  :: z_IntLat_vz
      !(out) 緯度積分された 1 次元動径格子点データ

      real(8), dimension(0:km)  :: z_IntLatTMP
      integer :: j

      z_IntLat_vz = 0.0d0
      do j=1,jc
         z_IntLat_vz(:) = z_IntLat_vz(:) + vz(j,:) * v_Lat_Weight(j)
      enddo
      z_IntLatTmp=z_IntLat_vz
      CALL MPI_ALLREDUCE(z_IntLatTMP,z_IntLat_vz,km+1,MPI_REAL8, MPI_SUM,MPI_COMM_WORLD,IERR)

    end function z_IntLat_vz
z_IntLat_yz( yz ) result(z_IntLat_yz)
Function :
z_IntLat_yz :real(8), dimension(0:km)
: (out) 緯度積分された 1 次元動径格子点データ
yz :real(8), dimension(jm,0:km), intent(in)
: (in) 2 次元緯度動径(子午面)格子点データ

緯度積分

2 次元(YZ)格子点データの緯度方向域積分.

2 次元データ f(φ,r) に対して∫f(φ,r) cosφ dφ を計算する.

Original external subprogram is wt_module#z_IntLat_yz

Function :
z_IntLonLat_xvz :real(8), dimension(0:km)
: (out) 緯度経度(水平, 球面)積分された 1 次元動径格子点データ
xvz :real(8), dimension(0:im-1,jc,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

緯度経度(水平)積分

3 次元格子点データの緯度経度(水平, 球面)積分

3 次元データ f(λ,φ,r) に対して

   ∫f(λ,φ,r) cosφ dλdφ

を計算する.

[Source]

    function z_IntLonLat_xvz(xvz)  ! 緯度経度(水平)積分
      !
      ! 3 次元格子点データの緯度経度(水平, 球面)積分
      ! 
      ! 3 次元データ f(λ,φ,r) に対して
      !
      !    ∫f(λ,φ,r) cosφ dλdφ 
      !
      ! を計算する.
      !
      real(8), dimension(0:im-1,jc,0:km), intent(in) :: xvz
      !(in) 3 次元経度緯度動径格子点データ

      real(8), dimension(0:km)     :: z_IntLonLat_xvz
      !(out) 緯度経度(水平, 球面)積分された 1 次元動径格子点データ

      real(8), dimension(0:km)     :: z_IntLonLatTMP
      integer :: i, j

      z_IntLonLat_xvz = 0
      do j=1,jc
         do i=0,im-1
            z_IntLonLat_xvz = z_IntLonLat_xvz + xvz(i,j,:) * x_Lon_Weight(i) * v_Lat_Weight(j)
         enddo
      enddo

      z_IntLonLatTmp=z_IntLonLat_xvz
      CALL MPI_ALLREDUCE(z_IntLonLatTMP,z_IntLonLat_xvz,km+1,MPI_REAL8, MPI_SUM,MPI_COMM_WORLD,IERR)

    end function z_IntLonLat_xvz
z_IntLonLat_xyz( xyz ) result(z_IntLonLat_xyz)
Function :
z_IntLonLat_xyz :real(8), dimension(0:km)
: (out) 緯度経度(水平, 球面)積分された 1 次元動径格子点データ
xyz :real(8), dimension(0:im-1,1:jm,0:km), intent(in)
: (in) 3 次元経度緯度動径格子点データ

緯度経度(水平)積分

3 次元格子点データの緯度経度(水平, 球面)積分

3 次元データ f(λ,φ,r) に対して

   ∫f(λ,φ,r) cosφ dλdφ

を計算する.

Original external subprogram is wt_module#z_IntLonLat_xyz

z_IntLon_xz( xz ) result(z_IntLon_xz)
Function :
z_IntLon_xz :real(8), dimension(0:km)
: (out) 経度積分された 1 次元動径格子点データ
xz :real(8), dimension(0:im-1,0:km), intent(in)
: (in) 2 次元緯度動径格子点データ

2 次元(XZ)格子点データの経度方向積分.

2 次元データ f(λ,r) に対して ∫f(λ,r)dλ を計算する.

Original external subprogram is wt_module#z_IntLon_xz

z_RAD
Variable :
g_X(:) :real(8), allocatable
: 格子点座標 km 次のチェビシェフ多項式の零点から定まる格子点

Original external subprogram is wt_module#z_RAD

z_RAD_WEIGHT
Variable :
g_X_Weight(:) :real(8), allocatable
: 格子点重み座標 各格子点における積分のための重みが格納してある

Original external subprogram is wt_module#z_RAD_WEIGHT