/[lmdze]/trunk/phylmd/physiq.f
ViewVC logotype

Contents of /trunk/phylmd/physiq.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 209 - (show annotations)
Wed Dec 7 17:37:21 2016 UTC (7 years, 5 months ago) by guez
Original Path: trunk/Sources/phylmd/physiq.f
File size: 40294 byte(s)
The program did not work with cycle_diurne set to false. mu0 in
physiq, which is supposed to be a cosine, was set to -999.999. So prmu
in swu had a value of the order of 1e3. So zrmum1 in sw2s had a value
of the order of 1e3. So zrayl in sw2s had a value of the order of
1e15. So ztray and ptauaz in swclr also had a large value. So zcorae
at line 138 in swclr had a large negative value, which resulted in
overflow at line 138.

This assignment of -999.999 to mu0 dates from somewhere between
revisions 348 and 524 of LMDZ. It was corrected in revision 1068 of
LMDZ with a call to angle which was present in revision 348. However,
procedure angle was removed from LMDZE in revision 22 because it was
not used. Hesitated to bring back angle but, finally, just removed the
option of having no diurnal cycle.

1 module physiq_m
2
3 IMPLICIT none
4
5 contains
6
7 SUBROUTINE physiq(lafin, dayvrai, time, paprs, play, pphi, pphis, u, v, t, &
8 qx, omega, d_u, d_v, d_t, d_qx)
9
10 ! From phylmd/physiq.F, version 1.22 2006/02/20 09:38:28
11 ! (subversion revision 678)
12
13 ! Author: Z. X. Li (LMD/CNRS) 1993
14
15 ! This is the main procedure for the "physics" part of the program.
16
17 use aaam_bud_m, only: aaam_bud
18 USE abort_gcm_m, ONLY: abort_gcm
19 use ajsec_m, only: ajsec
20 use calltherm_m, only: calltherm
21 USE clesphys, ONLY: cdhmax, cdmmax, ecrit_ins, ksta, ksta_ter, ok_kzmin, &
22 ok_instan
23 USE clesphys2, ONLY: conv_emanuel, nbapp_rad, new_oliq, ok_orodr, ok_orolf
24 USE clmain_m, ONLY: clmain
25 use clouds_gno_m, only: clouds_gno
26 use comconst, only: dtphys
27 USE comgeomphy, ONLY: airephy
28 USE concvl_m, ONLY: concvl
29 USE conf_gcm_m, ONLY: offline, lmt_pas
30 USE conf_phys_m, ONLY: conf_phys
31 use conflx_m, only: conflx
32 USE ctherm, ONLY: iflag_thermals, nsplit_thermals
33 use diagcld2_m, only: diagcld2
34 USE dimens_m, ONLY: llm, nqmx
35 USE dimphy, ONLY: klon
36 USE dimsoil, ONLY: nsoilmx
37 use drag_noro_m, only: drag_noro
38 use dynetat0_m, only: day_ref, annee_ref
39 USE fcttre, ONLY: foeew, qsatl, qsats
40 use fisrtilp_m, only: fisrtilp
41 USE hgardfou_m, ONLY: hgardfou
42 USE histsync_m, ONLY: histsync
43 USE histwrite_phy_m, ONLY: histwrite_phy
44 USE indicesol, ONLY: clnsurf, epsfra, is_lic, is_oce, is_sic, is_ter, &
45 nbsrf
46 USE ini_histins_m, ONLY: ini_histins, nid_ins
47 use netcdf95, only: NF95_CLOSE
48 use newmicro_m, only: newmicro
49 use nr_util, only: assert
50 use nuage_m, only: nuage
51 USE orbite_m, ONLY: orbite
52 USE ozonecm_m, ONLY: ozonecm
53 USE phyetat0_m, ONLY: phyetat0, rlat, rlon
54 USE phyredem_m, ONLY: phyredem
55 USE phyredem0_m, ONLY: phyredem0
56 USE phystokenc_m, ONLY: phystokenc
57 USE phytrac_m, ONLY: phytrac
58 use radlwsw_m, only: radlwsw
59 use yoegwd, only: sugwd
60 USE suphec_m, ONLY: rcpd, retv, rg, rlvtt, romega, rsigma, rtt
61 use time_phylmdz, only: itap, increment_itap
62 use transp_m, only: transp
63 use transp_lay_m, only: transp_lay
64 use unit_nml_m, only: unit_nml
65 USE ymds2ju_m, ONLY: ymds2ju
66 USE yoethf_m, ONLY: r2es, rvtmp2
67 use zenang_m, only: zenang
68
69 logical, intent(in):: lafin ! dernier passage
70
71 integer, intent(in):: dayvrai
72 ! current day number, based at value 1 on January 1st of annee_ref
73
74 REAL, intent(in):: time ! heure de la journ\'ee en fraction de jour
75
76 REAL, intent(in):: paprs(:, :) ! (klon, llm + 1)
77 ! pression pour chaque inter-couche, en Pa
78
79 REAL, intent(in):: play(:, :) ! (klon, llm)
80 ! pression pour le mileu de chaque couche (en Pa)
81
82 REAL, intent(in):: pphi(:, :) ! (klon, llm)
83 ! géopotentiel de chaque couche (référence sol)
84
85 REAL, intent(in):: pphis(:) ! (klon) géopotentiel du sol
86
87 REAL, intent(in):: u(:, :) ! (klon, llm)
88 ! vitesse dans la direction X (de O a E) en m / s
89
90 REAL, intent(in):: v(:, :) ! (klon, llm) vitesse Y (de S a N) en m / s
91 REAL, intent(in):: t(:, :) ! (klon, llm) temperature (K)
92
93 REAL, intent(in):: qx(:, :, :) ! (klon, llm, nqmx)
94 ! (humidit\'e sp\'ecifique et fractions massiques des autres traceurs)
95
96 REAL, intent(in):: omega(:, :) ! (klon, llm) vitesse verticale en Pa / s
97 REAL, intent(out):: d_u(:, :) ! (klon, llm) tendance physique de "u" (m s-2)
98 REAL, intent(out):: d_v(:, :) ! (klon, llm) tendance physique de "v" (m s-2)
99 REAL, intent(out):: d_t(:, :) ! (klon, llm) tendance physique de "t" (K / s)
100
101 REAL, intent(out):: d_qx(:, :, :) ! (klon, llm, nqmx)
102 ! tendance physique de "qx" (s-1)
103
104 ! Local:
105
106 LOGICAL:: firstcal = .true.
107
108 LOGICAL, PARAMETER:: ok_stratus = .FALSE.
109 ! Ajouter artificiellement les stratus
110
111 ! pour phystoke avec thermiques
112 REAL fm_therm(klon, llm + 1)
113 REAL entr_therm(klon, llm)
114 real, save:: q2(klon, llm + 1, nbsrf)
115
116 INTEGER, PARAMETER:: ivap = 1 ! indice de traceur pour vapeur d'eau
117 INTEGER, PARAMETER:: iliq = 2 ! indice de traceur pour eau liquide
118
119 REAL, save:: t_ancien(klon, llm), q_ancien(klon, llm)
120 LOGICAL, save:: ancien_ok
121
122 REAL d_t_dyn(klon, llm) ! tendance dynamique pour "t" (K / s)
123 REAL d_q_dyn(klon, llm) ! tendance dynamique pour "q" (kg / kg / s)
124
125 real da(klon, llm), phi(klon, llm, llm), mp(klon, llm)
126
127 REAL, save:: swdn0(klon, llm + 1), swdn(klon, llm + 1)
128 REAL, save:: swup0(klon, llm + 1), swup(klon, llm + 1)
129
130 REAL, save:: lwdn0(klon, llm + 1), lwdn(klon, llm + 1)
131 REAL, save:: lwup0(klon, llm + 1), lwup(klon, llm + 1)
132
133 ! prw: precipitable water
134 real prw(klon)
135
136 ! flwp, fiwp = Liquid Water Path & Ice Water Path (kg / m2)
137 ! flwc, fiwc = Liquid Water Content & Ice Water Content (kg / kg)
138 REAL flwp(klon), fiwp(klon)
139 REAL flwc(klon, llm), fiwc(klon, llm)
140
141 ! Variables propres a la physique
142
143 INTEGER, save:: radpas
144 ! Radiative transfer computations are made every "radpas" call to
145 ! "physiq".
146
147 REAL, save:: radsol(klon) ! bilan radiatif au sol calcule par code radiatif
148 REAL, save:: ftsol(klon, nbsrf) ! skin temperature of surface fraction
149
150 REAL, save:: ftsoil(klon, nsoilmx, nbsrf)
151 ! soil temperature of surface fraction
152
153 REAL, save:: fevap(klon, nbsrf) ! evaporation
154 REAL, save:: fluxlat(klon, nbsrf)
155
156 REAL, save:: fqsurf(klon, nbsrf)
157 ! humidite de l'air au contact de la surface
158
159 REAL, save:: qsol(klon)
160 ! column-density of water in soil, in kg m-2
161
162 REAL, save:: fsnow(klon, nbsrf) ! epaisseur neigeuse
163 REAL, save:: falbe(klon, nbsrf) ! albedo visible par type de surface
164
165 ! Param\`etres de l'orographie \`a l'\'echelle sous-maille (OESM) :
166 REAL, save:: zmea(klon) ! orographie moyenne
167 REAL, save:: zstd(klon) ! deviation standard de l'OESM
168 REAL, save:: zsig(klon) ! pente de l'OESM
169 REAL, save:: zgam(klon) ! anisotropie de l'OESM
170 REAL, save:: zthe(klon) ! orientation de l'OESM
171 REAL, save:: zpic(klon) ! Maximum de l'OESM
172 REAL, save:: zval(klon) ! Minimum de l'OESM
173 REAL, save:: rugoro(klon) ! longueur de rugosite de l'OESM
174 REAL zulow(klon), zvlow(klon)
175 INTEGER igwd, itest(klon)
176
177 REAL, save:: agesno(klon, nbsrf) ! age de la neige
178 REAL, save:: run_off_lic_0(klon)
179
180 ! Variables li\'ees \`a la convection d'Emanuel :
181 REAL, save:: Ma(klon, llm) ! undilute upward mass flux
182 REAL, save:: qcondc(klon, llm) ! in-cld water content from convect
183 REAL, save:: sig1(klon, llm), w01(klon, llm)
184
185 ! Variables pour la couche limite (Alain Lahellec) :
186 REAL cdragh(klon) ! drag coefficient pour T and Q
187 REAL cdragm(klon) ! drag coefficient pour vent
188
189 ! Pour phytrac :
190 REAL ycoefh(klon, llm) ! coef d'echange pour phytrac
191 REAL yu1(klon) ! vents dans la premiere couche U
192 REAL yv1(klon) ! vents dans la premiere couche V
193
194 REAL, save:: ffonte(klon, nbsrf)
195 ! flux thermique utilise pour fondre la neige
196
197 REAL, save:: fqcalving(klon, nbsrf)
198 ! flux d'eau "perdue" par la surface et necessaire pour limiter la
199 ! hauteur de neige, en kg / m2 / s
200
201 REAL zxffonte(klon), zxfqcalving(klon)
202
203 REAL, save:: pfrac_impa(klon, llm)! Produits des coefs lessivage impaction
204 REAL, save:: pfrac_nucl(klon, llm)! Produits des coefs lessivage nucleation
205
206 REAL, save:: pfrac_1nucl(klon, llm)
207 ! Produits des coefs lessi nucl (alpha = 1)
208
209 REAL frac_impa(klon, llm) ! fractions d'aerosols lessivees (impaction)
210 REAL frac_nucl(klon, llm) ! idem (nucleation)
211
212 REAL, save:: rain_fall(klon)
213 ! liquid water mass flux (kg / m2 / s), positive down
214
215 REAL, save:: snow_fall(klon)
216 ! solid water mass flux (kg / m2 / s), positive down
217
218 REAL rain_tiedtke(klon), snow_tiedtke(klon)
219
220 REAL evap(klon) ! flux d'\'evaporation au sol
221 real devap(klon) ! derivative of the evaporation flux at the surface
222 REAL sens(klon) ! flux de chaleur sensible au sol
223 real dsens(klon) ! derivee du flux de chaleur sensible au sol
224 REAL, save:: dlw(klon) ! derivee infra rouge
225 REAL bils(klon) ! bilan de chaleur au sol
226 REAL, save:: fder(klon) ! Derive de flux (sensible et latente)
227 REAL ve(klon) ! integr. verticale du transport meri. de l'energie
228 REAL vq(klon) ! integr. verticale du transport meri. de l'eau
229 REAL ue(klon) ! integr. verticale du transport zonal de l'energie
230 REAL uq(klon) ! integr. verticale du transport zonal de l'eau
231
232 REAL, save:: frugs(klon, nbsrf) ! longueur de rugosite
233 REAL zxrugs(klon) ! longueur de rugosite
234
235 ! Conditions aux limites
236
237 INTEGER julien
238 REAL, save:: pctsrf(klon, nbsrf) ! percentage of surface
239 REAL, save:: albsol(klon) ! albedo du sol total visible
240 REAL, SAVE:: wo(klon, llm) ! column density of ozone in a cell, in kDU
241
242 real, save:: clwcon(klon, llm), rnebcon(klon, llm)
243 real, save:: clwcon0(klon, llm), rnebcon0(klon, llm)
244
245 REAL rhcl(klon, llm) ! humiditi relative ciel clair
246 REAL dialiq(klon, llm) ! eau liquide nuageuse
247 REAL diafra(klon, llm) ! fraction nuageuse
248 REAL cldliq(klon, llm) ! eau liquide nuageuse
249 REAL cldfra(klon, llm) ! fraction nuageuse
250 REAL cldtau(klon, llm) ! epaisseur optique
251 REAL cldemi(klon, llm) ! emissivite infrarouge
252
253 REAL flux_q(klon, nbsrf) ! flux turbulent d'humidite à la surface
254 REAL flux_t(klon, nbsrf) ! flux turbulent de chaleur à la surface
255 REAL flux_u(klon, nbsrf) ! flux turbulent de vitesse u à la surface
256 REAL flux_v(klon, nbsrf) ! flux turbulent de vitesse v à la surface
257
258 ! Le rayonnement n'est pas calcul\'e tous les pas, il faut donc que
259 ! les variables soient r\'emanentes.
260 REAL, save:: heat(klon, llm) ! chauffage solaire
261 REAL, save:: heat0(klon, llm) ! chauffage solaire ciel clair
262 REAL, save:: cool(klon, llm) ! refroidissement infrarouge
263 REAL, save:: cool0(klon, llm) ! refroidissement infrarouge ciel clair
264 REAL, save:: topsw(klon), toplw(klon), solsw(klon)
265 REAL, save:: sollw(klon) ! rayonnement infrarouge montant \`a la surface
266 real, save:: sollwdown(klon) ! downward LW flux at surface
267 REAL, save:: topsw0(klon), toplw0(klon), solsw0(klon), sollw0(klon)
268 REAL, save:: albpla(klon)
269 REAL fsollw(klon, nbsrf) ! bilan flux IR pour chaque sous-surface
270 REAL fsolsw(klon, nbsrf) ! flux solaire absorb\'e pour chaque sous-surface
271
272 REAL conv_q(klon, llm) ! convergence de l'humidite (kg / kg / s)
273 REAL conv_t(klon, llm) ! convergence of temperature (K / s)
274
275 REAL cldl(klon), cldm(klon), cldh(klon) ! nuages bas, moyen et haut
276 REAL cldt(klon), cldq(klon) ! nuage total, eau liquide integree
277
278 REAL zxqsurf(klon), zxsnow(klon), zxfluxlat(klon)
279
280 REAL dist, mu0(klon), fract(klon)
281 real longi
282 REAL z_avant(klon), z_apres(klon), z_factor(klon)
283 REAL zb
284 REAL zx_t, zx_qs, zcor
285 real zqsat(klon, llm)
286 INTEGER i, k, iq, nsrf
287 REAL zphi(klon, llm)
288
289 ! cf. Anne Mathieu, variables pour la couche limite atmosphérique (hbtm)
290
291 REAL, SAVE:: pblh(klon, nbsrf) ! Hauteur de couche limite
292 REAL, SAVE:: plcl(klon, nbsrf) ! Niveau de condensation de la CLA
293 REAL, SAVE:: capCL(klon, nbsrf) ! CAPE de couche limite
294 REAL, SAVE:: oliqCL(klon, nbsrf) ! eau_liqu integree de couche limite
295 REAL, SAVE:: cteiCL(klon, nbsrf) ! cloud top instab. crit. couche limite
296 REAL, SAVE:: pblt(klon, nbsrf) ! T \`a la hauteur de couche limite
297 REAL, SAVE:: therm(klon, nbsrf)
298 REAL, SAVE:: trmb1(klon, nbsrf) ! deep_cape
299 REAL, SAVE:: trmb2(klon, nbsrf) ! inhibition
300 REAL, SAVE:: trmb3(klon, nbsrf) ! Point Omega
301 ! Grandeurs de sorties
302 REAL s_pblh(klon), s_lcl(klon), s_capCL(klon)
303 REAL s_oliqCL(klon), s_cteiCL(klon), s_pblt(klon)
304 REAL s_therm(klon), s_trmb1(klon), s_trmb2(klon)
305 REAL s_trmb3(klon)
306
307 ! Variables pour la convection de K. Emanuel :
308
309 REAL upwd(klon, llm) ! saturated updraft mass flux
310 REAL dnwd(klon, llm) ! saturated downdraft mass flux
311 REAL, save:: cape(klon)
312
313 INTEGER iflagctrl(klon) ! flag fonctionnement de convect
314
315 ! Variables du changement
316
317 ! con: convection
318 ! lsc: large scale condensation
319 ! ajs: ajustement sec
320 ! eva: \'evaporation de l'eau liquide nuageuse
321 ! vdf: vertical diffusion in boundary layer
322 REAL d_t_con(klon, llm), d_q_con(klon, llm)
323 REAL, save:: d_u_con(klon, llm), d_v_con(klon, llm)
324 REAL d_t_lsc(klon, llm), d_q_lsc(klon, llm), d_ql_lsc(klon, llm)
325 REAL d_t_ajs(klon, llm), d_q_ajs(klon, llm)
326 REAL d_u_ajs(klon, llm), d_v_ajs(klon, llm)
327 REAL rneb(klon, llm)
328
329 REAL mfu(klon, llm), mfd(klon, llm)
330 REAL pen_u(klon, llm), pen_d(klon, llm)
331 REAL pde_u(klon, llm), pde_d(klon, llm)
332 INTEGER kcbot(klon), kctop(klon), kdtop(klon)
333 REAL pmflxr(klon, llm + 1), pmflxs(klon, llm + 1)
334 REAL prfl(klon, llm + 1), psfl(klon, llm + 1)
335
336 INTEGER, save:: ibas_con(klon), itop_con(klon)
337 real ema_pct(klon) ! Emanuel pressure at cloud top, in Pa
338
339 REAL, save:: rain_con(klon)
340 real rain_lsc(klon)
341 REAL, save:: snow_con(klon) ! neige (mm / s)
342 real snow_lsc(klon)
343 REAL d_ts(klon, nbsrf)
344
345 REAL d_u_vdf(klon, llm), d_v_vdf(klon, llm)
346 REAL d_t_vdf(klon, llm), d_q_vdf(klon, llm)
347
348 REAL d_u_oro(klon, llm), d_v_oro(klon, llm)
349 REAL d_t_oro(klon, llm)
350 REAL d_u_lif(klon, llm), d_v_lif(klon, llm)
351 REAL d_t_lif(klon, llm)
352
353 REAL, save:: ratqs(klon, llm)
354 real ratqss(klon, llm), ratqsc(klon, llm)
355 real:: ratqsbas = 0.01, ratqshaut = 0.3
356
357 ! Parametres lies au nouveau schema de nuages (SB, PDF)
358 real:: fact_cldcon = 0.375
359 real:: facttemps = 1.e-4
360 logical:: ok_newmicro = .true.
361 real facteur
362
363 integer:: iflag_cldcon = 1
364 logical ptconv(klon, llm)
365
366 ! Variables pour effectuer les appels en s\'erie :
367
368 REAL t_seri(klon, llm), q_seri(klon, llm)
369 REAL ql_seri(klon, llm)
370 REAL u_seri(klon, llm), v_seri(klon, llm)
371 REAL tr_seri(klon, llm, nqmx - 2)
372
373 REAL zx_rh(klon, llm)
374
375 REAL zustrdr(klon), zvstrdr(klon)
376 REAL zustrli(klon), zvstrli(klon)
377 REAL zustrph(klon), zvstrph(klon)
378 REAL aam, torsfc
379
380 REAL ve_lay(klon, llm) ! transport meri. de l'energie a chaque niveau vert.
381 REAL vq_lay(klon, llm) ! transport meri. de l'eau a chaque niveau vert.
382 REAL ue_lay(klon, llm) ! transport zonal de l'energie a chaque niveau vert.
383 REAL uq_lay(klon, llm) ! transport zonal de l'eau a chaque niveau vert.
384
385 real date0
386 REAL ztsol(klon)
387
388 REAL d_t_ec(klon, llm)
389 ! tendance due \`a la conversion Ec en énergie thermique
390
391 REAL ZRCPD
392
393 REAL, save:: t2m(klon, nbsrf), q2m(klon, nbsrf)
394 ! temperature and humidity at 2 m
395
396 REAL, save:: u10m(klon, nbsrf), v10m(klon, nbsrf) ! vents a 10 m
397 REAL zt2m(klon), zq2m(klon) ! température, humidité 2 m moyenne sur 1 maille
398 REAL zu10m(klon), zv10m(klon) ! vents a 10 m moyennes sur 1 maille
399
400 ! Aerosol effects:
401
402 REAL sulfate(klon, llm) ! SO4 aerosol concentration (micro g / m3)
403
404 REAL, save:: sulfate_pi(klon, llm)
405 ! SO4 aerosol concentration, in \mu g / m3, pre-industrial value
406
407 REAL cldtaupi(klon, llm)
408 ! cloud optical thickness for pre-industrial aerosols
409
410 REAL re(klon, llm) ! Cloud droplet effective radius
411 REAL fl(klon, llm) ! denominator of re
412
413 ! Aerosol optical properties
414 REAL, save:: tau_ae(klon, llm, 2), piz_ae(klon, llm, 2)
415 REAL, save:: cg_ae(klon, llm, 2)
416
417 REAL, save:: topswad(klon), solswad(klon) ! aerosol direct effect
418 REAL, save:: topswai(klon), solswai(klon) ! aerosol indirect effect
419
420 LOGICAL:: ok_ade = .false. ! apply aerosol direct effect
421 LOGICAL:: ok_aie = .false. ! apply aerosol indirect effect
422
423 REAL:: bl95_b0 = 2., bl95_b1 = 0.2
424 ! Parameters in equation (D) of Boucher and Lohmann (1995, Tellus
425 ! B). They link cloud droplet number concentration to aerosol mass
426 ! concentration.
427
428 real zmasse(klon, llm)
429 ! (column-density of mass of air in a cell, in kg m-2)
430
431 integer, save:: ncid_startphy
432
433 namelist /physiq_nml/ fact_cldcon, facttemps, ok_newmicro, iflag_cldcon, &
434 ratqsbas, ratqshaut, ok_ade, ok_aie, bl95_b0, bl95_b1, &
435 iflag_thermals, nsplit_thermals
436
437 !----------------------------------------------------------------
438
439 IF (nqmx < 2) CALL abort_gcm('physiq', &
440 'eaux vapeur et liquide sont indispensables')
441
442 test_firstcal: IF (firstcal) THEN
443 ! initialiser
444 u10m = 0.
445 v10m = 0.
446 t2m = 0.
447 q2m = 0.
448 ffonte = 0.
449 fqcalving = 0.
450 piz_ae = 0.
451 tau_ae = 0.
452 cg_ae = 0.
453 rain_con = 0.
454 snow_con = 0.
455 topswai = 0.
456 topswad = 0.
457 solswai = 0.
458 solswad = 0.
459
460 d_u_con = 0.
461 d_v_con = 0.
462 rnebcon0 = 0.
463 clwcon0 = 0.
464 rnebcon = 0.
465 clwcon = 0.
466
467 pblh =0. ! Hauteur de couche limite
468 plcl =0. ! Niveau de condensation de la CLA
469 capCL =0. ! CAPE de couche limite
470 oliqCL =0. ! eau_liqu integree de couche limite
471 cteiCL =0. ! cloud top instab. crit. couche limite
472 pblt =0.
473 therm =0.
474 trmb1 =0. ! deep_cape
475 trmb2 =0. ! inhibition
476 trmb3 =0. ! Point Omega
477
478 iflag_thermals = 0
479 nsplit_thermals = 1
480 print *, "Enter namelist 'physiq_nml'."
481 read(unit=*, nml=physiq_nml)
482 write(unit_nml, nml=physiq_nml)
483
484 call conf_phys
485
486 ! Initialiser les compteurs:
487
488 frugs = 0.
489 CALL phyetat0(pctsrf, ftsol, ftsoil, fqsurf, qsol, fsnow, falbe, &
490 fevap, rain_fall, snow_fall, solsw, sollw, dlw, radsol, frugs, &
491 agesno, zmea, zstd, zsig, zgam, zthe, zpic, zval, t_ancien, &
492 q_ancien, ancien_ok, rnebcon, ratqs, clwcon, run_off_lic_0, sig1, &
493 w01, ncid_startphy)
494
495 ! ATTENTION : il faudra a terme relire q2 dans l'etat initial
496 q2 = 1e-8
497
498 radpas = lmt_pas / nbapp_rad
499 print *, "radpas = ", radpas
500
501 ! Initialisation pour le sch\'ema de convection d'Emanuel :
502 IF (conv_emanuel) THEN
503 ibas_con = 1
504 itop_con = 1
505 ENDIF
506
507 IF (ok_orodr) THEN
508 rugoro = MAX(1e-5, zstd * zsig / 2)
509 CALL SUGWD(paprs, play)
510 else
511 rugoro = 0.
512 ENDIF
513
514 ecrit_ins = NINT(ecrit_ins / dtphys)
515
516 ! Initialisation des sorties
517
518 call ini_histins(dtphys)
519 CALL ymds2ju(annee_ref, 1, day_ref, 0., date0)
520 ! Positionner date0 pour initialisation de ORCHIDEE
521 print *, 'physiq date0: ', date0
522 CALL phyredem0
523 ENDIF test_firstcal
524
525 ! We will modify variables *_seri and we will not touch variables
526 ! u, v, t, qx:
527 t_seri = t
528 u_seri = u
529 v_seri = v
530 q_seri = qx(:, :, ivap)
531 ql_seri = qx(:, :, iliq)
532 tr_seri = qx(:, :, 3:nqmx)
533
534 ztsol = sum(ftsol * pctsrf, dim = 2)
535
536 ! Diagnostic de la tendance dynamique :
537 IF (ancien_ok) THEN
538 DO k = 1, llm
539 DO i = 1, klon
540 d_t_dyn(i, k) = (t_seri(i, k) - t_ancien(i, k)) / dtphys
541 d_q_dyn(i, k) = (q_seri(i, k) - q_ancien(i, k)) / dtphys
542 ENDDO
543 ENDDO
544 ELSE
545 DO k = 1, llm
546 DO i = 1, klon
547 d_t_dyn(i, k) = 0.
548 d_q_dyn(i, k) = 0.
549 ENDDO
550 ENDDO
551 ancien_ok = .TRUE.
552 ENDIF
553
554 ! Ajouter le geopotentiel du sol:
555 DO k = 1, llm
556 DO i = 1, klon
557 zphi(i, k) = pphi(i, k) + pphis(i)
558 ENDDO
559 ENDDO
560
561 ! Check temperatures:
562 CALL hgardfou(t_seri, ftsol)
563
564 call increment_itap
565 julien = MOD(dayvrai, 360)
566 if (julien == 0) julien = 360
567
568 forall (k = 1: llm) zmasse(:, k) = (paprs(:, k) - paprs(:, k + 1)) / rg
569
570 ! Prescrire l'ozone :
571 wo = ozonecm(REAL(julien), paprs)
572
573 ! \'Evaporation de l'eau liquide nuageuse :
574 DO k = 1, llm
575 DO i = 1, klon
576 zb = MAX(0., ql_seri(i, k))
577 t_seri(i, k) = t_seri(i, k) &
578 - zb * RLVTT / RCPD / (1. + RVTMP2 * q_seri(i, k))
579 q_seri(i, k) = q_seri(i, k) + zb
580 ENDDO
581 ENDDO
582 ql_seri = 0.
583
584 frugs = MAX(frugs, 0.000015)
585 zxrugs = sum(frugs * pctsrf, dim = 2)
586
587 ! Calculs n\'ecessaires au calcul de l'albedo dans l'interface avec
588 ! la surface.
589
590 CALL orbite(REAL(julien), longi, dist)
591 CALL zenang(longi, time, dtphys * radpas, mu0, fract)
592
593 ! Calcul de l'abedo moyen par maille
594 albsol = sum(falbe * pctsrf, dim = 2)
595
596 ! R\'epartition sous maille des flux longwave et shortwave
597 ! R\'epartition du longwave par sous-surface lin\'earis\'ee
598
599 forall (nsrf = 1: nbsrf)
600 fsollw(:, nsrf) = sollw + 4. * RSIGMA * ztsol**3 &
601 * (ztsol - ftsol(:, nsrf))
602 fsolsw(:, nsrf) = solsw * (1. - falbe(:, nsrf)) / (1. - albsol)
603 END forall
604
605 fder = dlw
606
607 CALL clmain(dtphys, pctsrf, t_seri, q_seri, u_seri, v_seri, julien, mu0, &
608 ftsol, cdmmax, cdhmax, ksta, ksta_ter, ok_kzmin, ftsoil, qsol, &
609 paprs, play, fsnow, fqsurf, fevap, falbe, fluxlat, rain_fall, &
610 snow_fall, fsolsw, fsollw, fder, frugs, agesno, rugoro, d_t_vdf, &
611 d_q_vdf, d_u_vdf, d_v_vdf, d_ts, flux_t, flux_q, flux_u, flux_v, &
612 cdragh, cdragm, q2, dsens, devap, ycoefh, yu1, yv1, t2m, q2m, u10m, &
613 v10m, pblh, capCL, oliqCL, cteiCL, pblT, therm, trmb1, trmb2, trmb3, &
614 plcl, fqcalving, ffonte, run_off_lic_0)
615
616 ! Incr\'ementation des flux
617
618 sens = - sum(flux_t * pctsrf, dim = 2)
619 evap = - sum(flux_q * pctsrf, dim = 2)
620 fder = dlw + dsens + devap
621
622 DO k = 1, llm
623 DO i = 1, klon
624 t_seri(i, k) = t_seri(i, k) + d_t_vdf(i, k)
625 q_seri(i, k) = q_seri(i, k) + d_q_vdf(i, k)
626 u_seri(i, k) = u_seri(i, k) + d_u_vdf(i, k)
627 v_seri(i, k) = v_seri(i, k) + d_v_vdf(i, k)
628 ENDDO
629 ENDDO
630
631 ! Update surface temperature:
632
633 call assert(abs(sum(pctsrf, dim = 2) - 1.) <= EPSFRA, 'physiq: pctsrf')
634 ftsol = ftsol + d_ts
635 ztsol = sum(ftsol * pctsrf, dim = 2)
636 zxfluxlat = sum(fluxlat * pctsrf, dim = 2)
637 zt2m = sum(t2m * pctsrf, dim = 2)
638 zq2m = sum(q2m * pctsrf, dim = 2)
639 zu10m = sum(u10m * pctsrf, dim = 2)
640 zv10m = sum(v10m * pctsrf, dim = 2)
641 zxffonte = sum(ffonte * pctsrf, dim = 2)
642 zxfqcalving = sum(fqcalving * pctsrf, dim = 2)
643 s_pblh = sum(pblh * pctsrf, dim = 2)
644 s_lcl = sum(plcl * pctsrf, dim = 2)
645 s_capCL = sum(capCL * pctsrf, dim = 2)
646 s_oliqCL = sum(oliqCL * pctsrf, dim = 2)
647 s_cteiCL = sum(cteiCL * pctsrf, dim = 2)
648 s_pblT = sum(pblT * pctsrf, dim = 2)
649 s_therm = sum(therm * pctsrf, dim = 2)
650 s_trmb1 = sum(trmb1 * pctsrf, dim = 2)
651 s_trmb2 = sum(trmb2 * pctsrf, dim = 2)
652 s_trmb3 = sum(trmb3 * pctsrf, dim = 2)
653
654 ! Si une sous-fraction n'existe pas, elle prend la valeur moyenne :
655 DO nsrf = 1, nbsrf
656 DO i = 1, klon
657 IF (pctsrf(i, nsrf) < epsfra) then
658 ftsol(i, nsrf) = ztsol(i)
659 t2m(i, nsrf) = zt2m(i)
660 q2m(i, nsrf) = zq2m(i)
661 u10m(i, nsrf) = zu10m(i)
662 v10m(i, nsrf) = zv10m(i)
663 ffonte(i, nsrf) = zxffonte(i)
664 fqcalving(i, nsrf) = zxfqcalving(i)
665 pblh(i, nsrf) = s_pblh(i)
666 plcl(i, nsrf) = s_lcl(i)
667 capCL(i, nsrf) = s_capCL(i)
668 oliqCL(i, nsrf) = s_oliqCL(i)
669 cteiCL(i, nsrf) = s_cteiCL(i)
670 pblT(i, nsrf) = s_pblT(i)
671 therm(i, nsrf) = s_therm(i)
672 trmb1(i, nsrf) = s_trmb1(i)
673 trmb2(i, nsrf) = s_trmb2(i)
674 trmb3(i, nsrf) = s_trmb3(i)
675 end IF
676 ENDDO
677 ENDDO
678
679 ! Calculer la dérive du flux infrarouge
680
681 DO i = 1, klon
682 dlw(i) = - 4. * RSIGMA * ztsol(i)**3
683 ENDDO
684
685 ! Appeler la convection
686
687 if (conv_emanuel) then
688 CALL concvl(paprs, play, t_seri, q_seri, u_seri, v_seri, sig1, w01, &
689 d_t_con, d_q_con, d_u_con, d_v_con, rain_con, ibas_con, itop_con, &
690 upwd, dnwd, Ma, cape, iflagctrl, qcondc, pmflxr, da, phi, mp)
691 snow_con = 0.
692 clwcon0 = qcondc
693 mfu = upwd + dnwd
694
695 zqsat = MIN(0.5, r2es * FOEEW(t_seri, rtt >= t_seri) / play)
696 zqsat = zqsat / (1. - retv * zqsat)
697
698 ! Properties of convective clouds
699 clwcon0 = fact_cldcon * clwcon0
700 call clouds_gno(klon, llm, q_seri, zqsat, clwcon0, ptconv, ratqsc, &
701 rnebcon0)
702
703 forall (i = 1:klon) ema_pct(i) = paprs(i, itop_con(i) + 1)
704 mfd = 0.
705 pen_u = 0.
706 pen_d = 0.
707 pde_d = 0.
708 pde_u = 0.
709 else
710 conv_q = d_q_dyn + d_q_vdf / dtphys
711 conv_t = d_t_dyn + d_t_vdf / dtphys
712 z_avant = sum((q_seri + ql_seri) * zmasse, dim=2)
713 CALL conflx(dtphys, paprs, play, t_seri(:, llm:1:- 1), &
714 q_seri(:, llm:1:- 1), conv_t, conv_q, - evap, omega, &
715 d_t_con, d_q_con, rain_con, snow_con, mfu(:, llm:1:- 1), &
716 mfd(:, llm:1:- 1), pen_u, pde_u, pen_d, pde_d, kcbot, kctop, &
717 kdtop, pmflxr, pmflxs)
718 WHERE (rain_con < 0.) rain_con = 0.
719 WHERE (snow_con < 0.) snow_con = 0.
720 ibas_con = llm + 1 - kcbot
721 itop_con = llm + 1 - kctop
722 END if
723
724 DO k = 1, llm
725 DO i = 1, klon
726 t_seri(i, k) = t_seri(i, k) + d_t_con(i, k)
727 q_seri(i, k) = q_seri(i, k) + d_q_con(i, k)
728 u_seri(i, k) = u_seri(i, k) + d_u_con(i, k)
729 v_seri(i, k) = v_seri(i, k) + d_v_con(i, k)
730 ENDDO
731 ENDDO
732
733 IF (.not. conv_emanuel) THEN
734 z_apres = sum((q_seri + ql_seri) * zmasse, dim=2)
735 z_factor = (z_avant - (rain_con + snow_con) * dtphys) / z_apres
736 DO k = 1, llm
737 DO i = 1, klon
738 IF (z_factor(i) > 1. + 1E-8 .OR. z_factor(i) < 1. - 1E-8) THEN
739 q_seri(i, k) = q_seri(i, k) * z_factor(i)
740 ENDIF
741 ENDDO
742 ENDDO
743 ENDIF
744
745 ! Convection s\`eche (thermiques ou ajustement)
746
747 d_t_ajs = 0.
748 d_u_ajs = 0.
749 d_v_ajs = 0.
750 d_q_ajs = 0.
751 fm_therm = 0.
752 entr_therm = 0.
753
754 if (iflag_thermals == 0) then
755 ! Ajustement sec
756 CALL ajsec(paprs, play, t_seri, q_seri, d_t_ajs, d_q_ajs)
757 t_seri = t_seri + d_t_ajs
758 q_seri = q_seri + d_q_ajs
759 else
760 call calltherm(dtphys, play, paprs, pphi, u_seri, v_seri, t_seri, &
761 q_seri, d_u_ajs, d_v_ajs, d_t_ajs, d_q_ajs, fm_therm, entr_therm)
762 endif
763
764 ! Caclul des ratqs
765
766 ! ratqs convectifs \`a l'ancienne en fonction de (q(z = 0) - q) / q
767 ! on \'ecrase le tableau ratqsc calcul\'e par clouds_gno
768 if (iflag_cldcon == 1) then
769 do k = 1, llm
770 do i = 1, klon
771 if(ptconv(i, k)) then
772 ratqsc(i, k) = ratqsbas + fact_cldcon &
773 * (q_seri(i, 1) - q_seri(i, k)) / q_seri(i, k)
774 else
775 ratqsc(i, k) = 0.
776 endif
777 enddo
778 enddo
779 endif
780
781 ! ratqs stables
782 do k = 1, llm
783 do i = 1, klon
784 ratqss(i, k) = ratqsbas + (ratqshaut - ratqsbas) &
785 * min((paprs(i, 1) - play(i, k)) / (paprs(i, 1) - 3e4), 1.)
786 enddo
787 enddo
788
789 ! ratqs final
790 if (iflag_cldcon == 1 .or. iflag_cldcon == 2) then
791 ! les ratqs sont une conbinaison de ratqss et ratqsc
792 ! ratqs final
793 ! 1e4 (en gros 3 heures), en dur pour le moment, est le temps de
794 ! relaxation des ratqs
795 ratqs = max(ratqs * exp(- dtphys * facttemps), ratqss)
796 ratqs = max(ratqs, ratqsc)
797 else
798 ! on ne prend que le ratqs stable pour fisrtilp
799 ratqs = ratqss
800 endif
801
802 CALL fisrtilp(dtphys, paprs, play, t_seri, q_seri, ptconv, ratqs, &
803 d_t_lsc, d_q_lsc, d_ql_lsc, rneb, cldliq, rain_lsc, snow_lsc, &
804 pfrac_impa, pfrac_nucl, pfrac_1nucl, frac_impa, frac_nucl, prfl, &
805 psfl, rhcl)
806
807 WHERE (rain_lsc < 0) rain_lsc = 0.
808 WHERE (snow_lsc < 0) snow_lsc = 0.
809 DO k = 1, llm
810 DO i = 1, klon
811 t_seri(i, k) = t_seri(i, k) + d_t_lsc(i, k)
812 q_seri(i, k) = q_seri(i, k) + d_q_lsc(i, k)
813 ql_seri(i, k) = ql_seri(i, k) + d_ql_lsc(i, k)
814 cldfra(i, k) = rneb(i, k)
815 IF (.NOT.new_oliq) cldliq(i, k) = ql_seri(i, k)
816 ENDDO
817 ENDDO
818
819 ! PRESCRIPTION DES NUAGES POUR LE RAYONNEMENT
820
821 ! 1. NUAGES CONVECTIFS
822
823 IF (iflag_cldcon <= - 1) THEN
824 ! seulement pour Tiedtke
825 snow_tiedtke = 0.
826 if (iflag_cldcon == - 1) then
827 rain_tiedtke = rain_con
828 else
829 rain_tiedtke = 0.
830 do k = 1, llm
831 do i = 1, klon
832 if (d_q_con(i, k) < 0.) then
833 rain_tiedtke(i) = rain_tiedtke(i) - d_q_con(i, k) / dtphys &
834 * zmasse(i, k)
835 endif
836 enddo
837 enddo
838 endif
839
840 ! Nuages diagnostiques pour Tiedtke
841 CALL diagcld1(paprs, play, rain_tiedtke, snow_tiedtke, ibas_con, &
842 itop_con, diafra, dialiq)
843 DO k = 1, llm
844 DO i = 1, klon
845 IF (diafra(i, k) > cldfra(i, k)) THEN
846 cldliq(i, k) = dialiq(i, k)
847 cldfra(i, k) = diafra(i, k)
848 ENDIF
849 ENDDO
850 ENDDO
851 ELSE IF (iflag_cldcon == 3) THEN
852 ! On prend pour les nuages convectifs le maximum du calcul de
853 ! la convection et du calcul du pas de temps pr\'ec\'edent diminu\'e
854 ! d'un facteur facttemps.
855 facteur = dtphys * facttemps
856 do k = 1, llm
857 do i = 1, klon
858 rnebcon(i, k) = rnebcon(i, k) * facteur
859 if (rnebcon0(i, k) * clwcon0(i, k) &
860 > rnebcon(i, k) * clwcon(i, k)) then
861 rnebcon(i, k) = rnebcon0(i, k)
862 clwcon(i, k) = clwcon0(i, k)
863 endif
864 enddo
865 enddo
866
867 ! On prend la somme des fractions nuageuses et des contenus en eau
868 cldfra = min(max(cldfra, rnebcon), 1.)
869 cldliq = cldliq + rnebcon * clwcon
870 ENDIF
871
872 ! 2. Nuages stratiformes
873
874 IF (ok_stratus) THEN
875 CALL diagcld2(paprs, play, t_seri, q_seri, diafra, dialiq)
876 DO k = 1, llm
877 DO i = 1, klon
878 IF (diafra(i, k) > cldfra(i, k)) THEN
879 cldliq(i, k) = dialiq(i, k)
880 cldfra(i, k) = diafra(i, k)
881 ENDIF
882 ENDDO
883 ENDDO
884 ENDIF
885
886 ! Precipitation totale
887 DO i = 1, klon
888 rain_fall(i) = rain_con(i) + rain_lsc(i)
889 snow_fall(i) = snow_con(i) + snow_lsc(i)
890 ENDDO
891
892 ! Humidit\'e relative pour diagnostic :
893 DO k = 1, llm
894 DO i = 1, klon
895 zx_t = t_seri(i, k)
896 zx_qs = r2es * FOEEW(zx_t, rtt >= zx_t) / play(i, k)
897 zx_qs = MIN(0.5, zx_qs)
898 zcor = 1. / (1. - retv * zx_qs)
899 zx_qs = zx_qs * zcor
900 zx_rh(i, k) = q_seri(i, k) / zx_qs
901 zqsat(i, k) = zx_qs
902 ENDDO
903 ENDDO
904
905 ! Introduce the aerosol direct and first indirect radiative forcings:
906 tau_ae = 0.
907 piz_ae = 0.
908 cg_ae = 0.
909
910 ! Param\`etres optiques des nuages et quelques param\`etres pour
911 ! diagnostics :
912 if (ok_newmicro) then
913 CALL newmicro(paprs, play, t_seri, cldliq, cldfra, cldtau, cldemi, &
914 cldh, cldl, cldm, cldt, cldq, flwp, fiwp, flwc, fiwc, ok_aie, &
915 sulfate, sulfate_pi, bl95_b0, bl95_b1, cldtaupi, re, fl)
916 else
917 CALL nuage(paprs, play, t_seri, cldliq, cldfra, cldtau, cldemi, cldh, &
918 cldl, cldm, cldt, cldq, ok_aie, sulfate, sulfate_pi, bl95_b0, &
919 bl95_b1, cldtaupi, re, fl)
920 endif
921
922 IF (MOD(itap - 1, radpas) == 0) THEN
923 ! Appeler le rayonnement mais calculer tout d'abord l'albedo du sol.
924 ! Calcul de l'abedo moyen par maille
925 albsol = sum(falbe * pctsrf, dim = 2)
926
927 ! Rayonnement (compatible Arpege-IFS) :
928 CALL radlwsw(dist, mu0, fract, paprs, play, ztsol, albsol, t_seri, &
929 q_seri, wo, cldfra, cldemi, cldtau, heat, heat0, cool, cool0, &
930 radsol, albpla, topsw, toplw, solsw, sollw, sollwdown, topsw0, &
931 toplw0, solsw0, sollw0, lwdn0, lwdn, lwup0, lwup, swdn0, swdn, &
932 swup0, swup, ok_ade, ok_aie, tau_ae, piz_ae, cg_ae, topswad, &
933 solswad, cldtaupi, topswai, solswai)
934 ENDIF
935
936 ! Ajouter la tendance des rayonnements (tous les pas)
937 DO k = 1, llm
938 DO i = 1, klon
939 t_seri(i, k) = t_seri(i, k) + (heat(i, k) - cool(i, k)) * dtphys &
940 / 86400.
941 ENDDO
942 ENDDO
943
944 ! Calculer l'hydrologie de la surface
945 zxqsurf = sum(fqsurf * pctsrf, dim = 2)
946 zxsnow = sum(fsnow * pctsrf, dim = 2)
947
948 ! Calculer le bilan du sol et la d\'erive de temp\'erature (couplage)
949 DO i = 1, klon
950 bils(i) = radsol(i) - sens(i) + zxfluxlat(i)
951 ENDDO
952
953 ! Param\'etrisation de l'orographie \`a l'\'echelle sous-maille :
954
955 IF (ok_orodr) THEN
956 ! S\'election des points pour lesquels le sch\'ema est actif :
957 igwd = 0
958 DO i = 1, klon
959 itest(i) = 0
960 IF (zpic(i) - zmea(i) > 100. .AND. zstd(i) > 10.) THEN
961 itest(i) = 1
962 igwd = igwd + 1
963 ENDIF
964 ENDDO
965
966 CALL drag_noro(klon, llm, dtphys, paprs, play, zmea, zstd, zsig, zgam, &
967 zthe, zpic, zval, itest, t_seri, u_seri, v_seri, zulow, zvlow, &
968 zustrdr, zvstrdr, d_t_oro, d_u_oro, d_v_oro)
969
970 ! ajout des tendances
971 DO k = 1, llm
972 DO i = 1, klon
973 t_seri(i, k) = t_seri(i, k) + d_t_oro(i, k)
974 u_seri(i, k) = u_seri(i, k) + d_u_oro(i, k)
975 v_seri(i, k) = v_seri(i, k) + d_v_oro(i, k)
976 ENDDO
977 ENDDO
978 ENDIF
979
980 IF (ok_orolf) THEN
981 ! S\'election des points pour lesquels le sch\'ema est actif :
982 igwd = 0
983 DO i = 1, klon
984 itest(i) = 0
985 IF (zpic(i) - zmea(i) > 100.) THEN
986 itest(i) = 1
987 igwd = igwd + 1
988 ENDIF
989 ENDDO
990
991 CALL lift_noro(klon, llm, dtphys, paprs, play, rlat, zmea, zstd, zpic, &
992 itest, t_seri, u_seri, v_seri, zulow, zvlow, zustrli, zvstrli, &
993 d_t_lif, d_u_lif, d_v_lif)
994
995 ! Ajout des tendances :
996 DO k = 1, llm
997 DO i = 1, klon
998 t_seri(i, k) = t_seri(i, k) + d_t_lif(i, k)
999 u_seri(i, k) = u_seri(i, k) + d_u_lif(i, k)
1000 v_seri(i, k) = v_seri(i, k) + d_v_lif(i, k)
1001 ENDDO
1002 ENDDO
1003 ENDIF
1004
1005 ! Stress n\'ecessaires : toute la physique
1006
1007 DO i = 1, klon
1008 zustrph(i) = 0.
1009 zvstrph(i) = 0.
1010 ENDDO
1011 DO k = 1, llm
1012 DO i = 1, klon
1013 zustrph(i) = zustrph(i) + (u_seri(i, k) - u(i, k)) / dtphys &
1014 * zmasse(i, k)
1015 zvstrph(i) = zvstrph(i) + (v_seri(i, k) - v(i, k)) / dtphys &
1016 * zmasse(i, k)
1017 ENDDO
1018 ENDDO
1019
1020 CALL aaam_bud(rg, romega, rlat, rlon, pphis, zustrdr, zustrli, zustrph, &
1021 zvstrdr, zvstrli, zvstrph, paprs, u, v, aam, torsfc)
1022
1023 ! Calcul des tendances traceurs
1024 call phytrac(julien, time, firstcal, lafin, dtphys, t, paprs, play, mfu, &
1025 mfd, pde_u, pen_d, ycoefh, fm_therm, entr_therm, yu1, yv1, ftsol, &
1026 pctsrf, frac_impa, frac_nucl, da, phi, mp, upwd, dnwd, tr_seri, &
1027 zmasse, ncid_startphy)
1028
1029 IF (offline) call phystokenc(dtphys, t, mfu, mfd, pen_u, pde_u, pen_d, &
1030 pde_d, fm_therm, entr_therm, ycoefh, yu1, yv1, ftsol, pctsrf, &
1031 frac_impa, frac_nucl, pphis, airephy, dtphys)
1032
1033 ! Calculer le transport de l'eau et de l'energie (diagnostique)
1034 CALL transp(paprs, t_seri, q_seri, u_seri, v_seri, zphi, ve, vq, ue, uq)
1035
1036 ! diag. bilKP
1037
1038 CALL transp_lay(paprs, t_seri, q_seri, u_seri, v_seri, zphi, &
1039 ve_lay, vq_lay, ue_lay, uq_lay)
1040
1041 ! Accumuler les variables a stocker dans les fichiers histoire:
1042
1043 ! conversion Ec en énergie thermique
1044 DO k = 1, llm
1045 DO i = 1, klon
1046 ZRCPD = RCPD * (1. + RVTMP2 * q_seri(i, k))
1047 d_t_ec(i, k) = 0.5 / ZRCPD &
1048 * (u(i, k)**2 + v(i, k)**2 - u_seri(i, k)**2 - v_seri(i, k)**2)
1049 t_seri(i, k) = t_seri(i, k) + d_t_ec(i, k)
1050 d_t_ec(i, k) = d_t_ec(i, k) / dtphys
1051 END DO
1052 END DO
1053
1054 ! SORTIES
1055
1056 ! prw = eau precipitable
1057 DO i = 1, klon
1058 prw(i) = 0.
1059 DO k = 1, llm
1060 prw(i) = prw(i) + q_seri(i, k) * zmasse(i, k)
1061 ENDDO
1062 ENDDO
1063
1064 ! Convertir les incrementations en tendances
1065
1066 DO k = 1, llm
1067 DO i = 1, klon
1068 d_u(i, k) = (u_seri(i, k) - u(i, k)) / dtphys
1069 d_v(i, k) = (v_seri(i, k) - v(i, k)) / dtphys
1070 d_t(i, k) = (t_seri(i, k) - t(i, k)) / dtphys
1071 d_qx(i, k, ivap) = (q_seri(i, k) - qx(i, k, ivap)) / dtphys
1072 d_qx(i, k, iliq) = (ql_seri(i, k) - qx(i, k, iliq)) / dtphys
1073 ENDDO
1074 ENDDO
1075
1076 DO iq = 3, nqmx
1077 DO k = 1, llm
1078 DO i = 1, klon
1079 d_qx(i, k, iq) = (tr_seri(i, k, iq - 2) - qx(i, k, iq)) / dtphys
1080 ENDDO
1081 ENDDO
1082 ENDDO
1083
1084 ! Sauvegarder les valeurs de t et q a la fin de la physique:
1085 DO k = 1, llm
1086 DO i = 1, klon
1087 t_ancien(i, k) = t_seri(i, k)
1088 q_ancien(i, k) = q_seri(i, k)
1089 ENDDO
1090 ENDDO
1091
1092 CALL histwrite_phy("phis", pphis)
1093 CALL histwrite_phy("aire", airephy)
1094 CALL histwrite_phy("psol", paprs(:, 1))
1095 CALL histwrite_phy("precip", rain_fall + snow_fall)
1096 CALL histwrite_phy("plul", rain_lsc + snow_lsc)
1097 CALL histwrite_phy("pluc", rain_con + snow_con)
1098 CALL histwrite_phy("tsol", ztsol)
1099 CALL histwrite_phy("t2m", zt2m)
1100 CALL histwrite_phy("q2m", zq2m)
1101 CALL histwrite_phy("u10m", zu10m)
1102 CALL histwrite_phy("v10m", zv10m)
1103 CALL histwrite_phy("snow", snow_fall)
1104 CALL histwrite_phy("cdrm", cdragm)
1105 CALL histwrite_phy("cdrh", cdragh)
1106 CALL histwrite_phy("topl", toplw)
1107 CALL histwrite_phy("evap", evap)
1108 CALL histwrite_phy("sols", solsw)
1109 CALL histwrite_phy("soll", sollw)
1110 CALL histwrite_phy("solldown", sollwdown)
1111 CALL histwrite_phy("bils", bils)
1112 CALL histwrite_phy("sens", - sens)
1113 CALL histwrite_phy("fder", fder)
1114 CALL histwrite_phy("dtsvdfo", d_ts(:, is_oce))
1115 CALL histwrite_phy("dtsvdft", d_ts(:, is_ter))
1116 CALL histwrite_phy("dtsvdfg", d_ts(:, is_lic))
1117 CALL histwrite_phy("dtsvdfi", d_ts(:, is_sic))
1118
1119 DO nsrf = 1, nbsrf
1120 CALL histwrite_phy("pourc_"//clnsurf(nsrf), pctsrf(:, nsrf) * 100.)
1121 CALL histwrite_phy("fract_"//clnsurf(nsrf), pctsrf(:, nsrf))
1122 CALL histwrite_phy("sens_"//clnsurf(nsrf), flux_t(:, nsrf))
1123 CALL histwrite_phy("lat_"//clnsurf(nsrf), fluxlat(:, nsrf))
1124 CALL histwrite_phy("tsol_"//clnsurf(nsrf), ftsol(:, nsrf))
1125 CALL histwrite_phy("taux_"//clnsurf(nsrf), flux_u(:, nsrf))
1126 CALL histwrite_phy("tauy_"//clnsurf(nsrf), flux_v(:, nsrf))
1127 CALL histwrite_phy("rugs_"//clnsurf(nsrf), frugs(:, nsrf))
1128 CALL histwrite_phy("albe_"//clnsurf(nsrf), falbe(:, nsrf))
1129 END DO
1130
1131 CALL histwrite_phy("albs", albsol)
1132 CALL histwrite_phy("rugs", zxrugs)
1133 CALL histwrite_phy("s_pblh", s_pblh)
1134 CALL histwrite_phy("s_pblt", s_pblt)
1135 CALL histwrite_phy("s_lcl", s_lcl)
1136 CALL histwrite_phy("s_capCL", s_capCL)
1137 CALL histwrite_phy("s_oliqCL", s_oliqCL)
1138 CALL histwrite_phy("s_cteiCL", s_cteiCL)
1139 CALL histwrite_phy("s_therm", s_therm)
1140 CALL histwrite_phy("s_trmb1", s_trmb1)
1141 CALL histwrite_phy("s_trmb2", s_trmb2)
1142 CALL histwrite_phy("s_trmb3", s_trmb3)
1143
1144 if (conv_emanuel) then
1145 CALL histwrite_phy("ptop", ema_pct)
1146 CALL histwrite_phy("dnwd0", - mp)
1147 end if
1148
1149 CALL histwrite_phy("temp", t_seri)
1150 CALL histwrite_phy("vitu", u_seri)
1151 CALL histwrite_phy("vitv", v_seri)
1152 CALL histwrite_phy("geop", zphi)
1153 CALL histwrite_phy("pres", play)
1154 CALL histwrite_phy("dtvdf", d_t_vdf)
1155 CALL histwrite_phy("dqvdf", d_q_vdf)
1156 CALL histwrite_phy("rhum", zx_rh)
1157
1158 if (ok_instan) call histsync(nid_ins)
1159
1160 IF (lafin) then
1161 call NF95_CLOSE(ncid_startphy)
1162 CALL phyredem(pctsrf, ftsol, ftsoil, fqsurf, qsol, &
1163 fsnow, falbe, fevap, rain_fall, snow_fall, solsw, sollw, dlw, &
1164 radsol, frugs, agesno, zmea, zstd, zsig, zgam, zthe, zpic, zval, &
1165 t_ancien, q_ancien, rnebcon, ratqs, clwcon, run_off_lic_0, sig1, &
1166 w01)
1167 end IF
1168
1169 firstcal = .FALSE.
1170
1171 END SUBROUTINE physiq
1172
1173 end module physiq_m

  ViewVC Help
Powered by ViewVC 1.1.21