CMS 3D CMS Logo

CSCSegmentValidation.cc
Go to the documentation of this file.
2 #include <algorithm>
4 
5 
7 : CSCBaseValidation(inputTag),
8  theLayerHitsPerChamber(),
9  theChamberSegmentMap(),
10  theShowerThreshold(10)
11 {
13 }
14 
16 {
17 }
18 
20 {
21  theNPerEventPlot = iBooker.book1D("CSCSegmentsPerEvent", "Number of CSC segments per event", 100, 0, 50);
22  theNRecHitsPlot = iBooker.book1D("CSCRecHitsPerSegment", "Number of CSC rec hits per segment" , 8, 0, 7);
23  theNPerChamberTypePlot = iBooker.book1D("CSCSegmentsPerChamberType", "Number of CSC segments per chamber type", 11, 0, 10);
24  theTypePlot4HitsNoShower = iBooker.book1D("CSCSegments4HitsNoShower", "", 100, 0, 10);
25  theTypePlot4HitsNoShowerSeg = iBooker.book1D("CSCSegments4HitsNoShowerSeg", "", 100, 0, 10);
26  theTypePlot4HitsShower = iBooker.book1D("CSCSegments4HitsShower", "", 100, 0, 10);
27  theTypePlot4HitsShowerSeg = iBooker.book1D("CSCSegments4HitsShowerSeg", "", 100, 0, 10);
28  theTypePlot5HitsNoShower = iBooker.book1D("CSCSegments5HitsNoShower", "", 100, 0, 10);
29  theTypePlot5HitsNoShowerSeg = iBooker.book1D("CSCSegments5HitsNoShowerSeg", "", 100, 0, 10);
30  theTypePlot5HitsShower = iBooker.book1D("CSCSegments5HitsShower", "", 100, 0, 10);
31  theTypePlot5HitsShowerSeg = iBooker.book1D("CSCSegments5HitsShowerSeg", "", 100, 0, 10);
32  theTypePlot6HitsNoShower = iBooker.book1D("CSCSegments6HitsNoShower", "", 100, 0, 10);
33  theTypePlot6HitsNoShowerSeg = iBooker.book1D("CSCSegments6HitsNoShowerSeg", "", 100, 0, 10);
34  theTypePlot6HitsShower = iBooker.book1D("CSCSegments6HitsShower", "", 100, 0, 10);
35  theTypePlot6HitsShowerSeg = iBooker.book1D("CSCSegments6HitsShowerSeg", "", 100, 0, 10);
36  for(int i = 0; i < 10; ++i)
37  {
38  char title1[200], title2[200], title3[200], title4[200],
39  title5[200], title6[200], title7[200], title8[200];
40  sprintf(title1, "CSCSegmentRdPhiResolution%d", i+1);
41  sprintf(title2, "CSCSegmentRdPhiPull%d", i+1);
42  sprintf(title3, "CSCSegmentThetaResolution%d", i+1);
43  sprintf(title4, "CSCSegmentThetaPull%d", i+1);
44  sprintf(title5, "CSCSegmentdXdZResolution%d", i+1);
45  sprintf(title6, "CSCSegmentdXdZPull%d", i+1);
46  sprintf(title7, "CSCSegmentdYdZResolution%d", i+1);
47  sprintf(title8, "CSCSegmentdYdZPull%d", i+1);
48 
49  theRdPhiResolutionPlots[i] = iBooker.book1D(title1, title1, 100, -0.4, 0.4);
50  theRdPhiPullPlots[i] = iBooker.book1D(title2, title2, 100, -5, 5);
51  theThetaResolutionPlots[i] = iBooker.book1D(title3, title3, 100, -1, 1);
52  theThetaPullPlots[i] = iBooker.book1D(title4, title4, 100, -5, 5);
53  thedXdZResolutionPlots[i] = iBooker.book1D(title5, title5, 100, -1, 1);
54  thedXdZPullPlots[i] = iBooker.book1D(title6, title6, 100, -5, 5);
55  thedYdZResolutionPlots[i] = iBooker.book1D(title7, title7, 100, -1, 1);
56  thedYdZPullPlots[i] = iBooker.book1D(title8, title8, 100, -5, 5);
57  }
58 }
59 
61 {
62  // get the collection of CSCRecHsegmentItrD
64  e.getByToken(segments_Token_, hRecHits);
65  const CSCSegmentCollection * cscRecHits = hRecHits.product();
66 
67  theChamberSegmentMap.clear();
68  unsigned nPerEvent = 0;
69  for(CSCSegmentCollection::const_iterator segmentItr = cscRecHits->begin();
70  segmentItr != cscRecHits->end(); segmentItr++)
71  {
72  ++nPerEvent;
73  int detId = segmentItr->geographicalId().rawId();
74  int chamberType = whatChamberType(detId);
75 
76  theNRecHitsPlot->Fill(segmentItr->nRecHits());
77  theNPerChamberTypePlot->Fill(chamberType);
78  theChamberSegmentMap[detId].push_back(*segmentItr);
79 
80  // do the resolution plots
81  const PSimHit * hit = keyHit(detId);
82  if(hit != nullptr)
83  {
84  const CSCLayer * layer = findLayer(hit->detUnitId());
85  plotResolution(*hit, *segmentItr, layer, chamberType);
86  }
87  }
88 
89  theNPerEventPlot->Fill(nPerEvent);
90 
93 }
94 
96 {
97  // now plot efficiency by looping over all chambers with hits
98  for(ChamberHitMap::const_iterator mapItr = theLayerHitsPerChamber.begin(),
99  mapEnd = theLayerHitsPerChamber.end();
100  mapItr != mapEnd;
101  ++mapItr)
102  {
103  int chamberId = mapItr->first;
104  int nHitsInChamber = mapItr->second.size();
105  bool isShower = (nHitsInChamber > theShowerThreshold);
106  bool hasSeg = hasSegment(chamberId);
107  int chamberType = whatChamberType(chamberId);
108  // find how many layers were hit in this chamber
109  std::vector<int> v = mapItr->second;
110  std::sort(v.begin(), v.end());
111  // maybe can just count
112  v.erase(std::unique(v.begin(), v.end()), v.end());
113  int nLayersHit = v.size();
114 
115  if(nLayersHit == 4)
116  {
117 
118  if(isShower) theTypePlot4HitsShower->Fill(chamberType);
119  else theTypePlot4HitsNoShower->Fill(chamberType);
120 
121  if(hasSeg)
122  {
123  if(isShower) theTypePlot4HitsShowerSeg->Fill(chamberType);
124  else theTypePlot4HitsNoShowerSeg->Fill(chamberType);
125  }
126  }
127 
128  if(nLayersHit == 5)
129  {
130 
131  if(isShower) theTypePlot5HitsShower->Fill(chamberType);
132  else theTypePlot5HitsNoShower->Fill(chamberType);
133 
134  if(hasSeg)
135  {
136  if(isShower) theTypePlot5HitsShowerSeg->Fill(chamberType);
137  else theTypePlot5HitsNoShowerSeg->Fill(chamberType);
138  }
139  }
140 
141  if(nLayersHit == 6)
142  {
143 
144  if(isShower) theTypePlot6HitsShower->Fill(chamberType);
145  else theTypePlot6HitsNoShower->Fill(chamberType);
146 
147  if(hasSeg)
148  {
149  if(isShower) theTypePlot6HitsShowerSeg->Fill(chamberType);
150  else theTypePlot6HitsNoShowerSeg->Fill(chamberType);
151  }
152  }
153 
154 
155  }
156 }
157 
158 bool CSCSegmentValidation::hasSegment(int chamberId) const
159 {
160  return (theChamberSegmentMap.find(chamberId) != theChamberSegmentMap.end());
161 }
162 
163 
165 {
166  CSCDetId cscDetId(detId);
167  return CSCChamberSpecs::whatChamberType(cscDetId.station(), cscDetId.ring());
168 }
169 
170 
172  const CSCLayer * layer, int chamberType)
173 {
174  GlobalPoint simHitPos = layer->toGlobal(simHit.localPosition());
175  GlobalPoint segmentPos = layer->toGlobal(segment.localPosition());
176  LocalVector simHitDir = simHit.localDirection();
177  LocalVector segmentDir = segment.localDirection();
178 
179  double dphi = segmentPos.phi() - simHitPos.phi();
180  double rdphi = segmentPos.perp() * dphi;
181  double dtheta = segmentPos.theta() - simHitPos.theta();
182 
183  double sigmax = sqrt(segment.localPositionError().xx());
184  //double sigmay = sqrt(segment.localPositionError().yy());
185 
186  double ddxdz = segmentDir.x()/segmentDir.z() - simHitDir.x()/simHitDir.z();
187  double ddydz = segmentDir.y()/segmentDir.z() - simHitDir.y()/simHitDir.z();
188  double sigmadxdz = sqrt(segment.localDirectionError().xx());
189  double sigmadydz = sqrt(segment.localDirectionError().yy());
190 
191  theRdPhiResolutionPlots[chamberType-1]->Fill( rdphi );
192  theRdPhiPullPlots[chamberType-1]->Fill( rdphi/sigmax );
193  theThetaResolutionPlots[chamberType-1]->Fill( dtheta );
194  //theThetaPullPlots[chamberType-1]->Fill( dy/sigmay );
195 
196  thedXdZResolutionPlots[chamberType-1]->Fill( ddxdz );
197  thedXdZPullPlots[chamberType-1]->Fill( ddxdz/sigmadxdz );
198  thedYdZResolutionPlots[chamberType-1]->Fill( ddydz );
199  thedYdZPullPlots[chamberType-1]->Fill( ddydz/sigmadydz );
200 }
201 
202 
204 {
205  theLayerHitsPerChamber.clear();
206  std::vector<int> layersHit = theSimHitMap->detsWithHits();
207  for(std::vector<int>::const_iterator layerItr = layersHit.begin(),
208  layersHitEnd = layersHit.end();
209  layerItr != layersHitEnd;
210  ++layerItr)
211  {
212  CSCDetId layerId(*layerItr);
213  CSCDetId chamberId = layerId.chamberId();
214  int nhits = theSimHitMap->hits(*layerItr).size();
215  // multiple entries, so we can see showers
216  for(int i = 0; i < nhits; ++i) {
217  theLayerHitsPerChamber[chamberId.rawId()].push_back(*layerItr);
218  }
219  }
220 
221 }
222 
224  bool SimHitPabsLessThan(const PSimHit & p1, const PSimHit & p2)
225  {
226  return p1.pabs() < p2.pabs();
227  }
228 }
229 
230 
231 const PSimHit * CSCSegmentValidation::keyHit(int chamberId) const
232 {
233  const PSimHit * result = nullptr;
234  int layerId = chamberId + 3;
235  const edm::PSimHitContainer & layerHits = theSimHitMap->hits(layerId);
236 
237  if(!layerHits.empty())
238  {
239  // pick the hit with maximum energy
240  edm::PSimHitContainer::const_iterator hitItr = std::max_element(layerHits.begin(), layerHits.end(),
242  result = &(*hitItr);
243  }
244  return result;
245 }
MonitorElement * theNPerEventPlot
float xx() const
Definition: LocalError.h:24
void plotResolution(const PSimHit &simHit, const CSCSegment &recHit, const CSCLayer *layer, int chamberType)
MonitorElement * theTypePlot5HitsNoShowerSeg
T perp() const
Definition: PV3DBase.h:72
LocalVector localDirection() const override
Local direction.
Definition: CSCSegment.h:41
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
MonitorElement * theTypePlot4HitsShowerSeg
MonitorElement * theNPerChamberTypePlot
static int whatChamberType(int istation, int iring)
static int whatChamberType(int detId)
MonitorElement * theTypePlot6HitsShower
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:54
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
LocalError localDirectionError() const override
Error on the local direction.
Definition: CSCSegment.cc:51
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:47
T y() const
Definition: PV3DBase.h:63
MonitorElement * theTypePlot5HitsNoShower
edm::EDGetTokenT< CSCSegmentCollection > segments_Token_
MonitorElement * theNRecHitsPlot
MonitorElement * theTypePlot4HitsShower
MonitorElement * theRdPhiResolutionPlots[10]
MonitorElement * theTypePlot6HitsShowerSeg
void Fill(long long x)
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
MonitorElement * theTypePlot4HitsNoShower
float yy() const
Definition: LocalError.h:26
Local3DPoint localPosition() const
Definition: PSimHit.h:44
const PSimHitMap * theSimHitMap
T sqrt(T t)
Definition: SSEVec.h:18
def unique(seq, keepstr=True)
Definition: tier0.py:24
MonitorElement * theTypePlot6HitsNoShowerSeg
T z() const
Definition: PV3DBase.h:64
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
const edm::PSimHitContainer & hits(int detId) const
Definition: PSimHitMap.cc:22
LocalPoint localPosition() const override
Definition: CSCSegment.h:38
CSCDetId chamberId() const
Definition: CSCDetId.h:53
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
MonitorElement * thedYdZResolutionPlots[10]
std::vector< int > detsWithHits() const
Definition: PSimHitMap.cc:37
MonitorElement * theRdPhiPullPlots[10]
LocalVector localDirection() const
Obsolete. Same as momentumAtEntry().unit(), for backward compatibility.
Definition: PSimHit.h:52
int ring() const
Definition: CSCDetId.h:75
MonitorElement * theThetaPullPlots[10]
T const * product() const
Definition: Handle.h:81
const PSimHit * keyHit(int chamberId) const
MonitorElement * thedXdZPullPlots[10]
double p1[4]
Definition: TauolaWrapper.h:89
bool hasSegment(int chamberId) const
MonitorElement * theTypePlot5HitsShowerSeg
int station() const
Definition: CSCDetId.h:86
MonitorElement * theThetaResolutionPlots[10]
ChamberSegmentMap theChamberSegmentMap
std::vector< PSimHit > PSimHitContainer
MonitorElement * theTypePlot4HitsNoShowerSeg
LocalError localPositionError() const override
Definition: CSCSegment.cc:47
void bookHistograms(DQMStore::IBooker &)
bool SimHitPabsLessThan(const PSimHit &p1, const PSimHit &p2)
T x() const
Definition: PV3DBase.h:62
CSCSegmentValidation(const edm::InputTag &inputTag, edm::ConsumesCollector &&iC)
void analyze(const edm::Event &, const edm::EventSetup &) override
MonitorElement * theTypePlot5HitsShower
const CSCLayer * findLayer(int detId) const
unsigned int detUnitId() const
Definition: PSimHit.h:93
MonitorElement * thedXdZResolutionPlots[10]
MonitorElement * thedYdZPullPlots[10]