1 |
guez |
3 |
module phyetat0_m |
2 |
|
|
|
3 |
guez |
15 |
use dimphy, only: klon |
4 |
guez |
12 |
|
5 |
guez |
3 |
IMPLICIT none |
6 |
|
|
|
7 |
guez |
138 |
REAL, save:: rlat(klon), rlon(klon) |
8 |
|
|
! latitude and longitude of a point of the scalar grid identified |
9 |
|
|
! by a simple index, in degrees |
10 |
guez |
3 |
|
11 |
guez |
15 |
private klon |
12 |
guez |
3 |
|
13 |
|
|
contains |
14 |
|
|
|
15 |
guez |
99 |
SUBROUTINE phyetat0(pctsrf, tsol, tsoil, tslab, seaice, qsurf, qsol, & |
16 |
guez |
155 |
snow, albe, evap, rain_fall, snow_fall, solsw, sollw, fder, & |
17 |
guez |
99 |
radsol, frugs, agesno, zmea, zstd, zsig, zgam, zthe, zpic, zval, & |
18 |
|
|
t_ancien, q_ancien, ancien_ok, rnebcon, ratqs, clwcon, run_off_lic_0, & |
19 |
|
|
sig1, w01) |
20 |
guez |
3 |
|
21 |
guez |
49 |
! From phylmd/phyetat0.F, version 1.4 2005/06/03 10:03:07 |
22 |
|
|
! Author: Z.X. Li (LMD/CNRS) |
23 |
guez |
50 |
! Date: 1993/08/18 |
24 |
guez |
69 |
! Objet : lecture de l'état initial pour la physique |
25 |
guez |
3 |
|
26 |
guez |
69 |
use dimphy, only: zmasq, klev |
27 |
|
|
USE dimsoil, ONLY : nsoilmx |
28 |
guez |
12 |
USE indicesol, ONLY : epsfra, is_lic, is_oce, is_sic, is_ter, nbsrf |
29 |
guez |
101 |
use netcdf, only: nf90_global, nf90_inq_varid, NF90_NOERR, & |
30 |
guez |
50 |
NF90_NOWRITE |
31 |
guez |
101 |
use netcdf95, only: nf95_close, nf95_get_att, nf95_get_var, & |
32 |
|
|
nf95_inq_varid, nf95_inquire_variable, NF95_OPEN |
33 |
guez |
69 |
USE temps, ONLY : itau_phy |
34 |
guez |
12 |
|
35 |
guez |
72 |
REAL pctsrf(klon, nbsrf) |
36 |
guez |
49 |
REAL tsol(klon, nbsrf) |
37 |
|
|
REAL tsoil(klon, nsoilmx, nbsrf) |
38 |
guez |
3 |
REAL tslab(klon), seaice(klon) |
39 |
guez |
49 |
REAL qsurf(klon, nbsrf) |
40 |
guez |
99 |
REAL, intent(out):: qsol(:) ! (klon) |
41 |
guez |
49 |
REAL snow(klon, nbsrf) |
42 |
|
|
REAL albe(klon, nbsrf) |
43 |
|
|
REAL evap(klon, nbsrf) |
44 |
guez |
62 |
REAL, intent(out):: rain_fall(klon) |
45 |
guez |
3 |
REAL snow_fall(klon) |
46 |
|
|
real solsw(klon) |
47 |
guez |
72 |
REAL, intent(out):: sollw(klon) |
48 |
guez |
3 |
real fder(klon) |
49 |
guez |
72 |
REAL radsol(klon) |
50 |
guez |
49 |
REAL frugs(klon, nbsrf) |
51 |
|
|
REAL agesno(klon, nbsrf) |
52 |
guez |
3 |
REAL zmea(klon) |
53 |
guez |
13 |
REAL, intent(out):: zstd(klon) |
54 |
|
|
REAL, intent(out):: zsig(klon) |
55 |
guez |
3 |
REAL zgam(klon) |
56 |
|
|
REAL zthe(klon) |
57 |
|
|
REAL zpic(klon) |
58 |
|
|
REAL zval(klon) |
59 |
guez |
49 |
REAL t_ancien(klon, klev), q_ancien(klon, klev) |
60 |
|
|
LOGICAL, intent(out):: ancien_ok |
61 |
guez |
72 |
real rnebcon(klon, klev), ratqs(klon, klev), clwcon(klon, klev) |
62 |
|
|
REAL run_off_lic_0(klon) |
63 |
|
|
real, intent(out):: sig1(klon, klev) ! section adiabatic updraft |
64 |
guez |
3 |
|
65 |
guez |
72 |
real, intent(out):: w01(klon, klev) |
66 |
|
|
! vertical velocity within adiabatic updraft |
67 |
guez |
3 |
|
68 |
guez |
72 |
! Local: |
69 |
|
|
REAL fractint(klon) |
70 |
guez |
3 |
REAL xmin, xmax |
71 |
guez |
101 |
INTEGER ncid, varid, ndims |
72 |
guez |
156 |
INTEGER ierr, i |
73 |
guez |
3 |
|
74 |
|
|
!--------------------------------------------------------------- |
75 |
|
|
|
76 |
|
|
print *, "Call sequence information: phyetat0" |
77 |
|
|
|
78 |
guez |
72 |
! Fichier contenant l'état initial : |
79 |
guez |
99 |
call NF95_OPEN("startphy.nc", NF90_NOWRITE, ncid) |
80 |
guez |
3 |
|
81 |
guez |
101 |
call nf95_get_att(ncid, nf90_global, "itau_phy", itau_phy) |
82 |
guez |
3 |
|
83 |
|
|
! Lecture des latitudes (coordonnees): |
84 |
|
|
|
85 |
guez |
50 |
call NF95_INQ_VARID(ncid, "latitude", varid) |
86 |
|
|
call NF95_GET_VAR(ncid, varid, rlat) |
87 |
guez |
3 |
|
88 |
|
|
! Lecture des longitudes (coordonnees): |
89 |
|
|
|
90 |
guez |
50 |
call NF95_INQ_VARID(ncid, "longitude", varid) |
91 |
|
|
call NF95_GET_VAR(ncid, varid, rlon) |
92 |
guez |
3 |
|
93 |
|
|
! Lecture du masque terre mer |
94 |
|
|
|
95 |
guez |
101 |
call NF95_INQ_VARID(ncid, "masque", varid) |
96 |
|
|
call nf95_get_var(ncid, varid, zmasq) |
97 |
|
|
|
98 |
guez |
3 |
! Lecture des fractions pour chaque sous-surface |
99 |
|
|
|
100 |
|
|
! initialisation des sous-surfaces |
101 |
|
|
|
102 |
|
|
pctsrf = 0. |
103 |
|
|
|
104 |
|
|
! fraction de terre |
105 |
|
|
|
106 |
guez |
49 |
ierr = NF90_INQ_VARID(ncid, "FTER", varid) |
107 |
guez |
50 |
IF (ierr == NF90_NOERR) THEN |
108 |
|
|
call nf95_get_var(ncid, varid, pctsrf(:, is_ter)) |
109 |
guez |
3 |
else |
110 |
guez |
43 |
PRINT *, 'phyetat0: Le champ <FTER> est absent' |
111 |
guez |
3 |
ENDIF |
112 |
|
|
|
113 |
|
|
! fraction de glace de terre |
114 |
|
|
|
115 |
guez |
49 |
ierr = NF90_INQ_VARID(ncid, "FLIC", varid) |
116 |
guez |
50 |
IF (ierr == NF90_NOERR) THEN |
117 |
|
|
call nf95_get_var(ncid, varid, pctsrf(:, is_lic)) |
118 |
guez |
3 |
else |
119 |
guez |
43 |
PRINT *, 'phyetat0: Le champ <FLIC> est absent' |
120 |
guez |
3 |
ENDIF |
121 |
|
|
|
122 |
|
|
! fraction d'ocean |
123 |
|
|
|
124 |
guez |
49 |
ierr = NF90_INQ_VARID(ncid, "FOCE", varid) |
125 |
guez |
50 |
IF (ierr == NF90_NOERR) THEN |
126 |
|
|
call nf95_get_var(ncid, varid, pctsrf(:, is_oce)) |
127 |
guez |
3 |
else |
128 |
guez |
43 |
PRINT *, 'phyetat0: Le champ <FOCE> est absent' |
129 |
guez |
3 |
ENDIF |
130 |
|
|
|
131 |
|
|
! fraction glace de mer |
132 |
|
|
|
133 |
guez |
49 |
ierr = NF90_INQ_VARID(ncid, "FSIC", varid) |
134 |
guez |
50 |
IF (ierr == NF90_NOERR) THEN |
135 |
|
|
call nf95_get_var(ncid, varid, pctsrf(:, is_sic)) |
136 |
guez |
3 |
else |
137 |
guez |
43 |
PRINT *, 'phyetat0: Le champ <FSIC> est absent' |
138 |
guez |
3 |
ENDIF |
139 |
|
|
|
140 |
guez |
50 |
! Verification de l'adequation entre le masque et les sous-surfaces |
141 |
guez |
3 |
|
142 |
guez |
50 |
fractint = pctsrf(:, is_ter) + pctsrf(:, is_lic) |
143 |
guez |
3 |
DO i = 1 , klon |
144 |
guez |
50 |
IF ( abs(fractint(i) - zmasq(i) ) > EPSFRA ) THEN |
145 |
|
|
WRITE(*, *) 'phyetat0: attention fraction terre pas ', & |
146 |
guez |
3 |
'coherente ', i, zmasq(i), pctsrf(i, is_ter) & |
147 |
guez |
49 |
, pctsrf(i, is_lic) |
148 |
guez |
3 |
ENDIF |
149 |
|
|
END DO |
150 |
guez |
50 |
fractint = pctsrf(:, is_oce) + pctsrf(:, is_sic) |
151 |
guez |
3 |
DO i = 1 , klon |
152 |
guez |
50 |
IF ( abs( fractint(i) - (1. - zmasq(i))) > EPSFRA ) THEN |
153 |
|
|
WRITE(*, *) 'phyetat0 attention fraction ocean pas ', & |
154 |
guez |
3 |
'coherente ', i, zmasq(i) , pctsrf(i, is_oce) & |
155 |
guez |
49 |
, pctsrf(i, is_sic) |
156 |
guez |
3 |
ENDIF |
157 |
|
|
END DO |
158 |
|
|
|
159 |
|
|
! Lecture des temperatures du sol: |
160 |
guez |
101 |
call NF95_INQ_VARID(ncid, "TS", varid) |
161 |
|
|
call nf95_inquire_variable(ncid, varid, ndims = ndims) |
162 |
|
|
if (ndims == 2) then |
163 |
|
|
call NF95_GET_VAR(ncid, varid, tsol) |
164 |
|
|
else |
165 |
|
|
print *, "Found only one surface type for soil temperature." |
166 |
guez |
49 |
call nf95_get_var(ncid, varid, tsol(:, 1)) |
167 |
guez |
101 |
tsol(:, 2:nbsrf) = spread(tsol(:, 1), dim = 2, ncopies = nbsrf - 1) |
168 |
guez |
156 |
end if |
169 |
guez |
3 |
|
170 |
guez |
156 |
! Lecture des temperatures du sol profond: |
171 |
guez |
3 |
|
172 |
guez |
140 |
call NF95_INQ_VARID(ncid, 'Tsoil', varid) |
173 |
|
|
call NF95_GET_VAR(ncid, varid, tsoil) |
174 |
guez |
3 |
|
175 |
|
|
!IM "slab" ocean |
176 |
guez |
50 |
! Lecture de tslab (pour slab ocean seulement): |
177 |
guez |
99 |
tslab = 0. |
178 |
|
|
seaice = 0. |
179 |
guez |
3 |
|
180 |
|
|
! Lecture de l'humidite de l'air juste au dessus du sol: |
181 |
|
|
|
182 |
guez |
156 |
call NF95_INQ_VARID(ncid, "QS", varid) |
183 |
|
|
call nf95_get_var(ncid, varid, qsurf) |
184 |
|
|
xmin = 1.0E+20 |
185 |
|
|
xmax = -1.0E+20 |
186 |
|
|
DO i = 1, klon |
187 |
|
|
xmin = MIN(qsurf(i, 1), xmin) |
188 |
|
|
xmax = MAX(qsurf(i, 1), xmax) |
189 |
|
|
ENDDO |
190 |
|
|
PRINT *, 'Humidite pres du sol <QS>', xmin, xmax |
191 |
guez |
3 |
|
192 |
|
|
! Eau dans le sol (pour le modele de sol "bucket") |
193 |
|
|
|
194 |
guez |
49 |
ierr = NF90_INQ_VARID(ncid, "QSOL", varid) |
195 |
guez |
50 |
IF (ierr == NF90_NOERR) THEN |
196 |
guez |
49 |
call nf95_get_var(ncid, varid, qsol) |
197 |
guez |
3 |
else |
198 |
guez |
43 |
PRINT *, 'phyetat0: Le champ <QSOL> est absent' |
199 |
guez |
50 |
PRINT *, ' Valeur par defaut nulle' |
200 |
guez |
43 |
qsol = 0. |
201 |
guez |
3 |
ENDIF |
202 |
|
|
|
203 |
|
|
! Lecture de neige au sol: |
204 |
|
|
|
205 |
guez |
156 |
call NF95_INQ_VARID(ncid, "SNOW", varid) |
206 |
|
|
call nf95_get_var(ncid, varid, snow) |
207 |
|
|
xmin = 1.0E+20 |
208 |
|
|
xmax = -1.0E+20 |
209 |
|
|
DO i = 1, klon |
210 |
|
|
xmin = MIN(snow(i, 1), xmin) |
211 |
|
|
xmax = MAX(snow(i, 1), xmax) |
212 |
|
|
ENDDO |
213 |
|
|
PRINT *, 'Neige du sol <SNOW>', xmin, xmax |
214 |
guez |
3 |
|
215 |
|
|
! Lecture de albedo au sol: |
216 |
|
|
|
217 |
guez |
156 |
call NF95_INQ_VARID(ncid, "ALBE", varid) |
218 |
|
|
call nf95_get_var(ncid, varid, albe) |
219 |
|
|
xmin = 1.0E+20 |
220 |
|
|
xmax = -1.0E+20 |
221 |
|
|
DO i = 1, klon |
222 |
|
|
xmin = MIN(albe(i, 1), xmin) |
223 |
|
|
xmax = MAX(albe(i, 1), xmax) |
224 |
|
|
ENDDO |
225 |
|
|
PRINT *, 'Neige du sol <ALBE>', xmin, xmax |
226 |
guez |
3 |
|
227 |
guez |
50 |
! Lecture de evaporation: |
228 |
guez |
3 |
|
229 |
guez |
156 |
call NF95_INQ_VARID(ncid, "EVAP", varid) |
230 |
|
|
call nf95_get_var(ncid, varid, evap) |
231 |
|
|
xmin = 1.0E+20 |
232 |
|
|
xmax = -1.0E+20 |
233 |
|
|
DO i = 1, klon |
234 |
|
|
xmin = MIN(evap(i, 1), xmin) |
235 |
|
|
xmax = MAX(evap(i, 1), xmax) |
236 |
|
|
ENDDO |
237 |
|
|
PRINT *, 'Evap du sol <EVAP>', xmin, xmax |
238 |
guez |
3 |
|
239 |
|
|
! Lecture precipitation liquide: |
240 |
|
|
|
241 |
guez |
50 |
call NF95_INQ_VARID(ncid, "rain_f", varid) |
242 |
|
|
call NF95_GET_VAR(ncid, varid, rain_fall) |
243 |
guez |
3 |
|
244 |
|
|
! Lecture precipitation solide: |
245 |
|
|
|
246 |
guez |
50 |
call NF95_INQ_VARID(ncid, "snow_f", varid) |
247 |
|
|
call NF95_GET_VAR(ncid, varid, snow_fall) |
248 |
guez |
3 |
xmin = 1.0E+20 |
249 |
|
|
xmax = -1.0E+20 |
250 |
|
|
DO i = 1, klon |
251 |
guez |
49 |
xmin = MIN(snow_fall(i), xmin) |
252 |
|
|
xmax = MAX(snow_fall(i), xmax) |
253 |
guez |
3 |
ENDDO |
254 |
guez |
49 |
PRINT *, 'Precipitation solide snow_f:', xmin, xmax |
255 |
guez |
3 |
|
256 |
|
|
! Lecture rayonnement solaire au sol: |
257 |
|
|
|
258 |
guez |
49 |
ierr = NF90_INQ_VARID(ncid, "solsw", varid) |
259 |
|
|
IF (ierr /= NF90_NOERR) THEN |
260 |
guez |
43 |
PRINT *, 'phyetat0: Le champ <solsw> est absent' |
261 |
|
|
PRINT *, 'mis a zero' |
262 |
guez |
3 |
solsw = 0. |
263 |
|
|
ELSE |
264 |
guez |
49 |
call nf95_get_var(ncid, varid, solsw) |
265 |
guez |
3 |
ENDIF |
266 |
|
|
xmin = 1.0E+20 |
267 |
|
|
xmax = -1.0E+20 |
268 |
|
|
DO i = 1, klon |
269 |
guez |
49 |
xmin = MIN(solsw(i), xmin) |
270 |
|
|
xmax = MAX(solsw(i), xmax) |
271 |
guez |
3 |
ENDDO |
272 |
guez |
49 |
PRINT *, 'Rayonnement solaire au sol solsw:', xmin, xmax |
273 |
guez |
3 |
|
274 |
|
|
! Lecture rayonnement IF au sol: |
275 |
|
|
|
276 |
guez |
49 |
ierr = NF90_INQ_VARID(ncid, "sollw", varid) |
277 |
|
|
IF (ierr /= NF90_NOERR) THEN |
278 |
guez |
43 |
PRINT *, 'phyetat0: Le champ <sollw> est absent' |
279 |
|
|
PRINT *, 'mis a zero' |
280 |
guez |
3 |
sollw = 0. |
281 |
|
|
ELSE |
282 |
guez |
49 |
call nf95_get_var(ncid, varid, sollw) |
283 |
guez |
3 |
ENDIF |
284 |
guez |
72 |
PRINT *, 'Rayonnement IF au sol sollw:', minval(sollw), maxval(sollw) |
285 |
guez |
3 |
|
286 |
|
|
! Lecture derive des flux: |
287 |
|
|
|
288 |
guez |
49 |
ierr = NF90_INQ_VARID(ncid, "fder", varid) |
289 |
|
|
IF (ierr /= NF90_NOERR) THEN |
290 |
guez |
43 |
PRINT *, 'phyetat0: Le champ <fder> est absent' |
291 |
|
|
PRINT *, 'mis a zero' |
292 |
guez |
3 |
fder = 0. |
293 |
|
|
ELSE |
294 |
guez |
49 |
call nf95_get_var(ncid, varid, fder) |
295 |
guez |
3 |
ENDIF |
296 |
|
|
xmin = 1.0E+20 |
297 |
|
|
xmax = -1.0E+20 |
298 |
|
|
DO i = 1, klon |
299 |
guez |
49 |
xmin = MIN(fder(i), xmin) |
300 |
|
|
xmax = MAX(fder(i), xmax) |
301 |
guez |
3 |
ENDDO |
302 |
guez |
49 |
PRINT *, 'Derive des flux fder:', xmin, xmax |
303 |
guez |
3 |
|
304 |
|
|
! Lecture du rayonnement net au sol: |
305 |
|
|
|
306 |
guez |
50 |
call NF95_INQ_VARID(ncid, "RADS", varid) |
307 |
|
|
call NF95_GET_VAR(ncid, varid, radsol) |
308 |
guez |
3 |
xmin = 1.0E+20 |
309 |
|
|
xmax = -1.0E+20 |
310 |
|
|
DO i = 1, klon |
311 |
guez |
49 |
xmin = MIN(radsol(i), xmin) |
312 |
|
|
xmax = MAX(radsol(i), xmax) |
313 |
guez |
3 |
ENDDO |
314 |
guez |
49 |
PRINT *, 'Rayonnement net au sol radsol:', xmin, xmax |
315 |
guez |
3 |
|
316 |
|
|
! Lecture de la longueur de rugosite |
317 |
|
|
|
318 |
guez |
156 |
call NF95_INQ_VARID(ncid, "RUG", varid) |
319 |
|
|
call nf95_get_var(ncid, varid, frugs) |
320 |
|
|
xmin = 1.0E+20 |
321 |
|
|
xmax = -1.0E+20 |
322 |
|
|
DO i = 1, klon |
323 |
|
|
xmin = MIN(frugs(i, 1), xmin) |
324 |
|
|
xmax = MAX(frugs(i, 1), xmax) |
325 |
|
|
ENDDO |
326 |
|
|
PRINT *, 'rugosite <RUG>', xmin, xmax |
327 |
guez |
3 |
|
328 |
|
|
! Lecture de l'age de la neige: |
329 |
|
|
|
330 |
guez |
156 |
call NF95_INQ_VARID(ncid, "AGESNO", varid) |
331 |
|
|
call nf95_get_var(ncid, varid, agesno) |
332 |
|
|
xmin = 1.0E+20 |
333 |
|
|
xmax = -1.0E+20 |
334 |
|
|
DO i = 1, klon |
335 |
|
|
xmin = MIN(agesno(i, 1), xmin) |
336 |
|
|
xmax = MAX(agesno(i, 1), xmax) |
337 |
|
|
ENDDO |
338 |
|
|
PRINT *, 'Age de la neige <AGESNO>', xmin, xmax |
339 |
guez |
3 |
|
340 |
guez |
50 |
call NF95_INQ_VARID(ncid, "ZMEA", varid) |
341 |
|
|
call NF95_GET_VAR(ncid, varid, zmea) |
342 |
guez |
3 |
xmin = 1.0E+20 |
343 |
|
|
xmax = -1.0E+20 |
344 |
|
|
DO i = 1, klon |
345 |
guez |
49 |
xmin = MIN(zmea(i), xmin) |
346 |
|
|
xmax = MAX(zmea(i), xmax) |
347 |
guez |
3 |
ENDDO |
348 |
guez |
49 |
PRINT *, 'OROGRAPHIE SOUS-MAILLE zmea:', xmin, xmax |
349 |
guez |
3 |
|
350 |
guez |
50 |
call NF95_INQ_VARID(ncid, "ZSTD", varid) |
351 |
|
|
call NF95_GET_VAR(ncid, varid, zstd) |
352 |
guez |
3 |
xmin = 1.0E+20 |
353 |
|
|
xmax = -1.0E+20 |
354 |
|
|
DO i = 1, klon |
355 |
guez |
49 |
xmin = MIN(zstd(i), xmin) |
356 |
|
|
xmax = MAX(zstd(i), xmax) |
357 |
guez |
3 |
ENDDO |
358 |
guez |
49 |
PRINT *, 'OROGRAPHIE SOUS-MAILLE zstd:', xmin, xmax |
359 |
guez |
3 |
|
360 |
guez |
50 |
call NF95_INQ_VARID(ncid, "ZSIG", varid) |
361 |
|
|
call NF95_GET_VAR(ncid, varid, zsig) |
362 |
guez |
3 |
xmin = 1.0E+20 |
363 |
|
|
xmax = -1.0E+20 |
364 |
|
|
DO i = 1, klon |
365 |
guez |
49 |
xmin = MIN(zsig(i), xmin) |
366 |
|
|
xmax = MAX(zsig(i), xmax) |
367 |
guez |
3 |
ENDDO |
368 |
guez |
49 |
PRINT *, 'OROGRAPHIE SOUS-MAILLE zsig:', xmin, xmax |
369 |
guez |
3 |
|
370 |
guez |
50 |
call NF95_INQ_VARID(ncid, "ZGAM", varid) |
371 |
|
|
call NF95_GET_VAR(ncid, varid, zgam) |
372 |
guez |
3 |
xmin = 1.0E+20 |
373 |
|
|
xmax = -1.0E+20 |
374 |
|
|
DO i = 1, klon |
375 |
guez |
49 |
xmin = MIN(zgam(i), xmin) |
376 |
|
|
xmax = MAX(zgam(i), xmax) |
377 |
guez |
3 |
ENDDO |
378 |
guez |
49 |
PRINT *, 'OROGRAPHIE SOUS-MAILLE zgam:', xmin, xmax |
379 |
guez |
3 |
|
380 |
guez |
50 |
call NF95_INQ_VARID(ncid, "ZTHE", varid) |
381 |
|
|
call NF95_GET_VAR(ncid, varid, zthe) |
382 |
guez |
3 |
xmin = 1.0E+20 |
383 |
|
|
xmax = -1.0E+20 |
384 |
|
|
DO i = 1, klon |
385 |
guez |
49 |
xmin = MIN(zthe(i), xmin) |
386 |
|
|
xmax = MAX(zthe(i), xmax) |
387 |
guez |
3 |
ENDDO |
388 |
guez |
49 |
PRINT *, 'OROGRAPHIE SOUS-MAILLE zthe:', xmin, xmax |
389 |
guez |
3 |
|
390 |
guez |
50 |
call NF95_INQ_VARID(ncid, "ZPIC", varid) |
391 |
|
|
call NF95_GET_VAR(ncid, varid, zpic) |
392 |
guez |
3 |
xmin = 1.0E+20 |
393 |
|
|
xmax = -1.0E+20 |
394 |
|
|
DO i = 1, klon |
395 |
guez |
49 |
xmin = MIN(zpic(i), xmin) |
396 |
|
|
xmax = MAX(zpic(i), xmax) |
397 |
guez |
3 |
ENDDO |
398 |
guez |
49 |
PRINT *, 'OROGRAPHIE SOUS-MAILLE zpic:', xmin, xmax |
399 |
guez |
3 |
|
400 |
guez |
50 |
call NF95_INQ_VARID(ncid, "ZVAL", varid) |
401 |
|
|
call NF95_GET_VAR(ncid, varid, zval) |
402 |
guez |
3 |
xmin = 1.0E+20 |
403 |
|
|
xmax = -1.0E+20 |
404 |
|
|
DO i = 1, klon |
405 |
guez |
49 |
xmin = MIN(zval(i), xmin) |
406 |
|
|
xmax = MAX(zval(i), xmax) |
407 |
guez |
3 |
ENDDO |
408 |
guez |
49 |
PRINT *, 'OROGRAPHIE SOUS-MAILLE zval:', xmin, xmax |
409 |
guez |
3 |
|
410 |
|
|
ancien_ok = .TRUE. |
411 |
|
|
|
412 |
guez |
49 |
ierr = NF90_INQ_VARID(ncid, "TANCIEN", varid) |
413 |
|
|
IF (ierr /= NF90_NOERR) THEN |
414 |
guez |
43 |
PRINT *, "phyetat0: Le champ <TANCIEN> est absent" |
415 |
|
|
PRINT *, "Depart legerement fausse. Mais je continue" |
416 |
guez |
3 |
ancien_ok = .FALSE. |
417 |
|
|
ELSE |
418 |
guez |
49 |
call nf95_get_var(ncid, varid, t_ancien) |
419 |
guez |
3 |
ENDIF |
420 |
|
|
|
421 |
guez |
49 |
ierr = NF90_INQ_VARID(ncid, "QANCIEN", varid) |
422 |
|
|
IF (ierr /= NF90_NOERR) THEN |
423 |
guez |
43 |
PRINT *, "phyetat0: Le champ <QANCIEN> est absent" |
424 |
|
|
PRINT *, "Depart legerement fausse. Mais je continue" |
425 |
guez |
3 |
ancien_ok = .FALSE. |
426 |
|
|
ELSE |
427 |
guez |
49 |
call nf95_get_var(ncid, varid, q_ancien) |
428 |
guez |
3 |
ENDIF |
429 |
|
|
|
430 |
guez |
49 |
ierr = NF90_INQ_VARID(ncid, "CLWCON", varid) |
431 |
|
|
IF (ierr /= NF90_NOERR) THEN |
432 |
guez |
43 |
PRINT *, "phyetat0: Le champ CLWCON est absent" |
433 |
|
|
PRINT *, "Depart legerement fausse. Mais je continue" |
434 |
guez |
3 |
clwcon = 0. |
435 |
|
|
ELSE |
436 |
guez |
72 |
call nf95_get_var(ncid, varid, clwcon(:, 1)) |
437 |
|
|
clwcon(:, 2:) = 0. |
438 |
guez |
3 |
ENDIF |
439 |
|
|
xmin = 1.0E+20 |
440 |
|
|
xmax = -1.0E+20 |
441 |
|
|
xmin = MINval(clwcon) |
442 |
|
|
xmax = MAXval(clwcon) |
443 |
guez |
49 |
PRINT *, 'Eau liquide convective (ecart-type) clwcon:', xmin, xmax |
444 |
guez |
3 |
|
445 |
guez |
49 |
ierr = NF90_INQ_VARID(ncid, "RNEBCON", varid) |
446 |
|
|
IF (ierr /= NF90_NOERR) THEN |
447 |
guez |
43 |
PRINT *, "phyetat0: Le champ RNEBCON est absent" |
448 |
|
|
PRINT *, "Depart legerement fausse. Mais je continue" |
449 |
guez |
3 |
rnebcon = 0. |
450 |
|
|
ELSE |
451 |
guez |
72 |
call nf95_get_var(ncid, varid, rnebcon(:, 1)) |
452 |
|
|
rnebcon(:, 2:) = 0. |
453 |
guez |
3 |
ENDIF |
454 |
|
|
xmin = 1.0E+20 |
455 |
|
|
xmax = -1.0E+20 |
456 |
|
|
xmin = MINval(rnebcon) |
457 |
|
|
xmax = MAXval(rnebcon) |
458 |
guez |
49 |
PRINT *, 'Nebulosite convective (ecart-type) rnebcon:', xmin, xmax |
459 |
guez |
3 |
|
460 |
|
|
! Lecture ratqs |
461 |
|
|
|
462 |
guez |
49 |
ierr = NF90_INQ_VARID(ncid, "RATQS", varid) |
463 |
|
|
IF (ierr /= NF90_NOERR) THEN |
464 |
guez |
43 |
PRINT *, "phyetat0: Le champ <RATQS> est absent" |
465 |
|
|
PRINT *, "Depart legerement fausse. Mais je continue" |
466 |
guez |
3 |
ratqs = 0. |
467 |
|
|
ELSE |
468 |
guez |
72 |
call nf95_get_var(ncid, varid, ratqs(:, 1)) |
469 |
|
|
ratqs(:, 2:) = 0. |
470 |
guez |
3 |
ENDIF |
471 |
|
|
xmin = 1.0E+20 |
472 |
|
|
xmax = -1.0E+20 |
473 |
|
|
xmin = MINval(ratqs) |
474 |
|
|
xmax = MAXval(ratqs) |
475 |
guez |
49 |
PRINT *, '(ecart-type) ratqs:', xmin, xmax |
476 |
guez |
3 |
|
477 |
|
|
! Lecture run_off_lic_0 |
478 |
|
|
|
479 |
guez |
49 |
ierr = NF90_INQ_VARID(ncid, "RUNOFFLIC0", varid) |
480 |
|
|
IF (ierr /= NF90_NOERR) THEN |
481 |
guez |
43 |
PRINT *, "phyetat0: Le champ <RUNOFFLIC0> est absent" |
482 |
|
|
PRINT *, "Depart legerement fausse. Mais je continue" |
483 |
guez |
3 |
run_off_lic_0 = 0. |
484 |
|
|
ELSE |
485 |
guez |
49 |
call nf95_get_var(ncid, varid, run_off_lic_0) |
486 |
guez |
3 |
ENDIF |
487 |
|
|
xmin = 1.0E+20 |
488 |
|
|
xmax = -1.0E+20 |
489 |
|
|
xmin = MINval(run_off_lic_0) |
490 |
|
|
xmax = MAXval(run_off_lic_0) |
491 |
guez |
49 |
PRINT *, '(ecart-type) run_off_lic_0:', xmin, xmax |
492 |
guez |
3 |
|
493 |
guez |
72 |
call nf95_inq_varid(ncid, "sig1", varid) |
494 |
|
|
call nf95_get_var(ncid, varid, sig1) |
495 |
|
|
|
496 |
|
|
call nf95_inq_varid(ncid, "w01", varid) |
497 |
|
|
call nf95_get_var(ncid, varid, w01) |
498 |
|
|
|
499 |
guez |
49 |
call NF95_CLOSE(ncid) |
500 |
guez |
3 |
|
501 |
|
|
END SUBROUTINE phyetat0 |
502 |
|
|
|
503 |
|
|
end module phyetat0_m |