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

Annotation of /trunk/Sources/dyn3d/leapfrog.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 29 - (hide annotations)
Tue Mar 30 10:44:42 2010 UTC (14 years, 1 month ago) by guez
Original Path: trunk/libf/dyn3d/leapfrog.f90
File size: 8909 byte(s)
In "leapfrog", transformed some arrays with a single dimension for horizontal
position into arrays with two horizontal dimensions. Simplified the
computation of potential temperature and surface pressure at the
poles.

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

  ViewVC Help
Powered by ViewVC 1.1.21