CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PulseShapeFunctor.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <cmath>
3 #include <climits>
6 
7 namespace FitterFuncs {
8 
9  //Decalare the Pulse object take it in from Hcal and set some options
11  bool iPedestalConstraint,
12  bool iTimeConstraint,
13  bool iAddPulseJitter,
14  double iPulseJitter,
15  double iTimeMean,
16  double iPedMean,
17  unsigned nSamplesToFit)
18  : cntNANinfit(0),
19  acc25nsVec_(hcal::constants::maxPSshapeBin),
20  diff25nsItvlVec_(hcal::constants::maxPSshapeBin),
21  accVarLenIdxZEROVec_(hcal::constants::nsPerBX),
22  diffVarItvlIdxZEROVec_(hcal::constants::nsPerBX),
23  accVarLenIdxMinusOneVec_(hcal::constants::nsPerBX),
24  diffVarItvlIdxMinusOneVec_(hcal::constants::nsPerBX) {
25  //The raw pulse
26  for (int i = 0; i < hcal::constants::maxPSshapeBin; ++i)
27  pulse_hist[i] = pulse(i);
28  // Accumulate 25ns for each starting point of 0, 1, 2, 3...
29  for (int i = 0; i < hcal::constants::maxPSshapeBin; ++i) {
30  for (int j = i; j < i + hcal::constants::nsPerBX; ++j) { //sum over hcal::constants::nsPerBXns from point i
31  acc25nsVec_[i] +=
32  (j < hcal::constants::maxPSshapeBin ? pulse_hist[j] : pulse_hist[hcal::constants::maxPSshapeBin - 1]);
33  }
34  diff25nsItvlVec_[i] = (i + hcal::constants::nsPerBX < hcal::constants::maxPSshapeBin
36  : pulse_hist[hcal::constants::maxPSshapeBin - 1] - pulse_hist[i]);
37  }
38  // Accumulate different ns for starting point of index either 0 or -1
39  for (int i = 0; i < hcal::constants::nsPerBX; ++i) {
40  if (i == 0) {
43  } else {
45  accVarLenIdxMinusOneVec_[i] = accVarLenIdxMinusOneVec_[i - 1] + pulse_hist[i - 1];
46  }
48  diffVarItvlIdxMinusOneVec_[i] = pulse_hist[i] - pulse_hist[0];
49  }
50  for (int i = 0; i < hcal::constants::maxSamples; i++) {
51  psFit_x[i] = 0;
52  psFit_y[i] = 0;
53  psFit_erry[i] = 1.;
54  psFit_erry2[i] = 1.;
55  psFit_slew[i] = 0.;
56  }
57  //Constraints
58  pedestalConstraint_ = iPedestalConstraint;
59  timeConstraint_ = iTimeConstraint;
60  addPulseJitter_ = iAddPulseJitter;
61  pulseJitter_ = iPulseJitter * iPulseJitter;
62 
63  // for M2
64  timeMean_ = iTimeMean;
65  pedMean_ = iPedMean;
66  timeShift_ = 100.;
67  timeShift_ += 12.5; //we are trying to get BX
68 
69  nSamplesToFit_ = nSamplesToFit;
70  }
71 
72  void PulseShapeFunctor::funcShape(std::array<double, hcal::constants::maxSamples> &ntmpbin,
73  const double pulseTime,
74  const double pulseHeight,
75  const double slew,
76  bool scalePulse) {
77  // pulse shape components over a range of time 0 ns to 255 ns in 1 ns steps
78  constexpr int ns_per_bx = hcal::constants::nsPerBX;
79  //Get the starting time
80  int i_start = (-hcal::constants::iniTimeShift - pulseTime - slew > 0
81  ? 0
82  : (int)std::abs(-hcal::constants::iniTimeShift - pulseTime - slew) + 1);
83  double offset_start = i_start - hcal::constants::iniTimeShift - pulseTime -
84  slew; //-199-2*pars[0]-2.*slew (for pars[0] > 98.5) or just -98.5-pars[0]-slew;
85  // zeroing output binned pulse shape
86  ntmpbin.fill(0.0f);
87 
88  if (edm::isNotFinite(offset_start)) { //Check for nan
89  ++cntNANinfit;
90  } else {
91  if (offset_start == 1.0) {
92  offset_start = 0.;
93  i_start -= 1;
94  } //Deal with boundary
95 
96  const int bin_start = (int)offset_start; //bin off to integer
97  const int bin_0_start = (offset_start < bin_start + 0.5 ? bin_start - 1 : bin_start); //Round it
98  const int iTS_start = i_start / ns_per_bx; //Time Slice for time shift
99  const int distTo25ns_start = ns_per_bx - 1 - i_start % ns_per_bx; //Delta ns
100  const double factor = offset_start - bin_0_start - 0.5; //Small correction?
101 
102  //Build the new pulse
103  ntmpbin[iTS_start] =
104  (bin_0_start == -1
105  ? // Initial bin (I'm assuming this is ok)
106  accVarLenIdxMinusOneVec_[distTo25ns_start] + factor * diffVarItvlIdxMinusOneVec_[distTo25ns_start]
107  : accVarLenIdxZEROVec_[distTo25ns_start] + factor * diffVarItvlIdxZEROVec_[distTo25ns_start]);
108  //Fill the rest of the bins
109  for (int iTS = iTS_start + 1; iTS < hcal::constants::maxSamples; ++iTS) {
110  int bin_idx = distTo25ns_start + 1 + (iTS - iTS_start - 1) * ns_per_bx + bin_0_start;
111  ntmpbin[iTS] = acc25nsVec_[bin_idx] + factor * diff25nsItvlVec_[bin_idx];
112  }
113  //Scale the pulse
114  if (scalePulse) {
115  for (int i = iTS_start; i < hcal::constants::maxSamples; ++i) {
116  ntmpbin[i] *= pulseHeight;
117  }
118  }
119  }
120 
121  return;
122  }
123 
125 
126  void PulseShapeFunctor::EvalPulse(const float *pars) {
127  int time = (pars[0] + timeShift_ - timeMean_) * hcal::constants::invertnsPerBx;
128  float dummyPulseHeight = 0.f;
129  funcShape(pulse_shape_, pars[0], dummyPulseHeight, psFit_slew[time], false);
130  return;
131  }
132 
133  double PulseShapeFunctor::EvalPulseM2(const double *pars, const unsigned nPars) {
134  unsigned i = 0, j = 0;
135 
136  const double pedestal = pars[nPars - 1];
137 
138  //Stop crashes
139  for (i = 0; i < nPars; ++i)
140  if (edm::isNotFinite(pars[i])) {
141  ++cntNANinfit;
142  return 1e10;
143  }
144 
145  //calculate chisquare
146  double chisq = 0;
147  const unsigned parBy2 = (nPars - 1) / 2;
148  // std::array<float,hcal::constants::maxSamples> pulse_shape_;
149 
150  if (addPulseJitter_) {
151  int time = (pars[0] + timeShift_ - timeMean_) * hcal::constants::invertnsPerBx;
152  //Interpolate the fit (Quickly)
153  funcShape(pulse_shape_, pars[0], pars[1], psFit_slew[time], true);
154  for (j = 0; j < nSamplesToFit_; ++j) {
157  }
158 
159  for (i = 1; i < parBy2; ++i) {
160  time = (pars[i * 2] + timeShift_ - timeMean_) * hcal::constants::invertnsPerBx;
161  //Interpolate the fit (Quickly)
162  funcShape(pulse_shape_, pars[i * 2], pars[i * 2 + 1], psFit_slew[time], true);
163  // add an uncertainty from the pulse (currently noise * pulse height =>Ecal uses full cov)
165  for (j = 0; j < nSamplesToFit_; ++j) {
168  }
169  }
170  } else {
171  int time = (pars[0] + timeShift_ - timeMean_) * hcal::constants::invertnsPerBx;
172  //Interpolate the fit (Quickly)
173  funcShape(pulse_shape_, pars[0], pars[1], psFit_slew[time], true);
174  for (j = 0; j < nSamplesToFit_; ++j)
175  pulse_shape_sum_[j] = pulse_shape_[j] + pedestal;
176 
177  for (i = 1; i < parBy2; ++i) {
178  time = (pars[i * 2] + timeShift_ - timeMean_) * hcal::constants::invertnsPerBx;
179  //Interpolate the fit (Quickly)
180  funcShape(pulse_shape_, pars[i * 2], pars[i * 2 + 1], psFit_slew[time], true);
181  // add an uncertainty from the pulse (currently noise * pulse height =>Ecal uses full cov)
182  for (j = 0; j < nSamplesToFit_; ++j)
184  }
185  }
186 
187  for (i = 0; i < nSamplesToFit_; ++i) {
188  const double d = psFit_y[i] - pulse_shape_sum_[i];
189  chisq += d * d / psFit_erry2[i];
190  }
191 
192  if (pedestalConstraint_) {
193  //Add the pedestal Constraint to chi2
194  chisq += invertpedSig2_ * (pedestal - pedMean_) * (pedestal - pedMean_);
195  }
196  //Add the time Constraint to chi2
197  if (timeConstraint_) {
198  for (j = 0; j < parBy2; ++j) {
199  int time = (pars[j * 2] + timeShift_ - timeMean_) * (double)hcal::constants::invertnsPerBx;
200  double time1 = -100. + time * hcal::constants::nsPerBX;
201  chisq += inverttimeSig2_ * (pars[j * 2] - timeMean_ - time1) * (pars[j * 2] - timeMean_ - time1);
202  }
203  }
204  return chisq;
205  }
206 } // namespace FitterFuncs
double psFit_x[hcal::constants::maxSamples]
void EvalPulse(const float *pars)
double psFit_erry2[hcal::constants::maxSamples]
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
std::vector< float > accVarLenIdxMinusOneVec_
double psFit_slew[hcal::constants::maxSamples]
tuple d
Definition: ztail.py:151
double psFit_erry[hcal::constants::maxSamples]
constexpr int maxPSshapeBin
Definition: HcalConstants.h:7
PulseShapeFunctor(const HcalPulseShapes::Shape &pulse, bool iPedestalConstraint, bool iTimeConstraint, bool iAddPulseJitter, double iPulseJitter, double iTimeMean, double iPedMean, unsigned int nSamplesToFit)
std::vector< float > diff25nsItvlVec_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
constexpr int nsPerBX
Definition: HcalConstants.h:8
std::array< double, hcal::constants::maxSamples > pulse_shape_
std::vector< float > acc25nsVec_
std::array< float, hcal::constants::maxPSshapeBin > pulse_hist
double psFit_y[hcal::constants::maxSamples]
double pulse(double x, double y, double z, double t)
void funcShape(std::array< double, hcal::constants::maxSamples > &ntmpbin, const double pulseTime, const double pulseHeight, const double slew, bool scalePulse)
constexpr int maxSamples
Definition: HcalConstants.h:6
constexpr float iniTimeShift
Definition: HcalConstants.h:9
std::array< double, hcal::constants::maxSamples > pulse_shape_sum_
std::vector< float > accVarLenIdxZEROVec_
std::vector< float > diffVarItvlIdxZEROVec_
std::vector< float > diffVarItvlIdxMinusOneVec_
double EvalPulseM2(const double *pars, const unsigned nPar)
constexpr float invertnsPerBx
Definition: HcalConstants.h:10