!-------------------------------------- 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 SPGDAPAR(ZFM,KNI,KKSDIM,KNJ,ldstaggrid) 2,6
#if defined (DOC)
*
**s/r SPGDAPAR- Adjoint of the inverse spectral transform(PARALLEL LOOP)
*
*
*     Author  : P. Gauthier/P.Koclas *ARMA/AES and CMC
*     .                 - Nov 1997 Complete rewriting of the transform
*     .                 to have the most external loop to be on the zonal wavenumber
*     .                 and then to have a parallel treatment of the transform
*     Revision:
*       JM Belanger CMDA/SMC  Aug 2000
*                   . 32 bits conversion
*     . P. Koclas CMDA/SMC  Apr 2003
*     .             -removed commulti and added OpenMP fro IBM conversion
*     . P. Gauthier ARMA/MSC June 2003
*     .             Modification needed to perform the transform on a staggered grid in latitude
*     .             Correction of bug that would occur if an odd number of latitudes were treated
*     . L. Fillion ARMA/EC 11 May 2010 - Limit printout to processor 0.
*
*     Purpose: based on gdsppar to build the adjoint of the inverse
*     .        spectral transform (see Technical Notes)
Arguments:
*
*     ZFM()                 : Array of spectral coefficients
*     KNI                   : number of latitude circles on one hemisphere
*     KKSDIM                : truncation order (triangular)
*     KNJ                   : zonal wavenumber
*     LSTAGGRID             : logical key (.T. use the staggered grid in latitude)
*
*
#endif
*
      use modstag, only: nj_s, njlath_s, rwt_s
      USE procs_topo
      IMPLICIT NONE
*implicits
#include "pardim.cdk"
#include "comdim.cdk"
#include "comlun.cdk"
#include "comleg.cdk"
#include "comgem.cdk"
#include "comsp.cdk"
C
C     Arguments
C
      INTEGER :: KNI,KKSDIM,KNJ
      REAL*8  :: ZFM(KNI+2,KKSDIM, KNJ)
      logical :: ldstaggrid
C
C     Local variables
C
      INTEGER ISIZ
      INTEGER ILEN, JJ, JK, ILONR, ILONI, JM
     S     ,ILA, INM, JN, JM0, INS, JNS
      REAL*8 :: ZJM,DLRWT(NJBEG:NJEND)
      REAL*8 ::  DLALP(0:NTRUNCMX,NJLATH), DLDALP(0:NTRUNCMX, NJLATH)
      REAL*8 :: DLSP(0:NTRUNCMX,2,NKSDIM),DLSP2(0:NTRUNCMX,2,2*NFLEV)
      REAL*8 ::  ZFMS( NJLATH+1,2,NKSDIM), ZFMA( NJLATH+1,2,NKSDIM)
      integer :: i_njlath, i_nj

      INTEGER thdid,numthd,omp_get_thread_num,omp_get_num_threads
*     !
*     !    1. Set up according to the desired grid (staggered or not)
*     !  ---------------------
      IF(myid == 0) THEN
        write(nulout,fmt='(/,4x,A)')'SPGDAPAR- Adjoint of the inverse Legendre transform'
      endif
                                !
      if(ldstaggrid) then
         i_nj = nj_s
         i_njlath = njlath_s
         dlrwt(:) = rwt_s(:)
      else
         i_nj = nj
         i_njlath = njlath
         dlrwt(:) = rwt(:)
      end if
      if (mod(i_nj,2).ne.0) then
         dlrwt(i_njlath) = dlrwt(i_njlath)/2.d0
      end if
C
C**   2. Fourier transform all fields for all latitudes
C
!$OMP PARALLEL PRIVATE(DLALP,DLDALP,JM0,DLSP,DLSP2,ZFMS,ZFMA)
!$OMP+ PRIVATE(INM,ILA,JM,JN,JK,JJ,ZJM,JNS,INS,ILONR,ILONI)
!$OMP+ PRIVATE(thdid,numthd)
*
      thdid = omp_get_thread_num()
      numthd= omp_get_num_threads()
      DO  201 JM0 = thdid, NTRUNC/2,numthd
!      DO  201 JM0 = 0, NTRUNC/2
         INS=1
         IF(JM0.EQ.NTRUNC-JM0) INS = 0
         DO 202 JNS = 0,INS
            JM = (1-JNS)*JM0 + JNS*(NTRUNC - JM0)
            ILONR = 2 * JM + 1
            ILONI = ILONR + 1
            ZJM   = FLOAT(JM)
C
C**   .  2.1 Fetch the Legendre functions and their derivatives for this choice of "m"
C
            if(ldstaggrid) then
!     CRITICAL
               CALL GETALP_Stag(DLALP,DLDALP,I_NJLATH,NTRUNC,NTRUNCMX,JM)
!     END CRITICAL
            else
!     CRITICAL
               CALL GETALP(DLALP,DLDALP,I_NJLATH,NTRUNC,NTRUNCMX,JM)
!     END CRITICAL
            end if
C
C**   .  2.2  Build the symmetric and anti-symmetric Fourier coefficients including
C**   .       the appropriate quadrature weights (see scientific notes)
*
            DO 203 JJ = 1, I_NJLATH
C
C***  .       2.2.1  Coefficients for scalar fields
C
               DO JK = 2*NFLEV+1, NKSDIM
C     .                       SYMMETRIC COEFFICIENTS
                  ZFMS(JJ,1,JK) = DLRWT(JJ)*(ZFM(ILONR,JK,JJ)
     +                 + ZFM(ILONR,JK,1+I_NJ-JJ))
                  ZFMS(JJ,2,JK) = DLRWT(JJ)*(ZFM(ILONI,JK,JJ)
     +                 + ZFM(ILONI,JK,1+I_NJ-JJ))
C     .                        ANTISYMMETRIC COEFFICIENTS
                  ZFMA(JJ,1,JK) = DLRWT(JJ)*(ZFM(ILONR,JK,JJ)
     +                 - ZFM(ILONR,JK,1+I_NJ-JJ))
                  ZFMA(JJ,2,JK) = DLRWT(JJ)*(ZFM(ILONI,JK,JJ)
     +                 - ZFM(ILONI,JK,1+I_NJ-JJ))
               END DO
C
C***  .       2.2.2 Coefficients associated with the wind fields
C     .
               DO JK = 1, NFLEV
C     .                       VORTICITY: SYMMETRIC COEFFICIENTS
                  ZFMS(JJ,1,JK) = -ZJM*DLRWT(JJ)
     +                 *(ZFM(ILONI,JK+NFLEV,JJ)
     +                 +ZFM(ILONI,JK+NFLEV,1+I_NJ-JJ))
                  ZFMS(JJ,2,JK) = ZJM*DLRWT(JJ)
     +                 *(ZFM(ILONR,JK+NFLEV,JJ)
     +                 + ZFM(ILONR,JK+NFLEV,1+I_NJ-JJ))
C     .                       VORTICITY: ANTISYMMETRIC COEFFICIENTS
                  ZFMA(JJ,1,JK) = -ZJM*DLRWT(JJ)
     +                 *(ZFM(ILONI,JK+NFLEV,JJ)
     +                 - ZFM(ILONI,JK+NFLEV,1+I_NJ-JJ))
                  ZFMA(JJ,2,JK) = ZJM*DLRWT(JJ)
     +                 *(ZFM(ILONR,JK+NFLEV,JJ)
     +                 - ZFM(ILONR,JK+NFLEV,1+I_NJ-JJ))
C     .                        DIVERGENCE: SYMMETRIC COEFFICIENTS
                  ZFMS(JJ,1,JK+NFLEV) = -ZJM*DLRWT(JJ)
     +                 *(ZFM(ILONI,JK,JJ)+ ZFM(ILONI,JK,1+I_NJ-JJ))
                  ZFMS(JJ,2,JK+NFLEV) = ZJM*DLRWT(JJ)
     +                 *(ZFM(ILONR,JK,JJ)+ ZFM(ILONR,JK,1+I_NJ-JJ))
C     .                        DIVERGENCE: ANTISYMMETRIC COEFFICIENTS
                  ZFMA(JJ,1,JK+NFLEV) = -ZJM*DLRWT(JJ)
     +                 *(ZFM(ILONI,JK,JJ)- ZFM(ILONI,JK,1+I_NJ-JJ))
                  ZFMA(JJ,2,JK+NFLEV) = ZJM*DLRWT(JJ)
     +                 *(ZFM(ILONR,JK,JJ)- ZFM(ILONR,JK,1+I_NJ-JJ))
               END DO
 203        CONTINUE
C
C**   .    2.3 First one with ALP for all scalar fields and for half the terms
C**   .        required to define the divergence and vorticity
C
                                !
!      CRITICAL
            CALL LEGDIR3 (JM, ZFMS,ZFMA,DLSP,DLALP
     S           ,NKSDIM, I_NJLATH, NTRUNC, NTRUNCMX)
!     END CRITICAL
                                !
C
C**   .    2.4  Second transform with DALP to complete the construction of the
C**   .         vorticity and divergence fields
*
            DO 240 JJ = 1, I_NJLATH
               DO JK = 1, NFLEV
C     .                  Symmetric coefficients for zonal wind
                  ZFMS(JJ,1,JK) = DLRWT(JJ)*(ZFM(ILONR,JK,JJ)
     +                 + ZFM(ILONR,JK,1+I_NJ - JJ))
                  ZFMS(JJ,2,JK) = DLRWT(JJ)*(ZFM(ILONI,JK,JJ)
     +                 + ZFM(ILONI,JK,1+I_NJ-JJ))
C     .                   Antisymmetric coefficients for zonal wind
                  ZFMA(JJ,1,JK) = DLRWT(JJ)*(ZFM(ILONR,JK,JJ)
     +                 - ZFM(ILONR,JK,1+I_NJ-JJ))
                  ZFMA(JJ,2,JK) = DLRWT(JJ)*(ZFM(ILONI,JK,JJ)
     +                 - ZFM(ILONI,JK,1+I_NJ-JJ))
               END DO
               DO JK = NFLEV+1, 2*NFLEV
C     .                  Symmetric coefficients for meridional wind
                  ZFMS(JJ,1,JK) = -DLRWT(JJ)*(ZFM(ILONR,JK,JJ)
     +                 + ZFM(ILONR,JK,1+I_NJ-JJ))
                  ZFMS(JJ,2,JK) = -DLRWT(JJ)*(ZFM(ILONI,JK,JJ)
     +                 + ZFM(ILONI,JK,1+I_NJ-JJ))
C     .                   Antisymmetric coefficients for meridional wind
                  ZFMA(JJ,1,JK) = -DLRWT(JJ)*(ZFM(ILONR,JK,JJ)
     +                 - ZFM(ILONR,JK,1+I_NJ-JJ))
                  ZFMA(JJ,2,JK) = -DLRWT(JJ)*(ZFM(ILONI,JK,JJ)
     +                 - ZFM(ILONI,JK,1+I_NJ-JJ))
               END DO
 240        CONTINUE
!      CRITICAL
            CALL LEGDIR3 (JM, ZFMA,ZFMS,DLSP2,DLDALP
     S           ,2*NFLEV, i_NJLATH, NTRUNC, NTRUNCMX)
!     END CRITICAL
                                !
C
C**   .    2.5  Transfer the result in the global state
C
            DO JN = JM, NTRUNC
               ILA = NIND(JM) + JN - JM
               INM = JN - JM
C
               DO JK = 1, 2*NFLEV
                  SP(ILA,1,JK) = -R1SNP1(ILA)*(DLSP(INM,1,JK) + DLSP2(INM,1,JK))
                  SP(ILA,2,JK) = -R1SNP1(ILA)*(DLSP(INM,2,JK) + DLSP2(INM,2,JK))
               END DO
C
               DO JK = 2*NFLEV+1,NKSDIM
                  SP(ILA,1,JK) = DLSP(INM,1,JK)
                  SP(ILA,2,JK) = DLSP(INM,2,JK)
               END DO
            END DO
C
C     End of loop on zonal wavenumbers
C
                                !
 202     CONTINUE
 201  CONTINUE
!$OMP END PARALLEL
C
      END subroutine spgdapar