1 | MODULE nonoro_gwd_ran_mod |
---|
2 | |
---|
3 | IMPLICIT NONE |
---|
4 | |
---|
5 | REAL,allocatable,save :: du_nonoro_gwd(:,:) ! Zonal wind tendency due to GWD |
---|
6 | REAL,allocatable,save :: dv_nonoro_gwd(:,:) ! Meridional wind tendency due to GWD |
---|
7 | REAL,ALLOCATABLE,SAVE :: east_gwstress(:,:) ! Profile of eastward stress |
---|
8 | REAL,ALLOCATABLE,SAVE :: west_gwstress(:,:) ! Profile of westward stress |
---|
9 | |
---|
10 | !$OMP THREADPRIVATE(du_nonoro_gwd,dv_nonoro_gwd,east_gwstress,west_gwstress) |
---|
11 | |
---|
12 | CONTAINS |
---|
13 | |
---|
14 | SUBROUTINE NONORO_GWD_RAN(ngrid,nlayer,DTIME, cpnew, rnew, pp, & |
---|
15 | zmax_therm, pt, pu, pv, pdt, pdu, pdv, & |
---|
16 | zustr,zvstr,d_t, d_u, d_v) |
---|
17 | |
---|
18 | !-------------------------------------------------------------------------------- |
---|
19 | ! Parametrization of the momentum flux deposition due to a discrete |
---|
20 | ! number of gravity waves. |
---|
21 | ! F. Lott |
---|
22 | ! Version 14, Gaussian distribution of the source |
---|
23 | ! LMDz model online version |
---|
24 | ! ADAPTED FOR VENUS / F. LOTT + S. LEBONNOIS |
---|
25 | ! Version adapted on 03/04/2013: |
---|
26 | ! - input flux compensated in the deepest layers |
---|
27 | ! |
---|
28 | ! ADAPTED FOR MARS G.GILLI 02/2016 |
---|
29 | ! Revision with F.Forget 06/2016 Variable EP-flux according to |
---|
30 | ! PBL variation (max velocity thermals) |
---|
31 | ! UPDATED D.BARDET 01/2020 - reproductibility of the |
---|
32 | ! launching altitude calculation |
---|
33 | ! - wave characteristic |
---|
34 | ! calculation using MOD |
---|
35 | ! - adding east_gwstress and |
---|
36 | ! west_gwstress variables |
---|
37 | ! UPDATED J.LIU 12/2021 The rho (density) at the specific |
---|
38 | ! locations is introduced. The equation |
---|
39 | ! of EP-flux is corrected by adding the |
---|
40 | ! term of density at launch (source) |
---|
41 | ! altitude. |
---|
42 | ! |
---|
43 | !--------------------------------------------------------------------------------- |
---|
44 | |
---|
45 | use comcstfi_h, only: g, pi, r |
---|
46 | USE ioipsl_getin_p_mod, ONLY : getin_p |
---|
47 | use vertical_layers_mod, only : presnivs |
---|
48 | use geometry_mod, only: cell_area |
---|
49 | #ifdef CPP_XIOS |
---|
50 | use xios_output_mod, only: send_xios_field |
---|
51 | #endif |
---|
52 | implicit none |
---|
53 | include "yoegwd.h" |
---|
54 | include "callkeys.h" |
---|
55 | |
---|
56 | CHARACTER (LEN=20) :: modname='nonoro_gwd_ran' |
---|
57 | CHARACTER (LEN=80) :: abort_message |
---|
58 | |
---|
59 | |
---|
60 | ! 0. DECLARATIONS: |
---|
61 | |
---|
62 | ! 0.1 INPUTS |
---|
63 | INTEGER, intent(in):: ngrid ! number of atmospheric columns |
---|
64 | INTEGER, intent(in):: nlayer ! number of atmospheric layers |
---|
65 | REAL, intent(in):: DTIME ! Time step of the Physics(s) |
---|
66 | REAL, intent(in):: zmax_therm(ngrid) ! Altitude of max velocity thermals (m) |
---|
67 | REAL,INTENT(IN) :: cpnew(ngrid,nlayer)! Cp of the atmosphere |
---|
68 | REAL,INTENT(IN) :: rnew(ngrid,nlayer) ! R of the atmosphere |
---|
69 | REAL, intent(in):: pp(ngrid,nlayer) ! Pressure at full levels(Pa) |
---|
70 | REAL, intent(in):: pt(ngrid,nlayer) ! Temp at full levels(K) |
---|
71 | REAL, intent(in):: pu(ngrid,nlayer) ! Zonal wind at full levels(m/s) |
---|
72 | REAL, intent(in):: pv(ngrid,nlayer) ! Meridional winds at full levels(m/s) |
---|
73 | REAL,INTENT(in) :: pdt(ngrid,nlayer) ! Tendency on temperature (K/s) |
---|
74 | REAL,INTENT(in) :: pdu(ngrid,nlayer) ! Tendency on zonal wind (m/s/s) |
---|
75 | REAL,INTENT(in) :: pdv(ngrid,nlayer) ! Tendency on meridional wind (m/s/s) |
---|
76 | |
---|
77 | ! 0.2 OUTPUTS |
---|
78 | REAL, intent(out):: zustr(ngrid) ! Zonal surface stress |
---|
79 | REAL, intent(out):: zvstr(ngrid) ! Meridional surface stress |
---|
80 | REAL, intent(out):: d_t(ngrid, nlayer) ! Tendency on temperature (K/s) due to gravity waves (not used set to zero) |
---|
81 | REAL, intent(out):: d_u(ngrid, nlayer) ! Tendency on zonal wind (m/s/s) due to gravity waves |
---|
82 | REAL, intent(out):: d_v(ngrid, nlayer) ! Tendency on meridional wind (m/s/s) due to gravity waves |
---|
83 | |
---|
84 | ! 0.3 INTERNAL ARRAYS |
---|
85 | REAL :: TT(ngrid, nlayer) ! Temperature at full levels |
---|
86 | REAL :: RHO(ngrid, nlayer) ! Mass density at full levels |
---|
87 | REAL :: UU(ngrid, nlayer) ! Zonal wind at full levels |
---|
88 | REAL :: VV(ngrid, nlayer) ! Meridional winds at full levels |
---|
89 | REAL :: BVLOW(ngrid) ! initial N at each grid (not used) |
---|
90 | |
---|
91 | INTEGER II, JJ, LL |
---|
92 | |
---|
93 | ! 0.3.0 TIME SCALE OF THE LIFE CYCLE OF THE WAVES PARAMETERIZED |
---|
94 | REAL, parameter:: DELTAT = 24. * 3600. |
---|
95 | |
---|
96 | ! 0.3.1 GRAVITY-WAVES SPECIFICATIONS |
---|
97 | INTEGER, PARAMETER:: NK = 2 ! number of horizontal wavenumbers |
---|
98 | INTEGER, PARAMETER:: NP = 2 ! directions (eastward and westward) phase speed |
---|
99 | INTEGER, PARAMETER:: NO = 2 ! absolute values of phase speed |
---|
100 | INTEGER, PARAMETER:: NA = 5 ! number of realizations to get the phase speed |
---|
101 | INTEGER, PARAMETER:: NW = NK * NP * NO ! Total numbers of gravity waves |
---|
102 | INTEGER JK, JP, JO, JW ! Loop index for NK,NP,NO, and NW |
---|
103 | |
---|
104 | REAL, save :: kmax ! Max horizontal wavenumber (lambda min,lambda=2*pi/kmax),kmax=N/u, u(=30~50) zonal wind velocity at launch altitude |
---|
105 | !$OMP THREADPRIVATE(kmax) |
---|
106 | REAL, parameter:: kmin = 2.e-5 ! Min horizontal wavenumber (lambda max = 314 km,lambda=2*pi/kmin) |
---|
107 | REAL kstar ! Min horizontal wavenumber constrained by horizontal resolution |
---|
108 | REAL, save :: cmax ! Max horizontal absolute phase velocity (at the maginitide of zonal wind u at the launch altitude) |
---|
109 | !$OMP THREADPRIVATE(cmax) |
---|
110 | REAL, parameter:: cmin = 1. ! Min horizontal absolute phase velocity (not used) |
---|
111 | REAL CPHA ! absolute PHASE VELOCITY frequency |
---|
112 | REAL ZK(NW, ngrid) ! Horizontal wavenumber amplitude |
---|
113 | REAL ZP(NW, ngrid) ! Horizontal wavenumber angle |
---|
114 | REAL ZO(NW, ngrid) ! Absolute frequency |
---|
115 | |
---|
116 | REAL intr_freq_m(nw, ngrid) ! Waves Intr. freq. at the 1/2 lev below the full level (previous name: ZOM) |
---|
117 | REAL intr_freq_p(nw, ngrid) ! Waves Intr. freq. at the 1/2 lev above the full level (previous name: ZOP) |
---|
118 | REAL wwm(nw, ngrid) ! Wave EP-fluxes at the 1/2 level below the full level |
---|
119 | REAL wwp(nw, ngrid) ! Wave EP-fluxes at the 1/2 level above the full level |
---|
120 | REAL u_epflux_p(nw, ngrid) ! Partial zonal flux (=for each wave) at the 1/2 level above the full level (previous name: RUWP) |
---|
121 | REAL v_epflux_p(nw, ngrid) ! Partial meridional flux (=for each wave) at the 1/2 level above the full level (previous name: RVWP) |
---|
122 | REAL u_epflux_tot(ngrid, nlayer + 1) ! Total zonal flux (=for all waves (nw)) at the 1/2 level above the full level (3D) (previous name: RUW) |
---|
123 | REAL v_epflux_tot(ngrid, nlayer + 1) ! Total meridional flux (=for all waves (nw)) at the 1/2 level above the full level (3D) (previous name: RVW) |
---|
124 | REAL epflux_0(nw, ngrid) ! Fluxes at launching level (previous name: RUW0) |
---|
125 | REAL, save :: epflux_max ! Max EP flux value at launching altitude (previous name: RUWMAX, tunable) |
---|
126 | !$OMP THREADPRIVATE(epflux_max) |
---|
127 | INTEGER LAUNCH ! Launching altitude |
---|
128 | REAL, parameter:: xlaunch = 0.4 ! Control the launching altitude by pressure |
---|
129 | REAL, parameter:: zmaxth_top = 8000. ! Top of convective layer (approx. not used) |
---|
130 | |
---|
131 | |
---|
132 | ! 0.3.2 PARAMETERS OF WAVES DISSIPATIONS |
---|
133 | REAL, save :: sat ! saturation parameter(tunable) |
---|
134 | !$OMP THREADPRIVATE(sat) |
---|
135 | REAL, save :: rdiss ! dissipation coefficient (tunable) |
---|
136 | !$OMP THREADPRIVATE(rdiss) |
---|
137 | REAL, parameter:: zoisec = 1.e-10 ! security for intrinsic frequency |
---|
138 | |
---|
139 | ! 0.3.3 Background flow at 1/2 levels and vertical coordinate |
---|
140 | REAL H0bis(ngrid, nlayer) ! Characteristic Height of the atmosphere (specific locations) |
---|
141 | REAL, save:: H0 ! Characteristic Height of the atmosphere (constant) |
---|
142 | !$OMP THREADPRIVATE(H0) |
---|
143 | REAL, parameter:: pr = 250 ! Reference pressure [Pa] |
---|
144 | REAL, parameter:: tr = 220. ! Reference temperature [K] |
---|
145 | REAL ZH(ngrid, nlayer + 1) ! Log-pressure altitude (constant H0) |
---|
146 | REAL ZHbis(ngrid, nlayer + 1) ! Log-pressure altitude (varying H0bis) |
---|
147 | REAL UH(ngrid, nlayer + 1) ! zonal wind at 1/2 levels |
---|
148 | REAL VH(ngrid, nlayer + 1) ! meridional wind at 1/2 levels |
---|
149 | REAL PH(ngrid, nlayer + 1) ! Pressure at 1/2 levels |
---|
150 | REAL, parameter:: psec = 1.e-20 ! Security to avoid division by 0 pressure(!!IMPORTANT: should be lower than the topmost layer's pressure) |
---|
151 | REAL BV(ngrid, nlayer + 1) ! Brunt Vaisala freq. (N^2) at 1/2 levels |
---|
152 | REAL, parameter:: bvsec = 1.e-5 ! Security to avoid negative BV |
---|
153 | REAL HREF(nlayer + 1) ! Reference pressure for launching altitude |
---|
154 | |
---|
155 | |
---|
156 | REAL RAN_NUM_1,RAN_NUM_2,RAN_NUM_3 |
---|
157 | |
---|
158 | |
---|
159 | LOGICAL,SAVE :: firstcall = .true. |
---|
160 | !$OMP THREADPRIVATE(firstcall) |
---|
161 | |
---|
162 | !----------------------------------------------------------------------------------------------------------------------- |
---|
163 | ! 1. INITIALISATIONS |
---|
164 | !----------------------------------------------------------------------------------------------------------------------- |
---|
165 | IF (firstcall) THEN |
---|
166 | write(*,*) "nonoro_gwd_ran: FLott non-oro GW scheme is active!" |
---|
167 | epflux_max = 7.E-7 ! Mars' value !! |
---|
168 | call getin_p("nonoro_gwd_epflux_max", epflux_max) |
---|
169 | write(*,*) "nonoro_gwd_ran: epflux_max=", epflux_max |
---|
170 | sat = 1. ! default gravity waves saturation value !! |
---|
171 | call getin_p("nonoro_gwd_sat", sat) |
---|
172 | write(*,*) "nonoro_gwd_ran: sat=", sat |
---|
173 | cmax = 50. ! default gravity waves phase velocity value !! |
---|
174 | call getin_p("nonoro_gwd_cmax", cmax) |
---|
175 | write(*,*) "nonoro_gwd_ran: cmax=", cmax |
---|
176 | rdiss=0.01 ! default coefficient of dissipation !! |
---|
177 | call getin_p("nonoro_gwd_rdiss", rdiss) |
---|
178 | write(*,*) "nonoro_gwd_ran: rdiss=", rdiss |
---|
179 | kmax=7.E-4 ! default Max horizontal wavenumber !! |
---|
180 | call getin_p("nonoro_gwd_kmax", kmax) |
---|
181 | write(*,*) "nonoro_gwd_ran: kmax=", kmax |
---|
182 | |
---|
183 | ! Characteristic vertical scale height |
---|
184 | H0 = r * tr / g |
---|
185 | ! Control |
---|
186 | if (deltat .LT. dtime) THEN |
---|
187 | call abort_physic("nonoro_gwd_ran","gwd random: deltat lower than dtime!",1) |
---|
188 | endif |
---|
189 | if (nlayer .LT. nw) THEN |
---|
190 | call abort_physic("nonoro_gwd_ran","gwd random: nlayer lower than nw!",1) |
---|
191 | endif |
---|
192 | firstcall = .false. |
---|
193 | ENDIF |
---|
194 | |
---|
195 | ! Compute current values of temperature and winds |
---|
196 | tt(:,:)=pt(:,:)+dtime*pdt(:,:) |
---|
197 | uu(:,:)=pu(:,:)+dtime*pdu(:,:) |
---|
198 | vv(:,:)=pv(:,:)+dtime*pdv(:,:) |
---|
199 | ! Compute the real mass density by rho=p/R(T)T |
---|
200 | DO ll=1,nlayer |
---|
201 | DO ii=1,ngrid |
---|
202 | rho(ii,ll) = pp(ii,ll)/(rnew(ii,ll)*tt(ii,ll)) |
---|
203 | ENDDO |
---|
204 | ENDDO |
---|
205 | ! print*,'epflux_max just after firstcall:', epflux_max |
---|
206 | |
---|
207 | !----------------------------------------------------------------------------------------------------------------------- |
---|
208 | ! 2. EVALUATION OF THE BACKGROUND FLOW AT SEMI-LEVELS |
---|
209 | !----------------------------------------------------------------------------------------------------------------------- |
---|
210 | ! Pressure and inverse of pressure at 1/2 level |
---|
211 | DO LL = 2, nlayer |
---|
212 | PH(:, LL) = EXP((LOG(PP(:, LL)) + LOG(PP(:, LL - 1))) / 2.) |
---|
213 | end DO |
---|
214 | PH(:, nlayer + 1) = 0. |
---|
215 | PH(:, 1) = 2. * PP(:, 1) - PH(:, 2) |
---|
216 | |
---|
217 | call writediagfi(ngrid,'nonoro_pp','nonoro_pp', 'm',3,PP(:,1:nlayer)) |
---|
218 | call writediagfi(ngrid,'nonoro_ph','nonoro_ph', 'm',3,PH(:,1:nlayer)) |
---|
219 | |
---|
220 | ! Launching level for reproductible case |
---|
221 | !Pour revenir a la version non reproductible en changeant le nombre de |
---|
222 | !process |
---|
223 | ! Reprend la formule qui calcule PH en fonction de PP=play |
---|
224 | DO LL = 2, nlayer |
---|
225 | HREF(LL) = EXP((LOG(presnivs(LL))+ LOG(presnivs(LL - 1))) / 2.) |
---|
226 | end DO |
---|
227 | HREF(nlayer + 1) = 0. |
---|
228 | HREF(1) = 2. * presnivs(1) - HREF(2) |
---|
229 | |
---|
230 | LAUNCH=0 |
---|
231 | DO LL = 1, nlayer |
---|
232 | IF (HREF(LL) / HREF(1) > XLAUNCH) LAUNCH = LL |
---|
233 | ENDDO |
---|
234 | |
---|
235 | ! Log pressure vert. coordinate |
---|
236 | ZH(:,1) = 0. |
---|
237 | DO LL = 2, nlayer + 1 |
---|
238 | !ZH(:, LL) = H0 * LOG(PR / (PH(:, LL) + PSEC)) |
---|
239 | H0bis(:, LL-1) = r * TT(:, LL-1) / g |
---|
240 | ZH(:, LL) = ZH(:, LL-1) & |
---|
241 | + H0bis(:, LL-1)*(PH(:, LL-1)-PH(:,LL))/PP(:, LL-1) |
---|
242 | end DO |
---|
243 | ZH(:, 1) = H0 * LOG(PR / (PH(:, 1) + PSEC)) |
---|
244 | |
---|
245 | call writediagfi(ngrid,'nonoro_zh','nonoro_zh', 'm',3,ZH(:,2:nlayer+1)) |
---|
246 | |
---|
247 | ! Winds and Brunt Vaisala frequency |
---|
248 | DO LL = 2, nlayer |
---|
249 | UH(:, LL) = 0.5 * (UU(:, LL) + UU(:, LL - 1)) ! Zonal wind |
---|
250 | VH(:, LL) = 0.5 * (VV(:, LL) + VV(:, LL - 1)) ! Meridional wind |
---|
251 | ! Brunt Vaisala frequency (=g/T*[dT/dz + g/cp] ) |
---|
252 | BV(:, LL)= G/(0.5 * (TT(:, LL) + TT(:, LL - 1))) & |
---|
253 | *((TT(:, LL) - TT(:, LL - 1)) / (ZH(:, LL) - ZH(:, LL - 1))+ & |
---|
254 | G / cpnew(:,LL)) |
---|
255 | |
---|
256 | BV(:,LL) =MAX(1.E-12,BV(:,LL)) ! to ensure that BV is positive |
---|
257 | BV(:,LL) = MAX(BVSEC,SQRT(BV(:,LL))) ! make sure it is not too small |
---|
258 | end DO |
---|
259 | BV(:, 1) = BV(:, 2) |
---|
260 | UH(:, 1) = 0. |
---|
261 | VH(:, 1) = 0. |
---|
262 | BV(:, nlayer + 1) = BV(:, nlayer) |
---|
263 | UH(:, nlayer + 1) = UU(:, nlayer) |
---|
264 | VH(:, nlayer + 1) = VV(:, nlayer) |
---|
265 | |
---|
266 | call writediagfi(ngrid,'nonoro_bv','nonoro_bv', 'm',3,BV(:,2:nlayer+1)) |
---|
267 | |
---|
268 | !----------------------------------------------------------------------------------------------------------------------- |
---|
269 | ! 3. WAVES CHARACTERISTICS CHOSEN RANDOMLY |
---|
270 | !----------------------------------------------------------------------------------------------------------------------- |
---|
271 | ! The mod function of here a weird arguments are used to produce the waves characteristics in a stochastic way |
---|
272 | DO JW = 1, NW |
---|
273 | ! Angle |
---|
274 | DO II = 1, ngrid |
---|
275 | ! Angle (0 or PI so far) |
---|
276 | RAN_NUM_1=MOD(TT(II, JW) * 10., 1.) |
---|
277 | RAN_NUM_2= MOD(TT(II, JW) * 100., 1.) |
---|
278 | ZP(JW, II) = (SIGN(1., 0.5 - RAN_NUM_1) + 1.)* PI / 2. |
---|
279 | |
---|
280 | ! Horizontal wavenumber amplitude |
---|
281 | ! From Venus model: TN+GG April/June 2020 (rev2381) |
---|
282 | ! "Individual waves are not supposed to occupy the entire domain, but only a fraction of it" (Lott+2012) |
---|
283 | ! ZK(JW, II) = KMIN + (KMAX - KMIN) *RAN_NUM_2 |
---|
284 | KSTAR = PI/SQRT(cell_area(II)) |
---|
285 | ZK(JW, II) = MAX(KMIN,KSTAR) + (KMAX - MAX(KMIN,KSTAR)) *RAN_NUM_2 |
---|
286 | |
---|
287 | ! Horizontal phase speed |
---|
288 | ! this computation allows to get a gaussian distribution centered on 0 (right ?) |
---|
289 | ! then cmin is not useful, and it favors small values. |
---|
290 | CPHA = 0. |
---|
291 | DO JJ = 1, NA |
---|
292 | RAN_NUM_3=MOD(TT(II, JW+3*JJ)**2, 1.) |
---|
293 | CPHA = CPHA + 2.*CMAX* & |
---|
294 | (RAN_NUM_3 -0.5)* & |
---|
295 | SQRT(3.)/SQRT(NA*1.) |
---|
296 | END DO |
---|
297 | IF (CPHA.LT.0.) THEN |
---|
298 | CPHA = -1.*CPHA |
---|
299 | ZP(JW,II) = ZP(JW,II) + PI |
---|
300 | ENDIF |
---|
301 | ! otherwise, with the computation below, we get a uniform distribution between cmin and cmax. |
---|
302 | ! ran_num_3 = mod(tt(ii, jw)**2, 1.) |
---|
303 | ! cpha = cmin + (cmax - cmin) * ran_num_3 |
---|
304 | |
---|
305 | ! Intrinsic frequency |
---|
306 | ZO(JW, II) = CPHA * ZK(JW, II) |
---|
307 | ! Intrinsic frequency is imposed |
---|
308 | ZO(JW, II) = ZO(JW, II) & |
---|
309 | + ZK(JW, II) * COS(ZP(JW, II)) * UH(II, LAUNCH) & |
---|
310 | + ZK(JW, II) * SIN(ZP(JW, II)) * VH(II, LAUNCH) |
---|
311 | |
---|
312 | ! Momentum flux at launch level |
---|
313 | epflux_0(JW, II) = epflux_max & |
---|
314 | * MOD(100. * (UU(II, JW)**2 + VV(II, JW)**2), 1.) |
---|
315 | ENDDO |
---|
316 | end DO |
---|
317 | ! print*,'epflux_0 just after waves charac. ramdon:', epflux_0 |
---|
318 | |
---|
319 | !----------------------------------------------------------------------------------------------------------------------- |
---|
320 | ! 4. COMPUTE THE FLUXES |
---|
321 | !----------------------------------------------------------------------------------------------------------------------- |
---|
322 | ! 4.1 Vertical velocity at launching altitude to ensure the correct value to the imposed fluxes. |
---|
323 | !------------------------------------------------------ |
---|
324 | DO JW = 1, NW |
---|
325 | ! Evaluate intrinsic frequency at launching altitude: |
---|
326 | intr_freq_p(JW, :) = ZO(JW, :) & |
---|
327 | - ZK(JW, :) * COS(ZP(JW, :)) * UH(:, LAUNCH) & |
---|
328 | - ZK(JW, :) * SIN(ZP(JW, :)) * VH(:, LAUNCH) |
---|
329 | end DO |
---|
330 | |
---|
331 | ! VERSION WITHOUT CONVECTIVE SOURCE |
---|
332 | ! Vertical velocity at launch level, value to ensure the imposed |
---|
333 | ! mom flux: |
---|
334 | DO JW = 1, NW |
---|
335 | ! WW is directly a flux, here, not vertical velocity anymore |
---|
336 | WWP(JW, :) = epflux_0(JW,:) |
---|
337 | u_epflux_p(JW, :) = COS(ZP(JW, :)) * SIGN(1., intr_freq_p(JW, :)) * epflux_0(JW, :) |
---|
338 | v_epflux_p(JW, :) = SIN(ZP(JW, :)) * SIGN(1., intr_freq_p(JW, :)) * epflux_0(JW, :) |
---|
339 | end DO |
---|
340 | |
---|
341 | ! 4.2 Initial flux at launching altitude |
---|
342 | !------------------------------------------------------ |
---|
343 | u_epflux_tot(:, LAUNCH) = 0 |
---|
344 | v_epflux_tot(:, LAUNCH) = 0 |
---|
345 | DO JW = 1, NW |
---|
346 | u_epflux_tot(:, LAUNCH) = u_epflux_tot(:, LAUNCH) + u_epflux_p(JW, :) |
---|
347 | v_epflux_tot(:, LAUNCH) = v_epflux_tot(:, LAUNCH) + v_epflux_p(JW, :) |
---|
348 | end DO |
---|
349 | |
---|
350 | ! 4.3 Loop over altitudes, with passage from one level to the next done by: |
---|
351 | !---------------------------------------------------------------------------- |
---|
352 | ! i) conserving the EP flux, |
---|
353 | ! ii) dissipating a little, |
---|
354 | ! iii) testing critical levels, |
---|
355 | ! iv) testing the breaking. |
---|
356 | !---------------------------------------------------------------------------- |
---|
357 | DO LL = LAUNCH, nlayer - 1 |
---|
358 | ! W(KB)ARNING: ALL THE PHYSICS IS HERE (PASSAGE FROM ONE LEVEL TO THE NEXT) |
---|
359 | DO JW = 1, NW |
---|
360 | intr_freq_m(JW, :) = intr_freq_p(JW, :) |
---|
361 | WWM(JW, :) = WWP(JW, :) |
---|
362 | ! Intrinsic Frequency |
---|
363 | intr_freq_p(JW, :) = ZO(JW, :) - ZK(JW, :) * COS(ZP(JW,:)) * UH(:, LL + 1) & |
---|
364 | - ZK(JW, :) * SIN(ZP(JW,:)) * VH(:, LL + 1) |
---|
365 | ! Vertical velocity in flux formulation |
---|
366 | WWP(JW, :) = MIN( & |
---|
367 | ! No breaking (Eq.6) |
---|
368 | WWM(JW, :) & |
---|
369 | ! Dissipation (Eq. 8)(real rho used here rather than pressure |
---|
370 | ! parametration because the original code has a bug if the density of |
---|
371 | ! the planet at the launch altitude not approximate 1): ! |
---|
372 | * EXP(- RDISS*2./rho(:, LL + 1) & |
---|
373 | * ((BV(:, LL + 1) + BV(:, LL)) / 2.)**3 & |
---|
374 | / MAX(ABS(intr_freq_p(JW, :) + intr_freq_m(JW, :)) / 2., ZOISEC)**4 & |
---|
375 | * ZK(JW, :)**3 * (ZH(:, LL + 1) - ZH(:, LL))), & |
---|
376 | ! Critical levels (forced to zero if intrinsic frequency |
---|
377 | ! changes sign) |
---|
378 | MAX(0., SIGN(1., intr_freq_p(JW, :) * intr_freq_m(JW, :))) & |
---|
379 | ! Saturation (Eq. 12) (rho at launch altitude is imposed by J.Liu. |
---|
380 | ! Same reason with Eq. 8) |
---|
381 | * ABS(intr_freq_p(JW, :))**3 / (2.*BV(:, LL+1)) & |
---|
382 | * rho(:,launch)*exp(-zh(:, ll + 1) / H0) & |
---|
383 | * SAT**2 *KMIN**2 / ZK(JW, :)**4) |
---|
384 | end DO |
---|
385 | ! END OF W(KB)ARNING |
---|
386 | |
---|
387 | ! Evaluate EP-flux from Eq. 7 and give the right orientation to the stress |
---|
388 | DO JW = 1, NW |
---|
389 | u_epflux_p(JW, :) = SIGN(1.,intr_freq_p(JW, :)) * COS(ZP(JW, :)) * WWP(JW, :) |
---|
390 | v_epflux_p(JW, :) = SIGN(1.,intr_freq_p(JW, :)) * SIN(ZP(JW, :)) * WWP(JW, :) |
---|
391 | end DO |
---|
392 | |
---|
393 | u_epflux_tot(:, LL + 1) = 0. |
---|
394 | v_epflux_tot(:, LL + 1) = 0. |
---|
395 | DO JW = 1, NW |
---|
396 | u_epflux_tot(:, LL + 1) = u_epflux_tot(:, LL + 1) + u_epflux_p(JW, :) |
---|
397 | v_epflux_tot(:, LL + 1) = v_epflux_tot(:, LL + 1) + v_epflux_p(JW, :) |
---|
398 | EAST_GWSTRESS(:, LL)=EAST_GWSTRESS(:, LL)+MAX(0.,u_epflux_p(JW,:))/FLOAT(NW) |
---|
399 | WEST_GWSTRESS(:, LL)=WEST_GWSTRESS(:, LL)+MIN(0.,u_epflux_p(JW,:))/FLOAT(NW) |
---|
400 | end DO |
---|
401 | end DO ! DO LL = LAUNCH, nlayer - 1 |
---|
402 | |
---|
403 | !----------------------------------------------------------------------------------------------------------------------- |
---|
404 | ! 5. TENDENCY CALCULATIONS |
---|
405 | !----------------------------------------------------------------------------------------------------------------------- |
---|
406 | |
---|
407 | ! 5.1 Flow rectification at the top and in the low layers |
---|
408 | ! -------------------------------------------------------- |
---|
409 | ! Warning, this is the total on all GW |
---|
410 | u_epflux_tot(:, nlayer + 1) = 0. |
---|
411 | v_epflux_tot(:, nlayer + 1) = 0. |
---|
412 | |
---|
413 | ! Here, big change compared to FLott version: |
---|
414 | ! We compensate (u_epflux_tot(:, LAUNCH), ie total emitted upward flux |
---|
415 | ! over the layers max(1,LAUNCH-3) to LAUNCH-1 |
---|
416 | DO LL = 1, max(1,LAUNCH-3) |
---|
417 | u_epflux_tot(:, LL) = 0. |
---|
418 | v_epflux_tot(:, LL) = 0. |
---|
419 | end DO |
---|
420 | DO LL = max(2,LAUNCH-2), LAUNCH-1 |
---|
421 | u_epflux_tot(:, LL) = u_epflux_tot(:, LL - 1) + u_epflux_tot(:, LAUNCH) * & |
---|
422 | (PH(:,LL)-PH(:,LL-1)) / (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3))) |
---|
423 | v_epflux_tot(:, LL) = v_epflux_tot(:, LL - 1) + v_epflux_tot(:, LAUNCH) * & |
---|
424 | (PH(:,LL)-PH(:,LL-1)) / (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3))) |
---|
425 | EAST_GWSTRESS(:,LL) = EAST_GWSTRESS(:, LL - 1) + & |
---|
426 | EAST_GWSTRESS(:, LAUNCH) * (PH(:,LL)-PH(:,LL-1))/ & |
---|
427 | (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3))) |
---|
428 | WEST_GWSTRESS(:,LL) = WEST_GWSTRESS(:, LL - 1) + & |
---|
429 | WEST_GWSTRESS(:, LAUNCH) * (PH(:,LL)-PH(:,LL-1))/ & |
---|
430 | (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3))) |
---|
431 | end DO |
---|
432 | ! This way, the total flux from GW is zero, but there is a net transport |
---|
433 | ! (upward) that should be compensated by circulation |
---|
434 | ! and induce additional friction at the surface |
---|
435 | call writediagfi(ngrid,'nonoro_u_epflux_tot','nonoro_u_epflux_tot', '',3,u_epflux_tot(:,2:nlayer+1)) |
---|
436 | call writediagfi(ngrid,'nonoro_v_epflux_tot','nonoro_v_epflux_tot', '',3,v_epflux_tot(:,2:nlayer+1)) |
---|
437 | |
---|
438 | ! 5.2 AR-1 RECURSIVE FORMULA (13) IN VERSION 4 |
---|
439 | !--------------------------------------------- |
---|
440 | DO LL = 1, nlayer |
---|
441 | !Notice here the d_u and d_v are tendency (For Mars) but not increment(Venus). |
---|
442 | d_u(:, LL) = G * (u_epflux_tot(:, LL + 1) - u_epflux_tot(:, LL)) & |
---|
443 | / (PH(:, LL + 1) - PH(:, LL)) |
---|
444 | d_v(:, LL) = G * (v_epflux_tot(:, LL + 1) - v_epflux_tot(:, LL)) & |
---|
445 | / (PH(:, LL + 1) - PH(:, LL)) |
---|
446 | ENDDO |
---|
447 | d_t(:,:) = 0. |
---|
448 | call writediagfi(ngrid,'nonoro_d_u','nonoro_d_u', '',3,d_u) |
---|
449 | call writediagfi(ngrid,'nonoro_d_v','nonoro_d_v', '',3,d_v) |
---|
450 | |
---|
451 | ! 5.3 Update tendency of wind with the previous (and saved) values |
---|
452 | !----------------------------------------------------------------- |
---|
453 | d_u(:,:) = DTIME/DELTAT/REAL(NW) * d_u(:,:) & |
---|
454 | + (1.-DTIME/DELTAT) * du_nonoro_gwd(:,:) |
---|
455 | d_v(:,:) = DTIME/DELTAT/REAL(NW) * d_v(:,:) & |
---|
456 | + (1.-DTIME/DELTAT) * dv_nonoro_gwd(:,:) |
---|
457 | du_nonoro_gwd(:,:) = d_u(:,:) |
---|
458 | dv_nonoro_gwd(:,:) = d_v(:,:) |
---|
459 | |
---|
460 | call writediagfi(ngrid,'du_nonoro_gwd','du_nonoro_gwd', '',3,du_nonoro_gwd) |
---|
461 | call writediagfi(ngrid,'dv_nonoro_gwd','dv_nonoro_gwd', '',3,dv_nonoro_gwd) |
---|
462 | |
---|
463 | ! Cosmetic: evaluation of the cumulated stress |
---|
464 | ZUSTR(:) = 0. |
---|
465 | ZVSTR(:) = 0. |
---|
466 | DO LL = 1, nlayer |
---|
467 | ZUSTR(:) = ZUSTR(:) + D_U(:, LL) *rho(:, LL) * (ZH(:, LL + 1) - ZH(:, LL))*DTIME |
---|
468 | ZVSTR(:) = ZVSTR(:) + D_V(:, LL) *rho(:, LL) * (ZH(:, LL + 1) - ZH(:, LL))*DTIME |
---|
469 | ENDDO |
---|
470 | |
---|
471 | !#ifdef CPP_XIOS |
---|
472 | ! call send_xios_field("du_nonoro", d_u) |
---|
473 | ! call send_xios_field("dv_nonoro", d_v) |
---|
474 | ! call send_xios_field("east_gwstress", east_gwstress) |
---|
475 | ! call send_xios_field("west_gwstress", west_gwstress) |
---|
476 | !#endif |
---|
477 | |
---|
478 | END SUBROUTINE NONORO_GWD_RAN |
---|
479 | |
---|
480 | |
---|
481 | |
---|
482 | ! ======================================================== |
---|
483 | ! Subroutines used to allocate/deallocate module variables |
---|
484 | ! ======================================================== |
---|
485 | SUBROUTINE ini_nonoro_gwd_ran(ngrid,nlayer) |
---|
486 | |
---|
487 | IMPLICIT NONE |
---|
488 | |
---|
489 | INTEGER, INTENT (in) :: ngrid ! number of atmospheric columns |
---|
490 | INTEGER, INTENT (in) :: nlayer ! number of atmospheric layers |
---|
491 | |
---|
492 | allocate(du_nonoro_gwd(ngrid,nlayer)) |
---|
493 | allocate(dv_nonoro_gwd(ngrid,nlayer)) |
---|
494 | allocate(east_gwstress(ngrid,nlayer)) |
---|
495 | east_gwstress(:,:)=0 |
---|
496 | allocate(west_gwstress(ngrid,nlayer)) |
---|
497 | west_gwstress(:,:)=0 |
---|
498 | |
---|
499 | END SUBROUTINE ini_nonoro_gwd_ran |
---|
500 | ! ---------------------------------- |
---|
501 | SUBROUTINE end_nonoro_gwd_ran |
---|
502 | |
---|
503 | IMPLICIT NONE |
---|
504 | |
---|
505 | if (allocated(du_nonoro_gwd)) deallocate(du_nonoro_gwd) |
---|
506 | if (allocated(dv_nonoro_gwd)) deallocate(dv_nonoro_gwd) |
---|
507 | if (allocated(east_gwstress)) deallocate(east_gwstress) |
---|
508 | if (allocated(west_gwstress)) deallocate(west_gwstress) |
---|
509 | |
---|
510 | END SUBROUTINE end_nonoro_gwd_ran |
---|
511 | |
---|
512 | END MODULE nonoro_gwd_ran_mod |
---|