1 | SUBROUTINE inflow_turbulence |
---|
2 | |
---|
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 | ! |
---|
17 | ! Copyright 1997-2014 Leibniz Universitaet Hannover |
---|
18 | !--------------------------------------------------------------------------------! |
---|
19 | ! |
---|
20 | ! Current revisions: |
---|
21 | ! ----------------- |
---|
22 | ! ONLY-attribute added to USE-statements, |
---|
23 | ! kind-parameters added to all INTEGER and REAL declaration statements, |
---|
24 | ! kinds are defined in new module kinds, |
---|
25 | ! old module precision_kind is removed, |
---|
26 | ! revision history before 2012 removed, |
---|
27 | ! comment fields (!:) to be used for variable explanations added to |
---|
28 | ! all variable declaration statements |
---|
29 | ! |
---|
30 | ! module interfaces removed |
---|
31 | ! |
---|
32 | ! Former revisions: |
---|
33 | ! ----------------- |
---|
34 | ! $Id: inflow_turbulence.f90 1320 2014-03-20 08:40:49Z raasch $ |
---|
35 | ! |
---|
36 | ! 1092 2013-02-02 11:24:22Z raasch |
---|
37 | ! unused variables removed |
---|
38 | ! |
---|
39 | ! 1036 2012-10-22 13:43:42Z raasch |
---|
40 | ! code put under GPL (PALM 3.9) |
---|
41 | ! |
---|
42 | ! Initial version (2008/03/07) |
---|
43 | ! |
---|
44 | ! Description: |
---|
45 | ! ------------ |
---|
46 | ! Imposing turbulence at the respective inflow using the turbulence |
---|
47 | ! recycling method of Kataoka and Mizuno (2002). |
---|
48 | !------------------------------------------------------------------------------! |
---|
49 | |
---|
50 | USE arrays_3d, & |
---|
51 | ONLY: e, inflow_damping_factor, mean_inflow_profiles, pt, u, v, w |
---|
52 | |
---|
53 | USE control_parameters, & |
---|
54 | ONLY: recycling_plane |
---|
55 | |
---|
56 | USE cpulog, & |
---|
57 | ONLY: cpu_log, log_point |
---|
58 | |
---|
59 | USE grid_variables, & |
---|
60 | ONLY: |
---|
61 | |
---|
62 | USE indices, & |
---|
63 | ONLY: nbgp, nxl, ny, nyn, nys, nyng, nysg, nzb, nzt |
---|
64 | |
---|
65 | USE kinds |
---|
66 | |
---|
67 | USE pegrid |
---|
68 | |
---|
69 | |
---|
70 | IMPLICIT NONE |
---|
71 | |
---|
72 | INTEGER(iwp) :: i !: |
---|
73 | INTEGER(iwp) :: j !: |
---|
74 | INTEGER(iwp) :: k !: |
---|
75 | INTEGER(iwp) :: l !: |
---|
76 | INTEGER(iwp) :: ngp_ifd !: |
---|
77 | INTEGER(iwp) :: ngp_pr !: |
---|
78 | |
---|
79 | REAL(wp), DIMENSION(nzb:nzt+1,5,nbgp) :: & |
---|
80 | avpr, avpr_l !: |
---|
81 | REAL(wp), DIMENSION(nzb:nzt+1,nysg:nyng,5,nbgp) :: & |
---|
82 | inflow_dist !: |
---|
83 | |
---|
84 | CALL cpu_log( log_point(40), 'inflow_turbulence', 'start' ) |
---|
85 | |
---|
86 | ! |
---|
87 | !-- Carry out spanwise averaging in the recycling plane |
---|
88 | avpr_l = 0.0 |
---|
89 | ngp_pr = ( nzt - nzb + 2 ) * 5 * nbgp |
---|
90 | ngp_ifd = ngp_pr * ( nyn - nys + 1 + 2 * nbgp ) |
---|
91 | |
---|
92 | ! |
---|
93 | !-- First, local averaging within the recycling domain |
---|
94 | i = recycling_plane |
---|
95 | |
---|
96 | #if defined( __parallel ) |
---|
97 | IF ( myidx == id_recycling ) THEN |
---|
98 | |
---|
99 | DO l = 1, nbgp |
---|
100 | DO j = nys, nyn |
---|
101 | DO k = nzb, nzt + 1 |
---|
102 | |
---|
103 | avpr_l(k,1,l) = avpr_l(k,1,l) + u(k,j,i) |
---|
104 | avpr_l(k,2,l) = avpr_l(k,2,l) + v(k,j,i) |
---|
105 | avpr_l(k,3,l) = avpr_l(k,3,l) + w(k,j,i) |
---|
106 | avpr_l(k,4,l) = avpr_l(k,4,l) + pt(k,j,i) |
---|
107 | avpr_l(k,5,l) = avpr_l(k,5,l) + e(k,j,i) |
---|
108 | |
---|
109 | ENDDO |
---|
110 | ENDDO |
---|
111 | i = i + 1 |
---|
112 | ENDDO |
---|
113 | |
---|
114 | ENDIF |
---|
115 | ! |
---|
116 | !-- Now, averaging over all PEs |
---|
117 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
118 | CALL MPI_ALLREDUCE( avpr_l(nzb,1,1), avpr(nzb,1,1), ngp_pr, MPI_REAL, & |
---|
119 | MPI_SUM, comm2d, ierr ) |
---|
120 | |
---|
121 | #else |
---|
122 | DO l = 1, nbgp |
---|
123 | DO j = nys, nyn |
---|
124 | DO k = nzb, nzt + 1 |
---|
125 | |
---|
126 | avpr_l(k,1,l) = avpr_l(k,1,l) + u(k,j,i) |
---|
127 | avpr_l(k,2,l) = avpr_l(k,2,l) + v(k,j,i) |
---|
128 | avpr_l(k,3,l) = avpr_l(k,3,l) + w(k,j,i) |
---|
129 | avpr_l(k,4,l) = avpr_l(k,4,l) + pt(k,j,i) |
---|
130 | avpr_l(k,5,l) = avpr_l(k,5,l) + e(k,j,i) |
---|
131 | |
---|
132 | ENDDO |
---|
133 | ENDDO |
---|
134 | i = i + 1 |
---|
135 | ENDDO |
---|
136 | |
---|
137 | avpr = avpr_l |
---|
138 | #endif |
---|
139 | |
---|
140 | avpr = avpr / ( ny + 1 ) |
---|
141 | ! |
---|
142 | !-- Calculate the disturbances at the recycling plane |
---|
143 | i = recycling_plane |
---|
144 | |
---|
145 | #if defined( __parallel ) |
---|
146 | IF ( myidx == id_recycling ) THEN |
---|
147 | DO l = 1, nbgp |
---|
148 | DO j = nysg, nyng |
---|
149 | DO k = nzb, nzt + 1 |
---|
150 | |
---|
151 | inflow_dist(k,j,1,l) = u(k,j,i+1) - avpr(k,1,l) |
---|
152 | inflow_dist(k,j,2,l) = v(k,j,i) - avpr(k,2,l) |
---|
153 | inflow_dist(k,j,3,l) = w(k,j,i) - avpr(k,3,l) |
---|
154 | inflow_dist(k,j,4,l) = pt(k,j,i) - avpr(k,4,l) |
---|
155 | inflow_dist(k,j,5,l) = e(k,j,i) - avpr(k,5,l) |
---|
156 | |
---|
157 | ENDDO |
---|
158 | ENDDO |
---|
159 | i = i + 1 |
---|
160 | ENDDO |
---|
161 | |
---|
162 | ENDIF |
---|
163 | #else |
---|
164 | DO l = 1, nbgp |
---|
165 | DO j = nysg, nyng |
---|
166 | DO k = nzb, nzt+1 |
---|
167 | |
---|
168 | inflow_dist(k,j,1,l) = u(k,j,i+1) - avpr(k,1,l) |
---|
169 | inflow_dist(k,j,2,l) = v(k,j,i) - avpr(k,2,l) |
---|
170 | inflow_dist(k,j,3,l) = w(k,j,i) - avpr(k,3,l) |
---|
171 | inflow_dist(k,j,4,l) = pt(k,j,i) - avpr(k,4,l) |
---|
172 | inflow_dist(k,j,5,l) = e(k,j,i) - avpr(k,5,l) |
---|
173 | |
---|
174 | ENDDO |
---|
175 | ENDDO |
---|
176 | i = i + 1 |
---|
177 | ENDDO |
---|
178 | #endif |
---|
179 | |
---|
180 | ! |
---|
181 | !-- For parallel runs, send the disturbances to the respective inflow PE |
---|
182 | #if defined( __parallel ) |
---|
183 | IF ( myidx == id_recycling .AND. myidx /= id_inflow ) THEN |
---|
184 | |
---|
185 | CALL MPI_SEND( inflow_dist(nzb,nysg,1,1), ngp_ifd, MPI_REAL, & |
---|
186 | id_inflow, 1, comm1dx, ierr ) |
---|
187 | |
---|
188 | ELSEIF ( myidx /= id_recycling .AND. myidx == id_inflow ) THEN |
---|
189 | |
---|
190 | inflow_dist = 0.0 |
---|
191 | CALL MPI_RECV( inflow_dist(nzb,nysg,1,1), ngp_ifd, MPI_REAL, & |
---|
192 | id_recycling, 1, comm1dx, status, ierr ) |
---|
193 | |
---|
194 | ENDIF |
---|
195 | #endif |
---|
196 | |
---|
197 | ! |
---|
198 | !-- Add the disturbance at the inflow |
---|
199 | IF ( nxl == 0 ) THEN |
---|
200 | |
---|
201 | DO j = nysg, nyng |
---|
202 | DO k = nzb, nzt + 1 |
---|
203 | |
---|
204 | u(k,j,-nbgp+1:0) = mean_inflow_profiles(k,1) + & |
---|
205 | inflow_dist(k,j,1,1:nbgp) * inflow_damping_factor(k) |
---|
206 | v(k,j,-nbgp:-1) = mean_inflow_profiles(k,2) + & |
---|
207 | inflow_dist(k,j,2,1:nbgp) * inflow_damping_factor(k) |
---|
208 | w(k,j,-nbgp:-1) = & |
---|
209 | inflow_dist(k,j,3,1:nbgp) * inflow_damping_factor(k) |
---|
210 | pt(k,j,-nbgp:-1) = mean_inflow_profiles(k,4) + & |
---|
211 | inflow_dist(k,j,4,1:nbgp) * inflow_damping_factor(k) |
---|
212 | e(k,j,-nbgp:-1) = mean_inflow_profiles(k,5) + & |
---|
213 | inflow_dist(k,j,5,1:nbgp) * inflow_damping_factor(k) |
---|
214 | e(k,j,-nbgp:-1) = MAX( e(k,j,-nbgp:-1), 0.0 ) |
---|
215 | |
---|
216 | ENDDO |
---|
217 | ENDDO |
---|
218 | |
---|
219 | ENDIF |
---|
220 | |
---|
221 | CALL cpu_log( log_point(40), 'inflow_turbulence', 'stop' ) |
---|
222 | |
---|
223 | |
---|
224 | END SUBROUTINE inflow_turbulence |
---|