Mrweight.F90 6.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225
  1. ! First include the set of model-wide compiler flags
  2. #include "tm5.inc"
  3. module Mrweight
  4. !---------------------------------------------------------------------------
  5. ! Contains routine to compute the intensity at the top of the atmosphere
  6. ! for a Rayleigh scattering atmosphere with a Lambertian surface.
  7. !
  8. ! call rweight(wavelx,p_sfc,a_sfc,th0,th,dphix,refl)
  9. !
  10. ! Pepijn Veefkind, Henk Eskes, KNMI, 1999
  11. !---------------------------------------------------------------------------
  12. implicit none
  13. private
  14. public :: rweight
  15. contains
  16. subroutine csalbr(xtau,xalb)
  17. implicit none
  18. real,intent(in) :: xtau
  19. real,intent(out) :: xalb
  20. xalb=0
  21. xalb=(3*xtau-fintexp3(xtau)*(4+2*xtau)+2*exp(-xtau))
  22. xalb=xalb/(4.+3*xtau)
  23. end subroutine csalbr
  24. function fintexp3(xtau)
  25. implicit none
  26. real,intent(in) :: xtau
  27. real :: fintexp3
  28. real :: xx
  29. xx=(exp(-xtau)*(1.-xtau)+xtau*xtau*fint1exp(xtau))/2.
  30. fintexp3=xx
  31. end function fintexp3
  32. function fint1exp(xtau)
  33. ! accuracy 2e-07... for 0<xtau<1
  34. implicit none
  35. real,intent(in) :: xtau
  36. real :: fint1exp
  37. real :: xx,xftau
  38. real,dimension(0:5),parameter :: &
  39. a=(/-.57721566, 0.99999193,-0.24991055, &
  40. 0.05519968,-0.00976004, 0.00107857 /)
  41. integer :: i
  42. xx=a(0)
  43. xftau=1.
  44. do i=1,5
  45. xftau=xftau*xtau
  46. xx=xx+a(i)*xftau
  47. enddo
  48. fint1exp=xx-log(xtau)
  49. end function fint1exp
  50. subroutine chand (xphi,xmuv,xmus,xtau,xrray)
  51. !
  52. ! input parameters: xphi,xmus,xmuv,xtau
  53. ! xphi: azimuthal difference between sun and observation (xphi=0,
  54. ! in backscattering) and expressed in degree (0.:360.)
  55. ! xmus: cosine of the sun zenith angle
  56. ! xmuv: cosine of the observation zenith angle
  57. ! xtau: molecular optical depth
  58. ! output parameter:
  59. ! xrray : molecular reflectance (0.:1.)
  60. ! constant : xdep: depolarization factor (0.0279)
  61. !
  62. real,intent(in) :: xphi,xmuv,xmus,xtau
  63. real,intent(out) :: xrray
  64. !
  65. real :: xdep,pl(10)
  66. real :: fs0,fs1,fs2
  67. real :: as0(10),as1(2),as2(2)
  68. real :: fac,pi,phios,xcosf1,xcosf2
  69. real :: xcosf3,xbeta2,xfd,xph1,xph2,xph3,xitm, xp1, xp2, xp3
  70. real :: cfonc1,cfonc2,cfonc3,xlntau,xitot1,xitot2,xitot3
  71. integer :: i
  72. !
  73. as0 = (/.33243832,-6.777104e-02,.16285370 &
  74. ,1.577425e-03,-.30924818,-1.240906e-02,-.10324388 &
  75. ,3.241678e-02,.11493334,-3.503695e-02/)
  76. as1 = (/.19666292, -5.439061e-02/)
  77. as2 = (/.14545937, -2.910845e-02/)
  78. !
  79. pi=acos(-1.)
  80. fac=pi/180.
  81. phios=180.-xphi
  82. xcosf1=1.
  83. xcosf2=cos(phios*fac)
  84. xcosf3=cos(2*phios*fac)
  85. xbeta2=0.5
  86. xdep=0.0279
  87. xfd=xdep/(2-xdep)
  88. xfd=(1-xfd)/(1+2*xfd)
  89. xph1=1+(3*xmus*xmus-1)*(3*xmuv*xmuv-1)*xfd/8.
  90. xph2=-xmus*xmuv*sqrt(1-xmus*xmus)*sqrt(1-xmuv*xmuv)
  91. xph2=xph2*xfd*xbeta2*1.5
  92. xph3=(1-xmus*xmus)*(1-xmuv*xmuv)
  93. xph3=xph3*xfd*xbeta2*0.375
  94. xitm=(1-exp(-xtau*(1/xmus+1/xmuv)))*xmus/(4*(xmus+xmuv))
  95. xp1=xph1*xitm
  96. xp2=xph2*xitm
  97. xp3=xph3*xitm
  98. xitm=(1-exp(-xtau/xmus))*(1-exp(-xtau/xmuv))
  99. cfonc1=xph1*xitm
  100. cfonc2=xph2*xitm
  101. cfonc3=xph3*xitm
  102. xlntau=log(xtau)
  103. pl(1)=1.
  104. pl(2)=xlntau
  105. pl(3)=xmus+xmuv
  106. pl(4)=xlntau*pl(3)
  107. pl(5)=xmus*xmuv
  108. pl(6)=xlntau*pl(5)
  109. pl(7)=xmus*xmus+xmuv*xmuv
  110. pl(8)=xlntau*pl(7)
  111. pl(9)=xmus*xmus*xmuv*xmuv
  112. pl(10)=xlntau*pl(9)
  113. fs0=0.
  114. do i=1,10
  115. fs0=fs0+pl(i)*as0(i)
  116. enddo
  117. fs1=pl(1)*as1(1)+pl(2)*as1(2)
  118. fs2=pl(1)*as2(1)+pl(2)*as2(2)
  119. xitot1=xp1+cfonc1*fs0*xmus
  120. xitot2=xp2+cfonc2*fs1*xmus
  121. xitot3=xp3+cfonc3*fs2*xmus
  122. xrray=xitot1*xcosf1
  123. xrray=xrray+xitot2*xcosf2*2
  124. xrray=xrray+xitot3*xcosf3*2
  125. xrray=xrray/xmus
  126. !
  127. end subroutine chand
  128. subroutine rweight(wavelx,p_sfc,a_sfc,th0,th,dphix,refl)
  129. !-----------------------------------------------------------------------
  130. ! name :weight
  131. ! description:computes radiance weight as a function of cloud fraction
  132. ! date :18 August 1999
  133. ! version :0.0
  134. ! authors : J.P. Veefkind
  135. !-----------------------------------------------------------------------
  136. ! Variable declaration
  137. ! Input:
  138. ! Wavelength [nm] : wavel0
  139. ! Surface pressure [hPa] : p_sfc
  140. ! Surface albedo : a_sfc
  141. ! Solar zenith angle : th0
  142. ! Viewing zenith angle : th
  143. ! Relative azimuth angle : dphix
  144. ! Output:
  145. ! total reflectance : refl
  146. !-----------------------------------------------------------------------
  147. implicit none
  148. !
  149. real,intent(in) :: wavelx, p_sfc, a_sfc, th0, th, dphix
  150. real,intent(out) :: refl
  151. !
  152. real :: refl_a, t_mu, t_mu0, a_sph, a
  153. real :: mu, mu0, pi, dphi, wavel
  154. real :: tau_sfc, tau
  155. !
  156. pi=acos(-1.)
  157. !-----------------------------------------------------------------------
  158. ! Main
  159. !-----------------------------------------------------------------------
  160. ! Conversion of some parameters
  161. ! wavelength in micrometers:
  162. wavel=wavelx*1e-3
  163. ! compute mu and mu0
  164. mu0=cos(th0/180.*pi)
  165. mu =cos(th /180.*pi)
  166. !
  167. ! use symmetry of dphi
  168. dphi=abs(dphix)
  169. if( dphi > 180.0 ) dphi=360.0-dphi
  170. !
  171. ! redefine dphi
  172. dphi=180.-dphi
  173. !
  174. ! Compute the rayleigh optical depth for surface and cloud top pressure
  175. ! approximation formula by Hansen and Travis
  176. ! tau=0.008569 * wavel**(-4.) * (1. + 0.0113 * wavel**(-2.) + 0.00013 * wavel**(-4.))
  177. tau = (1. + 0.0113/(wavel*wavel) + 0.00013/(wavel*wavel*wavel*wavel))*0.008569/(wavel*wavel*wavel*wavel)
  178. tau_sfc= p_sfc / 1013.25 * tau
  179. tau=tau_sfc
  180. a=a_sfc
  181. ! compute the transmission, see Vermote and Tanre JQSRT 47, pp 305, 1992
  182. t_mu0=( (2./3.+mu0)+(2./3.-mu0)*exp(-tau/mu0) )/( (4./3.)+tau )
  183. t_mu =( (2./3.+mu )+(2./3.-mu )*exp(-tau/mu ) )/( (4./3.)+tau )
  184. ! compute sherical albedo
  185. call csalbr(tau, a_sph)
  186. ! compute atm reflectance
  187. call chand (dphi,mu,mu0,tau,refl_a)
  188. ! compute total reflectance
  189. refl=refl_a + t_mu0 * t_mu* a / (1. - a_sph * a)
  190. ! Write result of RT to screen
  191. ! print*,'dphix, dphi ',dphix, dphi
  192. ! print*,'mu0, mu ',mu0, mu
  193. ! write(*,'(a,f9.3)') 'Transmission mu0 ',t_mu0
  194. ! write(*,'(a,f9.3)') 'Transmission mu0 ',t_mu
  195. ! write(*,'(a,f9.3)') 'Spherical albedo ',a_sph
  196. ! write(*,'(a,f9.3)') 'Atmospheric reflectance ',refl_a
  197. ! write(*,'(a,f9.3)') 'Reflectance ',refl
  198. end subroutine rweight
  199. end module Mrweight