!-------------------------------------- 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 getdx(CDTLMADJ) 22,9
*
*     Author: Simon
*
*     Purpose: Get Adjoint values of perturbations from adjoint
*     model
*     Revisions:
*     M. Tanguay *RPN/SMC Feb. 2002
*     - Add adjust for change of norm
*     - Fix to adjoint image winds->true winds conversion
*     M. Buehner *ARMA/MSC August 2002
*     - Added character argument to allow use in SV calculation
*     P. Gauthier *ARMA/MSC August 2003
*     .     Adjustment for staggered winds
*     S. Pellerin, ARMA, August 2008
*         .Added calls to 'tmg_*' subroutines
*     L. Fillion *ARMA/SMC 16 May 2005 - Allow Limited-Area Analysis Option.
*     L. Fillion *ARMA/SMC 14 Dec 2005 - For LAM: Transfer only non-extended fields;
*                                        i.e.:  (mni_in,mnj_in)
*     L. Fillion *ARMA/EC 6 Apr 2006 - Introduce INMI option: TL.
*     L. Fillion *ARMA/EC 10 Apr 2007 - Change /conphy by *conima below for 'F' and 'I' options.
*     L. Fillion *ARMA/EC 10 May 2007 - Introduce INMI option: ADJ.
*     L. Fillion *ARMA/EC 14 Aug 2007 - Update to v_10_0_3.
*     L. Fillion *ARMA/EC 8 Jan 2009 - Update to v_10_1_2.

*
*     ------------------
*     Arguments:  CDTLMADJ =
*     'J' read adjoint state vector for INMI
*     'A' read adjoint state vector for 4Dvar or SV
*     'X' read adjoint state for SV without norm adjustment (RWT)
*     'F' read state from forward integration for SV
*     'I' read state from forward integration for INMI
*     ------------------
      use modstag, only: r1qm2_s, lstagwinds
      implicit none
*
      include 'prof_f.h'
#include "comlun.cdk"
#include "comdim.cdk"
#include "comcva.cdk"
#include "comstate.cdk"
#include "comgd0.cdk"
#include "pardim.cdk"
#include "comgem.cdk"
#include "comvfiles.cdk"
#include "comleg.cdk"
#include "comcst.cdk"
#include "comgrd_param.cdk"
#include "comgemla.cdk"
*
      character CDTLMADJ
*
*     Local variables
*
      integer ini,inj
      integer :: prof_rdrec,ihdl,istat,ievent,prof_wrrec
      integer :: jlev,icount,jlat, jlon,idim2d, jgl, imax, jstag, jj
      integer, dimension(2) :: idim
      character (len=128) :: clprof
      real*8, pointer, dimension(:,:) :: zbuff
      real*8, pointer, dimension(:) :: zbuff2d
*
      real*8 :: dlconphy_vv(njbeg:njend)
*
      nullify(zbuff,zbuff2d)
*     Opening Prof_File
*
      if(CDTLMADJ.eq.'A'.or.CDTLMADJ.eq.'J'.or.CDTLMADJ.eq.'L'.or.CDTLMADJ.eq.'X') THEN
         clprof = trim(CEXC4DV) // '/dwga.prof'
         write(nulout,*
     &        )'GETDX: Reading adjoint data from simulation no. ',nsim3d
      elseif(CDTLMADJ.eq.'F'.or.CDTLMADJ.eq.'I'.or.CDTLMADJ.eq.'K') THEN
         clprof = trim(CEXC4DV) // '/dwgf1.prof'
         write(nulout,*
     &        )'GETDX: Reading TLM data from simulation no. ',nsim3d
      endif
!
      if(grd_typ.eq.'LU') then
        ini = mni_in
        inj = mnj_in
      else
        ini = ni
        inj = nj
      endif
!
      write(nulout,*) 'BEFORE PROF_OPEN'
      call vflush(nulout)
      call tmg_start(82,'I/O_WAIT')
      ihdl = prof_open(clprof,'READ','FILE')
      call tmg_stop(82)
      write(nulout,*) 'AFTER PROF_OPEN'
      call vflush(nulout)
      call tmg_start(81,'PROF_R+W')
      istat = prof_rdrec(ihdl)
      call tmg_stop(81)
      istat = prof_gvar(ihdl,ievent,prm_evnt)
      if(ievent.eq.nsim3d) then
         write(nulout,*) 'GETDX: Data received for simulation no. '
     &        ,ievent

         if(ngexist(nguu).eq.1) then
            istat = prof_gvar(ihdl,zbuff,v3d_utru)
            idim  = ubound(zbuff)
            if (ini*inj.ne.idim(1)) then
               write(nulout,*) 'IN GETDX, BUFFER has ini*inj = '
     S              ,idim(1),' NFLEV = ',idim(2)
               write(nulout,*) 'Internal dimensions are ini = ', ini
     S              ,' inj = ',inj,' NFLEV = ',nflev
               call abort3d(nulout
     &              ,'GETDX- Inconsistent model states FOR UU')
            end if
            if(CDTLMADJ.eq.'A'.or.CDTLMADJ.eq.'J'.or.CDTLMADJ.eq.'L'.or.CDTLMADJ.eq.'X') THEN
               do jlev = 1, nflev
                  icount = 0
                  do jlat = 1, inj
                     do jlon = 1, ini
                        icount = icount + 1
                        ut0(jlon,jlev,jlat) = zbuff(icount,jlev)*conphy(jlat)
                     end do
                  end do
               end do
            elseif(CDTLMADJ.eq.'F'.or.CDTLMADJ.eq.'I'.or.CDTLMADJ.eq.'K') THEN
               do jlev = 1, nflev
                  icount = 0
                  do jlat = 1, inj
                     do jlon = 1, ini
                        icount = icount + 1
                        ut0(jlon,jlev,jlat) = zbuff(icount,jlev)*conima(jlat)
                     end do
                  end do
               end do
            endif
         end if
*
         if(ngexist(ngvv).eq.1) then
!
           if(lstagwinds) then
             if(grd_typ.ne.'LU') then
                 dlconphy_vv(:) = ra*r1qm2_s(:)
             else
               do jj=1,inj
                 dlconphy_vv(jj) = rrcos_an(1,jj) ! 1.0/cosh(theta)
               enddo
             endif
            else
               dlconphy_vv(:) = conphy(:)
            end if
!
            istat = prof_gvar(ihdl,zbuff,v3d_vtru)
            idim  = ubound(ZBUFF)
            if (ini*inj.ne.idim(1)) then
               write(nulout,*) 'IN GETDX, BUFFER has ini*inj = '
     S              ,idim(1),' NFLEV = ',idim(2)
               write(nulout,*) 'Internal dimensions are ini = ', ini
     S              ,' inj = ',inj,' NFLEV = ',nflev
               call abort3d(nulout
     &              ,'GETDX- Inconsistent model states FOR VV')
            end if
            if(CDTLMADJ.eq.'A'.or.CDTLMADJ.eq.'J'.or.CDTLMADJ.eq.'L'.or.CDTLMADJ.eq.'X') THEN
               if (lstagwinds) then
                 jstag = inj -1
               else
                 jstag = inj
               endif
               do jlev = 1, nflev
                  icount = 0
                  do jlat = 1, jstag
                     do jlon = 1, ini
                        icount = icount + 1
                        vt0(jlon,jlev,jlat) = zbuff(icount,jlev)*dlconphy_vv(jlat)
                     end do
                  end do
               end do
            elseif(CDTLMADJ.eq.'F'.or.CDTLMADJ.eq.'I'.or.CDTLMADJ.eq.'K') THEN
               if (lstagwinds) then
                 jstag = inj -1
               else
                 jstag = inj
               endif
               do jlev = 1, nflev
                  icount = 0
                  do jlat = 1, jstag
                     do jlon = 1, ini
                        icount = icount + 1
                        vt0(jlon,jlev,jlat) = zbuff(icount,jlev)/dlconphy_vv(jlat)
                     end do
                  end do
               end do
            endif
         end if
*
         if(ngexist(ngtt).eq.1) then
            ISTAT = PROF_GVAR(IHDL,ZBUFF,V3D_TEMP)
            IDIM  = UBOUND(ZBUFF)
            IF (ini*inj.NE.IDIM(1)) THEN
               WRITE(NULOUT,*) 'IN GETDX, BUFFER has ini*inj = '
     S              ,idim(1),' NFLEV = ',idim(2)
               write(nulout,*) 'Internal dimensions are ini = ', ini
     S              ,' inj = ',inj,' NFLEV = ',nflev
               call abort3d(nulout
     &              ,'GETDX- Inconsistent model states FOR TT')
            end if
            DO JLEV = 1, NFLEV
               ICOUNT = 0
               DO JLAT = 1, inj
                  DO JLON = 1, ini
                     ICOUNT = ICOUNT + 1
                     TT0(JLON,JLEV,JLAT) = ZBUFF(ICOUNT,JLEV)
                  END DO
               END DO
            END DO
         END IF
*
         IF(NGEXIST(NGQ).EQ.1) THEN
            ISTAT = PROF_GVAR(IHDL,ZBUFF,V3D_SPHU)
            IDIM  = UBOUND(ZBUFF)
            IF (ini*inj.NE.IDIM(1)) THEN
               WRITE(NULOUT,*) 'IN GETDX, BUFFER has ini*inj = '
     S              ,idim(1),' NFLEV = ',idim(2)
               write(nulout,*) 'Internal dimensions are ini = ', ini
     S              ,' inj = ',inj,' NFLEV = ',nflev
               call abort3d(nulout
     &              ,'GETDX- Inconsistent model states FOR Q')
            end if
            DO JLEV = 1, NFLEV
               ICOUNT = 0
               DO JLAT = 1, inj
                  DO JLON = 1, ini
                     ICOUNT = ICOUNT + 1
                     Q0(JLON,JLEV,JLAT) = ZBUFF(ICOUNT,JLEV)
                  END DO
               END DO
            END DO
         END IF
*
         deallocate(ZBUFF)

         call tmg_start(81,'PROF_R+W')
         istat = prof_rdrec(ihdl)
         call tmg_stop(81)
         istat = prof_gvar(ihdl,ievent,prm_evnt) + istat

         IF(NGEXIST(NGPS).EQ.1) THEN
            ISTAT = PROF_GVAR(IHDL,ZBUFF2D,V2D_PSUR)
            IDIM2D  = UBOUND(ZBUFF2D,1)
            IF (ini*inj.NE.IDIM2D) THEN
               WRITE(NULOUT,*) 'IN GETDX, BUFFER2D has ini*inj = '
     S              ,idim2d
               write(nulout,*) 'Internal dimensions should be ini = ', ini
     S              ,' inj = ',inj,' NFLEV = ',nflev
               call abort3d(nulout
     &              ,'GETDX- Inconsistent model states FOR PS')
            end if
*
            ICOUNT = 0
            DO JLAT = 1, inj
               DO JLON = 1, ini
                  ICOUNT = ICOUNT + 1
                  GPS0(JLON,1,JLAT)= ZBUFF2D(ICOUNT)
               END DO
            END DO
         END IF
*
         IF(NGEXIST(NGTG).EQ.1) THEN
c     ISTAT = PROF_GVAR(IHDL,ZBUFF2D,V2D_TGRN)
c     IDIM2D  = UBOUND(ZBUFF2D,1)
c     IF (ini*inj.NE.IDIM2D) THEN
c     WRITE(NULOUT,*) 'IN GETDX, BUFFER2D has ini*inj = '
c     S           ,idim2d
c     write(nulout,*) 'Internal dimensions should be ini = ', ini
c     S           ,' inj = ',inj,' NFLEV = ',nflev
c     call abort3d(nulout
c     &           ,'GETDX- Inconsistent model states for PS')
c     end if
c     *
c     ICOUNT = 0
c     DO JLAT = 1, inj
c     DO JLON = 1, ini
c     ICOUNT = ICOUNT + 1
c     GTG0(JLON,1,JLAT)= ZBUFF2D(ICOUNT)
c     END DO
c     END DO

            ICOUNT = 0
            DO JLAT = 1, inj
               DO JLON = 1, ini
                  ICOUNT = ICOUNT + 1
                  GTG0(JLON,1,JLAT)= 0.0
               END DO
            END DO
         END IF
*
         deallocate (ZBUFF2D)

         istat = prof_close(ihdl,.true.)
*
*     Add adjustment for change of norm as in BILINAD
*     -----------------------------------------------
         if(CDTLMADJ.eq.'A'.and.grd_typ.ne.'LU') THEN
            do jgl = 1, inj
               imax = nilon(jgl)
               do jlev = 1, nkgdim
                  do jlon = 1, imax
                     gd(jlon,jlev,jgl) = gd(jlon,jlev,jgl) *
     +                    nilon(jgl) / rwt(jgl)
                  enddo
               enddo
            enddo
         endif

      else
        write(nulout,*)
     &       'GETDX - ERROR - ABNORMAL TERMINATION : ievent = ',ievent
         clprof = 'dwgf.prof'
         IHDL = PROF_OPEN(clprof,'WRITE','FILE')
         ISTAT = PROF_PVAR(IHDL,evn_ferr,PRM_EVNT)
         call tmg_start(81,'PROF_R+W')
         istat = prof_wrrec(ihdl)
         istat = prof_close(ihdl)
         call tmg_stop(81)
         call abort3d(nulout,'GETDX')
      endif
*
      end subroutine getdx