source: LMDZ5/branches/LMDZ6_rc0/libf/phylmd/phys_output_write_mod.F90 @ 2272

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

Merged trunk changes -r2070:2158 into testing branch. Compilation problems introduced by revision r2155 have been corrected by hand

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