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  }
93 
95 
96  m_binary_collision_functor = CollisionFunctorType(collision_name, mypc, m_isSameSpecies);
97 
98  const amrex::ParmParse pp_collision_name(collision_name);
99  pp_collision_name.queryarr("product_species", m_product_species);
102  WARPX_ABORT_WITH_MESSAGE( "Binary collision " + collision_name +
103  " does not produce species. Thus, `product_species` should not be specified in the input script." );
104  }
105  m_copy_transform_functor = CopyTransformFunctorType(collision_name, mypc);
106  }
107 
108  ~BinaryCollision () override = default;
109 
110  BinaryCollision ( BinaryCollision const &) = default;
113  BinaryCollision& operator= ( BinaryCollision&& ) noexcept = default;
114 
122  void doCollisions (amrex::Real cur_time, amrex::Real dt, MultiParticleContainer* mypc) override
123  {
124  amrex::ignore_unused(cur_time);
125 
126  auto& species1 = mypc->GetParticleContainerFromName(m_species_names[0]);
127  auto& species2 = mypc->GetParticleContainerFromName(m_species_names[1]);
128 
129  // In case of particle creation, create the necessary vectors
130  const int n_product_species = m_product_species.size();
131  amrex::Vector<WarpXParticleContainer*> product_species_vector;
132  amrex::Vector<SmartCopyFactory> copy_factory_species1;
133  amrex::Vector<SmartCopyFactory> copy_factory_species2;
134  amrex::Vector<SmartCopy> copy_species1;
135  amrex::Vector<SmartCopy> copy_species2;
136  for (int i = 0; i < n_product_species; i++)
137  {
138  auto& product = mypc->GetParticleContainerFromName(m_product_species[i]);
139  product.defineAllParticleTiles();
140  product_species_vector.push_back(&product);
141  // Although the copy factories are not explicitly reused past this point, we need to
142  // store them in vectors so that the data that they own, which is used by the smart
143  // copy functors, does not go out of scope at the end of this for loop.
144  copy_factory_species1.push_back(SmartCopyFactory(species1, product));
145  copy_factory_species2.push_back(SmartCopyFactory(species2, product));
146  copy_species1.push_back(copy_factory_species1[i].getSmartCopy());
147  copy_species2.push_back(copy_factory_species2[i].getSmartCopy());
148  }
149 #ifdef AMREX_USE_GPU
150  amrex::Gpu::DeviceVector<SmartCopy> device_copy_species1(n_product_species);
151  amrex::Gpu::DeviceVector<SmartCopy> device_copy_species2(n_product_species);
152  amrex::Gpu::copyAsync(amrex::Gpu::hostToDevice, copy_species1.begin(),
153  copy_species1.end(), device_copy_species1.begin());
154  amrex::Gpu::copyAsync(amrex::Gpu::hostToDevice, copy_species2.begin(),
155  copy_species2.end(), device_copy_species2.begin());
157  auto copy_species1_data = device_copy_species1.data();
158  auto copy_species2_data = device_copy_species2.data();
159 #else
160  auto *copy_species1_data = copy_species1.data();
161  auto *copy_species2_data = copy_species2.data();
162 #endif
164  species1.defineAllParticleTiles();
165  species2.defineAllParticleTiles();
166  }
167 
168  // Enable tiling
169  amrex::MFItInfo info;
170  if (amrex::Gpu::notInLaunchRegion()) { info.EnableTiling(species1.tile_size); }
171 
172  // Loop over refinement levels
173  for (int lev = 0; lev <= species1.finestLevel(); ++lev){
174 
176 
177  // Loop over all grids/tiles at this level
178 #ifdef AMREX_USE_OMP
179  info.SetDynamic(true);
180 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
181 #endif
182  for (amrex::MFIter mfi = species1.MakeMFIter(lev, info); mfi.isValid(); ++mfi){
184  {
186  }
187  auto wt = static_cast<amrex::Real>(amrex::second());
188 
189  doCollisionsWithinTile( dt, lev, mfi, species1, species2, product_species_vector,
190  copy_species1_data, copy_species2_data);
191 
193  {
195  wt = static_cast<amrex::Real>(amrex::second()) - wt;
196  amrex::HostDevice::Atomic::Add( &(*cost)[mfi.index()], wt);
197  }
198  }
199  }
200  }
201 
215  amrex::Real dt, int const lev, amrex::MFIter const& mfi,
216  WarpXParticleContainer& species_1,
217  WarpXParticleContainer& species_2,
218  amrex::Vector<WarpXParticleContainer*> product_species_vector,
219  SmartCopy* copy_species1, SmartCopy* copy_species2)
220  {
221  using namespace ParticleUtils;
222  using namespace amrex::literals;
223 
224  CollisionFunctorType binary_collision_functor = m_binary_collision_functor;
225  const bool have_product_species = m_have_product_species;
226 
227  // Store product species data in vectors
228  const int n_product_species = m_product_species.size();
229  amrex::Vector<ParticleTileType*> tile_products;
230  amrex::Vector<GetParticlePosition<PIdx>> get_position_products;
231  amrex::Vector<index_type> products_np;
233  constexpr int getpos_offset = 0;
234  for (int i = 0; i < n_product_species; i++)
235  {
236  ParticleTileType& ptile_product = product_species_vector[i]->ParticlesAt(lev, mfi);
237  tile_products.push_back(&ptile_product);
238  get_position_products.push_back(GetParticlePosition<PIdx>(ptile_product,
239  getpos_offset));
240  products_np.push_back(ptile_product.numParticles());
241  products_mass.push_back(product_species_vector[i]->getMass());
242  }
243  auto *tile_products_data = tile_products.data();
244 
245  if ( m_isSameSpecies ) // species_1 == species_2
246  {
247  // Extract particles in the tile that `mfi` points to
248  ParticleTileType& ptile_1 = species_1.ParticlesAt(lev, mfi);
249 
250  // Find the particles that are in each cell of this tile
251  ParticleBins bins_1 = findParticlesInEachCell( lev, mfi, ptile_1 );
252 
253  // Loop over cells, and collide the particles in each cell
254 
255  // Extract low-level data
256  auto const n_cells = static_cast<int>(bins_1.numBins());
257  // - Species 1
258  const auto soa_1 = ptile_1.getParticleTileData();
259  index_type* AMREX_RESTRICT indices_1 = bins_1.permutationPtr();
260  index_type const* AMREX_RESTRICT cell_offsets_1 = bins_1.offsetsPtr();
261  const amrex::ParticleReal q1 = species_1.getCharge();
262  const amrex::ParticleReal m1 = species_1.getMass();
263  auto get_position_1 = GetParticlePosition<PIdx>(ptile_1, getpos_offset);
264  // Needed to access the particle id
266  particle_ptr_1 = ptile_1.GetArrayOfStructs()().data();
267 
268  amrex::Geometry const& geom = WarpX::GetInstance().Geom(lev);
269 #if defined WARPX_DIM_1D_Z
270  auto dV = geom.CellSize(0);
271 #elif defined WARPX_DIM_XZ
272  auto dV = geom.CellSize(0) * geom.CellSize(1);
273 #elif defined WARPX_DIM_RZ
274  amrex::Box const& cbx = mfi.tilebox(amrex::IntVect::TheZeroVector()); //Cell-centered box
275  const auto lo = lbound(cbx);
276  const auto hi = ubound(cbx);
277  int const nz = hi.y-lo.y+1;
278  auto dr = geom.CellSize(0);
279  auto dz = geom.CellSize(1);
280 #elif defined(WARPX_DIM_3D)
281  auto dV = geom.CellSize(0) * geom.CellSize(1) * geom.CellSize(2);
282 #endif
283 
284 
285  /*
286  The following calculations are only required when creating product particles
287  */
288  const int n_cells_products = have_product_species ? n_cells: 0;
289  amrex::Gpu::DeviceVector<index_type> n_pairs_in_each_cell(n_cells_products);
290  index_type* AMREX_RESTRICT p_n_pairs_in_each_cell = n_pairs_in_each_cell.dataPtr();
291 
292  // Compute how many pairs in each cell and store in n_pairs_in_each_cell array
293  // For a single species, the number of pair in a cell is half the number of particles
294  // in that cell, rounded up to the next higher integer.
295  amrex::ParallelFor( n_cells_products,
296  [=] AMREX_GPU_DEVICE (int i_cell) noexcept
297  {
298  const auto n_part_in_cell = cell_offsets_1[i_cell+1] - cell_offsets_1[i_cell];
299  // Particular case: if there's only 1 particle in a cell, then there's no pair
300  p_n_pairs_in_each_cell[i_cell] = (n_part_in_cell == 1)? 0: (n_part_in_cell+1)/2;
301  }
302  );
303 
304  // Start indices of the pairs in a cell. Will be used for particle creation.
305  amrex::Gpu::DeviceVector<index_type> pair_offsets(n_cells_products);
306  const index_type n_total_pairs = (n_cells_products == 0) ? 0:
307  amrex::Scan::ExclusiveSum(n_cells_products,
308  p_n_pairs_in_each_cell, pair_offsets.data());
309  index_type* AMREX_RESTRICT p_pair_offsets = pair_offsets.dataPtr();
310 
311  // mask: equal to 1 if particle creation occurs for a given pair, 0 otherwise
314  // Will be filled with the index of the first particle of a given pair
315  amrex::Gpu::DeviceVector<index_type> pair_indices_1(n_total_pairs);
316  index_type* AMREX_RESTRICT p_pair_indices_1 = pair_indices_1.dataPtr();
317  // Will be filled with the index of the second particle of a given pair
318  amrex::Gpu::DeviceVector<index_type> pair_indices_2(n_total_pairs);
319  index_type* AMREX_RESTRICT p_pair_indices_2 = pair_indices_2.dataPtr();
320  // How much weight should be given to the produced particles (and removed from the
321  // reacting particles)
322  amrex::Gpu::DeviceVector<amrex::ParticleReal> pair_reaction_weight(n_total_pairs);
323  amrex::ParticleReal* AMREX_RESTRICT p_pair_reaction_weight =
324  pair_reaction_weight.dataPtr();
325  /*
326  End of calculations only required when creating product particles
327  */
328 
329 
330  // Loop over cells
331  amrex::ParallelForRNG( n_cells,
332  [=] AMREX_GPU_DEVICE (int i_cell, amrex::RandomEngine const& engine) noexcept
333  {
334  // The particles from species1 that are in the cell `i_cell` are
335  // given by the `indices_1[cell_start_1:cell_stop_1]`
336  index_type const cell_start_1 = cell_offsets_1[i_cell];
337  index_type const cell_stop_1 = cell_offsets_1[i_cell+1];
338  index_type const cell_half_1 = (cell_start_1+cell_stop_1)/2;
339 
340  // Same but for the pairs
341  index_type const cell_start_pair = have_product_species?
342  p_pair_offsets[i_cell] : 0;
343 
344  // Do not collide if there is only one particle in the cell
345  if ( cell_stop_1 - cell_start_1 <= 1 ) { return; }
346 
347  // shuffle
349  indices_1, cell_start_1, cell_half_1, engine );
350 #if defined WARPX_DIM_RZ
351  const int ri = (i_cell - i_cell%nz) / nz;
352  auto dV = MathConst::pi*(2.0_prt*ri+1.0_prt)*dr*dr*dz;
353 #endif
354  // Call the function in order to perform collisions
355  // If there are product species, mask, p_pair_indices_1/2, and
356  // p_pair_reaction_weight are filled here
357  binary_collision_functor(
358  cell_start_1, cell_half_1,
359  cell_half_1, cell_stop_1,
360  indices_1, indices_1,
361  soa_1, soa_1, get_position_1, get_position_1,
362  q1, q1, m1, m1, dt, dV,
363  cell_start_pair, p_mask, p_pair_indices_1, p_pair_indices_2,
364  p_pair_reaction_weight, engine );
365  }
366  );
367 
368  // Create the new product particles and define their initial values
369  // num_added: how many particles of each product species have been created
370  const amrex::Vector<int> num_added = m_copy_transform_functor(n_total_pairs,
371  soa_1, soa_1,
372  product_species_vector,
373  tile_products_data,
374  particle_ptr_1, particle_ptr_1, m1, m1,
375  products_mass, p_mask, products_np,
376  copy_species1, copy_species2,
377  p_pair_indices_1, p_pair_indices_2,
378  p_pair_reaction_weight);
379 
380  for (int i = 0; i < n_product_species; i++)
381  {
382  setNewParticleIDs(*(tile_products_data[i]), static_cast<int>(products_np[i]), num_added[i]);
383  }
384  }
385  else // species_1 != species_2
386  {
387  // Extract particles in the tile that `mfi` points to
388  ParticleTileType& ptile_1 = species_1.ParticlesAt(lev, mfi);
389  ParticleTileType& ptile_2 = species_2.ParticlesAt(lev, mfi);
390 
391  // Find the particles that are in each cell of this tile
392  ParticleBins bins_1 = findParticlesInEachCell( lev, mfi, ptile_1 );
393  ParticleBins bins_2 = findParticlesInEachCell( lev, mfi, ptile_2 );
394 
395  // Loop over cells, and collide the particles in each cell
396 
397  // Extract low-level data
398  auto const n_cells = static_cast<int>(bins_1.numBins());
399  // - Species 1
400  const auto soa_1 = ptile_1.getParticleTileData();
401  index_type* AMREX_RESTRICT indices_1 = bins_1.permutationPtr();
402  index_type const* AMREX_RESTRICT cell_offsets_1 = bins_1.offsetsPtr();
403  const amrex::ParticleReal q1 = species_1.getCharge();
404  const amrex::ParticleReal m1 = species_1.getMass();
405  auto get_position_1 = GetParticlePosition<PIdx>(ptile_1, getpos_offset);
406  // Needed to access the particle id
408  particle_ptr_1 = ptile_1.GetArrayOfStructs()().data();
409  // - Species 2
410  const auto soa_2 = ptile_2.getParticleTileData();
411  index_type* AMREX_RESTRICT indices_2 = bins_2.permutationPtr();
412  index_type const* AMREX_RESTRICT cell_offsets_2 = bins_2.offsetsPtr();
413  const amrex::ParticleReal q2 = species_2.getCharge();
414  const amrex::ParticleReal m2 = species_2.getMass();
415  auto get_position_2 = GetParticlePosition<PIdx>(ptile_2, getpos_offset);
416  // Needed to access the particle id
418  particle_ptr_2 = ptile_2.GetArrayOfStructs()().data();
419 
420  amrex::Geometry const& geom = WarpX::GetInstance().Geom(lev);
421 #if defined WARPX_DIM_1D_Z
422  auto dV = geom.CellSize(0);
423 #elif defined WARPX_DIM_XZ
424  auto dV = geom.CellSize(0) * geom.CellSize(1);
425 #elif defined WARPX_DIM_RZ
426  amrex::Box const& cbx = mfi.tilebox(amrex::IntVect::TheZeroVector()); //Cell-centered box
427  const auto lo = lbound(cbx);
428  const auto hi = ubound(cbx);
429  int nz = hi.y-lo.y+1;
430  auto dr = geom.CellSize(0);
431  auto dz = geom.CellSize(1);
432 #elif defined(WARPX_DIM_3D)
433  auto dV = geom.CellSize(0) * geom.CellSize(1) * geom.CellSize(2);
434 #endif
435 
436 
437  /*
438  The following calculations are only required when creating product particles
439  */
440  const int n_cells_products = have_product_species ? n_cells: 0;
441  amrex::Gpu::DeviceVector<index_type> n_pairs_in_each_cell(n_cells_products);
442  index_type* AMREX_RESTRICT p_n_pairs_in_each_cell = n_pairs_in_each_cell.dataPtr();
443 
444  // Compute how many pairs in each cell and store in n_pairs_in_each_cell array
445  // For different species, the number of pairs in a cell is the number of particles of
446  // the species that has the most particles in that cell
447  amrex::ParallelFor( n_cells_products,
448  [=] AMREX_GPU_DEVICE (int i_cell) noexcept
449  {
450  const auto n_part_in_cell_1 = cell_offsets_1[i_cell+1] - cell_offsets_1[i_cell];
451  const auto n_part_in_cell_2 = cell_offsets_2[i_cell+1] - cell_offsets_2[i_cell];
452  // Particular case: no pair if a species has no particle in that cell
453  if (n_part_in_cell_1 == 0 || n_part_in_cell_2 == 0) {
454  p_n_pairs_in_each_cell[i_cell] = 0;
455  } else {
456  p_n_pairs_in_each_cell[i_cell] =
457  amrex::max(n_part_in_cell_1,n_part_in_cell_2);
458  }
459  }
460  );
461 
462  // Start indices of the pairs in a cell. Will be used for particle creation
463  amrex::Gpu::DeviceVector<index_type> pair_offsets(n_cells_products);
464  const index_type n_total_pairs = (n_cells_products == 0) ? 0:
465  amrex::Scan::ExclusiveSum(n_cells_products,
466  p_n_pairs_in_each_cell, pair_offsets.data());
467  index_type* AMREX_RESTRICT p_pair_offsets = pair_offsets.dataPtr();
468 
469  // mask: equal to 1 if particle creation occurs for a given pair, 0 otherwise
472  // Will be filled with the index of the first particle of a given pair
473  amrex::Gpu::DeviceVector<index_type> pair_indices_1(n_total_pairs);
474  index_type* AMREX_RESTRICT p_pair_indices_1 = pair_indices_1.dataPtr();
475  // Will be filled with the index of the second particle of a given pair
476  amrex::Gpu::DeviceVector<index_type> pair_indices_2(n_total_pairs);
477  index_type* AMREX_RESTRICT p_pair_indices_2 = pair_indices_2.dataPtr();
478  // How much weight should be given to the produced particles (and removed from the
479  // reacting particles)
480  amrex::Gpu::DeviceVector<amrex::ParticleReal> pair_reaction_weight(n_total_pairs);
481  amrex::ParticleReal* AMREX_RESTRICT p_pair_reaction_weight =
482  pair_reaction_weight.dataPtr();
483  /*
484  End of calculations only required when creating product particles
485  */
486 
487 
488  // Loop over cells
489  amrex::ParallelForRNG( n_cells,
490  [=] AMREX_GPU_DEVICE (int i_cell, amrex::RandomEngine const& engine) noexcept
491  {
492  // The particles from species1 that are in the cell `i_cell` are
493  // given by the `indices_1[cell_start_1:cell_stop_1]`
494  index_type const cell_start_1 = cell_offsets_1[i_cell];
495  index_type const cell_stop_1 = cell_offsets_1[i_cell+1];
496  // Same for species 2
497  index_type const cell_start_2 = cell_offsets_2[i_cell];
498  index_type const cell_stop_2 = cell_offsets_2[i_cell+1];
499  // Same but for the pairs
500  index_type const cell_start_pair = have_product_species?
501  p_pair_offsets[i_cell] : 0;
502 
503  // ux from species1 can be accessed like this:
504  // ux_1[ indices_1[i] ], where i is between
505  // cell_start_1 (inclusive) and cell_start_2 (exclusive)
506 
507  // Do not collide if one species is missing in the cell
508  if ( cell_stop_1 - cell_start_1 < 1 ||
509  cell_stop_2 - cell_start_2 < 1 ) { return; }
510 
511  // shuffle
512  ShuffleFisherYates(indices_1, cell_start_1, cell_stop_1, engine);
513  ShuffleFisherYates(indices_2, cell_start_2, cell_stop_2, engine);
514 #if defined WARPX_DIM_RZ
515  const int ri = (i_cell - i_cell%nz) / nz;
516  auto dV = MathConst::pi*(2.0_prt*ri+1.0_prt)*dr*dr*dz;
517 #endif
518  // Call the function in order to perform collisions
519  // If there are product species, p_mask, p_pair_indices_1/2, and
520  // p_pair_reaction_weight are filled here
521  binary_collision_functor(
522  cell_start_1, cell_stop_1, cell_start_2, cell_stop_2,
523  indices_1, indices_2,
524  soa_1, soa_2, get_position_1, get_position_2,
525  q1, q2, m1, m2, dt, dV,
526  cell_start_pair, p_mask, p_pair_indices_1, p_pair_indices_2,
527  p_pair_reaction_weight, engine );
528  }
529  );
530 
531 
532  // Create the new product particles and define their initial values
533  // num_added: how many particles of each product species have been created
534  const amrex::Vector<int> num_added = m_copy_transform_functor(n_total_pairs,
535  soa_1, soa_2,
536  product_species_vector,
537  tile_products_data,
538  particle_ptr_1, particle_ptr_2, m1, m2,
539  products_mass, p_mask, products_np,
540  copy_species1, copy_species2,
541  p_pair_indices_1, p_pair_indices_2,
542  p_pair_reaction_weight);
543 
544  for (int i = 0; i < n_product_species; i++)
545  {
546  setNewParticleIDs(*(tile_products_data[i]), static_cast<int>(products_np[i]), num_added[i]);
547  }
548 
549  } // end if ( m_isSameSpecies)
550 
551  }
552 
553 private:
554 
558  // functor that performs collisions within a cell
559  CollisionFunctorType m_binary_collision_functor;
560  // functor that creates new particles and initializes their parameters
561  CopyTransformFunctorType m_copy_transform_functor;
562 
563 };
564 
565 #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:561
BinaryCollision(BinaryCollision const &)=default
bool m_isSameSpecies
Definition: BinaryCollision.H:555
CollisionFunctorType m_binary_collision_functor
Definition: BinaryCollision.H:559
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:214
bool m_have_product_species
Definition: BinaryCollision.H:556
WarpXParticleContainer::ParticleType ParticleType
Definition: BinaryCollision.H:73
ParticleBins::index_type index_type
Definition: BinaryCollision.H:77
BinaryCollision(BinaryCollision &&) noexcept=default
amrex::Vector< std::string > m_product_species
Definition: BinaryCollision.H:557
void doCollisions(amrex::Real cur_time, amrex::Real dt, MultiParticleContainer *mypc) override
Definition: BinaryCollision.H:122
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:66
This class does nothing and is used as second template parameter for binary collisions that do not cr...
Definition: ParticleCreationFunc.H:297
A factory for creating SmartCopy functors.
Definition: SmartCopy.H:136
static WarpX & GetInstance()
Definition: WarpX.cpp:239
static amrex::LayoutData< amrex::Real > * getCosts(int lev)
Definition: WarpX.cpp:2996
static short load_balance_costs_update_algo
Definition: WarpX.H:178
Definition: WarpXParticleContainer.H:110
amrex::ParticleReal getCharge() const
Definition: WarpXParticleContainer.H:364
amrex::ParticleReal getMass() const
Definition: WarpXParticleContainer.H:366
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
default
Definition: run_alltests.py:113
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:133
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