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

Last change on this file since 3922 was 3922, checked in by mmaurice, 8 weeks ago

Generic PCM:

Add time-dependent dust scenario for (present-day) Mars from Montmessin
et al., 2004.

MM

  • Property svn:executable set to *
File size: 55.7 KB
Line 
1MODULE callcorrk_mod
2
3IMPLICIT NONE
4
5CONTAINS
6
7      subroutine callcorrk(ngrid,nlayer,pq,nq,qsurf,zls,       &
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,GSR_nu,                                &
14          int_dtaui,int_dtauv,                                 &
15          tau_col,cloudfrac,totcloudfrac,                      &
16          clearsky,p_var,frac_var,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
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, &
31                              iaero_back2lay, iaero_aurora,               &
32                              iaero_venus1,iaero_venus2,iaero_venus2p,    &
33                              iaero_venus3,iaero_venusUV
34      use aeropacity_mod, only: aeropacity
35      use aeroptproperties_mod, only: aeroptproperties
36      use tracer_h, only: igcm_h2o_ice, igcm_h2o_vap, igcm_co2_ice
37      use tracer_h, only: constants_epsi_generic
38      use comcstfi_mod, only: pi, mugaz, cpp
39      use callkeys_mod, only: varactive,diurnal,tracer,water,varfixed,satval, &
40                              diagdtau,kastprof,strictboundcorrk,specOLR, &
41                              CLFvarying,tplanckmin,tplanckmax,global1d, &
42                              generic_condensation, aerovenus, nvarlayer, varspec
43      use optcv_mod, only: optcv
44      use optci_mod, only: optci
45      use sfluxi_mod, only: sfluxi
46      use sfluxv_mod, only: sfluxv
47      use recombin_corrk_mod, only: corrk_recombin, call_recombin
48      use pindex_mod, only: pindex
49      use generic_cloud_common_h, only: Psat_generic, epsi_generic
50      use generic_tracer_index_mod, only: generic_tracer_index
51      use planetwide_mod, only: planetwide_maxval, planetwide_minval
52      implicit none
53
54!==================================================================
55!
56!     Purpose
57!     -------
58!     Solve the radiative transfer using the correlated-k method for
59!     the gaseous absorption and the Toon et al. (1989) method for
60!     scatttering due to aerosols.
61!
62!     Authors
63!     -------
64!     Emmanuel 01/2001, Forget 09/2001
65!     Robin Wordsworth (2009)
66!
67!==================================================================
68
69!-----------------------------------------------------------------------
70!     Declaration of the arguments (INPUT - OUTPUT) on the LMD GCM grid
71!     Layer #1 is the layer near the ground.
72!     Layer #nlayer is the layer at the top.
73!-----------------------------------------------------------------------
74
75
76      ! INPUT
77      INTEGER,INTENT(IN) :: ngrid                  ! Number of atmospheric columns.
78      INTEGER,INTENT(IN) :: nlayer                 ! Number of atmospheric layers.
79      REAL,INTENT(IN) :: pq(ngrid,nlayer,nq)       ! Tracers (kg/kg_of_air).
80      INTEGER,INTENT(IN) :: nq                     ! Number of tracers.
81      REAL,INTENT(IN) :: qsurf(ngrid,nq)           ! Tracers on surface (kg.m-2).
82      REAL,INTENT(IN) :: zls                       ! Stellar longitude (rad).
83      REAL,INTENT(IN) :: albedo(ngrid,L_NSPECTV)   ! Spectral Short Wavelengths Albedo. By MT2015
84      REAL,INTENT(IN) :: emis(ngrid)               ! Long Wave emissivity.
85      REAL,INTENT(IN) :: mu0(ngrid)                ! Cosine of sun incident angle.
86      REAL,INTENT(IN) :: pplev(ngrid,nlayer+1)     ! Inter-layer pressure (Pa).
87      REAL,INTENT(IN) :: pplay(ngrid,nlayer)       ! Mid-layer pressure (Pa).
88      REAL,INTENT(IN) :: pt(ngrid,nlayer)          ! Air temperature (K).
89      REAL,INTENT(IN) :: tsurf(ngrid)              ! Surface temperature (K).
90      REAL,INTENT(IN) :: fract(ngrid)              ! Fraction of day.
91      REAL,INTENT(IN) :: dist_star                 ! Distance star-planet (AU).
92      REAL,INTENT(IN) :: muvar(ngrid,nlayer+1)
93      REAL,INTENT(IN) :: cloudfrac(ngrid,nlayer)   ! Fraction of clouds (%).
94      REAL,INTENT(IN) :: frac_var(nvarlayer,ngasmx)! Variable molar fraction.
95      REAL,INTENT(IN) :: p_var(nvarlayer)          ! Pressure for frac_var interpolation (Pa)
96      logical,intent(in) :: clearsky
97      logical,intent(in) :: firstcall              ! Signals first call to physics.
98      logical,intent(in) :: lastcall               ! Signals last call to physics.
99     
100      ! OUTPUT
101      REAL,INTENT(OUT) :: aerosol(ngrid,nlayer,naerkind) ! Aerosol tau at reference wavelenght.
102      REAL,INTENT(OUT) :: dtlw(ngrid,nlayer)             ! Heating rate (K/s) due to LW radiation.
103      REAL,INTENT(OUT) :: dtsw(ngrid,nlayer)             ! Heating rate (K/s) due to SW radiation.
104      REAL,INTENT(OUT) :: fluxsurf_lw(ngrid)             ! Incident LW flux to surf (W/m2).
105      REAL,INTENT(OUT) :: fluxsurf_sw(ngrid)             ! Incident SW flux to surf (W/m2)
106      REAL,INTENT(OUT) :: fluxsurfabs_sw(ngrid)          ! Absorbed SW flux by the surface (W/m2). By MT2015.
107      REAL,INTENT(OUT) :: fluxtop_lw(ngrid)              ! Outgoing LW flux to space (W/m2).
108      REAL,INTENT(OUT) :: fluxabs_sw(ngrid)              ! SW flux absorbed by the planet (W/m2).
109      REAL,INTENT(OUT) :: fluxtop_dn(ngrid)              ! Incident top of atmosphere SW flux (W/m2).
110      REAL,INTENT(OUT) :: OLR_nu(ngrid,L_NSPECTI)        ! Outgoing LW radiation in each band (Normalized to the band width (W/m2/cm-1).
111      REAL,INTENT(OUT) :: OSR_nu(ngrid,L_NSPECTV)        ! Outgoing SW radiation in each band (Normalized to the band width (W/m2/cm-1).
112      REAL,INTENT(OUT) :: GSR_nu(ngrid,L_NSPECTV)        ! Surface SW radiation in each band (Normalized to the band width (W/m2/cm-1).
113      REAL,INTENT(OUT) :: tau_col(ngrid)                 ! Diagnostic from aeropacity.
114      REAL,INTENT(OUT) :: albedo_equivalent(ngrid)       ! Spectrally Integrated Albedo. For Diagnostic. By MT2015
115      REAL,INTENT(OUT) :: totcloudfrac(ngrid)            ! Column Fraction of clouds (%).
116      REAL,INTENT(OUT) :: int_dtaui(ngrid,nlayer,L_NSPECTI) ! VI optical thickness of layers within narrowbands for diags ().
117      REAL,INTENT(OUT) :: int_dtauv(ngrid,nlayer,L_NSPECTV) ! IR optical thickness of layers within narrowbands for diags ().
118     
119     
120     
121     
122
123! Globally varying aerosol optical properties on GCM grid ; not needed everywhere so not in radcommon_h.   
124! made "save" variables so they are allocated once in for all, not because
125! the values need be saved from a time step to the next
126      REAL,SAVE,ALLOCATABLE :: QVISsQREF3d(:,:,:,:)
127      REAL,SAVE,ALLOCATABLE :: omegaVIS3d(:,:,:,:)
128      REAL,SAVE,ALLOCATABLE :: gVIS3d(:,:,:,:)
129!$OMP THREADPRIVATE(QVISsQREF3d,omegaVIS3d,gVIS3d)
130      REAL,SAVE,ALLOCATABLE :: QIRsQREF3d(:,:,:,:)
131      REAL,SAVE,ALLOCATABLE :: omegaIR3d(:,:,:,:)
132      REAL,SAVE,ALLOCATABLE :: gIR3d(:,:,:,:)
133!$OMP THREADPRIVATE(QIRsQREF3d,omegaIR3d,gIR3d)
134
135!      REAL :: omegaREFvis3d(ngrid,nlayer,naerkind)
136!      REAL :: omegaREFir3d(ngrid,nlayer,naerkind) ! not sure of the point of these...
137
138      REAL,ALLOCATABLE,SAVE :: reffrad(:,:,:)  ! aerosol effective radius (m)
139      REAL,ALLOCATABLE,SAVE :: nueffrad(:,:,:) ! aerosol effective variance
140!$OMP THREADPRIVATE(reffrad,nueffrad)
141
142!-----------------------------------------------------------------------
143!     Declaration of the variables required by correlated-k subroutines
144!     Numbered from top to bottom (unlike in the GCM)
145!-----------------------------------------------------------------------
146
147      REAL*8 tmid(L_LEVELS),pmid(L_LEVELS)
148      REAL*8 tlevrad(L_LEVELS),plevrad(L_LEVELS)
149
150      ! Optical values for the optci/cv subroutines
151      REAL*8 stel(L_NSPECTV),stel_fract(L_NSPECTV)
152      ! NB: Arrays below are "save" to avoid reallocating them at every call
153      ! not because their content needs be reused from call to the next
154      REAL*8,allocatable,save :: dtaui(:,:,:)
155      REAL*8,allocatable,save :: dtauv(:,:,:)
156      REAL*8,allocatable,save :: cosbv(:,:,:)
157      REAL*8,allocatable,save :: cosbi(:,:,:)
158      REAL*8,allocatable,save :: wbari(:,:,:)
159      REAL*8,allocatable,save :: wbarv(:,:,:)
160!$OMP THREADPRIVATE(dtaui,dtauv,cosbv,cosbi,wbari,wbarv)
161      REAL*8,allocatable,save :: tauv(:,:,:)
162      REAL*8,allocatable,save :: taucumv(:,:,:)
163      REAL*8,allocatable,save :: taucumi(:,:,:)
164!$OMP THREADPRIVATE(tauv,taucumv,taucumi)
165      REAL*8,allocatable,save :: tauaero(:,:)
166!$OMP THREADPRIVATE(tauaero)
167      REAL*8 nfluxtopv,nfluxtopi,nfluxtop,fluxtopvdn
168      REAL*8 nfluxoutv_nu(L_NSPECTV)                 ! Outgoing band-resolved VI flux at TOA (W/m2).
169      REAL*8 nfluxtopi_nu(L_NSPECTI)                 ! Net band-resolved IR flux at TOA (W/m2).
170      REAL*8 fluxupi_nu(L_NLAYRAD,L_NSPECTI)         ! For 1D diagnostic.
171      REAL*8 fmneti(L_NLAYRAD),fmnetv(L_NLAYRAD)
172      REAL*8 fluxupv(L_NLAYRAD),fluxupi(L_NLAYRAD)
173      REAL*8 fluxdnv(L_NLAYRAD),fluxdni(L_NLAYRAD)
174      REAL*8 albi,acosz
175      REAL*8 albv(L_NSPECTV)                         ! Spectral Visible Albedo.
176
177      INTEGER ig,l,k,nw,iaer,iq
178
179      real*8,allocatable,save :: taugsurf(:,:)
180      real*8,allocatable,save :: taugsurfi(:,:)
181!$OMP THREADPRIVATE(taugsurf,taugsurfi)
182      real*8 qvar(L_LEVELS)   ! Mixing ratio of variable component (mol/mol). index 1 is the top of the atmosphere, index L_LEVELS is the bottom
183
184      ! Local aerosol optical properties for each column on RADIATIVE grid.
185      real*8,save,allocatable ::  QXVAER(:,:,:) ! Extinction coeff (QVISsQREF*QREFvis)
186      real*8,save,allocatable ::  QSVAER(:,:,:)
187      real*8,save,allocatable ::  GVAER(:,:,:)
188      real*8,save,allocatable ::  QXIAER(:,:,:) ! Extinction coeff (QIRsQREF*QREFir)
189      real*8,save,allocatable ::  QSIAER(:,:,:)
190      real*8,save,allocatable ::  GIAER(:,:,:)
191!$OMP THREADPRIVATE(QXVAER,QSVAER,GVAER,QXIAER,QSIAER,GIAER)
192      real, dimension(:,:,:), save, allocatable :: QREFvis3d
193      real, dimension(:,:,:), save, allocatable :: QREFir3d
194!$OMP THREADPRIVATE(QREFvis3d,QREFir3d)
195
196
197      ! Miscellaneous :
198      real*8  temp,temp1,temp2,pweight
199      character(len=10) :: tmp1
200      character(len=10) :: tmp2
201      character(len=100) :: message
202      character(len=10),parameter :: subname="callcorrk"
203
204      ! For fixed water vapour profiles.
205      integer i_var
206      real RH
207      real*8 pq_temp(nlayer)
208! real(KIND=r8) :: pq_temp(nlayer) ! better F90 way.. DOESNT PORT TO F77!!!
209      real psat,qsat
210
211      logical OLRz
212      real*8 NFLUXGNDV_nu(L_NSPECTV)
213
214      ! Included by RW for runaway greenhouse 1D study.
215      real vtmp(nlayer)
216      REAL*8 muvarrad(L_LEVELS)
217     
218      ! Included by MT for albedo calculations.     
219      REAL*8 albedo_temp(L_NSPECTV) ! For equivalent albedo calculation.
220      REAL*8 surface_stellar_flux   ! Stellar flux reaching the surface. Useful for equivalent albedo calculation.
221     
222      ! local variable
223      integer ok ! status (returned by NetCDF functions)
224
225      integer igcm_generic_vap, igcm_generic_ice! index of the vap and ice of generic_tracer
226      logical call_ice_vap_generic ! to call only one time the ice/vap pair of a tracer
227      real, save :: metallicity ! metallicity of planet --- is not used here, but necessary to call function Psat_generic
228!$OMP THREADPRIVATE(metallicity)
229      REAL, SAVE :: qvap_deep   ! deep mixing ratio of water vapor when simulating bottom less planets
230!$OMP THREADPRIVATE(qvap_deep)
231
232      REAL :: maxvalue,minvalue
233     
234      real :: frac_vari(L_LEVELS)
235      real :: fracvari(ngasmx,L_LEVELS)
236
237!===============================================================
238!           I.a Initialization on first call
239!===============================================================
240
241
242      if(firstcall) then
243
244        ! test on allocated necessary because of CLFvarying (two calls to callcorrk in physiq)
245        if(.not.allocated(QVISsQREF3d)) then
246          allocate(QVISsQREF3d(ngrid,nlayer,L_NSPECTV,naerkind))
247        endif
248        if(.not.allocated(omegaVIS3d)) then
249          allocate(omegaVIS3d(ngrid,nlayer,L_NSPECTV,naerkind))
250        endif
251        if(.not.allocated(gVIS3d)) then
252          allocate(gVIS3d(ngrid,nlayer,L_NSPECTV,naerkind))
253        endif
254        if (.not.allocated(QIRsQREF3d)) then
255          allocate(QIRsQREF3d(ngrid,nlayer,L_NSPECTI,naerkind))
256        endif
257        if (.not.allocated(omegaIR3d)) then
258          allocate(omegaIR3d(ngrid,nlayer,L_NSPECTI,naerkind))
259        endif
260        if (.not.allocated(gIR3d)) then
261          allocate(gIR3d(ngrid,nlayer,L_NSPECTI,naerkind))
262        endif
263        if (.not.allocated(tauaero)) then
264          allocate(tauaero(L_LEVELS,naerkind))
265        endif
266       
267        if(.not.allocated(QXVAER)) then
268          allocate(QXVAER(L_LEVELS,L_NSPECTV,naerkind), stat=ok)
269          if (ok /= 0) then
270             write(*,*) "memory allocation failed for QXVAER!"
271             call abort_physic(subname,'allocation failure for QXVAER',1)
272          endif
273        endif
274        if(.not.allocated(QSVAER)) then
275          allocate(QSVAER(L_LEVELS,L_NSPECTV,naerkind), stat=ok)
276          if (ok /= 0) then
277             write(*,*) "memory allocation failed for QSVAER!"
278             call abort_physic(subname,'allocation failure for QSVAER',1)
279          endif
280        endif
281        if(.not.allocated(GVAER)) then
282          allocate(GVAER(L_LEVELS,L_NSPECTV,naerkind), stat=ok)
283          if (ok /= 0) then
284             write(*,*) "memory allocation failed for GVAER!"
285             call abort_physic(subname,'allocation failure for GVAER',1)
286          endif
287        endif
288        if(.not.allocated(QXIAER)) then
289          allocate(QXIAER(L_LEVELS,L_NSPECTI,naerkind), stat=ok)
290          if (ok /= 0) then
291             write(*,*) "memory allocation failed for QXIAER!"
292             call abort_physic(subname,'allocation failure for QXIAER',1)
293          endif
294        endif
295        if(.not.allocated(QSIAER)) then
296          allocate(QSIAER(L_LEVELS,L_NSPECTI,naerkind), stat=ok)
297          if (ok /= 0) then
298             write(*,*) "memory allocation failed for QSIAER!"
299             call abort_physic(subname,'allocation failure for QSIAER',1)
300          endif
301        endif
302        if(.not.allocated(GIAER)) then
303          allocate(GIAER(L_LEVELS,L_NSPECTI,naerkind), stat=ok)
304          if (ok /= 0) then
305             write(*,*) "memory allocation failed for GIAER!"
306             call abort_physic(subname,'allocation failure for GIAER',1)
307          endif
308        endif
309
310         !!! ALLOCATED instances are necessary because of CLFvarying (strategy to call callcorrk twice in physiq...)
311         IF(.not.ALLOCATED(QREFvis3d))THEN
312           ALLOCATE(QREFvis3d(ngrid,nlayer,naerkind), stat=ok)
313           IF (ok/=0) THEN
314              write(*,*) "memory allocation failed for QREFvis3d!"
315              call abort_physic(subname,'allocation failure for QREFvis3d',1)
316           ENDIF
317         ENDIF
318         IF(.not.ALLOCATED(QREFir3d)) THEN
319           ALLOCATE(QREFir3d(ngrid,nlayer,naerkind), stat=ok)
320           IF (ok/=0) THEN
321              write(*,*) "memory allocation failed for QREFir3d!"
322              call abort_physic(subname,'allocation failure for QREFir3d',1)
323           ENDIF
324         ENDIF
325         ! Effective radius and variance of the aerosols
326         IF(.not.ALLOCATED(reffrad)) THEN
327           allocate(reffrad(ngrid,nlayer,naerkind), stat=ok)
328           IF (ok/=0) THEN
329              write(*,*) "memory allocation failed for reffrad!"
330              call abort_physic(subname,'allocation failure for reffrad',1)
331           ENDIF
332         ENDIF
333         IF(.not.ALLOCATED(nueffrad)) THEN
334           allocate(nueffrad(ngrid,nlayer,naerkind), stat=ok)
335           IF (ok/=0) THEN
336              write(*,*) "memory allocation failed for nueffrad!"
337              call abort_physic(subname,'allocation failure for nueffrad',1)
338           ENDIF
339         ENDIF
340
341#ifndef MESOSCALE
342         if (is_master) call system('rm -f surf_vals_long.out')
343#endif
344
345         call su_aer_radii(ngrid,nlayer,reffrad,nueffrad)
346         
347         
348!--------------------------------------------------
349!             Set up correlated k
350!--------------------------------------------------
351
352      !this block is now done at firstcall of physiq_mod
353         ! print*, "callcorrk: Correlated-k data base folder:",trim(datadir)
354         ! call getin_p("corrkdir",corrkdir)
355         ! print*, "corrkdir = ",corrkdir
356         ! write( tmp1, '(i3)' ) L_NSPECTI
357         ! write( tmp2, '(i3)' ) L_NSPECTV
358         ! banddir=trim(adjustl(tmp1))//'x'//trim(adjustl(tmp2))
359         ! banddir=trim(adjustl(corrkdir))//'/'//trim(adjustl(banddir))
360
361         ! call setspi            ! Basic infrared properties.
362         ! call setspv            ! Basic visible properties.
363         ! call sugas_corrk       ! Set up gaseous absorption properties.
364         ! call suaer_corrk       ! Set up aerosol optical properties.
365       
366
367         ! now that L_NGAUSS has been initialized (by sugas_corrk)
368         ! allocate related arrays
369         if(.not.allocated(dtaui)) then
370           ALLOCATE(dtaui(L_NLAYRAD,L_NSPECTI,L_NGAUSS), stat=ok)
371           if (ok/=0) then
372              write(*,*) "memory allocation failed for dtaui!"
373              call abort_physic(subname,'allocation failure for dtaui',1)
374           endif
375         endif
376         if(.not.allocated(dtauv)) then
377           ALLOCATE(dtauv(L_NLAYRAD,L_NSPECTV,L_NGAUSS), stat=ok)
378           if (ok/=0) then
379              write(*,*) "memory allocation failed for dtauv!"
380              call abort_physic(subname,'allocation failure for dtauv',1)
381           endif
382         endif
383         if(.not.allocated(cosbv)) then
384           ALLOCATE(cosbv(L_NLAYRAD,L_NSPECTV,L_NGAUSS), stat=ok)
385           if (ok/=0) then
386              write(*,*) "memory allocation failed for cosbv!"
387              call abort_physic(subname,'allocation failure for cobsv',1)
388           endif
389         endif
390         if(.not.allocated(cosbi)) then
391           ALLOCATE(cosbi(L_NLAYRAD,L_NSPECTI,L_NGAUSS), stat=ok)
392           if (ok/=0) then
393              write(*,*) "memory allocation failed for cosbi!"
394              call abort_physic(subname,'allocation failure for cobsi',1)
395           endif
396         endif
397         if(.not.allocated(wbari)) then
398           ALLOCATE(wbari(L_NLAYRAD,L_NSPECTI,L_NGAUSS), stat=ok)
399           if (ok/=0) then
400              write(*,*) "memory allocation failed for wbari!"
401              call abort_physic(subname,'allocation failure for wbari',1)
402           endif
403         endif
404         if(.not.allocated(wbarv)) then
405           ALLOCATE(wbarv(L_NLAYRAD,L_NSPECTV,L_NGAUSS), stat=ok)
406           if (ok/=0) then
407              write(*,*) "memory allocation failed for wbarv!"
408              call abort_physic(subname,'allocation failure for wbarv',1)
409           endif
410         endif
411         if(.not.allocated(tauv)) then
412           ALLOCATE(tauv(L_NLEVRAD,L_NSPECTV,L_NGAUSS), stat=ok)
413           if (ok/=0) then
414              write(*,*) "memory allocation failed for tauv!"
415              call abort_physic(subname,'allocation failure for tauv',1)
416           endif
417         endif
418         if(.not.allocated(taucumv)) then
419           ALLOCATE(taucumv(L_LEVELS,L_NSPECTV,L_NGAUSS), stat=ok)
420           if (ok/=0) then
421              write(*,*) "memory allocation failed for taucumv!"
422              call abort_physic(subname,'allocation failure for taucumv',1)
423           endif
424         endif
425         if(.not.allocated(taucumi)) then
426           ALLOCATE(taucumi(L_LEVELS,L_NSPECTI,L_NGAUSS), stat=ok)
427           if (ok/=0) then
428              write(*,*) "memory allocation failed for taucumi!"
429              call abort_physic(subname,'allocation failure for taucumi',1)
430           endif
431         endif
432         if(.not.allocated(taugsurf)) then
433           ALLOCATE(taugsurf(L_NSPECTV,L_NGAUSS-1), stat=ok)
434           if (ok/=0) then
435              write(*,*) "memory allocation failed for taugsurf!"
436              call abort_physic(subname,'allocation failure for taugsurf',1)
437           endif
438         endif
439         if(.not.allocated(taugsurfi)) then
440           ALLOCATE(taugsurfi(L_NSPECTI,L_NGAUSS-1), stat=ok)
441           if (ok/=0) then
442              write(*,*) "memory allocation failed for taugsurfi!"
443              call abort_physic(subname,'allocation failure for taugsurfi',1)
444           endif
445         endif
446
447         if((igcm_h2o_vap.eq.0) .and. varactive .and. water)then
448            message='varactive in callcorrk but no h2o_vap tracer.'
449            call abort_physic(subname,message,1)
450         endif
451
452         if(varfixed .and. generic_condensation .and. .not. water)then
453            write(*,*) "Deep generic tracer vapor mixing ratio ? (no effect if negative) "
454            qvap_deep=-1. ! default value
455            call getin_p("qvap_deep",qvap_deep)
456            write(*,*) " qvap_deep = ",qvap_deep
457
458            metallicity=0.0 ! default value --- is not used here but necessary to call function Psat_generic
459            call getin_p("metallicity",metallicity) ! --- is not used here but necessary to call function Psat_generic
460         endif
461
462      end if ! of if (firstcall)
463
464!=======================================================================
465!          I.b  Initialization on every call   
466!=======================================================================
467 
468      qxvaer(:,:,:)=0.0
469      qsvaer(:,:,:)=0.0
470      gvaer(:,:,:) =0.0
471
472      qxiaer(:,:,:)=0.0
473      qsiaer(:,:,:)=0.0
474      giaer(:,:,:) =0.0
475
476      OLR_nu(:,:) = 0.
477      OSR_nu(:,:) = 0.
478      GSR_nu(:,:) = 0.
479
480!--------------------------------------------------
481!     Effective radius and variance of the aerosols
482!--------------------------------------------------
483
484      do iaer=1,naerkind
485
486         if ((iaer.eq.iaero_co2).and.tracer.and.(igcm_co2_ice.gt.0)) then ! Treat condensed co2 particles.
487            call co2_reffrad(ngrid,nlayer,nq,pq,reffrad(1,1,iaero_co2))
488
489            call planetwide_maxval(reffrad(:,:,iaero_co2),maxvalue)
490            call planetwide_minval(reffrad(:,:,iaero_co2),minvalue)
491            if (is_master) then
492               print*,'Max. CO2 ice particle size = ',maxvalue/1.e-6,' um'
493               print*,'Min. CO2 ice particle size = ',minvalue/1.e-6,' um'
494            end if
495         end if
496         
497         if ((iaer.eq.iaero_h2o).and.water) then ! Treat condensed water particles. To be generalized for other aerosols ...
498            call h2o_reffrad(ngrid,nlayer,pq(1,1,igcm_h2o_ice),pt, &
499                             reffrad(1,1,iaero_h2o),nueffrad(1,1,iaero_h2o))
500           
501            call planetwide_maxval(reffrad(:,:,iaero_h2o),maxvalue)
502            call planetwide_minval(reffrad(:,:,iaero_h2o),minvalue)
503            if (is_master) then
504               print*,'Max. H2O cloud particle size = ',maxvalue/1.e-6,' um'
505               print*,'Min. H2O cloud particle size = ',minvalue/1.e-6,' um'
506            end if
507
508! Currently the variance is constant everywhere (see h2o_reffrad),
509! so no need to compute and print min/max
510!            call planetwide_maxval(nueffrad(:,:,iaero_h2o),maxvalue)
511!            call planetwide_minval(nueffrad(:,:,iaero_h2o),minvalue)
512!            if (is_master) then
513!               print*,'Max. H2O cloud particle variance = ',maxvalue
514!               print*,'Min. H2O cloud particle variance = ',minvalue
515!            end if
516         endif
517         
518         if(iaer.eq.iaero_dust)then
519            call dust_reffrad(ngrid,nlayer,reffrad(1,1,iaero_dust))
520            if (is_master) then
521               print*,'Dust particle size = ',reffrad(1,1,iaer)/1.e-6,' um'
522            end if
523         endif
524         
525         if(iaer.eq.iaero_h2so4)then
526            call h2so4_reffrad(ngrid,nlayer,reffrad(1,1,iaero_h2so4))
527            if (is_master) then
528               print*,'H2SO4 particle size =',reffrad(1,1,iaer)/1.e-6,' um'
529            end if
530         endif
531         
532          if(iaer.eq.iaero_back2lay)then
533            call back2lay_reffrad(ngrid,reffrad(1,1,iaero_back2lay),nlayer,pplev)
534         endif
535
536         !  For n-layer aerosol size set once for all at firstcall in su_aer_radii
537
538!         if(iaer.eq.iaero_aurora)then
539!           call aurora_reffrad(ngrid,nlayer,reffrad(1,1,iaero_aurora))
540!         endif
541       
542     end do !iaer=1,naerkind.
543
544
545      ! How much light do we get ?
546      do nw=1,L_NSPECTV
547         stel(nw)=stellarf(nw)/(dist_star**2)
548      end do
549
550      ! Get 3D aerosol optical properties.
551      call aeroptproperties(ngrid,nlayer,reffrad,nueffrad,         &
552           QVISsQREF3d,omegaVIS3d,gVIS3d,                          &
553           QIRsQREF3d,omegaIR3d,gIR3d,                             &
554           QREFvis3d,QREFir3d)                                     
555
556      ! Get aerosol optical depths.
557      call aeropacity(ngrid,nlayer,nq,pplay,pplev,pt,pq,zls,aerosol,      &
558           reffrad,nueffrad,QREFvis3d,QREFir3d,                             &
559           tau_col,cloudfrac,totcloudfrac,clearsky)               
560 
561!-----------------------------------------------------------------------   
562      do ig=1,ngrid ! Starting Big Loop over every GCM column
563!-----------------------------------------------------------------------
564
565
566!=======================================================================
567!              II.  Transformation of the GCM variables
568!=======================================================================
569
570
571!-----------------------------------------------------------------------
572!    Aerosol optical properties Qext, Qscat and g.
573!    The transformation in the vertical is the same as for temperature.
574!-----------------------------------------------------------------------
575           
576           
577            do iaer=1,naerkind
578               ! Shortwave.
579               do nw=1,L_NSPECTV
580               
581                  do l=1,nlayer
582
583                     temp1=QVISsQREF3d(ig,nlayer+1-l,nw,iaer)         &
584                         *QREFvis3d(ig,nlayer+1-l,iaer)
585
586                     temp2=QVISsQREF3d(ig,max(nlayer-l,1),nw,iaer)    &
587                         *QREFvis3d(ig,max(nlayer-l,1),iaer)
588
589                     qxvaer(2*l,nw,iaer)  = temp1
590                     qxvaer(2*l+1,nw,iaer)=(temp1+temp2)/2
591
592                     temp1=temp1*omegavis3d(ig,nlayer+1-l,nw,iaer)
593                     temp2=temp2*omegavis3d(ig,max(nlayer-l,1),nw,iaer)
594
595                     qsvaer(2*l,nw,iaer)  = temp1
596                     qsvaer(2*l+1,nw,iaer)=(temp1+temp2)/2
597
598                     temp1=gvis3d(ig,nlayer+1-l,nw,iaer)
599                     temp2=gvis3d(ig,max(nlayer-l,1),nw,iaer)
600
601                     gvaer(2*l,nw,iaer)  = temp1
602                     gvaer(2*l+1,nw,iaer)=(temp1+temp2)/2
603
604                  end do ! nlayer
605
606                  qxvaer(1,nw,iaer)=qxvaer(2,nw,iaer)
607                  qxvaer(2*nlayer+1,nw,iaer)=0.
608
609                  qsvaer(1,nw,iaer)=qsvaer(2,nw,iaer)
610                  qsvaer(2*nlayer+1,nw,iaer)=0.
611
612                  gvaer(1,nw,iaer)=gvaer(2,nw,iaer)
613                  gvaer(2*nlayer+1,nw,iaer)=0.
614
615               end do ! L_NSPECTV
616             
617               do nw=1,L_NSPECTI
618                  ! Longwave
619                  do l=1,nlayer
620
621                     temp1=QIRsQREF3d(ig,nlayer+1-l,nw,iaer)         &
622                          *QREFir3d(ig,nlayer+1-l,iaer)
623
624                     temp2=QIRsQREF3d(ig,max(nlayer-l,1),nw,iaer)    &
625                          *QREFir3d(ig,max(nlayer-l,1),iaer)
626
627                     qxiaer(2*l,nw,iaer)  = temp1
628                     qxiaer(2*l+1,nw,iaer)=(temp1+temp2)/2
629
630                     temp1=temp1*omegair3d(ig,nlayer+1-l,nw,iaer)
631                     temp2=temp2*omegair3d(ig,max(nlayer-l,1),nw,iaer)
632
633                     qsiaer(2*l,nw,iaer)  = temp1
634                     qsiaer(2*l+1,nw,iaer)=(temp1+temp2)/2
635
636                     temp1=gir3d(ig,nlayer+1-l,nw,iaer)
637                     temp2=gir3d(ig,max(nlayer-l,1),nw,iaer)
638
639                     giaer(2*l,nw,iaer)  = temp1
640                     giaer(2*l+1,nw,iaer)=(temp1+temp2)/2
641
642                  end do ! nlayer
643
644                  qxiaer(1,nw,iaer)=qxiaer(2,nw,iaer)
645                  qxiaer(2*nlayer+1,nw,iaer)=0.
646
647                  qsiaer(1,nw,iaer)=qsiaer(2,nw,iaer)
648                  qsiaer(2*nlayer+1,nw,iaer)=0.
649
650                  giaer(1,nw,iaer)=giaer(2,nw,iaer)
651                  giaer(2*nlayer+1,nw,iaer)=0.
652
653               end do ! L_NSPECTI
654               
655            end do ! naerkind
656
657            ! Test / Correct for freaky s. s. albedo values.
658            do iaer=1,naerkind
659               do k=1,L_LEVELS
660
661                  do nw=1,L_NSPECTV
662                     if(qsvaer(k,nw,iaer).gt.1.05*qxvaer(k,nw,iaer))then
663                        message='Serious problems with qsvaer values'
664                        call abort_physic(subname,message,1)
665                     endif
666                     if(qsvaer(k,nw,iaer).gt.qxvaer(k,nw,iaer))then
667                        qsvaer(k,nw,iaer)=qxvaer(k,nw,iaer)
668                     endif
669                  end do
670
671                  do nw=1,L_NSPECTI
672                     if(qsiaer(k,nw,iaer).gt.1.05*qxiaer(k,nw,iaer))then
673                        message='Serious problems with qsvaer values'
674                        call abort_physic(subname,message,1)
675                     endif
676                     if(qsiaer(k,nw,iaer).gt.qxiaer(k,nw,iaer))then
677                        qsiaer(k,nw,iaer)=qxiaer(k,nw,iaer)
678                     endif
679                  end do
680
681               end do ! L_LEVELS
682            end do ! naerkind
683
684!-----------------------------------------------------------------------
685!     Aerosol optical depths
686!-----------------------------------------------------------------------
687           
688         do iaer=1,naerkind     ! a bug was here           
689            do k=0,nlayer-1
690               
691               pweight=(pplay(ig,L_NLAYRAD-k)-pplev(ig,L_NLAYRAD-k+1))/   &
692                       (pplev(ig,L_NLAYRAD-k)-pplev(ig,L_NLAYRAD-k+1))
693               ! As 'aerosol' is at reference (visible) wavelenght we scale it as
694               ! it will be multplied by qxi/v in optci/v
695               temp=aerosol(ig,L_NLAYRAD-k,iaer)/QREFvis3d(ig,L_NLAYRAD-k,iaer)
696               tauaero(2*k+2,iaer)=max(temp*pweight,0.d0)
697               tauaero(2*k+3,iaer)=max(temp-tauaero(2*k+2,iaer),0.d0)
698
699            end do
700            ! boundary conditions
701            tauaero(1,iaer)          = tauaero(2,iaer)
702            !tauaero(1,iaer)          = 0.
703            !JL18 at time of testing, the two above conditions gave the same results bit for bit.
704           
705         end do ! naerkind
706
707         ! Albedo and Emissivity.
708         albi=1-emis(ig)   ! Long Wave.
709         DO nw=1,L_NSPECTV ! Short Wave loop.
710            albv(nw)=albedo(ig,nw)
711         ENDDO
712
713         acosz=mu0(ig) ! Cosine of sun incident angle : 3D simulations or local 1D simulations using latitude.
714
715!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
716!!! Note by JL13 : In the following, some indices were changed in the interpolations,
717!!!                so that the model results are less dependent on the number of layers !
718!!!
719!!!           ---  The older versions are commented with the comment !JL13index  ---
720!!!
721!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
722
723
724      !-----------------------------------------------------------------------
725      !     Water vapour (to be generalised for other gases eventually ...)
726      !-----------------------------------------------------------------------
727           
728      if (water) then
729         if(varactive)then
730
731            i_var=igcm_h2o_vap
732            do l=1,nlayer
733               qvar(2*l)   = pq(ig,nlayer+1-l,i_var)
734               qvar(2*l+1) = pq(ig,nlayer+1-l,i_var)   
735               !JL13index   qvar(2*l+1) = (pq(ig,nlayer+1-l,i_var)+pq(ig,max(nlayer-l,1),i_var))/2   
736               !JL13index   Average approximation as for temperature...
737            end do
738            qvar(1)=qvar(2)
739
740         elseif(varfixed)then
741
742            do l=1,nlayer ! Here we will assign fixed water vapour profiles globally.
743               RH = satval * ((pplay(ig,l)/pplev(ig,1) - 0.02) / 0.98)
744               if(RH.lt.0.0) RH=0.0
745               
746               call Psat_water(pt(ig,l),pplay(ig,l),psat,qsat)
747
748               !pq_temp(l) = qsat      ! fully saturated everywhere
749               pq_temp(l) = RH * qsat ! ~realistic profile (e.g. 80% saturation at ground)
750            end do
751           
752            do l=1,nlayer
753               qvar(2*l)   = pq_temp(nlayer+1-l)
754               qvar(2*l+1) = (pq_temp(nlayer+1-l)+pq_temp(max(nlayer-l,1)))/2
755            end do
756           
757            qvar(1)=qvar(2)
758
759            ! Lowest layer of atmosphere
760            RH = satval * (1 - 0.02) / 0.98
761            if(RH.lt.0.0) RH=0.0
762
763            qvar(2*nlayer+1)= RH * qsat ! ~realistic profile (e.g. 80% saturation at ground)
764   
765         else
766            do k=1,L_LEVELS
767               qvar(k) = 1.0D-7
768            end do
769         end if ! varactive/varfixed
770     
771      endif ! if (water)
772
773      !-----------------------------------------------------------------------
774      !  GCS (Generic Condensable Specie) Vapor
775      !  If you have GCS tracers and they are : variable & radiatively active
776      !
777      !  NC22
778      !-----------------------------------------------------------------------
779
780      if (generic_condensation .and. .not. water ) then
781
782         ! For now, only one GCS tracer can be both variable and radiatively active
783         ! If you set two GCS tracers, that are variable and radiatively active,
784         ! the last one in tracer.def will be chosen as the one that will be vadiatively active
785
786         do iq=1,nq
787
788            call generic_tracer_index(nq,iq,igcm_generic_vap,igcm_generic_ice,call_ice_vap_generic)
789           
790            if (call_ice_vap_generic) then ! to call only one time the ice/vap pair of a tracer
791
792               if(varactive)then
793
794                  i_var=igcm_generic_vap
795                  do l=1,nlayer
796                     qvar(2*l)   = pq(ig,nlayer+1-l,i_var)
797                     qvar(2*l+1) = pq(ig,nlayer+1-l,i_var)   
798                     !JL13index            qvar(2*l+1) = (pq(ig,nlayer+1-l,i_var)+pq(ig,max(nlayer-l,1),i_var))/2   
799                     !JL13index            ! Average approximation as for temperature...
800                  end do
801                  qvar(1)=qvar(2)
802
803               elseif(varfixed .and. (qvap_deep .ge. 0))then
804
805                  do l=1,nlayer ! Here we will assign fixed water vapour profiles globally.
806                                         
807                     call Psat_generic(pt(ig,l),pplay(ig,l),metallicity,psat,qsat)
808
809                     if (qsat .lt. qvap_deep) then
810                        pq_temp(l) = qsat      ! fully saturated everywhere
811                     else
812                        pq_temp(l) = qvap_deep
813                     end if
814
815                  end do
816                 
817                  do l=1,nlayer
818                     qvar(2*l)   = pq_temp(nlayer+1-l)
819                     qvar(2*l+1) = (pq_temp(nlayer+1-l)+pq_temp(max(nlayer-l,1)))/2
820                  end do
821                 
822                  qvar(1)=qvar(2)
823         
824               else
825                  do k=1,L_LEVELS
826                     qvar(k) = 1.0D-7
827                  end do
828               end if ! varactive/varfixed
829
830            endif
831
832         end do ! do iq=1,nq loop on tracers
833
834      end if ! if (generic_condensation .and. .not. water )
835
836      !-----------------------------------------------------------------------
837      !  No Water vapor and No GCS (Generic Condensable Specie) vapor
838      !-----------------------------------------------------------------------
839
840      if (.not. generic_condensation .and. .not. water ) then
841         do k=1,L_LEVELS
842            qvar(k) = 1.0D-7
843         end do
844      end if ! if (.not. generic_condensation .and. .not. water )
845
846
847      if(.not.kastprof)then
848         ! IMPORTANT: Now convert from kg/kg to mol/mol.
849         do k=1,L_LEVELS
850            if (water) then
851               qvar(k) = qvar(k)/(epsi+qvar(k)*(1.-epsi))
852            endif
853            if (generic_condensation .and. .not. water) then
854               do iq=1,nq
855                  call generic_tracer_index(nq,iq,igcm_generic_vap,igcm_generic_ice,call_ice_vap_generic)
856                  if (call_ice_vap_generic) then ! to call only one time the ice/vap pair of a tracer
857
858                     epsi_generic=constants_epsi_generic(iq)
859
860                     qvar(k) = qvar(k)/(epsi_generic+qvar(k)*(1.-epsi_generic))
861
862                  endif
863               end do ! do iq=1,nq loop on tracers
864            endif
865         end do
866      end if
867
868!-----------------------------------------------------------------------
869!     kcm mode only !
870!-----------------------------------------------------------------------
871
872      if(kastprof)then
873     
874         if(.not.global1d)then ! garde-fou/safeguard added by MT (to be removed in the future)
875           message='You have to fix mu0, the cosinus of the solar angle'
876           call abort_physic(subname,message,1)
877         endif
878         
879         ! Initial values equivalent to mugaz.
880         DO l=1,nlayer
881            muvarrad(2*l)   = mugaz
882            muvarrad(2*l+1) = mugaz
883         END DO
884
885         if(ngasmx.gt.1)then
886
887            DO l=1,nlayer
888               muvarrad(2*l)   =  muvar(ig,nlayer+2-l)
889               muvarrad(2*l+1) = (muvar(ig,nlayer+2-l) + &
890                                  muvar(ig,max(nlayer+1-l,1)))/2
891            END DO
892     
893            muvarrad(1) = muvarrad(2)
894            muvarrad(2*nlayer+1) = muvar(ig,1)
895
896            print*,'Recalculating qvar with VARIABLE epsi for kastprof'
897            print*,'Assumes that the variable gas is H2O!!!'
898            print*,'Assumes that there is only one tracer'
899           
900            !i_var=igcm_h2o_vap
901            i_var=1
902           
903            if(nq.gt.1)then
904               message='Need 1 tracer only to run kcm1d.e'
905               call abort_physic(subname,message,1)
906            endif
907           
908            do l=1,nlayer
909               vtmp(l)=pq(ig,l,i_var)/(epsi+pq(ig,l,i_var)*(1.-epsi))
910               !vtmp(l)=pq(ig,l,i_var)*muvar(ig,l+1)/mH2O !JL to be changed
911            end do
912
913            do l=1,nlayer
914               qvar(2*l)   = vtmp(nlayer+1-l)
915               qvar(2*l+1) = vtmp(nlayer+1-l)
916!               qvar(2*l+1) = ( vtmp(nlayer+1-l) + vtmp(max(nlayer-l,1)) )/2
917            end do
918            qvar(1)=qvar(2)
919
920            write(*,*)trim(subname),' :Warning: reducing qvar in callcorrk.'
921            write(*,*)trim(subname),' :Temperature profile no longer consistent ', &
922                   'with saturated H2O. qsat=',satval
923                   
924            do k=1,L_LEVELS
925               qvar(k) = qvar(k)*satval
926            end do
927
928         endif
929      else ! if kastprof
930         DO l=1,nlayer
931            muvarrad(2*l)   = muvar(ig,nlayer+2-l)
932            muvarrad(2*l+1) = (muvar(ig,nlayer+2-l)+muvar(ig,max(nlayer+1-l,1)))/2
933         END DO
934     
935         muvarrad(1) = muvarrad(2)
936         muvarrad(2*nlayer+1)=muvar(ig,1)         
937      endif ! if kastprof
938     
939      ! Keep values inside limits for which we have radiative transfer coefficients !!!
940      if(L_REFVAR.gt.1)then ! (there was a bug here)
941         do k=1,L_LEVELS
942            if(qvar(k).lt.wrefvar(1))then
943               qvar(k)=wrefvar(1)+1.0e-8
944            elseif(qvar(k).gt.wrefvar(L_REFVAR))then
945               qvar(k)=wrefvar(L_REFVAR)-1.0e-8
946            endif
947         end do
948      endif
949
950!-----------------------------------------------------------------------
951!     Pressure and temperature
952!-----------------------------------------------------------------------
953
954      DO l=1,nlayer
955         plevrad(2*l)   = pplay(ig,nlayer+1-l)/scalep
956         plevrad(2*l+1) = pplev(ig,nlayer+1-l)/scalep
957         tlevrad(2*l)   = pt(ig,nlayer+1-l)
958         tlevrad(2*l+1) = (pt(ig,nlayer+1-l)+pt(ig,max(nlayer-l,1)))/2
959      END DO
960     
961      plevrad(1) = 0.
962!      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.
963      if (aerovenus) then
964!!  GG19 modified below after SL routines
965        plevrad(2) = 0.
966      endif
967
968      tlevrad(1) = tlevrad(2)
969      tlevrad(2*nlayer+1)=tsurf(ig)
970     
971      pmid(1) = pplay(ig,nlayer)/scalep   
972      if (aerovenus) then
973!! GG19 modified below after SL routines
974        pmid(1) = max(pgasmin,0.0001*plevrad(3))
975      endif
976      pmid(2) =  pmid(1)
977
978      tmid(1) = tlevrad(2)
979      tmid(2) = tmid(1)
980   
981      DO l=1,L_NLAYRAD-1
982         tmid(2*l+1) = tlevrad(2*l+1)
983         tmid(2*l+2) = tlevrad(2*l+1)
984         pmid(2*l+1) = plevrad(2*l+1)
985         pmid(2*l+2) = plevrad(2*l+1)
986      END DO
987      pmid(L_LEVELS) = plevrad(L_LEVELS)
988      tmid(L_LEVELS) = tlevrad(L_LEVELS)
989
990!!Alternative interpolation:
991!         pmid(3) = pmid(1)
992!         pmid(4) = pmid(1)
993!         tmid(3) = tmid(1)
994!         tmid(4) = tmid(1)
995!      DO l=2,L_NLAYRAD-1
996!         tmid(2*l+1) = tlevrad(2*l)
997!         tmid(2*l+2) = tlevrad(2*l)
998!         pmid(2*l+1) = plevrad(2*l)
999!         pmid(2*l+2) = plevrad(2*l)
1000!      END DO
1001!      pmid(L_LEVELS) = plevrad(L_LEVELS-1)
1002!      tmid(L_LEVELS) = tlevrad(L_LEVELS-1)
1003
1004      ! Test for out-of-bounds pressure.
1005      if(plevrad(3).lt.pgasmin)then
1006         print*,'Minimum pressure is outside the radiative'
1007         print*,'transfer kmatrix bounds, exiting.'
1008         message="Minimum pressure outside of kmatrix bounds"
1009         call abort_physic(subname,message,1)
1010      elseif(plevrad(L_LEVELS).gt.pgasmax)then
1011         print*,'Maximum pressure is outside the radiative'
1012         print*,'transfer kmatrix bounds, exiting.'
1013         message="Minimum pressure outside of kmatrix bounds"
1014         call abort_physic(subname,message,1)
1015      endif
1016
1017      ! Test for out-of-bounds temperature.
1018      ! -- JVO 20 : Also add a sanity test checking that tlevrad is
1019      !             within Planck function temperature boundaries,
1020      !             which would cause gfluxi/sfluxi to crash.
1021      do k=1,L_LEVELS
1022
1023         if(tlevrad(k).lt.tgasmin)then
1024            print*,'Minimum temperature is outside the radiative'
1025            print*,'transfer kmatrix bounds'
1026            print*,"k=",k," tlevrad(k)=",tlevrad(k)
1027            print*,"tgasmin=",tgasmin
1028            if (strictboundcorrk) then
1029              message="Minimum temperature outside of kmatrix bounds"
1030              call abort_physic(subname,message,1)
1031            else
1032              print*,'***********************************************'
1033              print*,'we allow model to continue with tlevrad<tgasmin'
1034              print*,'  ... we assume we know what you are doing ... '
1035              print*,'  ... but do not let this happen too often ... '
1036              print*,'***********************************************'
1037              !tlevrad(k)=tgasmin ! Used in the source function !
1038            endif
1039         elseif(tlevrad(k).gt.tgasmax)then
1040            print*,'Maximum temperature is outside the radiative'
1041            print*,'transfer kmatrix bounds, exiting.'
1042            print*,"k=",k," tlevrad(k)=",tlevrad(k)
1043            print*,"tgasmax=",tgasmax
1044            if (strictboundcorrk) then
1045              message="Maximum temperature outside of kmatrix bounds"
1046              call abort_physic(subname,message,1)
1047            else
1048              print*,'***********************************************'
1049              print*,'we allow model to continue with tlevrad>tgasmax' 
1050              print*,'  ... we assume we know what you are doing ... '
1051              print*,'  ... but do not let this happen too often ... '
1052              print*,'***********************************************'
1053              !tlevrad(k)=tgasmax ! Used in the source function !
1054            endif
1055         endif
1056
1057         if (tlevrad(k).lt.tplanckmin) then
1058            print*,'Minimum temperature is outside the boundaries for'
1059            print*,'Planck function integration set in callphys.def, aborting.'
1060            print*,"k=",k," tlevrad(k)=",tlevrad(k)
1061            print*,"tplanckmin=",tplanckmin
1062            message="Minimum temperature outside Planck function bounds - Change tplanckmin in callphys.def"
1063            call abort_physic(subname,message,1)
1064          else if (tlevrad(k).gt.tplanckmax) then
1065            print*,'Maximum temperature is outside the boundaries for'
1066            print*,'Planck function integration set in callphys.def, aborting.'
1067            print*,"k=",k," tlevrad(k)=",tlevrad(k)
1068            print*,"tplanckmax=",tplanckmax
1069            message="Maximum temperature outside Planck function bounds - Change tplanckmax in callphys.def"
1070            call abort_physic(subname,message,1)
1071          endif
1072
1073      enddo
1074
1075      do k=1,L_NLAYRAD+1
1076         if(tmid(k).lt.tgasmin)then
1077            print*,'Minimum temperature is outside the radiative'
1078            print*,'transfer kmatrix bounds, exiting.'
1079            print*,"k=",k," tmid(k)=",tmid(k)
1080            print*,"tgasmin=",tgasmin
1081            if (strictboundcorrk) then
1082              message="Minimum temperature outside of kmatrix bounds"
1083              call abort_physic(subname,message,1)
1084            else
1085              print*,'***********************************************'
1086              print*,'we allow model to continue but with tmid=tgasmin'
1087              print*,'  ... we assume we know what you are doing ... '
1088              print*,'  ... but do not let this happen too often ... '
1089              print*,'***********************************************'
1090              tmid(k)=tgasmin
1091            endif
1092         elseif(tmid(k).gt.tgasmax)then
1093            print*,'Maximum temperature is outside the radiative'
1094            print*,'transfer kmatrix bounds, exiting.'
1095            print*,"k=",k," tmid(k)=",tmid(k)
1096            print*,"tgasmax=",tgasmax
1097            if (strictboundcorrk) then
1098              message="Maximum temperature outside of kmatrix bounds"
1099              call abort_physic(subname,message,1)
1100            else
1101              print*,'***********************************************'
1102              print*,'we allow model to continue but with tmid=tgasmax'
1103              print*,'  ... we assume we know what you are doing ... '
1104              print*,'  ... but do not let this happen too often ... '
1105              print*,'***********************************************'
1106              tmid(k)=tgasmax
1107            endif
1108         endif
1109      enddo
1110     
1111!-----------------------------------------------------------------------
1112!     Variation of molar fraction for CIAs
1113!-----------------------------------------------------------------------
1114
1115      if (varspec) then
1116        do k=1,ngasmx
1117          call pindex(p_var,frac_var(:,k),plevrad(:),nvarlayer,L_LEVELS,frac_vari)
1118          fracvari(k,:) = frac_vari
1119        enddo
1120      endif
1121
1122!=======================================================================
1123!          III. Calling the main radiative transfer subroutines
1124!=======================================================================
1125
1126! ----------------------------------------------------------------
1127! Recombine reference corrk tables if needed - Added by JVO, 2020.
1128         if (corrk_recombin) then
1129           call call_recombin(ig,nlayer,pq(ig,:,:),pplay(ig,:),pt(ig,:),qvar(:),tmid(:),pmid(:))
1130         endif
1131! ----------------------------------------------------------------
1132
1133         Cmk= 0.01 * 1.0 / (glat(ig) * mugaz * 1.672621e-27) ! q_main=1.0 assumed.
1134         glat_ig=glat(ig)
1135
1136!-----------------------------------------------------------------------
1137!        Short Wave Part
1138!-----------------------------------------------------------------------
1139
1140         if((fract(ig) .ge. 1.0e-4).or.(global1d)) then ! Only during daylight.
1141            if((ngrid.eq.1).and.(global1d))then
1142               do nw=1,L_NSPECTV
1143                  stel_fract(nw)= stel(nw)* 0.25 / acosz ! globally averaged = divide by 4, and we correct for solar zenith angle
1144               end do
1145            else
1146               do nw=1,L_NSPECTV
1147                  stel_fract(nw)= stel(nw) * fract(ig)
1148               end do
1149            endif
1150
1151            call optcv(dtauv,tauv,taucumv,plevrad,                 &
1152                 qxvaer,qsvaer,gvaer,wbarv,cosbv,tauaero,          &
1153                 tmid,pmid,taugsurf,qvar,muvarrad,fracvari)
1154
1155            call sfluxv(dtauv,tauv,taucumv,albv,dwnv,wbarv,cosbv,  &
1156                 acosz,stel_fract,                                 &
1157                 nfluxtopv,fluxtopvdn,nfluxoutv_nu,nfluxgndv_nu,   &
1158                 fmnetv,fluxupv,fluxdnv,fzerov,taugsurf)
1159
1160         else ! During the night, fluxes = 0.
1161            nfluxtopv       = 0.0d0
1162            fluxtopvdn      = 0.0d0
1163            nfluxoutv_nu(:) = 0.0d0
1164            nfluxgndv_nu(:) = 0.0d0
1165            do l=1,L_NLAYRAD
1166               fmnetv(l)=0.0d0
1167               fluxupv(l)=0.0d0
1168               fluxdnv(l)=0.0d0
1169            end do
1170         end if
1171
1172
1173         ! Equivalent Albedo Calculation (for OUTPUT). MT2015
1174         if((fract(ig) .ge. 1.0e-4).or.(global1d)) then ! equivalent albedo makes sense only during daylight.       
1175            surface_stellar_flux=sum(nfluxgndv_nu(1:L_NSPECTV))     
1176            if(surface_stellar_flux .gt. 1.0e-3) then ! equivalent albedo makes sense only if the stellar flux received by the surface is positive.
1177               DO nw=1,L_NSPECTV                 
1178                  albedo_temp(nw)=albedo(ig,nw)*nfluxgndv_nu(nw)
1179               ENDDO
1180               albedo_temp(1:L_NSPECTV)=albedo_temp(1:L_NSPECTV)/surface_stellar_flux
1181               albedo_equivalent(ig)=sum(albedo_temp(1:L_NSPECTV))
1182            else
1183               albedo_equivalent(ig)=0.0 ! Spectrally Integrated Albedo not defined for non-irradiated grid points. So we arbitrary set the equivalent albedo to 0.
1184            endif
1185         else
1186            albedo_equivalent(ig)=0.0 ! Spectrally Integrated Albedo not defined for non-irradiated grid points. So we arbitrary set the equivalent albedo to 0.
1187         endif
1188
1189
1190!-----------------------------------------------------------------------
1191!        Long Wave Part
1192!-----------------------------------------------------------------------
1193
1194         call optci(plevrad,tlevrad,dtaui,taucumi,                  &
1195              qxiaer,qsiaer,giaer,cosbi,wbari,tauaero,tmid,pmid,    &
1196              taugsurfi,qvar,muvarrad,fracvari)
1197
1198         call sfluxi(plevrad,tlevrad,dtaui,taucumi,ubari,albi,      &
1199              wnoi,dwni,cosbi,wbari,nfluxtopi,nfluxtopi_nu,         &
1200              fmneti,fluxupi,fluxdni,fluxupi_nu,fzeroi,taugsurfi)
1201
1202!-----------------------------------------------------------------------
1203!     Transformation of the correlated-k code outputs
1204!     (into dtlw, dtsw, fluxsurf_lw, fluxsurf_sw, fluxtop_lw, fluxtop_sw)
1205
1206!     Flux incident at the top of the atmosphere
1207         fluxtop_dn(ig)=fluxtopvdn
1208
1209         fluxtop_lw(ig)  = real(nfluxtopi)
1210         fluxabs_sw(ig)  = real(-nfluxtopv)
1211         fluxsurf_lw(ig) = real(fluxdni(L_NLAYRAD))
1212         fluxsurf_sw(ig) = real(fluxdnv(L_NLAYRAD))
1213         
1214!        Flux absorbed by the surface. By MT2015.         
1215         fluxsurfabs_sw(ig) = fluxsurf_sw(ig)*(1.-albedo_equivalent(ig))
1216
1217         if(fluxtop_dn(ig).lt.0.0)then
1218            print*,'Achtung! fluxtop_dn has lost the plot!'
1219            print*,'fluxtop_dn=',fluxtop_dn(ig)
1220            print*,'acosz=',acosz
1221            print*,'aerosol=',aerosol(ig,:,:)
1222            print*,'temp=   ',pt(ig,:)
1223            print*,'pplay=  ',pplay(ig,:)
1224            message="Achtung! fluxtop_dn has lost the plot!"
1225            call abort_physic(subname,message,1)
1226         endif
1227
1228!     Spectral output, for exoplanet observational comparison
1229         if(specOLR)then
1230            do nw=1,L_NSPECTI
1231               OLR_nu(ig,nw)=nfluxtopi_nu(nw)/DWNI(nw) !JL Normalize to the bandwidth
1232            end do
1233            do nw=1,L_NSPECTV
1234               GSR_nu(ig,nw)=nfluxgndv_nu(nw)/DWNV(nw)
1235               OSR_nu(ig,nw)=nfluxoutv_nu(nw)/DWNV(nw) !JL Normalize to the bandwidth
1236            end do
1237         endif
1238
1239!     Finally, the heating rates
1240
1241         DO l=2,L_NLAYRAD
1242            dtsw(ig,L_NLAYRAD+1-l)=(fmnetv(l)-fmnetv(l-1))  &
1243                *glat(ig)/(cpp*scalep*(plevrad(2*l+1)-plevrad(2*l-1)))
1244            dtlw(ig,L_NLAYRAD+1-l)=(fmneti(l)-fmneti(l-1))  &
1245                *glat(ig)/(cpp*scalep*(plevrad(2*l+1)-plevrad(2*l-1)))
1246         END DO     
1247
1248!     These are values at top of atmosphere
1249         dtsw(ig,L_NLAYRAD)=(fmnetv(1)-nfluxtopv)           &
1250             *glat(ig)/(cpp*scalep*(plevrad(3)-plevrad(2)))
1251         dtlw(ig,L_NLAYRAD)=(fmneti(1)-nfluxtopi)           &
1252             *glat(ig)/(cpp*scalep*(plevrad(3)-plevrad(2)))
1253
1254      !  Optical thickness diagnostics (added by JVO)
1255      if (diagdtau) then
1256        do l=1,L_NLAYRAD
1257          do nw=1,L_NSPECTV
1258            int_dtauv(ig,l,nw) = 0.0d0
1259             DO k=1,L_NGAUSS
1260              ! Output exp(-tau) because gweight ponderates exp and not tau itself
1261              int_dtauv(ig,l,nw)= int_dtauv(ig,l,nw) + exp(-dtauv(l,nw,k))*gweight(k)
1262             ENDDO
1263          enddo
1264          do nw=1,L_NSPECTI
1265           int_dtaui(ig,l,nw) = 0.0d0
1266             DO k=1,L_NGAUSS
1267              ! Output exp(-tau) because gweight ponderates exp and not tau itself
1268              int_dtaui(ig,l,nw)= int_dtaui(ig,l,nw) + exp(-dtaui(l,nw,k))*gweight(k)
1269             ENDDO
1270          enddo
1271        enddo
1272      endif       
1273
1274
1275!-----------------------------------------------------------------------   
1276      end do ! End of big loop over every GCM column.
1277!-----------------------------------------------------------------------
1278
1279
1280
1281!-----------------------------------------------------------------------
1282!     Additional diagnostics
1283!-----------------------------------------------------------------------
1284
1285      ! IR spectral output, for exoplanet observational comparison
1286      if(lastcall.and.(ngrid.eq.1))then  ! could disable the 1D output, they are in the diagfi and diagspec... JL12
1287
1288         print*,'Saving scalar quantities in surf_vals.out...'
1289         print*,'psurf = ', pplev(1,1),' Pa'
1290         open(116,file='surf_vals.out')
1291         write(116,*) tsurf(1),pplev(1,1),fluxtop_dn(1),         &
1292                      real(-nfluxtopv),real(nfluxtopi)
1293         close(116)
1294
1295
1296!          USEFUL COMMENT - Do Not Remove.
1297!
1298!           if(specOLR)then
1299!               open(117,file='OLRnu.out')
1300!               do nw=1,L_NSPECTI
1301!                  write(117,*) OLR_nu(1,nw)
1302!               enddo
1303!               close(117)
1304!
1305!               open(127,file='OSRnu.out')
1306!               do nw=1,L_NSPECTV
1307!                  write(127,*) OSR_nu(1,nw)
1308!               enddo
1309!               close(127)
1310!           endif
1311
1312           ! OLR vs altitude: do it as a .txt file.
1313         OLRz=.false.
1314         if(OLRz)then
1315            print*,'saving IR vertical flux for OLRz...'
1316            open(118,file='OLRz_plevs.out')
1317            open(119,file='OLRz.out')
1318            do l=1,L_NLAYRAD
1319               write(118,*) plevrad(2*l)
1320               do nw=1,L_NSPECTI
1321                  write(119,*) fluxupi_nu(l,nw)
1322               enddo
1323            enddo
1324            close(118)
1325            close(119)
1326         endif
1327
1328      endif
1329
1330      ! See physiq.F for explanations about CLFvarying. This is temporary.
1331      if (lastcall .and. .not.CLFvarying) then
1332        IF( ALLOCATED( gasi ) ) DEALLOCATE( gasi )
1333        IF( ALLOCATED( gasv ) ) DEALLOCATE( gasv )
1334!$OMP BARRIER
1335!$OMP MASTER
1336        IF( ALLOCATED( pgasref ) ) DEALLOCATE( pgasref )
1337        IF( ALLOCATED( tgasref ) ) DEALLOCATE( tgasref )
1338        IF( ALLOCATED( wrefvar ) ) DEALLOCATE( wrefvar )
1339        IF( ALLOCATED( pfgasref ) ) DEALLOCATE( pfgasref )
1340        IF( ALLOCATED( gweight ) ) DEALLOCATE( gweight )
1341!$OMP END MASTER
1342!$OMP BARRIER
1343        IF ( ALLOCATED(reffrad)) DEALLOCATE(reffrad)
1344        IF ( ALLOCATED(nueffrad)) DEALLOCATE(nueffrad)
1345      endif
1346
1347
1348    end subroutine callcorrk
1349
1350END MODULE callcorrk_mod
Note: See TracBrowser for help on using the repository browser.