source: trunk/LMDZ.GENERIC/libf/phystd/interpolateN2N2.F90 @ 716

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