source: palm/trunk/SOURCE/inflow_turbulence.f90 @ 1605

Last change on this file since 1605 was 1561, checked in by keck, 10 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 9.7 KB
RevLine 
[151]1 SUBROUTINE inflow_turbulence
2
[1036]3!--------------------------------------------------------------------------------!
4! This file is part of PALM.
5!
6! PALM is free software: you can redistribute it and/or modify it under the terms
7! of the GNU General Public License as published by the Free Software Foundation,
8! either version 3 of the License, or (at your option) any later 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!
[1310]17! Copyright 1997-2014 Leibniz Universitaet Hannover
[1036]18!--------------------------------------------------------------------------------!
19!
[484]20! Current revisions:
[151]21! -----------------
[1354]22!
[1561]23!
[151]24! Former revisions:
25! -----------------
26! $Id: inflow_turbulence.f90 1561 2015-03-06 11:00:41Z suehring $
27!
[1561]28! 1560 2015-03-06 10:48:54Z keck
29! Option recycling_yshift added. If this option is switched on, the turbulence
30! data, which is mapped from the recycling plane to the inflow, is shifted in
31! y direction (by ny * dy / 2 )
32!
[1354]33! 1353 2014-04-08 15:21:23Z heinze
34! REAL constants provided with KIND-attribute
35!
[1347]36! 1346 2014-03-27 13:18:20Z heinze
37! Bugfix: REAL constants provided with KIND-attribute especially in call of
38! intrinsic function like MAX, MIN, SIGN
39!
[1321]40! 1320 2014-03-20 08:40:49Z raasch
41! ONLY-attribute added to USE-statements,
42! kind-parameters added to all INTEGER and REAL declaration statements,
43! kinds are defined in new module kinds,
44! revision history before 2012 removed,
45! comment fields (!:) to be used for variable explanations added to
46! all variable declaration statements
47!
[1093]48! 1092 2013-02-02 11:24:22Z raasch
49! unused variables removed
50!
[1037]51! 1036 2012-10-22 13:43:42Z raasch
52! code put under GPL (PALM 3.9)
53!
[198]54! Initial version (2008/03/07)
[151]55!
56! Description:
57! ------------
58! Imposing turbulence at the respective inflow using the turbulence
59! recycling method of Kataoka and Mizuno (2002).
60!------------------------------------------------------------------------------!
61
[1320]62    USE arrays_3d,                                                             &
63        ONLY:  e, inflow_damping_factor, mean_inflow_profiles, pt, u, v, w
64       
65    USE control_parameters,                                                    &
[1560]66        ONLY:  recycling_plane, recycling_yshift
[1320]67       
68    USE cpulog,                                                                &
69        ONLY:  cpu_log, log_point
70       
71    USE grid_variables,                                                        &
72        ONLY: 
73       
74    USE indices,                                                               &
75        ONLY:  nbgp, nxl, ny, nyn, nys, nyng, nysg, nzb, nzt
76       
77    USE kinds
78   
[151]79    USE pegrid
80
81
82    IMPLICIT NONE
83
[1320]84    INTEGER(iwp) ::  i        !:
85    INTEGER(iwp) ::  j        !:
86    INTEGER(iwp) ::  k        !:
87    INTEGER(iwp) ::  l        !:
[1560]88    INTEGER(iwp) ::  next     !:
[1320]89    INTEGER(iwp) ::  ngp_ifd  !:
90    INTEGER(iwp) ::  ngp_pr   !:
[1560]91    INTEGER(iwp) ::  prev     !:
[151]92
[1320]93    REAL(wp), DIMENSION(nzb:nzt+1,5,nbgp)           ::                         &
94       avpr, avpr_l  !:
95    REAL(wp), DIMENSION(nzb:nzt+1,nysg:nyng,5,nbgp) ::                         &
[1560]96       inflow_dist, local_inflow_dist  !:
[151]97
98    CALL cpu_log( log_point(40), 'inflow_turbulence', 'start' )
99
100!
[667]101!-- Carry out spanwise averaging in the recycling plane
[1353]102    avpr_l = 0.0_wp
[667]103    ngp_pr = ( nzt - nzb + 2 ) * 5 * nbgp
104    ngp_ifd = ngp_pr * ( nyn - nys + 1 + 2 * nbgp )
[151]105
106!
107!-- First, local averaging within the recycling domain
[667]108    i = recycling_plane
[151]109
[667]110#if defined( __parallel )
111    IF ( myidx == id_recycling )  THEN
112       
113       DO  l = 1, nbgp
[151]114          DO  j = nys, nyn
[667]115             DO  k = nzb, nzt + 1
[151]116
[667]117                avpr_l(k,1,l) = avpr_l(k,1,l) + u(k,j,i)
118                avpr_l(k,2,l) = avpr_l(k,2,l) + v(k,j,i)
119                avpr_l(k,3,l) = avpr_l(k,3,l) + w(k,j,i)
120                avpr_l(k,4,l) = avpr_l(k,4,l) + pt(k,j,i)
121                avpr_l(k,5,l) = avpr_l(k,5,l) + e(k,j,i)
[151]122
123             ENDDO
124          ENDDO
[667]125          i = i + 1
[151]126       ENDDO
127
128    ENDIF
129!
130!-- Now, averaging over all PEs
[622]131    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
[709]132    CALL MPI_ALLREDUCE( avpr_l(nzb,1,1), avpr(nzb,1,1), ngp_pr, MPI_REAL, &
133                        MPI_SUM, comm2d, ierr )
[667]134
[151]135#else
[667]136    DO  l = 1, nbgp
137       DO  j = nys, nyn
138          DO  k = nzb, nzt + 1
139
140             avpr_l(k,1,l) = avpr_l(k,1,l) + u(k,j,i)
141             avpr_l(k,2,l) = avpr_l(k,2,l) + v(k,j,i)
142             avpr_l(k,3,l) = avpr_l(k,3,l) + w(k,j,i)
143             avpr_l(k,4,l) = avpr_l(k,4,l) + pt(k,j,i)
144             avpr_l(k,5,l) = avpr_l(k,5,l) + e(k,j,i)
145
146          ENDDO
147       ENDDO
148       i = i + 1 
149    ENDDO
150   
[151]151    avpr = avpr_l
152#endif
153
[667]154    avpr = avpr / ( ny + 1 )
[151]155!
156!-- Calculate the disturbances at the recycling plane
157    i = recycling_plane
158
[222]159#if defined( __parallel )
[163]160    IF ( myidx == id_recycling )  THEN
[667]161       DO  l = 1, nbgp
162          DO  j = nysg, nyng
163             DO  k = nzb, nzt + 1
[151]164
[667]165                inflow_dist(k,j,1,l) = u(k,j,i+1) - avpr(k,1,l)
166                inflow_dist(k,j,2,l) = v(k,j,i)   - avpr(k,2,l)
167                inflow_dist(k,j,3,l) = w(k,j,i)   - avpr(k,3,l)
168                inflow_dist(k,j,4,l) = pt(k,j,i)  - avpr(k,4,l)
169                inflow_dist(k,j,5,l) = e(k,j,i)   - avpr(k,5,l)
170             
171            ENDDO
[151]172          ENDDO
[667]173          i = i + 1
[151]174       ENDDO
175
176    ENDIF
[222]177#else
[667]178    DO  l = 1, nbgp
179       DO  j = nysg, nyng
180          DO  k = nzb, nzt+1
[151]181
[667]182             inflow_dist(k,j,1,l) = u(k,j,i+1) - avpr(k,1,l)
183             inflow_dist(k,j,2,l) = v(k,j,i)   - avpr(k,2,l)
184             inflow_dist(k,j,3,l) = w(k,j,i)   - avpr(k,3,l)
185             inflow_dist(k,j,4,l) = pt(k,j,i)  - avpr(k,4,l)
186             inflow_dist(k,j,5,l) = e(k,j,i)   - avpr(k,5,l)
187             
188          ENDDO
[222]189       ENDDO
[667]190       i = i + 1
[222]191    ENDDO
192#endif
193
[151]194!
195!-- For parallel runs, send the disturbances to the respective inflow PE
196#if defined( __parallel )
[163]197    IF ( myidx == id_recycling  .AND.  myidx /= id_inflow )  THEN
[151]198
[1560]199       CALL MPI_SEND( inflow_dist(nzb,nysg,1,1), ngp_ifd, MPI_REAL,            &
[151]200                      id_inflow, 1, comm1dx, ierr )
201
[163]202    ELSEIF ( myidx /= id_recycling  .AND.  myidx == id_inflow )  THEN
[151]203
[1353]204       inflow_dist = 0.0_wp
[1560]205       CALL MPI_RECV( inflow_dist(nzb,nysg,1,1), ngp_ifd, MPI_REAL,            &
[163]206                      id_recycling, 1, comm1dx, status, ierr )
[151]207
208    ENDIF
[1560]209
210   
211    IF ( recycling_yshift .AND. myidx == id_inflow ) THEN
212
213       IF ( pdims(2) >= 2 ) THEN
214 
215          IF ( myidy >= INT( pdims(2) / 2 ) ) THEN
216             prev = myidy - INT( pdims(2) / 2 )
217          ELSE
218             prev = pdims(2) - ( INT( pdims(2) / 2 ) - myidy )
219          ENDIF
220       
221          IF ( myidy < pdims(2) - INT( pdims(2) / 2 ) ) THEN
222             next = myidy + INT( pdims(2) / 2 )
223          ELSE
224             next = INT( pdims(2) / 2 ) - ( pdims(2) - myidy )
225          ENDIF
226
227       ENDIF
228
229       local_inflow_dist = 0.0_wp
230   
231       CALL MPI_SENDRECV( inflow_dist(nzb,nysg,1,1), ngp_ifd, MPI_REAL,        &
232                          next, 1, local_inflow_dist(nzb,nysg,1,1), ngp_ifd,   &
233                          MPI_REAL, prev, 1, comm1dy, status, ierr )
234       
235    ENDIF
236
[151]237#endif
238
239!
240!-- Add the disturbance at the inflow
241    IF ( nxl == 0 )  THEN
[1560]242       
243       IF ( recycling_yshift ) THEN       
[151]244
[1560]245          DO  j = nysg, nyng
246             DO  k = nzb, nzt + 1
[151]247
[1560]248                u(k,j,-nbgp+1:0) = mean_inflow_profiles(k,1) +                 &
249                                   local_inflow_dist(k,j,1,1:nbgp) *           &
250                                   inflow_damping_factor(k)
251                v(k,j,-nbgp:-1)  = mean_inflow_profiles(k,2) +                 &
252                                   local_inflow_dist(k,j,2,1:nbgp) *           &
253                                   inflow_damping_factor(k)
254                w(k,j,-nbgp:-1)  =                                             &
255                                   local_inflow_dist(k,j,3,1:nbgp) *           &
256                                   inflow_damping_factor(k)
257                pt(k,j,-nbgp:-1) = mean_inflow_profiles(k,4) +                 &
258                                   local_inflow_dist(k,j,4,1:nbgp) *           &
259                                   inflow_damping_factor(k)
260                e(k,j,-nbgp:-1)  = mean_inflow_profiles(k,5) +                 &
261                                   local_inflow_dist(k,j,5,1:nbgp) *           &
262                                   inflow_damping_factor(k)
263                e(k,j,-nbgp:-1)  = MAX( e(k,j,-nbgp:-1), 0.0_wp )
264
265             ENDDO
266          ENDDO
267
268       ELSE
269
270          DO  j = nysg, nyng
271             DO  k = nzb, nzt + 1
272 
273                 u(k,j,-nbgp+1:0) = mean_inflow_profiles(k,1) +                   &
274                              inflow_dist(k,j,1,1:nbgp) * inflow_damping_factor(k)
275                 v(k,j,-nbgp:-1)  = mean_inflow_profiles(k,2) +                   &
276                              inflow_dist(k,j,2,1:nbgp) * inflow_damping_factor(k)
277                 w(k,j,-nbgp:-1)  =                                               &
[709]278                           inflow_dist(k,j,3,1:nbgp) * inflow_damping_factor(k)
[1560]279                 pt(k,j,-nbgp:-1) = mean_inflow_profiles(k,4) +                   &
[667]280                           inflow_dist(k,j,4,1:nbgp) * inflow_damping_factor(k)
[1560]281                 e(k,j,-nbgp:-1)  = mean_inflow_profiles(k,5) +                   &
[667]282                           inflow_dist(k,j,5,1:nbgp) * inflow_damping_factor(k)
[1560]283                 e(k,j,-nbgp:-1)  = MAX( e(k,j,-nbgp:-1), 0.0_wp )
[151]284
[1560]285             ENDDO
[151]286          ENDDO
287
[1560]288       ENDIF
289   
[151]290    ENDIF
291
[1560]292
[151]293    CALL cpu_log( log_point(40), 'inflow_turbulence', 'stop' )
294
295
296 END SUBROUTINE inflow_turbulence
Note: See TracBrowser for help on using the repository browser.