CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_2_9_HLT1_bphpatch4/src/Validation/CSCRecHits/src/CSCSegmentValidation.cc

Go to the documentation of this file.
00001 #include "Validation/CSCRecHits/src/CSCSegmentValidation.h"
00002 #include "DataFormats/CSCRecHit/interface/CSCSegmentCollection.h"
00003 #include <algorithm>
00004 #include "DQMServices/Core/interface/DQMStore.h"
00005 
00006 
00007 
00008 CSCSegmentValidation::CSCSegmentValidation(DQMStore* dbe, const edm::InputTag & inputTag)
00009 : CSCBaseValidation(dbe, inputTag),
00010   theLayerHitsPerChamber(),
00011   theChamberSegmentMap(),
00012   theShowerThreshold(10),
00013   theNPerEventPlot( dbe_->book1D("CSCSegmentsPerEvent", "Number of CSC segments per event", 100, 0, 50) ),
00014   theNRecHitsPlot( dbe_->book1D("CSCRecHitsPerSegment", "Number of CSC rec hits per segment" , 8, 0, 7) ),
00015   theNPerChamberTypePlot( dbe_->book1D("CSCSegmentsPerChamberType", "Number of CSC segments per chamber type", 11, 0, 10) ),
00016   theTypePlot4HitsNoShower( dbe_->book1D("CSCSegments4HitsNoShower", "", 100, 0, 10) ),
00017   theTypePlot4HitsNoShowerSeg( dbe_->book1D("CSCSegments4HitsNoShowerSeg", "", 100, 0, 10) ),
00018   theTypePlot4HitsShower( dbe_->book1D("CSCSegments4HitsShower", "", 100, 0, 10) ),
00019   theTypePlot4HitsShowerSeg( dbe_->book1D("CSCSegments4HitsShowerSeg", "", 100, 0, 10) ),
00020   theTypePlot5HitsNoShower( dbe_->book1D("CSCSegments5HitsNoShower", "", 100, 0, 10) ),
00021   theTypePlot5HitsNoShowerSeg( dbe_->book1D("CSCSegments5HitsNoShowerSeg", "", 100, 0, 10) ),
00022   theTypePlot5HitsShower( dbe_->book1D("CSCSegments5HitsShower", "", 100, 0, 10) ),
00023   theTypePlot5HitsShowerSeg( dbe_->book1D("CSCSegments5HitsShowerSeg", "", 100, 0, 10) ),
00024   theTypePlot6HitsNoShower( dbe_->book1D("CSCSegments6HitsNoShower", "", 100, 0, 10) ),
00025   theTypePlot6HitsNoShowerSeg( dbe_->book1D("CSCSegments6HitsNoShowerSeg", "", 100, 0, 10) ),
00026   theTypePlot6HitsShower( dbe_->book1D("CSCSegments6HitsShower", "", 100, 0, 10) ),
00027   theTypePlot6HitsShowerSeg( dbe_->book1D("CSCSegments6HitsShowerSeg", "", 100, 0, 10) )
00028 {
00029    dbe_->setCurrentFolder("CSCRecHitsV/CSCRecHitTask");
00030    for(int i = 0; i < 10; ++i)
00031   {
00032     char title1[200], title2[200], title3[200], title4[200],  
00033          title5[200], title6[200], title7[200], title8[200];
00034     sprintf(title1, "CSCSegmentRdPhiResolution%d", i+1);
00035     sprintf(title2, "CSCSegmentRdPhiPull%d", i+1);
00036     sprintf(title3, "CSCSegmentThetaResolution%d", i+1);
00037     sprintf(title4, "CSCSegmentThetaPull%d", i+1);
00038     sprintf(title5, "CSCSegmentdXdZResolution%d", i+1);
00039     sprintf(title6, "CSCSegmentdXdZPull%d", i+1);
00040     sprintf(title7, "CSCSegmentdYdZResolution%d", i+1);
00041     sprintf(title8, "CSCSegmentdYdZPull%d", i+1);
00042 
00043 
00044     theRdPhiResolutionPlots[i] = dbe_->book1D(title1, title1, 100, -0.4, 0.4);
00045     theRdPhiPullPlots[i] = dbe_->book1D(title2, title2, 100, -5, 5);
00046     theThetaResolutionPlots[i] = dbe_->book1D(title3, title3, 100, -1, 1);
00047     theThetaPullPlots[i] = dbe_->book1D(title4, title4, 100, -5, 5);
00048     thedXdZResolutionPlots[i] = dbe_->book1D(title5, title5, 100, -1, 1);
00049     thedXdZPullPlots[i] = dbe_->book1D(title6, title6, 100, -5, 5);
00050     thedYdZResolutionPlots[i] = dbe_->book1D(title7, title7, 100, -1, 1);
00051     thedYdZPullPlots[i] = dbe_->book1D(title8, title8, 100, -5, 5);
00052 
00053   }
00054 }
00055 
00056 void CSCSegmentValidation::analyze(const edm::Event&e, const edm::EventSetup& eventSetup)
00057 {
00058   // get the collection of CSCRecHsegmentItrD
00059   edm::Handle<CSCSegmentCollection> hRecHits;
00060   e.getByLabel(theInputTag, hRecHits);
00061   const CSCSegmentCollection * cscRecHits = hRecHits.product();
00062 
00063   theChamberSegmentMap.clear();
00064   unsigned nPerEvent = 0;
00065   for(CSCSegmentCollection::const_iterator segmentItr = cscRecHits->begin(); 
00066       segmentItr != cscRecHits->end(); segmentItr++) 
00067   {
00068     ++nPerEvent;
00069     int detId = segmentItr->geographicalId().rawId();
00070     int chamberType = whatChamberType(detId);
00071 
00072     theNRecHitsPlot->Fill(segmentItr->nRecHits());
00073     theNPerChamberTypePlot->Fill(chamberType);
00074     theChamberSegmentMap[detId].push_back(*segmentItr);
00075 
00076     // do the resolution plots
00077     const PSimHit * hit = keyHit(detId);
00078     if(hit != 0) 
00079     {
00080       const CSCLayer * layer = findLayer(hit->detUnitId());
00081       plotResolution(*hit, *segmentItr, layer, chamberType);
00082     }  
00083   }
00084 
00085   theNPerEventPlot->Fill(nPerEvent);
00086 
00087   fillLayerHitsPerChamber();
00088   fillEfficiencyPlots();
00089 }
00090 
00091 
00092 void CSCSegmentValidation::fillEfficiencyPlots()
00093 {
00094     // now plot efficiency by looping over all chambers with hits
00095   for(ChamberHitMap::const_iterator mapItr = theLayerHitsPerChamber.begin(),
00096       mapEnd = theLayerHitsPerChamber.end();
00097       mapItr != mapEnd;
00098       ++mapItr)
00099   {
00100     int chamberId = mapItr->first;
00101     int nHitsInChamber = mapItr->second.size();
00102     bool isShower = (nHitsInChamber > theShowerThreshold);
00103     bool hasSeg = hasSegment(chamberId);
00104     int chamberType = whatChamberType(chamberId);
00105     // find how many layers were hit in this chamber
00106     std::vector<int> v = mapItr->second;
00107     std::sort(v.begin(), v.end());
00108     // maybe can just count
00109     v.erase(std::unique(v.begin(), v.end()), v.end());
00110     int nLayersHit = v.size();
00111 
00112     if(nLayersHit == 4)
00113     {
00114 
00115       if(isShower) theTypePlot4HitsShower->Fill(chamberType);
00116       else         theTypePlot4HitsNoShower->Fill(chamberType);
00117 
00118       if(hasSeg) 
00119       {
00120         if(isShower) theTypePlot4HitsShowerSeg->Fill(chamberType);
00121         else         theTypePlot4HitsNoShowerSeg->Fill(chamberType);
00122       } 
00123     }
00124 
00125     if(nLayersHit == 5)
00126     {
00127 
00128       if(isShower) theTypePlot5HitsShower->Fill(chamberType);
00129       else         theTypePlot5HitsNoShower->Fill(chamberType);
00130 
00131       if(hasSeg)
00132       {
00133         if(isShower) theTypePlot5HitsShowerSeg->Fill(chamberType);
00134         else         theTypePlot5HitsNoShowerSeg->Fill(chamberType);
00135       }
00136     }
00137 
00138     if(nLayersHit == 6)
00139     {
00140 
00141       if(isShower) theTypePlot6HitsShower->Fill(chamberType);
00142       else         theTypePlot6HitsNoShower->Fill(chamberType);
00143 
00144       if(hasSeg)
00145       {
00146         if(isShower) theTypePlot6HitsShowerSeg->Fill(chamberType);
00147         else         theTypePlot6HitsNoShowerSeg->Fill(chamberType);
00148       }
00149     }
00150 
00151 
00152   }
00153 }
00154 
00155 bool CSCSegmentValidation::hasSegment(int chamberId) const
00156 {
00157   return (theChamberSegmentMap.find(chamberId) != theChamberSegmentMap.end());
00158 }
00159 
00160 
00161 int CSCSegmentValidation::whatChamberType(int detId)
00162 {
00163   CSCDetId cscDetId(detId);
00164   return CSCChamberSpecs::whatChamberType(cscDetId.station(), cscDetId.ring());
00165 }
00166 
00167 
00168 void CSCSegmentValidation::plotResolution(const PSimHit & simHit, const CSCSegment & segment,
00169                                          const CSCLayer * layer, int chamberType)
00170 {
00171   GlobalPoint simHitPos = layer->toGlobal(simHit.localPosition());
00172   GlobalPoint segmentPos = layer->toGlobal(segment.localPosition());
00173   LocalVector simHitDir = simHit.localDirection();
00174   LocalVector segmentDir = segment.localDirection();
00175 
00176   double dphi = segmentPos.phi() - simHitPos.phi();
00177   double rdphi = segmentPos.perp() * dphi;
00178   double dtheta = segmentPos.theta() - simHitPos.theta();
00179 
00180   double sigmax = sqrt(segment.localPositionError().xx());
00181   //double sigmay = sqrt(segment.localPositionError().yy());
00182 
00183   double ddxdz = segmentDir.x()/segmentDir.z() - simHitDir.x()/simHitDir.z();
00184   double ddydz = segmentDir.y()/segmentDir.z() - simHitDir.y()/simHitDir.z();
00185   double sigmadxdz = sqrt(segment.localDirectionError().xx());
00186   double sigmadydz = sqrt(segment.localDirectionError().yy()); 
00187 
00188   theRdPhiResolutionPlots[chamberType-1]->Fill( rdphi );
00189   theRdPhiPullPlots[chamberType-1]->Fill( rdphi/sigmax );
00190   theThetaResolutionPlots[chamberType-1]->Fill( dtheta );
00191   //theThetaPullPlots[chamberType-1]->Fill( dy/sigmay );
00192 
00193   thedXdZResolutionPlots[chamberType-1]->Fill( ddxdz );
00194   thedXdZPullPlots[chamberType-1]->Fill( ddxdz/sigmadxdz );
00195   thedYdZResolutionPlots[chamberType-1]->Fill( ddydz );
00196   thedYdZPullPlots[chamberType-1]->Fill( ddydz/sigmadydz );
00197 }
00198 
00199 
00200 void CSCSegmentValidation::fillLayerHitsPerChamber()
00201 {
00202   theLayerHitsPerChamber.clear();
00203   std::vector<int> layersHit = theSimHitMap->detsWithHits();
00204   for(std::vector<int>::const_iterator layerItr = layersHit.begin(), 
00205       layersHitEnd = layersHit.end();
00206       layerItr != layersHitEnd;
00207       ++layerItr)
00208   {
00209     CSCDetId layerId(*layerItr);
00210     CSCDetId chamberId = layerId.chamberId();
00211     int nhits = theSimHitMap->hits(*layerItr).size();
00212     // multiple entries, so we can see showers
00213     for(int i = 0; i < nhits; ++i) {
00214       theLayerHitsPerChamber[chamberId.rawId()].push_back(*layerItr);
00215     }
00216   }
00217 
00218 }
00219 
00220 namespace CSCSegmentValidationUtils {
00221   bool SimHitPabsLessThan(const PSimHit & p1, const PSimHit & p2)
00222   {
00223     return p1.pabs() < p2.pabs();
00224   }
00225 }
00226 
00227 
00228 const PSimHit * CSCSegmentValidation::keyHit(int chamberId) const
00229 {
00230   const PSimHit * result = 0;
00231   int layerId = chamberId + 3;
00232   const edm::PSimHitContainer & layerHits = theSimHitMap->hits(layerId);
00233 
00234   if(!layerHits.empty())
00235   {
00236     // pick the hit with maximum energy
00237     edm::PSimHitContainer::const_iterator hitItr = std::max_element(layerHits.begin(), layerHits.end(),
00238                                                     CSCSegmentValidationUtils::SimHitPabsLessThan);
00239     result = &(*hitItr);
00240   }
00241   return result;
00242 }   
00243