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

Contents of /trunk/dyn3d/leapfrog.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 208 - (show annotations)
Wed Dec 7 16:44:53 2016 UTC (7 years, 5 months ago) by guez
Original Path: trunk/Sources/dyn3d/leapfrog.f
File size: 8633 byte(s)
Module academic was not used.

Useful values for iflag_phys were only 0 and 1 so changed type to logical.

Definition of fmagic was duplicated in procedures alboc and alboc_cd
so moved it up to interfsurf_hq and also moved multiplication by
fmagic (following LMDZ).

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: 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 enercin_m, only: enercin
31 USE exner_hyb_m, ONLY: exner_hyb
32 use filtreg_scal_m, only: filtreg_scal
33 use fluxstokenc_m, only: fluxstokenc
34 use geopot_m, only: geopot
35 USE guide_m, ONLY: guide
36 use inidissip_m, only: idissip
37 use integrd_m, only: integrd
38 use nr_util, only: assert
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 REAL p3d(iim + 1, jjm + 1, llm+1) ! pressure at layer interfaces, in Pa
102 ! ("p3d(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)",
103 ! for interface "l")
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)
116 pkf = pk
117 CALL filtreg_scal(pkf, direct = .true., intensive = .true.)
118
119 time_integration: do itau = 0, itaufin - 1
120 leapf = mod(itau, iperiod) /= 0
121 if (leapf) then
122 dt = 2 * dtvr
123 else
124 ! Matsuno
125 dt = dtvr
126 if (ok_guide) call guide(itau, ucov, vcov, teta, q(:, :, :, 1), ps)
127 vcovm1 = vcov
128 ucovm1 = ucov
129 tetam1 = teta
130 massem1 = masse
131 psm1 = ps
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, &
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, dt, leapf)
149
150 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
151 CALL exner_hyb(ps, p3d, pks, pk)
152 pkf = pk
153 CALL filtreg_scal(pkf, direct = .true., intensive = .true.)
154
155 if (.not. leapf) then
156 ! Matsuno backward
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, 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, dtvr, &
165 leapf=.false.)
166
167 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
168 CALL exner_hyb(ps, p3d, pks, pk)
169 pkf = pk
170 CALL filtreg_scal(pkf, direct = .true., intensive = .true.)
171 end if
172
173 IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys) THEN
174 CALL calfis(ucov, vcov, teta, q, p3d, pk, phis, phi, w, dufi, dvfi, &
175 dtetafi, dqfi, dayvrai = itau / day_step + day_ini, &
176 time = REAL(mod(itau, day_step)) / day_step, &
177 lafin = itau + 1 == itaufin)
178
179 CALL addfi(ucov, vcov, teta, q, dufi, dvfi, dtetafi, dqfi)
180 ENDIF
181
182 IF (MOD(itau + 1, idissip) == 0) THEN
183 ! Dissipation horizontale et verticale des petites \'echelles
184
185 ! calcul de l'\'energie cin\'etique avant dissipation
186 call covcont(llm, ucov, vcov, ucont, vcont)
187 call enercin(vcov, ucov, vcont, ucont, ecin0)
188
189 ! dissipation
190 CALL dissip(vcov, ucov, teta, p3d, dvdis, dudis, dtetadis)
191 ucov = ucov + dudis
192 vcov = vcov + dvdis
193
194 ! On ajoute la tendance due \`a la transformation \'energie
195 ! cin\'etique en \'energie thermique par la dissipation
196 call covcont(llm, ucov, vcov, ucont, vcont)
197 call enercin(vcov, ucov, vcont, ucont, ecin)
198 dtetadis = dtetadis + (ecin0 - ecin) / pk
199 teta = teta + dtetadis
200
201 ! Calcul de la valeur moyenne aux p\^oles :
202 forall (l = 1: llm)
203 teta(:, 1, l) = SUM(aire_2d(:iim, 1) * teta(:iim, 1, l)) &
204 / apoln
205 teta(:, jjm + 1, l) = SUM(aire_2d(:iim, jjm+1) &
206 * teta(:iim, jjm + 1, l)) / apols
207 END forall
208 END IF
209
210 IF (MOD(itau + 1, iperiod) == 0) THEN
211 ! \'Ecriture du fichier histoire moyenne:
212 CALL writedynav(vcov, ucov, teta, pk, phi, q, masse, ps, phis, &
213 time = itau + 1)
214 call bilan_dyn(ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, &
215 q(:, :, :, 1))
216 ENDIF
217
218 IF (MOD(itau + 1, iecri * day_step) == 0) THEN
219 CALL geopot(teta, pk, pks, phis, phi)
220 CALL writehist(itau, vcov, ucov, teta, phi, masse, ps)
221 END IF
222 end do time_integration
223
224 CALL dynredem1(vcov, ucov, teta, q, masse, ps, itau = itau_dyn + itaufin)
225
226 ! Calcul des tendances dynamiques:
227 CALL geopot(teta, pk, pks, phis, phi)
228 CALL caldyn(itaufin, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
229 dudyn, dv, dteta, dp, w, pbaru, pbarv, &
230 conser = MOD(itaufin, iconser) == 0)
231
232 END SUBROUTINE leapfrog
233
234 end module leapfrog_m

  ViewVC Help
Powered by ViewVC 1.1.21