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

Last change on this file since 1477 was 1397, checked in by milmd, 10 years ago

In LMDZ.GENERIC replacement of all phystd .h files by module files.

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