source: LMDZ.3.3/branches/rel-LF/libf/phylmd/physiq.F @ 230

Last change on this file since 230 was 230, checked in by lmdzadmin, 23 years ago

Merge de la physique avec la branche principale
LF

  • Property svn:eol-style set to native
  • Property svn:executable set to *
  • Property svn:keywords set to Author Date Id Revision
File size: 109.7 KB
Line 
1c
2c $Header$
3c
4      SUBROUTINE physiq (nlon,nlev,nqmax  ,
5     .            debut,lafin,rjourvrai,rjour_ecri,gmtime,pdtphys,
6     .            paprs,pplay,pphi,pphis,paire,presnivs,clesphy0,
7     .            u,v,t,qx,
8     .            omega, cufi, cvfi,
9     .            d_u, d_v, d_t, d_qx, d_ps)
10      USE ioipsl
11      USE histcom
12
13      IMPLICIT none
14c======================================================================
15c
16c Auteur(s) Z.X. Li (LMD/CNRS) date: 19930818
17c
18c Objet: Moniteur general de la physique du modele
19cAA      Modifications quant aux traceurs :
20cAA                  -  uniformisation des parametrisations ds phytrac
21cAA                  -  stockage des moyennes des champs necessaires
22cAA                     en mode traceur off-line
23c======================================================================
24c    modif   ( P. Le Van ,  12/10/98 )
25c
26c  Arguments:
27c
28c nlon----input-I-nombre de points horizontaux
29c nlev----input-I-nombre de couches verticales
30c nqmax---input-I-nombre de traceurs (y compris vapeur d'eau) = 1
31c debut---input-L-variable logique indiquant le premier passage
32c lafin---input-L-variable logique indiquant le dernier passage
33c rjour---input-R-numero du jour de l'experience
34c gmtime--input-R-temps universel dans la journee (0 a 86400 s)
35c pdtphys-input-R-pas d'integration pour la physique (seconde)
36c paprs---input-R-pression pour chaque inter-couche (en Pa)
37c pplay---input-R-pression pour le mileu de chaque couche (en Pa)
38c pphi----input-R-geopotentiel de chaque couche (g z) (reference sol)
39c pphis---input-R-geopotentiel du sol
40c paire---input-R-aire de chaque maille
41c presnivs-input_R_pressions approximat. des milieux couches ( en PA)
42c u-------input-R-vitesse dans la direction X (de O a E) en m/s
43c v-------input-R-vitesse Y (de S a N) en m/s
44c t-------input-R-temperature (K)
45c qx------input-R-humidite specifique (kg/kg) et d'autres traceurs
46c d_t_dyn-input-R-tendance dynamique pour "t" (K/s)
47c d_q_dyn-input-R-tendance dynamique pour "q" (kg/kg/s)
48c omega---input-R-vitesse verticale en Pa/s
49c cufi----input-R-resolution des mailles en x (m)
50c cvfi----input-R-resolution des mailles en y (m)
51c
52c d_u-----output-R-tendance physique de "u" (m/s/s)
53c d_v-----output-R-tendance physique de "v" (m/s/s)
54c d_t-----output-R-tendance physique de "t" (K/s)
55c d_qx----output-R-tendance physique de "qx" (kg/kg/s)
56c d_ps----output-R-tendance physique de la pression au sol
57c======================================================================
58#include "dimensions.h"
59      integer jjmp1
60      parameter (jjmp1=jjm+1-1/jjm)
61#include "dimphy.h"
62#include "regdim.h"
63#include "indicesol.h"
64#include "dimsoil.h"
65#include "clesphys.h"
66#include "control.h"
67#include "temps.h"
68c======================================================================
69      LOGICAL check ! Verifier la conservation du modele en eau
70      PARAMETER (check=.FALSE.)
71      LOGICAL ok_stratus ! Ajouter artificiellement les stratus
72      PARAMETER (ok_stratus=.FALSE.)
73c======================================================================
74c Parametres lies au coupleur OASIS:
75#include "oasis.h"
76      INTEGER,SAVE :: npas, nexca
77      logical rnpb
78      parameter(rnpb=.true.)
79c      PARAMETER (npas=1440)
80c      PARAMETER (nexca=48)
81      EXTERNAL fromcpl, intocpl, inicma
82c      ocean = type de modele ocean a utiliser: force, slab, couple
83      character*6 ocean
84      SAVE ocean
85
86c      parameter (ocean = 'force ')
87c     parameter (ocean = 'couple')
88      logical ok_ocean
89c======================================================================
90c Clef controlant l'activation du cycle diurne:
91ccc      LOGICAL cycle_diurne
92ccc      PARAMETER (cycle_diurne=.FALSE.)
93c======================================================================
94c Modele thermique du sol, a activer pour le cycle diurne:
95ccc      LOGICAL soil_model
96ccc      PARAMETER (soil_model=.FALSE.)
97      logical ok_veget
98      save ok_veget
99c     parameter (ok_veget = .true.)
100c      parameter (ok_veget = .false.)
101c======================================================================
102c Dans les versions precedentes, l'eau liquide nuageuse utilisee dans
103c le calcul du rayonnement est celle apres la precipitation des nuages.
104c Si cette cle new_oliq est activee, ce sera une valeur moyenne entre
105c la condensation et la precipitation. Cette cle augmente les impacts
106c radiatifs des nuages.
107ccc      LOGICAL new_oliq
108ccc      PARAMETER (new_oliq=.FALSE.)
109c======================================================================
110c Clefs controlant deux parametrisations de l'orographie:
111cc      LOGICAL ok_orodr
112ccc      PARAMETER (ok_orodr=.FALSE.)
113ccc      LOGICAL ok_orolf
114ccc      PARAMETER (ok_orolf=.FALSE.)
115c======================================================================
116      LOGICAL ok_journe ! sortir le fichier journalier
117      save ok_journe
118c      PARAMETER (ok_journe=.true.)
119c
120      LOGICAL ok_mensuel ! sortir le fichier mensuel
121      save ok_mensuel
122c      PARAMETER (ok_mensuel=.true.)
123c
124      LOGICAL ok_instan ! sortir le fichier instantane
125      save ok_instan
126c      PARAMETER (ok_instan=.true.)
127c
128      LOGICAL ok_region ! sortir le fichier regional
129      PARAMETER (ok_region=.FALSE.)
130c======================================================================
131c
132      INTEGER ivap          ! indice de traceurs pour vapeur d'eau
133      PARAMETER (ivap=1)
134      INTEGER iliq          ! indice de traceurs pour eau liquide
135      PARAMETER (iliq=2)
136
137      INTEGER nvm           ! nombre de vegetations
138      PARAMETER (nvm=8)
139      REAL veget(klon,nvm)  ! couverture vegetale
140      SAVE veget
141
142c
143c
144c Variables argument:
145c
146      INTEGER nlon
147      INTEGER nlev
148      INTEGER nqmax
149      REAL rjourvrai, rjour_ecri
150      REAL gmtime
151      REAL pdtphys
152      LOGICAL debut, lafin
153      REAL paprs(klon,klev+1)
154      REAL pplay(klon,klev)
155      REAL pphi(klon,klev)
156      REAL pphis(klon)
157      REAL paire(klon)
158      REAL presnivs(klev)
159      REAL znivsig(klev)
160      REAL zsurf(nbsrf)
161      real cufi(klon), cvfi(klon)
162
163      REAL u(klon,klev)
164      REAL v(klon,klev)
165      REAL t(klon,klev)
166      REAL qx(klon,klev,nqmax)
167
168      REAL t_ancien(klon,klev), q_ancien(klon,klev)
169      SAVE t_ancien, q_ancien
170      LOGICAL ancien_ok
171      SAVE ancien_ok
172
173      REAL d_t_dyn(klon,klev)
174      REAL d_q_dyn(klon,klev)
175
176      REAL omega(klon,klev)
177
178      REAL d_u(klon,klev)
179      REAL d_v(klon,klev)
180      REAL d_t(klon,klev)
181      REAL d_qx(klon,klev,nqmax)
182      REAL d_ps(klon)
183
184      INTEGER        longcles
185      PARAMETER    ( longcles = 20 )
186      REAL clesphy0( longcles      )
187c
188c Variables quasi-arguments
189c
190      REAL xjour
191      SAVE xjour
192c
193c
194c Variables propres a la physique
195c
196      REAL dtime
197      SAVE dtime                  ! pas temporel de la physique
198c
199      INTEGER radpas
200      SAVE radpas                 ! frequence d'appel rayonnement
201c
202      REAL radsol(klon)
203      SAVE radsol                 ! bilan radiatif au sol
204c
205      REAL rlat(klon)
206      SAVE rlat                   ! latitude pour chaque point
207c
208      REAL rlon(klon)
209      SAVE rlon                   ! longitude pour chaque point
210c
211cc      INTEGER iflag_con
212cc      SAVE iflag_con              ! indicateur de la convection
213c
214      INTEGER itap
215      SAVE itap                   ! compteur pour la physique
216c
217      REAL co2_ppm
218      SAVE co2_ppm                ! concentration du CO2
219c
220      REAL solaire
221      SAVE solaire                ! constante solaire
222c
223      REAL ftsol(klon,nbsrf)
224      SAVE ftsol                  ! temperature du sol
225c
226      REAL ftsoil(klon,nsoilmx,nbsrf)
227      SAVE ftsoil                 ! temperature dans le sol
228c
229      REAL fevap(klon,nbsrf)
230      SAVE fevap                 ! evaporation
231      REAL fluxlat(klon,nbsrf)
232      SAVE fluxlat
233c
234      REAL deltat(klon)
235      SAVE deltat                 ! ecart avec la SST de reference
236c
237      REAL fqsol(klon,nbsrf)
238      SAVE fqsol                  ! humidite du sol
239c
240      REAL fsnow(klon,nbsrf)
241      SAVE fsnow                  ! epaisseur neigeuse
242c
243      REAL falbe(klon,nbsrf)
244      SAVE falbe                  ! albedo par type de surface
245c
246c
247c  Parametres de l'Orographie a l'Echelle Sous-Maille (OESM):
248c
249      REAL zmea(klon)
250      SAVE zmea                   ! orographie moyenne
251c
252      REAL zstd(klon)
253      SAVE zstd                   ! deviation standard de l'OESM
254c
255      REAL zsig(klon)
256      SAVE zsig                   ! pente de l'OESM
257c
258      REAL zgam(klon)
259      save zgam                   ! anisotropie de l'OESM
260c
261      REAL zthe(klon)
262      SAVE zthe                   ! orientation de l'OESM
263c
264      REAL zpic(klon)
265      SAVE zpic                   ! Maximum de l'OESM
266c
267      REAL zval(klon)
268      SAVE zval                   ! Minimum de l'OESM
269c
270      REAL rugoro(klon)
271      SAVE rugoro                 ! longueur de rugosite de l'OESM
272c
273      REAL zulow(klon),zvlow(klon),zustr(klon), zvstr(klon)
274c
275      REAL zuthe(klon),zvthe(klon)
276      SAVE zuthe
277      SAVE zvthe
278      INTEGER igwd,idx(klon),itest(klon)
279c
280      REAL agesno(klon)
281      SAVE agesno                 ! age de la neige
282c
283      REAL alb_neig(klon)
284      SAVE alb_neig               ! albedo de la neige
285cKE43
286c Variables liees a la convection de K. Emanuel (sb):
287c
288      REAL ema_workcbmf(klon)   ! cloud base mass flux
289      SAVE ema_workcbmf
290
291      REAL ema_cbmf(klon)       ! cloud base mass flux
292      SAVE ema_cbmf
293
294      REAL ema_pcb(klon)        ! cloud base pressure
295      SAVE ema_pcb
296
297      REAL ema_pct(klon)        ! cloud top pressure
298      SAVE ema_pct
299
300      REAL bas, top             ! cloud base and top levels
301      SAVE bas
302      SAVE top
303
304      REAL Ma(klon,klev)        ! undilute upward mass flux
305      SAVE Ma
306      REAL ema_work1(klon, klev), ema_work2(klon, klev)
307      SAVE ema_work1, ema_work2
308      REAL wdn(klon), tdn(klon), qdn(klon)
309c Variables locales pour la couche limite (al1):
310c
311cAl1      REAL pblh(klon)           ! Hauteur de couche limite
312cAl1      SAVE pblh
313c34EK
314c
315c Variables locales:
316c
317      REAL cdragh(klon) ! drag coefficient pour T and Q
318      REAL cdragm(klon) ! drag coefficient pour vent
319cAA
320cAA  Pour phytrac
321cAA
322      REAL ycoefh(klon,klev)    ! coef d'echange pour phytrac
323      REAL yu1(klon)            ! vents dans la premiere couche U
324      REAL yv1(klon)            ! vents dans la premiere couche V
325      LOGICAL offline           ! Controle du stockage ds "physique"
326      PARAMETER (offline=.false.)
327      INTEGER physid
328      REAL pfrac_impa(klon,klev)! Produits des coefs lessivage impaction
329      save pfrac_impa
330      REAL pfrac_nucl(klon,klev)! Produits des coefs lessivage nucleation
331      save pfrac_nucl
332      REAL pfrac_1nucl(klon,klev)! Produits des coefs lessi nucl (alpha = 1)
333      save pfrac_1nucl
334      REAL frac_impa(klon,klev) ! fractions d'aerosols lessivees (impaction)
335      REAL frac_nucl(klon,klev) ! idem (nucleation)
336cAA
337      REAL rain_fall(klon) ! pluie
338      REAL snow_fall(klon) ! neige
339      save snow_fall, rain_fall
340      REAL evap(klon), devap(klon) ! evaporation et sa derivee
341      REAL sens(klon), dsens(klon) ! chaleur sensible et sa derivee
342      REAL bils(klon) ! bilan de chaleur au sol
343      REAL fder(klon) ! Derive de flux (sensible et latente)
344      save fder
345      REAL ve(klon) ! integr. verticale du transport meri. de l'energie
346      REAL vq(klon) ! integr. verticale du transport meri. de l'eau
347      REAL ue(klon) ! integr. verticale du transport zonal de l'energie
348      REAL uq(klon) ! integr. verticale du transport zonal de l'eau
349c
350      REAL frugs(klon,nbsrf) ! longueur de rugosite
351      save frugs
352      REAL zxrugs(klon) ! longueur de rugosite
353c
354c Conditions aux limites
355c
356      INTEGER julien
357      INTEGER idayvrai
358      SAVE idayvrai
359c
360      INTEGER lmt_pas
361      SAVE lmt_pas                ! frequence de mise a jour
362      REAL pctsrf(klon,nbsrf)
363      SAVE pctsrf                 ! sous-fraction du sol
364      REAL albsol(klon)
365      SAVE albsol                 ! albedo du sol total
366      REAL wo(klon,klev)
367      SAVE wo                     ! ozone
368c======================================================================
369c
370c Declaration des procedures appelees
371c
372      EXTERNAL angle     ! calculer angle zenithal du soleil
373      EXTERNAL alboc     ! calculer l'albedo sur ocean
374      EXTERNAL albsno    ! calculer albedo sur neige
375      EXTERNAL ajsec     ! ajustement sec
376      EXTERNAL clmain    ! couche limite
377      EXTERNAL condsurf  ! lire les conditions aux limites
378      EXTERNAL conlmd    ! convection (schema LMD)
379cKE43
380      EXTERNAL conema  ! convect4.3
381      EXTERNAL fisrtilp  ! schema de condensation a grande echelle (pluie)
382cAA
383      EXTERNAL fisrtilp_tr ! schema de condensation a grande echelle (pluie)
384c                          ! stockage des coefficients necessaires au
385c                          ! lessivage OFF-LINE et ON-LINE
386cAA
387      EXTERNAL hgardfou  ! verifier les temperatures
388      EXTERNAL nuage     ! calculer les proprietes radiatives
389      EXTERNAL o3cm      ! initialiser l'ozone
390      EXTERNAL orbite    ! calculer l'orbite terrestre
391      EXTERNAL ozonecm   ! prescrire l'ozone
392      EXTERNAL phyetat0  ! lire l'etat initial de la physique
393      EXTERNAL phyredem  ! ecrire l'etat de redemarrage de la physique
394      EXTERNAL radlwsw   ! rayonnements solaire et infrarouge
395      EXTERNAL suphec    ! initialiser certaines constantes
396      EXTERNAL transp    ! transport total de l'eau et de l'energie
397      EXTERNAL ecribina  ! ecrire le fichier binaire global
398      EXTERNAL ecribins  ! ecrire le fichier binaire global
399      EXTERNAL ecrirega  ! ecrire le fichier binaire regional
400      EXTERNAL ecriregs  ! ecrire le fichier binaire regional
401c
402c Variables locales
403c
404      REAL dialiq(klon,klev)  ! eau liquide nuageuse
405      REAL diafra(klon,klev)  ! fraction nuageuse
406      REAL cldliq(klon,klev)  ! eau liquide nuageuse
407      REAL cldfra(klon,klev)  ! fraction nuageuse
408      REAL cldtau(klon,klev)  ! epaisseur optique
409      REAL cldemi(klon,klev)  ! emissivite infrarouge
410c
411C§§§ PB
412      REAL fluxq(klon,klev, nbsrf)   ! flux turbulent d'humidite
413      REAL fluxt(klon,klev, nbsrf)   ! flux turbulent de chaleur
414      REAL fluxu(klon,klev, nbsrf)   ! flux turbulent de vitesse u
415      REAL fluxv(klon,klev, nbsrf)   ! flux turbulent de vitesse v
416c
417      REAL zxfluxt(klon, klev)
418      REAL zxfluxq(klon, klev)
419      REAL zxfluxu(klon, klev)
420      REAL zxfluxv(klon, klev)
421C§§§
422      REAL heat(klon,klev)    ! chauffage solaire
423      REAL heat0(klon,klev)   ! chauffage solaire ciel clair
424      REAL cool(klon,klev)    ! refroidissement infrarouge
425      REAL cool0(klon,klev)   ! refroidissement infrarouge ciel clair
426      REAL topsw(klon), toplw(klon), solsw(klon), sollw(klon)
427      real sollwdown(klon)    ! downward LW flux at surface
428      REAL topsw0(klon), toplw0(klon), solsw0(klon), sollw0(klon)
429      REAL albpla(klon)
430c Le rayonnement n'est pas calcule tous les pas, il faut donc
431c                      sauvegarder les sorties du rayonnement
432      SAVE  heat,cool,albpla,topsw,toplw,solsw,sollw,sollwdown
433      SAVE  topsw0,toplw0,solsw0,sollw0, heat0, cool0
434      INTEGER itaprad
435      SAVE itaprad
436c
437      REAL conv_q(klon,klev) ! convergence de l'humidite (kg/kg/s)
438      REAL conv_t(klon,klev) ! convergence de la temperature(K/s)
439c
440      REAL cldl(klon),cldm(klon),cldh(klon) !nuages bas, moyen et haut
441      REAL cldt(klon),cldq(klon) !nuage total, eau liquide integree
442c
443      REAL zxtsol(klon), zxqsol(klon), zxsnow(klon)
444c
445      REAL dist, rmu0(klon), fract(klon)
446      REAL zdtime, zlongi
447c
448      CHARACTER*2 str2
449c
450      REAL qcheck
451      REAL z_avant(klon), z_apres(klon), z_factor(klon)
452      LOGICAL zx_ajustq
453c
454      REAL za, zb
455      REAL zx_t, zx_qs, zdelta, zcor, zlvdcp, zlsdcp
456      INTEGER i, k, iq, nsrf, ll
457      REAL t_coup
458      PARAMETER (t_coup=234.0)
459c
460      REAL zphi(klon,klev)
461      REAL zx_tmp_x(iim), zx_tmp_yjjmp1
462      REAL zx_relief(iim,jjmp1)
463      REAL zx_aire(iim,jjmp1)
464cKE43
465c Variables locales pour la convection de K. Emanuel (sb):
466c
467      REAL upwd(klon,klev)      ! saturated updraft mass flux
468      REAL dnwd(klon,klev)      ! saturated downdraft mass flux
469      REAL dnwd0(klon,klev)     ! unsaturated downdraft mass flux
470      REAL tvp(klon,klev)       ! virtual temp of lifted parcel
471      REAL cape(klon)           ! CAPE
472      SAVE cape
473      REAL pbase(klon)          ! cloud base pressure
474      SAVE pbase
475      REAL bbase(klon)          ! cloud base buoyancy
476      SAVE bbase
477      REAL rflag(klon)          ! flag fonctionnement de convect
478c -- convect43:
479      INTEGER ntra              ! nb traceurs pour convect4.3
480      REAL pori_con(klon)    ! pressure at the origin level of lifted parcel
481      REAL plcl_con(klon),dtma_con(klon),dtlcl_con(klon)
482      REAL dtvpdt1(klon,klev), dtvpdq1(klon,klev)
483      REAL dplcldt(klon), dplcldr(klon)
484c?     .     condm_con(klon,klev),conda_con(klon,klev),
485c?     .     mr_con(klon,klev),ep_con(klon,klev)
486c?     .    ,sadiab(klon,klev),wadiab(klon,klev)
487c --
488c34EK
489c
490c Variables du changement
491c
492c con: convection
493c lsc: condensation a grande echelle (Large-Scale-Condensation)
494c ajs: ajustement sec
495c eva: evaporation de l'eau liquide nuageuse
496c vdf: couche limite (Vertical DiFfusion)
497      REAL d_t_con(klon,klev),d_q_con(klon,klev)
498      REAL d_u_con(klon,klev),d_v_con(klon,klev)
499      REAL d_t_lsc(klon,klev),d_q_lsc(klon,klev),d_ql_lsc(klon,klev)
500      REAL d_t_ajs(klon,klev), d_q_ajs(klon,klev)
501      REAL d_t_eva(klon,klev),d_q_eva(klon,klev)
502      REAL rneb(klon,klev)
503c
504      REAL pmfu(klon,klev), pmfd(klon,klev)
505      REAL pen_u(klon,klev), pen_d(klon,klev)
506      REAL pde_u(klon,klev), pde_d(klon,klev)
507      INTEGER kcbot(klon), kctop(klon), kdtop(klon)
508      REAL pmflxr(klon,klev+1), pmflxs(klon,klev+1)
509      REAL prfl(klon,klev+1), psfl(klon,klev+1)
510c
511      INTEGER ibas_con(klon), itop_con(klon)
512      REAL rain_con(klon), rain_lsc(klon)
513      REAL snow_con(klon), snow_lsc(klon)
514      REAL d_ts(klon,nbsrf)
515c
516      REAL d_u_vdf(klon,klev), d_v_vdf(klon,klev)
517      REAL d_t_vdf(klon,klev), d_q_vdf(klon,klev)
518c
519      REAL d_u_oro(klon,klev), d_v_oro(klon,klev)
520      REAL d_t_oro(klon,klev)
521      REAL d_u_lif(klon,klev), d_v_lif(klon,klev)
522      REAL d_t_lif(klon,klev)
523
524      REAL ratqs(klon,klev)
525      LOGICAL zpt_conv(klon,klev)
526
527c
528c Variables liees a l'ecriture de la bande histoire physique
529c
530      INTEGER ecrit_mth
531      SAVE ecrit_mth   ! frequence d'ecriture (fichier mensuel)
532c
533      INTEGER ecrit_day
534      SAVE ecrit_day   ! frequence d'ecriture (fichier journalier)
535c
536      INTEGER ecrit_ins
537      SAVE ecrit_ins   ! frequence d'ecriture (fichier instantane)
538c
539      INTEGER ecrit_reg
540      SAVE ecrit_reg   ! frequence d'ecriture
541c
542c
543c
544c Variables locales pour effectuer les appels en serie
545c
546      REAL t_seri(klon,klev), q_seri(klon,klev)
547      REAL ql_seri(klon,klev)
548      REAL u_seri(klon,klev), v_seri(klon,klev)
549c
550      REAL tr_seri(klon,klev,nbtr)
551
552      REAL zx_rh(klon,klev)
553
554      INTEGER        length
555      PARAMETER    ( length = 100 )
556      REAL tabcntr0( length       )
557c
558      INTEGER ndex2d(iim*jjmp1),ndex3d(iim*jjmp1*klev)
559      REAL zx_tmp_fi2d(klon)
560      REAL zx_tmp_2d(iim,jjmp1), zx_tmp_3d(iim,jjmp1,klev)
561      REAL zx_lon(iim,jjmp1), zx_lat(iim,jjmp1)
562c
563      INTEGER nid_day, nid_mth, nid_ins
564      SAVE nid_day, nid_mth, nid_ins
565c
566      INTEGER nhori, nvert
567      REAL zsto, zout, zjulian
568
569      character*20 modname
570      character*80 abort_message
571      logical ok_sync
572      real date0
573
574c
575c Declaration des constantes et des fonctions thermodynamiques
576c
577#include "YOMCST.h"
578#include "YOETHF.h"
579#include "FCTTRE.h"
580c======================================================================
581      modname = 'physiq'
582      ok_sync=.TRUE.
583      IF (nqmax .LT. 2) THEN
584         PRINT*, 'eaux vapeur et liquide sont indispensables'
585         CALL ABORT
586      ENDIF
587      IF (debut) THEN
588         CALL suphec ! initialiser constantes et parametres phys.
589      ENDIF
590
591
592c======================================================================
593      xjour = rjourvrai
594c
595c Si c'est le debut, il faut initialiser plusieurs choses
596c          ********
597c
598       IF (debut) THEN
599
600c
601c appel a la lecture du run.def physique
602c
603         call conf_phys(ocean, ok_veget, ok_journe, ok_mensuel,
604     .                  ok_instan)
605
606         DO k = 2, nvm          ! pas de vegetation
607            DO i = 1, klon
608               veget(i,k) = 0.0
609            ENDDO
610         ENDDO
611         DO i = 1, klon
612            veget(i,1) = 1.0    ! il n'y a que du desert
613         ENDDO
614         PRINT*, 'Pas de vegetation; desert partout'
615c
616c
617c Initialiser les compteurs:
618c
619
620         frugs = 0.
621         itap    = 0
622         itaprad = 0
623c
624         CALL phyetat0 ("startphy.nc",dtime,co2_ppm,solaire,
625     .       rlat,rlon,pctsrf, ftsol,ftsoil,deltat,fqsol,fsnow,
626     .       falbe, fevap, rain_fall,snow_fall,solsw, sollwdown,
627     .       fder,radsol,frugs,agesno,clesphy0,
628     .       zmea,zstd,zsig,zgam,zthe,zpic,zval,rugoro,tabcntr0,
629     .       t_ancien, q_ancien, ancien_ok )
630
631c
632         radpas = NINT( 86400./dtime/nbapp_rad)
633
634c
635         CALL printflag( tabcntr0,radpas,ok_ocean,ok_oasis ,ok_journe,
636     ,                    ok_instan, ok_region )
637c
638         IF (ABS(dtime-pdtphys).GT.0.001) THEN
639            PRINT*, 'Pas physique n est pas correcte',dtime,pdtphys
640            abort_message=' See above '
641            call abort_gcm(modname,abort_message,1)
642         ENDIF
643         IF (nlon .NE. klon) THEN
644            PRINT*, 'nlon et klon ne sont pas coherents', nlon, klon
645            abort_message=' See above '
646            call abort_gcm(modname,abort_message,1)
647         ENDIF
648         IF (nlev .NE. klev) THEN
649            PRINT*, 'nlev et klev ne sont pas coherents', nlev, klev
650            abort_message=' See above '
651            call abort_gcm(modname,abort_message,1)
652         ENDIF
653c
654         IF (dtime*FLOAT(radpas).GT.21600..AND.cycle_diurne) THEN
655           PRINT*, 'Nbre d appels au rayonnement insuffisant'
656           PRINT*, "Au minimum 4 appels par jour si cycle diurne"
657           abort_message=' See above '
658           call abort_gcm(modname,abort_message,1)
659         ENDIF
660         PRINT*, "Clef pour la convection, iflag_con=", iflag_con
661c
662cKE43
663c Initialisation pour la convection de K.E. (sb):
664         IF (iflag_con.EQ.4) THEN
665
666         PRINT*, "*** Convection de Kerry Emanuel 4.3  "
667         PRINT*, "On va utiliser le melange convectif des traceurs qui"
668         PRINT*, "est calcule dans convect4.3"
669         PRINT*, " !!! penser aux logical flags de phytrac"
670
671          DO i = 1, klon
672           ema_cbmf(i) = 0.
673           ema_pcb(i)  = 0.
674           ema_pct(i)  = 0.
675           ema_workcbmf(i) = 0.
676          ENDDO
677         ENDIF
678c34EK
679         IF (ok_orodr) THEN
680         DO i=1,klon
681         rugoro(i) = MAX(1.0e-05, zstd(i)*zsig(i)/2.0)
682         ENDDO
683         CALL SUGWD(klon,klev,paprs,pplay)
684         DO i=1,klon
685         zuthe(i)=0.
686         zvthe(i)=0.
687         if(zstd(i).gt.10.)then
688           zuthe(i)=(1.-zgam(i))*cos(zthe(i))
689           zvthe(i)=(1.-zgam(i))*sin(zthe(i))
690         endif
691         ENDDO
692         ENDIF
693c
694c
695         lmt_pas = NINT(86400./dtime * 1.0)   ! tous les jours
696         PRINT*,'La frequence de lecture surface est de ', lmt_pas
697c
698         ecrit_mth = NINT(86400./dtime *ecritphy)  ! tous les ecritphy jours
699         IF (ok_mensuel) THEN
700         PRINT*, 'La frequence de sortie mensuelle est de ', ecrit_mth
701         ENDIF
702         ecrit_day = NINT(86400./dtime *1.0)  ! tous les jours
703         IF (ok_journe) THEN
704         PRINT*, 'La frequence de sortie journaliere est de ',ecrit_day
705         ENDIF
706ccc         ecrit_ins = NINT(86400./dtime *0.5)  ! 2 fois par jour
707ccc         ecrit_ins = NINT(86400./dtime *0.25)  ! 4 fois par jour
708         ecrit_ins = NINT(86400./dtime/48.)  ! a chaque pas de temps
709         IF (ok_instan) THEN
710         PRINT*, 'La frequence de sortie instant. est de ', ecrit_ins
711         ENDIF
712         ecrit_reg = NINT(86400./dtime *0.25)  ! 4 fois par jour
713         IF (ok_region) THEN
714         PRINT*, 'La frequence de sortie region est de ', ecrit_reg
715         ENDIF
716
717c
718c Initialiser le couplage si necessaire
719c
720      npas = 0
721      nexca = 0
722      if (ocean == 'couple') then
723        npas = itaufin/ iphysiq
724        nexca = 86400 / dtime
725        write(*,*)' ##### Ocean couple #####'
726        write(*,*)' Valeurs des pas de temps'
727        write(*,*)' npas = ', npas
728        write(*,*)' nexca = ', nexca
729      endif       
730c
731c
732      IF (ok_journe) THEN
733c
734         CALL ymds2ju(anne_ini, 1, 1, 0.0, zjulian)
735         zjulian = zjulian + day_ini
736c
737         CALL gr_fi_ecrit(1,klon,iim,jjmp1,rlon,zx_lon)
738         DO i = 1, iim
739            zx_lon(i,1) = rlon(i+1)
740            zx_lon(i,jjmp1) = rlon(i+1)
741         ENDDO
742         DO ll=1,klev
743            znivsig(ll)=float(ll)
744         ENDDO
745         CALL gr_fi_ecrit(1,klon,iim,jjmp1,rlat,zx_lat)
746         CALL histbeg("histday", iim,zx_lon, jjmp1,zx_lat,
747     .                 1,iim,1,jjmp1, 0, zjulian, dtime,
748     .                 nhori, nid_day)
749         CALL histvert(nid_day, "presnivs", "Vertical levels", "mb",
750     .                 klev, presnivs, nvert)
751c        call histvert(nid_day, 'sig_s', 'Niveaux sigma','-',
752c    .              klev, znivsig, nvert)
753c
754         zsto = dtime
755         zout = dtime * FLOAT(ecrit_day)
756c
757         CALL histdef(nid_day, "phis", "Surface geop. height", "-",
758     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
759     .                "once", zsto,zout)
760c
761         CALL histdef(nid_day, "aire", "Grid area", "-",
762     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
763     .                "once", zsto,zout)
764c
765c Champs 2D:
766c
767         CALL histdef(nid_day, "tsol", "Surface Temperature", "K",
768     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
769     .                "ave(X)", zsto,zout)
770c
771         CALL histdef(nid_day, "tter", "Surface Temperature", "K",
772     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
773     .                "ave(X)", zsto,zout)
774c
775         CALL histdef(nid_day, "tlic", "Surface Temperature", "K",
776     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
777     .                "ave(X)", zsto,zout)
778c
779         CALL histdef(nid_day, "toce", "Surface Temperature", "K",
780     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
781     .                "ave(X)", zsto,zout)
782c
783         CALL histdef(nid_day, "tsic", "Surface Temperature", "K",
784     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
785     .                "ave(X)", zsto,zout)
786c
787         CALL histdef(nid_day, "psol", "Surface Pressure", "Pa",
788     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
789     .                "ave(X)", zsto,zout)
790c
791         CALL histdef(nid_day, "rain", "Precipitation", "mm/day",
792     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
793     .                "ave(X)", zsto,zout)
794c
795         CALL histdef(nid_day, "snow", "Snow fall", "mm/day",
796     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
797     .                "ave(X)", zsto,zout)
798c
799         CALL histdef(nid_day, "snow_cov", "Snow cover", "mm",
800     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
801     .                "ave(X)", zsto,zout)
802c
803         CALL histdef(nid_day, "evap", "Evaporation", "mm/day",
804     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
805     .                "ave(X)", zsto,zout)
806c
807         CALL histdef(nid_day, "tops", "Solar rad. at TOA", "W/m2",
808     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
809     .                "ave(X)", zsto,zout)
810c
811         CALL histdef(nid_day, "topl", "IR rad. at TOA", "W/m2",
812     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
813     .                "ave(X)", zsto,zout)
814c
815         CALL histdef(nid_day, "sols", "Solar rad. at surf.", "W/m2",
816     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
817     .                "ave(X)", zsto,zout)
818c
819         CALL histdef(nid_day, "soll", "IR rad. at surface", "W/m2",
820     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
821     .                "ave(X)", zsto,zout)
822c
823         CALL histdef(nid_day, "solldown", "Down. IR rad. at surface",
824     .                "W/m2", iim,jjmp1,nhori, 1,1,1, -99, 32,
825     .                "ave(X)", zsto,zout)
826c
827         CALL histdef(nid_day, "bils", "Surf. total heat flux", "W/m2",
828     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
829     .                "ave(X)", zsto,zout)
830c
831         CALL histdef(nid_day, "sens", "Sensible heat flux", "W/m2",
832     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
833     .                "ave(X)", zsto,zout)
834c
835         CALL histdef(nid_day, "fder", "Heat flux derivation", "W/m2",
836     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
837     .                "ave(X)", zsto,zout)
838c
839         CALL histdef(nid_day, "frtu", "Zonal wind stress", "Pa",
840     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
841     .                "ave(X)", zsto,zout)
842c
843         CALL histdef(nid_day, "frtv", "Meridional wind stress", "Pa",
844     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
845     .                "ave(X)", zsto,zout)
846c
847C §§§ PB flux pour chauqe sous surface
848C
849         DO nsrf = 1, nbsrf
850C
851           call histdef(nid_day, "pourc_"//clnsurf(nsrf),
852     $         "Fraction"//clnsurf(nsrf), "W/m2", 
853     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
854     $         "ave(X)", zsto,zout)
855C
856           call histdef(nid_day, "tsol_"//clnsurf(nsrf),
857     $         "Fraction"//clnsurf(nsrf), "W/m2", 
858     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
859     $         "ave(X)", zsto,zout)
860C
861           call histdef(nid_day, "sens_"//clnsurf(nsrf),
862     $         "Sensible heat flux "//clnsurf(nsrf), "W/m2", 
863     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
864     $         "ave(X)", zsto,zout)
865c
866           call histdef(nid_day, "lat_"//clnsurf(nsrf),
867     $         "Latent heat flux "//clnsurf(nsrf), "W/m2", 
868     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
869     $         "ave(X)", zsto,zout)
870C
871           call histdef(nid_day, "taux_"//clnsurf(nsrf),
872     $         "Zonal wind stress"//clnsurf(nsrf),"Pa",
873     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
874     $         "ave(X)", zsto,zout)
875
876           call histdef(nid_day, "tauy_"//clnsurf(nsrf),
877     $         "Meridional xind stress "//clnsurf(nsrf), "Pa", 
878     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
879     $         "ave(X)", zsto,zout)
880C
881           call histdef(nid_day, "albe_"//clnsurf(nsrf),
882     $         "Albedo surf. "//clnsurf(nsrf), "W/m2", 
883     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
884     $         "ave(X)", zsto,zout)
885C
886           call histdef(nid_day, "rugs_"//clnsurf(nsrf),
887     $         "Latent heat flux "//clnsurf(nsrf), "W/m2", 
888     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
889     $         "ave(X)", zsto,zout)
890
891C§§§
892         END DO
893           
894         CALL histdef(nid_day, "sicf", "Sea-ice fraction", "-",
895     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
896     .                "ave(X)", zsto,zout)
897c
898         CALL histdef(nid_day, "cldl", "Low-level cloudiness", "-",
899     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
900     .                "ave(X)", zsto,zout)
901c
902         CALL histdef(nid_day, "cldm", "Mid-level cloudiness", "-",
903     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
904     .                "ave(X)", zsto,zout)
905c
906         CALL histdef(nid_day, "cldh", "High-level cloudiness", "-",
907     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
908     .                "ave(X)", zsto,zout)
909c
910         CALL histdef(nid_day, "cldt", "Total cloudiness", "-",
911     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
912     .                "ave(X)", zsto,zout)
913c
914         CALL histdef(nid_day, "cldq", "Cloud liquid water path", "-",
915     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
916     .                "ave(X)", zsto,zout)
917c
918c Champs 3D:
919c
920         CALL histdef(nid_day, "temp", "Air temperature", "K",
921     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
922     .                "ave(X)", zsto,zout)
923c
924         CALL histdef(nid_day, "ovap", "Specific humidity", "Kg/Kg",
925     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
926     .                "ave(X)", zsto,zout)
927c
928         CALL histdef(nid_day, "geop", "Geopotential height", "m",
929     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
930     .                "ave(X)", zsto,zout)
931c
932         CALL histdef(nid_day, "vitu", "Zonal wind", "m/s",
933     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
934     .                "ave(X)", zsto,zout)
935c
936         CALL histdef(nid_day, "vitv", "Meridional wind", "m/s",
937     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
938     .                "ave(X)", zsto,zout)
939c
940         CALL histdef(nid_day, "vitw", "Vertical wind", "m/s",
941     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
942     .                "ave(X)", zsto,zout)
943c
944         CALL histdef(nid_day, "pres", "Air pressure", "Pa",
945     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
946     .                "ave(X)", zsto,zout)
947c
948         CALL histend(nid_day)
949c
950         ndex2d = 0
951         ndex3d = 0
952c
953      ENDIF ! fin de test sur ok_journe
954c
955      IF (ok_mensuel) THEN
956c
957         CALL ymds2ju(anne_ini, 1, 1, 0.0, zjulian)
958         zjulian = zjulian + day_ini
959c
960         CALL gr_fi_ecrit(1,klon,iim,jjmp1,rlon,zx_lon)
961         DO i = 1, iim
962            zx_lon(i,1) = rlon(i+1)
963            zx_lon(i,jjmp1) = rlon(i+1)
964         ENDDO
965         DO ll=1,klev
966            znivsig(ll)=float(ll)
967         ENDDO
968         CALL gr_fi_ecrit(1,klon,iim,jjmp1,rlat,zx_lat)
969         CALL histbeg("histmth", iim,zx_lon, jjmp1,zx_lat,
970     .                 1,iim,1,jjmp1, 0, zjulian, dtime,
971     .                 nhori, nid_mth)
972         CALL histvert(nid_mth, "presnivs", "Vertical levels", "mb",
973     .                 klev, presnivs, nvert)
974c        call histvert(nid_mth, 'sig_s', 'Niveaux sigma','-',
975c    .              klev, znivsig, nvert)
976c
977         zsto = dtime
978         zout = dtime * ecrit_mth
979c
980         CALL histdef(nid_mth, "phis", "Surface geop. height", "-",
981     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
982     .                "once",  zsto,zout)
983c
984         CALL histdef(nid_mth, "aire", "Grid area", "-",
985     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
986     .                "once",  zsto,zout)
987c
988c Champs 2D:
989c
990         CALL histdef(nid_mth, "tsol", "Surface Temperature", "K",
991     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
992     .                "ave(X)", zsto,zout)
993c
994         CALL histdef(nid_mth, "psol", "Surface Pressure", "Pa",
995     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
996     .                "ave(X)", zsto,zout)
997c
998         CALL histdef(nid_mth, "qsol", "Surface humidity", "mm",
999     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1000     .                "ave(X)", zsto,zout)
1001c
1002         CALL histdef(nid_mth, "rain", "Precipitation", "mm/day",
1003     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1004     .                "ave(X)", zsto,zout)
1005c
1006         CALL histdef(nid_mth, "plul", "Large-scale Precip.", "mm/day",
1007     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1008     .                "ave(X)", zsto,zout)
1009c
1010         CALL histdef(nid_mth, "pluc", "Convective Precip.", "mm/day",
1011     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1012     .                "ave(X)", zsto,zout)
1013c
1014         CALL histdef(nid_mth, "snow", "Snow fall", "mm/day",
1015     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1016     .                "ave(X)", zsto,zout)
1017c
1018         CALL histdef(nid_mth, "snow_cov", "Snow cover", "mm",
1019     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1020     .                "ave(X)", zsto,zout)
1021c
1022         CALL histdef(nid_mth, "ages", "Snow age", "day",
1023     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1024     .                "ave(X)", zsto,zout)
1025c
1026         CALL histdef(nid_mth, "evap", "Evaporation", "mm/day",
1027     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1028     .                "ave(X)", zsto,zout)
1029c
1030         CALL histdef(nid_mth, "tops", "Solar rad. at TOA", "W/m2",
1031     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1032     .                "ave(X)", zsto,zout)
1033c
1034         CALL histdef(nid_mth, "topl", "IR rad. at TOA", "W/m2",
1035     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1036     .                "ave(X)", zsto,zout)
1037c
1038         CALL histdef(nid_mth, "sols", "Solar rad. at surf.", "W/m2",
1039     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1040     .                "ave(X)", zsto,zout)
1041c
1042         CALL histdef(nid_mth, "soll", "IR rad. at surface", "W/m2",
1043     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1044     .                "ave(X)", zsto,zout)
1045c
1046         CALL histdef(nid_mth, "solldown", "Down. IR rad. at surface",
1047     .                "W/m2", iim,jjmp1,nhori, 1,1,1, -99, 32,
1048     .                "ave(X)", zsto,zout)
1049c
1050         CALL histdef(nid_mth, "tops0", "Solar rad. at TOA", "W/m2",
1051     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1052     .                "ave(X)", zsto,zout)
1053c
1054         CALL histdef(nid_mth, "topl0", "IR rad. at TOA", "W/m2",
1055     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1056     .                "ave(X)", zsto,zout)
1057c
1058         CALL histdef(nid_mth, "sols0", "Solar rad. at surf.", "W/m2",
1059     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1060     .                "ave(X)", zsto,zout)
1061c
1062         CALL histdef(nid_mth, "soll0", "IR rad. at surface", "W/m2",
1063     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1064     .                "ave(X)", zsto,zout)
1065c
1066         CALL histdef(nid_mth, "bils", "Surf. total heat flux", "W/m2",
1067     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1068     .                "ave(X)", zsto,zout)
1069c
1070         CALL histdef(nid_mth, "sens", "Sensible heat flux", "W/m2",
1071     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1072     .                "ave(X)", zsto,zout)
1073c
1074         CALL histdef(nid_mth, "fder", "Heat flux derivation", "W/m2",
1075     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1076     .                "ave(X)", zsto,zout)
1077c
1078         CALL histdef(nid_mth, "frtu", "Zonal wind stress", "Pa",
1079     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1080     .                "ave(X)", zsto,zout)
1081c
1082         CALL histdef(nid_mth, "frtv", "Meridional wind stress", "Pa",
1083     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1084     .                "ave(X)", zsto,zout)
1085c
1086         DO nsrf = 1, nbsrf
1087C
1088           call histdef(nid_mth, "pourc_"//clnsurf(nsrf),
1089     $         "Fraction "//clnsurf(nsrf), "W/m2", 
1090     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
1091     $         "ave(X)", zsto,zout)
1092C
1093           call histdef(nid_mth, "tsol_"//clnsurf(nsrf),
1094     $         "Fraction "//clnsurf(nsrf), "W/m2", 
1095     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
1096     $         "ave(X)", zsto,zout)
1097C
1098           call histdef(nid_mth, "sens_"//clnsurf(nsrf),
1099     $         "Sensible heat flux "//clnsurf(nsrf), "W/m2", 
1100     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
1101     $         "ave(X)", zsto,zout)
1102c
1103           call histdef(nid_mth, "lat_"//clnsurf(nsrf),
1104     $         "Latent heat flux "//clnsurf(nsrf), "W/m2", 
1105     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
1106     $         "ave(X)", zsto,zout)
1107C
1108           call histdef(nid_mth, "taux_"//clnsurf(nsrf),
1109     $         "Zonal wind stress"//clnsurf(nsrf), "Pa", 
1110     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
1111     $         "ave(X)", zsto,zout)
1112
1113           call histdef(nid_mth, "tauy_"//clnsurf(nsrf),
1114     $         "Meridional xind stress "//clnsurf(nsrf), "Pa", 
1115     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
1116     $         "ave(X)", zsto,zout)
1117c
1118           call histdef(nid_mth, "albe_"//clnsurf(nsrf),
1119     $         "Albedo surf. "//clnsurf(nsrf), "W/m2", 
1120     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
1121     $         "ave(X)", zsto,zout)
1122c
1123           call histdef(nid_mth, "rugs_"//clnsurf(nsrf),
1124     $         "Latent heat flux "//clnsurf(nsrf), "W/m2", 
1125     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
1126     $         "ave(X)", zsto,zout)
1127         END DO
1128C
1129         CALL histdef(nid_mth, "sicf", "Sea-ice fraction", "-",
1130     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1131     .                "ave(X)", zsto,zout)
1132c
1133         CALL histdef(nid_mth, "albs", "Surface albedo", "-",
1134     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1135     .                "ave(X)", zsto,zout)
1136c
1137         CALL histdef(nid_mth, "cdrm", "Momentum drag coef.", "-",
1138     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1139     .                "ave(X)", zsto,zout)
1140c
1141         CALL histdef(nid_mth, "cdrh", "Heat drag coef.", "-",
1142     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1143     .                "ave(X)", zsto,zout)
1144c
1145         CALL histdef(nid_mth, "cldl", "Low-level cloudiness", "-",
1146     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1147     .                "ave(X)", zsto,zout)
1148c
1149         CALL histdef(nid_mth, "cldm", "Mid-level cloudiness", "-",
1150     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1151     .                "ave(X)", zsto,zout)
1152c
1153         CALL histdef(nid_mth, "cldh", "High-level cloudiness", "-",
1154     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1155     .                "ave(X)", zsto,zout)
1156c
1157         CALL histdef(nid_mth, "cldt", "Total cloudiness", "-",
1158     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1159     .                "ave(X)", zsto,zout)
1160c
1161         CALL histdef(nid_mth, "cldq", "Cloud liquid water path", "-",
1162     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1163     .                "ave(X)", zsto,zout)
1164c
1165         CALL histdef(nid_mth, "ue", "Zonal energy transport", "-",
1166     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1167     .                "ave(X)", zsto,zout)
1168c
1169         CALL histdef(nid_mth, "ve", "Merid energy transport", "-",
1170     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1171     .                "ave(X)", zsto,zout)
1172c
1173         CALL histdef(nid_mth, "uq", "Zonal humidity transport", "-",
1174     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1175     .                "ave(X)", zsto,zout)
1176c
1177         CALL histdef(nid_mth, "vq", "Merid humidity transport", "-",
1178     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1179     .                "ave(X)", zsto,zout)
1180c
1181cKE43
1182      IF (iflag_con .EQ. 4) THEN ! sb
1183c
1184         CALL histdef(nid_mth, "cape", "Conv avlbl pot ener", "J/Kg",
1185     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1186     .                "ave(X)", zsto,zout)
1187c
1188         CALL histdef(nid_mth, "pbase", "Cld base pressure", "hPa",
1189     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1190     .                "ave(X)", zsto,zout)
1191c
1192         CALL histdef(nid_mth, "ptop", "Cld top pressure", "hPa",
1193     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1194     .                "ave(X)", zsto,zout)
1195c
1196         CALL histdef(nid_mth, "fbase", "Cld base mass flux", "Kg/m2/s",
1197     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1198     .                "ave(X)", zsto,zout)
1199c
1200c
1201      ENDIF
1202c34EK
1203c
1204c Champs 3D:
1205c
1206         CALL histdef(nid_mth, "temp", "Air temperature", "K",
1207     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1208     .                "ave(X)", zsto,zout)
1209c
1210         CALL histdef(nid_mth, "ovap", "Specific humidity", "Kg/Kg",
1211     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1212     .                "ave(X)", zsto,zout)
1213c
1214         CALL histdef(nid_mth, "geop", "Geopotential height", "m",
1215     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1216     .                "ave(X)", zsto,zout)
1217c
1218         CALL histdef(nid_mth, "vitu", "Zonal wind", "m/s",
1219     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1220     .                "ave(X)", zsto,zout)
1221c
1222         CALL histdef(nid_mth, "vitv", "Meridional wind", "m/s",
1223     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1224     .                "ave(X)", zsto,zout)
1225c
1226         CALL histdef(nid_mth, "vitw", "Vertical wind", "m/s",
1227     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1228     .                "ave(X)", zsto,zout)
1229c
1230         CALL histdef(nid_mth, "pres", "Air pressure", "Pa",
1231     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1232     .                "ave(X)", zsto,zout)
1233c
1234         CALL histdef(nid_mth, "rneb", "Cloud fraction", "-",
1235     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1236     .                "ave(X)", zsto,zout)
1237c
1238         CALL histdef(nid_mth, "rhum", "Relative humidity", "-",
1239     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1240     .                "ave(X)", zsto,zout)
1241c
1242         CALL histdef(nid_mth, "oliq", "Liquid water content", "kg/kg",
1243     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1244     .                "ave(X)", zsto,zout)
1245c
1246         CALL histdef(nid_mth, "dtdyn", "Dynamics dT", "K/s",
1247     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1248     .                "ave(X)", zsto,zout)
1249c
1250         CALL histdef(nid_mth, "dqdyn", "Dynamics dQ", "Kg/Kg/s",
1251     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1252     .                "ave(X)", zsto,zout)
1253c
1254         CALL histdef(nid_mth, "dtcon", "Convection dT", "K/s",
1255     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1256     .                "ave(X)", zsto,zout)
1257c
1258         CALL histdef(nid_mth, "dqcon", "Convection dQ", "Kg/Kg/s",
1259     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1260     .                "ave(X)", zsto,zout)
1261c
1262         CALL histdef(nid_mth, "dtlsc", "Condensation dT", "K/s",
1263     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1264     .                "ave(X)", zsto,zout)
1265c
1266         CALL histdef(nid_mth, "dqlsc", "Condensation dQ", "Kg/Kg/s",
1267     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1268     .                "ave(X)", zsto,zout)
1269c
1270         CALL histdef(nid_mth, "dtvdf", "Boundary-layer dT", "K/s",
1271     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1272     .                "ave(X)", zsto,zout)
1273c
1274         CALL histdef(nid_mth, "dqvdf", "Boundary-layer dQ", "Kg/Kg/s",
1275     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1276     .                "ave(X)", zsto,zout)
1277c
1278         CALL histdef(nid_mth, "dteva", "Reevaporation dT", "K/s",
1279     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1280     .                "ave(X)", zsto,zout)
1281c
1282         CALL histdef(nid_mth, "dqeva", "Reevaporation dQ", "Kg/Kg/s",
1283     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1284     .                "ave(X)", zsto,zout)
1285
1286         CALL histdef(nid_mth, "ptconv", "POINTS CONVECTIFS"," ",
1287     .                iim,jjm+1,nhori, klev,1,klev,nvert, 32,
1288     .                "ave(X)", zsto,zout)
1289
1290         CALL histdef(nid_mth, "ratqs", "RATQS"," ",
1291     .                iim,jjm+1,nhori, klev,1,klev,nvert, 32,
1292     .                "ave(X)", zsto,zout)
1293
1294c
1295         CALL histdef(nid_mth, "dtajs", "Dry adjust. dT", "K/s",
1296     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1297     .                "ave(X)", zsto,zout)
1298
1299         CALL histdef(nid_mth, "dqajs", "Dry adjust. dQ", "Kg/Kg/s",
1300     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1301     .                "ave(X)", zsto,zout)
1302c
1303         CALL histdef(nid_mth, "dtswr", "SW radiation dT", "K/s",
1304     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1305     .                "ave(X)", zsto,zout)
1306c
1307         CALL histdef(nid_mth, "dtsw0", "SW radiation dT", "K/s",
1308     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1309     .                "ave(X)", zsto,zout)
1310c
1311         CALL histdef(nid_mth, "dtlwr", "LW radiation dT", "K/s",
1312     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1313     .                "ave(X)", zsto,zout)
1314c
1315         CALL histdef(nid_mth, "dtlw0", "LW radiation dT", "K/s",
1316     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1317     .                "ave(X)", zsto,zout)
1318c
1319         CALL histdef(nid_mth, "duvdf", "Boundary-layer dU", "m/s2",
1320     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1321     .                "ave(X)", zsto,zout)
1322c
1323         CALL histdef(nid_mth, "dvvdf", "Boundary-layer dV", "m/s2",
1324     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1325     .                "ave(X)", zsto,zout)
1326c
1327         IF (ok_orodr) THEN
1328         CALL histdef(nid_mth, "duoro", "Orography dU", "m/s2",
1329     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1330     .                "ave(X)", zsto,zout)
1331c
1332         CALL histdef(nid_mth, "dvoro", "Orography dV", "m/s2",
1333     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1334     .                "ave(X)", zsto,zout)
1335c
1336         ENDIF
1337C
1338         IF (ok_orolf) THEN
1339         CALL histdef(nid_mth, "dulif", "Orography dU", "m/s2",
1340     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1341     .                "ave(X)", zsto,zout)
1342c
1343         CALL histdef(nid_mth, "dvlif", "Orography dV", "m/s2",
1344     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1345     .                "ave(X)", zsto,zout)
1346         ENDIF
1347C
1348         CALL histdef(nid_mth, "ozone", "Ozone concentration", "-",
1349     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1350     .                "ave(X)", zsto,zout)
1351c
1352         if (nqmax.GE.3) THEN
1353         DO iq=1,nqmax-2
1354         IF (iq.LE.99) THEN
1355         WRITE(str2,'(i2.2)') iq
1356         CALL histdef(nid_mth, "trac"//str2, "Tracer No."//str2, "-",
1357     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1358     .                "ave(X)", zsto,zout)
1359         ELSE
1360         PRINT*, "Trop de traceurs"
1361         CALL abort
1362         ENDIF
1363         ENDDO
1364         ENDIF
1365c
1366cKE43
1367      IF (iflag_con.EQ.4) THEN ! (sb)
1368c
1369         CALL histdef(nid_mth, "upwd", "saturated updraft", "Kg/m2/s",
1370     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1371     .                "ave(X)", zsto,zout)
1372c
1373         CALL histdef(nid_mth, "dnwd", "saturated downdraft","Kg/m2/s",
1374     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1375     .                "ave(X)", zsto,zout)
1376c
1377         CALL histdef(nid_mth, "dnwd0", "unsat. downdraft", "Kg/m2/s",
1378     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1379     .                "ave(X)", zsto,zout)
1380c
1381         CALL histdef(nid_mth,"Ma","undilute adiab updraft","Kg/m2/s",
1382     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1383     .                "ave(X)", zsto,zout)
1384c
1385c
1386      ENDIF
1387c34EK
1388         CALL histend(nid_mth)
1389c
1390         ndex2d = 0
1391         ndex3d = 0
1392c
1393      ENDIF ! fin de test sur ok_mensuel
1394c
1395c
1396      IF (ok_instan) THEN
1397c
1398         CALL ymds2ju(anne_ini, 1, 1, 0.0, zjulian)
1399         zjulian = zjulian + day_ini
1400c
1401         CALL gr_fi_ecrit(1,klon,iim,jjmp1,rlon,zx_lon)
1402         DO i = 1, iim
1403            zx_lon(i,1) = rlon(i+1)
1404            zx_lon(i,jjmp1) = rlon(i+1)
1405         ENDDO
1406         DO ll=1,klev
1407            znivsig(ll)=float(ll)
1408         ENDDO
1409         CALL gr_fi_ecrit(1,klon,iim,jjmp1,rlat,zx_lat)
1410         CALL histbeg("histins", iim,zx_lon, jjmp1,zx_lat,
1411     .                 1,iim,1,jjmp1, 0, zjulian, dtime,
1412     .                 nhori, nid_ins)
1413         CALL histvert(nid_ins, "presnivs", "Vertical levels", "mb",
1414     .                 klev, presnivs, nvert)
1415c        call histvert(nid_ins, 'sig_s', 'Niveaux sigma','-',
1416c    .              klev, znivsig, nvert)
1417c
1418c
1419         zsto = dtime * ecrit_ins
1420         zout = dtime * ecrit_ins
1421C
1422         CALL histdef(nid_ins, "phis", "Surface geop. height", "-",
1423     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1424     .                "once", zsto,zout)
1425c
1426         CALL histdef(nid_ins, "aire", "Grid area", "-",
1427     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1428     .                "once", zsto,zout)
1429c
1430c Champs 2D:
1431c
1432        CALL histdef(nid_ins, "tsol", "Surface Temperature", "K",
1433     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1434     .                "inst(X)", zsto,zout)
1435c
1436        CALL histdef(nid_ins, "psol", "Surface Pressure", "Pa",
1437     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1438     .                "inst(X)", zsto,zout)
1439c
1440         CALL histdef(nid_ins, "topl", "OLR", "W/m2",
1441     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1442     .                "inst(X)", zsto,zout)
1443c
1444         CALL histdef(nid_ins, "evap", "Evaporation", "mm/day",
1445     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1446     .                "inst(X)", zsto,zout)
1447c
1448         CALL histdef(nid_ins, "sols", "Solar rad. at surf.", "W/m2",
1449     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1450     .                "inst(X)", zsto,zout)
1451c
1452         CALL histdef(nid_ins, "soll", "IR rad. at surface", "W/m2",
1453     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1454     .                "inst(X)", zsto,zout)
1455c
1456         CALL histdef(nid_ins, "solldown", "Down. IR rad. at surface",
1457     .                "W/m2", iim,jjmp1,nhori, 1,1,1, -99, 32,
1458     .                "ave(X)", zsto,zout)
1459c
1460         CALL histdef(nid_ins, "bils", "Surf. total heat flux", "W/m2",
1461     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1462     .                "inst(X)", zsto,zout)
1463c
1464         CALL histdef(nid_ins, "sens", "Sensible heat flux", "W/m2",
1465     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1466     .                "inst(X)", zsto,zout)
1467c
1468         CALL histdef(nid_ins, "fder", "Heat flux derivation", "W/m2",
1469     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1470     .                "inst(X)", zsto,zout)
1471c
1472      CALL histdef(nid_ins, "dtsvdfo", "Boundary-layer dTs(o)", "K/s",
1473     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1474     .                "inst(X)", zsto,zout)
1475c
1476      CALL histdef(nid_ins, "dtsvdft", "Boundary-layer dTs(t)", "K/s",
1477     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1478     .                "inst(X)", zsto,zout)
1479c
1480      CALL histdef(nid_ins, "dtsvdfg", "Boundary-layer dTs(g)", "K/s",
1481     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1482     .                "inst(X)", zsto,zout)
1483c
1484      CALL histdef(nid_ins, "dtsvdfi", "Boundary-layer dTs(g)", "K/s",
1485     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1486     .                "inst(X)", zsto,zout)
1487
1488         DO nsrf = 1, nbsrf
1489C
1490           call histdef(nid_ins, "pourc_"//clnsurf(nsrf),
1491     $         "Fraction"//clnsurf(nsrf), "W/m2", 
1492     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
1493     $         "inst(X)", zsto,zout)
1494
1495           call histdef(nid_ins, "sens_"//clnsurf(nsrf),
1496     $         "Sensible heat flux "//clnsurf(nsrf), "W/m2", 
1497     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
1498     $         "inst(X)", zsto,zout)
1499c
1500           call histdef(nid_ins, "tsol_"//clnsurf(nsrf),
1501     $         "Surface Temperature"//clnsurf(nsrf), "W/m2", 
1502     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
1503     $         "inst(X)", zsto,zout)
1504c
1505           call histdef(nid_ins, "lat_"//clnsurf(nsrf),
1506     $         "Latent heat flux "//clnsurf(nsrf), "W/m2", 
1507     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
1508     $         "inst(X)", zsto,zout)
1509C
1510           call histdef(nid_ins, "taux_"//clnsurf(nsrf),
1511     $         "Zonal wind stress"//clnsurf(nsrf),"Pa",
1512     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
1513     $         "inst(X)", zsto,zout)
1514
1515           call histdef(nid_ins, "tauy_"//clnsurf(nsrf),
1516     $         "Meridional xind stress "//clnsurf(nsrf), "Pa", 
1517     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
1518     $         "inst(X)", zsto,zout)
1519c
1520           call histdef(nid_ins, "albe_"//clnsurf(nsrf),
1521     $         "Albedo surf. "//clnsurf(nsrf), "-", 
1522     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
1523     $         "inst(X)", zsto,zout)
1524c
1525           call histdef(nid_ins, "rugs_"//clnsurf(nsrf),
1526     $         "rugosite "//clnsurf(nsrf), "-", 
1527     $         iim,jjmp1,nhori, 1,1,1, -99, 32,
1528     $         "inst(X)", zsto,zout)
1529C§§§
1530         END DO
1531         CALL histdef(nid_ins, "rugs", "rugosity", "-",
1532     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1533     .                "inst(X)", zsto,zout)
1534
1535c
1536         CALL histdef(nid_ins, "albs", "Surface albedo", "-",
1537     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1538     .                "inst(X)", zsto,zout)
1539c
1540         CALL histdef(nid_ins, "snow_cov", "Snow cover", "mm",
1541     .                iim,jjmp1,nhori, 1,1,1, -99, 32,
1542     .                "inst(X)", zsto,zout)
1543c
1544c Champs 3D:
1545c
1546         CALL histdef(nid_ins, "temp", "Temperature", "K",
1547     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1548     .                "inst(X)", zsto,zout)
1549c
1550         CALL histdef(nid_ins, "vitu", "Zonal wind", "m/s",
1551     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1552     .                "inst(X)", zsto,zout)
1553c
1554         CALL histdef(nid_ins, "vitv", "Merid wind", "m/s",
1555     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1556     .                "inst(X)", zsto,zout)
1557c
1558         CALL histdef(nid_ins, "geop", "Geopotential height", "m",
1559     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1560     .                "inst(X)", zsto,zout)
1561c
1562         CALL histdef(nid_ins, "pres", "Air pressure", "Pa",
1563     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1564     .                "inst(X)", zsto,zout)
1565c
1566         CALL histdef(nid_ins, "dtvdf", "Boundary-layer dT", "K/s",
1567     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1568     .                "inst(X)", zsto,zout)
1569c
1570         CALL histdef(nid_ins, "dqvdf", "Boundary-layer dQ", "Kg/Kg/s",
1571     .                iim,jjmp1,nhori, klev,1,klev,nvert, 32,
1572     .                "inst(X)", zsto,zout)
1573c
1574
1575         CALL histend(nid_ins)
1576c
1577         ndex2d = 0
1578         ndex3d = 0
1579c
1580      ENDIF
1581c
1582c
1583c
1584c Prescrire l'ozone dans l'atmosphere
1585c
1586c
1587cc         DO i = 1, klon
1588cc         DO k = 1, klev
1589cc            CALL o3cm (paprs(i,k)/100.,paprs(i,k+1)/100., wo(i,k),20)
1590cc         ENDDO
1591cc         ENDDO
1592c
1593c
1594      ENDIF
1595c
1596c   ****************     Fin  de   IF ( debut  )   ***************
1597c
1598c
1599c Mettre a zero des variables de sortie (pour securite)
1600c
1601      DO i = 1, klon
1602         d_ps(i) = 0.0
1603      ENDDO
1604      DO k = 1, klev
1605      DO i = 1, klon
1606         d_t(i,k) = 0.0
1607         d_u(i,k) = 0.0
1608         d_v(i,k) = 0.0
1609      ENDDO
1610      ENDDO
1611      DO iq = 1, nqmax
1612      DO k = 1, klev
1613      DO i = 1, klon
1614         d_qx(i,k,iq) = 0.0
1615      ENDDO
1616      ENDDO
1617      ENDDO
1618c
1619c Ne pas affecter les valeurs entrees de u, v, h, et q
1620c
1621      DO k = 1, klev
1622      DO i = 1, klon
1623         t_seri(i,k)  = t(i,k)
1624         u_seri(i,k)  = u(i,k)
1625         v_seri(i,k)  = v(i,k)
1626         q_seri(i,k)  = qx(i,k,ivap)
1627         ql_seri(i,k) = qx(i,k,iliq)
1628      ENDDO
1629      ENDDO
1630      IF (nqmax.GE.3) THEN
1631      DO iq = 3, nqmax
1632      DO  k = 1, klev
1633      DO  i = 1, klon
1634         tr_seri(i,k,iq-2) = qx(i,k,iq)
1635      ENDDO
1636      ENDDO
1637      ENDDO
1638      ELSE
1639      DO k = 1, klev
1640      DO i = 1, klon
1641         tr_seri(i,k,1) = 0.0
1642      ENDDO
1643      ENDDO
1644      ENDIF
1645c
1646c Diagnostiquer la tendance dynamique
1647c
1648      IF (ancien_ok) THEN
1649         DO k = 1, klev
1650         DO i = 1, klon
1651            d_t_dyn(i,k) = (t_seri(i,k)-t_ancien(i,k))/dtime
1652            d_q_dyn(i,k) = (q_seri(i,k)-q_ancien(i,k))/dtime
1653         ENDDO
1654         ENDDO
1655      ELSE
1656         DO k = 1, klev
1657         DO i = 1, klon
1658            d_t_dyn(i,k) = 0.0
1659            d_q_dyn(i,k) = 0.0
1660         ENDDO
1661         ENDDO
1662         ancien_ok = .TRUE.
1663      ENDIF
1664c
1665c Ajouter le geopotentiel du sol:
1666c
1667      DO k = 1, klev
1668      DO i = 1, klon
1669         zphi(i,k) = pphi(i,k) + pphis(i)
1670      ENDDO
1671      ENDDO
1672c
1673c Verifier les temperatures
1674c
1675      CALL hgardfou(t_seri,ftsol,'debutphy')
1676c
1677c Incrementer le compteur de la physique
1678c
1679      itap   = itap + 1
1680      julien = MOD(NINT(xjour),360)
1681c
1682c Mettre en action les conditions aux limites (albedo, sst, etc.).
1683c Prescrire l'ozone et calculer l'albedo sur l'ocean.
1684c
1685      IF (MOD(itap-1,lmt_pas) .EQ. 0) THEN
1686         idayvrai = NINT(xjour)
1687         PRINT *,' PHYS cond  julien ',julien,idayvrai
1688         CALL ozonecm( FLOAT(julien), rlat, paprs, wo)
1689      ENDIF
1690c
1691c Re-evaporer l'eau liquide nuageuse
1692c
1693      DO k = 1, klev  ! re-evaporation de l'eau liquide nuageuse
1694      DO i = 1, klon
1695         zlvdcp=RLVTT/RCPD/(1.0+RVTMP2*q_seri(i,k))
1696         zlsdcp=RLSTT/RCPD/(1.0+RVTMP2*q_seri(i,k))
1697         zdelta = MAX(0.,SIGN(1.,RTT-t_seri(i,k)))
1698         zb = MAX(0.0,ql_seri(i,k))
1699         za = - MAX(0.0,ql_seri(i,k))
1700     .                  * (zlvdcp*(1.-zdelta)+zlsdcp*zdelta)
1701         t_seri(i,k) = t_seri(i,k) + za
1702         q_seri(i,k) = q_seri(i,k) + zb
1703         ql_seri(i,k) = 0.0
1704         d_t_eva(i,k) = za
1705         d_q_eva(i,k) = zb
1706      ENDDO
1707      ENDDO
1708c
1709c Appeler la diffusion verticale (programme de couche limite)
1710c
1711      DO i = 1, klon
1712c       if (.not. ok_veget) then
1713c          frugs(i,is_ter) = SQRT(frugs(i,is_ter)**2+rugoro(i)**2)
1714c       endif
1715c         frugs(i,is_lic) = rugoro(i)
1716c         frugs(i,is_oce) = rugmer(i)
1717c         frugs(i,is_sic) = 0.001
1718         zxrugs(i) = 0.0
1719      ENDDO
1720      DO nsrf = 1, nbsrf
1721      DO i = 1, klon
1722         frugs(i,nsrf) = MAX(frugs(i,nsrf),0.001)
1723      ENDDO
1724      ENDDO
1725      DO nsrf = 1, nbsrf
1726      DO i = 1, klon
1727            zxrugs(i) = zxrugs(i) + frugs(i,nsrf)*pctsrf(i,nsrf)
1728      ENDDO
1729      ENDDO
1730c
1731C calculs necessaires au calcul de l'albedo dans l'interface
1732c
1733      CALL orbite(FLOAT(julien),zlongi,dist)
1734      IF (cycle_diurne) THEN
1735        zdtime=dtime*FLOAT(radpas) ! pas de temps du rayonnement (s)
1736        CALL zenang(zlongi,gmtime,zdtime,rlat,rlon,rmu0,fract)
1737      ELSE
1738        rmu0 = -999.999
1739      ENDIF
1740
1741      fder = 0.
1742      date0 = day_ini
1743
1744      CALL clmain(dtime,itap,date0,pctsrf,
1745     e            t_seri,q_seri,u_seri,v_seri,
1746     e            julien, rmu0,
1747     e            ok_veget, ocean, npas, nexca, ftsol,
1748     $            soil_model,ftsoil,
1749     $            paprs,pplay,radsol, fsnow,fqsol,fevap,falbe,fluxlat,
1750     e            rain_fall, snow_fall, solsw, sollw, sollwdown, fder,
1751     e            rlon, rlat, cufi, cvfi, frugs,
1752     e            debut, lafin, agesno,rugoro ,
1753     s            d_t_vdf,d_q_vdf,d_u_vdf,d_v_vdf,d_ts,
1754     s            fluxt,fluxq,fluxu,fluxv,cdragh,cdragm,
1755     s            dsens, devap,
1756     s            ycoefh,yu1,yv1)
1757
1758c
1759C§§§ PB
1760C§§§ Incrementation des flux
1761C§§
1762      zxfluxt=0.
1763      zxfluxq=0.
1764      zxfluxu=0.
1765      zxfluxv=0.
1766      DO nsrf = 1, nbsrf
1767        DO k = 1, klev
1768          DO i = 1, klon
1769            zxfluxt(i,k) = zxfluxt(i,k) +
1770     $          fluxt(i,k,nsrf) * pctsrf( i, nsrf)
1771            zxfluxq(i,k) = zxfluxq(i,k) +
1772     $          fluxq(i,k,nsrf) * pctsrf( i, nsrf)
1773            zxfluxu(i,k) = zxfluxu(i,k) +
1774     $          fluxu(i,k,nsrf) * pctsrf( i, nsrf)
1775            zxfluxv(i,k) = zxfluxv(i,k) +
1776     $          fluxv(i,k,nsrf) * pctsrf( i, nsrf)
1777          END DO
1778        END DO
1779      END DO
1780      DO i = 1, klon
1781         sens(i) = - zxfluxt(i,1) ! flux de chaleur sensible au sol
1782c         evap(i) = - fluxq(i,1) ! flux d'evaporation au sol
1783         evap(i) = - zxfluxq(i,1) ! flux d'evaporation au sol
1784C LF test signe flux
1785         sens(i) = zxfluxt(i,1)
1786         evap(i) = zxfluxq(i,1)
1787         fder(i) = dsens(i) + devap(i)
1788      ENDDO
1789
1790      DO k = 1, klev
1791      DO i = 1, klon
1792         t_seri(i,k) = t_seri(i,k) + d_t_vdf(i,k)
1793         q_seri(i,k) = q_seri(i,k) + d_q_vdf(i,k)
1794         u_seri(i,k) = u_seri(i,k) + d_u_vdf(i,k)
1795         v_seri(i,k) = v_seri(i,k) + d_v_vdf(i,k)
1796      ENDDO
1797      ENDDO
1798c
1799c Incrementer la temperature du sol
1800c
1801      DO i = 1, klon
1802         zxtsol(i) = 0.0
1803         IF ( abs( pctsrf(i, is_ter) + pctsrf(i, is_lic) +
1804     $       pctsrf(i, is_oce) + pctsrf(i, is_sic)  - 1.) .GT. EPSFRA)
1805     $       THEN
1806             WRITE(*,*) 'physiq : pb sous surface au point ', i,
1807     $           pctsrf(i, 1 : nbsrf)
1808         ENDIF
1809      ENDDO
1810      DO nsrf = 1, nbsrf
1811      DO i = 1, klon
1812c$$$        IF (pctsrf(i,nsrf) .GE. EPSFRA) THEN
1813            ftsol(i,nsrf) = ftsol(i,nsrf) + d_ts(i,nsrf)
1814            zxtsol(i) = zxtsol(i) + ftsol(i,nsrf)*pctsrf(i,nsrf)
1815c$$$        ENDIF
1816      ENDDO
1817      ENDDO
1818
1819c
1820c Si une sous-fraction n'existe pas, elle prend la temp. moyenne
1821c
1822      DO nsrf = 1, nbsrf
1823        DO i = 1, klon
1824          IF (pctsrf(i,nsrf) .LT. epsfra) ftsol(i,nsrf) = zxtsol(i)
1825        ENDDO
1826      ENDDO
1827
1828c
1829c Calculer la derive du flux infrarouge
1830c
1831      DO nsrf = 1, nbsrf
1832      DO i = 1, klon
1833         fder(i) = fder(i) - 4.0*RSIGMA*zxtsol(i)**3 *
1834     .                       (ftsol(i,nsrf)-zxtsol(i))
1835     .                      *pctsrf(i,nsrf)
1836      ENDDO
1837      ENDDO
1838c
1839c Appeler la convection (au choix)
1840c
1841      DO k = 1, klev
1842      DO i = 1, klon
1843         conv_q(i,k) = d_q_dyn(i,k)
1844     .               + d_q_vdf(i,k)/dtime
1845         conv_t(i,k) = d_t_dyn(i,k)
1846     .               + d_t_vdf(i,k)/dtime
1847      ENDDO
1848      ENDDO
1849      IF (check) THEN
1850         za = qcheck(klon,klev,paprs,q_seri,ql_seri,paire)
1851         PRINT*, "avantcon=", za
1852      ENDIF
1853      zx_ajustq = .FALSE.
1854      IF (iflag_con.EQ.2) zx_ajustq=.TRUE.
1855      IF (zx_ajustq) THEN
1856         DO i = 1, klon
1857            z_avant(i) = 0.0
1858         ENDDO
1859         DO k = 1, klev
1860         DO i = 1, klon
1861            z_avant(i) = z_avant(i) + (q_seri(i,k)+ql_seri(i,k))
1862     .                        *(paprs(i,k)-paprs(i,k+1))/RG
1863         ENDDO
1864         ENDDO
1865      ENDIF
1866      IF (iflag_con.EQ.1) THEN
1867          stop'reactiver le call conlmd dans physiq.F'
1868c     CALL conlmd (dtime, paprs, pplay, t_seri, q_seri, conv_q,
1869c    .             d_t_con, d_q_con,
1870c    .             rain_con, snow_con, ibas_con, itop_con)
1871      ELSE IF (iflag_con.EQ.2) THEN
1872      CALL conflx(dtime, paprs, pplay, t_seri, q_seri,
1873     e            conv_t, conv_q, zxfluxq(1,1), omega,
1874     s            d_t_con, d_q_con, rain_con, snow_con,
1875     s            pmfu, pmfd, pen_u, pde_u, pen_d, pde_d,
1876     s            kcbot, kctop, kdtop, pmflxr, pmflxs)
1877      WHERE (rain_con < 0.) rain_con = 0.
1878      WHERE (snow_con < 0.) snow_con = 0.
1879      DO i = 1, klon
1880         ibas_con(i) = klev+1 - kcbot(i)
1881         itop_con(i) = klev+1 - kctop(i)
1882      ENDDO
1883      ELSE IF (iflag_con.EQ.3) THEN
1884          stop'reactiver le call conlmd dans physiq.F'
1885c     CALL conccm (dtime,paprs,pplay,t_seri,q_seri,conv_q,
1886c    s             d_t_con, d_q_con,
1887c    s             rain_con, snow_con, ibas_con, itop_con)
1888cKE43
1889      ELSE IF (iflag_con.EQ.4) THEN
1890c nb of tracers for the KE convection:
1891          if (nqmax .GE. 4) then
1892              ntra = nbtr
1893          else
1894              ntra = 1
1895          endif
1896cke43 (arguments inutiles enleves => des SAVE dans conema43?)
1897c$$$          CALL conema43(dtime,paprs,pplay,t_seri,q_seri,
1898c$$$     $        u_seri,v_seri,tr_seri,nbtr,
1899c$$$     .        ema_workcbmf,
1900c$$$     .        d_t_con,d_q_con,d_u_con,d_v_con,d_tr,
1901c$$$     .        wdn, tdn, qdn,
1902c$$$     .        rain_con, snow_con, ibas_con, itop_con,
1903c$$$     .        upwd,dnwd,dnwd0,bas,top,Ma,cape,tvp,rflag,
1904c$$$     .        pbase
1905c$$$     .        ,bbase,dtvpdt1,dtvpdq1,dplcldt,dplcldr,
1906c$$$     .        pori_con,plcl_con,dtma_con,dtlcl_con)
1907          CALL conema (dtime,paprs,pplay,t_seri,q_seri,
1908     $        u_seri,v_seri,tr_seri,nbtr,
1909     .        ema_work1,ema_work2,
1910     .        d_t_con,d_q_con,d_u_con,d_v_con,d_tr,
1911c$$$     .        wdn, tdn, qdn,
1912     .        rain_con, snow_con, ibas_con, itop_con,
1913     .        upwd,dnwd,dnwd0,bas,top,Ma,cape,tvp,rflag,
1914     .        pbase
1915     .        ,bbase,dtvpdt1,dtvpdq1,dplcldt,dplcldr)
1916c$$$     .        pori_con,plcl_con,dtma_con,dtlcl_con)
1917          DO i = 1, klon
1918            ema_pcb(i)  = pbase(i)
1919          ENDDO
1920          DO i = 1, klon
1921            ema_pct(i)  = paprs(i,itop_con(i))
1922          ENDDO
1923          DO i = 1, klon
1924            ema_cbmf(i) = ema_workcbmf(i)
1925          ENDDO     
1926      ELSE
1927          PRINT*, "iflag_con non-prevu", iflag_con
1928          CALL abort
1929      ENDIF
1930
1931      CALL homogene(paprs, q_seri, d_q_con, u_seri,v_seri,
1932     .              d_u_con, d_v_con)
1933      DO k = 1, klev
1934        DO i = 1, klon
1935         t_seri(i,k) = t_seri(i,k) + d_t_con(i,k)
1936         q_seri(i,k) = q_seri(i,k) + d_q_con(i,k)
1937         u_seri(i,k) = u_seri(i,k) + d_u_con(i,k)
1938         v_seri(i,k) = v_seri(i,k) + d_v_con(i,k)
1939        ENDDO
1940      ENDDO
1941      IF (check) THEN
1942          za = qcheck(klon,klev,paprs,q_seri,ql_seri,paire)
1943          PRINT*, "aprescon=", za
1944          zx_t = 0.0
1945          za = 0.0
1946          DO i = 1, klon
1947            za = za + paire(i)/FLOAT(klon)
1948            zx_t = zx_t + (rain_con(i)+snow_con(i))*paire(i)/FLOAT(klon)
1949          ENDDO
1950          zx_t = zx_t/za*dtime
1951          PRINT*, "Precip=", zx_t
1952      ENDIF
1953      IF (zx_ajustq) THEN
1954          DO i = 1, klon
1955            z_apres(i) = 0.0
1956          ENDDO
1957          DO k = 1, klev
1958            DO i = 1, klon
1959              z_apres(i) = z_apres(i) + (q_seri(i,k)+ql_seri(i,k))
1960     .            *(paprs(i,k)-paprs(i,k+1))/RG
1961            ENDDO
1962          ENDDO
1963          DO i = 1, klon
1964            z_factor(i) = (z_avant(i)-(rain_con(i)+snow_con(i))*dtime)
1965     .          /z_apres(i)
1966          ENDDO
1967          DO k = 1, klev
1968            DO i = 1, klon
1969              IF (z_factor(i).GT.(1.0+1.0E-08) .OR.
1970     .            z_factor(i).LT.(1.0-1.0E-08)) THEN
1971                  q_seri(i,k) = q_seri(i,k) * z_factor(i)
1972              ENDIF
1973            ENDDO
1974          ENDDO
1975      ENDIF
1976      zx_ajustq=.FALSE.
1977c
1978      IF (nqmax.GT.2) THEN !--melange convectif de traceurs
1979c
1980          IF (iflag_con .NE. 2 .AND.  iflag_con .NE. 4 ) THEN
1981              PRINT*, 'Pour l instant, seul conflx fonctionne ',
1982     $            'avec traceurs', iflag_con
1983              PRINT*,' Mettre iflag_con',
1984     $            ' = 2  ou 4 dans run.def et repasser'
1985              CALL abort
1986              ENDIF
1987c
1988      ENDIF !--nqmax.GT.2
1989c
1990c Appeler l'ajustement sec
1991c
1992      CALL ajsec(paprs, pplay, t_seri, q_seri, d_t_ajs, d_q_ajs)
1993      DO k = 1, klev
1994      DO i = 1, klon
1995         t_seri(i,k) = t_seri(i,k) + d_t_ajs(i,k)
1996         q_seri(i,k) = q_seri(i,k) + d_q_ajs(i,k)
1997      ENDDO
1998      ENDDO
1999
2000c   RATQS
2001      call calcratqs (
2002     I            paprs,pplay,q_seri,d_t_con,d_t_ajs
2003     O           ,ratqs,zpt_conv)
2004c
2005c Appeler le processus de condensation a grande echelle
2006c et le processus de precipitation
2007c
2008      CALL fisrtilp_tr(dtime,paprs,pplay,
2009     .           t_seri, q_seri,ratqs,
2010     .           d_t_lsc, d_q_lsc, d_ql_lsc, rneb, cldliq,
2011     .           rain_lsc, snow_lsc,
2012     .           pfrac_impa, pfrac_nucl, pfrac_1nucl,
2013     .           frac_impa, frac_nucl,
2014     .           prfl, psfl)
2015      WHERE (rain_lsc < 0) rain_lsc = 0.
2016      WHERE (snow_lsc < 0) snow_lsc = 0.
2017      DO k = 1, klev
2018      DO i = 1, klon
2019         t_seri(i,k) = t_seri(i,k) + d_t_lsc(i,k)
2020         q_seri(i,k) = q_seri(i,k) + d_q_lsc(i,k)
2021         ql_seri(i,k) = ql_seri(i,k) + d_ql_lsc(i,k)
2022         cldfra(i,k) = rneb(i,k)
2023         IF (.NOT.new_oliq) cldliq(i,k) = ql_seri(i,k)
2024      ENDDO
2025      ENDDO
2026      IF (check) THEN
2027         za = qcheck(klon,klev,paprs,q_seri,ql_seri,paire)
2028         PRINT*, "apresilp=", za
2029         zx_t = 0.0
2030         za = 0.0
2031         DO i = 1, klon
2032            za = za + paire(i)/FLOAT(klon)
2033            zx_t = zx_t + (rain_lsc(i)+snow_lsc(i))*paire(i)/FLOAT(klon)
2034        ENDDO
2035         zx_t = zx_t/za*dtime
2036         PRINT*, "Precip=", zx_t
2037      ENDIF
2038c
2039c Nuages diagnostiques:
2040c
2041      IF (iflag_con.EQ.2) THEN ! seulement pour Tiedtke
2042      CALL diagcld1(paprs,pplay,
2043     .             rain_con,snow_con,ibas_con,itop_con,
2044     .             diafra,dialiq)
2045      DO k = 1, klev
2046      DO i = 1, klon
2047      IF (diafra(i,k).GT.cldfra(i,k)) THEN
2048         cldliq(i,k) = dialiq(i,k)
2049         cldfra(i,k) = diafra(i,k)
2050      ENDIF
2051      ENDDO
2052      ENDDO
2053      ENDIF
2054c
2055c Nuages stratus artificiels:
2056c
2057      IF (ok_stratus) THEN
2058      CALL diagcld2(paprs,pplay,t_seri,q_seri, diafra,dialiq)
2059      DO k = 1, klev
2060      DO i = 1, klon
2061      IF (diafra(i,k).GT.cldfra(i,k)) THEN
2062         cldliq(i,k) = dialiq(i,k)
2063         cldfra(i,k) = diafra(i,k)
2064      ENDIF
2065      ENDDO
2066      ENDDO
2067      ENDIF
2068c
2069c Precipitation totale
2070c
2071      DO i = 1, klon
2072         rain_fall(i) = rain_con(i) + rain_lsc(i)
2073         snow_fall(i) = snow_con(i) + snow_lsc(i)
2074      ENDDO
2075c
2076c Calculer l'humidite relative pour diagnostique
2077c
2078      DO k = 1, klev
2079      DO i = 1, klon
2080         zx_t = t_seri(i,k)
2081         IF (thermcep) THEN
2082            zdelta = MAX(0.,SIGN(1.,rtt-zx_t))
2083            zx_qs  = r2es * FOEEW(zx_t,zdelta)/pplay(i,k)
2084            zx_qs  = MIN(0.5,zx_qs)
2085            zcor   = 1./(1.-retv*zx_qs)
2086            zx_qs  = zx_qs*zcor
2087         ELSE
2088           IF (zx_t.LT.t_coup) THEN
2089              zx_qs = qsats(zx_t)/pplay(i,k)
2090           ELSE
2091              zx_qs = qsatl(zx_t)/pplay(i,k)
2092           ENDIF
2093         ENDIF
2094         zx_rh(i,k) = q_seri(i,k)/zx_qs
2095      ENDDO
2096      ENDDO
2097c
2098c Calculer les parametres optiques des nuages et quelques
2099c parametres pour diagnostiques:
2100c
2101      CALL nuage (paprs, pplay,
2102     .            t_seri, cldliq, cldfra, cldtau, cldemi,
2103     .            cldh, cldl, cldm, cldt, cldq)
2104c
2105c Appeler le rayonnement mais calculer tout d'abord l'albedo du sol.
2106c
2107      IF (MOD(itaprad,radpas).EQ.0) THEN
2108      DO i = 1, klon
2109         albsol(i) = falbe(i,is_oce) * pctsrf(i,is_oce)
2110     .             + falbe(i,is_lic) * pctsrf(i,is_lic)
2111     .             + falbe(i,is_ter) * pctsrf(i,is_ter)
2112     .             + falbe(i,is_sic) * pctsrf(i,is_sic)
2113      ENDDO
2114      CALL radlwsw ! nouveau rayonnement (compatible Arpege-IFS)
2115     e            (dist, rmu0, fract, co2_ppm, solaire,
2116     e             paprs, pplay,zxtsol,albsol, t_seri,q_seri,wo,
2117     e             cldfra, cldemi, cldtau,
2118     s             heat,heat0,cool,cool0,radsol,albpla,
2119     s             topsw,toplw,solsw,sollw,
2120     s             sollwdown,
2121     s             topsw0,toplw0,solsw0,sollw0)
2122      itaprad = 0
2123      ENDIF
2124      itaprad = itaprad + 1
2125c
2126c Ajouter la tendance des rayonnements (tous les pas)
2127c
2128      DO k = 1, klev
2129      DO i = 1, klon
2130         t_seri(i,k) = t_seri(i,k)
2131     .               + (heat(i,k)-cool(i,k)) * dtime/86400.
2132      ENDDO
2133      ENDDO
2134c
2135c Calculer l'hydrologie de la surface
2136c
2137c      CALL hydrol(dtime,pctsrf,rain_fall, snow_fall, zxevap,
2138c     .            agesno, ftsol,fqsol,fsnow, ruis)
2139c
2140      DO i = 1, klon
2141         zxqsol(i) = 0.0
2142         zxsnow(i) = 0.0
2143      ENDDO
2144      DO nsrf = 1, nbsrf
2145      DO i = 1, klon
2146         zxqsol(i) = zxqsol(i) + fqsol(i,nsrf)*pctsrf(i,nsrf)
2147         zxsnow(i) = zxsnow(i) + fsnow(i,nsrf)*pctsrf(i,nsrf)
2148      ENDDO
2149      ENDDO
2150c
2151c Si une sous-fraction n'existe pas, elle prend la valeur moyenne
2152c
2153c$$$      DO nsrf = 1, nbsrf
2154c$$$      DO i = 1, klon
2155c$$$         IF (pctsrf(i,nsrf).LT.epsfra) THEN
2156c$$$            fqsol(i,nsrf) = zxqsol(i)
2157c$$$            fsnow(i,nsrf) = zxsnow(i)
2158c$$$         ENDIF
2159c$$$      ENDDO
2160c$$$      ENDDO
2161c
2162c Calculer le bilan du sol et la derive de temperature (couplage)
2163c
2164      DO i = 1, klon
2165         bils(i) = radsol(i) - sens(i) - evap(i)*RLVTT
2166      ENDDO
2167c
2168cmoddeblott(jan95)
2169c Appeler le programme de parametrisation de l'orographie
2170c a l'echelle sous-maille:
2171c
2172      IF (ok_orodr) THEN
2173c
2174c  selection des points pour lesquels le shema est actif:
2175        igwd=0
2176        DO i=1,klon
2177        itest(i)=0
2178c        IF ((zstd(i).gt.10.0)) THEN
2179        IF (((zpic(i)-zmea(i)).GT.100.).AND.(zstd(i).GT.10.0)) THEN
2180          itest(i)=1
2181          igwd=igwd+1
2182          idx(igwd)=i
2183        ENDIF
2184        ENDDO
2185c        igwdim=MAX(1,igwd)
2186c
2187        CALL drag_noro(klon,klev,dtime,paprs,pplay,
2188     e                   zmea,zstd, zsig, zgam, zthe,zpic,zval,
2189     e                   igwd,idx,itest,
2190     e                   t_seri, u_seri, v_seri,
2191     s                   zulow, zvlow, zustr, zvstr,
2192     s                   d_t_oro, d_u_oro, d_v_oro)
2193c
2194c  ajout des tendances
2195        DO k = 1, klev
2196        DO i = 1, klon
2197           t_seri(i,k) = t_seri(i,k) + d_t_oro(i,k)
2198           u_seri(i,k) = u_seri(i,k) + d_u_oro(i,k)
2199           v_seri(i,k) = v_seri(i,k) + d_v_oro(i,k)
2200        ENDDO
2201        ENDDO
2202c
2203      ENDIF ! fin de test sur ok_orodr
2204c
2205      IF (ok_orolf) THEN
2206c
2207c  selection des points pour lesquels le shema est actif:
2208        igwd=0
2209        DO i=1,klon
2210        itest(i)=0
2211        IF ((zpic(i)-zmea(i)).GT.100.) THEN
2212          itest(i)=1
2213          igwd=igwd+1
2214          idx(igwd)=i
2215        ENDIF
2216        ENDDO
2217c        igwdim=MAX(1,igwd)
2218c
2219        CALL lift_noro(klon,klev,dtime,paprs,pplay,
2220     e                   rlat,zmea,zstd,zpic,
2221     e                   itest,
2222     e                   t_seri, u_seri, v_seri,
2223     s                   zulow, zvlow, zustr, zvstr,
2224     s                   d_t_lif, d_u_lif, d_v_lif)
2225c
2226c  ajout des tendances
2227        DO k = 1, klev
2228        DO i = 1, klon
2229           t_seri(i,k) = t_seri(i,k) + d_t_lif(i,k)
2230           u_seri(i,k) = u_seri(i,k) + d_u_lif(i,k)
2231           v_seri(i,k) = v_seri(i,k) + d_v_lif(i,k)
2232        ENDDO
2233        ENDDO
2234c
2235      ENDIF ! fin de test sur ok_orolf
2236c
2237cAA
2238cAA Installation de l'interface online-offline pour traceurs
2239cAA
2240c====================================================================
2241c   Calcul  des tendances traceurs
2242c====================================================================
2243C Pascale : il faut quand meme apeller phytrac car il gere les sorties
2244cKE43       des traceurs => il faut donc mettre des flags a .false.
2245      IF (iflag_con.EQ.4) THEN
2246c           on ajoute les tendances calculees par KE43
2247        DO iq=1, nqmax-2 ! Sandrine a -3 ???
2248        DO k = 1, nlev
2249        DO i = 1, klon
2250          tr_seri(i,k,iq) = tr_seri(i,k,iq) + d_tr(i,k,iq)
2251        ENDDO
2252        ENDDO
2253        WRITE(iqn,'(i2.2)') iq
2254        CALL minmaxqfi(tr_seri(1,1,iq),0.,1.e33,'couche lim iq='//iqn)
2255        ENDDO
2256CMAF modif pour garder info du nombre de traceurs auxquels
2257C la physique s'applique
2258      ELSE
2259CMAF modif pour garder info du nombre de traceurs auxquels
2260C la physique s'applique
2261C
2262      call phytrac (rnpb,
2263     I                   debut,lafin,
2264     I                   nqmax-2,
2265     I                   nlon,nlev,dtime,
2266     I                   t,paprs,pplay,
2267     I                   pmfu, pmfd, pen_u, pde_u, pen_d, pde_d,
2268     I                   ycoefh,yu1,yv1,ftsol,pctsrf,rlat,
2269     I                   frac_impa, frac_nucl,
2270     I                   rlon,presnivs,paire,pphis,
2271     O                   tr_seri)
2272      ENDIF
2273
2274      IF (offline) THEN
2275
2276         call phystokenc (
2277     I                   nlon,nlev,pdtphys,rlon,rlat,
2278     I                   t,pmfu, pmfd, pen_u, pde_u, pen_d, pde_d,
2279     I                   ycoefh,yu1,yv1,ftsol,pctsrf,
2280     I                   frac_impa, frac_nucl,
2281     I                   pphis,paire,dtime,itap)
2282
2283
2284      ENDIF
2285
2286c
2287c Calculer le transport de l'eau et de l'energie (diagnostique)
2288c
2289      CALL transp (paprs,zxtsol,
2290     e                   t_seri, q_seri, u_seri, v_seri, zphi,
2291     s                   ve, vq, ue, uq)
2292c
2293c Accumuler les variables a stocker dans les fichiers histoire:
2294c
2295c
2296c
2297
2298      IF (ok_journe) THEN
2299c
2300      ndex2d = 0
2301      ndex3d = 0
2302c
2303c Champs 2D:
2304c
2305         zsto = dtime
2306         zout = dtime * FLOAT(ecrit_day)
2307
2308         i = NINT(zout/zsto)
2309         CALL gr_fi_ecrit(1,klon,iim,jjmp1,pphis,zx_tmp_2d)
2310         CALL histwrite(nid_day,"phis",i,zx_tmp_2d,iim*jjmp1,ndex2d)
2311c
2312         i = NINT(zout/zsto)
2313         CALL gr_fi_ecrit(1,klon,iim,jjmp1,paire,zx_tmp_2d)
2314         CALL histwrite(nid_day,"aire",i,zx_tmp_2d,iim*jjmp1,ndex2d)
2315C
2316      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zxtsol,zx_tmp_2d)
2317      CALL histwrite(nid_day,"tsol",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2318c
2319C
2320      zx_tmp_fi2d(1 : klon) = ftsol(1 : klon, is_ter)
2321      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d ,zx_tmp_2d)
2322      CALL histwrite(nid_day,"tter",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2323C
2324      zx_tmp_fi2d(1 : klon) = ftsol(1 : klon, is_lic)
2325      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zxtsol,zx_tmp_2d)
2326      CALL histwrite(nid_day,"tlic",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2327C
2328      zx_tmp_fi2d(1 : klon) = ftsol(1 : klon, is_oce)
2329      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zxtsol,zx_tmp_2d)
2330      CALL histwrite(nid_day,"toce",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2331C
2332      zx_tmp_fi2d(1 : klon) = ftsol(1 : klon, is_sic)
2333      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zxtsol,zx_tmp_2d)
2334      CALL histwrite(nid_day,"tsic",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2335C
2336      DO i = 1, klon
2337         zx_tmp_fi2d(i) = paprs(i,1)
2338      ENDDO
2339      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d,zx_tmp_2d)
2340      CALL histwrite(nid_day,"psol",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2341c
2342      DO i = 1, klon
2343         zx_tmp_fi2d(i) = rain_fall(i) + snow_fall(i)
2344      ENDDO
2345      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d,zx_tmp_2d)
2346      CALL histwrite(nid_day,"rain",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2347c
2348      CALL gr_fi_ecrit(1, klon,iim,jjmp1, snow_fall,zx_tmp_2d)
2349      CALL histwrite(nid_day,"snow",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2350c
2351      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zxsnow,zx_tmp_2d)
2352      CALL histwrite(nid_day,"snow_cov",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2353c
2354      CALL gr_fi_ecrit(1, klon,iim,jjmp1, evap,zx_tmp_2d)
2355      CALL histwrite(nid_day,"evap",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2356c
2357      CALL gr_fi_ecrit(1, klon,iim,jjmp1, topsw,zx_tmp_2d)
2358      CALL histwrite(nid_day,"tops",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2359c
2360      CALL gr_fi_ecrit(1, klon,iim,jjmp1, toplw,zx_tmp_2d)
2361      CALL histwrite(nid_day,"topl",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2362c
2363      CALL gr_fi_ecrit(1, klon,iim,jjmp1, solsw,zx_tmp_2d)
2364      CALL histwrite(nid_day,"sols",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2365c
2366      CALL gr_fi_ecrit(1, klon,iim,jjmp1, sollw,zx_tmp_2d)
2367      CALL histwrite(nid_day,"soll",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2368c
2369      CALL gr_fi_ecrit(1, klon,iim,jjmp1, sollwdown,zx_tmp_2d)
2370      CALL histwrite(nid_day,"solldown",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2371c
2372      CALL gr_fi_ecrit(1, klon,iim,jjmp1, bils,zx_tmp_2d)
2373      CALL histwrite(nid_day,"bils",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2374c
2375      CALL gr_fi_ecrit(1, klon,iim,jjmp1, sens,zx_tmp_2d)
2376      CALL histwrite(nid_day,"sens",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2377c
2378      CALL gr_fi_ecrit(1, klon,iim,jjmp1, fder,zx_tmp_2d)
2379      CALL histwrite(nid_day,"fder",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2380c
2381c
2382      DO nsrf = 1, nbsrf
2383C§§§
2384        zx_tmp_fi2d(1 : klon) = pctsrf( 1 : klon, nsrf)
2385        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2386        CALL histwrite(nid_day,"pourc_"//clnsurf(nsrf),itap,
2387     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2388C
2389        zx_tmp_fi2d(1 : klon) = ftsol( 1 : klon, nsrf)
2390        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2391        CALL histwrite(nid_day,"tsol_"//clnsurf(nsrf),itap,
2392     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2393C
2394        zx_tmp_fi2d(1 : klon) = fluxt( 1 : klon, 1, nsrf)
2395        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2396        CALL histwrite(nid_day,"sens_"//clnsurf(nsrf),itap,
2397     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2398C
2399        zx_tmp_fi2d(1 : klon) = fluxlat( 1 : klon, nsrf)
2400        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2401        CALL histwrite(nid_day,"lat_"//clnsurf(nsrf),itap,
2402     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2403C
2404        zx_tmp_fi2d(1 : klon) = fluxu( 1 : klon, 1, nsrf)
2405        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2406        CALL histwrite(nid_day,"taux_"//clnsurf(nsrf),itap,
2407     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2408C     
2409        zx_tmp_fi2d(1 : klon) = fluxv( 1 : klon, 1, nsrf)
2410        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2411        CALL histwrite(nid_day,"tauy_"//clnsurf(nsrf),itap,
2412     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2413C
2414        zx_tmp_fi2d(1 : klon) = falbe( 1 : klon, nsrf)
2415        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2416        CALL histwrite(nid_day,"albe_"//clnsurf(nsrf),itap,
2417     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2418C
2419        zx_tmp_fi2d(1 : klon) = frugs( 1 : klon, nsrf)
2420        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2421        CALL histwrite(nid_day,"rugs_"//clnsurf(nsrf),itap,
2422     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2423C
2424      END DO 
2425C
2426c$$$      DO i = 1, klon
2427c$$$         zx_tmp_fi2d(i) = pctsrf(i,is_sic)
2428c$$$      ENDDO
2429c$$$      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d,zx_tmp_2d)
2430c$$$      CALL histwrite(nid_day,"sicf",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2431c
2432      CALL gr_fi_ecrit(1, klon,iim,jjmp1, cldl,zx_tmp_2d)
2433      CALL histwrite(nid_day,"cldl",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2434c
2435      CALL gr_fi_ecrit(1, klon,iim,jjmp1, cldm,zx_tmp_2d)
2436      CALL histwrite(nid_day,"cldm",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2437c
2438      CALL gr_fi_ecrit(1, klon,iim,jjmp1, cldh,zx_tmp_2d)
2439      CALL histwrite(nid_day,"cldh",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2440c
2441      CALL gr_fi_ecrit(1, klon,iim,jjmp1, cldt,zx_tmp_2d)
2442      CALL histwrite(nid_day,"cldt",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2443c
2444      CALL gr_fi_ecrit(1, klon,iim,jjmp1, cldq,zx_tmp_2d)
2445      CALL histwrite(nid_day,"cldq",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2446c
2447c Champs 3D:
2448c
2449      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, t_seri, zx_tmp_3d)
2450      CALL histwrite(nid_day,"temp",itap,zx_tmp_3d,
2451     .                                   iim*jjmp1*klev,ndex3d)
2452c
2453      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, qx(1,1,ivap), zx_tmp_3d)
2454      CALL histwrite(nid_day,"ovap",itap,zx_tmp_3d,
2455     .                                   iim*jjmp1*klev,ndex3d)
2456c
2457      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, zphi, zx_tmp_3d)
2458      CALL histwrite(nid_day,"geop",itap,zx_tmp_3d,
2459     .                                   iim*jjmp1*klev,ndex3d)
2460c
2461      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, u_seri, zx_tmp_3d)
2462      CALL histwrite(nid_day,"vitu",itap,zx_tmp_3d,
2463     .                                   iim*jjmp1*klev,ndex3d)
2464c
2465      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, v_seri, zx_tmp_3d)
2466      CALL histwrite(nid_day,"vitv",itap,zx_tmp_3d,
2467     .                                   iim*jjmp1*klev,ndex3d)
2468c
2469      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, omega, zx_tmp_3d)
2470      CALL histwrite(nid_day,"vitw",itap,zx_tmp_3d,
2471     .                                   iim*jjmp1*klev,ndex3d)
2472c
2473      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, pplay, zx_tmp_3d)
2474      CALL histwrite(nid_day,"pres",itap,zx_tmp_3d,
2475     .                                   iim*jjmp1*klev,ndex3d)
2476c
2477      if (ok_sync) then
2478        call histsync(nid_day)
2479      endif
2480      ENDIF
2481C
2482      IF (ok_mensuel) THEN
2483c
2484      ndex2d = 0
2485      ndex3d = 0
2486c
2487c Champs 2D:
2488c
2489         zsto = dtime
2490         zout = dtime * ecrit_mth
2491
2492         i = NINT(zout/zsto)
2493         CALL gr_fi_ecrit(1,klon,iim,jjmp1,pphis,zx_tmp_2d)
2494         CALL histwrite(nid_mth,"phis",i,zx_tmp_2d,iim*jjmp1,ndex2d)
2495C
2496         i = NINT(zout/zsto)
2497         CALL gr_fi_ecrit(1,klon,iim,jjmp1,paire,zx_tmp_2d)
2498         CALL histwrite(nid_mth,"aire",i,zx_tmp_2d,iim*jjmp1,ndex2d)
2499
2500      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zxtsol,zx_tmp_2d)
2501      CALL histwrite(nid_mth,"tsol",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2502c
2503      DO i = 1, klon
2504         zx_tmp_fi2d(i) = paprs(i,1)
2505      ENDDO
2506      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d,zx_tmp_2d)
2507      CALL histwrite(nid_mth,"psol",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2508c
2509      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zxqsol,zx_tmp_2d)
2510      CALL histwrite(nid_mth,"qsol",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2511c
2512      DO i = 1, klon
2513         zx_tmp_fi2d(i) = rain_fall(i) + snow_fall(i)
2514      ENDDO
2515      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d,zx_tmp_2d)
2516      CALL histwrite(nid_mth,"rain",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2517c
2518      DO i = 1, klon
2519         zx_tmp_fi2d(i) = rain_lsc(i) + snow_lsc(i)
2520      ENDDO
2521      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d,zx_tmp_2d)
2522      CALL histwrite(nid_mth,"plul",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2523c
2524      DO i = 1, klon
2525         zx_tmp_fi2d(i) = rain_con(i) + snow_con(i)
2526      ENDDO
2527      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d,zx_tmp_2d)
2528      CALL histwrite(nid_mth,"pluc",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2529c
2530      CALL gr_fi_ecrit(1, klon,iim,jjmp1, snow_fall,zx_tmp_2d)
2531      CALL histwrite(nid_mth,"snow",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2532c
2533      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zxsnow,zx_tmp_2d)
2534      CALL histwrite(nid_mth,"snow_cov",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2535c
2536      CALL gr_fi_ecrit(1, klon,iim,jjmp1, agesno,zx_tmp_2d)
2537      CALL histwrite(nid_mth,"ages",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2538c
2539      CALL gr_fi_ecrit(1, klon,iim,jjmp1, evap,zx_tmp_2d)
2540      CALL histwrite(nid_mth,"evap",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2541c
2542      CALL gr_fi_ecrit(1, klon,iim,jjmp1, topsw,zx_tmp_2d)
2543      CALL histwrite(nid_mth,"tops",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2544c
2545      CALL gr_fi_ecrit(1, klon,iim,jjmp1, toplw,zx_tmp_2d)
2546      CALL histwrite(nid_mth,"topl",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2547c
2548      CALL gr_fi_ecrit(1, klon,iim,jjmp1, solsw,zx_tmp_2d)
2549      CALL histwrite(nid_mth,"sols",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2550c
2551      CALL gr_fi_ecrit(1, klon,iim,jjmp1, sollw,zx_tmp_2d)
2552      CALL histwrite(nid_mth,"soll",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2553c
2554      CALL gr_fi_ecrit(1, klon,iim,jjmp1, sollwdown,zx_tmp_2d)
2555      CALL histwrite(nid_mth,"solldown",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2556c
2557      CALL gr_fi_ecrit(1, klon,iim,jjmp1, topsw0,zx_tmp_2d)
2558      CALL histwrite(nid_mth,"tops0",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2559c
2560      CALL gr_fi_ecrit(1, klon,iim,jjmp1, toplw0,zx_tmp_2d)
2561      CALL histwrite(nid_mth,"topl0",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2562c
2563      CALL gr_fi_ecrit(1, klon,iim,jjmp1, solsw0,zx_tmp_2d)
2564      CALL histwrite(nid_mth,"sols0",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2565c
2566      CALL gr_fi_ecrit(1, klon,iim,jjmp1, sollw0,zx_tmp_2d)
2567      CALL histwrite(nid_mth,"soll0",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2568c
2569      CALL gr_fi_ecrit(1, klon,iim,jjmp1, bils,zx_tmp_2d)
2570      CALL histwrite(nid_mth,"bils",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2571c
2572      CALL gr_fi_ecrit(1, klon,iim,jjmp1, sens,zx_tmp_2d)
2573      CALL histwrite(nid_mth,"sens",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2574c
2575      CALL gr_fi_ecrit(1, klon,iim,jjmp1, fder,zx_tmp_2d)
2576      CALL histwrite(nid_mth,"fder",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2577c
2578c
2579c      DO i = 1, klon
2580c         zx_tmp_fi2d(i) = fluxu(i,1)
2581c      ENDDO
2582c      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d,zx_tmp_2d)
2583c      CALL histwrite(nid_mth,"frtu",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2584c
2585c      DO i = 1, klon
2586c         zx_tmp_fi2d(i) = fluxv(i,1)
2587c      ENDDO
2588c      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d,zx_tmp_2d)
2589c      CALL histwrite(nid_mth,"frtv",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2590c
2591      DO nsrf = 1, nbsrf
2592C§§§
2593        zx_tmp_fi2d(1 : klon) = pctsrf( 1 : klon, nsrf)
2594        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2595        CALL histwrite(nid_mth,"pourc_"//clnsurf(nsrf),itap,
2596     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2597C
2598        zx_tmp_fi2d(1 : klon) = ftsol( 1 : klon, nsrf)
2599        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2600        CALL histwrite(nid_mth,"tsol_"//clnsurf(nsrf),itap,
2601     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2602C
2603        zx_tmp_fi2d(1 : klon) = fluxt( 1 : klon, 1, nsrf)
2604        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2605        CALL histwrite(nid_mth,"sens_"//clnsurf(nsrf),itap,
2606     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2607C
2608        zx_tmp_fi2d(1 : klon) = fluxlat( 1 : klon, nsrf)
2609        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2610        CALL histwrite(nid_mth,"lat_"//clnsurf(nsrf),itap,
2611     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2612C
2613        zx_tmp_fi2d(1 : klon) = fluxu( 1 : klon, 1, nsrf)
2614        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2615        CALL histwrite(nid_mth,"taux_"//clnsurf(nsrf),itap,
2616     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2617C     
2618        zx_tmp_fi2d(1 : klon) = fluxv( 1 : klon, 1, nsrf)
2619        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2620        CALL histwrite(nid_mth,"tauy_"//clnsurf(nsrf),itap,
2621     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2622C
2623        zx_tmp_fi2d(1 : klon) = falbe( 1 : klon, nsrf)
2624        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2625        CALL histwrite(nid_mth,"albe_"//clnsurf(nsrf),itap,
2626     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2627C
2628        zx_tmp_fi2d(1 : klon) = frugs( 1 : klon, nsrf)
2629        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2630        CALL histwrite(nid_mth,"rugs_"//clnsurf(nsrf),itap,
2631     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2632
2633      END DO 
2634c$$$      DO i = 1, klon
2635c$$$         zx_tmp_fi2d(i) = pctsrf(i,is_sic)
2636c$$$      ENDDO
2637c$$$      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d,zx_tmp_2d)
2638c$$$      CALL histwrite(nid_mth,"sicf",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2639c
2640      CALL gr_fi_ecrit(1, klon,iim,jjmp1, albsol,zx_tmp_2d)
2641      CALL histwrite(nid_mth,"albs",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2642c
2643      CALL gr_fi_ecrit(1, klon,iim,jjmp1, cdragm,zx_tmp_2d)
2644      CALL histwrite(nid_mth,"cdrm",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2645c
2646      CALL gr_fi_ecrit(1, klon,iim,jjmp1, cdragh,zx_tmp_2d)
2647      CALL histwrite(nid_mth,"cdrh",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2648c
2649      CALL gr_fi_ecrit(1, klon,iim,jjmp1, cldl,zx_tmp_2d)
2650      CALL histwrite(nid_mth,"cldl",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2651c
2652      CALL gr_fi_ecrit(1, klon,iim,jjmp1, cldm,zx_tmp_2d)
2653      CALL histwrite(nid_mth,"cldm",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2654c
2655      CALL gr_fi_ecrit(1, klon,iim,jjmp1, cldh,zx_tmp_2d)
2656      CALL histwrite(nid_mth,"cldh",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2657c
2658      CALL gr_fi_ecrit(1, klon,iim,jjmp1, cldt,zx_tmp_2d)
2659      CALL histwrite(nid_mth,"cldt",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2660c
2661      CALL gr_fi_ecrit(1, klon,iim,jjmp1, cldq,zx_tmp_2d)
2662      CALL histwrite(nid_mth,"cldq",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2663c
2664      CALL gr_fi_ecrit(1, klon,iim,jjmp1, ue,zx_tmp_2d)
2665      CALL histwrite(nid_mth,"ue",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2666c
2667      CALL gr_fi_ecrit(1, klon,iim,jjmp1, ve,zx_tmp_2d)
2668      CALL histwrite(nid_mth,"ve",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2669c
2670      CALL gr_fi_ecrit(1, klon,iim,jjmp1, uq,zx_tmp_2d)
2671      CALL histwrite(nid_mth,"uq",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2672c
2673      CALL gr_fi_ecrit(1, klon,iim,jjmp1, vq,zx_tmp_2d)
2674      CALL histwrite(nid_mth,"vq",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2675cKE43
2676      IF (iflag_con .EQ. 4) THEN ! sb
2677c
2678      CALL gr_fi_ecrit(1, klon,iim,jjmp1, cape,zx_tmp_2d)
2679      CALL histwrite(nid_mth,"cape",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2680c
2681      CALL gr_fi_ecrit(1, klon,iim,jjmp1,pbase,zx_tmp_2d)
2682      CALL histwrite(nid_mth,"pbase",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2683c
2684      CALL gr_fi_ecrit(1, klon,iim,jjmp1,ema_pct,zx_tmp_2d)
2685      CALL histwrite(nid_mth,"ptop",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2686c
2687      CALL gr_fi_ecrit(1, klon,iim,jjmp1,ema_cbmf,zx_tmp_2d)
2688      CALL histwrite(nid_mth,"fbase",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2689c
2690c
2691      ENDIF
2692c34EK
2693c
2694c Champs 3D:
2695C
2696      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, t_seri, zx_tmp_3d)
2697      CALL histwrite(nid_mth,"temp",itap,zx_tmp_3d,
2698     .                                   iim*jjmp1*klev,ndex3d)
2699c
2700      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, qx(1,1,ivap), zx_tmp_3d)
2701      CALL histwrite(nid_mth,"ovap",itap,zx_tmp_3d,
2702     .                                   iim*jjmp1*klev,ndex3d)
2703c
2704      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, zphi, zx_tmp_3d)
2705      CALL histwrite(nid_mth,"geop",itap,zx_tmp_3d,
2706     .                                   iim*jjmp1*klev,ndex3d)
2707c
2708      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, u_seri, zx_tmp_3d)
2709      CALL histwrite(nid_mth,"vitu",itap,zx_tmp_3d,
2710     .                                   iim*jjmp1*klev,ndex3d)
2711c
2712      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, v_seri, zx_tmp_3d)
2713      CALL histwrite(nid_mth,"vitv",itap,zx_tmp_3d,
2714     .                                   iim*jjmp1*klev,ndex3d)
2715c
2716      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, omega, zx_tmp_3d)
2717      CALL histwrite(nid_mth,"vitw",itap,zx_tmp_3d,
2718     .                                   iim*jjmp1*klev,ndex3d)
2719c
2720      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, pplay, zx_tmp_3d)
2721      CALL histwrite(nid_mth,"pres",itap,zx_tmp_3d,
2722     .                                   iim*jjmp1*klev,ndex3d)
2723c
2724      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, cldfra, zx_tmp_3d)
2725      CALL histwrite(nid_mth,"rneb",itap,zx_tmp_3d,
2726     .                                   iim*jjmp1*klev,ndex3d)
2727c
2728      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, zx_rh, zx_tmp_3d)
2729      CALL histwrite(nid_mth,"rhum",itap,zx_tmp_3d,
2730     .                                   iim*jjmp1*klev,ndex3d)
2731c
2732      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, cldliq, zx_tmp_3d)
2733      CALL histwrite(nid_mth,"oliq",itap,zx_tmp_3d,
2734     .                                   iim*jjmp1*klev,ndex3d)
2735c
2736      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_t_dyn, zx_tmp_3d)
2737      CALL histwrite(nid_mth,"dtdyn",itap,zx_tmp_3d,
2738     .                                   iim*jjmp1*klev,ndex3d)
2739c
2740      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_q_dyn, zx_tmp_3d)
2741      CALL histwrite(nid_mth,"dqdyn",itap,zx_tmp_3d,
2742     .                                   iim*jjmp1*klev,ndex3d)
2743c
2744      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_t_con, zx_tmp_3d)
2745      CALL histwrite(nid_mth,"dtcon",itap,zx_tmp_3d,
2746     .                                   iim*jjmp1*klev,ndex3d)
2747c
2748      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_q_con, zx_tmp_3d)
2749      CALL histwrite(nid_mth,"dqcon",itap,zx_tmp_3d,
2750     .                                   iim*jjmp1*klev,ndex3d)
2751c
2752      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_t_lsc, zx_tmp_3d)
2753      CALL histwrite(nid_mth,"dtlsc",itap,zx_tmp_3d,
2754     .                                   iim*jjmp1*klev,ndex3d)
2755c
2756      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_q_lsc, zx_tmp_3d)
2757      CALL histwrite(nid_mth,"dqlsc",itap,zx_tmp_3d,
2758     .                                   iim*jjmp1*klev,ndex3d)
2759c
2760      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_t_vdf, zx_tmp_3d)
2761      CALL histwrite(nid_mth,"dtvdf",itap,zx_tmp_3d,
2762     .                                   iim*jjmp1*klev,ndex3d)
2763c
2764      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_q_vdf, zx_tmp_3d)
2765      CALL histwrite(nid_mth,"dqvdf",itap,zx_tmp_3d,
2766     .                                   iim*jjmp1*klev,ndex3d)
2767c
2768      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_t_eva, zx_tmp_3d)
2769      CALL histwrite(nid_mth,"dteva",itap,zx_tmp_3d,
2770     .                                   iim*jjmp1*klev,ndex3d)
2771c
2772      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_q_eva, zx_tmp_3d)
2773      CALL histwrite(nid_mth,"dqeva",itap,zx_tmp_3d,
2774     .                                   iim*jjmp1*klev,ndex3d)
2775c
2776      CALL gr_fi_ecrit(klev,klon,iim,jjm+1, zpt_conv, zx_tmp_3d)
2777      CALL histwrite(nid_mth,"ptconv",itap,zx_tmp_3d,
2778     .                                   iim*(jjm+1)*klev,ndex3d)
2779c
2780      CALL gr_fi_ecrit(klev,klon,iim,jjm+1, ratqs, zx_tmp_3d)
2781      CALL histwrite(nid_mth,"ratqs",itap,zx_tmp_3d,
2782     .                                   iim*(jjm+1)*klev,ndex3d)
2783c
2784      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_t_ajs, zx_tmp_3d)
2785      CALL histwrite(nid_mth,"dtajs",itap,zx_tmp_3d,
2786     .                                   iim*jjmp1*klev,ndex3d)
2787c
2788      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_q_ajs, zx_tmp_3d)
2789      CALL histwrite(nid_mth,"dqajs",itap,zx_tmp_3d,
2790     .                                   iim*jjmp1*klev,ndex3d)
2791c
2792      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, heat, zx_tmp_3d)
2793      CALL histwrite(nid_mth,"dtswr",itap,zx_tmp_3d,
2794     .                                   iim*jjmp1*klev,ndex3d)
2795c
2796      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, heat0, zx_tmp_3d)
2797      CALL histwrite(nid_mth,"dtsw0",itap,zx_tmp_3d,
2798     .                                   iim*jjmp1*klev,ndex3d)
2799c
2800      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, cool, zx_tmp_3d)
2801      CALL histwrite(nid_mth,"dtlwr",itap,zx_tmp_3d,
2802     .                                   iim*jjmp1*klev,ndex3d)
2803c
2804      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, cool0, zx_tmp_3d)
2805      CALL histwrite(nid_mth,"dtlw0",itap,zx_tmp_3d,
2806     .                                   iim*jjmp1*klev,ndex3d)
2807c
2808      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_u_vdf, zx_tmp_3d)
2809      CALL histwrite(nid_mth,"duvdf",itap,zx_tmp_3d,
2810     .                                   iim*jjmp1*klev,ndex3d)
2811c
2812      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_v_vdf, zx_tmp_3d)
2813      CALL histwrite(nid_mth,"dvvdf",itap,zx_tmp_3d,
2814     .                                   iim*jjmp1*klev,ndex3d)
2815c
2816      IF (ok_orodr) THEN
2817      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_u_oro, zx_tmp_3d)
2818      CALL histwrite(nid_mth,"duoro",itap,zx_tmp_3d,
2819     .                                   iim*jjmp1*klev,ndex3d)
2820c
2821      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_v_oro, zx_tmp_3d)
2822      CALL histwrite(nid_mth,"dvoro",itap,zx_tmp_3d,
2823     .                                   iim*jjmp1*klev,ndex3d)
2824c
2825      ENDIF
2826C
2827      IF (ok_orolf) THEN
2828      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_u_lif, zx_tmp_3d)
2829      CALL histwrite(nid_mth,"dulif",itap,zx_tmp_3d,
2830     .                                   iim*jjmp1*klev,ndex3d)
2831c
2832      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_v_lif, zx_tmp_3d)
2833      CALL histwrite(nid_mth,"dvlif",itap,zx_tmp_3d,
2834     .                                   iim*jjmp1*klev,ndex3d)
2835      ENDIF
2836C
2837      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, wo, zx_tmp_3d)
2838      CALL histwrite(nid_mth,"ozone",itap,zx_tmp_3d,
2839     .                                   iim*jjmp1*klev,ndex3d)
2840c
2841      IF (nqmax.GE.3) THEN
2842      DO iq=1,nqmax-2
2843      IF (iq.LE.99) THEN
2844         CALL gr_fi_ecrit(klev,klon,iim,jjmp1, qx(1,1,iq+2), zx_tmp_3d)
2845         WRITE(str2,'(i2.2)') iq
2846         CALL histwrite(nid_mth,"trac"//str2,itap,zx_tmp_3d,
2847     .                                   iim*jjmp1*klev,ndex3d)
2848      ELSE
2849         PRINT*, "Trop de traceurs"
2850         CALL abort
2851      ENDIF
2852      ENDDO
2853      ENDIF
2854cKE43
2855      IF (iflag_con.EQ.4) THEN ! (sb)
2856c
2857      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, upwd, zx_tmp_3d)
2858      CALL histwrite(nid_mth,"upwd",itap,zx_tmp_3d,
2859     .                                   iim*jjmp1*klev,ndex3d)
2860c
2861      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, dnwd, zx_tmp_3d)
2862      CALL histwrite(nid_mth,"dnwd",itap,zx_tmp_3d,
2863     .                                   iim*jjmp1*klev,ndex3d)
2864c
2865      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, dnwd0, zx_tmp_3d)
2866      CALL histwrite(nid_mth,"dnwd0",itap,zx_tmp_3d,
2867     .                                   iim*jjmp1*klev,ndex3d)
2868c
2869      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, Ma, zx_tmp_3d)
2870      CALL histwrite(nid_mth,"Ma",itap,zx_tmp_3d,
2871     .                                   iim*jjmp1*klev,ndex3d)
2872c
2873c
2874      ENDIF
2875c34EK
2876c
2877      if (ok_sync) then
2878        call histsync(nid_mth)
2879      endif
2880      ENDIF
2881c
2882      IF (ok_instan) THEN
2883c
2884      ndex2d = 0
2885      ndex3d = 0
2886c
2887c Champs 2D:
2888c
2889         zsto = dtime * ecrit_ins
2890         zout = dtime * ecrit_ins
2891
2892         i = NINT(zout/zsto)
2893         CALL gr_fi_ecrit(1,klon,iim,jjmp1,pphis,zx_tmp_2d)
2894         CALL histwrite(nid_ins,"phis",i,zx_tmp_2d,iim*jjmp1,ndex2d)
2895c
2896         i = NINT(zout/zsto)
2897         CALL gr_fi_ecrit(1,klon,iim,jjmp1,paire,zx_tmp_2d)
2898         CALL histwrite(nid_ins,"aire",i,zx_tmp_2d,iim*jjmp1,ndex2d)
2899
2900      DO i = 1, klon
2901         zx_tmp_fi2d(i) = paprs(i,1)
2902      ENDDO
2903      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d,zx_tmp_2d)
2904      CALL histwrite(nid_ins,"psol",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2905c
2906      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zxtsol,zx_tmp_2d)
2907      CALL histwrite(nid_ins,"tsol",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2908c
2909      CALL gr_fi_ecrit(1, klon,iim,jjmp1, toplw,zx_tmp_2d)
2910      CALL histwrite(nid_ins,"topl",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2911c
2912      CALL gr_fi_ecrit(1, klon,iim,jjmp1, evap,zx_tmp_2d)
2913      CALL histwrite(nid_ins,"evap",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2914c
2915      CALL gr_fi_ecrit(1, klon,iim,jjmp1, solsw,zx_tmp_2d)
2916      CALL histwrite(nid_ins,"sols",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2917c
2918      CALL gr_fi_ecrit(1, klon,iim,jjmp1, sollw,zx_tmp_2d)
2919      CALL histwrite(nid_ins,"soll",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2920c
2921      CALL gr_fi_ecrit(1, klon,iim,jjmp1, sollwdown,zx_tmp_2d)
2922      CALL histwrite(nid_ins,"solldown",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2923c
2924      CALL gr_fi_ecrit(1, klon,iim,jjmp1, bils,zx_tmp_2d)
2925      CALL histwrite(nid_ins,"bils",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2926c
2927      CALL gr_fi_ecrit(1, klon,iim,jjmp1, sens,zx_tmp_2d)
2928      CALL histwrite(nid_ins,"sens",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2929c
2930      CALL gr_fi_ecrit(1, klon,iim,jjmp1, fder,zx_tmp_2d)
2931      CALL histwrite(nid_ins,"fder",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2932c
2933      CALL gr_fi_ecrit(1, klon,iim,jjmp1, d_ts(1,is_oce),zx_tmp_2d)
2934      CALL histwrite(nid_ins,"dtsvdfo",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2935c
2936      CALL gr_fi_ecrit(1, klon,iim,jjmp1, d_ts(1,is_ter),zx_tmp_2d)
2937      CALL histwrite(nid_ins,"dtsvdft",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2938c
2939      CALL gr_fi_ecrit(1, klon,iim,jjmp1, d_ts(1,is_lic),zx_tmp_2d)
2940      CALL histwrite(nid_ins,"dtsvdfg",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2941c
2942      CALL gr_fi_ecrit(1, klon,iim,jjmp1, d_ts(1,is_sic),zx_tmp_2d)
2943      CALL histwrite(nid_ins,"dtsvdfi",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2944
2945      DO nsrf = 1, nbsrf
2946C§§§
2947        zx_tmp_fi2d(1 : klon) = pctsrf( 1 : klon, nsrf)
2948        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2949        CALL histwrite(nid_ins,"pourc_"//clnsurf(nsrf),itap,
2950     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2951C
2952        zx_tmp_fi2d(1 : klon) = fluxt( 1 : klon, 1, nsrf)
2953        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2954        CALL histwrite(nid_ins,"sens_"//clnsurf(nsrf),itap,
2955     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2956C
2957        zx_tmp_fi2d(1 : klon) = fluxlat( 1 : klon, nsrf)
2958        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2959        CALL histwrite(nid_ins,"lat_"//clnsurf(nsrf),itap,
2960     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2961C
2962        zx_tmp_fi2d(1 : klon) = ftsol( 1 : klon, nsrf)
2963        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2964        CALL histwrite(nid_ins,"tsol_"//clnsurf(nsrf),itap,
2965     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2966C
2967        zx_tmp_fi2d(1 : klon) = fluxu( 1 : klon, 1, nsrf)
2968        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2969        CALL histwrite(nid_ins,"taux_"//clnsurf(nsrf),itap,
2970     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2971C     
2972        zx_tmp_fi2d(1 : klon) = fluxv( 1 : klon, 1, nsrf)
2973        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2974        CALL histwrite(nid_ins,"tauy_"//clnsurf(nsrf),itap,
2975     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2976C
2977        zx_tmp_fi2d(1 : klon) = frugs( 1 : klon, nsrf)
2978        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2979        CALL histwrite(nid_ins,"rugs_"//clnsurf(nsrf),itap,
2980     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2981C
2982        zx_tmp_fi2d(1 : klon) = falbe( 1 : klon, nsrf)
2983        CALL gr_fi_ecrit(1, klon,iim,jjmp1, zx_tmp_fi2d , zx_tmp_2d)
2984        CALL histwrite(nid_ins,"albe_"//clnsurf(nsrf),itap,
2985     $      zx_tmp_2d,iim*jjmp1,ndex2d)
2986C
2987      END DO 
2988      CALL gr_fi_ecrit(1, klon,iim,jjmp1, albsol,zx_tmp_2d)
2989      CALL histwrite(nid_ins,"albs",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2990c
2991      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zxsnow,zx_tmp_2d)
2992      CALL histwrite(nid_ins,"snow_cov",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2993c
2994      CALL gr_fi_ecrit(1, klon,iim,jjmp1, zxrugs,zx_tmp_2d)
2995      CALL histwrite(nid_ins,"rugs",itap,zx_tmp_2d,iim*jjmp1,ndex2d)
2996c
2997c Champs 3D:
2998c
2999      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, t_seri, zx_tmp_3d)
3000      CALL histwrite(nid_ins,"temp",itap,zx_tmp_3d,
3001     .                                   iim*jjmp1*klev,ndex3d)
3002c
3003      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, u_seri, zx_tmp_3d)
3004      CALL histwrite(nid_ins,"vitu",itap,zx_tmp_3d,
3005     .                                   iim*jjmp1*klev,ndex3d)
3006c
3007      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, v_seri, zx_tmp_3d)
3008      CALL histwrite(nid_ins,"vitv",itap,zx_tmp_3d,
3009     .                                   iim*jjmp1*klev,ndex3d)
3010c
3011      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, zphi, zx_tmp_3d)
3012      CALL histwrite(nid_ins,"geop",itap,zx_tmp_3d,
3013     .                                   iim*jjmp1*klev,ndex3d)
3014c
3015      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, pplay, zx_tmp_3d)
3016      CALL histwrite(nid_ins,"pres",itap,zx_tmp_3d,
3017     .                                   iim*jjmp1*klev,ndex3d)
3018c
3019      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_t_vdf, zx_tmp_3d)
3020      CALL histwrite(nid_ins,"dtvdf",itap,zx_tmp_3d,
3021     .                                   iim*jjmp1*klev,ndex3d)
3022c
3023      CALL gr_fi_ecrit(klev,klon,iim,jjmp1, d_q_vdf, zx_tmp_3d)
3024      CALL histwrite(nid_ins,"dqvdf",itap,zx_tmp_3d,
3025     .                                   iim*jjmp1*klev,ndex3d)
3026
3027c
3028      if (ok_sync) then
3029        call histsync(nid_ins)
3030      endif
3031      ENDIF
3032c
3033c
3034c Ecrire la bande regionale (binaire grads)
3035      IF (ok_region .AND. mod(itap,ecrit_reg).eq.0) THEN
3036         CALL ecriregs(84,zxtsol)
3037         CALL ecriregs(84,paprs(1,1))
3038         CALL ecriregs(84,topsw)
3039         CALL ecriregs(84,toplw)
3040         CALL ecriregs(84,solsw)
3041         CALL ecriregs(84,sollw)
3042         CALL ecriregs(84,rain_fall)
3043         CALL ecriregs(84,snow_fall)
3044         CALL ecriregs(84,evap)
3045         CALL ecriregs(84,sens)
3046         CALL ecriregs(84,bils)
3047         CALL ecriregs(84,pctsrf(1,is_sic))
3048         CALL ecriregs(84,zxfluxu(1,1))
3049         CALL ecriregs(84,zxfluxv(1,1))
3050         CALL ecriregs(84,ue)
3051         CALL ecriregs(84,ve)
3052         CALL ecriregs(84,uq)
3053         CALL ecriregs(84,vq)
3054c
3055         CALL ecrirega(84,u_seri)
3056         CALL ecrirega(84,v_seri)
3057         CALL ecrirega(84,omega)
3058         CALL ecrirega(84,t_seri)
3059         CALL ecrirega(84,zphi)
3060         CALL ecrirega(84,q_seri)
3061         CALL ecrirega(84,cldfra)
3062         CALL ecrirega(84,cldliq)
3063         CALL ecrirega(84,pplay)
3064
3065
3066cc         CALL ecrirega(84,d_t_dyn)
3067cc         CALL ecrirega(84,d_q_dyn)
3068cc         CALL ecrirega(84,heat)
3069cc         CALL ecrirega(84,cool)
3070cc         CALL ecrirega(84,d_t_con)
3071cc         CALL ecrirega(84,d_q_con)
3072cc         CALL ecrirega(84,d_t_lsc)
3073cc         CALL ecrirega(84,d_q_lsc)
3074      ENDIF
3075c
3076c Convertir les incrementations en tendances
3077c
3078      DO k = 1, klev
3079      DO i = 1, klon
3080         d_u(i,k) = ( u_seri(i,k) - u(i,k) ) / dtime
3081         d_v(i,k) = ( v_seri(i,k) - v(i,k) ) / dtime
3082         d_t(i,k) = ( t_seri(i,k)-t(i,k) ) / dtime
3083         d_qx(i,k,ivap) = ( q_seri(i,k) - qx(i,k,ivap) ) / dtime
3084         d_qx(i,k,iliq) = ( ql_seri(i,k) - qx(i,k,iliq) ) / dtime
3085      ENDDO
3086      ENDDO
3087c
3088      IF (nqmax.GE.3) THEN
3089      DO iq = 3, nqmax
3090      DO  k = 1, klev
3091      DO  i = 1, klon
3092         d_qx(i,k,iq) = ( tr_seri(i,k,iq-2) - qx(i,k,iq) ) / dtime
3093      ENDDO
3094      ENDDO
3095      ENDDO
3096      ENDIF
3097c
3098c Sauvegarder les valeurs de t et q a la fin de la physique:
3099c
3100      DO k = 1, klev
3101      DO i = 1, klon
3102         t_ancien(i,k) = t_seri(i,k)
3103         q_ancien(i,k) = q_seri(i,k)
3104      ENDDO
3105      ENDDO
3106c
3107c====================================================================
3108c Si c'est la fin, il faut conserver l'etat de redemarrage
3109c====================================================================
3110c
3111      IF (lafin) THEN
3112ccc         IF (ok_oasis) CALL quitcpl
3113         CALL phyredem ("restartphy.nc",dtime,radpas,co2_ppm,solaire,
3114     .      rlat, rlon, pctsrf, ftsol, ftsoil, deltat, fqsol, fsnow,
3115     .      falbe, fevap, rain_fall, snow_fall,
3116     .      solsw, sollwdown,fder,
3117     .      radsol,frugs,agesno,
3118     .      zmea,zstd,zsig,zgam,zthe,zpic,zval,rugoro,
3119     .      t_ancien, q_ancien)
3120      ENDIF
3121
3122      RETURN
3123      END
3124      FUNCTION qcheck(klon,klev,paprs,q,ql,aire)
3125      IMPLICIT none
3126c
3127c Calculer et imprimer l'eau totale. A utiliser pour verifier
3128c la conservation de l'eau
3129c
3130#include "YOMCST.h"
3131      INTEGER klon,klev
3132      REAL paprs(klon,klev+1), q(klon,klev), ql(klon,klev)
3133      REAL aire(klon)
3134      REAL qtotal, zx, qcheck
3135      INTEGER i, k
3136c
3137      zx = 0.0
3138      DO i = 1, klon
3139         zx = zx + aire(i)
3140      ENDDO
3141      qtotal = 0.0
3142      DO k = 1, klev
3143      DO i = 1, klon
3144         qtotal = qtotal + (q(i,k)+ql(i,k)) * aire(i)
3145     .                     *(paprs(i,k)-paprs(i,k+1))/RG
3146      ENDDO
3147      ENDDO
3148c
3149      qcheck = qtotal/zx
3150c
3151      RETURN
3152      END
3153      SUBROUTINE gr_fi_ecrit(nfield,nlon,iim,jjmp1,fi,ecrit)
3154      IMPLICIT none
3155c
3156c Tranformer une variable de la grille physique a
3157c la grille d'ecriture
3158c
3159      INTEGER nfield,nlon,iim,jjmp1, jjm
3160      REAL fi(nlon,nfield), ecrit(iim*jjmp1,nfield)
3161c
3162      INTEGER i, n, ig
3163c
3164      jjm = jjmp1 - 1
3165      DO n = 1, nfield
3166         DO i=1,iim
3167            ecrit(i,n) = fi(1,n)
3168            ecrit(i+jjm*iim,n) = fi(nlon,n)
3169         ENDDO
3170         DO ig = 1, nlon - 2
3171           ecrit(iim+ig,n) = fi(1+ig,n)
3172         ENDDO
3173      ENDDO
3174      RETURN
3175      END
3176
Note: See TracBrowser for help on using the repository browser.