1 | SUBROUTINE prandtl_fluxes |
---|
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-2012 Leibniz University Hannover |
---|
18 | !--------------------------------------------------------------------------------! |
---|
19 | ! |
---|
20 | ! Current revisions: |
---|
21 | ! ----------------- |
---|
22 | ! |
---|
23 | ! |
---|
24 | ! Former revisions: |
---|
25 | ! ----------------- |
---|
26 | ! $Id: prandtl_fluxes.f90 1037 2012-10-22 14:10:22Z maronga $ |
---|
27 | ! |
---|
28 | ! 1036 2012-10-22 13:43:42Z raasch |
---|
29 | ! code put under GPL (PALM 3.9) |
---|
30 | ! |
---|
31 | ! 1015 2012-09-27 09:23:24Z raasch |
---|
32 | ! OpenACC statements added |
---|
33 | ! |
---|
34 | ! 978 2012-08-09 08:28:32Z fricke |
---|
35 | ! roughness length for scalar quantities z0h added |
---|
36 | ! |
---|
37 | ! 759 2011-09-15 13:58:31Z raasch |
---|
38 | ! Bugfix for ts limitation |
---|
39 | ! |
---|
40 | ! 709 2011-03-30 09:31:40Z raasch |
---|
41 | ! formatting adjustments |
---|
42 | ! |
---|
43 | ! 667 2010-12-23 12:06:00Z suehring/gryschka |
---|
44 | ! Changed surface boundary conditions for u and v from mirror to Dirichlet. |
---|
45 | ! Therefore u(uzb,:,:) and v(nzb,:,:) are now representative for height z0. |
---|
46 | ! nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng |
---|
47 | ! |
---|
48 | ! 315 2009-05-13 10:57:59Z raasch |
---|
49 | ! Saturation condition at (sea) surface is not used in precursor runs (only |
---|
50 | ! in the following coupled runs) |
---|
51 | ! Bugfix: qsws was calculated in case of constant heatflux = .FALSE. |
---|
52 | ! |
---|
53 | ! 187 2008-08-06 16:25:09Z letzel |
---|
54 | ! Bugfix: modification of the calculation of the vertical turbulent momentum |
---|
55 | ! fluxes u'w' and v'w' |
---|
56 | ! Bugfix: change definition of us_wall from 1D to 2D |
---|
57 | ! Change: modification of the integrated version of the profile function for |
---|
58 | ! momentum for unstable stratification (does not effect results) |
---|
59 | ! |
---|
60 | ! 108 2007-08-24 15:10:38Z letzel |
---|
61 | ! assume saturation at k=nzb_s_inner(j,i) for atmosphere coupled to ocean |
---|
62 | ! |
---|
63 | ! 75 2007-03-22 09:54:05Z raasch |
---|
64 | ! moisture renamed humidity |
---|
65 | ! |
---|
66 | ! RCS Log replace by Id keyword, revision history cleaned up |
---|
67 | ! |
---|
68 | ! Revision 1.19 2006/04/26 12:24:35 raasch |
---|
69 | ! +OpenMP directives and optimization (array assignments replaced by DO loops) |
---|
70 | ! |
---|
71 | ! Revision 1.1 1998/01/23 10:06:06 raasch |
---|
72 | ! Initial revision |
---|
73 | ! |
---|
74 | ! |
---|
75 | ! Description: |
---|
76 | ! ------------ |
---|
77 | ! Diagnostic computation of vertical fluxes in the Prandtl layer from the |
---|
78 | ! values of the variables at grid point k=1 |
---|
79 | !------------------------------------------------------------------------------! |
---|
80 | |
---|
81 | USE arrays_3d |
---|
82 | USE control_parameters |
---|
83 | USE grid_variables |
---|
84 | USE indices |
---|
85 | |
---|
86 | IMPLICIT NONE |
---|
87 | |
---|
88 | INTEGER :: i, j, k |
---|
89 | LOGICAL :: coupled_run |
---|
90 | REAL :: a, b, e_q, rifm, uv_total, z_p |
---|
91 | |
---|
92 | ! |
---|
93 | !-- Data information for accelerators |
---|
94 | !$acc data present( e, nzb_u_inner, nzb_v_inner, nzb_s_inner, pt, q, qs ) & |
---|
95 | !$acc present( qsws, rif, shf, ts, u, us, usws, v, vpt, vsws, zu, zw, z0, z0h ) |
---|
96 | ! |
---|
97 | !-- Compute theta* |
---|
98 | IF ( constant_heatflux ) THEN |
---|
99 | ! |
---|
100 | !-- For a given heat flux in the Prandtl layer: |
---|
101 | !-- for u* use the value from the previous time step |
---|
102 | !$OMP PARALLEL DO |
---|
103 | !$acc kernels do |
---|
104 | DO i = nxlg, nxrg |
---|
105 | DO j = nysg, nyng |
---|
106 | ts(j,i) = -shf(j,i) / ( us(j,i) + 1E-30 ) |
---|
107 | ! |
---|
108 | !-- ts must be limited, because otherwise overflow may occur in case of |
---|
109 | !-- us=0 when computing rif further below |
---|
110 | IF ( ts(j,i) < -1.05E5 ) ts(j,i) = -1.0E5 |
---|
111 | IF ( ts(j,i) > 1.0E5 ) ts(j,i) = 1.0E5 |
---|
112 | ENDDO |
---|
113 | ENDDO |
---|
114 | |
---|
115 | ELSE |
---|
116 | ! |
---|
117 | !-- For a given surface temperature: |
---|
118 | !-- (the Richardson number is still the one from the previous time step) |
---|
119 | !$OMP PARALLEL DO PRIVATE( a, b, k, z_p ) |
---|
120 | !$acc kernels do |
---|
121 | DO i = nxlg, nxrg |
---|
122 | DO j = nysg, nyng |
---|
123 | |
---|
124 | k = nzb_s_inner(j,i) |
---|
125 | z_p = zu(k+1) - zw(k) |
---|
126 | |
---|
127 | IF ( rif(j,i) >= 0.0 ) THEN |
---|
128 | ! |
---|
129 | !-- Stable stratification |
---|
130 | ts(j,i) = kappa * ( pt(k+1,j,i) - pt(k,j,i) ) / ( & |
---|
131 | LOG( z_p / z0h(j,i) ) + & |
---|
132 | 5.0 * rif(j,i) * ( z_p - z0h(j,i) ) / z_p & |
---|
133 | ) |
---|
134 | ELSE |
---|
135 | ! |
---|
136 | !-- Unstable stratification |
---|
137 | a = SQRT( 1.0 - 16.0 * rif(j,i) ) |
---|
138 | b = SQRT( 1.0 - 16.0 * rif(j,i) * z0h(j,i) / z_p ) |
---|
139 | |
---|
140 | ts(j,i) = kappa * ( pt(k+1,j,i) - pt(k,j,i) ) / ( & |
---|
141 | LOG( z_p / z0h(j,i) ) - & |
---|
142 | 2.0 * LOG( ( 1.0 + a ) / ( 1.0 + b ) ) ) |
---|
143 | ENDIF |
---|
144 | |
---|
145 | ENDDO |
---|
146 | ENDDO |
---|
147 | ENDIF |
---|
148 | |
---|
149 | ! |
---|
150 | !-- Compute z_p/L (corresponds to the Richardson-flux number) |
---|
151 | IF ( .NOT. humidity ) THEN |
---|
152 | !$OMP PARALLEL DO PRIVATE( k, z_p ) |
---|
153 | !$acc kernels do |
---|
154 | DO i = nxlg, nxrg |
---|
155 | DO j = nysg, nyng |
---|
156 | k = nzb_s_inner(j,i) |
---|
157 | z_p = zu(k+1) - zw(k) |
---|
158 | rif(j,i) = z_p * kappa * g * ts(j,i) / & |
---|
159 | ( pt(k+1,j,i) * ( us(j,i)**2 + 1E-30 ) ) |
---|
160 | ! |
---|
161 | !-- Limit the value range of the Richardson numbers. |
---|
162 | !-- This is necessary for very small velocities (u,v --> 0), because |
---|
163 | !-- the absolute value of rif can then become very large, which in |
---|
164 | !-- consequence would result in very large shear stresses and very |
---|
165 | !-- small momentum fluxes (both are generally unrealistic). |
---|
166 | IF ( rif(j,i) < rif_min ) rif(j,i) = rif_min |
---|
167 | IF ( rif(j,i) > rif_max ) rif(j,i) = rif_max |
---|
168 | ENDDO |
---|
169 | ENDDO |
---|
170 | ELSE |
---|
171 | !$OMP PARALLEL DO PRIVATE( k, z_p ) |
---|
172 | !$acc kernels do |
---|
173 | DO i = nxlg, nxrg |
---|
174 | DO j = nysg, nyng |
---|
175 | k = nzb_s_inner(j,i) |
---|
176 | z_p = zu(k+1) - zw(k) |
---|
177 | rif(j,i) = z_p * kappa * g * & |
---|
178 | ( ts(j,i) + 0.61 * pt(k+1,j,i) * qs(j,i) ) / & |
---|
179 | ( vpt(k+1,j,i) * ( us(j,i)**2 + 1E-30 ) ) |
---|
180 | ! |
---|
181 | !-- Limit the value range of the Richardson numbers. |
---|
182 | !-- This is necessary for very small velocities (u,v --> 0), because |
---|
183 | !-- the absolute value of rif can then become very large, which in |
---|
184 | !-- consequence would result in very large shear stresses and very |
---|
185 | !-- small momentum fluxes (both are generally unrealistic). |
---|
186 | IF ( rif(j,i) < rif_min ) rif(j,i) = rif_min |
---|
187 | IF ( rif(j,i) > rif_max ) rif(j,i) = rif_max |
---|
188 | ENDDO |
---|
189 | ENDDO |
---|
190 | ENDIF |
---|
191 | |
---|
192 | ! |
---|
193 | !-- Compute u* at the scalars' grid points |
---|
194 | !$OMP PARALLEL DO PRIVATE( a, b, k, uv_total, z_p ) |
---|
195 | !$acc kernels do |
---|
196 | DO i = nxl, nxr |
---|
197 | DO j = nys, nyn |
---|
198 | |
---|
199 | k = nzb_s_inner(j,i) |
---|
200 | z_p = zu(k+1) - zw(k) |
---|
201 | |
---|
202 | ! |
---|
203 | !-- Compute the absolute value of the horizontal velocity |
---|
204 | !-- (relative to the surface) |
---|
205 | uv_total = SQRT( ( 0.5 * ( u(k+1,j,i) + u(k+1,j,i+1) & |
---|
206 | - u(k,j,i) - u(k,j,i+1) ) )**2 + & |
---|
207 | ( 0.5 * ( v(k+1,j,i) + v(k+1,j+1,i) & |
---|
208 | - v(k,j,i) - v(k,j+1,i) ) )**2 ) |
---|
209 | |
---|
210 | |
---|
211 | IF ( rif(j,i) >= 0.0 ) THEN |
---|
212 | ! |
---|
213 | !-- Stable stratification |
---|
214 | us(j,i) = kappa * uv_total / ( & |
---|
215 | LOG( z_p / z0(j,i) ) + & |
---|
216 | 5.0 * rif(j,i) * ( z_p - z0(j,i) ) / z_p & |
---|
217 | ) |
---|
218 | ELSE |
---|
219 | ! |
---|
220 | !-- Unstable stratification |
---|
221 | a = SQRT( SQRT( 1.0 - 16.0 * rif(j,i) ) ) |
---|
222 | b = SQRT( SQRT( 1.0 - 16.0 * rif(j,i) / z_p * z0(j,i) ) ) |
---|
223 | |
---|
224 | us(j,i) = kappa * uv_total / ( & |
---|
225 | LOG( z_p / z0(j,i) ) - & |
---|
226 | LOG( ( 1.0 + a )**2 * ( 1.0 + a**2 ) / ( & |
---|
227 | ( 1.0 + b )**2 * ( 1.0 + b**2 ) ) ) + & |
---|
228 | 2.0 * ( ATAN( a ) - ATAN( b ) ) & |
---|
229 | ) |
---|
230 | ENDIF |
---|
231 | ENDDO |
---|
232 | ENDDO |
---|
233 | |
---|
234 | ! |
---|
235 | !-- Values of us at ghost point locations are needed for the evaluation of usws |
---|
236 | !-- and vsws. |
---|
237 | !$acc update host( us ) |
---|
238 | CALL exchange_horiz_2d( us ) |
---|
239 | !$acc update device( us ) |
---|
240 | |
---|
241 | ! |
---|
242 | !-- Compute u'w' for the total model domain. |
---|
243 | !-- First compute the corresponding component of u* and square it. |
---|
244 | !$OMP PARALLEL DO PRIVATE( a, b, k, rifm, z_p ) |
---|
245 | !$acc kernels do |
---|
246 | DO i = nxl, nxr |
---|
247 | DO j = nys, nyn |
---|
248 | |
---|
249 | k = nzb_u_inner(j,i) |
---|
250 | z_p = zu(k+1) - zw(k) |
---|
251 | |
---|
252 | ! |
---|
253 | !-- Compute Richardson-flux number for this point |
---|
254 | rifm = 0.5 * ( rif(j,i-1) + rif(j,i) ) |
---|
255 | IF ( rifm >= 0.0 ) THEN |
---|
256 | ! |
---|
257 | !-- Stable stratification |
---|
258 | usws(j,i) = kappa * ( u(k+1,j,i) - u(k,j,i) )/ ( & |
---|
259 | LOG( z_p / z0(j,i) ) + & |
---|
260 | 5.0 * rifm * ( z_p - z0(j,i) ) / z_p & |
---|
261 | ) |
---|
262 | ELSE |
---|
263 | ! |
---|
264 | !-- Unstable stratification |
---|
265 | a = SQRT( SQRT( 1.0 - 16.0 * rifm ) ) |
---|
266 | b = SQRT( SQRT( 1.0 - 16.0 * rifm / z_p * z0(j,i) ) ) |
---|
267 | |
---|
268 | usws(j,i) = kappa * ( u(k+1,j,i) - u(k,j,i) ) / ( & |
---|
269 | LOG( z_p / z0(j,i) ) - & |
---|
270 | LOG( (1.0 + a )**2 * ( 1.0 + a**2 ) / ( & |
---|
271 | (1.0 + b )**2 * ( 1.0 + b**2 ) ) ) + & |
---|
272 | 2.0 * ( ATAN( a ) - ATAN( b ) ) & |
---|
273 | ) |
---|
274 | ENDIF |
---|
275 | usws(j,i) = -usws(j,i) * 0.5 * ( us(j,i-1) + us(j,i) ) |
---|
276 | ENDDO |
---|
277 | ENDDO |
---|
278 | |
---|
279 | ! |
---|
280 | !-- Compute v'w' for the total model domain. |
---|
281 | !-- First compute the corresponding component of u* and square it. |
---|
282 | !$OMP PARALLEL DO PRIVATE( a, b, k, rifm, z_p ) |
---|
283 | !$acc kernels do |
---|
284 | DO i = nxl, nxr |
---|
285 | DO j = nys, nyn |
---|
286 | |
---|
287 | k = nzb_v_inner(j,i) |
---|
288 | z_p = zu(k+1) - zw(k) |
---|
289 | |
---|
290 | ! |
---|
291 | !-- Compute Richardson-flux number for this point |
---|
292 | rifm = 0.5 * ( rif(j-1,i) + rif(j,i) ) |
---|
293 | IF ( rifm >= 0.0 ) THEN |
---|
294 | ! |
---|
295 | !-- Stable stratification |
---|
296 | vsws(j,i) = kappa * ( v(k+1,j,i) - v(k,j,i) ) / ( & |
---|
297 | LOG( z_p / z0(j,i) ) + & |
---|
298 | 5.0 * rifm * ( z_p - z0(j,i) ) / z_p & |
---|
299 | ) |
---|
300 | ELSE |
---|
301 | ! |
---|
302 | !-- Unstable stratification |
---|
303 | a = SQRT( SQRT( 1.0 - 16.0 * rifm ) ) |
---|
304 | b = SQRT( SQRT( 1.0 - 16.0 * rifm / z_p * z0(j,i) ) ) |
---|
305 | |
---|
306 | vsws(j,i) = kappa * ( v(k+1,j,i) - v(k,j,i) ) / ( & |
---|
307 | LOG( z_p / z0(j,i) ) - & |
---|
308 | LOG( (1.0 + a )**2 * ( 1.0 + a**2 ) / ( & |
---|
309 | (1.0 + b )**2 * ( 1.0 + b**2 ) ) ) + & |
---|
310 | 2.0 * ( ATAN( a ) - ATAN( b ) ) & |
---|
311 | ) |
---|
312 | ENDIF |
---|
313 | vsws(j,i) = -vsws(j,i) * 0.5 * ( us(j-1,i) + us(j,i) ) |
---|
314 | ENDDO |
---|
315 | ENDDO |
---|
316 | |
---|
317 | ! |
---|
318 | !-- If required compute q* |
---|
319 | IF ( humidity .OR. passive_scalar ) THEN |
---|
320 | IF ( constant_waterflux ) THEN |
---|
321 | ! |
---|
322 | !-- For a given water flux in the Prandtl layer: |
---|
323 | !$OMP PARALLEL DO |
---|
324 | !$acc kernels do |
---|
325 | DO i = nxlg, nxrg |
---|
326 | DO j = nysg, nyng |
---|
327 | qs(j,i) = -qsws(j,i) / ( us(j,i) + 1E-30 ) |
---|
328 | ENDDO |
---|
329 | ENDDO |
---|
330 | |
---|
331 | ELSE |
---|
332 | coupled_run = ( coupling_mode == 'atmosphere_to_ocean' .AND. run_coupled ) |
---|
333 | !$OMP PARALLEL DO PRIVATE( a, b, k, z_p ) |
---|
334 | !$acc kernels do |
---|
335 | DO i = nxlg, nxrg |
---|
336 | DO j = nysg, nyng |
---|
337 | |
---|
338 | k = nzb_s_inner(j,i) |
---|
339 | z_p = zu(k+1) - zw(k) |
---|
340 | |
---|
341 | ! |
---|
342 | !-- Assume saturation for atmosphere coupled to ocean (but not |
---|
343 | !-- in case of precursor runs) |
---|
344 | IF ( coupled_run ) THEN |
---|
345 | e_q = 6.1 * & |
---|
346 | EXP( 0.07 * ( MIN(pt(0,j,i),pt(1,j,i)) - 273.15 ) ) |
---|
347 | q(k,j,i) = 0.622 * e_q / ( surface_pressure - e_q ) |
---|
348 | ENDIF |
---|
349 | IF ( rif(j,i) >= 0.0 ) THEN |
---|
350 | ! |
---|
351 | !-- Stable stratification |
---|
352 | qs(j,i) = kappa * ( q(k+1,j,i) - q(k,j,i) ) / ( & |
---|
353 | LOG( z_p / z0h(j,i) ) + & |
---|
354 | 5.0 * rif(j,i) * ( z_p - z0h(j,i) ) / z_p & |
---|
355 | ) |
---|
356 | ELSE |
---|
357 | ! |
---|
358 | !-- Unstable stratification |
---|
359 | a = SQRT( 1.0 - 16.0 * rif(j,i) ) |
---|
360 | b = SQRT( 1.0 - 16.0 * rif(j,i) * z0h(j,i) / z_p ) |
---|
361 | |
---|
362 | qs(j,i) = kappa * ( q(k+1,j,i) - q(k,j,i) ) / ( & |
---|
363 | LOG( z_p / z0h(j,i) ) - & |
---|
364 | 2.0 * LOG( (1.0 + a ) / ( 1.0 + b ) ) ) |
---|
365 | ENDIF |
---|
366 | |
---|
367 | ENDDO |
---|
368 | ENDDO |
---|
369 | ENDIF |
---|
370 | ENDIF |
---|
371 | |
---|
372 | ! |
---|
373 | !-- Exchange the boundaries for the momentum fluxes (only for sake of |
---|
374 | !-- completeness) |
---|
375 | !$acc update host( usws, vsws ) |
---|
376 | CALL exchange_horiz_2d( usws ) |
---|
377 | CALL exchange_horiz_2d( vsws ) |
---|
378 | !$acc update device( usws, vsws ) |
---|
379 | IF ( humidity .OR. passive_scalar ) THEN |
---|
380 | !$acc update host( qsws ) |
---|
381 | CALL exchange_horiz_2d( qsws ) |
---|
382 | !$acc update device( qsws ) |
---|
383 | ENDIF |
---|
384 | |
---|
385 | ! |
---|
386 | !-- Compute the vertical kinematic heat flux |
---|
387 | IF ( .NOT. constant_heatflux ) THEN |
---|
388 | !$OMP PARALLEL DO |
---|
389 | !$acc kernels do |
---|
390 | DO i = nxlg, nxrg |
---|
391 | DO j = nysg, nyng |
---|
392 | shf(j,i) = -ts(j,i) * us(j,i) |
---|
393 | ENDDO |
---|
394 | ENDDO |
---|
395 | ENDIF |
---|
396 | |
---|
397 | ! |
---|
398 | !-- Compute the vertical water/scalar flux |
---|
399 | IF ( .NOT. constant_waterflux .AND. ( humidity .OR. passive_scalar ) ) THEN |
---|
400 | !$OMP PARALLEL DO |
---|
401 | !$acc kernels do |
---|
402 | DO i = nxlg, nxrg |
---|
403 | DO j = nysg, nyng |
---|
404 | qsws(j,i) = -qs(j,i) * us(j,i) |
---|
405 | ENDDO |
---|
406 | ENDDO |
---|
407 | ENDIF |
---|
408 | |
---|
409 | ! |
---|
410 | !-- Bottom boundary condition for the TKE |
---|
411 | IF ( ibc_e_b == 2 ) THEN |
---|
412 | !$OMP PARALLEL DO |
---|
413 | !$acc kernels do |
---|
414 | DO i = nxlg, nxrg |
---|
415 | DO j = nysg, nyng |
---|
416 | e(nzb_s_inner(j,i)+1,j,i) = ( us(j,i) / 0.1 )**2 |
---|
417 | ! |
---|
418 | !-- As a test: cm = 0.4 |
---|
419 | ! e(nzb_s_inner(j,i)+1,j,i) = ( us(j,i) / 0.4 )**2 |
---|
420 | e(nzb_s_inner(j,i),j,i) = e(nzb_s_inner(j,i)+1,j,i) |
---|
421 | ENDDO |
---|
422 | ENDDO |
---|
423 | ENDIF |
---|
424 | |
---|
425 | !$acc end data |
---|
426 | |
---|
427 | END SUBROUTINE prandtl_fluxes |
---|