!-------------------------------------- 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 GDSPAPAR 1,3
#if defined (DOC)
*
C**s/r GDSPAPAR - Adjoint of Spectral transform( PARALLEL LOOP)
*
*
*Author : P.Gauthier/P.Koclas - *ARMA and CMC November 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 Jul 2000
* . 32 bits conversion
*
* P. Koclas CMDA/SMC Apr 2003
* -removed commulti
* -conversion to openmp
* P. Gauthier ARMA#MSC July 2003
* . Correction to handle an odd number of latitude circles
*
* Purpose: builds the grid point model state from the spectral
* state.
*
#endif
IMPLICIT NONE
*implicits
#include "comdim.cdk"
#include "comleg.cdk"
#include "comgd0.cdk"
#include "comsp.cdk"
*
#include <rpnmacros_f.h>
C
INTEGER ISIZ
INTEGER ILEN, JJ, JLATN, JLATS
S ,JM, ILONR,ILONIM,JK
S ,JN, ILA, INM, JM0, JNS, INS
C
REAL*8 ZJM,ZMONE
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)
! --------------------------
200 CONTINUE
*
!$OMP PARALLEL PRIVATE(DLALP,DLDALP,JM0,DLSP,ZFMS,ZFMA,JLATN,JLATS)
!$OMP DO PRIVATE(INM,ILA,JM,JN,JK,JJ,ZJM,JNS,INS,ILONR,ILONIM)
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, NKSDIM
DLSP(INM,1,JK) = SP(ILA,1,JK)
DLSP(INM,2,JK) = SP(ILA,2,JK)
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
*
CALL GETALP
(DLALP,DLDALP,NJLATH,NTRUNC,NTRUNCMX,JM)
*
C** . 2.3 Perform the inverse Legendre transform for all fields
*
CALL LEGINV3
(JM, ZFMS,ZFMA,DLSP,DLALP
S ,NKSDIM, NJLATH, NTRUNC, NTRUNCMX)
*
C** . 2.4 Passage to Fourier space
*
ILONR = 2*JM + 1
ILONIM = 2*JM + 2
ZMONE = -1.0d0
ZJM=FLOAT(JM)
DO 240 JJ = 1, NJLATH
JLATN = JJ
JLATS = 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, NJLATH, NTRUNC, NTRUNCMX)
C
ILONR = 2*JM + 1
ILONIM = 2*JM + 2
ZJM = FLOAT(JM)
DO JJ = 1, NJLATH
JLATN = JJ
JLATS = NJ - JLATN + 1
if(jlatn.ne.jlats) then
DO JK = 1, NFLEV
UT0(ILONR,JK,JLATN) = ZMONE*(UT0(ILONR,JK,JLATN)
S -(ZFMS(JJ,1,JK) + ZFMA(JJ,1,JK)))
UT0(ILONIM,JK,JLATN) = ZMONE*(UT0(ILONIM,JK,JLATN)
S -(ZFMS(JJ,2,JK)+ ZFMA(JJ,2,JK)))
VT0(ILONR,JK,JLATN) = ZMONE*(VT0(ILONR,JK,JLATN)
S +(ZFMS(JJ,1,JK+NFLEV) + ZFMA(JJ,1,JK+NFLEV)))
VT0(ILONIM,JK,JLATN) = ZMONE*(VT0(ILONIM,JK,JLATN)
S +(ZFMS(JJ,2,JK+NFLEV)+ ZFMA(JJ,2,JK+NFLEV)))
!
UT0(ILONR,JK,JLATS) = ZMONE*(UT0(ILONR,JK,JLATS)
S -(ZFMS(JJ,1,JK)- ZFMA(JJ,1,JK)))
UT0(ILONIM,JK,JLATS) = ZMONE*(UT0(ILONIM,JK,JLATS)
S -(ZFMS(JJ,2,JK)- ZFMA(JJ,2,JK)))
C
VT0(ILONR,JK,JLATS) = ZMONE*(VT0(ILONR,JK,JLATS)
S +(ZFMS(JJ,1,JK+NFLEV)- ZFMA(JJ,1,JK+NFLEV)))
VT0(ILONIM,JK,JLATS) = ZMONE*(VT0(ILONIM,JK,JLATS)
S +(ZFMS(JJ,2,JK+NFLEV)- ZFMA(JJ,2,JK+NFLEV)))
END DO
else
DO JK = 1, NFLEV
UT0(ILONR,JK,JLATN) = ZMONE*(UT0(ILONR,JK,JLATN)
S -(ZFMS(JJ,1,JK) + ZFMA(JJ,1,JK)))
UT0(ILONIM,JK,JLATN) = ZMONE*(UT0(ILONIM,JK,JLATN)
S -(ZFMS(JJ,2,JK)+ ZFMA(JJ,2,JK)))
VT0(ILONR,JK,JLATN) = ZMONE*(VT0(ILONR,JK,JLATN)
S +(ZFMS(JJ,1,JK+NFLEV) + ZFMA(JJ,1,JK+NFLEV)))
VT0(ILONIM,JK,JLATN) = ZMONE*(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 DO
!$OMP END PARALLEL
C
C** 3. DEALLOCATE LOCAL ARRAYS
C . -----------------------
300 CONTINUE
C
RETURN
END