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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 178 - (show annotations)
Fri Mar 11 18:47:26 2016 UTC (8 years, 2 months ago) by guez
File size: 8605 byte(s)
Moved variables date0, deltat, datasz_max, ncvar_ids, point, buff_pos,
buffer, regular from module histcom_var to modules where they are
defined.

Removed procedure ioipslmpp, useless for a sequential program.

Added argument datasz_max to histwrite_real (to avoid circular
dependency with histwrite).

Removed useless variables and computations everywhere.

Changed real litteral constants from default kind to double precision
in lwb, lwu, lwvn, sw1s, swtt, swtt1, swu.

Removed unused arguments: paer of sw, sw1s, sw2s, swclr; pcldsw of
sw1s, sw2s; pdsig, prayl of swr; co2_ppm of clmain, clqh; tsol of
transp_lay; nsrf of screenp; kcrit and kknu of gwstress; pstd of
orosetup.

Added output of relative humidity.

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

  ViewVC Help
Powered by ViewVC 1.1.21