WarpX
BinaryCollision.H
Go to the documentation of this file.
1 /* Copyright 2020-2021 Yinjian Zhao, David Grote, Neil Zaim
2  *
3  * This file is part of WarpX.
4  *
5  * License: BSD-3-Clause-LBNL
6  */
7 #ifndef WARPX_PARTICLES_COLLISION_BINARYCOLLISION_H_
8 #define WARPX_PARTICLES_COLLISION_BINARYCOLLISION_H_
9 
20 #include "Utils/ParticleUtils.H"
21 #include "Utils/TextMsg.H"
23 #include "WarpX.H"
24 
27 
28 #include <AMReX.H>
29 #include <AMReX_Algorithm.H>
30 #include <AMReX_BLassert.H>
31 #include <AMReX_Config.H>
32 #include <AMReX_DenseBins.H>
33 #include <AMReX_Extension.H>
34 #include <AMReX_Geometry.H>
35 #include <AMReX_GpuAtomic.H>
36 #include <AMReX_GpuContainers.H>
37 #include <AMReX_GpuControl.H>
38 #include <AMReX_GpuDevice.H>
39 #include <AMReX_GpuLaunch.H>
40 #include <AMReX_GpuQualifiers.H>
41 #include <AMReX_LayoutData.H>
42 #include <AMReX_MFIter.H>
43 #include <AMReX_PODVector.H>
44 #include <AMReX_ParmParse.H>
45 #include <AMReX_Particles.H>
46 #include <AMReX_ParticleTile.H>
47 #include <AMReX_Random.H>
48 #include <AMReX_REAL.H>
49 #include <AMReX_Scan.H>
50 #include <AMReX_Utility.H>
51 #include <AMReX_Vector.H>
52 
53 #include <AMReX_BaseFwd.H>
54 
55 #include <cmath>
56 #include <string>
57 
67 template <typename CollisionFunctorType,
68  typename CopyTransformFunctorType = NoParticleCreationFunc>
69 class BinaryCollision final
70  : public CollisionBase
71 {
72  // Define shortcuts for frequently-used type names
78 
79 public:
87  BinaryCollision (std::string collision_name, MultiParticleContainer const * const mypc)
88  : CollisionBase(collision_name)
89  {
90  if(m_species_names.size() != 2)
91  WARPX_ABORT_WITH_MESSAGE("Binary collision " + collision_name + " must have exactly two species.");
92 
94 
95  m_binary_collision_functor = CollisionFunctorType(collision_name, mypc, m_isSameSpecies);
96 
97  const amrex::ParmParse pp_collision_name(collision_name);
98  pp_collision_name.queryarr("product_species", m_product_species);
101  WARPX_ABORT_WITH_MESSAGE( "Binary collision " + collision_name +
102  " does not produce species. Thus, `product_species` should not be specified in the input script." );
103  }
104  m_copy_transform_functor = CopyTransformFunctorType(collision_name, mypc);
105  }
106 
107  ~BinaryCollision () override = default;
108 
109  BinaryCollision ( BinaryCollision const &) = default;
113 
121  void doCollisions (amrex::Real cur_time, amrex::Real dt, MultiParticleContainer* mypc) override
122  {
123  amrex::ignore_unused(cur_time);
124 
126  auto& species2 = mypc->GetParticleContainerFromName(m_species_names[1]);
127 
128  // In case of particle creation, create the necessary vectors
129  const int n_product_species = m_product_species.size();
130  amrex::Vector<WarpXParticleContainer*> product_species_vector;
131  amrex::Vector<SmartCopyFactory> copy_factory_species1;
132  amrex::Vector<SmartCopyFactory> copy_factory_species2;
133  amrex::Vector<SmartCopy> copy_species1;
134  amrex::Vector<SmartCopy> copy_species2;
135  for (int i = 0; i < n_product_species; i++)
136  {
137  auto& product = mypc->GetParticleContainerFromName(m_product_species[i]);
138  product.defineAllParticleTiles();
139  product_species_vector.push_back(&product);
140  // Although the copy factories are not explicitly reused past this point, we need to
141  // store them in vectors so that the data that they own, which is used by the smart
142  // copy functors, does not go out of scope at the end of this for loop.
143  copy_factory_species1.push_back(SmartCopyFactory(species1, product));
144  copy_factory_species2.push_back(SmartCopyFactory(species2, product));
145  copy_species1.push_back(copy_factory_species1[i].getSmartCopy());
146  copy_species2.push_back(copy_factory_species2[i].getSmartCopy());
147  }
148 #ifdef AMREX_USE_GPU
149  amrex::Gpu::DeviceVector<SmartCopy> device_copy_species1(n_product_species);
150  amrex::Gpu::DeviceVector<SmartCopy> device_copy_species2(n_product_species);
151  amrex::Gpu::copyAsync(amrex::Gpu::hostToDevice, copy_species1.begin(),
152  copy_species1.end(), device_copy_species1.begin());
153  amrex::Gpu::copyAsync(amrex::Gpu::hostToDevice, copy_species2.begin(),
154  copy_species2.end(), device_copy_species2.begin());
156  auto copy_species1_data = device_copy_species1.data();
157  auto copy_species2_data = device_copy_species2.data();
158 #else
159  auto copy_species1_data = copy_species1.data();
160  auto copy_species2_data = copy_species2.data();
161 #endif
163  species1.defineAllParticleTiles();
164  species2.defineAllParticleTiles();
165  }
166 
167  // Enable tiling
168  amrex::MFItInfo info;
169  if (amrex::Gpu::notInLaunchRegion()) info.EnableTiling(species1.tile_size);
170 
171  // Loop over refinement levels
172  for (int lev = 0; lev <= species1.finestLevel(); ++lev){
173 
175 
176  // Loop over all grids/tiles at this level
177 #ifdef AMREX_USE_OMP
178  info.SetDynamic(true);
179 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
180 #endif
181  for (amrex::MFIter mfi = species1.MakeMFIter(lev, info); mfi.isValid(); ++mfi){
183  {
185  }
186  auto wt = static_cast<amrex::Real>(amrex::second());
187 
188  doCollisionsWithinTile( dt, lev, mfi, species1, species2, product_species_vector,
189  copy_species1_data, copy_species2_data);
190 
192  {
194  wt = static_cast<amrex::Real>(amrex::second()) - wt;
195  amrex::HostDevice::Atomic::Add( &(*cost)[mfi.index()], wt);
196  }
197  }
198  }
199  }
200 
214  amrex::Real dt, int const lev, amrex::MFIter const& mfi,
215  WarpXParticleContainer& species_1,
216  WarpXParticleContainer& species_2,
217  amrex::Vector<WarpXParticleContainer*> product_species_vector,
218  SmartCopy* copy_species1, SmartCopy* copy_species2)
219  {
220  using namespace ParticleUtils;
221  using namespace amrex::literals;
222 
223  CollisionFunctorType binary_collision_functor = m_binary_collision_functor;
224  const bool have_product_species = m_have_product_species;
225 
226  // Store product species data in vectors
227  const int n_product_species = m_product_species.size();
228  amrex::Vector<ParticleTileType*> tile_products;
229  amrex::Vector<GetParticlePosition<PIdx>> get_position_products;
230  amrex::Vector<index_type> products_np;
232  constexpr int getpos_offset = 0;
233  for (int i = 0; i < n_product_species; i++)
234  {
235  ParticleTileType& ptile_product = product_species_vector[i]->ParticlesAt(lev, mfi);
236  tile_products.push_back(&ptile_product);
237  get_position_products.push_back(GetParticlePosition<PIdx>(ptile_product,
238  getpos_offset));
239  products_np.push_back(ptile_product.numParticles());
240  products_mass.push_back(product_species_vector[i]->getMass());
241  }
242  auto tile_products_data = tile_products.data();
243 
244  if ( m_isSameSpecies ) // species_1 == species_2
245  {
246  // Extract particles in the tile that `mfi` points to
247  ParticleTileType& ptile_1 = species_1.ParticlesAt(lev, mfi);
248 
249  // Find the particles that are in each cell of this tile
250  ParticleBins bins_1 = findParticlesInEachCell( lev, mfi, ptile_1 );
251 
252  // Loop over cells, and collide the particles in each cell
253 
254  // Extract low-level data
255  auto const n_cells = static_cast<int>(bins_1.numBins());
256  // - Species 1
257  const auto soa_1 = ptile_1.getParticleTileData();
258  index_type* AMREX_RESTRICT indices_1 = bins_1.permutationPtr();
259  index_type const* AMREX_RESTRICT cell_offsets_1 = bins_1.offsetsPtr();
260  const amrex::ParticleReal q1 = species_1.getCharge();
261  const amrex::ParticleReal m1 = species_1.getMass();
262  auto get_position_1 = GetParticlePosition<PIdx>(ptile_1, getpos_offset);
263  // Needed to access the particle id
265  particle_ptr_1 = ptile_1.GetArrayOfStructs()().data();
266 
267  amrex::Geometry const& geom = WarpX::GetInstance().Geom(lev);
268 #if defined WARPX_DIM_1D_Z
269  auto dV = geom.CellSize(0);
270 #elif defined WARPX_DIM_XZ
271  auto dV = geom.CellSize(0) * geom.CellSize(1);
272 #elif defined WARPX_DIM_RZ
273  amrex::Box const& cbx = mfi.tilebox(amrex::IntVect::TheZeroVector()); //Cell-centered box
274  const auto lo = lbound(cbx);
275  const auto hi = ubound(cbx);
276  int const nz = hi.y-lo.y+1;
277  auto dr = geom.CellSize(0);
278  auto dz = geom.CellSize(1);
279 #elif defined(WARPX_DIM_3D)
280  auto dV = geom.CellSize(0) * geom.CellSize(1) * geom.CellSize(2);
281 #endif
282 
283 
284  /*
285  The following calculations are only required when creating product particles
286  */
287  const int n_cells_products = have_product_species ? n_cells: 0;
288  amrex::Gpu::DeviceVector<index_type> n_pairs_in_each_cell(n_cells_products);
289  index_type* AMREX_RESTRICT p_n_pairs_in_each_cell = n_pairs_in_each_cell.dataPtr();
290 
291  // Compute how many pairs in each cell and store in n_pairs_in_each_cell array
292  // For a single species, the number of pair in a cell is half the number of particles
293  // in that cell, rounded up to the next higher integer.
294  amrex::ParallelFor( n_cells_products,
295  [=] AMREX_GPU_DEVICE (int i_cell) noexcept
296  {
297  const auto n_part_in_cell = cell_offsets_1[i_cell+1] - cell_offsets_1[i_cell];
298  // Particular case: if there's only 1 particle in a cell, then there's no pair
299  p_n_pairs_in_each_cell[i_cell] = (n_part_in_cell == 1)? 0: (n_part_in_cell+1)/2;
300  }
301  );
302 
303  // Start indices of the pairs in a cell. Will be used for particle creation.
304  amrex::Gpu::DeviceVector<index_type> pair_offsets(n_cells_products);
305  const index_type n_total_pairs = (n_cells_products == 0) ? 0:
306  amrex::Scan::ExclusiveSum(n_cells_products,
307  p_n_pairs_in_each_cell, pair_offsets.data());
308  index_type* AMREX_RESTRICT p_pair_offsets = pair_offsets.dataPtr();
309 
310  // mask: equal to 1 if particle creation occurs for a given pair, 0 otherwise
313  // Will be filled with the index of the first particle of a given pair
314  amrex::Gpu::DeviceVector<index_type> pair_indices_1(n_total_pairs);
315  index_type* AMREX_RESTRICT p_pair_indices_1 = pair_indices_1.dataPtr();
316  // Will be filled with the index of the second particle of a given pair
317  amrex::Gpu::DeviceVector<index_type> pair_indices_2(n_total_pairs);
318  index_type* AMREX_RESTRICT p_pair_indices_2 = pair_indices_2.dataPtr();
319  // How much weight should be given to the produced particles (and removed from the
320  // reacting particles)
321  amrex::Gpu::DeviceVector<amrex::ParticleReal> pair_reaction_weight(n_total_pairs);
322  amrex::ParticleReal* AMREX_RESTRICT p_pair_reaction_weight =
323  pair_reaction_weight.dataPtr();
324  /*
325  End of calculations only required when creating product particles
326  */
327 
328 
329  // Loop over cells
330  amrex::ParallelForRNG( n_cells,
331  [=] AMREX_GPU_DEVICE (int i_cell, amrex::RandomEngine const& engine) noexcept
332  {
333  // The particles from species1 that are in the cell `i_cell` are
334  // given by the `indices_1[cell_start_1:cell_stop_1]`
335  index_type const cell_start_1 = cell_offsets_1[i_cell];
336  index_type const cell_stop_1 = cell_offsets_1[i_cell+1];
337  index_type const cell_half_1 = (cell_start_1+cell_stop_1)/2;
338 
339  // Same but for the pairs
340  index_type const cell_start_pair = have_product_species?
341  p_pair_offsets[i_cell] : 0;
342 
343  // Do not collide if there is only one particle in the cell
344  if ( cell_stop_1 - cell_start_1 <= 1 ) return;
345 
346  // shuffle
348  indices_1, cell_start_1, cell_half_1, engine );
349 #if defined WARPX_DIM_RZ
350  const int ri = (i_cell - i_cell%nz) / nz;
351  auto dV = MathConst::pi*(2.0_prt*ri+1.0_prt)*dr*dr*dz;
352 #endif
353  // Call the function in order to perform collisions
354  // If there are product species, mask, p_pair_indices_1/2, and
355  // p_pair_reaction_weight are filled here
356  binary_collision_functor(
357  cell_start_1, cell_half_1,
358  cell_half_1, cell_stop_1,
359  indices_1, indices_1,
360  soa_1, soa_1, get_position_1, get_position_1,
361  q1, q1, m1, m1, dt, dV,
362  cell_start_pair, p_mask, p_pair_indices_1, p_pair_indices_2,
363  p_pair_reaction_weight, engine );
364  }
365  );
366 
367  // Create the new product particles and define their initial values
368  // num_added: how many particles of each product species have been created
369  const amrex::Vector<int> num_added = m_copy_transform_functor(n_total_pairs,
370  soa_1, soa_1, tile_products_data,
371  particle_ptr_1, particle_ptr_1, m1, m1,
372  products_mass, p_mask, products_np,
373  copy_species1, copy_species2,
374  p_pair_indices_1, p_pair_indices_2,
375  p_pair_reaction_weight);
376 
377  for (int i = 0; i < n_product_species; i++)
378  {
379  setNewParticleIDs(*(tile_products_data[i]), static_cast<int>(products_np[i]), num_added[i]);
380  }
381  }
382  else // species_1 != species_2
383  {
384  // Extract particles in the tile that `mfi` points to
385  ParticleTileType& ptile_1 = species_1.ParticlesAt(lev, mfi);
386  ParticleTileType& ptile_2 = species_2.ParticlesAt(lev, mfi);
387 
388  // Find the particles that are in each cell of this tile
389  ParticleBins bins_1 = findParticlesInEachCell( lev, mfi, ptile_1 );
390  ParticleBins bins_2 = findParticlesInEachCell( lev, mfi, ptile_2 );
391 
392  // Loop over cells, and collide the particles in each cell
393 
394  // Extract low-level data
395  auto const n_cells = static_cast<int>(bins_1.numBins());
396  // - Species 1
397  const auto soa_1 = ptile_1.getParticleTileData();
398  index_type* AMREX_RESTRICT indices_1 = bins_1.permutationPtr();
399  index_type const* AMREX_RESTRICT cell_offsets_1 = bins_1.offsetsPtr();
400  const amrex::ParticleReal q1 = species_1.getCharge();
401  const amrex::ParticleReal m1 = species_1.getMass();
402  auto get_position_1 = GetParticlePosition<PIdx>(ptile_1, getpos_offset);
403  // Needed to access the particle id
405  particle_ptr_1 = ptile_1.GetArrayOfStructs()().data();
406  // - Species 2
407  const auto soa_2 = ptile_2.getParticleTileData();
408  index_type* AMREX_RESTRICT indices_2 = bins_2.permutationPtr();
409  index_type const* AMREX_RESTRICT cell_offsets_2 = bins_2.offsetsPtr();
410  const amrex::ParticleReal q2 = species_2.getCharge();
411  const amrex::ParticleReal m2 = species_2.getMass();
412  auto get_position_2 = GetParticlePosition<PIdx>(ptile_2, getpos_offset);
413  // Needed to access the particle id
415  particle_ptr_2 = ptile_2.GetArrayOfStructs()().data();
416 
417  amrex::Geometry const& geom = WarpX::GetInstance().Geom(lev);
418 #if defined WARPX_DIM_1D_Z
419  auto dV = geom.CellSize(0);
420 #elif defined WARPX_DIM_XZ
421  auto dV = geom.CellSize(0) * geom.CellSize(1);
422 #elif defined WARPX_DIM_RZ
423  amrex::Box const& cbx = mfi.tilebox(amrex::IntVect::TheZeroVector()); //Cell-centered box
424  const auto lo = lbound(cbx);
425  const auto hi = ubound(cbx);
426  int nz = hi.y-lo.y+1;
427  auto dr = geom.CellSize(0);
428  auto dz = geom.CellSize(1);
429 #elif defined(WARPX_DIM_3D)
430  auto dV = geom.CellSize(0) * geom.CellSize(1) * geom.CellSize(2);
431 #endif
432 
433 
434  /*
435  The following calculations are only required when creating product particles
436  */
437  const int n_cells_products = have_product_species ? n_cells: 0;
438  amrex::Gpu::DeviceVector<index_type> n_pairs_in_each_cell(n_cells_products);
439  index_type* AMREX_RESTRICT p_n_pairs_in_each_cell = n_pairs_in_each_cell.dataPtr();
440 
441  // Compute how many pairs in each cell and store in n_pairs_in_each_cell array
442  // For different species, the number of pairs in a cell is the number of particles of
443  // the species that has the most particles in that cell
444  amrex::ParallelFor( n_cells_products,
445  [=] AMREX_GPU_DEVICE (int i_cell) noexcept
446  {
447  const auto n_part_in_cell_1 = cell_offsets_1[i_cell+1] - cell_offsets_1[i_cell];
448  const auto n_part_in_cell_2 = cell_offsets_2[i_cell+1] - cell_offsets_2[i_cell];
449  // Particular case: no pair if a species has no particle in that cell
450  if (n_part_in_cell_1 == 0 || n_part_in_cell_2 == 0)
451  p_n_pairs_in_each_cell[i_cell] = 0;
452  else
453  p_n_pairs_in_each_cell[i_cell] =
454  amrex::max(n_part_in_cell_1,n_part_in_cell_2);
455  }
456  );
457 
458  // Start indices of the pairs in a cell. Will be used for particle creation
459  amrex::Gpu::DeviceVector<index_type> pair_offsets(n_cells_products);
460  const index_type n_total_pairs = (n_cells_products == 0) ? 0:
461  amrex::Scan::ExclusiveSum(n_cells_products,
462  p_n_pairs_in_each_cell, pair_offsets.data());
463  index_type* AMREX_RESTRICT p_pair_offsets = pair_offsets.dataPtr();
464 
465  // mask: equal to 1 if particle creation occurs for a given pair, 0 otherwise
468  // Will be filled with the index of the first particle of a given pair
469  amrex::Gpu::DeviceVector<index_type> pair_indices_1(n_total_pairs);
470  index_type* AMREX_RESTRICT p_pair_indices_1 = pair_indices_1.dataPtr();
471  // Will be filled with the index of the second particle of a given pair
472  amrex::Gpu::DeviceVector<index_type> pair_indices_2(n_total_pairs);
473  index_type* AMREX_RESTRICT p_pair_indices_2 = pair_indices_2.dataPtr();
474  // How much weight should be given to the produced particles (and removed from the
475  // reacting particles)
476  amrex::Gpu::DeviceVector<amrex::ParticleReal> pair_reaction_weight(n_total_pairs);
477  amrex::ParticleReal* AMREX_RESTRICT p_pair_reaction_weight =
478  pair_reaction_weight.dataPtr();
479  /*
480  End of calculations only required when creating product particles
481  */
482 
483 
484  // Loop over cells
485  amrex::ParallelForRNG( n_cells,
486  [=] AMREX_GPU_DEVICE (int i_cell, amrex::RandomEngine const& engine) noexcept
487  {
488  // The particles from species1 that are in the cell `i_cell` are
489  // given by the `indices_1[cell_start_1:cell_stop_1]`
490  index_type const cell_start_1 = cell_offsets_1[i_cell];
491  index_type const cell_stop_1 = cell_offsets_1[i_cell+1];
492  // Same for species 2
493  index_type const cell_start_2 = cell_offsets_2[i_cell];
494  index_type const cell_stop_2 = cell_offsets_2[i_cell+1];
495  // Same but for the pairs
496  index_type const cell_start_pair = have_product_species?
497  p_pair_offsets[i_cell] : 0;
498 
499  // ux from species1 can be accessed like this:
500  // ux_1[ indices_1[i] ], where i is between
501  // cell_start_1 (inclusive) and cell_start_2 (exclusive)
502 
503  // Do not collide if one species is missing in the cell
504  if ( cell_stop_1 - cell_start_1 < 1 ||
505  cell_stop_2 - cell_start_2 < 1 ) return;
506 
507  // shuffle
508  ShuffleFisherYates(indices_1, cell_start_1, cell_stop_1, engine);
509  ShuffleFisherYates(indices_2, cell_start_2, cell_stop_2, engine);
510 #if defined WARPX_DIM_RZ
511  const int ri = (i_cell - i_cell%nz) / nz;
512  auto dV = MathConst::pi*(2.0_prt*ri+1.0_prt)*dr*dr*dz;
513 #endif
514  // Call the function in order to perform collisions
515  // If there are product species, p_mask, p_pair_indices_1/2, and
516  // p_pair_reaction_weight are filled here
517  binary_collision_functor(
518  cell_start_1, cell_stop_1, cell_start_2, cell_stop_2,
519  indices_1, indices_2,
520  soa_1, soa_2, get_position_1, get_position_2,
521  q1, q2, m1, m2, dt, dV,
522  cell_start_pair, p_mask, p_pair_indices_1, p_pair_indices_2,
523  p_pair_reaction_weight, engine );
524  }
525  );
526 
527 
528  // Create the new product particles and define their initial values
529  // num_added: how many particles of each product species have been created
530  const amrex::Vector<int> num_added = m_copy_transform_functor(n_total_pairs,
531  soa_1, soa_2, tile_products_data,
532  particle_ptr_1, particle_ptr_2, m1, m2,
533  products_mass, p_mask, products_np,
534  copy_species1, copy_species2,
535  p_pair_indices_1, p_pair_indices_2,
536  p_pair_reaction_weight);
537 
538  for (int i = 0; i < n_product_species; i++)
539  {
540  setNewParticleIDs(*(tile_products_data[i]), static_cast<int>(products_np[i]), num_added[i]);
541  }
542 
543  } // end if ( m_isSameSpecies)
544 
545  }
546 
547 private:
548 
552  // functor that performs collisions within a cell
553  CollisionFunctorType m_binary_collision_functor;
554  // functor that creates new particles and initializes their parameters
555  CopyTransformFunctorType m_copy_transform_functor;
556 
557 };
558 
559 #endif // WARPX_PARTICLES_COLLISION_BINARYCOLLISION_H_
#define AMREX_RESTRICT
#define AMREX_GPU_DEVICE
Array4< int const > mask
AMREX_GPU_HOST_DEVICE AMREX_INLINE void ShuffleFisherYates(T_index *array, T_index const is, T_index const ie, amrex::RandomEngine const &engine)
Definition: ShuffleFisherYates.H:20
void setNewParticleIDs(PTile &ptile, int old_size, int num_added)
Sets the ids of newly created particles to the next values.
Definition: SmartUtils.H:51
#define WARPX_ABORT_WITH_MESSAGE(MSG)
Definition: TextMsg.H:15
This class performs generic binary collisions.
Definition: BinaryCollision.H:71
CopyTransformFunctorType m_copy_transform_functor
Definition: BinaryCollision.H:555
BinaryCollision(BinaryCollision &&)=default
BinaryCollision(BinaryCollision const &)=default
bool m_isSameSpecies
Definition: BinaryCollision.H:549
CollisionFunctorType m_binary_collision_functor
Definition: BinaryCollision.H:553
void doCollisionsWithinTile(amrex::Real dt, int const lev, amrex::MFIter const &mfi, WarpXParticleContainer &species_1, WarpXParticleContainer &species_2, amrex::Vector< WarpXParticleContainer * > product_species_vector, SmartCopy *copy_species1, SmartCopy *copy_species2)
Definition: BinaryCollision.H:213
bool m_have_product_species
Definition: BinaryCollision.H:550
WarpXParticleContainer::ParticleType ParticleType
Definition: BinaryCollision.H:73
ParticleBins::index_type index_type
Definition: BinaryCollision.H:77
amrex::Vector< std::string > m_product_species
Definition: BinaryCollision.H:551
void doCollisions(amrex::Real cur_time, amrex::Real dt, MultiParticleContainer *mypc) override
Definition: BinaryCollision.H:121
BinaryCollision(std::string collision_name, MultiParticleContainer const *const mypc)
Constructor of the BinaryCollision class.
Definition: BinaryCollision.H:87
~BinaryCollision() override=default
BinaryCollision & operator=(BinaryCollision const &)=default
Definition: CollisionBase.H:18
amrex::Vector< std::string > m_species_names
Definition: CollisionBase.H:36
Definition: MultiParticleContainer.H:65
WarpXParticleContainer & GetParticleContainerFromName(const std::string &name) const
Definition: MultiParticleContainer.cpp:410
This class does nothing and is used as second template parameter for binary collisions that do not cr...
Definition: ParticleCreationFunc.H:275
A factory for creating SmartCopy functors.
Definition: SmartCopy.H:132
static WarpX & GetInstance()
Definition: WarpX.cpp:243
static amrex::LayoutData< amrex::Real > * getCosts(int lev)
Definition: WarpX.cpp:2941
static short load_balance_costs_update_algo
Definition: WarpX.H:179
Definition: WarpXParticleContainer.H:104
void defineAllParticleTiles() noexcept
Definition: WarpXParticleContainer.cpp:1293
amrex::ParticleReal getCharge() const
Definition: WarpXParticleContainer.H:353
amrex::ParticleReal getMass() const
Definition: WarpXParticleContainer.H:355
const Vector< Geometry > & Geom() const noexcept
const Real * CellSize() const noexcept
index_type * offsetsPtr() noexcept
index_type * permutationPtr() noexcept
unsigned int index_type
Long numBins() const noexcept
AMREX_GPU_HOST_DEVICE static constexpr AMREX_FORCE_INLINE IntVect TheZeroVector() noexcept
Box tilebox() const noexcept
T * data() noexcept
iterator begin() noexcept
T * dataPtr() noexcept
int queryarr(const char *name, std::vector< int > &ref, int start_ix=FIRST, int num_val=ALL) const
ParticleTile< ParticleType, NArrayReal, NArrayInt, Allocator > ParticleTileType
const ParticleTileType & ParticlesAt(int lev, int grid, int tile) const
Definition: ParticleUtils.cpp:25
ParticleBins findParticlesInEachCell(int const lev, MFIter const &mfi, ParticleTileType const &ptile)
Find the particles and count the particles that are in each cell. More specifically this function ret...
Definition: ParticleUtils.cpp:37
void synchronize() noexcept
void copyAsync(HostToDevice, InIter begin, InIter end, OutIter result) noexcept
static constexpr HostToDevice hostToDevice
void streamSynchronize() noexcept
bool notInLaunchRegion() noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void Add(T *const sum, T const value) noexcept
T ExclusiveSum(N n, T const *in, T *out, RetSum a_ret_sum=retSum)
void ParallelForRNG(T n, L &&f) noexcept
AMREX_GPU_HOST_DEVICE constexpr AMREX_FORCE_INLINE const T & max(const T &a, const T &b) noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Dim3 ubound(Array4< T > const &a) noexcept
double second() noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Dim3 lbound(Array4< T > const &a) noexcept
std::enable_if_t< std::is_integral< T >::value > ParallelFor(TypeList< CTOs... >, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)
i
Definition: check_interp_points_and_weights.py:174
int dz
Definition: compute_domain.py:36
data
Definition: run_alltests_1node.py:325
float dt
Definition: stencil.py:442
value
Definition: updateAMReX.py:141
string species1
Definition: video_yt.py:35
@ Timers
load balance according to in-code timer-based weights (i.e., with costs)
Definition: WarpXAlgorithmSelection.H:122
This is a functor for performing a "smart copy" that works in both host and device code.
Definition: SmartCopy.H:34
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int const * dataPtr() const noexcept
MFItInfo & EnableTiling(const IntVect &ts=FabArrayBase::mfiter_tile_size) noexcept
MFItInfo & SetDynamic(bool f) noexcept
ParticleTileData< StorageParticleType, NArrayReal, NArrayInt > ParticleTileDataType
ParticleTileDataType getParticleTileData()
int numParticles() const