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

Last change on this file since 537 was 471, checked in by aslmd, 13 years ago

LMDZ.GENERIC

13/12/2011 == AS

  • Same spirit as previous commit, but for ngasmx which is now read in gases.def -- before arrays w/ dim ngasmx are allocated dynamically
  • Allocation is done in su_gases.F90 which is called in inifis
  • Outside su_gases.F90, very few modifications to the code : the new module "gases_h.F90" simply replaces the old common "gases.h" !
  • Compiles fine. Tested with debugging options through pgdbg. Runs fine. Exact same results in Early Mars test case.
File size: 4.8 KB
Line 
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 CO2, N2 and
10!     H2 are possible).
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
27      use radcommon_h, only: WAVEV, BWNV, DWNV, tstellar, tauray, scalep
28      use gases_h
29
30      implicit none
31
32#include "comcstfi.h"
33#include "callkeys.h"
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      ! tau0/p0=tau/p (Hansen 1974)
48      ! Then calculate tau0 = (8*pi^3*p_1atm)/(3*N0^2) * 4*delta^2/(g*mugaz*lambda^4)
49      ! where delta=n-1 and N0 is an amagat
50
51      typeknown=.false.
52      do igas=1,ngasmx
53         if(igas.eq.vgas)then
54            print*,'Ignoring ',gnom(igas),' in Rayleigh scattering '// &
55            'as it is variable.'
56         elseif(gfrac(igas).lt.5.e-2)then
57            print*,'Ignoring ',gnom(igas),' in Rayleigh scattering '// &
58            'as its mixing ratio is less than 0.05.'
59            ! ignore variable gas in Rayleigh calculation
60            ! ignore gases of mixing ratio < 0.05 in Rayleigh calculation
61            tauconsti(igas) = 0.0
62         else
63            if(gnom(igas).eq.'CO2')then
64               tauconsti(igas) = (8.7/g)*1.527*scalep/P0
65            elseif(gnom(igas).eq.'N2_')then
66               tauconsti(igas) = (9.81/g)*8.569E-3*scalep/(P0/93.0)
67            elseif(gnom(igas).eq.'H2O')then
68               tauconsti(igas) = (10.0/g)*9.22E-3*scalep/101325.0
69            elseif(gnom(igas).eq.'H2_')then
70               tauconsti(igas) = (10.0/g)*0.0241*scalep/101325.0
71               !tauconsti(igas) = (10.0/g)*0.0487*scalep/(101325.0)
72               ! uses n=1.000132 from Optics, Fourth Edition.
73               ! was out by a factor 2!
74            elseif(gnom(igas).eq.'He_')then
75               print*,'Helium not ready yet!'
76            else
77               print*,'Error in calc_rayleigh: Gas species not recognised!'
78               call abort
79            endif
80
81            if(gfrac(igas).eq.1.0)then
82               print*,'Rayleigh scattering is for a pure ',gnom(igas),' atmosphere.'
83               typeknown=.true.
84            endif
85
86         endif
87      enddo
88
89      if(.not.typeknown)then
90         print*,'Rayleigh scattering is for a mixed gas atmosphere.'
91         typeknown=.true.
92      endif
93
94 
95      do N=1,L_NSPECTV
96         
97         tausum = 0.0
98         tauwei = 0.0
99         bstart = 10000.0/BWNV(N+1)
100         bwidth = (10000.0/BWNV(N)) - (10000.0/BWNV(N+1))
101         do ifine=1,Nfine
102            wl=bstart+dble(ifine)*bwidth/Nfine
103
104            tauvar=0.0
105            do igas=1,ngasmx
106               if((igas.eq.vgas).or.(gfrac(igas).lt.5.e-2))then
107                  ! ignore variable gas in Rayleigh calculation
108                  ! ignore gases of mixing ratio < 0.05 in Rayleigh calculation
109                  tauvari(igas)   = 0.0
110               else
111                  if(gnom(igas).eq.'CO2')then
112                     tauvari(igas) = (1.0+0.013/wl**2)/wl**4
113                  elseif(gnom(igas).eq.'N2_')then
114                     tauvari(igas) = (1.0+0.0113/wl**2+0.00013/wl**4)/wl**4
115                  elseif(gnom(igas).eq.'H2O')then
116                     tauvari(igas) = 1.0/wl**4 ! to be improved...
117                  elseif(gnom(igas).eq.'H2_')then
118                     tauvari(igas) = 1.0/wl**4
119                  else
120                     print*,'Error in calc_rayleigh: Gas species not recognised!'
121                     call abort
122                  endif
123               endif
124
125               tauvar=tauvar+tauconsti(igas)*tauvari(igas)*gfrac(igas)
126
127            enddo
128
129            call blackl(dble(wl*1e-6),dble(tstellar),df)
130            df=df*bwidth/Nfine
131            tauwei=tauwei+df
132            tausum=tausum+tauvar*df
133         
134         enddo
135         TAURAY(N)=tausum/tauwei
136         ! we multiply by scalep here (100) because plev, which is used in optcv,
137         ! is in units of mBar, so we need to convert.
138
139      end do
140
141      print*,'At 1 atm and lambda = ',WAVEV(L_NSPECTV-6),' um'
142      print*,'tau_R = ',TAURAY(L_NSPECTV-6)*1013.25
143      print*,'sig_R = ',TAURAY(L_NSPECTV-6)*g*mugaz*1.67e-27*100, &
144             'cm^2 molecule^-1'
145
146    end subroutine calc_rayleigh
Note: See TracBrowser for help on using the repository browser.