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

Contents of /trunk/phylmd/phyredem.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 101 - (show annotations)
Mon Jul 7 17:45:21 2014 UTC (9 years, 10 months ago) by guez
File size: 14161 byte(s)
Removed unused files "interfoce_slab.f" and "gath2cpl.f". Removed
unused variables coastalflow and riverflow of module
interface_surf. Removed unused arguments cal, radsol, dif_grnd,
fluxlat, fluxsens, dflux_s, dflux_l of procedure fonte_neige. Removed
unused arguments tslab, seaice of procedure interfsurf_hq and
clqh. Removed unused arguments seaice of procedure clmain.

In interfsurf_hq, used variable soil_model of module clesphys2 instead
of cascading it as an argument from physiq.

In phyetat0, stop if masque not found.

Variable TS instead of "TS[0-9][0-9]" in "(re)startphy.nc", with
additional dimension nbsrf.

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

  ViewVC Help
Powered by ViewVC 1.1.21