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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 61 - (show annotations)
Fri Apr 20 14:58:43 2012 UTC (12 years ago) by guez
File size: 9031 byte(s)
No more included file in LMDZE, not even "netcdf.inc".

Created a variable containing the list of common source files in
GNUmakefile. So we now also see clearly files that are specific to
each program.

Split module "histcom". Assembled resulting files in directory
"Histcom".

Removed aliasing in calls to "laplacien".

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

  ViewVC Help
Powered by ViewVC 1.1.21