source: LMDZ5/trunk/libf/dyn3dmem/leapfrog_loc.F @ 2597

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

Cleanup in the dynamics: get rid of comconst.h, make it a module comconst_mod.
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
File size: 55.8 KB
Line 
1!
2! $Id$
3!
4c
5c
6#define DEBUG_IO
7#undef DEBUG_IO
8
9
10      SUBROUTINE leapfrog_loc(ucov0,vcov0,teta0,ps0,
11     &                        masse0,phis0,q0,time_0)
12
13       USE misc_mod
14       USE parallel_lmdz
15       USE times
16       USE mod_hallo
17       USE Bands
18       USE Write_Field
19       USE Write_Field_p
20       USE vampir
21       USE timer_filtre, ONLY : print_filtre_timer
22       USE infotrac
23       USE guide_loc_mod, ONLY : guide_main
24       USE getparam
25       USE control_mod
26       USE mod_filtreg_p
27       USE write_field_loc
28       USE allocate_field_mod
29       USE call_dissip_mod, ONLY : call_dissip
30       USE call_calfis_mod, ONLY : call_calfis
31       USE leapfrog_mod
32       use exner_hyb_loc_m, only: exner_hyb_loc
33       use exner_milieu_loc_m, only: exner_milieu_loc
34       USE comconst_mod, ONLY: cpp, dtvr, ihf
35      IMPLICIT NONE
36
37c      ......   Version  du 10/01/98    ..........
38
39c             avec  coordonnees  verticales hybrides
40c   avec nouveaux operat. dissipation * ( gradiv2,divgrad2,nxgraro2 )
41
42c=======================================================================
43c
44c   Auteur:  P. Le Van /L. Fairhead/F.Hourdin
45c   -------
46c
47c   Objet:
48c   ------
49c
50c   GCM LMD nouvelle grille
51c
52c=======================================================================
53c
54c  ... Dans inigeom , nouveaux calculs pour les elongations  cu , cv
55c      et possibilite d'appeler une fonction f(y)  a derivee tangente
56c      hyperbolique a la  place de la fonction a derivee sinusoidale.
57
58c  ... Possibilite de choisir le shema pour l'advection de
59c        q  , en modifiant iadv dans traceur.def  (10/02) .
60c
61c      Pour Van-Leer + Vapeur d'eau saturee, iadv(1)=4. (F.Codron,10/99)
62c      Pour Van-Leer iadv=10
63c
64c-----------------------------------------------------------------------
65c   Declarations:
66c   -------------
67
68      include "dimensions.h"
69      include "paramet.h"
70      include "comdissnew.h"
71      include "comvert.h"
72      include "comgeom.h"
73      include "logic.h"
74      include "temps.h"
75      include "ener.h"
76      include "description.h"
77      include "serre.h"
78      include "iniprint.h"
79      include "academic.h"
80     
81      REAL,INTENT(IN) :: time_0 ! not used
82
83c   dynamical variables:
84      REAL,INTENT(IN) :: ucov0(ijb_u:ije_u,llm)    ! zonal covariant wind
85      REAL,INTENT(IN) :: vcov0(ijb_v:ije_v,llm)    ! meridional covariant wind
86      REAL,INTENT(IN) :: teta0(ijb_u:ije_u,llm)    ! potential temperature
87      REAL,INTENT(IN) :: q0(ijb_u:ije_u,llm,nqtot) ! advected tracers
88      REAL,INTENT(IN) :: ps0(ijb_u:ije_u)          ! surface pressure (Pa)
89      REAL,INTENT(IN) :: masse0(ijb_u:ije_u,llm)   ! air mass
90      REAL,INTENT(IN) :: phis0(ijb_u:ije_u)        ! geopotentiat at the surface
91
92      real zqmin,zqmax
93
94!      REAL,SAVE,ALLOCATABLE :: p (:,:  )               ! pression aux interfac.des couches
95!      REAL,SAVE,ALLOCATABLE :: pks(:)                      ! exner au  sol
96!      REAL,SAVE,ALLOCATABLE :: pk(:,:)                   ! exner au milieu des couches
97!      REAL,SAVE,ALLOCATABLE :: pkf(:,:)                  ! exner filt.au milieu des couches
98!      REAL,SAVE,ALLOCATABLE :: phi(:,:)                  ! geopotentiel
99!      REAL,SAVE,ALLOCATABLE :: w(:,:)                    ! vitesse verticale
100
101c variables dynamiques intermediaire pour le transport
102!      REAL,SAVE,ALLOCATABLE :: pbaru(:,:),pbarv(:,:) !flux de masse
103
104c   variables dynamiques au pas -1
105!      REAL,SAVE,ALLOCATABLE :: vcovm1(:,:),ucovm1(:,:)
106!      REAL,SAVE,ALLOCATABLE :: tetam1(:,:),psm1(:)
107!      REAL,SAVE,ALLOCATABLE :: massem1(:,:)
108
109c   tendances dynamiques
110!      REAL,SAVE,ALLOCATABLE :: dv(:,:),du(:,:)
111!      REAL,SAVE,ALLOCATABLE :: dteta(:,:),dp(:)
112!      REAL,DIMENSION(:,:,:), ALLOCATABLE, SAVE :: dq
113
114c   tendances de la dissipation
115!      REAL,SAVE,ALLOCATABLE :: dvdis(:,:),dudis(:,:)
116!      REAL,SAVE,ALLOCATABLE :: dtetadis(:,:)
117
118c   tendances physiques
119      REAL,SAVE,ALLOCATABLE :: dvfi(:,:),dufi(:,:)
120      REAL,SAVE,ALLOCATABLE :: dtetafi(:,:)
121      REAL,SAVE,ALLOCATABLE :: dpfi(:)
122      REAL,DIMENSION(:,:,:),ALLOCATABLE,SAVE :: dqfi
123
124c   variables pour le fichier histoire
125      REAL dtav      ! intervalle de temps elementaire
126
127      REAL tppn(iim),tpps(iim),tpn,tps
128c
129      INTEGER itau,itaufinp1,iav
130!      INTEGER  iday ! jour julien
131      REAL       time
132
133      REAL  SSUM
134!      REAL,SAVE,ALLOCATABLE :: finvmaold(:,:)
135
136cym      LOGICAL  lafin
137      LOGICAL :: lafin
138      INTEGER ij,iq,l
139      INTEGER ik
140
141      real time_step, t_wrt, t_ops
142
143! jD_cur: jour julien courant
144! jH_cur: heure julienne courante
145      REAL :: jD_cur, jH_cur
146      INTEGER :: an, mois, jour
147      REAL :: secondes
148
149      logical :: physic
150      LOGICAL first,callinigrads
151
152      data callinigrads/.true./
153      character*10 string10
154
155!      REAL,SAVE,ALLOCATABLE :: flxw(:,:) ! flux de masse verticale
156
157c+jld variables test conservation energie
158!      REAL,SAVE,ALLOCATABLE :: ecin(:,:),ecin0(:,:)
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,SAVE,ALLOCATABLE :: dtetaecdt(:,:)
163!      REAL,SAVE,ALLOCATABLE :: vcont(:,:),ucont(:,:)
164!      REAL,SAVE,ALLOCATABLE :: vnat(:,:),unat(:,:)
165      REAL      d_h_vcol, d_qt, d_qw, d_ql, d_ec
166      CHARACTER*15 ztit
167!!      INTEGER   ip_ebil_dyn  ! PRINT level for energy conserv. diag.
168!      SAVE      ip_ebil_dyn
169!      DATA      ip_ebil_dyn/0/
170c-jld
171
172      character*80 dynhist_file, dynhistave_file
173      character(len=*),parameter :: modname="leapfrog_loc"
174      character*80 abort_message
175
176
177      logical,PARAMETER :: dissip_conservative=.TRUE.
178 
179      INTEGER testita
180      PARAMETER (testita = 9)
181
182      logical , parameter :: flag_verif = .false.
183     
184c declaration liees au parallelisme
185      INTEGER :: ierr
186      LOGICAL :: FirstCaldyn
187      LOGICAL :: FirstPhysic
188      INTEGER :: ijb,ije,j,i
189      type(Request) :: TestRequest
190      type(Request) :: Request_Dissip
191      type(Request) :: Request_physic
192
193      INTEGER :: true_itau
194      INTEGER :: iapptrac
195      INTEGER :: AdjustCount
196!      INTEGER :: var_time
197      LOGICAL :: ok_start_timer=.FALSE.
198      LOGICAL, SAVE :: firstcall=.TRUE.
199      TYPE(distrib),SAVE :: new_dist
200
201      if (ok_iso_verif) then
202         call check_isotopes(q0,ijb_u,ije_u,'leapfrog204: debut')
203      endif !if (ok_iso_verif) then
204     
205c$OMP MASTER
206      ItCount=0
207c$OMP END MASTER     
208      true_itau=0
209      FirstCaldyn=.TRUE.
210      FirstPhysic=.TRUE.
211      iapptrac=0
212      AdjustCount = 0
213      lafin=.false.
214     
215      if (nday>=0) then
216         itaufin   = nday*day_step
217      else
218         itaufin   = -nday
219      endif
220
221      itaufinp1 = itaufin +1
222
223      if (ok_iso_verif) then
224        call check_isotopes(q0,ijb_u,ije_u,'leapfrog 226')
225      endif !if (ok_iso_verif) then
226
227      itau = 0
228      physic=.true.
229      if (iflag_phys==0.or.iflag_phys==2) physic=.false.
230      CALL init_nan
231      CALL leapfrog_allocate
232      ucov=ucov0
233      vcov=vcov0
234      teta=teta0
235      ps=ps0
236      masse=masse0
237      phis=phis0
238      q=q0
239
240      if (ok_iso_verif) then
241        call check_isotopes(q,ijb_u,ije_u,'leapfrog 239')
242      endif !if (ok_iso_verif) then
243     
244!      iday = day_ini+itau/day_step
245!      time = REAL(itau-(iday-day_ini)*day_step)/day_step+time_0
246!         IF(time.GT.1.) THEN
247!          time = time-1.
248!          iday = iday+1
249!         ENDIF
250
251c Allocate variables depending on dynamic variable nqtot
252!$OMP MASTER
253      if (firstcall) then
254!     
255!      ALLOCATE(p(ijb_u:ije_u,llmp1))
256!      ALLOCATE(pks(ijb_u:ije_u))
257!      ALLOCATE(pk(ijb_u:ije_u,llm))
258!      ALLOCATE(pkf(ijb_u:ije_u,llm))
259!      ALLOCATE(phi(ijb_u:ije_u,llm))
260!      ALLOCATE(w(ijb_u:ije_u,llm))
261!      ALLOCATE(pbaru(ip1jmp1,llm),pbarv(ip1jm,llm))
262!      ALLOCATE(vcovm1(ijb_v:ije_v,llm),ucovm1(ijb_u:ije_u,llm))
263!      ALLOCATE(tetam1(ijb_u:ije_u,llm),psm1(ijb_u:ije_u))
264!      ALLOCATE(massem1(ijb_u:ije_u,llm))
265!      ALLOCATE(dv(ijb_v:ije_v,llm),du(ijb_u:ije_u,llm))
266!      ALLOCATE(dteta(ijb_u:ije_u,llm),dp(ijb_u:ije_u))     
267!      ALLOCATE(dvdis(ijb_v:ije_v,llm),dudis(ijb_u:ije_u,llm))
268!      ALLOCATE(dtetadis(ijb_u:ije_u,llm))
269      ALLOCATE(dvfi(ijb_v:ije_v,llm),dufi(ijb_u:ije_u,llm))
270      ALLOCATE(dtetafi(ijb_u:ije_u,llm))
271      ALLOCATE(dpfi(ijb_u:ije_u))
272!      ALLOCATE(dq(ijb_u:ije_u,llm,nqtot))
273      ALLOCATE(dqfi(ijb_u:ije_u,llm,nqtot))
274!      ALLOCATE(dqfi_tmp(iip1,llm,nqtot))
275!      ALLOCATE(finvmaold(ijb_u:ije_u,llm))
276!      ALLOCATE(flxw(ijb_u:ije_u,llm))
277!      ALLOCATE(ecin(ijb_u:ije_u,llm),ecin0(ijb_u:ije_u,llm))
278!      ALLOCATE(dtetaecdt(ijb_u:ije_u,llm))
279!      ALLOCATE(vcont(ijb_v:ije_v,llm),ucont(ijb_u:ije_u,llm))
280!      ALLOCATE(vnat(ijb_v:ije_v,llm),unat(ijb_u:ije_u,llm))
281      endif
282!$OMP END MASTER     
283!$OMP BARRIER
284
285!                CALL dynredem1_loc("restart.nc",0.0,
286!     &                           vcov,ucov,teta,q,masse,ps)
287
288
289c-----------------------------------------------------------------------
290c   On initialise la pression et la fonction d'Exner :
291c   --------------------------------------------------
292
293c$OMP MASTER
294      dq(:,:,:)=0.
295      CALL pression ( ijnb_u, ap, bp, ps, p       )
296c$OMP END MASTER
297      if (pressure_exner) then
298      CALL exner_hyb_loc( ijnb_u, ps, p, pks, pk, pkf)
299      else
300        CALL exner_milieu_loc( ijnb_u, ps, p, pks, pk, pkf )
301      endif
302c-----------------------------------------------------------------------
303c   Debut de l'integration temporelle:
304c   ----------------------------------
305c et du parallelisme !!
306
307   1  CONTINUE ! Matsuno Forward step begins here
308
309c   date: (NB: date remains unchanged for Backward step)
310c   -----
311
312      jD_cur = jD_ref + day_ini - day_ref +                             &
313     &          (itau+1)/day_step
314      jH_cur = jH_ref + start_time +                                    &
315     &         mod(itau+1,day_step)/float(day_step)
316      if (jH_cur > 1.0 ) then
317        jD_cur = jD_cur +1.
318        jH_cur = jH_cur -1.
319      endif
320
321        if (ok_iso_verif) then
322           call check_isotopes(q,ijb_u,ije_u,'leapfrog 321')
323        endif !if (ok_iso_verif) then
324
325#ifdef CPP_IOIPSL
326      if (ok_guide) then
327        call guide_main(itau,ucov,vcov,teta,q,masse,ps)
328!$OMP BARRIER
329      endif
330#endif
331
332
333c
334c     IF( MOD( itau, 10* day_step ).EQ.0 )  THEN
335c       CALL  test_period ( ucov,vcov,teta,q,p,phis )
336c       PRINT *,' ----   Test_period apres continue   OK ! -----', itau
337c     ENDIF
338c
339cym      CALL SCOPY( ijmllm ,vcov , 1, vcovm1 , 1 )
340cym      CALL SCOPY( ijp1llm,ucov , 1, ucovm1 , 1 )
341cym      CALL SCOPY( ijp1llm,teta , 1, tetam1 , 1 )
342cym      CALL SCOPY( ijp1llm,masse, 1, massem1, 1 )
343cym      CALL SCOPY( ip1jmp1, ps  , 1,   psm1 , 1 )
344
345       if (FirstCaldyn) then
346c$OMP MASTER
347         ucovm1=ucov
348         vcovm1=vcov
349         tetam1= teta
350         massem1= masse
351         psm1= ps
352         
353! Ehouarn: finvmaold is actually not used       
354!         finvmaold = masse
355c$OMP END MASTER
356c$OMP BARRIER
357!         CALL filtreg_p ( finvmaold ,jjb_u,jje_u,jjb_u,jje_u,jjp1, llm,
358!     &                    -2,2, .TRUE., 1 )
359       else
360! Save fields obtained at previous time step as '...m1'
361         ijb=ij_begin
362         ije=ij_end
363
364c$OMP MASTER           
365         psm1     (ijb:ije) = ps    (ijb:ije)
366c$OMP END MASTER
367
368c$OMP DO SCHEDULE(STATIC,OMP_CHUNK)         
369         DO l=1,llm     
370           ije=ij_end
371           ucovm1   (ijb:ije,l) = ucov  (ijb:ije,l)
372           tetam1   (ijb:ije,l) = teta  (ijb:ije,l)
373           massem1  (ijb:ije,l) = masse (ijb:ije,l)
374!           finvmaold(ijb:ije,l)=masse(ijb:ije,l)
375                 
376           if (pole_sud) ije=ij_end-iip1
377           vcovm1(ijb:ije,l) = vcov  (ijb:ije,l)
378       
379
380         ENDDO
381c$OMP ENDDO 
382
383
384! Ehouarn: finvmaold not used
385!          CALL filtreg_p(finvmaold ,jjb_u,jje_u,jj_begin,jj_end,jjp1,
386!     .                    llm, -2,2, .TRUE., 1 )
387
388       endif ! of if (FirstCaldyn)
389       
390      forward = .TRUE.
391      leapf   = .FALSE.
392      dt      =  dtvr
393
394c   ...    P.Le Van .26/04/94  ....
395
396cym      CALL SCOPY   ( ijp1llm,   masse, 1, finvmaold,     1 )
397cym      CALL filtreg ( finvmaold ,jjp1, llm, -2,2, .TRUE., 1 )
398
399cym  ne sert a rien
400cym      call minmax(ijp1llm,q(:,:,3),zqmin,zqmax)
401
402
403        if (ok_iso_verif) then
404           call check_isotopes(q,ijb_u,ije_u,'leapfrog 400')
405        endif !if (ok_iso_verif) then
406
407   2  CONTINUE ! Matsuno backward or leapfrog step begins here
408
409
410        if (ok_iso_verif) then
411           call check_isotopes(q,ijb_u,ije_u,'leapfrog 402')
412        endif !if (ok_iso_verif) then
413
414c$OMP MASTER
415      ItCount=ItCount+1
416      if (MOD(ItCount,1)==1) then
417        debug=.true.
418      else
419        debug=.false.
420      endif
421c$OMP END MASTER
422c-----------------------------------------------------------------------
423
424c   date: (NB: only leapfrog step requires recomputing date)
425c   -----
426
427      IF (leapf) THEN
428        jD_cur = jD_ref + day_ini - day_ref +
429     &          (itau+1)/day_step
430        jH_cur = jH_ref + start_time +
431     &         mod(itau+1,day_step)/float(day_step)
432        if (jH_cur > 1.0 ) then
433          jD_cur = jD_cur +1.
434          jH_cur = jH_cur -1.
435        endif
436      ENDIF
437
438c   gestion des appels de la physique et des dissipations:
439c   ------------------------------------------------------
440c
441c   ...    P.Le Van  ( 6/02/95 )  ....
442
443      apphys = .FALSE.
444      statcl = .FALSE.
445      conser = .FALSE.
446      apdiss = .FALSE.
447
448      IF( purmats ) THEN
449      ! Purely Matsuno time stepping
450         IF( MOD(itau,iconser) .EQ.0.AND.  forward    ) conser = .TRUE.
451         IF( MOD(itau,dissip_period ).EQ.0.AND..NOT.forward )
452     s        apdiss = .TRUE.
453         IF( MOD(itau,iphysiq ).EQ.0.AND..NOT.forward
454     s          .and. physic                        ) apphys = .TRUE.
455      ELSE
456      ! Leapfrog/Matsuno time stepping
457         IF( MOD(itau   ,iconser) .EQ. 0              ) conser = .TRUE.
458         IF( MOD(itau+1,dissip_period).EQ.0 .AND. .NOT. forward )
459     s        apdiss = .TRUE.
460         IF( MOD(itau+1,iphysiq).EQ.0.AND.physic) apphys=.TRUE.
461      END IF
462
463! Ehouarn: for Shallow Water case (ie: 1 vertical layer),
464!          supress dissipation step
465      if (llm.eq.1) then
466        apdiss=.false.
467      endif
468
469cym    ---> Pour le moment     
470cym      apphys = .FALSE.
471      statcl = .FALSE.
472      conser = .FALSE. ! ie: no output of control variables to stdout in //
473     
474      if (firstCaldyn) then
475c$OMP MASTER
476          call Set_Distrib(distrib_caldyn)
477c$OMP END MASTER
478c$OMP BARRIER
479          firstCaldyn=.FALSE.
480cym          call InitTime
481c$OMP MASTER
482          call Init_timer
483c$OMP END MASTER
484      endif
485
486c$OMP MASTER     
487      IF (ok_start_timer) THEN
488        CALL InitTime
489        ok_start_timer=.FALSE.
490      ENDIF     
491c$OMP END MASTER     
492
493
494        if (ok_iso_verif) then
495           call check_isotopes(q,ijb_u,ije_u,'leapfrog 471')
496        endif !if (ok_iso_verif) then
497
498!ym  PAS D'AJUSTEMENT POUR LE MOMENT     
499      if (Adjust) then
500        AdjustCount=AdjustCount+1
501!        if (iapptrac==iapp_tracvl .and. (forward. OR . leapf)
502!     &         .and. itau/iphysiq>2 .and. Adjustcount>30) then
503        if (Adjustcount>1) then
504           AdjustCount=0
505c$OMP MASTER
506           call allgather_timer_average
507
508        if (prt_level > 9) then
509       
510        print *,'*********************************'
511        print *,'******    TIMER CALDYN     ******'
512        do i=0,mpi_size-1
513          print *,'proc',i,' :   Nb Bandes  :',jj_nb_caldyn(i),
514     &            '  : temps moyen :',
515     &             timer_average(jj_nb_caldyn(i),timer_caldyn,i),
516     &            '+-',timer_delta(jj_nb_caldyn(i),timer_caldyn,i)
517        enddo
518     
519        print *,'*********************************'
520        print *,'******    TIMER VANLEER    ******'
521        do i=0,mpi_size-1
522          print *,'proc',i,' :   Nb Bandes  :',jj_nb_vanleer(i),
523     &            '  : temps moyen :',
524     &             timer_average(jj_nb_vanleer(i),timer_vanleer,i),
525     &            '+-',timer_delta(jj_nb_vanleer(i),timer_vanleer,i)
526        enddo
527     
528        print *,'*********************************'
529        print *,'******    TIMER DISSIP    ******'
530        do i=0,mpi_size-1
531          print *,'proc',i,' :   Nb Bandes  :',jj_nb_dissip(i),
532     &            '  : temps moyen :',
533     &             timer_average(jj_nb_dissip(i),timer_dissip,i),
534     &             '+-',timer_delta(jj_nb_dissip(i),timer_dissip,i)
535        enddo
536       
537!        if (mpi_rank==0) call WriteBands
538       
539       endif
540       
541         call AdjustBands_caldyn(new_dist)
542!$OMP END MASTER
543!$OMP BARRIER
544         CALL leapfrog_switch_caldyn(new_dist)
545!$OMP BARRIER
546
547
548!$OMP MASTER
549         distrib_caldyn=new_dist
550         CALL set_distrib(distrib_caldyn)
551!$OMP END MASTER
552!$OMP BARRIER
553!         call Register_SwapFieldHallo(ucov,ucov,ip1jmp1,llm,
554!     &                                jj_Nb_caldyn,0,0,TestRequest)
555!         call Register_SwapFieldHallo(ucovm1,ucovm1,ip1jmp1,llm,
556!     &                                jj_Nb_caldyn,0,0,TestRequest)
557!         call Register_SwapFieldHallo(vcov,vcov,ip1jm,llm,
558!     &                                jj_Nb_caldyn,0,0,TestRequest)
559!         call Register_SwapFieldHallo(vcovm1,vcovm1,ip1jm,llm,
560!     &                                jj_Nb_caldyn,0,0,TestRequest)
561!         call Register_SwapFieldHallo(teta,teta,ip1jmp1,llm,
562!     &                                jj_Nb_caldyn,0,0,TestRequest)
563!         call Register_SwapFieldHallo(tetam1,tetam1,ip1jmp1,llm,
564!     &                                jj_Nb_caldyn,0,0,TestRequest)
565!         call Register_SwapFieldHallo(masse,masse,ip1jmp1,llm,
566!     &                                jj_Nb_caldyn,0,0,TestRequest)
567!         call Register_SwapFieldHallo(massem1,massem1,ip1jmp1,llm,
568!     &                                jj_Nb_caldyn,0,0,TestRequest)
569!         call Register_SwapFieldHallo(ps,ps,ip1jmp1,1,
570!     &                                jj_Nb_caldyn,0,0,TestRequest)
571!         call Register_SwapFieldHallo(psm1,psm1,ip1jmp1,1,
572!     &                                jj_Nb_caldyn,0,0,TestRequest)
573!         call Register_SwapFieldHallo(pkf,pkf,ip1jmp1,llm,
574!     &                                jj_Nb_caldyn,0,0,TestRequest)
575!         call Register_SwapFieldHallo(pk,pk,ip1jmp1,llm,
576!     &                                jj_Nb_caldyn,0,0,TestRequest)
577!         call Register_SwapFieldHallo(pks,pks,ip1jmp1,1,
578!     &                                jj_Nb_caldyn,0,0,TestRequest)
579!         call Register_SwapFieldHallo(phis,phis,ip1jmp1,1,
580!     &                                jj_Nb_caldyn,0,0,TestRequest)
581!         call Register_SwapFieldHallo(phi,phi,ip1jmp1,llm,
582!     &                                jj_Nb_caldyn,0,0,TestRequest)
583!         call Register_SwapFieldHallo(finvmaold,finvmaold,ip1jmp1,llm,
584!     &                                jj_Nb_caldyn,0,0,TestRequest)
585!
586!        do j=1,nqtot
587!         call Register_SwapFieldHallo(q(:,:,j),q(:,:,j),ip1jmp1,llm,
588!     &                                jj_nb_caldyn,0,0,TestRequest)
589!        enddo
590!
591!         call Set_Distrib(distrib_caldyn)
592!         call SendRequest(TestRequest)
593!         call WaitRequest(TestRequest)
594         
595!$OMP MASTER
596        call AdjustBands_dissip(new_dist)
597!$OMP END MASTER
598!$OMP BARRIER
599        CALL leapfrog_switch_dissip(new_dist)
600!$OMP BARRIER
601!$OMP MASTER
602        distrib_dissip=new_dist
603!$OMP END MASTER
604!$OMP BARRIER
605!        call AdjustBands_physic
606
607c$OMP MASTER 
608        if (mpi_rank==0) call WriteBands
609c$OMP END MASTER 
610
611
612      endif
613      endif       
614     
615     
616        if (ok_iso_verif) then
617           call check_isotopes(q,ijb_u,ije_u,'leapfrog 589')
618        endif !if (ok_iso_verif) then
619     
620c-----------------------------------------------------------------------
621c   calcul des tendances dynamiques:
622c   --------------------------------
623c$OMP BARRIER
624c$OMP MASTER
625       call VTb(VThallo)
626c$OMP END MASTER
627
628       call Register_Hallo_u(ucov,llm,1,1,1,1,TestRequest)
629       call Register_Hallo_v(vcov,llm,1,1,1,1,TestRequest)
630       call Register_Hallo_u(teta,llm,1,1,1,1,TestRequest)
631       call Register_Hallo_u(ps,1,1,2,2,1,TestRequest)
632       call Register_Hallo_u(pkf,llm,1,1,1,1,TestRequest)
633       call Register_Hallo_u(pk,llm,1,1,1,1,TestRequest)
634       call Register_Hallo_u(pks,1,1,1,1,1,TestRequest)
635       call Register_Hallo_u(p,llmp1,1,1,1,1,TestRequest)
636       
637c       do j=1,nqtot
638c         call Register_Hallo(q(1,1,j),ip1jmp1,llm,1,1,1,1,
639c     *                       TestRequest)
640c        enddo
641
642       call SendRequest(TestRequest)
643c$OMP BARRIER
644       call WaitRequest(TestRequest)
645
646c$OMP MASTER
647       call VTe(VThallo)
648c$OMP END MASTER
649c$OMP BARRIER
650     
651      if (debug) then       
652        call WriteField_u('ucov',ucov)
653        call WriteField_v('vcov',vcov)
654        call WriteField_u('teta',teta)
655        call WriteField_u('ps',ps)
656        call WriteField_u('masse',masse)
657        call WriteField_u('pk',pk)
658        call WriteField_u('pks',pks)
659        call WriteField_u('pkf',pkf)
660        call WriteField_u('phis',phis)
661        do j=1,nqtot
662          call WriteField_u('q'//trim(int2str(j)),
663     .                q(:,:,j))
664        enddo
665      endif
666
667     
668      True_itau=True_itau+1
669
670c$OMP MASTER
671      IF (prt_level>9) THEN
672        WRITE(lunout,*)"leapfrog_p: Iteration No",True_itau
673      ENDIF
674
675
676      call start_timer(timer_caldyn)
677
678      ! compute geopotential phi()
679      CALL geopot_loc  ( ip1jmp1, teta  , pk , pks,  phis  , phi   )
680       
681        if (ok_iso_verif) then
682           call check_isotopes(q,ijb_u,ije_u,'leapfrog 651')
683        endif !if (ok_iso_verif) then
684     
685      call VTb(VTcaldyn)
686c$OMP END MASTER
687!      var_time=time+iday-day_ini
688
689c$OMP BARRIER
690!      CALL FTRACE_REGION_BEGIN("caldyn")
691      time = jD_cur + jH_cur
692
693      CALL caldyn_loc
694     $  ( itau,ucov,vcov,teta,ps,masse,pk,pkf,phis ,
695     $    phi,conser,du,dv,dteta,dp,w, pbaru,pbarv, time )
696
697!      CALL FTRACE_REGION_END("caldyn")
698
699c$OMP MASTER
700      call VTe(VTcaldyn)
701c$OMP END MASTER     
702
703#ifdef DEBUG_IO   
704      call WriteField_u('du',du)
705      call WriteField_v('dv',dv)
706      call WriteField_u('dteta',dteta)
707      call WriteField_u('dp',dp)
708      call WriteField_u('w',w)
709      call WriteField_u('pbaru',pbaru)
710      call WriteField_v('pbarv',pbarv)
711      call WriteField_u('p',p)
712      call WriteField_u('masse',masse)
713      call WriteField_u('pk',pk)
714#endif
715c-----------------------------------------------------------------------
716c   calcul des tendances advection des traceurs (dont l'humidite)
717c   -------------------------------------------------------------
718
719        if (ok_iso_verif) then
720           call check_isotopes(q,ijb_u,ije_u,
721     &           'leapfrog 686: avant caladvtrac')
722        endif !if (ok_iso_verif) then
723     
724      IF( forward. OR . leapf )  THEN
725! Ehouarn: NB: fields sent to advtrac are those at the beginning of the time step
726        !write(*,*) 'leapfrog 679: avant CALL caladvtrac_loc'
727         CALL caladvtrac_loc(q,pbaru,pbarv,
728     *        p, masse, dq,  teta,
729     .        flxw,pk, iapptrac)
730
731         !write(*,*) 'leapfrog 719'
732         if (ok_iso_verif) then
733           call check_isotopes(q,ijb_u,ije_u,
734     &           'leapfrog 698: apres caladvtrac')
735         endif !if (ok_iso_verif) then
736
737!      do j=1,nqtot
738!        call WriteField_u('qadv'//trim(int2str(j)),q(:,:,j))
739!      enddo
740
741! Ehouarn: Storage of mass flux for off-line tracers... not implemented...
742
743      ENDIF ! of IF( forward. OR . leapf )
744
745
746c-----------------------------------------------------------------------
747c   integrations dynamique et traceurs:
748c   ----------------------------------
749
750c$OMP MASTER
751       call VTb(VTintegre)
752c$OMP END MASTER
753#ifdef DEBUG_IO   
754      if (true_itau>20) then
755      call WriteField_u('ucovm1',ucovm1)
756      call WriteField_v('vcovm1',vcovm1)
757      call WriteField_u('tetam1',tetam1)
758      call WriteField_u('psm1',psm1)
759      call WriteField_u('ucov_int',ucov)
760      call WriteField_v('vcov_int',vcov)
761      call WriteField_u('teta_int',teta)
762      call WriteField_u('ps_int',ps)
763      endif
764#endif
765c$OMP BARRIER
766!       CALL FTRACE_REGION_BEGIN("integrd")
767
768       !write(*,*) 'leapfrog 720'
769        if (ok_iso_verif) then
770           call check_isotopes(q,ijb_u,ije_u,'leapfrog 756')
771        endif !if (ok_iso_verif) then
772
773       ! CRisi: pourquoi aller jusqu'à 2 et non pas jusqu'à nqtot??
774       CALL integrd_loc ( nqtot,vcovm1,ucovm1,tetam1,psm1,massem1 ,
775     $         dv,du,dteta,dq,dp,vcov,ucov,teta,q,ps,masse,phis)
776!     $              finvmaold                                    )
777
778       !write(*,*) 'leapfrog 724'       
779        if (ok_iso_verif) then
780           call check_isotopes(q,ijb_u,ije_u,'leapfrog 762')
781        endif !if (ok_iso_verif) then
782 
783!       CALL FTRACE_REGION_END("integrd")
784c$OMP BARRIER
785#ifdef DEBUG_IO   
786      call WriteField_u('ucovm1',ucovm1)
787      call WriteField_v('vcovm1',vcovm1)
788      call WriteField_u('tetam1',tetam1)
789      call WriteField_u('psm1',psm1)
790      call WriteField_u('ucov_int',ucov)
791      call WriteField_v('vcov_int',vcov)
792      call WriteField_u('teta_int',teta)
793      call WriteField_u('ps_int',ps)
794#endif   
795
796        if (ok_iso_verif) then
797           call check_isotopes(q,ijb_u,ije_u,'leapfrog 775')
798        endif !if (ok_iso_verif) then
799
800c      do j=1,nqtot
801c        call WriteField_p('q'//trim(int2str(j)),
802c     .                reshape(q(:,:,j),(/iip1,jmp1,llm/)))
803c        call WriteField_p('dq'//trim(int2str(j)),
804c     .                reshape(dq(:,:,j),(/iip1,jmp1,llm/)))
805c      enddo
806
807
808c$OMP MASTER
809       call VTe(VTintegre)
810c$OMP END MASTER
811c .P.Le Van (26/04/94  ajout de  finvpold dans l'appel d'integrd)
812c
813c-----------------------------------------------------------------------
814c   calcul des tendances physiques:
815c   -------------------------------
816c    ########   P.Le Van ( Modif le  6/02/95 )   ###########
817c
818       IF( purmats )  THEN
819          IF( itau.EQ.itaufin.AND..NOT.forward ) lafin = .TRUE.
820       ELSE
821          IF( itau+1. EQ. itaufin )              lafin = .TRUE.
822       ENDIF
823
824cc$OMP END PARALLEL
825
826c
827c
828       IF( apphys )  THEN
829       
830         CALL call_calfis(itau,lafin,ucov,vcov,teta,masse,ps, 
831     &                     phis,q,flxw)
832! #ifdef DEBUG_IO   
833!         call WriteField_u('ucovfi',ucov)
834!         call WriteField_v('vcovfi',vcov)
835!         call WriteField_u('tetafi',teta)
836!         call WriteField_u('pfi',p)
837!         call WriteField_u('pkfi',pk)
838!         do j=1,nqtot
839!           call WriteField_u('qfi'//trim(int2str(j)),q(:,:,j))
840!         enddo
841! #endif
842! c
843! c     .......   Ajout   P.Le Van ( 17/04/96 )   ...........
844! c
845! cc$OMP PARALLEL DEFAULT(SHARED)
846! cc$OMP+         PRIVATE(rdaym_ini,rdayvrai,ijb,ije)
847
848! c$OMP MASTER
849!          call suspend_timer(timer_caldyn)
850
851!          write(lunout,*)
852!      &   'leapfrog_p: Entree dans la physique : Iteration No ',true_itau
853! c$OMP END MASTER
854
855!          CALL pression_loc (  ip1jmp1, ap, bp, ps,  p      )
856
857! c$OMP BARRIER
858!          CALL exner_hyb_loc(  ip1jmp1, ps, p,pks, pk, pkf )
859! c$OMP BARRIER
860!            jD_cur = jD_ref + day_ini - day_ref
861!      $        + int (itau * dtvr / daysec)
862!            jH_cur = jH_ref +                                            &
863!      &              (itau * dtvr / daysec - int(itau * dtvr / daysec))
864! !         call ju2ymds(jD_cur+jH_cur, an, mois, jour, secondes)
865
866! c rajout debug
867! c       lafin = .true.
868
869
870! c   Inbterface avec les routines de phylmd (phymars ... )
871! c   -----------------------------------------------------
872
873! c+jld
874
875! c  Diagnostique de conservation de l'energie : initialisation
876
877! c-jld
878! c$OMP BARRIER
879! c$OMP MASTER
880!         call VTb(VThallo)
881! c$OMP END MASTER
882
883! #ifdef DEBUG_IO   
884!         call WriteField_u('ucovfi',ucov)
885!         call WriteField_v('vcovfi',vcov)
886!         call WriteField_u('tetafi',teta)
887!         call WriteField_u('pfi',p)
888!         call WriteField_u('pkfi',pk)
889! #endif
890!         call SetTag(Request_physic,800)
891!         
892!         call Register_SwapField_u(ucov,ucov,distrib_physic,
893!      *                            Request_physic,up=2,down=2)
894!         
895!         call Register_SwapField_v(vcov,vcov,distrib_physic,
896!      *                            Request_physic,up=2,down=2)
897
898!         call Register_SwapField_u(teta,teta,distrib_physic,
899!      *                            Request_physic,up=2,down=2)
900!         
901!         call Register_SwapField_u(masse,masse,distrib_physic,
902!      *                            Request_physic,up=1,down=2)
903
904!         call Register_SwapField_u(p,p,distrib_physic,
905!      *                            Request_physic,up=2,down=2)
906!         
907!         call Register_SwapField_u(pk,pk,distrib_physic,
908!      *                            Request_physic,up=2,down=2)
909!         
910!         call Register_SwapField_u(phis,phis,distrib_physic,
911!      *                            Request_physic,up=2,down=2)
912!         
913!         call Register_SwapField_u(phi,phi,distrib_physic,
914!      *                            Request_physic,up=2,down=2)
915!         
916!         call Register_SwapField_u(w,w,distrib_physic,
917!      *                            Request_physic,up=2,down=2)
918!         
919!         call Register_SwapField_u(q,q,distrib_physic,
920!      *                            Request_physic,up=2,down=2)
921
922!         call Register_SwapField_u(flxw,flxw,distrib_physic,
923!      *                            Request_physic,up=2,down=2)
924!         
925!         call SendRequest(Request_Physic)
926! c$OMP BARRIER
927!         call WaitRequest(Request_Physic)       
928
929! c$OMP BARRIER
930! c$OMP MASTER
931!         call Set_Distrib(distrib_Physic)
932!         call VTe(VThallo)
933!         
934!         call VTb(VTphysiq)
935! c$OMP END MASTER
936! c$OMP BARRIER
937
938! #ifdef DEBUG_IO   
939!       call WriteField_u('ucovfi',ucov)
940!       call WriteField_v('vcovfi',vcov)
941!       call WriteField_u('tetafi',teta)
942!       call WriteField_u('pfi',p)
943!       call WriteField_u('pkfi',pk)
944!       do j=1,nqtot
945!         call WriteField_u('qfi'//trim(int2str(j)),q(:,:,j))
946!       enddo
947! #endif
948!        STOP
949! c$OMP BARRIER
950! !        CALL FTRACE_REGION_BEGIN("calfis")
951!         CALL calfis_loc(lafin ,jD_cur, jH_cur,
952!      $               ucov,vcov,teta,q,masse,ps,p,pk,phis,phi ,
953!      $               du,dv,dteta,dq,
954!      $               flxw,
955!      $               dufi,dvfi,dtetafi,dqfi,dpfi  )
956! !        CALL FTRACE_REGION_END("calfis")
957! !        ijb=ij_begin
958! !        ije=ij_end 
959! !        if ( .not. pole_nord) then
960! !c$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
961! !          DO l=1,llm
962! !          dufi_tmp(1:iip1,l)   = dufi(ijb:ijb+iim,l)
963! !          dvfi_tmp(1:iip1,l)   = dvfi(ijb:ijb+iim,l) 
964! !          dtetafi_tmp(1:iip1,l)= dtetafi(ijb:ijb+iim,l) 
965! !          dqfi_tmp(1:iip1,l,:) = dqfi(ijb:ijb+iim,l,:) 
966! !          ENDDO
967! !c$OMP END DO NOWAIT
968! !
969! !c$OMP MASTER
970! !          dpfi_tmp(1:iip1)     = dpfi(ijb:ijb+iim) 
971! !c$OMP END MASTER
972! !        endif ! of if ( .not. pole_nord)
973
974! !c$OMP BARRIER
975! !c$OMP MASTER
976! !        call Set_Distrib(distrib_physic_bis)
977
978! !        call VTb(VThallo)
979! !c$OMP END MASTER
980! !c$OMP BARRIER
981! !
982! !        call Register_Hallo_u(dufi,llm,
983! !     *                      1,0,0,1,Request_physic)
984! !       
985! !        call Register_Hallo_v(dvfi,llm,
986! !     *                      1,0,0,1,Request_physic)
987! !       
988! !        call Register_Hallo_u(dtetafi,llm,
989! !     *                      1,0,0,1,Request_physic)
990! !
991! !        call Register_Hallo_u(dpfi,1,
992! !     *                      1,0,0,1,Request_physic)
993! !
994! !        do j=1,nqtot
995! !          call Register_Hallo_u(dqfi(ijb_u,1,j),llm,
996! !     *                        1,0,0,1,Request_physic)
997! !        enddo
998! !       
999! !        call SendRequest(Request_Physic)
1000! !c$OMP BARRIER
1001! !        call WaitRequest(Request_Physic)
1002! !             
1003! !c$OMP BARRIER
1004! !c$OMP MASTER
1005! !        call VTe(VThallo)
1006! !
1007! !        call set_Distrib(distrib_Physic)
1008! !c$OMP END MASTER
1009! !c$OMP BARRIER       
1010! !                ijb=ij_begin
1011! !        if (.not. pole_nord) then
1012! !       
1013! !c$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
1014! !          DO l=1,llm
1015! !            dufi(ijb:ijb+iim,l) = dufi(ijb:ijb+iim,l)+dufi_tmp(1:iip1,l)
1016! !            dvfi(ijb:ijb+iim,l) = dvfi(ijb:ijb+iim,l)+dvfi_tmp(1:iip1,l)
1017! !            dtetafi(ijb:ijb+iim,l) = dtetafi(ijb:ijb+iim,l)
1018! !     &                              +dtetafi_tmp(1:iip1,l)
1019! !            dqfi(ijb:ijb+iim,l,:) = dqfi(ijb:ijb+iim,l,:)
1020! !     &                              + dqfi_tmp(1:iip1,l,:)
1021! !          ENDDO
1022! !c$OMP END DO NOWAIT
1023! !
1024! !c$OMP MASTER
1025! !          dpfi(ijb:ijb+iim)   = dpfi(ijb:ijb+iim)+ dpfi_tmp(1:iip1)
1026! !c$OMP END MASTER
1027! !         
1028! !        endif ! of if (.not. pole_nord)
1029
1030! #ifdef DEBUG_IO           
1031!         call WriteField_u('dufi',dufi)
1032!         call WriteField_v('dvfi',dvfi)
1033!         call WriteField_u('dtetafi',dtetafi)
1034!         call WriteField_u('dpfi',dpfi)
1035!         do j=1,nqtot
1036!           call WriteField_u('dqfi'//trim(int2str(j)),dqfi(:,:,j))
1037!        enddo
1038! #endif
1039
1040! c$OMP BARRIER
1041
1042! c      ajout des tendances physiques:
1043! c      ------------------------------
1044! #ifdef DEBUG_IO   
1045!         call WriteField_u('ucovfi',ucov)
1046!         call WriteField_v('vcovfi',vcov)
1047!         call WriteField_u('tetafi',teta)
1048!         call WriteField_u('psfi',ps)
1049!         do j=1,nqtot
1050!           call WriteField_u('qfi'//trim(int2str(j)),q(:,:,j))
1051!        enddo
1052! #endif
1053
1054!          IF (ok_strato) THEN
1055!            CALL top_bound_loc( vcov,ucov,teta,masse,dufi,dvfi,dtetafi)
1056!          ENDIF
1057
1058! #ifdef DEBUG_IO           
1059!         call WriteField_u('ucovfi',ucov)
1060!         call WriteField_v('vcovfi',vcov)
1061!         call WriteField_u('tetafi',teta)
1062!         call WriteField_u('psfi',ps)
1063!         do j=1,nqtot
1064!           call WriteField_u('qfi'//trim(int2str(j)),q(:,:,j))
1065!        enddo
1066! #endif
1067
1068!           CALL addfi_loc( dtphys, leapf, forward   ,
1069!      $                  ucov, vcov, teta , q   ,ps ,
1070!      $                 dufi, dvfi, dtetafi , dqfi ,dpfi  )
1071
1072! #ifdef DEBUG_IO   
1073!         call WriteField_u('ucovfi',ucov)
1074!         call WriteField_v('vcovfi',vcov)
1075!         call WriteField_u('tetafi',teta)
1076!         call WriteField_u('psfi',ps)
1077!         do j=1,nqtot
1078!           call WriteField_u('qfi'//trim(int2str(j)),q(:,:,j))
1079!        enddo
1080! #endif
1081
1082! c$OMP BARRIER
1083! c$OMP MASTER
1084!         call VTe(VTphysiq)
1085
1086!         call VTb(VThallo)
1087! c$OMP END MASTER
1088
1089!         call SetTag(Request_physic,800)
1090!         call Register_SwapField_u(ucov,ucov,
1091!      *                               distrib_caldyn,Request_physic)
1092!         
1093!         call Register_SwapField_v(vcov,vcov,
1094!      *                               distrib_caldyn,Request_physic)
1095!         
1096!         call Register_SwapField_u(teta,teta,
1097!      *                               distrib_caldyn,Request_physic)
1098!         
1099!         call Register_SwapField_u(masse,masse,
1100!      *                               distrib_caldyn,Request_physic)
1101
1102!         call Register_SwapField_u(p,p,
1103!      *                               distrib_caldyn,Request_physic)
1104!         
1105!         call Register_SwapField_u(pk,pk,
1106!      *                               distrib_caldyn,Request_physic)
1107!         
1108!         call Register_SwapField_u(phis,phis,
1109!      *                               distrib_caldyn,Request_physic)
1110!         
1111!         call Register_SwapField_u(phi,phi,
1112!      *                               distrib_caldyn,Request_physic)
1113!         
1114!         call Register_SwapField_u(w,w,
1115!      *                               distrib_caldyn,Request_physic)
1116
1117!         call Register_SwapField_u(q,q,
1118!      *                               distrib_caldyn,Request_physic)
1119!         
1120!         call SendRequest(Request_Physic)
1121! c$OMP BARRIER
1122!         call WaitRequest(Request_Physic)     
1123
1124! c$OMP BARRIER
1125! c$OMP MASTER
1126!        call VTe(VThallo)
1127!        call set_distrib(distrib_caldyn)
1128! c$OMP END MASTER
1129! c$OMP BARRIER
1130! c
1131! c  Diagnostique de conservation de l'energie : difference
1132!       IF (ip_ebil_dyn.ge.1 ) THEN
1133!           ztit='bil phys'
1134!           CALL diagedyn(ztit,2,1,1,dtphys
1135!      e  , ucov    , vcov , ps, p ,pk , teta , q(:,:,1), q(:,:,2))
1136!       ENDIF
1137
1138! #ifdef DEBUG_IO   
1139!         call WriteField_u('ucovfi',ucov)
1140!         call WriteField_v('vcovfi',vcov)
1141!         call WriteField_u('tetafi',teta)
1142!         call WriteField_u('psfi',ps)
1143!         do j=1,nqtot
1144!           call WriteField_u('qfi'//trim(int2str(j)),q(:,:,j))
1145!        enddo
1146! #endif
1147
1148
1149! c-jld
1150c$OMP MASTER
1151         if (FirstPhysic) then
1152           ok_start_timer=.TRUE.
1153           FirstPhysic=.false.
1154         endif
1155c$OMP END MASTER
1156       ENDIF ! of IF( apphys )
1157
1158        if (ok_iso_verif) then
1159           call check_isotopes(q,ijb_u,ije_u,'leapfrog 1132')
1160        endif !if (ok_iso_verif) then
1161        !write(*,*) 'leapfrog 1134: iflag_phys=',iflag_phys
1162
1163      IF(iflag_phys.EQ.2) THEN ! "Newtonian" case
1164c$OMP MASTER
1165         if (FirstPhysic) then
1166           ok_start_timer=.TRUE.
1167           FirstPhysic=.false.
1168         endif
1169c$OMP END MASTER
1170
1171
1172c   Calcul academique de la physique = Rappel Newtonien + fritcion
1173c   --------------------------------------------------------------
1174cym       teta(:,:)=teta(:,:)
1175cym     s  -iphysiq*dtvr*(teta(:,:)-tetarappel(:,:))/taurappel
1176       ijb=ij_begin
1177       ije=ij_end
1178!LF       teta(ijb:ije,:)=teta(ijb:ije,:)
1179!LF     s  -iphysiq*dtvr*(teta(ijb:ije,:)-tetarappel(ijb:ije,:))/taurappel
1180!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
1181       do l=1,llm
1182       teta(ijb:ije,l)=teta(ijb:ije,l) -dtvr*
1183     &        (teta(ijb:ije,l)-tetarappel(ijb:ije,l))*
1184     &                 (knewt_g+knewt_t(l)*clat4(ijb:ije))       
1185       enddo
1186!$OMP END DO
1187
1188!$OMP MASTER
1189       if (planet_type.eq."giant") then
1190         ! add an intrinsic heat flux at the base of the atmosphere
1191         teta(ijb:ije,1) = teta(ijb:ije,1)
1192     &        + dtvr * aire(ijb:ije) * ihf / cpp / masse(ijb:ije,1)
1193       endif
1194!$OMP END MASTER
1195!$OMP BARRIER
1196
1197
1198       call Register_Hallo_u(ucov,llm,0,1,1,0,Request_Physic)
1199       call Register_Hallo_v(vcov,llm,1,1,1,1,Request_Physic)
1200       call SendRequest(Request_Physic)
1201c$OMP BARRIER
1202       call WaitRequest(Request_Physic)     
1203c$OMP BARRIER
1204       call friction_loc(ucov,vcov,dtvr)
1205!$OMP BARRIER
1206
1207        ! Sponge layer (if any)
1208        IF (ok_strato) THEN
1209          CALL top_bound_loc(vcov,ucov,teta,masse,dtvr)
1210!$OMP BARRIER
1211        ENDIF ! of IF (ok_strato)
1212      ENDIF ! of IF(iflag_phys.EQ.2)
1213
1214
1215        CALL pression_loc ( ip1jmp1, ap, bp, ps, p                  )
1216c$OMP BARRIER
1217        if (pressure_exner) then
1218        CALL exner_hyb_loc( ijnb_u, ps, p, pks, pk, pkf )
1219        else
1220          CALL exner_milieu_loc( ijnb_u, ps, p, pks, pk, pkf )
1221        endif
1222c$OMP BARRIER
1223        CALL massdair_loc(p,masse)
1224c$OMP BARRIER
1225
1226cc$OMP END PARALLEL
1227        if (ok_iso_verif) then
1228           call check_isotopes(q,ijb_u,ije_u,'leapfrog 1196')
1229        endif !if (ok_iso_verif) then
1230
1231c-----------------------------------------------------------------------
1232c   dissipation horizontale et verticale  des petites echelles:
1233c   ----------------------------------------------------------
1234      !write(*,*) 'leapfrog 1163: apdiss=',apdiss
1235      IF(apdiss) THEN
1236     
1237        CALL call_dissip(ucov,vcov,teta,p,pk,ps)
1238!cc$OMP  PARALLEL DEFAULT(SHARED)
1239!cc$OMP+          PRIVATE(ijb,ije,tppn,tpn,tpps,tps)
1240!c$OMP MASTER
1241!        call suspend_timer(timer_caldyn)
1242!       
1243!c       print*,'Entree dans la dissipation : Iteration No ',true_itau
1244!c   calcul de l'energie cinetique avant dissipation
1245!c       print *,'Passage dans la dissipation'
1246
1247!        call VTb(VThallo)
1248!c$OMP END MASTER
1249
1250!c$OMP BARRIER
1251
1252!        call Register_SwapField_u(ucov,ucov,distrib_dissip,
1253!     *                            Request_dissip,up=1,down=1)
1254
1255!        call Register_SwapField_v(vcov,vcov,distrib_dissip,
1256!     *                            Request_dissip,up=1,down=1)
1257
1258!        call Register_SwapField_u(teta,teta,distrib_dissip,
1259!     *                            Request_dissip)
1260
1261!        call Register_SwapField_u(p,p,distrib_dissip,
1262!     *                            Request_dissip)
1263
1264!        call Register_SwapField_u(pk,pk,distrib_dissip,
1265!     *                            Request_dissip)
1266
1267!        call SendRequest(Request_dissip)       
1268!c$OMP BARRIER
1269!        call WaitRequest(Request_dissip)       
1270
1271!c$OMP BARRIER
1272!c$OMP MASTER
1273!        call set_distrib(distrib_dissip)
1274!        call VTe(VThallo)
1275!        call VTb(VTdissipation)
1276!        call start_timer(timer_dissip)
1277!c$OMP END MASTER
1278!c$OMP BARRIER
1279
1280!        call covcont_loc(llm,ucov,vcov,ucont,vcont)
1281!        call enercin_loc(vcov,ucov,vcont,ucont,ecin0)
1282
1283!c   dissipation
1284
1285!!        CALL FTRACE_REGION_BEGIN("dissip")
1286!        CALL dissip_loc(vcov,ucov,teta,p,dvdis,dudis,dtetadis)
1287
1288!#ifdef DEBUG_IO   
1289!        call WriteField_u('dudis',dudis)
1290!        call WriteField_v('dvdis',dvdis)
1291!        call WriteField_u('dtetadis',dtetadis)
1292!#endif
1293!
1294!!      CALL FTRACE_REGION_END("dissip")
1295!         
1296!        ijb=ij_begin
1297!        ije=ij_end
1298!c$OMP DO SCHEDULE(STATIC,OMP_CHUNK)       
1299!        DO l=1,llm
1300!          ucov(ijb:ije,l)=ucov(ijb:ije,l)+dudis(ijb:ije,l)
1301!        ENDDO
1302!c$OMP END DO NOWAIT       
1303!        if (pole_sud) ije=ije-iip1
1304!c$OMP DO SCHEDULE(STATIC,OMP_CHUNK)       
1305!        DO l=1,llm
1306!          vcov(ijb:ije,l)=vcov(ijb:ije,l)+dvdis(ijb:ije,l)
1307!        ENDDO
1308!c$OMP END DO NOWAIT       
1309
1310!c       teta=teta+dtetadis
1311
1312
1313!c------------------------------------------------------------------------
1314!        if (dissip_conservative) then
1315!C       On rajoute la tendance due a la transform. Ec -> E therm. cree
1316!C       lors de la dissipation
1317!c$OMP BARRIER
1318!c$OMP MASTER
1319!            call suspend_timer(timer_dissip)
1320!            call VTb(VThallo)
1321!c$OMP END MASTER
1322!            call Register_Hallo_u(ucov,llm,1,1,1,1,Request_Dissip)
1323!            call Register_Hallo_v(vcov,llm,1,1,1,1,Request_Dissip)
1324!            call SendRequest(Request_Dissip)
1325!c$OMP BARRIER
1326!            call WaitRequest(Request_Dissip)
1327!c$OMP MASTER
1328!            call VTe(VThallo)
1329!            call resume_timer(timer_dissip)
1330!c$OMP END MASTER
1331!c$OMP BARRIER           
1332!            call covcont_loc(llm,ucov,vcov,ucont,vcont)
1333!            call enercin_loc(vcov,ucov,vcont,ucont,ecin)
1334!           
1335!            ijb=ij_begin
1336!            ije=ij_end
1337!c$OMP DO SCHEDULE(STATIC,OMP_CHUNK)           
1338!            do l=1,llm
1339!              do ij=ijb,ije
1340!                dtetaecdt(ij,l)= (ecin0(ij,l)-ecin(ij,l))/ pk(ij,l)
1341!                dtetadis(ij,l)=dtetadis(ij,l)+dtetaecdt(ij,l)
1342!              enddo
1343!            enddo
1344!c$OMP END DO NOWAIT           
1345!       endif
1346
1347!       ijb=ij_begin
1348!       ije=ij_end
1349!c$OMP DO SCHEDULE(STATIC,OMP_CHUNK)           
1350!         do l=1,llm
1351!           do ij=ijb,ije
1352!              teta(ij,l)=teta(ij,l)+dtetadis(ij,l)
1353!           enddo
1354!         enddo
1355!c$OMP END DO NOWAIT         
1356!c------------------------------------------------------------------------
1357
1358
1359!c    .......        P. Le Van (  ajout  le 17/04/96  )   ...........
1360!c   ...      Calcul de la valeur moyenne, unique de h aux poles  .....
1361!c
1362
1363!        ijb=ij_begin
1364!        ije=ij_end
1365!         
1366!        if (pole_nord) then
1367!c$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
1368!          DO l  =  1, llm
1369!            DO ij =  1,iim
1370!             tppn(ij)  = aire(  ij    ) * teta(  ij    ,l)
1371!            ENDDO
1372!             tpn  = SSUM(iim,tppn,1)/apoln
1373
1374!            DO ij = 1, iip1
1375!             teta(  ij    ,l) = tpn
1376!            ENDDO
1377!          ENDDO
1378!c$OMP END DO NOWAIT
1379
1380!c$OMP MASTER               
1381!          DO ij =  1,iim
1382!            tppn(ij)  = aire(  ij    ) * ps (  ij    )
1383!          ENDDO
1384!            tpn  = SSUM(iim,tppn,1)/apoln
1385
1386!          DO ij = 1, iip1
1387!            ps(  ij    ) = tpn
1388!          ENDDO
1389!c$OMP END MASTER
1390!        endif
1391!       
1392!        if (pole_sud) then
1393!c$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
1394!          DO l  =  1, llm
1395!            DO ij =  1,iim
1396!             tpps(ij)  = aire(ij+ip1jm) * teta(ij+ip1jm,l)
1397!            ENDDO
1398!             tps  = SSUM(iim,tpps,1)/apols
1399
1400!            DO ij = 1, iip1
1401!             teta(ij+ip1jm,l) = tps
1402!            ENDDO
1403!          ENDDO
1404!c$OMP END DO NOWAIT
1405
1406!c$OMP MASTER               
1407!          DO ij =  1,iim
1408!            tpps(ij)  = aire(ij+ip1jm) * ps (ij+ip1jm)
1409!          ENDDO
1410!            tps  = SSUM(iim,tpps,1)/apols
1411
1412!          DO ij = 1, iip1
1413!            ps(ij+ip1jm) = tps
1414!          ENDDO
1415!c$OMP END MASTER
1416!        endif
1417
1418
1419!c$OMP BARRIER
1420!c$OMP MASTER
1421!        call VTe(VTdissipation)
1422
1423!        call stop_timer(timer_dissip)
1424!       
1425!        call VTb(VThallo)
1426!c$OMP END MASTER
1427!        call Register_SwapField_u(ucov,ucov,distrib_caldyn,
1428!     *                            Request_dissip)
1429
1430!        call Register_SwapField_v(vcov,vcov,distrib_caldyn,
1431!     *                            Request_dissip)
1432
1433!        call Register_SwapField_u(teta,teta,distrib_caldyn,
1434!     *                            Request_dissip)
1435
1436!        call Register_SwapField_u(p,p,distrib_caldyn,
1437!     *                            Request_dissip)
1438
1439!        call Register_SwapField_u(pk,pk,distrib_caldyn,
1440!     *                            Request_dissip)
1441
1442!        call SendRequest(Request_dissip)       
1443!c$OMP BARRIER
1444!        call WaitRequest(Request_dissip)       
1445
1446!c$OMP BARRIER
1447!c$OMP MASTER
1448!        call set_distrib(distrib_caldyn)
1449!        call VTe(VThallo)
1450!        call resume_timer(timer_caldyn)
1451!c        print *,'fin dissipation'
1452!c$OMP END MASTER
1453!c$OMP BARRIER
1454       END IF ! of IF(apdiss)
1455
1456cc$OMP END PARALLEL
1457
1458c ajout debug
1459c              IF( lafin ) then 
1460c                abort_message = 'Simulation finished'
1461c                call abort_gcm(modname,abort_message,0)
1462c              ENDIF
1463
1464        if (ok_iso_verif) then
1465           call check_isotopes(q,ijb_u,ije_u,'leapfrog 1430')
1466        endif !if (ok_iso_verif) then     
1467 
1468c   ********************************************************************
1469c   ********************************************************************
1470c   .... fin de l'integration dynamique  et physique pour le pas itau ..
1471c   ********************************************************************
1472c   ********************************************************************
1473
1474c   preparation du pas d'integration suivant  ......
1475cym      call WriteField('ucov',reshape(ucov,(/iip1,jmp1,llm/)))
1476cym      call WriteField('vcov',reshape(vcov,(/iip1,jjm,llm/)))
1477c$OMP MASTER     
1478      call stop_timer(timer_caldyn)
1479c$OMP END MASTER
1480      IF (itau==itaumax) then
1481c$OMP MASTER
1482         call allgather_timer_average
1483         call barrier
1484         if (mpi_rank==0) then
1485           
1486            print *,'*********************************'
1487            print *,'******    TIMER CALDYN     ******'
1488            do i=0,mpi_size-1
1489               print *,'proc',i,' :   Nb Bandes  :',jj_nb_caldyn(i),
1490     &              '  : temps moyen :',
1491     &              timer_average(jj_nb_caldyn(i),timer_caldyn,i)
1492            enddo
1493           
1494            print *,'*********************************'
1495            print *,'******    TIMER VANLEER    ******'
1496            do i=0,mpi_size-1
1497               print *,'proc',i,' :   Nb Bandes  :',jj_nb_vanleer(i),
1498     &              '  : temps moyen :',
1499     &              timer_average(jj_nb_vanleer(i),timer_vanleer,i)
1500            enddo
1501           
1502            print *,'*********************************'
1503            print *,'******    TIMER DISSIP    ******'
1504            do i=0,mpi_size-1
1505               print *,'proc',i,' :   Nb Bandes  :',jj_nb_dissip(i),
1506     &              '  : temps moyen :',
1507     &              timer_average(jj_nb_dissip(i),timer_dissip,i)
1508            enddo
1509           
1510            print *,'*********************************'
1511            print *,'******    TIMER PHYSIC    ******'
1512            do i=0,mpi_size-1
1513               print *,'proc',i,' :   Nb Bandes  :',jj_nb_physic(i),
1514     &              '  : temps moyen :',
1515     &              timer_average(jj_nb_physic(i),timer_physic,i)
1516            enddo
1517           
1518         endif 
1519         CALL barrier
1520         print *,'Taille du Buffer MPI (REAL*8)',MaxBufferSize
1521      print *,'Taille du Buffer MPI utilise (REAL*8)',MaxBufferSize_Used
1522       print *, 'Temps total ecoule sur la parallelisation :',DiffTime()
1523      print *, 'Temps CPU ecoule sur la parallelisation :',DiffCpuTime()
1524         CALL print_filtre_timer
1525c$OMP END MASTER
1526         CALL dynredem1_loc("restart.nc",0.0,
1527     .        vcov,ucov,teta,q,masse,ps)
1528c$OMP MASTER
1529         call fin_getparam
1530c$OMP END MASTER
1531
1532#ifdef INCA
1533         call finalize_inca
1534#endif
1535
1536c$OMP MASTER
1537         call finalize_parallel
1538c$OMP END MASTER
1539c$OMP BARRIER
1540         RETURN
1541      ENDIF
1542     
1543        if (ok_iso_verif) then
1544           call check_isotopes(q,ijb_u,ije_u,'leapfrog 1509')
1545        endif !if (ok_iso_verif) then
1546
1547      IF ( .NOT.purmats ) THEN
1548c       ........................................................
1549c       ..............  schema matsuno + leapfrog  ..............
1550c       ........................................................
1551
1552            IF(forward. OR. leapf) THEN
1553              itau= itau + 1
1554!              iday= day_ini+itau/day_step
1555!              time= REAL(itau-(iday-day_ini)*day_step)/day_step+time_0
1556!                IF(time.GT.1.) THEN
1557!                  time = time-1.
1558!                  iday = iday+1
1559!                ENDIF
1560            ENDIF
1561
1562
1563            IF( itau. EQ. itaufinp1 ) then
1564
1565              if (flag_verif) then
1566                write(79,*) 'ucov',ucov
1567                write(80,*) 'vcov',vcov
1568                write(81,*) 'teta',teta
1569                write(82,*) 'ps',ps
1570                write(83,*) 'q',q
1571                WRITE(85,*) 'q1 = ',q(:,:,1)
1572                WRITE(86,*) 'q3 = ',q(:,:,3)
1573              endif
1574 
1575
1576c$OMP MASTER
1577              call fin_getparam
1578c$OMP END MASTER
1579
1580#ifdef INCA
1581              call finalize_inca
1582#endif
1583
1584c$OMP MASTER
1585              call finalize_parallel
1586c$OMP END MASTER
1587              abort_message = 'Simulation finished'
1588              call abort_gcm(modname,abort_message,0)
1589              RETURN
1590            ENDIF
1591c-----------------------------------------------------------------------
1592c   ecriture du fichier histoire moyenne:
1593c   -------------------------------------
1594
1595            IF(MOD(itau,iperiod).EQ.0 .OR. itau.EQ.itaufin) THEN
1596c$OMP BARRIER
1597               IF(itau.EQ.itaufin) THEN
1598                  iav=1
1599               ELSE
1600                  iav=0
1601               ENDIF
1602
1603              ! Ehouarn: re-compute geopotential for outputs
1604c$OMP BARRIER
1605c$OMP MASTER
1606              CALL geopot_loc(ip1jmp1,teta,pk,pks,phis,phi)
1607c$OMP END MASTER
1608c$OMP BARRIER
1609
1610#ifdef CPP_IOIPSL
1611             IF (ok_dynzon) THEN
1612
1613              CALL bilan_dyn_loc(2,dtvr*iperiod,dtvr*day_step*periodav,
1614     ,             ps,masse,pk,pbaru,pbarv,teta,phi,ucov,vcov,q)
1615
1616              ENDIF !ok_dynzon
1617
1618              IF (ok_dyn_ave) THEN
1619                 CALL writedynav_loc(itau,vcov,
1620     &                 ucov,teta,pk,phi,q,masse,ps,phis)
1621              ENDIF
1622#endif
1623            ENDIF
1624
1625        if (ok_iso_verif) then
1626           call check_isotopes(q,ijb_u,ije_u,'leapfrog 1584')
1627        endif !if (ok_iso_verif) then
1628
1629c-----------------------------------------------------------------------
1630c   ecriture de la bande histoire:
1631c   ------------------------------
1632
1633            IF( MOD(itau,iecri).EQ.0) THEN
1634             ! Ehouarn: output only during LF or Backward Matsuno
1635             if (leapf.or.(.not.leapf.and.(.not.forward))) then
1636
1637c$OMP BARRIER
1638c$OMP MASTER
1639              CALL geopot_loc(ip1jmp1,teta,pk,pks,phis,phi)
1640c$OMP END MASTER
1641c$OMP BARRIER
1642       
1643#ifdef CPP_IOIPSL
1644             if (ok_dyn_ins) then
1645                 CALL writehist_loc(itau,vcov,ucov,teta,pk,phi,q,
1646     &                              masse,ps,phis)
1647             endif
1648#endif
1649            endif ! of if (leapf.or.(.not.leapf.and.(.not.forward)))
1650           ENDIF ! of IF(MOD(itau,iecri).EQ.0)
1651
1652            IF(itau.EQ.itaufin) THEN
1653
1654c$OMP BARRIER
1655
1656!              if (planet_type.eq."earth") then
1657! Write an Earth-format restart file
1658                CALL dynredem1_loc("restart.nc",0.0,
1659     &                           vcov,ucov,teta,q,masse,ps)
1660!              endif ! of if (planet_type.eq."earth")
1661
1662!              CLOSE(99)
1663            ENDIF ! of IF (itau.EQ.itaufin)
1664
1665        if (ok_iso_verif) then
1666           call check_isotopes(q,ijb_u,ije_u,'leapfrog 1624')
1667        endif !if (ok_iso_verif) then
1668
1669c-----------------------------------------------------------------------
1670c   gestion de l'integration temporelle:
1671c   ------------------------------------
1672
1673            IF( MOD(itau,iperiod).EQ.0 )    THEN
1674                    GO TO 1
1675            ELSE IF ( MOD(itau-1,iperiod). EQ. 0 ) THEN
1676
1677                   IF( forward )  THEN
1678c      fin du pas forward et debut du pas backward
1679
1680                      forward = .FALSE.
1681                        leapf = .FALSE.
1682                           GO TO 2
1683
1684                   ELSE
1685c      fin du pas backward et debut du premier pas leapfrog
1686
1687                        leapf =  .TRUE.
1688                        dt  =  2.*dtvr
1689                        GO TO 2
1690                   END IF
1691            ELSE
1692
1693c      ......   pas leapfrog  .....
1694
1695                 leapf = .TRUE.
1696                 dt  = 2.*dtvr
1697                 GO TO 2
1698            END IF ! of IF (MOD(itau,iperiod).EQ.0)
1699                   !    ELSEIF (MOD(itau-1,iperiod).EQ.0)
1700
1701
1702      ELSE ! of IF (.not.purmats)
1703
1704
1705        if (ok_iso_verif) then
1706           call check_isotopes(q,ijb_u,ije_u,'leapfrog 1664')
1707        endif !if (ok_iso_verif) then
1708
1709c       ........................................................
1710c       ..............       schema  matsuno        ...............
1711c       ........................................................
1712            IF( forward )  THEN
1713
1714             itau =  itau + 1
1715!             iday = day_ini+itau/day_step
1716!             time = REAL(itau-(iday-day_ini)*day_step)/day_step+time_0
1717!
1718!                  IF(time.GT.1.) THEN
1719!                   time = time-1.
1720!                   iday = iday+1
1721!                  ENDIF
1722
1723               forward =  .FALSE.
1724               IF( itau. EQ. itaufinp1 ) then 
1725c$OMP MASTER
1726                 call fin_getparam
1727c$OMP END MASTER
1728
1729#ifdef INCA
1730                 call finalize_inca
1731#endif
1732
1733c$OMP MASTER
1734                 call finalize_parallel
1735c$OMP END MASTER
1736                 abort_message = 'Simulation finished'
1737                 call abort_gcm(modname,abort_message,0)
1738                 RETURN
1739               ENDIF
1740               GO TO 2
1741
1742            ELSE ! of IF(forward) i.e. backward step
1743
1744             
1745        if (ok_iso_verif) then
1746           call check_isotopes(q,ijb_u,ije_u,'leapfrog 1698')
1747        endif !if (ok_iso_verif) then 
1748
1749              IF(MOD(itau,iperiod).EQ.0 .OR. itau.EQ.itaufin) THEN
1750               IF(itau.EQ.itaufin) THEN
1751                  iav=1
1752               ELSE
1753                  iav=0
1754               ENDIF
1755
1756#ifdef CPP_IOIPSL
1757              ! Ehouarn: re-compute geopotential for outputs
1758c$OMP BARRIER
1759c$OMP MASTER
1760              CALL geopot_loc(ip1jmp1,teta,pk,pks,phis,phi)
1761c$OMP END MASTER
1762c$OMP BARRIER
1763               
1764               IF (ok_dynzon) THEN
1765               CALL bilan_dyn_loc(2,dtvr*iperiod,dtvr*day_step*periodav,
1766     ,           ps,masse,pk,pbaru,pbarv,teta,phi,ucov,vcov,q)
1767               ENDIF
1768             
1769               IF (ok_dyn_ave) THEN
1770                 CALL writedynav_loc(itau,vcov,
1771     &                 ucov,teta,pk,phi,q,masse,ps,phis)
1772               ENDIF
1773#endif
1774              ENDIF ! of IF(MOD(itau,iperiod).EQ.0 .OR. itau.EQ.itaufin)
1775
1776
1777               IF(MOD(itau,iecri         ).EQ.0) THEN
1778
1779c$OMP BARRIER
1780c$OMP MASTER
1781              CALL geopot_loc(ip1jmp1,teta,pk,pks,phis,phi)
1782c$OMP END MASTER
1783c$OMP BARRIER
1784
1785
1786#ifdef CPP_IOIPSL
1787              if (ok_dyn_ins) then
1788                 CALL writehist_loc(itau,vcov,ucov,teta,pk,phi,q,
1789     &                              masse,ps,phis)
1790              endif ! of if (ok_dyn_ins)
1791#endif
1792              ENDIF ! of IF(MOD(itau,iecri).EQ.0)
1793             
1794
1795              IF(itau.EQ.itaufin) THEN
1796!                if (planet_type.eq."earth") then
1797                   CALL dynredem1_loc("restart.nc",0.0,
1798     .                               vcov,ucov,teta,q,masse,ps)
1799!               endif ! of if (planet_type.eq."earth")
1800              ENDIF ! of IF(itau.EQ.itaufin)
1801
1802              forward = .TRUE.
1803              GO TO  1
1804
1805            ENDIF ! of IF (forward)
1806
1807
1808        if (ok_iso_verif) then
1809           call check_isotopes(q,ijb_u,ije_u,'leapfrog 1750')
1810        endif !if (ok_iso_verif) then
1811
1812      END IF ! of IF(.not.purmats)
1813c$OMP MASTER
1814      call fin_getparam
1815c$OMP END MASTER
1816
1817#ifdef INCA
1818      call finalize_inca
1819#endif
1820
1821c$OMP MASTER
1822      call finalize_parallel
1823c$OMP END MASTER
1824      abort_message = 'Simulation finished'
1825      call abort_gcm(modname,abort_message,0)
1826      RETURN
1827      END
Note: See TracBrowser for help on using the repository browser.