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

Last change on this file since 960 was 960, checked in by lsce, 16 years ago
  • Ajoute du parametre config_inca dans conf_gcm.F config_inca='none'(sans INCA, par defaut) config_inca='chem'(avec INCA config chemie) config_inca='aero'(avec INCA config aerosol)
  • Menage parmis les cles CPP INCA
  • Enleve le calcul d'omega dans calfis.F et active le calcul correspondant dans advtrac.F(avant uniquement pour INCA).

JG

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