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

Last change on this file since 2106 was 2103, checked in by musat, 10 years ago

Ajout variables t2m_min_mon et t2m_max_mon dans les fichier
histmth. Ce sont les tempe©rature a minimales et
maximales journalie¨re) moyeneées sur le mo. Le calcul
est fait "manuelleÃment" dans calcul_divers.h.
Teste en sequentiel (en local) avec calendrier a 360 et 365 jours
et en parallele MPI avec calendrier 360 jours (sur adapp).
IM

  • This line, and those below, will be ignored--

M phys_output_ctrlout_mod.F90
M phys_local_var_mod.F90
M phys_output_write_mod.F90
M calcul_divers.h

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