CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/Validation/MuonCSCDigis/src/CSCWireDigiValidation.cc

Go to the documentation of this file.
00001 #include "Validation/MuonCSCDigis/src/CSCWireDigiValidation.h"
00002 #include "DataFormats/Common/interface/Handle.h"
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 #include "DataFormats/CSCDigi/interface/CSCWireDigiCollection.h"
00005 
00006 #include "Geometry/CSCGeometry/interface/CSCLayerGeometry.h"
00007 #include "Geometry/CSCGeometry/interface/CSCGeometry.h"
00008 #include "DQMServices/Core/interface/DQMStore.h"
00009 
00010 
00011 
00012 CSCWireDigiValidation::CSCWireDigiValidation(DQMStore* dbe, const edm::InputTag & inputTag, bool doSim)
00013 : CSCBaseValidation(dbe, inputTag),
00014   theDoSimFlag(doSim),
00015   theTimeBinPlots(),
00016   theNDigisPerLayerPlots(),
00017   theNDigisPerEventPlot( dbe_->book1D("CSCWireDigisPerEvent", "CSC Wire Digis per event", 100, 0, 100) )
00018 {
00019   for(int i = 0; i < 10; ++i)
00020   {
00021     char title1[200], title2[200], title3[200];
00022     sprintf(title1, "CSCWireDigiTimeType%d", i+1);
00023     sprintf(title2, "CSCWireDigisPerLayerType%d", i+1);
00024     sprintf(title3, "CSCWireDigiResolution%d", i+1);
00025     theTimeBinPlots[i] = dbe_->book1D(title1, title1, 9, 0, 8);
00026     theNDigisPerLayerPlots[i] = dbe_->book1D(title2, title2, 100, 0, 20);
00027     theResolutionPlots[i] = dbe_->book1D(title3, title3, 100, -10, 10);
00028   }
00029 }
00030 
00031 
00032 
00033 CSCWireDigiValidation::~CSCWireDigiValidation()
00034 {
00035 //   for(int i = 0; i < 10; ++i)
00036 //   {
00037 //     edm::LogInfo("CSCDigiValidation") << "Mean of " << theTimeBinPlots[i]->getName() 
00038 //       << " is " << theTimeBinPlots[i]->getMean() 
00039 //       << " +/- " << theTimeBinPlots[i]->getRMS();
00040 //     edm::LogInfo("CSCDigiValidation") << "RMS of " << theResolutionPlots[i]->getName() 
00041 //       << " is " << theResolutionPlots[i]->getRMS();
00042 //   }
00043 }
00044 
00045 
00046 void CSCWireDigiValidation::analyze(const edm::Event&e, const edm::EventSetup&)
00047 {
00048   edm::Handle<CSCWireDigiCollection> wires;
00049 
00050   e.getByLabel(theInputTag, wires);
00051   if (!wires.isValid()) {
00052     edm::LogError("CSCDigiDump") << "Cannot get wires by label " << theInputTag.encode();
00053   }
00054 
00055   unsigned nDigisPerEvent = 0;
00056 
00057   for (CSCWireDigiCollection::DigiRangeIterator j=wires->begin(); j!=wires->end(); j++) {
00058     std::vector<CSCWireDigi>::const_iterator beginDigi = (*j).second.first;
00059     std::vector<CSCWireDigi>::const_iterator endDigi = (*j).second.second;
00060     int detId = (*j).first.rawId();
00061     
00062     const CSCLayer * layer = findLayer(detId);
00063     int chamberType = layer->chamber()->specs()->chamberType();
00064     int nDigis = endDigi-beginDigi;
00065     nDigisPerEvent += nDigis;
00066     theNDigisPerLayerPlots[chamberType-1]->Fill(nDigis);
00067 
00068     for( std::vector<CSCWireDigi>::const_iterator digiItr = beginDigi;
00069          digiItr != endDigi; ++digiItr) 
00070     {
00071       theTimeBinPlots[chamberType-1]->Fill(digiItr->getTimeBin());
00072     }
00073 
00074     if(theDoSimFlag) 
00075     {
00076       const edm::PSimHitContainer simHits = theSimHitMap->hits(detId);
00077       if(nDigis == 1 && simHits.size() == 1)
00078       {
00079         plotResolution(simHits[0], *beginDigi, layer, chamberType);
00080       }
00081     }
00082   }
00083 
00084   theNDigisPerEventPlot->Fill(nDigisPerEvent);
00085 }
00086 
00087 
00088 void CSCWireDigiValidation::plotResolution(const PSimHit & hit, 
00089                                            const CSCWireDigi & digi, 
00090                                            const CSCLayer * layer, 
00091                                            int chamberType)
00092 {
00093   double hitX = hit.localPosition().x();
00094   double hitY = hit.localPosition().y();
00095   double digiY = layer->geometry()->yOfWireGroup(digi.getWireGroup(), hitX);
00096   theResolutionPlots[chamberType-1]->Fill(digiY - hitY);
00097 }
00098