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 {
79 
80  Eigen::initParallel();
81 
82 }
83 
85 
86 }
87 
88 bool PulseChiSqSNNLS::DoFit(const SampleVector &samples, const SampleMatrix &samplecor, double pederr, const BXVector &bxs, const FullSampleVector &fullpulse, const FullSampleMatrix &fullpulsecov) {
89 
90  //const unsigned int nsample = SampleVector::RowsAtCompileTime;
91  _npulsetot = bxs.rows();
92  //const unsigned int npulse = bxs.rows();
93 
94  _sampvec = samples;
95  _bxs = bxs;
96 
97  //_pulsemat = SamplePulseMatrix::Zero(nsample,npulse);
98  _pulsemat.resize(Eigen::NoChange,_npulsetot);
99  _ampvec = PulseVector::Zero(_npulsetot);
100  _errvec = PulseVector::Zero(_npulsetot);
101  _nP = 0;
102  _chisq = 0.;
103 
104  aTamat.resize(_npulsetot,_npulsetot);
105  wvec.resize(_npulsetot);
106 
107  //initialize pulse template matrix
108  for (unsigned int ipulse=0; ipulse<_npulsetot; ++ipulse) {
109  int bx = _bxs.coeff(ipulse);
110  //int firstsamplet = std::max(0,bx + 3);
111  int offset = 7-3-bx;
112 
113  //const unsigned int nsamplepulse = nsample-firstsamplet;
114  //_pulsemat.col(ipulse).segment(firstsamplet,nsamplepulse) = fullpulse.segment(firstsamplet+offset,nsamplepulse);
115 
116  _pulsemat.col(ipulse) = fullpulse.segment<SampleVector::RowsAtCompileTime>(offset);
117  }
118 
119  //do the actual fit
120  bool status = Minimize(samplecor,pederr,fullpulsecov);
122  _bxsmin = _bxs;
123 
124  if (!status) return status;
125 
126  if(!_computeErrors) return status;
127 
128  //compute MINOS-like uncertainties for in-time amplitude
129  bool foundintime = false;
130  unsigned int ipulseintime = 0;
131  for (unsigned int ipulse=0; ipulse<_npulsetot; ++ipulse) {
132  if (_bxs.coeff(ipulse)==0) {
133  ipulseintime = ipulse;
134  foundintime = true;
135  break;
136  }
137  }
138  if (!foundintime) return status;
139 
140  const unsigned int ipulseintimemin = ipulseintime;
141 
142  double approxerr = ComputeApproxUncertainty(ipulseintime);
143  double chisq0 = _chisq;
144  double x0 = _ampvecmin[ipulseintime];
145 
146  //move in time pulse first to active set if necessary
147  if (ipulseintime<_nP) {
148  _pulsemat.col(_nP-1).swap(_pulsemat.col(ipulseintime));
149  std::swap(_ampvec.coeffRef(_nP-1),_ampvec.coeffRef(ipulseintime));
150  std::swap(_bxs.coeffRef(_nP-1),_bxs.coeffRef(ipulseintime));
151  ipulseintime = _nP - 1;
152  --_nP;
153  }
154 
155 
156  SampleVector pulseintime = _pulsemat.col(ipulseintime);
157  _pulsemat.col(ipulseintime).setZero();
158 
159  //two point interpolation for upper uncertainty when amplitude is away from boundary
160  double xplus100 = x0 + approxerr;
161  _ampvec.coeffRef(ipulseintime) = xplus100;
162  _sampvec = samples - _ampvec.coeff(ipulseintime)*pulseintime;
163  status &= Minimize(samplecor,pederr,fullpulsecov);
164  if (!status) return status;
165  double chisqplus100 = ComputeChiSq();
166 
167  double sigmaplus = std::abs(xplus100-x0)/sqrt(chisqplus100-chisq0);
168 
169  //if amplitude is sufficiently far from the boundary, compute also the lower uncertainty and average them
170  if ( (x0/sigmaplus) > 0.5 ) {
171  for (unsigned int ipulse=0; ipulse<_npulsetot; ++ipulse) {
172  if (_bxs.coeff(ipulse)==0) {
173  ipulseintime = ipulse;
174  break;
175  }
176  }
177  double xminus100 = std::max(0.,x0-approxerr);
178  _ampvec.coeffRef(ipulseintime) = xminus100;
179  _sampvec = samples - _ampvec.coeff(ipulseintime)*pulseintime;
180  status &= Minimize(samplecor,pederr,fullpulsecov);
181  if (!status) return status;
182  double chisqminus100 = ComputeChiSq();
183 
184  double sigmaminus = std::abs(xminus100-x0)/sqrt(chisqminus100-chisq0);
185  _errvec[ipulseintimemin] = 0.5*(sigmaplus + sigmaminus);
186 
187  }
188  else {
189  _errvec[ipulseintimemin] = sigmaplus;
190  }
191 
192  _chisq = chisq0;
193 
194  return status;
195 
196 }
197 
198 bool PulseChiSqSNNLS::Minimize(const SampleMatrix &samplecor, double pederr, const FullSampleMatrix &fullpulsecov) {
199 
200 
201  const int maxiter = 50;
202  int iter = 0;
203  bool status = false;
204  while (true) {
205 
206  if (iter>=maxiter) {
207  edm::LogWarning("PulseChiSqSNNLS::Minimize") << "Max Iterations reached at iter " << iter << std::endl;
208  break;
209  }
210 
211  status = updateCov(samplecor,pederr,fullpulsecov);
212  if (!status) break;
213  status = NNLS();
214  if (!status) break;
215 
216  double chisqnow = ComputeChiSq();
217  double deltachisq = chisqnow-_chisq;
218 
219  _chisq = chisqnow;
220  if (std::abs(deltachisq)<1e-3) {
221  break;
222  }
223  ++iter;
224  }
225 
226  return status;
227 
228 }
229 
230 bool PulseChiSqSNNLS::updateCov(const SampleMatrix &samplecor, double pederr, const FullSampleMatrix &fullpulsecov) {
231 
232  const unsigned int nsample = SampleVector::RowsAtCompileTime;
233  const unsigned int npulse = _bxs.rows();
234 
235  const double pederr2 = pederr*pederr;
236  _invcov = pederr2*samplecor; //
237 
238  for (unsigned int ipulse=0; ipulse<npulse; ++ipulse) {
239  if (_ampvec.coeff(ipulse)==0.) continue;
240  int bx = _bxs.coeff(ipulse);
241  int firstsamplet = std::max(0,bx + 3);
242  int offset = 7-3-bx;
243 
244  const double ampveccoef = _ampvec.coeff(ipulse);
245  const double ampsq = ampveccoef*ampveccoef;
246 
247  const unsigned int nsamplepulse = nsample-firstsamplet;
248  _invcov.block(firstsamplet,firstsamplet,nsamplepulse,nsamplepulse) +=
249  ampsq*fullpulsecov.block(firstsamplet+offset,firstsamplet+offset,nsamplepulse,nsamplepulse);
250  }
251 
252  _covdecomp.compute(_invcov);
253 
254  bool status = true;
255  return status;
256 
257 }
258 
260 
261 // SampleVector resvec = _pulsemat*_ampvec - _sampvec;
262 // return resvec.transpose()*_covdecomp.solve(resvec);
263 
264  return _covdecomp.matrixL().solve(_pulsemat*_ampvec - _sampvec).squaredNorm();
265 
266 }
267 
268 double PulseChiSqSNNLS::ComputeApproxUncertainty(unsigned int ipulse) {
269  //compute approximate uncertainties
270  //(using 1/second derivative since full Hessian is not meaningful in
271  //presence of positive amplitude boundaries.)
272 
273  return 1./_covdecomp.matrixL().solve(_pulsemat.col(ipulse)).norm();
274 
275 }
276 
278 
279  //Fast NNLS (fnnls) algorithm as per http://citeseerx.ist.psu.edu/viewdoc/download?doi=10.1.1.157.9203&rep=rep1&type=pdf
280 
281  const unsigned int npulse = _bxs.rows();
282 
283  invcovp = _covdecomp.matrixL().solve(_pulsemat);
284  aTamat = invcovp.transpose()*invcovp; //.triangularView<Eigen::Lower>()
285  //aTamat = aTamat.selfadjointView<Eigen::Lower>();
286  aTbvec = invcovp.transpose()*_covdecomp.matrixL().solve(_sampvec);
287 
288  int iter = 0;
289  Index idxwmax;
290  double wmax = 0.0;
291  //work = PulseVector::zeros();
292  while (true) {
293  //can only perform this step if solution is guaranteed viable
294  if (iter>0 || _nP==0) {
295  if ( _nP==npulse ) break;
296 
297  const unsigned int nActive = npulse - _nP;
298 
300  wmax = updatework.tail(nActive).maxCoeff(&idxwmax);
301 
302  //convergence
303  if (wmax<1e-11) break;
304 
305  //unconstrain parameter
306  Index idxp = _nP + idxwmax;
307  //printf("adding index %i, orig index %i\n",int(idxp),int(_bxs.coeff(idxp)));
308  aTamat.col(_nP).swap(aTamat.col(idxp));
309  aTamat.row(_nP).swap(aTamat.row(idxp));
310  _pulsemat.col(_nP).swap(_pulsemat.col(idxp));
311  std::swap(aTbvec.coeffRef(_nP),aTbvec.coeffRef(idxp));
312  std::swap(_ampvec.coeffRef(_nP),_ampvec.coeffRef(idxp));
313  std::swap(_bxs.coeffRef(_nP),_bxs.coeffRef(idxp));
314 
315  // update now that we are done doing work
316  wvec.tail(nActive) = updatework.tail(nActive);
317  ++_nP;
318  }
319 
320 
321  while (true) {
322  //printf("iter in, idxsP = %i\n",int(_idxsP.size()));
323 
324  if (_nP==0) break;
325 
327 
328  //solve for unconstrained parameters
329  //need to have specialized function to call optimized versions
330  // of matrix solver... this is truly amazing...
332 
333  //check solution
334  auto ampvecpermhead = ampvecpermtest.head(_nP);
335  if ( ampvecpermhead.minCoeff()>0. ) {
336  _ampvec.head(_nP) = ampvecpermhead.head(_nP);
337  break;
338  }
339 
340  //update parameter vector
341  Index minratioidx=0;
342 
343  // no realizable optimization here (because it autovectorizes!)
344  double minratio = std::numeric_limits<double>::max();
345  for (unsigned int ipulse=0; ipulse<_nP; ++ipulse) {
346  if (ampvecpermtest.coeff(ipulse)<=0.) {
347  const double c_ampvec = _ampvec.coeff(ipulse);
348  const double ratio = c_ampvec/(c_ampvec-ampvecpermtest.coeff(ipulse));
349  if (ratio<minratio) {
350  minratio = ratio;
351  minratioidx = ipulse;
352  }
353  }
354  }
355 
356  _ampvec.head(_nP) += minratio*(ampvecpermhead - _ampvec.head(_nP));
357 
358  //avoid numerical problems with later ==0. check
359  _ampvec.coeffRef(minratioidx) = 0.;
360 
361  //printf("removing index %i, orig idx %i\n",int(minratioidx),int(_bxs.coeff(minratioidx)));
362  aTamat.col(_nP-1).swap(aTamat.col(minratioidx));
363  aTamat.row(_nP-1).swap(aTamat.row(minratioidx));
364  _pulsemat.col(_nP-1).swap(_pulsemat.col(minratioidx));
365  std::swap(aTbvec.coeffRef(_nP-1),aTbvec.coeffRef(minratioidx));
366  std::swap(_ampvec.coeffRef(_nP-1),_ampvec.coeffRef(minratioidx));
367  std::swap(_bxs.coeffRef(_nP-1),_bxs.coeffRef(minratioidx));
368  --_nP;
369  }
370  ++iter;
371  }
372 
373  return true;
374 
375 
376 }
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
const T & max(const T &a, const T &b)
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:48
SampleDecompLLT _covdecomp
PulseVector wvec
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static const MaxIter maxiter
unsigned int offset(bool)
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
tuple status
Definition: ntuplemaker.py:245
SamplePulseMatrix invcovp
PulseVector _ampvec
bool DoFit(const SampleVector &samples, const SampleMatrix &samplecor, double pederr, const BXVector &bxs, const FullSampleVector &fullpulse, const FullSampleMatrix &fullpulsecov)