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 N2, H2 |
---|
10 | ! and CH4 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 | ! B. de Batz de Trenquelleon (2023) : Added CH4 rayleigh. |
---|
17 | ! |
---|
18 | ! Called by |
---|
19 | ! --------- |
---|
20 | ! setspv.F |
---|
21 | ! |
---|
22 | ! Calls |
---|
23 | ! ----- |
---|
24 | ! none |
---|
25 | ! |
---|
26 | !================================================================== |
---|
27 | |
---|
28 | use radinc_h, only: L_NSPECTV |
---|
29 | use radcommon_h, only: WAVEV, BWNV, DWNV, tstellar, tauray, scalep |
---|
30 | use gases_h |
---|
31 | use comcstfi_mod, only: g, mugaz, pi |
---|
32 | |
---|
33 | implicit none |
---|
34 | |
---|
35 | real*8 wl |
---|
36 | integer N,Nfine,ifine,igas |
---|
37 | parameter(Nfine=500.0) |
---|
38 | real*8 :: P0 = 9.423D+6 ! Rayleigh scattering reference pressure in pascals. |
---|
39 | |
---|
40 | logical typeknown |
---|
41 | real*8 tauvar,tausum,tauwei,bwidth,bstart |
---|
42 | double precision df |
---|
43 | |
---|
44 | real*8 tauconsti(ngasmx) |
---|
45 | real*8 tauvari(ngasmx) |
---|
46 | |
---|
47 | integer icantbewrong |
---|
48 | |
---|
49 | ! we multiply by scalep here (100) because plev, which is used in optcv, |
---|
50 | ! is in units of mBar, so we need to convert. |
---|
51 | ! we calculate here TAURAY which is in m2/mBar |
---|
52 | |
---|
53 | ! tau0/p0=tau/p (Hansen 1974 : equation (2.31) ) |
---|
54 | ! Then calculate tau0 = (8*pi^3*p_1atm)/(3*N0^2) * 4*delta^2/(g*mugaz*lambda^4) (Hansen 1974 : equation (2.29) ) |
---|
55 | ! where delta=n-1 and N0 is an amagat |
---|
56 | |
---|
57 | ! In exo_k : (8*pi^3)/(3*N0^2) * 4 * delta^2 ------- is written : (24*pi^3)/(N0^2) * (4/9) * delta^2 |
---|
58 | |
---|
59 | ! (tauconsti * tauvari) = scalep * cross_section_molecule_exo_k / ( gravity * molecular_mass ) |
---|
60 | ! (tauconsti * tauvari) in m2/mBar because of scalep factor |
---|
61 | ! cross_section_molecule in m2/molecule |
---|
62 | ! molecular_mass is the mass of th considered molecule |
---|
63 | |
---|
64 | typeknown=.false. |
---|
65 | do igas=1,ngasmx |
---|
66 | if(gfrac(igas,nivref).lt.1.e-2)then |
---|
67 | print*,'Ignoring ',trim(gnom(igas)),' in Rayleigh scattering '// & |
---|
68 | 'as its mixing ratio is less than 0.01.' |
---|
69 | ! ignore variable gas in Rayleigh calculation |
---|
70 | ! ignore gases of mixing ratio < 0.01 in Rayleigh calculation |
---|
71 | tauconsti(igas) = 0.0 |
---|
72 | else |
---|
73 | if(igas.eq.igas_N2)then |
---|
74 | tauconsti(igas) = (9.81/g)*8.569E-3*scalep/(P0/93.0) |
---|
75 | elseif(igas.eq.igas_H2)then |
---|
76 | tauconsti(igas) = (10.0/g)*0.0241*scalep/101325.0 |
---|
77 | !tauconsti(igas) = (10.0/g)*0.0487*scalep/(101325.0) |
---|
78 | ! uses n=1.000132 from Optics, Fourth Edition. |
---|
79 | ! was out by a factor 2! |
---|
80 | elseif(igas.eq.igas_CH4)then |
---|
81 | ! For CH4 we use the formalism of exo_k (developed by Jeremy Leconte) |
---|
82 | ! the relation between exo_k formalism and LMDZ formalism is : (tauconsti*tauvari)=scalep*sigma_exok/(gravity*molecular_mass) |
---|
83 | ! Values are taken from Caldas (2019), equation 12 & appendix D |
---|
84 | ! 24.*pi**3/(1E5/(1.380648813E-23*273.15))**2 is a constant |
---|
85 | ! 4/9=(2/3)**2 is approximately ((n+1)/(n**2+2))**2 |
---|
86 | ! 1E24 = (1E6)**4 because wavelength is in micron |
---|
87 | ! 16.04*1.6726*1E-27 is CH4 molecular mass |
---|
88 | tauconsti(igas) = 24.*pi**3/(1E5/(1.380648813E-23*273.15))**2 * 4./9. * 1E24 * (1/( g *16.04*1.6726*1E-27)) * scalep |
---|
89 | else |
---|
90 | print*,'Error in calc_rayleigh: Gas species not recognised!' |
---|
91 | call abort |
---|
92 | endif |
---|
93 | |
---|
94 | if(gfrac(igas,nivref).eq.1.0)then |
---|
95 | print*,'Rayleigh scattering is for a pure ',trim(gnom(igas)),' atmosphere.' |
---|
96 | typeknown=.true. |
---|
97 | endif |
---|
98 | |
---|
99 | endif |
---|
100 | enddo |
---|
101 | |
---|
102 | if(.not.typeknown)then |
---|
103 | print*,'Rayleigh scattering is for a mixed gas atmosphere.' |
---|
104 | typeknown=.true. |
---|
105 | endif |
---|
106 | |
---|
107 | |
---|
108 | do N=1,L_NSPECTV |
---|
109 | |
---|
110 | tausum = 0.0 |
---|
111 | tauwei = 0.0 |
---|
112 | bstart = 10000.0/BWNV(N+1) |
---|
113 | bwidth = (10000.0/BWNV(N)) - (10000.0/BWNV(N+1)) |
---|
114 | do ifine=1,Nfine |
---|
115 | wl=bstart+dble(ifine)*bwidth/Nfine |
---|
116 | |
---|
117 | tauvar=0.0 |
---|
118 | do igas=1,ngasmx |
---|
119 | if(gfrac(igas,nivref).lt.1.e-2)then |
---|
120 | ! ignore variable gas in Rayleigh calculation |
---|
121 | ! ignore gases of mixing ratio < 0.01 in Rayleigh calculation |
---|
122 | tauvari(igas) = 0.0 |
---|
123 | else |
---|
124 | if(igas.eq.igas_N2)then |
---|
125 | tauvari(igas) = (1.0+0.0113/wl**2+0.00013/wl**4)/wl**4 |
---|
126 | elseif(igas.eq.igas_H2)then |
---|
127 | tauvari(igas) = 1.0/wl**4 |
---|
128 | elseif(igas.eq.igas_CH4)then |
---|
129 | tauvari(igas) = (4.6662E-4+4.02E-6/wl**2)**2 / wl**4 |
---|
130 | else |
---|
131 | print*,'Error in calc_rayleigh: Gas species not recognised!' |
---|
132 | call abort |
---|
133 | endif |
---|
134 | endif |
---|
135 | |
---|
136 | tauvar=tauvar+tauconsti(igas)*tauvari(igas)*gfrac(igas,nivref) |
---|
137 | |
---|
138 | enddo |
---|
139 | |
---|
140 | call blackl(dble(wl*1e-6),dble(tstellar),df) |
---|
141 | df=df*bwidth/Nfine |
---|
142 | tauwei=tauwei+df |
---|
143 | tausum=tausum+tauvar*df |
---|
144 | |
---|
145 | enddo |
---|
146 | TAURAY(N)=tausum/tauwei |
---|
147 | ! we multiply by scalep here (100) because plev, which is used in optcv, |
---|
148 | ! is in units of mBar, so we need to convert. |
---|
149 | |
---|
150 | end do |
---|
151 | |
---|
152 | IF (L_NSPECTV > 6) THEN |
---|
153 | icantbewrong = L_NSPECTV-6 |
---|
154 | print*,'At 1 atm and lambda = ',WAVEV(icantbewrong),' um' |
---|
155 | print*,'tau_R = ',TAURAY(icantbewrong)*1013.25 |
---|
156 | print*,'sig_R = ',TAURAY(icantbewrong)*g*mugaz*1.67e-27*100, & |
---|
157 | 'cm^2 molecule^-1' |
---|
158 | ENDIF |
---|
159 | |
---|
160 | end subroutine calc_rayleigh |
---|