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

Contents of /trunk/dyn3d/leapfrog.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 265 - (show annotations)
Tue Mar 20 09:35:59 2018 UTC (6 years, 1 month ago) by guez
File size: 8249 byte(s)
Rename module dimens_m to dimensions.
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 dimensions, 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 writehist_m, only: writehist
40
41 ! Variables dynamiques:
42 REAL, intent(inout):: ucov(:, :, :) ! (iim + 1, jjm + 1, llm) vent covariant
43 REAL, intent(inout):: vcov(:, :, :) ! (iim + 1, jjm, llm) ! vent covariant
44
45 REAL, intent(inout):: teta(:, :, :) ! (iim + 1, jjm + 1, llm)
46 ! potential temperature
47
48 REAL, intent(inout):: ps(:, :) ! (iim + 1, jjm + 1) pression au sol, en Pa
49 REAL, intent(inout):: masse(:, :, :) ! (iim + 1, jjm + 1, llm) masse d'air
50 REAL, intent(in):: phis(:, :) ! (iim + 1, jjm + 1) surface geopotential
51
52 REAL, intent(inout):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx)
53 ! mass fractions of advected fields
54
55 ! Local:
56
57 ! Variables dynamiques:
58
59 REAL pks(iim + 1, jjm + 1) ! exner au sol
60 REAL pk(iim + 1, jjm + 1, llm) ! exner au milieu des couches
61 REAL pkf(iim + 1, jjm + 1, llm) ! exner filtr\'e au milieu des couches
62 REAL phi(iim + 1, jjm + 1, llm) ! geopotential
63 REAL w(iim + 1, jjm + 1, llm) ! vitesse verticale
64
65 ! Variables dynamiques interm\'ediaires pour le transport
66 ! Flux de masse :
67 REAL pbaru(iim + 1, jjm + 1, llm), pbarv(iim + 1, jjm, llm)
68
69 ! Variables dynamiques au pas - 1
70 REAL vcovm1(iim + 1, jjm, llm), ucovm1(iim + 1, jjm + 1, llm)
71 REAL tetam1(iim + 1, jjm + 1, llm), psm1(iim + 1, jjm + 1)
72 REAL massem1(iim + 1, jjm + 1, llm)
73
74 ! Tendances dynamiques
75 REAL dv((iim + 1) * jjm, llm), du(iim + 1, jjm + 1, llm)
76 REAL dteta(iim + 1, jjm + 1, llm)
77 real dp(iim + 1, jjm + 1)
78
79 ! Tendances de la dissipation :
80 REAL dvdis(iim + 1, jjm, llm), dudis(iim + 1, jjm + 1, llm)
81 REAL dtetadis(iim + 1, jjm + 1, llm)
82
83 ! Tendances physiques
84 REAL dvfi(iim + 1, jjm, llm), dufi(iim + 1, jjm + 1, llm)
85 REAL dtetafi(iim + 1, jjm + 1, llm), dqfi(iim + 1, jjm + 1, llm, nqmx)
86
87 ! Variables pour le fichier histoire
88 INTEGER itau ! index of the time step of the dynamics, starts at 0
89 INTEGER itaufin
90 INTEGER l
91
92 ! Variables test conservation \'energie
93 REAL ecin(iim + 1, jjm + 1, llm), ecin0(iim + 1, jjm + 1, llm)
94
95 REAL vcont((iim + 1) * jjm, llm), ucont((iim + 1) * (jjm + 1), llm)
96 logical leapf
97 real dt ! time step, in s
98
99 REAL p3d(iim + 1, jjm + 1, llm + 1) ! pressure at layer interfaces, in Pa
100 ! ("p3d(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)",
101 ! for interface "l")
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 so is "itaufin".
110
111 ! On initialise la pression et la fonction d'Exner :
112 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
113 CALL exner_hyb(ps, p3d, pks, pk)
114 pkf = pk
115 CALL filtreg_scal(pkf, direct = .true., intensive = .true.)
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) call guide(itau, ucov, vcov, teta, q(:, :, :, 1), ps)
125 vcovm1 = vcov
126 ucovm1 = ucov
127 tetam1 = teta
128 massem1 = masse
129 psm1 = ps
130 end if
131
132 ! Calcul des tendances dynamiques:
133 CALL geopot(teta, pk, pks, phis, phi)
134 CALL caldyn(itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
135 du, dv, dteta, dp, w, pbaru, pbarv, &
136 conser = MOD(itau, iconser) == 0)
137
138 CALL caladvtrac(q, pbaru, pbarv, p3d, masse, teta, pk)
139
140 ! Int\'egrations dynamique et traceurs:
141 CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, du, dteta, &
142 dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, dt, leapf)
143
144 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
145 CALL exner_hyb(ps, p3d, pks, pk)
146 pkf = pk
147 CALL filtreg_scal(pkf, direct = .true., intensive = .true.)
148
149 if (.not. leapf) then
150 ! Matsuno backward
151 ! Calcul des tendances dynamiques:
152 CALL geopot(teta, pk, pks, phis, phi)
153 CALL caldyn(itau + 1, ucov, vcov, teta, ps, masse, pk, pkf, phis, &
154 phi, du, dv, dteta, dp, w, pbaru, pbarv, conser = .false.)
155
156 ! integrations dynamique et traceurs:
157 CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, du, &
158 dteta, dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, dtvr, &
159 leapf=.false.)
160
161 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
162 CALL exner_hyb(ps, p3d, pks, pk)
163 pkf = pk
164 CALL filtreg_scal(pkf, direct = .true., intensive = .true.)
165 end if
166
167 IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys) THEN
168 CALL calfis(ucov, vcov, teta, q, p3d, pk, phis, phi, w, dufi, dvfi, &
169 dtetafi, dqfi, dayvrai = itau / day_step + day_ini, &
170 time = REAL(mod(itau, day_step)) / day_step, &
171 lafin = itau + 1 == itaufin)
172
173 CALL addfi(ucov, vcov, teta, q, dufi, dvfi, dtetafi, dqfi)
174 ENDIF
175
176 IF (MOD(itau + 1, idissip) == 0) THEN
177 ! Dissipation horizontale et verticale des petites \'echelles
178
179 ! calcul de l'\'energie cin\'etique avant dissipation
180 call covcont(llm, ucov, vcov, ucont, vcont)
181 call enercin(vcov, ucov, vcont, ucont, ecin0)
182
183 ! dissipation
184 CALL dissip(vcov, ucov, teta, p3d, dvdis, dudis, dtetadis)
185 ucov = ucov + dudis
186 vcov = vcov + dvdis
187
188 ! On ajoute la tendance due \`a la transformation \'energie
189 ! cin\'etique en \'energie thermique par la dissipation
190 call covcont(llm, ucov, vcov, ucont, vcont)
191 call enercin(vcov, ucov, vcont, ucont, ecin)
192 dtetadis = dtetadis + (ecin0 - ecin) / pk
193 teta = teta + dtetadis
194
195 ! Calcul de la valeur moyenne aux p\^oles :
196 forall (l = 1: llm)
197 teta(:, 1, l) = SUM(aire_2d(:iim, 1) * teta(:iim, 1, l)) &
198 / apoln
199 teta(:, jjm + 1, l) = SUM(aire_2d(:iim, jjm + 1) &
200 * teta(:iim, jjm + 1, l)) / apols
201 END forall
202 END IF
203
204 IF (MOD(itau + 1, iperiod) == 0) THEN
205 call bilan_dyn(ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, &
206 q(:, :, :, 1))
207 ENDIF
208
209 IF (MOD(itau + 1, iecri) == 0) THEN
210 CALL geopot(teta, pk, pks, phis, phi)
211 CALL writehist(vcov, ucov, teta, pk, phi, q, masse, ps, &
212 itau_w = itau_dyn + itau + 1)
213 END IF
214 end do time_integration
215
216 CALL dynredem1(vcov, ucov, teta, q, masse, ps, itau = itau_dyn + itaufin)
217
218 ! Calcul des tendances dynamiques:
219 CALL geopot(teta, pk, pks, phis, phi)
220 CALL caldyn(itaufin, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
221 du, dv, dteta, dp, w, pbaru, pbarv, &
222 conser = MOD(itaufin, iconser) == 0)
223
224 END SUBROUTINE leapfrog
225
226 end module leapfrog_m

  ViewVC Help
Powered by ViewVC 1.1.21