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

Last change on this file since 2284 was 2284, checked in by jyg, 9 years ago

Bug fixing concerning aerosol scavenging (bugs in
cvltr_scav and in lsc_scav) and aerosol
concentration output (the values at beginning,
instead of end, of time step were output).

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