Class | w_base_mpi_module_sjpack |
In: |
libsrc/w_mpi_module_sjpack/w_base_mpi_module_sjpack.f90
|
spml/w_base_mpi_module_sjpack モジュールは球面上での 2 次元流体運動を 球面調和函数を用いたスペクトル法と MPI によって数値計算するための モジュール w_mpi_module_sjpack の下部モジュールであり, スペクトル法の 基本的なな Fortran90 関数を提供する. 内部で ISPACK の SJPACK-MPI の Fortran77 サブルーチンを呼んでいる. スペクトルデータおよび格子点データの格納方法や変換の詳しい計算法に ついては ISPACK/SJPACK のマニュアルを参照されたい.
Variable : | |||
mm =21 : | integer
|
Original external subprogram is w_base_module_sjpack#mm
Variable : | |||
nm =21 : | integer
|
Original external subprogram is w_base_module_sjpack#nm
Variable : | |||
nn =22 : | integer
|
Original external subprogram is w_base_module_sjpack#nn
Variable : | |||
np =1 : | integer
|
Original external subprogram is w_base_module_sjpack#np
Variable : | |||
openmp =.false. : | logical
|
Original external subprogram is w_base_module_sjpack#openmp
Subroutine : | |||
w_Psi((mm+1)*(mm+1)) : | real(8), intent(in)
| ||
w_Chi((mm+1)*(mm+1)) : | real(8), intent(in)
| ||
xv_U(0:im-1,1:jc) : | real(8), intent(out)
| ||
xv_V(0:im-1,1:jc) : | real(8), intent(out)
|
流線・ポテンシャル(スペクトルデータ)から速度場(格子データ)に (逆)変換する(1 層用, MPI)
スペクトル変換を用いず微分を計算するために, 変換回数が 2 回ですむ.
u cosφ = ∂χ/∂λ - cosφ∂ψ/∂φ, v cosφ = cosφ∂χ/∂φ + ∂ψ/∂λ
subroutine w_StreamPotential2VectorMPI(w_Psi, w_Chi, xv_U, xv_V) ! ! 流線・ポテンシャル(スペクトルデータ)から速度場(格子データ)に ! (逆)変換する(1 層用, MPI) ! ! スペクトル変換を用いず微分を計算するために, 変換回数が 2 回ですむ. ! ! u cosφ = ∂χ/∂λ - cosφ∂ψ/∂φ, ! v cosφ = cosφ∂χ/∂φ + ∂ψ/∂λ ! real(8), intent(in) :: w_Psi((mm+1)*(mm+1)) !(in) 流線関数 real(8), intent(in) :: w_Chi((mm+1)*(mm+1)) !(in) 速度ポテンシャル real(8), intent(out) :: xv_U(0:im-1,1:jc) !(out) 速度経度成分 real(8), intent(out) :: xv_V(0:im-1,1:jc) !(out) 速度緯度成分 real(8) :: w_Rdata((mm+4)*mm+2) ! 作業用スペクトルデータ(SJTS2G 出力用) real(8) :: w_Xdata((mm+1)*(mm+1)) ! 作業用スペクトルデータ(SJCS2X 出力用) real(8) :: w_Ydata((mm+4)*mm+2) ! 作業用スペクトルデータ(SJCS2Y 出力用) real(8) :: q(jm/2*7*np) ! 変換用作業配列 real(8) :: ws2(2*(nm+1)*np) ! 変換用作業配列 real(8) :: wg((im+2)*jm) ! 変換用作業配列 real(8) :: w((jm+1)*im) ! 変換用作業配列 logical :: first=.true. ! 初回判定スイッチ save first if ( .not. w_base_initialize ) then call MessageNotify('E','w_StreamPotential2VectorMPI', 'w_base_module not initialize yet.') endif if ( openmp .and. first ) then call MessageNotify('M','w_StreamPotential2Vector', 'OpenMP routine SJTSOG/SNPACK-MPI is used for spherical harmonic transformation.') endif first = .false. ! ! u cosφ = ∂χ/∂λ - cosφ∂ψ/∂φ の計算 ! call sjcs2x(mm,w_Chi,w_Xdata) call sjcs2y(mm,w_Psi,w_Ydata,c) call sjcrup(mm,nm,w_Xdata,w_Rdata) w_Rdata = w_Rdata - w_Ydata ! ! u の計算 ! if ( openmp ) then call sjtsog(mm,nm,nm,im,jc,w_Rdata,xv_U, it,t,p,q,r,ws2,wg,w,1) else call sjts2g(mm,nm,nm,im,jc,w_Rdata,xv_U, it,t,p,q,r,ws2,wg,w,1) endif ! ! v cosφ = cosφ∂χ/∂φ + ∂ψ/∂λ の計算 ! call sjcs2y(mm,w_Chi,w_Ydata,c) call sjcs2x(mm,w_Psi,w_Xdata) call sjcrup(mm,nm,w_Xdata,w_Rdata) w_Rdata= w_Rdata + w_Ydata ! ! v の計算 ! if ( openmp ) then call sjtsog(mm,nm,nm,im,jc,w_Rdata,xv_V, it,t,p,q,r,ws2,wg,w,1) else call sjts2g(mm,nm,nm,im,jc,w_Rdata,xv_V, it,t,p,q,r,ws2,wg,w,1) endif end subroutine w_StreamPotential2VectorMPI
Subroutine : | |||
xv_U(0:im-1,1:jc) : | real(8), intent(in)
| ||
xv_V(0:im-1,1:jc) : | real(8), intent(in)
| ||
w_Vor((mm+1)*(mm+1)) : | real(8), intent(out)
| ||
w_Div((mm+1)*(mm+1)) : | real(8), intent(out)
|
速度場(格子データ)から渦度・発散(スペクトルデータ)に (正)変換する(1 層用, MPI)
スペクトル変換を用いず微分を計算するために, 変換回数が 2 回ですむ.
ζ = 1/cosφ∂v/∂λ - 1/cosφ ∂(u cosφ)/∂φ D = 1/cosφ∂u/∂λ + 1/cosφ ∂(v cosφ)/∂φ
subroutine w_Vector2VorDivMPI(xv_U, xv_V, w_Vor, w_Div) ! ! 速度場(格子データ)から渦度・発散(スペクトルデータ)に ! (正)変換する(1 層用, MPI) ! ! スペクトル変換を用いず微分を計算するために, 変換回数が 2 回ですむ. ! ! ζ = 1/cosφ∂v/∂λ - 1/cosφ ∂(u cosφ)/∂φ ! D = 1/cosφ∂u/∂λ + 1/cosφ ∂(v cosφ)/∂φ ! real(8), intent(in) :: xv_U(0:im-1,1:jc) !(in) 速度経度成分 real(8), intent(in) :: xv_V(0:im-1,1:jc) !(in) 速度緯度成分 real(8), intent(out) :: w_Vor((mm+1)*(mm+1)) !(out) 流線関数 real(8), intent(out) :: w_Div((mm+1)*(mm+1)) !(out) 速度ポテンシャル real(8) :: w_Xdata((mm+1)*(mm+1)) ! 作業用スペクトルデータ(SJCS2X 出力用) real(8) :: w_Ydata((mm+4)*nm+2) ! 作業用スペクトルデータ(SJCY2S 出力用) real(8) :: w_Data1((mm+1)*(mm+1)) real(8) :: w_Data2((mm+1)*(mm+1)) real(8) :: q(jm/2*7*np) ! 変換用作業配列 real(8) :: ws2(2*(nm+1)*np+(2*NM+1-MM)*MM+NM+1) ! 変換用作業配列 real(8) :: wg((im+2)*jm) ! 変換用作業配列 real(8) :: w((jm+1)*im) ! 変換用作業配列 logical :: first=.true. ! 初回判定スイッチ save first if ( .not. w_base_initialize ) then call MessageNotify('E','w_Vector2VorDivMPI', 'w_base_module not initialize yet.') endif if ( openmp .and. first ) then call MessageNotify('M','w_Vector2VorDiv', 'OpenMP routine SJPGOS/SNPACK-MPI is used for spherical harmonic transformation.') endif first = .false. ! ! 1/cosφ∂u/∂λ, 1/cosφ ∂(u cosφ)/∂φ の計算 ! if ( openmp ) then call sjpgos(mm,nm,nm,im,jc,w_Ydata,xv_U, it,t,p,q,r,ws2,wg,w,1) else call sjpg2s(mm,nm,nm,im,jc,w_Ydata,xv_U, it,t,p,q,r,ws2,wg,w,1) endif call sjcrdn(mm,nm,w_Ydata,w_Xdata) call sjcs2x(mm,w_Xdata,w_Div) call sjcy2s(mm,w_Ydata,w_Data1,c) ! ! 1/cosφ∂v/∂λ, 1/cosφ ∂(v cosφ)/∂φ の計算 ! if ( openmp ) then call sjpgos(mm,nm,nm,im,jc,w_Ydata,xv_V, it,t,p,q,r,ws2,wg,w,1) else call sjpg2s(mm,nm,nm,im,jc,w_Ydata,xv_V, it,t,p,q,r,ws2,wg,w,1) endif call sjcrdn(mm,nm,w_Ydata,w_Xdata) call sjcs2x(mm,w_Xdata,w_Vor) call sjcy2s(mm,w_Ydata,w_Data2,c) ! ! 渦度・発散の計算 ! ! ζ = 1/cosφ∂v/∂λ - 1/cosφ ∂(u cosφ)/∂φ ! D = 1/cosφ∂u/∂λ + 1/cosφ ∂(v cosφ)/∂φ ! w_Vor = w_Vor - w_Data1 w_Div = w_Div + w_Data2 end subroutine w_Vector2VorDivMPI
Subroutine : | |||
xv_UCosLat(0:im-1,1:jc) : | real(8), intent(in)
| ||
xv_VCosLat(0:im-1,1:jc) : | real(8), intent(in)
| ||
w_Vor((mm+1)*(mm+1)) : | real(8), intent(out)
| ||
w_Div((mm+1)*(mm+1)) : | real(8), intent(out)
|
速度場(格子データ)から渦度・発散(スペクトルデータ)に (正)変換する(1 層用, MPI)
スペクトル変換を用いず微分を計算するために, 変換回数が 2 回ですむ.
ζ = 1/cosφ∂v/∂λ - 1/cosφ ∂(u cosφ)/∂φ D = 1/cosφ∂u/∂λ + 1/cosφ ∂(v cosφ)/∂φ
subroutine w_VectorCosLat2VorDivMPI(xv_UCosLat, xv_VCosLat, w_Vor, w_Div) ! ! 速度場(格子データ)から渦度・発散(スペクトルデータ)に ! (正)変換する(1 層用, MPI) ! ! スペクトル変換を用いず微分を計算するために, 変換回数が 2 回ですむ. ! ! ζ = 1/cosφ∂v/∂λ - 1/cosφ ∂(u cosφ)/∂φ ! D = 1/cosφ∂u/∂λ + 1/cosφ ∂(v cosφ)/∂φ ! real(8), intent(in) :: xv_UCosLat(0:im-1,1:jc) !(in) 速度経度成分 * cos(lat) real(8), intent(in) :: xv_VCosLat(0:im-1,1:jc) !(in) 速度緯度成分 * cos(lat) real(8), intent(out) :: w_Vor((mm+1)*(mm+1)) !(out) 流線関数 real(8), intent(out) :: w_Div((mm+1)*(mm+1)) !(out) 速度ポテンシャル real(8) :: w_Xdata((mm+1)*(mm+1)) ! 作業用スペクトルデータ(SJCS2X 出力用) real(8) :: w_Ydata((mm+4)*nm+2) ! 作業用スペクトルデータ(SJCY2S 出力用) real(8) :: w_Data1((mm+1)*(mm+1)) real(8) :: w_Data2((mm+1)*(mm+1)) real(8) :: q(jm/2*7*np) ! 変換用作業配列 real(8) :: ws2(2*(nm+1)*np+(2*NM+1-MM)*MM+NM+1) ! 変換用作業配列 real(8) :: wg((im+2)*jm) ! 変換用作業配列 real(8) :: w((jm+1)*im) ! 変換用作業配列 logical :: first=.true. ! 初回判定スイッチ save first if ( .not. w_base_initialize ) then call MessageNotify('E','w_VectorCosLat2VorDivMPI', 'w_base_module not initialize yet.') endif if ( openmp .and. first ) then call MessageNotify('M','w_VectorCosLat2VorDiv', 'OpenMP routine SJPGOS/SNPACK-MPI is used for spherical harmonic transformation.') endif first = .false. ! ! 1/cosφ∂u/∂λ, 1/cosφ ∂(u cosφ)/∂φ の計算 ! if ( openmp ) then call sjpgos(mm,nm,nm,im,jc,w_Ydata,xv_UCosLat, it,t,p,q,r,ws2,wg,w,2) else call sjpg2s(mm,nm,nm,im,jc,w_Ydata,xv_UCosLat, it,t,p,q,r,ws2,wg,w,2) endif call sjcrdn(mm,nm,w_Ydata,w_Xdata) call sjcs2x(mm,w_Xdata,w_Div) call sjcy2s(mm,w_Ydata,w_Data1,c) ! ! 1/cosφ∂v/∂λ, 1/cosφ ∂(v cosφ)/∂φ の計算 ! if ( openmp ) then call sjpgos(mm,nm,nm,im,jc,w_Ydata,xv_VCosLat, it,t,p,q,r,ws2,wg,w,2) else call sjpg2s(mm,nm,nm,im,jc,w_Ydata,xv_VCosLat, it,t,p,q,r,ws2,wg,w,2) endif call sjcrdn(mm,nm,w_Ydata,w_Xdata) call sjcs2x(mm,w_Xdata,w_Vor) call sjcy2s(mm,w_Ydata,w_Data2,c) ! ! 渦度・発散の計算 ! ! ζ = 1/cosφ∂v/∂λ - 1/cosφ ∂(u cosφ)/∂φ ! D = 1/cosφ∂u/∂λ + 1/cosφ ∂(v cosφ)/∂φ ! w_Vor = w_Vor - w_Data1 w_Div = w_Div + w_Data2 end subroutine w_VectorCosLat2VorDivMPI
Subroutine : |
モジュールの終了処理(割り付け配列の解放)をおこなう.
実際の使用には上位サブルーチン w_Finalize を用いること.
subroutine w_base_mpi_Finalize ! ! モジュールの終了処理(割り付け配列の解放)をおこなう. ! ! 実際の使用には上位サブルーチン w_Finalize を用いること. ! if ( .not. w_base_initialize ) then call MessageNotify('W','w_base_mpi_Finalize', 'w_base_mpi_module_sjpack not initialized yet') return endif deallocate(p) ! 変換用配列 deallocate(r) ! 変換用配列 deallocate(t) ! 変換用配列 deallocate(c) ! 変換用作業配列 deallocate(v_Lat,v_Lat_Weight) ! 格子点座標格納配列 deallocate(xv_Lon,xv_Lat) ! 格子点座標格納配列 deallocate(y) w_base_initialize = .false. call MessageNotify('M','w_base_mpi_Finalize', 'w_base_mpi_module_sjpack (2013/02/23) is finalized') end subroutine w_base_mpi_Finalize
Subroutine : |
スペクトル変換の格子点数, 波数および OPENMP 使用時の 最大スレッド数を設定する.
実際の使用には上位サブルーチン w_Initial を用いること.
subroutine w_base_mpi_Initial ! ! スペクトル変換の格子点数, 波数および OPENMP 使用時の ! 最大スレッド数を設定する. ! ! 実際の使用には上位サブルーチン w_Initial を用いること. ! integer :: i, j allocate(p(jm/2*(mm+4))) ! 変換用配列 allocate(r((mm+1)*(2*nm-mm-1)+1)) ! 変換用配列 allocate(t(im*6)) ! 変換用配列 allocate(c((mm+1)*(mm+1))) ! 変換用作業配列 ! 注意 : 別ルーチンによって w_base_Initial が呼んであることを仮定 call sjpini(mm,nm,jm,jc,im,p,r,it,t) call sjinic(mm,c) allocate(v_Lat(jc),v_Lat_Weight(jc)) ! 格子点座標格納配列 allocate(xv_Lon(0:im-1,jc),xv_Lat(0:im-1,jc)) ! 格子点座標格納配列 allocate(y(jc/2,mm+4)) y = reshape(p(1:jc/2*(mm+4)),(/jc/2,mm+4/)) do j=1,jc/2 v_Lat(jc/2+j) = asin(y(j,1)) ! 緯度座標 v_Lat(jc/2-j+1) = -asin(y(j,1)) ! 緯度座標 v_Lat_Weight(jc/2+j) = 2*y(j,2) ! 緯度重み(Gauss grid) v_Lat_Weight(jc/2-j+1) = 2*y(j,2) ! 緯度重み(Gauss grid) enddo do j=1,jc xv_Lon(:,j) = x_Lon enddo do i=0,im-1 xv_Lat(i,:) = v_Lat enddo w_base_initialize = .true. call MessageNotify('M','w_base_mpi_initial', 'w_base_mpi_module_sjpack (2013/02/23) is initialized') end subroutine w_base_mpi_Initial
Function : | |||
w_xv((mm+1)*(mm+1)) : | real(8)
| ||
xv_data(0:im-1,1:jc) : | real(8), intent(in)
| ||
ipow : | integer, intent(in), optional
| ||
iflag : | integer, intent(in), optional
|
格子データからスペクトルデータへ(正)変換する(1 層用).
function w_xv(xv_data,ipow,iflag) ! ! 格子データからスペクトルデータへ(正)変換する(1 層用). ! real(8) :: w_xv((mm+1)*(mm+1)) !(out) スペクトルデータ real(8), intent(in) :: xv_data(0:im-1,1:jc) !(in) 格子点データ integer, intent(in), optional :: ipow !(in) 変換時に同時に作用させる 1/cosφ の次数. 省略時は 0. integer, intent(in), optional :: iflag ! 変換の種類 ! 0 : 通常の正変換 ! -1 : 経度微分を作用させた正変換 ! 1 : 緯度微分 1/cosφ・∂(f cos^2φ)/∂φ を作用させた正変換 ! 2 : sinφを作用させた正変換 ! 省略時は 0. integer, parameter :: ipow_default = 0 ! スイッチデフォルト値 integer, parameter :: iflag_default = 0 ! スイッチデフォルト値 integer ipval, ifval real(8) :: w_Rdata((2*nn+1-mm)*mm+nn+1) ! 作業用スペクトルデータ(SJTS2G 出力用) real(8) :: w_Xdata((mm+1)*(mm+1)) ! 作業用スペクトルデータ(SJCS2X 出力用) real(8) :: w_Ydata((mm+4)*nm+2) ! 作業用スペクトルデータ(SJCY2S 出力用) real(8) :: q(jm/2*7*np) ! 変換用作業配列 real(8) :: ws(2*(nn+1)*np+(2*NN+1-MM)*MM+NN+1) ! 変換用作業配列 real(8) :: ws2(2*(nm+1)*np+(2*NM+1-MM)*MM+NM+1) ! 変換用作業配列 real(8) :: wg((im+2)*jm) ! 変換用作業配列 real(8) :: w((jm+1)*im) ! 変換用作業配列 logical :: first=.true. ! 初回判定スイッチ save first if ( .not. w_base_initialize ) then call MessageNotify('E','xv_w', 'w_base_mpi_module_sjpack not initialize yet.') endif if (present(ipow)) then ipval = ipow else ipval = ipow_default endif if (present(iflag)) then ifval = iflag else ifval = iflag_default endif if ( openmp .and. first ) then call MessageNotify('M','w_xv', 'OpenMP routine SJPGOS/SJPACK-MPI is used for spherical harmonic transformation.') endif if ( ifval == 0 ) then if ( openmp ) then call sjpgos(mm,nm,nn,im,jc,w_Rdata,xv_data, it,t,p,q,r,ws,wg,w,ipval) else call sjpg2s(mm,nm,nn,im,jc,w_Rdata,xv_data, it,t,p,q,r,ws,wg,w,ipval) endif call sjcrdn(mm,nn,w_Rdata,w_xv) else if ( ifval == -1 ) then if ( openmp ) then call sjpgos(mm,nm,nn,im,jc,w_Rdata,xv_data, it,t,p,q,r,ws,wg,w,ipval) else call sjpg2s(mm,nm,nn,im,jc,w_Rdata,xv_data, it,t,p,q,r,ws,wg,w,ipval) endif call sjcrdn(mm,nn,w_Rdata,w_Xdata) call sjcs2x(mm,w_Xdata,w_xv) else if ( ifval == 1 ) then if ( openmp ) then call sjpgos(mm,nm,nm,im,jc,w_Ydata,xv_data, it,t,p,q,r,ws2,wg,w,ipval) else call sjpg2s(mm,nm,nm,im,jc,w_Ydata,xv_data, it,t,p,q,r,ws2,wg,w,ipval) endif call sjcy2s(mm,w_Ydata,w_xv,c) else if ( ifval == 2 ) then if ( openmp ) then call sjpgos(mm,nm,nn,im,jc,w_Rdata,xv_data*sin(xv_Lat), it,t,p,q,r,ws,wg,w,ipval) else call sjpg2s(mm,nm,nn,im,jc,w_Rdata,xv_data*sin(xv_Lat), it,t,p,q,r,ws,wg,w,ipval) endif call sjcrdn(mm,nn,w_Rdata,w_xv) end if first = .false. end function w_xv
Function : | |||
xv_w(0:im-1,1:jc) : | real(8)
| ||
w_data((mm+1)*(mm+1)) : | real(8), intent(in)
| ||
ipow : | integer, intent(in), optional
| ||
iflag : | integer, intent(in), optional
|
スペクトルデータから格子データへ変換する(1 層用).
function xv_w(w_data,ipow,iflag) ! ! スペクトルデータから格子データへ変換する(1 層用). ! real(8) :: xv_w(0:im-1,1:jc) !(out) 格子点データ real(8), intent(in) :: w_data((mm+1)*(mm+1)) !(in) スペクトルデータ integer, intent(in), optional :: ipow !(in) 作用させる 1/cosφ の次数. 省略時は 0. integer, intent(in), optional :: iflag !(in) 変換の種類 ! 0 : 通常の正変換 ! -1 : 経度微分を作用させた逆変換 ! 1 : 緯度微分 cosφ・∂/∂φ を作用させた逆変換 ! 2 : sinφを作用させた逆変換 ! 省略時は 0. ! integer, parameter :: ipow_default = 0 integer, parameter :: iflag_default = 0 integer ipval, ifval real(8) :: w_Rdata((2*nn+1-mm)*mm+nn+1) ! 作業用スペクトルデータ(SJTS2G 出力用) real(8) :: w_Xdata((mm+1)*(mm+1)) ! 作業用スペクトルデータ(SJCS2X 出力用) real(8) :: w_Ydata((mm+4)*mm+2) ! 作業用スペクトルデータ(SJCS2Y 出力用) real(8) :: q(jm/2*7*np) ! 変換用作業配列 real(8) :: ws(2*(nn+1)*np) ! 変換用作業配列 real(8) :: ws2(2*(nm+1)*np) ! 変換用作業配列 real(8) :: wg((im+2)*jm) ! 変換用作業配列 real(8) :: w((jm+1)*im) ! 変換用作業配列 logical :: first=.true. ! 初回判定スイッチ save first if ( .not. w_base_initialize ) then call MessageNotify('E','xv_w', 'w_base_mpi_module_sjpack not initialize yet.') endif if (present(ipow)) then ipval = ipow else ipval = ipow_default endif if (present(iflag)) then ifval = iflag else ifval = iflag_default endif if ( openmp .and. first ) then call MessageNotify('M','xy_w', 'OpenMP routine SJTSOG/SJPACK-MPI is used for spherical harmonic transformation.') endif if ( ifval==0 ) then call sjcrup(mm,nn,w_data,w_Rdata) if ( openmp ) then call sjtsog(mm,nm,nn,im,jc,w_Rdata,xv_w, it,t,p,q,r,ws,wg,w,ipval) else call sjts2g(mm,nm,nn,im,jc,w_Rdata,xv_w, it,t,p,q,r,ws,wg,w,ipval) endif else if( ifval==-1 ) then call sjcs2x(mm,w_data,w_Xdata) call sjcrup(mm,nn,w_Xdata,w_Rdata) if ( openmp ) then call sjtsog(mm,nm,nn,im,jc,w_Rdata,xv_w, it,t,p,q,r,ws,wg,w,ipval) else call sjts2g(mm,nm,nn,im,jc,w_Rdata,xv_w, it,t,p,q,r,ws,wg,w,ipval) endif else if( ifval==1 ) then call sjcs2y(mm,w_data,w_Ydata,c) if ( openmp ) then call sjtsog(mm,nm,nm,im,jc,w_Ydata,xv_w, it,t,p,q,r,ws2,wg,w,ipval) else call sjts2g(mm,nm,nm,im,jc,w_Ydata,xv_w, it,t,p,q,r,ws2,wg,w,ipval) endif else if( ifval==2 ) then call sjcrup(mm,nn,w_data,w_Rdata) if ( openmp ) then call sjtsog(mm,nm,nn,im,jc,w_Rdata,xv_w, it,t,p,q,r,ws,wg,w,ipval) else call sjts2g(mm,nm,nn,im,jc,w_Rdata,xv_w, it,t,p,q,r,ws,wg,w,ipval) endif xv_w = xv_w * sin(xv_Lat) else call MessageNotify('E','xv_w','invalid value of iflag') endif first = .false. end function xv_w