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"
18 
19 #ifdef WARPX_QED
22 #endif
24 
25 #include <AMReX_Array.H>
26 #include <AMReX_FArrayBox.H>
27 #include <AMReX_GpuContainers.H>
28 #include <AMReX_INT.H>
29 #include <AMReX_ParIter.H>
30 #include <AMReX_Particles.H>
31 #include <AMReX_REAL.H>
32 #include <AMReX_StructOfArrays.H>
33 #include <AMReX_Vector.H>
34 
35 #include <AMReX_BaseFwd.H>
36 #include <AMReX_AmrCoreFwd.H>
37 
38 #include <array>
39 #include <iosfwd>
40 #include <map>
41 #include <memory>
42 #include <string>
43 #include <utility>
44 
46 {
47  const std::map<std::string, int> to_index = {
48  {"w", PIdx::w },
49  {"ux", PIdx::ux },
50  {"uy", PIdx::uy },
51  {"uz", PIdx::uz },
52 #ifdef WARPX_DIM_RZ
53  {"theta", PIdx::theta}
54 #endif
55  };
56 }
57 
59  : public amrex::ParIter<0,0,PIdx::nattribs>
60 {
61 public:
62  using amrex::ParIter<0,0,PIdx::nattribs>::ParIter;
63 
64  WarpXParIter (ContainerType& pc, int level);
65 
66  WarpXParIter (ContainerType& pc, int level, amrex::MFItInfo& info);
67 
68  const std::array<RealVector, PIdx::nattribs>& GetAttribs () const {
69  return GetStructOfArrays().GetRealData();
70  }
71 
72  std::array<RealVector, PIdx::nattribs>& GetAttribs () {
73  return GetStructOfArrays().GetRealData();
74  }
75 
76  const RealVector& GetAttribs (int comp) const {
77  return GetStructOfArrays().GetRealData(comp);
78  }
79 
80  RealVector& GetAttribs (int comp) {
81  return GetStructOfArrays().GetRealData(comp);
82  }
83 
84  IntVector& GetiAttribs (int comp) {
85  return GetStructOfArrays().GetIntData(comp);
86  }
87 };
88 
111  : public amrex::ParticleContainer<0,0,PIdx::nattribs>
112 {
113 public:
115 
116  // amrex::StructOfArrays with DiagIdx::nattribs amrex::ParticleReal components
117  // and 0 int components for the particle data.
118  using DiagnosticParticleData = amrex::StructOfArrays<DiagIdx::nattribs, 0>;
119  // DiagnosticParticles is a vector, with one element per MR level.
120  // DiagnosticParticles[lev] is typically a key-value pair where the key is
121  // a pair [grid_index, tile_index], and the value is the corresponding
122  // DiagnosticParticleData (see above) on this tile.
123  using DiagnosticParticles = amrex::Vector<std::map<std::pair<int, int>, DiagnosticParticleData> >;
124 
125  WarpXParticleContainer (amrex::AmrCore* amr_core, int ispecies);
127 
128  virtual void InitData () = 0;
129 
135  virtual void Evolve (int lev,
136  const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
137  const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz,
138  amrex::MultiFab& jx, amrex::MultiFab& jy, amrex::MultiFab& jz,
139  amrex::MultiFab* cjx, amrex::MultiFab* cjy, amrex::MultiFab* cjz,
140  amrex::MultiFab* rho, amrex::MultiFab* crho,
141  const amrex::MultiFab* cEx, const amrex::MultiFab* cEy, const amrex::MultiFab* cEz,
142  const amrex::MultiFab* cBx, const amrex::MultiFab* cBy, const amrex::MultiFab* cBz,
143  amrex::Real t, amrex::Real dt, DtType a_dt_type=DtType::Full, bool skip_deposition=false) = 0;
144 
145  virtual void PostRestart () = 0;
146 
147  virtual void GetParticleSlice(const int /*direction*/, const amrex::Real /*z_old*/,
148  const amrex::Real /*z_new*/, const amrex::Real /*t_boost*/,
149  const amrex::Real /*t_lab*/, const amrex::Real /*dt*/,
150  DiagnosticParticles& /*diagnostic_particles*/) {}
151 
152  void AllocData ();
153 
159  void PushX ( amrex::Real dt);
160  void PushX (int lev, amrex::Real dt);
161 
165  virtual void PushP (int lev, amrex::Real dt,
166  const amrex::MultiFab& Ex,
167  const amrex::MultiFab& Ey,
168  const amrex::MultiFab& Ez,
169  const amrex::MultiFab& Bx,
170  const amrex::MultiFab& By,
171  const amrex::MultiFab& Bz) = 0;
172 
187  void DepositCurrent (amrex::Vector<std::array< std::unique_ptr<amrex::MultiFab>, 3 > >& J,
188  const amrex::Real dt, const amrex::Real relative_t);
189 
200  void DepositCharge (amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
201  const bool local = false, const bool reset = false,
202  const bool do_rz_volume_scaling = false,
203  const bool interpolate_across_levels = true,
204  const int icomp = 0);
205 
206  std::unique_ptr<amrex::MultiFab> GetChargeDensity(int lev, bool local = false);
207 
208  virtual void DepositCharge (WarpXParIter& pti,
209  RealVector& wp,
210  const int * const ion_lev,
211  amrex::MultiFab* rho,
212  int icomp,
213  const long offset,
214  const long np_to_depose,
215  int thread_num,
216  int lev,
217  int depos_lev);
218 
219  virtual void DepositCurrent (WarpXParIter& pti,
220  RealVector const & wp,
221  RealVector const & uxp,
222  RealVector const & uyp,
223  RealVector const & uzp,
224  int const * const ion_lev,
225  amrex::MultiFab* const jx,
226  amrex::MultiFab* const jy,
227  amrex::MultiFab* const jz,
228  long const offset,
229  long const np_to_depose,
230  int const thread_num,
231  int const lev,
232  int const depos_lev,
233  amrex::Real const dt,
234  amrex::Real const relative_time);
235 
236  // If particles start outside of the domain, ContinuousInjection
237  // makes sure that they are initialized when they enter the domain, and
238  // NOT before. Virtual function, overriden by derived classes.
239  // Current status:
240  // PhysicalParticleContainer: implemented.
241  // LaserParticleContainer: implemented.
242  // RigidInjectedParticleContainer: not implemented.
243  virtual void ContinuousInjection(const amrex::RealBox& /*injection_box*/) {}
244  // Update optional sub-class-specific injection location.
245  virtual void UpdateContinuousInjectionPosition(amrex::Real /*dt*/) {}
246 
247  // Inject a continuous flux of particles from a defined plane
248  virtual void ContinuousFluxInjection(amrex::Real /*dt*/) {}
249 
254  amrex::Real sumParticleCharge(bool local = false);
255 
256  std::array<amrex::Real, 3> meanParticleVelocity(bool local = false);
257 
258  amrex::Real maxParticleVelocity(bool local = false);
259 
260  void AddNParticles (int lev,
261  int n, const amrex::ParticleReal* x, const amrex::ParticleReal* y, const amrex::ParticleReal* z,
262  const amrex::ParticleReal* vx, const amrex::ParticleReal* vy, const amrex::ParticleReal* vz,
263  int nattr, const amrex::ParticleReal* attr, int uniqueparticles, amrex::Long id=-1);
264 
265  virtual void ReadHeader (std::istream& is) = 0;
266 
267  virtual void WriteHeader (std::ostream& os) const = 0;
268 
269  virtual void ConvertUnits (ConvertDirection /*convert_dir*/){}
270 
271  static void ReadParameters ();
272 
273  static void BackwardCompatibility ();
274 
280  void ApplyBoundaryConditions ();
281 
282  bool do_splitting = false;
283  bool initialize_self_fields = false;
284  amrex::Real self_fields_required_precision =
285  amrex::Real(1.e-11);
286  int self_fields_max_iters = 200;
287  int self_fields_verbosity = 2;
288 
289  // split along diagonals (0) or axes (1)
290  int split_type = 0;
291 
292  using amrex::ParticleContainer<0, 0, PIdx::nattribs>::AddRealComp;
293  using amrex::ParticleContainer<0, 0, PIdx::nattribs>::AddIntComp;
294 
295  void AddRealComp (const std::string& name, bool comm=true)
296  {
297  particle_comps[name] = NumRealComps();
298  particle_runtime_comps[name] = NumRealComps() - PIdx::nattribs;
299  AddRealComp(comm);
300  }
301 
302  void AddIntComp (const std::string& name, bool comm=true)
303  {
304  particle_icomps[name] = NumIntComps();
305  particle_runtime_icomps[name] = NumIntComps() - 0;
306  AddIntComp(comm);
307  }
308 
309  int doBackTransformedDiagnostics () const { return do_back_transformed_diagnostics; }
310 
311  std::map<std::string, int> getParticleComps () const noexcept { return particle_comps;}
312  std::map<std::string, int> getParticleiComps () const noexcept { return particle_icomps;}
313  std::map<std::string, int> getParticleRuntimeComps () const noexcept { return particle_runtime_comps;}
314  std::map<std::string, int> getParticleRuntimeiComps () const noexcept { return particle_runtime_icomps;}
315 
316  //amrex::Real getCharge () {return charge;}
317  amrex::ParticleReal getCharge () const {return charge;}
318  //amrex::Real getMass () {return mass;}
319  amrex::ParticleReal getMass () const {return mass;}
320 
321  int DoFieldIonization() const { return do_field_ionization; }
322 
323 #ifdef WARPX_QED
324  //Species for which QED effects are relevant should override these methods
325  virtual bool has_quantum_sync() const {return false;}
326  virtual bool has_breit_wheeler() const {return false;}
327 
328  int DoQED() const { return has_quantum_sync() || has_breit_wheeler(); }
329 #else
330  int DoQED() const { return false; }
331 #endif
332 
333  /* \brief This function tests if the current species
334  * is of a given PhysicalSpecies (specified as a template parameter).
335  * @tparam PhysSpec the PhysicalSpecies to test against
336  * @return the result of the test
337  */
338  template<PhysicalSpecies PhysSpec>
339  bool AmIA () const noexcept {return (physical_species == PhysSpec);}
340 
341  amrex::Array<amrex::Real,3> get_v_galilean () {return m_v_galilean;}
342 
349  virtual void resample (const int /*timestep*/) {}
350 
357  void defineAllParticleTiles () noexcept;
358 
359 protected:
360  amrex::Array<amrex::Real,3> m_v_galilean = {{0}};
361  std::map<std::string, int> particle_comps;
362  std::map<std::string, int> particle_icomps;
363  std::map<std::string, int> particle_runtime_comps;
364  std::map<std::string, int> particle_runtime_icomps;
365 
367 
368  amrex::Real charge;
369  amrex::Real mass;
371 
372  // Controls boundaries for particles exiting the domain
374 
376  bool m_deposit_on_main_grid = false;
377 
379  bool m_gather_from_main_grid = false;
380 
381  int do_not_push = 0;
382  int do_not_deposit = 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;
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 
404  int do_back_transformed_diagnostics = 1;
405 
406 #ifdef WARPX_QED
407  //Species can receive a shared pointer to a QED engine (species for
408  //which this is relevant should override these functions)
409  virtual void
410  set_breit_wheeler_engine_ptr(std::shared_ptr<BreitWheelerEngine>){}
411  virtual void
412  set_quantum_sync_engine_ptr(std::shared_ptr<QuantumSynchrotronEngine>){}
413 
420 
421 #endif
422  amrex::Vector<amrex::FArrayBox> local_rho;
423  amrex::Vector<amrex::FArrayBox> local_jx;
424  amrex::Vector<amrex::FArrayBox> local_jy;
425  amrex::Vector<amrex::FArrayBox> local_jz;
426 
427 public:
428  using PairIndex = std::pair<int, int>;
429  using TmpParticleTile = std::array<amrex::Gpu::DeviceVector<amrex::ParticleReal>,
431  using TmpParticles = amrex::Vector<std::map<PairIndex, TmpParticleTile> >;
432 
433 protected:
435 
436 private:
437  virtual void particlePostLocate(ParticleType& p, const amrex::ParticleLocData& pld,
438  const int lev) override;
439 
440 };
441 
442 #endif
const RealVector & GetAttribs(int comp) const
Definition: WarpXParticleContainer.H:76
virtual void UpdateContinuousInjectionPosition(amrex::Real)
Definition: WarpXParticleContainer.H:245
DtType
Definition: WarpXDtType.H:10
std::array< amrex::Gpu::DeviceVector< amrex::ParticleReal >, TmpIdx::nattribs > TmpParticleTile
Definition: WarpXParticleContainer.H:430
virtual void GetParticleSlice(const int, const amrex::Real, const amrex::Real, const amrex::Real, const amrex::Real, const amrex::Real, DiagnosticParticles &)
Definition: WarpXParticleContainer.H:147
virtual void set_breit_wheeler_engine_ptr(std::shared_ptr< BreitWheelerEngine >)
Definition: WarpXParticleContainer.H:410
int m_qed_breit_wheeler_ele_product
Definition: WarpXParticleContainer.H:414
amrex::Vector< amrex::FArrayBox > local_rho
Definition: WarpXParticleContainer.H:422
Definition: WarpXParticleContainer_fwd.H:28
Definition: WarpXParticleContainer_fwd.H:30
amrex::Gpu::DeviceVector< amrex::Real > adk_power
Definition: WarpXParticleContainer.H:397
int m_qed_breit_wheeler_pos_product
Definition: WarpXParticleContainer.H:416
void AddRealComp(const std::string &name, bool comm=true)
Definition: WarpXParticleContainer.H:295
std::map< std::string, int > getParticleRuntimeiComps() const noexcept
Definition: WarpXParticleContainer.H:314
RealVector & GetAttribs(int comp)
Definition: WarpXParticleContainer.H:80
std::map< std::string, int > particle_icomps
Definition: WarpXParticleContainer.H:362
IntVector & GetiAttribs(int comp)
Definition: WarpXParticleContainer.H:84
std::map< std::string, int > particle_runtime_icomps
Definition: WarpXParticleContainer.H:364
def x
Definition: read_lab_particles.py:25
int DoFieldIonization() const
Definition: WarpXParticleContainer.H:321
Definition: ParticleBoundaries.H:19
PhysicalSpecies
Definition: SpeciesPhysicalProperties.H:19
const std::array< RealVector, PIdx::nattribs > & GetAttribs() const
Definition: WarpXParticleContainer.H:68
Definition: WarpXParticleContainer_fwd.H:28
virtual void set_quantum_sync_engine_ptr(std::shared_ptr< QuantumSynchrotronEngine >)
Definition: WarpXParticleContainer.H:412
bool AmIA() const noexcept
Definition: WarpXParticleContainer.H:339
int DoQED() const
Definition: WarpXParticleContainer.H:328
Definition: WarpXParticleContainer_fwd.H:32
virtual void ContinuousInjection(const amrex::RealBox &)
Definition: WarpXParticleContainer.H:243
virtual void ConvertUnits(ConvertDirection)
Definition: WarpXParticleContainer.H:269
amrex::ParticleReal getCharge() const
Definition: WarpXParticleContainer.H:317
amrex::Gpu::DeviceVector< amrex::Real > adk_prefactor
Definition: WarpXParticleContainer.H:398
def z
Definition: read_lab_particles.py:26
amrex::ParticleReal getMass() const
Definition: WarpXParticleContainer.H:319
ConvertDirection
Definition: WarpXParticleContainer_fwd.H:22
PhysicalSpecies physical_species
Definition: WarpXParticleContainer.H:370
std::map< std::string, int > getParticleComps() const noexcept
Definition: WarpXParticleContainer.H:311
amrex::Vector< amrex::FArrayBox > local_jz
Definition: WarpXParticleContainer.H:425
amrex::Vector< std::map< std::pair< int, int >, DiagnosticParticleData > > DiagnosticParticles
Definition: WarpXParticleContainer.H:123
virtual ~WarpXParticleContainer()
Definition: WarpXParticleContainer.H:126
amrex::StructOfArrays< DiagIdx::nattribs, 0 > DiagnosticParticleData
Definition: WarpXParticleContainer.H:118
Definition: WarpXParticleContainer_fwd.H:28
int doBackTransformedDiagnostics() const
Definition: WarpXParticleContainer.H:309
std::map< std::string, int > particle_comps
Definition: WarpXParticleContainer.H:361
amrex::Real charge
Definition: WarpXParticleContainer.H:368
std::string ionization_product_name
Definition: WarpXParticleContainer.H:393
amrex::Gpu::DeviceVector< amrex::Real > ionization_energies
Definition: WarpXParticleContainer.H:396
TmpParticles tmp_particle_data
Definition: WarpXParticleContainer.H:434
int n
Definition: run_libensemble_on_warpx.py:68
void AddIntComp(const std::string &name, bool comm=true)
Definition: WarpXParticleContainer.H:302
virtual void ContinuousFluxInjection(amrex::Real)
Definition: WarpXParticleContainer.H:248
int species_id
Definition: WarpXParticleContainer.H:366
amrex::Gpu::DeviceVector< amrex::Real > adk_exp_prefactor
Definition: WarpXParticleContainer.H:399
virtual bool has_quantum_sync() const
Definition: WarpXParticleContainer.H:325
amrex::Array< amrex::Real, 3 > get_v_galilean()
Definition: WarpXParticleContainer.H:341
amrex::Vector< std::map< PairIndex, TmpParticleTile > > TmpParticles
Definition: WarpXParticleContainer.H:431
Definition: WarpXParticleContainer.H:45
std::string physical_element
Definition: WarpXParticleContainer.H:400
std::map< std::string, int > getParticleRuntimeComps() const noexcept
Definition: WarpXParticleContainer.H:313
int ionization_product
Definition: WarpXParticleContainer.H:392
name
Definition: run_automated.py:204
WarpXParticleContainer::ParticleType ParticleType
Definition: ParticleUtils.cpp:29
std::map< std::string, int > particle_runtime_comps
Definition: WarpXParticleContainer.H:363
std::string m_qed_quantum_sync_phot_product_name
Definition: WarpXParticleContainer.H:419
virtual void resample(const int)
Virtual method to resample the species. Overriden by PhysicalParticleContainer only. Empty body is here because making the method purely virtual would mean that we need to override the method for every derived class. Note that in practice this function is never called because resample() is only called for PhysicalParticleContainers.
Definition: WarpXParticleContainer.H:349
virtual bool has_breit_wheeler() const
Definition: WarpXParticleContainer.H:326
std::string m_qed_breit_wheeler_ele_product_name
Definition: WarpXParticleContainer.H:415
const std::map< std::string, int > to_index
Definition: WarpXParticleContainer.H:47
std::pair< int, int > PairIndex
Definition: WarpXParticleContainer.H:428
amrex::Real mass
Definition: WarpXParticleContainer.H:369
amrex::Vector< amrex::FArrayBox > local_jx
Definition: WarpXParticleContainer.H:423
std::string m_qed_breit_wheeler_pos_product_name
Definition: WarpXParticleContainer.H:417
Definition: WarpXParticleContainer_fwd.H:50
int ion_atomic_number
Definition: WarpXParticleContainer.H:394
Definition: WarpXParticleContainer_fwd.H:27
std::map< std::string, int > getParticleiComps() const noexcept
Definition: WarpXParticleContainer.H:312
Definition: WarpXParticleContainer.H:58
std::array< RealVector, PIdx::nattribs > & GetAttribs()
Definition: WarpXParticleContainer.H:72
int m_qed_quantum_sync_phot_product
Definition: WarpXParticleContainer.H:418
friend MultiParticleContainer
Definition: WarpXParticleContainer.H:114
Definition: WarpXParticleContainer.H:110
amrex::Vector< amrex::FArrayBox > local_jy
Definition: WarpXParticleContainer.H:424
ParticleBoundaries m_boundary_conditions
Definition: WarpXParticleContainer.H:373