*=*=*=*= physiq.html =*=*=*=*
SUBROUTINE physiq

SUBROUTINE physiq


      SUBROUTINE physiq(ecri,ngrid,nlayer,nq,
     $            firstcall,lastcall,
     $            pday,ptime,ptimestep,
     $            pplev,pplay,pphi,
     $            pu,pv,pt,pq,
     $            pdudyn,pdvdyn,pdtdyn,pdqdyn,
     $            pw,
     $            pdu,pdv,pdt,pdq,pdpsrf)

      implicit none

#define   undim

c=======================================================================
c
c   subject:
c   --------
c
c   Organisation of the physical parametrisations of the LMD
c   martian atmospheric general circulation modele.
c   It includes:
c   raditive transfer (long and shortwave) for CO2 and dust.
c   vertical turbulent mixing
c   convective adjsutment
c   condensation and sublimation of carbon dioxide.
c
c   physical part history file writing
c
c
c   author: Frederic Hourdin	15/10/93
c   -------
c           Francois Forget		1994
c
c           Christophe Hourdin	02/1997
c               Cutting of physiq call in order to reduce the size
c               of the program (NDLO2 and ndomainsz in dimradmars.h)
c               Only for Radite and Gravity wave drag
c
c
c   arguments:
c   ----------
c
c   input:
c   ------
c
c    ngrid                 Size of the horizontal grid.
c                          All internal loops are performed on that grid.
c    nlayer                Number of vertical layers.
c    nq                    Number of advected fields
c    firstcall             True at the first call
c    lastcall              True at the last call
c    pday                  Number of days counted from the North. Spring
c                          equinoxe.
c    ptime                 hour (s)
c    ptimestep             timestep (s)
c    pplay(ngrid,nlayer+1) Pressure at the middle of the layers (Pa)
c    pplev(ngrid,nlayer+1) intermediate pressure levels (pa)
c    pphi(ngrid,nlayer)    Geopotential at the middle of the layers (m2s-2)
c    pu(ngrid,nlayer)      u component of the wind (ms-1)
c    pv(ngrid,nlayer)      v component of the wind (ms-1)
c    pt(ngrid,nlayer)      Temperature (K)
c    pq(ngrid,nlayer,nq)   Advected fields
c    pudyn(ngrid,nlayer)    \
c    pvdyn(ngrid,nlayer)     \ Dynamical temporal derivative for the
c    ptdyn(ngrid,nlayer)     / corresponding variables
c    pqdyn(ngrid,nlayer,nq) /
c    pw(ngrid,?)           vertical velocity
c
c   output:
c   -------
c
c    pdu(ngrid,nlayermx)        \
c    pdv(ngrid,nlayermx)         \  Temporal derivative of the corresponding
c    pdt(ngrid,nlayermx)         /  variables due to physical processes.
c    pdq(ngrid,nlayermx)        /
c    pdpsrf(ngrid)             /
c
c=======================================================================
c
c    0.  Declarations :
c    ------------------

#include "dimensions.h"
#include "dimphys.h"
#include "comgeomfi.h"
#include "surfdat.h"
#include "comdiurn.h"
#include "callkeys.h"
#include "comcstfi.h"
#include "planete.h"
#include "drsdef.h"
#include "comsaison.h"
#include "control.h"
#include "dimradmars.h"
#include "yomaer.h"
#include "yomrdu.h"
#include "comg1d.h"

#include "paramet.h"
#include "comvert.h"

c Arguments :
c -----------

c   inputs:
c   -------
      INTEGER ngrid,nlayer,nq

      REAL ecri
      REAL ptimestep
      REAL pplev(ngrid,nlayer+1),pplay(ngrid,nlayer)
      REAL pphi(ngrid,nlayer)
      REAL pu(ngrid,nlayer),pv(ngrid,nlayer)
      REAL pt(ngrid,nlayer),pq(ngrid,nlayer,nq)
      REAL pdu(ngrid,nlayer),pdv(ngrid,nlayer)

      CHARACTER*80 fichier
      data      fichier /'startfi'/

c  dynamical tendencies
      REAL pdtdyn(ngrid,nlayer),pdqdyn(ngrid,nlayer,nq)
      REAL pdudyn(ngrid,nlayer),pdvdyn(ngrid,nlayer)
      REAL pw(ngridmx,nlayer)

c  Time
      INTEGER pday
      REAL ptime

c   outputs:
c   --------

c  physical tendencies
      REAL pdt(ngrid,nlayer),pdq(ngrid,nlayer,nq)
      REAL pdpsrf(ngrid)
      LOGICAL firstcall,lastcall

c Local variables :
c -----------------
      CHARACTER*80 nomvar

      INTEGER l,ig,ierr,aslun,nlevel,igout,cllun,iq,itra
      integer ig1d
      INTEGER tapphys

      INTEGER*4 day_ini
      REAL time,zday
      REAL zh(ngridmx,nlayermx),z1,z2
      REAL zzlev(ngridmx,nlayermx+1),zzlay(ngridmx,nlayermx)
      REAL ztlev(ngridmx,nlayermx+1)
      REAL zdvfr(ngridmx,nlayermx),zdufr(ngridmx,nlayermx)
      REAL zdhfr(ngridmx,nlayermx),zdtsrf(ngridmx),zdtsrfr(ngridmx)
      REAL zdtc(ngridmx,nlayermx),zdtsrfc(ngridmx)
      REAL zdvc(ngridmx,nlayermx),zduc(ngridmx,nlayermx)
      REAL zflubid(ngridmx)
      REAL zplanck(ngridmx),zpopsk(ngridmx,nlayermx)
      REAL zrad(ngridmx,nlayermx),zradc(ngridmx,nlayermx)
      REAL zdum1(ngridmx,nlayermx)
      REAL zdum2(ngridmx,nlayermx)
      REAL zdum3(ngridmx,nlayermx)
      REAL zdum4(ngridmx,nlayermx)
      REAL zdum5(ngridmx,nlayermx)
      REAL stephan
      REAL ztim1,ztim2,ztim3
      REAL zco2,zp
      REAL zmu(ngridmx)
      REAL zvar,zls,zaervar
      REAL zdtlw(ngridmx,nlayermx),zdtsw(ngridmx,nlayermx)
      REAL zdtlwcl(ngridmx,nlayermx),zdtswcl(ngridmx,nlayermx)
      REAL zflux(ngridmx,6)
      REAL ztime_fin
      REAL zdaerosolc(ngridmx,nlayermx)
      REAL q2(ngridmx,nlayermx+1)
      REAL zulow(ngridmx),zvlow(ngridmx)
      REAL zustr(ngridmx),zvstr(ngridmx)
      REAL co2heat0,p0nonlte
      REAL zlsconst

      REAL sigtest(nlayermx+1)
      INTEGER igwd,igwdim,itest(ngridmx)
      logical ll

c Local saved variables:
c ----------------------

      save day_ini

      INTEGER icount
      SAVE icount

      REAL aerosol(ngridmx,nlayermx,5),cst_aer,pview(ngridmx)
      REAL tsurf(ngridmx),tsoil(ngridmx,nsoilmx)
      REAL co2ice(ngridmx),albedo(ngridmx,2),emis(ngridmx)
      REAL dtrad(ngridmx,nlayermx),fluxrad(ngridmx)
      REAL capcal(ngridmx),fluxgrd(ngridmx)
      REAL topdust0(ngridmx) , topdust
      SAVE aerosol,cst_aer,pview
      SAVE tsurf,tsoil
      SAVE co2ice,albedo,emis
      SAVE q2
      SAVE capcal,fluxgrd,dtrad,fluxrad
      SAVE topdust0

      LOGICAL ldrs,startdrs
      SAVE ldrs,startdrs
      SAVE stephan

      EXTERNAL vdif,convadj
      EXTERNAL orbite,mucorr
      EXTERNAL ismin,ismax

      DATA stephan/5.67e-08/

c Local variable for cutting of physiq call
c -----------------------------------------
      INTEGER ndomain
      parameter (ndomain = (ngridmx-1) / ndomainsz + 1)

c  cutting of radite call
      REAL zplev(ndomainsz,nlayermx+1)
      REAL zztlev(ndomainsz,nlayermx+1)
      REAL zplay(ndomainsz,nlayermx)
      REAL zt(ndomainsz,nlayermx)
      REAL zaerosol(ndomainsz,nlayermx,5)
      REAL zalbedo(ndomainsz,2)

      REAL zzrad(ndomainsz,nlayermx)
      REAL zzradc(ndomainsz,nlayermx)
      REAL zzdtlw(ndomainsz,nlayermx)
      REAL zzdtsw(ndomainsz,nlayermx)
      REAL zzdtlwcl(ndomainsz,nlayermx)
      REAL zzdtswcl(ndomainsz,nlayermx)
      REAL zzflux(ndomainsz,6)

      INTEGER jd,j,ig0,nd

c  cutting of Gravity wave drag call
      REAL zu(ndomainsz,nlayermx)
      REAL zv(ndomainsz,nlayermx)
      INTEGER zidx(ndomainsz)
      REAL zzdhfr(ndomainsz,nlayermx)
      REAL zzdufr(ndomainsz,nlayermx)
      REAL zzdvfr(ndomainsz,nlayermx)

c local variables only used for diagnostic (diagfi)
c -------------------------------------------------
      REAL ps(ngridmx)
      real ztmp_cl(ngridmx,nlayermx)
      real ztmp_cv(ngridmx,nlayermx)
      real ztmp_gw(ngridmx,nlayermx)

      real zu_cl(ngridmx,nlayermx)
      real zu_cv(ngridmx,nlayermx)
      real zu_gw(ngridmx,nlayermx)

      real zv_cl(ngridmx,nlayermx)
      real zv_cv(ngridmx,nlayermx)
      real zv_gw(ngridmx,nlayermx)

c pour etude des voisins proches
c-------------------------------

      real  netrad (ngridmx,nlayermx)            ! radiative budget (W/m2)
      real  znetrad(ndomainsz,nlayermx)

c=======================================================================
c
c Fonction Intrinseque
c --------------------
c     To be used for seasonal  varying dust:
      zaervar(zls)=0.7+.3*cos(zls+80.*pi/180.)

c-----------------------------------------------------------------------
c 1. Initialisation:
c -----------------
      zday=pday+ptime
      ldrs=.true.
c     igout=ngrid/2+1
      igout=60

c     Initialisation only at first call
c     --------------------------------------
      IF(firstcall) THEN

c        variables set to 0
c        """"""""""""""""""
         do itra=1,5
            do l=1,nlayer
               do ig=1,ngrid
                  aerosol(ig,l,itra)=0.
               enddo
            enddo
         enddo

         DO l=1,nlayer
            DO ig=1,ngrid
               dtrad(ig,l)=0.
            ENDDO
         ENDDO
         DO ig=1,ngrid
            fluxrad(ig)=0.
         ENDDO

c        reading startfi
c        """""""""""""""
         CALL readfi(fichier,0,0,ngrid,nsoilmx,startdrs,
     .      day_ini,time,co2ice,tsurf,tsoil,emis,q2)
cccccccccccccccccccccccccccccccccccccccccccccccccc
Chris met emis a 1   (il faut en + mettre callcond a FALSE dans callphys.def)
c        DO ig=1,ngrid
c           emis(ig)=1.
c        ENDDO
cccccccccccccccccccccccccccccccccccccccccccccccccc
         write (*,*) 'In physic day_ini =', day_ini
         write (*,*) 'emis au premier appel:',emis

         CALL surfini(ngrid,co2ice,albedo)
         CALL iniorbit(aphelie,periheli,year_day,peri_day,obliquit)

c        initialisation radiative transfer and dust
c        """"""""""""""""""""""""""""""""""""""""""
c        Total optical depth of "tauvis" for Ps=700Pa:
         cst_aer = tauvis/700.
c        altitude of the top of the aerosol layer (km) at Ls=2.76rad:
         DO ig=1,ngrid
             topdust0(ig)=60. -22.*SIN(lati(ig))**2
         END DO
         DO ig=1,ngrid
            pview(ig)=1.66     ! cosecant of viewing angle
         ENDDO
         CALL sucst(66,19900101,8000000,g,rad,mugaz,rcp,daysec)
         CALL surad(6)

c        initialisation soil
c        """""""""""""""""""
         IF(callsoil) THEN
            CALL soil(ngrid,nsoilmx,firstcall,inertiedat,
     s          ptimestep,tsurf,tsoil,capcal,fluxgrd)
         ELSE
            PRINT*,'WARNING!!! Thermal conduction in the soil
     s      turned off'
            DO ig=1,ngrid
               capcal(ig)=1.e5
               fluxgrd(ig)=0.
            ENDDO
         ENDIF
         icount=1

         IF (ngrid.EQ.1) GOTO 10002

c        On cree le fichier bande histoire histfi
c        """"""""""""""""""""""""""""""""""""""""
c        ATTENTION : On ecrase le "histphy" cree par gcm.F !
         tapphys  = 12
         if (tapphys.ne.12) stop 'Pb d''indice du fichier histfi'
         close (tapphys,STATUS='delete')
         close (tapphys+1,STATUS='delete')
         ierr = cllun (tapphys)

         ierr = aslun(tapphys,'histfi.dic',
     .             tapphys+1,'histfi.dat',IDRS_CREATE)
         IF (ierr.NE.0) THEN
            WRITE(6,*)' Pb d''ouverture du fichier histfi'
            WRITE(6,*)' ierr = ', ierr
            CALL EXIT(1)
         ENDIF
         CALL iniwritefi(tapphys,ldrs,pday)
c        Verification que la synchro des fichiers histoire dyn et phy est
c        possible (F.F)
         IF( MOD(iecri*day_step,iphysiq).ne.0) then
              write (*,*) 'pb de synchro a l ecriture fichier histoire'
              write (*,*) 'dyn et phy. Il faut iecri*day_step '
              write (*,*) 'multiple de iphysiq '
              stop
         ENDIF


10002    CONTINUE

c        initialisation Gravity Wave & orgraphy friction scheme
c        """""""""""""""""""""""""""""""""""""""""""""""""""""""
         if(calllott) then
            do l=1,nlayermx+1
              sigtest(l)=pplev(1,l)/pplev(1,1)
            enddo
            call sugwd(nlayermx,sigtest)
            print*,'WARNING!!! sigma levels are assumed'
         endif
c!-*-

c        initialisation pour les statistiques
c        """"""""""""""""""""""""""""""""""""

         IF (callstats) THEN
            CALL instat
         ENDIF

      ENDIF ! 	(du "if firstcall")

c    Initialisations at every physical timestep:
c    ------------------------------------------

      CALL solarlong(zday,zls)
c     print*,'zday = ',zday,' Ls= ',zls*180./pi

c     Temporal variation of the dust optical depth as a function
c     of solar longitude.
      IF(laervar) THEN
            zvar=zaervar(zls)
            WRITE(37,'(3e15.5)') zls,zvar,zday
      ELSE
            zvar=1.
      ENDIF


      IF (ngrid.NE.ngridmx) THEN
         PRINT*,'STOP in inifis'
         PRINT*,'Probleme de dimenesions :'
         PRINT*,'ngrid     = ',ngrid
         PRINT*,'ngridmx   = ',ngridmx
         STOP
      ENDIF

      DO l=1,nlayer
         DO ig=1,ngrid
            pdv(ig,l)=0.
            pdu(ig,l)=0.
            pdt(ig,l)=0.
         ENDDO
      ENDDO

      DO ig=1,ngrid
         pdpsrf(ig)=0.
         zflubid(ig)=0.
         zdtsrf(ig)=0.
      ENDDO
      CALL zerophys(ngrid*nlayer,zdum1)
      CALL zerophys(ngrid*nlayer,zdum2)
      CALL zerophys(ngrid*nlayer,zdum3)
      CALL zerophys(ngrid*nlayer,zdum4)
      CALL zerophys(ngrid*nlayer,zdum5)

c   calcul du geopotentiel aux niveaux intercouches
c   -----------------------------------------------
c   ponderation des altitudes au niveau des couches en dp/p

      DO l=1,nlayer
         DO ig=1,ngrid
            zzlay(ig,l)=pphi(ig,l)/g
         ENDDO
      ENDDO
      DO ig=1,ngrid
         zzlev(ig,1)=0.

      ENDDO
      DO l=2,nlayer
         DO ig=1,ngrid
            z1=(pplay(ig,l-1)+pplev(ig,l))/(pplay(ig,l-1)-pplev(ig,l))
            z2=(pplev(ig,l)+pplay(ig,l))/(pplev(ig,l)-pplay(ig,l))
            zzlev(ig,l)=(z1*zzlay(ig,l-1)+z2*zzlay(ig,l))/(z1+z2)
         ENDDO
      ENDDO


c   Transformation de la temperature en temperature potentielle
c   -----------------------------------------------------------
      DO l=1,nlayer
         DO ig=1,ngrid
            zpopsk(ig,l)=(pplay(ig,l)/pplev(ig,1))**rcp
         ENDDO
      ENDDO
      DO l=1,nlayer
         DO ig=1,ngrid
            zh(ig,l)=pt(ig,l)/zpopsk(ig,l)
         ENDDO
      ENDDO

c-----------------------------------------------------------------------
c    Stockage des statistiques :
c    -------------------------

         IF (callstats) THEN
            CALL stats(ngrid, nlayer, co2ice, emis, tsurf
     &                ,pplev, pplay, pt, pu, pv, q2)

         ENDIF

c-----------------------------------------------------------------------
      print*,''
      print*,''
      print*,'CALL of PHYSIQ ',icount

c-----------------------------------------------------------------------
c    2. Calcul of the radiative tendencies :
c    ---------------------------------------

      IF(callrad) THEN
         IF( MOD(icount-1,iradia).EQ.0) THEN

       print*, 'CALL of RADITE'
c      print*, 'icount =',icount,' et iradia=',iradia


c    2.1 Useful calculations to prepare radiative transfer
c    ------------------------------------------------------

c    Solar angle
c    ~~~~~~~~~~~

            CALL orbite(zls,dist_sol,declin)
            IF(diurnal) THEN
               ztim1=SIN(declin)
               ztim2=COS(declin)*COS(2.*pi*(zday-.5))
               ztim3=-COS(declin)*SIN(2.*pi*(zday-.5))
               CALL solang(ngrid,sinlon,coslon,sinlat,coslat,
     s         ztim1,ztim2,ztim3,
     s         mu0,fract)
               IF(lwrite) THEN
                  PRINT*,'day, declin, sinlon,coslon,sinlat,coslat'
                  PRINT*,zday, declin,
     s            sinlon(igout+1),coslon(igout+1),
     s            sinlat(igout),coslat(igout)
               ENDIF
            ELSE
               CALL mucorr(ngrid,declin, lati, mu0, fract,10000.,rad)
            ENDIF

c        Aerosol distribution:
c        ~~~~~~~~~~~~~~~~~~~~~

         if(esadust) then
c           new dust vertical distribution function:

            zlsconst=18.*SIN(zls-2.76)
            DO l=1,nlayer
               DO ig=1,ngrid
                 if(laervarlat) then
                   topdust=topdust0(ig)+zlsconst    ! viking year
                 else
                   topdust=topdustref               ! like before 03/96
                 endif
                 zp=(700./pplay(ig,l))**(70./topdust)
                 aerosol(ig,l,1)=zvar*cst_aer*
     s           (pplev(ig,l)-pplev(ig,l+1))
     s           *max( exp(.007*(1.-max(zp,1.))) , 1.E-3 )
               ENDDO
            ENDDO

            print*,'Ls= ',zls*180./pi,'   TAU=',zvar*tauvis
     &            ,'    esadust   D=',dist_sol
         else
c           old dust vertical distribution function ("conrath"):

            print*,'Ls= ',zls*180./pi,'   TAU=',zvar*tauvis
            DO l=1,nlayer
               DO ig=1,ngrid

                  zp=700./pplay(ig,l)
                  aerosol(ig,l,1)=
     s            zvar*cst_aer*(pplev(ig,l)-pplev(ig,l+1))
     s            *max( exp(.03*(1.-max(zp,1.))) , 1.E-3 )

               ENDDO
            ENDDO
         endif



c        Cutting variables for radite (and save memory)
c        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
      if (ngridmx .EQ. 1) then
        if (ndomainsz .NE. 1) then
          print*
          print*,'ATTENTION !!!'
          print*,'pour tourner en 1D, '
          print*,'fixer ndomainsz=1 dans phymars/dimradmars.h'
          print*
          call exit(1)
        endif
      endif

      IF(firstcall) THEN
      print*
      print*
      print*,'Radite cutting: ngridmx,ngrid,ndomainsz,ndomain',
     .		ngridmx,ngrid,ndomainsz,ndomain
      ENDIF

      DO jd=1,ndomain
        ig0=(jd-1)*ndomainsz
        if (jd.eq.ndomain) then
      	 nd=ngridmx-ig0
      	else	
      	 nd=ndomainsz
      	endif

      IF(firstcall) THEN
      print*,'domain number',jd, '   from ig0+1 = ',ig0+1,
     .  ' to ig0+nd = ',ig0+nd,' (',nd,')'
      ENDIF

c    Cutting of entries
      	do l=1,nlayer+1
      	 do ig = 1,nd
      	  zplev(ig,l) = pplev(ig0+ig,l)
      	 enddo
      	enddo

      	do l=1,nlayer
      	 do ig = 1,nd
      	  zplay(ig,l) = pplay(ig0+ig,l)
      	  zt(ig,l) = pt(ig0+ig,l)
      	 enddo
      	enddo

      	do j=1,5
      	 do l=1,nlayer
      	  do ig = 1,nd
      	   zaerosol(ig,l,j) = aerosol(ig0+ig,l,j)
      	  enddo
      	 enddo
      	enddo

      	do j=1,2
      	  do ig = 1,nd
      	   zalbedo(ig,j) = albedo(ig0+ig,j)
      	  enddo
      	enddo

c Intermediate  levels: (computing ztlev)
c ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
c Extrapolation for the air temperature above the surface
 	    DO ig=1,nd
              zztlev(ig,1)=zt(ig,1)+
     s        (zplev(ig,1)-zplay(ig,1))*
     s        (zt(ig,1)-zt(ig,2))/(zplay(ig,1)-zplay(ig,2))

c  Average with tsurf
              zztlev(ig,1)=0.5*(zztlev(ig,1)+tsurf(ig0+ig))

        ENDDO

        DO l=2,nlayer
 	     DO ig=1,nd
               zztlev(ig,l)=0.5*(zt(ig,l-1)+zt(ig,l))
         ENDDO
        ENDDO

c Value for the air temperature on the top of atmosphere
        nlevel=nlayer+1
 	     DO ig=1,nd
               zztlev(ig,nlevel)=zt(ig,nlayer)
         ENDDO

c    2.2 Calcul of the radiative tendencies and fluxes:
c    --------------------------------------------------

            zco2=.95
            CALL radite(ig0,icount,nd,nlayer,0,6,1,1,
     $         zaerosol,zalbedo,zco2,zdum4,emis(ig0+1),
     $         mu0(ig0+1),zdum5,
     $         zplev,zplay,zdum1,zdum2,zdum3,zztlev,tsurf(ig0+1),zt,
     $		   pview(ig0+1),
     $         zzdtlw,zzdtsw,zzdtlwcl,zzdtswcl,zzflux,zzrad,zzradc,
     $         fract(ig0+1),dist_sol
     .        ,znetrad)


c    2.3 total flux and tendencies:
c    ------------------------------

c    2.3.1 fluxes
c    ~~~~~~~~~~~~

        do l=1,nlayer
         do ig = 1,nd
          zrad(ig0+ig,l) = zzrad(ig,l)
          zradc(ig0+ig,l) = zzradc(ig,l)
          zdtlw(ig0+ig,l) = zzdtlw(ig,l)

cccccccccccccccccccccccccccccccccccccccccccccccccc
Chris coupe le sw
          zdtsw(ig0+ig,l) = zzdtsw(ig,l)
c         zdtsw(ig0+ig,l) = 0
cccccccccccccccccccccccccccccccccccccccccccccccccc

          zdtlwcl(ig0+ig,l) = zzdtlwcl(ig,l)
          zdtswcl(ig0+ig,l) = zzdtswcl(ig,l)
          netrad(ig0+ig,l) = znetrad(ig,l)
         enddo
        enddo

        do l=1,nlayer+1
         do ig = 1,nd
          ztlev(ig0+ig,l) = zztlev(ig,l)
         enddo
        enddo

        do j=1,6
         do ig = 1,nd
          zflux(ig0+ig,j) = zzflux(ig,j)
         enddo
        enddo

      ENDDO			! 	(boucle jd=1, ndomain)


            DO ig=1,ngrid

c si on veut on enleve ici le lw dans le bilan du sol avec zflux(ig,4)=0
c              fluxrad(ig)= 0
               fluxrad(ig)=zflux(ig,4)

     $         +zflux(ig,5)*(1.-albedo(ig,1))
     $         +zflux(ig,6)*(1.-albedo(ig,2))
               zplanck(ig)=tsurf(ig)*tsurf(ig)
               zplanck(ig)=emis(ig)*
     $         stephan*zplanck(ig)*zplanck(ig)
               fluxrad(ig)=fluxrad(ig)-zplanck(ig)
           ENDDO

c    2.4.2 temperature tendencies
c    ~~~~~~~~~~~~~~~~~~~~~~~~~~~
            DO l=1,nlayer
               DO ig=1,ngrid
                  dtrad(ig,l)=(zdtsw(ig,l)+zdtlwcl(ig,l))/daysec
c._.                  dtrad(ig,l)=(zdtsw(ig,l)+zdtlw(ig,l))/daysec
c             print*,'zdtlw, zdtsw',l,zdtlwcl(ig,l),zdtsw(ig,l)
               ENDDO
            ENDDO

cc   SORTIES SPECIALES DEBUG!
cc   ~~~~~~~~~~~~~~~~~~~~~~~~
cc      tests pour debuguer en cas de valeurs NaN sur fluxrad
c      do ig=1,ngrid
c      if(.not.fluxrad(ig).le.500.) then
cc     if(1.eq.0) then
c        jjj=(ig-2)/iim+2
c        iii=ig-1-iim*(jjj-2)
c        print*,'Sortie speciale debug physiq'
c        print*,'Ca plante au point ',ig,iii,jjj
c        print*,'mu0=',mu0(ig),'  fract=',fract(ig)
c        print*,'Alb.=',albedo(ig,1),albedo(ig,2),' emis',emis(ig)
c        print*,'Flux rad : ',fluxrad(ig)
c        print*,'l    P   P     T      T      dtlw     dtsw '
c        do l=1,nlayer
c           print*,l,pplev(ig,l),pplay(ig,l),ztlev(ig,l)
c    ,       ,pt(ig,l),zdtlwcl(ig,l),zdtsw(ig,l),aerosol(ig,l,1)
c        enddo
cc   appel du rayonnement avec diagnostiques au point de plantage
c        nimp=0
c        jlimprad=ig
c           CALL radite(ngrid,nlayer,0,6,1,1,
c    $         aerosol,albedo,zco2,zdum4,emis,
c    $         mu0,zdum5,
c    $         pplev,pplay,zdum1,zdum2,zdum3,ztlev,tsurf,pt,pview,
c    $         zdtlw,zdtsw,zdtlwcl,zdtswcl,zflux,zrad,zradc,
c    $         fract,dist_sol)
c
c        do l=1,nlayer
c           print*,l,pplev(ig,l),pplay(ig,l),ztlev(ig,l)
c    ,       ,pt(ig,l),zdtlwcl(ig,l),zdtsw(ig,l),aerosol(ig,l,1)
c        enddo
c        stop
c      endif
c      enddo


c    2.4 CO2 near infrared absorption
c    -------------------------------

            if(callco2v) then
c co2heat is the heating by CO2 at 700Pa for a zero zenithal angle.
c it is 1.3 K/day for the mean distance to sun (1.52UA).
               co2heat0=1.3*(1.52/dist_sol)**2/daysec
c p0noonlte is a pressure below which non LTE effects are significant.
c  Modif Fred H., 08 08 96
               p0nonlte=7.5e-3
               do ig=1,ngrid
                  zmu(ig)=sqrt(1224.*mu0(ig)*mu0(ig)+1.)/35.
               enddo
               if(callnlte) then
                 do l=1,nlayer
                    do ig=1,ngrid
                     if(fract(ig).gt.0.) dtrad(ig,l)=dtrad(ig,l)
     s               + co2heat0*sqrt((700.*zmu(ig))/pplay(ig,l))
     s               / (1+p0nonlte/pplay(ig,l))
                    enddo
                 enddo
               else
                 do l=1,nlayer
                    do ig=1,ngrid
                     if(fract(ig).gt.0.) dtrad(ig,l)=dtrad(ig,l)
     s               + co2heat0/sqrt((700.*zmu(ig))/pplay(ig,l))
                    enddo
                 enddo
               endif
            endif


        ENDIF ! mod(icount-1,iradia).eq.0


c    2.5 Transformation of the radiative tendencies:
c    -----------------------------------------------

         DO l=1,nlayer
            DO ig=1,ngrid
               pdt(ig,l)=pdt(ig,l)+dtrad(ig,l)
c            print*,'pdt',l,pdt(ig,l)
            ENDDO
         ENDDO

         IF(lwrite) THEN
            PRINT*,'Diagnotique for the radiation'
            PRINT*,'albedo, emissiv, mu0,fract,pview,Frad,Planck'
            PRINT*,albedo(igout,1),emis(igout),mu0(igout),
     s           fract(igout),pview(igout),
     s           fluxrad(igout),zplanck(igout)
            PRINT*,'Tlay Tlev Play Plev dT/dt SW dT/dt LW (K/day)'
            PRINT*,'daysec',daysec
            DO l=1,nlayer
               PRINT*,pt(igout,l),ztlev(igout,l),
     s         pplay(igout,l),pplev(igout,l),
     s         zdtsw(igout,l),zdtlw(igout,l)
            ENDDO
         ENDIF


      ENDIF

c-----------------------------------------------------------------------
c    3. Gravity wave drag :
c    ----------------------
c
c Appeler le programme de parametrisation de l'orographie
c a l'echelle sous-maille:
c

      IF(calllott)THEN

c        Cutting variables for gravity wave drag (and save memory)
c        ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
c     write(*,*)
c    .	'Gravity wave drag cutting: ngridmx,ngrid,ndomainsz,ndomain',
c    .	ngridmx,ngrid,ndomainsz,ndomain

      DO jd=1,ndomain
        ig0=(jd-1)*ndomainsz
        if (jd.eq.ndomain) then
         nd=ngridmx-ig0
        else
         nd=ndomainsz
        endif

c  selection des points pour lesquels le shema est actif:
        igwd=0
        DO ig=ig0+1,ig0+nd
        itest(ig)=0
cmod        ll=zstd(ig).gt.10.0
        ll=zstd(ig).gt.50.0
        IF(ll) then
          itest(ig)=1
          igwd=igwd+1
          zidx(igwd)=ig - ig0
        ENDIF
        ENDDO
        IGWDIM=MAX(1,IGWD)

c Cutting of entries
        do l=1,nlayer+1
         do ig = 1,nd
          zplev(ig,l) = pplev(ig0+ig,l)
         enddo
        enddo

        do l=1,nlayer
         do ig = 1,nd
          zplay(ig,l) = pplay(ig0+ig,l)
          zt(ig,l) = pt(ig0+ig,l)
          zu(ig,l) = pu(ig0+ig,l)
          zv(ig,l) = pv(ig0+ig,l)
         enddo
        enddo

c       Call gravity wave drag
c       ~~~~~~~~~~~~~~~~~~~~~~
        call drag_noro (nd,nlayer,ptimestep,zplay,zplev,
     e            zstd(ig0+1),zsig(ig0+1),zgam(ig0+1),zthe(ig0+1),
     e            igwd,igwdim,zidx,itest(ig0+1),
     e            zt, zu, zv,
     s            zulow(ig0+1),zvlow(ig0+1),zustr(ig0+1),zvstr(ig0+1),
     s            zzdhfr,zzdufr,zzdvfr)


        do l=1,nlayer
         do ig = 1,nd
          zdhfr(ig0+ig,l) = zzdhfr(ig,l)
          zdufr(ig0+ig,l) = zzdufr(ig,l)
          zdvfr(ig0+ig,l) = zzdvfr(ig,l)
         enddo
        enddo

      ENDDO         !   (boucle jd=1, ndomain)


c  ajout des tendances gravity wave drag
c  ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
        DO l=1,nlayer
          DO ig=1,ngrid
            pdv(ig,l)=pdv(ig,l)+zdvfr(ig,l)/ptimestep
            pdu(ig,l)=pdu(ig,l)+zdufr(ig,l)/ptimestep
            pdt(ig,l)=pdt(ig,l)+zdhfr(ig,l)/ptimestep

            ztmp_gw(ig,l)=zdhfr(ig,l)/ptimestep ! for diagnostic only
            zu_gw(ig,l)=zdufr(ig,l)/ptimestep ! for diagnostic only
            zv_gw(ig,l)=zdvfr(ig,l)/ptimestep ! for diagnostic only

          ENDDO
        ENDDO
      ENDIF

c-----------------------------------------------------------------------
c    4. Vertical diffusion (turbulent mixing):
c    -----------------------------------------
c
      IF(calldifv) THEN

         DO ig=1,ngrid
            zflubid(ig)=fluxrad(ig)+fluxgrd(ig)
         ENDDO

         CALL zerophys(ngrid*nlayer,zdum1)
         CALL zerophys(ngrid*nlayer,zdum2)
         CALL zerophys(ngrid*nlayer,zdum3)
         do l=1,nlayer
            do ig=1,ngrid
               zdum3(ig,l)=pdt(ig,l)/zpopsk(ig,l)
            enddo
         enddo

c        Appel vdif (version martienne AVEC condensation)
c        ~~~~~~~~~
         CALL vdifc(ngrid,nlayer,co2ice, zpopsk,
     $        ptimestep,capcal,
     $        pplay,pplev,zzlay,zzlev,
     $        z0,
     $        pu,pv,zh,tsurf,emis,
     $        zdum1,zdum2,zdum3,zflubid,
     $        zdufr,zdvfr,zdhfr,zdtsrfr,q2,
     $        lwrite)

         DO l=1,nlayer
            DO ig=1,ngrid
               pdv(ig,l)=pdv(ig,l)+zdvfr(ig,l)
               pdu(ig,l)=pdu(ig,l)+zdufr(ig,l)
               pdt(ig,l)=pdt(ig,l)+zdhfr(ig,l)*zpopsk(ig,l)

               ztmp_cl(ig,l)=zdhfr(ig,l)*zpopsk(ig,l) ! for diagnostic only
               zu_cl(ig,l)=zdufr(ig,l)                ! for diagnostic only
               zv_cl(ig,l)=zdvfr(ig,l)                ! for diagnostic only

            ENDDO
         ENDDO

         DO ig=1,ngrid
            zdtsrf(ig)=zdtsrf(ig)+zdtsrfr(ig)
         ENDDO

      ELSE
         DO ig=1,ngrid
            zdtsrf(ig)=zdtsrf(ig)+
     s      (fluxrad(ig)+fluxgrd(ig))/capcal(ig)
         ENDDO
      ENDIF
c
c
c
c-----------------------------------------------------------------------
c   4. Dry convective adjustment:
c   -----------------------------

      IF(calladj) THEN

         DO l=1,nlayer
            DO ig=1,ngrid
               zdum1(ig,l)=pdt(ig,l)/zpopsk(ig,l)
            ENDDO
         ENDDO
         CALL zerophys(ngrid*nlayer,zdufr)
         CALL zerophys(ngrid*nlayer,zdvfr)
         CALL zerophys(ngrid*nlayer,zdhfr)
         CALL convadj(ngrid,nlayer,ptimestep,
     $                pplay,pplev,zpopsk,
     $                pu,pv,zh,
     $                pdu,pdv,zdum1,
     $                zdufr,zdvfr,zdhfr)

         DO l=1,nlayer
            DO ig=1,ngrid
               pdu(ig,l)=pdu(ig,l)+zdufr(ig,l)
               pdv(ig,l)=pdv(ig,l)+zdvfr(ig,l)
               pdt(ig,l)=pdt(ig,l)+zdhfr(ig,l)*zpopsk(ig,l)

               ztmp_cv(ig,l)=zdhfr(ig,l)*zpopsk(ig,l) ! for diagnostic only
               zu_cv(ig,l)=zdufr(ig,l)                ! for diagnostic only
               zv_cv(ig,l)=zdvfr(ig,l)                ! for diagnostic only

            ENDDO
         ENDDO



      ENDIF
c-----------------------------------------------------------------------
c   5. Carbon dioxide condensation-sublimation:
c   -------------------------------------------

      IF(callcond) THEN
         CALL newcondens(ngrid,nlayer,ptimestep,
     $              capcal,pplay,pplev,tsurf,pt,
     $              pphi,pdt,pdu,pdv,zdtsrf,pu,pv,aerosol,
     $              co2ice,albedo,emis,
     $              zdtc,zdtsrfc,pdpsrf,zduc,zdvc,
     $              zdaerosolc)


         DO l=1,nlayer
           DO ig=1,ngrid
             pdt(ig,l)=pdt(ig,l)+zdtc(ig,l)
             pdv(ig,l)=pdv(ig,l)+zdvc(ig,l)
             pdu(ig,l)=pdu(ig,l)+zduc(ig,l)
         ENDDO
         ENDDO
         DO ig=1,ngrid
            zdtsrf(ig) = zdtsrf(ig) + zdtsrfc(ig)
         ENDDO
      ENDIF
c-----------------------------------------------------------------------
c   6. Traitement des traceurs dans la physique
c:   -------------------------------------------

      DO iq=1,nq
        DO l=1,nlayer
            DO ig=1,ngrid
               pdq(ig,l,iq) = 0.
            END DO
        END DO
      END DO



c**********************************************************************
#ifdef undim

Chris ecrit les sorties 1D avant d'ajouter les tendances

      ig1d = ngridmx/2 + 1
      print*,'ig1d pour grads1d',ig1d

      g1d_tmp1='temp'
      g1d_tmp2='temperature'
      CALL writeg1d(1,nlayer,pt(ig1d,1),g1d_tmp1,g1d_tmp2)

      if (1)  then
      g1d_tmp1='plev'
      g1d_tmp2='plev '
      CALL writeg1d(1,nlayer,pplev(ig1d,1),g1d_tmp1,g1d_tmp2)

      g1d_tmp1='ptop'
      g1d_tmp2='plev pour interface nlayer+1'
      CALL writeg1d(1,1,pplev(ig1d,nlayer+1),g1d_tmp1,g1d_tmp2)

      g1d_tmp1='ds'
      g1d_tmp2='ds '
      CALL writeg1d(1,nlayer,ds(1),g1d_tmp1,g1d_tmp2)

      g1d_tmp1='dsig'
      g1d_tmp2='dsig '
      CALL writeg1d(1,nlayer,dsig(1),g1d_tmp1,g1d_tmp2)

      g1d_tmp1='u'
      g1d_tmp2='vent : composante vers l est'
      CALL writeg1d(1,nlayer,pu(ig1d,1),g1d_tmp1,g1d_tmp2)

      g1d_tmp1='v'
      g1d_tmp2='vent : composante vers le nord'
      CALL writeg1d(1,nlayer,pv(ig1d,1),g1d_tmp1,g1d_tmp2)

      g1d_tmp1='p'
      g1d_tmp2='pression (play)'
      CALL writeg1d(1,nlayer,pplay(ig1d,1),g1d_tmp1,g1d_tmp2)

      g1d_tmp1='psurf'
      g1d_tmp2='pression au sol'
      CALL writeg1d(1,1,pplev(ig1d,1),g1d_tmp1,g1d_tmp2)

      g1d_tmp1='tsurf'
      g1d_tmp2='temperature au sol'
      CALL writeg1d(1,1,tsurf(ig1d),g1d_tmp1,g1d_tmp2)

      do l = 1, nsoilmx
      if (l .LE. 9) then
        write(g1d_tmp1,'(a2,i1)') 'tg',l
      else
        write(g1d_tmp1,'(a2,i2)') 'tg',l
      endif
      write(g1d_tmp2,*) 'temperature du sol couche ',l
      CALL writeg1d(1,1,tsoil(ig1d,l),g1d_tmp1,g1d_tmp2)
      enddo

      g1d_tmp1='dtlw'
      g1d_tmp2='cooling rate lw (K/j)'
      CALL writeg1d(1,nlayer,zdtlwcl(ig1d,1),g1d_tmp1,g1d_tmp2)

      g1d_tmp1='dtsw'
      g1d_tmp2='cooling rate sw (K/j)'
      CALL writeg1d(1,nlayer,zdtsw(ig1d,1),g1d_tmp1,g1d_tmp2)

      g1d_tmp1='net'
      g1d_tmp2='radiative budget lw (W/m2)'
      CALL writeg1d(1,nlayer,netrad(ig1d,1),g1d_tmp1,g1d_tmp2)

      g1d_tmp1='fmlw'
      g1d_tmp2='flux montant lw'
      CALL writeg1d(1,1,zplanck(ig1d),g1d_tmp1,g1d_tmp2)

      g1d_tmp1='fdlw'
      g1d_tmp2='flux descendant lw'
      CALL writeg1d(1,1,zflux(ig1d,4),g1d_tmp1,g1d_tmp2)

      g1d_tmp1='fmsw'
      g1d_tmp2='flux montant sw'
      CALL writeg1d(1,1,zflux(ig1d,5),g1d_tmp1,g1d_tmp2)

      g1d_tmp1='fdsw'
      g1d_tmp2='flux descendant sw'
      CALL writeg1d(1,1,zflux(ig1d,6),g1d_tmp1,g1d_tmp2)

      g1d_tmp1='aer'
      g1d_tmp2='aerosol'
      CALL writeg1d(1,nlayer,aerosol(ig1d,1,1),g1d_tmp1,g1d_tmp2)

      ENDIF

#endif
c**********************************************************************


c-----------------------------------------------------------------------
c   On incremente les tendances physiques de la temperature du sol:
c   ---------------------------------------------------------------

      DO ig=1,ngrid
         tsurf(ig)=tsurf(ig)+ptimestep*zdtsrf(ig)
      ENDDO
c
c-----------------------------------------------------------------------
c   soil temperatures:
c   --------------------


      IF (callsoil) THEN
         CALL soil(ngrid,nsoilmx,.false.,inertiedat,
     s          ptimestep,tsurf,tsoil,capcal,fluxgrd)
         IF(lwrite) THEN
            print*,'capacite du sol ',capcal(ngrid/2+1)
            PRINT*,'Surface Heat capacity,conduction Flux, Ts,
     s          dTs, dt'
            PRINT*,capcal(igout),fluxgrd(igout),tsurf(igout),
     s          zdtsrf(igout),ptimestep
         ENDIF

      ENDIF

      IF (lwrite) THEN
         PRINT*,'Global diagnostics for the physics'
         PRINT*,'Variables and their increments q and dq/dt * dt'
         WRITE(*,'(a6,a10,2a15)') 'Ts','dTs','ps','dps'
         WRITE(*,'(2f10.5,2f15.5)')
     s   tsurf(igout),zdtsrf(igout)*ptimestep,
     s   pplev(igout,1),pdpsrf(igout)*ptimestep
         WRITE(*,'(a4,a6,5a10)') 'l','u','du','v','dv','T','dT'
         WRITE(*,'(i4,6f10.5)') (l,
     s   pu(igout,l),pdu(igout,l)*ptimestep,
     s   pv(igout,l),pdv(igout,l)*ptimestep,
     s   pt(igout,l),pdt(igout,l)*ptimestep,
     s   l=1,nlayer)
      ENDIF
c-----------------------------------------------------------------------
c  Debut des ecritures
c  -------------------

      IF (ngrid.NE.1) THEN              !  (si pas un modele 1d)

c-----------------------------------------------------------------------
c   ecriture de la bande histoire physique (ajout 03/06/94 F.Forget)
c   --------------------------------------
c   La physique est sauvee dans histfi juste avant que
c   la dynamique ne le soit dans histoire.
c   mais entre maintenant et l'ecriture de histoire,
c   on aura itau = itau +1 et remise a jour de time.
c   Donc on stocke lorsque mod(itau+1,iecri*day_step) =0
c   Et avec temps = zday + dtvr/day_step - day_ini
c (en fait  dtvr = ptimestep/iphysiq)

c     write (*,*) 'ds la physique, itau =', iphysiq*icount -1

        IF( MOD(iphysiq*icount,iecri*day_step).EQ.0) THEN
         ztime_fin = zday - float(day_ini)
     &               + ptimestep/(float(iphysiq)*daysec)
         IF(ldrs) THEN
            CALL WRITEFI(ngrid,nsoilmx,
     s    12,ldrs,   ztime_fin
     &    ,1,co2ice,tsurf,tsoil,emis,q2)
         ENDIF
        ENDIF

c-----------------------------------------------------------------------

c    Ecriture du fichier de reinitialisation de la physique : restartfi'
c   -------------------------------------------------------------------
c   Remarque : On  stocke restarfi
c   juste avant que la dynamique ne le soit dans restart.
c   entre maintenant et l'ecriture de restart,
c   on aura itau = itau +1 et remise a jour de time.
c   (lastcall = .true. lorsque itau+1 = itaufin)
c   Donc on stocke  avec temps = time + dtvr

      IF(lastcall) THEN
         PRINT*,'Ecriture du fichier de reinitialisation de la physique'
         IF(startdrs) THEN
            ierr = aslun(92,'restartfi.dic',
     .                93,'restartfi.dat',IDRS_CREATE)
         ELSE
            OPEN(92,file='restartfi',form='unformatted',status='new',
     .      iostat=ierr)
         ENDIF
         if (ierr.ne.0) then
           write(6,*)' Pb d''ouverture du fichier restartfi'
           write(6,*)' ierr = ', ierr
           call exit(1)
         endif
         ztime_fin = ptime + ptimestep/(float(iphysiq)*daysec)
         CALL iniwritefi(92,startdrs,pday+NINT(ztime_fin))
         ztime_fin = ztime_fin - nint(ztime_fin)
         write(*,*)'pour writefi ztime_fin =',ztime_fin
         CALL WRITEFI(ngrid,nsoilmx,
     s    92,startdrs,ztime_fin,1,co2ice,tsurf,tsoil,emis,q2)
         IF(startdrs) THEN
            ierr=cllun(92)
         ELSE
            CLOSE(92)
         ENDIF
      ENDIF

c ***************************************************************
c SORTIES DIVERSES
c ****************
c     T20m = 999.
c     emism = 99.
c     do ig=1,ngrid
c         T20m = min(T20m,Tsurf(ig) -43*(1-emis(ig)))
c         emism = min(emism,emis(ig))
c     end do
c     write(*,*)'T20min =' , T20m
c     write(*,*)'presmin =' , presm

c  SORTIE DIAGFI  (Variable Diagnostique sous format DRS)
c **************
      do ig=1,ngrid
         ps(ig) = pplev(ig,1)
      end do


c     call WRITEDRSFI(ngridmx,44,'Ls','Lon.solaire','deg',0,zls*57.296)

      call WRITEDRSFI(ngridmx,44,'t','Temperature','K',3,pt)
      call WRITEDRSFI(ngridmx,44,'u','Vent zonal','m.s-1',3,pu)
      call WRITEDRSFI(ngridmx,44,'v','Vent merid','m.s-1',3,pv)
c     call WRITEDRSFI(ngridmx,44,'W','Vent Vertical','m.s-1',3,pw)
      call WRITEDRSFI(ngridmx,44,'q2','wind variance','m2.s-2',3,q2)

      call WRITEDRSFI(ngridmx,44,'tsurf','Surf T','K',2,tsurf)
      call WRITEDRSFI(ngridmx,44,'ps','Psurf','Pa',2,ps)
      call WRITEDRSFI(ngridmx,44,'co2ice','couche de glace co2',
     &  'kg/m2',2,co2ice)
      call WRITEDRSFI(ngridmx,44,'emis','grd emis',' ',2,emis)

      do l = 1, nsoilmx
        if (l .LE. 9) then
          write(nomvar,'(a2,i1)') 'tg',l
        else
          write(nomvar,'(a2,i2)') 'tg',l
        endif
        call WRITEDRSFI(ngridmx,44,nomvar,nomvar,'K',2,tsoil(1,l))
      enddo

      call WRITEDRSFI(ngridmx,44,'dtlw','dT/dt LWc','K',3,zdtlwcl)
      call WRITEDRSFI(ngridmx,44,'dtsw','dT/dt SW','K',3,zdtsw)
c     call WRITEDRSFI(ngridmx,44,'dtrad','dT/dt ','K',3,dtrad)

c     call WRITEDRSFI(ngridmx,44,'zdtc','dt condens','K.s-1',3,zdtc)
c     call WRITEDRSFI(ngridmx,44,'zduc','du condens','',3,zduc)
c     call WRITEDRSFI(ngridmx,44,'zdvc','dv condens','',3,zdvc)

c     call WRITEDRSFI(ngridmx,44,'ztmpgw','dt gw','K.s-1',3,ztmp_gw)
c     call WRITEDRSFI(ngridmx,44,'zugw','du gw','',3,zu_gw)
c     call WRITEDRSFI(ngridmx,44,'zvgw','dv gw','',3,zv_gw)

c     call WRITEDRSFI(ngridmx,44,'ztmpcl','dt vdiff','K.s-1',3,ztmp_cl)
c     call WRITEDRSFI(ngridmx,44,'zucl','du vdiff','',3,zu_cl)
c     call WRITEDRSFI(ngridmx,44,'zvcl','dv vdiff','',3,zv_cl)

c     call WRITEDRSFI(ngridmx,44,'ztmpcv','dt conv','K.s-1',3,ztmp_cv)
c     call WRITEDRSFI(ngridmx,44,'zucv','du conv','',3,zu_cv)
c     call WRITEDRSFI(ngridmx,44,'zvcv','dv conv','',3,zv_cv)

      call WRITEDRSFI(ngridmx,44,'dt','dt','',3,pdt)
      call WRITEDRSFI(ngridmx,44,'du','du','',3,pdu)
      call WRITEDRSFI(ngridmx,44,'dv','dv','',3,pdv)


      ENDIF

      icount=icount+1

      RETURN
      END