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

Annotation of /trunk/phylmd/physiq.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 97 - (hide annotations)
Fri Apr 25 14:58:31 2014 UTC (10 years ago) by guez
File size: 65089 byte(s)
Module pressure_var is now only used in gcm. Created local variables
pls and p3d in etat0, added argument p3d to regr_pr_o3.

In leapfrog, moved computation of p3d and exner function immediately
after integrd, for clarity (does not change the execution).

Removed unused arguments: ntra, tra1 and tra of cv3_compress; ntra,
tra and traent of cv3_mixing; ntra, ftra, ftra1 of cv3_uncompress;
ntra, tra, trap of cv3_unsat; ntra, tra, trap, traent, ftra of
cv3_yield; tra, tvp, pbase, bbase, dtvpdt1, dtvpdq1, dplcldt,
dplcldr, ntra of concvl; ndp1, ntra, tra1 of cv_driver

Removed argument d_tra and computation of d_tra in concvl. Removed
argument ftra1 and computation of ftra1 in cv_driver. ftra1 was just
set to 0 in cv_driver, associated to d_tra in concvl, and set again to
zero in concvl.

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

  ViewVC Help
Powered by ViewVC 1.1.21