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

Last change on this file since 4089 was 4089, checked in by fhourdin, 2 years ago

Reecriture des thermiques

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