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

Contents of /trunk/dyn3d/leapfrog.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 43 - (show annotations)
Fri Apr 8 12:43:31 2011 UTC (13 years, 1 month ago) by guez
Original Path: trunk/libf/dyn3d/leapfrog.f90
File size: 8833 byte(s)
"start_init_phys" is now called directly by "etat0" instead of through
"start_init_dyn". "qsol_2d" is no longer a variable of module
"start_init_phys_m", it is an argument of
"start_init_phys". "start_init_dyn" now receives "tsol_2d" from
"etat0".

Split file "vlspltqs.f" into "vlspltqs.f90", "vlxqs.f90" and
""vlyqs.f90".

In "start_init_orog", replaced calls to "flin*" by calls to NetCDF95.

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 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 logic, ONLY: iflag_phys, ok_guide
34 USE paramet_m, ONLY: ip1jmp1
35 USE pressure_var, ONLY: p3d
36 USE temps, ONLY: itau_dyn
37
38 ! Variables dynamiques:
39 REAL, intent(inout):: ucov(ip1jmp1, 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(ip1jmp1, llm) ! masse d'air
47 REAL phis(ip1jmp1) ! 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(ip1jmp1) ! exner au sol
59 REAL pk(iim + 1, jjm + 1, llm) ! exner au milieu des couches
60 REAL pkf(ip1jmp1, llm) ! exner filt.au milieu des couches
61 REAL phi(ip1jmp1, llm) ! geopotential
62 REAL w(ip1jmp1, llm) ! vitesse verticale
63
64 ! variables dynamiques intermediaire pour le transport
65 REAL pbaru(ip1jmp1, llm), pbarv((iim + 1) * jjm, llm) !flux de masse
66
67 ! variables dynamiques au pas - 1
68 REAL vcovm1((iim + 1) * jjm, llm), ucovm1(ip1jmp1, llm)
69 REAL tetam1(iim + 1, jjm + 1, llm), psm1(iim + 1, jjm + 1)
70 REAL massem1(ip1jmp1, llm)
71
72 ! tendances dynamiques
73 REAL dv((iim + 1) * jjm, llm), du(ip1jmp1, llm)
74 REAL dteta(ip1jmp1, llm), dq(ip1jmp1, llm, nqmx), dp(ip1jmp1)
75
76 ! tendances de la dissipation
77 REAL dvdis((iim + 1) * jjm, llm), dudis(ip1jmp1, llm)
78 REAL dtetadis(iim + 1, jjm + 1, llm)
79
80 ! tendances physiques
81 REAL dvfi((iim + 1) * jjm, llm), dufi(ip1jmp1, llm)
82 REAL dtetafi(ip1jmp1, llm), dqfi(ip1jmp1, llm, nqmx), dpfi(ip1jmp1)
83
84 ! variables pour le fichier histoire
85
86 INTEGER itau ! index of the time step of the dynamics, starts at 0
87 INTEGER itaufin
88 REAL time ! time of day, as a fraction of day length
89 real finvmaold(ip1jmp1, llm)
90 INTEGER l
91 REAL rdayvrai, rdaym_ini
92
93 ! Variables test conservation energie
94 REAL ecin(iim + 1, jjm + 1, llm), ecin0(iim + 1, jjm + 1, llm)
95
96 REAL dtetaecdt(iim + 1, jjm + 1, llm)
97 ! tendance de la température potentielle due à la tansformation
98 ! d'énergie cinétique en énergie thermique créée par la dissipation
99
100 REAL vcont((iim + 1) * jjm, llm), ucont(ip1jmp1, llm)
101 logical leapf
102 real dt
103
104 !---------------------------------------------------
105
106 print *, "Call sequence information: 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(ip1jmp1, teta, pk, pks, phis, phi)
137 CALL caldyn(itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
138 MOD(itau, iconser) == 0, du, dv, dteta, dp, w, pbaru, pbarv, &
139 time_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, du, dteta, dp, &
150 vcov, ucov, teta, q(:, :, :, :2), ps, masse, finvmaold, dt, leapf)
151
152 if (.not. leapf) then
153 ! Matsuno backward
154 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
155 CALL exner_hyb(ps, p3d, pks, pk, pkf)
156
157 ! Calcul des tendances dynamiques:
158 CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)
159 CALL caldyn(itau + 1, ucov, vcov, teta, ps, masse, pk, pkf, phis, &
160 phi, .false., du, dv, dteta, dp, w, pbaru, pbarv, time_0)
161
162 ! integrations dynamique et traceurs:
163 CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, du, dteta, &
164 dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, finvmaold, &
165 dtvr, leapf=.false.)
166 end if
167
168 IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys /= 0) THEN
169 ! calcul des tendances physiques:
170
171 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
172 CALL exner_hyb(ps, p3d, pks, pk, pkf)
173
174 rdaym_ini = itau * dtvr / daysec
175 rdayvrai = rdaym_ini + day_ini
176 time = REAL(mod(itau, day_step)) / day_step + time_0
177 IF (time > 1.) time = time - 1.
178
179 CALL calfis(rdayvrai, time, ucov, vcov, teta, q, masse, ps, pk, &
180 phis, phi, du, dv, dteta, dq, w, dufi, dvfi, dtetafi, dqfi, &
181 dpfi, lafin=itau+1==itaufin)
182
183 ! ajout des tendances physiques:
184 CALL addfi(nqmx, dtphys, ucov, vcov, teta, q, ps, dufi, dvfi, &
185 dtetafi, dqfi, dpfi)
186 ENDIF
187
188 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
189 CALL exner_hyb(ps, p3d, pks, pk, pkf)
190
191 IF (MOD(itau + 1, idissip) == 0) THEN
192 ! dissipation horizontale et verticale des petites echelles:
193
194 ! calcul de l'energie cinetique avant dissipation
195 call covcont(llm, ucov, vcov, ucont, vcont)
196 call enercin(vcov, ucov, vcont, ucont, ecin0)
197
198 ! dissipation
199 CALL dissip(vcov, ucov, teta, p3d, dvdis, dudis, dtetadis)
200 ucov=ucov + dudis
201 vcov=vcov + dvdis
202
203 ! On rajoute la tendance due à la transformation Ec -> E
204 ! thermique créée lors de la dissipation
205 call covcont(llm, ucov, vcov, ucont, vcont)
206 call enercin(vcov, ucov, vcont, ucont, ecin)
207 dtetaecdt= (ecin0 - ecin) / pk
208 dtetadis=dtetadis + dtetaecdt
209 teta=teta + dtetadis
210
211 ! Calcul de la valeur moyenne aux pôles :
212 forall (l = 1: llm)
213 teta(:, 1, l) = SUM(aire_2d(:iim, 1) * teta(:iim, 1, l)) &
214 / apoln
215 teta(:, jjm + 1, l) = SUM(aire_2d(:iim, jjm+1) &
216 * teta(:iim, jjm + 1, l)) / apols
217 END forall
218
219 ps(:, 1) = SUM(aire_2d(:iim, 1) * ps(:iim, 1)) / apoln
220 ps(:, jjm + 1) = SUM(aire_2d(:iim, jjm+1) * ps(:iim, jjm + 1)) &
221 / apols
222 END IF
223
224 IF (MOD(itau + 1, iperiod) == 0) THEN
225 ! Écriture du fichier histoire moyenne:
226 CALL writedynav(histaveid, nqmx, itau + 1, vcov, ucov, teta, pk, &
227 phi, q, masse, ps, phis)
228 call bilan_dyn(ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, &
229 q(:, :, :, 1), dt_app = dtvr * iperiod, &
230 dt_cum = dtvr * day_step * periodav)
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(ip1jmp1, teta, pk, pks, phis, phi)
239 CALL caldyn(itaufin, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
240 MOD(itaufin, iconser) == 0, du, dv, dteta, dp, w, pbaru, pbarv, &
241 time_0)
242
243 END SUBROUTINE leapfrog
244
245 end module leapfrog_m

  ViewVC Help
Powered by ViewVC 1.1.21