1 | subroutine sugas_corrk |
---|
2 | |
---|
3 | !================================================================== |
---|
4 | ! |
---|
5 | ! Purpose |
---|
6 | ! ------- |
---|
7 | ! Set up gaseous absorption parameters used by the radiation code. |
---|
8 | ! This subroutine is a replacement for the old 'setrad', which contained |
---|
9 | ! both absorption and scattering data. |
---|
10 | ! |
---|
11 | ! Authors |
---|
12 | ! ------- |
---|
13 | ! Adapted and generalised from the NASA Ames code by Robin Wordsworth (2009) |
---|
14 | ! |
---|
15 | ! Summary |
---|
16 | ! ------- |
---|
17 | ! |
---|
18 | !================================================================== |
---|
19 | |
---|
20 | use radinc_h |
---|
21 | use radcommon_h, only : pgasref,pfgasref,pgasmin,pgasmax |
---|
22 | use radcommon_h, only : tgasref,tgasmin,tgasmax |
---|
23 | use radcommon_h, only : gasv,gasi,FZEROI,FZEROV,gweight |
---|
24 | use radcommon_h, only : wrefvar,gastype |
---|
25 | use datafile_mod, only: datadir |
---|
26 | implicit none |
---|
27 | |
---|
28 | #include "callkeys.h" |
---|
29 | #include "gases.h" |
---|
30 | |
---|
31 | !================================================================== |
---|
32 | |
---|
33 | logical file_ok |
---|
34 | |
---|
35 | integer n, nt, np, nh, ng, nw, m, i |
---|
36 | integer L_NGAUSScheck, L_NPREFcheck, L_NTREFcheck, L_REFVARcheck |
---|
37 | |
---|
38 | character(len=100) :: file_id |
---|
39 | character(len=300) :: file_path |
---|
40 | |
---|
41 | real*8 gasi8(L_NTREF,L_NPREF,L_REFVAR,L_NSPECTI,L_NGAUSS) |
---|
42 | real*8 gasv8(L_NTREF,L_NPREF,L_REFVAR,L_NSPECTV,L_NGAUSS) |
---|
43 | |
---|
44 | real*8 x, xi(4), yi(4), ans, p |
---|
45 | |
---|
46 | integer ngas, igas |
---|
47 | |
---|
48 | double precision testcont ! for continuum absorption initialisation |
---|
49 | |
---|
50 | !======================================================================= |
---|
51 | ! Load variable species data, exit if we have wrong database |
---|
52 | file_id='/corrk_data/' // TRIM(corrkdir) // '/Q.dat' |
---|
53 | file_path=TRIM(datadir)//TRIM(file_id) |
---|
54 | |
---|
55 | |
---|
56 | ! check that the file exists |
---|
57 | inquire(FILE=file_path,EXIST=file_ok) |
---|
58 | if(.not.file_ok) then |
---|
59 | write(*,*)'The file ',TRIM(file_path) |
---|
60 | write(*,*)'was not found by sugas_corrk.F90, exiting.' |
---|
61 | write(*,*)'Check that your path to datagcm:',trim(datadir) |
---|
62 | write(*,*)' is correct. You can change it in callphys.def with:' |
---|
63 | write(*,*)' datadir = /absolute/path/to/datagcm' |
---|
64 | write(*,*)'Also check that the corrkdir you chose in callphys.def exists.' |
---|
65 | call abort |
---|
66 | endif |
---|
67 | |
---|
68 | ! check that database matches varactive toggle |
---|
69 | open(111,file=TRIM(file_path),form='formatted') |
---|
70 | read(111,*) ngas |
---|
71 | |
---|
72 | if(ngas.ne.ngasmx)then |
---|
73 | print*,'Number of gases in radiative transfer data (',ngas,') does not', & |
---|
74 | 'match that in gases.def (',ngasmx,'), exiting.' |
---|
75 | call abort |
---|
76 | endif |
---|
77 | |
---|
78 | if(ngas.eq.1 .and. (varactive.or.varfixed))then |
---|
79 | print*,'You have varactive/fixed=.true. but the database [', & |
---|
80 | corrkdir(1:LEN_TRIM(corrkdir)), & |
---|
81 | '] has no variable species, exiting.' |
---|
82 | call abort |
---|
83 | elseif(ngas.eq.2 .and. (.not.varactive) .and. (.not.varfixed))then |
---|
84 | print*,'You have varactive and varfixed=.false. and the database [', & |
---|
85 | corrkdir(1:LEN_TRIM(corrkdir)), & |
---|
86 | '] has a variable species.' |
---|
87 | call abort |
---|
88 | elseif(ngas.gt.4 .or. ngas.lt.1)then |
---|
89 | print*,ngas,' species in database [', & |
---|
90 | corrkdir(1:LEN_TRIM(corrkdir)), & |
---|
91 | '], radiative code cannot handle this.' |
---|
92 | call abort |
---|
93 | endif |
---|
94 | |
---|
95 | if(ngas.gt.3)then |
---|
96 | print*,'ngas>3, EXPERIMENTAL!' |
---|
97 | endif |
---|
98 | |
---|
99 | |
---|
100 | do n=1,ngas |
---|
101 | read(111,*) gastype(n) |
---|
102 | print*,'Gas ',n,' is ',gastype(n) |
---|
103 | enddo |
---|
104 | |
---|
105 | ! check the array size is correct, load the coefficients |
---|
106 | open(111,file=file_path(1:LEN_TRIM(file_path)),form='formatted') |
---|
107 | read(111,*) L_REFVARcheck |
---|
108 | if(.not.(L_REFVARcheck.eq.L_REFVAR)) then |
---|
109 | print*, L_REFVARcheck |
---|
110 | print*, L_REFVAR |
---|
111 | print*,'The size of your radiative transfer mixing ratio array does ' |
---|
112 | print*,'not match the value given in Q.dat, exiting.' |
---|
113 | print*,'Check the value of L_NREFVAR in radinc_h.F90.' |
---|
114 | call abort |
---|
115 | endif |
---|
116 | read(111,*) wrefvar |
---|
117 | close(111) |
---|
118 | |
---|
119 | ! Check that gastype and gnom match |
---|
120 | do n=1,ngas |
---|
121 | print*,'Gas ',n,' is ',gnom(n) |
---|
122 | if(gnom(n).ne.gastype(n))then |
---|
123 | print*,'Name of a gas in radiative transfer data (',gastype(n),') does not ', & |
---|
124 | 'match that in gases.def (',gnom(n),'), exiting. You should compare ', & |
---|
125 | 'gases.def with Q.dat in your radiative transfer directory.' |
---|
126 | call abort |
---|
127 | endif |
---|
128 | enddo |
---|
129 | print*,'Confirmed gas match in radiative transfer and gases.def!' |
---|
130 | |
---|
131 | ! display the values |
---|
132 | print*,'Variable gas mixing ratios:' |
---|
133 | do n=1,L_REFVAR |
---|
134 | !print*,n,'.',wrefvar(n),' kg/kg' ! pay attention! |
---|
135 | print*,n,'.',wrefvar(n),' mol/mol' |
---|
136 | end do |
---|
137 | print*,'' |
---|
138 | |
---|
139 | !======================================================================= |
---|
140 | ! Set the weighting in g-space |
---|
141 | |
---|
142 | file_id='/corrk_data/' // TRIM(corrkdir) // '/g.dat' |
---|
143 | file_path=TRIM(datadir)//TRIM(file_id) |
---|
144 | |
---|
145 | ! check that the file exists |
---|
146 | inquire(FILE=file_path,EXIST=file_ok) |
---|
147 | if(.not.file_ok) then |
---|
148 | write(*,*)'The file ',TRIM(file_path) |
---|
149 | write(*,*)'was not found by sugas_corrk.F90, exiting.' |
---|
150 | write(*,*)'Check that your path to datagcm:',trim(datadir) |
---|
151 | write(*,*)' is correct. You can change it in callphys.def with:' |
---|
152 | write(*,*)' datadir = /absolute/path/to/datagcm' |
---|
153 | write(*,*)'Also check that the corrkdir you chose in callphys.def exists.' |
---|
154 | call abort |
---|
155 | endif |
---|
156 | |
---|
157 | ! check the array size is correct, load the coefficients |
---|
158 | open(111,file=TRIM(file_path),form='formatted') |
---|
159 | read(111,*) L_NGAUSScheck |
---|
160 | if(.not.(L_NGAUSScheck.eq.L_NGAUSS)) then |
---|
161 | print*,'The size of your radiative transfer g-space array does ' |
---|
162 | print*,'not match the value given in g.dat, exiting.' |
---|
163 | call abort |
---|
164 | endif |
---|
165 | read(111,*) gweight |
---|
166 | close(111) |
---|
167 | |
---|
168 | ! display the values |
---|
169 | print*,'Correlated-k g-space grid:' |
---|
170 | do n=1,L_NGAUSS |
---|
171 | print*,n,'.',gweight(n) |
---|
172 | end do |
---|
173 | print*,'' |
---|
174 | |
---|
175 | !======================================================================= |
---|
176 | ! Set the reference pressure and temperature arrays. These are |
---|
177 | ! the pressures and temperatures at which we have k-coefficients. |
---|
178 | |
---|
179 | !----------------------------------------------------------------------- |
---|
180 | ! pressure |
---|
181 | |
---|
182 | file_id='/corrk_data/' // TRIM(corrkdir) // '/p.dat' |
---|
183 | file_path=TRIM(datadir)//TRIM(file_id) |
---|
184 | |
---|
185 | ! check that the file exists |
---|
186 | inquire(FILE=file_path,EXIST=file_ok) |
---|
187 | if(.not.file_ok) then |
---|
188 | write(*,*)'The file ',TRIM(file_path) |
---|
189 | write(*,*)'was not found by sugas_corrk.F90, exiting.' |
---|
190 | write(*,*)'Check that your path to datagcm:',trim(datadir) |
---|
191 | write(*,*)' is correct. You can change it in callphys.def with:' |
---|
192 | write(*,*)' datadir = /absolute/path/to/datagcm' |
---|
193 | write(*,*)'Also check that the corrkdir you chose in callphys.def exists.' |
---|
194 | call abort |
---|
195 | endif |
---|
196 | |
---|
197 | ! check the array size is correct, load the coefficients |
---|
198 | open(111,file=TRIM(file_path),form='formatted') |
---|
199 | read(111,*) L_NPREFcheck |
---|
200 | if(.not.(L_NPREFcheck.eq.L_NPREF)) then |
---|
201 | print*,'The size of your radiative transfer pressure array does ' |
---|
202 | print*,'not match the value given in p.dat, exiting.' |
---|
203 | call abort |
---|
204 | endif |
---|
205 | read(111,*) pgasref |
---|
206 | close(111) |
---|
207 | |
---|
208 | ! display the values |
---|
209 | print*,'Correlated-k pressure grid (mBar):' |
---|
210 | do n=1,L_NPREF |
---|
211 | print*,n,'. 1 x 10^',pgasref(n),' mBar' |
---|
212 | end do |
---|
213 | print*,'' |
---|
214 | |
---|
215 | ! save the min / max matrix values |
---|
216 | pgasmin = 10.0**pgasref(1) |
---|
217 | pgasmax = 10.0**pgasref(L_NPREF) |
---|
218 | |
---|
219 | ! interpolate to finer grid |
---|
220 | do n=1,L_NPREF-1 |
---|
221 | do m=1,5 |
---|
222 | pfgasref((n-1)*5+m) = pgasref(n)+(m-1)*0.2 |
---|
223 | end do |
---|
224 | end do |
---|
225 | pfgasref(L_PINT) = pgasref(L_NPREF) |
---|
226 | print*,'Warning, pfgasref needs generalisation to uneven grids!!' |
---|
227 | |
---|
228 | !----------------------------------------------------------------------- |
---|
229 | ! temperature |
---|
230 | |
---|
231 | file_id='/corrk_data/' // TRIM(corrkdir) // '/T.dat' |
---|
232 | file_path=TRIM(datadir)//TRIM(file_id) |
---|
233 | |
---|
234 | ! check that the file exists |
---|
235 | inquire(FILE=file_path,EXIST=file_ok) |
---|
236 | if(.not.file_ok) then |
---|
237 | write(*,*)'The file ',TRIM(file_path) |
---|
238 | write(*,*)'was not found by sugas_corrk.F90, exiting.' |
---|
239 | write(*,*)'Check that your path to datagcm:',trim(datadir) |
---|
240 | write(*,*)' is correct. You can change it in callphys.def with:' |
---|
241 | write(*,*)' datadir = /absolute/path/to/datagcm' |
---|
242 | write(*,*)'Also check that the corrkdir you chose in callphys.def exists.' |
---|
243 | call abort |
---|
244 | endif |
---|
245 | |
---|
246 | ! check the array size is correct, load the coefficients |
---|
247 | open(111,file=TRIM(file_path),form='formatted') |
---|
248 | read(111,*) L_NTREFcheck |
---|
249 | if(.not.(L_NTREFcheck.eq.L_NTREF)) then |
---|
250 | print*,'The size of your radiative transfer temperature array does ' |
---|
251 | print*,'not match the value given in T.dat, exiting.' |
---|
252 | call abort |
---|
253 | endif |
---|
254 | read(111,*) tgasref |
---|
255 | close(111) |
---|
256 | |
---|
257 | ! display the values |
---|
258 | print*,'Correlated-k temperature grid:' |
---|
259 | do n=1,L_NTREF |
---|
260 | print*,n,'.',tgasref(n),' K' |
---|
261 | end do |
---|
262 | |
---|
263 | ! save the min / max matrix values |
---|
264 | tgasmin = tgasref(1) |
---|
265 | tgasmax = tgasref(L_NTREF) |
---|
266 | |
---|
267 | !======================================================================= |
---|
268 | ! Get gaseous k-coefficients and interpolate onto finer pressure grid |
---|
269 | |
---|
270 | ! VISIBLE |
---|
271 | if (callgasvis.and..not.TRIM(corrkdir).eq.'null') then |
---|
272 | file_id='/corrk_data/'//trim(adjustl(banddir))//'/corrk_gcm_VI.dat' |
---|
273 | file_path=TRIM(datadir)//TRIM(file_id) |
---|
274 | |
---|
275 | ! check that the file exists |
---|
276 | inquire(FILE=file_path,EXIST=file_ok) |
---|
277 | if(.not.file_ok) then |
---|
278 | write(*,*)'The file ',TRIM(file_path) |
---|
279 | write(*,*)'was not found by sugas_corrk.F90.' |
---|
280 | write(*,*)'Are you sure you have absorption data for these bands?' |
---|
281 | call abort |
---|
282 | endif |
---|
283 | |
---|
284 | open(111,file=TRIM(file_path),form='formatted') |
---|
285 | read(111,*) gasv8 |
---|
286 | close(111) |
---|
287 | |
---|
288 | else |
---|
289 | print*,'Visible corrk gaseous absorption is set to zero.' |
---|
290 | gasv8(:,:,:,:,:)=0.0 |
---|
291 | endif |
---|
292 | |
---|
293 | ! INFRA-RED |
---|
294 | if (.not.TRIM(corrkdir).eq.'null') then |
---|
295 | file_id='/corrk_data/'//trim(adjustl(banddir))//'/corrk_gcm_IR.dat' |
---|
296 | file_path=TRIM(datadir)//TRIM(file_id) |
---|
297 | |
---|
298 | ! check that the file exists |
---|
299 | inquire(FILE=file_path,EXIST=file_ok) |
---|
300 | if(.not.file_ok) then |
---|
301 | write(*,*)'The file ',TRIM(file_path) |
---|
302 | write(*,*)'was not found by sugas_corrk.F90.' |
---|
303 | write(*,*)'Are you sure you have absorption data for these bands?' |
---|
304 | call abort |
---|
305 | endif |
---|
306 | |
---|
307 | open(111,file=TRIM(file_path),form='formatted') |
---|
308 | read(111,*) gasi8 |
---|
309 | close(111) |
---|
310 | |
---|
311 | do nw=1,L_NSPECTI |
---|
312 | fzeroi(nw) = 0. |
---|
313 | ! do nt=1,L_NTREF |
---|
314 | ! do np=1,L_NPREF |
---|
315 | ! do nh=1,L_REFVAR |
---|
316 | ! do ng = 1,L_NGAUSS |
---|
317 | ! if(gasi8(nt,np,nh,nw,ng).lt.1.0e-25)then |
---|
318 | ! fzeroi(nw)=fzeroi(nw)+1. |
---|
319 | ! endif |
---|
320 | ! end do |
---|
321 | ! end do |
---|
322 | ! end do |
---|
323 | ! end do |
---|
324 | ! fzeroi(nw)=fzeroi(nw)/dble(L_NTREF*L_NPREF*L_REFVAR*L_NGAUSS) |
---|
325 | end do |
---|
326 | |
---|
327 | do nw=1,L_NSPECTV |
---|
328 | fzerov(nw) = 0. |
---|
329 | ! do nt=1,L_NTREF |
---|
330 | ! do np=1,L_NPREF |
---|
331 | ! do nh=1,L_REFVAR |
---|
332 | ! do ng = 1,L_NGAUSS |
---|
333 | ! if(gasv8(nt,np,nh,nw,ng).lt.1.0e-25)then |
---|
334 | ! fzerov(nw)=fzerov(nw)+1. |
---|
335 | ! endif |
---|
336 | ! end do |
---|
337 | ! end do |
---|
338 | ! end do |
---|
339 | ! end do |
---|
340 | ! fzerov(nw)=fzerov(nw)/dble(L_NTREF*L_NPREF*L_REFVAR*L_NGAUSS) |
---|
341 | end do |
---|
342 | |
---|
343 | |
---|
344 | do nw=1,L_NSPECTV |
---|
345 | fzerov(nw) = 0. |
---|
346 | end do |
---|
347 | |
---|
348 | else |
---|
349 | print*,'Infrared corrk gaseous absorption is set to zero.' |
---|
350 | gasi8(:,:,:,:,:)=0.0 |
---|
351 | endif |
---|
352 | |
---|
353 | ! Take log10 of the values - this is what we will interpolate. |
---|
354 | ! Smallest value is 1.0E-200. |
---|
355 | |
---|
356 | do nt=1,L_NTREF |
---|
357 | do np=1,L_NPREF |
---|
358 | do nh=1,L_REFVAR |
---|
359 | do ng = 1,L_NGAUSS |
---|
360 | |
---|
361 | do nw=1,L_NSPECTV |
---|
362 | if(gasv8(nt,np,nh,nw,ng).gt.1.0d-200) then |
---|
363 | gasv8(nt,np,nh,nw,ng) = & |
---|
364 | log10(gasv8(nt,np,nh,nw,ng)) |
---|
365 | else |
---|
366 | gasv8(nt,np,nh,nw,ng) = -200.0 |
---|
367 | end if |
---|
368 | end do |
---|
369 | |
---|
370 | do nw=1,L_NSPECTI |
---|
371 | if(gasi8(nt,np,nh,nw,ng).gt.1.0d-200) then |
---|
372 | gasi8(nt,np,nh,nw,ng) = & |
---|
373 | log10(gasi8(nt,np,nh,nw,ng)) |
---|
374 | else |
---|
375 | gasi8(nt,np,nh,nw,ng) = -200.0 |
---|
376 | end if |
---|
377 | end do |
---|
378 | |
---|
379 | end do |
---|
380 | end do |
---|
381 | end do |
---|
382 | end do |
---|
383 | |
---|
384 | ! Interpolate the values: first the longwave |
---|
385 | |
---|
386 | do nt=1,L_NTREF |
---|
387 | do nh=1,L_REFVAR |
---|
388 | do nw=1,L_NSPECTI |
---|
389 | do ng=1,L_NGAUSS |
---|
390 | |
---|
391 | ! First, the initial interval |
---|
392 | |
---|
393 | n = 1 |
---|
394 | do m=1,5 |
---|
395 | x = pfgasref(m) |
---|
396 | xi(1) = pgasref(n) |
---|
397 | xi(2) = pgasref(n+1) |
---|
398 | xi(3) = pgasref(n+2) |
---|
399 | xi(4) = pgasref(n+3) |
---|
400 | yi(1) = gasi8(nt,n,nh,nw,ng) |
---|
401 | yi(2) = gasi8(nt,n+1,nh,nw,ng) |
---|
402 | yi(3) = gasi8(nt,n+2,nh,nw,ng) |
---|
403 | yi(4) = gasi8(nt,n+3,nh,nw,ng) |
---|
404 | call lagrange(x,xi,yi,ans) |
---|
405 | gasi(nt,m,nh,nw,ng) = 10.0**ans |
---|
406 | end do |
---|
407 | |
---|
408 | do n=2,L_NPREF-2 |
---|
409 | do m=1,5 |
---|
410 | i = (n-1)*5+m |
---|
411 | x = pfgasref(i) |
---|
412 | xi(1) = pgasref(n-1) |
---|
413 | xi(2) = pgasref(n) |
---|
414 | xi(3) = pgasref(n+1) |
---|
415 | xi(4) = pgasref(n+2) |
---|
416 | yi(1) = gasi8(nt,n-1,nh,nw,ng) |
---|
417 | yi(2) = gasi8(nt,n,nh,nw,ng) |
---|
418 | yi(3) = gasi8(nt,n+1,nh,nw,ng) |
---|
419 | yi(4) = gasi8(nt,n+2,nh,nw,ng) |
---|
420 | call lagrange(x,xi,yi,ans) |
---|
421 | gasi(nt,i,nh,nw,ng) = 10.0**ans |
---|
422 | end do |
---|
423 | end do |
---|
424 | |
---|
425 | ! Now, get the last interval |
---|
426 | |
---|
427 | n = L_NPREF-1 |
---|
428 | do m=1,5 |
---|
429 | i = (n-1)*5+m |
---|
430 | x = pfgasref(i) |
---|
431 | xi(1) = pgasref(n-2) |
---|
432 | xi(2) = pgasref(n-1) |
---|
433 | xi(3) = pgasref(n) |
---|
434 | xi(4) = pgasref(n+1) |
---|
435 | yi(1) = gasi8(nt,n-2,nh,nw,ng) |
---|
436 | yi(2) = gasi8(nt,n-1,nh,nw,ng) |
---|
437 | yi(3) = gasi8(nt,n,nh,nw,ng) |
---|
438 | yi(4) = gasi8(nt,n+1,nh,nw,ng) |
---|
439 | call lagrange(x,xi,yi,ans) |
---|
440 | gasi(nt,i,nh,nw,ng) = 10.0**ans |
---|
441 | end do |
---|
442 | |
---|
443 | ! Fill the last pressure point |
---|
444 | |
---|
445 | gasi(nt,L_PINT,nh,nw,ng) = & |
---|
446 | 10.0**gasi8(nt,L_NPREF,nh,nw,ng) |
---|
447 | |
---|
448 | end do |
---|
449 | end do |
---|
450 | end do |
---|
451 | end do |
---|
452 | |
---|
453 | ! Interpolate the values: now the shortwave |
---|
454 | |
---|
455 | do nt=1,L_NTREF |
---|
456 | do nh=1,L_REFVAR |
---|
457 | do nw=1,L_NSPECTV |
---|
458 | do ng=1,L_NGAUSS |
---|
459 | |
---|
460 | ! First, the initial interval |
---|
461 | |
---|
462 | n = 1 |
---|
463 | do m=1,5 |
---|
464 | x = pfgasref(m) |
---|
465 | xi(1) = pgasref(n) |
---|
466 | xi(2) = pgasref(n+1) |
---|
467 | xi(3) = pgasref(n+2) |
---|
468 | xi(4) = pgasref(n+3) |
---|
469 | yi(1) = gasv8(nt,n,nh,nw,ng) |
---|
470 | yi(2) = gasv8(nt,n+1,nh,nw,ng) |
---|
471 | yi(3) = gasv8(nt,n+2,nh,nw,ng) |
---|
472 | yi(4) = gasv8(nt,n+3,nh,nw,ng) |
---|
473 | call lagrange(x,xi,yi,ans) |
---|
474 | gasv(nt,m,nh,nw,ng) = 10.0**ans |
---|
475 | end do |
---|
476 | |
---|
477 | do n=2,L_NPREF-2 |
---|
478 | do m=1,5 |
---|
479 | i = (n-1)*5+m |
---|
480 | x = pfgasref(i) |
---|
481 | xi(1) = pgasref(n-1) |
---|
482 | xi(2) = pgasref(n) |
---|
483 | xi(3) = pgasref(n+1) |
---|
484 | xi(4) = pgasref(n+2) |
---|
485 | yi(1) = gasv8(nt,n-1,nh,nw,ng) |
---|
486 | yi(2) = gasv8(nt,n,nh,nw,ng) |
---|
487 | yi(3) = gasv8(nt,n+1,nh,nw,ng) |
---|
488 | yi(4) = gasv8(nt,n+2,nh,nw,ng) |
---|
489 | call lagrange(x,xi,yi,ans) |
---|
490 | gasv(nt,i,nh,nw,ng) = 10.0**ans |
---|
491 | end do |
---|
492 | end do |
---|
493 | |
---|
494 | ! Now, get the last interval |
---|
495 | |
---|
496 | n = L_NPREF-1 |
---|
497 | do m=1,5 |
---|
498 | i = (n-1)*5+m |
---|
499 | x = pfgasref(i) |
---|
500 | xi(1) = pgasref(n-2) |
---|
501 | xi(2) = pgasref(n-1) |
---|
502 | xi(3) = pgasref(n) |
---|
503 | xi(4) = pgasref(n+1) |
---|
504 | yi(1) = gasv8(nt,n-2,nh,nw,ng) |
---|
505 | yi(2) = gasv8(nt,n-1,nh,nw,ng) |
---|
506 | yi(3) = gasv8(nt,n,nh,nw,ng) |
---|
507 | yi(4) = gasv8(nt,n+1,nh,nw,ng) |
---|
508 | call lagrange(x,xi,yi,ans) |
---|
509 | gasv(nt,i,nh,nw,ng) = 10.0**ans |
---|
510 | end do |
---|
511 | |
---|
512 | ! Fill the last pressure point |
---|
513 | |
---|
514 | gasv(nt,L_PINT,nh,nw,ng) = & |
---|
515 | 10.0**gasv8(nt,L_NPREF,nh,nw,ng) |
---|
516 | |
---|
517 | end do |
---|
518 | end do |
---|
519 | end do |
---|
520 | end do |
---|
521 | |
---|
522 | |
---|
523 | |
---|
524 | do igas=1,ngasmx |
---|
525 | if(gnom(igas).eq.'H2_')then |
---|
526 | call interpolateH2H2(500.D+0,250.D+0,17500.D+0,testcont,.true.) |
---|
527 | elseif(gnom(igas).eq.'H2O')then |
---|
528 | call interpolateH2Ocont(990.D+0,296.D+0,683.2D+0*2,0.D+0,testcont,.true.) |
---|
529 | endif |
---|
530 | enddo |
---|
531 | |
---|
532 | return |
---|
533 | end subroutine sugas_corrk |
---|