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

Annotation of /trunk/dyn3d/leapfrog.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 88 - (hide annotations)
Tue Mar 11 15:09:02 2014 UTC (10 years, 2 months ago) by guez
File size: 8944 byte(s)
Removed useless argument mode of subroutine read_reanalyse.

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

  ViewVC Help
Powered by ViewVC 1.1.21