/[lmdze]/trunk/libf/phylmd/ini_hist.f90
ViewVC logotype

Contents of /trunk/libf/phylmd/ini_hist.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 15 - (show annotations)
Fri Aug 1 15:24:12 2008 UTC (15 years, 9 months ago) by guez
File size: 16561 byte(s)
-- Minor modification of input/output:

Added variable "Sigma_O3_Royer" to "histday.nc". "ecrit_day" is not
modified in "physiq". Removed variables "pyu1", "pyv1", "ftsol1",
"ftsol2", "ftsol3", "ftsol4", "psrf1", "psrf2", "psrf3", "psrf4"
"mfu", "mfd", "en_u", "en_d", "de_d", "de_u", "coefh" from
"histrac.nc".

Variable "raz_date" of module "conf_gcm_m" has logical type instead of
integer type.

-- Should not change any result at run time:

Modified calls to "IOIPSL_Lionel" procedures because the interfaces of
these procedures have been simplified.

Changed name of variable in module "start_init_orog_m": "masque" to
"mask".

Created a module containing procedure "phyredem".

Removed arguments "punjours", "pdayref" and "ptimestep" of procedure
"iniphysiq".

Renamed procedure "gr_phy_write" to "gr_phy_write_2d". Created
procedure "gr_phy_write_3d".

Removed procedures "ini_undefstd", "moy_undefSTD", "calcul_STDlev",
"calcul_divers".

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, presnivs, 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 ioipsl, only: ymds2ju, histbeg_totreg, histvert, histend
17 use phyetat0_m, only: rlon, rlat
18
19 REAL, intent(in):: dtime ! pas temporel de la physique (s)
20 real, intent(in):: presnivs(:)
21 integer, intent(out):: nid_hf, nid_hf3d
22
23 REAL zx_lon(iim, jjm + 1), zx_lat(iim, jjm + 1)
24 integer idayref
25 real zjulian
26 integer i, nhori, nvert
27
28 !-----------------------------------------------
29
30 idayref = day_ref
31 CALL ymds2ju(annee_ref, 1, idayref, 0.0, zjulian)
32
33 CALL gr_fi_ecrit(1, klon, iim, (jjm + 1), rlon, zx_lon)
34 DO i = 1, iim
35 zx_lon(i, 1) = rlon(i+1)
36 zx_lon(i, (jjm + 1)) = rlon(i+1)
37 ENDDO
38
39 CALL gr_fi_ecrit(1, klon, iim, (jjm + 1), rlat, zx_lat)
40 CALL histbeg_totreg("histhf", zx_lon(:, 1), zx_lat(1, :), 1, iim, 1, &
41 (jjm + 1), itau_phy, zjulian, dtime, nhori, nid_hf)
42
43 CALL histvert(nid_hf, "presnivs", "Vertical levels", "mb", &
44 llm, presnivs/100., nvert)
45
46 call ini_histhf3d(dtime, presnivs, nid_hf3d)
47 CALL histend(nid_hf)
48
49 end subroutine ini_histhf
50
51 !******************************************************************
52
53 subroutine ini_histhf3d(dtime, presnivs, nid_hf3d)
54
55 ! From phylmd/ini_histhf3d.h, v 1.2 2005/05/25 13:10:09
56
57 ! sorties hf 3d
58
59 use dimens_m, only: iim, jjm, llm
60 use dimphy, only: klon, nbtr
61 use temps, only: itau_phy, day_ref, annee_ref
62 use clesphys, only: ecrit_hf
63 use phyetat0_m, only: rlon, rlat
64 USE ioipsl, only: ymds2ju, histbeg_totreg, histvert, histend, histdef
65
66 REAL, intent(in):: dtime ! pas temporel de la physique (s)
67 real, intent(in):: presnivs(:)
68 integer, intent(out):: nid_hf3d
69
70 real zstohf, zout
71 REAL zx_lon(iim, jjm + 1), zx_lat(iim, jjm + 1)
72 real zjulian
73 integer i, nhori, nvert, idayref
74
75 !------------------------------------------
76
77 zstohf = dtime * REAL(ecrit_hf)
78 zout = dtime * REAL(ecrit_hf)
79
80 idayref = day_ref
81 CALL ymds2ju(annee_ref, 1, idayref, 0.0, zjulian)
82
83 CALL gr_fi_ecrit(1, klon, iim, (jjm + 1), rlon, zx_lon)
84 DO i = 1, iim
85 zx_lon(i, 1) = rlon(i+1)
86 zx_lon(i, (jjm + 1)) = rlon(i+1)
87 ENDDO
88
89 CALL gr_fi_ecrit(1, klon, iim, (jjm + 1), rlat, zx_lat)
90 CALL histbeg_totreg("histhf3d", zx_lon(:, 1), zx_lat(1, :), 1, iim, 1, &
91 (jjm + 1), itau_phy, zjulian, dtime, nhori, nid_hf3d)
92
93 CALL histvert(nid_hf3d, "presnivs", "Vertical levels", "mb", &
94 llm, presnivs/100., nvert)
95
96 ! Champs 3D:
97
98 CALL histdef(nid_hf3d, "temp", "Air temperature", "K", &
99 iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
100 "ave(X)", zstohf, zout)
101
102 CALL histdef(nid_hf3d, "ovap", "Specific humidity", "kg/kg", &
103 iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
104 "ave(X)", zstohf, zout)
105
106 CALL histdef(nid_hf3d, "vitu", "Zonal wind", "m/s", &
107 iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
108 "ave(X)", zstohf, zout)
109
110 CALL histdef(nid_hf3d, "vitv", "Meridional wind", "m/s", &
111 iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
112 "ave(X)", zstohf, zout)
113
114 if (nbtr >= 3) then
115 CALL histdef(nid_hf3d, "O3", "Ozone mass fraction", "?", iim, &
116 (jjm + 1), nhori, llm, 1, llm, nvert, "ave(X)", zstohf, &
117 zout)
118 end if
119
120 CALL histend(nid_hf3d)
121
122 end subroutine ini_histhf3d
123
124 !******************************************************************
125
126 subroutine ini_histday(dtime, presnivs, ok_journe, nid_day)
127
128 ! From phylmd/ini_histday.h, v 1.3 2005/05/25 13:10:09
129
130 use dimens_m, only: iim, jjm, llm
131 use dimphy, only: klon
132 use temps, only: itau_phy, day_ref, annee_ref
133 USE ioipsl, only: ymds2ju, histbeg_totreg, histvert, histend, histdef
134 use phyetat0_m, only: rlon, rlat
135 use clesphys, only: ecrit_day
136
137 REAL, intent(in):: dtime ! pas temporel de la physique (s)
138 real, intent(in):: presnivs(:)
139 logical, intent(in):: ok_journe
140 integer, intent(out):: nid_day
141
142 REAL zx_lon(iim, jjm + 1), zx_lat(iim, jjm + 1)
143 integer i, nhori, nvert, idayref
144 real zjulian
145
146 !--------------------------------
147
148 IF (ok_journe) THEN
149 idayref = day_ref
150 CALL ymds2ju(annee_ref, 1, idayref, 0.0, zjulian)
151
152 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, rlon, zx_lon)
153 DO i = 1, iim
154 zx_lon(i, 1) = rlon(i+1)
155 zx_lon(i, jjm + 1) = rlon(i+1)
156 ENDDO
157 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, rlat, zx_lat)
158 CALL histbeg_totreg("histday", zx_lon(:, 1), zx_lat(1, :), 1, iim, 1, &
159 jjm + 1, itau_phy, zjulian, dtime, nhori, nid_day)
160 CALL histvert(nid_day, "presnivs", "Vertical levels", "mb", &
161 llm, presnivs/100., nvert)
162 call histdef(nid_day, "Sigma_O3_Royer", &
163 "column-density of ozone, in a cell, from Royer", "DU", &
164 pxsize=iim, pysize=jjm+1, phoriid=nhori, pzsize=llm, par_oriz=1, &
165 par_szz=llm, pzid=nvert, popp="ave(X)", pfreq_opp=dtime, &
166 pfreq_wrt=real(ecrit_day))
167 CALL histend(nid_day)
168 ENDIF
169
170 end subroutine ini_histday
171
172 !****************************************************
173
174 subroutine ini_histins(dtime, presnivs, 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 ioipsl, only: ymds2ju, histbeg_totreg, histvert, histend, histdef
184 use phyetat0_m, only: rlon, rlat
185
186 REAL, intent(in):: dtime ! pas temporel de la physique (s)
187 real, intent(in):: presnivs(:)
188 logical, intent(in):: ok_instan
189 integer, intent(out):: nid_ins
190
191 REAL zx_lon(iim, jjm + 1), zx_lat(iim, jjm + 1)
192 real zjulian, zsto, zout
193 integer i, nhori, nvert, idayref, nsrf
194
195 !-------------------------------------------------------------------
196
197 IF (ok_instan) THEN
198
199 zsto = dtime * ecrit_ins
200 zout = dtime * ecrit_ins
201
202 idayref = day_ref
203 CALL ymds2ju(annee_ref, 1, idayref, 0.0, zjulian)
204
205 CALL gr_fi_ecrit(1, klon, iim, (jjm + 1), rlon, zx_lon)
206 DO i = 1, iim
207 zx_lon(i, 1) = rlon(i+1)
208 zx_lon(i, (jjm + 1)) = rlon(i+1)
209 ENDDO
210 CALL gr_fi_ecrit(1, klon, iim, (jjm + 1), rlat, zx_lat)
211 CALL histbeg_totreg("histins", zx_lon(:, 1), zx_lat(1, :), 1, iim, 1, &
212 jjm + 1, itau_phy, zjulian, dtime, nhori, nid_ins)
213 write(*, *)'Inst ', itau_phy, zjulian
214 CALL histvert(nid_ins, "presnivs", "Vertical levels", "mb", &
215 llm, presnivs/100., nvert)
216
217 CALL histdef(nid_ins, "phis", "Surface geop. height", "-", &
218 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
219 "once", zsto, zout)
220
221 CALL histdef(nid_ins, "aire", "Grid area", "-", &
222 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
223 "once", zsto, zout)
224
225 ! Champs 2D:
226
227 CALL histdef(nid_ins, "tsol", "Surface Temperature", "K", &
228 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
229 "inst(X)", zsto, zout)
230
231 CALL histdef(nid_ins, "t2m", "Temperature 2m", "K", &
232 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
233 "inst(X)", zsto, zout)
234
235 CALL histdef(nid_ins, "q2m", "Specific humidity 2m", "Kg/Kg", &
236 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
237 "inst(X)", zsto, zout)
238
239 CALL histdef(nid_ins, "u10m", "Vent zonal 10m", "m/s", &
240 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
241 "inst(X)", zsto, zout)
242
243 CALL histdef(nid_ins, "v10m", "Vent meridien 10m", "m/s", &
244 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
245 "inst(X)", zsto, zout)
246
247 CALL histdef(nid_ins, "psol", "Surface Pressure", "Pa", &
248 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
249 "inst(X)", zsto, zout)
250
251 CALL histdef(nid_ins, "plul", "Large-scale Precip.", "mm/day", &
252 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
253 "inst(X)", zsto, zout)
254
255 CALL histdef(nid_ins, "pluc", "Convective Precip.", "mm/day", &
256 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
257 "inst(X)", zsto, zout)
258
259 CALL histdef(nid_ins, "cdrm", "Momentum drag coef.", "-", &
260 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
261 "inst(X)", zsto, zout)
262
263 CALL histdef(nid_ins, "cdrh", "Heat drag coef.", "-", &
264 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
265 "inst(X)", zsto, zout)
266
267 CALL histdef(nid_ins, "precip", "Precipitation Totale liq+sol", &
268 "kg/(s*m2)", &
269 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
270 "inst(X)", zsto, zout)
271
272 CALL histdef(nid_ins, "snow", "Snow fall", "kg/(s*m2)", &
273 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
274 "inst(X)", zsto, zout)
275
276 ! CALL histdef(nid_ins, "snow_mass", "Snow Mass", "kg/m2",
277 ! . iim, (jjm + 1), nhori, 1, 1, 1, -99,
278 ! . "inst(X)", zsto, zout)
279
280 CALL histdef(nid_ins, "topl", "OLR", "W/m2", &
281 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
282 "inst(X)", zsto, zout)
283
284 CALL histdef(nid_ins, "evap", "Evaporation", "kg/(s*m2)", &
285 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
286 "inst(X)", zsto, zout)
287
288 CALL histdef(nid_ins, "sols", "Solar rad. at surf.", "W/m2", &
289 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
290 "inst(X)", zsto, zout)
291
292 CALL histdef(nid_ins, "soll", "IR rad. at surface", "W/m2", &
293 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
294 "inst(X)", zsto, zout)
295
296 CALL histdef(nid_ins, "solldown", "Down. IR rad. at surface", &
297 "W/m2", iim, (jjm + 1), nhori, 1, 1, 1, -99, &
298 "inst(X)", zsto, zout)
299
300 CALL histdef(nid_ins, "bils", "Surf. total heat flux", "W/m2", &
301 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
302 "inst(X)", zsto, zout)
303
304 CALL histdef(nid_ins, "sens", "Sensible heat flux", "W/m2", &
305 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
306 "inst(X)", zsto, zout)
307
308 CALL histdef(nid_ins, "fder", "Heat flux derivation", "W/m2", &
309 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
310 "inst(X)", zsto, zout)
311
312 CALL histdef(nid_ins, "dtsvdfo", "Boundary-layer dTs(o)", "K/s", &
313 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
314 "inst(X)", zsto, zout)
315
316 CALL histdef(nid_ins, "dtsvdft", "Boundary-layer dTs(t)", "K/s", &
317 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
318 "inst(X)", zsto, zout)
319
320 CALL histdef(nid_ins, "dtsvdfg", "Boundary-layer dTs(g)", "K/s", &
321 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
322 "inst(X)", zsto, zout)
323
324 CALL histdef(nid_ins, "dtsvdfi", "Boundary-layer dTs(g)", "K/s", &
325 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
326 "inst(X)", zsto, zout)
327
328 DO nsrf = 1, nbsrf
329
330 call histdef(nid_ins, "pourc_"//clnsurf(nsrf), &
331 "% "//clnsurf(nsrf), "%", &
332 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
333 "inst(X)", zsto, zout)
334
335 call histdef(nid_ins, "fract_"//clnsurf(nsrf), &
336 "Fraction "//clnsurf(nsrf), "1", &
337 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
338 "inst(X)", zsto, zout)
339
340 call histdef(nid_ins, "sens_"//clnsurf(nsrf), &
341 "Sensible heat flux "//clnsurf(nsrf), "W/m2", &
342 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
343 "inst(X)", zsto, zout)
344
345 call histdef(nid_ins, "tsol_"//clnsurf(nsrf), &
346 "Surface Temperature"//clnsurf(nsrf), "W/m2", &
347 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
348 "inst(X)", zsto, zout)
349
350 call histdef(nid_ins, "lat_"//clnsurf(nsrf), &
351 "Latent heat flux "//clnsurf(nsrf), "W/m2", &
352 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
353 "inst(X)", zsto, zout)
354
355 call histdef(nid_ins, "taux_"//clnsurf(nsrf), &
356 "Zonal wind stress"//clnsurf(nsrf), "Pa", &
357 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
358 "inst(X)", zsto, zout)
359
360 call histdef(nid_ins, "tauy_"//clnsurf(nsrf), &
361 "Meridional xind stress "//clnsurf(nsrf), "Pa", &
362 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
363 "inst(X)", zsto, zout)
364
365 call histdef(nid_ins, "albe_"//clnsurf(nsrf), &
366 "Albedo "//clnsurf(nsrf), "-", &
367 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
368 "inst(X)", zsto, zout)
369
370 call histdef(nid_ins, "rugs_"//clnsurf(nsrf), &
371 "rugosite "//clnsurf(nsrf), "-", &
372 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
373 "inst(X)", zsto, zout)
374 !XXX
375 END DO
376 CALL histdef(nid_ins, "rugs", "rugosity", "-", &
377 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
378 "inst(X)", zsto, zout)
379
380 CALL histdef(nid_ins, "albs", "Surface albedo", "-", &
381 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
382 "inst(X)", zsto, zout)
383 CALL histdef(nid_ins, "albslw", "Surface albedo LW", "-", &
384 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
385 "inst(X)", zsto, zout)
386
387 !IM cf. AM 081204 BEG
388 ! HBTM2
389 CALL histdef(nid_ins, "s_pblh", "Boundary Layer Height", "m", &
390 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
391 "inst(X)", zsto, zout)
392
393 CALL histdef(nid_ins, "s_pblt", "T at Boundary Layer Height", &
394 "K", &
395 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
396 "inst(X)", zsto, zout)
397
398 CALL histdef(nid_ins, "s_lcl", "Condensation level", "m", &
399 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
400 "inst(X)", zsto, zout)
401
402 CALL histdef(nid_ins, "s_capCL", "Conv avlbl pot ener for ABL", "J/m2", &
403 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
404 "inst(X)", zsto, zout)
405
406 CALL histdef(nid_ins, "s_oliqCL", "Liq Water in BL", "kg/m2", &
407 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
408 "inst(X)", zsto, zout)
409
410 CALL histdef(nid_ins, "s_cteiCL", "Instability criteria (ABL)", "K", &
411 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
412 "inst(X)", zsto, zout)
413
414 CALL histdef(nid_ins, "s_therm", "Exces du thermique", "K", &
415 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
416 "inst(X)", zsto, zout)
417
418 CALL histdef(nid_ins, "s_trmb1", "deep_cape(HBTM2)", "J/m2", &
419 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
420 "inst(X)", zsto, zout)
421
422 CALL histdef(nid_ins, "s_trmb2", "inhibition (HBTM2)", "J/m2", &
423 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
424 "inst(X)", zsto, zout)
425
426 CALL histdef(nid_ins, "s_trmb3", "Point Omega (HBTM2)", "m", &
427 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
428 "inst(X)", zsto, zout)
429
430 !IM cf. AM 081204 END
431
432 ! Champs 3D:
433
434 CALL histdef(nid_ins, "temp", "Temperature", "K", &
435 iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
436 "inst(X)", zsto, zout)
437
438 CALL histdef(nid_ins, "vitu", "Zonal wind", "m/s", &
439 iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
440 "inst(X)", zsto, zout)
441
442 CALL histdef(nid_ins, "vitv", "Merid wind", "m/s", &
443 iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
444 "inst(X)", zsto, zout)
445
446 CALL histdef(nid_ins, "geop", "Geopotential height", "m", &
447 iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
448 "inst(X)", zsto, zout)
449
450 CALL histdef(nid_ins, "pres", "Air pressure", "Pa", &
451 iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
452 "inst(X)", zsto, zout)
453
454 CALL histdef(nid_ins, "dtvdf", "Boundary-layer dT", "K/s", &
455 iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
456 "inst(X)", zsto, zout)
457
458 CALL histdef(nid_ins, "dqvdf", "Boundary-layer dQ", "Kg/Kg/s", &
459 iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
460 "inst(X)", zsto, zout)
461
462 CALL histend(nid_ins)
463 ENDIF
464
465 end subroutine ini_histins
466
467 end module ini_hist

  ViewVC Help
Powered by ViewVC 1.1.21