!-------------------------------------- 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 readstd3d 2,13
#if defined (DOC)
*
***s/r READSTD3D - Read in 3D background Std Dev fields
*
*Author  : M. Buehner *ARMA/MSC  March, 2002
*Revision:
*
* Arguments:
*
#endif
      implicit none
*implicits
*
*     Global variables
*
#include "comct0.cdk"
#include "pardim.cdk"
#include "comdim.cdk"
#include "comgem.cdk"
#include "comlun.cdk"
#include "comspg.cdk"
#include "rpnstd.cdk"
#include "comleg.cdk"
#include "comsp.cdk"
#include "compost.cdk"
#include "cominterp.cdk"
#include "comstate.cdk"
#include "compstat.cdk"
#include "comcva.cdk"
*
*     Arguments
*
*     local variables
*
      integer inbrvar3d,inbrvar2d
      parameter(inbrvar3d=4,inbrvar2d=2)
      integer ikey, ilen,jk,jj,ji,jlev,jvar,jn,jm,ila,iulssf,jlon,jlat
      real*8 zgr(ni,nj),zgr2(nj,nflev)
      character*2 clvar3d(inbrvar3d),clvar2d(inbrvar2d)
c
      integer vfstlir,vfstecr
      external vfstlir,vfstecr
c
c      data clvar3d/'PP','UC','UT','LQ'/
c      data clvar2d/'UP','TG'/
*---------------------------------------------------------------------
*
      if(nconf.eq.800.or.nconf.eq.801) then
        clvar3d(1)='PP'
        clvar3d(2)='CC'
        clvar3d(3)='TT'
        clvar3d(4)='**'
        clvar2d(1)='P0'
        clvar2d(2)='**'
      else
        clvar3d(1)='PP'
        clvar3d(2)='CC'
        clvar3d(3)='TT'
        clvar3d(4)='LQ'
        clvar2d(1)='P0'
        clvar2d(2)='TG'
c        clvar3d(1)='PP'
c        clvar3d(2)='UC'
c        clvar3d(3)='UT'
c        clvar3d(4)='LQ'
c        clvar2d(1)='UP'
c        clvar2d(2)='TG'
      endif


*
*     Initializing all the variances fields to zero
*
      write(nulout,*)'Zeroing RGSIG3D'
      ilen = (niend -nibeg +1)*(njend -njbeg +1)*nkgdim
      call zero(ilen,rgsig3d)
*
*     2. Reading the data
*
*     .  2.1 Background error standard deviations
*
      idate(1) = -1
      ip1      = -1
      ip2      = -1
      ip3      = -1
*
      cletiket = 'STDDEV3D'
      cltypvar =' '
***********************************************************************
      write(nulout,*) 'Reading 3D variables'
      do jvar = 1, inbrvar3d
        clnomvar = clvar3d(jvar)
        if(clnomvar.ne.'**') then
        write(nulout,*)'Reading ',clnomvar
c
        do jk = 1,nflev
          ip1 = nip1(jk)
          ikey = vfstlir(zgr,nulbgst,ini,inj,ink,idate(1)
     s         ,cletiket,ip1,ip2,ip3,cltypvar,clnomvar)
*
          if(ikey .lt.0 ) then
            call abort3d(nulout
     &           ,'READSTD3D: Problem with background stat file')
          endif
c
          if(clnomvar.eq.'PP') then
          do jj = 1, nj
            do ji = 1, ni
              rgsiguu3d(ji,jk,nj-jj+1) = zgr(ji,jj)
            enddo
          enddo
          elseif(clnomvar.eq.'UC'.or.clnomvar.eq.'CC') then
          do jj = 1, nj
            do ji = 1, ni
              rgsigvv3d(ji,jk,nj-jj+1) = zgr(ji,jj)
            enddo
          enddo
          elseif(clnomvar.eq.'UT'.or.clnomvar.eq.'TT') then
          do jj = 1, nj
            do ji = 1, ni
              rgsigtt3d(ji,jk,nj-jj+1) = zgr(ji,jj)
            enddo
          enddo
          elseif(clnomvar.eq.'LQ') then
          do jj = 1, nj
            do ji = 1, ni
              rgsigq3d(ji,jk,nj-jj+1) = zgr(ji,jj)
            enddo
          enddo
          endif
c
        enddo
c
c TEMPORARALLY READ IN OLD LQ STATS
c          write(nulout,*) 'Reading in 2D old LQ STDDEV'
c          ikey = vfstlir(zgr2,nulbgst,ini,inj,ink,-1
c     s         ,'STDDEV  ',-1,-1,-1,' ',clnomvar)
c          write(nulout,*) 'IKEY=',ikey
c          do jk = 1, nflev
c            do jj = 1, nj
c              do ji = 1, ni
c                rgsigq3d(ji,jk,nj-jj+1) = zgr2(jj,jk)
c              enddo
c            enddo
c          enddo
c
        endif
      enddo
*
      write(nulout,*) 'Reading 2D variables'
      do jvar = 1, inbrvar2d
c
        clnomvar = clvar2d(jvar)
        if(clnomvar.ne.'**') then
        if((clnomvar.eq.'TG'.and.nsexist(nstg).eq.1).or.
     +      clnomvar.ne.'TG') then
          write(nulout,*)'Reading ',clnomvar
          jk = jvar
          ip1 = 0
          ikey = vfstlir(zgr,nulbgst,ini,inj,ink,-1
     s       ,cletiket,-1,-1,-1,cltypvar,clnomvar)
*
          if(ikey .lt.0 ) then
            call abort3d(nulout
     &           ,'READSTD3D: Problem with background stat file')
          endif
c
          if(clnomvar.eq.'UP'.or.clnomvar.eq.'P0') then
          do jj = 1, nj
            do ji = 1, ni
              rgsigps3d(ji,1,nj-jj+1) = zgr(ji,jj)*1.0D2
            enddo
          enddo
          elseif(clnomvar.eq.'TG') then
          do jj = 1, nj
            do ji = 1, ni
c              rgsigtg3d(ji,1,nj-jj+1) = zgr(ji,jj)
              rgsigtg3d(ji,1,jj) = tgstdbg(ji,jj) ! USE STDDEV SPECIFIED IN SUTG
            enddo
          enddo
          endif
c
        endif
        endif
c
      enddo

cbue temporarily use horizontally constant variances
c        do jj = 1, nj
c          do jk = 1, nkgdim
c            do ji = 1, ni
c              rgsig3d(ji,jk,jj) = rgsig3d(1,jk,1)
c            enddo
c          enddo
c        enddo
cbue
cbue temporarily copy 6th level to levels 1 to 5
c      do jk = 1, 5
c        do jj = 1, nj
c          do ji = 1, ni
c            rgsiguu3d(ji,jk,jj) = rgsiguu3d(ji,6,jj)
c            rgsigvv3d(ji,jk,jj) = rgsigvv3d(ji,6,jj)
c          enddo
c        enddo
c      enddo
cbue

c
c filter the 3d variances fields to smooth out horizontal gradients
c
      if(ntrunc3dstd.gt.0) then
      CALL REESPE(NKSDIM,SP,RGSIG3D
     S             ,NLA,NIBEG,NIEND,NJBEG,NJEND,NKSDIM)
      DO JK=1,NKSDIM
        DO JN=(NTRUNC3DSTD+1),NTRUNC
          DO JM=0,JN
            ILA    = NIND(JM)   +JN-JM
            SP(ILA,1,JK)= 0.0
            SP(ILA,2,JK)= 0.0
          ENDDO
        ENDDO
      ENDDO
c copy level 3 to levels 1 and 2
c      DO JN=0,NTRUNC3DSTD
c        DO JK=1,2
c          DO JM=0,JN
c            ILA    = NIND(JM)   +JN-JM
c            SPVOR(ILA,1,JK)=SPVOR(ILA,1,3)
c            SPVOR(ILA,2,JK)=SPVOR(ILA,2,3)
c            SPDIV(ILA,1,JK)=SPDIV(ILA,1,3)
c            SPDIV(ILA,2,JK)=SPDIV(ILA,2,3)
c            SPTT(ILA,1,JK)=SPTT(ILA,1,3)
c            SPTT(ILA,2,JK)=SPTT(ILA,2,3)
c            SPQ(ILA,1,JK)=SPQ(ILA,1,3)
c            SPQ(ILA,2,JK)=SPQ(ILA,2,3)
c          ENDDO
c        ENDDO
c      ENDDO
      CALL SPEREE(NKSDIM,SP,RGSIG3D
     S             ,NLA,NIBEG,NIEND,NJBEG,NJEND,NKSDIM)
c
c     2. Write out 3D fields to file (in MKS units)
c
      IULSSF=0
      IERR = FNOM(IULSSF,'std_3d_filt.fst','RND',0)
      IERR =  FSTOUV(IULSSF,'RND')
C
CCC   Define a valid date stamp to satisfy fstecr
C                         YYYYMMDD HHMMSShh
      CALL NEWDATE(IDATEO,19991125,12000000,3)
C
      IDEET  = 0
      INPAS  = 0
      INK    = 1
      INJ    = NJ
      INI    = NI
      IP2    = 0
      IP3    = 0
C*    Parameters obtained via (compost)
C*
      CLGRTYP  =  CGRTYP
      IG1      =  NIG1
      IG2      =  NIG2
      IG3      =  NIG3
      IG4      =  NIG4
      IDATYP   =  NIDATYP
      IPAK     =  NPAK
      CLETIKET = 'STDDEV3D'
C
      CLTYPVAR = 'R'
C
      CLNOMVAR = 'LQ'
      do jk=1,NFLEV
        ip1=nip1(jk)
        DO jlon = 1, NI
          DO jlat = 1, NJ
            ZGR(JLON,NJ-JLAT+1) = RGSIGQ3D(JLON,JK,JLAT)
          END DO
        END DO
        IERR = VFSTECR(ZGR,ZGR,IPAK,IULSSF,IDATEO
     S       ,IDEET,IPAS,INI,INJ,INK,IP1,IP2,IP3,CLTYPVAR
     S       ,CLNOMVAR,CLETIKET,CLGRTYP,IG1,IG2,IG3,IG4,IDATYP,.TRUE.)
      enddo
c
      CLNOMVAR = 'PP'
      do jk=1,NFLEV
        ip1=nip1(jk)
        DO jlon = 1, NI
          DO jlat = 1, NJ
            ZGR(JLON,NJ-JLAT+1) = RGSIGUU3D(JLON,JK,JLAT)
          END DO
        END DO
        IERR = VFSTECR(ZGR,ZGR,IPAK,IULSSF,IDATEO
     S       ,IDEET,IPAS,INI,INJ,INK,IP1,IP2,IP3,CLTYPVAR
     S       ,CLNOMVAR,CLETIKET,CLGRTYP,IG1,IG2,IG3,IG4,IDATYP,.TRUE.)
      enddo
c
      CLNOMVAR = 'UC'
      do jk=1,NFLEV
        ip1=nip1(jk)
        DO jlon = 1, NI
          DO jlat = 1, NJ
            ZGR(JLON,NJ-JLAT+1) = RGSIGVV3D(JLON,JK,JLAT)
          END DO
        END DO
        IERR = VFSTECR(ZGR,ZGR,IPAK,IULSSF,IDATEO
     S       ,IDEET,IPAS,INI,INJ,INK,IP1,IP2,IP3,CLTYPVAR
     S       ,CLNOMVAR,CLETIKET,CLGRTYP,IG1,IG2,IG3,IG4,IDATYP,.TRUE.)
      enddo
c
      CLNOMVAR = 'UT'
c      CLNOMVAR = 'TT'
      do jk=1,NFLEV
        ip1=nip1(jk)
        DO jlon = 1, NI
          DO jlat = 1, NJ
            ZGR(JLON,NJ-JLAT+1) = RGSIGTT3D(JLON,JK,JLAT)
          END DO
        END DO
C
        IERR = VFSTECR(ZGR,ZGR,IPAK,IULSSF,IDATEO
     S       ,IDEET,IPAS,INI,INJ,INK,IP1,IP2,IP3,CLTYPVAR
     S       ,CLNOMVAR,CLETIKET,CLGRTYP,IG1,IG2,IG3,IG4,IDATYP,.TRUE.)
      enddo
c
      CLNOMVAR = 'UP'
        ip1=0
        DO jlon = 1, NI
          DO jlat = 1, NJ
            ZGR(JLON,NJ-JLAT+1) = RGSIGPS3D(JLON,1,JLAT)*1.0D-2
          END DO
        END DO
        IERR = VFSTECR(ZGR,ZGR,IPAK,IULSSF,IDATEO
     S       ,IDEET,IPAS,INI,INJ,INK,IP1,IP2,IP3,CLTYPVAR
     S       ,CLNOMVAR,CLETIKET,CLGRTYP,IG1,IG2,IG3,IG4,IDATYP,.TRUE.)
c
      CLNOMVAR = 'TG'
        ip1=0
        DO jlon = 1, NI
          DO jlat = 1, NJ
            ZGR(JLON,NJ-JLAT+1) = RGSIGTG3D(JLON,1,JLAT)
          END DO
        END DO
        IERR = VFSTECR(ZGR,ZGR,IPAK,IULSSF,IDATEO
     S       ,IDEET,IPAS,INI,INJ,INK,IP1,IP2,IP3,CLTYPVAR
     S       ,CLNOMVAR,CLETIKET,CLGRTYP,IG1,IG2,IG3,IG4,IDATYP,.TRUE.)
c
      IERR =  FSTFRM (IULSSF)
      IERR =  FCLOS  (IULSSF)

      endif
*
      WRITE(nulout,*)'DONE in READSTD3D'
*
      RETURN
      END