source: LMDZ6/branches/Portage_acc/libf/phylmd/physiq_mod.F90

Last change on this file was 4743, checked in by Laurent Fairhead, 7 months ago

Merge of ACC branch with 4740 revision from trunk

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