SUBROUTINE lpm_advec !--------------------------------------------------------------------------------! ! This file is part of PALM. ! ! PALM is free software: you can redistribute it and/or modify it under the terms ! of the GNU General Public License as published by the Free Software Foundation, ! either version 3 of the License, or (at your option) any later version. ! ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. ! ! You should have received a copy of the GNU General Public License along with ! PALM. If not, see . ! ! Copyright 1997-2012 Leibniz University Hannover !--------------------------------------------------------------------------------! ! ! Current revisions: ! ------------------ ! ! ! Former revisions: ! ----------------- ! $Id: lpm_advec.f90 1037 2012-10-22 14:10:22Z Giersch $ ! ! 1036 2012-10-22 13:43:42Z raasch ! code put under GPL (PALM 3.9) ! ! 849 2012-03-15 10:35:09Z raasch ! initial revision (former part of advec_particles) ! ! ! Description: ! ------------ ! Calculation of new particle positions due to advection using a simple Euler ! scheme. Particles may feel inertia effects. SGS transport can be included ! using the stochastic model of Weil et al. (2004, JAS, 61, 2877-2887). !------------------------------------------------------------------------------! USE arrays_3d USE control_parameters USE grid_variables USE indices USE particle_attributes USE statistics IMPLICIT NONE INTEGER :: i, j, k, n REAL :: aa, bb, cc, dd, dens_ratio, exp_arg, exp_term, gg, u_int, & u_int_l, u_int_u, v_int, v_int_l, v_int_u, w_int, w_int_l, & w_int_u, x, y INTEGER :: agp, kw, num_gp INTEGER :: gp_outside_of_building(1:8) REAL :: d_sum, de_dx_int, de_dx_int_l, de_dx_int_u, de_dy_int, & de_dy_int_l, de_dy_int_u, de_dt, de_dt_min, de_dz_int, & de_dz_int_l, de_dz_int_u, diss_int, diss_int_l, diss_int_u, & dt_gap, dt_particle, dt_particle_m, e_int, e_int_l, e_int_u, & e_mean_int, fs_int, lagr_timescale, random_gauss, vv_int REAL :: location(1:30,1:3) REAL, DIMENSION(1:30) :: de_dxi, de_dyi, de_dzi, dissi, d_gp_pl, ei DO n = 1, number_of_particles ! !-- Move particle only if the LES timestep has not (approximately) been !-- reached IF ( ( dt_3d - particles(n)%dt_sum ) < 1E-8 ) CYCLE ! !-- Interpolate u velocity-component, determine left, front, bottom !-- index of u-array i = ( particles(n)%x + 0.5 * dx ) * ddx j = particles(n)%y * ddy k = ( particles(n)%z + 0.5 * dz * atmos_ocean_sign ) / dz & + offset_ocean_nzt ! only exact if equidistant ! !-- Interpolation of the velocity components in the xy-plane x = particles(n)%x + ( 0.5 - i ) * dx y = particles(n)%y - j * dy aa = x**2 + y**2 bb = ( dx - x )**2 + y**2 cc = x**2 + ( dy - y )**2 dd = ( dx - x )**2 + ( dy - y )**2 gg = aa + bb + cc + dd u_int_l = ( ( gg - aa ) * u(k,j,i) + ( gg - bb ) * u(k,j,i+1) & + ( gg - cc ) * u(k,j+1,i) + ( gg - dd ) * u(k,j+1,i+1) & ) / ( 3.0 * gg ) - u_gtrans IF ( k+1 == nzt+1 ) THEN u_int = u_int_l ELSE u_int_u = ( ( gg-aa ) * u(k+1,j,i) + ( gg-bb ) * u(k+1,j,i+1) & + ( gg-cc ) * u(k+1,j+1,i) + ( gg-dd ) * u(k+1,j+1,i+1) & ) / ( 3.0 * gg ) - u_gtrans u_int = u_int_l + ( particles(n)%z - zu(k) ) / dz * & ( u_int_u - u_int_l ) ENDIF ! !-- Same procedure for interpolation of the v velocity-component (adopt !-- index k from u velocity-component) i = particles(n)%x * ddx j = ( particles(n)%y + 0.5 * dy ) * ddy x = particles(n)%x - i * dx y = particles(n)%y + ( 0.5 - j ) * dy aa = x**2 + y**2 bb = ( dx - x )**2 + y**2 cc = x**2 + ( dy - y )**2 dd = ( dx - x )**2 + ( dy - y )**2 gg = aa + bb + cc + dd v_int_l = ( ( gg - aa ) * v(k,j,i) + ( gg - bb ) * v(k,j,i+1) & + ( gg - cc ) * v(k,j+1,i) + ( gg - dd ) * v(k,j+1,i+1) & ) / ( 3.0 * gg ) - v_gtrans IF ( k+1 == nzt+1 ) THEN v_int = v_int_l ELSE v_int_u = ( ( gg-aa ) * v(k+1,j,i) + ( gg-bb ) * v(k+1,j,i+1) & + ( gg-cc ) * v(k+1,j+1,i) + ( gg-dd ) * v(k+1,j+1,i+1) & ) / ( 3.0 * gg ) - v_gtrans v_int = v_int_l + ( particles(n)%z - zu(k) ) / dz * & ( v_int_u - v_int_l ) ENDIF ! !-- Same procedure for interpolation of the w velocity-component (adopt !-- index i from v velocity-component) IF ( vertical_particle_advection(particles(n)%group) ) THEN j = particles(n)%y * ddy k = particles(n)%z / dz + offset_ocean_nzt_m1 x = particles(n)%x - i * dx y = particles(n)%y - j * dy aa = x**2 + y**2 bb = ( dx - x )**2 + y**2 cc = x**2 + ( dy - y )**2 dd = ( dx - x )**2 + ( dy - y )**2 gg = aa + bb + cc + dd w_int_l = ( ( gg - aa ) * w(k,j,i) + ( gg - bb ) * w(k,j,i+1) & + ( gg - cc ) * w(k,j+1,i) + ( gg - dd ) * w(k,j+1,i+1) & ) / ( 3.0 * gg ) IF ( k+1 == nzt+1 ) THEN w_int = w_int_l ELSE w_int_u = ( ( gg-aa ) * w(k+1,j,i) + & ( gg-bb ) * w(k+1,j,i+1) + & ( gg-cc ) * w(k+1,j+1,i) + & ( gg-dd ) * w(k+1,j+1,i+1) & ) / ( 3.0 * gg ) w_int = w_int_l + ( particles(n)%z - zw(k) ) / dz * & ( w_int_u - w_int_l ) ENDIF ELSE w_int = 0.0 ENDIF ! !-- Interpolate and calculate quantities needed for calculating the SGS !-- velocities IF ( use_sgs_for_particles ) THEN ! !-- Interpolate TKE i = particles(n)%x * ddx j = particles(n)%y * ddy k = ( particles(n)%z + 0.5 * dz * atmos_ocean_sign ) / dz & + offset_ocean_nzt ! only exact if eq.dist IF ( topography == 'flat' ) THEN x = particles(n)%x - i * dx y = particles(n)%y - j * dy aa = x**2 + y**2 bb = ( dx - x )**2 + y**2 cc = x**2 + ( dy - y )**2 dd = ( dx - x )**2 + ( dy - y )**2 gg = aa + bb + cc + dd e_int_l = ( ( gg-aa ) * e(k,j,i) + ( gg-bb ) * e(k,j,i+1) & + ( gg-cc ) * e(k,j+1,i) + ( gg-dd ) * e(k,j+1,i+1) & ) / ( 3.0 * gg ) IF ( k+1 == nzt+1 ) THEN e_int = e_int_l ELSE e_int_u = ( ( gg - aa ) * e(k+1,j,i) + & ( gg - bb ) * e(k+1,j,i+1) + & ( gg - cc ) * e(k+1,j+1,i) + & ( gg - dd ) * e(k+1,j+1,i+1) & ) / ( 3.0 * gg ) e_int = e_int_l + ( particles(n)%z - zu(k) ) / dz * & ( e_int_u - e_int_l ) ENDIF ! !-- Interpolate the TKE gradient along x (adopt incides i,j,k and !-- all position variables from above (TKE)) de_dx_int_l = ( ( gg - aa ) * de_dx(k,j,i) + & ( gg - bb ) * de_dx(k,j,i+1) + & ( gg - cc ) * de_dx(k,j+1,i) + & ( gg - dd ) * de_dx(k,j+1,i+1) & ) / ( 3.0 * gg ) IF ( ( k+1 == nzt+1 ) .OR. ( k == nzb ) ) THEN de_dx_int = de_dx_int_l ELSE de_dx_int_u = ( ( gg - aa ) * de_dx(k+1,j,i) + & ( gg - bb ) * de_dx(k+1,j,i+1) + & ( gg - cc ) * de_dx(k+1,j+1,i) + & ( gg - dd ) * de_dx(k+1,j+1,i+1) & ) / ( 3.0 * gg ) de_dx_int = de_dx_int_l + ( particles(n)%z - zu(k) ) / dz * & ( de_dx_int_u - de_dx_int_l ) ENDIF ! !-- Interpolate the TKE gradient along y de_dy_int_l = ( ( gg - aa ) * de_dy(k,j,i) + & ( gg - bb ) * de_dy(k,j,i+1) + & ( gg - cc ) * de_dy(k,j+1,i) + & ( gg - dd ) * de_dy(k,j+1,i+1) & ) / ( 3.0 * gg ) IF ( ( k+1 == nzt+1 ) .OR. ( k == nzb ) ) THEN de_dy_int = de_dy_int_l ELSE de_dy_int_u = ( ( gg - aa ) * de_dy(k+1,j,i) + & ( gg - bb ) * de_dy(k+1,j,i+1) + & ( gg - cc ) * de_dy(k+1,j+1,i) + & ( gg - dd ) * de_dy(k+1,j+1,i+1) & ) / ( 3.0 * gg ) de_dy_int = de_dy_int_l + ( particles(n)%z - zu(k) ) / dz * & ( de_dy_int_u - de_dy_int_l ) ENDIF ! !-- Interpolate the TKE gradient along z IF ( particles(n)%z < 0.5 * dz ) THEN de_dz_int = 0.0 ELSE de_dz_int_l = ( ( gg - aa ) * de_dz(k,j,i) + & ( gg - bb ) * de_dz(k,j,i+1) + & ( gg - cc ) * de_dz(k,j+1,i) + & ( gg - dd ) * de_dz(k,j+1,i+1) & ) / ( 3.0 * gg ) IF ( ( k+1 == nzt+1 ) .OR. ( k == nzb ) ) THEN de_dz_int = de_dz_int_l ELSE de_dz_int_u = ( ( gg - aa ) * de_dz(k+1,j,i) + & ( gg - bb ) * de_dz(k+1,j,i+1) + & ( gg - cc ) * de_dz(k+1,j+1,i) + & ( gg - dd ) * de_dz(k+1,j+1,i+1) & ) / ( 3.0 * gg ) de_dz_int = de_dz_int_l + ( particles(n)%z - zu(k) ) / dz * & ( de_dz_int_u - de_dz_int_l ) ENDIF ENDIF ! !-- Interpolate the dissipation of TKE diss_int_l = ( ( gg - aa ) * diss(k,j,i) + & ( gg - bb ) * diss(k,j,i+1) + & ( gg - cc ) * diss(k,j+1,i) + & ( gg - dd ) * diss(k,j+1,i+1) & ) / ( 3.0 * gg ) IF ( k+1 == nzt+1 ) THEN diss_int = diss_int_l ELSE diss_int_u = ( ( gg - aa ) * diss(k+1,j,i) + & ( gg - bb ) * diss(k+1,j,i+1) + & ( gg - cc ) * diss(k+1,j+1,i) + & ( gg - dd ) * diss(k+1,j+1,i+1) & ) / ( 3.0 * gg ) diss_int = diss_int_l + ( particles(n)%z - zu(k) ) / dz * & ( diss_int_u - diss_int_l ) ENDIF ELSE ! !-- In case that there are buildings it has to be determined !-- how many of the gridpoints defining the particle box are !-- situated within a building !-- gp_outside_of_building(1): i,j,k !-- gp_outside_of_building(2): i,j+1,k !-- gp_outside_of_building(3): i,j,k+1 !-- gp_outside_of_building(4): i,j+1,k+1 !-- gp_outside_of_building(5): i+1,j,k !-- gp_outside_of_building(6): i+1,j+1,k !-- gp_outside_of_building(7): i+1,j,k+1 !-- gp_outside_of_building(8): i+1,j+1,k+1 gp_outside_of_building = 0 location = 0.0 num_gp = 0 IF ( k > nzb_s_inner(j,i) .OR. nzb_s_inner(j,i) == 0 ) THEN num_gp = num_gp + 1 gp_outside_of_building(1) = 1 location(num_gp,1) = i * dx location(num_gp,2) = j * dy location(num_gp,3) = k * dz - 0.5 * dz ei(num_gp) = e(k,j,i) dissi(num_gp) = diss(k,j,i) de_dxi(num_gp) = de_dx(k,j,i) de_dyi(num_gp) = de_dy(k,j,i) de_dzi(num_gp) = de_dz(k,j,i) ENDIF IF ( k > nzb_s_inner(j+1,i) .OR. nzb_s_inner(j+1,i) == 0 ) & THEN num_gp = num_gp + 1 gp_outside_of_building(2) = 1 location(num_gp,1) = i * dx location(num_gp,2) = (j+1) * dy location(num_gp,3) = k * dz - 0.5 * dz ei(num_gp) = e(k,j+1,i) dissi(num_gp) = diss(k,j+1,i) de_dxi(num_gp) = de_dx(k,j+1,i) de_dyi(num_gp) = de_dy(k,j+1,i) de_dzi(num_gp) = de_dz(k,j+1,i) ENDIF IF ( k+1 > nzb_s_inner(j,i) .OR. nzb_s_inner(j,i) == 0 ) THEN num_gp = num_gp + 1 gp_outside_of_building(3) = 1 location(num_gp,1) = i * dx location(num_gp,2) = j * dy location(num_gp,3) = (k+1) * dz - 0.5 * dz ei(num_gp) = e(k+1,j,i) dissi(num_gp) = diss(k+1,j,i) de_dxi(num_gp) = de_dx(k+1,j,i) de_dyi(num_gp) = de_dy(k+1,j,i) de_dzi(num_gp) = de_dz(k+1,j,i) ENDIF IF ( k+1 > nzb_s_inner(j+1,i) .OR. nzb_s_inner(j+1,i) == 0 ) & THEN num_gp = num_gp + 1 gp_outside_of_building(4) = 1 location(num_gp,1) = i * dx location(num_gp,2) = (j+1) * dy location(num_gp,3) = (k+1) * dz - 0.5 * dz ei(num_gp) = e(k+1,j+1,i) dissi(num_gp) = diss(k+1,j+1,i) de_dxi(num_gp) = de_dx(k+1,j+1,i) de_dyi(num_gp) = de_dy(k+1,j+1,i) de_dzi(num_gp) = de_dz(k+1,j+1,i) ENDIF IF ( k > nzb_s_inner(j,i+1) .OR. nzb_s_inner(j,i+1) == 0 ) & THEN num_gp = num_gp + 1 gp_outside_of_building(5) = 1 location(num_gp,1) = (i+1) * dx location(num_gp,2) = j * dy location(num_gp,3) = k * dz - 0.5 * dz ei(num_gp) = e(k,j,i+1) dissi(num_gp) = diss(k,j,i+1) de_dxi(num_gp) = de_dx(k,j,i+1) de_dyi(num_gp) = de_dy(k,j,i+1) de_dzi(num_gp) = de_dz(k,j,i+1) ENDIF IF ( k > nzb_s_inner(j+1,i+1) .OR. nzb_s_inner(j+1,i+1) == 0 ) & THEN num_gp = num_gp + 1 gp_outside_of_building(6) = 1 location(num_gp,1) = (i+1) * dx location(num_gp,2) = (j+1) * dy location(num_gp,3) = k * dz - 0.5 * dz ei(num_gp) = e(k,j+1,i+1) dissi(num_gp) = diss(k,j+1,i+1) de_dxi(num_gp) = de_dx(k,j+1,i+1) de_dyi(num_gp) = de_dy(k,j+1,i+1) de_dzi(num_gp) = de_dz(k,j+1,i+1) ENDIF IF ( k+1 > nzb_s_inner(j,i+1) .OR. nzb_s_inner(j,i+1) == 0 ) & THEN num_gp = num_gp + 1 gp_outside_of_building(7) = 1 location(num_gp,1) = (i+1) * dx location(num_gp,2) = j * dy location(num_gp,3) = (k+1) * dz - 0.5 * dz ei(num_gp) = e(k+1,j,i+1) dissi(num_gp) = diss(k+1,j,i+1) de_dxi(num_gp) = de_dx(k+1,j,i+1) de_dyi(num_gp) = de_dy(k+1,j,i+1) de_dzi(num_gp) = de_dz(k+1,j,i+1) ENDIF IF ( k+1 > nzb_s_inner(j+1,i+1) .OR. nzb_s_inner(j+1,i+1) == 0)& THEN num_gp = num_gp + 1 gp_outside_of_building(8) = 1 location(num_gp,1) = (i+1) * dx location(num_gp,2) = (j+1) * dy location(num_gp,3) = (k+1) * dz - 0.5 * dz ei(num_gp) = e(k+1,j+1,i+1) dissi(num_gp) = diss(k+1,j+1,i+1) de_dxi(num_gp) = de_dx(k+1,j+1,i+1) de_dyi(num_gp) = de_dy(k+1,j+1,i+1) de_dzi(num_gp) = de_dz(k+1,j+1,i+1) ENDIF ! !-- If all gridpoints are situated outside of a building, then the !-- ordinary interpolation scheme can be used. IF ( num_gp == 8 ) THEN x = particles(n)%x - i * dx y = particles(n)%y - j * dy aa = x**2 + y**2 bb = ( dx - x )**2 + y**2 cc = x**2 + ( dy - y )**2 dd = ( dx - x )**2 + ( dy - y )**2 gg = aa + bb + cc + dd e_int_l = (( gg-aa ) * e(k,j,i) + ( gg-bb ) * e(k,j,i+1) & + ( gg-cc ) * e(k,j+1,i) + ( gg-dd ) * e(k,j+1,i+1)& ) / ( 3.0 * gg ) IF ( k+1 == nzt+1 ) THEN e_int = e_int_l ELSE e_int_u = ( ( gg - aa ) * e(k+1,j,i) + & ( gg - bb ) * e(k+1,j,i+1) + & ( gg - cc ) * e(k+1,j+1,i) + & ( gg - dd ) * e(k+1,j+1,i+1) & ) / ( 3.0 * gg ) e_int = e_int_l + ( particles(n)%z - zu(k) ) / dz * & ( e_int_u - e_int_l ) ENDIF ! !-- Interpolate the TKE gradient along x (adopt incides i,j,k !-- and all position variables from above (TKE)) de_dx_int_l = ( ( gg - aa ) * de_dx(k,j,i) + & ( gg - bb ) * de_dx(k,j,i+1) + & ( gg - cc ) * de_dx(k,j+1,i) + & ( gg - dd ) * de_dx(k,j+1,i+1) & ) / ( 3.0 * gg ) IF ( ( k+1 == nzt+1 ) .OR. ( k == nzb ) ) THEN de_dx_int = de_dx_int_l ELSE de_dx_int_u = ( ( gg - aa ) * de_dx(k+1,j,i) + & ( gg - bb ) * de_dx(k+1,j,i+1) + & ( gg - cc ) * de_dx(k+1,j+1,i) + & ( gg - dd ) * de_dx(k+1,j+1,i+1) & ) / ( 3.0 * gg ) de_dx_int = de_dx_int_l + ( particles(n)%z - zu(k) ) / & dz * ( de_dx_int_u - de_dx_int_l ) ENDIF ! !-- Interpolate the TKE gradient along y de_dy_int_l = ( ( gg - aa ) * de_dy(k,j,i) + & ( gg - bb ) * de_dy(k,j,i+1) + & ( gg - cc ) * de_dy(k,j+1,i) + & ( gg - dd ) * de_dy(k,j+1,i+1) & ) / ( 3.0 * gg ) IF ( ( k+1 == nzt+1 ) .OR. ( k == nzb ) ) THEN de_dy_int = de_dy_int_l ELSE de_dy_int_u = ( ( gg - aa ) * de_dy(k+1,j,i) + & ( gg - bb ) * de_dy(k+1,j,i+1) + & ( gg - cc ) * de_dy(k+1,j+1,i) + & ( gg - dd ) * de_dy(k+1,j+1,i+1) & ) / ( 3.0 * gg ) de_dy_int = de_dy_int_l + ( particles(n)%z - zu(k) ) / & dz * ( de_dy_int_u - de_dy_int_l ) ENDIF ! !-- Interpolate the TKE gradient along z IF ( particles(n)%z < 0.5 * dz ) THEN de_dz_int = 0.0 ELSE de_dz_int_l = ( ( gg - aa ) * de_dz(k,j,i) + & ( gg - bb ) * de_dz(k,j,i+1) + & ( gg - cc ) * de_dz(k,j+1,i) + & ( gg - dd ) * de_dz(k,j+1,i+1) & ) / ( 3.0 * gg ) IF ( ( k+1 == nzt+1 ) .OR. ( k == nzb ) ) THEN de_dz_int = de_dz_int_l ELSE de_dz_int_u = ( ( gg - aa ) * de_dz(k+1,j,i) + & ( gg - bb ) * de_dz(k+1,j,i+1) + & ( gg - cc ) * de_dz(k+1,j+1,i) + & ( gg - dd ) * de_dz(k+1,j+1,i+1) & ) / ( 3.0 * gg ) de_dz_int = de_dz_int_l + ( particles(n)%z - zu(k) ) /& dz * ( de_dz_int_u - de_dz_int_l ) ENDIF ENDIF ! !-- Interpolate the dissipation of TKE diss_int_l = ( ( gg - aa ) * diss(k,j,i) + & ( gg - bb ) * diss(k,j,i+1) + & ( gg - cc ) * diss(k,j+1,i) + & ( gg - dd ) * diss(k,j+1,i+1) & ) / ( 3.0 * gg ) IF ( k+1 == nzt+1 ) THEN diss_int = diss_int_l ELSE diss_int_u = ( ( gg - aa ) * diss(k+1,j,i) + & ( gg - bb ) * diss(k+1,j,i+1) + & ( gg - cc ) * diss(k+1,j+1,i) + & ( gg - dd ) * diss(k+1,j+1,i+1) & ) / ( 3.0 * gg ) diss_int = diss_int_l + ( particles(n)%z - zu(k) ) / dz *& ( diss_int_u - diss_int_l ) ENDIF ELSE ! !-- If wall between gridpoint 1 and gridpoint 5, then !-- Neumann boundary condition has to be applied IF ( gp_outside_of_building(1) == 1 .AND. & gp_outside_of_building(5) == 0 ) THEN num_gp = num_gp + 1 location(num_gp,1) = i * dx + 0.5 * dx location(num_gp,2) = j * dy location(num_gp,3) = k * dz - 0.5 * dz ei(num_gp) = e(k,j,i) dissi(num_gp) = diss(k,j,i) de_dxi(num_gp) = 0.0 de_dyi(num_gp) = de_dy(k,j,i) de_dzi(num_gp) = de_dz(k,j,i) ENDIF IF ( gp_outside_of_building(5) == 1 .AND. & gp_outside_of_building(1) == 0 ) THEN num_gp = num_gp + 1 location(num_gp,1) = i * dx + 0.5 * dx location(num_gp,2) = j * dy location(num_gp,3) = k * dz - 0.5 * dz ei(num_gp) = e(k,j,i+1) dissi(num_gp) = diss(k,j,i+1) de_dxi(num_gp) = 0.0 de_dyi(num_gp) = de_dy(k,j,i+1) de_dzi(num_gp) = de_dz(k,j,i+1) ENDIF ! !-- If wall between gridpoint 5 and gridpoint 6, then !-- then Neumann boundary condition has to be applied IF ( gp_outside_of_building(5) == 1 .AND. & gp_outside_of_building(6) == 0 ) THEN num_gp = num_gp + 1 location(num_gp,1) = (i+1) * dx location(num_gp,2) = j * dy + 0.5 * dy location(num_gp,3) = k * dz - 0.5 * dz ei(num_gp) = e(k,j,i+1) dissi(num_gp) = diss(k,j,i+1) de_dxi(num_gp) = de_dx(k,j,i+1) de_dyi(num_gp) = 0.0 de_dzi(num_gp) = de_dz(k,j,i+1) ENDIF IF ( gp_outside_of_building(6) == 1 .AND. & gp_outside_of_building(5) == 0 ) THEN num_gp = num_gp + 1 location(num_gp,1) = (i+1) * dx location(num_gp,2) = j * dy + 0.5 * dy location(num_gp,3) = k * dz - 0.5 * dz ei(num_gp) = e(k,j+1,i+1) dissi(num_gp) = diss(k,j+1,i+1) de_dxi(num_gp) = de_dx(k,j+1,i+1) de_dyi(num_gp) = 0.0 de_dzi(num_gp) = de_dz(k,j+1,i+1) ENDIF ! !-- If wall between gridpoint 2 and gridpoint 6, then !-- Neumann boundary condition has to be applied IF ( gp_outside_of_building(2) == 1 .AND. & gp_outside_of_building(6) == 0 ) THEN num_gp = num_gp + 1 location(num_gp,1) = i * dx + 0.5 * dx location(num_gp,2) = (j+1) * dy location(num_gp,3) = k * dz - 0.5 * dz ei(num_gp) = e(k,j+1,i) dissi(num_gp) = diss(k,j+1,i) de_dxi(num_gp) = 0.0 de_dyi(num_gp) = de_dy(k,j+1,i) de_dzi(num_gp) = de_dz(k,j+1,i) ENDIF IF ( gp_outside_of_building(6) == 1 .AND. & gp_outside_of_building(2) == 0 ) THEN num_gp = num_gp + 1 location(num_gp,1) = i * dx + 0.5 * dx location(num_gp,2) = (j+1) * dy location(num_gp,3) = k * dz - 0.5 * dz ei(num_gp) = e(k,j+1,i+1) dissi(num_gp) = diss(k,j+1,i+1) de_dxi(num_gp) = 0.0 de_dyi(num_gp) = de_dy(k,j+1,i+1) de_dzi(num_gp) = de_dz(k,j+1,i+1) ENDIF ! !-- If wall between gridpoint 1 and gridpoint 2, then !-- Neumann boundary condition has to be applied IF ( gp_outside_of_building(1) == 1 .AND. & gp_outside_of_building(2) == 0 ) THEN num_gp = num_gp + 1 location(num_gp,1) = i * dx location(num_gp,2) = j * dy + 0.5 * dy location(num_gp,3) = k * dz - 0.5 * dz ei(num_gp) = e(k,j,i) dissi(num_gp) = diss(k,j,i) de_dxi(num_gp) = de_dx(k,j,i) de_dyi(num_gp) = 0.0 de_dzi(num_gp) = de_dz(k,j,i) ENDIF IF ( gp_outside_of_building(2) == 1 .AND. & gp_outside_of_building(1) == 0 ) THEN num_gp = num_gp + 1 location(num_gp,1) = i * dx location(num_gp,2) = j * dy + 0.5 * dy location(num_gp,3) = k * dz - 0.5 * dz ei(num_gp) = e(k,j+1,i) dissi(num_gp) = diss(k,j+1,i) de_dxi(num_gp) = de_dx(k,j+1,i) de_dyi(num_gp) = 0.0 de_dzi(num_gp) = de_dz(k,j+1,i) ENDIF ! !-- If wall between gridpoint 3 and gridpoint 7, then !-- Neumann boundary condition has to be applied IF ( gp_outside_of_building(3) == 1 .AND. & gp_outside_of_building(7) == 0 ) THEN num_gp = num_gp + 1 location(num_gp,1) = i * dx + 0.5 * dx location(num_gp,2) = j * dy location(num_gp,3) = k * dz + 0.5 * dz ei(num_gp) = e(k+1,j,i) dissi(num_gp) = diss(k+1,j,i) de_dxi(num_gp) = 0.0 de_dyi(num_gp) = de_dy(k+1,j,i) de_dzi(num_gp) = de_dz(k+1,j,i) ENDIF IF ( gp_outside_of_building(7) == 1 .AND. & gp_outside_of_building(3) == 0 ) THEN num_gp = num_gp + 1 location(num_gp,1) = i * dx + 0.5 * dx location(num_gp,2) = j * dy location(num_gp,3) = k * dz + 0.5 * dz ei(num_gp) = e(k+1,j,i+1) dissi(num_gp) = diss(k+1,j,i+1) de_dxi(num_gp) = 0.0 de_dyi(num_gp) = de_dy(k+1,j,i+1) de_dzi(num_gp) = de_dz(k+1,j,i+1) ENDIF ! !-- If wall between gridpoint 7 and gridpoint 8, then !-- Neumann boundary condition has to be applied IF ( gp_outside_of_building(7) == 1 .AND. & gp_outside_of_building(8) == 0 ) THEN num_gp = num_gp + 1 location(num_gp,1) = (i+1) * dx location(num_gp,2) = j * dy + 0.5 * dy location(num_gp,3) = k * dz + 0.5 * dz ei(num_gp) = e(k+1,j,i+1) dissi(num_gp) = diss(k+1,j,i+1) de_dxi(num_gp) = de_dx(k+1,j,i+1) de_dyi(num_gp) = 0.0 de_dzi(num_gp) = de_dz(k+1,j,i+1) ENDIF IF ( gp_outside_of_building(8) == 1 .AND. & gp_outside_of_building(7) == 0 ) THEN num_gp = num_gp + 1 location(num_gp,1) = (i+1) * dx location(num_gp,2) = j * dy + 0.5 * dy location(num_gp,3) = k * dz + 0.5 * dz ei(num_gp) = e(k+1,j+1,i+1) dissi(num_gp) = diss(k+1,j+1,i+1) de_dxi(num_gp) = de_dx(k+1,j+1,i+1) de_dyi(num_gp) = 0.0 de_dzi(num_gp) = de_dz(k+1,j+1,i+1) ENDIF ! !-- If wall between gridpoint 4 and gridpoint 8, then !-- Neumann boundary condition has to be applied IF ( gp_outside_of_building(4) == 1 .AND. & gp_outside_of_building(8) == 0 ) THEN num_gp = num_gp + 1 location(num_gp,1) = i * dx + 0.5 * dx location(num_gp,2) = (j+1) * dy location(num_gp,3) = k * dz + 0.5 * dz ei(num_gp) = e(k+1,j+1,i) dissi(num_gp) = diss(k+1,j+1,i) de_dxi(num_gp) = 0.0 de_dyi(num_gp) = de_dy(k+1,j+1,i) de_dzi(num_gp) = de_dz(k+1,j+1,i) ENDIF IF ( gp_outside_of_building(8) == 1 .AND. & gp_outside_of_building(4) == 0 ) THEN num_gp = num_gp + 1 location(num_gp,1) = i * dx + 0.5 * dx location(num_gp,2) = (j+1) * dy location(num_gp,3) = k * dz + 0.5 * dz ei(num_gp) = e(k+1,j+1,i+1) dissi(num_gp) = diss(k+1,j+1,i+1) de_dxi(num_gp) = 0.0 de_dyi(num_gp) = de_dy(k+1,j+1,i+1) de_dzi(num_gp) = de_dz(k+1,j+1,i+1) ENDIF ! !-- If wall between gridpoint 3 and gridpoint 4, then !-- Neumann boundary condition has to be applied IF ( gp_outside_of_building(3) == 1 .AND. & gp_outside_of_building(4) == 0 ) THEN num_gp = num_gp + 1 location(num_gp,1) = i * dx location(num_gp,2) = j * dy + 0.5 * dy location(num_gp,3) = k * dz + 0.5 * dz ei(num_gp) = e(k+1,j,i) dissi(num_gp) = diss(k+1,j,i) de_dxi(num_gp) = de_dx(k+1,j,i) de_dyi(num_gp) = 0.0 de_dzi(num_gp) = de_dz(k+1,j,i) ENDIF IF ( gp_outside_of_building(4) == 1 .AND. & gp_outside_of_building(3) == 0 ) THEN num_gp = num_gp + 1 location(num_gp,1) = i * dx location(num_gp,2) = j * dy + 0.5 * dy location(num_gp,3) = k * dz + 0.5 * dz ei(num_gp) = e(k+1,j+1,i) dissi(num_gp) = diss(k+1,j+1,i) de_dxi(num_gp) = de_dx(k+1,j+1,i) de_dyi(num_gp) = 0.0 de_dzi(num_gp) = de_dz(k+1,j+1,i) ENDIF ! !-- If wall between gridpoint 1 and gridpoint 3, then !-- Neumann boundary condition has to be applied !-- (only one case as only building beneath is possible) IF ( gp_outside_of_building(1) == 0 .AND. & gp_outside_of_building(3) == 1 ) THEN num_gp = num_gp + 1 location(num_gp,1) = i * dx location(num_gp,2) = j * dy location(num_gp,3) = k * dz ei(num_gp) = e(k+1,j,i) dissi(num_gp) = diss(k+1,j,i) de_dxi(num_gp) = de_dx(k+1,j,i) de_dyi(num_gp) = de_dy(k+1,j,i) de_dzi(num_gp) = 0.0 ENDIF ! !-- If wall between gridpoint 5 and gridpoint 7, then !-- Neumann boundary condition has to be applied !-- (only one case as only building beneath is possible) IF ( gp_outside_of_building(5) == 0 .AND. & gp_outside_of_building(7) == 1 ) THEN num_gp = num_gp + 1 location(num_gp,1) = (i+1) * dx location(num_gp,2) = j * dy location(num_gp,3) = k * dz ei(num_gp) = e(k+1,j,i+1) dissi(num_gp) = diss(k+1,j,i+1) de_dxi(num_gp) = de_dx(k+1,j,i+1) de_dyi(num_gp) = de_dy(k+1,j,i+1) de_dzi(num_gp) = 0.0 ENDIF ! !-- If wall between gridpoint 2 and gridpoint 4, then !-- Neumann boundary condition has to be applied !-- (only one case as only building beneath is possible) IF ( gp_outside_of_building(2) == 0 .AND. & gp_outside_of_building(4) == 1 ) THEN num_gp = num_gp + 1 location(num_gp,1) = i * dx location(num_gp,2) = (j+1) * dy location(num_gp,3) = k * dz ei(num_gp) = e(k+1,j+1,i) dissi(num_gp) = diss(k+1,j+1,i) de_dxi(num_gp) = de_dx(k+1,j+1,i) de_dyi(num_gp) = de_dy(k+1,j+1,i) de_dzi(num_gp) = 0.0 ENDIF ! !-- If wall between gridpoint 6 and gridpoint 8, then !-- Neumann boundary condition has to be applied !-- (only one case as only building beneath is possible) IF ( gp_outside_of_building(6) == 0 .AND. & gp_outside_of_building(8) == 1 ) THEN num_gp = num_gp + 1 location(num_gp,1) = (i+1) * dx location(num_gp,2) = (j+1) * dy location(num_gp,3) = k * dz ei(num_gp) = e(k+1,j+1,i+1) dissi(num_gp) = diss(k+1,j+1,i+1) de_dxi(num_gp) = de_dx(k+1,j+1,i+1) de_dyi(num_gp) = de_dy(k+1,j+1,i+1) de_dzi(num_gp) = 0.0 ENDIF ! !-- Carry out the interpolation IF ( num_gp == 1 ) THEN ! !-- If only one of the gridpoints is situated outside of the !-- building, it follows that the values at the particle !-- location are the same as the gridpoint values e_int = ei(num_gp) diss_int = dissi(num_gp) de_dx_int = de_dxi(num_gp) de_dy_int = de_dyi(num_gp) de_dz_int = de_dzi(num_gp) ELSE IF ( num_gp > 1 ) THEN d_sum = 0.0 ! !-- Evaluation of the distances between the gridpoints !-- contributing to the interpolated values, and the particle !-- location DO agp = 1, num_gp d_gp_pl(agp) = ( particles(n)%x-location(agp,1) )**2 & + ( particles(n)%y-location(agp,2) )**2 & + ( particles(n)%z-location(agp,3) )**2 d_sum = d_sum + d_gp_pl(agp) ENDDO ! !-- Finally the interpolation can be carried out e_int = 0.0 diss_int = 0.0 de_dx_int = 0.0 de_dy_int = 0.0 de_dz_int = 0.0 DO agp = 1, num_gp e_int = e_int + ( d_sum - d_gp_pl(agp) ) * & ei(agp) / ( (num_gp-1) * d_sum ) diss_int = diss_int + ( d_sum - d_gp_pl(agp) ) * & dissi(agp) / ( (num_gp-1) * d_sum ) de_dx_int = de_dx_int + ( d_sum - d_gp_pl(agp) ) * & de_dxi(agp) / ( (num_gp-1) * d_sum ) de_dy_int = de_dy_int + ( d_sum - d_gp_pl(agp) ) * & de_dyi(agp) / ( (num_gp-1) * d_sum ) de_dz_int = de_dz_int + ( d_sum - d_gp_pl(agp) ) * & de_dzi(agp) / ( (num_gp-1) * d_sum ) ENDDO ENDIF ENDIF ENDIF ! !-- Vertically interpolate the horizontally averaged SGS TKE and !-- resolved-scale velocity variances and use the interpolated values !-- to calculate the coefficient fs, which is a measure of the ratio !-- of the subgrid-scale turbulent kinetic energy to the total amount !-- of turbulent kinetic energy. IF ( k == 0 ) THEN e_mean_int = hom(0,1,8,0) ELSE e_mean_int = hom(k,1,8,0) + & ( hom(k+1,1,8,0) - hom(k,1,8,0) ) / & ( zu(k+1) - zu(k) ) * & ( particles(n)%z - zu(k) ) ENDIF kw = particles(n)%z / dz IF ( k == 0 ) THEN aa = hom(k+1,1,30,0) * ( particles(n)%z / & ( 0.5 * ( zu(k+1) - zu(k) ) ) ) bb = hom(k+1,1,31,0) * ( particles(n)%z / & ( 0.5 * ( zu(k+1) - zu(k) ) ) ) cc = hom(kw+1,1,32,0) * ( particles(n)%z / & ( 1.0 * ( zw(kw+1) - zw(kw) ) ) ) ELSE aa = hom(k,1,30,0) + ( hom(k+1,1,30,0) - hom(k,1,30,0) ) * & ( ( particles(n)%z - zu(k) ) / ( zu(k+1) - zu(k) ) ) bb = hom(k,1,31,0) + ( hom(k+1,1,31,0) - hom(k,1,31,0) ) * & ( ( particles(n)%z - zu(k) ) / ( zu(k+1) - zu(k) ) ) cc = hom(kw,1,32,0) + ( hom(kw+1,1,32,0)-hom(kw,1,32,0) ) *& ( ( particles(n)%z - zw(kw) ) / ( zw(kw+1)-zw(kw) ) ) ENDIF vv_int = ( 1.0 / 3.0 ) * ( aa + bb + cc ) fs_int = ( 2.0 / 3.0 ) * e_mean_int / & ( vv_int + ( 2.0 / 3.0 ) * e_mean_int ) ! !-- Calculate the Lagrangian timescale according to Weil et al. (2004). lagr_timescale = ( 4.0 * e_int ) / & ( 3.0 * fs_int * c_0 * diss_int ) ! !-- Calculate the next particle timestep. dt_gap is the time needed to !-- complete the current LES timestep. dt_gap = dt_3d - particles(n)%dt_sum dt_particle = MIN( dt_3d, 0.025 * lagr_timescale, dt_gap ) ! !-- The particle timestep should not be too small in order to prevent !-- the number of particle timesteps of getting too large IF ( dt_particle < dt_min_part .AND. dt_min_part < dt_gap ) & THEN dt_particle = dt_min_part ENDIF ! !-- Calculate the SGS velocity components IF ( particles(n)%age == 0.0 ) THEN ! !-- For new particles the SGS components are derived from the SGS !-- TKE. Limit the Gaussian random number to the interval !-- [-5.0*sigma, 5.0*sigma] in order to prevent the SGS velocities !-- from becoming unrealistically large. particles(n)%rvar1 = SQRT( 2.0 * sgs_wfu_part * e_int ) * & ( random_gauss( iran_part, 5.0 ) - 1.0 ) particles(n)%rvar2 = SQRT( 2.0 * sgs_wfv_part * e_int ) * & ( random_gauss( iran_part, 5.0 ) - 1.0 ) particles(n)%rvar3 = SQRT( 2.0 * sgs_wfw_part * e_int ) * & ( random_gauss( iran_part, 5.0 ) - 1.0 ) ELSE ! !-- Restriction of the size of the new timestep: compared to the !-- previous timestep the increase must not exceed 200% dt_particle_m = particles(n)%age - particles(n)%age_m IF ( dt_particle > 2.0 * dt_particle_m ) THEN dt_particle = 2.0 * dt_particle_m ENDIF ! !-- For old particles the SGS components are correlated with the !-- values from the previous timestep. Random numbers have also to !-- be limited (see above). !-- As negative values for the subgrid TKE are not allowed, the !-- change of the subgrid TKE with time cannot be smaller than !-- -e_int/dt_particle. This value is used as a lower boundary !-- value for the change of TKE de_dt_min = - e_int / dt_particle de_dt = ( e_int - particles(n)%e_m ) / dt_particle_m IF ( de_dt < de_dt_min ) THEN de_dt = de_dt_min ENDIF particles(n)%rvar1 = particles(n)%rvar1 - fs_int * c_0 * & diss_int * particles(n)%rvar1 * dt_particle /& ( 4.0 * sgs_wfu_part * e_int ) + & ( 2.0 * sgs_wfu_part * de_dt * & particles(n)%rvar1 / & ( 2.0 * sgs_wfu_part * e_int ) + de_dx_int & ) * dt_particle / 2.0 + & SQRT( fs_int * c_0 * diss_int ) * & ( random_gauss( iran_part, 5.0 ) - 1.0 ) * & SQRT( dt_particle ) particles(n)%rvar2 = particles(n)%rvar2 - fs_int * c_0 * & diss_int * particles(n)%rvar2 * dt_particle /& ( 4.0 * sgs_wfv_part * e_int ) + & ( 2.0 * sgs_wfv_part * de_dt * & particles(n)%rvar2 / & ( 2.0 * sgs_wfv_part * e_int ) + de_dy_int & ) * dt_particle / 2.0 + & SQRT( fs_int * c_0 * diss_int ) * & ( random_gauss( iran_part, 5.0 ) - 1.0 ) * & SQRT( dt_particle ) particles(n)%rvar3 = particles(n)%rvar3 - fs_int * c_0 * & diss_int * particles(n)%rvar3 * dt_particle /& ( 4.0 * sgs_wfw_part * e_int ) + & ( 2.0 * sgs_wfw_part * de_dt * & particles(n)%rvar3 / & ( 2.0 * sgs_wfw_part * e_int ) + de_dz_int & ) * dt_particle / 2.0 + & SQRT( fs_int * c_0 * diss_int ) * & ( random_gauss( iran_part, 5.0 ) - 1.0 ) * & SQRT( dt_particle ) ENDIF u_int = u_int + particles(n)%rvar1 v_int = v_int + particles(n)%rvar2 w_int = w_int + particles(n)%rvar3 ! !-- Store the SGS TKE of the current timelevel which is needed for !-- for calculating the SGS particle velocities at the next timestep particles(n)%e_m = e_int ELSE ! !-- If no SGS velocities are used, only the particle timestep has to !-- be set dt_particle = dt_3d ENDIF ! !-- Store the old age of the particle ( needed to prevent that a !-- particle crosses several PEs during one timestep, and for the !-- evaluation of the subgrid particle velocity fluctuations ) particles(n)%age_m = particles(n)%age ! !-- Particle advection IF ( particle_groups(particles(n)%group)%density_ratio == 0.0 ) THEN ! !-- Pure passive transport (without particle inertia) particles(n)%x = particles(n)%x + u_int * dt_particle particles(n)%y = particles(n)%y + v_int * dt_particle particles(n)%z = particles(n)%z + w_int * dt_particle particles(n)%speed_x = u_int particles(n)%speed_y = v_int particles(n)%speed_z = w_int ELSE ! !-- Transport of particles with inertia particles(n)%x = particles(n)%x + particles(n)%speed_x * & dt_particle particles(n)%y = particles(n)%y + particles(n)%speed_y * & dt_particle particles(n)%z = particles(n)%z + particles(n)%speed_z * & dt_particle ! !-- Update of the particle velocity dens_ratio = particle_groups(particles(n)%group)%density_ratio IF ( cloud_droplets ) THEN exp_arg = 4.5 * dens_ratio * molecular_viscosity / & ( particles(n)%radius )**2 * & ( 1.0 + 0.15 * ( 2.0 * particles(n)%radius * & SQRT( ( u_int - particles(n)%speed_x )**2 + & ( v_int - particles(n)%speed_y )**2 + & ( w_int - particles(n)%speed_z )**2 ) / & molecular_viscosity )**0.687 & ) exp_term = EXP( -exp_arg * dt_particle ) ELSEIF ( use_sgs_for_particles ) THEN exp_arg = particle_groups(particles(n)%group)%exp_arg exp_term = EXP( -exp_arg * dt_particle ) ELSE exp_arg = particle_groups(particles(n)%group)%exp_arg exp_term = particle_groups(particles(n)%group)%exp_term ENDIF particles(n)%speed_x = particles(n)%speed_x * exp_term + & u_int * ( 1.0 - exp_term ) particles(n)%speed_y = particles(n)%speed_y * exp_term + & v_int * ( 1.0 - exp_term ) particles(n)%speed_z = particles(n)%speed_z * exp_term + & ( w_int - ( 1.0 - dens_ratio ) * g / exp_arg )& * ( 1.0 - exp_term ) ENDIF ! !-- Increment the particle age and the total time that the particle !-- has advanced within the particle timestep procedure particles(n)%age = particles(n)%age + dt_particle particles(n)%dt_sum = particles(n)%dt_sum + dt_particle ! !-- Check whether there is still a particle that has not yet completed !-- the total LES timestep IF ( ( dt_3d - particles(n)%dt_sum ) > 1E-8 ) THEN dt_3d_reached_l = .FALSE. ENDIF ENDDO END SUBROUTINE lpm_advec