!-------------------------------------- LICENCE BEGIN ------------------------------------
!Environment Canada - Atmospheric Science and Technology License/Disclaimer,
! version 3; Last Modified: May 7, 2008.
!This is free but copyrighted software; you can use/redistribute/modify it under the terms
!of the Environment Canada - Atmospheric Science and Technology License/Disclaimer
!version 3 or (at your option) any later version that should be found at:
!http://collaboration.cmc.ec.gc.ca/science/rpn.comm/license.html
!
!This software is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
!without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
!See the above mentioned License/Disclaimer for more details.
!You should have received a copy of the License/Disclaimer along with this software;
!if not, you can write to: EC-RPN COMM Group, 2121 TransCanada, suite 500, Dorval (Quebec),
!CANADA, H9P 1J3; or send e-mail to service.rpn@ec.gc.ca
!-------------------------------------- LICENCE END --------------------------------------
!
SUBROUTINE spgda_uv,5
#if defined (DOC)
*
C*s/r spgda_uv- Adjoint of the U,V inverse spectral transform
*
*
* Author : Luc Fillion - ARMA/EC - 25 Jun 2010.
*
* Revision:
*
#endif
*
use modstag
, only: level2_staggrid, nj_s, njlath_s
IMPLICIT NONE
*implicits
#include "comgrd_param.cdk"
#include "pardim.cdk"
#include "comdim.cdk"
#include "comsp.cdk"
#include "comgd0.cdk"
#include "comfftla.cdk"
#include "comgem.cdk"
#include "comcva.cdk"
C
INTEGER :: ILEN, JK, JLON, JLAT, INCPU, ISIZ, I_NJ
integer ji,jj,jla
C
integer idum1,idum2,idum3,idum4,injdim
integer inksdim
real*8 zmin,zmax
REAL*8 :: ZFM(NI+2,2*nflev,NJ)
real*8 zfmla(ni,2*nflev,nj), zsp(nla,2)
real*8 zgdxy(ni,nj)
C
C
C 1. Fourier transform all fields for all latitudes
C ----------------------------------------------
C
inksdim = 2*nflev
zfm(:,:,:) = 0.d0
if(level2_staggrid) then
i_nj = nj_s
else
i_nj = nj
end if
!
if (grd_typ.ne.'LU') then
DO JK = 1, inksdim
DO JLAT = 1, I_NJ
DO JLON = 1, NI
ZFM(JLON,JK,JLAT) = GD(JLON,JK,JLAT)
END DO
END DO
END DO
endif
C
injdim = i_nj
if(lcva_hsp) injdim = njlath
if (grd_typ.eq.'GU') then
CALL FFT3DVAR
(ZFM,1,NI+2,inksdim,1,nj,inksdim,-1)
else if (grd_typ.eq.'LU') then
do jk = 1, inksdim
do jj=1,nj
do ji=1,ni
zfmla(ji,jk,jj) = gd(ji,jk,jj)
zgdxy(ji,jj) = gd(ji,jk,jj)
enddo
enddo
if(lrpnfft) then
call aidft2dr
(zsp,zgdxy)
else
! call aidft2d(zsp,zfmla,sdft2d,nindxy,mlen2d, update...
! & jk,nla,ni,nj,inksdim)
endif
do jla = 1, nla
sp(jla,1,jk) = zsp(jla,1)
sp(jla,2,jk) = zsp(jla,2)
enddo
enddo
endif
C
C 2. Direct Legendre transform including wind transformations
C --------------------------------------------------------
C
if (grd_typ.eq.'GU') then
if (lcva_hsp) then
CALL SPGDAPAR_hem
(ZFM,NI,inksdim,njlath,level2_staggrid)
else
CALL SPGDAPAR
(ZFM,NI,inksdim,NJ,level2_staggrid)
endif
endif
C
END subroutine spgda_uv