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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 138 - (hide annotations)
Fri May 22 23:13:19 2015 UTC (8 years, 11 months ago) by guez
File size: 8626 byte(s)
Moved variable nb_files from module histcom_var to module
histbeg_totreg_m.

Removed unused argument q of writehist.

No history file is created in program ce0l so there is no need to call
histclo in etat0.

In phyredem, access variables rlat and rlon directly from module
phyetat0_m instead of having them as arguments. This is clearer for
the program gcm. There are bad side effects for the program ce0l: we
have to modify the module variables rlat and rlon in procedure etat0,
and we need the additional file phyetat0.f to compile ce0l.

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

  ViewVC Help
Powered by ViewVC 1.1.21