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