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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 19 - (show annotations)
Thu Aug 7 15:46:20 2008 UTC (15 years, 9 months ago) by guez
Original Path: trunk/libf/phylmd/phyredem.f90
File size: 14032 byte(s)
Inlined procedures "regr_pr_av" and "regr_pr_int" in "regr_pr_o3",
"regr_pr_av_coefoz" and "regr_pr_int_coefoz".

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

  ViewVC Help
Powered by ViewVC 1.1.21