CMS 3D CMS Logo

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, bool iTimeConstraint,bool iAddPulseJitter,
12  double iPulseJitter,double iTimeMean,double iPedMean,
13  unsigned nSamplesToFit) :
14  cntNANinfit(0),
15  acc25nsVec(HcalConst::maxPSshapeBin), diff25nsItvlVec(HcalConst::maxPSshapeBin),
16  accVarLenIdxZEROVec(HcalConst::nsPerBX), diffVarItvlIdxZEROVec(HcalConst::nsPerBX),
17  accVarLenIdxMinusOneVec(HcalConst::nsPerBX), diffVarItvlIdxMinusOneVec(HcalConst::nsPerBX) {
18  //The raw pulse
19  for(int i=0;i<HcalConst::maxPSshapeBin;++i) pulse_hist[i] = pulse(i);
20  // Accumulate 25ns for each starting point of 0, 1, 2, 3...
21  for(int i=0; i<HcalConst::maxPSshapeBin; ++i){
22  for(int j=i; j<i+HcalConst::nsPerBX; ++j){ //sum over HcalConst::nsPerBXns from point i
23  acc25nsVec[i] += ( j < HcalConst::maxPSshapeBin? pulse_hist[j] : pulse_hist[HcalConst::maxPSshapeBin-1]);
24  }
25  diff25nsItvlVec[i] = ( i+HcalConst::nsPerBX < HcalConst::maxPSshapeBin? pulse_hist[i+HcalConst::nsPerBX] - pulse_hist[i] : pulse_hist[HcalConst::maxPSshapeBin-1] - pulse_hist[i]);
26  }
27  // Accumulate different ns for starting point of index either 0 or -1
28  for(int i=0; i<HcalConst::nsPerBX; ++i){
29  if( i==0 ){
32  } else{
34  accVarLenIdxMinusOneVec[i] = accVarLenIdxMinusOneVec[i-1] + pulse_hist[i-1];
35  }
37  diffVarItvlIdxMinusOneVec[i] = pulse_hist[i] - pulse_hist[0];
38  }
39  for(int i = 0; i < HcalConst::maxSamples; i++) {
40  psFit_x[i] = 0;
41  psFit_y[i] = 0;
42  psFit_erry[i] = 1.;
43  psFit_erry2[i] = 1.;
44  psFit_slew [i] = 0.;
45  }
46  //Constraints
47  pedestalConstraint_ = iPedestalConstraint;
48  timeConstraint_ = iTimeConstraint;
49  addPulseJitter_ = iAddPulseJitter;
50  pulseJitter_ = iPulseJitter*iPulseJitter;
51 
52  // for M2
53  timeMean_ = iTimeMean;
54  pedMean_ = iPedMean;
55  timeShift_ = 100.;
56  timeShift_ += 12.5;//we are trying to get BX
57 
58  nSamplesToFit_ = nSamplesToFit;
59 
60  }
61 
62  void PulseShapeFunctor::funcShape(std::array<double,HcalConst::maxSamples> & ntmpbin, const double pulseTime, const double pulseHeight,const double slew) {
63  // pulse shape components over a range of time 0 ns to 255 ns in 1 ns steps
64  constexpr int ns_per_bx = HcalConst::nsPerBX;
66  constexpr int num_bx = num_ns/ns_per_bx;
67  //Get the starting time
68  int i_start = ( -HcalConst::iniTimeShift - pulseTime - slew >0 ? 0 : (int)std::abs(-HcalConst::iniTimeShift-pulseTime-slew) + 1);
69  double offset_start = i_start - HcalConst::iniTimeShift - pulseTime - slew; //-199-2*pars[0]-2.*slew (for pars[0] > 98.5) or just -98.5-pars[0]-slew;
70  // zeroing output binned pulse shape
71  ntmpbin.fill(0.0f);
72 
73  if( edm::isNotFinite(offset_start) ){ //Check for nan
74  ++ cntNANinfit;
75  }else{
76  if( offset_start == 1.0 ){ offset_start = 0.; i_start-=1; } //Deal with boundary
77 
78  const int bin_start = (int) offset_start; //bin off to integer
79  const int bin_0_start = ( offset_start < bin_start + 0.5 ? bin_start -1 : bin_start ); //Round it
80  const int iTS_start = i_start/ns_per_bx; //Time Slice for time shift
81  const int distTo25ns_start = HcalConst::nsPerBX - 1 - i_start%ns_per_bx; //Delta ns
82  const double factor = offset_start - bin_0_start - 0.5; //Small correction?
83 
84  //Build the new pulse
85  ntmpbin[iTS_start] = (bin_0_start == -1 ? // Initial bin (I'm assuming this is ok)
86  accVarLenIdxMinusOneVec[distTo25ns_start] + factor * diffVarItvlIdxMinusOneVec[distTo25ns_start]
87  : accVarLenIdxZEROVec [distTo25ns_start] + factor * diffVarItvlIdxZEROVec [distTo25ns_start]);
88  //Fill the rest of the bins
89  for(int iTS = iTS_start+1; iTS < num_bx; ++iTS){
90  int bin_idx = distTo25ns_start + 1 + (iTS-iTS_start-1)*ns_per_bx + bin_0_start;
91  ntmpbin[iTS] = acc25nsVec[bin_idx] + factor * diff25nsItvlVec[bin_idx];
92  }
93  //Scale the pulse
94  for(int i=iTS_start; i < num_bx; ++i) {
95  ntmpbin[i] *= pulseHeight;
96  }
97 
98  }
99 
100  return;
101  }
102 
104  }
105 
106  double PulseShapeFunctor::EvalPulse(const double *pars, const unsigned nPars) {
107 
108  unsigned i =0, j=0;
109 
110  const double pedestal=pars[nPars-1];
111 
112  //Stop crashes
113  for(i =0; i < nPars; ++i ) if( edm::isNotFinite(pars[i]) ){ ++ cntNANinfit; return 1e10; }
114 
115  //calculate chisquare
116  double chisq = 0;
117  const unsigned parBy2=(nPars-1)/2;
118  // std::array<float,HcalConst::maxSamples> pulse_shape_;
119 
120  if(addPulseJitter_) {
122  //Interpolate the fit (Quickly)
123  funcShape(pulse_shape_, pars[0],pars[1],psFit_slew[time]);
124  for (j=0; j<nSamplesToFit_; ++j) {
127  }
128 
129  for (i=1; i<parBy2; ++i) {
130  time = (pars[i*2]+timeShift_-timeMean_)*HcalConst::invertnsPerBx;
131  //Interpolate the fit (Quickly)
132  funcShape(pulse_shape_, pars[i*2],pars[i*2+1],psFit_slew[time]);
133  // add an uncertainty from the pulse (currently noise * pulse height =>Ecal uses full cov)
135  for (j=0; j<nSamplesToFit_; ++j) {
138  }
139  }
140  }
141  else{
143  //Interpolate the fit (Quickly)
144  funcShape(pulse_shape_, pars[0],pars[1],psFit_slew[time]);
145  for(j=0; j<nSamplesToFit_; ++j)
146  pulse_shape_sum_[j] = pulse_shape_[j] + pedestal;
147 
148  for (i=1; i<parBy2; ++i) {
149  time = (pars[i*2]+timeShift_-timeMean_)*HcalConst::invertnsPerBx;
150  //Interpolate the fit (Quickly)
151  funcShape(pulse_shape_, pars[i*2],pars[i*2+1],psFit_slew[time]);
152  // add an uncertainty from the pulse (currently noise * pulse height =>Ecal uses full cov)
153  for(j=0; j<nSamplesToFit_; ++j)
155  }
156  }
157 
158  for (i=0;i<nSamplesToFit_; ++i)
159  {
160  const double d = psFit_y[i]- pulse_shape_sum_[i];
161  chisq += d*d/psFit_erry2[i];
162  }
163 
164  if(pedestalConstraint_) {
165  //Add the pedestal Constraint to chi2
166  chisq += invertpedSig2_*(pedestal - pedMean_)*(pedestal - pedMean_);
167  }
168  //Add the time Constraint to chi2
169  if(timeConstraint_) {
170  for(j=0; j< parBy2; ++j ){
171  int time = (pars[j*2]+timeShift_-timeMean_)*(double)HcalConst::invertnsPerBx;
172  double time1 = -100.+time*HcalConst::nsPerBX;
173  chisq += inverttimeSig2_*(pars[j*2] - timeMean_ - time1)*(pars[j*2] - timeMean_ - time1);
174  }
175  }
176  return chisq;
177  }
178 
179  double PulseShapeFunctor::singlePulseShapeFunc( const double *x ) {
180  return EvalPulse(x,3);
181  }
182 
183  double PulseShapeFunctor::doublePulseShapeFunc( const double *x ) {
184  return EvalPulse(x,5);
185  }
186 
187  double PulseShapeFunctor::triplePulseShapeFunc( const double *x ) {
188  return EvalPulse(x,7);
189  }
190 
191 }
double psFit_erry[HcalConst::maxSamples]
double singlePulseShapeFunc(const double *x)
double EvalPulse(const double *pars, const unsigned nPar)
std::array< double, HcalConst::maxSamples > pulse_shape_
std::vector< float > acc25nsVec
#define constexpr
float iniTimeShift
double psFit_y[HcalConst::maxSamples]
std::vector< float > diffVarItvlIdxZEROVec
bool isNotFinite(T x)
Definition: isFinite.h:10
PulseShapeFunctor(const HcalPulseShapes::Shape &pulse, bool iPedestalConstraint, bool iTimeConstraint, bool iAddPulseJitter, double iPulseJitter, double iTimeMean, double iPedMean, unsigned int nSamplesToFit)
void funcShape(std::array< double, HcalConst::maxSamples > &ntmpbin, const double pulseTime, const double pulseHeight, const double slew)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::array< float, HcalConst::maxPSshapeBin > pulse_hist
double f[11][100]
std::vector< float > diff25nsItvlVec
std::vector< float > accVarLenIdxMinusOneVec
std::vector< float > diffVarItvlIdxMinusOneVec
double invertnsPerBx
double pulse(double x, double y, double z, double t)
double triplePulseShapeFunc(const double *x)
double doublePulseShapeFunc(const double *x)
double psFit_erry2[HcalConst::maxSamples]
double psFit_x[HcalConst::maxSamples]
std::array< double, HcalConst::maxSamples > pulse_shape_sum_
double psFit_slew[HcalConst::maxSamples]
std::vector< float > accVarLenIdxZEROVec