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