70 const InitFunctor& create_functor, ParticleListType& particle_list,
71 const int particles_per_cell, LocalGridType& local_grid,
72 const std::size_t previous_num_particles = 0,
73 const bool shrink_to_fit =
true,
const uint64_t seed = 123456,
78 using memory_space =
typename ParticleListType::memory_space;
84 const auto& global_grid = local_grid.globalGrid();
87 auto owned_cells = local_grid.indexSpace(
Own(),
Cell(),
Local() );
90 const auto local_seed = seed + global_grid.blockId();
91 using rnd_type = Kokkos::Random_XorShift64_Pool<ExecutionSpace>;
93#if ( KOKKOS_VERSION < 40700 )
95 pool.init( local_seed, owned_cells.size() );
97 rnd_type pool( local_seed, owned_cells.size() );
101 auto& aosoa = particle_list.aosoa();
105 int num_particles = particles_per_cell * owned_cells.size();
106 aosoa.resize( previous_num_particles + num_particles );
109 auto count = Kokkos::View<int*, memory_space>(
"particle_count", 1 );
110 Kokkos::deep_copy( count, previous_num_particles );
114 "Cabana::Grid::ParticleInit::Random", exec_space, owned_cells,
115 KOKKOS_LAMBDA(
const int i,
const int j,
const int k ) {
117 int i_own = i - owned_cells.min( Dim::I );
118 int j_own = j - owned_cells.min( Dim::J );
119 int k_own = k - owned_cells.min( Dim::K );
121 i_own + owned_cells.extent( Dim::I ) *
122 ( j_own + k_own * owned_cells.extent( Dim::J ) );
125 int low_node[3] = { i, j, k };
126 double low_coords[3];
127 local_mesh.coordinates(
Node(), low_node, low_coords );
130 int high_node[3] = { i + 1, j + 1, k + 1 };
131 double high_coords[3];
132 local_mesh.coordinates(
Node(), high_node, high_coords );
135 auto rand = pool.get_state( cell_id );
142 local_mesh.measure(
Cell(), low_node ) / particles_per_cell;
145 for (
int p = 0; p < particles_per_cell; ++p )
149 previous_num_particles + cell_id * particles_per_cell + p;
153 for (
int d = 0; d < 3; ++d )
155 px[d] = Kokkos::rand<
decltype( rand ),
double>::draw(
156 rand, low_coords[d], high_coords[d] );
160 auto particle = particle_list.getParticle( pid );
161 bool created = create_functor( pid, px, pv, particle );
166 auto c = Kokkos::atomic_fetch_add( &count( 0 ), 1 );
167 particle_list.setParticle( particle, c );
174 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), count );
175 aosoa.resize( count_host( 0 ) );
178 return count_host( 0 );
203 ParticleListType& particle_list,
const int particles_per_cell,
204 LocalGridType& local_grid,
const std::size_t previous_num_particles = 0,
205 const bool shrink_to_fit =
true,
const uint64_t seed = 123456,
209 using exec_space =
typename ParticleListType::memory_space::execution_space;
210 return createParticles( tag, exec_space{}, create_functor, particle_list,
211 particles_per_cell, local_grid,
212 previous_num_particles, shrink_to_fit, seed );
231 PositionType& positions,
const int particles_per_cell,
232 LocalGridType& local_grid,
const std::size_t previous_num_particles = 0,
233 const uint64_t seed = 123456,
235 Kokkos::is_view<PositionType>::value ),
239 using memory_space =
typename PositionType::memory_space;
245 const auto& global_grid = local_grid.globalGrid();
248 auto owned_cells = local_grid.indexSpace(
Own(),
Cell(),
Local() );
251 const auto local_seed = seed + global_grid.blockId();
252 using rnd_type = Kokkos::Random_XorShift64_Pool<ExecutionSpace>;
254#if ( KOKKOS_VERSION < 40700 )
256 pool.init( local_seed, owned_cells.size() );
258 rnd_type pool( local_seed, owned_cells.size() );
262 assert( positions.size() ==
static_cast<std::size_t
>( particles_per_cell *
263 owned_cells.size() ) +
264 previous_num_particles );
268 "Cabana::Grid::ParticleInit::Random", exec_space, owned_cells,
269 KOKKOS_LAMBDA(
const int i,
const int j,
const int k ) {
271 int i_own = i - owned_cells.min( Dim::I );
272 int j_own = j - owned_cells.min( Dim::J );
273 int k_own = k - owned_cells.min( Dim::K );
275 i_own + owned_cells.extent( Dim::I ) *
276 ( j_own + k_own * owned_cells.extent( Dim::J ) );
279 int low_node[3] = { i, j, k };
280 double low_coords[3];
281 local_mesh.coordinates(
Node(), low_node, low_coords );
284 int high_node[3] = { i + 1, j + 1, k + 1 };
285 double high_coords[3];
286 local_mesh.coordinates(
Node(), high_node, high_coords );
289 auto rand = pool.get_state( cell_id );
292 for (
int p = 0; p < particles_per_cell; ++p )
296 previous_num_particles + cell_id * particles_per_cell + p;
300 for (
int d = 0; d < 3; ++d )
302 positions( pid, d ) =
303 Kokkos::rand<
decltype( rand ),
double>::draw(
304 rand, low_coords[d], high_coords[d] );
325 const int particles_per_cell, LocalGridType& local_grid,
326 const std::size_t previous_num_particles = 0,
const uint64_t seed = 123456,
328 Kokkos::is_view<PositionType>::value ),
331 using exec_space =
typename PositionType::execution_space;
333 local_grid, previous_num_particles, seed );
362 const InitFunctor& create_functor, ParticleListType& particle_list,
363 const int particles_per_cell_dim, LocalGridType& local_grid,
364 const std::size_t previous_num_particles = 0,
365 const bool shrink_to_fit =
true,
370 using memory_space =
typename ParticleListType::memory_space;
376 auto owned_cells = local_grid.indexSpace(
Own(),
Cell(),
Local() );
379 auto& aosoa = particle_list.aosoa();
382 int particles_per_cell = particles_per_cell_dim * particles_per_cell_dim *
383 particles_per_cell_dim;
384 int num_particles = particles_per_cell * owned_cells.size();
385 aosoa.resize( previous_num_particles + num_particles );
388 auto count = Kokkos::View<int*, memory_space>(
"particle_count", 1 );
389 Kokkos::deep_copy( count, previous_num_particles );
393 "Cabana::Grid::ParticleInit::Uniform", exec_space, owned_cells,
394 KOKKOS_LAMBDA(
const int i,
const int j,
const int k ) {
396 int i_own = i - owned_cells.min( Dim::I );
397 int j_own = j - owned_cells.min( Dim::J );
398 int k_own = k - owned_cells.min( Dim::K );
400 i_own + owned_cells.extent( Dim::I ) *
401 ( j_own + k_own * owned_cells.extent( Dim::J ) );
404 int low_node[3] = { i, j, k };
405 double low_coords[3];
406 local_mesh.coordinates(
Node(), low_node, low_coords );
409 int high_node[3] = { i + 1, j + 1, k + 1 };
410 double high_coords[3];
411 local_mesh.coordinates(
Node(), high_node, high_coords );
414 double spacing[3] = { ( high_coords[Dim::I] - low_coords[Dim::I] ) /
415 particles_per_cell_dim,
416 ( high_coords[Dim::J] - low_coords[Dim::J] ) /
417 particles_per_cell_dim,
418 ( high_coords[Dim::K] - low_coords[Dim::K] ) /
419 particles_per_cell_dim };
426 local_mesh.measure(
Cell(), low_node ) / particles_per_cell;
429 for (
int ip = 0; ip < particles_per_cell_dim; ++ip )
430 for (
int jp = 0; jp < particles_per_cell_dim; ++jp )
431 for (
int kp = 0; kp < particles_per_cell_dim; ++kp )
434 int pid = previous_num_particles +
435 cell_id * particles_per_cell + ip +
436 particles_per_cell_dim *
437 ( jp + particles_per_cell_dim * kp );
440 px[Dim::I] = 0.5 * spacing[Dim::I] +
441 ip * spacing[Dim::I] + low_coords[Dim::I];
442 px[Dim::J] = 0.5 * spacing[Dim::J] +
443 jp * spacing[Dim::J] + low_coords[Dim::J];
444 px[Dim::K] = 0.5 * spacing[Dim::K] +
445 kp * spacing[Dim::K] + low_coords[Dim::K];
449 auto particle = particle_list.getParticle( pid );
450 bool created = create_functor( pid, px, pv, particle );
457 auto c = Kokkos::atomic_fetch_add( &count( 0 ), 1 );
458 particle_list.setParticle( particle, c );
465 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), count );
466 aosoa.resize( count_host( 0 ) );
469 return count_host( 0 );
494 ParticleListType& particle_list,
const int particles_per_cell_dim,
495 LocalGridType& local_grid,
const std::size_t previous_num_particles = 0,
496 const bool shrink_to_fit =
true,
500 using exec_space =
typename ParticleListType::memory_space::execution_space;
501 return createParticles( tag, exec_space{}, create_functor, particle_list,
502 particles_per_cell_dim, local_grid,
503 previous_num_particles, shrink_to_fit );
522 PositionType& positions,
const int particles_per_cell_dim,
523 LocalGridType& local_grid,
const std::size_t previous_num_particles = 0,
525 Kokkos::is_view<PositionType>::value ),
529 using memory_space =
typename PositionType::memory_space;
535 auto owned_cells = local_grid.indexSpace(
Own(),
Cell(),
Local() );
537 int particles_per_cell = particles_per_cell_dim * particles_per_cell_dim *
538 particles_per_cell_dim;
541 assert( positions.size() ==
static_cast<std::size_t
>( particles_per_cell *
542 owned_cells.size() ) +
543 previous_num_particles );
547 "Cabana::Grid::ParticleInit::Uniform", exec_space, owned_cells,
548 KOKKOS_LAMBDA(
const int i,
const int j,
const int k ) {
550 int i_own = i - owned_cells.min( Dim::I );
551 int j_own = j - owned_cells.min( Dim::J );
552 int k_own = k - owned_cells.min( Dim::K );
554 i_own + owned_cells.extent( Dim::I ) *
555 ( j_own + k_own * owned_cells.extent( Dim::J ) );
558 int low_node[3] = { i, j, k };
559 double low_coords[3];
560 local_mesh.coordinates(
Node(), low_node, low_coords );
563 int high_node[3] = { i + 1, j + 1, k + 1 };
564 double high_coords[3];
565 local_mesh.coordinates(
Node(), high_node, high_coords );
568 double spacing[3] = { ( high_coords[Dim::I] - low_coords[Dim::I] ) /
569 particles_per_cell_dim,
570 ( high_coords[Dim::J] - low_coords[Dim::J] ) /
571 particles_per_cell_dim,
572 ( high_coords[Dim::K] - low_coords[Dim::K] ) /
573 particles_per_cell_dim };
576 for (
int ip = 0; ip < particles_per_cell_dim; ++ip )
577 for (
int jp = 0; jp < particles_per_cell_dim; ++jp )
578 for (
int kp = 0; kp < particles_per_cell_dim; ++kp )
581 int pid = previous_num_particles +
582 cell_id * particles_per_cell + ip +
583 particles_per_cell_dim *
584 ( jp + particles_per_cell_dim * kp );
587 positions( pid, 0 ) = 0.5 * spacing[Dim::I] +
588 ip * spacing[Dim::I] +
590 positions( pid, 1 ) = 0.5 * spacing[Dim::J] +
591 jp * spacing[Dim::J] +
593 positions( pid, 2 ) = 0.5 * spacing[Dim::K] +
594 kp * spacing[Dim::K] +
616 const int particles_per_cell_dim, LocalGridType& local_grid,
617 const std::size_t previous_num_particles = 0,
619 Kokkos::is_view<PositionType>::value ),
622 using exec_space =
typename PositionType::execution_space;
624 local_grid, previous_num_particles );
void grid_parallel_for(const std::string &label, const ExecutionSpace &exec_space, const IndexSpace< N > &index_space, const FunctorType &functor)
Execute a functor in parallel with a multidimensional execution policy specified by the given index s...
Definition Cabana_Grid_Parallel.hpp:52
int createParticles(Cabana::InitRandom, const ExecutionSpace &exec_space, const InitFunctor &create_functor, ParticleListType &particle_list, const int particles_per_cell, LocalGridType &local_grid, const std::size_t previous_num_particles=0, const bool shrink_to_fit=true, const uint64_t seed=123456, typename std::enable_if< is_particle_list< ParticleListType >::value, int >::type *=0)
Initialize a random number of particles in each cell given an initialization functor.
Definition Cabana_Grid_ParticleInit.hpp:68
Application-level particle/mesh storage and single particle access.
Particle creation utilities.
ParticleList static type checker.
Definition Cabana_Grid_ParticleList.hpp:119