#include <Validation/CSCRecHits/src/CSCRecHit2DValidation.h>
Public Member Functions | |
virtual void | analyze (const edm::Event &, const edm::EventSetup &) |
CSCRecHit2DValidation (DQMStore *dbe, const edm::InputTag &inputTag) | |
virtual | ~CSCRecHit2DValidation () |
Private Member Functions | |
void | plotResolution (const PSimHit &simHit, const CSCRecHit2D &recHit, const CSCLayer *layer, int chamberType) |
Private Attributes | |
MonitorElement * | theNPerEventPlot |
MonitorElement * | thePullPlots [10] |
MonitorElement * | theRecHitPosInStrip [10] |
MonitorElement * | theResolutionPlots [10] |
MonitorElement * | theScatterPlots [10] |
MonitorElement * | theSimHitPosInStrip [10] |
MonitorElement * | theSimHitScatterPlots [10] |
MonitorElement * | theYPullPlots [10] |
MonitorElement * | theYResolutionPlots [10] |
Definition at line 10 of file CSCRecHit2DValidation.h.
CSCRecHit2DValidation::CSCRecHit2DValidation | ( | DQMStore * | dbe, | |
const edm::InputTag & | inputTag | |||
) |
Definition at line 7 of file CSCRecHit2DValidation.cc.
References DQMStore::book1D(), DQMStore::book2D(), CSCBaseValidation::dbe_, i, DQMStore::setCurrentFolder(), thePullPlots, theRecHitPosInStrip, theResolutionPlots, theScatterPlots, theSimHitPosInStrip, theSimHitScatterPlots, theYPullPlots, and theYResolutionPlots.
00008 : CSCBaseValidation(dbe, inputTag), 00009 theNPerEventPlot( dbe_->book1D("CSCRecHitsPerEvent", "Number of CSC Rec Hits per event", 100, 0, 500) ) 00010 { 00011 dbe_->setCurrentFolder("CSCRecHitsV/CSCRecHitTask"); 00012 00013 for(int i = 0; i < 10; ++i) 00014 { 00015 char title1[200], title2[200], title3[200], title4[200], title5[200], title6[200], title7[200], title8[200]; 00016 sprintf(title1, "CSCRecHitResolution%d", i+1); 00017 sprintf(title2, "CSCRecHitPull%d", i+1); 00018 sprintf(title3, "CSCRecHitYResolution%d", i+1); 00019 sprintf(title4, "CSCRecHitYPull%d", i+1); 00020 sprintf(title5, "CSCRecHitPosInStrip%d", i+1); 00021 sprintf(title6, "CSCSimHitPosInStrip%d", i+1); 00022 sprintf(title7, "CSCRecHit%d", i+1); 00023 sprintf(title8, "CSCSimHit%d", i+1); 00024 00025 00026 theResolutionPlots[i] = dbe_->book1D(title1, title1, 100, -0.2, 0.2); 00027 thePullPlots[i] = dbe_->book1D(title2, title2, 100, -3, 3); 00028 theYResolutionPlots[i] = dbe_->book1D(title3, title3, 100, -5, 5); 00029 theYPullPlots[i] = dbe_->book1D(title4, title4, 100, -3, 3); 00030 theRecHitPosInStrip[i] = dbe_->book1D(title5, title5, 100, -2, 2); 00031 theSimHitPosInStrip[i] = dbe_->book1D(title6, title6, 100, -2, 2); 00032 theScatterPlots[i] = dbe->book2D(title7, title7, 200, -20, 20, 200, -250, 250); 00033 theSimHitScatterPlots[i] = dbe->book2D(title8, title8, 200, -20, 20, 200, -250, 250); 00034 } 00035 00036 }
CSCRecHit2DValidation::~CSCRecHit2DValidation | ( | ) | [virtual] |
Definition at line 38 of file CSCRecHit2DValidation.cc.
References GenMuonPlsPt100GeV_cfg::cout, lat::endl(), MonitorElement::getName(), MonitorElement::getRMS(), i, and theResolutionPlots.
00039 { 00040 for(int i = 0; i < 10; ++i) 00041 { 00042 std::cout << "Resolution of " << theResolutionPlots[i]->getName() << " is " << theResolutionPlots[i]->getRMS() << std::endl; 00043 } 00044 }
void CSCRecHit2DValidation::analyze | ( | const edm::Event & | e, | |
const edm::EventSetup & | eventSetup | |||
) | [virtual] |
Implements CSCBaseValidation.
Definition at line 47 of file CSCRecHit2DValidation.cc.
References CSCLayer::chamber(), CSCChamberSpecs::chamberType(), detId, PSimHitMap::detsWithHits(), MonitorElement::Fill(), CSCBaseValidation::findLayer(), edm::Event::getByLabel(), PSimHitMap::hits(), i, PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), plotResolution(), edm::Handle< T >::product(), DetId::rawId(), trackerHits::simHits, CSCChamber::specs(), CSCBaseValidation::theInputTag, theNPerEventPlot, theScatterPlots, CSCBaseValidation::theSimHitMap, theSimHitScatterPlots, and GeomDet::toGlobal().
Referenced by CSCRecHitValidation::analyze().
00048 { 00049 // get the collection of CSCRecHrecHitItrD 00050 edm::Handle<CSCRecHit2DCollection> hRecHits; 00051 e.getByLabel(theInputTag, hRecHits); 00052 const CSCRecHit2DCollection * cscRecHits = hRecHits.product(); 00053 00054 unsigned nPerEvent = 0; 00055 00056 for(CSCRecHit2DCollection::const_iterator recHitItr = cscRecHits->begin(); 00057 recHitItr != cscRecHits->end(); recHitItr++) 00058 { 00059 ++nPerEvent; 00060 int detId = (*recHitItr).cscDetId().rawId(); 00061 edm::PSimHitContainer simHits = theSimHitMap->hits(detId); 00062 const CSCLayer * layer = findLayer(detId); 00063 int chamberType = layer->chamber()->specs()->chamberType(); 00064 00065 00066 if(simHits.size() == 1) 00067 { 00068 plotResolution(simHits[0], *recHitItr, layer, chamberType); 00069 } 00070 float localX = recHitItr->localPosition().x(); 00071 float localY = recHitItr->localPosition().y(); 00072 //theYPlots[chamberType-1]->Fill(localY); 00073 // find a local phi 00074 float globalR = layer->toGlobal(LocalPoint(0.,0.,0.)).perp(); 00075 GlobalPoint axisThruChamber(globalR+localY, localX, 0.); 00076 float localPhi = axisThruChamber.phi().degrees(); 00077 //thePhiPlots[chamberType-1]->Fill(axisThruChamber.phi().degrees()); 00078 theScatterPlots[chamberType-1]->Fill( localPhi, localY); 00079 } 00080 theNPerEventPlot->Fill(nPerEvent); 00081 return; 00082 // fill sim hits 00083 std::vector<int> layersWithSimHits = theSimHitMap->detsWithHits(); 00084 for(int i = 0; i < layersWithSimHits.size(); ++i) 00085 { 00086 edm::PSimHitContainer simHits = theSimHitMap->hits(layersWithSimHits[i]); 00087 for(edm::PSimHitContainer::const_iterator hitItr = simHits.begin(); hitItr != simHits.end(); ++hitItr) 00088 { 00089 const CSCLayer * layer = findLayer(layersWithSimHits[i]); 00090 int chamberType = layer->chamber()->specs()->chamberType(); 00091 float localX = hitItr->localPosition().x(); 00092 float localY = hitItr->localPosition().y(); 00093 //theYPlots[chamberType-1]->Fill(localY); 00094 // find a local phi 00095 float globalR = layer->toGlobal(LocalPoint(0.,0.,0.)).perp(); 00096 GlobalPoint axisThruChamber(globalR+localY, localX, 0.); 00097 float localPhi = axisThruChamber.phi().degrees(); 00098 //thePhiPlots[chamberType-1]->Fill(axisThruChamber.phi().degrees()); 00099 theSimHitScatterPlots[chamberType-1]->Fill( localPhi, localY); 00100 } 00101 } 00102 00103 }
void CSCRecHit2DValidation::plotResolution | ( | const PSimHit & | simHit, | |
const CSCRecHit2D & | recHit, | |||
const CSCLayer * | layer, | |||
int | chamberType | |||
) | [private] |
Definition at line 106 of file CSCRecHit2DValidation.cc.
References MonitorElement::Fill(), CSCLayer::geometry(), PSimHit::localPosition(), CSCRecHit2D::localPosition(), CSCRecHit2D::localPositionError(), PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), funct::sqrt(), CSCLayerGeometry::strip(), thePullPlots, theRecHitPosInStrip, theResolutionPlots, theSimHitPosInStrip, theYPullPlots, theYResolutionPlots, GeomDet::toGlobal(), LocalError::xx(), PV3DBase< T, PVType, FrameType >::y(), and LocalError::yy().
Referenced by analyze().
00108 { 00109 GlobalPoint simHitPos = layer->toGlobal(simHit.localPosition()); 00110 GlobalPoint recHitPos = layer->toGlobal(recHit.localPosition()); 00111 00112 double dphi = recHitPos.phi() - simHitPos.phi(); 00113 double rdphi = recHitPos.perp() * dphi; 00114 theResolutionPlots[chamberType-1]->Fill( rdphi ); 00115 thePullPlots[chamberType-1]->Fill( rdphi/ sqrt(recHit.localPositionError().xx()) ); 00116 double dy = recHit.localPosition().y() - simHit.localPosition().y(); 00117 theYResolutionPlots[chamberType-1]->Fill( dy ); 00118 theYPullPlots[chamberType-1]->Fill( dy/ sqrt(recHit.localPositionError().yy()) ); 00119 00120 const CSCLayerGeometry * layerGeometry = layer->geometry(); 00121 float recStrip = layerGeometry->strip(recHit.localPosition()); 00122 float simStrip = layerGeometry->strip(simHit.localPosition()); 00123 theRecHitPosInStrip[chamberType-1]->Fill( recStrip - int(recStrip) ); 00124 theSimHitPosInStrip[chamberType-1]->Fill( simStrip - int(simStrip) ); 00125 }
MonitorElement* CSCRecHit2DValidation::thePullPlots[10] [private] |
Definition at line 25 of file CSCRecHit2DValidation.h.
Referenced by CSCRecHit2DValidation(), and plotResolution().
MonitorElement* CSCRecHit2DValidation::theRecHitPosInStrip[10] [private] |
Definition at line 30 of file CSCRecHit2DValidation.h.
Referenced by CSCRecHit2DValidation(), and plotResolution().
MonitorElement* CSCRecHit2DValidation::theResolutionPlots[10] [private] |
Definition at line 24 of file CSCRecHit2DValidation.h.
Referenced by CSCRecHit2DValidation(), plotResolution(), and ~CSCRecHit2DValidation().
MonitorElement* CSCRecHit2DValidation::theScatterPlots[10] [private] |
Definition at line 28 of file CSCRecHit2DValidation.h.
Referenced by analyze(), and CSCRecHit2DValidation().
MonitorElement* CSCRecHit2DValidation::theSimHitPosInStrip[10] [private] |
Definition at line 31 of file CSCRecHit2DValidation.h.
Referenced by CSCRecHit2DValidation(), and plotResolution().
MonitorElement* CSCRecHit2DValidation::theSimHitScatterPlots[10] [private] |
Definition at line 29 of file CSCRecHit2DValidation.h.
Referenced by analyze(), and CSCRecHit2DValidation().
MonitorElement* CSCRecHit2DValidation::theYPullPlots[10] [private] |
Definition at line 27 of file CSCRecHit2DValidation.h.
Referenced by CSCRecHit2DValidation(), and plotResolution().
MonitorElement* CSCRecHit2DValidation::theYResolutionPlots[10] [private] |
Definition at line 26 of file CSCRecHit2DValidation.h.
Referenced by CSCRecHit2DValidation(), and plotResolution().