!-------------------------------------- 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 cv2gd_gl,4
#if defined (DOC)
*
***s/r cv2gd_gl - For grd_typ = 'GU, lcva_hemis = .true.: Convert analysis variables to model state variables
* .
*Author : Luc Fillion - 10 Jul 2009
* .
*Revision:
#endif
C
IMPLICIT NONE
#include "pardim.cdk"
#include "comdim.cdk"
#include "comlun.cdk"
#include "comcst.cdk"
#include "comgem.cdk"
#include "comleg.cdk"
#include "comcva.cdk"
#include "comsp.cdk"
#include "comspg.cdk"
#include "comgd0.cdk"
#include "compstat.cdk"
#include "comcorr.cdk"
#include "comcse1.cdk"
#include "comstate.cdk"
INTEGER ILEN
INTEGER JN,JM,JK,ILA,jk1,jk2
! REAL*8 ZSP(NKSDIM,2,0:NTRUNC),ZSP2(NKSDIM,2,0:NTRUNC),SQ2
! POINTER(PXZSP,ZSP), (PXZSP2,ZSP2)
REAL*8 SQ2
REAL*8 , ALLOCATABLE,DIMENSION(:,:,:) :: ZSP, ZSP2
INTEGER ILENSP, ILENGD, IERR
S , ILON, JLEV, JLON, JLAT, JLA, IULOUT
REAL*8 Z1MNU2, ZSQRTNU2, ZCORIOLIS
REAL*8 ZGDPSI(NIBEG:NIEND,NFLEV,NJBEG:NJEND)
S ,ZGDCHI(NIBEG:NIEND,NFLEV,NJBEG:NJEND)
REAL*8 DLZSOMME, DLDSOMME, DLA2, DL1SA2, DLFACT
S ,DLNORMPSI, DLNORMCHI
!
INTEGER :: JN0,INS,JNS,NTRUNCHF
INTEGER :: thdid,numthd,omp_get_thread_num,omp_get_num_threads
real*8 two
data two /2.0D0/
*-----------------------------------------------------------------
!
!*1. Build gridpoint PSI,CHI
! -----------------------
!
DLA2 = DBLE(RA)*DBLE(RA)
DL1SA2 = 1.D0/DLA2
!
CALL SPEREE
(NKSDIM,SP,GD
S ,NLA,NIBEG,NIEND,NJBEG,NJEND,NKSDIM)
!
DO JLAT = 1, NJ
ILON = NILON(JLAT)
DO JLEV = 1, NFLEV
DO JLON = 1, ILON
ZGDPSI(JLON,JLEV,JLAT)= UT0(JLON,JLEV,JLAT)
ZGDCHI(JLON,JLEV,JLAT)= VT0(JLON,JLEV,JLAT)
END DO
END DO
END DO
!
!*2. Denormalize
! -----------
!
IF(LUSE3DSTD) THEN
write(nulout,*) '*****USING STD3D*****'
DO JLAT = 1, NJ
DO JLEV = 1, NFLEV
ILON = NILON(JLAT)
DO JLON = 1, ILON
ZGDPSI(JLON,JLEV,JLAT)=
& ZGDPSI(JLON,JLEV,JLAT)*RGSIGUU3D(JLON,JLEV,JLAT)
ZGDCHI(JLON,JLEV,JLAT)=
& ZGDCHI(JLON,JLEV,JLAT)*RGSIGVV3D(JLON,JLEV,JLAT)
TT0(JLON,JLEV,JLAT)=
& TT0(JLON,JLEV,JLAT)*RGSIGTT3D(JLON,JLEV,JLAT)
Q0(JLON,JLEV,JLAT)=
& Q0(JLON,JLEV,JLAT)*RGSIGQ3D(JLON,JLEV,JLAT)
END DO
END DO
END DO
DO JLAT = 1, NJ
ILON = NILON(JLAT)
DO JLON = 1, ILON
GPS0(JLON,1,JLAT)=
& GPS0(JLON,1,JLAT)*RGSIGPS3D(JLON,1,JLAT)
END DO
END DO
IF(nsexist(nstg).eq.1) THEN
DO JLAT = 1, NJ
ILON = NILON(JLAT)
DO JLON = 1, ILON
GTG0(JLON,1,JLAT)=
& GTG0(JLON,1,JLAT)*RGSIGTG3D(JLON,1,JLAT)
END DO
END DO
ENDIF
ELSE
CALL FGERR
('M')
cpik
!$OMP PARALLEL DO PRIVATE(DLNORMPSI,DLNORMCHI,ILON,JLEV,JLAT,JLON)
DO JLEV = 1, NFLEV
DO JLAT = 1, NJ
DLNORMPSI = RGSIGUU(JLAT,JLEV)
DLNORMCHI = RGSIGVV(JLAT,JLEV)
ILON = NILON(JLAT)
DO JLON = 1, ILON
ZGDPSI(JLON,JLEV,JLAT)
S = ZGDPSI(JLON,JLEV,JLAT)*DLNORMPSI
ZGDCHI(JLON,JLEV,JLAT)
S = ZGDCHI(JLON,JLEV,JLAT)*DLNORMCHI
END DO
END DO
END DO
!$OMP END PARALLEL DO
ENDIF
!
!*3. Apply INMI balancing of increment
! ---------------------------------
!
if(linmi) then
call lbalgl
(zgdpsi,zgdchi)
endif
!
!*4. Spectral transform all fields
! -----------------------------
!
DO JLAT = 1, NJ
ILON = NILON(JLAT)
DO JLEV = 1, NFLEV
DO JLON = 1, ILON
UT0(JLON,JLEV,JLAT) = ZGDPSI(JLON,JLEV,JLAT)
VT0(JLON,JLEV,JLAT) = ZGDCHI(JLON,JLEV,JLAT)
END DO
END DO
END DO
!
!*5. Apply 3D amplification factor after constructing full variables
! ---------------------------------------------------------------
!
DO JLAT = 1, NJ
DO JLEV = 1, NKGDIM
DO JLON = 1, NI
GD(JLON,JLEV,JLAT)=
& GD(JLON,JLEV,JLAT)*damplibg(JLON,JLEV,JLAT)
END DO
END DO
END DO
!
CALL REESPE
(NKSDIM,SP,GD
& ,NLA,NIBEG,NIEND,NJBEG,NJEND,NKSDIM)
!
!*6. Rederive the vorticity and divergence from PSI and CHI
! ------------------------------------------------------
!
!$OMP PARALLEL DO
DO JLEV = 1, NFLEV
DO JLA = 1, NLA
SPVOR(JLA,1,JLEV) = SPVOR(JLA,1,JLEV)*DL1SA2*RNNP1(JLA)
SPVOR(JLA,2,JLEV) = SPVOR(JLA,2,JLEV)*DL1SA2*RNNP1(JLA)
SPDIV(JLA,1,JLEV) = SPDIV(JLA,1,JLEV)*DL1SA2*RNNP1(JLA)
SPDIV(JLA,2,JLEV) = SPDIV(JLA,2,JLEV)*DL1SA2*RNNP1(JLA)
END DO
END DO
!$OMP END PARALLEL DO
!
RETURN
END