WarpX
MultiParticleContainer.H
Go to the documentation of this file.
1 /* Copyright 2019-2020 Andrew Myers, Ann Almgren, Axel Huebl
2  * David Grote, Jean-Luc Vay, Junmin Gu
3  * Luca Fedeli, Mathieu Lobet, Maxence Thevenet
4  * Remi Lehe, Revathi Jambunathan, Weiqun Zhang
5  * Yinjian Zhao
6  *
7  * This file is part of WarpX.
8  *
9  * License: BSD-3-Clause-LBNL
10  */
11 #ifndef WARPX_ParticleContainer_H_
12 #define WARPX_ParticleContainer_H_
13 
14 #include "WarpXParticleContainer.H"
25 
26 #ifdef WARPX_QED
30 #endif
31 
32 #include <AMReX_Particles.H>
33 
34 #include <memory>
35 #include <map>
36 #include <string>
37 #include <algorithm>
38 
54 {
55 
56 public:
57 
58  MultiParticleContainer (amrex::AmrCore* amr_core);
59 
61 
63  GetParticleContainer (int ispecies) const {return *allcontainers[ispecies];}
64 
66  GetParticleContainerPtr (int ispecies) const {return allcontainers[ispecies].get();}
67 
68 #ifdef WARPX_USE_OPENPMD
69  std::unique_ptr<WarpXParticleContainer>& GetUniqueContainer(int ispecies) {
70  return allcontainers[ispecies];
71  }
72 #endif
73  std::array<amrex::Real, 3> meanParticleVelocity(int ispecies) {
74  return allcontainers[ispecies]->meanParticleVelocity();
75  }
76 
77  void AllocData ();
78 
79  void InitData ();
80 
86  void Evolve (int lev,
87  const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
88  const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz,
89  const amrex::MultiFab& Ex_avg, const amrex::MultiFab& Ey_avg, const amrex::MultiFab& Ez_avg,
90  const amrex::MultiFab& Bx_avg, const amrex::MultiFab& By_avg, const amrex::MultiFab& Bz_avg,
91  amrex::MultiFab& jx, amrex::MultiFab& jy, amrex::MultiFab& jz,
92  amrex::MultiFab* cjx, amrex::MultiFab* cjy, amrex::MultiFab* cjz,
93  amrex::MultiFab* rho, amrex::MultiFab* crho,
94  const amrex::MultiFab* cEx, const amrex::MultiFab* cEy, const amrex::MultiFab* cEz,
95  const amrex::MultiFab* cBx, const amrex::MultiFab* cBy, const amrex::MultiFab* cBz,
96  amrex::Real t, amrex::Real dt, DtType a_dt_type=DtType::Full);
97 
103  void PushX (amrex::Real dt);
104 
111  void PushP (int lev, amrex::Real dt,
112  const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
113  const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz);
114 
119  std::unique_ptr<amrex::MultiFab> GetChargeDensity(int lev, bool local = false);
120 
121  void doFieldIonization (int lev,
122  const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
123  const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz);
124 
125  void doCoulombCollisions ();
126 
132  void doResampling (const int timestep);
133 
134 #ifdef WARPX_QED
135 
142  void doQEDSchwinger ();
143 #endif
144 
145  void Restart (const std::string& dir);
146 
147  void PostRestart ();
148 
149  void ReadHeader (std::istream& is);
150 
151  void WriteHeader (std::ostream& os) const;
152 
153  void SortParticlesByBin (amrex::IntVect bin_size);
154 
155  void Redistribute ();
156 
157  void RedistributeLocal (const int num_ghost);
158 
161  void ApplyBoundaryConditions ();
162 
163  amrex::Vector<long> NumberOfParticlesInGrid(int lev) const;
164 
165  void Increment (amrex::MultiFab& mf, int lev);
166 
167  void SetParticleBoxArray (int lev, amrex::BoxArray& new_ba);
168  void SetParticleDistributionMap (int lev, amrex::DistributionMapping& new_dm);
169 
170  int nSpecies() const {return species_names.size();}
171 
175 
177  bool const onMainGrid = true;
178  auto const & v = m_deposit_on_main_grid;
179  return std::count( v.begin(), v.end(), onMainGrid );
180  }
181 
183  bool const fromMainGrid = true;
184  auto const & v = m_gather_from_main_grid;
185  return std::count( v.begin(), v.end(), fromMainGrid );
186  }
187 
188  void GetLabFrameData(const std::string& snapshot_name,
189  const int i_lab, const int direction,
190  const amrex::Real z_old, const amrex::Real z_new,
191  const amrex::Real t_boost, const amrex::Real t_lab, const amrex::Real dt,
192  amrex::Vector<WarpXParticleContainer::DiagnosticParticleData>& parts) const;
193 
194  // Inject particles during the simulation (for particles entering the
195  // simulation domain after some iterations, due to flowing plasma and/or
196  // moving window).
197  void ContinuousInjection(const amrex::RealBox& injection_box) const;
198  // Update injection position for continuously-injected species.
199  void UpdateContinuousInjectionPosition(amrex::Real dt) const;
200  int doContinuousInjection() const;
201 
202  std::vector<std::string> GetSpeciesNames() const { return species_names; }
203 
205 
206  std::string m_B_ext_particle_s = "default";
207  std::string m_E_ext_particle_s = "default";
208  // External fields added to particle fields.
209  amrex::Vector<amrex::Real> m_B_external_particle;
210  amrex::Vector<amrex::Real> m_E_external_particle;
211  // ParserWrapper for B_external on the particle
212  std::unique_ptr<ParserWrapper<4> > m_Bx_particle_parser;
213  std::unique_ptr<ParserWrapper<4> > m_By_particle_parser;
214  std::unique_ptr<ParserWrapper<4> > m_Bz_particle_parser;
215  // ParserWrapper for E_external on the particle
216  std::unique_ptr<ParserWrapper<4> > m_Ex_particle_parser;
217  std::unique_ptr<ParserWrapper<4> > m_Ey_particle_parser;
218  std::unique_ptr<ParserWrapper<4> > m_Ez_particle_parser;
219 
220 #ifdef WARPX_QED
221 
224  void doQedEvents (int lev,
225  const amrex::MultiFab& Ex,
226  const amrex::MultiFab& Ey,
227  const amrex::MultiFab& Ez,
228  const amrex::MultiFab& Bx,
229  const amrex::MultiFab& By,
230  const amrex::MultiFab& Bz);
231 #endif
232 
233  int getSpeciesID (std::string product_str) const;
234 
235 protected:
236 
237 #ifdef WARPX_QED
238 
241  void doQedBreitWheeler (int lev,
242  const amrex::MultiFab& Ex,
243  const amrex::MultiFab& Ey,
244  const amrex::MultiFab& Ez,
245  const amrex::MultiFab& Bx,
246  const amrex::MultiFab& By,
247  const amrex::MultiFab& Bz);
248 
252  void doQedQuantumSync (int lev,
253  const amrex::MultiFab& Ex,
254  const amrex::MultiFab& Ey,
255  const amrex::MultiFab& Ez,
256  const amrex::MultiFab& Bx,
257  const amrex::MultiFab& By,
258  const amrex::MultiFab& Bz);
259 #endif
260 
261  // Particle container types
263 
264  std::vector<std::string> species_names;
265 
266  std::vector<std::string> lasers_names;
267 
268  std::vector<std::string> collision_names;
269 
270  amrex::Vector<std::unique_ptr<CollisionType> > allcollisions;
271 
273  std::vector<bool> m_deposit_on_main_grid;
274 
276  std::vector<bool> m_gather_from_main_grid;
277 
278  std::vector<PCTypes> species_types;
279 
282 
283  template<typename ...Args>
284  amrex::MFItInfo getMFItInfo (const WarpXParticleContainer& pc_src,
285  Args const&... pc_dsts) const noexcept
286  {
287  amrex::MFItInfo info;
288 
289  MFItInfoCheckTiling(pc_src, pc_dsts...);
290 
291  if (pc_src.do_tiling && amrex::Gpu::notInLaunchRegion()) {
292  info.EnableTiling(pc_src.tile_size);
293  }
294 
295 #ifdef _OPENMP
296  info.SetDynamic(true);
297 #endif
298 
299  return info;
300  }
301 
302 
303 #ifdef WARPX_QED
304  // The QED engines
305  std::shared_ptr<BreitWheelerEngine> m_shr_p_bw_engine;
306  std::shared_ptr<QuantumSynchrotronEngine> m_shr_p_qs_engine;
307  //_______________________________
308 
313  void InitQED ();
314 
315  //Variables to store how many species need a QED process
318  //________
319 
321  static_cast<amrex::ParticleReal>(
322  2.0 * PhysConst::m_e * PhysConst::c * PhysConst::c );
332 
337 
341  void InitQuantumSync ();
342 
346  void InitBreitWheeler ();
347 
353 
359 
377 
378 #endif
379 
380 private:
381 
382  // physical particles (+ laser)
383  amrex::Vector<std::unique_ptr<WarpXParticleContainer>> allcontainers;
384  // Temporary particle container, used e.g. for particle splitting.
385  std::unique_ptr<PhysicalParticleContainer> pc_tmp;
386 
387  void ReadParameters ();
388 
389  void mapSpeciesProduct ();
390 
391  // Number of species dumped in BackTransformedDiagnostics
393  // map_species_back_transformed_diagnostics[i] is the species ID in
394  // MultiParticleContainer for 0<i<nspecies_back_transformed_diagnostics
397 
398  void MFItInfoCheckTiling(const WarpXParticleContainer& /*pc_src*/) const noexcept
399  {
400  return;
401  }
402 
403  template<typename First, typename ...Args>
405  First const& pc_dst, Args const&... others) const noexcept
406  {
407  if (pc_src.do_tiling && amrex::Gpu::notInLaunchRegion()) {
408  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(pc_dst.do_tiling,
409  "For particle creation processes, either all or none of the "
410  "particle species must use tiling.");
411  }
412 
413  MFItInfoCheckTiling(pc_src, others...);
414  }
415 
422 
423 #ifdef WARPX_QED
424 
429  void CheckQEDProductSpecies();
430 #endif
431 
432 
433 };
434 #endif /*WARPX_ParticleContainer_H_*/
int getSpeciesID(std::string product_str) const
Definition: MultiParticleContainer.cpp:608
void doQedQuantumSync(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)
Performs QED photon emission for the species for which it is enabled.
Definition: MultiParticleContainer.cpp:1222
ParticleBC
Definition: WarpXParticleContainer.H:29
DtType
Definition: WarpXDtType.H:10
std::unique_ptr< amrex::MultiFab > GetChargeDensity(int lev, bool local=false)
Definition: MultiParticleContainer.cpp:360
void Increment(amrex::MultiFab &mf, int lev)
Definition: MultiParticleContainer.cpp:422
std::vector< int > map_species_back_transformed_diagnostics
Definition: MultiParticleContainer.H:395
void InitQED()
Definition: MultiParticleContainer.cpp:741
void doQEDSchwinger()
Definition: MultiParticleContainer.cpp:1043
std::unique_ptr< ParserWrapper< 4 > > m_Ey_particle_parser
Definition: MultiParticleContainer.H:217
std::string m_B_ext_particle_s
Definition: MultiParticleContainer.H:206
void ApplyBoundaryConditions()
Definition: MultiParticleContainer.cpp:399
int NSpeciesQuantumSync() const
Definition: MultiParticleContainer.H:331
std::vector< std::string > lasers_names
Definition: MultiParticleContainer.H:266
std::vector< bool > m_deposit_on_main_grid
instead of depositing (current, charge) on the finest patch level, deposit to the coarsest grid ...
Definition: MultiParticleContainer.H:273
std::vector< std::string > species_names
Definition: MultiParticleContainer.H:264
void doQedEvents(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)
Performs QED events (Breit-Wheeler process and photon emission)
Definition: MultiParticleContainer.cpp:1141
ParticleBC m_boundary_conditions
Definition: MultiParticleContainer.H:281
Definition: PhysicalParticleContainer.H:39
Definition: MultiParticleContainer.H:53
int nSpecies() const
Definition: MultiParticleContainer.H:170
int m_qed_schwinger_pos_product
Definition: MultiParticleContainer.H:369
PCTypes
Definition: MultiParticleContainer.H:262
void GetLabFrameData(const std::string &snapshot_name, const int i_lab, const int direction, const amrex::Real z_old, const amrex::Real z_new, const amrex::Real t_boost, const amrex::Real t_lab, const amrex::Real dt, amrex::Vector< WarpXParticleContainer::DiagnosticParticleData > &parts) const
Definition: MultiParticleContainer.cpp:456
std::shared_ptr< QuantumSynchrotronEngine > m_shr_p_qs_engine
Definition: MultiParticleContainer.H:306
std::unique_ptr< WarpXParticleContainer > & GetUniqueContainer(int ispecies)
Definition: MultiParticleContainer.H:69
void QuantumSyncGenerateTable()
Definition: MultiParticleContainer.cpp:885
void mapSpeciesProduct()
Definition: MultiParticleContainer.cpp:563
std::vector< std::string > collision_names
Definition: MultiParticleContainer.H:268
int count
Definition: run_alltests.py:318
void doResampling(const int timestep)
This function loops over all species and performs resampling if appropriate.
Definition: MultiParticleContainer.cpp:716
void ContinuousInjection(const amrex::RealBox &injection_box) const
Definition: MultiParticleContainer.cpp:522
bool m_do_qed_schwinger
Definition: MultiParticleContainer.H:361
~MultiParticleContainer()
Definition: MultiParticleContainer.H:60
void SetParticleBoxArray(int lev, amrex::BoxArray &new_ba)
Definition: MultiParticleContainer.cpp:430
void ReadParameters()
Definition: MultiParticleContainer.cpp:85
direction
Definition: AnyFFT.H:60
void InitData()
Definition: MultiParticleContainer.cpp:295
std::unique_ptr< ParserWrapper< 4 > > m_Ex_particle_parser
Definition: MultiParticleContainer.H:216
int m_nspecies_quantum_sync
Definition: MultiParticleContainer.H:316
int mapSpeciesBackTransformedDiagnostics(int i) const
Definition: MultiParticleContainer.H:173
std::vector< PCTypes > species_types
Definition: MultiParticleContainer.H:278
std::string m_qed_schwinger_ele_product_name
Definition: MultiParticleContainer.H:363
int doContinuousInjection() const
Definition: MultiParticleContainer.cpp:547
int m_qed_schwinger_ele_product
Definition: MultiParticleContainer.H:367
std::string m_qed_schwinger_pos_product_name
Definition: MultiParticleContainer.H:365
void Restart(const std::string &dir)
Definition: ParticleIO.cpp:78
PhysicalParticleContainer & GetPCtmp()
Definition: MultiParticleContainer.H:204
void CheckQEDProductSpecies()
Definition: MultiParticleContainer.cpp:1287
void MFItInfoCheckTiling(const WarpXParticleContainer &pc_src, First const &pc_dst, Args const &... others) const noexcept
Definition: MultiParticleContainer.H:404
std::string m_E_ext_particle_s
Definition: MultiParticleContainer.H:207
WarpXParticleContainer * GetParticleContainerPtr(int ispecies) const
Definition: MultiParticleContainer.H:66
i
Definition: check_interp_points_and_weights.py:171
void PushX(amrex::Real dt)
Definition: MultiParticleContainer.cpp:342
void AllocData()
Definition: MultiParticleContainer.cpp:286
int nSpeciesGatherFromMainGrid() const
Definition: MultiParticleContainer.H:182
void InitBreitWheeler()
Definition: MultiParticleContainer.cpp:833
void doCoulombCollisions()
Definition: MultiParticleContainer.cpp:683
std::unique_ptr< ParserWrapper< 4 > > m_By_particle_parser
Definition: MultiParticleContainer.H:213
void BreitWheelerGenerateTable()
Definition: MultiParticleContainer.cpp:966
std::unique_ptr< ParserWrapper< 4 > > m_Bz_particle_parser
Definition: MultiParticleContainer.H:214
std::vector< std::string > GetSpeciesNames() const
Definition: MultiParticleContainer.H:202
MultiParticleContainer(amrex::AmrCore *amr_core)
Definition: MultiParticleContainer.cpp:31
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, const amrex::MultiFab &Ex_avg, const amrex::MultiFab &Ey_avg, const amrex::MultiFab &Ez_avg, const amrex::MultiFab &Bx_avg, const amrex::MultiFab &By_avg, const amrex::MultiFab &Bz_avg, 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)
Definition: MultiParticleContainer.cpp:315
int m_qed_schwinger_threshold_poisson_gaussian
Definition: MultiParticleContainer.H:376
WarpXParticleContainer & GetParticleContainer(int ispecies) const
Definition: MultiParticleContainer.H:63
int nSpeciesDepositOnMainGrid() const
Definition: MultiParticleContainer.H:176
std::array< amrex::Real, 3 > meanParticleVelocity(int ispecies)
Definition: MultiParticleContainer.H:73
std::unique_ptr< PhysicalParticleContainer > pc_tmp
Definition: MultiParticleContainer.H:385
void CheckIonizationProductSpecies()
Definition: MultiParticleContainer.cpp:729
int doBackTransformedDiagnostics() const
Definition: MultiParticleContainer.H:174
std::vector< bool > m_gather_from_main_grid
instead of gathering fields from the finest patch level, gather from the coarsest ...
Definition: MultiParticleContainer.H:276
void UpdateContinuousInjectionPosition(amrex::Real dt) const
Definition: MultiParticleContainer.cpp:537
amrex::Vector< amrex::Real > m_B_external_particle
Definition: MultiParticleContainer.H:209
void MFItInfoCheckTiling(const WarpXParticleContainer &) const noexcept
Definition: MultiParticleContainer.H:398
amrex::Real m_qed_schwinger_y_size
Definition: MultiParticleContainer.H:371
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)
Definition: MultiParticleContainer.cpp:350
std::unique_ptr< ParserWrapper< 4 > > m_Ez_particle_parser
Definition: MultiParticleContainer.H:218
int nspecies_back_transformed_diagnostics
Definition: MultiParticleContainer.H:392
int m_nspecies_breit_wheeler
Definition: MultiParticleContainer.H:317
amrex::MFItInfo getMFItInfo(const WarpXParticleContainer &pc_src, Args const &... pc_dsts) const noexcept
Definition: MultiParticleContainer.H:284
int nSpeciesBackTransformedDiagnostics() const
Definition: MultiParticleContainer.H:172
static constexpr auto m_default_quantum_sync_photon_creation_energy_threshold
Definition: MultiParticleContainer.H:320
void ReadHeader(std::istream &is)
Definition: ParticleIO.cpp:86
void doQedBreitWheeler(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)
Performs Breit-Wheeler process for the species for which it is enabled.
Definition: MultiParticleContainer.cpp:1155
void InitQuantumSync()
Definition: MultiParticleContainer.cpp:770
void SortParticlesByBin(amrex::IntVect bin_size)
Definition: MultiParticleContainer.cpp:375
int NSpeciesBreitWheeler() const
Definition: MultiParticleContainer.H:336
std::unique_ptr< ParserWrapper< 4 > > m_Bx_particle_parser
Definition: MultiParticleContainer.H:212
void WriteHeader(std::ostream &os) const
Definition: ParticleIO.cpp:94
void doFieldIonization(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)
Definition: MultiParticleContainer.cpp:632
void RedistributeLocal(const int num_ghost)
Definition: MultiParticleContainer.cpp:391
std::shared_ptr< BreitWheelerEngine > m_shr_p_bw_engine
Definition: MultiParticleContainer.H:305
void Redistribute()
Definition: MultiParticleContainer.cpp:383
amrex::Vector< std::unique_ptr< CollisionType > > allcollisions
Definition: MultiParticleContainer.H:270
amrex::Vector< long > NumberOfParticlesInGrid(int lev) const
Definition: MultiParticleContainer.cpp:407
amrex::Vector< amrex::Real > m_E_external_particle
Definition: MultiParticleContainer.H:210
void SetParticleDistributionMap(int lev, amrex::DistributionMapping &new_dm)
Definition: MultiParticleContainer.cpp:438
int do_back_transformed_diagnostics
Definition: MultiParticleContainer.H:396
amrex::ParticleReal m_quantum_sync_photon_creation_energy_threshold
Definition: MultiParticleContainer.H:325
Definition: WarpXParticleContainer.H:131
amrex::Vector< std::unique_ptr< WarpXParticleContainer > > allcontainers
Definition: MultiParticleContainer.H:383
void PostRestart()
Definition: MultiParticleContainer.cpp:446