source: trunk/LMDZ.GENERIC/libf/phystd/interpolateH2Ocont.F90 @ 527

Last change on this file since 527 was 526, checked in by jleconte, 13 years ago

13/02/2012 == JL + AS

  • All outputs are now in netCDF format. Even in 1D (No more G1D)
  • Clean up of the call to callcorrk when CLFvarying=true
  • Corrects a bug in writediagspecIR/VI. Output are now in W/m2/cm-1 as a function of the wavenumber in cm-1
  • Enable writediagspecIR/V to work in the CLFvarying=true case (output now done in Physiq after writediagfi)
  • Add a simple treatment for the supersaturation of CO2 (see forget et al 2012)
  • corrects a small bug when no clouds are present in aeropacity
  • Property svn:executable set to *
File size: 8.0 KB
Line 
1     subroutine interpolateH2Ocont(wn,temp,presS,presF,abcoef,firstcall)
2
3!==================================================================
4!     
5!     Purpose
6!     -------
7!     Calculates the H2O continuum opacity, using a lookup table from
8!     Clough (2005)
9!
10!     Authors
11!     -------
12!     R. Wordsworth (2011)
13!     
14!==================================================================
15
16      use datafile_mod, only: datadir
17      implicit none
18
19      ! input
20      double precision wn                 ! wavenumber             (cm^-1)
21      double precision temp               ! temperature            (Kelvin)
22      double precision presS              ! self-pressure          (Pascals)
23      double precision presF              ! foreign (air) pressure (Pascals)
24
25      ! output
26      double precision abcoef             ! absorption coefficient (m^-1)
27
28      integer nS,nT
29      parameter(nS=1001)
30      parameter(nT=11)
31
32      double precision kB
33      parameter(kB=1.3806488e-23)
34
35      double precision amagatS, amagatF, abcoefS, abcoefF, Nmolec
36      double precision wn_arr(nS)
37      double precision temp_arr(nT)
38      double precision abs_arrS(nS,nT)
39      double precision abs_arrF(nS,nT)
40      double precision data_tmp(nT)
41
42      integer k
43      logical firstcall
44
45      save wn_arr, temp_arr, abs_arrS, abs_arrF
46
47      character*100 dt_file
48      integer strlen,ios
49
50      amagatS=(273.15/temp)*(presS/101325.0)
51      amagatF=(273.15/temp)*(presF/101325.0)
52
53      if(firstcall)then
54
55!     1.1 Open the ASCII files
56
57         ! nu array
58         dt_file=TRIM(datadir)//'/continuum_data/H2O_CONT_NU.dat'
59         open(33,file=dt_file,form='formatted',status='old',iostat=ios)
60         if (ios.ne.0) then        ! file not found
61           write(*,*) 'Error from interpolateH2O_cont'
62           write(*,*) 'Data file ',trim(dt_file),' not found.'
63           write(*,*)'Check that your path to datagcm:',trim(datadir)
64           write(*,*)' is correct. You can change it in callphys.def with:'
65           write(*,*)' datadir = /absolute/path/to/datagcm'
66           write(*,*)' Also check that there is a continuum_data/H2O_CONT_NU.dat there.'
67           call abort
68         else
69            do k=1,nS
70               read(33,*) wn_arr(k)
71            enddo
72         endif
73         close(33)
74
75         ! self broadening
76         dt_file=TRIM(datadir)//'/continuum_data/H2O_CONT_SELF.dat'
77         open(34,file=dt_file,form='formatted',status='old',iostat=ios)
78         if (ios.ne.0) then        ! file not found
79           write(*,*) 'Error from interpolateH2O_cont'
80           write(*,*) 'Data file ',trim(dt_file),' not found.'
81           write(*,*)'Check that your path to datagcm:',trim(datadir)
82           write(*,*)' is correct. You can change it in callphys.def with:'
83           write(*,*)' datadir = /absolute/path/to/datagcm'
84           write(*,*)' Also check that there is a continuum_data/H2O_CONT_SELF.dat there.'
85           call abort
86         else
87            do k=1,nS
88               read(34,*) data_tmp
89               abs_arrS(k,1:nT)=data_tmp(1:nT)
90            end do
91         endif
92         close(34)
93
94         ! foreign (N2+O2+Ar) broadening
95         dt_file=TRIM(datadir)//'/continuum_data/H2O_CONT_FOREIGN.dat'
96         open(35,file=dt_file,form='formatted',status='old',iostat=ios)
97         if (ios.ne.0) then        ! file not found
98           write(*,*) 'Error from interpolateH2O_cont'
99           write(*,*) 'Data file ',trim(dt_file),' not found.'
100           write(*,*)'Check that your path to datagcm:',trim(datadir)
101           write(*,*)' is correct. You can change it in callphys.def with:'
102           write(*,*)' datadir = /absolute/path/to/datagcm'
103           write(*,*)' Also check that there is a continuum_data/H2O_CONT_FOREIGN.dat there.'
104           call abort
105         else
106            do k=1,nS
107               read(35,*) data_tmp
108               abs_arrF(k,1:nT)=data_tmp(1:nT)
109            end do
110         endif
111         close(35)
112
113         temp_arr(1)  = 200.
114         temp_arr(2)  = 250.
115         temp_arr(3)  = 300.
116         temp_arr(4)  = 350.
117         temp_arr(5)  = 400.
118         temp_arr(6)  = 450.
119         temp_arr(7)  = 500.
120         temp_arr(8)  = 550.
121         temp_arr(9)  = 600.
122         temp_arr(10) = 650.
123         temp_arr(11) = 700.
124
125         print*,'interpolateH2Ocont: At wavenumber ',wn,' cm^-1'
126         print*,'   temperature ',temp,' K'
127         print*,'   H2O pressure ',presS,' Pa'
128         print*,'   air pressure ',presF,' Pa'
129
130         call bilinearH2Ocont(wn_arr,temp_arr,abs_arrS,wn,temp,abcoefS)
131         print*,'the self absorption is ',abcoefS,' cm^2 molecule^-1'
132
133         call bilinearH2Ocont(wn_arr,temp_arr,abs_arrF,wn,temp,abcoefF)
134         print*,'the foreign absorption is ',abcoefF,' cm^2 molecule^-1'
135
136         print*,'We have ',amagatS,' amagats of H2O vapour'
137         print*,'and ',amagatF,' amagats of air'
138
139         abcoef = abcoefS*amagatS + abcoefF*amagatF ! Eq. (15) in Clough (1989)
140
141         Nmolec = (presS+presF)/(kB*temp)
142
143         print*,'Total number of molecules per m^3 is',Nmolec
144
145         abcoef = abcoef*Nmolec/(100.0**2) ! convert to m^-1
146         abcoef = abcoef*(presS/(presF+presS))     ! take H2O mixing ratio into account
147
148         print*,'So the total absorption is ',abcoef,' m^-1'
149         print*,'And optical depth / km : ',1000.0*abcoef
150
151      else
152
153         call bilinearH2Ocont(wn_arr,temp_arr,abs_arrS,wn,temp,abcoefS)
154         call bilinearH2Ocont(wn_arr,temp_arr,abs_arrF,wn,temp,abcoefF)
155
156         abcoef = abcoefS*amagatS + abcoefF*amagatF ! Eq. (15) in Clough (1989)
157
158         !abcoef = (presS/100000.0)*3.0e-24
159         !print*,'Matsui TEST'
160
161         Nmolec = (presS+presF)/(kB*temp)  ! sure this is correct??
162         abcoef = abcoef*Nmolec/(100.0**2) ! convert to m^-1
163         abcoef = abcoef*(presS/(presF+presS))     ! take H2O mixing ratio into account
164
165
166         ! unlike for Rayleigh scattering, we do not currently weight by the BB function
167         ! however our bands are normally thin, so this is no big deal.
168
169
170      endif
171
172      return
173    end subroutine interpolateH2Ocont
174
175
176
177
178!-------------------------------------------------------------------------
179      subroutine bilinearH2Ocont(x_arr,y_arr,f2d_arr,x_in,y_in,f)
180!     Necessary for interpolation of continuum data
181
182      implicit none
183
184      integer nX,nY,i,j,a,b
185      parameter(nX=1001)
186      parameter(nY=11)
187
188      real*8 x_in,y_in,x,y,x1,x2,y1,y2
189      real*8 f,f11,f12,f21,f22,fA,fB
190      real*8 x_arr(nX)
191      real*8 y_arr(nY)
192      real*8 f2d_arr(nX,nY)
193     
194      integer strlen
195      character*100 label
196      label='subroutine bilinear'
197
198      x=x_in
199      y=y_in
200
201!     1st check we're within the wavenumber range
202      if ((x.lt.x_arr(2)).or.(x.gt.x_arr(nX-2))) then
203         f=0.0D+0
204         return
205      else
206         
207!     in the x (wavenumber) direction 1st
208         i=1
209 10      if (i.lt.(nX+1)) then
210            if (x_arr(i).gt.x) then
211               x1=x_arr(i-1)
212               x2=x_arr(i)
213               a=i-1
214               i=9999
215            endif
216            i=i+1
217            goto 10
218         endif
219      endif
220     
221      if ((y.lt.y_arr(1)).or.(y.gt.y_arr(nY))) then
222         write(*,*) 'Warning from bilinearH2Ocont:'
223         write(*,*) 'Outside continuum temperature range!'
224         if(y.lt.y_arr(1))then
225            y=y_arr(1)+0.01
226         endif
227         if(y.gt.y_arr(nY))then
228            y=y_arr(nY)-0.01
229         endif
230      else
231
232!     in the y (temperature) direction 2nd
233         j=1
234 20      if (j.lt.(nY+1)) then
235            if (y_arr(j).gt.y) then
236               y1=y_arr(j-1)
237               y2=y_arr(j)
238               b=j-1
239               j=9999
240            endif
241            j=j+1
242            goto 20
243         endif
244      endif
245     
246      f11=f2d_arr(a,b)
247      f21=f2d_arr(a+1,b)
248      f12=f2d_arr(a,b+1)
249      f22=f2d_arr(a+1,b+1)
250     
251!     1st in x-direction
252      fA=f11*(x2-x)/(x2-x1)+f21*(x-x1)/(x2-x1)
253      fB=f12*(x2-x)/(x2-x1)+f22*(x-x1)/(x2-x1)
254     
255!     then in y-direction
256      f=fA*(y2-y)/(y2-y1)+fB*(y-y1)/(y2-y1)
257
258      return
259    end subroutine bilinearH2Ocont
Note: See TracBrowser for help on using the repository browser.