CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalUncalibRecHitWorkerMultiFit.cc
Go to the documentation of this file.
2 
8 
17 
18 
21  noisecorEBg12(SampleMatrix::Zero()), noisecorEEg12(SampleMatrix::Zero()),
22  noisecorEBg6(SampleMatrix::Zero()), noisecorEEg6(SampleMatrix::Zero()),
23  noisecorEBg1(SampleMatrix::Zero()), noisecorEEg1(SampleMatrix::Zero()),
24  fullpulseEB(FullSampleVector::Zero()),fullpulseEE(FullSampleVector::Zero()),
25  fullpulsecovEB(FullSampleMatrix::Zero()),fullpulsecovEE(FullSampleMatrix::Zero()) {
26 
27  // get the pulse shape, amplitude covariances and noise correlations
28  EcalPulseShapeParameters_ = ps.getParameter<edm::ParameterSet>("EcalPulseShapeParameters");
30 
31  // get the BX for the pulses to be activated
32  std::vector<int32_t> activeBXs = ps.getParameter< std::vector<int32_t> >("activeBXs");
33  activeBX.resize(activeBXs.size());
34  for (unsigned int ibx=0; ibx<activeBXs.size(); ++ibx) {
35  activeBX.coeffRef(ibx) = activeBXs[ibx];
36  }
37 
38  // uncertainty calculation (CPU intensive)
39  ampErrorCalculation_ = ps.getParameter<bool>("ampErrorCalculation");
40  useLumiInfoRunHeader_ = ps.getParameter<bool>("useLumiInfoRunHeader");
41 
42  if (useLumiInfoRunHeader_) {
43  bunchSpacing_ = c.consumes<int>(edm::InputTag("addPileupInfo","bunchSpacing"));
44  }
45 
46  doPrefitEB_ = ps.getParameter<bool>("doPrefitEB");
47  doPrefitEE_ = ps.getParameter<bool>("doPrefitEE");
48 
49  prefitMaxChiSqEB_ = ps.getParameter<double>("prefitMaxChiSqEB");
50  prefitMaxChiSqEE_ = ps.getParameter<double>("prefitMaxChiSqEE");
51 
52  // algorithm to be used for timing
53  timealgo_ = ps.getParameter<std::string>("timealgo");
54 
55  // ratio method parameters
56  EBtimeFitParameters_ = ps.getParameter<std::vector<double> >("EBtimeFitParameters");
57  EEtimeFitParameters_ = ps.getParameter<std::vector<double> >("EEtimeFitParameters");
58  EBamplitudeFitParameters_ = ps.getParameter<std::vector<double> >("EBamplitudeFitParameters");
59  EEamplitudeFitParameters_ = ps.getParameter<std::vector<double> >("EEamplitudeFitParameters");
60  EBtimeFitLimits_.first = ps.getParameter<double>("EBtimeFitLimits_Lower");
61  EBtimeFitLimits_.second = ps.getParameter<double>("EBtimeFitLimits_Upper");
62  EEtimeFitLimits_.first = ps.getParameter<double>("EEtimeFitLimits_Lower");
63  EEtimeFitLimits_.second = ps.getParameter<double>("EEtimeFitLimits_Upper");
64  EBtimeConstantTerm_=ps.getParameter<double>("EBtimeConstantTerm");
65  EEtimeConstantTerm_=ps.getParameter<double>("EEtimeConstantTerm");
66 
67  // leading edge parameters
68  ebPulseShape_ = ps.getParameter<std::vector<double> >("ebPulseShape");
69  eePulseShape_ = ps.getParameter<std::vector<double> >("eePulseShape");
70 
71  // chi2 parameters for flags determination
72  kPoorRecoFlagEB_ = ps.getParameter<bool>("kPoorRecoFlagEB");
73  kPoorRecoFlagEE_ = ps.getParameter<bool>("kPoorRecoFlagEE");;
74  chi2ThreshEB_=ps.getParameter<double>("chi2ThreshEB_");
75  chi2ThreshEE_=ps.getParameter<double>("chi2ThreshEE_");
76 
77 }
78 
79 
80 
81 // EcalUncalibRecHitWorkerMultiFit::EcalUncalibRecHitWorkerMultiFit(const edm::ParameterSet&ps) :
82 // EcalUncalibRecHitWorkerBaseClass(ps)
83 // {
84 // // ratio method parameters
85 // EBtimeFitParameters_ = ps.getParameter<std::vector<double> >("EBtimeFitParameters");
86 // EEtimeFitParameters_ = ps.getParameter<std::vector<double> >("EEtimeFitParameters");
87 // EBamplitudeFitParameters_ = ps.getParameter<std::vector<double> >("EBamplitudeFitParameters");
88 // EEamplitudeFitParameters_ = ps.getParameter<std::vector<double> >("EEamplitudeFitParameters");
89 // EBtimeFitLimits_.first = ps.getParameter<double>("EBtimeFitLimits_Lower");
90 // EBtimeFitLimits_.second = ps.getParameter<double>("EBtimeFitLimits_Upper");
91 // EEtimeFitLimits_.first = ps.getParameter<double>("EEtimeFitLimits_Lower");
92 // EEtimeFitLimits_.second = ps.getParameter<double>("EEtimeFitLimits_Upper");
93 // EBtimeConstantTerm_=ps.getParameter<double>("EBtimeConstantTerm");
94 // EEtimeConstantTerm_=ps.getParameter<double>("EEtimeConstantTerm");
95 //
96 // // leading edge parameters
97 // ebPulseShape_ = ps.getParameter<std::vector<double> >("ebPulseShape");
98 // eePulseShape_ = ps.getParameter<std::vector<double> >("eePulseShape");
99 //
100 // }
101 
102 
103 
104 
105 
106 
107 
108 
109 
110 
111 
112 
113 
114 void
116 {
117 
118  // common setup
119  es.get<EcalGainRatiosRcd>().get(gains);
120  es.get<EcalPedestalsRcd>().get(peds);
121 
122  // for the multifit method
124 
125  // weights parameters for the time
126  es.get<EcalWeightXtalGroupsRcd>().get(grps);
127  es.get<EcalTBWeightsRcd>().get(wgts);
128 
129  // which of the samples need be used
131 
132  // for the ratio method
133 
134  // for the leading edge method
137 
138  // for the time correction methods
140 }
141 
142 void
144 {
145 
146  if (useLumiInfoRunHeader_) {
147 
148  int bunchspacing = 450;
149 
150  if (evt.isRealData()) {
151  edm::RunNumber_t run = evt.run();
152  if (run == 178003 ||
153  run == 178004 ||
154  run == 209089 ||
155  run == 209106 ||
156  run == 209109 ||
157  run == 209146 ||
158  run == 209148 ||
159  run == 209151) {
160  bunchspacing = 25;
161  }
162  else {
163  bunchspacing = 50;
164  }
165  }
166  else {
167  edm::Handle<int> bunchSpacingH;
168  evt.getByToken(bunchSpacing_,bunchSpacingH);
169  bunchspacing = *bunchSpacingH;
170  }
171 
172  if (bunchspacing == 25) {
173  activeBX.resize(10);
174  activeBX << -5,-4,-3,-2,-1,0,1,2,3,4;
175  }
176  else {
177  //50ns configuration otherwise (also for no pileup)
178  activeBX.resize(5);
179  activeBX << -4,-2,0,2,4;
180  }
181  }
182 
183 }
184 
195  float ampli,
196  const std::vector<float>& amplitudeBins,
197  const std::vector<float>& shiftBins) {
198 
199  // computed initially in ns. Than turned in the BX's, as
200  // EcalUncalibratedRecHit need be.
201  double theCorrection = 0;
202 
203  // sanity check for arrays
204  if (amplitudeBins.size() == 0) {
205  edm::LogError("EcalRecHitError")
206  << "timeCorrAmplitudeBins is empty, forcing no time bias corrections.";
207 
208  return 0;
209  }
210 
211  if (amplitudeBins.size() != shiftBins.size()) {
212  edm::LogError("EcalRecHitError")
213  << "Size of timeCorrAmplitudeBins different from "
214  "timeCorrShiftBins. Forcing no time bias corrections. ";
215 
216  return 0;
217  }
218 
219  int myBin = -1;
220  for (int bin = 0; bin < (int) amplitudeBins.size(); bin++) {
221  if (ampli > amplitudeBins.at(bin)) {
222  myBin = bin;
223  } else {
224  break;
225  }
226  }
227 
228  if (myBin == -1) {
229  theCorrection = shiftBins.at(0);
230  } else if (myBin == ((int)(amplitudeBins.size() - 1))) {
231  theCorrection = shiftBins.at(myBin);
232  } else if (-1 < myBin && myBin < ((int) amplitudeBins.size() - 1)) {
233  // interpolate linearly between two assingned points
234  theCorrection = (shiftBins.at(myBin + 1) - shiftBins.at(myBin));
235  theCorrection *= (((double) ampli) - amplitudeBins.at(myBin)) /
236  (amplitudeBins.at(myBin + 1) - amplitudeBins.at(myBin));
237  theCorrection += shiftBins.at(myBin);
238  } else {
239  edm::LogError("EcalRecHitError")
240  << "Assigning time correction impossible. Setting it to 0 ";
241  theCorrection = 0.;
242  }
243 
244  // convert ns into clocks
245  return theCorrection / 25.;
246 }
247 
248 
249 
250 bool
254 {
255  DetId detid(itdg->id());
256 
257  const EcalSampleMask *sampleMask_ = sampleMaskHand_.product();
258 
259  // intelligence for recHit computation
260  EcalUncalibratedRecHit uncalibRecHit;
261 
262 
263  const EcalPedestals::Item * aped = 0;
264  const EcalMGPAGainRatio * aGain = 0;
265  const EcalXtalGroupId * gid = 0;
266 
267  if (detid.subdetId()==EcalEndcap) {
268  unsigned int hashedIndex = EEDetId(detid).hashedIndex();
269  aped = &peds->endcap(hashedIndex);
270  aGain = &gains->endcap(hashedIndex);
271  gid = &grps->endcap(hashedIndex);
274  } else {
275  unsigned int hashedIndex = EBDetId(detid).hashedIndex();
276  aped = &peds->barrel(hashedIndex);
277  aGain = &gains->barrel(hashedIndex);
278  gid = &grps->barrel(hashedIndex);
281  }
282 
283  pedVec[0] = aped->mean_x12;
284  pedVec[1] = aped->mean_x6;
285  pedVec[2] = aped->mean_x1;
286  pedRMSVec[0] = aped->rms_x12;
287  pedRMSVec[1] = aped->rms_x6;
288  pedRMSVec[2] = aped->rms_x1;
289  gainRatios[0] = 1.;
290  gainRatios[1] = aGain->gain12Over6();
291  gainRatios[2] = aGain->gain6Over1()*aGain->gain12Over6();
292 
293 
294  // === amplitude computation ===
295  int leadingSample = ((EcalDataFrame)(*itdg)).lastUnsaturatedSample();
296 
297  if ( leadingSample >= 0 ) { // saturation
298  if ( leadingSample != 4 ) {
299  // all samples different from the fifth are not reliable for the amplitude estimation
300  // put by default the energy at the saturation threshold and flag as saturated
301  float sratio = 1;
302  if ( detid.subdetId()==EcalBarrel) {
303  sratio = ebPulseShape_[5] / ebPulseShape_[4];
304  } else {
305  sratio = eePulseShape_[5] / eePulseShape_[4];
306  }
307  uncalibRecHit = EcalUncalibratedRecHit( (*itdg).id(), 4095*12*sratio, 0, 0, 0);
309  } else {
310  // float clockToNsConstant = 25.;
311  // reconstruct the rechit
312  if (detid.subdetId()==EcalEndcap) {
314  // float mult = (float)eePulseShape_.size() / (float)(*itdg).size();
315  // bin (or some analogous mapping) will be used instead of the leadingSample
316  //int bin = (int)(( (mult * leadingSample + mult/2) * clockToNsConstant + itimeconst ) / clockToNsConstant);
317  // bin is not uset for the moment
319  uncalibRecHit = leadingEdgeMethod_endcap_.makeRecHit(*itdg, pedVec, gainRatios, 0, 0);
322  } else {
324  // float mult = (float)ebPulseShape_.size() / (float)(*itdg).size();
325  // bin (or some analogous mapping) will be used instead of the leadingSample
326  //int bin = (int)(( (mult * leadingSample + mult/2) * clockToNsConstant + itimeconst ) / clockToNsConstant);
327  // bin is not uset for the moment
329  uncalibRecHit = leadingEdgeMethod_barrel_.makeRecHit(*itdg, pedVec, gainRatios, 0, 0);
332  }
333  }
334  // do not propagate the default chi2 = -1 value to the calib rechit (mapped to 64), set it to 0 when saturation
335  uncalibRecHit.setChi2(0);
336  } else {
337  // multifit
338  bool barrel = detid.subdetId()==EcalBarrel;
339  int gain = 12;
340  if (((EcalDataFrame)(*itdg)).hasSwitchToGain6()) {
341  gain = 6;
342  }
343  if (((EcalDataFrame)(*itdg)).hasSwitchToGain1()) {
344  gain = 1;
345  }
346  const SampleMatrix &noisecormat = noisecor(barrel,gain);
347  const FullSampleVector &fullpulse = barrel ? fullpulseEB : fullpulseEE;
348  const FullSampleMatrix &fullpulsecov = barrel ? fullpulsecovEB : fullpulsecovEE;
349 
350  uncalibRecHit = multiFitMethod_.makeRecHit(*itdg, aped, aGain, noisecormat,fullpulse,fullpulsecov,activeBX);
351 
352  // === time computation ===
353  if(timealgo_.compare("RatioMethod")==0) {
354  // ratio method
355  float const clockToNsConstant = 25.;
356  if (detid.subdetId()==EcalEndcap) {
357  ratioMethod_endcap_.init( *itdg, *sampleMask_, pedVec, pedRMSVec, gainRatios );
361  double theTimeCorrectionEE = timeCorrection(uncalibRecHit.amplitude(),
362  timeCorrBias_->EETimeCorrAmplitudeBins, timeCorrBias_->EETimeCorrShiftBins);
363 
364  uncalibRecHit.setJitter( crh.timeMax - 5 + theTimeCorrectionEE);
365  uncalibRecHit.setJitterError( std::sqrt(pow(crh.timeError,2) + std::pow(EEtimeConstantTerm_,2)/std::pow(clockToNsConstant,2)) );
366 
367  } else {
368  ratioMethod_barrel_.init( *itdg, *sampleMask_, pedVec, pedRMSVec, gainRatios );
373 
374  double theTimeCorrectionEB = timeCorrection(uncalibRecHit.amplitude(),
375  timeCorrBias_->EBTimeCorrAmplitudeBins, timeCorrBias_->EBTimeCorrShiftBins);
376 
377  uncalibRecHit.setJitter( crh.timeMax - 5 + theTimeCorrectionEB);
378 
379  uncalibRecHit.setJitterError( std::sqrt(std::pow(crh.timeError,2) + std::pow(EBtimeConstantTerm_,2)/std::pow(clockToNsConstant,2)) );
380  }
381  } else if(timealgo_.compare("WeightsMethod")==0) {
382  // weights method on the PU subtracted pulse shape
383  std::vector<double> amplitudes;
384  for(unsigned int ibx=0; ibx<activeBX.size(); ++ibx) amplitudes.push_back(uncalibRecHit.outOfTimeAmplitude(ibx));
385 
386  EcalTBWeights::EcalTDCId tdcid(1);
387  EcalTBWeights::EcalTBWeightMap const & wgtsMap = wgts->getMap();
388  EcalTBWeights::EcalTBWeightMap::const_iterator wit;
389  wit = wgtsMap.find( std::make_pair(*gid,tdcid) );
390  if( wit == wgtsMap.end() ) {
391  edm::LogError("EcalUncalibRecHitError") << "No weights found for EcalGroupId: "
392  << gid->id() << " and EcalTDCId: " << tdcid
393  << "\n skipping digi with id: " << detid.rawId();
394 
395  return false;
396  }
397  const EcalWeightSet& wset = wit->second; // this is the EcalWeightSet
398 
401 
402  weights[0] = &mat1;
403  weights[1] = &mat2;
404 
405  double timerh;
406  if (detid.subdetId()==EcalEndcap) {
407  timerh = weightsMethod_endcap_.time( *itdg, amplitudes, aped, aGain, fullpulse, weights);
408  } else {
409  timerh = weightsMethod_barrel_.time( *itdg, amplitudes, aped, aGain, fullpulse, weights);
410  }
411  uncalibRecHit.setJitter( timerh );
412  uncalibRecHit.setJitterError( 0. ); // not computed with weights
413  } else if(timealgo_.compare("None")==0) {
414  uncalibRecHit.setJitter( 0. );
415  uncalibRecHit.setJitterError( 0. );
416  } else {
417  edm::LogError("EcalUncalibRecHitError") << "No time estimation algorithm called "
418  << timealgo_
419  << "\n setting jitter to 0. and jitter uncertainty to 10000. ";
420 
421  uncalibRecHit.setJitter( 0. );
422  uncalibRecHit.setJitterError( 10000. );
423  }
424  }
425 
426  // set flags if gain switch has occurred
427  if( ((EcalDataFrame)(*itdg)).hasSwitchToGain6() ) uncalibRecHit.setFlagBit( EcalUncalibratedRecHit::kHasSwitchToGain6 );
428  if( ((EcalDataFrame)(*itdg)).hasSwitchToGain1() ) uncalibRecHit.setFlagBit( EcalUncalibratedRecHit::kHasSwitchToGain1 );
429 
430  // set quality flags based on chi2 of the the fit
431  /*
432  if(detid.subdetId()==EcalEndcap) {
433  if(kPoorRecoFlagEE_ && uncalibRecHit.chi2()>chi2ThreshEE_) {
434  bool samplesok = true;
435  for (int sample =0; sample < EcalDataFrame::MAXSAMPLES; ++sample) {
436  if (!sampleMask_->useSampleEE(sample)) {
437  samplesok = false;
438  break;
439  }
440  }
441  if (samplesok) uncalibRecHit.setFlagBit(EcalUncalibratedRecHit::kPoorReco);
442  }
443  } else {
444  if(kPoorRecoFlagEB_ && uncalibRecHit.chi2()>chi2ThreshEB_) {
445  bool samplesok = true;
446  for (int sample =0; sample < EcalDataFrame::MAXSAMPLES; ++sample) {
447  if (!sampleMask_->useSampleEB(sample)) {
448  samplesok = false;
449  break;
450  }
451  }
452  if (samplesok) uncalibRecHit.setFlagBit(EcalUncalibratedRecHit::kPoorReco);
453  }
454  }
455  */
456 
457  // put the recHit in the collection
458  if (detid.subdetId()==EcalEndcap) {
459  result.push_back( uncalibRecHit );
460  } else {
461  result.push_back( uncalibRecHit );
462  }
463 
464  return true;
465 }
466 
467 
469  if (barrel) {
470  if (gain==6) {
471  return noisecorEBg6;
472  }
473  else if (gain==1) {
474  return noisecorEBg1;
475  }
476  else {
477  return noisecorEBg12;
478  }
479  }
480  else {
481  if (gain==6) {
482  return noisecorEEg6;
483  }
484  else if (gain==1) {
485  return noisecorEEg1;
486  }
487  else {
488  return noisecorEEg12;
489  }
490  }
491 
492  return noisecorEBg12;
493 
494 }
495 
497 
498  const std::vector<double> ebCorMatG12 = params.getParameter< std::vector<double> >("EBCorrNoiseMatrixG12");
499  const std::vector<double> eeCorMatG12 = params.getParameter< std::vector<double> >("EECorrNoiseMatrixG12");
500  const std::vector<double> ebCorMatG06 = params.getParameter< std::vector<double> >("EBCorrNoiseMatrixG06");
501  const std::vector<double> eeCorMatG06 = params.getParameter< std::vector<double> >("EECorrNoiseMatrixG06");
502  const std::vector<double> ebCorMatG01 = params.getParameter< std::vector<double> >("EBCorrNoiseMatrixG01");
503  const std::vector<double> eeCorMatG01 = params.getParameter< std::vector<double> >("EECorrNoiseMatrixG01");
504 
505  int nnoise = ebCorMatG12.size();
506 
507  // fill correlation matrices: noise (HF (+) LF)
508  for (int i=0; i<nnoise; ++i) {
509  for (int j=0; j<nnoise; ++j) {
510  int vidx = std::abs(j-i);
511  noisecorEBg12(i,j) = ebCorMatG12[vidx];
512  noisecorEEg12(i,j) = eeCorMatG12[vidx];
513  noisecorEBg6(i,j) = ebCorMatG06[vidx];
514  noisecorEEg6(i,j) = eeCorMatG06[vidx];
515  noisecorEBg1(i,j) = ebCorMatG01[vidx];
516  noisecorEEg1(i,j) = eeCorMatG01[vidx];
517  }
518  }
519 
520  // fill shape: from simulation for samples 3-9, from alpha/beta shape for 10-14
521  const std::vector<double> ebPulse = params.getParameter< std::vector<double> >("EBPulseShapeTemplate");
522  const std::vector<double> eePulse = params.getParameter< std::vector<double> >("EEPulseShapeTemplate");
523  int nShapeSamples = ebPulse.size();
524  for (int i=0; i<nShapeSamples; ++i) {
525  fullpulseEB(i+7) = ebPulse[i];
526  fullpulseEE(i+7) = eePulse[i];
527  }
528 
529  const std::vector<double> ebPulseCov = params.getParameter< std::vector<double> >("EBPulseShapeCovariance");
530  const std::vector<double> eePulseCov = params.getParameter< std::vector<double> >("EEPulseShapeCovariance");
531  for(int k=0; k<std::pow(nShapeSamples,2); ++k) {
532  int i = k/nShapeSamples;
533  int j = k%nShapeSamples;
534  fullpulsecovEB(i+7,j+7) = ebPulseCov[k];
535  fullpulsecovEE(i+7,j+7) = eePulseCov[k];
536  }
537 
538 }
539 
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
std::map< std::pair< EcalXtalGroupId, EcalTDCId >, EcalWeightSet > EcalTBWeightMap
Definition: EcalTBWeights.h:21
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
int hashedIndex() const
get a compact index for arrays
Definition: EBDetId.h:86
void computeAmplitude(std::vector< double > &amplitudeFitParameters)
EcalUncalibRecHitWorkerMultiFit(const edm::ParameterSet &, edm::ConsumesCollector &c)
unsigned size(int bx) const
EcalUncalibRecHitMultiFitAlgo multiFitMethod_
boost::transform_iterator< IterHelp, boost::counting_iterator< int > > const_iterator
const SampleMatrix & noisecor(bool barrel, int gain) const
void setJitter(float jitter)
EcalUncalibRecHitTimeWeightsAlgo< EEDataFrame > weightsMethod_endcap_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
EcalUncalibRecHitLeadingEdgeAlgo< EBDataFrame > leadingEdgeMethod_barrel_
double time(const C &dataFrame, const std::vector< double > &amplitudes, const EcalPedestals::Item *aped, const EcalMGPAGainRatio *aGain, const FullSampleVector &fullpulse, const EcalWeightSet::EcalWeightMatrix **weights)
Compute time.
EcalUncalibratedRecHit makeRecHit(const EcalDataFrame &dataFrame, const EcalPedestals::Item *aped, const EcalMGPAGainRatio *aGain, const SampleMatrix &noisecor, const FullSampleVector &fullpulse, const FullSampleMatrix &fullpulsecov, const BXVector &activeBX)
compute rechits
void computeTime(std::vector< double > &timeFitParameters, std::pair< double, double > &timeFitLimits, std::vector< double > &amplitudeFitParameters)
edm::ESHandle< EcalWeightXtalGroups > grps
void push_back(T const &t)
void setJitterError(float jitterErr)
Eigen::Matrix< double, 19, 19 > FullSampleMatrix
bool isRealData() const
Definition: EventBase.h:60
EcalUncalibRecHitRatioMethodAlgo< EBDataFrame > ratioMethod_barrel_
void setPulseShape(std::vector< double > &shape)
const unsigned int id() const
bool run(const edm::Event &evt, const EcalDigiCollection::const_iterator &digi, EcalUncalibratedRecHitCollection &result) override
int hashedIndex(int ieta, int iphi)
Definition: EcalPyUtils.cc:42
float outOfTimeAmplitude(int bx) const
EcalUncalibRecHitRatioMethodAlgo< EEDataFrame > ratioMethod_endcap_
T sqrt(T t)
Definition: SSEVec.h:48
tuple result
Definition: query.py:137
RunNumber_t run() const
Definition: Event.h:85
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
edm::ESHandle< EcalTimeCalibConstants > itime
double timeCorrection(float ampli, const std::vector< float > &amplitudeBins, const std::vector< float > &shiftBins)
float gain6Over1() const
Eigen::Matrix< double, 19, 1 > FullSampleVector
EcalWeightMatrix & getWeightsAfterGainSwitch()
Definition: EcalWeightSet.h:30
virtual EcalUncalibratedRecHit makeRecHit(const C &dataFrame, const double *pedestals, const double *gainRatios, const EcalWeightSet::EcalWeightMatrix **weights, const EcalWeightSet::EcalChi2WeightMatrix **chi2Matrix)
Compute parameters.
edm::ESHandle< EcalTimeBiasCorrections > timeCorrBias_
Definition: DetId.h:18
EcalUncalibRecHitLeadingEdgeAlgo< EEDataFrame > leadingEdgeMethod_endcap_
void set(const edm::EventSetup &es) override
int hashedIndex() const
Definition: EEDetId.h:182
EcalWeightMatrix & getWeightsBeforeGainSwitch()
Definition: EcalWeightSet.h:29
const T & get() const
Definition: EventSetup.h:55
void resize(int bx, unsigned size)
T const * product() const
Definition: ESHandle.h:86
float gain12Over6() const
EcalUncalibRecHitTimeWeightsAlgo< EBDataFrame > weightsMethod_barrel_
edm::ESHandle< EcalTimeOffsetConstant > offtime
math::Matrix< 3, 10 >::type EcalWeightMatrix
Definition: EcalWeightSet.h:22
#define DEFINE_EDM_PLUGIN(factory, type, name)
unsigned int RunNumber_t
edm::ESHandle< EcalSampleMask > sampleMaskHand_
const EcalWeightSet::EcalWeightMatrix * weights[2]
Eigen::Matrix< double, 10, 10 > SampleMatrix
void init(const C &dataFrame, const EcalSampleMask &sampleMask, const double *pedestals, const double *pedestalRMSes, const double *gainRatios)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void fillInputs(const edm::ParameterSet &params)