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