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

Last change on this file since 2470 was 2470, checked in by yjaziri, 4 years ago

Generic GCM:
global1d and szangle for 1D simulation moved from callcorrk to callkeys
to defined a consistent 1D sza in physiq_mod used also in chemistry
+ typo
YJ

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