CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/Validation/EcalDigis/src/EcalPreshowerNoiseDistrib.cc

Go to the documentation of this file.
00001 /*
00002  * \file EcalPreshowerNoiseDistrib.cc
00003  *
00004 */
00005 
00006 #include <Validation/EcalDigis/interface/EcalPreshowerNoiseDistrib.h>
00007 #include "DQMServices/Core/interface/DQMStore.h"
00008 using namespace cms;
00009 using namespace edm;
00010 using namespace std;
00011 
00012 EcalPreshowerNoiseDistrib::EcalPreshowerNoiseDistrib(const ParameterSet& ps):
00013   ESdigiCollection_(ps.getParameter<edm::InputTag>("ESdigiCollection"))
00014 {
00015   
00016   // verbosity switch
00017   verbose_ = ps.getUntrackedParameter<bool>("verbose", false);
00018     
00019   dbe_ = 0;
00020   
00021   // get hold of back-end interface
00022   dbe_ = Service<DQMStore>().operator->();
00023   
00024   if ( dbe_ ) {
00025     if ( verbose_ ) {
00026       dbe_->setVerbose(1);
00027     } else {
00028       dbe_->setVerbose(0);
00029     }
00030   }
00031   
00032   if ( dbe_ ) {
00033     if ( verbose_ ) dbe_->showDirStructure();
00034   }
00035   
00036   // histos
00037   meESDigiMultiplicity_=0;
00038   for (int ii=0; ii<3; ii++ ) { meESDigiADC_[ii] = 0; }
00039 
00040   Char_t histo[200];
00041   if ( dbe_ ) {
00042     sprintf (histo, "multiplicity" ) ;
00043     meESDigiMultiplicity_ = dbe_->book1D(histo, histo, 1000, 0., 137728);
00044     
00045     for ( int ii = 0; ii < 3 ; ii++ ) {
00046       sprintf (histo, "esRefHistos%02d", ii) ;
00047       meESDigiADC_[ii] = dbe_->book1D(histo, histo, 35, 983.5, 1018.5) ;
00048     }
00049     
00050     for ( int ii = 0; ii < 3 ; ii++ ) {
00051       sprintf (histo, "esRefHistosCorr%02d", ii) ;
00052       meESDigiCorr_[ii] = dbe_->book2D(histo, histo, 35, 983.5, 1018.5, 35, 983.5, 1018.5) ;
00053     }
00054         
00055     meESDigi3D_ = dbe_->book3D("meESDigi3D_", "meESDigi3D_", 35, 983.5, 1018.5, 35, 983.5, 1018.5, 35, 983.5, 1018.5) ;
00056   }
00057 }
00058 
00059 
00060 void EcalPreshowerNoiseDistrib::analyze(const Event& e, const EventSetup& c){
00061 
00062   Handle<ESDigiCollection> EcalDigiES;
00063   
00064   e.getByLabel( ESdigiCollection_ , EcalDigiES );
00065 
00066   // retrun if no data
00067   if( !EcalDigiES.isValid() ) return;
00068   
00069   // loop over Digis
00070   const ESDigiCollection * preshowerDigi = EcalDigiES.product () ;
00071   
00072   std::vector<double> esADCCounts ;
00073   esADCCounts.reserve(ESDataFrame::MAXSAMPLES);
00074   
00075   int nDigis = 0;
00076 
00077   for (unsigned int digis=0; digis<EcalDigiES->size(); ++digis) {
00078     nDigis++;
00079     ESDataFrame esdf=(*preshowerDigi)[digis];
00080     int nrSamples=esdf.size();
00081     for (int sample = 0 ; sample < nrSamples; ++sample) {
00082       ESSample mySample = esdf[sample];
00083       if (meESDigiADC_[sample]) { meESDigiADC_[sample] ->Fill(mySample.adc()); }
00084     }
00085 
00086     // to study correlations
00087     if(meESDigiCorr_[0]){ meESDigiCorr_[0]->Fill(esdf[0].adc(),esdf[1].adc()); } 
00088     if(meESDigiCorr_[1]){ meESDigiCorr_[1]->Fill(esdf[0].adc(),esdf[2].adc()); } 
00089     if(meESDigiCorr_[2]){ meESDigiCorr_[2]->Fill(esdf[1].adc(),esdf[2].adc()); } 
00090 
00091     // reference histo: sample0, sample1, sample2
00092     if ( meESDigi3D_ ) meESDigi3D_ -> Fill(esdf[0].adc(),esdf[1].adc(),esdf[2].adc());
00093   }
00094   
00095   if ( meESDigiMultiplicity_ ) meESDigiMultiplicity_->Fill(nDigis);
00096   
00097 }
00098