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 
11 #include "Evolve/WarpXDtType.H"
14 
16 
17 #include <AMReX_REAL.H>
18 #include <AMReX_Vector.H>
19 
20 #include <AMReX_BaseFwd.H>
21 #include <AMReX_AmrCoreFwd.H>
22 
23 #include <iosfwd>
24 #include <string>
25 
46 {
47 public:
49  int ispecies,
50  const std::string& name);
51  ~RigidInjectedParticleContainer () override = default;
52 
57 
58 
59  void InitData() override;
60 
61  virtual void RemapParticles();
62 
63  void Evolve (int lev,
64  const amrex::MultiFab& Ex,
65  const amrex::MultiFab& Ey,
66  const amrex::MultiFab& Ez,
67  const amrex::MultiFab& Bx,
68  const amrex::MultiFab& By,
69  const amrex::MultiFab& Bz,
70  amrex::MultiFab& jx,
71  amrex::MultiFab& jy,
72  amrex::MultiFab& jz,
73  amrex::MultiFab* cjx,
74  amrex::MultiFab* cjy,
75  amrex::MultiFab* cjz,
76  amrex::MultiFab* rho,
77  amrex::MultiFab* crho,
78  const amrex::MultiFab* cEx,
79  const amrex::MultiFab* cEy,
80  const amrex::MultiFab* cEz,
81  const amrex::MultiFab* cBx,
82  const amrex::MultiFab* cBy,
83  const amrex::MultiFab* cBz,
84  amrex::Real t,
85  amrex::Real dt,
86  DtType a_dt_type=DtType::Full,
87  bool skip_deposition=false ) override;
88 
89  void PushPX (WarpXParIter& pti,
90  amrex::FArrayBox const * exfab,
91  amrex::FArrayBox const * eyfab,
92  amrex::FArrayBox const * ezfab,
93  amrex::FArrayBox const * bxfab,
94  amrex::FArrayBox const * byfab,
95  amrex::FArrayBox const * bzfab,
96  amrex::IntVect ngEB, int /*e_is_nodal*/,
97  long offset,
98  long np_to_push,
99  int lev, int gather_lev,
100  amrex::Real dt, ScaleFields scaleFields,
101  DtType a_dt_type=DtType::Full) override;
102 
103  void PushP (int lev, amrex::Real dt,
104  const amrex::MultiFab& Ex,
105  const amrex::MultiFab& Ey,
106  const amrex::MultiFab& Ez,
107  const amrex::MultiFab& Bx,
108  const amrex::MultiFab& By,
109  const amrex::MultiFab& Bz) override;
110 
111  void ReadHeader (std::istream& is) override;
112 
113  void WriteHeader (std::ostream& os) const override;
114 
115 private:
116 
117  // User input quantities
118  amrex::ParticleReal zinject_plane = 0.;
119  bool rigid_advance = true; // When true, particles are advance with vzbar before injection
120 
121  amrex::ParticleReal vzbeam_ave_boosted;
122 
124 
125  // Temporary quantites
126  amrex::ParticleReal zinject_plane_lev;
127  amrex::ParticleReal zinject_plane_lev_previous;
129 
130 };
131 
132 #endif
DtType
Definition: WarpXDtType.H:11
Definition: PhysicalParticleContainer.H:46
Definition: RigidInjectedParticleContainer.H:46
void WriteHeader(std::ostream &os) const override
Definition: ParticleIO.cpp:92
void InitData() override
Definition: RigidInjectedParticleContainer.cpp:72
amrex::ParticleReal zinject_plane_lev_previous
Definition: RigidInjectedParticleContainer.H:127
amrex::Vector< amrex::ParticleReal > zinject_plane_levels
Definition: RigidInjectedParticleContainer.H:123
virtual void RemapParticles()
Definition: RigidInjectedParticleContainer.cpp:89
bool done_injecting_lev
Definition: RigidInjectedParticleContainer.H:128
amrex::ParticleReal zinject_plane
Definition: RigidInjectedParticleContainer.H:118
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:331
~RigidInjectedParticleContainer() override=default
RigidInjectedParticleContainer & operator=(RigidInjectedParticleContainer const &)=delete
void ReadHeader(std::istream &is) override
Definition: ParticleIO.cpp:68
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) override
Evolve is the central function PhysicalParticleContainer that advances plasma particles for a time dt...
Definition: RigidInjectedParticleContainer.cpp:294
amrex::ParticleReal vzbeam_ave_boosted
Definition: RigidInjectedParticleContainer.H:121
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, amrex::IntVect ngEB, int, long offset, long np_to_push, int lev, int gather_lev, amrex::Real dt, ScaleFields scaleFields, DtType a_dt_type=DtType::Full) override
Definition: RigidInjectedParticleContainer.cpp:158
RigidInjectedParticleContainer(amrex::AmrCore *amr_core, int ispecies, const std::string &name)
Definition: RigidInjectedParticleContainer.cpp:59
RigidInjectedParticleContainer(RigidInjectedParticleContainer &&)=default
bool rigid_advance
Definition: RigidInjectedParticleContainer.H:119
amrex::ParticleReal zinject_plane_lev
Definition: RigidInjectedParticleContainer.H:126
RigidInjectedParticleContainer(RigidInjectedParticleContainer const &)=delete
Definition: WarpXParticleContainer.H:52
name
Definition: run_automated.py:229
float dt
Definition: stencil.py:442
Functor that scales E and B by a factor before pushing the particles. This is used for rigid injectio...
Definition: ScaleFields.H:14