[135] | 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. |
---|
[253] | 9 | ! Works for an arbitrary mix of gases (currently CO2, N2 and |
---|
| 10 | ! H2 are possible). |
---|
[135] | 11 | ! |
---|
| 12 | ! Authors |
---|
| 13 | ! ------- |
---|
| 14 | ! Robin Wordsworth (2010) |
---|
[1016] | 15 | ! Jeremy Leconte (2012): Added option for variable gas. Improved water rayleigh (Bucholtz 1995). |
---|
[135] | 16 | ! |
---|
| 17 | ! Called by |
---|
| 18 | ! --------- |
---|
| 19 | ! setspv.F |
---|
| 20 | ! |
---|
| 21 | ! Calls |
---|
| 22 | ! ----- |
---|
| 23 | ! none |
---|
| 24 | ! |
---|
| 25 | !================================================================== |
---|
| 26 | |
---|
| 27 | use radinc_h, only: L_NSPECTV |
---|
[1016] | 28 | use radcommon_h, only: WAVEV, BWNV, DWNV, tstellar, tauray, taurayvar, scalep |
---|
[471] | 29 | use gases_h |
---|
[1384] | 30 | use comcstfi_mod, only: g, mugaz |
---|
[135] | 31 | |
---|
| 32 | implicit none |
---|
| 33 | |
---|
| 34 | real*8 wl |
---|
[253] | 35 | integer N,Nfine,ifine,igas |
---|
[135] | 36 | parameter(Nfine=500.0) |
---|
| 37 | real*8 :: P0 = 9.423D+6 ! Rayleigh scattering reference pressure in pascals. |
---|
| 38 | |
---|
[253] | 39 | logical typeknown |
---|
[1016] | 40 | real*8 tauvar,tauvarvar,tausum,tausumvar,tauwei,tauweivar,bwidth,bstart |
---|
[135] | 41 | double precision df |
---|
| 42 | |
---|
[253] | 43 | real*8 tauconsti(ngasmx) |
---|
| 44 | real*8 tauvari(ngasmx) |
---|
[135] | 45 | |
---|
[861] | 46 | integer icantbewrong |
---|
| 47 | |
---|
[253] | 48 | ! tau0/p0=tau/p (Hansen 1974) |
---|
[305] | 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 |
---|
[253] | 51 | |
---|
| 52 | typeknown=.false. |
---|
| 53 | do igas=1,ngasmx |
---|
| 54 | if(igas.eq.vgas)then |
---|
[1016] | 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 |
---|
[869] | 58 | print*,'Ignoring ',trim(gnom(igas)),' in Rayleigh scattering '// & |
---|
[253] | 59 | 'as its mixing ratio is less than 0.05.' |
---|
| 60 | ! ignore variable gas in Rayleigh calculation |
---|
[305] | 61 | ! ignore gases of mixing ratio < 0.05 in Rayleigh calculation |
---|
[253] | 62 | tauconsti(igas) = 0.0 |
---|
[135] | 63 | else |
---|
[869] | 64 | if(igas.eq.igas_CO2) then |
---|
[253] | 65 | tauconsti(igas) = (8.7/g)*1.527*scalep/P0 |
---|
[869] | 66 | elseif(igas.eq.igas_N2)then |
---|
[253] | 67 | tauconsti(igas) = (9.81/g)*8.569E-3*scalep/(P0/93.0) |
---|
[869] | 68 | elseif(igas.eq.igas_H2O)then |
---|
[1016] | 69 | ! tauconsti(igas) = (10.0/g)*9.22E-3*scalep/101325.0 |
---|
| 70 | tauconsti(igas) = 1.98017E-10/(g*mugaz*100.) |
---|
[869] | 71 | elseif(igas.eq.igas_H2)then |
---|
[305] | 72 | tauconsti(igas) = (10.0/g)*0.0241*scalep/101325.0 |
---|
| 73 | !tauconsti(igas) = (10.0/g)*0.0487*scalep/(101325.0) |
---|
[253] | 74 | ! uses n=1.000132 from Optics, Fourth Edition. |
---|
[305] | 75 | ! was out by a factor 2! |
---|
[869] | 76 | elseif(igas.eq.igas_He)then |
---|
[1025] | 77 | tauconsti(igas) = (10.0/g)*0.00086*scalep/101325.0 |
---|
[253] | 78 | else |
---|
| 79 | print*,'Error in calc_rayleigh: Gas species not recognised!' |
---|
| 80 | call abort |
---|
| 81 | endif |
---|
| 82 | |
---|
[1016] | 83 | if((gfrac(igas).eq.1.0).and.(vgas.eq.0))then |
---|
[869] | 84 | print*,'Rayleigh scattering is for a pure ',trim(gnom(igas)),' atmosphere.' |
---|
[253] | 85 | typeknown=.true. |
---|
| 86 | endif |
---|
| 87 | |
---|
[135] | 88 | endif |
---|
[253] | 89 | enddo |
---|
| 90 | |
---|
| 91 | if(.not.typeknown)then |
---|
| 92 | print*,'Rayleigh scattering is for a mixed gas atmosphere.' |
---|
| 93 | typeknown=.true. |
---|
[135] | 94 | endif |
---|
[253] | 95 | |
---|
[135] | 96 | |
---|
| 97 | do N=1,L_NSPECTV |
---|
| 98 | |
---|
[253] | 99 | tausum = 0.0 |
---|
| 100 | tauwei = 0.0 |
---|
[1016] | 101 | tausumvar = 0.0 |
---|
| 102 | tauweivar = 0.0 |
---|
[253] | 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 |
---|
[135] | 107 | |
---|
[253] | 108 | tauvar=0.0 |
---|
[1016] | 109 | tauvarvar=0.0 |
---|
[253] | 110 | do igas=1,ngasmx |
---|
[1016] | 111 | if((igas/=vgas).and.(gfrac(igas).lt.5.e-2))then |
---|
[253] | 112 | ! ignore variable gas in Rayleigh calculation |
---|
[305] | 113 | ! ignore gases of mixing ratio < 0.05 in Rayleigh calculation |
---|
[253] | 114 | tauvari(igas) = 0.0 |
---|
| 115 | else |
---|
[869] | 116 | if(igas.eq.igas_CO2)then |
---|
[253] | 117 | tauvari(igas) = (1.0+0.013/wl**2)/wl**4 |
---|
[869] | 118 | elseif(igas.eq.igas_N2)then |
---|
[253] | 119 | tauvari(igas) = (1.0+0.0113/wl**2+0.00013/wl**4)/wl**4 |
---|
[869] | 120 | elseif(igas.eq.igas_H2O)then |
---|
[1016] | 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 |
---|
[869] | 123 | elseif(igas.eq.igas_H2)then |
---|
[305] | 124 | tauvari(igas) = 1.0/wl**4 |
---|
[869] | 125 | elseif(igas.eq.igas_He)then |
---|
[1025] | 126 | tauvari(igas) = 1.0/wl**4 |
---|
[253] | 127 | else |
---|
| 128 | print*,'Error in calc_rayleigh: Gas species not recognised!' |
---|
| 129 | call abort |
---|
| 130 | endif |
---|
| 131 | endif |
---|
[135] | 132 | |
---|
[1016] | 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 |
---|
[253] | 139 | |
---|
[135] | 140 | enddo |
---|
| 141 | |
---|
[253] | 142 | call blackl(dble(wl*1e-6),dble(tstellar),df) |
---|
| 143 | df=df*bwidth/Nfine |
---|
| 144 | tauwei=tauwei+df |
---|
| 145 | tausum=tausum+tauvar*df |
---|
[1016] | 146 | tauweivar=tauweivar+df |
---|
| 147 | tausumvar=tausumvar+tauvarvar*df |
---|
[253] | 148 | |
---|
| 149 | enddo |
---|
| 150 | TAURAY(N)=tausum/tauwei |
---|
[1016] | 151 | TAURAYVAR(N)=tausumvar/tauweivar |
---|
[135] | 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. |
---|
[253] | 154 | |
---|
[135] | 155 | end do |
---|
| 156 | |
---|
[848] | 157 | IF (L_NSPECTV > 6) THEN |
---|
[861] | 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, & |
---|
[848] | 162 | 'cm^2 molecule^-1' |
---|
| 163 | ENDIF |
---|
[253] | 164 | |
---|
[135] | 165 | end subroutine calc_rayleigh |
---|