CMS 3D CMS Logo

PulseChiSqSNNLS.cc
Go to the documentation of this file.
2 #include <cmath>
4 #include <iostream>
5 
6 void eigen_solve_submatrix(PulseMatrix &mat, PulseVector &invec, PulseVector &outvec, unsigned NP) {
7  using namespace Eigen;
8  switch (NP) { // pulse matrix is always square.
9  case 10: {
10  Matrix<double, 10, 10> temp = mat.topLeftCorner<10, 10>();
11  outvec.head<10>() = temp.ldlt().solve(invec.head<10>());
12  } break;
13  case 9: {
14  Matrix<double, 9, 9> temp = mat.topLeftCorner<9, 9>();
15  outvec.head<9>() = temp.ldlt().solve(invec.head<9>());
16  } break;
17  case 8: {
18  Matrix<double, 8, 8> temp = mat.topLeftCorner<8, 8>();
19  outvec.head<8>() = temp.ldlt().solve(invec.head<8>());
20  } break;
21  case 7: {
22  Matrix<double, 7, 7> temp = mat.topLeftCorner<7, 7>();
23  outvec.head<7>() = temp.ldlt().solve(invec.head<7>());
24  } break;
25  case 6: {
26  Matrix<double, 6, 6> temp = mat.topLeftCorner<6, 6>();
27  outvec.head<6>() = temp.ldlt().solve(invec.head<6>());
28  } break;
29  case 5: {
30  Matrix<double, 5, 5> temp = mat.topLeftCorner<5, 5>();
31  outvec.head<5>() = temp.ldlt().solve(invec.head<5>());
32  } break;
33  case 4: {
34  Matrix<double, 4, 4> temp = mat.topLeftCorner<4, 4>();
35  outvec.head<4>() = temp.ldlt().solve(invec.head<4>());
36  } break;
37  case 3: {
38  Matrix<double, 3, 3> temp = mat.topLeftCorner<3, 3>();
39  outvec.head<3>() = temp.ldlt().solve(invec.head<3>());
40  } break;
41  case 2: {
42  Matrix<double, 2, 2> temp = mat.topLeftCorner<2, 2>();
43  outvec.head<2>() = temp.ldlt().solve(invec.head<2>());
44  } break;
45  case 1: {
46  Matrix<double, 1, 1> temp = mat.topLeftCorner<1, 1>();
47  outvec.head<1>() = temp.ldlt().solve(invec.head<1>());
48  } break;
49  default:
50  throw cms::Exception("MultFitWeirdState")
51  << "Weird number of pulses encountered in multifit, module is configured incorrectly!";
52  }
53 }
54 
55 PulseChiSqSNNLS::PulseChiSqSNNLS() : _chisq(0.), _computeErrors(true), _maxiters(50), _maxiterwarnings(true) {}
56 
58 
60  const SampleMatrix &samplecov,
61  const BXVector &bxs,
62  const FullSampleVector &fullpulse,
63  const FullSampleMatrix &fullpulsecov,
64  const SampleGainVector &gains,
65  const SampleGainVector &badSamples) {
66  int npulse = bxs.rows();
67 
68  _sampvec = samples;
69  _bxs = bxs;
70  _pulsemat.resize(Eigen::NoChange, npulse);
71 
72  //construct dynamic pedestals if applicable
73  int ngains = gains.maxCoeff() + 1;
74  int nPedestals = 0;
75  for (int gainidx = 0; gainidx < ngains; ++gainidx) {
76  SampleGainVector mask = gainidx * SampleGainVector::Ones();
77  SampleVector pedestal = (gains.array() == mask.array()).cast<SampleVector::value_type>();
78  if (pedestal.maxCoeff() > 0.) {
79  ++nPedestals;
80  _bxs.resize(npulse + nPedestals);
81  _bxs[npulse + nPedestals - 1] = 100 + gainidx; //bx values >=100 indicate dynamic pedestals
82  _pulsemat.resize(Eigen::NoChange, npulse + nPedestals);
83  _pulsemat.col(npulse + nPedestals - 1) = pedestal;
84  }
85  }
86 
87  //construct negative step functions for saturated or potentially slew-rate-limited samples
88  for (int isample = 0; isample < SampleVector::RowsAtCompileTime; ++isample) {
89  if (badSamples.coeff(isample) > 0) {
90  SampleVector step = SampleVector::Zero();
91  //step correction has negative sign for saturated or slew-limited samples which have been forced to zero
92  step[isample] = -1.;
93 
94  ++nPedestals;
95  _bxs.resize(npulse + nPedestals);
96  _bxs[npulse + nPedestals - 1] =
97  -100 - isample; //bx values <=-100 indicate step corrections for saturated or slew-limited samples
98  _pulsemat.resize(Eigen::NoChange, npulse + nPedestals);
99  _pulsemat.col(npulse + nPedestals - 1) = step;
100  }
101  }
102 
103  _npulsetot = npulse + nPedestals;
104 
105  _ampvec = PulseVector::Zero(_npulsetot);
106  _errvec = PulseVector::Zero(_npulsetot);
107  _nP = 0;
108  _chisq = 0.;
109 
110  if (_bxs.rows() == 1 && std::abs(_bxs.coeff(0)) < 100) {
111  _ampvec.coeffRef(0) = _sampvec.coeff(_bxs.coeff(0) + 5);
112  }
113 
114  aTamat.resize(_npulsetot, _npulsetot);
115 
116  //initialize pulse template matrix
117  for (int ipulse = 0; ipulse < npulse; ++ipulse) {
118  int bx = _bxs.coeff(ipulse);
119  int offset = 7 - 3 - bx;
120  _pulsemat.col(ipulse) = fullpulse.segment<SampleVector::RowsAtCompileTime>(offset);
121  }
122 
123  //unconstrain pedestals already for first iteration since they should always be non-zero
124  if (nPedestals > 0) {
125  for (int i = 0; i < _bxs.rows(); ++i) {
126  int bx = _bxs.coeff(i);
127  if (bx >= 100) {
129  }
130  }
131  }
132 
133  //do the actual fit
134  bool status = Minimize(samplecov, fullpulsecov);
136  _bxsmin = _bxs;
137 
138  if (!status)
139  return status;
140 
141  if (!_computeErrors)
142  return status;
143 
144  //compute MINOS-like uncertainties for in-time amplitude
145  bool foundintime = false;
146  unsigned int ipulseintime = 0;
147  for (unsigned int ipulse = 0; ipulse < _npulsetot; ++ipulse) {
148  if (_bxs.coeff(ipulse) == 0) {
149  ipulseintime = ipulse;
150  foundintime = true;
151  break;
152  }
153  }
154  if (!foundintime)
155  return status;
156 
157  const unsigned int ipulseintimemin = ipulseintime;
158 
159  double approxerr = ComputeApproxUncertainty(ipulseintime);
160  double chisq0 = _chisq;
161  double x0 = _ampvecmin[ipulseintime];
162 
163  //move in time pulse first to active set if necessary
164  if (ipulseintime < _nP) {
165  _pulsemat.col(_nP - 1).swap(_pulsemat.col(ipulseintime));
166  std::swap(_ampvec.coeffRef(_nP - 1), _ampvec.coeffRef(ipulseintime));
167  std::swap(_bxs.coeffRef(_nP - 1), _bxs.coeffRef(ipulseintime));
168  ipulseintime = _nP - 1;
169  --_nP;
170  }
171 
172  SampleVector pulseintime = _pulsemat.col(ipulseintime);
173  _pulsemat.col(ipulseintime).setZero();
174 
175  //two point interpolation for upper uncertainty when amplitude is away from boundary
176  double xplus100 = x0 + approxerr;
177  _ampvec.coeffRef(ipulseintime) = xplus100;
178  _sampvec = samples - _ampvec.coeff(ipulseintime) * pulseintime;
179  status &= Minimize(samplecov, fullpulsecov);
180  if (!status)
181  return status;
182  double chisqplus100 = ComputeChiSq();
183 
184  double sigmaplus = std::abs(xplus100 - x0) / sqrt(chisqplus100 - chisq0);
185 
186  //if amplitude is sufficiently far from the boundary, compute also the lower uncertainty and average them
187  if ((x0 / sigmaplus) > 0.5) {
188  for (unsigned int ipulse = 0; ipulse < _npulsetot; ++ipulse) {
189  if (_bxs.coeff(ipulse) == 0) {
190  ipulseintime = ipulse;
191  break;
192  }
193  }
194  double xminus100 = std::max(0., x0 - approxerr);
195  _ampvec.coeffRef(ipulseintime) = xminus100;
196  _sampvec = samples - _ampvec.coeff(ipulseintime) * pulseintime;
197  status &= Minimize(samplecov, fullpulsecov);
198  if (!status)
199  return status;
200  double chisqminus100 = ComputeChiSq();
201 
202  double sigmaminus = std::abs(xminus100 - x0) / sqrt(chisqminus100 - chisq0);
203  _errvec[ipulseintimemin] = 0.5 * (sigmaplus + sigmaminus);
204 
205  } else {
206  _errvec[ipulseintimemin] = sigmaplus;
207  }
208 
209  _chisq = chisq0;
210 
211  return status;
212 }
213 
214 bool PulseChiSqSNNLS::Minimize(const SampleMatrix &samplecov, const FullSampleMatrix &fullpulsecov) {
215  const unsigned int npulse = _bxs.rows();
216 
217  int iter = 0;
218  bool status = false;
219  while (true) {
220  if (iter >= _maxiters) {
221  if (_maxiterwarnings) {
222  LogDebug("PulseChiSqSNNLS::Minimize") << "Max Iterations reached at iter " << iter;
223  }
224  break;
225  }
226 
227  status = updateCov(samplecov, fullpulsecov);
228  if (!status)
229  break;
230  if (npulse > 1) {
231  status = NNLS();
232  } else {
233  //special case for one pulse fit (performance optimized)
235  }
236  if (!status)
237  break;
238 
239  double chisqnow = ComputeChiSq();
240  double deltachisq = chisqnow - _chisq;
241 
242  _chisq = chisqnow;
243  if (std::abs(deltachisq) < 1e-3) {
244  break;
245  }
246  ++iter;
247  }
248 
249  return status;
250 }
251 
252 bool PulseChiSqSNNLS::updateCov(const SampleMatrix &samplecov, const FullSampleMatrix &fullpulsecov) {
253  const unsigned int nsample = SampleVector::RowsAtCompileTime;
254  const unsigned int npulse = _bxs.rows();
255 
256  _invcov = samplecov; //
257 
258  for (unsigned int ipulse = 0; ipulse < npulse; ++ipulse) {
259  if (_ampvec.coeff(ipulse) == 0.)
260  continue;
261  int bx = _bxs.coeff(ipulse);
262  if (std::abs(bx) >= 100)
263  continue; //no contribution to covariance from pedestal or saturation/slew step correction
264 
265  int firstsamplet = std::max(0, bx + 3);
266  int offset = 7 - 3 - bx;
267 
268  const double ampveccoef = _ampvec.coeff(ipulse);
269  const double ampsq = ampveccoef * ampveccoef;
270 
271  const unsigned int nsamplepulse = nsample - firstsamplet;
272  _invcov.block(firstsamplet, firstsamplet, nsamplepulse, nsamplepulse) +=
273  ampsq * fullpulsecov.block(firstsamplet + offset, firstsamplet + offset, nsamplepulse, nsamplepulse);
274  }
275 
276  _covdecomp.compute(_invcov);
277 
278  bool status = true;
279  return status;
280 }
281 
283  // SampleVector resvec = _pulsemat*_ampvec - _sampvec;
284  // return resvec.transpose()*_covdecomp.solve(resvec);
285 
286  return _covdecomp.matrixL().solve(_pulsemat * _ampvec - _sampvec).squaredNorm();
287 }
288 
289 double PulseChiSqSNNLS::ComputeApproxUncertainty(unsigned int ipulse) {
290  //compute approximate uncertainties
291  //(using 1/second derivative since full Hessian is not meaningful in
292  //presence of positive amplitude boundaries.)
293 
294  return 1. / _covdecomp.matrixL().solve(_pulsemat.col(ipulse)).norm();
295 }
296 
298  //Fast NNLS (fnnls) algorithm as per http://citeseerx.ist.psu.edu/viewdoc/download?doi=10.1.1.157.9203&rep=rep1&type=pdf
299 
300  const unsigned int npulse = _bxs.rows();
301  constexpr unsigned int nsamples = SampleVector::RowsAtCompileTime;
302 
303  invcovp = _covdecomp.matrixL().solve(_pulsemat);
304  aTamat.noalias() = invcovp.transpose().lazyProduct(invcovp);
305  aTbvec.noalias() = invcovp.transpose().lazyProduct(_covdecomp.matrixL().solve(_sampvec));
306 
307  int iter = 0;
308  Index idxwmax = 0;
309  double wmax = 0.0;
310  double threshold = 1e-11;
311  while (true) {
312  //can only perform this step if solution is guaranteed viable
313  if (iter > 0 || _nP == 0) {
314  if (_nP == std::min(npulse, nsamples))
315  break;
316 
317  const unsigned int nActive = npulse - _nP;
318 
320  Index idxwmaxprev = idxwmax;
321  double wmaxprev = wmax;
322  wmax = updatework.tail(nActive).maxCoeff(&idxwmax);
323 
324  //convergence
325  if (wmax < threshold || (idxwmax == idxwmaxprev && wmax == wmaxprev))
326  break;
327 
328  //worst case protection
329  if (iter >= 500) {
330  LogDebug("PulseChiSqSNNLS::NNLS()") << "Max Iterations reached at iter " << iter;
331  break;
332  }
333 
334  //unconstrain parameter
335  Index idxp = _nP + idxwmax;
337  }
338 
339  while (true) {
340  //printf("iter in, idxsP = %i\n",int(_idxsP.size()));
341 
342  if (_nP == 0)
343  break;
344 
346 
347  //solve for unconstrained parameters
348  //need to have specialized function to call optimized versions
349  // of matrix solver... this is truly amazing...
351 
352  //check solution
353  bool positive = true;
354  for (unsigned int i = 0; i < _nP; ++i)
355  positive &= (ampvecpermtest(i) > 0);
356  if (positive) {
357  _ampvec.head(_nP) = ampvecpermtest.head(_nP);
358  break;
359  }
360 
361  //update parameter vector
362  Index minratioidx = 0;
363 
364  // no realizable optimization here (because it autovectorizes!)
365  double minratio = std::numeric_limits<double>::max();
366  for (unsigned int ipulse = 0; ipulse < _nP; ++ipulse) {
367  if (ampvecpermtest.coeff(ipulse) <= 0.) {
368  const double c_ampvec = _ampvec.coeff(ipulse);
369  const double ratio = c_ampvec / (c_ampvec - ampvecpermtest.coeff(ipulse));
370  if (ratio < minratio) {
371  minratio = ratio;
372  minratioidx = ipulse;
373  }
374  }
375  }
376 
377  _ampvec.head(_nP) += minratio * (ampvecpermtest.head(_nP) - _ampvec.head(_nP));
378 
379  //avoid numerical problems with later ==0. check
380  _ampvec.coeffRef(minratioidx) = 0.;
381 
382  //printf("removing index %i, orig idx %i\n",int(minratioidx),int(_bxs.coeff(minratioidx)));
383  NNLSConstrainParameter(minratioidx);
384  }
385  ++iter;
386 
387  //adaptive convergence threshold to avoid infinite loops but still
388  //ensure best value is used
389  if (iter % 16 == 0) {
390  threshold *= 2;
391  }
392  }
393 
394  return true;
395 }
396 
398  aTamat.col(_nP).swap(aTamat.col(idxp));
399  aTamat.row(_nP).swap(aTamat.row(idxp));
400  _pulsemat.col(_nP).swap(_pulsemat.col(idxp));
401  std::swap(aTbvec.coeffRef(_nP), aTbvec.coeffRef(idxp));
402  std::swap(_ampvec.coeffRef(_nP), _ampvec.coeffRef(idxp));
403  std::swap(_bxs.coeffRef(_nP), _bxs.coeffRef(idxp));
404  ++_nP;
405 }
406 
408  aTamat.col(_nP - 1).swap(aTamat.col(minratioidx));
409  aTamat.row(_nP - 1).swap(aTamat.row(minratioidx));
410  _pulsemat.col(_nP - 1).swap(_pulsemat.col(minratioidx));
411  std::swap(aTbvec.coeffRef(_nP - 1), aTbvec.coeffRef(minratioidx));
412  std::swap(_ampvec.coeffRef(_nP - 1), _ampvec.coeffRef(minratioidx));
413  std::swap(_bxs.coeffRef(_nP - 1), _bxs.coeffRef(minratioidx));
414  --_nP;
415 }
416 
418  //Fast NNLS (fnnls) algorithm as per http://citeseerx.ist.psu.edu/viewdoc/download?doi=10.1.1.157.9203&rep=rep1&type=pdf
419 
420  // const unsigned int npulse = 1;
421 
422  invcovp = _covdecomp.matrixL().solve(_pulsemat);
423  // aTamat = invcovp.transpose()*invcovp;
424  // aTbvec = invcovp.transpose()*_covdecomp.matrixL().solve(_sampvec);
425 
426  SingleMatrix aTamatval = invcovp.transpose() * invcovp;
427  SingleVector aTbvecval = invcovp.transpose() * _covdecomp.matrixL().solve(_sampvec);
428  _ampvec.coeffRef(0) = std::max(0., aTbvecval.coeff(0) / aTamatval.coeff(0));
429 
430  return true;
431 }
PulseMatrix aTamat
void NNLSUnconstrainParameter(Index idxp)
unsigned int _npulsetot
Eigen::Matrix< double, 1, 1 > SingleVector
unsigned int _nP
Eigen::Matrix< double, FullSampleVectorSize, FullSampleVectorSize > FullSampleMatrix
PulseVector _ampvecmin
Eigen::Matrix< double, FullSampleVectorSize, 1 > FullSampleVector
PulseVector aTbvec
PulseVector ampvecpermtest
bool DoFit(const SampleVector &samples, const SampleMatrix &samplecov, const BXVector &bxs, const FullSampleVector &fullpulse, const FullSampleMatrix &fullpulsecov, const SampleGainVector &gains=-1 *SampleGainVector::Ones(), const SampleGainVector &badSamples=SampleGainVector::Zero())
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, PulseVectorSize, 1 > PulseVector
constexpr uint32_t mask
Definition: gpuClustering.h:24
SamplePulseMatrix _pulsemat
double ComputeApproxUncertainty(unsigned int ipulse)
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
T sqrt(T t)
Definition: SSEVec.h:19
SampleDecompLLT _covdecomp
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool updateCov(const SampleMatrix &samplecov, const FullSampleMatrix &fullpulsecov)
bool Minimize(const SampleMatrix &samplecov, const FullSampleMatrix &fullpulsecov)
constexpr float gains[NGAINS]
Definition: EcalConstants.h:11
PulseVector updatework
SampleVector _sampvec
void resize(int bx, unsigned size)
Eigen::Matrix< double, SampleVectorSize, 1 > SampleVector
Eigen::Matrix< double, SampleVectorSize, SampleVectorSize > SampleMatrix
SampleMatrix _invcov
Eigen::Matrix< char, SampleVectorSize, 1 > SampleGainVector
void eigen_solve_submatrix(PulseMatrix &mat, PulseVector &invec, PulseVector &outvec, unsigned NP)
PulseVector _errvec
Eigen::Matrix< double, 1, 1 > SingleMatrix
BXVector::Index Index
step
Definition: StallMonitor.cc:98
void NNLSConstrainParameter(Index minratioidx)
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, PulseVectorSize, PulseVectorSize > PulseMatrix
SamplePulseMatrix invcovp
PulseVector _ampvec
#define LogDebug(id)
__host__ __device__ V V wmax