[1682] | 1 | !> @file lpm_init.f90 |
---|
[1036] | 2 | !--------------------------------------------------------------------------------! |
---|
| 3 | ! This file is part of PALM. |
---|
| 4 | ! |
---|
| 5 | ! PALM is free software: you can redistribute it and/or modify it under the terms |
---|
| 6 | ! of the GNU General Public License as published by the Free Software Foundation, |
---|
| 7 | ! either version 3 of the License, or (at your option) any later version. |
---|
| 8 | ! |
---|
| 9 | ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY |
---|
| 10 | ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR |
---|
| 11 | ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. |
---|
| 12 | ! |
---|
| 13 | ! You should have received a copy of the GNU General Public License along with |
---|
| 14 | ! PALM. If not, see <http://www.gnu.org/licenses/>. |
---|
| 15 | ! |
---|
[1310] | 16 | ! Copyright 1997-2014 Leibniz Universitaet Hannover |
---|
[1036] | 17 | !--------------------------------------------------------------------------------! |
---|
| 18 | ! |
---|
[254] | 19 | ! Current revisions: |
---|
[1] | 20 | ! ----------------- |
---|
[1682] | 21 | ! Code annotations made doxygen readable |
---|
[1576] | 22 | ! |
---|
[1321] | 23 | ! Former revisions: |
---|
| 24 | ! ----------------- |
---|
| 25 | ! $Id: lpm_init.f90 1682 2015-10-07 23:56:08Z knoop $ |
---|
| 26 | ! |
---|
[1576] | 27 | ! 1575 2015-03-27 09:56:27Z raasch |
---|
| 28 | ! initial vertical particle position is allowed to follow the topography |
---|
| 29 | ! |
---|
[1360] | 30 | ! 1359 2014-04-11 17:15:14Z hoffmann |
---|
| 31 | ! New particle structure integrated. |
---|
| 32 | ! Kind definition added to all floating point numbers. |
---|
| 33 | ! lpm_init changed form a subroutine to a module. |
---|
| 34 | ! |
---|
[1329] | 35 | ! 1327 2014-03-21 11:00:16Z raasch |
---|
| 36 | ! -netcdf_output |
---|
| 37 | ! |
---|
[1323] | 38 | ! 1322 2014-03-20 16:38:49Z raasch |
---|
| 39 | ! REAL functions provided with KIND-attribute |
---|
| 40 | ! |
---|
[1321] | 41 | ! 1320 2014-03-20 08:40:49Z raasch |
---|
[1320] | 42 | ! ONLY-attribute added to USE-statements, |
---|
| 43 | ! kind-parameters added to all INTEGER and REAL declaration statements, |
---|
| 44 | ! kinds are defined in new module kinds, |
---|
| 45 | ! revision history before 2012 removed, |
---|
| 46 | ! comment fields (!:) to be used for variable explanations added to |
---|
| 47 | ! all variable declaration statements |
---|
| 48 | ! bugfix: #if defined( __parallel ) added |
---|
[850] | 49 | ! |
---|
[1315] | 50 | ! 1314 2014-03-14 18:25:17Z suehring |
---|
| 51 | ! Vertical logarithmic interpolation of horizontal particle speed for particles |
---|
| 52 | ! between roughness height and first vertical grid level. |
---|
| 53 | ! |
---|
[1093] | 54 | ! 1092 2013-02-02 11:24:22Z raasch |
---|
| 55 | ! unused variables removed |
---|
| 56 | ! |
---|
[1037] | 57 | ! 1036 2012-10-22 13:43:42Z raasch |
---|
| 58 | ! code put under GPL (PALM 3.9) |
---|
| 59 | ! |
---|
[850] | 60 | ! 849 2012-03-15 10:35:09Z raasch |
---|
[849] | 61 | ! routine renamed: init_particles -> lpm_init |
---|
| 62 | ! de_dx, de_dy, de_dz are allocated here (instead of automatic arrays in |
---|
| 63 | ! advec_particles), |
---|
| 64 | ! sort_particles renamed lpm_sort_arrays, user_init_particles renamed lpm_init |
---|
[392] | 65 | ! |
---|
[829] | 66 | ! 828 2012-02-21 12:00:36Z raasch |
---|
| 67 | ! call of init_kernels, particle feature color renamed class |
---|
| 68 | ! |
---|
[826] | 69 | ! 824 2012-02-17 09:09:57Z raasch |
---|
| 70 | ! particle attributes speed_x|y|z_sgs renamed rvar1|2|3, |
---|
| 71 | ! array particles implemented as pointer |
---|
| 72 | ! |
---|
[668] | 73 | ! 667 2010-12-23 12:06:00Z suehring/gryschka |
---|
| 74 | ! nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng for allocation |
---|
| 75 | ! of arrays. |
---|
| 76 | ! |
---|
[1] | 77 | ! Revision 1.1 1999/11/25 16:22:38 raasch |
---|
| 78 | ! Initial revision |
---|
| 79 | ! |
---|
| 80 | ! |
---|
| 81 | ! Description: |
---|
| 82 | ! ------------ |
---|
[1682] | 83 | !> This routine initializes a set of particles and their attributes (position, |
---|
| 84 | !> radius, ..) which are used by the Lagrangian particle model (see lpm). |
---|
[1] | 85 | !------------------------------------------------------------------------------! |
---|
[1682] | 86 | MODULE lpm_init_mod |
---|
| 87 | |
---|
[1] | 88 | |
---|
[1320] | 89 | USE arrays_3d, & |
---|
| 90 | ONLY: de_dx, de_dy, de_dz, zu, zw, z0 |
---|
| 91 | |
---|
| 92 | USE cloud_parameters, & |
---|
| 93 | ONLY: curvature_solution_effects |
---|
| 94 | |
---|
| 95 | USE control_parameters, & |
---|
[1359] | 96 | ONLY: cloud_droplets, current_timestep_number, dz, & |
---|
| 97 | initializing_actions, message_string, netcdf_data_format, & |
---|
| 98 | ocean, prandtl_layer, simulated_time |
---|
[1320] | 99 | |
---|
| 100 | USE dvrp_variables, & |
---|
| 101 | ONLY: particle_color, particle_dvrpsize |
---|
| 102 | |
---|
| 103 | USE grid_variables, & |
---|
[1359] | 104 | ONLY: ddx, dx, ddy, dy |
---|
[1320] | 105 | |
---|
| 106 | USE indices, & |
---|
[1575] | 107 | ONLY: nx, nxl, nxlg, nxrg, nxr, ny, nyn, nys, nyng, nysg, nz, nzb, & |
---|
| 108 | nzb_w_inner, nzt |
---|
[1320] | 109 | |
---|
| 110 | USE kinds |
---|
| 111 | |
---|
| 112 | USE lpm_collision_kernels_mod, & |
---|
| 113 | ONLY: init_kernels |
---|
| 114 | |
---|
| 115 | USE particle_attributes, & |
---|
[1359] | 116 | ONLY: alloc_factor, bc_par_b, bc_par_lr, bc_par_ns, bc_par_t, & |
---|
| 117 | block_offset, block_offset_def, collision_kernel, & |
---|
| 118 | density_ratio, dvrp_psize, grid_particles, & |
---|
| 119 | initial_weighting_factor, ibc_par_b, ibc_par_lr, ibc_par_ns, & |
---|
| 120 | ibc_par_t, iran_part, log_z_z0, & |
---|
| 121 | max_number_of_particle_groups, maximum_number_of_particles, & |
---|
| 122 | maximum_number_of_tailpoints, maximum_number_of_tails, & |
---|
| 123 | minimum_tailpoint_distance, min_nr_particle, & |
---|
| 124 | mpi_particle_type, new_tail_id, & |
---|
[1320] | 125 | number_of_initial_tails, number_of_particles, & |
---|
| 126 | number_of_particle_groups, number_of_sublayers, & |
---|
[1359] | 127 | number_of_tails, offset_ocean_nzt, offset_ocean_nzt_m1, & |
---|
| 128 | particles, particle_advection_start, particle_groups, & |
---|
| 129 | particle_groups_type, particles_per_point, & |
---|
[1320] | 130 | particle_tail_coordinates, particle_type, pdx, pdy, pdz, & |
---|
[1359] | 131 | prt_count, psb, psl, psn, psr, pss, pst, & |
---|
[1320] | 132 | radius, random_start_position, read_particles_from_restartfile,& |
---|
[1575] | 133 | seed_follows_topography, skip_particles_for_tail, sort_count, & |
---|
[1359] | 134 | tail_mask, total_number_of_particles, total_number_of_tails, & |
---|
[1320] | 135 | use_particle_tails, use_sgs_for_particles, & |
---|
[1359] | 136 | write_particle_statistics, uniform_particles, zero_particle, & |
---|
| 137 | z0_av_global |
---|
[1320] | 138 | |
---|
[1] | 139 | USE pegrid |
---|
| 140 | |
---|
[1320] | 141 | USE random_function_mod, & |
---|
| 142 | ONLY: random_function |
---|
[1] | 143 | |
---|
[1359] | 144 | IMPLICIT NONE |
---|
[1320] | 145 | |
---|
[1359] | 146 | PRIVATE |
---|
| 147 | |
---|
[1682] | 148 | INTEGER(iwp), PARAMETER :: PHASE_INIT = 1 !< |
---|
| 149 | INTEGER(iwp), PARAMETER, PUBLIC :: PHASE_RELEASE = 2 !< |
---|
[1359] | 150 | |
---|
| 151 | INTERFACE lpm_init |
---|
| 152 | MODULE PROCEDURE lpm_init |
---|
| 153 | END INTERFACE lpm_init |
---|
| 154 | |
---|
| 155 | INTERFACE lpm_create_particle |
---|
| 156 | MODULE PROCEDURE lpm_create_particle |
---|
| 157 | END INTERFACE lpm_create_particle |
---|
| 158 | |
---|
| 159 | PUBLIC lpm_init, lpm_create_particle |
---|
| 160 | |
---|
| 161 | CONTAINS |
---|
| 162 | |
---|
[1682] | 163 | !------------------------------------------------------------------------------! |
---|
| 164 | ! Description: |
---|
| 165 | ! ------------ |
---|
| 166 | !> @todo Missing subroutine description. |
---|
| 167 | !------------------------------------------------------------------------------! |
---|
[1359] | 168 | SUBROUTINE lpm_init |
---|
| 169 | |
---|
| 170 | USE lpm_collision_kernels_mod, & |
---|
| 171 | ONLY: init_kernels |
---|
| 172 | |
---|
[1] | 173 | IMPLICIT NONE |
---|
| 174 | |
---|
[1682] | 175 | INTEGER(iwp) :: i !< |
---|
| 176 | INTEGER(iwp) :: ip !< |
---|
| 177 | INTEGER(iwp) :: j !< |
---|
| 178 | INTEGER(iwp) :: jp !< |
---|
| 179 | INTEGER(iwp) :: k !< |
---|
| 180 | INTEGER(iwp) :: kp !< |
---|
| 181 | INTEGER(iwp) :: n !< |
---|
| 182 | INTEGER(iwp) :: nn !< |
---|
[1320] | 183 | |
---|
[1] | 184 | #if defined( __parallel ) |
---|
[1682] | 185 | INTEGER(iwp), DIMENSION(3) :: blocklengths !< |
---|
| 186 | INTEGER(iwp), DIMENSION(3) :: displacements !< |
---|
| 187 | INTEGER(iwp), DIMENSION(3) :: types !< |
---|
[1] | 188 | #endif |
---|
[1682] | 189 | LOGICAL :: uniform_particles_l !< |
---|
[1] | 190 | |
---|
[1682] | 191 | REAL(wp) :: height_int !< |
---|
| 192 | REAL(wp) :: height_p !< |
---|
| 193 | REAL(wp) :: z_p !< |
---|
| 194 | REAL(wp) :: z0_av_local !< |
---|
[1] | 195 | |
---|
| 196 | #if defined( __parallel ) |
---|
| 197 | ! |
---|
| 198 | !-- Define MPI derived datatype for FORTRAN datatype particle_type (see module |
---|
[82] | 199 | !-- particle_attributes). Integer length is 4 byte, Real is 8 byte |
---|
[1359] | 200 | #if defined( __twocachelines ) |
---|
| 201 | blocklengths(1) = 7; blocklengths(2) = 18; blocklengths(3) = 1 |
---|
| 202 | displacements(1) = 0; displacements(2) = 64; displacements(3) = 128 |
---|
[82] | 203 | |
---|
[1359] | 204 | types(1) = MPI_REAL ! 64 bit words |
---|
| 205 | types(2) = MPI_INTEGER ! 32 Bit words |
---|
| 206 | types(3) = MPI_UB |
---|
| 207 | #else |
---|
| 208 | blocklengths(1) = 19; blocklengths(2) = 6; blocklengths(3) = 1 |
---|
| 209 | displacements(1) = 0; displacements(2) = 152; displacements(3) = 176 |
---|
| 210 | |
---|
[1] | 211 | types(1) = MPI_REAL |
---|
| 212 | types(2) = MPI_INTEGER |
---|
| 213 | types(3) = MPI_UB |
---|
[1359] | 214 | #endif |
---|
[1] | 215 | CALL MPI_TYPE_STRUCT( 3, blocklengths, displacements, types, & |
---|
| 216 | mpi_particle_type, ierr ) |
---|
| 217 | CALL MPI_TYPE_COMMIT( mpi_particle_type, ierr ) |
---|
| 218 | #endif |
---|
| 219 | |
---|
| 220 | ! |
---|
[150] | 221 | !-- In case of oceans runs, the vertical index calculations need an offset, |
---|
| 222 | !-- because otherwise the k indices will become negative |
---|
| 223 | IF ( ocean ) THEN |
---|
| 224 | offset_ocean_nzt = nzt |
---|
| 225 | offset_ocean_nzt_m1 = nzt - 1 |
---|
| 226 | ENDIF |
---|
| 227 | |
---|
[1359] | 228 | ! |
---|
| 229 | !-- Define block offsets for dividing a gridcell in 8 sub cells |
---|
[150] | 230 | |
---|
[1359] | 231 | block_offset(0) = block_offset_def (-1,-1,-1) |
---|
| 232 | block_offset(1) = block_offset_def (-1,-1, 0) |
---|
| 233 | block_offset(2) = block_offset_def (-1, 0,-1) |
---|
| 234 | block_offset(3) = block_offset_def (-1, 0, 0) |
---|
| 235 | block_offset(4) = block_offset_def ( 0,-1,-1) |
---|
| 236 | block_offset(5) = block_offset_def ( 0,-1, 0) |
---|
| 237 | block_offset(6) = block_offset_def ( 0, 0,-1) |
---|
| 238 | block_offset(7) = block_offset_def ( 0, 0, 0) |
---|
[150] | 239 | ! |
---|
[1] | 240 | !-- Check the number of particle groups. |
---|
| 241 | IF ( number_of_particle_groups > max_number_of_particle_groups ) THEN |
---|
[274] | 242 | WRITE( message_string, * ) 'max_number_of_particle_groups =', & |
---|
| 243 | max_number_of_particle_groups , & |
---|
[254] | 244 | '&number_of_particle_groups reset to ', & |
---|
| 245 | max_number_of_particle_groups |
---|
[849] | 246 | CALL message( 'lpm_init', 'PA0213', 0, 1, 0, 6, 0 ) |
---|
[1] | 247 | number_of_particle_groups = max_number_of_particle_groups |
---|
| 248 | ENDIF |
---|
| 249 | |
---|
| 250 | ! |
---|
| 251 | !-- Set default start positions, if necessary |
---|
[1359] | 252 | IF ( psl(1) == 9999999.9_wp ) psl(1) = -0.5_wp * dx |
---|
| 253 | IF ( psr(1) == 9999999.9_wp ) psr(1) = ( nx + 0.5_wp ) * dx |
---|
| 254 | IF ( pss(1) == 9999999.9_wp ) pss(1) = -0.5_wp * dy |
---|
| 255 | IF ( psn(1) == 9999999.9_wp ) psn(1) = ( ny + 0.5_wp ) * dy |
---|
| 256 | IF ( psb(1) == 9999999.9_wp ) psb(1) = zu(nz/2) |
---|
| 257 | IF ( pst(1) == 9999999.9_wp ) pst(1) = psb(1) |
---|
[1] | 258 | |
---|
[1359] | 259 | IF ( pdx(1) == 9999999.9_wp .OR. pdx(1) == 0.0_wp ) pdx(1) = dx |
---|
| 260 | IF ( pdy(1) == 9999999.9_wp .OR. pdy(1) == 0.0_wp ) pdy(1) = dy |
---|
| 261 | IF ( pdz(1) == 9999999.9_wp .OR. pdz(1) == 0.0_wp ) pdz(1) = zu(2) - zu(1) |
---|
[1] | 262 | |
---|
| 263 | DO j = 2, number_of_particle_groups |
---|
[1359] | 264 | IF ( psl(j) == 9999999.9_wp ) psl(j) = psl(j-1) |
---|
| 265 | IF ( psr(j) == 9999999.9_wp ) psr(j) = psr(j-1) |
---|
| 266 | IF ( pss(j) == 9999999.9_wp ) pss(j) = pss(j-1) |
---|
| 267 | IF ( psn(j) == 9999999.9_wp ) psn(j) = psn(j-1) |
---|
| 268 | IF ( psb(j) == 9999999.9_wp ) psb(j) = psb(j-1) |
---|
| 269 | IF ( pst(j) == 9999999.9_wp ) pst(j) = pst(j-1) |
---|
| 270 | IF ( pdx(j) == 9999999.9_wp .OR. pdx(j) == 0.0_wp ) pdx(j) = pdx(j-1) |
---|
| 271 | IF ( pdy(j) == 9999999.9_wp .OR. pdy(j) == 0.0_wp ) pdy(j) = pdy(j-1) |
---|
| 272 | IF ( pdz(j) == 9999999.9_wp .OR. pdz(j) == 0.0_wp ) pdz(j) = pdz(j-1) |
---|
[1] | 273 | ENDDO |
---|
| 274 | |
---|
| 275 | ! |
---|
[849] | 276 | !-- Allocate arrays required for calculating particle SGS velocities |
---|
| 277 | IF ( use_sgs_for_particles ) THEN |
---|
| 278 | ALLOCATE( de_dx(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 279 | de_dy(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 280 | de_dz(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
| 281 | ENDIF |
---|
| 282 | |
---|
| 283 | ! |
---|
[1314] | 284 | !-- Allocate array required for logarithmic vertical interpolation of |
---|
| 285 | !-- horizontal particle velocities between the surface and the first vertical |
---|
| 286 | !-- grid level. In order to avoid repeated CPU cost-intensive CALLS of |
---|
| 287 | !-- intrinsic FORTRAN procedure LOG(z/z0), LOG(z/z0) is precalculated for |
---|
| 288 | !-- several heights. Splitting into 20 sublayers turned out to be sufficient. |
---|
| 289 | !-- To obtain exact height levels of particles, linear interpolation is applied |
---|
| 290 | !-- (see lpm_advec.f90). |
---|
| 291 | IF ( prandtl_layer ) THEN |
---|
| 292 | |
---|
| 293 | ALLOCATE ( log_z_z0(0:number_of_sublayers) ) |
---|
| 294 | z_p = zu(nzb+1) - zw(nzb) |
---|
| 295 | |
---|
| 296 | ! |
---|
| 297 | !-- Calculate horizontal mean value of z0 used for logartihmic |
---|
| 298 | !-- interpolation. Note: this is not exact for heterogeneous z0. |
---|
| 299 | !-- However, sensitivity studies showed that the effect is |
---|
| 300 | !-- negligible. |
---|
| 301 | z0_av_local = SUM( z0(nys:nyn,nxl:nxr) ) |
---|
[1359] | 302 | z0_av_global = 0.0_wp |
---|
[1314] | 303 | |
---|
[1320] | 304 | #if defined( __parallel ) |
---|
[1314] | 305 | CALL MPI_ALLREDUCE(z0_av_local, z0_av_global, 1, MPI_REAL, MPI_SUM, & |
---|
| 306 | comm2d, ierr ) |
---|
[1320] | 307 | #else |
---|
| 308 | z0_av_global = z0_av_local |
---|
| 309 | #endif |
---|
[1314] | 310 | |
---|
| 311 | z0_av_global = z0_av_global / ( ( ny + 1 ) * ( nx + 1 ) ) |
---|
| 312 | ! |
---|
| 313 | !-- Horizontal wind speed is zero below and at z0 |
---|
[1359] | 314 | log_z_z0(0) = 0.0_wp |
---|
[1314] | 315 | ! |
---|
| 316 | !-- Calculate vertical depth of the sublayers |
---|
[1322] | 317 | height_int = ( z_p - z0_av_global ) / REAL( number_of_sublayers, KIND=wp ) |
---|
[1314] | 318 | ! |
---|
| 319 | !-- Precalculate LOG(z/z0) |
---|
[1359] | 320 | height_p = 0.0_wp |
---|
[1314] | 321 | DO k = 1, number_of_sublayers |
---|
| 322 | |
---|
| 323 | height_p = height_p + height_int |
---|
| 324 | log_z_z0(k) = LOG( height_p / z0_av_global ) |
---|
| 325 | |
---|
| 326 | ENDDO |
---|
| 327 | |
---|
| 328 | |
---|
| 329 | ENDIF |
---|
| 330 | |
---|
| 331 | ! |
---|
[1359] | 332 | !-- Check boundary condition and set internal variables |
---|
| 333 | SELECT CASE ( bc_par_b ) |
---|
| 334 | |
---|
| 335 | CASE ( 'absorb' ) |
---|
| 336 | ibc_par_b = 1 |
---|
| 337 | |
---|
| 338 | CASE ( 'reflect' ) |
---|
| 339 | ibc_par_b = 2 |
---|
| 340 | |
---|
| 341 | CASE DEFAULT |
---|
| 342 | WRITE( message_string, * ) 'unknown boundary condition ', & |
---|
| 343 | 'bc_par_b = "', TRIM( bc_par_b ), '"' |
---|
| 344 | CALL message( 'lpm_init', 'PA0217', 1, 2, 0, 6, 0 ) |
---|
| 345 | |
---|
| 346 | END SELECT |
---|
| 347 | SELECT CASE ( bc_par_t ) |
---|
| 348 | |
---|
| 349 | CASE ( 'absorb' ) |
---|
| 350 | ibc_par_t = 1 |
---|
| 351 | |
---|
| 352 | CASE ( 'reflect' ) |
---|
| 353 | ibc_par_t = 2 |
---|
| 354 | |
---|
| 355 | CASE DEFAULT |
---|
| 356 | WRITE( message_string, * ) 'unknown boundary condition ', & |
---|
| 357 | 'bc_par_t = "', TRIM( bc_par_t ), '"' |
---|
| 358 | CALL message( 'lpm_init', 'PA0218', 1, 2, 0, 6, 0 ) |
---|
| 359 | |
---|
| 360 | END SELECT |
---|
| 361 | SELECT CASE ( bc_par_lr ) |
---|
| 362 | |
---|
| 363 | CASE ( 'cyclic' ) |
---|
| 364 | ibc_par_lr = 0 |
---|
| 365 | |
---|
| 366 | CASE ( 'absorb' ) |
---|
| 367 | ibc_par_lr = 1 |
---|
| 368 | |
---|
| 369 | CASE ( 'reflect' ) |
---|
| 370 | ibc_par_lr = 2 |
---|
| 371 | |
---|
| 372 | CASE DEFAULT |
---|
| 373 | WRITE( message_string, * ) 'unknown boundary condition ', & |
---|
| 374 | 'bc_par_lr = "', TRIM( bc_par_lr ), '"' |
---|
| 375 | CALL message( 'lpm_init', 'PA0219', 1, 2, 0, 6, 0 ) |
---|
| 376 | |
---|
| 377 | END SELECT |
---|
| 378 | SELECT CASE ( bc_par_ns ) |
---|
| 379 | |
---|
| 380 | CASE ( 'cyclic' ) |
---|
| 381 | ibc_par_ns = 0 |
---|
| 382 | |
---|
| 383 | CASE ( 'absorb' ) |
---|
| 384 | ibc_par_ns = 1 |
---|
| 385 | |
---|
| 386 | CASE ( 'reflect' ) |
---|
| 387 | ibc_par_ns = 2 |
---|
| 388 | |
---|
| 389 | CASE DEFAULT |
---|
| 390 | WRITE( message_string, * ) 'unknown boundary condition ', & |
---|
| 391 | 'bc_par_ns = "', TRIM( bc_par_ns ), '"' |
---|
| 392 | CALL message( 'lpm_init', 'PA0220', 1, 2, 0, 6, 0 ) |
---|
| 393 | |
---|
| 394 | END SELECT |
---|
| 395 | |
---|
| 396 | ! |
---|
[828] | 397 | !-- Initialize collision kernels |
---|
| 398 | IF ( collision_kernel /= 'none' ) CALL init_kernels |
---|
| 399 | |
---|
| 400 | ! |
---|
[1] | 401 | !-- For the first model run of a possible job chain initialize the |
---|
[849] | 402 | !-- particles, otherwise read the particle data from restart file. |
---|
[1] | 403 | IF ( TRIM( initializing_actions ) == 'read_restart_data' & |
---|
| 404 | .AND. read_particles_from_restartfile ) THEN |
---|
| 405 | |
---|
[849] | 406 | CALL lpm_read_restart_file |
---|
[1] | 407 | |
---|
| 408 | ELSE |
---|
| 409 | |
---|
| 410 | ! |
---|
| 411 | !-- Allocate particle arrays and set attributes of the initial set of |
---|
| 412 | !-- particles, which can be also periodically released at later times. |
---|
| 413 | !-- Also allocate array for particle tail coordinates, if needed. |
---|
[1359] | 414 | ALLOCATE( prt_count(nzb:nzt+1,nysg:nyng,nxlg:nxrg), & |
---|
| 415 | grid_particles(nzb+1:nzt,nys:nyn,nxl:nxr) ) |
---|
[1] | 416 | |
---|
[1359] | 417 | maximum_number_of_particles = 0 |
---|
| 418 | number_of_particles = 0 |
---|
[792] | 419 | |
---|
| 420 | sort_count = 0 |
---|
[1359] | 421 | prt_count = 0 |
---|
[792] | 422 | |
---|
[1] | 423 | ! |
---|
| 424 | !-- Initialize all particles with dummy values (otherwise errors may |
---|
| 425 | !-- occur within restart runs). The reason for this is still not clear |
---|
| 426 | !-- and may be presumably caused by errors in the respective user-interface. |
---|
[1359] | 427 | #if defined( __twocachelines ) |
---|
| 428 | zero_particle = particle_type( 0.0_wp, 0.0_sp, 0.0_sp, 0.0_sp, 0.0_sp, & |
---|
| 429 | 0.0_sp, 0.0_sp, 0.0_wp, 0.0_wp, 0.0_wp, & |
---|
| 430 | 0, .FALSE., 0.0_wp, 0.0_wp, 0.0_wp, & |
---|
| 431 | 0.0_sp, 0.0_sp, 0.0_sp, 0.0_sp, 0.0_sp, & |
---|
| 432 | 0.0_sp, 0, 0, 0, -1) |
---|
| 433 | #else |
---|
| 434 | zero_particle = particle_type( 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, & |
---|
| 435 | 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, & |
---|
| 436 | 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, & |
---|
| 437 | 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0, 0, 0, & |
---|
| 438 | 0, .FALSE., -1) |
---|
| 439 | #endif |
---|
| 440 | particle_groups = particle_groups_type( 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp ) |
---|
[1] | 441 | |
---|
| 442 | ! |
---|
| 443 | !-- Set the default particle size used for dvrp plots |
---|
[1359] | 444 | IF ( dvrp_psize == 9999999.9_wp ) dvrp_psize = 0.2_wp * dx |
---|
[1] | 445 | |
---|
| 446 | ! |
---|
| 447 | !-- Set values for the density ratio and radius for all particle |
---|
| 448 | !-- groups, if necessary |
---|
[1359] | 449 | IF ( density_ratio(1) == 9999999.9_wp ) density_ratio(1) = 0.0_wp |
---|
| 450 | IF ( radius(1) == 9999999.9_wp ) radius(1) = 0.0_wp |
---|
[1] | 451 | DO i = 2, number_of_particle_groups |
---|
[1359] | 452 | IF ( density_ratio(i) == 9999999.9_wp ) THEN |
---|
[1] | 453 | density_ratio(i) = density_ratio(i-1) |
---|
| 454 | ENDIF |
---|
[1359] | 455 | IF ( radius(i) == 9999999.9_wp ) radius(i) = radius(i-1) |
---|
[1] | 456 | ENDDO |
---|
| 457 | |
---|
| 458 | DO i = 1, number_of_particle_groups |
---|
[1359] | 459 | IF ( density_ratio(i) /= 0.0_wp .AND. radius(i) == 0 ) THEN |
---|
[254] | 460 | WRITE( message_string, * ) 'particle group #', i, 'has a', & |
---|
| 461 | 'density ratio /= 0 but radius = 0' |
---|
[849] | 462 | CALL message( 'lpm_init', 'PA0215', 1, 2, 0, 6, 0 ) |
---|
[1] | 463 | ENDIF |
---|
| 464 | particle_groups(i)%density_ratio = density_ratio(i) |
---|
| 465 | particle_groups(i)%radius = radius(i) |
---|
| 466 | ENDDO |
---|
| 467 | |
---|
[1359] | 468 | CALL lpm_create_particle (PHASE_INIT) |
---|
[1] | 469 | ! |
---|
[1359] | 470 | !-- Set a seed value for the random number generator to be exclusively |
---|
| 471 | !-- used for the particle code. The generated random numbers should be |
---|
| 472 | !-- different on the different PEs. |
---|
| 473 | iran_part = iran_part + myid |
---|
| 474 | |
---|
| 475 | ! |
---|
| 476 | !-- User modification of initial particles |
---|
| 477 | CALL user_lpm_init |
---|
| 478 | |
---|
| 479 | ! |
---|
| 480 | !-- Open file for statistical informations about particle conditions |
---|
| 481 | IF ( write_particle_statistics ) THEN |
---|
| 482 | CALL check_open( 80 ) |
---|
| 483 | WRITE ( 80, 8000 ) current_timestep_number, simulated_time, & |
---|
| 484 | number_of_particles, & |
---|
| 485 | maximum_number_of_particles |
---|
| 486 | CALL close_file( 80 ) |
---|
| 487 | ENDIF |
---|
| 488 | |
---|
| 489 | ! |
---|
| 490 | !-- Check if particles are really uniform in color and radius (dvrp_size) |
---|
| 491 | !-- (uniform_particles is preset TRUE) |
---|
| 492 | IF ( uniform_particles ) THEN |
---|
| 493 | DO ip = nxl, nxr |
---|
| 494 | DO jp = nys, nyn |
---|
| 495 | DO kp = nzb+1, nzt |
---|
| 496 | |
---|
| 497 | n = prt_count(kp,jp,ip) |
---|
| 498 | IF ( MINVAL( grid_particles(kp,jp,ip)%particles(1:n)%dvrp_psize ) == & |
---|
| 499 | MAXVAL( grid_particles(kp,jp,ip)%particles(1:n)%dvrp_psize ) .AND. & |
---|
| 500 | MINVAL( grid_particles(kp,jp,ip)%particles(1:n)%class ) == & |
---|
| 501 | MAXVAL( grid_particles(kp,jp,ip)%particles(1:n)%class ) ) THEN |
---|
| 502 | uniform_particles_l = .TRUE. |
---|
| 503 | ELSE |
---|
| 504 | uniform_particles_l = .FALSE. |
---|
| 505 | ENDIF |
---|
| 506 | |
---|
| 507 | ENDDO |
---|
| 508 | ENDDO |
---|
| 509 | ENDDO |
---|
| 510 | |
---|
| 511 | #if defined( __parallel ) |
---|
| 512 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 513 | CALL MPI_ALLREDUCE( uniform_particles_l, uniform_particles, 1, & |
---|
| 514 | MPI_LOGICAL, MPI_LAND, comm2d, ierr ) |
---|
| 515 | #else |
---|
| 516 | uniform_particles = uniform_particles_l |
---|
| 517 | #endif |
---|
| 518 | |
---|
| 519 | ENDIF |
---|
| 520 | |
---|
| 521 | ! |
---|
| 522 | !-- Particles will probably become none-uniform, if their size and color |
---|
| 523 | !-- will be determined by flow variables |
---|
| 524 | IF ( particle_color /= 'none' .OR. particle_dvrpsize /= 'none' ) THEN |
---|
| 525 | uniform_particles = .FALSE. |
---|
| 526 | ENDIF |
---|
| 527 | |
---|
| 528 | ! !kk Not implemented aft individual particle array fort every gridcell |
---|
| 529 | ! ! |
---|
| 530 | ! !-- Set the beginning of the particle tails and their age |
---|
| 531 | ! IF ( use_particle_tails ) THEN |
---|
| 532 | ! ! |
---|
| 533 | ! !-- Choose the maximum number of tails with respect to the maximum number |
---|
| 534 | ! !-- of particles and skip_particles_for_tail |
---|
| 535 | ! maximum_number_of_tails = maximum_number_of_particles / & |
---|
| 536 | ! skip_particles_for_tail |
---|
| 537 | ! |
---|
| 538 | ! ! |
---|
| 539 | ! !-- Create a minimum number of tails in case that there is no tail |
---|
| 540 | ! !-- initially (otherwise, index errors will occur when adressing the |
---|
| 541 | ! !-- arrays below) |
---|
| 542 | ! IF ( maximum_number_of_tails == 0 ) maximum_number_of_tails = 10 |
---|
| 543 | ! |
---|
| 544 | ! ALLOCATE( particle_tail_coordinates(maximum_number_of_tailpoints,5, & |
---|
| 545 | ! maximum_number_of_tails), & |
---|
| 546 | ! new_tail_id(maximum_number_of_tails), & |
---|
| 547 | ! tail_mask(maximum_number_of_tails) ) |
---|
| 548 | ! |
---|
| 549 | ! particle_tail_coordinates = 0.0_wp |
---|
| 550 | ! minimum_tailpoint_distance = minimum_tailpoint_distance**2 |
---|
| 551 | ! number_of_initial_tails = number_of_tails |
---|
| 552 | ! |
---|
| 553 | ! nn = 0 |
---|
| 554 | ! DO n = 1, number_of_particles |
---|
| 555 | ! ! |
---|
| 556 | ! !-- Only for those particles marked above with a provisional tail_id |
---|
| 557 | ! !-- tails will be created. Particles now get their final tail_id. |
---|
| 558 | ! IF ( particles(n)%tail_id /= 0 ) THEN |
---|
| 559 | ! |
---|
| 560 | ! nn = nn + 1 |
---|
| 561 | ! particles(n)%tail_id = nn |
---|
| 562 | ! |
---|
| 563 | ! particle_tail_coordinates(1,1,nn) = particles(n)%x |
---|
| 564 | ! particle_tail_coordinates(1,2,nn) = particles(n)%y |
---|
| 565 | ! particle_tail_coordinates(1,3,nn) = particles(n)%z |
---|
| 566 | ! particle_tail_coordinates(1,4,nn) = particles(n)%class |
---|
| 567 | ! particles(n)%tailpoints = 1 |
---|
| 568 | ! IF ( minimum_tailpoint_distance /= 0.0_wp ) THEN |
---|
| 569 | ! particle_tail_coordinates(2,1,nn) = particles(n)%x |
---|
| 570 | ! particle_tail_coordinates(2,2,nn) = particles(n)%y |
---|
| 571 | ! particle_tail_coordinates(2,3,nn) = particles(n)%z |
---|
| 572 | ! particle_tail_coordinates(2,4,nn) = particles(n)%class |
---|
| 573 | ! particle_tail_coordinates(1:2,5,nn) = 0.0_wp |
---|
| 574 | ! particles(n)%tailpoints = 2 |
---|
| 575 | ! ENDIF |
---|
| 576 | ! |
---|
| 577 | ! ENDIF |
---|
| 578 | ! ENDDO |
---|
| 579 | ! ENDIF |
---|
| 580 | ! |
---|
| 581 | ! ! |
---|
| 582 | ! !-- Plot initial positions of particles (only if particle advection is |
---|
| 583 | ! !-- switched on from the beginning of the simulation (t=0)) |
---|
| 584 | ! IF ( particle_advection_start == 0.0_wp ) CALL data_output_dvrp |
---|
| 585 | |
---|
| 586 | ENDIF |
---|
| 587 | |
---|
| 588 | ! |
---|
| 589 | !-- To avoid programm abort, assign particles array to the local version of |
---|
| 590 | !-- first grid cell |
---|
| 591 | number_of_particles = prt_count(nzb+1,nys,nxl) |
---|
| 592 | particles => grid_particles(nzb+1,nys,nxl)%particles(1:number_of_particles) |
---|
| 593 | |
---|
| 594 | ! |
---|
| 595 | !-- Formats |
---|
| 596 | 8000 FORMAT (I6,1X,F7.2,4X,I10,71X,I10) |
---|
| 597 | |
---|
| 598 | END SUBROUTINE lpm_init |
---|
| 599 | |
---|
[1682] | 600 | !------------------------------------------------------------------------------! |
---|
| 601 | ! Description: |
---|
| 602 | ! ------------ |
---|
| 603 | !> @todo Missing subroutine description. |
---|
| 604 | !------------------------------------------------------------------------------! |
---|
[1359] | 605 | SUBROUTINE lpm_create_particle (phase) |
---|
| 606 | |
---|
| 607 | USE lpm_exchange_horiz_mod, & |
---|
| 608 | ONLY: lpm_exchange_horiz, lpm_move_particle, realloc_particles_array |
---|
| 609 | |
---|
| 610 | USE lpm_pack_arrays_mod, & |
---|
| 611 | ONLY: lpm_pack_all_arrays |
---|
| 612 | |
---|
| 613 | IMPLICIT NONE |
---|
| 614 | |
---|
[1682] | 615 | INTEGER(iwp) :: alloc_size !< |
---|
| 616 | INTEGER(iwp) :: i !< |
---|
| 617 | INTEGER(iwp) :: ip !< |
---|
| 618 | INTEGER(iwp) :: j !< |
---|
| 619 | INTEGER(iwp) :: jp !< |
---|
| 620 | INTEGER(iwp) :: kp !< |
---|
| 621 | INTEGER(iwp) :: loop_stride !< |
---|
| 622 | INTEGER(iwp) :: n !< |
---|
| 623 | INTEGER(iwp) :: new_size !< |
---|
| 624 | INTEGER(iwp) :: nn !< |
---|
[1359] | 625 | |
---|
[1682] | 626 | INTEGER(iwp), INTENT(IN) :: phase !< |
---|
[1359] | 627 | |
---|
[1682] | 628 | INTEGER(iwp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) :: local_count !< |
---|
| 629 | INTEGER(iwp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) :: local_start !< |
---|
[1359] | 630 | |
---|
[1682] | 631 | LOGICAL :: first_stride !< |
---|
[1359] | 632 | |
---|
[1682] | 633 | REAL(wp) :: pos_x !< |
---|
| 634 | REAL(wp) :: pos_y !< |
---|
| 635 | REAL(wp) :: pos_z !< |
---|
[1359] | 636 | |
---|
[1682] | 637 | TYPE(particle_type),TARGET :: tmp_particle !< |
---|
[1359] | 638 | |
---|
| 639 | ! |
---|
| 640 | !-- Calculate particle positions and store particle attributes, if |
---|
| 641 | !-- particle is situated on this PE |
---|
| 642 | DO loop_stride = 1, 2 |
---|
| 643 | first_stride = (loop_stride == 1) |
---|
| 644 | IF ( first_stride ) THEN |
---|
| 645 | local_count = 0 ! count number of particles |
---|
| 646 | ELSE |
---|
| 647 | local_count = prt_count ! Start address of new particles |
---|
| 648 | ENDIF |
---|
| 649 | |
---|
[1] | 650 | n = 0 |
---|
| 651 | DO i = 1, number_of_particle_groups |
---|
| 652 | |
---|
| 653 | pos_z = psb(i) |
---|
| 654 | |
---|
| 655 | DO WHILE ( pos_z <= pst(i) ) |
---|
| 656 | |
---|
| 657 | pos_y = pss(i) |
---|
| 658 | |
---|
| 659 | DO WHILE ( pos_y <= psn(i) ) |
---|
| 660 | |
---|
[1359] | 661 | IF ( pos_y >= ( nys - 0.5_wp ) * dy .AND. & |
---|
| 662 | pos_y < ( nyn + 0.5_wp ) * dy ) THEN |
---|
[1] | 663 | |
---|
| 664 | pos_x = psl(i) |
---|
| 665 | |
---|
[1575] | 666 | xloop: DO WHILE ( pos_x <= psr(i) ) |
---|
[1] | 667 | |
---|
[1359] | 668 | IF ( pos_x >= ( nxl - 0.5_wp ) * dx .AND. & |
---|
| 669 | pos_x < ( nxr + 0.5_wp ) * dx ) THEN |
---|
[1] | 670 | |
---|
| 671 | DO j = 1, particles_per_point |
---|
| 672 | |
---|
| 673 | n = n + 1 |
---|
[1359] | 674 | #if defined( __twocachelines ) |
---|
| 675 | tmp_particle%x = pos_x |
---|
| 676 | tmp_particle%y = pos_y |
---|
| 677 | tmp_particle%z = pos_z |
---|
| 678 | tmp_particle%age = 0.0_sp |
---|
| 679 | tmp_particle%age_m = 0.0_sp |
---|
| 680 | tmp_particle%dt_sum = 0.0_wp |
---|
| 681 | tmp_particle%dvrp_psize = dvrp_psize |
---|
| 682 | tmp_particle%e_m = 0.0_sp |
---|
| 683 | IF ( curvature_solution_effects ) THEN |
---|
| 684 | ! |
---|
| 685 | !-- Initial values (internal timesteps, derivative) |
---|
| 686 | !-- for Rosenbrock method |
---|
| 687 | tmp_particle%rvar1 = 1.0E-12_wp |
---|
| 688 | tmp_particle%rvar2 = 1.0E-3_wp |
---|
| 689 | tmp_particle%rvar3 = -9999999.9_wp |
---|
| 690 | ELSE |
---|
| 691 | ! |
---|
| 692 | !-- Initial values for SGS velocities |
---|
| 693 | tmp_particle%rvar1 = 0.0_wp |
---|
| 694 | tmp_particle%rvar2 = 0.0_wp |
---|
| 695 | tmp_particle%rvar3 = 0.0_wp |
---|
[1] | 696 | ENDIF |
---|
[1359] | 697 | tmp_particle%speed_x = 0.0_sp |
---|
| 698 | tmp_particle%speed_y = 0.0_sp |
---|
| 699 | tmp_particle%speed_z = 0.0_sp |
---|
| 700 | tmp_particle%origin_x = pos_x |
---|
| 701 | tmp_particle%origin_y = pos_y |
---|
| 702 | tmp_particle%origin_z = pos_z |
---|
| 703 | tmp_particle%radius = particle_groups(i)%radius |
---|
| 704 | tmp_particle%weight_factor = initial_weighting_factor |
---|
| 705 | tmp_particle%class = 1 |
---|
| 706 | tmp_particle%group = i |
---|
| 707 | tmp_particle%tailpoints = 0 |
---|
| 708 | tmp_particle%particle_mask = .TRUE. |
---|
| 709 | #else |
---|
| 710 | tmp_particle%x = pos_x |
---|
| 711 | tmp_particle%y = pos_y |
---|
| 712 | tmp_particle%z = pos_z |
---|
| 713 | tmp_particle%age = 0.0_wp |
---|
| 714 | tmp_particle%age_m = 0.0_wp |
---|
| 715 | tmp_particle%dt_sum = 0.0_wp |
---|
| 716 | tmp_particle%dvrp_psize = dvrp_psize |
---|
| 717 | tmp_particle%e_m = 0.0_wp |
---|
[824] | 718 | IF ( curvature_solution_effects ) THEN |
---|
| 719 | ! |
---|
| 720 | !-- Initial values (internal timesteps, derivative) |
---|
| 721 | !-- for Rosenbrock method |
---|
[1359] | 722 | tmp_particle%rvar1 = 1.0E-12_wp |
---|
| 723 | tmp_particle%rvar2 = 1.0E-3_wp |
---|
| 724 | tmp_particle%rvar3 = -9999999.9_wp |
---|
[824] | 725 | ELSE |
---|
| 726 | ! |
---|
| 727 | !-- Initial values for SGS velocities |
---|
[1359] | 728 | tmp_particle%rvar1 = 0.0_wp |
---|
| 729 | tmp_particle%rvar2 = 0.0_wp |
---|
| 730 | tmp_particle%rvar3 = 0.0_wp |
---|
[824] | 731 | ENDIF |
---|
[1359] | 732 | tmp_particle%speed_x = 0.0_wp |
---|
| 733 | tmp_particle%speed_y = 0.0_wp |
---|
| 734 | tmp_particle%speed_z = 0.0_wp |
---|
| 735 | tmp_particle%origin_x = pos_x |
---|
| 736 | tmp_particle%origin_y = pos_y |
---|
| 737 | tmp_particle%origin_z = pos_z |
---|
| 738 | tmp_particle%radius = particle_groups(i)%radius |
---|
| 739 | tmp_particle%weight_factor = initial_weighting_factor |
---|
| 740 | tmp_particle%class = 1 |
---|
| 741 | tmp_particle%group = i |
---|
| 742 | tmp_particle%tailpoints = 0 |
---|
| 743 | tmp_particle%particle_mask = .TRUE. |
---|
| 744 | #endif |
---|
[1] | 745 | IF ( use_particle_tails .AND. & |
---|
| 746 | MOD( n, skip_particles_for_tail ) == 0 ) THEN |
---|
| 747 | number_of_tails = number_of_tails + 1 |
---|
| 748 | ! |
---|
| 749 | !-- This is a temporary provisional setting (see |
---|
| 750 | !-- further below!) |
---|
[1359] | 751 | tmp_particle%tail_id = number_of_tails |
---|
[1] | 752 | ELSE |
---|
[1359] | 753 | tmp_particle%tail_id = 0 |
---|
[1] | 754 | ENDIF |
---|
[1575] | 755 | ! |
---|
| 756 | !-- Determine the grid indices of the particle position |
---|
[1359] | 757 | ip = ( tmp_particle%x + 0.5_wp * dx ) * ddx |
---|
| 758 | jp = ( tmp_particle%y + 0.5_wp * dy ) * ddy |
---|
| 759 | kp = tmp_particle%z / dz + 1 + offset_ocean_nzt_m1 |
---|
[1] | 760 | |
---|
[1575] | 761 | IF ( seed_follows_topography ) THEN |
---|
| 762 | ! |
---|
| 763 | !-- Particle height is given relative to topography |
---|
| 764 | kp = kp + nzb_w_inner(jp,ip) |
---|
| 765 | tmp_particle%z = tmp_particle%z + zw(kp) |
---|
| 766 | IF ( kp > nzt ) THEN |
---|
| 767 | pos_x = pos_x + pdx(i) |
---|
| 768 | CYCLE xloop |
---|
| 769 | ENDIF |
---|
| 770 | ENDIF |
---|
| 771 | |
---|
[1359] | 772 | local_count(kp,jp,ip) = local_count(kp,jp,ip) + 1 |
---|
| 773 | IF ( .NOT. first_stride ) THEN |
---|
| 774 | IF ( ip < nxl .OR. jp < nys .OR. kp < nzb+1 ) THEN |
---|
| 775 | write(6,*) 'xl ',ip,jp,kp,nxl,nys,nzb+1 |
---|
| 776 | ENDIF |
---|
| 777 | IF ( ip > nxr .OR. jp > nyn .OR. kp > nzt ) THEN |
---|
| 778 | write(6,*) 'xu ',ip,jp,kp,nxr,nyn,nzt |
---|
| 779 | ENDIF |
---|
| 780 | grid_particles(kp,jp,ip)%particles(local_count(kp,jp,ip)) = tmp_particle |
---|
| 781 | ENDIF |
---|
[1] | 782 | ENDDO |
---|
| 783 | |
---|
| 784 | ENDIF |
---|
| 785 | |
---|
| 786 | pos_x = pos_x + pdx(i) |
---|
| 787 | |
---|
[1575] | 788 | ENDDO xloop |
---|
[1] | 789 | |
---|
| 790 | ENDIF |
---|
| 791 | |
---|
| 792 | pos_y = pos_y + pdy(i) |
---|
| 793 | |
---|
| 794 | ENDDO |
---|
| 795 | |
---|
| 796 | pos_z = pos_z + pdz(i) |
---|
| 797 | |
---|
| 798 | ENDDO |
---|
| 799 | |
---|
| 800 | ENDDO |
---|
| 801 | |
---|
[1359] | 802 | IF ( first_stride ) THEN |
---|
| 803 | DO ip = nxl, nxr |
---|
| 804 | DO jp = nys, nyn |
---|
| 805 | DO kp = nzb+1, nzt |
---|
| 806 | IF ( phase == PHASE_INIT ) THEN |
---|
| 807 | IF ( local_count(kp,jp,ip) > 0 ) THEN |
---|
| 808 | alloc_size = MAX( INT( local_count(kp,jp,ip) * & |
---|
| 809 | ( 1.0_wp + alloc_factor / 100.0_wp ) ), & |
---|
| 810 | min_nr_particle ) |
---|
| 811 | ELSE |
---|
| 812 | alloc_size = min_nr_particle |
---|
| 813 | ENDIF |
---|
| 814 | ALLOCATE(grid_particles(kp,jp,ip)%particles(1:alloc_size)) |
---|
| 815 | DO n = 1, alloc_size |
---|
| 816 | grid_particles(kp,jp,ip)%particles(n) = zero_particle |
---|
| 817 | ENDDO |
---|
| 818 | ELSEIF ( phase == PHASE_RELEASE ) THEN |
---|
| 819 | IF ( local_count(kp,jp,ip) > 0 ) THEN |
---|
| 820 | new_size = local_count(kp,jp,ip) + prt_count(kp,jp,ip) |
---|
| 821 | alloc_size = MAX( INT( new_size * ( 1.0_wp + & |
---|
| 822 | alloc_factor / 100.0_wp ) ), min_nr_particle ) |
---|
| 823 | IF( alloc_size > SIZE( grid_particles(kp,jp,ip)%particles) ) THEN |
---|
| 824 | CALL realloc_particles_array(ip,jp,kp,alloc_size) |
---|
| 825 | ENDIF |
---|
| 826 | ENDIF |
---|
| 827 | ENDIF |
---|
| 828 | ENDDO |
---|
| 829 | ENDDO |
---|
| 830 | ENDDO |
---|
| 831 | ENDIF |
---|
| 832 | ENDDO |
---|
[1] | 833 | |
---|
[1359] | 834 | local_start = prt_count+1 |
---|
| 835 | prt_count = local_count |
---|
[1] | 836 | ! |
---|
[1359] | 837 | !-- Add random fluctuation to particle positions |
---|
| 838 | IF ( random_start_position ) THEN |
---|
| 839 | DO ip = nxl, nxr |
---|
| 840 | DO jp = nys, nyn |
---|
| 841 | DO kp = nzb+1, nzt |
---|
| 842 | number_of_particles = prt_count(kp,jp,ip) |
---|
| 843 | IF ( number_of_particles <= 0 ) CYCLE |
---|
| 844 | particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles) |
---|
[1] | 845 | |
---|
[1359] | 846 | DO n = local_start(kp,jp,ip), number_of_particles !Move only new particles |
---|
| 847 | IF ( psl(particles(n)%group) /= psr(particles(n)%group) ) THEN |
---|
| 848 | particles(n)%x = particles(n)%x + & |
---|
| 849 | ( random_function( iran_part ) - 0.5_wp ) * & |
---|
| 850 | pdx(particles(n)%group) |
---|
| 851 | ENDIF |
---|
| 852 | IF ( pss(particles(n)%group) /= psn(particles(n)%group) ) THEN |
---|
| 853 | particles(n)%y = particles(n)%y + & |
---|
| 854 | ( random_function( iran_part ) - 0.5_wp ) * & |
---|
| 855 | pdy(particles(n)%group) |
---|
| 856 | ENDIF |
---|
| 857 | IF ( psb(particles(n)%group) /= pst(particles(n)%group) ) THEN |
---|
| 858 | particles(n)%z = particles(n)%z + & |
---|
| 859 | ( random_function( iran_part ) - 0.5_wp ) * & |
---|
| 860 | pdz(particles(n)%group) |
---|
| 861 | ENDIF |
---|
| 862 | ENDDO |
---|
[1] | 863 | ! |
---|
[1359] | 864 | !-- Identify particles located outside the model domain |
---|
| 865 | CALL lpm_boundary_conds( 'bottom/top' ) |
---|
| 866 | ENDDO |
---|
| 867 | ENDDO |
---|
| 868 | ENDDO |
---|
[1] | 869 | ! |
---|
[1359] | 870 | !-- Exchange particles between grid cells and processors |
---|
| 871 | CALL lpm_move_particle |
---|
| 872 | CALL lpm_exchange_horiz |
---|
[1] | 873 | |
---|
[1359] | 874 | ENDIF |
---|
[1] | 875 | ! |
---|
[1359] | 876 | !-- In case of random_start_position, delete particles identified by |
---|
| 877 | !-- lpm_exchange_horiz and lpm_boundary_conds. Then sort particles into blocks, |
---|
| 878 | !-- which is needed for a fast interpolation of the LES fields on the particle |
---|
| 879 | !-- position. |
---|
| 880 | CALL lpm_pack_all_arrays |
---|
[1] | 881 | |
---|
| 882 | ! |
---|
[1359] | 883 | !-- Determine maximum number of particles (i.e., all possible particles that |
---|
| 884 | !-- have been allocated) and the current number of particles |
---|
| 885 | DO ip = nxl, nxr |
---|
| 886 | DO jp = nys, nyn |
---|
| 887 | DO kp = nzb+1, nzt |
---|
| 888 | maximum_number_of_particles = maximum_number_of_particles & |
---|
| 889 | + SIZE(grid_particles(kp,jp,ip)%particles) |
---|
| 890 | number_of_particles = number_of_particles & |
---|
| 891 | + prt_count(kp,jp,ip) |
---|
[1] | 892 | ENDDO |
---|
[1359] | 893 | ENDDO |
---|
| 894 | ENDDO |
---|
[1] | 895 | ! |
---|
[1359] | 896 | !-- Calculate the number of particles and tails of the total domain |
---|
[1] | 897 | #if defined( __parallel ) |
---|
[1359] | 898 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 899 | CALL MPI_ALLREDUCE( number_of_particles, total_number_of_particles, 1, & |
---|
| 900 | MPI_INTEGER, MPI_SUM, comm2d, ierr ) |
---|
| 901 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 902 | CALL MPI_ALLREDUCE( number_of_tails, total_number_of_tails, 1, & |
---|
| 903 | MPI_INTEGER, MPI_SUM, comm2d, ierr ) |
---|
[1] | 904 | #else |
---|
[1359] | 905 | total_number_of_particles = number_of_particles |
---|
| 906 | total_number_of_tails = number_of_tails |
---|
[1] | 907 | #endif |
---|
| 908 | |
---|
[1359] | 909 | RETURN |
---|
[1] | 910 | |
---|
[1359] | 911 | END SUBROUTINE lpm_create_particle |
---|
[336] | 912 | |
---|
[1359] | 913 | END MODULE lpm_init_mod |
---|