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

Last change on this file since 2379 was 2377, checked in by musat, 9 years ago

Correct zfull and zhalf CMIP6' diagnostics :
zhalf correspond to the interface levels (and paprs) and include

both the top of the model atmosphere and surface levels

zfull is the actual height above mean sea level and corresponds to the

geopotential height (and pplay levels)

ASima/AI/IM

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