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