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)
 
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

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 10 of file CSCSegmentValidation.h.

Member Typedef Documentation

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

Definition at line 27 of file CSCSegmentValidation.h.

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

Definition at line 34 of file CSCSegmentValidation.h.

Constructor & Destructor Documentation

CSCSegmentValidation::CSCSegmentValidation ( DQMStore dbe,
const edm::InputTag inputTag 
)

Definition at line 8 of file CSCSegmentValidation.cc.

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

9 : CSCBaseValidation(dbe, inputTag),
13  theNPerEventPlot( dbe_->book1D("CSCSegmentsPerEvent", "Number of CSC segments per event", 100, 0, 50) ),
14  theNRecHitsPlot( dbe_->book1D("CSCRecHitsPerSegment", "Number of CSC rec hits per segment" , 8, 0, 7) ),
15  theNPerChamberTypePlot( dbe_->book1D("CSCSegmentsPerChamberType", "Number of CSC segments per chamber type", 11, 0, 10) ),
16  theTypePlot4HitsNoShower( dbe_->book1D("CSCSegments4HitsNoShower", "", 100, 0, 10) ),
17  theTypePlot4HitsNoShowerSeg( dbe_->book1D("CSCSegments4HitsNoShowerSeg", "", 100, 0, 10) ),
18  theTypePlot4HitsShower( dbe_->book1D("CSCSegments4HitsShower", "", 100, 0, 10) ),
19  theTypePlot4HitsShowerSeg( dbe_->book1D("CSCSegments4HitsShowerSeg", "", 100, 0, 10) ),
20  theTypePlot5HitsNoShower( dbe_->book1D("CSCSegments5HitsNoShower", "", 100, 0, 10) ),
21  theTypePlot5HitsNoShowerSeg( dbe_->book1D("CSCSegments5HitsNoShowerSeg", "", 100, 0, 10) ),
22  theTypePlot5HitsShower( dbe_->book1D("CSCSegments5HitsShower", "", 100, 0, 10) ),
23  theTypePlot5HitsShowerSeg( dbe_->book1D("CSCSegments5HitsShowerSeg", "", 100, 0, 10) ),
24  theTypePlot6HitsNoShower( dbe_->book1D("CSCSegments6HitsNoShower", "", 100, 0, 10) ),
25  theTypePlot6HitsNoShowerSeg( dbe_->book1D("CSCSegments6HitsNoShowerSeg", "", 100, 0, 10) ),
26  theTypePlot6HitsShower( dbe_->book1D("CSCSegments6HitsShower", "", 100, 0, 10) ),
27  theTypePlot6HitsShowerSeg( dbe_->book1D("CSCSegments6HitsShowerSeg", "", 100, 0, 10) )
28 {
29  dbe_->setCurrentFolder("CSCRecHitsV/CSCRecHitTask");
30  for(int i = 0; i < 10; ++i)
31  {
32  char title1[200], title2[200], title3[200], title4[200],
33  title5[200], title6[200], title7[200], title8[200];
34  sprintf(title1, "CSCSegmentRdPhiResolution%d", i+1);
35  sprintf(title2, "CSCSegmentRdPhiPull%d", i+1);
36  sprintf(title3, "CSCSegmentThetaResolution%d", i+1);
37  sprintf(title4, "CSCSegmentThetaPull%d", i+1);
38  sprintf(title5, "CSCSegmentdXdZResolution%d", i+1);
39  sprintf(title6, "CSCSegmentdXdZPull%d", i+1);
40  sprintf(title7, "CSCSegmentdYdZResolution%d", i+1);
41  sprintf(title8, "CSCSegmentdYdZPull%d", i+1);
42 
43 
44  theRdPhiResolutionPlots[i] = dbe_->book1D(title1, title1, 100, -0.4, 0.4);
45  theRdPhiPullPlots[i] = dbe_->book1D(title2, title2, 100, -5, 5);
46  theThetaResolutionPlots[i] = dbe_->book1D(title3, title3, 100, -1, 1);
47  theThetaPullPlots[i] = dbe_->book1D(title4, title4, 100, -5, 5);
48  thedXdZResolutionPlots[i] = dbe_->book1D(title5, title5, 100, -1, 1);
49  thedXdZPullPlots[i] = dbe_->book1D(title6, title6, 100, -5, 5);
50  thedYdZResolutionPlots[i] = dbe_->book1D(title7, title7, 100, -1, 1);
51  thedYdZPullPlots[i] = dbe_->book1D(title8, title8, 100, -5, 5);
52 
53  }
54 }
MonitorElement * theNPerEventPlot
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:717
MonitorElement * theTypePlot4HitsShowerSeg
MonitorElement * theNPerChamberTypePlot
MonitorElement * theTypePlot6HitsShower
MonitorElement * theTypePlot5HitsNoShower
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:429
MonitorElement * thedXdZResolutionPlots[10]
MonitorElement * thedYdZPullPlots[10]
virtual CSCSegmentValidation::~CSCSegmentValidation ( )
inlinevirtual

Definition at line 15 of file CSCSegmentValidation.h.

15 {}

Member Function Documentation

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

Implements CSCBaseValidation.

Definition at line 56 of file CSCSegmentValidation.cc.

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

Referenced by CSCRecHitValidation::analyze().

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

Definition at line 92 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, v, and whatChamberType().

Referenced by analyze().

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

Definition at line 200 of file CSCSegmentValidation.cc.

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

Referenced by analyze().

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

Definition at line 155 of file CSCSegmentValidation.cc.

References theChamberSegmentMap.

Referenced by fillEfficiencyPlots().

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

Definition at line 228 of file CSCSegmentValidation.cc.

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

Referenced by analyze().

229 {
230  const PSimHit * result = 0;
231  int layerId = chamberId + 3;
232  const edm::PSimHitContainer & layerHits = theSimHitMap->hits(layerId);
233 
234  if(!layerHits.empty())
235  {
236  // pick the hit with maximum energy
237  edm::PSimHitContainer::const_iterator hitItr = std::max_element(layerHits.begin(), layerHits.end(),
239  result = &(*hitItr);
240  }
241  return result;
242 }
const PSimHitMap * theSimHitMap
tuple result
Definition: query.py:137
const edm::PSimHitContainer & hits(int detId) const
Definition: PSimHitMap.cc:32
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 168 of file CSCSegmentValidation.cc.

References MonitorElement::Fill(), CSCSegment::localDirection(), PSimHit::localDirection(), CSCSegment::localDirectionError(), CSCSegment::localPosition(), PSimHit::localPosition(), CSCSegment::localPositionError(), 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(), LocalError::xx(), PV3DBase< T, PVType, FrameType >::y(), LocalError::yy(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by analyze().

170 {
171  GlobalPoint simHitPos = layer->toGlobal(simHit.localPosition());
172  GlobalPoint segmentPos = layer->toGlobal(segment.localPosition());
173  LocalVector simHitDir = simHit.localDirection();
174  LocalVector segmentDir = segment.localDirection();
175 
176  double dphi = segmentPos.phi() - simHitPos.phi();
177  double rdphi = segmentPos.perp() * dphi;
178  double dtheta = segmentPos.theta() - simHitPos.theta();
179 
180  double sigmax = sqrt(segment.localPositionError().xx());
181  //double sigmay = sqrt(segment.localPositionError().yy());
182 
183  double ddxdz = segmentDir.x()/segmentDir.z() - simHitDir.x()/simHitDir.z();
184  double ddydz = segmentDir.y()/segmentDir.z() - simHitDir.y()/simHitDir.z();
185  double sigmadxdz = sqrt(segment.localDirectionError().xx());
186  double sigmadydz = sqrt(segment.localDirectionError().yy());
187 
188  theRdPhiResolutionPlots[chamberType-1]->Fill( rdphi );
189  theRdPhiPullPlots[chamberType-1]->Fill( rdphi/sigmax );
190  theThetaResolutionPlots[chamberType-1]->Fill( dtheta );
191  //theThetaPullPlots[chamberType-1]->Fill( dy/sigmay );
192 
193  thedXdZResolutionPlots[chamberType-1]->Fill( ddxdz );
194  thedXdZPullPlots[chamberType-1]->Fill( ddxdz/sigmadxdz );
195  thedYdZResolutionPlots[chamberType-1]->Fill( ddydz );
196  thedYdZPullPlots[chamberType-1]->Fill( ddydz/sigmadydz );
197 }
T perp() const
Definition: PV3DBase.h:71
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:68
T y() const
Definition: PV3DBase.h:62
MonitorElement * theRdPhiResolutionPlots[10]
void Fill(long long x)
Geom::Theta< T > theta() const
Definition: PV3DBase.h:74
Local3DPoint localPosition() const
Definition: PSimHit.h:44
T sqrt(T t)
Definition: SSEVec.h:46
T z() const
Definition: PV3DBase.h:63
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:61
MonitorElement * thedXdZResolutionPlots[10]
MonitorElement * thedYdZPullPlots[10]
int CSCSegmentValidation::whatChamberType ( int  detId)
staticprivate

Definition at line 161 of file CSCSegmentValidation.cc.

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

Referenced by analyze(), and fillEfficiencyPlots().

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

Member Data Documentation

ChamberSegmentMap CSCSegmentValidation::theChamberSegmentMap
private

Definition at line 35 of file CSCSegmentValidation.h.

Referenced by analyze(), and hasSegment().

MonitorElement* CSCSegmentValidation::thedXdZPullPlots[10]
private

Definition at line 47 of file CSCSegmentValidation.h.

Referenced by CSCSegmentValidation(), and plotResolution().

MonitorElement* CSCSegmentValidation::thedXdZResolutionPlots[10]
private

Definition at line 46 of file CSCSegmentValidation.h.

Referenced by CSCSegmentValidation(), and plotResolution().

MonitorElement* CSCSegmentValidation::thedYdZPullPlots[10]
private

Definition at line 49 of file CSCSegmentValidation.h.

Referenced by CSCSegmentValidation(), and plotResolution().

MonitorElement* CSCSegmentValidation::thedYdZResolutionPlots[10]
private

Definition at line 48 of file CSCSegmentValidation.h.

Referenced by CSCSegmentValidation(), and plotResolution().

ChamberHitMap CSCSegmentValidation::theLayerHitsPerChamber
private

Definition at line 28 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots(), and fillLayerHitsPerChamber().

MonitorElement* CSCSegmentValidation::theNPerChamberTypePlot
private

Definition at line 41 of file CSCSegmentValidation.h.

Referenced by analyze().

MonitorElement* CSCSegmentValidation::theNPerEventPlot
private

Definition at line 39 of file CSCSegmentValidation.h.

Referenced by analyze().

MonitorElement* CSCSegmentValidation::theNRecHitsPlot
private

Definition at line 40 of file CSCSegmentValidation.h.

Referenced by analyze().

MonitorElement* CSCSegmentValidation::theRdPhiPullPlots[10]
private

Definition at line 43 of file CSCSegmentValidation.h.

Referenced by CSCSegmentValidation(), and plotResolution().

MonitorElement* CSCSegmentValidation::theRdPhiResolutionPlots[10]
private

Definition at line 42 of file CSCSegmentValidation.h.

Referenced by CSCSegmentValidation(), and plotResolution().

int CSCSegmentValidation::theShowerThreshold
private

Definition at line 37 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theThetaPullPlots[10]
private

Definition at line 45 of file CSCSegmentValidation.h.

Referenced by CSCSegmentValidation().

MonitorElement* CSCSegmentValidation::theThetaResolutionPlots[10]
private

Definition at line 44 of file CSCSegmentValidation.h.

Referenced by CSCSegmentValidation(), and plotResolution().

MonitorElement* CSCSegmentValidation::theTypePlot4HitsNoShower
private

Definition at line 53 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot4HitsNoShowerSeg
private

Definition at line 54 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot4HitsShower
private

Definition at line 55 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot4HitsShowerSeg
private

Definition at line 56 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot5HitsNoShower
private

Definition at line 57 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot5HitsNoShowerSeg
private

Definition at line 58 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot5HitsShower
private

Definition at line 59 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot5HitsShowerSeg
private

Definition at line 60 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot6HitsNoShower
private

Definition at line 61 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot6HitsNoShowerSeg
private

Definition at line 62 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot6HitsShower
private

Definition at line 63 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().

MonitorElement* CSCSegmentValidation::theTypePlot6HitsShowerSeg
private

Definition at line 64 of file CSCSegmentValidation.h.

Referenced by fillEfficiencyPlots().