!--
!----------------------------------------------------------------------
!     Copyright (c) 2016 Shin-ichi Takehiro. All rights reserved.
!----------------------------------------------------------------------
!表題  w_base_module_svpack
!
!  spml/w_base_module_svpack モジュールは球面上での 2 次元流体運動を
!  球面調和函数を用いたスペクトル法によって数値計算するためのモジュール 
!  w_module_svpack の下部モジュールであり, スペクトル計算の基本的な 
!  Fortran90 関数を提供する.
!
!  内部で ISPACK2 の SVPACK の Fortran77 サブルーチンを呼んでいる. 
!  スペクトルデータおよび格子点データの格納方法や変換の詳しい計算
!  法については ISPACK2/SVPACK のマニュアルを参照されたい.
!
!== 履歴
!
!      2016/03/09  竹広真一  
!
!      制限
!         ・変換する格子点データ, スペクトルデータの配列の大きさは決めうち
!         ・波数切断の仕方は三角波数切断に決めうち. 
!
!++
module w_base_module_svpack
  !
  != w_base_module_svpack
  !
  ! Authors:: Shin-ichi Takehiro, Youhei SASAKI
  ! Version:: $Id:$
  ! Copyright&License:: See COPYRIGHT[link:../COPYRIGHT]
  !
  !== 概要.
  !
  !  spml/w_base_module_svpack モジュールは球面上での 2 次元流体運動を
  !  球面調和函数を用いたスペクトル法によって数値計算するためのモジュール 
  !  w_module_svpack の下部モジュールであり, スペクトル計算の基本的な 
  !  Fortran90 関数を提供する.
  !
  !  内部で ISPACK2 の SVPACK の Fortran77 サブルーチンを呼んでいる. 
  !  スペクトルデータおよび格子点データの格納方法や変換の詳しい計算
  !  法については ISPACK2/SVPACK のマニュアルを参照されたい.
  !
  use dc_message, only : MessageNotify
  use iso_c_binding
  implicit none

  integer               :: im=64            ! 格子点の設定(東西)
  integer               :: jm=32            ! 格子点の設定(南北)
  integer               :: nm=21            ! 計算する最大の全波数の設定
  integer               :: nn=22            ! 切断波数(全波数)の設定
  integer               :: mm=21            ! 切断波数(東西波数)の設定
  integer               :: np=1             ! OPENMP 最大スレッド数

  integer, allocatable  :: it(:)            ! 変換用配列
  real(8), allocatable  :: t(:)             ! 変換用配列
  real(8), allocatable  :: p(:,:), r(:)     ! 変換用配列
  integer, allocatable  :: jc(:)            ! 作業配列
  real(8), allocatable  :: wp(:)            ! 変換用作業配列

  real(8), allocatable  :: c(:)             ! 作業配列
  
  real(8), allocatable  :: x_Lon(:), y_Lat(:)                ! 緯度経度
  real(8), allocatable  :: x_Lon_Weight(:), y_Lat_Weight(:)  ! 座標重み
  real(8), allocatable  :: xy_Lon(:,:), xy_Lat(:,:)

  real(8), pointer:: w(:),xy_pdata(:,:)
  type(C_PTR)     :: pw, pg
  
  logical               :: w_base_initialize=.false.   ! 初期化フラッグ

  real(8), parameter    :: pi=3.1415926535897932385D0

  private

  public im, jm, nn, mm, nm                   ! 格子点数, 切断波数
  public it, t, p, r, jc                      ! 変換用作業配列

  public w_base_Initial                       ! 初期化サブルーチン
  public w_base_Finalize                      ! 終了サブルーチン

  public x_Lon, y_Lat                         ! 格子座標
  public x_Lon_Weight, y_Lat_Weight           ! 格子座標重み
  public xy_Lon, xy_Lat                       ! 格子座標(im,jm)
  public l_nm, nm_l                           ! 波数格納位置
  public xy_w, w_xy                           ! 変換関数
  public w_StreamPotential2Vector             ! 流線ポテンシャルから速度場計算
  public w_Vector2VorDiv                      ! 速度場から渦度発散を計算
  public w_VectorCosLat2VorDiv                ! 速度場から渦度発散を計算

  interface l_nm
     module procedure l_nm_array00
     module procedure l_nm_array01
     module procedure l_nm_array10
     module procedure l_nm_array11
  end interface

  interface nm_l
     module procedure nm_l_int
     module procedure nm_l_array
  end interface

  save im, jm, nm, mm, nn                     ! 格子点数, 切断波数を記憶
  save it, t, p, r, jc                        ! 変換用配列を記憶
  save c                                      ! 変換用配列の大きさ
  save w, xy_pdata, pw, pg
  save w_base_initialize                      ! 初期化フラグ

  contains
  !--------------- 初期化 -----------------
    subroutine w_base_Initial(n_in,i_in,j_in,np_in)
      !
      ! スペクトル変換の格子点数, 波数および OPENMP 使用時の
      ! 最大スレッド数を設定する.
      !
      ! 実際の使用には上位サブルーチン w_Initial を用いること.
      !
      integer,intent(in) :: i_in              !(in) 格子点数(東西), 2の巾乗(<=2048)
      integer,intent(in) :: j_in              !(in) 格子点数(南北), 4 の倍数
      integer,intent(in) :: n_in              !(in) 切断全波数
      integer,intent(in), optional :: np_in   !(in) 互換性のためのダミー変数

      integer :: i, j

      im = i_in   ; jm = j_in
      nn = n_in   ; nm = n_in+1 ;  mm = n_in      ! default は三角波数切断

      if ( present(np_in) )then
         call MessageNotify('M','w_base_Initial', &
              'Optional argument "NP" is dummy in w_base_module_svpack')
      endif

      allocate(it(im/2))                      ! 変換用配列
      allocate(t(im))                         ! 変換用配列
      allocate(p(jm/2,2*mm+5))                ! 変換用配列
      allocate(r(((mm+1)*(2*nm-mm-1)+1)/4*3+(2*nm-mm)*(mm+1)/2))  ! 変換用配列
      allocate(jc((mm*(2*nm-mm-1)/8+mm)))     ! 変換用配列
      allocate(wp(jm/2*mm))                   ! 変換用作業配列

      allocate(c((mm+1)*(mm+1)))              ! 変換用配列
      
      allocate(x_Lon(0:im-1))                 ! 格子点座標格納配列(経度)
      allocate(x_Lon_Weight(0:im-1))
      allocate(xy_Lon(0:im-1,1:jm))
      allocate(y_Lat(1:jm))
      allocate(y_Lat_Weight(1:jm))             ! 格子点座標格納配列
      allocate(xy_Lat(0:im-1,1:jm))        ! 格子点座標格納配列

      call svinit(mm,nm,im,jm,it,t,p,r,jc,wp)
      
      call svinic(mm,c)

      do i=0,im-1
         x_Lon(i)  = 2*pi/im*i               ! 経度座標
         x_Lon_Weight(i) = 2*pi/im           ! 経度座標重み
      enddo

      do j=1,jm/2
         y_Lat(jm/2+j)   =  asin(p(j,1))        ! 緯度座標
         y_Lat(jm/2-j+1) = -asin(p(j,1))        ! 緯度座標
         y_Lat_Weight(jm/2+j)   = 2*p(j,2)      ! 緯度重み(Gauss grid)
         y_Lat_Weight(jm/2-j+1) = 2*p(j,2)      ! 緯度重み(Gauss grid)
      enddo

      do j=1,jm
         xy_Lon(:,j) = x_Lon
      enddo

      do i=0,im-1
         xy_Lat(i,:) = y_Lat
      enddo

      call mlallc(pg,jm*im)
      call mlallc(pw,jm*im*2)
      call c_f_pointer(pg, xy_pdata, [IM,JM])
      call c_f_pointer(pw, w, [jm*im*2])

      w_base_initialize = .true.

      call MessageNotify('M','w_base_initial',&
           'w_base_module_svpack (2016/03/09) for AVX/FMA is initialized')

    end subroutine w_base_Initial

  !--------------- 基本関数 -----------------

    function l_nm_array00(n,m)
      !
      ! 全波数(n)と東西波数(m)からそのスペクトルデータの格納位置を返す.
      ! 
      ! 引数 n,m がともに整数値の場合, 整数値を返す. 
      !
      integer               :: l_nm_array00   
      !(out) スペクトルデータの格納位置 

      integer, intent(in)   :: n     !(in) 全波数
      integer, intent(in)   :: m     !(in) 帯状波数  

      if ( .not. w_base_initialize ) then
         call MessageNotify('E','l_nm_array00',&
              'w_base_module_svpack not initialize yet. Use svnm2l routine in ISPACK2 directly.')
      else
         call svnm2l(nn,n,m,l_nm_array00)
      endif

    end function l_nm_array00

    function l_nm_array01(n,marray)           ! スペクトルデータの格納位置 
      !
      ! 全波数(n)と東西波数(m)からそのスペクトルデータの格納位置を返す.
      ! 
      ! 第 1 引数 n が整数, 第 2 引数 marray が整数 1 次元配列の場合, 
      ! marray と同じ大きさの 1 次元整数配列を返す. 
      !
      integer, intent(in)  :: n               !(in) 全波数
      integer, intent(in)  :: marray(:)       !(in) 帯状波数
      integer              :: l_nm_array01(size(marray))
      !(out) スペクトルデータ位置

      integer              :: i 

      do i=1, size(marray)
         l_nm_array01(i) = l_nm_array00(n,marray(i))
      enddo
    end function l_nm_array01

    function l_nm_array10(narray,m)
      !
      ! 全波数(n)と東西波数(m)からそのスペクトルデータの格納位置を返す.
      ! 
      ! 第 1 引数 narray が整数 1 次元配列, 第 2 引数  m が整数の場合, 
      ! narray と同じ大きさの 1 次元整数配列を返す. 
      !
      integer, intent(in)  :: narray(:)           !(in) 全波数  
      integer, intent(in)  :: m                   !(in) 帯状波数
      integer              :: l_nm_array10(size(narray))
      !(out) スペクトルデータ位置

      integer              :: i 

      do i=1, size(narray)
         l_nm_array10(i) = l_nm_array00(narray(i),m)
      enddo
    end function l_nm_array10

    function l_nm_array11(narray,marray)
      !
      ! 全波数(n)と東西波数(m)からそのスペクトルデータの格納位置を返す.
      ! 
      ! 第 1,2 引数 narray, marray がともに整数 1 次元配列の場合, 
      ! narray, marray と同じ大きさの 1 次元整数配列を返す. 
      ! narray, marray は同じ大きさでなければならない. 
      !
      integer, intent(in)  :: narray(:)          !(in) 全波数  
      integer, intent(in)  :: marray(:)          !(in) 帯状波数
      integer              :: l_nm_array11(size(narray))
      !(out) スペクトルデータ位置

      integer              :: i 

      if ( size(narray) .ne. size(marray) ) then
         call MessageNotify('E','l_nm_array11',&
              'dimensions of input arrays  n and m are different.')
      endif

      do i=1, size(narray)
         l_nm_array11(i) = l_nm_array00(narray(i),marray(i))
      enddo
    end function l_nm_array11

    function nm_l_int(l)
      ! 
      ! スペクトルデータの格納位置(l)から全波数(n)と東西波数(m)を返す.
      !
      ! 引数 l が整数値の場合, 対応する全波数と帯状波数を
      ! 長さ 2 の 1 次元整数値を返す. 
      ! nm_l(1) が全波数, nm_l(2) が帯状波数である. 
      !
      integer               :: nm_l_int(2)  !(out) 全波数, 帯状波数
      integer, intent(in)   :: l            !(in) スペクトルデータの格納位置
      
      if ( .not. w_base_initialize ) then
         call MessageNotify('E','nm_l_int',&
              'w_base_module_svpack not initialize yet. Use svl2nm routine in ISPACK2 directly.')
      else
         call svl2nm(nn,l,nm_l_int(1),nm_l_int(2))
      endif

    end function nm_l_int

    function nm_l_array(larray)
      ! 
      ! スペクトルデータの格納位置(l)から全波数(n)と東西波数(m)を返す.
      !
      ! 引数 larray が整数 1 次元配列の場合, 
      ! larray に対応する n, m を格納した 2 次元整数配列を返す. 
      ! nm_l_array(:,1) が全波数, nm_l_array(:,2) が帯状波数である. 
      !
      integer, intent(in)  :: larray(:)
      !(out) 全波数, 帯状波数

      integer              :: nm_l_array(size(larray),2)
      !(in) スペクトルデータの格納位置

      integer              :: i

      do i=1, size(larray)
         nm_l_array(i,:) = nm_l_int(larray(i))
      enddo
    end function nm_l_array

  !--------------- 基本変換 -----------------

    function xy_w(w_data,ipow,iflag)
      !
      ! スペクトルデータから格子データへ変換する(1 層用).
      !
      real(8)               :: xy_w(0:im-1,1:jm)
      !(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_Xdata((mm+1)*(mm+1))
      ! 作業用スペクトルデータ(SVCS2X 出力用)
      real(8)             :: w_Ydata((mm+4)*mm+2)
      ! 作業用スペクトルデータ(SVCS2Y 出力用)

      real(8)  :: ws((nn+1)*2*(mm+1))           ! 変換用作業配列

      if ( .not. w_base_initialize ) then
         call MessageNotify('E','xy_w',&
              'w_base_module_svpack 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 ( ifval==0 ) then
         call svts2g(mm,nm,nn,im,jm,w_data,xy_pdata,it,t,p,r,jc,ws,w,ipval)
         xy_w = xy_pdata
      else if( ifval==-1 ) then
         call svcs2x(mm,w_data,w_Xdata)
         call svts2g(mm,nm,nn,im,jm,w_Xdata,xy_pdata,it,t,p,r,jc,ws,w,ipval)
         xy_w = xy_pdata
      else if( ifval==1 ) then
         call svcs2y(mm,w_data,w_Ydata,c)
         call svts2g(mm,nm,nm,im,jm,w_Ydata,xy_pdata,it,t,p,r,jc,ws,w,ipval)
         xy_w = xy_pdata
      else if( ifval==2 ) then
         call svts2g(mm,nm,nn,im,jm,w_data,xy_pdata,it,t,p,r,jc,ws,w,ipval)
         xy_w = xy_pdata * sin(xy_Lat)
      else
         call MessageNotify('E','xy_w','invalid value of iflag')
      endif

    end function xy_w

    function w_xy(xy_data,ipow,iflag)
      !
      ! 格子データからスペクトルデータへ(正)変換する(1 層用).
      !
      real(8)               :: w_xy((mm+1)*(mm+1))
      !(out) スペクトルデータ

      real(8), intent(in)   :: xy_data(0:im-1,1:jm)
      !(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_Xdata((mm+1)*(mm+1))
      ! 作業用スペクトルデータ(SVCS2X 出力用)
      real(8)             :: w_Ydata((mm+4)*nm+2)
      ! 作業用スペクトルデータ(SVCY2S 出力用)

      real(8)  :: ws((nn+1)*2*(mm+1))           ! 変換用作業配列

      if ( .not. w_base_initialize ) then
         call MessageNotify('E','xy_w',&
              'w_base_module_svpack 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 ( ifval == 0 ) then
         xy_pdata = xy_data
         call svtg2s(mm,nm,nn,im,jm,w_xy,xy_pdata,it,t,p,r,jc,ws,w,ipval)
      else if ( ifval == -1 ) then
         xy_pdata = xy_data
         call svtg2s(mm,nm,nn,im,jm,w_Xdata,xy_pdata,it,t,p,r,jc,ws,w,ipval)
         call svcs2x(mm,w_Xdata,w_xy)
      else if ( ifval == 1 ) then
         xy_pdata = xy_data         
         call svtg2s(mm,nm,nm,im,jm,w_Ydata,xy_pdata,it,t,p,r,jc,ws,w,ipval)
         call svcy2s(mm,w_Ydata,w_xy,c)
      else if ( ifval == 2 ) then
         xy_pdata = xy_data*sin(xy_Lat)
         call svtg2s(mm,nm,nn,im,jm,w_xy,xy_pdata,it,t,p,r,jc,ws,w,ipval)
      else
         call MessageNotify('E','w_xy','invalid value of iflag')
      endif

    end function w_xy

  !----------- 速度, 渦度・発散, 流線・ポテンシャル計算 -------------

    subroutine w_StreamPotential2Vector(w_Psi, w_Chi, xy_U, xy_V)
      !
      ! 流線・ポテンシャル(スペクトルデータ)から速度場(格子データ)に
      ! (逆)変換する(1 層用)
      !
      ! スペクトル変換を用いず微分を計算するために, 変換回数が 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)   :: xy_U(0:im-1,1:jm)
      !(out) 速度経度成分
      real(8), intent(out)   :: xy_V(0:im-1,1:jm)
      !(out) 速度緯度成分

      real(8)             :: w_Rdata((mm+4)*mm+2)
      ! 作業用スペクトルデータ(SVTS2Y 出力用)
       real(8)             :: w_Xdata((mm+1)*(mm+1))
      ! 作業用スペクトルデータ(SVCS2X 出力用)
      real(8)             :: w_Ydata((mm+4)*mm+2)
      ! 作業用スペクトルデータ(SVCS2Y 出力用)

      real(8)  :: ws((nn+1)*2*(mm+1))           ! 変換用作業配列

      if ( .not. w_base_initialize ) then
         call MessageNotify('E','w_StreamPotential2Vector',&
              'w_base_module_svpack not initialize yet.')
      endif
      !
      !   u cosφ = ∂χ/∂λ - cosφ∂ψ/∂φ の計算
      !
      call svcs2x(mm,w_Chi,w_Xdata)
      call svcs2y(mm,w_Psi,w_Ydata,c)
      call svcrup(mm,nm,w_Xdata,w_Rdata)
      w_Rdata = w_Rdata - w_Ydata
      !
      !   u の計算
      !
      call svts2g(mm,nm,nm,im,jm,w_Rdata,xy_pdata,it,t,p,r,jc,ws,w,1)
      xy_U = xy_pdata
      !
      !   v cosφ = cosφ∂χ/∂φ + ∂ψ/∂λ の計算
      !
      call svcs2y(mm,w_Chi,w_Ydata,c)
      call svcs2x(mm,w_Psi,w_Xdata)
      call svcrup(mm,nm,w_Xdata,w_Rdata)
      w_Rdata= w_Rdata + w_Ydata
      !
      !   v の計算
      !
      call svts2g(mm,nm,nm,im,jm,w_Rdata,xy_pdata,it,t,p,r,jc,ws,w,1)
      xy_V = xy_pdata

    end subroutine w_StreamPotential2Vector

    subroutine w_Vector2VorDiv(xy_U, xy_V, w_Vor, w_Div)
      !
      ! 速度場(格子データ)から渦度・発散(スペクトルデータ)に
      ! (正)変換する(1 層用)
      ! 
      ! スペクトル変換を用いず微分を計算するために, 変換回数が 2 回ですむ. 
      !
      !   ζ = 1/cosφ∂v/∂λ - 1/cosφ ∂(u cosφ)/∂φ 
      !    D = 1/cosφ∂u/∂λ + 1/cosφ ∂(v cosφ)/∂φ
      !
      real(8), intent(in)   :: xy_U(0:im-1,1:jm)
      !(in) 速度経度成分
      real(8), intent(in)   :: xy_V(0:im-1,1:jm)
      !(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))
      ! 作業用スペクトルデータ(SVCS2X 出力用)
      real(8)             :: w_Ydata((mm+4)*nm+2)
      ! 作業用スペクトルデータ(SVCY2S 出力用)

      real(8)  :: ws((nn+1)*2*(mm+1))           ! 変換用作業配列

      real(8)  :: w_Data1((mm+1)*(mm+1))
      real(8)  :: w_Data2((mm+1)*(mm+1))

      if ( .not. w_base_initialize ) then
         call MessageNotify('E','w_Vector2VorDiv',&
              'w_base_module_svpack not initialize yet.')
      endif
      !
      ! 1/cosφ∂u/∂λ, 1/cosφ ∂(u cosφ)/∂φ の計算
      !
      xy_pdata = xy_U
      call svtg2s(mm,nm,nm,im,jm,w_Ydata,xy_pdata,it,t,p,r,jc,ws,w,1)
      call svcrdn(mm,nm,w_Ydata,w_Xdata)
      call svcs2x(mm,w_Xdata,w_Div)
      call svcy2s(mm,w_Ydata,w_Data1,c)
      !
      ! 1/cosφ∂v/∂λ, 1/cosφ ∂(v cosφ)/∂φ の計算
      !
      xy_pdata = xy_V
      call svtg2s(mm,nm,nm,im,jm,w_Ydata,xy_pdata,it,t,p,r,jc,ws,w,1)
      call svcrdn(mm,nm,w_Ydata,w_Xdata)
      call svcs2x(mm,w_Xdata,w_Vor)
      call svcy2s(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_Vector2VorDiv

    subroutine w_VectorCosLat2VorDiv(xy_UCosLat, xy_VCosLat, w_Vor, w_Div)
      !
      ! 速度場(格子データ,cos(lat) 重み)から渦度・発散(スペクトルデータ)に
      ! (正)変換する(1 層用)
      ! 
      ! スペクトル変換を用いず微分を計算するために, 変換回数が 2 回ですむ. 
      !
      !   ζ = 1/cosφ∂v/∂λ - 1/cosφ ∂(u cosφ)/∂φ 
      !    D = 1/cosφ∂u/∂λ + 1/cosφ ∂(v cosφ)/∂φ
      !
      real(8), intent(in)   :: xy_UCosLat(0:im-1,1:jm)
      !(in) 速度経度成分 * cos(lat)
      real(8), intent(in)   :: xy_VCosLat(0:im-1,1:jm)
      !(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))
      ! 作業用スペクトルデータ(SVCS2X 出力用)
      real(8)             :: w_Ydata((mm+4)*nm+2)
      ! 作業用スペクトルデータ(SVCY2S 出力用)

      real(8)  :: ws((nn+1)*2*(mm+1))           ! 変換用作業配列

      real(8)  :: w_Data1((mm+1)*(mm+1))
      real(8)  :: w_Data2((mm+1)*(mm+1))

      if ( .not. w_base_initialize ) then
         call MessageNotify('E','w_VectorCosLat2VorDiv',&
              'w_base_module_svpack not initialize yet.')
      endif
      !
      ! 1/cosφ∂u/∂λ, 1/cosφ ∂(u cosφ)/∂φ の計算
      !
      xy_pdata = xy_UCosLat
      call svtg2s(mm,nm,nm,im,jm,w_Ydata,xy_pdata,it,t,p,r,jc,ws,w,2)
      call svcrdn(mm,nm,w_Ydata,w_Xdata)
      call svcs2x(mm,w_Xdata,w_Div)
      call svcy2s(mm,w_Ydata,w_Data1,c)
      !
      ! 1/cosφ∂v/∂λ, 1/cosφ ∂(v cosφ)/∂φ の計算
      !
      xy_pdata = xy_VCosLat
      call svtg2s(mm,nm,nm,im,jm,w_Ydata,xy_pdata,it,t,p,r,jc,ws,w,2)
      call svcrdn(mm,nm,w_Ydata,w_Xdata)
      call svcs2x(mm,w_Xdata,w_Vor)
      call svcy2s(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_VectorCosLat2VorDiv
     

  !--------------- 終了処理 -----------------
    subroutine w_base_Finalize
      !
      ! モジュールの終了処理(割り付け配列の解放)をおこなう. 
      !
      ! 実際の使用には上位サブルーチン w_Finalize を用いること.
      !
      if ( .not. w_base_initialize ) then
         call MessageNotify('W','w_base_Finalize',&
              'w_base_module not initialized yet')
         return
      endif

      deallocate(it)                 ! 変換用配列
      deallocate(t)                  ! 変換用配列
      deallocate(p)                  ! 変換用配列
      deallocate(r)                  ! 変換用配列
      deallocate(jc)                 ! 変換用作業配列

      deallocate(wp)                 ! 変換用作業配列

      deallocate(c)                  ! 変換用作業配列

      nullify(w)                     ! 変換用作業配列
      nullify(xy_pdata)              ! 変換用作業配列
      call mlfree(pg)                ! 変換用作業配列
      call mlfree(pw)                ! 変換用作業配列

      deallocate(x_Lon)              ! 格子点座標格納配列(経度)
      deallocate(x_Lon_Weight)
      deallocate(xy_Lon)
      deallocate(y_Lat)
      deallocate(y_Lat_Weight)       ! 格子点座標格納配列
      deallocate(xy_Lat)             ! 格子点座標格納配列

      w_base_initialize = .false.

      call MessageNotify('M','w_base_Finalize',&
           'w_base_module_svpack (2016/03/09) for AVX/FMA is finalized')

    end subroutine w_base_Finalize

end module w_base_module_svpack
