CMS 3D CMS Logo

EcalDigisValidation.cc
Go to the documentation of this file.
1 /*
2  * \file EcalDigisValidation.cc
3  *
4  * \author F. Cossutti
5  *
6 */
7 
13 
15  HepMCToken_( consumes<edm::HepMCProduct>( edm::InputTag( ps.getParameter<std::string>( "moduleLabelMC" ) ) ) ),
16  g4TkInfoToken_( consumes<edm::SimTrackContainer>( edm::InputTag( ps.getParameter<std::string>( "moduleLabelG4" ) ) ) ),
17  g4VtxInfoToken_( consumes<edm::SimVertexContainer>( edm::InputTag( ps.getParameter<std::string>( "moduleLabelG4" ) ) ) ),
18  EBdigiCollectionToken_( consumes<EBDigiCollection>( ps.getParameter<edm::InputTag>( "EBdigiCollection" ) ) ),
19  EEdigiCollectionToken_( consumes<EEDigiCollection>( ps.getParameter<edm::InputTag>( "EEdigiCollection" ) ) ),
20  ESdigiCollectionToken_( consumes<ESDigiCollection>( ps.getParameter<edm::InputTag>( "ESdigiCollection" ) ) ),
21  crossingFramePCaloHitEBToken_( consumes< CrossingFrame<PCaloHit> >( edm::InputTag( std::string( "mix" )
22  , ps.getParameter<std::string>( "moduleLabelG4" ) + std::string( "EcalHitsEB" )
23  )
24  )
25  ),
26  crossingFramePCaloHitEEToken_( consumes< CrossingFrame<PCaloHit> >( edm::InputTag( std::string( "mix" )
27  , ps.getParameter<std::string>( "moduleLabelG4" ) + std::string( "EcalHitsEE" )
28  )
29  )
30  ),
31  crossingFramePCaloHitESToken_( consumes< CrossingFrame<PCaloHit> >( edm::InputTag( std::string( "mix" )
32  , ps.getParameter<std::string>( "moduleLabelG4" ) + std::string( "EcalHitsES" )
33  )
34  )
35  ) {
36 
37 
38  // DQM ROOT output
39  outputFile_ = ps.getUntrackedParameter<std::string>("outputFile", "");
40 
41  if ( outputFile_.size() != 0 ) {
42  edm::LogInfo("OutputInfo") << " Ecal Digi Task histograms will be saved to '" << outputFile_.c_str() << "'";
43  } else {
44  edm::LogInfo("OutputInfo") << " Ecal Digi Task histograms will NOT be saved";
45  }
46 
47  // verbosity switch
48  verbose_ = ps.getUntrackedParameter<bool>("verbose", false);
49 
50  gainConv_[1] = 1.;
51  gainConv_[2] = 2.;
52  gainConv_[3] = 12.;
53  gainConv_[0] = 12.; // saturated channels
54  barrelADCtoGeV_ = 0.035;
55  endcapADCtoGeV_ = 0.06;
56 
57  meGunEnergy_ = 0;
58  meGunEta_ = 0;
59  meGunPhi_ = 0;
60 
63 
66 
69 
70 }
71 
73 
74 }
75 
77 
79 
80 }
81 
83 
84  Char_t histo[200];
85 
86  ibooker.setCurrentFolder("EcalDigisV/EcalDigiTask");
87 
88  sprintf (histo, "EcalDigiTask Gun Momentum" ) ;
89  meGunEnergy_ = ibooker.book1D(histo, histo, 100, 0., 1000.);
90 
91  sprintf (histo, "EcalDigiTask Gun Eta" ) ;
92  meGunEta_ = ibooker.book1D(histo, histo, 700, -3.5, 3.5);
93 
94  sprintf (histo, "EcalDigiTask Gun Phi" ) ;
95  meGunPhi_ = ibooker.book1D(histo, histo, 360, 0., 360.);
96 
97  sprintf (histo, "EcalDigiTask Barrel maximum Digi over Sim ratio" ) ;
98  meEBDigiSimRatio_ = ibooker.book1D(histo, histo, 100, 0., 2.) ;
99 
100  sprintf (histo, "EcalDigiTask Endcap maximum Digi over Sim ratio" ) ;
101  meEEDigiSimRatio_ = ibooker.book1D(histo, histo, 100, 0., 2.) ;
102 
103  sprintf (histo, "EcalDigiTask Barrel maximum Digi over Sim ratio gt 10 ADC" ) ;
104  meEBDigiSimRatiogt10ADC_ = ibooker.book1D(histo, histo, 100, 0., 2.) ;
105 
106  sprintf (histo, "EcalDigiTask Endcap maximum Digi over Sim ratio gt 20 ADC" ) ;
107  meEEDigiSimRatiogt20ADC_ = ibooker.book1D(histo, histo, 100, 0., 2.) ;
108 
109  sprintf (histo, "EcalDigiTask Barrel maximum Digi over Sim ratio gt 100 ADC" ) ;
110  meEBDigiSimRatiogt100ADC_ = ibooker.book1D(histo, histo, 100, 0., 2.) ;
111 
112  sprintf (histo, "EcalDigiTask Endcap maximum Digi over Sim ratio gt 100 ADC" ) ;
113  meEEDigiSimRatiogt100ADC_ = ibooker.book1D(histo, histo, 100, 0., 2.) ;
114 
115 }
116 
118 
119  edm::LogInfo("EventInfo") << " Run = " << e.id().run() << " Event = " << e.id().event();
120 
121  std::vector<SimTrack> theSimTracks;
122  std::vector<SimVertex> theSimVertexes;
123 
127  edm::Handle<CrossingFrame<PCaloHit> > crossingFrame;
131 
132  bool skipMC = false;
133  e.getByToken( HepMCToken_, MCEvt );
134  if (!MCEvt.isValid()) { skipMC = true; }
135  e.getByToken( g4TkInfoToken_, SimTk );
136  e.getByToken( g4VtxInfoToken_, SimVtx );
137 
138  const EBDigiCollection* EBdigis =0;
139  const EEDigiCollection* EEdigis =0;
140  const ESDigiCollection* ESdigis =0;
141 
142  bool isBarrel = true;
143  e.getByToken( EBdigiCollectionToken_, EcalDigiEB );
144  if (EcalDigiEB.isValid()) {
145  EBdigis = EcalDigiEB.product();
146  LogDebug("DigiInfo") << "total # EBdigis: " << EBdigis->size() ;
147  if ( EBdigis->size() == 0 ) isBarrel = false;
148  } else {
149  isBarrel = false;
150  }
151 
152  bool isEndcap = true;
153  e.getByToken( EEdigiCollectionToken_, EcalDigiEE );
154  if (EcalDigiEE.isValid()) {
155  EEdigis = EcalDigiEE.product();
156  LogDebug("DigiInfo") << "total # EEdigis: " << EEdigis->size() ;
157  if ( EEdigis->size() == 0 ) isEndcap = false;
158  } else {
159  isEndcap = false;
160  }
161 
162  bool isPreshower = true;
163  e.getByToken( ESdigiCollectionToken_, EcalDigiES );
164  if (EcalDigiES.isValid()) {
165  ESdigis = EcalDigiES.product();
166  LogDebug("DigiInfo") << "total # ESdigis: " << ESdigis->size() ;
167  if ( ESdigis->size() == 0 ) isPreshower = false;
168  } else {
169  isPreshower = false;
170  }
171 
172  theSimTracks.insert(theSimTracks.end(),SimTk->begin(),SimTk->end());
173  theSimVertexes.insert(theSimVertexes.end(),SimVtx->begin(),SimVtx->end());
174 
175  if ( ! skipMC ) {
176  double theGunEnergy = 0.;
177  for ( HepMC::GenEvent::particle_const_iterator p = MCEvt->GetEvent()->particles_begin();
178  p != MCEvt->GetEvent()->particles_end(); ++p ) {
179 
180  theGunEnergy = (*p)->momentum().e();
181  double htheta = (*p)->momentum().theta();
182  double heta = -log(tan(htheta * 0.5));
183  double hphi = (*p)->momentum().phi();
184  hphi = (hphi>=0) ? hphi : hphi+2*M_PI;
185  hphi = hphi / M_PI * 180.;
186  LogDebug("EventInfo") << "Particle gun type form MC = " << abs((*p)->pdg_id()) << "\n" << "Energy = "<< (*p)->momentum().e() << " Eta = " << heta << " Phi = " << hphi;
187 
188  if (meGunEnergy_) meGunEnergy_->Fill(theGunEnergy);
189  if (meGunEta_) meGunEta_->Fill(heta);
190  if (meGunPhi_) meGunPhi_->Fill(hphi);
191 
192  }
193  }
194 
195  int nvtx = 0;
196  for (std::vector<SimVertex>::iterator isimvtx = theSimVertexes.begin();
197  isimvtx != theSimVertexes.end(); ++isimvtx){
198  LogDebug("EventInfo") <<" Vertex index = " << nvtx << " event Id = " << isimvtx->eventId().rawId() << "\n" << " vertex dump: " << *isimvtx ;
199  ++nvtx;
200  }
201 
202  int ntrk = 0;
203  for (std::vector<SimTrack>::iterator isimtrk = theSimTracks.begin();
204  isimtrk != theSimTracks.end(); ++isimtrk){
205  LogDebug("EventInfo") <<" Track index = " << ntrk << " track Id = " << isimtrk->trackId() << " event Id = " << isimtrk->eventId().rawId() << "\n" << " track dump: " << *isimtrk ;
206  ++ntrk;
207  }
208 
209  // BARREL
210 
211  // loop over simHits
212 
213  if ( isBarrel ) {
214 
215  e.getByToken( crossingFramePCaloHitEBToken_, crossingFrame );
216  const MixCollection<PCaloHit> barrelHits(crossingFrame.product());
217 
218  MapType ebSimMap;
219  for ( auto const & iHit : barrelHits ) {
220 
221  EBDetId ebid = EBDetId(iHit.id()) ;
222 
223  LogDebug("HitInfo")
224  << " CaloHit " << iHit.getName() << "\n"
225  << " DetID = " << iHit.id()<< " EBDetId = " << ebid.ieta() << " " << ebid.iphi() << "\n"
226  << " Time = " << iHit.time() << " Event id. = " << iHit.eventId().rawId() << "\n"
227  << " Track Id = " << iHit.geantTrackId() << "\n"
228  << " Energy = " << iHit.energy();
229 
230  uint32_t crystid = ebid.rawId();
231  ebSimMap[crystid] += iHit.energy();
232 
233  }
234 
235  // loop over Digis
236 
237  const EBDigiCollection * barrelDigi = EcalDigiEB.product () ;
238 
239  std::vector<double> ebAnalogSignal ;
240  std::vector<double> ebADCCounts ;
241  std::vector<double> ebADCGains ;
242  ebAnalogSignal.reserve(EBDataFrame::MAXSAMPLES);
243  ebADCCounts.reserve(EBDataFrame::MAXSAMPLES);
244  ebADCGains.reserve(EBDataFrame::MAXSAMPLES);
245 
246  for (unsigned int digis=0; digis<EcalDigiEB->size(); ++digis) {
247 
248  EBDataFrame ebdf=(*barrelDigi)[digis];
249  int nrSamples=ebdf.size();
250 
251  EBDetId ebid = ebdf.id () ;
252 
253  double Emax = 0. ;
254  int Pmax = 0 ;
255  double pedestalPreSample = 0.;
256  double pedestalPreSampleAnalog = 0.;
257 
258  for (int sample = 0 ; sample < nrSamples; ++sample) {
259  ebAnalogSignal[sample] = 0.;
260  ebADCCounts[sample] = 0.;
261  ebADCGains[sample] = -1.;
262  }
263 
264  for (int sample = 0 ; sample < nrSamples; ++sample) {
265 
267 
268  ebADCCounts[sample] = (mySample.adc()) ;
269  ebADCGains[sample] = (mySample.gainId()) ;
270  ebAnalogSignal[sample] = (ebADCCounts[sample]*gainConv_[(int)ebADCGains[sample]]*barrelADCtoGeV_);
271  if (Emax < ebAnalogSignal[sample] ) {
272  Emax = ebAnalogSignal[sample] ;
273  Pmax = sample ;
274  }
275  if ( sample < 3 ) {
276  pedestalPreSample += ebADCCounts[sample] ;
277  pedestalPreSampleAnalog += ebADCCounts[sample]*gainConv_[(int)ebADCGains[sample]]*barrelADCtoGeV_ ;
278  }
279  LogDebug("DigiInfo") << "EB sample " << sample << " ADC counts = " << ebADCCounts[sample] << " Gain Id = " << ebADCGains[sample] << " Analog eq = " << ebAnalogSignal[sample];
280  }
281 
282  pedestalPreSample /= 3. ;
283  pedestalPreSampleAnalog /= 3. ;
284  double Erec = Emax - pedestalPreSampleAnalog*gainConv_[(int)ebADCGains[Pmax]];
285 
286  if ( ebSimMap[ebid.rawId()] != 0. ) {
287  LogDebug("DigiInfo") << " Digi / Hit = " << Erec << " / " << ebSimMap[ebid.rawId()] << " gainConv " << gainConv_[(int)ebADCGains[Pmax]];
288  if ( meEBDigiSimRatio_ ) meEBDigiSimRatio_->Fill( Erec/ebSimMap[ebid.rawId()] ) ;
289  if ( Erec > 10.*barrelADCtoGeV_ && meEBDigiSimRatiogt10ADC_ ) meEBDigiSimRatiogt10ADC_->Fill( Erec/ebSimMap[ebid.rawId()] );
290  if ( Erec > 100.*barrelADCtoGeV_ && meEBDigiSimRatiogt100ADC_ ) meEBDigiSimRatiogt100ADC_->Fill( Erec/ebSimMap[ebid.rawId()] );
291 
292  }
293 
294  }
295 
296  }
297 
298  // ENDCAP
299 
300  // loop over simHits
301 
302  if ( isEndcap ) {
303 
304  e.getByToken( crossingFramePCaloHitEEToken_, crossingFrame );
305  const MixCollection<PCaloHit> endcapHits(crossingFrame.product());
306 
307  MapType eeSimMap;
308  for ( auto const & iHit : endcapHits ) {
309 
310  EEDetId eeid = EEDetId(iHit.id()) ;
311 
312  LogDebug("HitInfo")
313  << " CaloHit " << iHit.getName() << "\n"
314  << " DetID = "<<iHit.id()<< " EEDetId side = " << eeid.zside() << " = " << eeid.ix() << " " << eeid.iy() << "\n"
315  << " Time = " << iHit.time() << " Event id. = " << iHit.eventId().rawId() << "\n"
316  << " Track Id = " << iHit.geantTrackId() << "\n"
317  << " Energy = " << iHit.energy();
318 
319  uint32_t crystid = eeid.rawId();
320  eeSimMap[crystid] += iHit.energy();
321 
322  }
323 
324  // loop over Digis
325 
326  const EEDigiCollection * endcapDigi = EcalDigiEE.product () ;
327 
328  std::vector<double> eeAnalogSignal ;
329  std::vector<double> eeADCCounts ;
330  std::vector<double> eeADCGains ;
331  eeAnalogSignal.reserve(EEDataFrame::MAXSAMPLES);
332  eeADCCounts.reserve(EEDataFrame::MAXSAMPLES);
333  eeADCGains.reserve(EEDataFrame::MAXSAMPLES);
334 
335  for (unsigned int digis=0; digis<EcalDigiEE->size(); ++digis) {
336 
337  EEDataFrame eedf=(*endcapDigi)[digis];
338  int nrSamples=eedf.size();
339 
340  EEDetId eeid = eedf.id () ;
341 
342  double Emax = 0. ;
343  int Pmax = 0 ;
344  double pedestalPreSample = 0.;
345  double pedestalPreSampleAnalog = 0.;
346 
347  for (int sample = 0 ; sample < nrSamples; ++sample) {
348  eeAnalogSignal[sample] = 0.;
349  eeADCCounts[sample] = 0.;
350  eeADCGains[sample] = -1.;
351  }
352 
353  for (int sample = 0 ; sample < nrSamples; ++sample) {
354 
356 
357  eeADCCounts[sample] = (mySample.adc()) ;
358  eeADCGains[sample] = (mySample.gainId()) ;
359  eeAnalogSignal[sample] = (eeADCCounts[sample]*gainConv_[(int)eeADCGains[sample]]*endcapADCtoGeV_);
360  if (Emax < eeAnalogSignal[sample] ) {
361  Emax = eeAnalogSignal[sample] ;
362  Pmax = sample ;
363  }
364  if ( sample < 3 ) {
365  pedestalPreSample += eeADCCounts[sample] ;
366  pedestalPreSampleAnalog += eeADCCounts[sample]*gainConv_[(int)eeADCGains[sample]]*endcapADCtoGeV_ ;
367  }
368  LogDebug("DigiInfo") << "EE sample " << sample << " ADC counts = " << eeADCCounts[sample] << " Gain Id = " << eeADCGains[sample] << " Analog eq = " << eeAnalogSignal[sample];
369  }
370  pedestalPreSample /= 3. ;
371  pedestalPreSampleAnalog /= 3. ;
372  double Erec = Emax - pedestalPreSampleAnalog*gainConv_[(int)eeADCGains[Pmax]];
373 
374  if (eeSimMap[eeid.rawId()] != 0. ) {
375  LogDebug("DigiInfo") << " Digi / Hit = " << Erec << " / " << eeSimMap[eeid.rawId()] << " gainConv " << gainConv_[(int)eeADCGains[Pmax]];
376  if ( meEEDigiSimRatio_) meEEDigiSimRatio_->Fill( Erec/eeSimMap[eeid.rawId()] ) ;
377  if ( Erec > 20.*endcapADCtoGeV_ && meEEDigiSimRatiogt20ADC_ ) meEEDigiSimRatiogt20ADC_->Fill( Erec/eeSimMap[eeid.rawId()] );
378  if ( Erec > 100.*endcapADCtoGeV_ && meEEDigiSimRatiogt100ADC_ ) meEEDigiSimRatiogt100ADC_->Fill( Erec/eeSimMap[eeid.rawId()] );
379  }
380 
381  }
382 
383  }
384 
385  if ( isPreshower) {
386 
387  e.getByToken( crossingFramePCaloHitESToken_, crossingFrame );
388  const MixCollection<PCaloHit> preshowerHits(crossingFrame.product());
389  for ( auto const &iHit : preshowerHits ) {
390 
391  ESDetId esid = ESDetId(iHit.id()) ;
392 
393  LogDebug("HitInfo")
394  << " CaloHit " << iHit.getName() << "\n"
395  << " DetID = " << iHit.id()<< "ESDetId: z side " << esid.zside() << " plane " << esid.plane() << esid.six() << ',' << esid.siy() << ':' << esid.strip() << "\n"
396  << " Time = " << iHit.time() << " Event id. = " << iHit.eventId().rawId() << "\n"
397  << " Track Id = " << iHit.geantTrackId() << "\n"
398  << " Energy = " << iHit.energy();
399 
400  }
401 
402  }
403 
404 }
405 
407 {
408 
409  // ADC -> GeV Scale
411  eventSetup.get<EcalADCToGeVConstantRcd>().get(pAgc);
412  const EcalADCToGeVConstant* agc = pAgc.product();
413 
414  EcalMGPAGainRatio * defaultRatios = new EcalMGPAGainRatio();
415 
416  gainConv_[1] = 1.;
417  gainConv_[2] = defaultRatios->gain12Over6() ;
418  gainConv_[3] = gainConv_[2]*(defaultRatios->gain6Over1()) ;
419  gainConv_[0] = gainConv_[2]*(defaultRatios->gain6Over1()) ; // saturated channels
420 
421  LogDebug("EcalDigi") << " Gains conversions: " << "\n" << " g1 = " << gainConv_[1] << "\n" << " g2 = " << gainConv_[2] << "\n" << " g3 = " << gainConv_[3];
422  LogDebug("EcalDigi") << " Gains conversions: " << "\n" << " saturation = " << gainConv_[0];
423 
424  delete defaultRatios;
425 
426  const double barrelADCtoGeV_ = agc->getEBValue();
427  LogDebug("EcalDigi") << " Barrel GeV/ADC = " << barrelADCtoGeV_;
428  const double endcapADCtoGeV_ = agc->getEEValue();
429  LogDebug("EcalDigi") << " Endcap GeV/ADC = " << endcapADCtoGeV_;
430 
431 }
#define LogDebug(id)
RunNumber_t run() const
Definition: EventID.h:39
edm::EDGetTokenT< EEDigiCollection > EEdigiCollectionToken_
EventNumber_t event() const
Definition: EventID.h:41
T getUntrackedParameter(std::string const &, T const &) const
int strip() const
Definition: ESDetId.h:53
int ix() const
Definition: EEDetId.h:76
key_type id() const
Definition: EBDataFrame.h:31
edm::EDGetTokenT< EBDigiCollection > EBdigiCollectionToken_
std::map< uint32_t, float, std::less< uint32_t > > MapType
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
~EcalDigisValidation()
Destructor.
edm::EDGetTokenT< edm::SimVertexContainer > g4VtxInfoToken_
int six() const
Definition: ESDetId.h:49
void dqmBeginRun(edm::Run const &, edm::EventSetup const &) override
MonitorElement * meEEDigiSimRatiogt100ADC_
std::map< int, double, std::less< int > > gainConv_
int gainId() const
get the gainId (2 bits)
int size() const
Definition: EcalDataFrame.h:26
MonitorElement * meGunEnergy_
void Fill(long long x)
int iphi() const
get the crystal iphi
Definition: EBDetId.h:53
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
MonitorElement * meEBDigiSimRatiogt10ADC_
void analyze(edm::Event const &e, edm::EventSetup const &c) override
Analyze.
MonitorElement * meEBDigiSimRatio_
int siy() const
Definition: ESDetId.h:51
int zside() const
Definition: EEDetId.h:70
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< CrossingFrame< PCaloHit > > crossingFramePCaloHitESToken_
int iy() const
Definition: EEDetId.h:82
float gain6Over1() const
int ieta() const
get the crystal ieta
Definition: EBDetId.h:51
int zside() const
Definition: ESDetId.h:45
bool isValid() const
Definition: HandleBase.h:74
MonitorElement * meGunEta_
bool isEndcap(GeomDetEnumerators::SubDetector m)
key_type id() const
Definition: EEDataFrame.h:28
#define M_PI
void reserve(size_t isize)
edm::EDGetTokenT< CrossingFrame< PCaloHit > > crossingFramePCaloHitEBToken_
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
T const * product() const
Definition: Handle.h:81
const T & get() const
Definition: EventSetup.h:55
std::vector< SimVertex > SimVertexContainer
float gain12Over6() const
EcalDigisValidation(const edm::ParameterSet &ps)
Constructor.
MonitorElement * meEEDigiSimRatio_
edm::EventID id() const
Definition: EventBase.h:60
HLT enums.
MonitorElement * meGunPhi_
int plane() const
Definition: ESDetId.h:47
edm::EDGetTokenT< ESDigiCollection > ESdigiCollectionToken_
MonitorElement * meEEDigiSimRatiogt20ADC_
edm::EDGetTokenT< edm::HepMCProduct > HepMCToken_
void checkCalibrations(edm::EventSetup const &c)
edm::EDGetTokenT< CrossingFrame< PCaloHit > > crossingFramePCaloHitEEToken_
void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
T const * product() const
Definition: ESHandle.h:86
std::vector< SimTrack > SimTrackContainer
MonitorElement * meEBDigiSimRatiogt100ADC_
static const int MAXSAMPLES
Definition: EcalDataFrame.h:48
edm::EDGetTokenT< edm::SimTrackContainer > g4TkInfoToken_
Definition: Run.h:42
int adc() const
get the ADC sample (12 bits)