ATLAS Offline Software
RecoMuonTrackPlots.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 
7 RecoMuonTrackPlots::RecoMuonTrackPlots(PlotBase* pParent, const std::string& sDir):PlotBase(pParent, sDir),
8 m_oAllPlots(this, "/", "Reco Muon"),
9 m_oImpactPlots(this, "/"),
10 m_oTrkRecoInfoPlots(this, "/"),
11 m_pt_broad(nullptr),
12 m_eta_phi_broad(nullptr)
13 {}
14 
16  //booking histograms
17  m_pt_broad = Book1D("_pt_broad", "High p_{T} Distribution", 70, 100, 1500);
18  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);
19 }
20 
21 //when the plot function called with a Muon Container
22 //loop through each muon, get the corresponding link and fill it
23 
24 //when the plot function called with a Muon
25 //get's the corresponding link and fill it
26 void RecoMuonTrackPlots::fill(const xAOD::Muon& mu, int component){
27 
28  if (component == 0 ) {
29  const ElementLink<xAOD::TrackParticleContainer>& Mu_MStrack = mu.muonSpectrometerTrackParticleLink();
30  if(Mu_MStrack.isValid()){
31  const xAOD::TrackParticle* trk = *Mu_MStrack;
32  fill(*trk);
33  }
34  }
35 
36  if (component == 1 ){
37 
38  const ElementLink<xAOD::TrackParticleContainer>& Mu_metrack = mu.trackParticleLink(xAOD::Muon::TrackParticleType::ExtrapolatedMuonSpectrometerTrackParticle);
39  if(Mu_metrack.isValid()){
40  const xAOD::TrackParticle* trk = *Mu_metrack;
41  fill(*trk);
42  }
43  }
44 }
45 
46 
48  m_oAllPlots.fill(muTP);
49  m_oImpactPlots.fill(muTP);
51 
52  if (muTP.pt()/1000.0 > 100) {//ony for high pt muons
53  m_pt_broad->Fill(muTP.pt()/1000.0);
54  m_eta_phi_broad->Fill(muTP.eta(), muTP.phi());
55  }
56 
57 }
58 
59 
Trk::ParamPlots::fill
void fill(const xAOD::IParticle &prt, float weight=1.0)
Definition: Tracking/TrkValidation/TrkValHistUtils/Root/ParamPlots.cxx:40
xAOD::TrackParticle_v1::pt
virtual double pt() const override final
The transverse momentum ( ) of the particle.
Definition: TrackParticle_v1.cxx:73
RecoMuonTrackPlots::m_pt_broad
TH1 * m_pt_broad
Definition: RecoMuonTrackPlots.h:34
PlotBase
Definition: PlotBase.h:34
RecoMuonTrackPlots::m_oImpactPlots
Trk::ImpactPlots m_oImpactPlots
Definition: RecoMuonTrackPlots.h:31
python.copyTCTOutput.sDir
sDir
Definition: copyTCTOutput.py:60
xAOD::TrackParticle_v1::eta
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
Definition: TrackParticle_v1.cxx:77
RecoMuonTrackPlots::m_oTrkRecoInfoPlots
Trk::RecoInfoPlots m_oTrkRecoInfoPlots
Definition: RecoMuonTrackPlots.h:32
PlotBase::Book2D
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
xAOD::Muon_v1
Class describing a Muon.
Definition: Muon_v1.h:38
PlotBase::Book1D
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
RecoMuonTrackPlots.h
RecoMuonTrackPlots::m_oAllPlots
Trk::ParamPlots m_oAllPlots
Definition: RecoMuonTrackPlots.h:30
RecoMuonTrackPlots::m_eta_phi_broad
TH2 * m_eta_phi_broad
Definition: RecoMuonTrackPlots.h:35
RecoMuonTrackPlots::initializePlots
void initializePlots()
Definition: RecoMuonTrackPlots.cxx:15
Trk::RecoInfoPlots::fill
void fill(const xAOD::TrackParticle &trkprt, float weight=1.0)
Definition: Tracking/TrkValidation/TrkValHistUtils/Root/RecoInfoPlots.cxx:28
Trk::ImpactPlots::fill
void fill(const xAOD::TrackParticle &trkprt, float weight=1.0)
Definition: ImpactPlots.cxx:27
RecoMuonTrackPlots::RecoMuonTrackPlots
RecoMuonTrackPlots(PlotBase *pParent, const std::string &sDir)
Definition: RecoMuonTrackPlots.cxx:7
RecoMuonTrackPlots::fill
void fill(const xAOD::Muon &mu, int component)
Definition: RecoMuonTrackPlots.cxx:26
xAOD::TrackParticle_v1
Class describing a TrackParticle.
Definition: TrackParticle_v1.h:43
CaloNoise_fillDB.mu
mu
Definition: CaloNoise_fillDB.py:53
xAOD::TrackParticle_v1::phi
virtual double phi() const override final
The azimuthal angle ( ) of the particle (has range to .)