BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | Friends
BDSBunch Class Reference

The base class for bunch distribution generators. More...

#include <BDSBunch.hh>

Inheritance diagram for BDSBunch:
Inheritance graph
Collaboration diagram for BDSBunch:
Collaboration graph

Public Member Functions

 BDSBunch (const G4String &nameIn)
 
virtual void SetOptions (const BDSParticleDefinition *beamParticle, const GMAD::Beam &beam, const BDSBunchType &distrType, G4Transform3D beamlineTransformIn=G4Transform3D::Identity, const G4double beamlineS=0)
 
virtual void CheckParameters ()
 
virtual void Initialise ()
 Any initialisation - to be used after SetOptions, then CheckParameters.
 
BDSParticleCoordsFullGlobal GetNextParticle ()
 
virtual G4bool ExpectChangingParticleType () const
 A hint of whether we expect to require and extended particle set (ie pions, kaons, muons).
 
virtual BDSParticleCoordsFullGlobal GetNextParticleValid (G4int maxTries=100)
 
virtual void BeginOfRunAction (G4int numberOfEvents, G4bool batchMode)
 
virtual const BDSParticleDefinitionParticleDefinition () const
 Access the beam particle definition.
 
virtual void SetGeneratePrimariesOnly (G4bool generatePrimariesOnlyIn)
 
virtual BDSParticleCoordsFull GetNextParticleLocal ()
 
G4bool UseCurvilinearTransform () const
 Access whether there's a finite S offset and therefore we're using a CL transform.
 
virtual void RecreateAdvanceToEvent (G4int eventOffset)
 
G4bool BeamParticleIsAnIon () const
 Access whether the beam particle is an ion or not.
 
virtual void UpdateIonDefinition ()
 
G4bool ParticleDefinitionHasBeenUpdated () const
 
G4String Name () const
 Distribution name.
 
G4int CurrentBunchIndex () const
 Get the current bunch index for writing to output.
 
void CalculateBunchIndex (G4int eventIndex)
 Calculate which bunch index we should be at given an event index.
 

Static Public Member Functions

static G4double CalculateZp (G4double xp, G4double yp, G4double Zp0)
 Calculate zp safely based on other components.
 
static void SetEmittances (const BDSParticleDefinition *beamParticle, const GMAD::Beam &beam, G4double &emittGeometricX, G4double &emittGeometricY, G4double &emittNormalisedX, G4double &emittNormalisedY)
 

Protected Member Functions

BDSParticleCoordsFullGlobal ApplyTransform (const BDSParticleCoordsFull &localIn) const
 
void ApplyTilt (BDSParticleCoordsFull &localIn) const
 Apply a rotation about unitZ for the local coordinates according to member variable tilt.
 
void ApplyBunchTiming (BDSParticleCoordsFullGlobal &localIn) const
 Add on the offset in T for the current bunch number (i*bunchPeriod).
 
BDSParticleCoordsFullGlobal ApplyCurvilinearTransform (const BDSParticleCoordsFull &localIn) const
 Calculate the global coordinates from curvilinear coordinates of a beam line.
 

Protected Attributes

G4String name
 Name of distribution.
 
G4bool useCurvilinear
 Whether to ignore z and use s and transform for curvilinear coordinates.
 
BDSParticleDefinitionparticleDefinition
 Particle definition for bunch - this class owns it.
 
G4bool particleDefinitionHasBeenUpdated
 
G4bool finiteTilt
 
G4bool generatePrimariesOnly
 
G4double X0
 Centre of distributions.
 
G4double Y0
 Centre of distributions.
 
G4double Z0
 Centre of distributions.
 
G4double S0
 Centre of distributions.
 
G4double T0
 Centre of distributions.
 
G4double Xp0
 Centre of distributions.
 
G4double Yp0
 Centre of distributions.
 
G4double Zp0
 Centre of distributions.
 
G4double E0
 Centre of distributions.
 
G4double P0
 central momentum
 
G4double tilt
 Centre of distributions.
 
G4double sigmaT
 Centre of distributions.
 
G4double sigmaP
 Centre of distributions.
 
G4double sigmaE
 Centre of distributions.
 
G4double sigmaEk
 Centre of distributions.
 
bool useBunchTiming
 Bunch offset in time parameters.
 
G4int currentBunchIndex
 Bunch offset in time parameters.
 
G4int eventsPerBunch
 Bunch offset in time parameters.
 
G4double bunchPeriod
 Bunch offset in time parameters.
 
G4bool finiteSigmaE
 Flags to ignore random number generator in case of no finite E or T.
 
G4bool finiteSigmaT
 Flags to ignore random number generator in case of no finite E or T.
 

Private Attributes

G4Transform3D beamlineTransform
 Transform that beam line starts with that will also be applied to coordinates.
 
G4double beamlineS
 Beamline initial S position.
 
G4double mass2
 Cache of mass squared as required to convert from p to E.
 
const BDSBeamlinebeamline
 A reference to the fully constructed beamline that's lazily instantiated.
 

Friends

class BDSPrimaryGeneratorFileHEPMC
 Make BDSPrimaryGeneratorFileHEPMC a friend so it can use the protected ApplyTransform function.
 
class BDSPrimaryGeneratorFileSampler
 

Detailed Description

The base class for bunch distribution generators.

After instantiation of a class, use SetOptions(), then CheckParmeters(), then Initialise().

Author
Stewart Boogert

Definition at line 46 of file BDSBunch.hh.

Constructor & Destructor Documentation

◆ BDSBunch() [1/2]

BDSBunch::BDSBunch ( )

Definition at line 50 of file BDSBunch.cc.

◆ BDSBunch() [2/2]

BDSBunch::BDSBunch ( const G4String &  nameIn)
explicit

Definition at line 54 of file BDSBunch.cc.

◆ ~BDSBunch()

BDSBunch::~BDSBunch ( )
virtual

Definition at line 77 of file BDSBunch.cc.

Member Function Documentation

◆ ApplyBunchTiming()

void BDSBunch::ApplyBunchTiming ( BDSParticleCoordsFullGlobal localIn) const
protected

Add on the offset in T for the current bunch number (i*bunchPeriod).

Definition at line 316 of file BDSBunch.cc.

References bunchPeriod, and currentBunchIndex.

Referenced by ApplyTransform().

Here is the caller graph for this function:

◆ ApplyCurvilinearTransform()

BDSParticleCoordsFullGlobal BDSBunch::ApplyCurvilinearTransform ( const BDSParticleCoordsFull localIn) const
protected

Calculate the global coordinates from curvilinear coordinates of a beam line.

Definition at line 322 of file BDSBunch.cc.

References beamline, BDSParticleCoordsFullGlobal::beamlineIndex, BDSAcceleratorModel::BeamlineMain(), generatePrimariesOnly, BDSBeamline::GetGlobalEuclideanTransform(), and S0.

Referenced by ApplyTransform().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ApplyTilt()

void BDSBunch::ApplyTilt ( BDSParticleCoordsFull localIn) const
protected

Apply a rotation about unitZ for the local coordinates according to member variable tilt.

Definition at line 304 of file BDSBunch.cc.

References tilt.

Referenced by GetNextParticle().

Here is the caller graph for this function:

◆ ApplyTransform()

BDSParticleCoordsFullGlobal BDSBunch::ApplyTransform ( const BDSParticleCoordsFull localIn) const
protected

Apply either the curvilinear transform if we're using curvilinear coordinates or just apply the general beam line offset in global coordinates to the 'local' curvilinear coordinates.

Definition at line 292 of file BDSBunch.cc.

References ApplyBunchTiming(), ApplyCurvilinearTransform(), BDSParticleCoords::ApplyTransform(), beamlineTransform, useBunchTiming, and useCurvilinear.

Referenced by GetNextParticle(), BDSPrimaryGeneratorFileHEPMC::HepMC2G4(), and BDSPrimaryGeneratorFileSampler::ReadSingleEvent().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ BeamParticleIsAnIon()

G4bool BDSBunch::BeamParticleIsAnIon ( ) const
inline

Access whether the beam particle is an ion or not.

Definition at line 119 of file BDSBunch.hh.

References BDSParticleDefinition::IsAnIon(), and particleDefinition.

Referenced by BDSPrimaryGeneratorAction::BDSPrimaryGeneratorAction().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ BeginOfRunAction()

void BDSBunch::BeginOfRunAction ( G4int  numberOfEvents,
G4bool  batchMode 
)
virtual

An action that is called at the beginning of a run when we know the number of events that'll be generated. By default this is nothing, but can be used to calculate sample mean offsets in some derived classes.

Reimplemented in BDSBunchFileBased, and BDSBunchGaussian.

Definition at line 285 of file BDSBunch.cc.

Referenced by BDSIM::GeneratePrimariesOnly().

Here is the caller graph for this function:

◆ CalculateBunchIndex()

void BDSBunch::CalculateBunchIndex ( G4int  eventIndex)

Calculate which bunch index we should be at given an event index.

Definition at line 215 of file BDSBunch.cc.

References currentBunchIndex, eventsPerBunch, and useBunchTiming.

Referenced by BDSPrimaryGeneratorAction::GeneratePrimaries(), and RecreateAdvanceToEvent().

Here is the caller graph for this function:

◆ CalculateZp()

G4double BDSBunch::CalculateZp ( G4double  xp,
G4double  yp,
G4double  Zp0 
)
static

◆ CheckParameters()

void BDSBunch::CheckParameters ( )
virtual

◆ CurrentBunchIndex()

G4int BDSBunch::CurrentBunchIndex ( ) const
inline

Get the current bunch index for writing to output.

Definition at line 148 of file BDSBunch.hh.

References currentBunchIndex.

Referenced by BDSPrimaryGeneratorAction::GeneratePrimaries().

Here is the caller graph for this function:

◆ ExpectChangingParticleType()

virtual G4bool BDSBunch::ExpectChangingParticleType ( ) const
inlinevirtual

A hint of whether we expect to require and extended particle set (ie pions, kaons, muons).

Reimplemented in BDSBunchUserFile< T >.

Definition at line 80 of file BDSBunch.hh.

Referenced by BDSIM::Initialise().

Here is the caller graph for this function:

◆ GetNextParticle()

BDSParticleCoordsFullGlobal BDSBunch::GetNextParticle ( )

Main interface. Calls GetNextParticleLocal() and then applies the curvilinear transform if required.

Definition at line 262 of file BDSBunch.cc.

References ApplyTilt(), ApplyTransform(), finiteTilt, GetNextParticleLocal(), and particleDefinitionHasBeenUpdated.

Referenced by BDSBunchUserFile< T >::GetNextParticleValid(), and GetNextParticleValid().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetNextParticleLocal()

BDSParticleCoordsFull BDSBunch::GetNextParticleLocal ( )
virtual

Each derived class can override this default method of reference position. If S0 > 0 or derived class changes member bool 'curvilinear' z0 will be treated as S and the global z0 be calculated.

Reimplemented in BDSBunchBox, BDSBunchCircle, BDSBunchComposite, BDSBunchCompositeSDE, BDSBunchEShell, BDSBunchGaussian, BDSBunchHalo, BDSBunchHaloFlatSigma, BDSBunchPtc, BDSBunchRing, BDSBunchSixTrack, BDSBunchSixTrackLink, BDSBunchSphere, BDSBunchSquare, and BDSBunchUserFile< T >.

Definition at line 272 of file BDSBunch.cc.

References E0, S0, T0, X0, Xp0, Y0, Yp0, Z0, and Zp0.

Referenced by BDSLinkPrimaryGeneratorAction::GeneratePrimaries(), GetNextParticle(), BDSBunchComposite::GetNextParticleLocal(), BDSBunchCompositeSDE::GetNextParticleLocal(), BDSPrimaryGeneratorFileHEPMC::HepMC2G4(), and BDSPrimaryGeneratorFileSampler::ReadSingleEvent().

Here is the caller graph for this function:

◆ GetNextParticleValid()

BDSParticleCoordsFullGlobal BDSBunch::GetNextParticleValid ( G4int  maxTries = 100)
virtual

Interface to allow multiple calls until a safe particle is generated. This will repeatedly call GetNextParticle() until a particle is generated where the total energy is greater than the rest mass. This can be overridden by derived classes if they do not wish this behaviour. This may throw an exception if the maximum number of attempts (maxTries) is exceeded.

Reimplemented in BDSBunchUserFile< T >.

Definition at line 234 of file BDSBunch.cc.

References GetNextParticle(), BDSParticleDefinition::Mass(), particleDefinition, and particleDefinitionHasBeenUpdated.

Referenced by BDSPrimaryGeneratorAction::GeneratePrimaries(), and BDSIM::GeneratePrimariesOnly().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ Initialise()

void BDSBunch::Initialise ( )
virtual

Any initialisation - to be used after SetOptions, then CheckParameters.

Reimplemented in BDSBunchPtc, and BDSBunchUserFile< T >.

Definition at line 231 of file BDSBunch.cc.

◆ Name()

G4String BDSBunch::Name ( ) const
inline

Distribution name.

Definition at line 145 of file BDSBunch.hh.

References name.

Referenced by BDSIM::Initialise(), BDSBunchComposite::SetOptions(), and BDSBunchCompositeSDE::SetOptions().

Here is the caller graph for this function:

◆ ParticleDefinition()

virtual const BDSParticleDefinition * BDSBunch::ParticleDefinition ( ) const
inlinevirtual

Access the beam particle definition.

Reimplemented in BDSBunchSixTrackLink.

Definition at line 96 of file BDSBunch.hh.

References particleDefinition.

Referenced by BDSLinkPrimaryGeneratorAction::GeneratePrimaries(), BDSPrimaryGeneratorAction::GeneratePrimaries(), BDSIM::GeneratePrimariesOnly(), BDSIM::Initialise(), and BDSIMLink::Initialise().

Here is the caller graph for this function:

◆ ParticleDefinitionHasBeenUpdated()

G4bool BDSBunch::ParticleDefinitionHasBeenUpdated ( ) const
inline

Whether the particle definition has been updated since the last call to GetNextParticle() or GetNextParticleValid().

Definition at line 129 of file BDSBunch.hh.

References particleDefinitionHasBeenUpdated.

Referenced by BDSBunchComposite::GetNextParticleLocal(), and BDSBunchCompositeSDE::GetNextParticleLocal().

Here is the caller graph for this function:

◆ RecreateAdvanceToEvent()

void BDSBunch::RecreateAdvanceToEvent ( G4int  eventOffset)
virtual

When recreating events, it's possible that setting the seed state may not be sufficient for the bunch to get the right distribution. This is true when the bunch coordinates are based on an external source of data i.e. user bunch file. This default method allows such a distribution to advance to the correct event number. Updates the bunch index by default.

Reimplemented in BDSBunchPtc, and BDSBunchUserFile< T >.

Definition at line 280 of file BDSBunch.cc.

References CalculateBunchIndex().

Referenced by BDSPrimaryGeneratorAction::BDSPrimaryGeneratorAction(), and BDSBunchUserFile< T >::RecreateAdvanceToEvent().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetEmittances()

void BDSBunch::SetEmittances ( const BDSParticleDefinition beamParticle,
const GMAD::Beam beam,
G4double &  emittGeometricX,
G4double &  emittGeometricY,
G4double &  emittNormalisedX,
G4double &  emittNormalisedY 
)
static

Work out whether either the geometric or normalised emittances are set and update the variables by reference with the values. Can throw exception if more than one set in either x and y.

Definition at line 175 of file BDSBunch.cc.

References BDS::ConflictingParametersSet(), GMAD::BeamBase::emitNX, GMAD::BeamBase::emitNY, GMAD::BeamBase::emitx, GMAD::BeamBase::emity, BDSParticleDefinition::Gamma(), BDS::IsFinite(), and BDS::NBeamParametersSet().

Referenced by BDSBunchHalo::SetOptions(), BDSBunchHaloFlatSigma::SetOptions(), and BDSBunchTwiss::SetOptions().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetGeneratePrimariesOnly()

void BDSBunch::SetGeneratePrimariesOnly ( G4bool  generatePrimariesOnlyIn)
virtual

Set the flag to whether we're only generating primaries only. This sets the member variable generatePrimariesOnly which skips trying to perform a curvilinear transform if true. i.e. if we're generating primaries only, there's no beam line.

Reimplemented in BDSBunchComposite, and BDSBunchCompositeSDE.

Definition at line 289 of file BDSBunch.cc.

References generatePrimariesOnly.

Referenced by BDSBunchComposite::SetGeneratePrimariesOnly(), and BDSBunchCompositeSDE::SetGeneratePrimariesOnly().

Here is the caller graph for this function:

◆ SetOptions()

void BDSBunch::SetOptions ( const BDSParticleDefinition beamParticle,
const GMAD::Beam beam,
const BDSBunchType distrType,
G4Transform3D  beamlineTransformIn = G4Transform3D::Identity,
const G4double  beamlineS = 0 
)
virtual

Extract and set the relevant options from the beam definition. The distribution type is explicitly required as this function may be used inside a nested bunch distribution. This argument is for the most part ignored, but there's no way to have a default for it. Also, some classes can cover multiple input distributions so need to know which one they're meant to be.

Reimplemented in BDSBunchBox, BDSBunchComposite, BDSBunchCompositeSDE, BDSBunchEShell, BDSBunchFileBased, BDSBunchGaussian, BDSBunchHalo, BDSBunchHaloFlatSigma, BDSBunchPtc, BDSBunchRing, BDSBunchSigmaMatrix, BDSBunchSixTrack, BDSBunchSquare, BDSBunchTwiss, BDSBunchUserFile< T >, BDSBunchCircle, and BDSBunchEventGenerator.

Definition at line 82 of file BDSBunch.cc.

References beamlineS, beamlineTransform, BDSParticleDefinition::Beta(), GMAD::BeamBase::bunchFrequency, bunchPeriod, GMAD::BeamBase::bunchPeriod, CalculateZp(), BDS::ConflictingParametersSet(), E0, eventsPerBunch, GMAD::BeamBase::eventsPerBunch, finiteSigmaE, finiteSigmaT, finiteTilt, BDSGlobalConstants::Instance(), BDS::IsFinite(), BDSParticleDefinition::KineticEnergy(), BDSParticleDefinition::Momentum(), BDS::NBeamParametersSet(), P0, particleDefinition, S0, GMAD::BeamBase::S0, sigmaE, GMAD::BeamBase::sigmaE, sigmaEk, sigmaP, sigmaT, GMAD::BeamBase::sigmaT, T0, GMAD::BeamBase::T0, tilt, GMAD::BeamBase::tilt, BDSParticleDefinition::TotalEnergy(), useBunchTiming, useCurvilinear, X0, GMAD::BeamBase::X0, Xp0, GMAD::BeamBase::Xp0, Y0, GMAD::BeamBase::Y0, Yp0, GMAD::BeamBase::Yp0, Z0, GMAD::BeamBase::Z0, Zp0, and GMAD::BeamBase::Zp0.

Referenced by BDSIMLink::Initialise(), BDSBunchBox::SetOptions(), BDSBunchComposite::SetOptions(), BDSBunchCompositeSDE::SetOptions(), BDSBunchEShell::SetOptions(), BDSBunchFileBased::SetOptions(), BDSBunchGaussian::SetOptions(), BDSBunchHalo::SetOptions(), BDSBunchHaloFlatSigma::SetOptions(), BDSBunchRing::SetOptions(), BDSBunchSixTrack::SetOptions(), BDSBunchSquare::SetOptions(), and BDSBunchCircle::SetOptions().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ UpdateIonDefinition()

void BDSBunch::UpdateIonDefinition ( )
virtual

Update the Geant4 Particle Definition inside particleDefinition member in the case it's an ion. This can only be done later one once the run has started. Since the particle definition is kept here in the bunch this interface allows control over it being updated.

Reimplemented in BDSBunchSixTrackLink.

Definition at line 382 of file BDSBunch.cc.

References BDSIonDefinition::A(), BDSIonDefinition::ExcitationEnergy(), BDS::FixGeant105ThreshholdsForParticle(), BDSParticleDefinition::IonDefinition(), BDSParticleDefinition::IsAnIon(), particleDefinition, BDSParticleDefinition::UpdateG4ParticleDefinition(), and BDSIonDefinition::Z().

Referenced by BDSPrimaryGeneratorAction::GeneratePrimaries().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ UseCurvilinearTransform()

G4bool BDSBunch::UseCurvilinearTransform ( ) const
inline

Access whether there's a finite S offset and therefore we're using a CL transform.

Definition at line 109 of file BDSBunch.hh.

References useCurvilinear.

Referenced by BDSPrimaryGeneratorAction::GeneratePrimaries().

Here is the caller graph for this function:

Friends And Related Function Documentation

◆ BDSPrimaryGeneratorFileHEPMC

friend class BDSPrimaryGeneratorFileHEPMC
friend

Make BDSPrimaryGeneratorFileHEPMC a friend so it can use the protected ApplyTransform function.

Definition at line 54 of file BDSBunch.hh.

◆ BDSPrimaryGeneratorFileSampler

friend class BDSPrimaryGeneratorFileSampler
friend

Definition at line 55 of file BDSBunch.hh.

Field Documentation

◆ beamline

const BDSBeamline* BDSBunch::beamline
mutableprivate

A reference to the fully constructed beamline that's lazily instantiated.

Definition at line 223 of file BDSBunch.hh.

Referenced by ApplyCurvilinearTransform().

◆ beamlineS

G4double BDSBunch::beamlineS
private

Beamline initial S position.

Definition at line 219 of file BDSBunch.hh.

Referenced by SetOptions().

◆ beamlineTransform

G4Transform3D BDSBunch::beamlineTransform
private

Transform that beam line starts with that will also be applied to coordinates.

Definition at line 217 of file BDSBunch.hh.

Referenced by ApplyTransform(), and SetOptions().

◆ bunchPeriod

G4double BDSBunch::bunchPeriod
protected

Bunch offset in time parameters.

Definition at line 192 of file BDSBunch.hh.

Referenced by ApplyBunchTiming(), and SetOptions().

◆ currentBunchIndex

G4int BDSBunch::currentBunchIndex
protected

Bunch offset in time parameters.

Definition at line 190 of file BDSBunch.hh.

Referenced by ApplyBunchTiming(), CalculateBunchIndex(), and CurrentBunchIndex().

◆ E0

G4double BDSBunch::E0
protected

◆ eventsPerBunch

G4int BDSBunch::eventsPerBunch
protected

Bunch offset in time parameters.

Definition at line 191 of file BDSBunch.hh.

Referenced by CalculateBunchIndex(), and SetOptions().

◆ finiteSigmaE

G4bool BDSBunch::finiteSigmaE
protected

Flags to ignore random number generator in case of no finite E or T.

Definition at line 207 of file BDSBunch.hh.

Referenced by BDSBunchGaussian::GetNextParticleLocalCoords(), SetOptions(), BDSBunchGaussian::SetOptions(), and BDSBunchSigmaMatrix::SetOptions().

◆ finiteSigmaT

G4bool BDSBunch::finiteSigmaT
protected

Flags to ignore random number generator in case of no finite E or T.

Definition at line 208 of file BDSBunch.hh.

Referenced by BDSBunchGaussian::GetNextParticleLocalCoords(), SetOptions(), BDSBunchGaussian::SetOptions(), and BDSBunchSigmaMatrix::SetOptions().

◆ finiteTilt

G4bool BDSBunch::finiteTilt
protected

Flag of whether to apply beam rotation.

Definition at line 205 of file BDSBunch.hh.

Referenced by GetNextParticle(), and SetOptions().

◆ generatePrimariesOnly

G4bool BDSBunch::generatePrimariesOnly
protected

Whether we're only generating the primaries only and therefore there's no beamline to get a transform from and we shouldn't try.

Definition at line 213 of file BDSBunch.hh.

Referenced by ApplyCurvilinearTransform(), and SetGeneratePrimariesOnly().

◆ mass2

G4double BDSBunch::mass2
private

Cache of mass squared as required to convert from p to E.

Definition at line 220 of file BDSBunch.hh.

◆ name

G4String BDSBunch::name
protected

Name of distribution.

Definition at line 168 of file BDSBunch.hh.

Referenced by Name(), BDSBunchComposite::SetOptions(), and BDSBunchCompositeSDE::SetOptions().

◆ P0

G4double BDSBunch::P0
protected

central momentum

Definition at line 180 of file BDSBunch.hh.

Referenced by SetOptions().

◆ particleDefinition

BDSParticleDefinition* BDSBunch::particleDefinition
protected

◆ particleDefinitionHasBeenUpdated

G4bool BDSBunch::particleDefinitionHasBeenUpdated
protected

◆ S0

G4double BDSBunch::S0
protected

◆ sigmaE

G4double BDSBunch::sigmaE
protected

◆ sigmaEk

G4double BDSBunch::sigmaEk
protected

Centre of distributions.

Definition at line 185 of file BDSBunch.hh.

Referenced by SetOptions().

◆ sigmaP

G4double BDSBunch::sigmaP
protected

Centre of distributions.

Definition at line 183 of file BDSBunch.hh.

Referenced by SetOptions(), and BDSBunchTwiss::SetOptions().

◆ sigmaT

G4double BDSBunch::sigmaT
protected

◆ T0

G4double BDSBunch::T0
protected

◆ tilt

G4double BDSBunch::tilt
protected

Centre of distributions.

Definition at line 181 of file BDSBunch.hh.

Referenced by ApplyTilt(), and SetOptions().

◆ useBunchTiming

bool BDSBunch::useBunchTiming
protected

Bunch offset in time parameters.

Definition at line 189 of file BDSBunch.hh.

Referenced by ApplyTransform(), CalculateBunchIndex(), and SetOptions().

◆ useCurvilinear

G4bool BDSBunch::useCurvilinear
protected

Whether to ignore z and use s and transform for curvilinear coordinates.

Definition at line 196 of file BDSBunch.hh.

Referenced by ApplyTransform(), BDSBunchUserFile< T >::ParseFileFormat(), SetOptions(), and UseCurvilinearTransform().

◆ X0

G4double BDSBunch::X0
protected

◆ Xp0

G4double BDSBunch::Xp0
protected

◆ Y0

G4double BDSBunch::Y0
protected

◆ Yp0

G4double BDSBunch::Yp0
protected

◆ Z0

G4double BDSBunch::Z0
protected

◆ Zp0

G4double BDSBunch::Zp0
protected

The documentation for this class was generated from the following files: