CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/src/Validation/MuonCSCDigis/src/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 //   for(int i = 0; i < 10; ++i)
00039 //   {
00040 //     edm::LogInfo("CSCDigiValidation") << "Mean of " << theTimeBinPlots[i]->getName()
00041 //       << " is " << theTimeBinPlots[i]->getMean()
00042 //       << " +/- " << theTimeBinPlots[i]->getRMS();
00043 //     edm::LogInfo("CSCDigiValidation") << "Mean charge of " << the3StripPlots[i]->getName()
00044 //       << " is " << the3StripPlots[i]->getMean();
00045 //   }
00046 }
00047 
00048 
00049 void CSCComparatorDigiValidation::analyze(const edm::Event&e, const edm::EventSetup&)
00050 {
00051   edm::Handle<CSCComparatorDigiCollection> comparators;
00052   edm::Handle<CSCStripDigiCollection> stripDigis;
00053 
00054   e.getByLabel(theInputTag, comparators);
00055   if (!comparators.isValid()) {
00056     edm::LogError("CSCDigiDump") << "Cannot get comparators by label " << theInputTag.encode();
00057   }
00058   e.getByLabel(theStripDigiInputTag, stripDigis);
00059   if (!stripDigis.isValid()) {
00060     edm::LogError("CSCDigiDump") << "Cannot get comparators by label " << theInputTag.encode();
00061   }
00062   
00063   unsigned nDigisPerEvent = 0;
00064 
00065   for (CSCComparatorDigiCollection::DigiRangeIterator j=comparators->begin(); j!=comparators->end(); j++) {
00066     std::vector<CSCComparatorDigi>::const_iterator digiItr = (*j).second.first;
00067     std::vector<CSCComparatorDigi>::const_iterator last = (*j).second.second;
00068 
00069     CSCDetId detId((*j).first);
00070     const CSCLayer * layer = findLayer(detId.rawId());
00071     int chamberType = layer->chamber()->specs()->chamberType();
00072 
00073     CSCStripDigiCollection::Range stripRange = stripDigis->get(detId);
00074 
00075     theNDigisPerLayerPlots[chamberType-1]->Fill(last-digiItr);
00076 
00077     for( ; digiItr != last; ++digiItr) {
00078       ++nDigisPerEvent;
00079       theTimeBinPlots[chamberType-1]->Fill(digiItr->getTimeBin());
00080 
00081       int strip = digiItr->getStrip();
00082       for(std::vector<CSCStripDigi>::const_iterator stripItr = stripRange.first;
00083           stripItr != stripRange.second; ++stripItr)
00084       {
00085         if(stripItr->getStrip() == strip)
00086         {
00087           std::vector<int> adc = stripItr->getADCCounts();
00088           float pedc = 0.5*(adc[0]+adc[1]);
00089           float amp = adc[4] - pedc;
00090           theStripDigiPlots[chamberType-1]->Fill(amp);
00091           // check neighbors
00092           if(stripItr != stripRange.first && stripItr != stripRange.second-1)
00093           {
00094             std::vector<int> adcl = (stripItr-1)->getADCCounts();
00095             std::vector<int> adcr = (stripItr+1)->getADCCounts();
00096             float pedl = 0.5*(adcl[0]+adcl[1]);
00097             float pedr = 0.5*(adcr[0]+adcr[1]);
00098             float three = adcl[4]-pedl
00099                         + adcr[4]-pedr
00100                         + amp;
00101             the3StripPlots[chamberType-1]->Fill(three);
00102           }
00103         }
00104       }
00105     }
00106 
00107   }
00108 
00109   theNDigisPerEventPlot->Fill(nDigisPerEvent);
00110 }
00111