ATLAS Offline Software
Loading...
Searching...
No Matches
TGCHitAnalysis.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 "TGCHitAnalysis.h"
6
7// Section of includes for TGC of the Muon Spectrometer tests
9
11#include "CLHEP/Vector/LorentzVector.h"
12
13#include "TH1.h"
14#include "TTree.h"
15#include "TString.h"
16
17#include <algorithm>
18#include <math.h>
19#include <functional>
20#include <iostream>
21#include <stdio.h>
22
24 ATH_MSG_DEBUG( "Initializing TGCHitAnalysis" );
25
26 // Grab the Ntuple and histogramming service for the tree
27 ATH_CHECK(m_readKey.initialize());
28
30 m_h_hits_x = new TH1D("h_hits_tgc_x","hits_x", 100,-5000, 5000);
31 m_h_hits_x->StatOverflows();
32 ATH_CHECK(histSvc()->regHist(m_path + m_h_hits_x->GetName(), m_h_hits_x));
33
34 m_h_hits_y = new TH1D("h_hits_tgc_y", "hits_y", 100,-5000,5000);
35 m_h_hits_y->StatOverflows();
36 ATH_CHECK(histSvc()->regHist(m_path + m_h_hits_y->GetName(), m_h_hits_y));
37
38 m_h_hits_z = new TH1D("h_hits_tgc_z", "hits_z", 100,-12000,12000);
39 m_h_hits_z->StatOverflows();
40 ATH_CHECK(histSvc()->regHist(m_path + m_h_hits_z->GetName(), m_h_hits_z));
41
42 m_h_hits_r = new TH1D("h_hits_tgc_r", "hits_r", 100,2000,10000);
43 m_h_hits_r->StatOverflows();
44 ATH_CHECK(histSvc()->regHist(m_path + m_h_hits_r->GetName(), m_h_hits_r));
45
46 m_h_xy = new TH2D("h_tgc_xy", "xy", 100,-5000.,5000.,100, -5000., 5000.);
47 m_h_xy->StatOverflows();
48 ATH_CHECK(histSvc()->regHist( m_path+m_h_xy->GetName(), m_h_xy));
49
50 m_h_rz = new TH2D("h_tgc_rz", "rz", 100,2000.,10000.,100, -12000., 12000.);
51 m_h_rz->StatOverflows();
52 ATH_CHECK(histSvc()->regHist( m_path+m_h_rz->GetName(), m_h_rz));
53
54 m_h_hits_eta = new TH1D("h_hits_tgc_eta", "hits_eta", 100,-10.0,10.0);
55 m_h_hits_eta->StatOverflows();
56 ATH_CHECK(histSvc()->regHist(m_path + m_h_hits_eta->GetName(), m_h_hits_eta));
57
58 m_h_hits_phi = new TH1D("h_hits_tgc_phi", "hits_phi", 100,-3.2,3.2);
59 m_h_hits_phi->StatOverflows();
60 ATH_CHECK(histSvc()->regHist(m_path + m_h_hits_phi->GetName(), m_h_hits_phi));
61
62 m_h_hits_lx = new TH1D("h_hits_tgc_lx","hits_lx", 100,-800, 800);
63 m_h_hits_lx->StatOverflows();
64 ATH_CHECK(histSvc()->regHist(m_path + m_h_hits_lx->GetName(), m_h_hits_lx));
65
66 m_h_hits_ly = new TH1D("h_hits_tgc_ly", "hits_ly", 100,-800,800);
67 m_h_hits_ly->StatOverflows();
68 ATH_CHECK(histSvc()->regHist(m_path + m_h_hits_ly->GetName(), m_h_hits_ly));
69
70 m_h_hits_lz = new TH1D("h_hits_tgc_lz", "hits_lz", 100,-800,800);
71 m_h_hits_lz->StatOverflows();
72 ATH_CHECK(histSvc()->regHist(m_path + m_h_hits_lz->GetName(), m_h_hits_lz));
73
74 m_h_hits_dcx = new TH1D("h_hits_tgc_dcx","hits_dcx", 100,-1, 1);
75 m_h_hits_dcx->StatOverflows();
76 ATH_CHECK(histSvc()->regHist(m_path + m_h_hits_dcx->GetName(), m_h_hits_dcx));
77
78 m_h_hits_dcy = new TH1D("h_hits_tgc_dcy", "hits_dcy", 100,-1,1);
79 m_h_hits_dcy->StatOverflows();
80 ATH_CHECK(histSvc()->regHist(m_path + m_h_hits_dcy->GetName(), m_h_hits_dcy));
81
82 m_h_hits_dcz = new TH1D("h_hits_tgc_dcz", "hits_dcz", 100,-1,1);
83 m_h_hits_dcz->StatOverflows();
84 ATH_CHECK(histSvc()->regHist(m_path + m_h_hits_dcz->GetName(), m_h_hits_dcz));
85
86 m_h_hits_time = new TH1D("h_hits_tgc_time","hits_time", 100,0, 250);
87 m_h_hits_time->StatOverflows();
88 ATH_CHECK(histSvc()->regHist(m_path + m_h_hits_time->GetName(), m_h_hits_time));
89
90 m_h_hits_edep = new TH1D("h_hits_tgc_edep", "hits_edep", 100,0,0.5);
91 m_h_hits_edep->StatOverflows();
92 ATH_CHECK(histSvc()->regHist(m_path + m_h_hits_edep->GetName(), m_h_hits_edep));
93
94 m_h_hits_kine = new TH1D("h_hits_tgc_kine", "hits_kine", 100,0,1000);
95 m_h_hits_kine->StatOverflows();
96 ATH_CHECK(histSvc()->regHist(m_path + m_h_hits_kine->GetName(), m_h_hits_kine));
97
98 m_h_hits_step = new TH1D("h_hits_tgc_step", "hits_step", 100,0,50);
99 m_h_hits_step->StatOverflows();
100 ATH_CHECK(histSvc()->regHist(m_path + m_h_hits_step->GetName(), m_h_hits_step));
101
102
103 return StatusCode::SUCCESS;
104}
105
106
108 ATH_MSG_DEBUG( "In TGCHitAnalysis::execute()" );
109
110 const EventContext& ctx{Gaudi::Hive::currentContext()};
111 const TGCSimHitCollection* tgc_container{nullptr};
112 ATH_CHECK(SG::get(tgc_container, m_readKey, ctx));
113 for (TGCSimHitCollection::const_iterator i_hit = tgc_container->begin(); i_hit != tgc_container->end(); ++i_hit) {
114 //TGCSimHitCollection::const_iterator i_hit;
115 //for(auto i_hit : *tgc_container){
116 GeoTGCHit ghit(*i_hit);
117 if (!ghit) continue;
118
120 m_h_hits_x->Fill(p.x());
121 m_h_hits_y->Fill(p.y());
122 m_h_hits_z->Fill(p.z());
123 m_h_hits_r->Fill(p.perp());
124 m_h_xy->Fill(p.x(), p.y());
125 m_h_rz->Fill(p.perp(), p.z());
126 m_h_hits_eta->Fill(p.eta());
127 m_h_hits_phi->Fill(p.phi());
128 m_h_hits_lx->Fill((*i_hit).localPosition().x());
129 m_h_hits_ly->Fill((*i_hit).localPosition().y());
130 m_h_hits_lz->Fill((*i_hit).localPosition().z());
131 m_h_hits_dcx->Fill((*i_hit).localDireCos().x());
132 m_h_hits_dcy->Fill((*i_hit).localDireCos().y());
133 m_h_hits_dcz->Fill((*i_hit).localDireCos().z());
134 m_h_hits_edep->Fill((*i_hit).energyDeposit());
135 m_h_hits_time->Fill((*i_hit).globalTime());
136 m_h_hits_step->Fill((*i_hit).stepLength());
137 m_h_hits_kine->Fill((*i_hit).kineticEnergy());
138
139 }
140
141 return StatusCode::SUCCESS;
142}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_DEBUG(x)
AtlasHitsVector< TGCSimHit > TGCSimHitCollection
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
const_iterator begin() const
const_iterator end() const
Amg::Vector3D getGlobalPosition() const
SG::ReadHandleKey< TGCSimHitCollection > m_readKey
TH1 * m_h_hits_x
Some variables.
virtual StatusCode initialize() override
Gaudi::Property< std::string > m_path
virtual StatusCode execute() override
Eigen::Matrix< double, 3, 1 > Vector3D
const T * get(const ReadCondHandleKey< T > &key, const EventContext &ctx)
Convenience function to retrieve an object given a ReadCondHandleKey.