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

Contents of /trunk/phylmd/physiq.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 90 - (show annotations)
Wed Mar 12 21:16:36 2014 UTC (10 years, 2 months ago) by guez
File size: 66329 byte(s)
Removed procedures ini_histday, ini_histhf, write_histday and
write_histhf.

Divided file regr_pr_coefoz.f into regr_pr_av.f and
regr_pr_int.f. (Following LMDZ.) Divided module regr_pr_coefoz into
modules regr_pr_av_m and regr_pr_int_m. Renamed regr_pr_av_coefoz to
regr_pr_av and regr_pr_int_coefoz to regr_pr_int. The idea is that
those procedures are more general than Mobidic.

Removed argument dudyn of calfis and physiq. dudyn is not used either
in LMDZ. Removed computation in calfis of unused variable zpsrf (not
used either in LMDZ). Removed useless computation of dqfi in calfis
(part 62): the results were overwritten. (Same in LMDZ.)

1 module physiq_m
2
3 IMPLICIT none
4
5 contains
6
7 SUBROUTINE physiq(lafin, rdayvrai, time, dtphys, paprs, play, pphi, pphis, &
8 u, v, t, qx, omega, d_u, d_v, d_t, d_qx, d_ps)
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 aeropt_m, only: aeropt
20 use ajsec_m, only: ajsec
21 USE calendar, ONLY: ymds2ju
22 use calltherm_m, only: calltherm
23 USE clesphys, ONLY: cdhmax, cdmmax, co2_ppm, ecrit_hf, ecrit_ins, &
24 ecrit_mth, ecrit_reg, ecrit_tra, ksta, ksta_ter, ok_kzmin
25 USE clesphys2, ONLY: cycle_diurne, iflag_con, nbapp_rad, new_oliq, &
26 ok_orodr, ok_orolf, soil_model
27 USE clmain_m, ONLY: clmain
28 use clouds_gno_m, only: clouds_gno
29 USE comgeomphy, ONLY: airephy, cuphy, cvphy
30 USE concvl_m, ONLY: concvl
31 USE conf_gcm_m, ONLY: offline, raz_date
32 USE conf_phys_m, ONLY: conf_phys
33 use conflx_m, only: conflx
34 USE ctherm, ONLY: iflag_thermals, nsplit_thermals
35 use diagcld2_m, only: diagcld2
36 use diagetpq_m, only: diagetpq
37 use diagphy_m, only: diagphy
38 USE dimens_m, ONLY: llm, nqmx
39 USE dimphy, ONLY: klon, nbtr
40 USE dimsoil, ONLY: nsoilmx
41 use drag_noro_m, only: drag_noro
42 USE fcttre, ONLY: foeew, qsatl, qsats, thermcep
43 use fisrtilp_m, only: fisrtilp
44 USE hgardfou_m, ONLY: hgardfou
45 USE indicesol, ONLY: clnsurf, epsfra, is_lic, is_oce, is_sic, is_ter, &
46 nbsrf
47 USE ini_histins_m, ONLY: ini_histins
48 use newmicro_m, only: newmicro
49 USE oasis_m, ONLY: ok_oasis
50 USE orbite_m, ONLY: orbite, zenang
51 USE ozonecm_m, ONLY: ozonecm
52 USE phyetat0_m, ONLY: phyetat0, rlat, rlon
53 USE phyredem_m, ONLY: phyredem
54 USE phystokenc_m, ONLY: phystokenc
55 USE phytrac_m, ONLY: phytrac
56 USE qcheck_m, ONLY: qcheck
57 use radlwsw_m, only: radlwsw
58 use readsulfate_m, only: readsulfate
59 use sugwd_m, only: sugwd
60 USE suphec_m, ONLY: ra, rcpd, retv, rg, rlvtt, romega, rsigma, rtt
61 USE temps, ONLY: annee_ref, day_ref, itau_phy
62 use unit_nml_m, only: unit_nml
63 USE yoethf_m, ONLY: r2es, rvtmp2
64
65 ! Arguments:
66
67 REAL, intent(in):: rdayvrai
68 ! (elapsed time since January 1st 0h of the starting year, in days)
69
70 REAL, intent(in):: time ! heure de la journ\'ee en fraction de jour
71 REAL, intent(in):: dtphys ! pas d'integration pour la physique (seconde)
72 logical, intent(in):: lafin ! dernier passage
73
74 REAL, intent(in):: paprs(klon, llm + 1)
75 ! (pression pour chaque inter-couche, en Pa)
76
77 REAL, intent(in):: play(klon, llm)
78 ! (input pression pour le mileu de chaque couche (en Pa))
79
80 REAL, intent(in):: pphi(klon, llm)
81 ! (input geopotentiel de chaque couche (g z) (reference sol))
82
83 REAL, intent(in):: pphis(klon) ! input geopotentiel du sol
84
85 REAL, intent(in):: u(klon, llm)
86 ! vitesse dans la direction X (de O a E) en m/s
87
88 REAL, intent(in):: v(klon, llm) ! vitesse Y (de S a N) en m/s
89 REAL, intent(in):: t(klon, llm) ! input temperature (K)
90
91 REAL, intent(in):: qx(klon, llm, nqmx)
92 ! (humidit\'e sp\'ecifique et fractions massiques des autres traceurs)
93
94 REAL omega(klon, llm) ! input vitesse verticale en Pa/s
95 REAL, intent(out):: d_u(klon, llm) ! tendance physique de "u" (m/s/s)
96 REAL, intent(out):: d_v(klon, llm) ! tendance physique de "v" (m/s/s)
97 REAL, intent(out):: d_t(klon, llm) ! tendance physique de "t" (K/s)
98 REAL d_qx(klon, llm, nqmx) ! output tendance physique de "qx" (kg/kg/s)
99 REAL d_ps(klon) ! output tendance physique de la pression au sol
100
101 LOGICAL:: firstcal = .true.
102
103 INTEGER nbteta
104 PARAMETER(nbteta = 3)
105
106 LOGICAL ok_gust ! pour activer l'effet des gust sur flux surface
107 PARAMETER (ok_gust = .FALSE.)
108
109 LOGICAL check ! Verifier la conservation du modele en eau
110 PARAMETER (check = .FALSE.)
111
112 LOGICAL, PARAMETER:: ok_stratus = .FALSE.
113 ! Ajouter artificiellement les stratus
114
115 ! Parametres lies au coupleur OASIS:
116 INTEGER, SAVE:: npas, nexca
117 logical rnpb
118 parameter(rnpb = .true.)
119
120 character(len = 6):: ocean = 'force '
121 ! (type de mod\`ele oc\'ean \`a utiliser: "force" ou "slab" mais pas "couple")
122
123 ! "slab" ocean
124 REAL, save:: tslab(klon) ! temperature of ocean slab
125 REAL, save:: seaice(klon) ! glace de mer (kg/m2)
126 REAL fluxo(klon) ! flux turbulents ocean-glace de mer
127 REAL fluxg(klon) ! flux turbulents ocean-atmosphere
128
129 ! Modele thermique du sol, a activer pour le cycle diurne:
130 logical:: ok_veget = .false. ! type de modele de vegetation utilise
131
132 logical:: ok_journe = .false., ok_mensuel = .true., ok_instan = .false.
133 ! sorties journalieres, mensuelles et instantanees dans les
134 ! fichiers histday, histmth et histins
135
136 LOGICAL ok_region ! sortir le fichier regional
137 PARAMETER (ok_region = .FALSE.)
138
139 ! pour phsystoke avec thermiques
140 REAL fm_therm(klon, llm + 1)
141 REAL entr_therm(klon, llm)
142 real, save:: q2(klon, llm + 1, nbsrf)
143
144 INTEGER ivap ! indice de traceurs pour vapeur d'eau
145 PARAMETER (ivap = 1)
146 INTEGER iliq ! indice de traceurs pour eau liquide
147 PARAMETER (iliq = 2)
148
149 REAL, save:: t_ancien(klon, llm), q_ancien(klon, llm)
150 LOGICAL, save:: ancien_ok
151
152 REAL d_t_dyn(klon, llm) ! tendance dynamique pour "t" (K/s)
153 REAL d_q_dyn(klon, llm) ! tendance dynamique pour "q" (kg/kg/s)
154
155 real da(klon, llm), phi(klon, llm, llm), mp(klon, llm)
156
157 !IM Amip2 PV a theta constante
158
159 CHARACTER(LEN = 3) ctetaSTD(nbteta)
160 DATA ctetaSTD/'350', '380', '405'/
161 REAL rtetaSTD(nbteta)
162 DATA rtetaSTD/350., 380., 405./
163
164 !MI Amip2 PV a theta constante
165
166 REAL swdn0(klon, llm + 1), swdn(klon, llm + 1)
167 REAL swup0(klon, llm + 1), swup(klon, llm + 1)
168 SAVE swdn0, swdn, swup0, swup
169
170 REAL lwdn0(klon, llm + 1), lwdn(klon, llm + 1)
171 REAL lwup0(klon, llm + 1), lwup(klon, llm + 1)
172 SAVE lwdn0, lwdn, lwup0, lwup
173
174 !IM Amip2
175 ! variables a une pression donnee
176
177 integer nlevSTD
178 PARAMETER(nlevSTD = 17)
179 real rlevSTD(nlevSTD)
180 DATA rlevSTD/100000., 92500., 85000., 70000., &
181 60000., 50000., 40000., 30000., 25000., 20000., &
182 15000., 10000., 7000., 5000., 3000., 2000., 1000./
183 CHARACTER(LEN = 4) clevSTD(nlevSTD)
184 DATA clevSTD/'1000', '925 ', '850 ', '700 ', '600 ', &
185 '500 ', '400 ', '300 ', '250 ', '200 ', '150 ', '100 ', &
186 '70 ', '50 ', '30 ', '20 ', '10 '/
187
188 ! prw: precipitable water
189 real prw(klon)
190
191 ! flwp, fiwp = Liquid Water Path & Ice Water Path (kg/m2)
192 ! flwc, fiwc = Liquid Water Content & Ice Water Content (kg/kg)
193 REAL flwp(klon), fiwp(klon)
194 REAL flwc(klon, llm), fiwc(klon, llm)
195
196 INTEGER kmax, lmax
197 PARAMETER(kmax = 8, lmax = 8)
198 INTEGER kmaxm1, lmaxm1
199 PARAMETER(kmaxm1 = kmax-1, lmaxm1 = lmax-1)
200
201 REAL zx_tau(kmaxm1), zx_pc(lmaxm1)
202 DATA zx_tau/0., 0.3, 1.3, 3.6, 9.4, 23., 60./
203 DATA zx_pc/50., 180., 310., 440., 560., 680., 800./
204
205 ! cldtopres pression au sommet des nuages
206 REAL cldtopres(lmaxm1)
207 DATA cldtopres/50., 180., 310., 440., 560., 680., 800./
208
209 ! taulev: numero du niveau de tau dans les sorties ISCCP
210 CHARACTER(LEN = 4) taulev(kmaxm1)
211
212 DATA taulev/'tau0', 'tau1', 'tau2', 'tau3', 'tau4', 'tau5', 'tau6'/
213 CHARACTER(LEN = 3) pclev(lmaxm1)
214 DATA pclev/'pc1', 'pc2', 'pc3', 'pc4', 'pc5', 'pc6', 'pc7'/
215
216 CHARACTER(LEN = 28) cnameisccp(lmaxm1, kmaxm1)
217 DATA cnameisccp/'pc< 50hPa, tau< 0.3', 'pc= 50-180hPa, tau< 0.3', &
218 'pc= 180-310hPa, tau< 0.3', 'pc= 310-440hPa, tau< 0.3', &
219 'pc= 440-560hPa, tau< 0.3', 'pc= 560-680hPa, tau< 0.3', &
220 'pc= 680-800hPa, tau< 0.3', 'pc< 50hPa, tau= 0.3-1.3', &
221 'pc= 50-180hPa, tau= 0.3-1.3', 'pc= 180-310hPa, tau= 0.3-1.3', &
222 'pc= 310-440hPa, tau= 0.3-1.3', 'pc= 440-560hPa, tau= 0.3-1.3', &
223 'pc= 560-680hPa, tau= 0.3-1.3', 'pc= 680-800hPa, tau= 0.3-1.3', &
224 'pc< 50hPa, tau= 1.3-3.6', 'pc= 50-180hPa, tau= 1.3-3.6', &
225 'pc= 180-310hPa, tau= 1.3-3.6', 'pc= 310-440hPa, tau= 1.3-3.6', &
226 'pc= 440-560hPa, tau= 1.3-3.6', 'pc= 560-680hPa, tau= 1.3-3.6', &
227 'pc= 680-800hPa, tau= 1.3-3.6', 'pc< 50hPa, tau= 3.6-9.4', &
228 'pc= 50-180hPa, tau= 3.6-9.4', 'pc= 180-310hPa, tau= 3.6-9.4', &
229 'pc= 310-440hPa, tau= 3.6-9.4', 'pc= 440-560hPa, tau= 3.6-9.4', &
230 'pc= 560-680hPa, tau= 3.6-9.4', 'pc= 680-800hPa, tau= 3.6-9.4', &
231 'pc< 50hPa, tau= 9.4-23', 'pc= 50-180hPa, tau= 9.4-23', &
232 'pc= 180-310hPa, tau= 9.4-23', 'pc= 310-440hPa, tau= 9.4-23', &
233 'pc= 440-560hPa, tau= 9.4-23', 'pc= 560-680hPa, tau= 9.4-23', &
234 'pc= 680-800hPa, tau= 9.4-23', 'pc< 50hPa, tau= 23-60', &
235 'pc= 50-180hPa, tau= 23-60', 'pc= 180-310hPa, tau= 23-60', &
236 'pc= 310-440hPa, tau= 23-60', 'pc= 440-560hPa, tau= 23-60', &
237 'pc= 560-680hPa, tau= 23-60', 'pc= 680-800hPa, tau= 23-60', &
238 'pc< 50hPa, tau> 60.', 'pc= 50-180hPa, tau> 60.', &
239 'pc= 180-310hPa, tau> 60.', 'pc= 310-440hPa, tau> 60.', &
240 'pc= 440-560hPa, tau> 60.', 'pc= 560-680hPa, tau> 60.', &
241 'pc= 680-800hPa, tau> 60.'/
242
243 !IM ISCCP simulator v3.4
244
245 integer nid_hf, nid_hf3d
246 save nid_hf, nid_hf3d
247
248 ! Variables propres a la physique
249
250 INTEGER, save:: radpas
251 ! (Radiative transfer computations are made every "radpas" call to
252 ! "physiq".)
253
254 REAL radsol(klon)
255 SAVE radsol ! bilan radiatif au sol calcule par code radiatif
256
257 INTEGER, SAVE:: itap ! number of calls to "physiq"
258
259 REAL, save:: ftsol(klon, nbsrf) ! skin temperature of surface fraction
260
261 REAL, save:: ftsoil(klon, nsoilmx, nbsrf)
262 ! soil temperature of surface fraction
263
264 REAL, save:: fevap(klon, nbsrf) ! evaporation
265 REAL fluxlat(klon, nbsrf)
266 SAVE fluxlat
267
268 REAL fqsurf(klon, nbsrf)
269 SAVE fqsurf ! humidite de l'air au contact de la surface
270
271 REAL, save:: qsol(klon) ! hauteur d'eau dans le sol
272
273 REAL fsnow(klon, nbsrf)
274 SAVE fsnow ! epaisseur neigeuse
275
276 REAL falbe(klon, nbsrf)
277 SAVE falbe ! albedo par type de surface
278 REAL falblw(klon, nbsrf)
279 SAVE falblw ! albedo par type de surface
280
281 ! Param\`etres de l'orographie \`a l'\'echelle sous-maille (OESM) :
282 REAL, save:: zmea(klon) ! orographie moyenne
283 REAL, save:: zstd(klon) ! deviation standard de l'OESM
284 REAL, save:: zsig(klon) ! pente de l'OESM
285 REAL, save:: zgam(klon) ! anisotropie de l'OESM
286 REAL, save:: zthe(klon) ! orientation de l'OESM
287 REAL, save:: zpic(klon) ! Maximum de l'OESM
288 REAL, save:: zval(klon) ! Minimum de l'OESM
289 REAL, save:: rugoro(klon) ! longueur de rugosite de l'OESM
290
291 REAL zulow(klon), zvlow(klon)
292
293 INTEGER igwd, idx(klon), itest(klon)
294
295 REAL agesno(klon, nbsrf)
296 SAVE agesno ! age de la neige
297
298 REAL run_off_lic_0(klon)
299 SAVE run_off_lic_0
300 !KE43
301 ! Variables liees a la convection de K. Emanuel (sb):
302
303 REAL bas, top ! cloud base and top levels
304 SAVE bas
305 SAVE top
306
307 REAL Ma(klon, llm) ! undilute upward mass flux
308 SAVE Ma
309 REAL qcondc(klon, llm) ! in-cld water content from convect
310 SAVE qcondc
311 REAL, save:: sig1(klon, llm), w01(klon, llm)
312 REAL, save:: wd(klon)
313
314 ! Variables locales pour la couche limite (al1):
315
316 ! Variables locales:
317
318 REAL cdragh(klon) ! drag coefficient pour T and Q
319 REAL cdragm(klon) ! drag coefficient pour vent
320
321 ! Pour phytrac :
322 REAL ycoefh(klon, llm) ! coef d'echange pour phytrac
323 REAL yu1(klon) ! vents dans la premiere couche U
324 REAL yv1(klon) ! vents dans la premiere couche V
325 REAL ffonte(klon, nbsrf) !Flux thermique utilise pour fondre la neige
326 REAL fqcalving(klon, nbsrf) !Flux d'eau "perdue" par la surface
327 ! !et necessaire pour limiter la
328 ! !hauteur de neige, en kg/m2/s
329 REAL zxffonte(klon), zxfqcalving(klon)
330
331 REAL pfrac_impa(klon, llm)! Produits des coefs lessivage impaction
332 save pfrac_impa
333 REAL pfrac_nucl(klon, llm)! Produits des coefs lessivage nucleation
334 save pfrac_nucl
335 REAL pfrac_1nucl(klon, llm)! Produits des coefs lessi nucl (alpha = 1)
336 save pfrac_1nucl
337 REAL frac_impa(klon, llm) ! fractions d'aerosols lessivees (impaction)
338 REAL frac_nucl(klon, llm) ! idem (nucleation)
339
340 REAL, save:: rain_fall(klon) ! pluie
341 REAL, save:: snow_fall(klon) ! neige
342
343 REAL rain_tiedtke(klon), snow_tiedtke(klon)
344
345 REAL evap(klon), devap(klon) ! evaporation and its derivative
346 REAL sens(klon), dsens(klon) ! chaleur sensible et sa derivee
347 REAL dlw(klon) ! derivee infra rouge
348 SAVE dlw
349 REAL bils(klon) ! bilan de chaleur au sol
350 REAL fder(klon) ! Derive de flux (sensible et latente)
351 save fder
352 REAL ve(klon) ! integr. verticale du transport meri. de l'energie
353 REAL vq(klon) ! integr. verticale du transport meri. de l'eau
354 REAL ue(klon) ! integr. verticale du transport zonal de l'energie
355 REAL uq(klon) ! integr. verticale du transport zonal de l'eau
356
357 REAL frugs(klon, nbsrf) ! longueur de rugosite
358 save frugs
359 REAL zxrugs(klon) ! longueur de rugosite
360
361 ! Conditions aux limites
362
363 INTEGER julien
364
365 INTEGER, SAVE:: lmt_pas ! number of time steps of "physics" per day
366 REAL, save:: pctsrf(klon, nbsrf) ! percentage of surface
367 REAL pctsrf_new(klon, nbsrf) ! pourcentage surfaces issus d'ORCHIDEE
368
369 REAL albsol(klon)
370 SAVE albsol ! albedo du sol total
371 REAL albsollw(klon)
372 SAVE albsollw ! albedo du sol total
373
374 REAL, SAVE:: wo(klon, llm) ! column density of ozone in a cell, in kDU
375
376 ! Declaration des procedures appelees
377
378 EXTERNAL alboc ! calculer l'albedo sur ocean
379 !KE43
380 EXTERNAL conema3 ! convect4.3
381 EXTERNAL nuage ! calculer les proprietes radiatives
382 EXTERNAL transp ! transport total de l'eau et de l'energie
383
384 ! Variables locales
385
386 real, save:: clwcon(klon, llm), rnebcon(klon, llm)
387 real, save:: clwcon0(klon, llm), rnebcon0(klon, llm)
388
389 REAL rhcl(klon, llm) ! humiditi relative ciel clair
390 REAL dialiq(klon, llm) ! eau liquide nuageuse
391 REAL diafra(klon, llm) ! fraction nuageuse
392 REAL cldliq(klon, llm) ! eau liquide nuageuse
393 REAL cldfra(klon, llm) ! fraction nuageuse
394 REAL cldtau(klon, llm) ! epaisseur optique
395 REAL cldemi(klon, llm) ! emissivite infrarouge
396
397 REAL fluxq(klon, llm, nbsrf) ! flux turbulent d'humidite
398 REAL fluxt(klon, llm, nbsrf) ! flux turbulent de chaleur
399 REAL fluxu(klon, llm, nbsrf) ! flux turbulent de vitesse u
400 REAL fluxv(klon, llm, nbsrf) ! flux turbulent de vitesse v
401
402 REAL zxfluxt(klon, llm)
403 REAL zxfluxq(klon, llm)
404 REAL zxfluxu(klon, llm)
405 REAL zxfluxv(klon, llm)
406
407 ! Le rayonnement n'est pas calcul\'e tous les pas, il faut donc que
408 ! les variables soient r\'emanentes.
409 REAL, save:: heat(klon, llm) ! chauffage solaire
410 REAL heat0(klon, llm) ! chauffage solaire ciel clair
411 REAL, save:: cool(klon, llm) ! refroidissement infrarouge
412 REAL cool0(klon, llm) ! refroidissement infrarouge ciel clair
413 REAL, save:: topsw(klon), toplw(klon), solsw(klon)
414 REAL, save:: sollw(klon) ! rayonnement infrarouge montant \`a la surface
415 real, save:: sollwdown(klon) ! downward LW flux at surface
416 REAL, save:: topsw0(klon), toplw0(klon), solsw0(klon), sollw0(klon)
417 REAL albpla(klon)
418 REAL fsollw(klon, nbsrf) ! bilan flux IR pour chaque sous surface
419 REAL fsolsw(klon, nbsrf) ! flux solaire absorb. pour chaque sous surface
420 SAVE albpla
421 SAVE heat0, cool0
422
423 INTEGER itaprad
424 SAVE itaprad
425
426 REAL conv_q(klon, llm) ! convergence de l'humidite (kg/kg/s)
427 REAL conv_t(klon, llm) ! convergence of temperature (K/s)
428
429 REAL cldl(klon), cldm(klon), cldh(klon) !nuages bas, moyen et haut
430 REAL cldt(klon), cldq(klon) !nuage total, eau liquide integree
431
432 REAL zxtsol(klon), zxqsurf(klon), zxsnow(klon), zxfluxlat(klon)
433
434 REAL dist, rmu0(klon), fract(klon)
435 REAL zdtime ! pas de temps du rayonnement (s)
436 real zlongi
437 REAL z_avant(klon), z_apres(klon), z_factor(klon)
438 REAL za, zb
439 REAL zx_t, zx_qs, zdelta, zcor
440 real zqsat(klon, llm)
441 INTEGER i, k, iq, nsrf
442 REAL, PARAMETER:: t_coup = 234.
443 REAL zphi(klon, llm)
444
445 !IM cf. AM Variables locales pour la CLA (hbtm2)
446
447 REAL, SAVE:: pblh(klon, nbsrf) ! Hauteur de couche limite
448 REAL, SAVE:: plcl(klon, nbsrf) ! Niveau de condensation de la CLA
449 REAL, SAVE:: capCL(klon, nbsrf) ! CAPE de couche limite
450 REAL, SAVE:: oliqCL(klon, nbsrf) ! eau_liqu integree de couche limite
451 REAL, SAVE:: cteiCL(klon, nbsrf) ! cloud top instab. crit. couche limite
452 REAL, SAVE:: pblt(klon, nbsrf) ! T a la Hauteur de couche limite
453 REAL, SAVE:: therm(klon, nbsrf)
454 REAL, SAVE:: trmb1(klon, nbsrf) ! deep_cape
455 REAL, SAVE:: trmb2(klon, nbsrf) ! inhibition
456 REAL, SAVE:: trmb3(klon, nbsrf) ! Point Omega
457 ! Grdeurs de sorties
458 REAL s_pblh(klon), s_lcl(klon), s_capCL(klon)
459 REAL s_oliqCL(klon), s_cteiCL(klon), s_pblt(klon)
460 REAL s_therm(klon), s_trmb1(klon), s_trmb2(klon)
461 REAL s_trmb3(klon)
462
463 ! Variables locales pour la convection de K. Emanuel :
464
465 REAL upwd(klon, llm) ! saturated updraft mass flux
466 REAL dnwd(klon, llm) ! saturated downdraft mass flux
467 REAL dnwd0(klon, llm) ! unsaturated downdraft mass flux
468 REAL tvp(klon, llm) ! virtual temp of lifted parcel
469 REAL cape(klon) ! CAPE
470 SAVE cape
471
472 REAL pbase(klon) ! cloud base pressure
473 SAVE pbase
474 REAL bbase(klon) ! cloud base buoyancy
475 SAVE bbase
476 REAL rflag(klon) ! flag fonctionnement de convect
477 INTEGER iflagctrl(klon) ! flag fonctionnement de convect
478 ! -- convect43:
479 REAL dtvpdt1(klon, llm), dtvpdq1(klon, llm)
480 REAL dplcldt(klon), dplcldr(klon)
481
482 ! Variables du changement
483
484 ! con: convection
485 ! lsc: large scale condensation
486 ! ajs: ajustement sec
487 ! eva: \'evaporation de l'eau liquide nuageuse
488 ! vdf: vertical diffusion in boundary layer
489 REAL d_t_con(klon, llm), d_q_con(klon, llm)
490 REAL d_u_con(klon, llm), d_v_con(klon, llm)
491 REAL d_t_lsc(klon, llm), d_q_lsc(klon, llm), d_ql_lsc(klon, llm)
492 REAL d_t_ajs(klon, llm), d_q_ajs(klon, llm)
493 REAL d_u_ajs(klon, llm), d_v_ajs(klon, llm)
494 REAL rneb(klon, llm)
495
496 REAL mfu(klon, llm), mfd(klon, llm)
497 REAL pen_u(klon, llm), pen_d(klon, llm)
498 REAL pde_u(klon, llm), pde_d(klon, llm)
499 INTEGER kcbot(klon), kctop(klon), kdtop(klon)
500 REAL pmflxr(klon, llm + 1), pmflxs(klon, llm + 1)
501 REAL prfl(klon, llm + 1), psfl(klon, llm + 1)
502
503 INTEGER, save:: ibas_con(klon), itop_con(klon)
504
505 REAL rain_con(klon), rain_lsc(klon)
506 REAL snow_con(klon), snow_lsc(klon)
507 REAL d_ts(klon, nbsrf)
508
509 REAL d_u_vdf(klon, llm), d_v_vdf(klon, llm)
510 REAL d_t_vdf(klon, llm), d_q_vdf(klon, llm)
511
512 REAL d_u_oro(klon, llm), d_v_oro(klon, llm)
513 REAL d_t_oro(klon, llm)
514 REAL d_u_lif(klon, llm), d_v_lif(klon, llm)
515 REAL d_t_lif(klon, llm)
516
517 REAL, save:: ratqs(klon, llm)
518 real ratqss(klon, llm), ratqsc(klon, llm)
519 real:: ratqsbas = 0.01, ratqshaut = 0.3
520
521 ! Parametres lies au nouveau schema de nuages (SB, PDF)
522 real:: fact_cldcon = 0.375
523 real:: facttemps = 1.e-4
524 logical:: ok_newmicro = .true.
525 real facteur
526
527 integer:: iflag_cldcon = 1
528 logical ptconv(klon, llm)
529
530 ! Variables locales pour effectuer les appels en s\'erie :
531
532 REAL t_seri(klon, llm), q_seri(klon, llm)
533 REAL ql_seri(klon, llm), qs_seri(klon, llm)
534 REAL u_seri(klon, llm), v_seri(klon, llm)
535
536 REAL tr_seri(klon, llm, nbtr)
537 REAL d_tr(klon, llm, nbtr)
538
539 REAL zx_rh(klon, llm)
540
541 REAL zustrdr(klon), zvstrdr(klon)
542 REAL zustrli(klon), zvstrli(klon)
543 REAL zustrph(klon), zvstrph(klon)
544 REAL aam, torsfc
545
546 REAL zx_tmp_fi2d(klon) ! variable temporaire grille physique
547
548 INTEGER, SAVE:: nid_day, nid_ins
549
550 REAL ve_lay(klon, llm) ! transport meri. de l'energie a chaque niveau vert.
551 REAL vq_lay(klon, llm) ! transport meri. de l'eau a chaque niveau vert.
552 REAL ue_lay(klon, llm) ! transport zonal de l'energie a chaque niveau vert.
553 REAL uq_lay(klon, llm) ! transport zonal de l'eau a chaque niveau vert.
554
555 REAL zsto
556
557 logical ok_sync
558 real date0
559
560 ! Variables li\'ees au bilan d'\'energie et d'enthalpie :
561 REAL ztsol(klon)
562 REAL d_h_vcol, d_qt, d_qw, d_ql, d_qs, d_ec
563 REAL, SAVE:: d_h_vcol_phy
564 REAL fs_bound, fq_bound
565 REAL zero_v(klon)
566 CHARACTER(LEN = 15) tit
567 INTEGER:: ip_ebil = 0 ! print level for energy conservation diagnostics
568 INTEGER:: if_ebil = 0 ! verbosity for diagnostics of energy conservation
569
570 REAL d_t_ec(klon, llm) ! tendance due \`a la conversion Ec -> E thermique
571 REAL ZRCPD
572
573 REAL t2m(klon, nbsrf), q2m(klon, nbsrf) ! temperature and humidity at 2 m
574 REAL u10m(klon, nbsrf), v10m(klon, nbsrf) ! vents a 10 m
575 REAL zt2m(klon), zq2m(klon) ! temp., hum. 2 m moyenne s/ 1 maille
576 REAL zu10m(klon), zv10m(klon) ! vents a 10 m moyennes s/1 maille
577
578 ! Aerosol effects:
579
580 REAL sulfate(klon, llm) ! SO4 aerosol concentration (micro g/m3)
581
582 REAL, save:: sulfate_pi(klon, llm)
583 ! SO4 aerosol concentration, in micro g/m3, pre-industrial value
584
585 REAL cldtaupi(klon, llm)
586 ! cloud optical thickness for pre-industrial (pi) aerosols
587
588 REAL re(klon, llm) ! Cloud droplet effective radius
589 REAL fl(klon, llm) ! denominator of re
590
591 ! Aerosol optical properties
592 REAL, save:: tau_ae(klon, llm, 2), piz_ae(klon, llm, 2)
593 REAL, save:: cg_ae(klon, llm, 2)
594
595 REAL topswad(klon), solswad(klon) ! aerosol direct effect
596 REAL topswai(klon), solswai(klon) ! aerosol indirect effect
597
598 REAL aerindex(klon) ! POLDER aerosol index
599
600 LOGICAL:: ok_ade = .false. ! apply aerosol direct effect
601 LOGICAL:: ok_aie = .false. ! apply aerosol indirect effect
602
603 REAL:: bl95_b0 = 2., bl95_b1 = 0.2
604 ! Parameters in equation (D) of Boucher and Lohmann (1995, Tellus
605 ! B). They link cloud droplet number concentration to aerosol mass
606 ! concentration.
607
608 SAVE u10m
609 SAVE v10m
610 SAVE t2m
611 SAVE q2m
612 SAVE ffonte
613 SAVE fqcalving
614 SAVE rain_con
615 SAVE snow_con
616 SAVE topswai
617 SAVE topswad
618 SAVE solswai
619 SAVE solswad
620 SAVE d_u_con
621 SAVE d_v_con
622
623 real zmasse(klon, llm)
624 ! (column-density of mass of air in a cell, in kg m-2)
625
626 real, parameter:: dobson_u = 2.1415e-05 ! Dobson unit, in kg m-2
627
628 namelist /physiq_nml/ ocean, ok_veget, ok_journe, ok_mensuel, ok_instan, &
629 fact_cldcon, facttemps, ok_newmicro, iflag_cldcon, ratqsbas, &
630 ratqshaut, if_ebil, ok_ade, ok_aie, bl95_b0, bl95_b1, iflag_thermals, &
631 nsplit_thermals
632
633 !----------------------------------------------------------------
634
635 IF (if_ebil >= 1) zero_v = 0.
636 ok_sync = .TRUE.
637 IF (nqmx < 2) CALL abort_gcm('physiq', &
638 'eaux vapeur et liquide sont indispensables', 1)
639
640 test_firstcal: IF (firstcal) THEN
641 ! initialiser
642 u10m = 0.
643 v10m = 0.
644 t2m = 0.
645 q2m = 0.
646 ffonte = 0.
647 fqcalving = 0.
648 piz_ae = 0.
649 tau_ae = 0.
650 cg_ae = 0.
651 rain_con(:) = 0.
652 snow_con(:) = 0.
653 topswai(:) = 0.
654 topswad(:) = 0.
655 solswai(:) = 0.
656 solswad(:) = 0.
657
658 d_u_con = 0.
659 d_v_con = 0.
660 rnebcon0 = 0.
661 clwcon0 = 0.
662 rnebcon = 0.
663 clwcon = 0.
664
665 pblh =0. ! Hauteur de couche limite
666 plcl =0. ! Niveau de condensation de la CLA
667 capCL =0. ! CAPE de couche limite
668 oliqCL =0. ! eau_liqu integree de couche limite
669 cteiCL =0. ! cloud top instab. crit. couche limite
670 pblt =0. ! T a la Hauteur de couche limite
671 therm =0.
672 trmb1 =0. ! deep_cape
673 trmb2 =0. ! inhibition
674 trmb3 =0. ! Point Omega
675
676 IF (if_ebil >= 1) d_h_vcol_phy = 0.
677
678 iflag_thermals = 0
679 nsplit_thermals = 1
680 print *, "Enter namelist 'physiq_nml'."
681 read(unit=*, nml=physiq_nml)
682 write(unit_nml, nml=physiq_nml)
683
684 call conf_phys
685
686 ! Initialiser les compteurs:
687
688 frugs = 0.
689 itap = 0
690 itaprad = 0
691 CALL phyetat0("startphy.nc", pctsrf, ftsol, ftsoil, ocean, tslab, &
692 seaice, fqsurf, qsol, fsnow, falbe, falblw, fevap, rain_fall, &
693 snow_fall, solsw, sollw, dlw, radsol, frugs, agesno, zmea, &
694 zstd, zsig, zgam, zthe, zpic, zval, t_ancien, q_ancien, &
695 ancien_ok, rnebcon, ratqs, clwcon, run_off_lic_0, sig1, w01)
696
697 ! ATTENTION : il faudra a terme relire q2 dans l'etat initial
698 q2 = 1e-8
699
700 radpas = NINT(86400. / dtphys / nbapp_rad)
701
702 ! on remet le calendrier a zero
703 IF (raz_date) itau_phy = 0
704
705 PRINT *, 'cycle_diurne = ', cycle_diurne
706 CALL printflag(radpas, ocean /= 'force', ok_oasis, ok_journe, &
707 ok_instan, ok_region)
708
709 IF (dtphys * REAL(radpas) > 21600. .AND. cycle_diurne) THEN
710 print *, "Au minimum 4 appels par jour si cycle diurne"
711 call abort_gcm('physiq', &
712 "Nombre d'appels au rayonnement insuffisant", 1)
713 ENDIF
714
715 ! Initialisation pour le sch\'ema de convection d'Emanuel :
716 IF (iflag_con >= 3) THEN
717 ibas_con = 1
718 itop_con = 1
719 ENDIF
720
721 IF (ok_orodr) THEN
722 rugoro = MAX(1e-5, zstd * zsig / 2)
723 CALL SUGWD(paprs, play)
724 else
725 rugoro = 0.
726 ENDIF
727
728 lmt_pas = NINT(86400. / dtphys) ! tous les jours
729 print *, 'Number of time steps of "physics" per day: ', lmt_pas
730
731 ecrit_ins = NINT(ecrit_ins/dtphys)
732 ecrit_hf = NINT(ecrit_hf/dtphys)
733 ecrit_mth = NINT(ecrit_mth/dtphys)
734 ecrit_tra = NINT(86400.*ecrit_tra/dtphys)
735 ecrit_reg = NINT(ecrit_reg/dtphys)
736
737 ! Initialiser le couplage si necessaire
738
739 npas = 0
740 nexca = 0
741
742 ! Initialisation des sorties
743
744 call ini_histins(dtphys, ok_instan, nid_ins)
745 CALL ymds2ju(annee_ref, 1, int(day_ref), 0., date0)
746 ! Positionner date0 pour initialisation de ORCHIDEE
747 print *, 'physiq date0: ', date0
748 ENDIF test_firstcal
749
750 ! Mettre a zero des variables de sortie (pour securite)
751
752 DO i = 1, klon
753 d_ps(i) = 0.
754 ENDDO
755 DO iq = 1, nqmx
756 DO k = 1, llm
757 DO i = 1, klon
758 d_qx(i, k, iq) = 0.
759 ENDDO
760 ENDDO
761 ENDDO
762 da = 0.
763 mp = 0.
764 phi = 0.
765
766 ! Ne pas affecter les valeurs entr\'ees de u, v, h, et q :
767
768 DO k = 1, llm
769 DO i = 1, klon
770 t_seri(i, k) = t(i, k)
771 u_seri(i, k) = u(i, k)
772 v_seri(i, k) = v(i, k)
773 q_seri(i, k) = qx(i, k, ivap)
774 ql_seri(i, k) = qx(i, k, iliq)
775 qs_seri(i, k) = 0.
776 ENDDO
777 ENDDO
778 IF (nqmx >= 3) THEN
779 tr_seri(:, :, :nqmx-2) = qx(:, :, 3:nqmx)
780 ELSE
781 tr_seri(:, :, 1) = 0.
782 ENDIF
783
784 DO i = 1, klon
785 ztsol(i) = 0.
786 ENDDO
787 DO nsrf = 1, nbsrf
788 DO i = 1, klon
789 ztsol(i) = ztsol(i) + ftsol(i, nsrf)*pctsrf(i, nsrf)
790 ENDDO
791 ENDDO
792
793 IF (if_ebil >= 1) THEN
794 tit = 'after dynamics'
795 CALL diagetpq(airephy, tit, ip_ebil, 1, 1, dtphys, t_seri, q_seri, &
796 ql_seri, qs_seri, u_seri, v_seri, paprs, d_h_vcol, d_qt, d_qw, &
797 d_ql, d_qs, d_ec)
798 ! Comme les tendances de la physique sont ajout\'es dans la
799 ! dynamique, la variation d'enthalpie par la dynamique devrait
800 ! \^etre \'egale \`a la variation de la physique au pas de temps
801 ! pr\'ec\'edent. Donc la somme de ces 2 variations devrait \^etre
802 ! nulle.
803 call diagphy(airephy, tit, ip_ebil, zero_v, zero_v, zero_v, zero_v, &
804 zero_v, zero_v, zero_v, zero_v, ztsol, d_h_vcol + d_h_vcol_phy, &
805 d_qt, 0., fs_bound, fq_bound)
806 END IF
807
808 ! Diagnostic de la tendance dynamique :
809 IF (ancien_ok) THEN
810 DO k = 1, llm
811 DO i = 1, klon
812 d_t_dyn(i, k) = (t_seri(i, k) - t_ancien(i, k)) / dtphys
813 d_q_dyn(i, k) = (q_seri(i, k) - q_ancien(i, k)) / dtphys
814 ENDDO
815 ENDDO
816 ELSE
817 DO k = 1, llm
818 DO i = 1, klon
819 d_t_dyn(i, k) = 0.
820 d_q_dyn(i, k) = 0.
821 ENDDO
822 ENDDO
823 ancien_ok = .TRUE.
824 ENDIF
825
826 ! Ajouter le geopotentiel du sol:
827 DO k = 1, llm
828 DO i = 1, klon
829 zphi(i, k) = pphi(i, k) + pphis(i)
830 ENDDO
831 ENDDO
832
833 ! Check temperatures:
834 CALL hgardfou(t_seri, ftsol)
835
836 ! Incrementer le compteur de la physique
837 itap = itap + 1
838 julien = MOD(NINT(rdayvrai), 360)
839 if (julien == 0) julien = 360
840
841 forall (k = 1: llm) zmasse(:, k) = (paprs(:, k)-paprs(:, k + 1)) / rg
842
843 ! Mettre en action les conditions aux limites (albedo, sst etc.).
844
845 ! Prescrire l'ozone et calculer l'albedo sur l'ocean.
846 wo = ozonecm(REAL(julien), paprs)
847
848 ! \'Evaporation de l'eau liquide nuageuse :
849 DO k = 1, llm
850 DO i = 1, klon
851 zb = MAX(0., ql_seri(i, k))
852 t_seri(i, k) = t_seri(i, k) &
853 - zb * RLVTT / RCPD / (1. + RVTMP2 * q_seri(i, k))
854 q_seri(i, k) = q_seri(i, k) + zb
855 ENDDO
856 ENDDO
857 ql_seri = 0.
858
859 IF (if_ebil >= 2) THEN
860 tit = 'after reevap'
861 CALL diagetpq(airephy, tit, ip_ebil, 2, 1, dtphys, t_seri, q_seri, &
862 ql_seri, qs_seri, u_seri, v_seri, paprs, d_h_vcol, d_qt, d_qw, &
863 d_ql, d_qs, d_ec)
864 call diagphy(airephy, tit, ip_ebil, zero_v, zero_v, zero_v, zero_v, &
865 zero_v, zero_v, zero_v, zero_v, ztsol, d_h_vcol, d_qt, d_ec, &
866 fs_bound, fq_bound)
867
868 END IF
869
870 ! Appeler la diffusion verticale (programme de couche limite)
871
872 DO i = 1, klon
873 zxrugs(i) = 0.
874 ENDDO
875 DO nsrf = 1, nbsrf
876 DO i = 1, klon
877 frugs(i, nsrf) = MAX(frugs(i, nsrf), 0.000015)
878 ENDDO
879 ENDDO
880 DO nsrf = 1, nbsrf
881 DO i = 1, klon
882 zxrugs(i) = zxrugs(i) + frugs(i, nsrf)*pctsrf(i, nsrf)
883 ENDDO
884 ENDDO
885
886 ! calculs necessaires au calcul de l'albedo dans l'interface
887
888 CALL orbite(REAL(julien), zlongi, dist)
889 IF (cycle_diurne) THEN
890 zdtime = dtphys * REAL(radpas)
891 CALL zenang(zlongi, time, zdtime, rmu0, fract)
892 ELSE
893 rmu0 = -999.999
894 ENDIF
895
896 ! Calcul de l'abedo moyen par maille
897 albsol(:) = 0.
898 albsollw(:) = 0.
899 DO nsrf = 1, nbsrf
900 DO i = 1, klon
901 albsol(i) = albsol(i) + falbe(i, nsrf) * pctsrf(i, nsrf)
902 albsollw(i) = albsollw(i) + falblw(i, nsrf) * pctsrf(i, nsrf)
903 ENDDO
904 ENDDO
905
906 ! R\'epartition sous maille des flux longwave et shortwave
907 ! R\'epartition du longwave par sous-surface lin\'earis\'ee
908
909 DO nsrf = 1, nbsrf
910 DO i = 1, klon
911 fsollw(i, nsrf) = sollw(i) &
912 + 4. * RSIGMA * ztsol(i)**3 * (ztsol(i) - ftsol(i, nsrf))
913 fsolsw(i, nsrf) = solsw(i) * (1. - falbe(i, nsrf)) / (1. - albsol(i))
914 ENDDO
915 ENDDO
916
917 fder = dlw
918
919 ! Couche limite:
920
921 CALL clmain(dtphys, itap, pctsrf, pctsrf_new, t_seri, q_seri, &
922 u_seri, v_seri, julien, rmu0, co2_ppm, ok_veget, ocean, &
923 ftsol, soil_model, cdmmax, cdhmax, ksta, ksta_ter, ok_kzmin, ftsoil, &
924 qsol, paprs, play, fsnow, fqsurf, fevap, falbe, falblw, fluxlat, &
925 rain_fall, snow_fall, fsolsw, fsollw, fder, rlon, rlat, &
926 frugs, firstcal, agesno, rugoro, d_t_vdf, &
927 d_q_vdf, d_u_vdf, d_v_vdf, d_ts, fluxt, fluxq, fluxu, fluxv, cdragh, &
928 cdragm, q2, dsens, devap, ycoefh, yu1, yv1, t2m, q2m, u10m, v10m, &
929 pblh, capCL, oliqCL, cteiCL, pblT, therm, trmb1, trmb2, trmb3, plcl, &
930 fqcalving, ffonte, run_off_lic_0, fluxo, fluxg, tslab, seaice)
931
932 ! Incr\'ementation des flux
933
934 zxfluxt = 0.
935 zxfluxq = 0.
936 zxfluxu = 0.
937 zxfluxv = 0.
938 DO nsrf = 1, nbsrf
939 DO k = 1, llm
940 DO i = 1, klon
941 zxfluxt(i, k) = zxfluxt(i, k) + fluxt(i, k, nsrf) * pctsrf(i, nsrf)
942 zxfluxq(i, k) = zxfluxq(i, k) + fluxq(i, k, nsrf) * pctsrf(i, nsrf)
943 zxfluxu(i, k) = zxfluxu(i, k) + fluxu(i, k, nsrf) * pctsrf(i, nsrf)
944 zxfluxv(i, k) = zxfluxv(i, k) + fluxv(i, k, nsrf) * pctsrf(i, nsrf)
945 END DO
946 END DO
947 END DO
948 DO i = 1, klon
949 sens(i) = - zxfluxt(i, 1) ! flux de chaleur sensible au sol
950 evap(i) = - zxfluxq(i, 1) ! flux d'\'evaporation au sol
951 fder(i) = dlw(i) + dsens(i) + devap(i)
952 ENDDO
953
954 DO k = 1, llm
955 DO i = 1, klon
956 t_seri(i, k) = t_seri(i, k) + d_t_vdf(i, k)
957 q_seri(i, k) = q_seri(i, k) + d_q_vdf(i, k)
958 u_seri(i, k) = u_seri(i, k) + d_u_vdf(i, k)
959 v_seri(i, k) = v_seri(i, k) + d_v_vdf(i, k)
960 ENDDO
961 ENDDO
962
963 IF (if_ebil >= 2) THEN
964 tit = 'after clmain'
965 CALL diagetpq(airephy, tit, ip_ebil, 2, 2, dtphys, t_seri, q_seri, &
966 ql_seri, qs_seri, u_seri, v_seri, paprs, d_h_vcol, d_qt, d_qw, &
967 d_ql, d_qs, d_ec)
968 call diagphy(airephy, tit, ip_ebil, zero_v, zero_v, zero_v, zero_v, &
969 sens, evap, zero_v, zero_v, ztsol, d_h_vcol, d_qt, d_ec, &
970 fs_bound, fq_bound)
971 END IF
972
973 ! Update surface temperature:
974
975 DO i = 1, klon
976 zxtsol(i) = 0.
977 zxfluxlat(i) = 0.
978
979 zt2m(i) = 0.
980 zq2m(i) = 0.
981 zu10m(i) = 0.
982 zv10m(i) = 0.
983 zxffonte(i) = 0.
984 zxfqcalving(i) = 0.
985
986 s_pblh(i) = 0.
987 s_lcl(i) = 0.
988 s_capCL(i) = 0.
989 s_oliqCL(i) = 0.
990 s_cteiCL(i) = 0.
991 s_pblT(i) = 0.
992 s_therm(i) = 0.
993 s_trmb1(i) = 0.
994 s_trmb2(i) = 0.
995 s_trmb3(i) = 0.
996
997 IF (abs(pctsrf(i, is_ter) + pctsrf(i, is_lic) + pctsrf(i, is_oce) &
998 + pctsrf(i, is_sic) - 1.) > EPSFRA) print *, &
999 'physiq : probl\`eme sous surface au point ', i, pctsrf(i, 1 : nbsrf)
1000 ENDDO
1001 DO nsrf = 1, nbsrf
1002 DO i = 1, klon
1003 ftsol(i, nsrf) = ftsol(i, nsrf) + d_ts(i, nsrf)
1004 zxtsol(i) = zxtsol(i) + ftsol(i, nsrf)*pctsrf(i, nsrf)
1005 zxfluxlat(i) = zxfluxlat(i) + fluxlat(i, nsrf)*pctsrf(i, nsrf)
1006
1007 zt2m(i) = zt2m(i) + t2m(i, nsrf)*pctsrf(i, nsrf)
1008 zq2m(i) = zq2m(i) + q2m(i, nsrf)*pctsrf(i, nsrf)
1009 zu10m(i) = zu10m(i) + u10m(i, nsrf)*pctsrf(i, nsrf)
1010 zv10m(i) = zv10m(i) + v10m(i, nsrf)*pctsrf(i, nsrf)
1011 zxffonte(i) = zxffonte(i) + ffonte(i, nsrf)*pctsrf(i, nsrf)
1012 zxfqcalving(i) = zxfqcalving(i) + &
1013 fqcalving(i, nsrf)*pctsrf(i, nsrf)
1014 s_pblh(i) = s_pblh(i) + pblh(i, nsrf)*pctsrf(i, nsrf)
1015 s_lcl(i) = s_lcl(i) + plcl(i, nsrf)*pctsrf(i, nsrf)
1016 s_capCL(i) = s_capCL(i) + capCL(i, nsrf) *pctsrf(i, nsrf)
1017 s_oliqCL(i) = s_oliqCL(i) + oliqCL(i, nsrf) *pctsrf(i, nsrf)
1018 s_cteiCL(i) = s_cteiCL(i) + cteiCL(i, nsrf) *pctsrf(i, nsrf)
1019 s_pblT(i) = s_pblT(i) + pblT(i, nsrf) *pctsrf(i, nsrf)
1020 s_therm(i) = s_therm(i) + therm(i, nsrf) *pctsrf(i, nsrf)
1021 s_trmb1(i) = s_trmb1(i) + trmb1(i, nsrf) *pctsrf(i, nsrf)
1022 s_trmb2(i) = s_trmb2(i) + trmb2(i, nsrf) *pctsrf(i, nsrf)
1023 s_trmb3(i) = s_trmb3(i) + trmb3(i, nsrf) *pctsrf(i, nsrf)
1024 ENDDO
1025 ENDDO
1026
1027 ! Si une sous-fraction n'existe pas, elle prend la temp. moyenne
1028
1029 DO nsrf = 1, nbsrf
1030 DO i = 1, klon
1031 IF (pctsrf(i, nsrf) < epsfra) ftsol(i, nsrf) = zxtsol(i)
1032
1033 IF (pctsrf(i, nsrf) < epsfra) t2m(i, nsrf) = zt2m(i)
1034 IF (pctsrf(i, nsrf) < epsfra) q2m(i, nsrf) = zq2m(i)
1035 IF (pctsrf(i, nsrf) < epsfra) u10m(i, nsrf) = zu10m(i)
1036 IF (pctsrf(i, nsrf) < epsfra) v10m(i, nsrf) = zv10m(i)
1037 IF (pctsrf(i, nsrf) < epsfra) ffonte(i, nsrf) = zxffonte(i)
1038 IF (pctsrf(i, nsrf) < epsfra) &
1039 fqcalving(i, nsrf) = zxfqcalving(i)
1040 IF (pctsrf(i, nsrf) < epsfra) pblh(i, nsrf) = s_pblh(i)
1041 IF (pctsrf(i, nsrf) < epsfra) plcl(i, nsrf) = s_lcl(i)
1042 IF (pctsrf(i, nsrf) < epsfra) capCL(i, nsrf) = s_capCL(i)
1043 IF (pctsrf(i, nsrf) < epsfra) oliqCL(i, nsrf) = s_oliqCL(i)
1044 IF (pctsrf(i, nsrf) < epsfra) cteiCL(i, nsrf) = s_cteiCL(i)
1045 IF (pctsrf(i, nsrf) < epsfra) pblT(i, nsrf) = s_pblT(i)
1046 IF (pctsrf(i, nsrf) < epsfra) therm(i, nsrf) = s_therm(i)
1047 IF (pctsrf(i, nsrf) < epsfra) trmb1(i, nsrf) = s_trmb1(i)
1048 IF (pctsrf(i, nsrf) < epsfra) trmb2(i, nsrf) = s_trmb2(i)
1049 IF (pctsrf(i, nsrf) < epsfra) trmb3(i, nsrf) = s_trmb3(i)
1050 ENDDO
1051 ENDDO
1052
1053 ! Calculer la derive du flux infrarouge
1054
1055 DO i = 1, klon
1056 dlw(i) = - 4. * RSIGMA * zxtsol(i)**3
1057 ENDDO
1058
1059 ! Appeler la convection (au choix)
1060
1061 DO k = 1, llm
1062 DO i = 1, klon
1063 conv_q(i, k) = d_q_dyn(i, k) + d_q_vdf(i, k)/dtphys
1064 conv_t(i, k) = d_t_dyn(i, k) + d_t_vdf(i, k)/dtphys
1065 ENDDO
1066 ENDDO
1067
1068 IF (check) THEN
1069 za = qcheck(klon, llm, paprs, q_seri, ql_seri, airephy)
1070 print *, "avantcon = ", za
1071 ENDIF
1072
1073 if (iflag_con == 2) then
1074 z_avant = sum((q_seri + ql_seri) * zmasse, dim=2)
1075 CALL conflx(dtphys, paprs, play, t_seri(:, llm:1:-1), &
1076 q_seri(:, llm:1:-1), conv_t, conv_q, zxfluxq(:, 1), omega, &
1077 d_t_con, d_q_con, rain_con, snow_con, mfu(:, llm:1:-1), &
1078 mfd(:, llm:1:-1), pen_u, pde_u, pen_d, pde_d, kcbot, kctop, &
1079 kdtop, pmflxr, pmflxs)
1080 WHERE (rain_con < 0.) rain_con = 0.
1081 WHERE (snow_con < 0.) snow_con = 0.
1082 ibas_con = llm + 1 - kcbot
1083 itop_con = llm + 1 - kctop
1084 else
1085 ! iflag_con >= 3
1086
1087 CALL concvl(dtphys, paprs, play, t_seri, q_seri, u_seri, &
1088 v_seri, tr_seri, sig1, w01, d_t_con, d_q_con, &
1089 d_u_con, d_v_con, d_tr, rain_con, snow_con, ibas_con, &
1090 itop_con, upwd, dnwd, dnwd0, Ma, cape, tvp, iflagctrl, &
1091 pbase, bbase, dtvpdt1, dtvpdq1, dplcldt, dplcldr, qcondc, &
1092 wd, pmflxr, pmflxs, da, phi, mp, ntra=1)
1093 ! (number of tracers for the convection scheme of Kerry Emanuel:
1094 ! la partie traceurs est faite dans phytrac
1095 ! on met ntra = 1 pour limiter les appels mais on peut
1096 ! supprimer les calculs / ftra.)
1097
1098 clwcon0 = qcondc
1099 mfu = upwd + dnwd
1100 IF (.NOT. ok_gust) wd = 0.
1101
1102 ! Calcul des propri\'et\'es des nuages convectifs
1103
1104 DO k = 1, llm
1105 DO i = 1, klon
1106 zx_t = t_seri(i, k)
1107 IF (thermcep) THEN
1108 zdelta = MAX(0., SIGN(1., rtt-zx_t))
1109 zx_qs = r2es * FOEEW(zx_t, zdelta) / play(i, k)
1110 zx_qs = MIN(0.5, zx_qs)
1111 zcor = 1./(1.-retv*zx_qs)
1112 zx_qs = zx_qs*zcor
1113 ELSE
1114 IF (zx_t < t_coup) THEN
1115 zx_qs = qsats(zx_t)/play(i, k)
1116 ELSE
1117 zx_qs = qsatl(zx_t)/play(i, k)
1118 ENDIF
1119 ENDIF
1120 zqsat(i, k) = zx_qs
1121 ENDDO
1122 ENDDO
1123
1124 ! calcul des proprietes des nuages convectifs
1125 clwcon0 = fact_cldcon * clwcon0
1126 call clouds_gno(klon, llm, q_seri, zqsat, clwcon0, ptconv, ratqsc, &
1127 rnebcon0)
1128
1129 mfd = 0.
1130 pen_u = 0.
1131 pen_d = 0.
1132 pde_d = 0.
1133 pde_u = 0.
1134 END if
1135
1136 DO k = 1, llm
1137 DO i = 1, klon
1138 t_seri(i, k) = t_seri(i, k) + d_t_con(i, k)
1139 q_seri(i, k) = q_seri(i, k) + d_q_con(i, k)
1140 u_seri(i, k) = u_seri(i, k) + d_u_con(i, k)
1141 v_seri(i, k) = v_seri(i, k) + d_v_con(i, k)
1142 ENDDO
1143 ENDDO
1144
1145 IF (if_ebil >= 2) THEN
1146 tit = 'after convect'
1147 CALL diagetpq(airephy, tit, ip_ebil, 2, 2, dtphys, t_seri, q_seri, &
1148 ql_seri, qs_seri, u_seri, v_seri, paprs, d_h_vcol, d_qt, d_qw, &
1149 d_ql, d_qs, d_ec)
1150 call diagphy(airephy, tit, ip_ebil, zero_v, zero_v, zero_v, zero_v, &
1151 zero_v, zero_v, rain_con, snow_con, ztsol, d_h_vcol, d_qt, d_ec, &
1152 fs_bound, fq_bound)
1153 END IF
1154
1155 IF (check) THEN
1156 za = qcheck(klon, llm, paprs, q_seri, ql_seri, airephy)
1157 print *, "aprescon = ", za
1158 zx_t = 0.
1159 za = 0.
1160 DO i = 1, klon
1161 za = za + airephy(i)/REAL(klon)
1162 zx_t = zx_t + (rain_con(i)+ &
1163 snow_con(i))*airephy(i)/REAL(klon)
1164 ENDDO
1165 zx_t = zx_t/za*dtphys
1166 print *, "Precip = ", zx_t
1167 ENDIF
1168
1169 IF (iflag_con == 2) THEN
1170 z_apres = sum((q_seri + ql_seri) * zmasse, dim=2)
1171 z_factor = (z_avant - (rain_con + snow_con) * dtphys) / z_apres
1172 DO k = 1, llm
1173 DO i = 1, klon
1174 IF (z_factor(i) > 1. + 1E-8 .OR. z_factor(i) < 1. - 1E-8) THEN
1175 q_seri(i, k) = q_seri(i, k) * z_factor(i)
1176 ENDIF
1177 ENDDO
1178 ENDDO
1179 ENDIF
1180
1181 ! Convection s\`eche (thermiques ou ajustement)
1182
1183 d_t_ajs = 0.
1184 d_u_ajs = 0.
1185 d_v_ajs = 0.
1186 d_q_ajs = 0.
1187 fm_therm = 0.
1188 entr_therm = 0.
1189
1190 if (iflag_thermals == 0) then
1191 ! Ajustement sec
1192 CALL ajsec(paprs, play, t_seri, q_seri, d_t_ajs, d_q_ajs)
1193 t_seri = t_seri + d_t_ajs
1194 q_seri = q_seri + d_q_ajs
1195 else
1196 ! Thermiques
1197 call calltherm(dtphys, play, paprs, pphi, u_seri, v_seri, t_seri, &
1198 q_seri, d_u_ajs, d_v_ajs, d_t_ajs, d_q_ajs, fm_therm, entr_therm)
1199 endif
1200
1201 IF (if_ebil >= 2) THEN
1202 tit = 'after dry_adjust'
1203 CALL diagetpq(airephy, tit, ip_ebil, 2, 2, dtphys, t_seri, q_seri, &
1204 ql_seri, qs_seri, u_seri, v_seri, paprs, d_h_vcol, d_qt, d_qw, &
1205 d_ql, d_qs, d_ec)
1206 END IF
1207
1208 ! Caclul des ratqs
1209
1210 ! ratqs convectifs \`a l'ancienne en fonction de (q(z = 0) - q) / q
1211 ! on \'ecrase le tableau ratqsc calcul\'e par clouds_gno
1212 if (iflag_cldcon == 1) then
1213 do k = 1, llm
1214 do i = 1, klon
1215 if(ptconv(i, k)) then
1216 ratqsc(i, k) = ratqsbas + fact_cldcon &
1217 * (q_seri(i, 1) - q_seri(i, k)) / q_seri(i, k)
1218 else
1219 ratqsc(i, k) = 0.
1220 endif
1221 enddo
1222 enddo
1223 endif
1224
1225 ! ratqs stables
1226 do k = 1, llm
1227 do i = 1, klon
1228 ratqss(i, k) = ratqsbas + (ratqshaut - ratqsbas) &
1229 * min((paprs(i, 1) - play(i, k)) / (paprs(i, 1) - 3e4), 1.)
1230 enddo
1231 enddo
1232
1233 ! ratqs final
1234 if (iflag_cldcon == 1 .or. iflag_cldcon == 2) then
1235 ! les ratqs sont une conbinaison de ratqss et ratqsc
1236 ! ratqs final
1237 ! 1e4 (en gros 3 heures), en dur pour le moment, est le temps de
1238 ! relaxation des ratqs
1239 ratqs = max(ratqs * exp(- dtphys * facttemps), ratqss)
1240 ratqs = max(ratqs, ratqsc)
1241 else
1242 ! on ne prend que le ratqs stable pour fisrtilp
1243 ratqs = ratqss
1244 endif
1245
1246 CALL fisrtilp(dtphys, paprs, play, t_seri, q_seri, ptconv, ratqs, &
1247 d_t_lsc, d_q_lsc, d_ql_lsc, rneb, cldliq, rain_lsc, snow_lsc, &
1248 pfrac_impa, pfrac_nucl, pfrac_1nucl, frac_impa, frac_nucl, prfl, &
1249 psfl, rhcl)
1250
1251 WHERE (rain_lsc < 0) rain_lsc = 0.
1252 WHERE (snow_lsc < 0) snow_lsc = 0.
1253 DO k = 1, llm
1254 DO i = 1, klon
1255 t_seri(i, k) = t_seri(i, k) + d_t_lsc(i, k)
1256 q_seri(i, k) = q_seri(i, k) + d_q_lsc(i, k)
1257 ql_seri(i, k) = ql_seri(i, k) + d_ql_lsc(i, k)
1258 cldfra(i, k) = rneb(i, k)
1259 IF (.NOT.new_oliq) cldliq(i, k) = ql_seri(i, k)
1260 ENDDO
1261 ENDDO
1262 IF (check) THEN
1263 za = qcheck(klon, llm, paprs, q_seri, ql_seri, airephy)
1264 print *, "apresilp = ", za
1265 zx_t = 0.
1266 za = 0.
1267 DO i = 1, klon
1268 za = za + airephy(i)/REAL(klon)
1269 zx_t = zx_t + (rain_lsc(i) &
1270 + snow_lsc(i))*airephy(i)/REAL(klon)
1271 ENDDO
1272 zx_t = zx_t/za*dtphys
1273 print *, "Precip = ", zx_t
1274 ENDIF
1275
1276 IF (if_ebil >= 2) THEN
1277 tit = 'after fisrt'
1278 CALL diagetpq(airephy, tit, ip_ebil, 2, 2, dtphys, t_seri, q_seri, &
1279 ql_seri, qs_seri, u_seri, v_seri, paprs, d_h_vcol, d_qt, d_qw, &
1280 d_ql, d_qs, d_ec)
1281 call diagphy(airephy, tit, ip_ebil, zero_v, zero_v, zero_v, zero_v, &
1282 zero_v, zero_v, rain_lsc, snow_lsc, ztsol, d_h_vcol, d_qt, d_ec, &
1283 fs_bound, fq_bound)
1284 END IF
1285
1286 ! PRESCRIPTION DES NUAGES POUR LE RAYONNEMENT
1287
1288 ! 1. NUAGES CONVECTIFS
1289
1290 IF (iflag_cldcon <= -1) THEN
1291 ! seulement pour Tiedtke
1292 snow_tiedtke = 0.
1293 if (iflag_cldcon == -1) then
1294 rain_tiedtke = rain_con
1295 else
1296 rain_tiedtke = 0.
1297 do k = 1, llm
1298 do i = 1, klon
1299 if (d_q_con(i, k) < 0.) then
1300 rain_tiedtke(i) = rain_tiedtke(i)-d_q_con(i, k)/dtphys &
1301 *zmasse(i, k)
1302 endif
1303 enddo
1304 enddo
1305 endif
1306
1307 ! Nuages diagnostiques pour Tiedtke
1308 CALL diagcld1(paprs, play, rain_tiedtke, snow_tiedtke, ibas_con, &
1309 itop_con, diafra, dialiq)
1310 DO k = 1, llm
1311 DO i = 1, klon
1312 IF (diafra(i, k) > cldfra(i, k)) THEN
1313 cldliq(i, k) = dialiq(i, k)
1314 cldfra(i, k) = diafra(i, k)
1315 ENDIF
1316 ENDDO
1317 ENDDO
1318 ELSE IF (iflag_cldcon == 3) THEN
1319 ! On prend pour les nuages convectifs le maximum du calcul de
1320 ! la convection et du calcul du pas de temps pr\'ec\'edent diminu\'e
1321 ! d'un facteur facttemps.
1322 facteur = dtphys * facttemps
1323 do k = 1, llm
1324 do i = 1, klon
1325 rnebcon(i, k) = rnebcon(i, k) * facteur
1326 if (rnebcon0(i, k) * clwcon0(i, k) &
1327 > rnebcon(i, k) * clwcon(i, k)) then
1328 rnebcon(i, k) = rnebcon0(i, k)
1329 clwcon(i, k) = clwcon0(i, k)
1330 endif
1331 enddo
1332 enddo
1333
1334 ! On prend la somme des fractions nuageuses et des contenus en eau
1335 cldfra = min(max(cldfra, rnebcon), 1.)
1336 cldliq = cldliq + rnebcon*clwcon
1337 ENDIF
1338
1339 ! 2. Nuages stratiformes
1340
1341 IF (ok_stratus) THEN
1342 CALL diagcld2(paprs, play, t_seri, q_seri, diafra, dialiq)
1343 DO k = 1, llm
1344 DO i = 1, klon
1345 IF (diafra(i, k) > cldfra(i, k)) THEN
1346 cldliq(i, k) = dialiq(i, k)
1347 cldfra(i, k) = diafra(i, k)
1348 ENDIF
1349 ENDDO
1350 ENDDO
1351 ENDIF
1352
1353 ! Precipitation totale
1354 DO i = 1, klon
1355 rain_fall(i) = rain_con(i) + rain_lsc(i)
1356 snow_fall(i) = snow_con(i) + snow_lsc(i)
1357 ENDDO
1358
1359 IF (if_ebil >= 2) CALL diagetpq(airephy, "after diagcld", ip_ebil, 2, 2, &
1360 dtphys, t_seri, q_seri, ql_seri, qs_seri, u_seri, v_seri, paprs, &
1361 d_h_vcol, d_qt, d_qw, d_ql, d_qs, d_ec)
1362
1363 ! Humidit\'e relative pour diagnostic :
1364 DO k = 1, llm
1365 DO i = 1, klon
1366 zx_t = t_seri(i, k)
1367 IF (thermcep) THEN
1368 zdelta = MAX(0., SIGN(1., rtt-zx_t))
1369 zx_qs = r2es * FOEEW(zx_t, zdelta)/play(i, k)
1370 zx_qs = MIN(0.5, zx_qs)
1371 zcor = 1./(1.-retv*zx_qs)
1372 zx_qs = zx_qs*zcor
1373 ELSE
1374 IF (zx_t < t_coup) THEN
1375 zx_qs = qsats(zx_t)/play(i, k)
1376 ELSE
1377 zx_qs = qsatl(zx_t)/play(i, k)
1378 ENDIF
1379 ENDIF
1380 zx_rh(i, k) = q_seri(i, k)/zx_qs
1381 zqsat(i, k) = zx_qs
1382 ENDDO
1383 ENDDO
1384
1385 ! Introduce the aerosol direct and first indirect radiative forcings:
1386 IF (ok_ade .OR. ok_aie) THEN
1387 ! Get sulfate aerosol distribution :
1388 CALL readsulfate(rdayvrai, firstcal, sulfate)
1389 CALL readsulfate_preind(rdayvrai, firstcal, sulfate_pi)
1390
1391 CALL aeropt(play, paprs, t_seri, sulfate, rhcl, tau_ae, piz_ae, cg_ae, &
1392 aerindex)
1393 ELSE
1394 tau_ae = 0.
1395 piz_ae = 0.
1396 cg_ae = 0.
1397 ENDIF
1398
1399 ! Param\`etres optiques des nuages et quelques param\`etres pour diagnostics :
1400 if (ok_newmicro) then
1401 CALL newmicro(paprs, play, t_seri, cldliq, cldfra, cldtau, cldemi, &
1402 cldh, cldl, cldm, cldt, cldq, flwp, fiwp, flwc, fiwc, ok_aie, &
1403 sulfate, sulfate_pi, bl95_b0, bl95_b1, cldtaupi, re, fl)
1404 else
1405 CALL nuage(paprs, play, t_seri, cldliq, cldfra, cldtau, cldemi, cldh, &
1406 cldl, cldm, cldt, cldq, ok_aie, sulfate, sulfate_pi, bl95_b0, &
1407 bl95_b1, cldtaupi, re, fl)
1408 endif
1409
1410 ! Appeler le rayonnement mais calculer tout d'abord l'albedo du sol.
1411 IF (MOD(itaprad, radpas) == 0) THEN
1412 DO i = 1, klon
1413 albsol(i) = falbe(i, is_oce) * pctsrf(i, is_oce) &
1414 + falbe(i, is_lic) * pctsrf(i, is_lic) &
1415 + falbe(i, is_ter) * pctsrf(i, is_ter) &
1416 + falbe(i, is_sic) * pctsrf(i, is_sic)
1417 albsollw(i) = falblw(i, is_oce) * pctsrf(i, is_oce) &
1418 + falblw(i, is_lic) * pctsrf(i, is_lic) &
1419 + falblw(i, is_ter) * pctsrf(i, is_ter) &
1420 + falblw(i, is_sic) * pctsrf(i, is_sic)
1421 ENDDO
1422 ! Rayonnement (compatible Arpege-IFS) :
1423 CALL radlwsw(dist, rmu0, fract, paprs, play, zxtsol, albsol, &
1424 albsollw, t_seri, q_seri, wo, cldfra, cldemi, cldtau, heat, &
1425 heat0, cool, cool0, radsol, albpla, topsw, toplw, solsw, sollw, &
1426 sollwdown, topsw0, toplw0, solsw0, sollw0, lwdn0, lwdn, lwup0, &
1427 lwup, swdn0, swdn, swup0, swup, ok_ade, ok_aie, tau_ae, piz_ae, &
1428 cg_ae, topswad, solswad, cldtaupi, topswai, solswai)
1429 itaprad = 0
1430 ENDIF
1431 itaprad = itaprad + 1
1432
1433 ! Ajouter la tendance des rayonnements (tous les pas)
1434
1435 DO k = 1, llm
1436 DO i = 1, klon
1437 t_seri(i, k) = t_seri(i, k) + (heat(i, k)-cool(i, k)) * dtphys/86400.
1438 ENDDO
1439 ENDDO
1440
1441 IF (if_ebil >= 2) THEN
1442 tit = 'after rad'
1443 CALL diagetpq(airephy, tit, ip_ebil, 2, 2, dtphys, t_seri, q_seri, &
1444 ql_seri, qs_seri, u_seri, v_seri, paprs, d_h_vcol, d_qt, d_qw, &
1445 d_ql, d_qs, d_ec)
1446 call diagphy(airephy, tit, ip_ebil, topsw, toplw, solsw, sollw, &
1447 zero_v, zero_v, zero_v, zero_v, ztsol, d_h_vcol, d_qt, d_ec, &
1448 fs_bound, fq_bound)
1449 END IF
1450
1451 ! Calculer l'hydrologie de la surface
1452 DO i = 1, klon
1453 zxqsurf(i) = 0.
1454 zxsnow(i) = 0.
1455 ENDDO
1456 DO nsrf = 1, nbsrf
1457 DO i = 1, klon
1458 zxqsurf(i) = zxqsurf(i) + fqsurf(i, nsrf)*pctsrf(i, nsrf)
1459 zxsnow(i) = zxsnow(i) + fsnow(i, nsrf)*pctsrf(i, nsrf)
1460 ENDDO
1461 ENDDO
1462
1463 ! Calculer le bilan du sol et la d\'erive de temp\'erature (couplage)
1464
1465 DO i = 1, klon
1466 bils(i) = radsol(i) - sens(i) + zxfluxlat(i)
1467 ENDDO
1468
1469 ! Param\'etrisation de l'orographie \`a l'\'echelle sous-maille :
1470
1471 IF (ok_orodr) THEN
1472 ! selection des points pour lesquels le shema est actif:
1473 igwd = 0
1474 DO i = 1, klon
1475 itest(i) = 0
1476 IF (((zpic(i)-zmea(i)) > 100.).AND.(zstd(i) > 10.)) THEN
1477 itest(i) = 1
1478 igwd = igwd + 1
1479 idx(igwd) = i
1480 ENDIF
1481 ENDDO
1482
1483 CALL drag_noro(klon, llm, dtphys, paprs, play, zmea, zstd, zsig, zgam, &
1484 zthe, zpic, zval, igwd, idx, itest, t_seri, u_seri, v_seri, &
1485 zulow, zvlow, zustrdr, zvstrdr, d_t_oro, d_u_oro, d_v_oro)
1486
1487 ! ajout des tendances
1488 DO k = 1, llm
1489 DO i = 1, klon
1490 t_seri(i, k) = t_seri(i, k) + d_t_oro(i, k)
1491 u_seri(i, k) = u_seri(i, k) + d_u_oro(i, k)
1492 v_seri(i, k) = v_seri(i, k) + d_v_oro(i, k)
1493 ENDDO
1494 ENDDO
1495 ENDIF
1496
1497 IF (ok_orolf) THEN
1498 ! S\'election des points pour lesquels le sch\'ema est actif :
1499 igwd = 0
1500 DO i = 1, klon
1501 itest(i) = 0
1502 IF ((zpic(i) - zmea(i)) > 100.) THEN
1503 itest(i) = 1
1504 igwd = igwd + 1
1505 idx(igwd) = i
1506 ENDIF
1507 ENDDO
1508
1509 CALL lift_noro(klon, llm, dtphys, paprs, play, rlat, zmea, zstd, zpic, &
1510 itest, t_seri, u_seri, v_seri, zulow, zvlow, zustrli, zvstrli, &
1511 d_t_lif, d_u_lif, d_v_lif)
1512
1513 ! Ajout des tendances :
1514 DO k = 1, llm
1515 DO i = 1, klon
1516 t_seri(i, k) = t_seri(i, k) + d_t_lif(i, k)
1517 u_seri(i, k) = u_seri(i, k) + d_u_lif(i, k)
1518 v_seri(i, k) = v_seri(i, k) + d_v_lif(i, k)
1519 ENDDO
1520 ENDDO
1521 ENDIF
1522
1523 ! Stress n\'ecessaires : toute la physique
1524
1525 DO i = 1, klon
1526 zustrph(i) = 0.
1527 zvstrph(i) = 0.
1528 ENDDO
1529 DO k = 1, llm
1530 DO i = 1, klon
1531 zustrph(i) = zustrph(i) + (u_seri(i, k) - u(i, k)) / dtphys &
1532 * zmasse(i, k)
1533 zvstrph(i) = zvstrph(i) + (v_seri(i, k) - v(i, k)) / dtphys &
1534 * zmasse(i, k)
1535 ENDDO
1536 ENDDO
1537
1538 CALL aaam_bud(ra, rg, romega, rlat, rlon, pphis, zustrdr, zustrli, &
1539 zustrph, zvstrdr, zvstrli, zvstrph, paprs, u, v, aam, torsfc)
1540
1541 IF (if_ebil >= 2) CALL diagetpq(airephy, 'after orography', ip_ebil, 2, &
1542 2, dtphys, t_seri, q_seri, ql_seri, qs_seri, u_seri, v_seri, paprs, &
1543 d_h_vcol, d_qt, d_qw, d_ql, d_qs, d_ec)
1544
1545 ! Calcul des tendances traceurs
1546 call phytrac(rnpb, itap, lmt_pas, julien, time, firstcal, lafin, nqmx-2, &
1547 dtphys, u, t, paprs, play, mfu, mfd, pde_u, pen_d, ycoefh, fm_therm, &
1548 entr_therm, yu1, yv1, ftsol, pctsrf, frac_impa, frac_nucl, pphis, &
1549 albsol, rhcl, cldfra, rneb, diafra, cldliq, pmflxr, pmflxs, prfl, &
1550 psfl, da, phi, mp, upwd, dnwd, tr_seri, zmasse)
1551
1552 IF (offline) call phystokenc(dtphys, rlon, rlat, t, mfu, mfd, pen_u, &
1553 pde_u, pen_d, pde_d, fm_therm, entr_therm, ycoefh, yu1, yv1, ftsol, &
1554 pctsrf, frac_impa, frac_nucl, pphis, airephy, dtphys, itap)
1555
1556 ! Calculer le transport de l'eau et de l'energie (diagnostique)
1557 CALL transp(paprs, zxtsol, t_seri, q_seri, u_seri, v_seri, zphi, ve, vq, &
1558 ue, uq)
1559
1560 ! diag. bilKP
1561
1562 CALL transp_lay(paprs, zxtsol, t_seri, q_seri, u_seri, v_seri, zphi, &
1563 ve_lay, vq_lay, ue_lay, uq_lay)
1564
1565 ! Accumuler les variables a stocker dans les fichiers histoire:
1566
1567 ! conversion Ec -> E thermique
1568 DO k = 1, llm
1569 DO i = 1, klon
1570 ZRCPD = RCPD * (1. + RVTMP2 * q_seri(i, k))
1571 d_t_ec(i, k) = 0.5 / ZRCPD &
1572 * (u(i, k)**2 + v(i, k)**2 - u_seri(i, k)**2 - v_seri(i, k)**2)
1573 t_seri(i, k) = t_seri(i, k) + d_t_ec(i, k)
1574 d_t_ec(i, k) = d_t_ec(i, k) / dtphys
1575 END DO
1576 END DO
1577
1578 IF (if_ebil >= 1) THEN
1579 tit = 'after physic'
1580 CALL diagetpq(airephy, tit, ip_ebil, 1, 1, dtphys, t_seri, q_seri, &
1581 ql_seri, qs_seri, u_seri, v_seri, paprs, d_h_vcol, d_qt, d_qw, &
1582 d_ql, d_qs, d_ec)
1583 ! Comme les tendances de la physique sont ajoute dans la dynamique,
1584 ! on devrait avoir que la variation d'entalpie par la dynamique
1585 ! est egale a la variation de la physique au pas de temps precedent.
1586 ! Donc la somme de ces 2 variations devrait etre nulle.
1587 call diagphy(airephy, tit, ip_ebil, topsw, toplw, solsw, sollw, sens, &
1588 evap, rain_fall, snow_fall, ztsol, d_h_vcol, d_qt, d_ec, &
1589 fs_bound, fq_bound)
1590
1591 d_h_vcol_phy = d_h_vcol
1592
1593 END IF
1594
1595 ! SORTIES
1596
1597 ! prw = eau precipitable
1598 DO i = 1, klon
1599 prw(i) = 0.
1600 DO k = 1, llm
1601 prw(i) = prw(i) + q_seri(i, k)*zmasse(i, k)
1602 ENDDO
1603 ENDDO
1604
1605 ! Convertir les incrementations en tendances
1606
1607 DO k = 1, llm
1608 DO i = 1, klon
1609 d_u(i, k) = (u_seri(i, k) - u(i, k)) / dtphys
1610 d_v(i, k) = (v_seri(i, k) - v(i, k)) / dtphys
1611 d_t(i, k) = (t_seri(i, k) - t(i, k)) / dtphys
1612 d_qx(i, k, ivap) = (q_seri(i, k) - qx(i, k, ivap)) / dtphys
1613 d_qx(i, k, iliq) = (ql_seri(i, k) - qx(i, k, iliq)) / dtphys
1614 ENDDO
1615 ENDDO
1616
1617 IF (nqmx >= 3) THEN
1618 DO iq = 3, nqmx
1619 DO k = 1, llm
1620 DO i = 1, klon
1621 d_qx(i, k, iq) = (tr_seri(i, k, iq-2) - qx(i, k, iq)) / dtphys
1622 ENDDO
1623 ENDDO
1624 ENDDO
1625 ENDIF
1626
1627 ! Sauvegarder les valeurs de t et q a la fin de la physique:
1628 DO k = 1, llm
1629 DO i = 1, klon
1630 t_ancien(i, k) = t_seri(i, k)
1631 q_ancien(i, k) = q_seri(i, k)
1632 ENDDO
1633 ENDDO
1634
1635 ! Ecriture des sorties
1636 call write_histins
1637
1638 ! Si c'est la fin, il faut conserver l'etat de redemarrage
1639 IF (lafin) THEN
1640 itau_phy = itau_phy + itap
1641 CALL phyredem("restartphy.nc", rlat, rlon, pctsrf, ftsol, ftsoil, &
1642 tslab, seaice, fqsurf, qsol, fsnow, falbe, falblw, fevap, &
1643 rain_fall, snow_fall, solsw, sollw, dlw, radsol, frugs, &
1644 agesno, zmea, zstd, zsig, zgam, zthe, zpic, zval, t_ancien, &
1645 q_ancien, rnebcon, ratqs, clwcon, run_off_lic_0, sig1, w01)
1646 ENDIF
1647
1648 firstcal = .FALSE.
1649
1650 contains
1651
1652 subroutine write_histins
1653
1654 ! From phylmd/write_histins.h, version 1.2 2005/05/25 13:10:09
1655
1656 use dimens_m, only: iim, jjm
1657 USE histsync_m, ONLY: histsync
1658 USE histwrite_m, ONLY: histwrite
1659
1660 real zout
1661 integer itau_w ! pas de temps ecriture
1662 REAL zx_tmp_2d(iim, jjm + 1), zx_tmp_3d(iim, jjm + 1, llm)
1663
1664 !--------------------------------------------------
1665
1666 IF (ok_instan) THEN
1667 ! Champs 2D:
1668
1669 zsto = dtphys * ecrit_ins
1670 zout = dtphys * ecrit_ins
1671 itau_w = itau_phy + itap
1672
1673 i = NINT(zout/zsto)
1674 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, pphis, zx_tmp_2d)
1675 CALL histwrite(nid_ins, "phis", itau_w, zx_tmp_2d)
1676
1677 i = NINT(zout/zsto)
1678 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, airephy, zx_tmp_2d)
1679 CALL histwrite(nid_ins, "aire", itau_w, zx_tmp_2d)
1680
1681 DO i = 1, klon
1682 zx_tmp_fi2d(i) = paprs(i, 1)
1683 ENDDO
1684 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zx_tmp_fi2d, zx_tmp_2d)
1685 CALL histwrite(nid_ins, "psol", itau_w, zx_tmp_2d)
1686
1687 DO i = 1, klon
1688 zx_tmp_fi2d(i) = rain_fall(i) + snow_fall(i)
1689 ENDDO
1690 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zx_tmp_fi2d, zx_tmp_2d)
1691 CALL histwrite(nid_ins, "precip", itau_w, zx_tmp_2d)
1692
1693 DO i = 1, klon
1694 zx_tmp_fi2d(i) = rain_lsc(i) + snow_lsc(i)
1695 ENDDO
1696 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zx_tmp_fi2d, zx_tmp_2d)
1697 CALL histwrite(nid_ins, "plul", itau_w, zx_tmp_2d)
1698
1699 DO i = 1, klon
1700 zx_tmp_fi2d(i) = rain_con(i) + snow_con(i)
1701 ENDDO
1702 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zx_tmp_fi2d, zx_tmp_2d)
1703 CALL histwrite(nid_ins, "pluc", itau_w, zx_tmp_2d)
1704
1705 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zxtsol, zx_tmp_2d)
1706 CALL histwrite(nid_ins, "tsol", itau_w, zx_tmp_2d)
1707 !ccIM
1708 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zt2m, zx_tmp_2d)
1709 CALL histwrite(nid_ins, "t2m", itau_w, zx_tmp_2d)
1710
1711 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zq2m, zx_tmp_2d)
1712 CALL histwrite(nid_ins, "q2m", itau_w, zx_tmp_2d)
1713
1714 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zu10m, zx_tmp_2d)
1715 CALL histwrite(nid_ins, "u10m", itau_w, zx_tmp_2d)
1716
1717 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zv10m, zx_tmp_2d)
1718 CALL histwrite(nid_ins, "v10m", itau_w, zx_tmp_2d)
1719
1720 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, snow_fall, zx_tmp_2d)
1721 CALL histwrite(nid_ins, "snow", itau_w, zx_tmp_2d)
1722
1723 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, cdragm, zx_tmp_2d)
1724 CALL histwrite(nid_ins, "cdrm", itau_w, zx_tmp_2d)
1725
1726 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, cdragh, zx_tmp_2d)
1727 CALL histwrite(nid_ins, "cdrh", itau_w, zx_tmp_2d)
1728
1729 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, toplw, zx_tmp_2d)
1730 CALL histwrite(nid_ins, "topl", itau_w, zx_tmp_2d)
1731
1732 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, evap, zx_tmp_2d)
1733 CALL histwrite(nid_ins, "evap", itau_w, zx_tmp_2d)
1734
1735 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, solsw, zx_tmp_2d)
1736 CALL histwrite(nid_ins, "sols", itau_w, zx_tmp_2d)
1737
1738 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, sollw, zx_tmp_2d)
1739 CALL histwrite(nid_ins, "soll", itau_w, zx_tmp_2d)
1740
1741 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, sollwdown, zx_tmp_2d)
1742 CALL histwrite(nid_ins, "solldown", itau_w, zx_tmp_2d)
1743
1744 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, bils, zx_tmp_2d)
1745 CALL histwrite(nid_ins, "bils", itau_w, zx_tmp_2d)
1746
1747 zx_tmp_fi2d(1:klon) = -1*sens(1:klon)
1748 ! CALL gr_fi_ecrit(1, klon, iim, jjm + 1, sens, zx_tmp_2d)
1749 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zx_tmp_fi2d, zx_tmp_2d)
1750 CALL histwrite(nid_ins, "sens", itau_w, zx_tmp_2d)
1751
1752 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, fder, zx_tmp_2d)
1753 CALL histwrite(nid_ins, "fder", itau_w, zx_tmp_2d)
1754
1755 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, d_ts(1, is_oce), zx_tmp_2d)
1756 CALL histwrite(nid_ins, "dtsvdfo", itau_w, zx_tmp_2d)
1757
1758 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, d_ts(1, is_ter), zx_tmp_2d)
1759 CALL histwrite(nid_ins, "dtsvdft", itau_w, zx_tmp_2d)
1760
1761 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, d_ts(1, is_lic), zx_tmp_2d)
1762 CALL histwrite(nid_ins, "dtsvdfg", itau_w, zx_tmp_2d)
1763
1764 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, d_ts(1, is_sic), zx_tmp_2d)
1765 CALL histwrite(nid_ins, "dtsvdfi", itau_w, zx_tmp_2d)
1766
1767 DO nsrf = 1, nbsrf
1768 !XXX
1769 zx_tmp_fi2d(1 : klon) = pctsrf(1 : klon, nsrf)*100.
1770 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zx_tmp_fi2d, zx_tmp_2d)
1771 CALL histwrite(nid_ins, "pourc_"//clnsurf(nsrf), itau_w, &
1772 zx_tmp_2d)
1773
1774 zx_tmp_fi2d(1 : klon) = pctsrf(1 : klon, nsrf)
1775 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zx_tmp_fi2d, zx_tmp_2d)
1776 CALL histwrite(nid_ins, "fract_"//clnsurf(nsrf), itau_w, &
1777 zx_tmp_2d)
1778
1779 zx_tmp_fi2d(1 : klon) = fluxt(1 : klon, 1, nsrf)
1780 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zx_tmp_fi2d, zx_tmp_2d)
1781 CALL histwrite(nid_ins, "sens_"//clnsurf(nsrf), itau_w, &
1782 zx_tmp_2d)
1783
1784 zx_tmp_fi2d(1 : klon) = fluxlat(1 : klon, nsrf)
1785 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zx_tmp_fi2d, zx_tmp_2d)
1786 CALL histwrite(nid_ins, "lat_"//clnsurf(nsrf), itau_w, &
1787 zx_tmp_2d)
1788
1789 zx_tmp_fi2d(1 : klon) = ftsol(1 : klon, nsrf)
1790 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zx_tmp_fi2d, zx_tmp_2d)
1791 CALL histwrite(nid_ins, "tsol_"//clnsurf(nsrf), itau_w, &
1792 zx_tmp_2d)
1793
1794 zx_tmp_fi2d(1 : klon) = fluxu(1 : klon, 1, nsrf)
1795 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zx_tmp_fi2d, zx_tmp_2d)
1796 CALL histwrite(nid_ins, "taux_"//clnsurf(nsrf), itau_w, &
1797 zx_tmp_2d)
1798
1799 zx_tmp_fi2d(1 : klon) = fluxv(1 : klon, 1, nsrf)
1800 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zx_tmp_fi2d, zx_tmp_2d)
1801 CALL histwrite(nid_ins, "tauy_"//clnsurf(nsrf), itau_w, &
1802 zx_tmp_2d)
1803
1804 zx_tmp_fi2d(1 : klon) = frugs(1 : klon, nsrf)
1805 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zx_tmp_fi2d, zx_tmp_2d)
1806 CALL histwrite(nid_ins, "rugs_"//clnsurf(nsrf), itau_w, &
1807 zx_tmp_2d)
1808
1809 zx_tmp_fi2d(1 : klon) = falbe(1 : klon, nsrf)
1810 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zx_tmp_fi2d, zx_tmp_2d)
1811 CALL histwrite(nid_ins, "albe_"//clnsurf(nsrf), itau_w, &
1812 zx_tmp_2d)
1813
1814 END DO
1815 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, albsol, zx_tmp_2d)
1816 CALL histwrite(nid_ins, "albs", itau_w, zx_tmp_2d)
1817 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, albsollw, zx_tmp_2d)
1818 CALL histwrite(nid_ins, "albslw", itau_w, zx_tmp_2d)
1819
1820 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, zxrugs, zx_tmp_2d)
1821 CALL histwrite(nid_ins, "rugs", itau_w, zx_tmp_2d)
1822
1823 !HBTM2
1824
1825 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, s_pblh, zx_tmp_2d)
1826 CALL histwrite(nid_ins, "s_pblh", itau_w, zx_tmp_2d)
1827
1828 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, s_pblt, zx_tmp_2d)
1829 CALL histwrite(nid_ins, "s_pblt", itau_w, zx_tmp_2d)
1830
1831 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, s_lcl, zx_tmp_2d)
1832 CALL histwrite(nid_ins, "s_lcl", itau_w, zx_tmp_2d)
1833
1834 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, s_capCL, zx_tmp_2d)
1835 CALL histwrite(nid_ins, "s_capCL", itau_w, zx_tmp_2d)
1836
1837 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, s_oliqCL, zx_tmp_2d)
1838 CALL histwrite(nid_ins, "s_oliqCL", itau_w, zx_tmp_2d)
1839
1840 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, s_cteiCL, zx_tmp_2d)
1841 CALL histwrite(nid_ins, "s_cteiCL", itau_w, zx_tmp_2d)
1842
1843 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, s_therm, zx_tmp_2d)
1844 CALL histwrite(nid_ins, "s_therm", itau_w, zx_tmp_2d)
1845
1846 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, s_trmb1, zx_tmp_2d)
1847 CALL histwrite(nid_ins, "s_trmb1", itau_w, zx_tmp_2d)
1848
1849 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, s_trmb2, zx_tmp_2d)
1850 CALL histwrite(nid_ins, "s_trmb2", itau_w, zx_tmp_2d)
1851
1852 CALL gr_fi_ecrit(1, klon, iim, jjm + 1, s_trmb3, zx_tmp_2d)
1853 CALL histwrite(nid_ins, "s_trmb3", itau_w, zx_tmp_2d)
1854
1855 ! Champs 3D:
1856
1857 CALL gr_fi_ecrit(llm, klon, iim, jjm + 1, t_seri, zx_tmp_3d)
1858 CALL histwrite(nid_ins, "temp", itau_w, zx_tmp_3d)
1859
1860 CALL gr_fi_ecrit(llm, klon, iim, jjm + 1, u_seri, zx_tmp_3d)
1861 CALL histwrite(nid_ins, "vitu", itau_w, zx_tmp_3d)
1862
1863 CALL gr_fi_ecrit(llm, klon, iim, jjm + 1, v_seri, zx_tmp_3d)
1864 CALL histwrite(nid_ins, "vitv", itau_w, zx_tmp_3d)
1865
1866 CALL gr_fi_ecrit(llm, klon, iim, jjm + 1, zphi, zx_tmp_3d)
1867 CALL histwrite(nid_ins, "geop", itau_w, zx_tmp_3d)
1868
1869 CALL gr_fi_ecrit(llm, klon, iim, jjm + 1, play, zx_tmp_3d)
1870 CALL histwrite(nid_ins, "pres", itau_w, zx_tmp_3d)
1871
1872 CALL gr_fi_ecrit(llm, klon, iim, jjm + 1, d_t_vdf, zx_tmp_3d)
1873 CALL histwrite(nid_ins, "dtvdf", itau_w, zx_tmp_3d)
1874
1875 CALL gr_fi_ecrit(llm, klon, iim, jjm + 1, d_q_vdf, zx_tmp_3d)
1876 CALL histwrite(nid_ins, "dqvdf", itau_w, zx_tmp_3d)
1877
1878 if (ok_sync) then
1879 call histsync(nid_ins)
1880 endif
1881 ENDIF
1882
1883 end subroutine write_histins
1884
1885 END SUBROUTINE physiq
1886
1887 end module physiq_m

  ViewVC Help
Powered by ViewVC 1.1.21