/[lmdze]/trunk/phylmd/phyetat0.f
ViewVC logotype

Contents of /trunk/phylmd/phyetat0.f

Parent Directory Parent Directory | Revision Log Revision Log


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

  ViewVC Help
Powered by ViewVC 1.1.21