source: LMDZ5/trunk/libf/dyn3dpar/leapfrog_p.F @ 2038

Last change on this file since 2038 was 2038, checked in by fhourdin, 10 years ago

Possibilité d'imposer le nombre de pas de temps de la simulation avec nday<0
If nday<0, nday is the total number of time step of the simulation

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