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 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):: fsnow(klon, nbsrf) |
44 |
REAL, intent(out):: albe(klon, nbsrf) |
45 |
REAL, intent(out):: rain_fall(klon) |
46 |
REAL, intent(out):: snow_fall(klon) |
47 |
real, intent(out):: solsw(klon) |
48 |
REAL, intent(out):: sollw(klon) |
49 |
real, intent(out):: fder(klon) |
50 |
REAL, intent(out):: radsol(klon) |
51 |
REAL, intent(out):: frugs(klon, nbsrf) |
52 |
REAL, intent(out):: agesno(klon, nbsrf) |
53 |
REAL, intent(out):: zmea(klon) |
54 |
REAL, intent(out):: zstd(klon) |
55 |
REAL, intent(out):: zsig(klon) |
56 |
REAL, intent(out):: zgam(klon) |
57 |
REAL, intent(out):: zthe(klon) |
58 |
REAL, intent(out):: zpic(klon) |
59 |
REAL, intent(out):: zval(klon) |
60 |
REAL, intent(out):: t_ancien(klon, klev), q_ancien(klon, klev) |
61 |
LOGICAL, intent(out):: ancien_ok |
62 |
real, intent(out):: rnebcon(klon, klev), ratqs(klon, klev) |
63 |
REAL, intent(out):: clwcon(klon, klev), run_off_lic_0(klon) |
64 |
real, intent(out):: sig1(klon, klev) ! section adiabatic updraft |
65 |
|
66 |
real, intent(out):: w01(klon, klev) |
67 |
! vertical velocity within adiabatic updraft |
68 |
|
69 |
integer, intent(out):: ncid_startphy |
70 |
|
71 |
! Local: |
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 |
IF (raz_date) then |
83 |
itau_phy = 0 |
84 |
else |
85 |
call nf95_get_att(ncid_startphy, nf90_global, "itau_phy", itau_phy) |
86 |
end IF |
87 |
|
88 |
! Lecture des latitudes (coordonnees): |
89 |
|
90 |
call NF95_INQ_VARID(ncid_startphy, "latitude", varid) |
91 |
call NF95_GET_VAR(ncid_startphy, varid, rlat) |
92 |
|
93 |
! Lecture des longitudes (coordonnees): |
94 |
|
95 |
call NF95_INQ_VARID(ncid_startphy, "longitude", varid) |
96 |
call NF95_GET_VAR(ncid_startphy, varid, rlon) |
97 |
|
98 |
! Lecture du masque terre mer |
99 |
|
100 |
call NF95_INQ_VARID(ncid_startphy, "masque", varid) |
101 |
call nf95_get_var(ncid_startphy, varid, masque) |
102 |
|
103 |
! Lecture des fractions pour chaque sous-surface |
104 |
|
105 |
! initialisation des sous-surfaces |
106 |
|
107 |
call NF95_INQ_VARID(ncid_startphy, "pctsrf", varid) |
108 |
call nf95_get_var(ncid_startphy, varid, pctsrf) |
109 |
|
110 |
! Verification de l'adequation entre le masque et les sous-surfaces |
111 |
|
112 |
DO i = 1 , klon |
113 |
IF (abs(pctsrf(i, is_ter) + pctsrf(i, is_lic) - masque(i)) > EPSFRA) THEN |
114 |
print *, & |
115 |
'phyetat0: pctsrf does not agree with masque for continents', & |
116 |
i, masque(i), pctsrf(i, is_ter), pctsrf(i, is_lic) |
117 |
ENDIF |
118 |
END DO |
119 |
|
120 |
DO i = 1 , klon |
121 |
IF (abs(pctsrf(i, is_oce) + pctsrf(i, is_sic) - (1. - masque(i))) & |
122 |
> EPSFRA) THEN |
123 |
print *, 'phyetat0: pctsrf does not agree with masque for ocean ', & |
124 |
'and sea-ice', i, masque(i) , pctsrf(i, is_oce), & |
125 |
pctsrf(i, is_sic) |
126 |
ENDIF |
127 |
END DO |
128 |
|
129 |
! Lecture des temperatures du sol: |
130 |
call NF95_INQ_VARID(ncid_startphy, "TS", varid) |
131 |
call nf95_inquire_variable(ncid_startphy, varid, ndims = ndims) |
132 |
if (ndims == 2) then |
133 |
call NF95_GET_VAR(ncid_startphy, varid, ftsol) |
134 |
else |
135 |
print *, "Found only one surface type for soil temperature." |
136 |
call nf95_get_var(ncid_startphy, varid, ftsol(:, 1)) |
137 |
ftsol(:, 2:nbsrf) = spread(ftsol(:, 1), dim = 2, ncopies = nbsrf - 1) |
138 |
end if |
139 |
|
140 |
! Lecture des temperatures du sol profond: |
141 |
|
142 |
call NF95_INQ_VARID(ncid_startphy, 'Tsoil', varid) |
143 |
call NF95_GET_VAR(ncid_startphy, varid, ftsoil) |
144 |
|
145 |
! Lecture de l'humidite de l'air juste au dessus du sol: |
146 |
|
147 |
call NF95_INQ_VARID(ncid_startphy, "QS", varid) |
148 |
call nf95_get_var(ncid_startphy, varid, qsurf) |
149 |
|
150 |
ierr = NF90_INQ_VARID(ncid_startphy, "QSOL", varid) |
151 |
IF (ierr == NF90_NOERR) THEN |
152 |
call nf95_get_var(ncid_startphy, varid, qsol) |
153 |
else |
154 |
PRINT *, 'phyetat0: Le champ <QSOL> est absent' |
155 |
PRINT *, ' Valeur par defaut nulle' |
156 |
qsol = 0. |
157 |
ENDIF |
158 |
|
159 |
! Lecture de neige au sol: |
160 |
|
161 |
call NF95_INQ_VARID(ncid_startphy, "SNOW", varid) |
162 |
call nf95_get_var(ncid_startphy, varid, fsnow) |
163 |
|
164 |
! Lecture de albedo au sol: |
165 |
|
166 |
call NF95_INQ_VARID(ncid_startphy, "ALBE", varid) |
167 |
call nf95_get_var(ncid_startphy, varid, albe) |
168 |
|
169 |
! Lecture precipitation liquide: |
170 |
|
171 |
call NF95_INQ_VARID(ncid_startphy, "rain_f", varid) |
172 |
call NF95_GET_VAR(ncid_startphy, varid, rain_fall) |
173 |
|
174 |
! Lecture precipitation solide: |
175 |
|
176 |
call NF95_INQ_VARID(ncid_startphy, "snow_f", varid) |
177 |
call NF95_GET_VAR(ncid_startphy, varid, snow_fall) |
178 |
|
179 |
! Lecture rayonnement solaire au sol: |
180 |
|
181 |
ierr = NF90_INQ_VARID(ncid_startphy, "solsw", varid) |
182 |
IF (ierr /= NF90_NOERR) THEN |
183 |
PRINT *, 'phyetat0: Le champ <solsw> est absent' |
184 |
PRINT *, 'mis a zero' |
185 |
solsw = 0. |
186 |
ELSE |
187 |
call nf95_get_var(ncid_startphy, varid, solsw) |
188 |
ENDIF |
189 |
|
190 |
! Lecture rayonnement IF au sol: |
191 |
|
192 |
ierr = NF90_INQ_VARID(ncid_startphy, "sollw", varid) |
193 |
IF (ierr /= NF90_NOERR) THEN |
194 |
PRINT *, 'phyetat0: Le champ <sollw> est absent' |
195 |
PRINT *, 'mis a zero' |
196 |
sollw = 0. |
197 |
ELSE |
198 |
call nf95_get_var(ncid_startphy, varid, sollw) |
199 |
ENDIF |
200 |
|
201 |
! Lecture derive des flux: |
202 |
|
203 |
ierr = NF90_INQ_VARID(ncid_startphy, "fder", varid) |
204 |
IF (ierr /= NF90_NOERR) THEN |
205 |
PRINT *, 'phyetat0: Le champ <fder> est absent' |
206 |
PRINT *, 'mis a zero' |
207 |
fder = 0. |
208 |
ELSE |
209 |
call nf95_get_var(ncid_startphy, varid, fder) |
210 |
ENDIF |
211 |
|
212 |
! Lecture du rayonnement net au sol: |
213 |
|
214 |
call NF95_INQ_VARID(ncid_startphy, "RADS", varid) |
215 |
call NF95_GET_VAR(ncid_startphy, varid, radsol) |
216 |
|
217 |
! Lecture de la longueur de rugosite |
218 |
|
219 |
call NF95_INQ_VARID(ncid_startphy, "RUG", varid) |
220 |
call nf95_get_var(ncid_startphy, varid, frugs) |
221 |
|
222 |
! Lecture de l'age de la neige: |
223 |
|
224 |
call NF95_INQ_VARID(ncid_startphy, "AGESNO", varid) |
225 |
call nf95_get_var(ncid_startphy, varid, agesno) |
226 |
|
227 |
call NF95_INQ_VARID(ncid_startphy, "ZMEA", varid) |
228 |
call NF95_GET_VAR(ncid_startphy, varid, zmea) |
229 |
|
230 |
call NF95_INQ_VARID(ncid_startphy, "ZSTD", varid) |
231 |
call NF95_GET_VAR(ncid_startphy, varid, zstd) |
232 |
|
233 |
call NF95_INQ_VARID(ncid_startphy, "ZSIG", varid) |
234 |
call NF95_GET_VAR(ncid_startphy, varid, zsig) |
235 |
|
236 |
call NF95_INQ_VARID(ncid_startphy, "ZGAM", varid) |
237 |
call NF95_GET_VAR(ncid_startphy, varid, zgam) |
238 |
|
239 |
call NF95_INQ_VARID(ncid_startphy, "ZTHE", varid) |
240 |
call NF95_GET_VAR(ncid_startphy, varid, zthe) |
241 |
|
242 |
call NF95_INQ_VARID(ncid_startphy, "ZPIC", varid) |
243 |
call NF95_GET_VAR(ncid_startphy, varid, zpic) |
244 |
|
245 |
call NF95_INQ_VARID(ncid_startphy, "ZVAL", varid) |
246 |
call NF95_GET_VAR(ncid_startphy, varid, zval) |
247 |
|
248 |
ancien_ok = .TRUE. |
249 |
|
250 |
ierr = NF90_INQ_VARID(ncid_startphy, "TANCIEN", varid) |
251 |
IF (ierr /= NF90_NOERR) THEN |
252 |
PRINT *, "phyetat0: Le champ <TANCIEN> est absent" |
253 |
PRINT *, "Depart legerement fausse. Mais je continue" |
254 |
ancien_ok = .FALSE. |
255 |
ELSE |
256 |
call nf95_get_var(ncid_startphy, varid, t_ancien) |
257 |
ENDIF |
258 |
|
259 |
ierr = NF90_INQ_VARID(ncid_startphy, "QANCIEN", varid) |
260 |
IF (ierr /= NF90_NOERR) THEN |
261 |
PRINT *, "phyetat0: Le champ <QANCIEN> est absent" |
262 |
PRINT *, "Depart legerement fausse. Mais je continue" |
263 |
ancien_ok = .FALSE. |
264 |
ELSE |
265 |
call nf95_get_var(ncid_startphy, varid, q_ancien) |
266 |
ENDIF |
267 |
|
268 |
ierr = NF90_INQ_VARID(ncid_startphy, "CLWCON", varid) |
269 |
IF (ierr /= NF90_NOERR) THEN |
270 |
PRINT *, "phyetat0: Le champ CLWCON est absent" |
271 |
PRINT *, "Depart legerement fausse. Mais je continue" |
272 |
clwcon = 0. |
273 |
ELSE |
274 |
call nf95_get_var(ncid_startphy, varid, clwcon(:, 1)) |
275 |
clwcon(:, 2:) = 0. |
276 |
ENDIF |
277 |
|
278 |
ierr = NF90_INQ_VARID(ncid_startphy, "RNEBCON", varid) |
279 |
|
280 |
IF (ierr /= NF90_NOERR) THEN |
281 |
PRINT *, "phyetat0: le champ RNEBCON est absent." |
282 |
PRINT *, "D\'epart l\'eg\`erement fauss\'e. 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 |
|
343 |
end subroutine phyetat0_new |
344 |
|
345 |
end module phyetat0_m |