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

Annotation of /trunk/Sources/phylmd/phyredem.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 72 - (hide annotations)
Tue Jul 23 13:00:07 2013 UTC (10 years, 10 months ago) by guez
Original Path: trunk/libf/phylmd/phyredem.f90
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 guez 15 module phyredem_m
2 guez 3
3 guez 15 IMPLICIT NONE
4 guez 3
5 guez 15 contains
6 guez 3
7 guez 15 SUBROUTINE phyredem(fichnom, rlat, rlon, pctsrf, tsol, tsoil, tslab, &
8 guez 62 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 guez 72 run_off_lic_0, sig1, w01)
12 guez 3
13 guez 72 ! 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 guez 12
18 guez 72 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 guez 21 nf95_def_var, nf95_enddef, nf95_redef, nf95_put_var, nf95_close
24 guez 72 USE temps, ONLY: itau_phy
25 guez 12
26 guez 15 CHARACTER(len=*) fichnom
27 guez 62 REAL, INTENT(IN):: rlat(klon), rlon(klon)
28 guez 72 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 guez 43 REAL, intent(in):: qsol(klon)
34 guez 72 REAL snow(klon, nbsrf)
35     REAL albedo(klon, nbsrf)
36     REAL alblw(klon, nbsrf)
37     REAL evap(klon, nbsrf)
38 guez 62 REAL, INTENT(IN):: rain_fall(klon)
39 guez 72 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 guez 15 REAL, intent(in):: zstd(klon)
48     REAL, intent(in):: zsig(klon)
49 guez 72 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 guez 12
58 guez 72 real, intent(in):: w01(klon, klev)
59     ! vertical velocity within adiabatic updraft
60 guez 12
61 guez 72 ! Local:
62 guez 12
63 guez 72 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 guez 15 !------------------------------------------------------------
70 guez 12
71 guez 15 PRINT *, 'Call sequence information: phyredem'
72 guez 72 CALL nf95_create(fichnom, nf90_clobber, ncid)
73 guez 12
74 guez 72 call nf95_put_att(ncid, nf90_global, 'title', &
75 guez 19 'Fichier redémarrage physique')
76 guez 72 call nf95_put_att(ncid, nf90_global, "itau_phy", itau_phy)
77 guez 12
78 guez 72 call nf95_def_dim(ncid, 'points_physiques', klon, idim2)
79     call nf95_def_dim(ncid, 'klev', klev, idim3)
80 guez 12
81 guez 72 call nf95_def_var(ncid, 'longitude', nf90_float, idim2, varid)
82     call nf95_put_att(ncid, varid, 'title', &
83 guez 15 'Longitudes de la grille physique')
84 guez 72 call nf95_enddef(ncid)
85     call nf95_put_var(ncid, varid, rlon)
86 guez 12
87 guez 72 call nf95_redef(ncid)
88     call nf95_def_var(ncid, 'latitude', nf90_float, idim2, varid)
89     call nf95_put_att(ncid, varid, 'title', &
90 guez 15 'Latitudes de la grille physique')
91 guez 72 call nf95_enddef(ncid)
92     call nf95_put_var(ncid, varid, rlat)
93 guez 12
94 guez 15 ! PB ajout du masque terre/mer
95 guez 12
96 guez 72 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 guez 15 ! BP ajout des fraction de chaque sous-surface
102 guez 12
103 guez 15 ! 1. fraction de terre
104 guez 12
105 guez 72 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 guez 12
111 guez 15 ! 2. Fraction de glace de terre
112 guez 12
113 guez 72 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 guez 12
119 guez 15 ! 3. fraction ocean
120 guez 12
121 guez 72 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 guez 12
127 guez 15 ! 4. Fraction glace de mer
128 guez 3
129 guez 72 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 guez 12
135 guez 15 DO nsrf = 1, nbsrf
136     IF (nsrf<=99) THEN
137     WRITE (str2, '(i2.2)') nsrf
138 guez 72 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 guez 15 'Temperature de surface No.'//str2)
142 guez 72 call nf95_enddef(ncid)
143 guez 15 ELSE
144     PRINT *, 'Trop de sous-mailles'
145     STOP 1
146     END IF
147 guez 72 call nf95_put_var(ncid, varid, tsol(:, nsrf))
148 guez 15 END DO
149 guez 12
150 guez 15 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 guez 72 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 guez 15 'Temperature du sol No.'//str7)
158 guez 72 call nf95_enddef(ncid)
159 guez 15 ELSE
160     PRINT *, 'Trop de couches'
161     STOP 1
162     END IF
163 guez 72 call nf95_put_var(ncid, varid, tsoil(:, isoil, nsrf))
164 guez 15 END DO
165     END DO
166 guez 12
167 guez 15 !IM "slab" ocean
168 guez 72 call nf95_redef(ncid)
169     call nf95_def_var(ncid, 'TSLAB', nf90_float, idim2, varid)
170     call nf95_put_att(ncid, varid, 'title', &
171 guez 15 'Ecart de la SST (pour slab-ocean)')
172 guez 72 call nf95_enddef(ncid)
173     call nf95_put_var(ncid, varid, tslab)
174 guez 12
175 guez 72 call nf95_redef(ncid)
176     call nf95_def_var(ncid, 'SEAICE', nf90_float, idim2, varid)
177     call nf95_put_att(ncid, varid, 'title', &
178 guez 15 'Glace de mer kg/m2 (pour slab-ocean)')
179 guez 72 call nf95_enddef(ncid)
180     call nf95_put_var(ncid, varid, seaice)
181 guez 12
182 guez 15 DO nsrf = 1, nbsrf
183     IF (nsrf<=99) THEN
184     WRITE (str2, '(i2.2)') nsrf
185 guez 72 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 guez 15 'Humidite de surface No.'//str2)
189 guez 72 call nf95_enddef(ncid)
190 guez 15 ELSE
191     PRINT *, 'Trop de sous-mailles'
192     STOP 1
193     END IF
194 guez 72 call nf95_put_var(ncid, varid, qsurf(:, nsrf))
195 guez 15 END DO
196 guez 12
197 guez 72 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 guez 12
203 guez 15 DO nsrf = 1, nbsrf
204     IF (nsrf<=99) THEN
205     WRITE (str2, '(i2.2)') nsrf
206 guez 72 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 guez 15 'albedo de surface No.'//str2)
210 guez 72 call nf95_enddef(ncid)
211 guez 15 ELSE
212     PRINT *, 'Trop de sous-mailles'
213     STOP 1
214     END IF
215 guez 72 call nf95_put_var(ncid, varid, albedo(:, nsrf))
216 guez 15 END DO
217 guez 12
218 guez 15 !IM BEG albedo LW
219     DO nsrf = 1, nbsrf
220     IF (nsrf<=99) THEN
221     WRITE (str2, '(i2.2)') nsrf
222 guez 72 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 guez 15 'albedo LW de surface No.'//str2)
226 guez 72 call nf95_enddef(ncid)
227 guez 15 ELSE
228     PRINT *, 'Trop de sous-mailles'
229     STOP 1
230     END IF
231 guez 72 call nf95_put_var(ncid, varid, alblw(:, nsrf))
232 guez 15 END DO
233     !IM END albedo LW
234 guez 12
235 guez 15 DO nsrf = 1, nbsrf
236     IF (nsrf<=99) THEN
237     WRITE (str2, '(i2.2)') nsrf
238 guez 72 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 guez 15 'Evaporation de surface No.'//str2)
242 guez 72 call nf95_enddef(ncid)
243 guez 15 ELSE
244     PRINT *, 'Trop de sous-mailles'
245     STOP 1
246     END IF
247 guez 72 call nf95_put_var(ncid, varid, evap(:, nsrf))
248 guez 15 END DO
249 guez 12
250 guez 15 DO nsrf = 1, nbsrf
251     IF (nsrf<=99) THEN
252     WRITE (str2, '(i2.2)') nsrf
253 guez 72 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 guez 15 'Neige de surface No.'//str2)
257 guez 72 call nf95_enddef(ncid)
258 guez 15 ELSE
259     PRINT *, 'Trop de sous-mailles'
260     STOP 1
261     END IF
262 guez 72 call nf95_put_var(ncid, varid, snow(:, nsrf))
263 guez 15 END DO
264 guez 12
265 guez 72 call nf95_redef(ncid)
266     call nf95_def_var(ncid, 'RADS', nf90_float, idim2, varid)
267     call nf95_put_att(ncid, varid, 'title', &
268 guez 15 'Rayonnement net a la surface')
269 guez 72 call nf95_enddef(ncid)
270     call nf95_put_var(ncid, varid, radsol)
271 guez 12
272 guez 72 call nf95_redef(ncid)
273     call nf95_def_var(ncid, 'solsw', nf90_float, idim2, varid)
274     call nf95_put_att(ncid, varid, 'title', &
275 guez 15 'Rayonnement solaire a la surface')
276 guez 72 call nf95_enddef(ncid)
277     call nf95_put_var(ncid, varid, solsw)
278 guez 12
279 guez 72 call nf95_redef(ncid)
280     call nf95_def_var(ncid, 'sollw', nf90_float, idim2, varid)
281     call nf95_put_att(ncid, varid, 'title', &
282 guez 15 'Rayonnement IF a la surface')
283 guez 72 call nf95_enddef(ncid)
284     call nf95_put_var(ncid, varid, sollw)
285 guez 12
286 guez 72 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 guez 12
292 guez 72 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 guez 12
298 guez 72 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 guez 12
304 guez 15 DO nsrf = 1, nbsrf
305     IF (nsrf<=99) THEN
306     WRITE (str2, '(i2.2)') nsrf
307 guez 72 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 guez 15 'rugosite de surface No.'//str2)
311 guez 72 call nf95_enddef(ncid)
312 guez 15 ELSE
313     PRINT *, 'Trop de sous-mailles'
314     STOP 1
315     END IF
316 guez 72 call nf95_put_var(ncid, varid, frugs(:, nsrf))
317 guez 15 END DO
318 guez 12
319 guez 15 DO nsrf = 1, nbsrf
320     IF (nsrf<=99) THEN
321     WRITE (str2, '(i2.2)') nsrf
322 guez 72 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 guez 15 'Age de la neige surface No.'//str2)
326 guez 72 call nf95_enddef(ncid)
327 guez 15 ELSE
328     PRINT *, 'Trop de sous-mailles'
329     STOP 1
330     END IF
331 guez 72 call nf95_put_var(ncid, varid, agesno(:, nsrf))
332 guez 15 END DO
333 guez 12
334 guez 72 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 guez 12
339 guez 72 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 guez 12
364 guez 72 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 guez 12
369 guez 72 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 guez 12
374 guez 72 call nf95_redef(ncid)
375     call nf95_def_var(ncid, 'RUGMER', nf90_float, idim2, varid)
376     call nf95_put_att(ncid, varid, 'title', &
377 guez 15 'Longueur de rugosite sur mer')
378 guez 72 call nf95_enddef(ncid)
379     call nf95_put_var(ncid, varid, frugs(:, is_oce))
380 guez 12
381 guez 72 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 guez 12
387 guez 72 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 guez 12
393 guez 72 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 guez 15
399 guez 72 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 guez 15
404 guez 72 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 guez 15
408 guez 72 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 guez 15
412 guez 72 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 guez 15 END SUBROUTINE phyredem
421    
422     end module phyredem_m

  ViewVC Help
Powered by ViewVC 1.1.21