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

Annotation of /trunk/phylmd/phyetat0.f90

Parent Directory Parent Directory | Revision Log Revision Log


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

  ViewVC Help
Powered by ViewVC 1.1.21