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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 138 - (hide annotations)
Fri May 22 23:13:19 2015 UTC (9 years ago) by guez
File size: 14159 byte(s)
Moved variable nb_files from module histcom_var to module
histbeg_totreg_m.

Removed unused argument q of writehist.

No history file is created in program ce0l so there is no need to call
histclo in etat0.

In phyredem, access variables rlat and rlon directly from module
phyetat0_m instead of having them as arguments. This is clearer for
the program gcm. There are bad side effects for the program ce0l: we
have to modify the module variables rlat and rlon in procedure etat0,
and we need the additional file phyetat0.f to compile ce0l.

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

  ViewVC Help
Powered by ViewVC 1.1.21