WarpX
RigidInjectedParticleContainer.H
Go to the documentation of this file.
1 /* Copyright 2019 Andrew Myers, David Grote, Maxence Thevenet
2  * Weiqun Zhang
3  *
4  * This file is part of WarpX.
5  *
6  * License: BSD-3-Clause-LBNL
7  */
8 #ifndef WARPX_RigidInjectedParticleContainer_H_
9 #define WARPX_RigidInjectedParticleContainer_H_
10 
12 #include <AMReX_Vector.H>
13 
34 {
35 public:
36  RigidInjectedParticleContainer (amrex::AmrCore* amr_core,
37  int ispecies,
38  const std::string& name);
40 
41  virtual void InitData() override;
42 
43  virtual void RemapParticles();
44  virtual void BoostandRemapParticles();
45 
46  virtual void Evolve (int lev,
47  const amrex::MultiFab& Ex,
48  const amrex::MultiFab& Ey,
49  const amrex::MultiFab& Ez,
50  const amrex::MultiFab& Bx,
51  const amrex::MultiFab& By,
52  const amrex::MultiFab& Bz,
53  const amrex::MultiFab& Ex_avg,
54  const amrex::MultiFab& Ey_avg,
55  const amrex::MultiFab& Ez_avg,
56  const amrex::MultiFab& Bx_avg,
57  const amrex::MultiFab& By_avg,
58  const amrex::MultiFab& Bz_avg,
59  amrex::MultiFab& jx,
60  amrex::MultiFab& jy,
61  amrex::MultiFab& jz,
62  amrex::MultiFab* cjx,
63  amrex::MultiFab* cjy,
64  amrex::MultiFab* cjz,
65  amrex::MultiFab* rho,
66  amrex::MultiFab* crho,
67  const amrex::MultiFab* cEx,
68  const amrex::MultiFab* cEy,
69  const amrex::MultiFab* cEz,
70  const amrex::MultiFab* cBx,
71  const amrex::MultiFab* cBy,
72  const amrex::MultiFab* cBz,
73  amrex::Real t,
74  amrex::Real dt,
75  DtType a_dt_type=DtType::Full) override;
76 
77  virtual void PushPX (WarpXParIter& pti,
78  amrex::FArrayBox const * exfab,
79  amrex::FArrayBox const * eyfab,
80  amrex::FArrayBox const * ezfab,
81  amrex::FArrayBox const * bxfab,
82  amrex::FArrayBox const * byfab,
83  amrex::FArrayBox const * bzfab,
84  const int ngE, const int /*e_is_nodal*/,
85  const long offset,
86  const long np_to_push,
87  int lev, int gather_lev,
88  amrex::Real dt, ScaleFields scaleFields,
89  DtType a_dt_type=DtType::Full) override;
90 
91  virtual void PushP (int lev, amrex::Real dt,
92  const amrex::MultiFab& Ex,
93  const amrex::MultiFab& Ey,
94  const amrex::MultiFab& Ez,
95  const amrex::MultiFab& Bx,
96  const amrex::MultiFab& By,
97  const amrex::MultiFab& Bz) override;
98 
99  virtual void ReadHeader (std::istream& is) override;
100 
101  virtual void WriteHeader (std::ostream& os) const override;
102 
103 private:
104 
105  // User input quantities
106  amrex::Real zinject_plane = 0.;
107  bool projected = true; // When true, particle transverse positions are directly projected (without adjusment)
108  bool focused = false; // When true, particle transverse positions are adjusted to account for distance between zinject and z=0
109  bool rigid_advance = true; // When true, particles are advance with vzbar before injection
110 
111  amrex::Real vzbeam_ave_boosted;
112 
113  amrex::Vector<int> done_injecting;
114  amrex::Vector<amrex::Real> zinject_plane_levels;
115 
116  // Temporary quantites
117  amrex::Real zinject_plane_lev;
120 
121 };
122 
123 #endif
bool done_injecting_lev
Definition: RigidInjectedParticleContainer.H:119
DtType
Definition: WarpXDtType.H:10
RigidInjectedParticleContainer(amrex::AmrCore *amr_core, int ispecies, const std::string &name)
Definition: RigidInjectedParticleContainer.cpp:33
amrex::Real zinject_plane_lev_previous
Definition: RigidInjectedParticleContainer.H:118
Functor that scales E and B by a factor before pushing the particles. This is used for rigid injectio...
Definition: ScaleFields.H:13
Definition: PhysicalParticleContainer.H:39
Definition: RigidInjectedParticleContainer.H:32
virtual void WriteHeader(std::ostream &os) const override
Definition: ParticleIO.cpp:47
amrex::Real zinject_plane_lev
Definition: RigidInjectedParticleContainer.H:117
amrex::Vector< amrex::Real > zinject_plane_levels
Definition: RigidInjectedParticleContainer.H:114
bool focused
Definition: RigidInjectedParticleContainer.H:108
virtual void PushPX(WarpXParIter &pti, amrex::FArrayBox const *exfab, amrex::FArrayBox const *eyfab, amrex::FArrayBox const *ezfab, amrex::FArrayBox const *bxfab, amrex::FArrayBox const *byfab, amrex::FArrayBox const *bzfab, const int ngE, const int, const long offset, const long np_to_push, int lev, int gather_lev, amrex::Real dt, ScaleFields scaleFields, DtType a_dt_type=DtType::Full) override
Definition: RigidInjectedParticleContainer.cpp:229
bool rigid_advance
Definition: RigidInjectedParticleContainer.H:109
virtual void InitData() override
Definition: RigidInjectedParticleContainer.cpp:47
amrex::Real vzbeam_ave_boosted
Definition: RigidInjectedParticleContainer.H:111
virtual ~RigidInjectedParticleContainer()
Definition: RigidInjectedParticleContainer.H:39
virtual void RemapParticles()
Definition: RigidInjectedParticleContainer.cpp:61
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, 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) override
Evolve is the central function PhysicalParticleContainer that advances plasma particles for a time dt...
Definition: RigidInjectedParticleContainer.cpp:340
virtual void BoostandRemapParticles()
Definition: RigidInjectedParticleContainer.cpp:134
name
Definition: run_automated.py:204
bool projected
Definition: RigidInjectedParticleContainer.H:107
amrex::Vector< int > done_injecting
Definition: RigidInjectedParticleContainer.H:113
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) override
Definition: RigidInjectedParticleContainer.cpp:383
amrex::Real zinject_plane
Definition: RigidInjectedParticleContainer.H:106
Definition: WarpXParticleContainer.H:76
virtual void ReadHeader(std::istream &is) override
Definition: ParticleIO.cpp:17