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

Contents of /trunk/phylmd/physiq.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 69 - (show annotations)
Mon Feb 18 16:33:12 2013 UTC (11 years, 3 months ago) by guez
Original Path: trunk/libf/phylmd/physiq.f90
File size: 68855 byte(s)
Deleted files cvparam3.f90 and nuagecom.f90. Moved variables from
module cvparam3 to module cv3_param_m. Moved variables rad_chau1 and
rad_chau2 from module nuagecom to module conf_phys_m.

Read clesphys2_nml from conf_phys instead of gcm.

Removed argument iflag_con from several procedures. Access module
variable instead.

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

  ViewVC Help
Powered by ViewVC 1.1.21