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

Contents of /trunk/dyn3d/leapfrog.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 161 - (show annotations)
Fri Jul 24 14:27:59 2015 UTC (8 years, 10 months ago) by guez
Original Path: trunk/Sources/dyn3d/leapfrog.f
File size: 8467 byte(s)
rlon[uv] and rlat[uv] are already in start.nc.

Just encapsulated covcont in a module.

finvmaold was not used in leapfrog. Downgraded it from dummy argument
to local variable of SUBROUTINE integrd.

Simplified handling of mass in integrd: down from five 3-dimensional
arrays (masse, massem1, finvmaold, massescr and finvmasse) to three
(masse, massem1, finvmaold).

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

  ViewVC Help
Powered by ViewVC 1.1.21