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

Annotation of /trunk/dyn3d/leapfrog.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 41 - (hide annotations)
Tue Feb 22 15:09:57 2011 UTC (13 years, 2 months ago) by guez
Original Path: trunk/libf/dyn3d/leapfrog.f90
File size: 8765 byte(s)
Moved "Test_ozonecm" out of SVN.
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 27 ! From dyn3d/leapfrog.F, version 1.6, 2005/04/13 08:58:34
10     ! 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 26 USE calfis_m, ONLY: calfis
17     USE com_io_dyn, ONLY: histaveid
18     USE comconst, ONLY: daysec, dtphys, dtvr
19 guez 29 USE comgeom, ONLY: aire_2d, apoln, apols
20 guez 26 USE comvert, ONLY: ap, bp
21 guez 27 USE conf_gcm_m, ONLY: day_step, iconser, iperiod, iphysiq, nday, offline, &
22     periodav
23 guez 29 USE dimens_m, ONLY: iim, jjm, llm, nqmx
24 guez 26 USE dynetat0_m, ONLY: day_ini
25 guez 27 use dynredem1_m, only: dynredem1
26 guez 26 USE exner_hyb_m, ONLY: exner_hyb
27 guez 27 use filtreg_m, only: filtreg
28 guez 26 USE guide_m, ONLY: guide
29     use inidissip_m, only: idissip
30 guez 32 use integrd_m, only: integrd
31 guez 26 USE logic, ONLY: iflag_phys, ok_guide
32 guez 29 USE paramet_m, ONLY: ip1jmp1
33 guez 26 USE pressure_var, ONLY: p3d
34 guez 28 USE temps, ONLY: itau_dyn
35 guez 3
36 guez 10 ! Variables dynamiques:
37 guez 36 REAL, intent(inout):: ucov(ip1jmp1, llm) ! vent covariant
38 guez 32 REAL, intent(inout):: vcov((iim + 1) * jjm, llm) ! vent covariant
39 guez 29 REAL, intent(inout):: teta(iim + 1, jjm + 1, llm) ! potential temperature
40 guez 41 REAL, intent(inout):: ps(iim + 1, jjm + 1) ! pression au sol, en Pa
41 guez 10 REAL masse(ip1jmp1, llm) ! masse d'air
42     REAL phis(ip1jmp1) ! geopotentiel au sol
43 guez 40
44     REAL, intent(inout):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx)
45     ! mass fractions of advected fields
46    
47 guez 24 REAL, intent(in):: time_0
48 guez 10
49     ! Variables local to the procedure:
50    
51     ! Variables dynamiques:
52    
53 guez 3 REAL pks(ip1jmp1) ! exner au sol
54 guez 29 REAL pk(iim + 1, jjm + 1, llm) ! exner au milieu des couches
55 guez 3 REAL pkf(ip1jmp1, llm) ! exner filt.au milieu des couches
56     REAL phi(ip1jmp1, llm) ! geopotential
57     REAL w(ip1jmp1, llm) ! vitesse verticale
58    
59     ! variables dynamiques intermediaire pour le transport
60 guez 29 REAL pbaru(ip1jmp1, llm), pbarv((iim + 1) * jjm, llm) !flux de masse
61 guez 3
62     ! variables dynamiques au pas - 1
63 guez 29 REAL vcovm1((iim + 1) * jjm, llm), ucovm1(ip1jmp1, llm)
64     REAL tetam1(iim + 1, jjm + 1, llm), psm1(iim + 1, jjm + 1)
65 guez 3 REAL massem1(ip1jmp1, llm)
66    
67     ! tendances dynamiques
68 guez 29 REAL dv((iim + 1) * jjm, llm), du(ip1jmp1, llm)
69 guez 3 REAL dteta(ip1jmp1, llm), dq(ip1jmp1, llm, nqmx), dp(ip1jmp1)
70    
71     ! tendances de la dissipation
72 guez 29 REAL dvdis((iim + 1) * jjm, llm), dudis(ip1jmp1, llm)
73     REAL dtetadis(iim + 1, jjm + 1, llm)
74 guez 3
75     ! tendances physiques
76 guez 29 REAL dvfi((iim + 1) * jjm, llm), dufi(ip1jmp1, llm)
77 guez 3 REAL dtetafi(ip1jmp1, llm), dqfi(ip1jmp1, llm, nqmx), dpfi(ip1jmp1)
78    
79     ! variables pour le fichier histoire
80    
81 guez 22 INTEGER itau ! index of the time step of the dynamics, starts at 0
82 guez 27 INTEGER itaufin
83 guez 20 REAL time ! time of day, as a fraction of day length
84 guez 10 real finvmaold(ip1jmp1, llm)
85 guez 33 INTEGER l
86 guez 3 REAL rdayvrai, rdaym_ini
87    
88 guez 24 ! Variables test conservation energie
89 guez 29 REAL ecin(iim + 1, jjm + 1, llm), ecin0(iim + 1, jjm + 1, llm)
90 guez 3 ! Tendance de la temp. potentiel d (theta) / d t due a la
91     ! tansformation d'energie cinetique en energie thermique
92     ! cree par la dissipation
93 guez 29 REAL dtetaecdt(iim + 1, jjm + 1, llm)
94     REAL vcont((iim + 1) * jjm, llm), ucont(ip1jmp1, llm)
95 guez 33 logical leapf
96     real dt
97 guez 3
98     !---------------------------------------------------
99    
100     print *, "Call sequence information: leapfrog"
101    
102     itaufin = nday * day_step
103 guez 30 ! "day_step" is a multiple of "iperiod", therefore "itaufin" is one too
104    
105 guez 24 dq = 0.
106 guez 33
107 guez 3 ! On initialise la pression et la fonction d'Exner :
108 guez 37 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
109 guez 10 CALL exner_hyb(ps, p3d, pks, pk, pkf)
110 guez 3
111 guez 40 time_integration: do itau = 0, itaufin - 1
112 guez 33 leapf = mod(itau, iperiod) /= 0
113     if (leapf) then
114     dt = 2 * dtvr
115     else
116     ! Matsuno
117     dt = dtvr
118     if (ok_guide .and. (itaufin - itau - 1) * dtvr > 21600.) &
119     call guide(itau, ucov, vcov, teta, q, masse, ps)
120     vcovm1 = vcov
121     ucovm1 = ucov
122     tetam1 = teta
123     massem1 = masse
124     psm1 = ps
125     finvmaold = masse
126     CALL filtreg(finvmaold, jjm + 1, llm, - 2, 2, .TRUE., 1)
127     end if
128 guez 30
129     ! Calcul des tendances dynamiques:
130     CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)
131     CALL caldyn(itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
132     MOD(itau, iconser) == 0, du, dv, dteta, dp, w, pbaru, pbarv, &
133 guez 33 time_0)
134 guez 30
135     ! Calcul des tendances advection des traceurs (dont l'humidité)
136     CALL caladvtrac(q, pbaru, pbarv, p3d, masse, dq, teta, pk)
137 guez 33
138 guez 30 ! Stokage du flux de masse pour traceurs offline:
139     IF (offline) CALL fluxstokenc(pbaru, pbarv, masse, teta, phi, phis, &
140     dtvr, itau)
141    
142     ! integrations dynamique et traceurs:
143 guez 39 CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, du, dteta, dp, &
144 guez 40 vcov, ucov, teta, q(:, :, :, :2), ps, masse, finvmaold, dt, leapf)
145 guez 30
146 guez 33 if (.not. leapf) then
147     ! Matsuno backward
148 guez 37 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
149 guez 33 CALL exner_hyb(ps, p3d, pks, pk, pkf)
150 guez 30
151 guez 27 ! Calcul des tendances dynamiques:
152 guez 3 CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)
153 guez 33 CALL caldyn(itau + 1, ucov, vcov, teta, ps, masse, pk, pkf, phis, &
154     phi, .false., du, dv, dteta, dp, w, pbaru, pbarv, time_0)
155 guez 3
156     ! integrations dynamique et traceurs:
157 guez 39 CALL integrd(vcovm1, ucovm1, tetam1, psm1, massem1, dv, du, dteta, &
158 guez 40 dp, vcov, ucov, teta, q(:, :, :, :2), ps, masse, finvmaold, &
159     dtvr, leapf=.false.)
160 guez 33 end if
161 guez 3
162 guez 33 IF (MOD(itau + 1, iphysiq) == 0 .AND. iflag_phys /= 0) THEN
163     ! calcul des tendances physiques:
164 guez 3
165 guez 37 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
166 guez 33 CALL exner_hyb(ps, p3d, pks, pk, pkf)
167 guez 3
168 guez 33 rdaym_ini = itau * dtvr / daysec
169     rdayvrai = rdaym_ini + day_ini
170     time = REAL(mod(itau, day_step)) / day_step + time_0
171     IF (time > 1.) time = time - 1.
172 guez 3
173 guez 37 CALL calfis(rdayvrai, time, ucov, vcov, teta, q, masse, ps, pk, &
174     phis, phi, du, dv, dteta, dq, w, dufi, dvfi, dtetafi, dqfi, &
175     dpfi, lafin=itau+1==itaufin)
176 guez 3
177 guez 33 ! ajout des tendances physiques:
178     CALL addfi(nqmx, dtphys, ucov, vcov, teta, q, ps, dufi, dvfi, &
179     dtetafi, dqfi, dpfi)
180     ENDIF
181 guez 3
182 guez 37 forall (l = 1: llm + 1) p3d(:, :, l) = ap(l) + bp(l) * ps
183 guez 33 CALL exner_hyb(ps, p3d, pks, pk, pkf)
184 guez 3
185 guez 33 IF (MOD(itau + 1, idissip) == 0) THEN
186     ! dissipation horizontale et verticale des petites echelles:
187 guez 3
188 guez 33 ! calcul de l'energie cinetique avant dissipation
189     call covcont(llm, ucov, vcov, ucont, vcont)
190     call enercin(vcov, ucov, vcont, ucont, ecin0)
191 guez 3
192 guez 33 ! dissipation
193     CALL dissip(vcov, ucov, teta, p3d, dvdis, dudis, dtetadis)
194     ucov=ucov + dudis
195     vcov=vcov + dvdis
196 guez 3
197 guez 33 ! On rajoute la tendance due à la transformation Ec -> E
198     ! thermique créée lors de la dissipation
199     call covcont(llm, ucov, vcov, ucont, vcont)
200     call enercin(vcov, ucov, vcont, ucont, ecin)
201     dtetaecdt= (ecin0 - ecin) / pk
202     dtetadis=dtetadis + dtetaecdt
203     teta=teta + dtetadis
204 guez 3
205 guez 33 ! Calcul de la valeur moyenne aux pôles :
206     forall (l = 1: llm)
207     teta(:, 1, l) = SUM(aire_2d(:iim, 1) * teta(:iim, 1, l)) &
208     / apoln
209     teta(:, jjm + 1, l) = SUM(aire_2d(:iim, jjm+1) &
210     * teta(:iim, jjm + 1, l)) / apols
211     END forall
212 guez 3
213 guez 33 ps(:, 1) = SUM(aire_2d(:iim, 1) * ps(:iim, 1)) / apoln
214     ps(:, jjm + 1) = SUM(aire_2d(:iim, jjm+1) * ps(:iim, jjm + 1)) &
215     / apols
216     END IF
217 guez 3
218 guez 33 IF (MOD(itau + 1, iperiod) == 0) THEN
219 guez 40 ! Écriture du fichier histoire moyenne:
220 guez 33 CALL writedynav(histaveid, nqmx, itau + 1, vcov, ucov, teta, pk, &
221     phi, q, masse, ps, phis)
222 guez 40 call bilan_dyn(ps, masse, pk, pbaru, pbarv, teta, phi, ucov, vcov, &
223     q(:, :, :, 1), dt_app = dtvr * iperiod, &
224     dt_cum = dtvr * day_step * periodav)
225 guez 33 ENDIF
226 guez 40 end do time_integration
227 guez 3
228 guez 30 CALL dynredem1("restart.nc", vcov, ucov, teta, q, masse, ps, &
229     itau=itau_dyn+itaufin)
230    
231     ! Calcul des tendances dynamiques:
232     CALL geopot(ip1jmp1, teta, pk, pks, phis, phi)
233     CALL caldyn(itaufin, ucov, vcov, teta, ps, masse, pk, pkf, phis, phi, &
234     MOD(itaufin, iconser) == 0, du, dv, dteta, dp, w, pbaru, pbarv, &
235 guez 33 time_0)
236 guez 30
237 guez 3 END SUBROUTINE leapfrog
238    
239     end module leapfrog_m

  ViewVC Help
Powered by ViewVC 1.1.21