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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 37 - (show annotations)
Tue Dec 21 15:45:48 2010 UTC (13 years, 4 months ago) by guez
File size: 8562 byte(s)
Inlined procedure "pression".

Split "guide.f90" into "guide.f90" and "tau2alpha.f90". Split
"read_reanalyse.f" into single-procedure files in directory
"Read_reanalyse".

Useless copy of variables in "iniphysiq". Directly define module
variables in "gcm" and remove procedure "iniphysiq".

Added "pressure-altitude" in "test_disvert".

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

  ViewVC Help
Powered by ViewVC 1.1.21