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

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

Further modifications to enforce physics/dynamics separation:

  • moved iniprint.h and misc_mod back to dyn3d_common, as these should only be used by dynamics.
  • created print_control_mod in the physics to store flags prt_level, lunout, debug to be local to physics (should be used rather than iniprint.h)
  • created abort_physic.F90 , which does the same job as abort_gcm() did, but should be used instead when in physics.
  • reactivated inifis (turned it into a module, inifis_mod.F90) to initialize physical constants and print_control_mod flags.

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