source: LMDZ4/branches/LMDZ4-dev-20091210/libf/phylmd/phys_output_mod.F90 @ 2063

Last change on this file since 2063 was 1266, checked in by Laurent Fairhead, 15 years ago

Correction d'un bug dans aeropt_5wv.F90 : une variable était utilisée sans
être calculée
Ajout de tests sur l'écriture de variables dans les fichiers de sortie
Changement des arguments d'appels à INCA
ACo

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 71.1 KB
Line 
1!
2! $Id: phys_output_mod.F90 1266 2009-11-20 16:27:40Z fairhead $
3!
4! Abderrahmane 12 2007
5!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
6!!! Ecreture des Sorties du modele dans les fichiers Netcdf :
7! histmth.nc : moyennes mensuelles
8! histday.nc : moyennes journalieres
9! histhf.nc  : moyennes toutes les 3 heures
10! histins.nc : valeurs instantanees
11!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
12
13MODULE phys_output_mod
14
15  IMPLICIT NONE
16
17  private histdef2d, histdef3d, conf_physoutputs
18
19
20   integer, parameter                           :: nfiles = 5
21   logical, dimension(nfiles), save             :: clef_files
22   integer, dimension(nfiles), save             :: lev_files
23   integer, dimension(nfiles), save             :: nid_files
24!!$OMP THREADPRIVATE(clef_files, lev_files,nid_files)
25 
26   integer, dimension(nfiles), private, save :: nhorim, nvertm
27   integer, dimension(nfiles), private, save :: nvertap, nvertbp, nvertAlt
28!   integer, dimension(nfiles), private, save :: nvertp0
29   real, dimension(nfiles), private, save                :: zoutm
30   real,                    private, save                :: zdtime
31   CHARACTER(len=20), dimension(nfiles), private, save   :: type_ecri
32!$OMP THREADPRIVATE(nhorim, nvertm, zoutm,zdtime,type_ecri)
33
34!   integer, save                     :: nid_hf3d
35
36!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
37!! Definition pour chaque variable du niveau d ecriture dans chaque fichier
38!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!/ histmth, histday, histhf, histins /),'!!!!!!!!!!!!
39!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
40
41  integer, private:: levmin(nfiles) = 1
42  integer, private:: levmax(nfiles)
43
44  TYPE ctrl_out
45   integer,dimension(5) :: flag
46   character(len=20)     :: name
47  END TYPE ctrl_out
48
49!!! Comosentes de la coordonnee sigma-hybride
50!!! Ap et Bp
51  type(ctrl_out),save :: o_Ahyb         = ctrl_out((/ 1, 1, 1, 1, 1 /), 'Ap')
52  type(ctrl_out),save :: o_Bhyb         = ctrl_out((/ 1, 1, 1, 1, 1 /), 'Bp')
53  type(ctrl_out),save :: o_Alt          = ctrl_out((/ 1, 1, 1, 1, 1 /), 'Alt')
54
55!!! 1D
56  type(ctrl_out),save :: o_phis         = ctrl_out((/ 1, 1, 10, 1, 1 /), 'phis')
57  type(ctrl_out),save :: o_aire         = ctrl_out((/ 1, 1, 10,  1, 1 /),'aire')
58  type(ctrl_out),save :: o_contfracATM  = ctrl_out((/ 10, 1,  1, 10, 10 /),'contfracATM')
59  type(ctrl_out),save :: o_contfracOR   = ctrl_out((/ 10, 1,  1, 10, 10 /),'contfracOR')
60  type(ctrl_out),save :: o_aireTER      = ctrl_out((/ 10, 10, 1, 10, 10 /),'aireTER')
61 
62!!! 2D
63  type(ctrl_out),save :: o_flat         = ctrl_out((/ 10, 1, 10, 10, 1 /),'flat')
64  type(ctrl_out),save :: o_slp          = ctrl_out((/ 1, 1, 1, 10, 1 /),'slp')
65  type(ctrl_out),save :: o_tsol         = ctrl_out((/ 1, 1, 1, 1, 1 /),'tsol')
66  type(ctrl_out),save :: o_t2m          = ctrl_out((/ 1, 1, 1, 1, 1 /),'t2m')
67  type(ctrl_out),save :: o_t2m_min      = ctrl_out((/ 1, 1, 10, 10, 10 /),'t2m_min')
68  type(ctrl_out),save :: o_t2m_max      = ctrl_out((/ 1, 1, 10, 10, 10 /),'t2m_max')
69  type(ctrl_out),save,dimension(4) :: o_t2m_srf      = (/ ctrl_out((/ 10, 4, 10, 10, 10 /),'t2m_ter'), &
70                                                 ctrl_out((/ 10, 4, 10, 10, 10 /),'t2m_lic'), &
71                                                 ctrl_out((/ 10, 4, 10, 10, 10 /),'t2m_oce'), &
72                                                 ctrl_out((/ 10, 4, 10, 10, 10 /),'t2m_sic') /)
73
74  type(ctrl_out),save :: o_wind10m      = ctrl_out((/ 1, 1, 1, 10, 10 /),'wind10m')
75  type(ctrl_out),save :: o_wind10max    = ctrl_out((/ 10, 1, 10, 10, 10 /),'wind10max')
76  type(ctrl_out),save :: o_sicf         = ctrl_out((/ 1, 1, 10, 10, 10 /),'sicf')
77  type(ctrl_out),save :: o_q2m          = ctrl_out((/ 1, 1, 1, 1, 1 /),'q2m')
78  type(ctrl_out),save :: o_u10m         = ctrl_out((/ 1, 1, 1, 1, 1 /),'u10m')
79  type(ctrl_out),save :: o_v10m         = ctrl_out((/ 1, 1, 1, 1, 1 /),'v10m')
80  type(ctrl_out),save :: o_psol         = ctrl_out((/ 1, 1, 1, 1, 1 /),'psol')
81  type(ctrl_out),save :: o_qsurf        = ctrl_out((/ 1, 10, 10, 10, 10 /),'qsurf')
82
83  type(ctrl_out),save,dimension(4) :: o_u10m_srf     = (/ ctrl_out((/ 10, 4, 10, 10, 10 /),'u10m_ter'), &
84                                              ctrl_out((/ 10, 4, 10, 10, 10 /),'u10m_lic'), &
85                                              ctrl_out((/ 10, 4, 10, 10, 10 /),'u10m_oce'), &
86                                              ctrl_out((/ 10, 4, 10, 10, 10 /),'u10m_sic') /)
87
88  type(ctrl_out),save,dimension(4) :: o_v10m_srf     = (/ ctrl_out((/ 10, 4, 10, 10, 10 /),'v10m_ter'), &
89                                              ctrl_out((/ 10, 4, 10, 10, 10 /),'v10m_lic'), &
90                                              ctrl_out((/ 10, 4, 10, 10, 10 /),'v10m_oce'), &
91                                              ctrl_out((/ 10, 4, 10, 10, 10 /),'v10m_sic') /)
92
93  type(ctrl_out),save :: o_qsol         = ctrl_out((/ 1, 10, 10, 1, 1 /),'qsol')
94
95  type(ctrl_out),save :: o_ndayrain     = ctrl_out((/ 1, 10, 10, 10, 10 /),'ndayrain')
96  type(ctrl_out),save :: o_precip       = ctrl_out((/ 1, 1, 1, 1, 1 /),'precip')
97  type(ctrl_out),save :: o_plul         = ctrl_out((/ 1, 1, 1, 1, 10 /),'plul')
98
99  type(ctrl_out),save :: o_pluc         = ctrl_out((/ 1, 1, 1, 1, 10 /),'pluc')
100  type(ctrl_out),save :: o_snow         = ctrl_out((/ 1, 1, 10, 1, 10 /),'snow')
101  type(ctrl_out),save :: o_evap         = ctrl_out((/ 1, 1, 10, 1, 10 /),'evap')
102  type(ctrl_out),save :: o_tops         = ctrl_out((/ 1, 1, 10, 10, 10 /),'tops')
103  type(ctrl_out),save :: o_tops0        = ctrl_out((/ 1, 5, 10, 10, 10 /),'tops0')
104  type(ctrl_out),save :: o_topl         = ctrl_out((/ 1, 1, 10, 1, 10 /),'topl')
105  type(ctrl_out),save :: o_topl0        = ctrl_out((/ 1, 5, 10, 10, 10 /),'topl0')
106  type(ctrl_out),save :: o_SWupTOA      = ctrl_out((/ 1, 4, 10, 10, 10 /),'SWupTOA')
107  type(ctrl_out),save :: o_SWupTOAclr   = ctrl_out((/ 1, 4, 10, 10, 10 /),'SWupTOAclr')
108  type(ctrl_out),save :: o_SWdnTOA      = ctrl_out((/ 1, 4, 10, 10, 10 /),'SWdnTOA')
109  type(ctrl_out),save :: o_SWdnTOAclr   = ctrl_out((/ 1, 4, 10, 10, 10 /),'SWdnTOAclr')
110  type(ctrl_out),save :: o_SWup200      = ctrl_out((/ 1, 10, 10, 10, 10 /),'SWup200')
111  type(ctrl_out),save :: o_SWup200clr   = ctrl_out((/ 10, 1, 10, 10, 10 /),'SWup200clr')
112  type(ctrl_out),save :: o_SWdn200      = ctrl_out((/ 1, 10, 10, 10, 10 /),'SWdn200')
113  type(ctrl_out),save :: o_SWdn200clr   = ctrl_out((/ 10, 1, 10, 10, 10 /),'SWdn200clr')
114
115! arajouter
116!  type(ctrl_out),save :: o_LWupTOA     = ctrl_out((/ 1, 4, 10, 10, 10 /),'LWupTOA')
117!  type(ctrl_out),save :: o_LWupTOAclr  = ctrl_out((/ 1, 4, 10, 10, 10 /),'LWupTOAclr')
118!  type(ctrl_out),save :: o_LWdnTOA     = ctrl_out((/ 1, 4, 10, 10, 10 /),'LWdnTOA')
119!  type(ctrl_out),save :: o_LWdnTOAclr  = ctrl_out((/ 1, 4, 10, 10, 10 /),'LWdnTOAclr')
120
121  type(ctrl_out),save :: o_LWup200      = ctrl_out((/ 1, 10, 10, 10, 10 /),'LWup200')
122  type(ctrl_out),save :: o_LWup200clr   = ctrl_out((/ 1, 10, 10, 10, 10 /),'LWup200clr')
123  type(ctrl_out),save :: o_LWdn200      = ctrl_out((/ 1, 10, 10, 10, 10 /),'LWdn200')
124  type(ctrl_out),save :: o_LWdn200clr   = ctrl_out((/ 1, 10, 10, 10, 10 /),'LWdn200clr')
125  type(ctrl_out),save :: o_sols         = ctrl_out((/ 1, 1, 10, 1, 10 /),'sols')
126  type(ctrl_out),save :: o_sols0        = ctrl_out((/ 1, 5, 10, 10, 10 /),'sols0')
127  type(ctrl_out),save :: o_soll         = ctrl_out((/ 1, 1, 10, 1, 10 /),'soll')
128  type(ctrl_out),save :: o_soll0        = ctrl_out((/ 1, 5, 10, 10, 10 /),'soll0')
129  type(ctrl_out),save :: o_radsol       = ctrl_out((/ 1, 1, 10, 10, 10 /),'radsol')
130  type(ctrl_out),save :: o_SWupSFC      = ctrl_out((/ 1, 4, 10, 10, 10 /),'SWupSFC')
131  type(ctrl_out),save :: o_SWupSFCclr   = ctrl_out((/ 1, 4, 10, 10, 10 /),'SWupSFCclr')
132  type(ctrl_out),save :: o_SWdnSFC      = ctrl_out((/ 1, 1, 10, 10, 10 /),'SWdnSFC')
133  type(ctrl_out),save :: o_SWdnSFCclr   = ctrl_out((/ 1, 4, 10, 10, 10 /),'SWdnSFCclr')
134  type(ctrl_out),save :: o_LWupSFC      = ctrl_out((/ 1, 4, 10, 10, 10 /),'LWupSFC')
135  type(ctrl_out),save :: o_LWupSFCclr   = ctrl_out((/ 1, 4, 10, 10, 10 /),'LWupSFCclr')
136  type(ctrl_out),save :: o_LWdnSFC      = ctrl_out((/ 1, 4, 10, 10, 10 /),'LWdnSFC')
137  type(ctrl_out),save :: o_LWdnSFCclr   = ctrl_out((/ 1, 4, 10, 10, 10 /),'LWdnSFCclr')
138  type(ctrl_out),save :: o_bils         = ctrl_out((/ 1, 2, 10, 1, 10 /),'bils')
139  type(ctrl_out),save :: o_sens         = ctrl_out((/ 1, 1, 10, 1, 1 /),'sens')
140  type(ctrl_out),save :: o_fder         = ctrl_out((/ 1, 2, 10, 1, 10 /),'fder')
141  type(ctrl_out),save :: o_ffonte       = ctrl_out((/ 1, 10, 10, 10, 10 /),'ffonte')
142  type(ctrl_out),save :: o_fqcalving    = ctrl_out((/ 1, 10, 10, 10, 10 /),'fqcalving')
143  type(ctrl_out),save :: o_fqfonte      = ctrl_out((/ 1, 10, 10, 10, 10 /),'fqfonte')
144
145  type(ctrl_out),save,dimension(4) :: o_taux_srf     = (/ ctrl_out((/ 1, 4, 10, 1, 10 /),'taux_ter'), &
146                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'taux_lic'), &
147                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'taux_oce'), &
148                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'taux_sic') /)
149
150  type(ctrl_out),save,dimension(4) :: o_tauy_srf     = (/ ctrl_out((/ 1, 4, 10, 1, 10 /),'tauy_ter'), &
151                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'tauy_lic'), &
152                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'tauy_oce'), &
153                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'tauy_sic') /)
154
155
156  type(ctrl_out),save,dimension(4) :: o_pourc_srf    = (/ ctrl_out((/ 1, 4, 10, 1, 10 /),'pourc_ter'), &
157                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'pourc_lic'), &
158                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'pourc_oce'), &
159                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'pourc_sic') /)     
160
161  type(ctrl_out),save,dimension(4) :: o_fract_srf    = (/ ctrl_out((/ 1, 4, 10, 1, 10 /),'fract_ter'), &
162                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'fract_lic'), &
163                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'fract_oce'), &
164                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'fract_sic') /)
165
166  type(ctrl_out),save,dimension(4) :: o_tsol_srf     = (/ ctrl_out((/ 1, 4, 10, 1, 10 /),'tsol_ter'), &
167                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'tsol_lic'), &
168                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'tsol_oce'), &
169                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'tsol_sic') /)
170
171  type(ctrl_out),save,dimension(4) :: o_sens_srf     = (/ ctrl_out((/ 1, 4, 10, 1, 10 /),'sens_ter'), &
172                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'sens_lic'), &
173                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'sens_oce'), &
174                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'sens_sic') /)
175
176  type(ctrl_out),save,dimension(4) :: o_lat_srf      = (/ ctrl_out((/ 1, 4, 10, 1, 10 /),'lat_ter'), &
177                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'lat_lic'), &
178                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'lat_oce'), &
179                                                 ctrl_out((/ 1, 4, 10, 1, 10 /),'lat_sic') /)
180
181  type(ctrl_out),save,dimension(4) :: o_flw_srf      = (/ ctrl_out((/ 1, 10, 10, 10, 10 /),'flw_ter'), &
182                                                 ctrl_out((/ 1, 10, 10, 10, 10 /),'flw_lic'), &
183                                                 ctrl_out((/ 1, 10, 10, 10, 10 /),'flw_oce'), &
184                                                 ctrl_out((/ 1, 10, 10, 10, 10 /),'flw_sic') /)
185                                                 
186  type(ctrl_out),save,dimension(4) :: o_fsw_srf      = (/ ctrl_out((/ 1, 10, 10, 10, 10 /),'fsw_ter'), &
187                                                  ctrl_out((/ 1, 10, 10, 10, 10 /),'fsw_lic'), &
188                                                  ctrl_out((/ 1, 10, 10, 10, 10 /),'fsw_oce'), &
189                                                  ctrl_out((/ 1, 10, 10, 10, 10 /),'fsw_sic') /)
190
191  type(ctrl_out),save,dimension(4) :: o_wbils_srf    = (/ ctrl_out((/ 1, 10, 10, 10, 10 /),'wbils_ter'), &
192                                                 ctrl_out((/ 1, 10, 10, 10, 10 /),'wbils_lic'), &
193                                                 ctrl_out((/ 1, 10, 10, 10, 10 /),'wbils_oce'), &
194                                                 ctrl_out((/ 1, 10, 10, 10, 10 /),'wbils_sic') /)
195
196  type(ctrl_out),save,dimension(4) :: o_wbilo_srf    = (/ ctrl_out((/ 1, 10, 10, 10, 10 /),'wbilo_ter'), &
197                                                     ctrl_out((/ 1, 10, 10, 10, 10 /),'wbilo_lic'), &
198                                                 ctrl_out((/ 1, 10, 10, 10, 10 /),'wbilo_oce'), &
199                                                 ctrl_out((/ 1, 10, 10, 10, 10 /),'wbilo_sic') /)
200
201
202  type(ctrl_out),save :: o_cdrm         = ctrl_out((/ 1, 10, 10, 1, 10 /),'cdrm')
203  type(ctrl_out),save :: o_cdrh         = ctrl_out((/ 1, 10, 10, 1, 10 /),'cdrh')
204  type(ctrl_out),save :: o_cldl         = ctrl_out((/ 1, 1, 10, 10, 10 /),'cldl')
205  type(ctrl_out),save :: o_cldm         = ctrl_out((/ 1, 1, 10, 10, 10 /),'cldm')
206  type(ctrl_out),save :: o_cldh         = ctrl_out((/ 1, 1, 10, 10, 10 /),'cldh')
207  type(ctrl_out),save :: o_cldt         = ctrl_out((/ 1, 1, 2, 10, 10 /),'cldt')
208  type(ctrl_out),save :: o_cldq         = ctrl_out((/ 1, 1, 10, 10, 10 /),'cldq')
209  type(ctrl_out),save :: o_lwp          = ctrl_out((/ 1, 5, 10, 10, 10 /),'lwp')
210  type(ctrl_out),save :: o_iwp          = ctrl_out((/ 1, 5, 10, 10, 10 /),'iwp')
211  type(ctrl_out),save :: o_ue           = ctrl_out((/ 1, 10, 10, 10, 10 /),'ue')
212  type(ctrl_out),save :: o_ve           = ctrl_out((/ 1, 10, 10, 10, 10 /),'ve')
213  type(ctrl_out),save :: o_uq           = ctrl_out((/ 1, 10, 10, 10, 10 /),'uq')
214  type(ctrl_out),save :: o_vq           = ctrl_out((/ 1, 10, 10, 10, 10 /),'vq')
215 
216  type(ctrl_out),save :: o_cape         = ctrl_out((/ 1, 10, 10, 10, 10 /),'cape')
217  type(ctrl_out),save :: o_pbase        = ctrl_out((/ 1, 10, 10, 10, 10 /),'pbase')
218  type(ctrl_out),save :: o_ptop         = ctrl_out((/ 1, 4, 10, 10, 10 /),'ptop')
219  type(ctrl_out),save :: o_fbase        = ctrl_out((/ 1, 10, 10, 10, 10 /),'fbase')
220  type(ctrl_out),save :: o_prw          = ctrl_out((/ 1, 1, 10, 10, 10 /),'prw')
221
222  type(ctrl_out),save :: o_s_pblh       = ctrl_out((/ 1, 10, 10, 1, 1 /),'s_pblh')
223  type(ctrl_out),save :: o_s_pblt       = ctrl_out((/ 1, 10, 10, 1, 1 /),'s_pblt')
224  type(ctrl_out),save :: o_s_lcl        = ctrl_out((/ 1, 10, 10, 1, 10 /),'s_lcl')
225  type(ctrl_out),save :: o_s_capCL      = ctrl_out((/ 1, 10, 10, 1, 10 /),'s_capCL')
226  type(ctrl_out),save :: o_s_oliqCL     = ctrl_out((/ 1, 10, 10, 1, 10 /),'s_oliqCL')
227  type(ctrl_out),save :: o_s_cteiCL     = ctrl_out((/ 1, 10, 10, 1, 1 /),'s_cteiCL')
228  type(ctrl_out),save :: o_s_therm      = ctrl_out((/ 1, 10, 10, 1, 1 /),'s_therm')
229  type(ctrl_out),save :: o_s_trmb1      = ctrl_out((/ 1, 10, 10, 1, 10 /),'s_trmb1')
230  type(ctrl_out),save :: o_s_trmb2      = ctrl_out((/ 1, 10, 10, 1, 10 /),'s_trmb2')
231  type(ctrl_out),save :: o_s_trmb3      = ctrl_out((/ 1, 10, 10, 1, 10 /),'s_trmb3')
232
233  type(ctrl_out),save :: o_slab_bils    = ctrl_out((/ 1, 1, 10, 10, 10 /),'slab_bils_oce')
234
235  type(ctrl_out),save :: o_ale_bl       = ctrl_out((/ 1, 1, 1, 1, 10 /),'ale_bl')
236  type(ctrl_out),save :: o_alp_bl       = ctrl_out((/ 1, 1, 1, 1, 10 /),'alp_bl')
237  type(ctrl_out),save :: o_ale_wk       = ctrl_out((/ 1, 1, 1, 1, 10 /),'ale_wk')
238  type(ctrl_out),save :: o_alp_wk       = ctrl_out((/ 1, 1, 1, 1, 10 /),'alp_wk')
239
240  type(ctrl_out),save :: o_ale          = ctrl_out((/ 1, 1, 1, 1, 10 /),'ale')
241  type(ctrl_out),save :: o_alp          = ctrl_out((/ 1, 1, 1, 1, 10 /),'alp')
242  type(ctrl_out),save :: o_cin          = ctrl_out((/ 1, 1, 1, 1, 10 /),'cin')
243  type(ctrl_out),save :: o_wape         = ctrl_out((/ 1, 1, 1, 1, 10 /),'wape')
244
245
246! Champs interpolles sur des niveaux de pression ??? a faire correctement
247                                             
248  type(ctrl_out),save,dimension(6) :: o_uSTDlevs     = (/ ctrl_out((/ 1, 1, 3, 10, 10 /),'u850'), &
249                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'u700'), &
250                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'u500'), &
251                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'u200'), &
252                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'u50'), &
253                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'u10') /)
254                                                     
255
256  type(ctrl_out),save,dimension(6) :: o_vSTDlevs     = (/ ctrl_out((/ 1, 1, 3, 10, 10 /),'v850'), &
257                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'v700'), &
258                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'v500'), &
259                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'v200'), &
260                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'v50'), &
261                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'v10') /)
262
263  type(ctrl_out),save,dimension(6) :: o_wSTDlevs     = (/ ctrl_out((/ 1, 1, 3, 10, 10 /),'w850'), &
264                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'w700'), &
265                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'w500'), &
266                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'w200'), &
267                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'w50'), &
268                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'w10') /)
269
270  type(ctrl_out),save,dimension(6) :: o_tSTDlevs     = (/ ctrl_out((/ 1, 1, 3, 10, 10 /),'t850'), &
271                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'t700'), &
272                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'t500'), &
273                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'t200'), &
274                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'t50'), &
275                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'t10') /)
276
277  type(ctrl_out),save,dimension(6) :: o_qSTDlevs     = (/ ctrl_out((/ 1, 1, 3, 10, 10 /),'q850'), &
278                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'q700'), &
279                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'q500'), &
280                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'q200'), &
281                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'q50'), &
282                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'q10') /)
283
284  type(ctrl_out),save,dimension(6) :: o_phiSTDlevs   = (/ ctrl_out((/ 1, 1, 3, 10, 10 /),'phi850'), &
285                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'phi700'), &
286                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'phi500'), &
287                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'phi200'), &
288                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'phi50'), &
289                                                     ctrl_out((/ 1, 1, 3, 10, 10 /),'phi10') /)
290
291
292  type(ctrl_out),save :: o_t_oce_sic    = ctrl_out((/ 1, 10, 10, 10, 10 /),'t_oce_sic')
293
294  type(ctrl_out),save :: o_weakinv      = ctrl_out((/ 10, 1, 10, 10, 10 /),'weakinv')
295  type(ctrl_out),save :: o_dthmin       = ctrl_out((/ 10, 1, 10, 10, 10 /),'dthmin')
296  type(ctrl_out),save,dimension(4) :: o_u10_srf      = (/ ctrl_out((/ 10, 4, 10, 10, 10 /),'u10_ter'), &
297                                                     ctrl_out((/ 10, 4, 10, 10, 10 /),'u10_lic'), &
298                                                     ctrl_out((/ 10, 4, 10, 10, 10 /),'u10_oce'), &
299                                                     ctrl_out((/ 10, 4, 10, 10, 10 /),'u10_sic') /)
300
301  type(ctrl_out),save,dimension(4) :: o_v10_srf      = (/ ctrl_out((/ 10, 4, 10, 10, 10 /),'v10_ter'), &
302                                                     ctrl_out((/ 10, 4, 10, 10, 10 /),'v10_lic'), &
303                                                     ctrl_out((/ 10, 4, 10, 10, 10 /),'v10_oce'), &
304                                                     ctrl_out((/ 10, 4, 10, 10, 10 /),'v10_sic') /)
305
306  type(ctrl_out),save :: o_cldtau       = ctrl_out((/ 10, 5, 10, 10, 10 /),'cldtau')                     
307  type(ctrl_out),save :: o_cldemi       = ctrl_out((/ 10, 5, 10, 10, 10 /),'cldemi')
308  type(ctrl_out),save :: o_rh2m         = ctrl_out((/ 10, 5, 10, 10, 10 /),'rh2m')
309  type(ctrl_out),save :: o_qsat2m       = ctrl_out((/ 10, 5, 10, 10, 10 /),'qsat2m')
310  type(ctrl_out),save :: o_tpot         = ctrl_out((/ 10, 5, 10, 10, 10 /),'tpot')
311  type(ctrl_out),save :: o_tpote        = ctrl_out((/ 10, 5, 10, 10, 10 /),'tpote')
312  type(ctrl_out),save :: o_tke          = ctrl_out((/ 4, 10, 10, 10, 10 /),'tke ')
313  type(ctrl_out),save :: o_tke_max      = ctrl_out((/ 4, 10, 10, 10, 10 /),'tke_max')
314
315  type(ctrl_out),save,dimension(4) :: o_tke_srf      = (/ ctrl_out((/ 10, 4, 10, 10, 10 /),'tke_ter'), &
316                                                     ctrl_out((/ 10, 4, 10, 10, 10 /),'tke_lic'), &
317                                                     ctrl_out((/ 10, 4, 10, 10, 10 /),'tke_oce'), &
318                                                     ctrl_out((/ 10, 4, 10, 10, 10 /),'tke_sic') /)
319
320  type(ctrl_out),save,dimension(4) :: o_tke_max_srf  = (/ ctrl_out((/ 10, 4, 10, 10, 10 /),'tke_max_ter'), &
321                                                     ctrl_out((/ 10, 4, 10, 10, 10 /),'tke_max_lic'), &
322                                                     ctrl_out((/ 10, 4, 10, 10, 10 /),'tke_max_oce'), &
323                                                     ctrl_out((/ 10, 4, 10, 10, 10 /),'tke_max_sic') /)
324
325  type(ctrl_out),save :: o_kz           = ctrl_out((/ 4, 10, 10, 10, 10 /),'kz')
326  type(ctrl_out),save :: o_kz_max       = ctrl_out((/ 4, 10, 10, 10, 10 /),'kz_max')
327  type(ctrl_out),save :: o_SWnetOR      = ctrl_out((/ 10, 10, 2, 10, 10 /),'SWnetOR')
328  type(ctrl_out),save :: o_SWdownOR     = ctrl_out((/ 10, 10, 2, 10, 10 /),'SWdownOR')
329  type(ctrl_out),save :: o_LWdownOR     = ctrl_out((/ 10, 10, 2, 10, 10 /),'LWdownOR')
330
331  type(ctrl_out),save :: o_snowl        = ctrl_out((/ 10, 1, 10, 10, 10 /),'snowl')
332  type(ctrl_out),save :: o_cape_max     = ctrl_out((/ 10, 1, 10, 10, 10 /),'cape_max')
333  type(ctrl_out),save :: o_solldown     = ctrl_out((/ 10, 1, 10, 1, 10 /),'solldown')
334
335  type(ctrl_out),save :: o_dtsvdfo      = ctrl_out((/ 10, 10, 10, 1, 10 /),'dtsvdfo')
336  type(ctrl_out),save :: o_dtsvdft      = ctrl_out((/ 10, 10, 10, 1, 10 /),'dtsvdft')
337  type(ctrl_out),save :: o_dtsvdfg      = ctrl_out((/ 10, 10, 10, 1, 10 /),'dtsvdfg')
338  type(ctrl_out),save :: o_dtsvdfi      = ctrl_out((/ 10, 10, 10, 1, 10 /),'dtsvdfi')
339  type(ctrl_out),save :: o_rugs         = ctrl_out((/ 10, 10, 10, 1, 1 /),'rugs')
340
341  type(ctrl_out),save :: o_topswad      = ctrl_out((/ 4, 10, 10, 10, 10 /),'topswad')
342  type(ctrl_out),save :: o_topswai      = ctrl_out((/ 4, 10, 10, 10, 10 /),'topswai')
343  type(ctrl_out),save :: o_solswad      = ctrl_out((/ 4, 10, 10, 10, 10 /),'solswad')
344  type(ctrl_out),save :: o_solswai      = ctrl_out((/ 4, 10, 10, 10, 10 /),'solswai')
345
346  type(ctrl_out),save,dimension(10) :: o_tausumaero  = (/ ctrl_out((/ 4, 10, 10, 10, 10 /),'OD550_ASBCM'), &
347                                                     ctrl_out((/ 4, 10, 10, 10, 10 /),'OD550_ASPOMM'), &
348                                                     ctrl_out((/ 4, 10, 10, 10, 10 /),'OD550_ASSO4M'), &
349                                                     ctrl_out((/ 4, 10, 10, 10, 10 /),'OD550_CSSO4M'), &
350                                                     ctrl_out((/ 4, 10, 10, 10, 10 /),'OD550_SSSSM'), &
351                                                     ctrl_out((/ 4, 10, 10, 10, 10 /),'OD550_ASSSM'), &
352                                                     ctrl_out((/ 4, 10, 10, 10, 10 /),'OD550_CSSSM'), &
353                                                     ctrl_out((/ 4, 10, 10, 10, 10 /),'OD550_CIDUSTM'), &
354                                                     ctrl_out((/ 4, 10, 10, 10, 10 /),'OD550_AIBCM'), &
355                                                     ctrl_out((/ 4, 10, 10, 10, 10 /),'OD550_AIPOMM') /)
356
357  type(ctrl_out),save :: o_swtoaas_nat      = ctrl_out((/ 4, 10, 10, 10, 10 /),'swtoaas_nat')
358  type(ctrl_out),save :: o_swsrfas_nat      = ctrl_out((/ 4, 10, 10, 10, 10 /),'swsrfas_nat')
359  type(ctrl_out),save :: o_swtoacs_nat      = ctrl_out((/ 4, 10, 10, 10, 10 /),'swtoacs_nat')
360  type(ctrl_out),save :: o_swsrfcs_nat      = ctrl_out((/ 4, 10, 10, 10, 10 /),'swsrfcs_nat')
361
362  type(ctrl_out),save :: o_swtoaas_ant      = ctrl_out((/ 4, 10, 10, 10, 10 /),'swtoaas_ant')
363  type(ctrl_out),save :: o_swsrfas_ant      = ctrl_out((/ 4, 10, 10, 10, 10 /),'swsrfas_ant')
364  type(ctrl_out),save :: o_swtoacs_ant      = ctrl_out((/ 4, 10, 10, 10, 10 /),'swtoacs_ant')
365  type(ctrl_out),save :: o_swsrfcs_ant      = ctrl_out((/ 4, 10, 10, 10, 10 /),'swsrfcs_ant')
366
367  type(ctrl_out),save :: o_swtoacf_nat      = ctrl_out((/ 4, 10, 10, 10, 10 /),'swtoacf_nat')
368  type(ctrl_out),save :: o_swsrfcf_nat      = ctrl_out((/ 4, 10, 10, 10, 10 /),'swsrfcf_nat')
369  type(ctrl_out),save :: o_swtoacf_ant      = ctrl_out((/ 4, 10, 10, 10, 10 /),'swtoacf_ant')
370  type(ctrl_out),save :: o_swsrfcf_ant      = ctrl_out((/ 4, 10, 10, 10, 10 /),'swsrfcf_ant')
371  type(ctrl_out),save :: o_swtoacf_zero      = ctrl_out((/ 4, 10, 10, 10, 10 /),'swtoacf_zero')
372  type(ctrl_out),save :: o_swsrfcf_zero      = ctrl_out((/ 4, 10, 10, 10, 10 /),'swsrfcf_zero')
373
374
375!!!!!!!!!!!!!!!!!!!!!! 3D !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
376  type(ctrl_out),save :: o_lwcon        = ctrl_out((/ 2, 5, 10, 10, 1 /),'lwcon')
377  type(ctrl_out),save :: o_iwcon        = ctrl_out((/ 2, 5, 10, 10, 10 /),'iwcon')
378  type(ctrl_out),save :: o_temp         = ctrl_out((/ 2, 3, 4, 1, 1 /),'temp')
379  type(ctrl_out),save :: o_theta        = ctrl_out((/ 2, 3, 4, 1, 1 /),'theta')
380  type(ctrl_out),save :: o_ovap         = ctrl_out((/ 2, 3, 4, 1, 1 /),'ovap')
381  type(ctrl_out),save :: o_ovapinit         = ctrl_out((/ 2, 3, 4, 1, 1 /),'ovapinit')
382  type(ctrl_out),save :: o_wvapp        = ctrl_out((/ 2, 10, 10, 10, 10 /),'wvapp')
383  type(ctrl_out),save :: o_geop         = ctrl_out((/ 2, 3, 10, 1, 1 /),'geop')
384  type(ctrl_out),save :: o_vitu         = ctrl_out((/ 2, 3, 4, 1, 1 /),'vitu')
385  type(ctrl_out),save :: o_vitv         = ctrl_out((/ 2, 3, 4, 1, 1 /),'vitv')
386  type(ctrl_out),save :: o_vitw         = ctrl_out((/ 2, 3, 10, 10, 1 /),'vitw')
387  type(ctrl_out),save :: o_pres         = ctrl_out((/ 2, 3, 10, 1, 1 /),'pres')
388  type(ctrl_out),save :: o_rneb         = ctrl_out((/ 2, 5, 10, 10, 1 /),'rneb')
389  type(ctrl_out),save :: o_rnebcon      = ctrl_out((/ 2, 5, 10, 10, 1 /),'rnebcon')
390  type(ctrl_out),save :: o_rhum         = ctrl_out((/ 2, 10, 10, 10, 10 /),'rhum')
391  type(ctrl_out),save :: o_ozone        = ctrl_out((/ 2, 10, 10, 10, 10 /),'ozone')
392  type(ctrl_out),save :: o_ozone_light        = ctrl_out((/ 2, 10, 10, 10, 10 /),'ozone_daylight')
393  type(ctrl_out),save :: o_upwd         = ctrl_out((/ 2, 10, 10, 10, 10 /),'upwd')
394  type(ctrl_out),save :: o_dtphy        = ctrl_out((/ 2, 10, 10, 10, 1 /),'dtphy')
395  type(ctrl_out),save :: o_dqphy        = ctrl_out((/ 2, 10, 10, 10, 1 /),'dqphy')
396  type(ctrl_out),save :: o_pr_con_l     = ctrl_out((/ 2, 10, 10, 10, 10 /),'pr_con_l')
397  type(ctrl_out),save :: o_pr_con_i     = ctrl_out((/ 2, 10, 10, 10, 10 /),'pr_con_i')
398  type(ctrl_out),save :: o_pr_lsc_l     = ctrl_out((/ 2, 10, 10, 10, 10 /),'pr_lsc_l')
399  type(ctrl_out),save :: o_pr_lsc_i     = ctrl_out((/ 2, 10, 10, 10, 10 /),'pr_lsc_i')
400!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
401
402  type(ctrl_out),save,dimension(4) :: o_albe_srf     = (/ ctrl_out((/ 3, 4, 10, 1, 10 /),'albe_ter'), &
403                                                     ctrl_out((/ 3, 4, 10, 1, 10 /),'albe_lic'), &
404                                                     ctrl_out((/ 3, 4, 10, 1, 10 /),'albe_oce'), &
405                                                     ctrl_out((/ 3, 4, 10, 1, 10 /),'albe_sic') /)
406
407  type(ctrl_out),save,dimension(4) :: o_ages_srf     = (/ ctrl_out((/ 3, 10, 10, 10, 10 /),'ages_ter'), &
408                                                     ctrl_out((/ 3, 10, 10, 10, 10 /),'ages_lic'), &
409                                                     ctrl_out((/ 3, 10, 10, 10, 10 /),'ages_oce'), &
410                                                     ctrl_out((/ 3, 10, 10, 10, 10 /),'ages_sic') /)
411
412  type(ctrl_out),save,dimension(4) :: o_rugs_srf     = (/ ctrl_out((/ 3, 4, 10, 1, 10 /),'rugs_ter'), &
413                                                     ctrl_out((/ 3, 4, 10, 1, 10 /),'rugs_lic'), &
414                                                     ctrl_out((/ 3, 4, 10, 1, 10 /),'rugs_oce'), &
415                                                     ctrl_out((/ 3, 4, 10, 1, 10 /),'rugs_sic') /)
416
417  type(ctrl_out),save :: o_albs         = ctrl_out((/ 3, 10, 10, 1, 10 /),'albs')
418  type(ctrl_out),save :: o_albslw       = ctrl_out((/ 3, 10, 10, 1, 10 /),'albslw')
419
420  type(ctrl_out),save :: o_clwcon       = ctrl_out((/ 4, 10, 10, 10, 10 /),'clwcon')
421  type(ctrl_out),save :: o_Ma           = ctrl_out((/ 4, 10, 10, 10, 10 /),'Ma')
422  type(ctrl_out),save :: o_dnwd         = ctrl_out((/ 4, 10, 10, 10, 10 /),'dnwd')
423  type(ctrl_out),save :: o_dnwd0        = ctrl_out((/ 4, 10, 10, 10, 10 /),'dnwd0')
424  type(ctrl_out),save :: o_dtdyn        = ctrl_out((/ 4, 10, 10, 10, 1 /),'dtdyn')
425  type(ctrl_out),save :: o_dqdyn        = ctrl_out((/ 4, 10, 10, 10, 1 /),'dqdyn')
426  type(ctrl_out),save :: o_dudyn        = ctrl_out((/ 4, 10, 10, 10, 1 /),'dudyn')  !AXC
427  type(ctrl_out),save :: o_dvdyn        = ctrl_out((/ 4, 10, 10, 10, 1 /),'dvdyn')  !AXC
428  type(ctrl_out),save :: o_dtcon        = ctrl_out((/ 4, 5, 10, 10, 10 /),'dtcon')
429  type(ctrl_out),save :: o_ducon        = ctrl_out((/ 4, 10, 10, 10, 10 /),'ducon')
430  type(ctrl_out),save :: o_dqcon        = ctrl_out((/ 4, 5, 10, 10, 10 /),'dqcon')
431  type(ctrl_out),save :: o_dtwak        = ctrl_out((/ 4, 5, 10, 10, 10 /),'dtwak')
432  type(ctrl_out),save :: o_dqwak        = ctrl_out((/ 4, 5, 10, 10, 10 /),'dqwak')
433  type(ctrl_out),save :: o_wake_h       = ctrl_out((/ 4, 5, 10, 10, 10 /),'wake_h')
434  type(ctrl_out),save :: o_wake_s       = ctrl_out((/ 4, 5, 10, 10, 10 /),'wake_s')
435  type(ctrl_out),save :: o_wake_deltat  = ctrl_out((/ 4, 5, 10, 10, 10 /),'wake_deltat')
436  type(ctrl_out),save :: o_wake_deltaq  = ctrl_out((/ 4, 5, 10, 10, 10 /),'wake_deltaq')
437  type(ctrl_out),save :: o_wake_omg     = ctrl_out((/ 4, 5, 10, 10, 10 /),'wake_omg')
438  type(ctrl_out),save :: o_Vprecip      = ctrl_out((/ 10, 10, 10, 10, 10 /),'Vprecip')
439  type(ctrl_out),save :: o_ftd          = ctrl_out((/ 4, 5, 10, 10, 10 /),'ftd')
440  type(ctrl_out),save :: o_fqd          = ctrl_out((/ 4, 5, 10, 10, 10 /),'fqd')
441  type(ctrl_out),save :: o_dtlsc        = ctrl_out((/ 4, 10, 10, 10, 10 /),'dtlsc')
442  type(ctrl_out),save :: o_dtlschr      = ctrl_out((/ 4, 10, 10, 10, 10 /),'dtlschr')
443  type(ctrl_out),save :: o_dqlsc        = ctrl_out((/ 4, 10, 10, 10, 10 /),'dqlsc')
444  type(ctrl_out),save :: o_dtvdf        = ctrl_out((/ 4, 10, 10, 1, 10 /),'dtvdf')
445  type(ctrl_out),save :: o_dqvdf        = ctrl_out((/ 4, 10, 10, 1, 10 /),'dqvdf')
446  type(ctrl_out),save :: o_dteva        = ctrl_out((/ 4, 10, 10, 10, 10 /),'dteva')
447  type(ctrl_out),save :: o_dqeva        = ctrl_out((/ 4, 10, 10, 10, 10 /),'dqeva')
448  type(ctrl_out),save :: o_ptconv       = ctrl_out((/ 4, 10, 10, 10, 10 /),'ptconv')
449  type(ctrl_out),save :: o_ratqs        = ctrl_out((/ 4, 10, 10, 10, 10 /),'ratqs')
450  type(ctrl_out),save :: o_dtthe        = ctrl_out((/ 4, 10, 10, 10, 10 /),'dtthe')
451  type(ctrl_out),save :: o_f_th         = ctrl_out((/ 4, 10, 10, 10, 10 /),'f_th')
452  type(ctrl_out),save :: o_e_th         = ctrl_out((/ 4, 10, 10, 10, 10 /),'e_th')
453  type(ctrl_out),save :: o_w_th         = ctrl_out((/ 4, 10, 10, 10, 10 /),'w_th')
454  type(ctrl_out),save :: o_lambda_th    = ctrl_out((/ 10, 10, 10, 10, 10 /),'lambda_th')
455  type(ctrl_out),save :: o_q_th         = ctrl_out((/ 4, 10, 10, 10, 10 /),'q_th')
456  type(ctrl_out),save :: o_a_th         = ctrl_out((/ 4, 10, 10, 10, 10 /),'a_th')
457  type(ctrl_out),save :: o_d_th         = ctrl_out((/ 4, 10, 10, 10, 10 /),'d_th')
458  type(ctrl_out),save :: o_f0_th        = ctrl_out((/ 4, 10, 10, 10, 10 /),'f0_th')
459  type(ctrl_out),save :: o_zmax_th      = ctrl_out((/ 4, 10, 10, 10, 10 /),'zmax_th')
460  type(ctrl_out),save :: o_dqthe        = ctrl_out((/ 4, 10, 10, 10, 1 /),'dqthe')
461  type(ctrl_out),save :: o_dtajs        = ctrl_out((/ 4, 10, 10, 10, 10 /),'dtajs')
462  type(ctrl_out),save :: o_dqajs        = ctrl_out((/ 4, 10, 10, 10, 10 /),'dqajs')
463  type(ctrl_out),save :: o_dtswr        = ctrl_out((/ 4, 10, 10, 10, 1 /),'dtswr')
464  type(ctrl_out),save :: o_dtsw0        = ctrl_out((/ 4, 10, 10, 10, 10 /),'dtsw0')
465  type(ctrl_out),save :: o_dtlwr        = ctrl_out((/ 4, 10, 10, 10, 1 /),'dtlwr')
466  type(ctrl_out),save :: o_dtlw0        = ctrl_out((/ 4, 10, 10, 10, 10 /),'dtlw0')
467  type(ctrl_out),save :: o_dtec         = ctrl_out((/ 4, 10, 10, 10, 10 /),'dtec')
468  type(ctrl_out),save :: o_duvdf        = ctrl_out((/ 4, 10, 10, 10, 10 /),'duvdf')
469  type(ctrl_out),save :: o_dvvdf        = ctrl_out((/ 4, 10, 10, 10, 10 /),'dvvdf')
470  type(ctrl_out),save :: o_duoro        = ctrl_out((/ 4, 10, 10, 10, 10 /),'duoro')
471  type(ctrl_out),save :: o_dvoro        = ctrl_out((/ 4, 10, 10, 10, 10 /),'dvoro')
472  type(ctrl_out),save :: o_dulif        = ctrl_out((/ 4, 10, 10, 10, 10 /),'dulif')
473  type(ctrl_out),save :: o_dvlif        = ctrl_out((/ 4, 10, 10, 10, 10 /),'dvlif')
474
475! Attention a refaire correctement
476  type(ctrl_out),save,dimension(2) :: o_trac         = (/ ctrl_out((/ 4, 10, 10, 10, 10 /),'trac01'), &
477                                                     ctrl_out((/ 4, 10, 10, 10, 10 /),'trac02') /)
478    CONTAINS
479
480!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
481!!!!!!!!! Ouverture des fichier et definition des variable de sortie !!!!!!!!
482!! histbeg, histvert et histdef
483!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
484 
485  SUBROUTINE phys_output_open(jjmp1,nlevSTD,clevSTD,nbteta, &
486       ctetaSTD,dtime, ok_veget, &
487       type_ocean, iflag_pbl,ok_mensuel,ok_journe, &
488       ok_hf,ok_instan,ok_LES,ok_ade,ok_aie, read_climoz, &
489       new_aod, aerosol_couple)   
490
491
492  USE iophy
493  USE dimphy
494  USE infotrac
495  USE ioipsl
496  USE mod_phys_lmdz_para
497  USE aero_mod, only : naero_spc,name_aero
498
499  IMPLICIT NONE
500  include "dimensions.h"
501  include "temps.h"
502  include "indicesol.h"
503  include "clesphys.h"
504  include "thermcell.h"
505  include "comvert.h"
506
507  integer                               :: jjmp1
508  integer                               :: nbteta, nlevSTD, radpas
509  logical                               :: ok_mensuel, ok_journe, ok_hf, ok_instan
510  logical                               :: ok_LES,ok_ade,ok_aie
511  logical                               :: new_aod, aerosol_couple
512  integer, intent(in)::  read_climoz ! read ozone climatology
513  !     Allowed values are 0, 1 and 2
514  !     0: do not read an ozone climatology
515  !     1: read a single ozone climatology that will be used day and night
516  !     2: read two ozone climatologies, the average day and night
517  !     climatology and the daylight climatology
518
519  real                                  :: dtime
520  integer                               :: idayref
521  real                                  :: zjulian
522  real, dimension(klev)                 :: Ahyb, Bhyb, Alt
523  character(len=4), dimension(nlevSTD)  :: clevSTD
524  integer                               :: nsrf, k, iq, iiq, iff, i, j, ilev
525  integer                               :: naero
526  logical                               :: ok_veget
527  integer                               :: iflag_pbl
528  CHARACTER(len=4)                      :: bb2
529  CHARACTER(len=2)                      :: bb3
530  character(len=6)                      :: type_ocean
531  CHARACTER(len=3)                      :: ctetaSTD(nbteta)
532  real, dimension(nfiles)               :: ecrit_files
533  CHARACTER(len=20), dimension(nfiles)  :: phys_out_filenames
534  INTEGER, dimension(iim*jjmp1)         ::  ndex2d
535  INTEGER, dimension(iim*jjmp1*klev)    :: ndex3d
536  integer                               :: imin_ins, imax_ins
537  integer                               :: jmin_ins, jmax_ins
538  integer, dimension(nfiles)            :: phys_out_levmin, phys_out_levmax
539  integer, dimension(nfiles)            :: phys_out_filelevels
540  CHARACTER(len=20), dimension(nfiles)  :: type_ecri_files, phys_out_filetypes
541  character(len=20), dimension(nfiles)  :: chtimestep   = (/ 'DefFreq', 'DefFreq','DefFreq', 'DefFreq', 'DefFreq' /)
542  logical, dimension(nfiles)            :: phys_out_filekeys
543
544!!!!!!!!!! stockage dans une region limitee pour chaque fichier !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
545!                 entre [phys_out_lonmin,phys_out_lonmax] et [phys_out_latmin,phys_out_latmax]
546
547  logical, dimension(nfiles), save  :: phys_out_regfkey       = (/ .false., .false., .false., .false., .false. /)
548  real, dimension(nfiles), save     :: phys_out_lonmin        = (/ -180., -180., -180., -180., -180. /)
549  real, dimension(nfiles), save     :: phys_out_lonmax        = (/ 180., 180., 180., 180., 180. /)
550  real, dimension(nfiles), save     :: phys_out_latmin        = (/ -90., -90., -90., -90., -90. /)
551  real, dimension(nfiles), save     :: phys_out_latmax        = (/ 90., 90., 90., 90., 90. /)
552 
553 
554
555!
556   print*,'Debut phys_output_mod.F90'
557! Initialisations (Valeurs par defaut
558   levmax = (/ klev, klev, klev, klev, klev /)
559
560   phys_out_filenames(1) = 'histmth'
561   phys_out_filenames(2) = 'histday'
562   phys_out_filenames(3) = 'histhf'
563   phys_out_filenames(4) = 'histins'
564   phys_out_filenames(5) = 'histLES'
565
566   type_ecri(1) = 'ave(X)'
567   type_ecri(2) = 'ave(X)'
568   type_ecri(3) = 'ave(X)'
569   type_ecri(4) = 'inst(X)'
570   type_ecri(5) = 'inst(X)'
571
572   clef_files(1) = ok_mensuel
573   clef_files(2) = ok_journe
574   clef_files(3) = ok_hf
575   clef_files(4) = ok_instan
576   clef_files(5) = ok_LES
577
578   lev_files(1) = lev_histmth
579   lev_files(2) = lev_histday
580   lev_files(3) = lev_histhf
581   lev_files(4) = lev_histins
582   lev_files(5) = lev_histLES
583
584
585   ecrit_files(1) = ecrit_mth
586   ecrit_files(2) = ecrit_day
587   ecrit_files(3) = ecrit_hf
588   ecrit_files(4) = ecrit_ins
589   ecrit_files(5) = ecrit_LES
590 
591!! Lectures des parametres de sorties dans physiq.def
592
593   call getin('phys_out_regfkey',phys_out_regfkey)
594   call getin('phys_out_lonmin',phys_out_lonmin)
595   call getin('phys_out_lonmax',phys_out_lonmax)
596   call getin('phys_out_latmin',phys_out_latmin)
597   call getin('phys_out_latmax',phys_out_latmax)
598     phys_out_levmin(:)=levmin(:)
599   call getin('phys_out_levmin',levmin)
600     phys_out_levmax(:)=levmax(:)
601   call getin('phys_out_levmax',levmax)
602   call getin('phys_out_filenames',phys_out_filenames)
603     phys_out_filekeys(:)=clef_files(:)
604   call getin('phys_out_filekeys',clef_files)
605     phys_out_filelevels(:)=lev_files(:)
606   call getin('phys_out_filelevels',lev_files)
607   call getin('phys_out_filetimesteps',chtimestep)
608     phys_out_filetypes(:)=type_ecri(:)
609   call getin('phys_out_filetypes',type_ecri)
610
611   type_ecri_files(:)=type_ecri(:)
612
613   print*,'phys_out_lonmin=',phys_out_lonmin
614   print*,'phys_out_lonmax=',phys_out_lonmax
615   print*,'phys_out_latmin=',phys_out_latmin
616   print*,'phys_out_latmax=',phys_out_latmax
617   print*,'phys_out_filenames=',phys_out_filenames
618   print*,'phys_out_filetypes=',type_ecri
619   print*,'phys_out_filekeys=',clef_files
620   print*,'phys_out_filelevels=',lev_files
621
622!!!!!!!!!!!!!!!!!!!!!!! Boucle sur les fichiers !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
623! Appel de histbeg et histvert pour creer le fichier et les niveaux verticaux !!
624! Appel des histbeg pour definir les variables (nom, moy ou inst, freq de sortie ..
625!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
626
627 zdtime = dtime         ! Frequence ou l on moyenne
628
629! Calcul des Ahyb, Bhyb et Alt
630         do k=1,klev
631          Ahyb(k)=(ap(k)+ap(k+1))/2.
632          Bhyb(k)=(bp(k)+bp(k+1))/2.
633          Alt(k)=log(preff/presnivs(k))*8.
634         enddo
635!          if(prt_level.ge.1) then
636           print*,'Ap Hybrid = ',Ahyb(1:klev)
637           print*,'Bp Hybrid = ',Bhyb(1:klev)
638           print*,'Alt approx des couches pour une haut d echelle de 8km = ',Alt(1:klev)
639!          endif
640 DO iff=1,nfiles
641
642    IF (clef_files(iff)) THEN
643
644      if ( chtimestep(iff).eq.'DefFreq' ) then
645! Par defaut ecrit_files = (ecrit_mensuel ecrit_jour ecrit_hf ...)*86400.
646        ecrit_files(iff)=ecrit_files(iff)*86400.
647      else
648        call convers_timesteps(chtimestep(iff),ecrit_files(iff))
649      endif
650       print*,'ecrit_files(',iff,')= ',ecrit_files(iff)
651
652      zoutm(iff) = ecrit_files(iff) ! Frequence ou l on ecrit en seconde
653
654      idayref = day_ref
655      CALL ymds2ju(annee_ref, 1, idayref, 0.0, zjulian)
656
657!!!!!!!!!!!!!!!!! Traitement dans le cas ou l'on veut stocker sur un domaine limite !!
658!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
659     if (phys_out_regfkey(iff)) then
660
661        imin_ins=1
662        imax_ins=iim
663        jmin_ins=1
664        jmax_ins=jjmp1
665
666! correction abderr       
667        do i=1,iim
668           print*,'io_lon(i)=',io_lon(i)
669           if (io_lon(i).le.phys_out_lonmin(iff)) imin_ins=i
670           if (io_lon(i).le.phys_out_lonmax(iff)) imax_ins=i+1
671        enddo
672
673        do j=1,jjmp1
674            print*,'io_lat(j)=',io_lat(j)
675            if (io_lat(j).ge.phys_out_latmin(iff)) jmax_ins=j+1
676            if (io_lat(j).ge.phys_out_latmax(iff)) jmin_ins=j
677        enddo
678
679        print*,'On stoke le fichier histoire numero ',iff,' sur ', &
680         imin_ins,imax_ins,jmin_ins,jmax_ins
681         print*,'longitudes : ', &
682         io_lon(imin_ins),io_lon(imax_ins), &
683         'latitudes : ', &
684         io_lat(jmax_ins),io_lat(jmin_ins)
685
686 CALL histbeg(phys_out_filenames(iff),iim,io_lon,jjmp1,io_lat, &
687              imin_ins,imax_ins-imin_ins+1, &
688              jmin_ins,jmax_ins-jmin_ins+1, &
689              itau_phy,zjulian,dtime,nhorim(iff),nid_files(iff))
690!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
691       else
692 CALL histbeg_phy(phys_out_filenames(iff),itau_phy,zjulian,dtime,nhorim(iff),nid_files(iff))
693       endif
694 
695      CALL histvert(nid_files(iff), "presnivs", "Vertical levels", "mb", &
696           levmax(iff) - levmin(iff) + 1, &
697           presnivs(levmin(iff):levmax(iff))/100., nvertm(iff))
698
699!!!!!!!!!!!!! Traitement des champs 3D pour histhf !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
700!!!!!!!!!!!!!!! A Revoir plus tard !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
701!          IF (iff.eq.3.and.lev_files(iff).ge.4) THEN
702!          CALL histbeg_phy("histhf3d",itau_phy, &
703!     &                     zjulian, dtime, &
704!     &                     nhorim, nid_hf3d)
705
706!         CALL histvert(nid_hf3d, "presnivs", &
707!     &                 "Vertical levels", "mb", &
708!     &                 klev, presnivs/100., nvertm)
709!          ENDIF
710!
711!!!! Composentes de la coordonnee sigma-hybride
712   CALL histvert(nid_files(iff), "Ahyb","Ahyb comp of Hyb Cord ", "Pa", &
713                 levmax(iff) - levmin(iff) + 1,Ahyb,nvertap(iff))
714
715   CALL histvert(nid_files(iff), "Bhyb","Bhyb comp of Hyb Cord", " ", &
716                 levmax(iff) - levmin(iff) + 1,Bhyb,nvertbp(iff))
717
718   CALL histvert(nid_files(iff), "Alt","Height approx for scale heigh of 8km at levels", "Km", &
719                 levmax(iff) - levmin(iff) + 1,Alt,nvertAlt(iff))
720
721!   CALL histvert(nid_files(iff), "preff","Reference pressure", "Pa", &
722!                 1,preff,nvertp0(iff))
723!!! Champs 1D !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
724 CALL histdef2d(iff,o_phis%flag,o_phis%name,"Surface geop.height", "m2/s2")
725   type_ecri(1) = 'once'
726   type_ecri(2) = 'once'
727   type_ecri(3) = 'once'
728   type_ecri(4) = 'once'
729   type_ecri(5) = 'once'
730 CALL histdef2d(iff,o_aire%flag,o_aire%name,"Grid area", "-")
731 CALL histdef2d(iff,o_contfracATM%flag,o_contfracATM%name,"% sfce ter+lic", "-")
732   type_ecri(:) = type_ecri_files(:)
733
734!!! Champs 2D !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
735 CALL histdef2d(iff,o_contfracOR%flag,o_contfracOR%name,"% sfce terre OR", "-" )
736 CALL histdef2d(iff,o_aireTER%flag,o_aireTER%name,"Grid area CONT", "-" )
737 CALL histdef2d(iff,o_flat%flag,o_flat%name, "Latent heat flux", "W/m2")
738 CALL histdef2d(iff,o_slp%flag,o_slp%name, "Sea Level Pressure", "Pa" )
739 CALL histdef2d(iff,o_tsol%flag,o_tsol%name, "Surface Temperature", "K")
740 CALL histdef2d(iff,o_t2m%flag,o_t2m%name, "Temperature 2m", "K" )
741   type_ecri(1) = 't_min(X)'
742   type_ecri(2) = 't_min(X)'
743   type_ecri(3) = 't_min(X)'
744   type_ecri(4) = 't_min(X)'
745   type_ecri(5) = 't_min(X)'
746 CALL histdef2d(iff,o_t2m_min%flag,o_t2m_min%name, "Temp 2m min", "K" )
747   type_ecri(1) = 't_max(X)'
748   type_ecri(2) = 't_max(X)'
749   type_ecri(3) = 't_max(X)'
750   type_ecri(4) = 't_max(X)'
751   type_ecri(5) = 't_max(X)'
752 CALL histdef2d(iff,o_t2m_max%flag,o_t2m_max%name, "Temp 2m max", "K" )
753   type_ecri(:) = type_ecri_files(:)
754 CALL histdef2d(iff,o_wind10m%flag,o_wind10m%name, "10-m wind speed", "m/s")
755 CALL histdef2d(iff,o_wind10max%flag,o_wind10max%name, "10m wind speed max", "m/s")
756 CALL histdef2d(iff,o_sicf%flag,o_sicf%name, "Sea-ice fraction", "-" )
757 CALL histdef2d(iff,o_q2m%flag,o_q2m%name, "Specific humidity 2m", "kg/kg")
758 CALL histdef2d(iff,o_u10m%flag,o_u10m%name, "Vent zonal 10m", "m/s" )
759 CALL histdef2d(iff,o_v10m%flag,o_v10m%name, "Vent meridien 10m", "m/s")
760 CALL histdef2d(iff,o_psol%flag,o_psol%name, "Surface Pressure", "Pa" )
761 CALL histdef2d(iff,o_qsurf%flag,o_qsurf%name, "Surface Air humidity", "kg/kg")
762
763  if (.not. ok_veget) then
764 CALL histdef2d(iff,o_qsol%flag,o_qsol%name, "Soil watter content", "mm" )
765  endif
766
767 CALL histdef2d(iff,o_ndayrain%flag,o_ndayrain%name, "Number of dayrain(liq+sol)", "-")
768 CALL histdef2d(iff,o_precip%flag,o_precip%name, "Precip Totale liq+sol", "kg/(s*m2)" )
769 CALL histdef2d(iff,o_plul%flag,o_plul%name, "Large-scale Precip.", "kg/(s*m2)")
770 CALL histdef2d(iff,o_pluc%flag,o_pluc%name, "Convective Precip.", "kg/(s*m2)")
771 CALL histdef2d(iff,o_snow%flag,o_snow%name, "Snow fall", "kg/(s*m2)" )
772 CALL histdef2d(iff,o_evap%flag,o_evap%name, "Evaporat", "kg/(s*m2)" )
773 CALL histdef2d(iff,o_tops%flag,o_tops%name, "Solar rad. at TOA", "W/m2")
774 CALL histdef2d(iff,o_tops0%flag,o_tops0%name, "CS Solar rad. at TOA", "W/m2")
775 CALL histdef2d(iff,o_topl%flag,o_topl%name, "IR rad. at TOA", "W/m2" )
776 CALL histdef2d(iff,o_topl0%flag,o_topl0%name, "IR rad. at TOA", "W/m2")
777 CALL histdef2d(iff,o_SWupTOA%flag,o_SWupTOA%name, "SWup at TOA", "W/m2")
778 CALL histdef2d(iff,o_SWupTOAclr%flag,o_SWupTOAclr%name, "SWup clear sky at TOA", "W/m2")
779 CALL histdef2d(iff,o_SWdnTOA%flag,o_SWdnTOA%name, "SWdn at TOA", "W/m2" )
780 CALL histdef2d(iff,o_SWdnTOAclr%flag,o_SWdnTOAclr%name, "SWdn clear sky at TOA", "W/m2")
781 CALL histdef2d(iff,o_SWup200%flag,o_SWup200%name, "SWup at 200mb", "W/m2" )
782 CALL histdef2d(iff,o_SWup200clr%flag,o_SWup200clr%name, "SWup clear sky at 200mb", "W/m2")
783 CALL histdef2d(iff,o_SWdn200%flag,o_SWdn200%name, "SWdn at 200mb", "W/m2" )
784 CALL histdef2d(iff,o_SWdn200clr%flag,o_SWdn200clr%name, "SWdn clear sky at 200mb", "W/m2")
785 CALL histdef2d(iff,o_LWup200%flag,o_LWup200%name, "LWup at 200mb", "W/m2")
786 CALL histdef2d(iff,o_LWup200clr%flag,o_LWup200clr%name, "LWup clear sky at 200mb", "W/m2")
787 CALL histdef2d(iff,o_LWdn200%flag,o_LWdn200%name, "LWdn at 200mb", "W/m2")
788 CALL histdef2d(iff,o_LWdn200clr%flag,o_LWdn200clr%name, "LWdn clear sky at 200mb", "W/m2")
789 CALL histdef2d(iff,o_sols%flag,o_sols%name, "Solar rad. at surf.", "W/m2")
790 CALL histdef2d(iff,o_sols0%flag,o_sols0%name, "Solar rad. at surf.", "W/m2")
791 CALL histdef2d(iff,o_soll%flag,o_soll%name, "IR rad. at surface", "W/m2") 
792 CALL histdef2d(iff,o_radsol%flag,o_radsol%name, "Rayonnement au sol", "W/m2")
793 CALL histdef2d(iff,o_soll0%flag,o_soll0%name, "IR rad. at surface", "W/m2")
794 CALL histdef2d(iff,o_SWupSFC%flag,o_SWupSFC%name, "SWup at surface", "W/m2")
795 CALL histdef2d(iff,o_SWupSFCclr%flag,o_SWupSFCclr%name, "SWup clear sky at surface", "W/m2")
796 CALL histdef2d(iff,o_SWdnSFC%flag,o_SWdnSFC%name, "SWdn at surface", "W/m2")
797 CALL histdef2d(iff,o_SWdnSFCclr%flag,o_SWdnSFCclr%name, "SWdn clear sky at surface", "W/m2")
798 CALL histdef2d(iff,o_LWupSFC%flag,o_LWupSFC%name, "Upwd. IR rad. at surface", "W/m2")
799 CALL histdef2d(iff,o_LWdnSFC%flag,o_LWdnSFC%name, "Down. IR rad. at surface", "W/m2")
800 CALL histdef2d(iff,o_LWupSFCclr%flag,o_LWupSFCclr%name, "CS Upwd. IR rad. at surface", "W/m2")
801 CALL histdef2d(iff,o_LWdnSFCclr%flag,o_LWdnSFCclr%name, "Down. CS IR rad. at surface", "W/m2")
802 CALL histdef2d(iff,o_bils%flag,o_bils%name, "Surf. total heat flux", "W/m2")
803 CALL histdef2d(iff,o_sens%flag,o_sens%name, "Sensible heat flux", "W/m2")
804 CALL histdef2d(iff,o_fder%flag,o_fder%name, "Heat flux derivation", "W/m2")
805 CALL histdef2d(iff,o_ffonte%flag,o_ffonte%name, "Thermal flux for snow melting", "W/m2")
806 CALL histdef2d(iff,o_fqcalving%flag,o_fqcalving%name, "Ice Calving", "kg/m2/s")
807 CALL histdef2d(iff,o_fqfonte%flag,o_fqfonte%name, "Land ice melt", "kg/m2/s")
808
809     DO nsrf = 1, nbsrf
810 CALL histdef2d(iff,o_pourc_srf(nsrf)%flag,o_pourc_srf(nsrf)%name,"% "//clnsurf(nsrf),"%")
811 CALL histdef2d(iff,o_fract_srf(nsrf)%flag,o_fract_srf(nsrf)%name,"Fraction "//clnsurf(nsrf),"1")
812 CALL histdef2d(iff,o_taux_srf(nsrf)%flag,o_taux_srf(nsrf)%name,"Zonal wind stress"//clnsurf(nsrf),"Pa")
813 CALL histdef2d(iff,o_tauy_srf(nsrf)%flag,o_tauy_srf(nsrf)%name,"Meridional wind stress "//clnsurf(nsrf),"Pa")
814 CALL histdef2d(iff,o_tsol_srf(nsrf)%flag,o_tsol_srf(nsrf)%name,"Temperature "//clnsurf(nsrf),"K")
815 CALL histdef2d(iff,o_u10m_srf(nsrf)%flag,o_u10m_srf(nsrf)%name,"Vent Zonal 10m "//clnsurf(nsrf),"m/s")
816 CALL histdef2d(iff,o_v10m_srf(nsrf)%flag,o_v10m_srf(nsrf)%name,"Vent meredien 10m "//clnsurf(nsrf),"m/s")
817 CALL histdef2d(iff,o_t2m_srf(nsrf)%flag,o_t2m_srf(nsrf)%name,"Temp 2m "//clnsurf(nsrf),"K")
818 CALL histdef2d(iff,o_sens_srf(nsrf)%flag,o_sens_srf(nsrf)%name,"Sensible heat flux "//clnsurf(nsrf),"W/m2")
819 CALL histdef2d(iff,o_lat_srf(nsrf)%flag,o_lat_srf(nsrf)%name,"Latent heat flux "//clnsurf(nsrf),"W/m2")
820 CALL histdef2d(iff,o_flw_srf(nsrf)%flag,o_flw_srf(nsrf)%name,"LW "//clnsurf(nsrf),"W/m2")
821 CALL histdef2d(iff,o_fsw_srf(nsrf)%flag,o_fsw_srf(nsrf)%name,"SW "//clnsurf(nsrf),"W/m2")
822 CALL histdef2d(iff,o_wbils_srf(nsrf)%flag,o_wbils_srf(nsrf)%name,"Bilan sol "//clnsurf(nsrf),"W/m2" )
823 CALL histdef2d(iff,o_wbilo_srf(nsrf)%flag,o_wbilo_srf(nsrf)%name,"Bilan eau "//clnsurf(nsrf),"kg/(m2*s)")
824  if (iflag_pbl>1 .and. lev_files(iff).gt.10 ) then
825 CALL histdef2d(iff,o_tke_srf(nsrf)%flag,o_tke_srf(nsrf)%name,"Max Turb. Kinetic Energy "//clnsurf(nsrf),"-")
826   type_ecri(1) = 't_max(X)'
827   type_ecri(2) = 't_max(X)'
828   type_ecri(3) = 't_max(X)'
829   type_ecri(4) = 't_max(X)'
830   type_ecri(5) = 't_max(X)'
831 CALL histdef2d(iff,o_tke_max_srf(nsrf)%flag,o_tke_max_srf(nsrf)%name,"Max Turb. Kinetic Energy "//clnsurf(nsrf),"-")
832   type_ecri(:) = type_ecri_files(:)
833  endif
834 CALL histdef2d(iff,o_albe_srf(nsrf)%flag,o_albe_srf(nsrf)%name,"Albedo surf. "//clnsurf(nsrf),"-")
835 CALL histdef2d(iff,o_rugs_srf(nsrf)%flag,o_rugs_srf(nsrf)%name,"Latent heat flux "//clnsurf(nsrf),"W/m2")
836 CALL histdef2d(iff,o_ages_srf(nsrf)%flag,o_ages_srf(nsrf)%name,"Snow age", "day")
837END DO
838
839IF (new_aod .AND. (.NOT. aerosol_couple)) THEN
840  DO naero = 1, naero_spc
841  CALL histdef2d(iff,o_tausumaero(naero)%flag,o_tausumaero(naero)%name,"Aerosol Optical depth at 550 nm "//name_aero(naero),"1")
842  END DO
843ENDIF
844
845 IF (ok_ade) THEN
846  CALL histdef2d(iff,o_topswad%flag,o_topswad%name, "ADE at TOA", "W/m2")
847  CALL histdef2d(iff,o_solswad%flag,o_solswad%name, "ADE at SRF", "W/m2")
848
849 CALL histdef2d(iff,o_swtoaas_nat%flag,o_swtoaas_nat%name, "Natural aerosol radiative forcing all-sky at TOA", "W/m2")
850 CALL histdef2d(iff,o_swsrfas_nat%flag,o_swsrfas_nat%name, "Natural aerosol radiative forcing all-sky at SRF", "W/m2")
851 CALL histdef2d(iff,o_swtoacs_nat%flag,o_swtoacs_nat%name, "Natural aerosol radiative forcing clear-sky at TOA", "W/m2")
852 CALL histdef2d(iff,o_swsrfcs_nat%flag,o_swsrfcs_nat%name, "Natural aerosol radiative forcing clear-sky at SRF", "W/m2")
853
854 CALL histdef2d(iff,o_swtoaas_ant%flag,o_swtoaas_ant%name, "Anthropogenic aerosol radiative forcing all-sky at TOA", "W/m2")
855 CALL histdef2d(iff,o_swsrfas_ant%flag,o_swsrfas_ant%name, "Anthropogenic aerosol radiative forcing all-sky at SRF", "W/m2")
856 CALL histdef2d(iff,o_swtoacs_ant%flag,o_swtoacs_ant%name, "Anthropogenic aerosol radiative forcing clear-sky at TOA", "W/m2")
857 CALL histdef2d(iff,o_swsrfcs_ant%flag,o_swsrfcs_ant%name, "Anthropogenic aerosol radiative forcing clear-sky at SRF", "W/m2")
858
859 IF (.NOT. aerosol_couple) THEN
860 CALL histdef2d(iff,o_swtoacf_nat%flag,o_swtoacf_nat%name, "Natural aerosol impact on cloud radiative forcing at TOA", "W/m2")
861 CALL histdef2d(iff,o_swsrfcf_nat%flag,o_swsrfcf_nat%name, "Natural aerosol impact on cloud radiative forcing  at SRF", "W/m2")
862 CALL histdef2d(iff,o_swtoacf_ant%flag,o_swtoacf_ant%name, "Anthropogenic aerosol impact on cloud radiative forcing at TOA", "W/m2")
863 CALL histdef2d(iff,o_swsrfcf_ant%flag,o_swsrfcf_ant%name, "Anthropogenic aerosol impact on cloud radiative forcing at SRF", "W/m2")
864 CALL histdef2d(iff,o_swtoacf_zero%flag,o_swtoacf_zero%name, "Cloud radiative forcing (allsky-clearsky fluxes) at TOA", "W/m2")
865 CALL histdef2d(iff,o_swsrfcf_zero%flag,o_swsrfcf_zero%name, "Cloud radiative forcing (allsky-clearsky fluxes) at SRF", "W/m2")
866 ENDIF
867
868 ENDIF
869
870 IF (ok_aie) THEN
871  CALL histdef2d(iff,o_topswai%flag,o_topswai%name, "AIE at TOA", "W/m2")
872  CALL histdef2d(iff,o_solswai%flag,o_solswai%name, "AIE at SFR", "W/m2")
873 ENDIF
874
875
876 CALL histdef2d(iff,o_albs%flag,o_albs%name, "Surface albedo", "-")
877 CALL histdef2d(iff,o_albslw%flag,o_albslw%name, "Surface albedo LW", "-")
878 CALL histdef2d(iff,o_cdrm%flag,o_cdrm%name, "Momentum drag coef.", "-")
879 CALL histdef2d(iff,o_cdrh%flag,o_cdrh%name, "Heat drag coef.", "-" )
880 CALL histdef2d(iff,o_cldl%flag,o_cldl%name, "Low-level cloudiness", "-")
881 CALL histdef2d(iff,o_cldm%flag,o_cldm%name, "Mid-level cloudiness", "-")
882 CALL histdef2d(iff,o_cldh%flag,o_cldh%name, "High-level cloudiness", "-")
883 CALL histdef2d(iff,o_cldt%flag,o_cldt%name, "Total cloudiness", "%")
884 CALL histdef2d(iff,o_cldq%flag,o_cldq%name, "Cloud liquid water path", "kg/m2")
885 CALL histdef2d(iff,o_lwp%flag,o_lwp%name, "Cloud water path", "kg/m2")
886 CALL histdef2d(iff,o_iwp%flag,o_iwp%name, "Cloud ice water path", "kg/m2" )
887 CALL histdef2d(iff,o_ue%flag,o_ue%name, "Zonal energy transport", "-")
888 CALL histdef2d(iff,o_ve%flag,o_ve%name, "Merid energy transport", "-")
889 CALL histdef2d(iff,o_uq%flag,o_uq%name, "Zonal humidity transport", "-")
890 CALL histdef2d(iff,o_vq%flag,o_vq%name, "Merid humidity transport", "-")
891
892     IF(iflag_con.GE.3) THEN ! sb
893 CALL histdef2d(iff,o_cape%flag,o_cape%name, "Conv avlbl pot ener", "J/kg")
894 CALL histdef2d(iff,o_pbase%flag,o_pbase%name, "Cld base pressure", "mb")
895 CALL histdef2d(iff,o_ptop%flag,o_ptop%name, "Cld top pressure", "mb")
896 CALL histdef2d(iff,o_fbase%flag,o_fbase%name, "Cld base mass flux", "kg/m2/s")
897 CALL histdef2d(iff,o_prw%flag,o_prw%name, "Precipitable water", "kg/m2")
898   type_ecri(1) = 't_max(X)'
899   type_ecri(2) = 't_max(X)'
900   type_ecri(3) = 't_max(X)'
901   type_ecri(4) = 't_max(X)'
902   type_ecri(5) = 't_max(X)'
903 CALL histdef2d(iff,o_cape_max%flag,o_cape_max%name, "CAPE max.", "J/kg")
904   type_ecri(:) = type_ecri_files(:)
905 CALL histdef3d(iff,o_upwd%flag,o_upwd%name, "saturated updraft", "kg/m2/s")
906 CALL histdef3d(iff,o_Ma%flag,o_Ma%name, "undilute adiab updraft", "kg/m2/s")
907 CALL histdef3d(iff,o_dnwd%flag,o_dnwd%name, "saturated downdraft", "kg/m2/s")
908 CALL histdef3d(iff,o_dnwd0%flag,o_dnwd0%name, "unsat. downdraft", "kg/m2/s")
909     ENDIF !iflag_con .GE. 3
910
911 CALL histdef2d(iff,o_s_pblh%flag,o_s_pblh%name, "Boundary Layer Height", "m")
912 CALL histdef2d(iff,o_s_pblt%flag,o_s_pblt%name, "t at Boundary Layer Height", "K")
913 CALL histdef2d(iff,o_s_lcl%flag,o_s_lcl%name, "Condensation level", "m")
914 CALL histdef2d(iff,o_s_capCL%flag,o_s_capCL%name, "Conv avlbl pot enerfor ABL", "J/m2" )
915 CALL histdef2d(iff,o_s_oliqCL%flag,o_s_oliqCL%name, "Liq Water in BL", "kg/m2")
916 CALL histdef2d(iff,o_s_cteiCL%flag,o_s_cteiCL%name, "Instability criteria(ABL)", "K")
917 CALL histdef2d(iff,o_s_therm%flag,o_s_therm%name, "Exces du thermique", "K")
918 CALL histdef2d(iff,o_s_trmb1%flag,o_s_trmb1%name, "deep_cape(HBTM2)", "J/m2")
919 CALL histdef2d(iff,o_s_trmb2%flag,o_s_trmb2%name, "inhibition (HBTM2)", "J/m2")
920 CALL histdef2d(iff,o_s_trmb3%flag,o_s_trmb3%name, "Point Omega (HBTM2)", "m")
921
922! Champs interpolles sur des niveaux de pression
923
924   type_ecri(1) = 'inst(X)'
925   type_ecri(2) = 'inst(X)'
926   type_ecri(3) = 'inst(X)'
927   type_ecri(4) = 'inst(X)'
928   type_ecri(5) = 'inst(X)'
929
930! Attention a reverifier
931
932        ilev=0       
933        DO k=1, nlevSTD
934!     IF(k.GE.2.AND.k.LE.12) bb2=clevSTD(k)
935     bb2=clevSTD(k)
936     IF(bb2.EQ."850".OR.bb2.EQ."700".OR.bb2.EQ."500".OR.bb2.EQ."200".OR.bb2.EQ."50".OR.bb2.EQ."10")THEN
937      ilev=ilev+1
938      print*,'ilev k bb2 flag name ',ilev,k, bb2,o_uSTDlevs(ilev)%flag,o_uSTDlevs(ilev)%name
939 CALL histdef2d(iff,o_uSTDlevs(ilev)%flag,o_uSTDlevs(ilev)%name,"Zonal wind "//bb2//"mb", "m/s")
940 CALL histdef2d(iff,o_vSTDlevs(ilev)%flag,o_vSTDlevs(ilev)%name,"Meridional wind "//bb2//"mb", "m/s")
941 CALL histdef2d(iff,o_wSTDlevs(ilev)%flag,o_wSTDlevs(ilev)%name,"Vertical wind "//bb2//"mb", "Pa/s")
942 CALL histdef2d(iff,o_phiSTDlevs(ilev)%flag,o_phiSTDlevs(ilev)%name,"Geopotential "//bb2//"mb", "m")
943 CALL histdef2d(iff,o_qSTDlevs(ilev)%flag,o_qSTDlevs(ilev)%name,"Specific humidity "//bb2//"mb", "kg/kg" )
944 CALL histdef2d(iff,o_tSTDlevs(ilev)%flag,o_tSTDlevs(ilev)%name,"Temperature "//bb2//"mb", "K")
945     ENDIF !(bb2.EQ."850".OR.bb2.EQ."700".OR."500".OR.bb2.EQ."200".OR.bb2.EQ."50".OR.bb2.EQ."10")
946       ENDDO
947   type_ecri(:) = type_ecri_files(:)
948
949 CALL histdef2d(iff,o_t_oce_sic%flag,o_t_oce_sic%name, "Temp mixte oce-sic", "K")
950
951 IF (type_ocean=='slab') &
952     CALL histdef2d(iff,o_slab_bils%flag, o_slab_bils%name,"Bilan au sol sur ocean slab", "W/m2")
953
954! Couplage conv-CL
955 IF (iflag_con.GE.3) THEN
956    IF (iflag_coupl.EQ.1) THEN
957 CALL histdef2d(iff,o_ale_bl%flag,o_ale_bl%name, "ALE BL", "m2/s2")
958 CALL histdef2d(iff,o_alp_bl%flag,o_alp_bl%name, "ALP BL", "m2/s2")
959    ENDIF
960 ENDIF !(iflag_con.GE.3)
961
962 CALL histdef2d(iff,o_weakinv%flag,o_weakinv%name, "Weak inversion", "-")
963 CALL histdef2d(iff,o_dthmin%flag,o_dthmin%name, "dTheta mini", "K/m")
964 CALL histdef2d(iff,o_rh2m%flag,o_rh2m%name, "Relative humidity at 2m", "%" )
965 CALL histdef2d(iff,o_qsat2m%flag,o_qsat2m%name, "Saturant humidity at 2m", "%")
966 CALL histdef2d(iff,o_tpot%flag,o_tpot%name, "Surface air potential temperature", "K")
967 CALL histdef2d(iff,o_tpote%flag,o_tpote%name, "Surface air equivalent potential temperature", "K")
968 CALL histdef2d(iff,o_SWnetOR%flag,o_SWnetOR%name, "Sfce net SW radiation OR", "W/m2")
969 CALL histdef2d(iff,o_SWdownOR%flag,o_SWdownOR%name, "Sfce incident SW radiation OR", "W/m2")
970 CALL histdef2d(iff,o_LWdownOR%flag,o_LWdownOR%name, "Sfce incident LW radiation OR", "W/m2")
971 CALL histdef2d(iff,o_snowl%flag,o_snowl%name, "Solid Large-scale Precip.", "kg/(m2*s)")
972
973 CALL histdef2d(iff,o_solldown%flag,o_solldown%name, "Down. IR rad. at surface", "W/m2")
974 CALL histdef2d(iff,o_dtsvdfo%flag,o_dtsvdfo%name, "Boundary-layer dTs(o)", "K/s")
975 CALL histdef2d(iff,o_dtsvdft%flag,o_dtsvdft%name, "Boundary-layer dTs(t)", "K/s")
976 CALL histdef2d(iff,o_dtsvdfg%flag,o_dtsvdfg%name, "Boundary-layer dTs(g)", "K/s")
977 CALL histdef2d(iff,o_dtsvdfi%flag,o_dtsvdfi%name, "Boundary-layer dTs(g)", "K/s")
978 CALL histdef2d(iff,o_rugs%flag,o_rugs%name, "rugosity", "-" )
979
980! Champs 3D:
981 CALL histdef3d(iff,o_lwcon%flag,o_lwcon%name, "Cloud liquid water content", "kg/kg")
982 CALL histdef3d(iff,o_iwcon%flag,o_iwcon%name, "Cloud ice water content", "kg/kg")
983 CALL histdef3d(iff,o_temp%flag,o_temp%name, "Air temperature", "K" )
984 CALL histdef3d(iff,o_theta%flag,o_theta%name, "Potential air temperature", "K" )
985 CALL histdef3d(iff,o_ovap%flag,o_ovap%name, "Specific humidity + dqphy", "kg/kg" )
986 CALL histdef3d(iff,o_ovapinit%flag,o_ovapinit%name, "Specific humidity", "kg/kg" )
987 CALL histdef3d(iff,o_geop%flag,o_geop%name, "Geopotential height", "m2/s2")
988 CALL histdef3d(iff,o_vitu%flag,o_vitu%name, "Zonal wind", "m/s" )
989 CALL histdef3d(iff,o_vitv%flag,o_vitv%name, "Meridional wind", "m/s" )
990 CALL histdef3d(iff,o_vitw%flag,o_vitw%name, "Vertical wind", "Pa/s" )
991 CALL histdef3d(iff,o_pres%flag,o_pres%name, "Air pressure", "Pa" )
992 CALL histdef3d(iff,o_rneb%flag,o_rneb%name, "Cloud fraction", "-")
993 CALL histdef3d(iff,o_rnebcon%flag,o_rnebcon%name, "Convective Cloud Fraction", "-")
994 CALL histdef3d(iff,o_rhum%flag,o_rhum%name, "Relative humidity", "-")
995 CALL histdef3d(iff,o_ozone%flag,o_ozone%name, "Ozone mole fraction", "-")
996 if (read_climoz == 2) &
997      CALL histdef3d(iff,o_ozone_light%flag,o_ozone_light%name, &
998      "Daylight ozone mole fraction", "-")
999 CALL histdef3d(iff,o_dtphy%flag,o_dtphy%name, "Physics dT", "K/s")
1000 CALL histdef3d(iff,o_dqphy%flag,o_dqphy%name, "Physics dQ", "(kg/kg)/s")
1001 CALL histdef3d(iff,o_cldtau%flag,o_cldtau%name, "Cloud optical thickness", "1")
1002 CALL histdef3d(iff,o_cldemi%flag,o_cldemi%name, "Cloud optical emissivity", "1")
1003!IM: bug ?? dimensionnement variables (klon,klev+1) pmflxr, pmflxs, prfl, psfl
1004 CALL histdef3d(iff,o_pr_con_l%flag,o_pr_con_l%name, "Convective precipitation lic", " ")
1005 CALL histdef3d(iff,o_pr_con_i%flag,o_pr_con_i%name, "Convective precipitation ice", " ")
1006 CALL histdef3d(iff,o_pr_lsc_l%flag,o_pr_lsc_l%name, "Large scale precipitation lic", " ")
1007 CALL histdef3d(iff,o_pr_lsc_i%flag,o_pr_lsc_i%name, "Large scale precipitation ice", " ")
1008!FH Sorties pour la couche limite
1009     if (iflag_pbl>1) then
1010 CALL histdef3d(iff,o_tke%flag,o_tke%name, "TKE", "m2/s2")
1011   type_ecri(1) = 't_max(X)'
1012   type_ecri(2) = 't_max(X)'
1013   type_ecri(3) = 't_max(X)'
1014   type_ecri(4) = 't_max(X)'
1015   type_ecri(5) = 't_max(X)'
1016 CALL histdef3d(iff,o_tke_max%flag,o_tke_max%name, "TKE max", "m2/s2")
1017   type_ecri(:) = type_ecri_files(:)
1018     endif
1019
1020 CALL histdef3d(iff,o_kz%flag,o_kz%name, "Kz melange", "m2/s")
1021   type_ecri(1) = 't_max(X)'
1022   type_ecri(2) = 't_max(X)'
1023   type_ecri(3) = 't_max(X)'
1024   type_ecri(4) = 't_max(X)'
1025   type_ecri(5) = 't_max(X)'
1026 CALL histdef3d(iff,o_kz_max%flag,o_kz_max%name, "Kz melange max", "m2/s" )
1027   type_ecri(:) = type_ecri_files(:)
1028 CALL histdef3d(iff,o_clwcon%flag,o_clwcon%name, "Convective Cloud Liquid water content", "kg/kg")
1029 CALL histdef3d(iff,o_dtdyn%flag,o_dtdyn%name, "Dynamics dT", "K/s")
1030 CALL histdef3d(iff,o_dqdyn%flag,o_dqdyn%name, "Dynamics dQ", "(kg/kg)/s")
1031 CALL histdef3d(iff,o_dudyn%flag,o_dudyn%name, "Dynamics dU", "m/s2")
1032 CALL histdef3d(iff,o_dvdyn%flag,o_dvdyn%name, "Dynamics dV", "m/s2")
1033 CALL histdef3d(iff,o_dtcon%flag,o_dtcon%name, "Convection dT", "K/s")
1034 CALL histdef3d(iff,o_ducon%flag,o_ducon%name, "Convection du", "m/s2")
1035 CALL histdef3d(iff,o_dqcon%flag,o_dqcon%name, "Convection dQ", "(kg/kg)/s")
1036
1037! Wakes
1038 IF(iflag_con.EQ.3) THEN
1039 IF (iflag_wake == 1) THEN
1040   CALL histdef2d(iff,o_ale_wk%flag,o_ale_wk%name, "ALE WK", "m2/s2")
1041   CALL histdef2d(iff,o_alp_wk%flag,o_alp_wk%name, "ALP WK", "m2/s2")
1042   CALL histdef2d(iff,o_ale%flag,o_ale%name, "ALE", "m2/s2")
1043   CALL histdef2d(iff,o_alp%flag,o_alp%name, "ALP", "W/m2")
1044   CALL histdef2d(iff,o_cin%flag,o_cin%name, "Convective INhibition", "m2/s2")
1045   CALL histdef2d(iff,o_wape%flag,o_WAPE%name, "WAPE", "m2/s2")
1046   CALL histdef2d(iff,o_wake_h%flag,o_wake_h%name, "wake_h", "-")
1047   CALL histdef2d(iff,o_wake_s%flag,o_wake_s%name, "wake_s", "-")
1048   CALL histdef3d(iff,o_dtwak%flag,o_dtwak%name, "Wake dT", "K/s")
1049   CALL histdef3d(iff,o_dqwak%flag,o_dqwak%name, "Wake dQ", "(kg/kg)/s")
1050   CALL histdef3d(iff,o_wake_deltat%flag,o_wake_deltat%name, "wake_deltat", " ")
1051   CALL histdef3d(iff,o_wake_deltaq%flag,o_wake_deltaq%name, "wake_deltaq", " ")
1052   CALL histdef3d(iff,o_wake_omg%flag,o_wake_omg%name, "wake_omg", " ")
1053 ENDIF
1054   CALL histdef3d(iff,o_Vprecip%flag,o_Vprecip%name, "precipitation vertical profile", "-")
1055   CALL histdef3d(iff,o_ftd%flag,o_ftd%name, "tend temp due aux descentes precip", "-")
1056   CALL histdef3d(iff,o_fqd%flag,o_fqd%name,"tend vap eau due aux descentes precip", "-")
1057 ENDIF !(iflag_con.EQ.3)
1058
1059 CALL histdef3d(iff,o_dtlsc%flag,o_dtlsc%name, "Condensation dT", "K/s")
1060 CALL histdef3d(iff,o_dtlschr%flag,o_dtlschr%name,"Large-scale condensational heating rate","K/s")
1061 CALL histdef3d(iff,o_dqlsc%flag,o_dqlsc%name, "Condensation dQ", "(kg/kg)/s")
1062 CALL histdef3d(iff,o_dtvdf%flag,o_dtvdf%name, "Boundary-layer dT", "K/s")
1063 CALL histdef3d(iff,o_dqvdf%flag,o_dqvdf%name, "Boundary-layer dQ", "(kg/kg)/s")
1064 CALL histdef3d(iff,o_dteva%flag,o_dteva%name, "Reevaporation dT", "K/s")
1065 CALL histdef3d(iff,o_dqeva%flag,o_dqeva%name, "Reevaporation dQ", "(kg/kg)/s")
1066 CALL histdef3d(iff,o_ptconv%flag,o_ptconv%name, "POINTS CONVECTIFS", " ")
1067 CALL histdef3d(iff,o_ratqs%flag,o_ratqs%name, "RATQS", " ")
1068 CALL histdef3d(iff,o_dtthe%flag,o_dtthe%name, "Dry adjust. dT", "K/s")
1069
1070if(iflag_thermals.gt.1) THEN
1071 CALL histdef3d(iff,o_f_th%flag,o_f_th%name, "Thermal plume mass flux", "K/s")
1072 CALL histdef3d(iff,o_e_th%flag,o_e_th%name,"Thermal plume entrainment","K/s")
1073 CALL histdef3d(iff,o_w_th%flag,o_w_th%name,"Thermal plume vertical velocity","m/s")
1074 CALL histdef3d(iff,o_lambda_th%flag,o_lambda_th%name,"Thermal plume vertical velocity","m/s")
1075 CALL histdef3d(iff,o_q_th%flag,o_q_th%name, "Thermal plume total humidity", "kg/kg")
1076 CALL histdef3d(iff,o_a_th%flag,o_a_th%name, "Thermal plume fraction", "")
1077 CALL histdef3d(iff,o_d_th%flag,o_d_th%name, "Thermal plume detrainment", "K/s")
1078endif !iflag_thermals.gt.1
1079 CALL histdef2d(iff,o_f0_th%flag,o_f0_th%name, "Thermal closure mass flux", "K/s")
1080 CALL histdef2d(iff,o_zmax_th%flag,o_zmax_th%name, "Thermal plume height", "K/s")
1081 CALL histdef3d(iff,o_dqthe%flag,o_dqthe%name, "Dry adjust. dQ", "(kg/kg)/s")
1082 CALL histdef3d(iff,o_dtajs%flag,o_dtajs%name, "Dry adjust. dT", "K/s")
1083 CALL histdef3d(iff,o_dqajs%flag,o_dqajs%name, "Dry adjust. dQ", "(kg/kg)/s")
1084 CALL histdef3d(iff,o_dtswr%flag,o_dtswr%name, "SW radiation dT", "K/s")
1085 CALL histdef3d(iff,o_dtsw0%flag,o_dtsw0%name, "CS SW radiation dT", "K/s")
1086 CALL histdef3d(iff,o_dtlwr%flag,o_dtlwr%name, "LW radiation dT", "K/s")
1087 CALL histdef3d(iff,o_dtlw0%flag,o_dtlw0%name, "CS LW radiation dT", "K/s")
1088 CALL histdef3d(iff,o_dtec%flag,o_dtec%name, "Cinetic dissip dT", "K/s")
1089 CALL histdef3d(iff,o_duvdf%flag,o_duvdf%name, "Boundary-layer dU", "m/s2")
1090 CALL histdef3d(iff,o_dvvdf%flag,o_dvvdf%name, "Boundary-layer dV", "m/s2")
1091
1092     IF (ok_orodr) THEN
1093 CALL histdef3d(iff,o_duoro%flag,o_duoro%name, "Orography dU", "m/s2")
1094 CALL histdef3d(iff,o_dvoro%flag,o_dvoro%name, "Orography dV", "m/s2")
1095     ENDIF
1096
1097     IF (ok_orolf) THEN
1098 CALL histdef3d(iff,o_dulif%flag,o_dulif%name, "Orography dU", "m/s2")
1099 CALL histdef3d(iff,o_dvlif%flag,o_dvlif%name, "Orography dV", "m/s2")
1100     ENDIF
1101
1102      if (nqtot>=3) THEN
1103!Attention    DO iq=3,nqtot
1104    DO iq=3,4 
1105       iiq=niadv(iq)
1106! CALL histdef3d (iff, o_trac%flag,'o_'//tnom(iq)%name,ttext(iiq), "-" )
1107  CALL histdef3d (iff, o_trac(iq-2)%flag,o_trac(iq-2)%name,ttext(iiq), "-" )
1108    ENDDO
1109      endif
1110
1111        CALL histend(nid_files(iff))
1112
1113         ndex2d = 0
1114         ndex3d = 0
1115
1116         ENDIF ! clef_files
1117
1118         ENDDO !  iff
1119     print*,'Fin phys_output_mod.F90'
1120      end subroutine phys_output_open
1121
1122      SUBROUTINE histdef2d (iff,flag_var,nomvar,titrevar,unitvar)
1123     
1124       use ioipsl
1125       USE dimphy
1126       USE mod_phys_lmdz_para
1127
1128       IMPLICIT NONE
1129       
1130       include "dimensions.h"
1131       include "temps.h"
1132       include "indicesol.h"
1133       include "clesphys.h"
1134
1135       integer                          :: iff
1136       integer, dimension(nfiles)       :: flag_var
1137       character(len=20)                 :: nomvar
1138       character(len=*)                 :: titrevar
1139       character(len=*)                 :: unitvar
1140
1141       real zstophym
1142
1143       if (type_ecri(iff)=='inst(X)'.OR.type_ecri(iff)=='once') then
1144         zstophym=zoutm(iff)
1145       else
1146         zstophym=zdtime
1147       endif
1148
1149! Appel a la lecture des noms et niveau d'ecriture des variables dans output.def
1150       call conf_physoutputs(nomvar,flag_var)
1151       
1152       if ( flag_var(iff)<=lev_files(iff) ) then
1153 call histdef (nid_files(iff),nomvar,titrevar,unitvar, &
1154               iim,jj_nb,nhorim(iff), 1,1,1, -99, 32, &
1155               type_ecri(iff), zstophym,zoutm(iff))               
1156       endif                     
1157      end subroutine histdef2d
1158
1159      SUBROUTINE histdef3d (iff,flag_var,nomvar,titrevar,unitvar)
1160
1161       use ioipsl
1162       USE dimphy
1163       USE mod_phys_lmdz_para
1164
1165       IMPLICIT NONE
1166
1167       include "dimensions.h"
1168       include "temps.h"
1169       include "indicesol.h"
1170       include "clesphys.h"
1171
1172       integer                          :: iff
1173       integer, dimension(nfiles)       :: flag_var
1174       character(len=20)                 :: nomvar
1175       character(len=*)                 :: titrevar
1176       character(len=*)                 :: unitvar
1177
1178       real zstophym
1179
1180! Appel a la lecture des noms et niveau d'ecriture des variables dans output.def
1181       call conf_physoutputs(nomvar,flag_var)
1182
1183       if (type_ecri(iff)=='inst(X)'.OR.type_ecri(iff)=='once') then
1184         zstophym=zoutm(iff)
1185       else
1186         zstophym=zdtime
1187       endif
1188
1189       if ( flag_var(iff)<=lev_files(iff) ) then
1190          call histdef (nid_files(iff), nomvar, titrevar, unitvar, &
1191               iim, jj_nb, nhorim(iff), klev, levmin(iff), &
1192               levmax(iff)-levmin(iff)+1, nvertm(iff), 32, type_ecri(iff), &
1193               zstophym, zoutm(iff))
1194       endif
1195      end subroutine histdef3d
1196
1197      SUBROUTINE conf_physoutputs(nam_var,flag_var)
1198!!! Lecture des noms et niveau de sortie des variables dans output.def
1199!   en utilisant les routines getin de IOIPSL 
1200       use ioipsl
1201
1202       IMPLICIT NONE
1203
1204       include 'iniprint.h'
1205
1206       character(len=20)                :: nam_var
1207       integer, dimension(nfiles)      :: flag_var
1208
1209        IF(prt_level>10) WRITE(lunout,*)'Avant getin: nam_var flag_var ',nam_var,flag_var(:)
1210        call getin('flag_'//nam_var,flag_var)
1211        call getin('name_'//nam_var,nam_var)
1212        IF(prt_level>10) WRITE(lunout,*)'Apres getin: nam_var flag_var ',nam_var,flag_var(:)
1213
1214      END SUBROUTINE conf_physoutputs
1215
1216      SUBROUTINE convers_timesteps(str,timestep)
1217
1218        use ioipsl
1219
1220        IMPLICIT NONE
1221
1222        character(len=20)   :: str
1223        character(len=10)   :: type
1224        integer             :: ipos,il
1225        real                :: ttt,xxx,timestep,dayseconde
1226        parameter (dayseconde=86400.)
1227        include "temps.h"
1228        include "comconst.h"
1229
1230        ipos=scan(str,'0123456789.',.true.)
1231
1232        il=len_trim(str)
1233        print*,ipos,il
1234        read(str(1:ipos),*) ttt
1235        print*,ttt
1236        type=str(ipos+1:il)
1237
1238
1239        if ( il == ipos ) then
1240        type='day'
1241        endif
1242
1243        if ( type == 'day'.or.type == 'days'.or.type == 'jours'.or.type == 'jour' ) timestep = ttt * dayseconde
1244        if ( type == 'mounths'.or.type == 'mth'.or.type == 'mois' ) then
1245           print*,'annee_ref,day_ref mon_len',annee_ref,day_ref,ioget_mon_len(annee_ref,day_ref)
1246           timestep = ttt * dayseconde * ioget_mon_len(annee_ref,day_ref)
1247        endif
1248        if ( type == 'hours'.or.type == 'hr'.or.type == 'heurs') timestep = ttt * dayseconde / 24.
1249        if ( type == 'mn'.or.type == 'minutes'  ) timestep = ttt * 60.
1250        if ( type == 's'.or.type == 'sec'.or.type == 'secondes'   ) timestep = ttt
1251        if ( type == 'TS' ) timestep = dtphys
1252
1253        print*,'type =      ',type
1254        print*,'nb j/h/m =  ',ttt
1255        print*,'timestep(s)=',timestep
1256
1257        END SUBROUTINE convers_timesteps
1258
1259END MODULE phys_output_mod
1260
Note: See TracBrowser for help on using the repository browser.