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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 47 - (show annotations)
Fri Jul 1 15:00:48 2011 UTC (12 years, 10 months ago) by guez
File size: 8959 byte(s)
Split "thermcell.f" and "cv3_routines.f".
Removed copies of files that are now in "L_util".
Moved "mva9" and "diagetpq" to their own files.
Unified variable names across procedures.

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

  ViewVC Help
Powered by ViewVC 1.1.21