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

Last change on this file since 3693 was 3693, checked in by emillour, 3 months ago

Generic PCM:
Minor corrections:

  • about OpenMP in rad_common_h.F90 : (unclosed bracket in the THREADPRIVATE statement)
  • in interpolate_continuum.F90 : "filename" character string size must be specified as '*' (i.e. variable) and not a fixed number

Took the opportunity to also clean-up interpolate_continuum.F90 by
making it a module, adding some intent() statements, specifying
why saved variables are not threadprivate, get rid of useless "return"
statement at the end of routines, etc.
EM

File size: 34.2 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!     Modern traceur.def & corrk recombing scheme by J.Vatant d'Ollone (2020)
17!
18!     Summary
19!     -------
20!
21!==================================================================
22
23      use radinc_h, only: corrkdir, banddir, L_NSPECTI, L_NSPECTV, &
24                          L_NGAUSS, L_NPREF, L_NTREF, L_REFVAR, L_PINT
25      use radcommon_h, only : pgasref,pfgasref,pgasmin,pgasmax
26      use radcommon_h, only : tgasref,tgasmin,tgasmax
27      use radcommon_h, only : gasv,gasi,FZEROI,FZEROV,gweight
28      use radcommon_h, only : wrefvar,WNOI,WNOV
29      use datafile_mod, only: datadir
30      use comcstfi_mod, only: mugaz
31      use gases_h, only: gnom, ngasmx, &
32                         igas_H2, igas_H2O, igas_He, igas_N2, igas_CH4, &
33                         igas_CO2, igas_O2
34      use ioipsl_getin_p_mod, only: getin_p
35      use callkeys_mod, only: varactive,varfixed,graybody,callgasvis,&
36                continuum,generic_continuum_database
37      use tracer_h, only : nqtot, moderntracdef, is_recomb, noms
38      use recombin_corrk_mod, only: su_recombin,        &
39                corrk_recombin, use_premix, nrecomb_tot, rcb2tot_idx
40      use interpolate_continuum_mod, only: interpolate_continuum
41      implicit none
42
43!==================================================================
44
45      logical file_ok
46
47      integer n, nt, np, nh, ng, nw, m, i
48
49      character(len=200) :: file_id
50      character(len=500) :: file_path
51
52      ! ALLOCATABLE ARRAYS -- AS 12/2011
53      REAL*8, DIMENSION(:,:,:,:,:), ALLOCATABLE,SAVE :: gasi8, gasv8    !read by master
54      character*20,allocatable,DIMENSION(:),SAVE :: gastype ! for check with gnom, read by master
55
56      real*8 x, xi(4), yi(4), ans, p
57!     For gray case (JL12)
58      real kappa_IR, kappa_VI, IR_VI_wnlimit
59      integer nVI_limit,nIR_limit
60
61      integer ngas, igas, jgas
62
63      double precision testcont ! for continuum absorption initialisation
64
65      integer :: dummy
66
67      if (.not. moderntracdef) use_premix=.true. ! Added by JVO for compatibility with 'old' traceur.def
68     
69!$OMP MASTER
70      if (use_premix) then ! use_premix flag added by JVO, thus if pure recombining then premix is skipped
71
72!=======================================================================
73!     Load variable species data, exit if we have wrong database
74      file_id='/corrk_data/' // TRIM(corrkdir) // '/Q.dat'
75      file_path=TRIM(datadir)//TRIM(file_id)
76
77      ! check that the file exists
78      inquire(FILE=file_path,EXIST=file_ok)
79      if(.not.file_ok) then
80         write(*,*)'The file ',TRIM(file_path)
81         write(*,*)'was not found by sugas_corrk.F90, exiting.'
82         write(*,*)'Check that your path to datagcm:',trim(datadir)
83         write(*,*)' is correct. You can change it in callphys.def with:'
84         write(*,*)' datadir = /absolute/path/to/datagcm'
85         write(*,*)'Also check that the corrkdir you chose in callphys.def exists.'
86         call abort
87      endif
88
89      ! check that database matches varactive toggle
90      open(111,file=TRIM(file_path),form='formatted')
91      read(111,*) ngas
92
93      if(moderntracdef) then
94           ! JVO 20 - TODO : Sanity check with nspcrad !
95           print*, 'Warning : Sanity check on # of gases still not implemented !!'
96      else
97        if(ngas.ne.ngasmx)then
98           print*,'Number of gases in radiative transfer data (',ngas,') does not', &
99                  'match that in gases.def (',ngasmx,'), exiting.'
100           call abort
101        endif
102      endif
103
104      if(ngas.eq.1 .and. (varactive.or.varfixed))then
105         print*,'You have varactive/fixed=.true. but the database [',  &
106                     corrkdir(1:LEN_TRIM(corrkdir)),           &
107                '] has no variable species, exiting.'
108         call abort
109      elseif(ngas.gt.5 .or. ngas.lt.1)then
110         print*,ngas,' species in database [',               &
111                     corrkdir(1:LEN_TRIM(corrkdir)),           &
112                '], radiative code cannot handle this.'
113         call abort
114      endif
115
116      ! dynamically allocate gastype and read from Q.dat
117      IF ( .NOT. ALLOCATED( gastype ) ) ALLOCATE( gastype( ngas ) )
118
119      do igas=1,ngas
120         read(111,*) gastype(igas)
121         if(corrk_recombin) then
122            print*,'Premix gas ',igas,' is ',gastype(igas)
123         else
124           print*,'Gas ',igas,' is ',gastype(igas)
125         endif
126      enddo
127
128      ! get array size, load the coefficients
129      open(111,file=TRIM(file_path),form='formatted')
130      read(111,*) L_REFVAR
131      IF( .NOT. ALLOCATED( wrefvar ) ) ALLOCATE( WREFVAR(L_REFVAR) )
132      read(111,*) wrefvar
133      close(111)
134
135      if(L_REFVAR.gt.1 .and. (.not.varactive) .and. (.not.varfixed))then
136         print*,'You have varactive and varfixed=.false. and the database [', &
137                     corrkdir(1:LEN_TRIM(corrkdir)),           &
138                '] has a variable species.'
139         call abort
140      endif
141
142      if (moderntracdef) then
143          ! JVO 20 - TODO : Sanity check with nspcrad !
144          print*, 'Warning : Sanity check on name of gases still not implemented !!'
145      else
146        ! Check that gastype and gnom match
147        do igas=1,ngas
148           print*,'Gas ',igas,' is ',trim(gnom(igas))
149           if (trim(gnom(igas)).ne.trim(gastype(igas))) then
150              print*,'Name of a gas in radiative transfer data (',trim(gastype(igas)),') does not ', &
151                   'match that in gases.def (',trim(gnom(igas)),'), exiting. You should compare ', &
152                   'gases.def with Q.dat in your radiative transfer directory.'
153              call abort
154           endif
155        enddo
156        print*,'Confirmed gas match in radiative transfer and gases.def!'
157      endif
158
159      ! display the values
160      print*,'Variable gas volume mixing ratios:'
161      do n=1,L_REFVAR
162         !print*,n,'.',wrefvar(n),' kg/kg' ! pay attention!
163         print*,n,'.',wrefvar(n),' mol/mol'
164      end do
165      print*,''
166     
167      else
168        L_REFVAR = 1
169      endif ! use_premix
170
171!=======================================================================
172!     Set the weighting in g-space
173
174      file_id='/corrk_data/' // TRIM(corrkdir) // '/g.dat'
175      file_path=TRIM(datadir)//TRIM(file_id)
176
177      ! check that the file exists
178      inquire(FILE=file_path,EXIST=file_ok)
179      if(.not.file_ok) then
180         write(*,*)'The file ',TRIM(file_path)
181         write(*,*)'was not found by sugas_corrk.F90, exiting.'
182         write(*,*)'Check that your path to datagcm:',trim(datadir)
183         write(*,*)' is correct. You can change it in callphys.def with:'
184         write(*,*)' datadir = /absolute/path/to/datagcm'
185         write(*,*)'Also check that the corrkdir you chose in callphys.def exists.'
186         call abort
187      endif
188     
189      ! check the array size is correct, load the coefficients
190      open(111,file=TRIM(file_path),form='formatted')
191      read(111,*) L_NGAUSS
192      IF( .NOT. ALLOCATED( gweight ) ) ALLOCATE( GWEIGHT(L_NGAUSS) )
193      read(111,*) gweight
194      close(111)
195 
196      ! display the values
197      print*,'Correlated-k g-space grid:'
198      do n=1,L_NGAUSS
199         print*,n,'.',gweight(n)
200      end do
201      print*,''
202
203!=======================================================================
204!     Set the reference pressure and temperature arrays.  These are
205!     the pressures and temperatures at which we have k-coefficients.
206
207!-----------------------------------------------------------------------
208! pressure
209
210      file_id='/corrk_data/' // TRIM(corrkdir) // '/p.dat'
211      file_path=TRIM(datadir)//TRIM(file_id)
212
213      ! check that the file exists
214      inquire(FILE=file_path,EXIST=file_ok)
215      if(.not.file_ok) then
216         write(*,*)'The file ',TRIM(file_path)
217         write(*,*)'was not found by sugas_corrk.F90, exiting.'
218         write(*,*)'Check that your path to datagcm:',trim(datadir)
219         write(*,*)' is correct. You can change it in callphys.def with:'
220         write(*,*)' datadir = /absolute/path/to/datagcm'
221         write(*,*)'Also check that the corrkdir you chose in callphys.def exists.'
222         call abort
223      endif
224     
225      ! get array size, load the coefficients
226      open(111,file=TRIM(file_path),form='formatted')
227      read(111,*) L_NPREF
228      IF( .NOT. ALLOCATED( pgasref ) ) ALLOCATE( PGASREF(L_NPREF) )
229      read(111,*) pgasref
230      close(111)
231      L_PINT = (L_NPREF-1)*5+1
232      IF( .NOT. ALLOCATED( pfgasref ) ) ALLOCATE( PFGASREF(L_PINT) )
233
234      ! display the values
235      print*,'Correlated-k pressure grid (mBar):'
236      do n=1,L_NPREF
237         print*,n,'. 1 x 10^',pgasref(n),' mBar'
238      end do
239      print*,''
240
241      ! save the min / max matrix values
242      pgasmin = 10.0**pgasref(1)
243      pgasmax = 10.0**pgasref(L_NPREF)
244
245      ! interpolate to finer grid, adapted to uneven grids
246      do n=1,L_NPREF-1
247         do m=1,5
248            pfgasref((n-1)*5+m) = pgasref(n)+(m-1)*(pgasref(n+1) - pgasref(n))/5.
249         end do
250      end do
251      pfgasref(L_PINT) = pgasref(L_NPREF)
252
253!-----------------------------------------------------------------------
254! temperature
255
256      file_id='/corrk_data/' // TRIM(corrkdir) // '/T.dat'
257      file_path=TRIM(datadir)//TRIM(file_id)
258
259      ! check that the file exists
260      inquire(FILE=file_path,EXIST=file_ok)
261      if(.not.file_ok) then
262         write(*,*)'The file ',TRIM(file_path)
263         write(*,*)'was not found by sugas_corrk.F90, exiting.'
264         write(*,*)'Check that your path to datagcm:',trim(datadir)
265         write(*,*)' is correct. You can change it in callphys.def with:'
266         write(*,*)' datadir = /absolute/path/to/datagcm'
267         write(*,*)'Also check that the corrkdir you chose in callphys.def exists.'
268         call abort
269      endif
270
271      ! get array size, load the coefficients
272      open(111,file=TRIM(file_path),form='formatted')
273      read(111,*) L_NTREF
274      IF( .NOT. ALLOCATED( tgasref ) ) ALLOCATE( TGASREF(L_NTREF) )
275      read(111,*) tgasref
276      close(111)
277
278      ! display the values
279      print*,'Correlated-k temperature grid:'
280      do n=1,L_NTREF
281         print*,n,'.',tgasref(n),' K'
282      end do
283
284      ! save the min / max matrix values
285      tgasmin = tgasref(1)
286      tgasmax = tgasref(L_NTREF)
287
288      if(corrk_recombin) then ! even if no premix we keep a L_REFVAR=1 to store precombined at firstcall (see
289         IF(.NOT. ALLOCATED(gasi8)) ALLOCATE(gasi8(L_NTREF,L_NPREF,L_REFVAR+nrecomb_tot,L_NSPECTI,L_NGAUSS))
290         IF(.NOT. ALLOCATED(gasv8)) ALLOCATE(gasv8(L_NTREF,L_NPREF,L_REFVAR+nrecomb_tot,L_NSPECTV,L_NGAUSS))
291      else
292         IF(.NOT. ALLOCATED(gasi8)) ALLOCATE(gasi8(L_NTREF,L_NPREF,L_REFVAR,L_NSPECTI,L_NGAUSS))
293         IF(.NOT. ALLOCATED(gasv8)) ALLOCATE(gasv8(L_NTREF,L_NPREF,L_REFVAR,L_NSPECTV,L_NGAUSS))
294      endif
295!$OMP END MASTER
296!$OMP BARRIER
297
298! JVO 20 : In these gasi/gasi8 matrix  we stack in same dim. variable specie and species to recombine (to keep code small)
299
300!-----------------------------------------------------------------------
301! allocate the multidimensional arrays in radcommon_h
302     if(corrk_recombin) then
303        IF(.NOT. ALLOCATED(gasi)) ALLOCATE(gasi(L_NTREF,L_PINT,L_REFVAR+nrecomb_tot,L_NSPECTI,L_NGAUSS))
304        IF(.NOT. ALLOCATED(gasv)) ALLOCATE(gasv(L_NTREF,L_PINT,L_REFVAR+nrecomb_tot,L_NSPECTV,L_NGAUSS))
305        ! display the values
306        print*,''
307        print*,'Correlated-k matrix size:'
308        print*,'[',L_NTREF,',',L_NPREF,',',L_REFVAR+nrecomb_tot,' (',L_REFVAR,'+',nrecomb_tot,'),',L_NGAUSS,']'
309      else
310        IF(.NOT. ALLOCATED(gasi)) ALLOCATE(gasi(L_NTREF,L_PINT,L_REFVAR,L_NSPECTI,L_NGAUSS))
311        IF(.NOT. ALLOCATED(gasv)) ALLOCATE(gasv(L_NTREF,L_PINT,L_REFVAR,L_NSPECTV,L_NGAUSS))
312        ! display the values
313        print*,''
314        print*,'Correlated-k matrix size:'
315        print*,'[',L_NTREF,',',L_NPREF,',',L_REFVAR,',',L_NGAUSS,']'
316      endif
317
318!=======================================================================
319!     Get gaseous k-coefficients and interpolate onto finer pressure grid
320
321
322!        wavelength used to separate IR from VI in graybody. We will need that anyway
323         IR_VI_wnlimit=3000.
324         write(*,*)"graybody: Visible / Infrared separation set at",10000./IR_VI_wnlimit,"um"
325         
326         nVI_limit=0
327         Do nw=1,L_NSPECTV
328            if ((WNOV(nw).gt.IR_VI_wnlimit).and.(L_NSPECTV.gt.1)) then
329               nVI_limit=nw-1
330               exit
331            endif
332         End do
333         nIR_limit=L_NSPECTI
334         Do nw=1,L_NSPECTI
335            if ((WNOI(nw).gt.IR_VI_wnlimit).and.(L_NSPECTI.gt.1)) then
336               nIR_limit=nw-1
337               exit
338            endif
339         End do
340
341      if (graybody) then
342!        constant absorption coefficient in visible
343         write(*,*)"graybody: constant absorption coefficient in visible:"
344         kappa_VI=-100000.
345         call getin_p("kappa_VI",kappa_VI)
346         write(*,*)" kappa_VI = ",kappa_VI
347         kappa_VI=kappa_VI*1.e4* mugaz * 1.672621e-27    ! conversion from m^2/kg to cm^2/molecule         
348     
349!        constant absorption coefficient in IR
350         write(*,*)"graybody: constant absorption coefficient in InfraRed:"
351         kappa_IR=-100000.
352         call getin_p("kappa_IR",kappa_IR)
353         write(*,*)" kappa_IR = ",kappa_IR       
354         kappa_IR=kappa_IR*1.e4* mugaz * 1.672621e-27    ! conversion from m^2/kg to cm^2/molecule
355
356         write(*,*)"graybody: Visible / Infrared separation set at band: IR=",nIR_limit,", VI=",nVI_limit
357               
358      Else
359         kappa_VI=1.e-30     
360         kappa_IR=1.e-30       
361      End if
362
363!$OMP MASTER         
364
365      ! VISIBLE
366      if (callgasvis) then
367           
368        ! Looking for premixed corrk files corrk_gcm_VI.dat if needed
369        if (use_premix) then
370           ! print*,corrkdir(1:4)
371           if ((corrkdir(1:4).eq.'null'))then   !(TRIM(corrkdir).eq.'null_LowTeffStar')) then
372              gasv8(:,:,:,:,:)=0.0
373              print*,'using no corrk data'
374              print*,'Visible corrk gaseous absorption is set to zero if graybody=F'
375           else
376              file_id='/corrk_data/'//trim(adjustl(banddir))//'/corrk_gcm_VI.dat'
377              file_path=TRIM(datadir)//TRIM(file_id)
378             
379              ! check that the file exists
380              inquire(FILE=file_path,EXIST=file_ok)
381              if(.not.file_ok) then
382                 write(*,*)'The file ',TRIM(file_path)
383                 write(*,*)'was not found by sugas_corrk.F90.'
384                 write(*,*)'Are you sure you have absorption data for these bands?'
385                 call abort
386              endif
387           
388              open(111,file=TRIM(file_path),form='formatted')
389              read(111,*) gasv8(:,:,:L_REFVAR,:,:)
390              close(111)
391           end if
392        else
393           gasv8(:,:,1,:,:)=0.0 ! dummy but needs to be initialized
394        endif ! use_premix
395       
396        ! Looking for others files corrk_gcm_VI_XXX.dat if needed
397        if (corrk_recombin) then
398          do igas=1,nrecomb_tot
399
400            file_id='/corrk_data/'//trim(adjustl(banddir))//'/corrk_gcm_VI_'//trim(adjustl(noms(rcb2tot_idx(igas))))//'.dat'
401            file_path=TRIM(datadir)//TRIM(file_id)
402             
403            ! check that the file exists
404            inquire(FILE=file_path,EXIST=file_ok)
405            if(.not.file_ok) then
406               write(*,*)'The file ',TRIM(file_path)
407               write(*,*)'was not found by sugas_corrk.F90.'
408               write(*,*)'Are you sure you have absorption data for this specie at these bands?'
409               call abort
410            endif
411         
412            open(111,file=TRIM(file_path),form='formatted')
413            read(111,*) gasv8(:,:,L_REFVAR+igas,:,:)
414            close(111)
415          enddo                           
416        endif ! corrk_recombin
417
418        if(nVI_limit.eq.0) then
419         gasv8(:,:,:,:,:)= gasv8(:,:,:,:,:)+kappa_VI
420           else if (nVI_limit.eq.L_NSPECTV) then
421         gasv8(:,:,:,:,:)= gasv8(:,:,:,:,:)+kappa_IR
422      else
423         gasv8(:,:,:,1:nVI_limit,:)= gasv8(:,:,:,1:nVI_limit,:)+kappa_IR
424         gasv8(:,:,:,nVI_limit+1:L_NSPECTV,:)= gasv8(:,:,:,nVI_limit+1:L_NSPECTV,:)+kappa_VI
425      end if
426           
427         else
428            print*,'Visible corrk gaseous absorption is set to zero.'
429            gasv8(:,:,:,:,:)=0.0
430         endif ! callgasvis
431         
432!$OMP END MASTER
433!$OMP BARRIER
434
435      ! INFRA-RED
436     
437      ! Looking for premixed corrk files corrk_gcm_IR.dat if needed
438      if (use_premix) then
439        if ((corrkdir(1:4).eq.'null'))then       !.or.(TRIM(corrkdir).eq.'null_LowTeffStar')) then
440           print*,'Infrared corrk gaseous absorption is set to zero if graybody=F'
441!$OMP MASTER         
442           gasi8(:,:,:,:,:)=0.0
443!$OMP END MASTER
444!$OMP BARRIER
445        else
446           file_id='/corrk_data/'//trim(adjustl(banddir))//'/corrk_gcm_IR.dat'
447           file_path=TRIM(datadir)//TRIM(file_id)
448       
449           ! check that the file exists
450           inquire(FILE=file_path,EXIST=file_ok)
451           if(.not.file_ok) then
452              write(*,*)'The file ',TRIM(file_path)
453              write(*,*)'was not found by sugas_corrk.F90.'
454              write(*,*)'Are you sure you have absorption data for these bands?'
455              call abort
456           endif
457         
458!$OMP MASTER                   
459           open(111,file=TRIM(file_path),form='formatted')
460           read(111,*) gasi8(:,:,:L_REFVAR,:,:)
461           close(111)
462!$OMP END MASTER
463!$OMP BARRIER
464     
465           ! 'fzero' is a currently unused feature that allows optimisation
466           ! of the radiative transfer by neglecting bands where absorption
467           ! is close to zero. As it could be useful in the future, this
468           ! section of the code has been kept commented and not erased.
469           ! RW 7/3/12.
470
471           do nw=1,L_NSPECTI
472              fzeroi(nw) = 0.d0
473!              do nt=1,L_NTREF
474!                 do np=1,L_NPREF
475!                    do nh=1,L_REFVAR
476!                       do ng = 1,L_NGAUSS
477!                          if(gasi8(nt,np,nh,nw,ng).lt.1.0e-25)then
478!                             fzeroi(nw)=fzeroi(nw)+1.d0
479!                          endif
480!                       end do
481!                    end do
482!                 end do
483!              end do
484!              fzeroi(nw)=fzeroi(nw)/dble(L_NTREF*L_NPREF*L_REFVAR*L_NGAUSS)
485           end do
486
487           do nw=1,L_NSPECTV
488              fzerov(nw) = 0.d0
489!              do nt=1,L_NTREF
490!                 do np=1,L_NPREF
491!                    do nh=1,L_REFVAR
492!                       do ng = 1,L_NGAUSS
493!                          if(gasv8(nt,np,nh,nw,ng).lt.1.0e-25)then
494!                             fzerov(nw)=fzerov(nw)+1.d0
495!                          endif
496!                       end do
497!                    end do
498!                 end do
499!              end do
500!              fzerov(nw)=fzerov(nw)/dble(L_NTREF*L_NPREF*L_REFVAR*L_NGAUSS)
501           end do
502
503        endif ! if corrkdir=null
504
505      else
506         gasi8(:,:,1,:,:)=0.0 ! dummy but needs to be initialized
507      endif ! use_premix
508   
509      ! Looking for others files corrk_gcm_IR_XXX.dat if needed
510      if (corrk_recombin) then
511        do igas=1,nrecomb_tot
512
513          file_id='/corrk_data/'//trim(adjustl(banddir))//'/corrk_gcm_IR_'//trim(adjustl(noms(rcb2tot_idx(igas))))//'.dat'
514          file_path=TRIM(datadir)//TRIM(file_id)
515           
516          ! check that the file exists
517          inquire(FILE=file_path,EXIST=file_ok)
518          if(.not.file_ok) then
519             write(*,*)'The file ',TRIM(file_path)
520             write(*,*)'was not found by sugas_corrk.F90.'
521             write(*,*)'Are you sure you have absorption data for this specie at these bands?'
522             call abort
523          endif
524!$OMP MASTER           
525          open(111,file=TRIM(file_path),form='formatted')
526          read(111,*) gasi8(:,:,L_REFVAR+igas,:,:)
527          close(111)
528!$OMP END MASTER
529!$OMP BARRIER
530        enddo                           
531      endif ! corrk_recombin
532
533!$OMP MASTER                 
534      if(nIR_limit.eq.0) then
535         gasi8(:,:,:,:,:)= gasi8(:,:,:,:,:)+kappa_VI
536      else if (nIR_limit.eq.L_NSPECTI) then
537         gasi8(:,:,:,:,:)= gasi8(:,:,:,:,:)+kappa_IR
538      else
539         gasi8(:,:,:,1:nIR_limit,:)= gasi8(:,:,:,1:nIR_limit,:)+kappa_IR
540         gasi8(:,:,:,nIR_limit+1:L_NSPECTI,:)= gasi8(:,:,:,nIR_limit+1:L_NSPECTI,:)+kappa_VI
541      end if
542
543
544!     Take log10 of the values - this is what we will interpolate.
545!     Smallest value is 1.0E-200.
546
547      do nt=1,L_NTREF
548         do np=1,L_NPREF
549            do nh=1,L_REFVAR+nrecomb_tot
550               do ng = 1,L_NGAUSS
551
552                  do nw=1,L_NSPECTV
553                     if(gasv8(nt,np,nh,nw,ng).gt.1.0d-200) then
554                        gasv8(nt,np,nh,nw,ng) = log10(gasv8(nt,np,nh,nw,ng))
555                     else
556                        gasv8(nt,np,nh,nw,ng) = -200.0
557                     end if
558                  end do
559
560                  do nw=1,L_NSPECTI
561                     if(gasi8(nt,np,nh,nw,ng).gt.1.0d-200) then
562                        gasi8(nt,np,nh,nw,ng) = log10(gasi8(nt,np,nh,nw,ng))
563                     else
564                        gasi8(nt,np,nh,nw,ng) = -200.0
565                     end if
566                  end do
567                 
568               end do
569            end do
570         end do
571      end do
572!$OMP END MASTER
573!$OMP BARRIER
574
575!     Interpolate the values:  first the longwave
576
577      do nt=1,L_NTREF
578         do nh=1,L_REFVAR+nrecomb_tot
579            do nw=1,L_NSPECTI
580               do ng=1,L_NGAUSS
581
582!     First, the initial interval
583
584                  n = 1
585                  do m=1,5
586                     x     = pfgasref(m)
587                     xi(1) = pgasref(n)
588                     xi(2) = pgasref(n+1)
589                     xi(3) = pgasref(n+2)
590                     xi(4) = pgasref(n+3)
591                     yi(1) = gasi8(nt,n,nh,nw,ng)
592                     yi(2) = gasi8(nt,n+1,nh,nw,ng)
593                     yi(3) = gasi8(nt,n+2,nh,nw,ng)
594                     yi(4) = gasi8(nt,n+3,nh,nw,ng)
595                     call lagrange(x,xi,yi,ans)
596                     gasi(nt,m,nh,nw,ng) = 10.0**ans
597                  end do
598                 
599                  do n=2,L_NPREF-2
600                     do m=1,5
601                        i     = (n-1)*5+m
602                        x     = pfgasref(i)
603                        xi(1) = pgasref(n-1)
604                        xi(2) = pgasref(n)
605                        xi(3) = pgasref(n+1)
606                        xi(4) = pgasref(n+2)
607                        yi(1) = gasi8(nt,n-1,nh,nw,ng)
608                        yi(2) = gasi8(nt,n,nh,nw,ng)
609                        yi(3) = gasi8(nt,n+1,nh,nw,ng)
610                        yi(4) = gasi8(nt,n+2,nh,nw,ng)
611                        call lagrange(x,xi,yi,ans)
612                        gasi(nt,i,nh,nw,ng) = 10.0**ans
613                     end do
614                  end do
615
616!     Now, get the last interval
617
618                  n = L_NPREF-1                 
619                  do m=1,5
620                     i     = (n-1)*5+m
621                     x     = pfgasref(i)
622                     xi(1) = pgasref(n-2)
623                     xi(2) = pgasref(n-1)
624                     xi(3) = pgasref(n)
625                     xi(4) = pgasref(n+1)
626                     yi(1) = gasi8(nt,n-2,nh,nw,ng)
627                     yi(2) = gasi8(nt,n-1,nh,nw,ng)
628                     yi(3) = gasi8(nt,n,nh,nw,ng)
629                     yi(4) = gasi8(nt,n+1,nh,nw,ng)
630                     call lagrange(x,xi,yi,ans)
631                     gasi(nt,i,nh,nw,ng) = 10.0**ans
632                  end do 
633
634!     Fill the last pressure point
635
636                  gasi(nt,L_PINT,nh,nw,ng) = &
637                       10.0**gasi8(nt,L_NPREF,nh,nw,ng)
638
639               end do
640            end do
641         end do
642      end do
643
644!     Interpolate the values:  now the shortwave
645
646      do nt=1,L_NTREF
647         do nh=1,L_REFVAR+nrecomb_tot
648            do nw=1,L_NSPECTV
649               do ng=1,L_NGAUSS
650
651!     First, the initial interval
652
653                  n = 1
654                  do m=1,5
655                     x     = pfgasref(m)
656                     xi(1) = pgasref(n)
657                     xi(2) = pgasref(n+1)
658                     xi(3) = pgasref(n+2)
659                     xi(4) = pgasref(n+3)
660                     yi(1) = gasv8(nt,n,nh,nw,ng)
661                     yi(2) = gasv8(nt,n+1,nh,nw,ng)
662                     yi(3) = gasv8(nt,n+2,nh,nw,ng)
663                     yi(4) = gasv8(nt,n+3,nh,nw,ng)
664                     call lagrange(x,xi,yi,ans)
665                     gasv(nt,m,nh,nw,ng) = 10.0**ans
666                  end do
667                 
668                  do n=2,L_NPREF-2
669                     do m=1,5
670                        i     = (n-1)*5+m
671                        x     = pfgasref(i)
672                        xi(1) = pgasref(n-1)
673                        xi(2) = pgasref(n)
674                        xi(3) = pgasref(n+1)
675                        xi(4) = pgasref(n+2)
676                        yi(1) = gasv8(nt,n-1,nh,nw,ng)
677                        yi(2) = gasv8(nt,n,nh,nw,ng)
678                        yi(3) = gasv8(nt,n+1,nh,nw,ng)
679                        yi(4) = gasv8(nt,n+2,nh,nw,ng)
680                        call lagrange(x,xi,yi,ans)
681                        gasv(nt,i,nh,nw,ng) = 10.0**ans
682                     end do
683                  end do
684
685!     Now, get the last interval
686
687                  n = L_NPREF-1
688                  do m=1,5
689                     i     = (n-1)*5+m
690                     x     = pfgasref(i)
691                     xi(1) = pgasref(n-2)
692                     xi(2) = pgasref(n-1)
693                     xi(3) = pgasref(n)
694                     xi(4) = pgasref(n+1)
695                     yi(1) = gasv8(nt,n-2,nh,nw,ng)
696                     yi(2) = gasv8(nt,n-1,nh,nw,ng)
697                     yi(3) = gasv8(nt,n,nh,nw,ng)
698                     yi(4) = gasv8(nt,n+1,nh,nw,ng)
699                     call lagrange(x,xi,yi,ans)
700                     gasv(nt,i,nh,nw,ng) = 10.0**ans
701                  end do 
702
703!     Fill the last pressure point
704
705                  gasv(nt,L_PINT,nh,nw,ng) = &
706                      10.0**gasv8(nt,L_NPREF,nh,nw,ng)
707                 
708               end do
709            end do
710         end do
711      end do
712
713! Allocate and initialise arrays for corrk_recombin
714      if (corrk_recombin) call su_recombin
715
716!=======================================================================
717!     Initialise the continuum absorption data
718      if(continuum)then
719      if(generic_continuum_database)then
720       do igas=1,ngasmx ! we loop on all pairs of molecules that have data available
721       ! data can be downloaded from https://web.lmd.jussieu.fr/~lmdz/planets/LMDZ.GENERIC/datagcm/continuum_data/
722        if (igas .eq. igas_N2) then
723         file_id='/continuum_data/' // 'N2-N2_continuum_70-500K_2025.dat'
724         file_path=TRIM(datadir)//TRIM(file_id)
725         call interpolate_continuum(file_path,igas_N2,igas_N2,'IR',1,300.D+0,10000.D+0,20000.D+0,testcont,.true.)
726         do jgas=1,ngasmx
727          if (jgas .eq. igas_H2) then
728           file_id='/continuum_data/' // 'H2-N2_continuum_40-600K_2025.dat'
729           file_path=TRIM(datadir)//TRIM(file_id)
730           call interpolate_continuum(file_path,igas_N2,igas_H2,'IR',1,300.D+0,10000.D+0,20000.D+0,testcont,.true.)
731          elseif (jgas .eq. igas_O2) then
732           file_id='/continuum_data/' // 'O2-N2_continuum_100-500K_2025.dat'
733           file_path=TRIM(datadir)//TRIM(file_id)
734           call interpolate_continuum(file_path,igas_N2,igas_O2,'IR',1,300.D+0,10000.D+0,20000.D+0,testcont,.true.)
735          elseif (jgas .eq. igas_CH4) then
736           file_id='/continuum_data/' // 'CH4-N2_continuum_40-600K_2025.dat'
737           file_path=TRIM(datadir)//TRIM(file_id)
738           call interpolate_continuum(file_path,igas_N2,igas_CH4,'IR',1,300.D+0,10000.D+0,20000.D+0,testcont,.true.)
739          endif
740         enddo
741        elseif (igas .eq. igas_O2) then
742         file_id='/continuum_data/' // 'O2-O2_continuum_100-400K_2025.dat'
743         file_path=TRIM(datadir)//TRIM(file_id)
744         call interpolate_continuum(file_path,igas_O2,igas_O2,'IR',1,300.D+0,10000.D+0,20000.D+0,testcont,.true.)
745         do jgas=1,ngasmx
746          if (jgas .eq. igas_CO2) then
747           file_id='/continuum_data/' // 'CO2-O2_continuum_150-600K_2025.dat'
748           file_path=TRIM(datadir)//TRIM(file_id)
749           call interpolate_continuum(file_path,igas_CO2,igas_O2,'IR',1,300.D+0,10000.D+0,20000.D+0,testcont,.true.)
750          endif
751         enddo
752        elseif (igas .eq. igas_H2) then
753         file_id='/continuum_data/' // 'H2-H2_continuum_40-7000K_2025.dat'
754         file_path=TRIM(datadir)//TRIM(file_id)
755         call interpolate_continuum(file_path,igas_H2,igas_H2,'IR',1,300.D+0,10000.D+0,20000.D+0,testcont,.true.)
756         do jgas=1,ngasmx
757          if (jgas .eq. igas_CH4) then
758           file_id='/continuum_data/' // 'H2-CH4_continuum_40-600K_2025.dat'
759           file_path=TRIM(datadir)//TRIM(file_id)
760           call interpolate_continuum(file_path,igas_H2,igas_CH4,'IR',1,300.D+0,10000.D+0,20000.D+0,testcont,.true.)
761          elseif (jgas .eq. igas_He) then
762           file_id='/continuum_data/' // 'H2-He_continuum_40-5500K_2025.dat'
763           file_path=TRIM(datadir)//TRIM(file_id)
764           call interpolate_continuum(file_path,igas_H2,igas_He,'IR',1,300.D+0,10000.D+0,20000.D+0,testcont,.true.)
765          endif
766         enddo   
767        elseif (igas .eq. igas_CH4) then
768         file_id='/continuum_data/' // 'CH4-CH4_continuum_40-500K_2025.dat'
769         file_path=TRIM(datadir)//TRIM(file_id)
770         call interpolate_continuum(file_path,igas_CH4,igas_CH4,'IR',1,300.D+0,10000.D+0,20000.D+0,testcont,.true.)
771        elseif (igas .eq. igas_H2O) then
772         file_id='/continuum_data/' // 'H2O-H2O_continuum_100-2000K_2025.dat'
773         file_path=TRIM(datadir)//TRIM(file_id)
774         call interpolate_continuum(file_path,igas_H2O,igas_H2O,'IR',1,300.D+0,10000.D+0,20000.D+0,testcont,.true.)
775         do jgas=1,ngasmx
776          if (jgas .eq. igas_N2) then
777           file_id='/continuum_data/' // 'H2O-N2_continuum_100-2000K_2025.dat'
778           file_path=TRIM(datadir)//TRIM(file_id)
779           call interpolate_continuum(file_path,igas_H2O,igas_N2,'IR',1,300.D+0,10000.D+0,20000.D+0,testcont,.true.)
780          elseif (jgas .eq. igas_O2) then
781           file_id='/continuum_data/' // 'H2O-O2_continuum_100-2000K_2025.dat'
782           file_path=TRIM(datadir)//TRIM(file_id)
783           call interpolate_continuum(file_path,igas_H2O,igas_O2,'IR',1,300.D+0,10000.D+0,20000.D+0,testcont,.true.)
784          elseif (jgas .eq. igas_CO2) then
785           file_id='/continuum_data/' // 'H2O-CO2_continuum_100-800K_2025.dat'
786           file_path=TRIM(datadir)//TRIM(file_id)
787           call interpolate_continuum(file_path,igas_H2O,igas_CO2,'IR',1,300.D+0,10000.D+0,20000.D+0,testcont,.true.)
788          endif
789         enddo   
790        elseif (igas .eq. igas_CO2) then
791         file_id='/continuum_data/' // 'CO2-CO2_continuum_100-800K_2025.dat'
792         file_path=TRIM(datadir)//TRIM(file_id)
793         call interpolate_continuum(file_path,igas_CO2,igas_CO2,'IR',1,300.D+0,10000.D+0,20000.D+0,testcont,.true.)
794         do jgas=1,ngasmx
795          if (jgas .eq. igas_H2) then
796           file_id='/continuum_data/' // 'CO2-H2_continuum_100-800K_2025.dat'
797           file_path=TRIM(datadir)//TRIM(file_id)
798           call interpolate_continuum(file_path,igas_CO2,igas_H2,'IR',1,300.D+0,10000.D+0,20000.D+0,testcont,.true.)
799          elseif (jgas .eq. igas_CH4) then
800           file_id='/continuum_data/' // 'CO2-CH4_continuum_100-800K_2025.dat'
801           file_path=TRIM(datadir)//TRIM(file_id)
802           call interpolate_continuum(file_path,igas_CO2,igas_CH4,'IR',1,300.D+0,10000.D+0,20000.D+0,testcont,.true.)
803          endif
804         enddo
805        endif
806       enddo ! igas=1,ngasmx
807       
808      else ! generic_continuum_database tag
809        do igas=1,ngasmx
810         if (igas .eq. igas_N2) then
811
812            dummy = -9999
813            call interpolateN2N2(100.D+0,250.D+0,17500.D+0,testcont,.true.,dummy)
814
815         elseif (igas .eq. igas_H2) then
816
817            ! first do self-induced absorption
818            dummy = -9999
819            call interpolateH2H2(500.D+0,250.D+0,17500.D+0,testcont,.true.,dummy)
820            ! then cross-interactions with other gases
821            do jgas=1,ngasmx
822               if (jgas .eq. igas_N2) then
823                  dummy = -9999
824                  call interpolateN2H2(592.D+0,278.15D+0,200000.D+0,10000.D+0,testcont,.true.,dummy)
825               elseif (jgas .eq. igas_CO2) then
826                  dummy = -9999
827                  call interpolateCO2H2(592.D+0,278.15D+0,200000.D+0,10000.D+0,testcont,.true.,dummy)
828               elseif (jgas .eq. igas_He) then
829                  dummy = -9999
830                  call interpolateH2He(500.D+0,250.D+0,200000.D+0,10000.D+0,testcont,.true.,dummy)
831               endif
832            enddo
833               
834         
835         elseif (igas .eq. igas_CH4) then
836         
837            ! first do self-induced absorption
838            dummy = -9999
839            call interpolateCH4CH4(600.0,66.7,400.0,testcont,.true.,dummy)
840            ! then cross-interactions with other gases
841            do jgas=1,ngasmx
842               if (jgas .eq. igas_H2) then
843                  dummy = -9999
844                  call interpolateH2CH4(500.D+0,250.D+0,200000.D+0,10000.D+0,testcont,.true.,dummy)
845               elseif (jgas .eq. igas_CO2) then
846                  dummy = -9999
847                  call interpolateCO2CH4(592.D+0,278.15D+0,200000.D+0,10000.D+0,testcont,.true.,dummy)
848               elseif (jgas .eq. igas_He) then
849                  dummy = -9999
850                  call interpolateHeCH4(500.D+0,250.D+0,200000.D+0,10000.D+0,testcont,.true.,dummy)
851               endif
852            enddo
853           
854
855         elseif (igas .eq. igas_H2O) then
856
857            ! Compute self and foreign (with air) continuum of H2O
858            dummy = -9999
859            call interpolateH2O_self_foreign(990.D+0,296.D+0,683.2D+0*2,0.D+0,testcont,.true.,dummy) 
860
861         endif
862
863      enddo ! igas=1,ngasmx
864      endif ! generic_continuum_database tag
865      endif ! continuum flag
866
867      print*,'----------------------------------------------------'
868      print*,'And that`s all we have. It`s possible that other'
869      print*,'continuum absorption may be present, but if it is we'
870      print*,'don`t yet have data for it...'
871      print*,''
872
873!     Deallocate local arrays
874!$OMP BARRIER
875!$OMP MASTER
876      IF( ALLOCATED( gasi8 ) ) DEALLOCATE( gasi8 )
877      IF( ALLOCATED( gasv8 ) ) DEALLOCATE( gasv8 )
878      IF( ALLOCATED( pgasref ) ) DEALLOCATE( pgasref )
879      IF( ALLOCATED( gastype ) ) DEALLOCATE( gastype )
880!$OMP END MASTER
881!$OMP BARRIER
882
883    end subroutine sugas_corrk
Note: See TracBrowser for help on using the repository browser.