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

Contents of /trunk/phylmd/phyredem.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 3 - (show annotations)
Wed Feb 27 13:16:39 2008 UTC (16 years, 2 months ago) by guez
Original Path: trunk/libf/phylmd/phyredem.f
File size: 16282 byte(s)
Initial import
1 !
2 ! $Header: /home/cvsroot/LMDZ4/libf/phylmd/phyredem.F,v 1.3 2005/05/25 13:10:09 fairhead Exp $
3 !
4 c
5 SUBROUTINE phyredem (fichnom,dtime,radpas,
6 . rlat,rlon, pctsrf,tsol,tsoil,
7 cIM "slab" ocean
8 . tslab,seaice,
9 . qsurf,qsol,snow,
10 . albedo, alblw, evap, rain_fall, snow_fall,
11 . solsw, sollw,fder,
12 . radsol,frugs,agesno,
13 . zmea,zstd,zsig,zgam,zthe,zpic,zval,rugsrel,
14 . t_ancien, q_ancien, rnebcon, ratqs, clwcon,
15 . run_off_lic_0)
16 use dimens_m
17 use indicesol
18 use dimphy
19 use conf_gcm_m
20 use dimsoil
21 use temps
22 use clesphys
23 IMPLICIT none
24 c======================================================================
25 c Auteur(s) Z.X. Li (LMD/CNRS) date: 19930818
26 c Objet: Ecriture de l'etat de redemarrage pour la physique
27 c======================================================================
28 include "netcdf.inc"
29 c======================================================================
30 CHARACTER*(*) fichnom
31 REAL dtime
32 INTEGER radpas
33 REAL, intent(in):: rlat(klon), rlon(klon)
34 REAL tsol(klon,nbsrf)
35 REAL tsoil(klon,nsoilmx,nbsrf)
36 cIM "slab" ocean
37 REAL tslab(klon), seaice(klon)
38 REAL qsurf(klon,nbsrf)
39 REAL qsol(klon)
40 REAL snow(klon,nbsrf)
41 REAL albedo(klon,nbsrf)
42 cIM BEG
43 REAL alblw(klon,nbsrf)
44 cIM END
45 REAL evap(klon,nbsrf)
46 REAL rain_fall(klon)
47 REAL snow_fall(klon)
48 real solsw(klon)
49 real sollw(klon)
50 real fder(klon)
51 REAL radsol(klon)
52 REAL frugs(klon,nbsrf)
53 REAL agesno(klon,nbsrf)
54 REAL zmea(klon)
55 REAL zstd(klon)
56 REAL zsig(klon)
57 REAL zgam(klon)
58 REAL zthe(klon)
59 REAL zpic(klon)
60 REAL zval(klon)
61 REAL rugsrel(klon)
62 REAL pctsrf(klon, nbsrf)
63 REAL t_ancien(klon,klev), q_ancien(klon,klev)
64 real clwcon(klon,klev),rnebcon(klon,klev),ratqs(klon,klev)
65 REAL run_off_lic_0(klon)
66 c
67 INTEGER nid, nvarid, idim1, idim2, idim3
68 INTEGER ierr
69 INTEGER length
70 PARAMETER (length=100)
71 REAL tab_cntrl(length)
72 c
73 INTEGER isoil, nsrf
74 CHARACTER*7 str7
75 CHARACTER*2 str2
76 c
77 print *, "Call sequence information: phyredem"
78 ierr = NF_CREATE(fichnom, NF_CLOBBER, nid)
79 IF (ierr.NE.NF_NOERR) THEN
80 write(6,*)' Pb d''ouverture du fichier '//fichnom
81 write(6,*)' ierr = ', ierr
82 STOP 1
83 ENDIF
84 c
85 ierr = NF_PUT_ATT_TEXT (nid, NF_GLOBAL, "title", 28,
86 . "Fichier redemmarage physique")
87 c
88 ierr = NF_DEF_DIM (nid, "index", length, idim1)
89 ierr = NF_DEF_DIM (nid, "points_physiques", klon, idim2)
90 ierr = NF_DEF_DIM (nid, "horizon_vertical", klon*klev, idim3)
91 c
92 ierr = NF_ENDDEF(nid)
93 c
94 DO ierr = 1, length
95 tab_cntrl(ierr) = 0.0
96 ENDDO
97 tab_cntrl(1) = dtime
98 tab_cntrl(2) = radpas
99 tab_cntrl(3) = co2_ppm
100 tab_cntrl(4) = solaire
101 tab_cntrl(5) = iflag_con
102 tab_cntrl(6) = nbapp_rad
103
104 IF( cycle_diurne ) tab_cntrl( 7 ) = 1.
105 IF( soil_model ) tab_cntrl( 8 ) = 1.
106 IF( new_oliq ) tab_cntrl( 9 ) = 1.
107 IF( ok_orodr ) tab_cntrl(10 ) = 1.
108 IF( ok_orolf ) tab_cntrl(11 ) = 1.
109
110 tab_cntrl(13) = day_end
111 tab_cntrl(14) = annee_ref
112 tab_cntrl(15) = itau_phy
113 c
114 ierr = NF_REDEF (nid)
115 ierr = NF_DEF_VAR (nid, "controle", NF_FLOAT, 1, idim1,nvarid)
116 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 22,
117 . "Parametres de controle")
118 ierr = NF_ENDDEF(nid)
119 ierr = NF_PUT_VAR_REAL (nid,nvarid,tab_cntrl)
120 c
121 ierr = NF_REDEF (nid)
122 ierr = NF_DEF_VAR (nid, "longitude", NF_FLOAT, 1, idim2,nvarid)
123 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 32,
124 . "Longitudes de la grille physique")
125 ierr = NF_ENDDEF(nid)
126 ierr = NF_PUT_VAR_REAL (nid,nvarid,rlon)
127 c
128 ierr = NF_REDEF (nid)
129 ierr = NF_DEF_VAR (nid, "latitude", NF_FLOAT, 1, idim2,nvarid)
130 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 31,
131 . "Latitudes de la grille physique")
132 ierr = NF_ENDDEF(nid)
133 ierr = NF_PUT_VAR_REAL (nid,nvarid,rlat)
134 c
135 C PB ajout du masque terre/mer
136 C
137 ierr = NF_REDEF (nid)
138 ierr = NF_DEF_VAR (nid, "masque", NF_FLOAT, 1, idim2,nvarid)
139 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 16,
140 . "masque terre mer")
141 ierr = NF_ENDDEF(nid)
142 ierr = NF_PUT_VAR_REAL (nid,nvarid,zmasq)
143 c BP ajout des fraction de chaque sous-surface
144 C
145 C 1. fraction de terre
146 C
147 ierr = NF_REDEF (nid)
148 ierr = NF_DEF_VAR (nid, "FTER", NF_FLOAT, 1, idim2,nvarid)
149 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 21,
150 . "fraction de continent")
151 ierr = NF_ENDDEF(nid)
152 ierr = NF_PUT_VAR_REAL (nid,nvarid,pctsrf(1 : klon, is_ter))
153 C
154 C 2. Fraction de glace de terre
155 C
156 ierr = NF_REDEF (nid)
157 ierr = NF_DEF_VAR (nid, "FLIC", NF_FLOAT, 1, idim2,nvarid)
158 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 24,
159 . "fraction glace de terre")
160 ierr = NF_ENDDEF(nid)
161 ierr = NF_PUT_VAR_REAL (nid,nvarid,pctsrf(1 : klon, is_lic))
162 C
163 C 3. fraction ocean
164 C
165 ierr = NF_REDEF (nid)
166 ierr = NF_DEF_VAR (nid, "FOCE", NF_FLOAT, 1, idim2,nvarid)
167 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 14,
168 . "fraction ocean")
169 ierr = NF_ENDDEF(nid)
170 ierr = NF_PUT_VAR_REAL (nid,nvarid,pctsrf(1 : klon, is_oce))
171 C
172 C 4. Fraction glace de mer
173 C
174 ierr = NF_REDEF (nid)
175 ierr = NF_DEF_VAR (nid, "FSIC", NF_FLOAT, 1, idim2,nvarid)
176 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 18,
177 . "fraction glace mer")
178 ierr = NF_ENDDEF(nid)
179 ierr = NF_PUT_VAR_REAL (nid,nvarid,pctsrf(1 : klon, is_sic))
180 C
181 C
182 c
183 DO nsrf = 1, nbsrf
184 IF (nsrf.LE.99) THEN
185 WRITE(str2,'(i2.2)') nsrf
186 ierr = NF_REDEF (nid)
187 ierr = NF_DEF_VAR (nid, "TS"//str2, NF_FLOAT, 1, idim2,nvarid)
188 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 28,
189 . "Temperature de surface No."//str2)
190 ierr = NF_ENDDEF(nid)
191 ELSE
192 PRINT*, "Trop de sous-mailles"
193 stop 1
194 ENDIF
195 ierr = NF_PUT_VAR_REAL (nid,nvarid,tsol(1,nsrf))
196 ENDDO
197 c
198 DO nsrf = 1, nbsrf
199 DO isoil=1, nsoilmx
200 IF (isoil.LE.99 .AND. nsrf.LE.99) THEN
201 WRITE(str7,'(i2.2,"srf",i2.2)') isoil,nsrf
202 ierr = NF_REDEF (nid)
203 ierr = NF_DEF_VAR (nid, "Tsoil"//str7,NF_FLOAT,1,idim2,nvarid)
204 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 29,
205 . "Temperature du sol No."//str7)
206 ierr = NF_ENDDEF(nid)
207 ELSE
208 PRINT*, "Trop de couches"
209 stop 1
210 ENDIF
211 ierr = NF_PUT_VAR_REAL (nid,nvarid,tsoil(1,isoil,nsrf))
212 ENDDO
213 ENDDO
214 c
215 cIM "slab" ocean
216 ierr = NF_REDEF (nid)
217 ierr = NF_DEF_VAR (nid, "TSLAB", NF_FLOAT, 1, idim2,nvarid)
218 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 33,
219 . "Ecart de la SST (pour slab-ocean)")
220 ierr = NF_ENDDEF(nid)
221 ierr = NF_PUT_VAR_REAL (nid,nvarid,tslab)
222 c
223 ierr = NF_REDEF (nid)
224 ierr = NF_DEF_VAR (nid, "SEAICE", NF_FLOAT, 1, idim2,nvarid)
225 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 33,
226 . "Glace de mer kg/m2 (pour slab-ocean)")
227 ierr = NF_ENDDEF(nid)
228 ierr = NF_PUT_VAR_REAL (nid,nvarid,seaice)
229 c
230 DO nsrf = 1, nbsrf
231 IF (nsrf.LE.99) THEN
232 WRITE(str2,'(i2.2)') nsrf
233 ierr = NF_REDEF (nid)
234 ierr = NF_DEF_VAR (nid,"QS"//str2,NF_FLOAT,1,idim2,nvarid)
235 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 25,
236 . "Humidite de surface No."//str2)
237 ierr = NF_ENDDEF(nid)
238 ELSE
239 PRINT*, "Trop de sous-mailles"
240 stop 1
241 ENDIF
242 ierr = NF_PUT_VAR_REAL (nid,nvarid,qsurf(1,nsrf))
243 END DO
244 C
245 ierr = NF_REDEF (nid)
246 ierr = NF_DEF_VAR (nid,"QSOL",NF_FLOAT,1,idim2,nvarid)
247 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 20,
248 . "Eau dans le sol (mm)")
249 ierr = NF_ENDDEF(nid)
250 ierr = NF_PUT_VAR_REAL (nid,nvarid,qsol)
251 c
252 DO nsrf = 1, nbsrf
253 IF (nsrf.LE.99) THEN
254 WRITE(str2,'(i2.2)') nsrf
255 ierr = NF_REDEF (nid)
256 ierr = NF_DEF_VAR (nid,"ALBE"//str2,NF_FLOAT,1,idim2,nvarid)
257 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 23,
258 . "albedo de surface No."//str2)
259 ierr = NF_ENDDEF(nid)
260 ELSE
261 PRINT*, "Trop de sous-mailles"
262 stop 1
263 ENDIF
264 ierr = NF_PUT_VAR_REAL (nid,nvarid,albedo(1,nsrf))
265 ENDDO
266
267 cIM BEG albedo LW
268 DO nsrf = 1, nbsrf
269 IF (nsrf.LE.99) THEN
270 WRITE(str2,'(i2.2)') nsrf
271 ierr = NF_REDEF (nid)
272 ierr = NF_DEF_VAR (nid,"ALBLW"//str2,NF_FLOAT,1,idim2,nvarid)
273 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 23,
274 . "albedo LW de surface No."//str2)
275 ierr = NF_ENDDEF(nid)
276 ELSE
277 PRINT*, "Trop de sous-mailles"
278 stop 1
279 ENDIF
280 ierr = NF_PUT_VAR_REAL (nid,nvarid,alblw(1,nsrf))
281 ENDDO
282 cIM END albedo LW
283 c
284 DO nsrf = 1, nbsrf
285 IF (nsrf.LE.99) THEN
286 WRITE(str2,'(i2.2)') nsrf
287 ierr = NF_REDEF (nid)
288 ierr = NF_DEF_VAR (nid,"EVAP"//str2,NF_FLOAT,1,idim2,nvarid)
289 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 28,
290 . "Evaporation de surface No."//str2)
291 ierr = NF_ENDDEF(nid)
292 ELSE
293 PRINT*, "Trop de sous-mailles"
294 stop 1
295 ENDIF
296 ierr = NF_PUT_VAR_REAL (nid,nvarid,evap(1,nsrf))
297 ENDDO
298
299 c
300 DO nsrf = 1, nbsrf
301 IF (nsrf.LE.99) THEN
302 WRITE(str2,'(i2.2)') nsrf
303 ierr = NF_REDEF (nid)
304 ierr = NF_DEF_VAR (nid,"SNOW"//str2,NF_FLOAT,1,idim2,nvarid)
305 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 22,
306 . "Neige de surface No."//str2)
307 ierr = NF_ENDDEF(nid)
308 ELSE
309 PRINT*, "Trop de sous-mailles"
310 stop 1
311 ENDIF
312 ierr = NF_PUT_VAR_REAL (nid,nvarid,snow(1,nsrf))
313 ENDDO
314
315 c
316 ierr = NF_REDEF (nid)
317 ierr = NF_DEF_VAR (nid, "RADS", NF_FLOAT, 1, idim2,nvarid)
318 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 28,
319 . "Rayonnement net a la surface")
320 ierr = NF_ENDDEF(nid)
321 ierr = NF_PUT_VAR_REAL (nid,nvarid,radsol)
322 c
323 ierr = NF_REDEF (nid)
324 ierr = NF_DEF_VAR (nid, "solsw", NF_FLOAT, 1, idim2,nvarid)
325 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 32,
326 . "Rayonnement solaire a la surface")
327 ierr = NF_ENDDEF(nid)
328 ierr = NF_PUT_VAR_REAL (nid,nvarid,solsw)
329 c
330 ierr = NF_REDEF (nid)
331 ierr = NF_DEF_VAR (nid, "sollw", NF_FLOAT, 1, idim2,nvarid)
332 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 27,
333 . "Rayonnement IF a la surface")
334 ierr = NF_ENDDEF(nid)
335 ierr = NF_PUT_VAR_REAL (nid,nvarid,sollw)
336 c
337 ierr = NF_REDEF (nid)
338 ierr = NF_DEF_VAR (nid, "fder", NF_FLOAT, 1, idim2,nvarid)
339 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 14,
340 . "Derive de flux")
341 ierr = NF_ENDDEF(nid)
342 ierr = NF_PUT_VAR_REAL (nid,nvarid,fder)
343 c
344 ierr = NF_REDEF (nid)
345 ierr = NF_DEF_VAR (nid, "rain_f", NF_FLOAT, 1, idim2,nvarid)
346 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 21,
347 . "precipitation liquide")
348 ierr = NF_ENDDEF(nid)
349 ierr = NF_PUT_VAR_REAL (nid,nvarid,rain_fall)
350 c
351 ierr = NF_REDEF (nid)
352 ierr = NF_DEF_VAR (nid, "snow_f", NF_FLOAT, 1, idim2,nvarid)
353 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 20,
354 . "precipitation solide")
355 ierr = NF_ENDDEF(nid)
356 ierr = NF_PUT_VAR_REAL (nid,nvarid,snow_fall)
357 c
358 DO nsrf = 1, nbsrf
359 IF (nsrf.LE.99) THEN
360 WRITE(str2,'(i2.2)') nsrf
361 ierr = NF_REDEF (nid)
362 ierr = NF_DEF_VAR (nid,"RUG"//str2,NF_FLOAT,1,idim2,nvarid)
363 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 23,
364 . "rugosite de surface No."//str2)
365 ierr = NF_ENDDEF(nid)
366 ELSE
367 PRINT*, "Trop de sous-mailles"
368 stop 1
369 ENDIF
370 ierr = NF_PUT_VAR_REAL (nid,nvarid,frugs(1,nsrf))
371 ENDDO
372 c
373 DO nsrf = 1, nbsrf
374 IF (nsrf.LE.99) THEN
375 WRITE(str2,'(i2.2)') nsrf
376 ierr = NF_REDEF (nid)
377 ierr = NF_DEF_VAR (nid,"AGESNO"//str2,NF_FLOAT,1,idim2
378 $ ,nvarid)
379 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 15,
380 . "Age de la neige surface No."//str2)
381 ierr = NF_ENDDEF(nid)
382 ELSE
383 PRINT*, "Trop de sous-mailles"
384 stop 1
385 ENDIF
386 ierr = NF_PUT_VAR_REAL (nid,nvarid,agesno(1,nsrf))
387 ENDDO
388 c
389 ierr = NF_REDEF (nid)
390 ierr = NF_DEF_VAR (nid, "ZMEA", NF_FLOAT, 1, idim2,nvarid)
391 ierr = NF_ENDDEF(nid)
392 ierr = NF_PUT_VAR_REAL (nid,nvarid,zmea)
393 c
394 ierr = NF_REDEF (nid)
395 ierr = NF_DEF_VAR (nid, "ZSTD", NF_FLOAT, 1, idim2,nvarid)
396 ierr = NF_ENDDEF(nid)
397 ierr = NF_PUT_VAR_REAL (nid,nvarid,zstd)
398 ierr = NF_REDEF (nid)
399 ierr = NF_DEF_VAR (nid, "ZSIG", NF_FLOAT, 1, idim2,nvarid)
400 ierr = NF_ENDDEF(nid)
401 ierr = NF_PUT_VAR_REAL (nid,nvarid,zsig)
402 ierr = NF_REDEF (nid)
403 ierr = NF_DEF_VAR (nid, "ZGAM", NF_FLOAT, 1, idim2,nvarid)
404 ierr = NF_ENDDEF(nid)
405 ierr = NF_PUT_VAR_REAL (nid,nvarid,zgam)
406 ierr = NF_REDEF (nid)
407 ierr = NF_DEF_VAR (nid, "ZTHE", NF_FLOAT, 1, idim2,nvarid)
408 ierr = NF_ENDDEF(nid)
409 ierr = NF_PUT_VAR_REAL (nid,nvarid,zthe)
410 ierr = NF_REDEF (nid)
411 ierr = NF_DEF_VAR (nid, "ZPIC", NF_FLOAT, 1, idim2,nvarid)
412 ierr = NF_ENDDEF(nid)
413 ierr = NF_PUT_VAR_REAL (nid,nvarid,zpic)
414 ierr = NF_REDEF (nid)
415 ierr = NF_DEF_VAR (nid, "ZVAL", NF_FLOAT, 1, idim2,nvarid)
416 ierr = NF_ENDDEF(nid)
417 ierr = NF_PUT_VAR_REAL (nid,nvarid,zval)
418 ierr = NF_REDEF (nid)
419 ierr = NF_DEF_VAR (nid, "RUGSREL", NF_FLOAT, 1, idim2,nvarid)
420 ierr = NF_ENDDEF(nid)
421 ierr = NF_PUT_VAR_REAL (nid,nvarid,rugsrel)
422 c
423 ierr = NF_REDEF (nid)
424 ierr = NF_DEF_VAR (nid, "TANCIEN", NF_FLOAT, 1, idim3,nvarid)
425 ierr = NF_ENDDEF(nid)
426 ierr = NF_PUT_VAR_REAL (nid,nvarid,t_ancien)
427 c
428 ierr = NF_REDEF (nid)
429 ierr = NF_DEF_VAR (nid, "QANCIEN", NF_FLOAT, 1, idim3,nvarid)
430 ierr = NF_ENDDEF(nid)
431 ierr = NF_PUT_VAR_REAL (nid,nvarid,q_ancien)
432 c
433 ierr = NF_REDEF (nid)
434 ierr = NF_DEF_VAR (nid, "RUGMER", NF_FLOAT, 1, idim2,nvarid)
435 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 28,
436 . "Longueur de rugosite sur mer")
437 ierr = NF_ENDDEF(nid)
438 ierr = NF_PUT_VAR_REAL (nid,nvarid,frugs(1,is_oce))
439 c
440 ierr = NF_REDEF (nid)
441 ierr = NF_DEF_VAR (nid, "CLWCON", NF_FLOAT, 1, idim2,nvarid)
442 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 28,
443 . "Eau liquide convective")
444 ierr = NF_ENDDEF(nid)
445 ierr = NF_PUT_VAR_REAL (nid,nvarid,clwcon)
446 c
447 ierr = NF_REDEF (nid)
448 ierr = NF_DEF_VAR (nid, "RNEBCON", NF_FLOAT, 1, idim2,nvarid)
449 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 28,
450 . "Nebulosite convective")
451 ierr = NF_ENDDEF(nid)
452 ierr = NF_PUT_VAR_REAL (nid,nvarid,rnebcon)
453 c
454 ierr = NF_REDEF (nid)
455 ierr = NF_DEF_VAR (nid, "RATQS", NF_FLOAT, 1, idim2,nvarid)
456 ierr = NF_PUT_ATT_TEXT(nid,nvarid,"title", 5,
457 . "Ratqs")
458 ierr = NF_ENDDEF(nid)
459 ierr = NF_PUT_VAR_REAL (nid,nvarid,ratqs)
460 c
461 c run_off_lic_0
462 c
463 ierr = NF_REDEF (nid)
464 ierr=NF_DEF_VAR(nid,"RUNOFFLIC0",NF_FLOAT, 1,idim2,nvarid)
465 ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 10,
466 . "Runofflic0")
467 ierr = NF_ENDDEF(nid)
468 ierr = NF_PUT_VAR_REAL (nid,nvarid,run_off_lic_0)
469 c
470 c
471 ierr = NF_CLOSE(nid)
472 c
473 RETURN
474 END

  ViewVC Help
Powered by ViewVC 1.1.21