ATLAS Offline Software
Loading...
Searching...
No Matches
FlukaHaloGenerator Class Reference

A class to provide conversion from a FLUKA format ASCII input record into HepMC format, with or without sampling. More...

#include <FlukaHaloGenerator.h>

Inheritance diagram for FlukaHaloGenerator:
Collaboration diagram for FlukaHaloGenerator:

Public Member Functions

 FlukaHaloGenerator (int type, const HepPDT::ParticleDataTable *particleTable, const std::string &inputFile, const std::vector< std::string > &generatorSettings)
virtual ~FlukaHaloGenerator ()=default
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)
 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)
 A function to read one event in a simplified format.
virtual int readParticle (BeamHaloParticle *beamHaloParticle)
 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

bool m_sameEvent
bool m_firstEvent
FlukaParticle m_flukaParticle
FlukaParticle m_lastFlukaParticle
std::vector< std::string > m_generatorSettings
 A vector of generator settings in string command form.

Detailed Description

A class to provide conversion from a FLUKA format ASCII input record into HepMC format, with or without sampling.

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

Definition at line 18 of file FlukaHaloGenerator.h.

Member Enumeration Documentation

◆ enumCounters

enum BeamHaloGenerator::enumCounters
protectedinherited

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

◆ FlukaHaloGenerator()

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

Definition at line 15 of file FlukaHaloGenerator.cxx.

18 :
19 BeamHaloGenerator(particleTable, inputFile, generatorSettings),
20 m_sameEvent(true),
21 m_firstEvent(true),
22 m_flukaParticle(type),
24}
BeamHaloGenerator(const HepPDT::ParticleDataTable *particleTable, const std::string &inputFile, const std::vector< std::string > &generatorSettings)
FlukaParticle m_lastFlukaParticle
FlukaParticle m_flukaParticle

◆ ~FlukaHaloGenerator()

virtual FlukaHaloGenerator::~FlukaHaloGenerator ( )
virtualdefault

Member Function Documentation

◆ convertEvent()

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

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}
double m_interfacePlane
The position of the interface plane in mm.
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()

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

A function to create one event in HepMC format.

Implements BeamHaloGenerator.

Definition at line 48 of file FlukaHaloGenerator.cxx.

49 {
50 std::vector<BeamHaloParticle> beamHaloEvent;
51 int ret_val;
52
53 // Read one FLUKA event which passes the generator settings.
54 if((ret_val = readEvent(&beamHaloEvent, engine)) != 0) return ret_val;
55
56 // Convert the particles to GenParticles and attach them to the
57 // event. Flip the event if needed.
58 if((ret_val = BeamHaloGenerator::convertEvent(&beamHaloEvent, evt, engine)) != 0) return ret_val;
59
60 // Set the event number
61 evt->set_event_number(m_eventNumber);
62
63 // Set the signal process
65
66 // Increment the event number
68
69 return 0;
70}
long m_eventNumber
A data member to count the event number.
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.
virtual int readEvent(std::vector< BeamHaloParticle > *beamHaloEvent, CLHEP::HepRandomEngine *engine)
A function to read one event in a simplified format.
void set_signal_process_id(GenEvent *e, const int i)
Definition GenEvent.h:641

◆ flipEvent()

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

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}
float m_flipProbability
Flip probability.
bool m_enableFlip
Flag for flipping event.
bool m_debug
A flat to turn on or off debug print out.

◆ genFinalize()

int FlukaHaloGenerator::genFinalize ( )
virtual

A function to finalise the generator.

Reimplemented from BeamHaloGenerator.

Definition at line 74 of file FlukaHaloGenerator.cxx.

74 {
75 int ret_val;
76
77 m_asciiInput->close();
78
79 // Finalise base class
80 if((ret_val = BeamHaloGenerator::genFinalize()) != 0) return ret_val;
81
82 return 0;
83}
std::unique_ptr< AsciiInput > m_asciiInput
A pointer to an AsciiInput object, used to read data from the Ascii input file.
virtual int genFinalize()
A function to finalise the generator.

◆ genInitialize()

int FlukaHaloGenerator::genInitialize ( )
virtual

A function to initialise the generator.

Reimplemented from BeamHaloGenerator.

Definition at line 29 of file FlukaHaloGenerator.cxx.

29 {
30 int ret_val;
31
32 // Initialise base class
33 if((ret_val = BeamHaloGenerator::genInitialize()) != 0) return ret_val;
34
35 if(m_asciiInput->open() != 0) {
36 std::cout << "Error: Could not open ascii input file " << m_inputFile << std::endl;
37 return 1;
38 }
39
40 std::cout << "Info: Reading ascii input file " << m_inputFile << std::endl;
41
42 return 0;
43}
std::string m_inputFile
Input file name.
virtual int genInitialize()
A function to initialise the generator.

◆ readEvent()

int FlukaHaloGenerator::readEvent ( std::vector< BeamHaloParticle > * beamHaloEvent,
CLHEP::HepRandomEngine * engine )
protectedvirtual

A function to read one event in a simplified format.

Implements BeamHaloGenerator.

Definition at line 95 of file FlukaHaloGenerator.cxx.

96 {
97 BeamHaloParticle beamHaloParticle;
98
99 // Clear the event
100 beamHaloEvent->clear();
101
102 // If there was a last event.
103 if(!m_firstEvent) {
104
105 // If the last event caused the same event flag to be set to false
106 // copy the last particle into the vector of those in this event.
107 if(!m_sameEvent) {
108 // Fill the BeamHaloParticle with the data in the FlukaParticle
109 if(beamHaloParticle.fill(m_particleTable, &m_lastFlukaParticle)) {
110 std::cout << "Error: Conversion from FlukaParticle to BeamHaloParticle failed." << std::endl;
111 return 1;
112 }
113 // Append the BeamHalo particle to the event.
114 beamHaloEvent->push_back(beamHaloParticle);
115
116 // Set the same event flag to enter the while loop to read the
117 // rest of this event.
118 m_sameEvent = true;
119 }
120 }
121
122 // Loop over the ascii input and read each particle until a new
123 // event is found or there are no more particles.
124 std::vector<std::string> row;
125 bool endOfFile = false;
126 while(m_sameEvent && !endOfFile) {
127 row = m_asciiInput->readRow(); // Read one line of the ascii file.11
128
129 if(row.size() == 0) {
130 endOfFile = true;
131 continue;
132 }
133
134 if(m_flukaParticle.read(&row)) { // Fill the particle from the string vector
135 continue;
136 }
137
138 // Check if there was a last particle.
139 if(!m_firstEvent) {
140
141 m_sameEvent = false;
142
143 // Check if the event id of the last particle is the same as this particle.
145 && (m_lastFlukaParticle.eventId() == m_flukaParticle.eventId())
146 && (m_lastFlukaParticle.partGenNum() == m_flukaParticle.partGenNum())) m_sameEvent = true;
147
149 && (m_lastFlukaParticle.eventId() == m_flukaParticle.eventId())) m_sameEvent = true;
150 }
151 else {
152 // For the first event.
153 m_firstEvent = false;
154 m_sameEvent = true;
155 }
156
157 // If this is the same event copy the particle into the vector for
158 // this event.
159 if(m_sameEvent) {
160
161 // Fill the BeamHaloParticle with the data in the FlukaParticle
162 if(beamHaloParticle.fill(m_particleTable, &m_flukaParticle)) {
163 std::cout << "Error: Conversion from FlukaParticle to BeamHaloParticle failed." << std::endl;
164 return 1;
165 }
166
167 // Append the BeamHalo particle to the event.
168 beamHaloEvent->push_back(beamHaloParticle);
169 }
170
171 // Copy this particle into the last particle.
173
174 }
175
176 if(beamHaloEvent->size() == 0) {
177 std::cout << "Error: No particles read from " << m_inputFile << std::endl;
178 return 1;
179 }
180
182 m_wsums[TOT_READ]+= 1.0;
183
184 // Check if one of the particles in the event passes the generator settings.
185 std::vector<BeamHaloParticle>::iterator itr = beamHaloEvent->begin();
186 std::vector<BeamHaloParticle>::iterator itr_end = beamHaloEvent->end();
187 bool passed = false;
188 for(;itr!=itr_end;++itr) {
189 beamHaloParticle = *itr;
190 // Check the generator settings for this particle.
191 if(m_beamHaloGeneratorSettings->checkParticle(&beamHaloParticle)) {
192 passed = true;
193 }
194 else {
195 if(m_debug) std::cout << "Debug: Particle fails generator settings cuts." << std::endl;
196 }
197 }
198
200 m_wsums[TOT_AFTER]+= 1.0;
202 m_wsums[TOT_GEN]+= 1.0;
203
204 // If all of the particles from this event fail read another event.
205 // If there are no more events this function will exit with a
206 // WARNING.
207 if(!passed) {
208 int ret_val;
209 if((ret_val = readEvent(beamHaloEvent, engine)) != 0) return ret_val;
210 }
211
212 return 0;
213}
bool passed(DecisionID id, const DecisionIDContainer &)
checks if required decision ID is in the set of IDs in the container
double m_wsums[NUM_COUNTERS]
Sum of weights of particles or events of dimension enumCounters.
long m_counters[NUM_COUNTERS]
Particle or event counters of dimension enumCounters.
std::unique_ptr< BeamHaloGeneratorSettings > m_beamHaloGeneratorSettings
A pointer to a BeamHaloGeneratorSettings object used to filter particles.
const HepPDT::ParticleDataTable * m_particleTable
A pointer to the particle data table.
int fill(const HepPDT::ParticleDataTable *particleDataTable, MarsParticle *marsParticle)
A function to fill the data members from an input MarsParticle object.
row
Appending html table to final .html summary file.

◆ readParticle()

int FlukaHaloGenerator::readParticle ( BeamHaloParticle * beamHaloParticle)
protectedvirtual

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

Implements BeamHaloGenerator.

Definition at line 87 of file FlukaHaloGenerator.cxx.

87 {
88
89 std::cout << "Warning: FlukaHaloGenerator::readParticle is not yet available: "<< beamHaloParticle << std::endl;
90 return 0;
91}

◆ setBufferFileName()

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

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; }
std::string m_bufferFileName
The name of the binary buffer file, needed for sampling from a converted file.

◆ setDebugEnable()

void BeamHaloGenerator::setDebugEnable ( bool debug)
inlineinherited

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)
inlineinherited

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)
inlineinherited

Turn on or off sampling from the input ASCII file.

Definition at line 73 of file BeamHaloGenerator.h.

73{ m_enableSampling = enableSampling; }
bool m_enableSampling
Flag to enable or disable sampling.

◆ setFlipProbability()

void BeamHaloGenerator::setFlipProbability ( float flipProbability)
inlineinherited

Set probability of flipping an event about Z=0.

Definition at line 70 of file BeamHaloGenerator.h.

◆ setInterfacePlane()

void BeamHaloGenerator::setInterfacePlane ( double interfacePlane)
inlineinherited

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
protectedinherited

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
protectedinherited

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
protectedinherited

Binary particle buffer for caching converted events.

Definition at line 128 of file BeamHaloGenerator.h.

◆ m_bufferFileName

std::string BeamHaloGenerator::m_bufferFileName
protectedinherited

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]
protectedinherited

Particle or event counters of dimension enumCounters.

Definition at line 149 of file BeamHaloGenerator.h.

◆ m_debug

bool BeamHaloGenerator::m_debug
protectedinherited

A flat to turn on or off debug print out.

Definition at line 155 of file BeamHaloGenerator.h.

◆ m_enableFlip

bool BeamHaloGenerator::m_enableFlip
protectedinherited

Flag for flipping event.

Definition at line 115 of file BeamHaloGenerator.h.

◆ m_enableSampling

bool BeamHaloGenerator::m_enableSampling
protectedinherited

Flag to enable or disable sampling.

Definition at line 121 of file BeamHaloGenerator.h.

◆ m_eventNumber

long BeamHaloGenerator::m_eventNumber
protectedinherited

A data member to count the event number.

Definition at line 139 of file BeamHaloGenerator.h.

◆ m_firstEvent

bool FlukaHaloGenerator::m_firstEvent
private

Definition at line 50 of file FlukaHaloGenerator.h.

◆ m_flipProbability

float BeamHaloGenerator::m_flipProbability
protectedinherited

Flip probability.

Definition at line 118 of file BeamHaloGenerator.h.

◆ m_flukaParticle

FlukaParticle FlukaHaloGenerator::m_flukaParticle
private

Definition at line 51 of file FlukaHaloGenerator.h.

◆ m_generatorSettings

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

A vector of generator settings in string command form.

Definition at line 160 of file BeamHaloGenerator.h.

◆ m_inputFile

std::string BeamHaloGenerator::m_inputFile
protectedinherited

Input file name.

Definition at line 109 of file BeamHaloGenerator.h.

◆ m_interfacePlane

double BeamHaloGenerator::m_interfacePlane
protectedinherited

The position of the interface plane in mm.

Definition at line 112 of file BeamHaloGenerator.h.

◆ m_lastFlukaParticle

FlukaParticle FlukaHaloGenerator::m_lastFlukaParticle
private

Definition at line 52 of file FlukaHaloGenerator.h.

◆ m_particleTable

const HepPDT::ParticleDataTable* BeamHaloGenerator::m_particleTable
protectedinherited

A pointer to the particle data table.

Definition at line 106 of file BeamHaloGenerator.h.

◆ m_sameEvent

bool FlukaHaloGenerator::m_sameEvent
private

Definition at line 49 of file FlukaHaloGenerator.h.

◆ m_wsums

double BeamHaloGenerator::m_wsums[NUM_COUNTERS]
protectedinherited

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: