!-------------------------------------- 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 convert_avhrr(avhrr,sunzen) 1,4
! conversion des radiance IR en temperatures de brillance
! et des radiances visibles en "albedo"
use avhrr_var_mod
implicit none
type (avhrr_var) ,intent(inout) :: avhrr
real(8) ,intent(in) :: sunzen
integer :: ICL
REAL (8) :: tb(NIR),dtbsdrad(NIR)
REAL (8) :: FREQ(NIR),OFFSET(NIR),SLOPE(NIR)

DATA FREQ / 0.2687000000D+04 , 0.9272000000D+03 , 0.8377000000D+03 /
DATA OFFSET / 0.2066990000D+01 , 0.5512600000D+00 , 0.3471600000D+00 /
DATA SLOPE / 0.9965770000D+00  , 0.9985090000D+00  , 0.9989470000D+00 /

DO ICL=1,NCLASSAVHRR
   call calcbt(avhrr % radmoy(ICL,4:6), tb, dtbsdrad,freq,offset,slope)
   avhrr % tbmoy(ICL,4:6)=tb(1:3)
!   avhrr % tbmoy(ICL,5)=avhrr % tbmoy(ICL,5) + 0.3d0 ! application d'une correction de biais plate
!   avhrr % tbmoy(ICL,6)=avhrr % tbmoy(ICL,6) + 0.6d0 ! a la main (provisoire ?) 
!mise a jour des radiances pour rester coherent
!   avhrr %radmoy(ICL,5)=1000.d0*dplanck(freq(2), avhrr % tbmoy(ICL,5) ,offset(2),slope(2))
!   avhrr %radmoy(ICL,6)=1000.d0*dplanck(freq(3), avhrr % tbmoy(ICL,6) ,offset(3),slope(3))
   avhrr % tbstd(ICL,4:6)=avhrr % radstd(ICL,4:6) * dtbsdrad(1:3)

   call calcreflect(avhrr % radmoy(ICL,1:3) ,sunzen,avhrr % ALBEDMOY(ICL,1:3) )
   call calcreflect(avhrr % radstd(ICL,1:3) ,sunzen,avhrr % ALBEDSTD(ICL,1:3) )

ENDDO

contains

subroutine calcreflect(rad,sunzen,reflect) 2
implicit none

INTEGER ,parameter :: NVIS=3 ! Number of visible channels
REAL (8) , INTENT(IN) ,dimension(nvis) :: rad
REAL (8) , INTENT(IN) :: sunzen
REAL (8) , INTENT(out),dimension(nvis) :: reflect ! reflectivite en %
!************
REAL (8) ,DIMENSION(NVIS) :: SOLAR_FILTERED_IRRADIANCE
DATA SOLAR_FILTERED_IRRADIANCE /139.873215,232.919556,14.016470/
!# equivalent widths, integrated solar irradiance,  effective central wavelength
!0.084877,139.873215,0.632815
!0.229421,232.919556,0.841679
!0.056998,14.016470,1.606119
! pour la definition de l'albedo voir http://calval.cr.usgs.gov/PDF/Rao.CRN_IJRS.24.9.2003_Chander.pdf
REAL (8) :: RADB ! radiance en W/m2/str
!
REAL (8) :: PI
integer :: i
!**************************************************************

PI=ACOS(-1.0)

Do i = 1, nvis
   if (rad(i)>= 0.0 ) THEN
      radb=rad(i) / 1000.0
      reflect(i)=(PI*radb)/SOLAR_FILTERED_IRRADIANCE(I)
      IF (sunzen < 90.0 ) reflect(i)= reflect(i) / COS(sunzen*PI/180.0)
   else
      reflect(i)=-1
   end if
End Do

end subroutine calcreflect


subroutine calcbt(rad,tb,dtbsdrad,freq,offset,slope) 1
implicit none
INTEGER,parameter  :: nchan=3
Real(8) ,parameter :: c1= 1.19106590D-05   ! first planck constant
Real(8) ,parameter :: c2= 1.438833D0     ! second planck constant 
REAL (8) , INTENT(IN) ,dimension(nchan) :: rad,freq,offset,slope
REAL (8) , INTENT(out) ,dimension(nchan) :: tb,dtbsdrad
!************
integer :: i
REAL (8) ::  radtotal,tstore,planck1,planck2


Do i = 1, nchan
   if (rad(i)>1.d-20) THEN
      planck2= c2 * freq(I)
      planck1= c1 * ( freq(I) **3 ) 
      tstore = planck2 / Log( 1+planck1/rad(i) )
      tb(i) = ( tstore - offset(i) ) / slope(i)
     
      radtotal = rad(i)
   
      dtbsdrad(i) = planck1 * tstore**2 / ( planck2 * radtotal * ( radtotal + planck1 ) )
       
      dtbsdrad(i) = dtbsdrad(i) / slope(i)

      else
         tb(i) =0.d0
         dtbsdrad(i) = 0.d0
      end if

End Do

end subroutine calcbt


function dplanck(nu,t,offset,slope)
!    fonction de planck en double precision
!    nu en cm-1 t en Kelvin  planck en Watt / ( m2 strd cm-1 )
!    c en m/s  h en J.s  k en J/K
implicit none
real (8) :: nu,nu0,t,scale,offset,slope
real (8) :: dplanck,c,h,k,tt

c=299792458.D0
k=1.3806505D-23
h=6.62606876D-34
scale=100.d0
dplanck=-1.d0

if (t>0.d0) then
   nu0=nu*scale
   tt=t*slope+offset
   dplanck=scale*2.d0*h*c**2*nu0**3/(dexp(h*c/k*nu0/tt)-1.d0)
endif

return

end function dplanck


end subroutine convert_avhrr