1 | !> @file inflow_turbulence.f90 |
---|
2 | !------------------------------------------------------------------------------! |
---|
3 | ! This file is part of the PALM model system. |
---|
4 | ! |
---|
5 | ! PALM is free software: you can redistribute it and/or modify it under the |
---|
6 | ! terms of the GNU General Public License as published by the Free Software |
---|
7 | ! Foundation, either version 3 of the License, or (at your option) any later |
---|
8 | ! version. |
---|
9 | ! |
---|
10 | ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY |
---|
11 | ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR |
---|
12 | ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. |
---|
13 | ! |
---|
14 | ! You should have received a copy of the GNU General Public License along with |
---|
15 | ! PALM. If not, see <http://www.gnu.org/licenses/>. |
---|
16 | ! |
---|
17 | ! Copyright 1997-2017 Leibniz Universitaet Hannover |
---|
18 | !------------------------------------------------------------------------------! |
---|
19 | ! |
---|
20 | ! Current revisions: |
---|
21 | ! ----------------- |
---|
22 | ! |
---|
23 | ! |
---|
24 | ! Former revisions: |
---|
25 | ! ----------------- |
---|
26 | ! $Id: inflow_turbulence.f90 2696 2017-12-14 17:12:51Z maronga $ |
---|
27 | ! |
---|
28 | ! 2000 2016-08-20 18:09:15Z knoop |
---|
29 | ! Forced header and separation lines into 80 columns |
---|
30 | ! |
---|
31 | ! 1960 2016-07-12 16:34:24Z suehring |
---|
32 | ! Separate humidity and passive scalar |
---|
33 | ! |
---|
34 | ! 1806 2016-04-05 18:55:35Z gronemeier |
---|
35 | ! Added comments to variables and code segments. Removed code redundancies. |
---|
36 | ! |
---|
37 | ! 1682 2015-10-07 23:56:08Z knoop |
---|
38 | ! Code annotations made doxygen readable |
---|
39 | ! |
---|
40 | ! 1615 2015-07-08 18:49:19Z suehring |
---|
41 | ! Enable turbulent inflow for passive_scalar and humidity |
---|
42 | ! |
---|
43 | ! 1560 2015-03-06 10:48:54Z keck |
---|
44 | ! Option recycling_yshift added. If this option is switched on, the turbulence |
---|
45 | ! data, which is mapped from the recycling plane to the inflow, is shifted in |
---|
46 | ! y direction (by ny * dy / 2 ) |
---|
47 | ! |
---|
48 | ! 1353 2014-04-08 15:21:23Z heinze |
---|
49 | ! REAL constants provided with KIND-attribute |
---|
50 | ! |
---|
51 | ! 1346 2014-03-27 13:18:20Z heinze |
---|
52 | ! Bugfix: REAL constants provided with KIND-attribute especially in call of |
---|
53 | ! intrinsic function like MAX, MIN, SIGN |
---|
54 | ! |
---|
55 | ! 1320 2014-03-20 08:40:49Z raasch |
---|
56 | ! ONLY-attribute added to USE-statements, |
---|
57 | ! kind-parameters added to all INTEGER and REAL declaration statements, |
---|
58 | ! kinds are defined in new module kinds, |
---|
59 | ! revision history before 2012 removed, |
---|
60 | ! comment fields (!:) to be used for variable explanations added to |
---|
61 | ! all variable declaration statements |
---|
62 | ! |
---|
63 | ! 1092 2013-02-02 11:24:22Z raasch |
---|
64 | ! unused variables removed |
---|
65 | ! |
---|
66 | ! 1036 2012-10-22 13:43:42Z raasch |
---|
67 | ! code put under GPL (PALM 3.9) |
---|
68 | ! |
---|
69 | ! Initial version (2008/03/07) |
---|
70 | ! |
---|
71 | ! Description: |
---|
72 | ! ------------ |
---|
73 | !> Imposing turbulence at the respective inflow using the turbulence |
---|
74 | !> recycling method of Kataoka and Mizuno (2002). |
---|
75 | !------------------------------------------------------------------------------! |
---|
76 | SUBROUTINE inflow_turbulence |
---|
77 | |
---|
78 | |
---|
79 | USE arrays_3d, & |
---|
80 | ONLY: e, inflow_damping_factor, mean_inflow_profiles, pt, q, s, u, v, w |
---|
81 | |
---|
82 | USE control_parameters, & |
---|
83 | ONLY: humidity, passive_scalar, recycling_plane, recycling_yshift |
---|
84 | |
---|
85 | USE cpulog, & |
---|
86 | ONLY: cpu_log, log_point |
---|
87 | |
---|
88 | USE indices, & |
---|
89 | ONLY: nbgp, nxl, ny, nyn, nys, nyng, nysg, nzb, nzt |
---|
90 | |
---|
91 | USE kinds |
---|
92 | |
---|
93 | USE pegrid |
---|
94 | |
---|
95 | |
---|
96 | IMPLICIT NONE |
---|
97 | |
---|
98 | INTEGER(iwp) :: i !< loop index |
---|
99 | INTEGER(iwp) :: j !< loop index |
---|
100 | INTEGER(iwp) :: k !< loop index |
---|
101 | INTEGER(iwp) :: l !< loop index |
---|
102 | INTEGER(iwp) :: next !< ID of receiving PE for y-shift |
---|
103 | INTEGER(iwp) :: ngp_ifd !< number of grid points stored in avpr |
---|
104 | INTEGER(iwp) :: ngp_pr !< number of grid points stored in inflow_dist |
---|
105 | INTEGER(iwp) :: prev !< ID of sending PE for y-shift |
---|
106 | |
---|
107 | REAL(wp), DIMENSION(nzb:nzt+1,7,nbgp) :: & |
---|
108 | avpr !< stores averaged profiles at recycling plane |
---|
109 | REAL(wp), DIMENSION(nzb:nzt+1,7,nbgp) :: & |
---|
110 | avpr_l !< auxiliary variable to calculate avpr |
---|
111 | REAL(wp), DIMENSION(nzb:nzt+1,nysg:nyng,7,nbgp) :: & |
---|
112 | inflow_dist !< turbulence signal of vars, added at inflow boundary |
---|
113 | REAL(wp), DIMENSION(nzb:nzt+1,nysg:nyng,7,nbgp) :: & |
---|
114 | local_inflow_dist !< auxiliary variable for inflow_dist, used for yshift |
---|
115 | |
---|
116 | CALL cpu_log( log_point(40), 'inflow_turbulence', 'start' ) |
---|
117 | |
---|
118 | ! |
---|
119 | !-- Carry out spanwise averaging in the recycling plane |
---|
120 | avpr_l = 0.0_wp |
---|
121 | ngp_pr = ( nzt - nzb + 2 ) * 7 * nbgp |
---|
122 | ngp_ifd = ngp_pr * ( nyn - nys + 1 + 2 * nbgp ) |
---|
123 | |
---|
124 | ! |
---|
125 | !-- First, local averaging within the recycling domain |
---|
126 | i = recycling_plane |
---|
127 | |
---|
128 | #if defined( __parallel ) |
---|
129 | IF ( myidx == id_recycling ) THEN |
---|
130 | |
---|
131 | DO l = 1, nbgp |
---|
132 | DO j = nys, nyn |
---|
133 | DO k = nzb, nzt + 1 |
---|
134 | |
---|
135 | avpr_l(k,1,l) = avpr_l(k,1,l) + u(k,j,i) |
---|
136 | avpr_l(k,2,l) = avpr_l(k,2,l) + v(k,j,i) |
---|
137 | avpr_l(k,3,l) = avpr_l(k,3,l) + w(k,j,i) |
---|
138 | avpr_l(k,4,l) = avpr_l(k,4,l) + pt(k,j,i) |
---|
139 | avpr_l(k,5,l) = avpr_l(k,5,l) + e(k,j,i) |
---|
140 | IF ( humidity ) & |
---|
141 | avpr_l(k,6,l) = avpr_l(k,6,l) + q(k,j,i) |
---|
142 | IF ( passive_scalar ) & |
---|
143 | avpr_l(k,7,l) = avpr_l(k,7,l) + s(k,j,i) |
---|
144 | |
---|
145 | ENDDO |
---|
146 | ENDDO |
---|
147 | i = i + 1 |
---|
148 | ENDDO |
---|
149 | |
---|
150 | ENDIF |
---|
151 | ! |
---|
152 | !-- Now, averaging over all PEs |
---|
153 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
154 | CALL MPI_ALLREDUCE( avpr_l(nzb,1,1), avpr(nzb,1,1), ngp_pr, MPI_REAL, & |
---|
155 | MPI_SUM, comm2d, ierr ) |
---|
156 | |
---|
157 | #else |
---|
158 | DO l = 1, nbgp |
---|
159 | DO j = nys, nyn |
---|
160 | DO k = nzb, nzt + 1 |
---|
161 | |
---|
162 | avpr_l(k,1,l) = avpr_l(k,1,l) + u(k,j,i) |
---|
163 | avpr_l(k,2,l) = avpr_l(k,2,l) + v(k,j,i) |
---|
164 | avpr_l(k,3,l) = avpr_l(k,3,l) + w(k,j,i) |
---|
165 | avpr_l(k,4,l) = avpr_l(k,4,l) + pt(k,j,i) |
---|
166 | avpr_l(k,5,l) = avpr_l(k,5,l) + e(k,j,i) |
---|
167 | IF ( humidity ) & |
---|
168 | avpr_l(k,6,l) = avpr_l(k,6,l) + q(k,j,i) |
---|
169 | IF ( passive_scalar ) & |
---|
170 | avpr_l(k,7,l) = avpr_l(k,7,l) + s(k,j,i) |
---|
171 | |
---|
172 | ENDDO |
---|
173 | ENDDO |
---|
174 | i = i + 1 |
---|
175 | ENDDO |
---|
176 | |
---|
177 | avpr = avpr_l |
---|
178 | #endif |
---|
179 | |
---|
180 | avpr = avpr / ( ny + 1 ) |
---|
181 | ! |
---|
182 | !-- Calculate the disturbances at the recycling plane |
---|
183 | i = recycling_plane |
---|
184 | |
---|
185 | #if defined( __parallel ) |
---|
186 | IF ( myidx == id_recycling ) THEN |
---|
187 | DO l = 1, nbgp |
---|
188 | DO j = nysg, nyng |
---|
189 | DO k = nzb, nzt + 1 |
---|
190 | |
---|
191 | inflow_dist(k,j,1,l) = u(k,j,i+1) - avpr(k,1,l) |
---|
192 | inflow_dist(k,j,2,l) = v(k,j,i) - avpr(k,2,l) |
---|
193 | inflow_dist(k,j,3,l) = w(k,j,i) - avpr(k,3,l) |
---|
194 | inflow_dist(k,j,4,l) = pt(k,j,i) - avpr(k,4,l) |
---|
195 | inflow_dist(k,j,5,l) = e(k,j,i) - avpr(k,5,l) |
---|
196 | IF ( humidity ) & |
---|
197 | inflow_dist(k,j,6,l) = q(k,j,i) - avpr(k,6,l) |
---|
198 | IF ( passive_scalar ) & |
---|
199 | inflow_dist(k,j,7,l) = s(k,j,i) - avpr(k,7,l) |
---|
200 | ENDDO |
---|
201 | ENDDO |
---|
202 | i = i + 1 |
---|
203 | ENDDO |
---|
204 | |
---|
205 | ENDIF |
---|
206 | #else |
---|
207 | DO l = 1, nbgp |
---|
208 | DO j = nysg, nyng |
---|
209 | DO k = nzb, nzt+1 |
---|
210 | |
---|
211 | inflow_dist(k,j,1,l) = u(k,j,i+1) - avpr(k,1,l) |
---|
212 | inflow_dist(k,j,2,l) = v(k,j,i) - avpr(k,2,l) |
---|
213 | inflow_dist(k,j,3,l) = w(k,j,i) - avpr(k,3,l) |
---|
214 | inflow_dist(k,j,4,l) = pt(k,j,i) - avpr(k,4,l) |
---|
215 | inflow_dist(k,j,5,l) = e(k,j,i) - avpr(k,5,l) |
---|
216 | IF ( humidity ) & |
---|
217 | inflow_dist(k,j,6,l) = q(k,j,i) - avpr(k,6,l) |
---|
218 | IF ( passive_scalar ) & |
---|
219 | inflow_dist(k,j,7,l) = s(k,j,i) - avpr(k,7,l) |
---|
220 | |
---|
221 | ENDDO |
---|
222 | ENDDO |
---|
223 | i = i + 1 |
---|
224 | ENDDO |
---|
225 | #endif |
---|
226 | |
---|
227 | ! |
---|
228 | !-- For parallel runs, send the disturbances to the respective inflow PE |
---|
229 | #if defined( __parallel ) |
---|
230 | IF ( myidx == id_recycling .AND. myidx /= id_inflow ) THEN |
---|
231 | |
---|
232 | CALL MPI_SEND( inflow_dist(nzb,nysg,1,1), ngp_ifd, MPI_REAL, & |
---|
233 | id_inflow, 1, comm1dx, ierr ) |
---|
234 | |
---|
235 | ELSEIF ( myidx /= id_recycling .AND. myidx == id_inflow ) THEN |
---|
236 | |
---|
237 | inflow_dist = 0.0_wp |
---|
238 | CALL MPI_RECV( inflow_dist(nzb,nysg,1,1), ngp_ifd, MPI_REAL, & |
---|
239 | id_recycling, 1, comm1dx, status, ierr ) |
---|
240 | |
---|
241 | ENDIF |
---|
242 | |
---|
243 | ! |
---|
244 | !-- y-shift for inflow_dist |
---|
245 | !-- Shift inflow_dist in positive y direction by a distance of INT( npey / 2 ) |
---|
246 | IF ( recycling_yshift .AND. myidx == id_inflow ) THEN |
---|
247 | ! |
---|
248 | !-- Calculate the ID of the PE which sends data to this PE (prev) and of the |
---|
249 | !-- PE which receives data from this PE (next). |
---|
250 | IF ( myidy >= INT( pdims(2) / 2 ) ) THEN |
---|
251 | prev = myidy - INT( pdims(2) / 2 ) |
---|
252 | ELSE |
---|
253 | prev = pdims(2) - ( INT( pdims(2) / 2 ) - myidy ) |
---|
254 | ENDIF |
---|
255 | |
---|
256 | IF ( myidy < pdims(2) - INT( pdims(2) / 2 ) ) THEN |
---|
257 | next = myidy + INT( pdims(2) / 2 ) |
---|
258 | ELSE |
---|
259 | next = INT( pdims(2) / 2 ) - ( pdims(2) - myidy ) |
---|
260 | ENDIF |
---|
261 | |
---|
262 | local_inflow_dist = 0.0_wp |
---|
263 | |
---|
264 | CALL MPI_SENDRECV( inflow_dist(nzb,nysg,1,1), ngp_ifd, MPI_REAL, & |
---|
265 | next, 1, local_inflow_dist(nzb,nysg,1,1), ngp_ifd, & |
---|
266 | MPI_REAL, prev, 1, comm1dy, status, ierr ) |
---|
267 | |
---|
268 | inflow_dist = local_inflow_dist |
---|
269 | |
---|
270 | ENDIF |
---|
271 | |
---|
272 | #endif |
---|
273 | |
---|
274 | ! |
---|
275 | !-- Add the disturbance at the inflow |
---|
276 | IF ( nxl == 0 ) THEN |
---|
277 | |
---|
278 | DO j = nysg, nyng |
---|
279 | DO k = nzb, nzt + 1 |
---|
280 | |
---|
281 | u(k,j,-nbgp+1:0) = mean_inflow_profiles(k,1) + & |
---|
282 | inflow_dist(k,j,1,1:nbgp) * inflow_damping_factor(k) |
---|
283 | v(k,j,-nbgp:-1) = mean_inflow_profiles(k,2) + & |
---|
284 | inflow_dist(k,j,2,1:nbgp) * inflow_damping_factor(k) |
---|
285 | w(k,j,-nbgp:-1) = & |
---|
286 | inflow_dist(k,j,3,1:nbgp) * inflow_damping_factor(k) |
---|
287 | pt(k,j,-nbgp:-1) = mean_inflow_profiles(k,4) + & |
---|
288 | inflow_dist(k,j,4,1:nbgp) * inflow_damping_factor(k) |
---|
289 | e(k,j,-nbgp:-1) = mean_inflow_profiles(k,5) + & |
---|
290 | inflow_dist(k,j,5,1:nbgp) * inflow_damping_factor(k) |
---|
291 | e(k,j,-nbgp:-1) = MAX( e(k,j,-nbgp:-1), 0.0_wp ) |
---|
292 | |
---|
293 | IF ( humidity ) & |
---|
294 | q(k,j,-nbgp:-1) = mean_inflow_profiles(k,6) + & |
---|
295 | inflow_dist(k,j,6,1:nbgp) * inflow_damping_factor(k) |
---|
296 | IF ( passive_scalar ) & |
---|
297 | s(k,j,-nbgp:-1) = mean_inflow_profiles(k,7) + & |
---|
298 | inflow_dist(k,j,7,1:nbgp) * inflow_damping_factor(k) |
---|
299 | |
---|
300 | ENDDO |
---|
301 | ENDDO |
---|
302 | |
---|
303 | ENDIF |
---|
304 | |
---|
305 | |
---|
306 | CALL cpu_log( log_point(40), 'inflow_turbulence', 'stop' ) |
---|
307 | |
---|
308 | |
---|
309 | END SUBROUTINE inflow_turbulence |
---|