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

Contents of /trunk/dyn3d/leapfrog.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 266 - (show annotations)
Thu Apr 19 17:54:55 2018 UTC (6 years ago) by guez
File size: 8130 byte(s)
Define macros of the preprocessor CPP_IIM, CPP_JJM, CPP_LLM so we can
control the resolution from the compilation command, and automate
compilation for several resolutions.

In module yoethf_m, transform variables into named constants. So we do
not need procedure yoethf any longer.

Bug fix in program test_inter_barxy, missing calls to fyhyp and fxhyp,
and definition of rlatu.

Remove variable iecri of module conf_gcm_m. The files dyn_hist*.nc are
written every time step. We are simplifying the output system, pending
replacement by a whole new system.

Modify possible value of vert_sampling from "param" to
"strato_custom", following LMDZ. Default values of corresponding
namelist variables are now the values used for LMDZ CMIP6.

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

  ViewVC Help
Powered by ViewVC 1.1.21