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