!-------------------------------------- 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 SPA2SPAD 5,13
#if defined (DOC)
*
***s/r SPA2SPAD - Adjoint of the conversion of analysis variables to
* . model state increments
* .
* Purpose
* . To revert from the analysis variables to the model state variables
* . as defined in spectral space
* . IMPORTANT: after a call to SP2SPAD, the fields SPVOR and SPDIV
* . may contain PSI and CHI. This will eventually become
* . the default option
*Author : P. Gauthier *ARMA/AES March 25, 1996
* .
*Revision:
* . P. Gauthier *ARMA/AES August 5, 1996
* . Redefine the balance relationship and reformulate the
* . analysis variable in terms of PSI and CHI instead of
* . vorticity and divergence.
* . S. Pellerin *ARMA/AES Sept 97.
* - Control of the different model state of the 3Dvar
* through COMSTATE, COMSTATEC and COMSTNUM common
* blocks variables (comstate.cdk).
* . M. Buehner November 97
* . Introduce transformation for new definition of control
* . vector that reduces background covariances to the identity
* . matrix (for NANALVAR=3)
* . L. Fillion/M. Buehner 20 jul 98
* . Model's coordinate
* C. Charette *ARMA/AES Nov 1998
* Changed LDIVBAL to LBALDIV
* JM Belanger CMDA/SMC Aug 2000
* . 32 bits conversion
* (MXMA8, replace 2.0 by REAL*8 in SQRT)
* P. Koclas CMDA/SMC Apr 2003
* cnnert call to mxma8 to mxmaop
* M. Buehner *ARMA/SMC October 2004
* - changed damplibg to 3D and to be after balance operators
* L. Fillion *ARMA/EC 12 Jun 2009 - Update and validate option ldobal .false.
#endif
C
IMPLICIT NONE
#include "pardim.cdk"
#include "comdim.cdk"
#include "comcst.cdk"
#include "comgem.cdk"
#include "comleg.cdk"
#include "comcva.cdk"
#include "comsp.cdk"
#include "comspg.cdk"
#include "comgd0.cdk"
#include "compstat.cdk"
#include "comcorr.cdk"
#include "comcse1.cdk"
#include "comstate.cdk"
#include "comlun.cdk"
integer idum1,idum2,idum3,idum4
INTEGER ILEN
INTEGER JN,JM,JK,ILA
! REAL*8 ZSP(NKSDIM,2,0:NTRUNC),ZSP2(NKSDIM,2,0:NTRUNC),SQ2
! POINTER(PXZSP,ZSP), (PXZSP2,ZSP2)
real zmin,zmax
REAL*8 SQ2
REAL*8 ,ALLOCATABLE,DIMENSION(:,:,:) :: zsp,zsp2
INTEGER ILENGD, IERR, ILON, JLEV, JLON, JLAT, JLA
REAL*8 Z1SA2,ZNORMPSI, ZNORMCHI
! REAL*8 ZGDPSI(NIBEG:NIEND,NFLEV,NJBEG:NJEND)
! S ,ZGDCHI(NIBEG:NIEND,NFLEV,NJBEG:NJEND)
! POINTER (PXGDPSI,ZGDPSI),(PXGDCHI,ZGDCHI)
REAL*8,ALLOCATABLE,DIMENSION(:,:,:) :: ZGDPSI ,ZGDCHI
real*8 zgd(ni,nkgdim,nj)
!
INTEGER :: NTRUNCHF,JN0,INS,JNS,thdid,numthd
INTEGER :: omp_get_thread_num, omp_get_num_threads
real*8 two
*
two = 2.d0
*-------------------------------------------------------------------
C
C 1. Create local arrays for CHI and PSI if needed
C
100 CONTINUE
C
C . 1.1 Case where the analysis variable is X - Xb
C . (SPA2SP = Identity)
C
110 CONTINUE
IF(NANALVAR.EQ.0) THEN
RETURN
END IF
C
C . 1.2 Memory allocation
C
C
IF(NANALVAR.GE.2) THEN
ALLOCATE(ZGDPSI(NIBEG:NIEND,NFLEV,NJBEG:NJEND))
ALLOCATE(ZGDCHI(NIBEG:NIEND,NFLEV,NJBEG:NJEND))
END IF
C
C 2. Previous way of implementing the multivariate analysis
C . (based on the linear balance relationship)
C
200 CONTINUE
IF(NANALVAR.EQ.1) THEN
IF(NCNTVAR.EQ.2.AND.CFGERR.EQ.'G')THEN
CALL PROJIA
(.TRUE.)
CALL GDSPA
CALL FGERR
('M')
CALL SPGDA
CALL FGERR
('O')
END IF
IF(NCNTVAR.EQ.2.AND.CFGERR.EQ.'S')THEN
CALL PROJIA
(.TRUE.)
CALL FGERR
('M')
END IF
RETURN
END IF
C
C 3. Formulation using PSI and CHI
C
300 CONTINUE
C
C . 3.5 Rederive the vorticity and divergence from PSI and CHI
C
Z1SA2 = 1./(RA*RA)
!$OMP PARALLEL DO
DO JLEV = 1, NFLEV
DO JLA = 1, NLA
SPVOR(JLA,1,JLEV) = SPVOR(JLA,1,JLEV)*Z1SA2*RNNP1(JLA)
SPVOR(JLA,2,JLEV) = SPVOR(JLA,2,JLEV)*Z1SA2*RNNP1(JLA)
SPDIV(JLA,1,JLEV) = SPDIV(JLA,1,JLEV)*Z1SA2*RNNP1(JLA)
SPDIV(JLA,2,JLEV) = SPDIV(JLA,2,JLEV)*Z1SA2*RNNP1(JLA)
END DO
END DO
!$OMP END PARALLEL DO
C
C . 3.4 Spectral transform all fields
C
340 CONTINUE
CALL SPEREE
(NKSDIM,SP,GD
S ,NLA,NIBEG,NIEND,NJBEG,NJEND,NKSDIM)
C
C Apply 3D amplification factor
C
!cluc
! DO JLAT = 1, NJ
! DO JLEV = 1, NKGDIM
! DO JLON = 1, NI
! zgd(JLON,JLEV,JLAT)= damplibg(JLON,JLEV,JLAT)
! END DO
! END DO
! END DO
! call maxmin(zgd,ni,nj,nkgdim,zmin,zmax,
! & idum1,idum2,idum3,idum4,'spa2spad ',
! & 'ZGD')
! DO JLAT = 1, NJ
! DO JLEV = 1, NKGDIM
! DO JLON = 1, NI
! zgd(JLON,JLEV,JLAT)= gd(JLON,JLEV,JLAT)
! END DO
! END DO
! END DO
! call maxmin(zgd,ni,nj,nkgdim,zmin,zmax,
! & idum1,idum2,idum3,idum4,'spa2spad ',
! & 'ZGD')
!
DO JLAT = 1, NJ
DO JLEV = 1, NKGDIM
DO JLON = 1, NI
GD(JLON,JLEV,JLAT)=
+ GD(JLON,JLEV,JLAT)*damplibg(JLON,JLEV,JLAT)
END DO
END DO
END DO
C
C . 3.3 Obtain the full geopotential
C
DO JLAT = 1, NJ
ILON = NILON(JLAT)
DO JLEV = 1, NFLEV
DO JLON = 1, ILON
ZGDPSI(JLON,JLEV,JLAT) = UT0(JLON,JLEV,JLAT)
ZGDCHI(JLON,JLEV,JLAT) = VT0(JLON,JLEV,JLAT)
END DO
END DO
END DO
if(ldobal) then
IF(LBALDIV) THEN
CALL ADIVBAL
(zgdpsi,zgdchi)
ENDIF
call abmass
(zgdpsi,zgdchi)
endif
C
C . 3.2 Multiply by the background error variances
C
write(nulout,*) 'spa2spad: luse3dstd = ',LUSE3DSTD
IF(LUSE3DSTD) THEN
DO JLEV = 1, NFLEV
DO JLAT = 1, NJ
ILON = NILON(JLAT)
DO JLON = 1, ILON
UT0(JLON,JLEV,JLAT)=
+ ZGDPSI(JLON,JLEV,JLAT)*RGSIGUU3D(JLON,JLEV,JLAT)
VT0(JLON,JLEV,JLAT)=
+ ZGDCHI(JLON,JLEV,JLAT)*RGSIGVV3D(JLON,JLEV,JLAT)
TT0(JLON,JLEV,JLAT)=
+ TT0(JLON,JLEV,JLAT)*RGSIGTT3D(JLON,JLEV,JLAT)
Q0(JLON,JLEV,JLAT)=
+ Q0(JLON,JLEV,JLAT)*RGSIGQ3D(JLON,JLEV,JLAT)
END DO
END DO
END DO
DO JLAT = 1, NJ
ILON = NILON(JLAT)
DO JLON = 1, ILON
GPS0(JLON,1,JLAT)=
+ GPS0(JLON,1,JLAT)*RGSIGPS3D(JLON,1,JLAT)
END DO
END DO
IF(nsexist(nstg).eq.1) THEN
DO JLAT = 1, NJ
ILON = NILON(JLAT)
DO JLON = 1, ILON
GTG0(JLON,1,JLAT)=
+ GTG0(JLON,1,JLAT)*RGSIGTG3D(JLON,1,JLAT)
END DO
END DO
ENDIF
ELSE
320 CONTINUE
CALL FGERR
('M')
!$OMP PARALLEL DO PRIVATE(ILON,ZNORMPSI,ZNORMCHI)
DO JLAT = 1, NJ
ILON = NILON(JLAT)
DO JLEV = 1, NFLEV
ZNORMPSI = RGSIGUU(JLAT,JLEV)
ZNORMCHI = RGSIGVV(JLAT,JLEV)
DO JLON = 1, ILON
UT0(JLON,JLEV,JLAT)
S = ZGDPSI(JLON,JLEV,JLAT)*ZNORMPSI
VT0(JLON,JLEV,JLAT)
S = ZGDCHI(JLON,JLEV,JLAT)*ZNORMCHI
END DO
END DO
END DO
!$OMP END PARALLEL DO
ENDIF
C
C . 3.1 Transform PSI and CHI and the other variables
C . to physical space
C
310 CONTINUE
C
CALL REESPE
(NKSDIM,SP,GD
S ,NLA,NIBEG,NIEND,NJBEG,NJEND,NKSDIM)
C
C 3.0 Multiply by CORNS to add in vert corr (scale/rotate)
C
IF(NANALVAR.EQ.3) THEN
NTRUNCHF=NTRUNC/2
SQ2=sqrt(two)
!$OMP PARALLEL PRIVATE(thdid,numthd,JN,JM,JK,ILA,JN0,INS,JNS,zsp,zsp2)
thdid=omp_get_thread_num()
numthd=omp_get_num_threads()
ALLOCATE(zsp(NKSDIM,2,0:NTRUNC))
ALLOCATE(zsp2(NKSDIM,2,0:NTRUNC))
DO JN0 = thdid, NTRUNCHF,numthd
INS=1
IF(JN0 == (NTRUNC-JN0))INS=0
DO 222 JNS=0,INS
JN=(1-JNS)*JN0+JNS*(NTRUNC-JN0)
C
C* Transfer Dx from global to working array
C
DO JM = 0, JN
ILA = NIND(JM) + JN - JM
DO JK = 1, NKSDIM
ZSP(JK,1,JM) = SP(ILA,1,JK)
ZSP(JK,2,JM) = SP(ILA,2,JK)
END DO
END DO
C
C* Compute (CORNS(NKSDIM,NKSDIM,JN)' * ZSP)
C
CALL MXMAOP
(CORNS(1,1,JN,1),1,NKSDIM,ZSP(1,1,0),1,NKSDIM
+ ,ZSP2(1,1,0),1,NKSDIM,NKSDIM,NKSDIM,2*(JN+1))
C
C* Transfer from working array to global array
C
ILA = NIND(0) +JN
DO JK = 1, NKSDIM
SP(ILA,1,JK) = ZSP2(JK,1,0)*SQ2
SP(ILA,2,JK) = ZSP2(JK,2,0)*SQ2
END DO
C
DO JM = 1, JN
ILA = NIND(JM) +JN - JM
DO JK = 1, NKSDIM
SP(ILA,1,JK) = ZSP2(JK,1,JM)
SP(ILA,2,JK) = ZSP2(JK,2,JM)
END DO
END DO
C
222 CONTINUE
END DO
DEALLOCATE(zsp)
DEALLOCATE(zsp2)
!$OMP END PARALLEL
C
DO JK=1,NKSDIM
DO ILA=1,NTRUNC+1
SP(ILA,2,JK)=0.0
ENDDO
ENDDO
C
ENDIF
C
C 9. Deallocate local arrays
C
900 CONTINUE
C
C
IF(NANALVAR.GE.2) THEN
DEALLOCATE(ZGDPSI)
DEALLOCATE(ZGDCHI)
END IF
C
RETURN
END