/[lmdze]/trunk/Sources/phylmd/ini_histins.f
ViewVC logotype

Contents of /trunk/Sources/phylmd/ini_histins.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 215 - (show annotations)
Tue Mar 28 12:46:28 2017 UTC (7 years, 1 month ago) by guez
File size: 13113 byte(s)
size(snow) is now knon in interfsurf_hq.

Renamed snow to fsnow in clmain, same name as corresponding actual
argument. We can then rename ysnow to simply snow in clmain, same name
as corresponding dummy argument of clqh. No need to initialize local
snow to 0 since it is only used with indices 1:knon and already
initialized from fsnow for each type of surface. If there is no point
for a given type of surface, fsnow should be reset to 0 for this
type. We need to give a valid value to fsnow in this case even if it
will be multiplied by pctsrf = 0 in physiq.

In physiq, no need for intermediate zxsnow for output.

Removed unused arguments tsurf, p1lay, beta, coef1lay, ps, t1lay,
q1lay, u1lay, v1lay, petAcoef, peqAcoef, petBcoef, peqBcoef of
fonte_neige, with unused computations of zx_qs and zcor. (Same was
done in LMDZ.)

1 module ini_histins_m
2
3 implicit none
4
5 integer, save:: nid_ins
6
7 contains
8
9 subroutine ini_histins(dtime)
10
11 ! From phylmd/ini_histins.h, version 1.2, 2005/05/25 13:10:09
12
13 use clesphys, only: ecrit_ins, ok_instan
14 use clesphys2, only: conv_emanuel
15 use dimens_m, only: iim, jjm, llm, nqmx
16 use disvert_m, only: presnivs
17 use dynetat0_m, only: day_ref, annee_ref, rlatu, rlonv
18 USE histbeg_totreg_m, ONLY : histbeg_totreg
19 USE histdef_m, ONLY : histdef
20 USE histend_m, ONLY : histend
21 USE histvert_m, ONLY : histvert
22 use indicesol, only: nbsrf, clnsurf
23 use iniadvtrac_m, only: tname, ttext
24 use nr_util, only: pi
25 use phyetat0_m, only: itau_phy
26 USE ymds2ju_m, only: ymds2ju
27
28 REAL, intent(in):: dtime ! pas temporel de la physique (s)
29
30 ! Local:
31 real zjulian, zsto, zout
32 integer nhori, nvert, nsrf, iq, it
33
34 !-------------------------------------------------------------------
35
36 print *, 'Call sequence information: ini_histins'
37
38 test_ok_instan: IF (ok_instan) THEN
39 zsto = dtime * ecrit_ins
40 zout = dtime * ecrit_ins
41 CALL ymds2ju(annee_ref, 1, day_ref, 0.0, zjulian)
42 CALL histbeg_totreg("histins", rlonv(:iim) / pi * 180., &
43 rlatu / pi * 180., 1, iim, &
44 1, jjm + 1, itau_phy, zjulian, dtime, nhori, nid_ins)
45 print *, 'itau_phy = ', itau_phy
46 print *, "zjulian = ", zjulian
47 CALL histvert(nid_ins, "presnivs", "Vertical levels", "mb", &
48 presnivs/100., nvert)
49
50 CALL histdef(nid_ins, "phis", "Surface geop. height", "-", &
51 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
52 "once", zsto, zout)
53 CALL histdef(nid_ins, "aire", "Grid area", "-", &
54 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
55 "once", zsto, zout)
56
57 ! Champs 2D:
58
59 CALL histdef(nid_ins, "tsol", "Surface Temperature", "K", &
60 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
61 "inst(X)", zsto, zout)
62 CALL histdef(nid_ins, "t2m", "Temperature 2m", "K", &
63 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
64 "inst(X)", zsto, zout)
65 CALL histdef(nid_ins, "q2m", "Specific humidity 2m", "Kg/Kg", &
66 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
67 "inst(X)", zsto, zout)
68 CALL histdef(nid_ins, "u10m", "Vent zonal 10m", "m/s", &
69 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
70 "inst(X)", zsto, zout)
71 CALL histdef(nid_ins, "v10m", "Vent meridien 10m", "m/s", &
72 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
73 "inst(X)", zsto, zout)
74 CALL histdef(nid_ins, "psol", "Surface Pressure", "Pa", &
75 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
76 "inst(X)", zsto, zout)
77 CALL histdef(nid_ins, "plul", "Large-scale Precip.", "mm/day", &
78 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
79 "inst(X)", zsto, zout)
80 CALL histdef(nid_ins, "pluc", "Convective Precip.", "mm/day", &
81 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
82 "inst(X)", zsto, zout)
83 CALL histdef(nid_ins, "cdrm", "Momentum drag coef.", "-", &
84 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
85 "inst(X)", zsto, zout)
86 CALL histdef(nid_ins, "cdrh", "Heat drag coef.", "-", &
87 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
88 "inst(X)", zsto, zout)
89 CALL histdef(nid_ins, "precip", "Precipitation Totale liq+sol", &
90 "kg/(s*m2)", &
91 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
92 "inst(X)", zsto, zout)
93 CALL histdef(nid_ins, "snow", "Snow fall", "kg/(s*m2)", &
94 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
95 "inst(X)", zsto, zout)
96 CALL histdef(nid_ins, "topl", "OLR", "W/m2", &
97 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
98 "inst(X)", zsto, zout)
99 CALL histdef(nid_ins, "evap", "Evaporation", "kg/(s*m2)", &
100 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
101 "inst(X)", zsto, zout)
102 CALL histdef(nid_ins, "sols", "Solar rad. at surf.", "W/m2", &
103 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
104 "inst(X)", zsto, zout)
105 CALL histdef(nid_ins, "soll", "IR rad. at surface", "W/m2", &
106 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
107 "inst(X)", zsto, zout)
108 CALL histdef(nid_ins, "solldown", "Down. IR rad. at surface", &
109 "W/m2", iim, (jjm + 1), nhori, 1, 1, 1, -99, &
110 "inst(X)", zsto, zout)
111 CALL histdef(nid_ins, "bils", "Surf. total heat flux", "W/m2", &
112 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
113 "inst(X)", zsto, zout)
114 CALL histdef(nid_ins, "sens", "Sensible heat flux", "W/m2", &
115 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
116 "inst(X)", zsto, zout)
117 CALL histdef(nid_ins, "fder", "Heat flux derivation", "W/m2", &
118 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
119 "inst(X)", zsto, zout)
120 CALL histdef(nid_ins, "dtsvdfo", "Boundary-layer dTs(o)", "K/s", &
121 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
122 "inst(X)", zsto, zout)
123 CALL histdef(nid_ins, "dtsvdft", "Boundary-layer dTs(t)", "K/s", &
124 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
125 "inst(X)", zsto, zout)
126 CALL histdef(nid_ins, "dtsvdfg", "Boundary-layer dTs(g)", "K/s", &
127 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
128 "inst(X)", zsto, zout)
129 CALL histdef(nid_ins, "dtsvdfi", "Boundary-layer dTs(g)", "K/s", &
130 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
131 "inst(X)", zsto, zout)
132 CALL histdef(nid_ins, "msnow", "surface snow amount", "kg/m2", &
133 iim, jjm + 1, nhori, 1, 1, 1, -99, "inst(X)", zsto, zout)
134
135 DO nsrf = 1, nbsrf
136 call histdef(nid_ins, "pourc_"//clnsurf(nsrf), &
137 "% "//clnsurf(nsrf), "%", &
138 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
139 "inst(X)", zsto, zout)
140 call histdef(nid_ins, "fract_"//clnsurf(nsrf), &
141 "Fraction "//clnsurf(nsrf), "1", &
142 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
143 "inst(X)", zsto, zout)
144 call histdef(nid_ins, "sens_"//clnsurf(nsrf), &
145 "Sensible heat flux "//clnsurf(nsrf), "W/m2", &
146 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
147 "inst(X)", zsto, zout)
148 call histdef(nid_ins, "tsol_"//clnsurf(nsrf), &
149 "Surface Temperature"//clnsurf(nsrf), "W/m2", &
150 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
151 "inst(X)", zsto, zout)
152 call histdef(nid_ins, "lat_"//clnsurf(nsrf), &
153 "Latent heat flux "//clnsurf(nsrf), "W/m2", &
154 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
155 "inst(X)", zsto, zout)
156 call histdef(nid_ins, "taux_"//clnsurf(nsrf), &
157 "Zonal wind stress"//clnsurf(nsrf), "Pa", &
158 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
159 "inst(X)", zsto, zout)
160 call histdef(nid_ins, "tauy_"//clnsurf(nsrf), &
161 "Meridional xind stress "//clnsurf(nsrf), "Pa", &
162 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
163 "inst(X)", zsto, zout)
164 call histdef(nid_ins, "albe_"//clnsurf(nsrf), &
165 "Albedo "//clnsurf(nsrf), "-", &
166 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
167 "inst(X)", zsto, zout)
168 call histdef(nid_ins, "rugs_"//clnsurf(nsrf), &
169 "rugosite "//clnsurf(nsrf), "-", &
170 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
171 "inst(X)", zsto, zout)
172 END DO
173
174 CALL histdef(nid_ins, "rugs", "rugosity", "-", &
175 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
176 "inst(X)", zsto, zout)
177 CALL histdef(nid_ins, "albs", "Surface albedo", "-", &
178 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
179 "inst(X)", zsto, zout)
180 CALL histdef(nid_ins, "s_pblh", "Boundary Layer Height", "m", &
181 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
182 "inst(X)", zsto, zout)
183 CALL histdef(nid_ins, "s_pblt", "T at Boundary Layer Height", &
184 "K", &
185 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
186 "inst(X)", zsto, zout)
187 CALL histdef(nid_ins, "s_lcl", "Condensation level", "m", &
188 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
189 "inst(X)", zsto, zout)
190 CALL histdef(nid_ins, "s_capCL", "Conv avlbl pot ener for ABL", "J/m2", &
191 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
192 "inst(X)", zsto, zout)
193 CALL histdef(nid_ins, "s_oliqCL", "Liq Water in BL", "kg/m2", &
194 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
195 "inst(X)", zsto, zout)
196 CALL histdef(nid_ins, "s_cteiCL", "Instability criteria (ABL)", "K", &
197 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
198 "inst(X)", zsto, zout)
199 CALL histdef(nid_ins, "s_therm", "Exces du thermique", "K", &
200 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
201 "inst(X)", zsto, zout)
202 CALL histdef(nid_ins, "s_trmb1", "deep_cape(HBTM2)", "J/m2", &
203 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
204 "inst(X)", zsto, zout)
205 CALL histdef(nid_ins, "s_trmb2", "inhibition (HBTM2)", "J/m2", &
206 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
207 "inst(X)", zsto, zout)
208 CALL histdef(nid_ins, "s_trmb3", "Point Omega (HBTM2)", "m", &
209 iim, (jjm + 1), nhori, 1, 1, 1, -99, &
210 "inst(X)", zsto, zout)
211
212 if (conv_emanuel) then
213 CALL histdef(nid_ins, "ptop", "cloud top pressure", &
214 "Pa", iim, jjm + 1, nhori, 1, 1, 1, -99, "inst(X)", zsto, zout)
215 CALL histdef(nid_ins, "dnwd0", "unsaturated downdraft", &
216 "kg/m2/s", iim, jjm + 1, nhori, llm, 1, llm, nvert, "inst(X)", &
217 zsto, zout)
218 end if
219
220 ! Champs 3D:
221
222 CALL histdef(nid_ins, "tro3", "ozone mole fraction", "-", &
223 iim, jjm + 1, nhori, llm, 1, llm, nvert, "inst(X)", zsto, zout)
224 CALL histdef(nid_ins, "temp", "Temperature", "K", &
225 iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
226 "inst(X)", zsto, zout)
227 CALL histdef(nid_ins, "vitu", "Zonal wind", "m/s", &
228 iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
229 "inst(X)", zsto, zout)
230 CALL histdef(nid_ins, "vitv", "Merid wind", "m/s", &
231 iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
232 "inst(X)", zsto, zout)
233 CALL histdef(nid_ins, "geop", "Geopotential height", "m", &
234 iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
235 "inst(X)", zsto, zout)
236 CALL histdef(nid_ins, "pres", "Air pressure", "Pa", &
237 iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
238 "inst(X)", zsto, zout)
239 CALL histdef(nid_ins, "dtvdf", "Boundary-layer dT", "K/s", &
240 iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
241 "inst(X)", zsto, zout)
242 CALL histdef(nid_ins, "dqvdf", "Boundary-layer dQ", "Kg/Kg/s", &
243 iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
244 "inst(X)", zsto, zout)
245 CALL histdef(nid_ins, "zmasse", "column density of air in cell", &
246 "kg m-2", iim, jjm + 1, nhori, llm, 1, llm, nvert, "inst(X)", &
247 zsto, zout)
248 CALL histdef(nid_ins, "rhum", "Relative humidity", &
249 "", iim, jjm + 1, nhori, llm, 1, llm, nvert, "inst(X)", &
250 zsto, zout)
251 CALL histdef(nid_ins, "d_t_ec", "kinetic dissipation dT", &
252 "K/s", iim, jjm + 1, nhori, llm, 1, llm, nvert, "inst(X)", &
253 zsto, zout)
254 CALL histdef(nid_ins, "dtsw0", "CS SW radiation dT", &
255 "K/s", iim, jjm + 1, nhori, llm, 1, llm, nvert, "inst(X)", &
256 zsto, zout)
257 CALL histdef(nid_ins, "dtlw0", "CS LW radiation dT", &
258 "K/s", iim, jjm + 1, nhori, llm, 1, llm, nvert, "inst(X)", &
259 zsto, zout)
260
261 DO it = 1, nqmx - 2
262 ! champ 2D
263 iq=it+2
264 CALL histdef(nid_ins, tname(iq), ttext(iq), "U/kga", iim, jjm+1, &
265 nhori, llm, 1, llm, nvert, "inst(X)", zsto, zout)
266 CALL histdef(nid_ins, "fl"//tname(iq), "Flux "//ttext(iq), &
267 "U/m2/s", iim, jjm+1, nhori, llm, 1, llm, nvert, &
268 "inst(X)", zsto, zout)
269 CALL histdef(nid_ins, "d_tr_th_"//tname(iq), &
270 "tendance thermique"// ttext(iq), "?", &
271 iim, jjm+1, nhori, llm, 1, llm, nvert, &
272 "inst(X)", zsto, zout)
273 CALL histdef(nid_ins, "d_tr_cv_"//tname(iq), &
274 "tendance convection"// ttext(iq), "?", &
275 iim, jjm+1, nhori, llm, 1, llm, nvert, &
276 "inst(X)", zsto, zout)
277 CALL histdef(nid_ins, "d_tr_cl_"//tname(iq), &
278 "tendance couche limite"// ttext(iq), "?", &
279 iim, jjm+1, nhori, llm, 1, llm, nvert, &
280 "inst(X)", zsto, zout)
281 ENDDO
282
283 CALL histend(nid_ins)
284 ENDIF test_ok_instan
285
286 end subroutine ini_histins
287
288 end module ini_histins_m

  ViewVC Help
Powered by ViewVC 1.1.21