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

Last change on this file since 2603 was 2603, checked in by Ehouarn Millour, 8 years ago

Cleanup in the dynamics: turn logic.h into module logic_mod.F90
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: 28.6 KB
Line 
1!
2! $Id: leapfrog.F 2603 2016-07-25 09:31:56Z emillour $
3!
4c
5c
6      SUBROUTINE leapfrog(ucov,vcov,teta,ps,masse,phis,q,time_0)
7
8
9cIM : pour sortir les param. du modele dans un fis. netcdf 110106
10#ifdef CPP_IOIPSL
11      use IOIPSL
12#endif
13      USE infotrac, ONLY: nqtot,ok_iso_verif
14      USE guide_mod, ONLY : guide_main
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
20      use exner_hyb_m, only: exner_hyb
21      use exner_milieu_m, only: exner_milieu
22      USE comvert_mod, ONLY: ap,bp,pressure_exner,presnivs
23      USE comconst_mod, ONLY: cpp, dtphys, dtvr, pi, ihf
24      USE logic_mod, ONLY: iflag_phys,ok_guide,forward,leapf,apphys,
25     &                     statcl,conser,apdiss,purmats,ok_strato
26      USE temps_mod, ONLY: jD_ref,jH_ref,itaufin,day_ini,day_ref,
27     &                        start_time,dt
28
29      IMPLICIT NONE
30
31c      ......   Version  du 10/01/98    ..........
32
33c             avec  coordonnees  verticales hybrides
34c   avec nouveaux operat. dissipation * ( gradiv2,divgrad2,nxgraro2 )
35
36c=======================================================================
37c
38c   Auteur:  P. Le Van /L. Fairhead/F.Hourdin
39c   -------
40c
41c   Objet:
42c   ------
43c
44c   GCM LMD nouvelle grille
45c
46c=======================================================================
47c
48c  ... Dans inigeom , nouveaux calculs pour les elongations  cu , cv
49c      et possibilite d'appeler une fonction f(y)  a derivee tangente
50c      hyperbolique a la  place de la fonction a derivee sinusoidale.
51
52c  ... Possibilite de choisir le shema pour l'advection de
53c        q  , en modifiant iadv dans traceur.def  (10/02) .
54c
55c      Pour Van-Leer + Vapeur d'eau saturee, iadv(1)=4. (F.Codron,10/99)
56c      Pour Van-Leer iadv=10
57c
58c-----------------------------------------------------------------------
59c   Declarations:
60c   -------------
61
62      include "dimensions.h"
63      include "paramet.h"
64      include "comdissnew.h"
65      include "comgeom.h"
66      include "ener.h"
67      include "description.h"
68      include "iniprint.h"
69      include "academic.h"
70
71      REAL,INTENT(IN) :: time_0 ! not used
72
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
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)
101      REAL dteta(ip1jmp1,llm),dq(ip1jmp1,llm,nqtot),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,nqtot),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  iday ! jour julien
118      REAL       time
119
120      REAL  SSUM
121!     REAL 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! 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
137      LOGICAL first,callinigrads
138cIM : pour sortir les param. du modele dans un fis. netcdf 110106
139      save first
140      data first/.true./
141      real dt_cum
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.)
152      logical physic
153
154      data callinigrads/.true./
155      character*10 string10
156
157      REAL :: flxw(ip1jmp1,llm)  ! flux de masse verticale
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
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/
172c-jld
173
174      character*80 dynhist_file, dynhistave_file
175      character(len=*),parameter :: modname="leapfrog"
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
188      logical , parameter :: flag_verif = .false.
189     
190
191      integer itau_w   ! pas de temps ecriture = itap + itau_phy
192
193
194      if (nday>=0) then
195         itaufin   = nday*day_step
196      else
197         itaufin   = -nday
198      endif
199      itaufinp1 = itaufin +1
200      itau = 0
201      physic=.true.
202      if (iflag_phys==0.or.iflag_phys==2) physic=.false.
203
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
210
211
212c-----------------------------------------------------------------------
213c   On initialise la pression et la fonction d'Exner :
214c   --------------------------------------------------
215
216      dq(:,:,:)=0.
217      CALL pression ( ip1jmp1, ap, bp, ps, p       )
218      if (pressure_exner) then
219        CALL exner_hyb( ip1jmp1, ps, p, pks, pk, pkf )
220      else
221        CALL exner_milieu( ip1jmp1, ps, p, pks, pk, pkf )
222      endif
223
224c-----------------------------------------------------------------------
225c   Debut de l'integration temporelle:
226c   ----------------------------------
227
228   1  CONTINUE ! Matsuno Forward step begins here
229
230c   date: (NB: date remains unchanged for Backward step)
231c   -----
232
233      jD_cur = jD_ref + day_ini - day_ref +                             &
234     &          (itau+1)/day_step
235      jH_cur = jH_ref + start_time +                                    &
236     &          mod(itau+1,day_step)/float(day_step)
237      jD_cur = jD_cur + int(jH_cur)
238      jH_cur = jH_cur - int(jH_cur)
239
240        if (ok_iso_verif) then
241           call check_isotopes_seq(q,ip1jmp1,'leapfrog 321')
242        endif !if (ok_iso_verif) then
243
244#ifdef CPP_IOIPSL
245      if (ok_guide) then
246        call guide_main(itau,ucov,vcov,teta,q,masse,ps)
247      endif
248#endif
249
250
251c
252c     IF( MOD( itau, 10* day_step ).EQ.0 )  THEN
253c       CALL  test_period ( ucov,vcov,teta,q,p,phis )
254c       PRINT *,' ----   Test_period apres continue   OK ! -----', itau
255c     ENDIF
256c
257
258! Save fields obtained at previous time step as '...m1'
259      CALL SCOPY( ijmllm ,vcov , 1, vcovm1 , 1 )
260      CALL SCOPY( ijp1llm,ucov , 1, ucovm1 , 1 )
261      CALL SCOPY( ijp1llm,teta , 1, tetam1 , 1 )
262      CALL SCOPY( ijp1llm,masse, 1, massem1, 1 )
263      CALL SCOPY( ip1jmp1, ps  , 1,   psm1 , 1 )
264
265      forward = .TRUE.
266      leapf   = .FALSE.
267      dt      =  dtvr
268
269c   ...    P.Le Van .26/04/94  ....
270! Ehouarn: finvmaold is actually not used
271!      CALL SCOPY   ( ijp1llm,   masse, 1, finvmaold,     1 )
272!      CALL filtreg ( finvmaold ,jjp1, llm, -2,2, .TRUE., 1 )
273
274        if (ok_iso_verif) then
275           call check_isotopes_seq(q,ip1jmp1,'leapfrog 400')
276        endif !if (ok_iso_verif) then
277
278   2  CONTINUE ! Matsuno backward or leapfrog step begins here
279
280c-----------------------------------------------------------------------
281
282c   date: (NB: only leapfrog step requires recomputing date)
283c   -----
284
285      IF (leapf) THEN
286        jD_cur = jD_ref + day_ini - day_ref +
287     &            (itau+1)/day_step
288        jH_cur = jH_ref + start_time +
289     &            mod(itau+1,day_step)/float(day_step)
290        jD_cur = jD_cur + int(jH_cur)
291        jH_cur = jH_cur - int(jH_cur)
292      ENDIF
293
294
295c   gestion des appels de la physique et des dissipations:
296c   ------------------------------------------------------
297c
298c   ...    P.Le Van  ( 6/02/95 )  ....
299
300      apphys = .FALSE.
301      statcl = .FALSE.
302      conser = .FALSE.
303      apdiss = .FALSE.
304
305      IF( purmats ) THEN
306      ! Purely Matsuno time stepping
307         IF( MOD(itau,iconser) .EQ.0.AND.  forward    ) conser = .TRUE.
308         IF( MOD(itau,dissip_period ).EQ.0.AND..NOT.forward )
309     s        apdiss = .TRUE.
310         IF( MOD(itau,iphysiq ).EQ.0.AND..NOT.forward
311     s          .and. physic                        ) apphys = .TRUE.
312      ELSE
313      ! Leapfrog/Matsuno time stepping
314         IF( MOD(itau   ,iconser) .EQ. 0              ) conser = .TRUE.
315         IF( MOD(itau+1,dissip_period).EQ.0 .AND. .NOT. forward )
316     s        apdiss = .TRUE.
317         IF( MOD(itau+1,iphysiq).EQ.0.AND.physic       ) apphys=.TRUE.
318      END IF
319
320! Ehouarn: for Shallow Water case (ie: 1 vertical layer),
321!          supress dissipation step
322      if (llm.eq.1) then
323        apdiss=.false.
324      endif
325
326
327        if (ok_iso_verif) then
328           call check_isotopes_seq(q,ip1jmp1,'leapfrog 589')
329        endif !if (ok_iso_verif) then
330
331c-----------------------------------------------------------------------
332c   calcul des tendances dynamiques:
333c   --------------------------------
334
335      ! compute geopotential phi()
336      CALL geopot  ( ip1jmp1, teta  , pk , pks,  phis  , phi   )
337
338      time = jD_cur + jH_cur
339      CALL caldyn
340     $  ( itau,ucov,vcov,teta,ps,masse,pk,pkf,phis ,
341     $    phi,conser,du,dv,dteta,dp,w, pbaru,pbarv, time )
342
343
344c-----------------------------------------------------------------------
345c   calcul des tendances advection des traceurs (dont l'humidite)
346c   -------------------------------------------------------------
347
348        if (ok_iso_verif) then
349           call check_isotopes_seq(q,ip1jmp1,
350     &           'leapfrog 686: avant caladvtrac')
351        endif !if (ok_iso_verif) then
352
353      IF( forward. OR . leapf )  THEN
354! Ehouarn: NB: fields sent to advtrac are those at the beginning of the time step
355         CALL caladvtrac(q,pbaru,pbarv,
356     *        p, masse, dq,  teta,
357     .        flxw, pk)
358          !write(*,*) 'caladvtrac 346'
359
360         
361         IF (offline) THEN
362Cmaf stokage du flux de masse pour traceurs OFF-LINE
363
364#ifdef CPP_IOIPSL
365           CALL fluxstokenc(pbaru,pbarv,masse,teta,phi,phis,
366     .   dtvr, itau)
367#endif
368
369
370         ENDIF ! of IF (offline)
371c
372      ENDIF ! of IF( forward. OR . leapf )
373
374
375c-----------------------------------------------------------------------
376c   integrations dynamique et traceurs:
377c   ----------------------------------
378
379        if (ok_iso_verif) then
380           write(*,*) 'leapfrog 720'
381           call check_isotopes_seq(q,ip1jmp1,'leapfrog 756')
382        endif !if (ok_iso_verif) then
383       
384       CALL integrd ( nqtot,vcovm1,ucovm1,tetam1,psm1,massem1 ,
385     $         dv,du,dteta,dq,dp,vcov,ucov,teta,q,ps,masse,phis )
386!     $              finvmaold                                    )
387
388       if (ok_iso_verif) then
389          write(*,*) 'leapfrog 724'
390           call check_isotopes_seq(q,ip1jmp1,'leapfrog 762')
391        endif !if (ok_iso_verif) then
392
393c .P.Le Van (26/04/94  ajout de  finvpold dans l'appel d'integrd)
394c
395c-----------------------------------------------------------------------
396c   calcul des tendances physiques:
397c   -------------------------------
398c    ########   P.Le Van ( Modif le  6/02/95 )   ###########
399c
400       IF( purmats )  THEN
401          IF( itau.EQ.itaufin.AND..NOT.forward ) lafin = .TRUE.
402       ELSE
403          IF( itau+1. EQ. itaufin )              lafin = .TRUE.
404       ENDIF
405c
406c
407       IF( apphys )  THEN
408c
409c     .......   Ajout   P.Le Van ( 17/04/96 )   ...........
410c
411
412         CALL pression (  ip1jmp1, ap, bp, ps,  p      )
413         if (pressure_exner) then
414           CALL exner_hyb(  ip1jmp1, ps, p,pks, pk, pkf )
415         else
416           CALL exner_milieu( ip1jmp1, ps, p, pks, pk, pkf )
417         endif
418
419! Appel a geopot ajoute le 2014/05/08 pour garantir la convergence numerique
420! avec dyn3dmem
421         CALL geopot  ( ip1jmp1, teta  , pk , pks,  phis  , phi   )
422
423!           rdaym_ini  = itau * dtvr / daysec
424!           rdayvrai   = rdaym_ini  + day_ini
425!           jD_cur = jD_ref + day_ini - day_ref
426!     $        + int (itau * dtvr / daysec)
427!           jH_cur = jH_ref +                                            &
428!     &              (itau * dtvr / daysec - int(itau * dtvr / daysec))
429           jD_cur = jD_ref + day_ini - day_ref +                        &
430     &          (itau+1)/day_step
431
432           IF (planet_type .eq."generic") THEN
433              ! AS: we make jD_cur to be pday
434              jD_cur = int(day_ini + itau/day_step)
435           ENDIF
436
437           jH_cur = jH_ref + start_time +                               &
438     &              mod(itau+1,day_step)/float(day_step)
439           jD_cur = jD_cur + int(jH_cur)
440           jH_cur = jH_cur - int(jH_cur)
441!         write(lunout,*)'itau, jD_cur = ', itau, jD_cur, jH_cur
442!         call ju2ymds(jD_cur+jH_cur, an, mois, jour, secondes)
443!         write(lunout,*)'current date = ',an, mois, jour, secondes
444
445c rajout debug
446c       lafin = .true.
447
448
449c   Inbterface avec les routines de phylmd (phymars ... )
450c   -----------------------------------------------------
451
452c+jld
453
454c  Diagnostique de conservation de l'énergie : initialisation
455         IF (ip_ebil_dyn.ge.1 ) THEN
456          ztit='bil dyn'
457! Ehouarn: be careful, diagedyn is Earth-specific!
458           IF (planet_type.eq."earth") THEN
459            CALL diagedyn(ztit,2,1,1,dtphys
460     &    , ucov    , vcov , ps, p ,pk , teta , q(:,:,1), q(:,:,2))
461           ENDIF
462         ENDIF ! of IF (ip_ebil_dyn.ge.1 )
463c-jld
464#ifdef CPP_IOIPSL
465cIM decommenter les 6 lignes suivantes pour sortir quelques parametres dynamiques de LMDZ
466cIM uncomment next 6 lines to get some parameters for LMDZ dynamics
467c        IF (first) THEN
468c         first=.false.
469c#include "ini_paramLMDZ_dyn.h"
470c        ENDIF
471c
472c#include "write_paramLMDZ_dyn.h"
473c
474#endif
475! #endif of #ifdef CPP_IOIPSL
476#ifdef CPP_PHYS
477         CALL calfis( lafin , jD_cur, jH_cur,
478     $               ucov,vcov,teta,q,masse,ps,p,pk,phis,phi ,
479     $               du,dv,dteta,dq,
480     $               flxw,dufi,dvfi,dtetafi,dqfi,dpfi  )
481#endif
482c      ajout des tendances physiques:
483c      ------------------------------
484          CALL addfi( dtphys, leapf, forward   ,
485     $                  ucov, vcov, teta , q   ,ps ,
486     $                 dufi, dvfi, dtetafi , dqfi ,dpfi  )
487          ! since addfi updates ps(), also update p(), masse() and pk()
488          CALL pression (ip1jmp1,ap,bp,ps,p)
489          CALL massdair(p,masse)
490          if (pressure_exner) then
491            CALL exner_hyb(ip1jmp1,ps,p,pks,pk,pkf)
492          else
493            CALL exner_milieu(ip1jmp1,ps,p,pks,pk,pkf)
494          endif
495
496         IF (ok_strato) THEN
497           CALL top_bound( vcov,ucov,teta,masse,dtphys)
498         ENDIF
499       
500c
501c  Diagnostique de conservation de l'énergie : difference
502         IF (ip_ebil_dyn.ge.1 ) THEN
503          ztit='bil phys'
504          IF (planet_type.eq."earth") THEN
505           CALL diagedyn(ztit,2,1,1,dtphys
506     &     , ucov    , vcov , ps, p ,pk , teta , q(:,:,1), q(:,:,2))
507          ENDIF
508         ENDIF ! of IF (ip_ebil_dyn.ge.1 )
509
510       ENDIF ! of IF( apphys )
511
512      IF(iflag_phys.EQ.2) THEN ! "Newtonian" case
513!   Academic case : Simple friction and Newtonan relaxation
514!   -------------------------------------------------------
515        DO l=1,llm   
516          DO ij=1,ip1jmp1
517           teta(ij,l)=teta(ij,l)-dtvr*
518     &      (teta(ij,l)-tetarappel(ij,l))*(knewt_g+knewt_t(l)*clat4(ij))
519          ENDDO
520        ENDDO ! of DO l=1,llm
521       
522        if (planet_type.eq."giant") then
523          ! add an intrinsic heat flux at the base of the atmosphere
524          teta(:,1)=teta(:,1)+dtvr*aire(:)*ihf/cpp/masse(:,1)
525        endif
526
527        call friction(ucov,vcov,dtvr)
528       
529        ! Sponge layer (if any)
530        IF (ok_strato) THEN
531!          dufi(:,:)=0.
532!          dvfi(:,:)=0.
533!          dtetafi(:,:)=0.
534!          dqfi(:,:,:)=0.
535!          dpfi(:)=0.
536!          CALL top_bound(vcov,ucov,teta,masse,dufi,dvfi,dtetafi)
537           CALL top_bound( vcov,ucov,teta,masse,dtvr)
538!          CALL addfi( dtvr, leapf, forward   ,
539!     $                  ucov, vcov, teta , q   ,ps ,
540!     $                 dufi, dvfi, dtetafi , dqfi ,dpfi  )
541        ENDIF ! of IF (ok_strato)
542      ENDIF ! of IF (iflag_phys.EQ.2)
543
544
545c-jld
546
547        CALL pression ( ip1jmp1, ap, bp, ps, p                  )
548        if (pressure_exner) then
549          CALL exner_hyb( ip1jmp1, ps, p, pks, pk, pkf )
550        else
551          CALL exner_milieu( ip1jmp1, ps, p, pks, pk, pkf )
552        endif
553        CALL massdair(p,masse)
554
555        if (ok_iso_verif) then
556           call check_isotopes_seq(q,ip1jmp1,'leapfrog 1196')
557        endif !if (ok_iso_verif) then
558
559c-----------------------------------------------------------------------
560c   dissipation horizontale et verticale  des petites echelles:
561c   ----------------------------------------------------------
562
563      IF(apdiss) THEN
564
565
566c   calcul de l'energie cinetique avant dissipation
567        call covcont(llm,ucov,vcov,ucont,vcont)
568        call enercin(vcov,ucov,vcont,ucont,ecin0)
569
570c   dissipation
571        CALL dissip(vcov,ucov,teta,p,dvdis,dudis,dtetadis)
572        ucov=ucov+dudis
573        vcov=vcov+dvdis
574c       teta=teta+dtetadis
575
576
577c------------------------------------------------------------------------
578        if (dissip_conservative) then
579C       On rajoute la tendance due a la transform. Ec -> E therm. cree
580C       lors de la dissipation
581            call covcont(llm,ucov,vcov,ucont,vcont)
582            call enercin(vcov,ucov,vcont,ucont,ecin)
583            dtetaecdt= (ecin0-ecin)/ pk
584c           teta=teta+dtetaecdt
585            dtetadis=dtetadis+dtetaecdt
586        endif
587        teta=teta+dtetadis
588c------------------------------------------------------------------------
589
590
591c    .......        P. Le Van (  ajout  le 17/04/96  )   ...........
592c   ...      Calcul de la valeur moyenne, unique de h aux poles  .....
593c
594
595        DO l  =  1, llm
596          DO ij =  1,iim
597           tppn(ij)  = aire(  ij    ) * teta(  ij    ,l)
598           tpps(ij)  = aire(ij+ip1jm) * teta(ij+ip1jm,l)
599          ENDDO
600           tpn  = SSUM(iim,tppn,1)/apoln
601           tps  = SSUM(iim,tpps,1)/apols
602
603          DO ij = 1, iip1
604           teta(  ij    ,l) = tpn
605           teta(ij+ip1jm,l) = tps
606          ENDDO
607        ENDDO
608
609        if (1 == 0) then
610!!! Ehouarn: lines here 1) kill 1+1=2 in the dynamics
611!!!                     2) should probably not be here anyway
612!!! but are kept for those who would want to revert to previous behaviour
613           DO ij =  1,iim
614             tppn(ij)  = aire(  ij    ) * ps (  ij    )
615             tpps(ij)  = aire(ij+ip1jm) * ps (ij+ip1jm)
616           ENDDO
617             tpn  = SSUM(iim,tppn,1)/apoln
618             tps  = SSUM(iim,tpps,1)/apols
619
620           DO ij = 1, iip1
621             ps(  ij    ) = tpn
622             ps(ij+ip1jm) = tps
623           ENDDO
624        endif ! of if (1 == 0)
625
626      END IF ! of IF(apdiss)
627
628c ajout debug
629c              IF( lafin ) then 
630c                abort_message = 'Simulation finished'
631c                call abort_gcm(modname,abort_message,0)
632c              ENDIF
633       
634c   ********************************************************************
635c   ********************************************************************
636c   .... fin de l'integration dynamique  et physique pour le pas itau ..
637c   ********************************************************************
638c   ********************************************************************
639
640c   preparation du pas d'integration suivant  ......
641
642        if (ok_iso_verif) then
643           call check_isotopes_seq(q,ip1jmp1,'leapfrog 1509')
644        endif !if (ok_iso_verif) then
645
646      IF ( .NOT.purmats ) THEN
647c       ........................................................
648c       ..............  schema matsuno + leapfrog  ..............
649c       ........................................................
650
651            IF(forward. OR. leapf) THEN
652              itau= itau + 1
653c              iday= day_ini+itau/day_step
654c              time= REAL(itau-(iday-day_ini)*day_step)/day_step+time_0
655c                IF(time.GT.1.) THEN
656c                  time = time-1.
657c                  iday = iday+1
658c                ENDIF
659            ENDIF
660
661
662            IF( itau. EQ. itaufinp1 ) then 
663              if (flag_verif) then
664                write(79,*) 'ucov',ucov
665                write(80,*) 'vcov',vcov
666                write(81,*) 'teta',teta
667                write(82,*) 'ps',ps
668                write(83,*) 'q',q
669                WRITE(85,*) 'q1 = ',q(:,:,1)
670                WRITE(86,*) 'q3 = ',q(:,:,3)
671              endif
672
673              abort_message = 'Simulation finished'
674
675              call abort_gcm(modname,abort_message,0)
676            ENDIF
677c-----------------------------------------------------------------------
678c   ecriture du fichier histoire moyenne:
679c   -------------------------------------
680
681            IF(MOD(itau,iperiod).EQ.0 .OR. itau.EQ.itaufin) THEN
682               IF(itau.EQ.itaufin) THEN
683                  iav=1
684               ELSE
685                  iav=0
686               ENDIF
687               
688!              ! Ehouarn: re-compute geopotential for outputs
689               CALL geopot(ip1jmp1,teta,pk,pks,phis,phi)
690
691               IF (ok_dynzon) THEN
692#ifdef CPP_IOIPSL
693                 CALL bilan_dyn(2,dtvr*iperiod,dtvr*day_step*periodav,
694     &                 ps,masse,pk,pbaru,pbarv,teta,phi,ucov,vcov,q)
695#endif
696               END IF
697               IF (ok_dyn_ave) THEN
698#ifdef CPP_IOIPSL
699                 CALL writedynav(itau,vcov,
700     &                 ucov,teta,pk,phi,q,masse,ps,phis)
701#endif
702               ENDIF
703
704            ENDIF ! of IF((MOD(itau,iperiod).EQ.0).OR.(itau.EQ.itaufin))
705
706        if (ok_iso_verif) then
707           call check_isotopes_seq(q,ip1jmp1,'leapfrog 1584')
708        endif !if (ok_iso_verif) then
709
710c-----------------------------------------------------------------------
711c   ecriture de la bande histoire:
712c   ------------------------------
713
714            IF( MOD(itau,iecri).EQ.0) THEN
715             ! Ehouarn: output only during LF or Backward Matsuno
716             if (leapf.or.(.not.leapf.and.(.not.forward))) then
717              CALL geopot(ip1jmp1,teta,pk,pks,phis,phi)
718              unat=0.
719              do l=1,llm
720                unat(iip2:ip1jm,l)=ucov(iip2:ip1jm,l)/cu(iip2:ip1jm)
721                vnat(:,l)=vcov(:,l)/cv(:)
722              enddo
723#ifdef CPP_IOIPSL
724              if (ok_dyn_ins) then
725!               write(lunout,*) "leapfrog: call writehist, itau=",itau
726               CALL writehist(itau,vcov,ucov,teta,phi,q,masse,ps,phis)
727!               call WriteField('ucov',reshape(ucov,(/iip1,jmp1,llm/)))
728!               call WriteField('vcov',reshape(vcov,(/iip1,jjm,llm/)))
729!              call WriteField('teta',reshape(teta,(/iip1,jmp1,llm/)))
730!               call WriteField('ps',reshape(ps,(/iip1,jmp1/)))
731!               call WriteField('masse',reshape(masse,(/iip1,jmp1,llm/)))
732              endif ! of if (ok_dyn_ins)
733#endif
734! For some Grads outputs of fields
735              if (output_grads_dyn) then
736#include "write_grads_dyn.h"
737              endif
738             endif ! of if (leapf.or.(.not.leapf.and.(.not.forward)))
739            ENDIF ! of IF(MOD(itau,iecri).EQ.0)
740
741            IF(itau.EQ.itaufin) THEN
742
743
744!              if (planet_type.eq."earth") then
745! Write an Earth-format restart file
746                CALL dynredem1("restart.nc",start_time,
747     &                         vcov,ucov,teta,q,masse,ps)
748!              endif ! of if (planet_type.eq."earth")
749
750              CLOSE(99)
751              !!! Ehouarn: Why not stop here and now?
752            ENDIF ! of IF (itau.EQ.itaufin)
753
754c-----------------------------------------------------------------------
755c   gestion de l'integration temporelle:
756c   ------------------------------------
757
758            IF( MOD(itau,iperiod).EQ.0 )    THEN
759                    GO TO 1
760            ELSE IF ( MOD(itau-1,iperiod). EQ. 0 ) THEN
761
762                   IF( forward )  THEN
763c      fin du pas forward et debut du pas backward
764
765                      forward = .FALSE.
766                        leapf = .FALSE.
767                           GO TO 2
768
769                   ELSE
770c      fin du pas backward et debut du premier pas leapfrog
771
772                        leapf =  .TRUE.
773                        dt  =  2.*dtvr
774                        GO TO 2
775                   END IF ! of IF (forward)
776            ELSE
777
778c      ......   pas leapfrog  .....
779
780                 leapf = .TRUE.
781                 dt  = 2.*dtvr
782                 GO TO 2
783            END IF ! of IF (MOD(itau,iperiod).EQ.0)
784                   !    ELSEIF (MOD(itau-1,iperiod).EQ.0)
785
786      ELSE ! of IF (.not.purmats)
787
788        if (ok_iso_verif) then
789           call check_isotopes_seq(q,ip1jmp1,'leapfrog 1664')
790        endif !if (ok_iso_verif) then
791
792c       ........................................................
793c       ..............       schema  matsuno        ...............
794c       ........................................................
795            IF( forward )  THEN
796
797             itau =  itau + 1
798c             iday = day_ini+itau/day_step
799c             time = REAL(itau-(iday-day_ini)*day_step)/day_step+time_0
800c
801c                  IF(time.GT.1.) THEN
802c                   time = time-1.
803c                   iday = iday+1
804c                  ENDIF
805
806               forward =  .FALSE.
807               IF( itau. EQ. itaufinp1 ) then 
808                 abort_message = 'Simulation finished'
809                 call abort_gcm(modname,abort_message,0)
810               ENDIF
811               GO TO 2
812
813            ELSE ! of IF(forward) i.e. backward step
814 
815        if (ok_iso_verif) then
816           call check_isotopes_seq(q,ip1jmp1,'leapfrog 1698')
817        endif !if (ok_iso_verif) then 
818
819              IF(MOD(itau,iperiod).EQ.0 .OR. itau.EQ.itaufin) THEN
820               IF(itau.EQ.itaufin) THEN
821                  iav=1
822               ELSE
823                  iav=0
824               ENDIF
825
826!              ! Ehouarn: re-compute geopotential for outputs
827               CALL geopot(ip1jmp1,teta,pk,pks,phis,phi)
828
829               IF (ok_dynzon) THEN
830#ifdef CPP_IOIPSL
831                 CALL bilan_dyn(2,dtvr*iperiod,dtvr*day_step*periodav,
832     &                 ps,masse,pk,pbaru,pbarv,teta,phi,ucov,vcov,q)
833#endif
834               ENDIF
835               IF (ok_dyn_ave) THEN
836#ifdef CPP_IOIPSL
837                 CALL writedynav(itau,vcov,
838     &                 ucov,teta,pk,phi,q,masse,ps,phis)
839#endif
840               ENDIF
841
842              ENDIF ! of IF(MOD(itau,iperiod).EQ.0 .OR. itau.EQ.itaufin)
843
844              IF(MOD(itau,iecri         ).EQ.0) THEN
845c              IF(MOD(itau,iecri*day_step).EQ.0) THEN
846                CALL geopot(ip1jmp1,teta,pk,pks,phis,phi)
847                unat=0.
848                do l=1,llm
849                  unat(iip2:ip1jm,l)=ucov(iip2:ip1jm,l)/cu(iip2:ip1jm)
850                  vnat(:,l)=vcov(:,l)/cv(:)
851                enddo
852#ifdef CPP_IOIPSL
853              if (ok_dyn_ins) then
854!                write(lunout,*) "leapfrog: call writehist (b)",
855!     &                        itau,iecri
856                CALL writehist(itau,vcov,ucov,teta,phi,q,masse,ps,phis)
857              endif ! of if (ok_dyn_ins)
858#endif
859! For some Grads outputs
860                if (output_grads_dyn) then
861#include "write_grads_dyn.h"
862                endif
863
864              ENDIF ! of IF(MOD(itau,iecri         ).EQ.0)
865
866              IF(itau.EQ.itaufin) THEN
867!                if (planet_type.eq."earth") then
868                  CALL dynredem1("restart.nc",start_time,
869     &                           vcov,ucov,teta,q,masse,ps)
870!                endif ! of if (planet_type.eq."earth")
871              ENDIF ! of IF(itau.EQ.itaufin)
872
873              forward = .TRUE.
874              GO TO  1
875
876            ENDIF ! of IF (forward)
877
878      END IF ! of IF(.not.purmats)
879
880      STOP
881      END
Note: See TracBrowser for help on using the repository browser.