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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 99 - (hide annotations)
Wed Jul 2 18:39:15 2014 UTC (9 years, 10 months ago) by guez
Original Path: trunk/phylmd/phyredem.f
File size: 14278 byte(s)
Created procedure test_disvert (following LMDZ). Added procedures
hybrid and funcd in module disvert_m. Upgraded compute_ab from
internal procedure of disvert to module procedure. Added variables y,
ya in module disvert_m. Upgraded s from local variable of procedure
disvert to module variable.

Renamed allowed value of variable vert_sampling in procedure disvert
from "read" to "read_hybrid". Added possibility to read pressure
values, value "read_pressure". Replaced vertical distribution for
value "param" by the distribution "strato_correct" from LMDZ (but kept
the value "param"). In case "tropo", replaced 1 by dsigmin (following
LMDZ). In case "strato", replaced 0.3 by dsigmin (following LMDZ).

Changed computation of bp in procedure compute_ab.

Removed debugindex case in clmain. Removed useless argument rlon of
procedure clmain. Removed useless variables ytaux, ytauy of procedure
clmain.

Removed intermediary variables tsol, qsol, tsolsrf, tslab in procedure
etat0.

Removed variable ok_veget:. coupling with the model Orchid is not
possible. Removed variable ocean: modeling an ocean slab is not
possible.

Removed useless variables tmp_rriv and tmp_rcoa from module
interface_surf.

Moved initialization of variables da, mp, phi in procedure physiq to
to inside the test iflag_con >= 3.

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

  ViewVC Help
Powered by ViewVC 1.1.21