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

Annotation of /trunk/Sources/dyn3d/leapfrog.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 207 - (hide annotations)
Thu Sep 1 10:30:53 2016 UTC (7 years, 8 months ago) by guez
File size: 8638 byte(s)
New philosophy on compiler options.

Removed source code for thermcep = f. (Not used in LMDZ either.)

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

  ViewVC Help
Powered by ViewVC 1.1.21