source: trunk/LMDZ.COMMON/libf/dyn3d/leapfrog.F @ 1572

Last change on this file since 1572 was 1572, checked in by emillour, 8 years ago

All GCMs:
Further adaptations to keep up with changes in LMDZ5 concerning
physics/dynamics separation (up to rev r2500 of LMDZ5)

  • arch:
  • remove ifort debug option '-check all', replace it with '-check bounds,format,output_conversion,pointers,uninit' (i.e. get it to stop complaining about copying into temporary arrays)
  • dyn3d_common:
  • comconst_mod.F90 : add ngroup
  • dyn3d:
  • gcm.F90 : minor bug fix (arguments to a call_abort())
  • leapfrog.F90 : recompute geopotential for bilan_dyn outputs
  • conf_gcm.F90 : read "ngroup" from run.def
  • groupe.F , groupeun.F : ngroup no longer a local parameter
  • dyn3d_par:
  • conf_gcm.F90 : read "ngroup" from run.def
  • groupe_p.F , groupeun_p.F : ngroup no longer a local parameter
  • misc:
  • regr1_step_av_m.F90 : removed (not used)
  • phy_common:
  • mod_phys_lmdz_mpi_transfert.F90 , mod_phys_lmdz_mpi_data.F90 : change is_north_pole and is_south_pole to is_north_pole_dyn and is_south_pole_dyn
  • mod_phys_lmdz_omp_data.F90 : introduce is_nort_pole_phy and is_south_pole_phy
  • dynphy_lonlat:
  • mod_interface_dyn_phys.F90 : use is_north_pole_dyn and is_south_pole_dyn
  • calfis_p.F : use is_north_pole_dyn and is_south_pole_dyn
  • phyvenus:
  • physiq_mod , write_hist*.h : use is_north_pole_phy and is_south_pole_phy to correctly compute mesh area at poles to send to hist*nc files.
  • phytitan:
  • physiq_mod , write_hist*.h : use is_north_pole_phy and is_south_pole_phy to correctly compute mesh area at poles to send to hist*nc files.

EM

File size: 32.8 KB
Line 
1!
2! $Id: leapfrog.F 1446 2010-10-22 09:27:25Z emillour $
3!
4c
5c
6      SUBROUTINE leapfrog(ucov,vcov,teta,ps,masse,phis,q,
7     &                    time_0)
8
9
10cIM : pour sortir les param. du modele dans un fis. netcdf 110106
11#ifdef CPP_IOIPSL
12      use IOIPSL
13#endif
14      USE infotrac, ONLY: nqtot,ok_iso_verif
15      USE guide_mod, ONLY : guide_main
16      USE write_field, ONLY: writefield
17      USE control_mod, ONLY: planet_type,nday,day_step,iperiod,iphysiq,
18     &                       less1day,fractday,ndynstep,iconser,
19     &                       dissip_period,offline,ip_ebil_dyn,
20     &                       ok_dynzon,periodav,ok_dyn_ave,iecri,
21     &                       ok_dyn_ins,output_grads_dyn
22      use exner_hyb_m, only: exner_hyb
23      use exner_milieu_m, only: exner_milieu
24      use cpdet_mod, only: cpdet,tpot2t,t2tpot
25      use sponge_mod, only: callsponge,mode_sponge,sponge
26       use comuforc_h
27      USE comvert_mod, ONLY: ap,bp,pressure_exner,presnivs
28      USE comconst_mod, ONLY: daysec,dtvr,dtphys,dtdiss,
29     .                  cpp,ihf,iflag_top_bound,pi
30      USE logic_mod, ONLY: iflag_phys,ok_guide,forward,leapf,apphys,
31     .                  statcl,conser,apdiss,purmats,tidal,ok_strato
32      USE temps_mod, ONLY: jD_ref,jH_ref,itaufin,day_ini,day_ref,
33     .                  start_time,dt
34
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 "comgeom.h"
72!#include "com_io_dyn.h"
73#include "iniprint.h"
74#include "academic.h"
75
76! FH 2008/05/09 On elimine toutes les clefs physiques dans la dynamique
77! #include "clesphys.h"
78
79      REAL,INTENT(IN) :: time_0 ! not used
80
81c   dynamical variables:
82      REAL,INTENT(INOUT) :: ucov(ip1jmp1,llm)    ! zonal covariant wind
83      REAL,INTENT(INOUT) :: vcov(ip1jm,llm)      ! meridional covariant wind
84      REAL,INTENT(INOUT) :: teta(ip1jmp1,llm)    ! potential temperature
85      REAL,INTENT(INOUT) :: ps(ip1jmp1)          ! surface pressure (Pa)
86      REAL,INTENT(INOUT) :: masse(ip1jmp1,llm)   ! air mass
87      REAL,INTENT(INOUT) :: phis(ip1jmp1)        ! geopotentiat at the surface
88      REAL,INTENT(INOUT) :: q(ip1jmp1,llm,nqtot) ! advected tracers
89
90      REAL p (ip1jmp1,llmp1  )               ! interlayer pressure
91      REAL pks(ip1jmp1)                      ! exner at the surface
92      REAL pk(ip1jmp1,llm)                   ! exner at mid-layer
93      REAL pkf(ip1jmp1,llm)                  ! filtered exner at mid-layer
94      REAL phi(ip1jmp1,llm)                  ! geopotential
95      REAL w(ip1jmp1,llm)                    ! vertical velocity
96! ADAPTATION GCM POUR CP(T)
97      REAL temp(ip1jmp1,llm)                 ! temperature 
98      REAL tsurpk(ip1jmp1,llm)               ! cpp*T/pk 
99
100!      real zqmin,zqmax
101
102c variables dynamiques intermediaire pour le transport
103      REAL pbaru(ip1jmp1,llm),pbarv(ip1jm,llm) !flux de masse
104
105c   variables dynamiques au pas -1
106      REAL vcovm1(ip1jm,llm),ucovm1(ip1jmp1,llm)
107      REAL tetam1(ip1jmp1,llm),psm1(ip1jmp1)
108      REAL massem1(ip1jmp1,llm)
109
110c   tendances dynamiques en */s
111      REAL dv(ip1jm,llm),du(ip1jmp1,llm)
112      REAL dteta(ip1jmp1,llm),dq(ip1jmp1,llm,nqtot),dp(ip1jmp1)
113
114c   tendances de la dissipation en */s
115      REAL dvdis(ip1jm,llm),dudis(ip1jmp1,llm)
116      REAL dtetadis(ip1jmp1,llm)
117
118c   tendances de la couche superieure */s
119c      REAL dvtop(ip1jm,llm)
120      REAL dutop(ip1jmp1,llm)
121c      REAL dtetatop(ip1jmp1,llm)
122c      REAL dqtop(ip1jmp1,llm,nqtot),dptop(ip1jmp1)
123
124c   TITAN : tendances due au forces de marees */s
125      REAL dvtidal(ip1jm,llm),dutidal(ip1jmp1,llm)
126
127c   tendances physiques */s
128      REAL dvfi(ip1jm,llm),dufi(ip1jmp1,llm)
129      REAL dtetafi(ip1jmp1,llm),dqfi(ip1jmp1,llm,nqtot),dpfi(ip1jmp1)
130
131c   variables pour le fichier histoire
132!      REAL dtav      ! intervalle de temps elementaire
133
134      REAL tppn(iim),tpps(iim),tpn,tps
135c
136      INTEGER itau,itaufinp1,iav
137!      INTEGER  iday ! jour julien
138      REAL       time
139
140      REAL  SSUM
141!     REAL finvmaold(ip1jmp1,llm)
142
143cym      LOGICAL  lafin
144      LOGICAL :: lafin=.false.
145      INTEGER ij,iq,l
146!      INTEGER ik
147
148!      real time_step, t_wrt, t_ops
149
150      REAL rdaym_ini
151! jD_cur: jour julien courant
152! jH_cur: heure julienne courante
153      REAL :: jD_cur, jH_cur
154!      INTEGER :: an, mois, jour
155!      REAL :: secondes
156
157      LOGICAL first,callinigrads
158cIM : pour sortir les param. du modele dans un fis. netcdf 110106
159      save first
160      data first/.true./
161!      real dt_cum
162!      character*10 infile
163!      integer zan, tau0, thoriid
164!      integer nid_ctesGCM
165!      save nid_ctesGCM
166!      real degres
167!      real rlong(iip1), rlatg(jjp1)
168!      real zx_tmp_2d(iip1,jjp1)
169!      integer ndex2d(iip1*jjp1)
170      logical ok_sync
171      parameter (ok_sync = .true.)
172      logical physics
173
174      data callinigrads/.true./
175      character*10 string10
176
177!      REAL alpha(ip1jmp1,llm),beta(ip1jmp1,llm)
178      REAL :: flxw(ip1jmp1,llm)  ! flux de masse verticale
179
180c+jld variables test conservation energie
181      REAL ecin(ip1jmp1,llm),ecin0(ip1jmp1,llm)
182C     Tendance de la temp. potentiel d (theta)/ d t due a la
183C     tansformation d'energie cinetique en energie thermique
184C     cree par la dissipation
185      REAL dtetaecdt(ip1jmp1,llm)
186      REAL vcont(ip1jm,llm),ucont(ip1jmp1,llm)
187      REAL vnat(ip1jm,llm),unat(ip1jmp1,llm)
188!      REAL      d_h_vcol, d_qt, d_qw, d_ql, d_ec
189      CHARACTER*15 ztit
190!IM   INTEGER   ip_ebil_dyn  ! PRINT level for energy conserv. diag.
191!IM   SAVE      ip_ebil_dyn
192!IM   DATA      ip_ebil_dyn/0/
193c-jld
194
195!      integer :: itau_w ! for write_paramLMDZ_dyn.h
196
197!      character*80 dynhist_file, dynhistave_file
198      character(len=*),parameter :: modname="leapfrog"
199      character*80 abort_message
200
201      logical dissip_conservative
202      save dissip_conservative
203      data dissip_conservative/.true./
204
205      INTEGER testita
206      PARAMETER (testita = 9)
207
208      logical , parameter :: flag_verif = .false.
209     
210      ! for CP(T)
211      real :: dtec
212      real :: ztetaec(ip1jmp1,llm)
213
214      if (nday>=0) then
215         itaufin   = nday*day_step
216      else
217         ! to run a given (-nday) number of dynamical steps
218         itaufin   = -nday
219      endif
220      if (less1day) then
221c MODIF VENUS: to run less than one day:
222        itaufin   = int(fractday*day_step)
223      endif
224      if (ndynstep.gt.0) then
225        ! running a given number of dynamical steps
226        itaufin=ndynstep
227      endif
228      itaufinp1 = itaufin +1
229     
230c INITIALISATIONS
231        dudis(:,:)   =0.
232        dvdis(:,:)   =0.
233        dtetadis(:,:)=0.
234        dutop(:,:)   =0.
235c        dvtop(:,:)   =0.
236c        dtetatop(:,:)=0.
237c        dqtop(:,:,:) =0.
238c        dptop(:)     =0.
239        dufi(:,:)   =0.
240        dvfi(:,:)   =0.
241        dtetafi(:,:)=0.
242        dqfi(:,:,:) =0.
243        dpfi(:)     =0.
244
245      itau = 0
246      physics=.true.
247      if (iflag_phys==0.or.iflag_phys==2) physics=.false.
248
249c      iday = day_ini+itau/day_step
250c      time = REAL(itau-(iday-day_ini)*day_step)/day_step+time_0
251c         IF(time.GT.1.) THEN
252c          time = time-1.
253c          iday = iday+1
254c         ENDIF
255
256
257c-----------------------------------------------------------------------
258c   On initialise la pression et la fonction d'Exner :
259c   --------------------------------------------------
260
261      dq(:,:,:)=0.
262      CALL pression ( ip1jmp1, ap, bp, ps, p       )
263      if (pressure_exner) then
264        CALL exner_hyb( ip1jmp1, ps, p, pks, pk, pkf )
265      else
266        CALL exner_milieu( ip1jmp1, ps, p, pks, pk, pkf )
267      endif
268
269c------------------
270c TEST PK MONOTONE
271c------------------
272      write(*,*) "Test PK"
273      do ij=1,ip1jmp1
274        do l=2,llm
275          if(pk(ij,l).gt.pk(ij,l-1)) then
276c           write(*,*) ij,l,pk(ij,l)
277            abort_message = 'PK non strictement decroissante'
278            call abort_gcm(modname,abort_message,1)
279c           write(*,*) "ATTENTION, Test PK deconnecté..."
280          endif
281        enddo
282      enddo
283      write(*,*) "Fin Test PK"
284c     stop
285c------------------
286
287c-----------------------------------------------------------------------
288c   Debut de l'integration temporelle:
289c   ----------------------------------
290
291   1  CONTINUE ! Matsuno Forward step begins here
292
293c   date: (NB: date remains unchanged for Backward step)
294c   -----
295
296      jD_cur = jD_ref + day_ini - day_ref +                             &
297     &          (itau+1)/day_step
298      jH_cur = jH_ref + start_time +                                    &
299     &          mod(itau+1,day_step)/float(day_step)
300      jD_cur = jD_cur + int(jH_cur)
301      jH_cur = jH_cur - int(jH_cur)
302
303        if (ok_iso_verif) then
304           call check_isotopes_seq(q,ip1jmp1,'leapfrog 321')
305        endif !if (ok_iso_verif) then
306
307#ifdef CPP_IOIPSL
308      IF (planet_type.eq."earth") THEN
309      if (ok_guide) then
310        call guide_main(itau,ucov,vcov,teta,q,masse,ps)
311      endif
312      ENDIF
313#endif
314
315
316c
317c     IF( MOD( itau, 10* day_step ).EQ.0 )  THEN
318c       CALL  test_period ( ucov,vcov,teta,q,p,phis )
319c       PRINT *,' ----   Test_period apres continue   OK ! -----', itau
320c     ENDIF
321c
322
323! Save fields obtained at previous time step as '...m1'
324      CALL SCOPY( ijmllm ,vcov , 1, vcovm1 , 1 )
325      CALL SCOPY( ijp1llm,ucov , 1, ucovm1 , 1 )
326      CALL SCOPY( ijp1llm,teta , 1, tetam1 , 1 )
327      CALL SCOPY( ijp1llm,masse, 1, massem1, 1 )
328      CALL SCOPY( ip1jmp1, ps  , 1,   psm1 , 1 )
329
330      forward = .TRUE.
331      leapf   = .FALSE.
332      dt      =  dtvr
333
334c   ...    P.Le Van .26/04/94  ....
335! Ehouarn: finvmaold is actually not used
336!      CALL SCOPY   ( ijp1llm,   masse, 1, finvmaold,     1 )
337!      CALL filtreg ( finvmaold ,jjp1, llm, -2,2, .TRUE., 1 )
338
339        if (ok_iso_verif) then
340           call check_isotopes_seq(q,ip1jmp1,'leapfrog 400')
341        endif !if (ok_iso_verif) then
342
343   2  CONTINUE ! Matsuno backward or leapfrog step begins here
344
345c-----------------------------------------------------------------------
346
347c   date: (NB: only leapfrog step requires recomputing date)
348c   -----
349
350      IF (leapf) THEN
351        jD_cur = jD_ref + day_ini - day_ref +
352     &            (itau+1)/day_step
353        jH_cur = jH_ref + start_time +
354     &            mod(itau+1,day_step)/float(day_step)
355        jD_cur = jD_cur + int(jH_cur)
356        jH_cur = jH_cur - int(jH_cur)
357      ENDIF
358
359
360c   gestion des appels de la physique et des dissipations:
361c   ------------------------------------------------------
362c
363c   ...    P.Le Van  ( 6/02/95 )  ....
364
365      apphys = .FALSE.
366      statcl = .FALSE.
367      conser = .FALSE.
368      apdiss = .FALSE.
369
370      IF( purmats ) THEN
371      ! Purely Matsuno time stepping
372         IF( MOD(itau,iconser) .EQ.0.AND.  forward    ) conser = .TRUE.
373         IF( MOD(itau,dissip_period ).EQ.0.AND..NOT.forward )
374     s        apdiss = .TRUE.
375         IF( MOD(itau,iphysiq ).EQ.0.AND..NOT.forward
376     s          .and. physics                        ) apphys = .TRUE.
377      ELSE
378      ! Leapfrog/Matsuno time stepping
379         IF( MOD(itau   ,iconser) .EQ. 0              ) conser = .TRUE.
380         IF( MOD(itau+1,dissip_period).EQ.0 .AND. .NOT. forward )
381     s        apdiss = .TRUE.
382         IF( MOD(itau+1,iphysiq).EQ.0.AND.physics       ) apphys=.TRUE.
383      END IF
384
385! Ehouarn: for Shallow Water case (ie: 1 vertical layer),
386!          supress dissipation step
387      if (llm.eq.1) then
388        apdiss=.false.
389      endif
390
391#ifdef NODYN
392      apdiss=.false.
393#endif
394
395        if (ok_iso_verif) then
396           call check_isotopes_seq(q,ip1jmp1,'leapfrog 589')
397        endif !if (ok_iso_verif) then
398
399c-----------------------------------------------------------------------
400c   calcul des tendances dynamiques:
401c   --------------------------------
402
403! ADAPTATION GCM POUR CP(T)
404      call tpot2t(ijp1llm,teta,temp,pk)
405      tsurpk = cpp*temp/pk
406      ! compute geopotential phi()
407      CALL geopot  ( ip1jmp1, tsurpk  , pk , pks,  phis  , phi   )
408
409           rdaym_ini  = itau * dtvr / daysec
410
411      time = jD_cur + jH_cur
412
413#ifdef NODYN
414      dv(:,:) = 0.D+0
415      du(:,:) = 0.D+0
416      dteta(:,:) = 0.D+0
417      dq(:,:,:) = 0.D+0
418      dp(:) = 0.D+0
419#else
420      CALL caldyn
421     $  ( itau,ucov,vcov,teta,ps,masse,pk,pkf,tsurpk,phis ,
422     $    phi,conser,du,dv,dteta,dp,w, pbaru,pbarv, time )
423
424      ! Simple zonal wind nudging for generic planetary model
425      ! AS 09/2013
426      ! ---------------------------------------------------
427      if (planet_type.eq."generic") then
428       if (ok_guide) then
429         du(:,:) = du(:,:) + ((uforc(:,:)-ucov(:,:)) / facwind)
430       endif
431      endif
432
433c-----------------------------------------------------------------------
434c   calcul des tendances advection des traceurs (dont l'humidite)
435c   -------------------------------------------------------------
436
437        if (ok_iso_verif) then
438           call check_isotopes_seq(q,ip1jmp1,
439     &           'leapfrog 686: avant caladvtrac')
440        endif !if (ok_iso_verif) then
441
442      IF( forward. OR . leapf )  THEN
443! Ehouarn: NB: fields sent to advtrac are those at the beginning of the time step
444         CALL caladvtrac(q,pbaru,pbarv,
445     *        p, masse, dq,  teta,
446     .        flxw, pk)
447         
448         IF (offline) THEN
449Cmaf stokage du flux de masse pour traceurs OFF-LINE
450
451#ifdef CPP_IOIPSL
452           CALL fluxstokenc(pbaru,pbarv,masse,teta,phi,phis,
453     .   dtvr, itau)
454#endif
455
456
457         ENDIF ! of IF (offline)
458c
459      ENDIF ! of IF( forward. OR . leapf )
460
461
462c-----------------------------------------------------------------------
463c   integrations dynamique et traceurs:
464c   ----------------------------------
465
466        if (ok_iso_verif) then
467           write(*,*) 'leapfrog 720'
468           call check_isotopes_seq(q,ip1jmp1,'leapfrog 756')
469        endif !if (ok_iso_verif) then
470       
471       CALL integrd ( nqtot,vcovm1,ucovm1,tetam1,psm1,massem1 ,
472     $         dv,du,dteta,dq,dp,vcov,ucov,teta,q,ps,masse,phis )
473!     $              finvmaold                                    )
474
475       if (ok_iso_verif) then
476          write(*,*) 'leapfrog 724'
477           call check_isotopes_seq(q,ip1jmp1,'leapfrog 762')
478        endif !if (ok_iso_verif) then
479
480       IF ((planet_type.eq."titan").and.(tidal)) then
481c-----------------------------------------------------------------------
482c   Marées gravitationnelles causées par Saturne
483c   B. Charnay (28/10/2010)
484c   ----------------------------------------------------------
485            CALL tidal_forces(rdaym_ini, dutidal, dvtidal)
486            ucov=ucov+dutidal*dt
487            vcov=vcov+dvtidal*dt
488       ENDIF
489
490! NODYN precompiling flag
491#endif
492
493c .P.Le Van (26/04/94  ajout de  finvpold dans l'appel d'integrd)
494c
495c-----------------------------------------------------------------------
496c   calcul des tendances physiques:
497c   -------------------------------
498c    ########   P.Le Van ( Modif le  6/02/95 )   ###########
499c
500       IF( purmats )  THEN
501          IF( itau.EQ.itaufin.AND..NOT.forward ) lafin = .TRUE.
502       ELSE
503          IF( itau+1. EQ. itaufin )              lafin = .TRUE.
504       ENDIF
505c
506c
507       IF( apphys )  THEN
508c
509c     .......   Ajout   P.Le Van ( 17/04/96 )   ...........
510c
511
512         CALL pression (  ip1jmp1, ap, bp, ps,  p      )
513         if (pressure_exner) then
514           CALL exner_hyb(  ip1jmp1, ps, p,pks, pk, pkf )
515         else
516           CALL exner_milieu( ip1jmp1, ps, p, pks, pk, pkf )
517         endif
518
519! Compute geopotential (physics might need it)
520         CALL geopot  ( ip1jmp1, teta  , pk , pks,  phis  , phi   )
521
522           jD_cur = jD_ref + day_ini - day_ref +                        &
523     &          (itau+1)/day_step
524
525           IF ((planet_type .eq."generic").or.
526     &         (planet_type .eq."mars")) THEN
527              ! AS: we make jD_cur to be pday
528              jD_cur = int(day_ini + itau/day_step)
529           ENDIF
530
531           jH_cur = jH_ref + start_time +                               &
532     &          mod(itau+1,day_step)/float(day_step)
533           IF ((planet_type .eq."generic").or.
534     &         (planet_type .eq."mars")) THEN
535             jH_cur = jH_ref + start_time +                               &
536     &          mod(itau,day_step)/float(day_step)
537           ENDIF
538           jD_cur = jD_cur + int(jH_cur)
539           jH_cur = jH_cur - int(jH_cur)
540!         write(lunout,*)'itau, jD_cur = ', itau, jD_cur, jH_cur
541!         call ju2ymds(jD_cur+jH_cur, an, mois, jour, secondes)
542!         write(lunout,*)'current date = ',an, mois, jour, secondes
543
544c rajout debug
545c       lafin = .true.
546
547
548c   Interface avec les routines de phylmd (phymars ... )
549c   -----------------------------------------------------
550
551c+jld
552
553c  Diagnostique de conservation de l'énergie : initialisation
554         IF (ip_ebil_dyn.ge.1 ) THEN
555          ztit='bil dyn'
556! Ehouarn: be careful, diagedyn is Earth-specific!
557           IF (planet_type.eq."earth") THEN
558            CALL diagedyn(ztit,2,1,1,dtphys
559     &    , ucov    , vcov , ps, p ,pk , teta , q(:,:,1), q(:,:,2))
560           ENDIF
561         ENDIF ! of IF (ip_ebil_dyn.ge.1 )
562c-jld
563#ifdef CPP_IOIPSL
564cIM decommenter les 6 lignes suivantes pour sortir quelques parametres dynamiques de LMDZ
565cIM uncomment next 6 lines to get some parameters for LMDZ dynamics
566c        IF (first) THEN
567c         first=.false.
568c#include "ini_paramLMDZ_dyn.h"
569c        ENDIF
570c
571c#include "write_paramLMDZ_dyn.h"
572c
573#endif
574! #endif of #ifdef CPP_IOIPSL
575
576c          call WriteField('pfi',reshape(p,(/iip1,jmp1,llmp1/)))
577
578         CALL calfis( lafin , jD_cur, jH_cur,
579     $               ucov,vcov,teta,q,masse,ps,p,pk,phis,phi ,
580     $               du,dv,dteta,dq,
581     $               flxw,
582     $               dufi,dvfi,dtetafi,dqfi,dpfi  )
583
584c          call WriteField('dufi',reshape(dufi,(/iip1,jmp1,llm/)))
585c          call WriteField('dvfi',reshape(dvfi,(/iip1,jjm,llm/)))
586c          call WriteField('dtetafi',reshape(dtetafi,(/iip1,jmp1,llm/)))
587
588c      ajout des tendances physiques:
589c      ------------------------------
590          CALL addfi( dtphys, leapf, forward   ,
591     $                  ucov, vcov, teta , q   ,ps ,
592     $                 dufi, dvfi, dtetafi , dqfi ,dpfi  )
593          ! since addfi updates ps(), also update p(), masse() and pk()
594          CALL pression (ip1jmp1,ap,bp,ps,p)
595          CALL massdair(p,masse)
596          if (pressure_exner) then
597            CALL exner_hyb( ip1jmp1, ps, p, pks, pk, pkf )
598          else
599            CALL exner_milieu( ip1jmp1, ps, p, pks, pk, pkf )
600          endif
601         
602c      Couche superieure :
603c      -------------------
604         IF (iflag_top_bound > 0) THEN
605           CALL top_bound(vcov,ucov,teta,masse,dtphys,dutop)
606           dutop(:,:)=dutop(:,:)/dtphys   ! convert to a tendency in (m/s)/s
607         ENDIF
608
609c  Diagnostique de conservation de l'énergie : difference
610         IF (ip_ebil_dyn.ge.1 ) THEN
611          ztit='bil phys'
612          IF (planet_type.eq."earth") THEN
613           CALL diagedyn(ztit,2,1,1,dtphys
614     &     , ucov    , vcov , ps, p ,pk , teta , q(:,:,1), q(:,:,2))
615          ENDIF
616         ENDIF ! of IF (ip_ebil_dyn.ge.1 )
617
618       ENDIF ! of IF( apphys )
619
620      IF(iflag_phys.EQ.2) THEN ! "Newtonian" case
621!   Academic case : Simple friction and Newtonan relaxation
622!   -------------------------------------------------------
623        DO l=1,llm   
624          DO ij=1,ip1jmp1
625           teta(ij,l)=teta(ij,l)-dtvr*
626     &      (teta(ij,l)-tetarappel(ij,l))*(knewt_g+knewt_t(l)*clat4(ij))
627          ENDDO
628        ENDDO ! of DO l=1,llm
629   
630        if (planet_type.eq."giant") then
631          ! add an intrinsic heat flux at the base of the atmosphere
632          teta(:,1)=teta(:,1)+dtvr*aire(:)*ihf/cpp/masse(:,1)
633        endif
634
635        call friction(ucov,vcov,dtvr)
636
637   
638        ! Sponge layer (if any)
639        IF (ok_strato) THEN
640           CALL top_bound(vcov,ucov,teta,masse,dtvr,dutop)
641           dutop(:,:)=dutop(:,:)/dtvr   ! convert to a tendency in (m/s)/s
642        ENDIF ! of IF (ok_strato)
643      ENDIF ! of IF (iflag_phys.EQ.2)
644
645
646c-jld
647
648        CALL pression ( ip1jmp1, ap, bp, ps, p                  )
649        if (pressure_exner) then
650          CALL exner_hyb( ip1jmp1, ps, p, pks, pk, pkf )
651        else
652          CALL exner_milieu( ip1jmp1, ps, p, pks, pk, pkf )
653        endif
654        CALL massdair(p,masse)
655
656        if (ok_iso_verif) then
657           call check_isotopes_seq(q,ip1jmp1,'leapfrog 1196')
658        endif !if (ok_iso_verif) then
659
660c-----------------------------------------------------------------------
661c   dissipation horizontale et verticale  des petites echelles:
662c   ----------------------------------------------------------
663
664      IF(apdiss) THEN
665
666        ! sponge layer
667        if (callsponge) then
668          CALL sponge(ucov,vcov,teta,ps,dtdiss,mode_sponge)
669        endif
670
671c   calcul de l'energie cinetique avant dissipation
672        call covcont(llm,ucov,vcov,ucont,vcont)
673        call enercin(vcov,ucov,vcont,ucont,ecin0)
674
675c   dissipation
676! ADAPTATION GCM POUR CP(T)
677        call tpot2t(ijp1llm,teta,temp,pk)
678
679        CALL dissip(vcov,ucov,teta,p,dvdis,dudis,dtetadis)
680        ucov=ucov+dudis
681        vcov=vcov+dvdis
682        dudis=dudis/dtdiss   ! passage en (m/s)/s
683        dvdis=dvdis/dtdiss   ! passage en (m/s)/s
684
685c------------------------------------------------------------------------
686        if (dissip_conservative) then
687C       On rajoute la tendance due a la transform. Ec -> E therm. cree
688C       lors de la dissipation
689            call covcont(llm,ucov,vcov,ucont,vcont)
690            call enercin(vcov,ucov,vcont,ucont,ecin)
691! ADAPTATION GCM POUR CP(T)
692            do ij=1,ip1jmp1
693              do l=1,llm
694                dtec = (ecin0(ij,l)-ecin(ij,l))/cpdet(temp(ij,l))
695                temp(ij,l) = temp(ij,l) + dtec
696              enddo
697            enddo
698            call t2tpot(ijp1llm,temp,ztetaec,pk)
699            dtetaecdt=ztetaec-teta
700            dtetadis=dtetadis+dtetaecdt
701        endif
702        teta=teta+dtetadis
703        dtetadis=dtetadis/dtdiss   ! passage en K/s
704c------------------------------------------------------------------------
705
706
707c    .......        P. Le Van (  ajout  le 17/04/96  )   ...........
708c   ...      Calcul de la valeur moyenne, unique de h aux poles  .....
709c
710
711        DO l  =  1, llm
712          DO ij =  1,iim
713           tppn(ij)  = aire(  ij    ) * teta(  ij    ,l)
714           tpps(ij)  = aire(ij+ip1jm) * teta(ij+ip1jm,l)
715          ENDDO
716           tpn  = SSUM(iim,tppn,1)/apoln
717           tps  = SSUM(iim,tpps,1)/apols
718
719          DO ij = 1, iip1
720           teta(  ij    ,l) = tpn
721           teta(ij+ip1jm,l) = tps
722          ENDDO
723        ENDDO
724
725        if (1 == 0) then
726!!! Ehouarn: lines here 1) kill 1+1=2 in the dynamics
727!!!                     2) should probably not be here anyway
728!!! but are kept for those who would want to revert to previous behaviour
729           DO ij =  1,iim
730             tppn(ij)  = aire(  ij    ) * ps (  ij    )
731             tpps(ij)  = aire(ij+ip1jm) * ps (ij+ip1jm)
732           ENDDO
733             tpn  = SSUM(iim,tppn,1)/apoln
734             tps  = SSUM(iim,tpps,1)/apols
735
736           DO ij = 1, iip1
737             ps(  ij    ) = tpn
738             ps(ij+ip1jm) = tps
739           ENDDO
740        endif ! of if (1 == 0)
741
742      END IF ! of IF(apdiss)
743
744c ajout debug
745c              IF( lafin ) then 
746c                abort_message = 'Simulation finished'
747c                call abort_gcm(modname,abort_message,0)
748c              ENDIF
749       
750c   ********************************************************************
751c   ********************************************************************
752c   .... fin de l'integration dynamique  et physique pour le pas itau ..
753c   ********************************************************************
754c   ********************************************************************
755
756c   preparation du pas d'integration suivant  ......
757
758        if (ok_iso_verif) then
759           call check_isotopes_seq(q,ip1jmp1,'leapfrog 1509')
760        endif !if (ok_iso_verif) then
761
762      IF ( .NOT.purmats ) THEN
763c       ........................................................
764c       ..............  schema matsuno + leapfrog  ..............
765c       ........................................................
766
767            IF(forward. OR. leapf) THEN
768              itau= itau + 1
769c              iday= day_ini+itau/day_step
770c              time= REAL(itau-(iday-day_ini)*day_step)/day_step+time_0
771c                IF(time.GT.1.) THEN
772c                  time = time-1.
773c                  iday = iday+1
774c                ENDIF
775            ENDIF
776
777
778            IF( itau. EQ. itaufinp1 ) then 
779              if (flag_verif) then
780                write(79,*) 'ucov',ucov
781                write(80,*) 'vcov',vcov
782                write(81,*) 'teta',teta
783                write(82,*) 'ps',ps
784                write(83,*) 'q',q
785                WRITE(85,*) 'q1 = ',q(:,:,1)
786                WRITE(86,*) 'q3 = ',q(:,:,3)
787              endif
788
789              abort_message = 'Simulation finished'
790
791              call abort_gcm(modname,abort_message,0)
792            ENDIF
793c-----------------------------------------------------------------------
794c   ecriture du fichier histoire moyenne:
795c   -------------------------------------
796
797            IF(MOD(itau,iperiod).EQ.0 .OR. itau.EQ.itaufin) THEN
798               IF(itau.EQ.itaufin) THEN
799                  iav=1
800               ELSE
801                  iav=0
802               ENDIF
803               
804!              ! Ehouarn: re-compute geopotential for outputs
805               CALL geopot(ip1jmp1,teta,pk,pks,phis,phi)
806
807               IF (ok_dynzon) THEN
808#ifdef CPP_IOIPSL
809c les traceurs ne sont pas sortis, trop lourd.
810c Peut changer eventuellement si besoin.
811                 CALL bilan_dyn(dtvr*iperiod,dtvr*day_step*periodav,
812     &                 ps,masse,pk,pbaru,pbarv,teta,phi,ucov,vcov,
813     &                 du,dudis,dutop,dufi)
814#endif
815               END IF
816               IF (ok_dyn_ave) THEN
817#ifdef CPP_IOIPSL
818                 CALL writedynav(itau,vcov,
819     &                 ucov,teta,pk,phi,q,masse,ps,phis)
820#endif
821               ENDIF
822
823            ENDIF ! of IF((MOD(itau,iperiod).EQ.0).OR.(itau.EQ.itaufin))
824
825        if (ok_iso_verif) then
826           call check_isotopes_seq(q,ip1jmp1,'leapfrog 1584')
827        endif !if (ok_iso_verif) then
828
829c-----------------------------------------------------------------------
830c   ecriture de la bande histoire:
831c   ------------------------------
832
833            IF( MOD(itau,iecri).EQ.0) THEN
834             ! Ehouarn: output only during LF or Backward Matsuno
835             if (leapf.or.(.not.leapf.and.(.not.forward))) then
836! ADAPTATION GCM POUR CP(T)
837              call tpot2t(ijp1llm,teta,temp,pk)
838              tsurpk = cpp*temp/pk
839              CALL geopot(ip1jmp1,tsurpk,pk,pks,phis,phi)
840              unat=0.
841              do l=1,llm
842                unat(iip2:ip1jm,l)=ucov(iip2:ip1jm,l)/cu(iip2:ip1jm)
843                vnat(:,l)=vcov(:,l)/cv(:)
844              enddo
845#ifdef CPP_IOIPSL
846              if (ok_dyn_ins) then
847!               write(lunout,*) "leapfrog: call writehist, itau=",itau
848               CALL writehist(itau,vcov,ucov,teta,phi,q,masse,ps,phis)
849!               call WriteField('ucov',reshape(ucov,(/iip1,jmp1,llm/)))
850!               call WriteField('vcov',reshape(vcov,(/iip1,jjm,llm/)))
851!              call WriteField('teta',reshape(teta,(/iip1,jmp1,llm/)))
852!               call WriteField('ps',reshape(ps,(/iip1,jmp1/)))
853!               call WriteField('masse',reshape(masse,(/iip1,jmp1,llm/)))
854              endif ! of if (ok_dyn_ins)
855#endif
856! For some Grads outputs of fields
857              if (output_grads_dyn) then
858#include "write_grads_dyn.h"
859              endif
860             endif ! of if (leapf.or.(.not.leapf.and.(.not.forward)))
861            ENDIF ! of IF(MOD(itau,iecri).EQ.0)
862
863            IF(itau.EQ.itaufin) THEN
864
865              if (planet_type=="mars") then
866                CALL dynredem1("restart.nc",REAL(itau)/REAL(day_step),
867     &                         vcov,ucov,teta,q,masse,ps)
868              else
869                CALL dynredem1("restart.nc",start_time,
870     &                         vcov,ucov,teta,q,masse,ps)
871              endif
872              CLOSE(99)
873              !!! Ehouarn: Why not stop here and now?
874            ENDIF ! of IF (itau.EQ.itaufin)
875
876c-----------------------------------------------------------------------
877c   gestion de l'integration temporelle:
878c   ------------------------------------
879
880            IF( MOD(itau,iperiod).EQ.0 )    THEN
881                    GO TO 1
882            ELSE IF ( MOD(itau-1,iperiod). EQ. 0 ) THEN
883
884                   IF( forward )  THEN
885c      fin du pas forward et debut du pas backward
886
887                      forward = .FALSE.
888                        leapf = .FALSE.
889                           GO TO 2
890
891                   ELSE
892c      fin du pas backward et debut du premier pas leapfrog
893
894                        leapf =  .TRUE.
895                        dt  =  2.*dtvr
896                        GO TO 2
897                   END IF ! of IF (forward)
898            ELSE
899
900c      ......   pas leapfrog  .....
901
902                 leapf = .TRUE.
903                 dt  = 2.*dtvr
904                 GO TO 2
905            END IF ! of IF (MOD(itau,iperiod).EQ.0)
906                   !    ELSEIF (MOD(itau-1,iperiod).EQ.0)
907
908      ELSE ! of IF (.not.purmats)
909
910        if (ok_iso_verif) then
911           call check_isotopes_seq(q,ip1jmp1,'leapfrog 1664')
912        endif !if (ok_iso_verif) then
913
914c       ........................................................
915c       ..............       schema  matsuno        ...............
916c       ........................................................
917            IF( forward )  THEN
918
919             itau =  itau + 1
920c             iday = day_ini+itau/day_step
921c             time = REAL(itau-(iday-day_ini)*day_step)/day_step+time_0
922c
923c                  IF(time.GT.1.) THEN
924c                   time = time-1.
925c                   iday = iday+1
926c                  ENDIF
927
928               forward =  .FALSE.
929               IF( itau. EQ. itaufinp1 ) then 
930                 abort_message = 'Simulation finished'
931                 call abort_gcm(modname,abort_message,0)
932               ENDIF
933               GO TO 2
934
935            ELSE ! of IF(forward) i.e. backward step
936
937        if (ok_iso_verif) then
938           call check_isotopes_seq(q,ip1jmp1,'leapfrog 1698')
939        endif !if (ok_iso_verif) then 
940
941              IF(MOD(itau,iperiod).EQ.0 .OR. itau.EQ.itaufin) THEN
942               IF(itau.EQ.itaufin) THEN
943                  iav=1
944               ELSE
945                  iav=0
946               ENDIF
947
948!              ! Ehouarn: re-compute geopotential for outputs
949               CALL geopot(ip1jmp1,teta,pk,pks,phis,phi)
950
951               IF (ok_dynzon) THEN
952#ifdef CPP_IOIPSL
953c les traceurs ne sont pas sortis, trop lourd.
954c Peut changer eventuellement si besoin.
955                 CALL bilan_dyn(dtvr*iperiod,dtvr*day_step*periodav,
956     &                 ps,masse,pk,pbaru,pbarv,teta,phi,ucov,vcov,
957     &                 du,dudis,dutop,dufi)
958#endif
959               ENDIF
960               IF (ok_dyn_ave) THEN
961#ifdef CPP_IOIPSL
962                 CALL writedynav(itau,vcov,
963     &                 ucov,teta,pk,phi,q,masse,ps,phis)
964#endif
965               ENDIF
966
967              ENDIF ! of IF(MOD(itau,iperiod).EQ.0 .OR. itau.EQ.itaufin)
968
969              IF(MOD(itau,iecri         ).EQ.0) THEN
970c              IF(MOD(itau,iecri*day_step).EQ.0) THEN
971! ADAPTATION GCM POUR CP(T)
972                call tpot2t(ijp1llm,teta,temp,pk)
973                tsurpk = cpp*temp/pk
974                CALL geopot(ip1jmp1,tsurpk,pk,pks,phis,phi)
975                unat=0.
976                do l=1,llm
977                  unat(iip2:ip1jm,l)=ucov(iip2:ip1jm,l)/cu(iip2:ip1jm)
978                  vnat(:,l)=vcov(:,l)/cv(:)
979                enddo
980#ifdef CPP_IOIPSL
981              if (ok_dyn_ins) then
982!                write(lunout,*) "leapfrog: call writehist (b)",
983!     &                        itau,iecri
984                CALL writehist(itau,vcov,ucov,teta,phi,q,masse,ps,phis)
985              endif ! of if (ok_dyn_ins)
986#endif
987! For some Grads outputs
988                if (output_grads_dyn) then
989#include "write_grads_dyn.h"
990                endif
991
992              ENDIF ! of IF(MOD(itau,iecri         ).EQ.0)
993
994              IF(itau.EQ.itaufin) THEN
995                if (planet_type=="mars") then
996                  CALL dynredem1("restart.nc",REAL(itau)/REAL(day_step),
997     &                         vcov,ucov,teta,q,masse,ps)
998                else
999                  CALL dynredem1("restart.nc",start_time,
1000     &                         vcov,ucov,teta,q,masse,ps)
1001                endif
1002              ENDIF ! of IF(itau.EQ.itaufin)
1003
1004              forward = .TRUE.
1005              GO TO  1
1006
1007            ENDIF ! of IF (forward)
1008
1009      END IF ! of IF(.not.purmats)
1010
1011      STOP
1012      END
Note: See TracBrowser for help on using the repository browser.