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

Last change on this file since 277 was 253, checked in by emillour, 14 years ago

Generic GCM

  • Massive update to version 0.7

EM+RW

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