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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 190 - (hide annotations)
Thu Apr 14 15:15:56 2016 UTC (8 years, 1 month ago) by guez
File size: 52648 byte(s)
Created module cv_thermo_m around procedure cv_thermo. Moved variables
from module cvthermo to module cv_thermo_m, where they are defined.

In ini_histins and initphysto, using part of rlon and rlat from
phyetat0_m is pretending that we do not know about the dynamical grid,
while the way we extract zx_lon(:, 1) and zx_lat(1, :) depends on
ordering inside rlon and rlat. So we might as well simplify and
clarify by using directly rlonv and rlatu.

Removed intermediary variables in write_histins and phystokenc.

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 190 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 190 LOGICAL, PARAMETER:: check = .FALSE.
112 guez 98 ! 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 190 SAVE lwdn0, lwdn, lwup0, lwup
147 guez 3
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 190 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 190 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 guez 190 REAL, SAVE:: trmb2(klon, nbsrf) ! inhibition
323 guez 49 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 97 INTEGER, SAVE:: nid_ins
405 guez 3
406     REAL ve_lay(klon, llm) ! transport meri. de l'energie a chaque niveau vert.
407     REAL vq_lay(klon, llm) ! transport meri. de l'eau a chaque niveau vert.
408     REAL ue_lay(klon, llm) ! transport zonal de l'energie a chaque niveau vert.
409     REAL uq_lay(klon, llm) ! transport zonal de l'eau a chaque niveau vert.
410    
411     real date0
412    
413 guez 90 ! Variables li\'ees au bilan d'\'energie et d'enthalpie :
414 guez 3 REAL ztsol(klon)
415 guez 98 REAL d_h_vcol, d_qt, d_ec
416 guez 51 REAL, SAVE:: d_h_vcol_phy
417 guez 47 REAL zero_v(klon)
418 guez 97 CHARACTER(LEN = 20) tit
419 guez 51 INTEGER:: ip_ebil = 0 ! print level for energy conservation diagnostics
420 guez 190 INTEGER:: if_ebil = 0 ! verbosity for diagnostics of energy conservation
421 guez 51
422 guez 90 REAL d_t_ec(klon, llm) ! tendance due \`a la conversion Ec -> E thermique
423 guez 3 REAL ZRCPD
424 guez 51
425 guez 49 REAL t2m(klon, nbsrf), q2m(klon, nbsrf) ! temperature and humidity at 2 m
426 guez 69 REAL u10m(klon, nbsrf), v10m(klon, nbsrf) ! vents a 10 m
427     REAL zt2m(klon), zq2m(klon) ! temp., hum. 2 m moyenne s/ 1 maille
428     REAL zu10m(klon), zv10m(klon) ! vents a 10 m moyennes s/1 maille
429 guez 3
430 guez 69 ! Aerosol effects:
431    
432     REAL sulfate(klon, llm) ! SO4 aerosol concentration (micro g/m3)
433    
434 guez 49 REAL, save:: sulfate_pi(klon, llm)
435 guez 175 ! SO4 aerosol concentration, in \mu g/m3, pre-industrial value
436 guez 3
437     REAL cldtaupi(klon, llm)
438 guez 69 ! cloud optical thickness for pre-industrial (pi) aerosols
439 guez 3
440 guez 47 REAL re(klon, llm) ! Cloud droplet effective radius
441     REAL fl(klon, llm) ! denominator of re
442 guez 3
443     ! Aerosol optical properties
444 guez 68 REAL, save:: tau_ae(klon, llm, 2), piz_ae(klon, llm, 2)
445     REAL, save:: cg_ae(klon, llm, 2)
446 guez 3
447 guez 68 REAL topswad(klon), solswad(klon) ! aerosol direct effect
448 guez 62 REAL topswai(klon), solswai(klon) ! aerosol indirect effect
449 guez 3
450 guez 47 REAL aerindex(klon) ! POLDER aerosol index
451 guez 3
452 guez 68 LOGICAL:: ok_ade = .false. ! apply aerosol direct effect
453     LOGICAL:: ok_aie = .false. ! apply aerosol indirect effect
454 guez 3
455 guez 68 REAL:: bl95_b0 = 2., bl95_b1 = 0.2
456 guez 69 ! Parameters in equation (D) of Boucher and Lohmann (1995, Tellus
457     ! B). They link cloud droplet number concentration to aerosol mass
458     ! concentration.
459 guez 68
460 guez 3 SAVE u10m
461     SAVE v10m
462     SAVE t2m
463     SAVE q2m
464     SAVE ffonte
465     SAVE fqcalving
466     SAVE rain_con
467     SAVE topswai
468     SAVE topswad
469     SAVE solswai
470     SAVE solswad
471     SAVE d_u_con
472     SAVE d_v_con
473    
474 guez 190 real zmasse(klon, llm)
475 guez 17 ! (column-density of mass of air in a cell, in kg m-2)
476    
477 guez 175 integer, save:: ncid_startphy, itau_phy
478 guez 17
479 guez 99 namelist /physiq_nml/ ok_journe, ok_mensuel, ok_instan, fact_cldcon, &
480     facttemps, ok_newmicro, iflag_cldcon, ratqsbas, ratqshaut, if_ebil, &
481     ok_ade, ok_aie, bl95_b0, bl95_b1, iflag_thermals, nsplit_thermals
482 guez 68
483 guez 3 !----------------------------------------------------------------
484    
485 guez 69 IF (if_ebil >= 1) zero_v = 0.
486     IF (nqmx < 2) CALL abort_gcm('physiq', &
487 guez 171 'eaux vapeur et liquide sont indispensables')
488 guez 3
489 guez 7 test_firstcal: IF (firstcal) THEN
490 guez 47 ! initialiser
491 guez 51 u10m = 0.
492     v10m = 0.
493     t2m = 0.
494     q2m = 0.
495     ffonte = 0.
496     fqcalving = 0.
497     piz_ae = 0.
498     tau_ae = 0.
499     cg_ae = 0.
500 guez 98 rain_con = 0.
501     snow_con = 0.
502     topswai = 0.
503     topswad = 0.
504     solswai = 0.
505     solswad = 0.
506 guez 3
507 guez 72 d_u_con = 0.
508     d_v_con = 0.
509     rnebcon0 = 0.
510     clwcon0 = 0.
511     rnebcon = 0.
512     clwcon = 0.
513 guez 3
514 guez 47 pblh =0. ! Hauteur de couche limite
515     plcl =0. ! Niveau de condensation de la CLA
516     capCL =0. ! CAPE de couche limite
517     oliqCL =0. ! eau_liqu integree de couche limite
518     cteiCL =0. ! cloud top instab. crit. couche limite
519     pblt =0. ! T a la Hauteur de couche limite
520     therm =0.
521     trmb1 =0. ! deep_cape
522 guez 190 trmb2 =0. ! inhibition
523 guez 47 trmb3 =0. ! Point Omega
524 guez 3
525 guez 51 IF (if_ebil >= 1) d_h_vcol_phy = 0.
526 guez 3
527 guez 68 iflag_thermals = 0
528     nsplit_thermals = 1
529     print *, "Enter namelist 'physiq_nml'."
530     read(unit=*, nml=physiq_nml)
531     write(unit_nml, nml=physiq_nml)
532    
533     call conf_phys
534 guez 3
535     ! Initialiser les compteurs:
536    
537     frugs = 0.
538 guez 175 CALL phyetat0(pctsrf, ftsol, ftsoil, fqsurf, qsol, &
539 guez 157 fsnow, falbe, fevap, rain_fall, snow_fall, solsw, sollw, dlw, &
540     radsol, frugs, agesno, zmea, zstd, zsig, zgam, zthe, zpic, zval, &
541     t_ancien, q_ancien, ancien_ok, rnebcon, ratqs, clwcon, &
542 guez 175 run_off_lic_0, sig1, w01, ncid_startphy, itau_phy)
543 guez 3
544 guez 47 ! ATTENTION : il faudra a terme relire q2 dans l'etat initial
545 guez 69 q2 = 1e-8
546 guez 3
547 guez 154 lmt_pas = day_step / iphysiq
548     print *, 'Number of time steps of "physics" per day: ', lmt_pas
549 guez 3
550 guez 154 radpas = lmt_pas / nbapp_rad
551    
552     ! On remet le calendrier a zero
553 guez 15 IF (raz_date) itau_phy = 0
554 guez 3
555 guez 99 CALL printflag(radpas, ok_journe, ok_instan, ok_region)
556 guez 3
557 guez 90 ! Initialisation pour le sch\'ema de convection d'Emanuel :
558 guez 182 IF (conv_emanuel) THEN
559 guez 69 ibas_con = 1
560     itop_con = 1
561 guez 3 ENDIF
562    
563     IF (ok_orodr) THEN
564 guez 13 rugoro = MAX(1e-5, zstd * zsig / 2)
565 guez 54 CALL SUGWD(paprs, play)
566 guez 13 else
567     rugoro = 0.
568 guez 3 ENDIF
569    
570 guez 47 ecrit_ins = NINT(ecrit_ins/dtphys)
571     ecrit_hf = NINT(ecrit_hf/dtphys)
572     ecrit_mth = NINT(ecrit_mth/dtphys)
573     ecrit_tra = NINT(86400.*ecrit_tra/dtphys)
574     ecrit_reg = NINT(ecrit_reg/dtphys)
575 guez 3
576 guez 47 ! Initialisation des sorties
577 guez 3
578 guez 175 call ini_histins(dtphys, ok_instan, nid_ins, itau_phy)
579 guez 129 CALL ymds2ju(annee_ref, 1, day_ref, 0., date0)
580 guez 69 ! Positionner date0 pour initialisation de ORCHIDEE
581     print *, 'physiq date0: ', date0
582 guez 175 CALL phyredem0(lmt_pas, itau_phy)
583 guez 7 ENDIF test_firstcal
584 guez 3
585 guez 91 ! We will modify variables *_seri and we will not touch variables
586 guez 98 ! u, v, t, qx:
587     t_seri = t
588     u_seri = u
589     v_seri = v
590     q_seri = qx(:, :, ivap)
591     ql_seri = qx(:, :, iliq)
592 guez 157 tr_seri = qx(:, :, 3:nqmx)
593 guez 3
594 guez 98 ztsol = sum(ftsol * pctsrf, dim = 2)
595 guez 3
596 guez 190 IF (if_ebil >= 1) THEN
597 guez 62 tit = 'after dynamics'
598     CALL diagetpq(airephy, tit, ip_ebil, 1, 1, dtphys, t_seri, q_seri, &
599 guez 98 ql_seri, u_seri, v_seri, paprs, d_h_vcol, d_qt, d_ec)
600 guez 90 ! Comme les tendances de la physique sont ajout\'es dans la
601 guez 190 ! dynamique, la variation d'enthalpie par la dynamique devrait
602     ! \^etre \'egale \`a la variation de la physique au pas de temps
603     ! pr\'ec\'edent. Donc la somme de ces 2 variations devrait \^etre
604     ! nulle.
605 guez 62 call diagphy(airephy, tit, ip_ebil, zero_v, zero_v, zero_v, zero_v, &
606 guez 51 zero_v, zero_v, zero_v, zero_v, ztsol, d_h_vcol + d_h_vcol_phy, &
607 guez 98 d_qt, 0.)
608 guez 3 END IF
609    
610 guez 51 ! Diagnostic de la tendance dynamique :
611 guez 3 IF (ancien_ok) THEN
612     DO k = 1, llm
613     DO i = 1, klon
614 guez 49 d_t_dyn(i, k) = (t_seri(i, k) - t_ancien(i, k)) / dtphys
615     d_q_dyn(i, k) = (q_seri(i, k) - q_ancien(i, k)) / dtphys
616 guez 3 ENDDO
617     ENDDO
618     ELSE
619     DO k = 1, llm
620     DO i = 1, klon
621 guez 72 d_t_dyn(i, k) = 0.
622     d_q_dyn(i, k) = 0.
623 guez 3 ENDDO
624     ENDDO
625     ancien_ok = .TRUE.
626     ENDIF
627    
628     ! Ajouter le geopotentiel du sol:
629     DO k = 1, llm
630     DO i = 1, klon
631     zphi(i, k) = pphi(i, k) + pphis(i)
632     ENDDO
633     ENDDO
634    
635 guez 49 ! Check temperatures:
636 guez 3 CALL hgardfou(t_seri, ftsol)
637    
638 guez 98 ! Incrémenter le compteur de la physique
639 guez 7 itap = itap + 1
640 guez 130 julien = MOD(dayvrai, 360)
641 guez 3 if (julien == 0) julien = 360
642    
643 guez 103 forall (k = 1: llm) zmasse(:, k) = (paprs(:, k) - paprs(:, k + 1)) / rg
644 guez 17
645 guez 98 ! Prescrire l'ozone :
646 guez 57 wo = ozonecm(REAL(julien), paprs)
647 guez 3
648 guez 90 ! \'Evaporation de l'eau liquide nuageuse :
649 guez 51 DO k = 1, llm
650 guez 3 DO i = 1, klon
651 guez 51 zb = MAX(0., ql_seri(i, k))
652     t_seri(i, k) = t_seri(i, k) &
653     - zb * RLVTT / RCPD / (1. + RVTMP2 * q_seri(i, k))
654 guez 3 q_seri(i, k) = q_seri(i, k) + zb
655     ENDDO
656     ENDDO
657 guez 51 ql_seri = 0.
658 guez 3
659 guez 190 IF (if_ebil >= 2) THEN
660 guez 62 tit = 'after reevap'
661     CALL diagetpq(airephy, tit, ip_ebil, 2, 1, dtphys, t_seri, q_seri, &
662 guez 98 ql_seri, u_seri, v_seri, paprs, d_h_vcol, d_qt, d_ec)
663 guez 62 call diagphy(airephy, tit, ip_ebil, zero_v, zero_v, zero_v, zero_v, &
664 guez 98 zero_v, zero_v, zero_v, zero_v, ztsol, d_h_vcol, d_qt, d_ec)
665 guez 3 END IF
666    
667 guez 98 frugs = MAX(frugs, 0.000015)
668     zxrugs = sum(frugs * pctsrf, dim = 2)
669 guez 3
670 guez 118 ! Calculs nécessaires au calcul de l'albedo dans l'interface avec
671     ! la surface.
672 guez 3
673 guez 118 CALL orbite(REAL(julien), longi, dist)
674 guez 3 IF (cycle_diurne) THEN
675 guez 125 CALL zenang(longi, time, dtphys * radpas, mu0, fract)
676 guez 3 ELSE
677 guez 174 mu0 = - 999.999
678 guez 3 ENDIF
679    
680 guez 47 ! Calcul de l'abedo moyen par maille
681 guez 98 albsol = sum(falbe * pctsrf, dim = 2)
682 guez 3
683 guez 90 ! R\'epartition sous maille des flux longwave et shortwave
684     ! R\'epartition du longwave par sous-surface lin\'earis\'ee
685 guez 3
686 guez 98 forall (nsrf = 1: nbsrf)
687     fsollw(:, nsrf) = sollw + 4. * RSIGMA * ztsol**3 &
688     * (ztsol - ftsol(:, nsrf))
689     fsolsw(:, nsrf) = solsw * (1. - falbe(:, nsrf)) / (1. - albsol)
690     END forall
691 guez 3
692     fder = dlw
693    
694 guez 38 ! Couche limite:
695 guez 3
696 guez 98 CALL clmain(dtphys, itap, pctsrf, pctsrf_new, t_seri, q_seri, u_seri, &
697 guez 178 v_seri, julien, mu0, ftsol, cdmmax, cdhmax, ksta, ksta_ter, &
698 guez 154 ok_kzmin, ftsoil, qsol, paprs, play, fsnow, fqsurf, fevap, falbe, &
699 guez 155 fluxlat, rain_fall, snow_fall, fsolsw, fsollw, fder, rlat, frugs, &
700     firstcal, agesno, rugoro, d_t_vdf, d_q_vdf, d_u_vdf, d_v_vdf, d_ts, &
701     fluxt, fluxq, fluxu, fluxv, cdragh, cdragm, q2, dsens, devap, &
702 guez 154 ycoefh, yu1, yv1, t2m, q2m, u10m, v10m, pblh, capCL, oliqCL, cteiCL, &
703     pblT, therm, trmb1, trmb2, trmb3, plcl, fqcalving, ffonte, &
704 guez 175 run_off_lic_0)
705 guez 3
706 guez 90 ! Incr\'ementation des flux
707 guez 40
708 guez 51 zxfluxt = 0.
709     zxfluxq = 0.
710     zxfluxu = 0.
711     zxfluxv = 0.
712 guez 3 DO nsrf = 1, nbsrf
713     DO k = 1, llm
714     DO i = 1, klon
715 guez 70 zxfluxt(i, k) = zxfluxt(i, k) + fluxt(i, k, nsrf) * pctsrf(i, nsrf)
716     zxfluxq(i, k) = zxfluxq(i, k) + fluxq(i, k, nsrf) * pctsrf(i, nsrf)
717     zxfluxu(i, k) = zxfluxu(i, k) + fluxu(i, k, nsrf) * pctsrf(i, nsrf)
718     zxfluxv(i, k) = zxfluxv(i, k) + fluxv(i, k, nsrf) * pctsrf(i, nsrf)
719 guez 3 END DO
720     END DO
721     END DO
722     DO i = 1, klon
723     sens(i) = - zxfluxt(i, 1) ! flux de chaleur sensible au sol
724 guez 90 evap(i) = - zxfluxq(i, 1) ! flux d'\'evaporation au sol
725 guez 3 fder(i) = dlw(i) + dsens(i) + devap(i)
726     ENDDO
727    
728     DO k = 1, llm
729     DO i = 1, klon
730     t_seri(i, k) = t_seri(i, k) + d_t_vdf(i, k)
731     q_seri(i, k) = q_seri(i, k) + d_q_vdf(i, k)
732     u_seri(i, k) = u_seri(i, k) + d_u_vdf(i, k)
733     v_seri(i, k) = v_seri(i, k) + d_v_vdf(i, k)
734     ENDDO
735     ENDDO
736    
737 guez 190 IF (if_ebil >= 2) THEN
738 guez 62 tit = 'after clmain'
739     CALL diagetpq(airephy, tit, ip_ebil, 2, 2, dtphys, t_seri, q_seri, &
740 guez 98 ql_seri, u_seri, v_seri, paprs, d_h_vcol, d_qt, d_ec)
741 guez 62 call diagphy(airephy, tit, ip_ebil, zero_v, zero_v, zero_v, zero_v, &
742 guez 98 sens, evap, zero_v, zero_v, ztsol, d_h_vcol, d_qt, d_ec)
743 guez 3 END IF
744    
745 guez 49 ! Update surface temperature:
746 guez 3
747     DO i = 1, klon
748 guez 72 zxtsol(i) = 0.
749     zxfluxlat(i) = 0.
750 guez 3
751 guez 72 zt2m(i) = 0.
752     zq2m(i) = 0.
753     zu10m(i) = 0.
754     zv10m(i) = 0.
755     zxffonte(i) = 0.
756     zxfqcalving(i) = 0.
757 guez 3
758 guez 190 s_pblh(i) = 0.
759     s_lcl(i) = 0.
760 guez 72 s_capCL(i) = 0.
761     s_oliqCL(i) = 0.
762     s_cteiCL(i) = 0.
763     s_pblT(i) = 0.
764     s_therm(i) = 0.
765     s_trmb1(i) = 0.
766     s_trmb2(i) = 0.
767     s_trmb3(i) = 0.
768 guez 3
769 guez 69 IF (abs(pctsrf(i, is_ter) + pctsrf(i, is_lic) + pctsrf(i, is_oce) &
770 guez 190 + pctsrf(i, is_sic) - 1.) > EPSFRA) print *, &
771 guez 97 'physiq : probl\`eme sous surface au point ', i, &
772     pctsrf(i, 1 : nbsrf)
773 guez 3 ENDDO
774     DO nsrf = 1, nbsrf
775     DO i = 1, klon
776     ftsol(i, nsrf) = ftsol(i, nsrf) + d_ts(i, nsrf)
777     zxtsol(i) = zxtsol(i) + ftsol(i, nsrf)*pctsrf(i, nsrf)
778     zxfluxlat(i) = zxfluxlat(i) + fluxlat(i, nsrf)*pctsrf(i, nsrf)
779    
780     zt2m(i) = zt2m(i) + t2m(i, nsrf)*pctsrf(i, nsrf)
781     zq2m(i) = zq2m(i) + q2m(i, nsrf)*pctsrf(i, nsrf)
782     zu10m(i) = zu10m(i) + u10m(i, nsrf)*pctsrf(i, nsrf)
783     zv10m(i) = zv10m(i) + v10m(i, nsrf)*pctsrf(i, nsrf)
784     zxffonte(i) = zxffonte(i) + ffonte(i, nsrf)*pctsrf(i, nsrf)
785 guez 47 zxfqcalving(i) = zxfqcalving(i) + &
786 guez 3 fqcalving(i, nsrf)*pctsrf(i, nsrf)
787     s_pblh(i) = s_pblh(i) + pblh(i, nsrf)*pctsrf(i, nsrf)
788     s_lcl(i) = s_lcl(i) + plcl(i, nsrf)*pctsrf(i, nsrf)
789     s_capCL(i) = s_capCL(i) + capCL(i, nsrf) *pctsrf(i, nsrf)
790     s_oliqCL(i) = s_oliqCL(i) + oliqCL(i, nsrf) *pctsrf(i, nsrf)
791     s_cteiCL(i) = s_cteiCL(i) + cteiCL(i, nsrf) *pctsrf(i, nsrf)
792     s_pblT(i) = s_pblT(i) + pblT(i, nsrf) *pctsrf(i, nsrf)
793     s_therm(i) = s_therm(i) + therm(i, nsrf) *pctsrf(i, nsrf)
794     s_trmb1(i) = s_trmb1(i) + trmb1(i, nsrf) *pctsrf(i, nsrf)
795     s_trmb2(i) = s_trmb2(i) + trmb2(i, nsrf) *pctsrf(i, nsrf)
796     s_trmb3(i) = s_trmb3(i) + trmb3(i, nsrf) *pctsrf(i, nsrf)
797     ENDDO
798     ENDDO
799    
800 guez 97 ! Si une sous-fraction n'existe pas, elle prend la température moyenne :
801 guez 3 DO nsrf = 1, nbsrf
802     DO i = 1, klon
803 guez 47 IF (pctsrf(i, nsrf) < epsfra) ftsol(i, nsrf) = zxtsol(i)
804 guez 3
805 guez 47 IF (pctsrf(i, nsrf) < epsfra) t2m(i, nsrf) = zt2m(i)
806     IF (pctsrf(i, nsrf) < epsfra) q2m(i, nsrf) = zq2m(i)
807     IF (pctsrf(i, nsrf) < epsfra) u10m(i, nsrf) = zu10m(i)
808     IF (pctsrf(i, nsrf) < epsfra) v10m(i, nsrf) = zv10m(i)
809     IF (pctsrf(i, nsrf) < epsfra) ffonte(i, nsrf) = zxffonte(i)
810     IF (pctsrf(i, nsrf) < epsfra) &
811 guez 3 fqcalving(i, nsrf) = zxfqcalving(i)
812 guez 51 IF (pctsrf(i, nsrf) < epsfra) pblh(i, nsrf) = s_pblh(i)
813     IF (pctsrf(i, nsrf) < epsfra) plcl(i, nsrf) = s_lcl(i)
814     IF (pctsrf(i, nsrf) < epsfra) capCL(i, nsrf) = s_capCL(i)
815     IF (pctsrf(i, nsrf) < epsfra) oliqCL(i, nsrf) = s_oliqCL(i)
816     IF (pctsrf(i, nsrf) < epsfra) cteiCL(i, nsrf) = s_cteiCL(i)
817     IF (pctsrf(i, nsrf) < epsfra) pblT(i, nsrf) = s_pblT(i)
818     IF (pctsrf(i, nsrf) < epsfra) therm(i, nsrf) = s_therm(i)
819     IF (pctsrf(i, nsrf) < epsfra) trmb1(i, nsrf) = s_trmb1(i)
820     IF (pctsrf(i, nsrf) < epsfra) trmb2(i, nsrf) = s_trmb2(i)
821     IF (pctsrf(i, nsrf) < epsfra) trmb3(i, nsrf) = s_trmb3(i)
822 guez 3 ENDDO
823     ENDDO
824    
825 guez 98 ! Calculer la dérive du flux infrarouge
826 guez 3
827     DO i = 1, klon
828 guez 190 dlw(i) = - 4. * RSIGMA * zxtsol(i)**3
829 guez 3 ENDDO
830    
831 guez 103 IF (check) print *, "avantcon = ", qcheck(paprs, q_seri, ql_seri)
832    
833 guez 190 ! Appeler la convection
834 guez 3
835 guez 182 if (conv_emanuel) then
836 guez 99 da = 0.
837     mp = 0.
838     phi = 0.
839 guez 97 CALL concvl(dtphys, paprs, play, t_seri, q_seri, u_seri, v_seri, sig1, &
840 guez 183 w01, d_t_con, d_q_con, d_u_con, d_v_con, rain_con, ibas_con, &
841 guez 189 itop_con, upwd, dnwd, dnwd0, Ma, cape, iflagctrl, qcondc, pmflxr, &
842     da, phi, mp)
843 guez 183 snow_con = 0.
844 guez 62 clwcon0 = qcondc
845 guez 71 mfu = upwd + dnwd
846 guez 3
847 guez 103 IF (thermcep) THEN
848     zqsat = MIN(0.5, r2es * FOEEW(t_seri, rtt >= t_seri) / play)
849     zqsat = zqsat / (1. - retv * zqsat)
850     ELSE
851     zqsat = merge(qsats(t_seri), qsatl(t_seri), t_seri < t_coup) / play
852     ENDIF
853 guez 3
854 guez 103 ! Properties of convective clouds
855 guez 71 clwcon0 = fact_cldcon * clwcon0
856 guez 62 call clouds_gno(klon, llm, q_seri, zqsat, clwcon0, ptconv, ratqsc, &
857     rnebcon0)
858 guez 72
859 guez 190 forall (i = 1:klon) ema_pct(i) = paprs(i, itop_con(i) + 1)
860 guez 72 mfd = 0.
861     pen_u = 0.
862     pen_d = 0.
863     pde_d = 0.
864     pde_u = 0.
865 guez 182 else
866     conv_q = d_q_dyn + d_q_vdf / dtphys
867     conv_t = d_t_dyn + d_t_vdf / dtphys
868     z_avant = sum((q_seri + ql_seri) * zmasse, dim=2)
869     CALL conflx(dtphys, paprs, play, t_seri(:, llm:1:- 1), &
870     q_seri(:, llm:1:- 1), conv_t, conv_q, zxfluxq(:, 1), omega, &
871     d_t_con, d_q_con, rain_con, snow_con, mfu(:, llm:1:- 1), &
872     mfd(:, llm:1:- 1), pen_u, pde_u, pen_d, pde_d, kcbot, kctop, &
873     kdtop, pmflxr, pmflxs)
874     WHERE (rain_con < 0.) rain_con = 0.
875     WHERE (snow_con < 0.) snow_con = 0.
876     ibas_con = llm + 1 - kcbot
877     itop_con = llm + 1 - kctop
878 guez 69 END if
879 guez 3
880     DO k = 1, llm
881     DO i = 1, klon
882     t_seri(i, k) = t_seri(i, k) + d_t_con(i, k)
883     q_seri(i, k) = q_seri(i, k) + d_q_con(i, k)
884     u_seri(i, k) = u_seri(i, k) + d_u_con(i, k)
885     v_seri(i, k) = v_seri(i, k) + d_v_con(i, k)
886     ENDDO
887     ENDDO
888    
889 guez 190 IF (if_ebil >= 2) THEN
890 guez 62 tit = 'after convect'
891     CALL diagetpq(airephy, tit, ip_ebil, 2, 2, dtphys, t_seri, q_seri, &
892 guez 98 ql_seri, u_seri, v_seri, paprs, d_h_vcol, d_qt, d_ec)
893 guez 62 call diagphy(airephy, tit, ip_ebil, zero_v, zero_v, zero_v, zero_v, &
894 guez 98 zero_v, zero_v, rain_con, snow_con, ztsol, d_h_vcol, d_qt, d_ec)
895 guez 3 END IF
896    
897     IF (check) THEN
898 guez 98 za = qcheck(paprs, q_seri, ql_seri)
899 guez 62 print *, "aprescon = ", za
900 guez 72 zx_t = 0.
901     za = 0.
902 guez 3 DO i = 1, klon
903     za = za + airephy(i)/REAL(klon)
904     zx_t = zx_t + (rain_con(i)+ &
905     snow_con(i))*airephy(i)/REAL(klon)
906     ENDDO
907 guez 47 zx_t = zx_t/za*dtphys
908 guez 62 print *, "Precip = ", zx_t
909 guez 3 ENDIF
910 guez 69
911 guez 182 IF (.not. conv_emanuel) THEN
912 guez 69 z_apres = sum((q_seri + ql_seri) * zmasse, dim=2)
913     z_factor = (z_avant - (rain_con + snow_con) * dtphys) / z_apres
914 guez 3 DO k = 1, llm
915     DO i = 1, klon
916 guez 52 IF (z_factor(i) > 1. + 1E-8 .OR. z_factor(i) < 1. - 1E-8) THEN
917 guez 3 q_seri(i, k) = q_seri(i, k) * z_factor(i)
918     ENDIF
919     ENDDO
920     ENDDO
921     ENDIF
922    
923 guez 90 ! Convection s\`eche (thermiques ou ajustement)
924 guez 3
925 guez 51 d_t_ajs = 0.
926     d_u_ajs = 0.
927     d_v_ajs = 0.
928     d_q_ajs = 0.
929     fm_therm = 0.
930     entr_therm = 0.
931 guez 3
932 guez 47 if (iflag_thermals == 0) then
933     ! Ajustement sec
934     CALL ajsec(paprs, play, t_seri, q_seri, d_t_ajs, d_q_ajs)
935 guez 13 t_seri = t_seri + d_t_ajs
936     q_seri = q_seri + d_q_ajs
937 guez 3 else
938 guez 47 ! Thermiques
939     call calltherm(dtphys, play, paprs, pphi, u_seri, v_seri, t_seri, &
940     q_seri, d_u_ajs, d_v_ajs, d_t_ajs, d_q_ajs, fm_therm, entr_therm)
941 guez 3 endif
942    
943 guez 190 IF (if_ebil >= 2) THEN
944 guez 62 tit = 'after dry_adjust'
945     CALL diagetpq(airephy, tit, ip_ebil, 2, 2, dtphys, t_seri, q_seri, &
946 guez 98 ql_seri, u_seri, v_seri, paprs, d_h_vcol, d_qt, d_ec)
947 guez 3 END IF
948    
949 guez 47 ! Caclul des ratqs
950 guez 3
951 guez 90 ! ratqs convectifs \`a l'ancienne en fonction de (q(z = 0) - q) / q
952     ! on \'ecrase le tableau ratqsc calcul\'e par clouds_gno
953 guez 3 if (iflag_cldcon == 1) then
954 guez 51 do k = 1, llm
955     do i = 1, klon
956 guez 3 if(ptconv(i, k)) then
957 guez 70 ratqsc(i, k) = ratqsbas + fact_cldcon &
958     * (q_seri(i, 1) - q_seri(i, k)) / q_seri(i, k)
959 guez 3 else
960 guez 51 ratqsc(i, k) = 0.
961 guez 3 endif
962     enddo
963     enddo
964     endif
965    
966 guez 47 ! ratqs stables
967 guez 51 do k = 1, llm
968     do i = 1, klon
969 guez 70 ratqss(i, k) = ratqsbas + (ratqshaut - ratqsbas) &
970 guez 190 * min((paprs(i, 1) - play(i, k)) / (paprs(i, 1) - 3e4), 1.)
971 guez 3 enddo
972     enddo
973    
974 guez 47 ! ratqs final
975 guez 69 if (iflag_cldcon == 1 .or. iflag_cldcon == 2) then
976 guez 47 ! les ratqs sont une conbinaison de ratqss et ratqsc
977     ! ratqs final
978     ! 1e4 (en gros 3 heures), en dur pour le moment, est le temps de
979     ! relaxation des ratqs
980 guez 70 ratqs = max(ratqs * exp(- dtphys * facttemps), ratqss)
981 guez 51 ratqs = max(ratqs, ratqsc)
982 guez 3 else
983 guez 47 ! on ne prend que le ratqs stable pour fisrtilp
984 guez 51 ratqs = ratqss
985 guez 3 endif
986    
987 guez 51 CALL fisrtilp(dtphys, paprs, play, t_seri, q_seri, ptconv, ratqs, &
988     d_t_lsc, d_q_lsc, d_ql_lsc, rneb, cldliq, rain_lsc, snow_lsc, &
989     pfrac_impa, pfrac_nucl, pfrac_1nucl, frac_impa, frac_nucl, prfl, &
990     psfl, rhcl)
991 guez 3
992     WHERE (rain_lsc < 0) rain_lsc = 0.
993     WHERE (snow_lsc < 0) snow_lsc = 0.
994     DO k = 1, llm
995     DO i = 1, klon
996     t_seri(i, k) = t_seri(i, k) + d_t_lsc(i, k)
997     q_seri(i, k) = q_seri(i, k) + d_q_lsc(i, k)
998     ql_seri(i, k) = ql_seri(i, k) + d_ql_lsc(i, k)
999     cldfra(i, k) = rneb(i, k)
1000     IF (.NOT.new_oliq) cldliq(i, k) = ql_seri(i, k)
1001     ENDDO
1002     ENDDO
1003     IF (check) THEN
1004 guez 98 za = qcheck(paprs, q_seri, ql_seri)
1005 guez 62 print *, "apresilp = ", za
1006 guez 72 zx_t = 0.
1007     za = 0.
1008 guez 3 DO i = 1, klon
1009     za = za + airephy(i)/REAL(klon)
1010     zx_t = zx_t + (rain_lsc(i) &
1011     + snow_lsc(i))*airephy(i)/REAL(klon)
1012     ENDDO
1013 guez 47 zx_t = zx_t/za*dtphys
1014 guez 62 print *, "Precip = ", zx_t
1015 guez 3 ENDIF
1016    
1017 guez 190 IF (if_ebil >= 2) THEN
1018 guez 62 tit = 'after fisrt'
1019     CALL diagetpq(airephy, tit, ip_ebil, 2, 2, dtphys, t_seri, q_seri, &
1020 guez 98 ql_seri, u_seri, v_seri, paprs, d_h_vcol, d_qt, d_ec)
1021 guez 62 call diagphy(airephy, tit, ip_ebil, zero_v, zero_v, zero_v, zero_v, &
1022 guez 98 zero_v, zero_v, rain_lsc, snow_lsc, ztsol, d_h_vcol, d_qt, d_ec)
1023 guez 3 END IF
1024    
1025 guez 47 ! PRESCRIPTION DES NUAGES POUR LE RAYONNEMENT
1026 guez 3
1027     ! 1. NUAGES CONVECTIFS
1028    
1029 guez 174 IF (iflag_cldcon <= - 1) THEN
1030 guez 62 ! seulement pour Tiedtke
1031 guez 51 snow_tiedtke = 0.
1032 guez 174 if (iflag_cldcon == - 1) then
1033 guez 51 rain_tiedtke = rain_con
1034 guez 3 else
1035 guez 51 rain_tiedtke = 0.
1036     do k = 1, llm
1037     do i = 1, klon
1038 guez 7 if (d_q_con(i, k) < 0.) then
1039 guez 174 rain_tiedtke(i) = rain_tiedtke(i) - d_q_con(i, k)/dtphys &
1040 guez 17 *zmasse(i, k)
1041 guez 3 endif
1042     enddo
1043     enddo
1044     endif
1045    
1046     ! Nuages diagnostiques pour Tiedtke
1047 guez 69 CALL diagcld1(paprs, play, rain_tiedtke, snow_tiedtke, ibas_con, &
1048     itop_con, diafra, dialiq)
1049 guez 3 DO k = 1, llm
1050     DO i = 1, klon
1051 guez 51 IF (diafra(i, k) > cldfra(i, k)) THEN
1052 guez 3 cldliq(i, k) = dialiq(i, k)
1053     cldfra(i, k) = diafra(i, k)
1054     ENDIF
1055     ENDDO
1056     ENDDO
1057     ELSE IF (iflag_cldcon == 3) THEN
1058 guez 72 ! On prend pour les nuages convectifs le maximum du calcul de
1059 guez 90 ! la convection et du calcul du pas de temps pr\'ec\'edent diminu\'e
1060 guez 72 ! d'un facteur facttemps.
1061     facteur = dtphys * facttemps
1062 guez 51 do k = 1, llm
1063     do i = 1, klon
1064 guez 70 rnebcon(i, k) = rnebcon(i, k) * facteur
1065 guez 72 if (rnebcon0(i, k) * clwcon0(i, k) &
1066     > rnebcon(i, k) * clwcon(i, k)) then
1067 guez 51 rnebcon(i, k) = rnebcon0(i, k)
1068     clwcon(i, k) = clwcon0(i, k)
1069 guez 3 endif
1070     enddo
1071     enddo
1072    
1073 guez 47 ! On prend la somme des fractions nuageuses et des contenus en eau
1074 guez 51 cldfra = min(max(cldfra, rnebcon), 1.)
1075     cldliq = cldliq + rnebcon*clwcon
1076 guez 3 ENDIF
1077    
1078 guez 51 ! 2. Nuages stratiformes
1079 guez 3
1080     IF (ok_stratus) THEN
1081 guez 47 CALL diagcld2(paprs, play, t_seri, q_seri, diafra, dialiq)
1082 guez 3 DO k = 1, llm
1083     DO i = 1, klon
1084 guez 51 IF (diafra(i, k) > cldfra(i, k)) THEN
1085 guez 3 cldliq(i, k) = dialiq(i, k)
1086     cldfra(i, k) = diafra(i, k)
1087     ENDIF
1088     ENDDO
1089     ENDDO
1090     ENDIF
1091    
1092     ! Precipitation totale
1093     DO i = 1, klon
1094     rain_fall(i) = rain_con(i) + rain_lsc(i)
1095     snow_fall(i) = snow_con(i) + snow_lsc(i)
1096     ENDDO
1097    
1098 guez 62 IF (if_ebil >= 2) CALL diagetpq(airephy, "after diagcld", ip_ebil, 2, 2, &
1099 guez 98 dtphys, t_seri, q_seri, ql_seri, u_seri, v_seri, paprs, d_h_vcol, &
1100     d_qt, d_ec)
1101 guez 3
1102 guez 90 ! Humidit\'e relative pour diagnostic :
1103 guez 3 DO k = 1, llm
1104     DO i = 1, klon
1105     zx_t = t_seri(i, k)
1106     IF (thermcep) THEN
1107 guez 103 zx_qs = r2es * FOEEW(zx_t, rtt >= zx_t)/play(i, k)
1108 guez 47 zx_qs = MIN(0.5, zx_qs)
1109 guez 174 zcor = 1./(1. - retv*zx_qs)
1110 guez 47 zx_qs = zx_qs*zcor
1111 guez 3 ELSE
1112 guez 7 IF (zx_t < t_coup) THEN
1113 guez 47 zx_qs = qsats(zx_t)/play(i, k)
1114 guez 3 ELSE
1115 guez 47 zx_qs = qsatl(zx_t)/play(i, k)
1116 guez 3 ENDIF
1117     ENDIF
1118     zx_rh(i, k) = q_seri(i, k)/zx_qs
1119 guez 51 zqsat(i, k) = zx_qs
1120 guez 3 ENDDO
1121     ENDDO
1122 guez 52
1123     ! Introduce the aerosol direct and first indirect radiative forcings:
1124     IF (ok_ade .OR. ok_aie) THEN
1125 guez 68 ! Get sulfate aerosol distribution :
1126 guez 130 CALL readsulfate(dayvrai, time, firstcal, sulfate)
1127     CALL readsulfate_preind(dayvrai, time, firstcal, sulfate_pi)
1128 guez 3
1129 guez 52 CALL aeropt(play, paprs, t_seri, sulfate, rhcl, tau_ae, piz_ae, cg_ae, &
1130     aerindex)
1131 guez 3 ELSE
1132 guez 52 tau_ae = 0.
1133     piz_ae = 0.
1134     cg_ae = 0.
1135 guez 3 ENDIF
1136    
1137 guez 97 ! Param\`etres optiques des nuages et quelques param\`etres pour
1138     ! diagnostics :
1139 guez 3 if (ok_newmicro) then
1140 guez 69 CALL newmicro(paprs, play, t_seri, cldliq, cldfra, cldtau, cldemi, &
1141     cldh, cldl, cldm, cldt, cldq, flwp, fiwp, flwc, fiwc, ok_aie, &
1142     sulfate, sulfate_pi, bl95_b0, bl95_b1, cldtaupi, re, fl)
1143 guez 3 else
1144 guez 52 CALL nuage(paprs, play, t_seri, cldliq, cldfra, cldtau, cldemi, cldh, &
1145     cldl, cldm, cldt, cldq, ok_aie, sulfate, sulfate_pi, bl95_b0, &
1146     bl95_b1, cldtaupi, re, fl)
1147 guez 3 endif
1148    
1149 guez 154 IF (MOD(itap - 1, radpas) == 0) THEN
1150 guez 118 ! Appeler le rayonnement mais calculer tout d'abord l'albedo du sol.
1151 guez 155 ! Calcul de l'abedo moyen par maille
1152     albsol = sum(falbe * pctsrf, dim = 2)
1153    
1154 guez 62 ! Rayonnement (compatible Arpege-IFS) :
1155 guez 155 CALL radlwsw(dist, mu0, fract, paprs, play, zxtsol, albsol, t_seri, &
1156     q_seri, wo, cldfra, cldemi, cldtau, heat, heat0, cool, cool0, &
1157     radsol, albpla, topsw, toplw, solsw, sollw, sollwdown, topsw0, &
1158     toplw0, solsw0, sollw0, lwdn0, lwdn, lwup0, lwup, swdn0, swdn, &
1159     swup0, swup, ok_ade, ok_aie, tau_ae, piz_ae, cg_ae, topswad, &
1160     solswad, cldtaupi, topswai, solswai)
1161 guez 3 ENDIF
1162 guez 118
1163 guez 3 ! Ajouter la tendance des rayonnements (tous les pas)
1164    
1165     DO k = 1, llm
1166     DO i = 1, klon
1167 guez 174 t_seri(i, k) = t_seri(i, k) + (heat(i, k) - cool(i, k)) * dtphys/86400.
1168 guez 3 ENDDO
1169     ENDDO
1170    
1171 guez 190 IF (if_ebil >= 2) THEN
1172 guez 62 tit = 'after rad'
1173     CALL diagetpq(airephy, tit, ip_ebil, 2, 2, dtphys, t_seri, q_seri, &
1174 guez 98 ql_seri, u_seri, v_seri, paprs, d_h_vcol, d_qt, d_ec)
1175 guez 62 call diagphy(airephy, tit, ip_ebil, topsw, toplw, solsw, sollw, &
1176 guez 98 zero_v, zero_v, zero_v, zero_v, ztsol, d_h_vcol, d_qt, d_ec)
1177 guez 3 END IF
1178    
1179     ! Calculer l'hydrologie de la surface
1180     DO i = 1, klon
1181 guez 72 zxqsurf(i) = 0.
1182     zxsnow(i) = 0.
1183 guez 3 ENDDO
1184     DO nsrf = 1, nbsrf
1185     DO i = 1, klon
1186     zxqsurf(i) = zxqsurf(i) + fqsurf(i, nsrf)*pctsrf(i, nsrf)
1187     zxsnow(i) = zxsnow(i) + fsnow(i, nsrf)*pctsrf(i, nsrf)
1188     ENDDO
1189     ENDDO
1190    
1191 guez 90 ! Calculer le bilan du sol et la d\'erive de temp\'erature (couplage)
1192 guez 3
1193     DO i = 1, klon
1194     bils(i) = radsol(i) - sens(i) + zxfluxlat(i)
1195     ENDDO
1196    
1197 guez 90 ! Param\'etrisation de l'orographie \`a l'\'echelle sous-maille :
1198 guez 3
1199     IF (ok_orodr) THEN
1200 guez 174 ! S\'election des points pour lesquels le sch\'ema est actif :
1201 guez 51 igwd = 0
1202     DO i = 1, klon
1203     itest(i) = 0
1204 guez 174 IF (zpic(i) - zmea(i) > 100. .AND. zstd(i) > 10.) THEN
1205 guez 51 itest(i) = 1
1206     igwd = igwd + 1
1207 guez 3 ENDIF
1208     ENDDO
1209    
1210 guez 51 CALL drag_noro(klon, llm, dtphys, paprs, play, zmea, zstd, zsig, zgam, &
1211 guez 150 zthe, zpic, zval, itest, t_seri, u_seri, v_seri, zulow, zvlow, &
1212     zustrdr, zvstrdr, d_t_oro, d_u_oro, d_v_oro)
1213 guez 3
1214 guez 47 ! ajout des tendances
1215 guez 3 DO k = 1, llm
1216     DO i = 1, klon
1217     t_seri(i, k) = t_seri(i, k) + d_t_oro(i, k)
1218     u_seri(i, k) = u_seri(i, k) + d_u_oro(i, k)
1219     v_seri(i, k) = v_seri(i, k) + d_v_oro(i, k)
1220     ENDDO
1221     ENDDO
1222 guez 13 ENDIF
1223 guez 3
1224     IF (ok_orolf) THEN
1225 guez 90 ! S\'election des points pour lesquels le sch\'ema est actif :
1226 guez 51 igwd = 0
1227     DO i = 1, klon
1228     itest(i) = 0
1229 guez 174 IF (zpic(i) - zmea(i) > 100.) THEN
1230 guez 51 itest(i) = 1
1231     igwd = igwd + 1
1232 guez 3 ENDIF
1233     ENDDO
1234    
1235 guez 47 CALL lift_noro(klon, llm, dtphys, paprs, play, rlat, zmea, zstd, zpic, &
1236     itest, t_seri, u_seri, v_seri, zulow, zvlow, zustrli, zvstrli, &
1237 guez 3 d_t_lif, d_u_lif, d_v_lif)
1238    
1239 guez 51 ! Ajout des tendances :
1240 guez 3 DO k = 1, llm
1241     DO i = 1, klon
1242     t_seri(i, k) = t_seri(i, k) + d_t_lif(i, k)
1243     u_seri(i, k) = u_seri(i, k) + d_u_lif(i, k)
1244     v_seri(i, k) = v_seri(i, k) + d_v_lif(i, k)
1245     ENDDO
1246     ENDDO
1247 guez 49 ENDIF
1248 guez 3
1249 guez 90 ! Stress n\'ecessaires : toute la physique
1250 guez 3
1251     DO i = 1, klon
1252 guez 51 zustrph(i) = 0.
1253     zvstrph(i) = 0.
1254 guez 3 ENDDO
1255     DO k = 1, llm
1256     DO i = 1, klon
1257 guez 62 zustrph(i) = zustrph(i) + (u_seri(i, k) - u(i, k)) / dtphys &
1258     * zmasse(i, k)
1259     zvstrph(i) = zvstrph(i) + (v_seri(i, k) - v(i, k)) / dtphys &
1260     * zmasse(i, k)
1261 guez 3 ENDDO
1262     ENDDO
1263    
1264 guez 171 CALL aaam_bud(rg, romega, rlat, rlon, pphis, zustrdr, zustrli, zustrph, &
1265     zvstrdr, zvstrli, zvstrph, paprs, u, v, aam, torsfc)
1266 guez 3
1267 guez 62 IF (if_ebil >= 2) CALL diagetpq(airephy, 'after orography', ip_ebil, 2, &
1268 guez 98 2, dtphys, t_seri, q_seri, ql_seri, u_seri, v_seri, paprs, d_h_vcol, &
1269     d_qt, d_ec)
1270 guez 3
1271 guez 47 ! Calcul des tendances traceurs
1272 guez 118 call phytrac(itap, lmt_pas, julien, time, firstcal, lafin, dtphys, t, &
1273 guez 98 paprs, play, mfu, mfd, pde_u, pen_d, ycoefh, fm_therm, entr_therm, &
1274 guez 159 yu1, yv1, ftsol, pctsrf, frac_impa, frac_nucl, da, phi, mp, upwd, &
1275 guez 175 dnwd, tr_seri, zmasse, ncid_startphy, nid_ins, itau_phy)
1276 guez 3
1277 guez 190 IF (offline) call phystokenc(dtphys, t, mfu, mfd, pen_u, pde_u, pen_d, &
1278     pde_d, fm_therm, entr_therm, ycoefh, yu1, yv1, ftsol, pctsrf, &
1279     frac_impa, frac_nucl, pphis, airephy, dtphys, itap)
1280 guez 3
1281     ! Calculer le transport de l'eau et de l'energie (diagnostique)
1282 guez 171 CALL transp(paprs, t_seri, q_seri, u_seri, v_seri, zphi, ve, vq, ue, uq)
1283 guez 3
1284 guez 31 ! diag. bilKP
1285 guez 3
1286 guez 178 CALL transp_lay(paprs, t_seri, q_seri, u_seri, v_seri, zphi, &
1287 guez 3 ve_lay, vq_lay, ue_lay, uq_lay)
1288    
1289     ! Accumuler les variables a stocker dans les fichiers histoire:
1290    
1291 guez 51 ! conversion Ec -> E thermique
1292 guez 3 DO k = 1, llm
1293     DO i = 1, klon
1294 guez 51 ZRCPD = RCPD * (1. + RVTMP2 * q_seri(i, k))
1295     d_t_ec(i, k) = 0.5 / ZRCPD &
1296     * (u(i, k)**2 + v(i, k)**2 - u_seri(i, k)**2 - v_seri(i, k)**2)
1297     t_seri(i, k) = t_seri(i, k) + d_t_ec(i, k)
1298     d_t_ec(i, k) = d_t_ec(i, k) / dtphys
1299 guez 3 END DO
1300     END DO
1301 guez 51
1302 guez 190 IF (if_ebil >= 1) THEN
1303 guez 62 tit = 'after physic'
1304     CALL diagetpq(airephy, tit, ip_ebil, 1, 1, dtphys, t_seri, q_seri, &
1305 guez 98 ql_seri, u_seri, v_seri, paprs, d_h_vcol, d_qt, d_ec)
1306 guez 190 ! Comme les tendances de la physique sont ajoute dans la dynamique,
1307 guez 47 ! on devrait avoir que la variation d'entalpie par la dynamique
1308     ! est egale a la variation de la physique au pas de temps precedent.
1309     ! Donc la somme de ces 2 variations devrait etre nulle.
1310 guez 62 call diagphy(airephy, tit, ip_ebil, topsw, toplw, solsw, sollw, sens, &
1311 guez 98 evap, rain_fall, snow_fall, ztsol, d_h_vcol, d_qt, d_ec)
1312 guez 51 d_h_vcol_phy = d_h_vcol
1313 guez 3 END IF
1314    
1315 guez 47 ! SORTIES
1316 guez 3
1317 guez 69 ! prw = eau precipitable
1318 guez 3 DO i = 1, klon
1319     prw(i) = 0.
1320     DO k = 1, llm
1321 guez 17 prw(i) = prw(i) + q_seri(i, k)*zmasse(i, k)
1322 guez 3 ENDDO
1323     ENDDO
1324    
1325     ! Convertir les incrementations en tendances
1326    
1327     DO k = 1, llm
1328     DO i = 1, klon
1329 guez 49 d_u(i, k) = (u_seri(i, k) - u(i, k)) / dtphys
1330     d_v(i, k) = (v_seri(i, k) - v(i, k)) / dtphys
1331     d_t(i, k) = (t_seri(i, k) - t(i, k)) / dtphys
1332     d_qx(i, k, ivap) = (q_seri(i, k) - qx(i, k, ivap)) / dtphys
1333     d_qx(i, k, iliq) = (ql_seri(i, k) - qx(i, k, iliq)) / dtphys
1334 guez 3 ENDDO
1335     ENDDO
1336    
1337 guez 98 DO iq = 3, nqmx
1338     DO k = 1, llm
1339     DO i = 1, klon
1340 guez 174 d_qx(i, k, iq) = (tr_seri(i, k, iq - 2) - qx(i, k, iq)) / dtphys
1341 guez 3 ENDDO
1342     ENDDO
1343 guez 98 ENDDO
1344 guez 3
1345     ! Sauvegarder les valeurs de t et q a la fin de la physique:
1346     DO k = 1, llm
1347     DO i = 1, klon
1348     t_ancien(i, k) = t_seri(i, k)
1349     q_ancien(i, k) = q_seri(i, k)
1350     ENDDO
1351     ENDDO
1352    
1353     call write_histins
1354    
1355 guez 157 IF (lafin) then
1356     call NF95_CLOSE(ncid_startphy)
1357 guez 175 CALL phyredem(pctsrf, ftsol, ftsoil, fqsurf, qsol, &
1358 guez 157 fsnow, falbe, fevap, rain_fall, snow_fall, solsw, sollw, dlw, &
1359     radsol, frugs, agesno, zmea, zstd, zsig, zgam, zthe, zpic, zval, &
1360     t_ancien, q_ancien, rnebcon, ratqs, clwcon, run_off_lic_0, sig1, &
1361     w01)
1362     end IF
1363 guez 3
1364 guez 35 firstcal = .FALSE.
1365    
1366 guez 3 contains
1367    
1368     subroutine write_histins
1369    
1370 guez 47 ! From phylmd/write_histins.h, version 1.2 2005/05/25 13:10:09
1371 guez 3
1372 guez 157 ! Ecriture des sorties
1373    
1374 guez 189 use gr_phy_write_m, only: gr_phy_write
1375 guez 90 USE histsync_m, ONLY: histsync
1376     USE histwrite_m, ONLY: histwrite
1377    
1378 guez 190 integer itau_w ! pas de temps d'\'ecriture
1379 guez 3
1380     !--------------------------------------------------
1381    
1382     IF (ok_instan) THEN
1383     itau_w = itau_phy + itap
1384 guez 190 CALL histwrite(nid_ins, "phis", itau_w, gr_phy_write(pphis))
1385     CALL histwrite(nid_ins, "aire", itau_w, gr_phy_write(airephy))
1386     CALL histwrite(nid_ins, "psol", itau_w, gr_phy_write(paprs(:, 1)))
1387     CALL histwrite(nid_ins, "precip", itau_w, &
1388     gr_phy_write(rain_fall + snow_fall))
1389     CALL histwrite(nid_ins, "plul", itau_w, &
1390     gr_phy_write(rain_lsc + snow_lsc))
1391     CALL histwrite(nid_ins, "pluc", itau_w, &
1392     gr_phy_write(rain_con + snow_con))
1393     CALL histwrite(nid_ins, "tsol", itau_w, gr_phy_write(zxtsol))
1394     CALL histwrite(nid_ins, "t2m", itau_w, gr_phy_write(zt2m))
1395     CALL histwrite(nid_ins, "q2m", itau_w, gr_phy_write(zq2m))
1396     CALL histwrite(nid_ins, "u10m", itau_w, gr_phy_write(zu10m))
1397     CALL histwrite(nid_ins, "v10m", itau_w, gr_phy_write(zv10m))
1398     CALL histwrite(nid_ins, "snow", itau_w, gr_phy_write(snow_fall))
1399     CALL histwrite(nid_ins, "cdrm", itau_w, gr_phy_write(cdragm))
1400     CALL histwrite(nid_ins, "cdrh", itau_w, gr_phy_write(cdragh))
1401     CALL histwrite(nid_ins, "topl", itau_w, gr_phy_write(toplw))
1402     CALL histwrite(nid_ins, "evap", itau_w, gr_phy_write(evap))
1403     CALL histwrite(nid_ins, "sols", itau_w, gr_phy_write(solsw))
1404     CALL histwrite(nid_ins, "soll", itau_w, gr_phy_write(sollw))
1405     CALL histwrite(nid_ins, "solldown", itau_w, gr_phy_write(sollwdown))
1406     CALL histwrite(nid_ins, "bils", itau_w, gr_phy_write(bils))
1407     CALL histwrite(nid_ins, "sens", itau_w, gr_phy_write(- sens))
1408     CALL histwrite(nid_ins, "fder", itau_w, gr_phy_write(fder))
1409     CALL histwrite(nid_ins, "dtsvdfo", itau_w, &
1410     gr_phy_write(d_ts(:, is_oce)))
1411     CALL histwrite(nid_ins, "dtsvdft", itau_w, &
1412     gr_phy_write(d_ts(:, is_ter)))
1413     CALL histwrite(nid_ins, "dtsvdfg", itau_w, &
1414     gr_phy_write(d_ts(:, is_lic)))
1415     CALL histwrite(nid_ins, "dtsvdfi", itau_w, &
1416     gr_phy_write(d_ts(:, is_sic)))
1417 guez 3
1418     DO nsrf = 1, nbsrf
1419     CALL histwrite(nid_ins, "pourc_"//clnsurf(nsrf), itau_w, &
1420 guez 190 gr_phy_write(pctsrf(:, nsrf)*100.))
1421 guez 3 CALL histwrite(nid_ins, "fract_"//clnsurf(nsrf), itau_w, &
1422 guez 190 gr_phy_write(pctsrf(:, nsrf)))
1423 guez 3 CALL histwrite(nid_ins, "sens_"//clnsurf(nsrf), itau_w, &
1424 guez 190 gr_phy_write(fluxt(:, 1, nsrf)))
1425 guez 3 CALL histwrite(nid_ins, "lat_"//clnsurf(nsrf), itau_w, &
1426 guez 190 gr_phy_write(fluxlat(:, nsrf)))
1427 guez 3 CALL histwrite(nid_ins, "tsol_"//clnsurf(nsrf), itau_w, &
1428 guez 190 gr_phy_write(ftsol(:, nsrf)))
1429 guez 3 CALL histwrite(nid_ins, "taux_"//clnsurf(nsrf), itau_w, &
1430 guez 190 gr_phy_write(fluxu(:, 1, nsrf)))
1431 guez 3 CALL histwrite(nid_ins, "tauy_"//clnsurf(nsrf), itau_w, &
1432 guez 190 gr_phy_write(fluxv(:, 1, nsrf)))
1433 guez 3 CALL histwrite(nid_ins, "rugs_"//clnsurf(nsrf), itau_w, &
1434 guez 190 gr_phy_write(frugs(:, nsrf)))
1435 guez 3 CALL histwrite(nid_ins, "albe_"//clnsurf(nsrf), itau_w, &
1436 guez 190 gr_phy_write(falbe(:, nsrf)))
1437 guez 3 END DO
1438    
1439 guez 190 CALL histwrite(nid_ins, "albs", itau_w, gr_phy_write(albsol))
1440     CALL histwrite(nid_ins, "rugs", itau_w, gr_phy_write(zxrugs))
1441     CALL histwrite(nid_ins, "s_pblh", itau_w, gr_phy_write(s_pblh))
1442     CALL histwrite(nid_ins, "s_pblt", itau_w, gr_phy_write(s_pblt))
1443     CALL histwrite(nid_ins, "s_lcl", itau_w, gr_phy_write(s_lcl))
1444     CALL histwrite(nid_ins, "s_capCL", itau_w, gr_phy_write(s_capCL))
1445     CALL histwrite(nid_ins, "s_oliqCL", itau_w, gr_phy_write(s_oliqCL))
1446     CALL histwrite(nid_ins, "s_cteiCL", itau_w, gr_phy_write(s_cteiCL))
1447     CALL histwrite(nid_ins, "s_therm", itau_w, gr_phy_write(s_therm))
1448     CALL histwrite(nid_ins, "s_trmb1", itau_w, gr_phy_write(s_trmb1))
1449     CALL histwrite(nid_ins, "s_trmb2", itau_w, gr_phy_write(s_trmb2))
1450     CALL histwrite(nid_ins, "s_trmb3", itau_w, gr_phy_write(s_trmb3))
1451     if (conv_emanuel) CALL histwrite(nid_ins, "ptop", itau_w, &
1452     gr_phy_write(ema_pct))
1453     CALL histwrite(nid_ins, "temp", itau_w, gr_phy_write(t_seri))
1454     CALL histwrite(nid_ins, "vitu", itau_w, gr_phy_write(u_seri))
1455     CALL histwrite(nid_ins, "vitv", itau_w, gr_phy_write(v_seri))
1456     CALL histwrite(nid_ins, "geop", itau_w, gr_phy_write(zphi))
1457     CALL histwrite(nid_ins, "pres", itau_w, gr_phy_write(play))
1458     CALL histwrite(nid_ins, "dtvdf", itau_w, gr_phy_write(d_t_vdf))
1459     CALL histwrite(nid_ins, "dqvdf", itau_w, gr_phy_write(d_q_vdf))
1460     CALL histwrite(nid_ins, "rhum", itau_w, gr_phy_write(zx_rh))
1461 guez 91 call histsync(nid_ins)
1462 guez 3 ENDIF
1463    
1464     end subroutine write_histins
1465    
1466     END SUBROUTINE physiq
1467    
1468     end module physiq_m

  ViewVC Help
Powered by ViewVC 1.1.21