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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 103 - (hide annotations)
Fri Aug 29 13:00:05 2014 UTC (9 years, 8 months ago) by guez
Original Path: trunk/phylmd/physiq.f
File size: 61734 byte(s)
Renamed module cvparam to cv_param. Deleted procedure
cv_param. Changed variables of module cv_param into parameters.

In procedures cv_driver, cv_uncompress and cv3_uncompress, removed
some arguments giving dimensions and used module variables klon and
klev instead.

In procedures gradiv2, laplacien_gam and laplacien, changed
declarations of local variables because klevel is not always klev.

Removed code for nudging surface pressure.

Removed arguments pim and pjm of tau2alpha. Added assignment of false
to variable first.

Replaced real argument del of procedures foeew and FOEDE by logical
argument.

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

  ViewVC Help
Powered by ViewVC 1.1.21