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

Last change on this file since 1423 was 1423, checked in by sglmd, 10 years ago

Second part of solving the top-of-atmosphere problem in the radiative transfer: Return to a previous interpolation scheme

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