1 | module calc_rayleigh_mod |
---|
2 | |
---|
3 | implicit none |
---|
4 | |
---|
5 | contains |
---|
6 | |
---|
7 | subroutine calc_rayleigh(qvar,muvar,PMID,TMID,tauray) |
---|
8 | |
---|
9 | !================================================================== |
---|
10 | ! |
---|
11 | ! Purpose |
---|
12 | ! ------- |
---|
13 | ! Average the Rayleigh scattering in each band, weighting the |
---|
14 | ! average by the blackbody function at temperature tstellar. |
---|
15 | ! Works for an arbitrary mix of gases. |
---|
16 | ! |
---|
17 | ! Authors |
---|
18 | ! ------- |
---|
19 | ! Robin Wordsworth (2010) |
---|
20 | ! Jeremy Leconte (2012): Added option for variable gas. Improved water rayleigh (Bucholtz 1995). |
---|
21 | ! Noe Clement (2022) : Additionnal comments & Methane+CO Rayleigh |
---|
22 | ! Gwenael Milcareck (2025): Rewriting the code |
---|
23 | ! |
---|
24 | ! Called by |
---|
25 | ! --------- |
---|
26 | ! setspv.F |
---|
27 | ! |
---|
28 | ! Calls |
---|
29 | ! ----- |
---|
30 | ! none |
---|
31 | ! |
---|
32 | !================================================================== |
---|
33 | |
---|
34 | use radinc_h, only: L_NSPECTV, L_LEVELS |
---|
35 | use radcommon_h, only: WAVEV, BWNV, DWNV, tstellar, scalep |
---|
36 | use gases_h, only: ngasmx, vgas, gnom, gfrac, massmol, igas_CO2, igas_H2, & |
---|
37 | igas_H2O, igas_He, igas_N2, igas_CH4, igas_CO, igas_Ar, igas_O2 |
---|
38 | use comcstfi_mod, only: g, pi |
---|
39 | use callkeys_mod, only: strictboundrayleigh |
---|
40 | |
---|
41 | implicit none |
---|
42 | |
---|
43 | real, intent(in) :: qvar(L_LEVELS) ! mol/mol |
---|
44 | real, intent(in) :: muvar(L_LEVELS) ! g/mol |
---|
45 | real, intent(in) :: PMID(L_LEVELS) ! mbar |
---|
46 | real, intent(in) :: TMID(L_LEVELS) ! K |
---|
47 | real, intent(out) :: tauray(L_LEVELS,L_NSPECTV) |
---|
48 | real*8 wl,wn |
---|
49 | integer N,Nfine,ifine,igas,k |
---|
50 | parameter(Nfine=500.0) |
---|
51 | real*8 :: Fk ! King factor for the depolarization |
---|
52 | real*8 :: ng(L_LEVELS) ! real refractive index |
---|
53 | real*8 :: P0(L_LEVELS) ! reference pressure |
---|
54 | real*8 :: T0(L_LEVELS) ! reference temperature |
---|
55 | |
---|
56 | ! Parameters for H2O |
---|
57 | real*8 :: a0, a1, a2, a3, a4, a5, a6, a7 |
---|
58 | real*8 :: luv, lir |
---|
59 | real*8 :: rhor,tr,lr |
---|
60 | real*8 :: rho(L_LEVELS),rhos(L_LEVELS),ts(L_LEVELS) |
---|
61 | real*8 :: b(L_LEVELS) |
---|
62 | |
---|
63 | real*8 mass_frac(ngasmx,L_LEVELS) |
---|
64 | real*8 tauvar(L_LEVELS),tausum(L_LEVELS) |
---|
65 | real*8 tauwei,bwidth,bstart |
---|
66 | double precision df |
---|
67 | |
---|
68 | real*8 tauconsti(ngasmx,L_LEVELS) |
---|
69 | real*8 tauvari(ngasmx,L_LEVELS) |
---|
70 | |
---|
71 | ! Miscellaneous : |
---|
72 | character(len=200) :: message |
---|
73 | character(len=10),parameter :: subname="rayleigh" |
---|
74 | logical, save :: firstcall=.true. |
---|
75 | !$OMP THREADPRIVATE(firstcall) |
---|
76 | |
---|
77 | integer icantbewrong |
---|
78 | |
---|
79 | ! This module calculates the Rayleigh scattering (also known as the Cabannes peak) |
---|
80 | ! Rayleigh wings, Brillouin scattering and Raman scattering are not taken into account. |
---|
81 | |
---|
82 | ! we calculate here TAURAY which is in m2/mBar |
---|
83 | |
---|
84 | ! The cross section for ith particles of small size compared to the wavenumber |
---|
85 | ! and in the electric dipole approximation is: |
---|
86 | ! sigma_i = 24*pi**3*wn**4/N**2 * ((n_i(wn)**2 - 1)/(n_i(wn)**2 + 2))**2 * Fk_i(wn) |
---|
87 | ! nu is the wavenumber |
---|
88 | ! N is the number density of the gas (molecule/m3) |
---|
89 | ! n_i is the real refractive index of the ith gas |
---|
90 | ! Fk_i is the King factor of ith gas equals to (6+3*delta_i)/(6-7*delta_i) |
---|
91 | ! where delta_i is the depolarization factor of the ith gas |
---|
92 | |
---|
93 | ! The rayleigh opacity is expressed by: |
---|
94 | ! tau_r = P/(g*mu) * sum_{i=1}^Ntot [ x_i*sigma_i ] |
---|
95 | ! P is the pressure |
---|
96 | ! g is the standard gravity |
---|
97 | ! mu is the mean molecular weight |
---|
98 | ! x_i is the mass fraction of the ith gas |
---|
99 | ! The pressure P dependence is calculated in optcv.F90 |
---|
100 | |
---|
101 | if(firstcall) then |
---|
102 | |
---|
103 | if ((BWNV(L_NSPECTV+1).gt.60000.).and.(strictboundrayleigh)) then |
---|
104 | message="Rayleigh scattering is unknown for wn>60000 cm-1 - some singularities are present for higher wavenumber - if you know what you are doing, use strictboundrayleigh=.false." |
---|
105 | call abort_physic(subname,message,1) |
---|
106 | elseif ((BWNV(L_NSPECTV+1).gt.60000.).and..not.(strictboundrayleigh)) then |
---|
107 | print*,'**********************************************' |
---|
108 | print*,' we allow model to continue with wn>60000 cm-1' |
---|
109 | print*,' ... we assume we know what you are doing ... ' |
---|
110 | print*,' ... but do not let this happen too often ... ' |
---|
111 | print*,'**********************************************' |
---|
112 | endif |
---|
113 | firstcall = .false. |
---|
114 | endif |
---|
115 | |
---|
116 | |
---|
117 | do igas=1,ngasmx |
---|
118 | ! Convert qvar mol/mol -> kg/kg |
---|
119 | if((igas.eq.vgas).and.(maxval(QVAR(:)).ge.1.e-2))then |
---|
120 | ! print*,'variable gas is ',trim(gnom(igas)),' in Rayleigh scattering ' |
---|
121 | mass_frac(igas,:) = QVAR(:)*massmol(igas)/muvar(:) |
---|
122 | elseif((igas/=vgas).and.(gfrac(igas).ge.1.e-2))then |
---|
123 | mass_frac(igas,:) = gfrac(igas)*(1.-QVAR(:))*massmol(igas)/muvar(:) |
---|
124 | else |
---|
125 | ! print*,'Ignoring ',trim(gnom(igas)),' in Rayleigh scattering '// & |
---|
126 | ! 'as its mixing ratio is less than 0.01.' |
---|
127 | ! ignore variable gas in Rayleigh calculation |
---|
128 | ! ignore gases of mixing ratio < 0.01 in Rayleigh calculation |
---|
129 | mass_frac(igas,:) = 0.0 |
---|
130 | endif |
---|
131 | tauvari(igas,:) = 0. |
---|
132 | enddo |
---|
133 | |
---|
134 | ! WARNING, beyond 60000 cm-1, for all molecules, there are singularities due to the interpolation formula. |
---|
135 | |
---|
136 | |
---|
137 | do N=1,L_NSPECTV |
---|
138 | |
---|
139 | ! The refractive index depend on temperature and pressure |
---|
140 | ! It isn't the case here. Must be implemented in the future... |
---|
141 | ! But in the current scientific litterature (2024), it's difficult |
---|
142 | ! to find something that depends on temperature and pressure... |
---|
143 | ! except for H2O |
---|
144 | |
---|
145 | tausum = 0.0 |
---|
146 | tauwei = 0.0 |
---|
147 | bstart = 10000.0/BWNV(N+1) ! BWNV is in cm-1 so 10000.0/BWNV is in micron |
---|
148 | bwidth = (10000.0/BWNV(N)) - (10000.0/BWNV(N+1)) |
---|
149 | do ifine=1,Nfine |
---|
150 | wl=bstart+dble(ifine)*bwidth/Nfine |
---|
151 | wn=BWNV(N)+dble(ifine)*(BWNV(N+1)-BWNV(N))/Nfine |
---|
152 | |
---|
153 | tauvar(:)=0.0 |
---|
154 | do igas=1,ngasmx |
---|
155 | if (maxval(mass_frac(igas,:)).ge.1e-2) then |
---|
156 | |
---|
157 | if(igas.eq.igas_CO2)then |
---|
158 | ! Sneep et al, 2005 |
---|
159 | ! doi:10.1016/j.jqsrt.2004.07.025 |
---|
160 | T0(:) = 288.15 |
---|
161 | P0(:) = 1.01325e5 |
---|
162 | ! ng -> valid range of the measurements : 0.1807 - 1.8172 um |
---|
163 | ng(:) = 1. + 1.1427e3*(5799.25/(128908.9**2 - wn**2) + 120.05/(89223.8**2 - wn**2) + 5.3334/(75037.5**2 - wn**2) + 4.3244/(67837.7**2 - wn**2) + 0.1218145e-4/(2418.136**2 - wn**2)) ! there is an error on the paper 1.1427e6 -> 1.1427e3 |
---|
164 | Fk = 1.1364 + 25.3e-12*wn**2 |
---|
165 | tauvari(igas,:) = mass_frac(igas,:)*((ng(:)**2-1.)/(ng(:)**2+2.))**2 * Fk * (wn*100.)**4 ! wn*100 -> cm-1 to m-1 |
---|
166 | elseif(igas.eq.igas_N2)then |
---|
167 | ! Thalman et al, 2014 |
---|
168 | ! doi:10.1016/j.jqsrt.2014.05.030 |
---|
169 | T0(:) = 288.15 |
---|
170 | P0(:) = 1.01325e5 |
---|
171 | if(wn.gt.21360)then !between 21360 and 39370 cm-1. We extrapolate above. |
---|
172 | ng(:) = 1. + (5677.465 + 318.81874e13/(14.4e9 - wn**2))*1e-8 !there is an error on the paper e12 -> e13 |
---|
173 | else !between 4860 and 21360 cm-1. We extrapolate below. |
---|
174 | ng(:) = 1. + (6498.2 + 307.4335e13/(14.4e9 - wn**2))*1e-8 |
---|
175 | endif |
---|
176 | Fk = 1.034 + 3.17e-12*wn |
---|
177 | tauvari(igas,:) = mass_frac(igas,:)*((ng(:)**2-1.)/(ng(:)**2+2.))**2 * Fk * (wn*100.)**4 |
---|
178 | elseif(igas.eq.igas_H2O)then |
---|
179 | ! Harvey et al, 1998 |
---|
180 | ! doi:10.1063/1.556029 |
---|
181 | ! doi:10.1364/AO.35.001566 for wn<4840 cm-1 |
---|
182 | a0 = 0.244257733 |
---|
183 | a1 = 9.74634476e-3 |
---|
184 | a2 = -3.73234996e-3 |
---|
185 | a3 = 2.68678472e-4 |
---|
186 | a4 = 1.58920570e-3 |
---|
187 | a5 = 2.45934259e-3 |
---|
188 | a6 = 0.900704920 |
---|
189 | a7 = -1.66626219e-2 |
---|
190 | luv = 0.2292020 |
---|
191 | lir = 5.432937 |
---|
192 | Tr = 273.15 |
---|
193 | rhor = 1000. |
---|
194 | T0(:) = tmid(:) |
---|
195 | P0(:) = pmid(:)*scalep |
---|
196 | lr = 0.589 |
---|
197 | rho(:) = mass_frac(igas,:)*muvar(:)/massmol(igas)*P0(:)/(8.314463*T0(:)/(muvar(:)/1000.)) |
---|
198 | rhos(:) = rho(:)/rhor |
---|
199 | ts(:) = T0(:)/Tr |
---|
200 | |
---|
201 | b(:) = (a0 + a1*rhos(:) + a2*ts(:) + a3*ts(:)*(10000./wn/lr)**2 + a4/(10000./wn/lr)**2 + a5/((10000./wn/lr)**2 - luv**2) + a6/((10000./wn/lr)**2 - lir**2) + a7*rhos(:)**2)*rhos(:) |
---|
202 | ! ng -> valid range of the measurements : 0.2 - 1.1 um |
---|
203 | ng(:) = sqrt(2.*b(:)+1.)/sqrt(1.-b(:)) |
---|
204 | if(wn<4840.) then ! necessary to prevent a singularity at 3230 cm-1 |
---|
205 | ng(:) = 1. + 1.022e-8*(295.235 + 2.6422*(wn*1e-4)**2 - 0.032380*(wn*1e-4)**4 + 0.004028*(wn*1e-4)**6) |
---|
206 | endif |
---|
207 | Fk = (6.+3.*3e-4)/(6.-7.*3e-4) ! delta=3e-4 Murphy 1977 doi:10.1063/1.434794 |
---|
208 | tauvari(igas,:) = mass_frac(igas,:)*((ng(:)**2-1.)/(ng(:)**2+2.))**2 * Fk * (wn*100.)**4 |
---|
209 | elseif(igas.eq.igas_H2)then |
---|
210 | ! Peck and Hung, 1977 |
---|
211 | ! doi:10.1364/JOSA.67.001550 |
---|
212 | T0(:) = 273.15 |
---|
213 | P0(:) = 1.01325e5 |
---|
214 | ! ng -> valid range of the measurements : 0.1680 - 1.6945 um |
---|
215 | ng(:) = 1. + (14895.6/(180.7 - (wn*1e-4)**2) + 4903.7/(92.-(wn*1e-4)**2))*1e-6 |
---|
216 | Fk = (6.+3.*0.02)/(6.-7.*0.02) ! delta=0.02 Hansen 1974 |
---|
217 | tauvari(igas,:) = mass_frac(igas,:)*((ng(:)**2-1.)/(ng(:)**2+2.))**2 * Fk * (wn*100.)**4 |
---|
218 | elseif(igas.eq.igas_He)then |
---|
219 | ! Thalman et al, 2014 |
---|
220 | ! doi:10.1016/j.jqsrt.2014.05.030 |
---|
221 | T0(:) = 288.15 |
---|
222 | P0(:) = 1.01325e5 |
---|
223 | ! ng -> valid range of the measurements : 0.2753 - 20.5813 um |
---|
224 | ng(:) = 1. + (2283. + 1.8102e13/(1.5342e10 - wn**2))*1e-8 |
---|
225 | Fk = 1. |
---|
226 | tauvari(igas,:) = mass_frac(igas,:)*((ng(:)**2-1.)/(ng(:)**2+2.))**2 * Fk * (wn*100.)**4 |
---|
227 | elseif(igas.eq.igas_CH4)then |
---|
228 | ! Sneep et al, 2005 |
---|
229 | ! doi:10.1016/j.jqsrt.2004.07.025 |
---|
230 | T0(:) = 288.15 |
---|
231 | P0(:) = 1.01325e5 |
---|
232 | ! ng -> valid range of the measurements : 0.3251 - 0.6330 um |
---|
233 | ng(:) = 1. + 46662e-8 + 4.02e-14*wn**2 |
---|
234 | Fk = 1. |
---|
235 | tauvari(igas,:) = mass_frac(igas,:)*((ng(:)**2-1.)/(ng(:)**2+2.))**2 * Fk * (wn*100.)**4 |
---|
236 | elseif(igas.eq.igas_CO)then |
---|
237 | ! Sneep et al, 2005 |
---|
238 | ! doi:10.1016/j.jqsrt.2004.07.025 |
---|
239 | T0(:) = 288.15 |
---|
240 | P0(:) = 1.01325e5 |
---|
241 | ! ng -> valid range of the measurements : 0.168 - 0.288 um |
---|
242 | ng(:) = 1. + 22851e-8 + 0.456e4/(71427.**2 - wn**2) |
---|
243 | Fk = 1.016 |
---|
244 | tauvari(igas,:) = mass_frac(igas,:)*((ng(:)**2-1.)/(ng(:)**2+2.))**2 * Fk * (wn*100.)**4 |
---|
245 | elseif(igas.eq.igas_Ar)then |
---|
246 | ! Thalman et al, 2014 |
---|
247 | ! doi:10.1016/j.jqsrt.2014.05.030 |
---|
248 | T0(:) = 288.15 |
---|
249 | P0(:) = 1.01325e5 |
---|
250 | ! ng -> valid range of the measurements : 0.288 - 0.546 um |
---|
251 | ng(:) = 1. + (6432.135 + 286.06021e12/(14.4e9 - wn**2))*1e-8 |
---|
252 | Fk = 1. |
---|
253 | tauvari(igas,:) = mass_frac(igas,:)*((ng(:)**2-1.)/(ng(:)**2+2.))**2 * Fk * (wn*100.)**4 |
---|
254 | elseif(igas.eq.igas_O2)then |
---|
255 | ! Thalman et al, 2014 |
---|
256 | ! doi:10.1016/j.jqsrt.2014.05.030 |
---|
257 | T0(:) = 273.15 |
---|
258 | P0(:) = 1.01325e5 |
---|
259 | ! ng -> valid range of the measurements : 0.303 - 2.0 um |
---|
260 | ng(:) = 1. + (20564.8 + 2.480899e13/(4.09e9 - wn**2))*1e-8 |
---|
261 | Fk = 1.09 + 1.385e-11*wn**2 + 1.488e-20*wn**4 |
---|
262 | tauvari(igas,:) = mass_frac(igas,:)*((ng(:)**2-1.)/(ng(:)**2+2.))**2 * Fk * (wn*100.)**4 |
---|
263 | else |
---|
264 | print*,'No rayleigh scattering for ',trim(gnom(igas)),'. No data found.' |
---|
265 | endif |
---|
266 | |
---|
267 | tauconsti(igas,:) = 24.*pi**3 *6.022141E+023 / (g*(muvar(:)/1000.)*(P0(:)/(1.380649E-23*T0(:)))**2) |
---|
268 | ! N=P/(kB*T) |
---|
269 | ! pmid*scalep -> mbar to Pa |
---|
270 | ! muvar/1000 -> g/mol to kg/mol |
---|
271 | |
---|
272 | tauvar(:)=tauvar(:)+tauconsti(igas,:)*tauvari(igas,:) |
---|
273 | |
---|
274 | endif !greater than 0.01 |
---|
275 | |
---|
276 | enddo !ngasmx |
---|
277 | |
---|
278 | call blackn(dble(wn*1e2),dble(tstellar),df) |
---|
279 | df=df*(BWNV(N+1)-BWNV(N))/Nfine |
---|
280 | tauwei=tauwei+df |
---|
281 | tausum(:)=tausum(:)+tauvar(:)*df |
---|
282 | |
---|
283 | enddo !Nfine |
---|
284 | TAURAY(:,N)=tausum(:)/tauwei |
---|
285 | |
---|
286 | end do !L_NSPECTV |
---|
287 | |
---|
288 | |
---|
289 | end subroutine calc_rayleigh |
---|
290 | |
---|
291 | end module calc_rayleigh_mod |
---|