!-------------------------------------- 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 REESPEPAR(PSP,ZFM,KNI,KNJ,KLA,KDIM,KFIELD) 1,3
#if defined (DOC)
*
***s/r REESPEPAR -Spectral Transform for a number of scalar
*                 fields(PARALLEL LOOP)
*
*Author  : P. Gauthier/P.koclas *ARMA/AES CMC  Nov, 1997
*Revision:
*          JM Belanger CMDA/SMC  Jul 2000
*                   . 32 bits conversion
*          P. Koclas CMDA/SMC  Apr  2003
*           -remove commulti
*           -conversion to openmp
*          M. Buehner ARMA  May 2008
*           -use input argument KDIM, KFIELD instead of comdeck variables
*          L. Fillion ARMA/EC  11 May 2010 - Limit printout to processor 0.
*
*Arguments
*     Purpose:
*     i   KFIELD            : number of fields to be processed
*     o   PSP(KLA,2,KFIELD) :  array of spectral coefficients
*     i   ZFM(KNI+2,KFIELD, KNJ)
*     Following parameters are only necessary for the dimensioning
*     of PSP and PGD. Their equivalent global variables are used.
*     i   KLA               : total number of spectral coefficients
*     i   KDIM              : equivalent of NKSDIM
*     i   KNI               :
*     i   KNJ               :
*     .
#endif
*
      USE procs_topo
      IMPLICIT NONE
*implicits
#include "comdim.cdk"
#include "comlun.cdk"
#include "comleg.cdk"
#include <rpnmacros_f.h>
C
      INTEGER :: KFIELD
      INTEGER :: KLA,KDIM,KNI,KNJ
      REAL*8  :: PSP(KLA,2,KDIM)
      REAL*8  :: ZFM(KNI+2,KFIELD, KNJ)
C
      INTEGER :: ISIZ
      INTEGER :: ILEN, JJ, JK,JK2,ILONR, ILONI, JM,
     S     ILA, INM, JN, JM0, INS, JNS
      REAL*8 :: ZJM
      REAL*8 ::  DLALP(0:NTRUNCMX,NJLATH), DLDALP(0:NTRUNCMX, NJLATH)
      REAL*8 ::  DLSP(0:NTRUNCMX,2,KDIM)
      REAL*8 ::  ZFMS( NJLATH+1,2,KDIM), ZFMA( NJLATH+1,2,KDIM)
      real*8 :: dlrwt(njbeg:njend)

      INTEGER  thdid,numthd,omp_get_thread_num,omp_get_num_threads
C
C
C
                                !
                                !     1. Adjustment needed when an odd number of latitudes is considered
                                !     --------------------------------------------------------------------
      IF(myid == 0) THEN
        write(nulout,fmt='(/,4x,A)')'REESPEPAR- Direct Legendre transform'
      endif
!
      dlrwt(:)   = rwt  (:)
      if (mod(nj,2).ne.0) then
         dlrwt(njlath)   = dlrwt(njlath)/2.d0
      end if
*
*
!$OMP  PARALLEL PRIVATE(DLALP,DLDALP,JM0,DLSP,ZFMS,ZFMA,thdid,numthd)
!$OMP+ PRIVATE(INM,ILA,JM,JN,JK,JK2,JJ,ZJM,JNS,INS,ILONR,ILONI)
       thdid = omp_get_thread_num()
       numthd= omp_get_num_threads()
!      DO  201 JM0 = 0, NTRUNC/2
      DO  201 JM0 = thdid, NTRUNC/2,numthd
         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)
C
            DO 303 JJ = 1, NJLATH
C
C***  .       2.2.1  Coefficients for scalar fields
C
               DO JK = 1,KFIELD
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
 303        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,KFIELD, 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
            DO JK2 = 1,KFIELD
               PSP(ILA,1,JK2) = DLSP(INM,1,JK2)
               PSP(ILA,2,JK2) = DLSP(INM,2,JK2)
            END DO
         END DO
C
C     End of loop on zonal wavenumbers
C
 202     CONTINUE
 201  CONTINUE
!$OMP END PARALLEL
C
C     -----------------------
C
 300  CONTINUE
      RETURN
      END