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

Contents of /trunk/phylmd/phyetat0.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 330 - (show annotations)
Wed Jul 31 14:55:23 2019 UTC (4 years, 9 months ago) by guez
File size: 11003 byte(s)
Bug fix in `CMakeLists.txt`: default value of jjm.

Use the same variable name across procedures for column density of
snow at the surface: fsnow.

Remove useless intermediary variable bils in procedure physiq.

In procedure interfsurf_hq, no need to force zfra between 0 and 1: the
result of the computation is already between 0 and 1 since snow is
necessarily greater than or equal to 0.

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:: masque(KLON) ! fraction of land
13
14 private klon
15
16 contains
17
18 SUBROUTINE phyetat0(pctsrf, ftsol, ftsoil, qsurf, qsol, fsnow, 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):: fsnow(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 INTEGER varid, ndims
74 INTEGER ierr, i
75
76 !---------------------------------------------------------------
77
78 print *, "Call sequence information: phyetat0"
79
80 ! Fichier contenant l'état initial :
81 call NF95_OPEN("startphy.nc", NF90_NOWRITE, ncid_startphy)
82
83 IF (raz_date) then
84 itau_phy = 0
85 else
86 call nf95_get_att(ncid_startphy, nf90_global, "itau_phy", itau_phy)
87 end IF
88
89 ! Lecture des latitudes (coordonnees):
90
91 call NF95_INQ_VARID(ncid_startphy, "latitude", varid)
92 call NF95_GET_VAR(ncid_startphy, varid, rlat)
93
94 ! Lecture des longitudes (coordonnees):
95
96 call NF95_INQ_VARID(ncid_startphy, "longitude", varid)
97 call NF95_GET_VAR(ncid_startphy, varid, rlon)
98
99 ! Lecture du masque terre mer
100
101 call NF95_INQ_VARID(ncid_startphy, "masque", varid)
102 call nf95_get_var(ncid_startphy, varid, masque)
103
104 ! Lecture des fractions pour chaque sous-surface
105
106 ! initialisation des sous-surfaces
107
108 call NF95_INQ_VARID(ncid_startphy, "pctsrf", varid)
109 call nf95_get_var(ncid_startphy, varid, pctsrf)
110
111 ! Verification de l'adequation entre le masque et les sous-surfaces
112
113 DO i = 1 , klon
114 IF (abs(pctsrf(i, is_ter) + pctsrf(i, is_lic) - masque(i)) > EPSFRA) THEN
115 print *, &
116 'phyetat0: pctsrf does not agree with masque for continents', &
117 i, masque(i), pctsrf(i, is_ter), pctsrf(i, is_lic)
118 ENDIF
119 END DO
120
121 DO i = 1 , klon
122 IF (abs(pctsrf(i, is_oce) + pctsrf(i, is_sic) - (1. - masque(i))) &
123 > EPSFRA) THEN
124 print *, 'phyetat0: pctsrf does not agree with masque for ocean ', &
125 'and sea-ice', i, masque(i) , pctsrf(i, is_oce), &
126 pctsrf(i, is_sic)
127 ENDIF
128 END DO
129
130 ! Lecture des temperatures du sol:
131 call NF95_INQ_VARID(ncid_startphy, "TS", varid)
132 call nf95_inquire_variable(ncid_startphy, varid, ndims = ndims)
133 if (ndims == 2) then
134 call NF95_GET_VAR(ncid_startphy, varid, ftsol)
135 else
136 print *, "Found only one surface type for soil temperature."
137 call nf95_get_var(ncid_startphy, varid, ftsol(:, 1))
138 ftsol(:, 2:nbsrf) = spread(ftsol(:, 1), dim = 2, ncopies = nbsrf - 1)
139 end if
140
141 ! Lecture des temperatures du sol profond:
142
143 call NF95_INQ_VARID(ncid_startphy, 'Tsoil', varid)
144 call NF95_GET_VAR(ncid_startphy, varid, ftsoil)
145
146 ! Lecture de l'humidite de l'air juste au dessus du sol:
147
148 call NF95_INQ_VARID(ncid_startphy, "QS", varid)
149 call nf95_get_var(ncid_startphy, varid, qsurf)
150
151 ierr = NF90_INQ_VARID(ncid_startphy, "QSOL", varid)
152 IF (ierr == NF90_NOERR) THEN
153 call nf95_get_var(ncid_startphy, varid, qsol)
154 else
155 PRINT *, 'phyetat0: Le champ <QSOL> est absent'
156 PRINT *, ' Valeur par defaut nulle'
157 qsol = 0.
158 ENDIF
159
160 ! Lecture de neige au sol:
161
162 call NF95_INQ_VARID(ncid_startphy, "SNOW", varid)
163 call nf95_get_var(ncid_startphy, varid, fsnow)
164
165 ! Lecture de albedo au sol:
166
167 call NF95_INQ_VARID(ncid_startphy, "ALBE", varid)
168 call nf95_get_var(ncid_startphy, varid, albe)
169
170 ! Lecture precipitation liquide:
171
172 call NF95_INQ_VARID(ncid_startphy, "rain_f", varid)
173 call NF95_GET_VAR(ncid_startphy, varid, rain_fall)
174
175 ! Lecture precipitation solide:
176
177 call NF95_INQ_VARID(ncid_startphy, "snow_f", varid)
178 call NF95_GET_VAR(ncid_startphy, varid, snow_fall)
179
180 ! Lecture rayonnement solaire au sol:
181
182 ierr = NF90_INQ_VARID(ncid_startphy, "solsw", varid)
183 IF (ierr /= NF90_NOERR) THEN
184 PRINT *, 'phyetat0: Le champ <solsw> est absent'
185 PRINT *, 'mis a zero'
186 solsw = 0.
187 ELSE
188 call nf95_get_var(ncid_startphy, varid, solsw)
189 ENDIF
190
191 ! Lecture rayonnement IF au sol:
192
193 ierr = NF90_INQ_VARID(ncid_startphy, "sollw", varid)
194 IF (ierr /= NF90_NOERR) THEN
195 PRINT *, 'phyetat0: Le champ <sollw> est absent'
196 PRINT *, 'mis a zero'
197 sollw = 0.
198 ELSE
199 call nf95_get_var(ncid_startphy, varid, sollw)
200 ENDIF
201
202 ! Lecture derive des flux:
203
204 ierr = NF90_INQ_VARID(ncid_startphy, "fder", varid)
205 IF (ierr /= NF90_NOERR) THEN
206 PRINT *, 'phyetat0: Le champ <fder> est absent'
207 PRINT *, 'mis a zero'
208 fder = 0.
209 ELSE
210 call nf95_get_var(ncid_startphy, varid, fder)
211 ENDIF
212
213 ! Lecture du rayonnement net au sol:
214
215 call NF95_INQ_VARID(ncid_startphy, "RADS", varid)
216 call NF95_GET_VAR(ncid_startphy, varid, radsol)
217
218 ! Lecture de la longueur de rugosite
219
220 call NF95_INQ_VARID(ncid_startphy, "RUG", varid)
221 call nf95_get_var(ncid_startphy, varid, frugs)
222
223 ! Lecture de l'age de la neige:
224
225 call NF95_INQ_VARID(ncid_startphy, "AGESNO", varid)
226 call nf95_get_var(ncid_startphy, varid, agesno)
227
228 call NF95_INQ_VARID(ncid_startphy, "ZMEA", varid)
229 call NF95_GET_VAR(ncid_startphy, varid, zmea)
230
231 call NF95_INQ_VARID(ncid_startphy, "ZSTD", varid)
232 call NF95_GET_VAR(ncid_startphy, varid, zstd)
233
234 call NF95_INQ_VARID(ncid_startphy, "ZSIG", varid)
235 call NF95_GET_VAR(ncid_startphy, varid, zsig)
236
237 call NF95_INQ_VARID(ncid_startphy, "ZGAM", varid)
238 call NF95_GET_VAR(ncid_startphy, varid, zgam)
239
240 call NF95_INQ_VARID(ncid_startphy, "ZTHE", varid)
241 call NF95_GET_VAR(ncid_startphy, varid, zthe)
242
243 call NF95_INQ_VARID(ncid_startphy, "ZPIC", varid)
244 call NF95_GET_VAR(ncid_startphy, varid, zpic)
245
246 call NF95_INQ_VARID(ncid_startphy, "ZVAL", varid)
247 call NF95_GET_VAR(ncid_startphy, varid, zval)
248
249 ancien_ok = .TRUE.
250
251 ierr = NF90_INQ_VARID(ncid_startphy, "TANCIEN", varid)
252 IF (ierr /= NF90_NOERR) THEN
253 PRINT *, "phyetat0: Le champ <TANCIEN> est absent"
254 PRINT *, "Depart legerement fausse. Mais je continue"
255 ancien_ok = .FALSE.
256 ELSE
257 call nf95_get_var(ncid_startphy, varid, t_ancien)
258 ENDIF
259
260 ierr = NF90_INQ_VARID(ncid_startphy, "QANCIEN", varid)
261 IF (ierr /= NF90_NOERR) THEN
262 PRINT *, "phyetat0: Le champ <QANCIEN> est absent"
263 PRINT *, "Depart legerement fausse. Mais je continue"
264 ancien_ok = .FALSE.
265 ELSE
266 call nf95_get_var(ncid_startphy, varid, q_ancien)
267 ENDIF
268
269 ierr = NF90_INQ_VARID(ncid_startphy, "CLWCON", varid)
270 IF (ierr /= NF90_NOERR) THEN
271 PRINT *, "phyetat0: Le champ CLWCON est absent"
272 PRINT *, "Depart legerement fausse. Mais je continue"
273 clwcon = 0.
274 ELSE
275 call nf95_get_var(ncid_startphy, varid, clwcon(:, 1))
276 clwcon(:, 2:) = 0.
277 ENDIF
278
279 ierr = NF90_INQ_VARID(ncid_startphy, "RNEBCON", varid)
280 IF (ierr /= NF90_NOERR) THEN
281 PRINT *, "phyetat0: Le champ RNEBCON est absent"
282 PRINT *, "Depart legerement fausse. Mais je continue"
283 rnebcon = 0.
284 ELSE
285 call nf95_get_var(ncid_startphy, varid, rnebcon(:, 1))
286 rnebcon(:, 2:) = 0.
287 ENDIF
288
289 ! Lecture ratqs
290
291 ierr = NF90_INQ_VARID(ncid_startphy, "RATQS", varid)
292 IF (ierr /= NF90_NOERR) THEN
293 PRINT *, "phyetat0: Le champ <RATQS> est absent"
294 PRINT *, "Depart legerement fausse. Mais je continue"
295 ratqs = 0.
296 ELSE
297 call nf95_get_var(ncid_startphy, varid, ratqs(:, 1))
298 ratqs(:, 2:) = 0.
299 ENDIF
300
301 ! Lecture run_off_lic_0
302
303 ierr = NF90_INQ_VARID(ncid_startphy, "RUNOFFLIC0", varid)
304 IF (ierr /= NF90_NOERR) THEN
305 PRINT *, "phyetat0: Le champ <RUNOFFLIC0> est absent"
306 PRINT *, "Depart legerement fausse. Mais je continue"
307 run_off_lic_0 = 0.
308 ELSE
309 call nf95_get_var(ncid_startphy, varid, run_off_lic_0)
310 ENDIF
311
312 call nf95_inq_varid(ncid_startphy, "sig1", varid)
313 call nf95_get_var(ncid_startphy, varid, sig1)
314
315 call nf95_inq_varid(ncid_startphy, "w01", varid)
316 call nf95_get_var(ncid_startphy, varid, w01)
317
318 END SUBROUTINE phyetat0
319
320 !*********************************************************************
321
322 subroutine phyetat0_new
323
324 use nr_util, only: rad_to_deg
325
326 use dimensions, only: iim, jjm
327 use dynetat0_m, only: rlatu, rlonv
328 use grid_change, only: dyn_phy
329 USE start_init_orog_m, only: mask
330
331 !-------------------------------------------------------------------------
332
333 rlat(1) = 90.
334 rlat(2:klon-1) = pack(spread(rlatu(2:jjm), 1, iim), .true.) * rad_to_deg
335 rlat(klon) = - 90.
336
337 rlon(1) = 0.
338 rlon(2:klon-1) = pack(spread(rlonv(:iim), 2, jjm - 1), .true.) * rad_to_deg
339 rlon(klon) = 0.
340
341 masque = pack(mask, dyn_phy)
342 itau_phy = 0
343
344 end subroutine phyetat0_new
345
346 end module phyetat0_m

  ViewVC Help
Powered by ViewVC 1.1.21