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

Contents of /trunk/Sources/phylmd/phyetat0.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 157 - (show annotations)
Mon Jul 20 16:01:49 2015 UTC (8 years, 10 months ago) by guez
File size: 11176 byte(s)
Just encapsulated SUBROUTINE vlsplt in a module and cleaned it.

In procedure vlx, local variables dxqu and adxqu only need indices
iip2:ip1jm. Otherwise, just cleaned vlx.

Procedures dynredem0 and dynredem1 no longer have argument fichnom,
they just operate on a file named "restart.nc". The programming
guideline here is that gcm should not be more complex than it needs by
itself, other programs (ce0l etc.) just have to adapt to gcm. So ce0l
now creates files "restart.nc" and "restartphy.nc".

In order to facilitate decentralizing the writing of "restartphy.nc",
created a procedure phyredem0 out of phyredem. phyredem0 creates the
NetCDF header of "restartphy.nc" while phyredem writes the NetCDF
variables. As the global attribute itau_phy needs to be filled in
phyredem0, at the beginnig of the run, we must compute its value
instead of just using itap. So we have a dummy argument lmt_pas of
phyredem0. Also, the ncid of "startphy.nc" is upgraded from local
variable of phyetat0 to dummy argument. phyetat0 no longer closes
"startphy.nc".

Following the same decentralizing objective, the ncid of "restart.nc"
is upgraded from local variable of dynredem0 to module variable of
dynredem0_m. "restart.nc" is not closed at the end of dynredem0 nor
opened at the beginning of dynredem1.

In procedure etat0, instead of creating many vectors of size klon
which will be filled with zeroes, just create one array null_array.

In procedure phytrac, instead of writing trs(: 1) to a text file,
write it to "restartphy.nc" (following LMDZ). This is better because
now trs(: 1) is next to its coordinates. We can write to
"restartphy.nc" from phytrac directly, and not add trs(: 1) to the
long list of variables in physiq, thanks to the decentralizing of
"restartphy.nc".

In procedure phyetat0, we no longer write to standard output the
minimum and maximum values of read arrays. It is ok to check input and
abort on invalid values but just printing statistics on input seems too
much useless computation and out of place clutter.

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

  ViewVC Help
Powered by ViewVC 1.1.21