source: trunk/LMDZ.GENERIC/libf/phystd/interpolateH2He.F90 @ 848

Last change on this file since 848 was 716, checked in by rwordsworth, 12 years ago

Mainly updates to radiative transfer and gas management scheme.
Most CIA data now read from standard HITRAN datafiles. For the H2O
continuum, two options have been added: the standard CKD continuum,
and the empirical formula in PPC (Pierrehumbert 2010). Use the toggle
'H2Ocont_simple' in callphys.def to choose.

Note to Martians: I've changed the default values of 'sedimentation' and
'co2cond' in inifis to false. Both these are defined in the standard deftank
callphys.def file, so there shouldn't be any compatibility problems.

File size: 6.1 KB
Line 
1     subroutine interpolateH2He(wn,temp,presH2,presHe,abcoef,firstcall)
2
3!==================================================================
4!     
5!     Purpose
6!     -------
7!     Calculates the H2-He CIA opacity, using a lookup table from
8!     HITRAN (2011)
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 presH2             ! H2 partial pressure    (Pascals)
23      double precision presHe             ! He partial pressure    (Pascals)
24
25      ! output
26      double precision abcoef             ! absorption coefficient (m^-1)
27
28      integer nS,nT
29      parameter(nS=2428)
30      parameter(nT=10)
31
32      double precision, parameter :: losch = 2.6867774e19
33      ! Loschmit's number (molecule cm^-3 at STP)
34      ! converts cm^5 molecule^-2 --> cm^-1 amagat^-2
35      ! see Richard et al. 2011, JQSRT for details
36
37      double precision amagatH2
38      double precision amagatHe
39      double precision wn_arr(nS)
40      double precision temp_arr(nT)
41      double precision abs_arr(nS,nT)
42
43      integer k,iT
44      logical firstcall
45
46      save wn_arr, temp_arr, abs_arr
47
48      character*100 dt_file
49      integer strlen,ios
50
51      character(LEN=*), parameter :: fmat1 = "(A20,F10.3,F10.3,I7,F7.1,E10.3,F5.3)"
52
53      character*20 bleh
54      double precision blah, Ttemp
55      integer nres
56
57
58      if(temp.gt.400)then
59         print*,'Your temperatures are too high for this H2-He CIA dataset.'
60         print*,'Please run mixed H2-He atmospheres below T = 400 K.'     
61         stop
62      endif
63
64      amagatH2 = (273.15/temp)*(presH2/101325.0)
65      amagatHe = (273.15/temp)*(presHe/101325.0)
66
67      if(firstcall)then ! called by sugas_corrk only
68         print*,'----------------------------------------------------'
69         print*,'Initialising H2-He continuum from HITRAN database...'
70
71!     1.1 Open the ASCII files
72         dt_file=TRIM(datadir)//'/continuum_data/H2-He_norm_2011.cia'
73
74         open(33,file=dt_file,form='formatted',status='old',iostat=ios)
75         if (ios.ne.0) then        ! file not found
76           write(*,*) 'Error from interpolateH2He'
77           write(*,*) 'Data file ',trim(dt_file),' not found.'
78           write(*,*) 'Check that your path to datagcm:',trim(datadir)
79           write(*,*) 'is correct. You can change it in callphys.def with:'
80           write(*,*) 'datadir = /absolute/path/to/datagcm'
81           write(*,*) 'Also check that the continuum data continuum_data/H2-He_norm_2011.cia is there.'
82           call abort
83         else
84
85            do iT=1,nT
86
87               read(33,fmat1) bleh,blah,blah,nres,Ttemp
88               if(nS.ne.nres)then
89                  print*,'Resolution given in file: ',trim(dt_file)
90                  print*,'is ',nres,', which does not match nS.'
91                  print*,'Please adjust nS value in interpolateH2He.F90'
92                  stop
93               endif
94               temp_arr(iT)=Ttemp
95
96               do k=1,nS
97                  read(33,*) wn_arr(k),abs_arr(k,it)
98               end do
99
100            end do
101     
102         endif
103         close(33)
104
105         print*,'interpolateH2He: At wavenumber ',wn,' cm^-1'
106         print*,'   temperature                 ',temp,' K'
107         print*,'   H2 partial pressure         ',presH2,' Pa'
108         print*,'   and He partial pressure     ',presHe,' Pa'
109
110         call bilinearH2He(wn_arr,temp_arr,abs_arr,wn,temp,abcoef)
111
112         print*,'the absorption is ',abcoef,' cm^5 molecule^-2'
113         print*,'or ',abcoef*losch**2,' cm^-1 amagat^-2'
114
115         abcoef=abcoef*losch**2*100.0*amagatH2*amagatHe ! convert to m^-1
116
117         print*,'We have ',amagatH2,' amagats of H2'
118         print*,'and     ',amagatHe,' amagats of He'
119         print*,'So the absorption is ',abcoef,' m^-1'
120
121      else
122 
123         call bilinearH2He(wn_arr,temp_arr,abs_arr,wn,temp,abcoef)
124         abcoef=abcoef*losch**2*100.0*amagatH2*amagatHe ! convert to m^-1
125
126         ! unlike for Rayleigh scattering, we do not currently weight by the BB function
127         ! however our bands are normally thin, so this is no big deal.
128
129      endif
130
131      return
132    end subroutine interpolateH2He
133
134
135!-------------------------------------------------------------------------
136      subroutine bilinearH2He(x_arr,y_arr,f2d_arr,x_in,y_in,f)
137!     Necessary for interpolation of continuum data
138
139      implicit none
140
141      integer nX,nY,i,j,a,b
142      parameter(nX=2428)
143      parameter(nY=10)
144
145      real*8 x_in,y_in,x,y,x1,x2,y1,y2
146      real*8 f,f11,f12,f21,f22,fA,fB
147      real*8 x_arr(nX)
148      real*8 y_arr(nY)
149      real*8 f2d_arr(nX,nY)
150     
151      integer strlen
152      character*100 label
153      label='subroutine bilinear'
154
155      x=x_in
156      y=y_in
157
158!     1st check we're within the wavenumber range
159      if ((x.lt.x_arr(2)).or.(x.gt.x_arr(nX-2))) then
160         f=0.0D+0
161         return
162      else
163         
164!     in the x (wavenumber) direction 1st
165         i=1
166 10      if (i.lt.(nX+1)) then
167            if (x_arr(i).gt.x) then
168               x1=x_arr(i-1)
169               x2=x_arr(i)
170               a=i-1
171               i=9999
172            endif
173            i=i+1
174            goto 10
175         endif
176      endif
177     
178      if ((y.lt.y_arr(1)).or.(y.gt.y_arr(nY))) then
179         write(*,*) 'Warning from bilinearH2He:'
180         write(*,*) 'Outside continuum temperature range!'
181         if(y.lt.y_arr(1))then
182            y=y_arr(1)+0.01
183         endif
184         if(y.gt.y_arr(nY))then
185            y=y_arr(nY)-0.01
186         endif
187      else
188
189!     in the y (temperature) direction 2nd
190         j=1
191 20      if (j.lt.(nY+1)) then
192            if (y_arr(j).gt.y) then
193               y1=y_arr(j-1)
194               y2=y_arr(j)
195               b=j-1
196               j=9999
197            endif
198            j=j+1
199            goto 20
200         endif
201      endif
202     
203      f11=f2d_arr(a,b)
204      f21=f2d_arr(a+1,b)
205      f12=f2d_arr(a,b+1)
206      f22=f2d_arr(a+1,b+1)
207
208      call bilinear(f,f11,f21,f12,f22,x,x1,x2,y,y1,y2)
209
210     
211      return
212    end subroutine bilinearH2He
Note: See TracBrowser for help on using the repository browser.