!-------------------------------------- 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(pgomu,pgomv,pgomt,pgomgz,pgomq,pgomps,pgomtgr, & 4,13
 &                       prppobs,knkgdimo,knlev)
#if defined (doc)
!
!**s/r tovs_fill_profiles  - interpolation of obs. space to 43-level rttov space
!                      (adapted from part of code of dobstovs)
!
!
!author        : j. halle *cmda/aes  december 13, 2004
!
!revision 001  : a. beaulne *cmda/smc  june 2006
!                    -add ozone from climatology to all sensors
!                    -modifications for AIRS :
!                       + addition of geopotential field to subr argument
!                       + addition of latitude, longitude, height field,
!                         sun zenith angle and cloud fraction 
!                         to personnalized profile structure
!
!revision 002  : j. halle  *cmda/smc  march 2007
!                    -fix zvlev for hybrid coordinate
!
!    -------------------
!     purpose: fill profiles structure with info from obs. space
!
!arguments
!
!
#endif
  use mod_tovs
  use ozoneclim
 
  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 "partov.cdk"
#include "comtov.cdk"
#include "comtovst.cdk"

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

  integer, allocatable :: ksurf     (:) 
  integer, allocatable :: iptobs    (:) 
  integer, allocatable :: iptobscma (:) 
  integer, allocatable :: isatzen   (:)
  integer, allocatable :: isatazim  (:)
  integer, allocatable :: isunza    (:)
  integer, allocatable :: iclfr     (:)
  real*8,  allocatable :: zlat      (:)
  real*8,  allocatable :: zlon      (:)
  real*8,  allocatable :: ozo     (:,:)

  integer :: alloc_status(40)

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

  real*8 :: pgomu  (knkgdimo,*)
  real*8 :: pgomv  (knkgdimo,*)
  real*8 :: pgomt  (knkgdimo,*)
  real*8 :: pgomgz (knkgdimo,*)
  real*8 :: pgomq  (knkgdimo,*)
  real*8 :: pgomps (knkgdimo,*)
  real*8 :: pgomtgr(knkgdimo,*)
  real*8 :: prppobs(knlev   ,*)

  real*8, allocatable :: to    (:,:)
  real*8, allocatable :: lqo   (:,:)
  real*8, allocatable :: zho   (:,:)
  real*8, allocatable :: toext (:,:)
  real*8, allocatable :: qoext (:,:)
  real*8, allocatable :: zhoext(:,:)
  real*8, allocatable :: zvlev (:,:)
  real*8, allocatable :: zt    (:,:)
  real*8, allocatable :: zlq   (:,:)
  real*8, allocatable :: zht   (:,:)

  real*8 :: zptop, zptopmbs, ztodeg

  real*8, allocatable :: toto3obs(:)

  integer  isrcheq
  external isrcheq
!  external lintv2, extrap
  external intavg, extrap
  external lintv2
  external exthum4

  if(nobtov.eq.0) return    ! exit if there are no tovs data

  ztodeg = 180.0 / rpi

!     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 (ksurf     (jppf*nthreads)        ,stat= alloc_status(1) )
  allocate (iptobs    (jppf*nthreads)        ,stat= alloc_status(2) )
  allocate (iptobscma (jppf*nthreads)        ,stat= alloc_status(3) )
  allocate (isatzen   (jppf*nthreads)        ,stat= alloc_status(4) )
  allocate (isatazim  (jppf*nthreads)        ,stat= alloc_status(5) )
  allocate (isunza    (jppf*nthreads)        ,stat= alloc_status(6) )
  allocate (iclfr     (jppf*nthreads)        ,stat= alloc_status(7) )
  allocate (zlat      (jppf*nthreads)        ,stat= alloc_status(8) )
  allocate (zlon      (jppf*nthreads)        ,stat= alloc_status(9) )
  allocate (ozo       (jplev  ,jppf*nthreads),stat= alloc_status(10))
  allocate (to        (jpmolev,jppf*nthreads),stat= alloc_status(11))
  allocate (lqo       (jpmolev,jppf*nthreads),stat= alloc_status(12))
  allocate (zho       (jpmolev,jppf*nthreads),stat= alloc_status(13))
  allocate (toext     (jplev  ,jppf*nthreads),stat= alloc_status(14))
  allocate (qoext     (jplev  ,jppf*nthreads),stat= alloc_status(15))
  allocate (zhoext    (jplev  ,jppf*nthreads),stat= alloc_status(16))
  allocate (zvlev     (jpnflev,jppf*nthreads),stat= alloc_status(17))
  allocate (zt        (jpnflev,jppf*nthreads),stat= alloc_status(18))
  allocate (zlq       (jpnflev,jppf*nthreads),stat= alloc_status(19))
  allocate (zht       (jpnflev,jppf*nthreads),stat= alloc_status(20))

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

!     1.2   Read ozone climatology
!     .     ----------------------

  call readozoneclim 

!
!     2.  Fill profiles structure
!     .   -----------------------

  200  continue

! loop over all obs.

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

     knpf = knpf + 1
!    extract land/sea/sea-ice flag (0=land, 1=sea, 2=sea-ice)

     ilansea = mobhdr(ncmofl,jo)
     indxm = ilansea
     if (indxm .eq. 2 ) indxm = 0
     indxreg = isrcheq ( mlisreg, nregst, indxm )
     ksurf(knpf) = ilansea

!    extract satellite zenith and azimuth angle, 
!    sun zenith angle, cloud fraction, latitude and longitude

     isatzen(knpf) = mobhdr(ncmbox,jo)/10000
     isatazim(knpf) = mobhdr(ncmaza,jo)
     isunza(knpf) = mobhdr(ncmsun,jo)
     iclfr(knpf) = mobhdr(ncmclf,jo)
     zlat(knpf) = robhdr(ncmlat,jo) * ztodeg
     zlon(knpf) = robhdr(ncmlon,jo) * ztodeg

     do jl = 1, knlev
        zt   (jl,knpf) = pgomt(jl,jo)
        zlq  (jl,knpf) = pgomq(jl,jo)
        zvlev(jl,knpf) = prppobs(jl,jo) * rpatmb
        zht  (jl,knpf) = pgomgz(jl,jo) / rg
     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.1  Vertical interpolation of model temperature, logarithm of
!             specific humidity and height levels to pressure levels
!             required by tovs rt model
!     .       --------------------------------------------------------------

!$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 lintv2 (zvlev(1,jn),zt(1,jn),jpnflev,knlev,ilen, &
!                        jpmolev,xpres(jpmotop),to(1,jn))
!           call lintv2 (zvlev(1,jn),zlq(1,jn),jpnflev,knlev,ilen, &
!                        jpmolev,xpres(jpmotop),lqo(1,jn))
           call intavg (zvlev(1,jn),zt(1,jn),jpnflev,knlev,ilen, &
                        jpmolev,xpres(jpmotop),to(1,jn))
           call intavg (zvlev(1,jn),zlq(1,jn),jpnflev,knlev,ilen, &
                        jpmolev,xpres(jpmotop),lqo(1,jn))
           call lintv2 (zvlev(1,jn),zht(1,jn),jpnflev,knlev,ilen, &
                        jpmolev,xpres(jpmotop),zho(1,jn))
        endif
     enddo
!$omp end do
!$omp end parallel

!     .  2.2  Extrapolation of temperature profile above model top
!     .       ----------------------------------------------------

     toext(:,:) = 0.0
     call extrap (to,toext,jpmolev,jplev,knpf)

!     .  2.3  Extrapolation of height profile above model top
!     .       -----------------------------------------------

     zhoext(:,:) = 0.0
     call htextrap (zhoext,zho,xpres,jplev,jpmolev,jpmotop,knpf)

!     .  2.4  Extrapolation of humidity profile (kg/kg)
!             above rlimlvhu (normally 300mbs or 70mbs)
!     .       -----------------------------------------

     qoext(:,:) = 0.0
     do jk = 1, jpmolev
        do jn = 1, knpf
           qoext(jplev-jpmolev+jk,jn) = exp(lqo(jk,jn)) 
        enddo
     enddo
     if ( ldbgtov ) then
        do jn = 1, knpf
           write(nulout,*)'qoext*1000 avant exthum4    = '
           write(nulout,9263)(qoext(i,jn)*1000.,i=1,jplev)
           write(nulout,*)' '
        enddo
     endif
     call exthum4 (knpf,jplev,xpres,qoext)
     if ( ldbgtov ) then
        do jn = 1, knpf
           write(nulout,*)'qoext*1000 apres exthum4    = '
           write(nulout,9263)(qoext(i,jn)*1000.,i=1,jplev)
           write(nulout,*)' '
        enddo
     endif

!     .  2.5  Get ozone profiles (ppmv)
!     .       -------------------------

     allocate ( toto3obs(knpf) )     

     toto3obs(:) = 0.

     call new_ozone (ozo, fozo,totozo,toto3obs,zlat,  &
  &                       xpres,jplev,po3,nlato3,nlevo3,knpf)

     deallocate ( toto3obs )

!     .  2.6  Fill profiles structure
!     .       -----------------------

     no_id = 1
     do  j = 1 , knpf 
        profiles(iptobs(j)) % nlevels         = coef(no_id) % nlevels
        profiles(iptobs(j)) % ozone_data      = .true.
        profiles(iptobs(j)) % co2_data        = .false.
        profiles(iptobs(j)) % clw_data        = .false.
        profiles(iptobs(j)) % p(:)            = coef(no_id) % ref_prfl_p(:)
        profiles(iptobs(j)) % t(:)            = toext(:,j)
        profiles(iptobs(j)) % q(:)            = qoext(1:jplev,j) * q_mixratio_to_ppmv
        profiles(iptobs(j)) % o3(:)           = ozo(:,j)
        profiles(iptobs(j)) % ctp             = 500.0
        profiles(iptobs(j)) % cfraction       = 0.0
        profiles(iptobs(j)) % zenangle        = (isatzen(j)-9000)/100.0
        profiles(iptobs(j)) % azangle         = isatazim(j)/100.0 
        profiles(iptobs(j)) % skin % surftype = ksurf(j)
        profiles(iptobs(j)) % skin % t        = pgomtgr(1,iptobscma(j))
        profiles(iptobs(j)) % skin % fastem(:)= 0.0
        profiles(iptobs(j)) % s2m % t         = pgomt(ilowlvl,iptobscma(j))
        profiles(iptobs(j)) % s2m % q         = exp(pgomq(ilowlvl,iptobscma(j))) * q_mixratio_to_ppmv
        profiles(iptobs(j)) % s2m % p         = pgomps(1,iptobscma(j))*rpatmb
        profiles(iptobs(j)) % s2m % u         = pgomu(ilowlvl,iptobscma(j))
        profiles(iptobs(j)) % s2m % v         = pgomv(ilowlvl,iptobscma(j))
        profiles_qc(iptobs(j)) % lat          = zlat(j)
        profiles_qc(iptobs(j)) % lon          = zlon(j)
        profiles_qc(iptobs(j)) % z(:)         = zhoext(:,j)
        profiles_qc(iptobs(j)) % sunza        = (isunza(j)-9000)/100.0
        profiles_qc(iptobs(j)) % clfr         = iclfr(j)
     end do
 
!    next bunch !

     knpf = 0

  290 continue

  alloc_status(:) = 0
  deallocate (ksurf     ,stat= alloc_status(1) )
  deallocate (iptobs    ,stat= alloc_status(2) )
  deallocate (iptobscma ,stat= alloc_status(3) )
  deallocate (isatzen   ,stat= alloc_status(4) )
  deallocate (isatazim  ,stat= alloc_status(5) )
  deallocate (isunza    ,stat= alloc_status(6) )
  deallocate (iclfr     ,stat= alloc_status(7) )
  deallocate (zlat      ,stat= alloc_status(8) )
  deallocate (zlon      ,stat= alloc_status(9) )
  deallocate (ozo       ,stat= alloc_status(10))
  deallocate (to        ,stat= alloc_status(11))
  deallocate (lqo       ,stat= alloc_status(12))
  deallocate (zho       ,stat= alloc_status(13))
  deallocate (toext     ,stat= alloc_status(14))
  deallocate (qoext     ,stat= alloc_status(15))
  deallocate (zhoext    ,stat= alloc_status(16))
  deallocate (zvlev     ,stat= alloc_status(17))
  deallocate (zt        ,stat= alloc_status(18))
  deallocate (zlq       ,stat= alloc_status(19))
  deallocate (zht       ,stat= alloc_status(20))

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

  9263 format(1x,10f8.4)


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

  return

end subroutine tovs_fill_profiles