source: LMDZ5/trunk/libf/phylmd/phys_output_write_mod.F90 @ 2469

Last change on this file since 2469 was 2429, checked in by Laurent Fairhead, 8 years ago

Correction on the calculation of the surface of the grid at the poles (problem was introduced
in r2222).
Due to the different distribution of OMP tasks in the dynamics and the physics, had to
introduce 2 new logical variables, is_pole_north_phy and is_pole_south_phy, and so decided
to rename the old is_north_pole/is_south_pole to is_north_pole_dyn/is_south_pole_dyn to
stay coherent and, hopefully, clear things up a bit.
LF

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:keywords set to Id
File size: 67.1 KB
Line 
1!
2! $Id: phys_output_write_mod.F90 2429 2016-01-27 12:43:09Z lguez $
3!
4MODULE phys_output_write_mod
5
6  USE phytrac_mod, ONLY : d_tr_cl, d_tr_th, d_tr_cv, d_tr_lessi_impa, &
7       d_tr_lessi_nucl, d_tr_insc, d_tr_bcscav, d_tr_evapls, d_tr_ls,  &
8       d_tr_trsp, d_tr_sscav, d_tr_sat, d_tr_uscav
9
10  ! Author: Abderrahmane IDELKADI (original include file)
11  ! Author: Laurent FAIRHEAD (transformation to module/subroutine)
12  ! Author: Ulysse GERARD (effective implementation)
13
14CONTAINS
15
16  ! ug Routine pour définir (los du premier passageà) ET sortir les variables
17  SUBROUTINE phys_output_write(itap, pdtphys, paprs, pphis, &
18       pplay, lmax_th, aerosol_couple,         &
19       ok_ade, ok_aie, ivap, new_aod, ok_sync, &
20       ptconv, read_climoz, clevSTD, ptconvth, &
21       d_t, qx, d_qx, zmasse, flag_aerosol, flag_aerosol_strat, ok_cdnc)
22
23    ! This subroutine does the actual writing of diagnostics that were
24    ! defined and initialised in phys_output_mod.F90
25
26    USE dimphy, only: klon, klev, klevp1, nslay
27    USE mod_phys_lmdz_para, ONLY: is_north_pole_phy,is_south_pole_phy
28    USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat
29    USE time_phylmdz_mod, only: day_step_phy, start_time, itau_phy
30    USE phys_output_ctrlout_mod, only: o_phis, o_aire, is_ter, is_lic, is_oce, &
31         is_ave, is_sic, o_contfracATM, o_contfracOR, &
32         o_aireTER, o_flat, o_slp, o_ptstar, o_pt0, o_tsol, &
33         o_t2m, o_t2m_min, o_t2m_max, &
34         o_t2m_min_mon, o_t2m_max_mon, &
35         o_q2m, o_ustar, o_u10m, o_v10m, &
36         o_wind10m, o_wind10max, o_gusts, o_sicf, &
37         o_psol, o_mass, o_qsurf, o_qsol, &
38         o_precip, o_ndayrain, o_plul, o_pluc, &
39         o_snow, o_msnow, o_fsnow, o_evap, &
40         o_tops, o_tops0, o_topl, o_topl0, &
41         o_SWupTOA, o_SWupTOAclr, o_SWdnTOA, &
42         o_SWdnTOAclr, o_nettop, o_SWup200, &
43         o_SWup200clr, o_SWdn200, o_SWdn200clr, &
44         o_LWup200, o_LWup200clr, o_LWdn200, &
45         o_LWdn200clr, o_sols, o_sols0, &
46         o_soll, o_radsol, o_soll0, o_SWupSFC, &
47         o_SWupSFCclr, o_SWdnSFC, o_SWdnSFCclr, &
48         o_LWupSFC, o_LWdnSFC, o_LWupSFCclr, &
49         o_LWdnSFCclr, o_bils, o_bils_diss, &
50         o_bils_ec,o_bils_ech, o_bils_tke, o_bils_kinetic, &
51         o_bils_latent, o_bils_enthalp, o_sens, &
52         o_fder, o_ffonte, o_fqcalving, o_fqfonte, &
53         o_taux, o_tauy, o_snowsrf, o_qsnow, &
54         o_snowhgt, o_toice, o_sissnow, o_runoff, &
55         o_albslw3, o_pourc_srf, o_fract_srf, &
56         o_taux_srf, o_tauy_srf, o_tsol_srf, &
57         o_evappot_srf, o_ustar_srf, o_u10m_srf, &
58         o_v10m_srf, o_t2m_srf, o_evap_srf, &
59         o_sens_srf, o_lat_srf, o_flw_srf, &
60         o_fsw_srf, o_wbils_srf, o_wbilo_srf, &
61         o_tke_srf, o_tke_max_srf,o_dltpbltke_srf, o_wstar, &
62         o_cdrm, o_cdrh, o_cldl, o_cldm, o_cldh, &
63         o_cldt, o_JrNt, o_cldljn, o_cldmjn, &
64         o_cldhjn, o_cldtjn, o_cldq, o_lwp, o_iwp, &
65         o_ue, o_ve, o_uq, o_vq, o_cape, o_pbase, &
66         o_ptop, o_fbase, o_plcl, o_plfc, &
67         o_wbeff, o_cape_max, o_upwd, o_Ma, &
68         o_dnwd, o_dnwd0, o_ftime_con, o_mc, &
69         o_prw, o_s_pblh, o_s_pblt, o_s_lcl, &
70         o_s_therm, o_uSTDlevs, o_vSTDlevs, &
71         o_wSTDlevs, o_zSTDlevs, o_qSTDlevs, &
72         o_tSTDlevs, epsfra, o_t_oce_sic, &
73         o_ale_bl, o_alp_bl, o_ale_wk, o_alp_wk, &
74         o_dtvdf_x    , o_dtvdf_w    , o_dqvdf_x    , o_dqvdf_w    , &
75         o_sens_x     , o_sens_w     , o_flat_x     , o_flat_w     , &
76         o_delta_tsurf, &
77         o_cdragh_x   , o_cdragh_w   , o_cdragm_x   , o_cdragm_w   , &
78         o_kh         , o_kh_x       , o_kh_w       , &
79         o_ale, o_alp, o_cin, o_WAPE, o_wake_h, &
80         o_wake_s, o_wake_deltat, o_wake_deltaq, &
81         o_wake_omg, o_dtwak, o_dqwak, o_Vprecip, &
82         o_ftd, o_fqd, o_wdtrainA, o_wdtrainM, &
83         o_n2, o_s2, o_proba_notrig, &
84         o_random_notrig, o_ale_bl_stat, &
85         o_ale_bl_trig, o_alp_bl_det, &
86         o_alp_bl_fluct_m, o_alp_bl_fluct_tke, &
87         o_alp_bl_conv, o_alp_bl_stat, &
88         o_slab_qflux, o_tslab, o_slab_bils, &
89         o_slab_bilg, o_slab_sic, o_slab_tice, &
90         o_weakinv, o_dthmin, o_cldtau, &
91         o_cldemi, o_pr_con_l, o_pr_con_i, &
92         o_pr_lsc_l, o_pr_lsc_i, o_re, o_fl, &
93         o_rh2m, o_rh2m_min, o_rh2m_max, &
94         o_qsat2m, o_tpot, o_tpote, o_SWnetOR, &
95         o_SWdownOR, o_LWdownOR, o_snowl, &
96         o_solldown, o_dtsvdfo, o_dtsvdft, &
97         o_dtsvdfg, o_dtsvdfi, o_z0m, o_z0h, o_od550aer, &
98         o_od865aer, o_absvisaer, o_od550lt1aer, &
99         o_sconcso4, o_sconcno3, o_sconcoa, o_sconcbc, &
100         o_sconcss, o_sconcdust, o_concso4, o_concno3, &
101         o_concoa, o_concbc, o_concss, o_concdust, &
102         o_loadso4, o_loadoa, o_loadbc, o_loadss, &
103         o_loaddust, o_tausumaero, o_tausumaero_lw, &
104         o_topswad, o_topswad0, o_solswad, o_solswad0, &
105         o_toplwad, o_toplwad0, o_sollwad, o_sollwad0, &
106         o_swtoaas_nat, o_swsrfas_nat, &
107         o_swtoacs_nat, o_swtoaas_ant, &
108         o_swsrfas_ant, o_swtoacs_ant, &
109         o_swsrfcs_ant, o_swtoacf_nat, &
110         o_swsrfcf_nat, o_swtoacf_ant, &
111         o_swsrfcs_nat, o_swsrfcf_ant, &
112         o_swtoacf_zero, o_swsrfcf_zero, &
113         o_topswai, o_solswai, o_scdnc, &
114         o_cldncl, o_reffclws, o_reffclwc, &
115         o_cldnvi, o_lcc, o_lcc3d, o_lcc3dcon, &
116         o_lcc3dstra, o_reffclwtop, o_ec550aer, &
117         o_lwcon, o_iwcon, o_temp, o_theta, &
118         o_ovapinit, o_ovap, o_oliq, o_geop, &
119         o_vitu, o_vitv, o_vitw, o_pres, o_paprs, &
120         o_zfull, o_zhalf, o_rneb, o_rnebjn, o_rnebcon, &
121         o_rnebls, o_rhum, o_ozone, o_ozone_light, &
122         o_dtphy, o_dqphy, o_albe_srf, o_z0m_srf, o_z0h_srf, &
123         o_ages_srf, o_snow_srf, o_alb1, o_alb2, o_tke, &
124         o_tke_max, o_kz, o_kz_max, o_clwcon, &
125         o_dtdyn, o_dqdyn, o_dudyn, o_dvdyn, &
126         o_dtcon, o_tntc, o_ducon, o_dvcon, &
127         o_dqcon, o_tnhusc, o_tnhusc, o_dtlsc, &
128         o_dtlschr, o_dqlsc, o_beta_prec, &
129         o_dtlscth, o_dtlscst, o_dqlscth, &
130         o_dqlscst, o_plulth, o_plulst, &
131         o_ptconvth, o_lmaxth, o_dtvdf, &
132         o_dtdis, o_dqvdf, o_dteva, o_dqeva, &
133         o_ptconv, o_ratqs, o_dtthe, &
134         o_duthe, o_dvthe, o_ftime_th, &
135         o_f_th, o_e_th, o_w_th, o_q_th, &
136         o_a_th, o_d_th, o_f0_th, o_zmax_th, &
137         o_dqthe, o_dtajs, o_dqajs, o_dtswr, &
138         o_dtsw0, o_dtlwr, o_dtlw0, o_dtec, &
139         o_duvdf, o_dvvdf, o_duoro, o_dvoro, &
140         o_dtoro, o_dulif, o_dvlif, o_dtlif, &
141         o_du_gwd_hines, o_dv_gwd_hines, o_dthin, o_dqch4, o_rsu, &
142         o_du_gwd_front, o_dv_gwd_front, &
143         o_east_gwstress, o_west_gwstress, &
144         o_rsd, o_rlu, o_rld, o_rsucs, o_rsdcs, &
145         o_rlucs, o_rldcs, o_tnt, o_tntr, &
146         o_tntscpbl, o_tnhus, o_tnhusscpbl, &
147         o_evu, o_h2o, o_mcd, o_dmc, o_ref_liq, &
148         o_ref_ice, o_rsut4co2, o_rlut4co2, &
149         o_rsutcs4co2, o_rlutcs4co2, o_rsu4co2, &
150         o_rlu4co2, o_rsucs4co2, o_rlucs4co2, &
151         o_rsd4co2, o_rld4co2, o_rsdcs4co2, &
152         o_rldcs4co2, o_tnondef, o_ta, o_zg, &
153         o_hus, o_hur, o_ua, o_va, o_wap, &
154         o_psbg, o_tro3, o_tro3_daylight, &
155         o_uxv, o_vxq, o_vxT, o_wxq, o_vxphi, &
156         o_wxT, o_uxu, o_vxv, o_TxT, o_trac, &
157         o_dtr_vdf, o_dtr_the, o_dtr_con, &
158         o_dtr_lessi_impa, o_dtr_lessi_nucl, &
159         o_dtr_insc, o_dtr_bcscav, o_dtr_evapls, &
160         o_dtr_ls, o_dtr_trsp, o_dtr_sscav, &
161         o_dtr_sat, o_dtr_uscav, o_trac_cum, o_du_gwd_rando, o_dv_gwd_rando, &
162         o_ustr_gwd_hines,o_vstr_gwd_hines,o_ustr_gwd_rando,o_vstr_gwd_rando, &
163         o_ustr_gwd_front,o_vstr_gwd_front
164
165    USE phys_state_var_mod, only: pctsrf, paire_ter, rain_fall, snow_fall, &
166         qsol, z0m, z0h, fevap, agesno, &
167         nday_rain, rain_con, snow_con, &
168         topsw, toplw, toplw0, swup, swdn, &
169         topsw0, swup0, swdn0, SWup200, SWup200clr, &
170         SWdn200, SWdn200clr, LWup200, LWup200clr, &
171         LWdn200, LWdn200clr, solsw, solsw0, sollw, &
172         radsol, swradcorr, sollw0, sollwdown, sollw, gustiness, &
173         sollwdownclr, lwdn0, ftsol, ustar, u10m, &
174         v10m, pbl_tke, wake_delta_pbl_TKE, &
175         delta_tsurf, &
176         wstar, cape, ema_pcb, ema_pct, &
177         ema_cbmf, Ma, fm_therm, ale_bl, alp_bl, ale, &
178         alp, cin, wake_pe, wake_s, wake_deltat, &
179         wake_deltaq, ftd, fqd, ale_bl_trig, albsol1, &
180         rnebcon, wo, falb1, albsol2, coefh, clwcon0, &
181         ratqs, entr_therm, zqasc, detr_therm, f0, &
182         lwup, lwdn, lwup0, coefm, &
183         swupp, lwupp, swup0p, lwup0p, swdnp, lwdnp, &
184         swdn0p, lwdn0p, tnondef, O3sumSTD, uvsumSTD, &
185         vqsumSTD, vTsumSTD, O3daysumSTD, wqsumSTD, &
186         vphisumSTD, wTsumSTD, u2sumSTD, v2sumSTD, &
187         T2sumSTD, nlevSTD, du_gwd_rando, du_gwd_front, &
188         ulevSTD, vlevSTD, wlevSTD, philevSTD, qlevSTD, tlevSTD, &
189         rhlevSTD, O3STD, O3daySTD, uvSTD, vqSTD, vTSTD, wqSTD, &
190         vphiSTD, wTSTD, u2STD, v2STD, T2STD, missing_val_nf90
191
192    USE phys_local_var_mod, only: zxfluxlat, slp, ptstar, pt0, zxtsol, zt2m, &
193         t2m_min_mon, t2m_max_mon, evap, &
194         zu10m, zv10m, zq2m, zustar, zxqsurf, &
195         rain_lsc, snow_lsc, bils, sens, fder, &
196         zxffonte, zxfqcalving, zxfqfonte, fluxu, &
197         fluxv, zxsnow, qsnow, snowhgt, to_ice, &
198         sissnow, runoff, albsol3_lic, evap_pot, &
199         t2m, fluxt, fluxlat, fsollw, fsolsw, &
200         wfbils, wfbilo, cdragm, cdragh, cldl, cldm, &
201         cldh, cldt, JrNt, cldljn, cldmjn, cldhjn, &
202         cldtjn, cldq, flwp, fiwp, ue, ve, uq, vq, &
203         plcl, plfc, wbeff, upwd, dnwd, dnwd0, prw, &
204         s_pblh, s_pblt, s_lcl, s_therm, uwriteSTD, &
205         vwriteSTD, wwriteSTD, phiwriteSTD, qwriteSTD, &
206         twriteSTD, ale_wake, alp_wake, &
207         dtvdf_x    ,dtvdf_w    ,dqvdf_x    ,dqvdf_w    , &
208         sens_x     ,sens_w     ,zxfluxlat_x,zxfluxlat_w, &
209         cdragh_x   ,cdragh_w   ,cdragm_x   ,cdragm_w   , &
210         kh         ,kh_x       ,kh_w       , &
211         wake_h, &
212         wake_omg, d_t_wake, d_q_wake, Vprecip, &
213         wdtrainA, wdtrainM, n2, s2, proba_notrig, &
214         random_notrig, ale_bl_stat, &
215         alp_bl_det, alp_bl_fluct_m, alp_bl_conv, &
216         alp_bl_stat, alp_bl_fluct_tke, slab_wfbils, &
217         weak_inversion, dthmin, cldtau, cldemi, &
218         pmflxr, pmflxs, prfl, psfl, re, fl, rh2m, &
219         qsat2m, tpote, tpot, d_ts, od550aer, &
220         od865aer, absvisaer, od550lt1aer, sconcso4, sconcno3, &
221         sconcoa, sconcbc, sconcss, sconcdust, concso4, concno3, &
222         concoa, concbc, concss, concdust, loadso4, &
223         loadoa, loadbc, loadss, loaddust, tausum_aero, &
224         topswad_aero, topswad0_aero, solswad_aero, &
225         solswad0_aero, topsw_aero, solsw_aero, &
226         topsw0_aero, solsw0_aero, topswcf_aero, &
227         solswcf_aero, topswai_aero, solswai_aero, &
228         toplwad_aero, toplwad0_aero, sollwad_aero, &
229         sollwad0_aero, toplwai_aero, sollwai_aero, &
230         scdnc, cldncl, reffclws, reffclwc, cldnvi, &
231         lcc, lcc3d, lcc3dcon, lcc3dstra, reffclwtop, &
232         ec550aer, flwc, fiwc, t_seri, theta, q_seri, &
233!jyg<
234!!         ql_seri, zphi, u_seri, v_seri, omega, cldfra, &
235         ql_seri, tr_seri, &
236         zphi, u_seri, v_seri, omega, cldfra, &
237!>jyg
238         rneb, rnebjn, zx_rh, d_t_dyn, d_q_dyn, &
239         d_u_dyn, d_v_dyn, d_t_con, d_t_ajsb, d_t_ajs, &
240         d_u_ajs, d_v_ajs, &
241         d_u_con, d_v_con, d_q_con, d_q_ajs, d_t_lsc, &
242         d_t_lwr,d_t_lw0,d_t_swr,d_t_sw0, &
243         d_t_eva, d_q_lsc, beta_prec, d_t_lscth, &
244         d_t_lscst, d_q_lscth, d_q_lscst, plul_th, &
245         plul_st, d_t_vdf, d_t_diss, d_q_vdf, d_q_eva, &
246         zw2, fraca, zmax_th, d_q_ajsb, d_t_ec, d_u_vdf, &
247         d_v_vdf, d_u_oro, d_v_oro, d_t_oro, d_u_lif, &
248         d_v_lif, d_t_lif, du_gwd_hines, dv_gwd_hines, d_t_hin, &
249         dv_gwd_rando, dv_gwd_front, &
250         east_gwstress, west_gwstress, &
251         d_q_ch4, pmfd, pmfu, ref_liq, ref_ice, rhwriteSTD
252
253    USE phys_output_var_mod, only: vars_defined, snow_o, zfra_o, bils_diss, &
254         bils_ec,bils_ech, bils_tke, bils_kinetic, bils_latent, bils_enthalp, &
255         itau_con, nfiles, clef_files, nid_files, &
256         zustr_gwd_hines, zvstr_gwd_hines,zustr_gwd_rando, zvstr_gwd_rando, &
257         zustr_gwd_front, zvstr_gwd_front                                   
258    USE ocean_slab_mod, only: tslab, slab_bils, slab_bilg, tice, seaice
259    USE pbl_surface_mod, only: snow
260    USE indice_sol_mod, only: nbsrf
261    USE infotrac_phy, only: nqtot, nqo, type_trac
262    USE geometry_mod, only: cell_area
263    USE surface_data, only: type_ocean, version_ocean, ok_veget, ok_snow
264!    USE aero_mod, only: naero_spc
265    USE aero_mod, only: naero_tot, id_STRAT_phy
266    USE ioipsl, only: histend, histsync
267    USE iophy, only: set_itau_iophy, histwrite_phy
268    USE netcdf, only: nf90_fill_real
269    USE print_control_mod, ONLY: prt_level,lunout
270
271
272#ifdef CPP_XIOS
273    ! ug Pour les sorties XIOS
274    USE xios, ONLY: xios_update_calendar
275    USE wxios, only: wxios_closedef, missing_val
276#endif
277    USE phys_cal_mod, only : mth_len
278
279
280    IMPLICIT NONE
281
282
283    INCLUDE "clesphys.h"
284    INCLUDE "thermcell.h"
285    INCLUDE "compbl.h"
286    INCLUDE "YOMCST.h"
287
288    ! Input
289    INTEGER :: itap, ivap, read_climoz
290    INTEGER, DIMENSION(klon) :: lmax_th
291    LOGICAL :: aerosol_couple, ok_sync
292    LOGICAL :: ok_ade, ok_aie, new_aod
293    LOGICAL, DIMENSION(klon, klev) :: ptconv, ptconvth
294    REAL :: pdtphys
295    CHARACTER (LEN=4), DIMENSION(nlevSTD) :: clevSTD
296    REAL, DIMENSION(klon,nlevSTD) :: zx_tmp_fi3d_STD
297    REAL, DIMENSION(klon) :: pphis
298    REAL, DIMENSION(klon, klev) :: pplay, d_t
299    REAL, DIMENSION(klon, klev+1) :: paprs
300    REAL, DIMENSION(klon,klev,nqtot) :: qx, d_qx
301    REAL, DIMENSION(klon, klev) :: zmasse
302    LOGICAL :: flag_aerosol_strat
303    INTEGER :: flag_aerosol
304    LOGICAL :: ok_cdnc
305    REAL, DIMENSION(3) :: freq_moyNMC
306
307    ! Local
308    INTEGER :: itau_w
309    INTEGER :: i, iinit, iinitend=1, iff, iq, nsrf, k, ll, naero
310    REAL, DIMENSION (klon) :: zx_tmp_fi2d
311    REAL, DIMENSION (klon,klev) :: zx_tmp_fi3d, zpt_conv
312    REAL, DIMENSION (klon,klev+1) :: zx_tmp_fi3d1
313    CHARACTER (LEN=4)              :: bb2
314    INTEGER, DIMENSION(nbp_lon*nbp_lat)  :: ndex2d
315    INTEGER, DIMENSION(nbp_lon*nbp_lat*klev) :: ndex3d
316    REAL, PARAMETER :: dobson_u = 2.1415e-05 ! Dobson unit, in kg m-2
317!   REAL, PARAMETER :: missing_val=nf90_fill_real
318#ifndef CPP_XIOS
319    REAL :: missing_val
320#endif
321    REAL, PARAMETER :: un_jour=86400.
322
323    ! On calcul le nouveau tau:
324    itau_w = itau_phy + itap + start_time * day_step_phy
325    ! On le donne à iophy pour que les histwrite y aient accès:
326    CALL set_itau_iophy(itau_w)
327
328    IF(.NOT.vars_defined) THEN
329       iinitend = 2
330    ELSE
331       iinitend = 1
332    ENDIF
333
334    ! ug la boucle qui suit ne sert qu'une fois, pour l'initialisation, sinon il n'y a toujours qu'un seul passage:
335    DO iinit=1, iinitend
336#ifdef CPP_XIOS
337       !$OMP MASTER
338       IF (vars_defined) THEN
339          if (prt_level >= 10) then
340             write(lunout,*)"phys_output_write: call xios_update_calendar, itau_w=",itau_w
341          endif
342!          CALL xios_update_calendar(itau_w)
343          CALL xios_update_calendar(itap)
344       END IF
345       !$OMP END MASTER
346       !$OMP BARRIER
347#endif
348       ! On procède à l'écriture ou à la définition des nombreuses variables:
349!!! Champs 1D !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
350       CALL histwrite_phy(o_phis, pphis)
351
352       zx_tmp_fi2d = cell_area
353       if (is_north_pole_phy) then
354         zx_tmp_fi2d(1) = cell_area(1)/nbp_lon
355       endif
356       if (is_south_pole_phy) then
357         zx_tmp_fi2d(klon) = cell_area(klon)/nbp_lon
358       endif
359       CALL histwrite_phy(o_aire, zx_tmp_fi2d)
360
361       IF (vars_defined) THEN
362          DO i=1, klon
363             zx_tmp_fi2d(i)=pctsrf(i,is_ter)+pctsrf(i,is_lic)
364          ENDDO
365       ENDIF
366
367       CALL histwrite_phy(o_contfracATM, zx_tmp_fi2d)
368       CALL histwrite_phy(o_contfracOR, pctsrf(:,is_ter))
369       CALL histwrite_phy(o_aireTER, paire_ter)
370!!! Champs 2D !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
371       CALL histwrite_phy(o_flat, zxfluxlat)
372       CALL histwrite_phy(o_ptstar, ptstar)
373       CALL histwrite_phy(o_pt0, pt0)
374       CALL histwrite_phy(o_slp, slp)
375       CALL histwrite_phy(o_tsol, zxtsol)
376       CALL histwrite_phy(o_t2m, zt2m)
377       CALL histwrite_phy(o_t2m_min, zt2m)
378       CALL histwrite_phy(o_t2m_max, zt2m)
379       CALL histwrite_phy(o_t2m_max_mon, t2m_max_mon)
380       CALL histwrite_phy(o_t2m_min_mon, t2m_min_mon)
381
382       IF (vars_defined) THEN
383          DO i=1, klon
384             zx_tmp_fi2d(i)=SQRT(zu10m(i)*zu10m(i)+zv10m(i)*zv10m(i))
385          ENDDO
386       ENDIF
387       CALL histwrite_phy(o_wind10m, zx_tmp_fi2d)
388
389       IF (vars_defined) THEN
390          DO i=1, klon
391             zx_tmp_fi2d(i)=SQRT(zu10m(i)*zu10m(i)+zv10m(i)*zv10m(i))
392          ENDDO
393       ENDIF
394       CALL histwrite_phy(o_wind10max, zx_tmp_fi2d)
395
396       CALL histwrite_phy(o_gusts, gustiness)
397
398       IF (vars_defined) THEN
399          DO i = 1, klon
400             zx_tmp_fi2d(i) = pctsrf(i,is_sic)
401          ENDDO
402       ENDIF
403       CALL histwrite_phy(o_sicf, zx_tmp_fi2d)
404       CALL histwrite_phy(o_q2m, zq2m)
405       CALL histwrite_phy(o_ustar, zustar)
406       CALL histwrite_phy(o_u10m, zu10m)
407       CALL histwrite_phy(o_v10m, zv10m)
408
409       IF (vars_defined) THEN
410          DO i = 1, klon
411             zx_tmp_fi2d(i) = paprs(i,1)
412          ENDDO
413       ENDIF
414       CALL histwrite_phy(o_psol, zx_tmp_fi2d)
415       CALL histwrite_phy(o_mass, zmasse)
416       CALL histwrite_phy(o_qsurf, zxqsurf)
417
418       IF (.NOT. ok_veget) THEN
419          CALL histwrite_phy(o_qsol, qsol)
420       ENDIF
421
422       IF (vars_defined) THEN
423          DO i = 1, klon
424             zx_tmp_fi2d(i) = rain_fall(i) + snow_fall(i)
425          ENDDO
426       ENDIF
427
428       CALL histwrite_phy(o_precip, zx_tmp_fi2d)
429       CALL histwrite_phy(o_ndayrain, nday_rain)
430
431       IF (vars_defined) THEN
432          DO i = 1, klon
433             zx_tmp_fi2d(i) = rain_lsc(i) + snow_lsc(i)
434          ENDDO
435       ENDIF
436       CALL histwrite_phy(o_plul, zx_tmp_fi2d)
437
438       IF (vars_defined) THEN
439          DO i = 1, klon
440             zx_tmp_fi2d(i) = rain_con(i) + snow_con(i)
441          ENDDO
442       ENDIF
443       CALL histwrite_phy(o_pluc, zx_tmp_fi2d)
444       CALL histwrite_phy(o_snow, snow_fall)
445       CALL histwrite_phy(o_msnow, zxsnow)
446       CALL histwrite_phy(o_fsnow, zfra_o)
447       CALL histwrite_phy(o_evap, evap)
448       CALL histwrite_phy(o_tops, topsw*swradcorr)
449       CALL histwrite_phy(o_tops0, topsw0*swradcorr)
450       CALL histwrite_phy(o_topl, toplw)
451       CALL histwrite_phy(o_topl0, toplw0)
452
453       IF (vars_defined) THEN
454          zx_tmp_fi2d(:) = swup(:,klevp1)*swradcorr(:)
455       ENDIF
456       CALL histwrite_phy(o_SWupTOA, zx_tmp_fi2d)
457
458       IF (vars_defined) THEN
459          zx_tmp_fi2d(:) = swup0(:,klevp1)*swradcorr(:)
460       ENDIF
461       CALL histwrite_phy(o_SWupTOAclr, zx_tmp_fi2d)
462
463       IF (vars_defined) THEN
464          zx_tmp_fi2d(:) = swdn(:,klevp1)*swradcorr(:)
465       ENDIF
466       CALL histwrite_phy(o_SWdnTOA, zx_tmp_fi2d)
467
468       IF (vars_defined) THEN
469          zx_tmp_fi2d(:) = swdn0(:,klevp1)*swradcorr(:)
470       ENDIF
471       CALL histwrite_phy(o_SWdnTOAclr, zx_tmp_fi2d)
472
473       IF (vars_defined) THEN
474          zx_tmp_fi2d(:) = topsw(:)*swradcorr(:)-toplw(:)
475       ENDIF
476       CALL histwrite_phy(o_nettop, zx_tmp_fi2d)
477       CALL histwrite_phy(o_SWup200, SWup200*swradcorr)
478       CALL histwrite_phy(o_SWup200clr, SWup200clr*swradcorr)
479       CALL histwrite_phy(o_SWdn200, SWdn200*swradcorr)
480       CALL histwrite_phy(o_SWdn200clr, SWdn200clr*swradcorr)
481       CALL histwrite_phy(o_LWup200, LWup200)
482       CALL histwrite_phy(o_LWup200clr, LWup200clr)
483       CALL histwrite_phy(o_LWdn200, LWdn200)
484       CALL histwrite_phy(o_LWdn200clr, LWdn200clr)
485       CALL histwrite_phy(o_sols, solsw*swradcorr)
486       CALL histwrite_phy(o_sols0, solsw0*swradcorr)
487       CALL histwrite_phy(o_soll, sollw)
488       CALL histwrite_phy(o_soll0, sollw0)
489       CALL histwrite_phy(o_radsol, radsol)
490
491       IF (vars_defined) THEN
492          zx_tmp_fi2d(:) = swup(:,1)*swradcorr(:)
493       ENDIF
494       CALL histwrite_phy(o_SWupSFC, zx_tmp_fi2d)
495
496       IF (vars_defined) THEN
497          zx_tmp_fi2d(:) = swup0(:,1)*swradcorr(:)
498       ENDIF
499       CALL histwrite_phy(o_SWupSFCclr, zx_tmp_fi2d)
500
501       IF (vars_defined) THEN
502          zx_tmp_fi2d(:) = swdn(:,1)*swradcorr(:)
503       ENDIF
504       CALL histwrite_phy(o_SWdnSFC, zx_tmp_fi2d)
505
506       IF (vars_defined) THEN
507          zx_tmp_fi2d(:) = swdn0(:,1)*swradcorr(:)
508       ENDIF
509       CALL histwrite_phy(o_SWdnSFCclr, zx_tmp_fi2d)
510
511       IF (vars_defined) THEN
512          zx_tmp_fi2d(:)=sollwdown(:)-sollw(:)
513       ENDIF
514       CALL histwrite_phy(o_LWupSFC, zx_tmp_fi2d)
515       CALL histwrite_phy(o_LWdnSFC, sollwdown)
516
517       IF (vars_defined) THEN
518          sollwdownclr(1:klon) = -1.*lwdn0(1:klon,1)
519          zx_tmp_fi2d(1:klon)=sollwdownclr(1:klon)-sollw0(1:klon)
520       ENDIF
521       CALL histwrite_phy(o_LWupSFCclr, zx_tmp_fi2d)
522       CALL histwrite_phy(o_LWdnSFCclr, sollwdownclr)
523       CALL histwrite_phy(o_bils, bils)
524       CALL histwrite_phy(o_bils_diss, bils_diss)
525       CALL histwrite_phy(o_bils_ec, bils_ec)
526       IF (iflag_ener_conserv>=1) THEN
527         CALL histwrite_phy(o_bils_ech, bils_ech)
528       ENDIF
529       CALL histwrite_phy(o_bils_tke, bils_tke)
530       CALL histwrite_phy(o_bils_kinetic, bils_kinetic)
531       CALL histwrite_phy(o_bils_latent, bils_latent)
532       CALL histwrite_phy(o_bils_enthalp, bils_enthalp)
533
534       IF (vars_defined) THEN
535          zx_tmp_fi2d(1:klon)=-1*sens(1:klon)
536       ENDIF
537       CALL histwrite_phy(o_sens, zx_tmp_fi2d)
538       CALL histwrite_phy(o_fder, fder)
539       CALL histwrite_phy(o_ffonte, zxffonte)
540       CALL histwrite_phy(o_fqcalving, zxfqcalving)
541       CALL histwrite_phy(o_fqfonte, zxfqfonte)
542       IF (vars_defined) THEN
543          zx_tmp_fi2d=0.
544          DO nsrf=1,nbsrf
545             zx_tmp_fi2d(:)=zx_tmp_fi2d(:)+pctsrf(:,nsrf)*fluxu(:,1,nsrf)
546          ENDDO
547       ENDIF
548       CALL histwrite_phy(o_taux, zx_tmp_fi2d)
549
550       IF (vars_defined) THEN
551          zx_tmp_fi2d=0.
552          DO nsrf=1,nbsrf
553             zx_tmp_fi2d(:)=zx_tmp_fi2d(:)+pctsrf(:,nsrf)*fluxv(:,1,nsrf)
554          ENDDO
555       ENDIF
556       CALL histwrite_phy(o_tauy, zx_tmp_fi2d)
557
558       IF (ok_snow) THEN
559          CALL histwrite_phy(o_snowsrf, snow_o)
560          CALL histwrite_phy(o_qsnow, qsnow)
561          CALL histwrite_phy(o_snowhgt,snowhgt)
562          CALL histwrite_phy(o_toice,to_ice)
563          CALL histwrite_phy(o_sissnow,sissnow)
564          CALL histwrite_phy(o_runoff,runoff)
565          CALL histwrite_phy(o_albslw3,albsol3_lic)
566       ENDIF
567
568       DO nsrf = 1, nbsrf
569          IF (vars_defined)             zx_tmp_fi2d(1 : klon) = pctsrf( 1 : klon, nsrf)*100.
570          CALL histwrite_phy(o_pourc_srf(nsrf), zx_tmp_fi2d)
571          IF (vars_defined)           zx_tmp_fi2d(1 : klon) = pctsrf( 1 : klon, nsrf)
572          CALL histwrite_phy(o_fract_srf(nsrf), zx_tmp_fi2d)
573          IF (vars_defined)         zx_tmp_fi2d(1 : klon) = fluxu( 1 : klon, 1, nsrf)
574          CALL histwrite_phy(o_taux_srf(nsrf), zx_tmp_fi2d)
575          IF (vars_defined)         zx_tmp_fi2d(1 : klon) = fluxv( 1 : klon, 1, nsrf)
576          CALL histwrite_phy(o_tauy_srf(nsrf), zx_tmp_fi2d)
577          IF (vars_defined)         zx_tmp_fi2d(1 : klon) = ftsol( 1 : klon, nsrf)
578          CALL histwrite_phy(o_tsol_srf(nsrf), zx_tmp_fi2d)
579          IF (vars_defined)         zx_tmp_fi2d(1 : klon) = evap_pot( 1 : klon, nsrf)
580          CALL histwrite_phy(o_evappot_srf(nsrf), zx_tmp_fi2d)
581          IF (vars_defined)       zx_tmp_fi2d(1 : klon) = ustar(1 : klon, nsrf)
582          CALL histwrite_phy(o_ustar_srf(nsrf), zx_tmp_fi2d)
583          IF (vars_defined)       zx_tmp_fi2d(1 : klon) = u10m(1 : klon, nsrf)
584          CALL histwrite_phy(o_u10m_srf(nsrf), zx_tmp_fi2d)
585          IF (vars_defined)       zx_tmp_fi2d(1 : klon) = v10m(1 : klon, nsrf)
586          CALL histwrite_phy(o_v10m_srf(nsrf), zx_tmp_fi2d)
587          IF (vars_defined)       zx_tmp_fi2d(1 : klon) = t2m(1 : klon, nsrf)
588          CALL histwrite_phy(o_t2m_srf(nsrf), zx_tmp_fi2d)
589          IF (vars_defined)       zx_tmp_fi2d(1 : klon) = fevap(1 : klon, nsrf)
590          CALL histwrite_phy(o_evap_srf(nsrf), zx_tmp_fi2d)
591          IF (vars_defined)        zx_tmp_fi2d(1 : klon) = fluxt( 1 : klon, 1, nsrf)
592          CALL histwrite_phy(o_sens_srf(nsrf), zx_tmp_fi2d)
593          IF (vars_defined)         zx_tmp_fi2d(1 : klon) = fluxlat( 1 : klon, nsrf)
594          CALL histwrite_phy(o_lat_srf(nsrf), zx_tmp_fi2d)
595          IF (vars_defined)         zx_tmp_fi2d(1 : klon) = fsollw( 1 : klon, nsrf)
596          CALL histwrite_phy(o_flw_srf(nsrf), zx_tmp_fi2d)
597          IF (vars_defined)         zx_tmp_fi2d(1 : klon) = fsolsw( 1 : klon, nsrf)
598          CALL histwrite_phy(o_fsw_srf(nsrf), zx_tmp_fi2d)
599          IF (vars_defined)         zx_tmp_fi2d(1 : klon) = wfbils( 1 : klon, nsrf)
600          CALL histwrite_phy(o_wbils_srf(nsrf), zx_tmp_fi2d)
601          IF (vars_defined)         zx_tmp_fi2d(1 : klon) = wfbilo( 1 : klon, nsrf)
602          CALL histwrite_phy(o_wbilo_srf(nsrf), zx_tmp_fi2d)
603
604          IF (iflag_pbl > 1) THEN
605             CALL histwrite_phy(o_tke_srf(nsrf),  pbl_tke(:,1:klev,nsrf))
606             CALL histwrite_phy(o_tke_max_srf(nsrf),  pbl_tke(:,1:klev,nsrf))
607          ENDIF
608!jyg<
609          IF (iflag_pbl > 1 .AND. iflag_wake>=1  .AND. iflag_pbl_split >=1) THEN
610             CALL histwrite_phy(o_dltpbltke_srf(nsrf), wake_delta_pbl_TKE(:,1:klev,nsrf))
611          ENDIF
612!>jyg
613
614       ENDDO
615       DO nsrf=1,nbsrf+1
616          CALL histwrite_phy(o_wstar(nsrf), wstar(1 : klon, nsrf))
617       ENDDO
618
619       CALL histwrite_phy(o_cdrm, cdragm)
620       CALL histwrite_phy(o_cdrh, cdragh)
621       CALL histwrite_phy(o_cldl, cldl)
622       CALL histwrite_phy(o_cldm, cldm)
623       CALL histwrite_phy(o_cldh, cldh)
624       CALL histwrite_phy(o_cldt, cldt)
625       CALL histwrite_phy(o_JrNt, JrNt)
626       CALL histwrite_phy(o_cldljn, cldl*JrNt)
627       CALL histwrite_phy(o_cldmjn, cldm*JrNt)
628       CALL histwrite_phy(o_cldhjn, cldh*JrNt)
629       CALL histwrite_phy(o_cldtjn, cldt*JrNt)
630       CALL histwrite_phy(o_cldq, cldq)
631       IF (vars_defined)       zx_tmp_fi2d(1:klon) = flwp(1:klon)
632       CALL histwrite_phy(o_lwp, zx_tmp_fi2d)
633       IF (vars_defined)       zx_tmp_fi2d(1:klon) = fiwp(1:klon)
634       CALL histwrite_phy(o_iwp, zx_tmp_fi2d)
635       CALL histwrite_phy(o_ue, ue)
636       CALL histwrite_phy(o_ve, ve)
637       CALL histwrite_phy(o_uq, uq)
638       CALL histwrite_phy(o_vq, vq)
639       IF(iflag_con.GE.3) THEN ! sb
640          CALL histwrite_phy(o_cape, cape)
641          CALL histwrite_phy(o_pbase, ema_pcb)
642          CALL histwrite_phy(o_ptop, ema_pct)
643          CALL histwrite_phy(o_fbase, ema_cbmf)
644          if (iflag_con /= 30) then
645             CALL histwrite_phy(o_plcl, plcl)
646             CALL histwrite_phy(o_plfc, plfc)
647             CALL histwrite_phy(o_wbeff, wbeff)
648          end if
649
650          CALL histwrite_phy(o_cape_max, cape)
651
652          CALL histwrite_phy(o_upwd, upwd)
653          CALL histwrite_phy(o_Ma, Ma)
654          CALL histwrite_phy(o_dnwd, dnwd)
655          CALL histwrite_phy(o_dnwd0, dnwd0)
656          IF (vars_defined)         zx_tmp_fi2d=float(itau_con)/float(itap)
657          CALL histwrite_phy(o_ftime_con, zx_tmp_fi2d)
658          IF (vars_defined) THEN
659             IF(iflag_thermals>=1)THEN
660                zx_tmp_fi3d=dnwd+dnwd0+upwd+fm_therm(:,1:klev)
661             ELSE
662                zx_tmp_fi3d=dnwd+dnwd0+upwd
663             ENDIF
664          ENDIF
665          CALL histwrite_phy(o_mc, zx_tmp_fi3d)
666       ENDIF !iflag_con .GE. 3
667       CALL histwrite_phy(o_prw, prw)
668       CALL histwrite_phy(o_s_pblh, s_pblh)
669       CALL histwrite_phy(o_s_pblt, s_pblt)
670       CALL histwrite_phy(o_s_lcl, s_lcl)
671       CALL histwrite_phy(o_s_therm, s_therm)
672       !IM : Les champs suivants (s_capCL, s_oliqCL, s_cteiCL, s_trmb1, s_trmb2, s_trmb3) ne sont pas definis dans HBTM.F
673       !       IF (o_s_capCL%flag(iff)<=lev_files(iff)) THEN
674       !     CALL histwrite_phy(nid_files(iff),clef_stations(iff),
675       !    $o_s_capCL%name,itau_w,s_capCL)
676       !       ENDIF
677       !       IF (o_s_oliqCL%flag(iff)<=lev_files(iff)) THEN
678       !     CALL histwrite_phy(nid_files(iff),clef_stations(iff),
679       !    $o_s_oliqCL%name,itau_w,s_oliqCL)
680       !       ENDIF
681       !       IF (o_s_cteiCL%flag(iff)<=lev_files(iff)) THEN
682       !     CALL histwrite_phy(nid_files(iff),clef_stations(iff),
683       !    $o_s_cteiCL%name,itau_w,s_cteiCL)
684       !       ENDIF
685       !       IF (o_s_trmb1%flag(iff)<=lev_files(iff)) THEN
686       !     CALL histwrite_phy(nid_files(iff),clef_stations(iff),
687       !    $o_s_trmb1%name,itau_w,s_trmb1)
688       !       ENDIF
689       !       IF (o_s_trmb2%flag(iff)<=lev_files(iff)) THEN
690       !     CALL histwrite_phy(nid_files(iff),clef_stations(iff),
691       !    $o_s_trmb2%name,itau_w,s_trmb2)
692       !       ENDIF
693       !       IF (o_s_trmb3%flag(iff)<=lev_files(iff)) THEN
694       !     CALL histwrite_phy(nid_files(iff),clef_stations(iff),
695       !    $o_s_trmb3%name,itau_w,s_trmb3)
696       !       ENDIF
697
698#ifdef CPP_IOIPSL
699#ifndef CPP_XIOS
700  IF (.NOT.ok_all_xml) THEN
701       ! ATTENTION, LES ANCIENS HISTWRITE ONT ETES CONSERVES EN ATTENDANT MIEUX:
702       ! Champs interpolles sur des niveaux de pression
703       missing_val=missing_val_nf90
704       DO iff=1, nfiles
705          ll=0
706          DO k=1, nlevSTD
707             bb2=clevSTD(k)
708             IF(bb2.EQ."850".OR.bb2.EQ."700".OR. &
709                  bb2.EQ."500".OR.bb2.EQ."200".OR. &
710                  bb2.EQ."100".OR. &
711                  bb2.EQ."50".OR.bb2.EQ."10") THEN
712
713                ! a refaire correctement !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
714                ll=ll+1
715                CALL histwrite_phy(o_uSTDlevs(ll),uwriteSTD(:,k,iff), iff)
716                CALL histwrite_phy(o_vSTDlevs(ll),vwriteSTD(:,k,iff), iff)
717                CALL histwrite_phy(o_wSTDlevs(ll),wwriteSTD(:,k,iff), iff)
718                CALL histwrite_phy(o_zSTDlevs(ll),phiwriteSTD(:,k,iff), iff)
719                CALL histwrite_phy(o_qSTDlevs(ll),qwriteSTD(:,k,iff), iff)
720                CALL histwrite_phy(o_tSTDlevs(ll),twriteSTD(:,k,iff), iff)
721
722             ENDIF !(bb2.EQ."850".OR.bb2.EQ."700".OR.
723          ENDDO
724       ENDDO
725  ENDIF
726#endif
727#endif
728#ifdef CPP_XIOS
729  IF(ok_all_xml) THEN
730!XIOS  CALL xios_get_field_attr("u850",default_value=missing_val)
731!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
732          ll=0
733          DO k=1, nlevSTD
734             bb2=clevSTD(k)
735             IF(bb2.EQ."850".OR.bb2.EQ."700".OR. &
736                bb2.EQ."500".OR.bb2.EQ."200".OR. &
737                bb2.EQ."100".OR. &
738                bb2.EQ."50".OR.bb2.EQ."10") THEN
739                ll=ll+1
740                CALL histwrite_phy(o_uSTDlevs(ll),ulevSTD(:,k))
741                CALL histwrite_phy(o_vSTDlevs(ll),vlevSTD(:,k))
742                CALL histwrite_phy(o_wSTDlevs(ll),wlevSTD(:,k))
743                CALL histwrite_phy(o_zSTDlevs(ll),philevSTD(:,k))
744                CALL histwrite_phy(o_qSTDlevs(ll),qlevSTD(:,k))
745                CALL histwrite_phy(o_tSTDlevs(ll),tlevSTD(:,k))
746             ENDIF !(bb2.EQ."850".OR.bb2.EQ."700".OR.
747          ENDDO
748  ENDIF
749#endif
750       IF (vars_defined) THEN
751          DO i=1, klon
752             IF (pctsrf(i,is_oce).GT.epsfra.OR. &
753                  pctsrf(i,is_sic).GT.epsfra) THEN
754                zx_tmp_fi2d(i) = (ftsol(i, is_oce) * pctsrf(i,is_oce)+ &
755                     ftsol(i, is_sic) * pctsrf(i,is_sic))/ &
756                     (pctsrf(i,is_oce)+pctsrf(i,is_sic))
757             ELSE
758                zx_tmp_fi2d(i) = 273.15
759             ENDIF
760          ENDDO
761       ENDIF
762       CALL histwrite_phy(o_t_oce_sic, zx_tmp_fi2d)
763
764       ! Couplage convection-couche limite
765       IF (iflag_con.GE.3) THEN
766          IF (iflag_coupl>=1) THEN
767             CALL histwrite_phy(o_ale_bl, ale_bl)
768             CALL histwrite_phy(o_alp_bl, alp_bl)
769          ENDIF !iflag_coupl>=1
770       ENDIF !(iflag_con.GE.3)
771       ! Wakes
772       IF (iflag_con.EQ.3) THEN
773          IF (iflag_wake>=1) THEN
774             CALL histwrite_phy(o_ale_wk, ale_wake)
775             CALL histwrite_phy(o_alp_wk, alp_wake)
776             IF (iflag_pbl_split>=1) THEN
777               IF (vars_defined) zx_tmp_fi3d(1:klon,1:klev)=dtvdf_x(1:klon,1:klev)/pdtphys
778               CALL histwrite_phy(o_dtvdf_x    ,zx_tmp_fi3d)
779               IF (vars_defined) zx_tmp_fi3d(1:klon,1:klev)=dtvdf_w(1:klon,1:klev)/pdtphys
780               CALL histwrite_phy(o_dtvdf_w    ,zx_tmp_fi3d)
781               IF (vars_defined) zx_tmp_fi3d(1:klon,1:klev)=dqvdf_x(1:klon,1:klev)/pdtphys
782               CALL histwrite_phy(o_dqvdf_x    ,zx_tmp_fi3d)
783               IF (vars_defined) zx_tmp_fi3d(1:klon,1:klev)=dqvdf_w(1:klon,1:klev)/pdtphys
784               CALL histwrite_phy(o_dqvdf_w    ,zx_tmp_fi3d)
785               CALL histwrite_phy(o_sens_x     ,sens_x     )
786               CALL histwrite_phy(o_sens_w     ,sens_w     )
787               CALL histwrite_phy(o_flat_x     ,zxfluxlat_x)
788               CALL histwrite_phy(o_flat_w     ,zxfluxlat_w)
789               CALL histwrite_phy(o_delta_tsurf,delta_tsurf)
790               CALL histwrite_phy(o_cdragh_x   ,cdragh_x   )
791               CALL histwrite_phy(o_cdragh_w   ,cdragh_w   )
792               CALL histwrite_phy(o_cdragm_x   ,cdragm_x   )
793               CALL histwrite_phy(o_cdragm_w   ,cdragm_w   )
794               CALL histwrite_phy(o_kh         ,kh         )
795               CALL histwrite_phy(o_kh_x       ,kh_x       )
796               CALL histwrite_phy(o_kh_w       ,kh_w       )
797             ENDIF   ! (iflag_pbl_split>=1)
798             CALL histwrite_phy(o_ale, ale)
799             CALL histwrite_phy(o_alp, alp)
800             CALL histwrite_phy(o_cin, cin)
801             CALL histwrite_phy(o_WAPE, wake_pe)
802             CALL histwrite_phy(o_wake_h, wake_h)
803             CALL histwrite_phy(o_wake_s, wake_s)
804             CALL histwrite_phy(o_wake_deltat, wake_deltat)
805             CALL histwrite_phy(o_wake_deltaq, wake_deltaq)
806             CALL histwrite_phy(o_wake_omg, wake_omg)
807             IF (vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_t_wake(1:klon,1:klev) &
808                  /pdtphys
809             CALL histwrite_phy(o_dtwak, zx_tmp_fi3d)
810             IF (vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_q_wake(1:klon,1:klev)/pdtphys
811             CALL histwrite_phy(o_dqwak, zx_tmp_fi3d)
812          ENDIF ! iflag_wake>=1
813          CALL histwrite_phy(o_ftd, ftd)
814          CALL histwrite_phy(o_fqd, fqd)
815       ENDIF !(iflag_con.EQ.3)
816       IF (iflag_con.EQ.3.OR.iflag_con.EQ.30) THEN
817          ! sortie RomP convection descente insaturee iflag_con=30
818          ! etendue a iflag_con=3 (jyg)
819          CALL histwrite_phy(o_Vprecip, Vprecip)
820          CALL histwrite_phy(o_wdtrainA, wdtrainA)
821          CALL histwrite_phy(o_wdtrainM, wdtrainM)
822       ENDIF !(iflag_con.EQ.3.or.iflag_con.EQ.30)
823!!! nrlmd le 10/04/2012
824       IF (iflag_trig_bl>=1) THEN
825          CALL histwrite_phy(o_n2, n2)
826          CALL histwrite_phy(o_s2, s2)
827          CALL histwrite_phy(o_proba_notrig, proba_notrig)
828          CALL histwrite_phy(o_random_notrig, random_notrig)
829          CALL histwrite_phy(o_ale_bl_stat, ale_bl_stat)
830          CALL histwrite_phy(o_ale_bl_trig, ale_bl_trig)
831       ENDIF  !(iflag_trig_bl>=1)
832       IF (iflag_clos_bl>=1) THEN
833          CALL histwrite_phy(o_alp_bl_det, alp_bl_det)
834          CALL histwrite_phy(o_alp_bl_fluct_m, alp_bl_fluct_m)
835          CALL histwrite_phy(o_alp_bl_fluct_tke,  &
836               alp_bl_fluct_tke)
837          CALL histwrite_phy(o_alp_bl_conv, alp_bl_conv)
838          CALL histwrite_phy(o_alp_bl_stat, alp_bl_stat)
839       ENDIF  !(iflag_clos_bl>=1)
840!!! fin nrlmd le 10/04/2012
841       ! Output of slab ocean variables
842       IF (type_ocean=='slab ') THEN
843          CALL histwrite_phy(o_slab_qflux, slab_wfbils)
844          CALL histwrite_phy(o_slab_bils, slab_bils)
845          IF (nslay.EQ.1) THEN
846              zx_tmp_fi2d(:)=tslab(:,1)
847              CALL histwrite_phy(o_tslab, zx_tmp_fi2d)
848          ELSE
849              CALL histwrite_phy(o_tslab, tslab)
850          END IF
851          IF (version_ocean=='sicINT') THEN
852              CALL histwrite_phy(o_slab_bilg, slab_bilg)
853              CALL histwrite_phy(o_slab_tice, tice)
854              CALL histwrite_phy(o_slab_sic, seaice)
855          END IF
856       ENDIF !type_ocean == force/slab
857       CALL histwrite_phy(o_weakinv, weak_inversion)
858       CALL histwrite_phy(o_dthmin, dthmin)
859       CALL histwrite_phy(o_cldtau, cldtau)
860       CALL histwrite_phy(o_cldemi, cldemi)
861       CALL histwrite_phy(o_pr_con_l, pmflxr(:,1:klev))
862       CALL histwrite_phy(o_pr_con_i, pmflxs(:,1:klev))
863       CALL histwrite_phy(o_pr_lsc_l, prfl(:,1:klev))
864       CALL histwrite_phy(o_pr_lsc_i, psfl(:,1:klev))
865       CALL histwrite_phy(o_re, re)
866       CALL histwrite_phy(o_fl, fl)
867       IF (vars_defined) THEN
868          DO i=1, klon
869             zx_tmp_fi2d(i)=MIN(100.,rh2m(i)*100.)
870          ENDDO
871       ENDIF
872       CALL histwrite_phy(o_rh2m, zx_tmp_fi2d)
873
874       IF (vars_defined) THEN
875          DO i=1, klon
876             zx_tmp_fi2d(i)=MIN(100.,rh2m(i)*100.)
877          ENDDO
878       ENDIF
879       CALL histwrite_phy(o_rh2m_min, zx_tmp_fi2d)
880
881       IF (vars_defined) THEN
882          DO i=1, klon
883             zx_tmp_fi2d(i)=MIN(100.,rh2m(i)*100.)
884          ENDDO
885       ENDIF
886       CALL histwrite_phy(o_rh2m_max, zx_tmp_fi2d)
887
888       CALL histwrite_phy(o_qsat2m, qsat2m)
889       CALL histwrite_phy(o_tpot, tpot)
890       CALL histwrite_phy(o_tpote, tpote)
891       IF (vars_defined) zx_tmp_fi2d(1 : klon) = fsolsw( 1 : klon, is_ter)
892       CALL histwrite_phy(o_SWnetOR,  zx_tmp_fi2d)
893       IF (vars_defined) zx_tmp_fi2d(1:klon) = solsw(1:klon)/(1.-albsol1(1:klon))
894       CALL histwrite_phy(o_SWdownOR,  zx_tmp_fi2d)
895       CALL histwrite_phy(o_LWdownOR, sollwdown)
896       CALL histwrite_phy(o_snowl, snow_lsc)
897       CALL histwrite_phy(o_solldown, sollwdown)
898       CALL histwrite_phy(o_dtsvdfo, d_ts(:,is_oce))
899       CALL histwrite_phy(o_dtsvdft, d_ts(:,is_ter))
900       CALL histwrite_phy(o_dtsvdfg,  d_ts(:,is_lic))
901       CALL histwrite_phy(o_dtsvdfi, d_ts(:,is_sic))
902       CALL histwrite_phy(o_z0m, z0m(:,nbsrf+1))
903       CALL histwrite_phy(o_z0h, z0h(:,nbsrf+1))
904       ! OD550 per species
905!--OLIVIER
906!This is warranted by treating INCA aerosols as offline aerosols
907!       IF (new_aod .and. (.not. aerosol_couple)) THEN
908       IF (new_aod) THEN
909          IF (flag_aerosol.GT.0) THEN
910             CALL histwrite_phy(o_od550aer, od550aer)
911             CALL histwrite_phy(o_od865aer, od865aer)
912             CALL histwrite_phy(o_absvisaer, absvisaer)
913             CALL histwrite_phy(o_od550lt1aer, od550lt1aer)
914             CALL histwrite_phy(o_sconcso4, sconcso4)
915             CALL histwrite_phy(o_sconcno3, sconcno3)
916             CALL histwrite_phy(o_sconcoa, sconcoa)
917             CALL histwrite_phy(o_sconcbc, sconcbc)
918             CALL histwrite_phy(o_sconcss, sconcss)
919             CALL histwrite_phy(o_sconcdust, sconcdust)
920             CALL histwrite_phy(o_concso4, concso4)
921             CALL histwrite_phy(o_concno3, concno3)
922             CALL histwrite_phy(o_concoa, concoa)
923             CALL histwrite_phy(o_concbc, concbc)
924             CALL histwrite_phy(o_concss, concss)
925             CALL histwrite_phy(o_concdust, concdust)
926             CALL histwrite_phy(o_loadso4, loadso4)
927             CALL histwrite_phy(o_loadoa, loadoa)
928             CALL histwrite_phy(o_loadbc, loadbc)
929             CALL histwrite_phy(o_loadss, loadss)
930             CALL histwrite_phy(o_loaddust, loaddust)
931             !--STRAT AER
932          ENDIF
933          IF (flag_aerosol.GT.0.OR.flag_aerosol_strat) THEN
934!             DO naero = 1, naero_spc
935!--correction mini bug OB
936             DO naero = 1, naero_tot
937                CALL histwrite_phy(o_tausumaero(naero), &
938                     tausum_aero(:,2,naero) )
939             END DO
940          ENDIF
941          IF (flag_aerosol_strat) THEN
942             CALL histwrite_phy(o_tausumaero_lw, &
943                  tausum_aero(:,6,id_STRAT_phy) )
944          ENDIF
945       ENDIF
946       IF (ok_ade) THEN
947          CALL histwrite_phy(o_topswad, topswad_aero*swradcorr)
948          CALL histwrite_phy(o_topswad0, topswad0_aero*swradcorr)
949          CALL histwrite_phy(o_solswad, solswad_aero*swradcorr)
950          CALL histwrite_phy(o_solswad0, solswad0_aero*swradcorr)
951          CALL histwrite_phy(o_toplwad, toplwad_aero)
952          CALL histwrite_phy(o_toplwad0, toplwad0_aero)
953          CALL histwrite_phy(o_sollwad, sollwad_aero)
954          CALL histwrite_phy(o_sollwad0, sollwad0_aero)
955          !====MS forcing diagnostics
956          if (new_aod) then
957             zx_tmp_fi2d(:)=topsw_aero(:,1)*swradcorr(:)
958             CALL histwrite_phy(o_swtoaas_nat,zx_tmp_fi2d)
959             zx_tmp_fi2d(:)=solsw_aero(:,1)*swradcorr(:)
960             CALL histwrite_phy(o_swsrfas_nat,zx_tmp_fi2d)
961             zx_tmp_fi2d(:)=topsw0_aero(:,1)*swradcorr(:)
962             CALL histwrite_phy(o_swtoacs_nat,zx_tmp_fi2d)
963             zx_tmp_fi2d(:)=solsw0_aero(:,1)*swradcorr(:)
964             CALL histwrite_phy(o_swsrfcs_nat,zx_tmp_fi2d)
965             !ant
966             zx_tmp_fi2d(:)=topsw_aero(:,2)*swradcorr(:)
967             CALL histwrite_phy(o_swtoaas_ant,zx_tmp_fi2d)
968             zx_tmp_fi2d(:)=solsw_aero(:,2)*swradcorr(:)
969             CALL histwrite_phy(o_swsrfas_ant,zx_tmp_fi2d)
970             zx_tmp_fi2d(:)=topsw0_aero(:,2)*swradcorr(:)
971             CALL histwrite_phy(o_swtoacs_ant,zx_tmp_fi2d)
972             zx_tmp_fi2d(:)=solsw0_aero(:,2)*swradcorr(:)
973             CALL histwrite_phy(o_swsrfcs_ant,zx_tmp_fi2d)
974             !cf
975             if (.not. aerosol_couple) then
976                zx_tmp_fi2d(:)=topswcf_aero(:,1)*swradcorr(:)
977                CALL histwrite_phy(o_swtoacf_nat,zx_tmp_fi2d)
978                zx_tmp_fi2d(:)=solswcf_aero(:,1)*swradcorr(:)
979                CALL histwrite_phy(o_swsrfcf_nat,zx_tmp_fi2d)
980                zx_tmp_fi2d(:)=topswcf_aero(:,2)*swradcorr(:)
981                CALL histwrite_phy(o_swtoacf_ant,zx_tmp_fi2d)
982                zx_tmp_fi2d(:)=solswcf_aero(:,2)*swradcorr(:)
983                CALL histwrite_phy(o_swsrfcf_ant,zx_tmp_fi2d)
984                zx_tmp_fi2d(:)=topswcf_aero(:,3)*swradcorr(:)
985                CALL histwrite_phy(o_swtoacf_zero,zx_tmp_fi2d)
986                zx_tmp_fi2d(:)=solswcf_aero(:,3)*swradcorr(:)
987                CALL histwrite_phy(o_swsrfcf_zero,zx_tmp_fi2d)
988             endif
989          endif ! new_aod
990          !====MS forcing diagnostics
991       ENDIF
992       IF (ok_aie) THEN
993          CALL histwrite_phy(o_topswai, topswai_aero*swradcorr)
994          CALL histwrite_phy(o_solswai, solswai_aero*swradcorr)
995       ENDIF
996       IF (flag_aerosol.GT.0.AND.ok_cdnc) THEN
997          CALL histwrite_phy(o_scdnc, scdnc)
998          CALL histwrite_phy(o_cldncl, cldncl)
999          CALL histwrite_phy(o_reffclws, reffclws)
1000          CALL histwrite_phy(o_reffclwc, reffclwc)
1001          CALL histwrite_phy(o_cldnvi, cldnvi)
1002          CALL histwrite_phy(o_lcc, lcc)
1003          CALL histwrite_phy(o_lcc3d, lcc3d)
1004          CALL histwrite_phy(o_lcc3dcon, lcc3dcon)
1005          CALL histwrite_phy(o_lcc3dstra, lcc3dstra)
1006          CALL histwrite_phy(o_reffclwtop, reffclwtop)
1007       ENDIF
1008       ! Champs 3D:
1009       IF (ok_ade .OR. ok_aie) then
1010          CALL histwrite_phy(o_ec550aer, ec550aer)
1011       ENDIF
1012       CALL histwrite_phy(o_lwcon, flwc)
1013       CALL histwrite_phy(o_iwcon, fiwc)
1014       CALL histwrite_phy(o_temp, t_seri)
1015       CALL histwrite_phy(o_theta, theta)
1016       CALL histwrite_phy(o_ovapinit, qx(:,:,ivap))
1017       CALL histwrite_phy(o_ovap, q_seri)
1018       CALL histwrite_phy(o_oliq, ql_seri)
1019       CALL histwrite_phy(o_geop, zphi)
1020       CALL histwrite_phy(o_vitu, u_seri)
1021       CALL histwrite_phy(o_vitv, v_seri)
1022       CALL histwrite_phy(o_vitw, omega)
1023       CALL histwrite_phy(o_pres, pplay)
1024       CALL histwrite_phy(o_paprs, paprs(:,1:klev))
1025       CALL histwrite_phy(o_zfull,zphi/RG)
1026
1027       IF (vars_defined)  THEN
1028        zx_tmp_fi3d(:,1)= pphis(:)/RG
1029        DO k = 2, klev
1030         DO i = 1, klon
1031            zx_tmp_fi3d(i,k) = zphi(i,k)/RG + &
1032                          (zphi(i,k)-zphi(i,k-1))/RG * &
1033                          (paprs(i,k)-pplay(i,k))/(pplay(i,k)-pplay(i,k-1))
1034         ENDDO
1035        ENDDO
1036       ENDIF
1037       CALL histwrite_phy(o_zhalf, zx_tmp_fi3d)
1038       CALL histwrite_phy(o_rneb, cldfra)
1039       CALL histwrite_phy(o_rnebcon, rnebcon)
1040       CALL histwrite_phy(o_rnebls, rneb)
1041       IF (vars_defined)  THEN
1042          DO k=1, klev
1043             DO i=1, klon
1044                zx_tmp_fi3d(i,k)=cldfra(i,k)*JrNt(i)
1045             ENDDO
1046          ENDDO
1047       ENDIF
1048       CALL histwrite_phy(o_rnebjn, zx_tmp_fi3d)
1049       CALL histwrite_phy(o_rhum, zx_rh)
1050       CALL histwrite_phy(o_ozone, &
1051            wo(:, :, 1) * dobson_u * 1e3 / zmasse / rmo3 * rmd)
1052
1053       IF (read_climoz == 2) THEN
1054          CALL histwrite_phy(o_ozone_light, &
1055               wo(:, :, 2) * dobson_u * 1e3 / zmasse / rmo3 * rmd)
1056       ENDIF
1057
1058       CALL histwrite_phy(o_dtphy, d_t)
1059       CALL histwrite_phy(o_dqphy,  d_qx(:,:,ivap))
1060       DO nsrf=1, nbsrf
1061          IF (vars_defined) zx_tmp_fi2d(1 : klon) = falb1( 1 : klon, nsrf)
1062          CALL histwrite_phy(o_albe_srf(nsrf), zx_tmp_fi2d)
1063          IF (vars_defined) zx_tmp_fi2d(1 : klon) = z0m( 1 : klon, nsrf)
1064          CALL histwrite_phy(o_z0m_srf(nsrf), zx_tmp_fi2d)
1065          IF (vars_defined) zx_tmp_fi2d(1 : klon) = z0h( 1 : klon, nsrf)
1066          CALL histwrite_phy(o_z0h_srf(nsrf), zx_tmp_fi2d)
1067          IF (vars_defined) zx_tmp_fi2d(1 : klon) = agesno( 1 : klon, nsrf)
1068          CALL histwrite_phy(o_ages_srf(nsrf), zx_tmp_fi2d)
1069          IF (vars_defined) zx_tmp_fi2d(1 : klon) = snow( 1 : klon, nsrf)
1070          CALL histwrite_phy(o_snow_srf(nsrf), zx_tmp_fi2d)
1071       ENDDO !nsrf=1, nbsrf
1072       CALL histwrite_phy(o_alb1, albsol1)
1073       CALL histwrite_phy(o_alb2, albsol2)
1074       !FH Sorties pour la couche limite
1075       if (iflag_pbl>1) then
1076          zx_tmp_fi3d=0.
1077          IF (vars_defined) THEN
1078             do nsrf=1,nbsrf
1079                do k=1,klev
1080                   zx_tmp_fi3d(:,k)=zx_tmp_fi3d(:,k) &
1081                        +pctsrf(:,nsrf)*pbl_tke(:,k,nsrf)
1082                enddo
1083             enddo
1084          ENDIF
1085          CALL histwrite_phy(o_tke, zx_tmp_fi3d)
1086
1087          CALL histwrite_phy(o_tke_max, zx_tmp_fi3d)
1088       ENDIF
1089
1090       CALL histwrite_phy(o_kz, coefh(:,:,is_ave))
1091
1092       CALL histwrite_phy(o_kz_max, coefh(:,:,is_ave))
1093
1094       CALL histwrite_phy(o_clwcon, clwcon0)
1095       CALL histwrite_phy(o_dtdyn, d_t_dyn)
1096       CALL histwrite_phy(o_dqdyn, d_q_dyn)
1097       CALL histwrite_phy(o_dudyn, d_u_dyn)
1098       CALL histwrite_phy(o_dvdyn, d_v_dyn)
1099
1100       IF (vars_defined) THEN
1101          zx_tmp_fi3d(1:klon,1:klev)=d_t_con(1:klon,1:klev)/pdtphys
1102       ENDIF
1103       CALL histwrite_phy(o_dtcon, zx_tmp_fi3d)
1104       if(iflag_thermals.eq.0)then
1105          IF (vars_defined) THEN
1106             zx_tmp_fi3d(1:klon,1:klev)=d_t_con(1:klon,1:klev)/pdtphys + &
1107                  d_t_ajsb(1:klon,1:klev)/pdtphys
1108          ENDIF
1109          CALL histwrite_phy(o_tntc, zx_tmp_fi3d)
1110       else if(iflag_thermals.ge.1.and.iflag_wake.EQ.1)then
1111          IF (vars_defined) THEN
1112             zx_tmp_fi3d(1:klon,1:klev)=d_t_con(1:klon,1:klev)/pdtphys + &
1113                  d_t_ajs(1:klon,1:klev)/pdtphys + &
1114                  d_t_wake(1:klon,1:klev)/pdtphys
1115          ENDIF
1116          CALL histwrite_phy(o_tntc, zx_tmp_fi3d)
1117       endif
1118       IF (vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_u_con(1:klon,1:klev)/pdtphys
1119       CALL histwrite_phy(o_ducon, zx_tmp_fi3d)
1120       IF (vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_v_con(1:klon,1:klev)/pdtphys
1121       CALL histwrite_phy(o_dvcon, zx_tmp_fi3d)
1122       IF (vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_q_con(1:klon,1:klev)/pdtphys
1123       CALL histwrite_phy(o_dqcon, zx_tmp_fi3d)
1124
1125       IF(iflag_thermals.EQ.0) THEN
1126          IF (vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_q_con(1:klon,1:klev)/pdtphys
1127          CALL histwrite_phy(o_tnhusc, zx_tmp_fi3d)
1128       ELSE IF(iflag_thermals.GE.1.AND.iflag_wake.EQ.1) THEN
1129          IF (vars_defined) THEN
1130             zx_tmp_fi3d(1:klon,1:klev)=d_q_con(1:klon,1:klev)/pdtphys + &
1131                  d_q_ajs(1:klon,1:klev)/pdtphys + &
1132                  d_q_wake(1:klon,1:klev)/pdtphys
1133          ENDIF
1134          CALL histwrite_phy(o_tnhusc, zx_tmp_fi3d)
1135       ENDIF
1136
1137       IF (vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_t_lsc(1:klon,1:klev)/pdtphys
1138       CALL histwrite_phy(o_dtlsc, zx_tmp_fi3d)
1139       IF (vars_defined) zx_tmp_fi3d(1:klon, 1:klev)=(d_t_lsc(1:klon,1:klev)+ &
1140            d_t_eva(1:klon,1:klev))/pdtphys
1141       CALL histwrite_phy(o_dtlschr, zx_tmp_fi3d)
1142       IF (vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_q_lsc(1:klon,1:klev)/pdtphys
1143       CALL histwrite_phy(o_dqlsc, zx_tmp_fi3d)
1144       IF (vars_defined) zx_tmp_fi3d(1:klon,1:klev)=beta_prec(1:klon,1:klev)
1145       CALL histwrite_phy(o_beta_prec, zx_tmp_fi3d)
1146!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1147       ! Sorties specifiques a la separation thermiques/non thermiques
1148       if (iflag_thermals>=1) then
1149          IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_t_lscth(1:klon,1:klev)/pdtphys
1150          CALL histwrite_phy(o_dtlscth, zx_tmp_fi3d)
1151          IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_t_lscst(1:klon,1:klev)/pdtphys
1152          CALL histwrite_phy(o_dtlscst, zx_tmp_fi3d)
1153          IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_q_lscth(1:klon,1:klev)/pdtphys
1154          CALL histwrite_phy(o_dqlscth, zx_tmp_fi3d)
1155          IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_q_lscst(1:klon,1:klev)/pdtphys
1156          CALL histwrite_phy(o_dqlscst, zx_tmp_fi3d)
1157          CALL histwrite_phy(o_plulth, plul_th)
1158          CALL histwrite_phy(o_plulst, plul_st)
1159          IF (vars_defined) THEN
1160             do k=1,klev
1161                do i=1,klon
1162                   if (ptconvth(i,k)) then
1163                      zx_tmp_fi3d(i,k)=1.
1164                   else
1165                      zx_tmp_fi3d(i,k)=0.
1166                   endif
1167                enddo
1168             enddo
1169          ENDIF
1170          CALL histwrite_phy(o_ptconvth, zx_tmp_fi3d)
1171          IF (vars_defined) THEN
1172             do i=1,klon
1173                zx_tmp_fi2d(1:klon)=lmax_th(:)
1174             enddo
1175          ENDIF
1176          CALL histwrite_phy(o_lmaxth, zx_tmp_fi2d)
1177       endif ! iflag_thermals>=1
1178!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1179       IF (vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_t_vdf(1:klon,1:klev)/pdtphys
1180       CALL histwrite_phy(o_dtvdf, zx_tmp_fi3d)
1181       IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_t_diss(1:klon,1:klev)/pdtphys
1182       CALL histwrite_phy(o_dtdis, zx_tmp_fi3d)
1183       IF (vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_q_vdf(1:klon,1:klev)/pdtphys
1184       CALL histwrite_phy(o_dqvdf, zx_tmp_fi3d)
1185       IF (vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_t_eva(1:klon,1:klev)/pdtphys
1186       CALL histwrite_phy(o_dteva, zx_tmp_fi3d)
1187       IF (vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_q_eva(1:klon,1:klev)/pdtphys
1188       CALL histwrite_phy(o_dqeva, zx_tmp_fi3d)
1189       zpt_conv = 0.
1190       WHERE (ptconv) zpt_conv = 1.
1191       CALL histwrite_phy(o_ptconv, zpt_conv)
1192       CALL histwrite_phy(o_ratqs, ratqs)
1193       IF (vars_defined) THEN
1194          zx_tmp_fi3d(1:klon,1:klev)=d_t_ajs(1:klon,1:klev)/pdtphys - &
1195               d_t_ajsb(1:klon,1:klev)/pdtphys
1196       ENDIF
1197       CALL histwrite_phy(o_dtthe, zx_tmp_fi3d)
1198       IF (vars_defined) THEN
1199          zx_tmp_fi3d(1:klon,1:klev)=d_u_ajs(1:klon,1:klev)/pdtphys
1200       ENDIF
1201       CALL histwrite_phy(o_duthe, zx_tmp_fi3d)
1202       IF (vars_defined) THEN
1203          zx_tmp_fi3d(1:klon,1:klev)=d_v_ajs(1:klon,1:klev)/pdtphys
1204       ENDIF
1205       CALL histwrite_phy(o_dvthe, zx_tmp_fi3d)
1206
1207       IF (iflag_thermals>=1) THEN
1208          ! Pour l instant 0 a y reflichir pour les thermiques
1209          zx_tmp_fi2d=0.
1210          CALL histwrite_phy(o_ftime_th, zx_tmp_fi2d)
1211          CALL histwrite_phy(o_f_th, fm_therm)
1212          CALL histwrite_phy(o_e_th, entr_therm)
1213          CALL histwrite_phy(o_w_th, zw2)
1214          CALL histwrite_phy(o_q_th, zqasc)
1215          CALL histwrite_phy(o_a_th, fraca)
1216          CALL histwrite_phy(o_d_th, detr_therm)
1217          CALL histwrite_phy(o_f0_th, f0)
1218          CALL histwrite_phy(o_zmax_th, zmax_th)
1219          IF (vars_defined) THEN
1220             zx_tmp_fi3d(1:klon,1:klev)=d_q_ajs(1:klon,1:klev)/pdtphys - &
1221                  d_q_ajsb(1:klon,1:klev)/pdtphys
1222          ENDIF
1223          CALL histwrite_phy(o_dqthe, zx_tmp_fi3d)
1224       ENDIF !iflag_thermals
1225       IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_t_ajsb(1:klon,1:klev)/pdtphys
1226       CALL histwrite_phy(o_dtajs, zx_tmp_fi3d)
1227       IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_q_ajsb(1:klon,1:klev)/pdtphys
1228       CALL histwrite_phy(o_dqajs, zx_tmp_fi3d)
1229       IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_t_swr(1:klon,1:klev)/pdtphys
1230       CALL histwrite_phy(o_dtswr, zx_tmp_fi3d)
1231       IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_t_sw0(1:klon,1:klev)/pdtphys
1232       CALL histwrite_phy(o_dtsw0, zx_tmp_fi3d)
1233       IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_t_lwr(1:klon,1:klev)/pdtphys
1234       CALL histwrite_phy(o_dtlwr, zx_tmp_fi3d)
1235       IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_t_lw0(1:klon,1:klev)/pdtphys
1236       CALL histwrite_phy(o_dtlw0, zx_tmp_fi3d)
1237       IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_t_ec(1:klon,1:klev)/pdtphys
1238       CALL histwrite_phy(o_dtec, zx_tmp_fi3d)
1239       IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_u_vdf(1:klon,1:klev)/pdtphys
1240       CALL histwrite_phy(o_duvdf, zx_tmp_fi3d)
1241       IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_v_vdf(1:klon,1:klev)/pdtphys
1242       CALL histwrite_phy(o_dvvdf, zx_tmp_fi3d)
1243       IF (ok_orodr) THEN
1244          IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_u_oro(1:klon,1:klev)/pdtphys
1245          CALL histwrite_phy(o_duoro, zx_tmp_fi3d)
1246          IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_v_oro(1:klon,1:klev)/pdtphys
1247          CALL histwrite_phy(o_dvoro, zx_tmp_fi3d)
1248          IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_t_oro(1:klon,1:klev)/pdtphys
1249          CALL histwrite_phy(o_dtoro, zx_tmp_fi3d)
1250       ENDIF
1251       IF (ok_orolf) THEN
1252          IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_u_lif(1:klon,1:klev)/pdtphys
1253          CALL histwrite_phy(o_dulif, zx_tmp_fi3d)
1254
1255          IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_v_lif(1:klon,1:klev)/pdtphys
1256          CALL histwrite_phy(o_dvlif, zx_tmp_fi3d)
1257
1258          IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_t_lif(1:klon,1:klev)/pdtphys
1259          CALL histwrite_phy(o_dtlif, zx_tmp_fi3d)
1260       ENDIF
1261
1262       IF (ok_hines) THEN
1263          CALL histwrite_phy(o_du_gwd_hines, du_gwd_hines/pdtphys)
1264          CALL histwrite_phy(o_dv_gwd_hines, dv_gwd_hines/pdtphys)
1265          IF(vars_defined) zx_tmp_fi3d(1:klon,1:klev)=d_t_hin(1:klon,1:klev)/pdtphys
1266          CALL histwrite_phy(o_dthin, zx_tmp_fi3d)
1267          CALL histwrite_phy(o_ustr_gwd_hines, zustr_gwd_hines)
1268          CALL histwrite_phy(o_vstr_gwd_hines, zvstr_gwd_hines)
1269       end IF
1270
1271       if (.not. ok_hines .and. ok_gwd_rando) then
1272          CALL histwrite_phy(o_du_gwd_front, du_gwd_front / pdtphys)
1273          CALL histwrite_phy(o_dv_gwd_front, dv_gwd_front / pdtphys)
1274          CALL histwrite_phy(o_ustr_gwd_front, zustr_gwd_front)
1275          CALL histwrite_phy(o_vstr_gwd_front, zvstr_gwd_front)
1276       ENDIF
1277
1278       IF (ok_gwd_rando) then
1279          CALL histwrite_phy(o_du_gwd_rando, du_gwd_rando / pdtphys)
1280          CALL histwrite_phy(o_dv_gwd_rando, dv_gwd_rando / pdtphys)
1281          CALL histwrite_phy(o_ustr_gwd_rando, zustr_gwd_rando)
1282          CALL histwrite_phy(o_vstr_gwd_rando, zvstr_gwd_rando)
1283          CALL histwrite_phy(o_east_gwstress, east_gwstress )
1284          CALL histwrite_phy(o_west_gwstress, west_gwstress )
1285       end IF
1286
1287       IF (ok_qch4) then
1288          CALL histwrite_phy(o_dqch4, d_q_ch4 / pdtphys)
1289       ENDIF
1290
1291       DO k=1, klevp1
1292         zx_tmp_fi3d1(:,k)=swup(:,k)*swradcorr(:)
1293       ENDDO
1294       CALL histwrite_phy(o_rsu, zx_tmp_fi3d1)
1295       DO k=1, klevp1
1296         zx_tmp_fi3d1(:,k)=swdn(:,k)*swradcorr(:)
1297       ENDDO
1298       CALL histwrite_phy(o_rsd, zx_tmp_fi3d1)
1299       DO k=1, klevp1
1300         zx_tmp_fi3d1(:,k)=swup0(:,k)*swradcorr(:)
1301       ENDDO
1302       CALL histwrite_phy(o_rsucs, zx_tmp_fi3d1)
1303       DO k=1, klevp1
1304         zx_tmp_fi3d1(:,k)=swdn0(:,k)*swradcorr(:)
1305       ENDDO
1306       CALL histwrite_phy(o_rsdcs, zx_tmp_fi3d1)
1307
1308       CALL histwrite_phy(o_rlu, lwup)
1309       CALL histwrite_phy(o_rld, lwdn)
1310       CALL histwrite_phy(o_rlucs, lwup0)
1311       CALL histwrite_phy(o_rldcs, lwdn0)
1312
1313       IF(vars_defined) THEN
1314          zx_tmp_fi3d(1:klon,1:klev)=d_t(1:klon,1:klev)+ &
1315               d_t_dyn(1:klon,1:klev)
1316       ENDIF
1317       CALL histwrite_phy(o_tnt, zx_tmp_fi3d)
1318
1319       IF(vars_defined) THEN
1320          zx_tmp_fi3d(1:klon,1:klev)=d_t_swr(1:klon,1:klev)/pdtphys + &
1321               d_t_lwr(1:klon,1:klev)/pdtphys
1322       ENDIF
1323       CALL histwrite_phy(o_tntr, zx_tmp_fi3d)
1324       IF(vars_defined) THEN
1325          zx_tmp_fi3d(1:klon,1:klev)= (d_t_lsc(1:klon,1:klev)+ &
1326               d_t_eva(1:klon,1:klev)+ &
1327               d_t_vdf(1:klon,1:klev))/pdtphys
1328       ENDIF
1329       CALL histwrite_phy(o_tntscpbl, zx_tmp_fi3d)
1330       IF(vars_defined) THEN
1331          zx_tmp_fi3d(1:klon,1:klev)=d_qx(1:klon,1:klev,ivap)+ &
1332               d_q_dyn(1:klon,1:klev)
1333       ENDIF
1334       CALL histwrite_phy(o_tnhus, zx_tmp_fi3d)
1335       IF(vars_defined) THEN
1336          zx_tmp_fi3d(1:klon,1:klev)=d_q_lsc(1:klon,1:klev)/pdtphys+ &
1337               d_q_eva(1:klon,1:klev)/pdtphys
1338       ENDIF
1339       CALL histwrite_phy(o_tnhusscpbl, zx_tmp_fi3d)
1340       CALL histwrite_phy(o_evu, coefm(:,:,is_ave))
1341       IF(vars_defined) THEN
1342          zx_tmp_fi3d(1:klon,1:klev)=q_seri(1:klon,1:klev)+ &
1343               ql_seri(1:klon,1:klev)
1344       ENDIF
1345       CALL histwrite_phy(o_h2o, zx_tmp_fi3d)
1346       if (iflag_con >= 3) then
1347          IF(vars_defined) THEN
1348             zx_tmp_fi3d(1:klon,1:klev)=-1 * (dnwd(1:klon,1:klev)+ &
1349                  dnwd0(1:klon,1:klev))
1350          ENDIF
1351          CALL histwrite_phy(o_mcd, zx_tmp_fi3d)
1352          IF(vars_defined) THEN
1353             zx_tmp_fi3d(1:klon,1:klev)=upwd(1:klon,1:klev) + &
1354                  dnwd(1:klon,1:klev)+ dnwd0(1:klon,1:klev)
1355          ENDIF
1356          CALL histwrite_phy(o_dmc, zx_tmp_fi3d)
1357       else if (iflag_con == 2) then
1358          CALL histwrite_phy(o_mcd,  pmfd)
1359          CALL histwrite_phy(o_dmc,  pmfu + pmfd)
1360       end if
1361       CALL histwrite_phy(o_ref_liq, ref_liq)
1362       CALL histwrite_phy(o_ref_ice, ref_ice)
1363       if (RCO2_per.NE.RCO2_act.OR.RCH4_per.NE.RCH4_act.OR. &
1364            RN2O_per.NE.RN2O_act.OR.RCFC11_per.NE.RCFC11_act.OR. &
1365            RCFC12_per.NE.RCFC12_act) THEN
1366          IF(vars_defined) zx_tmp_fi2d(:) = swupp(:,klevp1)*swradcorr(:)
1367          CALL histwrite_phy(o_rsut4co2, zx_tmp_fi2d)
1368          IF(vars_defined) zx_tmp_fi2d(:) = lwupp(:,klevp1)
1369          CALL histwrite_phy(o_rlut4co2, zx_tmp_fi2d)
1370          IF(vars_defined) zx_tmp_fi2d(:) = swup0p(:,klevp1)*swradcorr(:)
1371          CALL histwrite_phy(o_rsutcs4co2, zx_tmp_fi2d)
1372          IF(vars_defined) zx_tmp_fi2d(:) = lwup0p(:,klevp1)
1373          CALL histwrite_phy(o_rlutcs4co2, zx_tmp_fi2d)
1374          DO k=1, klevp1
1375            zx_tmp_fi3d1(:,k)=swupp(:,k)*swradcorr(:)
1376          ENDDO
1377          CALL histwrite_phy(o_rsu4co2, zx_tmp_fi3d1)
1378          DO k=1, klevp1
1379            zx_tmp_fi3d1(:,k)=swup0p(:,k)*swradcorr(:)
1380          ENDDO
1381          CALL histwrite_phy(o_rsucs4co2, zx_tmp_fi3d1)
1382          DO k=1, klevp1
1383            zx_tmp_fi3d1(:,k)=swdnp(:,k)*swradcorr(:)
1384          ENDDO
1385          CALL histwrite_phy(o_rsd4co2, zx_tmp_fi3d1)
1386          DO k=1, klevp1
1387            zx_tmp_fi3d1(:,k)=swdn0p(:,k)*swradcorr(:)
1388          ENDDO
1389          CALL histwrite_phy(o_rsdcs4co2, zx_tmp_fi3d1)
1390          CALL histwrite_phy(o_rlu4co2, lwupp)
1391          CALL histwrite_phy(o_rlucs4co2, lwup0p)
1392          CALL histwrite_phy(o_rld4co2, lwdnp)
1393          CALL histwrite_phy(o_rldcs4co2, lwdn0p)
1394       ENDIF
1395!!!!!!!!!!!! Sorties niveaux de pression NMC !!!!!!!!!!!!!!!!!!!!
1396#ifdef CPP_IOIPSL
1397#ifndef CPP_XIOS
1398  IF (.NOT.ok_all_xml) THEN
1399       ! ATTENTION, LES ANCIENS HISTWRITE ONT ETES CONSERVES EN ATTENDANT MIEUX:
1400       ! Champs interpolles sur des niveaux de pression
1401       missing_val=missing_val_nf90
1402       DO iff=7, nfiles
1403
1404          CALL histwrite_phy(o_tnondef,tnondef(:,:,iff-6),iff)
1405          CALL histwrite_phy(o_ta,twriteSTD(:,:,iff-6),iff)
1406          CALL histwrite_phy(o_zg,phiwriteSTD(:,:,iff-6),iff)
1407          CALL histwrite_phy(o_hus,qwriteSTD(:,:,iff-6),iff)
1408          CALL histwrite_phy(o_hur,rhwriteSTD(:,:,iff-6),iff)
1409          CALL histwrite_phy(o_ua,uwriteSTD(:,:,iff-6),iff)
1410          CALL histwrite_phy(o_va,vwriteSTD(:,:,iff-6),iff)
1411          CALL histwrite_phy(o_wap,wwriteSTD(:,:,iff-6),iff)
1412          IF(vars_defined) THEN
1413             DO k=1, nlevSTD
1414                DO i=1, klon
1415                   IF(tnondef(i,k,iff-6).NE.missing_val) THEN
1416                      IF(freq_outNMC(iff-6).LT.0) THEN
1417                         freq_moyNMC(iff-6)=(mth_len*un_jour)/freq_calNMC(iff-6)
1418                      ELSE
1419                         freq_moyNMC(iff-6)=freq_outNMC(iff-6)/freq_calNMC(iff-6)
1420                      ENDIF
1421                      zx_tmp_fi3d_STD(i,k) = (100.*tnondef(i,k,iff-6))/freq_moyNMC(iff-6)
1422                   ELSE
1423                      zx_tmp_fi3d_STD(i,k) = missing_val
1424                   ENDIF
1425                ENDDO
1426             ENDDO
1427          ENDIF
1428          CALL histwrite_phy(o_psbg,zx_tmp_fi3d_STD,iff)
1429          IF(vars_defined) THEN
1430             DO k=1, nlevSTD
1431                DO i=1, klon
1432                   IF(O3sumSTD(i,k,iff-6).NE.missing_val) THEN
1433                      zx_tmp_fi3d_STD(i,k) = O3sumSTD(i,k,iff-6) * 1.e+9
1434                   ELSE
1435                      zx_tmp_fi3d_STD(i,k) = missing_val
1436                   ENDIF
1437                ENDDO
1438             ENDDO !k=1, nlevSTD
1439          ENDIF
1440          CALL histwrite_phy(o_tro3,zx_tmp_fi3d_STD,iff)
1441          if (read_climoz == 2) THEN
1442             IF(vars_defined) THEN
1443                DO k=1, nlevSTD
1444                   DO i=1, klon
1445                      IF(O3daysumSTD(i,k,iff-6).NE.missing_val) THEN
1446                         zx_tmp_fi3d_STD(i,k) = O3daysumSTD(i,k,iff-6) * 1.e+9
1447                      ELSE
1448                         zx_tmp_fi3d_STD(i,k) = missing_val
1449                      ENDIF
1450                   ENDDO
1451                ENDDO !k=1, nlevSTD
1452             ENDIF
1453             CALL histwrite_phy(o_tro3_daylight,zx_tmp_fi3d_STD,iff)
1454          endif
1455          CALL histwrite_phy(o_uxv,uvsumSTD(:,:,iff-6),iff)
1456          CALL histwrite_phy(o_vxq,vqsumSTD(:,:,iff-6),iff)
1457          CALL histwrite_phy(o_vxT,vTsumSTD(:,:,iff-6),iff)
1458          CALL histwrite_phy(o_wxq,wqsumSTD(:,:,iff-6),iff)
1459          CALL histwrite_phy(o_vxphi,vphisumSTD(:,:,iff-6),iff)
1460          CALL histwrite_phy(o_wxT,wTsumSTD(:,:,iff-6),iff)
1461          CALL histwrite_phy(o_uxu,u2sumSTD(:,:,iff-6),iff)
1462          CALL histwrite_phy(o_vxv,v2sumSTD(:,:,iff-6),iff)
1463          CALL histwrite_phy(o_TxT,T2sumSTD(:,:,iff-6),iff)
1464       ENDDO !nfiles
1465  ENDIF
1466#endif
1467#endif
1468#ifdef CPP_XIOS
1469  IF(ok_all_xml) THEN
1470!      DO iff=7, nfiles
1471
1472!         CALL histwrite_phy(o_tnondef,tnondef(:,:,3))
1473          CALL histwrite_phy(o_ta,tlevSTD(:,:))
1474          CALL histwrite_phy(o_zg,philevSTD(:,:))
1475          CALL histwrite_phy(o_hus,qlevSTD(:,:))
1476          CALL histwrite_phy(o_hur,rhlevSTD(:,:))
1477          CALL histwrite_phy(o_ua,ulevSTD(:,:))
1478          CALL histwrite_phy(o_va,vlevSTD(:,:))
1479          CALL histwrite_phy(o_wap,wlevSTD(:,:))
1480!         IF(vars_defined) THEN
1481!            DO k=1, nlevSTD
1482!               DO i=1, klon
1483!                  IF(tnondef(i,k,3).NE.missing_val) THEN
1484!                     IF(freq_outNMC(iff-6).LT.0) THEN
1485!                        freq_moyNMC(iff-6)=(mth_len*un_jour)/freq_calNMC(iff-6)
1486!                     ELSE
1487!                        freq_moyNMC(iff-6)=freq_outNMC(iff-6)/freq_calNMC(iff-6)
1488!                     ENDIF
1489!                     zx_tmp_fi3d_STD(i,k) = (100.*tnondef(i,k,3))/freq_moyNMC(iff-6)
1490!                  ELSE
1491!                     zx_tmp_fi3d_STD(i,k) = missing_val
1492!                  ENDIF
1493!               ENDDO
1494!            ENDDO
1495!         ENDIF
1496!         CALL histwrite_phy(o_psbg,zx_tmp_fi3d_STD)
1497          IF(vars_defined) THEN
1498             DO k=1, nlevSTD
1499                DO i=1, klon
1500                   IF(O3STD(i,k).NE.missing_val) THEN
1501                      zx_tmp_fi3d_STD(i,k) = O3STD(i,k) * 1.e+9
1502                   ELSE
1503                      zx_tmp_fi3d_STD(i,k) = missing_val
1504                   ENDIF
1505                ENDDO
1506             ENDDO !k=1, nlevSTD
1507          ENDIF
1508          CALL histwrite_phy(o_tro3,zx_tmp_fi3d_STD)
1509          if (read_climoz == 2) THEN
1510             IF(vars_defined) THEN
1511                DO k=1, nlevSTD
1512                   DO i=1, klon
1513                      IF(O3daySTD(i,k).NE.missing_val) THEN
1514                         zx_tmp_fi3d_STD(i,k) = O3daySTD(i,k) * 1.e+9
1515                      ELSE
1516                         zx_tmp_fi3d_STD(i,k) = missing_val
1517                      ENDIF
1518                   ENDDO
1519                ENDDO !k=1, nlevSTD
1520             ENDIF
1521             CALL histwrite_phy(o_tro3_daylight,zx_tmp_fi3d_STD)
1522          endif
1523          CALL histwrite_phy(o_uxv,uvSTD(:,:))
1524          CALL histwrite_phy(o_vxq,vqSTD(:,:))
1525          CALL histwrite_phy(o_vxT,vTSTD(:,:))
1526          CALL histwrite_phy(o_wxq,wqSTD(:,:))
1527          CALL histwrite_phy(o_vxphi,vphiSTD(:,:))
1528          CALL histwrite_phy(o_wxT,wTSTD(:,:))
1529          CALL histwrite_phy(o_uxu,u2STD(:,:))
1530          CALL histwrite_phy(o_vxv,v2STD(:,:))
1531          CALL histwrite_phy(o_TxT,T2STD(:,:))
1532!      ENDDO !nfiles
1533  ENDIF
1534#endif
1535!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1536        IF (nqtot.GE.nqo+1) THEN
1537            DO iq=nqo+1,nqtot
1538              IF (type_trac == 'lmdz' .OR. type_trac == 'repr') THEN
1539
1540!jyg<
1541!!             CALL histwrite_phy(o_trac(iq-nqo), qx(:,:,iq))
1542             CALL histwrite_phy(o_trac(iq-nqo), tr_seri(:,:,iq-nqo))
1543!>jyg
1544             CALL histwrite_phy(o_dtr_vdf(iq-nqo),d_tr_cl(:,:,iq-nqo))
1545             CALL histwrite_phy(o_dtr_the(iq-nqo),d_tr_th(:,:,iq-nqo))
1546             CALL histwrite_phy(o_dtr_con(iq-nqo),d_tr_cv(:,:,iq-nqo))
1547             CALL histwrite_phy(o_dtr_lessi_impa(iq-nqo),d_tr_lessi_impa(:,:,iq-nqo))
1548             CALL histwrite_phy(o_dtr_lessi_nucl(iq-nqo),d_tr_lessi_nucl(:,:,iq-nqo))
1549             CALL histwrite_phy(o_dtr_insc(iq-nqo),d_tr_insc(:,:,iq-nqo))
1550             CALL histwrite_phy(o_dtr_bcscav(iq-nqo),d_tr_bcscav(:,:,iq-nqo))
1551             CALL histwrite_phy(o_dtr_evapls(iq-nqo),d_tr_evapls(:,:,iq-nqo))
1552             CALL histwrite_phy(o_dtr_ls(iq-nqo),d_tr_ls(:,:,iq-nqo))
1553             CALL histwrite_phy(o_dtr_trsp(iq-nqo),d_tr_trsp(:,:,iq-nqo))
1554             CALL histwrite_phy(o_dtr_sscav(iq-nqo),d_tr_sscav(:,:,iq-nqo))
1555             CALL histwrite_phy(o_dtr_sat(iq-nqo),d_tr_sat(:,:,iq-nqo))
1556             CALL histwrite_phy(o_dtr_uscav(iq-nqo),d_tr_uscav(:,:,iq-nqo))
1557             zx_tmp_fi2d=0.
1558             IF(vars_defined) THEN
1559                DO k=1,klev
1560!jyg<
1561!!                   zx_tmp_fi2d(:)=zx_tmp_fi2d(:)+zmasse(:,k)*qx(:,k,iq)
1562                   zx_tmp_fi2d(:)=zx_tmp_fi2d(:)+zmasse(:,k)*tr_seri(:,k,iq-nqo)
1563!>jyg
1564                ENDDO
1565             ENDIF
1566             CALL histwrite_phy(o_trac_cum(iq-nqo), zx_tmp_fi2d)
1567             endif
1568          ENDDO
1569       ENDIF
1570
1571       IF(.NOT.vars_defined) THEN
1572          !$OMP MASTER
1573#ifndef CPP_IOIPSL_NO_OUTPUT
1574          DO iff=1,nfiles
1575             IF (clef_files(iff)) THEN
1576                CALL histend(nid_files(iff))
1577                ndex2d = 0
1578                ndex3d = 0
1579
1580             ENDIF ! clef_files
1581          ENDDO !  iff
1582#endif
1583#ifdef CPP_XIOS
1584          !On finalise l'initialisation:
1585          CALL wxios_closedef()
1586#endif
1587
1588          !$OMP END MASTER
1589          !$OMP BARRIER
1590          vars_defined = .TRUE.
1591
1592       END IF
1593
1594    END DO
1595
1596    IF(vars_defined) THEN
1597       ! On synchronise les fichiers pour IOIPSL
1598#ifndef CPP_IOIPSL_NO_OUTPUT
1599       !$OMP MASTER
1600       DO iff=1,nfiles
1601          IF (ok_sync .AND. clef_files(iff)) THEN
1602             CALL histsync(nid_files(iff))
1603          ENDIF
1604       END DO
1605       !$OMP END MASTER
1606#endif
1607    ENDIF
1608
1609  END SUBROUTINE phys_output_write
1610
1611END MODULE phys_output_write_mod
Note: See TracBrowser for help on using the repository browser.