CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Types | Private Member Functions | Static Private Member Functions | Private Attributes
CSCSegmentValidation Class Reference

#include <CSCSegmentValidation.h>

Inheritance diagram for CSCSegmentValidation:
CSCBaseValidation

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 
 CSCSegmentValidation (DQMStore *dbe, const edm::InputTag &inputTag, edm::ConsumesCollector &&iC)
 
virtual ~CSCSegmentValidation ()
 
- Public Member Functions inherited from CSCBaseValidation
 CSCBaseValidation (DQMStore *dbe, const edm::InputTag &inputTag)
 
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)
 

Static Private Member Functions

static int whatChamberType (int detId)
 

Private Attributes

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

Additional Inherited Members

- Protected Member Functions inherited from CSCBaseValidation
const CSCLayerfindLayer (int detId) const
 
- Protected Attributes inherited from CSCBaseValidation
DQMStoredbe_
 
const CSCGeometrytheCSCGeometry
 
edm::InputTag theInputTag
 
const PSimHitMaptheSimHitMap
 

Detailed Description

Definition at line 12 of file CSCSegmentValidation.h.

Member Typedef Documentation

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

Definition at line 31 of file CSCSegmentValidation.h.

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

Definition at line 38 of file CSCSegmentValidation.h.

Constructor & Destructor Documentation

CSCSegmentValidation::CSCSegmentValidation ( DQMStore dbe,
const edm::InputTag inputTag,
edm::ConsumesCollector &&  iC 
)

Definition at line 7 of file CSCSegmentValidation.cc.

References DQMStore::book1D(), CSCBaseValidation::dbe_, i, segments_Token_, DQMStore::setCurrentFolder(), thedXdZPullPlots, thedXdZResolutionPlots, thedYdZPullPlots, thedYdZResolutionPlots, theRdPhiPullPlots, theRdPhiResolutionPlots, theThetaPullPlots, and theThetaResolutionPlots.

8 : CSCBaseValidation(dbe, inputTag),
12  theNPerEventPlot( dbe_->book1D("CSCSegmentsPerEvent", "Number of CSC segments per event", 100, 0, 50) ),
13  theNRecHitsPlot( dbe_->book1D("CSCRecHitsPerSegment", "Number of CSC rec hits per segment" , 8, 0, 7) ),
14  theNPerChamberTypePlot( dbe_->book1D("CSCSegmentsPerChamberType", "Number of CSC segments per chamber type", 11, 0, 10) ),
15  theTypePlot4HitsNoShower( dbe_->book1D("CSCSegments4HitsNoShower", "", 100, 0, 10) ),
16  theTypePlot4HitsNoShowerSeg( dbe_->book1D("CSCSegments4HitsNoShowerSeg", "", 100, 0, 10) ),
17  theTypePlot4HitsShower( dbe_->book1D("CSCSegments4HitsShower", "", 100, 0, 10) ),
18  theTypePlot4HitsShowerSeg( dbe_->book1D("CSCSegments4HitsShowerSeg", "", 100, 0, 10) ),
19  theTypePlot5HitsNoShower( dbe_->book1D("CSCSegments5HitsNoShower", "", 100, 0, 10) ),
20  theTypePlot5HitsNoShowerSeg( dbe_->book1D("CSCSegments5HitsNoShowerSeg", "", 100, 0, 10) ),
21  theTypePlot5HitsShower( dbe_->book1D("CSCSegments5HitsShower", "", 100, 0, 10) ),
22  theTypePlot5HitsShowerSeg( dbe_->book1D("CSCSegments5HitsShowerSeg", "", 100, 0, 10) ),
23  theTypePlot6HitsNoShower( dbe_->book1D("CSCSegments6HitsNoShower", "", 100, 0, 10) ),
24  theTypePlot6HitsNoShowerSeg( dbe_->book1D("CSCSegments6HitsNoShowerSeg", "", 100, 0, 10) ),
25  theTypePlot6HitsShower( dbe_->book1D("CSCSegments6HitsShower", "", 100, 0, 10) ),
26  theTypePlot6HitsShowerSeg( dbe_->book1D("CSCSegments6HitsShowerSeg", "", 100, 0, 10) )
27 {
29 
30  dbe_->setCurrentFolder("CSCRecHitsV/CSCRecHitTask");
31  for(int i = 0; i < 10; ++i)
32  {
33  char title1[200], title2[200], title3[200], title4[200],
34  title5[200], title6[200], title7[200], title8[200];
35  sprintf(title1, "CSCSegmentRdPhiResolution%d", i+1);
36  sprintf(title2, "CSCSegmentRdPhiPull%d", i+1);
37  sprintf(title3, "CSCSegmentThetaResolution%d", i+1);
38  sprintf(title4, "CSCSegmentThetaPull%d", i+1);
39  sprintf(title5, "CSCSegmentdXdZResolution%d", i+1);
40  sprintf(title6, "CSCSegmentdXdZPull%d", i+1);
41  sprintf(title7, "CSCSegmentdYdZResolution%d", i+1);
42  sprintf(title8, "CSCSegmentdYdZPull%d", i+1);
43 
44 
45  theRdPhiResolutionPlots[i] = dbe_->book1D(title1, title1, 100, -0.4, 0.4);
46  theRdPhiPullPlots[i] = dbe_->book1D(title2, title2, 100, -5, 5);
47  theThetaResolutionPlots[i] = dbe_->book1D(title3, title3, 100, -1, 1);
48  theThetaPullPlots[i] = dbe_->book1D(title4, title4, 100, -5, 5);
49  thedXdZResolutionPlots[i] = dbe_->book1D(title5, title5, 100, -1, 1);
50  thedXdZPullPlots[i] = dbe_->book1D(title6, title6, 100, -5, 5);
51  thedYdZResolutionPlots[i] = dbe_->book1D(title7, title7, 100, -1, 1);
52  thedYdZPullPlots[i] = dbe_->book1D(title8, title8, 100, -5, 5);
53 
54  }
55 }
MonitorElement * theNPerEventPlot
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int i
Definition: DBlmapReader.cc:9
MonitorElement * theTypePlot5HitsNoShowerSeg
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:928
MonitorElement * theTypePlot4HitsShowerSeg
MonitorElement * theNPerChamberTypePlot
MonitorElement * theTypePlot6HitsShower
MonitorElement * theTypePlot5HitsNoShower
edm::EDGetTokenT< CSCSegmentCollection > segments_Token_
MonitorElement * theNRecHitsPlot
MonitorElement * theTypePlot4HitsShower
MonitorElement * theRdPhiResolutionPlots[10]
MonitorElement * theTypePlot6HitsShowerSeg
MonitorElement * theTypePlot4HitsNoShower
MonitorElement * theTypePlot6HitsNoShowerSeg
MonitorElement * theTypePlot6HitsNoShower
ChamberHitMap theLayerHitsPerChamber
MonitorElement * thedYdZResolutionPlots[10]
CSCBaseValidation(DQMStore *dbe, const edm::InputTag &inputTag)
MonitorElement * theRdPhiPullPlots[10]
MonitorElement * theThetaPullPlots[10]
MonitorElement * thedXdZPullPlots[10]
MonitorElement * theTypePlot5HitsShowerSeg
MonitorElement * theThetaResolutionPlots[10]
ChamberSegmentMap theChamberSegmentMap
MonitorElement * theTypePlot4HitsNoShowerSeg
MonitorElement * theTypePlot5HitsShower
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:639
MonitorElement * thedXdZResolutionPlots[10]
MonitorElement * thedYdZPullPlots[10]
virtual CSCSegmentValidation::~CSCSegmentValidation ( )
inlinevirtual

Definition at line 17 of file CSCSegmentValidation.h.

17 {}

Member Function Documentation

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

Implements CSCBaseValidation.

Definition at line 57 of file CSCSegmentValidation.cc.

References PSimHit::detUnitId(), MonitorElement::Fill(), fillEfficiencyPlots(), fillLayerHitsPerChamber(), CSCBaseValidation::findLayer(), edm::Event::getByToken(), keyHit(), plotResolution(), edm::Handle< T >::product(), segments_Token_, theChamberSegmentMap, theNPerChamberTypePlot, theNPerEventPlot, theNRecHitsPlot, and whatChamberType().

Referenced by CSCRecHitValidation::analyze().

58 {
59  // get the collection of CSCRecHsegmentItrD
61  e.getByToken(segments_Token_, hRecHits);
62  const CSCSegmentCollection * cscRecHits = hRecHits.product();
63 
64  theChamberSegmentMap.clear();
65  unsigned nPerEvent = 0;
66  for(CSCSegmentCollection::const_iterator segmentItr = cscRecHits->begin();
67  segmentItr != cscRecHits->end(); segmentItr++)
68  {
69  ++nPerEvent;
70  int detId = segmentItr->geographicalId().rawId();
71  int chamberType = whatChamberType(detId);
72 
73  theNRecHitsPlot->Fill(segmentItr->nRecHits());
74  theNPerChamberTypePlot->Fill(chamberType);
75  theChamberSegmentMap[detId].push_back(*segmentItr);
76 
77  // do the resolution plots
78  const PSimHit * hit = keyHit(detId);
79  if(hit != 0)
80  {
81  const CSCLayer * layer = findLayer(hit->detUnitId());
82  plotResolution(*hit, *segmentItr, layer, chamberType);
83  }
84  }
85 
86  theNPerEventPlot->Fill(nPerEvent);
87 
90 }
MonitorElement * theNPerEventPlot
void plotResolution(const PSimHit &simHit, const CSCSegment &recHit, const CSCLayer *layer, int chamberType)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
MonitorElement * theNPerChamberTypePlot
static int whatChamberType(int detId)
edm::EDGetTokenT< CSCSegmentCollection > segments_Token_
MonitorElement * theNRecHitsPlot
void Fill(long long x)
const PSimHit * keyHit(int chamberId) const
T const * product() const
Definition: Handle.h:81
ChamberSegmentMap theChamberSegmentMap
const CSCLayer * findLayer(int detId) const
unsigned int detUnitId() const
Definition: PSimHit.h:93
void CSCSegmentValidation::fillEfficiencyPlots ( )
private

Definition at line 93 of file CSCSegmentValidation.cc.

References MonitorElement::Fill(), hasSegment(), python.multivaluedict::sort(), theLayerHitsPerChamber, theShowerThreshold, theTypePlot4HitsNoShower, theTypePlot4HitsNoShowerSeg, theTypePlot4HitsShower, theTypePlot4HitsShowerSeg, theTypePlot5HitsNoShower, theTypePlot5HitsNoShowerSeg, theTypePlot5HitsShower, theTypePlot5HitsShowerSeg, theTypePlot6HitsNoShower, theTypePlot6HitsNoShowerSeg, theTypePlot6HitsShower, theTypePlot6HitsShowerSeg, findQualityFiles::v, and whatChamberType().

Referenced by analyze().

94 {
95  // now plot efficiency by looping over all chambers with hits
96  for(ChamberHitMap::const_iterator mapItr = theLayerHitsPerChamber.begin(),
97  mapEnd = theLayerHitsPerChamber.end();
98  mapItr != mapEnd;
99  ++mapItr)
100  {
101  int chamberId = mapItr->first;
102  int nHitsInChamber = mapItr->second.size();
103  bool isShower = (nHitsInChamber > theShowerThreshold);
104  bool hasSeg = hasSegment(chamberId);
105  int chamberType = whatChamberType(chamberId);
106  // find how many layers were hit in this chamber
107  std::vector<int> v = mapItr->second;
108  std::sort(v.begin(), v.end());
109  // maybe can just count
110  v.erase(std::unique(v.begin(), v.end()), v.end());
111  int nLayersHit = v.size();
112 
113  if(nLayersHit == 4)
114  {
115 
116  if(isShower) theTypePlot4HitsShower->Fill(chamberType);
117  else theTypePlot4HitsNoShower->Fill(chamberType);
118 
119  if(hasSeg)
120  {
121  if(isShower) theTypePlot4HitsShowerSeg->Fill(chamberType);
122  else theTypePlot4HitsNoShowerSeg->Fill(chamberType);
123  }
124  }
125 
126  if(nLayersHit == 5)
127  {
128 
129  if(isShower) theTypePlot5HitsShower->Fill(chamberType);
130  else theTypePlot5HitsNoShower->Fill(chamberType);
131 
132  if(hasSeg)
133  {
134  if(isShower) theTypePlot5HitsShowerSeg->Fill(chamberType);
135  else theTypePlot5HitsNoShowerSeg->Fill(chamberType);
136  }
137  }
138 
139  if(nLayersHit == 6)
140  {
141 
142  if(isShower) theTypePlot6HitsShower->Fill(chamberType);
143  else theTypePlot6HitsNoShower->Fill(chamberType);
144 
145  if(hasSeg)
146  {
147  if(isShower) theTypePlot6HitsShowerSeg->Fill(chamberType);
148  else theTypePlot6HitsNoShowerSeg->Fill(chamberType);
149  }
150  }
151 
152 
153  }
154 }
MonitorElement * theTypePlot5HitsNoShowerSeg
MonitorElement * theTypePlot4HitsShowerSeg
static int whatChamberType(int detId)
MonitorElement * theTypePlot6HitsShower
MonitorElement * theTypePlot5HitsNoShower
MonitorElement * theTypePlot4HitsShower
MonitorElement * theTypePlot6HitsShowerSeg
void Fill(long long x)
MonitorElement * theTypePlot4HitsNoShower
MonitorElement * theTypePlot6HitsNoShowerSeg
MonitorElement * theTypePlot6HitsNoShower
ChamberHitMap theLayerHitsPerChamber
bool hasSegment(int chamberId) const
MonitorElement * theTypePlot5HitsShowerSeg
MonitorElement * theTypePlot4HitsNoShowerSeg
MonitorElement * theTypePlot5HitsShower
void CSCSegmentValidation::fillLayerHitsPerChamber ( )
private

Definition at line 201 of file CSCSegmentValidation.cc.

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

Referenced by analyze().

202 {
203  theLayerHitsPerChamber.clear();
204  std::vector<int> layersHit = theSimHitMap->detsWithHits();
205  for(std::vector<int>::const_iterator layerItr = layersHit.begin(),
206  layersHitEnd = layersHit.end();
207  layerItr != layersHitEnd;
208  ++layerItr)
209  {
210  CSCDetId layerId(*layerItr);
211  CSCDetId chamberId = layerId.chamberId();
212  int nhits = theSimHitMap->hits(*layerItr).size();
213  // multiple entries, so we can see showers
214  for(int i = 0; i < nhits; ++i) {
215  theLayerHitsPerChamber[chamberId.rawId()].push_back(*layerItr);
216  }
217  }
218 
219 }
int i
Definition: DBlmapReader.cc:9
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
const PSimHitMap * theSimHitMap
const edm::PSimHitContainer & hits(int detId) const
Definition: PSimHitMap.cc:22
CSCDetId chamberId() const
Definition: CSCDetId.h:66
ChamberHitMap theLayerHitsPerChamber
std::vector< int > detsWithHits() const
Definition: PSimHitMap.cc:37
bool CSCSegmentValidation::hasSegment ( int  chamberId) const
private

Definition at line 156 of file CSCSegmentValidation.cc.

References theChamberSegmentMap.

Referenced by fillEfficiencyPlots().

157 {
158  return (theChamberSegmentMap.find(chamberId) != theChamberSegmentMap.end());
159 }
ChamberSegmentMap theChamberSegmentMap
const PSimHit * CSCSegmentValidation::keyHit ( int  chamberId) const
private

Definition at line 229 of file CSCSegmentValidation.cc.

References PSimHitMap::hits(), query::result, CSCSegmentValidationUtils::SimHitPabsLessThan(), and CSCBaseValidation::theSimHitMap.

Referenced by analyze().

230 {
231  const PSimHit * result = 0;
232  int layerId = chamberId + 3;
233  const edm::PSimHitContainer & layerHits = theSimHitMap->hits(layerId);
234 
235  if(!layerHits.empty())
236  {
237  // pick the hit with maximum energy
238  edm::PSimHitContainer::const_iterator hitItr = std::max_element(layerHits.begin(), layerHits.end(),
240  result = &(*hitItr);
241  }
242  return result;
243 }
const PSimHitMap * theSimHitMap
tuple result
Definition: query.py:137
const edm::PSimHitContainer & hits(int detId) const
Definition: PSimHitMap.cc:22
std::vector< PSimHit > PSimHitContainer
bool SimHitPabsLessThan(const PSimHit &p1, const PSimHit &p2)
void CSCSegmentValidation::plotResolution ( const PSimHit simHit,
const CSCSegment recHit,
const CSCLayer layer,
int  chamberType 
)
private

Definition at line 169 of file CSCSegmentValidation.cc.

References MonitorElement::Fill(), PSimHit::localDirection(), PSimHit::localPosition(), PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), mathSSE::sqrt(), thedXdZPullPlots, thedXdZResolutionPlots, thedYdZPullPlots, thedYdZResolutionPlots, theRdPhiPullPlots, theRdPhiResolutionPlots, PV3DBase< T, PVType, FrameType >::theta(), theThetaResolutionPlots, GeomDet::toGlobal(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by analyze().

171 {
172  GlobalPoint simHitPos = layer->toGlobal(simHit.localPosition());
173  GlobalPoint segmentPos = layer->toGlobal(segment.localPosition());
174  LocalVector simHitDir = simHit.localDirection();
175  LocalVector segmentDir = segment.localDirection();
176 
177  double dphi = segmentPos.phi() - simHitPos.phi();
178  double rdphi = segmentPos.perp() * dphi;
179  double dtheta = segmentPos.theta() - simHitPos.theta();
180 
181  double sigmax = sqrt(segment.localPositionError().xx());
182  //double sigmay = sqrt(segment.localPositionError().yy());
183 
184  double ddxdz = segmentDir.x()/segmentDir.z() - simHitDir.x()/simHitDir.z();
185  double ddydz = segmentDir.y()/segmentDir.z() - simHitDir.y()/simHitDir.z();
186  double sigmadxdz = sqrt(segment.localDirectionError().xx());
187  double sigmadydz = sqrt(segment.localDirectionError().yy());
188 
189  theRdPhiResolutionPlots[chamberType-1]->Fill( rdphi );
190  theRdPhiPullPlots[chamberType-1]->Fill( rdphi/sigmax );
191  theThetaResolutionPlots[chamberType-1]->Fill( dtheta );
192  //theThetaPullPlots[chamberType-1]->Fill( dy/sigmay );
193 
194  thedXdZResolutionPlots[chamberType-1]->Fill( ddxdz );
195  thedXdZPullPlots[chamberType-1]->Fill( ddxdz/sigmadxdz );
196  thedYdZResolutionPlots[chamberType-1]->Fill( ddydz );
197  thedYdZPullPlots[chamberType-1]->Fill( ddydz/sigmadydz );
198 }
T perp() const
Definition: PV3DBase.h:72
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:47
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
MonitorElement * theRdPhiResolutionPlots[10]
void Fill(long long x)
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
Local3DPoint localPosition() const
Definition: PSimHit.h:44
T sqrt(T t)
Definition: SSEVec.h:48
T z() const
Definition: PV3DBase.h:64
MonitorElement * thedYdZResolutionPlots[10]
MonitorElement * theRdPhiPullPlots[10]
LocalVector localDirection() const
Obsolete. Same as momentumAtEntry().unit(), for backward compatibility.
Definition: PSimHit.h:52
MonitorElement * thedXdZPullPlots[10]
MonitorElement * theThetaResolutionPlots[10]
T x() const
Definition: PV3DBase.h:62
MonitorElement * thedXdZResolutionPlots[10]
MonitorElement * thedYdZPullPlots[10]
int CSCSegmentValidation::whatChamberType ( int  detId)
staticprivate

Definition at line 162 of file CSCSegmentValidation.cc.

References CSCDetId::ring(), CSCDetId::station(), and CSCChamberSpecs::whatChamberType().

Referenced by analyze(), and fillEfficiencyPlots().

163 {
164  CSCDetId cscDetId(detId);
165  return CSCChamberSpecs::whatChamberType(cscDetId.station(), cscDetId.ring());
166 }
static int whatChamberType(int istation, int iring)

Member Data Documentation

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

Definition at line 27 of file CSCSegmentValidation.h.

Referenced by analyze(), and CSCSegmentValidation().

ChamberSegmentMap CSCSegmentValidation::theChamberSegmentMap
private

Definition at line 39 of file CSCSegmentValidation.h.

Referenced by analyze(), and hasSegment().

MonitorElement* CSCSegmentValidation::thedXdZPullPlots[10]
private

Definition at line 51 of file CSCSegmentValidation.h.

Referenced by CSCSegmentValidation(), and plotResolution().

MonitorElement* CSCSegmentValidation::thedXdZResolutionPlots[10]
private

Definition at line 50 of file CSCSegmentValidation.h.

Referenced by CSCSegmentValidation(), and plotResolution().

MonitorElement* CSCSegmentValidation::thedYdZPullPlots[10]
private

Definition at line 53 of file CSCSegmentValidation.h.

Referenced by CSCSegmentValidation(), and plotResolution().

MonitorElement* CSCSegmentValidation::thedYdZResolutionPlots[10]
private

Definition at line 52 of file CSCSegmentValidation.h.

Referenced by CSCSegmentValidation(), and plotResolution().

ChamberHitMap CSCSegmentValidation::theLayerHitsPerChamber
private

Definition at line 32 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots(), and fillLayerHitsPerChamber().

MonitorElement* CSCSegmentValidation::theNPerChamberTypePlot
private

Definition at line 45 of file CSCSegmentValidation.h.

Referenced by analyze().

MonitorElement* CSCSegmentValidation::theNPerEventPlot
private

Definition at line 43 of file CSCSegmentValidation.h.

Referenced by analyze().

MonitorElement* CSCSegmentValidation::theNRecHitsPlot
private

Definition at line 44 of file CSCSegmentValidation.h.

Referenced by analyze().

MonitorElement* CSCSegmentValidation::theRdPhiPullPlots[10]
private

Definition at line 47 of file CSCSegmentValidation.h.

Referenced by CSCSegmentValidation(), and plotResolution().

MonitorElement* CSCSegmentValidation::theRdPhiResolutionPlots[10]
private

Definition at line 46 of file CSCSegmentValidation.h.

Referenced by CSCSegmentValidation(), and plotResolution().

int CSCSegmentValidation::theShowerThreshold
private

Definition at line 41 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theThetaPullPlots[10]
private

Definition at line 49 of file CSCSegmentValidation.h.

Referenced by CSCSegmentValidation().

MonitorElement* CSCSegmentValidation::theThetaResolutionPlots[10]
private

Definition at line 48 of file CSCSegmentValidation.h.

Referenced by CSCSegmentValidation(), and plotResolution().

MonitorElement* CSCSegmentValidation::theTypePlot4HitsNoShower
private

Definition at line 57 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot4HitsNoShowerSeg
private

Definition at line 58 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot4HitsShower
private

Definition at line 59 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot4HitsShowerSeg
private

Definition at line 60 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot5HitsNoShower
private

Definition at line 61 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot5HitsNoShowerSeg
private

Definition at line 62 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot5HitsShower
private

Definition at line 63 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot5HitsShowerSeg
private

Definition at line 64 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot6HitsNoShower
private

Definition at line 65 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot6HitsNoShowerSeg
private

Definition at line 66 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot6HitsShower
private

Definition at line 67 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot6HitsShowerSeg
private

Definition at line 68 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().