source: trunk/LMDZ.GENERIC/libf/phystd/calc_rayleigh.F90 @ 1243

Last change on this file since 1243 was 1025, checked in by sglmd, 11 years ago

Helium case added

File size: 5.7 KB
Line 
1      subroutine calc_rayleigh
2
3!==================================================================
4!     
5!     Purpose
6!     -------
7!     Average the Rayleigh scattering in each band, weighting the
8!     average by the blackbody function at temperature tstellar.
9!     Works for an arbitrary mix of gases (currently CO2, N2 and
10!     H2 are possible).
11!     
12!     Authors
13!     -------
14!     Robin Wordsworth (2010)
15!     Jeremy Leconte (2012): Added option for variable gas. Improved water rayleigh (Bucholtz 1995).
16!
17!     Called by
18!     ---------
19!     setspv.F
20!     
21!     Calls
22!     -----
23!     none
24!     
25!==================================================================
26
27      use radinc_h, only: L_NSPECTV
28      use radcommon_h, only: WAVEV, BWNV, DWNV, tstellar, tauray, taurayvar, scalep
29      use gases_h
30
31      implicit none
32
33#include "comcstfi.h"
34#include "callkeys.h"
35
36      real*8 wl
37      integer N,Nfine,ifine,igas
38      parameter(Nfine=500.0)
39      real*8 :: P0 = 9.423D+6   ! Rayleigh scattering reference pressure in pascals.   
40
41      logical typeknown
42      real*8 tauvar,tauvarvar,tausum,tausumvar,tauwei,tauweivar,bwidth,bstart
43      double precision df
44
45      real*8 tauconsti(ngasmx)
46      real*8 tauvari(ngasmx)
47
48      integer icantbewrong
49
50      ! tau0/p0=tau/p (Hansen 1974)
51      ! Then calculate tau0 = (8*pi^3*p_1atm)/(3*N0^2) * 4*delta^2/(g*mugaz*lambda^4)
52      ! where delta=n-1 and N0 is an amagat
53
54      typeknown=.false.
55      do igas=1,ngasmx
56         if(igas.eq.vgas)then
57            print*,'variable gas is ',trim(gnom(igas)),' in Rayleigh scattering '
58         endif
59         if((igas/=vgas).and.(gfrac(igas).lt.5.e-2))then
60            print*,'Ignoring ',trim(gnom(igas)),' in Rayleigh scattering '// &
61            'as its mixing ratio is less than 0.05.'
62            ! ignore variable gas in Rayleigh calculation
63            ! ignore gases of mixing ratio < 0.05 in Rayleigh calculation
64            tauconsti(igas) = 0.0
65         else
66            if(igas.eq.igas_CO2) then
67               tauconsti(igas) = (8.7/g)*1.527*scalep/P0
68            elseif(igas.eq.igas_N2)then
69               tauconsti(igas) = (9.81/g)*8.569E-3*scalep/(P0/93.0)
70            elseif(igas.eq.igas_H2O)then
71!               tauconsti(igas) = (10.0/g)*9.22E-3*scalep/101325.0
72               tauconsti(igas) = 1.98017E-10/(g*mugaz*100.)
73            elseif(igas.eq.igas_H2)then
74               tauconsti(igas) = (10.0/g)*0.0241*scalep/101325.0
75               !tauconsti(igas) = (10.0/g)*0.0487*scalep/(101325.0)
76               ! uses n=1.000132 from Optics, Fourth Edition.
77               ! was out by a factor 2!
78            elseif(igas.eq.igas_He)then
79               tauconsti(igas) = (10.0/g)*0.00086*scalep/101325.0
80            else
81               print*,'Error in calc_rayleigh: Gas species not recognised!'
82               call abort
83            endif
84
85            if((gfrac(igas).eq.1.0).and.(vgas.eq.0))then
86               print*,'Rayleigh scattering is for a pure ',trim(gnom(igas)),' atmosphere.'
87               typeknown=.true.
88            endif
89
90         endif
91      enddo
92
93      if(.not.typeknown)then
94         print*,'Rayleigh scattering is for a mixed gas atmosphere.'
95         typeknown=.true.
96      endif
97
98 
99      do N=1,L_NSPECTV
100         
101         tausum = 0.0
102         tauwei = 0.0
103         tausumvar = 0.0
104         tauweivar = 0.0
105         bstart = 10000.0/BWNV(N+1)
106         bwidth = (10000.0/BWNV(N)) - (10000.0/BWNV(N+1))
107         do ifine=1,Nfine
108            wl=bstart+dble(ifine)*bwidth/Nfine
109
110            tauvar=0.0
111            tauvarvar=0.0
112            do igas=1,ngasmx
113               if((igas/=vgas).and.(gfrac(igas).lt.5.e-2))then
114                  ! ignore variable gas in Rayleigh calculation
115                  ! ignore gases of mixing ratio < 0.05 in Rayleigh calculation
116                  tauvari(igas)   = 0.0
117               else
118                  if(igas.eq.igas_CO2)then
119                     tauvari(igas) = (1.0+0.013/wl**2)/wl**4
120                  elseif(igas.eq.igas_N2)then
121                     tauvari(igas) = (1.0+0.0113/wl**2+0.00013/wl**4)/wl**4
122                  elseif(igas.eq.igas_H2O)then
123!                     tauvari(igas) = 1.0/wl**4 ! to be improved...
124                     tauvari(igas) = (5.7918E6/(2.38E2-1/wl**2)+1.679E5/(57.36E0-1/wl**2))**2/wl**4
125                  elseif(igas.eq.igas_H2)then
126                     tauvari(igas) = 1.0/wl**4
127                  elseif(igas.eq.igas_He)then
128                     tauvari(igas) = 1.0/wl**4
129                  else
130                     print*,'Error in calc_rayleigh: Gas species not recognised!'
131                     call abort
132                  endif
133               endif
134
135               if(igas.eq.vgas) then
136                  tauvarvar=tauvarvar+tauconsti(igas)*tauvari(igas)
137                  tauvar=tauvar+0.0*0.0*gfrac(igas)
138               else
139                  tauvar=tauvar+tauconsti(igas)*tauvari(igas)*gfrac(igas)
140               endif
141
142            enddo
143
144            call blackl(dble(wl*1e-6),dble(tstellar),df)
145            df=df*bwidth/Nfine
146            tauwei=tauwei+df
147            tausum=tausum+tauvar*df
148            tauweivar=tauweivar+df
149            tausumvar=tausumvar+tauvarvar*df
150         
151         enddo
152         TAURAY(N)=tausum/tauwei
153         TAURAYVAR(N)=tausumvar/tauweivar
154         ! we multiply by scalep here (100) because plev, which is used in optcv,
155         ! is in units of mBar, so we need to convert.
156
157      end do
158
159      IF (L_NSPECTV > 6) THEN
160        icantbewrong = L_NSPECTV-6
161        print*,'At 1 atm and lambda = ',WAVEV(icantbewrong),' um'
162        print*,'tau_R = ',TAURAY(icantbewrong)*1013.25
163        print*,'sig_R = ',TAURAY(icantbewrong)*g*mugaz*1.67e-27*100, &
164               'cm^2 molecule^-1'
165      ENDIF
166
167    end subroutine calc_rayleigh
Note: See TracBrowser for help on using the repository browser.