ATLAS Offline Software
Loading...
Searching...
No Matches
BeamHaloGenerator.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
10#include "AtlasHepMC/GenEvent.h"
11#include "CLHEP/Random/RandFlat.h"
12#include "CLHEP/Units/PhysicalConstants.h"
13#include "TMath.h"
14#include <cmath>
15
16BeamHaloGenerator::BeamHaloGenerator(const HepPDT::ParticleDataTable* particleTable,
17 const std::string& inputFile,
18 const std::vector<std::string>& generatorSettings):
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}
34
35//------------------------------------------------------------------
36
39
40//------------------------------------------------------------------
41
43 m_asciiInput = std::make_unique<AsciiInput>(m_inputFile);
44 m_beamHaloGeneratorSettings = std::make_unique< BeamHaloGeneratorSettings>(m_generatorSettings);
45
46 return 0;
47}
48
49//------------------------------------------------------------------
50
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}
71
72//------------------------------------------------------------------
73
74bool BeamHaloGenerator::flipEvent(CLHEP::HepRandomEngine* engine) {
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}
85
86//------------------------------------------------------------------
87
88int BeamHaloGenerator::convertEvent(std::vector<BeamHaloParticle>* beamHaloEvent,
89 HepMC::GenEvent* evt,
90 CLHEP::HepRandomEngine* engine) {
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}
148
149//------------------------------------------------------------------
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.
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.
std::unique_ptr< AsciiInput > m_asciiInput
A pointer to an AsciiInput object, used to read data from the Ascii input file.
std::string m_inputFile
Input file name.
bool m_enableFlip
Flag for flipping event.
virtual int genInitialize()
A function to initialise the generator.
BeamHaloGenerator(const HepPDT::ParticleDataTable *particleTable, const std::string &inputFile, const std::vector< std::string > &generatorSettings)
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.
virtual int genFinalize()
A function to finalise the generator.
bool m_debug
A flat to turn on or off debug print out.
std::unique_ptr< BeamHaloGeneratorSettings > m_beamHaloGeneratorSettings
A pointer to a BeamHaloGeneratorSettings object used to filter particles.
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.
bool flipEvent(CLHEP::HepRandomEngine *engine)
A member function to check if the event should be flipped.
std::string m_bufferFileName
The name of the binary buffer file, needed for sampling from a converted file.
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