CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ESRecHitFitAlgo.cc
Go to the documentation of this file.
3 
4 #include "TMath.h"
5 #include "TGraph.h"
6 
7 #include <iostream>
8 #include <cmath>
9 
10 // fit function
11 Double_t fitf(Double_t *x, Double_t *par) {
12 
13  Double_t wc = 0.07291;
14  Double_t n = 1.798; // n-1 (in fact)
15  Double_t v1 = pow(wc/n*(x[0]-par[1]), n);
16  Double_t v2 = TMath::Exp(n-wc*(x[0]-par[1]));
17  Double_t v = par[0]*v1*v2;
18 
19  if (x[0] < par[1]) v = 0;
20 
21  return v;
22 }
23 
25 
26  fit_ = new TF1("fit", fitf, -200, 200, 2);
27  fit_->SetParameters(50, 10);
28 
29 }
30 
32  delete fit_;
33 }
34 
35 double* ESRecHitFitAlgo::EvalAmplitude(const ESDataFrame& digi, double ped) const {
36 
37  double *fitresults = new double[3];
38  double adc[3];
39  double tx[3];
40  tx[0] = -5.;
41  tx[1] = 20.;
42  tx[2] = 45.;
43 
44  for (int i=0; i<digi.size(); i++)
45  adc[i] = digi.sample(i).adc() - ped;
46 
47  double status = 0;
48  if (adc[0] > 20) status = 14;
49  if (adc[1] < 0 || adc[2] < 0) status = 10;
50  if (adc[0] > adc[1] && adc[0] > adc[2]) status = 8;
51  if (adc[2] > adc[1] && adc[2] > adc[0]) status = 9;
52  double r12 = (adc[1] != 0) ? adc[0]/adc[1] : 99999;
53  double r23 = (adc[2] != 0) ? adc[1]/adc[2] : 99999;
54  if (r12 > ratioCuts_->getR12High()) status = 5;
55  if (r23 > ratioCuts_->getR23High()) status = 6;
56  if (r23 < ratioCuts_->getR23Low()) status = 7;
57 
58  if (int(status) == 0) {
59  double para[10];
60  TGraph *gr = new TGraph(3, tx, adc);
61  fit_->SetParameters(50, 10);
62  gr->Fit("fit", "MQ");
63  fit_->GetParameters(para);
64  fitresults[0] = para[0];
65  fitresults[1] = para[1];
66  delete gr;
67 
68  if (adc[1] > 2800 && adc[2] > 2800) status = 11;
69  else if (adc[1] > 2800) status = 12;
70  else if (adc[2] > 2800) status = 13;
71 
72  } else {
73  fitresults[0] = 0;
74  fitresults[1] = -999;
75  }
76 
77  fitresults[2] = status;
78 
79  return fitresults;
80 }
81 
83 
84  ESPedestals::const_iterator it_ped = peds_->find(digi.id());
85 
88 
90 
91  double* results;
92 
93  results = EvalAmplitude(digi, it_ped->getMean());
94 
95  double energy = results[0];
96  double t0 = results[1];
97  int status = (int) results[2];
98  delete[] results;
99 
100  double mipCalib = (std::fabs(cos(*it_ang)) != 0.) ? (*it_mip)/fabs(cos(*it_ang)) : 0.;
101  energy *= (mipCalib != 0.) ? MIPGeV_/mipCalib : 0.;
102 
103  LogDebug("ESRecHitFitAlgo") << "ESRecHitFitAlgo : reconstructed energy "<<energy<<" T0 "<<t0;
104 
105  EcalRecHit rechit(digi.id(), energy, t0);
106 
107  if (it_status->getStatusCode() == 1) {
109  } else {
110  if (status == 0)
111  rechit.setFlag(EcalRecHit::kESGood);
112  else if (status == 5)
113  rechit.setFlag(EcalRecHit::kESBadRatioFor12);
114  else if (status == 6)
115  rechit.setFlag(EcalRecHit::kESBadRatioFor23Upper);
116  else if (status == 7)
117  rechit.setFlag(EcalRecHit::kESBadRatioFor23Lower);
118  else if (status == 8)
119  rechit.setFlag(EcalRecHit::kESTS1Largest);
120  else if (status == 9)
121  rechit.setFlag(EcalRecHit::kESTS3Largest);
122  else if (status == 10)
123  rechit.setFlag(EcalRecHit::kESTS3Negative);
124  else if (status == 11)
125  rechit.setFlag(EcalRecHit::kESSaturated);
126  else if (status == 12)
127  rechit.setFlag(EcalRecHit::kESTS2Saturated);
128  else if (status == 13)
129  rechit.setFlag(EcalRecHit::kESTS3Saturated);
130  else if (status == 14)
131  rechit.setFlag(EcalRecHit::kESTS13Sigmas);
132  }
133 
134  return rechit;
135 }
136 
#define LogDebug(id)
int adc(sample_type sample)
get the ADC sample (12 bits)
int i
Definition: DBlmapReader.cc:9
const ESChannelStatus * channelStatus_
double fitf(double *x, double *par)
Definition: ESTimingTask.cc:30
const ESDetId & id() const
Definition: ESDataFrame.h:21
const self & getMap() const
void setFlag(int flag)
set the flags (from Flags or ESFlags)
Definition: EcalRecHit.h:99
int size() const
Definition: ESDataFrame.h:23
const ESRecHitRatioCuts * ratioCuts_
const ESAngleCorrectionFactors * ang_
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
const_iterator find(uint32_t rawId) const
const ESSample & sample(int i) const
Definition: ESDataFrame.h:26
EcalRecHit reconstruct(const ESDataFrame &digi) const
std::vector< Item >::const_iterator const_iterator
float getR12High() const
float getR23High() const
double * EvalAmplitude(const ESDataFrame &digi, double ped) const
int adc() const
get the ADC sample (singed 16 bits)
Definition: ESSample.h:18
Definition: DDAxes.h:10
tuple status
Definition: ntuplemaker.py:245
const ESIntercalibConstants * mips_
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
const ESPedestals * peds_