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

Annotation of /trunk/libf/dyn3d/leapfrog.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 34 - (hide annotations)
Wed Jun 2 11:01:12 2010 UTC (13 years, 11 months ago) by guez
File size: 8496 byte(s)
Split "ini_hist.f90" into single-procedure files.

In "calfis" and "physiq", removed dummy argument "nq" since "nq" must
be equal to "nqmx".

In "calfis", renamed dummy argument "pq" to "q", same name as actual
argument in "leapfrog". Renamed local variable "zqfi" to "qx", same
name as dummy argument in "physiq".

Removed arguments "itop_con" and "ibas_con" of "phytrac", which were
not used.

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

  ViewVC Help
Powered by ViewVC 1.1.21