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

Last change on this file since 2176 was 2032, checked in by emillour, 6 years ago

Generic GCM:

  • correct a bug introduced in revision 2026; now that L_NGAUSS is a parameter read in via sugas_corrk (called at first call by callcorrk), automatic arrays of size L_NGAUSS cannot be declared in callcorrk and must be allocated once the value of L_NGAUSS has been set.
  • turned optci, optcv and callcorrk into modules in the process.

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, only: corrkdir, banddir, L_NSPECTI, L_NSPECTV, &
23                          L_NGAUSS, L_NPREF, L_NTREF, L_REFVAR, L_PINT
24      use radcommon_h, only : pgasref,pfgasref,pgasmin,pgasmax
25      use radcommon_h, only : tgasref,tgasmin,tgasmax
26      use radcommon_h, only : gasv,gasi,FZEROI,FZEROV,gweight
27      use radcommon_h, only : wrefvar,WNOI,WNOV
28      use datafile_mod, only: datadir
29      use comcstfi_mod, only: mugaz
30      use gases_h, only: gnom, ngasmx, &
31                         igas_H2, igas_H2O, igas_He, igas_N2
32      use ioipsl_getin_p_mod, only: getin_p
33      use callkeys_mod, only: varactive,varfixed,graybody,callgasvis,&
34                continuum,H2Ocont_simple
35      implicit none
36
37!==================================================================
38
39      logical file_ok
40
41      integer n, nt, np, nh, ng, nw, m, i
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_NGAUSS
164      IF( .NOT. ALLOCATED( gweight ) ) ALLOCATE( GWEIGHT(L_NGAUSS) )
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      ! get array size, load the coefficients
198      open(111,file=TRIM(file_path),form='formatted')
199      read(111,*) L_NPREF
200      IF( .NOT. ALLOCATED( pgasref ) ) ALLOCATE( PGASREF(L_NPREF) )
201      read(111,*) pgasref
202      close(111)
203      L_PINT = (L_NPREF-1)*5+1
204      IF( .NOT. ALLOCATED( pfgasref ) ) ALLOCATE( PFGASREF(L_PINT) )
205
206      ! display the values
207      print*,'Correlated-k pressure grid (mBar):'
208      do n=1,L_NPREF
209         print*,n,'. 1 x 10^',pgasref(n),' mBar'
210      end do
211      print*,''
212
213      ! save the min / max matrix values
214      pgasmin = 10.0**pgasref(1)
215      pgasmax = 10.0**pgasref(L_NPREF)
216
217      ! interpolate to finer grid, adapted to uneven grids
218      do n=1,L_NPREF-1
219         do m=1,5
220            pfgasref((n-1)*5+m) = pgasref(n)+(m-1)*(pgasref(n+1) - pgasref(n))/5.
221         end do
222      end do
223      pfgasref(L_PINT) = pgasref(L_NPREF)
224
225!-----------------------------------------------------------------------
226! temperature
227
228      file_id='/corrk_data/' // TRIM(corrkdir) // '/T.dat'
229      file_path=TRIM(datadir)//TRIM(file_id)
230
231      ! check that the file exists
232      inquire(FILE=file_path,EXIST=file_ok)
233      if(.not.file_ok) then
234         write(*,*)'The file ',TRIM(file_path)
235         write(*,*)'was not found by sugas_corrk.F90, exiting.'
236         write(*,*)'Check that your path to datagcm:',trim(datadir)
237         write(*,*)' is correct. You can change it in callphys.def with:'
238         write(*,*)' datadir = /absolute/path/to/datagcm'
239         write(*,*)'Also check that the corrkdir you chose in callphys.def exists.'
240         call abort
241      endif
242
243      ! get array size, load the coefficients
244      open(111,file=TRIM(file_path),form='formatted')
245      read(111,*) L_NTREF
246      IF( .NOT. ALLOCATED( tgasref ) ) ALLOCATE( TGASREF(L_NTREF) )
247      read(111,*) tgasref
248      close(111)
249
250      ! display the values
251      print*,'Correlated-k temperature grid:'
252      do n=1,L_NTREF
253         print*,n,'.',tgasref(n),' K'
254      end do
255
256      ! save the min / max matrix values
257      tgasmin = tgasref(1)
258      tgasmax = tgasref(L_NTREF)
259
260      IF( .NOT. ALLOCATED( gasi8 ) ) ALLOCATE( gasi8(L_NTREF,L_NPREF,L_REFVAR,L_NSPECTI,L_NGAUSS) )
261      IF( .NOT. ALLOCATED( gasv8 ) ) ALLOCATE( gasv8(L_NTREF,L_NPREF,L_REFVAR,L_NSPECTV,L_NGAUSS) )
262!$OMP END MASTER
263!$OMP BARRIER
264
265!-----------------------------------------------------------------------
266! allocate the multidimensional arrays in radcommon_h
267      IF( .NOT. ALLOCATED( gasi ) ) ALLOCATE( gasi(L_NTREF,L_PINT,L_REFVAR,L_NSPECTI,L_NGAUSS) )
268      IF( .NOT. ALLOCATED( gasv ) ) ALLOCATE( gasv(L_NTREF,L_PINT,L_REFVAR,L_NSPECTV,L_NGAUSS) )
269
270      ! display the values
271      print*,''
272      print*,'Correlated-k matrix size:'
273      print*,'[',L_NTREF,',',L_NPREF,',',L_REFVAR,',',L_NGAUSS,']'
274
275!=======================================================================
276!     Get gaseous k-coefficients and interpolate onto finer pressure grid
277
278
279!        wavelength used to separate IR from VI in graybody. We will need that anyway
280         IR_VI_wnlimit=3000.
281         write(*,*)"graybody: Visible / Infrared separation set at",10000./IR_VI_wnlimit,"um"
282         
283         nVI_limit=0
284         Do nw=1,L_NSPECTV
285            if ((WNOV(nw).gt.IR_VI_wnlimit).and.(L_NSPECTV.gt.1)) then
286               nVI_limit=nw-1
287               exit
288            endif
289         End do
290         nIR_limit=L_NSPECTI
291         Do nw=1,L_NSPECTI
292            if ((WNOI(nw).gt.IR_VI_wnlimit).and.(L_NSPECTI.gt.1)) then
293               nIR_limit=nw-1
294               exit
295            endif
296         End do
297
298      if (graybody) then
299!        constant absorption coefficient in visible
300         write(*,*)"graybody: constant absorption coefficient in visible:"
301         kappa_VI=-100000.
302         call getin_p("kappa_VI",kappa_VI)
303         write(*,*)" kappa_VI = ",kappa_VI
304         kappa_VI=kappa_VI*1.e4* mugaz * 1.672621e-27    ! conversion from m^2/kg to cm^2/molecule         
305     
306!        constant absorption coefficient in IR
307         write(*,*)"graybody: constant absorption coefficient in InfraRed:"
308         kappa_IR=-100000.
309         call getin_p("kappa_IR",kappa_IR)
310         write(*,*)" kappa_IR = ",kappa_IR       
311         kappa_IR=kappa_IR*1.e4* mugaz * 1.672621e-27    ! conversion from m^2/kg to cm^2/molecule
312
313         write(*,*)"graybody: Visible / Infrared separation set at band: IR=",nIR_limit,", VI=",nVI_limit
314               
315      Else
316         kappa_VI=1.e-30     
317         kappa_IR=1.e-30       
318      End if
319
320!$OMP MASTER         
321!      print*,corrkdir(1:4)
322      ! VISIBLE
323      if (callgasvis) then
324         if ((corrkdir(1:4).eq.'null'))then   !(TRIM(corrkdir).eq.'null_LowTeffStar')) then
325            gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTV,1:L_NGAUSS)=0.0
326            print*,'using no corrk data'
327            print*,'Visible corrk gaseous absorption is set to zero if graybody=F'
328         else
329            file_id='/corrk_data/'//trim(adjustl(banddir))//'/corrk_gcm_VI.dat'
330            file_path=TRIM(datadir)//TRIM(file_id)
331           
332            ! check that the file exists
333            inquire(FILE=file_path,EXIST=file_ok)
334            if(.not.file_ok) then
335               write(*,*)'The file ',TRIM(file_path)
336               write(*,*)'was not found by sugas_corrk.F90.'
337               write(*,*)'Are you sure you have absorption data for these bands?'
338               call abort
339            endif
340         
341            open(111,file=TRIM(file_path),form='formatted')
342            read(111,*) gasv8
343            close(111)
344         end if
345
346         if(nVI_limit.eq.0) then
347            gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTV,1:L_NGAUSS)=   &
348                  gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTV,1:L_NGAUSS)+kappa_VI
349         else if (nVI_limit.eq.L_NSPECTV) then
350            gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTV,1:L_NGAUSS)=   &
351                  gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTV,1:L_NGAUSS)+kappa_IR
352         else
353            gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:nVI_limit,1:L_NGAUSS)=   &
354                  gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:nVI_limit,1:L_NGAUSS)+kappa_IR
355            gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,nVI_limit+1:L_NSPECTV,1:L_NGAUSS)=   &
356                  gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,nVI_limit+1:L_NSPECTV,1:L_NGAUSS)+kappa_VI
357         end if
358      else
359         print*,'Visible corrk gaseous absorption is set to zero.'
360         gasv8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTV,1:L_NGAUSS)=0.0
361      endif
362!$OMP END MASTER
363!$OMP BARRIER
364
365      ! INFRA-RED
366      if ((corrkdir(1:4).eq.'null'))then       !.or.(TRIM(corrkdir).eq.'null_LowTeffStar')) then
367         print*,'Infrared corrk gaseous absorption is set to zero if graybody=F'
368!$OMP MASTER         
369         gasi8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTI,1:L_NGAUSS)=0.0
370!$OMP END MASTER
371!$OMP BARRIER
372      else
373         file_id='/corrk_data/'//trim(adjustl(banddir))//'/corrk_gcm_IR.dat'
374         file_path=TRIM(datadir)//TRIM(file_id)
375     
376         ! check that the file exists
377         inquire(FILE=file_path,EXIST=file_ok)
378         if(.not.file_ok) then
379            write(*,*)'The file ',TRIM(file_path)
380            write(*,*)'was not found by sugas_corrk.F90.'
381            write(*,*)'Are you sure you have absorption data for these bands?'
382            call abort
383         endif
384         
385!$OMP MASTER                   
386         open(111,file=TRIM(file_path),form='formatted')
387         read(111,*) gasi8
388         close(111)
389!$OMP END MASTER
390!$OMP BARRIER
391     
392         ! 'fzero' is a currently unused feature that allows optimisation
393         ! of the radiative transfer by neglecting bands where absorption
394         ! is close to zero. As it could be useful in the future, this
395         ! section of the code has been kept commented and not erased.
396         ! RW 7/3/12.
397
398         do nw=1,L_NSPECTI
399            fzeroi(nw) = 0.d0
400!            do nt=1,L_NTREF
401!               do np=1,L_NPREF
402!                  do nh=1,L_REFVAR
403!                     do ng = 1,L_NGAUSS
404!                        if(gasi8(nt,np,nh,nw,ng).lt.1.0e-25)then
405!                           fzeroi(nw)=fzeroi(nw)+1.d0
406!                        endif
407!                     end do
408!                  end do
409!               end do
410!            end do
411!            fzeroi(nw)=fzeroi(nw)/dble(L_NTREF*L_NPREF*L_REFVAR*L_NGAUSS)
412         end do
413
414         do nw=1,L_NSPECTV
415            fzerov(nw) = 0.d0
416!            do nt=1,L_NTREF
417!               do np=1,L_NPREF
418!                  do nh=1,L_REFVAR
419!                     do ng = 1,L_NGAUSS
420!                        if(gasv8(nt,np,nh,nw,ng).lt.1.0e-25)then
421!                           fzerov(nw)=fzerov(nw)+1.d0
422!                        endif
423!                     end do
424!                  end do
425!               end do
426!            end do
427!            fzerov(nw)=fzerov(nw)/dble(L_NTREF*L_NPREF*L_REFVAR*L_NGAUSS)
428         end do
429
430      endif
431
432!$OMP MASTER                 
433      if(nIR_limit.eq.0) then
434         gasi8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTI,1:L_NGAUSS)=   &
435                  gasi8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTI,1:L_NGAUSS)+kappa_VI
436      else if (nIR_limit.eq.L_NSPECTI) then
437         gasi8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTI,1:L_NGAUSS)=   &
438                  gasi8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:L_NSPECTI,1:L_NGAUSS)+kappa_IR
439      else
440         gasi8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:nIR_limit,1:L_NGAUSS)=   &
441                  gasi8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,1:nIR_limit,1:L_NGAUSS)+kappa_IR
442         gasi8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,nIR_limit+1:L_NSPECTI,1:L_NGAUSS)=   &
443                  gasi8(1:L_NTREF,1:L_NPREF,1:L_REFVAR,nIR_limit+1:L_NSPECTI,1:L_NGAUSS)+kappa_VI
444      end if
445
446
447!     Take log10 of the values - this is what we will interpolate.
448!     Smallest value is 1.0E-200.
449
450      do nt=1,L_NTREF
451         do np=1,L_NPREF
452            do nh=1,L_REFVAR
453               do ng = 1,L_NGAUSS
454
455                  do nw=1,L_NSPECTV
456                     if(gasv8(nt,np,nh,nw,ng).gt.1.0d-200) then
457                        gasv8(nt,np,nh,nw,ng) = log10(gasv8(nt,np,nh,nw,ng))
458                     else
459                        gasv8(nt,np,nh,nw,ng) = -200.0
460                     end if
461                  end do
462
463                  do nw=1,L_NSPECTI
464                     if(gasi8(nt,np,nh,nw,ng).gt.1.0d-200) then
465                        gasi8(nt,np,nh,nw,ng) = log10(gasi8(nt,np,nh,nw,ng))
466                     else
467                        gasi8(nt,np,nh,nw,ng) = -200.0
468                     end if
469                  end do
470                 
471               end do
472            end do
473         end do
474      end do
475!$OMP END MASTER
476!$OMP BARRIER
477
478!     Interpolate the values:  first the longwave
479
480      do nt=1,L_NTREF
481         do nh=1,L_REFVAR
482            do nw=1,L_NSPECTI
483               do ng=1,L_NGAUSS
484
485!     First, the initial interval
486
487                  n = 1
488                  do m=1,5
489                     x     = pfgasref(m)
490                     xi(1) = pgasref(n)
491                     xi(2) = pgasref(n+1)
492                     xi(3) = pgasref(n+2)
493                     xi(4) = pgasref(n+3)
494                     yi(1) = gasi8(nt,n,nh,nw,ng)
495                     yi(2) = gasi8(nt,n+1,nh,nw,ng)
496                     yi(3) = gasi8(nt,n+2,nh,nw,ng)
497                     yi(4) = gasi8(nt,n+3,nh,nw,ng)
498                     call lagrange(x,xi,yi,ans)
499                     gasi(nt,m,nh,nw,ng) = 10.0**ans
500                  end do
501                 
502                  do n=2,L_NPREF-2
503                     do m=1,5
504                        i     = (n-1)*5+m
505                        x     = pfgasref(i)
506                        xi(1) = pgasref(n-1)
507                        xi(2) = pgasref(n)
508                        xi(3) = pgasref(n+1)
509                        xi(4) = pgasref(n+2)
510                        yi(1) = gasi8(nt,n-1,nh,nw,ng)
511                        yi(2) = gasi8(nt,n,nh,nw,ng)
512                        yi(3) = gasi8(nt,n+1,nh,nw,ng)
513                        yi(4) = gasi8(nt,n+2,nh,nw,ng)
514                        call lagrange(x,xi,yi,ans)
515                        gasi(nt,i,nh,nw,ng) = 10.0**ans
516                     end do
517                  end do
518
519!     Now, get the last interval
520
521                  n = L_NPREF-1                 
522                  do m=1,5
523                     i     = (n-1)*5+m
524                     x     = pfgasref(i)
525                     xi(1) = pgasref(n-2)
526                     xi(2) = pgasref(n-1)
527                     xi(3) = pgasref(n)
528                     xi(4) = pgasref(n+1)
529                     yi(1) = gasi8(nt,n-2,nh,nw,ng)
530                     yi(2) = gasi8(nt,n-1,nh,nw,ng)
531                     yi(3) = gasi8(nt,n,nh,nw,ng)
532                     yi(4) = gasi8(nt,n+1,nh,nw,ng)
533                     call lagrange(x,xi,yi,ans)
534                     gasi(nt,i,nh,nw,ng) = 10.0**ans
535                  end do 
536
537!     Fill the last pressure point
538
539                  gasi(nt,L_PINT,nh,nw,ng) = &
540                       10.0**gasi8(nt,L_NPREF,nh,nw,ng)
541
542               end do
543            end do
544         end do
545      end do
546
547!     Interpolate the values:  now the shortwave
548
549      do nt=1,L_NTREF
550         do nh=1,L_REFVAR
551            do nw=1,L_NSPECTV
552               do ng=1,L_NGAUSS
553
554!     First, the initial interval
555
556                  n = 1
557                  do m=1,5
558                     x     = pfgasref(m)
559                     xi(1) = pgasref(n)
560                     xi(2) = pgasref(n+1)
561                     xi(3) = pgasref(n+2)
562                     xi(4) = pgasref(n+3)
563                     yi(1) = gasv8(nt,n,nh,nw,ng)
564                     yi(2) = gasv8(nt,n+1,nh,nw,ng)
565                     yi(3) = gasv8(nt,n+2,nh,nw,ng)
566                     yi(4) = gasv8(nt,n+3,nh,nw,ng)
567                     call lagrange(x,xi,yi,ans)
568                     gasv(nt,m,nh,nw,ng) = 10.0**ans
569                  end do
570                 
571                  do n=2,L_NPREF-2
572                     do m=1,5
573                        i     = (n-1)*5+m
574                        x     = pfgasref(i)
575                        xi(1) = pgasref(n-1)
576                        xi(2) = pgasref(n)
577                        xi(3) = pgasref(n+1)
578                        xi(4) = pgasref(n+2)
579                        yi(1) = gasv8(nt,n-1,nh,nw,ng)
580                        yi(2) = gasv8(nt,n,nh,nw,ng)
581                        yi(3) = gasv8(nt,n+1,nh,nw,ng)
582                        yi(4) = gasv8(nt,n+2,nh,nw,ng)
583                        call lagrange(x,xi,yi,ans)
584                        gasv(nt,i,nh,nw,ng) = 10.0**ans
585                     end do
586                  end do
587
588!     Now, get the last interval
589
590                  n = L_NPREF-1
591                  do m=1,5
592                     i     = (n-1)*5+m
593                     x     = pfgasref(i)
594                     xi(1) = pgasref(n-2)
595                     xi(2) = pgasref(n-1)
596                     xi(3) = pgasref(n)
597                     xi(4) = pgasref(n+1)
598                     yi(1) = gasv8(nt,n-2,nh,nw,ng)
599                     yi(2) = gasv8(nt,n-1,nh,nw,ng)
600                     yi(3) = gasv8(nt,n,nh,nw,ng)
601                     yi(4) = gasv8(nt,n+1,nh,nw,ng)
602                     call lagrange(x,xi,yi,ans)
603                     gasv(nt,i,nh,nw,ng) = 10.0**ans
604                  end do 
605
606!     Fill the last pressure point
607
608                  gasv(nt,L_PINT,nh,nw,ng) = &
609                      10.0**gasv8(nt,L_NPREF,nh,nw,ng)
610                 
611               end do
612            end do
613         end do
614      end do
615
616
617!=======================================================================
618!     Initialise the continuum absorption data
619      if(continuum)then
620      do igas=1,ngasmx
621
622         if (igas .eq. igas_N2) then
623
624            dummy = -9999
625            call interpolateN2N2(100.D+0,250.D+0,17500.D+0,testcont,.true.,dummy)
626
627         elseif (igas .eq. igas_H2) then
628
629            ! first do self-induced absorption
630            dummy = -9999
631            call interpolateH2H2(500.D+0,250.D+0,17500.D+0,testcont,.true.,dummy)
632            ! then cross-interactions with other gases
633            do jgas=1,ngasmx
634               if (jgas .eq. igas_N2) then
635                  dummy = -9999
636                  call interpolateN2H2(592.D+0,278.15D+0,200000.D+0,10000.D+0,testcont,.true.,dummy)
637               elseif (jgas .eq. igas_He) then
638                  dummy = -9999
639                  call interpolateH2He(500.D+0,250.D+0,200000.D+0,10000.D+0,testcont,.true.,dummy)
640               endif
641            enddo
642
643         elseif (igas .eq. igas_H2O) then
644
645            ! H2O is special
646            if(H2Ocont_simple)then
647               call interpolateH2Ocont_PPC(990.D+0,296.D+0,683.2D+0*2,0.D+0,testcont,.true.)
648            else
649               dummy = -9999
650               call interpolateH2Ocont_CKD(990.D+0,296.D+0,683.2D+0*2,0.D+0,testcont,.true.,dummy)
651            endif
652
653         endif
654
655      enddo
656      endif
657
658      print*,'----------------------------------------------------'
659      print*,'And that`s all we have. It`s possible that other'
660      print*,'continuum absorption may be present, but if it is we'
661      print*,'don`t yet have data for it...'
662      print*,''
663
664!     Deallocate local arrays
665!$OMP BARRIER
666!$OMP MASTER
667      IF( ALLOCATED( gasi8 ) ) DEALLOCATE( gasi8 )
668      IF( ALLOCATED( gasv8 ) ) DEALLOCATE( gasv8 )
669      IF( ALLOCATED( pgasref ) ) DEALLOCATE( pgasref )
670      IF( ALLOCATED( gastype ) ) DEALLOCATE( gastype )
671!$OMP END MASTER
672!$OMP BARRIER
673
674      return
675    end subroutine sugas_corrk
Note: See TracBrowser for help on using the repository browser.