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

Last change on this file since 2114 was 2114, checked in by Laurent Fairhead, 10 years ago

Suite (et fin?) des modifications pour permettre un controle 'tout xml'
des fichiers de sorties XIOS.
Un nouveau paramètre logique est introduit, ok_all_xml, false par défaut, et lu
dans le run.def qui permet de faire du 'tout xml'


Follow-up modifications to ensure total xlm control over the output files
from XIOS.
A new logical parameter, ok_all_wml, is introduced. False by default, it is
read from the run.def file and, if true, will give over control to the
XIOS xml files

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