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

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