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
BDSPrimaryGeneratorFile Class Referenceabstract

Common interface for any primary generators that are file based. More...

#include <BDSPrimaryGeneratorFile.hh>

Inheritance diagram for BDSPrimaryGeneratorFile:
Inheritance graph
Collaboration diagram for BDSPrimaryGeneratorFile:
Collaboration graph

Public Member Functions

 BDSPrimaryGeneratorFile (G4bool loopFileIn, BDSBunchEventGenerator *bunchIn)
 
G4bool GeneratePrimaryVertexSafe (G4Event *event)
 Return false if not able to generate a primary vertex.
 
virtual void RecreateAdvanceToEvent (G4int eventOffset)=0
 Advance into the file as required.
 
G4long NEventsInFile () const
 
G4long NEventsLeftInFile () const
 
G4bool OKToLoopFile () const
 
G4long NEventsReadThatPassedFilters () const
 Accessor.
 
G4long NEventsSkipped () const
 Return the offset into the file if any.
 
G4bool DistributionIsFinished () const
 Report whether the distribution is finished generating.
 
void ThrowExceptionIfRecreateOffsetTooHigh (G4long eventOffset) const
 

Static Public Member Functions

static BDSPrimaryGeneratorFileConstructGenerator (const GMAD::Beam &beam, BDSBunch *bunchIn, G4bool recreate, G4int eventOffset, G4bool batchMode)
 

Protected Member Functions

G4bool VertexInsideWorld (const G4ThreeVector &pos) const
 Utility function for derived classes to check a position is inside the world.
 

Protected Attributes

G4bool loopFile
 
BDSBunchEventGeneratorbunch
 
G4bool endOfFileReached
 
G4bool vertexGeneratedSuccessfully
 
G4long currentFileEventIndex
 
G4long nEventsInFile
 
G4long nEventsReadThatPassedFilters
 
G4long nEventsSkipped
 
G4VSolid * worldSolid
 

Detailed Description

Common interface for any primary generators that are file based.

Any derived class must implement G4VPrimaryGeneratorAction::GeneratePrimaryVertex(event). In this, it must update the protected members in this class of endOfFileReached and vertexGeneratedSuccessfully.

Matching the file length is done externally but by using the DistributionIsFinished() method in this class. Whether looping is permitted is known about in this class.

Author
Laurie Nevay

Definition at line 49 of file BDSPrimaryGeneratorFile.hh.

Constructor & Destructor Documentation

◆ BDSPrimaryGeneratorFile()

BDSPrimaryGeneratorFile::BDSPrimaryGeneratorFile ( G4bool  loopFileIn,
BDSBunchEventGenerator bunchIn 
)

Definition at line 47 of file BDSPrimaryGeneratorFile.cc.

◆ ~BDSPrimaryGeneratorFile()

BDSPrimaryGeneratorFile::~BDSPrimaryGeneratorFile ( )
virtual

Definition at line 60 of file BDSPrimaryGeneratorFile.cc.

Member Function Documentation

◆ ConstructGenerator()

BDSPrimaryGeneratorFile * BDSPrimaryGeneratorFile::ConstructGenerator ( const GMAD::Beam beam,
BDSBunch bunchIn,
G4bool  recreate,
G4int  eventOffset,
G4bool  batchMode 
)
static

◆ DistributionIsFinished()

G4bool BDSPrimaryGeneratorFile::DistributionIsFinished ( ) const
inline

Report whether the distribution is finished generating.

Definition at line 90 of file BDSPrimaryGeneratorFile.hh.

Referenced by BDSPrimaryGeneratorAction::GeneratePrimariesFromFile().

Here is the caller graph for this function:

◆ GeneratePrimaryVertexSafe()

G4bool BDSPrimaryGeneratorFile::GeneratePrimaryVertexSafe ( G4Event *  event)

Return false if not able to generate a primary vertex.

Definition at line 73 of file BDSPrimaryGeneratorFile.cc.

Referenced by BDSPrimaryGeneratorAction::GeneratePrimariesFromFile().

Here is the caller graph for this function:

◆ NEventsInFile()

G4long BDSPrimaryGeneratorFile::NEventsInFile ( ) const
inline

Return the number of events in the file - not necessarily the number that match the filters but that are there in total.

Definition at line 73 of file BDSPrimaryGeneratorFile.hh.

◆ NEventsLeftInFile()

G4long BDSPrimaryGeneratorFile::NEventsLeftInFile ( ) const

Return number of available events (excluding any filters in derived classes) in the file. This is nominally nEventsInFile - currentFileEventIndex.

Definition at line 82 of file BDSPrimaryGeneratorFile.cc.

Referenced by ConstructGenerator().

Here is the caller graph for this function:

◆ NEventsReadThatPassedFilters()

G4long BDSPrimaryGeneratorFile::NEventsReadThatPassedFilters ( ) const
inline

Accessor.

Definition at line 84 of file BDSPrimaryGeneratorFile.hh.

Referenced by BDSPrimaryGeneratorAction::GeneratePrimariesFromFile().

Here is the caller graph for this function:

◆ NEventsSkipped()

G4long BDSPrimaryGeneratorFile::NEventsSkipped ( ) const
inline

Return the offset into the file if any.

Definition at line 87 of file BDSPrimaryGeneratorFile.hh.

◆ OKToLoopFile()

G4bool BDSPrimaryGeneratorFile::OKToLoopFile ( ) const

Return whether at least 1 event has passed a filter in the file. If we have no events in the file that pass then we will loop infinitely to find one.

Definition at line 87 of file BDSPrimaryGeneratorFile.cc.

Referenced by BDSPrimaryGeneratorFileSampler::GeneratePrimaryVertex(), and BDSPrimaryGeneratorFileHEPMC::ReadSingleEvent().

Here is the caller graph for this function:

◆ RecreateAdvanceToEvent()

virtual void BDSPrimaryGeneratorFile::RecreateAdvanceToEvent ( G4int  eventOffset)
pure virtual

Advance into the file as required.

Implemented in BDSPrimaryGeneratorFileHEPMC, and BDSPrimaryGeneratorFileSampler.

Referenced by ConstructGenerator().

Here is the caller graph for this function:

◆ ThrowExceptionIfRecreateOffsetTooHigh()

void BDSPrimaryGeneratorFile::ThrowExceptionIfRecreateOffsetTooHigh ( G4long  eventOffset) const

Utility function to check eventOffset < nEventsInFile. This can only possibly happen If the recreation is done with the same filename but different contents or indeed somehow a different file.

Definition at line 63 of file BDSPrimaryGeneratorFile.cc.

References BDSBunchFileBased::DistrFileLoopNTimes().

Referenced by BDSPrimaryGeneratorFileHEPMC::RecreateAdvanceToEvent(), and BDSPrimaryGeneratorFileSampler::RecreateAdvanceToEvent().

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

◆ VertexInsideWorld()

G4bool BDSPrimaryGeneratorFile::VertexInsideWorld ( const G4ThreeVector &  pos) const
protected

Utility function for derived classes to check a position is inside the world.

Definition at line 156 of file BDSPrimaryGeneratorFile.cc.

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

Here is the caller graph for this function:

Field Documentation

◆ bunch

BDSBunchEventGenerator* BDSPrimaryGeneratorFile::bunch
protected

Definition at line 102 of file BDSPrimaryGeneratorFile.hh.

◆ currentFileEventIndex

G4long BDSPrimaryGeneratorFile::currentFileEventIndex
protected

Definition at line 105 of file BDSPrimaryGeneratorFile.hh.

◆ endOfFileReached

G4bool BDSPrimaryGeneratorFile::endOfFileReached
protected

Definition at line 103 of file BDSPrimaryGeneratorFile.hh.

◆ loopFile

G4bool BDSPrimaryGeneratorFile::loopFile
protected

Definition at line 101 of file BDSPrimaryGeneratorFile.hh.

◆ nEventsInFile

G4long BDSPrimaryGeneratorFile::nEventsInFile
protected

Definition at line 106 of file BDSPrimaryGeneratorFile.hh.

◆ nEventsReadThatPassedFilters

G4long BDSPrimaryGeneratorFile::nEventsReadThatPassedFilters
protected

Definition at line 107 of file BDSPrimaryGeneratorFile.hh.

◆ nEventsSkipped

G4long BDSPrimaryGeneratorFile::nEventsSkipped
protected

Definition at line 108 of file BDSPrimaryGeneratorFile.hh.

◆ vertexGeneratedSuccessfully

G4bool BDSPrimaryGeneratorFile::vertexGeneratedSuccessfully
protected

Definition at line 104 of file BDSPrimaryGeneratorFile.hh.

◆ worldSolid

G4VSolid* BDSPrimaryGeneratorFile::worldSolid
mutableprotected

Definition at line 109 of file BDSPrimaryGeneratorFile.hh.


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