CMS 3D CMS Logo

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

#include <SimpleHBHEPhase1Algo.h>

Inheritance diagram for SimpleHBHEPhase1Algo:
AbsHBHEPhase1Algo

Public Member Functions

void beginRun (const edm::Run &, const edm::EventSetup &) override
 
void endRun () override
 
int getFirstSampleShift () const
 
float getPhaseNS () const
 
int getRunNumber () const
 
int getSamplesToAdd () const
 
float getTimeShift () const
 
bool isConfigurable () const override
 
bool isCorrectingForPhaseContainment () const
 
HBHERecHit reconstruct (const HBHEChannelInfo &info, const HcalRecoParam *params, const HcalCalibrations &calibs, bool isRealData) override
 
 SimpleHBHEPhase1Algo (int firstSampleShift, int samplesToAdd, float phaseNS, float timeShift, bool correctForPhaseContainment, std::unique_ptr< PulseShapeFitOOTPileupCorrection > m2, std::unique_ptr< HcalDeterministicFit > detFit, std::unique_ptr< MahiFit > mahi)
 
 ~SimpleHBHEPhase1Algo () override
 
- Public Member Functions inherited from AbsHBHEPhase1Algo
virtual bool configure (const AbsHcalAlgoData *)
 
virtual ~AbsHBHEPhase1Algo ()
 

Public Attributes

const HcalTimeSlewhcalTimeSlew_delay_
 

Protected Member Functions

float hbminusCorrectionFactor (const HcalDetId &cell, float energy, bool isRealData) const
 
float m0Energy (const HBHEChannelInfo &info, double reconstructedCharge, bool applyContainmentCorrection, double phaseNS, int nSamplesToAdd)
 
float m0Time (const HBHEChannelInfo &info, double reconstructedCharge, const HcalCalibrations &calibs, int nSamplesToExamine) const
 

Private Attributes

bool corrFPC_
 
int firstSampleShift_
 
std::unique_ptr< HcalDeterministicFithltOOTpuCorr_
 
std::unique_ptr< MahiFitmahiOOTpuCorr_
 
float phaseNS_
 
std::unique_ptr< PulseShapeFitOOTPileupCorrectionpsFitOOTpuCorr_
 
HcalPulseContainmentManager pulseCorr_
 
int runnum_
 
int samplesToAdd_
 
HcalPulseShapes theHcalPulseShapes_
 
float timeShift_
 

Detailed Description

Definition at line 18 of file SimpleHBHEPhase1Algo.h.

Constructor & Destructor Documentation

SimpleHBHEPhase1Algo::SimpleHBHEPhase1Algo ( int  firstSampleShift,
int  samplesToAdd,
float  phaseNS,
float  timeShift,
bool  correctForPhaseContainment,
std::unique_ptr< PulseShapeFitOOTPileupCorrection m2,
std::unique_ptr< HcalDeterministicFit detFit,
std::unique_ptr< MahiFit mahi 
)

Definition at line 20 of file SimpleHBHEPhase1Algo.cc.

References hcalTimeSlew_delay_.

32  phaseNS_(phaseNS),
34  runnum_(0),
37  hltOOTpuCorr_(std::move(detFit)),
39 {
40  hcalTimeSlew_delay_ = nullptr;
41 }
const HcalTimeSlew * hcalTimeSlew_delay_
std::unique_ptr< HcalDeterministicFit > hltOOTpuCorr_
float PulseContainmentFractionalError
HcalPulseContainmentManager pulseCorr_
std::unique_ptr< MahiFit > mahiOOTpuCorr_
std::unique_ptr< PulseShapeFitOOTPileupCorrection > psFitOOTpuCorr_
def move(src, dest)
Definition: eostools.py:510
SimpleHBHEPhase1Algo::~SimpleHBHEPhase1Algo ( )
inlineoverride

Definition at line 52 of file SimpleHBHEPhase1Algo.h.

References beginRun(), and endRun().

52 {}

Member Function Documentation

void SimpleHBHEPhase1Algo::beginRun ( const edm::Run r,
const edm::EventSetup es 
)
overridevirtual

Reimplemented from AbsHBHEPhase1Algo.

Definition at line 43 of file SimpleHBHEPhase1Algo.cc.

References HcalPulseContainmentManager::beginRun(), edm::EventSetup::get(), hcalTimeSlew_delay_, pulseCorr_, edm::RunBase::run(), and runnum_.

Referenced by ~SimpleHBHEPhase1Algo().

45 {
47  es.get<HcalTimeSlewRecord>().get("HBHE", delay);
48  hcalTimeSlew_delay_ = &*delay;
49 
50  runnum_ = r.run();
51  pulseCorr_.beginRun(es);
52 }
const HcalTimeSlew * hcalTimeSlew_delay_
RunNumber_t run() const
Definition: RunBase.h:40
HcalPulseContainmentManager pulseCorr_
const T & get() const
Definition: EventSetup.h:59
void beginRun(edm::EventSetup const &es)
void SimpleHBHEPhase1Algo::endRun ( )
overridevirtual

Reimplemented from AbsHBHEPhase1Algo.

Definition at line 54 of file SimpleHBHEPhase1Algo.cc.

References runnum_.

Referenced by ~SimpleHBHEPhase1Algo().

55 {
56  runnum_ = 0;
57 }
int SimpleHBHEPhase1Algo::getFirstSampleShift ( ) const
inline

Definition at line 65 of file SimpleHBHEPhase1Algo.h.

References firstSampleShift_.

float SimpleHBHEPhase1Algo::getPhaseNS ( ) const
inline

Definition at line 67 of file SimpleHBHEPhase1Algo.h.

References phaseNS_.

67 {return phaseNS_;}
int SimpleHBHEPhase1Algo::getRunNumber ( ) const
inline

Definition at line 70 of file SimpleHBHEPhase1Algo.h.

References runnum_.

int SimpleHBHEPhase1Algo::getSamplesToAdd ( ) const
inline

Definition at line 66 of file SimpleHBHEPhase1Algo.h.

References samplesToAdd_.

float SimpleHBHEPhase1Algo::getTimeShift ( ) const
inline

Definition at line 68 of file SimpleHBHEPhase1Algo.h.

References timeShift_.

float SimpleHBHEPhase1Algo::hbminusCorrectionFactor ( const HcalDetId cell,
float  energy,
bool  isRealData 
) const
protected

Definition at line 159 of file SimpleHBHEPhase1Algo.cc.

References corr, hbminus_special_ecorr(), HcalBarrel, HcalDetId::ieta(), HcalDetId::iphi(), runnum_, and HcalDetId::subdet().

Referenced by reconstruct().

162 {
163  float corr = 1.f;
164  if (isRealData && runnum_ > 0)
165  if (cell.subdet() == HcalBarrel)
166  {
167  const int ieta = cell.ieta();
168  const int iphi = cell.iphi();
169  corr = hbminus_special_ecorr(ieta, iphi, energy, runnum_);
170  }
171  return corr;
172 }
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:49
float hbminus_special_ecorr(int ieta, int iphi, double energy, int runnum)
Special energy correction for some HB- cells.
int ieta() const
get the cell ieta
Definition: HcalDetId.h:56
JetCorrectorParameters corr
Definition: classes.h:5
int iphi() const
get the cell iphi
Definition: HcalDetId.cc:103
bool SimpleHBHEPhase1Algo::isConfigurable ( ) const
inlineoverridevirtual

Implements AbsHBHEPhase1Algo.

Definition at line 58 of file SimpleHBHEPhase1Algo.h.

References info(), and reconstruct().

58 {return false;}
bool SimpleHBHEPhase1Algo::isCorrectingForPhaseContainment ( ) const
inline

Definition at line 69 of file SimpleHBHEPhase1Algo.h.

References corrFPC_.

float SimpleHBHEPhase1Algo::m0Energy ( const HBHEChannelInfo info,
double  reconstructedCharge,
bool  applyContainmentCorrection,
double  phaseNS,
int  nSamplesToAdd 
)
protected

Definition at line 174 of file SimpleHBHEPhase1Algo.cc.

References MillePedeFileConverter_cfg::e, HBHEChannelInfo::energyInWindow(), firstSampleShift_, HcalPulseContainmentManager::get(), HBHEChannelInfo::id(), pulseCorr_, and HBHEChannelInfo::soi().

Referenced by reconstruct().

179 {
180  int ibeg = static_cast<int>(info.soi()) + firstSampleShift_;
181  if (ibeg < 0)
182  ibeg = 0;
183  double e = info.energyInWindow(ibeg, ibeg + nSamplesToAdd);
184 
185  // Pulse containment correction
186  {
187  double corrFactor = 1.0;
188  if (applyContainmentCorrection)
189  corrFactor = pulseCorr_.get(info.id(), nSamplesToAdd, phaseNs)->getCorrection(fc_ampl);
190  e *= corrFactor;
191  }
192 
193  return e;
194 }
const HcalPulseContainmentCorrection * get(const HcalDetId &detId, int toAdd, float fixedphase_ns)
unsigned soi() const
HcalPulseContainmentManager pulseCorr_
HcalDetId id() const
double energyInWindow(const unsigned begin, const unsigned end) const
float SimpleHBHEPhase1Algo::m0Time ( const HBHEChannelInfo info,
double  reconstructedCharge,
const HcalCalibrations calibs,
int  nSamplesToExamine 
) const
protected

Definition at line 196 of file SimpleHBHEPhase1Algo.cc.

References HcalTimeSlew::delay(), f, firstSampleShift_, hcalTimeSlew_delay_, SiStripPI::max, findQualityFiles::maxI, HBHEChannelInfo::MAXSAMPLES, HcalTimeSlew::Medium, hgc_digi::nSamples, HBHEChannelInfo::nSamples(), HBHEChannelInfo::peakEnergyTS(), HBHEChannelInfo::soi(), cscNeutronWriter_cfi::t0, reco::t2, ntuplemaker::time, HcalCalibrations::timecorr(), timeshift_ns_hbheho(), HBHEChannelInfo::tsEnergy(), and mitigatedMETSequence_cff::U.

Referenced by reconstruct().

200 {
201  float time = -9999.f; // historic value
202 
203  const unsigned nSamples = info.nSamples();
204  if (nSamples > 2U)
205  {
206  const int soi = info.soi();
207  int ibeg = soi + firstSampleShift_;
208  if (ibeg < 0)
209  ibeg = 0;
210  const int iend = ibeg + nSamplesToExamine;
211  unsigned maxI = info.peakEnergyTS(ibeg, iend);
212  if (maxI < HBHEChannelInfo::MAXSAMPLES)
213  {
214  if (!maxI)
215  maxI = 1U;
216  else if (maxI >= nSamples - 1U)
217  maxI = nSamples - 2U;
218 
219  // The remaining code in this scope emulates
220  // the historic algorithm
221  float t0 = info.tsEnergy(maxI - 1U);
222  float maxA = info.tsEnergy(maxI);
223  float t2 = info.tsEnergy(maxI + 1U);
224 
225  // Handle negative excursions by moving "zero"
226  float minA = t0;
227  if (maxA < minA) minA = maxA;
228  if (t2 < minA) minA=t2;
229  if (minA < 0.f) { maxA-=minA; t0-=minA; t2-=minA; }
230  float wpksamp = (t0 + maxA + t2);
231  if (wpksamp) wpksamp = (maxA + 2.f*t2) / wpksamp;
232  time = (maxI - soi)*25.f + timeshift_ns_hbheho(wpksamp);
233 
234  // Legacy QIE8 timing correction
235  time -= hcalTimeSlew_delay_->delay(std::max(1.0, fc_ampl), HcalTimeSlew::Medium);
236  // Time calibration
237  time -= calibs.timecorr();
238  }
239  }
240  return time;
241 }
const HcalTimeSlew * hcalTimeSlew_delay_
float timeshift_ns_hbheho(float wpksamp)
double delay(double fC, BiasSetting bias=Medium) const
Returns the amount (ns) by which a pulse of the given number of fC will be delayed by the timeslew ef...
Definition: HcalTimeSlew.cc:14
unsigned soi() const
static const unsigned MAXSAMPLES
auto const T2 &decltype(t1.eta()) t2
Definition: deltaR.h:16
double f[11][100]
constexpr size_t nSamples
double tsEnergy(const unsigned ts) const
double timecorr() const
get time correction factor
unsigned nSamples() const
unsigned peakEnergyTS(const unsigned begin, const unsigned end) const
HBHERecHit SimpleHBHEPhase1Algo::reconstruct ( const HBHEChannelInfo info,
const HcalRecoParam params,
const HcalCalibrations calibs,
bool  isRealData 
)
overridevirtual

Implements AbsHBHEPhase1Algo.

Definition at line 59 of file SimpleHBHEPhase1Algo.cc.

References HBHEChannelInfo::chargeInWindow(), vertices_cff::chi2, HcalRecoParam::correctForPhaseContainment(), HcalRecoParam::correctionPhaseNS(), corrFPC_, firstSampleShift_, HcalPulseShapes::getShape(), HBHEChannelInfo::hasTimeInfo(), HcalPhase1FlagLabels::HBHEPulseFitBit, hbminusCorrectionFactor(), hcalTimeSlew_delay_, hltOOTpuCorr_, HBHEChannelInfo::id(), HcalSpecialTimes::isSpecial(), m0Energy(), m0Time(), mahiOOTpuCorr_, HBHEChannelInfo::nSamples(), HcalDeterministicFit::phase1Apply(), PulseShapeFitOOTPileupCorrection::phase1Apply(), MahiFit::phase1Apply(), phaseNS_, psFitOOTpuCorr_, HBHEChannelInfo::recoShape(), HcalRecoParam::samplesToAdd(), samplesToAdd_, HBHERecHitAuxSetter::setAux(), HBHERecHit::setAuxEnergy(), HBHERecHit::setChiSquared(), CaloRecHit::setFlagField(), HBHERecHit::setRawEnergy(), HBHEChannelInfo::soi(), HBHEChannelInfo::soiRiseTime(), theHcalPulseShapes_, and timeShift_.

Referenced by isConfigurable().

63 {
64  HBHERecHit rh;
65 
66  const HcalDetId channelId(info.id());
67 
68  // Calculate "Method 0" quantities
69  float m0t = 0.f, m0E = 0.f;
70  {
71  int ibeg = static_cast<int>(info.soi()) + firstSampleShift_;
72  if (ibeg < 0)
73  ibeg = 0;
74  const int nSamplesToAdd = params ? params->samplesToAdd() : samplesToAdd_;
75  const double fc_ampl = info.chargeInWindow(ibeg, ibeg + nSamplesToAdd);
76  const bool applyContainment = params ? params->correctForPhaseContainment() : corrFPC_;
77  const float phasens = params ? params->correctionPhaseNS() : phaseNS_;
78  m0E = m0Energy(info, fc_ampl, applyContainment, phasens, nSamplesToAdd);
79  m0E *= hbminusCorrectionFactor(channelId, m0E, isData);
80  m0t = m0Time(info, fc_ampl, calibs, nSamplesToAdd);
81  }
82 
83  // Run "Method 2"
84  float m2t = 0.f, m2E = 0.f, chi2 = -1.f;
85  bool useTriple = false;
87  if (method2)
88  {
89  psFitOOTpuCorr_->setPulseShapeTemplate(theHcalPulseShapes_.getShape(info.recoShape()),
91  // "phase1Apply" call below sets m2E, m2t, useTriple, and chi2.
92  // These parameters are pased by non-const reference.
93  method2->phase1Apply(info, m2E, m2t, useTriple, chi2);
94  m2E *= hbminusCorrectionFactor(channelId, m2E, isData);
95  }
96 
97  // Run "Method 3"
98  float m3t = 0.f, m3E = 0.f;
99  const HcalDeterministicFit* method3 = hltOOTpuCorr_.get();
100  if (method3)
101  {
102  // "phase1Apply" sets m3E and m3t (pased by non-const reference)
103  method3->phase1Apply(info, m3E, m3t, hcalTimeSlew_delay_);
104  m3E *= hbminusCorrectionFactor(channelId, m3E, isData);
105  }
106 
107  // Run Mahi
108  float m4E = 0.f, m4chi2 = -1.f;
109  float m4T = 0.f;
110  bool m4UseTriple=false;
111 
112  const MahiFit* mahi = mahiOOTpuCorr_.get();
113 
114  if (mahi) {
116  mahi->phase1Apply(info,m4E,m4T,m4UseTriple,m4chi2);
117  m4E *= hbminusCorrectionFactor(channelId, m4E, isData);
118  }
119 
120  // Finally, construct the rechit
121  float rhE = m0E;
122  float rht = m0t;
123  float rhX = -1.f;
124  if (mahi)
125  {
126  rhE = m4E;
127  rht = m4T;
128  rhX = m4chi2;
129  }
130  else if (method2)
131  {
132  rhE = m2E;
133  rht = m2t;
134  rhX = chi2;
135  }
136  else if (method3)
137  {
138  rhE = m3E;
139  rht = m3t;
140  }
141  float tdcTime = info.soiRiseTime();
142  if (!HcalSpecialTimes::isSpecial(tdcTime))
143  tdcTime += timeShift_;
144  rh = HBHERecHit(channelId, rhE, rht, tdcTime);
145  rh.setRawEnergy(m0E);
146  rh.setAuxEnergy(m3E);
147  rh.setChiSquared(rhX);
148 
149  // Set rechit aux words
150  HBHERecHitAuxSetter::setAux(info, &rh);
151 
152  // Set some rechit flags (here, for Method 2/Mahi)
153  if (useTriple || m4UseTriple)
155 
156  return rh;
157 }
const Shape & getShape(int shapeType) const
const HcalTimeSlew * hcalTimeSlew_delay_
void phase1Apply(const HBHEChannelInfo &channelData, float &reconstructedEnergy, float &reconstructedTime, bool &useTriple, float &chi2) const
Definition: MahiFit.cc:40
bool hasTimeInfo() const
std::unique_ptr< HcalDeterministicFit > hltOOTpuCorr_
bool correctForPhaseContainment() const
Definition: HcalRecoParam.h:29
void setFlagField(uint32_t value, int base, int width=1)
Definition: CaloRecHit.cc:20
unsigned soi() const
void phase1Apply(const HBHEChannelInfo &channelData, float &reconstructedEnergy, float &reconstructedTime, bool &useTriple, float &chi2) const
HcalDetId id() const
bool isSpecial(const float t)
float m0Time(const HBHEChannelInfo &info, double reconstructedCharge, const HcalCalibrations &calibs, int nSamplesToExamine) const
float soiRiseTime() const
double chargeInWindow(const unsigned begin, const unsigned end) const
int recoShape() const
unsigned int samplesToAdd() const
Definition: HcalRecoParam.h:33
float hbminusCorrectionFactor(const HcalDetId &cell, float energy, bool isRealData) const
std::unique_ptr< MahiFit > mahiOOTpuCorr_
float correctionPhaseNS() const
Definition: HcalRecoParam.h:31
HcalPulseShapes theHcalPulseShapes_
float m0Energy(const HBHEChannelInfo &info, double reconstructedCharge, bool applyContainmentCorrection, double phaseNS, int nSamplesToAdd)
static void setAux(const HBHEChannelInfo &info, HBHERecHit *rechit)
void setAuxEnergy(const float en)
Definition: HBHERecHit.h:33
void setChiSquared(const float chi2)
Definition: HBHERecHit.h:27
unsigned nSamples() const
std::unique_ptr< PulseShapeFitOOTPileupCorrection > psFitOOTpuCorr_
void setRawEnergy(const float en)
Definition: HBHERecHit.h:30
void phase1Apply(const HBHEChannelInfo &channelData, float &reconstructedEnergy, float &reconstructedTime, const HcalTimeSlew *hcalTimeSlew_delay) const

Member Data Documentation

bool SimpleHBHEPhase1Algo::corrFPC_
private

Definition at line 100 of file SimpleHBHEPhase1Algo.h.

Referenced by isCorrectingForPhaseContainment(), and reconstruct().

int SimpleHBHEPhase1Algo::firstSampleShift_
private

Definition at line 95 of file SimpleHBHEPhase1Algo.h.

Referenced by getFirstSampleShift(), m0Energy(), m0Time(), and reconstruct().

const HcalTimeSlew* SimpleHBHEPhase1Algo::hcalTimeSlew_delay_

Definition at line 72 of file SimpleHBHEPhase1Algo.h.

Referenced by beginRun(), m0Time(), reconstruct(), and SimpleHBHEPhase1Algo().

std::unique_ptr<HcalDeterministicFit> SimpleHBHEPhase1Algo::hltOOTpuCorr_
private

Definition at line 106 of file SimpleHBHEPhase1Algo.h.

Referenced by reconstruct().

std::unique_ptr<MahiFit> SimpleHBHEPhase1Algo::mahiOOTpuCorr_
private

Definition at line 109 of file SimpleHBHEPhase1Algo.h.

Referenced by reconstruct().

float SimpleHBHEPhase1Algo::phaseNS_
private

Definition at line 97 of file SimpleHBHEPhase1Algo.h.

Referenced by getPhaseNS(), and reconstruct().

std::unique_ptr<PulseShapeFitOOTPileupCorrection> SimpleHBHEPhase1Algo::psFitOOTpuCorr_
private

Definition at line 103 of file SimpleHBHEPhase1Algo.h.

Referenced by reconstruct().

HcalPulseContainmentManager SimpleHBHEPhase1Algo::pulseCorr_
private

Definition at line 93 of file SimpleHBHEPhase1Algo.h.

Referenced by beginRun(), and m0Energy().

int SimpleHBHEPhase1Algo::runnum_
private

Definition at line 99 of file SimpleHBHEPhase1Algo.h.

Referenced by beginRun(), endRun(), getRunNumber(), and hbminusCorrectionFactor().

int SimpleHBHEPhase1Algo::samplesToAdd_
private

Definition at line 96 of file SimpleHBHEPhase1Algo.h.

Referenced by getSamplesToAdd(), and reconstruct().

HcalPulseShapes SimpleHBHEPhase1Algo::theHcalPulseShapes_
private

Definition at line 111 of file SimpleHBHEPhase1Algo.h.

Referenced by reconstruct().

float SimpleHBHEPhase1Algo::timeShift_
private

Definition at line 98 of file SimpleHBHEPhase1Algo.h.

Referenced by getTimeShift(), and reconstruct().