source: LMDZ6/trunk/libf/phylmd/physiq_mod.F90 @ 5813

Last change on this file since 5813 was 5813, checked in by musat, 2 months ago

Ajout iflag_tropo pour choisir le calcul de la tropopause
pour l'ozone:
iflag_tropo=0 : tropopause "dynamique"=> valeur par defaut pour

retrocompatibilite

iflag_tropo=1 : tropopause "lapse-rate" comme pour les aerosols

stratospheriques ==> valeur recommandee

Ionela Musat

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