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

Contents of /trunk/phylmd/phyetat0.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 304 - (show annotations)
Thu Sep 6 15:51:09 2018 UTC (5 years, 8 months ago) by guez
File size: 11997 byte(s)
Variable fevap of physiq is not used. Remove it from physiq and from
the restart file. Remove the corresponding argument evap of
pbl_surface.

Use directly yqsurf instead of qairsol in pbl_surface.

1 module phyetat0_m
2
3 use dimphy, only: klon
4
5 IMPLICIT none
6
7 REAL, save, protected:: 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, protected:: itau_phy
12 REAL, save, protected:: zmasq(KLON) ! fraction of land
13
14 private klon
15
16 contains
17
18 SUBROUTINE phyetat0(pctsrf, ftsol, ftsoil, qsurf, qsol, snow, albe, &
19 rain_fall, snow_fall, solsw, sollw, fder, radsol, frugs, agesno, zmea, &
20 zstd, zsig, zgam, zthe, zpic, zval, t_ancien, q_ancien, ancien_ok, &
21 rnebcon, ratqs, clwcon, run_off_lic_0, sig1, w01, ncid_startphy)
22
23 ! From phylmd/phyetat0.F, version 1.4 2005/06/03 10:03:07
24 ! Author: Z.X. Li (LMD/CNRS)
25 ! Date: 1993/08/18
26 ! Objet : lecture de l'état initial pour la physique
27
28 USE conf_gcm_m, ONLY: raz_date
29 use dimphy, only: klev
30 USE dimsoil, ONLY : nsoilmx
31 USE indicesol, ONLY : epsfra, is_lic, is_oce, is_sic, is_ter, nbsrf
32 use netcdf, only: nf90_global, nf90_inq_varid, NF90_NOERR, NF90_NOWRITE
33 use netcdf95, only: nf95_get_att, nf95_get_var, nf95_inq_varid, &
34 nf95_inquire_variable, NF95_OPEN
35
36 REAL, intent(out):: pctsrf(klon, nbsrf)
37 REAL, intent(out):: ftsol(klon, nbsrf)
38 REAL, intent(out):: ftsoil(klon, nsoilmx, nbsrf)
39 REAL, intent(out):: qsurf(klon, nbsrf)
40
41 REAL, intent(out):: qsol(:)
42 ! (klon) column-density of water in soil, in kg m-2
43
44 REAL, intent(out):: snow(klon, nbsrf)
45 REAL, intent(out):: albe(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 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 !*********************************************************************
355
356 subroutine phyetat0_new
357
358 use nr_util, only: pi
359
360 use dimensions, only: iim, jjm
361 use dynetat0_m, only: rlatu, rlonv
362 use grid_change, only: dyn_phy
363 USE start_init_orog_m, only: mask
364
365 !-------------------------------------------------------------------------
366
367 rlat(1) = 90.
368 rlat(2:klon-1) = pack(spread(rlatu(2:jjm), 1, iim), .true.) * 180. / pi
369 ! (with conversion to degrees)
370 rlat(klon) = - 90.
371
372 rlon(1) = 0.
373 rlon(2:klon-1) = pack(spread(rlonv(:iim), 2, jjm - 1), .true.) * 180. / pi
374 ! (with conversion to degrees)
375 rlon(klon) = 0.
376
377 zmasq = pack(mask, dyn_phy)
378 itau_phy = 0
379
380 end subroutine phyetat0_new
381
382 end module phyetat0_m

  ViewVC Help
Powered by ViewVC 1.1.21