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

Diff of /trunk/dyn3d/leapfrog.f

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

trunk/libf/dyn3d/leapfrog.f90 revision 47 by guez, Fri Jul 1 15:00:48 2011 UTC trunk/Sources/dyn3d/leapfrog.f revision 224 by guez, Fri Apr 28 13:40:59 2017 UTC
# Line 4  module leapfrog_m Line 4  module leapfrog_m
4    
5  contains  contains
6    
7    SUBROUTINE leapfrog(ucov, vcov, teta, ps, masse, phis, q, time_0)    SUBROUTINE leapfrog(ucov, vcov, teta, ps, masse, phis, q)
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.  
12        ! Intégration temporelle du modèle : Matsuno-leapfrog scheme.
13    
14      use addfi_m, only: addfi      use addfi_m, only: addfi
15      use bilan_dyn_m, only: bilan_dyn      use bilan_dyn_m, only: bilan_dyn
16      use caladvtrac_m, only: caladvtrac      use caladvtrac_m, only: caladvtrac
17      use caldyn_m, only: caldyn      use caldyn_m, only: caldyn
18      USE calfis_m, ONLY: calfis      USE calfis_m, ONLY: calfis
19      USE com_io_dyn, ONLY: histaveid      USE comconst, ONLY: dtvr
     USE comconst, ONLY: daysec, dtphys, dtvr  
20      USE comgeom, ONLY: aire_2d, apoln, apols      USE comgeom, ONLY: aire_2d, apoln, apols
21      USE comvert, ONLY: ap, bp      use covcont_m, only: covcont
22      USE conf_gcm_m, ONLY: day_step, iconser, iperiod, iphysiq, nday, offline, &      USE disvert_m, ONLY: ap, bp
23           periodav      USE conf_gcm_m, ONLY: day_step, iconser, iperiod, iphysiq, nday, &
24             iflag_phys, iecri
25        USE conf_guide_m, ONLY: ok_guide
26      USE dimens_m, ONLY: iim, jjm, llm, nqmx      USE dimens_m, ONLY: iim, jjm, llm, nqmx
27      use dissip_m, only: dissip      use dissip_m, only: dissip
28      USE dynetat0_m, ONLY: day_ini      USE dynetat0_m, ONLY: day_ini
29      use dynredem1_m, only: dynredem1      use dynredem1_m, only: dynredem1
30        use enercin_m, only: enercin
31      USE exner_hyb_m, ONLY: exner_hyb      USE exner_hyb_m, ONLY: exner_hyb
32      use filtreg_m, only: filtreg      use filtreg_scal_m, only: filtreg_scal
33      use geopot_m, only: geopot      use geopot_m, only: geopot
34      USE guide_m, ONLY: guide      USE guide_m, ONLY: guide
35      use inidissip_m, only: idissip      use inidissip_m, only: idissip
36      use integrd_m, only: integrd      use integrd_m, only: integrd
37      USE logic, ONLY: iflag_phys, ok_guide      use nr_util, only: assert
     USE paramet_m, ONLY: ip1jmp1  
     USE pressure_var, ONLY: p3d  
38      USE temps, ONLY: itau_dyn      USE temps, ONLY: itau_dyn
39        use writedynav_m, only: writedynav
40        use writehist_m, only: writehist
41    
42      ! Variables dynamiques:      ! Variables dynamiques:
43      REAL, intent(inout):: ucov(ip1jmp1, llm) ! vent covariant      REAL, intent(inout):: ucov(:, :, :) ! (iim + 1, jjm + 1, llm) vent covariant
44      REAL, intent(inout):: vcov((iim + 1) * jjm, llm) ! vent covariant      REAL, intent(inout):: vcov(:, :, :) ! (iim + 1, jjm, llm) ! vent covariant
45    
46      REAL, intent(inout):: teta(:, :, :) ! (iim + 1, jjm + 1, llm)      REAL, intent(inout):: teta(:, :, :) ! (iim + 1, jjm + 1, llm)
47      ! potential temperature      ! potential temperature
48    
49      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
50      REAL masse(ip1jmp1, llm) ! masse d'air      REAL, intent(inout):: masse(:, :, :) ! (iim + 1, jjm + 1, llm) masse d'air
51      REAL phis(ip1jmp1) ! geopotentiel au sol      REAL, intent(in):: phis(:, :) ! (iim + 1, jjm + 1) surface geopotential
52    
53      REAL, intent(inout):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx)      REAL, intent(inout):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx)
54      ! mass fractions of advected fields      ! mass fractions of advected fields
55    
56      REAL, intent(in):: time_0      ! Local:
   
     ! Variables local to the procedure:  
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\'e 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)
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)
   
     ! variables pour le fichier histoire  
87    
88        ! Variables pour le fichier histoire
89      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
90      INTEGER itaufin      INTEGER itaufin
     REAL time ! time of day, as a fraction of day length  
     real finvmaold(ip1jmp1, llm)  
91      INTEGER l      INTEGER l
     REAL rdayvrai, rdaym_ini  
92    
93      ! Variables test conservation energie      ! Variables test conservation \'energie
94      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)
95    
96      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)  
97      logical leapf      logical leapf
98      real dt      real dt ! time step, in s
99    
100        REAL p3d(iim + 1, jjm + 1, llm + 1) ! pressure at layer interfaces, in Pa
101        ! ("p3d(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)",
102        ! for interface "l")
103    
104      !---------------------------------------------------      !---------------------------------------------------
105    
106      print *, "Call sequence information: leapfrog"      print *, "Call sequence information: leapfrog"
107        call assert(shape(ucov) == (/iim + 1, jjm + 1, llm/), "leapfrog")
108    
109      itaufin = nday * day_step      itaufin = nday * day_step
110      ! "day_step" is a multiple of "iperiod", therefore "itaufin" is one too      ! "day_step" is a multiple of "iperiod", therefore so is "itaufin".
   
     dq = 0.  
111    
112      ! On initialise la pression et la fonction d'Exner :      ! On initialise la pression et la fonction d'Exner :
113      forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps      forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
114      CALL exner_hyb(ps, p3d, pks, pk, pkf)      CALL exner_hyb(ps, p3d, pks, pk)
115        pkf = pk
116        CALL filtreg_scal(pkf, direct = .true., intensive = .true.)
117    
118      time_integration: do itau = 0, itaufin - 1      time_integration: do itau = 0, itaufin - 1
119         leapf = mod(itau, iperiod) /= 0         leapf = mod(itau, iperiod) /= 0
# Line 122  contains Line 122  contains
122         else         else
123            ! Matsuno            ! Matsuno
124            dt = dtvr            dt = dtvr
125            if (ok_guide .and. (itaufin - itau - 1) * dtvr > 21600.) &            if (ok_guide) call guide(itau, ucov, vcov, teta, q(:, :, :, 1), ps)
                call guide(itau, ucov, vcov, teta, q, masse, ps)  
126            vcovm1 = vcov            vcovm1 = vcov
127            ucovm1 = ucov            ucovm1 = ucov
128            tetam1 = teta            tetam1 = teta
129            massem1 = masse            massem1 = masse
130            psm1 = ps            psm1 = ps
           finvmaold = masse  
           CALL filtreg(finvmaold, jjm + 1, llm, - 2, 2, .TRUE., 1)  
131         end if         end if
132    
133         ! Calcul des tendances dynamiques:         ! Calcul des tendances dynamiques:
134         CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)         CALL geopot(teta, pk, pks, phis, phi)
135         CALL caldyn(itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &         CALL caldyn(itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
136              dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &              dudyn, dv, dteta, dp, w, pbaru, pbarv, &
137              conser=MOD(itau, iconser)==0)              conser = MOD(itau, iconser) == 0)
138    
139         ! Calcul des tendances advection des traceurs (dont l'humidité)         CALL caladvtrac(q, pbaru, pbarv, p3d, masse, teta, pk)
        CALL caladvtrac(q, pbaru, pbarv, p3d, masse, dq, teta, pk)  
140    
141         ! Stokage du flux de masse pour traceurs offline:         ! Int\'egrations dynamique et traceurs:
        IF (offline) CALL fluxstokenc(pbaru, pbarv, masse, teta, phi, phis, &  
             dtvr, itau)  
   
        ! Integrations dynamique et traceurs:  
142         CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, dteta, &         CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, dteta, &
143              dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, finvmaold, dt, &              dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, dt, leapf)
144              leapf)  
145           forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
146           CALL exner_hyb(ps, p3d, pks, pk)
147           pkf = pk
148           CALL filtreg_scal(pkf, direct = .true., intensive = .true.)
149    
150         if (.not. leapf) then         if (.not. leapf) then
151            ! Matsuno backward            ! Matsuno backward
           forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps  
           CALL exner_hyb(ps, p3d, pks, pk, pkf)  
   
152            ! Calcul des tendances dynamiques:            ! Calcul des tendances dynamiques:
153            CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)            CALL geopot(teta, pk, pks, phis, phi)
154            CALL caldyn(itau + 1, ucov, vcov, teta, ps, masse, pk, pkf, phis, &            CALL caldyn(itau + 1, ucov, vcov, teta, ps, masse, pk, pkf, phis, &
155                 phi, dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &                 phi, dudyn, dv, dteta, dp, w, pbaru, pbarv, conser = .false.)
                conser=.false.)  
156    
157            ! integrations dynamique et traceurs:            ! integrations dynamique et traceurs:
158            CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, &            CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, &
159                 dteta, dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, &                 dteta, dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, dtvr, &
160                 finvmaold, dtvr, leapf=.false.)                 leapf=.false.)
        end if  
   
        IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys /= 0) THEN  
           ! calcul des tendances physiques:  
161    
162            forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps            forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
163            CALL exner_hyb(ps, p3d, pks, pk, pkf)            CALL exner_hyb(ps, p3d, pks, pk)
164              pkf = pk
165              CALL filtreg_scal(pkf, direct = .true., intensive = .true.)
166           end if
167    
168            rdaym_ini = itau * dtvr / daysec         IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys) THEN
169            rdayvrai = rdaym_ini + day_ini            CALL calfis(ucov, vcov, teta, q, p3d, pk, phis, phi, w, dufi, dvfi, &
170            time = REAL(mod(itau, day_step)) / day_step + time_0                 dtetafi, dqfi, dayvrai = itau / day_step + day_ini, &
171            IF (time > 1.) time = time - 1.                 time = REAL(mod(itau, day_step)) / day_step, &
172                   lafin = itau + 1 == itaufin)
           CALL calfis(rdayvrai, time, ucov, vcov, teta, q, masse, ps, pk, &  
                phis, phi, dudyn, dv, dq, w, dufi, dvfi, dtetafi, dqfi, dpfi, &  
                lafin=itau+1==itaufin)  
   
           ! ajout des tendances physiques:  
           CALL addfi(nqmx, dtphys, ucov, vcov, teta, q, ps, dufi, dvfi, &  
                dtetafi, dqfi, dpfi)  
        ENDIF  
173    
174         forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps            CALL addfi(ucov, vcov, teta, q, dufi, dvfi, dtetafi, dqfi)
175         CALL exner_hyb(ps, p3d, pks, pk, pkf)         ENDIF
176    
177         IF (MOD(itau + 1, idissip) == 0) THEN         IF (MOD(itau + 1, idissip) == 0) THEN
178            ! dissipation horizontale et verticale des petites echelles:            ! Dissipation horizontale et verticale des petites \'echelles
179    
180            ! calcul de l'energie cinetique avant dissipation            ! calcul de l'\'energie cin\'etique avant dissipation
181            call covcont(llm, ucov, vcov, ucont, vcont)            call covcont(llm, ucov, vcov, ucont, vcont)
182            call enercin(vcov, ucov, vcont, ucont, ecin0)            call enercin(vcov, ucov, vcont, ucont, ecin0)
183    
184            ! dissipation            ! dissipation
185            CALL dissip(vcov, ucov, teta, p3d, dvdis, dudis, dtetadis)            CALL dissip(vcov, ucov, teta, p3d, dvdis, dudis, dtetadis)
186            ucov=ucov + dudis            ucov = ucov + dudis
187            vcov=vcov + dvdis            vcov = vcov + dvdis
188    
189            ! On rajoute la tendance due à la transformation Ec -> E            ! On ajoute la tendance due \`a la transformation \'energie
190            ! thermique créée lors de la dissipation            ! cin\'etique en \'energie thermique par la dissipation
191            call covcont(llm, ucov, vcov, ucont, vcont)            call covcont(llm, ucov, vcov, ucont, vcont)
192            call enercin(vcov, ucov, vcont, ucont, ecin)            call enercin(vcov, ucov, vcont, ucont, ecin)
193            dtetaecdt= (ecin0 - ecin) / pk            dtetadis = dtetadis + (ecin0 - ecin) / pk
194            dtetadis=dtetadis + dtetaecdt            teta = teta + dtetadis
           teta=teta + dtetadis  
195    
196            ! Calcul de la valeur moyenne aux pôles :            ! Calcul de la valeur moyenne aux p\^oles :
197            forall (l = 1: llm)            forall (l = 1: llm)
198               teta(:, 1, l) = SUM(aire_2d(:iim, 1) * teta(:iim, 1, l)) &               teta(:, 1, l) = SUM(aire_2d(:iim, 1) * teta(:iim, 1, l)) &
199                    / apoln                    / apoln
200               teta(:, jjm + 1, l) = SUM(aire_2d(:iim, jjm+1) &               teta(:, jjm + 1, l) = SUM(aire_2d(:iim, jjm + 1) &
201                    * teta(:iim, jjm + 1, l)) / apols                    * teta(:iim, jjm + 1, l)) / apols
202            END forall            END forall
   
           ps(:, 1) = SUM(aire_2d(:iim, 1) * ps(:iim, 1)) / apoln  
           ps(:, jjm + 1) = SUM(aire_2d(:iim, jjm+1) * ps(:iim, jjm + 1)) &  
                / apols  
203         END IF         END IF
204    
205         IF (MOD(itau + 1, iperiod) == 0) THEN         IF (MOD(itau + 1, iperiod) == 0) THEN
206            ! Écriture du fichier histoire moyenne:            ! \'Ecriture du fichier histoire moyenne:
207            CALL writedynav(histaveid, nqmx, itau + 1, vcov, ucov, teta, pk, &            CALL writedynav(vcov, ucov, teta, pk, phi, q, masse, ps, phis, &
208                 phi, q, masse, ps, phis)                 time = itau + 1)
209            call bilan_dyn(ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, &            call bilan_dyn(ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, &
210                 q(:, :, :, 1), dt_app = dtvr * iperiod, &                 q(:, :, :, 1))
                dt_cum = dtvr * day_step * periodav)  
211         ENDIF         ENDIF
212    
213           IF (MOD(itau + 1, iecri * day_step) == 0) THEN
214              CALL geopot(teta, pk, pks, phis, phi)
215              CALL writehist(itau, vcov, ucov, teta, phi, masse, ps)
216           END IF
217      end do time_integration      end do time_integration
218    
219      CALL dynredem1("restart.nc", vcov, ucov, teta, q, masse, ps, &      CALL dynredem1(vcov, ucov, teta, q, masse, ps, itau = itau_dyn + itaufin)
          itau=itau_dyn+itaufin)  
220    
221      ! Calcul des tendances dynamiques:      ! Calcul des tendances dynamiques:
222      CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)      CALL geopot(teta, pk, pks, phis, phi)
223      CALL caldyn(itaufin, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &      CALL caldyn(itaufin, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
224           dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &           dudyn, dv, dteta, dp, w, pbaru, pbarv, &
225           conser=MOD(itaufin, iconser)==0)           conser = MOD(itaufin, iconser) == 0)
226    
227    END SUBROUTINE leapfrog    END SUBROUTINE leapfrog
228    
229  end module leapfrog_m  end module leapfrog_m

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

  ViewVC Help
Powered by ViewVC 1.1.21