ATLAS Offline Software
Loading...
Searching...
No Matches
IDHitSummaryPlots.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2017 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include <utility>
6
9
10namespace Muon{
11
13HitFracTypePlots::HitFracTypePlots(PlotBase* pParent,std::string sHitType, std::string sHitLabel):
14 PlotBase(pParent, ""), m_sHitType(std::move(sHitType)), m_sHitLabel(std::move(sHitLabel))
15{}
17{
18 fracHits = Book1D(m_sHitType, m_sHitLabel + ";" + m_sHitLabel + ";Entries", 30 ,0., 0.6);
19 fracHitsVsEta = Book2D(m_sHitType + "vsEta", m_sHitLabel + ";" + m_sHitLabel + ";#eta;Entries", 30 , 0., 0.3, 54, -2.7, 2.7);
20}
21 void HitFracTypePlots::fill(float fHits, float fEta, float weight)
22{
23 fracHits->Fill(fHits, weight);
24 // unclear how to weight this
25 fracHitsVsEta->Fill(fHits, fEta);
26}
27
28
29IDHitSummaryPlots::IDHitSummaryPlots(PlotBase* pParent, const std::string& sDir):
30 PlotBase(pParent, sDir)
31 , nBLayerHitsIfExpected(this, "nBLayerHitsIfExpected", "B-Layer clusters",0,4)
32 , nPixelHitsPlusDead(this, "nPixelHitsPlusDead","Pixel clusters",0,9)
33 , nSCTHitsPlusDead(this, "nSCTHitsPlusDead","SCT clusters",0,20)
34 , nTRTHitsPlusDead(this, "nTRTHitsPlusDead","TRT clusters",0,60)
35 , nTRTHitsPlusOutliers(this, "nTRTHitsPlusOutliers","TRT hits+Outliers",0,60)
36 , nPixSCTHoles(this, "nPixSCTHoles","Pix+SCT holes",0,5)
37 , fPixelOutliers(this, "fPixelOutliers","Fraction of Pix Outliers")
38 , fSCTOutliers(this, "fSCTOutliers","Fraction of SCT outliers")
39 , fTRTOutliers(this, "fTRTOutliers","Fraction of TRT Outliers")
40{}
41
42 void IDHitSummaryPlots::fill(const xAOD::TrackParticle& trk, float weight)
43{
44 float eta=trk.eta();
45 float phi=trk.phi();
46
47 uint8_t iBLayerHits(0),bExpectBLayerHit(0);
50 nBLayerHitsIfExpected.fill(bExpectBLayerHit*iBLayerHits,eta,phi,weight);
51
52 uint8_t iPixHits(0),iPixDead(0);
55 nPixelHitsPlusDead.fill(iPixHits+iPixDead,eta,phi,weight);
56
57 uint8_t iSCTHits(0),iSCTDead(0);
60 nSCTHitsPlusDead.fill(iSCTHits+iSCTDead,eta,phi,weight);
61
62 uint8_t iTRTHits(0),iTRTDead(0),iTRTOutliers(0);
66 nTRTHitsPlusDead.fill(iTRTHits+iTRTDead,eta,phi,weight);
67 nTRTHitsPlusOutliers.fill(iTRTHits+iTRTOutliers,eta,phi,weight);
68
69 uint8_t iPixHoles(0),iSCTHoles(0);
72 nPixSCTHoles.fill(iPixHoles+iSCTHoles,eta,phi,weight);
73
74 uint8_t iPixelOutliers(0);
75 trk.summaryValue(iPixelOutliers,xAOD::numberOfPixelOutliers);
76 float nPix=(int)iPixelOutliers+(int)iPixHits;
77 float fracPixOutliers = (nPix>0)? 1.*((int)iPixelOutliers)/nPix : 0;
78 fPixelOutliers.fill(fracPixOutliers,eta,weight);
79
80 uint8_t iSCTOutliers(0);
82 float nSCT=(int)iSCTOutliers+(int)iSCTHits;
83 float fracSCTOutliers = (nSCT>0)? 1.*((int)iSCTOutliers)/nSCT : 0;
84 fSCTOutliers.fill(fracSCTOutliers,eta,weight);
85
86 float nTRT=((int)iTRTOutliers)+((int)iTRTHits);
87 float fracTRTOutliers = (nTRT>0)? ((int)iTRTOutliers)/nTRT : 0;
88 fTRTOutliers.fill(fracTRTOutliers,eta,weight);
89}
90
91
92
93}
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
void fill(float hitval, float trketa, float weight=1.0)
HitFracTypePlots(PlotBase *pParent, std::string sHitType, std::string sHitLabel)
Trk::HitTypePlots nSCTHitsPlusDead
void fill(const xAOD::TrackParticle &trk, float weight=1.0)
Trk::HitTypePlots nPixSCTHoles
IDHitSummaryPlots(PlotBase *pParent, const std::string &sDir)
Trk::HitTypePlots nTRTHitsPlusDead
HitFracTypePlots fTRTOutliers
HitFracTypePlots fSCTOutliers
Trk::HitTypePlots nBLayerHitsIfExpected
Trk::HitTypePlots nPixelHitsPlusDead
HitFracTypePlots fPixelOutliers
Trk::HitTypePlots nTRTHitsPlusOutliers
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
virtual double phi() const override final
The azimuthal angle ( ) of the particle (has range to .)
bool summaryValue(uint8_t &value, const SummaryType &information) const
Accessor for TrackSummary values.
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
NRpcCablingAlg reads raw condition data and writes derived condition data to the condition store.
STL namespace.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
@ expectInnermostPixelLayerHit
Do we expect a 0th-layer barrel hit for this track?
@ numberOfTRTDeadStraws
number of dead TRT straws crossed [unit8_t].
@ numberOfPixelHoles
number of pixel layers on track with absence of hits [unit8_t].
@ numberOfPixelOutliers
these are the pixel outliers, including the b-layer [unit8_t].
@ numberOfTRTHits
number of TRT hits [unit8_t].
@ numberOfSCTDeadSensors
number of dead SCT sensors crossed [unit8_t].
@ numberOfSCTHits
number of hits in SCT [unit8_t].
@ numberOfSCTOutliers
number of SCT outliers [unit8_t].
@ numberOfInnermostPixelLayerHits
these are the hits in the 0th pixel barrel layer
@ numberOfPixelHits
these are the pixel hits, including the b-layer [unit8_t].
@ numberOfTRTOutliers
number of TRT outliers [unit8_t].
@ numberOfPixelDeadSensors
number of dead pixel sensors crossed [unit8_t].
@ numberOfSCTHoles
number of SCT holes [unit8_t].