Particle database class for two dimensions. More...
#include <particledatabase.hpp>
 
  
 | Public Member Functions | |
| ParticleDataBase2D (const Geometry &geom) | |
| Constructor.  More... | |
| ParticleDataBase2D (const ParticleDataBase2D &pdb) | |
| Copy constructor.  More... | |
| ParticleDataBase2D (std::istream &s, const Geometry &geom) | |
| Constructor for loading particle statistics from a file.  More... | |
| ~ParticleDataBase2D () | |
| Destructor.  More... | |
| const ParticleDataBase2D & | operator= (const ParticleDataBase2D &pdb) | 
| Assignment.  More... | |
| virtual Particle2D & | particle (uint32_t i) | 
| Returns a reference to particle i.  More... | |
| virtual const Particle2D & | particle (uint32_t i) const | 
| Returns a const reference to particle i.  More... | |
| virtual const ParticleP2D & | trajectory_point (uint32_t i, uint32_t j) const | 
| Gets the particle i trajectory point j as particle point.  More... | |
| void | trajectories_at_plane (std::vector< Particle2D > &tdata, coordinate_axis_e axis, double val) const | 
| Gets trajectory data at plane axis = val as particles in vector tdata.  More... | |
| void | add_particle (double IQ, double q, double m, const ParticleP2D &x) | 
| Add one particle.  More... | |
| void | add_particle (const Particle2D &p) | 
| Add one particle.  More... | |
| void | add_2d_beam_with_energy (uint32_t N, double J, double q, double m, double E, double Tp, double Tt, double x1, double y1, double x2, double y2) | 
| Add a 2d beam with energies.  More... | |
| void | add_2d_beam_with_velocity (uint32_t N, double J, double q, double m, double v, double dvp, double dvt, double x1, double y1, double x2, double y2) | 
| Add a 2d beam with velocities.  More... | |
| void | add_2d_KV_beam_with_emittance (uint32_t N, double I, double q, double m, double a, double b, double e, double Ex, double x0, double y0) | 
| Add a 2d beam with defined KV emittance.  More... | |
| void | add_2d_gaussian_beam_with_emittance (uint32_t N, double I, double q, double m, double a, double b, double e, double Ex, double x0, double y0) | 
| Add a 2d beam with defined gaussian emittance.  More... | |
| virtual void | save (const std::string &filename) const | 
| Saves data to a new file filename.  More... | |
| virtual void | save (std::ostream &s) const | 
| Saves data to stream.  More... | |
| virtual void | debug_print (std::ostream &os) const | 
| Print debugging information to os.  More... | |
|  Public Member Functions inherited from ParticleDataBase | |
| virtual | ~ParticleDataBase () | 
| Virtual destructor.  More... | |
| void | set_thread_count (uint32_t threadcount) | 
| Set the number of threads used for calculation.  More... | |
| void | set_accuracy (double epsabs, double epsrel) | 
| Set the accuracy requirement for calculation.  More... | |
| void | set_bfield_suppression (const CallbackFunctorD_V *functor) | 
| Set magnetic field suppression location depedent callback functor.  More... | |
| void | set_trajectory_handler_callback (TrajectoryHandlerCallback *thand_cb) | 
| Set trajectory handler callback.  More... | |
| void | set_trajectory_end_callback (TrajectoryEndCallback *tend_cb) | 
| Set trajectory end callback.  More... | |
| void | set_trajectory_surface_collision_callback (TrajectorySurfaceCollisionCallback *tsur_cb) | 
| Set trajectory surface collision callback.  More... | |
| void | set_relativistic (bool enable) | 
| Set relativistic particle iteration.  More... | |
| void | set_surface_collision (bool surface_collision) | 
| Set surface collision model to be used.  More... | |
| void | set_polyint (bool polyint) | 
| Set the interpolation type to polynomial(true) or linear(false).  More... | |
| bool | get_polyint (void) const | 
| Get current interpolation type.  More... | |
| void | set_trajectory_interpolation (trajectory_interpolation_e intrp) | 
| Set trajectory interpolation type.  More... | |
| trajectory_interpolation_e | get_trajectory_interpolation (void) const | 
| Get trajectory interpolation type.  More... | |
| void | set_scharge_deposition (scharge_deposition_e type) | 
| Set space charge deposition type.  More... | |
| scharge_deposition_e | get_scharge_deposition (void) const | 
| Get space charge deposition type.  More... | |
| void | set_max_steps (uint32_t maxsteps) | 
| Set maximum number of steps to iterate.  More... | |
| void | set_max_time (double maxt) | 
| Set maximum lifetime of particle in simulation.  More... | |
| void | set_save_all_points (bool save_points) | 
| Set trajectory saving.  More... | |
| void | set_save_trajectories (uint32_t div) | 
| Set trajectory saving.  More... | |
| uint32_t | get_save_trajectories (void) const | 
| Get trajectory saving.  More... | |
| void | set_mirror (const bool mirror[6]) | 
| Set particle mirroring on boundaries.  More... | |
| void | get_mirror (bool mirror[6]) const | 
| Get particle mirroring on boundaries.  More... | |
| int | get_iteration_number (void) const | 
| Get the number of iteration rounds done.  More... | |
| double | get_rhosum (void) const | 
| Return sum of defined beam space charge density.  More... | |
| void | set_rhosum (double rhosum) | 
| Set sum of defined beam space charge density.  More... | |
| const ParticleStatistics & | get_statistics (void) const | 
| Get particle iterator statistics.  More... | |
| geom_mode_e | geom_mode () const | 
| Return geometry mode.  More... | |
| size_t | size (void) const | 
| Returns particle count.  More... | |
| double | traj_length (uint32_t i) const | 
| Returns the length of trajectory i.  More... | |
| size_t | traj_size (uint32_t i) const | 
| Returns number of trajectory points for particle i.  More... | |
| void | trajectory_point (double &t, Vec3D &loc, Vec3D &vel, uint32_t i, uint32_t j) const | 
| Gets the particle i trajectory point j into vel, loc and t.  More... | |
| void | trajectories_at_plane (TrajectoryDiagnosticData &tdata, coordinate_axis_e axis, double val, const std::vector< trajectory_diagnostic_e > &diagnostics) const | 
| Gets trajectory diagnostic diagnostics at plane axis = val in trajectory diagnostic data object tdata.  More... | |
| void | build_trajectory_density_field (MeshScalarField &tdens) const | 
| Build trajectory density field.  More... | |
| void | clear (void) | 
| Clears the particle database of all particles.  More... | |
| void | clear_trajectories (void) | 
| Clears the particle trajectory database.  More... | |
| void | clear_trajectory (size_t a) | 
| Clears particle a in the particle trajectory database.  More... | |
| void | reserve (size_t size) | 
| Reserve memory for size particles.  More... | |
| void | iterate_trajectories (MeshScalarField &scharge, const VectorField &efield, const VectorField &bfield) | 
| Iterate particles through the geometry.  More... | |
| void | step_particles (MeshScalarField &scharge, const VectorField &efield, const VectorField &bfield, double dt) | 
| Step particles forward by time step dt.  More... | |
| Additional Inherited Members | |
|  Protected Member Functions inherited from ParticleDataBase | |
| ParticleDataBase () | |
| Constructor.  More... | |
| ParticleDataBase (const ParticleDataBase &pdb) | |
| Copy constructor.  More... | |
| const ParticleDataBase & | operator= (const ParticleDataBase &pdb) | 
| Assignment.  More... | |
| void | set_implementation_pointer (class ParticleDataBaseImp *imp) | 
| Set particle database implementation pointer.  More... | |
Detailed Description
Particle database class for two dimensions.
Particle database holds the definitions of particles and possibly the trajectories of particles it the particle iterator has saved them. ParticleDataBase2D provides a variety of convenience functions for defining particle beams.
Particles are always stored in the database in the order they are defined. When reading back the simulation results, the order can be used to identify the particles.
Constructor & Destructor Documentation
| ParticleDataBase2D::ParticleDataBase2D | ( | const Geometry & | geom | ) | 
Constructor.
| ParticleDataBase2D::ParticleDataBase2D | ( | const ParticleDataBase2D & | pdb | ) | 
Copy constructor.
| ParticleDataBase2D::ParticleDataBase2D | ( | std::istream & | s, | 
| const Geometry & | geom | ||
| ) | 
Constructor for loading particle statistics from a file.
| ParticleDataBase2D::~ParticleDataBase2D | ( | ) | 
Destructor.
Member Function Documentation
| void ParticleDataBase2D::add_2d_beam_with_energy | ( | uint32_t | N, | 
| double | J, | ||
| double | q, | ||
| double | m, | ||
| double | E, | ||
| double | Tp, | ||
| double | Tt, | ||
| double | x1, | ||
| double | y1, | ||
| double | x2, | ||
| double | y2 | ||
| ) | 
Add a 2d beam with energies.
Adds a beam consisting of N particles. The beam current density is J (A/m^2), charge of beam particle is q (in multiples of e), mass is m (u). The beam is defined on a line from (x1, y1) to (x2, y2). The beam propagates into a direction 90 degrees clockwise from the direction of vector pointing from (x1, y1) to (x2, y2) with a mean energy E (eV). The beam has parallel temperature Tp (eV) and transverse temperature Tt (eV).
The particle speeds of the beam in direction i are sampled from a gaussian distribution with standard deviation dv_i = sqrt(T_i*e/m), where T_i is the beam temperature in direction (eV), e is electron charge (C) and m is the mass of the ion (kg).
Space charge J/v is added to the rhosum variable.
| void ParticleDataBase2D::add_2d_beam_with_velocity | ( | uint32_t | N, | 
| double | J, | ||
| double | q, | ||
| double | m, | ||
| double | v, | ||
| double | dvp, | ||
| double | dvt, | ||
| double | x1, | ||
| double | y1, | ||
| double | x2, | ||
| double | y2 | ||
| ) | 
Add a 2d beam with velocities.
Adds a beam consisting of N particles. The beam current density is J (A/m^2), charge of beam particle is q (in multiples of e), mass is m (u). The beam is defined on a line from (x1, y1) to (x2, y2). The beam propagates into a direction 90 degrees clockwise from the direction of vector pointing from (x1, y1) to (x2, y2) with a mean velocity v (m/s). The beam has parallel gaussian velocity distribution with standard deviation dvp (m/s) and transverse gaussian velocity distribution with standard deviation dvt (m/s).
Space charge J/v is added to the rhosum variable.
| void ParticleDataBase2D::add_2d_gaussian_beam_with_emittance | ( | uint32_t | N, | 
| double | I, | ||
| double | q, | ||
| double | m, | ||
| double | a, | ||
| double | b, | ||
| double | e, | ||
| double | Ex, | ||
| double | x0, | ||
| double | y0 | ||
| ) | 
Add a 2d beam with defined gaussian emittance.
Adds a beam consisting of N particles and a total current of I (A). The particles are defined to have equal currents, charge q (in multiples of e) and mass m (u). The starting energy of the beam is Ex (eV) and the starting location (x0,y0) (center point). The beam propagates to positive x-direction. The beam is made to match Twiss parameters  (a) and
 (a) and  (b) in projectional direction (y,y'). The rms-emittance of the beam is made to match
 (b) in projectional direction (y,y'). The rms-emittance of the beam is made to match  (e).
 (e).
The beam spread in the projectional space is made according to Gaussian distribution.
| void ParticleDataBase2D::add_2d_KV_beam_with_emittance | ( | uint32_t | N, | 
| double | I, | ||
| double | q, | ||
| double | m, | ||
| double | a, | ||
| double | b, | ||
| double | e, | ||
| double | Ex, | ||
| double | x0, | ||
| double | y0 | ||
| ) | 
Add a 2d beam with defined KV emittance.
Adds a beam consisting of N particles and a total current of I (A). The particles are defined to have equal currents, charge q (in multiples of e) and mass m (u). The starting energy of the beam is Ex (eV) and the starting location (x0,y0) (center point). The beam propagates to positive x-direction. The beam is made to match Twiss parameters  (a) and
 (a) and  (b) in projectional direction (y,y'). The rms-emittance of the beam is made to match
 (b) in projectional direction (y,y'). The rms-emittance of the beam is made to match  (e).
 (e).
The beam spread in the projectional space is made according to KV/hard-edged (Kapchinsky-Vladimirsky) distribution.
| void ParticleDataBase2D::add_particle | ( | double | IQ, | 
| double | q, | ||
| double | m, | ||
| const ParticleP2D & | x | ||
| ) | 
Add one particle.
Adds one particle to particle database. Particle properties are: IQ is the current (A) in time-independent or charge (C) in time-dependent simulations carried by the particle cloud that the simulated particle represents, q is the charge state of the microscopic particle (in multiples of e), m is the mass of the microscopic particle (u) and x contains the time, position (m) and velocity (m/s) of the particle.
| void ParticleDataBase2D::add_particle | ( | const Particle2D & | p | ) | 
Add one particle.
Adds one particle to database.
| 
 | virtual | 
Print debugging information to os.
Implements ParticleDataBase.
| const ParticleDataBase2D & ParticleDataBase2D::operator= | ( | const ParticleDataBase2D & | pdb | ) | 
Assignment.
| 
 | virtual | 
Returns a reference to particle i.
Implements ParticleDataBase.
| 
 | virtual | 
Returns a const reference to particle i.
Implements ParticleDataBase.
| 
 | virtual | 
Saves data to a new file filename.
Implements ParticleDataBase.
| 
 | virtual | 
Saves data to stream.
Implements ParticleDataBase.
| void ParticleDataBase2D::trajectories_at_plane | ( | std::vector< Particle2D > & | tdata, | 
| coordinate_axis_e | axis, | ||
| double | val | ||
| ) | const | 
Gets trajectory data at plane axis = val as particles in vector tdata.
| 
 | virtual | 
Gets the particle i trajectory point j as particle point.
Implements ParticleDataBase.
The documentation for this class was generated from the following files:
