!-------------------------------------- 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 meancvgd(kulstat,koutfile,ldhelm,ldfplane,ldtb_psi, 1,97
     &                    ldcoriol_one,pmfpp,pmfcu,pmftu,pmflq,pmfpsu,pmftg)
#if defined (DOC)
*
* * s/r meancvgd: Compute mean and variances of gridpoint control-vector 
*                 variables from forecast error samples for the 
*                 LAM4D configuration. Then, store on file.
*
* Author:   L. Fillion - ARMA/MSC -  1 Oct 2005.
* Revision: L. Fillion - ARMA/EC  - 27 Sep 2006. 
*                        Introduce filtering of s-dev option.
*           JF Caron   - ARMA/EC  - 15 Fev 2008. 
*                        Computation of stats for extra variables (Tb, Psb)
*                        and ratio of TT-var explained by Pb
*           L. Fillion - ARMA/EC  - Correct bug: Include unbiasing step on GD prior to initgdla.ftn (as for ptotla.ftn)
*           L. Fillion - ARMA/EC  - 12 Dec 2008 - Get rid of comstdd.cdk and include temporary treatment of TG.
*           L. Fillion - 18 Nov 2008 - Introduce pmfpp,pmfcu,pmftu,pmflq,pmfpsu,pmftg fields in case of MC case.
*           L. Fillion - 1st May 2009 - Finalize the vort,div approach.
!           L. Fillion - 27 Jan 2010 - Reduce considerably the amount of arrays needed so as to allow very high res. analysis.
*
* Arguments: 
*     KULSTAT  logical unit number
*
#endif
      IMPLICIT NONE
*implicits
#include "taglam4d.cdk"
#include "pardim.cdk"
#include "comdim.cdk"
#include "comlun.cdk"
#include "comct0.cdk"
#include "comcst.cdk"
#include "comgem.cdk"
#include "comsp.cdk"
#include "comgd0.cdk"
#include "comgdpar.cdk"
#include "comcse1.cdk"
#include "comgrd_param.cdk"
#include "comgrd.cdk"
#include "comgemla.cdk"
#include "comfftla.cdk"
#include "comcva.cdk"
#include "namgdpar.cdk"
*
      logical ldhelm,ldfplane,ldtb_psi,llzdpc,ldcoriol_one
      integer iflag
      INTEGER KULSTAT,koutfile
      real*8 pmfpp(ni,nflev,nj)
      real*8 pmfcu(ni,nflev,nj)
      real*8 pmftu(ni,nflev,nj)
      real*8 pmflq(ni,nflev,nj)
      real*8 pmfpsu(ni,nj)
      real*8 pmftg(ni,nj)
!
      logical llfiltersdev,llfilt,llvproj,llvfilt
      INTEGER JENS, IENS, jk1, IERR, JFILE, iensemble
      integer ji,jj,jk,jvar,inbvar,icase,itrunc,inip1,injp1
      parameter(inbvar=20)
      character*1 clpart,clflt,clgrid
      character*2 clvar(inbvar)
!
      CHARACTER*8 cletik
      INTEGER FNOM, FSTOUV, FSTFRM, FCLOS, FSTPRM, FSTINF, FSTINL
      INTEGER VFSTECR, vfstlir
      integer newdate
!
!     RPN Standard files parameters
!
      INTEGER INI,INJ,INK, INPAS, INBITS, IDATYP, IDEET
     +     ,IP1,IP2,IP3,IG1,IG2,IG3,IG4,ISWA,ILENGTH,IDLTF
     +     ,IUBC,IEXTR1,IEXTR2,IEXTR3
      INTEGER ILISTE(100),IDATE(100), IDATV(100), IDIMAX, INFON, IFSTRUN, IHH
      real*8 HEURES
      CHARACTER*1 CLTYPVAR,CLGRTYP
      CHARACTER*2 CLNOMVAR
      CHARACTER*8 CLETIKET
      character (len = 1) :: etk1
      character (len = 2) :: etk2
      character (len = 2) :: etk
!
      integer idum1,idum2,idum3,idum4
      real*8 zthreshold,zni
      REAL*8 DLA2
      REAL*8 ZFACT,zmin,zmax
      INTEGER IPAK, IDATEO
      CHARACTER*128 CLFLFILE
!
      real*8 zx8(ni),zy8(nj)
!
      real*8 bufwk3dijk(ni,nj,nflev)
!
      real*8 zmpp(ni,nflev,nj)
      real*8 zmcu(ni,nflev,nj)
      real*8 zmtu(ni,nflev,nj)
      real*8 zmlq(ni,nflev,nj)
      real*8 zmpsu(ni,nj)
      real*8 zmtg(ni,nj)
!
      real*8 zspp(ni,nflev,nj)
      real*8 zscu(ni,nflev,nj)
      real*8 zstu(ni,nflev,nj)
      real*8 zslq(ni,nflev,nj)
      real*8 zspsu(ni,nj)
      real*8 zstg(ni,nj)
!
      real*8 zmeanpsu(ni,nj)
      real*8 zmeantg(ni,nj)
      real*8 zmeanpb(ni,nj)
!
      real*8 zstdpsu(ni,nj)
      real*8 zstdtg(ni,nj)
      real*8 zstdpb(ni,nj)
!
      real*8 z2d(ni,nj)
      real*8 z2din(mni_in,mnj_in)
      real*8 zwork2d(ni,nj)
!
      real*8 zgdpsi(ni,nflev,nj)
      real*8 zgdchi(ni,nflev,nj)
      real*8 zvort(ni,nflev,nj)
      real*8 zdiv(ni,nflev,nj)
      real*8 zpb(ni,nflev,nj)
      real*8 zt0(ni,nflev,nj)
      real*8 ztb(ni,nflev,nj)
      real*8 ztu(ni,nflev,nj)
!
      real*8 zpsb(ni,nj)
      real*8 zpsu(ni,nj)
!
      real*8 zxmtb(ni,nflev,nj)
      real*8 zstb(ni,nflev,nj)
      real*8 zxmpsb(ni,nj)
      real*8 zspsb(ni,nj)
!  
!      real*8 zwrkflt(ni,nflev,nj)
!
!!
      llfiltersdev = .false.
      llfilt = .false.
      if(lflt_low) then
        llfilt = .true.
        clflt = 'L'
      else if(lflt_high) then
        llfilt = .true.
        clflt = 'H'
      endif
      llvfilt = .false.
      llvproj = .false.
!
      clvar(1) = 'PP' 
      clvar(2) = 'CC' 
      if(.not.ldhelm) then
        clvar(1) = 'QQ' 
        clvar(2) = 'DD' 
      endif
      llzdpc = .false.
      if(cfstvar(1).eq.'UU'.or.cfstvar(1).eq.'QQ'.and.ldhelm) llzdpc = .true.
      write(nulout,*) 'meancvgd: cfstvar(1) = ',cfstvar(1)
      write(nulout,*) 'meancvgd: llzdpc = ',llzdpc
      clvar(3) = 'TT' 
      clvar(4) = 'LQ' 
      clvar(5) = 'TB' 
      clvar(6) = 'PU' 
      clvar(7) = 'PB' 
!
!
      inip1 = ni+1
      injp1 = nj+1
!
      zthreshold = 1.e-15
      DLA2 = DBLE(RA) * DBLE(RA)
      IDIMAX = 100
      write(nulout,*) 'meancvgd: BEGIN '
      write(nulout,*) 'meancvgd: NFLSTAT = ',NFLSTAT

!     -----------------------
!*1   Initialize accumulators
!     -----------------------

      zmpp(:,:,:) = 0.0
      zmcu(:,:,:) = 0.0
      zmtu(:,:,:) = 0.0
      zmlq(:,:,:) = 0.0
      zxmtb(:,:,:) = 0.0
      zmpsu(:,:) = 0.0
      zxmpsb(:,:) = 0.0
      zmtg(:,:) = 0.0
!
      zspp(:,:,:) = 0.0
      zscu(:,:,:) = 0.0
      zstu(:,:,:) = 0.0
      zstb(:,:,:) = 0.0
      zslq(:,:,:) = 0.0
      zspsb(:,:) = 0.0
      zspsu(:,:) = 0.0
      zstg(:,:) = 0.0
!
!cluc      call getmeangdla(nulbgst)  ! read from file: precomputed mean of fields
!
!     -------------------------------------------------------------
!*2   Access the increments from a set of files (loop on the files)
!     -------------------------------------------------------------
!
      iensemble = 0

      DO 201 JFILE = 1, NFLSTAT
!
         call openinc(kulstat,jfile)
!
!*    .  2.1 Find how many cases there are to be treated
!
         IP1 = -1
         IP2 = -1
         IP3 = -1
         CLNOMVAR = 'P0'
         write(NULOUT,*)
         IERR = FSTINL (KULSTAT,INI,INJ,INK
     &        ,-1,CETIKETERR,IP1,IP2,IP3,' '
     &        ,CLNOMVAR,ILISTE,INFON,IDIMAX)
!
!        Ensure namelist parameters mextendx, mextendy are feasible
!        w.r.t. basic-state GEM fields
!        ----------------------------------------------------------
!
         if(mni_in.gt.ini) then
           write(nulout,*) 'meancvgd: mni_in,INI=',mni_in,INI
           call abort3d(nulout,'meancvgd: mni_in.gt.INI')
           write(nulout,*) 'meancvgd: verify CETIKETERR in NAMCSE1 as compared to input file etiket...'
         else if(mnj_in.gt.inj) then
           write(nulout,*) 'meancvgd: mnj_in,INJ=',mnj_in,INJ
           call abort3d(nulout,'meancvgd: mnj_in.gt.INJ')
         endif
!
         WRITE(NULOUT,9210)INFON
 9210    FORMAT(//,4X,"meancvgd: Ensemble of ",I4," increments")
         IF(INFON.EQ.0) THEN
            WRITE(NULOUT,*)' THIS FILE IS EMPTY. 
     $           CHECK THE SELECTION CRITERIA'
            CALL ABORT3D(NULOUT,'meancvgd: problem with FSTINL')
         END IF
         IENS = INFON
!
!*       2.2  Get all the dates at which increments are available
!
         DO JENS = 1, IENS
            IERR = FSTPRM(ILISTE(JENS),IDATE(JENS),IDEET,INPAS
     &           ,INI,INJ,INK, INBITS, IDATYP
     &           ,IP1,IP2,IP3,CLTYPVAR,CLNOMVAR,CLETIKET,CLGRTYP
     &           ,IG1,IG2,IG3,IG4,ISWA,ILENGTH,IDLTF
     &           ,IUBC,IEXTR1,IEXTR2,IEXTR3)
!
!        Ensure namelist parameters mextendx, mextendy are feasible
!        w.r.t. basic-state GEM fields
!        ----------------------------------------------------------
!
            if(grd_typ.eq.'LU') then
              if(mni_in.gt.ini) then
                write(nulout,*) 'MEANCVGD: INI, mni_in = ',INI,mni_in
                call abort3d(nulout,'MEANCVGD: mni_in.gt.INI')
              else if(mnj_in.gt.inj) then
                write(nulout,*) 'MEANCVGD: INJ, mnj_in = ',INJ,mnj_in
                call abort3d(nulout,'MEANCVGD: mnj_in.gt.INJ')
              endif
            endif
!
            heures = real(INPAS*IDEET/3600)
!
            CALL INCDATR(IDATV(JENS),IDATE(JENS),heures)
            ierr= NEWDATE(IDATV(JENS),IFSTRUN,IHH,-3)
            WRITE(NULOUT,9320)JENS, IFSTRUN,IHH
         END DO
 9320    FORMAT(5X,"Case No. ",I3,5x,"Date and time: ",I10,5x,I8)
!
         IF(iensemble.EQ.0) THEN
            NDATESTAT = IDATE(1)
         END IF
!
!        2.3 Loop on the ensemble contained in the current file of differences
!
         DO 231 JENS = 1, IENS
!
!          Get the increment in grid-point form
!
           if(lmcstats) then
             NSTAMPN = -1
             icase=jens       ! i.e. will use IP3 as a search parameter and ignore the date in current file since all same
           else
             NSTAMPN = IDATE(JENS) ! i.e. will use the current date of validity of the current error sample in standard file.
             icase = -1  ! ignore IP3 as a search parameter in vfstlir
           endif
!
! 2.3.1 
           call geterr(kulstat,'G','E',icase)
!
           do ji=1,ni
           do jj=1,nj
             zpsb(ji,jj)=gps0(ji,1,jj)
           enddo
           enddo
!
           call maxmin(zpsb,ni,nj,1,zmin,zmax,
     &            idum1,idum2,idum3,idum4,'meancvgd    ',
     &            'P0 ')
!
           if(cfstvar(1).eq.'PP'.or.cfstvar(1).eq.'QQ') then 
             iflag = 0
             do jk=1,nflev
               do jj=1,nj
                 do ji=1,ni
                   zgdpsi(ji,jk,jj) = ut0(ji,jk,jj)
                   zgdchi(ji,jk,jj) = vt0(ji,jk,jj)
                   if(ut0(ji,jk,jj).ne.0.) iflag = 1
                 enddo
               enddo
             enddo
           endif 
!
           if(iflag.eq.0) then
             write(nulout,*) 'meancvgd: Variable',cfstvar(1)
             call abort3d(nulout,'meancvgd: field identically zero') 
           endif
!
!2.3.4     Subtract forecast mean if in MC mode
!
           write(nulout,*) 'meancvgd: lmcstats = ',lmcstats
           if(lmcstats) then
             write(nulout,*) 'meancvgd: Subtract fsct ensemble mean from current fcst sample'
             do jj=1,nj
               do ji=1,ni
                 do jk=1,nflev
                   zgdpsi(ji,jk,jj) = zgdpsi(ji,jk,jj) - pmfpp(ji,jk,jj)
                   zgdchi(ji,jk,jj) = zgdchi(ji,jk,jj) - pmfcu(ji,jk,jj)
                   tt0(ji,jk,jj) = tt0(ji,jk,jj) - pmftu(ji,jk,jj)
                   q0(ji,jk,jj) = q0(ji,jk,jj) - pmflq(ji,jk,jj)
                 enddo
                 gps0(ji,1,jj) = gps0(ji,1,jj) - pmfpsu(ji,jj)
                 gtg0(ji,1,jj) = gtg0(ji,1,jj) - pmftg(ji,jj)
               enddo
             enddo
           endif
!
           ztb(:,:,:) = 0.0
           zpb(:,:,:) = 0.0
           zpsb(:,:) = 0.0
           zt0(:,:,:) = tt0(:,:,:)
!
           write(nulout,*) ' '
           write(nulout,*) 'meancvgd: ****************************************'
           write(nulout,*) 'meancvgd: Uses Balance operators of Order = ',mbal_order
           write(nulout,*) 'meancvgd: ****************************************'
           write(nulout,*) ' '
!
           if(mbal_order.eq.1) then
!
!            Linear-Geostrophic
!            ------------------
!
             write(nulout,*) 'meancvgd: ldcoriol_one =',ldcoriol_one
             if(ldtb_psi) then
                clgrid = 'P'
                if (lpsifromglb) clgrid = 'S'
                call linbal_la(zpb,zgdpsi,ldfplane,ldcoriol_one,clgrid)
             else
               do jk=1,nflev
               do jj=1,nj
               do ji=1,ni
                 zpb(ji,jk,jj)=zvort(ji,jk,jj)
                enddo
                enddo
                enddo
             endif
             if(llvfilt) call vfilt(zpb,5,'L')
             if(llvproj) call vproj(zpb,zpb,nflev)
             call p2tpsb(ztb,zpsb,zpb,cptot,.false.,.true.)  ! build Tb, psb from psi
             if(llvproj) call viproj(ztb,ztb,nflev)
!
           else if(mbal_order.eq.2) then
!
!            Tangent-Linear  2nd-Order Balance
!            ---------------------------------
!
           endif
!
            if(llvfilt) call vfilt(zt0,5,'L')
            call maxmin(zt0,ni,nj,nflev,zmin,zmax,
     &           idum1,idum2,idum3,idum4,'meancvgd    ',
     &           'T0 ')
!
           do ji=1,ni
           do jj=1,nj
             do jk=1,nflev
               ut0(ji,jk,jj) = zgdpsi(ji,jk,jj)
               vt0(ji,jk,jj) = zgdchi(ji,jk,jj)
               ztu(ji,jk,jj) = zt0(ji,jk,jj)-ztb(ji,jk,jj)
!               zsavet0(ji,jj,jk) = zt0(ji,jk,jj)
!               zsavetb(ji,jj,jk) = ztb(ji,jk,jj)
!               zsavepsi(ji,jj,jk) = zgdpsi(ji,jk,jj)
!               zsavevort(ji,jj,jk) = zvort(ji,jk,jj)
             enddo
           enddo
           enddo
!   
!          ensure T_u, Ps_u are bi-periodic
!
           if(mbal_order.ne.0) then
             call mach3(ztu,ni,nj,nflev,inip1,injp1)
             call mach3(ztb,ni,nj,nflev,inip1,injp1)
!
!            Transfer Tu to tt0
!
             do ji=1,ni
             do jj=1,nj
             do jk=1,nflev
                tt0(ji,jk,jj) = ztu(ji,jk,jj)
             enddo
             enddo
             enddo
!
             call maxmin(ztb,ni,nj,nflev,zmin,zmax,
     &            idum1,idum2,idum3,idum4,'meancvgd    ',
     &            'TB ')
!
             call maxmin(ztu,ni,nj,nflev,zmin,zmax,
     &            idum1,idum2,idum3,idum4,'meancvgd    ',
     &            'TU ')
!
             do ji=1,ni
             do jj=1,nj
                zpsu(ji,jj)=gps0(ji,1,jj)-zpsb(ji,jj)
             enddo
             enddo
!
             write(nulout,*) 'meancvgd: Balanced Surface Pressure'
             call maxmin(zpsb,ni,nj,1,zmin,zmax,
     &            idum1,idum2,idum3,idum4,'meancvgd    ',
     &            'P0B')
!
!            ensure ps_u is bi-periodic
!
             call mach2(zpsu,ni,nj,inip1,injp1)
!
             call maxmin(zpsu,ni,nj,1,zmin,zmax,
     &              idum1,idum2,idum3,idum4,'meancvgd    ',
     &              'P0U')
             write(nulout,*) 'meancvgd: point 2: P0, (min,max)=',zmin,zmax
!
             do ji=1,ni
             do jj=1,nj
                gps0(ji,1,jj) = zpsu(ji,jj)
             enddo
             enddo
           endif
!
!          ACCUMULATE SUM OF ELEMENTS AND SUM OF SQUARED ELEMENTS
!
!    ldhelm:       Loop for Psi, Chi, Tu, lq, ps_u, TG  (control variables)
!    .not.ldhelm:  Loop for Vort ,  Div, Tu, lq, ps_u, TG  (control variables)
!
           do jj = 1, nj
             do ji = 1, ni
               do jk1 = 1, nflev
                 zmpp(ji,jk1,jj) = zmpp(ji,jk1,jj) + ut0(ji,jk1,jj)
                 zmcu(ji,jk1,jj) = zmcu(ji,jk1,jj) + vt0(ji,jk1,jj)
                 zmtu(ji,jk1,jj) = zmtu(ji,jk1,jj) + tt0(ji,jk1,jj)
                 zmlq(ji,jk1,jj) = zmlq(ji,jk1,jj) + q0(ji,jk1,jj)
!
                 zspp(ji,jk1,jj) = zspp(ji,jk1,jj)
     &               + ut0(ji,jk1,jj)*ut0(ji,jk1,jj)
                 zscu(ji,jk1,jj) = zscu(ji,jk1,jj)
     &               + vt0(ji,jk1,jj)*vt0(ji,jk1,jj)
                 zstu(ji,jk1,jj) = zstu(ji,jk1,jj)
     &               + tt0(ji,jk1,jj)*tt0(ji,jk1,jj)
                 zslq(ji,jk1,jj) = zslq(ji,jk1,jj)
     &               + q0(ji,jk1,jj)*q0(ji,jk1,jj)
               enddo
               zmpsu(ji,jj) = zmpsu(ji,jj) + gps0(ji,1,jj)
               zmtg(ji,jj) = zmtg(ji,jj) + tt0(ji,nflev,jj)  ! Use Tu at surface until TG is abord.... cluc
!
               zspsu(ji,jj) = zspsu(ji,jj)
     &               + gps0(ji,1,jj)*gps0(ji,1,jj)
               zstg(ji,jj) = zstg(ji,jj)
     &               + tt0(ji,nflev,jj)*tt0(ji,nflev,jj) ! Use Tu at surface until TG is abord.... cluc
             enddo
           enddo

!          Loop for Tb, Psb (extra variables)
!
           do jj = 1, nj
             do jk1 = 1, NFLEV
               do ji = 1, ni
                 zxmtb(ji,jk1,jj) = zxmtb(ji,jk1,jj) + ztb(ji,jk1,jj)
                 zstb(ji,jk1,jj) = zstb(ji,jk1,jj) +
     $                ztb(ji,jk1,jj)*ztb(ji,jk1,jj)
               enddo
             enddo
           enddo
           do jj = 1, nj
              do ji = 1, ni
                 zxmpsb(ji,jj) = zxmpsb(ji,jj) + zpsb(ji,jj)
                 zspsb(ji,jj) = zspsb(ji,jj) +
     $                zpsb(ji,jj)*zpsb(ji,jj)
              enddo
           enddo
!
 231     continue ! end loop on the ensemble within current file
!
         iensemble = iensemble + IENS
         IERR =  FSTFRM (KULSTAT)
         IERR =  FCLOS  (KULSTAT)
 201  CONTINUE   ! end loop on jfile
      write(nulout,*) 'meancvgd: iensemble = ',iensemble
!
      write(nulout,*) 'meancvgd: zxmtb apres loop sur les echantillons'
      if(mbal_order.gt.0) then
        call maxmin(zxmtb,ni,nj,nflev,zmin,zmax,
     &            idum1,idum2,idum3,idum4,'meancvgd    ',
     &            'TB ')
      endif
!
!     ----------------------------------------
!*3.  COMPUTE VARIANCES OF GD FOR THE ENSEMBLE
!     ----------------------------------------

!     Loop for Psi, Chi, Tu, Psu, lq (control variables)
      DO jj = 1,nj
        DO ji = 1,ni
          do jk1 = 1, nflev
            zspp(ji,jk1,jj) = ( zspp(ji,jk1,jj) -
     +      ((zmpp(ji,jk1,jj)*zmpp(ji,jk1,jj)) / iensemble )) /
     +       (iensemble - 1)
            if (zspp(ji,jk1,jj).eq.0.0) zspp(ji,jk1,jj)=zthreshold
!
            zscu(ji,jk1,jj) = ( zscu(ji,jk1,jj) -
     +      ((zmcu(ji,jk1,jj)*zmcu(ji,jk1,jj)) / iensemble )) /
     +       (iensemble - 1)
            if (zscu(ji,jk1,jj).eq.0.0) zscu(ji,jk1,jj)=zthreshold
!
            zstu(ji,jk1,jj) = ( zstu(ji,jk1,jj) -
     +      ((zmtu(ji,jk1,jj)*zmtu(ji,jk1,jj)) / iensemble )) /
     +       (iensemble - 1)
            if (zstu(ji,jk1,jj).eq.0.0) zstu(ji,jk1,jj)=zthreshold
!
            zslq(ji,jk1,jj) = ( zslq(ji,jk1,jj) -
     +      ((zmlq(ji,jk1,jj)*zmlq(ji,jk1,jj)) / iensemble )) /
     +       (iensemble - 1)
            if (zslq(ji,jk1,jj).eq.0.0) zslq(ji,jk1,jj)=zthreshold
          ENDDO
          zspsu(ji,jj) = ( zspsu(ji,jj) -
     +    ((zmpsu(ji,jj)*zmpsu(ji,jj)) / iensemble )) /
     +     (iensemble - 1)
          if (zspsu(ji,jj).eq.0.0) zspsu(ji,jj)=zthreshold
          zstg(ji,jj) = ( zstg(ji,jj) -
     +    ((zmtg(ji,jj)*zmtg(ji,jj)) / iensemble )) /
     +     (iensemble - 1)
          if (zstg(ji,jj).eq.0.0) zstg(ji,jj)=zthreshold
        ENDDO
      ENDDO

!     Loop for Tb, Psb (extra variables)
      DO jj = 1,nj
        do jk1 = 1, NFLEV
          DO ji = 1,ni
            ZSTB(ji,jk1,jj) = ( ZSTB(ji,jk1,jj) -
     +      ((ZXMTB(ji,jk1,jj)*ZXMTB(ji,jk1,jj)) / iensemble )) /
     +       (iensemble - 1)
            if (zstb(ji,jk1,jj).eq.0.0) zstb(ji,jk1,jj)=zthreshold
          ENDDO
        ENDDO
      ENDDO
      DO jj = 1,nj
          DO ji = 1,ni
            ZSPSB(ji,jj) = ( ZSPSB(ji,jj) -
     +      ((ZXMPSB(ji,jj)*ZXMPSB(ji,jj)) / iensemble )) /
     +       (iensemble - 1)
            if (zspsb(ji,jj).eq.0.0) zspsb(ji,jj)=zthreshold
          ENDDO
      ENDDO
!     ---------------------------------------
!*4.  COMPUTE THE MEAN OF GD FOR THE ENSEMBLE
!     ---------------------------------------

c     Loop for Psi, Chi, Tu, Psu, lq (control variables)
      DO jj = 1,nj
        DO ji = 1,ni
          do jk1 = 1, nflev
            zmpp(ji,jk1,jj) = zmpp(ji,jk1,jj) / iensemble
            zmcu(ji,jk1,jj) = zmcu(ji,jk1,jj) / iensemble
            zmtu(ji,jk1,jj) = zmtu(ji,jk1,jj) / iensemble
            zmlq(ji,jk1,jj) = zmlq(ji,jk1,jj) / iensemble
          ENDDO
          zmpsu(ji,jj) = zmpsu(ji,jj) / iensemble
          zmtg(ji,jj) = zmtg(ji,jj) / iensemble
        ENDDO
      ENDDO

c     Loop for Tb, Psb (extra variables)
      DO jj = 1,nj
        do jk1 = 1, NFLEV
          DO ji = 1,ni
            ZXMTB(ji,jk1,jj) = ZXMTB(ji,jk1,jj) / iensemble
          ENDDO
        ENDDO
      ENDDO
      DO jj = 1,nj
          DO ji = 1,ni
            ZXMPSB(ji,jj) = ZXMPSB(ji,jj) / iensemble
          ENDDO
      ENDDO

!     Transfer data
      do jj=1,nj
        do ji=1,ni
          zmeanpsu(ji,jj) = zmpsu(ji,jj)
          zmeantg(ji,jj) = zmtg(ji,jj)
          zmeanpb(ji,jj) = zxmpsb(ji,jj)
!          do jk=1,nflev
!            zmeanpp(ji,jj,jk) = zmpp(ji,jk,jj)
!            zmeancu(ji,jj,jk) = zmcu(ji,jk,jj)
!            zmeantu(ji,jj,jk) = zmtu(ji,jk,jj)
!            zmeanlq(ji,jj,jk) = zmlq(ji,jk,jj)
!            zmeantb(ji,jj,jk) = zxmtb(ji,jk,jj)
!          enddo
        enddo
      enddo
!
      call maxmin(zmeanpsu,ni,nj,1,zmin,zmax,
     &            idum1,idum2,idum3,idum4,'meancvgd    ',
     &            'PUM')
      call maxmin(zmeantg,ni,nj,1,zmin,zmax,
     &            idum1,idum2,idum3,idum4,'meancvgd    ',
     &            'TGM')
!
!     -------
!*5   STD DEV
!     -------

      do jj=1,nj
        do ji=1,ni
!         Psu
          if (zspsu(ji,jj).le.0.0) then
           write(nulout,*) 'meancvgd: ji,jj,zspsu=',
     &            ji,jj,zspsu(ji,jj)
           call abort3d(nulout,'meancvgd: ps_u st-dev problem!')
          endif
          zspsu(ji,jj) = sqrt(zspsu(ji,jj))
          zstdpsu(ji,jj) = zspsu(ji,jj)
!         TG
          if (zstg(ji,jj).le.0.0) then
           write(nulout,*) 'meancvgd: ji,jj,zstg=',
     &            ji,jj,zstg(ji,jj)
           call abort3d(nulout,'meancvgd: TG st-dev problem!')
          endif
          zstg(ji,jj) = sqrt(zstg(ji,jj))
          zstdtg(ji,jj) = zstg(ji,jj)

!         Psb
          if (zspsb(ji,jj).le.0.0) then
           write(nulout,*) 'meancvgd: ji,jj,zspspb=',
     &            ji,jj,zspsb(ji,jj)
           call abort3d(nulout,'meancvgd: ps_b st-dev problem!')
          endif
          zstdpb(ji,jj) = sqrt(zspsb(ji,jj))

          do jk=1,nflev
!
            if(zspp(ji,jk,jj).le.0.0) then        ! Psi
              write(nulout,*) 'meancvgd: ji,jk,jj,zspp=',
     &                         ji,jk,jj,zspp(ji,jk,jj)
              call abort3d(nulout,'meancvgd: PSI st-dev problem!')
            else if(zscu(ji,jk,jj).le.0.0) then   ! Chi
              write(nulout,*) 'meancvgd: ji,jk,jj,zspp=',
     &              ji,jk,jj,zspp(ji,jk,jj)
              call abort3d(nulout,'meancvgd: CHI_u st-dev problem!')
            else if(zstu(ji,jk,jj).le.0.0) then   ! Tu
              write(nulout,*) 'meancvgd: ji,jk,jj,zstu=',
     &              ji,jk,jj,zstu(ji,jk,jj)
              call abort3d(nulout,'meancvgd: T_u st-dev problem!')
            else if(zslq(ji,jk,jj).le.0.0) then   ! lq
              write(nulout,*) 'meancvgd: ji,jk,jj,zslq=',
     &              ji,jk,jj,zslq(ji,jk,jj)
              call abort3d(nulout,'meancvgd: Lnq st-dev problem!')
            else if(zstb(ji,jk,jj).le.0.0) then   ! Tb
              write(nulout,*) 'meancvgd: ji,jk,jj,zstb=',
     &             ji,jk,jj,zstb(ji,jk,jj)
              call abort3d(nulout,'meancvgd: T_b st-dev problem!')
            endif
!
            zspp(ji,jk,jj) = sqrt(zspp(ji,jk,jj))
            zscu(ji,jk,jj) = sqrt(zscu(ji,jk,jj))
            zstu(ji,jk,jj) = sqrt(zstu(ji,jk,jj))
            zslq(ji,jk,jj) = sqrt(zslq(ji,jk,jj))
            zstb(ji,jk,jj) = sqrt(zstb(ji,jk,jj))
!
!            zstdpp(ji,jj,jk) = zspp(ji,jk,jj)
!            zstdcu(ji,jj,jk) = zscu(ji,jk,jj)
!            zstdtu(ji,jj,jk) = zstu(ji,jk,jj)
!            zstdlq(ji,jj,jk) = zslq(ji,jk,jj)
!            zstdtb(ji,jj,jk) = zstb(ji,jk,jj)
!
!            zonalpp(jj,jk) = zonalpp(jj,jk)+zstdpp(ji,jj,jk)
!            zonalcu(jj,jk) = zonalcu(jj,jk)+zstdcu(ji,jj,jk)
!            zonaltu(jj,jk) = zonaltu(jj,jk)+zstdtu(ji,jj,jk)
!            zonallq(jj,jk) = zonallq(jj,jk)+zstdlq(ji,jj,jk)
!            zonaltb(jj,jk) = zonaltb(jj,jk)+zstdtb(ji,jj,jk)
          enddo
!          zonalpsu(jj) = zonalpsu(jj)+zstdpsu(ji,jj)
!          zonaltg(jj) = zonaltg(jj)+zstdtg(ji,jj)
!          zonalpb(jj) = zonalpb(jj)+zstdpb(ji,jj)
        enddo
      enddo
!
!      zni=real(ni)
!      do jj=1,nj
!        do jk=1,nflev
!          zonalpp(jj,jk) = zonalpp(jj,jk)/zni
!          zonalcu(jj,jk) = zonalcu(jj,jk)/zni
!          zonaltu(jj,jk) = zonaltu(jj,jk)/zni
!          zonallq(jj,jk) = zonallq(jj,jk)/zni
!          zonaltb(jj,jk) = zonaltb(jj,jk)/zni
!        enddo
!        zonalpsu(jj) = zonalpsu(jj)/zni
!        zonaltg(jj) = zonaltg(jj)/zni
!        zonalpb(jj) = zonalpb(jj)/zni
!      enddo

!     ------------------------------------------
!     Apply spectral filter to st-dev if desired
!     ------------------------------------------
      if(llfiltersdev) then
!
!* Psi 
!
        call gdtruncr(zspp,'T',itrunc,clpart,nflev)
!
!* Chi_u
!
        call gdtruncr(zscu,'T',itrunc,clpart,nflev)
!
!* T_u
!
        call gdtruncr(zstu,'T',itrunc,clpart,nflev)
!
!* lq
!
        call gdtruncr(zslq,'T',itrunc,clpart,nflev)
!
!* T_b
!
        call gdtruncr(zstb,'T',itrunc,clpart,nflev)
!
!* ps_u
!
        call gdtruncr(zstdpsu,'T',itrunc,clpart,1)
!
!* TG
!
        call gdtruncr(zstdtg,'T',itrunc,clpart,1)
!
!* ps_b
!
        call gdtruncr(zstdpb,'T',itrunc,clpart,1)

      endif
!
!     -----------
!*6.  OUTPUT MEAN
!     -----------

      IPAK = -32
      IDATYP = 5
      IP1 = 0
      IP2 = 0
      IP3 = iensemble
      IDATEO = NDATESTAT
!
      cletik = 'CVGDMEAN'
!
!      write(nulout,*) 'meancvgd: print_mean'
!      write(nulout,*) 'meancvgd: PP MEAN'
!      call maxmin(zmeanpp,ni,nflev,nj,zmin,zmax,
!     &            idum1,idum2,idum3,idum4,'meancvgd    ',
!     &            'PPM')
!      write(nulout,*) 'meancvgd: CU MEAN'
!      call maxmin(zmeancu,ni,nflev,nj,zmin,zmax,
!     &            idum1,idum2,idum3,idum4,'meancvgd    ',
!     &            'CUM')
!      write(nulout,*) 'meancvgd: TU MEAN'
!      call maxmin(zmeantu,ni,nflev,nj,zmin,zmax,
!     &            idum1,idum2,idum3,idum4,'meancvgd    ',
!     &            'TUM')
!      write(nulout,*) 'meancvgd: LQ MEAN'
!      call maxmin(zmeanlq,ni,nflev,nj,zmin,zmax,
!     &            idum1,idum2,idum3,idum4,'meancvgd    ',
!     &            'LQM')
!      write(nulout,*) 'meancvgd: TB MEAN'
!      call maxmin(zmeantb,ni,nflev,nj,zmin,zmax,
!     &            idum1,idum2,idum3,idum4,'meancvgd    ',
!     &            'TBM')
!
      do jvar=1,inbvar
        do jk=1,nflev
          IP1      =  NIP1(jk)
          if(clvar(jvar).eq.'UU') then
            do jj=1,nj
            do ji=1,ni
              z2d(ji,jj) = zmpp(ji,jk,jj)
            enddo
            enddo
            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
     &            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),cletik,
     &            'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
          else if(clvar(jvar).eq.'VV') then
            do jj=1,nj
            do ji=1,ni
              z2d(ji,jj) = zmcu(ji,jk,jj)
            enddo
            enddo
            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
     &            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),cletik,
     &            'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
          else if((clvar(jvar).eq.'PP').or.(clvar(jvar).eq.'QQ')) then
            do jj=1,nj
            do ji=1,ni
              z2d(ji,jj) = zmpp(ji,jk,jj)
            enddo
            enddo
            call innerfld(z2din,z2d,ni,nj,mni_in,mnj_in)
            IERR = VFSTECR(z2din,z2din,IPAK,koutfile,IDATEO,
     &           0,0,mni_in,mnj_in,1,IP1,IP2,IP3,'E',clvar(jvar),'MEANGDIN',
     &                     'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
!
            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
     &            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),cletik,
     &            'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
          else if((clvar(jvar).eq.'CC').or.(clvar(jvar).eq.'DD')) then
            do jj=1,nj
            do ji=1,ni
              z2d(ji,jj) = zmcu(ji,jk,jj)
            enddo
            enddo
            call innerfld(z2din,z2d,ni,nj,mni_in,mnj_in)
            IERR = VFSTECR(z2din,z2din,IPAK,koutfile,IDATEO,
     &             0,0,mni_in,mnj_in,1,IP1,IP2,IP3,'E',clvar(jvar),'MEANGDIN',
     &                     'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)

            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
     &            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),cletik,
     &            'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
          else if(clvar(jvar).eq.'QQ') then
            do jj=1,nj
            do ji=1,ni
              z2d(ji,jj) = zmpp(ji,jk,jj)
            enddo
            enddo
            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
     &            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),cletik,
     &            'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
          else if(clvar(jvar).eq.'DD') then
            do jj=1,nj
            do ji=1,ni
              z2d(ji,jj) = zmcu(ji,jk,jj)
            enddo
            enddo
            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
     $            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),cletik,
     &            'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
          else if(clvar(jvar).eq.'TT') then
            do jj=1,nj
            do ji=1,ni
              z2d(ji,jj) = zmtu(ji,jk,jj)
            enddo
            enddo
            call innerfld(z2din,z2d,ni,nj,mni_in,mnj_in)
            IERR = VFSTECR(z2din,z2din,IPAK,koutfile,IDATEO,
     &            0,0,mni_in,mnj_in,1,IP1,IP2,IP3,'E',clvar(jvar),'MEANGDIN',
     &                     'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
!
            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
     $            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),cletik,
     &            'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
          else if(clvar(jvar).eq.'TB') then
            do jj=1,nj
            do ji=1,ni
              z2d(ji,jj) = zxmtb(ji,jk,jj)
            enddo
            enddo
            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
     $            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),
     $            cletik,'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
!            IERR = VFSTECR(zsavet0(1,1,jk),zsavet0(1,1,jk),IPAK,koutfile,
!     $            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),
!     $           'TT0     ','Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
!            IERR = VFSTECR(zsavetb(1,1,jk),zsavetb(1,1,jk),IPAK,koutfile,
!     $            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),
!     $           'TTB     ','Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
!            IERR = VFSTECR(zsavepsi(1,1,jk),zsavepsi(1,1,jk),IPAK,koutfile,
!     $            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E','PP',
!     $           'ERR-LU  ','Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
!            IERR = VFSTECR(zsavevort(1,1,jk),zsavevort(1,1,jk),IPAK,koutfile,
!     $            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E','QQ',
!     $           'ERR-LU  ','Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
          else if(clvar(jvar).eq.'LQ') then
            do jj=1,nj
            do ji=1,ni
              z2d(ji,jj) = zmlq(ji,jk,jj)
            enddo
            enddo
            call innerfld(z2din,z2d,ni,nj,mni_in,mnj_in)
            IERR = VFSTECR(z2din,z2din,IPAK,koutfile,IDATEO,
     &          0,0,mni_in,mnj_in,1,IP1,IP2,IP3,'E',clvar(jvar),'MEANGDIN',
     &                     'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
!
            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
     $            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),cletik,
     &            'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
          endif
        enddo
      enddo

!     Convert Pa to hPa
      do jj=1,nj
        do ji=1,ni
          zmeanpsu(ji,jj) = zmeanpsu(ji,jj) * 0.01
          zmeanpb(ji,jj) = zmeanpb(ji,jj) * 0.01
        enddo
      enddo
!
      IP1=0
      IP2=0
!
      write(nulout,*) 'meancvgd: Pb MEAN'
      call maxmin(zmeanpb,ni,nj,1,zmin,zmax,
     &            idum1,idum2,idum3,idum4,'meancvgd    ',
     &            'Pbm')
      IERR = VFSTECR(zmeanpb,zmeanpb,IPAK,koutfile,IDATEO,0,0,ni,nj,1,
     $     IP1,IP2,IP3,'E','PB',cletik,'Z',mig1flda,mig2flda,mig3flda,
     $     0,IDATYP,.TRUE.)
      write(nulout,*) 'meancvgd: Pu MEAN'
      call maxmin(zmeanpsu,ni,nj,1,zmin,zmax,
     &            idum1,idum2,idum3,idum4,'meancvgd    ',
     &            'Pum')
      call innerfld(z2din,zmeanpsu,ni,nj,mni_in,mnj_in) ! on inner grid for visualization
      IERR = VFSTECR(z2din,z2din,IPAK,koutfile,IDATEO,0,0,mni_in,mnj_in,
     &                   1  ,IP1,IP2,IP3,'E','PU','MEANGDIN',
     &                   'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
!
      IERR = VFSTECR(zmeanpsu,zmeanpsu,IPAK,koutfile,IDATEO,0,0,ni,nj,1,
     $     IP1,IP2,IP3,'E','PU',cletik,'Z',mig1flda,mig2flda,mig3flda,
     $     0,IDATYP,.TRUE.)
      write(nulout,*) 'meancvgd: Tg MEAN'
      call maxmin(zmeantg,ni,nj,1,zmin,zmax,
     &            idum1,idum2,idum3,idum4,'meancvgd    ',
     &            'Tgm')
      call innerfld(z2din,zmeantg,ni,nj,mni_in,mnj_in) ! on inner grid for visualization
      IERR = VFSTECR(z2din,z2din,IPAK,koutfile,IDATEO,0,0,mni_in,mnj_in,
     &                   1  ,IP1,IP2,IP3,'E','TG','MEANGDIN',
     &                   'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
      IERR = VFSTECR(zmeantg,zmeantg,IPAK,koutfile,IDATEO,0,0,ni,nj,1,
     $     IP1,IP2,IP3,'E','TG',cletik,'Z',mig1flda,mig2flda,mig3flda,
     $     0,IDATYP,.TRUE.)
!
      IERR = VFSTECR(zpsb,zpsb,IPAK,koutfile,IDATEO,0,0,ni,nj,1,
     $     IP1,IP2,IP3,'E','P0','PSB     ','Z',mig1flda,mig2flda,mig3flda,
     $     0,IDATYP,.TRUE.)
!
!     ---------------------------------------------
!*7.  OUTPUT STD DEV 
!     ---------------------------------------------

      IPAK = -32
      IDATYP = 5
      IP1 = 0
      IP2 = 0
      IP3 = iensemble
      IDATEO = NDATESTAT
!
      cletik = 'CVGDSDEV'
!
      write(nulout,*) 'meancvgd: print_stdev'
      write(nulout,*) 'meancvgd: PP ST-DEV'
!      call maxmin(zstdpp,ni,nflev,nj,zmin,zmax,  ! attention aux dimensions ici....
!     &            idum1,idum2,idum3,idum4,'meancvgd    ',
!     &            'PP ')
!      write(nulout,*) 'meancvgd: CU ST-DEV'
!      call maxmin(zstdcu,ni,nflev,nj,zmin,zmax,  ! attention aux dimensions ici....
!     &            idum1,idum2,idum3,idum4,'meancvgd    ',
!     &            'CU ')
!      write(nulout,*) 'meancvgd: TU ST-DEV'
!      call maxmin(zstdtu,ni,nflev,nj,zmin,zmax,  ! attention aux dimensions ici....
!     &            idum1,idum2,idum3,idum4,'meancvgd    ',
!     &            'TU ')
!      write(nulout,*) 'meancvgd: LQ ST-DEV'
!      call maxmin(zstdlq,ni,nflev,nj,zmin,zmax,  ! attention aux dimensions ici....
!     &            idum1,idum2,idum3,idum4,'meancvgd    ',
!     &            'LQ ')
!      write(nulout,*) 'meancvgd: TB ST-DEV'
!      call maxmin(zstdtb,ni,nflev,nj,zmin,zmax,  ! attention aux dimensions ici....
!     &            idum1,idum2,idum3,idum4,'meancvgd    ',
!     &            'TB ')
!
      do jvar=1,inbvar
        do jk=1,nflev
          IP1      =  NIP1(jk)
          if(clvar(jvar).eq.'UU') then
            do jj=1,nj
            do ji=1,ni
              z2d(ji,jj) = zspp(ji,jk,jj)
            enddo
            enddo
            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
     $            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),cletik,
     &            'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
          else if(clvar(jvar).eq.'VV') then
            do jj=1,nj
            do ji=1,ni
              z2d(ji,jj) = zscu(ji,jk,jj)
            enddo
            enddo
            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
     $            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),cletik,
     &            'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
          else if((clvar(jvar).eq.'PP').or.(clvar(jvar).eq.'QQ')) then
!
            do jj=1,nj
            do ji=1,ni
              z2d(ji,jj) = zspp(ji,jk,jj)
            enddo
            enddo
             call innerfld(z2din,z2d,ni,nj,mni_in,mnj_in)
             IERR = VFSTECR(z2din,z2din,IPAK,koutfile,IDATEO,
     &        0,0,mni_in,mnj_in,1,IP1,IP2,IP3,'E',clvar(jvar),'STDDEVIN',
     &       'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
!
            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
     $            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),cletik,
     &            'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
          else if((clvar(jvar).eq.'CC').or.(clvar(jvar).eq.'DD')) then
            do jj=1,nj
            do ji=1,ni
              z2d(ji,jj) = zscu(ji,jk,jj)
            enddo
            enddo
            call innerfld(z2din,z2d,ni,nj,mni_in,mnj_in)
            IERR = VFSTECR(z2din,z2din,IPAK,koutfile,IDATEO,
     &         0,0,mni_in,mnj_in,1,IP1,IP2,IP3,'E',clvar(jvar),'STDDEVIN',
     &                     'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
!
            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
     $            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),cletik,
     &            'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
          else if(clvar(jvar).eq.'TT') then
c     Tu
            do jj=1,nj
            do ji=1,ni
              z2d(ji,jj) = zstu(ji,jk,jj)
            enddo
            enddo
            call innerfld(z2din,z2d,ni,nj,mni_in,mnj_in)
            IERR = VFSTECR(z2din,z2din,IPAK,koutfile,IDATEO,
     &             0,0,mni_in,mnj_in,1,IP1,IP2,IP3,'E',clvar(jvar),'STDDEVIN',
     &                     'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
!
            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
     $            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),cletik,
     &            'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
            do ji=1,ni
               do jj=1,nj
                     z2d(ji,jj) = z2d(ji,jj)*z2d(ji,jj)
               enddo
            enddo
            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
     $            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),
     $           'TTU_VAR ','Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
          else if(clvar(jvar).eq.'TB') then
c     Tb
            do jj=1,nj
            do ji=1,ni
              z2d(ji,jj) = zstb(ji,jk,jj)
            enddo
            enddo
            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
     $            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),
     $            cletik,'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
            do ji=1,ni
               do jj=1,nj
                     z2d(ji,jj) = z2d(ji,jj)*z2d(ji,jj)
               enddo
            enddo
            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
     $            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),
     $           'TTB_VAR ','Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
c     TT
!            do ji=1,ni
!               do jj=1,nj
!                     z2d(ji,jj) = zstt(ji,jk,jj)*zstt(ji,jk,jj)
!               enddo
!            enddo
!            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
!     $            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),
!     $           'TTT_VAR ','Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
          else if(clvar(jvar).eq.'LQ') then
            do jj=1,nj
            do ji=1,ni
              z2d(ji,jj) = zslq(ji,jk,jj)
            enddo
            enddo
            call innerfld(z2din,z2d,ni,nj,mni_in,mnj_in)
            IERR = VFSTECR(z2din,z2din,IPAK,koutfile,IDATEO,
     &            0,0,mni_in,mnj_in,1,IP1,IP2,IP3,'E',clvar(jvar),'STDDEVIN',
     &                     'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
!
            IERR = VFSTECR(z2d,z2d,IPAK,koutfile,
     $            IDATEO,0,0,ni,nj,1,IP1,IP2,IP3,'E',clvar(jvar),cletik,
     &            'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
          endif
        enddo
      enddo

!     Convert Pa to hPa
      do jj=1,nj
        do ji=1,ni
          zstdpsu(ji,jj) = zstdpsu(ji,jj) * 0.01
          zstdpb(ji,jj) = zstdpb(ji,jj) * 0.01
        enddo
      enddo
      write(nulout,*) 'meancvgd: Pb ST-DEV'
      call maxmin(zstdpb,ni,nj,1,zmin,zmax,
     &            idum1,idum2,idum3,idum4,'meancvgd    ',
     &            'Pb')
!
      IP1=0
      IP2=0
      call innerfld(z2din,zstdpsu,ni,nj,mni_in,mnj_in)
      IERR = VFSTECR(z2din,z2din,IPAK,koutfile,IDATEO,0,0,mni_in,mnj_in,
     &               1  ,IP1,IP2,IP3,'E','PU','STDDEVIN',
     &               'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
!
      IERR = VFSTECR(zstdpsu,zstdpsu,IPAK,koutfile,IDATEO,0,0,ni,nj,1,
     $     IP1,IP2,IP3,'E','PU',cletik,'Z',mig1flda,mig2flda,mig3flda,
     $     0,IDATYP,.TRUE.)
      IERR = VFSTECR(zstdtg,zstdtg,IPAK,koutfile,IDATEO,0,0,ni,nj,1,
     $     IP1,IP2,IP3,'E','TG',cletik,'Z',mig1flda,mig2flda,mig3flda,
     $     0,IDATYP,.TRUE.)
      IERR = VFSTECR(zstdpb,zstdpb,IPAK,koutfile,IDATEO,0,0,ni,nj,1,
     $     IP1,IP2,IP3,'E','PB',cletik,'Z',mig1flda,mig2flda,mig3flda,
     $     0,IDATYP,.TRUE.)

!
!     Zonal
!

!      do jj=1,nj
!        zonalpsu(jj) = zonalpsu(jj) * 0.01
!        zonalpb(jj) = zonalpb(jj) * 0.01
!      enddo
!
!      IERR = VFSTECR(zonalpsu,zonalpsu,IPAK,koutfile,IDATEO,0,0,1,nj,
!     &               1  ,IP1,IP2,IP3,'E','P0','PU-ZONAL',
!     &               'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!      IERR = VFSTECR(zonalpb,zonalpb,IPAK,koutfile,IDATEO,0,0,1,nj,
!     &               1  ,IP1,IP2,IP3,'E','P0','PB-ZONAL',
!     &               'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!      if(ldhelm) then
!        clnomvar = 'PP'
!      else
!        clnomvar = 'QQ'
!      endif
!      IERR = VFSTECR(zonalpp,zonalpp,IPAK,koutfile,IDATEO,
!     &               0,0,nj,nflev,1,IP1,IP2,IP3,'E',clnomvar,'PP-ZONAL',
!     &               'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!      if(ldhelm) then
!        clnomvar = 'CC'
!      else
!        clnomvar = 'DD'
!      endif
!      IERR = VFSTECR(zonalcu,zonalcu,IPAK,koutfile,IDATEO,
!     &               0,0,nj,nflev,1,IP1,IP2,IP3,'E',clnomvar,'CU-ZONAL',
!     &               'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!      IERR = VFSTECR(zonaltu,zonaltu,IPAK,koutfile,IDATEO,
!     &               0,0,nj,nflev,1,IP1,IP2,IP3,'E','TT','TU-ZONAL',
!     &               'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!      IERR = VFSTECR(zonaltb,zonaltb,IPAK,koutfile,IDATEO,
!     &               0,0,nj,nflev,1,IP1,IP2,IP3,'E','TT','TB-ZONAL',
!     &               'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!      IERR = VFSTECR(zonallq,zonallq,IPAK,koutfile,IDATEO,
!     &               0,0,nj,nflev,1,IP1,IP2,IP3,'E','LQ','LQ-ZONAL',
!     &               'Z',mig1flda,mig2flda,mig3flda,0,IDATYP,.TRUE.)
!
!     -----------------------------
!     Writing positional parameters
!     -----------------------------

      cletiket = ' '
      cltypvar = ' '
!
      do ji=1,ni
        zx8(ji)=grd_x_8(ji)
      enddo
      do jj=1,nj
        zy8(jj)=grd_y_8(jj)
      enddo
!
      ierr = vfstecr(zx8,zx8,ipak,koutfile,idateo,
     &          0,0,ni,1,1,mig1flda,mig2flda,mig3flda,cltypvar
     &          ,'>>',cletiket,'E',mig1tic,mig2tic,mig3tic,mig4tic,idatyp
     &          ,.true.)
!
      ierr = vfstecr(zy8, zy8, ipak, koutfile, idateo,
     &          0,0, 1, nj, 1, mig1flda,mig2flda,mig3flda,cltypvar
     &          ,'^^',cletiket,'E',mig1tic,mig2tic,mig3tic,mig4tic,idatyp
     &          ,.true.)
!
!     2nd: Tictic/tactac of the inner(GEM) domain (i.e. excluding extension zones)
!
      ierr = vfstecr(zx8,zx8,ipak,koutfile,idateo,
     &          0,0,mni_in,1,1,mig1in,mig2in,mig3in,cltypvar
     &          ,'>>',cletiket,'E',mig1tic,mig2tic,mig3tic,mig4tic,idatyp
     &          ,.true.)
!
      ierr = vfstecr(zy8, zy8, ipak,koutfile,idateo,
     &          0,0, 1, mnj_in, 1, mig1in,mig2in,mig3in,cltypvar
     &          ,'^^',cletiket,'E',mig1tic,mig2tic,mig3tic,mig4tic,idatyp
     &          ,.true.)
!
      write(nulout,*) 'meancvgd: END '
      nensemble = iensemble
!
      return
      end