!-------------------------------------- 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 spereepar_hem(KFIELD,PSP,PGD, 1,2
     S      KLA, KIBEG, KIEND, KJBEG, KJEND, KDIM)
#if defined (DOC)
*
**s/r spereepar_hem -Hemispheric version of spereepar.ftn
*                fields(PARALLEL LOOP)
*
*Author  : L. Fillion - ARMA/EC - 13 May 2010.
*Revision:
*Arguments
*     i   KFIELD            : number of fields to be processed
*     i   PSP(KLA,2,KDIM)   :  array of spectral coefficients
*     o   PGD(KIBEG:KIEND,KDIM,KJBEG:KJEND)
*     .                     : grid-point field
*     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   KIBEG, KIEND      : equivalent of NIBEG and NIEND
*     i   KJBEG, KJEND      : equivalent of NJBEG and NJEND
*
#endif
      IMPLICIT NONE
*implicits
#include "comdim.cdk"
#include "comleg.cdk"
#include <rpnmacros_f.h>
C
      INTEGER KFIELD
      INTEGER KLA, KIBEG, KIEND, KJBEG, KJEND, KDIM
      REAL*8 PSP(KLA,2,KDIM), PGD(KIBEG:KIEND,KDIM,KJBEG:KJEND)
C
      INTEGER ISIZ
      INTEGER ILEN, JJ, JLATN, JLATS
     S     , JM, ILONR,ILONIM,JK
     S     , JN, ILA, INM, JM0, JNS, INS
C
      REAL*8 ZJM
      REAL*8 DLALP(0:NTRUNCMX,NJLATH), DLDALP(0:NTRUNCMX,NJLATH)
      REAL*8 ZFMS(NJLATH+1,2,KFIELD), ZFMA(NJLATH+1,2,KFIELD)
      REAL*8 DLSP(0:NTRUNCMX,2,KFIELD)

      INTEGER  thdid,numthd,omp_get_thread_num,omp_get_num_threads
C
C
C
C        --------------------------
 200     CONTINUE
C
!$OMP  PARALLEL PRIVATE(DLALP,DLDALP,DLSP,ZFMS,ZFMA,JLATN,JLATS,
!$OMP+  INM,ILA,JM,JN,JK,JJ,ZJM,JNS,INS,ILONR,ILONIM,JM0)
         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)
C
C**   .   2.1 Copy global spectral state into local spectral state
C
               DO JN = JM, NTRUNC
                  ILA = NIND(JM) + JN - JM
                  INM = JN - JM
                  DO JK =1,KFIELD
                     DLSP(INM,1,JK) = PSP(ILA,1,JK)
                     DLSP(INM,2,JK) = PSP(ILA,2,JK)
                  END DO
               END DO
C
C     .   2.2  Get Legendre polynomial (and its derivative) for all latitudes
C     .        but for the chosen value of "m" from the global array
C
               CALL GETALP(DLALP,DLDALP,NJLATH,NTRUNC,NTRUNCMX,JM)
C
C     .   2.3  Perform the inverse Legendre transform for all fields
C
               CALL LEGINV3_hem (JM, ZFMS,ZFMA,DLSP,DLALP
     S              ,KFIELD, NJLATH, NTRUNC, NTRUNCMX)
*
C     .   2.4 Passage to Fourier space
*
               ILONR  = 2*JM + 1
               ILONIM = 2*JM + 2
               ZJM = FLOAT(JM)
               DO 240 JJ = 1, NJLATH
                  JLATN = JJ
                  JLATS = NJ - JLATN + 1
                  DO  JK = 1,KFIELD
                     PGD(ILONR,JK,JLATN) = ZFMS(JJ,1,JK)
                     PGD(ILONIM,JK,JLATN) = ZFMS(JJ,2,JK)
                     PGD(ILONR,JK,JLATS) = ZFMS(JJ,1,JK)
                     PGD(ILONIM,JK,JLATS) = ZFMS(JJ,2,JK)
                  END DO
C
 240           CONTINUE
C
C     End of loop on zonal wavenumbers
C
 202     CONTINUE
 201  CONTINUE
!$OMP END PARALLEL
C
C*    3. DEALLOCATE LOCAL ARRAYS
C     -----------------------
C
 300  CONTINUE
C
      RETURN
      END