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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 91 - (hide annotations)
Wed Mar 26 17:18:58 2014 UTC (10 years, 1 month ago) by guez
Original Path: trunk/phylmd/physiq.f
File size: 66008 byte(s)
Removed unused variables lock_startdate and time_stamp of module
calendar.

Noticed that physiq does not change the surface pressure. So removed
arguments ps and dpfi of subroutine addfi. dpfi was always 0. The
computation of ps in addfi included some averaging at the poles. In
principle, this does not change ps but in practice it does because of
finite numerical precision. So the results of the simulation are
changed. Removed arguments ps and dpfi of calfis. Removed argument
d_ps of physiq.

du at the poles is not computed by dudv1, so declare only the
corresponding latitudes in dudv1. caldyn passes only a section of the
array dudyn as argument.

Removed variable niadv of module iniadvtrac_m.

Declared arguments of exner_hyb as assumed-shape arrays and made all
other horizontal sizes in exner_hyb dynamic. This allows the external
program test_disvert to use exner_hyb at a single horizontal position.

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

  ViewVC Help
Powered by ViewVC 1.1.21