CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CSCSegmentValidation.cc
Go to the documentation of this file.
3 #include <algorithm>
5 
6 
7 
9 : CSCBaseValidation(dbe, inputTag),
10  theLayerHitsPerChamber(),
11  theChamberSegmentMap(),
12  theShowerThreshold(10),
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 }
55 
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 }
90 
91 
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 }
154 
155 bool CSCSegmentValidation::hasSegment(int chamberId) const
156 {
157  return (theChamberSegmentMap.find(chamberId) != theChamberSegmentMap.end());
158 }
159 
160 
162 {
163  CSCDetId cscDetId(detId);
164  return CSCChamberSpecs::whatChamberType(cscDetId.station(), cscDetId.ring());
165 }
166 
167 
168 void CSCSegmentValidation::plotResolution(const PSimHit & simHit, const CSCSegment & segment,
169  const CSCLayer * layer, int chamberType)
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 }
198 
199 
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 }
219 
220 namespace CSCSegmentValidationUtils {
221  bool SimHitPabsLessThan(const PSimHit & p1, const PSimHit & p2)
222  {
223  return p1.pabs() < p2.pabs();
224  }
225 }
226 
227 
228 const PSimHit * CSCSegmentValidation::keyHit(int chamberId) const
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 }
243 
MonitorElement * theNPerEventPlot
int i
Definition: DBlmapReader.cc:9
edm::InputTag theInputTag
void plotResolution(const PSimHit &simHit, const CSCSegment &recHit, const CSCLayer *layer, int chamberType)
MonitorElement * theTypePlot5HitsNoShowerSeg
T perp() const
Definition: PV3DBase.h:72
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:873
MonitorElement * theTypePlot4HitsShowerSeg
MonitorElement * theNPerChamberTypePlot
static int whatChamberType(int istation, int iring)
static int whatChamberType(int detId)
MonitorElement * theTypePlot6HitsShower
CSCSegmentValidation(DQMStore *dbe, const edm::InputTag &inputTag)
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 * theTypePlot5HitsNoShower
MonitorElement * theNRecHitsPlot
MonitorElement * theTypePlot4HitsShower
MonitorElement * theRdPhiResolutionPlots[10]
MonitorElement * theTypePlot6HitsShowerSeg
void Fill(long long x)
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
MonitorElement * theTypePlot4HitsNoShower
Local3DPoint localPosition() const
Definition: PSimHit.h:44
const PSimHitMap * theSimHitMap
T sqrt(T t)
Definition: SSEVec.h:48
MonitorElement * theTypePlot6HitsNoShowerSeg
T z() const
Definition: PV3DBase.h:64
tuple result
Definition: query.py:137
const edm::PSimHitContainer & hits(int detId) const
Definition: PSimHitMap.cc:32
CSCDetId chamberId() const
Definition: CSCDetId.h:55
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
Definition: PSimHit.h:63
MonitorElement * theTypePlot6HitsNoShower
double p2[4]
Definition: TauolaWrapper.h:90
ChamberHitMap theLayerHitsPerChamber
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
DQMStore * dbe_
MonitorElement * thedYdZResolutionPlots[10]
std::vector< int > detsWithHits() const
Definition: PSimHitMap.cc:47
MonitorElement * theRdPhiPullPlots[10]
LocalVector localDirection() const
Obsolete. Same as momentumAtEntry().unit(), for backward compatibility.
Definition: PSimHit.h:52
int ring() const
Definition: CSCDetId.h:77
virtual void analyze(const edm::Event &, const edm::EventSetup &)
MonitorElement * theThetaPullPlots[10]
const PSimHit * keyHit(int chamberId) const
MonitorElement * thedXdZPullPlots[10]
T const * product() const
Definition: Handle.h:81
double p1[4]
Definition: TauolaWrapper.h:89
bool hasSegment(int chamberId) const
MonitorElement * theTypePlot5HitsShowerSeg
int station() const
Definition: CSCDetId.h:88
MonitorElement * theThetaResolutionPlots[10]
ChamberSegmentMap theChamberSegmentMap
std::vector< PSimHit > PSimHitContainer
MonitorElement * theTypePlot4HitsNoShowerSeg
bool SimHitPabsLessThan(const PSimHit &p1, const PSimHit &p2)
T x() const
Definition: PV3DBase.h:62
MonitorElement * theTypePlot5HitsShower
const CSCLayer * findLayer(int detId) const
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:585
unsigned int detUnitId() const
Definition: PSimHit.h:93
MonitorElement * thedXdZResolutionPlots[10]
MonitorElement * thedYdZPullPlots[10]