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 using namespace amrex::literals;
50 
52  : public amrex::ParIter<0,0,PIdx::nattribs>
53 {
54 public:
56 
57  WarpXParIter (ContainerType& pc, int level);
58 
59  WarpXParIter (ContainerType& pc, int level, amrex::MFItInfo& info);
60 
61  const std::array<RealVector, PIdx::nattribs>& GetAttribs () const {
62  return GetStructOfArrays().GetRealData();
63  }
64 
65  std::array<RealVector, PIdx::nattribs>& GetAttribs () {
66  return GetStructOfArrays().GetRealData();
67  }
68 
69  const RealVector& GetAttribs (int comp) const {
70  return GetStructOfArrays().GetRealData(comp);
71  }
72 
73  RealVector& GetAttribs (int comp) {
74  return GetStructOfArrays().GetRealData(comp);
75  }
76 
77  IntVector& GetiAttribs (int comp) {
78  return GetStructOfArrays().GetIntData(comp);
79  }
80 };
81 
104  : public NamedComponentParticleContainer<amrex::DefaultAllocator>
105 {
106 public:
108 
109  // amrex::StructOfArrays with DiagIdx::nattribs amrex::ParticleReal components
110  // and 0 int components for the particle data.
112  // DiagnosticParticles is a vector, with one element per MR level.
113  // DiagnosticParticles[lev] is typically a key-value pair where the key is
114  // a pair [grid_index, tile_index], and the value is the corresponding
115  // DiagnosticParticleData (see above) on this tile.
117 
118  WarpXParticleContainer (amrex::AmrCore* amr_core, int ispecies);
120 
121  virtual void InitData () = 0;
122 
123  virtual void InitIonizationModule () {}
124 
125  /*
126  * \brief Virtual function that returns a pointer to the plasma injector,
127  * for derived classes that define one (PhysicalParticleContainer).
128  */
129  virtual PlasmaInjector* GetPlasmaInjector () { return nullptr; }
130 
136  virtual void Evolve (int lev,
137  const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
138  const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz,
141  amrex::MultiFab* rho, amrex::MultiFab* crho,
142  const amrex::MultiFab* cEx, const amrex::MultiFab* cEy, const amrex::MultiFab* cEz,
143  const amrex::MultiFab* cBx, const amrex::MultiFab* cBy, const amrex::MultiFab* cBz,
144  amrex::Real t, amrex::Real dt, DtType a_dt_type=DtType::Full, bool skip_deposition=false) = 0;
145 
146  virtual void PostRestart () = 0;
147 
148  void AllocData ();
149 
156  NArrayReal, NArrayInt,
157  amrex::PinnedArenaAllocator>& pinned_tile,
158  const int n_external_attr_real,
159  const int n_external_attr_int,
160  const amrex::RandomEngine& engine) = 0;
161 
167  void PushX ( amrex::Real dt);
168  void PushX (int lev, amrex::Real dt);
169 
173  virtual void PushP (int lev, amrex::Real dt,
174  const amrex::MultiFab& Ex,
175  const amrex::MultiFab& Ey,
176  const amrex::MultiFab& Ez,
177  const amrex::MultiFab& Bx,
178  const amrex::MultiFab& By,
179  const amrex::MultiFab& Bz) = 0;
180 
192  void DepositCurrent (amrex::Vector<std::array< std::unique_ptr<amrex::MultiFab>, 3 > >& J,
193  const amrex::Real dt, const amrex::Real relative_time);
194 
205  void DepositCharge (amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
206  const bool local = false, const bool reset = false,
207  const bool apply_boundary_and_scale_volume = false,
208  const bool interpolate_across_levels = true,
209  const int icomp = 0);
210  void DepositCharge (std::unique_ptr<amrex::MultiFab>& rho, const int lev,
211  const bool local = false, const bool reset = false,
212  const bool apply_boundary_and_scale_volume = false,
213  const int icomp = 0);
214 
215  std::unique_ptr<amrex::MultiFab> GetChargeDensity(int lev, bool local = false);
216 
217  virtual void DepositCharge (WarpXParIter& pti,
218  RealVector const & wp,
219  const int * const ion_lev,
220  amrex::MultiFab* rho,
221  const int icomp,
222  const long offset,
223  const long np_to_depose,
224  const int thread_num,
225  const int lev,
226  const int depos_lev);
227 
228  virtual void DepositCurrent (WarpXParIter& pti,
229  RealVector const & wp,
230  RealVector const & uxp,
231  RealVector const & uyp,
232  RealVector const & uzp,
233  int const * const ion_lev,
234  amrex::MultiFab* const jx,
235  amrex::MultiFab* const jy,
236  amrex::MultiFab* const jz,
237  long const offset,
238  long const np_to_depose,
239  int const thread_num,
240  int const lev,
241  int const depos_lev,
242  amrex::Real const dt,
243  amrex::Real const relative_time);
244 
245  // If particles start outside of the domain, ContinuousInjection
246  // makes sure that they are initialized when they enter the domain, and
247  // NOT before. Virtual function, overriden by derived classes.
248  // Current status:
249  // PhysicalParticleContainer: implemented.
250  // LaserParticleContainer: implemented.
251  // RigidInjectedParticleContainer: not implemented.
252  virtual void ContinuousInjection(const amrex::RealBox& /*injection_box*/) {}
253  // Update optional sub-class-specific injection location.
254  virtual void UpdateContinuousInjectionPosition(amrex::Real /*dt*/) {}
255  bool doContinuousInjection() const {return do_continuous_injection;}
256 
257  // Inject a continuous flux of particles from a defined plane
258  virtual void ContinuousFluxInjection(amrex::Real /*t*/, amrex::Real /*dt*/) {}
259 
260  int getSpeciesId() const {return species_id;}
261 
266  amrex::ParticleReal sumParticleCharge(bool local = false);
267 
268  std::array<amrex::ParticleReal, 3> meanParticleVelocity(bool local = false);
269 
270  amrex::ParticleReal maxParticleVelocity(bool local = false);
271 
296  void AddNParticles (int lev,
297  int n, const amrex::ParticleReal* x, const amrex::ParticleReal* y, const amrex::ParticleReal* z,
298  const amrex::ParticleReal* ux, const amrex::ParticleReal* uy, const amrex::ParticleReal* uz,
299  const int nattr_real, const amrex::ParticleReal* attr_real,
300  const int nattr_int, const int* attr_int,
301  int uniqueparticles, amrex::Long id=-1);
302 
303  virtual void ReadHeader (std::istream& is) = 0;
304 
305  virtual void WriteHeader (std::ostream& os) const = 0;
306 
307  static void ReadParameters ();
308 
309  static void BackwardCompatibility ();
310 
313  void ApplyBoundaryConditions ();
314 
315  bool do_splitting = false;
316  int do_not_deposit = 0;
317  bool initialize_self_fields = false;
318  amrex::Real self_fields_required_precision = amrex::Real(1.e-11);
319  amrex::Real self_fields_absolute_tolerance = amrex::Real(0.0);
320  int self_fields_max_iters = 200;
321  int self_fields_verbosity = 2;
322 
325 
326  // split along diagonals (0) or axes (1)
327  int split_type = 0;
328 
333  void SetDoBackTransformedParticles(const bool do_back_transformed_particles) {
334  m_do_back_transformed_particles = do_back_transformed_particles;
335  }
336 
337  //amrex::Real getCharge () {return charge;}
338  amrex::ParticleReal getCharge () const {return charge;}
339  //amrex::Real getMass () {return mass;}
340  amrex::ParticleReal getMass () const {return mass;}
341 
342  int DoFieldIonization() const { return do_field_ionization; }
343 
344 #ifdef WARPX_QED
345  //Species for which QED effects are relevant should override these methods
346  virtual bool has_quantum_sync() const {return false;}
347  virtual bool has_breit_wheeler() const {return false;}
348 
349  int DoQED() const { return has_quantum_sync() || has_breit_wheeler(); }
350 #else
351  int DoQED() const { return false; }
352 #endif
353 
354  /* \brief This function tests if the current species
355  * is of a given PhysicalSpecies (specified as a template parameter).
356  * @tparam PhysSpec the PhysicalSpecies to test against
357  * @return the result of the test
358  */
359  template<PhysicalSpecies PhysSpec>
360  bool AmIA () const noexcept {return (physical_species == PhysSpec);}
361 
365  std::string getSpeciesTypeName () const {return species::get_name(physical_species);}
366 
373  virtual void resample (const int /*timestep*/) {}
374 
381  void defineAllParticleTiles () noexcept;
382 
383 protected:
384  int species_id;
385 
386  amrex::ParticleReal charge;
387  amrex::ParticleReal mass;
388  PhysicalSpecies physical_species;
389 
390  // Controls boundaries for particles exiting the domain
391  ParticleBoundaries m_boundary_conditions;
392 
394  bool m_deposit_on_main_grid = false;
395 
397  bool m_gather_from_main_grid = false;
398 
399  int do_not_push = 0;
400  int do_not_gather = 0;
401 
402  // Whether to allow particles outside of the simulation domain to be
403  // initialized when they enter the domain.
404  // This is currently required because continuous injection does not
405  // support all features allowed by direct injection.
406  int do_continuous_injection = 0;
407 
408  int do_field_ionization = 0;
409  int ionization_product;
410  std::string ionization_product_name;
411  int ion_atomic_number;
412  int ionization_initial_level = 0;
413  amrex::Gpu::DeviceVector<amrex::Real> ionization_energies;
414  amrex::Gpu::DeviceVector<amrex::Real> adk_power;
415  amrex::Gpu::DeviceVector<amrex::Real> adk_prefactor;
416  amrex::Gpu::DeviceVector<amrex::Real> adk_exp_prefactor;
417  std::string physical_element;
418 
419  int do_resampling = 0;
420 
422  bool m_do_back_transformed_particles = false;
423 
424 #ifdef WARPX_QED
425  //Species can receive a shared pointer to a QED engine (species for
426  //which this is relevant should override these functions)
427  virtual void
428  set_breit_wheeler_engine_ptr(std::shared_ptr<BreitWheelerEngine>){}
429  virtual void
430  set_quantum_sync_engine_ptr(std::shared_ptr<QuantumSynchrotronEngine>){}
431 
438 
439 #endif
444 
445 public:
446  using PairIndex = std::pair<int, int>;
447  using TmpParticleTile = std::array<amrex::Gpu::DeviceVector<amrex::ParticleReal>,
450 
451  TmpParticles getTmpParticleData () const noexcept {return tmp_particle_data;}
452 protected:
454 
455 private:
456  virtual void particlePostLocate(ParticleType& p, const amrex::ParticleLocData& pld,
457  const int lev) override;
458 
459 };
460 
461 #endif
PhysicalSpecies
Definition: SpeciesPhysicalProperties.H:16
DtType
Definition: WarpXDtType.H:11
Definition: NamedComponentParticleContainer.H:48
Definition: PlasmaInjector.H:41
Definition: WarpXParticleContainer.H:53
std::array< RealVector, PIdx::nattribs > & GetAttribs()
Definition: WarpXParticleContainer.H:65
const RealVector & GetAttribs(int comp) const
Definition: WarpXParticleContainer.H:69
IntVector & GetiAttribs(int comp)
Definition: WarpXParticleContainer.H:77
RealVector & GetAttribs(int comp)
Definition: WarpXParticleContainer.H:73
const std::array< RealVector, PIdx::nattribs > & GetAttribs() const
Definition: WarpXParticleContainer.H:61
Definition: WarpXParticleContainer.H:105
bool AmIA() const noexcept
Definition: WarpXParticleContainer.H:360
virtual void PostRestart()=0
virtual void ReadHeader(std::istream &is)=0
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:324
std::string m_qed_breit_wheeler_ele_product_name
Definition: WarpXParticleContainer.H:433
virtual void WriteHeader(std::ostream &os) const =0
virtual void set_quantum_sync_engine_ptr(std::shared_ptr< QuantumSynchrotronEngine >)
Definition: WarpXParticleContainer.H:430
virtual void ContinuousFluxInjection(amrex::Real, amrex::Real)
Definition: WarpXParticleContainer.H:258
std::string getSpeciesTypeName() const
This function returns a string containing the name of the species type.
Definition: WarpXParticleContainer.H:365
int m_qed_quantum_sync_phot_product
Definition: WarpXParticleContainer.H:436
virtual void InitIonizationModule()
Definition: WarpXParticleContainer.H:123
virtual void InitData()=0
void SetDoBackTransformedParticles(const bool do_back_transformed_particles)
Definition: WarpXParticleContainer.H:333
amrex::Vector< amrex::FArrayBox > local_jx
Definition: WarpXParticleContainer.H:441
bool doContinuousInjection() const
Definition: WarpXParticleContainer.H:255
amrex::Vector< amrex::FArrayBox > local_jy
Definition: WarpXParticleContainer.H:442
std::array< amrex::Gpu::DeviceVector< amrex::ParticleReal >, TmpIdx::nattribs > TmpParticleTile
Definition: WarpXParticleContainer.H:448
virtual ~WarpXParticleContainer()
Definition: WarpXParticleContainer.H:119
amrex::Vector< amrex::FArrayBox > local_jz
Definition: WarpXParticleContainer.H:443
virtual void resample(const int)
Virtual method to resample the species. Overriden by PhysicalParticleContainer only....
Definition: WarpXParticleContainer.H:373
amrex::ParticleReal getCharge() const
Definition: WarpXParticleContainer.H:338
TmpParticles tmp_particle_data
Definition: WarpXParticleContainer.H:453
virtual void set_breit_wheeler_engine_ptr(std::shared_ptr< BreitWheelerEngine >)
Definition: WarpXParticleContainer.H:428
virtual bool has_breit_wheeler() const
Definition: WarpXParticleContainer.H:347
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
virtual void UpdateContinuousInjectionPosition(amrex::Real)
Definition: WarpXParticleContainer.H:254
std::string m_qed_breit_wheeler_pos_product_name
Definition: WarpXParticleContainer.H:435
std::string m_qed_quantum_sync_phot_product_name
Definition: WarpXParticleContainer.H:437
virtual void ContinuousInjection(const amrex::RealBox &)
Definition: WarpXParticleContainer.H:252
int m_qed_breit_wheeler_pos_product
Definition: WarpXParticleContainer.H:434
int m_qed_breit_wheeler_ele_product
Definition: WarpXParticleContainer.H:432
virtual void DefaultInitializeRuntimeAttributes(amrex::ParticleTile< amrex::Particle< NStructReal, NStructInt >, NArrayReal, NArrayInt, amrex::PinnedArenaAllocator > &pinned_tile, const int n_external_attr_real, const int n_external_attr_int, const amrex::RandomEngine &engine)=0
Virtual method to initialize runtime attributes. Must be overriden by each derived class.
virtual bool has_quantum_sync() const
Definition: WarpXParticleContainer.H:346
TmpParticles getTmpParticleData() const noexcept
Definition: WarpXParticleContainer.H:451
friend MultiParticleContainer
Definition: WarpXParticleContainer.H:107
virtual PlasmaInjector * GetPlasmaInjector()
Definition: WarpXParticleContainer.H:129
int getSpeciesId() const
Definition: WarpXParticleContainer.H:260
int DoFieldIonization() const
Definition: WarpXParticleContainer.H:342
amrex::ParticleReal getMass() const
Definition: WarpXParticleContainer.H:340
int DoQED() const
Definition: WarpXParticleContainer.H:349
static void BackwardCompatibility()
amrex::Vector< amrex::FArrayBox > local_rho
Definition: WarpXParticleContainer.H:440
std::pair< int, int > PairIndex
Definition: WarpXParticleContainer.H:446
typename SoA::RealVector RealVector
typename SoA::IntVector IntVector
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:440
Definition: ParticleBoundaries.H:19
@ nattribs
Definition: WarpXParticleContainer_fwd.H:39