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

Annotation of /trunk/dyn3d/leapfrog.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 339 - (hide annotations)
Thu Sep 26 17:08:42 2019 UTC (4 years, 7 months ago) by guez
File size: 8065 byte(s)
Simplify newmicro and rename dummy arguments

Rename dummy argument ucov of procedure advect to ang_3d. Rename dummy
argument radliq of procedure fisrtilp to cldliq. Rename dummy argument
qlwp of procedure newmicro to cldliq. Motivation: same name across
procedures.

Remove useless intermediary local variable rel in procedure
newmicro. The value of rad_chaud can be used instead of 1 in the
computation of cldtau if zflwp is 0: it does not matter.

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

  ViewVC Help
Powered by ViewVC 1.1.21