/[lmdze]/trunk/libf/dyn3d/leapfrog.f90
ViewVC logotype

Diff of /trunk/libf/dyn3d/leapfrog.f90

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

revision 22 by guez, Fri Jul 31 15:18:47 2009 UTC revision 55 by guez, Mon Dec 12 13:25:01 2011 UTC
# Line 1  Line 1 
1  module leapfrog_m  module leapfrog_m
2    
   ! This module is clean: no C preprocessor directive, no include line.  
   
3    IMPLICIT NONE    IMPLICIT NONE
4    
5  contains  contains
6    
7    SUBROUTINE leapfrog(ucov, vcov, teta, ps, masse, phis, nq, q, time_0)    SUBROUTINE leapfrog(ucov, vcov, teta, ps, masse, phis, q, time_0)
8    
9      ! From dyn3d/leapfrog.F, version 1.6 2005/04/13 08:58:34      ! From dyn3d/leapfrog.F, version 1.6, 2005/04/13 08:58:34
10        ! Authors: P. Le Van, L. Fairhead, F. Hourdin
11        ! Matsuno-leapfrog scheme.
12    
13        use addfi_m, only: addfi
14        use bilan_dyn_m, only: bilan_dyn
15        use caladvtrac_m, only: caladvtrac
16        use caldyn_m, only: caldyn
17        USE calfis_m, ONLY: calfis
18        USE com_io_dyn, ONLY: histaveid
19        USE comconst, ONLY: daysec, dtphys, dtvr
20        USE comgeom, ONLY: aire_2d, apoln, apols
21        USE comvert, ONLY: ap, bp
22        USE conf_gcm_m, ONLY: day_step, iconser, iperiod, iphysiq, nday, offline, &
23             periodav
24        USE dimens_m, ONLY: iim, jjm, llm, nqmx
25        use dissip_m, only: dissip
26        USE dynetat0_m, ONLY: day_ini
27        use dynredem1_m, only: dynredem1
28        USE exner_hyb_m, ONLY: exner_hyb
29        use filtreg_m, only: filtreg
30        use geopot_m, only: geopot
31        USE guide_m, ONLY: guide
32        use inidissip_m, only: idissip
33        use integrd_m, only: integrd
34        USE logic, ONLY: iflag_phys, ok_guide
35        use nr_util, only: assert
36        USE pressure_var, ONLY: p3d
37        USE temps, ONLY: itau_dyn
38    
39      ! Version du 10/01/98, avec coordonnees verticales hybrides, avec      ! Variables dynamiques:
40      ! nouveaux operat. dissipation * (gradiv2, divgrad2, nxgraro2)      REAL, intent(inout):: ucov(:, :, :) ! (iim + 1, jjm + 1, llm) vent covariant
41        REAL, intent(inout):: vcov(:, :, :) ! (iim + 1, jjm, llm) ! vent covariant
42    
43      ! Auteur: P. Le Van /L. Fairhead/F.Hourdin      REAL, intent(inout):: teta(:, :, :) ! (iim + 1, jjm + 1, llm)
44      ! Objet:      ! potential temperature
     ! GCM LMD nouvelle grille  
   
     ! ... Dans inigeom, nouveaux calculs pour les elongations cu, cv  
     ! et possibilite d'appeler une fonction f(y) a derivee tangente  
     ! hyperbolique a la place de la fonction a derivee sinusoidale.  
   
     ! ... Possibilité de choisir le schéma pour l'advection de  
     ! q, en modifiant iadv dans "traceur.def" (10/02) .  
   
     ! Pour Van-Leer + Vapeur d'eau saturee, iadv(1)=4. (F.Codron, 10/99)  
     ! Pour Van-Leer iadv=10  
   
     use dimens_m, only: iim, jjm, llm, nqmx  
     use paramet_m, only: ip1jmp1, ip1jm, ijmllm, ijp1llm, jjp1, iip1, iip2  
     use comconst, only: dtvr, daysec, dtphys  
     use comvert, only: ap, bp  
     use conf_gcm_m, only: day_step, iconser, idissip, iphysiq, iperiod, nday, &  
          offline, periodav  
     use logic, only: ok_guide, iflag_phys  
     use comgeom  
     use serre  
     use temps, only: itaufin, day_ini, dt  
     use iniprint, only: prt_level  
     use com_io_dyn  
     use ener  
     use calfis_m, only: calfis  
     use exner_hyb_m, only: exner_hyb  
     use guide_m, only: guide  
     use pression_m, only: pression  
     use pressure_var, only: p3d  
45    
46      integer, intent(in):: nq      REAL, intent(inout):: ps(:, :) ! (iim + 1, jjm + 1) pression au sol, en Pa
47        REAL masse((iim + 1) * (jjm + 1), llm) ! masse d'air
48        REAL phis((iim + 1) * (jjm + 1)) ! geopotentiel au sol
49    
50      ! Variables dynamiques:      REAL, intent(inout):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx)
51      REAL vcov(ip1jm, llm), ucov(ip1jmp1, llm) ! vents covariants      ! mass fractions of advected fields
     REAL teta(ip1jmp1, llm) ! temperature potentielle  
     REAL q(ip1jmp1, llm, nqmx) ! mass fractions of advected fields  
     REAL ps(ip1jmp1) ! pression au sol, en Pa  
     REAL masse(ip1jmp1, llm) ! masse d'air  
     REAL phis(ip1jmp1) ! geopotentiel au sol  
52    
53      REAL time_0      REAL, intent(in):: time_0
54    
55      ! Variables local to the procedure:      ! Variables local to the procedure:
56    
57      ! Variables dynamiques:      ! Variables dynamiques:
58    
59      REAL pks(ip1jmp1) ! exner au sol      REAL pks((iim + 1) * (jjm + 1)) ! exner au sol
60      REAL pk(ip1jmp1, llm) ! exner au milieu des couches      REAL pk(iim + 1, jjm + 1, llm) ! exner au milieu des couches
61      REAL pkf(ip1jmp1, llm) ! exner filt.au milieu des couches      REAL pkf((iim + 1) * (jjm + 1), llm) ! exner filt.au milieu des couches
62      REAL phi(ip1jmp1, llm) ! geopotential      REAL phi(iim + 1, jjm + 1, llm) ! geopotential
63      REAL w(ip1jmp1, llm) ! vitesse verticale      REAL w((iim + 1) * (jjm + 1), llm) ! vitesse verticale
64    
65      ! variables dynamiques intermediaire pour le transport      ! Variables dynamiques intermediaire pour le transport
66      REAL pbaru(ip1jmp1, llm), pbarv(ip1jm, llm) !flux de masse      ! Flux de masse :
67        REAL pbaru((iim + 1) * (jjm + 1), llm), pbarv((iim + 1) * jjm, llm)
68    
69      ! variables dynamiques au pas - 1      ! variables dynamiques au pas - 1
70      REAL vcovm1(ip1jm, llm), ucovm1(ip1jmp1, llm)      REAL vcovm1(iim + 1, jjm, llm), ucovm1(iim + 1, jjm + 1, llm)
71      REAL tetam1(ip1jmp1, llm), psm1(ip1jmp1)      REAL tetam1(iim + 1, jjm + 1, llm), psm1(iim + 1, jjm + 1)
72      REAL massem1(ip1jmp1, llm)      REAL massem1((iim + 1) * (jjm + 1), llm)
73    
74      ! tendances dynamiques      ! tendances dynamiques
75      REAL dv(ip1jm, llm), du(ip1jmp1, llm)      REAL dv((iim + 1) * jjm, llm), dudyn((iim + 1) * (jjm + 1), llm)
76      REAL dteta(ip1jmp1, llm), dq(ip1jmp1, llm, nqmx), dp(ip1jmp1)      REAL dteta(iim + 1, jjm + 1, llm), dq((iim + 1) * (jjm + 1), llm, nqmx)
77        real dp((iim + 1) * (jjm + 1))
78    
79      ! tendances de la dissipation      ! tendances de la dissipation
80      REAL dvdis(ip1jm, llm), dudis(ip1jmp1, llm)      REAL dvdis(iim + 1, jjm, llm), dudis(iim + 1, jjm + 1, llm)
81      REAL dtetadis(ip1jmp1, llm)      REAL dtetadis(iim + 1, jjm + 1, llm)
82    
83      ! tendances physiques      ! tendances physiques
84      REAL dvfi(ip1jm, llm), dufi(ip1jmp1, llm)      REAL dvfi((iim + 1) * jjm, llm), dufi((iim + 1) * (jjm + 1), llm)
85      REAL dtetafi(ip1jmp1, llm), dqfi(ip1jmp1, llm, nqmx), dpfi(ip1jmp1)      REAL dtetafi(iim + 1, jjm + 1, llm), dqfi((iim + 1) * (jjm + 1), llm, nqmx)
86        real dpfi((iim + 1) * (jjm + 1))
87    
88      ! variables pour le fichier histoire      ! variables pour le fichier histoire
89    
     REAL tppn(iim), tpps(iim), tpn, tps  
   
90      INTEGER itau ! index of the time step of the dynamics, starts at 0      INTEGER itau ! index of the time step of the dynamics, starts at 0
91      integer itaufinp1      INTEGER itaufin
     INTEGER iday ! jour julien  
92      REAL time ! time of day, as a fraction of day length      REAL time ! time of day, as a fraction of day length
93        real finvmaold((iim + 1) * (jjm + 1), llm)
94      REAL SSUM      INTEGER l
     real finvmaold(ip1jmp1, llm)  
   
     LOGICAL :: lafin=.false.  
     INTEGER ij, l  
   
95      REAL rdayvrai, rdaym_ini      REAL rdayvrai, rdaym_ini
     LOGICAL:: callinigrads = .true.  
96    
97      !+jld variables test conservation energie      ! Variables test conservation energie
98      REAL ecin(ip1jmp1, llm), ecin0(ip1jmp1, llm)      REAL ecin(iim + 1, jjm + 1, llm), ecin0(iim + 1, jjm + 1, llm)
99      ! Tendance de la temp. potentiel d (theta) / d t due a la  
100      ! tansformation d'energie cinetique en energie thermique      REAL dtetaecdt(iim + 1, jjm + 1, llm)
101      ! cree par la dissipation      ! tendance de la température potentielle due à la transformation
102      REAL dtetaecdt(ip1jmp1, llm)      ! d'énergie cinétique en énergie thermique par la dissipation
103      REAL vcont(ip1jm, llm), ucont(ip1jmp1, llm)  
104      CHARACTER*15 ztit      REAL vcont((iim + 1) * jjm, llm), ucont((iim + 1) * (jjm + 1), llm)
105      INTEGER:: ip_ebil_dyn = 0 ! PRINT level for energy conserv. diag.      logical leapf
106        real dt
     logical:: dissip_conservative = .true.  
     LOGICAL:: prem = .true.  
     logical forward, leapf, apphys, conser, apdiss  
107    
108      !---------------------------------------------------      !---------------------------------------------------
109    
110      print *, "Call sequence information: leapfrog"      print *, "Call sequence information: leapfrog"
111        call assert(shape(ucov) == (/iim + 1, jjm + 1, llm/), "leapfrog")
112    
113      itaufin = nday * day_step      itaufin = nday * day_step
114      itaufinp1 = itaufin + 1      ! "day_step" is a multiple of "iperiod", therefore "itaufin" is one too
115    
116      itau = 0      dq = 0.
     iday = day_ini  
     time = time_0  
     IF (time > 1.) THEN  
        time = time - 1.  
        iday = iday + 1  
     ENDIF  
117    
118      ! On initialise la pression et la fonction d'Exner :      ! On initialise la pression et la fonction d'Exner :
119      dq=0.      forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
     CALL pression(ip1jmp1, ap, bp, ps, p3d)  
120      CALL exner_hyb(ps, p3d, pks, pk, pkf)      CALL exner_hyb(ps, p3d, pks, pk, pkf)
121    
122      ! Debut de l'integration temporelle:      time_integration: do itau = 0, itaufin - 1
123      outer_loop:do         leapf = mod(itau, iperiod) /= 0
124         if (ok_guide.and.(itaufin - itau - 1) * dtvr > 21600) then         if (leapf) then
125            call guide(itau, ucov, vcov, teta, q, masse, ps)            dt = 2 * dtvr
126         else         else
127            IF (prt_level > 9) print *, &            ! Matsuno
128                 'Attention : on ne guide pas les 6 dernieres heures.'            dt = dtvr
129         endif            if (ok_guide .and. (itaufin - itau - 1) * dtvr > 21600.) &
130                   call guide(itau, ucov, vcov, teta, q, masse, ps)
131         CALL SCOPY(ijmllm, vcov, 1, vcovm1, 1)            vcovm1 = vcov
132         CALL SCOPY(ijp1llm, ucov, 1, ucovm1, 1)            ucovm1 = ucov
133         CALL SCOPY(ijp1llm, teta, 1, tetam1, 1)            tetam1 = teta
134         CALL SCOPY(ijp1llm, masse, 1, massem1, 1)            massem1 = masse
135         CALL SCOPY(ip1jmp1, ps, 1, psm1, 1)            psm1 = ps
136              finvmaold = masse
137         forward = .TRUE.            CALL filtreg(finvmaold, jjm + 1, llm, - 2, 2, .TRUE., 1)
138         leapf = .FALSE.         end if
139         dt = dtvr  
140           ! Calcul des tendances dynamiques:
141         CALL SCOPY(ijp1llm, masse, 1, finvmaold, 1)         CALL geopot((iim + 1) * (jjm + 1), teta, pk, pks, phis, phi)
142         CALL filtreg(finvmaold, jjp1, llm, - 2, 2, .TRUE., 1)         CALL caldyn(itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
143                dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &
144         do              conser=MOD(itau, iconser)==0)
145            ! gestion des appels de la physique et des dissipations:  
146           ! Calcul des tendances advection des traceurs (dont l'humidité)
147            apphys = .FALSE.         CALL caladvtrac(q, pbaru, pbarv, p3d, masse, dq, teta, pk)
148            conser = .FALSE.  
149            apdiss = .FALSE.         ! Stokage du flux de masse pour traceurs offline:
150           IF (offline) CALL fluxstokenc(pbaru, pbarv, masse, teta, phi, phis, &
151            IF (MOD(itau, iconser) == 0) conser = .TRUE.              dtvr, itau)
152            IF (MOD(itau + 1, idissip) == 0) apdiss = .TRUE.  
153            IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys /= 0) apphys=.TRUE.         ! Integrations dynamique et traceurs:
154           CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, dteta, &
155            ! calcul des tendances dynamiques:              dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, finvmaold, dt, &
156                leapf)
157            CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)  
158           if (.not. leapf) then
159            CALL caldyn(itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &            ! Matsuno backward
160                 conser, du, dv, dteta, dp, w, pbaru, pbarv, &            forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
161                 time + iday - day_ini)            CALL exner_hyb(ps, p3d, pks, pk, pkf)
162    
163            ! calcul des tendances advection des traceurs (dont l'humidite)            ! Calcul des tendances dynamiques:
164              CALL geopot((iim + 1) * (jjm + 1), teta, pk, pks, phis, phi)
165            IF (forward .OR. leapf) THEN            CALL caldyn(itau + 1, ucov, vcov, teta, ps, masse, pk, pkf, phis, &
166               CALL caladvtrac(q, pbaru, pbarv, p3d, masse, dq, teta, pk)                 phi, dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &
167               IF (offline) THEN                 conser=.false.)
                 !maf stokage du flux de masse pour traceurs OFF-LINE  
                 CALL fluxstokenc(pbaru, pbarv, masse, teta, phi, phis, dtvr, &  
                      itau)  
              ENDIF  
           ENDIF  
168    
169            ! integrations dynamique et traceurs:            ! integrations dynamique et traceurs:
170            CALL integrd(2, vcovm1, ucovm1, tetam1, psm1, massem1, dv, du, &            CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, &
171                 dteta, dq, dp, vcov, ucov, teta, q, ps, masse, phis, &                 dteta, dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, &
172                 finvmaold, leapf)                 finvmaold, dtvr, leapf=.false.)
173           end if
174    
175           IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys /= 0) THEN
176            ! calcul des tendances physiques:            ! calcul des tendances physiques:
177    
178            IF (apphys) THEN            forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
              IF (itau + 1 == itaufin) lafin = .TRUE.  
   
              CALL pression(ip1jmp1, ap, bp, ps, p3d)  
              CALL exner_hyb(ps, p3d, pks, pk, pkf)  
   
              rdaym_ini = itau * dtvr / daysec  
              rdayvrai = rdaym_ini + day_ini  
   
              ! Interface avec les routines de phylmd (phymars ...)  
   
              ! Diagnostique de conservation de l'énergie : initialisation  
              IF (ip_ebil_dyn >= 1) THEN  
                 ztit='bil dyn'  
                 CALL diagedyn(ztit, 2, 1, 1, dtphys, ucov, vcov, ps, p3d, pk, &  
                      teta, q(:, :, 1), q(:, :, 2))  
              ENDIF  
   
              CALL calfis(nq, lafin, rdayvrai, time, ucov, vcov, teta, q, &  
                   masse, ps, pk, phis, phi, du, dv, dteta, dq, w, &  
                   dufi, dvfi, dtetafi, dqfi, dpfi)  
   
              ! ajout des tendances physiques:  
              CALL addfi(nqmx, dtphys, &  
                   ucov, vcov, teta, q, ps, &  
                   dufi, dvfi, dtetafi, dqfi, dpfi)  
   
              ! Diagnostique de conservation de l'énergie : difference  
              IF (ip_ebil_dyn >= 1) THEN  
                 ztit = 'bil phys'  
                 CALL diagedyn(ztit, 2, 1, 1, dtphys, ucov, vcov, ps, p3d, pk, &  
                      teta, q(:, :, 1), q(:, :, 2))  
              ENDIF  
           ENDIF  
   
           CALL pression(ip1jmp1, ap, bp, ps, p3d)  
179            CALL exner_hyb(ps, p3d, pks, pk, pkf)            CALL exner_hyb(ps, p3d, pks, pk, pkf)
180    
181            ! dissipation horizontale et verticale des petites echelles:            rdaym_ini = itau * dtvr / daysec
182              rdayvrai = rdaym_ini + day_ini
183            IF (apdiss) THEN            time = REAL(mod(itau, day_step)) / day_step + time_0
184               ! calcul de l'energie cinetique avant dissipation            IF (time > 1.) time = time - 1.
185               call covcont(llm, ucov, vcov, ucont, vcont)  
186               call enercin(vcov, ucov, vcont, ucont, ecin0)            CALL calfis(rdayvrai, time, ucov, vcov, teta, q, masse, ps, pk, &
187                   phis, phi, dudyn, dv, dq, w, dufi, dvfi, dtetafi, dqfi, dpfi, &
188               ! dissipation                 lafin=itau+1==itaufin)
189               CALL dissip(vcov, ucov, teta, p3d, dvdis, dudis, dtetadis)  
190               ucov=ucov + dudis            ! ajout des tendances physiques:
191               vcov=vcov + dvdis            CALL addfi(nqmx, dtphys, ucov, vcov, teta, q, ps, dufi, dvfi, &
192                   dtetafi, dqfi, dpfi)
193               if (dissip_conservative) then         ENDIF
194                  ! On rajoute la tendance due a la transform. Ec -> E  
195                  ! therm. cree lors de la dissipation         forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
196                  call covcont(llm, ucov, vcov, ucont, vcont)         CALL exner_hyb(ps, p3d, pks, pk, pkf)
197                  call enercin(vcov, ucov, vcont, ucont, ecin)  
198                  dtetaecdt= (ecin0 - ecin) / pk         IF (MOD(itau + 1, idissip) == 0) THEN
199                  dtetadis=dtetadis + dtetaecdt            ! Dissipation horizontale et verticale des petites échelles
200               endif  
201               teta=teta + dtetadis            ! calcul de l'énergie cinétique avant dissipation
202              call covcont(llm, ucov, vcov, ucont, vcont)
203               ! Calcul de la valeur moyenne, unique de h aux poles .....            call enercin(vcov, ucov, vcont, ucont, ecin0)
204    
205               DO l = 1, llm            ! dissipation
206                  DO ij = 1, iim            CALL dissip(vcov, ucov, teta, p3d, dvdis, dudis, dtetadis)
207                     tppn(ij) = aire(ij) * teta(ij, l)            ucov = ucov + dudis
208                     tpps(ij) = aire(ij + ip1jm) * teta(ij + ip1jm, l)            vcov = vcov + dvdis
209                  ENDDO  
210                  tpn = SSUM(iim, tppn, 1) / apoln            ! On ajoute la tendance due à la transformation énergie
211                  tps = SSUM(iim, tpps, 1) / apols            ! cinétique en énergie thermique par la dissipation
212              call covcont(llm, ucov, vcov, ucont, vcont)
213                  DO ij = 1, iip1            call enercin(vcov, ucov, vcont, ucont, ecin)
214                     teta(ij, l) = tpn            dtetaecdt= (ecin0 - ecin) / pk
215                     teta(ij + ip1jm, l) = tps            dtetadis = dtetadis + dtetaecdt
216                  ENDDO            teta = teta + dtetadis
217               ENDDO  
218              ! Calcul de la valeur moyenne aux pôles :
219               DO ij = 1, iim            forall (l = 1: llm)
220                  tppn(ij) = aire(ij) * ps(ij)               teta(:, 1, l) = SUM(aire_2d(:iim, 1) * teta(:iim, 1, l)) &
221                  tpps(ij) = aire(ij + ip1jm) * ps(ij + ip1jm)                    / apoln
222               ENDDO               teta(:, jjm + 1, l) = SUM(aire_2d(:iim, jjm+1) &
223               tpn = SSUM(iim, tppn, 1) / apoln                    * teta(:iim, jjm + 1, l)) / apols
224               tps = SSUM(iim, tpps, 1) / apols            END forall
225    
226               DO ij = 1, iip1            ps(:, 1) = SUM(aire_2d(:iim, 1) * ps(:iim, 1)) / apoln
227                  ps(ij) = tpn            ps(:, jjm + 1) = SUM(aire_2d(:iim, jjm+1) * ps(:iim, jjm + 1)) &
228                  ps(ij + ip1jm) = tps                 / apols
229               ENDDO         END IF
230    
231            END IF         IF (MOD(itau + 1, iperiod) == 0) THEN
232              ! Écriture du fichier histoire moyenne:
233            ! fin de l'intégration dynamique et physique pour le pas "itau"            CALL writedynav(histaveid, nqmx, itau + 1, vcov, ucov, teta, pk, &
234            ! préparation du pas d'intégration suivant                 phi, q, masse, ps, phis)
235              call bilan_dyn(ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, &
236            ! schema matsuno + leapfrog                 q(:, :, :, 1), dt_app = dtvr * iperiod, &
237            IF (forward .OR. leapf) THEN                 dt_cum = dtvr * day_step * periodav)
238               itau = itau + 1         ENDIF
239               iday = day_ini + itau / day_step      end do time_integration
240               time = REAL(itau - (iday - day_ini) * day_step) / day_step &  
241                    + time_0      CALL dynredem1("restart.nc", vcov, ucov, teta, q, masse, ps, &
242               IF (time > 1.) THEN           itau=itau_dyn+itaufin)
243                  time = time - 1.  
244                  iday = iday + 1      ! Calcul des tendances dynamiques:
245               ENDIF      CALL geopot((iim + 1) * (jjm + 1), teta, pk, pks, phis, phi)
246            ENDIF      CALL caldyn(itaufin, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
247             dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &
248            IF (itau == itaufinp1) exit outer_loop           conser=MOD(itaufin, iconser)==0)
   
           ! ecriture du fichier histoire moyenne:  
   
           ! Comment out the following calls when you do not want the output  
           ! files "dyn_hist_ave.nc" and "dynzon.nc"  
           IF (MOD(itau, iperiod) == 0 .OR. itau == itaufin) THEN  
              CALL writedynav(histaveid, nqmx, itau, vcov, &  
                   ucov, teta, pk, phi, q, masse, ps, phis)  
              call bilan_dyn(2, dtvr * iperiod, dtvr * day_step * periodav, &  
                   ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, q)  
           ENDIF  
   
           IF (itau == itaufin) THEN  
              CALL dynredem1("restart.nc", 0., vcov, ucov, teta, q, masse, ps)  
              CLOSE(99)  
           ENDIF  
   
           ! gestion de l'integration temporelle:  
   
           IF (MOD(itau, iperiod) == 0) exit  
           IF (MOD(itau - 1, iperiod) == 0) THEN  
              IF (forward) THEN  
                 ! fin du pas forward et debut du pas backward  
                 forward = .FALSE.  
                 leapf = .FALSE.  
              ELSE  
                 ! fin du pas backward et debut du premier pas leapfrog  
                 leapf = .TRUE.  
                 dt = 2. * dtvr  
              END IF  
           ELSE  
              ! pas leapfrog  
              leapf = .TRUE.  
              dt = 2. * dtvr  
           END IF  
        end do  
     end do outer_loop  
   
249    END SUBROUTINE leapfrog    END SUBROUTINE leapfrog
250    
251  end module leapfrog_m  end module leapfrog_m

Legend:
Removed from v.22  
changed lines
  Added in v.55

  ViewVC Help
Powered by ViewVC 1.1.21