!-------------------------------------- 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 --------------------------------------
***s/r v4d_getdx - Read increments (or adjoint increments)
*                  from 3D-Var and prepare them for GEM 
*
#include "model_macros_f.h"
*

      subroutine v4d_getdx(kstatus) 5,25
*
      use v4d_prof, only: Pr_mode_S, Pr_llfrm_L, Pr_dsnooze_8, Pr_nsim4d
      use v4dz,     only: V4dzgauss_ni,V4dzgauss_nj
*
      implicit none
*
      integer, intent(inout):: kstatus
*
*author
*     P. Gauthier
*
*revision
* v3_00 - P. Gauthier        - initial MPI version
* v3_00 - M. Tanguay         - add v4d_gauss2gem_ad/Simon's exchange
* v3_01 - Tanguay/Buehner    - introduce gem2gauss for singular vectors
* v3_02 - Tanguay M.         - locate HU in tracers 
* v3_02 - Buehner M.         - read simulation no. and verify consistency
* v3_11 - Tanguay M.         - Remove V4dg_ga_eq_ge_L
* v3_30 - Fillion/Tanguay    - Adapt diagnostics for LAM
*
*object
*     -------------------------
*     If V4dg_di_L or V4dg_tl_L
*     -------------------------
*     1) Proc0: Read increments from 3D-Var
*     2) Proc0: Transfert from Gaussian grid to GEM scalar grid
*     3) All processors: Conversion from 3D-Var units to GEM units and Staggering 
*
*     ------------
*     If V4dg_ad_L
*     ------------
*     1) Proc0: Read adjoint increments from 3D-Var
*     2) Proc0: Adjoint of [Transfert from GEM scalar grid to Gaussian grid]  
*     3) All processors: Adjoint of [Conversion from GEM units and Staggering to 3D-Var units] 
*
*arguments
* Name         I/O                 Description
*----------------------------------------------------------------
* kstatus      I                   Status of the job
*----------------------------------------------------------------
*
*implicits
#include "glb_ld.cdk"
#include "grd.cdk"
#include "lun.cdk"
#include "v4dg.cdk"
#include "vt1.cdk"
#include "vt1m.cdk"
#include "ptopo.cdk"
#include "tr3d.cdk"
#include "path.cdk"
#include <clib_interface.cdk>
#include <prof_f.h>
*
*     Local variables
*     ---------------
      integer idim(2),istat,ihdlin,icount,nvar,jlon,jlat,jlev,i,j,k,
     %        ierr,pnlkey1(9),err,pnerr,nigauss,njgauss,status,ininj,
     %        nsim3d,nstag,inn
*
      real*8,pointer    :: dlbuff_8(:,:),dlbuff2d_8(:)
      real,  allocatable:: zbuff(:,:,:),zbuff2d(:,:)
      real,  allocatable::  gut1(:,:,:),gvt1(:,:,:),gtpt1(:,:,:),
     %                     ghut1(:,:,:),gst1(:,:)
*
      character*256 pathdwgf_S,pathdwga_S
*     
      integer  vmmlod,vmmget,vmmuld,fnom,prof_rdrec
      external vmmlod,vmmget,vmmuld,fnom,prof_rdrec
*
      integer key1(Tr3d_ntr), key1_, key1m(Tr3d_ntr), key1m_, n 
      real hut1, hut1m
      pointer (pahu1, hut1(LDIST_SHAPE,*)), (pahu1m, hut1m(LDIST_SHAPE,*))
*
      logical plpr_L
*
*     ______________________________________________________
*
      if ( V4dg_di_L     ) call gem_stop('v4d_getdx',-1)
*     ______________________________________________________
*
      write(Lun_out,2000)
*
*     Nullify pointers for prof_gvar
*     ------------------------------
      nullify (dlbuff_8,dlbuff2d_8)
*
*     Flag for diagnostics
*     --------------------
      plpr_L=.false.
*
*     Get fields in memory
*     --------------------
      pnlkey1(1)= VMM_KEY( ut1)
      pnlkey1(2)= VMM_KEY( vt1)
      pnlkey1(3)= VMM_KEY( tpt1)
      pnlkey1(4)= VMM_KEY( st1)
      nvar= 4 
*
      err = vmmlod(pnlkey1,nvar)
*
      err = VMM_GET_VAR( ut1)
      err = VMM_GET_VAR( vt1)
      err = VMM_GET_VAR( tpt1)
      err = VMM_GET_VAR( st1)
*
*     Load humidity field
*     -------------------
      key1_ = VMM_KEY (trt1)
      do n=1,Tr3d_ntr
         key1(n) = key1_ + n
      end do
      err = vmmlod(key1,Tr3d_ntr)
      do n=1,Tr3d_ntr
      if (Tr3d_name_S(n).eq.'HU') err = vmmget(key1(n),pahu1,hut1) 
      end do
*     
      kstatus = 0
*
      if(Ptopo_myproc.eq.0) then
*
      if(V4dg_di_L.or.V4dg_tl_L) then
*     ---------------------------
*     Read increments from 3D-Var
*     ---------------------------
*
*     1A.Open dwgf PROF file
*        -------------------
         write(Lun_out,*) 'Opening file dwgf PROF file'
*
         pathdwgf_S = trim(Path_xchg_S)//'/dwgf.prof'
         ihdlin = prof_open(pathdwgf_S,'READ',Pr_mode_S,Pr_dsnooze_8)
*
         if(ihdlin.le.0) then
            write(Lun_out,*) 'Problem opening dwgf PROF file'
            kstatus = -99
            goto 1001
         endif
*
      elseif(V4dg_ad_L) then
*     -----------------------------------
*     Read adjoint increments from 3D-Var
*     -----------------------------------
*
*     1B.Open dwga PROF file
*        -------------------
         write(Lun_out,*) 'Opening file dwga PROF file'
*
         pathdwga_S = trim(Path_xchg_S)//'/dwga1.prof'
         ihdlin = prof_open(pathdwga_S,'READ',Pr_mode_S,Pr_dsnooze_8)
*
         if(ihdlin.le.0) then
            write(Lun_out,*) 'Problem opening dwga PROF file'
            kstatus = -99
            goto 1001
         end if
*
      end if
*
*     2. Get dimensions NI NJ of Gaussian grid 
*        -------------------------------------
         write(Lun_out,*) 'Accessing first record with 3D fields...'
*
         istat = prof_rdrec(ihdlin)
*
         if(istat.ne.0) then
            write(Lun_out,*) 'Problem accessing first record with 3D fields'
            kstatus = -99 
            goto 1001
         endif
*
         write(Lun_out,*) 'Reading dimensions NI NJ of Gaussian grid...'
*
         istat = prof_gvar(ihdlin,ininj,PRM_NINJ)
*
         if(istat.ne.0) then
            write(Lun_out,*) 'Problem reading dimensions NI NJ of Gaussian grid'
            kstatus = -99
            goto 1001
         endif
*
         call mvbits(ininj,16,16,V4dzgauss_ni,0)
         call mvbits(ininj, 0,16,V4dzgauss_nj,0)
*
*        Get simulation number
*        ---------------------
         istat = prof_gvar(ihdlin,nsim3d,PRM_EVNT)
         if(nsim3d.ne.Pr_nsim4d) then
            write(Lun_out,*) 'WRONG SIMULATION NUMBER NSIM3D = ',nsim3d,' NSIM4D = ',Pr_nsim4d
            call gem_stop('v4d_getdx',-1)
         else
            write(Lun_out,*) 'NSIM3D IS THE RIGHT SIMULATION NUMBER = ',nsim3d
         endif
*
*        Get prm_stag parameter 
*        ----------------------
         istat = prof_gvar(ihdlin,nstag,PRM_STAG)
*
*        If nstag = 1, the staggering was already done by 3D-Var
*        If nstag = 0, the staggering     will be done by GEM     
*        -------------------------------------------------------
         if(nstag.ne.0.and.nstag.ne.1) then
*
            write(Lun_out,*) 'PRM_STAG IS WRONG = ',nstag
            kstatus = -99
            goto 1001
*
         endif
*
      end if
*
 1001 call rpn_comm_bcast(kstatus,     1,"MPI_INTEGER",0,"GRID",ierr)
*
      if(kstatus.ne.0) return 
*
      call rpn_comm_bcast(V4dzgauss_ni,1,"MPI_INTEGER",0,"GRID",ierr)
      call rpn_comm_bcast(V4dzgauss_nj,1,"MPI_INTEGER",0,"GRID",ierr)
*
      call rpn_comm_bcast(nstag,1,"MPI_INTEGER",0,"GRID",ierr)
*
*     If nstag = 1, the staggering was already done by 3D-Var
*     If nstag = 0, the staggering     will be done by GEM
*     -------------------------------------------------------
      if(nstag.eq.1) then
*
         write(Lun_out,*) 'WIND STAGGERING ALREADY DONE BY 3D-VAR: PRM_STAG = ',nstag
         V4dg_ustag_L = .true.
         V4dg_vstag_L = .true.
*
      elseif(nstag.eq.0) then
*
         write(Lun_out,*) 'WIND STAGGERING WILL BE DONE BY GEM-DM: PRM_STAG = ',nstag
         V4dg_ustag_L = .false.
         V4dg_vstag_L = .false.
*
      endif
*
*     Store dimensions of 3D-Var Gaussian grid 
*     ----------------------------------------
      nigauss = V4dzgauss_ni 
      njgauss = V4dzgauss_nj 
*
      if(Ptopo_myproc.eq.0) then
         if(.not.allocated(gut1 )) allocate(gut1 (nigauss,njgauss,G_nk), STAT=status)
         if(.not.allocated(gvt1 )) allocate(gvt1 (nigauss,njgauss,G_nk), STAT=status)
         if(.not.allocated(gtpt1)) allocate(gtpt1(nigauss,njgauss,G_nk), STAT=status)
         if(.not.allocated(ghut1)) allocate(ghut1(nigauss,njgauss,G_nk), STAT=status)
         if(.not.allocated(gst1 )) allocate(gst1 (nigauss,njgauss)     , STAT=status)
      endif
*
      if(Ptopo_myproc.eq.0) then
*
*     3. Get all 3D dynamical fields
*        ---------------------------
*
*           U component of winds
*           --------------------
            call v4d_getfld('UU',kstatus)
*
            if(kstatus.ne.0) goto 1002
*
            gut1(:,:,:) = zbuff(:,:,:) 
*
*           V component of winds
*           --------------------
            call v4d_getfld('VV',kstatus)
*
            if(kstatus.ne.0) goto 1002
*     
            gvt1(:,:,:) = zbuff(:,:,:) 
*     
*           Temperature
*           -----------
            call v4d_getfld('TT',kstatus) 
*
            if(kstatus.ne.0) goto 1002
*     
            gtpt1(:,:,:) = zbuff(:,:,:)
*     
*           Specific humidity
*           -----------------
            call v4d_getfld('HU',kstatus) 
*
            if(kstatus.ne.0) goto 1002
*     
            ghut1(:,:,:) = zbuff(:,:,:) 
*
*     4. Get all 2D dynamical fields
*        ---------------------------
*     
*           Read second record with 2D fields
*           --------------------------------- 
            write(Lun_out,*) 'Accessing second record with 2D fields...'
*
            istat = prof_rdrec(ihdlin)
*
            if(istat.ne.0) then
               write(Lun_out,*) 'Problem accessing second record with 2D fields'
               kstatus = -99 
               goto 1002
            endif
*
*           Surface pressure
*           ----------------
            call v4d_getfld('PS',kstatus) 
*
            if(kstatus.ne.0) goto 1002
*     
            gst1(:,:) = zbuff2d(:,:) 
*
            if(associated(dlbuff_8 )) deallocate(dlbuff_8  )
            if(allocated ( zbuff   )) deallocate( zbuff    )
*
            if(associated(dlbuff2d_8))deallocate(dlbuff2d_8)
            if(allocated ( zbuff2d  ))deallocate( zbuff2d  )
*
      endif
*
 1002 call rpn_comm_bcast(kstatus,     1,"MPI_INTEGER",0,"GRID",ierr)
*
      if(kstatus.ne.0) return
*
      if(Ptopo_myproc.eq.0) then
*
      if(V4dg_di_L.or.V4dg_tl_L) then
*
*     5A.Close dwgf PROF file
*        --------------------
         write(Lun_out,*) 'Closing file dwgf PROF file'
*
         istat = prof_close(ihdlin,Pr_llfrm_L)
*
         if(istat.ne.0) then
            write(Lun_out,*) 'Problem closing file dwgf PROF file'
            kstatus = -99
            goto 1003
         endif
*
      elseif(V4dg_ad_L) then
*
*     5B.Close dwga PROF file
*        --------------------
         write(Lun_out,*) 'Closing file dwga PROF file'
*
         istat = prof_close(ihdlin,Pr_llfrm_L)
*
         if(istat.ne.0) then
            write(Lun_out,*) 'Problem closing file dwga PROF file'
            kstatus = -99
            goto 1003
         endif
*
      end if
*
      end if
*
 1003 call rpn_comm_bcast(kstatus,1,"MPI_INTEGER",0,"GRID",ierr)
*
      if(kstatus.ne.0) return 
*
      if(V4dg_di_L.or.V4dg_tl_L) then
*     -----------------------------------------------
*     Transfert from Gaussian grid to GEM scalar grid
*     -----------------------------------------------
      call v4d_gauss2gem( ut1, vt1, tpt1, hut1, st1, LDIST_DIM,
     %                   gut1,gvt1,gtpt1,ghut1,gst1,nigauss,njgauss,G_nk)
*
      elseif(V4dg_ad_L) then
*     ------------------------------------------------------------
*     Adjoint of [Transfert from GEM scalar grid to Gaussian grid] 
*     ------------------------------------------------------------
      call v4d_gem2gauss_ad( ut1, vt1, tpt1, hut1, st1, LDIST_DIM,
     %                      gut1,gvt1,gtpt1,ghut1,gst1,nigauss,njgauss,G_nk)
*
      endif
*
      if(plpr_L) then
         inn= 0
         if (G_lam) then
             inn=1
         endif
         write(Lun_out,*) 'AFTER GAUSS2GEM or AFTER GEM2GAUSS_AD'
         call glbstat(ut1 ,'UU',LDIST_DIM,G_nk,1,G_ni-inn,1,G_nj,  1,G_nk)
         call glbstat(vt1 ,'VV',LDIST_DIM,G_nk,1,G_ni,1,G_nj-1,1,G_nk)
         call glbstat(tpt1,'TP',LDIST_DIM,G_nk,1,G_ni,1,G_nj,  1,G_nk)
         call glbstat(st1 ,'4S',LDIST_DIM,   1,1,G_ni,1,G_nj,  1,   1)
         call glbstat(hut1,'HU',LDIST_DIM,G_nk,1,G_ni,1,G_nj,  1,G_nk)
         write(Lun_out,*) '-----------------------'
      endif
*
      if(Ptopo_myproc.eq.0) then
         if(allocated (gut1      )) deallocate(gut1    )
         if(allocated (gvt1      )) deallocate(gvt1    )
         if(allocated (gtpt1     )) deallocate(gtpt1   )
         if(allocated (ghut1     )) deallocate(ghut1   )
         if(allocated (gst1      )) deallocate(gst1    )
      endif
*
*     Load TRAJECTORY fields needed in VARCONV or VARCONV_AD
*     ------------------------------------------------------
      if(V4dg_tl_L.or.V4dg_ad_L) then
*
         pnlkey1(1)= VMM_KEY(tpt1m)
         pnlkey1(2)= VMM_KEY( st1m)
         nvar= 2
         err = vmmlod(pnlkey1,nvar)
*
         err = VMM_GET_VAR( tpt1m)
         err = VMM_GET_VAR( st1m)
*
*        Load TRAJ humidity field
*        ------------------------
         key1m_ = VMM_KEY (trt1m)
         do n=1,Tr3d_ntr
            key1m(n) = key1m_ + n
         end do
         err = vmmlod(key1m,Tr3d_ntr)
         do n=1,Tr3d_ntr
         if (Tr3d_name_S(n).eq.'HU') err = vmmget(key1m(n),pahu1m,hut1m) 
         end do
*
      endif
*
      if(V4dg_di_L.or.V4dg_tl_L) then
*     --------------------------------------------------------
*     Conversion from 3D-Var units to GEM units and Staggering
*     --------------------------------------------------------
*
*     Direct (nonlinear)
*     ------------------ 
      if(V4dg_di_L) then 
         call v4d_varconv(ut1,vt1,tpt1,hut1,st1,LDIST_DIM,l_nk,.TRUE.)
*     
*     TLM 
*     ---
      elseif(V4dg_tl_L) then
         call v4d_varconv_tl(ut1,vt1,tpt1,hut1,st1,
     $                       tpt1m,hut1m,st1m,LDIST_DIM,l_nk,.TRUE.)
      end if
*
      elseif(V4dg_ad_L) then
*     ---------------------------------------------------------------------
*     Adjoint of [Conversion from GEM units and Staggering to 3D-Var units] 
*     ---------------------------------------------------------------------
*
         call v4d_varconv_ad(ut1,vt1,tpt1,hut1,st1,
     $                       tpt1m,hut1m,st1m,LDIST_DIM,l_nk,.FALSE.)
      end if
*
      if(plpr_L) then
         write(Lun_out,*) 'AFTER VARCONV or VARCONV_AD'
         if(G_lam) then
           call glbstat(ut1 ,'UU',LDIST_DIM,G_nk,1,G_ni-1,1,G_nj,  1,G_nk)
         else
           call glbstat(ut1 ,'UU',LDIST_DIM,G_nk,1,G_ni,  1,G_nj,  1,G_nk)
         endif
         call glbstat(vt1 ,'VV',LDIST_DIM,G_nk,1,G_ni,1,G_nj-1,1,G_nk)
         call glbstat(tpt1,'TP',LDIST_DIM,G_nk,1,G_ni,1,G_nj,  1,G_nk)
         call glbstat(st1 ,'4S',LDIST_DIM,   1,1,G_ni,1,G_nj,  1,   1)
         call glbstat(hut1,'HU',LDIST_DIM,G_nk,1,G_ni,1,G_nj,  1,G_nk)
         write(Lun_out,*) '-----------------------'
      endif
*
      pnerr = vmmuld(-1,0)
*
      write(Lun_out,2001) kstatus
*
 2000 format(/,'V4D_GETDX: Read Model state sent by 3D-Var',
     +       /,'==========================================')
 2001 format(/,'V4D_GETDX: Model state received from 3D-Var --- Status = ',I8,
     +       /,'======================================================')
*
      return
*
*     Host subroutine
*     ---------------
      contains

      subroutine v4d_getfld(cdvar,kstatus) 5
*
      implicit none
*
      character*2, intent(in):: cdvar
      integer, intent(inout) :: kstatus
*
*author
*     P. Gauthier
*
*revision
* v3_00 - P. Gauthier        - initial MPI version
* v3_01 - M. Tanguay         - introduce gem2gauss for singular vectors
* v3_11 - P. Gauthier        - Adjust latitude reversing when V4dg_vstag_L 
* v3_30 - Fillion/Tanguay    - Allow Limited-Area option 
*
*object
*
*arguments
* Name         I/O                 Description
*----------------------------------------------------------------
* cdvar        I                   Type of profile
* kstatus      I                   Status of the job
*----------------------------------------------------------------
*
*implicits
#include "v4dg.cdk"
*
      integer njx
*
*     Read 3D field if request 
*     ------------------------
      select case(cdvar)
      case('UU')
         istat = prof_gvar(ihdlin,dlbuff_8,V3D_UTRU)
      case('VV')
         istat = prof_gvar(ihdlin,dlbuff_8,V3D_VTRU)
      case('TT')
         istat = prof_gvar(ihdlin,dlbuff_8,V3D_TEMP)
      case('HU')
         istat = prof_gvar(ihdlin,dlbuff_8,V3D_SPHU)
*
*     Read 2D field if request 
*     ------------------------
      case('PS')
         istat = prof_gvar(ihdlin,dlbuff2d_8,V2D_PSUR)
      end select
*     
*     Change accuracy and reverse latitude if 3D field   
*     ------------------------------------------------
      select case(cdvar)
      case('UU','VV','TT','HU')
         idim = ubound(dlbuff_8)
         write(Lun_out,*)' - Dimension of ',cdvar,' field 3D ',idim
         if(istat.ne.0) then
            write(Lun_out,*)'Problem in getting ',cdvar,' field 3D'
            kstatus = -99
         else
            if(idim(1).ne.nigauss*njgauss.or.idim(2).ne.G_nk) then
               kstatus = -99
            else
*     
*              Transfer real*8 to real
*              -----------------------
               if(.not.allocated(zbuff))allocate(zbuff(nigauss,njgauss,G_nk))
               zbuff(:,:,:) = 0.
*
               if(.not.G_lam) then
                 njx = njgauss
                 if(cdvar.eq.'VV'.and.V4dg_vstag_L) njx = njgauss -1
                 do jlev = 1, G_nk
                    icount = 0
                    do jlat = 1,njx
                       do jlon = 1,nigauss
                          icount = icount+1
                          zbuff(jlon,njx -jlat+1,jlev) = dlbuff_8(icount,jlev) 
                       end do
                    end do
                 end do
               else
                 njx = njgauss
                 if(cdvar.eq.'VV'.and.V4dg_vstag_L) njx = njgauss -1
                 do jlev = 1, G_nk
                    icount = 0
                    do jlat = 1,njx
                       do jlon = 1,nigauss
                          icount = icount+1
                          zbuff(jlon,jlat,jlev) = dlbuff_8(icount,jlev) 
                       end do
                    end do
                 end do
               endif
            end if
         end if
*
*     Change accuracy and reverse latitude if 2D field   
*     ------------------------------------------------
      case('PS')
         idim(1) = size ( dlbuff2d_8, 1 ) 
         write(Lun_out,*)' - Dimension of ',cdvar,' field 2D ',idim(1)
         if (istat.ne.0) then
            write(Lun_out,*)'Problem in getting ',cdvar,' field 2D'
            kstatus = -99
         else
            if(idim(1).ne.nigauss*njgauss) then
               kstatus = -99
            else
*
*              Transfer real*8 to real
*              -----------------------
               if(.not.allocated(zbuff2d)) allocate(zbuff2d(nigauss,njgauss))
               zbuff2d(:,:) = 0.
               icount = 0
               if(.not.G_lam) then
                  do jlat = 1,njgauss
                     do jlon = 1,nigauss
                        icount = icount+1
                        zbuff2d(jlon,njgauss -jlat+1) = dlbuff2d_8(icount) 
                     end do
                  end do
               else
                  do jlat = 1,njgauss
                     do jlon = 1,nigauss
                        icount = icount+1
                        zbuff2d(jlon,jlat) = dlbuff2d_8(icount) 
                     end do
                  end do
               endif
            end if
         end if
      end select
*
      end subroutine v4d_getfld
      end subroutine v4d_getdx