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"
16 #include "Evolve/WarpXPushType.H"
20 
21 #ifdef WARPX_QED
24 #endif
27 
28 #include <AMReX_Array.H>
29 #include <AMReX_FArrayBox.H>
30 #include <AMReX_GpuAllocators.H>
31 #include <AMReX_GpuContainers.H>
32 #include <AMReX_INT.H>
33 #include <AMReX_ParIter.H>
34 #include <AMReX_Particles.H>
35 #include <AMReX_Random.H>
36 #include <AMReX_REAL.H>
37 #include <AMReX_StructOfArrays.H>
38 #include <AMReX_Vector.H>
39 
40 #include <AMReX_BaseFwd.H>
41 #include <AMReX_AmrCoreFwd.H>
42 
43 #include <array>
44 #include <iosfwd>
45 #include <map>
46 #include <memory>
47 #include <string>
48 #include <utility>
49 
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  [[nodiscard]] const std::array<RealVector, PIdx::nattribs>& GetAttribs () const
62  {
63  return GetStructOfArrays().GetRealData();
64  }
65 
66  [[nodiscard]] std::array<RealVector, PIdx::nattribs>& GetAttribs ()
67  {
68  return GetStructOfArrays().GetRealData();
69  }
70 
71  [[nodiscard]] const RealVector& GetAttribs (int comp) const
72  {
73  return GetStructOfArrays().GetRealData(comp);
74  }
75 
76  [[nodiscard]] RealVector& GetAttribs (int comp)
77  {
78  return GetStructOfArrays().GetRealData(comp);
79  }
80 
81  [[nodiscard]] IntVector& GetiAttribs (int comp)
82  {
83  return GetStructOfArrays().GetIntData(comp);
84  }
85 };
86 
109  : public NamedComponentParticleContainer<amrex::DefaultAllocator>
110 {
111 public:
113 
114  // amrex::StructOfArrays with DiagIdx::nattribs amrex::ParticleReal components
115  // and 0 int components for the particle data.
117  // DiagnosticParticles is a vector, with one element per MR level.
118  // DiagnosticParticles[lev] is typically a key-value pair where the key is
119  // a pair [grid_index, tile_index], and the value is the corresponding
120  // DiagnosticParticleData (see above) on this tile.
122 
123  WarpXParticleContainer (amrex::AmrCore* amr_core, int ispecies);
124  ~WarpXParticleContainer() override = default;
125 
126  // Move and copy operations
131 
132  virtual void InitData () = 0;
133 
134  virtual void InitIonizationModule () {}
135 
136  /*
137  * \brief Virtual function that returns a pointer to the plasma injector,
138  * for derived classes that define one (PhysicalParticleContainer).
139  */
140  virtual PlasmaInjector* GetPlasmaInjector (const int /*i*/) { return nullptr; }
141 
147  virtual void Evolve (int lev,
148  const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
149  const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz,
152  amrex::MultiFab* rho, amrex::MultiFab* crho,
153  const amrex::MultiFab* cEx, const amrex::MultiFab* cEy, const amrex::MultiFab* cEz,
154  const amrex::MultiFab* cBx, const amrex::MultiFab* cBy, const amrex::MultiFab* cBz,
155  amrex::Real t, amrex::Real dt, DtType a_dt_type=DtType::Full, bool skip_deposition=false,
156  PushType push_type=PushType::Explicit) = 0;
157 
158  virtual void PostRestart () = 0;
159 
160  void AllocData ();
161 
169  amrex::PinnedArenaAllocator>& pinned_tile,
170  int n_external_attr_real,
171  int n_external_attr_int) = 0;
172 
178  void PushX ( amrex::Real dt);
179  void PushX (int lev, amrex::Real dt);
180 
184  virtual void PushP (int lev, amrex::Real dt,
185  const amrex::MultiFab& Ex,
186  const amrex::MultiFab& Ey,
187  const amrex::MultiFab& Ez,
188  const amrex::MultiFab& Bx,
189  const amrex::MultiFab& By,
190  const amrex::MultiFab& Bz) = 0;
191 
203  void DepositCurrent (amrex::Vector<std::array< std::unique_ptr<amrex::MultiFab>, 3 > >& J,
204  amrex::Real dt, amrex::Real relative_time);
205 
216  void DepositCharge (amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
217  bool local = false, bool reset = false,
218  bool apply_boundary_and_scale_volume = false,
219  bool interpolate_across_levels = true,
220  int icomp = 0);
221  void DepositCharge (std::unique_ptr<amrex::MultiFab>& rho, int lev,
222  bool local = false, bool reset = false,
223  bool apply_boundary_and_scale_volume = false,
224  int icomp = 0);
225 
226  std::unique_ptr<amrex::MultiFab> GetChargeDensity(int lev, bool local = false);
227 
228  virtual void DepositCharge (WarpXParIter& pti,
229  RealVector const & wp,
230  const int* ion_lev,
231  amrex::MultiFab* rho,
232  int icomp,
233  long offset,
234  long np_to_deposit,
235  int thread_num,
236  int lev,
237  int depos_lev);
238 
239  virtual void DepositCurrent (WarpXParIter& pti,
240  RealVector const & wp,
241  RealVector const & uxp,
242  RealVector const & uyp,
243  RealVector const & uzp,
244  int const* ion_lev,
245  amrex::MultiFab* jx,
246  amrex::MultiFab* jy,
247  amrex::MultiFab* jz,
248  long offset,
249  long np_to_deposit,
250  int thread_num,
251  int lev,
252  int depos_lev,
253  amrex::Real dt,
254  amrex::Real relative_time,
255  PushType push_type);
256 
257  // If particles start outside of the domain, ContinuousInjection
258  // makes sure that they are initialized when they enter the domain, and
259  // NOT before. Virtual function, overriden by derived classes.
260  // Current status:
261  // PhysicalParticleContainer: implemented.
262  // LaserParticleContainer: implemented.
263  // RigidInjectedParticleContainer: not implemented.
264  virtual void ContinuousInjection(const amrex::RealBox& /*injection_box*/) {}
265 
270  virtual void UpdateAntennaPosition(const amrex::Real /*dt*/) {}
271 
273 
274  // Inject a continuous flux of particles from a defined plane
275  virtual void ContinuousFluxInjection(amrex::Real /*t*/, amrex::Real /*dt*/) {}
276 
277  int getSpeciesId() const {return species_id;}
278 
283  amrex::ParticleReal sumParticleCharge(bool local = false);
284 
285  std::array<amrex::ParticleReal, 3> meanParticleVelocity(bool local = false);
286 
287  amrex::ParticleReal maxParticleVelocity(bool local = false);
288 
313  void AddNParticles (int lev, long n,
320  int nattr_real,
322  int nattr_int, amrex::Vector<amrex::Vector<int>> const & attr_int,
323  int uniqueparticles, amrex::Long id=-1);
324 
325  virtual void ReadHeader (std::istream& is) = 0;
326 
327  virtual void WriteHeader (std::ostream& os) const = 0;
328 
329  static void ReadParameters ();
330 
331  static void BackwardCompatibility ();
332 
335  void ApplyBoundaryConditions ();
336 
337  bool do_splitting = false;
338  int do_not_deposit = 0;
340  amrex::Real self_fields_required_precision = amrex::Real(1.e-11);
341  amrex::Real self_fields_absolute_tolerance = amrex::Real(0.0);
344 
345  // External fields added to particle fields.
348 
351 
352  // split along diagonals (0) or axes (1)
353  int split_type = 0;
354 
359  void SetDoBackTransformedParticles(const bool do_back_transformed_particles) {
360  m_do_back_transformed_particles = do_back_transformed_particles;
361  }
362 
363  //amrex::Real getCharge () {return charge;}
364  amrex::ParticleReal getCharge () const {return charge;}
365  //amrex::Real getMass () {return mass;}
366  amrex::ParticleReal getMass () const {return mass;}
367 
368  int DoFieldIonization() const { return do_field_ionization; }
369 
370 #ifdef WARPX_QED
371  //Species for which QED effects are relevant should override these methods
372  virtual bool has_quantum_sync() const {return false;}
373  virtual bool has_breit_wheeler() const {return false;}
374 
375  int DoQED() const { return has_quantum_sync() || has_breit_wheeler(); }
376 #else
377  int DoQED() const { return false; }
378 #endif
379 
380  /* \brief This function tests if the current species
381  * is of a given PhysicalSpecies (specified as a template parameter).
382  * @tparam PhysSpec the PhysicalSpecies to test against
383  * @return the result of the test
384  */
385  template<PhysicalSpecies PhysSpec>
386  bool AmIA () const noexcept {return (physical_species == PhysSpec);}
387 
391  std::string getSpeciesTypeName () const {return species::get_name(physical_species);}
392 
399  virtual void resample (const int /*timestep*/, bool /*verbose*/) {}
400 
407  void defineAllParticleTiles () noexcept;
408 
409  virtual std::vector<std::string> getUserIntAttribs () const { return std::vector<std::string>{}; }
410 
411  virtual std::vector<std::string> getUserRealAttribs () const { return std::vector<std::string>{}; }
412 
414 
416 
417 #ifdef WARPX_QED
418  virtual BreitWheelerEngine* get_breit_wheeler_engine_ptr () const {return nullptr;}
419  virtual QuantumSynchrotronEngine* get_quantum_sync_engine_ptr () const {return nullptr;}
420 #endif
421 
422 protected:
424 
425  amrex::ParticleReal charge;
426  amrex::ParticleReal mass;
428 
429  // Controls boundaries for particles exiting the domain
431 
434 
437 
438  int do_not_push = 0;
439  int do_not_gather = 0;
440 
441  // Whether to allow particles outside of the simulation domain to be
442  // initialized when they enter the domain.
443  // This is currently required because continuous injection does not
444  // support all features allowed by direct injection.
446 
456  std::string physical_element;
457 
458  int do_resampling = 0;
459 
462 
463 #ifdef WARPX_QED
464  //Species can receive a shared pointer to a QED engine (species for
465  //which this is relevant should override these functions)
466  virtual void
467  set_breit_wheeler_engine_ptr(std::shared_ptr<BreitWheelerEngine>){}
468  virtual void
469  set_quantum_sync_engine_ptr(std::shared_ptr<QuantumSynchrotronEngine>){}
470 
477 
478 #endif
483 
484 public:
485  using PairIndex = std::pair<int, int>;
486  using TmpParticleTile = std::array<amrex::Gpu::DeviceVector<amrex::ParticleReal>,
489 
491 
492  int getIonizationInitialLevel () const noexcept {return ionization_initial_level;}
493 
494 protected:
496 
497 private:
498  void particlePostLocate(ParticleType& p, const amrex::ParticleLocData& pld, int lev) override;
499 
500 };
501 
502 #endif
PhysicalSpecies
Definition: SpeciesPhysicalProperties.H:16
DtType
Definition: WarpXDtType.H:11
PushType
Definition: WarpXPushType.H:12
Definition: BreitWheelerEngineWrapper.H:294
Definition: NamedComponentParticleContainer.H:49
Definition: PlasmaInjector.H:42
Definition: QuantumSyncEngineWrapper.H:273
Definition: WarpXParticleContainer.H:53
WarpXParIter(ContainerType &pc, int level)
Definition: WarpXParticleContainer.cpp:74
std::array< RealVector, PIdx::nattribs > & GetAttribs()
Definition: WarpXParticleContainer.H:66
const RealVector & GetAttribs(int comp) const
Definition: WarpXParticleContainer.H:71
IntVector & GetiAttribs(int comp)
Definition: WarpXParticleContainer.H:81
RealVector & GetAttribs(int comp)
Definition: WarpXParticleContainer.H:76
const std::array< RealVector, PIdx::nattribs > & GetAttribs() const
Definition: WarpXParticleContainer.H:61
Definition: WarpXParticleContainer.H:110
amrex::Real self_fields_required_precision
Definition: WarpXParticleContainer.H:340
int do_continuous_injection
Definition: WarpXParticleContainer.H:445
virtual std::vector< std::string > getUserIntAttribs() const
Definition: WarpXParticleContainer.H:409
~WarpXParticleContainer() override=default
bool AmIA() const noexcept
Definition: WarpXParticleContainer.H:386
amrex::Vector< amrex::ParticleReal > m_B_external_particle
Definition: WarpXParticleContainer.H:346
bool m_gather_from_main_grid
instead of gathering fields from the finest patch level, gather from the coarsest
Definition: WarpXParticleContainer.H:436
void AllocData()
Definition: WarpXParticleContainer.cpp:149
virtual PlasmaInjector * GetPlasmaInjector(const int)
Definition: WarpXParticleContainer.H:140
void defineAllParticleTiles() noexcept
Definition: WarpXParticleContainer.cpp:1394
virtual void PostRestart()=0
void particlePostLocate(ParticleType &p, const amrex::ParticleLocData &pld, int lev) override
Definition: WarpXParticleContainer.cpp:1411
virtual void ReadHeader(std::istream &is)=0
int self_fields_max_iters
Definition: WarpXParticleContainer.H:342
amrex::Real m_current_injection_position
Current injection position.
Definition: WarpXParticleContainer.H:350
int species_id
Definition: WarpXParticleContainer.H:423
std::string m_qed_breit_wheeler_ele_product_name
Definition: WarpXParticleContainer.H:472
virtual void WriteHeader(std::ostream &os) const =0
amrex::Gpu::DeviceVector< amrex::Real > adk_exp_prefactor
Definition: WarpXParticleContainer.H:455
amrex::Gpu::DeviceVector< amrex::Real > adk_power
Definition: WarpXParticleContainer.H:453
PhysicalSpecies physical_species
Definition: WarpXParticleContainer.H:427
virtual void set_quantum_sync_engine_ptr(std::shared_ptr< QuantumSynchrotronEngine >)
Definition: WarpXParticleContainer.H:469
virtual void ContinuousFluxInjection(amrex::Real, amrex::Real)
Definition: WarpXParticleContainer.H:275
std::string getSpeciesTypeName() const
This function returns a string containing the name of the species type.
Definition: WarpXParticleContainer.H:391
int m_qed_quantum_sync_phot_product
Definition: WarpXParticleContainer.H:475
amrex::ParticleReal maxParticleVelocity(bool local=false)
Definition: WarpXParticleContainer.cpp:1300
virtual void InitIonizationModule()
Definition: WarpXParticleContainer.H:134
virtual void InitData()=0
amrex::Gpu::DeviceVector< amrex::Real > adk_prefactor
Definition: WarpXParticleContainer.H:454
void SetDoBackTransformedParticles(const bool do_back_transformed_particles)
Definition: WarpXParticleContainer.H:359
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, PushType push_type=PushType::Explicit)=0
amrex::ParticleReal sumParticleCharge(bool local=false)
Definition: WarpXParticleContainer.cpp:1185
int ionization_initial_level
Definition: WarpXParticleContainer.H:451
amrex::Real self_fields_absolute_tolerance
Definition: WarpXParticleContainer.H:341
amrex::Vector< amrex::FArrayBox > local_jx
Definition: WarpXParticleContainer.H:480
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:270
amrex::ParticleReal mass
Definition: WarpXParticleContainer.H:426
bool doContinuousInjection() const
Definition: WarpXParticleContainer.H:272
amrex::Vector< amrex::FArrayBox > local_jy
Definition: WarpXParticleContainer.H:481
std::unique_ptr< amrex::MultiFab > GetChargeDensity(int lev, bool local=false)
Definition: WarpXParticleContainer.cpp:1160
static void ReadParameters()
Definition: WarpXParticleContainer.cpp:137
int getIonizationInitialLevel() const noexcept
Definition: WarpXParticleContainer.H:492
std::array< amrex::Gpu::DeviceVector< amrex::ParticleReal >, TmpIdx::nattribs > TmpParticleTile
Definition: WarpXParticleContainer.H:487
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)=0
Virtual method to initialize runtime attributes. Must be overriden by each derived class.
int ion_atomic_number
Definition: WarpXParticleContainer.H:450
bool m_do_back_transformed_particles
Definition: WarpXParticleContainer.H:461
amrex::Vector< amrex::FArrayBox > local_jz
Definition: WarpXParticleContainer.H:482
amrex::Gpu::DeviceVector< amrex::Real > ionization_energies
Definition: WarpXParticleContainer.H:452
amrex::ParticleReal getCharge() const
Definition: WarpXParticleContainer.H:364
TmpParticles tmp_particle_data
Definition: WarpXParticleContainer.H:495
virtual amrex::Vector< amrex::Parser * > getUserRealAttribParser() const
Definition: WarpXParticleContainer.H:415
void ApplyBoundaryConditions()
Apply particle BC.
Definition: WarpXParticleContainer.cpp:1433
virtual void set_breit_wheeler_engine_ptr(std::shared_ptr< BreitWheelerEngine >)
Definition: WarpXParticleContainer.H:467
virtual bool has_breit_wheeler() const
Definition: WarpXParticleContainer.H:373
virtual amrex::Vector< amrex::Parser * > getUserIntAttribParser() const
Definition: WarpXParticleContainer.H:413
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:474
std::string m_qed_quantum_sync_phot_product_name
Definition: WarpXParticleContainer.H:476
virtual void ContinuousInjection(const amrex::RealBox &)
Definition: WarpXParticleContainer.H:264
WarpXParticleContainer & operator=(const WarpXParticleContainer &)=delete
ParticleBoundaries m_boundary_conditions
Definition: WarpXParticleContainer.H:430
int m_qed_breit_wheeler_pos_product
Definition: WarpXParticleContainer.H:473
int ionization_product
Definition: WarpXParticleContainer.H:448
virtual BreitWheelerEngine * get_breit_wheeler_engine_ptr() const
Definition: WarpXParticleContainer.H:418
int m_qed_breit_wheeler_ele_product
Definition: WarpXParticleContainer.H:471
bool do_splitting
Definition: WarpXParticleContainer.H:337
std::array< amrex::ParticleReal, 3 > meanParticleVelocity(bool local=false)
Definition: WarpXParticleContainer.cpp:1214
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:712
bool initialize_self_fields
Definition: WarpXParticleContainer.H:339
int self_fields_verbosity
Definition: WarpXParticleContainer.H:343
virtual std::vector< std::string > getUserRealAttribs() const
Definition: WarpXParticleContainer.H:411
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:1057
virtual bool has_quantum_sync() const
Definition: WarpXParticleContainer.H:372
WarpXParticleContainer(amrex::AmrCore *amr_core, int ispecies)
Definition: WarpXParticleContainer.cpp:86
std::string physical_element
Definition: WarpXParticleContainer.H:456
virtual QuantumSynchrotronEngine * get_quantum_sync_engine_ptr() const
Definition: WarpXParticleContainer.H:419
int split_type
Definition: WarpXParticleContainer.H:353
void PushX(amrex::Real dt)
Definition: WarpXParticleContainer.cpp:1327
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:158
TmpParticles getTmpParticleData() const noexcept
Definition: WarpXParticleContainer.H:490
friend MultiParticleContainer
Definition: WarpXParticleContainer.H:112
amrex::Vector< amrex::ParticleReal > m_E_external_particle
Definition: WarpXParticleContainer.H:347
std::string ionization_product_name
Definition: WarpXParticleContainer.H:449
int getSpeciesId() const
Definition: WarpXParticleContainer.H:277
virtual void resample(const int, bool)
Virtual method to resample the species. Overriden by PhysicalParticleContainer only....
Definition: WarpXParticleContainer.H:399
int DoFieldIonization() const
Definition: WarpXParticleContainer.H:368
int do_not_deposit
Definition: WarpXParticleContainer.H:338
WarpXParticleContainer(const WarpXParticleContainer &)=delete
amrex::ParticleReal getMass() const
Definition: WarpXParticleContainer.H:366
int do_not_push
Definition: WarpXParticleContainer.H:438
int DoQED() const
Definition: WarpXParticleContainer.H:375
int do_resampling
Definition: WarpXParticleContainer.H:458
int do_not_gather
Definition: WarpXParticleContainer.H:439
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:433
amrex::Vector< amrex::FArrayBox > local_rho
Definition: WarpXParticleContainer.H:479
std::pair< int, int > PairIndex
Definition: WarpXParticleContainer.H:485
int do_field_ionization
Definition: WarpXParticleContainer.H:447
amrex::ParticleReal charge
Definition: WarpXParticleContainer.H:425
typename SoA::RealVector RealVector
typename SoA::IntVector IntVector
static constexpr int NArrayInt
static constexpr int NArrayReal
typename SoA::RealVector RealVector
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:295
float dt
Definition: stencil.py:442
Definition: ParticleBoundaries.H:19
@ nattribs
Definition: WarpXParticleContainer_fwd.H:39