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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 15 - (hide annotations)
Fri Aug 1 15:24:12 2008 UTC (15 years, 9 months ago) by guez
Original Path: trunk/libf/phylmd/phyredem.f90
File size: 14021 byte(s)
-- Minor modification of input/output:

Added variable "Sigma_O3_Royer" to "histday.nc". "ecrit_day" is not
modified in "physiq". Removed variables "pyu1", "pyv1", "ftsol1",
"ftsol2", "ftsol3", "ftsol4", "psrf1", "psrf2", "psrf3", "psrf4"
"mfu", "mfd", "en_u", "en_d", "de_d", "de_u", "coefh" from
"histrac.nc".

Variable "raz_date" of module "conf_gcm_m" has logical type instead of
integer type.

-- Should not change any result at run time:

Modified calls to "IOIPSL_Lionel" procedures because the interfaces of
these procedures have been simplified.

Changed name of variable in module "start_init_orog_m": "masque" to
"mask".

Created a module containing procedure "phyredem".

Removed arguments "punjours", "pdayref" and "ptimestep" of procedure
"iniphysiq".

Renamed procedure "gr_phy_write" to "gr_phy_write_2d". Created
procedure "gr_phy_write_3d".

Removed procedures "ini_undefstd", "moy_undefSTD", "calcul_STDlev",
"calcul_divers".

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

  ViewVC Help
Powered by ViewVC 1.1.21