!-------------------------------------- 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 glb2lam_1(kulstat,kulstdev) 1,58
#if defined (DOC)
*
***s/r glb2lam_1:  Process files of forecast error samples and generate/add target LAM interpolated fields.
*
*Author:  L. Fillion *ARMA/EC 09 Oct 2008
*Revision: L. Fillion *ARMA/EC 26 Nov 2008 - Extend to output LAM fields PP,CC,TT,GZ,LQ,P0 ready for LAM stats.
*  L. Fillion *ARMA/EC 18 Dec 2008 - Introduce resolution dependent horizontal interpolation via cldegint.
*
*Purpose: To feed LAM4D background error statistics in terms of global estimates of PSI/CHI and other regular error samples.
*         This avoids the problem of biperiodicization in a lam4d context which is problematical.
*
*Arguments   KULSTAT  logical unit number
*
#endif
      IMPLICIT NONE
*
*implicits
#include "pardim.cdk"
#include "comdim.cdk"
#include "comlun.cdk"
#include "comct0.cdk"
#include "comcva.cdk"
#include "comcst.cdk"
#include "comgem.cdk"
#include "comsp.cdk"
#include "comgd0.cdk"
#include "comgdpar.cdk"
#include "comcse1.cdk"
#include "comstdd.cdk"
#include "comgrd_param.cdk"
#include "comgrd.cdk"
#include "comgemla.cdk"
#include "comode.cdk"
*
      INTEGER KULSTAT,kulstdev
C
      CHARACTER*1 clflt
      logical llfilt
      integer ifois
      data ifois/0/
      INTEGER JENS, IENS, JK1, IERR, JFILE, ji,jj,JK, jla,JLAT, JLON
      integer idum1,idum2,idum3,idum4
C
      INTEGER FNOM, FSTOUV, FSTFRM, FCLOS, FSTPRM, FSTINL
      INTEGER VFSTECR,ezqkdef
C
C*    RPN Standard files parameters
C
      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
!
      logical llqd,lloutglb,llvfilt
      integer jvar,jlev, itrlgid,ind
      integer iip1s(jpnflev),iip1,iip2,iip3,itrlnlev
      integer ipmode,ipkind,ip1_pak_trl,ip1_vco_trl
      integer :: k,koutmpg  !  the unit which has the selected records.
      real    zlev(jpnflev)
      character*1 clstring
      character*(6) cldegint
      REAL*8 DHEURES
      CHARACTER*1 CLTYPVAR,CLGRTYP
      CHARACTER*2 CLNOMVAR,clmoist
      CHARACTER*8 CLETIKET
C
      REAL*8 DLA2
      REAL*8 ZFACT
      INTEGER IPAK, IDATEO, IKULFILE, iflag,ikount
      CHARACTER*128 CLFLFILE
!
      integer igridid_glb,itargetid,ezgdef_fmem
!
      real ax(mni_in),ay(mnj_in)
      real*8 zglb_dx,zlam_dx
      real*8 zmin,zmax
      real*8 zwork
      real*8 zvortglb(ni,nj,nflev)
      real*8 zdivglb(ni,nj,nflev)
      real*8 zpsiglb(ni,nj,nflev)
      real*8 zchiglb(ni,nj,nflev)
      real*8 zttglb(ni,nj,nflev)
      real*8 zgzglb(ni,nj,nflev)
      real*8 zlqglb(ni,nj,nflev)
      real*8 zpsglb(ni,nj)
!
      real*8 z3d(nila,nflev,njla)
      real*8 z3d2(nila,nflev,njla)
      real*8 z2d(nila,njla)
!
      real*8 zvortlam(mni_in,mnj_in,nflev)
      real*8 zdivlam(mni_in,mnj_in,nflev)
      real*8 zpsilam(mni_in,mnj_in,nflev)
      real*8 zchilam(mni_in,mnj_in,nflev)
      real*8 zttlam(mni_in,mnj_in,nflev)
      real*8 zgzlam(mni_in,mnj_in,nflev)
      real*8 zgzlam2(mni_in,mnj_in,nflev)
      real*8 zlqlam(mni_in,mnj_in,nflev)
      real*8 zpslam(mni_in,mnj_in)
*
      real*8 ztg(ni,nflev,nj), zqg(ni,nflev,nj), zesg(ni,nflev,nj)
      real*8 zug(ni,nflev,nj), zvg(ni,nflev,nj), zgzg(ni,nflev,nj)
      real*8 zpsg(ni,nj)
!
      real*8 zwrksp(nla,2,nflev)
!
!!
      if(grd_typ.ne.'GU'.or.multi_grd.lt.1) then
        call abort3d(nulout,'glb2lam_1: grd_typ, multi_grd bad')
      endif
!
      llfilt = .false.
      if(lflt_low) then
        llfilt = .true.
        clflt = 'L'
      else if(lflt_high) then
        llfilt = .true.
        clflt = 'H'
      endif
!
      llqd = .true.
      lloutglb = .false.
      DLA2 = DBLE(RA) * DBLE(RA)
      IKULFILE = kulstat
!
!     Check which moisture input filed is used
!     
      NFSTVAR = 0
      CFSTVAR(1) = ' '
      CFSTVAR(2) = ' '
      CFSTVAR(3) = ' '
      CFSTVAR(4) = ' '
      CFSTVAR(5) = ' '
      CFSTVAR(6) = ' '
      CFSTVAR(7) = ' '
!
      CALL READNML('NAMGDPAR',IFLAG)
 
      ikount = 0
      do jvar = 1, nfstvar
        if(cfstvar(jvar).eq.'HU'.or.cfstvar(jvar).eq.'LQ') then
          clmoist = cfstvar(jvar)
          ikount = ikount +1
        endif
      enddo
      if(ikount.gt.1) then
        CALL ABORT3D(NULOUT,'glb2lam_1: moisture cfstvar...Only one allowed!')
      else
        write(nulout,*) 'glb2lam_1: clmoist = ',clmoist
      endif
!
! 1.  Set target analysis grid onto which interpolation is done
!     ---------------------------------------------------------
!
      do ji=1,mni_in
        ax(ji)=grd_x_8(ji)
      enddo
      do jj=1,mnj_in
        ay(jj)=grd_y_8(jj)
      enddo
!
      itargetid= ezgdef_fmem(mni_in,mnj_in,'Z','E',mig1tic, mig2tic,mig3tic,mig4tic,
     &                       ax,ay)  ! tic tac same as extended grid
!
!*    2. Access the increments from a set of files
!     .  (loop on the files)
!
      IDIMAX = 100
      DO 201 JFILE = 1, NFLSTAT
!
         CALL GETINCR(KULSTAT,JFILE)
!
!*    .  2.1 Find how many cases there are to be treated
!
         IP1 = -1
         IP2 = -1
         IP3 = -1
         CLNOMVAR = 'P0'
         IP1 =0
!
         write(NULOUT,*)
         IERR = FSTINL (KULSTAT,INI,INJ,INK
     &        ,-1,CETIKETN,IP1,IP2,IP3,' '
     &        ,CLNOMVAR,ILISTE,INFON,IDIMAX)
         WRITE(NULOUT,9210)INFON
 9210    FORMAT(//,4X,"Ensemble of ",I4," increments")
         IF(INFON.EQ.0) THEN
            WRITE(NULOUT,*)' THIS FILE IS EMPTY. CHECK THE SELECTION CRITERIA'
            CALL ABORT3D(NULOUT,'glb2lam_1: 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)
!
            DHEURES = DBLE(INPAS*IDEET/3600)
!
            CALL INCDATR(IDATV(JENS),IDATE(JENS),SNGL(DHEURES))
            CALL 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(NENSEMBLE.EQ.0) THEN
            NDATESTAT = IDATE(1)
         END IF
!
         CTYPVARN = ' '
         CETIKETN = CLETIKET
!
         write(nulout,*) 'glb2lam_1: INI,INJ = ',INI,INJ
         if(ini.ne.ni.or.inj.ne.nj) then
           call abort3d(nulout,'glb2lam_1: ini.ne.ni.or.inj.ne.nj')
         endif
!
         igridid_glb = ezqkdef(INI, INJ, 'G', 0,0,0,0,0)
!
!        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
!
            CALL NEWDATE(IDATV(JENS),IFSTRUN,IHH,-3)
            WRITE(NULOUT,9310)JENS, IFSTRUN,IHH
 9310       FORMAT(///,5X,"--- Case No. ",I3,5x,"Date and time: ",I10,5x
     &           ,I8)
            NSTAMPN = IDATE(JENS)
!            if(lmcstats) then
!              call getens(ztg,zgzg,zqg,zug,zvg,zesg,zpsg,kulstat)
!            else
              CALL GETFST(KULSTAT,'G','N',IP3)
!            endif
!
            CALL GDSP
!
            if(llqd) then  ! will produce both QQ,DD and PP,CC onto target LAM grid
              CALL SPEREE(NKSDIM,SP,GD,NLA,NIBEG,NIEND,NJBEG,NJEND,NKSDIM)
              do jk = 1, nflev
                do jj = 1, nj
                  ind = nj-jj+1
                  do ji = 1, ni
                    zvortglb(ji,ind,jk) = ut0(ji,jk,jj)
                    zdivglb(ji,ind,jk) = vt0(ji,jk,jj)
                  enddo
                enddo
              enddo
            endif
!
!           CONVERT FROM VORT/DIV TO PSI/CHI
!
            DO JK = 1, NFLEV
              DO jla = 1, NLA
                SPVOR(jla,1,JK) = SPVOR(jla,1,JK) * DLA2*R1SNP1(jla)
                SPVOR(jla,2,JK) = SPVOR(jla,2,JK) * DLA2*R1SNP1(jla)
                SPDIV(jla,1,JK) = SPDIV(jla,1,JK) * DLA2*R1SNP1(jla)
                SPDIV(jla,2,JK) = SPDIV(jla,2,JK) * DLA2*R1SNP1(jla)
              ENDDO
            ENDDO
!
            CALL SPEREE(NKSDIM,SP,GD,NLA,NIBEG,NIEND,NJBEG,NJEND,NKSDIM)
!
!           Set fields ready for horizontal interpolation 
!           and in proper output units onto RPN standard file
!
            do jk = 1, nflev
              do jj = 1, nj
                ind = nj-jj+1
                do ji = 1, ni
                  zpsiglb(ji,ind,jk) = ut0(ji,jk,jj)
                  zchiglb(ji,ind,jk) = vt0(ji,jk,jj)
                  zttglb(ji,ind,jk) = tt0(ji,jk,jj)
                  zgzglb(ji,ind,jk) = RDECAM*gz0(ji,jk,jj)
                  if(clmoist.eq.'LQ') then
                    zlqglb(ji,ind,jk) = q0(ji,jk,jj)
                  else if(clmoist.eq.'HU') then
                    zlqglb(ji,ind,jk) = dlog(max(q0(ji,jk,jj),1.d-6))
                  endif
                enddo
              enddo
            enddo
!
            do jj = 1, nj
              ind = nj-jj+1
              do ji = 1, ni
                zpsglb(ji,ind) = RPATMB*gps0(ji,1,jj)
              enddo
            enddo
!
            if(ifois.eq.0) then
              if(llqd) then
                call maxmin(zvortglb,ni,nj,nflev,zmin,zmax,
     &                    idum1,idum2,idum3,idum4,'glb2lam_1   ',
     &                    'VOR')
                call maxmin(zdivglb,ni,nj,nflev,zmin,zmax,
     &                    idum1,idum2,idum3,idum4,'glb2lam_1   ',
     &                    'DIV')
              endif
              call maxmin(zpsiglb,ni,nj,nflev,zmin,zmax,
     &                  idum1,idum2,idum3,idum4,'glb2lam_1   ',
     &                  'PSI')
              call maxmin(zchiglb,ni,nj,nflev,zmin,zmax,
     &                  idum1,idum2,idum3,idum4,'glb2lam_1   ',
     &                  'CHI')
              call maxmin(zttglb,ni,nj,nflev,zmin,zmax,
     &                  idum1,idum2,idum3,idum4,'glb2lam_1   ',
     &                  'T0 ')
              call maxmin(zgzglb,ni,nj,nflev,zmin,zmax,
     &                  idum1,idum2,idum3,idum4,'glb2lam_1   ',
     &                  'GZ ')
              call maxmin(zlqglb,ni,nj,nflev,zmin,zmax,
     &                  idum1,idum2,idum3,idum4,'glb2lam_1   ',
     &                  'LQ ')
              call maxmin(zpsglb,ni,nj,1,zmin,zmax,
     &                  idum1,idum2,idum3,idum4,'glb2lam_1   ',
     &                  'P0 ')
              ifois = 1
            endif
!
!           Take a decision on order of interpolation to be used upon source and target grids
!
            zglb_dx = 40.e6/real(ni)  ! resolution at the equator of the Gaussian grid
            zlam_dx = dxlam(1,njla/2) ! resolution at the equator of the computational LAM4D grid.
            cldegint = 'CUBIC '
            if(zglb_dx.lt.zlam_dx) cldegint = 'LINEAR'
            write(nulout,*) 'glb2lam_1: Order of Horizontal Interpolation = ',cldegint
!
            if(llqd) then
              call hintscal(zvortglb,ni*nj,igridid_glb,
     &                      zvortlam,mni_in*mnj_in,itargetid,nflev,cldegint)
              call hintscal(zdivglb,ni*nj,igridid_glb,
     &                      zdivlam,mni_in*mnj_in,itargetid,nflev,cldegint)
            endif
            call hintscal(zpsiglb,ni*nj,igridid_glb,
     &                    zpsilam,mni_in*mnj_in,itargetid,nflev,cldegint)
            call hintscal(zchiglb,ni*nj,igridid_glb,
     &                    zchilam,mni_in*mnj_in,itargetid,nflev,cldegint)
            call hintscal(zttglb,ni*nj,igridid_glb,
     &                    zttlam,mni_in*mnj_in,itargetid,nflev,cldegint)
            call hintscal(zgzglb,ni*nj,igridid_glb,
     &                    zgzlam,mni_in*mnj_in,itargetid,nflev,cldegint)
            call hintscal(zlqglb,ni*nj,igridid_glb,
     &                    zlqlam,mni_in*mnj_in,itargetid,nflev,cldegint)
            call hintscal(zpsglb,ni*nj,igridid_glb,
     &                    zpslam,mni_in*mnj_in,itargetid,1,cldegint)
!
!*3.  Apply spatial filtering
!     -----------------------
!
      if(llfilt) then
        write(nulout,*) 'glb2lam_1: **************************************'
        write(nulout,*) 'glb2lam_1: Error sample is spectrally filtered at wvnb =',mflt_trunc
        write(nulout,*) 'glb2lam_1: Filtering part: ',clflt
        write(nulout,*) 'glb2lam_1: **************************************'
!PP
        z3d(:,:,:) = 0.0
        do jj=1,mnj_in
        do jk=1,nflev
        do ji=1,mni_in
          z3d(ji,jk,jj) = zpsilam(ji,jj,jk)
        enddo
        enddo
        enddo
        call initgdla1(z3d,nila,njla,nflev)
        call gdtruncr(z3d,zwrksp,'T',mflt_trunc,clflt,.false.,nflev)
        do jj=1,mnj_in
        do jk=1,nflev
        do ji=1,mni_in
          zpsilam(ji,jj,jk)=z3d(ji,jk,jj)
        enddo
        enddo
        enddo
!CC
        z3d(:,:,:) = 0.0
        do jj=1,mnj_in
        do jk=1,nflev
        do ji=1,mni_in
          z3d(ji,jk,jj) = zchilam(ji,jj,jk)
        enddo
        enddo
        enddo
        call initgdla1(z3d,nila,njla,nflev)
        call gdtruncr(z3d,zwrksp,'T',mflt_trunc,clflt,.false.,nflev)
        do jj=1,mnj_in
        do jk=1,nflev
        do ji=1,mni_in
          zchilam(ji,jj,jk)=z3d(ji,jk,jj)
        enddo
        enddo
        enddo
!TT
        z3d(:,:,:) = 0.0
        do jj=1,mnj_in
        do jk=1,nflev
        do ji=1,mni_in
          z3d(ji,jk,jj) = zttlam(ji,jj,jk)
        enddo
        enddo
        enddo
        call initgdla1(z3d,nila,njla,nflev)
        call gdtruncr(z3d,zwrksp,'T',mflt_trunc,clflt,.false.,nflev)
        do jj=1,mnj_in
        do jk=1,nflev
        do ji=1,mni_in
          zttlam(ji,jj,jk)=z3d(ji,jk,jj)
        enddo
        enddo
        enddo
!Lnq
        z3d(:,:,:) = 0.0
        do jj=1,mnj_in
        do jk=1,nflev
        do ji=1,mni_in
          z3d(ji,jk,jj) = zlqlam(ji,jj,jk)
        enddo
        enddo
        enddo
        call initgdla1(z3d,nila,njla,nflev)
        call gdtruncr(z3d,zwrksp,'T',mflt_trunc,clflt,.false.,nflev)
        do jj=1,mnj_in
        do jk=1,nflev
        do ji=1,mni_in
          zlqlam(ji,jj,jk)=z3d(ji,jk,jj)
        enddo
        enddo
        enddo
!P0
        z2d(:,:) = 0.0
        do jj=1,mnj_in
        do ji=1,mni_in
          z2d(ji,jj) = zpslam(ji,jj)
        enddo
        enddo
        call initgdla1(z2d,nila,njla,1)
        call gdtruncr(z2d,zwrksp,'T',mflt_trunc,clflt,.false.,1)
        do jj=1,mnj_in
        do ji=1,mni_in
          zpslam(ji,jj)=z2d(ji,jj)
        enddo
        enddo
      else
        write(nulout,*) 'glb2lam_1: NO FILTER applied to Error sample'
      endif
!
!*3.  Build Generalized Geopotential P
!     --------------------------------
!
! Build GZ first
!
        z3d(:,:,:) = 0.0
        z3d2(:,:,:) = 0.0
        write(nulout,*) 'glb2lam_1: mni_in,mnj_in,nila,njla=',
     & mni_in,mnj_in,nila,njla
!
        do jj=1,mnj_in
        do jlev=1,nflev
        do ji=1,mni_in
          z3d(ji,jlev,jj) = zttlam(ji,jj,jlev)
        enddo
        enddo
        enddo
        call initgdla1(z3d,nila,njla,nflev)
        call maxmin(z3d,nila,nflev,njla,zmin,zmax,
     &            idum1,idum2,idum3,idum4,'glb2lam_1   ',
     &            'T0')
!
        call ltt2phi_inmi(z3d2,z3d,nila,njla,nflev)
        call initgdla1(z3d2,nila,njla,nflev)
        call maxmin(z3d2,nila,nflev,njla,zmin,zmax,
     &            idum1,idum2,idum3,idum4,'glb2lam_1   ',
     &            'GZ')
        do jj=1,mnj_in
        do jlev=1,nflev
        do ji=1,mni_in
          zgzlam2(ji,jj,jlev) = rdecam*z3d2(ji,jlev,jj)
        enddo
        enddo
        enddo
!
! Build ln_ps increment part
!
        z2d(:,:) = 0.0
!
        do jj = 1,mnj_in
          do ji = 1,mni_in
            z2d(ji,jj) = zpslam(ji,jj)/1.e5
!            z2d(ji,jj) = zpslam(ji,jj)/gpsg(ji,1,jj) ! quand on aura disponible le champ de base....
          enddo
        enddo
        call initgdla1(z2d,nila,njla,1)
!
! Build P = GZ + R Tstar * del(lnps): Resuult in zgz array
!
        do jk=1,nflev
          do ji=1,mni_in
            do jj=1,mnj_in
              zgzlam2(ji,jj,jk)=zgzlam2(ji,jj,jk)+rd*tstar*z2d(ji,jj)
            enddo
          enddo
        enddo
!
! Filter vertically and horizontally
!
            if(llvfilt) call vfilt(z3d,5,'L')  ! vertical normal mode filtering

            IPAK = -32
            iDATYP = 5
            ip2 = 0
            ip3 = jfile
!
            if(lloutglb) then
              do jk = 1, nflev
                IP1      =  NIP1(jk)
                if(llqd) then
                  IERR = VFSTECR(zvortglb(1,1,jk),zvortglb(1,1,jk),IPAK,
     &                         IKULFILE,NDATESTAT,0,0,ni,nj,1,
     &                         ip1,ip2,ip3,'E','QQ',CETIKETN,
     &                         'G',0,0,0,0,IDATYP,.TRUE.)
                  IERR = VFSTECR(zdivglb(1,1,jk),zdivglb(1,1,jk),IPAK,
     &                         IKULFILE,NDATESTAT,0,0,ni,nj,1,
     &                         ip1,ip2,ip3,'E','DD',CETIKETN,
     &                         'G',0,0,0,0,IDATYP,.TRUE.)
                endif
                IERR = VFSTECR(zpsiglb(1,1,jk),zwork,IPAK,
     &                       IKULFILE,NDATESTAT,0,0,ni,nj,1,
     &                       ip1,ip2,ip3,'E','PP',CETIKETN,
     &                       'G',0,0,0,0,IDATYP,.TRUE.)
                IERR = VFSTECR(zchiglb(1,1,jk),zwork,IPAK,
     &                       IKULFILE,NDATESTAT,0,0,ni,nj,1,
     &                       ip1,ip2,ip3,'E','CC',CETIKETN,
     &                       'G',0,0,0,0,IDATYP,.TRUE.)
              enddo
            endif
!
            cletiket = cetiketerr
!
            do jk = 1, nflev
              IP1      =  NIP1(jk)
              if(llqd) then
                IERR = VFSTECR(zvortlam(1,1,jk),zwork,IPAK,IKULFILE,NDATESTAT,
     &                 0,0,mni_in,mnj_in,1,IP1,IP2,IP3,'E','QQ',cletiket,
     &                 'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
!
                IERR = VFSTECR(zdivlam(1,1,jk),zwork,IPAK,IKULFILE,NDATESTAT,
     &                 0,0,mni_in,mnj_in,1,IP1,IP2,IP3,'E','DD',cletiket,
     &                 'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
              endif
              IERR = VFSTECR(zpsilam(1,1,jk),zwork,IPAK,IKULFILE,NDATESTAT,
     &               0,0,mni_in,mnj_in,1,IP1,IP2,IP3,'E','PP',cletiket,
     &               'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
!
              IERR = VFSTECR(zchilam(1,1,jk),zwork,IPAK,IKULFILE,NDATESTAT,
     &               0,0,mni_in,mnj_in,1,IP1,IP2,IP3,'E','CC',cletiket,
     &               'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
!
              IERR = VFSTECR(zttlam(1,1,jk),zwork,IPAK,IKULFILE,NDATESTAT,
     &               0,0,mni_in,mnj_in,1,IP1,IP2,IP3,'E','TT',cletiket,
     &               'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
!
              IERR = VFSTECR(zgzlam(1,1,jk),zwork,IPAK,IKULFILE,NDATESTAT,
     &               0,0,mni_in,mnj_in,1,IP1,IP2,IP3,'E','GZ',cletiket,
     &               'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
!
              IERR = VFSTECR(zgzlam2(1,1,jk),zwork,IPAK,IKULFILE,NDATESTAT,
     &               0,0,mni_in,mnj_in,1,IP1,IP2,IP3,'E','GZ','LTT2PHI ',
     &               'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
!
              IERR = VFSTECR(zlqlam(1,1,jk),zwork,IPAK,IKULFILE,NDATESTAT,
     &               0,0,mni_in,mnj_in,1,IP1,IP2,IP3,'E','LQ',cletiket,
     &               'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
            enddo
            IP1      =  NIP1(NFLEV)
            IERR = VFSTECR(zpslam(1,1),zwork,IPAK,IKULFILE,NDATESTAT,
     &             0,0,mni_in,mnj_in,1,IP1,IP2,IP3,'E','P0',cletiket,
     &             'Z',mig1in,mig2in,mig3in,0,IDATYP,.TRUE.)
!
 231     CONTINUE  ! end loop on the ensemble within current file
!
!
!        Writing positional parameters
!        -----------------------------
!
         IPAK = -32
         iDATYP = 5
!
! Write >>:
         cltypvar = 'X'
         clnomvar = '>>'
         clgrtyp = 'E'
         ierr = vfstecr(grd_x_8,zwork,IPAK,IKULFILE,IDATE(1),
     &          ideet,inpas,mni_in,1,1,mig1in,mig2in,mig3in,cltypvar
     &          ,clnomvar,cetikinc,clgrtyp,mig1tic,mig2tic,mig3tic,mig4tic
     &          ,iDATYP,.true.)
!
! Write ^^
         cltypvar = 'X'
         clnomvar = '^^'
         clgrtyp = 'E'
         ierr = vfstecr(grd_y_8,zwork,IPAK,IKULFILE,IDATE(1),
     &          ideet,inpas,1,mnj_in,1,mig1in,mig2in,mig3in,cltypvar
     &          ,clnomvar,cetikinc,clgrtyp,mig1tic,mig2tic,mig3tic,mig4tic
     &          ,iDATYP,.true.)
!
         IERR =  FSTFRM (KULSTAT)
         IERR =  FCLOS  (KULSTAT)
 201  CONTINUE ! end loop on jfile
!
      write(nulout,*) 'glb2lam_1: END'
!
      RETURN
      END