CMS 3D CMS Logo

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