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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 140 - (hide annotations)
Fri Jun 5 18:58:06 2015 UTC (8 years, 11 months ago) by guez
File size: 13820 byte(s)
Changed unit of variables lat_min_guide and lat_max_guide from module
conf_guide_m from degrees to rad. Then we do not have to convert the
whole array rlat from rad to degrees in SUBROUTINE tau2alpha.

Removed some useless computations in inigeom.

Removed module coefils. Moved variables sddv, unsddv, sddu, unsddu,
eignfnu, eignfnv of module coefils to module inifgn_m. Downgraded
variables coefilu, coefilu2, coefilv, coefilv2, modfrstu, modfrstv of
module coefils to local variables of SUBROUTINE inifilr.

Write and read a 3-dimensional variable Tsoil in restartphy.nc and
startphy.nc instead of multiple variables for the different
subs-urfaces and soil layers. This does not allow any longer to
provide only the surface value in startphy.nc and spread it to other
layers. Instead, if necessary, pre-process the file startphy.nc to
spread the surface value.

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

  ViewVC Help
Powered by ViewVC 1.1.21