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.
2 #include <algorithm>
4 
5 
6 
8 : CSCBaseValidation(dbe, inputTag),
9  theLayerHitsPerChamber(),
10  theChamberSegmentMap(),
11  theShowerThreshold(10),
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 {
28  segments_Token_ = iC.consumes<CSCSegmentCollection>(inputTag);
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 }
56 
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 }
91 
92 
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 }
155 
156 bool CSCSegmentValidation::hasSegment(int chamberId) const
157 {
158  return (theChamberSegmentMap.find(chamberId) != theChamberSegmentMap.end());
159 }
160 
161 
163 {
164  CSCDetId cscDetId(detId);
165  return CSCChamberSpecs::whatChamberType(cscDetId.station(), cscDetId.ring());
166 }
167 
168 
169 void CSCSegmentValidation::plotResolution(const PSimHit & simHit, const CSCSegment & segment,
170  const CSCLayer * layer, int chamberType)
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 }
199 
200 
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 }
220 
221 namespace CSCSegmentValidationUtils {
222  bool SimHitPabsLessThan(const PSimHit & p1, const PSimHit & p2)
223  {
224  return p1.pabs() < p2.pabs();
225  }
226 }
227 
228 
229 const PSimHit * CSCSegmentValidation::keyHit(int chamberId) const
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 }
244 
MonitorElement * theNPerEventPlot
int i
Definition: DBlmapReader.cc:9
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
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:954
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
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:47
LocalError localPositionError() const
Definition: CSCSegment.cc:47
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
MonitorElement * theTypePlot5HitsNoShower
edm::EDGetTokenT< CSCSegmentCollection > segments_Token_
MonitorElement * theNRecHitsPlot
MonitorElement * theTypePlot4HitsShower
MonitorElement * theRdPhiResolutionPlots[10]
CSCSegmentValidation(DQMStore *dbe, const edm::InputTag &inputTag, edm::ConsumesCollector &&iC)
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
LocalPoint localPosition() const
Definition: CSCSegment.h:38
LocalVector localDirection() const
Local direction.
Definition: CSCSegment.h:41
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: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:22
CSCDetId chamberId() const
Definition: CSCDetId.h:66
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
DQMStore * dbe_
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:88
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:99
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:667
LocalError localDirectionError() const
Error on the local direction.
Definition: CSCSegment.cc:51
unsigned int detUnitId() const
Definition: PSimHit.h:93
MonitorElement * thedXdZResolutionPlots[10]
MonitorElement * thedYdZPullPlots[10]