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

Last change on this file since 2236 was 2209, checked in by Ehouarn Millour, 9 years ago

Update of the slab ocean by Francis Codron. There are now 3 possibilities for the "version_ocean" slab type:
sicOBS = prescribed ice fraction. Water temperature nearby is set to -1.8°C and cannot become lower.
sicNO = ignore sea ice. One can prescribe a fraction, but the nearby ocean evolves freely, depending on surface fluxes: temperature can go below freezing point or above...
sicINT = interactive sea ice.
EM

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