CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalPreshowerNoiseDistrib.cc
Go to the documentation of this file.
1 /*
2  * \file EcalPreshowerNoiseDistrib.cc
3  *
4 */
5 
7 
9  ESdigiCollectionToken_( consumes<ESDigiCollection>( ps.getParameter<edm::InputTag>( "ESdigiCollection" ) ) )
10 {
11 
12  // verbosity switch
13  verbose_ = ps.getUntrackedParameter<bool>("verbose", false);
14 
15  dbe_ = 0;
16 
17  // get hold of back-end interface
18  dbe_ = edm::Service<DQMStore>().operator->();
19 
20  if ( dbe_ ) {
21  if ( verbose_ ) {
22  dbe_->setVerbose(1);
23  } else {
24  dbe_->setVerbose(0);
25  }
26  }
27 
28  if ( dbe_ ) {
29  if ( verbose_ ) dbe_->showDirStructure();
30  }
31 
32  // histos
34  for (int ii=0; ii<3; ii++ ) { meESDigiADC_[ii] = 0; }
35 
36  Char_t histo[200];
37  if ( dbe_ ) {
38  sprintf (histo, "multiplicity" ) ;
39  meESDigiMultiplicity_ = dbe_->book1D(histo, histo, 1000, 0., 137728);
40 
41  for ( int ii = 0; ii < 3 ; ii++ ) {
42  sprintf (histo, "esRefHistos%02d", ii) ;
43  meESDigiADC_[ii] = dbe_->book1D(histo, histo, 35, 983.5, 1018.5) ;
44  }
45 
46  for ( int ii = 0; ii < 3 ; ii++ ) {
47  sprintf (histo, "esRefHistosCorr%02d", ii) ;
48  meESDigiCorr_[ii] = dbe_->book2D(histo, histo, 35, 983.5, 1018.5, 35, 983.5, 1018.5) ;
49  }
50 
51  meESDigi3D_ = dbe_->book3D("meESDigi3D_", "meESDigi3D_", 35, 983.5, 1018.5, 35, 983.5, 1018.5, 35, 983.5, 1018.5) ;
52  }
53 }
54 
55 
57 
59 
60  e.getByToken( ESdigiCollectionToken_ , EcalDigiES );
61 
62  // retrun if no data
63  if( !EcalDigiES.isValid() ) return;
64 
65  // loop over Digis
66  const ESDigiCollection * preshowerDigi = EcalDigiES.product () ;
67 
68  std::vector<double> esADCCounts ;
69  esADCCounts.reserve(ESDataFrame::MAXSAMPLES);
70 
71  int nDigis = 0;
72 
73  for (unsigned int digis=0; digis<EcalDigiES->size(); ++digis) {
74  nDigis++;
75  ESDataFrame esdf=(*preshowerDigi)[digis];
76  int nrSamples=esdf.size();
77  for (int sample = 0 ; sample < nrSamples; ++sample) {
78  ESSample mySample = esdf[sample];
79  if (meESDigiADC_[sample]) { meESDigiADC_[sample] ->Fill(mySample.adc()); }
80  }
81 
82  // to study correlations
83  if(meESDigiCorr_[0]){ meESDigiCorr_[0]->Fill(esdf[0].adc(),esdf[1].adc()); }
84  if(meESDigiCorr_[1]){ meESDigiCorr_[1]->Fill(esdf[0].adc(),esdf[2].adc()); }
85  if(meESDigiCorr_[2]){ meESDigiCorr_[2]->Fill(esdf[1].adc(),esdf[2].adc()); }
86 
87  // reference histo: sample0, sample1, sample2
88  if ( meESDigi3D_ ) meESDigi3D_ -> Fill(esdf[0].adc(),esdf[1].adc(),esdf[2].adc());
89  }
90 
92 
93 }
94 
int adc(sample_type sample)
get the ADC sample (12 bits)
T getUntrackedParameter(std::string const &, T const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
int size() const
Definition: ESDataFrame.h:23
int ii
Definition: cuy.py:588
EcalPreshowerNoiseDistrib(const edm::ParameterSet &ps)
Constructor.
void Fill(long long x)
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
static const int MAXSAMPLES
Definition: ESDataFrame.h:32
bool isValid() const
Definition: HandleBase.h:75
T const * product() const
Definition: Handle.h:81
edm::EDGetTokenT< ESDigiCollection > ESdigiCollectionToken_
void analyze(const edm::Event &e, const edm::EventSetup &c)
Analyze.
int adc() const
get the ADC sample (singed 16 bits)
Definition: ESSample.h:18