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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 91 - (show annotations)
Wed Mar 26 17:18:58 2014 UTC (10 years, 1 month ago) by guez
Original Path: trunk/dyn3d/leapfrog.f
File size: 8852 byte(s)
Removed unused variables lock_startdate and time_stamp of module
calendar.

Noticed that physiq does not change the surface pressure. So removed
arguments ps and dpfi of subroutine addfi. dpfi was always 0. The
computation of ps in addfi included some averaging at the poles. In
principle, this does not change ps but in practice it does because of
finite numerical precision. So the results of the simulation are
changed. Removed arguments ps and dpfi of calfis. Removed argument
d_ps of physiq.

du at the poles is not computed by dudv1, so declare only the
corresponding latitudes in dudv1. caldyn passes only a section of the
array dudyn as argument.

Removed variable niadv of module iniadvtrac_m.

Declared arguments of exner_hyb as assumed-shape arrays and made all
other horizontal sizes in exner_hyb dynamic. This allows the external
program test_disvert to use exner_hyb at a single horizontal position.

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

  ViewVC Help
Powered by ViewVC 1.1.21