ATLAS Offline Software
Loading...
Searching...
No Matches
ALFA_BeamTransport.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6// ALFA_BeamTransport.cxx, (c) ATLAS Detector software
8
11
12#include "GaudiKernel/ServiceHandle.h"
13
14// FrameWork includes
15#include "GaudiKernel/ITHistSvc.h"
16#include "Gaudi/Property.h"
17
18
19#include "AtlasHepMC/GenEvent.h"
23
24//ROOT headers
25#include "TFile.h"
26#include "TH1.h"
27#include "TString.h"
28
29
30#include <iostream>
31#include <vector>
32
33//================ Constructor =================================================
34
35ALFA_BeamTransport::ALFA_BeamTransport(const std::string& name, ISvcLocator* pSvcLocator)
36 :
37 AthAlgorithm(name,pSvcLocator)
38{
39 declareProperty("ConfDir", m_FPConfig.ConfDir="./config");
40 declareProperty("UseALFA", m_FPConfig.UseALFA=true);
41 declareProperty("Debug", m_WriteDebugOutput=false);
42 declareProperty("IP",m_FPConfig.IP = 1);
43 declareProperty("UseAper",m_FPConfig.useaper = true);
44 declareProperty("apermb",m_FPConfig.apermb = 0.);
45 declareProperty("xcol1",m_FPConfig.xcol1 = 999.);//999e-4 but we need open collimators
46 declareProperty("xcol2",m_FPConfig.xcol2 = 999.);//999e-4
47 declareProperty("BeamEnergy",m_FPConfig.pbeam0 = 7000.);//beam energy
48 declareProperty("RPDistance",m_FPConfig.RPDistance = 237.398);
49 declareProperty("absZMagMax",m_FPConfig.absZMagMax = 500.);//dont read magnets after this value
50
51 //Cuts for production
52 declareProperty("EtaCut",m_EtaCut = 8.);
53 declareProperty("XiCut",m_XiCut = 0.2);
54
55 declareProperty("FPOutputBeam1", m_FPOutputBeam1="Beam1");
56 declareProperty("FPOutputBeam2", m_FPOutputBeam2="Beam2");
57
58}
59
61
62
64{
65 ATH_MSG_INFO("Initializing " << name() << "...");
66 ATH_CHECK(m_eventInfoKey.initialize());
67 ATH_CHECK(m_MCKey.initialize());
68 // open files for beamtracking output
69 // string in char
70 char* cBeam1Name = new char[m_FPOutputBeam1.length()];
71 strcpy(cBeam1Name, m_FPOutputBeam1.c_str());
72 char* cBeam2Name = new char[m_FPOutputBeam2.length()];
73 strcpy(cBeam2Name, m_FPOutputBeam2.c_str());
74
76 m_FileBeam1.open(cBeam1Name);
77 m_FileBeam2.open(cBeam2Name);
78 }
79
80 //delete char arrays
81 delete[] cBeam1Name;
82 delete[] cBeam2Name;
83
84 //Initialisation of FPTracker
85 m_BeamTracker.ALFA_BeamTrack::initialize(m_FPConfig);
86
87 //-----------------------------------------------------------------------------------------
88 ATH_MSG_INFO( name() << " initialize()" );
89 ATH_MSG_INFO( "initialize() successful in " << name() );
90 return StatusCode::SUCCESS;
91}
92
93//================ Finalisation =================================================
94
96{
97 // close outputfiles for FPTracker
99 m_FileBeam1.close();
100 m_FileBeam2.close();
101 }
102 // Code entered here will be executed once at the end of the program run.
103 return StatusCode::SUCCESS;
104}
105
106//================ Execution ====================================================
107
109{
110 // Code entered here will be executed once per event
111 ATH_MSG_DEBUG ("Executing " << name() << "...");
112 //-----------------------------------------------------------------------------------------
113 int run_number;
114 uint64_t evt_number;
115
116 //Set particles counter to zero!
117 //counter for particles marked as outgoing in event record
118 m_pcount=0;
119 //counter for particles marked as incomming in event record
120 m_pint=0;
121
122 //Load event info
123 SG::ReadHandle<xAOD::EventInfo> eventInfo (m_eventInfoKey,getContext());
124 if(!eventInfo.isValid()) {
125 ATH_MSG_ERROR("Could not retrieve EventInfo");
126 return StatusCode::FAILURE;
127 }
128 else {
129 run_number=eventInfo->runNumber();
130 evt_number=eventInfo->eventNumber();
131
132 ATH_MSG_DEBUG("run: " << run_number << " event: " << evt_number);
133 }
134
135 SG::ReadHandle<McEventCollection> mcColl(m_MCKey, getContext());
136 if (!mcColl.isValid()) {
137 ATH_MSG_WARNING("Could not retrieve McEventCollection: " << m_MCKey);
138 return StatusCode::SUCCESS;
139 }
140 // Loop over all events in McEventCollectio
141 ATH_MSG_INFO("successful load of HepMC info");
142 for (const HepMC::GenEvent* itr : *mcColl){
143 HepMC::GenEvent evt = (*itr);
144 MC::MeVToGeV(&evt); //Only scales momenta and masses
145 HepMC::Print::line(std::cout, evt);
146
147 // Select final state particle from event generator
148 // set event status !=1 (final state)
149 // fill member variables with particle data
150 // ALFA_BeamTransport::SelectParticles(evt);
151
152 // There should be a check if the first particle is really the
153 // particle for beam 1!!!!
154
155 // run Funktion which does the beamtracking
156 // ALFA_BeamTransport::DoBeamTracking(evt_number);
157
158
159 TransportSelectedParticle(evt, evt_number);
160
161 // Print new data collection on screen
162 HepMC::Print::line(std::cout, evt);
163 }
164
165 ATH_MSG_INFO("after running the process function");
166
167 //-----------------------------------------------------------------------------------------
168
169 return StatusCode::SUCCESS;
170}
171
172
173
175{
176 //do particle tracking for beam 1
177
178 //tracking funktion
179 m_BeamTracker.ALFA_BeamTrack::CalculatePosRP(m_Particle1);//calculates position and momentum at RP1
180 //Position at RP
181 m_PosRP1=m_BeamTracker.ALFA_BeamTrack::PosRP();
182 //Momentum at RP
183 m_MomRP1=m_BeamTracker.ALFA_BeamTrack::MomRP();
184
185 //do particle tracking for beam 2
186
187 //tracking funktion
188 m_BeamTracker.ALFA_BeamTrack::CalculatePosRP(m_Particle2);//gives position and momentum at RP3
189 //Position at RP
190 m_PosRP3=m_BeamTracker.ALFA_BeamTrack::PosRP();
191 //Momentum at RP
192 m_MomRP3=m_BeamTracker.ALFA_BeamTrack::MomRP();
193
194 //Write Output
196 {
197 m_FileBeam1 << evt_number << "\t" << std::setprecision(17)
198 << m_PosRP1.x() << std::setprecision(17) << "\t"
199 << m_PosRP1.y() << "\t" << std::setprecision(17)
200 << m_PosRP1.z() << "\t" << std::setprecision(17)
201 << m_MomRP1.x() << std::setprecision(17) << "\t"
202 << m_MomRP1.y() << "\t" << std::setprecision(17)
203 << m_MomRP1.z() << "\t" << std::setprecision(17)
204 << m_EnergyRP1 << std::endl;
205 m_FileBeam2 << evt_number << "\t" << std::setprecision(17)
206 << m_PosRP3.x() << std::setprecision(17) << "\t"
207 << m_PosRP3.y() << "\t" << std::setprecision(17)
208 << m_PosRP3.z() << "\t" << std::setprecision(17)
209 << m_MomRP3.x() << std::setprecision(17) << "\t"
210 << m_MomRP3.y() << "\t" << std::setprecision(17)
211 << m_MomRP3.z() << "\t" << std::setprecision(17)
212 << m_EnergyRP3 << std::endl;
213 }
214
215 return true;
216}
217
218int ALFA_BeamTransport::TransportSelectedParticle(HepMC::GenEvent& evt, int evt_number){
219 HepMC::GenParticlePtr p1{nullptr};
220 HepMC::GenParticlePtr p2{nullptr};
221
222 std::vector<FPTracker::Point> PosAtRP1;
223 std::vector<FPTracker::Point> PosAtRP3;
224 std::vector<FPTracker::Point> MomAtPR1;
225 std::vector<FPTracker::Point> MomAtPR3;
226 std::vector<double> EnergyRP1;
227 std::vector<double> EnergyRP3;
228
229
230 double mom=0.;
231 double eta=0.;
232 double theta=0.;
233
234 // First we have to select the final state particles from the MC
235 for (const HepMC::GenParticlePtr& p : evt) {
236
237 // Simple Eta Pt cut to remove particles from BeamTransportation which
238 // have no chance to reach RP plane
239 mom = std::sqrt(std::pow(p->momentum().px(), 2) +
240 std::pow(p->momentum().py(), 2) +
241 std::pow(p->momentum().pz(), 2));
242 theta = std::acos(std::abs(p->momentum().pz()) / mom);
243 eta = -std::log(std::tan(theta / 2));
244
245 if (MC::isStable(p) && (!p->end_vertex())) {
246
247 int pid = p->pdg_id();
248 if (eta > m_EtaCut && 1 - std::abs(mom / m_FPConfig.pbeam0) < m_XiCut) {
249
250 // save a copy of the particles which passed the cut
251 HepMC::FourVector Position = p->production_vertex()->position();
252 HepMC::FourVector Momentum = p->momentum();
253 HepMC::GenVertexPtr Vertex = HepMC::newGenVertexPtr(Position); // copy of the vertex
254 HepMC::GenParticlePtr Particle = HepMC::newGenParticlePtr(Momentum, pid, 1);
255
256 Vertex->add_particle_out(std::move(Particle));
257 evt.add_vertex(std::move(Vertex));
258
259 // select direction of particle
260 if (p->momentum().pz() > 0. && pid == 2212) { // Beam1 TODO Tracking only works for protons!!!!
261 p1 = p;
262 // now we want to track the final particle if it's a protons
263 // Positions are given in mm FPTracker needs them in meter
265 p1->production_vertex()->position().x() / 1000.,
266 p1->production_vertex()->position().y() / 1000.,
267 p1->production_vertex()->position().z() / 1000.,
268 p1->momentum().px(),
269 p1->momentum().py(),
270 p1->momentum().pz());
271
272 // do particle tracking for beam 1
273
274 // tracking funktion
275 m_BeamTracker.ALFA_BeamTrack::CalculatePosRP(m_Particle1); // calculates position and momentum at RP1
276 // Position at RP
277 m_PosRP1 = m_BeamTracker.ALFA_BeamTrack::PosRP();
278 PosAtRP1.push_back(m_PosRP1);
279 // Momentum at RP
280 m_MomRP1 = m_BeamTracker.ALFA_BeamTrack::MomRP();
281 MomAtPR1.push_back(m_MomRP1);
282 // no Energy change between IP and RP
283 m_EnergyRP1 = p1->momentum().e();
284 EnergyRP1.push_back(m_EnergyRP1);
285 if (m_WriteDebugOutput) {
286 m_FileBeam1 << evt_number << "\t" << std::setprecision(17)
287 << m_PosRP1.x() << std::setprecision(17) << "\t"
288 << m_PosRP1.y() << "\t" << std::setprecision(17)
289 << m_PosRP1.z() << "\t" << std::setprecision(17)
290 << m_MomRP1.x() << std::setprecision(17) << "\t"
291 << m_MomRP1.y() << "\t" << std::setprecision(17)
292 << m_MomRP1.z() << "\t" << std::setprecision(17)
293 << m_EnergyRP1 << std::endl;
294 }
295
296 } else if (p->momentum().pz() < 0. && pid == 2212) { // beam 2
297 p2 = p;
298
300 p2->production_vertex()->position().x() / 1000.,
301 p2->production_vertex()->position().y() / 1000.,
302 p2->production_vertex()->position().z() / 1000.,
303 p2->momentum().px(),
304 p2->momentum().py(),
305 p2->momentum().pz());
306 // tracking funktion
307 m_BeamTracker.ALFA_BeamTrack::CalculatePosRP(m_Particle2); // gives position and momentum at RP3
308 // Position at RP
309 m_PosRP3 = m_BeamTracker.ALFA_BeamTrack::PosRP();
310 PosAtRP3.push_back(m_PosRP3);
311 // Momentum at RP
312 m_MomRP3 = m_BeamTracker.ALFA_BeamTrack::MomRP();
313 MomAtPR3.push_back(m_MomRP3);
314
315 m_EnergyRP3 = p2->momentum().e();
316 EnergyRP3.push_back(m_EnergyRP3);
317 // Write Output
318 if (m_WriteDebugOutput) {
319 m_FileBeam2 << evt_number << "\t" << std::setprecision(17)
320 << m_PosRP3.x() << std::setprecision(17) << "\t"
321 << m_PosRP3.y() << "\t" << std::setprecision(17)
322 << m_PosRP3.z() << "\t" << std::setprecision(17)
323 << m_MomRP3.x() << std::setprecision(17) << "\t"
324 << m_MomRP3.y() << "\t" << std::setprecision(17)
325 << m_MomRP3.z() << "\t" << std::setprecision(17)
326 << m_EnergyRP3 << std::endl;
327 }
328
329 } else {
330 ATH_MSG_ERROR("Strange: Particle rests at IP");
331 }
332
333 if (pid == 2212) { // Find the protons
334 m_pcount++;
335 if (m_pcount > 2) {
336 ATH_MSG_ERROR("Strange: More than two protons in this event!");
337 }
338 }
339 }
340 }
341 }
342
343 //from here the transported particles are added to HepMC
344
345 ATH_MSG_INFO("Add transproted particle into HepMC event record Beam 1");
346 //Add Data for HepMC Collection
347
348 for(int i=0;i<(int)PosAtRP1.size();i++){//Beam1
349 //The factor 1000 comes from the fact that HepMC saves length in mm
350 HepMC::FourVector PositionVectorRP1 = HepMC::FourVector(PosAtRP1.at(i).x()*1000.,PosAtRP1.at(i).y()*1000.,PosAtRP1.at(i).z()*1000.,0.*1000.);
351 HepMC::FourVector MomentumVectorRP1 = HepMC::FourVector(MomAtPR1.at(i).x(),MomAtPR1.at(i).y(),MomAtPR1.at(i).z(),EnergyRP1.at(i));
352
353 HepMC::GenVertexPtr VertexRP1 = HepMC::newGenVertexPtr(PositionVectorRP1);
354 HepMC::GenParticlePtr ParticleRP1 = HepMC::newGenParticlePtr(MomentumVectorRP1,2212,1); //save the transported particle with status code 1 (added 120124) preview was 201
355 //Add particle to vertex
356 VertexRP1->add_particle_out(std::move(ParticleRP1));
357 //add new vertex to HepMC event record
358 if(m_PosRP1.x()!=-99){ // add vertex to event record if the particle in the first beam was not lost
359 evt.add_vertex(std::move(VertexRP1));
360 }
361 }
362 ATH_MSG_INFO ("Add transproted particle into HepMC event record Beam 2" );
363 for(int i=0;i<(int)PosAtRP3.size();i++){
364 HepMC::FourVector PositionVectorRP3 = HepMC::FourVector(PosAtRP3.at(i).x()*1000.,PosAtRP3.at(i).y()*1000.,PosAtRP3.at(i).z()*1000.,0.*1000.);
365 HepMC::FourVector MomentumVectorRP3 = HepMC::FourVector(MomAtPR3.at(i).x(),MomAtPR3.at(i).y(),MomAtPR3.at(i).z(),EnergyRP3.at(i));
366 HepMC::GenVertexPtr VertexRP3 = HepMC::newGenVertexPtr(PositionVectorRP3);
367 HepMC::GenParticlePtr ParticleRP3 = HepMC::newGenParticlePtr(MomentumVectorRP3,2212,1);//save the transported particle with status code 1 (added 120124) preview was 201
368 VertexRP3->add_particle_out(std::move(ParticleRP3));
369 if (m_PosRP3.x() != -99) { // add vertex to event record if the // particle in the second beam was not lost
370 evt.add_vertex(std::move(VertexRP3));
371 }
372 }
373
374 MC::GeVToMeV(&evt); //Only scales momenta and masses
375 return true;
376}
Scalar eta() const
pseudorapidity method
Scalar theta() const
theta method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
boost::graph_traits< boost::adjacency_list< boost::vecS, boost::vecS, boost::bidirectionalS > >::vertex_descriptor Vertex
ATLAS-specific HepMC functions.
int TransportSelectedParticle(HepMC::GenEvent &evt, int evt_number)
This function is event selection, tracking and HepMC save ine one function.
FPTracker::Point m_PosRP3
ALFA_BeamTrack m_BeamTracker
FPTracker::Particle m_Particle1
StatusCode finalize()
standard Athena-Algorithm method
StatusCode initialize()
standard Athena-Algorithm method
StatusCode execute()
standard Athena-Algorithm method
FPTracker::Point m_PosRP1
ALFA_BeamTransport(const std::string &name, ISvcLocator *pSvcLocator)
Standard Athena-Algorithm Constructor.
SG::ReadHandleKey< McEventCollection > m_MCKey
SG::ReadHandleKey< xAOD::EventInfo > m_eventInfoKey
int DoBeamTracking(int evt_number)
Function which calls BeamTrack class to calcualte Position at RPs.
FPTracker::Particle m_Particle2
FPTracker::Point m_MomRP3
~ALFA_BeamTransport()
Default Destructor.
FPTracker::Point m_MomRP1
AthAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor with parameters:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
virtual bool isValid() override final
Can the handle be successfully dereferenced?
void line(std::ostream &os, const GenEvent &e)
Definition GenEvent.h:676
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
void MeVToGeV(HepMC::GenEvent *evt)
bool isStable(const T &p)
Identify if the particle is stable, i.e. has not decayed.
void GeVToMeV(HepMC::GenEvent *evt)