!-------------------------------------- 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 GDSPPAR(ZFM,KNI,KKSDIM,KNJ) 1,3
#if defined (DOC)
*
***s/r GDSPPAR  - Spectral Transform(PARALLEL LOOP)
*
*     Purpose: spectral transform of the model state including conversion of
*     .        wind images to vorticity and divergence
*
*     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:
*     P. Koclas CMC Apr 2003
*     .         - removed commulti and added openmp for ibm conversion
*     P. Gauthier ARMA/MSC July 2003
*     .         - modifications to handle an odd number of latitudes circles
*Arguments:
*     ZFM()                 : Array of spectral coefficients
*     KNI                   : number of latitude circles on one hemisphere
*     KKSDIM                : truncation order (triangular)
*     KNJ                   : zonal wavenumber
*
#endif
*
      IMPLICIT NONE
*implicits
#include "comdim.cdk"
#include "comlun.cdk"
#include "comleg.cdk"
#include "comsp.cdk"
#include <rpnmacros_f.h>
C     
      INTEGER :: KNI,KKSDIM,KNJ
      REAL*8  :: ZFM(KNI+2,KKSDIM, KNJ)
C
      INTEGER :: ISIZ
      INTEGER :: ILEN, JLA, JJ, JK, ILONR, ILONI, JM
     S     ,ILA, INM, JN, JM0, INS, JNS
      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)
      REAL*8  :: ZJM, DLRWT(NJBEG:NJEND), dlrwocs(NJBEG:NJEND)
C     
                                !     
                                !     1. Adjustment needed when an odd number of latitudes is considered
                                !     --------------------------------------------------------------------
      write(nulout,fmt='(/,4x,A)')'GDSPPAR- Direct Legendre transform'
      dlrwt(:)   = rwt  (:)
      dlrwocs(:) = rwocs(:)
      if (mod(nj,2).ne.0) then
         dlrwt(njlath)   = dlrwt(njlath)/2.d0
         dlrwocs(njlath) = dlrwocs(njlath)/2.d0
      end if
*
!$OMP PARALLEL PRIVATE(DLALP,DLDALP,JM0,DLSP,DLSP2,ZFMS,ZFMA)
!$OMP DO PRIVATE(INM,ILA,JM,JN,JK,JJ,ZJM,JNS,INS,ILONR,ILONI)
      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     
            CALL GETALP(DLALP,DLDALP,NJLATH,NTRUNC,NTRUNCMX,JM)
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, 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+NJ-JJ))
                  ZFMS(JJ,2,JK) = DLRWT(JJ)*(ZFM(ILONI,JK,JJ)
     +                 + ZFM(ILONI,JK,1+NJ-JJ))
C     .                        ANTISYMMETRIC COEFFICIENTS
                  ZFMA(JJ,1,JK) = DLRWT(JJ)*(ZFM(ILONR,JK,JJ)
     +                 - ZFM(ILONR,JK,1+NJ-JJ))
                  ZFMA(JJ,2,JK) = DLRWT(JJ)*(ZFM(ILONI,JK,JJ)
     +                 - ZFM(ILONI,JK,1+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*DLRWOCS(JJ)
     +                 *(ZFM(ILONI,JK+NFLEV,JJ)
     +                 +ZFM(ILONI,JK+NFLEV,1+NJ-JJ))
                  ZFMS(JJ,2,JK) = ZJM*DLRWOCS(JJ)
     +                 *(ZFM(ILONR,JK+NFLEV,JJ)
     +                 + ZFM(ILONR,JK+NFLEV,1+NJ-JJ))
C     .                       VORTICITY: ANTISYMMETRIC COEFFICIENTS
                  ZFMA(JJ,1,JK) = -ZJM*DLRWOCS(JJ)
     +                 *(ZFM(ILONI,JK+NFLEV,JJ)
     +                 - ZFM(ILONI,JK+NFLEV,1+NJ-JJ))
                  ZFMA(JJ,2,JK) = ZJM*DLRWOCS(JJ)
     +                 *(ZFM(ILONR,JK+NFLEV,JJ)
     +                 - ZFM(ILONR,JK+NFLEV,1+NJ-JJ))
C     .                        DIVERGENCE: SYMMETRIC COEFFICIENTS
                  ZFMS(JJ,1,JK+NFLEV) = -ZJM*DLRWOCS(JJ)
     +                 *(ZFM(ILONI,JK,JJ)+ ZFM(ILONI,JK,1+NJ-JJ))
                  ZFMS(JJ,2,JK+NFLEV) = ZJM*DLRWOCS(JJ)
     +                 *(ZFM(ILONR,JK,JJ)+ ZFM(ILONR,JK,1+NJ-JJ))
C     .                        DIVERGENCE: ANTISYMMETRIC COEFFICIENTS
                  ZFMA(JJ,1,JK+NFLEV) = -ZJM*DLRWOCS(JJ)
     +                 *(ZFM(ILONI,JK,JJ)- ZFM(ILONI,JK,1+NJ-JJ))
                  ZFMA(JJ,2,JK+NFLEV) = ZJM*DLRWOCS(JJ)
     +                 *(ZFM(ILONR,JK,JJ)- ZFM(ILONR,JK,1+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     
            CALL LEGDIR3 (JM, ZFMS,ZFMA,DLSP
     S           ,DLALP,NKSDIM, NJLATH, NTRUNC, NTRUNCMX)

C     
C**   .    2.4  Second transform with DALP to complete the construction of the
C**   .         vorticity and divergence fields
*     
            DO 241 JJ = 1, NJLATH
               DO JK = 1, NFLEV
C     .                  Symmetric coefficients for zonal wind
                  ZFMS(JJ,1,JK) = DLRWOCS(JJ)*(ZFM(ILONR,JK,JJ)
     +                 + ZFM(ILONR,JK,1+NJ - JJ))
                  ZFMS(JJ,2,JK) = DLRWOCS(JJ)*(ZFM(ILONI,JK,JJ)
     +                 + ZFM(ILONI,JK,1+NJ-JJ))
C     .                   Antisymmetric coefficients for zonal wind
                  ZFMA(JJ,1,JK) = DLRWOCS(JJ)*(ZFM(ILONR,JK,JJ)
     +                 - ZFM(ILONR,JK,1+NJ-JJ))
                  ZFMA(JJ,2,JK) = DLRWOCS(JJ)*(ZFM(ILONI,JK,JJ)
     +                 - ZFM(ILONI,JK,1+NJ-JJ))
               END DO
               DO JK = NFLEV+1, 2*NFLEV
C     .                  Symmetric coefficients for zonal wind
                  ZFMS(JJ,1,JK) = -DLRWOCS(JJ)*(ZFM(ILONR,JK,JJ)
     +                 + ZFM(ILONR,JK,1+NJ-JJ))
                  ZFMS(JJ,2,JK) = -DLRWOCS(JJ)*(ZFM(ILONI,JK,JJ)
     +                 + ZFM(ILONI,JK,1+NJ-JJ))
C     .                   Antisymmetric coefficients for zonal wind
                  ZFMA(JJ,1,JK) = -DLRWOCS(JJ)*(ZFM(ILONR,JK,JJ)
     +                 - ZFM(ILONR,JK,1+NJ-JJ))
                  ZFMA(JJ,2,JK) = -DLRWOCS(JJ)*(ZFM(ILONI,JK,JJ)
     +                 - ZFM(ILONI,JK,1+NJ-JJ))
               END DO
 241     CONTINUE
         CALL LEGDIR3 (JM, ZFMA,ZFMS,DLSP2,
     S         DLDALP,2*NFLEV, NJLATH, NTRUNC, NTRUNCMX)
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) = DLSP(INM,1,JK) + DLSP2(INM,1,JK)
               SP(ILA,2,JK) = 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 DO
!$OMP END PARALLEL
C
      RETURN
      END