CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Chisq_Constrainer.cc
Go to the documentation of this file.
1 //
2 // $Id: Chisq_Constrainer.cc,v 1.2 2011/10/13 09:49:48 snaumann Exp $
3 //
4 // File: src/Chisq_Constrainer.cc
5 // Purpose: Minimize a chisq subject to a set of constraints.
6 // Based on the SQUAW algorithm.
7 // Created: Jul, 2000, sss, based on run 1 mass analysis code.
8 //
9 // CMSSW File : src/Chisq_Constrainer.cc
10 // Original Author : Scott Stuart Snyder <snyder@bnl.gov> for D0
11 // Imported to CMSSW by Haryo Sumowidagdo <Suharyo.Sumowidagdo@cern.ch>
12 //
13 
41 #include <cmath>
42 #include <cassert>
43 #include <iostream>
44 #include <iomanip>
45 
46 using std::abs;
47 using std::cout;
48 using std::fixed;
49 using std::ios_base;
50 using std::ostream;
51 using std::resetiosflags;
52 using std::setiosflags;
53 using std::sqrt;
54 
55 
56 //*************************************************************************
57 
58 
59 namespace hitfit {
60 
61 
63 //
64 // Purpose: Constructor.
65 //
66 // Inputs:
67 // defs - The Defaults instance from which to initialize.
68 //
69  : _base_constrainer_args (defs)
70 {
71  _use_G = defs.get_bool ("use_G");
72  _printfit = defs.get_bool ("printfit");
73  _constraint_sum_eps = defs.get_float ("constraint_sum_eps");
74  _chisq_diff_eps = defs.get_float ("chisq_diff_eps");
75  _maxit = defs.get_int ("maxit");
76  _max_cut = defs.get_int ("maxcut");
77  _cutsize = defs.get_float ("cutsize");
78  _min_tot_cutsize = defs.get_float ("min_tot_cutsize");
79  _chisq_test_eps = defs.get_float ("chisq_test_eps");
80 }
81 
82 
84 //
85 // Purpose: Return the printfit parameter.
86 // See the header for documentation.
87 //
88 {
89  return _printfit;
90 }
91 
92 
94 //
95 // Purpose: Return the use_G parameter.
96 // See the header for documentation.
97 //
98 {
99  return _use_G;
100 }
101 
102 
104 //
105 // Purpose: Return the constraint_sum_eps parameter.
106 // See the header for documentation.
107 //
108 {
109  return _constraint_sum_eps;
110 }
111 
112 
114 //
115 // Purpose: Return the chisq_diff_eps parameter.
116 // See the header for documentation.
117 //
118 {
119  return _chisq_diff_eps;
120 }
121 
122 
124 //
125 // Purpose: Return the maxit parameter.
126 // See the header for documentation.
127 //
128 {
129  return _maxit;
130 }
131 
132 
134 //
135 // Purpose: Return the max_cut parameter.
136 // See the header for documentation.
137 //
138 {
139  return _max_cut;
140 }
141 
142 
144 //
145 // Purpose: Return the cutsize parameter.
146 // See the header for documentation.
147 //
148 {
149  return _cutsize;
150 }
151 
152 
154 //
155 // Purpose: Return the min_tot_cutsize parameter.
156 // See the header for documentation.
157 //
158 {
159  return _min_tot_cutsize;
160 }
161 
162 
164 //
165 // Purpose: Return the chisq_test_eps parameter.
166 // See the header for documentation.
167 //
168 {
169  return _chisq_test_eps;
170 }
171 
172 
175 //
176 // Purpose: Return the contained subobject parameters.
177 // See the header for documentation.
178 //
179 {
180  return _base_constrainer_args;
181 }
182 
183 
184 } // namespace hitfit
185 
186 
187 //*************************************************************************
188 
189 
190 namespace {
191 
192 
193 using namespace hitfit;
194 
195 
212 bool solve_linear_system (const Matrix& H,
213  const Diagonal_Matrix& Y,
214  const Matrix& By,
215  const Row_Vector& r,
217  Column_Vector& d,
218  Matrix& W,
219  Matrix& U,
220  Matrix& V)
221 //
222 // Purpose: Solve the system
223 //
224 // [ -H B.t ] [ alpha ] [ r ]
225 // [ ] * [ ] = [ ]
226 // [ B Y ] [ d ] [ 0 ]
227 //
228 // for alpha and d.
229 //
230 // Also returns the inverse matrices:
231 //
232 // [ W V.t ] [ -H B.t ]
233 // [ ] = [ ] ^ -1
234 // [ V U ] [ B Y ]
235 //
236 // Returns true if successful, false if not.
237 //
238 {
239  int nconstraints = H.num_row();
240  int nbadvars = Y.num_row();
241 
242  // Form the matrix on the LHS from H, By, and Y.
243  Matrix A (nconstraints+nbadvars, nconstraints+nbadvars);
244  A.sub (1, 1, -H);
245  if (nbadvars > 0) {
246  A.sub (nconstraints+1, nconstraints+1, Y);
247  A.sub (1, nconstraints+1, By.T());
248  A.sub (nconstraints+1, 1, By);
249  }
250 
251  // Form the RHS vector from r.
252  Column_Vector yy(nconstraints + nbadvars, 0);
253  yy.sub (1, r.T());
254 
255  // Invert the matrix.
256  // Try to handle singularities correctly.
257  Matrix Ai;
258  int ierr = 0;
259  do {
260  Ai = A.inverse (ierr);
261  if (ierr) {
262  int allzero = 0;
263  for (int i=1; i<=nconstraints; i++) {
264  allzero = 1;
265  for (int j=1; j<=nconstraints; j++) {
266  if (A(i,j) != 0) {
267  allzero = 0;
268  break;
269  }
270  }
271  if (allzero) {
272  A(i,i) = 1;
273  break;
274  }
275  }
276  if (!allzero) return false;
277  }
278  } while (ierr != 0);
279 
280  // Solve the system of equations.
281  Column_Vector xx = Ai * yy;
282 
283  // Extract the needed pieces from the inverted matrix
284  // and the solution vector.
285  W = Ai.sub (1, nconstraints, 1, nconstraints);
286  if (nbadvars > 0) {
287  U = Ai.sub (nconstraints+1, nconstraints+nbadvars,
288  nconstraints+1, nconstraints+nbadvars);
289  V = Ai.sub (nconstraints+1, nconstraints+nbadvars, 1, nconstraints);
290  d = xx.sub (nconstraints+1, nconstraints+nbadvars);
291  }
292 
293  alpha = xx.sub (1, nconstraints);
294 
295  return true;
296 }
297 
298 
299 } // unnamed namespace
300 
301 
302 namespace hitfit {
303 
304 
305 //*************************************************************************
306 
307 
309 //
310 // Purpose: Constructor.
311 //
312 // Inputs:
313 // args - The parameter settings for this instance.
314 //
315  : Base_Constrainer (args.base_constrainer_args()),
316  _args (args)
317 {
318 }
319 
320 
321 double Chisq_Constrainer::fit (Constraint_Calculator& constraint_calculator,
322  const Column_Vector& xm,
323  Column_Vector& x,
324  const Column_Vector& ym,
325  Column_Vector& y,
326  const Matrix& G_i,
327  const Diagonal_Matrix& Y,
328  Column_Vector& pullx,
329  Column_Vector& pully,
330  Matrix& Q,
331  Matrix& R,
332  Matrix& S)
333 //
334 // Purpose: Do a constrained fit.
335 //
336 // Call the number of well-measured variables Nw, the number of
337 // poorly-measured variables Np, and the number of constraints Nc.
338 //
339 // Inputs:
340 // constraint_calculator - The object that will be used to evaluate
341 // the constraints.
342 // xm(Nw) - The measured values of the well-measured variables.
343 // ym(Np) - The measured values of the poorly-measured variables.
344 // x(Nw) - The starting values for the well-measured variables.
345 // y(Np) - The starting values for the poorly-measured variables.
346 // G_i(Nw,Nw) - The error matrix for the well-measured variables.
347 // Y(Np,Np) - The inverse error matrix for the poorly-measured variables.
348 //
349 // Outputs:
350 // x(Nw) - The fit values of the well-measured variables.
351 // y(Np) - The fit values of the poorly-measured variables.
352 // pullx(Nw) - The pull quantities for the well-measured variables.
353 // pully(Nw) - The pull quantities for the poorly-measured variables.
354 // Q(Nw,Nw) - The final error matrix for the well-measured variables.
355 // R(Np,Np) - The final error matrix for the poorly-measured variables.
356 // S(Nw,Np) - The final cross error matrix for the two sets of variables.
357 //
358 // Returns:
359 // The minimum chisq satisfying the constraints.
360 // Returns a value < 0 if the fit failed to converge.
361 //
362 {
363  // Check that the various matrices we've been passed have consistent
364  // dimensionalities.
365  int nvars = x.num_row();
366  assert (nvars == G_i.num_col());
367  assert (nvars == xm.num_row());
368 
369  int nbadvars = y.num_row();
370  assert (nbadvars == Y.num_col());
371  assert (nbadvars == ym.num_row());
372 
373  // If we're going to check the chisq calculation by explicitly using G,
374  // calculate it now from its inverse G_i.
375  Matrix G (nvars, nvars);
376  if (_args.use_G()) {
377  int ierr = 0;
378  G = G_i.inverse (ierr);
379  assert (!ierr);
380  }
381 
382  int nconstraints = constraint_calculator.nconstraints ();
383 
384  // Results of the constraint evaluation function.
385  Row_Vector F (nconstraints); // Constraint vector.
386  Matrix Bx (nvars, nconstraints); // Gradients wrt x
387  Matrix By (nbadvars, nconstraints); // Gradients wrt y
388 
389  // (2) Evaluate the constraints at the starting point.
390  // If the starting point is rejected as invalid,
391  // give up and return an error.
392  if (! call_constraint_fcn (constraint_calculator, x, y, F, Bx, By)) {
393  // cout << "Bad initial values!";
394  // return -1000;
395  return -999.0;
396  }
397 
398  // (3) Initialize variables for the fitting loop.
399  double constraint_sum_last = -1000;
400  double chisq_last = -1000;
401  bool near_convergence = false;
402  double last_step_cutsize = 1;
403 
404  unsigned nit = 0;
405 
406  // Initialize the displacement vectors c and d.
407  Column_Vector c = x - xm;
408  Column_Vector d = y - ym;
409 
410  Matrix E (nvars, nconstraints);
411  Matrix W (nconstraints, nconstraints);
412  Matrix U (nbadvars, nbadvars);
413  Matrix V (nbadvars, nconstraints);
414 
415  // (4) Fitting loop:
416  do {
417  // (5) Calculate E, H, and r.
418  E = G_i * Bx;
419  Matrix H = E.T() * Bx;
420  Row_Vector r = c.T() * Bx + d.T() * By - F;
421 
422  // (6) Solve the linearized system for the new values
423  // of the Lagrange multipliers
424  // $\alpha$ and the new value for the displacements d.
425  Column_Vector alpha (nvars);
426  Column_Vector d1 (nbadvars);
427  if (!solve_linear_system (H, Y, By, r,
428  alpha, d1, W, U, V)) {
430  // return -1000;
431  return -998.0;
432  }
433 
434  // (7) Compute the new values for the displacements c and the chisq.
435  Column_Vector c1 = -E * alpha;
436  double chisq = - scalar (r * alpha);
437 
438  double psi_cut = 0;
439 
440  // (8) Find where this step is going to be taking us.
441  x = c1 + xm;
442  y = d1 + ym;
443 
444  // (9) Set up for cutting this step, should we have to.
445  Matrix save_By = By;
446  Row_Vector save_negF = - F;
447  double this_step_cutsize = 1;
448  double constraint_sum = -1;
449  unsigned ncut = 0;
450 
451  // (10) Evaluate the constraints at the new point.
452  // If the point is rejected, we have to try to cut the step.
453  // We accept the step if:
454  // The constraint sum is below the convergence threshold
455  // constraint_sum_eps, or
456  // This is the first iteration, or
457  // The constraint sum has decreased since the last iteration.
458  // Otherwise, the constraints have gotten worse, and we
459  // try to cut the step.
460  while (! call_constraint_fcn (constraint_calculator, x, y, F, Bx, By) ||
461  ((constraint_sum = norm_infinity (F))
463  nit > 0 &&
464  constraint_sum > constraint_sum_last))
465  {
466 
467  // Doing step cutting...
468  if (nit > 0 && _args.printfit() && ncut == 0) {
469  cout << "(" << chisq << " " << chisq_last << ") ";
470  }
471 
472  // (10a) If this is the first time we've tried to cut this step,
473  // test to see if the chisq is stationary. If it hasn't changed
474  // since the last iteration, try a directed step.
475  if (ncut == 0 &&
476  abs (chisq - chisq_last) < _args.chisq_diff_eps()) {
477 
478  // Trying a directed step now.
479  // Try to make the smallest step which satisfies the
480  // (linearized) constraints.
481  if (_args.printfit())
482  cout << " directed step ";
483 
484  // (10a.i) Solve the linearized system for $\beta$ and
485  // the y-displacement vector $\delta$.
486  Column_Vector beta (nconstraints);
487  Column_Vector delta (nbadvars);
488  solve_linear_system (H, Y, save_By, save_negF,
489  beta, delta, W, U, V);
490 
491  // (10a.ii) Get the x-displacement vector $\gamma$.
492  Column_Vector gamma = -E * beta;
493 
494  // (10a.iii) Find the destination of the directed step.
495  x = c + xm + gamma;
496  y = d + ym + delta;
497 
498  // (10a.iv) Accept this point if it's not rejected by the constraint
499  // function, and the constraints improve.
500  if (call_constraint_fcn (constraint_calculator, x, y, F, Bx, By) &&
501  (constraint_sum = norm_infinity (F)) > 0 &&
502  (constraint_sum < constraint_sum_last)) {
503 
504  // Accept this step. Calculate the chisq and new displacement
505  // vectors.
506  chisq = chisq_last - scalar ((-save_negF + r*2) * beta);
507  c1 = x - xm;
508  d1 = y - ym;
509 
510  // Exit from step cutting loop.
511  break;
512  }
513  }
514 
515  // If this is the first time we're cutting the step,
516  // initialize $\psi$.
517  if (ncut == 0)
518  psi_cut = scalar ((save_negF - r) * alpha);
519 
520  // (10b) Give up if we've tried to cut this step too many times.
521  if (++ncut > _args.max_cut()) {
522  // cout << " Too many cut steps ";
523  // return -1000;
524  return -997.0;
525  }
526 
527  // (10c) Set up the size by which we're going to cut this step.
528  // Normally, this is cutsize. But if this is the first time we're
529  // cutting this step and the last step was also cut, set the cut
530  // size to twice the final cut size from the last step (provided
531  // that it is less than cutsize).
532  double this_cutsize = _args.cutsize();
533  if (ncut == 1 && last_step_cutsize < 1) {
534  this_cutsize = 2 * last_step_cutsize;
535  if (this_cutsize > _args.cutsize())
536  this_cutsize = _args.cutsize();
537  }
538 
539  // (10d) Keep track of the total amount by which we've cut this step.
540  this_step_cutsize *= this_cutsize;
541 
542  // If it falls below min_tot_cutsize, give up.
543  if (this_step_cutsize < _args.min_tot_cutsize()) {
544  // cout << "Cut size underflow ";
545  // return -1000;
546  return -996.0;
547  }
548 
549  // (10e) Cut the step: calculate the new displacement vectors.
550  double cutleft = 1 - this_cutsize;
551  c1 = c1 * this_cutsize + c * cutleft;
552  d1 = d1 * this_cutsize + d * cutleft;
553 
554  // (10f) Calculate the new chisq.
555  if (chisq_last >= 0) {
556  chisq = this_cutsize*this_cutsize * chisq +
557  cutleft*cutleft * chisq_last +
558  2*this_cutsize*cutleft * psi_cut;
559  psi_cut = this_cutsize * psi_cut + cutleft * chisq_last;
560  }
561  else
562  chisq = chisq_last;
563 
564  // Log what we've done.
565  if (_args.printfit()) {
566  cout << constraint_sum << " cut " << ncut << " size "
567  << setiosflags (ios_base::scientific)
568  << this_cutsize << " tot size " << this_step_cutsize
569  << resetiosflags (ios_base::scientific)
570  << " " << chisq << "\n";
571  }
572 
573  // Find the new step destination.
574  x = c1 + xm;
575  y = d1 + ym;
576 
577  // Now, go and test the step again for acceptability.
578  }
579 
580  // (11) At this point, we have an acceptable step.
581  // Shuffle things around to prepare for the next step.
582  last_step_cutsize = this_step_cutsize;
583 
584  // If requested, calculate the chisq using G to test for
585  // possible loss of precision.
586  double chisq_b = 0;
587  if (_args.use_G()) {
588  chisq_b = scalar (c1.T() * G * c1) + scalar (d1.T() * Y * d1);
589  if (chisq >= 0 &&
590  abs ((chisq - chisq_b) / chisq) > _args.chisq_test_eps()) {
591  cout << chisq << " " << chisq_b
592  << "lost precision?\n";
593  abort ();
594  }
595  }
596 
597  // Log what we're doing.
598  if (_args.printfit()) {
599  cout << chisq << " ";
600  if (_args.use_G())
601  cout << chisq_b << " ";
602  }
603 
604  double z2 = abs (chisq - chisq_last);
605 
606  if (_args.printfit()) {
607  cout << constraint_sum << " " << z2 << "\n";
608  }
609 
610  c = c1;
611  d = d1;
612  chisq_last = chisq;
613  constraint_sum_last = constraint_sum;
614 
615  // (12) Test for convergence. The conditions must be satisfied
616  // for two iterations in a row.
617  if (chisq >= 0 && constraint_sum < _args.constraint_sum_eps() &&
618  z2 < _args.chisq_diff_eps())
619  {
620  if (near_convergence) break; // Converged! Exit loop.
621  near_convergence = true;
622  }
623  else
624  near_convergence = false;
625 
626  // (13) Give up if we've done this too many times.
627  if (++nit > _args.maxit()) {
628  // cout << "too many iterations";
629  // return -1000;
630  return -995.0;
631  }
632 
633  } while (1);
634 
635  // (15) Ok, we have a successful fit!
636 
637 
638  // Calculate the error matrices.
639  Q = E * W * E.T();
640  S = - E * V.T();
641  R = U;
642 
643  // And the vectors of pull functions.
644  pullx = Column_Vector (nvars);
645  for (int i=1; i<=nvars; i++) {
646  double a = Q(i,i);
647  if (a < 0)
648  pullx(i) = c(i) / sqrt (-a);
649  else {
650  pullx(i) = 0;
651  // cout << " bad pull fcn for var " << i << " (" << a << ") ";
652  }
653  }
654 
655  pully = Column_Vector (nbadvars);
656  for (int i=1; i<=nbadvars; i++) {
657  double a = 1 - Y(i,i)*R(i,i);
658  if (a > 0)
659  pully(i) = d(i) * sqrt (Y(i,i) / a);
660  else {
661  pully(i) = 0;
662  // cout << " bad pull fcn for badvar " << i << " ";
663  }
664  }
665 
666  // Finish calculation of Q.
667  Q = Q + G_i;
668 
669  // Return the final chisq.
670  return chisq_last;
671 }
672 
673 
674 std::ostream& Chisq_Constrainer::print (std:: ostream& s) const
675 //
676 // Purpose: Print our state.
677 //
678 // Inputs:
679 // s - The stream to which to write.
680 //
681 // Returns:
682 // The stream S.
683 //
684 {
686  s << " printfit: " << _args.printfit()
687  << " use_G: " << _args.use_G() << "\n";
688  s << " constraint_sum_eps: " << _args.constraint_sum_eps()
689  << " chisq_diff_eps: " << _args.chisq_diff_eps()
690  << " chisq_test_eps: " << _args.chisq_test_eps() << "\n";
691  s << " maxit: " << _args.maxit()
692  << " max_cut: " << _args.max_cut()
693  << " min_tot_cutsize: " << _args.min_tot_cutsize()
694  << " cutsize: " << _args.cutsize() << "\n";
695  return s;
696 }
697 
698 
699 } // namespace hitfit
const double beta
dbl * delta
Definition: mlp_gen.cc:36
int i
Definition: DBlmapReader.cc:9
float alpha
Definition: AMPTWrapper.h:95
const Base_Constrainer_Args & base_constrainer_args() const
Define an abstract interface for getting parameter settings.
CLHEP::HepVector Column_Vector
Definition: matutil.h:67
#define abs(x)
Definition: mlp_lapack.h:159
Abstract base class for evaluating constraints. Users derive from this class and implement the eval()...
CLHEP::HepMatrix Matrix
Definition: matutil.h:66
virtual bool get_bool(std::string name) const =0
Chisq_Constrainer(const Chisq_Constrainer_Args &args)
Base class for constrained fitter.
Chisq_Constrainer_Args(const Defaults &defs)
Constructor, creates an instance of Chisq_Constrainer_Args from a Defaults object.
T sqrt(T t)
Definition: SSEVec.h:48
const Base_Constrainer_Args _base_constrainer_args
Minimize a subject to a set of constraints, based on the SQUAW algorithm.
virtual std::ostream & print(std::ostream &s) const
Print the state of this instance of Chisq_Constrainer.
int j
Definition: DBlmapReader.cc:9
bool call_constraint_fcn(Constraint_Calculator &constraint_calculator, const Column_Vector &x, const Column_Vector &y, Row_Vector &F, Matrix &Bx, Matrix &By) const
Helper function to evaluate constraints. This takes care of checking what the user function returns a...
Hold on to parameters for the Chisq_Constrainer class.
virtual double fit(Constraint_Calculator &constraint_calculator, const Column_Vector &xm, Column_Vector &x, const Column_Vector &ym, Column_Vector &y, const Matrix &G_i, const Diagonal_Matrix &Y, Column_Vector &pullx, Column_Vector &pully, Matrix &Q, Matrix &R, Matrix &S)
virtual double get_float(std::string name) const =0
Hold on to parameters for the Base_Constrainer class.
CLHEP::HepDiagMatrix Diagonal_Matrix
Definition: matutil.h:68
string const
Definition: compareJSON.py:14
dictionary args
Define an interface for getting parameter settings.
Definition: Defaults.h:62
double a
Definition: hdecay.h:121
const Chisq_Constrainer_Args _args
tuple cout
Definition: gather_cfg.py:121
Definition: DDAxes.h:10
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:281
virtual int get_int(std::string name) const =0
virtual std::ostream & print(std::ostream &s) const
Print out internal state to output stream.
Row-vector class. CLHEP doesn&#39;t have a row-vector class, so HitFit uses its own. This is only a simpl...
Definition: matutil.h:84
double scalar(const CLHEP::HepGenMatrix &m)
Return the matrix as a scalar. Raise an assertion if the matris is not .
Definition: matutil.cc:184