source: LMDZ4/trunk/libf/dyn3d/leapfrog.F @ 724

Last change on this file since 724 was 692, checked in by lmdzadmin, 19 years ago

Ajout fichier contenant les flags/parametres de la dynamique paramLMDZ_dyn.nc
IM

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