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

Annotation of /trunk/dyn3d/leapfrog.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 224 - (hide annotations)
Fri Apr 28 13:40:59 2017 UTC (7 years ago) by guez
Original Path: trunk/Sources/dyn3d/leapfrog.f
File size: 8428 byte(s)
offline allows to save variables for another run with offline
transport, but there is no provision for this other run with offline
transport in LMDZE. The procedures read_pstoke and read_pstoke0 in
LMDZ are never called. So removed the possibility offline = T.

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

  ViewVC Help
Powered by ViewVC 1.1.21