CMS 3D CMS Logo

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

#include <ESRecHitFitAlgo.h>

Public Member Functions

 ESRecHitFitAlgo ()
 
double * EvalAmplitude (const ESDataFrame &digi, double ped) const
 
EcalRecHit reconstruct (const ESDataFrame &digi) const
 
void setAngleCorrectionFactors (const ESAngleCorrectionFactors *ang)
 
void setChannelStatus (const ESChannelStatus *status)
 
void setESGain (const double &value)
 
void setIntercalibConstants (const ESIntercalibConstants *mips)
 
void setMIPGeV (const double &value)
 
void setPedestals (const ESPedestals *peds)
 
void setRatioCuts (const ESRecHitRatioCuts *ratioCuts)
 
 ~ESRecHitFitAlgo ()
 

Private Attributes

const ESAngleCorrectionFactorsang_
 
const ESChannelStatuschannelStatus_
 
TF1 * fit_
 
double gain_
 
double MIPGeV_
 
const ESIntercalibConstantsmips_
 
const ESPedestalspeds_
 
const ESRecHitRatioCutsratioCuts_
 

Detailed Description

Definition at line 14 of file ESRecHitFitAlgo.h.

Constructor & Destructor Documentation

◆ ESRecHitFitAlgo()

ESRecHitFitAlgo::ESRecHitFitAlgo ( )

Definition at line 24 of file ESRecHitFitAlgo.cc.

24  {
25  fit_ = new TF1("fit", fitf, -200, 200, 2);
26  fit_->SetParameters(50, 10);
27 }

References fit_, and fitf().

◆ ~ESRecHitFitAlgo()

ESRecHitFitAlgo::~ESRecHitFitAlgo ( )

Definition at line 29 of file ESRecHitFitAlgo.cc.

29 { delete fit_; }

References fit_.

Member Function Documentation

◆ EvalAmplitude()

double * ESRecHitFitAlgo::EvalAmplitude ( const ESDataFrame digi,
double  ped 
) const

Definition at line 31 of file ESRecHitFitAlgo.cc.

31  {
32  double* fitresults = new double[3]{};
33  double adc[3]{};
34  double tx[3]{};
35  tx[0] = -5.;
36  tx[1] = 20.;
37  tx[2] = 45.;
38 
39  for (int i = 0; i < digi.size(); i++)
40  adc[i] = digi.sample(i).adc() - ped;
41 
42  double status = 0;
43  if (adc[0] > 20)
44  status = 14;
45  if (adc[1] < 0 || adc[2] < 0)
46  status = 10;
47  if (adc[0] > adc[1] && adc[0] > adc[2])
48  status = 8;
49  if (adc[2] > adc[1] && adc[2] > adc[0])
50  status = 9;
51  double r12 = (adc[1] != 0) ? adc[0] / adc[1] : 99999;
52  double r23 = (adc[2] != 0) ? adc[1] / adc[2] : 99999;
53  if (r12 > ratioCuts_->getR12High())
54  status = 5;
55  if (r23 > ratioCuts_->getR23High())
56  status = 6;
57  if (r23 < ratioCuts_->getR23Low())
58  status = 7;
59 
60  if (int(status) == 0) {
61  double para[10]{};
62  TGraph* gr = new TGraph(3, tx, adc);
63  fit_->SetParameters(50, 10);
64  gr->Fit(fit_, "MQ");
65  fit_->GetParameters(para);
66  fitresults[0] = para[0];
67  fitresults[1] = para[1];
68  delete gr;
69 
70  if (adc[1] > 2800 && adc[2] > 2800)
71  status = 11;
72  else if (adc[1] > 2800)
73  status = 12;
74  else if (adc[2] > 2800)
75  status = 13;
76 
77  } else {
78  fitresults[0] = 0;
79  fitresults[1] = -999;
80  }
81 
82  fitresults[2] = status;
83 
84  return fitresults;
85 }

References ESSample::adc(), gpuClustering::adc, fit_, ESRecHitRatioCuts::getR12High(), ESRecHitRatioCuts::getR23High(), mps_fire::i, ratioCuts_, ESDataFrame::sample(), ESDataFrame::size(), and mps_update::status.

Referenced by reconstruct().

◆ reconstruct()

EcalRecHit ESRecHitFitAlgo::reconstruct ( const ESDataFrame digi) const

Definition at line 87 of file ESRecHitFitAlgo.cc.

87  {
88  ESPedestals::const_iterator it_ped = peds_->find(digi.id());
89 
92 
94 
95  double* results;
96 
97  results = EvalAmplitude(digi, it_ped->getMean());
98 
99  double energy = results[0];
100  double t0 = results[1];
101  int status = (int)results[2];
102  delete[] results;
103 
104  double mipCalib = (std::fabs(cos(*it_ang)) != 0.) ? (*it_mip) / fabs(cos(*it_ang)) : 0.;
105  energy *= (mipCalib != 0.) ? MIPGeV_ / mipCalib : 0.;
106 
107  LogDebug("ESRecHitFitAlgo") << "ESRecHitFitAlgo : reconstructed energy " << energy << " T0 " << t0;
108 
109  EcalRecHit rechit(digi.id(), energy, t0);
110 
111  if (it_status->getStatusCode() == 1) {
113  } else {
114  if (status == 0)
115  rechit.setFlag(EcalRecHit::kESGood);
116  else if (status == 5)
117  rechit.setFlag(EcalRecHit::kESBadRatioFor12);
118  else if (status == 6)
119  rechit.setFlag(EcalRecHit::kESBadRatioFor23Upper);
120  else if (status == 7)
121  rechit.setFlag(EcalRecHit::kESBadRatioFor23Lower);
122  else if (status == 8)
123  rechit.setFlag(EcalRecHit::kESTS1Largest);
124  else if (status == 9)
125  rechit.setFlag(EcalRecHit::kESTS3Largest);
126  else if (status == 10)
127  rechit.setFlag(EcalRecHit::kESTS3Negative);
128  else if (status == 11)
129  rechit.setFlag(EcalRecHit::kESSaturated);
130  else if (status == 12)
131  rechit.setFlag(EcalRecHit::kESTS2Saturated);
132  else if (status == 13)
133  rechit.setFlag(EcalRecHit::kESTS3Saturated);
134  else if (status == 14)
135  rechit.setFlag(EcalRecHit::kESTS13Sigmas);
136  }
137 
138  return rechit;
139 }

References ang_, channelStatus_, funct::cos(), HCALHighEnergyHPDFilter_cfi::energy, EvalAmplitude(), ESCondObjectContainer< T >::find(), ESCondObjectContainer< T >::getMap(), ESDataFrame::id(), createfilelist::int, EcalRecHit::kESBadRatioFor12, EcalRecHit::kESBadRatioFor23Lower, EcalRecHit::kESBadRatioFor23Upper, EcalRecHit::kESDead, EcalRecHit::kESGood, EcalRecHit::kESSaturated, EcalRecHit::kESTS13Sigmas, EcalRecHit::kESTS1Largest, EcalRecHit::kESTS2Saturated, EcalRecHit::kESTS3Largest, EcalRecHit::kESTS3Negative, EcalRecHit::kESTS3Saturated, LogDebug, MIPGeV_, mips_, peds_, bookConverter::results, EcalRecHit::setFlag(), mps_update::status, and FrontierCondition_GT_autoExpress_cfi::t0.

Referenced by ESRecHitWorker::run().

◆ setAngleCorrectionFactors()

void ESRecHitFitAlgo::setAngleCorrectionFactors ( const ESAngleCorrectionFactors ang)
inline

Definition at line 25 of file ESRecHitFitAlgo.h.

25 { ang_ = ang; }

References ang_.

Referenced by ESRecHitWorker::set().

◆ setChannelStatus()

void ESRecHitFitAlgo::setChannelStatus ( const ESChannelStatus status)
inline

Definition at line 23 of file ESRecHitFitAlgo.h.

References channelStatus_, and mps_update::status.

Referenced by ESRecHitWorker::set().

◆ setESGain()

void ESRecHitFitAlgo::setESGain ( const double &  value)
inline

Definition at line 19 of file ESRecHitFitAlgo.h.

19 { gain_ = value; }

References gain_, and relativeConstraints::value.

Referenced by ESRecHitWorker::set().

◆ setIntercalibConstants()

void ESRecHitFitAlgo::setIntercalibConstants ( const ESIntercalibConstants mips)
inline

Definition at line 22 of file ESRecHitFitAlgo.h.

22 { mips_ = mips; }

References mips_.

Referenced by ESRecHitWorker::set().

◆ setMIPGeV()

void ESRecHitFitAlgo::setMIPGeV ( const double &  value)
inline

Definition at line 20 of file ESRecHitFitAlgo.h.

20 { MIPGeV_ = value; }

References MIPGeV_, and relativeConstraints::value.

Referenced by ESRecHitWorker::set().

◆ setPedestals()

void ESRecHitFitAlgo::setPedestals ( const ESPedestals peds)
inline

Definition at line 21 of file ESRecHitFitAlgo.h.

21 { peds_ = peds; }

References peds_.

Referenced by ESRecHitWorker::set().

◆ setRatioCuts()

void ESRecHitFitAlgo::setRatioCuts ( const ESRecHitRatioCuts ratioCuts)
inline

Definition at line 24 of file ESRecHitFitAlgo.h.

24 { ratioCuts_ = ratioCuts; }

References ratioCuts_.

Referenced by ESRecHitWorker::set().

Member Data Documentation

◆ ang_

const ESAngleCorrectionFactors* ESRecHitFitAlgo::ang_
private

Definition at line 36 of file ESRecHitFitAlgo.h.

Referenced by reconstruct(), and setAngleCorrectionFactors().

◆ channelStatus_

const ESChannelStatus* ESRecHitFitAlgo::channelStatus_
private

Definition at line 34 of file ESRecHitFitAlgo.h.

Referenced by reconstruct(), and setChannelStatus().

◆ fit_

TF1* ESRecHitFitAlgo::fit_
private

Definition at line 30 of file ESRecHitFitAlgo.h.

Referenced by ESRecHitFitAlgo(), EvalAmplitude(), and ~ESRecHitFitAlgo().

◆ gain_

double ESRecHitFitAlgo::gain_
private

Definition at line 31 of file ESRecHitFitAlgo.h.

Referenced by setESGain().

◆ MIPGeV_

double ESRecHitFitAlgo::MIPGeV_
private

Definition at line 37 of file ESRecHitFitAlgo.h.

Referenced by reconstruct(), and setMIPGeV().

◆ mips_

const ESIntercalibConstants* ESRecHitFitAlgo::mips_
private

Definition at line 33 of file ESRecHitFitAlgo.h.

Referenced by reconstruct(), and setIntercalibConstants().

◆ peds_

const ESPedestals* ESRecHitFitAlgo::peds_
private

Definition at line 32 of file ESRecHitFitAlgo.h.

Referenced by reconstruct(), and setPedestals().

◆ ratioCuts_

const ESRecHitRatioCuts* ESRecHitFitAlgo::ratioCuts_
private

Definition at line 35 of file ESRecHitFitAlgo.h.

Referenced by EvalAmplitude(), and setRatioCuts().

ESRecHitFitAlgo::fit_
TF1 * fit_
Definition: ESRecHitFitAlgo.h:30
ESRecHitFitAlgo::mips_
const ESIntercalibConstants * mips_
Definition: ESRecHitFitAlgo.h:33
EcalRecHit
Definition: EcalRecHit.h:15
mps_fire.i
i
Definition: mps_fire.py:428
EcalRecHit::kESBadRatioFor23Lower
Definition: EcalRecHit.h:52
ESSample::adc
int adc() const
get the ADC sample (singed 16 bits)
Definition: ESSample.h:16
mps_update.status
status
Definition: mps_update.py:68
gpuClustering::adc
uint16_t *__restrict__ uint16_t const *__restrict__ adc
Definition: gpuClusterChargeCut.h:20
ESCondObjectContainer::find
const_iterator find(uint32_t rawId) const
Definition: ESCondObjectContainer.h:33
ESRecHitFitAlgo::channelStatus_
const ESChannelStatus * channelStatus_
Definition: ESRecHitFitAlgo.h:34
ESDataFrame::size
int size() const
Definition: ESDataFrame.h:21
EcalRecHit::kESTS3Saturated
Definition: EcalRecHit.h:58
bookConverter.results
results
Definition: bookConverter.py:144
ESRecHitRatioCuts::getR12High
float getR12High() const
Definition: ESRecHitRatioCuts.h:19
ESRecHitFitAlgo::gain_
double gain_
Definition: ESRecHitFitAlgo.h:31
ESDataFrame::sample
const ESSample & sample(int i) const
Definition: ESDataFrame.h:24
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
FrontierCondition_GT_autoExpress_cfi.t0
t0
Definition: FrontierCondition_GT_autoExpress_cfi.py:149
ESCondObjectContainer< ESPedestal >::const_iterator
std::vector< Item >::const_iterator const_iterator
Definition: ESCondObjectContainer.h:17
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
ESRecHitFitAlgo::ang_
const ESAngleCorrectionFactors * ang_
Definition: ESRecHitFitAlgo.h:36
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:233
EcalRecHit::kESSaturated
Definition: EcalRecHit.h:56
EcalRecHit::kESTS1Largest
Definition: EcalRecHit.h:53
EcalRecHit::kESTS3Negative
Definition: EcalRecHit.h:55
EcalRecHit::kESTS13Sigmas
Definition: EcalRecHit.h:59
createfilelist.int
int
Definition: createfilelist.py:10
EcalRecHit::kESBadRatioFor12
Definition: EcalRecHit.h:50
ESRecHitFitAlgo::peds_
const ESPedestals * peds_
Definition: ESRecHitFitAlgo.h:32
ESRecHitFitAlgo::EvalAmplitude
double * EvalAmplitude(const ESDataFrame &digi, double ped) const
Definition: ESRecHitFitAlgo.cc:31
ESCondObjectContainer::getMap
const self & getMap() const
Definition: ESCondObjectContainer.h:41
ESRecHitRatioCuts::getR23High
float getR23High() const
Definition: ESRecHitRatioCuts.h:21
relativeConstraints.value
value
Definition: relativeConstraints.py:53
ESRecHitFitAlgo::MIPGeV_
double MIPGeV_
Definition: ESRecHitFitAlgo.h:37
ESRecHitFitAlgo::ratioCuts_
const ESRecHitRatioCuts * ratioCuts_
Definition: ESRecHitFitAlgo.h:35
EcalRecHit::setFlag
void setFlag(int flag)
set the flags (from Flags or ESFlags)
Definition: EcalRecHit.h:183
EcalRecHit::kESTS2Saturated
Definition: EcalRecHit.h:57
fitf
Double_t fitf(Double_t *x, Double_t *par)
Definition: ESRecHitFitAlgo.cc:11
EcalRecHit::kESBadRatioFor23Upper
Definition: EcalRecHit.h:51
EcalRecHit::kESDead
Definition: EcalRecHit.h:46
EcalRecHit::kESGood
Definition: EcalRecHit.h:45
ESDataFrame::id
const ESDetId & id() const
Definition: ESDataFrame.h:19
EcalRecHit::kESTS3Largest
Definition: EcalRecHit.h:54