source: LMDZ4/branches/LMDZ4-dev/libf/dyn3d/leapfrog.F @ 1143

Last change on this file since 1143 was 1143, checked in by jghattas, 15 years ago

Recuperation des developpements fait uniquement sur la branche LMDZ4_V4_patches :

  • ajoute de la nouvelle flag ok_dynzon
  • ajoute du parametre aer_type
  • optimisation : isccp_cloud_types.F

+ bug pour le slab dans conf_phys.F90

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 21.4 KB
Line 
1!
2c
3c
4      SUBROUTINE leapfrog(ucov,vcov,teta,ps,masse,phis,q,clesphy0,
5     &                    time_0)
6
7
8cIM : pour sortir les param. du modele dans un fis. netcdf 110106
9#ifdef CPP_IOIPSL
10      use IOIPSL
11#endif
12      USE infotrac
13
14      IMPLICIT NONE
15
16c      ......   Version  du 10/01/98    ..........
17
18c             avec  coordonnees  verticales hybrides
19c   avec nouveaux operat. dissipation * ( gradiv2,divgrad2,nxgraro2 )
20
21c=======================================================================
22c
23c   Auteur:  P. Le Van /L. Fairhead/F.Hourdin
24c   -------
25c
26c   Objet:
27c   ------
28c
29c   GCM LMD nouvelle grille
30c
31c=======================================================================
32c
33c  ... Dans inigeom , nouveaux calculs pour les elongations  cu , cv
34c      et possibilite d'appeler une fonction f(y)  a derivee tangente
35c      hyperbolique a la  place de la fonction a derivee sinusoidale.
36
37c  ... Possibilite de choisir le shema pour l'advection de
38c        q  , en modifiant iadv dans traceur.def  (10/02) .
39c
40c      Pour Van-Leer + Vapeur d'eau saturee, iadv(1)=4. (F.Codron,10/99)
41c      Pour Van-Leer iadv=10
42c
43c-----------------------------------------------------------------------
44c   Declarations:
45c   -------------
46
47#include "dimensions.h"
48#include "paramet.h"
49#include "comconst.h"
50#include "comdissnew.h"
51#include "comvert.h"
52#include "comgeom.h"
53#include "logic.h"
54#include "temps.h"
55#include "control.h"
56#include "ener.h"
57#include "description.h"
58#include "serre.h"
59#include "com_io_dyn.h"
60#include "iniprint.h"
61#include "academic.h"
62
63! FH 2008/05/09 On elimine toutes les clefs physiques dans la dynamique
64! #include "clesphys.h"
65
66      INTEGER         longcles
67      PARAMETER     ( longcles = 20 )
68      REAL  clesphy0( longcles )
69
70      real zqmin,zqmax
71      INTEGER nbetatmoy, nbetatdem,nbetat
72
73c   variables dynamiques
74      REAL vcov(ip1jm,llm),ucov(ip1jmp1,llm) ! vents covariants
75      REAL teta(ip1jmp1,llm)                 ! temperature potentielle
76      REAL q(ip1jmp1,llm,nqtot)               ! champs advectes
77      REAL ps(ip1jmp1)                       ! pression  au sol
78      REAL p (ip1jmp1,llmp1  )               ! pression aux interfac.des couches
79      REAL pks(ip1jmp1)                      ! exner au  sol
80      REAL pk(ip1jmp1,llm)                   ! exner au milieu des couches
81      REAL pkf(ip1jmp1,llm)                  ! exner filt.au milieu des couches
82      REAL masse(ip1jmp1,llm)                ! masse d'air
83      REAL phis(ip1jmp1)                     ! geopotentiel au sol
84      REAL phi(ip1jmp1,llm)                  ! geopotentiel
85      REAL w(ip1jmp1,llm)                    ! vitesse verticale
86
87c variables dynamiques intermediaire pour le transport
88      REAL pbaru(ip1jmp1,llm),pbarv(ip1jm,llm) !flux de masse
89
90c   variables dynamiques au pas -1
91      REAL vcovm1(ip1jm,llm),ucovm1(ip1jmp1,llm)
92      REAL tetam1(ip1jmp1,llm),psm1(ip1jmp1)
93      REAL massem1(ip1jmp1,llm)
94
95c   tendances dynamiques
96      REAL dv(ip1jm,llm),du(ip1jmp1,llm)
97      REAL dteta(ip1jmp1,llm),dq(ip1jmp1,llm,nqtot),dp(ip1jmp1)
98
99c   tendances de la dissipation
100      REAL dvdis(ip1jm,llm),dudis(ip1jmp1,llm)
101      REAL dtetadis(ip1jmp1,llm)
102
103c   tendances physiques
104      REAL dvfi(ip1jm,llm),dufi(ip1jmp1,llm)
105      REAL dtetafi(ip1jmp1,llm),dqfi(ip1jmp1,llm,nqtot),dpfi(ip1jmp1)
106
107c   variables pour le fichier histoire
108      REAL dtav      ! intervalle de temps elementaire
109
110      REAL tppn(iim),tpps(iim),tpn,tps
111c
112      INTEGER itau,itaufinp1,iav
113      INTEGER*4  iday ! jour julien
114      REAL       time ! Heure de la journee en fraction d'1 jour
115
116      REAL  SSUM
117      REAL time_0 , finvmaold(ip1jmp1,llm)
118
119cym      LOGICAL  lafin
120      LOGICAL :: lafin=.false.
121      INTEGER ij,iq,l
122      INTEGER ik
123
124      real time_step, t_wrt, t_ops
125
126      REAL rdayvrai,rdaym_ini
127      LOGICAL first,callinigrads
128cIM : pour sortir les param. du modele dans un fis. netcdf 110106
129      save first
130      data first/.true./
131      real dt_cum, zjulian
132      character*10 infile
133      integer zan, tau0, thoriid
134      integer nid_ctesGCM
135      save nid_ctesGCM
136      real degres
137      real rlong(iip1), rlatg(jjp1)
138      real zx_tmp_2d(iip1,jjp1)
139      integer ndex2d(iip1*jjp1)
140      logical ok_sync
141      parameter (ok_sync = .true.)
142
143      data callinigrads/.true./
144      character*10 string10
145
146      REAL alpha(ip1jmp1,llm),beta(ip1jmp1,llm)
147      REAL :: flxw(ip1jmp1,llm)  ! flux de masse verticale
148
149c+jld variables test conservation energie
150      REAL ecin(ip1jmp1,llm),ecin0(ip1jmp1,llm)
151C     Tendance de la temp. potentiel d (theta)/ d t due a la
152C     tansformation d'energie cinetique en energie thermique
153C     cree par la dissipation
154      REAL dtetaecdt(ip1jmp1,llm)
155      REAL vcont(ip1jm,llm),ucont(ip1jmp1,llm)
156      REAL vnat(ip1jm,llm),unat(ip1jmp1,llm)
157      REAL      d_h_vcol, d_qt, d_qw, d_ql, d_ec
158      CHARACTER*15 ztit
159!IM   INTEGER   ip_ebil_dyn  ! PRINT level for energy conserv. diag.
160!IM   SAVE      ip_ebil_dyn
161!IM   DATA      ip_ebil_dyn/0/
162c-jld
163
164      character*80 dynhist_file, dynhistave_file
165      character(len=20) :: modname
166      character*80 abort_message
167
168C Calendrier
169      LOGICAL true_calendar
170      PARAMETER (true_calendar = .false.)
171
172      logical dissip_conservative
173      save dissip_conservative
174      data dissip_conservative/.true./
175
176      LOGICAL prem
177      save prem
178      DATA prem/.true./
179      INTEGER testita
180      PARAMETER (testita = 9)
181
182      logical , parameter :: flag_verif = .false.
183     
184
185      integer itau_w   ! pas de temps ecriture = itap + itau_phy
186
187
188      itaufin   = nday*day_step
189      itaufinp1 = itaufin +1
190      modname="leapfrog"
191     
192
193      itau = 0
194      iday = day_ini+itau/day_step
195      time = FLOAT(itau-(iday-day_ini)*day_step)/day_step+time_0
196         IF(time.GT.1.) THEN
197          time = time-1.
198          iday = iday+1
199         ENDIF
200
201
202c-----------------------------------------------------------------------
203c   On initialise la pression et la fonction d'Exner :
204c   --------------------------------------------------
205
206      dq=0.
207      CALL pression ( ip1jmp1, ap, bp, ps, p       )
208      CALL exner_hyb( ip1jmp1, ps, p,alpha,beta, pks, pk, pkf )
209
210c-----------------------------------------------------------------------
211c   Debut de l'integration temporelle:
212c   ----------------------------------
213
214   1  CONTINUE
215
216
217#ifdef CPP_IOIPSL
218      if (ok_guide.and.(itaufin-itau-1)*dtvr.gt.21600) then
219        call guide(itau,ucov,vcov,teta,q,masse,ps)
220      else
221        IF(prt_level>9)WRITE(lunout,*)'leapfrog: attention on ne ',
222     .    'guide pas les 6 dernieres heures'
223      endif
224#endif
225c
226c     IF( MOD( itau, 10* day_step ).EQ.0 )  THEN
227c       CALL  test_period ( ucov,vcov,teta,q,p,phis )
228c       PRINT *,' ----   Test_period apres continue   OK ! -----', itau
229c     ENDIF
230c
231
232! Save fields obtained at previous time step as '...m1'
233      CALL SCOPY( ijmllm ,vcov , 1, vcovm1 , 1 )
234      CALL SCOPY( ijp1llm,ucov , 1, ucovm1 , 1 )
235      CALL SCOPY( ijp1llm,teta , 1, tetam1 , 1 )
236      CALL SCOPY( ijp1llm,masse, 1, massem1, 1 )
237      CALL SCOPY( ip1jmp1, ps  , 1,   psm1 , 1 )
238
239      forward = .TRUE.
240      leapf   = .FALSE.
241      dt      =  dtvr
242
243c   ...    P.Le Van .26/04/94  ....
244
245      CALL SCOPY   ( ijp1llm,   masse, 1, finvmaold,     1 )
246      CALL filtreg ( finvmaold ,jjp1, llm, -2,2, .TRUE., 1 )
247
248! Ehouarn: what is this for? zqmin & zqmax are not used anyway ...
249!      call minmax(ijp1llm,q(:,:,3),zqmin,zqmax)
250
251   2  CONTINUE
252
253c-----------------------------------------------------------------------
254
255c   date:
256c   -----
257
258
259c   gestion des appels de la physique et des dissipations:
260c   ------------------------------------------------------
261c
262c   ...    P.Le Van  ( 6/02/95 )  ....
263
264      apphys = .FALSE.
265      statcl = .FALSE.
266      conser = .FALSE.
267      apdiss = .FALSE.
268
269      IF( purmats ) THEN
270         IF( MOD(itau,iconser) .EQ.0.AND.  forward    ) conser = .TRUE.
271         IF( MOD(itau,idissip ).EQ.0.AND..NOT.forward ) apdiss = .TRUE.
272         IF( MOD(itau,iphysiq ).EQ.0.AND..NOT.forward
273     s          .and. iflag_phys.EQ.1                 ) apphys = .TRUE.
274      ELSE
275         IF( MOD(itau   ,iconser) .EQ. 0              ) conser = .TRUE.
276         IF( MOD(itau+1,idissip)  .EQ. 0              ) apdiss = .TRUE.
277         IF( MOD(itau+1,iphysiq).EQ.0.AND.iflag_phys.EQ.1) apphys=.TRUE.
278      END IF
279
280c-----------------------------------------------------------------------
281c   calcul des tendances dynamiques:
282c   --------------------------------
283
284      CALL geopot  ( ip1jmp1, teta  , pk , pks,  phis  , phi   )
285
286      CALL caldyn
287     $  ( itau,ucov,vcov,teta,ps,masse,pk,pkf,phis ,
288     $    phi,conser,du,dv,dteta,dp,w, pbaru,pbarv, time+iday-day_ini )
289
290c-----------------------------------------------------------------------
291c   calcul des tendances advection des traceurs (dont l'humidite)
292c   -------------------------------------------------------------
293
294      IF( forward. OR . leapf )  THEN
295
296         CALL caladvtrac(q,pbaru,pbarv,
297     *        p, masse, dq,  teta,
298     .        flxw, pk)
299         
300         IF (offline) THEN
301Cmaf stokage du flux de masse pour traceurs OFF-LINE
302
303#ifdef CPP_IOIPSL
304           CALL fluxstokenc(pbaru,pbarv,masse,teta,phi,phis,
305     .   dtvr, itau)
306#endif
307
308
309         ENDIF ! of IF (offline)
310c
311      ENDIF ! of IF( forward. OR . leapf )
312
313
314c-----------------------------------------------------------------------
315c   integrations dynamique et traceurs:
316c   ----------------------------------
317
318
319       CALL integrd ( 2,vcovm1,ucovm1,tetam1,psm1,massem1 ,
320     $         dv,du,dteta,dq,dp,vcov,ucov,teta,q,ps,masse,phis ,
321     $              finvmaold                                    )
322
323
324c .P.Le Van (26/04/94  ajout de  finvpold dans l'appel d'integrd)
325c
326c-----------------------------------------------------------------------
327c   calcul des tendances physiques:
328c   -------------------------------
329c    ########   P.Le Van ( Modif le  6/02/95 )   ###########
330c
331       IF( purmats )  THEN
332          IF( itau.EQ.itaufin.AND..NOT.forward ) lafin = .TRUE.
333       ELSE
334          IF( itau+1. EQ. itaufin )              lafin = .TRUE.
335       ENDIF
336c
337c
338       IF( apphys )  THEN
339c
340c     .......   Ajout   P.Le Van ( 17/04/96 )   ...........
341c
342
343         CALL pression (  ip1jmp1, ap, bp, ps,  p      )
344         CALL exner_hyb(  ip1jmp1, ps, p,alpha,beta,pks, pk, pkf )
345
346           rdaym_ini  = itau * dtvr / daysec
347           rdayvrai   = rdaym_ini  + day_ini
348
349
350c rajout debug
351c       lafin = .true.
352
353
354c   Inbterface avec les routines de phylmd (phymars ... )
355c   -----------------------------------------------------
356
357c+jld
358
359c  Diagnostique de conservation de l'énergie : initialisation
360         IF (ip_ebil_dyn.ge.1 ) THEN
361          ztit='bil dyn'
362! Ehouarn: be careful, diagedyn is Earth-specific (includes ../phylmd/..)!
363           IF (planet_type.eq."earth") THEN
364            CALL diagedyn(ztit,2,1,1,dtphys
365     &    , ucov    , vcov , ps, p ,pk , teta , q(:,:,1), q(:,:,2))
366           ENDIF
367         ENDIF ! of IF (ip_ebil_dyn.ge.1 )
368c-jld
369#ifdef CPP_IOIPSL
370cIM : pour sortir les param. du modele dans un fis. netcdf 110106
371         IF (first) THEN
372          first=.false.
373#include "ini_paramLMDZ_dyn.h"
374         ENDIF
375c
376#include "write_paramLMDZ_dyn.h"
377c
378#endif
379! #endif of #ifdef CPP_IOIPSL
380         CALL calfis( lafin ,rdayvrai,time  ,
381     $               ucov,vcov,teta,q,masse,ps,p,pk,phis,phi ,
382     $               du,dv,dteta,dq,
383     $               flxw,
384     $               clesphy0, dufi,dvfi,dtetafi,dqfi,dpfi  )
385
386         IF (ok_strato) THEN
387           CALL top_bound( vcov,ucov,teta, dufi,dvfi,dtetafi)
388         ENDIF
389       
390c      ajout des tendances physiques:
391c      ------------------------------
392          CALL addfi( dtphys, leapf, forward   ,
393     $                  ucov, vcov, teta , q   ,ps ,
394     $                 dufi, dvfi, dtetafi , dqfi ,dpfi  )
395c
396c  Diagnostique de conservation de l'énergie : difference
397         IF (ip_ebil_dyn.ge.1 ) THEN
398          ztit='bil phys'
399          IF (planet_type.eq."earth") THEN
400           CALL diagedyn(ztit,2,1,1,dtphys
401     &     , ucov    , vcov , ps, p ,pk , teta , q(:,:,1), q(:,:,2))
402          ENDIF
403         ENDIF ! of IF (ip_ebil_dyn.ge.1 )
404
405       ENDIF ! of IF( apphys )
406
407      IF(iflag_phys.EQ.2) THEN ! "Newtonian" case
408c   Calcul academique de la physique = Rappel Newtonien + friction
409c   --------------------------------------------------------------
410       teta(:,:)=teta(:,:)
411     s  -iphysiq*dtvr*(teta(:,:)-tetarappel(:,:))/taurappel
412       call friction(ucov,vcov,iphysiq*dtvr)
413      ENDIF
414
415
416c-jld
417
418        CALL pression ( ip1jmp1, ap, bp, ps, p                  )
419        CALL exner_hyb( ip1jmp1, ps, p,alpha,beta, pks, pk, pkf )
420
421
422c-----------------------------------------------------------------------
423c   dissipation horizontale et verticale  des petites echelles:
424c   ----------------------------------------------------------
425
426      IF(apdiss) THEN
427
428
429c   calcul de l'energie cinetique avant dissipation
430        call covcont(llm,ucov,vcov,ucont,vcont)
431        call enercin(vcov,ucov,vcont,ucont,ecin0)
432
433c   dissipation
434        CALL dissip(vcov,ucov,teta,p,dvdis,dudis,dtetadis)
435        ucov=ucov+dudis
436        vcov=vcov+dvdis
437c       teta=teta+dtetadis
438
439
440c------------------------------------------------------------------------
441        if (dissip_conservative) then
442C       On rajoute la tendance due a la transform. Ec -> E therm. cree
443C       lors de la dissipation
444            call covcont(llm,ucov,vcov,ucont,vcont)
445            call enercin(vcov,ucov,vcont,ucont,ecin)
446            dtetaecdt= (ecin0-ecin)/ pk
447c           teta=teta+dtetaecdt
448            dtetadis=dtetadis+dtetaecdt
449        endif
450        teta=teta+dtetadis
451c------------------------------------------------------------------------
452
453
454c    .......        P. Le Van (  ajout  le 17/04/96  )   ...........
455c   ...      Calcul de la valeur moyenne, unique de h aux poles  .....
456c
457
458        DO l  =  1, llm
459          DO ij =  1,iim
460           tppn(ij)  = aire(  ij    ) * teta(  ij    ,l)
461           tpps(ij)  = aire(ij+ip1jm) * teta(ij+ip1jm,l)
462          ENDDO
463           tpn  = SSUM(iim,tppn,1)/apoln
464           tps  = SSUM(iim,tpps,1)/apols
465
466          DO ij = 1, iip1
467           teta(  ij    ,l) = tpn
468           teta(ij+ip1jm,l) = tps
469          ENDDO
470        ENDDO
471
472        DO ij =  1,iim
473          tppn(ij)  = aire(  ij    ) * ps (  ij    )
474          tpps(ij)  = aire(ij+ip1jm) * ps (ij+ip1jm)
475        ENDDO
476          tpn  = SSUM(iim,tppn,1)/apoln
477          tps  = SSUM(iim,tpps,1)/apols
478
479        DO ij = 1, iip1
480          ps(  ij    ) = tpn
481          ps(ij+ip1jm) = tps
482        ENDDO
483
484
485      END IF ! of IF(apdiss)
486
487c ajout debug
488c              IF( lafin ) then 
489c                abort_message = 'Simulation finished'
490c                call abort_gcm(modname,abort_message,0)
491c              ENDIF
492       
493c   ********************************************************************
494c   ********************************************************************
495c   .... fin de l'integration dynamique  et physique pour le pas itau ..
496c   ********************************************************************
497c   ********************************************************************
498
499c   preparation du pas d'integration suivant  ......
500
501      IF ( .NOT.purmats ) THEN
502c       ........................................................
503c       ..............  schema matsuno + leapfrog  ..............
504c       ........................................................
505
506            IF(forward. OR. leapf) THEN
507              itau= itau + 1
508              iday= day_ini+itau/day_step
509              time= FLOAT(itau-(iday-day_ini)*day_step)/day_step+time_0
510                IF(time.GT.1.) THEN
511                  time = time-1.
512                  iday = iday+1
513                ENDIF
514            ENDIF
515
516
517            IF( itau. EQ. itaufinp1 ) then 
518              if (flag_verif) then
519                write(79,*) 'ucov',ucov
520                write(80,*) 'vcov',vcov
521                write(81,*) 'teta',teta
522                write(82,*) 'ps',ps
523                write(83,*) 'q',q
524                WRITE(85,*) 'q1 = ',q(:,:,1)
525                WRITE(86,*) 'q3 = ',q(:,:,3)
526              endif
527
528              abort_message = 'Simulation finished'
529
530              call abort_gcm(modname,abort_message,0)
531            ENDIF
532c-----------------------------------------------------------------------
533c   ecriture du fichier histoire moyenne:
534c   -------------------------------------
535
536            IF(MOD(itau,iperiod).EQ.0 .OR. itau.EQ.itaufin) THEN
537               IF(itau.EQ.itaufin) THEN
538                  iav=1
539               ELSE
540                  iav=0
541               ENDIF
542               
543               IF (ok_dynzon) THEN
544#ifdef CPP_IOIPSL
545                  CALL writedynav(histaveid, itau,vcov ,
546     ,                 ucov,teta,pk,phi,q,masse,ps,phis)
547                  CALL bilan_dyn (2,dtvr*iperiod,dtvr*day_step*periodav,
548     ,                 ps,masse,pk,pbaru,pbarv,teta,phi,ucov,vcov,q)
549#endif
550               END IF
551
552            ENDIF
553
554c-----------------------------------------------------------------------
555c   ecriture de la bande histoire:
556c   ------------------------------
557
558            IF( MOD(itau,iecri         ).EQ.0) THEN
559c           IF( MOD(itau,iecri*day_step).EQ.0) THEN
560
561              nbetat = nbetatdem
562              CALL geopot(ip1jmp1,teta,pk,pks,phis,phi)
563              unat=0.
564              do l=1,llm
565                unat(iip2:ip1jm,l)=ucov(iip2:ip1jm,l)/cu(iip2:ip1jm)
566                vnat(:,l)=vcov(:,l)/cv(:)
567              enddo
568#ifdef CPP_IOIPSL
569c             CALL writehist(histid,histvid,itau,vcov,
570c     &                      ucov,teta,phi,q,masse,ps,phis)
571#endif
572! For some Grads outputs of fields
573             if (output_grads_dyn) then
574#include "write_grads_dyn.h"
575             endif
576
577            ENDIF ! of IF(MOD(itau,iecri).EQ.0)
578
579            IF(itau.EQ.itaufin) THEN
580
581
582              if (planet_type.eq."earth") then
583#ifdef CPP_EARTH
584! Write an Earth-format restart file
585                CALL dynredem1("restart.nc",0.0,
586     &                         vcov,ucov,teta,q,masse,ps)
587#endif
588              endif ! of if (planet_type.eq."earth")
589
590              CLOSE(99)
591            ENDIF ! of IF (itau.EQ.itaufin)
592
593c-----------------------------------------------------------------------
594c   gestion de l'integration temporelle:
595c   ------------------------------------
596
597            IF( MOD(itau,iperiod).EQ.0 )    THEN
598                    GO TO 1
599            ELSE IF ( MOD(itau-1,iperiod). EQ. 0 ) THEN
600
601                   IF( forward )  THEN
602c      fin du pas forward et debut du pas backward
603
604                      forward = .FALSE.
605                        leapf = .FALSE.
606                           GO TO 2
607
608                   ELSE
609c      fin du pas backward et debut du premier pas leapfrog
610
611                        leapf =  .TRUE.
612                        dt  =  2.*dtvr
613                        GO TO 2
614                   END IF ! of IF (forward)
615            ELSE
616
617c      ......   pas leapfrog  .....
618
619                 leapf = .TRUE.
620                 dt  = 2.*dtvr
621                 GO TO 2
622            END IF ! of IF (MOD(itau,iperiod).EQ.0)
623                   !    ELSEIF (MOD(itau-1,iperiod).EQ.0)
624
625      ELSE ! of IF (.not.purmats)
626
627c       ........................................................
628c       ..............       schema  matsuno        ...............
629c       ........................................................
630            IF( forward )  THEN
631
632             itau =  itau + 1
633             iday = day_ini+itau/day_step
634             time = FLOAT(itau-(iday-day_ini)*day_step)/day_step+time_0
635
636                  IF(time.GT.1.) THEN
637                   time = time-1.
638                   iday = iday+1
639                  ENDIF
640
641               forward =  .FALSE.
642               IF( itau. EQ. itaufinp1 ) then 
643                 abort_message = 'Simulation finished'
644                 call abort_gcm(modname,abort_message,0)
645               ENDIF
646               GO TO 2
647
648            ELSE ! of IF(forward)
649
650              IF(MOD(itau,iperiod).EQ.0 .OR. itau.EQ.itaufin) THEN
651               IF(itau.EQ.itaufin) THEN
652                  iav=1
653               ELSE
654                  iav=0
655               ENDIF
656
657               IF (ok_dynzon) THEN
658#ifdef CPP_IOIPSL
659                  CALL writedynav(histaveid, itau,vcov ,
660     ,                 ucov,teta,pk,phi,q,masse,ps,phis)
661                  CALL bilan_dyn (2,dtvr*iperiod,dtvr*day_step*periodav,
662     ,                 ps,masse,pk,pbaru,pbarv,teta,phi,ucov,vcov,q)
663#endif
664               END IF
665
666              ENDIF ! of IF(MOD(itau,iperiod).EQ.0 .OR. itau.EQ.itaufin)
667
668              IF(MOD(itau,iecri         ).EQ.0) THEN
669c              IF(MOD(itau,iecri*day_step).EQ.0) THEN
670                nbetat = nbetatdem
671                CALL geopot(ip1jmp1,teta,pk,pks,phis,phi)
672                unat=0.
673                do l=1,llm
674                  unat(iip2:ip1jm,l)=ucov(iip2:ip1jm,l)/cu(iip2:ip1jm)
675                  vnat(:,l)=vcov(:,l)/cv(:)
676                enddo
677#ifdef CPP_IOIPSL
678c               CALL writehist( histid, histvid, itau,vcov ,
679c    &                           ucov,teta,phi,q,masse,ps,phis)
680#endif
681! For some Grads outputs
682                if (output_grads_dyn) then
683#include "write_grads_dyn.h"
684                endif
685
686              ENDIF ! of IF(MOD(itau,iecri         ).EQ.0)
687
688              IF(itau.EQ.itaufin) THEN
689                if (planet_type.eq."earth") then
690#ifdef CPP_EARTH
691                  CALL dynredem1("restart.nc",0.0,
692     &                           vcov,ucov,teta,q,masse,ps)
693#endif
694                endif ! of if (planet_type.eq."earth")
695              ENDIF ! of IF(itau.EQ.itaufin)
696
697              forward = .TRUE.
698              GO TO  1
699
700            ENDIF ! of IF (forward)
701
702      END IF ! of IF(.not.purmats)
703
704      STOP
705      END
Note: See TracBrowser for help on using the repository browser.