CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ESRecHitSimAlgo.cc
Go to the documentation of this file.
3 
4 #include <cmath>
5 #include <vdt/vdtMath.h>
6 
7 #include<cstdlib>
8 #include<cstring>
9 #include<cassert>
10 #include <iostream>
11 
13 
14  float energy = 0;
15  float adc[3];
16  float pw[3];
17  pw[0] = w0_;
18  pw[1] = w1_;
19  pw[2] = w2_;
20 
21  for (int i=0; i<digi.size(); i++) {
22  energy += pw[i]*(digi.sample(i).adc()-ped);
23  LogDebug("ESRecHitSimAlgo") << "ESRecHitSimAlgo : Digi "<<i<<" ADC counts "<<digi.sample(i).adc()<<" Ped "<<ped;
24  //std::cout<<i<<" "<<digi.sample(i).adc()<<" "<<ped<<" "<<pw[i]<<std::endl;
25  adc[i] = digi.sample(i).adc() - ped;
26  }
27 
29  if (adc[0] > 20.f) status = EcalRecHit::kESTS13Sigmas; // 14;
30  if (adc[1] <= 0 || adc[2] <= 0) status = EcalRecHit::kESTS3Negative; // 10;
31  if (adc[0] > adc[1] && adc[0] > adc[2]) status = EcalRecHit::kESTS1Largest; // 8;
32  if (adc[2] > adc[1] && adc[2] > adc[0]) status = EcalRecHit::kESTS3Largest; // 9;
33  auto r12 = (adc[1] != 0) ? adc[0]/adc[1] : 99.f;
34  auto r23 = (adc[2] != 0) ? adc[1]/adc[2] : 99.f;
35  if (r12 > ratioCuts_->getR12High()) status = EcalRecHit::kESBadRatioFor12;// 5;
36  if (r23 > ratioCuts_->getR23High()) status = EcalRecHit::kESBadRatioFor23Upper; // 6;
37  if (r23 < ratioCuts_->getR23Low()) status = EcalRecHit::kESBadRatioFor23Lower; // 7
38 
39  auto A1 = adc[1];
40  auto A2 = adc[2];
41 
42  // t0 from analytical formula:
43  constexpr float n = 1.798;
44  constexpr float w = 0.07291;
45  constexpr float DeltaT = 25.;
46  auto aaa = (A2 > 0 && A1 > 0) ? std::log(A2/A1)/n : 20.f; // if A1=0, t0=20
47  constexpr float bbb = w/n*DeltaT;
48  auto ccc= std::exp(aaa+bbb);
49 
50  auto t0 = (2.f-ccc)/(1.f-ccc) * DeltaT - 5.f;
51 
52  // A from analytical formula:
53  constexpr float t1 = 20.;
54  constexpr float A_1 = 1./( std::pow(w/n*(t1),n) * std::exp(n-w*(t1)) );
55  auto AA1 = A1 * A_1 ;
56 
57  if (adc[1] > 2800.f && adc[2] > 2800.f) status = EcalRecHit::kESSaturated;
58  else if (adc[1] > 2800.f) status = EcalRecHit::kESTS2Saturated;
59  else if (adc[2] > 2800.f) status = EcalRecHit::kESTS3Saturated;
60 
61  results[0] = energy; // energy with weight method
62  results[1] = t0; // timing
63  results[2] = AA1; // energy with analytic method
64 
65  return status;
66 
67 }
68 
70 
71 
72  auto ind = digi.id().hashedIndex();
73 
74  auto const & ped = peds_->preshower(ind);
75  auto const & mip = mips_->getMap().preshower(ind);
76  auto const & ang = ang_->getMap().preshower(ind);
77  auto const & statusCh = channelStatus_->getMap().preshower(ind);
78 
79  float results[3];
80 
81  auto status = evalAmplitude(results, digi, ped.getMean());
82 
83  auto energy = results[0];
84  auto t0 = results[1];
85  auto otenergy = results[2] * 1000000.f; // set out-of-time energy to keV
86 
87 
88  auto mipCalib = (mip != 0.f) ? MIPGeV_*std::abs(vdt::fast_cosf(ang))/(mip) : 0.f;
89  energy *= mipCalib;
90  otenergy *= mipCalib;
91 
92  LogDebug("ESRecHitSimAlgo") << "ESRecHitSimAlgo : reconstructed energy "<<energy;
93 
94  EcalRecHit rechit(digi.id(), energy, t0);
95  rechit.setOutOfTimeEnergy(otenergy);
96 
97  rechit.setFlag(statusCh.getStatusCode() == 1 ? EcalRecHit::kESDead : status);
98 
99  return rechit;
100 
101 }
102 
103 /*
104 
105  auto oldHit = oldreconstruct(digi);
106 
107  assert(rechit.recoFlag()==oldHit.recoFlag());
108 
109  if (oldHit.energy()>0)
110  std::cout << "ESd " << digi.id() <<" : "
111  << rechit.energy()<<"," << oldHit.energy() << " "
112  << rechit.time()<<"," << oldHit.time() << " "
113  << rechit.outOfTimeEnergy()<<"," << oldHit.outOfTimeEnergy() << std::endl;
114 
115 
116  auto bitdiff = [](float a, float b)->int {
117  int ia, ib;
118  memcpy(&ia,&a,4);
119  memcpy(&ib,&b,4);
120  return std::abs(ia-ib);
121  };
122 
123  auto d0 = bitdiff( rechit.energy(), rechit.energy() );
124  auto d1 = bitdiff( rechit.time(), rechit.time() );
125  auto d2 = bitdiff(rechit.outOfTimeEnergy(), oldHit.outOfTimeEnergy());
126 
127  struct aMax {
128  ~aMax() { std::cout << "\n\nmax deviation " << m << " " << mp << std::endl; }
129  int m=0;
130  int mp=0;
131  };
132 
133  static aMax am;
134  am.m = std::max(am.m,std::max(d0,std::max(d1,d2)));
135  if (oldHit.energy()>0) am.mp = std::max(am.m,std::max(d0,std::max(d1,d2)));
136 
137 
138  return rechit;
139 }
140 */
141 
142 //
143 
144 double* ESRecHitSimAlgo::oldEvalAmplitude(const ESDataFrame& digi, const double& ped, const double& w0, const double& w1, const double& w2) const {
145 
146  double *results = new double[4];
147  float energy = 0;
148  double adc[3];
149  float pw[3];
150  pw[0] = w0;
151  pw[1] = w1;
152  pw[2] = w2;
153 
154  for (int i=0; i<digi.size(); i++) {
155  energy += pw[i]*(digi.sample(i).adc()-ped);
156  LogDebug("ESRecHitSimAlgo") << "ESRecHitSimAlgo : Digi "<<i<<" ADC counts "<<digi.sample(i).adc()<<" Ped "<<ped;
157  //std::cout<<i<<" "<<digi.sample(i).adc()<<" "<<ped<<" "<<pw[i]<<std::endl;
158  adc[i] = digi.sample(i).adc() - ped;
159  }
160 
161  double status = 0;
162  if (adc[0] > 20) status = 14;
163  if (adc[1] <= 0 || adc[2] <= 0) status = 10;
164  if (adc[0] > adc[1] && adc[0] > adc[2]) status = 8;
165  if (adc[2] > adc[1] && adc[2] > adc[0]) status = 9;
166  double r12 = (adc[1] != 0) ? adc[0]/adc[1] : 99;
167  double r23 = (adc[2] != 0) ? adc[1]/adc[2] : 99;
168  if (r12 > ratioCuts_->getR12High()) status = 5;
169  if (r23 > ratioCuts_->getR23High()) status = 6;
170  if (r23 < ratioCuts_->getR23Low()) status = 7;
171 
172  double A1 = adc[1];
173  double A2 = adc[2];
174 
175  // t0 from analytical formula:
176  double n = 1.798;
177  double w = 0.07291;
178  double DeltaT = 25.;
179  double aaa = (A2 > 0 && A1 > 0) ? log(A2/A1)/n : 20.; // if A1=0, t0=20
180  double bbb = w/n*DeltaT;
181  double ccc= exp(aaa+bbb);
182 
183  double t0 = (2.-ccc)/(1.-ccc) * DeltaT - 5;
184 
185  // A from analytical formula:
186  double t1 = 20.;
187  double A_1 = pow(w/n*(t1),n) * exp(n-w*(t1));
188  double AA1 = (A_1 != 0.) ? A1 / A_1 : 0.;
189 
190  if (adc[1] > 2800 && adc[2] > 2800) status = 11;
191  else if (adc[1] > 2800) status = 12;
192  else if (adc[2] > 2800) status = 13;
193 
194  results[0] = energy; // energy with weight method
195  results[1] = t0; // timing
196  results[2] = status; // hit status
197  results[3] = AA1; // energy with analytic method
198 
199  return results;
200 }
201 
203 
204  ESPedestals::const_iterator it_ped = peds_->find(digi.id());
205 
208 
210 
211  double* results;
212 
213  results = oldEvalAmplitude(digi, it_ped->getMean(), w0_, w1_, w2_);
214 
215  double energy = results[0];
216  double t0 = results[1];
217  int status = (int) results[2];
218  double otenergy = results[3] * 1000000.; // set out-of-time energy to keV
219  delete[] results;
220 
221  double mipCalib = (fabs(cos(*it_ang)) != 0.) ? (*it_mip)/fabs(cos(*it_ang)) : 0.;
222  energy *= (mipCalib != 0.) ? MIPGeV_/mipCalib : 0.;
223  otenergy *= (mipCalib != 0.) ? MIPGeV_/mipCalib : 0.;
224 
225  LogDebug("ESRecHitSimAlgo") << "ESRecHitSimAlgo : reconstructed energy "<<energy;
226 
227  EcalRecHit rechit(digi.id(), energy, t0);
228  rechit.setOutOfTimeEnergy(otenergy);
229 
230  if (it_status->getStatusCode() == 1) {
231  rechit.setFlag(EcalRecHit::kESDead);
232  } else {
233  if (status == 0)
234  rechit.setFlag(EcalRecHit::kESGood);
235  else if (status == 5)
236  rechit.setFlag(EcalRecHit::kESBadRatioFor12);
237  else if (status == 6)
238  rechit.setFlag(EcalRecHit::kESBadRatioFor23Upper);
239  else if (status == 7)
240  rechit.setFlag(EcalRecHit::kESBadRatioFor23Lower);
241  else if (status == 8)
242  rechit.setFlag(EcalRecHit::kESTS1Largest);
243  else if (status == 9)
244  rechit.setFlag(EcalRecHit::kESTS3Largest);
245  else if (status == 10)
246  rechit.setFlag(EcalRecHit::kESTS3Negative);
247  else if (status == 11)
248  rechit.setFlag(EcalRecHit::kESSaturated);
249  else if (status == 12)
250  rechit.setFlag(EcalRecHit::kESTS2Saturated);
251  else if (status == 13)
252  rechit.setFlag(EcalRecHit::kESTS3Saturated);
253  else if (status == 14)
254  rechit.setFlag(EcalRecHit::kESTS13Sigmas);
255  }
256 
257  return rechit;
258 }
259 
#define LogDebug(id)
int adc(sample_type sample)
get the ADC sample (12 bits)
int i
Definition: DBlmapReader.cc:9
EcalRecHit reconstruct(const ESDataFrame &digi) const
common ppss p3p6s2 common epss epspn46 common const1 w2
Definition: inclppp.h:1
const double w
Definition: UKUtility.cc:23
const ESDetId & id() const
Definition: ESDataFrame.h:21
const self & getMap() const
int size() const
Definition: ESDataFrame.h:23
const ESPedestals * peds_
#define constexpr
const ESIntercalibConstants * mips_
const ESChannelStatus * channelStatus_
EcalRecHit oldreconstruct(const ESDataFrame &digi) const
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
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
const Item & preshower(size_t hashedIndex) const
double * oldEvalAmplitude(const ESDataFrame &digi, const double &ped, const double &w0, const double &w1, const double &w2) const
const ESRecHitRatioCuts * ratioCuts_
const ESAngleCorrectionFactors * ang_
std::vector< Item >::const_iterator const_iterator
float getR12High() const
float getR23High() const
void setOutOfTimeEnergy(float energy)
Definition: EcalRecHit.h:137
EcalRecHit::ESFlags evalAmplitude(float *result, const ESDataFrame &digi, float ped) const
int adc() const
get the ADC sample (singed 16 bits)
Definition: ESSample.h:18
tuple status
Definition: ntuplemaker.py:245
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
int hashedIndex() const
get a compact index for arrays [TODO: NEEDS WORK]
Definition: ESDetId.cc:69