1 |
module phyetat0_m |
2 |
|
3 |
use dimphy, only: klon |
4 |
|
5 |
IMPLICIT none |
6 |
|
7 |
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 |
|
11 |
private klon |
12 |
|
13 |
contains |
14 |
|
15 |
SUBROUTINE phyetat0(pctsrf, tsol, tsoil, qsurf, qsol, snow, albe, evap, & |
16 |
rain_fall, snow_fall, solsw, sollw, fder, radsol, frugs, agesno, zmea, & |
17 |
zstd, zsig, zgam, zthe, zpic, zval, t_ancien, q_ancien, ancien_ok, & |
18 |
rnebcon, ratqs, clwcon, run_off_lic_0, sig1, w01, ncid_startphy, & |
19 |
itau_phy) |
20 |
|
21 |
! From phylmd/phyetat0.F, version 1.4 2005/06/03 10:03:07 |
22 |
! Author: Z.X. Li (LMD/CNRS) |
23 |
! Date: 1993/08/18 |
24 |
! Objet : lecture de l'état initial pour la physique |
25 |
|
26 |
use dimphy, only: zmasq, klev |
27 |
USE dimsoil, ONLY : nsoilmx |
28 |
USE indicesol, ONLY : epsfra, is_lic, is_oce, is_sic, is_ter, nbsrf |
29 |
use netcdf, only: nf90_global, nf90_inq_varid, NF90_NOERR, NF90_NOWRITE |
30 |
use netcdf95, only: nf95_close, nf95_get_att, nf95_get_var, & |
31 |
nf95_inq_varid, nf95_inquire_variable, NF95_OPEN |
32 |
|
33 |
REAL, intent(out):: pctsrf(klon, nbsrf) |
34 |
REAL, intent(out):: tsol(klon, nbsrf) |
35 |
REAL, intent(out):: tsoil(klon, nsoilmx, nbsrf) |
36 |
REAL, intent(out):: qsurf(klon, nbsrf) |
37 |
REAL, intent(out):: qsol(:) ! (klon) |
38 |
REAL, intent(out):: snow(klon, nbsrf) |
39 |
REAL, intent(out):: albe(klon, nbsrf) |
40 |
REAL, intent(out):: evap(klon, nbsrf) |
41 |
REAL, intent(out):: rain_fall(klon) |
42 |
REAL, intent(out):: snow_fall(klon) |
43 |
real, intent(out):: solsw(klon) |
44 |
REAL, intent(out):: sollw(klon) |
45 |
real, intent(out):: fder(klon) |
46 |
REAL, intent(out):: radsol(klon) |
47 |
REAL, intent(out):: frugs(klon, nbsrf) |
48 |
REAL, intent(out):: agesno(klon, nbsrf) |
49 |
REAL, intent(out):: zmea(klon) |
50 |
REAL, intent(out):: zstd(klon) |
51 |
REAL, intent(out):: zsig(klon) |
52 |
REAL, intent(out):: zgam(klon) |
53 |
REAL, intent(out):: zthe(klon) |
54 |
REAL, intent(out):: zpic(klon) |
55 |
REAL, intent(out):: zval(klon) |
56 |
REAL, intent(out):: t_ancien(klon, klev), q_ancien(klon, klev) |
57 |
LOGICAL, intent(out):: ancien_ok |
58 |
real, intent(out):: rnebcon(klon, klev), ratqs(klon, klev) |
59 |
REAL, intent(out):: clwcon(klon, klev), run_off_lic_0(klon) |
60 |
real, intent(out):: sig1(klon, klev) ! section adiabatic updraft |
61 |
|
62 |
real, intent(out):: w01(klon, klev) |
63 |
! vertical velocity within adiabatic updraft |
64 |
|
65 |
integer, intent(out):: ncid_startphy, itau_phy |
66 |
|
67 |
! Local: |
68 |
REAL fractint(klon) |
69 |
INTEGER varid, ndims |
70 |
INTEGER ierr, i |
71 |
|
72 |
!--------------------------------------------------------------- |
73 |
|
74 |
print *, "Call sequence information: phyetat0" |
75 |
|
76 |
! Fichier contenant l'état initial : |
77 |
call NF95_OPEN("startphy.nc", NF90_NOWRITE, ncid_startphy) |
78 |
|
79 |
call nf95_get_att(ncid_startphy, nf90_global, "itau_phy", itau_phy) |
80 |
|
81 |
! Lecture des latitudes (coordonnees): |
82 |
|
83 |
call NF95_INQ_VARID(ncid_startphy, "latitude", varid) |
84 |
call NF95_GET_VAR(ncid_startphy, varid, rlat) |
85 |
|
86 |
! Lecture des longitudes (coordonnees): |
87 |
|
88 |
call NF95_INQ_VARID(ncid_startphy, "longitude", varid) |
89 |
call NF95_GET_VAR(ncid_startphy, varid, rlon) |
90 |
|
91 |
! Lecture du masque terre mer |
92 |
|
93 |
call NF95_INQ_VARID(ncid_startphy, "masque", varid) |
94 |
call nf95_get_var(ncid_startphy, varid, zmasq) |
95 |
|
96 |
! Lecture des fractions pour chaque sous-surface |
97 |
|
98 |
! initialisation des sous-surfaces |
99 |
|
100 |
pctsrf = 0. |
101 |
|
102 |
! fraction de terre |
103 |
|
104 |
ierr = NF90_INQ_VARID(ncid_startphy, "FTER", varid) |
105 |
IF (ierr == NF90_NOERR) THEN |
106 |
call nf95_get_var(ncid_startphy, varid, pctsrf(:, is_ter)) |
107 |
else |
108 |
PRINT *, 'phyetat0: Le champ <FTER> est absent' |
109 |
ENDIF |
110 |
|
111 |
! fraction de glace de terre |
112 |
|
113 |
ierr = NF90_INQ_VARID(ncid_startphy, "FLIC", varid) |
114 |
IF (ierr == NF90_NOERR) THEN |
115 |
call nf95_get_var(ncid_startphy, varid, pctsrf(:, is_lic)) |
116 |
else |
117 |
PRINT *, 'phyetat0: Le champ <FLIC> est absent' |
118 |
ENDIF |
119 |
|
120 |
! fraction d'ocean |
121 |
|
122 |
ierr = NF90_INQ_VARID(ncid_startphy, "FOCE", varid) |
123 |
IF (ierr == NF90_NOERR) THEN |
124 |
call nf95_get_var(ncid_startphy, varid, pctsrf(:, is_oce)) |
125 |
else |
126 |
PRINT *, 'phyetat0: Le champ <FOCE> est absent' |
127 |
ENDIF |
128 |
|
129 |
! fraction glace de mer |
130 |
|
131 |
ierr = NF90_INQ_VARID(ncid_startphy, "FSIC", varid) |
132 |
IF (ierr == NF90_NOERR) THEN |
133 |
call nf95_get_var(ncid_startphy, varid, pctsrf(:, is_sic)) |
134 |
else |
135 |
PRINT *, 'phyetat0: Le champ <FSIC> est absent' |
136 |
ENDIF |
137 |
|
138 |
! Verification de l'adequation entre le masque et les sous-surfaces |
139 |
|
140 |
fractint = pctsrf(:, is_ter) + pctsrf(:, is_lic) |
141 |
DO i = 1 , klon |
142 |
IF ( abs(fractint(i) - zmasq(i) ) > EPSFRA ) THEN |
143 |
WRITE(*, *) 'phyetat0: attention fraction terre pas ', & |
144 |
'coherente ', i, zmasq(i), pctsrf(i, is_ter) & |
145 |
, pctsrf(i, is_lic) |
146 |
ENDIF |
147 |
END DO |
148 |
fractint = pctsrf(:, is_oce) + pctsrf(:, is_sic) |
149 |
DO i = 1 , klon |
150 |
IF ( abs( fractint(i) - (1. - zmasq(i))) > EPSFRA ) THEN |
151 |
WRITE(*, *) 'phyetat0 attention fraction ocean pas ', & |
152 |
'coherente ', i, zmasq(i) , pctsrf(i, is_oce) & |
153 |
, pctsrf(i, is_sic) |
154 |
ENDIF |
155 |
END DO |
156 |
|
157 |
! Lecture des temperatures du sol: |
158 |
call NF95_INQ_VARID(ncid_startphy, "TS", varid) |
159 |
call nf95_inquire_variable(ncid_startphy, varid, ndims = ndims) |
160 |
if (ndims == 2) then |
161 |
call NF95_GET_VAR(ncid_startphy, varid, tsol) |
162 |
else |
163 |
print *, "Found only one surface type for soil temperature." |
164 |
call nf95_get_var(ncid_startphy, varid, tsol(:, 1)) |
165 |
tsol(:, 2:nbsrf) = spread(tsol(:, 1), dim = 2, ncopies = nbsrf - 1) |
166 |
end if |
167 |
|
168 |
! Lecture des temperatures du sol profond: |
169 |
|
170 |
call NF95_INQ_VARID(ncid_startphy, 'Tsoil', varid) |
171 |
call NF95_GET_VAR(ncid_startphy, varid, tsoil) |
172 |
|
173 |
! Lecture de l'humidite de l'air juste au dessus du sol: |
174 |
|
175 |
call NF95_INQ_VARID(ncid_startphy, "QS", varid) |
176 |
call nf95_get_var(ncid_startphy, varid, qsurf) |
177 |
|
178 |
! Eau dans le sol (pour le modele de sol "bucket") |
179 |
|
180 |
ierr = NF90_INQ_VARID(ncid_startphy, "QSOL", varid) |
181 |
IF (ierr == NF90_NOERR) THEN |
182 |
call nf95_get_var(ncid_startphy, varid, qsol) |
183 |
else |
184 |
PRINT *, 'phyetat0: Le champ <QSOL> est absent' |
185 |
PRINT *, ' Valeur par defaut nulle' |
186 |
qsol = 0. |
187 |
ENDIF |
188 |
|
189 |
! Lecture de neige au sol: |
190 |
|
191 |
call NF95_INQ_VARID(ncid_startphy, "SNOW", varid) |
192 |
call nf95_get_var(ncid_startphy, varid, snow) |
193 |
|
194 |
! Lecture de albedo au sol: |
195 |
|
196 |
call NF95_INQ_VARID(ncid_startphy, "ALBE", varid) |
197 |
call nf95_get_var(ncid_startphy, varid, albe) |
198 |
|
199 |
! Lecture de evaporation: |
200 |
|
201 |
call NF95_INQ_VARID(ncid_startphy, "EVAP", varid) |
202 |
call nf95_get_var(ncid_startphy, varid, evap) |
203 |
|
204 |
! Lecture precipitation liquide: |
205 |
|
206 |
call NF95_INQ_VARID(ncid_startphy, "rain_f", varid) |
207 |
call NF95_GET_VAR(ncid_startphy, varid, rain_fall) |
208 |
|
209 |
! Lecture precipitation solide: |
210 |
|
211 |
call NF95_INQ_VARID(ncid_startphy, "snow_f", varid) |
212 |
call NF95_GET_VAR(ncid_startphy, varid, snow_fall) |
213 |
|
214 |
! Lecture rayonnement solaire au sol: |
215 |
|
216 |
ierr = NF90_INQ_VARID(ncid_startphy, "solsw", varid) |
217 |
IF (ierr /= NF90_NOERR) THEN |
218 |
PRINT *, 'phyetat0: Le champ <solsw> est absent' |
219 |
PRINT *, 'mis a zero' |
220 |
solsw = 0. |
221 |
ELSE |
222 |
call nf95_get_var(ncid_startphy, varid, solsw) |
223 |
ENDIF |
224 |
|
225 |
! Lecture rayonnement IF au sol: |
226 |
|
227 |
ierr = NF90_INQ_VARID(ncid_startphy, "sollw", varid) |
228 |
IF (ierr /= NF90_NOERR) THEN |
229 |
PRINT *, 'phyetat0: Le champ <sollw> est absent' |
230 |
PRINT *, 'mis a zero' |
231 |
sollw = 0. |
232 |
ELSE |
233 |
call nf95_get_var(ncid_startphy, varid, sollw) |
234 |
ENDIF |
235 |
|
236 |
! Lecture derive des flux: |
237 |
|
238 |
ierr = NF90_INQ_VARID(ncid_startphy, "fder", varid) |
239 |
IF (ierr /= NF90_NOERR) THEN |
240 |
PRINT *, 'phyetat0: Le champ <fder> est absent' |
241 |
PRINT *, 'mis a zero' |
242 |
fder = 0. |
243 |
ELSE |
244 |
call nf95_get_var(ncid_startphy, varid, fder) |
245 |
ENDIF |
246 |
|
247 |
! Lecture du rayonnement net au sol: |
248 |
|
249 |
call NF95_INQ_VARID(ncid_startphy, "RADS", varid) |
250 |
call NF95_GET_VAR(ncid_startphy, varid, radsol) |
251 |
|
252 |
! Lecture de la longueur de rugosite |
253 |
|
254 |
call NF95_INQ_VARID(ncid_startphy, "RUG", varid) |
255 |
call nf95_get_var(ncid_startphy, varid, frugs) |
256 |
|
257 |
! Lecture de l'age de la neige: |
258 |
|
259 |
call NF95_INQ_VARID(ncid_startphy, "AGESNO", varid) |
260 |
call nf95_get_var(ncid_startphy, varid, agesno) |
261 |
|
262 |
call NF95_INQ_VARID(ncid_startphy, "ZMEA", varid) |
263 |
call NF95_GET_VAR(ncid_startphy, varid, zmea) |
264 |
|
265 |
call NF95_INQ_VARID(ncid_startphy, "ZSTD", varid) |
266 |
call NF95_GET_VAR(ncid_startphy, varid, zstd) |
267 |
|
268 |
call NF95_INQ_VARID(ncid_startphy, "ZSIG", varid) |
269 |
call NF95_GET_VAR(ncid_startphy, varid, zsig) |
270 |
|
271 |
call NF95_INQ_VARID(ncid_startphy, "ZGAM", varid) |
272 |
call NF95_GET_VAR(ncid_startphy, varid, zgam) |
273 |
|
274 |
call NF95_INQ_VARID(ncid_startphy, "ZTHE", varid) |
275 |
call NF95_GET_VAR(ncid_startphy, varid, zthe) |
276 |
|
277 |
call NF95_INQ_VARID(ncid_startphy, "ZPIC", varid) |
278 |
call NF95_GET_VAR(ncid_startphy, varid, zpic) |
279 |
|
280 |
call NF95_INQ_VARID(ncid_startphy, "ZVAL", varid) |
281 |
call NF95_GET_VAR(ncid_startphy, varid, zval) |
282 |
|
283 |
ancien_ok = .TRUE. |
284 |
|
285 |
ierr = NF90_INQ_VARID(ncid_startphy, "TANCIEN", varid) |
286 |
IF (ierr /= NF90_NOERR) THEN |
287 |
PRINT *, "phyetat0: Le champ <TANCIEN> est absent" |
288 |
PRINT *, "Depart legerement fausse. Mais je continue" |
289 |
ancien_ok = .FALSE. |
290 |
ELSE |
291 |
call nf95_get_var(ncid_startphy, varid, t_ancien) |
292 |
ENDIF |
293 |
|
294 |
ierr = NF90_INQ_VARID(ncid_startphy, "QANCIEN", varid) |
295 |
IF (ierr /= NF90_NOERR) THEN |
296 |
PRINT *, "phyetat0: Le champ <QANCIEN> est absent" |
297 |
PRINT *, "Depart legerement fausse. Mais je continue" |
298 |
ancien_ok = .FALSE. |
299 |
ELSE |
300 |
call nf95_get_var(ncid_startphy, varid, q_ancien) |
301 |
ENDIF |
302 |
|
303 |
ierr = NF90_INQ_VARID(ncid_startphy, "CLWCON", varid) |
304 |
IF (ierr /= NF90_NOERR) THEN |
305 |
PRINT *, "phyetat0: Le champ CLWCON est absent" |
306 |
PRINT *, "Depart legerement fausse. Mais je continue" |
307 |
clwcon = 0. |
308 |
ELSE |
309 |
call nf95_get_var(ncid_startphy, varid, clwcon(:, 1)) |
310 |
clwcon(:, 2:) = 0. |
311 |
ENDIF |
312 |
|
313 |
ierr = NF90_INQ_VARID(ncid_startphy, "RNEBCON", varid) |
314 |
IF (ierr /= NF90_NOERR) THEN |
315 |
PRINT *, "phyetat0: Le champ RNEBCON est absent" |
316 |
PRINT *, "Depart legerement fausse. Mais je continue" |
317 |
rnebcon = 0. |
318 |
ELSE |
319 |
call nf95_get_var(ncid_startphy, varid, rnebcon(:, 1)) |
320 |
rnebcon(:, 2:) = 0. |
321 |
ENDIF |
322 |
|
323 |
! Lecture ratqs |
324 |
|
325 |
ierr = NF90_INQ_VARID(ncid_startphy, "RATQS", varid) |
326 |
IF (ierr /= NF90_NOERR) THEN |
327 |
PRINT *, "phyetat0: Le champ <RATQS> est absent" |
328 |
PRINT *, "Depart legerement fausse. Mais je continue" |
329 |
ratqs = 0. |
330 |
ELSE |
331 |
call nf95_get_var(ncid_startphy, varid, ratqs(:, 1)) |
332 |
ratqs(:, 2:) = 0. |
333 |
ENDIF |
334 |
|
335 |
! Lecture run_off_lic_0 |
336 |
|
337 |
ierr = NF90_INQ_VARID(ncid_startphy, "RUNOFFLIC0", varid) |
338 |
IF (ierr /= NF90_NOERR) THEN |
339 |
PRINT *, "phyetat0: Le champ <RUNOFFLIC0> est absent" |
340 |
PRINT *, "Depart legerement fausse. Mais je continue" |
341 |
run_off_lic_0 = 0. |
342 |
ELSE |
343 |
call nf95_get_var(ncid_startphy, varid, run_off_lic_0) |
344 |
ENDIF |
345 |
|
346 |
call nf95_inq_varid(ncid_startphy, "sig1", varid) |
347 |
call nf95_get_var(ncid_startphy, varid, sig1) |
348 |
|
349 |
call nf95_inq_varid(ncid_startphy, "w01", varid) |
350 |
call nf95_get_var(ncid_startphy, varid, w01) |
351 |
|
352 |
END SUBROUTINE phyetat0 |
353 |
|
354 |
end module phyetat0_m |