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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 54 - (hide annotations)
Tue Dec 6 15:07:04 2011 UTC (12 years, 5 months ago) by guez
Original Path: trunk/libf/phylmd/phyredem.f90
File size: 13599 byte(s)
Removed Numerical Recipes procedure "ran1". Replaced calls to "ran1"
in "inidissip" by calls to intrinsic procedures.

Split file "interface_surf.f90" into a file with a module containing
only variables, "interface_surf", and single-procedure files. Gathered
files into directory "Interface_surf".

Added argument "cdivu" to "gradiv" and "gradiv2", "cdivh" to
"divgrad2" and "divgrad", and "crot" to "nxgraro2" and
"nxgrarot". "dissip" now uses variables "cdivu", "cdivh" and "crot"
from module "inidissip_m", so it can pass them to "gradiv2",
etc. Thanks to this modification, we avoid a circular dependency
betwwen "inidissip.f90" and "gradiv2.f90", etc. The value -1. used by
"gradiv2", for instance, during computation of eigenvalues is not the
value "cdivu" computed by "inidissip".

Extracted procedure "start_inter_3d" from module "startdyn", to its
own module.

In "inidissip", unrolled loop on "ii". I find it clearer now.

Moved variables "matriceun", "matriceus", "matricevn", "matricevs",
"matrinvn" and "matrinvs" from module "parafilt" to module
"inifilr_m". Moved variables "jfiltnu", "jfiltnv", "jfiltsu",
"jfiltsv" from module "coefils" to module "inifilr_m".

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

  ViewVC Help
Powered by ViewVC 1.1.21