source: LMDZ6/trunk/libf/phylmdiso/physiq_mod.F90 @ 4613

Last change on this file since 4613 was 4613, checked in by fhourdin, 12 months ago

Integration/replay_isation travail Louis (ratqs)

File size: 253.7 KB
Line 
1!
2! $Id: physiq_mod.F90 3908 2021-05-20 07:11:13Z idelkadi $
3!
4!#define IO_DEBUG
5MODULE physiq_mod
6
7  IMPLICIT NONE
8
9CONTAINS
10
11  SUBROUTINE physiq (nlon,nlev, &
12       debut,lafin,pdtphys_, &
13       paprs,pplay,pphi,pphis,presnivs, &
14       u,v,rot,t,qx, &
15       flxmass_w, &
16       d_u, d_v, d_t, d_qx, d_ps)
17
18! For clarity, the "USE" section is now arranged in alphabetical order,
19! with a separate section for CPP keys
20! PLEASE try to follow this rule
21
22    USE ACAMA_GWD_rando_m, only: ACAMA_GWD_rando
23    USE aero_mod
24    USE add_phys_tend_mod, only : add_pbl_tend, add_phys_tend, diag_phys_tend, prt_enerbil, &
25  &      fl_ebil, fl_cor_ebil
26    USE assert_m, only: assert
27    USE change_srf_frac_mod
28    USE conf_phys_m, only: conf_phys
29    USE carbon_cycle_mod, ONLY : infocfields_init, RCO2_glo, carbon_cycle_rad
30    USE CFMIP_point_locations   ! IM stations CFMIP
31    USE cmp_seri_mod
32    USE dimphy
33    USE etat0_limit_unstruct_mod
34    USE FLOTT_GWD_rando_m, only: FLOTT_GWD_rando
35    USE fonte_neige_mod, ONLY  : fonte_neige_get_vars
36    USE geometry_mod, ONLY: cell_area, latitude_deg, longitude_deg
37    USE ioipsl, only: histbeg, histvert, histdef, histend, histsync, &
38         histwrite, ju2ymds, ymds2ju, getin
39    USE ioipsl_getin_p_mod, ONLY : getin_p
40    USE indice_sol_mod
41    USE infotrac_phy, ONLY: nqtot, nbtr, nqo, tracers, type_trac
42    USE readTracFiles_mod, ONLY: addPhase
43    USE strings_mod,  ONLY: strIdx
44    USE iophy
45    USE limit_read_mod, ONLY : init_limit_read
46    USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat, nbp_lev, klon_glo, grid1dTo2d_glo, grid_type, unstructured
47    USE mod_phys_lmdz_mpi_data, only: is_mpi_root
48    USE mod_phys_lmdz_para
49    USE netcdf95, only: nf95_close
50    USE netcdf, only: nf90_fill_real     ! IM for NMC files
51    USE open_climoz_m, only: open_climoz ! ozone climatology from a file
52    USE ozonecm_m, only: ozonecm ! ozone of J.-F. Royer
53    USE pbl_surface_mod, ONLY : pbl_surface
54    USE phyaqua_mod, only: zenang_an
55    USE phyetat0_mod, only: phyetat0
56    USE phystokenc_mod, ONLY: offline, phystokenc
57    USE phys_cal_mod, only: year_len, mth_len, days_elapsed, jh_1jan, &
58         year_cur, mth_cur,jD_cur, jH_cur, jD_ref, day_cur, hour, calend
59!!  USE phys_local_var_mod, ONLY : a long list of variables
60!!              ==> see below, after "CPP Keys" section
61    USE phys_state_var_mod ! Variables sauvegardees de la physique
62    USE phys_output_mod
63    USE phys_output_ctrlout_mod
64    USE print_control_mod, ONLY: mydebug=>debug , lunout, prt_level, &
65         alert_first_call, call_alert, prt_alerte
66    USE readaerosol_mod, ONLY : init_aero_fromfile
67    USE readaerosolstrato_m, ONLY : init_readaerosolstrato
68    USE radlwsw_m, only: radlwsw
69    USE regr_horiz_time_climoz_m, ONLY: regr_horiz_time_climoz
70    USE regr_pr_time_av_m, only: regr_pr_time_av
71    USE surface_data,     ONLY : type_ocean, ok_veget
72    USE time_phylmdz_mod, only: current_time, itau_phy, pdtphys, raz_date, update_time
73    USE tracinca_mod, ONLY: config_inca
74    USE tropopause_m,     ONLY: dyn_tropopause
75    USE ice_sursat_mod,  ONLY: flight_init, airplane
76    USE vampir
77    USE write_field_phy
78#ifdef CPP_XIOS
79    USE wxios, ONLY: g_ctx, wxios_set_context
80#endif
81    USE lscp_mod, ONLY : lscp
82    USE lmdz_wake_ini, ONLY : wake_ini
83    USE yamada_ini_mod, ONLY : yamada_ini
84    USE atke_turbulence_ini_mod, ONLY : atke_ini
85    USE lmdz_thermcell_ini, ONLY : thermcell_ini
86    USE lmdz_thermcell_dtke, ONLY : thermcell_dtke
87    USE blowing_snow_ini_mod, ONLY : blowing_snow_ini , qbst_bs
88    USE lscp_ini_mod, ONLY : lscp_ini
89    USE lmdz_ratqs_main, ONLY : ratqs_main
90    USE lmdz_ratqs_ini, ONLY : ratqs_ini
91
92    !USE cmp_seri_mod
93!    USE add_phys_tend_mod, only : add_pbl_tend, add_phys_tend, diag_phys_tend, prt_enerbil, &
94!  &      fl_ebil, fl_cor_ebil
95
96!!!!!!!!!!!!!!!!!! "USE" section for CPP keys !!!!!!!!!!!!!!!!!!!!!!!!
97!
98!
99#ifdef CPP_Dust
100    USE phytracr_spl_mod, ONLY: phytracr_spl, phytracr_spl_out_init
101    USE phys_output_write_spl_mod
102#else
103    USE phytrac_mod, ONLY : phytrac_init, phytrac
104    USE phys_output_write_mod
105#endif
106
107
108#ifdef INCA
109    USE geometry_mod,      ONLY: longitude, latitude, boundslon, boundslat, ind_cell_glo
110    USE time_phylmdz_mod,  ONLY: ndays
111    USE infotrac_phy,      ONLY: nqCO2
112#endif
113#ifdef REPROBUS
114    USE chem_rep, ONLY: Init_chem_rep_xjour, d_q_rep, d_ql_rep, d_qi_rep, &
115                        ptrop, ttrop, ztrop, gravit, itroprep, Z1, Z2, fac, B
116#endif
117#if defined INCA || defined REPROBUS
118    USE time_phylmdz_mod,    ONLY: annee_ref, day_ini, day_ref, start_time
119    USE vertical_layers_mod, ONLY: aps, bps, ap, bp
120#endif
121
122
123#ifdef CPP_RRTM
124    USE YOERAD, ONLY : NRADLP
125!    USE YOESW, ONLY : RSUN
126#endif
127
128
129#ifdef CPP_StratAer
130    USE strataer_local_var_mod
131    USE strataer_nuc_mod, ONLY: strataer_nuc_init
132    USE strataer_emiss_mod, ONLY: strataer_emiss_init
133#endif
134
135
136#ifdef CPP_XIOS
137    USE xios, ONLY: xios_update_calendar, xios_context_finalize
138    USE xios, ONLY: xios_get_field_attr, xios_field_is_active, xios_context
139    USE xios, ONLY: xios_set_current_context
140    USE wxios, ONLY: missing_val, missing_val_omp
141#endif
142#ifndef CPP_XIOS
143    USE paramLMDZ_phy_mod
144#endif
145
146
147#ifdef ISO
148    USE infotrac_phy, ONLY: iqIsoPha,niso, ntraciso=>ntiso, nzone
149    USE isotopes_mod, ONLY: iso_eau,iso_HDO,iso_O18,iso_O17,iso_HTO, &
150        & bidouille_anti_divergence,ok_bidouille_wake, &
151        & modif_ratqs,essai_convergence,iso_init,ridicule_rain,tnat, &
152        & ridicule,ridicule_snow
153    USE isotopes_routines_mod, ONLY: iso_tritium
154#ifdef ISOVERIF
155    USE isotopes_verif_mod, ONLY: errmax,errmaxrel, &
156        & deltalim,deltalim_snow, deltaD,deltaO,nlevmaxO17, &
157        & iso_verif_positif_nostop,iso_verif_egalite_vect2D, &
158        & iso_verif_aberrant_enc_vect2D,iso_verif_noNaN_vect2D, &
159        & iso_verif_noNaN,iso_verif_aberrant_vect2D, &
160        & iso_verif_o18_aberrant,iso_verif_aberrant_o17, &
161        & iso_verif_aberrant,iso_verif_egalite_choix, &
162        & iso_verif_aberrant_encadre,iso_verif_egalite, &
163        & iso_verif_aberrant_choix,iso_verif_positif, &
164        & iso_verif_positif_choix_vect,iso_verif_o18_aberrant_nostop, &
165        & iso_verif_init,iso_verif_aberrant_enc_choix_nostop,&
166        & iso_verif_positif_strict_nostop,iso_verif_O18_aberrant_enc_vect2D
167#endif
168#ifdef ISOTRAC
169    USE isotrac_mod, ONLY: option_traceurs,iso_traceurs_init, &
170&       option_tmin,nzone_temp
171    USE isotrac_routines_mod, only: initialise_bassins_boites, &
172        & isotrac_recolorise_general,isotrac_recolorise_conv, &
173        & iso_verif_traceur_jbid_vect,iso_verif_traceur_jbid_pos_vect
174#ifdef ISOVERIF
175    USE isotrac_routines_mod, only: iso_verif_traceur_pbidouille
176    USE isotopes_verif_mod, ONLY: iso_verif_traceur, &
177&       iso_verif_traceur_justmass,iso_verif_traceur_vect, &
178&       iso_verif_trac17_q_deltad,iso_verif_trac_masse_vect, &
179&       iso_verif_tag17_q_deltaD_vect, iso_verif_tracpos_choix_nostop
180#endif
181#endif
182#endif
183
184!
185!
186!!!!!!!!!!!!!!!!!!  END "USE" for CPP keys !!!!!!!!!!!!!!!!!!!!!!
187
188USE phys_local_var_mod, ONLY: phys_local_var_init, phys_local_var_end, &
189       ! [Variables internes non sauvegardees de la physique]
190       ! Variables locales pour effectuer les appels en serie
191       t_seri,q_seri,ql_seri,qs_seri,qbs_seri,u_seri,v_seri,tr_seri,rneb_seri, &
192       ! Dynamic tendencies (diagnostics)
193       d_t_dyn,d_q_dyn,d_ql_dyn,d_qs_dyn,d_qbs_dyn,d_u_dyn,d_v_dyn,d_tr_dyn,d_rneb_dyn, &
194       d_q_dyn2d,d_ql_dyn2d,d_qs_dyn2d,d_qbs_dyn2d, &
195       ! Physic tendencies
196       d_t_con,d_q_con,d_u_con,d_v_con, &
197       d_tr, &                              !! to be removed?? (jyg)
198       d_t_wake,d_q_wake, &
199       d_t_lwr,d_t_lw0,d_t_swr,d_t_sw0, &
200       d_t_ajsb,d_q_ajsb, &
201       d_t_ajs,d_q_ajs,d_u_ajs,d_v_ajs, &
202!       d_t_ajs_w,d_q_ajs_w, &
203!       d_t_ajs_x,d_q_ajs_x, &
204       !
205       d_t_eva,d_q_eva,d_ql_eva,d_qi_eva, &
206       d_t_lsc,d_q_lsc,d_ql_lsc,d_qi_lsc, &
207       d_t_lscst,d_q_lscst, &
208       d_t_lscth,d_q_lscth, &
209       plul_st,plul_th, &
210       !
211       d_t_vdf,d_q_vdf, d_qbs_vdf, d_u_vdf,d_v_vdf,d_t_diss, &
212       d_t_vdf_x, d_t_vdf_w, &
213       d_q_vdf_x, d_q_vdf_w, &
214       d_ts, &
215       !
216       d_t_bs,d_q_bs,d_qbs_bs, &
217       !
218!       d_t_oli,d_u_oli,d_v_oli, &
219       d_t_oro,d_u_oro,d_v_oro, &
220       d_t_oro_gw,d_u_oro_gw,d_v_oro_gw, &
221       d_t_lif,d_u_lif,d_v_lif, &
222       d_t_ec, &
223       !
224       du_gwd_hines,dv_gwd_hines,d_t_hin, &
225       dv_gwd_rando,dv_gwd_front, &
226       east_gwstress,west_gwstress, &
227       d_q_ch4, &
228       !  Special RRTM
229       ZLWFT0_i,ZSWFT0_i,ZFLDN0,  &
230       ZFLUP0,ZFSDN0,ZFSUP0,      &
231       !
232       topswad_aero,solswad_aero,   &
233       topswai_aero,solswai_aero,   &
234       topswad0_aero,solswad0_aero, &
235       !LW additional
236       toplwad_aero,sollwad_aero,   &
237       toplwai_aero,sollwai_aero,   &
238       toplwad0_aero,sollwad0_aero, &
239       !
240       topsw_aero,solsw_aero,       &
241       topsw0_aero,solsw0_aero,     &
242       topswcf_aero,solswcf_aero,   &
243       tausum_aero,tau3d_aero,      &
244       drytausum_aero,              &
245       !
246       !variables CFMIP2/CMIP5
247       topswad_aerop, solswad_aerop,   &
248       topswai_aerop, solswai_aerop,   &
249       topswad0_aerop, solswad0_aerop, &
250       topsw_aerop, topsw0_aerop,      &
251       solsw_aerop, solsw0_aerop,      &
252       topswcf_aerop, solswcf_aerop,   &
253       !LW diagnostics
254       toplwad_aerop, sollwad_aerop,   &
255       toplwai_aerop, sollwai_aerop,   &
256       toplwad0_aerop, sollwad0_aerop, &
257       !
258       ptstar, pt0, slp, &
259       !
260       bils, &
261       !
262       cldh, cldl,cldm, cldq, cldt,      &
263       JrNt,                             &
264       dthmin, evap, snowerosion,fder, plcl, plfc,   &
265       prw, prlw, prsw, prbsw,                  &
266       s_lcl, s_pblh, s_pblt, s_therm,   &
267       cdragm, cdragh,                   &
268       zustar, zu10m, zv10m, rh2m, qsat2m, &
269       zq2m, zt2m, zn2mout, weak_inversion, &
270       zt2m_min_mon, zt2m_max_mon,   &         ! pour calcul_divers.h
271       t2m_min_mon, t2m_max_mon,  &            ! pour calcul_divers.h
272       !
273       s_pblh_x, s_pblh_w, &
274       s_lcl_x, s_lcl_w,   &
275       !
276       slab_wfbils, tpot, tpote,               &
277       ue, uq, ve, vq, zxffonte,               &
278       uwat, vwat,                             &
279       zxfqcalving, zxfluxlat,                 &
280       zxrunofflic,                            &
281       zxtsol, snow_lsc, zxfqfonte, zxqsurf,   &
282       delta_qsurf,                            &
283       rain_lsc, rain_num,                     &
284       !
285       sens_x, sens_w, &
286       zxfluxlat_x, zxfluxlat_w, &
287       !
288       pbl_tke_input, tke_dissip, l_mix, wprime,&
289       t_therm, q_therm, u_therm, v_therm, &
290       cdragh_x, cdragh_w, &
291       cdragm_x, cdragm_w, &
292       kh, kh_x, kh_w, &
293       !
294       wake_k, &
295       alp_wake, &
296       wake_h, wake_omg, &
297                       ! tendencies of delta T and delta q:
298       d_deltat_wk, d_deltaq_wk, &         ! due to wakes
299       d_deltat_wk_gw, d_deltaq_wk_gw, &   ! due to wake induced gravity waves
300       d_deltat_vdf, d_deltaq_vdf, &       ! due to vertical diffusion
301       d_deltat_the, d_deltaq_the, &       ! due to thermals
302       d_deltat_ajs_cv, d_deltaq_ajs_cv, & ! due to dry adjustment of (w) before convection
303                       ! tendencies of wake fractional area and wake number per unit area:
304       d_s_wk,  d_dens_a_wk,  d_dens_wk, &  ! due to wakes
305!!!       d_s_vdf, d_dens_a_vdf, d_dens_vdf, & ! due to vertical diffusion
306!!!       d_s_the, d_dens_a_the, d_dens_the, & ! due to thermals
307       !                                 
308       ptconv, ratqsc, &
309       wbeff, convoccur, zmax_th, &
310       sens, flwp, fiwp,  &
311       alp_bl_conv,alp_bl_det,  &
312       alp_bl_fluct_m,alp_bl_fluct_tke,  &
313       alp_bl_stat, n2, s2,  &
314       proba_notrig, random_notrig,  &
315!!       cv_gen,  &  !moved to phys_state_var_mod
316       !
317       dnwd0,  &
318       omega,  &
319       epmax_diag,  &
320       !    Deep convective variables used in phytrac
321       pmflxr, pmflxs,  &
322       wdtrainA, wdtrainS, wdtrainM,  &
323       upwd, dnwd, &
324       ep,  &
325       da, mp, &
326       phi, &
327       wght_cvfd, &
328       phi2, &
329       d1a, dam, &
330       ev, &
331       elij, &
332       qtaa, &
333       clw, &
334       epmlmMm, eplaMm, &
335       sij, &
336       !
337       rneblsvol, &
338       pfraclr,pfracld, &
339       distcltop, &
340       zqsatl, zqsats, &
341       qclr, qcld, qss, qvc, rnebclr, rnebss, gamma_ss, &
342       Tcontr, qcontr, qcontr2, fcontrN, fcontrP, &
343       cldemi,  &
344       cldfra, cldtau, fiwc,  &
345       fl, re, flwc,  &
346       ref_liq, ref_ice, theta,  &
347       ref_liq_pi, ref_ice_pi,  &
348       zphi, zx_rh, zx_rhl, zx_rhi,  &
349       pmfd, pmfu,  &
350       !
351       t2m, fluxlat,  &
352       fsollw, evap_pot,  &
353       fsolsw, wfbils, wfbilo,  &
354       wfevap, wfrain, wfsnow,  & 
355       prfl, psfl,bsfl, fraca, Vprecip,  &
356       zw2,  &
357       !
358       fluxu, fluxv,  &
359       fluxt,  &
360       !
361       uwriteSTD, vwriteSTD, &                !pour calcul_STDlev.h
362       wwriteSTD, phiwriteSTD, &              !pour calcul_STDlev.h
363       qwriteSTD, twriteSTD, rhwriteSTD, &    !pour calcul_STDlev.h
364       !
365       beta_prec,  &
366       rneb,  &
367       zxsnow,snowhgt,qsnow,to_ice,sissnow,runoff,albsol3_lic, &
368       zxfluxt,zxfluxq
369
370
371#ifdef ISO
372       USE phys_local_var_mod, ONLY: xt_seri,xtl_seri,xts_seri, &
373       d_xt_eva,d_xtl_eva,d_xti_eva,           &
374       d_xt_dyn,d_xtl_dyn,d_xts_dyn,           &
375       d_xt_con, d_xt_wake,                    &
376       d_xt_ajsb, d_xt_ajs,                    &
377       d_xt_ajs_w, d_xt_ajs_x,                 &
378       d_xt_vdf, d_xt_vdf_w, d_xt_vdf_x,       &
379       d_xt_eva, d_xt_lsc, d_xtl_lsc, d_xti_lsc, &
380       d_xt_ch4,d_xt_prod_nucl,d_xt_cosmo,d_xt_decroiss, &
381       zxfxtcalving, xtevap,xtprw,    &
382       xt_therm, dxtvdf_w, dxtvdf_x,           &
383       zxfxtcalving, zxfxtfonte,       &
384       zxxtrunofflic, & ! pas besoin du runoff, c'est seulement quand on couple à sisvat
385       h1_diag,runoff_diag,xtrunoff_diag,zxxtsnow,xtVprecip,xtVprecipi, &
386       xtrain_lsc,xtsnow_lsc,pxtrfl,pxtsfl,xtclw,d_deltaxt_wk, &
387       d_deltaxt_the,d_deltaxt_vdf,d_deltaxt_ajs_cv, &
388       xtwdtrainA,xtev
389#ifdef DIAGISO
390       USE phys_local_var_mod, ONLY: fq_detrainement,fq_ddft,fq_fluxmasse,fq_evapprecip, &
391          fxt_detrainement,fxt_ddft,fxt_fluxmasse,fxt_evapprecip, &
392          f_detrainement,q_detrainement,xt_detrainement, &
393          qlp,xtlp,qvp,xtvp,q_the,xt_the
394#endif
395#endif
396       !
397
398    IMPLICIT NONE
399    !>======================================================================
400    !!
401    !! Auteur(s) Z.X. Li (LMD/CNRS) date: 19930818
402    !!
403    !! Objet: Moniteur general de la physique du modele
404    !!AA      Modifications quant aux traceurs :
405    !!AA                  -  uniformisation des parametrisations ds phytrac
406    !!AA                  -  stockage des moyennes des champs necessaires
407    !!AA                     en mode traceur off-line
408    !!======================================================================
409    !!   CLEFS CPP POUR LES IO
410    !!   =====================
411#define histNMC
412    !!======================================================================
413    !!    modif   ( P. Le Van ,  12/10/98 )
414    !!
415    !!  Arguments:
416    !!
417    !! nlon----input-I-nombre de points horizontaux
418    !! nlev----input-I-nombre de couches verticales, doit etre egale a klev
419    !! debut---input-L-variable logique indiquant le premier passage
420    !! lafin---input-L-variable logique indiquant le dernier passage
421    !! jD_cur       -R-jour courant a l'appel de la physique (jour julien)
422    !! jH_cur       -R-heure courante a l'appel de la physique (jour julien)
423    !! pdtphys-input-R-pas d'integration pour la physique (seconde)
424    !! paprs---input-R-pression pour chaque inter-couche (en Pa)
425    !! pplay---input-R-pression pour le mileu de chaque couche (en Pa)
426    !! pphi----input-R-geopotentiel de chaque couche (g z) (reference sol)
427    !! pphis---input-R-geopotentiel du sol
428    !! presnivs-input_R_pressions approximat. des milieux couches ( en PA)
429    !! u-------input-R-vitesse dans la direction X (de O a E) en m/s
430    !! v-------input-R-vitesse Y (de S a N) en m/s
431    !! t-------input-R-temperature (K)
432    !! qx------input-R-humidite specifique (kg/kg) et d'autres traceurs
433    !! d_t_dyn-input-R-tendance dynamique pour "t" (K/s)
434    !! d_q_dyn-input-R-tendance dynamique pour "q" (kg/kg/s)
435    !! d_ql_dyn-input-R-tendance dynamique pour "ql" (kg/kg/s)
436    !! d_qs_dyn-input-R-tendance dynamique pour "qs" (kg/kg/s)
437    !! flxmass_w -input-R- flux de masse verticale
438    !! d_u-----output-R-tendance physique de "u" (m/s/s)
439    !! d_v-----output-R-tendance physique de "v" (m/s/s)
440    !! d_t-----output-R-tendance physique de "t" (K/s)
441    !! d_qx----output-R-tendance physique de "qx" (kg/kg/s)
442    !! d_ps----output-R-tendance physique de la pression au sol
443    !!======================================================================
444    integer jjmp1
445    !  parameter (jjmp1=jjm+1-1/jjm) ! => (jjmp1=nbp_lat-1/(nbp_lat-1))
446    !  integer iip1
447    !  parameter (iip1=iim+1)
448
449    include "regdim.h"
450    include "dimsoil.h"
451    include "clesphys.h"
452    include "alpale.h"
453    include "dimpft.h"
454    !======================================================================
455    LOGICAL, SAVE :: ok_volcan ! pour activer les diagnostics volcaniques
456    !$OMP THREADPRIVATE(ok_volcan)
457    INTEGER, SAVE :: flag_volc_surfstrat ! pour imposer le cool/heat rate à la surf/strato
458    !$OMP THREADPRIVATE(flag_volc_surfstrat)
459    LOGICAL ok_cvl  ! pour activer le nouveau driver pour convection KE
460    PARAMETER (ok_cvl=.TRUE.)
461    LOGICAL ok_gust ! pour activer l'effet des gust sur flux surface
462    PARAMETER (ok_gust=.FALSE.)
463    INTEGER, SAVE :: iflag_radia     ! active ou non le rayonnement (MPL)
464    !$OMP THREADPRIVATE(iflag_radia)
465    !======================================================================
466    LOGICAL check ! Verifier la conservation du modele en eau
467    PARAMETER (check=.FALSE.)
468    LOGICAL ok_stratus ! Ajouter artificiellement les stratus
469    PARAMETER (ok_stratus=.FALSE.)
470    !======================================================================
471    REAL amn, amx
472    INTEGER igout
473    !======================================================================
474    ! Clef iflag_cycle_diurne controlant l'activation du cycle diurne:
475    ! en attente du codage des cles par Fred
476    ! iflag_cycle_diurne est initialise par conf_phys et se trouve
477    ! dans clesphys.h (IM)
478    !======================================================================
479    ! Modele thermique du sol, a activer pour le cycle diurne:
480    !cc      LOGICAL soil_model
481    !cc      PARAMETER (soil_model=.FALSE.)
482    !======================================================================
483    ! Dans les versions precedentes, l'eau liquide nuageuse utilisee dans
484    ! le calcul du rayonnement est celle apres la precipitation des nuages.
485    ! Si cette cle new_oliq est activee, ce sera une valeur moyenne entre
486    ! la condensation et la precipitation. Cette cle augmente les impacts
487    ! radiatifs des nuages.
488    !cc      LOGICAL new_oliq
489    !cc      PARAMETER (new_oliq=.FALSE.)
490    !======================================================================
491    ! Clefs controlant deux parametrisations de l'orographie:
492    !c      LOGICAL ok_orodr
493    !cc      PARAMETER (ok_orodr=.FALSE.)
494    !cc      LOGICAL ok_orolf
495    !cc      PARAMETER (ok_orolf=.FALSE.)
496    !======================================================================
497    LOGICAL ok_journe ! sortir le fichier journalier
498    SAVE ok_journe
499    !$OMP THREADPRIVATE(ok_journe)
500    !
501    LOGICAL ok_mensuel ! sortir le fichier mensuel
502    SAVE ok_mensuel
503    !$OMP THREADPRIVATE(ok_mensuel)
504    !
505    LOGICAL ok_instan ! sortir le fichier instantane
506    SAVE ok_instan
507    !$OMP THREADPRIVATE(ok_instan)
508    !
509    LOGICAL ok_LES ! sortir le fichier LES
510    SAVE ok_LES                           
511    !$OMP THREADPRIVATE(ok_LES)                 
512    !
513    LOGICAL callstats ! sortir le fichier stats
514    SAVE callstats                           
515    !$OMP THREADPRIVATE(callstats)                 
516    !
517    LOGICAL ok_region ! sortir le fichier regional
518    PARAMETER (ok_region=.FALSE.)
519    !======================================================================
520    REAL seuil_inversion
521    SAVE seuil_inversion
522    !$OMP THREADPRIVATE(seuil_inversion)
523   
524   
525   
526    real facteur
527
528    REAL wmax_th(klon)
529    REAL tau_overturning_th(klon)
530
531    INTEGER lmax_th(klon)
532    INTEGER limbas(klon)
533    REAL ratqscth(klon,klev)
534    REAL ratqsdiff(klon,klev)
535    REAL zqsatth(klon,klev)
536
537    !======================================================================
538    !
539    ! indices de traceurs eau vapeur, liquide, glace, fraction nuageuse LS (optional), blowing snow (optional)
540    INTEGER,SAVE :: ivap, iliq, isol, irneb, ibs
541!$OMP THREADPRIVATE(ivap, iliq, isol, irneb, ibs)
542    !
543    !
544    ! Variables argument:
545    !
546    INTEGER nlon
547    INTEGER nlev
548    REAL,INTENT(IN) :: pdtphys_
549    ! NB: pdtphys to be used in physics is in time_phylmdz_mod
550    LOGICAL debut, lafin
551    REAL paprs(klon,klev+1)
552    REAL pplay(klon,klev)
553    REAL pphi(klon,klev)
554    REAL pphis(klon)
555    REAL presnivs(klev)
556!JLD    REAL znivsig(klev)
557!JLD    real pir
558
559    REAL u(klon,klev)
560    REAL v(klon,klev)
561
562    REAL, intent(in):: rot(klon, klev)
563    ! relative vorticity, in s-1, needed for frontal waves
564
565    REAL t(klon,klev),thetal(klon,klev)
566    ! thetal: ligne suivante a decommenter si vous avez les fichiers
567    !     MPL 20130625
568    ! fth_fonctions.F90 et parkind1.F90
569    ! sinon thetal=theta
570    !     REAL fth_thetae,fth_thetav,fth_thetal
571    REAL qx(klon,klev,nqtot)
572    REAL flxmass_w(klon,klev)
573    REAL d_u(klon,klev)
574    REAL d_v(klon,klev)
575    REAL d_t(klon,klev)
576    REAL d_qx(klon,klev,nqtot)
577    REAL d_ps(klon)
578  ! variables pour tend_to_tke
579    REAL duadd(klon,klev)
580    REAL dvadd(klon,klev)
581    REAL dtadd(klon,klev)
582
583#ifndef CPP_XIOS
584    REAL, SAVE :: missing_val=nf90_fill_real
585#endif
586!!   Variables moved to phys_local_var_mod
587!!    ! Variables pour le transport convectif
588!!    real da(klon,klev),phi(klon,klev,klev),mp(klon,klev)
589!!    real wght_cvfd(klon,klev)
590!!    ! Variables pour le lessivage convectif
591!!    ! RomP >>>
592!!    real phi2(klon,klev,klev)
593!!    real d1a(klon,klev),dam(klon,klev)
594!!    real ev(klon,klev)
595!!    real clw(klon,klev),elij(klon,klev,klev)
596!!    real epmlmMm(klon,klev,klev),eplaMm(klon,klev)
597!!    ! RomP <<<
598    !IM definition dynamique o_trac dans phys_output_open
599    !      type(ctrl_out) :: o_trac(nqtot)
600
601    ! variables a une pression donnee
602    !
603    include "declare_STDlev.h"
604    !
605    !
606    include "radopt.h"
607    !
608    !
609    INTEGER n
610    !ym      INTEGER npoints
611    !ym      PARAMETER(npoints=klon)
612    !
613    INTEGER nregISCtot
614    PARAMETER(nregISCtot=1)
615    !
616    ! imin_debut, nbpti, jmin_debut, nbptj : parametres pour sorties
617    ! sur 1 region rectangulaire y compris pour 1 point
618    ! imin_debut : indice minimum de i; nbpti : nombre de points en
619    ! direction i (longitude)
620    ! jmin_debut : indice minimum de j; nbptj : nombre de points en
621    ! direction j (latitude)
622!JLD    INTEGER imin_debut, nbpti
623!JLD    INTEGER jmin_debut, nbptj
624    !IM: region='3d' <==> sorties en global
625    CHARACTER*3 region
626    PARAMETER(region='3d')
627    LOGICAL ok_hf
628    !
629    SAVE ok_hf
630    !$OMP THREADPRIVATE(ok_hf)
631
632    INTEGER, PARAMETER :: longcles=20
633    REAL, SAVE :: clesphy0(longcles)
634    !$OMP THREADPRIVATE(clesphy0)
635    !
636    ! Variables propres a la physique
637    INTEGER, SAVE :: itap         ! compteur pour la physique
638    !$OMP THREADPRIVATE(itap)
639
640    INTEGER, SAVE :: abortphy=0   ! Reprere si on doit arreter en fin de phys
641    !$OMP THREADPRIVATE(abortphy)
642    !
643    REAL,SAVE ::  solarlong0
644    !$OMP THREADPRIVATE(solarlong0)
645
646    !
647    !  Parametres de l'Orographie a l'Echelle Sous-Maille (OESM):
648    !
649    !IM 141004     REAL zulow(klon),zvlow(klon),zustr(klon), zvstr(klon)
650    REAL zulow(klon),zvlow(klon)
651    !
652    INTEGER igwd,idx(klon),itest(klon)
653    !
654    !      REAL,allocatable,save :: run_off_lic_0(:)
655    ! !$OMP THREADPRIVATE(run_off_lic_0)
656    !ym      SAVE run_off_lic_0
657    !KE43
658    ! Variables liees a la convection de K. Emanuel (sb):
659    !
660    REAL, SAVE :: bas, top             ! cloud base and top levels
661    !$OMP THREADPRIVATE(bas, top)
662    !------------------------------------------------------------------
663    ! Upmost level reached by deep convection and related variable (jyg)
664    !
665!    INTEGER izero
666    INTEGER k_upper_cv
667    !------------------------------------------------------------------
668    ! Compteur de l'occurence de cvpas=1
669    INTEGER Ncvpaseq1
670    SAVE Ncvpaseq1
671    !$OMP THREADPRIVATE(Ncvpaseq1)
672    !
673    !==========================================================================
674    !CR04.12.07: on ajoute les nouvelles variables du nouveau schema
675    !de convection avec poches froides
676    ! Variables li\'ees \`a la poche froide (jyg)
677
678!!    REAL mipsh(klon,klev)  ! mass flux shed by the adiab ascent at each level
679!!      Moved to phys_state_var_mod
680    !
681    REAL wape_prescr, fip_prescr
682    INTEGER it_wape_prescr
683    SAVE wape_prescr, fip_prescr, it_wape_prescr
684    !$OMP THREADPRIVATE(wape_prescr, fip_prescr, it_wape_prescr)
685    !
686    ! variables supplementaires de concvl
687    REAL Tconv(klon,klev)
688!!    variable moved to phys_local_var_mod
689!!    REAL sij(klon,klev,klev)
690!!    !
691!!    ! variables pour tester la conservation de l'energie dans concvl
692!!    REAL, DIMENSION(klon,klev)     :: d_t_con_sat
693!!    REAL, DIMENSION(klon,klev)     :: d_q_con_sat
694!!    REAL, DIMENSION(klon,klev)     :: dql_sat
695
696    REAL, SAVE :: alp_bl_prescr=0.
697    REAL, SAVE :: ale_bl_prescr=0.
698    REAL, SAVE :: wake_s_min_lsp=0.1
699    !$OMP THREADPRIVATE(alp_bl_prescr,ale_bl_prescr)
700    !$OMP THREADPRIVATE(wake_s_min_lsp)
701
702    REAL ok_wk_lsp(klon)
703
704    !RC
705    ! Variables li\'ees \`a la poche froide (jyg et rr)
706
707    INTEGER,  SAVE               :: iflag_wake_tend  ! wake: if =0, then wake state variables are
708                                                     ! updated within calwake
709    !$OMP THREADPRIVATE(iflag_wake_tend)
710    INTEGER,  SAVE               :: iflag_alp_wk_cond=0 ! wake: if =0, then Alp_wk is the average lifting
711                                                        ! power provided by the wakes; else, Alp_wk is the
712                                                        ! lifting power conditionned on the presence of a
713                                                        ! gust-front in the grid cell.
714    !$OMP THREADPRIVATE(iflag_alp_wk_cond)
715
716    REAL t_w(klon,klev),q_w(klon,klev) ! temperature and moisture profiles in the wake region
717    REAL t_x(klon,klev),q_x(klon,klev) ! temperature and moisture profiles in the off-wake region
718
719    REAL wake_dth(klon,klev)        ! wake : temp pot difference
720
721    REAL wake_omgbdth(klon,klev)    ! Wake : flux of Delta_Theta
722    ! transported by LS omega
723    REAL wake_dp_omgb(klon,klev)    ! Wake : vertical gradient of
724    ! large scale omega
725    REAL wake_dtKE(klon,klev)       ! Wake : differential heating
726    ! (wake - unpertubed) CONV
727    REAL wake_dqKE(klon,klev)       ! Wake : differential moistening
728    ! (wake - unpertubed) CONV
729    REAL wake_dp_deltomg(klon,klev) ! Wake : gradient vertical de wake_omg
730    REAL wake_spread(klon,klev)     ! spreading term in wake_delt
731    !
732    !pourquoi y'a pas de save??
733    !
734!!!    INTEGER, SAVE, DIMENSION(klon)   :: wake_k
735!!!    !$OMP THREADPRIVATE(wake_k)
736#ifdef ISO
737    REAL xt_x(ntraciso,klon,klev)
738    REAL xt_w(ntraciso,klon,klev)
739#endif
740    !
741    !jyg<
742    !cc      REAL wake_pe(klon)              ! Wake potential energy - WAPE
743    !>jyg
744
745    REAL wake_fip_0(klon)           ! Average Front Incoming Power (unconditionned)
746    REAL wake_gfl(klon)             ! Gust Front Length
747!!!    REAL wake_dens(klon)         ! moved to phys_state_var_mod
748    !
749    !
750    REAL dt_dwn(klon,klev)
751    REAL dq_dwn(klon,klev)
752    REAL M_dwn(klon,klev)
753    REAL M_up(klon,klev)
754    REAL dt_a(klon,klev)
755    REAL dq_a(klon,klev)
756    REAL d_t_adjwk(klon,klev)                !jyg
757    REAL d_q_adjwk(klon,klev)                !jyg
758    LOGICAL,SAVE :: ok_adjwk=.FALSE.
759    !$OMP THREADPRIVATE(ok_adjwk)
760    INTEGER,SAVE :: iflag_adjwk=0            !jyg
761    !$OMP THREADPRIVATE(iflag_adjwk)         !jyg
762    REAL,SAVE :: oliqmax=999.,oicemax=999.
763    !$OMP THREADPRIVATE(oliqmax,oicemax)
764    REAL, SAVE :: alp_offset
765    !$OMP THREADPRIVATE(alp_offset)
766    REAL, SAVE :: dtcon_multistep_max=1.e6
767    !$OMP THREADPRIVATE(dtcon_multistep_max)
768    REAL, SAVE :: dqcon_multistep_max=1.e6
769    !$OMP THREADPRIVATE(dqcon_multistep_max)
770#ifdef ISO
771    REAL dxt_dwn(ntraciso,klon,klev)
772    REAL dxt_a(ntraciso,klon,klev) 
773    REAL d_xt_adjwk(ntraciso,klon,klev)
774#endif
775 
776    !
777    !RR:fin declarations poches froides
778    !==========================================================================
779
780    REAL ztv(klon,klev),ztva(klon,klev)
781    REAL zpspsk(klon,klev)
782    REAL ztla(klon,klev),zqla(klon,klev)
783    REAL zthl(klon,klev)
784
785    !cc nrlmd le 10/04/2012
786
787    !--------Stochastic Boundary Layer Triggering: ALE_BL--------
788    !---Propri\'et\'es du thermiques au LCL
789    real zlcl_th(klon)          ! Altitude du LCL calcul\'e
790    ! continument (pcon dans
791    ! thermcell_main.F90)
792    real fraca0(klon)           ! Fraction des thermiques au LCL
793    real w0(klon)               ! Vitesse des thermiques au LCL
794    real w_conv(klon)           ! Vitesse verticale de grande \'echelle au LCL
795    real tke0(klon,klev+1)      ! TKE au d\'ebut du pas de temps
796    real therm_tke_max0(klon)   ! TKE dans les thermiques au LCL
797    real env_tke_max0(klon)     ! TKE dans l'environnement au LCL
798
799!JLD    !---D\'eclenchement stochastique
800!JLD    integer :: tau_trig(klon)
801
802    REAL,SAVE :: random_notrig_max=1.
803    !$OMP THREADPRIVATE(random_notrig_max)
804
805    !--------Statistical Boundary Layer Closure: ALP_BL--------
806    !---Profils de TKE dans et hors du thermique
807    real therm_tke_max(klon,klev)   ! Profil de TKE dans les thermiques
808    real env_tke_max(klon,klev)     ! Profil de TKE dans l'environnement
809
810    !-------Activer les tendances de TKE due a l'orograp??ie---------
811     INTEGER, SAVE :: addtkeoro
812    !$OMP THREADPRIVATE(addtkeoro)
813     REAL, SAVE :: alphatkeoro
814    !$OMP THREADPRIVATE(alphatkeoro)
815     LOGICAL, SAVE :: smallscales_tkeoro
816    !$OMP THREADPRIVATE(smallscales_tkeoro)
817
818
819
820    !cc fin nrlmd le 10/04/2012
821
822    ! Variables locales pour la couche limite (al1):
823    !
824    !Al1      REAL pblh(klon)           ! Hauteur de couche limite
825    !Al1      SAVE pblh
826    !34EK
827    !
828    ! Variables locales:
829    !
830    !AA
831    !AA  Pour phytrac
832    REAL u1(klon)             ! vents dans la premiere couche U
833    REAL v1(klon)             ! vents dans la premiere couche V
834
835    !@$$      LOGICAL offline           ! Controle du stockage ds "physique"
836    !@$$      PARAMETER (offline=.false.)
837    !@$$      INTEGER physid
838    REAL frac_impa(klon,klev) ! fractions d'aerosols lessivees (impaction)
839    REAL frac_nucl(klon,klev) ! idem (nucleation)
840    ! RomP >>>
841    REAL beta_prec_fisrt(klon,klev) ! taux de conv de l'eau cond (fisrt)
842    ! RomP <<<
843
844    !IM cf FH pour Tiedtke 080604
845    REAL rain_tiedtke(klon),snow_tiedtke(klon)
846    !
847    !IM 050204 END
848    REAL devap(klon) ! evaporation et sa derivee
849    REAL dsens(klon) ! chaleur sensible et sa derivee
850#ifdef ISO
851    REAL dxtevap(ntraciso,klon)
852#endif
853    !
854    ! Conditions aux limites
855    !
856    !
857    REAL :: day_since_equinox
858    ! Date de l'equinoxe de printemps
859    INTEGER, parameter :: mth_eq=3, day_eq=21
860    REAL :: jD_eq
861
862    LOGICAL, parameter :: new_orbit = .TRUE.
863
864    !
865    INTEGER lmt_pas
866    SAVE lmt_pas                ! frequence de mise a jour
867    !$OMP THREADPRIVATE(lmt_pas)
868    real zmasse(klon, nbp_lev),exner(klon, nbp_lev)
869    !     (column-density of mass of air in a cell, in kg m-2)
870    real, parameter:: dobson_u = 2.1415e-05 ! Dobson unit, in kg m-2
871
872    !IM sorties
873    REAL un_jour
874    PARAMETER(un_jour=86400.)
875    INTEGER itapm1 !pas de temps de la physique du(es) mois precedents
876    SAVE itapm1    !mis a jour le dernier pas de temps du mois en cours
877    !$OMP THREADPRIVATE(itapm1)
878    !======================================================================
879    !
880    ! Declaration des procedures appelees
881    !
882    EXTERNAL angle     ! calculer angle zenithal du soleil
883    EXTERNAL alboc     ! calculer l'albedo sur ocean
884    EXTERNAL ajsec     ! ajustement sec
885    EXTERNAL conlmd    ! convection (schema LMD)
886    !KE43
887    EXTERNAL conema3  ! convect4.3
888    EXTERNAL fisrtilp  ! schema de condensation a grande echelle (pluie)
889    !AA
890    ! JBM (3/14) fisrtilp_tr not loaded
891    ! EXTERNAL fisrtilp_tr ! schema de condensation a grande echelle (pluie)
892    !                          ! stockage des coefficients necessaires au
893    !                          ! lessivage OFF-LINE et ON-LINE
894    EXTERNAL hgardfou  ! verifier les temperatures
895    EXTERNAL nuage     ! calculer les proprietes radiatives
896    !C      EXTERNAL o3cm      ! initialiser l'ozone
897    EXTERNAL orbite    ! calculer l'orbite terrestre
898    EXTERNAL phyredem  ! ecrire l'etat de redemarrage de la physique
899    EXTERNAL suphel    ! initialiser certaines constantes
900    EXTERNAL transp    ! transport total de l'eau et de l'energie
901    !IM
902    EXTERNAL haut2bas  !variables de haut en bas
903    EXTERNAL ini_undefSTD  !initialise a 0 une variable a 1 niveau de pression
904    EXTERNAL undefSTD !somme les valeurs definies d'1 var a 1 niveau de pression
905    !     EXTERNAL moy_undefSTD  !moyenne d'1 var a 1 niveau de pression
906    ! EXTERNAL moyglo_aire
907    ! moyenne globale d'1 var ponderee par l'aire de la maille (moyglo_pondaire)
908    ! par la masse/airetot (moyglo_pondaima) et la vraie masse (moyglo_pondmass)
909    !
910    !
911    ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
912    ! Local variables
913    ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
914    !
915    REAL rhcl(klon,klev)    ! humiditi relative ciel clair
916    REAL dialiq(klon,klev)  ! eau liquide nuageuse
917    REAL diafra(klon,klev)  ! fraction nuageuse
918    REAL radocond(klon,klev)  ! eau condensee nuageuse
919    !
920    !XXX PB
921    REAL fluxq(klon,klev, nbsrf)   ! flux turbulent d'humidite
922    REAL fluxqbs(klon,klev, nbsrf)   ! flux turbulent de neige soufflee
923    !
924    !FC    REAL zxfluxt(klon, klev)
925    !FC    REAL zxfluxq(klon, klev)
926    REAL zxfluxqbs(klon,klev)
927    REAL zxfluxu(klon, klev)
928    REAL zxfluxv(klon, klev)
929#ifdef ISO
930    real fluxxt(ntraciso,klon,klev, nbsrf)
931    REAL zxfluxxt(ntraciso,klon, klev)
932#endif
933
934    ! Le rayonnement n'est pas calcule tous les pas, il faut donc
935    !                      sauvegarder les sorties du rayonnement
936    !ym      SAVE  heat,cool,albpla,topsw,toplw,solsw,sollw,sollwdown
937    !ym      SAVE  sollwdownclr, toplwdown, toplwdownclr
938    !ym      SAVE  topsw0,toplw0,solsw0,sollw0, heat0, cool0
939    !
940    INTEGER itaprad
941    SAVE itaprad
942    !$OMP THREADPRIVATE(itaprad)
943    !
944    REAL conv_q(klon,klev) ! convergence de l'humidite (kg/kg/s)
945    REAL conv_t(klon,klev) ! convergence de la temperature(K/s)
946    !
947#ifdef INCA
948    REAL zxsnow_dummy(klon)
949#endif
950    REAL zsav_tsol(klon)
951    !
952    REAL dist, rmu0(klon), fract(klon)
953    REAL zrmu0(klon), zfract(klon)
954    REAL zdtime, zdtime1, zdtime2, zlongi
955    !
956    REAL qcheck
957    REAL z_avant(klon), z_apres(klon), z_factor(klon)
958    LOGICAL zx_ajustq
959    !
960    REAL za
961    REAL zx_t, zx_qs, zdelta, zcor
962    real zqsat(klon,klev)
963    !
964    INTEGER i, k, iq, nsrf, l, itr
965#ifdef ISO
966    real zxt_apres(ntraciso,klon)
967    REAL xt_iso(klon,klev),xtl_iso(klon,klev)
968    REAL zxt_factor(ntraciso,klon),zxt_avant(ntraciso,klon)
969    INTEGER ixt
970    real wake_deltaq_prec(klon,klev) 
971#endif
972    !
973    REAL t_coup
974    PARAMETER (t_coup=234.0)
975
976    !ym A voir plus tard !!
977    !ym      REAL zx_relief(iim,jjmp1)
978    !ym      REAL zx_aire(iim,jjmp1)
979    !
980    ! Grandeurs de sorties
981    REAL s_capCL(klon)
982    REAL s_oliqCL(klon), s_cteiCL(klon)
983    REAL s_trmb1(klon), s_trmb2(klon)
984    REAL s_trmb3(klon)
985
986    ! La convection n'est pas calculee tous les pas, il faut donc
987    !                      sauvegarder les sorties de la convection
988    !ym      SAVE 
989    !ym      SAVE 
990    !ym      SAVE 
991    !
992    INTEGER itapcv, itapwk
993    SAVE itapcv, itapwk
994    !$OMP THREADPRIVATE(itapcv, itapwk)
995
996    !KE43
997    ! Variables locales pour la convection de K. Emanuel (sb):
998
999    REAL tvp(klon,klev)       ! virtual temp of lifted parcel
1000    CHARACTER*40 capemaxcels  !max(CAPE)
1001
1002    REAL rflag(klon)          ! flag fonctionnement de convect
1003    INTEGER iflagctrl(klon)          ! flag fonctionnement de convect
1004
1005    ! -- convect43:
1006    INTEGER ntra              ! nb traceurs pour convect4.3
1007    REAL dtvpdt1(klon,klev), dtvpdq1(klon,klev)
1008    REAL dplcldt(klon), dplcldr(klon)
1009    !?     .     condm_con(klon,klev),conda_con(klon,klev),
1010    !?     .     mr_con(klon,klev),ep_con(klon,klev)
1011    !?     .    ,sadiab(klon,klev),wadiab(klon,klev)
1012    ! --
1013    !34EK
1014    !
1015    ! Variables du changement
1016    !
1017    ! con: convection
1018    ! lsc: condensation a grande echelle (Large-Scale-Condensation)
1019    ! ajs: ajustement sec
1020    ! eva: evaporation de l'eau liquide nuageuse
1021    ! vdf: couche limite (Vertical DiFfusion)
1022    !
1023    ! tendance nulles
1024    REAL, dimension(klon,klev):: du0, dv0, dt0, dq0, dql0, dqi0, dqbs0
1025#ifdef ISO
1026    REAL, dimension(ntraciso,klon,klev):: dxt0, dxtl0, dxti0
1027#endif
1028    REAL, dimension(klon)     :: dsig0, ddens0
1029    INTEGER, dimension(klon)  :: wkoccur1
1030    ! tendance buffer pour appel de add_phys_tend
1031    REAL, DIMENSION(klon,klev)  :: d_q_ch4_dtime
1032#ifdef ISO
1033    REAL, DIMENSION(ntraciso,klon,klev)  :: d_xt_ch4_dtime
1034#endif
1035    !
1036    ! Flag pour pouvoir ne pas ajouter les tendances.
1037    ! Par defaut, les tendances doivente etre ajoutees et
1038    ! flag_inhib_tend = 0
1039    ! flag_inhib_tend > 0 : tendances non ajoutees, avec un nombre
1040    ! croissant de print quand la valeur du flag augmente
1041    !!! attention, ce flag doit etre change avec prudence !!!
1042    INTEGER :: flag_inhib_tend = 0 !  0 is the default value
1043!!    INTEGER :: flag_inhib_tend = 2
1044    !
1045    ! Logical switch to a bug : reseting to 0 convective variables at the
1046    ! begining of physiq.
1047    LOGICAL, SAVE :: ok_bug_cv_trac = .TRUE.
1048    !$OMP THREADPRIVATE(ok_bug_cv_trac)
1049    !
1050    ! Logical switch to a bug : changing wake_deltat when thermals are active
1051    ! even when there are no wakes.
1052    LOGICAL, SAVE :: ok_bug_split_th = .TRUE.
1053    !$OMP THREADPRIVATE(ok_bug_split_th)
1054
1055    !
1056    !********************************************************
1057    !     declarations
1058
1059    !********************************************************
1060    !IM 081204 END
1061    !
1062    REAL pen_u(klon,klev), pen_d(klon,klev)
1063    REAL pde_u(klon,klev), pde_d(klon,klev)
1064    INTEGER kcbot(klon), kctop(klon), kdtop(klon)
1065    !
1066    REAL ratqsbas,ratqshaut,tau_ratqs
1067    SAVE ratqsbas,ratqshaut,tau_ratqs
1068    !$OMP THREADPRIVATE(ratqsbas,ratqshaut,tau_ratqs)
1069    REAL, SAVE :: ratqsp0=50000., ratqsdp=20000.
1070    !$OMP THREADPRIVATE(ratqsp0, ratqsdp)
1071
1072    ! Parametres lies au nouveau schema de nuages (SB, PDF)
1073    REAL, SAVE :: fact_cldcon
1074    REAL, SAVE :: facttemps
1075    !$OMP THREADPRIVATE(fact_cldcon,facttemps)
1076    LOGICAL, SAVE :: ok_newmicro
1077    !$OMP THREADPRIVATE(ok_newmicro)
1078
1079    INTEGER, SAVE :: iflag_cld_th
1080    !$OMP THREADPRIVATE(iflag_cld_th)
1081!IM logical ptconv(klon,klev)  !passe dans phys_local_var_mod
1082    !IM cf. AM 081204 BEG
1083    LOGICAL ptconvth(klon,klev)
1084
1085    REAL picefra(klon,klev)
1086    REAL zrel_oro(klon)
1087    !IM cf. AM 081204 END
1088    !
1089    ! Variables liees a l'ecriture de la bande histoire physique
1090    !
1091    !======================================================================
1092    !
1093    !
1094!JLD    integer itau_w   ! pas de temps ecriture = itap + itau_phy
1095    !
1096    !
1097    ! Variables locales pour effectuer les appels en serie
1098    !
1099    !IM RH a 2m (la surface)
1100    REAL Lheat
1101
1102    INTEGER        length
1103    PARAMETER    ( length = 100 )
1104    REAL tabcntr0( length       )
1105    !
1106!JLD    INTEGER ndex2d(nbp_lon*nbp_lat)
1107    !IM
1108    !
1109    !IM AMIP2 BEG
1110!JLD    REAL moyglo, mountor
1111    !IM 141004 BEG
1112    REAL zustrdr(klon), zvstrdr(klon)
1113    REAL zustrli(klon), zvstrli(klon)
1114    REAL zustrph(klon), zvstrph(klon)
1115    REAL aam, torsfc
1116    !IM 141004 END
1117    !IM 190504 BEG
1118    !  INTEGER imp1jmp1
1119    !  PARAMETER(imp1jmp1=(iim+1)*jjmp1)
1120    !ym A voir plus tard
1121    !  REAL zx_tmp((nbp_lon+1)*nbp_lat)
1122    !  REAL airedyn(nbp_lon+1,nbp_lat)
1123    !IM 190504 END
1124!JLD    LOGICAL ok_msk
1125!JLD    REAL msk(klon)
1126    !ym A voir plus tard
1127    !ym      REAL zm_wo(jjmp1, klev)
1128    !IM AMIP2 END
1129    !
1130    REAL zx_tmp_fi2d(klon)      ! variable temporaire grille physique
1131    REAL zx_tmp_fi3d(klon,klev) ! variable temporaire pour champs 3D
1132!JLD    REAL zx_tmp_2d(nbp_lon,nbp_lat)
1133!JLD    REAL zx_lon(nbp_lon,nbp_lat)
1134!JLD    REAL zx_lat(nbp_lon,nbp_lat)
1135    !
1136    INTEGER nid_ctesGCM
1137    SAVE nid_ctesGCM
1138    !$OMP THREADPRIVATE(nid_ctesGCM)
1139    !
1140    !IM 280405 BEG
1141    !  INTEGER nid_bilKPins, nid_bilKPave
1142    !  SAVE nid_bilKPins, nid_bilKPave
1143    !  !$OMP THREADPRIVATE(nid_bilKPins, nid_bilKPave)
1144    !
1145    REAL ve_lay(klon,klev) ! transport meri. de l'energie a chaque niveau vert.
1146    REAL vq_lay(klon,klev) ! transport meri. de l'eau a chaque niveau vert.
1147    REAL ue_lay(klon,klev) ! transport zonal de l'energie a chaque niveau vert.
1148    REAL uq_lay(klon,klev) ! transport zonal de l'eau a chaque niveau vert.
1149    !
1150!JLD    REAL zjulian
1151!JLD    SAVE zjulian
1152!JLD!$OMP THREADPRIVATE(zjulian)
1153
1154!JLD    INTEGER nhori, nvert
1155!JLD    REAL zsto
1156!JLD    REAL zstophy, zout
1157
1158    CHARACTER (LEN=20) :: modname='physiq_mod'
1159    CHARACTER*80 abort_message
1160    LOGICAL, SAVE ::  ok_sync, ok_sync_omp
1161    !$OMP THREADPRIVATE(ok_sync)
1162    REAL date0
1163
1164    ! essai writephys
1165    INTEGER fid_day, fid_mth, fid_ins
1166    PARAMETER (fid_ins = 1, fid_day = 2, fid_mth = 3)
1167    INTEGER prof2d_on, prof3d_on, prof2d_av, prof3d_av
1168    PARAMETER (prof2d_on = 1, prof3d_on = 2, prof2d_av = 3, prof3d_av = 4)
1169    REAL ztsol(klon)
1170    REAL q2m(klon,nbsrf)  ! humidite a 2m
1171    REAL fsnowerosion(klon,nbsrf) ! blowing snow flux at surface
1172    REAL qbsfra  ! blowing snow fraction
1173#ifdef ISO
1174    REAL d_xtw(ntraciso),d_xtl(ntraciso), d_xts(ntraciso)
1175#endif
1176
1177    !IM: t2m, q2m, ustar, u10m, v10m et t2mincels, t2maxcels
1178    CHARACTER*40 t2mincels, t2maxcels       !t2m min., t2m max
1179    CHARACTER*40 tinst, tave
1180    REAL cldtaupi(klon,klev) ! Cloud optical thickness for
1181    ! pre-industrial (pi) aerosols
1182
1183    INTEGER :: naero
1184    ! Aerosol optical properties
1185    CHARACTER*4, DIMENSION(naero_grp) :: rfname
1186    REAL, DIMENSION(klon,klev)     :: mass_solu_aero ! total mass
1187    ! concentration
1188    ! for all soluble
1189    ! aerosols[ug/m3]
1190    REAL, DIMENSION(klon,klev)     :: mass_solu_aero_pi
1191    ! - " - (pre-industrial value)
1192
1193    ! Parameters
1194    LOGICAL ok_ade, ok_aie    ! Apply aerosol (in)direct effects or not
1195    LOGICAL ok_alw            ! Apply aerosol LW effect or not
1196    LOGICAL ok_cdnc ! ok cloud droplet number concentration (O. Boucher 01-2013)
1197    REAL bl95_b0, bl95_b1   ! Parameter in Boucher and Lohmann (1995)
1198    SAVE ok_ade, ok_aie, ok_alw, ok_cdnc, bl95_b0, bl95_b1
1199    !$OMP THREADPRIVATE(ok_ade, ok_aie, ok_alw, ok_cdnc, bl95_b0, bl95_b1)
1200    LOGICAL, SAVE :: aerosol_couple ! true  : calcul des aerosols dans INCA
1201    ! false : lecture des aerosol dans un fichier
1202    !$OMP THREADPRIVATE(aerosol_couple)   
1203    LOGICAL, SAVE :: chemistry_couple ! true  : use INCA chemistry O3
1204    ! false : use offline chemistry O3
1205    !$OMP THREADPRIVATE(chemistry_couple)   
1206    INTEGER, SAVE :: flag_aerosol
1207    !$OMP THREADPRIVATE(flag_aerosol)
1208    LOGICAL, SAVE :: flag_bc_internal_mixture
1209    !$OMP THREADPRIVATE(flag_bc_internal_mixture)
1210    !
1211    !--STRAT AEROSOL
1212    INTEGER, SAVE :: flag_aerosol_strat
1213    !$OMP THREADPRIVATE(flag_aerosol_strat)
1214    !
1215    !--INTERACTIVE AEROSOL FEEDBACK ON RADIATION
1216    LOGICAL, SAVE :: flag_aer_feedback
1217    !$OMP THREADPRIVATE(flag_aer_feedback)
1218
1219    !c-fin STRAT AEROSOL
1220    !
1221    ! Declaration des constantes et des fonctions thermodynamiques
1222    !
1223    LOGICAL,SAVE :: first=.TRUE.
1224    !$OMP THREADPRIVATE(first)
1225
1226    ! VARIABLES RELATED TO OZONE CLIMATOLOGIES ; all are OpenMP shared
1227    ! Note that pressure vectors are in Pa and in stricly ascending order
1228    INTEGER,SAVE :: read_climoz                ! Read ozone climatology
1229    !     (let it keep the default OpenMP shared attribute)
1230    !     Allowed values are 0, 1 and 2
1231    !     0: do not read an ozone climatology
1232    !     1: read a single ozone climatology that will be used day and night
1233    !     2: read two ozone climatologies, the average day and night
1234    !     climatology and the daylight climatology
1235    INTEGER,SAVE :: ncid_climoz                ! NetCDF file identifier
1236    REAL, ALLOCATABLE, SAVE :: press_cen_climoz(:) ! Pressure levels
1237    REAL, ALLOCATABLE, SAVE :: press_edg_climoz(:) ! Edges of pressure intervals
1238    REAL, ALLOCATABLE, SAVE :: time_climoz(:)      ! Time vector
1239    CHARACTER(LEN=13), PARAMETER :: vars_climoz(2) &
1240                                  = ["tro3         ","tro3_daylight"]
1241    ! vars_climoz(1:read_climoz): variables names in climoz file.
1242    ! vars_climoz(1:read_climoz-2) if read_climoz>2 (temporary)
1243    REAL :: ro3i ! 0<=ro3i<=360 ; required time index in NetCDF file for
1244                 ! the ozone fields, old method.
1245
1246    include "YOMCST.h"
1247    include "YOETHF.h"
1248    include "FCTTRE.h"
1249    !IM 100106 BEG : pouvoir sortir les ctes de la physique
1250    include "conema3.h"
1251    include "nuage.h"
1252    include "compbl.h"
1253    !IM 100106 END : pouvoir sortir les ctes de la physique
1254    !
1255    ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1256    ! Declarations pour Simulateur COSP
1257    !============================================================
1258    real :: mr_ozone(klon,klev), phicosp(klon,klev)
1259
1260    !IM stations CFMIP
1261    INTEGER, SAVE :: nCFMIP
1262    !$OMP THREADPRIVATE(nCFMIP)
1263    INTEGER, PARAMETER :: npCFMIP=120
1264    INTEGER, ALLOCATABLE, SAVE :: tabCFMIP(:)
1265    REAL, ALLOCATABLE, SAVE :: lonCFMIP(:), latCFMIP(:)
1266    !$OMP THREADPRIVATE(tabCFMIP, lonCFMIP, latCFMIP)
1267    INTEGER, ALLOCATABLE, SAVE :: tabijGCM(:)
1268    REAL, ALLOCATABLE, SAVE :: lonGCM(:), latGCM(:)
1269    !$OMP THREADPRIVATE(tabijGCM, lonGCM, latGCM)
1270    INTEGER, ALLOCATABLE, SAVE :: iGCM(:), jGCM(:)
1271    !$OMP THREADPRIVATE(iGCM, jGCM)
1272    logical, dimension(nfiles)            :: phys_out_filestations
1273    logical, parameter :: lNMC=.FALSE.
1274
1275    !IM betaCRF
1276    REAL, SAVE :: pfree, beta_pbl, beta_free
1277    !$OMP THREADPRIVATE(pfree, beta_pbl, beta_free)
1278    REAL, SAVE :: lon1_beta,  lon2_beta, lat1_beta, lat2_beta
1279    !$OMP THREADPRIVATE(lon1_beta,  lon2_beta, lat1_beta, lat2_beta)
1280    LOGICAL, SAVE :: mskocean_beta
1281    !$OMP THREADPRIVATE(mskocean_beta)
1282    REAL, dimension(klon, klev) :: beta ! facteur sur cldtaurad et
1283    ! cldemirad pour evaluer les
1284    ! retros liees aux CRF
1285    REAL, dimension(klon, klev) :: cldtaurad   ! epaisseur optique
1286    ! pour radlwsw pour
1287    ! tester "CRF off"
1288    REAL, dimension(klon, klev) :: cldtaupirad ! epaisseur optique
1289    ! pour radlwsw pour
1290    ! tester "CRF off"
1291    REAL, dimension(klon, klev) :: cldemirad   ! emissivite pour
1292    ! radlwsw pour tester
1293    ! "CRF off"
1294    REAL, dimension(klon, klev) :: cldfrarad   ! fraction nuageuse
1295
1296#ifdef INCA
1297    REAL :: calday, zxsnow_dummy(klon)
1298    ! set de variables utilisees pour l'initialisation des valeurs provenant de INCA
1299    REAL, DIMENSION(klon,klev,naero_grp,nbands) :: init_tauinca
1300    REAL, DIMENSION(klon,klev,naero_grp,nbands) :: init_pizinca
1301    REAL, DIMENSION(klon,klev,naero_grp,nbands) :: init_cginca
1302    REAL, DIMENSION(klon,klev,nbands) :: init_ccminca
1303#endif
1304    REAL, DIMENSION(klon,nbtr) :: init_source
1305
1306    !lwoff=y : offset LW CRE for radiation code and other schemes
1307    REAL, SAVE :: betalwoff
1308    !OMP THREADPRIVATE(betalwoff)
1309!
1310    INTEGER :: nbtr_tmp ! Number of tracer inside concvl
1311    REAL, dimension(klon,klev) :: sh_in ! Specific humidity entering in phytrac
1312    REAL, dimension(klon,klev) :: ch_in ! Condensed humidity entering in phytrac (eau liquide)
1313    integer iostat
1314
1315    REAL, dimension(klon,klev+1) :: tke_dissip_ave, l_mix_ave, wprime_ave
1316    REAL zzz
1317    !albedo SB >>>
1318    REAL,DIMENSION(6), SAVE :: SFRWL
1319!$OMP THREADPRIVATE(SFRWL)
1320    !albedo SB <<<
1321
1322    !--OB variables for mass fixer (hard coded for now)
1323    LOGICAL, PARAMETER :: mass_fixer=.FALSE.
1324    REAL qql1(klon),qql2(klon),corrqql
1325#ifdef ISO
1326    real xtql1(ntraciso,klon),xtql2(ntraciso,klon),corrxtql(ntraciso)
1327#endif
1328
1329    REAL pi
1330
1331    pi = 4. * ATAN(1.)
1332
1333    ! set-up call to alerte function
1334    call_alert = (alert_first_call .AND. is_master)
1335   
1336    ! Ehouarn: set value of jjmp1 since it is no longer a "fixed parameter"
1337    jjmp1=nbp_lat
1338
1339    !======================================================================
1340    ! Gestion calendrier : mise a jour du module phys_cal_mod
1341    !
1342    pdtphys=pdtphys_
1343    CALL update_time(pdtphys)
1344    phys_tstep=NINT(pdtphys)
1345#ifdef CPP_XIOS
1346! switch to XIOS LMDZ physics context
1347    IF (.NOT. debut .AND. is_omp_master) THEN
1348       CALL wxios_set_context()
1349       CALL xios_update_calendar(itap+1)
1350    ENDIF
1351#endif
1352
1353    !======================================================================
1354    ! Ecriture eventuelle d'un profil verticale en entree de la physique.
1355    ! Utilise notamment en 1D mais peut etre active egalement en 3D
1356    ! en imposant la valeur de igout.
1357    !======================================================================d
1358    IF (prt_level.ge.1) THEN
1359       igout=klon/2+1/klon
1360       write(lunout,*) 'DEBUT DE PHYSIQ !!!!!!!!!!!!!!!!!!!!'
1361       write(lunout,*) 'igout, lat, lon ',igout, latitude_deg(igout), &
1362            longitude_deg(igout)
1363       write(lunout,*) &
1364            'nlon,klev,nqtot,debut,lafin, jD_cur, jH_cur,pdtphys'
1365       write(lunout,*) &
1366            nlon,klev,nqtot,debut,lafin, jD_cur, jH_cur,pdtphys
1367
1368       write(lunout,*) 'paprs, play, phi, u, v, t'
1369       DO k=1,klev
1370          write(lunout,*) paprs(igout,k),pplay(igout,k),pphi(igout,k), &
1371               u(igout,k),v(igout,k),t(igout,k)
1372       ENDDO
1373       write(lunout,*) 'ovap (g/kg),  oliq (g/kg)'
1374       DO k=1,klev
1375          write(lunout,*) qx(igout,k,1)*1000,qx(igout,k,2)*1000.
1376       ENDDO
1377    ENDIF
1378
1379    ! Quick check on pressure levels:
1380    CALL assert(paprs(:, nbp_lev + 1) < paprs(:, nbp_lev), &
1381            "physiq_mod paprs bad order")
1382
1383    IF (first) THEN
1384       ivap = strIdx(tracers(:)%name, addPhase('H2O', 'g'))
1385       iliq = strIdx(tracers(:)%name, addPhase('H2O', 'l'))
1386       isol = strIdx(tracers(:)%name, addPhase('H2O', 's'))
1387       irneb= strIdx(tracers(:)%name, addPhase('H2O', 'r'))
1388       ibs  = strIdx(tracers(:)%name, addPhase('H2O', 'b'))
1389       CALL init_etat0_limit_unstruct
1390       IF (.NOT. create_etat0_limit) CALL init_limit_read(days_elapsed)
1391       !CR:nvelles variables convection/poches froides
1392
1393       WRITE(lunout,*) '================================================='
1394       WRITE(lunout,*) 'Allocation des variables locales et sauvegardees'
1395       WRITE(lunout,*) '================================================='
1396       CALL phys_local_var_init
1397       !
1398       !     appel a la lecture du run.def physique
1399       CALL conf_phys(ok_journe, ok_mensuel, &
1400            ok_instan, ok_hf, &
1401            ok_LES, &
1402            callstats, &
1403            solarlong0,seuil_inversion, &
1404            fact_cldcon, facttemps,ok_newmicro,iflag_radia, &
1405            iflag_cld_th,ratqsbas,ratqshaut,tau_ratqs, &
1406            ok_ade, ok_aie, ok_alw, ok_cdnc, ok_volcan, flag_volc_surfstrat, aerosol_couple, &
1407            chemistry_couple, flag_aerosol, flag_aerosol_strat, flag_aer_feedback, &
1408            flag_bc_internal_mixture, bl95_b0, bl95_b1, &
1409                                ! nv flags pour la convection et les
1410                                ! poches froides
1411            read_climoz, &
1412            alp_offset)
1413       CALL phys_state_var_init(read_climoz)
1414       CALL phys_output_var_init
1415       IF (read_climoz>=1 .AND. create_etat0_limit .AND. grid_type==unstructured) &
1416          CALL regr_horiz_time_climoz(read_climoz,ok_daily_climoz)
1417
1418#ifdef CPP_StratAer
1419       CALL strataer_init
1420       CALL strataer_nuc_init
1421       CALL strataer_emiss_init
1422#endif
1423
1424       print*, '================================================='
1425       !
1426       !CR: check sur le nb de traceurs de l eau
1427       IF ((iflag_ice_thermo.gt.0).and.(nqo==2)) THEN
1428          WRITE (lunout, *) ' iflag_ice_thermo==1 requires 3 H2O tracers ', &
1429               '(H2O_g, H2O_l, H2O_s) but nqo=', nqo, '. Might as well stop here.'
1430          abort_message='see above'
1431          CALL abort_physic(modname,abort_message,1)
1432       ENDIF
1433
1434       IF (ok_ice_sursat.AND.(iflag_ice_thermo.EQ.0)) THEN
1435          WRITE (lunout, *) ' ok_ice_sursat=y requires iflag_ice_thermo=1 as well'
1436          abort_message='see above'
1437          CALL abort_physic(modname,abort_message,1)
1438       ENDIF
1439
1440       IF (ok_ice_sursat.AND.(nqo.LT.4)) THEN
1441          WRITE (lunout, *) ' ok_ice_sursat=y requires 4 H2O tracers ', &
1442               '(H2O_g, H2O_l, H2O_s, H2O_r) but nqo=', nqo, '. Might as well stop here.'
1443          abort_message='see above'
1444          CALL abort_physic(modname,abort_message,1)
1445       ENDIF
1446
1447       IF (ok_plane_h2o.AND..NOT.ok_ice_sursat) THEN
1448          WRITE (lunout, *) ' ok_plane_h2o=y requires ok_ice_sursat=y '
1449          abort_message='see above'
1450          CALL abort_physic(modname,abort_message,1)
1451       ENDIF
1452
1453       IF (ok_plane_contrail.AND..NOT.ok_ice_sursat) THEN
1454          WRITE (lunout, *) ' ok_plane_contrail=y requires ok_ice_sursat=y '
1455          abort_message='see above'
1456          CALL abort_physic(modname,abort_message,1)
1457       ENDIF
1458
1459        IF (ok_bs) THEN
1460          abort_message='blowing snow cannot be activated with water isotopes yet'
1461          CALL abort_physic(modname,abort_message, 1)
1462         IF ((ok_ice_sursat.AND.nqo .LT.5).OR.(.NOT.ok_ice_sursat.AND.nqo.LT.4)) THEN
1463             WRITE (lunout, *) 'activation of blowing snow needs a specific H2O tracer', &
1464                               'but nqo=', nqo
1465             abort_message='see above'
1466             CALL abort_physic(modname,abort_message, 1)
1467         ENDIF
1468        ENDIF
1469       Ncvpaseq1 = 0
1470       dnwd0=0.0
1471       ftd=0.0
1472       fqd=0.0
1473       cin=0.
1474       !ym Attention pbase pas initialise dans concvl !!!!
1475       pbase=0
1476       !IM 180608
1477
1478       itau_con=0
1479       first=.FALSE.
1480
1481    ENDIF  ! first
1482
1483    !ym => necessaire pour iflag_con != 2   
1484    pmfd(:,:) = 0.
1485    pen_u(:,:) = 0.
1486    pen_d(:,:) = 0.
1487    pde_d(:,:) = 0.
1488    pde_u(:,:) = 0.
1489    aam=0.
1490    d_t_adjwk(:,:)=0
1491    d_q_adjwk(:,:)=0
1492#ifdef ISO
1493    d_xt_adjwk(:,:,:)=0
1494#endif
1495    alp_bl_conv(:)=0.
1496
1497    torsfc=0.
1498    forall (k=1: nbp_lev) zmasse(:, k) = (paprs(:, k)-paprs(:, k+1)) / rg
1499
1500
1501    IF (debut) THEN
1502       CALL suphel ! initialiser constantes et parametres phys.
1503! tau_gl : constante de rappel de la temperature a la surface de la glace - en
1504       tau_gl=5.
1505       CALL getin_p('tau_gl', tau_gl)
1506! tau_gl : constante de rappel de la temperature a la surface de la glace - en
1507! secondes
1508       tau_gl=86400.*tau_gl
1509       WRITE(lunout,*) 'debut physiq_mod tau_gl=',tau_gl
1510
1511       CALL getin_p('iflag_alp_wk_cond', iflag_alp_wk_cond)
1512       CALL getin_p('random_notrig_max',random_notrig_max)
1513       CALL getin_p('ok_adjwk',ok_adjwk)
1514       IF (ok_adjwk) iflag_adjwk=2  ! for compatibility with older versions
1515       ! iflag_adjwk: ! 0 = Default: no convective adjustment of w-region
1516                      ! 1 => convective adjustment but state variables are unchanged
1517                      ! 2 => convective adjustment and state variables are changed
1518       CALL getin_p('iflag_adjwk',iflag_adjwk)
1519       CALL getin_p('dtcon_multistep_max',dtcon_multistep_max)
1520       CALL getin_p('dqcon_multistep_max',dqcon_multistep_max)
1521       CALL getin_p('oliqmax',oliqmax)
1522       CALL getin_p('oicemax',oicemax)
1523       CALL getin_p('ratqsp0',ratqsp0)
1524       CALL getin_p('ratqsdp',ratqsdp)
1525       iflag_wake_tend = 0
1526       CALL getin_p('iflag_wake_tend',iflag_wake_tend)
1527       ok_bad_ecmwf_thermo=.TRUE. ! By default thermodynamical constants are set
1528                                  ! in rrtm/suphec.F90 (and rvtmp2 is set to 0).
1529       CALL getin_p('ok_bad_ecmwf_thermo',ok_bad_ecmwf_thermo)
1530       CALL getin_p('ok_bug_cv_trac',ok_bug_cv_trac)
1531       CALL getin_p('ok_bug_split_th',ok_bug_split_th)
1532       fl_ebil = 0 ! by default, conservation diagnostics are desactivated
1533       CALL getin_p('fl_ebil',fl_ebil)
1534       fl_cor_ebil = 0 ! by default, no correction to ensure energy conservation
1535       CALL getin_p('fl_cor_ebil',fl_cor_ebil)
1536       iflag_phytrac = 1 ! by default we do want to call phytrac
1537       CALL getin_p('iflag_phytrac',iflag_phytrac)
1538#ifdef CPP_Dust
1539       IF (iflag_phytrac.EQ.0) THEN
1540         WRITE(lunout,*) 'In order to run with SPLA, iflag_phytrac will be forced to 1'
1541         iflag_phytrac = 1
1542       ENDIF
1543#endif
1544       nvm_lmdz = 13
1545       CALL getin_p('NVM',nvm_lmdz)
1546
1547       WRITE(lunout,*) 'iflag_alp_wk_cond=',  iflag_alp_wk_cond
1548       WRITE(lunout,*) 'random_ntrig_max=',   random_notrig_max
1549       WRITE(lunout,*) 'ok_adjwk=',           ok_adjwk
1550       WRITE(lunout,*) 'iflag_adjwk=',        iflag_adjwk
1551       WRITE(lunout,*) 'qtcon_multistep_max=',dtcon_multistep_max
1552       WRITE(lunout,*) 'qdcon_multistep_max=',dqcon_multistep_max
1553       WRITE(lunout,*) 'ratqsp0=',            ratqsp0
1554       WRITE(lunout,*) 'ratqsdp=',            ratqsdp
1555       WRITE(lunout,*) 'iflag_wake_tend=',    iflag_wake_tend
1556       WRITE(lunout,*) 'ok_bad_ecmwf_thermo=',ok_bad_ecmwf_thermo
1557       WRITE(lunout,*) 'ok_bug_cv_trac=',     ok_bug_cv_trac
1558       WRITE(lunout,*) 'ok_bug_split_th=',    ok_bug_split_th
1559       WRITE(lunout,*) 'fl_ebil=',            fl_ebil
1560       WRITE(lunout,*) 'fl_cor_ebil=',        fl_cor_ebil
1561       WRITE(lunout,*) 'iflag_phytrac=',      iflag_phytrac
1562       WRITE(lunout,*) 'NVM=',                nvm_lmdz
1563
1564       !--PC: defining fields to be exchanged between LMDz, ORCHIDEE and NEMO
1565       WRITE(lunout,*) 'Call to infocfields from physiq'
1566       CALL infocfields_init
1567
1568    ENDIF
1569
1570    IF (prt_level.ge.1) print *,'CONVERGENCE PHYSIQUE THERM 1 '
1571
1572    !======================================================================
1573    ! Gestion calendrier : mise a jour du module phys_cal_mod
1574    !
1575    !     CALL phys_cal_update(jD_cur,jH_cur)
1576
1577    !
1578    ! Si c'est le debut, il faut initialiser plusieurs choses
1579    !          ********
1580    !
1581    IF (debut) THEN
1582       !rv CRinitialisation de wght_th et lalim_conv pour la
1583       !definition de la couche alimentation de la convection a partir
1584       !des caracteristiques du thermique
1585       wght_th(:,:)=1.
1586       lalim_conv(:)=1
1587       !RC
1588       ustar(:,:)=0.
1589!       u10m(:,:)=0.
1590!       v10m(:,:)=0.
1591       rain_con(:)=0.
1592       snow_con(:)=0.
1593       topswai(:)=0.
1594       topswad(:)=0.
1595       solswai(:)=0.
1596       solswad(:)=0.
1597
1598       wmax_th(:)=0.
1599       tau_overturning_th(:)=0.
1600
1601       IF (ANY(type_trac == ['inca','inco'])) THEN
1602          ! jg : initialisation jusqu'au ces variables sont dans restart
1603          ccm(:,:,:) = 0.
1604          tau_aero(:,:,:,:) = 0.
1605          piz_aero(:,:,:,:) = 0.
1606          cg_aero(:,:,:,:) = 0.
1607          d_q_ch4(:,:) = 0.
1608
1609          config_inca='none' ! default
1610          CALL getin_p('config_inca',config_inca)
1611
1612       ELSE
1613          config_inca='none' ! default
1614       ENDIF
1615
1616       tau_aero(:,:,:,:) = 1.e-15
1617       piz_aero(:,:,:,:) = 1.
1618       cg_aero(:,:,:,:)  = 0.
1619       d_q_ch4(:,:) = 0.
1620
1621       IF (aerosol_couple .AND. (config_inca /= "aero" &
1622            .AND. config_inca /= "aeNP ")) THEN
1623          abort_message &
1624               = 'if aerosol_couple is activated, config_inca need to be ' &
1625               // 'aero or aeNP'
1626          CALL abort_physic (modname,abort_message,1)
1627       ENDIF
1628
1629       rnebcon0(:,:) = 0.0
1630       clwcon0(:,:) = 0.0
1631       rnebcon(:,:) = 0.0
1632       clwcon(:,:) = 0.0
1633
1634       !
1635       print*,'iflag_coupl,iflag_clos,iflag_wake', &
1636            iflag_coupl,iflag_clos,iflag_wake
1637       print*,'iflag_cycle_diurne', iflag_cycle_diurne
1638       !
1639       IF (iflag_con.EQ.2.AND.iflag_cld_th.GT.-1) THEN
1640          abort_message = 'Tiedtke needs iflag_cld_th=-2 or -1'
1641          CALL abort_physic (modname,abort_message,1)
1642       ENDIF
1643       !
1644       !
1645       ! Initialiser les compteurs:
1646       !
1647       itap    = 0
1648       itaprad = 0
1649       itapcv = 0
1650       itapwk = 0
1651
1652! C Risi: vérifier compatibilité des options isotopiques entre
1653! dyn3dmem et physiq
1654#ifdef ISO
1655    write(*,*) 'physiq 1846a: ok_isotopes,ntraciso,niso=',niso>0,ntraciso,niso
1656    IF(niso  <= 0) CALL abort_gcm('physiq 1756','options iso incompatibles',1)
1657#ifdef ISOTRAC
1658    IF(nzone <= 0) CALL abort_gcm('physiq 1758','options isotrac incompatibles',1)
1659#else
1660    IF(nzone  > 0) CALL abort_gcm('physiq 1762','options isotrac incompatibles',1)
1661#endif
1662#else
1663    if(niso   > 0) CALL abort_gcm('physiq 1772','options iso incompatibles',1)
1664#endif
1665
1666#ifdef ISO
1667        ! initialisations isotopiques
1668#ifdef ISOVERIF
1669           write(*,*) 'physiq 1366: call iso_init'
1670           write(*,*) 'ok_isotopes=',niso > 0
1671#endif
1672        if (niso > 0) call iso_init()
1673#ifdef ISOTRAC
1674IF(nzone > 0) then
1675        write(*,*) 'physiq 1416: call iso_traceurs_init'
1676        call iso_traceurs_init()
1677endif
1678#endif
1679!write(*,*) 'gcm 265: ntraciso=',ntraciso
1680#ifdef ISOVERIF
1681        write(*,*) 'physiq 1421: call iso_verif_init'
1682        call iso_verif_init()
1683#endif
1684#endif
1685
1686
1687       ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1688       !! Un petit travail \`a faire ici.
1689       ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1690
1691       IF (iflag_pbl>1) THEN
1692          PRINT*, "Using method MELLOR&YAMADA"
1693       ENDIF
1694
1695       ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1696       ! FH 2008/05/02 changement lie a la lecture de nbapp_rad dans
1697       ! phylmd plutot que dyn3d
1698       ! Attention : la version precedente n'etait pas tres propre.
1699       ! Il se peut qu'il faille prendre une valeur differente de nbapp_rad
1700       ! pour obtenir le meme resultat.
1701!jyg for fh<
1702       WRITE(lunout,*) 'Pas de temps phys_tstep pdtphys ',phys_tstep,pdtphys
1703       IF (abs(phys_tstep-pdtphys)>1.e-10) THEN
1704          abort_message='pas de temps doit etre entier en seconde pour orchidee et XIOS'
1705          CALL abort_physic(modname,abort_message,1)
1706       ENDIF
1707!>jyg
1708       IF (MOD(NINT(86400./phys_tstep),nbapp_rad).EQ.0) THEN
1709          radpas = NINT( 86400./phys_tstep)/nbapp_rad
1710       ELSE
1711          WRITE(lunout,*) 'le nombre de pas de temps physique doit etre un ', &
1712               'multiple de nbapp_rad'
1713          WRITE(lunout,*) 'changer nbapp_rad ou alors commenter ce test ', &
1714               'mais 1+1<>2'
1715          abort_message='nbre de pas de temps physique n est pas multiple ' &
1716               // 'de nbapp_rad'
1717          CALL abort_physic(modname,abort_message,1)
1718       ENDIF
1719       IF (nbapp_cv .EQ. 0) nbapp_cv=86400./phys_tstep
1720       IF (nbapp_wk .EQ. 0) nbapp_wk=86400./phys_tstep
1721       print *,'physiq, nbapp_cv, nbapp_wk ',nbapp_cv,nbapp_wk
1722       IF (MOD(NINT(86400./phys_tstep),nbapp_cv).EQ.0) THEN
1723          cvpas_0 = NINT( 86400./phys_tstep)/nbapp_cv
1724          cvpas = cvpas_0
1725       print *,'physiq, cvpas ',cvpas
1726       ELSE
1727          WRITE(lunout,*) 'le nombre de pas de temps physique doit etre un ', &
1728               'multiple de nbapp_cv'
1729          WRITE(lunout,*) 'changer nbapp_cv ou alors commenter ce test ', &
1730               'mais 1+1<>2'
1731          abort_message='nbre de pas de temps physique n est pas multiple ' &
1732               // 'de nbapp_cv'
1733          CALL abort_physic(modname,abort_message,1)
1734       ENDIF
1735       IF (MOD(NINT(86400./phys_tstep),nbapp_wk).EQ.0) THEN
1736          wkpas = NINT( 86400./phys_tstep)/nbapp_wk
1737!       print *,'physiq, wkpas ',wkpas
1738       ELSE
1739          WRITE(lunout,*) 'le nombre de pas de temps physique doit etre un ', &
1740               'multiple de nbapp_wk'
1741          WRITE(lunout,*) 'changer nbapp_wk ou alors commenter ce test ', &
1742               'mais 1+1<>2'
1743          abort_message='nbre de pas de temps physique n est pas multiple ' &
1744               // 'de nbapp_wk'
1745          CALL abort_physic(modname,abort_message,1)
1746       ENDIF
1747       ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1748       CALL init_iophy_new(latitude_deg,longitude_deg)
1749
1750          !===================================================================
1751          !IM stations CFMIP
1752          nCFMIP=npCFMIP
1753          OPEN(98,file='npCFMIP_param.data',status='old', &
1754               form='formatted',iostat=iostat)
1755          IF (iostat == 0) THEN
1756             READ(98,*,end=998) nCFMIP
1757998          CONTINUE
1758             CLOSE(98)
1759             CONTINUE
1760             IF(nCFMIP.GT.npCFMIP) THEN
1761                print*,'nCFMIP > npCFMIP : augmenter npCFMIP et recompiler'
1762                CALL abort_physic("physiq", "", 1)
1763             ELSE
1764                print*,'physiq npCFMIP=',npCFMIP,'nCFMIP=',nCFMIP
1765             ENDIF
1766
1767             !
1768             ALLOCATE(tabCFMIP(nCFMIP))
1769             ALLOCATE(lonCFMIP(nCFMIP), latCFMIP(nCFMIP))
1770             ALLOCATE(tabijGCM(nCFMIP))
1771             ALLOCATE(lonGCM(nCFMIP), latGCM(nCFMIP))
1772             ALLOCATE(iGCM(nCFMIP), jGCM(nCFMIP))
1773             !
1774             ! lecture des nCFMIP stations CFMIP, de leur numero
1775             ! et des coordonnees geographiques lonCFMIP, latCFMIP
1776             !
1777             CALL read_CFMIP_point_locations(nCFMIP, tabCFMIP,  &
1778                  lonCFMIP, latCFMIP)
1779             !
1780             ! identification des
1781             ! 1) coordonnees lonGCM, latGCM des points CFMIP dans la
1782             ! grille de LMDZ
1783             ! 2) indices points tabijGCM de la grille physique 1d sur
1784             ! klon points
1785             ! 3) indices iGCM, jGCM de la grille physique 2d
1786             !
1787             CALL LMDZ_CFMIP_point_locations(nCFMIP, lonCFMIP, latCFMIP, &
1788                  tabijGCM, lonGCM, latGCM, iGCM, jGCM)
1789             !
1790          ELSE
1791             ALLOCATE(tabijGCM(0))
1792             ALLOCATE(lonGCM(0), latGCM(0))
1793             ALLOCATE(iGCM(0), jGCM(0))
1794          ENDIF
1795
1796#ifdef CPP_IOIPSL
1797
1798       !$OMP MASTER
1799       ! FH : if ok_sync=.true. , the time axis is written at each time step
1800       ! in the output files. Only at the end in the opposite case
1801       ok_sync_omp=.FALSE.
1802       CALL getin('ok_sync',ok_sync_omp)
1803       CALL phys_output_open(longitude_deg,latitude_deg,nCFMIP,tabijGCM, &
1804            iGCM,jGCM,lonGCM,latGCM, &
1805            jjmp1,nlevSTD,clevSTD,rlevSTD, phys_tstep,ok_veget, &
1806            type_ocean,iflag_pbl,iflag_pbl_split,ok_mensuel,ok_journe, &
1807            ok_hf,ok_instan,ok_LES,ok_ade,ok_aie, &
1808            read_climoz, phys_out_filestations, &
1809            aerosol_couple, &
1810            flag_aerosol_strat, pdtphys, paprs, pphis,  &
1811            pplay, lmax_th, ptconv, ptconvth, ivap,  &
1812            d_u, d_t, qx, d_qx, zmasse, ok_sync_omp)
1813       !$OMP END MASTER
1814       !$OMP BARRIER
1815       ok_sync=ok_sync_omp
1816
1817       freq_outNMC(1) = ecrit_files(7)
1818       freq_outNMC(2) = ecrit_files(8)
1819       freq_outNMC(3) = ecrit_files(9)
1820       WRITE(lunout,*)'OK freq_outNMC(1)=',freq_outNMC(1)
1821       WRITE(lunout,*)'OK freq_outNMC(2)=',freq_outNMC(2)
1822       WRITE(lunout,*)'OK freq_outNMC(3)=',freq_outNMC(3)
1823
1824#ifndef CPP_XIOS
1825       CALL ini_paramLMDZ_phy(phys_tstep,nid_ctesGCM)
1826#endif
1827
1828#endif
1829       ecrit_reg = ecrit_reg * un_jour
1830       ecrit_tra = ecrit_tra * un_jour
1831
1832       !XXXPB Positionner date0 pour initialisation de ORCHIDEE
1833       date0 = jD_ref
1834       WRITE(*,*) 'physiq date0 : ',date0
1835       !
1836
1837!       CALL create_climoz(read_climoz)
1838      IF (.NOT. create_etat0_limit) CALL init_aero_fromfile(flag_aerosol)  !! initialise aero from file for XIOS interpolation (unstructured_grid)
1839      IF (.NOT. create_etat0_limit) CALL init_readaerosolstrato(flag_aerosol_strat)  !! initialise aero strato from file for XIOS interpolation (unstructured_grid)
1840
1841#ifdef CPP_COSP
1842      IF (ok_cosp) THEN
1843!           DO k = 1, klev
1844!             DO i = 1, klon
1845!               phicosp(i,k) = pphi(i,k) + pphis(i)
1846!             ENDDO
1847!           ENDDO
1848        CALL phys_cosp(itap,phys_tstep,freq_cosp, &
1849               ok_mensuelCOSP,ok_journeCOSP,ok_hfCOSP, &
1850               ecrit_mth,ecrit_day,ecrit_hf, ok_all_xml, missing_val, &
1851               klon,klev,longitude_deg,latitude_deg,presnivs,overlap, &
1852               JrNt,ref_liq,ref_ice, &
1853               pctsrf(:,is_ter)+pctsrf(:,is_lic), &
1854               zu10m,zv10m,pphis, &
1855               zphi,paprs(:,1:klev),pplay,zxtsol,t_seri, &
1856               qx(:,:,ivap),zx_rh,cldfra,rnebcon,flwc,fiwc, &
1857               prfl(:,1:klev),psfl(:,1:klev), &
1858               pmflxr(:,1:klev),pmflxs(:,1:klev), &
1859               mr_ozone,cldtau, cldemi)
1860      ENDIF
1861#endif
1862
1863#ifdef CPP_COSP2
1864        IF (ok_cosp) THEN
1865!           DO k = 1, klev
1866!             DO i = 1, klon
1867!               phicosp(i,k) = pphi(i,k) + pphis(i)
1868!             ENDDO
1869!           ENDDO
1870          CALL phys_cosp2(itap,phys_tstep,freq_cosp, &
1871               ok_mensuelCOSP,ok_journeCOSP,ok_hfCOSP, &
1872               ecrit_mth,ecrit_day,ecrit_hf, ok_all_xml, missing_val, &
1873               klon,klev,longitude_deg,latitude_deg,presnivs,overlap, &
1874               JrNt,ref_liq,ref_ice, &
1875               pctsrf(:,is_ter)+pctsrf(:,is_lic), &
1876               zu10m,zv10m,pphis, &
1877               zphi,paprs(:,1:klev),pplay,zxtsol,t_seri, &
1878               qx(:,:,ivap),zx_rh,cldfra,rnebcon,flwc,fiwc, &
1879               prfl(:,1:klev),psfl(:,1:klev), &
1880               pmflxr(:,1:klev),pmflxs(:,1:klev), &
1881               mr_ozone,cldtau, cldemi)
1882       ENDIF
1883#endif
1884
1885#ifdef CPP_COSPV2
1886        IF (ok_cosp) THEN
1887           DO k = 1, klev
1888             DO i = 1, klon
1889               phicosp(i,k) = pphi(i,k) + pphis(i)
1890             ENDDO
1891           ENDDO
1892          CALL lmdz_cosp_interface(itap,phys_tstep,freq_cosp, &
1893               ok_mensuelCOSP,ok_journeCOSP,ok_hfCOSP, &
1894               ecrit_mth,ecrit_day,ecrit_hf, ok_all_xml, missing_val, &
1895               klon,klev,longitude_deg,latitude_deg,presnivs,overlap, &
1896               JrNt,ref_liq,ref_ice, &
1897               pctsrf(:,is_ter)+pctsrf(:,is_lic), &
1898               zu10m,zv10m,pphis, &
1899               phicosp,paprs(:,1:klev),pplay,zxtsol,t_seri, &
1900               qx(:,:,ivap),zx_rh,cldfra,rnebcon,flwc,fiwc, &
1901               prfl(:,1:klev),psfl(:,1:klev), &
1902               pmflxr(:,1:klev),pmflxs(:,1:klev), &
1903               mr_ozone,cldtau, cldemi)
1904       ENDIF
1905#endif
1906
1907       !
1908       !
1909!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1910       ! Nouvelle initialisation pour le rayonnement RRTM
1911!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1912
1913       CALL iniradia(klon,klev,paprs(1,1:klev+1))
1914
1915!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1916       CALL wake_ini(rg,rd,rv,prt_level)
1917       CALL yamada_ini(klon,lunout,prt_level)
1918       CALL atke_ini(prt_level, lunout, RG, RD, RPI, RCPD)
1919       CALL thermcell_ini(iflag_thermals,prt_level,tau_thermals,lunout, &
1920   &    RG,RD,RCPD,RKAPPA,RLVTT,RETV)
1921       CALL ratqs_ini(klon,klev,iflag_thermals,lunout,nbsrf,is_lic,is_ter,RG,RV,RD,RCPD,RLSTT,RLVTT,RTT)
1922       CALL lscp_ini(pdtphys,ok_ice_sursat,RCPD, RLSTT, RLVTT, RLMLT, RVTMP2, RTT,RD,RG)
1923       CALL blowing_snow_ini(prt_level,lunout, &
1924                             RCPD, RLSTT, RLVTT, RLMLT, &
1925                             RVTMP2, RTT,RD,RG)
1926
1927!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1928
1929       !
1930!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1931       ! Initialisation des champs dans phytrac* qui sont utilises par phys_output_write*
1932       !
1933!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1934
1935#ifdef CPP_Dust
1936       ! Quand on utilise SPLA, on force iflag_phytrac=1
1937       CALL phytracr_spl_out_init()
1938       CALL phys_output_write_spl(itap, pdtphys, paprs, pphis,                  &
1939                                pplay, lmax_th, aerosol_couple,                 &
1940                                ok_ade, ok_aie, ivap, ok_sync,                  &
1941                                ptconv, read_climoz, clevSTD,                   &
1942                                ptconvth, d_t, qx, d_qx, d_tr_dyn, zmasse,      &
1943                                flag_aerosol, flag_aerosol_strat, ok_cdnc)
1944#else
1945       ! phys_output_write écrit des variables traceurs seulement si iflag_phytrac == 1
1946       ! donc seulement dans ce cas on doit appeler phytrac_init()
1947       IF (iflag_phytrac == 1 ) THEN
1948          CALL phytrac_init()
1949       ENDIF
1950       CALL phys_output_write(itap, pdtphys, paprs, pphis,                    &
1951                              pplay, lmax_th, aerosol_couple,                 &
1952                              ok_ade, ok_aie, ok_volcan, ivap, iliq, isol, ibs,  ok_sync,&
1953                              ptconv, read_climoz, clevSTD,                   &
1954                              ptconvth, d_u, d_t, qx, d_qx, zmasse,           &
1955                              flag_aerosol, flag_aerosol_strat, ok_cdnc,t,u1,v1)
1956#endif
1957
1958
1959#ifdef CPP_XIOS
1960       IF (is_omp_master) CALL xios_update_calendar(1)
1961#endif
1962       IF(read_climoz>=1 .AND. create_etat0_limit) CALL regr_horiz_time_climoz(read_climoz,ok_daily_climoz)
1963       CALL create_etat0_limit_unstruct
1964       CALL phyetat0 ("startphy.nc",clesphy0,tabcntr0)
1965
1966!jyg<
1967       IF (iflag_pbl<=1) THEN
1968          ! No TKE for Standard Physics
1969          pbl_tke(:,:,:)=0.
1970
1971       ELSE IF (klon_glo==1) THEN
1972          pbl_tke(:,:,is_ave) = 0.
1973          DO nsrf=1,nbsrf
1974            DO k = 1,klev+1
1975                 pbl_tke(:,k,is_ave) = pbl_tke(:,k,is_ave) &
1976                     +pctsrf(:,nsrf)*pbl_tke(:,k,nsrf)
1977            ENDDO
1978          ENDDO
1979       ELSE
1980          pbl_tke(:,:,is_ave) = 0. !ym missing init : maybe must be initialized in the same way that for klon_glo==1 ??
1981!>jyg
1982       ENDIF
1983       !IM begin
1984       print*,'physiq: clwcon rnebcon ratqs',clwcon(1,1),rnebcon(1,1) &
1985            ,ratqs(1,1)
1986       !IM end
1987
1988#ifdef ISO
1989#ifdef ISOTRAC
1990! on initialise les cartes des bassins pour les traceurs si besoin:
1991       call initialise_bassins_boites(presnivs)
1992#endif
1993#endif
1994
1995       ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1996       !
1997       ! on remet le calendrier a zero
1998       !
1999       IF (raz_date .eq. 1) THEN
2000          itau_phy = 0
2001       ENDIF
2002
2003!       IF (ABS(phys_tstep-pdtphys).GT.0.001) THEN
2004!          WRITE(lunout,*) 'Pas physique n est pas correct',phys_tstep, &
2005!               pdtphys
2006!          abort_message='Pas physique n est pas correct '
2007!          !           call abort_physic(modname,abort_message,1)
2008!          phys_tstep=pdtphys
2009!       ENDIF
2010       IF (nlon .NE. klon) THEN
2011          WRITE(lunout,*)'nlon et klon ne sont pas coherents', nlon,  &
2012               klon
2013          abort_message='nlon et klon ne sont pas coherents'
2014          CALL abort_physic(modname,abort_message,1)
2015       ENDIF
2016       IF (nlev .NE. klev) THEN
2017          WRITE(lunout,*)'nlev et klev ne sont pas coherents', nlev, &
2018               klev
2019          abort_message='nlev et klev ne sont pas coherents'
2020          CALL abort_physic(modname,abort_message,1)
2021       ENDIF
2022       !
2023       IF (phys_tstep*REAL(radpas).GT.21600..AND.iflag_cycle_diurne.GE.1) THEN
2024          WRITE(lunout,*)'Nbre d appels au rayonnement insuffisant'
2025          WRITE(lunout,*)"Au minimum 4 appels par jour si cycle diurne"
2026          abort_message='Nbre d appels au rayonnement insuffisant'
2027          CALL abort_physic(modname,abort_message,1)
2028       ENDIF
2029
2030!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2031       ! Initialisation pour la convection de K.E. et pour les poches froides
2032       !
2033!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2034
2035       WRITE(lunout,*)"Clef pour la convection, iflag_con=", iflag_con
2036       WRITE(lunout,*)"Clef pour le driver de la convection, ok_cvl=", ok_cvl
2037       !
2038       !KE43
2039       ! Initialisation pour la convection de K.E. (sb):
2040       IF (iflag_con.GE.3) THEN
2041
2042          WRITE(lunout,*)"*** Convection de Kerry Emanuel 4.3  "
2043          WRITE(lunout,*) &
2044               "On va utiliser le melange convectif des traceurs qui"
2045          WRITE(lunout,*)"est calcule dans convect4.3"
2046          WRITE(lunout,*)" !!! penser aux logical flags de phytrac"
2047
2048          DO i = 1, klon
2049             ema_cbmf(i) = 0.
2050             ema_pcb(i)  = 0.
2051             ema_pct(i)  = 0.
2052             !          ema_workcbmf(i) = 0.
2053          ENDDO
2054          !IM15/11/02 rajout initialisation ibas_con,itop_con cf. SB =>BEG
2055          DO i = 1, klon
2056             ibas_con(i) = 1
2057             itop_con(i) = 1
2058          ENDDO
2059          !IM15/11/02 rajout initialisation ibas_con,itop_con cf. SB =>END
2060          !================================================================
2061          !CR:04.12.07: initialisations poches froides
2062          ! Controle de ALE et ALP pour la fermeture convective (jyg)
2063          IF (iflag_wake>=1) THEN
2064             CALL ini_wake(0.,0.,it_wape_prescr,wape_prescr,fip_prescr &
2065                  ,alp_bl_prescr, ale_bl_prescr)
2066             ! 11/09/06 rajout initialisation ALE et ALP du wake et PBL(YU)
2067             !        print*,'apres ini_wake iflag_cld_th=', iflag_cld_th
2068             !
2069             ! Initialize tendencies of wake state variables (for some flag values
2070             ! they are not computed).
2071             d_deltat_wk(:,:) = 0.
2072             d_deltaq_wk(:,:) = 0.
2073             d_deltat_wk_gw(:,:) = 0.
2074             d_deltaq_wk_gw(:,:) = 0.
2075             d_deltat_vdf(:,:) = 0.
2076             d_deltaq_vdf(:,:) = 0.
2077             d_deltat_the(:,:) = 0.
2078             d_deltaq_the(:,:) = 0.
2079             d_deltat_ajs_cv(:,:) = 0.
2080             d_deltaq_ajs_cv(:,:) = 0.
2081             d_s_wk(:) = 0.
2082             d_dens_wk(:) = 0.
2083#ifdef ISO
2084             d_deltaxt_wk(:,:,:) = 0.
2085             d_deltaxt_vdf(:,:,:) = 0.
2086             d_deltaxt_the(:,:,:) = 0.
2087             d_deltaxt_ajs_cv(:,:,:) = 0.
2088#endif
2089          ENDIF  !  (iflag_wake>=1)
2090
2091          !        do i = 1,klon
2092          !           Ale_bl(i)=0.
2093          !           Alp_bl(i)=0.
2094          !        enddo
2095
2096       !ELSE
2097       !   ALLOCATE(tabijGCM(0))
2098       !   ALLOCATE(lonGCM(0), latGCM(0))
2099       !   ALLOCATE(iGCM(0), jGCM(0))
2100       ENDIF  !  (iflag_con.GE.3)
2101       !
2102       DO i=1,klon
2103          rugoro(i) = f_rugoro * MAX(1.0e-05, zstd(i)*zsig(i)/2.0)
2104       ENDDO
2105
2106       !34EK
2107       IF (ok_orodr) THEN
2108
2109          ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2110          ! FH sans doute a enlever de finitivement ou, si on le
2111          ! garde, l'activer justement quand ok_orodr = false.
2112          ! ce rugoro est utilise par la couche limite et fait double emploi
2113          ! avec les param\'etrisations sp\'ecifiques de Francois Lott.
2114          !           DO i=1,klon
2115          !             rugoro(i) = MAX(1.0e-05, zstd(i)*zsig(i)/2.0)
2116          !           ENDDO
2117          ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2118          IF (ok_strato) THEN
2119             CALL SUGWD_strato(klon,klev,paprs,pplay)
2120          ELSE
2121             CALL SUGWD(klon,klev,paprs,pplay)
2122          ENDIF
2123
2124          DO i=1,klon
2125             zuthe(i)=0.
2126             zvthe(i)=0.
2127             IF (zstd(i).gt.10.) THEN
2128                zuthe(i)=(1.-zgam(i))*cos(zthe(i))
2129                zvthe(i)=(1.-zgam(i))*sin(zthe(i))
2130             ENDIF
2131          ENDDO
2132       ENDIF
2133       !
2134       !
2135       lmt_pas = NINT(86400./phys_tstep * 1.0)   ! tous les jours
2136       WRITE(lunout,*)'La frequence de lecture surface est de ',  &
2137            lmt_pas
2138       !
2139       capemaxcels = 't_max(X)'
2140       t2mincels = 't_min(X)'
2141       t2maxcels = 't_max(X)'
2142       tinst = 'inst(X)'
2143       tave = 'ave(X)'
2144       !IM cf. AM 081204 BEG
2145       write(lunout,*)'AVANT HIST IFLAG_CON=',iflag_con
2146       !IM cf. AM 081204 END
2147       !
2148       !=============================================================
2149       !   Initialisation des sorties
2150       !=============================================================
2151
2152#ifdef CPP_XIOS
2153       ! Get "missing_val" value from XML files (from temperature variable)
2154       !$OMP MASTER
2155       CALL xios_get_field_attr("temp",default_value=missing_val_omp)
2156       !$OMP END MASTER
2157       !$OMP BARRIER
2158       missing_val=missing_val_omp
2159#endif
2160
2161#ifdef CPP_XIOS
2162! Need to put this initialisation after phyetat0 as in the coupled model the XIOS context is only
2163! initialised at that moment
2164       ! Get "missing_val" value from XML files (from temperature variable)
2165       !$OMP MASTER
2166       CALL xios_get_field_attr("temp",default_value=missing_val_omp)
2167       !$OMP END MASTER
2168       !$OMP BARRIER
2169       missing_val=missing_val_omp
2170       !
2171       ! Now we activate some double radiation call flags only if some
2172       ! diagnostics are requested, otherwise there is no point in doing this
2173       IF (is_master) THEN
2174         !--setting up swaero_diag to TRUE in XIOS case
2175         IF (xios_field_is_active("topswad").OR.xios_field_is_active("topswad0").OR. &
2176            xios_field_is_active("solswad").OR.xios_field_is_active("solswad0").OR. &
2177            xios_field_is_active("topswai").OR.xios_field_is_active("solswai").OR.  &
2178              (iflag_rrtm==1.AND.(xios_field_is_active("toplwad").OR.xios_field_is_active("toplwad0").OR. &
2179                                  xios_field_is_active("sollwad").OR.xios_field_is_active("sollwad0"))))  &
2180            !!!--for now these fields are not in the XML files so they are omitted
2181            !!!  xios_field_is_active("toplwai").OR.xios_field_is_active("sollwai") !))) &
2182            swaero_diag=.TRUE.
2183 
2184         !--setting up swaerofree_diag to TRUE in XIOS case
2185         IF (xios_field_is_active("SWdnSFCcleanclr").OR.xios_field_is_active("SWupSFCcleanclr").OR. &
2186            xios_field_is_active("SWupTOAcleanclr").OR.xios_field_is_active("rsucsaf").OR.   &
2187            xios_field_is_active("rsdcsaf") .OR. xios_field_is_active("LWdnSFCcleanclr").OR. &
2188            xios_field_is_active("LWupTOAcleanclr")) &
2189            swaerofree_diag=.TRUE.
2190 
2191         !--setting up dryaod_diag to TRUE in XIOS case
2192         DO naero = 1, naero_tot-1
2193          IF (xios_field_is_active("dryod550_"//name_aero_tau(naero))) dryaod_diag=.TRUE.
2194         ENDDO
2195         !
2196         !--setting up ok_4xCO2atm to TRUE in XIOS case
2197         IF (xios_field_is_active("rsut4co2").OR.xios_field_is_active("rlut4co2").OR. &
2198            xios_field_is_active("rsutcs4co2").OR.xios_field_is_active("rlutcs4co2").OR. &
2199            xios_field_is_active("rsu4co2").OR.xios_field_is_active("rsucs4co2").OR. &
2200            xios_field_is_active("rsd4co2").OR.xios_field_is_active("rsdcs4co2").OR. &
2201            xios_field_is_active("rlu4co2").OR.xios_field_is_active("rlucs4co2").OR. &
2202            xios_field_is_active("rld4co2").OR.xios_field_is_active("rldcs4co2")) &
2203            ok_4xCO2atm=.TRUE.
2204       ENDIF
2205       !$OMP BARRIER
2206       CALL bcast(swaero_diag)
2207       CALL bcast(swaerofree_diag)
2208       CALL bcast(dryaod_diag)
2209       CALL bcast(ok_4xCO2atm)
2210#endif
2211       !
2212       CALL printflag( tabcntr0,radpas,ok_journe, &
2213            ok_instan, ok_region )
2214       !
2215       !
2216       ! Prescrire l'ozone dans l'atmosphere
2217       !
2218       !c         DO i = 1, klon
2219       !c         DO k = 1, klev
2220       !c            CALL o3cm (paprs(i,k)/100.,paprs(i,k+1)/100., wo(i,k),20)
2221       !c         ENDDO
2222       !c         ENDDO
2223       !
2224       IF (ANY(type_trac == ['inca','inco'])) THEN ! ModThL
2225#ifdef INCA
2226          CALL VTe(VTphysiq)
2227          CALL VTb(VTinca)
2228          calday = REAL(days_elapsed) + jH_cur
2229          WRITE(lunout,*) 'initial time chemini', days_elapsed, calday
2230
2231          call init_const_lmdz( &
2232          ndays, nbsrf, is_oce,is_sic, is_ter,is_lic, calend, &
2233          config_inca)
2234
2235          CALL init_inca_geometry( &
2236               longitude, latitude, &
2237               boundslon, boundslat, &
2238               cell_area, ind_cell_glo)
2239
2240          if (grid_type==unstructured) THEN
2241             CALL chemini(  pplay, &
2242                  nbp_lon, nbp_lat, &
2243                  latitude_deg, &
2244                  longitude_deg, &
2245                  presnivs, &
2246                  calday, &
2247                  klon, &
2248                  nqtot, &
2249                  nqo+nqCO2, &
2250                  pdtphys, &
2251                  annee_ref, &
2252                  year_cur, &
2253                  day_ref,  &
2254                  day_ini, &
2255                  start_time, &
2256                  itau_phy, &
2257                  date0, &
2258                  chemistry_couple, &
2259                  init_source, &
2260                  init_tauinca, &
2261                  init_pizinca, &
2262                  init_cginca, &
2263                  init_ccminca)
2264          ELSE
2265             CALL chemini(  pplay, &
2266                  nbp_lon, nbp_lat, &
2267                  latitude_deg, &
2268                  longitude_deg, &
2269                  presnivs, &
2270                  calday, &
2271                  klon, &
2272                  nqtot, &
2273                  nqo+nqCO2, &
2274                  pdtphys, &
2275                  annee_ref, &
2276                  year_cur, &
2277                  day_ref,  &
2278                  day_ini, &
2279                  start_time, &
2280                  itau_phy, &
2281                  date0, &
2282                  chemistry_couple, &
2283                  init_source, &
2284                  init_tauinca, &
2285                  init_pizinca, &
2286                  init_cginca, &
2287                  init_ccminca, &
2288                  io_lon, &
2289                  io_lat)
2290          ENDIF
2291
2292
2293          ! initialisation des variables depuis le restart de inca
2294          ccm(:,:,:) = init_ccminca
2295          tau_aero(:,:,:,:) = init_tauinca
2296          piz_aero(:,:,:,:) = init_pizinca
2297          cg_aero(:,:,:,:) = init_cginca
2298!         
2299
2300
2301          CALL VTe(VTinca)
2302          CALL VTb(VTphysiq)
2303#endif
2304       ENDIF
2305       !
2306       IF (type_trac == 'repr') THEN
2307#ifdef REPROBUS
2308          CALL chemini_rep(  &
2309               presnivs, &
2310               pdtphys, &
2311               annee_ref, &
2312               day_ref,  &
2313               day_ini, &
2314               start_time, &
2315               itau_phy, &
2316               io_lon, &
2317               io_lat)
2318#endif
2319       ENDIF
2320
2321       !$omp single
2322       IF (read_climoz >= 1) CALL open_climoz(ncid_climoz, press_cen_climoz,   &
2323           press_edg_climoz, time_climoz, ok_daily_climoz, adjust_tropopause)
2324       !$omp end single
2325       !
2326       !IM betaCRF
2327       pfree=70000. !Pa
2328       beta_pbl=1.
2329       beta_free=1.
2330       lon1_beta=-180.
2331       lon2_beta=+180.
2332       lat1_beta=90.
2333       lat2_beta=-90.
2334       mskocean_beta=.FALSE.
2335
2336       !albedo SB >>>
2337       SELECT CASE(nsw)
2338       CASE(2)
2339          SFRWL(1)=0.45538747
2340          SFRWL(2)=0.54461211
2341       CASE(4)
2342          SFRWL(1)=0.45538747
2343          SFRWL(2)=0.32870591
2344          SFRWL(3)=0.18568763
2345          SFRWL(4)=3.02191470E-02
2346       CASE(6)
2347          SFRWL(1)=1.28432794E-03
2348          SFRWL(2)=0.12304168
2349          SFRWL(3)=0.33106142
2350          SFRWL(4)=0.32870591
2351          SFRWL(5)=0.18568763
2352          SFRWL(6)=3.02191470E-02
2353       END SELECT
2354       !albedo SB <<<
2355
2356       OPEN(99,file='beta_crf.data',status='old', &
2357            form='formatted',err=9999)
2358       READ(99,*,end=9998) pfree
2359       READ(99,*,end=9998) beta_pbl
2360       READ(99,*,end=9998) beta_free
2361       READ(99,*,end=9998) lon1_beta
2362       READ(99,*,end=9998) lon2_beta
2363       READ(99,*,end=9998) lat1_beta
2364       READ(99,*,end=9998) lat2_beta
2365       READ(99,*,end=9998) mskocean_beta
23669998   Continue
2367       CLOSE(99)
23689999   Continue
2369       WRITE(*,*)'pfree=',pfree
2370       WRITE(*,*)'beta_pbl=',beta_pbl
2371       WRITE(*,*)'beta_free=',beta_free
2372       WRITE(*,*)'lon1_beta=',lon1_beta
2373       WRITE(*,*)'lon2_beta=',lon2_beta
2374       WRITE(*,*)'lat1_beta=',lat1_beta
2375       WRITE(*,*)'lat2_beta=',lat2_beta
2376       WRITE(*,*)'mskocean_beta=',mskocean_beta
2377
2378      !lwoff=y : offset LW CRE for radiation code and other schemes
2379      !lwoff=y : betalwoff=1.
2380      betalwoff=0.
2381      IF (ok_lwoff) THEN
2382         betalwoff=1.
2383      ENDIF
2384      WRITE(*,*)'ok_lwoff=',ok_lwoff
2385      !
2386      !lwoff=y to begin only sollw and sollwdown are set up to CS values
2387      sollw = sollw + betalwoff * (sollw0 - sollw)
2388      sollwdown(:)= sollwdown(:) + betalwoff *(-1.*ZFLDN0(:,1) - &
2389                    sollwdown(:))
2390
2391
2392
2393    ENDIF
2394    !
2395    !   ****************     Fin  de   IF ( debut  )   ***************
2396    !
2397    !
2398    ! Incrementer le compteur de la physique
2399    !
2400    itap   = itap + 1
2401    IF (is_master .OR. prt_level > 9) THEN
2402      IF (prt_level > 5 .or. MOD(itap,5) == 0) THEN
2403         WRITE(LUNOUT,*)'Entering physics elapsed seconds since start ', current_time
2404         WRITE(LUNOUT,100)year_cur,mth_cur,day_cur,hour/3600.
2405 100     FORMAT('Date = ',i4.4,' / ',i2.2, ' / ',i2.2,' : ',f20.17)
2406      ENDIF
2407    ENDIF
2408    !
2409    !
2410    ! Update fraction of the sub-surfaces (pctsrf) and
2411    ! initialize, where a new fraction has appeared, all variables depending
2412    ! on the surface fraction.
2413    !
2414    CALL change_srf_frac(itap, phys_tstep, days_elapsed+1,  &
2415         pctsrf, fevap, z0m, z0h, agesno,              &
2416         falb_dir, falb_dif, ftsol, ustar, u10m, v10m, pbl_tke &
2417#ifdef ISO
2418     ,fxtevap  &
2419#endif
2420&       )
2421
2422    ! Update time and other variables in Reprobus
2423    IF (type_trac == 'repr') THEN
2424#ifdef REPROBUS
2425       CALL Init_chem_rep_xjour(jD_cur-jD_ref+day_ref)
2426       print*,'xjour equivalent rjourvrai',jD_cur-jD_ref+day_ref
2427       CALL Rtime(debut)
2428#endif
2429    ENDIF
2430
2431    ! Tendances bidons pour les processus qui n'affectent pas certaines
2432    ! variables.
2433    du0(:,:)=0.
2434    dv0(:,:)=0.
2435    dt0 = 0.
2436    dq0(:,:)=0.
2437    dql0(:,:)=0.
2438    dqi0(:,:)=0.
2439    dqbs0(:,:)=0.
2440#ifdef ISO
2441      dxt0(:,:,:)=0.
2442      dxtl0(:,:,:)=0.
2443      dxti0(:,:,:)=0.
2444#endif
2445    dsig0(:) = 0.
2446    ddens0(:) = 0.
2447    wkoccur1(:)=1
2448    !
2449    ! Mettre a zero des variables de sortie (pour securite)
2450    !
2451    DO i = 1, klon
2452       d_ps(i) = 0.0
2453    ENDDO
2454    DO k = 1, klev
2455       DO i = 1, klon
2456          d_t(i,k) = 0.0
2457          d_u(i,k) = 0.0
2458          d_v(i,k) = 0.0
2459       ENDDO
2460    ENDDO
2461    DO iq = 1, nqtot
2462       DO k = 1, klev
2463          DO i = 1, klon
2464             d_qx(i,k,iq) = 0.0
2465          ENDDO
2466       ENDDO
2467    ENDDO
2468    beta_prec_fisrt(:,:)=0.
2469    beta_prec(:,:)=0.
2470    !
2471    !   Output variables from the convective scheme should not be set to 0
2472    !   since convection is not always called at every time step.
2473    IF (ok_bug_cv_trac) THEN
2474      da(:,:)=0.
2475      mp(:,:)=0.
2476      phi(:,:,:)=0.
2477      ! RomP >>>
2478      phi2(:,:,:)=0.
2479      epmlmMm(:,:,:)=0.
2480      eplaMm(:,:)=0.
2481      d1a(:,:)=0.
2482      dam(:,:)=0.
2483      pmflxr(:,:)=0.
2484      pmflxs(:,:)=0.
2485      ! RomP <<<
2486    ENDIF
2487    !
2488    ! Ne pas affecter les valeurs entrees de u, v, h, et q
2489    !
2490    DO k = 1, klev
2491       DO i = 1, klon
2492          t_seri(i,k)  = t(i,k)
2493          u_seri(i,k)  = u(i,k)
2494          v_seri(i,k)  = v(i,k)
2495          q_seri(i,k)  = qx(i,k,ivap)
2496          ql_seri(i,k) = qx(i,k,iliq)
2497          qbs_seri(i,k) = 0.
2498          !CR: ATTENTION, on rajoute la variable glace
2499          IF (nqo.EQ.2) THEN             !--vapour and liquid only
2500             qs_seri(i,k) = 0.
2501             rneb_seri(i,k) = 0.
2502          ELSE IF (nqo.EQ.3) THEN        !--vapour, liquid and ice
2503             qs_seri(i,k) = qx(i,k,isol)
2504             rneb_seri(i,k) = 0.
2505          ELSE IF (nqo.GE.4) THEN        !--vapour, liquid, ice and rneb and blowing snow
2506             qs_seri(i,k) = qx(i,k,isol)
2507             IF (ok_ice_sursat) THEN
2508             rneb_seri(i,k) = qx(i,k,irneb)
2509             ENDIF
2510             IF (ok_bs) THEN
2511             qbs_seri(i,k)= qx(i,k,ibs)
2512             ENDIF
2513
2514          ENDIF
2515
2516
2517       ENDDO
2518    ENDDO
2519
2520    ! C Risi: dispatcher les isotopes dans les xt_seri
2521#ifdef ISO
2522#ifdef ISOVERIF
2523    write(*,*) 'physiq 1847: qx(1,1,:)=',qx(1,1,:)
2524    write(*,*) 'iqIsoPha(:,ivap)=',iqIsoPha(:,ivap)
2525    write(*,*) 'physiq 1846b: ok_isotopes,ntraciso,niso=',niso>0,ntraciso,niso
2526#endif
2527    do ixt=1,ntraciso
2528#ifdef ISOVERIF
2529      write(*,*) 'physiq tmp 1762a: ixt,iqiso_vap=',ixt,iqIsoPha(ixt,ivap)
2530      write(*,*) 'physiq tmp 1762b: ixt,iqiso_liq=',ixt,iqIsoPha(ixt,iliq)
2531      if (nqo.eq.3) then 
2532        write(*,*) 'physiq tmp 1762c: ixt,iqiso_liq=',ixt,iqIsoPha(ixt,iliq)
2533      endif !if (nqo.eq.3) then
2534#endif
2535      if (ixt.gt.niso) write(*,*) 'izone=',tracers(iqIsoPha(ixt,ivap))%iso_iZone
2536      DO k = 1, klev
2537       DO i = 1, klon
2538          xt_seri(ixt,i,k)  = qx(i,k,iqIsoPha(ixt,ivap))
2539          xtl_seri(ixt,i,k) = qx(i,k,iqIsoPha(ixt,iliq))
2540          if (nqo.eq.2) then
2541             xts_seri(ixt,i,k) = 0.
2542          else if (nqo.eq.3) then
2543             xts_seri(ixt,i,k) = qx(i,k,iqIsoPha(ixt,isol))
2544          endif
2545       enddo !DO i = 1, klon
2546      enddo ! DO k = 1, klev
2547      !write(*,*) 'physiq tmp 1897: dispatch termine pour ixt=',ixt
2548    enddo !do ixt=1,niso
2549!    write(*,*) 'physiq tmp 1898: dispatch termine pour tous'
2550#endif
2551! #ifdef ISO
2552    !
2553    !--OB mass fixer
2554    IF (mass_fixer) THEN
2555    !--store initial water burden
2556    qql1(:)=0.0
2557    DO k = 1, klev
2558      qql1(:)=qql1(:)+(q_seri(:,k)+ql_seri(:,k)+qs_seri(:,k)+qbs_seri(:,k))*zmasse(:,k)
2559    ENDDO
2560#ifdef ISO
2561#ifdef ISOVERIF
2562        write(*,*) 'physiq tmp 1913'
2563#endif
2564    do ixt=1,ntraciso
2565    xtql1(ixt,:)=0.0
2566    DO k = 1, klev
2567      xtql1(ixt,:)=xtql1(ixt,:)+(xt_seri(ixt,:,k)+xtl_seri(ixt,:,k)+xts_seri(ixt,:,k))*zmasse(:,k)
2568    ENDDO
2569    enddo !do ixt=1,ntraciso
2570#endif
2571    ENDIF
2572    !--fin mass fixer
2573
2574    tke0(:,:)=pbl_tke(:,:,is_ave)
2575    IF (nqtot > nqo) THEN
2576       ! water isotopes are not included in tr_seri
2577       itr = 0
2578       DO iq = 1, nqtot
2579         IF(.NOT.tracers(iq)%isInPhysics) CYCLE
2580         itr = itr+1
2581!#ifdef ISOVERIF
2582!         write(*,*) 'physiq 1973: itr,iq=',itr,iq
2583!         write(*,*) 'qx(1,1,iq)=',qx(1,1,iq)
2584!#endif
2585         DO  k = 1, klev
2586             DO  i = 1, klon
2587                tr_seri(i,k,itr) = qx(i,k,iq)
2588             ENDDO
2589          ENDDO
2590       ENDDO
2591    ELSE
2592! DC: make sure the final "1" index was meant for 1st H2O phase (vapor) !!!
2593       tr_seri(:,:,strIdx(tracers(:)%name,addPhase('H2O','g'))) = 0.0
2594    ENDIF
2595!
2596! Temporary solutions adressing ticket #104 and the non initialisation of tr_ancien
2597! LF
2598    IF (debut) THEN
2599      WRITE(lunout,*)' WARNING: tr_ancien initialised to tr_seri'
2600       itr = 0
2601       do iq = 1, nqtot
2602         IF(.NOT.tracers(iq)%isInPhysics) CYCLE
2603         itr = itr+1
2604         tr_ancien(:,:,itr)=tr_seri(:,:,itr)       
2605       enddo
2606    ENDIF
2607    !
2608    DO i = 1, klon
2609       ztsol(i) = 0.
2610    ENDDO
2611    DO nsrf = 1, nbsrf
2612       DO i = 1, klon
2613          ztsol(i) = ztsol(i) + ftsol(i,nsrf)*pctsrf(i,nsrf)
2614       ENDDO
2615    ENDDO
2616    ! Initialize variables used for diagnostic purpose
2617    IF (flag_inhib_tend .ne. 0) CALL init_cmp_seri
2618
2619    ! Diagnostiquer la tendance dynamique
2620#ifdef ISOVERIF
2621      if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then 
2622        do i=1,klon
2623         do k=1,klev
2624            if (q_seri(i,k).gt.ridicule) then 
2625               if (iso_verif_o18_aberrant_nostop( &
2626     &              xt_seri(iso_HDO,i,k)/q_seri(i,k), &
2627     &              xt_seri(iso_O18,i,k)/q_seri(i,k), &
2628     &              'physiq 2099 qv').eq.1) then
2629                  write(*,*) 'physic 2444: i,k,q_seri(i,k)=',i,k,q_seri(i,k)
2630                  write(*,*) 'xt_seri(:,i,k)=',xt_seri(:,i,k)
2631                  stop
2632              endif !  if (iso_verif_o18_aberrant_nostop
2633            endif !if (q_seri(i,k).gt.errmax) then 
2634            if (ql_seri(i,k).gt.ridicule) then 
2635               if (iso_verif_o18_aberrant_nostop( &
2636     &              xtl_seri(iso_HDO,i,k)/ql_seri(i,k), &
2637     &              xtl_seri(iso_O18,i,k)/ql_seri(i,k), &
2638     &              'physiq 2099 ql').eq.1) then
2639                  write(*,*) 'i,k,ql_seri(i,k)=',i,k,ql_seri(i,k)
2640                  !stop
2641              endif !  if (iso_verif_o18_aberrant_nostop
2642            endif !if (q_seri(i,k).gt.errmax) then 
2643            if (qs_seri(i,k).gt.ridicule) then 
2644               if (iso_verif_o18_aberrant_nostop( &
2645     &              xts_seri(iso_HDO,i,k)/qs_seri(i,k), &
2646     &              xts_seri(iso_O18,i,k)/qs_seri(i,k), &
2647     &              'physiq 2099 qs').eq.1) then
2648                  write(*,*) 'i,k,qs_seri(i,k)=',i,k,qs_seri(i,k)
2649                  !stop
2650              endif !  if (iso_verif_o18_aberrant_nostop
2651            endif !if (q_seri(i,k).gt.errmax) then
2652          enddo !k=1,klev
2653         enddo  !i=1,klon
2654        endif !if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then         
2655#ifdef ISOTRAC 
2656         DO k = 1, klev
2657           DO i = 1, klon   
2658             call iso_verif_traceur(xt_seri(1,i,k),'physiq 2620a')
2659             call iso_verif_traceur(xtl_seri(1,i,k),'physiq 2620b')
2660             call iso_verif_traceur(xts_seri(1,i,k),'physiq 2620c')
2661           enddo
2662         enddo
2663#endif
2664#endif
2665    !
2666    IF (ancien_ok) THEN
2667    !
2668       d_u_dyn(:,:)  = (u_seri(:,:)-u_ancien(:,:))/phys_tstep
2669       d_v_dyn(:,:)  = (v_seri(:,:)-v_ancien(:,:))/phys_tstep
2670       d_t_dyn(:,:)  = (t_seri(:,:)-t_ancien(:,:))/phys_tstep
2671       d_q_dyn(:,:)  = (q_seri(:,:)-q_ancien(:,:))/phys_tstep
2672       d_ql_dyn(:,:) = (ql_seri(:,:)-ql_ancien(:,:))/phys_tstep
2673       d_qs_dyn(:,:) = (qs_seri(:,:)-qs_ancien(:,:))/phys_tstep
2674       d_qbs_dyn(:,:) = (qbs_seri(:,:)-qbs_ancien(:,:))/phys_tstep
2675       CALL water_int(klon,klev,q_seri,zmasse,zx_tmp_fi2d)
2676       d_q_dyn2d(:)=(zx_tmp_fi2d(:)-prw_ancien(:))/phys_tstep
2677       CALL water_int(klon,klev,ql_seri,zmasse,zx_tmp_fi2d)
2678       d_ql_dyn2d(:)=(zx_tmp_fi2d(:)-prlw_ancien(:))/phys_tstep
2679       CALL water_int(klon,klev,qs_seri,zmasse,zx_tmp_fi2d)
2680       d_qs_dyn2d(:)=(zx_tmp_fi2d(:)-prsw_ancien(:))/phys_tstep
2681       CALL water_int(klon,klev,qbs_seri,zmasse,zx_tmp_fi2d)
2682       d_qbs_dyn2d(:)=(zx_tmp_fi2d(:)-prbsw_ancien(:))/phys_tstep
2683       ! !! RomP >>>   td dyn traceur
2684       IF (nqtot > nqo) d_tr_dyn(:,:,:)=(tr_seri(:,:,:)-tr_ancien(:,:,:))/phys_tstep
2685       ! !! RomP <<<
2686       !!d_rneb_dyn(:,:)=(rneb_seri(:,:)-rneb_ancien(:,:))/phys_tstep
2687       d_rneb_dyn(:,:)=0.0
2688
2689#ifdef ISO
2690!       write(*,*) 'physiq tmp 2112'
2691         DO k = 1, klev
2692         DO i = 1, klon
2693            do ixt=1,ntraciso
2694              d_xt_dyn(ixt,i,k) =  &
2695     &           (xt_seri(ixt,i,k)-xt_ancien(ixt,i,k))/phys_tstep
2696              d_xtl_dyn(ixt,i,k) =  &
2697     &           (xtl_seri(ixt,i,k)-xtl_ancien(ixt,i,k))/phys_tstep
2698              d_xts_dyn(ixt,i,k) =  &
2699     &           (xts_seri(ixt,i,k)-xts_ancien(ixt,i,k))/phys_tstep
2700            enddo ! do ixt=1,ntraciso
2701         ENDDO
2702         ENDDO   
2703#ifdef ISOVERIF
2704         DO k = 1, klev
2705         DO i = 1, klon
2706           do ixt=1,ntraciso
2707           call iso_verif_noNaN(xt_seri(ixt,i,k),'physiq 2220')
2708           call iso_verif_noNaN(xtl_seri(ixt,i,k),'physiq 2220b')
2709           call iso_verif_noNaN(xts_seri(ixt,i,k),'physiq 2220b')
2710           call iso_verif_noNaN(d_xt_dyn(ixt,i,k),'physiq 2220c')
2711           call iso_verif_noNaN(d_xtl_dyn(ixt,i,k),'physiq 2220d')
2712           call iso_verif_noNaN(d_xts_dyn(ixt,i,k),'physiq 2220e')
2713           enddo ! do ixt=1,ntraciso
2714         enddo
2715         enddo
2716#endif   
2717#ifdef ISOVERIF     
2718         if (iso_eau.gt.0) then
2719         DO k = 1, klev
2720         DO i = 1, klon   
2721           call iso_verif_egalite_choix(xt_seri(iso_eau,i,k), &
2722     &            q_seri(i,k),'physiq 2033a',errmax,errmaxrel)
2723           call iso_verif_egalite_choix(xtl_seri(iso_eau,i,k), &
2724     &            ql_seri(i,k),'physiq 2033al',errmax,errmaxrel)
2725           call iso_verif_egalite_choix(xts_seri(iso_eau,i,k), &
2726     &            qs_seri(i,k),'physiq 2033as',errmax,errmaxrel)
2727           call iso_verif_egalite_choix(xt_ancien(iso_eau,i,k), &
2728     &            q_ancien(i,k),'physiq 2033b',errmax,errmaxrel)
2729           call iso_verif_egalite_choix(d_xt_dyn(iso_eau,i,k), &
2730     &            d_q_dyn(i,k),'physiq 2033c',errmax,errmaxrel)
2731         enddo
2732         enddo   
2733         endif
2734         if (iso_HDO.gt.0) then
2735            DO k = 1, klev
2736             DO i = 1, klon
2737              if (q_seri(i,k).gt.3e-3) then
2738              call iso_verif_positif(deltaD(xt_seri(iso_HDO,i,k) &
2739     &           /q_seri(i,k))+400.0,'physiq 2045a')
2740              call iso_verif_positif(deltaD(xt_ancien(iso_HDO,i,k) &
2741     &            /q_ancien(i,k))+400.0,'physiq 2045b')
2742!              call iso_verif_egalite_choix(d_xt_dyn(iso_hdo,i,k) &
2743!     &          /tnat(iso_hdo),d_q_dyn(i,k),'physiq 2045c',1e-7,0.4)
2744              endif
2745             enddo
2746            enddo
2747         endif           
2748#ifdef ISOTRAC 
2749         DO k = 1, klev
2750           DO i = 1, klon   
2751             call iso_verif_traceur(xt_seri(1,i,k),'physiq 2065')
2752             call iso_verif_traceur(xt_ancien(1,i,k),'physiq 2066')
2753           enddo
2754         enddo
2755#endif           
2756#endif         
2757#endif
2758
2759    ELSE ! ancien_OK
2760       d_u_dyn(:,:)  = 0.0
2761       d_v_dyn(:,:)  = 0.0
2762       d_t_dyn(:,:)  = 0.0
2763       d_q_dyn(:,:)  = 0.0
2764       d_ql_dyn(:,:) = 0.0
2765       d_qs_dyn(:,:) = 0.0
2766       d_q_dyn2d(:)  = 0.0
2767       d_ql_dyn2d(:) = 0.0
2768       d_qs_dyn2d(:) = 0.0
2769       d_qbs_dyn2d(:)= 0.0
2770
2771#ifdef ISO
2772!        write(*,*) 'physiq tmp 2115'
2773!        i=496
2774!        k=18
2775!        write(*,*) 'physiq 2105: q_seri(i,k),ql_seri(i,k),qs_seri(i,k),ridicule=', &
2776!                q_seri(i,k),ql_seri(i,k),qs_seri(i,k),ridicule
2777         DO k = 1, klev
2778          DO i = 1, klon
2779           do ixt=1,ntraciso
2780            d_xt_dyn(ixt,i,k)= 0.0
2781            d_xtl_dyn(ixt,i,k)= 0.0
2782            d_xts_dyn(ixt,i,k)= 0.0
2783           enddo
2784          enddo
2785         enddo
2786!        write(*,*) 'physiq tmp 2125'
2787#endif
2788
2789       ! !! RomP >>>   td dyn traceur
2790       IF (nqtot > nqo) d_tr_dyn(:,:,:)= 0.0
2791       ! !! RomP <<<
2792       d_rneb_dyn(:,:)=0.0
2793       d_qbs_dyn(:,:)=0.0
2794       ancien_ok = .TRUE.
2795    ENDIF
2796    !
2797    ! Ajouter le geopotentiel du sol:
2798    !
2799    DO k = 1, klev
2800       DO i = 1, klon
2801          zphi(i,k) = pphi(i,k) + pphis(i)
2802       ENDDO
2803    ENDDO
2804    !
2805    ! Verifier les temperatures
2806    !
2807    !IM BEG
2808    IF (check) THEN
2809       amn=MIN(ftsol(1,is_ter),1000.)
2810       amx=MAX(ftsol(1,is_ter),-1000.)
2811       DO i=2, klon
2812          amn=MIN(ftsol(i,is_ter),amn)
2813          amx=MAX(ftsol(i,is_ter),amx)
2814       ENDDO
2815       !
2816       PRINT*,' debut avant hgardfou min max ftsol',itap,amn,amx
2817    ENDIF !(check) THEN
2818    !IM END
2819    !
2820    CALL hgardfou(t_seri,ftsol,'debutphy',abortphy)
2821    IF (abortphy==1) Print*,'ERROR ABORT hgardfou debutphy'
2822
2823    !
2824    !IM BEG
2825    IF (check) THEN
2826       amn=MIN(ftsol(1,is_ter),1000.)
2827       amx=MAX(ftsol(1,is_ter),-1000.)
2828       DO i=2, klon
2829          amn=MIN(ftsol(i,is_ter),amn)
2830          amx=MAX(ftsol(i,is_ter),amx)
2831       ENDDO
2832       !
2833       PRINT*,' debut apres hgardfou min max ftsol',itap,amn,amx
2834    ENDIF !(check) THEN
2835    !IM END
2836    !
2837    ! Mettre en action les conditions aux limites (albedo, sst, etc.).
2838    ! Prescrire l'ozone et calculer l'albedo sur l'ocean.
2839    !
2840    ! Update ozone if day change
2841    IF (MOD(itap-1,lmt_pas) == 0) THEN
2842       IF (read_climoz <= 0) THEN
2843          ! Once per day, update ozone from Royer:
2844          IF (solarlong0<-999.) then
2845             ! Generic case with evolvoing season
2846             zzz=real(days_elapsed+1)
2847          ELSE IF (abs(solarlong0-1000.)<1.e-4) then
2848             ! Particular case with annual mean insolation
2849             zzz=real(90) ! could be revisited
2850             IF (read_climoz/=-1) THEN
2851                abort_message ='read_climoz=-1 is recommended when ' &
2852                     // 'solarlong0=1000.'
2853                CALL abort_physic (modname,abort_message,1)
2854             ENDIF
2855          ELSE
2856             ! Case where the season is imposed with solarlong0
2857             zzz=real(90) ! could be revisited
2858          ENDIF
2859
2860          wo(:,:,1)=ozonecm(latitude_deg, paprs,read_climoz,rjour=zzz)
2861#ifdef REPROBUS
2862          ptrop=dyn_tropopause(t_seri, ztsol, paprs, pplay, rot)/100.
2863          DO i = 1, klon
2864             Z1=t_seri(i,itroprep(i)+1)
2865             Z2=t_seri(i,itroprep(i))
2866             fac=(Z1-Z2)/alog(pplay(i,itroprep(i)+1)/pplay(i,itroprep(i)))
2867             B=Z2-fac*alog(pplay(i,itroprep(i)))
2868             ttrop(i)= fac*alog(ptrop(i))+B
2869!       
2870             Z1= 1.e-3 * ( pphi(i,itroprep(i)+1)+pphis(i) ) / gravit
2871             Z2= 1.e-3 * ( pphi(i,itroprep(i))  +pphis(i) ) / gravit
2872             fac=(Z1-Z2)/alog(pplay(i,itroprep(i)+1)/pplay(i,itroprep(i)))
2873             B=Z2-fac*alog(pplay(i,itroprep(i)))
2874             ztrop(i)=fac*alog(ptrop(i))+B
2875          ENDDO
2876#endif
2877       ELSE
2878          !--- ro3i = elapsed days number since current year 1st january, 0h
2879          ro3i=days_elapsed+jh_cur-jh_1jan
2880          !--- scaling for old style files (360 records)
2881          IF(SIZE(time_climoz)==360.AND..NOT.ok_daily_climoz) ro3i=ro3i*360./year_len
2882          IF(adjust_tropopause) THEN
2883             CALL regr_pr_time_av(ncid_climoz, vars_climoz(1:read_climoz),   &
2884                      ro3i, 'C', press_cen_climoz, pplay, wo, paprs(:,1),    &
2885                      time_climoz ,  longitude_deg,   latitude_deg,          &
2886                      dyn_tropopause(t_seri, ztsol, paprs, pplay, rot))
2887          ELSE
2888             CALL regr_pr_time_av(ncid_climoz,  vars_climoz(1:read_climoz),  &
2889                      ro3i, 'C', press_cen_climoz, pplay, wo, paprs(:,1),    &
2890                      time_climoz )
2891          ENDIF
2892          ! Convert from mole fraction of ozone to column density of ozone in a
2893          ! cell, in kDU:
2894          FORALL (l = 1: read_climoz) wo(:, :, l) = wo(:, :, l) * rmo3 / rmd &
2895               * zmasse / dobson_u / 1e3
2896          ! (By regridding ozone values for LMDZ only once a day, we
2897          ! have already neglected the variation of pressure in one
2898          ! day. So do not recompute "wo" at each time step even if
2899          ! "zmasse" changes a little.)
2900       ENDIF
2901    ENDIF
2902
2903#ifdef ISO
2904#ifdef ISOVERIF
2905#ifdef ISOTRAC
2906         DO k = 1, klev
2907           DO i = 1, klon         
2908             call iso_verif_traceur(xt_seri(1,i,k), &
2909     &           'physiq 2856: avant reevp')
2910           enddo
2911         enddo
2912#endif
2913#endif
2914#endif
2915    !
2916    ! Re-evaporer l'eau liquide nuageuse
2917    !
2918     CALL reevap (klon,klev,iflag_ice_thermo,t_seri,q_seri,ql_seri,qs_seri, &
2919   &         d_t_eva,d_q_eva,d_ql_eva,d_qi_eva &
2920#ifdef ISO
2921             ,xt_seri,xtl_seri,xts_seri,d_xt_eva,d_xtl_eva,d_xti_eva &
2922#endif
2923   &     )
2924
2925     CALL add_phys_tend &
2926            (du0,dv0,d_t_eva,d_q_eva,d_ql_eva,d_qi_eva,dqbs0,paprs,&
2927               'eva',abortphy,flag_inhib_tend,itap,0 &
2928#ifdef ISO
2929     &    ,d_xt_eva,d_xtl_eva,d_xti_eva &
2930#endif     
2931     &   )
2932    CALL prt_enerbil('eva',itap)
2933
2934    !=========================================================================
2935    ! Calculs de l'orbite.
2936    ! Necessaires pour le rayonnement et la surface (calcul de l'albedo).
2937    ! doit donc etre plac\'e avant radlwsw et pbl_surface
2938
2939    ! !!   jyg 17 Sep 2010 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2940    CALL ymds2ju(year_cur, mth_eq, day_eq,0., jD_eq)
2941    day_since_equinox = (jD_cur + jH_cur) - jD_eq
2942    !
2943    !   choix entre calcul de la longitude solaire vraie ou valeur fixee a
2944    !   solarlong0
2945    IF (solarlong0<-999.) THEN
2946       IF (new_orbit) THEN
2947          ! calcul selon la routine utilisee pour les planetes
2948          CALL solarlong(day_since_equinox, zlongi, dist)
2949       ELSE
2950          ! calcul selon la routine utilisee pour l'AR4
2951          CALL orbite(REAL(days_elapsed+1),zlongi,dist)
2952       ENDIF
2953    ELSE
2954       zlongi=solarlong0  ! longitude solaire vraie
2955       dist=1.            ! distance au soleil / moyenne
2956    ENDIF
2957
2958    IF (prt_level.ge.1) write(lunout,*)'Longitude solaire ',zlongi,solarlong0,dist
2959
2960
2961    ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2962    ! Calcul de l'ensoleillement :
2963    ! ============================
2964    ! Pour une solarlong0=1000., on calcule un ensoleillement moyen sur
2965    ! l'annee a partir d'une formule analytique.
2966    ! Cet ensoleillement est sym\'etrique autour de l'\'equateur et
2967    ! non nul aux poles.
2968    IF (abs(solarlong0-1000.)<1.e-4) THEN
2969       CALL zenang_an(iflag_cycle_diurne.GE.1,jH_cur, &
2970            latitude_deg,longitude_deg,rmu0,fract)
2971       swradcorr(:) = 1.0
2972       JrNt(:) = 1.0
2973       zrmu0(:) = rmu0(:)
2974    ELSE
2975       ! recode par Olivier Boucher en sept 2015
2976       SELECT CASE (iflag_cycle_diurne)
2977       CASE(0) 
2978          !  Sans cycle diurne
2979          CALL angle(zlongi, latitude_deg, fract, rmu0)
2980          swradcorr = 1.0
2981          JrNt = 1.0
2982          zrmu0 = rmu0
2983       CASE(1) 
2984          !  Avec cycle diurne sans application des poids
2985          !  bit comparable a l ancienne formulation cycle_diurne=true
2986          !  on integre entre gmtime et gmtime+radpas
2987          zdtime=phys_tstep*REAL(radpas) ! pas de temps du rayonnement (s)
2988          CALL zenang(zlongi,jH_cur,0.0,zdtime, &
2989               latitude_deg,longitude_deg,rmu0,fract)
2990          zrmu0 = rmu0
2991          swradcorr = 1.0
2992          ! Calcul du flag jour-nuit
2993          JrNt = 0.0
2994          WHERE (fract.GT.0.0) JrNt = 1.0
2995       CASE(2) 
2996          !  Avec cycle diurne sans application des poids
2997          !  On integre entre gmtime-pdtphys et gmtime+pdtphys*(radpas-1)
2998          !  Comme cette routine est appele a tous les pas de temps de
2999          !  la physique meme si le rayonnement n'est pas appele je
3000          !  remonte en arriere les radpas-1 pas de temps
3001          !  suivant. Petite ruse avec MOD pour prendre en compte le
3002          !  premier pas de temps de la physique pendant lequel
3003          !  itaprad=0
3004          zdtime1=phys_tstep*REAL(-MOD(itaprad,radpas)-1)     
3005          zdtime2=phys_tstep*REAL(radpas-MOD(itaprad,radpas)-1)
3006          CALL zenang(zlongi,jH_cur,zdtime1,zdtime2, &
3007               latitude_deg,longitude_deg,rmu0,fract)
3008          !
3009          ! Calcul des poids
3010          !
3011          zdtime1=-phys_tstep !--on corrige le rayonnement pour representer le
3012          zdtime2=0.0    !--pas de temps de la physique qui se termine
3013          CALL zenang(zlongi,jH_cur,zdtime1,zdtime2, &
3014               latitude_deg,longitude_deg,zrmu0,zfract)
3015          swradcorr = 0.0
3016          WHERE (rmu0.GE.1.e-10 .OR. fract.GE.1.e-10) &
3017               swradcorr=zfract/fract*zrmu0/rmu0
3018          ! Calcul du flag jour-nuit
3019          JrNt = 0.0
3020          WHERE (zfract.GT.0.0) JrNt = 1.0
3021       END SELECT
3022    ENDIF
3023    sza_o = ACOS (rmu0) *180./pi
3024
3025    IF (mydebug) THEN
3026       CALL writefield_phy('u_seri',u_seri,nbp_lev)
3027       CALL writefield_phy('v_seri',v_seri,nbp_lev)
3028       CALL writefield_phy('t_seri',t_seri,nbp_lev)
3029       CALL writefield_phy('q_seri',q_seri,nbp_lev)
3030    ENDIF
3031
3032    !cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
3033    ! Appel au pbl_surface : Planetary Boudary Layer et Surface
3034    ! Cela implique tous les interactions des sous-surfaces et la
3035    ! partie diffusion turbulent du couche limit.
3036    !
3037    ! Certains varibales de sorties de pbl_surface sont utiliser que pour
3038    ! ecriture des fihiers hist_XXXX.nc, ces sont :
3039    !   qsol,      zq2m,      s_pblh,  s_lcl,
3040    !   s_capCL,   s_oliqCL,  s_cteiCL,s_pblT,
3041    !   s_therm,   s_trmb1,   s_trmb2, s_trmb3,
3042    !   zu10m,     zv10m,   fder,
3043    !   zxqsurf,   delta_qsurf,
3044    !   rh2m,      zxfluxu, zxfluxv,
3045    !   frugs,     agesno,    fsollw,  fsolsw,
3046    !   d_ts,      fevap,     fluxlat, t2m,
3047    !   wfbils,    wfbilo,    fluxt,   fluxu, fluxv,
3048    !
3049    ! Certains ne sont pas utiliser du tout :
3050    !   dsens, devap, zxsnow, zxfluxt, zxfluxq, q2m, fluxq
3051    !
3052
3053    ! Calcul de l'humidite de saturation au niveau du sol
3054
3055
3056
3057    IF (iflag_pbl/=0) THEN
3058
3059       !jyg+nrlmd<
3060!!jyg       IF (prt_level .ge. 2 .and. mod(iflag_pbl_split,2) .eq. 1) THEN
3061       IF (prt_level .ge. 2 .and. mod(iflag_pbl_split,10) .ge. 1) THEN
3062          print *,'debut du splitting de la PBL, wake_s = ', wake_s(:)
3063          print *,'debut du splitting de la PBL, wake_deltat = ', wake_deltat(:,1)
3064          print *,'debut du splitting de la PBL, wake_deltaq = ', wake_deltaq(:,1)
3065       ENDIF
3066       ! !!
3067       !>jyg+nrlmd
3068       !
3069       !-------gustiness calculation-------!
3070       !ym : Warning gustiness non inialized for iflag_gusts=2 & iflag_gusts=3
3071       gustiness=0  !ym missing init
3072       
3073       IF (iflag_gusts==0) THEN
3074          gustiness(1:klon)=0
3075       ELSE IF (iflag_gusts==1) THEN
3076          gustiness(1:klon)=f_gust_bl*ale_bl(1:klon)+f_gust_wk*ale_wake(1:klon)
3077       ELSE IF (iflag_gusts==2) THEN
3078          gustiness(1:klon)=f_gust_bl*ale_bl_stat(1:klon)+f_gust_wk*ale_wake(1:klon)
3079          ! ELSE IF (iflag_gusts==2) THEN
3080          !    do i = 1, klon
3081          !       gustiness(i)=f_gust_bl*ale_bl(i)+sigma_wk(i)*f_gust_wk&
3082          !           *ale_wake(i) !! need to make sigma_wk accessible here
3083          !    enddo
3084          ! ELSE IF (iflag_gusts==3) THEN
3085          !    do i = 1, klon
3086          !       gustiness(i)=f_gust_bl*alp_bl(i)+f_gust_wk*alp_wake(i)
3087          !    enddo
3088       ENDIF
3089
3090       CALL pbl_surface(  &
3091            phys_tstep,     date0,     itap,    days_elapsed+1, &
3092            debut,     lafin, &
3093            longitude_deg, latitude_deg, rugoro,  zrmu0,      &
3094            sollwdown,    cldt,      &
3095            rain_fall, snow_fall, bs_fall, solsw,   solswfdiff, sollw,     &
3096            gustiness,                                &
3097            t_seri,    q_seri,   qbs_seri,  u_seri,  v_seri,    &
3098                                !nrlmd+jyg<
3099            wake_deltat, wake_deltaq, wake_cstar, wake_s, &
3100                                !>nrlmd+jyg
3101            pplay,     paprs,     pctsrf,             &
3102            ftsol,SFRWL,falb_dir,falb_dif,ustar,u10m,v10m,wstar, &
3103                                !albedo SB <<<
3104            cdragh,    cdragm,  u1,    v1,            &
3105            beta_aridity, &
3106                                !albedo SB >>>
3107                                ! albsol1,   albsol2,   sens,    evap,      &
3108            albsol_dir,   albsol_dif,   sens,    evap, snowerosion, & 
3109                                !albedo SB <<<
3110            albsol3_lic,runoff,   snowhgt,   qsnow, to_ice, sissnow, &
3111            zxtsol,    zxfluxlat, zt2m,    qsat2m,  zn2mout, &
3112            d_t_vdf,   d_q_vdf, d_qbs_vdf,  d_u_vdf, d_v_vdf, d_t_diss, &
3113                                !nrlmd<
3114                                !jyg<
3115            d_t_vdf_w, d_q_vdf_w, &
3116            d_t_vdf_x, d_q_vdf_x, &
3117            sens_x, zxfluxlat_x, sens_w, zxfluxlat_w, &
3118                                !>jyg
3119            delta_tsurf,wake_dens, &
3120            cdragh_x,cdragh_w,cdragm_x,cdragm_w, &
3121            kh,kh_x,kh_w, &
3122                                !>nrlmd
3123            coefh(1:klon,1:klev,1:nbsrf+1), coefm(1:klon,1:klev,1:nbsrf+1), &
3124            slab_wfbils,                 &
3125            qsol,      zq2m,      s_pblh,  s_lcl, &
3126                                !jyg<
3127            s_pblh_x, s_lcl_x, s_pblh_w, s_lcl_w, &
3128                                !>jyg
3129            s_capCL,   s_oliqCL,  s_cteiCL,s_pblT, &
3130            s_therm,   s_trmb1,   s_trmb2, s_trmb3, &
3131            zustar, zu10m,     zv10m,   fder, &
3132            zxqsurf, delta_qsurf,   rh2m,      zxfluxu, zxfluxv, &
3133            z0m, z0h,     agesno,    fsollw,  fsolsw, &
3134            d_ts,      fevap,     fluxlat, t2m, &
3135            wfbils, wfbilo, wfevap, wfrain, wfsnow, &
3136            fluxt,   fluxu,  fluxv, &
3137            dsens,     devap,     zxsnow, &
3138            zxfluxt,   zxfluxq,  zxfluxqbs,  q2m, fluxq, fluxqbs, pbl_tke, &
3139                                !nrlmd+jyg<
3140            wake_delta_pbl_TKE, &
3141                                !>nrlmd+jyg
3142             treedrg &           
3143#ifdef ISO
3144     &   ,xtrain_fall, xtsnow_fall,xt_seri, &
3145     &   wake_deltaxt,xtevap,fxtevap, &
3146     &   d_xt_vdf,d_xt_vdf_w,d_xt_vdf_x, &
3147     &   xtsol,dxtevap,zxxtsnow,zxfluxxt,fluxxt, &
3148     &   h1_diag,runoff_diag,xtrunoff_diag &
3149#endif     
3150     &   )
3151!FC
3152
3153#ifdef ISO
3154       ! write(*,*) 'physiq 2402: apres pbl_surface'
3155#ifdef ISOVERIF
3156      do i=1,klon
3157       do k=1,klev
3158        do ixt=1,ntraciso
3159            call iso_verif_noNaN(d_xt_vdf(ixt,i,k),'physiq 1993a')
3160            call iso_verif_noNaN(xt_seri(ixt,i,k),'physiq 1993b')
3161        enddo !do ixt=1,ntraciso   
3162       enddo
3163      enddo
3164#endif 
3165#ifdef ISOVERIF   
3166      do i=1,klon
3167       do k=1,klev     
3168#ifdef ISOTRAC     
3169        call iso_verif_traceur_justmass(d_xt_vdf(1,i,k),'physiq 2443')
3170#endif           
3171       enddo
3172      enddo
3173       
3174      ! verif iso_eau
3175      !write(*,*) 'physiq tmp 2748: iso_eau=',iso_eau
3176      !write(*,*) 'iso_eau.gt.0=',iso_eau.gt.0
3177      !write(*,*) 'd_xt_vdf(iso_eau,1,1),d_q_vdf(1,1)=',d_xt_vdf(iso_eau,1,1),d_q_vdf(1,1)
3178      !write(*,*) 'bidouille_anti_divergence=',bidouille_anti_divergence
3179      if (iso_eau.gt.0) then
3180!        write(*,*) 'physiq 2665: d_q_vdf,d_xt_vdf(iso_eau,554,19)=',d_q_vdf(554,19),d_xt_vdf(iso_eau,554,19)
3181!        write(*,*) 'd_q_vdf,d_xt_vdf(iso_eau,2,1)=',d_q_vdf(2,1),d_xt_vdf(iso_eau,2,1)
3182      do i=1,klon
3183!        write(*,*) 'physiq 2667: i,k=',i,k
3184!        write(*,*) 'd_q_vdf,d_xt_vdf(iso_eau,554,19)=',d_q_vdf(554,19),d_xt_vdf(iso_eau,554,19)
3185!        write(*,*) 'd_q_vdf,d_xt_vdf(iso_eau,2,1)=',d_q_vdf(2,1),d_xt_vdf(iso_eau,2,1)
3186        do k=1,klev
3187!          write(*,*) 'physiq 2670: i,k=',i,k
3188!        write(*,*) 'd_q_vdf,d_xt_vdf(iso_eau,554,19)=',d_q_vdf(554,19),d_xt_vdf(iso_eau,554,19)
3189!        write(*,*) 'd_q_vdf,d_xt_vdf(iso_eau,2,1)=',d_q_vdf(2,1),d_xt_vdf(iso_eau,2,1)
3190!        write(*,*) 'd_q_vdf,d_xt_vdf(iso_eau,i,k)=',d_q_vdf(i,k),d_xt_vdf(iso_eau,i,k)
3191          call iso_verif_egalite_choix( &
3192     &           d_xt_vdf(iso_eau,i,k),d_q_vdf(i,k), &
3193     &           'physiq 1984',errmax,errmaxrel)
3194         call iso_verif_egalite_choix( &
3195     &           xt_seri(iso_eau,i,k),q_seri(i,k), &
3196     &           'physiq 1985',errmax,errmaxrel)
3197          do nsrf=1,nbsrf
3198              call iso_verif_egalite_choix(fluxxt(iso_eau,i,k,nsrf), &
3199     &           fluxq(i,k,nsrf),'physiq 1991',errmax,errmaxrel)
3200          enddo !do nsrf=1,nbsrf
3201        enddo !k=1,klev
3202      enddo  !i=1,klon
3203      endif !if (iso_eau.gt.0) then
3204      if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then 
3205        do i=1,klon
3206         do k=1,klev
3207            if (q_seri(i,k).gt.ridicule) then 
3208               if (iso_verif_o18_aberrant_nostop( &
3209     &              xt_seri(iso_HDO,i,k)/q_seri(i,k), &
3210     &              xt_seri(iso_O18,i,k)/q_seri(i,k), &
3211     &              'physiq 4177, apres pbl_surface').eq.1) then
3212                  write(*,*) 'i,k,q_seri(i,k)=',i,k,q_seri(i,k)                       
3213                  stop
3214              endif !  if (iso_verif_o18_aberrant_nostop
3215            endif !if (q_seri(i,k).gt.errmax) then   
3216          enddo !k=1,klev
3217         enddo  !i=1,klon
3218        endif !if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then
3219        !write(*,*) 'physiq tmp 2689'
3220#endif
3221!#ifdef ISOVERIF
3222#endif
3223!#ifdef ISO
3224       !
3225       !  Add turbulent diffusion tendency to the wake difference variables
3226!!jyg       IF (mod(iflag_pbl_split,2) .NE. 0) THEN
3227       IF (mod(iflag_pbl_split,10) .NE. 0) THEN
3228!jyg<
3229          d_deltat_vdf(:,:) = d_t_vdf_w(:,:)-d_t_vdf_x(:,:)
3230          d_deltaq_vdf(:,:) = d_q_vdf_w(:,:)-d_q_vdf_x(:,:)
3231          CALL add_wake_tend &
3232             (d_deltat_vdf, d_deltaq_vdf, dsig0, ddens0, ddens0, wkoccur1, 'vdf', abortphy &
3233#ifdef ISO
3234               ,d_deltaxt_vdf &
3235#endif
3236                )
3237       ELSE
3238          d_deltat_vdf(:,:) = 0.
3239          d_deltaq_vdf(:,:) = 0.
3240!>jyg
3241#ifdef ISO
3242          d_deltaxt_vdf(:,:,:) = 0.
3243#endif
3244       ENDIF
3245
3246       !---------------------------------------------------------------------
3247       ! ajout des tendances de la diffusion turbulente
3248       IF (klon_glo==1) THEN
3249          CALL add_pbl_tend &
3250               (d_u_vdf,d_v_vdf,d_t_vdf+d_t_diss,d_q_vdf,dql0,dqi0,d_qbs_vdf,paprs,&
3251               'vdf',abortphy,flag_inhib_tend,itap &
3252#ifdef ISO
3253     &    ,d_xt_vdf,dxtl0,dxti0 &
3254#endif     
3255     &   )
3256       ELSE
3257          CALL add_phys_tend &
3258               (d_u_vdf,d_v_vdf,d_t_vdf+d_t_diss,d_q_vdf,dql0,dqi0,d_qbs_vdf,paprs,&
3259               'vdf',abortphy,flag_inhib_tend,itap,0 &
3260#ifdef ISO
3261     &    ,d_xt_vdf,dxtl0,dxti0 &
3262#endif     
3263     &   )
3264       ENDIF
3265#ifdef ISOVERIF
3266        write(*,*) 'physiq tmp 2736'
3267#endif
3268
3269       CALL prt_enerbil('vdf',itap)
3270       !--------------------------------------------------------------------
3271
3272       IF (mydebug) THEN
3273          CALL writefield_phy('u_seri',u_seri,nbp_lev)
3274          CALL writefield_phy('v_seri',v_seri,nbp_lev)
3275          CALL writefield_phy('t_seri',t_seri,nbp_lev)
3276          CALL writefield_phy('q_seri',q_seri,nbp_lev)
3277       ENDIF
3278
3279       !albedo SB >>>
3280       albsol1=0.
3281       albsol2=0.
3282       falb1=0.
3283       falb2=0.
3284       SELECT CASE(nsw)
3285       CASE(2)
3286          albsol1=albsol_dir(:,1)
3287          albsol2=albsol_dir(:,2)
3288          falb1=falb_dir(:,1,:)
3289          falb2=falb_dir(:,2,:)
3290       CASE(4)
3291          albsol1=albsol_dir(:,1)
3292          albsol2=albsol_dir(:,2)*SFRWL(2)+albsol_dir(:,3)*SFRWL(3) &
3293               +albsol_dir(:,4)*SFRWL(4)
3294          albsol2=albsol2/(SFRWL(2)+SFRWL(3)+SFRWL(4))
3295          falb1=falb_dir(:,1,:)
3296          falb2=falb_dir(:,2,:)*SFRWL(2)+falb_dir(:,3,:)*SFRWL(3) &
3297               +falb_dir(:,4,:)*SFRWL(4)
3298          falb2=falb2/(SFRWL(2)+SFRWL(3)+SFRWL(4))
3299       CASE(6)
3300          albsol1=albsol_dir(:,1)*SFRWL(1)+albsol_dir(:,2)*SFRWL(2) &
3301               +albsol_dir(:,3)*SFRWL(3)
3302          albsol1=albsol1/(SFRWL(1)+SFRWL(2)+SFRWL(3))
3303          albsol2=albsol_dir(:,4)*SFRWL(4)+albsol_dir(:,5)*SFRWL(5) &
3304               +albsol_dir(:,6)*SFRWL(6)
3305          albsol2=albsol2/(SFRWL(4)+SFRWL(5)+SFRWL(6))
3306          falb1=falb_dir(:,1,:)*SFRWL(1)+falb_dir(:,2,:)*SFRWL(2) &
3307               +falb_dir(:,3,:)*SFRWL(3)
3308          falb1=falb1/(SFRWL(1)+SFRWL(2)+SFRWL(3))
3309          falb2=falb_dir(:,4,:)*SFRWL(4)+falb_dir(:,5,:)*SFRWL(5) &
3310               +falb_dir(:,6,:)*SFRWL(6)
3311          falb2=falb2/(SFRWL(4)+SFRWL(5)+SFRWL(6))
3312       END SELECt
3313       !albedo SB <<<
3314
3315
3316       CALL evappot(klon,nbsrf,ftsol,pplay(:,1),cdragh, &
3317            t_seri(:,1),q_seri(:,1),u_seri(:,1),v_seri(:,1),evap_pot)
3318
3319    ENDIF
3320
3321    ! ==================================================================
3322    ! Blowing snow sublimation and sedimentation
3323
3324    d_t_bs(:,:)=0.
3325    d_q_bs(:,:)=0.
3326    d_qbs_bs(:,:)=0.
3327    bsfl(:,:)=0.
3328    bs_fall(:)=0.
3329    IF (ok_bs) THEN
3330
3331     CALL call_blowing_snow_sublim_sedim(klon,klev,phys_tstep,t_seri,q_seri,qbs_seri,pplay,paprs, &
3332                                        d_t_bs,d_q_bs,d_qbs_bs,bsfl,bs_fall)
3333
3334     CALL add_phys_tend &
3335               (du0,dv0,d_t_bs,d_q_bs,dql0,dqi0,d_qbs_bs,paprs,&
3336               'bs',abortphy,flag_inhib_tend,itap,0  &
3337#ifdef ISO                                       
3338       &    ,dxt0,dxtl0,dxti0 &                   
3339#endif                                       
3340       &   )
3341
3342    ENDIF
3343
3344    ! =================================================================== c
3345    !   Calcul de Qsat
3346
3347    DO k = 1, klev
3348       DO i = 1, klon
3349          zx_t = t_seri(i,k)
3350          IF (thermcep) THEN
3351             zdelta = MAX(0.,SIGN(1.,rtt-zx_t))
3352             zx_qs  = r2es * FOEEW(zx_t,zdelta)/pplay(i,k)
3353             zx_qs  = MIN(0.5,zx_qs)
3354             zcor   = 1./(1.-retv*zx_qs)
3355             zx_qs  = zx_qs*zcor
3356          ELSE
3357             !!           IF (zx_t.LT.t_coup) THEN             !jyg
3358             IF (zx_t.LT.rtt) THEN                  !jyg
3359                zx_qs = qsats(zx_t)/pplay(i,k)
3360             ELSE
3361                zx_qs = qsatl(zx_t)/pplay(i,k)
3362             ENDIF
3363          ENDIF
3364          zqsat(i,k)=zx_qs
3365       ENDDO
3366    ENDDO
3367
3368    IF (prt_level.ge.1) THEN
3369       write(lunout,*) 'L   qsat (g/kg) avant clouds_gno'
3370       write(lunout,'(i4,f15.4)') (k,1000.*zqsat(igout,k),k=1,klev)
3371    ENDIF
3372    !
3373    ! Appeler la convection (au choix)
3374    !
3375    DO k = 1, klev
3376       DO i = 1, klon
3377          conv_q(i,k) = d_q_dyn(i,k)  &
3378               + d_q_vdf(i,k)/phys_tstep
3379          conv_t(i,k) = d_t_dyn(i,k)  &
3380               + d_t_vdf(i,k)/phys_tstep
3381       ENDDO
3382    ENDDO
3383    IF (check) THEN
3384       za = qcheck(klon,klev,paprs,q_seri,ql_seri,cell_area)
3385       WRITE(lunout,*) "avantcon=", za
3386    ENDIF
3387    zx_ajustq = .FALSE.
3388    IF (iflag_con.EQ.2) zx_ajustq=.TRUE.
3389    IF (zx_ajustq) THEN
3390       DO i = 1, klon
3391          z_avant(i) = 0.0
3392       ENDDO
3393       DO k = 1, klev
3394          DO i = 1, klon
3395             z_avant(i) = z_avant(i) + (q_seri(i,k)+ql_seri(i,k)) &
3396                  *(paprs(i,k)-paprs(i,k+1))/RG
3397          ENDDO
3398       ENDDO
3399    ENDIF
3400
3401    ! Calcule de vitesse verticale a partir de flux de masse verticale
3402    DO k = 1, klev
3403       DO i = 1, klon
3404          omega(i,k) = RG*flxmass_w(i,k) / cell_area(i)
3405       ENDDO
3406    ENDDO
3407
3408    IF (prt_level.ge.1) write(lunout,*) 'omega(igout, :) = ', &
3409         omega(igout, :)
3410    !
3411    ! Appel de la convection tous les "cvpas"
3412    !
3413!!jyg    IF (MOD(itapcv,cvpas).EQ.0) THEN
3414!!    print *,' physiq : itapcv, cvpas, itap-1, cvpas_0 ', &
3415!!                       itapcv, cvpas, itap-1, cvpas_0
3416    IF (MOD(itapcv,cvpas).EQ.0 .OR. MOD(itap-1,cvpas_0).EQ.0) THEN
3417
3418    !
3419    ! Mettre a zero des variables de sortie (pour securite)
3420    !
3421    pmflxr(:,:) = 0.
3422    pmflxs(:,:) = 0.
3423    wdtrainA(:,:) = 0.
3424    wdtrainS(:,:) = 0.
3425    wdtrainM(:,:) = 0.
3426    upwd(:,:) = 0.
3427    dnwd(:,:) = 0.
3428    ep(:,:) = 0.
3429    da(:,:)=0.
3430    mp(:,:)=0.
3431    wght_cvfd(:,:)=0.
3432    phi(:,:,:)=0.
3433    phi2(:,:,:)=0.
3434    epmlmMm(:,:,:)=0.
3435    eplaMm(:,:)=0.
3436    d1a(:,:)=0.
3437    dam(:,:)=0.
3438    elij(:,:,:)=0.
3439    ev(:,:)=0.
3440    qtaa(:,:)=0.
3441    clw(:,:)=0.
3442    sij(:,:,:)=0.
3443#ifdef ISO
3444    xtclw(:,:,:)=0.
3445    xtev(:,:,:)=0.
3446    xtwdtrainA(:,:,:) = 0.
3447#endif
3448    !
3449    IF (iflag_con.EQ.1) THEN
3450       abort_message ='reactiver le call conlmd dans physiq.F'
3451       CALL abort_physic (modname,abort_message,1)
3452       !     CALL conlmd (phys_tstep, paprs, pplay, t_seri, q_seri, conv_q,
3453       !    .             d_t_con, d_q_con,
3454       !    .             rain_con, snow_con, ibas_con, itop_con)
3455    ELSE IF (iflag_con.EQ.2) THEN
3456#ifdef ISO
3457      CALL abort_gcm('physiq 2770','isos pas prevus ici',1)
3458#endif
3459       CALL conflx(phys_tstep, paprs, pplay, t_seri, q_seri, &
3460            conv_t, conv_q, -evap, omega, &
3461            d_t_con, d_q_con, rain_con, snow_con, &
3462            pmfu, pmfd, pen_u, pde_u, pen_d, pde_d, &
3463            kcbot, kctop, kdtop, pmflxr, pmflxs)
3464       d_u_con = 0.
3465       d_v_con = 0.
3466
3467       WHERE (rain_con < 0.) rain_con = 0.
3468       WHERE (snow_con < 0.) snow_con = 0.
3469       DO i = 1, klon
3470          ibas_con(i) = klev+1 - kcbot(i)
3471          itop_con(i) = klev+1 - kctop(i)
3472       ENDDO
3473    ELSE IF (iflag_con.GE.3) THEN
3474       ! nb of tracers for the KE convection:
3475       ! MAF la partie traceurs est faite dans phytrac
3476       ! on met ntra=1 pour limiter les appels mais on peut
3477       ! supprimer les calculs / ftra.
3478
3479         
3480#ifdef ISOVERIF
3481       do k = 1, klev
3482        do i = 1, klon
3483          call iso_verif_positif(q_seri(i,k),'physic 2929')
3484        enddo
3485       enddo
3486#endif
3487#ifdef ISO
3488#ifdef ISOVERIF
3489        write(*,*) 'physiq 2877'
3490      call iso_verif_noNaN_vect2D( &
3491     &           wake_deltaxt, &
3492     &           'physiq 2704c, wake_deltaxt',ntraciso,klon,klev)     
3493      do k = 1, klev
3494        do i = 1, klon
3495           do ixt=1,ntraciso
3496             call iso_verif_noNaN(xt_seri(ixt,i,k),'physiq 2757')
3497           enddo ! do ixt=1,ntraciso
3498        enddo !do i = 1, klon
3499       enddo    !do k = 1, klev
3500      if (iso_eau.gt.0) then
3501      call iso_verif_egalite_vect2D( &
3502     &           xt_seri,q_seri, &
3503     &           'physiq 2704a, avant calcul undi',ntraciso,klon,klev)
3504      call iso_verif_egalite_vect2D( &
3505     &           wake_deltaxt,wake_deltaq, &
3506     &           'physiq 2704b, wake_deltaq',ntraciso,klon,klev)
3507      endif
3508      if (iso_HDO.gt.0) then
3509      call iso_verif_aberrant_enc_vect2D( &
3510     &           xt_seri,q_seri, &
3511     &           'physiq 2709, avant calcul undi',ntraciso,klon,klev)
3512      endif
3513      if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then 
3514        do k = 1, klev
3515        do i = 1, klon
3516            if (q_seri(i,k).gt.ridicule) then 
3517               if (iso_verif_o18_aberrant_nostop( &
3518     &              xt_seri(iso_HDO,i,k)/q_seri(i,k), &
3519     &              xt_seri(iso_O18,i,k)/q_seri(i,k), &
3520     &              'physiq 2947, avant calcul undi').eq.1) then
3521                  write(*,*) 'physic 3315: i,k,q_seri(i,k)=',i,k,q_seri(i,k)                       
3522                  stop
3523              endif !  if (iso_verif_o18_aberrant_nostop
3524            endif !if (q_seri(i,k).gt.errmax) then   
3525        enddo !do i = 1, klon
3526        enddo !do k = 1, klev
3527        endif !if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then 
3528#endif
3529#endif   
3530
3531       ntra = 1
3532
3533       !=======================================================================
3534       !ajout pour la parametrisation des poches froides: calcul de
3535       !t_w et t_x: si pas de poches froides, t_w=t_x=t_seri
3536       IF (iflag_wake>=1) THEN
3537         DO k=1,klev
3538            DO i=1,klon
3539                t_w(i,k) = t_seri(i,k) + (1-wake_s(i))*wake_deltat(i,k)
3540                q_w(i,k) = q_seri(i,k) + (1-wake_s(i))*wake_deltaq(i,k)
3541                t_x(i,k) = t_seri(i,k) - wake_s(i)*wake_deltat(i,k)
3542                q_x(i,k) = q_seri(i,k) - wake_s(i)*wake_deltaq(i,k)
3543#ifdef ISO
3544                do ixt=1,ntraciso
3545                  xt_w(ixt,i,k) = xt_seri(ixt,i,k) + (1-wake_s(i))*wake_deltaxt(ixt,i,k)
3546                  xt_x(ixt,i,k) = xt_seri(ixt,i,k) - wake_s(i)*wake_deltaxt(ixt,i,k)
3547                enddo
3548#ifdef ISOVERIF
3549                if (iso_eau.gt.0) then
3550                  call iso_verif_egalite(wake_deltaq(i,k),wake_deltaxt(iso_eau,i,k),'physiq 3337')
3551                  call iso_verif_egalite(q_w(i,k),xt_w(iso_eau,i,k),'physiq 3338')
3552                endif
3553                if (iso_HDO.gt.0) then
3554                  if ((iso_verif_aberrant_enc_choix_nostop(xt_x(iso_hdo,i,k),q_x(i,k), &
3555                        ridicule,deltalim,'physic 3462a xt_x').eq.1).or. &
3556                        (iso_verif_aberrant_enc_choix_nostop(xt_w(iso_hdo,i,k),q_w(i,k), &
3557                        ridicule,deltalim,'physic 3462b xt_w').eq.1)) then
3558                     write(*,*) 'i,k=',i,k
3559                     write(*,*) 'q_x(i,k),q_seri(i,k),wake_s(i),wake_deltaq(i,k)=', &
3560                                q_x(i,k),q_seri(i,k),wake_s(i),wake_deltaq(i,k)
3561                     write(*,*) 'xt_x(iso_hdo,i,k),xt_seri(iso_hdo,i,k),wake_s(i),wake_deltaxt(iso_hdo,i,k)=', &
3562                                xt_x(iso_hdo,i,k),xt_seri(iso_hdo,i,k),wake_s(i),wake_deltaxt(iso_hdo,i,k)
3563                     write(*,*) 'deltaD_seri,wake=',deltaD(xt_seri(iso_hdo,i,k)/q_seri(i,k)), &
3564                                deltaD(wake_deltaxt(iso_hdo,i,k)/wake_deltaq(i,k))
3565                     write(*,*) 'deltaD_x,deltaD_w=',deltaD(xt_x(iso_hdo,i,k)/q_x(i,k)),deltaD(xt_w(iso_hdo,i,k)/q_w(i,k))
3566                     stop
3567                  endif
3568                endif
3569#endif
3570#endif
3571            ENDDO
3572         ENDDO
3573       ELSE
3574                t_w(:,:) = t_seri(:,:)
3575                q_w(:,:) = q_seri(:,:)
3576                t_x(:,:) = t_seri(:,:)
3577                q_x(:,:) = q_seri(:,:)
3578#ifdef ISO
3579                do ixt=1,ntraciso
3580                  xt_w(ixt,:,:) = xt_seri(ixt,:,:)
3581                  xt_x(ixt,:,:) = xt_seri(ixt,:,:)
3582                enddo
3583#endif
3584       ENDIF !IF (iflag_wake>=1) THEN
3585
3586#ifdef ISOVERIF
3587         if (iso_eau.gt.0) then
3588           DO k=1,klev
3589            DO i=1,klon
3590                  call iso_verif_egalite(q_w(i,k),xt_w(iso_eau,i,k),'physiq 3358')
3591            enddo
3592          enddo
3593        endif ! if (iso_eau.gt.0) then
3594#endif
3595       !
3596       !jyg<
3597       ! Perform dry adiabatic adjustment on wake profile
3598       ! The corresponding tendencies are added to the convective tendencies
3599       ! after the call to the convective scheme.
3600       IF (iflag_wake>=1) then
3601          IF (iflag_adjwk >= 1) THEN
3602             limbas(:) = 1
3603             CALL ajsec(paprs, pplay, t_w, q_w, limbas, &
3604                  d_t_adjwk, d_q_adjwk &
3605#ifdef ISO
3606     &      ,xt_w,d_xt_adjwk &     
3607#endif         
3608     &   )
3609             !
3610             DO k=1,klev
3611                DO i=1,klon
3612                   IF (wake_s(i) .GT. 1.e-3) THEN
3613                      t_w(i,k) = t_w(i,k) + d_t_adjwk(i,k)
3614                      q_w(i,k) = q_w(i,k) + d_q_adjwk(i,k)
3615                      d_deltat_ajs_cv(i,k) = d_t_adjwk(i,k)
3616                      d_deltaq_ajs_cv(i,k) = d_q_adjwk(i,k)
3617#ifdef ISO
3618                      do ixt=1,ntraciso
3619                       xt_w(ixt,i,k) = xt_w(ixt,i,k) + d_xt_adjwk(ixt,i,k)
3620                       d_deltaxt_ajs_cv(ixt,i,k) = d_xt_adjwk(ixt,i,k)
3621                      enddo
3622#endif
3623                   ELSE
3624                      d_deltat_ajs_cv(i,k) = 0.
3625                      d_deltaq_ajs_cv(i,k) = 0.
3626#ifdef ISO
3627                      do ixt=1,ntraciso
3628                       d_deltaxt_ajs_cv(ixt,i,k) = 0.
3629                      enddo
3630#endif
3631                   ENDIF
3632                ENDDO
3633             ENDDO
3634             IF (iflag_adjwk == 2) THEN
3635               CALL add_wake_tend &
3636                 (d_deltat_ajs_cv, d_deltaq_ajs_cv, dsig0, ddens0, ddens0, wkoccur1, 'ajs_cv', abortphy &
3637#ifdef ISO
3638                      ,d_deltaxt_ajs_cv  &
3639#endif
3640                 )
3641             ENDIF  ! (iflag_adjwk == 2)
3642          ENDIF  ! (iflag_adjwk >= 1)
3643       ENDIF ! (iflag_wake>=1)
3644       !>jyg
3645       !
3646       
3647!!      print *,'physiq. q_w(1,k), q_x(1,k) ', &
3648!!             (k, q_w(1,k), q_x(1,k),k=1,25)
3649
3650!jyg<
3651       CALL alpale( debut, itap, phys_tstep, paprs, omega, t_seri,   &
3652                    alp_offset, it_wape_prescr,  wape_prescr, fip_prescr, &
3653                    ale_bl_prescr, alp_bl_prescr, &
3654                    wake_pe, wake_fip,  &
3655                    Ale_bl, Ale_bl_trig, Alp_bl, &
3656                    Ale, Alp , Ale_wake, Alp_wake)
3657!>jyg
3658!
3659       ! sb, oct02:
3660       ! Schema de convection modularise et vectorise:
3661       ! (driver commun aux versions 3 et 4)
3662       !
3663       IF (ok_cvl) THEN ! new driver for convectL
3664          !
3665       
3666#ifdef ISO
3667#ifdef ISOVERIF
3668       write(*,*) 'physic 2553: avant appel concvl'
3669       do k = 1, klev
3670        do i = 1, klon
3671           do ixt=1,ntraciso
3672             call iso_verif_noNaN(xt_seri(ixt,i,k),'physiq 2925a')
3673             call iso_verif_noNaN(xt_x(ixt,i,k),'physiq 2925b')
3674           enddo ! do ixt=1,ntraciso
3675           call iso_verif_positif(q_seri(i,k),'physic 3133a')
3676           if (iso_eau.gt.0) then
3677             call iso_verif_egalite_choix( &
3678     &           xt_seri(iso_eau,i,k),q_seri(i,k), &
3679     &          'physic 2559a',errmax,errmaxrel)
3680             call iso_verif_egalite_choix( &
3681     &           xt_x(iso_eau,i,k),q_x(i,k), &
3682     &          'physic 2559b',errmax,errmaxrel)
3683             call iso_verif_egalite_choix( &
3684     &           xt_w(iso_eau,i,k),q_w(i,k), &
3685     &          'physic 2559c',errmax,errmaxrel)
3686           endif !if (iso_eau.gt.0) then
3687           if (iso_HDO.gt.0) then
3688             if (q_seri(i,k).gt.ridicule) then
3689              call iso_verif_aberrant_encadre( &
3690     &           xt_seri(iso_hdo,i,k)/q_seri(i,k), &
3691     &          'physic 2657a')
3692             endif !if (q_seri(i,k).gt.ridicule) then
3693             if (q_x(i,k).gt.ridicule) then
3694              call iso_verif_aberrant_encadre( &
3695     &           xt_x(iso_hdo,i,k)/q_x(i,k), &
3696     &          'physic 2657b')
3697             endif !if (q_x(i,k).gt.ridicule) then
3698             if (q_w(i,k).gt.ridicule) then
3699              call iso_verif_aberrant_encadre( &
3700     &           xt_w(iso_hdo,i,k)/q_w(i,k), &
3701     &          'physic 2657c')
3702             endif !if (q_x(i,k).gt.ridicule) then
3703           endif !if (iso_HDO.gt.0) then
3704           if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then
3705            if (q_seri(i,k).gt.ridicule) then
3706             call iso_verif_O18_aberrant( &
3707     &              xt_seri(iso_HDO,i,k)/q_seri(i,k), &
3708     &              xt_seri(iso_O18,i,k)/q_seri(i,k), &
3709     &              'physiq 3285')
3710            endif ! if (q_seri(i,k).gt.ridicule) then
3711           endif ! if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then
3712           if ((iso_O17.gt.0).and.(iso_O18.gt.0)) then
3713             if ((q_seri(i,k).gt.ridicule).and.(k.lt.nlevmaxO17)) then
3714              call iso_verif_aberrant_o17(xt_seri(iso_o17,i,k) &
3715     &           /q_seri(i,k),xt_seri(iso_o18,i,k) &
3716     &           /q_seri(i,k),'physiq 2781: avat clmain')
3717             endif !if (q_seri(i,k).gt.ridicule) then
3718           endif ! if ((iso_O17.gt.0).and.(iso_O18.gt.0)) then
3719#ifdef ISOTRAC
3720           call iso_verif_traceur(xt_seri(1,i,k),'physiq 2975a')
3721           call iso_verif_traceur(xt_x(1,i,k),'physiq 2975b')
3722           call iso_verif_traceur(xt_w(1,i,k),'physiq 2975b')
3723           if (iso_verif_tracpos_choix_nostop(xt_seri(1,i,k), &
3724     &           'physiq 3060, avantconcvl',1e-5) &
3725     &           .eq.1) then
3726              write(*,*) 'i,k=',i,k
3727              stop
3728            endif !if (iso_verif_tracpos_choix_nostop(xt_seri(1,i,k),
3729            if (option_tmin.eq.1) then
3730               call iso_verif_trac17_q_deltaD(xt_seri(1,i,k), &
3731     &           'physiq 1456, avant concvl')
3732            endif   
3733#endif     
3734      enddo !do k=1,nlev
3735   enddo  !do i=1,klon   
3736#endif
3737!ISOVERIF         
3738          if ((bidouille_anti_divergence).and. &
3739     &           (iso_eau.gt.0)) then
3740           do k=1,klev
3741            do i=1,klon 
3742             xt_seri(iso_eau,i,k)= q_seri(i,k)
3743             xt_x(iso_eau,i,k)= q_x(i,k)
3744             xt_w(iso_eau,i,k)= q_w(i,k)
3745            enddo !do i=1,klon
3746           enddo !do k=1,klev
3747          endif !if ((bidouille_anti_divergence).and. &       
3748#endif
3749          !jyg<
3750          ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3751          ! Calculate the upmost level of deep convection loops: k_upper_cv
3752          !  (near 22 km)
3753          k_upper_cv = klev
3754          !izero = klon/2+1/klon
3755          !DO k = klev,1,-1
3756          !   IF (pphi(izero,k) > 22.e4) k_upper_cv = k
3757          !ENDDO
3758          ! FH : nouveau calcul base sur un profil global sans quoi
3759          ! le modele etait sensible au decoupage de domaines
3760          DO k = klev,1,-1
3761             IF (-7*log(presnivs(k)/presnivs(1)) > 25.) k_upper_cv = k
3762          ENDDO
3763          IF (prt_level .ge. 5) THEN
3764             Print *, 'upmost level of deep convection loops: k_upper_cv = ', &
3765                  k_upper_cv
3766          ENDIF
3767          !
3768          !>jyg
3769          IF (type_trac == 'repr') THEN
3770             nbtr_tmp=ntra
3771          ELSE
3772             nbtr_tmp=nbtr
3773          ENDIF
3774          !jyg   iflag_con est dans clesphys
3775          !c          CALL concvl (iflag_con,iflag_clos,
3776          CALL concvl (iflag_clos, &
3777               phys_tstep, paprs, pplay, k_upper_cv, t_x,q_x, &
3778               t_w,q_w,wake_s, &
3779               u_seri,v_seri,tr_seri,nbtr_tmp, &
3780               ALE,ALP, &
3781               sig1,w01, &
3782               d_t_con,d_q_con,fqcomp,d_u_con,d_v_con,d_tr, &
3783               rain_con, snow_con, ibas_con, itop_con, sigd, &
3784               ema_cbmf,plcl,plfc,wbeff,convoccur,upwd,dnwd,dnwd0, &
3785               Ma,mipsh,Vprecip,cape,cin,tvp,Tconv,iflagctrl, &
3786               pbase,bbase,dtvpdt1,dtvpdq1,dplcldt,dplcldr,qcondc,wd, &
3787                                ! RomP >>>
3788                                !!     .        pmflxr,pmflxs,da,phi,mp,
3789                                !!     .        ftd,fqd,lalim_conv,wght_th)
3790               pmflxr,pmflxs,da,phi,mp,phi2,d1a,dam,sij,qtaa,clw,elij, &
3791               ftd,fqd,lalim_conv,wght_th, &
3792               ev, ep,epmlmMm,eplaMm, &
3793               wdtrainA, wdtrainS, wdtrainM,wght_cvfd,qtc_cv,sigt_cv,detrain_cv, &
3794               tau_cld_cv,coefw_cld_cv,epmax_diag &
3795#ifdef ISO
3796     &       ,xt_x,xt_w,d_xt_con,xtrain_con,xtsnow_con &
3797     &       ,xtVprecip,xtVprecipi &
3798     &       ,xtclw,fxtd,xtev,xtwdtrainA &
3799#ifdef DIAGISO
3800     &       , qlp,xtlp,qvp,xtvp & ! juste diagnostique 
3801     &       , fq_detrainement,fq_ddft,fq_fluxmasse,fq_evapprecip &
3802     &       , fxt_detrainement,fxt_ddft,fxt_fluxmasse,fxt_evapprecip &
3803     &       , f_detrainement,q_detrainement,xt_detrainement &
3804#endif         
3805#endif
3806     &        )
3807          ! RomP <<<
3808
3809#ifdef ISO
3810#ifdef ISOVERIF
3811!          write(*,*) 'q_detrainement(1,:)=',q_detrainement(1,:)
3812      call iso_verif_noNaN_vect2D(d_xt_con, &
3813     &           'physiq 3203a apres conv',ntraciso,klon,klev)
3814      call iso_verif_noNaN_vect2D(xt_seri, &
3815     &           'physiq 3203b apres conv',ntraciso,klon,klev)
3816      if (iso_HDO.gt.0) then
3817      call iso_verif_aberrant_enc_vect2D( &
3818     &           xt_seri,q_seri, &
3819     &           'physiq 3619a',ntraciso,klon,klev)
3820      endif
3821      if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then
3822      call iso_verif_O18_aberrant_enc_vect2D( &
3823     &           xt_seri,q_seri, &
3824     &           'physiq 3619b',ntraciso,klon,klev)
3825      endif         
3826#endif
3827#ifdef ISOVERIF
3828#ifdef ISOTRAC
3829        call iso_verif_trac_masse_vect(d_xt_con,klon,klev, &
3830     &           'physiq 3003 apres tend concvl', &
3831     &            errmax,errmaxrel)
3832        call iso_verif_traceur_vect(xt_seri,klon,klev, &
3833     &           'physiq 3005 apres tend concvl')
3834#endif
3835#endif
3836#endif 
3837
3838          !IM begin
3839          !       print*,'physiq: cin pbase dnwd0 ftd fqd ',cin(1),pbase(1),
3840          !    .dnwd0(1,1),ftd(1,1),fqd(1,1)
3841          !IM end
3842          !IM cf. FH
3843          clwcon0=qcondc
3844          pmfu(:,:)=upwd(:,:)+dnwd(:,:)
3845          fm_cv(:,:)=upwd(:,:)+dnwd(:,:)+dnwd0(:,:)
3846          !
3847          !jyg<
3848          ! If convective tendencies are too large, then call convection
3849          !  every time step
3850          cvpas = cvpas_0
3851          DO k=1,k_upper_cv
3852             DO i=1,klon
3853               IF (d_t_con(i,k) > 6.721 .AND. d_t_con(i,k) < 6.722 .AND.&
3854                   d_q_con(i,k) > -.0002171 .AND. d_q_con(i,k) < -.0002170) THEN
3855                     dtcon_multistep_max = 3.
3856                     dqcon_multistep_max = 0.02
3857               ENDIF
3858             ENDDO
3859          ENDDO
3860!
3861          DO k=1,k_upper_cv
3862             DO i=1,klon
3863!!               IF (abs(d_t_con(i,k)) > 0.24 .OR. &
3864!!                   abs(d_q_con(i,k)) > 2.e-2) THEN
3865               IF (abs(d_t_con(i,k)) > dtcon_multistep_max .OR. &
3866                   abs(d_q_con(i,k)) > dqcon_multistep_max) THEN
3867                 cvpas = 1
3868!!                 print *,'physiq1, i,k,d_t_con(i,k),d_q_con(i,k) ', &
3869!!                                   i,k,d_t_con(i,k),d_q_con(i,k)
3870               ENDIF
3871             ENDDO
3872          ENDDO
3873!!!   Ligne a ne surtout pas remettre sans avoir murement reflechi (jyg)
3874!!!          call bcast(cvpas)
3875!!!   ------------------------------------------------------------
3876          !>jyg
3877          !
3878          DO i = 1, klon
3879             IF (iflagctrl(i).le.1) itau_con(i)=itau_con(i)+cvpas
3880          ENDDO
3881          !
3882          !jyg<
3883          !    Add the tendency due to the dry adjustment of the wake profile
3884          IF (iflag_wake>=1) THEN
3885            IF (iflag_adjwk == 2) THEN
3886              DO k=1,klev
3887                 DO i=1,klon
3888                    ftd(i,k) = ftd(i,k) + wake_s(i)*d_t_adjwk(i,k)/phys_tstep
3889                    fqd(i,k) = fqd(i,k) + wake_s(i)*d_q_adjwk(i,k)/phys_tstep
3890                    d_t_con(i,k) = d_t_con(i,k) + wake_s(i)*d_t_adjwk(i,k)
3891                    d_q_con(i,k) = d_q_con(i,k) + wake_s(i)*d_q_adjwk(i,k)
3892#ifdef ISO
3893                   do ixt=1,ntraciso
3894                     fxtd(ixt,i,k) = fxtd(ixt,i,k) + wake_s(i)*d_xt_adjwk(ixt,i,k)/phys_tstep
3895                     d_xt_con(ixt,i,k) = d_xt_con(ixt,i,k) + wake_s(i)*d_xt_adjwk(ixt,i,k)
3896                   enddo
3897#endif
3898                 ENDDO
3899              ENDDO
3900            ENDIF  ! (iflag_adjwk = 2)
3901          ENDIF   ! (iflag_wake>=1)
3902          !>jyg
3903          !
3904       ELSE ! ok_cvl
3905
3906          ! MAF conema3 ne contient pas les traceurs
3907          CALL conema3 (phys_tstep, &
3908               paprs,pplay,t_seri,q_seri, &
3909               u_seri,v_seri,tr_seri,ntra, &
3910               sig1,w01, &
3911               d_t_con,d_q_con,d_u_con,d_v_con,d_tr, &
3912               rain_con, snow_con, ibas_con, itop_con, &
3913               upwd,dnwd,dnwd0,bas,top, &
3914               Ma,cape,tvp,rflag, &
3915               pbase &
3916               ,bbase,dtvpdt1,dtvpdq1,dplcldt,dplcldr &
3917               ,clwcon0)
3918
3919       ENDIF ! ok_cvl
3920
3921       !
3922       ! Correction precip
3923       rain_con = rain_con * cvl_corr
3924       snow_con = snow_con * cvl_corr
3925#ifdef ISO
3926          xtrain_con = xtrain_con * cvl_corr
3927          xtsnow_con = xtsnow_con * cvl_corr
3928#endif
3929       !
3930
3931       IF (.NOT. ok_gust) THEN
3932          do i = 1, klon
3933             wd(i)=0.0
3934          enddo
3935       ENDIF
3936
3937       ! =================================================================== c
3938       ! Calcul des proprietes des nuages convectifs
3939       !
3940
3941       !   calcul des proprietes des nuages convectifs
3942       clwcon0(:,:)=fact_cldcon*clwcon0(:,:)
3943       IF (iflag_cld_cv == 0) THEN
3944          CALL clouds_gno &
3945               (klon,klev,q_seri,zqsat,clwcon0,ptconv,ratqsc,rnebcon0)
3946       ELSE
3947          CALL clouds_bigauss &
3948               (klon,klev,q_seri,zqsat,qtc_cv,sigt_cv,ptconv,ratqsc,rnebcon0)
3949       ENDIF
3950
3951
3952       ! =================================================================== c
3953
3954       DO i = 1, klon
3955          itop_con(i) = min(max(itop_con(i),1),klev)
3956          ibas_con(i) = min(max(ibas_con(i),1),itop_con(i))
3957       ENDDO
3958
3959       DO i = 1, klon
3960          ! C Risi modif: pour éviter pb de dépassement d'indice dans les cas
3961          ! où i n'est pas un point convectif et donc ibas_con(i)=0
3962          ! c'est un pb indépendant des isotopes
3963          if (ibas_con(i) > 0) then
3964             ema_pcb(i)  = paprs(i,ibas_con(i))
3965          else
3966             ema_pcb(i)  = 0.0
3967          endif
3968       ENDDO
3969       DO i = 1, klon
3970          ! L'idicage de itop_con peut cacher un pb potentiel
3971          ! FH sous la dictee de JYG, CR
3972          ema_pct(i)  = paprs(i,itop_con(i)+1)
3973
3974          IF (itop_con(i).gt.klev-3) THEN
3975             IF (prt_level >= 9) THEN
3976                write(lunout,*)'La convection monte trop haut '
3977                write(lunout,*)'itop_con(,',i,',)=',itop_con(i)
3978             ENDIF
3979          ENDIF
3980       ENDDO
3981    ELSE IF (iflag_con.eq.0) THEN
3982       write(lunout,*) 'On n appelle pas la convection'
3983       clwcon0=0.
3984       rnebcon0=0.
3985       d_t_con=0.
3986       d_q_con=0.
3987       d_u_con=0.
3988       d_v_con=0.
3989       rain_con=0.
3990       snow_con=0.
3991       bas=1
3992       top=1
3993#ifdef ISO
3994       d_xt_con=0.
3995       xtrain_con=0.
3996       xtsnow_con=0.
3997#endif
3998    ELSE
3999       WRITE(lunout,*) "iflag_con non-prevu", iflag_con
4000       CALL abort_physic("physiq", "", 1)
4001    ENDIF
4002
4003    !     CALL homogene(paprs, q_seri, d_q_con, u_seri,v_seri,
4004    !    .              d_u_con, d_v_con)
4005
4006!jyg    Reinitialize proba_notrig and itapcv when convection has been called
4007    proba_notrig(:) = 1.
4008    itapcv = 0
4009    ENDIF !  (MOD(itapcv,cvpas).EQ.0 .OR. MOD(itapcv,cvpas_0).EQ.0)
4010!
4011    itapcv = itapcv+1
4012    !
4013    ! Compter les steps ou cvpas=1
4014    IF (cvpas == 1) THEN
4015      Ncvpaseq1 = Ncvpaseq1+1
4016    ENDIF
4017    IF (mod(itap,1000) == 0) THEN
4018      print *,' physiq, nombre de steps ou cvpas = 1 : ', Ncvpaseq1
4019    ENDIF
4020
4021!!!jyg  Appel diagnostique a add_phys_tend pour tester la conservation de
4022!!!     l'energie dans les courants satures.
4023!!    d_t_con_sat(:,:) = d_t_con(:,:) - ftd(:,:)*dtime
4024!!    d_q_con_sat(:,:) = d_q_con(:,:) - fqd(:,:)*dtime
4025!!    dql_sat(:,:) = (wdtrainA(:,:)+wdtrainM(:,:))*dtime/zmasse(:,:)
4026!!    CALL add_phys_tend(d_u_con, d_v_con, d_t_con_sat, d_q_con_sat, dql_sat,   &
4027!!                     dqi0, paprs, 'convection_sat', abortphy, flag_inhib_tend,& 
4028!!                     itap, 1)
4029!!    call prt_enerbil('convection_sat',itap)
4030!!
4031!!
4032    CALL add_phys_tend(d_u_con, d_v_con, d_t_con, d_q_con, dql0, dqi0, dqbs0, paprs, &
4033         'convection',abortphy,flag_inhib_tend,itap,0 &
4034#ifdef ISO
4035     &    ,d_xt_con,dxtl0,dxti0 &
4036#endif     
4037     &   )
4038    CALL prt_enerbil('convection',itap)
4039
4040    !-------------------------------------------------------------------------
4041
4042    IF (mydebug) THEN
4043       CALL writefield_phy('u_seri',u_seri,nbp_lev)
4044       CALL writefield_phy('v_seri',v_seri,nbp_lev)
4045       CALL writefield_phy('t_seri',t_seri,nbp_lev)
4046       CALL writefield_phy('q_seri',q_seri,nbp_lev)
4047    ENDIF
4048
4049    IF (check) THEN
4050       za = qcheck(klon,klev,paprs,q_seri,ql_seri,cell_area)
4051       WRITE(lunout,*)"aprescon=", za
4052       zx_t = 0.0
4053       za = 0.0
4054       DO i = 1, klon
4055          za = za + cell_area(i)/REAL(klon)
4056          zx_t = zx_t + (rain_con(i)+ &
4057               snow_con(i))*cell_area(i)/REAL(klon)
4058       ENDDO
4059       zx_t = zx_t/za*phys_tstep
4060       WRITE(lunout,*)"Precip=", zx_t
4061    ENDIF
4062    IF (zx_ajustq) THEN
4063       DO i = 1, klon
4064          z_apres(i) = 0.0
4065#ifdef ISO
4066        do ixt=1,ntraciso
4067          zxt_apres(ixt,i) = 0.0
4068        enddo !do ixt=1,ntraciso
4069#endif
4070       ENDDO
4071       DO k = 1, klev
4072          DO i = 1, klon
4073             z_apres(i) = z_apres(i) + (q_seri(i,k)+ql_seri(i,k)) &
4074                  *(paprs(i,k)-paprs(i,k+1))/RG
4075          ENDDO
4076       ENDDO
4077#ifdef ISO
4078          DO k = 1, klev
4079            DO i = 1, klon
4080              do ixt=1,ntraciso
4081                zxt_apres(ixt,i) = zxt_apres(ixt,i)  &
4082     &            + (xt_seri(ixt,i,k)+xtl_seri(ixt,i,k)) &
4083     &            *(paprs(i,k)-paprs(i,k+1))/RG
4084              enddo ! do ixt=1,ntraciso
4085            ENDDO
4086          ENDDO 
4087#endif
4088       DO i = 1, klon
4089          z_factor(i) = (z_avant(i)-(rain_con(i)+snow_con(i))*phys_tstep) &
4090               /z_apres(i)
4091       ENDDO
4092#ifdef ISO
4093         DO i = 1, klon
4094            do ixt=1,ntraciso
4095              zxt_factor(ixt,i) = (zxt_avant(ixt,i)-(xtrain_con(ixt,i) &
4096     &          +xtsnow_con(ixt,i))*phys_tstep)/zxt_apres(ixt,i)
4097            enddo ! do ixt=1,ntraciso
4098         ENDDO   
4099#endif
4100       DO k = 1, klev
4101          DO i = 1, klon
4102             IF (z_factor(i).GT.(1.0+1.0E-08) .OR. &
4103                  z_factor(i).LT.(1.0-1.0E-08)) THEN
4104                q_seri(i,k) = q_seri(i,k) * z_factor(i)
4105#ifdef ISO
4106                do ixt=1,ntraciso
4107                  xt_seri(ixt,i,k)=xt_seri(ixt,i,k)*zxt_factor(ixt,i) 
4108                enddo ! do ixt=1,ntraciso
4109#endif
4110             ENDIF
4111          ENDDO
4112       ENDDO
4113    ENDIF
4114    zx_ajustq=.FALSE.
4115
4116#ifdef ISO
4117#ifdef ISOVERIF
4118     write(*,*) 'physiq 3425: apres convection'
4119      if (iso_HDO.gt.0) then
4120      call iso_verif_aberrant_enc_vect2D( &
4121     &           xt_seri,q_seri, &
4122     &           'physiq 3691a',ntraciso,klon,klev)
4123      endif
4124      if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then
4125      call iso_verif_O18_aberrant_enc_vect2D( &
4126     &           xt_seri,q_seri, &
4127     &           'physiq 3691b',ntraciso,klon,klev)
4128      endif
4129#ifdef ISOTRAC
4130        call iso_verif_traceur_vect(xt_seri,klon,klev, &
4131     &           'physiq 3386 avant wakes')
4132#endif
4133#endif
4134#endif
4135
4136    !
4137    !==========================================================================
4138    !RR:Evolution de la poche froide: on ne fait pas de separation wake/env
4139    !pour la couche limite diffuse pour l instant
4140    !
4141    !
4142    ! nrlmd le 22/03/2011---Si on met les poches hors des thermiques
4143    ! il faut rajouter cette tendance calcul\'ee hors des poches
4144    ! froides
4145    !
4146    IF (iflag_wake>=1) THEN
4147       !
4148       !
4149       ! Call wakes every "wkpas" step
4150       !
4151       IF (MOD(itapwk,wkpas).EQ.0) THEN
4152          !
4153          DO k=1,klev
4154             DO i=1,klon
4155                dt_dwn(i,k)  = ftd(i,k)
4156                dq_dwn(i,k)  = fqd(i,k)
4157                M_dwn(i,k)   = dnwd0(i,k)
4158                M_up(i,k)    = upwd(i,k)
4159                dt_a(i,k)    = d_t_con(i,k)/phys_tstep - ftd(i,k)
4160                dq_a(i,k)    = d_q_con(i,k)/phys_tstep - fqd(i,k)
4161#ifdef ISO
4162                do ixt=1,ntraciso
4163                dxt_dwn(ixt,i,k)  = fxtd(ixt,i,k) 
4164                dxt_a(ixt,i,k)    = d_xt_con(ixt,i,k)/phys_tstep - fxtd(ixt,i,k)
4165                enddo !do ixt=1,ntraciso
4166#endif
4167             ENDDO
4168          ENDDO
4169         
4170          IF (iflag_wake==2) THEN
4171             ok_wk_lsp(:)=max(sign(1.,wake_s(:)-wake_s_min_lsp),0.)
4172             DO k = 1,klev
4173                dt_dwn(:,k)= dt_dwn(:,k)+ &
4174                     ok_wk_lsp(:)*(d_t_eva(:,k)+d_t_lsc(:,k))/phys_tstep
4175                dq_dwn(:,k)= dq_dwn(:,k)+ &
4176                     ok_wk_lsp(:)*(d_q_eva(:,k)+d_q_lsc(:,k))/phys_tstep
4177#ifdef ISO
4178             do ixt=1,ntraciso
4179                dxt_dwn(:,k,ixt)= dxt_dwn(:,k,ixt)+ &
4180                     ok_wk_lsp(:)*(d_xt_eva(:,k,ixt)+d_xt_lsc(:,k,ixt))/phys_tstep
4181             enddo
4182#endif
4183             ENDDO
4184          ELSEIF (iflag_wake==3) THEN
4185             ok_wk_lsp(:)=max(sign(1.,wake_s(:)-wake_s_min_lsp),0.)
4186             DO k = 1,klev
4187                DO i=1,klon
4188                   IF (rneb(i,k)==0.) THEN
4189                      ! On ne tient compte des tendances qu'en dehors des
4190                      ! nuages (c'est-\`a-dire a priri dans une region ou
4191                      ! l'eau se reevapore).
4192                      dt_dwn(i,k)= dt_dwn(i,k)+ &
4193                           ok_wk_lsp(i)*d_t_lsc(i,k)/phys_tstep
4194                      dq_dwn(i,k)= dq_dwn(i,k)+ &
4195                           ok_wk_lsp(i)*d_q_lsc(i,k)/phys_tstep
4196#ifdef ISO
4197                      do ixt=1,ntraciso
4198                        dxt_dwn(:,k,ixt)= dxt_dwn(:,k,ixt)+ &
4199                           ok_wk_lsp(:)*d_xt_lsc(:,k,ixt)/phys_tstep
4200                      enddo
4201#endif
4202                   ENDIF
4203                ENDDO
4204             ENDDO
4205          ENDIF
4206         
4207#ifdef ISOVERIF
4208        write(*,*) 'physiq 3977: verif des inputs de calwake'
4209        DO k = 1,klev
4210          DO i=1,klon
4211            if (iso_eau.gt.0) then
4212              call iso_verif_egalite(q_seri(i,k),xt_seri(iso_eau,i,k),'physiq 3981a')
4213              call iso_verif_egalite(dq_dwn(i,k),dxt_dwn(iso_eau,i,k),'physiq 3981b')
4214              call iso_verif_egalite(dq_a(i,k),dxt_a(iso_eau,i,k),'physiq 3981c')
4215              call iso_verif_egalite(wake_deltaq(i,k),wake_deltaxt(iso_eau,i,k),'physiq 3981d')
4216            endif
4217          enddo
4218        enddo       
4219#endif
4220          !
4221          !calcul caracteristiques de la poche froide
4222          CALL calWAKE (iflag_wake_tend, paprs, pplay, phys_tstep, &
4223               t_seri, q_seri, omega,  &
4224               dt_dwn, dq_dwn, M_dwn, M_up,  &
4225               dt_a, dq_a, cv_gen,  &
4226               sigd, cin,  &
4227               wake_deltat, wake_deltaq, wake_s, awake_dens, wake_dens,  &
4228               wake_dth, wake_h,  &
4229!!               wake_pe, wake_fip, wake_gfl,  &
4230               wake_pe, wake_fip_0, wake_gfl,  &   !! jyg
4231               d_t_wake, d_q_wake,  &
4232               wake_k, t_x, q_x,  &
4233               wake_omgbdth, wake_dp_omgb,  &
4234               wake_dtKE, wake_dqKE,  &
4235               wake_omg, wake_dp_deltomg,  &
4236               wake_spread, wake_Cstar, d_deltat_wk_gw,  &
4237               d_deltat_wk, d_deltaq_wk, d_s_wk, d_dens_a_wk, d_dens_wk &
4238#ifdef ISO
4239     &               ,xt_seri,dxt_dwn,dxt_a &
4240     &               ,wake_deltaxt,d_xt_wake,xt_x,d_deltaxt_wk &
4241#endif     
4242     &   )
4243          !
4244          !jyg    Reinitialize itapwk when wakes have been called
4245          itapwk = 0
4246       ENDIF !  (MOD(itapwk,wkpas).EQ.0)
4247       !
4248       itapwk = itapwk+1
4249#ifdef ISOVERIF
4250        write(*,*) 'physiq 4020: verif des outputs de calwake'
4251        DO k = 1,klev
4252          DO i=1,klon
4253            if (iso_eau.gt.0) then
4254              call iso_verif_egalite(d_q_wake(i,k),d_xt_wake(iso_eau,i,k),'physiq 4023a')
4255              call iso_verif_egalite(q_x(i,k),xt_x(iso_eau,i,k),'physiq 4023b')
4256              call iso_verif_egalite(d_deltaq_wk(i,k),d_deltaxt_wk(iso_eau,i,k),'physiq 4023c')
4257            endif
4258          enddo
4259        enddo
4260#endif
4261       !
4262       !-----------------------------------------------------------------------
4263       ! ajout des tendances des poches froides
4264       CALL add_phys_tend(du0,dv0,d_t_wake,d_q_wake,dql0,dqi0,dqbs0,paprs,'wake', &
4265            abortphy,flag_inhib_tend,itap,0 &
4266#ifdef ISO
4267     &    ,d_xt_wake,dxtl0,dxti0 &
4268#endif     
4269     &   )
4270       CALL prt_enerbil('wake',itap)
4271       !------------------------------------------------------------------------
4272
4273       ! Increment Wake state variables
4274       IF (iflag_wake_tend .GT. 0.) THEN
4275
4276         CALL add_wake_tend &
4277            (d_deltat_wk, d_deltaq_wk, d_s_wk, d_dens_a_wk, d_dens_wk, wake_k, &
4278             'wake', abortphy &
4279#ifdef ISO
4280                ,d_deltaxt_wk &
4281#endif
4282         )
4283          CALL prt_enerbil('wake',itap)
4284       ENDIF   ! (iflag_wake_tend .GT. 0.)
4285       !
4286       IF (prt_level .GE. 10) THEN
4287         print *,' physiq, after calwake, wake_s: ',wake_s(:)
4288         print *,' physiq, after calwake, wake_deltat: ',wake_deltat(:,1)
4289         print *,' physiq, after calwake, wake_deltaq: ',wake_deltaq(:,1)
4290       ENDIF
4291
4292       IF (iflag_alp_wk_cond .GT. 0.) THEN
4293
4294         CALL alpale_wk(phys_tstep, cell_area, wake_k, wake_s, wake_dens, wake_fip_0, &
4295                        wake_fip)
4296       ELSE
4297         wake_fip(:) = wake_fip_0(:)
4298       ENDIF   ! (iflag_alp_wk_cond .GT. 0.)
4299
4300    ENDIF  ! (iflag_wake>=1)
4301    !
4302    !===================================================================
4303    ! Convection seche (thermiques ou ajustement)
4304    !===================================================================
4305    !
4306    CALL stratocu_if(klon,klev,pctsrf,paprs, pplay,t_seri &
4307         ,seuil_inversion,weak_inversion,dthmin)
4308
4309
4310
4311    d_t_ajsb(:,:)=0.
4312    d_q_ajsb(:,:)=0.
4313    d_t_ajs(:,:)=0.
4314    d_u_ajs(:,:)=0.
4315    d_v_ajs(:,:)=0.
4316    d_q_ajs(:,:)=0.
4317    clwcon0th(:,:)=0.
4318#ifdef ISO
4319      d_xt_ajs(:,:,:)=0.0
4320      d_xt_ajsb(:,:,:)=0.0
4321#endif
4322    !
4323    !      fm_therm(:,:)=0.
4324    !      entr_therm(:,:)=0.
4325    !      detr_therm(:,:)=0.
4326    !
4327    IF (prt_level>9) WRITE(lunout,*) &
4328         'AVANT LA CONVECTION SECHE , iflag_thermals=' &
4329         ,iflag_thermals,'   nsplit_thermals=',nsplit_thermals
4330    IF (iflag_thermals<0) THEN
4331       !  Rien
4332       !  ====
4333       IF (prt_level>9) WRITE(lunout,*)'pas de convection seche'
4334
4335
4336#ifdef ISO
4337#ifdef ISOVERIF
4338      call iso_verif_noNaN_vect2D(xt_seri,'physiq 3971', &
4339     &            ntraciso,klon,klev)
4340#endif         
4341#ifdef ISOVERIF         
4342      write(*,*) 'physiq 3570'
4343         do k=1,klev
4344           do i=1,klon
4345             if (iso_eau.gt.0) then
4346                call iso_verif_egalite_choix(xt_seri(iso_eau,i,k), &
4347     &                 q_seri(i,k),'physiq 2765',errmax,errmaxrel)
4348             endif !if (iso_eau.gt.0) then
4349             if (iso_HDO.gt.0) then
4350               if (q_seri(i,k).gt.ridicule) then
4351                  call iso_verif_aberrant_encadre(xt_seri(iso_HDO,i,k)/ &
4352     &                 q_seri(i,k),'physiq 2770')
4353               endif !if (abs(d_q_ajs(i,k)).gt.ridicule) then               
4354             endif !if (iso_HDO.gt.0) then
4355             if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then     
4356               if (q_seri(i,k).gt.ridicule) then 
4357                 if (iso_verif_o18_aberrant_nostop( &
4358     &              xt_seri(iso_HDO,i,k)/q_seri(i,k), &
4359     &              xt_seri(iso_O18,i,k)/q_seri(i,k), &
4360     &              'physiq 3978, avant calltherm').eq.1) then
4361                  write(*,*) 'physic 4115: i,k,q_seri(i,k)=',i,k,q_seri(i,k)                       
4362                  stop
4363                 endif !  if (iso_verif_o18_aberrant_nostop
4364               endif !if (q_seri(i,k).gt.errmax) then   
4365             endif !if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then   
4366             if ((iso_O17.gt.0).and.(iso_O18.gt.0)) then
4367               if ((q_seri(i,k).gt.ridicule).and.(k.lt.nlevmaxO17)) then
4368                 call iso_verif_aberrant_o17(xt_seri(iso_o17,i,k) &
4369     &           /q_seri(i,k),xt_seri(iso_o18,i,k) &
4370     &           /q_seri(i,k),'physiq 3178: avant clmain')
4371               endif !if (q_seri(i,k).gt.ridicule) then
4372             endif !if (iso_O17.gt.0) then
4373#ifdef ISOTRAC     
4374        call iso_verif_traceur(xt_seri(1,i,k),'physiq 3389')
4375        if (iso_verif_tracpos_choix_nostop(xt_seri(1,i,k), &
4376     &           'physiq 3481, avant ajsec',1e-5) &
4377     &           .eq.1) then
4378              write(*,*) 'i,k=',i,k
4379#ifdef ISOVERIF
4380                  stop
4381#endif
4382!#ifdef ISOVERIF     
4383        endif !if (iso_verif_tracpos_choix_nostop(xt_seri(1,i,k),             
4384#endif   
4385!#ifdef ISOTRAC                     
4386            enddo !do i=1,klon
4387         enddo !do k=1,klev
4388#endif
4389!#ifdef ISOVERIF         
4390       if ((iso_eau.gt.0).and.(bidouille_anti_divergence)) then
4391        do k=1,klev   
4392        do i=1,klon
4393            xt_seri(iso_eau,i,k)=q_seri(i,k)
4394        enddo !do i=1,klon
4395        enddo !do k=1,klev 
4396#ifdef ISOTRAC       
4397        call iso_verif_traceur_jbid_vect(xt_seri, &
4398     &            klon,klev)   
4399#endif       
4400      endif !if ((iso_eau.gt.0).and.(bidouille_anti_divergence)) then 
4401#endif
4402!#ifdef ISO     
4403
4404
4405    ELSE
4406
4407       !  Thermiques
4408       !  ==========
4409       IF (prt_level>9) WRITE(lunout,*)'JUSTE AVANT , iflag_thermals=' &
4410            ,iflag_thermals,'   nsplit_thermals=',nsplit_thermals
4411
4412
4413       !cc nrlmd le 10/04/2012
4414       DO k=1,klev+1
4415          DO i=1,klon
4416             pbl_tke_input(i,k,is_oce)=pbl_tke(i,k,is_oce)
4417             pbl_tke_input(i,k,is_ter)=pbl_tke(i,k,is_ter)
4418             pbl_tke_input(i,k,is_lic)=pbl_tke(i,k,is_lic)
4419             pbl_tke_input(i,k,is_sic)=pbl_tke(i,k,is_sic)
4420          ENDDO
4421       ENDDO
4422       !cc fin nrlmd le 10/04/2012
4423
4424       IF (iflag_thermals>=1) THEN
4425          !jyg<
4426!!       IF (mod(iflag_pbl_split/2,2) .EQ. 1) THEN
4427       IF (mod(iflag_pbl_split/10,10) .GE. 1) THEN
4428             !  Appel des thermiques avec les profils exterieurs aux poches
4429             DO k=1,klev
4430                DO i=1,klon
4431                   t_therm(i,k) = t_seri(i,k) - wake_s(i)*wake_deltat(i,k)
4432                   q_therm(i,k) = q_seri(i,k) - wake_s(i)*wake_deltaq(i,k)
4433                   u_therm(i,k) = u_seri(i,k)
4434                   v_therm(i,k) = v_seri(i,k)
4435#ifdef ISO
4436                   do ixt=1,ntraciso
4437                     xt_therm(ixt,i,k) = xt_seri(ixt,i,k) - wake_s(i)*wake_deltaxt(ixt,i,k)
4438                   enddo !do ixt=1,ntraciso
4439#endif
4440                ENDDO
4441             ENDDO
4442          ELSE
4443             !  Appel des thermiques avec les profils moyens
4444             DO k=1,klev
4445                DO i=1,klon
4446                   t_therm(i,k) = t_seri(i,k)
4447                   q_therm(i,k) = q_seri(i,k)
4448                   u_therm(i,k) = u_seri(i,k)
4449                   v_therm(i,k) = v_seri(i,k)
4450#ifdef ISO
4451                   do ixt=1,ntraciso
4452                     xt_therm(ixt,i,k) = xt_seri(ixt,i,k)
4453                   enddo !do ixt=1,ntraciso
4454#endif
4455                ENDDO
4456             ENDDO
4457          ENDIF
4458
4459
4460#ifdef ISO
4461#ifdef ISOTRAC   
4462        ! ajout C Risi juillet 2020: pas sure de ce que je fais...   
4463        call iso_verif_traceur_jbid_pos_vect(klon,klev,xt_seri)   
4464        call iso_verif_traceur_jbid_pos_vect(klon,klev,xt_therm) 
4465#endif
4466#endif
4467          !>jyg
4468          CALL calltherm(pdtphys &
4469               ,pplay,paprs,pphi,weak_inversion &
4470                        ! ,u_seri,v_seri,t_seri,q_seri,zqsat,debut & !jyg
4471               ,u_therm,v_therm,t_therm,q_therm,zqsat,debut &  !jyg
4472               ,d_u_ajs,d_v_ajs,d_t_ajs,d_q_ajs &
4473               ,fm_therm,entr_therm,detr_therm &
4474               ,zqasc,clwcon0th,lmax_th,ratqscth &
4475               ,ratqsdiff,zqsatth &
4476                                !on rajoute ale et alp, et les
4477                                !caracteristiques de la couche alim
4478               ,Ale_bl,Alp_bl,lalim_conv,wght_th, zmax0, f0, zw2,fraca &
4479               ,ztv,zpspsk,ztla,zthl &
4480                                !cc nrlmd le 10/04/2012
4481               ,pbl_tke_input,pctsrf,omega,cell_area &
4482               ,zlcl_th,fraca0,w0,w_conv,therm_tke_max0,env_tke_max0 &
4483               ,n2,s2,ale_bl_stat &
4484               ,therm_tke_max,env_tke_max &
4485               ,alp_bl_det,alp_bl_fluct_m,alp_bl_fluct_tke &
4486               ,alp_bl_conv,alp_bl_stat &
4487                                !cc fin nrlmd le 10/04/2012
4488               ,zqla,ztva &
4489#ifdef ISO         
4490     &      ,xt_therm,d_xt_ajs &
4491#ifdef DIAGISO         
4492     &      ,q_the,xt_the &
4493#endif         
4494#endif         
4495     &   )
4496          !
4497          !jyg<
4498!!jyg          IF (mod(iflag_pbl_split/2,2) .EQ. 1) THEN
4499          IF (mod(iflag_pbl_split/10,10) .GE. 1) THEN
4500             !  Si les thermiques ne sont presents que hors des
4501             !  poches, la tendance moyenne associ\'ee doit etre
4502             !  multipliee par la fraction surfacique qu'ils couvrent.
4503             DO k=1,klev
4504                DO i=1,klon
4505                   !
4506                   d_deltat_the(i,k) = - d_t_ajs(i,k)
4507                   d_deltaq_the(i,k) = - d_q_ajs(i,k)
4508                   !
4509                   d_u_ajs(i,k) = d_u_ajs(i,k)*(1.-wake_s(i))
4510                   d_v_ajs(i,k) = d_v_ajs(i,k)*(1.-wake_s(i))
4511                   d_t_ajs(i,k) = d_t_ajs(i,k)*(1.-wake_s(i))
4512                   d_q_ajs(i,k) = d_q_ajs(i,k)*(1.-wake_s(i))
4513#ifdef ISO
4514                   do ixt=1,ntraciso
4515                     d_deltaxt_the(ixt,i,k) = - d_xt_ajs(ixt,i,k)
4516                     d_xt_ajs(ixt,i,k) = d_xt_ajs(ixt,i,k)*(1.-wake_s(i))
4517                   enddo
4518#endif
4519                   !
4520                ENDDO
4521             ENDDO
4522          !
4523             IF (ok_bug_split_th) THEN
4524               CALL add_wake_tend &
4525                   (d_deltat_the, d_deltaq_the, dsig0, ddens0, ddens0, wkoccur1, 'the', abortphy &
4526#ifdef ISO
4527                   ,d_deltaxt_the &
4528#endif
4529        &)
4530             ELSE
4531               CALL add_wake_tend &
4532                   (d_deltat_the, d_deltaq_the, dsig0, ddens0, ddens0, wake_k, 'the', abortphy &
4533#ifdef ISO
4534                   ,d_deltaxt_the &
4535#endif
4536        &)
4537             ENDIF
4538             CALL prt_enerbil('the',itap)
4539          !
4540          ENDIF  ! (mod(iflag_pbl_split/10,10) .GE. 1)
4541          !
4542          CALL add_phys_tend(d_u_ajs,d_v_ajs,d_t_ajs,d_q_ajs,  &
4543                             dql0,dqi0,dqbs0,paprs,'thermals', abortphy,flag_inhib_tend,itap,0 &
4544#ifdef ISO
4545     &    ,d_xt_ajs,dxtl0,dxti0 &
4546#endif     
4547     &   )
4548          CALL prt_enerbil('thermals',itap)
4549          !
4550!
4551          CALL alpale_th( phys_tstep, lmax_th, t_seri, cell_area,  &
4552                          cin, s2, n2,  &
4553                          ale_bl_trig, ale_bl_stat, ale_bl,  &
4554                          alp_bl, alp_bl_stat, &
4555                          proba_notrig, random_notrig, cv_gen)
4556          !>jyg
4557
4558          ! ------------------------------------------------------------------
4559          ! Transport de la TKE par les panaches thermiques.
4560          ! FH : 2010/02/01
4561          !     if (iflag_pbl.eq.10) then
4562          !     call thermcell_dtke(klon,klev,nbsrf,pdtphys,fm_therm,entr_therm,
4563          !    s           rg,paprs,pbl_tke)
4564          !     endif
4565          ! -------------------------------------------------------------------
4566
4567          DO i=1,klon
4568             !           zmax_th(i)=pphi(i,lmax_th(i))/rg
4569             !CR:04/05/12:correction calcul zmax
4570             zmax_th(i)=zmax0(i)
4571          ENDDO
4572
4573       ENDIF
4574
4575       !  Ajustement sec
4576       !  ==============
4577
4578       ! Dans le cas o\`u on active les thermiques, on fait partir l'ajustement
4579       ! a partir du sommet des thermiques.
4580       ! Dans le cas contraire, on demarre au niveau 1.
4581
4582       IF (iflag_thermals>=13.or.iflag_thermals<=0) THEN
4583
4584          IF (iflag_thermals.eq.0) THEN
4585             IF (prt_level>9) WRITE(lunout,*)'ajsec'
4586             limbas(:)=1
4587          ELSE
4588             limbas(:)=lmax_th(:)
4589          ENDIF
4590
4591          ! Attention : le call ajsec_convV2 n'est maintenu que momentanneement
4592          ! pour des test de convergence numerique.
4593          ! Le nouveau ajsec est a priori mieux, meme pour le cas
4594          ! iflag_thermals = 0 (l'ancienne version peut faire des tendances
4595          ! non nulles numeriquement pour des mailles non concernees.
4596
4597#ifdef ISO
4598#ifdef ISOVERIF
4599      call iso_verif_noNaN_vect2D(xt_seri,'physiq 4112', &
4600     &            ntraciso,klon,klev)
4601#endif         
4602#ifdef ISOVERIF
4603      write(*,*) 'physiq 3691c: avant call ajsec'
4604      if (iso_eau.gt.0) then
4605      call iso_verif_egalite_vect2D( &
4606     &           xt_seri,q_seri, &
4607     &           'physiq 3727, avant ajsec',ntraciso,klon,klev)
4608      endif
4609      if (iso_HDO.gt.0) then
4610      call iso_verif_aberrant_enc_vect2D( &
4611     &           xt_seri,q_seri, &
4612     &           'physiq 3732, avant ajsec',ntraciso,klon,klev)
4613      endif   
4614      if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then     
4615        do k = 1, klev
4616        do i = 1, klon
4617            if (q_seri(i,k).gt.ridicule) then 
4618               if (iso_verif_o18_aberrant_nostop( &
4619     &              xt_seri(iso_HDO,i,k)/q_seri(i,k), &
4620     &              xt_seri(iso_O18,i,k)/q_seri(i,k), &
4621     &              'physiq 4051, avant ajsec').eq.1) then
4622                  write(*,*) 'physiq 4367: i,k,q_seri(i,k)=',i,k,q_seri(i,k)                       
4623                  stop
4624              endif !  if (iso_verif_o18_aberrant_nostop
4625            endif !if (q_seri(i,k).gt.errmax) then   
4626        enddo !do i = 1, klon
4627        enddo !do k = 1, klev
4628        endif !if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then     
4629#ifdef ISOTRAC
4630        call iso_verif_traceur_vect(xt_seri,klon,klev, &
4631     &           'physiq 3600 avant ajsec')
4632#endif
4633        !#ifdef ISOTRAC
4634#endif
4635        !#ifdef ISOVERIF
4636#endif 
4637        !#ifdef ISO
4638
4639          IF (iflag_thermals==0) THEN
4640             ! Calling adjustment alone (but not the thermal plume model)
4641             CALL ajsec_convV2(paprs, pplay, t_seri,q_seri &
4642                  , d_t_ajsb, d_q_ajsb &
4643#ifdef ISO
4644     &      ,xt_seri,d_xt_ajsb      &
4645#endif         
4646     &   )
4647          ELSE IF (iflag_thermals>0) THEN
4648             ! Calling adjustment above the top of thermal plumes
4649             CALL ajsec(paprs, pplay, t_seri,q_seri,limbas &
4650                  , d_t_ajsb, d_q_ajsb &
4651#ifdef ISO
4652     &      ,xt_seri,d_xt_ajsb &     
4653#endif         
4654     &   )
4655          ENDIF
4656
4657          !--------------------------------------------------------------------
4658          ! ajout des tendances de l'ajustement sec ou des thermiques
4659          CALL add_phys_tend(du0,dv0,d_t_ajsb,d_q_ajsb,dql0,dqi0,dqbs0,paprs, &
4660               'ajsb',abortphy,flag_inhib_tend,itap,0 &
4661#ifdef ISO
4662     &    ,d_xt_ajsb,dxtl0,dxti0 &
4663#endif     
4664     &   )
4665          CALL prt_enerbil('ajsb',itap)
4666          d_t_ajs(:,:)=d_t_ajs(:,:)+d_t_ajsb(:,:)
4667          d_q_ajs(:,:)=d_q_ajs(:,:)+d_q_ajsb(:,:)
4668#ifdef ISO
4669         d_xt_ajs(:,:,:)=d_xt_ajs(:,:,:)+d_xt_ajsb(:,:,:)
4670#endif
4671
4672          !---------------------------------------------------------------------
4673
4674       ENDIF
4675
4676    ENDIF
4677    !
4678    !===================================================================
4679    ! Computation of ratqs, the width (normalized) of the subrid scale
4680    ! water distribution
4681
4682    tke_dissip_ave(:,:)=0.
4683    l_mix_ave(:,:)=0.
4684    wprime_ave(:,:)=0.
4685
4686    DO nsrf = 1, nbsrf
4687       DO i = 1, klon
4688          tke_dissip_ave(i,:) = tke_dissip_ave(i,:) + tke_dissip(i,:,nsrf)*pctsrf(i,nsrf)
4689          l_mix_ave(i,:) = l_mix_ave(i,:) + l_mix(i,:,nsrf)*pctsrf(i,nsrf)
4690          wprime_ave(i,:) = wprime_ave(i,:) + wprime(i,:,nsrf)*pctsrf(i,nsrf)
4691       ENDDO
4692    ENDDO
4693
4694    CALL ratqs_main(klon,klev,nbsrf,prt_level,lunout,        &
4695         iflag_ratqs,iflag_con,iflag_cld_th,pdtphys,  &
4696         ratqsbas,ratqshaut,ratqsp0, ratqsdp, &
4697         pctsrf,s_pblh,zstd, &
4698         tau_ratqs,fact_cldcon,wake_s, wake_deltaq,   &
4699         ptconv,ptconvth,clwcon0th, rnebcon0th,     &
4700         paprs,pplay,t_seri,q_seri, &
4701         qtc_cv, sigt_cv,detrain_cv,fm_cv,fqd,fqcomp,sigd,zqsat, &
4702         omega,pbl_tke(:,:,is_ave),tke_dissip_ave,l_mix_ave,wprime_ave, &
4703         t2m,q2m,fm_therm,entr_therm,detr_therm,cell_area, &
4704         ratqs,ratqsc,ratqs_inter_)
4705
4706    !
4707    ! Appeler le processus de condensation a grande echelle
4708    ! et le processus de precipitation
4709    !-------------------------------------------------------------------------
4710    IF (prt_level .GE.10) THEN
4711       print *,'itap, ->fisrtilp ',itap
4712    ENDIF
4713    !
4714#ifdef ISO   
4715#ifdef ISOVERIF
4716      write(*,*) 'physiq 3837: verif avant ilp'
4717#endif         
4718#ifdef ISOVERIF
4719      if (iso_eau.gt.0) then
4720      call iso_verif_egalite_vect2D( &
4721     &           xt_seri,q_seri, &
4722     &           'physiq 3709, avant ilp',ntraciso,klon,klev)
4723      endif
4724      if (iso_HDO.gt.0) then
4725      call iso_verif_aberrant_enc_vect2D( &
4726     &           xt_seri,q_seri, &
4727     &           'physiq 3714, avant ilp',ntraciso,klon,klev)
4728      endif
4729#endif         
4730#ifdef ISOVERIF
4731       do k=1,klev
4732        do i=1,klon
4733           if ((iso_O17.gt.0).and.(iso_O18.gt.0)) then
4734            if ((q_seri(i,k).gt.ridicule).and.(k.lt.nlevmaxO17)) then
4735             call iso_verif_aberrant_o17(xt_seri(iso_o17,i,k) &
4736     &           /q_seri(i,k),xt_seri(iso_o18,i,k) &
4737     &           /q_seri(i,k),'physiq 3389')
4738            endif !if (q_seri(i,k).gt.ridicule) then
4739           endif !if (iso_O17.gt.0) then
4740           if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then   
4741            if (q_seri(i,k).gt.ridicule) then 
4742               if (iso_verif_o18_aberrant_nostop( &
4743     &              xt_seri(iso_HDO,i,k)/q_seri(i,k), &
4744     &              xt_seri(iso_O18,i,k)/q_seri(i,k), &
4745     &              'physiq 4177, avant il pleut').eq.1) then
4746                  write(*,*) 'physic 4475: i,k,q_seri(i,k)=',i,k,q_seri(i,k)                       
4747                  stop
4748              endif !  if (iso_verif_o18_aberrant_nostop
4749            endif !if (q_seri(i,k).gt.errmax) then   
4750        endif !if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then 
4751#ifdef ISOTRAC     
4752        call iso_verif_traceur(xt_seri(1,i,k),'physiq 3609')
4753        call iso_verif_traceur_pbidouille(xt_seri(1,i,k), &
4754     &           'physiq 3609')
4755        if (iso_verif_tracpos_choix_nostop(xt_seri(1,i,k), &
4756     &           'physiq 3707, avant fisrtilp',1e-5) &
4757     &           .eq.1) then
4758              write(*,*) 'i,k=',i,k
4759#ifdef ISOVERIF
4760                  stop
4761#endif
4762         endif
4763         if (option_tmin.eq.1) then
4764               call iso_verif_trac17_q_deltaD(xt_seri(1,i,k), &
4765     &           'physiq 3905, avant ilp')
4766         endif
4767#endif
4768!ISOTRAC         
4769         enddo !do i=1,klon
4770        enddo !do k=1,klev       
4771
4772        ! verif température
4773        do k=1,klev
4774           do i=1,klon
4775             call iso_verif_positif(370.0-t_seri(i,k), &
4776     &          'physiq 3535, avant il pleut')
4777             call iso_verif_positif(t_seri(i,k)-100.0, &
4778     &          'physiq 3537, avant il pleut')
4779           enddo
4780        enddo
4781      write(*,*) 'physiq 3113: appel de fisrtilp'
4782#endif     
4783!ISOVERIF
4784          if ((bidouille_anti_divergence).and. &
4785     &           (iso_eau.gt.0)) then
4786           do k=1,klev
4787            do i=1,klon 
4788             xt_seri(iso_eau,i,k)= q_seri(i,k)
4789            enddo !do i=1,klon
4790           enddo !do k=1,klev
4791          endif !if ((bidouille_anti_divergence).and. &   
4792#endif
4793!ISO
4794
4795    picefra(:,:)=0.
4796
4797    IF (ok_new_lscp) THEN
4798
4799    !--mise à jour de flight_m et flight_h2o dans leur module
4800    IF (ok_plane_h2o .OR. ok_plane_contrail) THEN
4801      CALL airplane(debut,pphis,pplay,paprs,t_seri)
4802    ENDIF
4803
4804    CALL lscp(klon,klev,phys_tstep,missing_val,paprs,pplay, &
4805         t_seri, q_seri,ptconv,ratqs, &
4806         d_t_lsc, d_q_lsc, d_ql_lsc, d_qi_lsc, rneb, rneblsvol, rneb_seri, &
4807         pfraclr,pfracld, &
4808         radocond, picefra, rain_lsc, snow_lsc, &
4809         frac_impa, frac_nucl, beta_prec_fisrt, &
4810         prfl, psfl, rhcl,  &
4811         zqasc, fraca,ztv,zpspsk,ztla,zthl,iflag_cld_th, &
4812         iflag_ice_thermo, ok_ice_sursat, zqsatl, zqsats, distcltop, &
4813         qclr, qcld, qss, qvc, rnebclr, rnebss, gamma_ss, &
4814         Tcontr, qcontr, qcontr2, fcontrN, fcontrP )
4815
4816    ELSE
4817
4818    CALL fisrtilp(phys_tstep,paprs,pplay, &
4819         t_seri, q_seri,ptconv,ratqs, &
4820         d_t_lsc, d_q_lsc, d_ql_lsc, d_qi_lsc, rneb, radocond, &
4821         rain_lsc, snow_lsc, &
4822         pfrac_impa, pfrac_nucl, pfrac_1nucl, &
4823         frac_impa, frac_nucl, beta_prec_fisrt, &
4824         prfl, psfl, rhcl,  &
4825         zqasc, fraca,ztv,zpspsk,ztla,zthl,iflag_cld_th, &
4826         iflag_ice_thermo &
4827#ifdef ISO     
4828     &          ,xt_seri,xtrain_lsc,xtsnow_lsc &
4829     &          ,d_xt_lsc,d_xtl_lsc,d_xti_lsc,pxtrfl,pxtsfl &
4830#endif     
4831     &   )
4832    ENDIF
4833    !
4834    WHERE (rain_lsc < 0) rain_lsc = 0.
4835    WHERE (snow_lsc < 0) snow_lsc = 0.
4836#ifdef ISO
4837      WHERE (xtrain_lsc < 0) xtrain_lsc = 0.
4838      WHERE (xtsnow_lsc < 0) xtsnow_lsc = 0.
4839#endif
4840
4841#ifdef ISO
4842#ifdef ISOVERIF
4843      DO k = 1, klev
4844        do i=1,klon
4845        if (iso_O18.gt.0.and.iso_HDO.gt.0) then
4846        if (ql_seri(i,k).gt.ridicule) then
4847               call iso_verif_aberrant( &
4848     &           xtl_seri(iso_HDO,i,k)/ql_seri(i,k),'physiq 4330')
4849               if (iso_O18.gt.0) then 
4850                 call iso_verif_o18_aberrant( &
4851     &                  xtl_seri(iso_HDO,i,k)/ql_seri(i,k), &
4852     &                  xtl_seri(iso_O18,i,k)/ql_seri(i,k), &
4853     &                  'physiq 4335, avant il pleut')
4854               endif ! if (iso_O18.gt.0) then 
4855             endif !if (ql_seri(i,k).gt.errmax) then 
4856         endif !if (iso_HDO.gt.0) then
4857       ENDDO ! i=1,klon
4858      ENDDO !DO k = 1, klev
4859#endif
4860#endif
4861
4862!+JLD
4863!    write(*,9000) 'phys lsc',"enerbil: bil_q, bil_e,",rain_lsc+snow_lsc &
4864!        & ,((rcw-rcpd)*rain_lsc + (rcs-rcpd)*snow_lsc)*t_seri(1,1)-rlvtt*rain_lsc+rlstt*snow_lsc &
4865!        & ,rain_lsc,snow_lsc
4866!    write(*,9000) "rcpv","rcw",rcpv,rcw,rcs,t_seri(1,1)
4867!-JLD
4868    CALL add_phys_tend(du0,dv0,d_t_lsc,d_q_lsc,d_ql_lsc,d_qi_lsc,dqbs0,paprs, &
4869         'lsc',abortphy,flag_inhib_tend,itap,0 &
4870#ifdef ISO
4871     &    ,d_xt_lsc,d_xtl_lsc,d_xti_lsc &
4872#endif     
4873     &   )
4874    CALL prt_enerbil('lsc',itap)
4875    rain_num(:)=0.
4876    DO k = 1, klev
4877       DO i = 1, klon
4878          IF (ql_seri(i,k)>oliqmax) THEN
4879             rain_num(i)=rain_num(i)+(ql_seri(i,k)-oliqmax)*zmasse(i,k)/pdtphys
4880#ifdef ISO
4881             do ixt=1,ntraciso
4882                xtl_seri(ixt,i,k)=xtl_seri(ixt,i,k)/ql_seri(i,k)*oliqmax
4883             enddo
4884#endif
4885             ql_seri(i,k)=oliqmax
4886          ENDIF
4887       ENDDO
4888    ENDDO
4889    IF (nqo >= 3) THEN
4890    DO k = 1, klev
4891       DO i = 1, klon
4892          IF (qs_seri(i,k)>oicemax) THEN
4893             rain_num(i)=rain_num(i)+(qs_seri(i,k)-oicemax)*zmasse(i,k)/pdtphys
4894#ifdef ISO
4895             do ixt=1,ntraciso
4896                xts_seri(ixt,i,k)=xts_seri(ixt,i,k)/qs_seri(i,k)*oliqmax
4897             enddo
4898#endif
4899             qs_seri(i,k)=oicemax
4900          ENDIF
4901       ENDDO
4902    ENDDO
4903    ENDIF
4904
4905
4906!---------------------------------------------------------------------------
4907    DO k = 1, klev
4908       DO i = 1, klon
4909          cldfra(i,k) = rneb(i,k)
4910          !CR: a quoi ca sert? Faut-il ajouter qs_seri?
4911          !EV: en effet etrange, j'ajouterais aussi qs_seri
4912          !    plus largement, je nettoierais (enleverrais) ces lignes
4913          IF (.NOT.new_oliq) radocond(i,k) = ql_seri(i,k)
4914       ENDDO
4915    ENDDO
4916
4917
4918    ! Option to activate the radiative effect of blowing snow (ok_rad_bs)
4919    ! makes sense only if the new large scale condensation scheme is active
4920    ! with the ok_icefra_lscp flag active as well
4921
4922    IF (ok_bs .AND. ok_rad_bs) THEN
4923       IF (ok_new_lscp .AND. ok_icefra_lscp) THEN
4924           DO k=1,klev
4925             DO i=1,klon
4926                radocond(i,k)=radocond(i,k)+qbs_seri(i,k)
4927                picefra(i,k)=(radocond(i,k)*picefra(i,k)+qbs_seri(i,k))/(radocond(i,k))
4928                qbsfra=min(qbs_seri(i,k)/qbst_bs,1.0)
4929                cldfra(i,k)=max(cldfra(i,k),qbsfra)
4930             ENDDO
4931           ENDDO
4932       ELSE
4933          WRITE(lunout,*)"PAY ATTENTION, you try to activate the radiative effect of blowing snow"
4934          WRITE(lunout,*)"with ok_new_lscp=false and/or ok_icefra_lscp=false"
4935          abort_message='inconsistency in cloud phase for blowing snow'
4936          CALL abort_physic(modname,abort_message,1)
4937       ENDIF
4938
4939    ENDIF
4940#ifdef ISO     
4941!#ifdef ISOVERIF
4942!      write(*,*) 'physiq 4030: verif apres ajout lsc'
4943!#endif     
4944#ifdef ISOVERIF
4945      call iso_verif_noNaN_vect2D(xt_seri,'physiq 4565', &
4946     &            ntraciso,klon,klev)
4947#endif
4948#ifdef ISOVERIF
4949      if (iso_eau.gt.0) then
4950      call iso_verif_egalite_vect2D( &
4951     &           xt_seri,q_seri, &
4952     &           'physiq 2599, apres ajout lsc',ntraciso,klon,klev)
4953      endif
4954      if (iso_HDO.gt.0) then
4955      call iso_verif_aberrant_enc_vect2D( &
4956     &           xt_seri,q_seri, &
4957     &           'physiq 2866, apres ajout lsc',ntraciso,klon,klev)
4958      endif
4959#endif     
4960#ifdef ISOVERIF
4961      DO k = 1, klev
4962      DO i = 1, klon
4963          if (iso_eau.gt.0) then
4964           call iso_verif_egalite_choix( &
4965     &           xtl_seri(iso_eau,i,k),ql_seri(i,k), &
4966     &            'physiq 2601',errmax,errmaxrel)
4967          endif !if (iso_eau.gt.0) then 
4968          if (iso_HDO.gt.0) then   
4969            if (q_seri(i,k).gt.ridicule) then 
4970              if (iso_O18.gt.0) then 
4971               if (iso_verif_o18_aberrant_nostop( &
4972     &              xt_seri(iso_HDO,i,k)/q_seri(i,k), &
4973     &              xt_seri(iso_O18,i,k)/q_seri(i,k), &
4974     &              'physiq 2863, after il pleut').eq.1) then
4975                  write(*,*) 'physic 4656: i,k,q_seri(i,k)=',i,k,q_seri(i,k)                       
4976                  stop
4977                endif
4978              endif ! if (iso_O18.gt.0) then 
4979            endif !if (q_seri(i,k).gt.errmax) then
4980            if (ql_seri(i,k).gt.ridicule) then
4981               call iso_verif_aberrant_choix(xtl_seri(iso_HDO,i,k),ql_seri(i,k), &
4982                        ridicule,deltalim_snow,'physiq 2871')
4983               if (iso_O18.gt.0) then 
4984                 if (iso_verif_o18_aberrant_nostop( &
4985     &                  xtl_seri(iso_HDO,i,k)/ql_seri(i,k), &
4986     &                  xtl_seri(iso_O18,i,k)/ql_seri(i,k), &
4987     &                  'physiq 2872, after il pleut').eq.1) then
4988                        write(*,*) 'i,k,ql_seri(i,k)=',i,k,ql_seri(i,k)
4989                        write(*,*) 'd_ql_lsc(i,k)=',d_ql_lsc(i,k)
4990                        write(*,*) 'deltaD(d_ql_lsc(i,k))=',deltaD( &
4991     &                           d_xtl_lsc(iso_HDO,i,k)/d_ql_lsc(i,k))
4992                        write(*,*) 'deltaO(d_ql_lsc(i,k))=',deltaO( &
4993     &                           d_xtl_lsc(iso_O18,i,k)/d_ql_lsc(i,k))
4994                        !stop
4995                 endif
4996               endif ! if (iso_O18.gt.0) then 
4997             endif !if (ql_seri(i,k).gt.errmax) then 
4998         endif !if (iso_HDO.gt.0) then
4999       ENDDO ! i=1,klon
5000      ENDDO !DO k = 1, klev
5001      if ((iso_O17.gt.0).and.(iso_O18.gt.0)) then
5002        do i=1,klon
5003          do k=1,nlev
5004           if ((q_seri(i,k).gt.ridicule).and.(k.lt.nlevmaxO17)) then
5005            call iso_verif_aberrant_o17(xt_seri(iso_o17,i,k) &
5006     &           /q_seri(i,k),xt_seri(iso_o18,i,k) &
5007     &           /q_seri(i,k),'physiq 3549')
5008           endif !if (q_seri(i,k).gt.ridicule) then
5009          enddo !do k=1,nlev
5010        enddo  !do i=1,klon       
5011      endif !if (iso_O17.gt.0) then
5012#endif 
5013!#ifdef ISOVERIF
5014#ifdef ISOTRAC
5015#ifdef ISOVERIF       
5016      call iso_verif_traceur_vect(xt_seri,klon,klev,'physiq 3794')
5017      if (option_tmin.eq.1) then
5018      if (nzone_temp.ge.5) then
5019            call iso_verif_tag17_q_deltaD_vect(xt_seri,klon,klev, &
5020     &          'physiq 3385: apres ilp')
5021      endif ! if (nzone_temp.ge.5) then
5022      endif !if (option_tmin.eq.1) then
5023#endif   
5024!#ifdef ISOVERIF               
5025      if (bidouille_anti_divergence) then
5026          call iso_verif_traceur_jbid_vect(xt_seri, &
5027     &            klon,klev) 
5028      endif
5029
5030        ! si tag 22, on recolorise les zones très convectives dans les
5031        ! tropiques
5032        if (option_traceurs.eq.22) then
5033          call isotrac_recolorise_conv(xt_seri,latitude_deg,presnivs,rain_con)
5034        endif
5035
5036#endif
5037!#ifdef ISOTRAC     
5038#endif   
5039!#ifdef ISO   
5040
5041    IF (check) THEN
5042       za = qcheck(klon,klev,paprs,q_seri,ql_seri,cell_area)
5043       WRITE(lunout,*)"apresilp=", za
5044       zx_t = 0.0
5045       za = 0.0
5046       DO i = 1, klon
5047          za = za + cell_area(i)/REAL(klon)
5048          zx_t = zx_t + (rain_lsc(i) &
5049               + snow_lsc(i))*cell_area(i)/REAL(klon)
5050       ENDDO
5051       zx_t = zx_t/za*phys_tstep
5052       WRITE(lunout,*)"Precip=", zx_t
5053    ENDIF
5054
5055    IF (mydebug) THEN
5056       CALL writefield_phy('u_seri',u_seri,nbp_lev)
5057       CALL writefield_phy('v_seri',v_seri,nbp_lev)
5058       CALL writefield_phy('t_seri',t_seri,nbp_lev)
5059       CALL writefield_phy('q_seri',q_seri,nbp_lev)
5060    ENDIF
5061
5062    !
5063    !-------------------------------------------------------------------
5064    !  PRESCRIPTION DES NUAGES POUR LE RAYONNEMENT
5065    !-------------------------------------------------------------------
5066
5067    ! 1. NUAGES CONVECTIFS
5068    !
5069    !IM cf FH
5070    !     IF (iflag_cld_th.eq.-1) THEN ! seulement pour Tiedtke
5071    IF (iflag_cld_th.le.-1) THEN ! seulement pour Tiedtke
5072       snow_tiedtke=0.
5073       !     print*,'avant calcul de la pseudo precip '
5074       !     print*,'iflag_cld_th',iflag_cld_th
5075       IF (iflag_cld_th.eq.-1) THEN
5076          rain_tiedtke=rain_con
5077       ELSE
5078          !       print*,'calcul de la pseudo precip '
5079          rain_tiedtke=0.
5080          !         print*,'calcul de la pseudo precip 0'
5081          DO k=1,klev
5082             DO i=1,klon
5083                IF (d_q_con(i,k).lt.0.) THEN
5084                   rain_tiedtke(i)=rain_tiedtke(i)-d_q_con(i,k)/pdtphys &
5085                        *(paprs(i,k)-paprs(i,k+1))/rg
5086                ENDIF
5087             ENDDO
5088          ENDDO
5089       ENDIF
5090       !
5091       !     call dump2d(iim,jjm,rain_tiedtke(2:klon-1),'PSEUDO PRECIP ')
5092       !
5093
5094       ! Nuages diagnostiques pour Tiedtke
5095       CALL diagcld1(paprs,pplay, &
5096                                !IM cf FH. rain_con,snow_con,ibas_con,itop_con,
5097            rain_tiedtke,snow_tiedtke,ibas_con,itop_con, &
5098            diafra,dialiq)
5099       DO k = 1, klev
5100          DO i = 1, klon
5101             IF (diafra(i,k).GT.cldfra(i,k)) THEN
5102                radocond(i,k) = dialiq(i,k)
5103                cldfra(i,k) = diafra(i,k)
5104             ENDIF
5105          ENDDO
5106       ENDDO
5107
5108    ELSE IF (iflag_cld_th.ge.3) THEN
5109       !  On prend pour les nuages convectifs le max du calcul de la
5110       !  convection et du calcul du pas de temps precedent diminue d'un facteur
5111       !  facttemps
5112       facteur = pdtphys *facttemps
5113       DO k=1,klev
5114          DO i=1,klon
5115             rnebcon(i,k)=rnebcon(i,k)*facteur
5116             IF (rnebcon0(i,k)*clwcon0(i,k).GT.rnebcon(i,k)*clwcon(i,k)) THEN
5117                rnebcon(i,k)=rnebcon0(i,k)
5118                clwcon(i,k)=clwcon0(i,k)
5119             ENDIF
5120          ENDDO
5121       ENDDO
5122
5123       !   On prend la somme des fractions nuageuses et des contenus en eau
5124
5125       IF (iflag_cld_th>=5) THEN
5126
5127          DO k=1,klev
5128             ptconvth(:,k)=fm_therm(:,k+1)>0.
5129          ENDDO
5130
5131          IF (iflag_coupl==4) THEN
5132
5133             ! Dans le cas iflag_coupl==4, on prend la somme des convertures
5134             ! convectives et lsc dans la partie des thermiques
5135             ! Le controle par iflag_coupl est peut etre provisoire.
5136             DO k=1,klev
5137                DO i=1,klon
5138                   IF (ptconv(i,k).AND.ptconvth(i,k)) THEN
5139                      radocond(i,k)=radocond(i,k)+rnebcon(i,k)*clwcon(i,k)
5140                      cldfra(i,k)=min(cldfra(i,k)+rnebcon(i,k),1.)
5141                   ELSE IF (ptconv(i,k)) THEN
5142                      cldfra(i,k)=rnebcon(i,k)
5143                      radocond(i,k)=rnebcon(i,k)*clwcon(i,k)
5144                   ENDIF
5145                ENDDO
5146             ENDDO
5147
5148          ELSE IF (iflag_coupl==5) THEN
5149             DO k=1,klev
5150                DO i=1,klon
5151                   cldfra(i,k)=min(cldfra(i,k)+rnebcon(i,k),1.)
5152                   radocond(i,k)=radocond(i,k)+rnebcon(i,k)*clwcon(i,k)
5153                ENDDO
5154             ENDDO
5155
5156          ELSE
5157
5158             ! Si on est sur un point touche par la convection
5159             ! profonde et pas par les thermiques, on prend la
5160             ! couverture nuageuse et l'eau nuageuse de la convection
5161             ! profonde.
5162
5163             !IM/FH: 2011/02/23
5164             ! definition des points sur lesquels ls thermiques sont actifs
5165
5166             DO k=1,klev
5167                DO i=1,klon
5168                   IF (ptconv(i,k).AND. .NOT.ptconvth(i,k)) THEN
5169                      cldfra(i,k)=rnebcon(i,k)
5170                      radocond(i,k)=rnebcon(i,k)*clwcon(i,k)
5171                   ENDIF
5172                ENDDO
5173             ENDDO
5174
5175          ENDIF
5176
5177       ELSE
5178
5179          ! Ancienne version
5180          cldfra(:,:)=min(max(cldfra(:,:),rnebcon(:,:)),1.)
5181          radocond(:,:)=radocond(:,:)+rnebcon(:,:)*clwcon(:,:)
5182       ENDIF
5183
5184    ENDIF
5185
5186    !     plulsc(:)=0.
5187    !     do k=1,klev,-1
5188    !        do i=1,klon
5189    !              zzz=prfl(:,k)+psfl(:,k)
5190    !           if (.not.ptconvth.zzz.gt.0.)
5191    !        enddo prfl, psfl,
5192    !     enddo
5193    !
5194    ! 2. NUAGES STARTIFORMES
5195    !
5196    IF (ok_stratus) THEN
5197       CALL diagcld2(paprs,pplay,t_seri,q_seri, diafra,dialiq)
5198       DO k = 1, klev
5199          DO i = 1, klon
5200             IF (diafra(i,k).GT.cldfra(i,k)) THEN
5201                radocond(i,k) = dialiq(i,k)
5202                cldfra(i,k) = diafra(i,k)
5203             ENDIF
5204          ENDDO
5205       ENDDO
5206    ENDIF
5207    !
5208    ! Precipitation totale
5209    !
5210    DO i = 1, klon
5211       rain_fall(i) = rain_con(i) + rain_lsc(i)
5212       snow_fall(i) = snow_con(i) + snow_lsc(i)
5213    ENDDO
5214
5215
5216
5217#ifdef ISO
5218      DO i = 1, klon
5219         do ixt=1,ntraciso
5220           xtrain_fall(ixt,i)=xtrain_con(ixt,i) + xtrain_lsc(ixt,i)
5221           xtsnow_fall(ixt,i)=xtsnow_con(ixt,i) + xtsnow_lsc(ixt,i)
5222         enddo
5223      ENDDO !DO i = 1, klon
5224#ifdef ISOVERIF
5225      if (iso_eau.gt.0) then
5226        DO i = 1, klon         
5227          call iso_verif_egalite_choix(xtsnow_fall(iso_eau,i), &
5228     &           snow_fall(i),'physiq 3387', &
5229     &           errmax,errmaxrel)
5230          if (iso_verif_positif_nostop(xtrain_lsc(iso_eau,i), &
5231     &           'physiq 4298').eq.1) then
5232            if (rain_lsc(i).ge.0.0) then
5233              write(*,*) 'rain_fall(i)=',rain_lsc(i)
5234              stop
5235            endif !if (rain_con(i).ge.0.0) then
5236          endif  !if (iso_verif_positif_nostop
5237          if (iso_verif_positif_nostop(xtrain_fall(iso_eau,i), &
5238     &           'physiq 3405').eq.1) then
5239            if (rain_fall(i).ge.0.0) then
5240              write(*,*) 'rain_fall(i)=',rain_fall(i)
5241              stop
5242            endif !if (rain_con(i).ge.0.0) then
5243          endif  !if (iso_verif_positif_nostop
5244        ENDDO !DO i = 1, klon   
5245      endif !if (iso_eau.gt.0) then
5246      if (iso_HDO.gt.0) then
5247          DO i = 1, klon
5248          call iso_verif_aberrant_choix(xtsnow_fall(iso_hdo,i), &
5249     &           snow_fall(i),ridicule_snow, &
5250     &           deltalim_snow,'physiq 3856')
5251          enddo !DO i = 1, klon
5252      endif
5253#ifdef ISOTRAC 
5254        DO i = 1, klon   
5255          call iso_verif_traceur(xtrain_fall(1,i),'physiq 4012')
5256          call iso_verif_traceur(xtsnow_fall(1,i),'physiq 4013')
5257        enddo
5258#endif     
5259#endif
5260#ifdef ISOVERIF   
5261       do i=1,klon   
5262         do ixt=1,ntraciso           
5263           call iso_verif_noNaN(xtsnow_con(ixt,i), &
5264     &             'physiq 4942')   
5265           call iso_verif_noNaN(xtsnow_lsc(ixt,i), &
5266     &             'physiq 4942')
5267           call iso_verif_noNaN(xtsnow_fall(ixt,i), &
5268     &             'physiq 4942')
5269         enddo  !do ixt=1,ntraciso
5270       enddo ! do klon
5271#endif     
5272#endif
5273#ifdef ISO
5274        if ((iso_eau.gt.0).and.(bidouille_anti_divergence)) then
5275        do i=1,klon
5276            xtrain_fall(iso_eau,i)=rain_fall(i)
5277        enddo !do i=1,klon
5278        endif
5279#endif   
5280    !
5281    ! Calculer l'humidite relative pour diagnostique
5282    !
5283    DO k = 1, klev
5284       DO i = 1, klon
5285          zx_t = t_seri(i,k)
5286          IF (thermcep) THEN
5287             !!           if (iflag_ice_thermo.eq.0) then                 !jyg
5288             zdelta = MAX(0.,SIGN(1.,rtt-zx_t))
5289             !!           else                                            !jyg
5290             !!           zdelta = MAX(0.,SIGN(1.,t_glace_min-zx_t))      !jyg
5291             !!           endif                                           !jyg
5292             zx_qs  = r2es * FOEEW(zx_t,zdelta)/pplay(i,k)
5293             zx_qs  = MIN(0.5,zx_qs)
5294             zcor   = 1./(1.-retv*zx_qs)
5295             zx_qs  = zx_qs*zcor
5296          ELSE
5297             !!           IF (zx_t.LT.t_coup) THEN             !jyg
5298             IF (zx_t.LT.rtt) THEN                  !jyg
5299                zx_qs = qsats(zx_t)/pplay(i,k)
5300             ELSE
5301                zx_qs = qsatl(zx_t)/pplay(i,k)
5302             ENDIF
5303          ENDIF
5304          zx_rh(i,k) = q_seri(i,k)/zx_qs
5305            IF (iflag_ice_thermo .GT. 0) THEN
5306          zx_rhl(i,k) = q_seri(i,k)/(qsatl(zx_t)/pplay(i,k))
5307          zx_rhi(i,k) = q_seri(i,k)/(qsats(zx_t)/pplay(i,k))
5308            ENDIF
5309          zqsat(i,k)=zx_qs
5310       ENDDO
5311    ENDDO
5312
5313    !IM Calcul temp.potentielle a 2m (tpot) et temp. potentielle
5314    !   equivalente a 2m (tpote) pour diagnostique
5315    !
5316    DO i = 1, klon
5317       tpot(i)=zt2m(i)*(100000./paprs(i,1))**RKAPPA
5318       IF (thermcep) THEN
5319          IF(zt2m(i).LT.RTT) then
5320             Lheat=RLSTT
5321          ELSE
5322             Lheat=RLVTT
5323          ENDIF
5324       ELSE
5325          IF (zt2m(i).LT.RTT) THEN
5326             Lheat=RLSTT
5327          ELSE
5328             Lheat=RLVTT
5329          ENDIF
5330       ENDIF
5331       tpote(i) = tpot(i)*      &
5332            EXP((Lheat *qsat2m(i))/(RCPD*zt2m(i)))
5333    ENDDO
5334
5335    IF (ANY(type_trac == ['inca','inco'])) THEN ! ModThL
5336#ifdef INCA
5337       CALL VTe(VTphysiq)
5338       CALL VTb(VTinca)
5339       calday = REAL(days_elapsed + 1) + jH_cur
5340
5341       CALL chemtime(itap+itau_phy-1, date0, phys_tstep, itap)
5342       CALL AEROSOL_METEO_CALC( &
5343            calday,pdtphys,pplay,paprs,t,pmflxr,pmflxs, &
5344            prfl,psfl,pctsrf,cell_area, &
5345            latitude_deg,longitude_deg,u10m,v10m)
5346
5347       zxsnow_dummy(:) = 0.0
5348
5349       CALL chemhook_begin (calday, &
5350            days_elapsed+1, &
5351            jH_cur, &
5352            pctsrf(1,1), &
5353            latitude_deg, &
5354            longitude_deg, &
5355            cell_area, &
5356            paprs, &
5357            pplay, &
5358            coefh(1:klon,1:klev,is_ave), &
5359            pphi, &
5360            t_seri, &
5361            u, &
5362            v, &
5363            rot, &
5364            wo(:, :, 1), &
5365            q_seri, &
5366            zxtsol, &
5367            zt2m, &
5368            zxsnow_dummy, &
5369            solsw, &
5370            albsol1, &
5371            rain_fall, &
5372            snow_fall, &
5373            itop_con, &
5374            ibas_con, &
5375            cldfra, &
5376            nbp_lon, &
5377            nbp_lat-1, &
5378            tr_seri(:,:,1+nqCO2:nbtr), &
5379            ftsol, &
5380            paprs, &
5381            cdragh, &
5382            cdragm, &
5383            pctsrf, &
5384            pdtphys, &
5385            itap)
5386
5387       CALL VTe(VTinca)
5388       CALL VTb(VTphysiq)
5389#endif
5390    ENDIF !type_trac = inca or inco
5391    IF (type_trac == 'repr') THEN
5392#ifdef REPROBUS
5393    !CALL chemtime_rep(itap+itau_phy-1, date0, dtime, itap)
5394    CALL chemtime_rep(itap+itau_phy-1, date0, phys_tstep, itap)
5395#endif
5396    ENDIF
5397
5398    !
5399    ! Appeler le rayonnement mais calculer tout d'abord l'albedo du sol.
5400    !
5401    IF (MOD(itaprad,radpas).EQ.0) THEN
5402
5403       !
5404       !jq - introduce the aerosol direct and first indirect radiative forcings
5405       !jq - Johannes Quaas, 27/11/2003 (quaas@lmd.jussieu.fr)
5406       IF (flag_aerosol .GT. 0) THEN
5407          IF (iflag_rrtm .EQ. 0) THEN !--old radiation
5408             IF (.NOT. aerosol_couple) THEN
5409                !
5410                CALL readaerosol_optic( &
5411                     debut, flag_aerosol, itap, jD_cur-jD_ref, &
5412                     pdtphys, pplay, paprs, t_seri, rhcl, presnivs,  &
5413                     mass_solu_aero, mass_solu_aero_pi,  &
5414                     tau_aero, piz_aero, cg_aero,  &
5415                     tausum_aero, tau3d_aero)
5416             ENDIF
5417          ELSE                       ! RRTM radiation
5418             IF (aerosol_couple .AND. config_inca == 'aero' ) THEN
5419                abort_message='config_inca=aero et rrtm=1 impossible'
5420                CALL abort_physic(modname,abort_message,1)
5421             ELSE
5422                !
5423#ifdef CPP_RRTM
5424                IF (NSW.EQ.6) THEN
5425                   !--new aerosol properties SW and LW
5426                   !
5427#ifdef CPP_Dust
5428                   !--SPL aerosol model
5429                   CALL splaerosol_optic_rrtm( ok_alw, pplay, paprs, t_seri, rhcl, &
5430                        tr_seri, mass_solu_aero, mass_solu_aero_pi,  &
5431                        tau_aero_sw_rrtm, piz_aero_sw_rrtm, cg_aero_sw_rrtm,  &
5432                        tausum_aero, tau3d_aero)
5433#else
5434                   !--climatologies or INCA aerosols
5435                   CALL readaerosol_optic_rrtm( debut, aerosol_couple, ok_alw, ok_volcan, &
5436                        flag_aerosol, flag_bc_internal_mixture, itap, jD_cur-jD_ref, &
5437                        pdtphys, pplay, paprs, t_seri, rhcl, presnivs,  &
5438                        tr_seri, mass_solu_aero, mass_solu_aero_pi,  &
5439                        tau_aero_sw_rrtm, piz_aero_sw_rrtm, cg_aero_sw_rrtm,  &
5440                        tausum_aero, drytausum_aero, tau3d_aero)
5441#endif
5442
5443                   IF (flag_aerosol .EQ. 7) THEN
5444                      CALL MACv2SP(pphis,pplay,paprs,longitude_deg,latitude_deg,  &
5445                                   tau_aero_sw_rrtm,piz_aero_sw_rrtm,cg_aero_sw_rrtm)
5446                   ENDIF
5447
5448                   !
5449                ELSE IF (NSW.EQ.2) THEN
5450                   !--for now we use the old aerosol properties
5451                   !
5452                   CALL readaerosol_optic( &
5453                        debut, flag_aerosol, itap, jD_cur-jD_ref, &
5454                        pdtphys, pplay, paprs, t_seri, rhcl, presnivs,  &
5455                        mass_solu_aero, mass_solu_aero_pi,  &
5456                        tau_aero, piz_aero, cg_aero,  &
5457                        tausum_aero, tau3d_aero)
5458                   !
5459                   !--natural aerosols
5460                   tau_aero_sw_rrtm(:,:,1,:)=tau_aero(:,:,3,:)
5461                   piz_aero_sw_rrtm(:,:,1,:)=piz_aero(:,:,3,:)
5462                   cg_aero_sw_rrtm (:,:,1,:)=cg_aero (:,:,3,:)
5463                   !--all aerosols
5464                   tau_aero_sw_rrtm(:,:,2,:)=tau_aero(:,:,2,:)
5465                   piz_aero_sw_rrtm(:,:,2,:)=piz_aero(:,:,2,:)
5466                   cg_aero_sw_rrtm (:,:,2,:)=cg_aero (:,:,2,:)
5467                   !
5468                   !--no LW optics
5469                   tau_aero_lw_rrtm(:,:,:,:) = 1.e-15
5470                   !
5471                ELSE
5472                   abort_message='Only NSW=2 or 6 are possible with ' &
5473                        // 'aerosols and iflag_rrtm=1'
5474                   CALL abort_physic(modname,abort_message,1)
5475                ENDIF
5476#else
5477                abort_message='You should compile with -rrtm if running ' &
5478                     // 'with iflag_rrtm=1'
5479                CALL abort_physic(modname,abort_message,1)
5480#endif
5481                !
5482             ENDIF
5483          ENDIF
5484       ELSE   !--flag_aerosol = 0
5485          tausum_aero(:,:,:) = 0.
5486          drytausum_aero(:,:) = 0.
5487          mass_solu_aero(:,:) = 0.
5488          mass_solu_aero_pi(:,:) = 0.
5489          IF (iflag_rrtm .EQ. 0) THEN !--old radiation
5490             tau_aero(:,:,:,:) = 1.e-15
5491             piz_aero(:,:,:,:) = 1.
5492             cg_aero(:,:,:,:)  = 0.
5493          ELSE
5494             tau_aero_sw_rrtm(:,:,:,:) = 1.e-15
5495             tau_aero_lw_rrtm(:,:,:,:) = 1.e-15
5496             piz_aero_sw_rrtm(:,:,:,:) = 1.0
5497             cg_aero_sw_rrtm(:,:,:,:)  = 0.0
5498          ENDIF
5499       ENDIF
5500       !
5501       !--WMO criterion to determine tropopause
5502       CALL stratosphere_mask(missing_val, pphis, t_seri, pplay, latitude_deg)
5503       !
5504       !--STRAT AEROSOL
5505       !--updates tausum_aero,tau_aero,piz_aero,cg_aero
5506       IF (flag_aerosol_strat.GT.0) THEN
5507          IF (prt_level .GE.10) THEN
5508             PRINT *,'appel a readaerosolstrat', mth_cur
5509          ENDIF
5510          IF (iflag_rrtm.EQ.0) THEN
5511           IF (flag_aerosol_strat.EQ.1) THEN
5512             CALL readaerosolstrato(debut)
5513           ELSE
5514             abort_message='flag_aerosol_strat must equal 1 for rrtm=0'
5515             CALL abort_physic(modname,abort_message,1)
5516           ENDIF
5517          ELSE
5518#ifdef CPP_RRTM
5519#ifndef CPP_StratAer
5520          !--prescribed strat aerosols
5521          !--only in the case of non-interactive strat aerosols
5522            IF (flag_aerosol_strat.EQ.1) THEN
5523             CALL readaerosolstrato1_rrtm(debut)
5524            ELSEIF (flag_aerosol_strat.EQ.2) THEN
5525             CALL readaerosolstrato2_rrtm(debut, ok_volcan)
5526            ELSE
5527             abort_message='flag_aerosol_strat must equal 1 or 2 for rrtm=1'
5528             CALL abort_physic(modname,abort_message,1)
5529            ENDIF
5530#endif
5531#else
5532             abort_message='You should compile with -rrtm if running ' &
5533                  // 'with iflag_rrtm=1'
5534             CALL abort_physic(modname,abort_message,1)
5535#endif
5536          ENDIF
5537       ELSE
5538          tausum_aero(:,:,id_STRAT_phy) = 0.
5539       ENDIF
5540!
5541#ifdef CPP_RRTM
5542#ifdef CPP_StratAer
5543       !--compute stratospheric mask
5544       CALL stratosphere_mask(missing_val, pphis, t_seri, pplay, latitude_deg)
5545       !--interactive strat aerosols
5546       CALL calcaerosolstrato_rrtm(pplay,t_seri,paprs,debut)
5547#endif
5548#endif
5549       !--fin STRAT AEROSOL
5550       !     
5551
5552       ! Calculer les parametres optiques des nuages et quelques
5553       ! parametres pour diagnostiques:
5554       !
5555       IF (aerosol_couple.AND.config_inca=='aero') THEN
5556          mass_solu_aero(:,:)    = ccm(:,:,1)
5557          mass_solu_aero_pi(:,:) = ccm(:,:,2)
5558       ENDIF
5559
5560       IF (ok_newmicro) then
5561! AI          IF (iflag_rrtm.NE.0) THEN
5562          IF (iflag_rrtm.EQ.1) THEN
5563#ifdef CPP_RRTM
5564             IF (ok_cdnc.AND.NRADLP.NE.3) THEN
5565             abort_message='RRTM choix incoherent NRADLP doit etre egal a 3 ' &
5566                  // 'pour ok_cdnc'
5567             CALL abort_physic(modname,abort_message,1)
5568             ENDIF
5569#else
5570
5571             abort_message='You should compile with -rrtm if running with '//'iflag_rrtm=1'
5572             CALL abort_physic(modname,abort_message,1)
5573#endif
5574          ENDIF
5575          CALL newmicro (flag_aerosol, ok_cdnc, bl95_b0, bl95_b1, &
5576               paprs, pplay, t_seri, radocond, picefra, cldfra, &
5577               cldtau, cldemi, cldh, cldl, cldm, cldt, cldq, &
5578               flwp, fiwp, flwc, fiwc, &
5579               mass_solu_aero, mass_solu_aero_pi, &
5580               cldtaupi, latitude_deg, distcltop, re, fl, ref_liq, ref_ice, &
5581               ref_liq_pi, ref_ice_pi)
5582       ELSE
5583          CALL nuage (paprs, pplay, &
5584               t_seri, radocond, picefra, cldfra, cldtau, cldemi, &
5585               cldh, cldl, cldm, cldt, cldq, &
5586               ok_aie, &
5587               mass_solu_aero, mass_solu_aero_pi, &
5588               bl95_b0, bl95_b1, distcltop, &
5589               cldtaupi, re, fl)
5590       ENDIF
5591       !
5592       !IM betaCRF
5593       !
5594       cldtaurad   = cldtau
5595       cldtaupirad = cldtaupi
5596       cldemirad   = cldemi
5597       cldfrarad   = cldfra
5598
5599       !
5600       IF (lon1_beta.EQ.-180..AND.lon2_beta.EQ.180..AND. &
5601           lat1_beta.EQ.90..AND.lat2_beta.EQ.-90.) THEN
5602          !
5603          ! global
5604          !
5605!IM 251017 begin
5606!               print*,'physiq betaCRF global zdtime=',zdtime
5607!IM 251017 end
5608          DO k=1, klev
5609             DO i=1, klon
5610                IF (pplay(i,k).GE.pfree) THEN
5611                   beta(i,k) = beta_pbl
5612                ELSE
5613                   beta(i,k) = beta_free
5614                ENDIF
5615                IF (mskocean_beta) THEN
5616                   beta(i,k) = beta(i,k) * pctsrf(i,is_oce)
5617                ENDIF
5618                cldtaurad(i,k)   = cldtau(i,k) * beta(i,k)
5619                cldtaupirad(i,k) = cldtaupi(i,k) * beta(i,k)
5620                cldemirad(i,k)   = cldemi(i,k) * beta(i,k)
5621                cldfrarad(i,k)   = cldfra(i,k) * beta(i,k)
5622             ENDDO
5623          ENDDO
5624          !
5625       ELSE
5626          !
5627          ! regional
5628          !
5629          DO k=1, klev
5630             DO i=1,klon
5631                !
5632                IF (longitude_deg(i).ge.lon1_beta.AND. &
5633                    longitude_deg(i).le.lon2_beta.AND. &
5634                    latitude_deg(i).le.lat1_beta.AND.  &
5635                    latitude_deg(i).ge.lat2_beta) THEN
5636                   IF (pplay(i,k).GE.pfree) THEN
5637                      beta(i,k) = beta_pbl
5638                   ELSE
5639                      beta(i,k) = beta_free
5640                   ENDIF
5641                   IF (mskocean_beta) THEN
5642                      beta(i,k) = beta(i,k) * pctsrf(i,is_oce)
5643                   ENDIF
5644                   cldtaurad(i,k)   = cldtau(i,k) * beta(i,k)
5645                   cldtaupirad(i,k) = cldtaupi(i,k) * beta(i,k)
5646                   cldemirad(i,k)   = cldemi(i,k) * beta(i,k)
5647                   cldfrarad(i,k)   = cldfra(i,k) * beta(i,k)
5648                ENDIF
5649             !
5650             ENDDO
5651          ENDDO
5652       !
5653       ENDIF
5654
5655       !lecture de la chlorophylle pour le nouvel albedo de Sunghye Baek
5656       IF (ok_chlorophyll) THEN
5657          print*,"-- reading chlorophyll"
5658          CALL readchlorophyll(debut)
5659       ENDIF
5660
5661!--if ok_suntime_rrtm we use ancillay data for RSUN
5662!--previous values are therefore overwritten
5663!--this is needed for CMIP6 runs
5664!--and only possible for new radiation scheme
5665       IF (iflag_rrtm.EQ.1.AND.ok_suntime_rrtm) THEN
5666#ifdef CPP_RRTM
5667         CALL read_rsun_rrtm(debut)
5668#endif
5669       ENDIF
5670
5671       IF (mydebug) THEN
5672          CALL writefield_phy('u_seri',u_seri,nbp_lev)
5673          CALL writefield_phy('v_seri',v_seri,nbp_lev)
5674          CALL writefield_phy('t_seri',t_seri,nbp_lev)
5675          CALL writefield_phy('q_seri',q_seri,nbp_lev)
5676       ENDIF
5677
5678       !
5679       !sonia : If Iflag_radia >=2, pertubation of some variables
5680       !input to radiation (DICE)
5681       !
5682       IF (iflag_radia .ge. 2) THEN
5683          zsav_tsol (:) = zxtsol(:)
5684          CALL perturb_radlwsw(zxtsol,iflag_radia)
5685       ENDIF
5686
5687       IF (aerosol_couple.AND.config_inca=='aero') THEN
5688#ifdef INCA
5689          CALL radlwsw_inca  &
5690               (chemistry_couple, kdlon,kflev,dist, rmu0, fract, solaire, &
5691               paprs, pplay,zxtsol,albsol1, albsol2, t_seri,q_seri, &
5692               size(wo,3), wo, &
5693               cldfrarad, cldemirad, cldtaurad, &
5694               heat,heat0,cool,cool0,albpla, &
5695               topsw,toplw,solsw,sollw, &
5696               sollwdown, &
5697               topsw0,toplw0,solsw0,sollw0, &
5698               lwdn0, lwdn, lwup0, lwup,  &
5699               swdn0, swdn, swup0, swup, &
5700               ok_ade, ok_aie, &
5701               tau_aero, piz_aero, cg_aero, &
5702               topswad_aero, solswad_aero, &
5703               topswad0_aero, solswad0_aero, &
5704               topsw_aero, topsw0_aero, &
5705               solsw_aero, solsw0_aero, &
5706               cldtaupirad, &
5707               topswai_aero, solswai_aero)
5708#endif
5709       ELSE
5710          !
5711          !IM calcul radiatif pour le cas actuel
5712          !
5713          RCO2 = RCO2_act
5714          RCH4 = RCH4_act
5715          RN2O = RN2O_act
5716          RCFC11 = RCFC11_act
5717          RCFC12 = RCFC12_act
5718          !
5719          !--interactive CO2 in ppm from carbon cycle
5720          IF (carbon_cycle_rad) RCO2=RCO2_glo
5721          !
5722          IF (prt_level .GE.10) THEN
5723             print *,' ->radlwsw, number 1 '
5724          ENDIF
5725          !
5726          CALL radlwsw &
5727               (dist, rmu0, fract,  &
5728                                !albedo SB >>>
5729                                !      paprs, pplay,zxtsol,albsol1, albsol2,  &
5730               paprs, pplay,zxtsol,SFRWL,albsol_dir, albsol_dif,  &
5731                                !albedo SB <<<
5732               t_seri,q_seri,wo, &
5733               cldfrarad, cldemirad, cldtaurad, &
5734               ok_ade.OR.flag_aerosol_strat.GT.0, ok_aie,  ok_volcan, flag_volc_surfstrat, &
5735               flag_aerosol, flag_aerosol_strat, flag_aer_feedback, &
5736               tau_aero, piz_aero, cg_aero, &
5737               tau_aero_sw_rrtm, piz_aero_sw_rrtm, cg_aero_sw_rrtm, &
5738               ! Rajoute par OB pour RRTM
5739               tau_aero_lw_rrtm, &
5740               cldtaupirad, &
5741!              zqsat, flwcrad, fiwcrad, &
5742               zqsat, flwc, fiwc, &
5743               ref_liq, ref_ice, ref_liq_pi, ref_ice_pi, &
5744               heat,heat0,cool,cool0,albpla, &
5745               heat_volc,cool_volc, &
5746               topsw,toplw,solsw,solswfdiff,sollw, &
5747               sollwdown, &
5748               topsw0,toplw0,solsw0,sollw0, &
5749               lwdnc0, lwdn0, lwdn, lwupc0, lwup0, lwup,  &
5750               swdnc0, swdn0, swdn, swupc0, swup0, swup, &
5751               topswad_aero, solswad_aero, &
5752               topswai_aero, solswai_aero, &
5753               topswad0_aero, solswad0_aero, &
5754               topsw_aero, topsw0_aero, &
5755               solsw_aero, solsw0_aero, &
5756               topswcf_aero, solswcf_aero, &
5757                                !-C. Kleinschmitt for LW diagnostics
5758               toplwad_aero, sollwad_aero,&
5759               toplwai_aero, sollwai_aero, &
5760               toplwad0_aero, sollwad0_aero,&
5761                                !-end
5762               ZLWFT0_i, ZFLDN0, ZFLUP0, &
5763               ZSWFT0_i, ZFSDN0, ZFSUP0)
5764
5765          !lwoff=y, betalwoff=1. : offset LW CRE for radiation code and other
5766          !schemes
5767          toplw = toplw + betalwoff * (toplw0 - toplw)
5768          sollw = sollw + betalwoff * (sollw0 - sollw)
5769          lwdn = lwdn + betalwoff * (lwdn0 - lwdn)
5770          lwup = lwup + betalwoff * (lwup0 - lwup)
5771          sollwdown(:)= sollwdown(:) + betalwoff *(-1.*ZFLDN0(:,1) - &
5772                        sollwdown(:))
5773          cool = cool + betalwoff * (cool0 - cool)
5774 
5775#ifndef CPP_XIOS
5776          !
5777          !IM 2eme calcul radiatif pour le cas perturbe ou au moins un
5778          !IM des taux doit etre different du taux actuel
5779          !IM Par defaut on a les taux perturbes egaux aux taux actuels
5780          !
5781          IF (RCO2_per.NE.RCO2_act.OR. &
5782              RCH4_per.NE.RCH4_act.OR. &
5783              RN2O_per.NE.RN2O_act.OR. &
5784              RCFC11_per.NE.RCFC11_act.OR. &
5785              RCFC12_per.NE.RCFC12_act) ok_4xCO2atm =.TRUE.
5786#endif
5787   !
5788          IF (ok_4xCO2atm) THEN
5789                !
5790                RCO2 = RCO2_per
5791                RCH4 = RCH4_per
5792                RN2O = RN2O_per
5793                RCFC11 = RCFC11_per
5794                RCFC12 = RCFC12_per
5795                !
5796                IF (prt_level .GE.10) THEN
5797                   print *,' ->radlwsw, number 2 '
5798                ENDIF
5799                !
5800                CALL radlwsw &
5801                     (dist, rmu0, fract,  &
5802                                !albedo SB >>>
5803                                !      paprs, pplay,zxtsol,albsol1, albsol2,  &
5804                     paprs, pplay,zxtsol,SFRWL,albsol_dir, albsol_dif, &
5805                                !albedo SB <<<
5806                     t_seri,q_seri,wo, &
5807                     cldfrarad, cldemirad, cldtaurad, &
5808                     ok_ade.OR.flag_aerosol_strat.GT.0, ok_aie,  ok_volcan, flag_volc_surfstrat, &
5809                     flag_aerosol, flag_aerosol_strat, flag_aer_feedback, &
5810                     tau_aero, piz_aero, cg_aero, &
5811                     tau_aero_sw_rrtm, piz_aero_sw_rrtm, cg_aero_sw_rrtm, &
5812                                ! Rajoute par OB pour RRTM
5813                     tau_aero_lw_rrtm, &
5814                     cldtaupi, &
5815!                    zqsat, flwcrad, fiwcrad, &
5816                     zqsat, flwc, fiwc, &
5817                     ref_liq, ref_ice, ref_liq_pi, ref_ice_pi, &
5818                     heatp,heat0p,coolp,cool0p,albplap, &
5819                     heat_volc,cool_volc, &
5820                     topswp,toplwp,solswp,solswfdiffp,sollwp, &
5821                     sollwdownp, &
5822                     topsw0p,toplw0p,solsw0p,sollw0p, &
5823                     lwdnc0p, lwdn0p, lwdnp, lwupc0p, lwup0p, lwupp,  &
5824                     swdnc0p, swdn0p, swdnp, swupc0p, swup0p, swupp, &
5825                     topswad_aerop, solswad_aerop, &
5826                     topswai_aerop, solswai_aerop, &
5827                     topswad0_aerop, solswad0_aerop, &
5828                     topsw_aerop, topsw0_aerop, &
5829                     solsw_aerop, solsw0_aerop, &
5830                     topswcf_aerop, solswcf_aerop, &
5831                                !-C. Kleinschmitt for LW diagnostics
5832                     toplwad_aerop, sollwad_aerop,&
5833                     toplwai_aerop, sollwai_aerop, &
5834                     toplwad0_aerop, sollwad0_aerop,&
5835                                !-end
5836                     ZLWFT0_i, ZFLDN0, ZFLUP0, &
5837                     ZSWFT0_i, ZFSDN0, ZFSUP0)
5838          ENDIF !ok_4xCO2atm
5839       ENDIF ! aerosol_couple
5840       itaprad = 0
5841       !
5842       !  If Iflag_radia >=2, reset pertubed variables
5843       !
5844       IF (iflag_radia .ge. 2) THEN
5845          zxtsol(:) = zsav_tsol (:)
5846       ENDIF
5847    ENDIF ! MOD(itaprad,radpas)
5848    itaprad = itaprad + 1
5849
5850    IF (iflag_radia.eq.0) THEN
5851       IF (prt_level.ge.9) THEN
5852          PRINT *,'--------------------------------------------------'
5853          PRINT *,'>>>> ATTENTION rayonnement desactive pour ce cas'
5854          PRINT *,'>>>>           heat et cool mis a zero '
5855          PRINT *,'--------------------------------------------------'
5856       ENDIF
5857       heat=0.
5858       cool=0.
5859       sollw=0.   ! MPL 01032011
5860       solsw=0.
5861       radsol=0.
5862       swup=0.    ! MPL 27102011 pour les fichiers AMMA_profiles et AMMA_scalars
5863       swup0=0.
5864       lwup=0.
5865       lwup0=0.
5866       lwdn=0.
5867       lwdn0=0.
5868    ENDIF
5869
5870    !
5871    ! Calculer radsol a l'exterieur de radlwsw
5872    ! pour prendre en compte le cycle diurne
5873    ! recode par Olivier Boucher en sept 2015
5874    !
5875    radsol=solsw*swradcorr+sollw
5876
5877    IF (ok_4xCO2atm) THEN
5878       radsolp=solswp*swradcorr+sollwp
5879    ENDIF
5880
5881    !
5882    ! Ajouter la tendance des rayonnements (tous les pas)
5883    ! avec une correction pour le cycle diurne dans le SW
5884    !
5885
5886    DO k=1, klev
5887       d_t_swr(:,k)=swradcorr(:)*heat(:,k)*phys_tstep/RDAY
5888       d_t_sw0(:,k)=swradcorr(:)*heat0(:,k)*phys_tstep/RDAY
5889       d_t_lwr(:,k)=-cool(:,k)*phys_tstep/RDAY
5890       d_t_lw0(:,k)=-cool0(:,k)*phys_tstep/RDAY
5891    ENDDO
5892
5893    CALL add_phys_tend(du0,dv0,d_t_swr,dq0,dql0,dqi0,dqbs0,paprs,'SW',abortphy,flag_inhib_tend,itap,0 &
5894#ifdef ISO
5895     &    ,dxt0,dxtl0,dxti0 &
5896#endif     
5897     &  )
5898    CALL prt_enerbil('SW',itap)
5899    CALL add_phys_tend(du0,dv0,d_t_lwr,dq0,dql0,dqi0,dqbs0,paprs,'LW',abortphy,flag_inhib_tend,itap,0 &
5900#ifdef ISO
5901     &    ,dxt0,dxtl0,dxti0 &
5902#endif     
5903     &  )
5904    CALL prt_enerbil('LW',itap)
5905
5906    !
5907    IF (mydebug) THEN
5908       CALL writefield_phy('u_seri',u_seri,nbp_lev)
5909       CALL writefield_phy('v_seri',v_seri,nbp_lev)
5910       CALL writefield_phy('t_seri',t_seri,nbp_lev)
5911       CALL writefield_phy('q_seri',q_seri,nbp_lev)
5912    ENDIF
5913
5914    ! Calculer l'hydrologie de la surface
5915    !
5916    !      CALL hydrol(dtime,pctsrf,rain_fall, snow_fall, zxevap,
5917    !     .            agesno, ftsol,fqsurf,fsnow, ruis)
5918    !
5919
5920    !
5921    ! Calculer le bilan du sol et la derive de temperature (couplage)
5922    !
5923    DO i = 1, klon
5924       !         bils(i) = radsol(i) - sens(i) - evap(i)*RLVTT
5925       ! a la demande de JLD
5926       bils(i) = radsol(i) - sens(i) + zxfluxlat(i)
5927    ENDDO
5928    !
5929    !moddeblott(jan95)
5930    ! Appeler le programme de parametrisation de l'orographie
5931    ! a l'echelle sous-maille:
5932    !
5933    IF (prt_level .GE.10) THEN
5934       print *,' call orography ? ', ok_orodr
5935    ENDIF
5936    !
5937    IF (ok_orodr) THEN
5938       !
5939       !  selection des points pour lesquels le shema est actif:
5940       igwd=0
5941       DO i=1,klon
5942          itest(i)=0
5943          zrel_oro(i)=zstd(i)/(max(zsig(i),1.E-8)*sqrt(cell_area(i)))
5944          !zrel_oro: relative mountain height wrt relief explained by mean slope
5945          ! -> condition on zrel_oro can deactivate the drag on tilted planar terrains
5946          !    such as ice sheets (work by V. Wiener)
5947          ! zpmm_orodr_t and zstd_orodr_t are activation thresholds set by F. Lott to
5948          ! earn computation time but they are not physical.
5949          IF (((zpic(i)-zmea(i)).GT.zpmm_orodr_t).AND.(zstd(i).GT.zstd_orodr_t).AND.(zrel_oro(i).LE.zrel_oro_t)) THEN
5950             itest(i)=1
5951             igwd=igwd+1
5952             idx(igwd)=i
5953          ENDIF
5954       ENDDO
5955       !        igwdim=MAX(1,igwd)
5956       !
5957       IF (ok_strato) THEN
5958
5959          CALL drag_noro_strato(0,klon,klev,phys_tstep,paprs,pplay, &
5960               zmea,zstd, zsig, zgam, zthe,zpic,zval, &
5961               igwd,idx,itest, &
5962               t_seri, u_seri, v_seri, &
5963               zulow, zvlow, zustrdr, zvstrdr, &
5964               d_t_oro, d_u_oro, d_v_oro)
5965
5966       ELSE
5967          CALL drag_noro(klon,klev,phys_tstep,paprs,pplay, &
5968               zmea,zstd, zsig, zgam, zthe,zpic,zval, &
5969               igwd,idx,itest, &
5970               t_seri, u_seri, v_seri, &
5971               zulow, zvlow, zustrdr, zvstrdr, &
5972               d_t_oro, d_u_oro, d_v_oro)
5973       ENDIF
5974       !
5975       !  ajout des tendances
5976       !-----------------------------------------------------------------------
5977       ! ajout des tendances de la trainee de l'orographie
5978       CALL add_phys_tend(d_u_oro,d_v_oro,d_t_oro,dq0,dql0,dqi0,dqbs0,paprs,'oro', &
5979            abortphy,flag_inhib_tend,itap,0 &
5980#ifdef ISO
5981     &    ,dxt0,dxtl0,dxti0 &
5982#endif     
5983     &   )
5984       CALL prt_enerbil('oro',itap)
5985       !----------------------------------------------------------------------
5986       !
5987    ENDIF ! fin de test sur ok_orodr
5988    !
5989    IF (mydebug) THEN
5990       CALL writefield_phy('u_seri',u_seri,nbp_lev)
5991       CALL writefield_phy('v_seri',v_seri,nbp_lev)
5992       CALL writefield_phy('t_seri',t_seri,nbp_lev)
5993       CALL writefield_phy('q_seri',q_seri,nbp_lev)
5994    ENDIF
5995
5996    IF (ok_orolf) THEN
5997       !
5998       !  selection des points pour lesquels le shema est actif:
5999       igwd=0
6000       DO i=1,klon
6001          itest(i)=0
6002          !zrel_oro: relative mountain height wrt relief explained by mean slope
6003          ! -> condition on zrel_oro can deactivate the lifting on tilted planar terrains
6004          !    such as ice sheets (work by V. Wiener)
6005          zrel_oro(i)=zstd(i)/(max(zsig(i),1.E-8)*sqrt(cell_area(i)))
6006          IF (((zpic(i)-zmea(i)).GT.zpmm_orolf_t).AND.(zrel_oro(i).LE.zrel_oro_t)) THEN
6007             itest(i)=1
6008             igwd=igwd+1
6009             idx(igwd)=i
6010          ENDIF
6011       ENDDO
6012       !        igwdim=MAX(1,igwd)
6013       !
6014       IF (ok_strato) THEN
6015
6016          CALL lift_noro_strato(klon,klev,phys_tstep,paprs,pplay, &
6017               latitude_deg,zmea,zstd,zpic,zgam,zthe,zpic,zval, &
6018               igwd,idx,itest, &
6019               t_seri, u_seri, v_seri, &
6020               zulow, zvlow, zustrli, zvstrli, &
6021               d_t_lif, d_u_lif, d_v_lif               )
6022
6023       ELSE
6024          CALL lift_noro(klon,klev,phys_tstep,paprs,pplay, &
6025               latitude_deg,zmea,zstd,zpic, &
6026               itest, &
6027               t_seri, u_seri, v_seri, &
6028               zulow, zvlow, zustrli, zvstrli, &
6029               d_t_lif, d_u_lif, d_v_lif)
6030       ENDIF
6031
6032       ! ajout des tendances de la portance de l'orographie
6033       CALL add_phys_tend(d_u_lif, d_v_lif, d_t_lif, dq0, dql0, dqi0, dqbs0,paprs, &
6034            'lif', abortphy,flag_inhib_tend,itap,0 &
6035#ifdef ISO
6036     &    ,dxt0,dxtl0,dxti0 &
6037#endif     
6038     &   )
6039       CALL prt_enerbil('lif',itap)
6040    ENDIF ! fin de test sur ok_orolf
6041
6042    IF (ok_hines) then
6043       !  HINES GWD PARAMETRIZATION
6044       east_gwstress=0.
6045       west_gwstress=0.
6046       du_gwd_hines=0.
6047       dv_gwd_hines=0.
6048       CALL hines_gwd(klon, klev, phys_tstep, paprs, pplay, latitude_deg, t_seri, &
6049            u_seri, v_seri, zustr_gwd_hines, zvstr_gwd_hines, d_t_hin, &
6050            du_gwd_hines, dv_gwd_hines)
6051       zustr_gwd_hines=0.
6052       zvstr_gwd_hines=0.
6053       DO k = 1, klev
6054          zustr_gwd_hines(:)=zustr_gwd_hines(:)+ du_gwd_hines(:, k)/phys_tstep &
6055               * (paprs(:, k)-paprs(:, k+1))/rg
6056          zvstr_gwd_hines(:)=zvstr_gwd_hines(:)+ dv_gwd_hines(:, k)/phys_tstep &
6057               * (paprs(:, k)-paprs(:, k+1))/rg
6058       ENDDO
6059
6060       d_t_hin(:, :)=0.
6061       CALL add_phys_tend(du_gwd_hines, dv_gwd_hines, d_t_hin, dq0, dql0, &
6062            dqi0, dqbs0,paprs, 'hin', abortphy,flag_inhib_tend,itap,0 &
6063#ifdef ISO
6064     &    ,dxt0,dxtl0,dxti0 &
6065#endif     
6066     &   )
6067       CALL prt_enerbil('hin',itap)
6068    ENDIF
6069
6070    IF (.not. ok_hines .and. ok_gwd_rando) then
6071       ! ym missing init for east_gwstress & west_gwstress -> added in phys_local_var_mod
6072       CALL acama_GWD_rando(PHYS_TSTEP, pplay, latitude_deg, t_seri, u_seri, &
6073            v_seri, rot, zustr_gwd_front, zvstr_gwd_front, du_gwd_front, &
6074            dv_gwd_front, east_gwstress, west_gwstress)
6075       zustr_gwd_front=0.
6076       zvstr_gwd_front=0.
6077       DO k = 1, klev
6078          zustr_gwd_front(:)=zustr_gwd_front(:)+ du_gwd_front(:, k)/phys_tstep &
6079               * (paprs(:, k)-paprs(:, k+1))/rg
6080          zvstr_gwd_front(:)=zvstr_gwd_front(:)+ dv_gwd_front(:, k)/phys_tstep &
6081               * (paprs(:, k)-paprs(:, k+1))/rg
6082       ENDDO
6083
6084       CALL add_phys_tend(du_gwd_front, dv_gwd_front, dt0, dq0, dql0, dqi0, dqbs0, &
6085            paprs, 'front_gwd_rando', abortphy,flag_inhib_tend,itap,0 &
6086#ifdef ISO
6087     &    ,dxt0,dxtl0,dxti0 &
6088#endif     
6089     &   )
6090       CALL prt_enerbil('front_gwd_rando',itap)
6091    ENDIF
6092
6093    IF (ok_gwd_rando) THEN
6094       CALL FLOTT_GWD_rando(PHYS_TSTEP, pplay, t_seri, u_seri, v_seri, &
6095            rain_fall + snow_fall, zustr_gwd_rando, zvstr_gwd_rando, &
6096            du_gwd_rando, dv_gwd_rando, east_gwstress, west_gwstress)
6097       CALL add_phys_tend(du_gwd_rando, dv_gwd_rando, dt0, dq0, dql0, dqi0, dqbs0, &
6098            paprs, 'flott_gwd_rando', abortphy,flag_inhib_tend,itap,0 &
6099#ifdef ISO
6100     &    ,dxt0,dxtl0,dxti0 &
6101#endif     
6102     &  )
6103       CALL prt_enerbil('flott_gwd_rando',itap)
6104       zustr_gwd_rando=0.
6105       zvstr_gwd_rando=0.
6106       DO k = 1, klev
6107          zustr_gwd_rando(:)=zustr_gwd_rando(:)+ du_gwd_rando(:, k)/phys_tstep &
6108               * (paprs(:, k)-paprs(:, k+1))/rg
6109          zvstr_gwd_rando(:)=zvstr_gwd_rando(:)+ dv_gwd_rando(:, k)/phys_tstep &
6110               * (paprs(:, k)-paprs(:, k+1))/rg
6111       ENDDO
6112    ENDIF
6113
6114    ! STRESS NECESSAIRES: TOUTE LA PHYSIQUE
6115
6116    IF (mydebug) THEN
6117       CALL writefield_phy('u_seri',u_seri,nbp_lev)
6118       CALL writefield_phy('v_seri',v_seri,nbp_lev)
6119       CALL writefield_phy('t_seri',t_seri,nbp_lev)
6120       CALL writefield_phy('q_seri',q_seri,nbp_lev)
6121    ENDIF
6122
6123    DO i = 1, klon
6124       zustrph(i)=0.
6125       zvstrph(i)=0.
6126    ENDDO
6127    DO k = 1, klev
6128       DO i = 1, klon
6129          zustrph(i)=zustrph(i)+(u_seri(i,k)-u(i,k))/phys_tstep* &
6130               (paprs(i,k)-paprs(i,k+1))/rg
6131          zvstrph(i)=zvstrph(i)+(v_seri(i,k)-v(i,k))/phys_tstep* &
6132               (paprs(i,k)-paprs(i,k+1))/rg
6133       ENDDO
6134    ENDDO
6135    !
6136    !IM calcul composantes axiales du moment angulaire et couple des montagnes
6137    !
6138    IF (is_sequential .and. ok_orodr) THEN
6139       CALL aaam_bud (27,klon,klev,jD_cur-jD_ref,jH_cur, &
6140            ra,rg,romega, &
6141            latitude_deg,longitude_deg,pphis, &
6142            zustrdr,zustrli,zustrph, &
6143            zvstrdr,zvstrli,zvstrph, &
6144            paprs,u,v, &
6145            aam, torsfc)
6146    ENDIF
6147    !IM cf. FLott END
6148    !DC Calcul de la tendance due au methane
6149    IF (ok_qch4) THEN
6150       CALL METHOX(1,klon,klon,klev,q_seri,d_q_ch4,pplay)
6151       ! ajout de la tendance d'humidite due au methane
6152       d_q_ch4_dtime(:,:) = d_q_ch4(:,:)*phys_tstep
6153#ifdef ISO
6154       d_xt_ch4_dtime(:,:,:) = d_xt_ch4(:,:,:)*phys_tstep
6155#endif
6156       CALL add_phys_tend(du0, dv0, dt0, d_q_ch4_dtime, dql0, dqi0, dqbs0, paprs, &
6157            'q_ch4', abortphy,flag_inhib_tend,itap,0 &
6158#ifdef ISO
6159     &    ,d_xt_ch4_dtime,dxtl0,dxti0 &
6160#endif     
6161     &   )
6162       d_q_ch4(:,:) = d_q_ch4_dtime(:,:)/phys_tstep
6163#ifdef ISO
6164       d_xt_ch4(:,:,:) = d_xt_ch4_dtime(:,:,:)/phys_tstep
6165#endif
6166    ENDIF
6167    !
6168    !
6169#ifdef CPP_StratAer
6170    IF (ok_qemiss) THEN
6171       flh2o=1
6172       IF(flag_verbose_strataer) THEN
6173          print *,'IN physiq_mod: ok_qemiss =yes (',ok_qemiss,'), flh2o=',flh2o
6174          print *,'IN physiq_mod: flag_emit=',flag_emit,', nErupt=',nErupt
6175          print *,'IN physiq_mod: nAerErupt=',nAerErupt
6176       ENDIF
6177       
6178       SELECT CASE(flag_emit)
6179       CASE(1) ! emission volc H2O dans LMDZ
6180          DO ieru=1, nErupt
6181             IF (year_cur==year_emit_vol(ieru).AND.&
6182                  mth_cur==mth_emit_vol(ieru).AND.&
6183                  day_cur>=day_emit_vol(ieru).AND.&
6184                  day_cur<(day_emit_vol(ieru)+injdur)) THEN
6185               
6186                IF(flag_verbose_strataer) print *,'IN physiq_mod: date=',year_cur,mth_cur,day_cur
6187                ! initialisation tendance q emission
6188                d_q_emiss(:,:)=0.
6189                ! daily injection mass emission - NL
6190                m_H2O_emiss_vol_daily = m_H2O_emiss_vol(ieru)/(REAL(injdur)&
6191                     *REAL(ponde_lonlat_vol(ieru)))
6192                !
6193                CALL STRATEMIT(pdtphys,pdtphys,latitude_deg,longitude_deg,t_seri,&
6194                     pplay,paprs,tr_seri,&
6195                     m_H2O_emiss_vol_daily,&
6196                     xlat_min_vol(ieru),xlat_max_vol(ieru),&
6197                     xlon_min_vol(ieru),xlon_max_vol(ieru),&
6198                     altemiss_vol(ieru),&
6199                     sigma_alt_vol(ieru),1,&
6200                     1,nAerErupt+1,0)
6201               
6202                IF(flag_verbose_strataer) print *,'IN physiq_mod: min max d_q_emiss=',&
6203                     minval(d_q_emiss),maxval(d_q_emiss)
6204               
6205                CALL add_phys_tend(du0, dv0, dt0, d_q_emiss, dql0, dqi0, dqbs0, paprs, &
6206                     'q_emiss',abortphy,flag_inhib_tend,itap,0)
6207                IF (abortphy==1) Print*,'ERROR ABORT TEND EMISS'
6208             ENDIF
6209          ENDDO
6210          flh2o=0
6211       END SELECT ! emission scenario (flag_emit)
6212    ENDIF
6213#endif
6214
6215
6216!===============================================================
6217!            Additional tendency of TKE due to orography
6218!===============================================================
6219!
6220! Inititialization
6221!------------------
6222
6223       addtkeoro=0   
6224       CALL getin_p('addtkeoro',addtkeoro)
6225     
6226       IF (prt_level.ge.5) &
6227            print*,'addtkeoro', addtkeoro
6228           
6229       alphatkeoro=1.   
6230       CALL getin_p('alphatkeoro',alphatkeoro)
6231       alphatkeoro=min(max(0.,alphatkeoro),1.)
6232
6233       smallscales_tkeoro=.FALSE.   
6234       CALL getin_p('smallscales_tkeoro',smallscales_tkeoro)
6235
6236
6237       dtadd(:,:)=0.
6238       duadd(:,:)=0.
6239       dvadd(:,:)=0.
6240
6241! Choices for addtkeoro:
6242!      ** 0 no TKE tendency from orography   
6243!      ** 1 we include a fraction alphatkeoro of the whole tendency duoro
6244!      ** 2 we include a fraction alphatkeoro of the gravity wave part of duoro
6245!
6246
6247       IF (addtkeoro .GT. 0 .AND. ok_orodr ) THEN
6248!      -------------------------------------------
6249
6250
6251       !  selection des points pour lesquels le schema est actif:
6252
6253
6254  IF (addtkeoro .EQ. 1 ) THEN
6255
6256            duadd(:,:)=alphatkeoro*d_u_oro(:,:)
6257            dvadd(:,:)=alphatkeoro*d_v_oro(:,:)
6258
6259  ELSE IF (addtkeoro .EQ. 2) THEN
6260
6261     IF (smallscales_tkeoro) THEN
6262       igwd=0
6263       DO i=1,klon
6264          itest(i)=0
6265! Etienne: ici je prends en compte plus de relief que la routine drag_noro_strato
6266! car on peut s'attendre a ce que les petites echelles produisent aussi de la TKE
6267! Mais attention, cela ne va pas dans le sens de la conservation de l'energie!
6268          IF ((zstd(i).GT.1.0) .AND.(zrel_oro(i).LE.zrel_oro_t)) THEN
6269             itest(i)=1
6270             igwd=igwd+1
6271             idx(igwd)=i
6272          ENDIF
6273       ENDDO
6274
6275     ELSE
6276
6277       igwd=0
6278       DO i=1,klon
6279          itest(i)=0
6280        IF (((zpic(i)-zmea(i)).GT.zpmm_orodr_t).AND.(zstd(i).GT.zstd_orodr_t).AND.(zrel_oro(i).LE.zrel_oro_t)) THEN
6281             itest(i)=1
6282             igwd=igwd+1
6283             idx(igwd)=i
6284        ENDIF
6285       ENDDO
6286
6287     ENDIF
6288
6289     CALL drag_noro_strato(addtkeoro,klon,klev,phys_tstep,paprs,pplay, &
6290               zmea,zstd, zsig, zgam, zthe,zpic,zval, &
6291               igwd,idx,itest, &
6292               t_seri, u_seri, v_seri, &
6293               zulow, zvlow, zustrdr, zvstrdr, &
6294               d_t_oro_gw, d_u_oro_gw, d_v_oro_gw)
6295
6296     zustrdr(:)=0.
6297     zvstrdr(:)=0.
6298     zulow(:)=0.
6299     zvlow(:)=0.
6300
6301     duadd(:,:)=alphatkeoro*d_u_oro_gw(:,:)
6302     dvadd(:,:)=alphatkeoro*d_v_oro_gw(:,:)
6303  ENDIF
6304
6305
6306   ! TKE update from subgrid temperature and wind tendencies
6307   !----------------------------------------------------------
6308    forall (k=1: nbp_lev) exner(:, k) = (pplay(:, k)/paprs(:,1))**RKAPPA
6309
6310
6311    CALL tend_to_tke(pdtphys,paprs,exner,t_seri,u_seri,v_seri,dtadd,duadd,dvadd,pctsrf,pbl_tke)
6312   !
6313   ! Prevent pbl_tke_w from becoming negative
6314    wake_delta_pbl_tke(:,:,:) = max(wake_delta_pbl_tke(:,:,:), -pbl_tke(:,:,:))
6315   !
6316
6317       ENDIF
6318!      -----
6319!===============================================================
6320
6321#ifdef ISO
6322#ifdef ISOVERIF
6323if (iso_HDO.gt.0) then
6324      call iso_verif_aberrant_enc_vect2D( &
6325     &           xt_seri,q_seri, &
6326     &           'physiq 5924, juste apres methox',ntraciso,klon,klev)
6327endif     
6328      if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then 
6329        do i=1,klon
6330         do k=1,klev
6331            if (q_seri(i,k).gt.ridicule) then 
6332               if (iso_verif_o18_aberrant_nostop( &
6333     &              xt_seri(iso_HDO,i,k)/q_seri(i,k), &
6334     &              xt_seri(iso_O18,i,k)/q_seri(i,k), &
6335     &              'physiq 5937a, juste apres methox, qv').eq.1) then
6336                  write(*,*) 'physic 2444: i,k,q_seri(i,k)=',i,k,q_seri(i,k)
6337                  write(*,*) 'xt_seri(:,i,k)=',xt_seri(:,i,k)
6338                  stop
6339              endif !  if (iso_verif_o18_aberrant_nostop
6340            endif !if (q_seri(i,k).gt.errmax) then 
6341            if (ql_seri(i,k).gt.ridicule) then 
6342               if (iso_verif_o18_aberrant_nostop( &
6343     &              xtl_seri(iso_HDO,i,k)/ql_seri(i,k), &
6344     &              xtl_seri(iso_O18,i,k)/ql_seri(i,k), &
6345     &              'physiq 5937b, juste apres methox, ql').eq.1) then
6346                  write(*,*) 'i,k,ql_seri(i,k)=',i,k,ql_seri(i,k)
6347                  !stop
6348              endif !  if (iso_verif_o18_aberrant_nostop
6349            endif !if (q_seri(i,k).gt.errmax) then 
6350            if (qs_seri(i,k).gt.ridicule) then 
6351               if (iso_verif_o18_aberrant_nostop( &
6352     &              xts_seri(iso_HDO,i,k)/qs_seri(i,k), &
6353     &              xts_seri(iso_O18,i,k)/qs_seri(i,k), &
6354     &              'physiq 5937c, juste apres methox, qs').eq.1) then
6355                  write(*,*) 'i,k,qs_seri(i,k)=',i,k,qs_seri(i,k)
6356                  !stop
6357              endif !  if (iso_verif_o18_aberrant_nostop
6358            endif !if (q_seri(i,k).gt.errmax) then
6359          enddo !k=1,klev
6360         enddo  !i=1,klon
6361        endif !if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then
6362#endif
6363#endif
6364
6365    !====================================================================
6366    ! Interface Simulateur COSP (Calipso, ISCCP, MISR, ..)
6367    !====================================================================
6368    ! Abderrahmane 24.08.09
6369
6370    IF (ok_cosp) THEN
6371       ! adeclarer
6372#ifdef CPP_COSP
6373       IF (itap.eq.1.or.MOD(itap,NINT(freq_cosp/phys_tstep)).EQ.0) THEN
6374
6375          IF (prt_level .GE.10) THEN
6376             print*,'freq_cosp',freq_cosp
6377          ENDIF
6378          mr_ozone=wo(:, :, 1) * dobson_u * 1e3 / zmasse
6379          !       print*,'Dans physiq.F avant appel cosp ref_liq,ref_ice=',
6380          !     s        ref_liq,ref_ice
6381          CALL phys_cosp(itap,phys_tstep,freq_cosp, &
6382               ok_mensuelCOSP,ok_journeCOSP,ok_hfCOSP, &
6383               ecrit_mth,ecrit_day,ecrit_hf, ok_all_xml, missing_val, &
6384               klon,klev,longitude_deg,latitude_deg,presnivs,overlap, &
6385               JrNt,ref_liq,ref_ice, &
6386               pctsrf(:,is_ter)+pctsrf(:,is_lic), &
6387               zu10m,zv10m,pphis, &
6388               zphi,paprs(:,1:klev),pplay,zxtsol,t_seri, &
6389               qx(:,:,ivap),zx_rh,cldfra,rnebcon,flwc,fiwc, &
6390               prfl(:,1:klev),psfl(:,1:klev), &
6391               pmflxr(:,1:klev),pmflxs(:,1:klev), &
6392               mr_ozone,cldtau, cldemi)
6393
6394          !     L         calipso2D,calipso3D,cfadlidar,parasolrefl,atb,betamol,
6395          !     L          cfaddbze,clcalipso2,dbze,cltlidarradar,
6396          !     M          clMISR,
6397          !     R          clisccp2,boxtauisccp,boxptopisccp,tclisccp,ctpisccp,
6398          !     I          tauisccp,albisccp,meantbisccp,meantbclrisccp)
6399
6400       ENDIF
6401#endif
6402
6403#ifdef CPP_COSP2
6404       IF (itap.eq.1.or.MOD(itap,NINT(freq_cosp/phys_tstep)).EQ.0) THEN
6405
6406          IF (prt_level .GE.10) THEN
6407             print*,'freq_cosp',freq_cosp
6408          ENDIF
6409          mr_ozone=wo(:, :, 1) * dobson_u * 1e3 / zmasse
6410                 print*,'Dans physiq.F avant appel '
6411          !     s        ref_liq,ref_ice
6412          CALL phys_cosp2(itap,phys_tstep,freq_cosp, &
6413               ok_mensuelCOSP,ok_journeCOSP,ok_hfCOSP, &
6414               ecrit_mth,ecrit_day,ecrit_hf, ok_all_xml, missing_val, &
6415               klon,klev,longitude_deg,latitude_deg,presnivs,overlap, &
6416               JrNt,ref_liq,ref_ice, &
6417               pctsrf(:,is_ter)+pctsrf(:,is_lic), &
6418               zu10m,zv10m,pphis, &
6419               zphi,paprs(:,1:klev),pplay,zxtsol,t_seri, &
6420               qx(:,:,ivap),zx_rh,cldfra,rnebcon,flwc,fiwc, &
6421               prfl(:,1:klev),psfl(:,1:klev), &
6422               pmflxr(:,1:klev),pmflxs(:,1:klev), &
6423               mr_ozone,cldtau, cldemi)
6424       ENDIF
6425#endif
6426
6427#ifdef CPP_COSPV2
6428       IF (itap.eq.1.or.MOD(itap,NINT(freq_cosp/phys_tstep)).EQ.0) THEN
6429!        IF (MOD(itap,NINT(freq_cosp/phys_tstep)).EQ.0) THEN
6430
6431          IF (prt_level .GE.10) THEN
6432             print*,'freq_cosp',freq_cosp
6433          ENDIF
6434           DO k = 1, klev
6435             DO i = 1, klon
6436               phicosp(i,k) = pphi(i,k) + pphis(i)
6437             ENDDO
6438           ENDDO
6439          mr_ozone=wo(:, :, 1) * dobson_u * 1e3 / zmasse
6440                 print*,'Dans physiq.F avant appel '
6441          !     s        ref_liq,ref_ice
6442          CALL lmdz_cosp_interface(itap,phys_tstep,freq_cosp, &
6443               ok_mensuelCOSP,ok_journeCOSP,ok_hfCOSP, &
6444               ecrit_mth,ecrit_day,ecrit_hf, ok_all_xml, missing_val, &
6445               klon,klev,longitude_deg,latitude_deg,presnivs,overlap, &
6446               JrNt,ref_liq,ref_ice, &
6447               pctsrf(:,is_ter)+pctsrf(:,is_lic), &
6448               zu10m,zv10m,pphis, &
6449               zphi,paprs(:,1:klev),pplay,zxtsol,t_seri, &
6450               qx(:,:,ivap),zx_rh,cldfra,rnebcon,flwc,fiwc, &
6451               prfl(:,1:klev),psfl(:,1:klev), &
6452               pmflxr(:,1:klev),pmflxs(:,1:klev), &
6453               mr_ozone,cldtau, cldemi)
6454       ENDIF
6455#endif
6456
6457    ENDIF  !ok_cosp
6458
6459!c====================================================================
6460!c   Ajout de la production de tritium (naturelle et essais nucléaires)
6461!c====================================================================
6462!C
6463#ifdef ISO
6464#ifdef ISOVERIF
6465      call iso_verif_noNaN_vect2D(xt_seri, &
6466     &     'physiq 5595: avant appel tritium',ntraciso,klon,klev)
6467#endif
6468        call iso_tritium(paprs,pplay, &
6469     &           zphi,phys_tstep, &
6470     &           d_xt_prod_nucl, &
6471     &           d_xt_cosmo, &
6472     &           d_xt_decroiss, &
6473     &           xt_seri)
6474#ifdef ISOVERIF
6475      call iso_verif_noNaN_vect2D(xt_seri, &
6476     &     'physiq 5607: apres appel tritium',ntraciso,klon,klev)
6477       if (iso_HTO.gt.0) then ! Tritium
6478       ixt=iso_HTO
6479       do i=1,klon
6480       do k=1,klev
6481          if (iso_verif_positif_strict_nostop(xt_seri(ixt,i,k), &
6482     &      'physiq 5620 : xt_seri(HTO) nul ou negatif').eq.1) then
6483          write(*,*) 'ixt,i,klon,k,klev=',ixt,i,klon,k,klev
6484          write(*,*) 'xt_seri(iso_HTO,i,k)=',xt_seri(ixt,i,k)
6485          stop
6486          endif
6487       enddo
6488       enddo
6489       endif
6490#endif
6491!      write(*,*)'itap=',itap
6492!      write(*,*)'itau_phy=',itau_phy
6493!      write(*,*)'jD_cur=',jD_cur
6494#endif
6495!ifdef ISO
6496
6497
6498! Marine
6499
6500  IF (ok_airs) then
6501
6502  IF (itap.eq.1.or.MOD(itap,NINT(freq_airs/phys_tstep)).EQ.0) THEN
6503     write(*,*) 'je vais appeler simu_airs, ok_airs, freq_airs=', ok_airs, freq_airs
6504     CALL simu_airs(itap,rneb, t_seri, cldemi, fiwc, ref_ice, pphi, pplay, paprs,&
6505        & map_prop_hc,map_prop_hist,&
6506        & map_emis_hc,map_iwp_hc,map_deltaz_hc,map_pcld_hc,map_tcld_hc,&
6507        & map_emis_Cb,map_pcld_Cb,map_tcld_Cb,&
6508        & map_emis_ThCi,map_pcld_ThCi,map_tcld_ThCi,&
6509        & map_emis_Anv,map_pcld_Anv,map_tcld_Anv,&
6510        & map_emis_hist,map_iwp_hist,map_deltaz_hist,map_rad_hist,&
6511        & map_ntot,map_hc,map_hist,&
6512        & map_Cb,map_ThCi,map_Anv,&
6513        & alt_tropo )
6514  ENDIF
6515
6516  ENDIF  ! ok_airs
6517
6518
6519    ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
6520    !AA
6521    !AA Installation de l'interface online-offline pour traceurs
6522    !AA
6523    !====================================================================
6524    !   Calcul  des tendances traceurs
6525    !====================================================================
6526    !
6527
6528    IF (type_trac == 'repr') THEN
6529!MM pas d'impact, car on recupere q_seri,tr_seri,t_seri via phys_local_var_mod
6530!MM                               dans Reprobus
6531       sh_in(:,:) = q_seri(:,:)
6532#ifdef REPROBUS
6533       d_q_rep(:,:) = 0.
6534       d_ql_rep(:,:) = 0.
6535       d_qi_rep(:,:) = 0.
6536#endif
6537    ELSE
6538       sh_in(:,:) = qx(:,:,ivap)
6539       IF (nqo >= 3) THEN
6540          ch_in(:,:) = qx(:,:,iliq) + qx(:,:,isol)
6541       ELSE
6542          ch_in(:,:) = qx(:,:,iliq)
6543       ENDIF
6544    ENDIF
6545
6546#ifdef CPP_Dust
6547    !  Avec SPLA, iflag_phytrac est forcé =1
6548    CALL       phytracr_spl ( debut,lafin , jD_cur,jH_cur,iflag_con,       &  ! I
6549                      pdtphys,ftsol,                                   &  ! I
6550                      t,q_seri,paprs,pplay,RHcl,                  &  ! I
6551                      pmfu, pmfd, pen_u, pde_u, pen_d, pde_d,          &  ! I
6552                      coefh(1:klon,1:klev,is_ave), cdragh, cdragm, u1, v1,                 &  ! I
6553                      u_seri, v_seri, latitude_deg, longitude_deg,  &
6554                      pphis,pctsrf,pmflxr,pmflxs,prfl,psfl,            &  ! I
6555                      da,phi,phi2,d1a,dam,mp,ep,sigd,sij,clw,elij,     &  ! I
6556                      epmlmMm,eplaMm,upwd,dnwd,itop_con,ibas_con,      &  ! I
6557                      ev,wdtrainA,  wdtrainM,wght_cvfd,              &  ! I
6558                      fm_therm, entr_therm, rneb,                      &  ! I
6559                      beta_prec_fisrt,beta_prec, & !I
6560                      zu10m,zv10m,wstar,ale_bl,ale_wake,               &  ! I
6561                      d_tr_dyn,tr_seri)
6562
6563#else
6564    IF (iflag_phytrac == 1 ) THEN
6565      CALL phytrac ( &
6566         itap,     days_elapsed+1,    jH_cur,   debut, &
6567         lafin,    phys_tstep,     u, v,     t, &
6568         paprs,    pplay,     pmfu,     pmfd, &
6569         pen_u,    pde_u,     pen_d,    pde_d, &
6570         cdragh,   coefh(1:klon,1:klev,is_ave),   fm_therm, entr_therm, &
6571         u1,       v1,        ftsol,    pctsrf, &
6572         zustar,   zu10m,     zv10m, &
6573         wstar(:,is_ave),    ale_bl,         ale_wake, &
6574         latitude_deg, longitude_deg, &
6575         frac_impa,frac_nucl, beta_prec_fisrt,beta_prec, &
6576         presnivs, pphis,     pphi,     albsol1, &
6577         sh_in,   ch_in,    rhcl,      cldfra,   rneb, &
6578         diafra,   radocond,    itop_con, ibas_con, &
6579         pmflxr,   pmflxs,    prfl,     psfl, &
6580         da,       phi,       mp,       upwd, &
6581         phi2,     d1a,       dam,      sij, wght_cvfd, &        !<<RomP+RL
6582         wdtrainA, wdtrainM,  sigd,     clw,elij, &   !<<RomP
6583         ev,       ep,        epmlmMm,  eplaMm, &     !<<RomP
6584         dnwd,     aerosol_couple,      flxmass_w, &
6585         tau_aero, piz_aero,  cg_aero,  ccm, &
6586         rfname, &
6587         d_tr_dyn, &                                 !<<RomP
6588         tr_seri, init_source)
6589#ifdef REPROBUS
6590
6591
6592          print*,'avt add phys rep',abortphy
6593
6594     CALL add_phys_tend &
6595            (du0,dv0,dt0,d_q_rep,d_ql_rep,d_qi_rep,dqbs0,paprs,&
6596             'rep',abortphy,flag_inhib_tend,itap,0)
6597        IF (abortphy==1) Print*,'ERROR ABORT REP'
6598
6599          print*,'apr add phys rep',abortphy
6600
6601#endif
6602    ENDIF    ! (iflag_phytrac=1)
6603
6604#endif
6605    !ENDIF    ! (iflag_phytrac=1)
6606
6607    IF (offline) THEN
6608
6609       IF (prt_level.ge.9) &
6610            print*,'Attention on met a 0 les thermiques pour phystoke'
6611       CALL phystokenc ( &
6612            nlon,klev,pdtphys,longitude_deg,latitude_deg, &
6613            t,pmfu, pmfd, pen_u, pde_u, pen_d, pde_d, &
6614            fm_therm,entr_therm, &
6615            cdragh,coefh(1:klon,1:klev,is_ave),u1,v1,ftsol,pctsrf, &
6616            frac_impa, frac_nucl, &
6617            pphis,cell_area,phys_tstep,itap, &
6618            qx(:,:,ivap),da,phi,mp,upwd,dnwd)
6619
6620
6621    ENDIF
6622
6623    !
6624    ! Calculer le transport de l'eau et de l'energie (diagnostique)
6625    !
6626    CALL transp (paprs,zxtsol, t_seri, q_seri, ql_seri, qs_seri, u_seri, v_seri, zphi, &
6627                 ue, ve, uq, vq, uwat, vwat)
6628    !
6629    !IM global posePB BEG
6630    IF(1.EQ.0) THEN
6631       !
6632       CALL transp_lay (paprs,zxtsol, t_seri, q_seri, u_seri, v_seri, zphi, &
6633            ve_lay, vq_lay, ue_lay, uq_lay)
6634       !
6635    ENDIF !(1.EQ.0) THEN
6636    !IM global posePB END
6637    !
6638    ! Accumuler les variables a stocker dans les fichiers histoire:
6639    !
6640
6641    !================================================================
6642    ! Conversion of kinetic and potential energy into heat, for
6643    ! parameterisation of subgrid-scale motions
6644    !================================================================
6645
6646    d_t_ec(:,:)=0.
6647    forall (k=1: nbp_lev) exner(:, k) = (pplay(:, k)/paprs(:,1))**RKAPPA
6648    CALL ener_conserv(klon,klev,pdtphys,u,v,t,qx,ivap,iliq,isol, &
6649         u_seri,v_seri,t_seri,q_seri,ql_seri,qs_seri,pbl_tke(:,:,is_ave)-tke0(:,:), &
6650         zmasse,exner,d_t_ec)
6651    t_seri(:,:)=t_seri(:,:)+d_t_ec(:,:)
6652
6653    !=======================================================================
6654    !   SORTIES
6655    !=======================================================================
6656    !
6657    !IM initialisation + calculs divers diag AMIP2
6658    !
6659    include "calcul_divers.h"
6660    !
6661    !IM Interpolation sur les niveaux de pression du NMC
6662    !   -------------------------------------------------
6663    !
6664    include "calcul_STDlev.h"
6665    !
6666    ! slp sea level pressure derived from Arpege-IFS : CALL ctstar + CALL pppmer
6667    CALL diag_slp(klon,t_seri,paprs,pplay,pphis,ptstar,pt0,slp)
6668    !
6669    !cc prw  = eau precipitable
6670    !   prlw = colonne eau liquide
6671    !   prlw = colonne eau solide
6672    !   prbsw = colonne neige soufflee
6673    prw(:) = 0.
6674    prlw(:) = 0.
6675    prsw(:) = 0.
6676    prbsw(:) = 0.
6677    DO k = 1, klev
6678       prw(:)  = prw(:)  + q_seri(:,k)*zmasse(:,k)
6679       prlw(:) = prlw(:) + ql_seri(:,k)*zmasse(:,k)
6680       prsw(:) = prsw(:) + qs_seri(:,k)*zmasse(:,k)
6681       prbsw(:)= prbsw(:) + qbs_seri(:,k)*zmasse(:,k)
6682    ENDDO
6683
6684#ifdef ISO
6685      DO i = 1, klon
6686      do ixt=1,ntraciso
6687       xtprw(ixt,i) = 0.
6688       DO k = 1, klev
6689        xtprw(ixt,i) = xtprw(ixt,i) + &
6690     &           xt_seri(ixt,i,k)*(paprs(i,k)-paprs(i,k+1))/RG
6691       ENDDO !DO k = 1, klev
6692      enddo !do ixt=1,ntraciso
6693      enddo !DO i = 1, klon
6694#endif
6695    !
6696    IF (ANY(type_trac == ['inca','inco'])) THEN
6697#ifdef INCA
6698       CALL VTe(VTphysiq)
6699       CALL VTb(VTinca)
6700
6701       CALL chemhook_end ( &
6702            phys_tstep, &
6703            pplay, &
6704            t_seri, &
6705            tr_seri(:,:,1+nqCO2:nbtr), &
6706            nbtr, &
6707            paprs, &
6708            q_seri, &
6709            cell_area, &
6710            pphi, &
6711            pphis, &
6712            zx_rh, &
6713            aps, bps, ap, bp, lafin)
6714
6715       CALL VTe(VTinca)
6716       CALL VTb(VTphysiq)
6717#endif
6718    ENDIF
6719
6720    IF (type_trac == 'repr') THEN
6721#ifdef REPROBUS
6722        CALL coord_hyb_rep(paprs, pplay, aps, bps, ap, bp, cell_area)
6723#endif
6724    ENDIF
6725
6726    !
6727    ! Convertir les incrementations en tendances
6728    !
6729    IF (prt_level .GE.10) THEN
6730       print *,'Convertir les incrementations en tendances '
6731    ENDIF
6732    !
6733    IF (mydebug) THEN
6734       CALL writefield_phy('u_seri',u_seri,nbp_lev)
6735       CALL writefield_phy('v_seri',v_seri,nbp_lev)
6736       CALL writefield_phy('t_seri',t_seri,nbp_lev)
6737       CALL writefield_phy('q_seri',q_seri,nbp_lev)
6738    ENDIF
6739
6740    DO k = 1, klev
6741       DO i = 1, klon
6742          d_u(i,k) = ( u_seri(i,k) - u(i,k) ) / phys_tstep
6743          d_v(i,k) = ( v_seri(i,k) - v(i,k) ) / phys_tstep
6744          d_t(i,k) = ( t_seri(i,k)-t(i,k) ) / phys_tstep
6745          d_qx(i,k,ivap) = ( q_seri(i,k) - qx(i,k,ivap) ) / phys_tstep
6746          d_qx(i,k,iliq) = ( ql_seri(i,k) - qx(i,k,iliq) ) / phys_tstep
6747          !CR: on ajoute le contenu en glace
6748          IF (nqo >= 3) THEN
6749             d_qx(i,k,isol) = ( qs_seri(i,k) - qx(i,k,isol) ) / phys_tstep
6750          ENDIF
6751          !--ice_sursat: nqo=4, on ajoute rneb
6752          IF (nqo.ge.4 .and. ok_ice_sursat) THEN
6753             d_qx(i,k,irneb) = ( rneb_seri(i,k) - qx(i,k,irneb) ) / phys_tstep
6754          ENDIF
6755
6756           IF (nqo.ge.4 .and. ok_bs) THEN
6757             d_qx(i,k,ibs) = ( qbs_seri(i,k) - qx(i,k,ibs) ) / phys_tstep
6758          ENDIF
6759
6760
6761       ENDDO
6762    ENDDO
6763
6764    ! C Risi: dispatcher les isotopes dans les xt_seri
6765#ifdef ISO
6766    do ixt=1,ntraciso
6767      DO k = 1, klev
6768       DO i = 1, klon
6769          iq=iqIsoPha(ixt,ivap)
6770          d_qx(i,k,iq) = ( xt_seri(ixt,i,k) - qx(i,k,iq) ) / phys_tstep
6771          iq=iqIsoPha(ixt,iliq)
6772          d_qx(i,k,iq) = ( xtl_seri(ixt,i,k) - qx(i,k,iq) ) / phys_tstep
6773          if (nqo.eq.3) then
6774             iq=iqIsoPha(ixt,isol)
6775             d_qx(i,k,iq) = ( xts_seri(ixt,i,k) - qx(i,k,iq) ) / phys_tstep
6776          endif
6777       enddo !DO i = 1, klon
6778      enddo ! DO k = 1, klev
6779    enddo !do ixt=1,ntraciso
6780#endif
6781    ! DC: All iterations are cycled if nqtot==nqo, so no nqtot>nqo condition required
6782    itr = 0
6783    DO iq = 1, nqtot
6784       IF(.NOT.tracers(iq)%isInPhysics) CYCLE
6785       itr = itr+1
6786       DO  k = 1, klev
6787          DO  i = 1, klon
6788             d_qx(i,k,iq) = ( tr_seri(i,k,itr) - qx(i,k,iq) ) / phys_tstep
6789          ENDDO
6790       ENDDO
6791    ENDDO
6792    !
6793    !IM rajout diagnostiques bilan KP pour analyse MJO par Jun-Ichi Yano
6794    !IM global posePB      include "write_bilKP_ins.h"
6795    !IM global posePB      include "write_bilKP_ave.h"
6796    !
6797
6798    !--OB mass fixer
6799    !--profile is corrected to force mass conservation of water
6800    IF (mass_fixer) THEN
6801    qql2(:)=0.0
6802    DO k = 1, klev
6803      qql2(:)=qql2(:)+(q_seri(:,k)+ql_seri(:,k)+qs_seri(:,k))*zmasse(:,k)
6804    ENDDO
6805   
6806#ifdef CPP_StratAer
6807    IF (ok_qemiss) THEN
6808       DO k = 1, klev
6809          qql1(:) = qql1(:)+d_q_emiss(:,k)*zmasse(:,k)
6810       ENDDO
6811    ENDIF
6812#endif
6813    IF (ok_qch4) THEN
6814       DO k = 1, klev
6815          qql1(:) = qql1(:)+d_q_ch4_dtime(:,k)*zmasse(:,k)
6816       ENDDO
6817    ENDIF
6818   
6819    DO i = 1, klon
6820      !--compute ratio of what q+ql should be with conservation to what it is
6821      corrqql=(qql1(i)+(evap(i)-rain_fall(i)-snow_fall(i))*pdtphys)/qql2(i)
6822      DO k = 1, klev
6823        q_seri(i,k) =q_seri(i,k)*corrqql
6824        ql_seri(i,k)=ql_seri(i,k)*corrqql
6825      ENDDO
6826    ENDDO
6827#ifdef ISO
6828    do ixt=1,ntraciso
6829    xtql2(ixt,:)=0.0
6830    DO k = 1, klev
6831      xtql2(ixt,:)=xtql2(ixt,:)+(xt_seri(ixt,:,k)+xtl_seri(ixt,:,k)+xts_seri(ixt,:,k))*zmasse(:,k)
6832    ENDDO
6833    DO i = 1, klon
6834      !--compute ratio of what q+ql should be with conservation to what it is
6835      corrxtql(ixt)=(xtql1(ixt,i)+(xtevap(ixt,i)-xtrain_fall(ixt,i)-xtsnow_fall(ixt,i))*pdtphys)/xtql2(ixt,i)
6836      DO k = 1, klev
6837        xt_seri(ixt,i,k) =xt_seri(ixt,i,k)*corrxtql(ixt)
6838        xtl_seri(ixt,i,k)=xtl_seri(ixt,i,k)*corrxtql(ixt)
6839      ENDDO
6840    ENDDO   
6841    enddo !do ixt=1,ntraciso
6842#endif
6843    ENDIF
6844    !--fin mass fixer
6845
6846    ! Sauvegarder les valeurs de t et q a la fin de la physique:
6847    !
6848    u_ancien(:,:)  = u_seri(:,:)
6849    v_ancien(:,:)  = v_seri(:,:)
6850    t_ancien(:,:)  = t_seri(:,:)
6851    q_ancien(:,:)  = q_seri(:,:)
6852    ql_ancien(:,:) = ql_seri(:,:)
6853    qs_ancien(:,:) = qs_seri(:,:)
6854    qbs_ancien(:,:) = qbs_seri(:,:)
6855    rneb_ancien(:,:) = rneb_seri(:,:)
6856#ifdef ISO
6857    xt_ancien(:,:,:)=xt_seri(:,:,:)
6858    xtl_ancien(:,:,:)=xtl_seri(:,:,:)
6859    xts_ancien(:,:,:)=xts_seri(:,:,:)
6860#endif
6861    CALL water_int(klon,klev,q_ancien,zmasse,prw_ancien)
6862    CALL water_int(klon,klev,ql_ancien,zmasse,prlw_ancien)
6863    CALL water_int(klon,klev,qs_ancien,zmasse,prsw_ancien)
6864    CALL water_int(klon,klev,qbs_ancien,zmasse,prbsw_ancien)
6865    ! !! RomP >>>
6866    IF (nqtot > nqo) tr_ancien(:,:,:) = tr_seri(:,:,:)
6867    ! !! RomP <<<
6868    !==========================================================================
6869    ! Sorties des tendances pour un point particulier
6870    ! a utiliser en 1D, avec igout=1 ou en 3D sur un point particulier
6871    ! pour le debug
6872    ! La valeur de igout est attribuee plus haut dans le programme
6873    !==========================================================================
6874
6875    IF (prt_level.ge.1) THEN
6876       write(lunout,*) 'FIN DE PHYSIQ !!!!!!!!!!!!!!!!!!!!'
6877       write(lunout,*) &
6878            'nlon,klev,nqtot,debut,lafin,jD_cur, jH_cur, pdtphys pct tlos'
6879       write(lunout,*) &
6880            nlon,klev,nqtot,debut,lafin, jD_cur, jH_cur ,pdtphys, &
6881            pctsrf(igout,is_ter), pctsrf(igout,is_lic),pctsrf(igout,is_oce), &
6882            pctsrf(igout,is_sic)
6883       write(lunout,*) 'd_t_dyn,d_t_con,d_t_lsc,d_t_ajsb,d_t_ajs,d_t_eva'
6884       DO k=1,klev
6885          write(lunout,*) d_t_dyn(igout,k),d_t_con(igout,k), &
6886               d_t_lsc(igout,k),d_t_ajsb(igout,k),d_t_ajs(igout,k), &
6887               d_t_eva(igout,k)
6888       ENDDO
6889       write(lunout,*) 'cool,heat'
6890       DO k=1,klev
6891          write(lunout,*) cool(igout,k),heat(igout,k)
6892       ENDDO
6893
6894       !jyg<     (En attendant de statuer sur le sort de d_t_oli)
6895       !jyg!     write(lunout,*) 'd_t_oli,d_t_vdf,d_t_oro,d_t_lif,d_t_ec'
6896       !jyg!     do k=1,klev
6897       !jyg!        write(lunout,*) d_t_oli(igout,k),d_t_vdf(igout,k), &
6898       !jyg!             d_t_oro(igout,k),d_t_lif(igout,k),d_t_ec(igout,k)
6899       !jyg!     enddo
6900       write(lunout,*) 'd_t_vdf,d_t_oro,d_t_lif,d_t_ec'
6901       DO k=1,klev
6902          write(lunout,*) d_t_vdf(igout,k), &
6903               d_t_oro(igout,k),d_t_lif(igout,k),d_t_ec(igout,k)
6904       ENDDO
6905       !>jyg
6906
6907       write(lunout,*) 'd_ps ',d_ps(igout)
6908       write(lunout,*) 'd_u, d_v, d_t, d_qx1, d_qx2 '
6909       DO k=1,klev
6910          write(lunout,*) d_u(igout,k),d_v(igout,k),d_t(igout,k), &
6911               d_qx(igout,k,1),d_qx(igout,k,2)
6912       ENDDO
6913    ENDIF
6914
6915    !============================================================
6916    !   Calcul de la temperature potentielle
6917    !============================================================
6918    DO k = 1, klev
6919       DO i = 1, klon
6920          !JYG/IM theta en debut du pas de temps
6921          !JYG/IM       theta(i,k)=t(i,k)*(100000./pplay(i,k))**(RD/RCPD)
6922          !JYG/IM theta en fin de pas de temps de physique
6923          theta(i,k)=t_seri(i,k)*(100000./pplay(i,k))**(RD/RCPD)
6924          ! thetal: 2 lignes suivantes a decommenter si vous avez les fichiers
6925          !     MPL 20130625
6926          ! fth_fonctions.F90 et parkind1.F90
6927          ! sinon thetal=theta
6928          !       thetal(i,k)=fth_thetal(pplay(i,k),t_seri(i,k),q_seri(i,k),
6929          !    :         ql_seri(i,k))
6930          thetal(i,k)=theta(i,k)
6931       ENDDO
6932    ENDDO
6933    !
6934
6935    ! 22.03.04 BEG
6936    !=============================================================
6937    !   Ecriture des sorties
6938    !=============================================================
6939#ifdef CPP_IOIPSL
6940
6941    ! Recupere des varibles calcule dans differents modules
6942    ! pour ecriture dans histxxx.nc
6943
6944    ! Get some variables from module fonte_neige_mod
6945    CALL fonte_neige_get_vars(pctsrf,  &
6946         zxfqcalving, zxfqfonte, zxffonte, zxrunofflic &
6947#ifdef ISO     
6948     &  ,zxfxtcalving, zxfxtfonte,zxxtrunofflic &
6949#endif     
6950     &  )
6951
6952
6953    !=============================================================
6954    ! Separation entre thermiques et non thermiques dans les sorties
6955    ! de fisrtilp
6956    !=============================================================
6957
6958    IF (iflag_thermals>=1) THEN
6959       d_t_lscth=0.
6960       d_t_lscst=0.
6961       d_q_lscth=0.
6962       d_q_lscst=0.
6963       DO k=1,klev
6964          DO i=1,klon
6965             IF (ptconvth(i,k)) THEN
6966                d_t_lscth(i,k)=d_t_eva(i,k)+d_t_lsc(i,k)
6967                d_q_lscth(i,k)=d_q_eva(i,k)+d_q_lsc(i,k)
6968             ELSE
6969                d_t_lscst(i,k)=d_t_eva(i,k)+d_t_lsc(i,k)
6970                d_q_lscst(i,k)=d_q_eva(i,k)+d_q_lsc(i,k)
6971             ENDIF
6972          ENDDO
6973       ENDDO
6974
6975       DO i=1,klon
6976          plul_st(i)=prfl(i,lmax_th(i)+1)+psfl(i,lmax_th(i)+1)
6977          plul_th(i)=prfl(i,1)+psfl(i,1)
6978       ENDDO
6979    ENDIF
6980
6981    !On effectue les sorties:
6982
6983#ifdef CPP_Dust
6984  CALL phys_output_write_spl(itap, pdtphys, paprs, pphis,  &
6985       pplay, lmax_th, aerosol_couple,                 &
6986       ok_ade, ok_aie, ivap, ok_sync,                  &
6987       ptconv, read_climoz, clevSTD,                   &
6988       ptconvth, d_t, qx, d_qx, d_tr_dyn, zmasse,      &
6989       flag_aerosol, flag_aerosol_strat, ok_cdnc)
6990#else
6991    CALL phys_output_write(itap, pdtphys, paprs, pphis,  &
6992         pplay, lmax_th, aerosol_couple,                 &
6993         ok_ade, ok_aie, ok_volcan, ivap, iliq, isol, ibs,   &
6994         ok_sync, ptconv, read_climoz, clevSTD,          &
6995         ptconvth, d_u, d_t, qx, d_qx, zmasse,           &
6996         flag_aerosol, flag_aerosol_strat, ok_cdnc,t,u1,v1)
6997#endif
6998
6999#ifndef CPP_XIOS
7000    CALL write_paramLMDZ_phy(itap,nid_ctesGCM,ok_sync)
7001#endif
7002
7003#endif
7004
7005! Pour XIOS : On remet des variables a .false. apres un premier appel
7006    IF (debut) THEN
7007#ifdef CPP_XIOS
7008      swaero_diag=.FALSE.
7009      swaerofree_diag=.FALSE.
7010      dryaod_diag=.FALSE.
7011      ok_4xCO2atm= .FALSE.
7012!      write (lunout,*)'ok_4xCO2atm= ',swaero_diag, swaerofree_diag, dryaod_diag, ok_4xCO2atm
7013
7014      IF (is_master) THEN
7015        !--setting up swaero_diag to TRUE in XIOS case
7016        IF (xios_field_is_active("topswad").OR.xios_field_is_active("topswad0").OR. &
7017           xios_field_is_active("solswad").OR.xios_field_is_active("solswad0").OR. &
7018           xios_field_is_active("topswai").OR.xios_field_is_active("solswai").OR.  &
7019             (iflag_rrtm==1.AND.(xios_field_is_active("toplwad").OR.xios_field_is_active("toplwad0").OR. &
7020                                 xios_field_is_active("sollwad").OR.xios_field_is_active("sollwad0"))))  &
7021           !!!--for now these fields are not in the XML files so they are omitted
7022           !!!  xios_field_is_active("toplwai").OR.xios_field_is_active("sollwai") !))) &
7023           swaero_diag=.TRUE.
7024
7025        !--setting up swaerofree_diag to TRUE in XIOS case
7026        IF (xios_field_is_active("SWdnSFCcleanclr").OR.xios_field_is_active("SWupSFCcleanclr").OR. &
7027           xios_field_is_active("SWupTOAcleanclr").OR.xios_field_is_active("rsucsaf").OR.   &
7028           xios_field_is_active("rsdcsaf") .OR. xios_field_is_active("LWdnSFCcleanclr").OR. &
7029           xios_field_is_active("LWupTOAcleanclr")) &
7030           swaerofree_diag=.TRUE.
7031
7032        !--setting up dryaod_diag to TRUE in XIOS case
7033        DO naero = 1, naero_tot-1
7034         IF (xios_field_is_active("dryod550_"//name_aero_tau(naero))) dryaod_diag=.TRUE.
7035        ENDDO
7036        !
7037        !--setting up ok_4xCO2atm to TRUE in XIOS case
7038        IF (xios_field_is_active("rsut4co2").OR.xios_field_is_active("rlut4co2").OR. &
7039           xios_field_is_active("rsutcs4co2").OR.xios_field_is_active("rlutcs4co2").OR. &
7040           xios_field_is_active("rsu4co2").OR.xios_field_is_active("rsucs4co2").OR. &
7041           xios_field_is_active("rsd4co2").OR.xios_field_is_active("rsdcs4co2").OR. &
7042           xios_field_is_active("rlu4co2").OR.xios_field_is_active("rlucs4co2").OR. &
7043           xios_field_is_active("rld4co2").OR.xios_field_is_active("rldcs4co2")) &
7044           ok_4xCO2atm=.TRUE.
7045      ENDIF
7046      !$OMP BARRIER
7047      CALL bcast(swaero_diag)
7048      CALL bcast(swaerofree_diag)
7049      CALL bcast(dryaod_diag)
7050      CALL bcast(ok_4xCO2atm)
7051!      write (lunout,*)'ok_4xCO2atm= ',swaero_diag, swaerofree_diag, dryaod_diag, ok_4xCO2atm
7052#endif
7053    ENDIF
7054
7055    !====================================================================
7056    ! Arret du modele apres hgardfou en cas de detection d'un
7057    ! plantage par hgardfou
7058    !====================================================================
7059
7060    IF (abortphy==1) THEN
7061       abort_message ='Plantage hgardfou'
7062       CALL abort_physic (modname,abort_message,1)
7063    ENDIF
7064
7065    ! 22.03.04 END
7066    !
7067    !====================================================================
7068    ! Si c'est la fin, il faut conserver l'etat de redemarrage
7069    !====================================================================
7070    !
7071#ifdef ISO   
7072#ifdef ISOVERIF
7073if (iso_HDO.gt.0) then
7074      call iso_verif_aberrant_enc_vect2D( &
7075     &           xt_ancien,q_ancien, &
7076     &           'physiq 6577',ntraciso,klon,klev)
7077endif
7078         write(*,*) 'physiq 3731: verif avant phyisoredem'   
7079         do k=1,klev
7080          do i=1,klon
7081            if (iso_eau.gt.0) then
7082               call iso_verif_egalite_choix(xt_ancien(iso_eau,i,k), &
7083     &           q_ancien(i,k),'physiq 3728: avant phyisoredem', &
7084     &           errmax,errmaxrel)
7085            endif ! if (iso_eau.gt.0) then
7086#ifdef ISOTRAC
7087            IF(nzone > 0) CALL iso_verif_traceur(xt_ancien(1,i,k),'physiq 4802')
7088#endif         
7089          enddo
7090         enddo !do k=1,klev
7091#endif         
7092#endif
7093! ISO
7094
7095    ! Disabling calls to the prt_alerte function
7096    alert_first_call = .FALSE.
7097   
7098    IF (lafin) THEN
7099       itau_phy = itau_phy + itap
7100       CALL phyredem ("restartphy.nc")
7101       !         open(97,form="unformatted",file="finbin")
7102       !         write(97) u_seri,v_seri,t_seri,q_seri
7103       !         close(97)
7104     
7105       IF (is_omp_master) THEN
7106       
7107         IF (read_climoz >= 1) THEN
7108           IF (is_mpi_root) CALL nf95_close(ncid_climoz)
7109            DEALLOCATE(press_edg_climoz) ! pointer
7110            DEALLOCATE(press_cen_climoz) ! pointer
7111         ENDIF
7112       
7113       ENDIF
7114#ifdef CPP_XIOS
7115       IF (is_omp_master) CALL xios_context_finalize
7116
7117#ifdef INCA
7118       if (type_trac == 'inca') then
7119          IF (is_omp_master .and. grid_type==unstructured) THEN
7120             CALL finalize_inca
7121          ENDIF
7122       endif
7123#endif
7124
7125#endif
7126       WRITE(lunout,*) ' physiq fin, nombre de steps ou cvpas = 1 : ', Ncvpaseq1
7127    ENDIF
7128
7129    !      first=.false.
7130
7131  END SUBROUTINE physiq
7132
7133END MODULE physiq_mod
Note: See TracBrowser for help on using the repository browser.