source: trunk/LMDZ.GENERIC/libf/phystd/callcorrk.F90 @ 704

Last change on this file since 704 was 650, checked in by jleconte, 13 years ago
  • Corrected the temperature used to differentiate sublimation and evaporation in watersat_grad
  • Minor name changes in watercommon
  • Better physical parametrization of the effective radius of liquid and icy water cloud particles in callcorrk

(for radfixed=true)

  • Added consistency check in inifis
  • Moved 1d water initialization from physiqu to rcm1d
  • Property svn:executable set to *
File size: 32.0 KB
Line 
1      subroutine callcorrk(ngrid,nlayer,pq,nq,qsurf,           &
2          albedo,emis,mu0,pplev,pplay,pt,                      &
3          tsurf,fract,dist_star,aerosol,muvar,           &
4          dtlw,dtsw,fluxsurf_lw,                               &
5          fluxsurf_sw,fluxtop_lw,fluxabs_sw,fluxtop_dn,        &
6          OLR_nu,OSR_nu,                                       &
7          reffrad,tau_col,cloudfrac,totcloudfrac,              &
8          clearsky,firstcall,lastcall)
9
10
11      use radinc_h
12      use radcommon_h
13      use watercommon_h
14      use datafile_mod, only: datadir
15      use ioipsl_getincom
16      use gases_h
17
18      implicit none
19
20!==================================================================
21!
22!     Purpose
23!     -------
24!     Solve the radiative transfer using the correlated-k method for
25!     the gaseous absorption and the Toon et al. (1989) method for
26!     scatttering due to aerosols.
27!
28!     Authors
29!     -------
30!     Emmanuel 01/2001, Forget 09/2001
31!     Robin Wordsworth (2009)
32!
33!==================================================================
34
35#include "dimphys.h"
36#include "comcstfi.h"
37#include "callkeys.h"
38#include "tracer.h"
39
40!-----------------------------------------------------------------------
41!     Declaration of the arguments (INPUT - OUTPUT) on the LMD GCM grid
42!     Layer #1 is the layer near the ground.
43!     Layer #nlayermx is the layer at the top.
44
45!     INPUT
46      INTEGER icount
47      INTEGER ngrid,nlayer
48      REAL aerosol(ngrid,nlayermx,naerkind) ! aerosol tau (kg/kg)
49      REAL albedo(ngrid)                    ! SW albedo
50      REAL emis(ngrid)                      ! LW emissivity
51      REAL pplay(ngrid,nlayermx)            ! pres. level in GCM mid of layer
52      REAL pplev(ngrid,nlayermx+1)          ! pres. level at GCM layer boundaries
53
54      REAL pt(ngrid,nlayermx)               ! air temperature (K)
55      REAL tsurf(ngrid)                     ! surface temperature (K)
56      REAL dist_star,mu0(ngrid)             ! distance star-planet (AU)
57      REAL fract(ngrid)                     ! fraction of day
58
59!     Globally varying aerosol optical properties on GCM grid
60!     Not needed everywhere so not in radcommon_h
61      REAL :: QVISsQREF3d(ngridmx,nlayermx,L_NSPECTV,naerkind)
62      REAL :: omegaVIS3d(ngridmx,nlayermx,L_NSPECTV,naerkind)
63      REAL :: gVIS3d(ngridmx,nlayermx,L_NSPECTV,naerkind)
64
65      REAL :: QIRsQREF3d(ngridmx,nlayermx,L_NSPECTI,naerkind)
66      REAL :: omegaIR3d(ngridmx,nlayermx,L_NSPECTI,naerkind)
67      REAL :: gIR3d(ngridmx,nlayermx,L_NSPECTI,naerkind)
68
69      REAL :: QREFvis3d(ngridmx,nlayermx,naerkind)
70      REAL :: QREFir3d(ngridmx,nlayermx,naerkind)
71
72!      REAL :: omegaREFvis3d(ngridmx,nlayermx,naerkind)
73!      REAL :: omegaREFir3d(ngridmx,nlayermx,naerkind) ! not sure of the point of these...
74
75      REAL reffrad(ngrid,nlayer,naerkind)
76      REAL nueffrad(ngrid,nlayer,naerkind)
77
78!     OUTPUT
79      REAL dtsw(ngridmx,nlayermx) ! heating rate (K/s) due to SW
80      REAL dtlw(ngridmx,nlayermx) ! heating rate (K/s) due to LW
81      REAL fluxsurf_lw(ngridmx)   ! incident LW flux to surf (W/m2)
82      REAL fluxtop_lw(ngridmx)    ! outgoing LW flux to space (W/m2)
83      REAL fluxsurf_sw(ngridmx)   ! incident SW flux to surf (W/m2)
84      REAL fluxabs_sw(ngridmx)    ! SW flux absorbed by planet (W/m2)
85      REAL fluxtop_dn(ngridmx)    ! incident top of atmosphere SW flux (W/m2)
86      REAL OLR_nu(ngrid,L_NSPECTI)! Outgoing LW radition in each band (Normalized to the band width (W/m2/cm-1)
87      REAL OSR_nu(ngrid,L_NSPECTV)! Outgoing SW radition in each band (Normalized to the band width (W/m2/cm-1)
88!-----------------------------------------------------------------------
89!     Declaration of the variables required by correlated-k subroutines
90!     Numbered from top to bottom unlike in the GCM!
91
92      REAL*8 tmid(L_LEVELS),pmid(L_LEVELS)
93      REAL*8 tlevrad(L_LEVELS),plevrad(L_LEVELS)
94
95!     Optical values for the optci/cv subroutines
96      REAL*8 stel(L_NSPECTV),stel_fract(L_NSPECTV)
97      REAL*8 dtaui(L_NLAYRAD,L_NSPECTI,L_NGAUSS)
98      REAL*8 dtauv(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
99      REAL*8 cosbv(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
100      REAL*8 cosbi(L_NLAYRAD,L_NSPECTI,L_NGAUSS)
101      REAL*8 wbari(L_NLAYRAD,L_NSPECTI,L_NGAUSS)
102      REAL*8 wbarv(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
103      REAL*8 tauv(L_NLEVRAD,L_NSPECTV,L_NGAUSS)
104      REAL*8 taucumv(L_LEVELS,L_NSPECTV,L_NGAUSS)
105      REAL*8 taucumi(L_LEVELS,L_NSPECTI,L_NGAUSS)
106
107      REAL*8 tauaero(L_LEVELS+1,naerkind)
108      REAL*8 nfluxtopv,nfluxtopi,nfluxtop
109      real*8 nfluxoutv_nu(L_NSPECTV) ! outgoing band-resolved VI flux at TOA (W/m2)
110      real*8 nfluxtopi_nu(L_NSPECTI) ! net band-resolved IR flux at TOA (W/m2)
111      real*8 fluxupi_nu(L_NLAYRAD,L_NSPECTI) ! for 1D diagnostic
112      REAL*8 fmneti(L_NLAYRAD),fmnetv(L_NLAYRAD)
113      REAL*8 fluxupv(L_NLAYRAD),fluxupi(L_NLAYRAD)
114      REAL*8 fluxdnv(L_NLAYRAD),fluxdni(L_NLAYRAD)
115      REAL*8 albi,albv,acosz
116
117      INTEGER ig,l,k,nw,iaer,irad
118
119      real fluxtoplanet
120      real szangle
121      logical global1d
122      save szangle,global1d
123      real*8 taugsurf(L_NSPECTV,L_NGAUSS-1)
124      real*8 taugsurfi(L_NSPECTI,L_NGAUSS-1)
125
126      real*8 qvar(L_LEVELS)          ! mixing ratio of variable component (mol/mol)
127      REAL pq(ngridmx,nlayer,nq)
128      REAL qsurf(ngridmx,nqmx)       ! tracer on surface (e.g. kg.m-2)
129      integer nq
130
131!     Local aerosol optical properties for each column on RADIATIVE grid
132      real*8  QXVAER(L_LEVELS+1,L_NSPECTV,naerkind)
133      real*8  QSVAER(L_LEVELS+1,L_NSPECTV,naerkind)
134      real*8  GVAER(L_LEVELS+1,L_NSPECTV,naerkind)
135      real*8  QXIAER(L_LEVELS+1,L_NSPECTI,naerkind)
136      real*8  QSIAER(L_LEVELS+1,L_NSPECTI,naerkind)
137      real*8  GIAER(L_LEVELS+1,L_NSPECTI,naerkind)
138
139      save qxvaer, qsvaer, gvaer
140      save qxiaer, qsiaer, giaer
141      save QREFvis3d, QREFir3d
142
143      REAL tau_col(ngrid) ! diagnostic from aeropacity
144
145!     Misc.
146      logical firstcall, lastcall, nantest
147      real*8  tempv(L_NSPECTV)
148      real*8  tempi(L_NSPECTI)
149      real*8  temp,temp1,temp2,pweight
150      character(len=10) :: tmp1
151      character(len=10) :: tmp2
152
153!     for fixed water vapour profiles
154      integer i_var
155      real RH
156      real*8 pq_temp(nlayer)
157      real ptemp, Ttemp, qsat
158
159!      real(KIND=r8) :: pq_temp(nlayer) ! better F90 way.. DOESNT PORT TO F77!!!
160
161      !real ptime, pday
162      logical OLRz
163      real*8 NFLUXGNDV_nu(L_NSPECTV)
164
165!     for H2O cloud fraction in aeropacity
166      real cloudfrac(ngridmx,nlayermx)
167      real totcloudfrac(ngridmx)
168      logical clearsky
169
170      ! for weird cloud test
171      real pqtest(ngridmx,nlayer,nq)
172
173!     are particle radii fixed?
174      logical, save ::  radfixed
175      real maxrad, minrad
176
177!     Local water cloud optical properties
178      real, parameter ::  rad_froid=35.0e-6
179      real, parameter ::  rad_chaud=10.0e-6
180      real, parameter ::  coef_chaud=0.13
181      real, parameter ::  coef_froid=0.09
182      real zfice
183     
184      real CBRT
185      external CBRT
186
187!     included by RW for runaway greenhouse 1D study
188      real muvar(ngridmx,nlayermx+1)
189      real vtmp(nlayermx)
190      REAL*8 muvarrad(L_LEVELS)
191
192!=======================================================================
193!     Initialization on first call
194
195      qxvaer(:,:,:)=0.0
196      qsvaer(:,:,:)=0.0
197      gvaer(:,:,:) =0.0
198
199      qxiaer(:,:,:)=0.0
200      qsiaer(:,:,:)=0.0
201      giaer(:,:,:) =0.0
202      radfixed=.false.
203
204
205      if(firstcall) then
206         if(kastprof)then
207            radfixed=.true.
208         endif 
209
210         call system('rm -f surf_vals_long.out')
211
212!--------------------------------------------------
213!     Effective radius and variance of the aerosols
214
215         do iaer=1,naerkind
216!     these values will change once the microphysics gets to work
217!     UNLESS tracer=.false., in which case we should be working with
218!     a fixed aerosol layer, and be able to define reffrad in a
219!     .def file. To be improved!
220
221            if(iaer.eq.1)then ! CO2 ice
222               do l=1,nlayer
223                  do ig=1,ngrid
224                     reffrad(ig,l,iaer)  = 1.e-4
225                     nueffrad(ig,l,iaer) = 0.1
226                  enddo
227               enddo
228            endif
229
230            if(iaer.eq.2)then ! H2O ice
231               do l=1,nlayer
232                  do ig=1,ngrid
233                     reffrad(ig,l,iaer)  = 1.e-5
234                     nueffrad(ig,l,iaer) = 0.1
235                  enddo
236               enddo
237            endif
238
239            if(iaer.eq.3)then ! dust
240               do l=1,nlayer
241                  do ig=1,ngrid
242                     reffrad(ig,l,iaer)  = 1.e-5
243                     nueffrad(ig,l,iaer) = 0.1
244                  enddo
245               enddo
246            endif
247 
248            if(iaer.gt.3)then
249               print*,'Error in callcorrk, naerkind is too high.'
250               print*,'The code still needs generalisation to arbitrary'
251               print*,'aerosol kinds and number.'
252               call abort
253            endif
254
255         enddo
256
257         print*, "callcorrk: Correlated-k data base folder:",trim(datadir)
258         call getin("corrkdir",corrkdir)
259         print*, "corrkdir = ",corrkdir
260         write( tmp1, '(i3)' ) L_NSPECTI
261         write( tmp2, '(i3)' ) L_NSPECTV
262         banddir=trim(adjustl(tmp1))//'x'//trim(adjustl(tmp2))
263         banddir=trim(adjustl(corrkdir))//'/'//trim(adjustl(banddir))
264
265         call sugas_corrk       ! set up gaseous absorption properties
266         call setspi            ! basic infrared properties
267         call setspv            ! basic visible properties
268         call suaer_corrk       ! set up aerosol optical properties
269
270         Cmk= 0.01 * 1.0 / (g * mugaz * 1.672621e-27) ! q_main=1.0 assumed
271
272         if((igcm_h2o_vap.eq.0) .and. varactive)then
273            print*,'varactive in callcorrk but no h2o_vap tracer.'
274            stop
275         endif
276
277         OLR_nu=0.
278         OSR_nu=0.
279
280         if (ngridmx.eq.1) then
281           PRINT*, 'Simulate global averaged conditions ?'
282           global1d = .false. ! default value
283           call getin("global1d",global1d)
284           write(*,*) "global1d = ",global1d
285           ! Test of incompatibility:
286           ! if global1d is true, there should not be any diurnal cycle
287           if (global1d.and.diurnal) then
288            print*,'if global1d is true, diurnal must be set to false'
289            stop
290           endif
291
292           if (global1d) then
293             PRINT *,'Solar Zenith angle (deg.) ?'
294             PRINT *,'(assumed for averaged solar flux S/4)'
295             szangle=60.0  ! default value
296             call getin("szangle",szangle)
297             write(*,*) "szangle = ",szangle
298           endif
299         endif
300
301         firstcall=.false.   
302
303      end if
304
305!=======================================================================
306!     Initialization on every call   
307
308      do l=1,nlayer
309         do ig=1,ngrid
310            do iaer=1,naerkind
311               nueffrad(ig,l,iaer) = 0.1 ! stays at 0.1
312            enddo
313         enddo
314      enddo
315
316
317
318      if(radfixed)then
319         do l=1,nlayer
320            do ig=1,ngrid
321               reffrad(ig,l,1)    = 5.e-5 ! CO2 ice
322            enddo
323         enddo
324         print*,'CO2 ice particle size = ',reffrad(1,1,1)/1.e-6,' um'
325         if(naerkind.ge.2)then
326            do l=1,nlayer
327               do ig=1,ngrid
328                  !reffrad(ig,l,2) = 2.e-5 ! H2O ice
329                   reffrad(ig,l,2) = rad_chaud ! H2O ice
330
331                   zfice = 1.0 - (pt(ig,l)-T_h2O_ice_clouds) / (T_h2O_ice_liq-T_h2O_ice_clouds)
332                   zfice = MIN(MAX(zfice,0.0),1.0)
333                   reffrad(ig,l,2)= rad_chaud * (1.-zfice) + rad_froid * zfice
334                   nueffrad(ig,l,2) = coef_chaud * (1.-zfice) + coef_froid * zfice
335               enddo
336            enddo
337            print*,'H2O ice particle size = ',reffrad(1,1,2)/1.e-6,' um'
338         endif
339         if(naerkind.eq.3)then
340            do l=1,nlayer
341               do ig=1,ngrid
342                  reffrad(ig,l,naerkind) = 2.e-6 ! dust
343               enddo
344            enddo
345            print*,'Dust particle size = ',reffrad(1,1,naerkind)/1.e-6,' um'
346         endif
347         if(naerkind.gt.3)then
348            print*,'Code not general enough to deal with naerkind > 3 yet.'
349            call abort
350         endif
351
352      else
353
354         maxrad=0.0
355         !averad=0.0
356         minrad=1.0
357         do l=1,nlayer
358
359            !masse = (pplev(ig,l) - pplev(ig,l+1))/g
360
361            do ig=1,ngrid
362               !if(tracer)then
363               if(tracer.and.igcm_co2_ice.gt.0)then
364
365                  if(igcm_co2_ice.lt.1)then
366                     print*,'Tracers but no CO2 ice still seems to be a problem...'
367                     print*,'Aborting in callcorrk.'
368                     stop
369                  endif
370
371                  reffrad(ig,l,1) = CBRT( 3*pq(ig,l,igcm_co2_ice)/ &
372                       (4*Nmix_co2*pi*rho_co2) )
373               endif
374               reffrad(ig,l,1) = max(reffrad(ig,l,1),1.e-6)
375               reffrad(ig,l,1) = min(reffrad(ig,l,1),500.e-6)
376
377               !averad = averad + reffrad(ig,l,1)*area(ig)
378               maxrad = max(reffrad(ig,l,1),maxrad)
379               minrad = min(reffrad(ig,l,1),minrad)
380            enddo
381         enddo
382         if(igcm_co2_ice.gt.0)then
383            print*,'Max. CO2 ice particle size = ',maxrad/1.e-6,' um'
384            print*,'Min. CO2 ice particle size = ',minrad/1.e-6,' um'
385         endif
386
387         if((naerkind.ge.2).and.water)then
388            maxrad=0.0
389            minrad=1.0
390            do l=1,nlayer
391               do ig=1,ngrid
392                  reffrad(ig,l,2) = CBRT( 3*pq(ig,l,igcm_h2o_ice)/ &
393                       (4*Nmix_h2o*pi*rho_ice) )
394                  reffrad(ig,l,2) = max(reffrad(ig,l,2),1.e-6)
395                  reffrad(ig,l,2) = min(reffrad(ig,l,2),100.e-6)
396
397                  maxrad = max(reffrad(ig,l,2),maxrad)
398                  minrad = min(reffrad(ig,l,2),minrad)
399               enddo
400            enddo
401            print*,'Max. H2O ice particle size = ',maxrad/1.e-6,' um'
402            print*,'Min. H2O ice particle size = ',minrad/1.e-6,' um'
403         endif
404
405         if(naerkind.eq.3)then
406            do l=1,nlayer
407               do ig=1,ngrid
408                  reffrad(ig,l,naerkind) = 2.e-6 ! dust
409               enddo
410            enddo
411         endif
412
413      endif
414
415
416!     how much light we get
417      fluxtoplanet=0
418      do nw=1,L_NSPECTV
419         stel(nw)=stellarf(nw)/(dist_star**2)
420         fluxtoplanet=fluxtoplanet + stel(nw)
421      end do
422
423      call aeroptproperties(ngrid,nlayer,reffrad,nueffrad,         &
424           QVISsQREF3d,omegaVIS3d,gVIS3d,                          &
425           QIRsQREF3d,omegaIR3d,gIR3d,                             &
426           QREFvis3d,QREFir3d)                                     ! get 3D aerosol optical properties
427
428      call aeropacity(ngrid,nlayer,nq,pplay,pplev,pq,aerosol,      &
429           reffrad,QREFvis3d,QREFir3d,                             &
430           tau_col,cloudfrac,totcloudfrac,clearsky)                ! get aerosol optical depths
431
432
433!-----------------------------------------------------------------------
434!     Starting Big Loop over every GCM column
435      do ig=1,ngridmx
436
437!=======================================================================
438!     Transformation of the GCM variables
439
440!-----------------------------------------------------------------------
441!     Aerosol optical properties Qext, Qscat and g
442!     The transformation in the vertical is the same as for temperature
443           
444!     shortwave
445            do iaer=1,naerkind
446               DO nw=1,L_NSPECTV
447                  do l=1,nlayermx
448
449                     temp1=QVISsQREF3d(ig,nlayermx+1-l,nw,iaer)         &
450                         *QREFvis3d(ig,nlayermx+1-l,iaer)
451
452                     temp2=QVISsQREF3d(ig,max(nlayermx-l,1),nw,iaer)    &
453                         *QREFvis3d(ig,max(nlayermx-l,1),iaer)
454
455                     qxvaer(2*l,nw,iaer)  = temp1
456                     qxvaer(2*l+1,nw,iaer)=(temp1+temp2)/2
457
458                     temp1=temp1*omegavis3d(ig,nlayermx+1-l,nw,iaer)
459                     temp2=temp2*omegavis3d(ig,max(nlayermx-l,1),nw,iaer)
460
461                     qsvaer(2*l,nw,iaer)  = temp1
462                     qsvaer(2*l+1,nw,iaer)=(temp1+temp2)/2
463
464                     temp1=gvis3d(ig,nlayermx+1-l,nw,iaer)
465                     temp2=gvis3d(ig,max(nlayermx-l,1),nw,iaer)
466
467                     gvaer(2*l,nw,iaer)  = temp1
468                     gvaer(2*l+1,nw,iaer)=(temp1+temp2)/2
469
470                  end do
471
472                  qxvaer(1,nw,iaer)=qxvaer(2,nw,iaer)
473                  qxvaer(2*nlayermx+1,nw,iaer)=0.
474
475                  qsvaer(1,nw,iaer)=qsvaer(2,nw,iaer)
476                  qsvaer(2*nlayermx+1,nw,iaer)=0.
477
478                  gvaer(1,nw,iaer)=gvaer(2,nw,iaer)
479                  gvaer(2*nlayermx+1,nw,iaer)=0.
480
481               end do
482
483!     longwave
484               DO nw=1,L_NSPECTI
485                  do l=1,nlayermx
486
487                     temp1=QIRsQREF3d(ig,nlayermx+1-l,nw,iaer)         &
488                          *QREFir3d(ig,nlayermx+1-l,iaer)
489
490                     temp2=QIRsQREF3d(ig,max(nlayermx-l,1),nw,iaer)    &
491                          *QREFir3d(ig,max(nlayermx-l,1),iaer)
492
493                     qxiaer(2*l,nw,iaer)  = temp1
494                     qxiaer(2*l+1,nw,iaer)=(temp1+temp2)/2
495
496                     temp1=temp1*omegair3d(ig,nlayermx+1-l,nw,iaer)
497                     temp2=temp2*omegair3d(ig,max(nlayermx-l,1),nw,iaer)
498
499                     qsiaer(2*l,nw,iaer)  = temp1
500                     qsiaer(2*l+1,nw,iaer)=(temp1+temp2)/2
501
502                     temp1=gir3d(ig,nlayermx+1-l,nw,iaer)
503                     temp2=gir3d(ig,max(nlayermx-l,1),nw,iaer)
504
505                     giaer(2*l,nw,iaer)  = temp1
506                     giaer(2*l+1,nw,iaer)=(temp1+temp2)/2
507
508                  end do
509
510                  qxiaer(1,nw,iaer)=qxiaer(2,nw,iaer)
511                  qxiaer(2*nlayermx+1,nw,iaer)=0.
512
513                  qsiaer(1,nw,iaer)=qsiaer(2,nw,iaer)
514                  qsiaer(2*nlayermx+1,nw,iaer)=0.
515
516                  giaer(1,nw,iaer)=giaer(2,nw,iaer)
517                  giaer(2*nlayermx+1,nw,iaer)=0.
518
519               end do
520            end do
521
522            ! test / correct for freaky s. s. albedo values
523            do iaer=1,naerkind
524               do k=1,L_LEVELS+1
525
526                  do nw=1,L_NSPECTV
527                     if(qsvaer(k,nw,iaer).gt.1.05*qxvaer(k,nw,iaer))then
528                        print*,'Serious problems with qsvaer values in callcorrk'
529                        call abort
530                     endif
531                     if(qsvaer(k,nw,iaer).gt.qxvaer(k,nw,iaer))then
532                        qsvaer(k,nw,iaer)=qxvaer(k,nw,iaer)
533                     endif
534                  end do
535
536                  do nw=1,L_NSPECTI
537                     if(qsiaer(k,nw,iaer).gt.1.05*qxiaer(k,nw,iaer))then
538                        print*,'Serious problems with qsiaer values in callcorrk'
539                        call abort
540                     endif
541                     if(qsiaer(k,nw,iaer).gt.qxiaer(k,nw,iaer))then
542                        qsiaer(k,nw,iaer)=qxiaer(k,nw,iaer)
543                     endif
544                  end do
545
546               end do
547            end do
548
549!-----------------------------------------------------------------------
550!     Aerosol optical depths
551           
552         do iaer=1,naerkind     ! a bug was here           
553            do k=0,nlayer-1
554               
555               pweight=(pplay(ig,L_NLAYRAD-k)-pplev(ig,L_NLAYRAD-k+1))/   &
556                        (pplev(ig,L_NLAYRAD-k)-pplev(ig,L_NLAYRAD-k+1))
557
558               temp=aerosol(ig,L_NLAYRAD-k,iaer)/QREFvis3d(ig,L_NLAYRAD-k,iaer)
559
560               tauaero(2*k+2,iaer)=max(temp*pweight,0.d0)
561               tauaero(2*k+3,iaer)=max(temp-tauaero(2*k+2,iaer),0.d0)
562!
563            end do
564            ! boundary conditions
565            tauaero(1,iaer)          = tauaero(2,iaer)
566            tauaero(L_LEVELS+1,iaer) = tauaero(L_LEVELS,iaer)
567            !tauaero(1,iaer)          = 0.
568            !tauaero(L_LEVELS+1,iaer) = 0.
569         end do
570
571!     Albedo and emissivity
572         albi=1-emis(ig)        ! longwave
573         albv=albedo(ig)        ! shortwave
574
575      if(noradsurf.and.(albv.gt.0.0))then
576         print*,'For open lower boundary in callcorrk must'
577         print*,'have surface albedo set to zero!'
578         call abort
579      endif
580
581      if ((ngridmx.eq.1).and.(global1d)) then       ! fixed zenith angle 'szangle' in 1D simulations w/ globally-averaged sunlight
582         acosz = cos(pi*szangle/180.0)
583         print*,'acosz=',acosz,', szangle=',szangle
584      else
585         acosz=mu0(ig)          ! cosine of sun incident angle : 3D simulations or local 1D simulations using latitude
586      endif
587
588!-----------------------------------------------------------------------
589!     Water vapour (to be generalised for other gases eventually)
590     
591      if(varactive)then
592
593         i_var=igcm_h2o_vap
594         do l=1,nlayer
595            qvar(2*l)   = pq(ig,nlayer+1-l,i_var)
596            qvar(2*l+1) = (pq(ig,nlayer+1-l,i_var)+pq(ig,max(nlayer-l,1),i_var))/2   
597            ! Average approximation as for temperature...
598         end do
599         qvar(1)=qvar(2)
600
601      elseif(varfixed)then
602
603         do l=1,nlayermx        ! here we will assign fixed water vapour profiles globally
604            RH = satval * ((pplay(ig,l)/pplev(ig,1) - 0.02) / 0.98)
605            if(RH.lt.0.0) RH=0.0
606           
607            ptemp=pplay(ig,l)
608            Ttemp=pt(ig,l)
609            call watersat(Ttemp,ptemp,qsat)
610
611            !pq_temp(l) = qsat      ! fully saturated everywhere
612            pq_temp(l) = RH * qsat ! ~realistic profile (e.g. 80% saturation at ground)
613         end do
614         
615         do l=1,nlayer
616            qvar(2*l)   = pq_temp(nlayer+1-l)
617            qvar(2*l+1) = (pq_temp(nlayer+1-l)+pq_temp(max(nlayer-l,1)))/2
618         end do
619         qvar(1)=qvar(2)
620
621         ! Lowest layer of atmosphere
622         RH = satval * (1 - 0.02) / 0.98
623         if(RH.lt.0.0) RH=0.0
624
625         ptemp = pplev(ig,1)
626         Ttemp = tsurf(ig)
627         call watersat(Ttemp,ptemp,qsat)
628
629         !qvar(2*nlayermx+1)=qsat      ! fully saturated everywhere
630         qvar(2*nlayermx+1)= RH * qsat ! ~realistic profile (e.g. 80% saturation at ground)
631         !qvar=0.005                   ! completely constant profile (JL)
632
633      else
634         do k=1,L_LEVELS
635            qvar(k) = 1.0D-7
636         end do
637      end if
638
639      if(.not.kastprof)then
640      ! IMPORTANT: Now convert from kg/kg to mol/mol
641      do k=1,L_LEVELS
642         qvar(k) = qvar(k)/epsi
643      end do
644      end if
645
646!-----------------------------------------------------------------------
647!     kcm mode only
648      if(kastprof)then
649
650         ! initial values equivalent to mugaz
651         DO l=1,nlayer
652            muvarrad(2*l)   = mugaz
653            muvarrad(2*l+1) = mugaz
654         END DO
655
656         !do k=1,L_LEVELS
657         !   qvar(k) = 0.0
658         !end do
659         !print*,'ASSUMING qH2O=0 EVERYWHERE IN CALLCORRK!'
660      endif
661
662
663      if(kastprof.and.(ngasmx.gt.1))then
664
665         DO l=1,nlayer
666            muvarrad(2*l)   = muvar(ig,nlayer+2-l)
667            muvarrad(2*l+1) = (muvar(ig,nlayer+2-l) + &
668                                muvar(ig,max(nlayer+1-l,1)))/2
669         END DO
670     
671         muvarrad(1) = muvarrad(2)
672         muvarrad(2*nlayermx+1)=muvar(ig,1)
673
674         print*,'Recalculating qvar with VARIABLE epsi for kastprof'
675         print*,'Assumes that the variable gas is H2O!!!'
676         print*,'Assumes that there is only one tracer'
677         !i_var=igcm_h2o_vap
678         i_var=1
679         if(nqmx.gt.1)then
680            print*,'Need 1 tracer only to run kcm1d.e'
681            stop
682         endif
683         do l=1,nlayer
684            vtmp(l)=pq(ig,l,i_var)*muvar(ig,l+1)/mH2O
685         end do
686
687         do l=1,nlayer
688            qvar(2*l)   = vtmp(nlayer+1-l)
689            qvar(2*l+1) = ( vtmp(nlayer+1-l) + vtmp(max(nlayer-l,1)) )/2
690         end do
691         qvar(1)=qvar(2)
692
693         print*,'Warning: reducing qvar in callcorrk.'
694         print*,'Temperature profile no longer consistent ', &
695                            'with saturated H2O.'
696         do k=1,L_LEVELS
697            qvar(k) = qvar(k)*satval
698         end do
699
700      endif
701
702      ! Keep values inside limits for which we have radiative transfer coefficients
703      if(L_REFVAR.gt.1)then ! there was a bug here!
704         do k=1,L_LEVELS
705            if(qvar(k).lt.wrefvar(1))then
706               qvar(k)=wrefvar(1)+1.0e-8
707            elseif(qvar(k).gt.wrefvar(L_REFVAR))then
708               qvar(k)=wrefvar(L_REFVAR)-1.0e-8
709            endif
710         end do
711      endif
712
713!-----------------------------------------------------------------------
714!     Pressure and temperature
715
716      DO l=1,nlayer
717         plevrad(2*l)   = pplay(ig,nlayer+1-l)/scalep
718         plevrad(2*l+1) = pplev(ig,nlayer+1-l)/scalep
719         tlevrad(2*l)   = pt(ig,nlayer+1-l)
720         tlevrad(2*l+1) = (pt(ig,nlayer+1-l)+pt(ig,max(nlayer-l,1)))/2
721      END DO
722     
723      plevrad(1) = 0.
724      plevrad(2) = max(pgasmin,0.0001*plevrad(3))
725
726      tlevrad(1) = tlevrad(2)
727      tlevrad(2*nlayermx+1)=tsurf(ig)
728     
729      tmid(1) = tlevrad(2)
730      tmid(2) = tlevrad(2)
731      pmid(1) = plevrad(2)
732      pmid(2) = plevrad(2)
733     
734      DO l=1,L_NLAYRAD-1
735         tmid(2*l+1) = tlevrad(2*l+1)
736         tmid(2*l+2) = tlevrad(2*l+1)
737         pmid(2*l+1) = plevrad(2*l+1)
738         pmid(2*l+2) = plevrad(2*l+1)
739      END DO
740      pmid(L_LEVELS) = plevrad(L_LEVELS)
741      tmid(L_LEVELS) = tlevrad(L_LEVELS)
742
743      ! test for out-of-bounds pressure
744      if(plevrad(3).lt.pgasmin)then
745         print*,'Minimum pressure is outside the radiative'
746         print*,'transfer kmatrix bounds, exiting.'
747         call abort
748      elseif(plevrad(L_LEVELS).gt.pgasmax)then
749         print*,'Maximum pressure is outside the radiative'
750         print*,'transfer kmatrix bounds, exiting.'
751         call abort
752      endif
753
754      ! test for out-of-bounds temperature
755      do k=1,L_LEVELS
756         if(tlevrad(k).lt.tgasmin)then
757            print*,'Minimum temperature is outside the radiative'
758            print*,'transfer kmatrix bounds, exiting.'
759            !print*,'WARNING, OVERRIDING FOR TEST'
760            call abort
761         elseif(tlevrad(k).gt.tgasmax)then
762            print*,'Maximum temperature is outside the radiative'
763            print*,'transfer kmatrix bounds, exiting.'
764            !print*,'WARNING, OVERRIDING FOR TEST'
765            call abort
766         endif
767      enddo
768
769!=======================================================================
770!     Calling the main radiative transfer subroutines
771
772
773!-----------------------------------------------------------------------
774!     Shortwave
775
776         if(fract(ig) .ge. 1.0e-4) then ! only during daylight!
777
778            fluxtoplanet=0.
779
780            if((ngridmx.eq.1).and.(global1d))then
781               do nw=1,L_NSPECTV
782                  stel_fract(nw)= stel(nw) * 0.25 / acosz
783                  fluxtoplanet=fluxtoplanet + stel_fract(nw)
784                                ! globally averaged = divide by 4
785                                ! but we correct for solar zenith angle
786               end do
787            else
788               do nw=1,L_NSPECTV
789                  stel_fract(nw)= stel(nw) * fract(ig)
790                  fluxtoplanet=fluxtoplanet + stel_fract(nw)
791               end do
792            endif
793
794            call optcv(dtauv,tauv,taucumv,plevrad,                 &
795                 qxvaer,qsvaer,gvaer,wbarv,cosbv,tauray,tauaero,   &
796                 tmid,pmid,taugsurf,qvar,muvarrad)
797
798            call sfluxv(dtauv,tauv,taucumv,albv,dwnv,wbarv,cosbv,  &
799                 acosz,stel_fract,gweight,                         &
800                 nfluxtopv,nfluxoutv_nu,nfluxgndv_nu,              &
801                 fmnetv,fluxupv,fluxdnv,fzerov,taugsurf)
802
803         else                          ! during the night, fluxes = 0
804            nfluxtopv       = 0.0
805            nfluxoutv_nu(:) = 0.0
806            nfluxgndv_nu(:) = 0.0
807            do l=1,L_NLAYRAD
808               fmnetv(l)=0.0
809               fluxupv(l)=0.0
810               fluxdnv(l)=0.0
811            end do
812         end if
813
814!-----------------------------------------------------------------------
815!     Longwave
816
817         call optci(plevrad,tlevrad,dtaui,taucumi,                  &
818              qxiaer,qsiaer,giaer,cosbi,wbari,tauaero,tmid,pmid,    &
819              taugsurfi,qvar,muvarrad)
820
821         call sfluxi(plevrad,tlevrad,dtaui,taucumi,ubari,albi,      &
822              wnoi,dwni,cosbi,wbari,gweight,nfluxtopi,nfluxtopi_nu, &
823              fmneti,fluxupi,fluxdni,fluxupi_nu,fzeroi,taugsurfi)
824
825!-----------------------------------------------------------------------
826!     Transformation of the correlated-k code outputs
827!     (into dtlw, dtsw, fluxsurf_lw, fluxsurf_sw, fluxtop_lw, fluxtop_sw)
828
829!     Flux incident at the top of the atmosphere
830         fluxtop_dn(ig)=fluxdnv(1)
831
832         fluxtop_lw(ig)  = real(nfluxtopi)
833         fluxabs_sw(ig)  = real(-nfluxtopv)
834         fluxsurf_lw(ig) = real(fluxdni(L_NLAYRAD))
835         fluxsurf_sw(ig) = real(fluxdnv(L_NLAYRAD))
836
837         if(fluxtop_dn(ig).lt.0.0)then
838            print*,'Achtung! fluxtop_dn has lost the plot!'
839            print*,'fluxtop_dn=',fluxtop_dn(ig)
840            print*,'acosz=',acosz
841            print*,'aerosol=',aerosol(ig,:,:)
842            print*,'temp=   ',pt(ig,:)
843            print*,'pplay=  ',pplay(ig,:)
844            call abort
845         endif
846
847!     Spectral output, for exoplanet observational comparison
848         if(specOLR)then
849            do nw=1,L_NSPECTI
850               OLR_nu(ig,nw)=nfluxtopi_nu(nw)/DWNI(nw) !JL Normalize to the bandwidth
851            end do
852            do nw=1,L_NSPECTV
853               !GSR_nu(ig,nw)=nfluxgndv_nu(nw)
854               OSR_nu(ig,nw)=nfluxoutv_nu(nw)/DWNV(nw) !JL Normalize to the bandwidth
855            end do
856         endif
857
858!     Finally, the heating rates
859
860         DO l=2,L_NLAYRAD
861            dtsw(ig,L_NLAYRAD+1-l)=(fmnetv(l)-fmnetv(l-1))  &
862                *g/(cpp*scalep*(plevrad(2*l+1)-plevrad(2*l-1)))
863            dtlw(ig,L_NLAYRAD+1-l)=(fmneti(l)-fmneti(l-1))  &
864                *g/(cpp*scalep*(plevrad(2*l+1)-plevrad(2*l-1)))
865         END DO     
866
867!     These are values at top of atmosphere
868         dtsw(ig,L_NLAYRAD)=(fmnetv(1)-nfluxtopv)           &
869             *g/(cpp*scalep*(plevrad(3)-plevrad(1)))
870         dtlw(ig,L_NLAYRAD)=(fmneti(1)-nfluxtopi)           &
871             *g/(cpp*scalep*(plevrad(3)-plevrad(1)))
872
873!     ---------------------------------------------------------------
874      end do                    ! end of big loop over every GCM column (ig = 1:ngrid)
875
876
877!-----------------------------------------------------------------------
878!     Additional diagnostics
879
880!     IR spectral output, for exoplanet observational comparison
881
882
883      if(lastcall.and.(ngrid.eq.1))then  ! could disable the 1D output, they are in the diagfi and diagspec... JL12
884
885           print*,'Saving scalar quantities in surf_vals.out...'
886           print*,'psurf = ', pplev(1,1),' Pa'
887           open(116,file='surf_vals.out')
888           write(116,*) tsurf(1),pplev(1,1),fluxtop_dn(1),         &
889                real(-nfluxtopv),real(nfluxtopi)
890           close(116)
891
892!          I am useful, please don`t remove me!
893!           if(specOLR)then
894!               open(117,file='OLRnu.out')
895!               do nw=1,L_NSPECTI
896!                  write(117,*) OLR_nu(1,nw)
897!               enddo
898!               close(117)
899!
900!               open(127,file='OSRnu.out')
901!               do nw=1,L_NSPECTV
902!                  write(127,*) OSR_nu(1,nw)
903!               enddo
904!               close(127)
905!           endif
906
907!     OLR vs altitude: do it as a .txt file
908           OLRz=.false.
909           if(OLRz)then
910              print*,'saving IR vertical flux for OLRz...'
911              open(118,file='OLRz_plevs.out')
912              open(119,file='OLRz.out')
913              do l=1,L_NLAYRAD
914                 write(118,*) plevrad(2*l)
915                 do nw=1,L_NSPECTI
916                     write(119,*) fluxupi_nu(l,nw)
917                  enddo
918              enddo
919              close(118)
920              close(119)
921           endif
922
923      endif
924
925      ! see physiq.F for explanations about CLFvarying. This is temporary.
926      if (lastcall .and. .not.CLFvarying) then
927        IF( ALLOCATED( gasi ) ) DEALLOCATE( gasi )
928        IF( ALLOCATED( gasv ) ) DEALLOCATE( gasv )
929        IF( ALLOCATED( pgasref ) ) DEALLOCATE( pgasref )
930        IF( ALLOCATED( tgasref ) ) DEALLOCATE( tgasref )
931        IF( ALLOCATED( wrefvar ) ) DEALLOCATE( wrefvar )
932        IF( ALLOCATED( pfgasref ) ) DEALLOCATE( pfgasref )
933      endif
934
935    end subroutine callcorrk
Note: See TracBrowser for help on using the repository browser.