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

Last change on this file since 868 was 861, checked in by aslmd, 13 years ago

LMDZ.GENERIC

  • Fixed an allocating bug which arises from previous modifications owing to the double use of callcorrk with CLF_varying
  • Fixed a small bug in a diagnostic in the end of calc_rayleigh. Some picky compilers complain.
  • Fixed a small bug with the array noms which is not allocated when tracer is false. But still need in physdem1.
File size: 5.1 KB
RevLine 
[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)
15!
16!     Called by
17!     ---------
18!     setspv.F
19!     
20!     Calls
21!     -----
22!     none
23!     
24!==================================================================
25
26      use radinc_h, only: L_NSPECTV
[253]27      use radcommon_h, only: WAVEV, BWNV, DWNV, tstellar, tauray, scalep
[471]28      use gases_h
[135]29
30      implicit none
31
32#include "comcstfi.h"
[253]33#include "callkeys.h"
[135]34
35      real*8 wl
[253]36      integer N,Nfine,ifine,igas
[135]37      parameter(Nfine=500.0)
38      real*8 :: P0 = 9.423D+6   ! Rayleigh scattering reference pressure in pascals.   
39
[253]40      logical typeknown
41      real*8 tauvar,tausum,tauwei,bwidth,bstart
[135]42      double precision df
43
[253]44      real*8 tauconsti(ngasmx)
45      real*8 tauvari(ngasmx)
[135]46
[861]47      integer icantbewrong
48
[253]49      ! tau0/p0=tau/p (Hansen 1974)
[305]50      ! Then calculate tau0 = (8*pi^3*p_1atm)/(3*N0^2) * 4*delta^2/(g*mugaz*lambda^4)
51      ! where delta=n-1 and N0 is an amagat
[253]52
53      typeknown=.false.
54      do igas=1,ngasmx
55         if(igas.eq.vgas)then
56            print*,'Ignoring ',gnom(igas),' in Rayleigh scattering '// &
57            'as it is variable.'
58         elseif(gfrac(igas).lt.5.e-2)then
59            print*,'Ignoring ',gnom(igas),' in Rayleigh scattering '// &
60            'as its mixing ratio is less than 0.05.'
61            ! ignore variable gas in Rayleigh calculation
[305]62            ! ignore gases of mixing ratio < 0.05 in Rayleigh calculation
[253]63            tauconsti(igas) = 0.0
[135]64         else
[253]65            if(gnom(igas).eq.'CO2')then
66               tauconsti(igas) = (8.7/g)*1.527*scalep/P0
67            elseif(gnom(igas).eq.'N2_')then
68               tauconsti(igas) = (9.81/g)*8.569E-3*scalep/(P0/93.0)
[305]69            elseif(gnom(igas).eq.'H2O')then
70               tauconsti(igas) = (10.0/g)*9.22E-3*scalep/101325.0
[253]71            elseif(gnom(igas).eq.'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!
[253]76            elseif(gnom(igas).eq.'He_')then
77               print*,'Helium not ready yet!'
[716]78               tauconsti(igas) = 0.0
79               call abort
[253]80            else
81               print*,'Error in calc_rayleigh: Gas species not recognised!'
82               call abort
83            endif
84
85            if(gfrac(igas).eq.1.0)then
86               print*,'Rayleigh scattering is for a pure ',gnom(igas),' atmosphere.'
87               typeknown=.true.
88            endif
89
[135]90         endif
[253]91      enddo
92
93      if(.not.typeknown)then
94         print*,'Rayleigh scattering is for a mixed gas atmosphere.'
95         typeknown=.true.
[135]96      endif
[253]97
[135]98 
99      do N=1,L_NSPECTV
100         
[253]101         tausum = 0.0
102         tauwei = 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
[135]107
[253]108            tauvar=0.0
109            do igas=1,ngasmx
110               if((igas.eq.vgas).or.(gfrac(igas).lt.5.e-2))then
111                  ! ignore variable gas in Rayleigh calculation
[305]112                  ! ignore gases of mixing ratio < 0.05 in Rayleigh calculation
[253]113                  tauvari(igas)   = 0.0
114               else
115                  if(gnom(igas).eq.'CO2')then
116                     tauvari(igas) = (1.0+0.013/wl**2)/wl**4
117                  elseif(gnom(igas).eq.'N2_')then
118                     tauvari(igas) = (1.0+0.0113/wl**2+0.00013/wl**4)/wl**4
[305]119                  elseif(gnom(igas).eq.'H2O')then
120                     tauvari(igas) = 1.0/wl**4 ! to be improved...
[253]121                  elseif(gnom(igas).eq.'H2_')then
[305]122                     tauvari(igas) = 1.0/wl**4
[716]123                  elseif(gnom(igas).eq.'He_')then
124                     print*,'Helium not ready yet!'
125                     tauvari(igas) = 0.0
[253]126                  else
127                     print*,'Error in calc_rayleigh: Gas species not recognised!'
128                     call abort
129                  endif
130               endif
[135]131
[253]132               tauvar=tauvar+tauconsti(igas)*tauvari(igas)*gfrac(igas)
133
[135]134            enddo
135
[253]136            call blackl(dble(wl*1e-6),dble(tstellar),df)
137            df=df*bwidth/Nfine
138            tauwei=tauwei+df
139            tausum=tausum+tauvar*df
140         
141         enddo
142         TAURAY(N)=tausum/tauwei
[135]143         ! we multiply by scalep here (100) because plev, which is used in optcv,
144         ! is in units of mBar, so we need to convert.
[253]145
[135]146      end do
147
[848]148      IF (L_NSPECTV > 6) THEN
[861]149        icantbewrong = L_NSPECTV-6
150        print*,'At 1 atm and lambda = ',WAVEV(icantbewrong),' um'
151        print*,'tau_R = ',TAURAY(icantbewrong)*1013.25
152        print*,'sig_R = ',TAURAY(icantbewrong)*g*mugaz*1.67e-27*100, &
[848]153               'cm^2 molecule^-1'
154      ENDIF
[253]155
[135]156    end subroutine calc_rayleigh
Note: See TracBrowser for help on using the repository browser.