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