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

Contents of /trunk/dyn3d/leapfrog.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 254 - (show annotations)
Mon Feb 5 10:39:38 2018 UTC (6 years, 3 months ago) by guez
File size: 8423 byte(s)
Move Sources/* to root directory.
1 module leapfrog_m
2
3 IMPLICIT NONE
4
5 contains
6
7 SUBROUTINE leapfrog(ucov, vcov, teta, ps, masse, phis, q)
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
12 ! Intégration temporelle du modèle : Matsuno-leapfrog scheme.
13
14 use addfi_m, only: addfi
15 use bilan_dyn_m, only: bilan_dyn
16 use caladvtrac_m, only: caladvtrac
17 use caldyn_m, only: caldyn
18 USE calfis_m, ONLY: calfis
19 USE comconst, ONLY: dtvr
20 USE comgeom, ONLY: aire_2d, apoln, apols
21 use covcont_m, only: covcont
22 USE disvert_m, ONLY: ap, bp
23 USE conf_gcm_m, ONLY: day_step, iconser, iperiod, iphysiq, nday, &
24 iflag_phys, iecri
25 USE conf_guide_m, ONLY: ok_guide
26 USE dimens_m, ONLY: iim, jjm, llm, nqmx
27 use dissip_m, only: dissip
28 USE dynetat0_m, ONLY: day_ini
29 use dynredem1_m, only: dynredem1
30 use enercin_m, only: enercin
31 USE exner_hyb_m, ONLY: exner_hyb
32 use filtreg_scal_m, only: filtreg_scal
33 use geopot_m, only: geopot
34 USE guide_m, ONLY: guide
35 use inidissip_m, only: idissip
36 use integrd_m, only: integrd
37 use nr_util, only: assert
38 USE temps, ONLY: itau_dyn
39 use writedynav_m, only: writedynav
40 use writehist_m, only: writehist
41
42 ! Variables dynamiques:
43 REAL, intent(inout):: ucov(:, :, :) ! (iim + 1, jjm + 1, llm) vent covariant
44 REAL, intent(inout):: vcov(:, :, :) ! (iim + 1, jjm, llm) ! vent covariant
45
46 REAL, intent(inout):: teta(:, :, :) ! (iim + 1, jjm + 1, llm)
47 ! potential temperature
48
49 REAL, intent(inout):: ps(:, :) ! (iim + 1, jjm + 1) pression au sol, en Pa
50 REAL, intent(inout):: masse(:, :, :) ! (iim + 1, jjm + 1, llm) masse d'air
51 REAL, intent(in):: phis(:, :) ! (iim + 1, jjm + 1) surface geopotential
52
53 REAL, intent(inout):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx)
54 ! mass fractions of advected fields
55
56 ! Local:
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\'e 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)
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
88 ! Variables pour le fichier histoire
89 INTEGER itau ! index of the time step of the dynamics, starts at 0
90 INTEGER itaufin
91 INTEGER l
92
93 ! Variables test conservation \'energie
94 REAL ecin(iim + 1, jjm + 1, llm), ecin0(iim + 1, jjm + 1, llm)
95
96 REAL vcont((iim + 1) * jjm, llm), ucont((iim + 1) * (jjm + 1), llm)
97 logical leapf
98 real dt ! time step, in s
99
100 REAL p3d(iim + 1, jjm + 1, llm + 1) ! pressure at layer interfaces, in Pa
101 ! ("p3d(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)",
102 ! for interface "l")
103
104 !---------------------------------------------------
105
106 print *, "Call sequence information: leapfrog"
107 call assert(shape(ucov) == (/iim + 1, jjm + 1, llm/), "leapfrog")
108
109 itaufin = nday * day_step
110 ! "day_step" is a multiple of "iperiod", therefore so is "itaufin".
111
112 ! On initialise la pression et la fonction d'Exner :
113 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
114 CALL exner_hyb(ps, p3d, pks, pk)
115 pkf = pk
116 CALL filtreg_scal(pkf, direct = .true., intensive = .true.)
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) call guide(itau, ucov, vcov, teta, q(:, :, :, 1), ps)
126 vcovm1 = vcov
127 ucovm1 = ucov
128 tetam1 = teta
129 massem1 = masse
130 psm1 = ps
131 end if
132
133 ! Calcul des tendances dynamiques:
134 CALL geopot(teta, pk, pks, phis, phi)
135 CALL caldyn(itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
136 dudyn, dv, dteta, dp, w, pbaru, pbarv, &
137 conser = MOD(itau, iconser) == 0)
138
139 CALL caladvtrac(q, pbaru, pbarv, p3d, masse, teta, pk)
140
141 ! Int\'egrations dynamique et traceurs:
142 CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, dteta, &
143 dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, dt, leapf)
144
145 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
146 CALL exner_hyb(ps, p3d, pks, pk)
147 pkf = pk
148 CALL filtreg_scal(pkf, direct = .true., intensive = .true.)
149
150 if (.not. leapf) then
151 ! Matsuno backward
152 ! Calcul des tendances dynamiques:
153 CALL geopot(teta, pk, pks, phis, phi)
154 CALL caldyn(itau + 1, ucov, vcov, teta, ps, masse, pk, pkf, phis, &
155 phi, dudyn, dv, dteta, dp, w, pbaru, pbarv, conser = .false.)
156
157 ! integrations dynamique et traceurs:
158 CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, dudyn, &
159 dteta, dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, dtvr, &
160 leapf=.false.)
161
162 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
163 CALL exner_hyb(ps, p3d, pks, pk)
164 pkf = pk
165 CALL filtreg_scal(pkf, direct = .true., intensive = .true.)
166 end if
167
168 IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys) THEN
169 CALL calfis(ucov, vcov, teta, q, p3d, pk, phis, phi, w, dufi, dvfi, &
170 dtetafi, dqfi, dayvrai = itau / day_step + day_ini, &
171 time = REAL(mod(itau, day_step)) / day_step, &
172 lafin = itau + 1 == itaufin)
173
174 CALL addfi(ucov, vcov, teta, q, dufi, dvfi, dtetafi, dqfi)
175 ENDIF
176
177 IF (MOD(itau + 1, idissip) == 0) THEN
178 ! Dissipation horizontale et verticale des petites \'echelles
179
180 ! calcul de l'\'energie cin\'etique avant dissipation
181 call covcont(llm, ucov, vcov, ucont, vcont)
182 call enercin(vcov, ucov, vcont, ucont, ecin0)
183
184 ! dissipation
185 CALL dissip(vcov, ucov, teta, p3d, dvdis, dudis, dtetadis)
186 ucov = ucov + dudis
187 vcov = vcov + dvdis
188
189 ! On ajoute la tendance due \`a la transformation \'energie
190 ! cin\'etique en \'energie thermique par la dissipation
191 call covcont(llm, ucov, vcov, ucont, vcont)
192 call enercin(vcov, ucov, vcont, ucont, ecin)
193 dtetadis = dtetadis + (ecin0 - ecin) / pk
194 teta = teta + dtetadis
195
196 ! Calcul de la valeur moyenne aux p\^oles :
197 forall (l = 1: llm)
198 teta(:, 1, l) = SUM(aire_2d(:iim, 1) * teta(:iim, 1, l)) &
199 / apoln
200 teta(:, jjm + 1, l) = SUM(aire_2d(:iim, jjm + 1) &
201 * teta(:iim, jjm + 1, l)) / apols
202 END forall
203 END IF
204
205 IF (MOD(itau + 1, iperiod) == 0) THEN
206 ! \'Ecriture du fichier histoire moyenne:
207 CALL writedynav(vcov, ucov, teta, pk, phi, q, masse, ps, phis, &
208 time = itau + 1)
209 call bilan_dyn(ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, &
210 q(:, :, :, 1))
211 ENDIF
212
213 IF (MOD(itau + 1, iecri * day_step) == 0) THEN
214 CALL geopot(teta, pk, pks, phis, phi)
215 CALL writehist(itau, vcov, ucov, teta, phi, masse, ps)
216 END IF
217 end do time_integration
218
219 CALL dynredem1(vcov, ucov, teta, q, masse, ps, itau = itau_dyn + itaufin)
220
221 ! Calcul des tendances dynamiques:
222 CALL geopot(teta, pk, pks, phis, phi)
223 CALL caldyn(itaufin, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
224 dudyn, dv, dteta, dp, w, pbaru, pbarv, &
225 conser = MOD(itaufin, iconser) == 0)
226
227 END SUBROUTINE leapfrog
228
229 end module leapfrog_m

  ViewVC Help
Powered by ViewVC 1.1.21