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