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

Annotation of /trunk/phylmd/phyetat0.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 310 - (hide annotations)
Thu Sep 27 16:29:06 2018 UTC (5 years, 7 months ago) by guez
Original Path: trunk/phylmd/phyetat0.f
File size: 10997 byte(s)
Read and write the whole pctsrf array in (re)startphy.nc, instead of
splitting it into FTER, FLIC, FOCE, FSIC.

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

  ViewVC Help
Powered by ViewVC 1.1.21