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

Last change on this file since 545 was 538, checked in by rwordsworth, 14 years ago

calc_cpp_mugaz changed to a check only in 3D

  • 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,                      &
[305]3          tsurf,fract,dist_star,aerosol,cpp3D,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
120      real*8 taugsurf(L_NSPECTV,L_NGAUSS-1)
121      real*8 taugsurfi(L_NSPECTI,L_NGAUSS-1)
122
[305]123      real*8 qvar(L_LEVELS)          ! mixing ratio of variable component (mol/mol)
[253]124      REAL pq(ngridmx,nlayer,nq)
125      REAL qsurf(ngridmx,nqmx)       ! tracer on surface (e.g. kg.m-2)
126      integer nq
127
128!     Local aerosol optical properties for each column on RADIATIVE grid
129      real*8  QXVAER(L_LEVELS+1,L_NSPECTV,naerkind)
130      real*8  QSVAER(L_LEVELS+1,L_NSPECTV,naerkind)
131      real*8  GVAER(L_LEVELS+1,L_NSPECTV,naerkind)
132      real*8  QXIAER(L_LEVELS+1,L_NSPECTI,naerkind)
133      real*8  QSIAER(L_LEVELS+1,L_NSPECTI,naerkind)
134      real*8  GIAER(L_LEVELS+1,L_NSPECTI,naerkind)
135
136      save qxvaer, qsvaer, gvaer
137      save qxiaer, qsiaer, giaer
138      save QREFvis3d, QREFir3d
139
140      REAL tau_col(ngrid) ! diagnostic from aeropacity
141
142!     Misc.
143      logical firstcall, lastcall, nantest
144      real*8  tempv(L_NSPECTV)
145      real*8  tempi(L_NSPECTI)
146      real*8  temp,temp1,temp2,pweight
147      character(len=10) :: tmp1
148      character(len=10) :: tmp2
149
150!     for fixed water vapour profiles
151      integer i_var
152      real RH
153      real*8 pq_temp(nlayer)
154      real ptemp, Ttemp, qsat
155
156!      real(KIND=r8) :: pq_temp(nlayer) ! better F90 way.. DOESNT PORT TO F77!!!
157
[305]158      !real ptime, pday
[253]159      logical OLRz
160      real*8 NFLUXGNDV_nu(L_NSPECTV)
161
162!     for H2O cloud fraction in aeropacity
163      real cloudfrac(ngridmx,nlayermx)
164      real totcloudfrac(ngridmx)
165      logical clearsky
166
167!     Allow variations in cp with location
168      real cpp3D(ngridmx,nlayermx)   ! specific heat capacity at const. pressure
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
[253]271         firstcall=.false.   
272
273      end if
274
275!=======================================================================
276!     Initialization on every call   
277
278      do l=1,nlayer
279         do ig=1,ngrid
280            do iaer=1,naerkind
281               nueffrad(ig,l,iaer) = 0.1 ! stays at 0.1
282            enddo
283         enddo
284      enddo
285
[305]286      if(kastprof)then
287         radfixed=.true.
288      endif
289
[253]290      if(radfixed)then
291         do l=1,nlayer
292            do ig=1,ngrid
293               reffrad(ig,l,1)    = 5.e-5 ! CO2 ice
294            enddo
295         enddo
296         print*,'CO2 ice particle size = ',reffrad(1,1,1)/1.e-6,' um'
297         if(naerkind.ge.2)then
298            do l=1,nlayer
299               do ig=1,ngrid
[305]300                  !reffrad(ig,l,2) = 2.e-5 ! H2O ice
301                  reffrad(ig,l,2) = 5.e-6 ! H2O ice
[253]302               enddo
303            enddo
304            print*,'H2O ice particle size = ',reffrad(1,1,2)/1.e-6,' um'
305         endif
306         if(naerkind.eq.3)then
307            do l=1,nlayer
308               do ig=1,ngrid
309                  reffrad(ig,l,naerkind) = 2.e-6 ! dust
310               enddo
311            enddo
312            print*,'Dust particle size = ',reffrad(1,1,naerkind)/1.e-6,' um'
313         endif
314         if(naerkind.gt.3)then
315            print*,'Code not general enough to deal with naerkind > 3 yet.'
316            call abort
317         endif
318
319      else
320
321         maxrad=0.0
[305]322         !averad=0.0
[253]323         minrad=1.0
324         do l=1,nlayer
[305]325
326            !masse = (pplev(ig,l) - pplev(ig,l+1))/g
327
[253]328            do ig=1,ngrid
[538]329               if(tracer)then
330               !if(tracer.and.igcm_co2_ice.gt.0)then
331
332                  if(igcm_co2_ice.lt.1)then
333                     print*,'Tracers but no CO2 ice still seems to be a problem...'
334                     print*,'Aborting in callcorrk.'
335                     stop
336                  endif
337
[253]338                  reffrad(ig,l,1) = CBRT( 3*pq(ig,l,igcm_co2_ice)/ &
339                       (4*Nmix_co2*pi*rho_co2) )
340               endif
341               reffrad(ig,l,1) = max(reffrad(ig,l,1),1.e-6)
342               reffrad(ig,l,1) = min(reffrad(ig,l,1),500.e-6)
343
[305]344               !averad = averad + reffrad(ig,l,1)*area(ig)
[253]345               maxrad = max(reffrad(ig,l,1),maxrad)
346               minrad = min(reffrad(ig,l,1),minrad)
347            enddo
348         enddo
349         if(igcm_co2_ice.gt.0)then
350            print*,'Max. CO2 ice particle size = ',maxrad/1.e-6,' um'
351            print*,'Min. CO2 ice particle size = ',minrad/1.e-6,' um'
352         endif
353
354         if((naerkind.ge.2).and.water)then
355            maxrad=0.0
356            minrad=1.0
357            do l=1,nlayer
358               do ig=1,ngrid
359                  reffrad(ig,l,2) = CBRT( 3*pq(ig,l,igcm_h2o_ice)/ &
360                       (4*Nmix_h2o*pi*rho_ice) )
361                  reffrad(ig,l,2) = max(reffrad(ig,l,2),1.e-6)
362                  reffrad(ig,l,2) = min(reffrad(ig,l,2),100.e-6)
363
364                  maxrad = max(reffrad(ig,l,2),maxrad)
365                  minrad = min(reffrad(ig,l,2),minrad)
366               enddo
367            enddo
368            print*,'Max. H2O ice particle size = ',maxrad/1.e-6,' um'
369            print*,'Min. H2O ice particle size = ',minrad/1.e-6,' um'
370         endif
371
372         if(naerkind.eq.3)then
373            do l=1,nlayer
374               do ig=1,ngrid
375                  reffrad(ig,l,naerkind) = 2.e-6 ! dust
376               enddo
377            enddo
378         endif
379
380      endif
381
[305]382
[253]383!     how much light we get
384      fluxtoplanet=0
385      do nw=1,L_NSPECTV
386         stel(nw)=stellarf(nw)/(dist_star**2)
387         fluxtoplanet=fluxtoplanet + stel(nw)
388      end do
389
390      call aeroptproperties(ngrid,nlayer,reffrad,nueffrad,         &
391           QVISsQREF3d,omegaVIS3d,gVIS3d,                          &
392           QIRsQREF3d,omegaIR3d,gIR3d,                             &
393           QREFvis3d,QREFir3d)                                     ! get 3D aerosol optical properties
394
395      call aeropacity(ngrid,nlayer,nq,pplay,pplev,pq,aerosol,      &
396           reffrad,QREFvis3d,QREFir3d,                             &
397           tau_col,cloudfrac,totcloudfrac,clearsky)                ! get aerosol optical depths
398
[305]399
[253]400!-----------------------------------------------------------------------
401!     Starting Big Loop over every GCM column
402      do ig=1,ngridmx
403
404!=======================================================================
405!     Transformation of the GCM variables
406
407!-----------------------------------------------------------------------
408!     Aerosol optical properties Qext, Qscat and g
409!     The transformation in the vertical is the same as for temperature
410           
411!     shortwave
412            do iaer=1,naerkind
413               DO nw=1,L_NSPECTV
414                  do l=1,nlayermx
415
416                     temp1=QVISsQREF3d(ig,nlayermx+1-l,nw,iaer)         &
417                         *QREFvis3d(ig,nlayermx+1-l,iaer)
418
419                     temp2=QVISsQREF3d(ig,max(nlayermx-l,1),nw,iaer)    &
420                         *QREFvis3d(ig,max(nlayermx-l,1),iaer)
421
422                     qxvaer(2*l,nw,iaer)  = temp1
423                     qxvaer(2*l+1,nw,iaer)=(temp1+temp2)/2
424
425                     temp1=temp1*omegavis3d(ig,nlayermx+1-l,nw,iaer)
426                     temp2=temp2*omegavis3d(ig,max(nlayermx-l,1),nw,iaer)
427
428                     qsvaer(2*l,nw,iaer)  = temp1
429                     qsvaer(2*l+1,nw,iaer)=(temp1+temp2)/2
430
431                     temp1=gvis3d(ig,nlayermx+1-l,nw,iaer)
432                     temp2=gvis3d(ig,max(nlayermx-l,1),nw,iaer)
433
434                     gvaer(2*l,nw,iaer)  = temp1
435                     gvaer(2*l+1,nw,iaer)=(temp1+temp2)/2
436
437                  end do
438
439                  qxvaer(1,nw,iaer)=qxvaer(2,nw,iaer)
440                  qxvaer(2*nlayermx+1,nw,iaer)=0.
441
442                  qsvaer(1,nw,iaer)=qsvaer(2,nw,iaer)
443                  qsvaer(2*nlayermx+1,nw,iaer)=0.
444
445                  gvaer(1,nw,iaer)=gvaer(2,nw,iaer)
446                  gvaer(2*nlayermx+1,nw,iaer)=0.
447
448               end do
449
450!     longwave
451               DO nw=1,L_NSPECTI
452                  do l=1,nlayermx
453
454                     temp1=QIRsQREF3d(ig,nlayermx+1-l,nw,iaer)         &
455                          *QREFir3d(ig,nlayermx+1-l,iaer)
456
457                     temp2=QIRsQREF3d(ig,max(nlayermx-l,1),nw,iaer)    &
458                          *QREFir3d(ig,max(nlayermx-l,1),iaer)
459
460                     qxiaer(2*l,nw,iaer)  = temp1
461                     qxiaer(2*l+1,nw,iaer)=(temp1+temp2)/2
462
463                     temp1=temp1*omegair3d(ig,nlayermx+1-l,nw,iaer)
464                     temp2=temp2*omegair3d(ig,max(nlayermx-l,1),nw,iaer)
465
466                     qsiaer(2*l,nw,iaer)  = temp1
467                     qsiaer(2*l+1,nw,iaer)=(temp1+temp2)/2
468
469                     temp1=gir3d(ig,nlayermx+1-l,nw,iaer)
470                     temp2=gir3d(ig,max(nlayermx-l,1),nw,iaer)
471
472                     giaer(2*l,nw,iaer)  = temp1
473                     giaer(2*l+1,nw,iaer)=(temp1+temp2)/2
474
475                  end do
476
477                  qxiaer(1,nw,iaer)=qxiaer(2,nw,iaer)
478                  qxiaer(2*nlayermx+1,nw,iaer)=0.
479
480                  qsiaer(1,nw,iaer)=qsiaer(2,nw,iaer)
481                  qsiaer(2*nlayermx+1,nw,iaer)=0.
482
483                  giaer(1,nw,iaer)=giaer(2,nw,iaer)
484                  giaer(2*nlayermx+1,nw,iaer)=0.
485
486               end do
487            end do
488
489            ! test / correct for freaky s. s. albedo values
490            do iaer=1,naerkind
491               do k=1,L_LEVELS+1
492
493                  do nw=1,L_NSPECTV
494                     if(qsvaer(k,nw,iaer).gt.1.05*qxvaer(k,nw,iaer))then
495                        print*,'Serious problems with qsvaer values in callcorrk'
496                        call abort
497                     endif
498                     if(qsvaer(k,nw,iaer).gt.qxvaer(k,nw,iaer))then
499                        qsvaer(k,nw,iaer)=qxvaer(k,nw,iaer)
500                     endif
501                  end do
502
503                  do nw=1,L_NSPECTI
504                     if(qsiaer(k,nw,iaer).gt.1.05*qxiaer(k,nw,iaer))then
505                        print*,'Serious problems with qsiaer values in callcorrk'
506                        call abort
507                     endif
508                     if(qsiaer(k,nw,iaer).gt.qxiaer(k,nw,iaer))then
509                        qsiaer(k,nw,iaer)=qxiaer(k,nw,iaer)
510                     endif
511                  end do
512
513               end do
514            end do
515
516!-----------------------------------------------------------------------
517!     Aerosol optical depths
518           
519         do iaer=1,naerkind     ! a bug was here           
520            do k=0,nlayer-1
521               
522               pweight=(pplay(ig,L_NLAYRAD-k)-pplev(ig,L_NLAYRAD-k+1))/   &
523                        (pplev(ig,L_NLAYRAD-k)-pplev(ig,L_NLAYRAD-k+1))
524
525               temp=aerosol(ig,L_NLAYRAD-k,iaer)/QREFvis3d(ig,L_NLAYRAD-k,iaer)
526
527               tauaero(2*k+2,iaer)=max(temp*pweight,0.0)
528               tauaero(2*k+3,iaer)=max(temp-tauaero(2*k+2,iaer),0.0)
529!
530            end do
531            ! boundary conditions
532            tauaero(1,iaer)          = tauaero(2,iaer)
533            tauaero(L_LEVELS+1,iaer) = tauaero(L_LEVELS,iaer)
534            !tauaero(1,iaer)          = 0.
535            !tauaero(L_LEVELS+1,iaer) = 0.
536         end do
537
538!     Albedo and emissivity
539         albi=1-emis(ig)        ! longwave
540         albv=albedo(ig)        ! shortwave
541
542      if(noradsurf.and.(albv.gt.0.0))then
543         print*,'For open lower boundary in callcorrk must'
544         print*,'have surface albedo set to zero!'
545         call abort
546      endif
547
[305]548      if(ngridmx.eq.1) then       ! fixed zenith angle 'szangle' in 1D
[253]549         acosz = cos(pi*szangle/180.0)
550         print*,'acosz=',acosz,', szangle=',szangle
551      else
552         acosz=mu0(ig)          ! cosine of sun incident angle
553      endif
554
555!-----------------------------------------------------------------------
[305]556!     Water vapour (to be generalised for other gases eventually)
[253]557     
[305]558      if(varactive)then
[253]559
560         i_var=igcm_h2o_vap
561         do l=1,nlayer
562            qvar(2*l)   = pq(ig,nlayer+1-l,i_var)
563            qvar(2*l+1) = (pq(ig,nlayer+1-l,i_var)+pq(ig,max(nlayer-l,1),i_var))/2   
564            ! Average approximation as for temperature...
565         end do
566         qvar(1)=qvar(2)
567
568      elseif(varfixed)then
569
570         do l=1,nlayermx        ! here we will assign fixed water vapour profiles globally
571            RH = satval * ((pplay(ig,l)/pplev(ig,1) - 0.02) / 0.98)
572            if(RH.lt.0.0) RH=0.0
573           
574            ptemp=pplay(ig,l)
575            Ttemp=pt(ig,l)
576            call watersat(Ttemp,ptemp,qsat)
577
578            !pq_temp(l) = qsat      ! fully saturated everywhere
579            pq_temp(l) = RH * qsat ! ~realistic profile (e.g. 80% saturation at ground)
580         end do
581         
582         do l=1,nlayer
583            qvar(2*l)   = pq_temp(nlayer+1-l)
584            qvar(2*l+1) = (pq_temp(nlayer+1-l)+pq_temp(max(nlayer-l,1)))/2
585         end do
586         qvar(1)=qvar(2)
587
588         ! Lowest layer of atmosphere
589         RH = satval * (1 - 0.02) / 0.98
590         if(RH.lt.0.0) RH=0.0
591
592         ptemp = pplev(ig,1)
593         Ttemp = tsurf(ig)
594         call watersat(Ttemp,ptemp,qsat)
595
596         !qvar(2*nlayermx+1)=qsat      ! fully saturated everywhere
597         qvar(2*nlayermx+1)= RH * qsat ! ~realistic profile (e.g. 80% saturation at ground)
[538]598         !qvar=0.005                   ! completely constant profile (JL)
[253]599
600      else
601         do k=1,L_LEVELS
602            qvar(k) = 1.0D-7
603         end do
604      end if
605
[538]606      if(.not.kastprof)then
[305]607      ! IMPORTANT: Now convert from kg/kg to mol/mol
[253]608      do k=1,L_LEVELS
609         qvar(k) = qvar(k)/epsi
610      end do
[538]611      end if
[253]612
[366]613!-----------------------------------------------------------------------
614!     kcm mode only
[305]615      if(kastprof)then
616
[486]617         ! initial values equivalent to mugaz
[305]618         DO l=1,nlayer
[366]619            muvarrad(2*l)   = mugaz
620            muvarrad(2*l+1) = mugaz
621         END DO
622
[486]623         !do k=1,L_LEVELS
624         !   qvar(k) = 0.0
625         !end do
626         !print*,'ASSUMING qH2O=0 EVERYWHERE IN CALLCORRK!'
[366]627      endif
628
629
630      if(kastprof.and.(ngasmx.gt.1))then
631
632         DO l=1,nlayer
[305]633            muvarrad(2*l)   = muvar(ig,nlayer+2-l)
634            muvarrad(2*l+1) = (muvar(ig,nlayer+2-l) + &
635                                muvar(ig,max(nlayer+1-l,1)))/2
636         END DO
637     
638         muvarrad(1) = muvarrad(2)
639         muvarrad(2*nlayermx+1)=muvar(ig,1)
640
641         print*,'Recalculating qvar with VARIABLE epsi for kastprof'
[538]642         print*,'Assumes that the variable gas is H2O!!!'
643         print*,'Assumes that there is only one tracer'
644         !i_var=igcm_h2o_vap
645         i_var=1
646         if(nqmx.gt.1)then
647            print*,'Need 1 tracer only to run kcm1d.e'
648            stop
649         endif
[305]650         do l=1,nlayer
651            vtmp(l)=pq(ig,l,i_var)*muvar(ig,l+1)/mH2O
652         end do
653
654         do l=1,nlayer
655            qvar(2*l)   = vtmp(nlayer+1-l)
656            qvar(2*l+1) = ( vtmp(nlayer+1-l) + vtmp(max(nlayer-l,1)) )/2
657         end do
658         qvar(1)=qvar(2)
659
[538]660         print*,'Warning: reducing qvar in callcorrk.'
661         print*,'Temperature profile no longer consistent ', &
662                            'with saturated H2O.'
663         do k=1,L_LEVELS
664            qvar(k) = qvar(k)*satval
665         end do
666
[305]667      endif
668
[253]669      ! Keep values inside limits for which we have radiative transfer coefficients
670      if(L_REFVAR.gt.1)then ! there was a bug here!
671         do k=1,L_LEVELS
672            if(qvar(k).lt.wrefvar(1))then
673               qvar(k)=wrefvar(1)+1.0e-8
674            elseif(qvar(k).gt.wrefvar(L_REFVAR))then
675               qvar(k)=wrefvar(L_REFVAR)-1.0e-8
676            endif
677         end do
678      endif
679
680!-----------------------------------------------------------------------
681!     Pressure and temperature
682
683      DO l=1,nlayer
684         plevrad(2*l)   = pplay(ig,nlayer+1-l)/scalep
685         plevrad(2*l+1) = pplev(ig,nlayer+1-l)/scalep
686         tlevrad(2*l)   = pt(ig,nlayer+1-l)
687         tlevrad(2*l+1) = (pt(ig,nlayer+1-l)+pt(ig,max(nlayer-l,1)))/2
688      END DO
689     
690      plevrad(1) = 0
691      plevrad(2) = max(pgasmin,0.0001*plevrad(3))
692
693      tlevrad(1) = tlevrad(2)
694      tlevrad(2*nlayermx+1)=tsurf(ig)
695     
696      tmid(1) = tlevrad(2)
697      tmid(2) = tlevrad(2)
698      pmid(1) = plevrad(2)
699      pmid(2) = plevrad(2)
700     
701      DO l=1,L_NLAYRAD-1
702         tmid(2*l+1) = tlevrad(2*l+1)
703         tmid(2*l+2) = tlevrad(2*l+1)
704         pmid(2*l+1) = plevrad(2*l+1)
705         pmid(2*l+2) = plevrad(2*l+1)
706      END DO
707      pmid(L_LEVELS) = plevrad(L_LEVELS)
708      tmid(L_LEVELS) = tlevrad(L_LEVELS)
709
710      ! test for out-of-bounds pressure
711      if(plevrad(3).lt.pgasmin)then
712         print*,'Minimum pressure is outside the radiative'
713         print*,'transfer kmatrix bounds, exiting.'
714         call abort
715      elseif(plevrad(L_LEVELS).gt.pgasmax)then
716         print*,'Maximum pressure is outside the radiative'
717         print*,'transfer kmatrix bounds, exiting.'
718         call abort
719      endif
720
721      ! test for out-of-bounds temperature
722      do k=1,L_LEVELS
723         if(tlevrad(k).lt.tgasmin)then
724            print*,'Minimum temperature is outside the radiative'
725            print*,'transfer kmatrix bounds, exiting.'
[486]726            !print*,'WARNING, OVERRIDING FOR TEST'
727            call abort
[253]728         elseif(tlevrad(k).gt.tgasmax)then
729            print*,'Maximum temperature is outside the radiative'
730            print*,'transfer kmatrix bounds, exiting.'
[486]731            !print*,'WARNING, OVERRIDING FOR TEST'
732            call abort
[253]733         endif
734      enddo
735
736!=======================================================================
737!     Calling the main radiative transfer subroutines
738
739
740!-----------------------------------------------------------------------
741!     Shortwave
742
743         if(fract(ig) .ge. 1.0e-4) then ! only during daylight!
744
745            fluxtoplanet=0.
746
747            if((ngridmx.eq.1).and.(.not.(diurnal.or.tlocked)))then
748               do nw=1,L_NSPECTV
749                  stel_fract(nw)= stel(nw) * 0.25 / acosz
750                  fluxtoplanet=fluxtoplanet + stel_fract(nw)
751                                ! globally averaged = divide by 4
752                                ! but we correct for solar zenith angle
753               end do
754            else
755               do nw=1,L_NSPECTV
756                  stel_fract(nw)= stel(nw) * fract(ig)
757                  fluxtoplanet=fluxtoplanet + stel_fract(nw)
758               end do
759            endif
760
761            call optcv(dtauv,tauv,taucumv,plevrad,                 &
762                 qxvaer,qsvaer,gvaer,wbarv,cosbv,tauray,tauaero,   &
[305]763                 tmid,pmid,taugsurf,qvar,muvarrad)
[253]764
765            call sfluxv(dtauv,tauv,taucumv,albv,dwnv,wbarv,cosbv,  &
[366]766                 acosz,stel_fract,gweight,                         &
767                 nfluxtopv,nfluxoutv_nu,nfluxgndv_nu,              &
[253]768                 fmnetv,fluxupv,fluxdnv,fzerov,taugsurf)
769
770         else                          ! during the night, fluxes = 0
[366]771            nfluxtopv       = 0.0
772            nfluxoutv_nu(:) = 0.0
773            nfluxgndv_nu(:) = 0.0
[253]774            do l=1,L_NLAYRAD
775               fmnetv(l)=0.0
776               fluxupv(l)=0.0
777               fluxdnv(l)=0.0
778            end do
779         end if
780
781!-----------------------------------------------------------------------
782!     Longwave
783
784         call optci(plevrad,tlevrad,dtaui,taucumi,                  &
785              qxiaer,qsiaer,giaer,cosbi,wbari,tauaero,tmid,pmid,    &
[305]786              taugsurfi,qvar,muvarrad)
[538]787
[253]788         call sfluxi(plevrad,tlevrad,dtaui,taucumi,ubari,albi,      &
789              wnoi,dwni,cosbi,wbari,gweight,nfluxtopi,nfluxtopi_nu, &
790              fmneti,fluxupi,fluxdni,fluxupi_nu,fzeroi,taugsurfi)
791
792!-----------------------------------------------------------------------
793!     Transformation of the correlated-k code outputs
794!     (into dtlw, dtsw, fluxsurf_lw, fluxsurf_sw, fluxtop_lw, fluxtop_sw)
795
796!     Flux incident at the top of the atmosphere
797         fluxtop_dn(ig)=fluxdnv(1)
798
799         fluxtop_lw(ig)  = real(nfluxtopi)
800         fluxabs_sw(ig)  = real(-nfluxtopv)
801         fluxsurf_lw(ig) = real(fluxdni(L_NLAYRAD))
802         fluxsurf_sw(ig) = real(fluxdnv(L_NLAYRAD))
803
[305]804
[253]805         if(fluxtop_dn(ig).lt.0.0)then
806            print*,'Achtung! fluxtop_dn has lost the plot!'
807            print*,'fluxtop_dn=',fluxtop_dn(ig)
808            print*,'acosz=',acosz
809            print*,'aerosol=',aerosol(ig,:,:)
810            print*,'temp=   ',pt(ig,:)
811            print*,'pplay=  ',pplay(ig,:)
812            call abort
813         endif
814
815!     Spectral output, for exoplanet observational comparison
816         if(specOLR)then
817            do nw=1,L_NSPECTI
[526]818               OLR_nu(ig,nw)=nfluxtopi_nu(nw)/DWNI(nw) !JL Normalize to the bandwidth
[253]819            end do
820            do nw=1,L_NSPECTV
[366]821               !GSR_nu(ig,nw)=nfluxgndv_nu(nw)
[526]822               OSR_nu(ig,nw)=nfluxoutv_nu(nw)/DWNV(nw) !JL Normalize to the bandwidth
[253]823            end do
824         endif
825
826!     Finally, the heating rates
827         if(nonideal)then
828
829            DO l=2,L_NLAYRAD
830               dtsw(ig,L_NLAYRAD+1-l)=(fmnetv(l)-fmnetv(l-1))  &
831                    *g/(cpp3D(ig,L_NLAYRAD+1-l)                &
832                    *scalep*(plevrad(2*l+1)-plevrad(2*l-1)))
833               dtlw(ig,L_NLAYRAD+1-l)=(fmneti(l)-fmneti(l-1))  &
834                    *g/(cpp3D(ig,L_NLAYRAD+1-l)                &
835                    *scalep*(plevrad(2*l+1)-plevrad(2*l-1)))           
836            END DO     
837
838!     These are values at top of atmosphere
839            dtsw(ig,L_NLAYRAD)=(fmnetv(1)-nfluxtopv)           &
840                 *g/(cpp3D(ig,L_NLAYRAD)*scalep*(plevrad(3)-plevrad(1)))
841            dtlw(ig,L_NLAYRAD)=(fmneti(1)-nfluxtopi)           &
842                 *g/(cpp3D(ig,L_NLAYRAD)*scalep*(plevrad(3)-plevrad(1)))
843
844         else
845
846            DO l=2,L_NLAYRAD
847               dtsw(ig,L_NLAYRAD+1-l)=(fmnetv(l)-fmnetv(l-1))  &
848                   *g/(cpp*scalep*(plevrad(2*l+1)-plevrad(2*l-1)))
849               dtlw(ig,L_NLAYRAD+1-l)=(fmneti(l)-fmneti(l-1))  &
850                   *g/(cpp*scalep*(plevrad(2*l+1)-plevrad(2*l-1)))
851            END DO     
852
853!     These are values at top of atmosphere
854            dtsw(ig,L_NLAYRAD)=(fmnetv(1)-nfluxtopv)           &
855                *g/(cpp*scalep*(plevrad(3)-plevrad(1)))
856            dtlw(ig,L_NLAYRAD)=(fmneti(1)-nfluxtopi)           &
857                *g/(cpp*scalep*(plevrad(3)-plevrad(1)))
858
859         endif
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.