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

Last change on this file since 997 was 962, checked in by jleconte, 12 years ago

finish correction ISR calculation in callcorrk

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