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

Last change on this file since 1714 was 1714, checked in by mturbet, 7 years ago

kcm1d like a phoenix rises from the ashes, part I

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