1 |
module ini_hist |
2 |
|
3 |
! This module is clean: no C preprocessor directive, no include line. |
4 |
|
5 |
IMPLICIT none |
6 |
|
7 |
contains |
8 |
|
9 |
subroutine ini_histhf(dtime, nid_hf, nid_hf3d) |
10 |
|
11 |
! From phylmd/ini_histhf.h, version 1.3 2005/05/25 13:10:09 |
12 |
|
13 |
use dimens_m, only: iim, jjm, llm |
14 |
use temps, only: day_ref, annee_ref, itau_phy |
15 |
use dimphy, only: klon |
16 |
USE calendar, only: ymds2ju |
17 |
use histcom, only: histbeg_totreg, histvert, histend |
18 |
use phyetat0_m, only: rlon, rlat |
19 |
use comvert, only: presnivs |
20 |
|
21 |
REAL, intent(in):: dtime ! pas temporel de la physique (s) |
22 |
integer, intent(out):: nid_hf, nid_hf3d |
23 |
|
24 |
REAL zx_lon(iim, jjm + 1), zx_lat(iim, jjm + 1) |
25 |
integer idayref |
26 |
real zjulian |
27 |
integer i, nhori, nvert |
28 |
|
29 |
!----------------------------------------------- |
30 |
|
31 |
idayref = day_ref |
32 |
CALL ymds2ju(annee_ref, 1, idayref, 0.0, zjulian) |
33 |
|
34 |
CALL gr_fi_ecrit(1, klon, iim, (jjm + 1), rlon, zx_lon) |
35 |
DO i = 1, iim |
36 |
zx_lon(i, 1) = rlon(i+1) |
37 |
zx_lon(i, (jjm + 1)) = rlon(i+1) |
38 |
ENDDO |
39 |
|
40 |
CALL gr_fi_ecrit(1, klon, iim, (jjm + 1), rlat, zx_lat) |
41 |
CALL histbeg_totreg("histhf", zx_lon(:, 1), zx_lat(1, :), 1, iim, 1, & |
42 |
(jjm + 1), itau_phy, zjulian, dtime, nhori, nid_hf) |
43 |
|
44 |
CALL histvert(nid_hf, "presnivs", "Vertical levels", "mb", & |
45 |
llm, presnivs/100., nvert) |
46 |
|
47 |
call ini_histhf3d(dtime, nid_hf3d) |
48 |
CALL histend(nid_hf) |
49 |
|
50 |
end subroutine ini_histhf |
51 |
|
52 |
!****************************************************************** |
53 |
|
54 |
subroutine ini_histhf3d(dtime, nid_hf3d) |
55 |
|
56 |
! From phylmd/ini_histhf3d.h, v 1.2 2005/05/25 13:10:09 |
57 |
|
58 |
! sorties hf 3d |
59 |
|
60 |
use dimens_m, only: iim, jjm, llm |
61 |
use dimphy, only: klon, nbtr |
62 |
use temps, only: itau_phy, day_ref, annee_ref |
63 |
use clesphys, only: ecrit_hf |
64 |
use phyetat0_m, only: rlon, rlat |
65 |
USE calendar, only: ymds2ju |
66 |
use histcom, only: histbeg_totreg, histvert, histend, histdef |
67 |
use comvert, only: presnivs |
68 |
|
69 |
REAL, intent(in):: dtime ! pas temporel de la physique (s) |
70 |
integer, intent(out):: nid_hf3d |
71 |
|
72 |
real zstohf, zout |
73 |
REAL zx_lon(iim, jjm + 1), zx_lat(iim, jjm + 1) |
74 |
real zjulian |
75 |
integer i, nhori, nvert, idayref |
76 |
|
77 |
!------------------------------------------ |
78 |
|
79 |
zstohf = dtime * REAL(ecrit_hf) |
80 |
zout = dtime * REAL(ecrit_hf) |
81 |
|
82 |
idayref = day_ref |
83 |
CALL ymds2ju(annee_ref, 1, idayref, 0.0, zjulian) |
84 |
|
85 |
CALL gr_fi_ecrit(1, klon, iim, (jjm + 1), rlon, zx_lon) |
86 |
DO i = 1, iim |
87 |
zx_lon(i, 1) = rlon(i+1) |
88 |
zx_lon(i, (jjm + 1)) = rlon(i+1) |
89 |
ENDDO |
90 |
|
91 |
CALL gr_fi_ecrit(1, klon, iim, (jjm + 1), rlat, zx_lat) |
92 |
CALL histbeg_totreg("histhf3d", zx_lon(:, 1), zx_lat(1, :), 1, iim, 1, & |
93 |
(jjm + 1), itau_phy, zjulian, dtime, nhori, nid_hf3d) |
94 |
|
95 |
CALL histvert(nid_hf3d, "presnivs", "Vertical levels", "mb", & |
96 |
llm, presnivs/100., nvert) |
97 |
|
98 |
! Champs 3D: |
99 |
|
100 |
CALL histdef(nid_hf3d, "temp", "Air temperature", "K", & |
101 |
iim, (jjm + 1), nhori, llm, 1, llm, nvert, & |
102 |
"ave(X)", zstohf, zout) |
103 |
|
104 |
CALL histdef(nid_hf3d, "ovap", "Specific humidity", "kg/kg", & |
105 |
iim, (jjm + 1), nhori, llm, 1, llm, nvert, & |
106 |
"ave(X)", zstohf, zout) |
107 |
|
108 |
CALL histdef(nid_hf3d, "vitu", "Zonal wind", "m/s", & |
109 |
iim, (jjm + 1), nhori, llm, 1, llm, nvert, & |
110 |
"ave(X)", zstohf, zout) |
111 |
|
112 |
CALL histdef(nid_hf3d, "vitv", "Meridional wind", "m/s", & |
113 |
iim, (jjm + 1), nhori, llm, 1, llm, nvert, & |
114 |
"ave(X)", zstohf, zout) |
115 |
|
116 |
if (nbtr >= 3) then |
117 |
CALL histdef(nid_hf3d, "O3", "Ozone mass fraction", "?", iim, & |
118 |
(jjm + 1), nhori, llm, 1, llm, nvert, "ave(X)", zstohf, & |
119 |
zout) |
120 |
end if |
121 |
|
122 |
CALL histend(nid_hf3d) |
123 |
|
124 |
end subroutine ini_histhf3d |
125 |
|
126 |
!****************************************************************** |
127 |
|
128 |
subroutine ini_histday(dtime, ok_journe, nid_day, nq) |
129 |
|
130 |
! From phylmd/ini_histday.h, v 1.3 2005/05/25 13:10:09 |
131 |
|
132 |
use dimens_m, only: iim, jjm, llm |
133 |
use temps, only: itau_phy, day_ref, annee_ref |
134 |
USE calendar, only: ymds2ju |
135 |
use histcom, only: histbeg_totreg, histvert, histend, histdef |
136 |
use phyetat0_m, only: rlon, rlat |
137 |
use clesphys, only: ecrit_day |
138 |
use grid_change, only: gr_phy_write_2d |
139 |
use comvert, only: presnivs |
140 |
|
141 |
REAL, intent(in):: dtime ! pas temporel de la physique (s) |
142 |
logical, intent(in):: ok_journe |
143 |
integer, intent(out):: nid_day |
144 |
INTEGER, intent(in):: nq ! nombre de traceurs (y compris vapeur d'eau) |
145 |
|
146 |
! Variables local to the procedure: |
147 |
REAL zx_lat(iim, jjm + 1) |
148 |
integer nhori, nvert |
149 |
real zjulian |
150 |
|
151 |
!-------------------------------- |
152 |
|
153 |
IF (ok_journe) THEN |
154 |
CALL ymds2ju(annee_ref, 1, day_ref, 0., zjulian) |
155 |
zx_lat = gr_phy_write_2d(rlat) |
156 |
CALL histbeg_totreg("histday", rlon(2: iim+1), zx_lat(1, :), 1, iim, & |
157 |
1, jjm + 1, itau_phy, zjulian, dtime, nhori, nid_day) |
158 |
CALL histvert(nid_day, "presnivs", "Vertical levels", "mb", & |
159 |
llm, presnivs/100., nvert) |
160 |
if (nq <= 4) then |
161 |
call histdef(nid_day, "Sigma_O3_Royer", & |
162 |
"column-density of ozone, in a cell, from Royer", "DU", & |
163 |
pxsize=iim, pysize=jjm+1, phoriid=nhori, pzsize=llm, & |
164 |
par_oriz=1, par_szz=llm, pzid=nvert, popp="ave(X)", & |
165 |
pfreq_opp=dtime, pfreq_wrt=real(ecrit_day)) |
166 |
end if |
167 |
CALL histend(nid_day) |
168 |
ENDIF |
169 |
|
170 |
end subroutine ini_histday |
171 |
|
172 |
!**************************************************** |
173 |
|
174 |
subroutine ini_histins(dtime, ok_instan, nid_ins) |
175 |
|
176 |
! From phylmd/ini_histins.h, v 1.2 2005/05/25 13:10:09 |
177 |
|
178 |
use dimens_m, only: iim, jjm, llm |
179 |
use dimphy, only: klon |
180 |
use temps, only: itau_phy, day_ref, annee_ref |
181 |
use clesphys, only: ecrit_ins |
182 |
use indicesol, only: nbsrf, clnsurf |
183 |
USE calendar, only: ymds2ju |
184 |
use histcom, only: histbeg_totreg, histvert, histend, histdef |
185 |
use phyetat0_m, only: rlon, rlat |
186 |
use comvert, only: presnivs |
187 |
|
188 |
REAL, intent(in):: dtime ! pas temporel de la physique (s) |
189 |
logical, intent(in):: ok_instan |
190 |
integer, intent(out):: nid_ins |
191 |
|
192 |
REAL zx_lon(iim, jjm + 1), zx_lat(iim, jjm + 1) |
193 |
real zjulian, zsto, zout |
194 |
integer i, nhori, nvert, idayref, nsrf |
195 |
|
196 |
!------------------------------------------------------------------- |
197 |
|
198 |
IF (ok_instan) THEN |
199 |
|
200 |
zsto = dtime * ecrit_ins |
201 |
zout = dtime * ecrit_ins |
202 |
|
203 |
idayref = day_ref |
204 |
CALL ymds2ju(annee_ref, 1, idayref, 0.0, zjulian) |
205 |
|
206 |
CALL gr_fi_ecrit(1, klon, iim, (jjm + 1), rlon, zx_lon) |
207 |
DO i = 1, iim |
208 |
zx_lon(i, 1) = rlon(i+1) |
209 |
zx_lon(i, (jjm + 1)) = rlon(i+1) |
210 |
ENDDO |
211 |
CALL gr_fi_ecrit(1, klon, iim, (jjm + 1), rlat, zx_lat) |
212 |
CALL histbeg_totreg("histins", zx_lon(:, 1), zx_lat(1, :), 1, iim, 1, & |
213 |
jjm + 1, itau_phy, zjulian, dtime, nhori, nid_ins) |
214 |
write(*, *)'Inst ', itau_phy, zjulian |
215 |
CALL histvert(nid_ins, "presnivs", "Vertical levels", "mb", & |
216 |
llm, presnivs/100., nvert) |
217 |
|
218 |
CALL histdef(nid_ins, "phis", "Surface geop. height", "-", & |
219 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
220 |
"once", zsto, zout) |
221 |
|
222 |
CALL histdef(nid_ins, "aire", "Grid area", "-", & |
223 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
224 |
"once", zsto, zout) |
225 |
|
226 |
! Champs 2D: |
227 |
|
228 |
CALL histdef(nid_ins, "tsol", "Surface Temperature", "K", & |
229 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
230 |
"inst(X)", zsto, zout) |
231 |
|
232 |
CALL histdef(nid_ins, "t2m", "Temperature 2m", "K", & |
233 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
234 |
"inst(X)", zsto, zout) |
235 |
|
236 |
CALL histdef(nid_ins, "q2m", "Specific humidity 2m", "Kg/Kg", & |
237 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
238 |
"inst(X)", zsto, zout) |
239 |
|
240 |
CALL histdef(nid_ins, "u10m", "Vent zonal 10m", "m/s", & |
241 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
242 |
"inst(X)", zsto, zout) |
243 |
|
244 |
CALL histdef(nid_ins, "v10m", "Vent meridien 10m", "m/s", & |
245 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
246 |
"inst(X)", zsto, zout) |
247 |
|
248 |
CALL histdef(nid_ins, "psol", "Surface Pressure", "Pa", & |
249 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
250 |
"inst(X)", zsto, zout) |
251 |
|
252 |
CALL histdef(nid_ins, "plul", "Large-scale Precip.", "mm/day", & |
253 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
254 |
"inst(X)", zsto, zout) |
255 |
|
256 |
CALL histdef(nid_ins, "pluc", "Convective Precip.", "mm/day", & |
257 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
258 |
"inst(X)", zsto, zout) |
259 |
|
260 |
CALL histdef(nid_ins, "cdrm", "Momentum drag coef.", "-", & |
261 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
262 |
"inst(X)", zsto, zout) |
263 |
|
264 |
CALL histdef(nid_ins, "cdrh", "Heat drag coef.", "-", & |
265 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
266 |
"inst(X)", zsto, zout) |
267 |
|
268 |
CALL histdef(nid_ins, "precip", "Precipitation Totale liq+sol", & |
269 |
"kg/(s*m2)", & |
270 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
271 |
"inst(X)", zsto, zout) |
272 |
|
273 |
CALL histdef(nid_ins, "snow", "Snow fall", "kg/(s*m2)", & |
274 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
275 |
"inst(X)", zsto, zout) |
276 |
|
277 |
! CALL histdef(nid_ins, "snow_mass", "Snow Mass", "kg/m2", |
278 |
! . iim, (jjm + 1), nhori, 1, 1, 1, -99, |
279 |
! . "inst(X)", zsto, zout) |
280 |
|
281 |
CALL histdef(nid_ins, "topl", "OLR", "W/m2", & |
282 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
283 |
"inst(X)", zsto, zout) |
284 |
|
285 |
CALL histdef(nid_ins, "evap", "Evaporation", "kg/(s*m2)", & |
286 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
287 |
"inst(X)", zsto, zout) |
288 |
|
289 |
CALL histdef(nid_ins, "sols", "Solar rad. at surf.", "W/m2", & |
290 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
291 |
"inst(X)", zsto, zout) |
292 |
|
293 |
CALL histdef(nid_ins, "soll", "IR rad. at surface", "W/m2", & |
294 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
295 |
"inst(X)", zsto, zout) |
296 |
|
297 |
CALL histdef(nid_ins, "solldown", "Down. IR rad. at surface", & |
298 |
"W/m2", iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
299 |
"inst(X)", zsto, zout) |
300 |
|
301 |
CALL histdef(nid_ins, "bils", "Surf. total heat flux", "W/m2", & |
302 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
303 |
"inst(X)", zsto, zout) |
304 |
|
305 |
CALL histdef(nid_ins, "sens", "Sensible heat flux", "W/m2", & |
306 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
307 |
"inst(X)", zsto, zout) |
308 |
|
309 |
CALL histdef(nid_ins, "fder", "Heat flux derivation", "W/m2", & |
310 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
311 |
"inst(X)", zsto, zout) |
312 |
|
313 |
CALL histdef(nid_ins, "dtsvdfo", "Boundary-layer dTs(o)", "K/s", & |
314 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
315 |
"inst(X)", zsto, zout) |
316 |
|
317 |
CALL histdef(nid_ins, "dtsvdft", "Boundary-layer dTs(t)", "K/s", & |
318 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
319 |
"inst(X)", zsto, zout) |
320 |
|
321 |
CALL histdef(nid_ins, "dtsvdfg", "Boundary-layer dTs(g)", "K/s", & |
322 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
323 |
"inst(X)", zsto, zout) |
324 |
|
325 |
CALL histdef(nid_ins, "dtsvdfi", "Boundary-layer dTs(g)", "K/s", & |
326 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
327 |
"inst(X)", zsto, zout) |
328 |
|
329 |
DO nsrf = 1, nbsrf |
330 |
|
331 |
call histdef(nid_ins, "pourc_"//clnsurf(nsrf), & |
332 |
"% "//clnsurf(nsrf), "%", & |
333 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
334 |
"inst(X)", zsto, zout) |
335 |
|
336 |
call histdef(nid_ins, "fract_"//clnsurf(nsrf), & |
337 |
"Fraction "//clnsurf(nsrf), "1", & |
338 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
339 |
"inst(X)", zsto, zout) |
340 |
|
341 |
call histdef(nid_ins, "sens_"//clnsurf(nsrf), & |
342 |
"Sensible heat flux "//clnsurf(nsrf), "W/m2", & |
343 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
344 |
"inst(X)", zsto, zout) |
345 |
|
346 |
call histdef(nid_ins, "tsol_"//clnsurf(nsrf), & |
347 |
"Surface Temperature"//clnsurf(nsrf), "W/m2", & |
348 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
349 |
"inst(X)", zsto, zout) |
350 |
|
351 |
call histdef(nid_ins, "lat_"//clnsurf(nsrf), & |
352 |
"Latent heat flux "//clnsurf(nsrf), "W/m2", & |
353 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
354 |
"inst(X)", zsto, zout) |
355 |
|
356 |
call histdef(nid_ins, "taux_"//clnsurf(nsrf), & |
357 |
"Zonal wind stress"//clnsurf(nsrf), "Pa", & |
358 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
359 |
"inst(X)", zsto, zout) |
360 |
|
361 |
call histdef(nid_ins, "tauy_"//clnsurf(nsrf), & |
362 |
"Meridional xind stress "//clnsurf(nsrf), "Pa", & |
363 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
364 |
"inst(X)", zsto, zout) |
365 |
|
366 |
call histdef(nid_ins, "albe_"//clnsurf(nsrf), & |
367 |
"Albedo "//clnsurf(nsrf), "-", & |
368 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
369 |
"inst(X)", zsto, zout) |
370 |
|
371 |
call histdef(nid_ins, "rugs_"//clnsurf(nsrf), & |
372 |
"rugosite "//clnsurf(nsrf), "-", & |
373 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
374 |
"inst(X)", zsto, zout) |
375 |
!XXX |
376 |
END DO |
377 |
CALL histdef(nid_ins, "rugs", "rugosity", "-", & |
378 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
379 |
"inst(X)", zsto, zout) |
380 |
|
381 |
CALL histdef(nid_ins, "albs", "Surface albedo", "-", & |
382 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
383 |
"inst(X)", zsto, zout) |
384 |
CALL histdef(nid_ins, "albslw", "Surface albedo LW", "-", & |
385 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
386 |
"inst(X)", zsto, zout) |
387 |
|
388 |
!IM cf. AM 081204 BEG |
389 |
! HBTM2 |
390 |
CALL histdef(nid_ins, "s_pblh", "Boundary Layer Height", "m", & |
391 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
392 |
"inst(X)", zsto, zout) |
393 |
|
394 |
CALL histdef(nid_ins, "s_pblt", "T at Boundary Layer Height", & |
395 |
"K", & |
396 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
397 |
"inst(X)", zsto, zout) |
398 |
|
399 |
CALL histdef(nid_ins, "s_lcl", "Condensation level", "m", & |
400 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
401 |
"inst(X)", zsto, zout) |
402 |
|
403 |
CALL histdef(nid_ins, "s_capCL", "Conv avlbl pot ener for ABL", "J/m2", & |
404 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
405 |
"inst(X)", zsto, zout) |
406 |
|
407 |
CALL histdef(nid_ins, "s_oliqCL", "Liq Water in BL", "kg/m2", & |
408 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
409 |
"inst(X)", zsto, zout) |
410 |
|
411 |
CALL histdef(nid_ins, "s_cteiCL", "Instability criteria (ABL)", "K", & |
412 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
413 |
"inst(X)", zsto, zout) |
414 |
|
415 |
CALL histdef(nid_ins, "s_therm", "Exces du thermique", "K", & |
416 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
417 |
"inst(X)", zsto, zout) |
418 |
|
419 |
CALL histdef(nid_ins, "s_trmb1", "deep_cape(HBTM2)", "J/m2", & |
420 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
421 |
"inst(X)", zsto, zout) |
422 |
|
423 |
CALL histdef(nid_ins, "s_trmb2", "inhibition (HBTM2)", "J/m2", & |
424 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
425 |
"inst(X)", zsto, zout) |
426 |
|
427 |
CALL histdef(nid_ins, "s_trmb3", "Point Omega (HBTM2)", "m", & |
428 |
iim, (jjm + 1), nhori, 1, 1, 1, -99, & |
429 |
"inst(X)", zsto, zout) |
430 |
|
431 |
!IM cf. AM 081204 END |
432 |
|
433 |
! Champs 3D: |
434 |
|
435 |
CALL histdef(nid_ins, "temp", "Temperature", "K", & |
436 |
iim, (jjm + 1), nhori, llm, 1, llm, nvert, & |
437 |
"inst(X)", zsto, zout) |
438 |
|
439 |
CALL histdef(nid_ins, "vitu", "Zonal wind", "m/s", & |
440 |
iim, (jjm + 1), nhori, llm, 1, llm, nvert, & |
441 |
"inst(X)", zsto, zout) |
442 |
|
443 |
CALL histdef(nid_ins, "vitv", "Merid wind", "m/s", & |
444 |
iim, (jjm + 1), nhori, llm, 1, llm, nvert, & |
445 |
"inst(X)", zsto, zout) |
446 |
|
447 |
CALL histdef(nid_ins, "geop", "Geopotential height", "m", & |
448 |
iim, (jjm + 1), nhori, llm, 1, llm, nvert, & |
449 |
"inst(X)", zsto, zout) |
450 |
|
451 |
CALL histdef(nid_ins, "pres", "Air pressure", "Pa", & |
452 |
iim, (jjm + 1), nhori, llm, 1, llm, nvert, & |
453 |
"inst(X)", zsto, zout) |
454 |
|
455 |
CALL histdef(nid_ins, "dtvdf", "Boundary-layer dT", "K/s", & |
456 |
iim, (jjm + 1), nhori, llm, 1, llm, nvert, & |
457 |
"inst(X)", zsto, zout) |
458 |
|
459 |
CALL histdef(nid_ins, "dqvdf", "Boundary-layer dQ", "Kg/Kg/s", & |
460 |
iim, (jjm + 1), nhori, llm, 1, llm, nvert, & |
461 |
"inst(X)", zsto, zout) |
462 |
|
463 |
CALL histend(nid_ins) |
464 |
ENDIF |
465 |
|
466 |
end subroutine ini_histins |
467 |
|
468 |
!************************************************* |
469 |
|
470 |
subroutine ini_histrac(nid_tra, pdtphys, nq_phys, lessivage) |
471 |
|
472 |
! From phylmd/ini_histrac.h, version 1.10 2006/02/21 08:08:30 |
473 |
|
474 |
use dimens_m, only: iim, jjm, llm |
475 |
USE calendar, only: ymds2ju |
476 |
use histcom, only: histbeg_totreg, histvert, histend, histdef |
477 |
use temps, only: annee_ref, day_ref, itau_phy |
478 |
use iniadvtrac_m, only: niadv, tnom, ttext |
479 |
use dimphy, only: klon |
480 |
use clesphys, only: ecrit_tra |
481 |
use grid_change, only: gr_phy_write_2d |
482 |
use phyetat0_m, only: rlon, rlat |
483 |
use comvert, only: presnivs |
484 |
|
485 |
INTEGER, intent(out):: nid_tra |
486 |
real, intent(in):: pdtphys ! pas d'integration pour la physique (s) |
487 |
|
488 |
integer, intent(in):: nq_phys |
489 |
! (nombre de traceurs auxquels on applique la physique) |
490 |
|
491 |
logical, intent(in):: lessivage |
492 |
|
493 |
! Variables local to the procedure: |
494 |
|
495 |
REAL zjulian |
496 |
REAL zx_lat(iim, jjm+1) |
497 |
INTEGER nhori, nvert |
498 |
REAL zsto, zout |
499 |
integer it, iq, iiq |
500 |
|
501 |
!--------------------------------------------------------- |
502 |
|
503 |
CALL ymds2ju(annee_ref, month=1, day=day_ref, sec=0.0, julian=zjulian) |
504 |
zx_lat(:, :) = gr_phy_write_2d(rlat) |
505 |
CALL histbeg_totreg("histrac", rlon(2:iim+1), zx_lat(1, :), & |
506 |
1, iim, 1, jjm+1, itau_phy, zjulian, pdtphys, nhori, nid_tra) |
507 |
CALL histvert(nid_tra, "presnivs", "Vertical levels", "mb", llm, & |
508 |
presnivs, nvert) |
509 |
|
510 |
zsto = pdtphys |
511 |
zout = pdtphys * REAL(ecrit_tra) |
512 |
|
513 |
CALL histdef(nid_tra, "phis", "Surface geop. height", "-", & |
514 |
iim, jjm+1, nhori, 1, 1, 1, -99, & |
515 |
"once", zsto, zout) |
516 |
CALL histdef(nid_tra, "aire", "Grid area", "-", & |
517 |
iim, jjm+1, nhori, 1, 1, 1, -99, & |
518 |
"once", zsto, zout) |
519 |
CALL histdef(nid_tra, "zmasse", "column density of air in cell", & |
520 |
"kg m-2", iim, jjm + 1, nhori, llm, 1, llm, nvert, "ave(X)", & |
521 |
zsto, zout) |
522 |
|
523 |
DO it = 1, nq_phys |
524 |
! champ 2D |
525 |
iq=it+2 |
526 |
iiq=niadv(iq) |
527 |
CALL histdef(nid_tra, tnom(iq), ttext(iiq), "U/kga", iim, jjm+1, & |
528 |
nhori, llm, 1, llm, nvert, "ave(X)", zsto, zout) |
529 |
if (lessivage) THEN |
530 |
CALL histdef(nid_tra, "fl"//tnom(iq), "Flux "//ttext(iiq), & |
531 |
"U/m2/s", iim, jjm+1, nhori, llm, 1, llm, nvert, & |
532 |
"ave(X)", zsto, zout) |
533 |
endif |
534 |
|
535 |
!---Ajout Olivia |
536 |
CALL histdef(nid_tra, "d_tr_th_"//tnom(iq), & |
537 |
"tendance thermique"// ttext(iiq), "?", & |
538 |
iim, jjm+1, nhori, llm, 1, llm, nvert, & |
539 |
"ave(X)", zsto, zout) |
540 |
CALL histdef(nid_tra, "d_tr_cv_"//tnom(iq), & |
541 |
"tendance convection"// ttext(iiq), "?", & |
542 |
iim, jjm+1, nhori, llm, 1, llm, nvert, & |
543 |
"ave(X)", zsto, zout) |
544 |
CALL histdef(nid_tra, "d_tr_cl_"//tnom(iq), & |
545 |
"tendance couche limite"// ttext(iiq), "?", & |
546 |
iim, jjm+1, nhori, llm, 1, llm, nvert, & |
547 |
"ave(X)", zsto, zout) |
548 |
!---fin Olivia |
549 |
|
550 |
ENDDO |
551 |
|
552 |
CALL histdef(nid_tra, "pplay", "", "-", & |
553 |
iim, jjm+1, nhori, llm, 1, llm, nvert, & |
554 |
"inst(X)", zout, zout) |
555 |
CALL histdef(nid_tra, "T", "temperature", "K", iim, jjm+1, nhori, llm, & |
556 |
1, llm, nvert, "inst(X)", zout, zout) |
557 |
|
558 |
CALL histend(nid_tra) |
559 |
|
560 |
end subroutine ini_histrac |
561 |
|
562 |
end module ini_hist |