CMS 3D CMS Logo

CSCComparatorDigiValidation.cc

Go to the documentation of this file.
00001 #include "Validation/MuonCSCDigis/src/CSCComparatorDigiValidation.h"
00002 #include "DataFormats/Common/interface/Handle.h"
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 #include "DataFormats/CSCDigi/interface/CSCComparatorDigiCollection.h"
00005 #include "DataFormats/CSCDigi/interface/CSCStripDigiCollection.h"
00006 #include "DQMServices/Core/interface/DQMStore.h"
00007 
00008 
00009 CSCComparatorDigiValidation::CSCComparatorDigiValidation(DQMStore* dbe, 
00010     const edm::InputTag & inputTag, const edm::InputTag & stripDigiInputTag)
00011 : CSCBaseValidation(dbe, inputTag),
00012   theStripDigiInputTag(stripDigiInputTag),
00013   theTimeBinPlots(),
00014   theNDigisPerLayerPlots(),
00015   theStripDigiPlots(),
00016   the3StripPlots(),
00017   theNDigisPerEventPlot( dbe_->book1D("CSCComparatorDigisPerEvent", "CSC Comparator Digis per event", 100, 0, 100) )
00018 {
00019   for(int i = 0; i < 10; ++i)
00020   {
00021     char title1[200], title2[200], title3[200], title4[200];
00022     sprintf(title1, "CSCComparatorDigiTimeType%d", i+1);
00023     sprintf(title2, "CSCComparatorDigisPerLayerType%d", i+1);
00024     sprintf(title3, "CSCComparatorStripAmplitudeType%d", i+1);
00025     sprintf(title4, "CSCComparator3StripAmplitudeType%d", i+1);
00026     theTimeBinPlots[i] = dbe_->book1D(title1, title1, 9, 0, 8);
00027     theNDigisPerLayerPlots[i] = dbe_->book1D(title2, title2, 100, 0, 20);
00028     theStripDigiPlots[i] = dbe_->book1D(title3, title3, 100, 0, 1000);
00029     the3StripPlots[i] = dbe_->book1D(title4, title4, 100, 0, 1000);
00030 
00031   }
00032 }
00033 
00034 
00035 
00036 CSCComparatorDigiValidation::~CSCComparatorDigiValidation()
00037 {
00038   /*
00039   for(int i = 0; i < 10; ++i)
00040   {
00041     std::cout << "Mean of " << theTimeBinPlots[i]->getName()
00042       << " is " << theTimeBinPlots[i]->getMean()
00043       << " +/- " << theTimeBinPlots[i]->getRMS() << std::endl;
00044     std::cout << "Mean charge of " << the3StripPlots[i]->getName()
00045       << " is " << the3StripPlots[i]->getMean() << std::endl;
00046   }
00047   */
00048 }
00049 
00050 
00051 void CSCComparatorDigiValidation::analyze(const edm::Event&e, const edm::EventSetup&)
00052 {
00053   edm::Handle<CSCComparatorDigiCollection> comparators;
00054   edm::Handle<CSCStripDigiCollection> stripDigis;
00055 
00056   e.getByLabel(theInputTag, comparators);
00057   if (!comparators.isValid()) {
00058     edm::LogError("CSCDigiDump") << "Cannot get comparators by label " << theInputTag.encode();
00059   }
00060   e.getByLabel(theStripDigiInputTag, stripDigis);
00061   if (!stripDigis.isValid()) {
00062     edm::LogError("CSCDigiDump") << "Cannot get comparators by label " << theInputTag.encode();
00063   }
00064   
00065   unsigned nDigisPerEvent = 0;
00066 
00067   for (CSCComparatorDigiCollection::DigiRangeIterator j=comparators->begin(); j!=comparators->end(); j++) {
00068     std::vector<CSCComparatorDigi>::const_iterator digiItr = (*j).second.first;
00069     std::vector<CSCComparatorDigi>::const_iterator last = (*j).second.second;
00070 
00071     CSCDetId detId((*j).first);
00072     const CSCLayer * layer = findLayer(detId.rawId());
00073     int chamberType = layer->chamber()->specs()->chamberType();
00074 
00075     CSCStripDigiCollection::Range stripRange = stripDigis->get(detId);
00076 
00077     theNDigisPerLayerPlots[chamberType-1]->Fill(last-digiItr);
00078 
00079     for( ; digiItr != last; ++digiItr) {
00080       ++nDigisPerEvent;
00081       theTimeBinPlots[chamberType-1]->Fill(digiItr->getTimeBin());
00082 
00083       int strip = digiItr->getStrip();
00084       for(std::vector<CSCStripDigi>::const_iterator stripItr = stripRange.first;
00085           stripItr != stripRange.second; ++stripItr)
00086       {
00087         if(stripItr->getStrip() == strip)
00088         {
00089           std::vector<int> adc = stripItr->getADCCounts();
00090           float pedc = 0.5*(adc[0]+adc[1]);
00091           float amp = adc[4] - pedc;
00092           theStripDigiPlots[chamberType-1]->Fill(amp);
00093           // check neighbors
00094           if(stripItr != stripRange.first && stripItr != stripRange.second-1)
00095           {
00096             std::vector<int> adcl = (stripItr-1)->getADCCounts();
00097             std::vector<int> adcr = (stripItr+1)->getADCCounts();
00098             float pedl = 0.5*(adcl[0]+adcl[1]);
00099             float pedr = 0.5*(adcr[0]+adcr[1]);
00100             float three = adcl[4]-pedl
00101                         + adcr[4]-pedr
00102                         + amp;
00103             the3StripPlots[chamberType-1]->Fill(three);
00104           }
00105         }
00106       }
00107     }
00108 
00109   }
00110 
00111   theNDigisPerEventPlot->Fill(nDigisPerEvent);
00112 }
00113 

Generated on Tue Jun 9 17:49:24 2009 for CMSSW by  doxygen 1.5.4