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

Last change on this file since 537 was 374, checked in by emillour, 13 years ago

Generic GCM: Upgrade: The location of the 'datagcm' directory can now be given in the callphys.def file ( datadir = /absolute/path/to/datagcm ). Changed "datafile.h" into a F90 module "datafile_mod.F90" and spread this change to all routines that used to use "datafile.h".
EM

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