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  const std::array<RealVector, PIdx::nattribs>& GetAttribs () const {
61  return GetStructOfArrays().GetRealData();
62  }
63 
64  std::array<RealVector, PIdx::nattribs>& GetAttribs () {
65  return GetStructOfArrays().GetRealData();
66  }
67 
68  const RealVector& GetAttribs (int comp) const {
69  return GetStructOfArrays().GetRealData(comp);
70  }
71 
72  RealVector& GetAttribs (int comp) {
73  return GetStructOfArrays().GetRealData(comp);
74  }
75 
76  IntVector& GetiAttribs (int comp) {
77  return GetStructOfArrays().GetIntData(comp);
78  }
79 };
80 
103  : public NamedComponentParticleContainer<amrex::DefaultAllocator>
104 {
105 public:
107 
108  // amrex::StructOfArrays with DiagIdx::nattribs amrex::ParticleReal components
109  // and 0 int components for the particle data.
111  // DiagnosticParticles is a vector, with one element per MR level.
112  // DiagnosticParticles[lev] is typically a key-value pair where the key is
113  // a pair [grid_index, tile_index], and the value is the corresponding
114  // DiagnosticParticleData (see above) on this tile.
116 
117  WarpXParticleContainer (amrex::AmrCore* amr_core, int ispecies);
119 
120  // Move and copy operations
125 
126  virtual void InitData () = 0;
127 
128  virtual void InitIonizationModule () {}
129 
130  /*
131  * \brief Virtual function that returns a pointer to the plasma injector,
132  * for derived classes that define one (PhysicalParticleContainer).
133  */
134  virtual PlasmaInjector* GetPlasmaInjector () { return nullptr; }
135 
141  virtual void Evolve (int lev,
142  const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
143  const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz,
146  amrex::MultiFab* rho, amrex::MultiFab* crho,
147  const amrex::MultiFab* cEx, const amrex::MultiFab* cEy, const amrex::MultiFab* cEz,
148  const amrex::MultiFab* cBx, const amrex::MultiFab* cBy, const amrex::MultiFab* cBz,
149  amrex::Real t, amrex::Real dt, DtType a_dt_type=DtType::Full, bool skip_deposition=false) = 0;
150 
151  virtual void PostRestart () = 0;
152 
153  void AllocData ();
154 
162  amrex::PinnedArenaAllocator>& pinned_tile,
163  int n_external_attr_real,
164  int n_external_attr_int,
165  const amrex::RandomEngine& engine) = 0;
166 
172  void PushX ( amrex::Real dt);
173  void PushX (int lev, amrex::Real dt);
174 
178  virtual void PushP (int lev, amrex::Real dt,
179  const amrex::MultiFab& Ex,
180  const amrex::MultiFab& Ey,
181  const amrex::MultiFab& Ez,
182  const amrex::MultiFab& Bx,
183  const amrex::MultiFab& By,
184  const amrex::MultiFab& Bz) = 0;
185 
197  void DepositCurrent (amrex::Vector<std::array< std::unique_ptr<amrex::MultiFab>, 3 > >& J,
198  amrex::Real dt, amrex::Real relative_time);
199 
210  void DepositCharge (amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
211  bool local = false, bool reset = false,
212  bool apply_boundary_and_scale_volume = false,
213  bool interpolate_across_levels = true,
214  int icomp = 0);
215  void DepositCharge (std::unique_ptr<amrex::MultiFab>& rho, int lev,
216  bool local = false, bool reset = false,
217  bool apply_boundary_and_scale_volume = false,
218  int icomp = 0);
219 
220  std::unique_ptr<amrex::MultiFab> GetChargeDensity(int lev, bool local = false);
221 
222  virtual void DepositCharge (WarpXParIter& pti,
223  RealVector const & wp,
224  const int* ion_lev,
225  amrex::MultiFab* rho,
226  int icomp,
227  long offset,
228  long np_to_depose,
229  int thread_num,
230  int lev,
231  int depos_lev);
232 
233  virtual void DepositCurrent (WarpXParIter& pti,
234  RealVector const & wp,
235  RealVector const & uxp,
236  RealVector const & uyp,
237  RealVector const & uzp,
238  int const* ion_lev,
239  amrex::MultiFab* jx,
240  amrex::MultiFab* jy,
241  amrex::MultiFab* jz,
242  long offset,
243  long np_to_depose,
244  int thread_num,
245  int lev,
246  int depos_lev,
247  amrex::Real dt,
248  amrex::Real relative_time);
249 
250  // If particles start outside of the domain, ContinuousInjection
251  // makes sure that they are initialized when they enter the domain, and
252  // NOT before. Virtual function, overriden by derived classes.
253  // Current status:
254  // PhysicalParticleContainer: implemented.
255  // LaserParticleContainer: implemented.
256  // RigidInjectedParticleContainer: not implemented.
257  virtual void ContinuousInjection(const amrex::RealBox& /*injection_box*/) {}
258 
263  virtual void UpdateAntennaPosition(const amrex::Real /*dt*/) {}
264 
266 
267  // Inject a continuous flux of particles from a defined plane
268  virtual void ContinuousFluxInjection(amrex::Real /*t*/, amrex::Real /*dt*/) {}
269 
270  int getSpeciesId() const {return species_id;}
271 
276  amrex::ParticleReal sumParticleCharge(bool local = false);
277 
278  std::array<amrex::ParticleReal, 3> meanParticleVelocity(bool local = false);
279 
280  amrex::ParticleReal maxParticleVelocity(bool local = false);
281 
306  void AddNParticles (int lev, long n,
313  int nattr_real,
315  int nattr_int, amrex::Vector<amrex::Vector<int>> const & attr_int,
316  int uniqueparticles, amrex::Long id=-1);
317 
318  virtual void ReadHeader (std::istream& is) = 0;
319 
320  virtual void WriteHeader (std::ostream& os) const = 0;
321 
322  static void ReadParameters ();
323 
324  static void BackwardCompatibility ();
325 
328  void ApplyBoundaryConditions ();
329 
330  bool do_splitting = false;
331  int do_not_deposit = 0;
333  amrex::Real self_fields_required_precision = amrex::Real(1.e-11);
334  amrex::Real self_fields_absolute_tolerance = amrex::Real(0.0);
337 
340 
341  // split along diagonals (0) or axes (1)
342  int split_type = 0;
343 
348  void SetDoBackTransformedParticles(const bool do_back_transformed_particles) {
349  m_do_back_transformed_particles = do_back_transformed_particles;
350  }
351 
352  //amrex::Real getCharge () {return charge;}
353  amrex::ParticleReal getCharge () const {return charge;}
354  //amrex::Real getMass () {return mass;}
355  amrex::ParticleReal getMass () const {return mass;}
356 
357  int DoFieldIonization() const { return do_field_ionization; }
358 
359 #ifdef WARPX_QED
360  //Species for which QED effects are relevant should override these methods
361  virtual bool has_quantum_sync() const {return false;}
362  virtual bool has_breit_wheeler() const {return false;}
363 
364  int DoQED() const { return has_quantum_sync() || has_breit_wheeler(); }
365 #else
366  int DoQED() const { return false; }
367 #endif
368 
369  /* \brief This function tests if the current species
370  * is of a given PhysicalSpecies (specified as a template parameter).
371  * @tparam PhysSpec the PhysicalSpecies to test against
372  * @return the result of the test
373  */
374  template<PhysicalSpecies PhysSpec>
375  bool AmIA () const noexcept {return (physical_species == PhysSpec);}
376 
380  std::string getSpeciesTypeName () const {return species::get_name(physical_species);}
381 
388  virtual void resample (const int /*timestep*/, bool /*verbose*/) {}
389 
396  void defineAllParticleTiles () noexcept;
397 
398 protected:
400 
401  amrex::ParticleReal charge;
402  amrex::ParticleReal mass;
404 
405  // Controls boundaries for particles exiting the domain
407 
410 
413 
414  int do_not_push = 0;
415  int do_not_gather = 0;
416 
417  // Whether to allow particles outside of the simulation domain to be
418  // initialized when they enter the domain.
419  // This is currently required because continuous injection does not
420  // support all features allowed by direct injection.
422 
432  std::string physical_element;
433 
434  int do_resampling = 0;
435 
438 
439 #ifdef WARPX_QED
440  //Species can receive a shared pointer to a QED engine (species for
441  //which this is relevant should override these functions)
442  virtual void
443  set_breit_wheeler_engine_ptr(std::shared_ptr<BreitWheelerEngine>){}
444  virtual void
445  set_quantum_sync_engine_ptr(std::shared_ptr<QuantumSynchrotronEngine>){}
446 
453 
454 #endif
459 
460 public:
461  using PairIndex = std::pair<int, int>;
462  using TmpParticleTile = std::array<amrex::Gpu::DeviceVector<amrex::ParticleReal>,
465 
467 protected:
469 
470 private:
471  void particlePostLocate(ParticleType& p, const amrex::ParticleLocData& pld, int lev) override;
472 
473 };
474 
475 #endif
PhysicalSpecies
Definition: SpeciesPhysicalProperties.H:16
DtType
Definition: WarpXDtType.H:11
Definition: NamedComponentParticleContainer.H:49
Definition: PlasmaInjector.H:41
Definition: WarpXParticleContainer.H:52
WarpXParIter(ContainerType &pc, int level)
Definition: WarpXParticleContainer.cpp:73
std::array< RealVector, PIdx::nattribs > & GetAttribs()
Definition: WarpXParticleContainer.H:64
const RealVector & GetAttribs(int comp) const
Definition: WarpXParticleContainer.H:68
IntVector & GetiAttribs(int comp)
Definition: WarpXParticleContainer.H:76
RealVector & GetAttribs(int comp)
Definition: WarpXParticleContainer.H:72
const std::array< RealVector, PIdx::nattribs > & GetAttribs() const
Definition: WarpXParticleContainer.H:60
Definition: WarpXParticleContainer.H:104
amrex::Real self_fields_required_precision
Definition: WarpXParticleContainer.H:333
int do_continuous_injection
Definition: WarpXParticleContainer.H:421
bool AmIA() const noexcept
Definition: WarpXParticleContainer.H:375
bool m_gather_from_main_grid
instead of gathering fields from the finest patch level, gather from the coarsest
Definition: WarpXParticleContainer.H:412
void AllocData()
Definition: WarpXParticleContainer.cpp:137
void defineAllParticleTiles() noexcept
Definition: WarpXParticleContainer.cpp:1293
virtual void PostRestart()=0
void particlePostLocate(ParticleType &p, const amrex::ParticleLocData &pld, int lev) override
Definition: WarpXParticleContainer.cpp:1310
virtual void ReadHeader(std::istream &is)=0
int self_fields_max_iters
Definition: WarpXParticleContainer.H:335
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:339
int species_id
Definition: WarpXParticleContainer.H:399
std::string m_qed_breit_wheeler_ele_product_name
Definition: WarpXParticleContainer.H:448
virtual void WriteHeader(std::ostream &os) const =0
amrex::Gpu::DeviceVector< amrex::Real > adk_exp_prefactor
Definition: WarpXParticleContainer.H:431
amrex::Gpu::DeviceVector< amrex::Real > adk_power
Definition: WarpXParticleContainer.H:429
PhysicalSpecies physical_species
Definition: WarpXParticleContainer.H:403
virtual void set_quantum_sync_engine_ptr(std::shared_ptr< QuantumSynchrotronEngine >)
Definition: WarpXParticleContainer.H:445
virtual void ContinuousFluxInjection(amrex::Real, amrex::Real)
Definition: WarpXParticleContainer.H:268
std::string getSpeciesTypeName() const
This function returns a string containing the name of the species type.
Definition: WarpXParticleContainer.H:380
int m_qed_quantum_sync_phot_product
Definition: WarpXParticleContainer.H:451
amrex::ParticleReal maxParticleVelocity(bool local=false)
Definition: WarpXParticleContainer.cpp:1199
virtual void InitIonizationModule()
Definition: WarpXParticleContainer.H:128
virtual void InitData()=0
amrex::Gpu::DeviceVector< amrex::Real > adk_prefactor
Definition: WarpXParticleContainer.H:430
void SetDoBackTransformedParticles(const bool do_back_transformed_particles)
Definition: WarpXParticleContainer.H:348
amrex::ParticleReal sumParticleCharge(bool local=false)
Definition: WarpXParticleContainer.cpp:1084
int ionization_initial_level
Definition: WarpXParticleContainer.H:427
amrex::Real self_fields_absolute_tolerance
Definition: WarpXParticleContainer.H:334
amrex::Vector< amrex::FArrayBox > local_jx
Definition: WarpXParticleContainer.H:456
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:263
amrex::ParticleReal mass
Definition: WarpXParticleContainer.H:402
bool doContinuousInjection() const
Definition: WarpXParticleContainer.H:265
amrex::Vector< amrex::FArrayBox > local_jy
Definition: WarpXParticleContainer.H:457
std::unique_ptr< amrex::MultiFab > GetChargeDensity(int lev, bool local=false)
Definition: WarpXParticleContainer.cpp:1060
static void ReadParameters()
Definition: WarpXParticleContainer.cpp:125
std::array< amrex::Gpu::DeviceVector< amrex::ParticleReal >, TmpIdx::nattribs > TmpParticleTile
Definition: WarpXParticleContainer.H:463
int ion_atomic_number
Definition: WarpXParticleContainer.H:426
bool m_do_back_transformed_particles
Definition: WarpXParticleContainer.H:437
amrex::Vector< amrex::FArrayBox > local_jz
Definition: WarpXParticleContainer.H:458
amrex::Gpu::DeviceVector< amrex::Real > ionization_energies
Definition: WarpXParticleContainer.H:428
amrex::ParticleReal getCharge() const
Definition: WarpXParticleContainer.H:353
TmpParticles tmp_particle_data
Definition: WarpXParticleContainer.H:468
void ApplyBoundaryConditions()
Apply particle BC.
Definition: WarpXParticleContainer.cpp:1332
virtual void set_breit_wheeler_engine_ptr(std::shared_ptr< BreitWheelerEngine >)
Definition: WarpXParticleContainer.H:443
virtual bool has_breit_wheeler() const
Definition: WarpXParticleContainer.H:362
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:450
std::string m_qed_quantum_sync_phot_product_name
Definition: WarpXParticleContainer.H:452
virtual void ContinuousInjection(const amrex::RealBox &)
Definition: WarpXParticleContainer.H:257
WarpXParticleContainer & operator=(const WarpXParticleContainer &)=delete
ParticleBoundaries m_boundary_conditions
Definition: WarpXParticleContainer.H:406
int m_qed_breit_wheeler_pos_product
Definition: WarpXParticleContainer.H:449
int ionization_product
Definition: WarpXParticleContainer.H:424
int m_qed_breit_wheeler_ele_product
Definition: WarpXParticleContainer.H:447
bool do_splitting
Definition: WarpXParticleContainer.H:330
std::array< amrex::ParticleReal, 3 > meanParticleVelocity(bool local=false)
Definition: WarpXParticleContainer.cpp:1113
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:612
bool initialize_self_fields
Definition: WarpXParticleContainer.H:332
int self_fields_verbosity
Definition: WarpXParticleContainer.H:336
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:957
virtual bool has_quantum_sync() const
Definition: WarpXParticleContainer.H:361
WarpXParticleContainer(amrex::AmrCore *amr_core, int ispecies)
Definition: WarpXParticleContainer.cpp:85
std::string physical_element
Definition: WarpXParticleContainer.H:432
int split_type
Definition: WarpXParticleContainer.H:342
void PushX(amrex::Real dt)
Definition: WarpXParticleContainer.cpp:1226
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:146
TmpParticles getTmpParticleData() const noexcept
Definition: WarpXParticleContainer.H:466
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:106
virtual PlasmaInjector * GetPlasmaInjector()
Definition: WarpXParticleContainer.H:134
std::string ionization_product_name
Definition: WarpXParticleContainer.H:425
int getSpeciesId() const
Definition: WarpXParticleContainer.H:270
virtual void resample(const int, bool)
Virtual method to resample the species. Overriden by PhysicalParticleContainer only....
Definition: WarpXParticleContainer.H:388
int DoFieldIonization() const
Definition: WarpXParticleContainer.H:357
int do_not_deposit
Definition: WarpXParticleContainer.H:331
WarpXParticleContainer(const WarpXParticleContainer &)=delete
amrex::ParticleReal getMass() const
Definition: WarpXParticleContainer.H:355
int do_not_push
Definition: WarpXParticleContainer.H:414
~WarpXParticleContainer() override
Definition: WarpXParticleContainer.H:118
int DoQED() const
Definition: WarpXParticleContainer.H:364
int do_resampling
Definition: WarpXParticleContainer.H:434
int do_not_gather
Definition: WarpXParticleContainer.H:415
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:409
amrex::Vector< amrex::FArrayBox > local_rho
Definition: WarpXParticleContainer.H:455
std::pair< int, int > PairIndex
Definition: WarpXParticleContainer.H:461
int do_field_ionization
Definition: WarpXParticleContainer.H:423
amrex::ParticleReal charge
Definition: WarpXParticleContainer.H:401
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