/[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 40 by guez, Tue Feb 22 13:49:36 2011 UTC trunk/Sources/dyn3d/leapfrog.f revision 161 by guez, Fri Jul 24 14:27:59 2015 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
18      USE calfis_m, ONLY: calfis      USE calfis_m, ONLY: calfis
19      USE com_io_dyn, ONLY: histaveid      USE comconst, ONLY: daysec, 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 disvert_m, ONLY: ap, bp
23      USE conf_gcm_m, ONLY: day_step, iconser, iperiod, iphysiq, nday, offline, &      USE conf_gcm_m, ONLY: day_step, iconser, iperiod, iphysiq, nday, offline, &
24           periodav           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
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 exner_hyb_m, ONLY: exner_hyb      USE exner_hyb_m, ONLY: exner_hyb
31      use filtreg_m, only: filtreg      use filtreg_scal_m, only: filtreg_scal
32        use fluxstokenc_m, only: fluxstokenc
33        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  
38      USE pressure_var, ONLY: p3d      USE pressure_var, ONLY: p3d
39      USE temps, ONLY: itau_dyn      USE temps, ONLY: itau_dyn
40        use writedynav_m, only: writedynav
41        use writehist_m, only: writehist
42    
43      ! Variables dynamiques:      ! Variables dynamiques:
44      REAL, intent(inout):: ucov(ip1jmp1, llm) ! vent covariant      REAL, intent(inout):: ucov(:, :, :) ! (iim + 1, jjm + 1, llm) vent covariant
45      REAL, intent(inout):: vcov((iim + 1) * jjm, llm) ! vent covariant      REAL, intent(inout):: vcov(:, :, :) ! (iim + 1, jjm, llm) ! vent covariant
46      REAL, intent(inout):: teta(iim + 1, jjm + 1, llm) ! potential temperature  
47      REAL ps(iim + 1, jjm + 1) ! pression au sol, en Pa      REAL, intent(inout):: teta(:, :, :) ! (iim + 1, jjm + 1, llm)
48      REAL masse(ip1jmp1, llm) ! masse d'air      ! potential temperature
49      REAL phis(ip1jmp1) ! geopotentiel au sol  
50        REAL, intent(inout):: ps(:, :) ! (iim + 1, jjm + 1) pression au sol, en Pa
51        REAL, intent(inout):: masse(:, :, :) ! (iim + 1, jjm + 1, llm) masse d'air
52        REAL, intent(in):: phis(:, :) ! (iim + 1, jjm + 1) surface geopotential
53    
54      REAL, intent(inout):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx)      REAL, intent(inout):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx)
55      ! mass fractions of advected fields      ! mass fractions of advected fields
56    
57      REAL, intent(in):: time_0      ! Local:
   
     ! Variables local to the procedure:  
58    
59      ! Variables dynamiques:      ! Variables dynamiques:
60    
61      REAL pks(ip1jmp1) ! exner au sol      REAL pks(iim + 1, jjm + 1) ! exner au sol
62      REAL pk(iim + 1, jjm + 1, llm) ! exner au milieu des couches      REAL pk(iim + 1, jjm + 1, llm) ! exner au milieu des couches
63      REAL pkf(ip1jmp1, llm) ! exner filt.au milieu des couches      REAL pkf(iim + 1, jjm + 1, llm) ! exner filtr\'e au milieu des couches
64      REAL phi(ip1jmp1, llm) ! geopotential      REAL phi(iim + 1, jjm + 1, llm) ! geopotential
65      REAL w(ip1jmp1, llm) ! vitesse verticale      REAL w(iim + 1, jjm + 1, llm) ! vitesse verticale
66    
67      ! variables dynamiques intermediaire pour le transport      ! Variables dynamiques intermediaire pour le transport
68      REAL pbaru(ip1jmp1, llm), pbarv((iim + 1) * jjm, llm) !flux de masse      ! Flux de masse :
69        REAL pbaru(iim + 1, jjm + 1, llm), pbarv(iim + 1, jjm, llm)
70    
71      ! variables dynamiques au pas - 1      ! Variables dynamiques au pas - 1
72      REAL vcovm1((iim + 1) * jjm, llm), ucovm1(ip1jmp1, llm)      REAL vcovm1(iim + 1, jjm, llm), ucovm1(iim + 1, jjm + 1, llm)
73      REAL tetam1(iim + 1, jjm + 1, llm), psm1(iim + 1, jjm + 1)      REAL tetam1(iim + 1, jjm + 1, llm), psm1(iim + 1, jjm + 1)
74      REAL massem1(ip1jmp1, llm)      REAL massem1(iim + 1, jjm + 1, llm)
75    
76      ! tendances dynamiques      ! Tendances dynamiques
77      REAL dv((iim + 1) * jjm, llm), du(ip1jmp1, llm)      REAL dv((iim + 1) * jjm, llm), dudyn(iim + 1, jjm + 1, llm)
78      REAL dteta(ip1jmp1, llm), dq(ip1jmp1, llm, nqmx), dp(ip1jmp1)      REAL dteta(iim + 1, jjm + 1, llm)
79        real dp((iim + 1) * (jjm + 1))
80    
81      ! tendances de la dissipation      ! Tendances de la dissipation :
82      REAL dvdis((iim + 1) * jjm, llm), dudis(ip1jmp1, llm)      REAL dvdis(iim + 1, jjm, llm), dudis(iim + 1, jjm + 1, llm)
83      REAL dtetadis(iim + 1, jjm + 1, llm)      REAL dtetadis(iim + 1, jjm + 1, llm)
84    
85      ! tendances physiques      ! Tendances physiques
86      REAL dvfi((iim + 1) * jjm, llm), dufi(ip1jmp1, llm)      REAL dvfi(iim + 1, jjm, llm), dufi(iim + 1, jjm + 1, llm)
87      REAL dtetafi(ip1jmp1, 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  
88    
89        ! Variables pour le fichier histoire
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 itaufin      INTEGER itaufin
     REAL time ! time of day, as a fraction of day length  
     real finvmaold(ip1jmp1, llm)  
92      INTEGER l      INTEGER l
     REAL rdayvrai, rdaym_ini  
93    
94      ! Variables test conservation energie      ! Variables test conservation \'energie
95      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)
96      ! Tendance de la temp. potentiel d (theta) / d t due a la  
97      ! tansformation d'energie cinetique en energie thermique      REAL vcont((iim + 1) * jjm, llm), ucont((iim + 1) * (jjm + 1), llm)
     ! cree par la dissipation  
     REAL dtetaecdt(iim + 1, jjm + 1, llm)  
     REAL vcont((iim + 1) * jjm, llm), ucont(ip1jmp1, llm)  
98      logical leapf      logical leapf
99      real dt      real dt ! time step, in s
100    
101      !---------------------------------------------------      !---------------------------------------------------
102    
103      print *, "Call sequence information: leapfrog"      print *, "Call sequence information: leapfrog"
104        call assert(shape(ucov) == (/iim + 1, jjm + 1, llm/), "leapfrog")
105    
106      itaufin = nday * day_step      itaufin = nday * day_step
107      ! "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.  
108    
109      ! On initialise la pression et la fonction d'Exner :      ! On initialise la pression et la fonction d'Exner :
110      forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps      forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
111      CALL exner_hyb(ps, p3d, pks, pk, pkf)      CALL exner_hyb(ps, p3d, pks, pk)
112        pkf = pk
113        CALL filtreg_scal(pkf, direct = .true., intensive = .true.)
114    
115      time_integration: do itau = 0, itaufin - 1      time_integration: do itau = 0, itaufin - 1
116         leapf = mod(itau, iperiod) /= 0         leapf = mod(itau, iperiod) /= 0
# Line 115  contains Line 119  contains
119         else         else
120            ! Matsuno            ! Matsuno
121            dt = dtvr            dt = dtvr
122            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)  
123            vcovm1 = vcov            vcovm1 = vcov
124            ucovm1 = ucov            ucovm1 = ucov
125            tetam1 = teta            tetam1 = teta
126            massem1 = masse            massem1 = masse
127            psm1 = ps            psm1 = ps
           finvmaold = masse  
           CALL filtreg(finvmaold, jjm + 1, llm, - 2, 2, .TRUE., 1)  
128         end if         end if
129    
130         ! Calcul des tendances dynamiques:         ! Calcul des tendances dynamiques:
131         CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)         CALL geopot(teta, pk, pks, phis, phi)
132         CALL caldyn(itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &         CALL caldyn(itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
133              MOD(itau, iconser) == 0, du, dv, dteta, dp, w, pbaru, pbarv, &              dudyn, dv, dteta, dp, w, pbaru, pbarv, &
134              time_0)              conser = MOD(itau, iconser) == 0)
135    
136         ! 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)  
137    
138         ! Stokage du flux de masse pour traceurs offline:         ! Stokage du flux de masse pour traceurs offline:
139         IF (offline) CALL fluxstokenc(pbaru, pbarv, masse, teta, phi, phis, &         IF (offline) CALL fluxstokenc(pbaru, pbarv, masse, teta, phi, phis, &
140              dtvr, itau)              dtvr, itau)
141    
142         ! integrations dynamique et traceurs:         ! Int\'egrations dynamique et traceurs:
143         CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, du, dteta, dp, &         CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, dteta, &
144              vcov, ucov, teta, q(:, :, :, :2), ps, masse, finvmaold, dt, leapf)              dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, dt, leapf)
145    
146           forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
147           CALL exner_hyb(ps, p3d, pks, pk)
148           pkf = pk
149           CALL filtreg_scal(pkf, direct = .true., intensive = .true.)
150    
151         if (.not. leapf) then         if (.not. leapf) then
152            ! Matsuno backward            ! Matsuno backward
           forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps  
           CALL exner_hyb(ps, p3d, pks, pk, pkf)  
   
153            ! Calcul des tendances dynamiques:            ! Calcul des tendances dynamiques:
154            CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)            CALL geopot(teta, pk, pks, phis, phi)
155            CALL caldyn(itau + 1, ucov, vcov, teta, ps, masse, pk, pkf, phis, &            CALL caldyn(itau + 1, ucov, vcov, teta, ps, masse, pk, pkf, phis, &
156                 phi, .false., du, dv, dteta, dp, w, pbaru, pbarv, time_0)                 phi, dudyn, dv, dteta, dp, w, pbaru, pbarv, conser = .false.)
157    
158            ! integrations dynamique et traceurs:            ! integrations dynamique et traceurs:
159            CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, du, dteta, &            CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, &
160                 dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, finvmaold, &                 dteta, dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, dtvr, &
161                 dtvr, leapf=.false.)                 leapf=.false.)
162    
163              forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
164              CALL exner_hyb(ps, p3d, pks, pk)
165              pkf = pk
166              CALL filtreg_scal(pkf, direct = .true., intensive = .true.)
167         end if         end if
168    
169         IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys /= 0) THEN         IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys /= 0) THEN
170            ! calcul des tendances physiques:            CALL calfis(ucov, vcov, teta, q, pk, phis, phi, w, dufi, dvfi, &
171                   dtetafi, dqfi, dayvrai = itau / day_step + day_ini, &
172            forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps                 time = REAL(mod(itau, day_step)) / day_step, &
173            CALL exner_hyb(ps, p3d, pks, pk, pkf)                 lafin = itau + 1 == itaufin)
174    
175            rdaym_ini = itau * dtvr / daysec            CALL addfi(ucov, vcov, teta, q, dufi, dvfi, dtetafi, dqfi)
           rdayvrai = rdaym_ini + day_ini  
           time = REAL(mod(itau, day_step)) / day_step + time_0  
           IF (time > 1.) time = time - 1.  
   
           CALL calfis(rdayvrai, time, ucov, vcov, teta, q, masse, ps, pk, &  
                phis, phi, du, dv, dteta, 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)  
176         ENDIF         ENDIF
177    
        forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps  
        CALL exner_hyb(ps, p3d, pks, pk, pkf)  
   
178         IF (MOD(itau + 1, idissip) == 0) THEN         IF (MOD(itau + 1, idissip) == 0) THEN
179            ! dissipation horizontale et verticale des petites echelles:            ! Dissipation horizontale et verticale des petites \'echelles
180    
181            ! calcul de l'energie cinetique avant dissipation            ! calcul de l'\'energie cin\'etique avant dissipation
182            call covcont(llm, ucov, vcov, ucont, vcont)            call covcont(llm, ucov, vcov, ucont, vcont)
183            call enercin(vcov, ucov, vcont, ucont, ecin0)            call enercin(vcov, ucov, vcont, ucont, ecin0)
184    
185            ! dissipation            ! dissipation
186            CALL dissip(vcov, ucov, teta, p3d, dvdis, dudis, dtetadis)            CALL dissip(vcov, ucov, teta, p3d, dvdis, dudis, dtetadis)
187            ucov=ucov + dudis            ucov = ucov + dudis
188            vcov=vcov + dvdis            vcov = vcov + dvdis
189    
190            ! On rajoute la tendance due à la transformation Ec -> E            ! On ajoute la tendance due \`a la transformation \'energie
191            ! thermique créée lors de la dissipation            ! cin\'etique en \'energie thermique par la dissipation
192            call covcont(llm, ucov, vcov, ucont, vcont)            call covcont(llm, ucov, vcov, ucont, vcont)
193            call enercin(vcov, ucov, vcont, ucont, ecin)            call enercin(vcov, ucov, vcont, ucont, ecin)
194            dtetaecdt= (ecin0 - ecin) / pk            dtetadis = dtetadis + (ecin0 - ecin) / pk
195            dtetadis=dtetadis + dtetaecdt            teta = teta + dtetadis
           teta=teta + dtetadis  
196    
197            ! Calcul de la valeur moyenne aux pôles :            ! Calcul de la valeur moyenne aux p\^oles :
198            forall (l = 1: llm)            forall (l = 1: llm)
199               teta(:, 1, l) = SUM(aire_2d(:iim, 1) * teta(:iim, 1, l)) &               teta(:, 1, l) = SUM(aire_2d(:iim, 1) * teta(:iim, 1, l)) &
200                    / apoln                    / apoln
201               teta(:, jjm + 1, l) = SUM(aire_2d(:iim, jjm+1) &               teta(:, jjm + 1, l) = SUM(aire_2d(:iim, jjm+1) &
202                    * teta(:iim, jjm + 1, l)) / apols                    * teta(:iim, jjm + 1, l)) / apols
203            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  
204         END IF         END IF
205    
206         IF (MOD(itau + 1, iperiod) == 0) THEN         IF (MOD(itau + 1, iperiod) == 0) THEN
207            ! Écriture du fichier histoire moyenne:            ! \'Ecriture du fichier histoire moyenne:
208            CALL writedynav(histaveid, nqmx, itau + 1, vcov, ucov, teta, pk, &            CALL writedynav(vcov, ucov, teta, pk, phi, q, masse, ps, phis, &
209                 phi, q, masse, ps, phis)                 time = itau + 1)
210            call bilan_dyn(ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, &            call bilan_dyn(ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, &
211                 q(:, :, :, 1), dt_app = dtvr * iperiod, &                 q(:, :, :, 1))
                dt_cum = dtvr * day_step * periodav)  
212         ENDIF         ENDIF
213    
214           IF (MOD(itau + 1, iecri * day_step) == 0) THEN
215              CALL geopot(teta, pk, pks, phis, phi)
216              CALL writehist(itau, vcov, ucov, teta, phi, masse, ps)
217           END IF
218      end do time_integration      end do time_integration
219    
220      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)  
221    
222      ! Calcul des tendances dynamiques:      ! Calcul des tendances dynamiques:
223      CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)      CALL geopot(teta, pk, pks, phis, phi)
224      CALL caldyn(itaufin, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &      CALL caldyn(itaufin, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
225           MOD(itaufin, iconser) == 0, du, dv, dteta, dp, w, pbaru, pbarv, &           dudyn, dv, dteta, dp, w, pbaru, pbarv, &
226           time_0)           conser = MOD(itaufin, iconser) == 0)
227    
228    END SUBROUTINE leapfrog    END SUBROUTINE leapfrog
229    

Legend:
Removed from v.40  
changed lines
  Added in v.161

  ViewVC Help
Powered by ViewVC 1.1.21