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

Annotation of /trunk/phylmd/phyredem.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 155 - (hide annotations)
Wed Jul 8 17:03:45 2015 UTC (8 years, 10 months ago) by guez
Original Path: trunk/Sources/phylmd/phyredem.f
File size: 13248 byte(s)
Do not write any longer to startphy.nc nor read from restartphy.nc the
NetCDF variable ALBLW: it was the same than ALBE. ALBE was for the
visible and ALBLW for the near infrared. In physiq, use only variables
falbe and albsol, removed falblw and albsollw. See revision 888 of
LMDZ.

Removed unused arguments pdp of SUBROUTINE lwbv, ptave of SUBROUTINE
lwv, kuaer of SUBROUTINE lwvd, nq of SUBROUTINE initphysto.

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

  ViewVC Help
Powered by ViewVC 1.1.21