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

Last change on this file since 2782 was 2736, checked in by aslmd, 4 years ago

Changed of callcorrk to physiq_mod

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