source: LMDZ4/branches/pre_V3/libf/dyn3dpar/leapfrog.F @ 1387

Last change on this file since 1387 was 630, checked in by Laurent Fairhead, 19 years ago

Import d'une version parallele de la dynamique YM
LF

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