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

Last change on this file since 323 was 305, checked in by rwordsworth, 13 years ago

Several new files added as part of the climate evolution model
(main program kcm.F90). Some general cleanup in physiq.F90 and
callcorrk.F90. Bugs in dust radiative transfer and H2 Rayleigh
scattering corrected.

  • Property svn:executable set to *
File size: 7.2 KB
RevLine 
[305]1     subroutine interpolateH2Ocont(wn,temp,presS,presF,abcoef,firstcall)
2      implicit none
3
4!==================================================================
5!     
6!     Purpose
7!     -------
8!     Calculates the H2O continuum opacity, using a lookup table from
9!     Clough (2005)
10!
11!     Authors
12!     -------
13!     R. Wordsworth (2011)
14!     
15!==================================================================
16
17#include "datafile.h"
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=datafile(1:LEN_TRIM(datafile))//'/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.for'
62            write(*,*) 'Data file could not be found:'
63            write(*,*) dt_file
64            call abort
65         else
66            do k=1,nS
67               read(33,*) wn_arr(k)
68            enddo
69         endif
70         close(33)
71
72         ! self broadening
73         dt_file=datafile(1:LEN_TRIM(datafile))//'/continuum_data/H2O_CONT_SELF.dat'
74         open(34,file=dt_file,form='formatted',status='old',iostat=ios)
75         if (ios.ne.0) then        ! file not found
76            write(*,*) 'Error from interpolateH2O_cont.for'
77            write(*,*) 'Data file could not be found:'
78            write(*,*) dt_file
79            call abort
80         else
81            do k=1,nS
82               read(34,*) data_tmp
83               abs_arrS(k,1:nT)=data_tmp(1:nT)
84            end do
85         endif
86         close(34)
87
88         ! foreign (N2+O2+Ar) broadening
89         dt_file=datafile(1:LEN_TRIM(datafile))//'/continuum_data/H2O_CONT_FOREIGN.dat'
90         open(35,file=dt_file,form='formatted',status='old',iostat=ios)
91         if (ios.ne.0) then        ! file not found
92            write(*,*) 'Error from interpolateH2O_cont.for'
93            write(*,*) 'Data file could not be found:'
94            write(*,*) dt_file
95            call abort
96         else
97            do k=1,nS
98               read(35,*) data_tmp
99               abs_arrF(k,1:nT)=data_tmp(1:nT)
100            end do
101         endif
102         close(35)
103
104         temp_arr(1)  = 200.
105         temp_arr(2)  = 250.
106         temp_arr(3)  = 300.
107         temp_arr(4)  = 350.
108         temp_arr(5)  = 400.
109         temp_arr(6)  = 450.
110         temp_arr(7)  = 500.
111         temp_arr(8)  = 550.
112         temp_arr(9)  = 600.
113         temp_arr(10) = 650.
114         temp_arr(11) = 700.
115
116         print*,'At wavenumber ',wn,' cm^-1'
117         print*,'   temperature ',temp,' K'
118         print*,'   H2O pressure ',presS,' Pa'
119         print*,'   air pressure ',presF,' Pa'
120
121         call bilinearH2Ocont(wn_arr,temp_arr,abs_arrS,wn,temp,abcoefS)
122         print*,'the self absorption is ',abcoefS,' cm^2 molecule^-1'
123
124         call bilinearH2Ocont(wn_arr,temp_arr,abs_arrF,wn,temp,abcoefF)
125         print*,'the foreign absorption is ',abcoefF,' cm^2 molecule^-1'
126
127         print*,'We have ',amagatS,' amagats of H2O vapour'
128         print*,'and ',amagatF,' amagats of air'
129
130         abcoef = abcoefS*amagatS + abcoefF*amagatF ! Eq. (15) in Clough (1989)
131
132         Nmolec = (presS+presF)/(kB*temp)
133
134         print*,'Total number of molecules per m^3 is',Nmolec
135
136         abcoef = abcoef*Nmolec/(100.0**2) ! convert to m^-1
137         abcoef = abcoef*(presS/(presF+presS))     ! take H2O mixing ratio into account
138
139         print*,'So the total absorption is ',abcoef,' m^-1'
140         print*,'And optical depth / km : ',1000.0*abcoef
141
142      else
143 
144         call bilinearH2Ocont(wn_arr,temp_arr,abs_arrS,wn,temp,abcoefS)
145         call bilinearH2Ocont(wn_arr,temp_arr,abs_arrF,wn,temp,abcoefF)
146
147         abcoef = abcoefS*amagatS + abcoefF*amagatF ! Eq. (15) in Clough (1989)
148
149         !abcoef = (presS/100000.0)*3.0e-24
150         !print*,'Matsui TEST'
151
152         Nmolec = (presS+presF)/(kB*temp)  ! sure this is correct??
153         abcoef = abcoef*Nmolec/(100.0**2) ! convert to m^-1
154         abcoef = abcoef*(presS/(presF+presS))     ! take H2O mixing ratio into account
155
156         ! unlike for Rayleigh scattering, we do not currently weight by the BB function
157         ! however our bands are normally thin, so this is no big deal.
158
159
160      endif
161
162      return
163    end subroutine interpolateH2Ocont
164
165
166
167
168!-------------------------------------------------------------------------
169      subroutine bilinearH2Ocont(x_arr,y_arr,f2d_arr,x_in,y_in,f)
170!     Necessary for interpolation of continuum data
171
172      implicit none
173
174      integer nX,nY,i,j,a,b
175      parameter(nX=1001)
176      parameter(nY=11)
177
178      real*8 x_in,y_in,x,y,x1,x2,y1,y2
179      real*8 f,f11,f12,f21,f22,fA,fB
180      real*8 x_arr(nX)
181      real*8 y_arr(nY)
182      real*8 f2d_arr(nX,nY)
183     
184      integer strlen
185      character*100 label
186      label='subroutine bilinear'
187
188      x=x_in
189      y=y_in
190
191!     1st check we're within the wavenumber range
192      if ((x.lt.x_arr(2)).or.(x.gt.x_arr(nX-2))) then
193         f=0.0D+0
194         return
195      else
196         
197!     in the x (wavenumber) direction 1st
198         i=1
199 10      if (i.lt.(nX+1)) then
200            if (x_arr(i).gt.x) then
201               x1=x_arr(i-1)
202               x2=x_arr(i)
203               a=i-1
204               i=9999
205            endif
206            i=i+1
207            goto 10
208         endif
209      endif
210     
211      if ((y.lt.y_arr(1)).or.(y.gt.y_arr(nY))) then
212         write(*,*) 'Warning from bilinear.for:'
213         write(*,*) 'Outside continuum temperature range!'
214         if(y.lt.y_arr(1))then
215            y=y_arr(1)+0.01
216         endif
217         if(y.gt.y_arr(nY))then
218            y=y_arr(nY)-0.01
219         endif
220      else
221
222!     in the y (temperature) direction 2nd
223         j=1
224 20      if (j.lt.(nY+1)) then
225            if (y_arr(j).gt.y) then
226               y1=y_arr(j-1)
227               y2=y_arr(j)
228               b=j-1
229               j=9999
230            endif
231            j=j+1
232            goto 20
233         endif
234      endif
235     
236      f11=f2d_arr(a,b)
237      f21=f2d_arr(a+1,b)
238      f12=f2d_arr(a,b+1)
239      f22=f2d_arr(a+1,b+1)
240     
241!     1st in x-direction
242      fA=f11*(x2-x)/(x2-x1)+f21*(x-x1)/(x2-x1)
243      fB=f12*(x2-x)/(x2-x1)+f22*(x-x1)/(x2-x1)
244     
245!     then in y-direction
246      f=fA*(y2-y)/(y2-y1)+fB*(y-y1)/(y2-y1)
247
248      return
249    end subroutine bilinearH2Ocont
Note: See TracBrowser for help on using the repository browser.