CMS 3D CMS Logo

CSCRecHit2DValidation.cc
Go to the documentation of this file.
4 
5 
6 
8  const edm::InputTag & inputTag,
10 : CSCBaseValidation(inputTag),
11  theNPerEventPlot(0)
12 {
14 }
15 
17 {
18  for(int i = 0; i < 10; ++i)
19  {
20  edm::LogInfo("CSCRecHitValidation") << "Resolution of " << theResolutionPlots[i]->getName() << " is " << theResolutionPlots[i]->getRMS();
21  edm::LogInfo("CSCRecHitValidation") << "Peak Time is " << theTPeaks[i]->getMean();
22  }
23 }
24 
26 {
27  theNPerEventPlot = iBooker.book1D("CSCRecHitsPerEvent", "Number of CSC Rec Hits per event", 100, 0, 500);
28  for(int i = 0; i < 10; ++i)
29  {
30  char title1[200], title2[200], title3[200], title4[200], title5[200], title6[200], title7[200], title8[200], title9[200];
31  sprintf(title1, "CSCRecHitResolution%d", i+1);
32  sprintf(title2, "CSCRecHitPull%d", i+1);
33  sprintf(title3, "CSCRecHitYResolution%d", i+1);
34  sprintf(title4, "CSCRecHitYPull%d", i+1);
35  sprintf(title5, "CSCRecHitPosInStrip%d", i+1);
36  sprintf(title6, "CSCSimHitPosInStrip%d", i+1);
37  sprintf(title7, "CSCRecHit%d", i+1);
38  sprintf(title8, "CSCSimHit%d", i+1);
39  sprintf(title9, "CSCTPeak%d", i+1);
40 
41  theResolutionPlots[i] = iBooker.book1D(title1, title1, 100, -0.2, 0.2);
42  thePullPlots[i] = iBooker.book1D(title2, title2, 100, -3, 3);
43  theYResolutionPlots[i] = iBooker.book1D(title3, title3, 100, -5, 5);
44  theYPullPlots[i] = iBooker.book1D(title4, title4, 100, -3, 3);
45  theRecHitPosInStrip[i] = iBooker.book1D(title5, title5, 100, -2, 2);
46  theSimHitPosInStrip[i] = iBooker.book1D(title6, title6, 100, -2, 2);
47  theScatterPlots[i] = iBooker.book2D(title7, title7, 200, -20, 20, 200, -250, 250);
48  theSimHitScatterPlots[i] = iBooker.book2D(title8, title8, 200, -20, 20, 200, -250, 250);
49  theTPeaks[i] = iBooker.book1D(title9, title9, 200, 0, 400);
50  }
51 }
52 
54 {
55  // get the collection of CSCRecHrecHitItrD
57  e.getByToken(rechits_Token_, hRecHits);
58  const CSCRecHit2DCollection * cscRecHits = hRecHits.product();
59 
60  unsigned nPerEvent = 0;
61 
62  for(CSCRecHit2DCollection::const_iterator recHitItr = cscRecHits->begin();
63  recHitItr != cscRecHits->end(); recHitItr++)
64  {
65  ++nPerEvent;
66  int detId = (*recHitItr).cscDetId().rawId();
68  const CSCLayer * layer = findLayer(detId);
69  int chamberType = layer->chamber()->specs()->chamberType();
70  theTPeaks[chamberType-1]->Fill(recHitItr->tpeak());
71  if(simHits.size() == 1)
72  {
73  plotResolution(simHits[0], *recHitItr, layer, chamberType);
74  }
75  float localX = recHitItr->localPosition().x();
76  float localY = recHitItr->localPosition().y();
77  //theYPlots[chamberType-1]->Fill(localY);
78  // find a local phi
79  float globalR = layer->toGlobal(LocalPoint(0.,0.,0.)).perp();
80  GlobalPoint axisThruChamber(globalR+localY, localX, 0.);
81  float localPhi = axisThruChamber.phi().degrees();
82  //thePhiPlots[chamberType-1]->Fill(axisThruChamber.phi().degrees());
83  theScatterPlots[chamberType-1]->Fill( localPhi, localY);
84  }
85  theNPerEventPlot->Fill(nPerEvent);
86 return;
87  // fill sim hits
88  std::vector<int> layersWithSimHits = theSimHitMap->detsWithHits();
89  for(unsigned i = 0; i < layersWithSimHits.size(); ++i)
90  {
91  edm::PSimHitContainer simHits = theSimHitMap->hits(layersWithSimHits[i]);
92  for(edm::PSimHitContainer::const_iterator hitItr = simHits.begin(); hitItr != simHits.end(); ++hitItr)
93  {
94  const CSCLayer * layer = findLayer(layersWithSimHits[i]);
95  int chamberType = layer->chamber()->specs()->chamberType();
96  float localX = hitItr->localPosition().x();
97  float localY = hitItr->localPosition().y();
98  //theYPlots[chamberType-1]->Fill(localY);
99  // find a local phi
100  float globalR = layer->toGlobal(LocalPoint(0.,0.,0.)).perp();
101  GlobalPoint axisThruChamber(globalR+localY, localX, 0.);
102  float localPhi = axisThruChamber.phi().degrees();
103  //thePhiPlots[chamberType-1]->Fill(axisThruChamber.phi().degrees());
104  theSimHitScatterPlots[chamberType-1]->Fill( localPhi, localY);
105  }
106  }
107 
108 }
109 
111  const CSCLayer * layer, int chamberType)
112 {
113  GlobalPoint simHitPos = layer->toGlobal(simHit.localPosition());
114  GlobalPoint recHitPos = layer->toGlobal(recHit.localPosition());
115 
116  double dphi = recHitPos.phi() - simHitPos.phi();
117  double rdphi = recHitPos.perp() * dphi;
118  theResolutionPlots[chamberType-1]->Fill( rdphi );
119  thePullPlots[chamberType-1]->Fill( rdphi/ sqrt(recHit.localPositionError().xx()) );
120  double dy = recHit.localPosition().y() - simHit.localPosition().y();
121  theYResolutionPlots[chamberType-1]->Fill( dy );
122  theYPullPlots[chamberType-1]->Fill( dy/ sqrt(recHit.localPositionError().yy()) );
123 
124  const CSCLayerGeometry * layerGeometry = layer->geometry();
125  float recStrip = layerGeometry->strip(recHit.localPosition());
126  float simStrip = layerGeometry->strip(simHit.localPosition());
127  theRecHitPosInStrip[chamberType-1]->Fill( recStrip - int(recStrip) );
128  theSimHitPosInStrip[chamberType-1]->Fill( simStrip - int(simStrip) );
129 }
130 
MonitorElement * theSimHitScatterPlots[10]
MonitorElement * theYResolutionPlots[10]
const std::string & getName(void) const
get name of ME
float xx() const
Definition: LocalError.h:24
MonitorElement * theScatterPlots[10]
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:32
T perp() const
Definition: PV3DBase.h:72
void plotResolution(const PSimHit &simHit, const CSCRecHit2D &recHit, const CSCLayer *layer, int chamberType)
edm::EDGetTokenT< CSCRecHit2DCollection > rechits_Token_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
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
T y() const
Definition: PV3DBase.h:63
void bookHistograms(DQMStore::IBooker &)
MonitorElement * theResolutionPlots[10]
double getMean(int axis=1) const
get mean value of histogram along x, y or z axis (axis=1, 2, 3 respectively)
CSCRecHit2DValidation(const edm::InputTag &inputTag, edm::ConsumesCollector &&iC)
void Fill(long long x)
MonitorElement * theSimHitPosInStrip[10]
MonitorElement * theYPullPlots[10]
const CSCChamberSpecs * specs() const
Definition: CSCChamber.h:42
LocalError localPositionError() const
Definition: CSCRecHit2D.h:51
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
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
const edm::PSimHitContainer & hits(int detId) const
Definition: PSimHitMap.cc:22
float strip(const LocalPoint &lp) const
virtual void analyze(const edm::Event &, const edm::EventSetup &)
std::vector< int > detsWithHits() const
Definition: PSimHitMap.cc:37
T const * product() const
Definition: Handle.h:81
int chamberType() const
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
MonitorElement * theNPerEventPlot
MonitorElement * thePullPlots[10]
MonitorElement * theRecHitPosInStrip[10]
double getRMS(int axis=1) const
get RMS of histogram along x, y or z axis (axis=1, 2, 3 respectively)
MonitorElement * theTPeaks[10]
std::vector< PSimHit > PSimHitContainer
LocalPoint localPosition() const
Definition: CSCRecHit2D.h:50
const CSCChamber * chamber() const
Definition: CSCLayer.h:52
const CSCLayerGeometry * geometry() const
Definition: CSCLayer.h:47
const CSCLayer * findLayer(int detId) const