ARTEMIS
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"
18 
19 #ifdef WARPX_QED
22 #endif
25 
26 #include <AMReX_Array.H>
27 #include <AMReX_FArrayBox.H>
28 #include <AMReX_GpuAllocators.H>
29 #include <AMReX_GpuContainers.H>
30 #include <AMReX_INT.H>
31 #include <AMReX_ParIter.H>
32 #include <AMReX_Particles.H>
33 #include <AMReX_Random.H>
34 #include <AMReX_REAL.H>
35 #include <AMReX_StructOfArrays.H>
36 #include <AMReX_Vector.H>
37 
38 #include <AMReX_BaseFwd.H>
39 #include <AMReX_AmrCoreFwd.H>
40 
41 #include <array>
42 #include <iosfwd>
43 #include <map>
44 #include <memory>
45 #include <string>
46 #include <utility>
47 
48 using namespace amrex::literals;
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  virtual void InitData () = 0;
121 
122  virtual void InitIonizationModule () {}
123 
129  virtual void Evolve (int lev,
130  const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
131  const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz,
134  amrex::MultiFab* rho, amrex::MultiFab* crho,
135  const amrex::MultiFab* cEx, const amrex::MultiFab* cEy, const amrex::MultiFab* cEz,
136  const amrex::MultiFab* cBx, const amrex::MultiFab* cBy, const amrex::MultiFab* cBz,
137  amrex::Real t, amrex::Real dt, DtType a_dt_type=DtType::Full, bool skip_deposition=false) = 0;
138 
139  virtual void PostRestart () = 0;
140 
141  void AllocData ();
142 
149  NArrayReal, NArrayInt,
150  amrex::PinnedArenaAllocator>& pinned_tile,
151  const int n_external_attr_real,
152  const int n_external_attr_int,
153  const amrex::RandomEngine& engine) = 0;
154 
160  void PushX ( amrex::Real dt);
161  void PushX (int lev, amrex::Real dt);
162 
166  virtual void PushP (int lev, amrex::Real dt,
167  const amrex::MultiFab& Ex,
168  const amrex::MultiFab& Ey,
169  const amrex::MultiFab& Ez,
170  const amrex::MultiFab& Bx,
171  const amrex::MultiFab& By,
172  const amrex::MultiFab& Bz) = 0;
173 
185  void DepositCurrent (amrex::Vector<std::array< std::unique_ptr<amrex::MultiFab>, 3 > >& J,
186  const amrex::Real dt, const amrex::Real relative_time);
187 
198  void DepositCharge (amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
199  const bool local = false, const bool reset = false,
200  const bool apply_boundary_and_scale_volume = false,
201  const bool interpolate_across_levels = true,
202  const int icomp = 0);
203 
204  std::unique_ptr<amrex::MultiFab> GetChargeDensity(int lev, bool local = false);
205 
206  virtual void DepositCharge (WarpXParIter& pti,
207  RealVector const & wp,
208  const int * const ion_lev,
209  amrex::MultiFab* rho,
210  const int icomp,
211  const long offset,
212  const long np_to_depose,
213  const int thread_num,
214  const int lev,
215  const int depos_lev);
216 
217  virtual void DepositCurrent (WarpXParIter& pti,
218  RealVector const & wp,
219  RealVector const & uxp,
220  RealVector const & uyp,
221  RealVector const & uzp,
222  int const * const ion_lev,
223  amrex::MultiFab* const jx,
224  amrex::MultiFab* const jy,
225  amrex::MultiFab* const jz,
226  long const offset,
227  long const np_to_depose,
228  int const thread_num,
229  int const lev,
230  int const depos_lev,
231  amrex::Real const dt,
232  amrex::Real const relative_time);
233 
234  // If particles start outside of the domain, ContinuousInjection
235  // makes sure that they are initialized when they enter the domain, and
236  // NOT before. Virtual function, overriden by derived classes.
237  // Current status:
238  // PhysicalParticleContainer: implemented.
239  // LaserParticleContainer: implemented.
240  // RigidInjectedParticleContainer: not implemented.
241  virtual void ContinuousInjection(const amrex::RealBox& /*injection_box*/) {}
242  // Update optional sub-class-specific injection location.
243  virtual void UpdateContinuousInjectionPosition(amrex::Real /*dt*/) {}
244 
245  // Inject a continuous flux of particles from a defined plane
246  virtual void ContinuousFluxInjection(amrex::Real /*t*/, amrex::Real /*dt*/) {}
247 
252  amrex::ParticleReal sumParticleCharge(bool local = false);
253 
254  std::array<amrex::ParticleReal, 3> meanParticleVelocity(bool local = false);
255 
256  amrex::ParticleReal maxParticleVelocity(bool local = false);
257 
282  void AddNParticles (int lev,
283  int n, const amrex::ParticleReal* x, const amrex::ParticleReal* y, const amrex::ParticleReal* z,
284  const amrex::ParticleReal* ux, const amrex::ParticleReal* uy, const amrex::ParticleReal* uz,
285  const int nattr_real, const amrex::ParticleReal* attr_real,
286  const int nattr_int, const int* attr_int,
287  int uniqueparticles, amrex::Long id=-1);
288 
289  virtual void ReadHeader (std::istream& is) = 0;
290 
291  virtual void WriteHeader (std::ostream& os) const = 0;
292 
293  static void ReadParameters ();
294 
295  static void BackwardCompatibility ();
296 
299  void ApplyBoundaryConditions ();
300 
301  bool do_splitting = false;
302  int do_not_deposit = 0;
303  bool initialize_self_fields = false;
304  amrex::Real self_fields_required_precision = amrex::Real(1.e-11);
305  amrex::Real self_fields_absolute_tolerance = amrex::Real(0.0);
306  int self_fields_max_iters = 200;
307  int self_fields_verbosity = 2;
308 
309  // split along diagonals (0) or axes (1)
310  int split_type = 0;
311 
316  void SetDoBackTransformedParticles(const bool do_back_transformed_particles) {
317  m_do_back_transformed_particles = do_back_transformed_particles;
318  }
319 
320  //amrex::Real getCharge () {return charge;}
321  amrex::ParticleReal getCharge () const {return charge;}
322  //amrex::Real getMass () {return mass;}
323  amrex::ParticleReal getMass () const {return mass;}
324 
325  int DoFieldIonization() const { return do_field_ionization; }
326 
327 #ifdef WARPX_QED
328  //Species for which QED effects are relevant should override these methods
329  virtual bool has_quantum_sync() const {return false;}
330  virtual bool has_breit_wheeler() const {return false;}
331 
332  int DoQED() const { return has_quantum_sync() || has_breit_wheeler(); }
333 #else
334  int DoQED() const { return false; }
335 #endif
336 
337  /* \brief This function tests if the current species
338  * is of a given PhysicalSpecies (specified as a template parameter).
339  * @tparam PhysSpec the PhysicalSpecies to test against
340  * @return the result of the test
341  */
342  template<PhysicalSpecies PhysSpec>
343  bool AmIA () const noexcept {return (physical_species == PhysSpec);}
344 
348  std::string getSpeciesTypeName () const {return species::get_name(physical_species);}
349 
356  virtual void resample (const int /*timestep*/) {}
357 
364  void defineAllParticleTiles () noexcept;
365 
366 protected:
367  int species_id;
368 
369  amrex::ParticleReal charge;
370  amrex::ParticleReal mass;
371  PhysicalSpecies physical_species;
372 
373  // Controls boundaries for particles exiting the domain
374  ParticleBoundaries m_boundary_conditions;
375 
377  bool m_deposit_on_main_grid = false;
378 
380  bool m_gather_from_main_grid = false;
381 
382  int do_not_push = 0;
383  int do_not_gather = 0;
384 
385  // Whether to allow particles outside of the simulation domain to be
386  // initialized when they enter the domain.
387  // This is currently required because continuous injection does not
388  // support all features allowed by direct injection.
389  int do_continuous_injection = 0;
390 
391  int do_field_ionization = 0;
392  int ionization_product;
393  std::string ionization_product_name;
394  int ion_atomic_number;
395  int ionization_initial_level = 0;
396  amrex::Gpu::DeviceVector<amrex::Real> ionization_energies;
397  amrex::Gpu::DeviceVector<amrex::Real> adk_power;
398  amrex::Gpu::DeviceVector<amrex::Real> adk_prefactor;
399  amrex::Gpu::DeviceVector<amrex::Real> adk_exp_prefactor;
400  std::string physical_element;
401 
402  int do_resampling = 0;
403 
405  bool m_do_back_transformed_particles = false;
406 
407 #ifdef WARPX_QED
408  //Species can receive a shared pointer to a QED engine (species for
409  //which this is relevant should override these functions)
410  virtual void
411  set_breit_wheeler_engine_ptr(std::shared_ptr<BreitWheelerEngine>){}
412  virtual void
413  set_quantum_sync_engine_ptr(std::shared_ptr<QuantumSynchrotronEngine>){}
414 
421 
422 #endif
427 
428 public:
429  using PairIndex = std::pair<int, int>;
430  using TmpParticleTile = std::array<amrex::Gpu::DeviceVector<amrex::ParticleReal>,
433 
434  TmpParticles getTmpParticleData () const noexcept {return tmp_particle_data;}
435 protected:
437 
438 private:
439  virtual void particlePostLocate(ParticleType& p, const amrex::ParticleLocData& pld,
440  const int lev) override;
441 
442 };
443 
444 #endif
PhysicalSpecies
Definition: SpeciesPhysicalProperties.H:16
DtType
Definition: WarpXDtType.H:11
Definition: NamedComponentParticleContainer.H:47
Definition: WarpXParticleContainer.H:52
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
bool AmIA() const noexcept
Definition: WarpXParticleContainer.H:343
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
std::string m_qed_breit_wheeler_ele_product_name
Definition: WarpXParticleContainer.H:416
virtual void WriteHeader(std::ostream &os) const =0
virtual void set_quantum_sync_engine_ptr(std::shared_ptr< QuantumSynchrotronEngine >)
Definition: WarpXParticleContainer.H:413
virtual void ContinuousFluxInjection(amrex::Real, amrex::Real)
Definition: WarpXParticleContainer.H:246
std::string getSpeciesTypeName() const
This function returns a string containing the name of the species type.
Definition: WarpXParticleContainer.H:348
int m_qed_quantum_sync_phot_product
Definition: WarpXParticleContainer.H:419
virtual void InitIonizationModule()
Definition: WarpXParticleContainer.H:122
virtual void InitData()=0
void SetDoBackTransformedParticles(const bool do_back_transformed_particles)
Definition: WarpXParticleContainer.H:316
amrex::Vector< amrex::FArrayBox > local_jx
Definition: WarpXParticleContainer.H:424
amrex::Vector< amrex::FArrayBox > local_jy
Definition: WarpXParticleContainer.H:425
std::array< amrex::Gpu::DeviceVector< amrex::ParticleReal >, TmpIdx::nattribs > TmpParticleTile
Definition: WarpXParticleContainer.H:431
virtual ~WarpXParticleContainer()
Definition: WarpXParticleContainer.H:118
amrex::Vector< amrex::FArrayBox > local_jz
Definition: WarpXParticleContainer.H:426
virtual void resample(const int)
Virtual method to resample the species. Overriden by PhysicalParticleContainer only....
Definition: WarpXParticleContainer.H:356
amrex::ParticleReal getCharge() const
Definition: WarpXParticleContainer.H:321
TmpParticles tmp_particle_data
Definition: WarpXParticleContainer.H:436
virtual void set_breit_wheeler_engine_ptr(std::shared_ptr< BreitWheelerEngine >)
Definition: WarpXParticleContainer.H:411
virtual bool has_breit_wheeler() const
Definition: WarpXParticleContainer.H:330
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:243
std::string m_qed_breit_wheeler_pos_product_name
Definition: WarpXParticleContainer.H:418
std::string m_qed_quantum_sync_phot_product_name
Definition: WarpXParticleContainer.H:420
virtual void ContinuousInjection(const amrex::RealBox &)
Definition: WarpXParticleContainer.H:241
int m_qed_breit_wheeler_pos_product
Definition: WarpXParticleContainer.H:417
int m_qed_breit_wheeler_ele_product
Definition: WarpXParticleContainer.H:415
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:329
TmpParticles getTmpParticleData() const noexcept
Definition: WarpXParticleContainer.H:434
friend MultiParticleContainer
Definition: WarpXParticleContainer.H:106
int DoFieldIonization() const
Definition: WarpXParticleContainer.H:325
amrex::ParticleReal getMass() const
Definition: WarpXParticleContainer.H:323
int DoQED() const
Definition: WarpXParticleContainer.H:332
static void BackwardCompatibility()
amrex::Vector< amrex::FArrayBox > local_rho
Definition: WarpXParticleContainer.H:423
std::pair< int, int > PairIndex
Definition: WarpXParticleContainer.H:429
typename SoA::RealVector RealVector
typename SoA::IntVector IntVector
typename SoA::RealVector RealVector
def y
Definition: Excitation_Flag_Generator.py:76
def z
Definition: Excitation_Flag_Generator.py:77
int n
Definition: run_libensemble_on_warpx.py:67
std::string get_name(const PhysicalSpecies &ps)
Returns the name associated to a PhysicalSpecies.
Definition: SpeciesPhysicalProperties.cpp:294
int dt
Definition: stencil.py:440
Definition: ParticleBoundaries.H:19
@ nattribs
Definition: WarpXParticleContainer_fwd.H:39