/[lmdze]/trunk/Sources/dyn3d/leapfrog.f
ViewVC logotype

Diff of /trunk/Sources/dyn3d/leapfrog.f

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

revision 47 by guez, Fri Jul 1 15:00:48 2011 UTC revision 70 by guez, Mon Jun 24 15:39:52 2013 UTC
# Line 6  contains Line 6  contains
6    
7    SUBROUTINE leapfrog(ucov, vcov, teta, ps, masse, phis, 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 revision 616
10      ! Authors: P. Le Van, L. Fairhead, F. Hourdin      ! Authors: P. Le Van, L. Fairhead, F. Hourdin
11      ! Matsuno-leapfrog scheme.      ! Matsuno-leapfrog scheme.
12    
# Line 15  contains Line 15  contains
15      use caladvtrac_m, only: caladvtrac      use caladvtrac_m, only: caladvtrac
16      use caldyn_m, only: caldyn      use caldyn_m, only: caldyn
17      USE calfis_m, ONLY: calfis      USE calfis_m, ONLY: calfis
     USE com_io_dyn, ONLY: histaveid  
18      USE comconst, ONLY: daysec, dtphys, dtvr      USE comconst, ONLY: daysec, dtphys, dtvr
19      USE comgeom, ONLY: aire_2d, apoln, apols      USE comgeom, ONLY: aire_2d, apoln, apols
20      USE comvert, ONLY: ap, bp      USE disvert_m, ONLY: ap, bp
21      USE conf_gcm_m, ONLY: day_step, iconser, iperiod, iphysiq, nday, offline, &      USE conf_gcm_m, ONLY: day_step, iconser, iperiod, iphysiq, nday, offline, &
22           periodav           iflag_phys, ok_guide, iecri
23      USE dimens_m, ONLY: iim, jjm, llm, nqmx      USE dimens_m, ONLY: iim, jjm, llm, nqmx
24      use dissip_m, only: dissip      use dissip_m, only: dissip
25      USE dynetat0_m, ONLY: day_ini      USE dynetat0_m, ONLY: day_ini
26      use dynredem1_m, only: dynredem1      use dynredem1_m, only: dynredem1
27      USE exner_hyb_m, ONLY: exner_hyb      USE exner_hyb_m, ONLY: exner_hyb
28      use filtreg_m, only: filtreg      use filtreg_m, only: filtreg
29        use fluxstokenc_m, only: fluxstokenc
30      use geopot_m, only: geopot      use geopot_m, only: geopot
31      USE guide_m, ONLY: guide      USE guide_m, ONLY: guide
32      use inidissip_m, only: idissip      use inidissip_m, only: idissip
33      use integrd_m, only: integrd      use integrd_m, only: integrd
34      USE logic, ONLY: iflag_phys, ok_guide      use nr_util, only: assert
     USE paramet_m, ONLY: ip1jmp1  
35      USE pressure_var, ONLY: p3d      USE pressure_var, ONLY: p3d
36      USE temps, ONLY: itau_dyn      USE temps, ONLY: itau_dyn
37        use writedynav_m, only: writedynav
38        use writehist_m, only: writehist
39    
40      ! Variables dynamiques:      ! Variables dynamiques:
41      REAL, intent(inout):: ucov(ip1jmp1, llm) ! vent covariant      REAL, intent(inout):: ucov(:, :, :) ! (iim + 1, jjm + 1, llm) vent covariant
42      REAL, intent(inout):: vcov((iim + 1) * jjm, llm) ! vent covariant      REAL, intent(inout):: vcov(:, :, :) ! (iim + 1, jjm, llm) ! vent covariant
43    
44      REAL, intent(inout):: teta(:, :, :) ! (iim + 1, jjm + 1, llm)      REAL, intent(inout):: teta(:, :, :) ! (iim + 1, jjm + 1, llm)
45      ! potential temperature      ! potential temperature
46    
47      REAL, intent(inout):: ps(:, :) ! (iim + 1, jjm + 1) pression au sol, en Pa      REAL, intent(inout):: ps(:, :) ! (iim + 1, jjm + 1) pression au sol, en Pa
48      REAL masse(ip1jmp1, llm) ! masse d'air      REAL, intent(inout):: masse(:, :, :) ! (iim + 1, jjm + 1, llm) masse d'air
49      REAL phis(ip1jmp1) ! geopotentiel au sol      REAL, intent(in):: phis(:, :) ! (iim + 1, jjm + 1) surface geopotential
50    
51      REAL, intent(inout):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx)      REAL, intent(inout):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx)
52      ! mass fractions of advected fields      ! mass fractions of advected fields
# Line 56  contains Line 57  contains
57    
58      ! Variables dynamiques:      ! Variables dynamiques:
59    
60      REAL pks(ip1jmp1) ! exner au sol      REAL pks(iim + 1, jjm + 1) ! exner au sol
61      REAL pk(iim + 1, jjm + 1, llm) ! exner au milieu des couches      REAL pk(iim + 1, jjm + 1, llm) ! exner au milieu des couches
62      REAL pkf(ip1jmp1, llm) ! exner filt.au milieu des couches      REAL pkf(iim + 1, jjm + 1, llm) ! exner filtré au milieu des couches
63      REAL phi(iim + 1, jjm + 1, llm) ! geopotential      REAL phi(iim + 1, jjm + 1, llm) ! geopotential
64      REAL w(ip1jmp1, llm) ! vitesse verticale      REAL w((iim + 1) * (jjm + 1), llm) ! vitesse verticale
65    
66      ! variables dynamiques intermediaire pour le transport      ! Variables dynamiques intermediaire pour le transport
67      REAL pbaru(ip1jmp1, llm), pbarv((iim + 1) * jjm, llm) !flux de masse      ! Flux de masse :
68        REAL pbaru((iim + 1) * (jjm + 1), llm), pbarv((iim + 1) * jjm, llm)
69    
70      ! variables dynamiques au pas - 1      ! Variables dynamiques au pas - 1
71      REAL vcovm1((iim + 1) * jjm, llm), ucovm1(ip1jmp1, llm)      REAL vcovm1(iim + 1, jjm, llm), ucovm1(iim + 1, jjm + 1, llm)
72      REAL tetam1(iim + 1, jjm + 1, llm), psm1(iim + 1, jjm + 1)      REAL tetam1(iim + 1, jjm + 1, llm), psm1(iim + 1, jjm + 1)
73      REAL massem1(ip1jmp1, llm)      REAL massem1(iim + 1, jjm + 1, llm)
74    
75      ! tendances dynamiques      ! Tendances dynamiques
76      REAL dv((iim + 1) * jjm, llm), dudyn(ip1jmp1, llm)      REAL dv((iim + 1) * jjm, llm), dudyn((iim + 1) * (jjm + 1), llm)
77      REAL dteta(iim + 1, jjm + 1, llm), dq(ip1jmp1, llm, nqmx), dp(ip1jmp1)      REAL dteta(iim + 1, jjm + 1, llm), dq((iim + 1) * (jjm + 1), llm, nqmx)
78        real dp((iim + 1) * (jjm + 1))
79    
80      ! tendances de la dissipation      ! Tendances de la dissipation :
81      REAL dvdis((iim + 1) * jjm, llm), dudis(ip1jmp1, llm)      REAL dvdis(iim + 1, jjm, llm), dudis(iim + 1, jjm + 1, llm)
82      REAL dtetadis(iim + 1, jjm + 1, llm)      REAL dtetadis(iim + 1, jjm + 1, llm)
83    
84      ! tendances physiques      ! Tendances physiques
85      REAL dvfi((iim + 1) * jjm, llm), dufi(ip1jmp1, llm)      REAL dvfi((iim + 1) * jjm, llm), dufi((iim + 1) * (jjm + 1), llm)
86      REAL dtetafi(iim + 1, jjm + 1, llm), dqfi(ip1jmp1, llm, nqmx), dpfi(ip1jmp1)      REAL dtetafi(iim + 1, jjm + 1, llm), dqfi((iim + 1) * (jjm + 1), llm, nqmx)
87        real dpfi((iim + 1) * (jjm + 1))
88    
89      ! variables pour le fichier histoire      ! Variables pour le fichier histoire
90    
91      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
92      INTEGER itaufin      INTEGER itaufin
93      REAL time ! time of day, as a fraction of day length      REAL time ! time of day, as a fraction of day length
94      real finvmaold(ip1jmp1, llm)      real finvmaold(iim + 1, jjm + 1, llm)
95      INTEGER l      INTEGER l
96      REAL rdayvrai, rdaym_ini      REAL rdayvrai, rdaym_ini
97    
98      ! Variables test conservation energie      ! Variables test conservation energie
99      REAL ecin(iim + 1, jjm + 1, llm), ecin0(iim + 1, jjm + 1, llm)      REAL ecin(iim + 1, jjm + 1, llm), ecin0(iim + 1, jjm + 1, llm)
100    
101      REAL dtetaecdt(iim + 1, jjm + 1, llm)      REAL vcont((iim + 1) * jjm, llm), ucont((iim + 1) * (jjm + 1), llm)
     ! tendance de la température potentielle due à la tansformation  
     ! d'énergie cinétique en énergie thermique créée par la dissipation  
   
     REAL vcont((iim + 1) * jjm, llm), ucont(ip1jmp1, llm)  
102      logical leapf      logical leapf
103      real dt      real dt
104    
105      !---------------------------------------------------      !---------------------------------------------------
106    
107      print *, "Call sequence information: leapfrog"      print *, "Call sequence information: leapfrog"
108        call assert(shape(ucov) == (/iim + 1, jjm + 1, llm/), "leapfrog")
109    
110      itaufin = nday * day_step      itaufin = nday * day_step
111      ! "day_step" is a multiple of "iperiod", therefore "itaufin" is one too      ! "day_step" is a multiple of "iperiod", therefore so is "itaufin".
112    
113      dq = 0.      dq = 0.
114    
# Line 130  contains Line 131  contains
131            massem1 = masse            massem1 = masse
132            psm1 = ps            psm1 = ps
133            finvmaold = masse            finvmaold = masse
134            CALL filtreg(finvmaold, jjm + 1, llm, - 2, 2, .TRUE., 1)            CALL filtreg(finvmaold, jjm + 1, llm, - 2, 2, .TRUE.)
135         end if         end if
136    
137         ! Calcul des tendances dynamiques:         ! Calcul des tendances dynamiques:
138         CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)         CALL geopot(teta, pk, pks, phis, phi)
139         CALL caldyn(itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &         CALL caldyn(itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
140              dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &              dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &
141              conser=MOD(itau, iconser)==0)              conser=MOD(itau, iconser)==0)
# Line 146  contains Line 147  contains
147         IF (offline) CALL fluxstokenc(pbaru, pbarv, masse, teta, phi, phis, &         IF (offline) CALL fluxstokenc(pbaru, pbarv, masse, teta, phi, phis, &
148              dtvr, itau)              dtvr, itau)
149    
150         ! Integrations dynamique et traceurs:         ! Intégrations dynamique et traceurs:
151         CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, dteta, &         CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, dteta, &
152              dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, finvmaold, dt, &              dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, finvmaold, dt, &
153              leapf)              leapf)
# Line 157  contains Line 158  contains
158            CALL exner_hyb(ps, p3d, pks, pk, pkf)            CALL exner_hyb(ps, p3d, pks, pk, pkf)
159    
160            ! Calcul des tendances dynamiques:            ! Calcul des tendances dynamiques:
161            CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)            CALL geopot(teta, pk, pks, phis, phi)
162            CALL caldyn(itau + 1, ucov, vcov, teta, ps, masse, pk, pkf, phis, &            CALL caldyn(itau + 1, ucov, vcov, teta, ps, masse, pk, pkf, phis, &
163                 phi, dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &                 phi, dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &
164                 conser=.false.)                 conser=.false.)
# Line 169  contains Line 170  contains
170         end if         end if
171    
172         IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys /= 0) THEN         IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys /= 0) THEN
173            ! calcul des tendances physiques:            ! Calcul des tendances physiques:
174    
175            forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps            forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
176            CALL exner_hyb(ps, p3d, pks, pk, pkf)            CALL exner_hyb(ps, p3d, pks, pk, pkf)
# Line 181  contains Line 182  contains
182    
183            CALL calfis(rdayvrai, time, ucov, vcov, teta, q, masse, ps, pk, &            CALL calfis(rdayvrai, time, ucov, vcov, teta, q, masse, ps, pk, &
184                 phis, phi, dudyn, dv, dq, w, dufi, dvfi, dtetafi, dqfi, dpfi, &                 phis, phi, dudyn, dv, dq, w, dufi, dvfi, dtetafi, dqfi, dpfi, &
185                 lafin=itau+1==itaufin)                 lafin = itau + 1 == itaufin)
186    
187            ! ajout des tendances physiques:            ! Ajout des tendances physiques:
188            CALL addfi(nqmx, dtphys, ucov, vcov, teta, q, ps, dufi, dvfi, &            CALL addfi(nqmx, ucov, vcov, teta, q, ps, dufi, dvfi, dtetafi, &
189                 dtetafi, dqfi, dpfi)                 dqfi, dpfi)
190         ENDIF         ENDIF
191    
192         forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps         forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
193         CALL exner_hyb(ps, p3d, pks, pk, pkf)         CALL exner_hyb(ps, p3d, pks, pk, pkf)
194    
195         IF (MOD(itau + 1, idissip) == 0) THEN         IF (MOD(itau + 1, idissip) == 0) THEN
196            ! dissipation horizontale et verticale des petites echelles:            ! Dissipation horizontale et verticale des petites échelles
197    
198            ! calcul de l'energie cinetique avant dissipation            ! calcul de l'énergie cinétique avant dissipation
199            call covcont(llm, ucov, vcov, ucont, vcont)            call covcont(llm, ucov, vcov, ucont, vcont)
200            call enercin(vcov, ucov, vcont, ucont, ecin0)            call enercin(vcov, ucov, vcont, ucont, ecin0)
201    
202            ! dissipation            ! dissipation
203            CALL dissip(vcov, ucov, teta, p3d, dvdis, dudis, dtetadis)            CALL dissip(vcov, ucov, teta, p3d, dvdis, dudis, dtetadis)
204            ucov=ucov + dudis            ucov = ucov + dudis
205            vcov=vcov + dvdis            vcov = vcov + dvdis
206    
207            ! On rajoute la tendance due à la transformation Ec -> E            ! On ajoute la tendance due à la transformation énergie
208            ! thermique créée lors de la dissipation            ! cinétique en énergie thermique par la dissipation
209            call covcont(llm, ucov, vcov, ucont, vcont)            call covcont(llm, ucov, vcov, ucont, vcont)
210            call enercin(vcov, ucov, vcont, ucont, ecin)            call enercin(vcov, ucov, vcont, ucont, ecin)
211            dtetaecdt= (ecin0 - ecin) / pk            dtetadis = dtetadis + (ecin0 - ecin) / pk
212            dtetadis=dtetadis + dtetaecdt            teta = teta + dtetadis
           teta=teta + dtetadis  
213    
214            ! Calcul de la valeur moyenne aux pôles :            ! Calcul de la valeur moyenne aux pôles :
215            forall (l = 1: llm)            forall (l = 1: llm)
# Line 226  contains Line 226  contains
226    
227         IF (MOD(itau + 1, iperiod) == 0) THEN         IF (MOD(itau + 1, iperiod) == 0) THEN
228            ! Écriture du fichier histoire moyenne:            ! Écriture du fichier histoire moyenne:
229            CALL writedynav(histaveid, nqmx, itau + 1, vcov, ucov, teta, pk, &            CALL writedynav(vcov, ucov, teta, pk, phi, q, masse, ps, phis, &
230                 phi, q, masse, ps, phis)                 time = itau + 1)
231            call bilan_dyn(ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, &            call bilan_dyn(ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, &
232                 q(:, :, :, 1), dt_app = dtvr * iperiod, &                 q(:, :, :, 1))
                dt_cum = dtvr * day_step * periodav)  
233         ENDIF         ENDIF
234    
235           IF (MOD(itau + 1, iecri * day_step) == 0) THEN
236              CALL geopot(teta, pk, pks, phis, phi)
237              CALL writehist(itau, vcov, ucov, teta, phi, q, masse, ps)
238           END IF
239      end do time_integration      end do time_integration
240    
241      CALL dynredem1("restart.nc", vcov, ucov, teta, q, masse, ps, &      CALL dynredem1("restart.nc", vcov, ucov, teta, q, masse, ps, &
242           itau=itau_dyn+itaufin)           itau = itau_dyn + itaufin)
243    
244      ! Calcul des tendances dynamiques:      ! Calcul des tendances dynamiques:
245      CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)      CALL geopot(teta, pk, pks, phis, phi)
246      CALL caldyn(itaufin, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &      CALL caldyn(itaufin, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
247           dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &           dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &
248           conser=MOD(itaufin, iconser)==0)           conser = MOD(itaufin, iconser) == 0)
249    
250    END SUBROUTINE leapfrog    END SUBROUTINE leapfrog
251    
252  end module leapfrog_m  end module leapfrog_m

Legend:
Removed from v.47  
changed lines
  Added in v.70

  ViewVC Help
Powered by ViewVC 1.1.21