CMS 3D CMS Logo

List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes
CSCSegmentValidation Class Reference

#include <CSCSegmentValidation.h>

Inheritance diagram for CSCSegmentValidation:
CSCBaseValidation

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void bookHistograms (DQMStore::IBooker &)
 
 CSCSegmentValidation (const edm::ParameterSet &, edm::ConsumesCollector &&iC)
 
 ~CSCSegmentValidation () override
 
- Public Member Functions inherited from CSCBaseValidation
 CSCBaseValidation (const edm::ParameterSet &ps)
 
void setGeometry (const CSCGeometry *geom)
 
void setSimHitMap (const PSimHitMap *simHitMap)
 
virtual ~CSCBaseValidation ()
 

Private Types

typedef std::map< int, std::vector< int > > ChamberHitMap
 
typedef std::map< int, std::vector< CSCSegment > > ChamberSegmentMap
 

Private Member Functions

void fillEfficiencyPlots ()
 
void fillLayerHitsPerChamber ()
 
bool hasSegment (int chamberId) const
 
const PSimHitkeyHit (int chamberId) const
 
void plotResolution (const PSimHit &simHit, const CSCSegment &recHit, const CSCLayer *layer, int chamberType)
 

Private Attributes

edm::InputTag inputTag_
 
edm::EDGetTokenT< CSCSegmentCollectionsegments_Token_
 
ChamberSegmentMap theChamberSegmentMap
 
MonitorElementthedXdZPullPlots [10]
 
MonitorElementthedXdZResolutionPlots [10]
 
MonitorElementthedYdZPullPlots [10]
 
MonitorElementthedYdZResolutionPlots [10]
 
ChamberHitMap theLayerHitsPerChamber
 
MonitorElementtheNPerChamberTypePlot
 
MonitorElementtheNPerEventPlot
 
MonitorElementtheNRecHitsPlot
 
MonitorElementtheRdPhiPullPlots [10]
 
MonitorElementtheRdPhiResolutionPlots [10]
 
int theShowerThreshold
 
MonitorElementtheThetaResolutionPlots [10]
 
MonitorElementtheTypePlot4HitsNoShower
 
MonitorElementtheTypePlot4HitsNoShowerSeg
 
MonitorElementtheTypePlot4HitsShower
 
MonitorElementtheTypePlot4HitsShowerSeg
 
MonitorElementtheTypePlot5HitsNoShower
 
MonitorElementtheTypePlot5HitsNoShowerSeg
 
MonitorElementtheTypePlot5HitsShower
 
MonitorElementtheTypePlot5HitsShowerSeg
 
MonitorElementtheTypePlot6HitsNoShower
 
MonitorElementtheTypePlot6HitsNoShowerSeg
 
MonitorElementtheTypePlot6HitsShower
 
MonitorElementtheTypePlot6HitsShowerSeg
 

Additional Inherited Members

- Public Types inherited from CSCBaseValidation
typedef dqm::legacy::DQMStore DQMStore
 
typedef dqm::legacy::MonitorElement MonitorElement
 
- Protected Member Functions inherited from CSCBaseValidation
const CSCLayerfindLayer (int detId) const
 
bool isSimTrackGood (const SimTrack &t) const
 
- Protected Attributes inherited from CSCBaseValidation
bool doSim_
 
double simTrackMaxEta_
 
double simTrackMinEta_
 
double simTrackMinPt_
 
const CSCGeometrytheCSCGeometry
 
const PSimHitMaptheSimHitMap
 

Detailed Description

Definition at line 8 of file CSCSegmentValidation.h.

Member Typedef Documentation

◆ ChamberHitMap

typedef std::map<int, std::vector<int> > CSCSegmentValidation::ChamberHitMap
private

Definition at line 26 of file CSCSegmentValidation.h.

◆ ChamberSegmentMap

typedef std::map<int, std::vector<CSCSegment> > CSCSegmentValidation::ChamberSegmentMap
private

Definition at line 33 of file CSCSegmentValidation.h.

Constructor & Destructor Documentation

◆ CSCSegmentValidation()

CSCSegmentValidation::CSCSegmentValidation ( const edm::ParameterSet ps,
edm::ConsumesCollector &&  iC 
)

Definition at line 5 of file CSCSegmentValidation.cc.

References edm::ParameterSet::getParameterSet(), inputTag_, muonDTDigis_cfi::pset, and segments_Token_.

7  const auto &pset = ps.getParameterSet("cscSegment");
8  inputTag_ = pset.getParameter<edm::InputTag>("inputTag");
10 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
ParameterSet const & getParameterSet(std::string const &) const
edm::EDGetTokenT< CSCSegmentCollection > segments_Token_
CSCBaseValidation(const edm::ParameterSet &ps)
ChamberHitMap theLayerHitsPerChamber
ChamberSegmentMap theChamberSegmentMap

◆ ~CSCSegmentValidation()

CSCSegmentValidation::~CSCSegmentValidation ( )
override

Definition at line 12 of file CSCSegmentValidation.cc.

12 {}

Member Function Documentation

◆ analyze()

void CSCSegmentValidation::analyze ( const edm::Event e,
const edm::EventSetup eventSetup 
)
overridevirtual

Implements CSCBaseValidation.

Definition at line 53 of file CSCSegmentValidation.cc.

References nano_mu_digi_cff::chamberType, hcalRecHitTable_cff::detId, MillePedeFileConverter_cfg::e, dqm::impl::MonitorElement::Fill(), fillEfficiencyPlots(), fillLayerHitsPerChamber(), CSCBaseValidation::findLayer(), keyHit(), nano_mu_digi_cff::layer, plotResolution(), edm::Handle< T >::product(), segments_Token_, theChamberSegmentMap, theNPerChamberTypePlot, theNPerEventPlot, and theNRecHitsPlot.

53  {
54  // get the collection of CSCRecHsegmentItrD
56  e.getByToken(segments_Token_, hRecHits);
57  const CSCSegmentCollection *cscRecHits = hRecHits.product();
58 
59  theChamberSegmentMap.clear();
60  unsigned nPerEvent = 0;
61  for (auto segmentItr = cscRecHits->begin(); segmentItr != cscRecHits->end(); segmentItr++) {
62  ++nPerEvent;
63  int detId = segmentItr->geographicalId().rawId();
64  int chamberType = segmentItr->cscDetId().iChamberType();
65 
66  theNRecHitsPlot->Fill(segmentItr->nRecHits());
68  theChamberSegmentMap[detId].push_back(*segmentItr);
69 
70  // do the resolution plots
71  const PSimHit *hit = keyHit(detId);
72  if (hit != nullptr) {
73  const CSCLayer *layer = findLayer(hit->detUnitId());
74  plotResolution(*hit, *segmentItr, layer, chamberType);
75  }
76  }
77 
78  theNPerEventPlot->Fill(nPerEvent);
79 
82 }
MonitorElement * theNPerEventPlot
void plotResolution(const PSimHit &simHit, const CSCSegment &recHit, const CSCLayer *layer, int chamberType)
MonitorElement * theNPerChamberTypePlot
T const * product() const
Definition: Handle.h:70
edm::EDGetTokenT< CSCSegmentCollection > segments_Token_
MonitorElement * theNRecHitsPlot
void Fill(long long x)
const CSCLayer * findLayer(int detId) const
const PSimHit * keyHit(int chamberId) const
ChamberSegmentMap theChamberSegmentMap

◆ bookHistograms()

void CSCSegmentValidation::bookHistograms ( DQMStore::IBooker iBooker)

Definition at line 14 of file CSCSegmentValidation.cc.

References dqm::implementation::IBooker::book1D(), CSCDetId::chamberName(), mps_fire::i, AlCaHLTBitMon_QueryRunRegistry::string, RandomServiceHelper::t1, RandomServiceHelper::t2, RandomServiceHelper::t3, thedXdZPullPlots, thedXdZResolutionPlots, thedYdZPullPlots, thedYdZResolutionPlots, theNPerChamberTypePlot, theNPerEventPlot, theNRecHitsPlot, theRdPhiPullPlots, theRdPhiResolutionPlots, theThetaResolutionPlots, theTypePlot4HitsNoShower, theTypePlot4HitsNoShowerSeg, theTypePlot4HitsShower, theTypePlot4HitsShowerSeg, theTypePlot5HitsNoShower, theTypePlot5HitsNoShowerSeg, theTypePlot5HitsShower, theTypePlot5HitsShowerSeg, theTypePlot6HitsNoShower, theTypePlot6HitsNoShowerSeg, theTypePlot6HitsShower, and theTypePlot6HitsShowerSeg.

14  {
15  theNPerEventPlot = iBooker.book1D("CSCSegmentsPerEvent", "Number of CSC segments per event", 100, 0, 50);
16  theNRecHitsPlot = iBooker.book1D("CSCRecHitsPerSegment", "Number of CSC rec hits per segment", 8, 0, 7);
18  iBooker.book1D("CSCSegmentsPerChamberType", "Number of CSC segments per chamber type", 11, 0, 10);
19  theTypePlot4HitsNoShower = iBooker.book1D("CSCSegments4HitsNoShower", "", 100, 0, 10);
20  theTypePlot4HitsNoShowerSeg = iBooker.book1D("CSCSegments4HitsNoShowerSeg", "", 100, 0, 10);
21  theTypePlot4HitsShower = iBooker.book1D("CSCSegments4HitsShower", "", 100, 0, 10);
22  theTypePlot4HitsShowerSeg = iBooker.book1D("CSCSegments4HitsShowerSeg", "", 100, 0, 10);
23  theTypePlot5HitsNoShower = iBooker.book1D("CSCSegments5HitsNoShower", "", 100, 0, 10);
24  theTypePlot5HitsNoShowerSeg = iBooker.book1D("CSCSegments5HitsNoShowerSeg", "", 100, 0, 10);
25  theTypePlot5HitsShower = iBooker.book1D("CSCSegments5HitsShower", "", 100, 0, 10);
26  theTypePlot5HitsShowerSeg = iBooker.book1D("CSCSegments5HitsShowerSeg", "", 100, 0, 10);
27  theTypePlot6HitsNoShower = iBooker.book1D("CSCSegments6HitsNoShower", "", 100, 0, 10);
28  theTypePlot6HitsNoShowerSeg = iBooker.book1D("CSCSegments6HitsNoShowerSeg", "", 100, 0, 10);
29  theTypePlot6HitsShower = iBooker.book1D("CSCSegments6HitsShower", "", 100, 0, 10);
30  theTypePlot6HitsShowerSeg = iBooker.book1D("CSCSegments6HitsShowerSeg", "", 100, 0, 10);
31 
32  for (int i = 1; i <= 10; ++i) {
34 
35  const std::string t1("CSCSegmentRdPhiResolution_" + cn);
36  const std::string t2("CSCSegmentRdPhiPull_" + cn);
37  const std::string t3("CSCSegmentThetaResolution_" + cn);
38  const std::string t5("CSCSegmentdXdZResolution_" + cn);
39  const std::string t6("CSCSegmentdXdZPull_" + cn);
40  const std::string t7("CSCSegmentdYdZResolution_" + cn);
41  const std::string t8("CSCSegmentdYdZPull_" + cn);
42 
43  theRdPhiResolutionPlots[i - 1] = iBooker.book1D(t1, t1, 100, -0.4, 0.4);
44  theRdPhiPullPlots[i - 1] = iBooker.book1D(t2, t2, 100, -5, 5);
45  theThetaResolutionPlots[i - 1] = iBooker.book1D(t3, t3, 100, -1, 1);
46  thedXdZResolutionPlots[i - 1] = iBooker.book1D(t5, t5, 100, -1, 1);
47  thedXdZPullPlots[i - 1] = iBooker.book1D(t6, t6, 100, -5, 5);
48  thedYdZResolutionPlots[i - 1] = iBooker.book1D(t7, t7, 100, -1, 1);
49  thedYdZPullPlots[i - 1] = iBooker.book1D(t8, t8, 100, -5, 5);
50  }
51 }
MonitorElement * theNPerEventPlot
MonitorElement * theTypePlot5HitsNoShowerSeg
MonitorElement * theTypePlot4HitsShowerSeg
MonitorElement * theNPerChamberTypePlot
MonitorElement * theTypePlot6HitsShower
MonitorElement * theTypePlot5HitsNoShower
MonitorElement * theNRecHitsPlot
MonitorElement * theTypePlot4HitsShower
MonitorElement * theRdPhiResolutionPlots[10]
MonitorElement * theTypePlot6HitsShowerSeg
MonitorElement * theTypePlot4HitsNoShower
std::string chamberName() const
Definition: CSCDetId.cc:92
MonitorElement * theTypePlot6HitsNoShowerSeg
MonitorElement * theTypePlot6HitsNoShower
MonitorElement * thedYdZResolutionPlots[10]
MonitorElement * theRdPhiPullPlots[10]
MonitorElement * thedXdZPullPlots[10]
MonitorElement * theTypePlot5HitsShowerSeg
MonitorElement * theThetaResolutionPlots[10]
MonitorElement * theTypePlot4HitsNoShowerSeg
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
MonitorElement * theTypePlot5HitsShower
MonitorElement * thedXdZResolutionPlots[10]
MonitorElement * thedYdZPullPlots[10]

◆ fillEfficiencyPlots()

void CSCSegmentValidation::fillEfficiencyPlots ( )
private

Definition at line 84 of file CSCSegmentValidation.cc.

References nano_mu_digi_cff::chamberType, dqm::impl::MonitorElement::Fill(), hasSegment(), CSCDetId::iChamberType(), jetUpdater_cfi::sort, theLayerHitsPerChamber, theShowerThreshold, theTypePlot4HitsNoShower, theTypePlot4HitsNoShowerSeg, theTypePlot4HitsShower, theTypePlot4HitsShowerSeg, theTypePlot5HitsNoShower, theTypePlot5HitsNoShowerSeg, theTypePlot5HitsShower, theTypePlot5HitsShowerSeg, theTypePlot6HitsNoShower, theTypePlot6HitsNoShowerSeg, theTypePlot6HitsShower, theTypePlot6HitsShowerSeg, tier0::unique(), and findQualityFiles::v.

Referenced by analyze().

84  {
85  // now plot efficiency by looping over all chambers with hits
86  for (auto mapItr = theLayerHitsPerChamber.begin(), mapEnd = theLayerHitsPerChamber.end(); mapItr != mapEnd;
87  ++mapItr) {
88  int chamberId = mapItr->first;
89  int nHitsInChamber = mapItr->second.size();
90  bool isShower = (nHitsInChamber > theShowerThreshold);
91  bool hasSeg = hasSegment(chamberId);
92  int chamberType = CSCDetId(chamberId).iChamberType();
93  // find how many layers were hit in this chamber
94  std::vector<int> v = mapItr->second;
95  std::sort(v.begin(), v.end());
96  // maybe can just count
97  v.erase(std::unique(v.begin(), v.end()), v.end());
98  int nLayersHit = v.size();
99 
100  if (nLayersHit == 4) {
101  if (isShower)
103  else
105 
106  if (hasSeg) {
107  if (isShower)
109  else
111  }
112  }
113 
114  if (nLayersHit == 5) {
115  if (isShower)
117  else
119 
120  if (hasSeg) {
121  if (isShower)
123  else
125  }
126  }
127 
128  if (nLayersHit == 6) {
129  if (isShower)
131  else
133 
134  if (hasSeg) {
135  if (isShower)
137  else
139  }
140  }
141  }
142 }
MonitorElement * theTypePlot5HitsNoShowerSeg
MonitorElement * theTypePlot4HitsShowerSeg
MonitorElement * theTypePlot6HitsShower
MonitorElement * theTypePlot5HitsNoShower
MonitorElement * theTypePlot4HitsShower
unsigned short iChamberType() const
Definition: CSCDetId.h:96
MonitorElement * theTypePlot6HitsShowerSeg
void Fill(long long x)
bool hasSegment(int chamberId) const
MonitorElement * theTypePlot4HitsNoShower
def unique(seq, keepstr=True)
Definition: tier0.py:24
MonitorElement * theTypePlot6HitsNoShowerSeg
MonitorElement * theTypePlot6HitsNoShower
ChamberHitMap theLayerHitsPerChamber
MonitorElement * theTypePlot5HitsShowerSeg
MonitorElement * theTypePlot4HitsNoShowerSeg
MonitorElement * theTypePlot5HitsShower

◆ fillLayerHitsPerChamber()

void CSCSegmentValidation::fillLayerHitsPerChamber ( )
private

Definition at line 178 of file CSCSegmentValidation.cc.

References CSCDetId::chamberId(), PSimHitMap::detsWithHits(), PSimHitMap::hits(), mps_fire::i, TrackingDataMCValidation_Standalone_cff::nhits, DetId::rawId(), theLayerHitsPerChamber, and CSCBaseValidation::theSimHitMap.

Referenced by analyze().

178  {
179  theLayerHitsPerChamber.clear();
180  std::vector<int> layersHit = theSimHitMap->detsWithHits();
181  for (auto layerItr = layersHit.begin(), layersHitEnd = layersHit.end(); layerItr != layersHitEnd; ++layerItr) {
182  CSCDetId layerId(*layerItr);
183  CSCDetId chamberId = layerId.chamberId();
184  int nhits = theSimHitMap->hits(*layerItr).size();
185  // multiple entries, so we can see showers
186  for (int i = 0; i < nhits; ++i) {
187  theLayerHitsPerChamber[chamberId.rawId()].push_back(*layerItr);
188  }
189  }
190 }
const PSimHitMap * theSimHitMap
const edm::PSimHitContainer & hits(int detId) const
Definition: PSimHitMap.cc:20
ChamberHitMap theLayerHitsPerChamber
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
CSCDetId chamberId() const
Definition: CSCDetId.h:47
std::vector< int > detsWithHits() const
Definition: PSimHitMap.cc:29

◆ hasSegment()

bool CSCSegmentValidation::hasSegment ( int  chamberId) const
private

Definition at line 144 of file CSCSegmentValidation.cc.

References theChamberSegmentMap.

Referenced by fillEfficiencyPlots().

144  {
145  return (theChamberSegmentMap.find(chamberId) != theChamberSegmentMap.end());
146 }
ChamberSegmentMap theChamberSegmentMap

◆ keyHit()

const PSimHit * CSCSegmentValidation::keyHit ( int  chamberId) const
private

Definition at line 192 of file CSCSegmentValidation.cc.

References PSimHitMap::hits(), LaserDQM_cfg::p1, SiStripOfflineCRack_cfg::p2, mps_fire::result, and CSCBaseValidation::theSimHitMap.

Referenced by analyze().

192  {
193  auto SimHitPabsLessThan = [](const PSimHit &p1, const PSimHit &p2) -> bool { return p1.pabs() < p2.pabs(); };
194 
195  const PSimHit *result = nullptr;
196  int layerId = chamberId + 3;
197  const auto &layerHits = theSimHitMap->hits(layerId);
198 
199  if (!layerHits.empty()) {
200  // pick the hit with maximum energy
201  auto hitItr = std::max_element(layerHits.begin(), layerHits.end(), SimHitPabsLessThan);
202  result = &(*hitItr);
203  }
204  return result;
205 }
const PSimHitMap * theSimHitMap
const edm::PSimHitContainer & hits(int detId) const
Definition: PSimHitMap.cc:20

◆ plotResolution()

void CSCSegmentValidation::plotResolution ( const PSimHit simHit,
const CSCSegment recHit,
const CSCLayer layer,
int  chamberType 
)
private

Definition at line 148 of file CSCSegmentValidation.cc.

References nano_mu_digi_cff::chamberType, dqm::impl::MonitorElement::Fill(), nano_mu_digi_cff::layer, CSCSegment::localDirection(), CSCSegment::localDirectionError(), CSCSegment::localPosition(), CSCSegment::localPositionError(), PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), rpcPointValidation_cfi::simHit, mathSSE::sqrt(), thedXdZPullPlots, thedXdZResolutionPlots, thedYdZPullPlots, thedYdZResolutionPlots, theRdPhiPullPlots, theRdPhiResolutionPlots, PV3DBase< T, PVType, FrameType >::theta(), theThetaResolutionPlots, PV3DBase< T, PVType, FrameType >::x(), LocalError::xx(), PV3DBase< T, PVType, FrameType >::y(), LocalError::yy(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by analyze().

151  {
152  GlobalPoint simHitPos = layer->toGlobal(simHit.localPosition());
153  GlobalPoint segmentPos = layer->toGlobal(segment.localPosition());
154  LocalVector simHitDir = simHit.localDirection();
155  LocalVector segmentDir = segment.localDirection();
156 
157  double dphi = segmentPos.phi() - simHitPos.phi();
158  double rdphi = segmentPos.perp() * dphi;
159  double dtheta = segmentPos.theta() - simHitPos.theta();
160 
161  double sigmax = sqrt(segment.localPositionError().xx());
162 
163  double ddxdz = segmentDir.x() / segmentDir.z() - simHitDir.x() / simHitDir.z();
164  double ddydz = segmentDir.y() / segmentDir.z() - simHitDir.y() / simHitDir.z();
165  double sigmadxdz = sqrt(segment.localDirectionError().xx());
166  double sigmadydz = sqrt(segment.localDirectionError().yy());
167 
169  theRdPhiPullPlots[chamberType - 1]->Fill(rdphi / sigmax);
171 
173  thedXdZPullPlots[chamberType - 1]->Fill(ddxdz / sigmadxdz);
175  thedYdZPullPlots[chamberType - 1]->Fill(ddydz / sigmadydz);
176 }
T perp() const
Definition: PV3DBase.h:69
T z() const
Definition: PV3DBase.h:61
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
MonitorElement * theRdPhiResolutionPlots[10]
void Fill(long long x)
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
T sqrt(T t)
Definition: SSEVec.h:19
MonitorElement * thedYdZResolutionPlots[10]
MonitorElement * theRdPhiPullPlots[10]
MonitorElement * thedXdZPullPlots[10]
MonitorElement * theThetaResolutionPlots[10]
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
MonitorElement * thedXdZResolutionPlots[10]
MonitorElement * thedYdZPullPlots[10]

Member Data Documentation

◆ inputTag_

edm::InputTag CSCSegmentValidation::inputTag_
private

Definition at line 21 of file CSCSegmentValidation.h.

Referenced by CSCSegmentValidation().

◆ segments_Token_

edm::EDGetTokenT<CSCSegmentCollection> CSCSegmentValidation::segments_Token_
private

Definition at line 20 of file CSCSegmentValidation.h.

Referenced by analyze(), and CSCSegmentValidation().

◆ theChamberSegmentMap

ChamberSegmentMap CSCSegmentValidation::theChamberSegmentMap
private

Definition at line 34 of file CSCSegmentValidation.h.

Referenced by analyze(), and hasSegment().

◆ thedXdZPullPlots

MonitorElement* CSCSegmentValidation::thedXdZPullPlots[10]
private

Definition at line 45 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and plotResolution().

◆ thedXdZResolutionPlots

MonitorElement* CSCSegmentValidation::thedXdZResolutionPlots[10]
private

Definition at line 44 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and plotResolution().

◆ thedYdZPullPlots

MonitorElement* CSCSegmentValidation::thedYdZPullPlots[10]
private

Definition at line 47 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and plotResolution().

◆ thedYdZResolutionPlots

MonitorElement* CSCSegmentValidation::thedYdZResolutionPlots[10]
private

Definition at line 46 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and plotResolution().

◆ theLayerHitsPerChamber

ChamberHitMap CSCSegmentValidation::theLayerHitsPerChamber
private

Definition at line 27 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots(), and fillLayerHitsPerChamber().

◆ theNPerChamberTypePlot

MonitorElement* CSCSegmentValidation::theNPerChamberTypePlot
private

Definition at line 40 of file CSCSegmentValidation.h.

Referenced by analyze(), and bookHistograms().

◆ theNPerEventPlot

MonitorElement* CSCSegmentValidation::theNPerEventPlot
private

Definition at line 38 of file CSCSegmentValidation.h.

Referenced by analyze(), and bookHistograms().

◆ theNRecHitsPlot

MonitorElement* CSCSegmentValidation::theNRecHitsPlot
private

Definition at line 39 of file CSCSegmentValidation.h.

Referenced by analyze(), and bookHistograms().

◆ theRdPhiPullPlots

MonitorElement* CSCSegmentValidation::theRdPhiPullPlots[10]
private

Definition at line 42 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and plotResolution().

◆ theRdPhiResolutionPlots

MonitorElement* CSCSegmentValidation::theRdPhiResolutionPlots[10]
private

Definition at line 41 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and plotResolution().

◆ theShowerThreshold

int CSCSegmentValidation::theShowerThreshold
private

Definition at line 36 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

◆ theThetaResolutionPlots

MonitorElement* CSCSegmentValidation::theThetaResolutionPlots[10]
private

Definition at line 43 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and plotResolution().

◆ theTypePlot4HitsNoShower

MonitorElement* CSCSegmentValidation::theTypePlot4HitsNoShower
private

Definition at line 49 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and fillEfficiencyPlots().

◆ theTypePlot4HitsNoShowerSeg

MonitorElement* CSCSegmentValidation::theTypePlot4HitsNoShowerSeg
private

Definition at line 50 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and fillEfficiencyPlots().

◆ theTypePlot4HitsShower

MonitorElement* CSCSegmentValidation::theTypePlot4HitsShower
private

Definition at line 51 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and fillEfficiencyPlots().

◆ theTypePlot4HitsShowerSeg

MonitorElement* CSCSegmentValidation::theTypePlot4HitsShowerSeg
private

Definition at line 52 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and fillEfficiencyPlots().

◆ theTypePlot5HitsNoShower

MonitorElement* CSCSegmentValidation::theTypePlot5HitsNoShower
private

Definition at line 53 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and fillEfficiencyPlots().

◆ theTypePlot5HitsNoShowerSeg

MonitorElement* CSCSegmentValidation::theTypePlot5HitsNoShowerSeg
private

Definition at line 54 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and fillEfficiencyPlots().

◆ theTypePlot5HitsShower

MonitorElement* CSCSegmentValidation::theTypePlot5HitsShower
private

Definition at line 55 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and fillEfficiencyPlots().

◆ theTypePlot5HitsShowerSeg

MonitorElement* CSCSegmentValidation::theTypePlot5HitsShowerSeg
private

Definition at line 56 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and fillEfficiencyPlots().

◆ theTypePlot6HitsNoShower

MonitorElement* CSCSegmentValidation::theTypePlot6HitsNoShower
private

Definition at line 57 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and fillEfficiencyPlots().

◆ theTypePlot6HitsNoShowerSeg

MonitorElement* CSCSegmentValidation::theTypePlot6HitsNoShowerSeg
private

Definition at line 58 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and fillEfficiencyPlots().

◆ theTypePlot6HitsShower

MonitorElement* CSCSegmentValidation::theTypePlot6HitsShower
private

Definition at line 59 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and fillEfficiencyPlots().

◆ theTypePlot6HitsShowerSeg

MonitorElement* CSCSegmentValidation::theTypePlot6HitsShowerSeg
private

Definition at line 60 of file CSCSegmentValidation.h.

Referenced by bookHistograms(), and fillEfficiencyPlots().