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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 69 - (show annotations)
Mon Feb 18 16:33:12 2013 UTC (11 years, 2 months ago) by guez
Original Path: trunk/libf/dyn3d/leapfrog.f90
File size: 9355 byte(s)
Deleted files cvparam3.f90 and nuagecom.f90. Moved variables from
module cvparam3 to module cv3_param_m. Moved variables rad_chau1 and
rad_chau2 from module nuagecom to module conf_phys_m.

Read clesphys2_nml from conf_phys instead of gcm.

Removed argument iflag_con from several procedures. Access module
variable instead.

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 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é 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), dq((iim + 1) * (jjm + 1), llm, nqmx)
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 dq = 0.
114
115 ! On initialise la pression et la fonction d'Exner :
116 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
117 CALL exner_hyb(ps, p3d, pks, pk, pkf)
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 .and. (itaufin - itau - 1) * dtvr > 21600.) &
127 call guide(itau, ucov, vcov, teta, q, masse, ps)
128 vcovm1 = vcov
129 ucovm1 = ucov
130 tetam1 = teta
131 massem1 = masse
132 psm1 = ps
133 finvmaold = masse
134 CALL filtreg(finvmaold, jjm + 1, llm, - 2, 2, .TRUE.)
135 end if
136
137 ! Calcul des tendances dynamiques:
138 CALL geopot((iim + 1) * (jjm + 1), teta, pk, pks, phis, phi)
139 CALL caldyn(itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
140 dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &
141 conser=MOD(itau, iconser)==0)
142
143 ! Calcul des tendances advection des traceurs (dont l'humidité)
144 CALL caladvtrac(q, pbaru, pbarv, p3d, masse, dq, teta, pk)
145
146 ! Stokage du flux de masse pour traceurs offline:
147 IF (offline) CALL fluxstokenc(pbaru, pbarv, masse, teta, phi, phis, &
148 dtvr, itau)
149
150 ! Integrations dynamique et traceurs:
151 CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, dteta, &
152 dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, finvmaold, dt, &
153 leapf)
154
155 if (.not. leapf) then
156 ! Matsuno backward
157 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
158 CALL exner_hyb(ps, p3d, pks, pk, pkf)
159
160 ! Calcul des tendances dynamiques:
161 CALL geopot((iim + 1) * (jjm + 1), teta, pk, pks, phis, phi)
162 CALL caldyn(itau + 1, ucov, vcov, teta, ps, masse, pk, pkf, phis, &
163 phi, dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &
164 conser=.false.)
165
166 ! integrations dynamique et traceurs:
167 CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, &
168 dteta, dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, &
169 finvmaold, dtvr, leapf=.false.)
170 end if
171
172 IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys /= 0) THEN
173 ! Calcul des tendances physiques:
174
175 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
176 CALL exner_hyb(ps, p3d, pks, pk, pkf)
177
178 rdaym_ini = itau * dtvr / daysec
179 rdayvrai = rdaym_ini + day_ini
180 time = REAL(mod(itau, day_step)) / day_step + time_0
181 IF (time > 1.) time = time - 1.
182
183 CALL calfis(rdayvrai, time, ucov, vcov, teta, q, masse, ps, pk, &
184 phis, phi, dudyn, dv, dq, w, dufi, dvfi, dtetafi, dqfi, dpfi, &
185 lafin = itau + 1 == itaufin)
186
187 ! Ajout des tendances physiques:
188 CALL addfi(nqmx, dtphys, ucov, vcov, teta, q, ps, dufi, dvfi, &
189 dtetafi, dqfi, dpfi)
190 ENDIF
191
192 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
193 CALL exner_hyb(ps, p3d, pks, pk, pkf)
194
195 IF (MOD(itau + 1, idissip) == 0) THEN
196 ! Dissipation horizontale et verticale des petites échelles
197
198 ! calcul de l'énergie cinétique avant dissipation
199 call covcont(llm, ucov, vcov, ucont, vcont)
200 call enercin(vcov, ucov, vcont, ucont, ecin0)
201
202 ! dissipation
203 CALL dissip(vcov, ucov, teta, p3d, dvdis, dudis, dtetadis)
204 ucov = ucov + dudis
205 vcov = vcov + dvdis
206
207 ! On ajoute la tendance due à la transformation énergie
208 ! cinétique en énergie thermique par la dissipation
209 call covcont(llm, ucov, vcov, ucont, vcont)
210 call enercin(vcov, ucov, vcont, ucont, ecin)
211 dtetadis = dtetadis + (ecin0 - ecin) / pk
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(vcov, ucov, teta, pk, phi, q, masse, ps, phis, &
230 time = itau + 1)
231 call bilan_dyn(ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, &
232 q(:, :, :, 1))
233 ENDIF
234
235 IF (MOD(itau + 1, iecri * day_step) == 0) THEN
236 CALL geopot((iim + 1) * (jjm + 1), teta, pk, pks, phis, phi)
237 CALL writehist(itau, vcov, ucov, teta, phi, q, masse, ps)
238 END IF
239 end do time_integration
240
241 CALL dynredem1("restart.nc", vcov, ucov, teta, q, masse, ps, &
242 itau = itau_dyn + itaufin)
243
244 ! Calcul des tendances dynamiques:
245 CALL geopot((iim + 1) * (jjm + 1), teta, pk, pks, phis, phi)
246 CALL caldyn(itaufin, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
247 dudyn, dv, dteta, dp, w, pbaru, pbarv, time_0, &
248 conser = MOD(itaufin, iconser) == 0)
249
250 END SUBROUTINE leapfrog
251
252 end module leapfrog_m

  ViewVC Help
Powered by ViewVC 1.1.21