!-------------------------------------- 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 tovs_fill_profiles_ad(pgomu_ad,pgomv_ad,pgomt_ad,pgomq_ad,pgomps_ad,pgomtgr_ad, & 1,7
 &                          prppobs,knkgdimo,knlev)
#if defined (doc)
!
!**s/r tovs_fill_profiles_ad  - adjoint of interpolation of obs.space to 43-level rttov space
!                          (adapted from part of code of avtov)
!
!
!author        : j. halle *cmda/aes  april 21, 2005
!
!revision 001  : j. halle  *cmda/smc  march 2007
!                    -fix zvlev for hybrid coordinate
!
!    -------------------
!     purpose: fill tangent linear profiles structure with info from obs. space
!
!arguments
!
!
#endif
  use mod_tovs

  implicit none
!implicits
#include "pardim.cdk"
#include "comdim.cdk"
#include "comgem.cdk"
#include "comdimo.cdk"
#include "comlun.cdk"
#include "comcst.cdk"
#include "comoabdy.cdk"
#include "comoahdr.cdk"
#include "comoba.cdk"
#include "commvog.cdk"
#include "commvo1.cdk"
#include "partov.cdk"
#include "comtov.cdk"
#include "comtovst.cdk"

  real*8, parameter :: q_mixratio_to_ppmv  = 1.60771704e+6

  integer, allocatable :: iptobs    (:) 
  integer, allocatable :: iptobscma (:) 

  integer :: alloc_status(40)

  integer :: omp_get_num_threads, nthreads
  integer :: istride, ilen, imodulo
  integer :: nlevels
  integer :: j, i, krtid, jf, iobs, jj
  integer :: ilowlvl, knpf, jo
  integer :: jk, jn, jl
  integer :: knkgdimo, knlev

  real*8 :: pgomu_ad  (knkgdimo,*)
  real*8 :: pgomv_ad  (knkgdimo,*)
  real*8 :: pgomt_ad  (knkgdimo,*)
  real*8 :: pgomq_ad  (knkgdimo,*)
  real*8 :: pgomps_ad (knkgdimo,*)
  real*8 :: pgomtgr_ad(knkgdimo,*)
  real*8 :: prppobs   (knlev   ,*)

  real*8, allocatable :: to_ad    (:,:)
  real*8, allocatable :: lqo_ad   (:,:)
  real*8, allocatable :: toext_ad (:,:)
  real*8, allocatable :: qoext_ad (:,:)
  real*8, allocatable :: zvlev    (:,:)
  real*8, allocatable :: zt_ad    (:,:)
  real*8, allocatable :: zlq_ad   (:,:)
  real*8, allocatable :: zt       (:,:)
  real*8, allocatable :: zlq      (:,:)
  real*8, allocatable :: qoext    (:,:)
  real*8, allocatable :: zps_ad   (:)

  real*8 :: zptop, zptopmbs

!  external alintv2, aextrap
  external intavgad, aextrap
  external aexthum4
         
  if(nobtov.eq.0) return    ! exit if there are not tovs data


!     1.    Set index for model's lowest level and model top
!     .     ------------------------------------------------
                 
  100   continue

  if (  prppobs(1,1) .lt. prppobs(knlev,1) ) then
     ilowlvl = knlev
  else
     ilowlvl = 1
  endif

! find model level top, within 0.000001 mbs.

  zptop    = prppobs(1,1)
  zptopmbs = zptop/100.
  zptopmbs = zptopmbs - 0.000001

  jpmotop = 1
  do jl = 2, jplev
     if ( zptopmbs .ge. xpres(jl-1) .and.       &
          zptopmbs .lt. xpres(jl)        ) then
        jpmotop = jl
        go to 110
     endif
  enddo

  110 continue
  jpmolev = (jplev-jpmotop+1)


!     1.1   Number of threads and memory allocation
!     .     ---------------------------------------

!$omp parallel 
  nthreads = omp_get_num_threads()
!$omp end parallel

  alloc_status(:) = 0
  allocate (iptobs   (jppf*nthreads)        ,stat= alloc_status(1) )
  allocate (iptobscma(jppf*nthreads)        ,stat= alloc_status(2) )
  allocate (to_ad    (jpmolev,jppf*nthreads),stat= alloc_status(3) )
  allocate (lqo_ad   (jpmolev,jppf*nthreads),stat= alloc_status(4) )
  allocate (toext_ad (jplev  ,jppf*nthreads),stat= alloc_status(5) )
  allocate (qoext_ad (jplev  ,jppf*nthreads),stat= alloc_status(6) )
  allocate (zvlev    (jpnflev,jppf*nthreads),stat= alloc_status(7) )
  allocate (zt_ad    (jpnflev,jppf*nthreads),stat= alloc_status(8) )
  allocate (zlq_ad   (jpnflev,jppf*nthreads),stat= alloc_status(9) )
  allocate (zt       (jpnflev,jppf*nthreads),stat= alloc_status(10))
  allocate (zlq      (jpnflev,jppf*nthreads),stat= alloc_status(11))
  allocate (qoext    (jplev  ,jppf*nthreads),stat= alloc_status(12))
  allocate (zps_ad   (jppf*nthreads)        ,stat= alloc_status(13))

  if( any(alloc_status /= 0) ) then
     write(nulout,*) ' tovs_fill_profiles_ad : memory allocation error'
     call abort3d(nulout,'tovs_fill_profiles_ad     ')
  end if
!
!     2.  Adjoint of fill profiles structure
!     .   ----------------------------------

  200  continue

! loop over all obs.

  knpf = 0
  do 290 iobs = 1, nobtov
     jo = lobsno(iobs)

     knpf = knpf + 1

     do jl = 1, knlev
        zt   (jl,knpf) = gomt1(jl,jo)
        zlq  (jl,knpf) = gomq1(jl,jo)
        zvlev(jl,knpf) = prppobs(jl,jo) * rpatmb
     enddo

     iptobs   (knpf) = iobs
     iptobscma(knpf) = jo
     if ( knpf .le. 0                                    ) go to 290
     if ( knpf .ne. jppf*nthreads .and. iobs .ne. nobtov ) go to 290

!     .  2.4  Adjoint of filling profiles_ad structure
!     .       ----------------------------------------

     do j = 1, knpf
        qoext(:,j) =  profiles(iptobs(j)) % q(:) / q_mixratio_to_ppmv
     enddo

     do  j = 1 , knpf 
        toext_ad(:,j)                  = profiles_ad(iptobs(j)) % t(:)
        qoext_ad(:,j)                  = profiles_ad(iptobs(j)) % q(:) * q_mixratio_to_ppmv
        pgomps_ad(1,iptobscma(j))      = profiles_ad(iptobs(j)) % s2m % p * rpatmb
        pgomtgr_ad(1,iptobscma(j))     = profiles_ad(iptobs(j)) % skin % t
        pgomt_ad(ilowlvl,iptobscma(j)) = profiles_ad(iptobs(j)) % s2m % t
!!!        pgomq_ad(ilowlvl,iptobscma(j)) = qoext(ilowlvl,j) * profiles_ad(iptobs(j)) % s2m % q * q_mixratio_to_ppmv
!!!        pgomq_ad(ilowlvl,iptobscma(j)) = qoext(ilowlvl,j) * profiles_ad(iptobs(j)) % s2m % q 
        pgomq_ad(ilowlvl,iptobscma(j)) = 0. 
        pgomu_ad(ilowlvl,iptobscma(j)) = profiles_ad(iptobs(j)) % s2m % u
        pgomv_ad(ilowlvl,iptobscma(j)) = profiles_ad(iptobs(j)) % s2m % v
     end do

!     .  2.3  Adjoint of extrapolation of humidity profile (kg/kg)
!             above rlimlvhu (normally 300mbs or 70mbs)
!     .  

     if ( ldbgtov ) then
        do jn = 1, knpf
           write(nulout,*)'qoext_ad*1000 avant aexthum4    = '
           write(nulout,9263)(qoext_ad(i,jn)*1000.,i=1,jplev)
           write(nulout,*)' '
        enddo
     endif
     call aexthum4 (knpf,jplev,xpres,qoext_ad,qoext)
     if ( ldbgtov ) then
        do jn = 1, knpf
           write(nulout,*)'qoext_ad*1000 apres aexthum4    = '
           write(nulout,9263)(qoext_ad(i,jn)*1000.,i=1,jplev)
           write(nulout,*)' '
        enddo
     endif

!    adjoint of conversion lnq --> q

     lqo_ad(:,:) = 0.0
     do jk = 1, jpmolev
        do jn = 1, knpf
           lqo_ad(jk,jn) = qoext_ad(jplev-jpmolev+jk,jn) * qoext(jplev-jpmolev+jk,jn)
        enddo
     enddo

!     .  2.2  Adjoint of extrapolation of temperature profile above 10mb
!     .       ----------------------------------------------------------

     to_ad(:,:) = 0.0
     call aextrap (to_ad,toext_ad,jpmolev,jplev,knpf)
 
!     .  2.1  Adjoint of vertical interpolation of model temperature and logarithm of
!             specific humidity to pressure levels required by tovs rt model
!     .       -----------------------------------------------------------------------

     zt_ad (:,:) = 0.0
     zlq_ad(:,:) = 0.0
     zps_ad(:)   = 0.0

!$omp parallel private(istride)
     imodulo = mod(knpf,nthreads)
     if ( imodulo .eq. 0 ) then
        istride = max(1,(knpf/nthreads))
     else
        istride = max(1,(knpf/nthreads)+1)
     endif
!$omp do private(jn,ilen)
     do jn=1,knpf,istride
        ilen = min (knpf-jn+1,istride)
        if ( ilen .gt. 0) then

!            call alintv2 (zvlev(1,jn),zt_ad(1,jn),zt(1,jn), &
!     &                    zps_ad(jn),jpnflev,knlev,ilen, &
!     &                    jpmolev,xpres(jpmotop),to_ad(1,jn))
!            call alintv2 (zvlev(1,jn),zlq_ad(1,jn),zlq(1,jn), &
!     &                    zps_ad(jn),jpnflev,knlev,ilen, &
!     &                    jpmolev,xpres(jpmotop),lqo_ad(1,jn))
            call intavgad (zvlev(1,jn),zt_ad(1,jn),zt(1,jn), &
     &                    zps_ad(jn),jpnflev,knlev,ilen, &
     &                    jpmolev,xpres(jpmotop),to_ad(1,jn))
            call intavgad (zvlev(1,jn),zlq_ad(1,jn),zlq(1,jn), &
     &                    zps_ad(jn),jpnflev,knlev,ilen, &
     &                    jpmolev,xpres(jpmotop),lqo_ad(1,jn))

        endif
     enddo
!$omp end do
!$omp end parallel
 
!     .  2.1  Store adjoints in commvo
!     .       ------------------------

     do  jn = 1 , knpf 
        pgomps_ad(1,iptobscma(jn)) = pgomps_ad(1,iptobscma(jn)) + zps_ad  (jn) * rpatmb
        do jl = 1, knlev
           pgomt_ad(jl,iptobscma(jn)) = pgomt_ad(jl,iptobscma(jn)) + zt_ad  (jl,jn)
           pgomq_ad(jl,iptobscma(jn)) = pgomq_ad(jl,iptobscma(jn)) + zlq_ad (jl,jn)
        enddo
     enddo
 
!    next bunch !

     knpf = 0

  290 continue

  alloc_status(:) = 0
  deallocate (iptobs   ,stat= alloc_status(1) )
  deallocate (iptobscma,stat= alloc_status(2) )
  deallocate (to_ad    ,stat= alloc_status(3) )
  deallocate (lqo_ad   ,stat= alloc_status(4) )
  deallocate (toext_ad ,stat= alloc_status(5) )
  deallocate (qoext_ad ,stat= alloc_status(6) )
  deallocate (zvlev    ,stat= alloc_status(7) )
  deallocate (zt_ad    ,stat= alloc_status(8) )
  deallocate (zlq_ad   ,stat= alloc_status(9) )
  deallocate (zt       ,stat= alloc_status(10))
  deallocate (zlq      ,stat= alloc_status(11))
  deallocate (qoext    ,stat= alloc_status(12))
  deallocate (zps_ad   ,stat= alloc_status(13))

  if( any(alloc_status /= 0) ) then
     write(nulout,*) ' tovs_fill_profiles_ad : memory deallocation error'
     call abort3d(nulout,'tovs_fill_profiles_ad        ')
  end if

  9263 format(1x,10f8.4)


!     3.  Close up
!     .   --------

  return

end subroutine tovs_fill_profiles_ad