source: trunk/LMDZ.GENERIC/libf/phystd/sugas_corrk.F90 @ 1477

Last change on this file since 1477 was 1450, checked in by emillour, 9 years ago

Generic GCM:
OpenMP fix: only master allocates/deallocates shared saved variables.
EM

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