/[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 217 - (show annotations)
Thu Mar 30 14:25:18 2017 UTC (7 years, 1 month ago) by guez
File size: 13584 byte(s)
run_off_lic downgraded from variable of module interface_surf to local
variable of fonte_neige.

Code could not work with ok_aie set to true, so removed this
possibility. tauae, piz_ae, cg_ae, topswai, solswai were then
0. cldtaupi was the same as cldtaupd.

In sw and procedures called by sw, flag_aer did not need to be double
precision, changed it to logical.

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

  ViewVC Help
Powered by ViewVC 1.1.21