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

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

Changement de nom de clef CPP:
CPP_NO_IOIPSL devient CPP_IOIPSL_NO_OUTPUT pour éviter la confusion. Elle
permet de ne pas sortir les fichiers IOIPSL "proprement"
L'option -io de makelmdz et makelmdz_fcm est changée:
avec la valeur ioipsl, on ne sort que les fichier IOIPSL

mix, on sort les fichiers IOIPSL et XIOS
xios, on ne sort que les fichiers XIOS


Change in the name of a CPP key:
CPP_NO_IOIPSL becomes CPP_IOIPSL_NO_OUTPUT. If defined, IOIPSL outputs are not
generated.
The -io option for makelmdz and makelmdz_fcm is changed as well:
with the value ioipsl, only IOIPSL files are output

mix, IOIPSL and XIOS files are output
xios, only XIOS files are output

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