ATLAS Offline Software
Loading...
Searching...
No Matches
AFPHitAnalysis.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include "AFPHitAnalysis.h"
6
9
10
11
12
14 ATH_MSG_DEBUG( "Initializing AFPHitAnalysis" );
15
16 ATH_CHECK(m_readKey.initialize());
18 m_h_hitID = new TH1D("h_hitID", "hitID",100, 0., 100.);
19 m_h_hitID->StatOverflows();
20 ATH_CHECK(histSvc()->regHist(m_path+m_h_hitID->GetName(), m_h_hitID));
21
22 m_h_pdgID = new TH1D("h_pdgID", "pdgID", 200, -100,100);
23 m_h_pdgID->StatOverflows();
24 ATH_CHECK(histSvc()->regHist(m_path+m_h_pdgID->GetName(), m_h_pdgID));
25
26 m_h_trackID = new TH1D("h_trackID", "trackID", 100, 0,100);
27 m_h_trackID->StatOverflows();
28 ATH_CHECK(histSvc()->regHist(m_path+m_h_trackID->GetName(), m_h_trackID));
29
30 m_h_kine = new TH1D("h_kine", "kine", 100, 0,1000);
31 m_h_kine->StatOverflows();
32 ATH_CHECK(histSvc()->regHist(m_path+m_h_kine->GetName(), m_h_kine));
33
34 m_h_edep = new TH1D("h_edep", "edep", 100, 0,1000);
35 m_h_edep->StatOverflows();
36 ATH_CHECK(histSvc()->regHist(m_path+m_h_edep->GetName(), m_h_edep));
37
38 m_h_stepX = new TH1D("h_stepX", "stepX", 100, 0,1000);
39 m_h_stepX->StatOverflows();
40 ATH_CHECK(histSvc()->regHist(m_path+m_h_stepX->GetName(), m_h_stepX));
41
42 m_h_stepY = new TH1D("h_stepY", "stepY", 100, 0,1000);
43 m_h_stepY->StatOverflows();
44 ATH_CHECK(histSvc()->regHist(m_path+m_h_stepY->GetName(), m_h_stepY));
45
46 m_h_stepZ = new TH1D("h_stepZ", "stepZ", 100, 0,1000);
47 m_h_stepZ->StatOverflows();
48 ATH_CHECK(histSvc()->regHist(m_path+m_h_stepZ->GetName(), m_h_stepZ));
49
50 m_h_stationID = new TH1D("h_stationID", "stationID", 50, 0,50);
51 m_h_stationID->StatOverflows();
52 ATH_CHECK(histSvc()->regHist(m_path+m_h_stationID->GetName(), m_h_stationID));
53
54 m_h_detID = new TH1D("h_detID", "detID", 50, 0,50);
55 m_h_detID->StatOverflows();
56 ATH_CHECK(histSvc()->regHist(m_path+m_h_detID->GetName(), m_h_detID));
57
58 m_h_pixelRow = new TH1D("h_pixelRow", "pixelRow", 20, 0,20);
59 m_h_pixelRow->StatOverflows();
60 ATH_CHECK(histSvc()->regHist(m_path+m_h_pixelRow->GetName(), m_h_pixelRow));
61
62 m_h_pixelCol = new TH1D("h_pixelCol", "pixelCol", 20, 0,20);
63 m_h_pixelCol->StatOverflows();
64 ATH_CHECK(histSvc()->regHist(m_path+m_h_pixelCol->GetName(), m_h_pixelCol));
65
67 m_tree = new TTree("AFP","AFP");
68 std::string fullNtupleName = "/" + m_ntupleFileName + "/";
69 ATH_CHECK(histSvc()->regTree(fullNtupleName,m_tree));
70
71 if (m_tree){
72 m_tree->Branch("hitID", &m_hitID);
73 m_tree->Branch("pdgID", &m_pdgID);
74 m_tree->Branch("trackID", &m_trackID);
75 m_tree->Branch("kine", &m_kine);
76 m_tree->Branch("edep", &m_edep);
77 m_tree->Branch("stepX", &m_stepX);
78 m_tree->Branch("stepY", &m_stepY);
79 m_tree->Branch("stepZ", &m_stepZ);
80 m_tree->Branch("stationID", &m_stationID);
81 m_tree->Branch("detID", &m_detID);
82 m_tree->Branch("pixelRow", &m_pixelRow);
83 m_tree->Branch("pixelCol", &m_pixelCol);
84 }
85
86 ATH_MSG_INFO("Exiting AFPHitAnalysis::initialize()");
87
88 return StatusCode::SUCCESS;
89}
90
91
93 ATH_MSG_DEBUG( "In AFPHitAnalysis::execute()" );
94
95 m_hitID->clear();
96 m_pdgID->clear();
97 m_trackID->clear();
98 m_kine->clear();
99 m_edep->clear();
100 m_stepX->clear();
101 m_stepY->clear();
102 m_stepZ->clear();
103 m_time->clear();
104 m_stationID->clear();
105 m_detID->clear();
106 m_pixelRow->clear();
107 m_pixelCol->clear();
108
109 ATH_MSG_INFO( "AFPHitAnalysis tree branches cleared" );
110
112
113 const EventContext& ctx{Gaudi::Hive::currentContext()};
114 const AFP_SIDSimHitCollection* iter{nullptr};
115 ATH_CHECK(SG::get(iter, m_readKey, ctx));
116
117 ATH_MSG_DEBUG( "AFP_SIDSSimHitCollection retrieved" );
118
119 for ( hi=(*iter).begin(); hi != (*iter).end(); ++hi ) {
120 AFP_SIDSimHit ghit(*hi);
121 m_h_hitID->Fill(ghit.m_nHitID);
122 m_h_pdgID->Fill(ghit.m_nParticleEncoding);
123 m_h_trackID->Fill(ghit.m_nTrackID);
124 m_h_kine->Fill(ghit.m_fKineticEnergy);
125 m_h_edep->Fill(ghit.m_fEnergyDeposit);
126 m_h_stepX->Fill(ghit.m_fPostStepX-ghit.m_fPreStepX);
127 m_h_stepY->Fill(ghit.m_fPostStepY-ghit.m_fPreStepY);
128 m_h_stepX->Fill(ghit.m_fPostStepZ-ghit.m_fPreStepZ);
129 m_h_time->Fill(ghit.m_fGlobalTime);
130 m_h_stationID->Fill(ghit.m_nStationID);
131 m_h_detID->Fill(ghit.m_nDetectorID);
132 m_h_pixelRow->Fill(ghit.m_nPixelRow);
133 m_h_pixelCol->Fill(ghit.m_nPixelCol);
134
135 m_hitID->push_back(ghit.m_nHitID);
136 m_pdgID->push_back(ghit.m_nParticleEncoding);
137 m_trackID->push_back(ghit.m_nTrackID);
138 m_kine->push_back(ghit.m_fKineticEnergy);
139 m_edep->push_back(ghit.m_fEnergyDeposit);
140 m_stepX->push_back(ghit.m_fPostStepX-ghit.m_fPreStepX);
141 m_stepY->push_back(ghit.m_fPostStepY-ghit.m_fPreStepY);
142 m_stepX->push_back(ghit.m_fPostStepZ-ghit.m_fPreStepZ);
143 m_time->push_back(ghit.m_fGlobalTime);
144 m_stationID->push_back(ghit.m_nStationID);
145 m_detID->push_back(ghit.m_nDetectorID);
146 m_pixelRow->push_back(ghit.m_nPixelRow);
147 m_pixelCol->push_back(ghit.m_nPixelCol);
148 }
149
150 if (m_tree) m_tree->Fill();
151
152 return StatusCode::SUCCESS;
153}
AtlasHitsVector< AFP_SIDSimHit >::const_iterator AFP_SIDSimHitConstIter
AtlasHitsVector< AFP_SIDSimHit > AFP_SIDSimHitCollection
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_INFO(x)
#define ATH_MSG_DEBUG(x)
Handle class for reading from StoreGate.
std::vector< float > * m_stepZ
std::vector< float > * m_edep
std::vector< float > * m_hitID
virtual StatusCode initialize() override final
std::vector< int > * m_pixelCol
Gaudi::Property< std::string > m_ntupleFileName
std::vector< float > * m_pdgID
SG::ReadHandleKey< AFP_SIDSimHitCollection > m_readKey
std::vector< int > * m_detID
std::vector< float > * m_stepY
virtual StatusCode execute() override final
std::vector< int > * m_stationID
std::vector< float > * m_stepX
Gaudi::Property< std::string > m_path
std::vector< float > * m_trackID
TH1 * m_h_hitID
Some histograms.
std::vector< float > * m_kine
std::vector< float > * m_time
std::vector< int > * m_pixelRow
float m_fKineticEnergy
float m_fEnergyDeposit
const ServiceHandle< ITHistSvc > & histSvc() const
The standard THistSvc (for writing histograms and TTrees and more to a root file) Returns (kind of) a...
const T * get(const ReadCondHandleKey< T > &key, const EventContext &ctx)
Convenience function to retrieve an object given a ReadCondHandleKey.