!-------------------------------------- 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 SPA2GDAD 2,11
      use modstag, only : lstagwinds
#if defined (DOC)
*
***s/r SPA2GDAD  - Adjoint of the conversion of analysis variables to
*     .            model state increments
*               ***Version of SPA2SPAD adapted only for NANALVAR=4
*               ***for new PtoT approach with localized Tb correlations
*               ***Created May 2008 by M. Buehner
*     .
* Purpose
*     .  To revert from the analysis variables to the model state variables
*     .  as defined in spectral space
*     .  IMPORTANT: after a call to SP2SPAD, the fields SPVOR and SPDIV
*     .             may contain PSI and CHI. This will eventually become
*     .             the default option
*Author  : P. Gauthier *ARMA/AES  March 25, 1996
*     .
*Revision:
*     .    P. Gauthier *ARMA/AES  August 5, 1996
*     .         Redefine the balance relationship and reformulate the
*     .         analysis variable in terms of PSI and CHI instead of
*     .         vorticity and divergence.
*     . S. Pellerin *ARMA/AES Sept 97.
*                   - Control of the different model state of the 3Dvar
*                     through COMSTATE, COMSTATEC and COMSTNUM common
*                     blocks variables (comstate.cdk).
*     . M. Buehner November 97
*     .         Introduce transformation for new definition of control
*     .         vector that reduces background covariances to the identity
*     .         matrix (for NANALVAR=3)
*     . L. Fillion/M. Buehner 20 jul 98
*     .         Model's coordinate
*       C. Charette *ARMA/AES Nov 1998
*               Changed LDIVBAL to LBALDIV
*       JM Belanger CMDA/SMC  Aug 2000
*                   . 32 bits conversion
*                    (MXMA8, replace 2.0 by REAL*8 in SQRT)
*       P. Koclas    CMDA/SMC  Apr 2003
*                    cnnert call to mxma8 to mxmaop
*       M. Buehner *ARMA/SMC October 2004
*          - changed damplibg to 3D and to be after balance operators
#endif
C
      IMPLICIT NONE
#include "pardim.cdk"
#include "comdim.cdk"
#include "comcst.cdk"
#include "comgem.cdk"
#include "comleg.cdk"
#include "comcva.cdk"
#include "comsp.cdk"
#include "comspg.cdk"
#include "comgd0.cdk"
#include "comgd1.cdk"
#include "compstat.cdk"
#include "comcorr.cdk"
#include "comcse1.cdk"
#include "comstate.cdk"
#include "comlun.cdk"
      INTEGER ILEN
      INTEGER JN,JM,JK,ILA,JLATBIN,JLATMIN,JLATMAX
!      REAL*8 ZSP(NKSDIM2,2,0:NTRUNC),ZSP2(NKSDIM2,2,0:NTRUNC),SQ2
!      POINTER(PXZSP,ZSP), (PXZSP2,ZSP2)
      REAL*8 SQ2
      REAL*8 ,ALLOCATABLE,DIMENSION(:,:,:) :: zsp,zsp2
      REAL*8 zp(NI,nflev,NJ)
      REAL*8 ZWINDOW(NJ)

      INTEGER ILENGD, IERR, ILON, JLEV, JLON, JLAT, JLA, KLATPTOT(NLATBIN)
      REAL*8 DL1SA2,DLA2,ZNORMPSI, ZNORMCHI, ZCORIOLIS, ZPSB(NI,NJ)
!      REAL*8 ZGDPSI(NIBEG:NIEND,NFLEV,NJBEG:NJEND)
!     S     ,ZGDCHI(NIBEG:NIEND,NFLEV,NJBEG:NJEND)
!      POINTER (PXGDPSI,ZGDPSI),(PXGDCHI,ZGDCHI)
      REAL*8,ALLOCATABLE,DIMENSION(:,:,:) :: ZGDPSI ,ZGDCHI
      REAL*8 ZGD(NI,NKGDIM,NJ)
!
      INTEGER :: NTRUNCHF,JN0,INS,JNS,thdid,numthd
      INTEGER :: omp_get_thread_num, omp_get_num_threads

      real*8 two
      REAL*8 ONE8,ZERO8
      data ONE8,ZERO8/1.D0,0.D0/
*
      two = 2.d0
      DLA2   = DBLE(RA)*DBLE(RA)
      DL1SA2 = 1.D0/DLA2
*-------------------------------------------------------------------
C
C     1. Create local arrays for CHI and PSI if needed
C
 100  CONTINUE
C
C     .  1.1 Case where the analysis variable is X - Xb
C     .      (SPA2SP = Identity)
C
 110  CONTINUE
      IF(NANALVAR.EQ.0) THEN
         RETURN
      END IF
C
C     .  1.2 Memory allocation
C
      IF(NANALVAR.GE.2) THEN
         ALLOCATE(ZGDPSI(NIBEG:NIEND,NFLEV,NJBEG:NJEND))
         ALLOCATE(ZGDCHI(NIBEG:NIEND,NFLEV,NJBEG:NJEND))
      END IF
C
      if(nlatbin.eq.1) then
        klatptot(1)=1
      elseif(nlatbin.eq.3) then
        klatptot(1)=1
        klatptot(2)=nj/2
        klatptot(3)=nj
      endif
C
C     3. Formulation using PSI and CHI
C
      DO JLAT = 1,NJ
         ILON = NILON(JLAT)
         DO JLEV = 1, NKSDIM
            DO JLON = 1, ILON
               ZGD(JLON,JLEV,JLAT)= GD(JLON,JLEV,JLAT)
            END DO
         END DO
      END DO
C
      DO JLATBIN=1,NLATBIN
C
C COMBINE LATITUDE BINS TO OBTAIN SINGLE GLOBAL 3D STATE INCREMENT
      ZWINDOW(:)=0.0d0
      CALL SUWINLATBIN(ZWINDOW,JLATBIN)
      call transfer('ZGD0')
      DO JLAT = 1,NJ
         ILON = NILON(JLAT)
         DO JLEV = 1, NKGDIM
            DO JLON = 1, ILON
               GD(JLON,JLEV,JLAT)= ZWINDOW(JLAT)*ZGD(JLON,JLEV,JLAT)
            END DO
         END DO
      END DO
C
      if (lstagwinds) then
         call stagwinds_ad(nulout)
      else
         CALL SPGDA
      end if
C
 300  CONTINUE
C
C     .  3.5 Rederive the vorticity and divergence from PSI and CHI
C
!$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
C
C     .  3.4 Spectral transform all fields
C
 340  CONTINUE
      CALL SPEREE(NKSDIM,SP,GD
     S     ,NLA,NIBEG,NIEND,NJBEG,NJEND,NKSDIM)
C
C Apply 3D amplification factor
C
        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
c
      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
C
C     .  3.2 Multiply by the background error variances
C
      DO JLAT = 1, NJ
        DO JLEV = 1, NFLEV
          ILON = NILON(JLAT)
          DO JLON = 1, ILON
            TB0(JLON,JLEV,JLAT)=TT0(JLON,JLEV,JLAT)
c            TB0(JLON,JLEV,JLAT)=0.0d0
            TB0(JLON,JLEV,JLAT)=
     +        TB0(JLON,JLEV,JLAT)*RGSIGTB(JLAT,JLEV)
          END DO
        END DO
      END DO
      DO JLAT = 1, NJ
        ILON = NILON(JLAT)
        DO JLON = 1, ILON
          ZPSB(JLON,JLAT)=GPS0(JLON,1,JLAT)
c          ZPSB(JLON,JLAT)=0.0
          ZPSB(JLON,JLAT)=
     +      ZPSB(JLON,JLAT)*RGSIGPSB(JLAT)
        END DO
      END DO
      CALL FGERR('M')

      IF(LBALDIV) THEN
        CALL ADIVBAL(zgdpsi,zgdchi)
      ENDIF

!$OMP PARALLEL DO PRIVATE(ILON,ZNORMPSI,ZNORMCHI)
      DO JLAT = 1, NJ
        ILON = NILON(JLAT)
        DO JLEV = 1, NFLEV
            ZNORMPSI = RGSIGUU(JLAT,JLEV)
            ZNORMCHI = RGSIGVV(JLAT,JLEV)
            DO JLON = 1, ILON
              ZGDPSI(JLON,JLEV,JLAT)
     S             = ZGDPSI(JLON,JLEV,JLAT)*ZNORMPSI
              ZGDCHI(JLON,JLEV,JLAT)
     S             = ZGDCHI(JLON,JLEV,JLAT)*ZNORMCHI
            END DO
        END DO
      END DO
!$OMP END PARALLEL DO

C
C     .  3.3 Obtain the full geopotential
C
      DO JLAT = 1, NJ
         ZCORIOLIS = 2.*ROMEGA*RMU(JLAT)
         ILON = NILON(JLAT)
         DO JLEV = 1, NFLEV
            DO JLON = 1, ILON
               TB0(JLON,JLEV,JLAT)= ZCORIOLIS*TB0(JLON,JLEV,JLAT)
            END DO
         END DO
      END DO

      DO JLAT = 1, NJ
        ILON = NILON(JLAT)
        DO JLON = 1, ILON
          do jlev=1,NFLEVPTOT
            zp(jlon,jlev,jlat) = PtoT(NFLEV+1,jlev,klatptot(jlatbin))*zpsb(jlon,jlat)
          END DO
        END DO
      ENDDO
      DO JLAT = 1, NJ
        ZCORIOLIS = 2.*ROMEGA*RMU(JLAT)
        DO JLEV = 1, NFLEVPTOT
          DO JLON = 1, ILON
            ZGDPSI(JLON,JLEV,JLAT) = ZCORIOLIS*zp(JLON,JLEV,JLAT)
     S             + ZGDPSI(JLON,JLEV,JLAT)
          END DO
        END DO
      END DO
c
      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
C
C
C     .  3.1 Transform PSI and CHI and the other variables
C     .      to physical space
C
 310  CONTINUE
C
      CALL REESPE(NKSDIM,SP,GD
     S     ,NLA,NIBEG,NIEND,NJBEG,NJEND,NKSDIM)
      CALL REESPE(NFLEV,SPTB,TB0
     S     ,NLA,NIBEG,NIEND,NJBEG,NJEND,NFLEV)
C
C        3.0 Multiply by CORNS to add in vert corr (scale/rotate)
C
       NTRUNCHF=NTRUNC/2
       SQ2=sqrt(two)
!$OMP PARALLEL PRIVATE(thdid,numthd,JN,JM,JK,ILA,JN0,INS,JNS,zsp,zsp2)
       thdid=omp_get_thread_num()
       numthd=omp_get_num_threads()

       ALLOCATE(zsp(NKSDIM,2,0:NTRUNC))
       ALLOCATE(zsp2(NKSDIM2,2,0:NTRUNC))

       DO JN0 = thdid, NTRUNCHF,numthd
         INS=1
         IF(JN0 == (NTRUNC-JN0))INS=0
         DO 222 JNS=0,INS
            JN=(1-JNS)*JN0+JNS*(NTRUNC-JN0)
C
C*    Transfer Dx from global to working array
C
         DO JM = 0, JN
            ILA = NIND(JM) + JN - JM
            DO JK = 1, NKSDIM
               ZSP2(JK,1,JM) = SP(ILA,1,JK)
               ZSP2(JK,2,JM) = SP(ILA,2,JK)
            END DO
         END DO
         DO JM = 0, JN
            ILA = NIND(JM) + JN - JM
            DO JK = 1, NFLEV
               ZSP2(JK+NKSDIM,1,JM) = SPTB(ILA,1,JK)
               ZSP2(JK+NKSDIM,2,JM) = SPTB(ILA,2,JK)
            END DO
         END DO
C
C*    Compute (CORNS(NKSDIM,NKSDIM,JN)' * ZSP)
C
c            CALL MXMAOP(CORNS(1,1,JN),NKSDIM2,1,ZSP2(1,1,0),1,NKSDIM2
c     +        ,ZSP(1,1,0),1,NKSDIM,NKSDIM,NKSDIM2,2*(JN+1))
            CALL DGEMM('T','N',NKSDIM,2*(JN+1),NKSDIM2,ONE8,CORNS(1,1,JN,JLATBIN),NKSDIM2,
     +                 ZSP2(1,1,0),NKSDIM2,ZERO8,ZSP(1,1,0),NKSDIM)
C
C*    Transfer from working array to global array
C
         ILA = NIND(0) +JN
         DO JK = 1, NKSDIM
            SPLAT(ILA,1,JK,JLATBIN) = ZSP(JK,1,0)*SQ2
            SPLAT(ILA,2,JK,JLATBIN) = ZSP(JK,2,0)*SQ2
         END DO
         ILA = NIND(0) +JN
C
         DO JM = 1, JN
            ILA = NIND(JM) +JN - JM
            DO JK = 1, NKSDIM
               SPLAT(ILA,1,JK,JLATBIN) = ZSP(JK,1,JM)
               SPLAT(ILA,2,JK,JLATBIN) = ZSP(JK,2,JM)
            END DO
         END DO
C
222    CONTINUE
       END DO
       DEALLOCATE(zsp)
       DEALLOCATE(zsp2)
!$OMP END PARALLEL

C
         DO JK=1,NKSDIM
            DO ILA=1,NTRUNC+1
               SPLAT(ILA,2,JK,JLATBIN)=0.0
            ENDDO
         ENDDO
C
C END LOOP ON JLATBIN
      ENDDO

C
C     9. Deallocate local arrays
C
 900  CONTINUE
C
C
      IF(NANALVAR.GE.2) THEN
         DEALLOCATE(ZGDPSI)
         DEALLOCATE(ZGDCHI)
      END IF
C
      write(nulout,*) 'end of spa2gdad!!!'
      call vflush(nulout)
      RETURN
      END