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

Contents of /trunk/dyn3d/leapfrog.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 108 - (show annotations)
Tue Sep 16 14:00:41 2014 UTC (9 years, 7 months ago) by guez
File size: 8677 byte(s)
Imported writefield from LMDZ. Close at the end of gcm the files which
were created by writefiled (not done in LMDZ).

Removed procedures for the output of Grads files. Removed calls to
dump2d. In guide, replaced calls to wrgrads by calls to writefield.

In vlspltqs, removed redundant programming of saturation
pressure. Call foeew from module FCTTRE instead.

Bug fix in interpre: size of w exceeding size of correponding actual
argument wg in advtrac.

In leapfrog, call guide until the end of the run, instead of six hours
before the end.

Bug fix in readsulfate_preind: type of arguments.

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

  ViewVC Help
Powered by ViewVC 1.1.21