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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 72 - (show annotations)
Tue Jul 23 13:00:07 2013 UTC (10 years, 9 months ago) by guez
File size: 14347 byte(s)
NaN to signalling NaN in gfortran_debug.mk.

Removed unused procedures in getincom and getincom2. In procedure
conf_interface, replaced call to getincom by new namelist. Moved
procedure conf_interface into module interface_surf.

Added variables sig1 and w01 to startphy.nc and restartphy.nc, for
procedure cv_driver. Renamed (ema_)?work1 and (ema_)?work2 to sig1 and
w01 in concvl and physiq.

Deleted unused arguments of clmain, clqh and intersurf_hq, among which
(y)?sollwdown. Following LMDZ, in physiq, read sollw instead of
sollwdown from startphy.nc, write sollw instead of sollwdown to
restartphy.nc.

In procedure sw, initialized zfs[ud][pn]a[di], for runs where ok_ade
and ok_aie are false. (Following LMDZ.)

Added dimension klev to startphy.nc and restartphy.nc, and deleted
dimension horizon_vertical. Made t_ancien and q_ancien two-dimensional
NetCDF variables. Bug fix: in phyetat0, define ratqs, clwcon and
rnebcon for vertical levels >=2.

Bug fix: set mfg, p[de]n_[ud] to 0. when iflag_con >= 3. (Following LMDZ.)

1 module phyredem_m
2
3 IMPLICIT NONE
4
5 contains
6
7 SUBROUTINE phyredem(fichnom, rlat, rlon, pctsrf, tsol, tsoil, tslab, &
8 seaice, qsurf, qsol, snow, albedo, alblw, evap, rain_fall, snow_fall, &
9 solsw, sollw, fder, radsol, frugs, agesno, zmea, zstd, zsig, zgam, &
10 zthe, zpic, zval, t_ancien, q_ancien, rnebcon, ratqs, clwcon, &
11 run_off_lic_0, sig1, w01)
12
13 ! From phylmd/phyredem.F, version 1.3 2005/05/25 13:10:09
14 ! Author: Z. X. Li (LMD/CNRS)
15 ! Date: 19930818
16 ! Objet : écriture de l'état de démarrage ou redémarrage pour la physique
17
18 USE dimphy, ONLY: klev, klon, zmasq
19 USE dimsoil, ONLY: nsoilmx
20 USE indicesol, ONLY: is_lic, is_oce, is_sic, is_ter, nbsrf
21 USE netcdf, ONLY: nf90_clobber, nf90_global, nf90_float
22 USE netcdf95, ONLY: nf95_create, nf95_put_att, nf95_def_dim, &
23 nf95_def_var, nf95_enddef, nf95_redef, nf95_put_var, nf95_close
24 USE temps, ONLY: itau_phy
25
26 CHARACTER(len=*) fichnom
27 REAL, INTENT(IN):: rlat(klon), rlon(klon)
28 REAL pctsrf(klon, nbsrf)
29 REAL tsol(klon, nbsrf)
30 REAL tsoil(klon, nsoilmx, nbsrf)
31 REAL tslab(klon), seaice(klon) !IM "slab" ocean
32 REAL qsurf(klon, nbsrf)
33 REAL, intent(in):: qsol(klon)
34 REAL snow(klon, nbsrf)
35 REAL albedo(klon, nbsrf)
36 REAL alblw(klon, nbsrf)
37 REAL evap(klon, nbsrf)
38 REAL, INTENT(IN):: rain_fall(klon)
39 REAL snow_fall(klon)
40 REAL solsw(klon)
41 REAL, INTENT(IN):: sollw(klon)
42 REAL fder(klon)
43 REAL radsol(klon)
44 REAL frugs(klon, nbsrf)
45 REAL agesno(klon, nbsrf)
46 REAL zmea(klon)
47 REAL, intent(in):: zstd(klon)
48 REAL, intent(in):: zsig(klon)
49 REAL zgam(klon)
50 REAL zthe(klon)
51 REAL zpic(klon)
52 REAL zval(klon)
53 REAL t_ancien(klon, klev), q_ancien(klon, klev)
54 REAL rnebcon(klon, klev), ratqs(klon, klev), clwcon(klon, klev)
55 REAL run_off_lic_0(klon)
56 real, intent(in):: sig1(klon, klev) ! section adiabatic updraft
57
58 real, intent(in):: w01(klon, klev)
59 ! vertical velocity within adiabatic updraft
60
61 ! Local:
62
63 INTEGER ncid, idim2, idim3
64 integer varid, varid_run_off_lic_0, varid_sig1, varid_w01
65 INTEGER isoil, nsrf
66 CHARACTER(len=7) str7
67 CHARACTER(len=2) str2
68
69 !------------------------------------------------------------
70
71 PRINT *, 'Call sequence information: phyredem'
72 CALL nf95_create(fichnom, nf90_clobber, ncid)
73
74 call nf95_put_att(ncid, nf90_global, 'title', &
75 'Fichier redémarrage physique')
76 call nf95_put_att(ncid, nf90_global, "itau_phy", itau_phy)
77
78 call nf95_def_dim(ncid, 'points_physiques', klon, idim2)
79 call nf95_def_dim(ncid, 'klev', klev, idim3)
80
81 call nf95_def_var(ncid, 'longitude', nf90_float, idim2, varid)
82 call nf95_put_att(ncid, varid, 'title', &
83 'Longitudes de la grille physique')
84 call nf95_enddef(ncid)
85 call nf95_put_var(ncid, varid, rlon)
86
87 call nf95_redef(ncid)
88 call nf95_def_var(ncid, 'latitude', nf90_float, idim2, varid)
89 call nf95_put_att(ncid, varid, 'title', &
90 'Latitudes de la grille physique')
91 call nf95_enddef(ncid)
92 call nf95_put_var(ncid, varid, rlat)
93
94 ! PB ajout du masque terre/mer
95
96 call nf95_redef(ncid)
97 call nf95_def_var(ncid, 'masque', nf90_float, idim2, varid)
98 call nf95_put_att(ncid, varid, 'title', 'masque terre mer')
99 call nf95_enddef(ncid)
100 call nf95_put_var(ncid, varid, zmasq)
101 ! BP ajout des fraction de chaque sous-surface
102
103 ! 1. fraction de terre
104
105 call nf95_redef(ncid)
106 call nf95_def_var(ncid, 'FTER', nf90_float, idim2, varid)
107 call nf95_put_att(ncid, varid, 'title', 'fraction de continent')
108 call nf95_enddef(ncid)
109 call nf95_put_var(ncid, varid, pctsrf(:, is_ter))
110
111 ! 2. Fraction de glace de terre
112
113 call nf95_redef(ncid)
114 call nf95_def_var(ncid, 'FLIC', nf90_float, idim2, varid)
115 call nf95_put_att(ncid, varid, 'title', 'fraction glace de terre')
116 call nf95_enddef(ncid)
117 call nf95_put_var(ncid, varid, pctsrf(:, is_lic))
118
119 ! 3. fraction ocean
120
121 call nf95_redef(ncid)
122 call nf95_def_var(ncid, 'FOCE', nf90_float, idim2, varid)
123 call nf95_put_att(ncid, varid, 'title', 'fraction ocean')
124 call nf95_enddef(ncid)
125 call nf95_put_var(ncid, varid, pctsrf(:, is_oce))
126
127 ! 4. Fraction glace de mer
128
129 call nf95_redef(ncid)
130 call nf95_def_var(ncid, 'FSIC', nf90_float, idim2, varid)
131 call nf95_put_att(ncid, varid, 'title', 'fraction glace mer')
132 call nf95_enddef(ncid)
133 call nf95_put_var(ncid, varid, pctsrf(:, is_sic))
134
135 DO nsrf = 1, nbsrf
136 IF (nsrf<=99) THEN
137 WRITE (str2, '(i2.2)') nsrf
138 call nf95_redef(ncid)
139 call nf95_def_var(ncid, 'TS'//str2, nf90_float, idim2, varid)
140 call nf95_put_att(ncid, varid, 'title', &
141 'Temperature de surface No.'//str2)
142 call nf95_enddef(ncid)
143 ELSE
144 PRINT *, 'Trop de sous-mailles'
145 STOP 1
146 END IF
147 call nf95_put_var(ncid, varid, tsol(:, nsrf))
148 END DO
149
150 DO nsrf = 1, nbsrf
151 DO isoil = 1, nsoilmx
152 IF (isoil<=99 .AND. nsrf<=99) THEN
153 WRITE (str7, '(i2.2, "srf", i2.2)') isoil, nsrf
154 call nf95_redef(ncid)
155 call nf95_def_var(ncid, 'Tsoil'//str7, nf90_float, idim2, varid)
156 call nf95_put_att(ncid, varid, 'title', &
157 'Temperature du sol No.'//str7)
158 call nf95_enddef(ncid)
159 ELSE
160 PRINT *, 'Trop de couches'
161 STOP 1
162 END IF
163 call nf95_put_var(ncid, varid, tsoil(:, isoil, nsrf))
164 END DO
165 END DO
166
167 !IM "slab" ocean
168 call nf95_redef(ncid)
169 call nf95_def_var(ncid, 'TSLAB', nf90_float, idim2, varid)
170 call nf95_put_att(ncid, varid, 'title', &
171 'Ecart de la SST (pour slab-ocean)')
172 call nf95_enddef(ncid)
173 call nf95_put_var(ncid, varid, tslab)
174
175 call nf95_redef(ncid)
176 call nf95_def_var(ncid, 'SEAICE', nf90_float, idim2, varid)
177 call nf95_put_att(ncid, varid, 'title', &
178 'Glace de mer kg/m2 (pour slab-ocean)')
179 call nf95_enddef(ncid)
180 call nf95_put_var(ncid, varid, seaice)
181
182 DO nsrf = 1, nbsrf
183 IF (nsrf<=99) THEN
184 WRITE (str2, '(i2.2)') nsrf
185 call nf95_redef(ncid)
186 call nf95_def_var(ncid, 'QS'//str2, nf90_float, idim2, varid)
187 call nf95_put_att(ncid, varid, 'title', &
188 'Humidite de surface No.'//str2)
189 call nf95_enddef(ncid)
190 ELSE
191 PRINT *, 'Trop de sous-mailles'
192 STOP 1
193 END IF
194 call nf95_put_var(ncid, varid, qsurf(:, nsrf))
195 END DO
196
197 call nf95_redef(ncid)
198 call nf95_def_var(ncid, 'QSOL', nf90_float, idim2, varid)
199 call nf95_put_att(ncid, varid, 'title', 'Eau dans le sol (mm)')
200 call nf95_enddef(ncid)
201 call nf95_put_var(ncid, varid, qsol)
202
203 DO nsrf = 1, nbsrf
204 IF (nsrf<=99) THEN
205 WRITE (str2, '(i2.2)') nsrf
206 call nf95_redef(ncid)
207 call nf95_def_var(ncid, 'ALBE'//str2, nf90_float, idim2, varid)
208 call nf95_put_att(ncid, varid, 'title', &
209 'albedo de surface No.'//str2)
210 call nf95_enddef(ncid)
211 ELSE
212 PRINT *, 'Trop de sous-mailles'
213 STOP 1
214 END IF
215 call nf95_put_var(ncid, varid, albedo(:, nsrf))
216 END DO
217
218 !IM BEG albedo LW
219 DO nsrf = 1, nbsrf
220 IF (nsrf<=99) THEN
221 WRITE (str2, '(i2.2)') nsrf
222 call nf95_redef(ncid)
223 call nf95_def_var(ncid, 'ALBLW'//str2, nf90_float, idim2, varid)
224 call nf95_put_att(ncid, varid, 'title', &
225 'albedo LW de surface No.'//str2)
226 call nf95_enddef(ncid)
227 ELSE
228 PRINT *, 'Trop de sous-mailles'
229 STOP 1
230 END IF
231 call nf95_put_var(ncid, varid, alblw(:, nsrf))
232 END DO
233 !IM END albedo LW
234
235 DO nsrf = 1, nbsrf
236 IF (nsrf<=99) THEN
237 WRITE (str2, '(i2.2)') nsrf
238 call nf95_redef(ncid)
239 call nf95_def_var(ncid, 'EVAP'//str2, nf90_float, idim2, varid)
240 call nf95_put_att(ncid, varid, 'title', &
241 'Evaporation de surface No.'//str2)
242 call nf95_enddef(ncid)
243 ELSE
244 PRINT *, 'Trop de sous-mailles'
245 STOP 1
246 END IF
247 call nf95_put_var(ncid, varid, evap(:, nsrf))
248 END DO
249
250 DO nsrf = 1, nbsrf
251 IF (nsrf<=99) THEN
252 WRITE (str2, '(i2.2)') nsrf
253 call nf95_redef(ncid)
254 call nf95_def_var(ncid, 'SNOW'//str2, nf90_float, idim2, varid)
255 call nf95_put_att(ncid, varid, 'title', &
256 'Neige de surface No.'//str2)
257 call nf95_enddef(ncid)
258 ELSE
259 PRINT *, 'Trop de sous-mailles'
260 STOP 1
261 END IF
262 call nf95_put_var(ncid, varid, snow(:, nsrf))
263 END DO
264
265 call nf95_redef(ncid)
266 call nf95_def_var(ncid, 'RADS', nf90_float, idim2, varid)
267 call nf95_put_att(ncid, varid, 'title', &
268 'Rayonnement net a la surface')
269 call nf95_enddef(ncid)
270 call nf95_put_var(ncid, varid, radsol)
271
272 call nf95_redef(ncid)
273 call nf95_def_var(ncid, 'solsw', nf90_float, idim2, varid)
274 call nf95_put_att(ncid, varid, 'title', &
275 'Rayonnement solaire a la surface')
276 call nf95_enddef(ncid)
277 call nf95_put_var(ncid, varid, solsw)
278
279 call nf95_redef(ncid)
280 call nf95_def_var(ncid, 'sollw', nf90_float, idim2, varid)
281 call nf95_put_att(ncid, varid, 'title', &
282 'Rayonnement IF a la surface')
283 call nf95_enddef(ncid)
284 call nf95_put_var(ncid, varid, sollw)
285
286 call nf95_redef(ncid)
287 call nf95_def_var(ncid, 'fder', nf90_float, idim2, varid)
288 call nf95_put_att(ncid, varid, 'title', 'Derive de flux')
289 call nf95_enddef(ncid)
290 call nf95_put_var(ncid, varid, fder)
291
292 call nf95_redef(ncid)
293 call nf95_def_var(ncid, 'rain_f', nf90_float, idim2, varid)
294 call nf95_put_att(ncid, varid, 'title', 'precipitation liquide')
295 call nf95_enddef(ncid)
296 call nf95_put_var(ncid, varid, rain_fall)
297
298 call nf95_redef(ncid)
299 call nf95_def_var(ncid, 'snow_f', nf90_float, idim2, varid)
300 call nf95_put_att(ncid, varid, 'title', 'precipitation solide')
301 call nf95_enddef(ncid)
302 call nf95_put_var(ncid, varid, snow_fall)
303
304 DO nsrf = 1, nbsrf
305 IF (nsrf<=99) THEN
306 WRITE (str2, '(i2.2)') nsrf
307 call nf95_redef(ncid)
308 call nf95_def_var(ncid, 'RUG'//str2, nf90_float, idim2, varid)
309 call nf95_put_att(ncid, varid, 'title', &
310 'rugosite de surface No.'//str2)
311 call nf95_enddef(ncid)
312 ELSE
313 PRINT *, 'Trop de sous-mailles'
314 STOP 1
315 END IF
316 call nf95_put_var(ncid, varid, frugs(:, nsrf))
317 END DO
318
319 DO nsrf = 1, nbsrf
320 IF (nsrf<=99) THEN
321 WRITE (str2, '(i2.2)') nsrf
322 call nf95_redef(ncid)
323 call nf95_def_var(ncid, 'AGESNO'//str2, nf90_float, idim2, varid)
324 call nf95_put_att(ncid, varid, 'title', &
325 'Age de la neige surface No.'//str2)
326 call nf95_enddef(ncid)
327 ELSE
328 PRINT *, 'Trop de sous-mailles'
329 STOP 1
330 END IF
331 call nf95_put_var(ncid, varid, agesno(:, nsrf))
332 END DO
333
334 call nf95_redef(ncid)
335 call nf95_def_var(ncid, 'ZMEA', nf90_float, idim2, varid)
336 call nf95_enddef(ncid)
337 call nf95_put_var(ncid, varid, zmea)
338
339 call nf95_redef(ncid)
340 call nf95_def_var(ncid, 'ZSTD', nf90_float, idim2, varid)
341 call nf95_enddef(ncid)
342 call nf95_put_var(ncid, varid, zstd)
343 call nf95_redef(ncid)
344 call nf95_def_var(ncid, 'ZSIG', nf90_float, idim2, varid)
345 call nf95_enddef(ncid)
346 call nf95_put_var(ncid, varid, zsig)
347 call nf95_redef(ncid)
348 call nf95_def_var(ncid, 'ZGAM', nf90_float, idim2, varid)
349 call nf95_enddef(ncid)
350 call nf95_put_var(ncid, varid, zgam)
351 call nf95_redef(ncid)
352 call nf95_def_var(ncid, 'ZTHE', nf90_float, idim2, varid)
353 call nf95_enddef(ncid)
354 call nf95_put_var(ncid, varid, zthe)
355 call nf95_redef(ncid)
356 call nf95_def_var(ncid, 'ZPIC', nf90_float, idim2, varid)
357 call nf95_enddef(ncid)
358 call nf95_put_var(ncid, varid, zpic)
359 call nf95_redef(ncid)
360 call nf95_def_var(ncid, 'ZVAL', nf90_float, idim2, varid)
361 call nf95_enddef(ncid)
362 call nf95_put_var(ncid, varid, zval)
363
364 call nf95_redef(ncid)
365 call nf95_def_var(ncid, 'TANCIEN', nf90_float, (/idim2, idim3/), varid)
366 call nf95_enddef(ncid)
367 call nf95_put_var(ncid, varid, t_ancien)
368
369 call nf95_redef(ncid)
370 call nf95_def_var(ncid, 'QANCIEN', nf90_float, (/idim2, idim3/), varid)
371 call nf95_enddef(ncid)
372 call nf95_put_var(ncid, varid, q_ancien)
373
374 call nf95_redef(ncid)
375 call nf95_def_var(ncid, 'RUGMER', nf90_float, idim2, varid)
376 call nf95_put_att(ncid, varid, 'title', &
377 'Longueur de rugosite sur mer')
378 call nf95_enddef(ncid)
379 call nf95_put_var(ncid, varid, frugs(:, is_oce))
380
381 call nf95_redef(ncid)
382 call nf95_def_var(ncid, 'CLWCON', nf90_float, idim2, varid)
383 call nf95_put_att(ncid, varid, 'title', 'Eau liquide convective')
384 call nf95_enddef(ncid)
385 call nf95_put_var(ncid, varid, clwcon(:, 1))
386
387 call nf95_redef(ncid)
388 call nf95_def_var(ncid, 'RNEBCON', nf90_float, idim2, varid)
389 call nf95_put_att(ncid, varid, 'title', 'Nebulosite convective')
390 call nf95_enddef(ncid)
391 call nf95_put_var(ncid, varid, rnebcon(:, 1))
392
393 call nf95_redef(ncid)
394 call nf95_def_var(ncid, 'RATQS', nf90_float, idim2, varid)
395 call nf95_put_att(ncid, varid, 'title', 'Ratqs')
396 call nf95_enddef(ncid)
397 call nf95_put_var(ncid, varid, ratqs(:, 1))
398
399 call nf95_redef(ncid)
400 call nf95_def_var(ncid, 'RUNOFFLIC0', nf90_float, idim2, &
401 varid_run_off_lic_0)
402 call nf95_put_att(ncid, varid_run_off_lic_0, 'title', 'Runofflic0')
403
404 call nf95_def_var(ncid, 'sig1', nf90_float, (/idim2, idim3/), varid_sig1)
405 call nf95_put_att(ncid, varid_sig1, 'long_name', &
406 'section adiabatic updraft')
407
408 call nf95_def_var(ncid, 'w01', nf90_float, (/idim2, idim3/), varid_w01)
409 call nf95_put_att(ncid, varid_w01, 'long_name', &
410 'vertical velocity within adiabatic updraft')
411
412 call nf95_enddef(ncid)
413
414 call nf95_put_var(ncid, varid_run_off_lic_0, run_off_lic_0)
415 call nf95_put_var(ncid, varid_sig1, sig1)
416 call nf95_put_var(ncid, varid_w01, w01)
417
418 call nf95_close(ncid)
419
420 END SUBROUTINE phyredem
421
422 end module phyredem_m

  ViewVC Help
Powered by ViewVC 1.1.21