WarpX
WarpXParticleContainer.H
Go to the documentation of this file.
1 /* Copyright 2019-2020 Andrew Myers, Axel Huebl, David Grote
2  * Jean-Luc Vay, Junmin Gu, Luca Fedeli
3  * Maxence Thevenet, Remi Lehe, Revathi Jambunathan
4  * Weiqun Zhang, Yinjian Zhao
5  *
6  * This file is part of WarpX.
7  *
8  * License: BSD-3-Clause-LBNL
9  */
10 #ifndef WARPX_WarpXParticleContainer_H_
11 #define WARPX_WarpXParticleContainer_H_
12 
14 
15 #include "Evolve/WarpXDtType.H"
19 
20 #ifdef WARPX_QED
23 #endif
26 
27 #include <AMReX_Array.H>
28 #include <AMReX_FArrayBox.H>
29 #include <AMReX_GpuAllocators.H>
30 #include <AMReX_GpuContainers.H>
31 #include <AMReX_INT.H>
32 #include <AMReX_ParIter.H>
33 #include <AMReX_Particles.H>
34 #include <AMReX_Random.H>
35 #include <AMReX_REAL.H>
36 #include <AMReX_StructOfArrays.H>
37 #include <AMReX_Vector.H>
38 
39 #include <AMReX_BaseFwd.H>
40 #include <AMReX_AmrCoreFwd.H>
41 
42 #include <array>
43 #include <iosfwd>
44 #include <map>
45 #include <memory>
46 #include <string>
47 #include <utility>
48 
49 
51  : public amrex::ParIter<0,0,PIdx::nattribs>
52 {
53 public:
55 
56  WarpXParIter (ContainerType& pc, int level);
57 
58  WarpXParIter (ContainerType& pc, int level, amrex::MFItInfo& info);
59 
60  [[nodiscard]] const std::array<RealVector, PIdx::nattribs>& GetAttribs () const
61  {
62  return GetStructOfArrays().GetRealData();
63  }
64 
65  [[nodiscard]] std::array<RealVector, PIdx::nattribs>& GetAttribs ()
66  {
67  return GetStructOfArrays().GetRealData();
68  }
69 
70  [[nodiscard]] const RealVector& GetAttribs (int comp) const
71  {
72  return GetStructOfArrays().GetRealData(comp);
73  }
74 
75  [[nodiscard]] RealVector& GetAttribs (int comp)
76  {
77  return GetStructOfArrays().GetRealData(comp);
78  }
79 
80  [[nodiscard]] IntVector& GetiAttribs (int comp)
81  {
82  return GetStructOfArrays().GetIntData(comp);
83  }
84 };
85 
108  : public NamedComponentParticleContainer<amrex::DefaultAllocator>
109 {
110 public:
112 
113  // amrex::StructOfArrays with DiagIdx::nattribs amrex::ParticleReal components
114  // and 0 int components for the particle data.
116  // DiagnosticParticles is a vector, with one element per MR level.
117  // DiagnosticParticles[lev] is typically a key-value pair where the key is
118  // a pair [grid_index, tile_index], and the value is the corresponding
119  // DiagnosticParticleData (see above) on this tile.
121 
122  WarpXParticleContainer (amrex::AmrCore* amr_core, int ispecies);
123  ~WarpXParticleContainer() override = default;
124 
125  // Move and copy operations
130 
131  virtual void InitData () = 0;
132 
133  virtual void InitIonizationModule () {}
134 
135  /*
136  * \brief Virtual function that returns a pointer to the plasma injector,
137  * for derived classes that define one (PhysicalParticleContainer).
138  */
139  virtual PlasmaInjector* GetPlasmaInjector (const int /*i*/) { return nullptr; }
140 
146  virtual void Evolve (int lev,
147  const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
148  const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz,
151  amrex::MultiFab* rho, amrex::MultiFab* crho,
152  const amrex::MultiFab* cEx, const amrex::MultiFab* cEy, const amrex::MultiFab* cEz,
153  const amrex::MultiFab* cBx, const amrex::MultiFab* cBy, const amrex::MultiFab* cBz,
154  amrex::Real t, amrex::Real dt, DtType a_dt_type=DtType::Full, bool skip_deposition=false) = 0;
155 
156  virtual void PostRestart () = 0;
157 
158  void AllocData ();
159 
167  amrex::PinnedArenaAllocator>& pinned_tile,
168  int n_external_attr_real,
169  int n_external_attr_int,
170  const amrex::RandomEngine& engine) = 0;
171 
177  void PushX ( amrex::Real dt);
178  void PushX (int lev, amrex::Real dt);
179 
183  virtual void PushP (int lev, amrex::Real dt,
184  const amrex::MultiFab& Ex,
185  const amrex::MultiFab& Ey,
186  const amrex::MultiFab& Ez,
187  const amrex::MultiFab& Bx,
188  const amrex::MultiFab& By,
189  const amrex::MultiFab& Bz) = 0;
190 
202  void DepositCurrent (amrex::Vector<std::array< std::unique_ptr<amrex::MultiFab>, 3 > >& J,
203  amrex::Real dt, amrex::Real relative_time);
204 
215  void DepositCharge (amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
216  bool local = false, bool reset = false,
217  bool apply_boundary_and_scale_volume = false,
218  bool interpolate_across_levels = true,
219  int icomp = 0);
220  void DepositCharge (std::unique_ptr<amrex::MultiFab>& rho, int lev,
221  bool local = false, bool reset = false,
222  bool apply_boundary_and_scale_volume = false,
223  int icomp = 0);
224 
225  std::unique_ptr<amrex::MultiFab> GetChargeDensity(int lev, bool local = false);
226 
227  virtual void DepositCharge (WarpXParIter& pti,
228  RealVector const & wp,
229  const int* ion_lev,
230  amrex::MultiFab* rho,
231  int icomp,
232  long offset,
233  long np_to_depose,
234  int thread_num,
235  int lev,
236  int depos_lev);
237 
238  virtual void DepositCurrent (WarpXParIter& pti,
239  RealVector const & wp,
240  RealVector const & uxp,
241  RealVector const & uyp,
242  RealVector const & uzp,
243  int const* ion_lev,
244  amrex::MultiFab* jx,
245  amrex::MultiFab* jy,
246  amrex::MultiFab* jz,
247  long offset,
248  long np_to_depose,
249  int thread_num,
250  int lev,
251  int depos_lev,
252  amrex::Real dt,
253  amrex::Real relative_time);
254 
255  // If particles start outside of the domain, ContinuousInjection
256  // makes sure that they are initialized when they enter the domain, and
257  // NOT before. Virtual function, overriden by derived classes.
258  // Current status:
259  // PhysicalParticleContainer: implemented.
260  // LaserParticleContainer: implemented.
261  // RigidInjectedParticleContainer: not implemented.
262  virtual void ContinuousInjection(const amrex::RealBox& /*injection_box*/) {}
263 
268  virtual void UpdateAntennaPosition(const amrex::Real /*dt*/) {}
269 
271 
272  // Inject a continuous flux of particles from a defined plane
273  virtual void ContinuousFluxInjection(amrex::Real /*t*/, amrex::Real /*dt*/) {}
274 
275  int getSpeciesId() const {return species_id;}
276 
281  amrex::ParticleReal sumParticleCharge(bool local = false);
282 
283  std::array<amrex::ParticleReal, 3> meanParticleVelocity(bool local = false);
284 
285  amrex::ParticleReal maxParticleVelocity(bool local = false);
286 
311  void AddNParticles (int lev, long n,
318  int nattr_real,
320  int nattr_int, amrex::Vector<amrex::Vector<int>> const & attr_int,
321  int uniqueparticles, amrex::Long id=-1);
322 
323  virtual void ReadHeader (std::istream& is) = 0;
324 
325  virtual void WriteHeader (std::ostream& os) const = 0;
326 
327  static void ReadParameters ();
328 
329  static void BackwardCompatibility ();
330 
333  void ApplyBoundaryConditions ();
334 
335  bool do_splitting = false;
336  int do_not_deposit = 0;
338  amrex::Real self_fields_required_precision = amrex::Real(1.e-11);
339  amrex::Real self_fields_absolute_tolerance = amrex::Real(0.0);
342 
343  // External fields added to particle fields.
346 
349 
350  // split along diagonals (0) or axes (1)
351  int split_type = 0;
352 
357  void SetDoBackTransformedParticles(const bool do_back_transformed_particles) {
358  m_do_back_transformed_particles = do_back_transformed_particles;
359  }
360 
361  //amrex::Real getCharge () {return charge;}
362  amrex::ParticleReal getCharge () const {return charge;}
363  //amrex::Real getMass () {return mass;}
364  amrex::ParticleReal getMass () const {return mass;}
365 
366  int DoFieldIonization() const { return do_field_ionization; }
367 
368 #ifdef WARPX_QED
369  //Species for which QED effects are relevant should override these methods
370  virtual bool has_quantum_sync() const {return false;}
371  virtual bool has_breit_wheeler() const {return false;}
372 
373  int DoQED() const { return has_quantum_sync() || has_breit_wheeler(); }
374 #else
375  int DoQED() const { return false; }
376 #endif
377 
378  /* \brief This function tests if the current species
379  * is of a given PhysicalSpecies (specified as a template parameter).
380  * @tparam PhysSpec the PhysicalSpecies to test against
381  * @return the result of the test
382  */
383  template<PhysicalSpecies PhysSpec>
384  bool AmIA () const noexcept {return (physical_species == PhysSpec);}
385 
389  std::string getSpeciesTypeName () const {return species::get_name(physical_species);}
390 
397  virtual void resample (const int /*timestep*/, bool /*verbose*/) {}
398 
405  void defineAllParticleTiles () noexcept;
406 
407 protected:
409 
410  amrex::ParticleReal charge;
411  amrex::ParticleReal mass;
413 
414  // Controls boundaries for particles exiting the domain
416 
419 
422 
423  int do_not_push = 0;
424  int do_not_gather = 0;
425 
426  // Whether to allow particles outside of the simulation domain to be
427  // initialized when they enter the domain.
428  // This is currently required because continuous injection does not
429  // support all features allowed by direct injection.
431 
441  std::string physical_element;
442 
443  int do_resampling = 0;
444 
447 
448 #ifdef WARPX_QED
449  //Species can receive a shared pointer to a QED engine (species for
450  //which this is relevant should override these functions)
451  virtual void
452  set_breit_wheeler_engine_ptr(std::shared_ptr<BreitWheelerEngine>){}
453  virtual void
454  set_quantum_sync_engine_ptr(std::shared_ptr<QuantumSynchrotronEngine>){}
455 
462 
463 #endif
468 
469 public:
470  using PairIndex = std::pair<int, int>;
471  using TmpParticleTile = std::array<amrex::Gpu::DeviceVector<amrex::ParticleReal>,
474 
476 protected:
478 
479 private:
480  void particlePostLocate(ParticleType& p, const amrex::ParticleLocData& pld, int lev) override;
481 
482 };
483 
484 #endif
PhysicalSpecies
Definition: SpeciesPhysicalProperties.H:16
DtType
Definition: WarpXDtType.H:11
Definition: NamedComponentParticleContainer.H:49
Definition: PlasmaInjector.H:42
Definition: WarpXParticleContainer.H:52
WarpXParIter(ContainerType &pc, int level)
Definition: WarpXParticleContainer.cpp:74
std::array< RealVector, PIdx::nattribs > & GetAttribs()
Definition: WarpXParticleContainer.H:65
const RealVector & GetAttribs(int comp) const
Definition: WarpXParticleContainer.H:70
IntVector & GetiAttribs(int comp)
Definition: WarpXParticleContainer.H:80
RealVector & GetAttribs(int comp)
Definition: WarpXParticleContainer.H:75
const std::array< RealVector, PIdx::nattribs > & GetAttribs() const
Definition: WarpXParticleContainer.H:60
Definition: WarpXParticleContainer.H:109
amrex::Real self_fields_required_precision
Definition: WarpXParticleContainer.H:338
int do_continuous_injection
Definition: WarpXParticleContainer.H:430
~WarpXParticleContainer() override=default
bool AmIA() const noexcept
Definition: WarpXParticleContainer.H:384
amrex::Vector< amrex::ParticleReal > m_B_external_particle
Definition: WarpXParticleContainer.H:344
bool m_gather_from_main_grid
instead of gathering fields from the finest patch level, gather from the coarsest
Definition: WarpXParticleContainer.H:421
void AllocData()
Definition: WarpXParticleContainer.cpp:149
virtual PlasmaInjector * GetPlasmaInjector(const int)
Definition: WarpXParticleContainer.H:139
void defineAllParticleTiles() noexcept
Definition: WarpXParticleContainer.cpp:1305
virtual void PostRestart()=0
void particlePostLocate(ParticleType &p, const amrex::ParticleLocData &pld, int lev) override
Definition: WarpXParticleContainer.cpp:1322
virtual void ReadHeader(std::istream &is)=0
int self_fields_max_iters
Definition: WarpXParticleContainer.H:340
virtual void Evolve(int lev, const amrex::MultiFab &Ex, const amrex::MultiFab &Ey, const amrex::MultiFab &Ez, const amrex::MultiFab &Bx, const amrex::MultiFab &By, const amrex::MultiFab &Bz, amrex::MultiFab &jx, amrex::MultiFab &jy, amrex::MultiFab &jz, amrex::MultiFab *cjx, amrex::MultiFab *cjy, amrex::MultiFab *cjz, amrex::MultiFab *rho, amrex::MultiFab *crho, const amrex::MultiFab *cEx, const amrex::MultiFab *cEy, const amrex::MultiFab *cEz, const amrex::MultiFab *cBx, const amrex::MultiFab *cBy, const amrex::MultiFab *cBz, amrex::Real t, amrex::Real dt, DtType a_dt_type=DtType::Full, bool skip_deposition=false)=0
amrex::Real m_current_injection_position
Current injection position.
Definition: WarpXParticleContainer.H:348
int species_id
Definition: WarpXParticleContainer.H:408
std::string m_qed_breit_wheeler_ele_product_name
Definition: WarpXParticleContainer.H:457
virtual void WriteHeader(std::ostream &os) const =0
amrex::Gpu::DeviceVector< amrex::Real > adk_exp_prefactor
Definition: WarpXParticleContainer.H:440
amrex::Gpu::DeviceVector< amrex::Real > adk_power
Definition: WarpXParticleContainer.H:438
PhysicalSpecies physical_species
Definition: WarpXParticleContainer.H:412
virtual void set_quantum_sync_engine_ptr(std::shared_ptr< QuantumSynchrotronEngine >)
Definition: WarpXParticleContainer.H:454
virtual void ContinuousFluxInjection(amrex::Real, amrex::Real)
Definition: WarpXParticleContainer.H:273
std::string getSpeciesTypeName() const
This function returns a string containing the name of the species type.
Definition: WarpXParticleContainer.H:389
int m_qed_quantum_sync_phot_product
Definition: WarpXParticleContainer.H:460
amrex::ParticleReal maxParticleVelocity(bool local=false)
Definition: WarpXParticleContainer.cpp:1211
virtual void InitIonizationModule()
Definition: WarpXParticleContainer.H:133
virtual void InitData()=0
amrex::Gpu::DeviceVector< amrex::Real > adk_prefactor
Definition: WarpXParticleContainer.H:439
void SetDoBackTransformedParticles(const bool do_back_transformed_particles)
Definition: WarpXParticleContainer.H:357
amrex::ParticleReal sumParticleCharge(bool local=false)
Definition: WarpXParticleContainer.cpp:1096
int ionization_initial_level
Definition: WarpXParticleContainer.H:436
amrex::Real self_fields_absolute_tolerance
Definition: WarpXParticleContainer.H:339
amrex::Vector< amrex::FArrayBox > local_jx
Definition: WarpXParticleContainer.H:465
virtual void UpdateAntennaPosition(const amrex::Real)
Update antenna position for continuous injection of lasers in a boosted frame. Empty function for con...
Definition: WarpXParticleContainer.H:268
amrex::ParticleReal mass
Definition: WarpXParticleContainer.H:411
bool doContinuousInjection() const
Definition: WarpXParticleContainer.H:270
amrex::Vector< amrex::FArrayBox > local_jy
Definition: WarpXParticleContainer.H:466
std::unique_ptr< amrex::MultiFab > GetChargeDensity(int lev, bool local=false)
Definition: WarpXParticleContainer.cpp:1072
static void ReadParameters()
Definition: WarpXParticleContainer.cpp:137
std::array< amrex::Gpu::DeviceVector< amrex::ParticleReal >, TmpIdx::nattribs > TmpParticleTile
Definition: WarpXParticleContainer.H:472
int ion_atomic_number
Definition: WarpXParticleContainer.H:435
bool m_do_back_transformed_particles
Definition: WarpXParticleContainer.H:446
amrex::Vector< amrex::FArrayBox > local_jz
Definition: WarpXParticleContainer.H:467
amrex::Gpu::DeviceVector< amrex::Real > ionization_energies
Definition: WarpXParticleContainer.H:437
amrex::ParticleReal getCharge() const
Definition: WarpXParticleContainer.H:362
TmpParticles tmp_particle_data
Definition: WarpXParticleContainer.H:477
void ApplyBoundaryConditions()
Apply particle BC.
Definition: WarpXParticleContainer.cpp:1344
virtual void set_breit_wheeler_engine_ptr(std::shared_ptr< BreitWheelerEngine >)
Definition: WarpXParticleContainer.H:452
virtual bool has_breit_wheeler() const
Definition: WarpXParticleContainer.H:371
virtual void PushP(int lev, amrex::Real dt, const amrex::MultiFab &Ex, const amrex::MultiFab &Ey, const amrex::MultiFab &Ez, const amrex::MultiFab &Bx, const amrex::MultiFab &By, const amrex::MultiFab &Bz)=0
std::string m_qed_breit_wheeler_pos_product_name
Definition: WarpXParticleContainer.H:459
std::string m_qed_quantum_sync_phot_product_name
Definition: WarpXParticleContainer.H:461
virtual void ContinuousInjection(const amrex::RealBox &)
Definition: WarpXParticleContainer.H:262
WarpXParticleContainer & operator=(const WarpXParticleContainer &)=delete
ParticleBoundaries m_boundary_conditions
Definition: WarpXParticleContainer.H:415
int m_qed_breit_wheeler_pos_product
Definition: WarpXParticleContainer.H:458
int ionization_product
Definition: WarpXParticleContainer.H:433
int m_qed_breit_wheeler_ele_product
Definition: WarpXParticleContainer.H:456
bool do_splitting
Definition: WarpXParticleContainer.H:335
std::array< amrex::ParticleReal, 3 > meanParticleVelocity(bool local=false)
Definition: WarpXParticleContainer.cpp:1125
WarpXParticleContainer & operator=(WarpXParticleContainer &&)=default
void DepositCurrent(amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > &J, amrex::Real dt, amrex::Real relative_time)
Deposit current density.
Definition: WarpXParticleContainer.cpp:624
bool initialize_self_fields
Definition: WarpXParticleContainer.H:337
int self_fields_verbosity
Definition: WarpXParticleContainer.H:341
void DepositCharge(amrex::Vector< std::unique_ptr< amrex::MultiFab > > &rho, bool local=false, bool reset=false, bool apply_boundary_and_scale_volume=false, bool interpolate_across_levels=true, int icomp=0)
Deposit charge density.
Definition: WarpXParticleContainer.cpp:969
virtual bool has_quantum_sync() const
Definition: WarpXParticleContainer.H:370
WarpXParticleContainer(amrex::AmrCore *amr_core, int ispecies)
Definition: WarpXParticleContainer.cpp:86
std::string physical_element
Definition: WarpXParticleContainer.H:441
int split_type
Definition: WarpXParticleContainer.H:351
void PushX(amrex::Real dt)
Definition: WarpXParticleContainer.cpp:1238
WarpXParticleContainer(WarpXParticleContainer &&)=default
void AddNParticles(int lev, long n, amrex::Vector< amrex::ParticleReal > const &x, amrex::Vector< amrex::ParticleReal > const &y, amrex::Vector< amrex::ParticleReal > const &z, amrex::Vector< amrex::ParticleReal > const &ux, amrex::Vector< amrex::ParticleReal > const &uy, amrex::Vector< amrex::ParticleReal > const &uz, int nattr_real, amrex::Vector< amrex::Vector< amrex::ParticleReal >> const &attr_real, int nattr_int, amrex::Vector< amrex::Vector< int >> const &attr_int, int uniqueparticles, amrex::Long id=-1)
Adds n particles to the simulation.
Definition: WarpXParticleContainer.cpp:158
TmpParticles getTmpParticleData() const noexcept
Definition: WarpXParticleContainer.H:475
virtual void DefaultInitializeRuntimeAttributes(amrex::ParticleTile< amrex::Particle< NStructReal, NStructInt >, NArrayReal, NArrayInt, amrex::PinnedArenaAllocator > &pinned_tile, int n_external_attr_real, int n_external_attr_int, const amrex::RandomEngine &engine)=0
Virtual method to initialize runtime attributes. Must be overriden by each derived class.
friend MultiParticleContainer
Definition: WarpXParticleContainer.H:111
amrex::Vector< amrex::ParticleReal > m_E_external_particle
Definition: WarpXParticleContainer.H:345
std::string ionization_product_name
Definition: WarpXParticleContainer.H:434
int getSpeciesId() const
Definition: WarpXParticleContainer.H:275
virtual void resample(const int, bool)
Virtual method to resample the species. Overriden by PhysicalParticleContainer only....
Definition: WarpXParticleContainer.H:397
int DoFieldIonization() const
Definition: WarpXParticleContainer.H:366
int do_not_deposit
Definition: WarpXParticleContainer.H:336
WarpXParticleContainer(const WarpXParticleContainer &)=delete
amrex::ParticleReal getMass() const
Definition: WarpXParticleContainer.H:364
int do_not_push
Definition: WarpXParticleContainer.H:423
int DoQED() const
Definition: WarpXParticleContainer.H:373
int do_resampling
Definition: WarpXParticleContainer.H:443
int do_not_gather
Definition: WarpXParticleContainer.H:424
static void BackwardCompatibility()
bool m_deposit_on_main_grid
instead of depositing (current, charge) on the finest patch level, deposit to the coarsest grid
Definition: WarpXParticleContainer.H:418
amrex::Vector< amrex::FArrayBox > local_rho
Definition: WarpXParticleContainer.H:464
std::pair< int, int > PairIndex
Definition: WarpXParticleContainer.H:470
int do_field_ionization
Definition: WarpXParticleContainer.H:432
amrex::ParticleReal charge
Definition: WarpXParticleContainer.H:410
typename SoA::RealVector RealVector
typename SoA::IntVector IntVector
static constexpr int NArrayInt
static constexpr int NArrayReal
typename SoA::RealVector RealVector
PODVector< T, ArenaAllocator< T > > DeviceVector
int n
Definition: run_libensemble_on_warpx.py:70
std::string get_name(const PhysicalSpecies &ps)
Returns the name associated to a PhysicalSpecies.
Definition: SpeciesPhysicalProperties.cpp:294
float dt
Definition: stencil.py:442
Definition: ParticleBoundaries.H:19
@ nattribs
Definition: WarpXParticleContainer_fwd.H:39