ATLAS Offline Software
Loading...
Searching...
No Matches
LucidHitAnalysis.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 "LucidHitAnalysis.h"
7
8
9
11 ATH_MSG_DEBUG( "Initializing LucidHitAnalysis" );
12
13 // Grab the Ntuple and histogramming service for the tree
14 ATH_CHECK(m_readKey.initialize());
16 m_h_hit_x = new TH1D("h_hit_x", "hit_x", 100,-150.,150.);
17 m_h_hit_x->StatOverflows();
18 ATH_CHECK(histSvc()->regHist( m_path+m_h_hit_x->GetName(), m_h_hit_x));
19
20 m_h_hit_y = new TH1D("h_hit_y", "hit_y", 100,-150.,150.);
21 m_h_hit_y->StatOverflows();
22 ATH_CHECK(histSvc()->regHist( m_path+m_h_hit_y->GetName(), m_h_hit_y));
23
24 m_h_hit_z = new TH1D("h_hit_z", "hit_z", 100,-20000.,20000.);
25 m_h_hit_z->StatOverflows();
26 ATH_CHECK(histSvc()->regHist( m_path+m_h_hit_z->GetName(), m_h_hit_z));
27
28 m_h_xy = new TH2D("h_xy", "hit_xy", 100,-150.,150.,100,-150,150);
29 m_h_xy->StatOverflows();
30 ATH_CHECK(histSvc()->regHist( m_path+m_h_xy->GetName(), m_h_xy));
31
32 m_h_zr = new TH2D("h_zr", "hit_zr", 100,-20000.,20000.,100,0,250);
33 m_h_zr->StatOverflows();
34 ATH_CHECK(histSvc()->regHist( m_path+m_h_zr->GetName(), m_h_zr));
35
36 m_h_hit_post_x = new TH1D("h_hit_post_x", "hit_post_x", 100,-150.,150.);
37 m_h_hit_post_x->StatOverflows();
38 ATH_CHECK(histSvc()->regHist( m_path+m_h_hit_post_x->GetName(), m_h_hit_post_x));
39
40 m_h_hit_post_y = new TH1D("h_hit_post_y", "hit_post_y", 100,-150,150.);
41 m_h_hit_post_y->StatOverflows();
42 ATH_CHECK(histSvc()->regHist( m_path+m_h_hit_post_y->GetName(), m_h_hit_post_y));
43
44 m_h_hit_post_z = new TH1D("h_hit_post_z", "hit_post_z", 100,-15000,15000.);
45 m_h_hit_post_z->StatOverflows();
46 ATH_CHECK(histSvc()->regHist( m_path+m_h_hit_post_z->GetName(), m_h_hit_post_z));
47
48 m_h_hit_edep = new TH1D("h_hit_edep", "hit_edep", 100,0.,20.);
49 m_h_hit_edep->StatOverflows();
50 ATH_CHECK(histSvc()->regHist( m_path+m_h_hit_edep->GetName(), m_h_hit_edep));
51
52 m_h_hit_pdgid = new TH1D("h_hit_pdgid", "hit_pdgid", 100,0.,7e6);
53 m_h_hit_pdgid->StatOverflows();
54 ATH_CHECK(histSvc()->regHist( m_path+m_h_hit_pdgid->GetName(), m_h_hit_pdgid));
55
56 m_h_hit_pretime = new TH1D("h_hit_pretime", "hit_pretime", 100,0.,100.);
57 m_h_hit_pretime->StatOverflows();
58 ATH_CHECK(histSvc()->regHist( m_path+m_h_hit_pretime->GetName(), m_h_hit_pretime));
59
60 m_h_hit_posttime = new TH1D("h_hit_posttime", "hit_posttime", 100,0.,100.);
61 m_h_hit_posttime->StatOverflows();
63
64 m_h_genvolume = new TH1D("h_genvolume", "genvolume", 20,0.,5.);
65 m_h_genvolume->StatOverflows();
66 ATH_CHECK(histSvc()->regHist( m_path+m_h_genvolume->GetName(), m_h_genvolume));
67
68 m_h_wavelength = new TH1D("m_wavelength", "wavelength", 150,0.,800.);
69 m_h_wavelength->StatOverflows();
70 ATH_CHECK(histSvc()->regHist( m_path+m_h_wavelength->GetName(), m_h_wavelength));
71
73 m_tree = new TTree("Lucid","Lucid");
74 std::string fullNtupleName = "/" + m_ntupleFileName + "/";
75 ATH_CHECK(histSvc()->regTree(fullNtupleName,m_tree));
76
77 m_tree->Branch("hit_x", &m_hit_x);
78 m_tree->Branch("hit_y", &m_hit_y);
79 m_tree->Branch("hit_z", &m_hit_z);
80 m_tree->Branch("hit_post_x", &m_hit_post_x);
81 m_tree->Branch("hit_post_y", &m_hit_post_y);
82 m_tree->Branch("hit_post_z", &m_hit_post_z);
83 m_tree->Branch("edep", &m_hit_edep);
84 m_tree->Branch("pdgid", &m_hit_pdgid);
85 m_tree->Branch("pretime", &m_hit_pretime);
86 m_tree->Branch("posttime", &m_hit_posttime);
87 m_tree->Branch("gen_volume", &m_gen_volume);
88 m_tree->Branch("wavelength", &m_wavelength);
89
90 return StatusCode::SUCCESS;
91}
92
93
95 ATH_MSG_DEBUG( "In LucidHitAnalysis::execute()" );
96
97 m_hit_x->clear();
98 m_hit_y->clear();
99 m_hit_z->clear();
100 m_hit_post_x->clear();
101 m_hit_post_y->clear();
102 m_hit_post_z->clear();
103 m_hit_edep->clear();
104 m_hit_pdgid->clear();
105 m_hit_pretime->clear();
106 m_hit_posttime->clear();
107 m_gen_volume->clear();
108 m_wavelength->clear();
109
110 const EventContext& ctx{Gaudi::Hive::currentContext()};
111 const LUCID_SimHitCollection* iter{nullptr};
112 ATH_CHECK(SG::get(iter, m_readKey, ctx));
113 for (LUCID_SimHitCollection::const_iterator i_hit = (*iter).begin(); i_hit != (*iter).end(); ++i_hit) {
114 double x = i_hit->GetX();
115 double y = i_hit->GetY();
116 double z = i_hit->GetZ();
117 double r = sqrt(x*x+y*y);
118
119 m_h_hit_x->Fill(i_hit->GetX());
120 m_h_hit_y->Fill(i_hit->GetY());
121 m_h_hit_z->Fill(i_hit->GetZ());
122 m_h_xy->Fill(x, y);
123 m_h_zr->Fill(z, r);
124 m_h_hit_post_x->Fill(i_hit->GetEPX());
125 m_h_hit_post_y->Fill(i_hit->GetEPY());
126 m_h_hit_post_z->Fill(i_hit->GetEPZ());
127 m_h_hit_edep->Fill(i_hit->GetEnergy());
128 m_h_hit_pdgid->Fill(i_hit->GetPdgCode());
129 m_h_hit_pretime->Fill(i_hit->GetPreStepTime());
130 m_h_hit_posttime->Fill(i_hit->GetPostStepTime());
131 m_h_genvolume->Fill(i_hit->GetGenVolume());
132 m_h_wavelength->Fill(i_hit->GetWavelength());
133
134 m_hit_x->push_back(i_hit->GetX());
135 m_hit_y->push_back(i_hit->GetY());
136 m_hit_z->push_back(i_hit->GetZ());
137 m_hit_post_x->push_back(i_hit->GetEPX());
138 m_hit_post_y->push_back(i_hit->GetEPY());
139 m_hit_post_z->push_back(i_hit->GetEPZ());
140 m_hit_edep->push_back(i_hit->GetEnergy());
141 m_hit_pdgid->push_back(i_hit->GetPdgCode());
142 m_hit_pretime->push_back(i_hit->GetPreStepTime());
143 m_hit_posttime->push_back(i_hit->GetPostStepTime());
144 m_gen_volume->push_back(i_hit->GetGenVolume());
145 m_wavelength->push_back(i_hit->GetWavelength());
146 }
147
148 if (m_tree) m_tree->Fill();
149
150 return StatusCode::SUCCESS;
151
152}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_DEBUG(x)
AtlasHitsVector< LUCID_SimHit > LUCID_SimHitCollection
Property holding a SG store/key/clid from which a ReadHandle is made.
#define y
#define x
#define z
const ServiceHandle< ITHistSvc > & histSvc() const
The standard THistSvc (for writing histograms and TTrees and more to a root file) Returns (kind of) a...
CONT::const_iterator const_iterator
StringProperty m_ntupleFileName
std::vector< float > * m_wavelength
std::vector< float > * m_hit_posttime
std::vector< float > * m_hit_z
SG::ReadHandleKey< LUCID_SimHitCollection > m_readKey
TH1 * m_h_hit_x
Some histograms.
StringProperty m_path
std::vector< float > * m_hit_x
virtual StatusCode initialize() override final
std::vector< float > * m_hit_pdgid
std::vector< float > * m_hit_pretime
std::vector< float > * m_gen_volume
std::vector< float > * m_hit_post_y
std::vector< float > * m_hit_post_x
std::vector< float > * m_hit_y
virtual StatusCode execute() override final
std::vector< float > * m_hit_post_z
std::vector< float > * m_hit_edep
int r
Definition globals.cxx:22
const T * get(const ReadCondHandleKey< T > &key, const EventContext &ctx)
Convenience function to retrieve an object given a ReadCondHandleKey.