source: codes/icosagcm/branches/SATURN_DYNAMICO/LMDZ.COMMON/libf/phystd/sugas_corrk.F90 @ 264

Last change on this file since 264 was 227, checked in by milmd, 10 years ago

Last LMDZ version (1315) with OpenMP directives and other stuff

File size: 24.0 KB
Line 
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!     Added double gray case by Jeremy Leconte (2012)
15!     New HITRAN continuum data section by RW (2012)
16!
17!     Summary
18!     -------
19!
20!==================================================================
21
22      use radinc_h
23      use radcommon_h, only : pgasref,pfgasref,pgasmin,pgasmax
24      use radcommon_h, only : tgasref,tgasmin,tgasmax
25      use radcommon_h, only : gasv,gasi,FZEROI,FZEROV,gweight
26      use radcommon_h, only : wrefvar,WNOI,WNOV
27      use datafile_mod, only: datadir
28
29      use gases_h
30!      use ioipsl_getincom
31      use ioipsl_getincom_p 
32      implicit none
33
34#include "callkeys.h"
35#include "comcstfi.h"
36
37!==================================================================
38
39      logical file_ok
40
41      integer n, nt, np, nh, ng, nw, m, i
42      integer L_NGAUSScheck
43
44      character(len=200) :: file_id
45      character(len=500) :: file_path
46
47      ! ALLOCATABLE ARRAYS -- AS 12/2011
48      REAL*8, DIMENSION(:,:,:,:,:), ALLOCATABLE,SAVE :: gasi8, gasv8    !read by master
49      character*20,allocatable,DIMENSION(:),SAVE :: gastype ! for check with gnom, read by master
50
51      real*8 x, xi(4), yi(4), ans, p
52!     For gray case (JL12)
53      real kappa_IR, kappa_VI, IR_VI_wnlimit
54      integer nVI_limit,nIR_limit
55
56      integer ngas, igas, jgas
57
58      double precision testcont ! for continuum absorption initialisation
59
60      integer :: dummy
61
62!=======================================================================
63!     Load variable species data, exit if we have wrong database
64      file_id='/corrk_data/' // TRIM(corrkdir) // '/Q.dat'
65      file_path=TRIM(datadir)//TRIM(file_id)
66
67      ! check that the file exists
68      inquire(FILE=file_path,EXIST=file_ok)
69      if(.not.file_ok) then
70         write(*,*)'The file ',TRIM(file_path)
71         write(*,*)'was not found by sugas_corrk.F90, exiting.'
72         write(*,*)'Check that your path to datagcm:',trim(datadir)
73         write(*,*)' is correct. You can change it in callphys.def with:'
74         write(*,*)' datadir = /absolute/path/to/datagcm'
75         write(*,*)'Also check that the corrkdir you chose in callphys.def exists.'
76         call abort
77      endif
78
79!$OMP MASTER
80      ! check that database matches varactive toggle
81      open(111,file=TRIM(file_path),form='formatted')
82      read(111,*) ngas
83
84      if(ngas.ne.ngasmx)then
85         print*,'Number of gases in radiative transfer data (',ngas,') does not', &
86                'match that in gases.def (',ngasmx,'), exiting.'
87         call abort
88      endif
89
90      if(ngas.eq.1 .and. (varactive.or.varfixed))then
91         print*,'You have varactive/fixed=.true. but the database [',  &
92                     corrkdir(1:LEN_TRIM(corrkdir)),           &
93                '] has no variable species, exiting.'
94         call abort
95      elseif(ngas.gt.5 .or. ngas.lt.1)then
96         print*,ngas,' species in database [',               &
97                     corrkdir(1:LEN_TRIM(corrkdir)),           &
98                '], radiative code cannot handle this.'
99         call abort
100      endif 
101
102      ! dynamically allocate gastype and read from Q.dat
103      IF ( .NOT. ALLOCATED( gastype ) ) ALLOCATE( gastype( ngas ) )
104
105      do igas=1,ngas
106         read(111,*) gastype(igas)
107         print*,'Gas ',igas,' is ',gastype(igas)
108      enddo
109
110      ! get array size, load the coefficients
111      open(111,file=TRIM(file_path),form='formatted')
112      read(111,*) L_REFVAR
113      IF( .NOT. ALLOCATED( wrefvar ) ) ALLOCATE( WREFVAR(L_REFVAR) )
114      read(111,*) wrefvar
115      close(111)
116
117      if(L_REFVAR.gt.1 .and. (.not.varactive) .and. (.not.varfixed))then
118         print*,'You have varactive and varfixed=.false. and the database [', &
119                     corrkdir(1:LEN_TRIM(corrkdir)),           &
120                '] has a variable species.'
121         call abort
122      endif
123
124      ! Check that gastype and gnom match
125      do igas=1,ngas
126         print*,'Gas ',igas,' is ',trim(gnom(igas))
127         if (trim(gnom(igas)).ne.trim(gastype(igas))) then
128            print*,'Name of a gas in radiative transfer data (',trim(gastype(igas)),') does not ', &
129                 'match that in gases.def (',trim(gnom(igas)),'), exiting. You should compare ', &
130                 'gases.def with Q.dat in your radiative transfer directory.' 
131            call abort
132         endif
133      enddo
134      print*,'Confirmed gas match in radiative transfer and gases.def!'
135
136      ! display the values
137      print*,'Variable gas volume mixing ratios:'
138      do n=1,L_REFVAR
139         !print*,n,'.',wrefvar(n),' kg/kg' ! pay attention!
140         print*,n,'.',wrefvar(n),' mol/mol'
141      end do
142      print*,''
143
144!=======================================================================
145!     Set the weighting in g-space
146
147      file_id='/corrk_data/' // TRIM(corrkdir) // '/g.dat'
148      file_path=TRIM(datadir)//TRIM(file_id)
149
150      ! check that the file exists
151      inquire(FILE=file_path,EXIST=file_ok)
152      if(.not.file_ok) then
153         write(*,*)'The file ',TRIM(file_path)
154         write(*,*)'was not found by sugas_corrk.F90, exiting.'
155         write(*,*)'Check that your path to datagcm:',trim(datadir)
156         write(*,*)' is correct. You can change it in callphys.def with:'
157         write(*,*)' datadir = /absolute/path/to/datagcm'
158         write(*,*)'Also check that the corrkdir you chose in callphys.def exists.'
159         call abort
160      endif
161     
162      ! check the array size is correct, load the coefficients
163      open(111,file=TRIM(file_path),form='formatted')
164      read(111,*) L_NGAUSScheck
165      if(.not.(L_NGAUSScheck.eq.L_NGAUSS)) then
166         print*,'The size of your radiative transfer g-space array does '
167         print*,'not match the value given in g.dat, exiting.'
168         call abort
169      endif
170      read(111,*) gweight
171      close(111)
172 
173      ! display the values
174      print*,'Correlated-k g-space grid:'
175      do n=1,L_NGAUSS
176         print*,n,'.',gweight(n)
177      end do
178      print*,''
179
180!=======================================================================
181!     Set the reference pressure and temperature arrays.  These are
182!     the pressures and temperatures at which we have k-coefficients.
183
184!-----------------------------------------------------------------------
185! pressure
186
187      file_id='/corrk_data/' // TRIM(corrkdir) // '/p.dat'
188      file_path=TRIM(datadir)//TRIM(file_id)
189
190      ! check that the file exists
191      inquire(FILE=file_path,EXIST=file_ok)
192      if(.not.file_ok) then
193         write(*,*)'The file ',TRIM(file_path)
194         write(*,*)'was not found by sugas_corrk.F90, exiting.'
195         write(*,*)'Check that your path to datagcm:',trim(datadir)
196         write(*,*)' is correct. You can change it in callphys.def with:'
197         write(*,*)' datadir = /absolute/path/to/datagcm'
198         write(*,*)'Also check that the corrkdir you chose in callphys.def exists.'
199         call abort
200      endif
201     
202      ! get array size, load the coefficients
203      open(111,file=TRIM(file_path),form='formatted')
204      read(111,*) L_NPREF
205      IF( .NOT. ALLOCATED( pgasref ) ) ALLOCATE( PGASREF(L_NPREF) )
206      read(111,*) pgasref
207      close(111)
208      L_PINT = (L_NPREF-1)*5+1
209      IF( .NOT. ALLOCATED( pfgasref ) ) ALLOCATE( PFGASREF(L_PINT) )
210
211      ! display the values
212      print*,'Correlated-k pressure grid (mBar):'
213      do n=1,L_NPREF
214         print*,n,'. 1 x 10^',pgasref(n),' mBar'
215      end do
216      print*,''
217
218      ! save the min / max matrix values
219      pgasmin = 10.0**pgasref(1)
220      pgasmax = 10.0**pgasref(L_NPREF)
221
222      ! interpolate to finer grid
223      do n=1,L_NPREF-1
224         do m=1,5
225            pfgasref((n-1)*5+m) = pgasref(n)+(m-1)*0.2
226         end do
227      end do
228      pfgasref(L_PINT) = pgasref(L_NPREF)
229      ! Warning, pfgasref needs generalisation to uneven grids!
230
231!-----------------------------------------------------------------------
232! temperature
233
234      file_id='/corrk_data/' // TRIM(corrkdir) // '/T.dat'
235      file_path=TRIM(datadir)//TRIM(file_id)
236
237      ! check that the file exists
238      inquire(FILE=file_path,EXIST=file_ok)
239      if(.not.file_ok) then
240         write(*,*)'The file ',TRIM(file_path)
241         write(*,*)'was not found by sugas_corrk.F90, exiting.'
242         write(*,*)'Check that your path to datagcm:',trim(datadir)
243         write(*,*)' is correct. You can change it in callphys.def with:'
244         write(*,*)' datadir = /absolute/path/to/datagcm'
245         write(*,*)'Also check that the corrkdir you chose in callphys.def exists.'
246         call abort
247      endif
248
249      ! get array size, load the coefficients
250      open(111,file=TRIM(file_path),form='formatted')
251      read(111,*) L_NTREF
252      IF( .NOT. ALLOCATED( tgasref ) ) ALLOCATE( TGASREF(L_NTREF) )
253      read(111,*) tgasref
254      close(111)
255
256      ! display the values
257      print*,'Correlated-k temperature grid:'
258      do n=1,L_NTREF
259         print*,n,'.',tgasref(n),' K'
260      end do
261
262      ! save the min / max matrix values
263      tgasmin = tgasref(1)
264      tgasmax = tgasref(L_NTREF)
265!$OMP END MASTER
266!$OMP BARRIER
267
268!-----------------------------------------------------------------------
269! allocate the multidimensional arrays in radcommon_h
270
271      IF( .NOT. ALLOCATED( gasi8 ) ) ALLOCATE( gasi8(L_NTREF,L_NPREF,L_REFVAR,L_NSPECTI,L_NGAUSS) )
272      IF( .NOT. ALLOCATED( gasv8 ) ) ALLOCATE( gasv8(L_NTREF,L_NPREF,L_REFVAR,L_NSPECTV,L_NGAUSS) )
273      IF( .NOT. ALLOCATED( gasi ) ) ALLOCATE( gasi(L_NTREF,L_PINT,L_REFVAR,L_NSPECTI,L_NGAUSS) )
274      IF( .NOT. ALLOCATED( gasv ) ) ALLOCATE( gasv(L_NTREF,L_PINT,L_REFVAR,L_NSPECTV,L_NGAUSS) )
275
276      ! display the values
277      print*,''
278      print*,'Correlated-k matrix size:' 
279      print*,'[',L_NTREF,',',L_NPREF,',',L_REFVAR,',',L_NGAUSS,']' 
280
281!=======================================================================
282!     Get gaseous k-coefficients and interpolate onto finer pressure grid
283
284
285!        wavelength used to separate IR from VI in graybody. We will need that anyway
286         IR_VI_wnlimit=3000.
287         write(*,*)"graybody: Visible / Infrared separation set at",10000./IR_VI_wnlimit,"um"
288         
289         nVI_limit=0
290         Do nw=1,L_NSPECTV
291            if ((WNOV(nw).gt.IR_VI_wnlimit).and.(L_NSPECTV.gt.1)) then
292               nVI_limit=nw-1
293               exit
294            endif
295         End do
296         nIR_limit=L_NSPECTI
297         Do nw=1,L_NSPECTI
298            if ((WNOI(nw).gt.IR_VI_wnlimit).and.(L_NSPECTI.gt.1)) then
299               nIR_limit=nw-1
300               exit
301            endif
302         End do
303
304      if (graybody) then
305!        constant absorption coefficient in visible
306         write(*,*)"graybody: constant absorption coefficient in visible:"
307         kappa_VI=-100000.
308         call getin_p("kappa_VI",kappa_VI)
309         write(*,*)" kappa_VI = ",kappa_VI
310         kappa_VI=kappa_VI*1.e4* mugaz * 1.672621e-27    ! conversion from m^2/kg to cm^2/molecule         
311     
312!        constant absorption coefficient in IR
313         write(*,*)"graybody: constant absorption coefficient in InfraRed:"
314         kappa_IR=-100000.
315         call getin_p("kappa_IR",kappa_IR)
316         write(*,*)" kappa_IR = ",kappa_IR       
317         kappa_IR=kappa_IR*1.e4* mugaz * 1.672621e-27    ! conversion from m^2/kg to cm^2/molecule
318
319         write(*,*)"graybody: Visible / Infrared separation set at band: IR=",nIR_limit,", VI=",nVI_limit
320               
321      Else
322         kappa_VI=1.e-30     
323         kappa_IR=1.e-30       
324      End if
325
326!$OMP MASTER         
327!      print*,corrkdir(1:4)
328      ! VISIBLE
329      if (callgasvis) then
330         if ((corrkdir(1:4).eq.'null'))then   !(TRIM(corrkdir).eq.'null_LowTeffStar')) then
331            gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTV,1:L_NGAUSS)=0.0
332            print*,'using no corrk data'
333            print*,'Visible corrk gaseous absorption is set to zero if graybody=F'
334         else
335            file_id='/corrk_data/'//trim(adjustl(banddir))//'/corrk_gcm_VI.dat' 
336            file_path=TRIM(datadir)//TRIM(file_id)
337           
338            ! check that the file exists
339            inquire(FILE=file_path,EXIST=file_ok)
340            if(.not.file_ok) then
341               write(*,*)'The file ',TRIM(file_path)
342               write(*,*)'was not found by sugas_corrk.F90.'
343               write(*,*)'Are you sure you have absorption data for these bands?'
344               call abort
345            endif
346         
347            open(111,file=TRIM(file_path),form='formatted')
348            read(111,*) gasv8
349            close(111)
350         end if
351
352         if(nVI_limit.eq.0) then
353            gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTV,1:L_NGAUSS)=   &
354                  gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTV,1:L_NGAUSS)+kappa_VI
355         else if (nVI_limit.eq.L_NSPECTV) then
356            gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTV,1:L_NGAUSS)=   &
357                  gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTV,1:L_NGAUSS)+kappa_IR
358         else
359            gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:nVI_limit,1:L_NGAUSS)=   &
360                  gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:nVI_limit,1:L_NGAUSS)+kappa_IR
361            gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,nVI_limit+1:L_NSPECTV,1:L_NGAUSS)=   &
362                  gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,nVI_limit+1:L_NSPECTV,1:L_NGAUSS)+kappa_VI
363         end if
364      else
365         print*,'Visible corrk gaseous absorption is set to zero.'
366         gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTV,1:L_NGAUSS)=0.0
367      endif
368!$OMP END MASTER
369!$OMP BARRIER
370
371      ! INFRA-RED
372      if ((corrkdir(1:4).eq.'null'))then       !.or.(TRIM(corrkdir).eq.'null_LowTeffStar')) then
373         print*,'Infrared corrk gaseous absorption is set to zero if graybody=F'
374!$OMP MASTER         
375         gasi8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTI,1:L_NGAUSS)=0.0
376!$OMP END MASTER
377!$OMP BARRIER
378      else
379         file_id='/corrk_data/'//trim(adjustl(banddir))//'/corrk_gcm_IR.dat' 
380         file_path=TRIM(datadir)//TRIM(file_id)
381     
382         ! check that the file exists
383         inquire(FILE=file_path,EXIST=file_ok)
384         if(.not.file_ok) then
385            write(*,*)'The file ',TRIM(file_path)
386            write(*,*)'was not found by sugas_corrk.F90.'
387            write(*,*)'Are you sure you have absorption data for these bands?'
388            call abort
389         endif
390         
391!$OMP MASTER                   
392         open(111,file=TRIM(file_path),form='formatted')
393         read(111,*) gasi8
394         close(111)
395!$OMP END MASTER
396!$OMP BARRIER
397     
398         ! 'fzero' is a currently unused feature that allows optimisation
399         ! of the radiative transfer by neglecting bands where absorption
400         ! is close to zero. As it could be useful in the future, this
401         ! section of the code has been kept commented and not erased.
402         ! RW 7/3/12.
403
404         do nw=1,L_NSPECTI
405            fzeroi(nw) = 0.d0
406!            do nt=1,L_NTREF
407!               do np=1,L_NPREF
408!                  do nh=1,L_REFVAR
409!                     do ng = 1,L_NGAUSS
410!                        if(gasi8(nt,np,nh,nw,ng).lt.1.0e-25)then
411!                           fzeroi(nw)=fzeroi(nw)+1.d0
412!                        endif
413!                     end do
414!                  end do
415!               end do
416!            end do
417!            fzeroi(nw)=fzeroi(nw)/dble(L_NTREF*L_NPREF*L_REFVAR*L_NGAUSS)
418         end do
419
420         do nw=1,L_NSPECTV
421            fzerov(nw) = 0.d0
422!            do nt=1,L_NTREF
423!               do np=1,L_NPREF
424!                  do nh=1,L_REFVAR
425!                     do ng = 1,L_NGAUSS
426!                        if(gasv8(nt,np,nh,nw,ng).lt.1.0e-25)then
427!                           fzerov(nw)=fzerov(nw)+1.d0
428!                        endif
429!                     end do
430!                  end do
431!               end do
432!            end do
433!            fzerov(nw)=fzerov(nw)/dble(L_NTREF*L_NPREF*L_REFVAR*L_NGAUSS)
434         end do
435
436      endif
437
438!$OMP MASTER                 
439      if(nIR_limit.eq.0) then
440         gasi8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTI,1:L_NGAUSS)=   &
441                  gasi8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTI,1:L_NGAUSS)+kappa_VI
442      else if (nIR_limit.eq.L_NSPECTI) then
443         gasi8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTI,1:L_NGAUSS)=   &
444                  gasi8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTI,1:L_NGAUSS)+kappa_IR
445      else
446         gasi8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:nIR_limit,1:L_NGAUSS)=   &
447                  gasi8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:nIR_limit,1:L_NGAUSS)+kappa_IR
448         gasi8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,nIR_limit+1:L_NSPECTI,1:L_NGAUSS)=   &
449                  gasi8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,nIR_limit+1:L_NSPECTI,1:L_NGAUSS)+kappa_VI
450      end if
451
452
453!     Take log10 of the values - this is what we will interpolate.
454!     Smallest value is 1.0E-200.
455
456      do nt=1,L_NTREF
457         do np=1,L_NPREF
458            do nh=1,L_REFVAR
459               do ng = 1,L_NGAUSS
460
461                  do nw=1,L_NSPECTV
462                     if(gasv8(nt,np,nh,nw,ng).gt.1.0d-200) then
463                        gasv8(nt,np,nh,nw,ng) = log10(gasv8(nt,np,nh,nw,ng))
464                     else
465                        gasv8(nt,np,nh,nw,ng) = -200.0
466                     end if
467                  end do
468
469                  do nw=1,L_NSPECTI
470                     if(gasi8(nt,np,nh,nw,ng).gt.1.0d-200) then
471                        gasi8(nt,np,nh,nw,ng) = log10(gasi8(nt,np,nh,nw,ng))
472                     else
473                        gasi8(nt,np,nh,nw,ng) = -200.0
474                     end if
475                  end do
476                 
477               end do
478            end do
479         end do
480      end do
481!$OMP END MASTER
482!$OMP BARRIER
483
484!     Interpolate the values:  first the longwave
485
486      do nt=1,L_NTREF
487         do nh=1,L_REFVAR
488            do nw=1,L_NSPECTI
489               do ng=1,L_NGAUSS
490
491!     First, the initial interval
492
493                  n = 1 
494                  do m=1,5
495                     x     = pfgasref(m)
496                     xi(1) = pgasref(n)
497                     xi(2) = pgasref(n+1)
498                     xi(3) = pgasref(n+2)
499                     xi(4) = pgasref(n+3)
500                     yi(1) = gasi8(nt,n,nh,nw,ng)
501                     yi(2) = gasi8(nt,n+1,nh,nw,ng)
502                     yi(3) = gasi8(nt,n+2,nh,nw,ng)
503                     yi(4) = gasi8(nt,n+3,nh,nw,ng)
504                     call lagrange(x,xi,yi,ans)
505                     gasi(nt,m,nh,nw,ng) = 10.0**ans
506                  end do
507                 
508                  do n=2,L_NPREF-2
509                     do m=1,5
510                        i     = (n-1)*5+m
511                        x     = pfgasref(i)
512                        xi(1) = pgasref(n-1)
513                        xi(2) = pgasref(n)
514                        xi(3) = pgasref(n+1)
515                        xi(4) = pgasref(n+2)
516                        yi(1) = gasi8(nt,n-1,nh,nw,ng)
517                        yi(2) = gasi8(nt,n,nh,nw,ng)
518                        yi(3) = gasi8(nt,n+1,nh,nw,ng)
519                        yi(4) = gasi8(nt,n+2,nh,nw,ng)
520                        call lagrange(x,xi,yi,ans)
521                        gasi(nt,i,nh,nw,ng) = 10.0**ans
522                     end do
523                  end do
524
525!     Now, get the last interval
526
527                  n = L_NPREF-1                 
528                  do m=1,5
529                     i     = (n-1)*5+m
530                     x     = pfgasref(i)
531                     xi(1) = pgasref(n-2)
532                     xi(2) = pgasref(n-1)
533                     xi(3) = pgasref(n)
534                     xi(4) = pgasref(n+1)
535                     yi(1) = gasi8(nt,n-2,nh,nw,ng)
536                     yi(2) = gasi8(nt,n-1,nh,nw,ng)
537                     yi(3) = gasi8(nt,n,nh,nw,ng)
538                     yi(4) = gasi8(nt,n+1,nh,nw,ng)
539                     call lagrange(x,xi,yi,ans)
540                     gasi(nt,i,nh,nw,ng) = 10.0**ans
541                  end do 
542
543!     Fill the last pressure point
544
545                  gasi(nt,L_PINT,nh,nw,ng) = &
546                       10.0**gasi8(nt,L_NPREF,nh,nw,ng)
547
548               end do
549            end do
550         end do
551      end do
552
553!     Interpolate the values:  now the shortwave
554
555      do nt=1,L_NTREF
556         do nh=1,L_REFVAR
557            do nw=1,L_NSPECTV
558               do ng=1,L_NGAUSS
559
560!     First, the initial interval
561
562                  n = 1 
563                  do m=1,5
564                     x     = pfgasref(m)
565                     xi(1) = pgasref(n)
566                     xi(2) = pgasref(n+1)
567                     xi(3) = pgasref(n+2)
568                     xi(4) = pgasref(n+3)
569                     yi(1) = gasv8(nt,n,nh,nw,ng)
570                     yi(2) = gasv8(nt,n+1,nh,nw,ng)
571                     yi(3) = gasv8(nt,n+2,nh,nw,ng)
572                     yi(4) = gasv8(nt,n+3,nh,nw,ng)
573                     call lagrange(x,xi,yi,ans)
574                     gasv(nt,m,nh,nw,ng) = 10.0**ans
575                  end do
576                 
577                  do n=2,L_NPREF-2
578                     do m=1,5
579                        i     = (n-1)*5+m
580                        x     = pfgasref(i)
581                        xi(1) = pgasref(n-1)
582                        xi(2) = pgasref(n)
583                        xi(3) = pgasref(n+1)
584                        xi(4) = pgasref(n+2)
585                        yi(1) = gasv8(nt,n-1,nh,nw,ng)
586                        yi(2) = gasv8(nt,n,nh,nw,ng)
587                        yi(3) = gasv8(nt,n+1,nh,nw,ng)
588                        yi(4) = gasv8(nt,n+2,nh,nw,ng)
589                        call lagrange(x,xi,yi,ans)
590                        gasv(nt,i,nh,nw,ng) = 10.0**ans
591                     end do
592                  end do
593
594!     Now, get the last interval
595
596                  n = L_NPREF-1
597                  do m=1,5
598                     i     = (n-1)*5+m
599                     x     = pfgasref(i)
600                     xi(1) = pgasref(n-2)
601                     xi(2) = pgasref(n-1)
602                     xi(3) = pgasref(n)
603                     xi(4) = pgasref(n+1)
604                     yi(1) = gasv8(nt,n-2,nh,nw,ng)
605                     yi(2) = gasv8(nt,n-1,nh,nw,ng)
606                     yi(3) = gasv8(nt,n,nh,nw,ng)
607                     yi(4) = gasv8(nt,n+1,nh,nw,ng)
608                     call lagrange(x,xi,yi,ans)
609                     gasv(nt,i,nh,nw,ng) = 10.0**ans
610                  end do 
611
612!     Fill the last pressure point
613
614                  gasv(nt,L_PINT,nh,nw,ng) = &
615                      10.0**gasv8(nt,L_NPREF,nh,nw,ng)
616                 
617               end do
618            end do
619         end do
620      end do
621
622
623!=======================================================================
624!     Initialise the continuum absorption data
625      if(continuum)then
626      do igas=1,ngasmx
627
628         if (igas .eq. igas_N2) then
629
630            dummy = -9999
631            call interpolateN2N2(100.D+0,250.D+0,17500.D+0,testcont,.true.,dummy)
632
633         elseif (igas .eq. igas_H2) then
634
635            ! first do self-induced absorption
636            dummy = -9999
637            call interpolateH2H2(500.D+0,250.D+0,17500.D+0,testcont,.true.,dummy)
638            ! then cross-interactions with other gases
639            do jgas=1,ngasmx
640               if (jgas .eq. igas_N2) then
641                  dummy = -9999
642                  call interpolateN2H2(592.D+0,278.15D+0,200000.D+0,10000.D+0,testcont,.true.,dummy)
643               elseif (jgas .eq. igas_He) then
644                  dummy = -9999
645                  call interpolateH2He(500.D+0,250.D+0,200000.D+0,10000.D+0,testcont,.true.,dummy)
646               endif
647            enddo
648
649         elseif (igas .eq. igas_H2O) then
650
651            ! H2O is special
652            if(H2Ocont_simple)then
653               call interpolateH2Ocont_PPC(990.D+0,296.D+0,683.2D+0*2,0.D+0,testcont,.true.)
654            else
655               dummy = -9999
656               call interpolateH2Ocont_CKD(990.D+0,296.D+0,683.2D+0*2,0.D+0,testcont,.true.,dummy)
657            endif
658
659         endif
660
661      enddo
662      endif
663
664      print*,'----------------------------------------------------'
665      print*,'And that`s all we have. It`s possible that other'
666      print*,'continuum absorption may be present, but if it is we'
667      print*,'don`t yet have data for it...'
668      print*,''
669
670!     Deallocate local arrays
671!$OMP BARRIER
672!$OMP MASTER
673      IF( ALLOCATED( gasi8 ) ) DEALLOCATE( gasi8 )
674      IF( ALLOCATED( gasv8 ) ) DEALLOCATE( gasv8 )
675      IF( ALLOCATED( pgasref ) ) DEALLOCATE( pgasref )
676      IF( ALLOCATED( gastype ) ) DEALLOCATE( gastype )
677!$OMP END MASTER
678!$OMP BARRIER
679
680      return
681    end subroutine sugas_corrk
Note: See TracBrowser for help on using the repository browser.