CMS 3D CMS Logo

List of all members | Public Member Functions | Public Attributes | Private Member Functions | Private Attributes
PulseShapeFitOOTPileupCorrection Class Reference

#include <PulseShapeFitOOTPileupCorrection.h>

Public Member Functions

void apply (const CaloSamples &cs, const std::vector< int > &capidvec, const HcalCalibrations &calibs, double &reconstructedEnergy, float &reconstructedTime, bool &useTriple, float &chi2) const
 
void phase1Apply (const HBHEChannelInfo &channelData, float &reconstructedEnergy, float &reconstructedTime, bool &useTriple, float &chi2) const
 
 PulseShapeFitOOTPileupCorrection ()
 
void resetPulseShapeTemplate (const HcalPulseShapes::Shape &ps)
 
void setChi2Term (bool isHPD)
 
void setPulseShapeTemplate (const HcalPulseShapes::Shape &ps, bool isHPD)
 
void setPUParams (bool iPedestalConstraint, bool iTimeConstraint, bool iAddPulseJitter, bool iApplyTimeSlew, double iTS4Min, const std::vector< double > &iTS4Max, double iPulseJitter, double iTimeMean, double iTimeSigHPD, double iTimeSigSiPM, double iPedMean, double iPedSigHPD, double iPedSigSiPM, double iNoiseHPD, double iNoiseSiPM, double iTMin, double iTMax, const std::vector< double > &its4Chi2, HcalTimeSlew::BiasSetting slewFlavor, int iFitTimes, bool iDCConstraint=false)
 
 ~PulseShapeFitOOTPileupCorrection ()
 

Public Attributes

const HcalPulseShapes::ShapecurrentPulseShape_ =nullptr
 

Private Member Functions

void fit (int iFit, float &timevalfit, float &chargevalfit, float &pedvalfit, float &chi2, bool &fitStatus, double &iTSMax, const double &iTSTOTen, double *iEnArr, int(&iBX)[3]) const
 
int pulseShapeFit (const double *energyArr, const double *pedenArr, const double *chargeArr, const double *pedArr, const double *gainArr, const double tsTOTen, std::vector< float > &fitParsVec, const double *ADCnoise) const
 

Private Attributes

bool addPulseJitter_
 
bool applyTimeSlew_
 
double chargeThreshold_
 
int cntsetPulseShape
 
bool dcConstraint_
 
std::unique_ptr< ROOT::Math::Functor > dpfunctor_
 
int fitTimes_
 
PSFitter::HybridMinimizerhybridfitter
 
std::array< double, HcalConst::maxSamplesiniTimesArr
 
bool isCurrentChannelHPD_
 
double noise_
 
double noiseHPD_
 
double noiseSiPM_
 
bool pedestalConstraint_
 
double pedMean_
 
double pedSig_
 
double pedSigHPD_
 
double pedSigSiPM_
 
std::unique_ptr< FitterFuncs::PulseShapeFunctorpsfPtr_
 
double pulseJitter_
 
HcalTimeSlew::BiasSetting slewFlavor_
 
std::unique_ptr< ROOT::Math::Functor > spfunctor_
 
bool timeConstraint_
 
double timeMean_
 
double timeSig_
 
double timeSigHPD_
 
double timeSigSiPM_
 
std::unique_ptr< ROOT::Math::Functor > tpfunctor_
 
double ts4Chi2_
 
double ts4Max_
 
double ts4Min_
 
int TSMax_
 
int TSMin_
 
bool unConstrainedFit_
 
std::vector< double > vts4Chi2_
 
std::vector< double > vts4Max_
 

Detailed Description

Definition at line 91 of file PulseShapeFitOOTPileupCorrection.h.

Constructor & Destructor Documentation

PulseShapeFitOOTPileupCorrection::PulseShapeFitOOTPileupCorrection ( )

Definition at line 210 of file PulseShapeFitOOTPileupCorrection.cc.

References hybridfitter, iniTimesArr, and PSFitter::HybridMinimizer::kMigrad.

210  : cntsetPulseShape(0),
211  psfPtr_(nullptr), spfunctor_(nullptr), dpfunctor_(nullptr), tpfunctor_(nullptr),
212  TSMin_(0), TSMax_(0), vts4Chi2_(0), pedestalConstraint_(false),
213  timeConstraint_(false), addPulseJitter_(false), applyTimeSlew_(false),
214  ts4Min_(0), vts4Max_(0), pulseJitter_(0), timeMean_(0), timeSig_(0), pedMean_(0), pedSig_(0),
215  noise_(0), dcConstraint_(false) {
217  iniTimesArr = { {-100,-75,-50,-25,0,25,50,75,100,125} };
218 }
std::array< double, HcalConst::maxSamples > iniTimesArr
std::unique_ptr< ROOT::Math::Functor > tpfunctor_
std::unique_ptr< ROOT::Math::Functor > dpfunctor_
std::unique_ptr< ROOT::Math::Functor > spfunctor_
std::unique_ptr< FitterFuncs::PulseShapeFunctor > psfPtr_
PulseShapeFitOOTPileupCorrection::~PulseShapeFitOOTPileupCorrection ( )

Definition at line 220 of file PulseShapeFitOOTPileupCorrection.cc.

References hybridfitter.

220  {
221  if(hybridfitter) delete hybridfitter;
222 }

Member Function Documentation

void PulseShapeFitOOTPileupCorrection::apply ( const CaloSamples cs,
const std::vector< int > &  capidvec,
const HcalCalibrations calibs,
double &  reconstructedEnergy,
float &  reconstructedTime,
bool &  useTriple,
float &  chi2 
) const

Definition at line 299 of file PulseShapeFitOOTPileupCorrection.cc.

References ALCARECOTkAlJpsiMuMu_cff::charge, muonCSCDigis_cfi::gain, HcalConst::maxSamples, noise_, HcalCalibrations::pedestal(), psfPtr_, pulseShapeFit(), HcalCalibrations::respcorrgain(), CaloSamples::size(), mathSSE::sqrt(), ts4Chi2_, ts4Max_, ts4Min_, vts4Chi2_, and vts4Max_.

Referenced by heavyIonTools.ConfigureHeavyIons::__call__(), editorTools.UserCodeTool::__call__(), HiCoreTools.RestrictInputToAOD::__call__(), coreTools.RunOnData::__call__(), trackTools.MakeAODTrackCandidates::__call__(), runJetUncertainties.RunJetUncertainties::__call__(), metTools.AddMETCollection::__call__(), heavyIonTools.ProductionDefaults::__call__(), editorTools.ChangeSource::__call__(), HiCoreTools.RemoveMCMatching::__call__(), coreTools.RemoveMCMatching::__call__(), trackTools.MakePATTrackCandidates::__call__(), trigTools.SwitchOnTrigger::__call__(), runMETCorrectionsAndUncertainties.RunMETCorrectionsAndUncertainties::__call__(), heavyIonTools.SelectionDefaults::__call__(), HiCoreTools.RemoveAllPATObjectsBut::__call__(), heavyIonTools.DisbaleMonteCarloDeps::__call__(), HiCoreTools.RemoveSpecificPATObjects::__call__(), trigTools.SwitchOnTriggerStandAlone::__call__(), trackTools.MakeTrackCandidates::__call__(), tauTools.AddTauCollection::__call__(), trigTools.SwitchOnTriggerMatching::__call__(), HiCoreTools.RemoveCleaning::__call__(), HiCoreTools.AddCleaning::__call__(), trigTools.SwitchOnTriggerMatchingStandAlone::__call__(), trigTools.SwitchOnTriggerMatchEmbedding::__call__(), jetTools.AddJetCollection::__call__(), jetTools.SwitchJetCollection::__call__(), jetTools.UpdateJetCollection::__call__(), jetTools.AddJetID::__call__(), jetTools.SetTagInfos::__call__(), HcalSimpleRecAlgoImpl::reco(), and HcalSimpleRecAlgoImpl::recoHBHE().

306 {
307  psfPtr_->setDefaultcntNANinfit();
308 
309  const unsigned int cssize = cs.size();
310 // initialize arrays to be zero
311  double chargeArr[HcalConst::maxSamples]={}, pedArr[HcalConst::maxSamples]={}, gainArr[HcalConst::maxSamples]={};
312  double energyArr[HcalConst::maxSamples]={}, pedenArr[HcalConst::maxSamples]={};
313  double noiseADCArr[HcalConst::maxSamples]={};
314  double noisePHArr[HcalConst::maxSamples]={};
315  double noiseArrSq[HcalConst::maxSamples]={};
316 
317  double tsTOT = 0, tstrig = 0; // in fC
318  double tsTOTen = 0; // in GeV
319  for(unsigned int ip=0; ip<cssize; ++ip){
320  if( ip >= (unsigned) HcalConst::maxSamples ) continue; // Too many samples than what we wanna fit (10 is enough...) -> skip them
321  const int capid = capidvec[ip];
322  double charge = cs[ip];
323  double ped = calibs.pedestal(capid);
324  double gain = calibs.respcorrgain(capid);
325 
326  double energy = charge*gain;
327  double peden = ped*gain;
328 
329  chargeArr[ip] = charge; pedArr[ip] = ped; gainArr[ip] = gain;
330  energyArr[ip] = energy; pedenArr[ip] = peden;
331 
332  // HPD fcByPE from
333  // https://github.com/cms-sw/cmssw/blob/CMSSW_8_1_0/SimCalorimetry/HcalSimProducers/python/hcalSimParameters_cfi.py#L62
334  if((charge-ped)>noise_) {
335  noisePHArr[ip] = sqrt((charge-ped)*0.3305); // Add photostatistics uncertainty
336  }
337 
338  noiseADCArr[ip] = psfPtr_->sigmaHPDQIE8(chargeArr[ip]); // Add Greg's channel discretization
339  noiseArrSq[ip] = noise_ * noise_ + noiseADCArr[ip] * noiseADCArr[ip] + noisePHArr[ip] * noisePHArr[ip];
340 
341  tsTOT += charge - ped;
342  tsTOTen += energy - peden;
343  if( ip ==4 || ip==5 ){
344  tstrig += charge - ped;
345  }
346  }
347 
348  ts4Max_=vts4Max_[0]; // HB and HE are identical for Run2
349  ts4Chi2_=vts4Chi2_[0]; // HB and HE are identical for Run2
350 
351  std::vector<float> fitParsVec;
352 
353  if(tstrig >= ts4Min_&& tsTOTen > 0.) { //Two sigma from 0
354  pulseShapeFit(energyArr, pedenArr, chargeArr, pedArr, gainArr, tsTOTen, fitParsVec, noiseArrSq);
355  } else {
356  fitParsVec.clear();
357  fitParsVec.push_back(0.); //charge
358  fitParsVec.push_back(-9999); // time
359  fitParsVec.push_back(0.); // ped
360  fitParsVec.push_back(-9999); // chi2
361  fitParsVec.push_back(false); // triple
362  }
363 
364  reconstructedEnergy=fitParsVec[0];
365  reconstructedTime=fitParsVec[1];
366  chi2 = fitParsVec[3];
367  useTriple=fitParsVec[4];
368 
369 }
double respcorrgain(int fCapId) const
get response corrected gain for capid=0..3
double pedestal(int fCapId) const
get pedestal for capid=0..3
T sqrt(T t)
Definition: SSEVec.h:18
int size() const
get the size
Definition: CaloSamples.h:24
std::unique_ptr< FitterFuncs::PulseShapeFunctor > psfPtr_
int pulseShapeFit(const double *energyArr, const double *pedenArr, const double *chargeArr, const double *pedArr, const double *gainArr, const double tsTOTen, std::vector< float > &fitParsVec, const double *ADCnoise) const
void PulseShapeFitOOTPileupCorrection::fit ( int  iFit,
float &  timevalfit,
float &  chargevalfit,
float &  pedvalfit,
float &  chi2,
bool &  fitStatus,
double &  iTSMax,
const double &  iTSTOTen,
double *  iEnArr,
int(&)  iBX[3] 
) const
private

Definition at line 437 of file PulseShapeFitOOTPileupCorrection.cc.

References funct::abs(), PSFitter::HybridMinimizer::Clear(), dpfunctor_, fitTimes_, hybridfitter, mps_fire::i, iniTimesArr, createfilelist::int, PSFitter::HybridMinimizer::kMigrad, PSFitter::HybridMinimizer::kScan, PSFitter::HybridMinimizer::Minimize(), PSFitter::HybridMinimizer::MinValue(), gen::n, pedMean_, pedSig_, mps_update::results, PSFitter::HybridMinimizer::SetFixedVariable(), PSFitter::HybridMinimizer::SetFunction(), PSFitter::HybridMinimizer::SetLimitedVariable(), PSFitter::HybridMinimizer::SetMinimizerType(), spfunctor_, timeMean_, timeSig_, tpfunctor_, TSMax_, TSMin_, varNames, and PSFitter::HybridMinimizer::X().

Referenced by trackingPlots.Iteration::modules(), and pulseShapeFit().

437  {
438  int n = 3;
439  if(iFit == 2) n = 5; //Two Pulse Fit
440  if(iFit == 3) n = 7; //Three Pulse Fit
441  //Step 1 Single Pulse fit
442  float pedMax = iTSMax; //=> max timeslice
443  float tMin = TSMin_; //Fitting Time Min
444  float tMax = TSMax_; //Fitting Time Max
445  //Checks to make sure fitting happens
446  if(pedMax < 1.) pedMax = 1.;
447  // Set starting values andf step sizes for parameters
448  double vstart[n];
449  for(int i = 0; i < int((n-1)/2); i++) {
450  vstart[2*i+0] = iniTimesArr[iBX[i]]+timeMean_;
451  vstart[2*i+1] = iEnArr[iBX[i]];
452  }
453  vstart[n-1] = pedMean_;
454 
455  double step[n];
456  for(int i = 0; i < n; i++) step[i] = 0.1;
457 
458  if(iFit == 1) hybridfitter->SetFunction(*spfunctor_);
459  if(iFit == 2) hybridfitter->SetFunction(*dpfunctor_);
460  if(iFit == 3) hybridfitter->SetFunction(*tpfunctor_);
461  hybridfitter->Clear();
462  //Times and amplitudes
463  for(int i = 0; i < int((n-1)/2); i++) {
464  hybridfitter->SetLimitedVariable(0+i*2, varNames[2*i+0] , vstart[0+i*2], step[0+i*2],iniTimesArr[iBX[i]]+tMin, iniTimesArr[ iBX[i] ]+tMax);
465  hybridfitter->SetLimitedVariable(1+i*2, varNames[2*i+1] , vstart[1+i*2], step[1+i*2], 0, 1.2*iTSTOTEn);
466  //Secret Option to fix the time
467  if(timeSig_ < 0) hybridfitter->SetFixedVariable(0+i*2, varNames[2*i+0],vstart[0+i*2]);
468  }
469  //Pedestal
470  if(vstart[n-1] > std::abs(pedMax)) vstart[n-1] = pedMax;
471  hybridfitter->SetLimitedVariable(n-1, varNames[n-1], vstart[n-1], step[n-1],-pedMax,pedMax);
472  //Secret Option to fix the pedestal
473  if(pedSig_ < 0) hybridfitter->SetFixedVariable(n-1,varNames[n-1],vstart[n-1]);
474  //a special number to label the initial condition
475  chi2=-1;
476  //3 fits why?!
477  const double *results = nullptr;
478  for(int tries=0; tries<=3;++tries){
479  if( fitTimes_ != 2 || tries !=1 ){
481  fitStatus = hybridfitter->Minimize();
482  }
483  double chi2valfit = hybridfitter->MinValue();
484  const double *newresults = hybridfitter->X();
485  if(chi2 == -1 || chi2>chi2valfit+0.01) {
486  results=newresults;
487  chi2=chi2valfit;
488  if( tries == 0 && fitTimes_ == 1 ) break;
489  if( tries == 1 && (fitTimes_ == 2 || fitTimes_ ==3 ) ) break;
490  if( tries == 2 && fitTimes_ == 4 ) break;
491  if( tries == 3 && fitTimes_ == 5 ) break;
492  //Secret option to speed up the fit => perhaps we should drop this
493  if(timeSig_ < 0 || pedSig_ < 0) break;
494  if(tries==0){
496  fitStatus = hybridfitter->Minimize();
497  } else if(tries==1){
498  hybridfitter->SetStrategy(1);
499  } else if(tries==2){
500  hybridfitter->SetStrategy(2);
501  }
502  } else {
503  break;
504  }
505  }
506  assert(results);
507 
508  timevalfit = results[0];
509  chargevalfit = results[1];
510  pedvalfit = results[n-1];
511 
512 }
double MinValue() const override
return minimum function value
constexpr char const * varNames[]
void SetFunction(const ROOT::Math::IMultiGenFunction &func) override
set the function to minimize
void SetMinimizerType(EMinimizerType type)
std::array< double, HcalConst::maxSamples > iniTimesArr
bool SetFixedVariable(unsigned int, const std::string &, double) override
set fixed variable (override if minimizer supports them )
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::unique_ptr< ROOT::Math::Functor > tpfunctor_
std::unique_ptr< ROOT::Math::Functor > dpfunctor_
std::unique_ptr< ROOT::Math::Functor > spfunctor_
const double * X() const override
return pointer to X values at the minimum
bool SetLimitedVariable(unsigned int ivar, const std::string &name, double val, double step, double, double) override
set upper/lower limited variable (override if minimizer supports them )
step
void PulseShapeFitOOTPileupCorrection::phase1Apply ( const HBHEChannelInfo channelData,
float &  reconstructedEnergy,
float &  reconstructedTime,
bool &  useTriple,
float &  chi2 
) const

Definition at line 514 of file PulseShapeFitOOTPileupCorrection.cc.

References ALCARECOTkAlJpsiMuMu_cff::charge, HBHEChannelInfo::darkCurrent(), dcConstraint_, HBHEChannelInfo::fcByPE(), muonCSCDigis_cfi::gain, HBHEChannelInfo::hasTimeInfo(), HBHEChannelInfo::lambda(), HcalConst::maxSamples, HBHEChannelInfo::nSamples(), psfPtr_, pulseShapeFit(), mathSSE::sqrt(), ts4Chi2_, ts4Max_, ts4Min_, HBHEChannelInfo::tsDFcPerADC(), HBHEChannelInfo::tsGain(), HBHEChannelInfo::tsPedestal(), HBHEChannelInfo::tsPedestalWidth(), HBHEChannelInfo::tsRawCharge(), vts4Chi2_, and vts4Max_.

Referenced by SimpleHBHEPhase1Algo::reconstruct().

519 {
520 
521  psfPtr_->setDefaultcntNANinfit();
522 
523  const unsigned cssize = channelData.nSamples();
524  // initialize arrays to be zero
525  double chargeArr[HcalConst::maxSamples]={}, pedArr[HcalConst::maxSamples]={}, gainArr[HcalConst::maxSamples]={};
526  double energyArr[HcalConst::maxSamples]={}, pedenArr[HcalConst::maxSamples]={};
527  double noiseADCArr[HcalConst::maxSamples]={};
528  double noiseDCArr[HcalConst::maxSamples]={};
529  double noiseArrSq[HcalConst::maxSamples]={};
530  double noisePHArr[HcalConst::maxSamples]={};
531  double tsTOT = 0, tstrig = 0; // in fC
532  double tsTOTen = 0; // in GeV
533  double sipmDarkCurrentWidth = psfPtr_->getSiPMDarkCurrent(channelData.darkCurrent(),channelData.fcByPE(),channelData.lambda());
534 
535  // go over the time slices
536  for(unsigned int ip=0; ip<cssize; ++ip){
537  if( ip >= (unsigned) HcalConst::maxSamples ) continue; // Too many samples than what we wanna fit (10 is enough...) -> skip them
538 
539  // const int capid = channelData.capid(); // not needed
540  double charge = channelData.tsRawCharge(ip);
541  double ped = channelData.tsPedestal(ip);
542  double gain = channelData.tsGain(ip);
543 
544  double energy = charge*gain;
545  double peden = ped*gain;
546 
547  chargeArr[ip] = charge; pedArr[ip] = ped; gainArr[ip] = gain;
548  energyArr[ip] = energy; pedenArr[ip] = peden;
549 
550  // quantization noise from the ADC (QIE8 or QIE10/11)
551  noiseADCArr[ip] = (1./sqrt(12))*channelData.tsDFcPerADC(ip);
552 
553  // dark current noise relevant for siPM
554  noiseDCArr[ip] = 0;
555  if(channelData.hasTimeInfo() && (charge-ped)>channelData.tsPedestalWidth(ip)) {
556  noiseDCArr[ip] = sipmDarkCurrentWidth;
557  }
558 
559  // Photo statistics uncertainties
560  // sigmaFC/FC = 1/sqrt(Ne);
561  // Note2. (from kPedro): the output number of photoelectrons after smearing is treated very differently for SiPMs: *each* pe is assigned a different time based on a random generation from the Y11 pulse plus the SimHit time. In HPDs, the overall pulse is shaped all at once using just the SimHit time.
562 
563  noisePHArr[ip] = 0;
564  if((charge-ped)>channelData.tsPedestalWidth(ip)) {
565  noisePHArr[ip] = sqrt((charge-ped)*channelData.fcByPE());
566  }
567 
568  // sum all in quadrature
569  noiseArrSq[ip]= noiseADCArr[ip]*noiseADCArr[ip] + noiseDCArr[ip]*noiseDCArr[ip] + channelData.tsPedestalWidth(ip)*channelData.tsPedestalWidth(ip) + noisePHArr[ip]*noisePHArr[ip];
570 
571  tsTOT += charge - ped;
572  tsTOTen += energy - peden;
573  if( ip ==4 || ip==5 ){
574  tstrig += charge - ped;
575  }
576  }
577 
578  double sipmDarkCurrentWidth2 = 0.;
579  if(dcConstraint_) sipmDarkCurrentWidth2 = sipmDarkCurrentWidth*sipmDarkCurrentWidth;
580  double averagePedSig2GeV=0.25*((channelData.tsPedestalWidth(0)*channelData.tsPedestalWidth(0)+sipmDarkCurrentWidth2)*channelData.tsGain(0)*channelData.tsGain(0) +
581  (channelData.tsPedestalWidth(1)*channelData.tsPedestalWidth(1)+sipmDarkCurrentWidth2)*channelData.tsGain(1)*channelData.tsGain(1) +
582  (channelData.tsPedestalWidth(2)*channelData.tsPedestalWidth(2)+sipmDarkCurrentWidth2)*channelData.tsGain(2)*channelData.tsGain(2) +
583  (channelData.tsPedestalWidth(3)*channelData.tsPedestalWidth(3)+sipmDarkCurrentWidth2)*channelData.tsGain(3)*channelData.tsGain(3));
584 
585  // redefine the invertpedSig2
586  psfPtr_->setinvertpedSig2(1./(averagePedSig2GeV));
587 
588  if(channelData.hasTimeInfo()) {
590  } else {
592  }
593 
594  std::vector<float> fitParsVec;
595  if(tstrig >= ts4Min_ && tsTOTen > 0.) { //Two sigma from 0
596  pulseShapeFit(energyArr, pedenArr, chargeArr, pedArr, gainArr, tsTOTen, fitParsVec,noiseArrSq);
597  } else {
598  fitParsVec.clear();
599  fitParsVec.push_back(0.); //charge
600  fitParsVec.push_back(-9999); // time
601  fitParsVec.push_back(0.); // ped
602  fitParsVec.push_back(-9999); // chi2
603  fitParsVec.push_back(false); // triple
604  }
605 
606 
607  reconstructedEnergy = fitParsVec[0];
608  reconstructedTime = fitParsVec[1];
609  chi2 = fitParsVec[3];
610  useTriple = fitParsVec[4];
611 
612 }
double lambda() const
bool hasTimeInfo() const
double tsGain(const unsigned ts) const
double tsPedestal(const unsigned ts) const
double tsRawCharge(const unsigned ts) const
T sqrt(T t)
Definition: SSEVec.h:18
double fcByPE() const
double darkCurrent() const
double tsPedestalWidth(const unsigned ts) const
std::unique_ptr< FitterFuncs::PulseShapeFunctor > psfPtr_
unsigned nSamples() const
int pulseShapeFit(const double *energyArr, const double *pedenArr, const double *chargeArr, const double *pedArr, const double *gainArr, const double tsTOTen, std::vector< float > &fitParsVec, const double *ADCnoise) const
float tsDFcPerADC(const unsigned ts) const
int PulseShapeFitOOTPileupCorrection::pulseShapeFit ( const double *  energyArr,
const double *  pedenArr,
const double *  chargeArr,
const double *  pedArr,
const double *  gainArr,
const double  tsTOTen,
std::vector< float > &  fitParsVec,
const double *  ADCnoise 
) const
private

Definition at line 373 of file PulseShapeFitOOTPileupCorrection.cc.

References funct::abs(), applyTimeSlew_, rpcdqm::BX, vertices_cff::chi2, HcalTimeSlew::delay(), fit(), mps_fire::i, SiStripPI::max, HcalConst::maxSamples, psfPtr_, slewFlavor_, mathSSE::sqrt(), ts4Chi2_, and ts4Max_.

Referenced by apply(), and phase1Apply().

373  {
374  double tsMAX=0;
376  double tstrig = 0; // in fC
377  for(int i=0;i<HcalConst::maxSamples;++i){
378  tmpx[i]=i;
379  tmpy[i]=energyArr[i]-pedenArr[i];
380  //Add Time Slew !!! does this need to be pedestal subtracted
381  tmpslew[i] = 0;
382  if(applyTimeSlew_) tmpslew[i] = HcalTimeSlew::delay(std::max(1.0,chargeArr[i]),slewFlavor_);
383  // add the noise components
384  tmperry2[i]=noiseArrSq[i];
385 
386  //Propagate it through
387  tmperry2[i]*=(gainArr[i]*gainArr[i]); //Convert from fC to GeV
388  tmperry [i]=sqrt(tmperry2[i]);
389 
390  if(std::abs(energyArr[i])>tsMAX) tsMAX=std::abs(tmpy[i]);
391  if( i ==4 || i ==5 ){
392  tstrig += chargeArr[i] - pedArr[i];
393  }
394  }
395  psfPtr_->setpsFitx (tmpx);
396  psfPtr_->setpsFity (tmpy);
397  psfPtr_->setpsFiterry (tmperry);
398  psfPtr_->setpsFiterry2(tmperry2);
399  psfPtr_->setpsFitslew (tmpslew);
400 
401  //Fit 1 single pulse
402  float timevalfit = 0;
403  float chargevalfit= 0;
404  float pedvalfit = 0;
405  float chi2 = 999; //cannot be zero
406  bool fitStatus = false;
407  bool useTriple = false;
408 
409  int BX[3] = {4,5,3};
410  if(ts4Chi2_ != 0) fit(1,timevalfit,chargevalfit,pedvalfit,chi2,fitStatus,tsMAX,tsTOTen,tmpy,BX);
411 // Based on the pulse shape ( 2. likely gives the same performance )
412  if(tmpy[2] > 3.*tmpy[3]) BX[2] = 2;
413 // Only do three-pulse fit when tstrig < ts4Max_, otherwise one-pulse fit is used (above)
414  if(chi2 > ts4Chi2_ && tstrig < ts4Max_) { //fails chi2 cut goes straight to 3 Pulse fit
415  fit(3,timevalfit,chargevalfit,pedvalfit,chi2,fitStatus,tsMAX,tsTOTen,tmpy,BX);
416  useTriple=true;
417  }
418 
419  /*
420  if(chi2 > ts345Chi2_) { //fails do two pulse chi2 for TS5
421  BX[1] = 5;
422  fit(3,timevalfit,chargevalfit,pedvalfit,chi2,fitStatus,tsMAX,tsTOTen,BX);
423  }
424  */
425  //Fix back the timeslew
426  //if(applyTimeSlew_) timevalfit+=HcalTimeSlew::delay(std::max(1.0,chargeArr[4]),slewFlavor_);
427  int outfitStatus = (fitStatus ? 1: 0 );
428  fitParsVec.clear();
429  fitParsVec.push_back(chargevalfit);
430  fitParsVec.push_back(timevalfit);
431  fitParsVec.push_back(pedvalfit);
432  fitParsVec.push_back(chi2);
433  fitParsVec.push_back(useTriple);
434  return outfitStatus;
435 }
void fit(int iFit, float &timevalfit, float &chargevalfit, float &pedvalfit, float &chi2, bool &fitStatus, double &iTSMax, const double &iTSTOTen, double *iEnArr, int(&iBX)[3]) const
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::unique_ptr< FitterFuncs::PulseShapeFunctor > psfPtr_
static double delay(double fC, BiasSetting bias=Medium)
Returns the amount (ns) by which a pulse of the given number of fC will be delayed by the timeslew ef...
Definition: HcalTimeSlew.cc:16
void PulseShapeFitOOTPileupCorrection::resetPulseShapeTemplate ( const HcalPulseShapes::Shape ps)

Definition at line 289 of file PulseShapeFitOOTPileupCorrection.cc.

References addPulseJitter_, applyTimeSlew_, cntsetPulseShape, FitterFuncs::PulseShapeFunctor::doublePulseShapeFunc(), dpfunctor_, noise_, pedestalConstraint_, pedMean_, pedSig_, psfPtr_, pulseJitter_, FitterFuncs::PulseShapeFunctor::singlePulseShapeFunc(), spfunctor_, timeConstraint_, timeMean_, timeSig_, tpfunctor_, and FitterFuncs::PulseShapeFunctor::triplePulseShapeFunc().

Referenced by setPulseShapeTemplate().

289  {
290  ++ cntsetPulseShape;
293  spfunctor_ = std::unique_ptr<ROOT::Math::Functor>( new ROOT::Math::Functor(psfPtr_.get(),&FitterFuncs::PulseShapeFunctor::singlePulseShapeFunc, 3) );
294  dpfunctor_ = std::unique_ptr<ROOT::Math::Functor>( new ROOT::Math::Functor(psfPtr_.get(),&FitterFuncs::PulseShapeFunctor::doublePulseShapeFunc, 5) );
295  tpfunctor_ = std::unique_ptr<ROOT::Math::Functor>( new ROOT::Math::Functor(psfPtr_.get(),&FitterFuncs::PulseShapeFunctor::triplePulseShapeFunc, 7) );
296 
297 }
std::unique_ptr< ROOT::Math::Functor > tpfunctor_
std::unique_ptr< ROOT::Math::Functor > dpfunctor_
std::unique_ptr< ROOT::Math::Functor > spfunctor_
std::unique_ptr< FitterFuncs::PulseShapeFunctor > psfPtr_
void PulseShapeFitOOTPileupCorrection::setChi2Term ( bool  isHPD)

Definition at line 224 of file PulseShapeFitOOTPileupCorrection.cc.

References noise_, noiseHPD_, noiseSiPM_, pedSig_, pedSigHPD_, pedSigSiPM_, timeSig_, timeSigHPD_, and timeSigSiPM_.

Referenced by setPulseShapeTemplate().

224  {
225 
226  if(isHPD) {
228  pedSig_ = pedSigHPD_; // unused in Phase1, take the average pedestal width from the DB
229  noise_ = noiseHPD_; // unused in Phase1, take the pedestal width from the DB
230  } else {
232  pedSig_ = pedSigSiPM_; // unused in Phase1, take the average pedestal width from the DB
233  noise_ = noiseSiPM_; // unused in Phase1, take the pedestal width from the DB
234  }
235 
236 }
void PulseShapeFitOOTPileupCorrection::setPulseShapeTemplate ( const HcalPulseShapes::Shape ps,
bool  isHPD 
)

Definition at line 277 of file PulseShapeFitOOTPileupCorrection.cc.

References currentPulseShape_, isCurrentChannelHPD_, resetPulseShapeTemplate(), and setChi2Term().

277  {
278  // initialize for every different channel types (HPD vs SiPM)
279 
280  if (!(&ps == currentPulseShape_ && isHPD == isCurrentChannelHPD_))
281  {
282  setChi2Term(isHPD);
284  currentPulseShape_ = &ps;
285  isCurrentChannelHPD_ = isHPD;
286  }
287 }
void resetPulseShapeTemplate(const HcalPulseShapes::Shape &ps)
const HcalPulseShapes::Shape * currentPulseShape_
void PulseShapeFitOOTPileupCorrection::setPUParams ( bool  iPedestalConstraint,
bool  iTimeConstraint,
bool  iAddPulseJitter,
bool  iApplyTimeSlew,
double  iTS4Min,
const std::vector< double > &  iTS4Max,
double  iPulseJitter,
double  iTimeMean,
double  iTimeSigHPD,
double  iTimeSigSiPM,
double  iPedMean,
double  iPedSigHPD,
double  iPedSigSiPM,
double  iNoiseHPD,
double  iNoiseSiPM,
double  iTMin,
double  iTMax,
const std::vector< double > &  its4Chi2,
HcalTimeSlew::BiasSetting  slewFlavor,
int  iFitTimes,
bool  iDCConstraint = false 
)

Definition at line 239 of file PulseShapeFitOOTPileupCorrection.cc.

References addPulseJitter_, applyTimeSlew_, dcConstraint_, fitTimes_, noiseHPD_, noiseSiPM_, pedestalConstraint_, pedMean_, pedSigHPD_, pedSigSiPM_, pulseJitter_, slewFlavor_, timeConstraint_, timeMean_, timeSigHPD_, timeSigSiPM_, ts4Min_, TSMax_, TSMin_, vts4Chi2_, and vts4Max_.

246  {
247 
248  TSMin_ = iTMin;
249  TSMax_ = iTMax;
250  // ts4Chi2_ = its4Chi2;
251  vts4Chi2_ = its4Chi2;
252  pedestalConstraint_ = iPedestalConstraint;
253  timeConstraint_ = iTimeConstraint;
254  addPulseJitter_ = iAddPulseJitter;
255  applyTimeSlew_ = iApplyTimeSlew;
256  ts4Min_ = iTS4Min;
257  // ts4Max_ = iTS4Max;
258  vts4Max_ = iTS4Max;
259  pulseJitter_ = iPulseJitter*iPulseJitter;
260  timeMean_ = iTimeMean;
261  // timeSig_ = iTimeSig;
262  timeSigHPD_ = iTimeSigHPD;
263  timeSigSiPM_ = iTimeSigSiPM;
264  pedMean_ = iPedMean;
265  // pedSig_ = iPedSig;
266  pedSigHPD_ = iPedSigHPD;
267  pedSigSiPM_ = iPedSigSiPM;
268  // noise_ = iNoise;
269  noiseHPD_ = iNoiseHPD;
270  noiseSiPM_ = iNoiseSiPM;
271  slewFlavor_ = slewFlavor;
272  fitTimes_ = iFitTimes;
273  dcConstraint_ = iDCConstraint;
274 
275 }

Member Data Documentation

bool PulseShapeFitOOTPileupCorrection::addPulseJitter_
private

Definition at line 148 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by resetPulseShapeTemplate(), and setPUParams().

bool PulseShapeFitOOTPileupCorrection::applyTimeSlew_
private
double PulseShapeFitOOTPileupCorrection::chargeThreshold_
private

Definition at line 135 of file PulseShapeFitOOTPileupCorrection.h.

int PulseShapeFitOOTPileupCorrection::cntsetPulseShape
private

Definition at line 133 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by resetPulseShapeTemplate().

const HcalPulseShapes::Shape* PulseShapeFitOOTPileupCorrection::currentPulseShape_ =nullptr

Definition at line 120 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by setPulseShapeTemplate().

bool PulseShapeFitOOTPileupCorrection::dcConstraint_
private

Definition at line 167 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by phase1Apply(), and setPUParams().

std::unique_ptr<ROOT::Math::Functor> PulseShapeFitOOTPileupCorrection::dpfunctor_
private

Definition at line 140 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by fit(), and resetPulseShapeTemplate().

int PulseShapeFitOOTPileupCorrection::fitTimes_
private

Definition at line 136 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by fit(), and setPUParams().

PSFitter::HybridMinimizer* PulseShapeFitOOTPileupCorrection::hybridfitter
private
std::array<double,HcalConst::maxSamples> PulseShapeFitOOTPileupCorrection::iniTimesArr
private

Definition at line 134 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by fit(), and PulseShapeFitOOTPileupCorrection().

bool PulseShapeFitOOTPileupCorrection::isCurrentChannelHPD_
private

Definition at line 169 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by setPulseShapeTemplate().

double PulseShapeFitOOTPileupCorrection::noise_
private

Definition at line 163 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by apply(), resetPulseShapeTemplate(), and setChi2Term().

double PulseShapeFitOOTPileupCorrection::noiseHPD_
private

Definition at line 164 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by setChi2Term(), and setPUParams().

double PulseShapeFitOOTPileupCorrection::noiseSiPM_
private

Definition at line 165 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by setChi2Term(), and setPUParams().

bool PulseShapeFitOOTPileupCorrection::pedestalConstraint_
private

Definition at line 146 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by resetPulseShapeTemplate(), and setPUParams().

double PulseShapeFitOOTPileupCorrection::pedMean_
private

Definition at line 159 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by fit(), resetPulseShapeTemplate(), and setPUParams().

double PulseShapeFitOOTPileupCorrection::pedSig_
private

Definition at line 160 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by fit(), resetPulseShapeTemplate(), and setChi2Term().

double PulseShapeFitOOTPileupCorrection::pedSigHPD_
private

Definition at line 161 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by setChi2Term(), and setPUParams().

double PulseShapeFitOOTPileupCorrection::pedSigSiPM_
private

Definition at line 162 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by setChi2Term(), and setPUParams().

std::unique_ptr<FitterFuncs::PulseShapeFunctor> PulseShapeFitOOTPileupCorrection::psfPtr_
private
double PulseShapeFitOOTPileupCorrection::pulseJitter_
private

Definition at line 154 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by resetPulseShapeTemplate(), and setPUParams().

HcalTimeSlew::BiasSetting PulseShapeFitOOTPileupCorrection::slewFlavor_
private

Definition at line 166 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by pulseShapeFit(), and setPUParams().

std::unique_ptr<ROOT::Math::Functor> PulseShapeFitOOTPileupCorrection::spfunctor_
private

Definition at line 139 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by fit(), and resetPulseShapeTemplate().

bool PulseShapeFitOOTPileupCorrection::timeConstraint_
private

Definition at line 147 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by resetPulseShapeTemplate(), and setPUParams().

double PulseShapeFitOOTPileupCorrection::timeMean_
private

Definition at line 155 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by fit(), resetPulseShapeTemplate(), and setPUParams().

double PulseShapeFitOOTPileupCorrection::timeSig_
private

Definition at line 156 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by fit(), resetPulseShapeTemplate(), and setChi2Term().

double PulseShapeFitOOTPileupCorrection::timeSigHPD_
private

Definition at line 157 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by setChi2Term(), and setPUParams().

double PulseShapeFitOOTPileupCorrection::timeSigSiPM_
private

Definition at line 158 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by setChi2Term(), and setPUParams().

std::unique_ptr<ROOT::Math::Functor> PulseShapeFitOOTPileupCorrection::tpfunctor_
private

Definition at line 141 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by fit(), and resetPulseShapeTemplate().

double PulseShapeFitOOTPileupCorrection::ts4Chi2_
mutableprivate

Definition at line 144 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by apply(), phase1Apply(), and pulseShapeFit().

double PulseShapeFitOOTPileupCorrection::ts4Max_
mutableprivate

Definition at line 152 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by apply(), phase1Apply(), and pulseShapeFit().

double PulseShapeFitOOTPileupCorrection::ts4Min_
private

Definition at line 151 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by apply(), phase1Apply(), and setPUParams().

int PulseShapeFitOOTPileupCorrection::TSMax_
private

Definition at line 143 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by fit(), and setPUParams().

int PulseShapeFitOOTPileupCorrection::TSMin_
private

Definition at line 142 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by fit(), and setPUParams().

bool PulseShapeFitOOTPileupCorrection::unConstrainedFit_
private

Definition at line 149 of file PulseShapeFitOOTPileupCorrection.h.

std::vector<double> PulseShapeFitOOTPileupCorrection::vts4Chi2_
private

Definition at line 145 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by apply(), phase1Apply(), and setPUParams().

std::vector<double> PulseShapeFitOOTPileupCorrection::vts4Max_
private

Definition at line 153 of file PulseShapeFitOOTPileupCorrection.h.

Referenced by apply(), phase1Apply(), and setPUParams().