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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 70 - (hide annotations)
Mon Jun 24 15:39:52 2013 UTC (10 years, 10 months ago) by guez
Original Path: trunk/libf/phylmd/physiq.f90
File size: 68729 byte(s)
In procedure, "addfi" access directly the module variable "dtphys"
instead of going through an argument.

In "conflx", do not create a local variable for temperature with
reversed order of vertical levels. Instead, give an actual argument
with reversed order in "physiq".

Changed names of variables "rmd" and "rmv" from module "suphec_m" to
"md" and "mv".

In "hgardfou", print only the first temperature out of range found.

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

  ViewVC Help
Powered by ViewVC 1.1.21