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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 138 - (show annotations)
Fri May 22 23:13:19 2015 UTC (9 years ago) by guez
File size: 14159 byte(s)
Moved variable nb_files from module histcom_var to module
histbeg_totreg_m.

Removed unused argument q of writehist.

No history file is created in program ce0l so there is no need to call
histclo in etat0.

In phyredem, access variables rlat and rlon directly from module
phyetat0_m instead of having them as arguments. This is clearer for
the program gcm. There are bad side effects for the program ce0l: we
have to modify the module variables rlat and rlon in procedure etat0,
and we need the additional file phyetat0.f to compile ce0l.

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

  ViewVC Help
Powered by ViewVC 1.1.21