CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch1/src/Validation/MuonCSCDigis/src/CSCStripDigiValidation.cc

Go to the documentation of this file.
00001 #include "Validation/MuonCSCDigis/src/CSCStripDigiValidation.h"
00002 #include "DataFormats/Common/interface/Handle.h"
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 #include "DataFormats/CSCDigi/interface/CSCStripDigiCollection.h"
00005 #include "DQMServices/Core/interface/DQMStore.h"
00006 
00007 
00008 CSCStripDigiValidation::CSCStripDigiValidation(DQMStore* dbe, 
00009                                                const edm::InputTag & inputTag,
00010                                                bool doSim)
00011 : CSCBaseValidation(dbe, inputTag),
00012   thePedestalSum(0),
00013   thePedestalCovarianceSum(0),
00014   thePedestalCount(0),
00015   theDoSimFlag(doSim),
00016   thePedestalPlot( dbe_->book1D("CSCPedestal", "CSC Pedestal ", 400, 550, 650) ),
00017   thePedestalTimeCorrelationPlot(0),
00018   thePedestalNeighborCorrelationPlot(0),
00019   theAmplitudePlot( dbe_->book1D("CSCStripAmplitude", "CSC Strip Amplitude", 200, 0, 2000) ),
00020   theRatio4to5Plot( dbe_->book1D("CSCStrip4to5", "CSC Strip Ratio tbin 4 to tbin 5", 100, 0, 1) ),
00021   theRatio6to5Plot( dbe_->book1D("CSCStrip6to5", "CSC Strip Ratio tbin 6 to tbin 5", 120, 0, 1.2) ),
00022   theNDigisPerLayerPlot( dbe_->book1D("CSCStripDigisPerLayer", "Number of CSC Strip Digis per layer", 48, 0, 48) ),
00023   theNDigisPerChamberPlot(0),
00024   theNDigisPerEventPlot( dbe_->book1D("CSCStripDigisPerEvent", "Number of CSC Strip Digis per event", 100, 0, 500) )
00025 {
00026   if(doSim) {
00027     for(int i = 0; i < 10; ++i)
00028     {
00029       char title1[200];
00030       sprintf(title1, "CSCStripDigiResolution%d", i+1);
00031       theResolutionPlots[i] = dbe_->book1D(title1, title1, 100, -5, 5);
00032     }
00033   }
00034 
00035 }
00036 
00037 
00038 CSCStripDigiValidation::~CSCStripDigiValidation() {
00039  
00040 //   edm::LogInfo("CSCDigiValidation") << "RATIO for strips 4 to 5 : " << theRatio4to5Plot->getMean();
00041 //   edm::LogInfo("CSCDigiValidation") << "RATIO for strips 6 to 5 : " << theRatio6to5Plot->getMean();
00042 //   edm::LogInfo("CSCDigiValidation") << "NDIGIS per event : " << theNDigisPerEventPlot->getMean();
00043 
00044 }
00045 
00046 
00047 void CSCStripDigiValidation::analyze(const edm::Event& e, 
00048                                      const edm::EventSetup&)
00049 {
00050   edm::Handle<CSCStripDigiCollection> strips;
00051   e.getByLabel(theInputTag, strips);
00052   if (!strips.isValid()) {
00053     edm::LogError("CSCDigiValidation") << "Cannot get strips by label " 
00054                                        << theInputTag.encode();
00055   }
00056 
00057  unsigned nDigisPerEvent = 0;
00058 
00059  for (CSCStripDigiCollection::DigiRangeIterator j=strips->begin(); j!=strips->end(); j++) {
00060     std::vector<CSCStripDigi>::const_iterator digiItr = (*j).second.first;
00061     std::vector<CSCStripDigi>::const_iterator last = (*j).second.second;
00062     int nDigis = last-digiItr;
00063     nDigisPerEvent += nDigis;
00064     theNDigisPerLayerPlot->Fill(nDigis);
00065 
00066     double maxAmplitude = 0.;
00067     // int maxStrip = 0;
00068 
00069     for( ; digiItr != last; ++digiItr) {
00070       // average up the pedestals
00071       std::vector<int> adcCounts = digiItr->getADCCounts();
00072       thePedestalSum += adcCounts[0];
00073       thePedestalSum += adcCounts[1];
00074       thePedestalCount += 2;
00075       float pedestal = thePedestalSum/thePedestalCount;
00076       if(adcCounts[4]-pedestal > maxAmplitude)
00077       {
00078       //  maxStrip = digiItr->getStrip();
00079         maxAmplitude = adcCounts[4]-pedestal;
00080       }
00081 
00082       // if we have enough pedestal statistics
00083       if(thePedestalCount > 100)
00084       {
00085         fillPedestalPlots(*digiItr);
00086        
00087         // see if it's big enough to count as "signal"
00088         if(adcCounts[5] > (thePedestalSum/thePedestalCount + 100))
00089         {
00090           fillSignalPlots(*digiItr);
00091         }
00092       }
00093     }
00094 
00095 /*
00096     int detId = (*j).first.rawId();
00097 
00098     edm::PSimHitContainer simHits = theSimHitMap->hits(detId);
00099 
00100     if(simHits.size() == 1)
00101     {
00102       const CSCLayer * layer = findLayer(detId);
00103       int chamberType = layer->chamber()->specs()->chamberType();
00104       plotResolution(simHits[0], maxStrip, layer, chamberType);
00105     }
00106 */
00107   } // loop over digis
00108 
00109   theNDigisPerEventPlot->Fill(nDigisPerEvent);
00110 }
00111 
00112 
00113 void CSCStripDigiValidation::fillPedestalPlots(const CSCStripDigi & digi)
00114 {
00115   std::vector<int> adcCounts = digi.getADCCounts();
00116   thePedestalPlot->Fill(adcCounts[0]);
00117   thePedestalPlot->Fill(adcCounts[1]);
00118 }
00119 
00120 
00121 
00122 void CSCStripDigiValidation::fillSignalPlots(const CSCStripDigi & digi)
00123 {
00124   std::vector<int> adcCounts = digi.getADCCounts();
00125   float pedestal = thePedestalSum/thePedestalCount;
00126   theAmplitudePlot->Fill(adcCounts[4] - pedestal);
00127   theRatio4to5Plot->Fill( (adcCounts[3]-pedestal) / (adcCounts[4]-pedestal) );
00128   theRatio6to5Plot->Fill( (adcCounts[5]-pedestal) / (adcCounts[4]-pedestal) );
00129 }
00130 
00131 
00132 void CSCStripDigiValidation::plotResolution(const PSimHit & hit, int strip,
00133                                            const CSCLayer * layer, int chamberType)
00134 {
00135   double hitX = hit.localPosition().x();
00136   double hitY = hit.localPosition().y();
00137   double digiX = layer->geometry()->xOfStrip(strip, hitY);
00138   theResolutionPlots[chamberType-1]->Fill(digiX - hitX);
00139 }
00140 
00141