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

Contents of /trunk/dyn3d/leapfrog.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 328 - (show annotations)
Thu Jun 13 14:40:06 2019 UTC (4 years, 11 months ago) by guez
File size: 8053 byte(s)
Change all `.f` suffixes to `.f90`. (The opposite was done in revision
82.)  Because of change of philosopy in GNUmakefile: we already had a
rewritten rule for `.f`, so it does not make the makefile longer to
replace it by a rule for `.f90`. And it spares us options of
makedepf90 and of the compiler. Also we prepare the way for a simpler
`CMakeLists.txt`.

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

  ViewVC Help
Powered by ViewVC 1.1.21