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

Last change on this file since 648 was 596, checked in by jleconte, 13 years ago
  • Added double gray case (if graybody=true in callphys.def):
    • opacities are set to a constant value in sugas_corrk.
    • the values are kappa_IR m2/kg in the infrared (to be read in callphys.def)

kappa_VI m2/kg in the visible (to be read in callphys.def)

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