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

Last change on this file since 177 was 177, checked in by lmdzadmin, 24 years ago

Lots of stuff, plus particulierement:

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