ATLAS Offline Software
Loading...
Searching...
No Matches
TruthElectronHistograms.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
6
8#include "GaudiKernel/ITHistSvc.h"
9
11
12#include "TH1D.h"
13#include "TH2D.h"
14
15using namespace egammaMonitoring;
16
18 return initializePlots (false);
19}
20
21StatusCode TruthElectronHistograms::initializePlots(bool reducedHistSet) {
22
23 const char* fN = m_name.c_str();
24
25 if (!reducedHistSet) {
26 histoMap["deltaPhi2"] = new TH1D(Form("%s_deltaPhi2",fN), ";deltaPhi2; Events", 40, -0.06, 0.06);
27 histoMap["deltaEta2"] = new TH1D(Form("%s_deltaEta2",fN), ";deltaEta2; Events", 40, -0.04, 0.04);
28 histoMap["deltaPhiRescaled2"] = new TH1D(Form("%s_deltaPhiRescaled2",fN), ";deltaPhiRescaled2; Events", 40, -0.04, 0.04);
29
30 histoMap["d0Oversigmad0"] = new TH1D(Form("%s_d0Oversigmad0",fN), "; d0Oversigmad0; Events", 40, -10, 10);
31 histoMap["qOverp_resolution"] = new TH1D(Form("%s_qOverp_resolution",fN), ";(q/P reco - q/P truth)/ q/p truth; Events", 60, -1, 1.5);
32
33 ATH_CHECK(m_rootHistSvc->regHist(m_folder+"deltaPhi2", histoMap["deltaPhi2"]));
34 ATH_CHECK(m_rootHistSvc->regHist(m_folder+"deltaEta2", histoMap["deltaEta2"]));
35 ATH_CHECK(m_rootHistSvc->regHist(m_folder+"deltaPhiRescaled2", histoMap["deltaPhiRescaled2"]));
36 ATH_CHECK(m_rootHistSvc->regHist(m_folder+"d0Oversigmad0", histoMap["d0Oversigmad0"]));
37 ATH_CHECK(m_rootHistSvc->regHist(m_folder+"qOverp_resolution", histoMap["qOverp_resolution"]));
38
39 // 2D only for truthPromptElectronWithRecoTrack (temporary)
40 if (m_name == "truthPromptElectronWithRecoTrack") {
41 histoMap2D["eta_deltaPhi2"] = new TH2D(Form("%s_eta_deltaPhi2",fN),
42 ";#eta;#Delta#phi_{2}; Events", 90, -4.5, 4.5, 40, -0.06, 0.06);
43 histoMap2D["eta_deltaEta2"] = new TH2D(Form("%s_eta_deltaEta2",fN),
44 ";#eta;#Delta#eta_{2}; Events", 90, -4.5, 4.5, 40, -0.06, 0.06);
45 histoMap2D["eta_deltaPhiRescaled2"] = new TH2D(Form("%s_eta_deltaPhiRescaled2",fN),
46 ";#eta;#Delta#phi_{2}^{Rescaled}; Events", 90, -4.5, 4.5, 40, -0.06, 0.06);
47 histoMap2D["eta_d0Oversigmad0"] = new TH2D(Form("%s_eta_d0Oversigmad0",fN),
48 ";#eta;d_{0}/#sigma_{d_{0}}; Events", 90, -4.5, 4.5, 40, -10, 10);
49 histoMap2D["eta_qOverp_resolution"] = new TH2D(Form("%s_eta_qOverp_resolution",fN),
50 ";#eta;(q/P_{reco})/(q/P_{truth}) -1; Events", 90, -4.5, 4.5, 60, -1, 1.5);
51
52 for (const auto& e : histoMap2D) {
53 ATH_CHECK(m_rootHistSvc->regHist(m_folder+e.first, e.second));
54 }
55 }
56 }
57
59
60 m_reducedHistSet = reducedHistSet;
61
62 return StatusCode::SUCCESS;
63
64}
65
66
68
70
71 if (!electron || m_reducedHistSet) return;
72
73 const xAOD::TrackParticle* track = electron->trackParticle();
74
75 // This can happen if we use it for forwardElectron
76 if (!track) return;
77
78 bool has2DHis = !histoMap2D.empty();
79
80 float dphires2(0.);
81 float dphi2(0.);
82 float deta2(0);
83
84 if (electron->trackCaloMatchValue(dphires2, xAOD::EgammaParameters::deltaPhiRescaled2)) {
85 histoMap["deltaPhiRescaled2"]->Fill(dphires2);
86 if (has2DHis) histoMap2D["eta_deltaPhiRescaled2"]->Fill(electron->eta(),dphires2);
87 }
88 if (electron->trackCaloMatchValue(dphi2, xAOD::EgammaParameters::deltaPhi2)) {
89 histoMap["deltaPhi2"]->Fill(dphi2);
90 if (has2DHis) histoMap2D["eta_deltaPhi2"]->Fill(electron->eta(),dphi2);
91 }
92 if (electron->trackCaloMatchValue(deta2, xAOD::EgammaParameters::deltaEta2)) {
93 histoMap["deltaEta2"]->Fill(deta2);
94 if (has2DHis) histoMap2D["eta_deltaEta2"]->Fill(electron->eta(),deta2);
95 }
96
97 float d0 = track->d0();
98 float reco_qp = track->qOverP();
99 float truth_qp = truth->charge()/(truth->pt()*cosh(truth->eta()));
100 float vard0 = track->definingParametersCovMatrix()(0, 0);
101
102 if (vard0 > 0) {
103 histoMap["d0Oversigmad0"]->Fill(d0/sqrtf(vard0));
104 if (has2DHis) histoMap2D["eta_d0Oversigmad0"]->Fill(electron->eta(),d0/sqrtf(vard0));
105 }
106
107 if (truth_qp > 0) {
108 histoMap["qOverp_resolution"]->Fill((reco_qp-truth_qp)/truth_qp);
109 if (has2DHis) histoMap2D["eta_qOverp_resolution"]->Fill(electron->eta(),(reco_qp-truth_qp)/truth_qp);
110 }
111
112
113}
114
#define ATH_CHECK
Evaluate an expression and check for errors.
SmartIF< ITHistSvc > m_rootHistSvc
Definition IHistograms.h:48
std::map< std::string, TH1D * > histoMap
Definition IHistograms.h:42
void fill(const xAOD::IParticle &egamma)
void fill(const xAOD::TruthParticle *truth, const xAOD::Electron *el=nullptr)
virtual double pt() const override final
The transverse momentum ( ) of the particle.
double charge() const
Physical charge.
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
@ deltaEta2
difference between the cluster eta (second sampling) and the eta of the track extrapolated to the sec...
@ deltaPhiRescaled2
difference between the cluster phi (second sampling) and the phi of the track extrapolated to the sec...
@ deltaPhi2
difference between the cluster phi (second sampling) and the phi of the track extrapolated to the sec...
TrackParticle_v1 TrackParticle
Reference the current persistent version:
TruthParticle_v1 TruthParticle
Typedef to implementation.
Electron_v1 Electron
Definition of the current "egamma version".