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

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

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