source: LMDZ5/trunk/libf/dyn3d/leapfrog.F @ 2242

Last change on this file since 2242 was 2239, checked in by Ehouarn Millour, 9 years ago

Reorganizing physics/dynamics interface:

  • what is related to dynamics-physics interface is now in a seperate directory: dynlmdz_phy* for physics in phy*
  • 1d model and related dependencies (including a couple from "dynamics", set up as symbolic links) is now in subdirectory "dyn1d" of phy*.
  • "bibio" directory is now "misc" and should only contain autonomous utilities.
  • "cosp" is now a subdirectory of phylmd.

EM

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 26.1 KB
RevLine 
[524]1!
[1279]2! $Id: leapfrog.F 2239 2015-03-23 07:27:30Z emillour $
3!
[524]4c
5c
[2221]6      SUBROUTINE leapfrog(ucov,vcov,teta,ps,masse,phis,q,time_0)
[524]7
8
[692]9cIM : pour sortir les param. du modele dans un fis. netcdf 110106
[1146]10#ifdef CPP_IOIPSL
11      use IOIPSL
12#endif
[1987]13      USE infotrac, ONLY: nqtot
[1279]14      USE guide_mod, ONLY : guide_main
[1987]15      USE write_field, ONLY: writefield
16      USE control_mod, ONLY: nday, day_step, planet_type, offline,
17     &                       iconser, iphysiq, iperiod, dissip_period,
18     &                       iecri, ip_ebil_dyn, ok_dynzon, ok_dyn_ins,
19     &                       periodav, ok_dyn_ave, output_grads_dyn
[2021]20      use exner_hyb_m, only: exner_hyb
21      use exner_milieu_m, only: exner_milieu
22
[524]23      IMPLICIT NONE
24
25c      ......   Version  du 10/01/98    ..........
26
27c             avec  coordonnees  verticales hybrides
28c   avec nouveaux operat. dissipation * ( gradiv2,divgrad2,nxgraro2 )
29
30c=======================================================================
31c
32c   Auteur:  P. Le Van /L. Fairhead/F.Hourdin
33c   -------
34c
35c   Objet:
36c   ------
37c
38c   GCM LMD nouvelle grille
39c
40c=======================================================================
41c
42c  ... Dans inigeom , nouveaux calculs pour les elongations  cu , cv
43c      et possibilite d'appeler une fonction f(y)  a derivee tangente
44c      hyperbolique a la  place de la fonction a derivee sinusoidale.
45
46c  ... Possibilite de choisir le shema pour l'advection de
47c        q  , en modifiant iadv dans traceur.def  (10/02) .
48c
49c      Pour Van-Leer + Vapeur d'eau saturee, iadv(1)=4. (F.Codron,10/99)
50c      Pour Van-Leer iadv=10
51c
52c-----------------------------------------------------------------------
53c   Declarations:
54c   -------------
55
56#include "dimensions.h"
57#include "paramet.h"
58#include "comconst.h"
59#include "comdissnew.h"
60#include "comvert.h"
61#include "comgeom.h"
62#include "logic.h"
63#include "temps.h"
64#include "ener.h"
65#include "description.h"
66#include "serre.h"
[1403]67!#include "com_io_dyn.h"
[524]68#include "iniprint.h"
69#include "academic.h"
70
[1987]71      REAL,INTENT(IN) :: time_0 ! not used
[524]72
[1987]73c   dynamical variables:
74      REAL,INTENT(INOUT) :: ucov(ip1jmp1,llm)    ! zonal covariant wind
75      REAL,INTENT(INOUT) :: vcov(ip1jm,llm)      ! meridional covariant wind
76      REAL,INTENT(INOUT) :: teta(ip1jmp1,llm)    ! potential temperature
77      REAL,INTENT(INOUT) :: ps(ip1jmp1)          ! surface pressure (Pa)
78      REAL,INTENT(INOUT) :: masse(ip1jmp1,llm)   ! air mass
79      REAL,INTENT(INOUT) :: phis(ip1jmp1)        ! geopotentiat at the surface
80      REAL,INTENT(INOUT) :: q(ip1jmp1,llm,nqtot) ! advected tracers
81
82      REAL p (ip1jmp1,llmp1  )               ! interlayer pressure
83      REAL pks(ip1jmp1)                      ! exner at the surface
84      REAL pk(ip1jmp1,llm)                   ! exner at mid-layer
85      REAL pkf(ip1jmp1,llm)                  ! filtered exner at mid-layer
86      REAL phi(ip1jmp1,llm)                  ! geopotential
87      REAL w(ip1jmp1,llm)                    ! vertical velocity
88
[524]89      real zqmin,zqmax
90
91c variables dynamiques intermediaire pour le transport
92      REAL pbaru(ip1jmp1,llm),pbarv(ip1jm,llm) !flux de masse
93
94c   variables dynamiques au pas -1
95      REAL vcovm1(ip1jm,llm),ucovm1(ip1jmp1,llm)
96      REAL tetam1(ip1jmp1,llm),psm1(ip1jmp1)
97      REAL massem1(ip1jmp1,llm)
98
99c   tendances dynamiques
100      REAL dv(ip1jm,llm),du(ip1jmp1,llm)
[1146]101      REAL dteta(ip1jmp1,llm),dq(ip1jmp1,llm,nqtot),dp(ip1jmp1)
[524]102
103c   tendances de la dissipation
104      REAL dvdis(ip1jm,llm),dudis(ip1jmp1,llm)
105      REAL dtetadis(ip1jmp1,llm)
106
107c   tendances physiques
108      REAL dvfi(ip1jm,llm),dufi(ip1jmp1,llm)
[1146]109      REAL dtetafi(ip1jmp1,llm),dqfi(ip1jmp1,llm,nqtot),dpfi(ip1jmp1)
[524]110
111c   variables pour le fichier histoire
112      REAL dtav      ! intervalle de temps elementaire
113
114      REAL tppn(iim),tpps(iim),tpn,tps
115c
116      INTEGER itau,itaufinp1,iav
[1279]117!      INTEGER  iday ! jour julien
118      REAL       time
[524]119
120      REAL  SSUM
[1616]121!     REAL finvmaold(ip1jmp1,llm)
[524]122
[566]123cym      LOGICAL  lafin
124      LOGICAL :: lafin=.false.
[524]125      INTEGER ij,iq,l
126      INTEGER ik
127
128      real time_step, t_wrt, t_ops
129
[1279]130!      REAL rdayvrai,rdaym_ini
131! jD_cur: jour julien courant
132! jH_cur: heure julienne courante
133      REAL :: jD_cur, jH_cur
134      INTEGER :: an, mois, jour
135      REAL :: secondes
136
[524]137      LOGICAL first,callinigrads
[692]138cIM : pour sortir les param. du modele dans un fis. netcdf 110106
139      save first
140      data first/.true./
[1279]141      real dt_cum
[692]142      character*10 infile
143      integer zan, tau0, thoriid
144      integer nid_ctesGCM
145      save nid_ctesGCM
146      real degres
147      real rlong(iip1), rlatg(jjp1)
148      real zx_tmp_2d(iip1,jjp1)
149      integer ndex2d(iip1*jjp1)
150      logical ok_sync
151      parameter (ok_sync = .true.)
[1529]152      logical physic
[524]153
154      data callinigrads/.true./
155      character*10 string10
156
[960]157      REAL :: flxw(ip1jmp1,llm)  ! flux de masse verticale
[524]158
159c+jld variables test conservation energie
160      REAL ecin(ip1jmp1,llm),ecin0(ip1jmp1,llm)
161C     Tendance de la temp. potentiel d (theta)/ d t due a la
162C     tansformation d'energie cinetique en energie thermique
163C     cree par la dissipation
164      REAL dtetaecdt(ip1jmp1,llm)
165      REAL vcont(ip1jm,llm),ucont(ip1jmp1,llm)
166      REAL vnat(ip1jm,llm),unat(ip1jmp1,llm)
167      REAL      d_h_vcol, d_qt, d_qw, d_ql, d_ec
168      CHARACTER*15 ztit
[692]169!IM   INTEGER   ip_ebil_dyn  ! PRINT level for energy conserv. diag.
170!IM   SAVE      ip_ebil_dyn
171!IM   DATA      ip_ebil_dyn/0/
[524]172c-jld
173
174      character*80 dynhist_file, dynhistave_file
[1520]175      character(len=*),parameter :: modname="leapfrog"
[524]176      character*80 abort_message
177
178      logical dissip_conservative
179      save dissip_conservative
180      data dissip_conservative/.true./
181
182      LOGICAL prem
183      save prem
184      DATA prem/.true./
185      INTEGER testita
186      PARAMETER (testita = 9)
187
[1286]188      logical , parameter :: flag_verif = .false.
[999]189     
190
[825]191      integer itau_w   ! pas de temps ecriture = itap + itau_phy
192
193
[2038]194      if (nday>=0) then
195         itaufin   = nday*day_step
196      else
197         itaufin   = -nday
198      endif
[524]199      itaufinp1 = itaufin +1
[1529]200      itau = 0
201      physic=.true.
202      if (iflag_phys==0.or.iflag_phys==2) physic=.false.
[524]203
[1403]204c      iday = day_ini+itau/day_step
205c      time = REAL(itau-(iday-day_ini)*day_step)/day_step+time_0
206c         IF(time.GT.1.) THEN
207c          time = time-1.
208c          iday = iday+1
209c         ENDIF
[524]210
211
212c-----------------------------------------------------------------------
213c   On initialise la pression et la fonction d'Exner :
214c   --------------------------------------------------
215
[1454]216      dq(:,:,:)=0.
[524]217      CALL pression ( ip1jmp1, ap, bp, ps, p       )
[1625]218      if (pressure_exner) then
[2021]219        CALL exner_hyb( ip1jmp1, ps, p, pks, pk, pkf )
[1625]220      else
[2021]221        CALL exner_milieu( ip1jmp1, ps, p, pks, pk, pkf )
[1520]222      endif
[524]223
224c-----------------------------------------------------------------------
225c   Debut de l'integration temporelle:
226c   ----------------------------------
227
[1614]228   1  CONTINUE ! Matsuno Forward step begins here
[524]229
[1577]230      jD_cur = jD_ref + day_ini - day_ref +                             &
[1626]231     &          itau/day_step
[1577]232      jH_cur = jH_ref + start_time +                                    &
[1626]233     &          mod(itau,day_step)/float(day_step)
[1577]234      jD_cur = jD_cur + int(jH_cur)
235      jH_cur = jH_cur - int(jH_cur)
[524]236
[1279]237
[524]238#ifdef CPP_IOIPSL
[1279]239      if (ok_guide) then
240        call guide_main(itau,ucov,vcov,teta,q,masse,ps)
[524]241      endif
242#endif
[1279]243
244
[524]245c
246c     IF( MOD( itau, 10* day_step ).EQ.0 )  THEN
247c       CALL  test_period ( ucov,vcov,teta,q,p,phis )
248c       PRINT *,' ----   Test_period apres continue   OK ! -----', itau
249c     ENDIF
250c
[1146]251
252! Save fields obtained at previous time step as '...m1'
[524]253      CALL SCOPY( ijmllm ,vcov , 1, vcovm1 , 1 )
254      CALL SCOPY( ijp1llm,ucov , 1, ucovm1 , 1 )
255      CALL SCOPY( ijp1llm,teta , 1, tetam1 , 1 )
256      CALL SCOPY( ijp1llm,masse, 1, massem1, 1 )
257      CALL SCOPY( ip1jmp1, ps  , 1,   psm1 , 1 )
258
259      forward = .TRUE.
260      leapf   = .FALSE.
261      dt      =  dtvr
262
263c   ...    P.Le Van .26/04/94  ....
[1616]264! Ehouarn: finvmaold is actually not used
265!      CALL SCOPY   ( ijp1llm,   masse, 1, finvmaold,     1 )
266!      CALL filtreg ( finvmaold ,jjp1, llm, -2,2, .TRUE., 1 )
[524]267
[1614]268   2  CONTINUE ! Matsuno backward or leapfrog step begins here
[524]269
270c-----------------------------------------------------------------------
271
272c   date:
273c   -----
274
275
276c   gestion des appels de la physique et des dissipations:
277c   ------------------------------------------------------
278c
279c   ...    P.Le Van  ( 6/02/95 )  ....
280
281      apphys = .FALSE.
282      statcl = .FALSE.
283      conser = .FALSE.
284      apdiss = .FALSE.
285
286      IF( purmats ) THEN
[1403]287      ! Purely Matsuno time stepping
[524]288         IF( MOD(itau,iconser) .EQ.0.AND.  forward    ) conser = .TRUE.
[1502]289         IF( MOD(itau,dissip_period ).EQ.0.AND..NOT.forward )
290     s        apdiss = .TRUE.
[524]291         IF( MOD(itau,iphysiq ).EQ.0.AND..NOT.forward
[1529]292     s          .and. physic                        ) apphys = .TRUE.
[524]293      ELSE
[1403]294      ! Leapfrog/Matsuno time stepping
[524]295         IF( MOD(itau   ,iconser) .EQ. 0              ) conser = .TRUE.
[1502]296         IF( MOD(itau+1,dissip_period).EQ.0 .AND. .NOT. forward )
297     s        apdiss = .TRUE.
[1529]298         IF( MOD(itau+1,iphysiq).EQ.0.AND.physic       ) apphys=.TRUE.
[524]299      END IF
300
[1403]301! Ehouarn: for Shallow Water case (ie: 1 vertical layer),
302!          supress dissipation step
303      if (llm.eq.1) then
304        apdiss=.false.
305      endif
306
[524]307c-----------------------------------------------------------------------
308c   calcul des tendances dynamiques:
309c   --------------------------------
310
[1614]311      ! compute geopotential phi()
[524]312      CALL geopot  ( ip1jmp1, teta  , pk , pks,  phis  , phi   )
313
[1279]314      time = jD_cur + jH_cur
[524]315      CALL caldyn
316     $  ( itau,ucov,vcov,teta,ps,masse,pk,pkf,phis ,
[1279]317     $    phi,conser,du,dv,dteta,dp,w, pbaru,pbarv, time )
[524]318
[1279]319
[524]320c-----------------------------------------------------------------------
321c   calcul des tendances advection des traceurs (dont l'humidite)
322c   -------------------------------------------------------------
323
324      IF( forward. OR . leapf )  THEN
[1987]325! Ehouarn: NB: fields sent to advtrac are those at the beginning of the time step
[960]326         CALL caladvtrac(q,pbaru,pbarv,
327     *        p, masse, dq,  teta,
328     .        flxw, pk)
329         
[524]330         IF (offline) THEN
331Cmaf stokage du flux de masse pour traceurs OFF-LINE
332
333#ifdef CPP_IOIPSL
[541]334           CALL fluxstokenc(pbaru,pbarv,masse,teta,phi,phis,
335     .   dtvr, itau)
[524]336#endif
337
338
[1146]339         ENDIF ! of IF (offline)
[524]340c
[1146]341      ENDIF ! of IF( forward. OR . leapf )
[524]342
343
344c-----------------------------------------------------------------------
345c   integrations dynamique et traceurs:
346c   ----------------------------------
347
348
349       CALL integrd ( 2,vcovm1,ucovm1,tetam1,psm1,massem1 ,
[1616]350     $         dv,du,dteta,dq,dp,vcov,ucov,teta,q,ps,masse,phis )
351!     $              finvmaold                                    )
[524]352
353
354c .P.Le Van (26/04/94  ajout de  finvpold dans l'appel d'integrd)
355c
356c-----------------------------------------------------------------------
357c   calcul des tendances physiques:
358c   -------------------------------
359c    ########   P.Le Van ( Modif le  6/02/95 )   ###########
360c
361       IF( purmats )  THEN
362          IF( itau.EQ.itaufin.AND..NOT.forward ) lafin = .TRUE.
363       ELSE
364          IF( itau+1. EQ. itaufin )              lafin = .TRUE.
365       ENDIF
366c
367c
368       IF( apphys )  THEN
369c
370c     .......   Ajout   P.Le Van ( 17/04/96 )   ...........
371c
372
373         CALL pression (  ip1jmp1, ap, bp, ps,  p      )
[1625]374         if (pressure_exner) then
[2021]375           CALL exner_hyb(  ip1jmp1, ps, p,pks, pk, pkf )
[1625]376         else
[2021]377           CALL exner_milieu( ip1jmp1, ps, p, pks, pk, pkf )
[1520]378         endif
[524]379
[2039]380! Appel a geopot ajoute le 2014/05/08 pour garantir la convergence numerique
381! avec dyn3dmem
382         CALL geopot  ( ip1jmp1, teta  , pk , pks,  phis  , phi   )
383
[1279]384!           rdaym_ini  = itau * dtvr / daysec
385!           rdayvrai   = rdaym_ini  + day_ini
[1577]386!           jD_cur = jD_ref + day_ini - day_ref
387!     $        + int (itau * dtvr / daysec)
388!           jH_cur = jH_ref +                                            &
389!     &              (itau * dtvr / daysec - int(itau * dtvr / daysec))
390           jD_cur = jD_ref + day_ini - day_ref +                        &
[1626]391     &          itau/day_step
[1676]392
393           IF (planet_type .eq."generic") THEN
394              ! AS: we make jD_cur to be pday
395              jD_cur = int(day_ini + itau/day_step)
396           ENDIF
397
[1577]398           jH_cur = jH_ref + start_time +                               &
[1626]399     &              mod(itau,day_step)/float(day_step)
[1577]400           jD_cur = jD_cur + int(jH_cur)
401           jH_cur = jH_cur - int(jH_cur)
[1279]402!         write(lunout,*)'itau, jD_cur = ', itau, jD_cur, jH_cur
403!         call ju2ymds(jD_cur+jH_cur, an, mois, jour, secondes)
404!         write(lunout,*)'current date = ',an, mois, jour, secondes
[524]405
406c rajout debug
407c       lafin = .true.
408
409
410c   Inbterface avec les routines de phylmd (phymars ... )
411c   -----------------------------------------------------
412
413c+jld
414
415c  Diagnostique de conservation de l'énergie : initialisation
[1146]416         IF (ip_ebil_dyn.ge.1 ) THEN
[524]417          ztit='bil dyn'
[1146]418! Ehouarn: be careful, diagedyn is Earth-specific (includes ../phylmd/..)!
419           IF (planet_type.eq."earth") THEN
[1615]420#ifdef CPP_EARTH
[1146]421            CALL diagedyn(ztit,2,1,1,dtphys
422     &    , ucov    , vcov , ps, p ,pk , teta , q(:,:,1), q(:,:,2))
[1615]423#endif
[1146]424           ENDIF
425         ENDIF ! of IF (ip_ebil_dyn.ge.1 )
[524]426c-jld
[1146]427#ifdef CPP_IOIPSL
[1656]428cIM decommenter les 6 lignes suivantes pour sortir quelques parametres dynamiques de LMDZ
429cIM uncomment next 6 lines to get some parameters for LMDZ dynamics
430c        IF (first) THEN
431c         first=.false.
432c#include "ini_paramLMDZ_dyn.h"
433c        ENDIF
[692]434c
[1656]435c#include "write_paramLMDZ_dyn.h"
[692]436c
[1146]437#endif
438! #endif of #ifdef CPP_IOIPSL
[2239]439#ifdef CPP_PHYS
[1279]440         CALL calfis( lafin , jD_cur, jH_cur,
[524]441     $               ucov,vcov,teta,q,masse,ps,p,pk,phis,phi ,
[960]442     $               du,dv,dteta,dq,
[2221]443     $               flxw,dufi,dvfi,dtetafi,dqfi,dpfi  )
[2239]444#endif
[524]445c      ajout des tendances physiques:
446c      ------------------------------
[1146]447          CALL addfi( dtphys, leapf, forward   ,
[524]448     $                  ucov, vcov, teta , q   ,ps ,
449     $                 dufi, dvfi, dtetafi , dqfi ,dpfi  )
[1987]450          ! since addfi updates ps(), also update p(), masse() and pk()
451          CALL pression (ip1jmp1,ap,bp,ps,p)
452          CALL massdair(p,masse)
453          if (pressure_exner) then
[2021]454            CALL exner_hyb(ip1jmp1,ps,p,pks,pk,pkf)
[1987]455          else
[2021]456            CALL exner_milieu(ip1jmp1,ps,p,pks,pk,pkf)
[1987]457          endif
[1793]458
459         IF (ok_strato) THEN
460           CALL top_bound( vcov,ucov,teta,masse,dtphys)
461         ENDIF
462       
[524]463c
464c  Diagnostique de conservation de l'énergie : difference
[1146]465         IF (ip_ebil_dyn.ge.1 ) THEN
[524]466          ztit='bil phys'
[1146]467          IF (planet_type.eq."earth") THEN
468           CALL diagedyn(ztit,2,1,1,dtphys
469     &     , ucov    , vcov , ps, p ,pk , teta , q(:,:,1), q(:,:,2))
470          ENDIF
471         ENDIF ! of IF (ip_ebil_dyn.ge.1 )
472
[1060]473       ENDIF ! of IF( apphys )
[524]474
[1146]475      IF(iflag_phys.EQ.2) THEN ! "Newtonian" case
[1454]476!   Academic case : Simple friction and Newtonan relaxation
477!   -------------------------------------------------------
478        DO l=1,llm   
479          DO ij=1,ip1jmp1
480           teta(ij,l)=teta(ij,l)-dtvr*
481     &      (teta(ij,l)-tetarappel(ij,l))*(knewt_g+knewt_t(l)*clat4(ij))
482          ENDDO
483        ENDDO ! of DO l=1,llm
484       
[1505]485        if (planet_type.eq."giant") then
486          ! add an intrinsic heat flux at the base of the atmosphere
487          teta(:,1)=teta(:,1)+dtvr*aire(:)*ihf/cpp/masse(:,1)
488        endif
489
490        call friction(ucov,vcov,dtvr)
491       
[1454]492        ! Sponge layer (if any)
493        IF (ok_strato) THEN
[1793]494!          dufi(:,:)=0.
495!          dvfi(:,:)=0.
496!          dtetafi(:,:)=0.
497!          dqfi(:,:,:)=0.
498!          dpfi(:)=0.
499!          CALL top_bound(vcov,ucov,teta,masse,dufi,dvfi,dtetafi)
500           CALL top_bound( vcov,ucov,teta,masse,dtvr)
501!          CALL addfi( dtvr, leapf, forward   ,
502!     $                  ucov, vcov, teta , q   ,ps ,
503!     $                 dufi, dvfi, dtetafi , dqfi ,dpfi  )
[1454]504        ENDIF ! of IF (ok_strato)
505      ENDIF ! of IF (iflag_phys.EQ.2)
[524]506
507
508c-jld
509
510        CALL pression ( ip1jmp1, ap, bp, ps, p                  )
[1625]511        if (pressure_exner) then
[2021]512          CALL exner_hyb( ip1jmp1, ps, p, pks, pk, pkf )
[1625]513        else
[2021]514          CALL exner_milieu( ip1jmp1, ps, p, pks, pk, pkf )
[1520]515        endif
[1987]516        CALL massdair(p,masse)
[524]517
518
519c-----------------------------------------------------------------------
520c   dissipation horizontale et verticale  des petites echelles:
521c   ----------------------------------------------------------
522
523      IF(apdiss) THEN
524
525
526c   calcul de l'energie cinetique avant dissipation
527        call covcont(llm,ucov,vcov,ucont,vcont)
528        call enercin(vcov,ucov,vcont,ucont,ecin0)
529
530c   dissipation
531        CALL dissip(vcov,ucov,teta,p,dvdis,dudis,dtetadis)
532        ucov=ucov+dudis
533        vcov=vcov+dvdis
534c       teta=teta+dtetadis
535
536
537c------------------------------------------------------------------------
538        if (dissip_conservative) then
539C       On rajoute la tendance due a la transform. Ec -> E therm. cree
540C       lors de la dissipation
541            call covcont(llm,ucov,vcov,ucont,vcont)
542            call enercin(vcov,ucov,vcont,ucont,ecin)
543            dtetaecdt= (ecin0-ecin)/ pk
544c           teta=teta+dtetaecdt
545            dtetadis=dtetadis+dtetaecdt
546        endif
547        teta=teta+dtetadis
548c------------------------------------------------------------------------
549
550
551c    .......        P. Le Van (  ajout  le 17/04/96  )   ...........
552c   ...      Calcul de la valeur moyenne, unique de h aux poles  .....
553c
554
555        DO l  =  1, llm
556          DO ij =  1,iim
557           tppn(ij)  = aire(  ij    ) * teta(  ij    ,l)
558           tpps(ij)  = aire(ij+ip1jm) * teta(ij+ip1jm,l)
559          ENDDO
560           tpn  = SSUM(iim,tppn,1)/apoln
561           tps  = SSUM(iim,tpps,1)/apols
562
563          DO ij = 1, iip1
564           teta(  ij    ,l) = tpn
565           teta(ij+ip1jm,l) = tps
566          ENDDO
567        ENDDO
568
[1614]569        if (1 == 0) then
570!!! Ehouarn: lines here 1) kill 1+1=2 in the dynamics
571!!!                     2) should probably not be here anyway
572!!! but are kept for those who would want to revert to previous behaviour
573           DO ij =  1,iim
574             tppn(ij)  = aire(  ij    ) * ps (  ij    )
575             tpps(ij)  = aire(ij+ip1jm) * ps (ij+ip1jm)
576           ENDDO
577             tpn  = SSUM(iim,tppn,1)/apoln
578             tps  = SSUM(iim,tpps,1)/apols
[524]579
[1614]580           DO ij = 1, iip1
581             ps(  ij    ) = tpn
582             ps(ij+ip1jm) = tps
583           ENDDO
584        endif ! of if (1 == 0)
[524]585
[1146]586      END IF ! of IF(apdiss)
[524]587
588c ajout debug
589c              IF( lafin ) then 
590c                abort_message = 'Simulation finished'
591c                call abort_gcm(modname,abort_message,0)
592c              ENDIF
593       
594c   ********************************************************************
595c   ********************************************************************
596c   .... fin de l'integration dynamique  et physique pour le pas itau ..
597c   ********************************************************************
598c   ********************************************************************
599
600c   preparation du pas d'integration suivant  ......
601
602      IF ( .NOT.purmats ) THEN
603c       ........................................................
604c       ..............  schema matsuno + leapfrog  ..............
605c       ........................................................
606
607            IF(forward. OR. leapf) THEN
608              itau= itau + 1
[1403]609c              iday= day_ini+itau/day_step
610c              time= REAL(itau-(iday-day_ini)*day_step)/day_step+time_0
611c                IF(time.GT.1.) THEN
612c                  time = time-1.
613c                  iday = iday+1
614c                ENDIF
[524]615            ENDIF
616
617
618            IF( itau. EQ. itaufinp1 ) then 
[999]619              if (flag_verif) then
[1279]620                write(79,*) 'ucov',ucov
621                write(80,*) 'vcov',vcov
622                write(81,*) 'teta',teta
623                write(82,*) 'ps',ps
624                write(83,*) 'q',q
[999]625                WRITE(85,*) 'q1 = ',q(:,:,1)
626                WRITE(86,*) 'q3 = ',q(:,:,3)
627              endif
[524]628
629              abort_message = 'Simulation finished'
630
631              call abort_gcm(modname,abort_message,0)
632            ENDIF
633c-----------------------------------------------------------------------
634c   ecriture du fichier histoire moyenne:
635c   -------------------------------------
636
637            IF(MOD(itau,iperiod).EQ.0 .OR. itau.EQ.itaufin) THEN
638               IF(itau.EQ.itaufin) THEN
639                  iav=1
640               ELSE
641                  iav=0
642               ENDIF
[1146]643               
644               IF (ok_dynzon) THEN
[524]645#ifdef CPP_IOIPSL
[1403]646                 CALL bilan_dyn(2,dtvr*iperiod,dtvr*day_step*periodav,
647     &                 ps,masse,pk,pbaru,pbarv,teta,phi,ucov,vcov,q)
[524]648#endif
[1146]649               END IF
[1403]650               IF (ok_dyn_ave) THEN
651#ifdef CPP_IOIPSL
652                 CALL writedynav(itau,vcov,
653     &                 ucov,teta,pk,phi,q,masse,ps,phis)
654#endif
655               ENDIF
[524]656
[1403]657            ENDIF ! of IF((MOD(itau,iperiod).EQ.0).OR.(itau.EQ.itaufin))
[524]658
659c-----------------------------------------------------------------------
660c   ecriture de la bande histoire:
661c   ------------------------------
662
[1403]663            IF( MOD(itau,iecri).EQ.0) THEN
664             ! Ehouarn: output only during LF or Backward Matsuno
665             if (leapf.or.(.not.leapf.and.(.not.forward))) then
[1146]666              CALL geopot(ip1jmp1,teta,pk,pks,phis,phi)
667              unat=0.
668              do l=1,llm
669                unat(iip2:ip1jm,l)=ucov(iip2:ip1jm,l)/cu(iip2:ip1jm)
670                vnat(:,l)=vcov(:,l)/cv(:)
671              enddo
[524]672#ifdef CPP_IOIPSL
[1403]673              if (ok_dyn_ins) then
674!               write(lunout,*) "leapfrog: call writehist, itau=",itau
675               CALL writehist(itau,vcov,ucov,teta,phi,q,masse,ps,phis)
676!               call WriteField('ucov',reshape(ucov,(/iip1,jmp1,llm/)))
677!               call WriteField('vcov',reshape(vcov,(/iip1,jjm,llm/)))
678!              call WriteField('teta',reshape(teta,(/iip1,jmp1,llm/)))
679!               call WriteField('ps',reshape(ps,(/iip1,jmp1/)))
680!               call WriteField('masse',reshape(masse,(/iip1,jmp1,llm/)))
681              endif ! of if (ok_dyn_ins)
[1146]682#endif
683! For some Grads outputs of fields
[1403]684              if (output_grads_dyn) then
[524]685#include "write_grads_dyn.h"
[1403]686              endif
687             endif ! of if (leapf.or.(.not.leapf.and.(.not.forward)))
[1146]688            ENDIF ! of IF(MOD(itau,iecri).EQ.0)
[524]689
690            IF(itau.EQ.itaufin) THEN
691
692
[1454]693!              if (planet_type.eq."earth") then
[1146]694! Write an Earth-format restart file
[1577]695                CALL dynredem1("restart.nc",start_time,
[1146]696     &                         vcov,ucov,teta,q,masse,ps)
[1454]697!              endif ! of if (planet_type.eq."earth")
[524]698
699              CLOSE(99)
[1614]700              !!! Ehouarn: Why not stop here and now?
[1146]701            ENDIF ! of IF (itau.EQ.itaufin)
[524]702
703c-----------------------------------------------------------------------
704c   gestion de l'integration temporelle:
705c   ------------------------------------
706
707            IF( MOD(itau,iperiod).EQ.0 )    THEN
708                    GO TO 1
709            ELSE IF ( MOD(itau-1,iperiod). EQ. 0 ) THEN
710
711                   IF( forward )  THEN
712c      fin du pas forward et debut du pas backward
713
714                      forward = .FALSE.
715                        leapf = .FALSE.
716                           GO TO 2
717
718                   ELSE
719c      fin du pas backward et debut du premier pas leapfrog
720
721                        leapf =  .TRUE.
722                        dt  =  2.*dtvr
[1146]723                        GO TO 2
724                   END IF ! of IF (forward)
[524]725            ELSE
726
727c      ......   pas leapfrog  .....
728
729                 leapf = .TRUE.
730                 dt  = 2.*dtvr
731                 GO TO 2
[1146]732            END IF ! of IF (MOD(itau,iperiod).EQ.0)
733                   !    ELSEIF (MOD(itau-1,iperiod).EQ.0)
[524]734
[1146]735      ELSE ! of IF (.not.purmats)
[524]736
737c       ........................................................
738c       ..............       schema  matsuno        ...............
739c       ........................................................
740            IF( forward )  THEN
741
742             itau =  itau + 1
[1403]743c             iday = day_ini+itau/day_step
744c             time = REAL(itau-(iday-day_ini)*day_step)/day_step+time_0
745c
746c                  IF(time.GT.1.) THEN
747c                   time = time-1.
748c                   iday = iday+1
749c                  ENDIF
[524]750
751               forward =  .FALSE.
752               IF( itau. EQ. itaufinp1 ) then 
753                 abort_message = 'Simulation finished'
754                 call abort_gcm(modname,abort_message,0)
755               ENDIF
756               GO TO 2
757
[1403]758            ELSE ! of IF(forward) i.e. backward step
[524]759
[1146]760              IF(MOD(itau,iperiod).EQ.0 .OR. itau.EQ.itaufin) THEN
[524]761               IF(itau.EQ.itaufin) THEN
762                  iav=1
763               ELSE
764                  iav=0
765               ENDIF
[1146]766
767               IF (ok_dynzon) THEN
[524]768#ifdef CPP_IOIPSL
[1403]769                 CALL bilan_dyn(2,dtvr*iperiod,dtvr*day_step*periodav,
770     &                 ps,masse,pk,pbaru,pbarv,teta,phi,ucov,vcov,q)
[524]771#endif
[1403]772               ENDIF
773               IF (ok_dyn_ave) THEN
774#ifdef CPP_IOIPSL
775                 CALL writedynav(itau,vcov,
776     &                 ucov,teta,pk,phi,q,masse,ps,phis)
777#endif
778               ENDIF
[524]779
[1146]780              ENDIF ! of IF(MOD(itau,iperiod).EQ.0 .OR. itau.EQ.itaufin)
[524]781
[1146]782              IF(MOD(itau,iecri         ).EQ.0) THEN
[524]783c              IF(MOD(itau,iecri*day_step).EQ.0) THEN
[1146]784                CALL geopot(ip1jmp1,teta,pk,pks,phis,phi)
785                unat=0.
786                do l=1,llm
787                  unat(iip2:ip1jm,l)=ucov(iip2:ip1jm,l)/cu(iip2:ip1jm)
788                  vnat(:,l)=vcov(:,l)/cv(:)
789                enddo
[524]790#ifdef CPP_IOIPSL
[1403]791              if (ok_dyn_ins) then
792!                write(lunout,*) "leapfrog: call writehist (b)",
793!     &                        itau,iecri
794                CALL writehist(itau,vcov,ucov,teta,phi,q,masse,ps,phis)
795              endif ! of if (ok_dyn_ins)
[1146]796#endif
797! For some Grads outputs
798                if (output_grads_dyn) then
[524]799#include "write_grads_dyn.h"
[1146]800                endif
[524]801
[1146]802              ENDIF ! of IF(MOD(itau,iecri         ).EQ.0)
[524]803
[1146]804              IF(itau.EQ.itaufin) THEN
[1454]805!                if (planet_type.eq."earth") then
[1577]806                  CALL dynredem1("restart.nc",start_time,
[1146]807     &                           vcov,ucov,teta,q,masse,ps)
[1454]808!                endif ! of if (planet_type.eq."earth")
[1146]809              ENDIF ! of IF(itau.EQ.itaufin)
[524]810
[1146]811              forward = .TRUE.
812              GO TO  1
[524]813
[1146]814            ENDIF ! of IF (forward)
[524]815
[1146]816      END IF ! of IF(.not.purmats)
[524]817
818      STOP
819      END
Note: See TracBrowser for help on using the repository browser.