ATLAS Offline Software
Loading...
Searching...
No Matches
RecoMuonIDTrackPlots.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2020 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7RecoMuonIDTrackPlots::RecoMuonIDTrackPlots(PlotBase* pParent, const std::string& sDir):PlotBase(pParent, sDir),
8m_oAllPlots(this, "/", "Reco Muon"),
9m_oImpactPlots(this, "/"),
10m_oTrkRecoInfoPlots(this, "/"),
11m_oIDHitPlots(this,"/"),
12m_pt_broad(nullptr),
13m_eta_phi_broad(nullptr)
14{}
15
17 //booking histograms
18 m_pt_broad = Book1D("_pt_broad", "High p_{T} Distribution", 70, 100, 1500);
19 m_eta_phi_broad = Book2D("_eta_phi_broad", "High p_{T} Muon #eta #phi Distribution;;#eta;#phi", 52, -2.6, 2.6, 64, -3.2, 3.2);
20}
21
22//when the plot function called with a Muon Container
23//loop through each muon, get the corresponding link and fill it
24
25//when the plot function called with a Muon
26//get's the corresponding link and fill it
27void RecoMuonIDTrackPlots::fill(const xAOD::Muon& mu, int component){
28 if (component == 2 ){
29 const ElementLink<xAOD::TrackParticleContainer>& Mu_idtrack = mu.inDetTrackParticleLink();
30 if(Mu_idtrack.isValid()){
31 const xAOD::TrackParticle* trk = *Mu_idtrack;
32 fill(*trk);
33 }
34 }
35}
36
37
39 m_oAllPlots.fill(muTP);
40 m_oImpactPlots.fill(muTP);
41 m_oTrkRecoInfoPlots.fill(muTP);
42 m_oIDHitPlots.fill(muTP);
43 if (muTP.pt()/1000.0 > 100) {//ony for high pt muons
44 m_pt_broad->Fill(muTP.pt()/1000.0);
45 m_eta_phi_broad->Fill(muTP.eta(), muTP.phi());
46 }
47
48}
49
50
TH1D * Book1D(const std::string &name, const std::string &labels, int nBins, float start, float end, bool prependDir=true)
Book a TH1D histogram.
Definition PlotBase.cxx:94
PlotBase(PlotBase *parent, const std::string &sDir)
Definition PlotBase.cxx:29
TH2F * Book2D(const std::string &name, const std::string &labels, int nBinsX, float startX, float endX, int nBinsY, float startY, float endY, bool prependDir=true)
Book a TH2F histogram.
Definition PlotBase.cxx:123
Trk::IDHitPlots m_oIDHitPlots
void fill(const xAOD::Muon &mu, int component)
RecoMuonIDTrackPlots(PlotBase *pParent, const std::string &sDir)
Trk::ImpactPlots m_oImpactPlots
Trk::RecoInfoPlots m_oTrkRecoInfoPlots
virtual double phi() const override final
The azimuthal angle ( ) of the particle (has range to .)
virtual double pt() const override final
The transverse momentum ( ) of the particle.
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
Muon_v1 Muon
Reference the current persistent version: