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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 90 - (show annotations)
Wed Mar 12 21:16:36 2014 UTC (10 years, 2 months ago) by guez
Original Path: trunk/dyn3d/leapfrog.f
File size: 8940 byte(s)
Removed procedures ini_histday, ini_histhf, write_histday and
write_histhf.

Divided file regr_pr_coefoz.f into regr_pr_av.f and
regr_pr_int.f. (Following LMDZ.) Divided module regr_pr_coefoz into
modules regr_pr_av_m and regr_pr_int_m. Renamed regr_pr_av_coefoz to
regr_pr_av and regr_pr_int_coefoz to regr_pr_int. The idea is that
those procedures are more general than Mobidic.

Removed argument dudyn of calfis and physiq. dudyn is not used either
in LMDZ. Removed computation in calfis of unused variable zpsrf (not
used either in LMDZ). Removed useless computation of dqfi in calfis
(part 62): the results were overwritten. (Same in LMDZ.)

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 ! Variables local to the procedure:
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 real dpfi((iim + 1) * (jjm + 1))
88
89 ! Variables pour le fichier histoire
90
91 INTEGER itau ! index of the time step of the dynamics, starts at 0
92 INTEGER itaufin
93 REAL time ! time of day, as a fraction of day length
94 real finvmaold(iim + 1, jjm + 1, llm)
95 INTEGER l
96 REAL rdayvrai, rdaym_ini
97
98 ! Variables test conservation \'energie
99 REAL ecin(iim + 1, jjm + 1, llm), ecin0(iim + 1, jjm + 1, llm)
100
101 REAL vcont((iim + 1) * jjm, llm), ucont((iim + 1) * (jjm + 1), llm)
102 logical leapf
103 real dt
104
105 !---------------------------------------------------
106
107 print *, "Call sequence information: leapfrog"
108 call assert(shape(ucov) == (/iim + 1, jjm + 1, llm/), "leapfrog")
109
110 itaufin = nday * day_step
111 ! "day_step" is a multiple of "iperiod", therefore so is "itaufin".
112
113 ! On initialise la pression et la fonction d'Exner :
114 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
115 CALL exner_hyb(ps, p3d, pks, pk, pkf)
116
117 time_integration: do itau = 0, itaufin - 1
118 leapf = mod(itau, iperiod) /= 0
119 if (leapf) then
120 dt = 2 * dtvr
121 else
122 ! Matsuno
123 dt = dtvr
124 if (ok_guide .and. (itaufin - itau - 1) * dtvr > 21600.) &
125 call guide(itau, ucov, vcov, teta, q, masse, ps)
126 vcovm1 = vcov
127 ucovm1 = ucov
128 tetam1 = teta
129 massem1 = masse
130 psm1 = ps
131 finvmaold = masse
132 CALL filtreg(finvmaold, jjm + 1, llm, - 2, 2, .TRUE.)
133 end if
134
135 ! Calcul des tendances dynamiques:
136 CALL geopot(teta, pk, pks, phis, phi)
137 CALL caldyn(itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
138 dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &
139 conser = MOD(itau, iconser) == 0)
140
141 CALL caladvtrac(q, pbaru, pbarv, p3d, masse, teta, pk)
142
143 ! Stokage du flux de masse pour traceurs offline:
144 IF (offline) CALL fluxstokenc(pbaru, pbarv, masse, teta, phi, phis, &
145 dtvr, itau)
146
147 ! Int\'egrations dynamique et traceurs:
148 CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, dteta, &
149 dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, finvmaold, dt, &
150 leapf)
151
152 if (.not. leapf) then
153 ! Matsuno backward
154 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
155 CALL exner_hyb(ps, p3d, pks, pk, pkf)
156
157 ! Calcul des tendances dynamiques:
158 CALL geopot(teta, pk, pks, phis, phi)
159 CALL caldyn(itau + 1, ucov, vcov, teta, ps, masse, pk, pkf, phis, &
160 phi, dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &
161 conser = .false.)
162
163 ! integrations dynamique et traceurs:
164 CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, &
165 dteta, dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, &
166 finvmaold, dtvr, leapf=.false.)
167 end if
168
169 IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys /= 0) THEN
170 ! Calcul des tendances physiques:
171
172 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
173 CALL exner_hyb(ps, p3d, pks, pk, pkf)
174
175 rdaym_ini = itau * dtvr / daysec
176 rdayvrai = rdaym_ini + day_ini
177 time = REAL(mod(itau, day_step)) / day_step + time_0
178 IF (time > 1.) time = time - 1.
179
180 CALL calfis(rdayvrai, time, ucov, vcov, teta, q, ps, pk, phis, phi, &
181 w, dufi, dvfi, dtetafi, dqfi, dpfi, lafin = itau + 1 == itaufin)
182
183 ! Ajout des tendances physiques:
184 CALL addfi(ucov, vcov, teta, q, ps, dufi, dvfi, dtetafi, dqfi, dpfi)
185 ENDIF
186
187 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
188 CALL exner_hyb(ps, p3d, pks, pk, pkf)
189
190 IF (MOD(itau + 1, idissip) == 0) THEN
191 ! Dissipation horizontale et verticale des petites \'echelles
192
193 ! calcul de l'\'energie cin\'etique avant dissipation
194 call covcont(llm, ucov, vcov, ucont, vcont)
195 call enercin(vcov, ucov, vcont, ucont, ecin0)
196
197 ! dissipation
198 CALL dissip(vcov, ucov, teta, p3d, dvdis, dudis, dtetadis)
199 ucov = ucov + dudis
200 vcov = vcov + dvdis
201
202 ! On ajoute la tendance due \`a la transformation \'energie
203 ! cin\'etique en \'energie thermique par la dissipation
204 call covcont(llm, ucov, vcov, ucont, vcont)
205 call enercin(vcov, ucov, vcont, ucont, ecin)
206 dtetadis = dtetadis + (ecin0 - ecin) / pk
207 teta = teta + dtetadis
208
209 ! Calcul de la valeur moyenne aux p\^oles :
210 forall (l = 1: llm)
211 teta(:, 1, l) = SUM(aire_2d(:iim, 1) * teta(:iim, 1, l)) &
212 / apoln
213 teta(:, jjm + 1, l) = SUM(aire_2d(:iim, jjm+1) &
214 * teta(:iim, jjm + 1, l)) / apols
215 END forall
216 END IF
217
218 IF (MOD(itau + 1, iperiod) == 0) THEN
219 ! \'Ecriture du fichier histoire moyenne:
220 CALL writedynav(vcov, ucov, teta, pk, phi, q, masse, ps, phis, &
221 time = itau + 1)
222 call bilan_dyn(ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, &
223 q(:, :, :, 1))
224 ENDIF
225
226 IF (MOD(itau + 1, iecri * day_step) == 0) THEN
227 CALL geopot(teta, pk, pks, phis, phi)
228 CALL writehist(itau, vcov, ucov, teta, phi, q, masse, ps)
229 END IF
230 end do time_integration
231
232 CALL dynredem1("restart.nc", vcov, ucov, teta, q, masse, ps, &
233 itau = itau_dyn + itaufin)
234
235 ! Calcul des tendances dynamiques:
236 CALL geopot(teta, pk, pks, phis, phi)
237 CALL caldyn(itaufin, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
238 dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &
239 conser = MOD(itaufin, iconser) == 0)
240
241 END SUBROUTINE leapfrog
242
243 end module leapfrog_m

  ViewVC Help
Powered by ViewVC 1.1.21