ATLAS Offline Software
Loading...
Searching...
No Matches
BeamHaloGenerator Class Referenceabstract

An abstract base class to provide generic beam halo generator functionality. More...

#include <BeamHaloGenerator.h>

Inheritance diagram for BeamHaloGenerator:
Collaboration diagram for BeamHaloGenerator:

Public Member Functions

 BeamHaloGenerator (const HepPDT::ParticleDataTable *particleTable, const std::string &inputFile, const std::vector< std::string > &generatorSettings)
virtual ~BeamHaloGenerator ()
virtual int genInitialize ()
 A function to initialise the generator.
virtual int genFinalize ()
 A function to finalise the generator.
virtual int fillEvt (HepMC::GenEvent *evt, CLHEP::HepRandomEngine *engine)=0
 A function to create one event in HepMC format.
void setInterfacePlane (double interfacePlane)
 Set the position of the interface plane in mm.
void setEnableFlip (bool enableFlip)
 Turn on or off the event flipping code.
void setFlipProbability (float flipProbability)
 Set probability of flipping an event about Z=0.
void setEnableSampling (bool enableSampling)
 Turn on or off sampling from the input ASCII file.
void setBufferFileName (const std::string &bufferFileName)
 Set the name of the binary buffer file, needed for sampling from a converted file.
void setDebugEnable (bool debug)
 A function to turn on or off debug print out.

Protected Types

enum  enumGeneratorProcessIds {
  MARS_READ , MARS_SINGLE , MARS_SHOWER , FLUKA_READ ,
  FLUKA_SINGLE , FLUKA_SHOWER
}
 An enum of generator process identifiers. More...
enum  enumCounters { TOT_READ , TOT_AFTER , TOT_GEN , NUM_COUNTERS }
 An enum of particle counter indices. More...

Protected Member Functions

virtual int readEvent (std::vector< BeamHaloParticle > *beamHaloEvent, CLHEP::HepRandomEngine *engine)=0
 A function to read one event in a simplified format.
virtual int readParticle (BeamHaloParticle *beamHaloParticle)=0
 A function to read one particle from the input ASCII file.
bool flipEvent (CLHEP::HepRandomEngine *engine)
 A member function to check if the event should be flipped.
int convertEvent (std::vector< BeamHaloParticle > *beamHaloEvent, HepMC::GenEvent *evt, CLHEP::HepRandomEngine *engine)
 A member function to convert a vector of beam halo particles into a GenEvent.

Protected Attributes

const HepPDT::ParticleDataTable * m_particleTable
 A pointer to the particle data table.
std::string m_inputFile
 Input file name.
double m_interfacePlane
 The position of the interface plane in mm.
bool m_enableFlip
 Flag for flipping event.
float m_flipProbability
 Flip probability.
bool m_enableSampling
 Flag to enable or disable sampling.
std::string m_bufferFileName
 The name of the binary buffer file, needed for sampling from a converted file.
std::unique_ptr< BeamHaloParticleBufferm_beamHaloParticleBuffer
 Binary particle buffer for caching converted events.
std::unique_ptr< AsciiInputm_asciiInput
 A pointer to an AsciiInput object, used to read data from the Ascii input file.
std::unique_ptr< BeamHaloGeneratorSettingsm_beamHaloGeneratorSettings
 A pointer to a BeamHaloGeneratorSettings object used to filter particles.
long m_eventNumber
 A data member to count the event number.
long m_counters [NUM_COUNTERS]
 Particle or event counters of dimension enumCounters.
double m_wsums [NUM_COUNTERS]
 Sum of weights of particles or events of dimension enumCounters.
bool m_debug
 A flat to turn on or off debug print out.

Private Attributes

std::vector< std::string > m_generatorSettings
 A vector of generator settings in string command form.

Detailed Description

An abstract base class to provide generic beam halo generator functionality.

Author
W. H. Bell W.Bel.nosp@m.l@ce.nosp@m.rn.ch

Definition at line 32 of file BeamHaloGenerator.h.

Member Enumeration Documentation

◆ enumCounters

An enum of particle counter indices.

Enumerator
TOT_READ 
TOT_AFTER 
TOT_GEN 
NUM_COUNTERS 

Definition at line 142 of file BeamHaloGenerator.h.

◆ enumGeneratorProcessIds

An enum of generator process identifiers.

_READ corresponds to no sampling and therefore simple conversion from input file to HepMC record.

Enumerator
MARS_READ 
MARS_SINGLE 
MARS_SHOWER 
FLUKA_READ 
FLUKA_SINGLE 
FLUKA_SHOWER 

Definition at line 87 of file BeamHaloGenerator.h.

Constructor & Destructor Documentation

◆ BeamHaloGenerator()

BeamHaloGenerator::BeamHaloGenerator ( const HepPDT::ParticleDataTable * particleTable,
const std::string & inputFile,
const std::vector< std::string > & generatorSettings )

Definition at line 16 of file BeamHaloGenerator.cxx.

18 :
19 m_particleTable(particleTable),
20 m_inputFile(inputFile),
22 m_enableFlip(false),
24 m_enableSampling(false),
25 m_bufferFileName("buffer.bin"),
27 m_debug(false),
28 m_generatorSettings(generatorSettings) {
29 for(int i=0;i<NUM_COUNTERS;i++) {
30 m_counters[i] = 0;
31 m_wsums[i] = 0.;
32 }
33}
long m_eventNumber
A data member to count the event number.
double m_wsums[NUM_COUNTERS]
Sum of weights of particles or events of dimension enumCounters.
float m_flipProbability
Flip probability.
std::string m_inputFile
Input file name.
bool m_enableFlip
Flag for flipping event.
long m_counters[NUM_COUNTERS]
Particle or event counters of dimension enumCounters.
std::vector< std::string > m_generatorSettings
A vector of generator settings in string command form.
bool m_debug
A flat to turn on or off debug print out.
bool m_enableSampling
Flag to enable or disable sampling.
const HepPDT::ParticleDataTable * m_particleTable
A pointer to the particle data table.
double m_interfacePlane
The position of the interface plane in mm.
std::string m_bufferFileName
The name of the binary buffer file, needed for sampling from a converted file.

◆ ~BeamHaloGenerator()

BeamHaloGenerator::~BeamHaloGenerator ( )
virtual

Definition at line 37 of file BeamHaloGenerator.cxx.

37 {
38}

Member Function Documentation

◆ convertEvent()

int BeamHaloGenerator::convertEvent ( std::vector< BeamHaloParticle > * beamHaloEvent,
HepMC::GenEvent * evt,
CLHEP::HepRandomEngine * engine )
protected

A member function to convert a vector of beam halo particles into a GenEvent.

Definition at line 88 of file BeamHaloGenerator.cxx.

90 {
91 double pz;
92 bool flipFlag = flipEvent(engine);
93 HepMC::GenParticlePtr genParticle;
94 HepMC::GenVertexPtr genVertex;
95 const double c = 2.99792458E+11; // speed of light in mm/s
96
97 std::vector<BeamHaloParticle>::iterator itr = beamHaloEvent->begin();
98 std::vector<BeamHaloParticle>::iterator itr_end = beamHaloEvent->end();
99
100 if(itr != itr_end) {
101 // Take the primary information from the first particle in the event. If the event is more than one
102 // particle, the particles come from the same primary proton.
103 evt->weights().push_back((*itr).weight()); // event weight
104 HepMC::FourVector primaryPosition = (*itr).positionAtPrimary();
105 evt->weights().push_back(primaryPosition.x()); // Starting position of primary particle: x position (cm)
106 evt->weights().push_back(primaryPosition.y()); // Starting position of primary particle: y position (cm)
107 evt->weights().push_back(primaryPosition.z()); // Starting position of primary particle: z position (cm)
108 evt->weights().push_back(primaryPosition.t()); // Starting position of primary particle: t (s)
109 }
110
111 // Append each particle to the GenEvent
112 for(;itr!=itr_end;++itr) {
113 HepMC::FourVector position = (*itr).positionAtScoringPlane();
114 HepMC::FourVector fourVector = (*itr).fourVector();
115
116 if(!flipFlag) {
117 genParticle = HepMC::newGenParticlePtr(fourVector,
118 (*itr).pdgId(),
119 1);
120
121 genVertex = HepMC::newGenVertexPtr(HepMC::FourVector(position.x(),
122 position.y(),
123 (position.z() + m_interfacePlane),
124 c*position.t() -1.0*std::fabs(m_interfacePlane)),
125 1);
126 }
127 else {
128 pz = fourVector.pz();
129 fourVector.setPz(-pz);
130
131 genParticle = HepMC::newGenParticlePtr(fourVector,
132 (*itr).pdgId(),
133 1);
134
135 genVertex = HepMC::newGenVertexPtr(HepMC::FourVector(position.x(),
136 position.y(),
137 -1.0*(position.z() + m_interfacePlane),
138 c*position.t() -1.0*std::fabs(m_interfacePlane)),
139 1);
140 }
141
142 genVertex->add_particle_out(genParticle);
143 evt->add_vertex(genVertex);
144 }
145
146 return 0;
147}
bool flipEvent(CLHEP::HepRandomEngine *engine)
A member function to check if the event should be flipped.
HepMC::GenVertex * GenVertexPtr
Definition GenVertex.h:59
GenVertexPtr newGenVertexPtr(const HepMC::FourVector &pos=HepMC::FourVector(0.0, 0.0, 0.0, 0.0), const int i=0)
Definition GenVertex.h:64
GenParticlePtr newGenParticlePtr(const HepMC::FourVector &mom=HepMC::FourVector(0.0, 0.0, 0.0, 0.0), int pid=0, int status=0)
Definition GenParticle.h:39
GenParticle * GenParticlePtr
Definition GenParticle.h:37

◆ fillEvt()

virtual int BeamHaloGenerator::fillEvt ( HepMC::GenEvent * evt,
CLHEP::HepRandomEngine * engine )
pure virtual

A function to create one event in HepMC format.

Implemented in FlukaHaloGenerator, and MarsHaloGenerator.

◆ flipEvent()

bool BeamHaloGenerator::flipEvent ( CLHEP::HepRandomEngine * engine)
protected

A member function to check if the event should be flipped.

Definition at line 74 of file BeamHaloGenerator.cxx.

74 {
75 if(!m_enableFlip) return false;
76
77 // Check to see if the event should be flipped or not
78 if(CLHEP::RandFlat::shoot(engine) <= m_flipProbability) {
79 if(m_debug) std::cout << "Debug: Flipping event" << std::endl;
80 return true;
81 }
82
83 return false;
84}

◆ genFinalize()

int BeamHaloGenerator::genFinalize ( )
virtual

A function to finalise the generator.

Reimplemented in FlukaHaloGenerator, and MarsHaloGenerator.

Definition at line 51 of file BeamHaloGenerator.cxx.

51 {
52
53 std::cout << "=================================================================" << std::endl;
54 std::cout << "| | Total Number | Total Weight |" << std::endl;
55 std::cout << "|---------------------------------------------------------------|" << std::endl;
56 std::cout << "| Particles read from input file | ";
57 std::cout << std::setw(12) << m_counters[TOT_READ] << " | ";
58 std::cout << std::setw(12) << std::setprecision(8) << m_wsums[TOT_READ] << " |" << std::endl;
59 std::cout << "| Particles after cuts | ";
60 std::cout << std::setw(12) << m_counters[TOT_AFTER] << " | ";
61 std::cout << std::setw(12) << std::setprecision(8) << m_wsums[TOT_AFTER] << " |" << std::endl;
62 std::cout << "| Particles generated | ";
63 std::cout << std::setw(12) << m_counters[TOT_GEN] << " | ";
64 std::cout << std::setw(12) << std::setprecision(8) << m_wsums[TOT_GEN] << " |" << std::endl;
65 std::cout << "=================================================================" << std::endl;
66
67 std::cout << std::setprecision(6);
68
69 return 0;
70}

◆ genInitialize()

int BeamHaloGenerator::genInitialize ( )
virtual

A function to initialise the generator.

Reimplemented in FlukaHaloGenerator, and MarsHaloGenerator.

Definition at line 42 of file BeamHaloGenerator.cxx.

42 {
43 m_asciiInput = std::make_unique<AsciiInput>(m_inputFile);
44 m_beamHaloGeneratorSettings = std::make_unique< BeamHaloGeneratorSettings>(m_generatorSettings);
45
46 return 0;
47}
std::unique_ptr< AsciiInput > m_asciiInput
A pointer to an AsciiInput object, used to read data from the Ascii input file.
std::unique_ptr< BeamHaloGeneratorSettings > m_beamHaloGeneratorSettings
A pointer to a BeamHaloGeneratorSettings object used to filter particles.

◆ readEvent()

virtual int BeamHaloGenerator::readEvent ( std::vector< BeamHaloParticle > * beamHaloEvent,
CLHEP::HepRandomEngine * engine )
protectedpure virtual

A function to read one event in a simplified format.

Implemented in FlukaHaloGenerator, and MarsHaloGenerator.

◆ readParticle()

virtual int BeamHaloGenerator::readParticle ( BeamHaloParticle * beamHaloParticle)
protectedpure virtual

A function to read one particle from the input ASCII file.

Implemented in FlukaHaloGenerator, and MarsHaloGenerator.

◆ setBufferFileName()

void BeamHaloGenerator::setBufferFileName ( const std::string & bufferFileName)
inline

Set the name of the binary buffer file, needed for sampling from a converted file.

Definition at line 77 of file BeamHaloGenerator.h.

77{ m_bufferFileName = bufferFileName; }

◆ setDebugEnable()

void BeamHaloGenerator::setDebugEnable ( bool debug)
inline

A function to turn on or off debug print out.

Definition at line 80 of file BeamHaloGenerator.h.

80{ m_debug = debug; }
const bool debug

◆ setEnableFlip()

void BeamHaloGenerator::setEnableFlip ( bool enableFlip)
inline

Turn on or off the event flipping code.

Definition at line 67 of file BeamHaloGenerator.h.

67{ m_enableFlip = enableFlip; }

◆ setEnableSampling()

void BeamHaloGenerator::setEnableSampling ( bool enableSampling)
inline

Turn on or off sampling from the input ASCII file.

Definition at line 73 of file BeamHaloGenerator.h.

73{ m_enableSampling = enableSampling; }

◆ setFlipProbability()

void BeamHaloGenerator::setFlipProbability ( float flipProbability)
inline

Set probability of flipping an event about Z=0.

Definition at line 70 of file BeamHaloGenerator.h.

◆ setInterfacePlane()

void BeamHaloGenerator::setInterfacePlane ( double interfacePlane)
inline

Set the position of the interface plane in mm.

Definition at line 64 of file BeamHaloGenerator.h.

Member Data Documentation

◆ m_asciiInput

std::unique_ptr<AsciiInput> BeamHaloGenerator::m_asciiInput
protected

A pointer to an AsciiInput object, used to read data from the Ascii input file.

Definition at line 132 of file BeamHaloGenerator.h.

◆ m_beamHaloGeneratorSettings

std::unique_ptr<BeamHaloGeneratorSettings> BeamHaloGenerator::m_beamHaloGeneratorSettings
protected

A pointer to a BeamHaloGeneratorSettings object used to filter particles.

Definition at line 136 of file BeamHaloGenerator.h.

◆ m_beamHaloParticleBuffer

std::unique_ptr<BeamHaloParticleBuffer> BeamHaloGenerator::m_beamHaloParticleBuffer
protected

Binary particle buffer for caching converted events.

Definition at line 128 of file BeamHaloGenerator.h.

◆ m_bufferFileName

std::string BeamHaloGenerator::m_bufferFileName
protected

The name of the binary buffer file, needed for sampling from a converted file.

Definition at line 125 of file BeamHaloGenerator.h.

◆ m_counters

long BeamHaloGenerator::m_counters[NUM_COUNTERS]
protected

Particle or event counters of dimension enumCounters.

Definition at line 149 of file BeamHaloGenerator.h.

◆ m_debug

bool BeamHaloGenerator::m_debug
protected

A flat to turn on or off debug print out.

Definition at line 155 of file BeamHaloGenerator.h.

◆ m_enableFlip

bool BeamHaloGenerator::m_enableFlip
protected

Flag for flipping event.

Definition at line 115 of file BeamHaloGenerator.h.

◆ m_enableSampling

bool BeamHaloGenerator::m_enableSampling
protected

Flag to enable or disable sampling.

Definition at line 121 of file BeamHaloGenerator.h.

◆ m_eventNumber

long BeamHaloGenerator::m_eventNumber
protected

A data member to count the event number.

Definition at line 139 of file BeamHaloGenerator.h.

◆ m_flipProbability

float BeamHaloGenerator::m_flipProbability
protected

Flip probability.

Definition at line 118 of file BeamHaloGenerator.h.

◆ m_generatorSettings

std::vector<std::string> BeamHaloGenerator::m_generatorSettings
private

A vector of generator settings in string command form.

Definition at line 160 of file BeamHaloGenerator.h.

◆ m_inputFile

std::string BeamHaloGenerator::m_inputFile
protected

Input file name.

Definition at line 109 of file BeamHaloGenerator.h.

◆ m_interfacePlane

double BeamHaloGenerator::m_interfacePlane
protected

The position of the interface plane in mm.

Definition at line 112 of file BeamHaloGenerator.h.

◆ m_particleTable

const HepPDT::ParticleDataTable* BeamHaloGenerator::m_particleTable
protected

A pointer to the particle data table.

Definition at line 106 of file BeamHaloGenerator.h.

◆ m_wsums

double BeamHaloGenerator::m_wsums[NUM_COUNTERS]
protected

Sum of weights of particles or events of dimension enumCounters.

Definition at line 152 of file BeamHaloGenerator.h.


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