ATLAS Offline Software
Loading...
Searching...
No Matches
egammaMonitoring::ClusterHistograms Class Reference

#include <ClusterHistograms.h>

Collaboration diagram for egammaMonitoring::ClusterHistograms:

Public Member Functions

 ClusterHistograms (std::string name, std::string title, std::string folder, SmartIF< ITHistSvc > rootHistSvc)
StatusCode initializePlots ()
void fill (const xAOD::Egamma &egamma)
void fill (const xAOD::Egamma &egamma, float mu)

Public Attributes

std::map< std::string, TH2D * > histo2DMap
std::map< std::string, TProfile * > profileMap

Protected Attributes

std::string m_name
std::string m_title
std::string m_folder
SmartIF< ITHistSvc > m_rootHistSvc

Detailed Description

Definition at line 21 of file ClusterHistograms.h.

Constructor & Destructor Documentation

◆ ClusterHistograms()

egammaMonitoring::ClusterHistograms::ClusterHistograms ( std::string name,
std::string title,
std::string folder,
SmartIF< ITHistSvc > rootHistSvc )
inline

Definition at line 25 of file ClusterHistograms.h.

29 :
30 m_name(std::move(name)),
31 m_title(std::move(title)),
32 m_folder(std::move(folder)),
33 m_rootHistSvc(std::move(rootHistSvc)) {}

Member Function Documentation

◆ fill() [1/2]

void ClusterHistograms::fill ( const xAOD::Egamma & egamma)

Definition at line 69 of file ClusterHistograms.cxx.

69 {
70 fill(egamma,0.);
71}
void fill(const xAOD::Egamma &egamma)

◆ fill() [2/2]

void ClusterHistograms::fill ( const xAOD::Egamma & egamma,
float mu = 0 )

Definition at line 73 of file ClusterHistograms.cxx.

73 {
74
75 const xAOD::CaloCluster *cluster = egamma.caloCluster();
76
78
79 if ( !truth_egamma ) return;
80
81 const auto Ereco = cluster->rawE();
82 const auto Etruth = truth_egamma->e();
83 const auto Eres = (Ereco - Etruth)/Etruth;
84
85 profileMap["Eraw_Etruth_vs_Etruth_profile"]->Fill(Etruth/1000,Ereco/Etruth);
86 profileMap["Eraw_Etruth_vs_eta_profile"]->Fill(truth_egamma->eta(),Ereco/Etruth);
87 histo2DMap["mu_energy_resolution_2D"]->Fill(mu,Eres);
88
89 const CaloClusterCellLink* cellLinks = cluster->getCellLinks();
90
91 if ( !cellLinks ) return;
92
93 std::map<int, int > cells_per_layer;
94 std::map<int, int > cells_per_sam;
95
96 for (const CaloCell* cell : *cellLinks) {
97 if (cell) {
98 int layer = cell->caloDDE()->getLayer();
99 cells_per_layer[layer]++;
100 int sam = cell->caloDDE()->getSampling();
101 cells_per_sam[sam]++;
102 }
103 }
104
105 auto associatedTopoCluster = xAOD::EgammaHelpers::getAssociatedTopoClusters(cluster);
106 profileMap["number_topocluster_vs_e_profile"]->Fill(Etruth/1000, associatedTopoCluster.size());
107 profileMap["number_topocluster_vs_eta_profile"]->Fill(truth_egamma->eta(), associatedTopoCluster.size());
108
109 for (auto const& x : cells_per_layer) {
110 profileMap["number_cell_in_layer"]->Fill(x.second,x.first);
111 profileMap["number_cells_vs_e_in_layer_"+std::to_string(x.first)+"_profile"]->Fill(Etruth/1000, x.second);
112 profileMap["number_cells_vs_eta_in_layer_"+std::to_string(x.first)+"_profile"]->Fill(truth_egamma->eta(), x.second);
113 }
114 for (auto const& x : cells_per_sam) {
115 profileMap[Form("number_cells_vs_eta_in_sampling_%i_profile",x.first)]->Fill(truth_egamma->eta(), x.second);
116 }
117
118}
#define x
std::map< std::string, TProfile * > profileMap
std::map< std::string, TH2D * > histo2DMap
flt_t rawE() const
const CaloClusterCellLink * getCellLinks() const
Get a pointer to the CaloClusterCellLink object (const version)
virtual double e() const override final
The total energy of the particle.
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
@ layer
Definition HitInfo.h:79
std::vector< const xAOD::CaloCluster * > getAssociatedTopoClusters(const xAOD::CaloCluster *cluster)
Return a vector of all the topo clusters associated with the egamma cluster.
const xAOD::TruthParticle * getTruthParticle(const xAOD::IParticle &p)
Return the truthParticle associated to the given IParticle (if any)
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
TruthParticle_v1 TruthParticle
Typedef to implementation.

◆ initializePlots()

StatusCode ClusterHistograms::initializePlots ( )

Definition at line 20 of file ClusterHistograms.cxx.

20 {
21
22 const char* fN = m_name.c_str();
23
24 for (int is = 0; is < 4; is++) {
25 TString pNm = Form("number_cells_vs_eta_in_layer_%i_profile",is);
26 TString pN = Form("%s_number_cells_vs_eta_in_layer_%i_profile",fN,is);
27 profileMap[pNm.Data()] = new TProfile(pN.Data(), "Number of cells;truth #eta", 90,-4.5,4.5, 0, 100);
28 pN = Form("%s%s",m_folder.c_str(),pNm.Data());
29 ATH_CHECK(m_rootHistSvc->regHist(pN.Data(), profileMap[pNm.Data()]));
30
31 pNm = Form("number_cells_vs_e_in_layer_%i_profile",is);
32 pN = Form("%s_number_cells_vs_e_in_layer_%i_profile",fN,is);
33 profileMap[pNm.Data()] = new TProfile(pN.Data(), "Number of cells;truth E [GeV]", 60, 0, 300, 0, 100);
34 pN = Form("%s%s",m_folder.c_str(),pNm.Data());
35 ATH_CHECK(m_rootHistSvc->regHist(pN.Data(), profileMap[pNm.Data()]));
36 }
37
38 for (int is = 0; is < 24; is++) {
39 TString pNm = Form("number_cells_vs_eta_in_sampling_%i_profile",is);
40 TString pN = Form("%s_number_cells_vs_eta_in_sampling_%i_profile",fN,is);
41 profileMap[pNm.Data()] = new TProfile(pN.Data(), "Number of cells;truth #eta", 90,-4.5,4.5, 0,100);
42 pN = Form("%s%s",m_folder.c_str(),pNm.Data());
43 ATH_CHECK(m_rootHistSvc->regHist(pN.Data(), profileMap[pNm.Data()]));
44 }
45
46 profileMap["Eraw_Etruth_vs_Etruth_profile"] = new TProfile(Form("%s_Eraw_Etruth_vs_Etruth_profile",fN), ";E^{truth};E^{raw}/E^{truth}", 100, 0., 200., 0.5, 1.5);
47 profileMap["Eraw_Etruth_vs_eta_profile"] = new TProfile(Form("%s_Eraw_Etruth_vs_eta_profile",fN), ";truth #eta;E^{raw}/E^{truth}", 90, -4.5, 4.5, 0.5, 1.5);
48
49 profileMap["number_topocluster_vs_e_profile"] = new TProfile(Form("%s_number_topocluster_vs_e_profile",fN), "Number of topocluster;truth E [GeV]", 60, 0, 300, -0.5, 14.5);
50 profileMap["number_topocluster_vs_eta_profile"] = new TProfile(Form("%s_number_topocluster_vs_eta_profile",fN), "Number of topocluster;truth #eta", 90,-4.5,4.5, -0.5, 14.5);
51
52 // do not undertstand this histo...
53 profileMap["number_cell_in_layer"] = new TProfile(Form("%s_number_cell_in_layer",fN), ";Number of cells;Layer", 100, 0, 200,0, 4);
54 histo2DMap["mu_energy_resolution_2D"] = new TH2D(Form("%s_mu_energy_resolution_2D",fN), ";<#mu>; Energy Resolution", 5, 0, 80, 20, -1, 1);
55
56 ATH_CHECK(m_rootHistSvc->regHist(m_folder+"Eraw_Etruth_vs_Etruth_profile", profileMap["Eraw_Etruth_vs_Etruth_profile"]));
57 ATH_CHECK(m_rootHistSvc->regHist(m_folder+"Eraw_Etruth_vs_eta_profile", profileMap["Eraw_Etruth_vs_eta_profile"]));
58
59 ATH_CHECK(m_rootHistSvc->regHist(m_folder+"number_topocluster_vs_e_profile", profileMap["number_topocluster_vs_e_profile"]));
60 ATH_CHECK(m_rootHistSvc->regHist(m_folder+"number_topocluster_vs_eta_profile", profileMap["number_topocluster_vs_eta_profile"]));
61
62 ATH_CHECK(m_rootHistSvc->regHist(m_folder+"number_cell_in_layer", profileMap["number_cell_in_layer"]));
63 ATH_CHECK(m_rootHistSvc->regHist(m_folder+"mu_energy_resolution_2D", histo2DMap["mu_energy_resolution_2D"]));
64
65 return StatusCode::SUCCESS;
66
67} // initializePlots
#define ATH_CHECK
Evaluate an expression and check for errors.

Member Data Documentation

◆ histo2DMap

std::map<std::string, TH2D* > egammaMonitoring::ClusterHistograms::histo2DMap

Definition at line 35 of file ClusterHistograms.h.

◆ m_folder

std::string egammaMonitoring::ClusterHistograms::m_folder
protected

Definition at line 45 of file ClusterHistograms.h.

◆ m_name

std::string egammaMonitoring::ClusterHistograms::m_name
protected

Definition at line 43 of file ClusterHistograms.h.

◆ m_rootHistSvc

SmartIF<ITHistSvc> egammaMonitoring::ClusterHistograms::m_rootHistSvc
protected

Definition at line 46 of file ClusterHistograms.h.

◆ m_title

std::string egammaMonitoring::ClusterHistograms::m_title
protected

Definition at line 44 of file ClusterHistograms.h.

◆ profileMap

std::map<std::string, TProfile* > egammaMonitoring::ClusterHistograms::profileMap

Definition at line 36 of file ClusterHistograms.h.


The documentation for this class was generated from the following files: