source: trunk/LMDZ.PLUTO/libf/phypluto/physiq_mod.F90 @ 3232

Last change on this file since 3232 was 3232, checked in by afalco, 9 months ago

Pluto PCM:
1D functional
AF

  • Property svn:executable set to *
File size: 78.1 KB
Line 
1      module physiq_mod
2
3      implicit none
4
5      contains
6
7      subroutine physiq(ngrid,nlayer,nq,   &
8                  firstcall,lastcall,      &
9                  pday,ptime,ptimestep,    &
10                  pplev,pplay,pphi,        &
11                  pu,pv,pt,pq,             &
12                  flxw,                    &
13                  pdu,pdv,pdt,pdq,pdpsrf)
14
15!!
16      use write_field_phy, only: Writefield_phy
17!!
18      use ioipsl_getin_p_mod, only: getin_p
19      use radinc_h, only : L_NSPECTI,L_NSPECTV,naerkind, corrkdir, banddir
20      use generic_cloud_common_h, only : epsi_generic, Psat_generic
21      use gases_h, only: gnom, gfrac
22      use radcommon_h, only: sigma, glat, grav, BWNV, WNOI, DWNI, DWNV, WNOV
23      use suaer_corrk_mod, only: suaer_corrk
24      use radii_mod, only: su_aer_radii,haze_reffrad_fix
25      use aerosol_mod, only: iaero_haze, i_haze, haze_prof
26      use surfdat_h, only: phisfi, zmea, zstd, zsig, zgam, zthe, &
27                           dryness
28      use comdiurn_h, only: coslat, sinlat, coslon, sinlon
29      use comsaison_h, only: mu0, fract, dist_star, declin, right_ascen
30      use comsoil_h, only: nsoilmx, layer, mlayer, inertiedat
31      use geometry_mod, only: latitude, longitude, cell_area
32      USE comgeomfi_h, only: totarea, totarea_planet
33      USE tracer_h, only: noms, mmol, radius, rho_q, qext, &
34                          igcm_n2,&
35                          alpha_lift, alpha_devil, qextrhor, &
36                          nesp, is_chim, is_condensable,constants_epsi_generic
37      use time_phylmdz_mod, only: ecritphy, iphysiq, nday
38      use phyetat0_mod, only: phyetat0
39      use wstats_mod, only: callstats, wstats, mkstats
40      use phyredem, only: physdem0, physdem1
41      use planetwide_mod, only: planetwide_minval,planetwide_maxval,planetwide_sumval
42      use mod_phys_lmdz_para, only : is_master
43      use planete_mod, only: apoastr, periastr, year_day, peri_day, &
44                            obliquit, nres, z0
45      use comcstfi_mod, only: pi, g, rcp, r, rad, mugaz, cpp
46      use time_phylmdz_mod, only: daysec
47      use callkeys_mod, only: fast,albedo_spectral_mode, calladj, calldifv, &
48                              callrad, callsoil, nosurf, &
49                              aerohaze, corrk, diagdtau,&
50                              diurnal, enertest, fat1au, &
51                              icetstep, intheat, iradia, kastprof, &
52                              lwrite, mass_redistrib, meanOLR, &
53                              n2cond,nearn2cond, noseason_day, &
54                              season, sedimentation,generic_condensation, &
55                              aerohaze, haze_proffix, &
56                              specOLR, &
57                              startphy_file, testradtimes, &
58                              tracer, UseTurbDiff, &
59                              global1d, szangle
60      use generic_tracer_index_mod, only: generic_tracer_index
61      use check_fields_mod, only: check_physics_fields
62      use conc_mod, only: rnew, cpnew, ini_conc_mod
63      use phys_state_var_mod
64      use callcorrk_mod, only: callcorrk
65      use vdifc_mod, only: vdifc
66      use turbdiff_mod, only: turbdiff
67      use turb_mod, only : q2,sensibFlux,turb_resolved
68      use mass_redistribution_mod, only: mass_redistribution
69      use condensation_generic_mod, only: condensation_generic
70      use datafile_mod, only: datadir
71#ifndef MESOSCALE
72      use vertical_layers_mod, only: presnivs, pseudoalt
73      use mod_phys_lmdz_omp_data, ONLY: is_omp_master
74#else
75      use comm_wrf, only : comm_HR_SW, comm_HR_LW,  &
76         comm_CLOUDFRAC,comm_TOTCLOUDFRAC, comm_RH, &
77         comm_HR_DYN,     &
78         comm_DQICE,comm_DQVAP,comm_ALBEQ,          &
79         comm_FLUXTOP_DN,comm_FLUXABS_SW,           &
80         comm_FLUXTOP_LW,comm_FLUXSURF_SW,          &
81         comm_FLUXSURF_LW,comm_FLXGRD,              &
82         comm_DTRAIN,comm_DTLSC,                    &
83         comm_LATENT_HF
84
85#endif
86
87#ifdef CPP_XIOS
88      use xios_output_mod, only: initialize_xios_output, &
89                                 update_xios_timestep, &
90                                 send_xios_field
91      use wxios, only: wxios_context_init, xios_context_finalize
92#endif
93
94      implicit none
95
96
97!==================================================================
98!
99!     Purpose
100!     -------
101!     Central subroutine for all the physics parameterisations in the
102!     universal model. Originally adapted from the Mars LMDZ model.
103!
104!     The model can be run without or with tracer transport
105!     depending on the value of "tracer" in file "callphys.def".
106!
107!
108!   It includes:
109!
110!      I. Initialization :
111!         I.1 Firstcall initializations.
112!         I.2 Initialization for every call to physiq.
113!
114!      II. Compute radiative transfer tendencies (longwave and shortwave) :
115!         II.a Option 1 : Call correlated-k radiative transfer scheme.
116!         II.b Option 2 : Call Newtonian cooling scheme.
117!         II.c Option 3 : Atmosphere has no radiative effect.
118!
119!      III. Vertical diffusion (turbulent mixing) :
120!
121!      IV. Convection :
122!         IV.a Thermal plume model
123!         IV.b Dry convective adjusment
124!
125!      V. Condensation and sublimation of gases (currently just N2).
126!
127!      VI. Tracers
128!         VI.1. Water and water ice.
129!         VI.2  Photochemistry
130!         VI.3. Aerosols and particles.
131!         VI.4. Updates (pressure variations, surface budget).
132!         VI.5. Slab Ocean.
133!         VI.6. Surface Tracer Update.
134!
135!      VII. Surface and sub-surface soil temperature.
136!
137!      VIII. Perform diagnostics and write output files.
138!
139!
140!   arguments
141!   ---------
142!
143!   INPUT
144!   -----
145!
146!    ngrid                 Size of the horizontal grid.
147!    nlayer                Number of vertical layers.
148!    nq                    Number of advected fields.
149!
150!    firstcall             True at the first call.
151!    lastcall              True at the last call.
152!
153!    pday                  Number of days counted from the North. Spring equinoxe.
154!    ptime                 Universal time (0<ptime<1): ptime=0.5 at 12:00 UT.
155!    ptimestep             timestep (s).
156!
157!    pplay(ngrid,nlayer)   Pressure at the middle of the layers (Pa).
158!    pplev(ngrid,nlayer+1) Intermediate pressure levels (Pa).
159!    pphi(ngrid,nlayer)    Geopotential at the middle of the layers (m2.s-2).
160!
161!    pu(ngrid,nlayer)      u, zonal component of the wind (ms-1).
162!    pv(ngrid,nlayer)      v, meridional component of the wind (ms-1).
163!
164!    pt(ngrid,nlayer)      Temperature (K).
165!
166!    pq(ngrid,nlayer,nq)   Advected fields.
167!
168!    pudyn(ngrid,nlayer)    \
169!    pvdyn(ngrid,nlayer)     \ Dynamical temporal derivative for the
170!    ptdyn(ngrid,nlayer)     / corresponding variables.
171!    pqdyn(ngrid,nlayer,nq) /
172!    flxw(ngrid,nlayer)      vertical mass flux (kg/s) at layer lower boundary
173!
174!   OUTPUT
175!   ------
176!
177!    pdu(ngrid,nlayer)        \
178!    pdv(ngrid,nlayer)         \  Temporal derivative of the corresponding
179!    pdt(ngrid,nlayer)         /  variables due to physical processes.
180!    pdq(ngrid,nlayer)        /
181!    pdpsrf(ngrid)             /
182!
183!
184!     Authors
185!     -------
186!           Frederic Hourdin        15/10/93
187!           Francois Forget        1994
188!           Christophe Hourdin        02/1997
189!           Subroutine completely rewritten by F. Forget (01/2000)
190!           Water ice clouds: Franck Montmessin (update 06/2003)
191!           Radiatively active tracers: J.-B. Madeleine (10/2008-06/2009)
192!           New correlated-k radiative scheme: R. Wordsworth (2009)
193!           Many specifically Martian subroutines removed: R. Wordsworth (2009)
194!           Improved water cycle: R. Wordsworth / B. Charnay (2010)
195!           To F90: R. Wordsworth (2010)
196!           New turbulent diffusion scheme: J. Leconte (2012)
197!           Loops converted to F90 matrix format: J. Leconte (2012)
198!           No more ngrid/nq, F90 commons and adaptation to parallel: A. Spiga (2012)
199!           Purge of the code : M. Turbet (2015)
200!           Photochemical core developped by F. Lefevre: B. Charnay (2017)
201!           Purge for Pluto model : A. Falco (2024)
202!==================================================================
203
204
205!    0.  Declarations :
206!    ------------------
207
208include "netcdf.inc"
209
210! Arguments :
211! -----------
212
213!   INPUTS:
214!   -------
215
216      integer,intent(in) :: ngrid             ! Number of atmospheric columns.
217      integer,intent(in) :: nlayer            ! Number of atmospheric layers.
218      integer,intent(in) :: nq                ! Number of tracers.
219
220      logical,intent(in) :: firstcall ! Signals first call to physics.
221      logical,intent(in) :: lastcall  ! Signals last call to physics.
222
223      real,intent(in) :: pday                  ! Number of elapsed sols since reference Ls=0.
224      real,intent(in) :: ptime                 ! "Universal time", given as fraction of sol (e.g.: 0.5 for noon).
225      real,intent(in) :: ptimestep             ! Physics timestep (s).
226      real,intent(in) :: pplev(ngrid,nlayer+1) ! Inter-layer pressure (Pa).
227      real,intent(in) :: pplay(ngrid,nlayer)   ! Mid-layer pressure (Pa).
228      real,intent(in) :: pphi(ngrid,nlayer)    ! Geopotential at mid-layer (m2s-2).
229      real,intent(in) :: pu(ngrid,nlayer)      ! Zonal wind component (m/s).
230      real,intent(in) :: pv(ngrid,nlayer)      ! Meridional wind component (m/s).
231      real,intent(in) :: pt(ngrid,nlayer)      ! Temperature (K).
232      real,intent(in) :: pq(ngrid,nlayer,nq)   ! Tracers (kg/kg_of_air).
233      real,intent(in) :: flxw(ngrid,nlayer)    ! Vertical mass flux (ks/s) at lower boundary of layer
234
235!   OUTPUTS:
236!   --------
237
238!     Physical tendencies :
239
240      real,intent(out) :: pdu(ngrid,nlayer)    ! Zonal wind tendencies (m/s/s).
241      real,intent(out) :: pdv(ngrid,nlayer)    ! Meridional wind tendencies (m/s/s).
242      real,intent(out) :: pdt(ngrid,nlayer)    ! Temperature tendencies (K/s).
243      real,intent(out) :: pdq(ngrid,nlayer,nq) ! Tracer tendencies (kg/kg_of_air/s).
244      real,intent(out) :: pdpsrf(ngrid)        ! Surface pressure tendency (Pa/s).
245
246! Local saved variables:
247! ----------------------
248      integer,save :: day_ini                                      ! Initial date of the run (sol since Ls=0).
249      integer,save :: icount                                       ! Counter of calls to physiq during the run.
250!$OMP THREADPRIVATE(day_ini,icount)
251
252! Local variables :
253! -----------------
254
255!     Aerosol (dust or ice) extinction optical depth  at reference wavelength
256!     for the "naerkind" optically active aerosols:
257
258      real,save,allocatable :: aerosol(:,:,:) ! Aerosols
259!$OMP THREADPRIVATE(aerosol)
260      real zh(ngrid,nlayer)               ! Potential temperature (K).
261      real pw(ngrid,nlayer)               ! Vertical velocity (m/s). (NOTE : >0 WHEN DOWNWARDS !!)
262      real omega(ngrid,nlayer)            ! omega velocity (Pa/s, >0 when downward)
263
264      integer l,ig,ierr,iq,nw,isoil,iesp, igcm_generic_vap, igcm_generic_ice
265      logical call_ice_vap_generic ! to call only one time the ice/vap pair of a tracer
266
267      real zls                       ! Solar longitude (radians).
268      real zlss                      ! Sub solar point longitude (radians).
269      real zday                      ! Date (time since Ls=0, calculated in sols).
270      real zzlay(ngrid,nlayer)       ! Altitude at the middle of the atmospheric layers.
271      real zzlev(ngrid,nlayer+1)     ! Altitude at the atmospheric layer boundaries.
272
273! VARIABLES for the thermal plume model (AF24: deleted)
274
275! TENDENCIES due to various processes :
276
277      ! For Surface Temperature : (K/s)
278      real zdtsurf(ngrid)     ! Cumulated tendencies.
279      real zdtsurfmr(ngrid)   ! Mass_redistribution routine.
280      real zdtsurfc(ngrid)    ! Condense_n2 routine.
281      real zdtsdif(ngrid)     ! Turbdiff/vdifc routines.
282      ! real zdtsurf_hyd(ngrid) ! Hydrol routine.
283
284      ! For Atmospheric Temperatures : (K/s)
285      real dtlscale(ngrid,nlayer)                             ! Largescale routine.
286      real dt_generic_condensation(ngrid,nlayer)              ! condensation_generic routine.
287      real zdtc(ngrid,nlayer)                                 ! Condense_n2 routine.
288      real zdtdif(ngrid,nlayer)                               ! Turbdiff/vdifc routines.
289      real zdtmr(ngrid,nlayer)                                ! Mass_redistribution routine.
290      real zdtsw1(ngrid,nlayer), zdtlw1(ngrid,nlayer)         ! Callcorrk routine.
291      real zdtchim(ngrid,nlayer)                              ! Calchim routine.
292
293      ! For Surface Tracers : (kg/m2/s)
294      real dqsurf(ngrid,nq)                 ! Cumulated tendencies.
295      !real zdqsurfc(ngrid)                  ! Condense_n2 routine.
296      REAL zdqsc(ngrid,nq)                  ! Condense_n2 routine.
297      real zdqsdif(ngrid,nq)                ! Turbdiff/vdifc routines.
298      real zdqssed(ngrid,nq)                ! Callsedim routine.
299      real zdqsurfmr(ngrid,nq)              ! Mass_redistribution routine.
300
301      ! For Tracers : (kg/kg_of_air/s)
302      real zdqc(ngrid,nlayer,nq)      ! Condense_n2 routine.
303      real zdqadj(ngrid,nlayer,nq)    ! Convadj routine.
304      real zdqdif(ngrid,nlayer,nq)    ! Turbdiff/vdifc routines.
305      real zdqevap(ngrid,nlayer)      ! Turbdiff routine.
306      real zdqsed(ngrid,nlayer,nq)    ! Callsedim routine.
307      real zdqmr(ngrid,nlayer,nq)     ! Mass_redistribution routine.
308      real dqvaplscale_generic(ngrid,nlayer,nq) ! condensation_generic routine.
309      real dqcldlscale_generic(ngrid,nlayer,nq) ! condensation_generic routine.
310      REAL,allocatable,save :: zdqchim(:,:,:) ! Calchim_asis routine
311      REAL,allocatable,save :: zdqschim(:,:)  ! Calchim_asis routine
312!$OMP THREADPRIVATE(zdqchim,zdqschim)
313
314
315      !! PLUTO variables
316      REAL zdqch4cloud(ngrid,nlayer,nq)
317      REAL zdqsch4cloud(ngrid,nq)
318      REAL zdtch4cloud(ngrid,nlayer)
319      REAL zdqcocloud(ngrid,nlayer,nq)
320      REAL zdqscocloud(ngrid,nq)
321      REAL zdtcocloud(ngrid,nlayer)
322      REAL rice_ch4(ngrid,nlayer)    ! Methane ice geometric mean radius (m)
323      REAL rice_co(ngrid,nlayer)     ! CO ice geometric mean radius (m)
324
325      REAL zdqsch4fast(ngrid)    ! used only for fast model nogcm
326      REAL zdqch4fast(ngrid)    ! used only for fast model nogcm
327      REAL zdqscofast(ngrid)    ! used only for fast model nogcm
328      REAL zdqcofast(ngrid)    ! used only for fast model nogcm
329      REAL zdqflow(ngrid,nq)
330
331      REAL zdteuv(ngrid,nlayer)    ! (K/s)
332      REAL zdtconduc(ngrid,nlayer) ! (K/s)
333      REAL zdumolvis(ngrid,nlayer)
334      REAL zdvmolvis(ngrid,nlayer)
335      real zdqmoldiff(ngrid,nlayer,nq)
336
337      ! Haze relatated tendancies
338      REAL zdqhaze(ngrid,nlayer,nq)
339      REAL zdqprodhaze(ngrid,nq)
340      REAL zdqsprodhaze(ngrid)
341      REAL zdqhaze_col(ngrid)
342      REAL zdqphot_prec(ngrid,nlayer)
343      REAL zdqphot_ch4(ngrid,nlayer)
344      REAL zdqconv_prec(ngrid,nlayer)
345      REAL zdq_source(ngrid,nlayer,nq)
346      ! Fast Haze relatated tendancies
347      REAL fluxbot(ngrid)
348      REAL gradflux(ngrid)
349      REAL fluxlym_sol_bot(ngrid)      ! Solar flux Lyman alpha ph.m-2.s-1 reaching the surface
350      REAL fluxlym_ipm_bot(ngrid)      ! IPM (Interplanetary) flux Lyman alpha ph.m-2.s-1 reaching the surface
351      REAL flym_sol(ngrid)      ! Incident Solar flux Lyman alpha ph.m-2.s-1
352      REAL flym_ipm(ngrid)      ! Incident IPM (Interplanetary) flux Lyman alpha ph.m-2.s-1
353
354
355
356      REAL array_zero1(ngrid)
357      REAL array_zero2(ngrid,nlayer)
358
359      ! For Winds : (m/s/s)
360      real zdvadj(ngrid,nlayer), zduadj(ngrid,nlayer)       ! Convadj routine.
361      real zdumr(ngrid,nlayer), zdvmr(ngrid,nlayer)         ! Mass_redistribution routine.
362      real zdvdif(ngrid,nlayer), zdudif(ngrid,nlayer)       ! Turbdiff/vdifc routines.
363      real zdhdif(ngrid,nlayer)                             ! Turbdiff/vdifc routines.
364      real zdhadj(ngrid,nlayer)                             ! Convadj routine.
365      REAL zdvc(ngrid,nlayer),zduc(ngrid,nlayer)            ! condense_n2 routine.
366
367      ! For Pressure and Mass :
368      real zdmassmr(ngrid,nlayer) ! Atmospheric Mass tendency for mass_redistribution (kg_of_air/m2/s).
369      real zdmassmr_col(ngrid)    ! Atmospheric Column Mass tendency for mass_redistribution (kg_of_air/m2/s).
370      real zdpsrfmr(ngrid)        ! Pressure tendency for mass_redistribution routine (Pa/s).
371
372
373
374! Local variables for LOCAL CALCULATIONS:
375! ---------------------------------------
376      real zflubid(ngrid)
377      real zplanck(ngrid),zpopsk(ngrid,nlayer)
378      real ztim1,ztim2,ztim3, z1,z2
379      real ztime_fin
380      real zdh(ngrid,nlayer)
381      real gmplanet
382      real taux(ngrid),tauy(ngrid)
383
384
385! local variables for DIAGNOSTICS : (diagfi & stat)
386! -------------------------------------------------
387      real ps(ngrid)                                     ! Surface Pressure.
388      real zt(ngrid,nlayer)                              ! Atmospheric Temperature.
389      real zu(ngrid,nlayer),zv(ngrid,nlayer)             ! Zonal and Meridional Winds.
390      real zq(ngrid,nlayer,nq)                           ! Atmospheric Tracers.
391      real zdtadj(ngrid,nlayer)                          ! Convadj Diagnostic.
392      real zdtdyn(ngrid,nlayer)                          ! Dynamical Heating (K/s).
393      real zdudyn(ngrid,nlayer)                          ! Dynamical Zonal Wind tendency (m.s-2).
394
395      real reff(ngrid,nlayer)                       ! Effective dust radius (used if doubleq=T).
396      real vmr(ngrid,nlayer)                        ! volume mixing ratio
397      real time_phys
398
399      real ISR,ASR,OLR,GND,DYN,GSR,Ts1,Ts2,Ts3,TsS ! for Diagnostic.
400
401      real qcol(ngrid,nq) ! Tracer Column Mass (kg/m2).
402
403      !     Pluto adding variables
404      real vmr_ch4(ngrid)  ! vmr ch4
405      real vmr_co(ngrid)  ! vmr co
406      real rho(ngrid,nlayer) ! density
407      real zrho_ch4(ngrid,nlayer) ! density methane kg.m-3
408      real zrho_co(ngrid,nlayer) ! density CO kg.m-3
409      real zrho_haze(ngrid,nlayer) ! density haze kg.m-3
410      real zdqrho_photprec(ngrid,nlayer) !photolysis rate kg.m-3.s-1
411      real zq1temp_ch4(ngrid) !
412      real qsat_ch4(ngrid) !
413      real qsat_ch4_l1(ngrid) !
414!      CHARACTER(LEN=20) :: txt ! to temporarly store text for eddy tracers
415      real profmmr(ngrid,nlayer) ! fixed profile of haze if haze_proffix
416      real sensiblehf1(ngrid) ! sensible heat flux
417      real sensiblehf2(ngrid) ! sensible heat flux
418
419!     included by RW for H2O Manabe scheme
420      real rneb_man(ngrid,nlayer) ! H2O cloud fraction (moistadj).
421      real rneb_lsc(ngrid,nlayer) ! H2O cloud fraction (large scale).
422
423
424!     to test energy conservation (RW+JL)
425      real mass(ngrid,nlayer),massarea(ngrid,nlayer)
426      real dEtot, dEtots, AtmToSurf_TurbFlux
427      real,save :: dEtotSW, dEtotsSW, dEtotLW, dEtotsLW
428!$OMP THREADPRIVATE(dEtotSW, dEtotsSW, dEtotLW, dEtotsLW)
429
430!JL12 conservation test for mean flow kinetic energy has been disabled temporarily
431
432      real dtmoist_max,dtmoist_min
433      real dItot, dItot_tmp, dVtot, dVtot_tmp
434
435
436      real dWtot, dWtot_tmp, dWtots, dWtots_tmp
437
438      real psat_tmp ! AF24: to remove?
439
440      real qsat_generic(ngrid,nlayer,nq) ! generic condensable tracers (GCS) specific concentration at saturation (kg/kg_of_air).
441      real RH_generic(ngrid,nlayer,nq)   ! generic condensable tracers (GCS) Relative humidity.
442      real rneb_generic(ngrid,nlayer,nq) ! GCS cloud fraction (generic condensation).
443      real psat_tmp_generic
444      real, save :: metallicity ! metallicity of planet --- is not used here, but necessary to call function Psat_generic
445!$OMP THREADPRIVATE(metallicity)
446
447      real reffrad_generic_zeros_for_wrf(ngrid,nlayer) !  !!! this is temporary, it is only a list of zeros, it will be replaced when a generic aerosol will be implemented
448
449      ! For Clear Sky Case. (AF24: deleted)
450
451      real nconsMAX, vdifcncons(ngrid), cadjncons(ngrid) ! Vdfic water conservation test. By RW
452
453      real muvar(ngrid,nlayer+1) ! For Runaway Greenhouse 1D study. By RW
454
455      real,save,allocatable :: reffcol(:,:)
456!$OMP THREADPRIVATE(reffcol)
457
458!  Non-oro GW tendencies
459      REAL d_u_hin(ngrid,nlayer), d_v_hin(ngrid,nlayer)
460      REAL d_t_hin(ngrid,nlayer)
461!  Diagnostics 2D of gw_nonoro
462      REAL zustrhi(ngrid), zvstrhi(ngrid)
463
464
465      real :: tsurf2(ngrid)
466!!      real :: flux_o(ngrid),flux_g(ngrid)
467      real :: flux_g(ngrid)
468      real :: flux_sens_lat(ngrid)
469      real :: qsurfint(ngrid,nq)
470#ifdef MESOSCALE
471      REAL :: lsf_dt(nlayer)
472      REAL :: lsf_dq(nlayer)
473#endif
474
475      ! flags to trigger extra sanity checks
476      logical, save ::  check_physics_inputs=.false.
477      logical, save ::  check_physics_outputs=.false.
478!$OPM THREADPRIVATE(check_physics_inputs,check_physics_outputs)
479
480      ! Misc
481      character*2 :: str2
482      character(len=10) :: tmp1
483      character(len=10) :: tmp2
484!==================================================================================================
485
486! -----------------
487! I. INITIALISATION
488! -----------------
489
490! --------------------------------
491! I.1   First Call Initialisation.
492! --------------------------------
493      if (firstcall) then
494         call getin_p("check_physics_inputs", check_physics_inputs)
495         call getin_p("check_physics_outputs", check_physics_outputs)
496
497         ! Allocate saved arrays (except for 1D model, where this has already
498         ! been done)
499#ifndef MESOSCALE
500         if (ngrid>1) call phys_state_var_init(nq)
501#endif
502
503!        Variables set to 0
504!        ~~~~~~~~~~~~~~~~~~
505         dtrad(:,:) = 0.0
506         fluxrad(:) = 0.0
507         tau_col(:) = 0.0
508         zdtsw(:,:) = 0.0
509         zdtlw(:,:) = 0.0
510
511!        Initialize tracer names, indexes and properties.
512!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
513         IF (.NOT.ALLOCATED(noms)) ALLOCATE(noms(nq)) ! (because noms is an argument of physdem1 whether or not tracer is on)
514         if (tracer) then
515            call initracer(ngrid,nq)
516            ! if(photochem) then !AF24: removed
517         endif
518!        Initialize aerosol indexes.
519!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~
520         ! call iniaerosol
521         ! allocate related local arrays
522         ! (need be allocated instead of automatic because of "naerkind")
523         allocate(aerosol(ngrid,nlayer,naerkind))
524         allocate(reffcol(ngrid,naerkind))
525
526#ifdef CPP_XIOS
527         ! Initialize XIOS context
528         write(*,*) "physiq: call wxios_context_init"
529         CALL wxios_context_init
530#endif
531
532!        Read 'startfi.nc' file.
533!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
534#ifndef MESOSCALE
535         call phyetat0(startphy_file,                                 &
536                       ngrid,nlayer,"startfi.nc",0,0,nsoilmx,nq,      &
537                       day_ini,time_phys,tsurf,tsoil,emis,q2,qsurf)
538
539#else
540
541         day_ini = pday
542#endif
543
544#ifndef MESOSCALE
545         if (.not.startphy_file) then
546           ! additionnal "academic" initialization of physics
547           if (is_master) write(*,*) "Physiq: initializing tsurf(:) to pt(:,1) !!"
548           tsurf(:)=pt(:,1)
549           if (is_master) write(*,*) "Physiq: initializing tsoil(:) to pt(:,1) !!"
550           do isoil=1,nsoilmx
551             tsoil(1:ngrid,isoil)=tsurf(1:ngrid)
552           enddo
553           if (is_master) write(*,*) "Physiq: initializing day_ini to pday !"
554           day_ini=pday
555         endif
556#endif
557         if (pday.ne.day_ini) then
558             write(*,*) "ERROR in physiq.F90:"
559             write(*,*) "bad synchronization between physics and dynamics"
560             write(*,*) "dynamics day: ",pday
561             write(*,*) "physics day:  ",day_ini
562             stop
563         endif
564
565         write (*,*) 'In physiq day_ini =', day_ini
566
567!        Initialize albedo calculation.
568!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
569         albedo(:,:)=0.0
570         albedo_bareground(:)=0.0
571         albedo_snow_SPECTV(:)=0.0
572         albedo_n2_ice_SPECTV(:)=0.0
573         call surfini(ngrid,nq,qsurf,albedo,albedo_bareground,albedo_snow_SPECTV,albedo_n2_ice_SPECTV)
574
575!        Initialize orbital calculation.
576!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
577         call iniorbit(apoastr,periastr,year_day,peri_day,obliquit)
578
579
580!        Initialize soil.
581!        ~~~~~~~~~~~~~~~~
582         if (callsoil) then
583            call soil(ngrid,nsoilmx,firstcall,lastcall,inertiedat, &
584                      ptimestep,tsurf,tsoil,capcal,fluxgrd)
585         else ! else of 'callsoil'.
586            print*,'WARNING! Thermal conduction in the soil turned off'
587            capcal(:)=1.e6
588            fluxgrd(:)=intheat
589            print*,'Flux from ground = ',intheat,' W m^-2'
590         endif ! end of 'callsoil'.
591
592         icount=1
593
594!        Initialize surface history variable.
595!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
596         qsurf_hist(:,:)=qsurf(:,:)
597
598!!         call WriteField_phy("post_qsurf_hist_qsurf",qsurf(1:ngrid,igcm_h2o_vap),1)
599
600!        Initialize variable for dynamical heating and zonal wind tendency diagnostic
601!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
602         ztprevious(:,:)=pt(:,:)
603         zuprevious(:,:)=pu(:,:)
604
605!        Set temperature just above condensation temperature (for Early Mars)
606!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
607         if(nearn2cond) then
608            write(*,*)' WARNING! Starting at Tcond+1K'
609            do l=1, nlayer
610               do ig=1,ngrid
611                  pdt(ig,l)= ((-3167.8)/(log(.01*pplay(ig,l))-23.23)+4     &
612                      -pt(ig,l)) / ptimestep
613               enddo
614            enddo
615         endif
616
617         if(meanOLR)then
618            call system('rm -f rad_bal.out') ! to record global radiative balance.
619            call system('rm -f tem_bal.out') ! to record global mean/max/min temperatures.
620            call system('rm -f h2o_bal.out') ! to record global hydrological balance.
621         endif
622
623
624         !! Initialize variables for water cycle ! AF24: removed
625
626!        Set metallicity for GCS
627!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
628         metallicity=0.0 ! default value --- is not used here but necessary to call function Psat_generic
629         call getin_p("metallicity",metallicity) ! --- is not used here but necessary to call function Psat_generic
630
631!        Set some parameters for the thermal plume model !AF24: removed
632!        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
633
634#ifndef MESOSCALE
635         if (ngrid.ne.1) then ! Note : no need to create a restart file in 1d.
636            call physdem0("restartfi.nc",longitude,latitude,nsoilmx,ngrid,nlayer,nq, &
637                         ptimestep,pday+nday,time_phys,cell_area,          &
638                         albedo_bareground,inertiedat,zmea,zstd,zsig,zgam,zthe)
639         endif
640
641!!         call WriteField_phy("post_physdem_qsurf",qsurf(1:ngrid,igcm_h2o_vap),1)
642#endif
643         if (corrk) then
644            ! We initialise the spectral grid here instead of
645            ! at firstcall of callcorrk so we can output XspecIR, XspecVI
646            ! when using Dynamico
647            print*, "physiq_mod: Correlated-k data base folder:",trim(datadir)
648            call getin_p("corrkdir",corrkdir)
649            print*,"corrkdir = ", corrkdir
650            write (tmp1, '(i4)') L_NSPECTI
651            write (tmp2, '(i4)') L_NSPECTV
652            banddir=trim(trim(adjustl(tmp1))//'x'//trim(adjustl(tmp2)))
653            banddir=trim(trim(adjustl(corrkdir))//'/'//trim(adjustl(banddir)))
654            call setspi !Basic infrared properties.
655            call setspv ! Basic visible properties.
656            call sugas_corrk       ! Set up gaseous absorption properties.
657            if (aerohaze) then
658               call suaer_corrk       ! Set up aerosol optical properties.
659            endif
660         endif
661
662!!         call WriteField_phy("post_corrk_firstcall_qsurf",qsurf(1:ngrid,igcm_h2o_vap),1)
663         ! XIOS outputs
664#ifdef CPP_XIOS
665
666         write(*,*) "physiq: call initialize_xios_output"
667         call initialize_xios_output(pday,ptime,ptimestep,daysec, &
668                                     year_day,presnivs,pseudoalt,WNOI,WNOV)
669#endif
670
671!!         call WriteField_phy("post_xios_qsurf",qsurf(1:ngrid,igcm_h2o_vap),1)
672
673         write(*,*) "physiq: end of firstcall"
674      endif ! end of 'firstcall'
675
676!!      call WriteField_phy("post_firstcall_qsurf",qsurf(1:ngrid,igcm_h2o_vap),1)
677!!      call writediagfi(ngrid,"firstcall_post_qsurf"," "," ",2,qsurf(1:ngrid,igcm_h2o_vap))
678
679      if (check_physics_inputs) then
680         !check the validity of input fields coming from the dynamics
681         call check_physics_fields("begin physiq:", pt, pu, pv, pplev, pq)
682      endif
683
684!      call writediagfi(ngrid,"pre_physical_rnat"," "," ",2,rnat)
685!      call writediagfi(ngrid,"pre_physical_capcal"," "," ",2,capcal)
686
687! ------------------------------------------------------
688! I.2   Initializations done at every physical timestep:
689! ------------------------------------------------------
690
691#ifdef CPP_XIOS
692      ! update XIOS time/calendar
693      call update_xios_timestep
694#endif
695
696      ! Initialize various variables
697      ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
698
699      if ( .not.nearn2cond ) then
700         pdt(1:ngrid,1:nlayer) = 0.0
701      endif
702      zdtsurf(1:ngrid)      = 0.0
703      pdq(1:ngrid,1:nlayer,1:nq) = 0.0
704      dqsurf(1:ngrid,1:nq)= 0.0
705      pdu(1:ngrid,1:nlayer) = 0.0
706      pdv(1:ngrid,1:nlayer) = 0.0
707      pdpsrf(1:ngrid)       = 0.0
708      zflubid(1:ngrid)      = 0.0
709      flux_sens_lat(1:ngrid) = 0.0
710      taux(1:ngrid) = 0.0
711      tauy(1:ngrid) = 0.0
712
713      zday=pday+ptime ! Compute time, in sols (and fraction thereof).
714
715      ! Compute Stellar Longitude (Ls), and orbital parameters.
716      ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
717      if (season) then
718         call stellarlong(zday,zls)
719      else
720         call stellarlong(noseason_day,zls)
721      end if
722
723      call orbite(zls,dist_star,declin,right_ascen)
724
725      if (diurnal) then
726              zlss=-2.*pi*(zday-.5)
727      else if(diurnal .eqv. .false.) then
728              zlss=9999.
729      endif
730
731
732      glat(:) = g !AF24: removed oblateness
733
734      ! Compute geopotential between layers.
735      ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
736      zzlay(1:ngrid,1:nlayer)=pphi(1:ngrid,1:nlayer)
737      do l=1,nlayer
738         zzlay(1:ngrid,l)= zzlay(1:ngrid,l)/glat(1:ngrid)
739      enddo
740
741      zzlev(1:ngrid,1)=0.
742
743      do l=2,nlayer
744         do ig=1,ngrid
745            z1=(pplay(ig,l-1)+pplev(ig,l))/(pplay(ig,l-1)-pplev(ig,l))
746            z2=(pplev(ig,l)+pplay(ig,l))/(pplev(ig,l)-pplay(ig,l))
747            zzlev(ig,l)=(z1*zzlay(ig,l-1)+z2*zzlay(ig,l))/(z1+z2)
748         enddo
749      enddo
750
751      !Altitude of top interface (nlayer+1), using the thicknesss of the level below the top one. LT22
752
753      zzlev(1:ngrid,nlayer+1) = 2*zzlev(1:ngrid,nlayer)-zzlev(1:ngrid,nlayer-1)
754
755      ! Compute potential temperature
756      ! Note : Potential temperature calculation may not be the same in physiq and dynamic...
757      ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
758      do l=1,nlayer
759         do ig=1,ngrid
760            zpopsk(ig,l)=(pplay(ig,l)/pplev(ig,1))**rcp
761            zh(ig,l)=pt(ig,l)/zpopsk(ig,l)
762            mass(ig,l)  = (pplev(ig,l) - pplev(ig,l+1))/glat(ig)
763            massarea(ig,l)=mass(ig,l)*cell_area(ig)
764         enddo
765      enddo
766
767     ! Compute vertical velocity (m/s) from vertical mass flux
768     ! w = F / (rho*area) and rho = P/(r*T)
769     ! But first linearly interpolate mass flux to mid-layers
770      if (.not.fast) then
771       do l=1,nlayer-1
772         pw(1:ngrid,l)=0.5*(flxw(1:ngrid,l)+flxw(1:ngrid,l+1))
773       enddo
774       pw(1:ngrid,nlayer)=0.5*flxw(1:ngrid,nlayer) ! since flxw(nlayer+1)=0
775       do l=1,nlayer
776         pw(1:ngrid,l)=(pw(1:ngrid,l)*r*pt(1:ngrid,l)) /  &
777                       (pplay(1:ngrid,l)*cell_area(1:ngrid))
778       enddo
779       ! omega in Pa/s
780       do l=1,nlayer-1
781         omega(1:ngrid,l)=0.5*(flxw(1:ngrid,l)+flxw(1:ngrid,l+1))
782       enddo
783       omega(1:ngrid,nlayer)=0.5*flxw(1:ngrid,nlayer) ! since flxw(nlayer+1)=0
784       do l=1,nlayer
785         omega(1:ngrid,l)=g*omega(1:ngrid,l)/cell_area(1:ngrid)
786       enddo
787      endif
788!---------------------------------
789! II. Compute radiative tendencies
790!---------------------------------
791
792      ! Compute local stellar zenith angles
793      ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
794      if (diurnal) then
795         ztim1=SIN(declin)
796         ztim2=COS(declin)*COS(2.*pi*(zday-.5))
797         ztim3=-COS(declin)*SIN(2.*pi*(zday-.5))
798
799         call stelang(ngrid,sinlon,coslon,sinlat,coslat,    &
800                        ztim1,ztim2,ztim3,mu0,fract)
801      else if(diurnal .eqv. .false.) then
802
803         call mucorr(ngrid,declin,latitude,mu0,fract,10000.,rad)
804         ! WARNING: this function appears not to work in 1D
805
806         if ((ngrid.eq.1).and.(global1d)) then ! Fixed zenith angle 'szangle' in 1D simulations w/ globally-averaged sunlight.
807            mu0 = cos(pi*szangle/180.0)
808         endif
809
810      endif
811
812      if (callrad) then
813         if( mod(icount-1,iradia).eq.0.or.lastcall) then
814
815            ! Eclipse incoming sunlight !AF24: removed
816
817!!            call writediagfi(ngrid,"corrk_pre_dqsurf"," "," ",2,dqsurf(1:ngrid,igcm_h2o_vap))
818!!            call writediagfi(ngrid,"corrk_pre_qsurf"," "," ",2,qsurf(1:ngrid,igcm_h2o_vap))
819
820
821            if (corrk) then
822
823! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
824! II.a Call correlated-k radiative transfer scheme
825! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
826               if(kastprof)then
827                  print*,'kastprof should not = true here'
828                  call abort
829               endif
830               ! if(water) then !AF24: removed
831
832               if(generic_condensation) then
833                  do iq=1,nq
834
835                     call generic_tracer_index(nq,iq,igcm_generic_vap,igcm_generic_ice,call_ice_vap_generic)
836
837                     if (call_ice_vap_generic) then ! to call only one time the ice/vap pair of a tracer
838
839                        epsi_generic=constants_epsi_generic(iq)
840
841                        muvar(1:ngrid,1:nlayer)=mugaz/(1.e0+(1.e0/epsi_generic-1.e0)*pq(1:ngrid,1:nlayer,igcm_generic_vap))
842                        muvar(1:ngrid,nlayer+1)=mugaz/(1.e0+(1.e0/epsi_generic-1.e0)*pq(1:ngrid,nlayer,igcm_generic_vap))
843
844                     endif
845                  end do ! do iq=1,nq loop on tracers
846               ! take into account generic condensable specie (GCS) effect on mean molecular weight
847
848               else
849                  muvar(1:ngrid,1:nlayer+1)=mugaz
850               endif
851
852               ! if(ok_slab_ocean) then !AF24: removed
853
854               ! standard callcorrk
855               ! clearsky=.false.
856               call callcorrk(ngrid,nlayer,pq,nq,qsurf,                           &
857                              albedo,albedo_equivalent,emis,mu0,pplev,pplay,pt,   &
858                              tsurf,fract,dist_star,aerosol,muvar,                &
859                              zdtlw,zdtsw,fluxsurf_lw,fluxsurf_sw,                &
860                              fluxsurfabs_sw,fluxtop_lw,                          &
861                              fluxabs_sw,fluxtop_dn,OLR_nu,OSR_nu,GSR_nu,         &
862                              int_dtaui,int_dtauv,                                &
863                              tau_col,cloudfrac,totcloudfrac,                     &
864                              .false.,firstcall,lastcall)
865
866                !GG (feb2021): Option to "artificially" decrease the raditive time scale in
867                !the deep atmosphere  press > 0.1 bar. Suggested by MT.
868                !! COEFF_RAD to be "tuned" to facilitate convergence of tendency
869
870                !coeff_rad=0.   ! 0 values, it doesn't accelerate the convergence
871                !coeff_rad=0.5
872                !coeff_rad=1.
873                !do l=1, nlayer
874                !  do ig=1,ngrid
875                !    if(pplay(ig,l).ge.1.d4) then
876                !      zdtsw(ig,l)=zdtsw(ig,l)*(pplay(ig,l)/1.d4)**coeff_rad
877                !      zdtlw(ig,l)=zdtlw(ig,l)*(pplay(ig,l)/1.d4)**coeff_rad
878                !    endif
879                !  enddo
880                !enddo
881
882                ! AF24: removed CLFvarying Option
883
884               ! if(ok_slab_ocean) then
885               !    tsurf(:)=tsurf2(:)
886               ! endif
887
888
889               ! Radiative flux from the sky absorbed by the surface (W.m-2).
890               GSR=0.0
891               fluxrad_sky(1:ngrid)=emis(1:ngrid)*fluxsurf_lw(1:ngrid)+fluxsurfabs_sw(1:ngrid)
892
893                            !if(noradsurf)then ! no lower surface; SW flux just disappears
894                            !   GSR = SUM(fluxsurf_sw(1:ngrid)*cell_area(1:ngrid))/totarea
895                            !   fluxrad_sky(1:ngrid)=emis(1:ngrid)*fluxsurf_lw(1:ngrid)
896                            !   print*,'SW lost in deep atmosphere = ',GSR,' W m^-2'
897                            !endif
898
899               ! Net atmospheric radiative heating rate (K.s-1)
900               dtrad(1:ngrid,1:nlayer)=zdtsw(1:ngrid,1:nlayer)+zdtlw(1:ngrid,1:nlayer)
901
902               ! Late initialization of the Ice Spectral Albedo. We needed the visible bands to do that !
903               if (firstcall .and. albedo_spectral_mode) then
904                  call spectral_albedo_calc(albedo_snow_SPECTV)
905               endif
906
907! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
908! II.b Call Newtonian cooling scheme !AF24: removed
909! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
910            else
911! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
912! II.c Atmosphere has no radiative effect
913! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
914               fluxtop_dn(1:ngrid)  = fract(1:ngrid)*mu0(1:ngrid)*Fat1AU/dist_star**2
915               if(ngrid.eq.1)then ! / by 4 globally in 1D case...
916                  fluxtop_dn(1)  = fract(1)*Fat1AU/dist_star**2/2.0
917               endif
918               fluxsurf_sw(1:ngrid) = fluxtop_dn(1:ngrid)
919               print*,'------------WARNING---WARNING------------' ! by MT2015.
920               print*,'You are in corrk=false mode, '
921               print*,'and the surface albedo is taken equal to the first visible spectral value'
922
923               albedo_equivalent(1:ngrid)=albedo(1:ngrid,1)
924               fluxsurfabs_sw(1:ngrid) = fluxtop_dn(1:ngrid)*(1.-albedo(1:ngrid,1))
925               fluxabs_sw(1:ngrid)=fluxsurfabs_sw(1:ngrid)
926               fluxrad_sky(1:ngrid)    = fluxsurfabs_sw(1:ngrid)
927               fluxtop_lw(1:ngrid)  = emis(1:ngrid)*sigma*tsurf(1:ngrid)**4
928
929               dtrad(1:ngrid,1:nlayer)=0.0 ! no atmospheric radiative heating
930
931            endif ! end of corrk
932
933         endif ! of if(mod(icount-1,iradia).eq.0)
934
935
936         ! Transformation of the radiative tendencies
937         ! ------------------------------------------
938         zplanck(1:ngrid)=tsurf(1:ngrid)*tsurf(1:ngrid)
939         zplanck(1:ngrid)=emis(1:ngrid)*sigma*zplanck(1:ngrid)*zplanck(1:ngrid)
940         fluxrad(1:ngrid)=fluxrad_sky(1:ngrid)-zplanck(1:ngrid)
941         pdt(1:ngrid,1:nlayer)=pdt(1:ngrid,1:nlayer)+dtrad(1:ngrid,1:nlayer)
942
943         ! Test of energy conservation
944         !----------------------------
945         if(enertest)then
946            call planetwide_sumval(cpp*massarea(:,:)*zdtsw(:,:)/totarea_planet,dEtotSW)
947            call planetwide_sumval(cpp*massarea(:,:)*zdtlw(:,:)/totarea_planet,dEtotLW)
948            !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
949            call planetwide_sumval(fluxsurfabs_sw(:)*cell_area(:)/totarea_planet,dEtotsSW) !JL13 carefull, albedo can have changed since the last time we called corrk
950            call planetwide_sumval((fluxsurf_lw(:)*emis(:)-zplanck(:))*cell_area(:)/totarea_planet,dEtotsLW)
951            dEzRadsw(:,:)=cpp*mass(:,:)*zdtsw(:,:)
952            dEzRadlw(:,:)=cpp*mass(:,:)*zdtlw(:,:)
953            if (is_master) then
954               print*,'---------------------------------------------------------------'
955               print*,'In corrk SW atmospheric heating       =',dEtotSW,' W m-2'
956               print*,'In corrk LW atmospheric heating       =',dEtotLW,' W m-2'
957               print*,'atmospheric net rad heating (SW+LW)   =',dEtotLW+dEtotSW,' W m-2'
958               print*,'In corrk SW surface heating           =',dEtotsSW,' W m-2'
959               print*,'In corrk LW surface heating           =',dEtotsLW,' W m-2'
960               print*,'surface net rad heating (SW+LW)       =',dEtotsLW+dEtotsSW,' W m-2'
961            endif
962         endif ! end of 'enertest'
963
964      endif ! of if (callrad)
965
966!!      call writediagfi(ngrid,"vdifc_pre_dqsurf"," "," ",2,dqsurf(1:ngrid,igcm_h2o_vap))
967!!      call writediagfi(ngrid,"vdifc_pre_qsurf"," "," ",2,qsurf(1:ngrid,igcm_h2o_vap))
968
969
970!  --------------------------------------------
971!  III. Vertical diffusion (turbulent mixing) :
972!  --------------------------------------------
973
974      if (calldifv) then
975
976         zflubid(1:ngrid)=fluxrad(1:ngrid)+fluxgrd(1:ngrid)
977
978         ! JL12 the following if test is temporarily there to allow us to compare the old vdifc with turbdiff.
979         if (UseTurbDiff) then
980
981            call turbdiff(ngrid,nlayer,nq,                  &
982                          ptimestep,capcal,                      &
983                          pplay,pplev,zzlay,zzlev,z0,            &
984                          pu,pv,pt,zpopsk,pq,tsurf,emis,qsurf,   &
985                          pdt,pdq,zflubid,                       &
986                          zdudif,zdvdif,zdtdif,zdtsdif,          &
987                          sensibFlux,q2,zdqdif,zdqevap,zdqsdif,  &
988                          taux,tauy)
989
990         else
991
992            zdh(1:ngrid,1:nlayer)=pdt(1:ngrid,1:nlayer)/zpopsk(1:ngrid,1:nlayer)
993
994            call vdifc(ngrid,nlayer,nq,zpopsk,           &
995                       ptimestep,capcal,lwrite,               &
996                       pplay,pplev,zzlay,zzlev,z0,            &
997                       pu,pv,zh,pq,tsurf,emis,qsurf,          &
998                       zdh,pdq,zflubid,                       &
999                       zdudif,zdvdif,zdhdif,zdtsdif,          &
1000                       sensibFlux,q2,zdqdif,zdqsdif)
1001
1002            zdtdif(1:ngrid,1:nlayer)=zdhdif(1:ngrid,1:nlayer)*zpopsk(1:ngrid,1:nlayer) ! for diagnostic only
1003            zdqevap(1:ngrid,1:nlayer)=0.
1004
1005         end if !end of 'UseTurbDiff'
1006
1007         zdtsurf(1:ngrid)=zdtsurf(1:ngrid)+zdtsdif(1:ngrid)
1008
1009         !!! this is always done, except for turbulence-resolving simulations
1010         if (.not. turb_resolved) then
1011           pdv(1:ngrid,1:nlayer)=pdv(1:ngrid,1:nlayer)+zdvdif(1:ngrid,1:nlayer)
1012           pdu(1:ngrid,1:nlayer)=pdu(1:ngrid,1:nlayer)+zdudif(1:ngrid,1:nlayer)
1013           pdt(1:ngrid,1:nlayer)=pdt(1:ngrid,1:nlayer)+zdtdif(1:ngrid,1:nlayer)
1014         endif
1015
1016         ! if(ok_slab_ocean)then !AF24: removed
1017         !    flux_sens_lat(1:ngrid)=(zdtsdif(1:ngrid)*capcal(1:ngrid)-fluxrad(1:ngrid))
1018         ! endif
1019
1020!!         call writediagfi(ngrid,"vdifc_post_zdqsdif"," "," ",2,zdqsdif(1:ngrid,igcm_h2o_vap))
1021
1022         if (tracer) then
1023           pdq(1:ngrid,1:nlayer,1:nq)=pdq(1:ngrid,1:nlayer,1:nq)+ zdqdif(1:ngrid,1:nlayer,1:nq)
1024           dqsurf(1:ngrid,1:nq)=dqsurf(1:ngrid,1:nq) + zdqsdif(1:ngrid,1:nq)
1025         end if ! of if (tracer)
1026
1027!!         call writediagfi(ngrid,"vdifc_post_dqsurf"," "," ",2,dqsurf(1:ngrid,igcm_h2o_vap))
1028!!         call writediagfi(ngrid,"vdifc_post_qsurf"," "," ",2,qsurf(1:ngrid,igcm_h2o_vap))
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. !AF24: removed
1063
1064      else ! calldifv
1065
1066         ! if(.not.newtonian)then
1067            zdtsurf(1:ngrid) = zdtsurf(1:ngrid) + (fluxrad(1:ngrid) + fluxgrd(1:ngrid))/capcal(1:ngrid)
1068
1069      endif ! end of 'calldifv'
1070
1071
1072!-------------------
1073!   IV. Convection :
1074!-------------------
1075
1076! ~~~~~~~~~~~~~~~~~~~~~~~~~~
1077! IV.a Thermal plume model : AF24: removed
1078! ~~~~~~~~~~~~~~~~~~~~~~~~~~
1079
1080! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1081! IV.b Dry convective adjustment :
1082! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1083
1084      if(calladj) then
1085
1086         zdh(1:ngrid,1:nlayer) = pdt(1:ngrid,1:nlayer)/zpopsk(1:ngrid,1:nlayer)
1087         zduadj(1:ngrid,1:nlayer)=0.0
1088         zdvadj(1:ngrid,1:nlayer)=0.0
1089         zdhadj(1:ngrid,1:nlayer)=0.0
1090
1091
1092         call convadj(ngrid,nlayer,nq,ptimestep,            &
1093                      pplay,pplev,zpopsk,                   &
1094                      pu,pv,zh,pq,                          &
1095                      pdu,pdv,zdh,pdq,                      &
1096                      zduadj,zdvadj,zdhadj,                 &
1097                      zdqadj)
1098
1099         pdu(1:ngrid,1:nlayer) = pdu(1:ngrid,1:nlayer) + zduadj(1:ngrid,1:nlayer)
1100         pdv(1:ngrid,1:nlayer) = pdv(1:ngrid,1:nlayer) + zdvadj(1:ngrid,1:nlayer)
1101         pdt(1:ngrid,1:nlayer)    = pdt(1:ngrid,1:nlayer) + zdhadj(1:ngrid,1:nlayer)*zpopsk(1:ngrid,1:nlayer)
1102         zdtadj(1:ngrid,1:nlayer) = zdhadj(1:ngrid,1:nlayer)*zpopsk(1:ngrid,1:nlayer) ! for diagnostic only
1103
1104         if(tracer) then
1105            pdq(1:ngrid,1:nlayer,1:nq) = pdq(1:ngrid,1:nlayer,1:nq) + zdqadj(1:ngrid,1:nlayer,1:nq)
1106         end if
1107
1108         ! Test energy conservation
1109         if(enertest)then
1110            call planetwide_sumval(cpp*massarea(:,:)*zdtadj(:,:)/totarea_planet,dEtot)
1111            if (is_master) print*,'In convadj atmospheric energy change  =',dEtot,' W m-2'
1112         endif
1113
1114         ! ! Test water conservation !AF24: removed
1115
1116      endif ! end of 'calladj'
1117!----------------------------------------------
1118!   Non orographic Gravity Waves: AF24: removed
1119!---------------------------------------------
1120
1121!-----------------------------------------------
1122!   V. Nitrogen condensation-sublimation :
1123!-----------------------------------------------
1124
1125      if (n2cond) then
1126
1127         if (.not.tracer) then
1128            print*,'We need a N2 ice tracer to condense N2'
1129            call abort
1130         endif
1131
1132         call condense_n2(ngrid,nlayer,nq,ptimestep,                    &
1133                           capcal,pplay,pplev,tsurf,pt,                 &
1134                           pphi,pdt,pdu,pdv,zdtsurf,pu,pv,pq,pdq,       &
1135                           qsurf(1,igcm_n2),albedo,emis,                &
1136                           zdtc,zdtsurfc,pdpsrf,zduc,zdvc,              &
1137                           zdqc,zdqsc(1,igcm_n2))
1138
1139         pdt(1:ngrid,1:nlayer) = pdt(1:ngrid,1:nlayer)+zdtc(1:ngrid,1:nlayer)
1140         pdv(1:ngrid,1:nlayer) = pdv(1:ngrid,1:nlayer)+zdvc(1:ngrid,1:nlayer)
1141         pdu(1:ngrid,1:nlayer) = pdu(1:ngrid,1:nlayer)+zduc(1:ngrid,1:nlayer)
1142         zdtsurf(1:ngrid)      = zdtsurf(1:ngrid) + zdtsurfc(1:ngrid)
1143
1144         pdq(1:ngrid,1:nlayer,1:nq)   = pdq(1:ngrid,1:nlayer,1:nq)+ zdqc(1:ngrid,1:nlayer,1:nq)
1145         dqsurf(1:ngrid,igcm_n2) = dqsurf(1:ngrid,igcm_n2) + zdqsc(1:ngrid,igcm_n2)
1146
1147!!         call writediagfi(ngrid,"condense_n2_post_dqsurf"," "," ",2,dqsurf(1:ngrid,igcm_h2o_vap))
1148!!         call writediagfi(ngrid,"condense_n2_post_qsurf"," "," ",2,qsurf(1:ngrid,igcm_h2o_vap))
1149
1150         ! test energy conservation
1151         if(enertest)then
1152            call planetwide_sumval(cpp*massarea(:,:)*zdtc(:,:)/totarea_planet,dEtot)
1153            call planetwide_sumval(capcal(:)*zdtsurfc(:)*cell_area(:)/totarea_planet,dEtots)
1154            if (is_master) then
1155               print*,'In n2cloud atmospheric energy change   =',dEtot,' W m-2'
1156               print*,'In n2cloud surface energy change       =',dEtots,' W m-2'
1157            endif
1158         endif
1159
1160      endif  ! end of 'n2cond'
1161
1162
1163!---------------------------------------------
1164!   VI. Specific parameterizations for tracers
1165!---------------------------------------------
1166
1167      if (tracer) then
1168
1169  ! ---------------------
1170  !   VI.1. Water and ice !AF24: removed
1171  ! ---------------------
1172  ! -------------------------
1173  !   VI.2. Photochemistry !AF24: removed
1174  ! -------------------------
1175  ! -------------------------
1176  !   VI.3. Aerosol particles
1177  ! -------------------------
1178
1179         !Generic Condensation
1180         if (generic_condensation) then
1181            call condensation_generic(ngrid,nlayer,nq,ptimestep,pplev,pplay,   &
1182                                          pt,pq,pdt,pdq,dt_generic_condensation, &
1183                                          dqvaplscale_generic,dqcldlscale_generic,rneb_generic)
1184            pdt(1:ngrid,1:nlayer) = pdt(1:ngrid,1:nlayer)+dt_generic_condensation(1:ngrid,1:nlayer)
1185            pdq(1:ngrid,1:nlayer,1:nq) = pdq(1:ngrid,1:nlayer,1:nq)+dqvaplscale_generic(1:ngrid,1:nlayer,1:nq)
1186            pdq(1:ngrid,1:nlayer,1:nq) = pdq(1:ngrid,1:nlayer,1:nq)+dqcldlscale_generic(1:ngrid,1:nlayer,1:nq)
1187
1188            if(enertest)then
1189               do ig=1,ngrid
1190                  genericconddE(ig) = cpp*SUM(mass(:,:)*dt_generic_condensation(:,:))
1191               enddo
1192
1193               call planetwide_sumval(cpp*massarea(:,:)*dt_generic_condensation(:,:)/totarea_planet,dEtot)
1194
1195               if (is_master) print*,'In generic condensation atmospheric energy change =',dEtot,' W m-2'
1196            end if
1197
1198            ! if (.not. water) then
1199               ! Compute GCS (Generic Condensable Specie) cloud fraction. For now we can not have both water cloud fraction and GCS cloud fraction
1200               ! Water is the priority
1201               ! If you have set water and generic_condensation, then cloudfrac will be water cloudfrac
1202               !
1203               ! If you have set generic_condensation (and not water) and you have set several GCS
1204               ! then cloudfrac will be only the cloudfrac of the last generic tracer
1205               ! (Because it is rewritten every tracer in the loop)
1206               !
1207               ! Maybe one should create a cloudfrac_generic(ngrid,nlayer,nq) with 3 dimensions, the last one for tracers
1208
1209               ! Let's loop on tracers
1210               cloudfrac(:,:)=0.0
1211               do iq=1,nq
1212                  call generic_tracer_index(nq,iq,igcm_generic_vap,igcm_generic_ice,call_ice_vap_generic)
1213                  if (call_ice_vap_generic) then ! to call only one time the ice/vap pair of a tracer
1214                     do l = 1, nlayer
1215                        do ig=1,ngrid
1216                           cloudfrac(ig,l)=rneb_generic(ig,l,iq)
1217                        enddo
1218                     enddo
1219                  endif
1220            end do ! do iq=1,nq loop on tracers
1221            ! endif ! .not. water
1222
1223         endif !generic_condensation
1224
1225         !Generic Rain !AF24: removed
1226
1227         ! Sedimentation.
1228         if (sedimentation) then
1229
1230            zdqsed(1:ngrid,1:nlayer,1:nq) = 0.0
1231            zdqssed(1:ngrid,1:nq)  = 0.0
1232
1233            ! if(watertest)then !AF24: removed
1234
1235            call callsedim(ngrid,nlayer,ptimestep,           &
1236                          pplev,zzlev,pt,pdt,pq,pdq,        &
1237                          zdqsed,zdqssed,nq)
1238
1239            ! if(watertest)then !AF24: removed
1240
1241            ! Whether it falls as rain or snow depends only on the surface temperature
1242            pdq(1:ngrid,1:nlayer,1:nq) = pdq(1:ngrid,1:nlayer,1:nq) + zdqsed(1:ngrid,1:nlayer,1:nq)
1243            dqsurf(1:ngrid,1:nq) = dqsurf(1:ngrid,1:nq) + zdqssed(1:ngrid,1:nq)
1244
1245!!            call writediagfi(ngrid,"callsedim_post_dqsurf"," "," ",2,dqsurf(1:ngrid,igcm_h2o_vap))
1246
1247            ! ! Test water conservation !AF24: removed
1248
1249         end if ! end of 'sedimentation'
1250
1251!!         call writediagfi(ngrid,"mass_redist_pre_dqsurf"," "," ",2,dqsurf(1:ngrid,igcm_h2o_vap))
1252!!         call writediagfi(ngrid,"mass_redist_pre_qsurf"," "," ",2,qsurf(1:ngrid,igcm_h2o_vap))
1253
1254  ! ---------------
1255  !   VI.4. Updates
1256  ! ---------------
1257
1258         ! Updating Atmospheric Mass and Tracers budgets.
1259         if(mass_redistrib) then
1260
1261            zdmassmr(1:ngrid,1:nlayer) = mass(1:ngrid,1:nlayer) * 0
1262            !    (   zdqevap(1:ngrid,1:nlayer)                          &
1263            !    !   + zdqrain(1:ngrid,1:nlayer,igcm_h2o_vap)             &
1264            !    !   + dqmoist(1:ngrid,1:nlayer,igcm_h2o_vap)             &
1265            !      + dqvaplscale(1:ngrid,1:nlayer) )
1266
1267            do ig = 1, ngrid
1268               zdmassmr_col(ig)=SUM(zdmassmr(ig,1:nlayer))
1269            enddo
1270
1271            ! call writediagfi(ngrid,"mass_evap","mass gain"," ",3,zdmassmr)
1272            ! call writediagfi(ngrid,"mass_evap_col","mass gain col"," ",2,zdmassmr_col)
1273            call writediagfi(ngrid,"mass","mass","kg/m2",3,mass)
1274
1275            call mass_redistribution(ngrid,nlayer,nq,ptimestep,                     &
1276                                     capcal,pplay,pplev,pt,tsurf,pq,qsurf,     &
1277                                     pu,pv,pdt,zdtsurf,pdq,pdu,pdv,zdmassmr,        &
1278                                     zdtmr,zdtsurfmr,zdpsrfmr,zdumr,zdvmr,zdqmr,zdqsurfmr)
1279
1280            pdq(1:ngrid,1:nlayer,1:nq) = pdq(1:ngrid,1:nlayer,1:nq) + zdqmr(1:ngrid,1:nlayer,1:nq)
1281            dqsurf(1:ngrid,1:nq)       = dqsurf(1:ngrid,1:nq) + zdqsurfmr(1:ngrid,1:nq)
1282            pdt(1:ngrid,1:nlayer)      = pdt(1:ngrid,1:nlayer) + zdtmr(1:ngrid,1:nlayer)
1283            pdu(1:ngrid,1:nlayer)      = pdu(1:ngrid,1:nlayer) + zdumr(1:ngrid,1:nlayer)
1284            pdv(1:ngrid,1:nlayer)      = pdv(1:ngrid,1:nlayer) + zdvmr(1:ngrid,1:nlayer)
1285            pdpsrf(1:ngrid)            = pdpsrf(1:ngrid) + zdpsrfmr(1:ngrid)
1286            zdtsurf(1:ngrid)           = zdtsurf(1:ngrid) + zdtsurfmr(1:ngrid)
1287
1288         endif
1289
1290!         call writediagfi(ngrid,"mass_redistribution_post_dqsurf"," "," ",2,dqsurf(1:ngrid,igcm_h2o_vap))
1291
1292!!         call writediagfi(ngrid,"slab_pre_dqsurf"," "," ",2,dqsurf(1:ngrid,igcm_h2o_vap))
1293!!         call writediagfi(ngrid,"slab_pre_qsurf"," "," ",2,qsurf(1:ngrid,igcm_h2o_vap))
1294
1295
1296! ------------------
1297!   VI.5. Slab Ocean !AF24: removed
1298! ------------------
1299
1300  ! -----------------------------
1301  !   VI.6. Surface Tracer Update
1302  ! -----------------------------
1303
1304         ! AF24: deleted hydrology
1305
1306         qsurf(1:ngrid,1:nq) = qsurf(1:ngrid,1:nq) + ptimestep*dqsurf(1:ngrid,1:nq)
1307
1308         ! Add qsurf to qsurf_hist, which is what we save in diagfi.nc. At the same time, we set the water
1309         ! content of ocean gridpoints back to zero, in order to avoid rounding errors in vdifc, rain.
1310         qsurf_hist(:,:) = qsurf(:,:)
1311
1312         ! if(ice_update)then
1313         !    ice_min(1:ngrid)=min(ice_min(1:ngrid),qsurf(1:ngrid,igcm_h2o_ice))
1314         ! endif
1315
1316      endif! end of if 'tracer'
1317
1318
1319!------------------------------------------------
1320!   VII. Surface and sub-surface soil temperature
1321!------------------------------------------------
1322
1323
1324      ! ! Increment surface temperature
1325      ! if(ok_slab_ocean)then  !AF24: removed
1326
1327      tsurf(1:ngrid)=tsurf(1:ngrid)+ptimestep*zdtsurf(1:ngrid)
1328      ! Compute soil temperatures and subsurface heat flux.
1329      if (callsoil) then
1330         call soil(ngrid,nsoilmx,.false.,lastcall,inertiedat,   &
1331                   ptimestep,tsurf,tsoil,capcal,fluxgrd)
1332      endif
1333
1334
1335!       if (ok_slab_ocean) then !AF24: removed
1336
1337      ! Test energy conservation
1338      if(enertest)then
1339         call planetwide_sumval(cell_area(:)*capcal(:)*zdtsurf(:)/totarea_planet,dEtots)
1340         if (is_master) print*,'Surface energy change                 =',dEtots,' W m-2'
1341      endif
1342
1343
1344!---------------------------------------------------
1345!   VIII. Perform diagnostics and write output files
1346!---------------------------------------------------
1347
1348      ! Note : For output only: the actual model integration is performed in the dynamics.
1349
1350
1351      ! Temperature, zonal and meridional winds.
1352      zt(1:ngrid,1:nlayer) = pt(1:ngrid,1:nlayer) + pdt(1:ngrid,1:nlayer)*ptimestep
1353      zu(1:ngrid,1:nlayer) = pu(1:ngrid,1:nlayer) + pdu(1:ngrid,1:nlayer)*ptimestep
1354      zv(1:ngrid,1:nlayer) = pv(1:ngrid,1:nlayer) + pdv(1:ngrid,1:nlayer)*ptimestep
1355
1356      !! Recast thermal plume vertical velocity array for outputs
1357      !! AF24: removed
1358
1359      ! Diagnostic.
1360      zdtdyn(1:ngrid,1:nlayer)     = (pt(1:ngrid,1:nlayer)-ztprevious(1:ngrid,1:nlayer)) / ptimestep
1361      ztprevious(1:ngrid,1:nlayer) = zt(1:ngrid,1:nlayer)
1362
1363      zdudyn(1:ngrid,1:nlayer)     = (pu(1:ngrid,1:nlayer)-zuprevious(1:ngrid,1:nlayer)) / ptimestep
1364      zuprevious(1:ngrid,1:nlayer) = zu(1:ngrid,1:nlayer)
1365
1366      if(firstcall)then
1367         zdtdyn(1:ngrid,1:nlayer)=0.0
1368         zdudyn(1:ngrid,1:nlayer)=0.0
1369      endif
1370
1371      ! Dynamical heating diagnostic.
1372      do ig=1,ngrid
1373         fluxdyn(ig)= SUM(zdtdyn(ig,:) *mass(ig,:))*cpp
1374      enddo
1375
1376      ! Tracers.
1377      zq(1:ngrid,1:nlayer,1:nq) = pq(1:ngrid,1:nlayer,1:nq) + pdq(1:ngrid,1:nlayer,1:nq)*ptimestep
1378
1379      ! Surface pressure.
1380      ps(1:ngrid) = pplev(1:ngrid,1) + pdpsrf(1:ngrid)*ptimestep
1381
1382
1383      ! Surface and soil temperature information
1384      call planetwide_sumval(cell_area(:)*tsurf(:)/totarea_planet,Ts1)
1385      call planetwide_minval(tsurf(:),Ts2)
1386      call planetwide_maxval(tsurf(:),Ts3)
1387      if(callsoil)then
1388         TsS = SUM(cell_area(:)*tsoil(:,nsoilmx))/totarea        ! mean temperature at bottom soil layer
1389         if (is_master) then
1390            print*,'          ave[Tsurf]             min[Tsurf]             max[Tsurf]             ave[Tdeep]'
1391            print*,Ts1,Ts2,Ts3,TsS
1392         end if
1393      else
1394         if (is_master) then
1395            print*,'          ave[Tsurf]             min[Tsurf]             max[Tsurf]'
1396            print*,Ts1,Ts2,Ts3
1397         endif
1398      end if
1399
1400
1401      ! Check the energy balance of the simulation during the run
1402      if(corrk)then
1403
1404         call planetwide_sumval(cell_area(:)*fluxtop_dn(:)/totarea_planet,ISR)
1405         call planetwide_sumval(cell_area(:)*fluxabs_sw(:)/totarea_planet,ASR)
1406         call planetwide_sumval(cell_area(:)*fluxtop_lw(:)/totarea_planet,OLR)
1407         call planetwide_sumval(cell_area(:)*fluxgrd(:)/totarea_planet,GND)
1408         call planetwide_sumval(cell_area(:)*fluxdyn(:)/totarea_planet,DYN)
1409         do ig=1,ngrid
1410            if(fluxtop_dn(ig).lt.0.0)then
1411               print*,'fluxtop_dn has gone crazy'
1412               print*,'fluxtop_dn=',fluxtop_dn(ig)
1413               print*,'tau_col=',tau_col(ig)
1414               print*,'aerosol=',aerosol(ig,:,:)
1415               print*,'temp=   ',pt(ig,:)
1416               print*,'pplay=  ',pplay(ig,:)
1417               call abort
1418            endif
1419         end do
1420
1421         if(ngrid.eq.1)then
1422            DYN=0.0
1423         endif
1424
1425         if (is_master) then
1426            print*,'  ISR            ASR            OLR            GND            DYN [W m^-2]'
1427            print*, ISR,ASR,OLR,GND,DYN
1428         endif
1429
1430         if(enertest .and. is_master)then
1431            print*,'SW flux/heating difference SW++ - ASR = ',dEtotSW+dEtotsSW-ASR,' W m-2'
1432            print*,'LW flux/heating difference LW++ - OLR = ',dEtotLW+dEtotsLW+OLR,' W m-2'
1433            print*,'LW energy balance LW++ + ASR = ',dEtotLW+dEtotsLW+ASR,' W m-2'
1434         endif
1435
1436         if(meanOLR .and. is_master)then
1437            if((ngrid.gt.1) .or. (mod(icount-1,ecritphy).eq.0))then
1438               ! to record global radiative balance
1439               open(92,file="rad_bal.out",form='formatted',position='append')
1440               write(92,*) zday,ISR,ASR,OLR
1441               close(92)
1442               open(93,file="tem_bal.out",form='formatted',position='append')
1443               if(callsoil)then
1444                         write(93,*) zday,Ts1,Ts2,Ts3,TsS
1445               else
1446                  write(93,*) zday,Ts1,Ts2,Ts3
1447               endif
1448               close(93)
1449            endif
1450         endif
1451
1452      endif ! end of 'corrk'
1453
1454
1455      ! Diagnostic to test radiative-convective timescales in code.
1456      if(testradtimes)then
1457         open(38,file="tau_phys.out",form='formatted',position='append')
1458         ig=1
1459         do l=1,nlayer
1460            write(38,*) -1./pdt(ig,l),pt(ig,l),pplay(ig,l)
1461         enddo
1462         close(38)
1463         print*,'As testradtimes enabled,'
1464         print*,'exiting physics on first call'
1465         call abort
1466      endif
1467
1468
1469      ! Compute column amounts (kg m-2) if tracers are enabled.
1470      if(tracer)then
1471         qcol(1:ngrid,1:nq)=0.0
1472         do iq=1,nq
1473            do ig=1,ngrid
1474               qcol(ig,iq) = SUM( zq(ig,1:nlayer,iq) * mass(ig,1:nlayer))
1475            enddo
1476         enddo
1477
1478         if (aerohaze) then
1479            ! Generalised for arbitrary aerosols now. By LK
1480         reffcol(1:ngrid,1:naerkind)=0.0
1481            ! call n2_reffrad(ngrid,nlayer,nq,zq,reffrad(1,1,iaero_haze))
1482         if (haze_proffix.and.i_haze.gt.0.) then
1483               call haze_prof(ngrid,nlayer,zzlay,pplay,pt, &
1484                                        reffrad,profmmr)
1485               zdqhaze(:,:,i_haze)=(profmmr(:,:)-pq(:,:,igcm_n2)) & ! AF: TODO: replace by igcm_haze?
1486                                                        /ptimestep
1487              else
1488               !! AF: TODO import from pluto.old?
1489               ! call hazecloud(ngrid,nlayer,nq,ptimestep,&
1490               !    pplay,pplev,pq,pdq,dist_sol,mu0,zfluxuv,zdqhaze,&
1491               !    zdqphot_prec,zdqphot_ch4,zdqconv_prec,declin)
1492              endif
1493
1494         DO iq=1, nq ! should be updated
1495            DO l=1,nlayer
1496              DO ig=1,ngrid
1497                   pdq(ig,l,iq)=pdq(ig,l,iq)+ zdqhaze(ig,l,iq)
1498              ENDDO
1499            ENDDO
1500         ENDDO
1501         endif ! end of aerohaze
1502      endif ! end of 'tracer'
1503
1504
1505      ! ! Test for water conservation. !AF24: removed
1506
1507      ! Calculate RH_generic (Generic Relative Humidity) for diagnostic.
1508      if(generic_condensation)then
1509         RH_generic(:,:,:)=0.0
1510         do iq=1,nq
1511
1512            call generic_tracer_index(nq,iq,igcm_generic_vap,igcm_generic_ice,call_ice_vap_generic)
1513
1514            if (call_ice_vap_generic) then ! to call only one time the ice/vap pair of a tracer
1515
1516               do l = 1, nlayer
1517                  do ig=1,ngrid
1518                     call Psat_generic(zt(ig,l),pplay(ig,l),metallicity,psat_tmp_generic,qsat_generic(ig,l,iq))
1519                     RH_generic(ig,l,iq) = zq(ig,l,igcm_generic_vap) / qsat_generic(ig,l,iq)
1520                  enddo
1521               enddo
1522
1523            end if
1524
1525         end do ! iq=1,nq
1526
1527      endif ! end of 'generic_condensation'
1528
1529
1530      if (is_master) print*,'--> Ls =',zls*180./pi
1531
1532
1533!----------------------------------------------------------------------
1534!        Writing NetCDF file  "RESTARTFI" at the end of the run
1535!----------------------------------------------------------------------
1536
1537!        Note: 'restartfi' is stored just before dynamics are stored
1538!              in 'restart'. Between now and the writting of 'restart',
1539!              there will have been the itau=itau+1 instruction and
1540!              a reset of 'time' (lastacll = .true. when itau+1= itaufin)
1541!              thus we store for time=time+dtvr
1542
1543
1544
1545      if(lastcall) then
1546         ztime_fin = ptime + ptimestep/(float(iphysiq)*daysec)
1547
1548         !! Update surface ice distribution to iterate to steady state if requested
1549         !! AF24: removed
1550
1551         ! endif
1552#ifndef MESOSCALE
1553
1554         if (ngrid.ne.1) then
1555            write(*,*)'PHYSIQ: for physdem ztime_fin =',ztime_fin
1556
1557            call physdem1("restartfi.nc",nsoilmx,ngrid,nlayer,nq, &
1558                          ptimestep,ztime_fin,                    &
1559                          tsurf,tsoil,emis,q2,qsurf_hist)
1560         endif
1561#endif
1562         ! if(ok_slab_ocean) then
1563         !    call ocean_slab_final!(tslab, seaice)
1564         ! end if
1565
1566    endif ! end of 'lastcall'
1567
1568
1569! -----------------------------------------------------------------
1570! WSTATS: Saving statistics
1571! -----------------------------------------------------------------
1572! ("stats" stores and accumulates key variables in file "stats.nc"
1573!  which can later be used to make the statistic files of the run:
1574!  if flag "callstats" from callphys.def is .true.)
1575
1576
1577         call wstats(ngrid,"ps","Surface pressure","Pa",2,ps)
1578         call wstats(ngrid,"tsurf","Surface temperature","K",2,tsurf)
1579         call wstats(ngrid,"fluxsurf_lw",                               &
1580                     "Thermal IR radiative flux to surface","W.m-2",2,  &
1581                     fluxsurf_lw)
1582         call wstats(ngrid,"fluxtop_lw",                                &
1583                     "Thermal IR radiative flux to space","W.m-2",2,    &
1584                     fluxtop_lw)
1585
1586!            call wstats(ngrid,"fluxsurf_sw",                               &
1587!                        "Solar radiative flux to surface","W.m-2",2,       &
1588!                         fluxsurf_sw_tot)
1589!            call wstats(ngrid,"fluxtop_sw",                                &
1590!                        "Solar radiative flux to space","W.m-2",2,         &
1591!                        fluxtop_sw_tot)
1592
1593
1594         call wstats(ngrid,"ISR","incoming stellar rad.","W m-2",2,fluxtop_dn)
1595         call wstats(ngrid,"ASR","absorbed stellar rad.","W m-2",2,fluxabs_sw)
1596         call wstats(ngrid,"OLR","outgoing longwave rad.","W m-2",2,fluxtop_lw)
1597         !call wstats(ngrid,"ALB","Surface albedo"," ",2,albedo_equivalent)
1598         !call wstats(ngrid,"ALB_1st","First Band Surface albedo"," ",2,albedo(:,1))
1599         call wstats(ngrid,"p","Pressure","Pa",3,pplay)
1600         call wstats(ngrid,"temp","Atmospheric temperature","K",3,zt)
1601         call wstats(ngrid,"u","Zonal (East-West) wind","m.s-1",3,zu)
1602         call wstats(ngrid,"v","Meridional (North-South) wind","m.s-1",3,zv)
1603         call wstats(ngrid,"w","Vertical (down-up) wind","m.s-1",3,pw)
1604         call wstats(ngrid,"q2","Boundary layer eddy kinetic energy","m2.s-2",3,q2)
1605
1606         if (tracer) then
1607            do iq=1,nq
1608               call wstats(ngrid,noms(iq),noms(iq),'kg/kg',3,zq(1,1,iq))
1609               call wstats(ngrid,trim(noms(iq))//'_surf',trim(noms(iq))//'_surf',  &
1610                           'kg m^-2',2,qsurf(1,iq) )
1611               call wstats(ngrid,trim(noms(iq))//'_col',trim(noms(iq))//'_col',    &
1612                          'kg m^-2',2,qcol(1,iq) )
1613
1614!              call wstats(ngrid,trim(noms(iq))//'_reff',                          &
1615!                          trim(noms(iq))//'_reff',                                   &
1616!                          'm',3,reffrad(1,1,iq))
1617
1618            end do
1619
1620         endif ! end of 'tracer'
1621
1622         !AF24: deleted slab ocean and water
1623
1624         if(lastcall.and.callstats) then
1625            write (*,*) "Writing stats..."
1626            call mkstats(ierr)
1627         endif
1628
1629
1630#ifndef MESOSCALE
1631
1632!-----------------------------------------------------------------------------------------------------
1633!           OUTPUT in netcdf file "DIAGFI.NC", containing any variable for diagnostic
1634!
1635!             Note 1 : output with  period "ecritphy", set in "run.def"
1636!
1637!             Note 2 : writediagfi can also be called from any other subroutine for any variable,
1638!                      but its preferable to keep all the calls in one place ...
1639!-----------------------------------------------------------------------------------------------------
1640
1641      call writediagfi(ngrid,"Ls","solar longitude","deg",0,zls*180./pi)
1642      call writediagfi(ngrid,"Lss","sub solar longitude","deg",0,zlss*180./pi)
1643      call writediagfi(ngrid,"RA","right ascension","deg",0,right_ascen*180./pi)
1644      call writediagfi(ngrid,"Declin","solar declination","deg",0,declin*180./pi)
1645      call writediagfi(ngrid,"tsurf","Surface temperature","K",2,tsurf)
1646      call writediagfi(ngrid,"ps","Surface pressure","Pa",2,ps)
1647
1648      if (.not.fast) then
1649       call writediagfi(ngrid,"temp","temperature","K",3,zt)
1650       call writediagfi(ngrid,"teta","potential temperature","K",3,zh)
1651       call writediagfi(ngrid,"u","Zonal wind","m.s-1",3,zu)
1652       call writediagfi(ngrid,"v","Meridional wind","m.s-1",3,zv)
1653       call writediagfi(ngrid,"w","Vertical wind","m.s-1",3,pw)
1654       call writediagfi(ngrid,"p","Pressure","Pa",3,pplay)
1655      endif
1656
1657!     Subsurface temperatures
1658!        call writediagsoil(ngrid,"tsurf","Surface temperature","K",2,tsurf)
1659!        call writediagsoil(ngrid,"temp","temperature","K",3,tsoil)
1660
1661!       ! Oceanic layers !AF24: removed
1662
1663      ! ! Thermal plume model  !AF24: removed
1664
1665!        GW non-oro outputs  !AF24: removed
1666
1667      ! Total energy balance diagnostics
1668      if(callrad)then
1669
1670         call writediagfi(ngrid,"ALB","Surface albedo"," ",2,albedo_equivalent)
1671         !call writediagfi(ngrid,"ALB_1st","First Band Surface albedo"," ",2,albedo(:,1))
1672         call writediagfi(ngrid,"ISR","incoming stellar rad.","W m-2",2,fluxtop_dn)
1673         call writediagfi(ngrid,"ASR","absorbed stellar rad.","W m-2",2,fluxabs_sw)
1674         call writediagfi(ngrid,"OLR","outgoing longwave rad.","W m-2",2,fluxtop_lw)
1675         call writediagfi(ngrid,"shad","rings"," ", 2, fract)
1676
1677!           call writediagfi(ngrid,"ASRcs","absorbed stellar rad (cs).","W m-2",2,fluxabs_sw1)
1678!           call writediagfi(ngrid,"OLRcs","outgoing longwave rad (cs).","W m-2",2,fluxtop_lw1)
1679!           call writediagfi(ngrid,"fluxsurfsw","sw surface flux.","W m-2",2,fluxsurf_sw)
1680!           call writediagfi(ngrid,"fluxsurflw","lw back radiation.","W m-2",2,fluxsurf_lw)
1681!           call writediagfi(ngrid,"fluxsurfswcs","sw surface flux (cs).","W m-2",2,fluxsurf_sw1)
1682!           call writediagfi(ngrid,"fluxsurflwcs","lw back radiation (cs).","W m-2",2,fluxsurf_lw1)
1683
1684         ! if(ok_slab_ocean) then
1685         !    call writediagfi(ngrid,"GND","heat flux from ground","W m-2",2,fluxgrdocean)
1686         ! else
1687         call writediagfi(ngrid,"GND","heat flux from ground","W m-2",2,fluxgrd)
1688         ! endif
1689
1690         call writediagfi(ngrid,"DYN","dynamical heat input","W m-2",2,fluxdyn)
1691
1692      endif ! end of 'callrad'
1693
1694      if(enertest) then
1695
1696         if (calldifv) then
1697
1698            call writediagfi(ngrid,"q2","turbulent kinetic energy","J.kg^-1",3,q2)
1699            call writediagfi(ngrid,"sensibFlux","sensible heat flux","w.m^-2",2,sensibFlux)
1700
1701!             call writediagfi(ngrid,"dEzdiff","turbulent diffusion heating (-sensible flux)","w.m^-2",3,dEzdiff)
1702!             call writediagfi(ngrid,"dEdiff","integrated turbulent diffusion heating (-sensible flux)","w.m^-2",2,dEdiff)
1703!             call writediagfi(ngrid,"dEdiffs","In TurbDiff (correc rad+latent heat) surf nrj change","w.m^-2",2,dEdiffs)
1704
1705         endif
1706
1707         if (corrk) then
1708            call writediagfi(ngrid,"dEzradsw","radiative heating","w.m^-2",3,dEzradsw)
1709            call writediagfi(ngrid,"dEzradlw","radiative heating","w.m^-2",3,dEzradlw)
1710         endif
1711
1712!          if(watercond) then  !AF24: removed
1713
1714         if (generic_condensation) then
1715
1716            call writediagfi(ngrid,"genericconddE","heat from generic condensation","W m-2",2,genericconddE)
1717            call writediagfi(ngrid,"dt_generic_condensation","heating from generic condensation","K s-1",3,dt_generic_condensation)
1718
1719         endif
1720
1721      endif ! end of 'enertest'
1722
1723        ! Diagnostics of optical thickness
1724        ! Warning this is exp(-tau), I let you postproc with -log to have tau itself - JVO 19
1725        if (diagdtau) then
1726          do nw=1,L_NSPECTV
1727            write(str2,'(i2.2)') nw
1728            call writediagfi(ngrid,'dtauv'//str2,'Layer optical thickness attenuation in VI band '//str2,'',1,int_dtauv(:,nlayer:1:-1,nw))
1729          enddo
1730          do nw=1,L_NSPECTI
1731            write(str2,'(i2.2)') nw
1732            call writediagfi(ngrid,'dtaui'//str2,'Layer optical thickness attenuation in IR band '//str2,'',1,int_dtaui(:,nlayer:1:-1,nw))
1733          enddo
1734        endif
1735
1736
1737        ! Temporary inclusions for heating diagnostics.
1738        call writediagfi(ngrid,"zdtsw","SW heating","T s-1",3,zdtsw)
1739        call writediagfi(ngrid,"zdtlw","LW heating","T s-1",3,zdtlw)
1740        call writediagfi(ngrid,"dtrad","radiative heating","K s-1",3,dtrad)
1741        call writediagfi(ngrid,"zdtdyn","Dyn. heating","T s-1",3,zdtdyn)
1742
1743        ! For Debugging.
1744        !call writediagfi(ngrid,'rnat','Terrain type',' ',2,real(rnat))
1745        !call writediagfi(ngrid,'pphi','Geopotential',' ',3,pphi)
1746
1747
1748      ! Output aerosols.!AF: TODO: write haze aerosols
1749      ! if (igcm_n2_ice.ne.0.and.iaero_haze.ne.0) &
1750      !    call writediagfi(ngrid,'N2ice_reff','N2ice_reff','m',3,reffrad(1,1,iaero_haze))
1751      ! if (igcm_n2_ice.ne.0.and.iaero_haze.ne.0) &
1752      !    call writediagfi(ngrid,'N2ice_reffcol','N2ice_reffcol','um kg m^-2',2,reffcol(1,iaero_haze))
1753      ! if (igcm_h2o_ice.ne.0.and.iaero_h2o.ne.0) &  !AF24: removed
1754
1755      ! Output tracers.
1756      if (tracer) then
1757
1758         do iq=1,nq
1759            call writediagfi(ngrid,noms(iq),noms(iq),'kg/kg',3,zq(1,1,iq))
1760            call writediagfi(ngrid,trim(noms(iq))//'_surf',trim(noms(iq))//'_surf',  &
1761                             'kg m^-2',2,qsurf_hist(1,iq) )
1762            call writediagfi(ngrid,trim(noms(iq))//'_col',trim(noms(iq))//'_col',    &
1763                            'kg m^-2',2,qcol(1,iq) )
1764!          call writediagfi(ngrid,trim(noms(iq))//'_surf',trim(noms(iq))//'_surf',  &
1765!                          'kg m^-2',2,qsurf(1,iq) )
1766
1767            ! if(watercond.or.CLFvarying)then  !AF24: removed
1768
1769            if(generic_condensation)then
1770               call writediagfi(ngrid,"rneb_generic","GCS cloud fraction (generic condensation)"," ",3,rneb_generic)
1771               call writediagfi(ngrid,"CLF","GCS cloud fraction"," ",3,cloudfrac)
1772               call writediagfi(ngrid,"RH_generic","GCS relative humidity"," ",3,RH_generic)
1773            endif
1774
1775            ! if(generic_rain)then  !AF24: removed
1776            ! if((hydrology).and.(.not.ok_slab_ocean))then  !AF24: removed
1777
1778            call writediagfi(ngrid,"tau_col","Total aerosol optical depth","[]",2,tau_col)
1779
1780         enddo ! end of 'nq' loop
1781
1782       endif ! end of 'tracer'
1783
1784
1785      ! Output spectrum.
1786      if(specOLR.and.corrk)then
1787         call writediagspecIR(ngrid,"OLR3D","OLR(lon,lat,band)","W/m^2/cm^-1",3,OLR_nu)
1788         call writediagspecVI(ngrid,"OSR3D","OSR(lon,lat,band)","W/m^2/cm^-1",3,OSR_nu)
1789         call writediagspecVI(ngrid,"GSR3D","GSR(lon,lat,band)","W/m^2/cm^-1",3,GSR_nu)
1790      endif
1791
1792#else
1793   comm_HR_SW(1:ngrid,1:nlayer) = zdtsw(1:ngrid,1:nlayer)
1794   comm_HR_LW(1:ngrid,1:nlayer) = zdtlw(1:ngrid,1:nlayer)
1795   comm_ALBEQ(1:ngrid)=albedo_equivalent(1:ngrid)
1796   if (.not.calldifv) comm_LATENT_HF(:)=0.0
1797   ! if ((tracer).and.(water)) then  !AF24: removed
1798
1799   if ((tracer).and.(generic_condensation)) then
1800   ! .and.(.not. water)
1801
1802      ! If you have set generic_condensation (and not water) and you have set several GCS
1803      ! then the outputs given to WRF will be only the ones for the last generic tracer
1804      ! (Because it is rewritten every tracer in the loop)
1805      ! WRF can take only one moist tracer
1806
1807      do iq=1,nq
1808         call generic_tracer_index(nq,iq,igcm_generic_vap,igcm_generic_ice,call_ice_vap_generic)
1809
1810         if (call_ice_vap_generic) then ! to call only one time the ice/vap pair of a tracer
1811
1812            reffrad_generic_zeros_for_wrf(:,:) = 1.
1813
1814            comm_CLOUDFRAC(1:ngrid,1:nlayer) = cloudfrac(1:ngrid,1:nlayer)
1815            comm_TOTCLOUDFRAC(1:ngrid) = totcloudfrac(1:ngrid) !??????
1816            comm_SURFRAIN(1:ngrid) = zdqsrain_generic(1:ngrid,iq)
1817            comm_DQVAP(1:ngrid,1:nlayer) = pdq(1:ngrid,1:nlayer,igcm_generic_vap)
1818            comm_DQICE(1:ngrid,1:nlayer)=pdq(1:ngrid,1:nlayer,igcm_generic_ice)
1819            ! comm_H2OICE_REFF(1:ngrid,1:nlayer) = reffrad_generic_zeros_for_wrf(1:ngrid,1:nlayer) ! for now zeros !
1820            !comm_H2OICE_REFF(1:ngrid,1:nlayer) = 0*zdtrain_generic(1:ngrid,1:nlayer) ! for now zeros !
1821            comm_REEVAP(1:ngrid) = reevap_precip_generic(1:ngrid,iq)
1822            comm_DTRAIN(1:ngrid,1:nlayer) = zdtrain_generic(1:ngrid,1:nlayer)
1823            comm_DTLSC(1:ngrid,1:nlayer) = dt_generic_condensation(1:ngrid,1:nlayer)
1824            comm_RH(1:ngrid,1:nlayer) = RH_generic(1:ngrid,1:nlayer,iq)
1825
1826         endif
1827      end do ! do iq=1,nq loop on tracers
1828
1829   else
1830      comm_CLOUDFRAC(1:ngrid,1:nlayer)=0.
1831      comm_TOTCLOUDFRAC(1:ngrid)=0.
1832      comm_SURFRAIN(1:ngrid)=0.
1833      comm_DQVAP(1:ngrid,1:nlayer)=0.
1834      comm_DQICE(1:ngrid,1:nlayer)=0.
1835      ! comm_H2OICE_REFF(1:ngrid,1:nlayer)=0.
1836      comm_REEVAP(1:ngrid)=0.
1837      comm_DTRAIN(1:ngrid,1:nlayer)=0.
1838      comm_DTLSC(1:ngrid,1:nlayer)=0.
1839      comm_RH(1:ngrid,1:nlayer)=0.
1840
1841   endif ! if water, if generic_condensation, else
1842
1843   comm_FLUXTOP_DN(1:ngrid)=fluxtop_dn(1:ngrid)
1844   comm_FLUXABS_SW(1:ngrid)=fluxabs_sw(1:ngrid)
1845   comm_FLUXTOP_LW(1:ngrid)=fluxtop_lw(1:ngrid)
1846   comm_FLUXSURF_SW(1:ngrid)=fluxsurf_sw(1:ngrid)
1847   comm_FLUXSURF_LW(1:ngrid)=fluxsurf_lw(1:ngrid)
1848   comm_FLXGRD(1:ngrid)=fluxgrd(1:ngrid)
1849   sensibFlux(1:ngrid) = zflubid(1:ngrid) - capcal(1:ngrid)*zdtsdif(1:ngrid) !!! ????
1850   comm_HR_DYN(1:ngrid,1:nlayer) = zdtdyn(1:ngrid,1:nlayer)
1851
1852!      if (turb_resolved) then
1853!        open(17,file='lsf.txt',form='formatted',status='old')
1854!        rewind(17)
1855!        DO l=1,nlayer
1856!          read(17,*) lsf_dt(l),lsf_dq(l)
1857!        ENDDO
1858!        close(17)
1859!        do ig=1,ngrid
1860!          if ((tracer).and.(water)) then
1861!           pdq(ig,:,igcm_h2o_vap) = pdq(ig,:,igcm_h2o_vap) + lsf_dq(:)
1862!          endif
1863!          pdt(ig,:) = pdt(ig,:) + lsf_dt(:)
1864!          comm_HR_DYN(ig,:) = lsf_dt(:)
1865!        enddo
1866!      endif
1867#endif
1868
1869! XIOS outputs
1870#ifdef CPP_XIOS
1871      ! Send fields to XIOS: (NB these fields must also be defined as
1872      ! <field id="..." /> in context_lmdz_physics.xml to be correctly used)
1873      CALL send_xios_field("ls",zls)
1874
1875      CALL send_xios_field("ps",ps)
1876      CALL send_xios_field("area",cell_area)
1877      CALL send_xios_field("p",pplay)
1878      CALL send_xios_field("temperature",zt)
1879      CALL send_xios_field("u",zu)
1880      CALL send_xios_field("v",zv)
1881      CALL send_xios_field("omega",omega)
1882
1883!       IF (calltherm) THEN  !AF24: removed
1884      ! IF (water) THEN  !AF24: removed
1885
1886      CALL send_xios_field("ISR",fluxtop_dn)
1887      CALL send_xios_field("OLR",fluxtop_lw)
1888      CALL send_xios_field("ASR",fluxabs_sw)
1889
1890      if (specOLR .and. corrk) then
1891         call send_xios_field("OLR3D",OLR_nu)
1892         call send_xios_field("IR_Bandwidth",DWNI)
1893         call send_xios_field("VI_Bandwidth",DWNV)
1894         call send_xios_field("OSR3D",OSR_nu)
1895         call send_xios_field("GSR3D",GSR_nu)
1896      endif
1897
1898      if (lastcall.and.is_omp_master) then
1899        write(*,*) "physiq: call xios_context_finalize"
1900        call xios_context_finalize
1901      endif
1902#endif
1903
1904      if (check_physics_outputs) then
1905         ! Check the validity of updated fields at the end of the physics step
1906         call check_physics_fields("end of physiq:", zt, zu, zv, pplev, zq)
1907      endif
1908
1909      icount=icount+1
1910
1911      end subroutine physiq
1912
1913   end module physiq_mod
Note: See TracBrowser for help on using the repository browser.