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 
20 
24 
27 {
28 
29  // get the BX for the pulses to be activated
30  std::vector<int32_t> activeBXs = ps.getParameter< std::vector<int32_t> >("activeBXs");
31  activeBX.resize(activeBXs.size());
32  for (unsigned int ibx=0; ibx<activeBXs.size(); ++ibx) {
33  activeBX.coeffRef(ibx) = activeBXs[ibx];
34  }
35 
36  // uncertainty calculation (CPU intensive)
37  ampErrorCalculation_ = ps.getParameter<bool>("ampErrorCalculation");
38  useLumiInfoRunHeader_ = ps.getParameter<bool>("useLumiInfoRunHeader");
39 
40  if (useLumiInfoRunHeader_) {
41  bunchSpacing_ = c.consumes<unsigned int>(edm::InputTag("bunchSpacingProducer"));
43  } else {
44  bunchSpacingManual_ = ps.getParameter<int>("bunchSpacing");
45  }
46 
47  doPrefitEB_ = ps.getParameter<bool>("doPrefitEB");
48  doPrefitEE_ = ps.getParameter<bool>("doPrefitEE");
49 
50  prefitMaxChiSqEB_ = ps.getParameter<double>("prefitMaxChiSqEB");
51  prefitMaxChiSqEE_ = ps.getParameter<double>("prefitMaxChiSqEE");
52 
53  dynamicPedestalsEB_ = ps.getParameter<bool>("dynamicPedestalsEB");
54  dynamicPedestalsEE_ = ps.getParameter<bool>("dynamicPedestalsEE");
55  mitigateBadSamplesEB_ = ps.getParameter<bool>("mitigateBadSamplesEB");
56  mitigateBadSamplesEE_ = ps.getParameter<bool>("mitigateBadSamplesEE");
57  gainSwitchUseMaxSampleEB_ = ps.getParameter<bool>("gainSwitchUseMaxSampleEB");
58  gainSwitchUseMaxSampleEE_ = ps.getParameter<bool>("gainSwitchUseMaxSampleEE");
59  selectiveBadSampleCriteriaEB_ = ps.getParameter<bool>("selectiveBadSampleCriteriaEB");
60  selectiveBadSampleCriteriaEE_ = ps.getParameter<bool>("selectiveBadSampleCriteriaEE");
61  addPedestalUncertaintyEB_ = ps.getParameter<double>("addPedestalUncertaintyEB");
62  addPedestalUncertaintyEE_ = ps.getParameter<double>("addPedestalUncertaintyEE");
63  simplifiedNoiseModelForGainSwitch_ = ps.getParameter<bool>("simplifiedNoiseModelForGainSwitch");
64 
65  // algorithm to be used for timing
66  auto const & timeAlgoName = ps.getParameter<std::string>("timealgo");
67  if(timeAlgoName=="RatioMethod") timealgo_=ratioMethod;
68  else if(timeAlgoName=="WeightsMethod") timealgo_=weightsMethod;
69  else if(timeAlgoName!="None")
70  edm::LogError("EcalUncalibRecHitError") << "No time estimation algorithm defined";
71 
72  // ratio method parameters
73  EBtimeFitParameters_ = ps.getParameter<std::vector<double> >("EBtimeFitParameters");
74  EEtimeFitParameters_ = ps.getParameter<std::vector<double> >("EEtimeFitParameters");
75  EBamplitudeFitParameters_ = ps.getParameter<std::vector<double> >("EBamplitudeFitParameters");
76  EEamplitudeFitParameters_ = ps.getParameter<std::vector<double> >("EEamplitudeFitParameters");
77  EBtimeFitLimits_.first = ps.getParameter<double>("EBtimeFitLimits_Lower");
78  EBtimeFitLimits_.second = ps.getParameter<double>("EBtimeFitLimits_Upper");
79  EEtimeFitLimits_.first = ps.getParameter<double>("EEtimeFitLimits_Lower");
80  EEtimeFitLimits_.second = ps.getParameter<double>("EEtimeFitLimits_Upper");
81  EBtimeConstantTerm_=ps.getParameter<double>("EBtimeConstantTerm");
82  EEtimeConstantTerm_=ps.getParameter<double>("EEtimeConstantTerm");
83  EBtimeNconst_=ps.getParameter<double>("EBtimeNconst");
84  EEtimeNconst_=ps.getParameter<double>("EEtimeNconst");
85  outOfTimeThreshG12pEB_ = ps.getParameter<double>("outOfTimeThresholdGain12pEB");
86  outOfTimeThreshG12mEB_ = ps.getParameter<double>("outOfTimeThresholdGain12mEB");
87  outOfTimeThreshG61pEB_ = ps.getParameter<double>("outOfTimeThresholdGain61pEB");
88  outOfTimeThreshG61mEB_ = ps.getParameter<double>("outOfTimeThresholdGain61mEB");
89  outOfTimeThreshG12pEE_ = ps.getParameter<double>("outOfTimeThresholdGain12pEE");
90  outOfTimeThreshG12mEE_ = ps.getParameter<double>("outOfTimeThresholdGain12mEE");
91  outOfTimeThreshG61pEE_ = ps.getParameter<double>("outOfTimeThresholdGain61pEE");
92  outOfTimeThreshG61mEE_ = ps.getParameter<double>("outOfTimeThresholdGain61mEE");
93  amplitudeThreshEB_ = ps.getParameter<double>("amplitudeThresholdEB");
94  amplitudeThreshEE_ = ps.getParameter<double>("amplitudeThresholdEE");
95 
96  // spike threshold
97  ebSpikeThresh_ = ps.getParameter<double>("ebSpikeThreshold");
98 
99  ebPulseShape_ = ps.getParameter<std::vector<double> >("ebPulseShape");
100  eePulseShape_ = ps.getParameter<std::vector<double> >("eePulseShape");
101 
102  // chi2 parameters for flags determination
103  kPoorRecoFlagEB_ = ps.getParameter<bool>("kPoorRecoFlagEB");
104  kPoorRecoFlagEE_ = ps.getParameter<bool>("kPoorRecoFlagEE");;
105  chi2ThreshEB_=ps.getParameter<double>("chi2ThreshEB_");
106  chi2ThreshEE_=ps.getParameter<double>("chi2ThreshEE_");
107 
108 }
109 
110 
111 
112 void
114 {
115 
116  // common setup
117  es.get<EcalGainRatiosRcd>().get(gains);
118  es.get<EcalPedestalsRcd>().get(peds);
119 
120  // for the multifit method
125 
126  // weights parameters for the time
127  es.get<EcalWeightXtalGroupsRcd>().get(grps);
128  es.get<EcalTBWeightsRcd>().get(wgts);
129 
130  // which of the samples need be used
132 
133  // for the ratio method
136 
137  // for the time correction methods
139 
140  int nnoise = SampleVector::RowsAtCompileTime;
141  SampleMatrix &noisecorEBg12 = noisecors_[1][0];
142  SampleMatrix &noisecorEBg6 = noisecors_[1][1];
143  SampleMatrix &noisecorEBg1 = noisecors_[1][2];
144  SampleMatrix &noisecorEEg12 = noisecors_[0][0];
145  SampleMatrix &noisecorEEg6 = noisecors_[0][1];
146  SampleMatrix &noisecorEEg1 = noisecors_[0][2];
147 
148  for (int i=0; i<nnoise; ++i) {
149  for (int j=0; j<nnoise; ++j) {
150  int vidx = std::abs(j-i);
151  noisecorEBg12(i,j) = noisecovariances->EBG12SamplesCorrelation[vidx];
152  noisecorEEg12(i,j) = noisecovariances->EEG12SamplesCorrelation[vidx];
153  noisecorEBg6(i,j) = noisecovariances->EBG6SamplesCorrelation[vidx];
154  noisecorEEg6(i,j) = noisecovariances->EEG6SamplesCorrelation[vidx];
155  noisecorEBg1(i,j) = noisecovariances->EBG1SamplesCorrelation[vidx];
156  noisecorEEg1(i,j) = noisecovariances->EEG1SamplesCorrelation[vidx];
157  }
158  }
159 }
160 
161 void
163 {
164 
165  unsigned int bunchspacing = 450;
166 
167  if (useLumiInfoRunHeader_) {
168 
169  edm::Handle<unsigned int> bunchSpacingH;
170  evt.getByToken(bunchSpacing_,bunchSpacingH);
171  bunchspacing = *bunchSpacingH;
172  }
173  else {
174  bunchspacing = bunchSpacingManual_;
175  }
176 
178  if (bunchspacing == 25) {
179  activeBX.resize(10);
180  activeBX << -5,-4,-3,-2,-1,0,1,2,3,4;
181  }
182  else {
183  //50ns configuration otherwise (also for no pileup)
184  activeBX.resize(5);
185  activeBX << -4,-2,0,2,4;
186  }
187  }
188 
189 }
190 
201  float ampli,
202  const std::vector<float>& amplitudeBins,
203  const std::vector<float>& shiftBins) {
204 
205  // computed initially in ns. Than turned in the BX's, as
206  // EcalUncalibratedRecHit need be.
207  double theCorrection = 0;
208 
209  // sanity check for arrays
210  if (amplitudeBins.size() == 0) {
211  edm::LogError("EcalRecHitError")
212  << "timeCorrAmplitudeBins is empty, forcing no time bias corrections.";
213 
214  return 0;
215  }
216 
217  if (amplitudeBins.size() != shiftBins.size()) {
218  edm::LogError("EcalRecHitError")
219  << "Size of timeCorrAmplitudeBins different from "
220  "timeCorrShiftBins. Forcing no time bias corrections. ";
221 
222  return 0;
223  }
224 
225  // FIXME? what about a binary search?
226  int myBin = -1;
227  for (int bin = 0; bin < (int) amplitudeBins.size(); bin++) {
228  if (ampli > amplitudeBins[bin]) {
229  myBin = bin;
230  } else {
231  break;
232  }
233  }
234 
235  if (myBin == -1) {
236  theCorrection = shiftBins[0];
237  } else if (myBin == ((int)(amplitudeBins.size() - 1))) {
238  theCorrection = shiftBins[myBin];
239  } else {
240  // interpolate linearly between two assingned points
241  theCorrection = (shiftBins[myBin + 1] - shiftBins[myBin]);
242  theCorrection *= (((double) ampli) - amplitudeBins[myBin]) /
243  (amplitudeBins[myBin + 1] - amplitudeBins[myBin]);
244  theCorrection += shiftBins[myBin];
245  }
246 
247  // convert ns into clocks
248  constexpr double inv25 = 1./25.;
249  return theCorrection * inv25;
250 }
251 
252 
253 
254 bool
258 {
259  DetId detid(itdg->id());
260 
261  const EcalSampleMask *sampleMask_ = sampleMaskHand_.product();
262 
263  // intelligence for recHit computation
264  EcalUncalibratedRecHit uncalibRecHit;
265  float offsetTime = 0;
266 
267  const EcalPedestals::Item * aped = 0;
268  const EcalMGPAGainRatio * aGain = 0;
269  const EcalXtalGroupId * gid = 0;
270  const EcalPulseShapes::Item * aPulse = 0;
271  const EcalPulseCovariances::Item * aPulseCov = 0;
272 
274 
275  if (detid.subdetId()==EcalEndcap) {
276  unsigned int hashedIndex = EEDetId(detid).hashedIndex();
277  aped = &peds->endcap(hashedIndex);
278  aGain = &gains->endcap(hashedIndex);
279  gid = &grps->endcap(hashedIndex);
280  aPulse = &pulseshapes->endcap(hashedIndex);
281  aPulseCov = &pulsecovariances->endcap(hashedIndex);
289  offsetTime = offtime->getEEValue();
290  } else {
291  unsigned int hashedIndex = EBDetId(detid).hashedIndex();
292  aped = &peds->barrel(hashedIndex);
293  aGain = &gains->barrel(hashedIndex);
294  gid = &grps->barrel(hashedIndex);
295  aPulse = &pulseshapes->barrel(hashedIndex);
296  aPulseCov = &pulsecovariances->barrel(hashedIndex);
304  offsetTime = offtime->getEBValue();
305  }
306 
307  double pedVec[3] = { aped->mean_x12, aped->mean_x6, aped->mean_x1 };
308  double pedRMSVec[3] = { aped->rms_x12, aped->rms_x6, aped->rms_x1};
309  double gainRatios[3] = { 1., aGain->gain12Over6(), aGain->gain6Over1()*aGain->gain12Over6()};
310 
311  FullSampleVector fullpulse(FullSampleVector::Zero());
312  FullSampleMatrix fullpulsecov(FullSampleMatrix::Zero());
313  for (int i=0; i<EcalPulseShape::TEMPLATESAMPLES; ++i)
314  fullpulse(i+7) = aPulse->pdfval[i];
315 
316  for(int i=0; i<EcalPulseShape::TEMPLATESAMPLES;i++)
317  for(int j=0; j<EcalPulseShape::TEMPLATESAMPLES;j++)
318  fullpulsecov(i+7,j+7) = aPulseCov->covval[i][j];
319 
320 
321  // compute the right bin of the pulse shape using time calibration constants
323  EcalTimeCalibConstant itimeconst = 0;
324  if( it != itime->end() ) {
325  itimeconst = (*it);
326  } else {
327  edm::LogError("EcalRecHitError") << "No time intercalib const found for xtal "
328  << detid.rawId()
329  << "! something wrong with EcalTimeCalibConstants in your DB? ";
330  }
331 
332  // === amplitude computation ===
333  int leadingSample = ((EcalDataFrame)(*itdg)).lastUnsaturatedSample();
334 
335  if ( leadingSample == 4 ) { // saturation on the expected max sample
336  uncalibRecHit = EcalUncalibratedRecHit( (*itdg).id(), 4095*12, 0, 0, 0);
338  // do not propagate the default chi2 = -1 value to the calib rechit (mapped to 64), set it to 0 when saturation
339  uncalibRecHit.setChi2(0);
340  } else if ( leadingSample >= 0 ) { // saturation on other samples: cannot extrapolate from the fourth one
341  int gainId = ((EcalDataFrame)(*itdg)).sample(5).gainId();
342  if (gainId==0) gainId=3;
343  auto pedestal = pedVec[gainId-1];
344  auto gainratio = gainRatios[gainId-1];
345  double amplitude = ((double)(((EcalDataFrame)(*itdg)).sample(5).adc()) - pedestal) * gainratio;
346  uncalibRecHit = EcalUncalibratedRecHit( (*itdg).id(), amplitude, 0, 0, 0);
348  // do not propagate the default chi2 = -1 value to the calib rechit (mapped to 64), set it to 0 when saturation
349  uncalibRecHit.setChi2(0);
350  } else {
351  // multifit
352  bool barrel = detid.subdetId()==EcalBarrel;
353  const SampleMatrixGainArray &noisecors = noisecor(barrel);
354 
355  uncalibRecHit = multiFitMethod_.makeRecHit(*itdg, aped, aGain, noisecors,fullpulse,fullpulsecov,activeBX);
356 
357  // === time computation ===
358  if(timealgo_==ratioMethod) {
359  // ratio method
360  constexpr float clockToNsConstant = 25.;
361  constexpr float invClockToNs = 1/clockToNsConstant;
362  if (detid.subdetId()==EcalEndcap) {
363  ratioMethod_endcap_.init( *itdg, *sampleMask_, pedVec, pedRMSVec, gainRatios );
367  double theTimeCorrectionEE = timeCorrection(uncalibRecHit.amplitude(),
368  timeCorrBias_->EETimeCorrAmplitudeBins, timeCorrBias_->EETimeCorrShiftBins);
369 
370  uncalibRecHit.setJitter( crh.timeMax - 5 + theTimeCorrectionEE);
371  uncalibRecHit.setJitterError( std::sqrt(std::pow(crh.timeError,2) + std::pow(EEtimeConstantTerm_*invClockToNs,2)) );
372 
373  // consider flagging as kOutOfTime only if above noise
374  if (uncalibRecHit.amplitude() > pedRMSVec[0] * amplitudeThreshEE_){
375  float outOfTimeThreshP = outOfTimeThreshG12pEE_;
376  float outOfTimeThreshM = outOfTimeThreshG12mEE_;
377  // determine if gain has switched away from gainId==1 (x12 gain)
378  // and determine cuts (number of 'sigmas') to ose for kOutOfTime
379  // >3k ADC is necessasry condition for gain switch to occur
380  if (uncalibRecHit.amplitude() > 3000.){
381  for (int iSample = 0; iSample < EEDataFrame::MAXSAMPLES; iSample++) {
382  int GainId = ((EcalDataFrame)(*itdg)).sample(iSample).gainId();
383  if (GainId!=1) {
384  outOfTimeThreshP = outOfTimeThreshG61pEE_;
385  outOfTimeThreshM = outOfTimeThreshG61mEE_;
386  break;
387  }
388  }}
389  float correctedTime = (crh.timeMax-5) * clockToNsConstant + itimeconst + offsetTime;
390  float cterm = EEtimeConstantTerm_;
391  float sigmaped = pedRMSVec[0]; // approx for lower gains
392  float nterm = EEtimeNconst_*sigmaped/uncalibRecHit.amplitude();
393  float sigmat = std::sqrt( nterm*nterm + cterm*cterm );
394  if ( ( correctedTime > sigmat*outOfTimeThreshP ) ||
395  ( correctedTime < -sigmat*outOfTimeThreshM) )
396  { uncalibRecHit.setFlagBit( EcalUncalibratedRecHit::kOutOfTime ); }
397  }
398 
399  } else {
400  ratioMethod_barrel_.init( *itdg, *sampleMask_, pedVec, pedRMSVec, gainRatios );
405 
406  double theTimeCorrectionEB = timeCorrection(uncalibRecHit.amplitude(),
407  timeCorrBias_->EBTimeCorrAmplitudeBins, timeCorrBias_->EBTimeCorrShiftBins);
408 
409  uncalibRecHit.setJitter( crh.timeMax - 5 + theTimeCorrectionEB);
410  uncalibRecHit.setJitterError( std::sqrt(std::pow(crh.timeError,2) + std::pow(EBtimeConstantTerm_,2)/std::pow(clockToNsConstant,2)) );
411 
412  // consider flagging as kOutOfTime only if above noise
413  if (uncalibRecHit.amplitude() > pedRMSVec[0] * amplitudeThreshEB_){
414  float outOfTimeThreshP = outOfTimeThreshG12pEB_;
415  float outOfTimeThreshM = outOfTimeThreshG12mEB_;
416  // determine if gain has switched away from gainId==1 (x12 gain)
417  // and determine cuts (number of 'sigmas') to ose for kOutOfTime
418  // >3k ADC is necessasry condition for gain switch to occur
419  if (uncalibRecHit.amplitude() > 3000.){
420  for (int iSample = 0; iSample < EBDataFrame::MAXSAMPLES; iSample++) {
421  int GainId = ((EcalDataFrame)(*itdg)).sample(iSample).gainId();
422  if (GainId!=1) {
423  outOfTimeThreshP = outOfTimeThreshG61pEB_;
424  outOfTimeThreshM = outOfTimeThreshG61mEB_;
425  break;}
426  } }
427  float correctedTime = (crh.timeMax-5) * clockToNsConstant + itimeconst + offsetTime;
428  float cterm = EBtimeConstantTerm_;
429  float sigmaped = pedRMSVec[0]; // approx for lower gains
430  float nterm = EBtimeNconst_*sigmaped/uncalibRecHit.amplitude();
431  float sigmat = std::sqrt( nterm*nterm + cterm*cterm );
432  if ( ( correctedTime > sigmat*outOfTimeThreshP ) ||
433  ( correctedTime < -sigmat*outOfTimeThreshM ))
434  { uncalibRecHit.setFlagBit( EcalUncalibratedRecHit::kOutOfTime ); }
435  }
436 
437  }
438  } else if(timealgo_==weightsMethod) {
439  // weights method on the PU subtracted pulse shape
440  std::vector<double> amplitudes;
441  for(unsigned int ibx=0; ibx<activeBX.size(); ++ibx) amplitudes.push_back(uncalibRecHit.outOfTimeAmplitude(ibx));
442 
443  EcalTBWeights::EcalTDCId tdcid(1);
444  EcalTBWeights::EcalTBWeightMap const & wgtsMap = wgts->getMap();
445  EcalTBWeights::EcalTBWeightMap::const_iterator wit;
446  wit = wgtsMap.find( std::make_pair(*gid,tdcid) );
447  if( wit == wgtsMap.end() ) {
448  edm::LogError("EcalUncalibRecHitError") << "No weights found for EcalGroupId: "
449  << gid->id() << " and EcalTDCId: " << tdcid
450  << "\n skipping digi with id: " << detid.rawId();
451 
452  return false;
453  }
454  const EcalWeightSet& wset = wit->second; // this is the EcalWeightSet
455 
458 
459  weights[0] = &mat1;
460  weights[1] = &mat2;
461 
462  double timerh;
463  if (detid.subdetId()==EcalEndcap) {
464  timerh = weightsMethod_endcap_.time( *itdg, amplitudes, aped, aGain, fullpulse, weights);
465  } else {
466  timerh = weightsMethod_barrel_.time( *itdg, amplitudes, aped, aGain, fullpulse, weights);
467  }
468  uncalibRecHit.setJitter( timerh );
469  uncalibRecHit.setJitterError( 0. ); // not computed with weights
470  } else { // no time method;
471  uncalibRecHit.setJitter( 0. );
472  uncalibRecHit.setJitterError( 0. );
473  }
474  }
475 
476  // set flags if gain switch has occurred
477  if( ((EcalDataFrame)(*itdg)).hasSwitchToGain6() ) uncalibRecHit.setFlagBit( EcalUncalibratedRecHit::kHasSwitchToGain6 );
478  if( ((EcalDataFrame)(*itdg)).hasSwitchToGain1() ) uncalibRecHit.setFlagBit( EcalUncalibratedRecHit::kHasSwitchToGain1 );
479 
480  // set quality flags based on chi2 of the the fit
481  /*
482  if(detid.subdetId()==EcalEndcap) {
483  if(kPoorRecoFlagEE_ && uncalibRecHit.chi2()>chi2ThreshEE_) {
484  bool samplesok = true;
485  for (int sample =0; sample < EcalDataFrame::MAXSAMPLES; ++sample) {
486  if (!sampleMask_->useSampleEE(sample)) {
487  samplesok = false;
488  break;
489  }
490  }
491  if (samplesok) uncalibRecHit.setFlagBit(EcalUncalibratedRecHit::kPoorReco);
492  }
493  } else {
494  if(kPoorRecoFlagEB_ && uncalibRecHit.chi2()>chi2ThreshEB_) {
495  bool samplesok = true;
496  for (int sample =0; sample < EcalDataFrame::MAXSAMPLES; ++sample) {
497  if (!sampleMask_->useSampleEB(sample)) {
498  samplesok = false;
499  break;
500  }
501  }
502  if (samplesok) uncalibRecHit.setFlagBit(EcalUncalibratedRecHit::kPoorReco);
503  }
504  }
505  */
506 
507  // put the recHit in the collection
508  if (detid.subdetId()==EcalEndcap) {
509  result.push_back( uncalibRecHit );
510  } else {
511  result.push_back( uncalibRecHit );
512  }
513 
514  return true;
515 }
516 
517 /*
518 const SampleMatrix &EcalUncalibRecHitWorkerMultiFit::noisecor(bool barrel, int gain) const {
519  if (barrel) {
520  if (gain==6) {
521  return noisecorEBg6;
522  }
523  else if (gain==1) {
524  return noisecorEBg1;
525  }
526  else {
527  return noisecorEBg12;
528  }
529  }
530  else {
531  if (gain==6) {
532  return noisecorEEg6;
533  }
534  else if (gain==1) {
535  return noisecorEEg1;
536  }
537  else {
538  return noisecorEEg12;
539  }
540  }
541 
542  return noisecorEBg12;
543 
544 }
545 */
546 
549 
551  psd0.addNode((edm::ParameterDescription<std::vector<double>>("EBPulseShapeTemplate", {1.13979e-02, 7.58151e-01, 1.00000e+00, 8.87744e-01, 6.73548e-01, 4.74332e-01, 3.19561e-01, 2.15144e-01, 1.47464e-01, 1.01087e-01, 6.93181e-02, 4.75044e-02}, true) and
552  edm::ParameterDescription<std::vector<double>>("EEPulseShapeTemplate", {1.16442e-01, 7.56246e-01, 1.00000e+00, 8.97182e-01, 6.86831e-01, 4.91506e-01, 3.44111e-01, 2.45731e-01, 1.74115e-01, 1.23361e-01, 8.74288e-02, 6.19570e-02}, true)));
553 
554  psd0.addNode((edm::ParameterDescription<std::string>("EEdigiCollection", "", true) and
555  edm::ParameterDescription<std::string>("EBdigiCollection", "", true) and
556  edm::ParameterDescription<std::string>("ESdigiCollection", "", true) and
557  edm::ParameterDescription<bool>("UseLCcorrection", true, false) and
558  edm::ParameterDescription<std::vector<double>>("EBCorrNoiseMatrixG12", {1.00000, 0.71073, 0.55721, 0.46089, 0.40449, 0.35931, 0.33924, 0.32439, 0.31581, 0.30481 }, true) and
559  edm::ParameterDescription<std::vector<double>>("EECorrNoiseMatrixG12", {1.00000, 0.71373, 0.44825, 0.30152, 0.21609, 0.14786, 0.11772, 0.10165, 0.09465, 0.08098 }, true) and
560  edm::ParameterDescription<std::vector<double>>("EBCorrNoiseMatrixG06", {1.00000, 0.70946, 0.58021, 0.49846, 0.45006, 0.41366, 0.39699, 0.38478, 0.37847, 0.37055 }, true) and
561  edm::ParameterDescription<std::vector<double>>("EECorrNoiseMatrixG06", {1.00000, 0.71217, 0.47464, 0.34056, 0.26282, 0.20287, 0.17734, 0.16256, 0.15618, 0.14443 }, true) and
562  edm::ParameterDescription<std::vector<double>>("EBCorrNoiseMatrixG01", {1.00000, 0.73354, 0.64442, 0.58851, 0.55425, 0.53082, 0.51916, 0.51097, 0.50732, 0.50409 }, true) and
563  edm::ParameterDescription<std::vector<double>>("EECorrNoiseMatrixG01", {1.00000, 0.72698, 0.62048, 0.55691, 0.51848, 0.49147, 0.47813, 0.47007, 0.46621, 0.46265 }, true) and
564  edm::ParameterDescription<bool>("EcalPreMixStage1", false, true) and
565  edm::ParameterDescription<bool>("EcalPreMixStage2", false, true)));
566 
567  psd0.addOptionalNode((edm::ParameterDescription<std::vector<double>>("EBPulseShapeCovariance", {3.001e-06, 1.233e-05, 0.000e+00, -4.416e-06, -4.571e-06, -3.614e-06, -2.636e-06, -1.286e-06, -8.410e-07, -5.296e-07, 0.000e+00, 0.000e+00,
568  1.233e-05, 6.154e-05, 0.000e+00, -2.200e-05, -2.309e-05, -1.838e-05, -1.373e-05, -7.334e-06, -5.088e-06, -3.745e-06, -2.428e-06, 0.000e+00,
569  0.000e+00, 0.000e+00, 0.000e+00, 0.000e+00, 0.000e+00, 0.000e+00, 0.000e+00, 0.000e+00, 0.000e+00, 0.000e+00, 0.000e+00, 0.000e+00,
570  -4.416e-06, -2.200e-05, 0.000e+00, 8.319e-06, 8.545e-06, 6.792e-06, 5.059e-06, 2.678e-06, 1.816e-06, 1.223e-06, 8.245e-07, 5.589e-07,
571  -4.571e-06, -2.309e-05, 0.000e+00, 8.545e-06, 9.182e-06, 7.219e-06, 5.388e-06, 2.853e-06, 1.944e-06, 1.324e-06, 9.083e-07, 6.335e-07,
572  -3.614e-06, -1.838e-05, 0.000e+00, 6.792e-06, 7.219e-06, 6.016e-06, 4.437e-06, 2.385e-06, 1.636e-06, 1.118e-06, 7.754e-07, 5.556e-07,
573  -2.636e-06, -1.373e-05, 0.000e+00, 5.059e-06, 5.388e-06, 4.437e-06, 3.602e-06, 1.917e-06, 1.322e-06, 9.079e-07, 6.529e-07, 4.752e-07,
574  -1.286e-06, -7.334e-06, 0.000e+00, 2.678e-06, 2.853e-06, 2.385e-06, 1.917e-06, 1.375e-06, 9.100e-07, 6.455e-07, 4.693e-07, 3.657e-07,
575  -8.410e-07, -5.088e-06, 0.000e+00, 1.816e-06, 1.944e-06, 1.636e-06, 1.322e-06, 9.100e-07, 9.115e-07, 6.062e-07, 4.436e-07, 3.422e-07,
576  -5.296e-07, -3.745e-06, 0.000e+00, 1.223e-06, 1.324e-06, 1.118e-06, 9.079e-07, 6.455e-07, 6.062e-07, 7.217e-07, 4.862e-07, 3.768e-07,
577  0.000e+00, -2.428e-06, 0.000e+00, 8.245e-07, 9.083e-07, 7.754e-07, 6.529e-07, 4.693e-07, 4.436e-07, 4.862e-07, 6.509e-07, 4.418e-07,
578  0.000e+00, 0.000e+00, 0.000e+00, 5.589e-07, 6.335e-07, 5.556e-07, 4.752e-07, 3.657e-07, 3.422e-07, 3.768e-07, 4.418e-07, 6.142e-07}, true) and
579  edm::ParameterDescription<std::vector<double>>("EEPulseShapeCovariance", {3.941e-05, 3.333e-05, 0.000e+00, -1.449e-05, -1.661e-05, -1.424e-05, -1.183e-05, -6.842e-06, -4.915e-06, -3.411e-06, 0.000e+00, 0.000e+00,
580  3.333e-05, 2.862e-05, 0.000e+00, -1.244e-05, -1.431e-05, -1.233e-05, -1.032e-05, -5.883e-06, -4.154e-06, -2.902e-06, -2.128e-06, 0.000e+00,
581  0.000e+00, 0.000e+00, 0.000e+00, 0.000e+00, 0.000e+00, 0.000e+00, 0.000e+00, 0.000e+00, 0.000e+00, 0.000e+00, 0.000e+00, 0.000e+00,
582  -1.449e-05, -1.244e-05, 0.000e+00, 5.840e-06, 6.649e-06, 5.720e-06, 4.812e-06, 2.708e-06, 1.869e-06, 1.330e-06, 9.186e-07, 6.446e-07,
583  -1.661e-05, -1.431e-05, 0.000e+00, 6.649e-06, 7.966e-06, 6.898e-06, 5.794e-06, 3.157e-06, 2.184e-06, 1.567e-06, 1.084e-06, 7.575e-07,
584  -1.424e-05, -1.233e-05, 0.000e+00, 5.720e-06, 6.898e-06, 6.341e-06, 5.347e-06, 2.859e-06, 1.991e-06, 1.431e-06, 9.839e-07, 6.886e-07,
585  -1.183e-05, -1.032e-05, 0.000e+00, 4.812e-06, 5.794e-06, 5.347e-06, 4.854e-06, 2.628e-06, 1.809e-06, 1.289e-06, 9.020e-07, 6.146e-07,
586  -6.842e-06, -5.883e-06, 0.000e+00, 2.708e-06, 3.157e-06, 2.859e-06, 2.628e-06, 1.863e-06, 1.296e-06, 8.882e-07, 6.108e-07, 4.283e-07,
587  -4.915e-06, -4.154e-06, 0.000e+00, 1.869e-06, 2.184e-06, 1.991e-06, 1.809e-06, 1.296e-06, 1.217e-06, 8.669e-07, 5.751e-07, 3.882e-07,
588  -3.411e-06, -2.902e-06, 0.000e+00, 1.330e-06, 1.567e-06, 1.431e-06, 1.289e-06, 8.882e-07, 8.669e-07, 9.522e-07, 6.717e-07, 4.293e-07,
589  0.000e+00, -2.128e-06, 0.000e+00, 9.186e-07, 1.084e-06, 9.839e-07, 9.020e-07, 6.108e-07, 5.751e-07, 6.717e-07, 7.911e-07, 5.493e-07,
590  0.000e+00, 0.000e+00, 0.000e+00, 6.446e-07, 7.575e-07, 6.886e-07, 6.146e-07, 4.283e-07, 3.882e-07, 4.293e-07, 5.493e-07, 7.027e-07}, true)), true);
591 
593  psd.addNode(edm::ParameterDescription<std::vector<int>>("activeBXs", {-5,-4,-3,-2,-1,0,1,2,3,4}, true) and
594  edm::ParameterDescription<bool>("ampErrorCalculation", true, true) and
595  edm::ParameterDescription<bool>("useLumiInfoRunHeader", true, true) and
596  edm::ParameterDescription<int>("bunchSpacing", 0, true) and
597  edm::ParameterDescription<bool>("doPrefitEB", false, true) and
598  edm::ParameterDescription<bool>("doPrefitEE", false, true) and
599  edm::ParameterDescription<double>("prefitMaxChiSqEB", 25., true) and
600  edm::ParameterDescription<double>("prefitMaxChiSqEE", 10., true) and
601  edm::ParameterDescription<bool>("dynamicPedestalsEB", false, true) and
602  edm::ParameterDescription<bool>("dynamicPedestalsEE", false, true) and
603  edm::ParameterDescription<bool>("mitigateBadSamplesEB", false, true) and
604  edm::ParameterDescription<bool>("mitigateBadSamplesEE", false, true) and
605  edm::ParameterDescription<bool>("gainSwitchUseMaxSampleEB", false, true) and
606  edm::ParameterDescription<bool>("gainSwitchUseMaxSampleEE", false, true) and
607  edm::ParameterDescription<bool>("selectiveBadSampleCriteriaEB", false, true) and
608  edm::ParameterDescription<bool>("selectiveBadSampleCriteriaEE", false, true) and
609  edm::ParameterDescription<double>("addPedestalUncertaintyEB", 0., true) and
610  edm::ParameterDescription<double>("addPedestalUncertaintyEE", 0., true) and
611  edm::ParameterDescription<bool>("simplifiedNoiseModelForGainSwitch", true, true) and
612  edm::ParameterDescription<std::string>("timealgo", "RatioMethod", true) and
613  edm::ParameterDescription<std::vector<double>>("EBtimeFitParameters", {-2.015452e+00, 3.130702e+00, -1.234730e+01, 4.188921e+01, -8.283944e+01, 9.101147e+01, -5.035761e+01, 1.105621e+01}, true) and
614  edm::ParameterDescription<std::vector<double>>("EEtimeFitParameters", {-2.390548e+00, 3.553628e+00, -1.762341e+01, 6.767538e+01, -1.332130e+02, 1.407432e+02, -7.541106e+01, 1.620277e+01}, true) and
615  edm::ParameterDescription<std::vector<double>>("EBamplitudeFitParameters", {1.138,1.652}, true) and
616  edm::ParameterDescription<std::vector<double>>("EEamplitudeFitParameters", {1.890,1.400}, true) and
617  edm::ParameterDescription<double>("EBtimeFitLimits_Lower", 0.2, true) and
618  edm::ParameterDescription<double>("EBtimeFitLimits_Upper", 1.4, true) and
619  edm::ParameterDescription<double>("EEtimeFitLimits_Lower", 0.2, true) and
620  edm::ParameterDescription<double>("EEtimeFitLimits_Upper", 1.4, true) and
621  edm::ParameterDescription<double>("EBtimeConstantTerm", .6, true) and
622  edm::ParameterDescription<double>("EEtimeConstantTerm", 1.0, true) and
623  edm::ParameterDescription<double>("EBtimeNconst", 28.5, true) and
624  edm::ParameterDescription<double>("EEtimeNconst", 31.8, true) and
625  edm::ParameterDescription<double>("outOfTimeThresholdGain12pEB", 5, true) and
626  edm::ParameterDescription<double>("outOfTimeThresholdGain12mEB", 5, true) and
627  edm::ParameterDescription<double>("outOfTimeThresholdGain61pEB", 5, true) and
628  edm::ParameterDescription<double>("outOfTimeThresholdGain61mEB", 5, true) and
629  edm::ParameterDescription<double>("outOfTimeThresholdGain12pEE", 1000, true) and
630  edm::ParameterDescription<double>("outOfTimeThresholdGain12mEE", 1000, true) and
631  edm::ParameterDescription<double>("outOfTimeThresholdGain61pEE", 1000, true) and
632  edm::ParameterDescription<double>("outOfTimeThresholdGain61mEE", 1000, true) and
633  edm::ParameterDescription<double>("amplitudeThresholdEB", 10, true) and
634  edm::ParameterDescription<double>("amplitudeThresholdEE", 10, true) and
635  edm::ParameterDescription<double>("ebSpikeThreshold", 1.042, true) and
636  edm::ParameterDescription<std::vector<double>>("ebPulseShape", {5.2e-05,-5.26e-05 , 6.66e-05, 0.1168, 0.7575, 1., 0.8876, 0.6732, 0.4741, 0.3194}, true) and
637  edm::ParameterDescription<std::vector<double>>("eePulseShape", {5.2e-05,-5.26e-05 , 6.66e-05, 0.1168, 0.7575, 1., 0.8876, 0.6732, 0.4741, 0.3194}, true) and
638  edm::ParameterDescription<bool>("kPoorRecoFlagEB", true, true) and
639  edm::ParameterDescription<bool>("kPoorRecoFlagEE", false, true) and
640  edm::ParameterDescription<double>("chi2ThreshEB_", 65.0, true) and
641  edm::ParameterDescription<double>("chi2ThreshEE_", 50.0, true) and
642  edm::ParameterDescription<edm::ParameterSetDescription>("EcalPulseShapeParameters", psd0, true));
643 
644  return psd;
645 }
646 
647 
648 
edm::ParameterSetDescription getAlgoDescription() override
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
std::array< SampleMatrixGainArray, 2 > noisecors_
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)
unsigned size(int bx) const
EcalUncalibRecHitMultiFitAlgo multiFitMethod_
boost::transform_iterator< IterHelp, boost::counting_iterator< int > > const_iterator
void setJitter(float jitter)
int gainId(sample_type sample)
get the gainId (2 bits)
EcalUncalibRecHitTimeWeightsAlgo< EEDataFrame > weightsMethod_endcap_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
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.
void computeTime(std::vector< double > &timeFitParameters, std::pair< double, double > &timeFitLimits, std::vector< double > &amplitudeFitParameters)
static const int TEMPLATESAMPLES
edm::ESHandle< EcalWeightXtalGroups > grps
void push_back(T const &t)
ParameterDescriptionNode * addNode(ParameterDescriptionNode const &node)
void setJitterError(float jitterErr)
Eigen::Matrix< double, SampleVectorSize, SampleVectorSize > SampleMatrix
edm::EDGetTokenT< unsigned int > bunchSpacing_
EcalUncalibRecHitRatioMethodAlgo< EBDataFrame > ratioMethod_barrel_
#define constexpr
const unsigned int id() const
bool run(const edm::Event &evt, const EcalDigiCollection::const_iterator &digi, EcalUncalibratedRecHitCollection &result) override
tuple result
Definition: mps_fire.py:84
int hashedIndex(int ieta, int iphi)
Definition: EcalPyUtils.cc:42
std::array< SampleMatrix, NGains > SampleMatrixGainArray
float outOfTimeAmplitude(int bx) const
EcalUncalibRecHitRatioMethodAlgo< EEDataFrame > ratioMethod_endcap_
T sqrt(T t)
Definition: SSEVec.h:18
edm::ESHandle< EcalPulseShapes > pulseshapes
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
ParameterDescriptionNode * addOptionalNode(ParameterDescriptionNode const &node, bool writeToCfi)
edm::ESHandle< EcalTimeCalibConstants > itime
double timeCorrection(float ampli, const std::vector< float > &amplitudeBins, const std::vector< float > &shiftBins)
float gain6Over1() const
EcalWeightMatrix & getWeightsAfterGainSwitch()
Definition: EcalWeightSet.h:30
float covval[EcalPulseShape::TEMPLATESAMPLES][EcalPulseShape::TEMPLATESAMPLES]
EcalUncalibratedRecHit makeRecHit(const EcalDataFrame &dataFrame, const EcalPedestals::Item *aped, const EcalMGPAGainRatio *aGain, const SampleMatrixGainArray &noisecors, const FullSampleVector &fullpulse, const FullSampleMatrix &fullpulsecov, const BXVector &activeBX)
compute rechits
edm::ESHandle< EcalTimeBiasCorrections > timeCorrBias_
Definition: DetId.h:18
void set(const edm::EventSetup &es) override
edm::ESHandle< EcalPulseCovariances > pulsecovariances
int hashedIndex() const
Definition: EEDetId.h:182
EcalWeightMatrix & getWeightsBeforeGainSwitch()
Definition: EcalWeightSet.h:29
const T & get() const
Definition: EventSetup.h:56
std::vector< Item >::const_iterator const_iterator
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
float pdfval[TEMPLATESAMPLES]
float EcalTimeCalibConstant
math::Matrix< 3, 10 >::type EcalWeightMatrix
Definition: EcalWeightSet.h:22
Eigen::Matrix< double, FullSampleVectorSize, 1 > FullSampleVector
#define DEFINE_EDM_PLUGIN(factory, type, name)
const SampleMatrix & noisecor(bool barrel, int gain) const
edm::ESHandle< EcalSampleMask > sampleMaskHand_
const EcalWeightSet::EcalWeightMatrix * weights[2]
Eigen::Matrix< double, FullSampleVectorSize, FullSampleVectorSize > FullSampleMatrix
edm::ESHandle< EcalSamplesCorrelation > noisecovariances
void init(const C &dataFrame, const EcalSampleMask &sampleMask, const double *pedestals, const double *pedestalRMSes, const double *gainRatios)
static const int MAXSAMPLES
Definition: EcalDataFrame.h:48
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40