source: trunk/LMDZ.GENERIC/libf/phystd/interpolateH2H2.F90 @ 344

Last change on this file since 344 was 305, checked in by rwordsworth, 14 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: 5.6 KB
Line 
1     subroutine interpolateH2H2(wn,temp,pres,abcoef,firstcall)
2      implicit none
3
4!==================================================================
5!     
6!     Purpose
7!     -------
8!     Calculates the H2-H2 CIA opacity, using a lookup table from
9!     Borysow (2002)
10!
11!     Authors
12!     -------
13!     R. Wordsworth (2009)
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 pres               ! pressure               (Pascals)
23
24      ! output
25      double precision abcoef             ! absorption coefficient (m^-1)
26
27      integer nS,nT
28      parameter(nS=1649)
29      parameter(nT=14)
30
31      double precision amagat
32      double precision wn_arr(nS)
33      double precision temp_arr(nT)
34      double precision abs_arr(nS,nT)
35      double precision data_tmp(nT/2+1)
36
37      integer k
38      logical firstcall
39
40      save wn_arr, temp_arr, abs_arr
41
42      character*100 dt_file
43      integer strlen,ios
44
45      amagat=(273.15/temp)*(pres/101325.0)
46
47      if(firstcall)then
48
49!     1.1 Open the ASCII files
50
51         ! cold
52         dt_file=datafile(1:LEN_TRIM(datafile))//'/continuum_data/H2H2_CIA_LT.dat'
53         open(33,file=dt_file,form='formatted',status='old',iostat=ios)
54         if (ios.ne.0) then        ! file not found
55            write(*,*) 'Error from interpolateH2H2.for'
56            write(*,*) 'Data file could not be found:'
57            write(*,*) dt_file
58            call abort
59         else
60            do k=1,nS
61               read(33,*) data_tmp
62               wn_arr(k)=data_tmp(1)
63               abs_arr(k,1:7)=data_tmp(2:8)
64            end do
65         endif
66         close(33)
67
68         ! hot
69         dt_file=datafile(1:LEN_TRIM(datafile))//'/continuum_data/H2H2_CIA_HT.dat'
70         open(34,file=dt_file,form='formatted',status='old',iostat=ios)
71         if (ios.ne.0) then        ! file not found
72            write(*,*) 'Error from interpolateH2H2.for'
73            write(*,*) 'Data file could not be found:'
74            write(*,*) dt_file
75            call abort
76         else
77            do k=1,nS
78               read(34,*) data_tmp
79               wn_arr(k)=data_tmp(1)
80               ! wn_arr is identical
81               abs_arr(k,8:14)=data_tmp(2:8)
82            end do
83         endif
84         close(34)
85
86         temp_arr(1)  = 60.
87         temp_arr(2)  = 100.
88         temp_arr(3)  = 150.
89         temp_arr(4)  = 200.
90         temp_arr(5)  = 250.
91         temp_arr(6)  = 300.
92         temp_arr(7)  = 350.
93         temp_arr(8)  = 400.
94         temp_arr(9)  = 500.
95         temp_arr(10) = 600.
96         temp_arr(11) = 700.
97         temp_arr(12) = 800.
98         temp_arr(13) = 900.
99         temp_arr(14) = 1000.
100
101
102         print*,'At wavenumber ',wn,' cm^-1'
103         print*,'   temperature ',temp,' K'
104         print*,'   pressure ',pres,' Pa'
105
106         call bilinearH2H2(wn_arr,temp_arr,abs_arr,wn,temp,abcoef)
107
108         print*,'the absorption is ',abcoef,' cm^-1 amg^-2'
109
110         abcoef=abcoef*100.0*amagat**2 ! convert to m^-1
111
112         print*,'We have ',amagat,' amagats'
113         print*,'So the absorption is ',abcoef,' m^-1'
114
115      else
116 
117         call bilinearH2H2(wn_arr,temp_arr,abs_arr,wn,temp,abcoef)
118         abcoef=abcoef*100.0*amagat**2 ! convert to m^-1
119         !print*,'The absorption is ',abcoef,' m^-1'
120
121         ! unlike for Rayleigh scattering, we do not currently weight by the BB function
122         ! however our bands are normally thin, so this is no big deal.
123
124      endif
125
126      return
127    end subroutine interpolateH2H2
128
129
130!-------------------------------------------------------------------------
131      subroutine bilinearH2H2(x_arr,y_arr,f2d_arr,x_in,y_in,f)
132!     Necessary for interpolation of continuum data
133
134      implicit none
135
136      integer nX,nY,i,j,a,b
137      parameter(nX=1649)
138      parameter(nY=14)
139
140      real*8 x_in,y_in,x,y,x1,x2,y1,y2
141      real*8 f,f11,f12,f21,f22,fA,fB
142      real*8 x_arr(nX)
143      real*8 y_arr(nY)
144      real*8 f2d_arr(nX,nY)
145     
146      integer strlen
147      character*100 label
148      label='subroutine bilinear'
149
150      x=x_in
151      y=y_in
152
153!     1st check we're within the wavenumber range
154      if ((x.lt.x_arr(2)).or.(x.gt.x_arr(nX-2))) then
155         f=0.0D+0
156         return
157      else
158         
159!     in the x (wavenumber) direction 1st
160         i=1
161 10      if (i.lt.(nX+1)) then
162            if (x_arr(i).gt.x) then
163               x1=x_arr(i-1)
164               x2=x_arr(i)
165               a=i-1
166               i=9999
167            endif
168            i=i+1
169            goto 10
170         endif
171      endif
172     
173      if ((y.lt.y_arr(1)).or.(y.gt.y_arr(nY))) then
174         write(*,*) 'Warning from bilinear.for:'
175         write(*,*) 'Outside continuum temperature range!'
176         if(y.lt.y_arr(1))then
177            y=y_arr(1)+0.01
178         endif
179         if(y.gt.y_arr(nY))then
180            y=y_arr(nY)-0.01
181         endif
182      else
183
184!     in the y (temperature) direction 2nd
185         j=1
186 20      if (j.lt.(nY+1)) then
187            if (y_arr(j).gt.y) then
188               y1=y_arr(j-1)
189               y2=y_arr(j)
190               b=j-1
191               j=9999
192            endif
193            j=j+1
194            goto 20
195         endif
196      endif
197     
198      f11=f2d_arr(a,b)
199      f21=f2d_arr(a+1,b)
200      f12=f2d_arr(a,b+1)
201      f22=f2d_arr(a+1,b+1)
202     
203!     1st in x-direction
204      fA=f11*(x2-x)/(x2-x1)+f21*(x-x1)/(x2-x1)
205      fB=f12*(x2-x)/(x2-x1)+f22*(x-x1)/(x2-x1)
206     
207!     then in y-direction
208      f=fA*(y2-y)/(y2-y1)+fB*(y-y1)/(y2-y1)
209     
210      return
211    end subroutine bilinearH2H2
Note: See TracBrowser for help on using the repository browser.