source: trunk/LMDZ.GENERIC/libf/phystd/physiq_mod.F90 @ 1859

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

add diagnostic of reevaporation of precipitation

  • Property svn:executable set to *
File size: 91.4 KB
Line 
1      module physiq_mod
2     
3      implicit none
4     
5      contains
6     
7      subroutine physiq(ngrid,nlayer,nq,   &
8                  nametrac,                &
9                  firstcall,lastcall,      &
10                  pday,ptime,ptimestep,    &
11                  pplev,pplay,pphi,        &
12                  pu,pv,pt,pq,             &
13                  flxw,                    &
14                  pdu,pdv,pdt,pdq,pdpsrf)
15 
16      use radinc_h, only : L_NSPECTI,L_NSPECTV,naerkind
17      use watercommon_h, only : RLVTT, Psat_water,epsi
18      use gases_h, only: gnom, gfrac
19      use radcommon_h, only: sigma, glat, grav, BWNV
20      use radii_mod, only: h2o_reffrad, co2_reffrad
21      use aerosol_mod, only: iaero_co2, iaero_h2o
22      use surfdat_h, only: phisfi, zmea, zstd, zsig, zgam, zthe, &
23                           dryness, watercaptag
24      use comdiurn_h, only: coslat, sinlat, coslon, sinlon
25      use comsaison_h, only: mu0, fract, dist_star, declin, right_ascen
26      use comsoil_h, only: nsoilmx, layer, mlayer, inertiedat
27      use geometry_mod, only: latitude, longitude, cell_area
28      USE comgeomfi_h, only: totarea, totarea_planet
29      USE tracer_h, only: noms, mmol, radius, rho_q, qext, &
30                          alpha_lift, alpha_devil, qextrhor, &
31                          igcm_h2o_ice, igcm_h2o_vap, igcm_dustbin, &
32                          igcm_co2_ice
33      use time_phylmdz_mod, only: ecritphy, iphysiq, nday
34      use phyetat0_mod, only: phyetat0
35      use phyredem, only: physdem0, physdem1
36      use slab_ice_h, only: capcalocean, capcalseaice,capcalsno, &
37                            noceanmx
38      use ocean_slab_mod, only :ocean_slab_init, ocean_slab_ice, &
39                                ini_surf_heat_transp_mod, &
40                                ocean_slab_get_vars,ocean_slab_final
41      use surf_heat_transp_mod,only :init_masquv
42      use planetwide_mod, only: planetwide_minval,planetwide_maxval,planetwide_sumval
43      use mod_phys_lmdz_para, only : is_master
44      use planete_mod, only: apoastr, periastr, year_day, peri_day, &
45                            obliquit, nres, z0
46      use comcstfi_mod, only: pi, g, rcp, r, rad, mugaz, cpp
47      use time_phylmdz_mod, only: daysec
48      use callkeys_mod
49      use conc_mod
50      use phys_state_var_mod
51      use turb_mod, only : q2,sensibFlux,turb_resolved
52#ifndef MESOSCALE
53      use vertical_layers_mod, only: presnivs, pseudoalt
54      use mod_phys_lmdz_omp_data, ONLY: is_omp_master
55#else
56      use comm_wrf, only : comm_HR_SW, comm_HR_LW, &
57                           comm_CLOUDFRAC,comm_TOTCLOUDFRAC,&
58                           comm_RAIN,comm_SNOW,comm_ALBEQ,&
59                           comm_FLUXTOP_DN,comm_FLUXABS_SW,&
60                           comm_FLUXTOP_LW,comm_FLUXSURF_SW,&
61                           comm_FLUXSURF_LW,comm_FLXGRD,&
62                           comm_LSCEZ,comm_H2OICE_REFF
63#endif
64
65#ifdef CPP_XIOS     
66      use xios_output_mod, only: initialize_xios_output, &
67                                 update_xios_timestep, &
68                                 send_xios_field
69      use wxios, only: wxios_context_init, xios_context_finalize
70#endif
71
72      implicit none
73
74
75!==================================================================
76!     
77!     Purpose
78!     -------
79!     Central subroutine for all the physics parameterisations in the
80!     universal model. Originally adapted from the Mars LMDZ model.
81!
82!     The model can be run without or with tracer transport
83!     depending on the value of "tracer" in file "callphys.def".
84!
85!
86!   It includes:
87!
88!      I. Initialization :
89!         I.1 Firstcall initializations.
90!         I.2 Initialization for every call to physiq.
91!
92!      II. Compute radiative transfer tendencies (longwave and shortwave) :
93!         II.a Option 1 : Call correlated-k radiative transfer scheme.
94!         II.b Option 2 : Call Newtonian cooling scheme.
95!         II.c Option 3 : Atmosphere has no radiative effect.
96!
97!      III. Vertical diffusion (turbulent mixing) :
98!
99!      IV. Dry Convective adjustment :
100!
101!      V. Condensation and sublimation of gases (currently just CO2).
102!
103!      VI. Tracers
104!         VI.1. Water and water ice.
105!         VI.2  Photochemistry
106!         VI.3. Aerosols and particles.
107!         VI.4. Updates (pressure variations, surface budget).
108!         VI.5. Slab Ocean.
109!         VI.6. Surface Tracer Update.
110!
111!      VII. Surface and sub-surface soil temperature.
112!
113!      VIII. Perform diagnostics and write output files.
114!
115!
116!   arguments
117!   ---------
118!
119!   INPUT
120!   -----
121!
122!    ngrid                 Size of the horizontal grid.
123!    nlayer                Number of vertical layers.
124!    nq                    Number of advected fields.
125!    nametrac              Name of corresponding advected fields.
126!
127!    firstcall             True at the first call.
128!    lastcall              True at the last call.
129!
130!    pday                  Number of days counted from the North. Spring equinoxe.
131!    ptime                 Universal time (0<ptime<1): ptime=0.5 at 12:00 UT.
132!    ptimestep             timestep (s).
133!
134!    pplay(ngrid,nlayer)   Pressure at the middle of the layers (Pa).
135!    pplev(ngrid,nlayer+1) Intermediate pressure levels (Pa).
136!    pphi(ngrid,nlayer)    Geopotential at the middle of the layers (m2.s-2).
137!
138!    pu(ngrid,nlayer)      u, zonal component of the wind (ms-1).
139!    pv(ngrid,nlayer)      v, meridional component of the wind (ms-1).
140!
141!    pt(ngrid,nlayer)      Temperature (K).
142!
143!    pq(ngrid,nlayer,nq)   Advected fields.
144!
145!    pudyn(ngrid,nlayer)    \
146!    pvdyn(ngrid,nlayer)     \ Dynamical temporal derivative for the
147!    ptdyn(ngrid,nlayer)     / corresponding variables.
148!    pqdyn(ngrid,nlayer,nq) /
149!    flxw(ngrid,nlayer)      vertical mass flux (kg/s) at layer lower boundary
150!
151!   OUTPUT
152!   ------
153!
154!    pdu(ngrid,nlayer)        \
155!    pdv(ngrid,nlayer)         \  Temporal derivative of the corresponding
156!    pdt(ngrid,nlayer)         /  variables due to physical processes.
157!    pdq(ngrid,nlayer)        /
158!    pdpsrf(ngrid)             /
159!
160!
161!     Authors
162!     -------
163!           Frederic Hourdin        15/10/93
164!           Francois Forget        1994
165!           Christophe Hourdin        02/1997
166!           Subroutine completely rewritten by F. Forget (01/2000)
167!           Water ice clouds: Franck Montmessin (update 06/2003)
168!           Radiatively active tracers: J.-B. Madeleine (10/2008-06/2009)
169!           New correlated-k radiative scheme: R. Wordsworth (2009)
170!           Many specifically Martian subroutines removed: R. Wordsworth (2009)       
171!           Improved water cycle: R. Wordsworth / B. Charnay (2010)
172!           To F90: R. Wordsworth (2010)
173!           New turbulent diffusion scheme: J. Leconte (2012)
174!           Loops converted to F90 matrix format: J. Leconte (2012)
175!           No more ngridmx/nqmx, F90 commons and adaptation to parallel: A. Spiga (2012)
176!           Purge of the code : M. Turbet (2015)
177!           Photochemical core developped by F. Lefevre: B. Charnay (2017)
178!==================================================================
179
180
181!    0.  Declarations :
182!    ------------------
183
184include "netcdf.inc"
185
186! Arguments :
187! -----------
188
189!   INPUTS:
190!   -------
191
192      integer,intent(in) :: ngrid             ! Number of atmospheric columns.
193      integer,intent(in) :: nlayer            ! Number of atmospheric layers.
194      integer,intent(in) :: nq                ! Number of tracers.
195      character*20,intent(in) :: nametrac(nq) ! Names of the tracers taken from dynamics.
196     
197      logical,intent(in) :: firstcall ! Signals first call to physics.
198      logical,intent(in) :: lastcall  ! Signals last call to physics.
199     
200      real,intent(in) :: pday                  ! Number of elapsed sols since reference Ls=0.
201      real,intent(in) :: ptime                 ! "Universal time", given as fraction of sol (e.g.: 0.5 for noon).
202      real,intent(in) :: ptimestep             ! Physics timestep (s).
203      real,intent(in) :: pplev(ngrid,nlayer+1) ! Inter-layer pressure (Pa).
204      real,intent(in) :: pplay(ngrid,nlayer)   ! Mid-layer pressure (Pa).
205      real,intent(in) :: pphi(ngrid,nlayer)    ! Geopotential at mid-layer (m2s-2).
206      real,intent(in) :: pu(ngrid,nlayer)      ! Zonal wind component (m/s).
207      real,intent(in) :: pv(ngrid,nlayer)      ! Meridional wind component (m/s).
208      real,intent(in) :: pt(ngrid,nlayer)      ! Temperature (K).
209      real,intent(in) :: pq(ngrid,nlayer,nq)   ! Tracers (kg/kg_of_air).
210      real,intent(in) :: flxw(ngrid,nlayer)    ! Vertical mass flux (ks/s) at lower boundary of layer
211
212!   OUTPUTS:
213!   --------
214
215!     Physical tendencies :
216
217      real,intent(out) :: pdu(ngrid,nlayer)    ! Zonal wind tendencies (m/s/s).
218      real,intent(out) :: pdv(ngrid,nlayer)    ! Meridional wind tendencies (m/s/s).
219      real,intent(out) :: pdt(ngrid,nlayer)    ! Temperature tendencies (K/s).
220      real,intent(out) :: pdq(ngrid,nlayer,nq) ! Tracer tendencies (kg/kg_of_air/s).
221      real,intent(out) :: pdpsrf(ngrid)        ! Surface pressure tendency (Pa/s).
222
223! Local saved variables:
224! ----------------------
225      integer,save :: day_ini                                      ! Initial date of the run (sol since Ls=0).
226      integer,save :: icount                                       ! Counter of calls to physiq during the run.
227!$OMP THREADPRIVATE(day_ini,icount)
228
229! Local variables :
230! -----------------
231
232!     Aerosol (dust or ice) extinction optical depth  at reference wavelength
233!     for the "naerkind" optically active aerosols:
234
235      real aerosol(ngrid,nlayer,naerkind) ! Aerosols.
236      real zh(ngrid,nlayer)               ! Potential temperature (K).
237      real pw(ngrid,nlayer)               ! Vertical velocity (m/s). (NOTE : >0 WHEN DOWNWARDS !!)
238
239      integer l,ig,ierr,iq,nw,isoil
240     
241       
242        !$OMP zdtlw,zdtsw,sensibFlux)
243
244      real zls                       ! Solar longitude (radians).
245      real zlss                      ! Sub solar point longitude (radians).
246      real zday                      ! Date (time since Ls=0, calculated in sols).
247      real zzlay(ngrid,nlayer)       ! Altitude at the middle of the atmospheric layers.
248      real zzlev(ngrid,nlayer+1)     ! Altitude at the atmospheric layer boundaries.
249
250! TENDENCIES due to various processes :
251
252      ! For Surface Temperature : (K/s)
253      real zdtsurf(ngrid)     ! Cumulated tendencies.
254      real zdtsurfmr(ngrid)   ! Mass_redistribution routine.
255      real zdtsurfc(ngrid)    ! Condense_co2 routine.
256      real zdtsdif(ngrid)     ! Turbdiff/vdifc routines.
257      real zdtsurf_hyd(ngrid) ! Hydrol routine.
258           
259      ! For Atmospheric Temperatures : (K/s)   
260      real dtlscale(ngrid,nlayer)                             ! Largescale routine.
261      real zdtc(ngrid,nlayer)                                 ! Condense_co2 routine.
262      real zdtdif(ngrid,nlayer)                                      ! Turbdiff/vdifc routines.
263      real zdtmr(ngrid,nlayer)                                ! Mass_redistribution routine.
264      real zdtrain(ngrid,nlayer)                              ! Rain routine.
265      real dtmoist(ngrid,nlayer)                              ! Moistadj routine.
266      real dt_ekman(ngrid,noceanmx), dt_hdiff(ngrid,noceanmx) ! Slab_ocean routine.
267      real zdtsw1(ngrid,nlayer), zdtlw1(ngrid,nlayer)         ! Callcorrk routine.
268                             
269      ! For Surface Tracers : (kg/m2/s)
270      real dqsurf(ngrid,nq)                 ! Cumulated tendencies.
271      real zdqsurfc(ngrid)                  ! Condense_co2 routine.
272      real zdqsdif(ngrid,nq)                ! Turbdiff/vdifc routines.
273      real zdqssed(ngrid,nq)                ! Callsedim routine.
274      real zdqsurfmr(ngrid,nq)              ! Mass_redistribution routine.
275      real zdqsrain(ngrid), zdqssnow(ngrid) ! Rain routine.
276      real dqs_hyd(ngrid,nq)                ! Hydrol routine.
277      real reevap_precip(ngrid)             ! re-evaporation flux of precipitation (integrated over the atmospheric column)
278                 
279      ! For Tracers : (kg/kg_of_air/s)
280      real zdqc(ngrid,nlayer,nq)      ! Condense_co2 routine.
281      real zdqadj(ngrid,nlayer,nq)    ! Convadj routine.
282      real zdqdif(ngrid,nlayer,nq)    ! Turbdiff/vdifc routines.
283      real zdqevap(ngrid,nlayer)      ! Turbdiff routine.
284      real zdqsed(ngrid,nlayer,nq)    ! Callsedim routine.
285      real zdqmr(ngrid,nlayer,nq)     ! Mass_redistribution routine.
286      real zdqrain(ngrid,nlayer,nq)   ! Rain routine.
287      real dqmoist(ngrid,nlayer,nq)   ! Moistadj routine.
288      real dqvaplscale(ngrid,nlayer)  ! Largescale routine.
289      real dqcldlscale(ngrid,nlayer)  ! Largescale routine.
290      REAL zdqchim(ngrid,nlayer,nq)   ! Calchim_asis routine
291      REAL zdqschim(ngrid,nq)         ! Calchim_asis routine
292
293      REAL array_zero1(ngrid)
294      REAL array_zero2(ngrid,nlayer)
295                       
296      ! For Winds : (m/s/s)
297      real zdvadj(ngrid,nlayer),zduadj(ngrid,nlayer) ! Convadj routine.
298      real zdumr(ngrid,nlayer),zdvmr(ngrid,nlayer)   ! Mass_redistribution routine.
299      real zdvdif(ngrid,nlayer),zdudif(ngrid,nlayer) ! Turbdiff/vdifc routines.
300      real zdhdif(ngrid,nlayer)                      ! Turbdiff/vdifc routines.
301      real zdhadj(ngrid,nlayer)                      ! Convadj routine.
302     
303      ! For Pressure and Mass :
304      real zdmassmr(ngrid,nlayer) ! Atmospheric Mass tendency for mass_redistribution (kg_of_air/m2/s).
305      real zdmassmr_col(ngrid)    ! Atmospheric Column Mass tendency for mass_redistribution (kg_of_air/m2/s).
306      real zdpsrfmr(ngrid)        ! Pressure tendency for mass_redistribution routine (Pa/s).
307
308     
309   
310! Local variables for LOCAL CALCULATIONS:
311! ---------------------------------------
312      real zflubid(ngrid)
313      real zplanck(ngrid),zpopsk(ngrid,nlayer)
314      real ztim1,ztim2,ztim3, z1,z2
315      real ztime_fin
316      real zdh(ngrid,nlayer)
317      real gmplanet
318      real taux(ngrid),tauy(ngrid)
319
320
321! local variables for DIAGNOSTICS : (diagfi & stat)
322! -------------------------------------------------
323      real ps(ngrid)                                     ! Surface Pressure.
324      real zt(ngrid,nlayer)                              ! Atmospheric Temperature.
325      real zu(ngrid,nlayer),zv(ngrid,nlayer)             ! Zonal and Meridional Winds.
326      real zq(ngrid,nlayer,nq)                           ! Atmospheric Tracers.
327      real zdtadj(ngrid,nlayer)                          ! Convadj Diagnostic.
328      real zdtdyn(ngrid,nlayer)                          ! Dynamical Heating (K/s).
329      real zdudyn(ngrid,nlayer)                          ! Dynamical Zonal Wind tendency (m.s-2).
330
331      real reff(ngrid,nlayer)                       ! Effective dust radius (used if doubleq=T).
332      real vmr(ngrid,nlayer)                        ! volume mixing ratio
333      real time_phys
334     
335      real ISR,ASR,OLR,GND,DYN,GSR,Ts1,Ts2,Ts3,TsS ! for Diagnostic.
336     
337      real qcol(ngrid,nq) ! Tracer Column Mass (kg/m2).
338
339!     included by RW for H2O Manabe scheme
340      real rneb_man(ngrid,nlayer) ! H2O cloud fraction (moistadj).
341      real rneb_lsc(ngrid,nlayer) ! H2O cloud fraction (large scale).
342
343
344!     to test energy conservation (RW+JL)
345      real mass(ngrid,nlayer),massarea(ngrid,nlayer)
346      real dEtot, dEtots, AtmToSurf_TurbFlux
347      real,save :: dEtotSW, dEtotsSW, dEtotLW, dEtotsLW
348!$OMP THREADPRIVATE(dEtotSW, dEtotsSW, dEtotLW, dEtotsLW)
349      real dEzRadsw(ngrid,nlayer),dEzRadlw(ngrid,nlayer),dEzdiff(ngrid,nlayer)
350      real dEdiffs(ngrid),dEdiff(ngrid)
351      real madjdE(ngrid), lscaledE(ngrid),madjdEz(ngrid,nlayer), lscaledEz(ngrid,nlayer)
352     
353!JL12 conservation test for mean flow kinetic energy has been disabled temporarily
354
355      real dtmoist_max,dtmoist_min     
356      real dItot, dItot_tmp, dVtot, dVtot_tmp
357
358
359      real h2otot                      ! Total Amount of water. For diagnostic.
360      real icesrf,liqsrf,icecol,vapcol ! Total Amounts of water (ice,liq,vap). For diagnostic.
361      real dWtot, dWtot_tmp, dWtots, dWtots_tmp
362      logical,save :: watertest
363!$OMP THREADPRIVATE(watertest)
364
365      real qsat(ngrid,nlayer) ! Water Vapor Volume Mixing Ratio at saturation (kg/kg_of_air).
366      real RH(ngrid,nlayer)   ! Relative humidity.
367      real H2Omaxcol(ngrid)   ! Maximum possible H2O column amount (at 100% saturation) (kg/m2).
368      real psat_tmp
369     
370      logical clearsky ! For double radiative transfer call. By BC
371     
372      ! For Clear Sky Case.
373      real fluxsurf_lw1(ngrid), fluxsurf_sw1(ngrid), fluxsurfabs_sw1(ngrid)  ! For SW/LW flux.
374      real fluxtop_lw1(ngrid), fluxabs_sw1(ngrid)                            ! For SW/LW flux.
375      real albedo_equivalent1(ngrid)                                         ! For Equivalent albedo calculation.
376      real tau_col1(ngrid)                                                   ! For aerosol optical depth diagnostic.
377      real OLR_nu1(ngrid,L_NSPECTI), OSR_nu1(ngrid,L_NSPECTV)                ! For Outgoing Radiation diagnostics.
378      real tf, ntf
379
380      real nconsMAX, vdifcncons(ngrid), cadjncons(ngrid) ! Vdfic water conservation test. By RW
381
382      real muvar(ngrid,nlayer+1) ! For Runaway Greenhouse 1D study. By RW
383
384      real reffcol(ngrid,naerkind)
385
386!     Sourceevol for 'accelerated ice evolution'. By RW
387      real delta_ice,ice_tot
388      integer num_run
389      logical,save :: ice_update
390
391
392      real :: tsurf2(ngrid)
393      real :: flux_o(ngrid),flux_g(ngrid),fluxgrdocean(ngrid)
394      real :: flux_sens_lat(ngrid)
395      real :: qsurfint(ngrid,nq)
396
397
398!==================================================================================================
399
400! -----------------
401! I. INITIALISATION
402! -----------------
403
404! --------------------------------
405! I.1   First Call Initialisation.
406! --------------------------------
407      if (firstcall) then
408        ! Allocate saved arrays.
409        call phys_state_var_init(nq)
410
411!        Variables set to 0
412!        ~~~~~~~~~~~~~~~~~~
413         dtrad(:,:) = 0.0
414         fluxrad(:) = 0.0
415         tau_col(:) = 0.0
416         zdtsw(:,:) = 0.0
417         zdtlw(:,:) = 0.0
418
419
420!        Initialize aerosol indexes.
421!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~
422         call iniaerosol()
423
424     
425!        Initialize tracer names, indexes and properties.
426!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
427         IF (.NOT.ALLOCATED(noms)) ALLOCATE(noms(nq)) ! (because noms is an argument of physdem1 whether or not tracer is on)
428         if (tracer) then
429            call initracer(ngrid,nq,nametrac)
430            if(photochem) then
431              call ini_conc_mod(ngrid,nlayer)
432            endif
433         endif
434
435#ifdef CPP_XIOS
436        ! Initialize XIOS context
437        write(*,*) "physiq: call wxios_context_init"
438        CALL wxios_context_init
439#endif
440
441!        Read 'startfi.nc' file.
442!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
443#ifndef MESOSCALE
444         call phyetat0(startphy_file,                                 &
445                       ngrid,nlayer,"startfi.nc",0,0,nsoilmx,nq,      &
446                       day_ini,time_phys,tsurf,tsoil,emis,q2,qsurf,   &
447                       cloudfrac,totcloudfrac,hice,                   &
448                       rnat,pctsrf_sic,tslab, tsea_ice,sea_ice)
449#else
450      emis(:)=0.0
451      q2(:,:)=0.0
452      qsurf(:,:)=0.0
453      day_ini = pday
454#endif
455
456#ifndef MESOSCALE
457         if (.not.startphy_file) then
458           ! additionnal "academic" initialization of physics
459           if (is_master) write(*,*) "Physiq: initializing tsurf(:) to pt(:,1) !!"
460           tsurf(:)=pt(:,1)
461           if (is_master) write(*,*) "Physiq: initializing tsoil(:) to pt(:,1) !!"
462           do isoil=1,nsoilmx
463             tsoil(1:ngrid,isoil)=tsurf(1:ngrid)
464           enddo
465           if (is_master) write(*,*) "Physiq: initializing day_ini to pdat !"
466           day_ini=pday
467         endif
468#endif
469
470         if (pday.ne.day_ini) then
471           write(*,*) "ERROR in physiq.F90:"
472           write(*,*) "bad synchronization between physics and dynamics"
473           write(*,*) "dynamics day: ",pday
474           write(*,*) "physics day:  ",day_ini
475           stop
476         endif
477
478         write (*,*) 'In physiq day_ini =', day_ini
479
480!        Initialize albedo calculation.
481!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
482         albedo(:,:)=0.0
483          albedo_bareground(:)=0.0
484          albedo_snow_SPECTV(:)=0.0
485          albedo_co2_ice_SPECTV(:)=0.0
486         call surfini(ngrid,nq,qsurf,albedo,albedo_bareground,albedo_snow_SPECTV,albedo_co2_ice_SPECTV)
487         
488!        Initialize orbital calculation.
489!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
490         call iniorbit(apoastr,periastr,year_day,peri_day,obliquit)
491
492
493         if(tlocked)then
494            print*,'Planet is tidally locked at resonance n=',nres
495            print*,'Make sure you have the right rotation rate!!!'
496         endif
497
498!        Initialize oceanic variables.
499!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
500
501         if (ok_slab_ocean)then
502
503           call ocean_slab_init(ngrid,ptimestep, tslab, &
504                                sea_ice, pctsrf_sic)
505
506           call ini_surf_heat_transp_mod()
507           
508           knindex(:) = 0
509
510           do ig=1,ngrid
511              zmasq(ig)=1
512              knindex(ig) = ig
513              if (nint(rnat(ig)).eq.0) then
514                 zmasq(ig)=0
515              endif
516           enddo
517
518           CALL init_masquv(ngrid,zmasq)
519
520         endif ! end of 'ok_slab_ocean'.
521
522
523!        Initialize soil.
524!        ~~~~~~~~~~~~~~~~
525         if (callsoil) then
526         
527            call soil(ngrid,nsoilmx,firstcall,lastcall,inertiedat, &
528                      ptimestep,tsurf,tsoil,capcal,fluxgrd)
529
530            if (ok_slab_ocean) then
531               do ig=1,ngrid
532                  if (nint(rnat(ig)).eq.2) then
533                     capcal(ig)=capcalocean
534                     if (pctsrf_sic(ig).gt.0.5) then
535                        capcal(ig)=capcalseaice
536                        if (qsurf(ig,igcm_h2o_ice).gt.0.) then
537                           capcal(ig)=capcalsno
538                        endif
539                     endif
540                  endif
541               enddo
542            endif ! end of 'ok_slab_ocean'.
543
544         else ! else of 'callsoil'.
545         
546            print*,'WARNING! Thermal conduction in the soil turned off'
547            capcal(:)=1.e6
548            fluxgrd(:)=intheat
549            print*,'Flux from ground = ',intheat,' W m^-2'
550           
551         endif ! end of 'callsoil'.
552         
553         icount=1
554
555!        Decide whether to update ice at end of run.
556!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
557         ice_update=.false.
558         if(sourceevol)then
559!$OMP MASTER
560            open(128,file='num_run',form='formatted', &
561                     status="old",iostat=ierr)
562            if (ierr.ne.0) then
563               write(*,*) "physiq: Error! No num_run file!"
564               write(*,*) " (which is needed for sourceevol option)"
565               stop
566            endif
567            read(128,*) num_run
568            close(128)
569!$OMP END MASTER
570!$OMP BARRIER
571       
572            if(num_run.ne.0.and.mod(num_run,2).eq.0)then
573               print*,'Updating ice at end of this year!'
574               ice_update=.true.
575               ice_min(:)=1.e4
576            endif
577           
578         endif ! end of 'sourceevol'.
579
580
581         ! Here is defined the type of the surface : Continent or Ocean.
582         ! BC2014 : This is now already done in newstart.
583         ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
584         if (.not.ok_slab_ocean) then
585         
586           rnat(:)=1.
587           do ig=1,ngrid
588              if(inertiedat(ig,1).gt.1.E4)then
589                 rnat(ig)=0
590              endif
591           enddo
592
593           print*,'WARNING! Surface type currently decided by surface inertia'
594           print*,'This should be improved e.g. in newstart.F'
595           
596         endif ! end of 'ok_slab_ocean'.
597
598
599!        Initialize surface history variable.
600!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
601         qsurf_hist(:,:)=qsurf(:,:)
602
603!        Initialize variable for dynamical heating and zonal wind tendency diagnostic
604!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
605         ztprevious(:,:)=pt(:,:)
606         zuprevious(:,:)=pu(:,:)
607
608!        Set temperature just above condensation temperature (for Early Mars)
609!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
610         if(nearco2cond) then
611            write(*,*)' WARNING! Starting at Tcond+1K'
612            do l=1, nlayer
613               do ig=1,ngrid
614                  pdt(ig,l)= ((-3167.8)/(log(.01*pplay(ig,l))-23.23)+4     &
615                      -pt(ig,l)) / ptimestep
616               enddo
617            enddo
618         endif
619
620         if(meanOLR)then         
621            call system('rm -f rad_bal.out') ! to record global radiative balance.           
622            call system('rm -f tem_bal.out') ! to record global mean/max/min temperatures.           
623            call system('rm -f h2o_bal.out') ! to record global hydrological balance.
624         endif
625
626
627         watertest=.false.       
628         if(water)then ! Initialize variables for water cycle
629           
630            if(enertest)then
631               watertest = .true.
632            endif
633
634            if(ice_update)then
635               ice_initial(:)=qsurf(:,igcm_h2o_ice)
636            endif
637
638         endif
639         
640         call su_watercycle ! even if we don't have a water cycle, we might
641                            ! need epsi for the wvp definitions in callcorrk.F
642#ifndef MESOSCALE
643         if (ngrid.ne.1) then ! Note : no need to create a restart file in 1d.
644            call physdem0("restartfi.nc",longitude,latitude,nsoilmx,ngrid,nlayer,nq, &
645                          ptimestep,pday+nday,time_phys,cell_area,          &
646                          albedo_bareground,inertiedat,zmea,zstd,zsig,zgam,zthe)
647         endif
648#endif         
649         
650         ! XIOS outputs
651#ifdef CPP_XIOS
652
653         write(*,*) "physiq: call initialize_xios_output"
654         call initialize_xios_output(pday,ptime,ptimestep,daysec, &
655                                     presnivs,pseudoalt)
656#endif
657         write(*,*) "physiq: end of firstcall"
658      endif ! end of 'firstcall'
659
660! ------------------------------------------------------
661! I.2   Initializations done at every physical timestep:
662! ------------------------------------------------------
663
664#ifdef CPP_XIOS     
665      ! update XIOS time/calendar
666      call update_xios_timestep
667#endif     
668
669      ! Initialize various variables
670      ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
671     
672      if ( .not.nearco2cond ) then
673         pdt(1:ngrid,1:nlayer) = 0.0
674      endif     
675      zdtsurf(1:ngrid)      = 0.0
676      pdq(1:ngrid,1:nlayer,1:nq) = 0.0
677      dqsurf(1:ngrid,1:nq)= 0.0
678      pdu(1:ngrid,1:nlayer) = 0.0
679      pdv(1:ngrid,1:nlayer) = 0.0
680      pdpsrf(1:ngrid)       = 0.0
681      zflubid(1:ngrid)      = 0.0     
682      flux_sens_lat(1:ngrid) = 0.0
683      taux(1:ngrid) = 0.0
684      tauy(1:ngrid) = 0.0
685
686      zday=pday+ptime ! Compute time, in sols (and fraction thereof).
687
688      ! Compute Stellar Longitude (Ls), and orbital parameters.
689      ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
690      if (season) then
691         call stellarlong(zday,zls)
692      else
693         call stellarlong(float(day_ini),zls)
694      end if
695
696      call orbite(zls,dist_star,declin,right_ascen)
697     
698      if (tlocked) then
699              zlss=Mod(-(2.*pi*(zday/year_day)*nres - right_ascen),2.*pi)
700      elseif (diurnal) then
701              zlss=-2.*pi*(zday-.5)
702      else if(diurnal .eqv. .false.) then
703              zlss=9999.
704      endif
705
706
707      ! Compute variations of g with latitude (oblate case).
708      ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
709      if (oblate .eqv. .false.) then     
710         glat(:) = g         
711      else if (flatten .eq. 0.0 .or. J2 .eq. 0.0 .or. Rmean .eq. 0.0 .or. MassPlanet .eq. 0.0) then     
712         print*,'I need values for flatten, J2, Rmean and MassPlanet to compute glat (else set oblate=.false.)'
713         call abort
714      else
715         gmplanet = MassPlanet*grav*1e24
716         do ig=1,ngrid
717            glat(ig)= gmplanet/(Rmean**2) * (1.D0 + 0.75 *J2 - 2.0*flatten/3. + (2.*flatten - 15./4.* J2) * cos(2. * (pi/2. - latitude(ig))))
718         end do
719      endif
720
721      ! Compute geopotential between layers.
722      ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
723      zzlay(1:ngrid,1:nlayer)=pphi(1:ngrid,1:nlayer)
724      do l=1,nlayer         
725         zzlay(1:ngrid,l)= zzlay(1:ngrid,l)/glat(1:ngrid)
726      enddo
727
728      zzlev(1:ngrid,1)=0.
729      zzlev(1:ngrid,nlayer+1)=1.e7 ! Dummy top of last layer above 10000 km...
730
731      do l=2,nlayer
732         do ig=1,ngrid
733            z1=(pplay(ig,l-1)+pplev(ig,l))/(pplay(ig,l-1)-pplev(ig,l))
734            z2=(pplev(ig,l)+pplay(ig,l))/(pplev(ig,l)-pplay(ig,l))
735            zzlev(ig,l)=(z1*zzlay(ig,l-1)+z2*zzlay(ig,l))/(z1+z2)
736         enddo
737      enddo     
738
739      ! Compute potential temperature
740      ! Note : Potential temperature calculation may not be the same in physiq and dynamic...
741      ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
742      do l=1,nlayer         
743         do ig=1,ngrid
744            zpopsk(ig,l)=(pplay(ig,l)/pplev(ig,1))**rcp
745            zh(ig,l)=pt(ig,l)/zpopsk(ig,l)
746            mass(ig,l)  = (pplev(ig,l) - pplev(ig,l+1))/glat(ig)
747            massarea(ig,l)=mass(ig,l)*cell_area(ig)
748         enddo
749      enddo
750
751     ! Compute vertical velocity (m/s) from vertical mass flux
752     ! w = F / (rho*area) and rho = P/(r*T)
753     ! But first linearly interpolate mass flux to mid-layers
754      do l=1,nlayer-1
755         pw(1:ngrid,l)=0.5*(flxw(1:ngrid,l)+flxw(1:ngrid,l+1))
756      enddo
757      pw(1:ngrid,nlayer)=0.5*flxw(1:ngrid,nlayer) ! since flxw(nlayer+1)=0
758      do l=1,nlayer
759         pw(1:ngrid,l)=(pw(1:ngrid,l)*r*pt(1:ngrid,l)) /  &
760                       (pplay(1:ngrid,l)*cell_area(1:ngrid))
761      enddo
762
763      ! ----------------------------------------------------------------
764      ! Compute mean mass, cp, and R
765      ! --------------------------------
766
767      if(photochem) then
768         call concentrations(ngrid,nlayer,nq,pplay,pt,pdt,pq,pdq,ptimestep)
769      endif
770
771!---------------------------------
772! II. Compute radiative tendencies
773!---------------------------------
774
775      if (callrad) then
776         if( mod(icount-1,iradia).eq.0.or.lastcall) then
777
778          ! Compute local stellar zenith angles
779          ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
780            if (tlocked) then
781            ! JL14 corrects tidally resonant (and inclined) cases. nres=omega_rot/omega_orb
782               ztim1=SIN(declin)
783               ztim2=COS(declin)*COS(zlss)
784               ztim3=COS(declin)*SIN(zlss)
785
786               call stelang(ngrid,sinlon,coslon,sinlat,coslat,    &
787                            ztim1,ztim2,ztim3,mu0,fract, flatten)
788
789            elseif (diurnal) then
790               ztim1=SIN(declin)
791               ztim2=COS(declin)*COS(2.*pi*(zday-.5))
792               ztim3=-COS(declin)*SIN(2.*pi*(zday-.5))
793
794               call stelang(ngrid,sinlon,coslon,sinlat,coslat,    &
795                            ztim1,ztim2,ztim3,mu0,fract, flatten)
796            else if(diurnal .eqv. .false.) then
797
798               call mucorr(ngrid,declin,latitude,mu0,fract,10000.,rad,flatten)
799               ! WARNING: this function appears not to work in 1D
800
801            endif
802           
803            ! Eclipse incoming sunlight (e.g. Saturn ring shadowing).       
804            if(rings_shadow) then
805                call call_rings(ngrid, ptime, pday, diurnal)
806            endif   
807
808
809            if (corrk) then
810           
811! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
812! II.a Call correlated-k radiative transfer scheme
813! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
814               if(kastprof)then
815                  print*,'kastprof should not = true here'
816                  call abort
817               endif
818               if(water) then
819                  muvar(1:ngrid,1:nlayer)=mugaz/(1.e0+(1.e0/epsi-1.e0)*pq(1:ngrid,1:nlayer,igcm_h2o_vap))
820                  muvar(1:ngrid,nlayer+1)=mugaz/(1.e0+(1.e0/epsi-1.e0)*pq(1:ngrid,nlayer,igcm_h2o_vap))
821                  ! take into account water effect on mean molecular weight
822               else
823                  muvar(1:ngrid,1:nlayer+1)=mugaz
824               endif
825     
826
827               if(ok_slab_ocean) then
828                  tsurf2(:)=tsurf(:)
829                  do ig=1,ngrid
830                     if (nint(rnat(ig))==0) then
831                        tsurf(ig)=((1.-pctsrf_sic(ig))*tslab(ig,1)**4+pctsrf_sic(ig)*tsea_ice(ig)**4)**0.25
832                     endif
833                  enddo
834               endif !(ok_slab_ocean)
835               
836               ! standard callcorrk
837               clearsky=.false.
838               call callcorrk(ngrid,nlayer,pq,nq,qsurf,                           &
839                              albedo,albedo_equivalent,emis,mu0,pplev,pplay,pt,   &
840                              tsurf,fract,dist_star,aerosol,muvar,                &
841                              zdtlw,zdtsw,fluxsurf_lw,fluxsurf_sw,                &
842                              fluxsurfabs_sw,fluxtop_lw,                          &
843                              fluxabs_sw,fluxtop_dn,OLR_nu,OSR_nu,                &
844                              tau_col,cloudfrac,totcloudfrac,                     &
845                              clearsky,firstcall,lastcall)     
846
847                ! Option to call scheme once more for clear regions 
848               if(CLFvarying)then
849
850                  ! ---> PROBLEMS WITH ALLOCATED ARRAYS : temporary solution in callcorrk: do not deallocate if CLFvarying ...
851                  clearsky=.true.
852                  call callcorrk(ngrid,nlayer,pq,nq,qsurf,                           &
853                                 albedo,albedo_equivalent1,emis,mu0,pplev,pplay,pt,   &
854                                 tsurf,fract,dist_star,aerosol,muvar,                &
855                                 zdtlw1,zdtsw1,fluxsurf_lw1,fluxsurf_sw1,            &
856                                 fluxsurfabs_sw1,fluxtop_lw1,                        &
857                                 fluxabs_sw1,fluxtop_dn,OLR_nu1,OSR_nu1,             &
858                                 tau_col1,cloudfrac,totcloudfrac,                    &
859                                 clearsky,firstcall,lastcall)
860                  clearsky = .false.  ! just in case.     
861
862                  ! Sum the fluxes and heating rates from cloudy/clear cases
863                  do ig=1,ngrid
864                     tf=totcloudfrac(ig)
865                     ntf=1.-tf                   
866                     fluxsurf_lw(ig)       = ntf*fluxsurf_lw1(ig)       + tf*fluxsurf_lw(ig)
867                     fluxsurf_sw(ig)       = ntf*fluxsurf_sw1(ig)       + tf*fluxsurf_sw(ig)
868                     albedo_equivalent(ig) = ntf*albedo_equivalent1(ig) + tf*albedo_equivalent(ig)
869                     fluxsurfabs_sw(ig)    = ntf*fluxsurfabs_sw1(ig)    + tf*fluxsurfabs_sw(ig)
870                     fluxtop_lw(ig)        = ntf*fluxtop_lw1(ig)        + tf*fluxtop_lw(ig)
871                     fluxabs_sw(ig)        = ntf*fluxabs_sw1(ig)        + tf*fluxabs_sw(ig)
872                     tau_col(ig)           = ntf*tau_col1(ig)           + tf*tau_col(ig)
873                   
874                     zdtlw(ig,1:nlayer) = ntf*zdtlw1(ig,1:nlayer) + tf*zdtlw(ig,1:nlayer)
875                     zdtsw(ig,1:nlayer) = ntf*zdtsw1(ig,1:nlayer) + tf*zdtsw(ig,1:nlayer)
876
877                     OSR_nu(ig,1:L_NSPECTV) = ntf*OSR_nu1(ig,1:L_NSPECTV) + tf*OSR_nu(ig,1:L_NSPECTV)                       
878                     OLR_nu(ig,1:L_NSPECTI) = ntf*OLR_nu1(ig,1:L_NSPECTI) + tf*OLR_nu(ig,1:L_NSPECTI)                       
879                  enddo                               
880
881               endif ! end of CLFvarying.
882
883               if(ok_slab_ocean) then
884                  tsurf(:)=tsurf2(:)
885               endif
886
887
888               ! Radiative flux from the sky absorbed by the surface (W.m-2).
889               GSR=0.0
890               fluxrad_sky(1:ngrid)=emis(1:ngrid)*fluxsurf_lw(1:ngrid)+fluxsurfabs_sw(1:ngrid)
891
892                            !if(noradsurf)then ! no lower surface; SW flux just disappears
893                            !   GSR = SUM(fluxsurf_sw(1:ngrid)*cell_area(1:ngrid))/totarea
894                            !   fluxrad_sky(1:ngrid)=emis(1:ngrid)*fluxsurf_lw(1:ngrid)
895                            !   print*,'SW lost in deep atmosphere = ',GSR,' W m^-2'
896                            !endif
897
898               ! Net atmospheric radiative heating rate (K.s-1)
899               dtrad(1:ngrid,1:nlayer)=zdtsw(1:ngrid,1:nlayer)+zdtlw(1:ngrid,1:nlayer)
900               
901               ! Late initialization of the Ice Spectral Albedo. We needed the visible bands to do that !
902               if (firstcall .and. albedo_spectral_mode) then
903                  call spectral_albedo_calc(albedo_snow_SPECTV)
904               endif
905
906            elseif(newtonian)then
907           
908! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
909! II.b Call Newtonian cooling scheme
910! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
911               call newtrelax(ngrid,nlayer,mu0,sinlat,zpopsk,pt,pplay,pplev,dtrad,firstcall)
912
913               zdtsurf(1:ngrid) = +(pt(1:ngrid,1)-tsurf(1:ngrid))/ptimestep
914               ! e.g. surface becomes proxy for 1st atmospheric layer ?
915
916            else
917           
918! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
919! II.c Atmosphere has no radiative effect
920! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
921               fluxtop_dn(1:ngrid)  = fract(1:ngrid)*mu0(1:ngrid)*Fat1AU/dist_star**2
922               if(ngrid.eq.1)then ! / by 4 globally in 1D case...
923                  fluxtop_dn(1)  = fract(1)*Fat1AU/dist_star**2/2.0
924               endif
925               fluxsurf_sw(1:ngrid) = fluxtop_dn(1:ngrid)
926               print*,'------------WARNING---WARNING------------' ! by MT2015.
927               print*,'You are in corrk=false mode, '
928               print*,'and the surface albedo is taken equal to the first visible spectral value'               
929               
930               fluxsurfabs_sw(1:ngrid) = fluxtop_dn(1:ngrid)*(1.-albedo(1:ngrid,1))
931               fluxrad_sky(1:ngrid)    = fluxsurfabs_sw(1:ngrid)
932               fluxtop_lw(1:ngrid)  = emis(1:ngrid)*sigma*tsurf(1:ngrid)**4
933
934               dtrad(1:ngrid,1:nlayer)=0.0 ! no atmospheric radiative heating
935
936            endif ! end of corrk
937
938         endif ! of if(mod(icount-1,iradia).eq.0)
939       
940
941         ! Transformation of the radiative tendencies
942         ! ------------------------------------------
943         zplanck(1:ngrid)=tsurf(1:ngrid)*tsurf(1:ngrid)
944         zplanck(1:ngrid)=emis(1:ngrid)*sigma*zplanck(1:ngrid)*zplanck(1:ngrid)
945         fluxrad(1:ngrid)=fluxrad_sky(1:ngrid)-zplanck(1:ngrid)
946         pdt(1:ngrid,1:nlayer)=pdt(1:ngrid,1:nlayer)+dtrad(1:ngrid,1:nlayer)
947         
948         ! Test of energy conservation
949         !----------------------------
950         if(enertest)then
951            call planetwide_sumval(cpp*massarea(:,:)*zdtsw(:,:)/totarea_planet,dEtotSW)
952            call planetwide_sumval(cpp*massarea(:,:)*zdtlw(:,:)/totarea_planet,dEtotLW)
953            !call planetwide_sumval(fluxsurf_sw(:)*(1.-albedo_equivalent(:))*cell_area(:)/totarea_planet,dEtotsSW) !JL13 carefull, albedo can have changed since the last time we called corrk
954            call planetwide_sumval(fluxsurfabs_sw(:)*cell_area(:)/totarea_planet,dEtotsSW) !JL13 carefull, albedo can have changed since the last time we called corrk
955            call planetwide_sumval((fluxsurf_lw(:)*emis(:)-zplanck(:))*cell_area(:)/totarea_planet,dEtotsLW)
956            dEzRadsw(:,:)=cpp*mass(:,:)*zdtsw(:,:)
957            dEzRadlw(:,:)=cpp*mass(:,:)*zdtlw(:,:)
958            if (is_master) then
959               print*,'---------------------------------------------------------------'
960               print*,'In corrk SW atmospheric heating       =',dEtotSW,' W m-2'
961               print*,'In corrk LW atmospheric heating       =',dEtotLW,' W m-2'
962               print*,'atmospheric net rad heating (SW+LW)   =',dEtotLW+dEtotSW,' W m-2'
963               print*,'In corrk SW surface heating           =',dEtotsSW,' W m-2'
964               print*,'In corrk LW surface heating           =',dEtotsLW,' W m-2'
965               print*,'surface net rad heating (SW+LW)       =',dEtotsLW+dEtotsSW,' W m-2'
966            endif
967         endif ! end of 'enertest'
968
969      endif ! of if (callrad)
970
971
972
973!  --------------------------------------------
974!  III. Vertical diffusion (turbulent mixing) :
975!  --------------------------------------------
976
977      if (calldifv) then
978     
979         zflubid(1:ngrid)=fluxrad(1:ngrid)+fluxgrd(1:ngrid)
980
981         ! JL12 the following if test is temporarily there to allow us to compare the old vdifc with turbdiff.
982         if (UseTurbDiff) then
983         
984            call turbdiff(ngrid,nlayer,nq,rnat,                  &
985                          ptimestep,capcal,lwrite,               &
986                          pplay,pplev,zzlay,zzlev,z0,            &
987                          pu,pv,pt,zpopsk,pq,tsurf,emis,qsurf,   &
988                          pdt,pdq,zflubid,                       &
989                          zdudif,zdvdif,zdtdif,zdtsdif,          &
990                          sensibFlux,q2,zdqdif,zdqevap,zdqsdif,  &
991                          taux,tauy,lastcall)
992
993         else
994         
995            zdh(1:ngrid,1:nlayer)=pdt(1:ngrid,1:nlayer)/zpopsk(1:ngrid,1:nlayer)
996 
997            call vdifc(ngrid,nlayer,nq,rnat,zpopsk,           &
998                       ptimestep,capcal,lwrite,               &
999                       pplay,pplev,zzlay,zzlev,z0,            &
1000                       pu,pv,zh,pq,tsurf,emis,qsurf,          &
1001                       zdh,pdq,zflubid,                       &
1002                       zdudif,zdvdif,zdhdif,zdtsdif,          &
1003                       sensibFlux,q2,zdqdif,zdqsdif,          &
1004                       taux,tauy,lastcall)
1005
1006            zdtdif(1:ngrid,1:nlayer)=zdhdif(1:ngrid,1:nlayer)*zpopsk(1:ngrid,1:nlayer) ! for diagnostic only
1007            zdqevap(1:ngrid,1:nlayer)=0.
1008
1009         end if !end of 'UseTurbDiff'
1010
1011         !!! this is always done, except for turbulence-resolving simulations
1012         if (.not. turb_resolved) then
1013           pdv(1:ngrid,1:nlayer)=pdv(1:ngrid,1:nlayer)+zdvdif(1:ngrid,1:nlayer)
1014           pdu(1:ngrid,1:nlayer)=pdu(1:ngrid,1:nlayer)+zdudif(1:ngrid,1:nlayer)
1015           pdt(1:ngrid,1:nlayer)=pdt(1:ngrid,1:nlayer)+zdtdif(1:ngrid,1:nlayer)
1016           zdtsurf(1:ngrid)=zdtsurf(1:ngrid)+zdtsdif(1:ngrid)
1017         endif
1018
1019         if(ok_slab_ocean)then
1020            flux_sens_lat(1:ngrid)=(zdtsdif(1:ngrid)*capcal(1:ngrid)-fluxrad(1:ngrid))
1021         endif
1022
1023
1024         if (tracer) then
1025           pdq(1:ngrid,1:nlayer,1:nq)=pdq(1:ngrid,1:nlayer,1:nq)+ zdqdif(1:ngrid,1:nlayer,1:nq)
1026           dqsurf(1:ngrid,1:nq)=dqsurf(1:ngrid,1:nq) + zdqsdif(1:ngrid,1:nq)
1027         end if ! of if (tracer)
1028
1029
1030         ! test energy conservation
1031         !-------------------------
1032         if(enertest)then
1033         
1034            dEzdiff(:,:)=cpp*mass(:,:)*zdtdif(:,:)
1035            do ig = 1, ngrid
1036               dEdiff(ig)=SUM(dEzdiff (ig,:))+ sensibFlux(ig)! subtract flux to the ground
1037               dEzdiff(ig,1)= dEzdiff(ig,1)+ sensibFlux(ig)! subtract flux to the ground
1038            enddo
1039           
1040            call planetwide_sumval(dEdiff(:)*cell_area(:)/totarea_planet,dEtot)
1041            dEdiffs(:)=capcal(:)*zdtsdif(:)-zflubid(:)-sensibFlux(:)
1042            call planetwide_sumval(dEdiffs(:)*cell_area(:)/totarea_planet,dEtots)
1043            call planetwide_sumval(sensibFlux(:)*cell_area(:)/totarea_planet,AtmToSurf_TurbFlux)
1044           
1045            if (is_master) then
1046           
1047               if (UseTurbDiff) then
1048                         print*,'In TurbDiff sensible flux (atm=>surf) =',AtmToSurf_TurbFlux,' W m-2'
1049                         print*,'In TurbDiff non-cons atm nrj change   =',dEtot,' W m-2'
1050                  print*,'In TurbDiff (correc rad+latent heat) surf nrj change =',dEtots,' W m-2'
1051               else
1052                         print*,'In vdifc sensible flux (atm=>surf)    =',AtmToSurf_TurbFlux,' W m-2'
1053                         print*,'In vdifc non-cons atm nrj change      =',dEtot,' W m-2'
1054                         print*,'In vdifc (correc rad+latent heat) surf nrj change =',dEtots,' W m-2'
1055               end if
1056            endif ! end of 'is_master'
1057           
1058         ! JL12 : note that the black body radiative flux emitted by the surface has been updated by the implicit scheme but not given back elsewhere.
1059         endif ! end of 'enertest'
1060
1061
1062         ! Test water conservation.
1063         if(watertest.and.water)then
1064         
1065            call planetwide_sumval(massarea(:,:)*zdqdif(:,:,igcm_h2o_vap)*ptimestep/totarea_planet,dWtot_tmp)
1066            call planetwide_sumval(zdqsdif(:,igcm_h2o_vap)*cell_area(:)*ptimestep/totarea_planet,dWtots_tmp)
1067            do ig = 1, ngrid
1068               vdifcncons(ig)=SUM(mass(ig,:)*zdqdif(ig,:,igcm_h2o_vap))
1069            enddo
1070            call planetwide_sumval(massarea(:,:)*zdqdif(:,:,igcm_h2o_ice)*ptimestep/totarea_planet,dWtot)
1071            call planetwide_sumval(zdqsdif(:,igcm_h2o_ice)*cell_area(:)*ptimestep/totarea_planet,dWtots)
1072            dWtot = dWtot + dWtot_tmp
1073            dWtots = dWtots + dWtots_tmp
1074            do ig = 1, ngrid
1075               vdifcncons(ig)=vdifcncons(ig) + SUM(mass(ig,:)*zdqdif(ig,:,igcm_h2o_ice))
1076            enddo           
1077            call planetwide_maxval(vdifcncons(:),nconsMAX)
1078
1079            if (is_master) then
1080               print*,'---------------------------------------------------------------'
1081               print*,'In difv atmospheric water change        =',dWtot,' kg m-2'
1082               print*,'In difv surface water change            =',dWtots,' kg m-2'
1083               print*,'In difv non-cons factor                 =',dWtot+dWtots,' kg m-2'
1084               print*,'In difv MAX non-cons factor             =',nconsMAX,' kg m-2 s-1'
1085            endif
1086
1087         endif ! end of 'watertest'
1088         !-------------------------
1089
1090      else ! calldifv
1091
1092         if(.not.newtonian)then
1093
1094            zdtsurf(1:ngrid) = zdtsurf(1:ngrid) + (fluxrad(1:ngrid) + fluxgrd(1:ngrid))/capcal(1:ngrid)
1095
1096         endif
1097
1098      endif ! end of 'calldifv'
1099
1100
1101!----------------------------------
1102!   IV. Dry convective adjustment :
1103!----------------------------------
1104
1105      if(calladj) then
1106
1107         zdh(1:ngrid,1:nlayer) = pdt(1:ngrid,1:nlayer)/zpopsk(1:ngrid,1:nlayer)
1108         zduadj(1:ngrid,1:nlayer)=0.0
1109         zdvadj(1:ngrid,1:nlayer)=0.0
1110         zdhadj(1:ngrid,1:nlayer)=0.0
1111
1112
1113         call convadj(ngrid,nlayer,nq,ptimestep,            &
1114                      pplay,pplev,zpopsk,                   &
1115                      pu,pv,zh,pq,                          &
1116                      pdu,pdv,zdh,pdq,                      &
1117                      zduadj,zdvadj,zdhadj,                 &
1118                      zdqadj)
1119
1120         pdu(1:ngrid,1:nlayer) = pdu(1:ngrid,1:nlayer) + zduadj(1:ngrid,1:nlayer)
1121         pdv(1:ngrid,1:nlayer) = pdv(1:ngrid,1:nlayer) + zdvadj(1:ngrid,1:nlayer)
1122         pdt(1:ngrid,1:nlayer)    = pdt(1:ngrid,1:nlayer) + zdhadj(1:ngrid,1:nlayer)*zpopsk(1:ngrid,1:nlayer)
1123         zdtadj(1:ngrid,1:nlayer) = zdhadj(1:ngrid,1:nlayer)*zpopsk(1:ngrid,1:nlayer) ! for diagnostic only
1124
1125         if(tracer) then
1126            pdq(1:ngrid,1:nlayer,1:nq) = pdq(1:ngrid,1:nlayer,1:nq) + zdqadj(1:ngrid,1:nlayer,1:nq)
1127         end if
1128
1129         ! Test energy conservation
1130         if(enertest)then
1131            call planetwide_sumval(cpp*massarea(:,:)*zdtadj(:,:)/totarea_planet,dEtot)
1132            if (is_master) print*,'In convadj atmospheric energy change  =',dEtot,' W m-2'
1133         endif
1134
1135         ! Test water conservation
1136         if(watertest)then
1137            call planetwide_sumval(massarea(:,:)*zdqadj(:,:,igcm_h2o_vap)*ptimestep/totarea_planet,dWtot_tmp)
1138            do ig = 1, ngrid
1139               cadjncons(ig)=SUM(mass(ig,:)*zdqadj(ig,:,igcm_h2o_vap))
1140            enddo
1141            call planetwide_sumval(massarea(:,:)*zdqadj(:,:,igcm_h2o_ice)*ptimestep/totarea_planet,dWtot)
1142            dWtot = dWtot + dWtot_tmp
1143            do ig = 1, ngrid
1144               cadjncons(ig)=cadjncons(ig) + SUM(mass(ig,:)*zdqadj(ig,:,igcm_h2o_ice))
1145            enddo           
1146            call planetwide_maxval(cadjncons(:),nconsMAX)
1147
1148            if (is_master) then
1149               print*,'In convadj atmospheric water change     =',dWtot,' kg m-2'
1150               print*,'In convadj MAX non-cons factor          =',nconsMAX,' kg m-2 s-1'
1151            endif
1152           
1153         endif ! end of 'watertest'
1154         
1155      endif ! end of 'calladj'
1156
1157!-----------------------------------------------
1158!   V. Carbon dioxide condensation-sublimation :
1159!-----------------------------------------------
1160
1161      if (co2cond) then
1162     
1163         if (.not.tracer) then
1164            print*,'We need a CO2 ice tracer to condense CO2'
1165            call abort
1166         endif
1167         call condense_co2(ngrid,nlayer,nq,ptimestep,                    &
1168                           capcal,pplay,pplev,tsurf,pt,                  &
1169                           pdt,zdtsurf,pq,pdq,                           &
1170                           qsurf,zdqsurfc,albedo,emis,                   &
1171                           albedo_bareground,albedo_co2_ice_SPECTV,      &
1172                           zdtc,zdtsurfc,pdpsrf,zdqc)
1173
1174         pdt(1:ngrid,1:nlayer) = pdt(1:ngrid,1:nlayer)+zdtc(1:ngrid,1:nlayer)
1175         zdtsurf(1:ngrid)      = zdtsurf(1:ngrid) + zdtsurfc(1:ngrid)
1176
1177         pdq(1:ngrid,1:nlayer,1:nq)   = pdq(1:ngrid,1:nlayer,1:nq)+ zdqc(1:ngrid,1:nlayer,1:nq)
1178         dqsurf(1:ngrid,igcm_co2_ice) = dqsurf(1:ngrid,igcm_co2_ice) + zdqsurfc(1:ngrid)
1179
1180         ! test energy conservation
1181         if(enertest)then
1182            call planetwide_sumval(cpp*massarea(:,:)*zdtc(:,:)/totarea_planet,dEtot)
1183            call planetwide_sumval(capcal(:)*zdtsurfc(:)*cell_area(:)/totarea_planet,dEtots)
1184            if (is_master) then
1185               print*,'In co2cloud atmospheric energy change   =',dEtot,' W m-2'
1186               print*,'In co2cloud surface energy change       =',dEtots,' W m-2'
1187            endif
1188         endif
1189
1190      endif  ! end of 'co2cond'
1191
1192
1193!---------------------------------------------
1194!   VI. Specific parameterizations for tracers
1195!---------------------------------------------
1196
1197      if (tracer) then
1198     
1199  ! ---------------------
1200  !   VI.1. Water and ice
1201  ! ---------------------
1202         if (water) then
1203
1204            ! Water ice condensation in the atmosphere
1205            if(watercond.and.(RLVTT.gt.1.e-8))then
1206
1207               dqmoist(1:ngrid,1:nlayer,1:nq)=0.
1208               dtmoist(1:ngrid,1:nlayer)=0.
1209               
1210               ! Moist Convective Adjustment.
1211               ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1212               call moistadj(ngrid,nlayer,nq,pt,pq,pdq,pplev,pplay,dtmoist,dqmoist,ptimestep,rneb_man)
1213
1214               pdq(1:ngrid,1:nlayer,igcm_h2o_vap) = pdq(1:ngrid,1:nlayer,igcm_h2o_vap)     &
1215                                                  + dqmoist(1:ngrid,1:nlayer,igcm_h2o_vap)
1216               pdq(1:ngrid,1:nlayer,igcm_h2o_ice) = pdq(1:ngrid,1:nlayer,igcm_h2o_ice)     &
1217                                                  + dqmoist(1:ngrid,1:nlayer,igcm_h2o_ice)
1218               pdt(1:ngrid,1:nlayer) = pdt(1:ngrid,1:nlayer)+dtmoist(1:ngrid,1:nlayer)
1219
1220               ! Test energy conservation.
1221               if(enertest)then
1222               
1223                  call planetwide_sumval(cpp*massarea(:,:)*dtmoist(:,:)/totarea_planet,dEtot)
1224                  call planetwide_maxval(dtmoist(:,:),dtmoist_max)
1225                  call planetwide_minval(dtmoist(:,:),dtmoist_min)
1226                  madjdEz(:,:)=cpp*mass(:,:)*dtmoist(:,:)
1227                  do ig=1,ngrid
1228                     madjdE(ig) = cpp*SUM(mass(:,:)*dtmoist(:,:))
1229                  enddo
1230                 
1231                  if (is_master) then
1232                     print*,'In moistadj atmospheric energy change   =',dEtot,' W m-2'
1233                     print*,'In moistadj MAX atmospheric energy change   =',dtmoist_max*ptimestep,'K/step'
1234                     print*,'In moistadj MIN atmospheric energy change   =',dtmoist_min*ptimestep,'K/step'
1235                  endif
1236                 
1237                  call planetwide_sumval(massarea(:,:)*dqmoist(:,:,igcm_h2o_vap)*ptimestep/totarea_planet+        &
1238                                           massarea(:,:)*dqmoist(:,:,igcm_h2o_ice)*ptimestep/totarea_planet,dWtot)
1239                  if (is_master) print*,'In moistadj atmospheric water change    =',dWtot,' kg m-2'
1240                 
1241               endif ! end of 'enertest'
1242               
1243
1244               ! Large scale condensation/evaporation.
1245               ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1246               call largescale(ngrid,nlayer,nq,ptimestep,pplev,pplay,pt,pq,pdt,pdq,dtlscale,dqvaplscale,dqcldlscale,rneb_lsc)
1247
1248               pdt(1:ngrid,1:nlayer) = pdt(1:ngrid,1:nlayer)+dtlscale(1:ngrid,1:nlayer)
1249               pdq(1:ngrid,1:nlayer,igcm_h2o_vap) = pdq(1:ngrid,1:nlayer,igcm_h2o_vap)+dqvaplscale(1:ngrid,1:nlayer)
1250               pdq(1:ngrid,1:nlayer,igcm_h2o_ice) = pdq(1:ngrid,1:nlayer,igcm_h2o_ice)+dqcldlscale(1:ngrid,1:nlayer)
1251
1252               ! Test energy conservation.
1253               if(enertest)then
1254                  lscaledEz(:,:) = cpp*mass(:,:)*dtlscale(:,:)
1255                  do ig=1,ngrid
1256                     lscaledE(ig) = cpp*SUM(mass(:,:)*dtlscale(:,:))
1257                  enddo
1258                  call planetwide_sumval(cpp*massarea(:,:)*dtlscale(:,:)/totarea_planet,dEtot)
1259
1260                  if (is_master) print*,'In largescale atmospheric energy change =',dEtot,' W m-2'
1261
1262                  ! Test water conservation.
1263                  call planetwide_sumval(massarea(:,:)*dqvaplscale(:,:)*ptimestep/totarea_planet+        &
1264                                           SUM(massarea(:,:)*dqcldlscale(:,:))*ptimestep/totarea_planet,dWtot)
1265                       
1266                  if (is_master) print*,'In largescale atmospheric water change  =',dWtot,' kg m-2'
1267               endif ! end of 'enertest'
1268
1269               ! Compute cloud fraction.
1270               do l = 1, nlayer
1271                  do ig=1,ngrid
1272                     cloudfrac(ig,l)=MAX(rneb_lsc(ig,l),rneb_man(ig,l))
1273                  enddo
1274               enddo
1275
1276            endif ! end of 'watercond'
1277           
1278            ! Water ice / liquid precipitation.
1279            ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1280            if(waterrain)then
1281
1282               zdqrain(1:ngrid,1:nlayer,1:nq) = 0.0
1283               zdqsrain(1:ngrid)    = 0.0
1284               zdqssnow(1:ngrid)    = 0.0
1285
1286               call rain(ngrid,nlayer,nq,ptimestep,pplev,pplay,pt,pdt,pq,pdq,            &
1287                         zdtrain,zdqrain,zdqsrain,zdqssnow,reevap_precip,cloudfrac)
1288
1289               pdq(1:ngrid,1:nlayer,igcm_h2o_vap) = pdq(1:ngrid,1:nlayer,igcm_h2o_vap) &
1290                                                  + zdqrain(1:ngrid,1:nlayer,igcm_h2o_vap)
1291               pdq(1:ngrid,1:nlayer,igcm_h2o_ice) = pdq(1:ngrid,1:nlayer,igcm_h2o_ice) &
1292                                                  + zdqrain(1:ngrid,1:nlayer,igcm_h2o_ice)
1293               pdt(1:ngrid,1:nlayer) = pdt(1:ngrid,1:nlayer)+zdtrain(1:ngrid,1:nlayer)
1294               
1295               dqsurf(1:ngrid,igcm_h2o_vap) = dqsurf(1:ngrid,igcm_h2o_vap)+zdqsrain(1:ngrid)
1296               dqsurf(1:ngrid,igcm_h2o_ice) = dqsurf(1:ngrid,igcm_h2o_ice)+zdqssnow(1:ngrid)
1297
1298               ! Test energy conservation.
1299               if(enertest)then
1300               
1301                  call planetwide_sumval(cpp*massarea(:,:)*zdtrain(:,:)/totarea_planet,dEtot)
1302                  if (is_master) print*,'In rain atmospheric T energy change       =',dEtot,' W m-2'
1303                  call planetwide_sumval(massarea(:,:)*zdqrain(:,:,igcm_h2o_ice)/totarea_planet*RLVTT/cpp,dItot_tmp)
1304                  call planetwide_sumval(cell_area(:)*zdqssnow(:)/totarea_planet*RLVTT/cpp,dItot)
1305                  dItot = dItot + dItot_tmp
1306                  call planetwide_sumval(massarea(:,:)*zdqrain(:,:,igcm_h2o_vap)*ptimestep/totarea_planet,dVtot_tmp)
1307                  call planetwide_sumval(cell_area(:)*zdqsrain(:)/totarea_planet*RLVTT/cpp,dVtot)
1308                  dVtot = dVtot + dVtot_tmp
1309                  dEtot = dItot + dVtot
1310                 
1311                  if (is_master) then
1312                     print*,'In rain dItot =',dItot,' W m-2'
1313                     print*,'In rain dVtot =',dVtot,' W m-2'
1314                     print*,'In rain atmospheric L energy change       =',dEtot,' W m-2'
1315                  endif
1316
1317                  ! Test water conservation
1318                  call planetwide_sumval(massarea(:,:)*zdqrain(:,:,igcm_h2o_vap)*ptimestep/totarea_planet+        &
1319                          massarea(:,:)*zdqrain(:,:,igcm_h2o_ice)*ptimestep/totarea_planet,dWtot)
1320                  call planetwide_sumval((zdqsrain(:)+zdqssnow(:))*cell_area(:)*ptimestep/totarea_planet,dWtots)
1321                 
1322                  if (is_master) then
1323                          print*,'In rain atmospheric water change        =',dWtot,' kg m-2'
1324                          print*,'In rain surface water change            =',dWtots,' kg m-2'
1325                          print*,'In rain non-cons factor                 =',dWtot+dWtots,' kg m-2'
1326                  endif
1327                 
1328               endif ! end of 'enertest'
1329
1330            end if ! enf of 'waterrain'
1331           
1332         end if ! end of 'water'
1333
1334 ! -------------------------
1335  !   VI.2. Photochemistry
1336  ! -------------------------
1337
1338         IF (photochem) then
1339
1340             DO ig=1,ngrid
1341               array_zero1(ig)=0.0
1342               DO l=1,nlayer
1343                 array_zero2(ig,l)=0.
1344               ENDDO
1345             ENDDO
1346
1347            call calchim_asis(ngrid,nlayer,nq,                          &
1348                        ptimestep,pplay,pplev,pt,pdt,dist_star,mu0,     &
1349                        fract,zzlev,zzlay,zday,pq,pdq,zdqchim,zdqschim, &
1350                        array_zero1,array_zero1,                        &
1351                        pu,pdu,pv,pdv,array_zero2,array_zero2)
1352
1353           ! increment values of tracers:
1354           DO iq=1,nq ! loop on all tracers; tendencies for non-chemistry
1355                      ! tracers is zero anyways
1356             DO l=1,nlayer
1357               DO ig=1,ngrid
1358                 pdq(ig,l,iq)=pdq(ig,l,iq)+zdqchim(ig,l,iq)
1359               ENDDO
1360             ENDDO
1361           ENDDO ! of DO iq=1,nq
1362
1363
1364           ! increment surface values of tracers:
1365           DO iq=1,nq ! loop on all tracers; tendencies for non-chemistry
1366                      ! tracers is zero anyways
1367             DO ig=1,ngrid
1368!               dqsurf(ig,iq)=dqsurf(ig,iq)+zdqschim(ig,iq)
1369             ENDDO
1370           ENDDO ! of DO iq=1,nq
1371
1372         END IF  ! of IF (photochem)
1373
1374
1375
1376  ! -------------------------
1377  !   VI.3. Aerosol particles
1378  ! -------------------------
1379
1380         ! Sedimentation.
1381         if (sedimentation) then
1382       
1383            zdqsed(1:ngrid,1:nlayer,1:nq) = 0.0
1384            zdqssed(1:ngrid,1:nq)  = 0.0
1385
1386            if(watertest)then
1387           
1388               iq=igcm_h2o_ice
1389               call planetwide_sumval(massarea(:,:)*pq(:,:,iq)*ptimestep/totarea_planet,dWtot)
1390               call planetwide_sumval(massarea(:,:)*pdq(:,:,iq)*ptimestep/totarea_planet,dWtots)
1391               if (is_master) then
1392                        print*,'Before sedim pq  =',dWtot,' kg m-2'
1393                  print*,'Before sedim pdq =',dWtots,' kg m-2'
1394               endif
1395            endif
1396           
1397            call callsedim(ngrid,nlayer,ptimestep,           &
1398                          pplev,zzlev,pt,pdt,pq,pdq,        &
1399                          zdqsed,zdqssed,nq)
1400
1401            if(watertest)then
1402               iq=igcm_h2o_ice
1403               call planetwide_sumval(massarea(:,:)*pq(:,:,iq)*ptimestep/totarea_planet,dWtot)
1404               call planetwide_sumval(massarea(:,:)*pdq(:,:,iq)*ptimestep/totarea_planet,dWtots)
1405               if (is_master) then
1406                        print*,'After sedim pq  =',dWtot,' kg m-2'
1407                        print*,'After sedim pdq =',dWtots,' kg m-2'
1408                 endif
1409            endif
1410
1411            ! Whether it falls as rain or snow depends only on the surface temperature
1412            pdq(1:ngrid,1:nlayer,1:nq) = pdq(1:ngrid,1:nlayer,1:nq) + zdqsed(1:ngrid,1:nlayer,1:nq)
1413            dqsurf(1:ngrid,1:nq) = dqsurf(1:ngrid,1:nq) + zdqssed(1:ngrid,1:nq)
1414
1415            ! Test water conservation
1416            if(watertest)then
1417               call planetwide_sumval(massarea(:,:)*(zdqsed(:,:,igcm_h2o_vap)+zdqsed(:,:,igcm_h2o_ice))*ptimestep/totarea_planet,dWtot)
1418               call planetwide_sumval((zdqssed(:,igcm_h2o_vap)+zdqssed(:,igcm_h2o_ice))*cell_area(:)*ptimestep/totarea_planet,dWtots)
1419               if (is_master) then
1420                        print*,'In sedim atmospheric ice change         =',dWtot,' kg m-2'
1421                        print*,'In sedim surface ice change             =',dWtots,' kg m-2'
1422                        print*,'In sedim non-cons factor                =',dWtot+dWtots,' kg m-2'
1423               endif
1424            endif
1425
1426         end if ! end of 'sedimentation'
1427
1428
1429  ! ---------------
1430  !   VI.4. Updates
1431  ! ---------------
1432
1433         ! Updating Atmospheric Mass and Tracers budgets.
1434         if(mass_redistrib) then
1435
1436            zdmassmr(1:ngrid,1:nlayer) = mass(1:ngrid,1:nlayer) *     &
1437               (   zdqevap(1:ngrid,1:nlayer)                          &
1438                 + zdqrain(1:ngrid,1:nlayer,igcm_h2o_vap)             &
1439                 + dqmoist(1:ngrid,1:nlayer,igcm_h2o_vap)             &
1440                 + dqvaplscale(1:ngrid,1:nlayer) )
1441
1442            do ig = 1, ngrid
1443               zdmassmr_col(ig)=SUM(zdmassmr(ig,1:nlayer))
1444            enddo
1445           
1446            call writediagfi(ngrid,"mass_evap","mass gain"," ",3,zdmassmr)
1447            call writediagfi(ngrid,"mass_evap_col","mass gain col"," ",2,zdmassmr_col)
1448            call writediagfi(ngrid,"mass","mass","kg/m2",3,mass)
1449
1450            call mass_redistribution(ngrid,nlayer,nq,ptimestep,                     &
1451                                     rnat,capcal,pplay,pplev,pt,tsurf,pq,qsurf,     &
1452                                     pu,pv,pdt,zdtsurf,pdq,pdu,pdv,zdmassmr,        &
1453                                     zdtmr,zdtsurfmr,zdpsrfmr,zdumr,zdvmr,zdqmr,zdqsurfmr)
1454         
1455            pdq(1:ngrid,1:nlayer,1:nq) = pdq(1:ngrid,1:nlayer,1:nq) + zdqmr(1:ngrid,1:nlayer,1:nq)
1456            dqsurf(1:ngrid,1:nq)       = dqsurf(1:ngrid,1:nq) + zdqsurfmr(1:ngrid,1:nq)
1457            pdt(1:ngrid,1:nlayer)      = pdt(1:ngrid,1:nlayer) + zdtmr(1:ngrid,1:nlayer)
1458            pdu(1:ngrid,1:nlayer)      = pdu(1:ngrid,1:nlayer) + zdumr(1:ngrid,1:nlayer)
1459            pdv(1:ngrid,1:nlayer)      = pdv(1:ngrid,1:nlayer) + zdvmr(1:ngrid,1:nlayer)
1460            pdpsrf(1:ngrid)            = pdpsrf(1:ngrid) + zdpsrfmr(1:ngrid)
1461            zdtsurf(1:ngrid)           = zdtsurf(1:ngrid) + zdtsurfmr(1:ngrid)
1462           
1463         endif
1464
1465  ! ------------------
1466  !   VI.5. Slab Ocean
1467  ! ------------------
1468
1469         if (ok_slab_ocean)then
1470
1471            do ig=1,ngrid
1472               qsurfint(:,igcm_h2o_ice)=qsurf(:,igcm_h2o_ice)
1473            enddo
1474
1475            call ocean_slab_ice(ptimestep,                          &
1476                                ngrid, knindex, tsea_ice, fluxrad,  &
1477                                zdqssnow, qsurf(:,igcm_h2o_ice),    &
1478                              - zdqsdif(:,igcm_h2o_vap),            &
1479                                flux_sens_lat,tsea_ice, pctsrf_sic, &
1480                                taux,tauy,icount)
1481
1482
1483            call ocean_slab_get_vars(ngrid,tslab,      &
1484                                     sea_ice, flux_o,  &
1485                                     flux_g, dt_hdiff, &
1486                                     dt_ekman)
1487   
1488            do ig=1,ngrid
1489               if (nint(rnat(ig)).eq.1)then
1490                  tslab(ig,1) = 0.
1491                  tslab(ig,2) = 0.
1492                  tsea_ice(ig) = 0.
1493                  sea_ice(ig) = 0.
1494                  pctsrf_sic(ig) = 0.
1495                  qsurf(ig,igcm_h2o_ice) = qsurfint(ig,igcm_h2o_ice)
1496               endif
1497            enddo
1498
1499         endif ! end of 'ok_slab_ocean'
1500
1501  ! -----------------------------
1502  !   VI.6. Surface Tracer Update
1503  ! -----------------------------
1504
1505         if(hydrology)then
1506           
1507            call hydrol(ngrid,nq,ptimestep,rnat,tsurf,qsurf,dqsurf,dqs_hyd, &
1508                        capcal,albedo,albedo_bareground,                    &
1509                        albedo_snow_SPECTV,albedo_co2_ice_SPECTV,           &
1510                        mu0,zdtsurf,zdtsurf_hyd,hice,pctsrf_sic,            &
1511                        sea_ice)
1512
1513            zdtsurf(1:ngrid)     = zdtsurf(1:ngrid) + zdtsurf_hyd(1:ngrid)
1514            dqsurf(1:ngrid,1:nq) = dqsurf(1:ngrid,1:nq) + dqs_hyd(1:ngrid,1:nq)
1515           
1516            qsurf(1:ngrid,1:nq) = qsurf(1:ngrid,1:nq) + ptimestep*dqsurf(1:ngrid,1:nq)
1517
1518            ! Test energy conservation
1519            if(enertest)then
1520               call planetwide_sumval(cell_area(:)*capcal(:)*zdtsurf_hyd(:)/totarea_planet,dEtots)
1521               if (is_master) print*,'In hydrol surface energy change     =',dEtots,' W m-2'
1522            endif
1523
1524            ! test water conservation
1525            if(watertest)then
1526               call planetwide_sumval(dqs_hyd(:,igcm_h2o_ice)*cell_area(:)*ptimestep/totarea_planet,dWtots)
1527               if (is_master) print*,'In hydrol surface ice change            =',dWtots,' kg m-2'
1528                  call planetwide_sumval(dqs_hyd(:,igcm_h2o_vap)*cell_area(:)*ptimestep/totarea_planet,dWtots)
1529               if (is_master) then
1530                  print*,'In hydrol surface water change          =',dWtots,' kg m-2'
1531                  print*,'---------------------------------------------------------------'
1532               endif
1533            endif
1534
1535         else ! of if (hydrology)
1536
1537            qsurf(1:ngrid,1:nq) = qsurf(1:ngrid,1:nq) + ptimestep*dqsurf(1:ngrid,1:nq)
1538
1539         end if ! of if (hydrology)
1540
1541         ! Add qsurf to qsurf_hist, which is what we save in diagfi.nc. At the same time, we set the water
1542         ! content of ocean gridpoints back to zero, in order to avoid rounding errors in vdifc, rain.
1543         qsurf_hist(:,:) = qsurf(:,:)
1544
1545         if(ice_update)then
1546            ice_min(1:ngrid)=min(ice_min(1:ngrid),qsurf(1:ngrid,igcm_h2o_ice))
1547         endif
1548
1549      endif! end of if 'tracer'
1550
1551
1552!------------------------------------------------
1553!   VII. Surface and sub-surface soil temperature
1554!------------------------------------------------
1555
1556
1557      ! Increment surface temperature
1558      if(ok_slab_ocean)then
1559         do ig=1,ngrid
1560           if (nint(rnat(ig)).eq.0)then
1561             zdtsurf(ig)= (tslab(ig,1)              &
1562             + pctsrf_sic(ig)*(tsea_ice(ig)-tslab(ig,1))-tsurf(ig))/ptimestep
1563           endif
1564           tsurf(ig)=tsurf(ig)+ptimestep*zdtsurf(ig)
1565         enddo
1566   
1567      else
1568        tsurf(1:ngrid)=tsurf(1:ngrid)+ptimestep*zdtsurf(1:ngrid)
1569      endif ! end of 'ok_slab_ocean'
1570
1571
1572      ! Compute soil temperatures and subsurface heat flux.
1573      if (callsoil) then
1574         call soil(ngrid,nsoilmx,.false.,lastcall,inertiedat,   &
1575                   ptimestep,tsurf,tsoil,capcal,fluxgrd)     
1576      endif
1577
1578
1579      if (ok_slab_ocean) then
1580     
1581         do ig=1,ngrid
1582         
1583            fluxgrdocean(ig)=fluxgrd(ig)
1584            if (nint(rnat(ig)).eq.0) then
1585               capcal(ig)=capcalocean
1586               fluxgrd(ig)=0.
1587               fluxgrdocean(ig)=pctsrf_sic(ig)*flux_g(ig)+(1-pctsrf_sic(ig))*(dt_hdiff(ig,1)+dt_ekman(ig,1))
1588               do iq=1,nsoilmx
1589                  tsoil(ig,iq)=tsurf(ig)
1590               enddo
1591               if (pctsrf_sic(ig).gt.0.5) then
1592                  capcal(ig)=capcalseaice
1593                  if (qsurf(ig,igcm_h2o_ice).gt.0.) then
1594                     capcal(ig)=capcalsno
1595                  endif
1596               endif               
1597            endif
1598           
1599         enddo
1600         
1601      endif !end of 'ok_slab_ocean'
1602
1603
1604      ! Test energy conservation
1605      if(enertest)then
1606         call planetwide_sumval(cell_area(:)*capcal(:)*zdtsurf(:)/totarea_planet,dEtots)         
1607         if (is_master) print*,'Surface energy change                 =',dEtots,' W m-2'
1608      endif
1609
1610
1611!---------------------------------------------------
1612!   VIII. Perform diagnostics and write output files
1613!---------------------------------------------------
1614
1615      ! Note : For output only: the actual model integration is performed in the dynamics.
1616
1617
1618 
1619      ! Temperature, zonal and meridional winds.
1620      zt(1:ngrid,1:nlayer) = pt(1:ngrid,1:nlayer) + pdt(1:ngrid,1:nlayer)*ptimestep
1621      zu(1:ngrid,1:nlayer) = pu(1:ngrid,1:nlayer) + pdu(1:ngrid,1:nlayer)*ptimestep
1622      zv(1:ngrid,1:nlayer) = pv(1:ngrid,1:nlayer) + pdv(1:ngrid,1:nlayer)*ptimestep
1623
1624      ! Diagnostic.
1625      zdtdyn(1:ngrid,1:nlayer)     = (pt(1:ngrid,1:nlayer)-ztprevious(1:ngrid,1:nlayer)) / ptimestep
1626      ztprevious(1:ngrid,1:nlayer) = zt(1:ngrid,1:nlayer)
1627
1628      zdudyn(1:ngrid,1:nlayer)     = (pu(1:ngrid,1:nlayer)-zuprevious(1:ngrid,1:nlayer)) / ptimestep
1629      zuprevious(1:ngrid,1:nlayer) = zu(1:ngrid,1:nlayer)
1630
1631      if(firstcall)then
1632         zdtdyn(1:ngrid,1:nlayer)=0.0
1633         zdudyn(1:ngrid,1:nlayer)=0.0
1634      endif
1635
1636      ! Dynamical heating diagnostic.
1637      do ig=1,ngrid
1638         fluxdyn(ig)= SUM(zdtdyn(ig,:) *mass(ig,:))*cpp
1639      enddo
1640
1641      ! Tracers.
1642      zq(1:ngrid,1:nlayer,1:nq) = pq(1:ngrid,1:nlayer,1:nq) + pdq(1:ngrid,1:nlayer,1:nq)*ptimestep
1643
1644      ! Surface pressure.
1645      ps(1:ngrid) = pplev(1:ngrid,1) + pdpsrf(1:ngrid)*ptimestep
1646
1647
1648
1649      ! Surface and soil temperature information
1650      call planetwide_sumval(cell_area(:)*tsurf(:)/totarea_planet,Ts1)
1651      call planetwide_minval(tsurf(:),Ts2)
1652      call planetwide_maxval(tsurf(:),Ts3)
1653      if(callsoil)then
1654         TsS = SUM(cell_area(:)*tsoil(:,nsoilmx))/totarea        ! mean temperature at bottom soil layer
1655         if (is_master) then
1656            print*,'          ave[Tsurf]             min[Tsurf]             max[Tsurf]             ave[Tdeep]'
1657            print*,Ts1,Ts2,Ts3,TsS
1658         end if
1659      else
1660         if (is_master) then
1661            print*,'          ave[Tsurf]             min[Tsurf]             max[Tsurf]'
1662            print*,Ts1,Ts2,Ts3
1663         endif
1664      end if
1665
1666
1667      ! Check the energy balance of the simulation during the run
1668      if(corrk)then
1669
1670         call planetwide_sumval(cell_area(:)*fluxtop_dn(:)/totarea_planet,ISR)
1671         call planetwide_sumval(cell_area(:)*fluxabs_sw(:)/totarea_planet,ASR)
1672         call planetwide_sumval(cell_area(:)*fluxtop_lw(:)/totarea_planet,OLR)
1673         call planetwide_sumval(cell_area(:)*fluxgrd(:)/totarea_planet,GND)
1674         call planetwide_sumval(cell_area(:)*fluxdyn(:)/totarea_planet,DYN)
1675         do ig=1,ngrid
1676            if(fluxtop_dn(ig).lt.0.0)then
1677               print*,'fluxtop_dn has gone crazy'
1678               print*,'fluxtop_dn=',fluxtop_dn(ig)
1679               print*,'tau_col=',tau_col(ig)
1680               print*,'aerosol=',aerosol(ig,:,:)
1681               print*,'temp=   ',pt(ig,:)
1682               print*,'pplay=  ',pplay(ig,:)
1683               call abort
1684            endif
1685         end do
1686                     
1687         if(ngrid.eq.1)then
1688            DYN=0.0
1689         endif
1690         
1691         if (is_master) then
1692            print*,'  ISR            ASR            OLR            GND            DYN [W m^-2]'
1693            print*, ISR,ASR,OLR,GND,DYN
1694         endif
1695
1696         if(enertest .and. is_master)then
1697            print*,'SW flux/heating difference SW++ - ASR = ',dEtotSW+dEtotsSW-ASR,' W m-2'
1698            print*,'LW flux/heating difference LW++ - OLR = ',dEtotLW+dEtotsLW+OLR,' W m-2'
1699            print*,'LW energy balance LW++ + ASR = ',dEtotLW+dEtotsLW+ASR,' W m-2'
1700         endif
1701
1702         if(meanOLR .and. is_master)then
1703            if((ngrid.gt.1) .or. (mod(icount-1,ecritphy).eq.0))then
1704               ! to record global radiative balance
1705               open(92,file="rad_bal.out",form='formatted',position='append')
1706               write(92,*) zday,ISR,ASR,OLR
1707               close(92)
1708               open(93,file="tem_bal.out",form='formatted',position='append')
1709               if(callsoil)then
1710                         write(93,*) zday,Ts1,Ts2,Ts3,TsS
1711               else
1712                  write(93,*) zday,Ts1,Ts2,Ts3
1713               endif
1714               close(93)
1715            endif
1716         endif
1717
1718      endif ! end of 'corrk'
1719
1720
1721      ! Diagnostic to test radiative-convective timescales in code.
1722      if(testradtimes)then
1723         open(38,file="tau_phys.out",form='formatted',position='append')
1724         ig=1
1725         do l=1,nlayer
1726            write(38,*) -1./pdt(ig,l),pt(ig,l),pplay(ig,l)
1727         enddo
1728         close(38)
1729         print*,'As testradtimes enabled,'
1730         print*,'exiting physics on first call'
1731         call abort
1732      endif
1733
1734
1735      ! Compute column amounts (kg m-2) if tracers are enabled.
1736      if(tracer)then
1737         qcol(1:ngrid,1:nq)=0.0
1738         do iq=1,nq
1739            do ig=1,ngrid
1740               qcol(ig,iq) = SUM( zq(ig,1:nlayer,iq) * mass(ig,1:nlayer))
1741            enddo
1742         enddo
1743
1744         ! Generalised for arbitrary aerosols now. By LK
1745         reffcol(1:ngrid,1:naerkind)=0.0
1746         if(co2cond.and.(iaero_co2.ne.0))then
1747            call co2_reffrad(ngrid,nlayer,nq,zq,reffrad(1,1,iaero_co2))
1748            do ig=1,ngrid
1749               reffcol(ig,iaero_co2) = SUM(zq(ig,1:nlayer,igcm_co2_ice)*reffrad(ig,1:nlayer,iaero_co2)*mass(ig,1:nlayer))
1750            enddo
1751         endif
1752         if(water.and.(iaero_h2o.ne.0))then
1753            call h2o_reffrad(ngrid,nlayer,zq(1,1,igcm_h2o_ice),zt, &
1754                             reffrad(1,1,iaero_h2o),nueffrad(1,1,iaero_h2o))
1755            do ig=1,ngrid
1756               reffcol(ig,iaero_h2o) = SUM(zq(ig,1:nlayer,igcm_h2o_ice)*reffrad(ig,1:nlayer,iaero_h2o)*mass(ig,1:nlayer))
1757            enddo
1758         endif
1759
1760      endif ! end of 'tracer'
1761
1762
1763      ! Test for water conservation.
1764      if(water)then
1765
1766         call planetwide_sumval(cell_area(:)*qsurf_hist(:,igcm_h2o_ice)/totarea_planet,icesrf)
1767         call planetwide_sumval(cell_area(:)*qsurf_hist(:,igcm_h2o_vap)/totarea_planet,liqsrf)
1768         call planetwide_sumval(cell_area(:)*qcol(:,igcm_h2o_ice)/totarea_planet,icecol)
1769         call planetwide_sumval(cell_area(:)*qcol(:,igcm_h2o_vap)/totarea_planet,vapcol)
1770
1771         h2otot = icesrf + liqsrf + icecol + vapcol
1772         
1773         if (is_master) then
1774            print*,' Total water amount [kg m^-2]: ',h2otot
1775            print*,' Surface ice    Surface liq.   Atmos. con.     Atmos. vap. [kg m^-2] '
1776            print*, icesrf,liqsrf,icecol,vapcol
1777         endif
1778
1779         if(meanOLR .and. is_master)then
1780            if((ngrid.gt.1) .or. (mod(icount-1,ecritphy).eq.0))then
1781               ! to record global water balance
1782               open(98,file="h2o_bal.out",form='formatted',position='append')
1783               write(98,*) zday,icesrf,liqsrf,icecol,vapcol
1784               close(98)
1785            endif
1786         endif
1787
1788      endif
1789
1790
1791      ! Calculate RH (Relative Humidity) for diagnostic.
1792      if(water)then
1793     
1794         do l = 1, nlayer
1795            do ig=1,ngrid
1796               call Psat_water(zt(ig,l),pplay(ig,l),psat_tmp,qsat(ig,l))
1797               RH(ig,l) = zq(ig,l,igcm_h2o_vap) / qsat(ig,l)
1798            enddo
1799         enddo
1800
1801         ! Compute maximum possible H2O column amount (100% saturation).
1802         do ig=1,ngrid
1803            H2Omaxcol(ig) = SUM( qsat(ig,:) * mass(ig,:))
1804         enddo
1805
1806      endif ! end of 'water'
1807
1808
1809      if (is_master) print*,'--> Ls =',zls*180./pi
1810     
1811     
1812!----------------------------------------------------------------------
1813!        Writing NetCDF file  "RESTARTFI" at the end of the run
1814!----------------------------------------------------------------------
1815
1816!        Note: 'restartfi' is stored just before dynamics are stored
1817!              in 'restart'. Between now and the writting of 'restart',
1818!              there will have been the itau=itau+1 instruction and
1819!              a reset of 'time' (lastacll = .true. when itau+1= itaufin)
1820!              thus we store for time=time+dtvr
1821
1822
1823
1824      if(lastcall) then
1825         ztime_fin = ptime + ptimestep/(float(iphysiq)*daysec)
1826
1827         ! Update surface ice distribution to iterate to steady state if requested
1828         if(ice_update)then
1829
1830            do ig=1,ngrid
1831
1832               delta_ice = (qsurf(ig,igcm_h2o_ice)-ice_initial(ig))
1833               
1834               ! add multiple years of evolution
1835               qsurf_hist(ig,igcm_h2o_ice) = qsurf_hist(ig,igcm_h2o_ice) + delta_ice*icetstep
1836
1837               ! if ice has gone -ve, set to zero
1838               if(qsurf_hist(ig,igcm_h2o_ice).lt.0.0)then
1839                  qsurf_hist(ig,igcm_h2o_ice) = 0.0
1840               endif
1841
1842               ! if ice is seasonal, set to zero (NEW)
1843               if(ice_min(ig).lt.0.01)then
1844                  qsurf_hist(ig,igcm_h2o_ice) = 0.0
1845               endif
1846
1847            enddo
1848
1849            ! enforce ice conservation
1850            ice_tot= SUM(qsurf_hist(:,igcm_h2o_ice)*cell_area(:) )/SUM(cell_area(:))
1851            qsurf_hist(:,igcm_h2o_ice) = qsurf_hist(:,igcm_h2o_ice)*(icesrf/ice_tot)
1852
1853         endif
1854#ifndef MESOSCALE
1855         
1856         if (ngrid.ne.1) then
1857            write(*,*)'PHYSIQ: for physdem ztime_fin =',ztime_fin
1858           
1859            call physdem1("restartfi.nc",nsoilmx,ngrid,nlayer,nq, &
1860                          ptimestep,ztime_fin,                    &
1861                          tsurf,tsoil,emis,q2,qsurf_hist,         &
1862                          cloudfrac,totcloudfrac,hice,            &
1863                          rnat,pctsrf_sic,tslab,tsea_ice,sea_ice)
1864         endif
1865#endif
1866         if(ok_slab_ocean) then
1867            call ocean_slab_final!(tslab, seaice)
1868         end if
1869
1870    endif ! end of 'lastcall'
1871
1872
1873!-----------------------------------
1874!        Saving statistics :
1875!-----------------------------------
1876
1877!    Note :("stats" stores and accumulates 8 key variables in file "stats.nc"
1878!           which can later be used to make the statistic files of the run:
1879!           "stats")          only possible in 3D runs !!!
1880
1881         
1882      if (callstats) then
1883
1884         call wstats(ngrid,"ps","Surface pressure","Pa",2,ps)
1885         call wstats(ngrid,"tsurf","Surface temperature","K",2,tsurf)
1886         call wstats(ngrid,"fluxsurf_lw",                               &
1887                     "Thermal IR radiative flux to surface","W.m-2",2,  &
1888                     fluxsurf_lw)
1889         call wstats(ngrid,"fluxtop_lw",                                &
1890                     "Thermal IR radiative flux to space","W.m-2",2,    &
1891                     fluxtop_lw)
1892                     
1893!            call wstats(ngrid,"fluxsurf_sw",                               &
1894!                        "Solar radiative flux to surface","W.m-2",2,       &
1895!                         fluxsurf_sw_tot)                     
1896!            call wstats(ngrid,"fluxtop_sw",                                &
1897!                        "Solar radiative flux to space","W.m-2",2,         &
1898!                        fluxtop_sw_tot)
1899
1900
1901         call wstats(ngrid,"ISR","incoming stellar rad.","W m-2",2,fluxtop_dn)
1902         call wstats(ngrid,"ASR","absorbed stellar rad.","W m-2",2,fluxabs_sw)
1903         call wstats(ngrid,"OLR","outgoing longwave rad.","W m-2",2,fluxtop_lw)
1904         !call wstats(ngrid,"ALB","Surface albedo"," ",2,albedo_equivalent)
1905         !call wstats(ngrid,"ALB_1st","First Band Surface albedo"," ",2,albedo(:,1))
1906         call wstats(ngrid,"p","Pressure","Pa",3,pplay)
1907         call wstats(ngrid,"temp","Atmospheric temperature","K",3,zt)
1908         call wstats(ngrid,"u","Zonal (East-West) wind","m.s-1",3,zu)
1909         call wstats(ngrid,"v","Meridional (North-South) wind","m.s-1",3,zv)
1910         call wstats(ngrid,"w","Vertical (down-up) wind","m.s-1",3,pw)
1911         call wstats(ngrid,"q2","Boundary layer eddy kinetic energy","m2.s-2",3,q2)
1912
1913         if (tracer) then
1914            do iq=1,nq
1915               call wstats(ngrid,noms(iq),noms(iq),'kg/kg',3,zq(1,1,iq))
1916               call wstats(ngrid,trim(noms(iq))//'_surf',trim(noms(iq))//'_surf',  &
1917                           'kg m^-2',2,qsurf(1,iq) )
1918               call wstats(ngrid,trim(noms(iq))//'_col',trim(noms(iq))//'_col',    &
1919                          'kg m^-2',2,qcol(1,iq) )
1920                         
1921!              call wstats(ngrid,trim(noms(iq))//'_reff',                          &
1922!                          trim(noms(iq))//'_reff',                                   &
1923!                          'm',3,reffrad(1,1,iq))
1924
1925            end do
1926           
1927            if (water) then
1928               vmr=zq(1:ngrid,1:nlayer,igcm_h2o_vap)*mugaz/mmol(igcm_h2o_vap)
1929               call wstats(ngrid,"vmr_h2ovapor",                             &
1930                           "H2O vapour volume mixing ratio","mol/mol",       &
1931                           3,vmr)
1932            endif
1933
1934         endif ! end of 'tracer'
1935
1936         if(watercond.and.CLFvarying)then
1937            call wstats(ngrid,"rneb_man","H2O cloud fraction (conv)"," ",3,rneb_man)
1938            call wstats(ngrid,"rneb_lsc","H2O cloud fraction (large scale)"," ",3,rneb_lsc)
1939            call wstats(ngrid,"CLF","H2O cloud fraction"," ",3,cloudfrac)
1940            call wstats(ngrid,"CLFt","H2O column cloud fraction"," ",2,totcloudfrac)
1941            call wstats(ngrid,"RH","relative humidity"," ",3,RH)
1942         endif
1943
1944         if (ok_slab_ocean) then
1945            call wstats(ngrid,"dt_hdiff1","dt_hdiff1","K/s",2,dt_hdiff(:,1))
1946            call wstats(ngrid,"dt_hdiff2","dt_hdiff2","K/s",2,dt_hdiff(:,2))
1947            call wstats(ngrid,"dt_ekman1","dt_ekman1","K/s",2,dt_ekman(:,1))
1948            call wstats(ngrid,"dt_ekman2","dt_ekman2","K/s",2,dt_ekman(:,2))
1949            call wstats(ngrid,"tslab1","tslab1","K",2,tslab(:,1))
1950            call wstats(ngrid,"tslab2","tslab2","K",2,tslab(:,2))
1951            call wstats(ngrid,"pctsrf_sic","pct ice/sea","",2,pctsrf_sic)
1952            call wstats(ngrid,"tsea_ice","tsea_ice","K",2,tsea_ice)
1953            call wstats(ngrid,"sea_ice","sea ice","kg/m2",2,sea_ice)
1954            call wstats(ngrid,"rnat","nature of the surface","",2,rnat)
1955         endif
1956
1957         if(lastcall) then
1958            write (*,*) "Writing stats..."
1959            call mkstats(ierr)
1960         endif
1961         
1962      endif ! end of 'callstats'
1963
1964#ifndef MESOSCALE
1965       
1966!-----------------------------------------------------------------------------------------------------
1967!           OUTPUT in netcdf file "DIAGFI.NC", containing any variable for diagnostic
1968!
1969!             Note 1 : output with  period "ecritphy", set in "run.def"
1970!
1971!             Note 2 : writediagfi can also be called from any other subroutine for any variable,
1972!                      but its preferable to keep all the calls in one place ...
1973!-----------------------------------------------------------------------------------------------------
1974
1975      call writediagfi(ngrid,"Ls","solar longitude","deg",0,zls*180./pi)
1976      call writediagfi(ngrid,"Lss","sub solar longitude","deg",0,zlss*180./pi)
1977      call writediagfi(ngrid,"RA","right ascension","deg",0,right_ascen*180./pi)
1978      call writediagfi(ngrid,"Declin","solar declination","deg",0,declin*180./pi)
1979      call writediagfi(ngrid,"tsurf","Surface temperature","K",2,tsurf)
1980      call writediagfi(ngrid,"ps","Surface pressure","Pa",2,ps)
1981      call writediagfi(ngrid,"temp","temperature","K",3,zt)
1982      call writediagfi(ngrid,"teta","potential temperature","K",3,zh)
1983      call writediagfi(ngrid,"u","Zonal wind","m.s-1",3,zu)
1984      call writediagfi(ngrid,"v","Meridional wind","m.s-1",3,zv)
1985      call writediagfi(ngrid,"w","Vertical wind","m.s-1",3,pw)
1986      call writediagfi(ngrid,"p","Pressure","Pa",3,pplay)
1987
1988!     Subsurface temperatures
1989!        call writediagsoil(ngrid,"tsurf","Surface temperature","K",2,tsurf)
1990!        call writediagsoil(ngrid,"temp","temperature","K",3,tsoil)
1991
1992      ! Oceanic layers
1993      if(ok_slab_ocean) then
1994         call writediagfi(ngrid,"pctsrf_sic","pct ice/sea","",2,pctsrf_sic)
1995         call writediagfi(ngrid,"tsea_ice","tsea_ice","K",2,tsea_ice)
1996         call writediagfi(ngrid,"sea_ice","sea ice","kg/m2",2,sea_ice)
1997         call writediagfi(ngrid,"tslab1","tslab1","K",2,tslab(:,1))
1998         call writediagfi(ngrid,"tslab2","tslab2","K",2,tslab(:,2))
1999         call writediagfi(ngrid,"dt_hdiff1","dt_hdiff1","K/s",2,dt_hdiff(:,1))
2000         call writediagfi(ngrid,"dt_hdiff2","dt_hdiff2","K/s",2,dt_hdiff(:,2))
2001         call writediagfi(ngrid,"dt_ekman1","dt_ekman1","K/s",2,dt_ekman(:,1))
2002         call writediagfi(ngrid,"dt_ekman2","dt_ekman2","K/s",2,dt_ekman(:,2))
2003         call writediagfi(ngrid,"rnat","nature of the surface","",2,rnat)
2004         call writediagfi(ngrid,"sensibFlux","sensible heat flux","w.m^-2",2,sensibFlux)
2005         call writediagfi(ngrid,"latentFlux","latent heat flux","w.m^-2",2,zdqsdif(:,igcm_h2o_vap)*RLVTT)
2006      endif
2007     
2008
2009      ! Total energy balance diagnostics
2010      if(callrad.and.(.not.newtonian))then
2011     
2012         !call writediagfi(ngrid,"ALB","Surface albedo"," ",2,albedo_equivalent)
2013         !call writediagfi(ngrid,"ALB_1st","First Band Surface albedo"," ",2,albedo(:,1))
2014         call writediagfi(ngrid,"ISR","incoming stellar rad.","W m-2",2,fluxtop_dn)
2015         call writediagfi(ngrid,"ASR","absorbed stellar rad.","W m-2",2,fluxabs_sw)
2016         call writediagfi(ngrid,"OLR","outgoing longwave rad.","W m-2",2,fluxtop_lw)
2017         call writediagfi(ngrid,"shad","rings"," ", 2, fract)
2018         
2019!           call writediagfi(ngrid,"ASRcs","absorbed stellar rad (cs).","W m-2",2,fluxabs_sw1)
2020!           call writediagfi(ngrid,"OLRcs","outgoing longwave rad (cs).","W m-2",2,fluxtop_lw1)
2021!           call writediagfi(ngrid,"fluxsurfsw","sw surface flux.","W m-2",2,fluxsurf_sw)
2022!           call writediagfi(ngrid,"fluxsurflw","lw back radiation.","W m-2",2,fluxsurf_lw)
2023!           call writediagfi(ngrid,"fluxsurfswcs","sw surface flux (cs).","W m-2",2,fluxsurf_sw1)
2024!           call writediagfi(ngrid,"fluxsurflwcs","lw back radiation (cs).","W m-2",2,fluxsurf_lw1)
2025
2026         if(ok_slab_ocean) then
2027            call writediagfi(ngrid,"GND","heat flux from ground","W m-2",2,fluxgrdocean)
2028         else
2029            call writediagfi(ngrid,"GND","heat flux from ground","W m-2",2,fluxgrd)
2030         endif
2031         
2032         call writediagfi(ngrid,"DYN","dynamical heat input","W m-2",2,fluxdyn)
2033         
2034      endif ! end of 'callrad'
2035       
2036      if(enertest) then
2037     
2038         if (calldifv) then
2039         
2040            call writediagfi(ngrid,"q2","turbulent kinetic energy","J.kg^-1",3,q2)
2041            call writediagfi(ngrid,"sensibFlux","sensible heat flux","w.m^-2",2,sensibFlux)
2042           
2043!             call writediagfi(ngrid,"dEzdiff","turbulent diffusion heating (-sensible flux)","w.m^-2",3,dEzdiff)
2044!             call writediagfi(ngrid,"dEdiff","integrated turbulent diffusion heating (-sensible flux)","w.m^-2",2,dEdiff)
2045!             call writediagfi(ngrid,"dEdiffs","In TurbDiff (correc rad+latent heat) surf nrj change","w.m^-2",2,dEdiffs)
2046             
2047         endif
2048         
2049         if (corrk) then
2050            call writediagfi(ngrid,"dEzradsw","radiative heating","w.m^-2",3,dEzradsw)
2051            call writediagfi(ngrid,"dEzradlw","radiative heating","w.m^-2",3,dEzradlw)
2052         endif
2053         
2054         if(watercond) then
2055
2056            call writediagfi(ngrid,"lscaledE","heat from largescale","W m-2",2,lscaledE)
2057            call writediagfi(ngrid,"madjdE","heat from moistadj","W m-2",2,madjdE)
2058            call writediagfi(ngrid,"qsatatm","atm qsat"," ",3,qsat)
2059             
2060!             call writediagfi(ngrid,"lscaledEz","heat from largescale","W m-2",3,lscaledEz)
2061!             call writediagfi(ngrid,"madjdEz","heat from moistadj","W m-2",3,madjdEz)             
2062!             call writediagfi(ngrid,"h2o_max_col","maximum H2O column amount","kg.m^-2",2,H2Omaxcol)
2063
2064         endif
2065         
2066      endif ! end of 'enertest'
2067
2068
2069        ! Temporary inclusions for heating diagnostics.
2070        !call writediagfi(ngrid,"zdtsw","SW heating","T s-1",3,zdtsw)
2071        !call writediagfi(ngrid,"zdtlw","LW heating","T s-1",3,zdtlw)
2072        !call writediagfi(ngrid,"dtrad","radiative heating","K s-1",3,dtrad)
2073        ! call writediagfi(ngrid,"zdtdyn","Dyn. heating","T s-1",3,zdtdyn)
2074       
2075        ! For Debugging.
2076        !call writediagfi(ngrid,'rnat','Terrain type',' ',2,real(rnat))
2077        !call writediagfi(ngrid,'pphi','Geopotential',' ',3,pphi)
2078       
2079
2080      ! Output aerosols.
2081      if (igcm_co2_ice.ne.0.and.iaero_co2.ne.0) &
2082         call writediagfi(ngrid,'CO2ice_reff','CO2ice_reff','m',3,reffrad(1,1,iaero_co2))
2083      if (igcm_h2o_ice.ne.0.and.iaero_h2o.ne.0) &
2084         call writediagfi(ngrid,'H2Oice_reff','H2Oice_reff','m',3,reffrad(:,:,iaero_h2o))
2085      if (igcm_co2_ice.ne.0.and.iaero_co2.ne.0) &
2086         call writediagfi(ngrid,'CO2ice_reffcol','CO2ice_reffcol','um kg m^-2',2,reffcol(1,iaero_co2))
2087      if (igcm_h2o_ice.ne.0.and.iaero_h2o.ne.0) &
2088         call writediagfi(ngrid,'H2Oice_reffcol','H2Oice_reffcol','um kg m^-2',2,reffcol(1,iaero_h2o))
2089
2090      ! Output tracers.
2091      if (tracer) then
2092     
2093         do iq=1,nq
2094            call writediagfi(ngrid,noms(iq),noms(iq),'kg/kg',3,zq(1,1,iq))
2095            call writediagfi(ngrid,trim(noms(iq))//'_surf',trim(noms(iq))//'_surf',  &
2096                             'kg m^-2',2,qsurf_hist(1,iq) )
2097            call writediagfi(ngrid,trim(noms(iq))//'_col',trim(noms(iq))//'_col',    &
2098                            'kg m^-2',2,qcol(1,iq) )
2099                         
2100!          call writediagfi(ngrid,trim(noms(iq))//'_surf',trim(noms(iq))//'_surf',  &
2101!                          'kg m^-2',2,qsurf(1,iq) )                         
2102
2103            if(watercond.or.CLFvarying)then
2104               call writediagfi(ngrid,"rneb_man","H2O cloud fraction (conv)"," ",3,rneb_man)
2105               call writediagfi(ngrid,"rneb_lsc","H2O cloud fraction (large scale)"," ",3,rneb_lsc)
2106               call writediagfi(ngrid,"CLF","H2O cloud fraction"," ",3,cloudfrac)
2107               call writediagfi(ngrid,"CLFt","H2O column cloud fraction"," ",2,totcloudfrac)
2108               call writediagfi(ngrid,"RH","relative humidity"," ",3,RH)
2109            endif
2110
2111            if(waterrain)then
2112               call writediagfi(ngrid,"rain","rainfall","kg m-2 s-1",2,zdqsrain)
2113               call writediagfi(ngrid,"snow","snowfall","kg m-2 s-1",2,zdqssnow)
2114               call writediagfi(ngrid,"reevap","reevaporation of precipitation","kg m-2 s-1",2,reevap_precip)
2115            endif
2116
2117            if((hydrology).and.(.not.ok_slab_ocean))then
2118               call writediagfi(ngrid,"hice","oceanic ice height","m",2,hice)
2119            endif
2120
2121            if(ice_update)then
2122               call writediagfi(ngrid,"ice_min","min annual ice","m",2,ice_min)
2123               call writediagfi(ngrid,"ice_ini","initial annual ice","m",2,ice_initial)
2124            endif
2125
2126            do ig=1,ngrid
2127               if(tau_col(ig).gt.1.e3)then
2128                  print*,'WARNING: tau_col=',tau_col(ig)
2129                  print*,'at ig=',ig,'in PHYSIQ'
2130               endif         
2131            end do
2132
2133            call writediagfi(ngrid,"tau_col","Total aerosol optical depth","[]",2,tau_col)
2134
2135         enddo ! end of 'nq' loop
2136         
2137       endif ! end of 'tracer'
2138
2139
2140      ! Output spectrum.
2141      if(specOLR.and.corrk)then     
2142         call writediagspecIR(ngrid,"OLR3D","OLR(lon,lat,band)","W/m^2/cm^-1",3,OLR_nu)
2143         call writediagspecVI(ngrid,"OSR3D","OSR(lon,lat,band)","W/m^2/cm^-1",3,OSR_nu)
2144      endif
2145
2146#else
2147      comm_HR_SW(1:ngrid,1:nlayer) = zdtsw(1:ngrid,1:nlayer)
2148      comm_HR_LW(1:ngrid,1:nlayer) = zdtlw(1:ngrid,1:nlayer)
2149      comm_CLOUDFRAC(1:ngrid,1:nlayer)=cloudfrac(1:ngrid,1:nlayer)
2150      comm_TOTCLOUDFRAC(1:ngrid)=totcloudfrac(1:ngrid)
2151      comm_RAIN(1:ngrid,1:nlayer)=zdqrain(1:ngrid,1:nlayer,igcm_h2o_vap)
2152      comm_SNOW(1:ngrid,1:nlayer)=zdqrain(1:ngrid,1:nlayer,igcm_h2o_ice)
2153      comm_ALBEQ(1:ngrid)=albedo_equivalent(1:ngrid)
2154      comm_FLUXTOP_DN(1:ngrid)=fluxtop_dn(1:ngrid)
2155      comm_FLUXABS_SW(1:ngrid)=fluxabs_sw(1:ngrid)
2156      comm_FLUXTOP_LW(1:ngrid)=fluxtop_lw(1:ngrid)
2157      comm_FLUXSURF_SW(1:ngrid)=fluxsurf_sw(1:ngrid)
2158      comm_FLUXSURF_LW(1:ngrid)=fluxsurf_lw(1:ngrid)
2159      comm_FLXGRD(1:ngrid)=fluxgrd(1:ngrid)
2160      comm_LSCEZ(1:ngrid,1:nlayer)=lscaledEz(1:ngrid,1:nlayer)
2161      comm_H2OICE_REFF(1:ngrid,1:nlayer)=reffrad(1:ngrid,1:nlayer,iaero_h2o)
2162      sensibFlux(1:ngrid) = zflubid(1:ngrid) - capcal(1:ngrid)*zdtsdif(1:ngrid) !!! ????
2163#endif
2164
2165! XIOS outputs
2166#ifdef CPP_XIOS     
2167      ! Send fields to XIOS: (NB these fields must also be defined as
2168      ! <field id="..." /> in context_lmdz_physics.xml to be correctly used)
2169      CALL send_xios_field("ls",zls)
2170     
2171      CALL send_xios_field("ps",ps)
2172      CALL send_xios_field("area",cell_area)
2173
2174      CALL send_xios_field("temperature",zt)
2175      CALL send_xios_field("u",zu)
2176      CALL send_xios_field("v",zv)
2177
2178      if (lastcall.and.is_omp_master) then
2179        write(*,*) "physiq: call xios_context_finalize"
2180        call xios_context_finalize
2181      endif
2182#endif
2183
2184      icount=icount+1
2185
2186    end subroutine physiq
2187   
2188    end module physiq_mod
Note: See TracBrowser for help on using the repository browser.