/[lmdze]/trunk/libf/dyn3d/leapfrog.f90
ViewVC logotype

Contents of /trunk/libf/dyn3d/leapfrog.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 41 - (show annotations)
Tue Feb 22 15:09:57 2011 UTC (13 years, 2 months ago) by guez
File size: 8765 byte(s)
Moved "Test_ozonecm" out of SVN.
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
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 calfis_m, ONLY: calfis
17 USE com_io_dyn, ONLY: histaveid
18 USE comconst, ONLY: daysec, dtphys, dtvr
19 USE comgeom, ONLY: aire_2d, apoln, apols
20 USE comvert, ONLY: ap, bp
21 USE conf_gcm_m, ONLY: day_step, iconser, iperiod, iphysiq, nday, offline, &
22 periodav
23 USE dimens_m, ONLY: iim, jjm, llm, nqmx
24 USE dynetat0_m, ONLY: day_ini
25 use dynredem1_m, only: dynredem1
26 USE exner_hyb_m, ONLY: exner_hyb
27 use filtreg_m, only: filtreg
28 USE guide_m, ONLY: guide
29 use inidissip_m, only: idissip
30 use integrd_m, only: integrd
31 USE logic, ONLY: iflag_phys, ok_guide
32 USE paramet_m, ONLY: ip1jmp1
33 USE pressure_var, ONLY: p3d
34 USE temps, ONLY: itau_dyn
35
36 ! Variables dynamiques:
37 REAL, intent(inout):: ucov(ip1jmp1, llm) ! vent covariant
38 REAL, intent(inout):: vcov((iim + 1) * jjm, llm) ! vent covariant
39 REAL, intent(inout):: teta(iim + 1, jjm + 1, llm) ! potential temperature
40 REAL, intent(inout):: ps(iim + 1, jjm + 1) ! pression au sol, en Pa
41 REAL masse(ip1jmp1, llm) ! masse d'air
42 REAL phis(ip1jmp1) ! geopotentiel au sol
43
44 REAL, intent(inout):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx)
45 ! mass fractions of advected fields
46
47 REAL, intent(in):: time_0
48
49 ! Variables local to the procedure:
50
51 ! Variables dynamiques:
52
53 REAL pks(ip1jmp1) ! exner au sol
54 REAL pk(iim + 1, jjm + 1, llm) ! exner au milieu des couches
55 REAL pkf(ip1jmp1, llm) ! exner filt.au milieu des couches
56 REAL phi(ip1jmp1, llm) ! geopotential
57 REAL w(ip1jmp1, llm) ! vitesse verticale
58
59 ! variables dynamiques intermediaire pour le transport
60 REAL pbaru(ip1jmp1, llm), pbarv((iim + 1) * jjm, llm) !flux de masse
61
62 ! variables dynamiques au pas - 1
63 REAL vcovm1((iim + 1) * jjm, llm), ucovm1(ip1jmp1, llm)
64 REAL tetam1(iim + 1, jjm + 1, llm), psm1(iim + 1, jjm + 1)
65 REAL massem1(ip1jmp1, llm)
66
67 ! tendances dynamiques
68 REAL dv((iim + 1) * jjm, llm), du(ip1jmp1, llm)
69 REAL dteta(ip1jmp1, llm), dq(ip1jmp1, llm, nqmx), dp(ip1jmp1)
70
71 ! tendances de la dissipation
72 REAL dvdis((iim + 1) * jjm, llm), dudis(ip1jmp1, llm)
73 REAL dtetadis(iim + 1, jjm + 1, llm)
74
75 ! tendances physiques
76 REAL dvfi((iim + 1) * jjm, llm), dufi(ip1jmp1, llm)
77 REAL dtetafi(ip1jmp1, llm), dqfi(ip1jmp1, llm, nqmx), dpfi(ip1jmp1)
78
79 ! variables pour le fichier histoire
80
81 INTEGER itau ! index of the time step of the dynamics, starts at 0
82 INTEGER itaufin
83 REAL time ! time of day, as a fraction of day length
84 real finvmaold(ip1jmp1, llm)
85 INTEGER l
86 REAL rdayvrai, rdaym_ini
87
88 ! Variables test conservation energie
89 REAL ecin(iim + 1, jjm + 1, llm), ecin0(iim + 1, jjm + 1, llm)
90 ! Tendance de la temp. potentiel d (theta) / d t due a la
91 ! tansformation d'energie cinetique en energie thermique
92 ! cree par la dissipation
93 REAL dtetaecdt(iim + 1, jjm + 1, llm)
94 REAL vcont((iim + 1) * jjm, llm), ucont(ip1jmp1, llm)
95 logical leapf
96 real dt
97
98 !---------------------------------------------------
99
100 print *, "Call sequence information: leapfrog"
101
102 itaufin = nday * day_step
103 ! "day_step" is a multiple of "iperiod", therefore "itaufin" is one too
104
105 dq = 0.
106
107 ! On initialise la pression et la fonction d'Exner :
108 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
109 CALL exner_hyb(ps, p3d, pks, pk, pkf)
110
111 time_integration: do itau = 0, itaufin - 1
112 leapf = mod(itau, iperiod) /= 0
113 if (leapf) then
114 dt = 2 * dtvr
115 else
116 ! Matsuno
117 dt = dtvr
118 if (ok_guide .and. (itaufin - itau - 1) * dtvr > 21600.) &
119 call guide(itau, ucov, vcov, teta, q, masse, ps)
120 vcovm1 = vcov
121 ucovm1 = ucov
122 tetam1 = teta
123 massem1 = masse
124 psm1 = ps
125 finvmaold = masse
126 CALL filtreg(finvmaold, jjm + 1, llm, - 2, 2, .TRUE., 1)
127 end if
128
129 ! Calcul des tendances dynamiques:
130 CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)
131 CALL caldyn(itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
132 MOD(itau, iconser) == 0, du, dv, dteta, dp, w, pbaru, pbarv, &
133 time_0)
134
135 ! Calcul des tendances advection des traceurs (dont l'humidité)
136 CALL caladvtrac(q, pbaru, pbarv, p3d, masse, dq, 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 ! integrations dynamique et traceurs:
143 CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, du, dteta, dp, &
144 vcov, ucov, teta, q(:, :, :, :2), ps, masse, finvmaold, dt, leapf)
145
146 if (.not. leapf) then
147 ! Matsuno backward
148 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
149 CALL exner_hyb(ps, p3d, pks, pk, pkf)
150
151 ! Calcul des tendances dynamiques:
152 CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)
153 CALL caldyn(itau + 1, ucov, vcov, teta, ps, masse, pk, pkf, phis, &
154 phi, .false., du, dv, dteta, dp, w, pbaru, pbarv, time_0)
155
156 ! integrations dynamique et traceurs:
157 CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, du, dteta, &
158 dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, finvmaold, &
159 dtvr, leapf=.false.)
160 end if
161
162 IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys /= 0) THEN
163 ! calcul des tendances physiques:
164
165 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
166 CALL exner_hyb(ps, p3d, pks, pk, pkf)
167
168 rdaym_ini = itau * dtvr / daysec
169 rdayvrai = rdaym_ini + day_ini
170 time = REAL(mod(itau, day_step)) / day_step + time_0
171 IF (time > 1.) time = time - 1.
172
173 CALL calfis(rdayvrai, time, ucov, vcov, teta, q, masse, ps, pk, &
174 phis, phi, du, dv, dteta, dq, w, dufi, dvfi, dtetafi, dqfi, &
175 dpfi, lafin=itau+1==itaufin)
176
177 ! ajout des tendances physiques:
178 CALL addfi(nqmx, dtphys, ucov, vcov, teta, q, ps, dufi, dvfi, &
179 dtetafi, dqfi, dpfi)
180 ENDIF
181
182 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
183 CALL exner_hyb(ps, p3d, pks, pk, pkf)
184
185 IF (MOD(itau + 1, idissip) == 0) THEN
186 ! dissipation horizontale et verticale des petites echelles:
187
188 ! calcul de l'energie cinetique avant dissipation
189 call covcont(llm, ucov, vcov, ucont, vcont)
190 call enercin(vcov, ucov, vcont, ucont, ecin0)
191
192 ! dissipation
193 CALL dissip(vcov, ucov, teta, p3d, dvdis, dudis, dtetadis)
194 ucov=ucov + dudis
195 vcov=vcov + dvdis
196
197 ! On rajoute la tendance due à la transformation Ec -> E
198 ! thermique créée lors de la dissipation
199 call covcont(llm, ucov, vcov, ucont, vcont)
200 call enercin(vcov, ucov, vcont, ucont, ecin)
201 dtetaecdt= (ecin0 - ecin) / pk
202 dtetadis=dtetadis + dtetaecdt
203 teta=teta + dtetadis
204
205 ! Calcul de la valeur moyenne aux pôles :
206 forall (l = 1: llm)
207 teta(:, 1, l) = SUM(aire_2d(:iim, 1) * teta(:iim, 1, l)) &
208 / apoln
209 teta(:, jjm + 1, l) = SUM(aire_2d(:iim, jjm+1) &
210 * teta(:iim, jjm + 1, l)) / apols
211 END forall
212
213 ps(:, 1) = SUM(aire_2d(:iim, 1) * ps(:iim, 1)) / apoln
214 ps(:, jjm + 1) = SUM(aire_2d(:iim, jjm+1) * ps(:iim, jjm + 1)) &
215 / apols
216 END IF
217
218 IF (MOD(itau + 1, iperiod) == 0) THEN
219 ! Écriture du fichier histoire moyenne:
220 CALL writedynav(histaveid, nqmx, itau + 1, vcov, ucov, teta, pk, &
221 phi, q, masse, ps, phis)
222 call bilan_dyn(ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, &
223 q(:, :, :, 1), dt_app = dtvr * iperiod, &
224 dt_cum = dtvr * day_step * periodav)
225 ENDIF
226 end do time_integration
227
228 CALL dynredem1("restart.nc", vcov, ucov, teta, q, masse, ps, &
229 itau=itau_dyn+itaufin)
230
231 ! Calcul des tendances dynamiques:
232 CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)
233 CALL caldyn(itaufin, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
234 MOD(itaufin, iconser) == 0, du, dv, dteta, dp, w, pbaru, pbarv, &
235 time_0)
236
237 END SUBROUTINE leapfrog
238
239 end module leapfrog_m

  ViewVC Help
Powered by ViewVC 1.1.21