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

Annotation of /trunk/Sources/phylmd/physiq.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 189 - (hide annotations)
Tue Mar 29 15:20:23 2016 UTC (8 years, 1 month ago) by guez
File size: 55611 byte(s)
There was a function gr_phy_write_3d in dyn3d and a function
gr_phy_write_2d in module grid_change. Moved them into a new module
gr_phy_write_m under a generic interface gr_phy_write. Replaced calls
to gr_fi_ecrit by calls to gr_phy_write.

Removed arguments len, nloc and nd of cv30_compress.

Removed arguments wd and wd1 of cv30_uncompress, wd of cv30_yield, wd
of concvl, wd1 of cv_driver. Was just filled with 0. Removed option
ok_gust in physiq, never used.

In cv30_unsat, cv30_yield and cv_driver, we only need to define b to
level nl - 1.

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

  ViewVC Help
Powered by ViewVC 1.1.21