BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
|
An intermediate layer for any bunch classes that are file based. More...
#include <BDSBunchFileBased.hh>
Public Member Functions | |
BDSBunchFileBased (const G4String &distributionName) | |
virtual void | BeginOfRunAction (G4int numberOfEvents, G4bool batchMode) |
void | SetNEventsInFile (unsigned long long int nEventsInFileIn) |
void | SetNOriginalEvents (unsigned long long int nOriginalEventsIn) |
void | IncrementNEventsInFileSkipped () |
void | IncrementNEventsInFileSkipped (unsigned long long int plus) |
virtual void | SetOptions (const BDSParticleDefinition *beamParticle, const GMAD::Beam &beam, const BDSBunchType &distrType, G4Transform3D beamlineTransformIn=G4Transform3D::Identity, const G4double beamlineS=0) |
Pull out the relevant options and then pass through to BDSBunch::SetOptions(). | |
BDSBunchFileBased & | operator= (const BDSBunchFileBased &)=delete |
Assignment and copy constructor not implemented nor used. | |
BDSBunchFileBased (BDSBunchFileBased &)=delete | |
Assignment and copy constructor not implemented nor used. | |
unsigned long long int | NOriginalEvents () const |
Accessor. | |
unsigned long long int | NEventsInFile () const |
Accessor. | |
unsigned long long int | NEventsInFileSkipped () const |
Accessor. | |
G4int | DistrFileLoopNTimes () const |
Accessor. | |
![]() | |
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 BDSParticleDefinition * | ParticleDefinition () 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. | |
Protected Attributes | |
unsigned long long int | nOriginalEvents |
nOriginalEvents from upstream file if skimmed - need to pass through. | |
unsigned long long int | nEventsInFile |
The number of entries in the file loaded. | |
unsigned long long int | nEventsInFileSkipped |
Number that are skipped as we go through the file due to filters. | |
G4bool | distrFileLoop |
G4int | distrFileLoopNTimes |
![]() | |
G4String | name |
Name of distribution. | |
G4bool | useCurvilinear |
Whether to ignore z and use s and transform for curvilinear coordinates. | |
BDSParticleDefinition * | particleDefinition |
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. | |
Additional Inherited Members | |
![]() | |
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) |
![]() | |
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. | |
An intermediate layer for any bunch classes that are file based.
This class is purely to hold a record of variables that will be passed to the end of run action. Even for 'distributions' like bdsimsampler and hepmc that don't really use the BDSBunch classes fully, they can publish information in the BDSBunchFileBased instance they get.
Definition at line 37 of file BDSBunchFileBased.hh.
|
explicit |
Definition at line 23 of file BDSBunchFileBased.cc.
|
virtual |
Definition at line 32 of file BDSBunchFileBased.cc.
|
virtual |
Customise the beginning of run action to allow file looping always when in interactive mode, i.e. not batch mode.
Reimplemented from BDSBunch.
Definition at line 46 of file BDSBunchFileBased.cc.
|
inline |
Accessor.
Definition at line 56 of file BDSBunchFileBased.hh.
Referenced by BDSPrimaryGeneratorFileHEPMC::SkipEvents(), BDSPrimaryGeneratorFileSampler::SkipEvents(), and BDSPrimaryGeneratorFile::ThrowExceptionIfRecreateOffsetTooHigh().
|
inline |
Definition at line 62 of file BDSBunchFileBased.hh.
|
inline |
Definition at line 63 of file BDSBunchFileBased.hh.
|
inline |
|
inline |
|
inline |
|
inline |
Definition at line 59 of file BDSBunchFileBased.hh.
|
inline |
Definition at line 60 of file BDSBunchFileBased.hh.
|
virtual |
Pull out the relevant options and then pass through to BDSBunch::SetOptions().
Reimplemented from BDSBunch.
Reimplemented in BDSBunchPtc, BDSBunchUserFile< T >, and BDSBunchEventGenerator.
Definition at line 35 of file BDSBunchFileBased.cc.
References GMAD::BeamBase::distrFileLoop, GMAD::BeamBase::distrFileLoopNTimes, and BDSBunch::SetOptions().
Referenced by BDSBunchPtc::SetOptions(), BDSBunchUserFile< T >::SetOptions(), and BDSBunchEventGenerator::SetOptions().
|
protected |
Definition at line 76 of file BDSBunchFileBased.hh.
|
protected |
Definition at line 77 of file BDSBunchFileBased.hh.
|
protected |
The number of entries in the file loaded.
Definition at line 74 of file BDSBunchFileBased.hh.
Referenced by BDSBunchPtc::Initialise(), BDSBunchUserFile< T >::Initialise(), and NEventsInFile().
|
protected |
Number that are skipped as we go through the file due to filters.
Definition at line 75 of file BDSBunchFileBased.hh.
Referenced by NEventsInFileSkipped().
|
protected |
nOriginalEvents from upstream file if skimmed - need to pass through.
Definition at line 73 of file BDSBunchFileBased.hh.
Referenced by NOriginalEvents().