ATLAS Offline Software
Loading...
Searching...
No Matches
BeamHaloGeneratorAlg.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
9
10#include "GaudiKernel/MsgStream.h"
11#include "AtlasHepMC/GenEvent.h"
13#include <cmath>
14
15//--------------------------------------------------------------------------
16
17BeamHaloGeneratorAlg::BeamHaloGeneratorAlg(const std::string& name, ISvcLocator* svcLocator)
18 : GenModule(name,svcLocator)
19{
20 for (int i=0;i<NPLOTS;i++) {
21 m_validationPlots[i] = nullptr;
22 }
23}
24
25//--------------------------------------------------------------------------
26
28 ATH_MSG_INFO("Initialising this BeamHaloGeneratorAlg: " << name());
29
30 ATH_MSG_INFO("================ Generator Settings =================");
31 ATH_MSG_INFO(" input type ------------------ " << m_inputTypeStr);
32 ATH_MSG_INFO(" input file ------------------ " << m_inputFile);
33 ATH_MSG_INFO(" interface plane ------------- " << m_interfacePlane << " (mm)");
34 ATH_MSG_INFO(" enable flip about z=0 ------- " << (m_enableFlip ? "True" : "False"));
35 ATH_MSG_INFO(" flip probability ----------- " << m_flipProbability);
36 ATH_MSG_INFO(" enable sampling ------------- " << (m_enableSampling ? "True" : "False"));
37 ATH_MSG_INFO(" binary buffer file ---------- " << m_bufferFileName);
38 std::vector<std::string>::iterator itr = m_generatorSettings.begin();
39 std::vector<std::string>::iterator itr_end = m_generatorSettings.end();
40 for(;itr!=itr_end;++itr) ATH_MSG_INFO((*itr));
41 ATH_MSG_INFO(" Produce monitoring plots ---- " << (m_doMonitoringPlots ? "True" : "False"));
42 ATH_MSG_INFO(" Random stream name ---------- " << m_randomStream);
43 ATH_MSG_INFO("=====================================================");
44
45 // Check the flip probability
46 if(m_flipProbability <= 0) {
47 ATH_MSG_INFO("Event flipping will be disabled.");
48 m_enableFlip = false;
49 }
50 else if(m_flipProbability > 1) {
51 ATH_MSG_WARNING("Flip probability " << m_flipProbability << " is out of range. Event flipping will be disabled.");
52 m_enableFlip = false;
53 }
54 else {
55 m_enableFlip = true;
56 }
57
58 // Retrieve pointer to THistSvc if monitoring plots are requested
61 ATH_CHECK(m_tHistSvc.retrieve());
62
63 // Create log10 bins for x-axis of E plots.
64 const Int_t nbins_E = 60;
65 double xmin = 1e-2;
66 double xmax = 3.5e3;
67 double logxmin = std::log10(xmin);
68 double logxmax = std::log10(xmax);
69 double binwidth = (logxmax-logxmin)/nbins_E;
70 Double_t xbins[nbins_E+1];
71 xbins[0] = xmin;
72 for (Int_t i=1;i<=nbins_E;i++) {
73 xbins[i] = xmin + std::pow(10,logxmin+i*binwidth);
74 }
75
76 // Create the monitoring plots
77 m_validationPlots[PRI_R] = new TH1F("primaryR",";Radius [m];Events/[m]",100,0., 0.2);
78 m_validationPlots[PRI_Z] = new TH1F("primaryZ",";z [m];Events/[m]",1100,0., 550.); // Beam gas and TCT
79 m_validationPlots[PRI_Z_TCT] = new TH1F("primaryZ_TCT",";z [m];Events/[m]",100,144., 149.); // TCT region
80 m_validationPlots[SP_R_ALL] = new TH1F("scoringPlaneR_all",";Radius [m];Particles/[m]",240,0., 24.);
81 m_validationPlots[SP_E_ALL] = new TH1F("scoringPlaneE_all",";Energy [GeV];Particles/[GeV]",nbins_E,xbins);
82 m_validationPlots[SP_PZ_ALL] = new TH1F("scoringPlanePz_all",";p_{z} [GeV];Particles/[GeV]",350,0., 3500.);
83 m_validationPlots[SP_PT_ALL] = new TH1F("scoringPlanePt_all",";p_{T} [GeV];Particles/[GeV]",500,0., 50.);
84 m_validationPlots[SP_R_PROTONS] = new TH1F("scoringPlaneR_protons",";Radius [m];Particles/[m]",240,0., 24.);
85 m_validationPlots[SP_E_PROTONS] = new TH1F("scoringPlaneE_protons",";Energy [GeV];Particles/[GeV]",nbins_E,xbins);
86 m_validationPlots[SP_PZ_PROTONS] = new TH1F("scoringPlanePz_protons",";p_{z} [GeV];Particles/[GeV]",350,0., 3500.);
87 m_validationPlots[SP_PT_PROTONS] = new TH1F("scoringPlanePt_protons",";p_{T} [GeV];Particles/[GeV]",500,0., 50.);
88 m_validationPlots[SP_R_MUONS] = new TH1F("scoringPlaneR_muons",";Radius [m];Particles/[m]",240,0., 24.);
89 m_validationPlots[SP_E_MUONS] = new TH1F("scoringPlaneE_muons",";Energy [GeV];Particles/[GeV]",nbins_E,xbins);
90 m_validationPlots[SP_PZ_MUONS] = new TH1F("scoringPlanePz_muons",";p_{z} [GeV];Particles/[GeV]",350,0., 3500.);
91 m_validationPlots[SP_PT_MUONS] = new TH1F("scoringPlanePt_muons",";p_{T} [GeV];Particles/[GeV]",500,0., 50.);
92
93 ATH_CHECK(m_tHistSvc->regHist("/BeamHalo/primaryR", m_validationPlots[PRI_R]));
94 ATH_CHECK(m_tHistSvc->regHist("/BeamHalo/primaryZ", m_validationPlots[PRI_Z]));
95 ATH_CHECK(m_tHistSvc->regHist("/BeamHalo/primaryZ_TCT", m_validationPlots[PRI_Z_TCT]));
96 ATH_CHECK(m_tHistSvc->regHist("/BeamHalo/scoringPlaneR_all", m_validationPlots[SP_R_ALL]));
97 ATH_CHECK(m_tHistSvc->regHist("/BeamHalo/scoringPlaneE_all", m_validationPlots[SP_E_ALL]));
98 ATH_CHECK(m_tHistSvc->regHist("/BeamHalo/scoringPlanePz_all", m_validationPlots[SP_PZ_ALL]));
99 ATH_CHECK(m_tHistSvc->regHist("/BeamHalo/scoringPlanePt_all", m_validationPlots[SP_PT_ALL]));
100 ATH_CHECK(m_tHistSvc->regHist("/BeamHalo/scoringPlaneR_protons", m_validationPlots[SP_R_PROTONS]));
101 ATH_CHECK(m_tHistSvc->regHist("/BeamHalo/scoringPlaneE_protons", m_validationPlots[SP_E_PROTONS]));
102 ATH_CHECK(m_tHistSvc->regHist("/BeamHalo/scoringPlanePz_protons", m_validationPlots[SP_PZ_PROTONS]));
103 ATH_CHECK(m_tHistSvc->regHist("/BeamHalo/scoringPlanePt_protons", m_validationPlots[SP_PT_PROTONS]));
104 ATH_CHECK(m_tHistSvc->regHist("/BeamHalo/scoringPlaneR_muons", m_validationPlots[SP_R_MUONS]));
105 ATH_CHECK(m_tHistSvc->regHist("/BeamHalo/scoringPlaneE_muons", m_validationPlots[SP_E_MUONS]));
106 ATH_CHECK(m_tHistSvc->regHist("/BeamHalo/scoringPlanePz_muons", m_validationPlots[SP_PZ_MUONS]));
107 ATH_CHECK(m_tHistSvc->regHist("/BeamHalo/scoringPlanePt_muons", m_validationPlots[SP_PT_MUONS]));
108 }
109
110 // Check the input type string
111 if (m_inputTypeStr == "MARS-NM") {
115 }
116 else if (m_inputTypeStr == "FLUKA-VT") {
120 }
121 else if (m_inputTypeStr == "FLUKA-RB") {
125 }
126 else {
128 ATH_MSG_FATAL("Input type " << m_inputTypeStr << " is not known. Available types are: MARS-NM or FLUKA-VT");
129 return StatusCode::FAILURE;
130 }
131
132 // Set the options using those specified to this algorithm
133 m_beamHaloGenerator->setInterfacePlane(m_interfacePlane);
134 m_beamHaloGenerator->setEnableFlip(m_enableFlip);
135 m_beamHaloGenerator->setFlipProbability(m_flipProbability);
136 m_beamHaloGenerator->setEnableSampling(m_enableSampling);
137 m_beamHaloGenerator->setBufferFileName(m_bufferFileName);
138 m_beamHaloGenerator->setDebugEnable(msgLvl(MSG::DEBUG));
139
140 // Initialise the generator.
141 if(m_beamHaloGenerator->genInitialize() != 0) return StatusCode::FAILURE;
142
143 return StatusCode::SUCCESS;
144}
145
146//---------------------------------------------------------------------------
147
149
150 const EventContext& ctx = Gaudi::Hive::currentContext();
151 CLHEP::HepRandomEngine* rndmEngine = this->getRandomEngine(m_randomStream, ctx);
152
153
154 // Clear the event ready for it to be filled with the next beam halo event.
155 m_evt.clear();
156
157 // Fill an event with particles from the converted ASCII input file.
158 // (If the end of file has already been reached this will return a
159 // non-zero value.)
160 if(m_beamHaloGenerator->fillEvt(&m_evt, rndmEngine) != 0) return StatusCode::FAILURE;
161
162 // Fill monitoring plots if requested
164 auto weightContainer = m_evt.weights();
165 if(weightContainer.size() != 5) {
166 ATH_MSG_WARNING("The number of weights for this event is not equal to 5.");
167 return StatusCode::SUCCESS;
168 }
169 double weight = weightContainer[0];
170 HepMC::FourVector primaryPosition(weightContainer[1],weightContainer[2],weightContainer[3],weightContainer[4]);
171 m_validationPlots[PRI_R]->Fill(primaryPosition.perp()/1000.0,weight);
172 m_validationPlots[PRI_Z]->Fill(weightContainer[3]/1000.0,weight);
173 m_validationPlots[PRI_Z_TCT]->Fill(std::fabs(weightContainer[3]/1000.0),weight);
174
175 double values[4];
176 int pdgId;
177 for (auto hepmc_part: m_evt) {
178 auto prodVertex = hepmc_part->production_vertex();
179 if(!prodVertex) continue;
180
181 // Store the values for use in the if conditions that follow
182 values[0] = prodVertex->position().perp()/1000.;
183 values[1] = hepmc_part->momentum().e()/1000.;
184 values[2] = hepmc_part->momentum().pz()/1000.;
185 values[3] = hepmc_part->momentum().perp()/1000.;
186
187 pdgId = hepmc_part->pdg_id(); if(pdgId<0) pdgId = -pdgId;
188 m_validationPlots[SP_R_ALL]->Fill(values[0],weight);
189 m_validationPlots[SP_E_ALL]->Fill(values[1],weight);
190 m_validationPlots[SP_PZ_ALL]->Fill(values[2],weight);
191 m_validationPlots[SP_PT_ALL]->Fill(values[3],weight);
192 if(pdgId == 2212) {
193 m_validationPlots[SP_R_PROTONS]->Fill(values[0],weight);
194 m_validationPlots[SP_E_PROTONS]->Fill(values[1],weight);
195 m_validationPlots[SP_PZ_PROTONS]->Fill(values[2],weight);
196 m_validationPlots[SP_PT_PROTONS]->Fill(values[3],weight);
197 }
198 else if(pdgId == 13) {
199 m_validationPlots[SP_R_MUONS]->Fill(values[0],weight);
200 m_validationPlots[SP_E_MUONS]->Fill(values[1],weight);
201 m_validationPlots[SP_PZ_MUONS]->Fill(values[2],weight);
202 m_validationPlots[SP_PT_MUONS]->Fill(values[3],weight);
203 }
204 }
205 }
206
207 return StatusCode::SUCCESS;
208}
209
210//---------------------------------------------------------------------------
211
212StatusCode BeamHaloGeneratorAlg::fillEvt(HepMC::GenEvent* event) {
213 *event = m_evt;
214 return StatusCode::SUCCESS;
215}
216
217//---------------------------------------------------------------------------
218
220 if(m_beamHaloGenerator->genFinalize() != 0) return StatusCode::FAILURE;
221 return StatusCode::SUCCESS;
222}
223
224//---------------------------------------------------------------------------
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
bool msgLvl(const MSG::Level lvl) const
virtual StatusCode genInitialize()
For initializing the generator, if required.
BooleanProperty m_enableFlip
Flag for flipping event.
ServiceHandle< ITHistSvc > m_tHistSvc
A pointer to the THist service for validation plots.
DoubleProperty m_interfacePlane
The position of the interface plane in mm.
virtual StatusCode genFinalize()
For finalising the generator, if required.
StringProperty m_inputFile
Input file name.
BooleanProperty m_doMonitoringPlots
A flag to allow monitoring plots to be turned on or off.
StringProperty m_inputTypeStr
Input file type and therefore associated beam halo generator that should be used.
virtual StatusCode fillEvt(HepMC::GenEvent *evt)
Fill the GenEvent pointer with the contents of the GenEvent cache.
TH1F * m_validationPlots[NPLOTS]
An array of TH1F pointers for validation plots.
StringProperty m_bufferFileName
The name of the binary buffer file, needed for sampling from a converted file.
virtual StatusCode callGenerator()
Read one event from the selected input and convert it into GenEvent format.
FloatProperty m_flipProbability
Flip probability.
StringProperty m_randomStream
Name of the random number stream.
BeamHaloGenerator * m_beamHaloGenerator
A pointer to the beam halo generator.
StringArrayProperty m_generatorSettings
A vector of strings defining generator settings.
HepMC::GenEvent m_evt
An empty GenEvent to cache the generate output between callGenerator and fillEvt.
BooleanProperty m_enableSampling
Flag to enable or disable sampling.
BeamHaloGeneratorAlg(const std::string &name, ISvcLocator *svcLocator)
A class to provide conversion from a FLUKA format ASCII input record into HepMC format,...
const HepPDT::ParticleDataTable & particleTable() const
Get a particle data table.
Definition GenBase.h:118
GenModule(const std::string &name, ISvcLocator *pSvcLocator)
Constructor.
Definition GenModule.cxx:14
CLHEP::HepRandomEngine * getRandomEngine(const std::string &streamName, const EventContext &ctx) const
Definition GenModule.cxx:34
A class to provide conversion from a MARS format ASCII input record into HepMC format,...
double xmax
Definition listroot.cxx:61
double xmin
Definition listroot.cxx:60
bool binwidth
Definition listroot.cxx:58