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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 128 - (show annotations)
Thu Feb 12 16:23:33 2015 UTC (9 years, 3 months ago) by guez
Original Path: trunk/dyn3d/leapfrog.f
File size: 8496 byte(s)
The variable temps of file restart.nc is always 0. So we remove the
possibility that it can be something else. So removed argument time_0
of caldyn, dynetat0, leapfrog.

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 ! 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, 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, iecri
23 USE conf_guide_m, ONLY: ok_guide
24 USE dimens_m, ONLY: iim, jjm, llm, nqmx
25 use dissip_m, only: dissip
26 USE dynetat0_m, ONLY: day_ini
27 use dynredem1_m, only: dynredem1
28 USE exner_hyb_m, ONLY: exner_hyb
29 use filtreg_m, only: filtreg
30 use fluxstokenc_m, only: fluxstokenc
31 use geopot_m, only: geopot
32 USE guide_m, ONLY: guide
33 use inidissip_m, only: idissip
34 use integrd_m, only: integrd
35 use nr_util, only: assert
36 USE pressure_var, ONLY: p3d
37 USE temps, ONLY: itau_dyn
38 use writedynav_m, only: writedynav
39 use writehist_m, only: writehist
40
41 ! Variables dynamiques:
42 REAL, intent(inout):: ucov(:, :, :) ! (iim + 1, jjm + 1, llm) vent covariant
43 REAL, intent(inout):: vcov(:, :, :) ! (iim + 1, jjm, llm) ! vent covariant
44
45 REAL, intent(inout):: teta(:, :, :) ! (iim + 1, jjm + 1, llm)
46 ! potential temperature
47
48 REAL, intent(inout):: ps(:, :) ! (iim + 1, jjm + 1) pression au sol, en Pa
49 REAL, intent(inout):: masse(:, :, :) ! (iim + 1, jjm + 1, llm) masse d'air
50 REAL, intent(in):: phis(:, :) ! (iim + 1, jjm + 1) surface geopotential
51
52 REAL, intent(inout):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx)
53 ! mass fractions of advected fields
54
55 ! Local:
56
57 ! Variables dynamiques:
58
59 REAL pks(iim + 1, jjm + 1) ! exner au sol
60 REAL pk(iim + 1, jjm + 1, llm) ! exner au milieu des couches
61 REAL pkf(iim + 1, jjm + 1, llm) ! exner filtr\'e au milieu des couches
62 REAL phi(iim + 1, jjm + 1, llm) ! geopotential
63 REAL w(iim + 1, jjm + 1, llm) ! vitesse verticale
64
65 ! Variables dynamiques intermediaire pour le transport
66 ! Flux de masse :
67 REAL pbaru(iim + 1, jjm + 1, llm), pbarv(iim + 1, jjm, llm)
68
69 ! Variables dynamiques au pas - 1
70 REAL vcovm1(iim + 1, jjm, llm), ucovm1(iim + 1, jjm + 1, llm)
71 REAL tetam1(iim + 1, jjm + 1, llm), psm1(iim + 1, jjm + 1)
72 REAL massem1(iim + 1, jjm + 1, llm)
73
74 ! Tendances dynamiques
75 REAL dv((iim + 1) * jjm, llm), dudyn(iim + 1, jjm + 1, llm)
76 REAL dteta(iim + 1, jjm + 1, llm)
77 real dp((iim + 1) * (jjm + 1))
78
79 ! Tendances de la dissipation :
80 REAL dvdis(iim + 1, jjm, llm), dudis(iim + 1, jjm + 1, llm)
81 REAL dtetadis(iim + 1, jjm + 1, llm)
82
83 ! Tendances physiques
84 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
87 ! Variables pour le fichier histoire
88 INTEGER itau ! index of the time step of the dynamics, starts at 0
89 INTEGER itaufin
90 REAL time ! time of day, as a fraction of day length
91 real finvmaold(iim + 1, jjm + 1, llm)
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, pkf)
112
113 time_integration: do itau = 0, itaufin - 1
114 leapf = mod(itau, iperiod) /= 0
115 if (leapf) then
116 dt = 2 * dtvr
117 else
118 ! Matsuno
119 dt = dtvr
120 if (ok_guide) call guide(itau, ucov, vcov, teta, q(:, :, :, 1), ps)
121 vcovm1 = vcov
122 ucovm1 = ucov
123 tetam1 = teta
124 massem1 = masse
125 psm1 = ps
126 finvmaold = masse
127 CALL filtreg(finvmaold, direct = .false., intensive = .false.)
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, finvmaold, dt, &
145 leapf)
146
147 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
148 CALL exner_hyb(ps, p3d, pks, pk, pkf)
149
150 if (.not. leapf) then
151 ! Matsuno backward
152 ! Calcul des tendances dynamiques:
153 CALL geopot(teta, pk, pks, phis, phi)
154 CALL caldyn(itau + 1, ucov, vcov, teta, ps, masse, pk, pkf, phis, &
155 phi, dudyn, dv, dteta, dp, w, pbaru, pbarv, conser = .false.)
156
157 ! integrations dynamique et traceurs:
158 CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, &
159 dteta, dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, &
160 finvmaold, dtvr, leapf=.false.)
161
162 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
163 CALL exner_hyb(ps, p3d, pks, pk, pkf)
164 end if
165
166 IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys /= 0) THEN
167 ! Calcul des tendances physiques :
168 time = REAL(mod(itau, day_step)) / day_step
169 IF (time > 1.) time = time - 1.
170 CALL calfis(itau * dtvr / daysec + day_ini, time, ucov, vcov, teta, &
171 q, pk, phis, phi, w, dufi, dvfi, dtetafi, dqfi, &
172 lafin = itau + 1 == itaufin)
173
174 CALL addfi(ucov, vcov, teta, q, dufi, dvfi, dtetafi, dqfi)
175 ENDIF
176
177 IF (MOD(itau + 1, idissip) == 0) THEN
178 ! Dissipation horizontale et verticale des petites \'echelles
179
180 ! calcul de l'\'energie cin\'etique avant dissipation
181 call covcont(llm, ucov, vcov, ucont, vcont)
182 call enercin(vcov, ucov, vcont, ucont, ecin0)
183
184 ! dissipation
185 CALL dissip(vcov, ucov, teta, p3d, dvdis, dudis, dtetadis)
186 ucov = ucov + dudis
187 vcov = vcov + dvdis
188
189 ! On ajoute la tendance due \`a la transformation \'energie
190 ! cin\'etique en \'energie thermique par la dissipation
191 call covcont(llm, ucov, vcov, ucont, vcont)
192 call enercin(vcov, ucov, vcont, ucont, ecin)
193 dtetadis = dtetadis + (ecin0 - ecin) / pk
194 teta = teta + dtetadis
195
196 ! Calcul de la valeur moyenne aux p\^oles :
197 forall (l = 1: llm)
198 teta(:, 1, l) = SUM(aire_2d(:iim, 1) * teta(:iim, 1, l)) &
199 / apoln
200 teta(:, jjm + 1, l) = SUM(aire_2d(:iim, jjm+1) &
201 * teta(:iim, jjm + 1, l)) / apols
202 END forall
203 END IF
204
205 IF (MOD(itau + 1, iperiod) == 0) THEN
206 ! \'Ecriture du fichier histoire moyenne:
207 CALL writedynav(vcov, ucov, teta, pk, phi, q, masse, ps, phis, &
208 time = itau + 1)
209 call bilan_dyn(ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, &
210 q(:, :, :, 1))
211 ENDIF
212
213 IF (MOD(itau + 1, iecri * day_step) == 0) THEN
214 CALL geopot(teta, pk, pks, phis, phi)
215 CALL writehist(itau, vcov, ucov, teta, phi, q, masse, ps)
216 END IF
217 end do time_integration
218
219 CALL dynredem1("restart.nc", vcov, ucov, teta, q, masse, ps, &
220 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