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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 55 - (hide annotations)
Mon Dec 12 13:25:01 2011 UTC (12 years, 5 months ago) by guez
Original Path: trunk/libf/dyn3d/leapfrog.f90
File size: 9353 byte(s)
-- In procedure "bilan_dyn", replaced average of "zvq" by integral of
"zvq", following a comment of Francis Codron :

Le calcul actuel donne des unités peu pratiques : transports de
chaleur en K m / s par exemple. C'est bien pour les sorties à 2
dimensions, latitude et pression, car alors le transport ne dépend pas
de l'espacement des niveaux, mieux pour comparer ou tracer en latitude
et pression. Par contre, quand on somme sur la verticale, on
préfèrerait avoir des transports d'énergie en watts, ou au moins an K
kg / s (à multiplier par "Cp" ou "L"). On doit pouvoir recalculer le
transport intégré à partir des fichiers de sortie, mais c'est embêtant
(calcul de "cv").

-- Gathered files in directory Dissipation.

1 guez 3 module leapfrog_m
2    
3     IMPLICIT NONE
4    
5     contains
6    
7 guez 23 SUBROUTINE leapfrog(ucov, vcov, teta, ps, masse, phis, q, time_0)
8 guez 3
9 guez 27 ! From dyn3d/leapfrog.F, version 1.6, 2005/04/13 08:58:34
10     ! Authors: P. Le Van, L. Fairhead, F. Hourdin
11 guez 36 ! Matsuno-leapfrog scheme.
12 guez 3
13 guez 37 use addfi_m, only: addfi
14 guez 40 use bilan_dyn_m, only: bilan_dyn
15     use caladvtrac_m, only: caladvtrac
16 guez 43 use caldyn_m, only: caldyn
17 guez 26 USE calfis_m, ONLY: calfis
18     USE com_io_dyn, ONLY: histaveid
19     USE comconst, ONLY: daysec, dtphys, dtvr
20 guez 29 USE comgeom, ONLY: aire_2d, apoln, apols
21 guez 26 USE comvert, ONLY: ap, bp
22 guez 27 USE conf_gcm_m, ONLY: day_step, iconser, iperiod, iphysiq, nday, offline, &
23     periodav
24 guez 29 USE dimens_m, ONLY: iim, jjm, llm, nqmx
25 guez 47 use dissip_m, only: dissip
26 guez 26 USE dynetat0_m, ONLY: day_ini
27 guez 27 use dynredem1_m, only: dynredem1
28 guez 26 USE exner_hyb_m, ONLY: exner_hyb
29 guez 27 use filtreg_m, only: filtreg
30 guez 43 use geopot_m, only: geopot
31 guez 26 USE guide_m, ONLY: guide
32     use inidissip_m, only: idissip
33 guez 32 use integrd_m, only: integrd
34 guez 26 USE logic, ONLY: iflag_phys, ok_guide
35 guez 55 use nr_util, only: assert
36 guez 26 USE pressure_var, ONLY: p3d
37 guez 28 USE temps, ONLY: itau_dyn
38 guez 3
39 guez 10 ! Variables dynamiques:
40 guez 55 REAL, intent(inout):: ucov(:, :, :) ! (iim + 1, jjm + 1, llm) vent covariant
41     REAL, intent(inout):: vcov(:, :, :) ! (iim + 1, jjm, llm) ! vent covariant
42 guez 43
43     REAL, intent(inout):: teta(:, :, :) ! (iim + 1, jjm + 1, llm)
44     ! potential temperature
45    
46 guez 45 REAL, intent(inout):: ps(:, :) ! (iim + 1, jjm + 1) pression au sol, en Pa
47 guez 55 REAL masse((iim + 1) * (jjm + 1), llm) ! masse d'air
48     REAL phis((iim + 1) * (jjm + 1)) ! geopotentiel au sol
49 guez 40
50     REAL, intent(inout):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx)
51     ! mass fractions of advected fields
52    
53 guez 24 REAL, intent(in):: time_0
54 guez 10
55     ! Variables local to the procedure:
56    
57     ! Variables dynamiques:
58    
59 guez 55 REAL pks((iim + 1) * (jjm + 1)) ! exner au sol
60 guez 29 REAL pk(iim + 1, jjm + 1, llm) ! exner au milieu des couches
61 guez 55 REAL pkf((iim + 1) * (jjm + 1), llm) ! exner filt.au milieu des couches
62 guez 47 REAL phi(iim + 1, jjm + 1, llm) ! geopotential
63 guez 55 REAL w((iim + 1) * (jjm + 1), llm) ! vitesse verticale
64 guez 3
65 guez 55 ! Variables dynamiques intermediaire pour le transport
66     ! Flux de masse :
67     REAL pbaru((iim + 1) * (jjm + 1), llm), pbarv((iim + 1) * jjm, llm)
68 guez 3
69     ! variables dynamiques au pas - 1
70 guez 55 REAL vcovm1(iim + 1, jjm, llm), ucovm1(iim + 1, jjm + 1, llm)
71 guez 29 REAL tetam1(iim + 1, jjm + 1, llm), psm1(iim + 1, jjm + 1)
72 guez 55 REAL massem1((iim + 1) * (jjm + 1), llm)
73 guez 3
74     ! tendances dynamiques
75 guez 55 REAL dv((iim + 1) * jjm, llm), dudyn((iim + 1) * (jjm + 1), llm)
76     REAL dteta(iim + 1, jjm + 1, llm), dq((iim + 1) * (jjm + 1), llm, nqmx)
77     real dp((iim + 1) * (jjm + 1))
78 guez 3
79     ! tendances de la dissipation
80 guez 55 REAL dvdis(iim + 1, jjm, llm), dudis(iim + 1, jjm + 1, llm)
81 guez 29 REAL dtetadis(iim + 1, jjm + 1, llm)
82 guez 3
83     ! tendances physiques
84 guez 55 REAL dvfi((iim + 1) * jjm, llm), dufi((iim + 1) * (jjm + 1), llm)
85     REAL dtetafi(iim + 1, jjm + 1, llm), dqfi((iim + 1) * (jjm + 1), llm, nqmx)
86     real dpfi((iim + 1) * (jjm + 1))
87 guez 3
88     ! variables pour le fichier histoire
89    
90 guez 22 INTEGER itau ! index of the time step of the dynamics, starts at 0
91 guez 27 INTEGER itaufin
92 guez 20 REAL time ! time of day, as a fraction of day length
93 guez 55 real finvmaold((iim + 1) * (jjm + 1), llm)
94 guez 33 INTEGER l
95 guez 3 REAL rdayvrai, rdaym_ini
96    
97 guez 24 ! Variables test conservation energie
98 guez 29 REAL ecin(iim + 1, jjm + 1, llm), ecin0(iim + 1, jjm + 1, llm)
99 guez 43
100 guez 29 REAL dtetaecdt(iim + 1, jjm + 1, llm)
101 guez 55 ! tendance de la température potentielle due à la transformation
102     ! d'énergie cinétique en énergie thermique par la dissipation
103 guez 43
104 guez 55 REAL vcont((iim + 1) * jjm, llm), ucont((iim + 1) * (jjm + 1), llm)
105 guez 33 logical leapf
106     real dt
107 guez 3
108     !---------------------------------------------------
109    
110     print *, "Call sequence information: leapfrog"
111 guez 55 call assert(shape(ucov) == (/iim + 1, jjm + 1, llm/), "leapfrog")
112 guez 3
113     itaufin = nday * day_step
114 guez 30 ! "day_step" is a multiple of "iperiod", therefore "itaufin" is one too
115    
116 guez 24 dq = 0.
117 guez 33
118 guez 3 ! On initialise la pression et la fonction d'Exner :
119 guez 37 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
120 guez 10 CALL exner_hyb(ps, p3d, pks, pk, pkf)
121 guez 3
122 guez 40 time_integration: do itau = 0, itaufin - 1
123 guez 33 leapf = mod(itau, iperiod) /= 0
124     if (leapf) then
125     dt = 2 * dtvr
126     else
127     ! Matsuno
128     dt = dtvr
129     if (ok_guide .and. (itaufin - itau - 1) * dtvr > 21600.) &
130     call guide(itau, ucov, vcov, teta, q, masse, ps)
131     vcovm1 = vcov
132     ucovm1 = ucov
133     tetam1 = teta
134     massem1 = masse
135     psm1 = ps
136     finvmaold = masse
137     CALL filtreg(finvmaold, jjm + 1, llm, - 2, 2, .TRUE., 1)
138     end if
139 guez 30
140     ! Calcul des tendances dynamiques:
141 guez 55 CALL geopot((iim + 1) * (jjm + 1), teta, pk, pks, phis, phi)
142 guez 30 CALL caldyn(itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
143 guez 47 dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &
144     conser=MOD(itau, iconser)==0)
145 guez 30
146     ! Calcul des tendances advection des traceurs (dont l'humidité)
147     CALL caladvtrac(q, pbaru, pbarv, p3d, masse, dq, teta, pk)
148 guez 33
149 guez 30 ! Stokage du flux de masse pour traceurs offline:
150     IF (offline) CALL fluxstokenc(pbaru, pbarv, masse, teta, phi, phis, &
151     dtvr, itau)
152    
153 guez 46 ! Integrations dynamique et traceurs:
154 guez 47 CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, dteta, &
155     dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, finvmaold, dt, &
156     leapf)
157 guez 30
158 guez 33 if (.not. leapf) then
159     ! Matsuno backward
160 guez 37 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
161 guez 33 CALL exner_hyb(ps, p3d, pks, pk, pkf)
162 guez 30
163 guez 27 ! Calcul des tendances dynamiques:
164 guez 55 CALL geopot((iim + 1) * (jjm + 1), teta, pk, pks, phis, phi)
165 guez 33 CALL caldyn(itau + 1, ucov, vcov, teta, ps, masse, pk, pkf, phis, &
166 guez 47 phi, dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &
167     conser=.false.)
168 guez 3
169     ! integrations dynamique et traceurs:
170 guez 47 CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, &
171     dteta, dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, &
172     finvmaold, dtvr, leapf=.false.)
173 guez 33 end if
174 guez 3
175 guez 33 IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys /= 0) THEN
176     ! calcul des tendances physiques:
177 guez 3
178 guez 37 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
179 guez 33 CALL exner_hyb(ps, p3d, pks, pk, pkf)
180 guez 3
181 guez 33 rdaym_ini = itau * dtvr / daysec
182     rdayvrai = rdaym_ini + day_ini
183     time = REAL(mod(itau, day_step)) / day_step + time_0
184     IF (time > 1.) time = time - 1.
185 guez 3
186 guez 37 CALL calfis(rdayvrai, time, ucov, vcov, teta, q, masse, ps, pk, &
187 guez 47 phis, phi, dudyn, dv, dq, w, dufi, dvfi, dtetafi, dqfi, dpfi, &
188 guez 44 lafin=itau+1==itaufin)
189 guez 3
190 guez 33 ! ajout des tendances physiques:
191     CALL addfi(nqmx, dtphys, ucov, vcov, teta, q, ps, dufi, dvfi, &
192     dtetafi, dqfi, dpfi)
193     ENDIF
194 guez 3
195 guez 37 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
196 guez 33 CALL exner_hyb(ps, p3d, pks, pk, pkf)
197 guez 3
198 guez 33 IF (MOD(itau + 1, idissip) == 0) THEN
199 guez 55 ! Dissipation horizontale et verticale des petites échelles
200 guez 3
201 guez 55 ! calcul de l'énergie cinétique avant dissipation
202 guez 33 call covcont(llm, ucov, vcov, ucont, vcont)
203     call enercin(vcov, ucov, vcont, ucont, ecin0)
204 guez 3
205 guez 33 ! dissipation
206     CALL dissip(vcov, ucov, teta, p3d, dvdis, dudis, dtetadis)
207 guez 55 ucov = ucov + dudis
208     vcov = vcov + dvdis
209 guez 3
210 guez 55 ! On ajoute la tendance due à la transformation énergie
211     ! cinétique en énergie thermique par la dissipation
212 guez 33 call covcont(llm, ucov, vcov, ucont, vcont)
213     call enercin(vcov, ucov, vcont, ucont, ecin)
214     dtetaecdt= (ecin0 - ecin) / pk
215 guez 55 dtetadis = dtetadis + dtetaecdt
216     teta = teta + dtetadis
217 guez 3
218 guez 33 ! Calcul de la valeur moyenne aux pôles :
219     forall (l = 1: llm)
220     teta(:, 1, l) = SUM(aire_2d(:iim, 1) * teta(:iim, 1, l)) &
221     / apoln
222     teta(:, jjm + 1, l) = SUM(aire_2d(:iim, jjm+1) &
223     * teta(:iim, jjm + 1, l)) / apols
224     END forall
225 guez 3
226 guez 33 ps(:, 1) = SUM(aire_2d(:iim, 1) * ps(:iim, 1)) / apoln
227     ps(:, jjm + 1) = SUM(aire_2d(:iim, jjm+1) * ps(:iim, jjm + 1)) &
228     / apols
229     END IF
230 guez 3
231 guez 33 IF (MOD(itau + 1, iperiod) == 0) THEN
232 guez 40 ! Écriture du fichier histoire moyenne:
233 guez 33 CALL writedynav(histaveid, nqmx, itau + 1, vcov, ucov, teta, pk, &
234     phi, q, masse, ps, phis)
235 guez 40 call bilan_dyn(ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, &
236     q(:, :, :, 1), dt_app = dtvr * iperiod, &
237     dt_cum = dtvr * day_step * periodav)
238 guez 33 ENDIF
239 guez 40 end do time_integration
240 guez 3
241 guez 30 CALL dynredem1("restart.nc", vcov, ucov, teta, q, masse, ps, &
242     itau=itau_dyn+itaufin)
243    
244     ! Calcul des tendances dynamiques:
245 guez 55 CALL geopot((iim + 1) * (jjm + 1), teta, pk, pks, phis, phi)
246 guez 30 CALL caldyn(itaufin, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
247 guez 47 dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &
248     conser=MOD(itaufin, iconser)==0)
249 guez 3 END SUBROUTINE leapfrog
250    
251     end module leapfrog_m

  ViewVC Help
Powered by ViewVC 1.1.21