!-------------------------------------- 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 SPPAR (lstaggrid) 1,5
#if defined (DOC)
*
C**s/r SPPAR  - Inverse spectral transform(PARALLEL LOOP)
*
*
*Author  : P.Gauthier/P.Koclas  - *ARMA and CMC  May 1997
*     .     Complete rewriting of the transform for a parallel treatment
*     .     based on having the outer loop on the zonal wavenumber
*     Revision:
*       JM Belanger CMDA/SMC  Aug 2000
*                   . 32 bits conversion
*       P. Koclas CMDA/SMC  Apr 2003
*                   . removed commulti and added openmp for ibm conversion
*     . P. Gauthier ARMA/MSC June 2003
*     .     Adaptation to produce winds on staggered grid in latitude
*     .     (see STAGWINDS for information on how to use this)
*
*     Purpose: builds the grid point model state from the spectral
*              state.A
*
*
** Method: the inverse transform is done by latitude bands and the
*          result is put in comgd0.cdk
*
*     Arguments:
*     LSTAGGRID: inverse legendre transform is done on the staggered grid for winds only
*
*
#endif
      use modstag, only: nj_s, njlath_s
      IMPLICIT NONE
*implicits
#include "pardim.cdk"
#include "comdim.cdk"
#include "comleg.cdk"
#include "comgem.cdk"
#include "comgd0.cdk"
#include "comsp.cdk"
                                !#include <rpnmacros_f.h>
*
*     Arguments
*
      logical :: lstaggrid
C
C     Local variables
C
      INTEGER ISIZ
      INTEGER ILEN, JJ, JLATN, JLATS
     S     ,JM, ILONR,ILONIM,JK
     S     ,JN, ILA, INM, JM0, JNS, INS
                                !
      integer :: i_njlath, i_nj
C
      REAL*8 ZJM
      REAL*8 DLALP(0:NTRUNCMX,NJLATH), DLDALP(0:NTRUNCMX,NJLATH)
      REAL*8 ZFMS(NJLATH+1,2,NKSDIM), ZFMA(NJLATH+1,2,NKSDIM)
      REAL*8 DLSP(0:NTRUNCMX,2,NKSDIM)

      INTEGER  thdid,numthd,omp_get_thread_num,omp_get_num_threads
C
                                !
                                !    1. Set up according to the desired grid (staggered or not)
                                !  ---------------------
      if(lstaggrid) then
         i_nj = nj_s
         i_njlath = njlath_s
      else
         i_nj = nj
         i_njlath = njlath
      end if
C
C     *    2. Inverse Legendre transform
C     --------------------------
 200  CONTINUE
*
!$OMP PARALLEL PRIVATE(DLALP,DLDALP,JM0,DLSP,ZFMS,ZFMA,JLATN,JLATS)
!$OMP+  PRIVATE(INM,ILA,JM,JN,JK,JJ,ZJM,JNS,INS,ILONR,ILONIM)
!$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)
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 = 2*NFLEV+1, NKSDIM
                  DLSP(INM,1,JK) = SP(ILA,1,JK)
                  DLSP(INM,2,JK) = SP(ILA,2,JK)
               END DO
               DO JK = 1, 2*NFLEV
                  DLSP(INM,1,JK) = SP(ILA,1,JK)*R1SNP1(ILA)
                  DLSP(INM,2,JK) = SP(ILA,2,JK)*R1SNP1(ILA)
               END DO
            END DO
*
C**   .   2.2  Get Legendre polynomial (and its derivative) for all latitudes
C**   .        but for the chosen value of "m" from the global array
*
            if (lstaggrid) then
               CALL GETALP_Stag(DLALP,DLDALP,I_NJLATH,NTRUNC,NTRUNCMX,JM)
            else
               CALL GETALP  (DLALP,DLDALP,I_NJLATH,NTRUNC,NTRUNCMX,JM)
            end if
*
**   .   2.3  Perform the inverse Legendre transform for all fields
*
            CALL LEGINV3 (JM, ZFMS,ZFMA,DLSP,DLALP
     S           ,NKSDIM, I_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, I_NJLATH
               JLATN = JJ
               JLATS = I_NJ - JLATN + 1
               DO  JK = 2*NFLEV+1, NKSDIM
                  GD(ILONR,JK,JLATN) = ZFMS(JJ,1,JK)
     S                 + ZFMA(JJ,1,JK)
                  GD(ILONIM,JK,JLATN) = ZFMS(JJ,2,JK)
     S                 + ZFMA(JJ,2,JK)
                  GD(ILONR,JK,JLATS) = ZFMS(JJ,1,JK)
     S                 - ZFMA(JJ,1,JK)
                  GD(ILONIM,JK,JLATS) = ZFMS(JJ,2,JK)
     S                 - ZFMA(JJ,2,JK)
               END DO
C
               DO JK = 1, NFLEV
                  UT0(ILONR,JK,JLATN)  = -ZJM*(ZFMS(JJ,2,JK+NFLEV)
     S                 + ZFMA(JJ,2,JK+NFLEV))
                  UT0(ILONIM,JK,JLATN) = ZJM*(ZFMS(JJ,1,JK+NFLEV)
     S                 + ZFMA(JJ,1,JK+NFLEV))
                  UT0(ILONR,JK,JLATS)  = -ZJM*(ZFMS(JJ,2,JK+NFLEV)
     S                 - ZFMA(JJ,2,JK+NFLEV))
                  UT0(ILONIM,JK,JLATS) = ZJM*(ZFMS(JJ,1,JK+NFLEV)
     S                 - ZFMA(JJ,1,JK+NFLEV))
*
                  VT0(ILONR,JK,JLATN)  = -ZJM*(ZFMS(JJ,2,JK)
     S                 + ZFMA(JJ,2,JK))
                  VT0(ILONIM,JK,JLATN) = ZJM*(ZFMS(JJ,1,JK)
     S                 + ZFMA(JJ,1,JK))
                  VT0(ILONR,JK,JLATS)  = -ZJM*(ZFMS(JJ,2,JK)
     S                 - ZFMA(JJ,2,JK))
                  VT0(ILONIM,JK,JLATS) = ZJM*(ZFMS(JJ,1,JK)
     S                 - ZFMA(JJ,1,JK))
               END DO
 240        CONTINUE
C
C**   .   2.5 Completion of the computation of the winds in Fourier space
C     .   -----------------------------------------------------------
C
            CALL LEGINV3 (JM, ZFMA,ZFMS,DLSP,DLDALP
     S           ,2*NFLEV, I_NJLATH, NTRUNC, NTRUNCMX)
C
            ILONR  = 2*JM + 1
            ILONIM = 2*JM + 2
            ZJM    = FLOAT(JM)
            DO JJ = 1, I_NJLATH
               JLATN = JJ
               JLATS = I_NJ - JLATN + 1
                                !
               if(jlatn.ne.jlats) then
                  do jk = 1,nflev
                                !
                                ! Northern latitudes
                                !
                     UT0(ILONR,JK,JLATN) = UT0(ILONR,JK,JLATN)
     S                    -(ZFMS(JJ,1,JK) + ZFMA(JJ,1,JK))
                     UT0(ILONIM,JK,JLATN) = UT0(ILONIM,JK,JLATN)
     S                    -(ZFMS(JJ,2,JK)+ ZFMA(JJ,2,JK))
C
                     VT0(ILONR,JK,JLATN) = VT0(ILONR,JK,JLATN)
     S                    +(ZFMS(JJ,1,JK+NFLEV) + ZFMA(JJ,1,JK+NFLEV))
                     VT0(ILONIM,JK,JLATN) = VT0(ILONIM,JK,JLATN)
     S                    +(ZFMS(JJ,2,JK+NFLEV)+ ZFMA(JJ,2,JK+NFLEV))
                                !
                                ! Southern latitudes
                                !
                     UT0(ILONR,JK,JLATS) = UT0(ILONR,JK,JLATS)
     S                    -(ZFMS(JJ,1,JK)- ZFMA(JJ,1,JK))
                     UT0(ILONIM,JK,JLATS) = UT0(ILONIM,JK,JLATS)
     S                    -(ZFMS(JJ,2,JK)- ZFMA(JJ,2,JK))
                     VT0(ILONR,JK,JLATS) = VT0(ILONR,JK,JLATS)
     S                    +(ZFMS(JJ,1,JK+NFLEV)- ZFMA(JJ,1,JK+NFLEV))
                     VT0(ILONIM,JK,JLATS) = VT0(ILONIM,JK,JLATS)
     S                    +(ZFMS(JJ,2,JK+NFLEV)- ZFMA(JJ,2,JK+NFLEV))
                  END DO
               else
                                !
                                ! Special case for the equator
                                !
                  DO JK = 1, NFLEV
                     UT0(ILONR,JK,JLATN) = UT0(ILONR,JK,JLATN)
     S                    -(ZFMS(JJ,1,JK) + ZFMA(JJ,1,JK))
                     UT0(ILONIM,JK,JLATN) = UT0(ILONIM,JK,JLATN)
     S                    -(ZFMS(JJ,2,JK)+ ZFMA(JJ,2,JK))
C
                     VT0(ILONR,JK,JLATN) = VT0(ILONR,JK,JLATN)
     S                    +(ZFMS(JJ,1,JK+NFLEV) + ZFMA(JJ,1,JK+NFLEV))
                     VT0(ILONIM,JK,JLATN) = VT0(ILONIM,JK,JLATN)
     S                    +(ZFMS(JJ,2,JK+NFLEV)+ ZFMA(JJ,2,JK+NFLEV))
                  end do
               end if
            END DO
C
 202     CONTINUE
C     End loop on m
 201  CONTINUE
!$OMP END PARALLEL
C
      END subroutine sppar