OPAL (Object Oriented Parallel Accelerator Library)
2.2.0
OPAL
|
Particle Bunch. More...
#include <PartBunch.h>
Public Types | |
enum | { Dim = Dimension } |
typedef IpplParticleBase < Layout_t > | pbase_t |
![]() | |
enum | UnitState_t |
typedef AbstractParticle < double, Dim >::ParticlePos_t | ParticlePos_t |
typedef AbstractParticle < double, Dim > ::ParticleIndex_t | ParticleIndex_t |
typedef AbstractParticle < double, Dim >::UpdateFlags | UpdateFlags_t |
typedef AbstractParticle < double, Dim >::Position_t | Position_t |
typedef std::pair< Vector_t, Vector_t > | VectorPair_t |
Public Member Functions | |
PartBunch (const PartData *ref) | |
Default constructor. More... | |
PartBunch ()=delete | |
PartBunch (const PartBunch &)=delete | |
PartBunch & | operator= (const PartBunch &)=delete |
~PartBunch () | |
void | runTests () |
void | initialize (FieldLayout_t *fLayout) |
void | do_binaryRepart () |
double | getRho (int x, int y, int z) |
const Mesh_t & | getMesh () const |
Mesh_t & | getMesh () |
FieldLayout_t & | getFieldLayout () |
void | setBCAllPeriodic () |
void | setBCAllOpen () |
void | setBCForDCBeam () |
VectorPair_t | getEExtrema () |
void | computeSelfFields () |
void | computeSelfFields (int b) |
void | computeSelfFields_cycl (double gamma) |
Calculates the self electric field from the charge density distribution for use in cyclotrons. More... | |
void | computeSelfFields_cycl (int b) |
Calculates the self electric field from the charge density distribution for use in cyclotrons. More... | |
void | resetInterpolationCache (bool clearCache=false) |
void | swap (unsigned int i, unsigned int j) |
Inform & | print (Inform &os) |
![]() | |
PartBunchBase (AbstractParticle< double, Dim > *pb) | |
PartBunchBase (AbstractParticle< double, Dim > *pb, const PartData *ref) | |
PartBunchBase (const PartBunchBase &rhs)=delete | |
virtual | ~PartBunchBase () |
bool | getIfBeamEmitting () |
int | getLastEmittedEnergyBin () |
size_t | getNumberOfEmissionSteps () |
int | getNumberOfEnergyBins () |
void | Rebin () |
void | setEnergyBins (int numberOfEnergyBins) |
bool | weHaveEnergyBins () |
void | switchToUnitlessPositions (bool use_dt_per_particle=false) |
void | switchOffUnitlessPositions (bool use_dt_per_particle=false) |
void | setDistribution (Distribution *d, std::vector< Distribution * > addedDistributions, size_t &np) |
bool | isGridFixed () |
bool | hasBinning () |
void | setTEmission (double t) |
double | getTEmission () |
bool | doEmission () |
bool | weHaveBins () const |
void | setPBins (PartBins *pbin) |
void | setPBins (PartBinsCyc *pbin) |
size_t | emitParticles (double eZ) |
Emit particles in the given bin i.e. copy the particles from the bin structure into the particle container. More... | |
void | updateNumTotal () |
void | rebin () |
int | getNumBins () |
int | getLastemittedBin () |
void | setLocalBinCount (size_t num, int bin) |
void | calcGammas () |
Compute the gammas of all bins. More... | |
void | calcGammas_cycl () |
double | getBinGamma (int bin) |
Get gamma of one bin. More... | |
virtual void | setBinCharge (int bin, double q) |
Set the charge of one bin to the value of q and all other to zero. More... | |
virtual void | setBinCharge (int bin) |
Set the charge of all other the ones in bin to zero. More... | |
size_t | calcNumPartsOutside (Vector_t x) |
returns the number of particles outside of a box defined by x More... | |
void | calcLineDensity (unsigned int nBins, std::vector< double > &lineDensity, std::pair< double, double > &meshInfo) |
void | setBeamFrequency (double v) |
virtual void | boundp () |
void | boundp_destroy () |
size_t | boundp_destroyT () |
size_t | destroyT () |
virtual double | getPx (int i) |
virtual double | getPy (int i) |
virtual double | getPz (int i) |
virtual double | getPx0 (int i) |
virtual double | getPy0 (int i) |
virtual double | getX (int i) |
virtual double | getY (int i) |
virtual double | getZ (int i) |
virtual double | getX0 (int i) |
virtual double | getY0 (int i) |
virtual void | setZ (int i, double zcoo) |
void | get_bounds (Vector_t &rmin, Vector_t &rmax) |
void | getLocalBounds (Vector_t &rmin, Vector_t &rmax) |
std::pair< Vector_t, double > | getBoundingSphere () |
std::pair< Vector_t, double > | getLocalBoundingSphere () |
void | push_back (OpalParticle p) |
void | set_part (FVector< double, 6 > z, int ii) |
void | set_part (OpalParticle p, int ii) |
OpalParticle | get_part (int ii) |
void | maximumAmplitudes (const FMatrix< double, 6, 6 > &D, double &axmax, double &aymax) |
Return maximum amplitudes. More... | |
void | setdT (double dt) |
double | getdT () const |
void | setT (double t) |
void | incrementT () |
double | getT () const |
double | get_sPos () |
void | set_sPos (double s) |
double | get_gamma () const |
double | get_meanKineticEnergy () const |
Vector_t | get_origin () const |
Vector_t | get_maxExtent () const |
Vector_t | get_centroid () const |
Vector_t | get_rrms () const |
Vector_t | get_rprms () const |
Vector_t | get_rmean () const |
Vector_t | get_prms () const |
Vector_t | get_pmean () const |
Vector_t | get_pmean_Distribution () const |
Vector_t | get_emit () const |
Vector_t | get_norm_emit () const |
Vector_t | get_halo () const |
virtual Vector_t | get_hr () const |
double | get_Dx () const |
double | get_Dy () const |
double | get_DDx () const |
double | get_DDy () const |
virtual void | set_meshEnlargement (double dh) |
void | gatherLoadBalanceStatistics () |
size_t | getLoadBalance (int p) const |
void | get_PBounds (Vector_t &min, Vector_t &max) const |
void | calcBeamParameters () |
void | calcBeamParametersInitial () |
double | getCouplingConstant () const |
void | setCouplingConstant (double c) |
void | setCharge (double q) |
void | setChargeZeroPart (double q) |
void | setMass (double mass) |
double | getEkin () const |
Need Ek for the Schottky effect calculation (eV) More... | |
double | getWorkFunctionRf () const |
Need the work function for the Schottky effect calculation (eV) More... | |
double | getLaserEnergy () const |
Need the laser energy for the Schottky effect calculation (eV) More... | |
double | getCharge () const |
get the total charge per simulation particle More... | |
double | getChargePerParticle () const |
get the macro particle charge More... | |
virtual void | setSolver (FieldSolver *fs) |
bool | hasFieldSolver () |
std::string | getFieldSolverType () const |
void | setStepsPerTurn (int n) |
int | getStepsPerTurn () const |
void | setGlobalTrackStep (long long n) |
step in multiple TRACK commands More... | |
long long | getGlobalTrackStep () const |
void | setLocalTrackStep (long long n) |
step in a TRACK command More... | |
void | incTrackSteps () |
long long | getLocalTrackStep () const |
void | setNumBunch (short n) |
short | getNumBunch () const |
void | setTotalNumPerBunch (size_t numpart, short n) |
size_t | getTotalNumPerBunch (short n) const |
void | setLocalNumPerBunch (size_t numpart, short n) |
size_t | getLocalNumPerBunch (short n) const |
void | countTotalNumPerBunch () |
void | setGlobalMeanR (Vector_t globalMeanR) |
Vector_t | getGlobalMeanR () |
void | setGlobalToLocalQuaternion (Quaternion_t globalToLocalQuaternion) |
Quaternion_t | getGlobalToLocalQuaternion () |
void | setSteptoLastInj (int n) |
int | getSteptoLastInj () |
double | calcMeanPhi () |
calculate average angle of longitudinal direction of bins More... | |
bool | resetPartBinID2 (const double eta) |
reset Bin[] for each particle according to the method given in paper PAST-AB(064402) by G. Fubiani et al. More... | |
bool | resetPartBinBunch () |
double | getdE () |
virtual double | getGamma (int i) |
virtual double | getBeta (int i) |
virtual void | actT () |
const PartData * | getReference () const |
double | getEmissionDeltaT () |
Quaternion_t | getQKs3D () |
void | setQKs3D (Quaternion_t q) |
Vector_t | getKs3DRefr () |
void | setKs3DRefr (Vector_t r) |
Vector_t | getKs3DRefp () |
void | setKs3DRefp (Vector_t p) |
void | iterateEmittedBin (int binNumber) |
void | calcEMean () |
void | correctEnergy (double avrgp) |
Inform & | print (Inform &os) |
size_t | getTotalNum () const |
void | setTotalNum (size_t n) |
void | setLocalNum (size_t n) |
size_t | getLocalNum () const |
size_t | getDestroyNum () const |
size_t | getGhostNum () const |
unsigned int | getMinimumNumberOfParticlesPerCore () const |
void | setMinimumNumberOfParticlesPerCore (unsigned int n) |
ParticleLayout< double, Dim > & | getLayout () |
const ParticleLayout< double, Dim > & | getLayout () const |
bool | getUpdateFlag (UpdateFlags_t f) const |
void | setUpdateFlag (UpdateFlags_t f, bool val) |
ParticleBConds< Position_t, Dimension > & | getBConds () |
void | setBConds (const ParticleBConds< Position_t, Dimension > &bc) |
bool | singleInitNode () const |
void | resetID () |
void | update () |
void | update (const ParticleAttrib< char > &canSwap) |
void | createWithID (unsigned id) |
void | create (size_t M) |
void | globalCreate (size_t np) |
void | destroy (size_t M, size_t I, bool doNow=false) |
void | performDestroy (bool updateLocalNum=false) |
void | ghostDestroy (size_t M, size_t I) |
FMatrix< double, 2 *Dim, 2 *Dim > | getSigmaMatrix () |
double | getQ () const |
Access to reference data. More... | |
double | getM () const |
double | getP () const |
double | getE () const |
ParticleType::type | getPType () const |
double | getInitialBeta () const |
double | getInitialGamma () const |
void | resetQ (double q) |
void | resetM (double m) |
void | setPType (ParticleType::type) |
Private Member Functions | |
void | updateDomainLength (Vektor< int, 3 > &grid) |
void | updateFields (const Vector_t &hr, const Vector_t &origin) |
void | resizeMesh () |
resize mesh to geometry specified More... | |
ParticleLayout< double, 3 > & | getLayout () |
const ParticleLayout< double, 3 > & | getLayout () const |
Private Attributes | |
BConds< double, 3, Mesh_t, Center_t > | bc_m |
for defining the boundary conditions More... | |
BConds< Vector_t, 3, Mesh_t, Center_t > | vbc_m |
bool | interpolationCacheSet_m |
ParticleAttrib< CacheDataCIC < double, 3U > > | interpolationCache_m |
Additional Inherited Members | |
![]() | |
static const unsigned | Dimension |
![]() | |
size_t | calcMoments () |
void | calcMomentsInitial () |
double | calculateAngle (double x, double y) |
angle range [0~2PI) degree More... | |
![]() | |
IpplTimings::TimerRef | boundpTimer_m |
IpplTimings::TimerRef | boundpBoundsTimer_m |
IpplTimings::TimerRef | boundpUpdateTimer_m |
IpplTimings::TimerRef | statParamTimer_m |
IpplTimings::TimerRef | histoTimer_m |
const PartData * | reference |
UnitState_t | unit_state_ |
UnitState_t | stateOfLastBoundP_ |
double | centroid_m [2 *Dim] |
holds the centroid of the beam More... | |
FMatrix< double, 2 *Dim, 2 *Dim > | moments_m |
6x6 matrix of the moments of the beam More... | |
double | dt_m |
holds the timestep in seconds More... | |
double | t_m |
holds the actual time of the integration More... | |
double | eKin_m |
mean energy of the bunch (MeV) More... | |
double | dE_m |
energy spread of the beam in MeV More... | |
double | spos_m |
the position along design trajectory More... | |
Vector_t | globalMeanR_m |
Quaternion_t | globalToLocalQuaternion_m |
Vector_t | rmax_m |
maximal extend of particles More... | |
Vector_t | rmin_m |
minimal extend of particles More... | |
Vector_t | rrms_m |
rms beam size (m) More... | |
Vector_t | prms_m |
rms momenta More... | |
Vector_t | rmean_m |
mean position (m) More... | |
Vector_t | pmean_m |
mean momenta More... | |
Vector_t | eps_m |
rms emittance (not normalized) More... | |
Vector_t | eps_norm_m |
rms normalized emittance More... | |
Vector_t | halo_m |
Vector_t | rprms_m |
rms correlation More... | |
double | Dx_m |
dispersion x & y More... | |
double | Dy_m |
double | DDx_m |
derivative of the dispersion More... | |
double | DDy_m |
Vector_t | hr_m |
meshspacing of cartesian mesh More... | |
Vektor< int, 3 > | nr_m |
meshsize of cartesian mesh More... | |
FieldSolver * | fs_m |
stores the used field solver More... | |
double | couplingConstant_m |
double | qi_m |
int | distDump_m |
counter to store the distribution dump More... | |
int | fieldDBGStep_m |
double | dh_m |
Mesh enlargement. More... | |
double | tEmission_m |
in % how much the mesh is enlarged More... | |
std::unique_ptr< double[]> | bingamma_m |
holds the gamma of the bin More... | |
std::unique_ptr< size_t[]> | binemitted_m |
int | stepsPerTurn_m |
steps per turn for OPAL-cycl More... | |
long long | localTrackStep_m |
step in a TRACK command More... | |
long long | globalTrackStep_m |
if multiple TRACK commands More... | |
short | numBunch_m |
current bunch number More... | |
std::vector< size_t > | bunchTotalNum_m |
number of particles per bunch More... | |
std::vector< size_t > | bunchLocalNum_m |
int | SteptoLastInj_m |
std::unique_ptr< size_t[]> | globalPartPerNode_m |
Distribution * | dist_m |
bool | dcBeam_m |
double | periodLength_m |
std::shared_ptr < AbstractParticle< double, Dim > > | pbase |
Particle Bunch.
Definition at line 30 of file PartBunch.h.
typedef IpplParticleBase<Layout_t> PartBunch::pbase_t |
Definition at line 33 of file PartBunch.h.
anonymous enum |
Enumerator | |
---|---|
Dim |
Definition at line 34 of file PartBunch.h.
|
explicit |
Default constructor.
Definition at line 54 of file PartBunch.cpp.
|
delete |
|
delete |
PartBunch::~PartBunch | ( | ) |
Definition at line 62 of file PartBunch.cpp.
|
virtual |
Magnetic field in x and y direction induced by the eletric field
\[ B_x = \gamma(B_x^{'} - \frac{beta}{c}E_y^{'}) = -\gamma \frac{beta}{c}E_y^{'} = -\frac{beta}{c}E_y \]
\[ B_y = \gamma(B_y^{'} - \frac{beta}{c}E_x^{'}) = +\gamma \frac{beta}{c}E_x^{'} = +\frac{beta}{c}E_x \]
\[ B_z = B_z^{'} = 0 \]
Implements PartBunchBase< double, 3 >.
Definition at line 463 of file PartBunch.cpp.
References PartBunchBase< double, 3 >::Bf, Physics::c, PoissonSolver::computePotential(), PartBunchBase< double, 3 >::dt, PartBunchBase< double, 3 >::Ef, eg_m, endl(), PartBunchBase< double, 3 >::fieldDBGStep_m, PartBunchBase< double, 3 >::fs_m, PartBunchBase< double, 3 >::getCouplingConstant(), PartBunchBase< double, 3 >::getdT(), getFieldLayout(), FieldSolver::getFieldSolverType(), OpalData::getInputBasename(), OpalData::getInstance(), FieldLayout< Dim >::getLocalNDIndex(), PartBunchBase< double, 3 >::getTotalNum(), gmsg, Grad(), FieldSolver::hasValidSolver(), PartBunchBase< double, 3 >::hr_m, INFOMSG, PartBunchBase< double, 3 >::nr_m, PartBunchBase< double, 3 >::P, PartBunchBase< double, 3 >::Q, PartBunchBase< double, 3 >::R, resizeMesh(), rho_m, ParticleAttrib< T >::scatter(), PartBunchBase< double, 3 >::selfFieldTimer_m, FieldSolver::solver_m, sqrt(), IpplTimings::startTimer(), IpplTimings::stopTimer(), and sum().
|
virtual |
/brief used for self fields with binned distribution
Set initial charge density to zero. Create image charge potential field.
Set initial E field to zero.
Mesh the whole domain
Scatter charge onto space charge grid.
Calculate mesh-scale factor and get gamma for this specific slice (bin).
Scale charge density to get charge density in real units. Account for Lorentz transformation in longitudinal direction.
Scale mesh spacing to real units (meters). Lorentz transform the longitudinal direction.
Find potential from charge in this bin (no image yet) using Poisson's equation (without coefficient: -1/(eps)).
Scale mesh back (to same units as particle locations.)
The scalar potential is given back in rho_m and must be converted to the right units.
IPPL Grad numerical computes gradient to find the electric field (in bin rest frame).
Scale field. Combine Lorentz transform with conversion to proper field units.
Interpolate electric field at particle positions. We reuse the cached information about where the particles are relative to the field, since the particles have not moved since this the most recent scatter operation.
Magnetic field in x and y direction induced by the electric field.
\[ B_x = \gamma(B_x^{'} - \frac{beta}{c}E_y^{'}) = -\gamma \frac{beta}{c}E_y^{'} = -\frac{beta}{c}E_y \]
\[ B_y = \gamma(B_y^{'} - \frac{beta}{c}E_x^{'}) = +\gamma \frac{beta}{c}E_x^{'} = +\frac{beta}{c}E_x \]
\[ B_z = B_z^{'} = 0 \]
Now compute field due to image charge. This is done separately as the image charge is moving to -infinity (instead of +infinity) so the Lorentz transform is different.
Find z shift for shifted Green's function.
Find potential from image charge in this bin using Poisson's equation (without coefficient: -1/(eps)).
Scale mesh back (to same units as particle locations.)
The scalar potential is given back in rho_m and must be converted to the right units.
IPPL Grad numerical computes gradient to find the electric field (in rest frame of this bin's image charge).
Scale field. Combine Lorentz transform with conversion to proper field units.
Interpolate electric field at particle positions. We reuse the cached information about where the particles are relative to the field, since the particles have not moved since this the most recent scatter operation.
Magnetic field in x and y direction induced by the image charge electric field. Note that beta will have the opposite sign from the bunch charge field, as the image charge is moving in the opposite direction.
\[ B_x = \gamma(B_x^{'} - \frac{beta}{c}E_y^{'}) = -\gamma \frac{beta}{c}E_y^{'} = -\frac{beta}{c}E_y \]
\[ B_y = \gamma(B_y^{'} - \frac{beta}{c}E_x^{'}) = +\gamma \frac{beta}{c}E_x^{'} = +\frac{beta}{c}E_x \]
\[ B_z = B_z^{'} = 0 \]
Implements PartBunchBase< double, 3 >.
Definition at line 119 of file PartBunch.cpp.
References PartBunchBase< double, 3 >::Bf, Physics::c, PoissonSolver::computePotential(), ParticleAttrib< T >::create(), ParticleAttrib< T >::destroy(), PartBunchBase< double, 3 >::dt, PartBunchBase< double, 3 >::Ef, PartBunchBase< double, 3 >::Eftmp, eg_m, endl(), PartBunchBase< double, 3 >::fieldDBGStep_m, PartBunchBase< double, 3 >::fs_m, PartBunchBase< double, 3 >::get_bounds(), Field< T, Dim, Mesh, Centering >::get_mesh(), PartBunchBase< double, 3 >::getBinGamma(), PartBunchBase< double, 3 >::getCouplingConstant(), FieldLayout< Dim >::getDomain(), PartBunchBase< double, 3 >::getdT(), getFieldLayout(), FieldSolver::getFieldSolverType(), OpalData::getInputBasename(), OpalData::getInstance(), FieldLayout< Dim >::getLocalNDIndex(), PartBunchBase< double, 3 >::getLocalNum(), IpplInfo::getNodes(), gmsg, Grad(), FieldSolver::hasValidSolver(), PartBunchBase< double, 3 >::hr_m, INFOMSG, interpolationCache_m, interpolationCacheSet_m, level1(), PartBunchBase< double, 3 >::nr_m, PartBunchBase< double, 3 >::Q, PartBunchBase< double, 3 >::R, resizeMesh(), rho_m, ParticleAttrib< T >::scatter(), PartBunchBase< double, 3 >::selfFieldTimer_m, ParticleAttrib< T >::size(), FieldSolver::solver_m, sqrt(), IpplTimings::startTimer(), and IpplTimings::stopTimer().
|
virtual |
Calculates the self electric field from the charge density distribution for use in cyclotrons.
Comments -DW: I have made some changes in here: -) Some refacturing to make more similar to computeSelfFields() -) Added meanR and quaternion to be handed to the function so that SAAMG solver knows how to rotate the boundary geometry correctly. -) Fixed an error where gamma was not taken into account correctly in direction of movement (y in cyclotron) -) Comment: There is no account for image charges in the cyclotron tracker (yet?)!
set initial charge density to zero.
set initial E field to zero
mesh the whole domain
scatter particles charge onto grid.
Lorentz transformation In particle rest frame, the longitudinal length (y for cyclotron) enlarged
from charge (C) to charge density (C/m^3).
now charge density is in rho_m calculate Possion equation (without coefficient: -1/(eps))
retrive coefficient: -1/(eps)
calculate electric field vectors from field potential
Back Lorentz transformation CAVE: y coordinate needs 1/gamma factor because IPPL function Grad() divides by hr_m which is not scaled appropriately with Lorentz contraction in y direction only hr_scaled is! -DW
interpolate electric field at particle positions.
calculate coefficient
calculate B field from E field
Implements PartBunchBase< double, 3 >.
Definition at line 643 of file PartBunch.cpp.
References PartBunchBase< double, 3 >::Bf, Physics::c, PoissonSolver::computePotential(), PartBunchBase< double, 3 >::Ef, eg_m, endl(), PartBunchBase< double, 3 >::fieldDBGStep_m, PartBunchBase< double, 3 >::fs_m, PartBunchBase< double, 3 >::getCouplingConstant(), getFieldLayout(), FieldSolver::getFieldSolverType(), OpalData::getInputBasename(), OpalData::getInstance(), FieldLayout< Dim >::getLocalNDIndex(), gmsg, Grad(), FieldSolver::hasValidSolver(), PartBunchBase< double, 3 >::hr_m, INFOMSG, PartBunchBase< double, 3 >::nr_m, PartBunchBase< double, 3 >::Q, PartBunchBase< double, 3 >::R, resizeMesh(), rho_m, ParticleAttrib< T >::scatter(), PartBunchBase< double, 3 >::selfFieldTimer_m, FieldSolver::solver_m, sqrt(), IpplTimings::startTimer(), and IpplTimings::stopTimer().
|
virtual |
Calculates the self electric field from the charge density distribution for use in cyclotrons.
Overloaded version for having multiple bins with separate gamma for each bin. This is necessary For multi-bunch mode.
Comments -DW: I have made some changes in here: -) Some refacturing to make more similar to computeSelfFields() -) Added meanR and quaternion to be handed to the function (TODO: fall back to meanR = 0 and unit quaternion if not specified) so that SAAMG solver knows how to rotate the boundary geometry correctly. -) Fixed an error where gamma was not taken into account correctly in direction of movement (y in cyclotron) -) Comment: There is no account for image charges in the cyclotron tracker (yet?)!
set initial charge dentsity to zero.
set initial E field to zero
get gamma of this bin
mesh the whole domain
scatter particles charge onto grid.
Lorentz transformation In particle rest frame, the longitudinal length (y for cyclotron) enlarged
from charge (C) to charge density (C/m^3).
now charge density is in rho_m calculate Possion equation (without coefficient: -1/(eps))
retrive coefficient: -1/(eps)
calculate electric field vectors from field potential
Back Lorentz transformation CAVE: y coordinate needs 1/gamma factor because IPPL function Grad() divides by hr_m which is not scaled appropriately with Lorentz contraction in y direction only hr_scaled is! -DW
Interpolate electric field at particle positions.
Calculate coefficient
Calculate B_bin field from E_bin field accumulate B and E field
Implements PartBunchBase< double, 3 >.
Definition at line 833 of file PartBunch.cpp.
References PartBunchBase< double, 3 >::Bf, Physics::c, PoissonSolver::computePotential(), PartBunchBase< double, 3 >::Ef, PartBunchBase< double, 3 >::Eftmp, eg_m, endl(), PartBunchBase< double, 3 >::fieldDBGStep_m, PartBunchBase< double, 3 >::fs_m, PartBunchBase< double, 3 >::getBinGamma(), PartBunchBase< double, 3 >::getCouplingConstant(), getFieldLayout(), FieldSolver::getFieldSolverType(), OpalData::getInputBasename(), OpalData::getInstance(), FieldLayout< Dim >::getLocalNDIndex(), gmsg, Grad(), FieldSolver::hasValidSolver(), PartBunchBase< double, 3 >::hr_m, INFOMSG, PartBunchBase< double, 3 >::nr_m, PartBunchBase< double, 3 >::Q, PartBunchBase< double, 3 >::R, resizeMesh(), rho_m, ParticleAttrib< T >::scatter(), PartBunchBase< double, 3 >::selfFieldTimer_m, FieldSolver::solver_m, sqrt(), IpplTimings::startTimer(), and IpplTimings::stopTimer().
|
virtual |
Reimplemented from PartBunchBase< double, 3 >.
Definition at line 106 of file PartBunch.cpp.
References BinaryRepartition(), PartBunchBase< double, 3 >::boundp(), PartBunchBase< double, 3 >::get_bounds(), PartBunchBase< double, 3 >::pbase, PartBunchBase< double, 3 >::rmax_m, PartBunchBase< double, 3 >::rmin_m, and PartBunchBase< double, 3 >::update().
|
inlinevirtual |
Here we emit particles from the cathode. All particles in a new simulation (not a restart) initially reside in the bin container "pbin_m" and are not part of the beam bunch (so they cannot "see" fields, space charge etc.). In pbin_m, particles are sorted into the bins of a time histogram that describes the longitudinal time distribution of the beam, where the number of bins is given by \(NBIN \times SBIN\). \(NBIN\) and \(SBIN\) are parameters given when defining the initial beam distribution. During emission, the time step of the simulation is set so that an integral number of these bins are emitted each step. Once all of the particles have been emitted, the simulation time step is reset to the value defined in the input file.
A typical integration time step, \(\Delta t\), is broken down into 3 sub-steps:
1) Drift particles for \(\frac{\Delta t}{2}\).
2) Calculate fields and advance momentum.
3) Drift particles for \(\frac{\Delta t}{2}\) at the new momentum to complete the full time step.
The difficulty for emission is that at the cathode position there is a step function discontinuity in the fields. If we apply the typical integration time step across this boundary, we get an artificial numerical bunching of the beam, especially at very high accelerating fields. This function takes the cathode position boundary into account in order to achieve smoother particle emission.
During an emission step, an integral number of time bins from the distribution histogram are emitted. However, each particle contained in those time bins will actually be emitted from the cathode at a different time, so will only spend some fraction of the total time step, \(\Delta t_{full-timestep}\), in the simulation. The trick to emission is to give each particle a unique time step, \(Delta t_{temp}\), that is equal to the actual time during the emission step that the particle exists in the simulation. For the next integration time step, the particle's time step is set back to the global time step, \(\Delta t_{full-timestep}\).
Implements PartBunchBase< double, 3 >.
Definition at line 1139 of file PartBunch.cpp.
References eg_m, max(), and min().
|
virtual |
Implements PartBunchBase< double, 3 >.
Definition at line 1026 of file PartBunch.cpp.
References RegionLayout< T, Dim, MeshType >::getFieldLayout(), ParticleSpatialLayout< T, Dim, Mesh, CachingPolicy >::getLayout(), and getLayout().
Referenced by computeSelfFields(), computeSelfFields_cycl(), resizeMesh(), runTests(), updateDomainLength(), and updateFields().
|
inlineprivate |
Definition at line 128 of file PartBunch.h.
References PartBunchBase< double, 3 >::pbase.
Referenced by getFieldLayout(), getMesh(), and initialize().
|
inlineprivate |
Definition at line 133 of file PartBunch.h.
References PartBunchBase< double, 3 >::pbase.
|
inline |
Definition at line 146 of file PartBunch.h.
References ParticleSpatialLayout< T, Dim, Mesh, CachingPolicy >::getLayout(), getLayout(), and RegionLayout< T, Dim, MeshType >::getMesh().
Referenced by resizeMesh(), runTests(), and updateFields().
|
inline |
Definition at line 152 of file PartBunch.h.
References ParticleSpatialLayout< T, Dim, Mesh, CachingPolicy >::getLayout(), getLayout(), and RegionLayout< T, Dim, MeshType >::getMesh().
|
inlinevirtual |
Implements PartBunchBase< double, 3 >.
Definition at line 141 of file PartBunch.h.
References rho_m.
|
virtual |
Implements PartBunchBase< double, 3 >.
Definition at line 71 of file PartBunch.cpp.
References RegionLayout< T, Dim, MeshType >::changeDomain(), ParticleSpatialLayout< T, Dim, Mesh, CachingPolicy >::getLayout(), and getLayout().
Definition at line 1166 of file PartBunch.cpp.
References PartBunchBase< T, Dim >::print().
Referenced by operator<<().
|
inlinevirtual |
Reimplemented from PartBunchBase< double, 3 >.
Definition at line 1149 of file PartBunch.cpp.
References ParticleAttrib< T >::destroy(), interpolationCache_m, interpolationCacheSet_m, and ParticleAttrib< T >::size().
|
private |
resize mesh to geometry specified
Definition at line 413 of file PartBunch.cpp.
References bc_m, PartBunchBase< double, 3 >::boundp(), PartBunchBase< double, 3 >::destroy(), eg_m, endl(), PartBunchBase< double, 3 >::fs_m, PartBunchBase< double, 3 >::get_bounds(), getFieldLayout(), PartBunchBase< double, 3 >::getLocalNum(), getMesh(), PoissonSolver::getXRangeMax(), PoissonSolver::getXRangeMin(), PoissonSolver::getYRangeMax(), PoissonSolver::getYRangeMin(), PoissonSolver::getZRangeMax(), PoissonSolver::getZRangeMin(), PartBunchBase< double, 3 >::hr_m, PartBunchBase< double, 3 >::ID, INFOMSG, Field< T, Dim, Mesh, Centering >::initialize(), level2(), Hypervolume::n, PartBunchBase< double, 3 >::nr_m, PartBunchBase< double, 3 >::R, rho_m, PartBunchBase< double, 3 >::rmax_m, PartBunchBase< double, 3 >::rmin_m, UniformCartesian< D, T >::set_meshSpacing(), UniformCartesian< D, T >::set_origin(), FieldSolver::solver_m, PartBunchBase< double, 3 >::update(), and vbc_m.
Referenced by computeSelfFields(), and computeSelfFields_cycl().
|
virtual |
Reimplemented from PartBunchBase< double, 3 >.
Definition at line 76 of file PartBunch.cpp.
References bc_m, PartBunchBase< double, 3 >::Dimension, eg_m, PartBunchBase< double, 3 >::fs_m, FieldLayout< Dim >::getDomain(), getFieldLayout(), getMesh(), PartBunchBase< double, 3 >::hr_m, Field< T, Dim, Mesh, Centering >::initialize(), PartBunchBase< double, 3 >::nr_m, rho_m, UniformCartesian< D, T >::set_meshSpacing(), UniformCartesian< D, T >::set_origin(), setBCAllPeriodic(), FieldSolver::solver_m, PoissonSolver::test(), and vbc_m.
|
virtual |
Reimplemented from PartBunchBase< double, 3 >.
Definition at line 1050 of file PartBunch.cpp.
References bc_m, PartBunchBase< double, 3 >::dcBeam_m, endl(), PartBunchBase< double, 3 >::getBConds(), INFOMSG, level3(), ParticleNoBCond(), and vbc_m.
|
virtual |
Reimplemented from PartBunchBase< double, 3 >.
Definition at line 1031 of file PartBunch.cpp.
References bc_m, PartBunchBase< double, 3 >::dcBeam_m, endl(), PartBunchBase< double, 3 >::getBConds(), IpplInfo::getNodes(), INFOMSG, level3(), ParticlePeriodicBCond(), and vbc_m.
Referenced by runTests().
|
virtual |
Reimplemented from PartBunchBase< double, 3 >.
Definition at line 1060 of file PartBunch.cpp.
References bc_m, PartBunchBase< double, 3 >::dcBeam_m, endl(), PartBunchBase< double, 3 >::getBConds(), IpplInfo::getNodes(), INFOMSG, level3(), ParticleNoBCond(), ParticlePeriodicBCond(), and vbc_m.
|
virtual |
Reimplemented from PartBunchBase< double, 3 >.
Definition at line 1156 of file PartBunch.cpp.
References interpolationCache_m, interpolationCacheSet_m, and PartBunchBase< T, Dim >::swap().
|
privatevirtual |
Implements PartBunchBase< double, 3 >.
Definition at line 1083 of file PartBunch.cpp.
References PartBunchBase< double, 3 >::Dimension, FieldLayout< Dim >::getDomain(), and getFieldLayout().
Reimplemented from PartBunchBase< double, 3 >.
Definition at line 1090 of file PartBunch.cpp.
References bc_m, eg_m, getFieldLayout(), getMesh(), PartBunchBase< double, 3 >::hr_m, Field< T, Dim, Mesh, Centering >::initialize(), rho_m, UniformCartesian< D, T >::set_meshSpacing(), UniformCartesian< D, T >::set_origin(), and vbc_m.
for defining the boundary conditions
Definition at line 119 of file PartBunch.h.
Referenced by resizeMesh(), runTests(), setBCAllOpen(), setBCAllPeriodic(), setBCForDCBeam(), and updateFields().
VField_t PartBunch::eg_m |
vector field on the grid
Definition at line 105 of file PartBunch.h.
Referenced by computeSelfFields(), computeSelfFields_cycl(), getEExtrema(), resizeMesh(), runTests(), and updateFields().
|
private |
Definition at line 125 of file PartBunch.h.
Referenced by computeSelfFields(), resetInterpolationCache(), and swap().
|
private |
Definition at line 123 of file PartBunch.h.
Referenced by computeSelfFields(), resetInterpolationCache(), and swap().
Field_t PartBunch::rho_m |
scalar potential
Definition at line 102 of file PartBunch.h.
Referenced by computeSelfFields(), computeSelfFields_cycl(), getRho(), resizeMesh(), runTests(), and updateFields().
Definition at line 120 of file PartBunch.h.
Referenced by resizeMesh(), runTests(), setBCAllOpen(), setBCAllPeriodic(), setBCForDCBeam(), and updateFields().