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

Last change on this file since 255 was 253, checked in by emillour, 14 years ago

Generic GCM

  • Massive update to version 0.7

EM+RW

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