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

Contents of /trunk/dyn3d/leapfrog.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 313 - (show annotations)
Mon Dec 10 15:54:30 2018 UTC (5 years, 5 months ago) by guez
File size: 8110 byte(s)
Remove module temps. Move variable itau_dyn from module temps to
module dynetat0_m, where it is defined.

Split module dynetat0_m into dynetat0_m and dynetat0_chosen_m. The
motivation is to create smaller modules. Procedures principal_cshift
and invert_zoomx had to stay in dynetat0_m because of circular
dependency. Now we will be able to move them away. Module variables
which are chosen by the user, not computed, in program ce0l go to
dynetat0_chosen_m: day_ref, annee_ref, clon, clat, grossismx,
grossismy, dzoomx, dzoomy, taux, tauy.

Move variable "pa" from module disvert_m to module
dynetat0_chosen_m. Define "pa" in dynetat0_chosen rather than etat0.

Define day_ref and annee_ref in procedure read_serre rather than
etat0.

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

  ViewVC Help
Powered by ViewVC 1.1.21