1 | subroutine interpolateN2H2(wn,temp,presN2,presH2,abcoef,firstcall) |
---|
2 | |
---|
3 | !================================================================== |
---|
4 | ! |
---|
5 | ! Purpose |
---|
6 | ! ------- |
---|
7 | ! Calculates the N2-H2 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 presN2 ! N2 partial pressure (Pascals) |
---|
23 | double precision presH2 ! H2 partial pressure (Pascals) |
---|
24 | |
---|
25 | ! output |
---|
26 | double precision abcoef ! absorption coefficient (m^-1) |
---|
27 | |
---|
28 | integer nS,nT |
---|
29 | parameter(nS=1914) |
---|
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 amagatN2 |
---|
38 | double precision amagatH2 |
---|
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 | if(temp.gt.400)then |
---|
58 | print*,'Your temperatures are too high for this N2-H2 CIA dataset.' |
---|
59 | print*,'Please run mixed N2-H2 atmospheres below T = 400 K.' |
---|
60 | stop |
---|
61 | endif |
---|
62 | |
---|
63 | amagatN2 = (273.15/temp)*(presN2/101325.0) |
---|
64 | amagatH2 = (273.15/temp)*(presH2/101325.0) |
---|
65 | |
---|
66 | if(firstcall)then ! called by sugas_corrk only |
---|
67 | print*,'----------------------------------------------------' |
---|
68 | print*,'Initialising N2-H2 continuum from HITRAN database...' |
---|
69 | |
---|
70 | ! 1.1 Open the ASCII files |
---|
71 | dt_file=TRIM(datadir)//'/continuum_data/N2-H2_2011.cia' |
---|
72 | |
---|
73 | open(33,file=dt_file,form='formatted',status='old',iostat=ios) |
---|
74 | if (ios.ne.0) then ! file not found |
---|
75 | write(*,*) 'Error from interpolateN2H2' |
---|
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 the continuum data continuum_data/N2-H2_2011.cia is there.' |
---|
81 | call abort |
---|
82 | else |
---|
83 | |
---|
84 | do iT=1,nT |
---|
85 | |
---|
86 | read(33,fmat1) bleh,blah,blah,nres,Ttemp |
---|
87 | if(nS.ne.nres)then |
---|
88 | print*,'Resolution given in file: ',trim(dt_file) |
---|
89 | print*,'is ',nres,', which does not match nS.' |
---|
90 | print*,'Please adjust nS value in interpolateN2H2.F90' |
---|
91 | stop |
---|
92 | endif |
---|
93 | temp_arr(iT)=Ttemp |
---|
94 | |
---|
95 | do k=1,nS |
---|
96 | read(33,*) wn_arr(k),abs_arr(k,it) |
---|
97 | end do |
---|
98 | |
---|
99 | end do |
---|
100 | |
---|
101 | endif |
---|
102 | close(33) |
---|
103 | |
---|
104 | print*,'interpolateN2H2: At wavenumber ',wn,' cm^-1' |
---|
105 | print*,' temperature ',temp,' K' |
---|
106 | print*,' N2 partial pressure ',presN2,' Pa' |
---|
107 | print*,' and H2 partial pressure ',presH2,' Pa' |
---|
108 | |
---|
109 | call bilinearN2H2(wn_arr,temp_arr,abs_arr,wn,temp,abcoef) |
---|
110 | |
---|
111 | print*,'the absorption is ',abcoef,' cm^5 molecule^-2' |
---|
112 | print*,'or ',abcoef*losch**2,' cm^-1 amagat^-2' |
---|
113 | |
---|
114 | abcoef=abcoef*losch**2*100.0*amagatN2*amagatH2 ! convert to m^-1 |
---|
115 | |
---|
116 | print*,'We have ',amagatN2,' amagats of N2' |
---|
117 | print*,'and ',amagatH2,' amagats of H2' |
---|
118 | print*,'So the absorption is ',abcoef,' m^-1' |
---|
119 | |
---|
120 | else |
---|
121 | |
---|
122 | call bilinearN2H2(wn_arr,temp_arr,abs_arr,wn,temp,abcoef) |
---|
123 | abcoef=abcoef*losch**2*100.0*amagatN2*amagatH2 ! convert to m^-1 |
---|
124 | |
---|
125 | ! unlike for Rayleigh scattering, we do not currently weight by the BB function |
---|
126 | ! however our bands are normally thin, so this is no big deal. |
---|
127 | |
---|
128 | endif |
---|
129 | |
---|
130 | return |
---|
131 | end subroutine interpolateN2H2 |
---|
132 | |
---|
133 | |
---|
134 | !------------------------------------------------------------------------- |
---|
135 | subroutine bilinearN2H2(x_arr,y_arr,f2d_arr,x_in,y_in,f) |
---|
136 | |
---|
137 | implicit none |
---|
138 | |
---|
139 | integer nX,nY,i,j,a,b |
---|
140 | parameter(nX=1914) |
---|
141 | parameter(nY=10) |
---|
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 bilinearN2H2:' |
---|
178 | write(*,*) 'Outside continuum 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 | call bilinear(f,f11,f21,f12,f22,x,x1,x2,y,y1,y2) |
---|
207 | |
---|
208 | return |
---|
209 | end subroutine bilinearN2H2 |
---|