CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalPreshowerDigisValidation.cc
Go to the documentation of this file.
1 /*
2  * \file EcalPreshowerDigisValidation.cc
3  *
4  * \author F. Cossutti
5  *
6 */
7 
9 
11  ESdigiCollectionToken_( consumes<ESDigiCollection>( ps.getParameter<edm::InputTag>( "ESdigiCollection" ) ) )
12 {
13 
14  // verbosity switch
15  verbose_ = ps.getUntrackedParameter<bool>("verbose", false);
16 
18 
19  for (int i = 0; i < 3 ; i++ ) {
20  meESDigiADC_[i] = 0;
21  }
22 
23 }
24 
26 
27  Char_t histo[200];
28 
29  ibooker.setCurrentFolder("EcalDigisV/EcalDigiTask");
30 
31  sprintf (histo, "EcalDigiTask Preshower digis multiplicity" ) ;
32  meESDigiMultiplicity_ = ibooker.book1D(histo, histo, 1000, 0., 137728);
33 
34  for ( int i = 0; i < 3 ; i++ ) {
35 
36  sprintf (histo, "EcalDigiTask Preshower ADC pulse %02d", i+1) ;
37  meESDigiADC_[i] = ibooker.book1D(histo, histo, 4096, -0.5, 4095.5) ;
38  }
39 }
40 
42 
43  //LogInfo("EventInfo") << " Run = " << e.id().run() << " Event = " << e.id().event();
44 
46 
47  e.getByToken( ESdigiCollectionToken_ , EcalDigiES );
48 
49  // Return if no preshower data
50  if( !EcalDigiES.isValid() ) return;
51 
52  // PRESHOWER
53 
54  // loop over Digis
55 
56  const ESDigiCollection * preshowerDigi = EcalDigiES.product () ;
57 
58  std::vector<double> esADCCounts ;
59  esADCCounts.reserve(ESDataFrame::MAXSAMPLES);
60 
61  int nDigis = 0;
62 
63  for (unsigned int digis=0; digis<EcalDigiES->size(); ++digis) {
64 
65  ESDataFrame esdf=(*preshowerDigi)[digis];
66  int nrSamples=esdf.size();
67 
68  ESDetId esid = esdf.id () ;
69 
70  nDigis++;
71 
72  for (int sample = 0 ; sample < nrSamples; ++sample) {
73  esADCCounts[sample] = 0.;
74  }
75 
76  for (int sample = 0 ; sample < nrSamples; ++sample) {
77  ESSample mySample = esdf[sample];
78  esADCCounts[sample] = (mySample.adc()) ;
79  }
80  if (verbose_) {
81  LogDebug("DigiInfo") << "Preshower Digi for ESDetId: z side " << esid.zside() << " plane " << esid.plane() << esid.six() << ',' << esid.siy() << ':' << esid.strip();
82  for ( int i = 0; i < 3 ; i++ ) {
83  LogDebug("DigiInfo") << "sample " << i << " ADC = " << esADCCounts[i];
84  }
85  }
86 
87  for ( int i = 0 ; i < 3 ; i++ ) {
88  if (meESDigiADC_[i]) meESDigiADC_[i]->Fill( esADCCounts[i] ) ;
89  }
90 
91  }
92 
94 
95 }
96 
97 
#define LogDebug(id)
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
int strip() const
Definition: ESDetId.h:52
void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
const ESDetId & id() const
Definition: ESDataFrame.h:21
int six() const
Definition: ESDetId.h:48
int size() const
Definition: ESDataFrame.h:23
void Fill(long long x)
void analyze(const edm::Event &e, const edm::EventSetup &c)
Analyze.
int siy() const
Definition: ESDetId.h:50
static const int MAXSAMPLES
Definition: ESDataFrame.h:32
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
int zside() const
Definition: ESDetId.h:44
bool isValid() const
Definition: HandleBase.h:75
edm::EDGetTokenT< ESDigiCollection > ESdigiCollectionToken_
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
T const * product() const
Definition: Handle.h:81
EcalPreshowerDigisValidation(const edm::ParameterSet &ps)
Constructor.
int plane() const
Definition: ESDetId.h:46
int adc() const
get the ADC sample (singed 16 bits)
Definition: ESSample.h:18
Definition: Run.h:43