| [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 |
|---|
| [1648] | 28 | use radcommon_h, only: WAVEV, BWNV, DWNV, tstellar, tauray, 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 |
|---|
| [1648] | 40 | real*8 tauvar,tausum,tauwei,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 |
|---|
| [1648] | 54 | if(gfrac(igas,nivref).lt.5.e-2)then |
|---|
| [869] | 55 | print*,'Ignoring ',trim(gnom(igas)),' in Rayleigh scattering '// & |
|---|
| [253] | 56 | 'as its mixing ratio is less than 0.05.' |
|---|
| 57 | ! ignore variable gas in Rayleigh calculation |
|---|
| [305] | 58 | ! ignore gases of mixing ratio < 0.05 in Rayleigh calculation |
|---|
| [253] | 59 | tauconsti(igas) = 0.0 |
|---|
| [135] | 60 | else |
|---|
| [1647] | 61 | if(igas.eq.igas_N2)then |
|---|
| [253] | 62 | tauconsti(igas) = (9.81/g)*8.569E-3*scalep/(P0/93.0) |
|---|
| [869] | 63 | elseif(igas.eq.igas_H2)then |
|---|
| [305] | 64 | tauconsti(igas) = (10.0/g)*0.0241*scalep/101325.0 |
|---|
| 65 | !tauconsti(igas) = (10.0/g)*0.0487*scalep/(101325.0) |
|---|
| [253] | 66 | ! uses n=1.000132 from Optics, Fourth Edition. |
|---|
| [305] | 67 | ! was out by a factor 2! |
|---|
| [253] | 68 | else |
|---|
| 69 | print*,'Error in calc_rayleigh: Gas species not recognised!' |
|---|
| 70 | call abort |
|---|
| 71 | endif |
|---|
| 72 | |
|---|
| [1648] | 73 | if(gfrac(igas,nivref).eq.1.0)then |
|---|
| [869] | 74 | print*,'Rayleigh scattering is for a pure ',trim(gnom(igas)),' atmosphere.' |
|---|
| [253] | 75 | typeknown=.true. |
|---|
| 76 | endif |
|---|
| 77 | |
|---|
| [135] | 78 | endif |
|---|
| [253] | 79 | enddo |
|---|
| 80 | |
|---|
| 81 | if(.not.typeknown)then |
|---|
| 82 | print*,'Rayleigh scattering is for a mixed gas atmosphere.' |
|---|
| 83 | typeknown=.true. |
|---|
| [135] | 84 | endif |
|---|
| [253] | 85 | |
|---|
| [135] | 86 | |
|---|
| 87 | do N=1,L_NSPECTV |
|---|
| 88 | |
|---|
| [253] | 89 | tausum = 0.0 |
|---|
| 90 | tauwei = 0.0 |
|---|
| 91 | bstart = 10000.0/BWNV(N+1) |
|---|
| 92 | bwidth = (10000.0/BWNV(N)) - (10000.0/BWNV(N+1)) |
|---|
| 93 | do ifine=1,Nfine |
|---|
| 94 | wl=bstart+dble(ifine)*bwidth/Nfine |
|---|
| [135] | 95 | |
|---|
| [253] | 96 | tauvar=0.0 |
|---|
| 97 | do igas=1,ngasmx |
|---|
| [1648] | 98 | if(gfrac(igas,nivref).lt.5.e-2)then |
|---|
| [253] | 99 | ! ignore variable gas in Rayleigh calculation |
|---|
| [305] | 100 | ! ignore gases of mixing ratio < 0.05 in Rayleigh calculation |
|---|
| [253] | 101 | tauvari(igas) = 0.0 |
|---|
| 102 | else |
|---|
| [1647] | 103 | if(igas.eq.igas_N2)then |
|---|
| [253] | 104 | tauvari(igas) = (1.0+0.0113/wl**2+0.00013/wl**4)/wl**4 |
|---|
| [869] | 105 | elseif(igas.eq.igas_H2)then |
|---|
| [305] | 106 | tauvari(igas) = 1.0/wl**4 |
|---|
| [253] | 107 | else |
|---|
| 108 | print*,'Error in calc_rayleigh: Gas species not recognised!' |
|---|
| 109 | call abort |
|---|
| 110 | endif |
|---|
| 111 | endif |
|---|
| [135] | 112 | |
|---|
| [1648] | 113 | tauvar=tauvar+tauconsti(igas)*tauvari(igas)*gfrac(igas,nivref) |
|---|
| [253] | 114 | |
|---|
| [135] | 115 | enddo |
|---|
| 116 | |
|---|
| [253] | 117 | call blackl(dble(wl*1e-6),dble(tstellar),df) |
|---|
| 118 | df=df*bwidth/Nfine |
|---|
| 119 | tauwei=tauwei+df |
|---|
| 120 | tausum=tausum+tauvar*df |
|---|
| 121 | |
|---|
| 122 | enddo |
|---|
| 123 | TAURAY(N)=tausum/tauwei |
|---|
| [135] | 124 | ! we multiply by scalep here (100) because plev, which is used in optcv, |
|---|
| 125 | ! is in units of mBar, so we need to convert. |
|---|
| [253] | 126 | |
|---|
| [135] | 127 | end do |
|---|
| 128 | |
|---|
| [848] | 129 | IF (L_NSPECTV > 6) THEN |
|---|
| [861] | 130 | icantbewrong = L_NSPECTV-6 |
|---|
| 131 | print*,'At 1 atm and lambda = ',WAVEV(icantbewrong),' um' |
|---|
| 132 | print*,'tau_R = ',TAURAY(icantbewrong)*1013.25 |
|---|
| 133 | print*,'sig_R = ',TAURAY(icantbewrong)*g*mugaz*1.67e-27*100, & |
|---|
| [848] | 134 | 'cm^2 molecule^-1' |
|---|
| 135 | ENDIF |
|---|
| [253] | 136 | |
|---|
| [135] | 137 | end subroutine calc_rayleigh |
|---|