ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TDecompBase.cxx
Go to the documentation of this file.
1 // @(#)root/matrix:$Id$
2 // Authors: Fons Rademakers, Eddy Offermann Dec 2003
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 ///////////////////////////////////////////////////////////////////////////
13 // //
14 // Decomposition Base class //
15 // //
16 // This class forms the base for all the decompositions methods in the //
17 // linear algebra package . //
18 // It or its derived classes have installed the methods to solve //
19 // equations,invert matrices and calculate determinants while monitoring //
20 // the accuracy. //
21 // //
22 // Each derived class has always the following methods available: //
23 // //
24 // Condition() : //
25 // In an iterative scheme the condition number for matrix inversion is //
26 // calculated . This number is of interest for estimating the accuracy //
27 // of x in the equation Ax=b //
28 // For example: //
29 // A is a (10x10) Hilbert matrix which looks deceivingly innocent //
30 // and simple, A(i,j) = 1/(i+j+1) //
31 // b(i) = Sum_j A(i,j), so a sum of a row in A //
32 // //
33 // the solution is x(i) = 1. i=0,.,9 //
34 // //
35 // However, //
36 // TMatrixD m....; TVectorD b..... //
37 // TDecompLU lu(m); lu.SetTol(1.0e-12); lu.Solve(b); b.Print() //
38 // gives, //
39 // //
40 // {1.000,1.000,1.000,1.000,0.998,1.000,0.993,1.001,0.996,1.000} //
41 // //
42 // Looking at the condition number, this is in line with expected the //
43 // accuracy . The condition number is 3.957e+12 . As a simple rule of //
44 // thumb, a condition number of 1.0e+n means that you lose up to n //
45 // digits of accuracy in a solution . Since doubles are stored with 15 //
46 // digits, we can expect the accuracy to be as small as 3 digits . //
47 // //
48 // Det(Double_t &d1,Double_t &d2) //
49 // The determinant is d1*TMath::Power(2.,d2) //
50 // Expressing the determinant this way makes under/over-flow very //
51 // unlikely . //
52 // //
53 // Decompose() //
54 // Here the actually decomposition is performed . One can change the //
55 // matrix A after the decomposition constructor has been called //
56 // without effecting the decomposition result //
57 // //
58 // Solve(TVectorD &b) //
59 // Solve A x = b . x is supplied through the argument and replaced with //
60 // the solution . //
61 // //
62 // TransSolve(TVectorD &b) //
63 // Solve A^T x = b . x is supplied through the argument and replaced //
64 // with the solution . //
65 // //
66 // MultiSolve(TMatrixD &B) //
67 // Solve A X = B . where X and are now matrices . X is supplied through //
68 // the argument and replaced with the solution . //
69 // //
70 // Invert(TMatrixD &inv) //
71 // This is of course just a call to MultiSolve with as input argument //
72 // the unit matrix . Note that for a matrix a(m,n) with m > n a //
73 // pseudo-inverse is calculated . //
74 // //
75 // Tolerances and Scaling //
76 // ---------------------- //
77 // The tolerance parameter (which is a member of this base class) plays //
78 // a crucial role in all operations of the decomposition classes . It //
79 // gives the user a powerful tool to monitor and steer the operations //
80 // Its default value is sqrt(epsilon) where 1+epsilon = 1 //
81 // //
82 // If you do not want to be bothered by the following considerations, //
83 // like in most other linear algebra packages, just set the tolerance //
84 // with SetTol to an arbitrary small number . //
85 // //
86 // The tolerance number is used by each decomposition method to decide //
87 // whether the matrix is near singular, except of course SVD which can //
88 // handle singular matrices . //
89 // For each decomposition this will be checked in a different way; in LU //
90 // the matrix is considered singular when, at some point in the //
91 // decomposition, a diagonal element < fTol . Therefore, we had to set in//
92 // the example above of the (10x10) Hilbert, which is near singular, the //
93 // tolerance on 10e-12 . (The fact that we have to set the tolerance < //
94 // sqrt(epsilon) is a clear indication that we are losing precision .) //
95 // //
96 // If the matrix is flagged as being singular, operations with the //
97 // decomposition will fail and will return matrices/vectors that are //
98 // invalid . //
99 // //
100 // The observant reader will notice that by scaling the complete matrix //
101 // by some small number the decomposition will detect a singular matrix .//
102 // In this case the user will have to reduce the tolerance number by this//
103 // factor . (For CPU time saving we decided not to make this an automatic//
104 // procedure) . //
105 // //
106 // Code for this could look as follows: //
107 // const Double_t max_abs = Abs(a).Max(); //
108 // const Double_t scale = TMath::Min(max_abs,1.); //
109 // a.SetTol(a.GetTol()*scale); //
110 // //
111 // For usage examples see $ROOTSYS/test/stressLinear.cxx //
112 ///////////////////////////////////////////////////////////////////////////
113 
114 #include "TDecompBase.h"
115 #include "TMath.h"
116 #include "TError.h"
117 
119 
120 ////////////////////////////////////////////////////////////////////////////////
121 /// Default constructor
122 
124 {
126  fDet1 = 0;
127  fDet2 = 0;
128  fCondition = -1.0;
129  fRowLwb = 0;
130  fColLwb = 0;
131 }
132 
133 ////////////////////////////////////////////////////////////////////////////////
134 /// Copy constructor
135 
136 TDecompBase::TDecompBase(const TDecompBase &another) : TObject(another)
137 {
138  *this = another;
139 }
140 
141 ////////////////////////////////////////////////////////////////////////////////
142 
144 {
145 // Estimates lower bound for norm1 of inverse of A. Returns norm
146 // estimate in est. iter sets the maximum number of iterations to be used.
147 // The return value indicates the number of iterations remaining on exit from
148 // loop, hence if this is non-zero the processed "converged".
149 // This routine uses Hager's Convex Optimisation Algorithm.
150 // See Applied Numerical Linear Algebra, p139 & SIAM J Sci Stat Comp 1984 pp 311-16
151 
152  est = -1.0;
153 
154  const TMatrixDBase &m = GetDecompMatrix();
155  if (!m.IsValid())
156  return iter;
157 
158  const Int_t n = m.GetNrows();
159 
160  TVectorD b(n); TVectorD y(n); TVectorD z(n);
161  b = Double_t(1.0/n);
162  Double_t inv_norm1 = 0.0;
163  Bool_t stop = kFALSE;
164  do {
165  y = b;
166  if (!Solve(y))
167  return iter;
168  const Double_t ynorm1 = y.Norm1();
169  if ( ynorm1 <= inv_norm1 ) {
170  stop = kTRUE;
171  } else {
172  inv_norm1 = ynorm1;
173  Int_t i;
174  for (i = 0; i < n; i++)
175  z(i) = ( y(i) >= 0.0 ? 1.0 : -1.0 );
176  if (!TransSolve(z))
177  return iter;
178  Int_t imax = 0;
179  Double_t maxz = TMath::Abs(z(0));
180  for (i = 1; i < n; i++) {
181  const Double_t absz = TMath::Abs(z(i));
182  if ( absz > maxz ) {
183  maxz = absz;
184  imax = i;
185  }
186  }
187  stop = (maxz <= b*z);
188  if (!stop) {
189  b = 0.0;
190  b(imax) = 1.0;
191  }
192  }
193  iter--;
194  } while (!stop && iter);
195  est = inv_norm1;
196 
197  return iter;
198 }
199 
200 ////////////////////////////////////////////////////////////////////////////////
201 
203 {
204 // Returns product of matrix diagonal elements in d1 and d2. d1 is a mantissa and d2
205 // an exponent for powers of 2. If matrix is in diagonal or triangular-matrix form this
206 // will be the determinant.
207 // Based on Bowler, Martin, Peters and Wilkinson in HACLA
208 
209  const Double_t zero = 0.0;
210  const Double_t one = 1.0;
211  const Double_t four = 4.0;
212  const Double_t sixteen = 16.0;
213  const Double_t sixteenth = 0.0625;
214 
215  const Int_t n = diag.GetNrows();
216 
217  Double_t t1 = 1.0;
218  Double_t t2 = 0.0;
219  Int_t niter2 =0;
220  Int_t niter3 =0;
221  for (Int_t i = 0; (((i < n) && (t1 !=zero ))); i++) {
222  if (TMath::Abs(diag(i)) > tol) {
223  t1 *= (Double_t) diag(i);
224  while ( TMath::Abs(t1) < one) {
225  t1 *= sixteenth;
226  t2 += four;
227  niter2++;
228  if ( niter2>100) break;
229  }
230  while ( TMath::Abs(t1) < sixteenth) {
231  t1 *= sixteen;
232  t2 -= four;
233  niter3++;
234  if (niter3>100) break;
235  }
236  } else {
237  t1 = zero;
238  t2 = zero;
239  }
240  }
241  d1 = t1;
242  d2 = t2;
243 
244  return;
245 }
246 
247 ////////////////////////////////////////////////////////////////////////////////
248 /// Matrix condition number
249 
251 {
252  if ( !TestBit(kCondition) ) {
253  fCondition = -1;
254  if (TestBit(kSingular))
255  return fCondition;
256  if ( !TestBit(kDecomposed) ) {
257  if (!Decompose())
258  return fCondition;
259  }
260  Double_t invNorm;
261  if (Hager(invNorm))
262  fCondition *= invNorm;
263  else // no convergence in Hager
264  Error("Condition()","Hager procedure did NOT converge");
266  }
267  return fCondition;
268 }
269 
270 ////////////////////////////////////////////////////////////////////////////////
271 /// Solve set of equations with RHS in columns of B
272 
274 {
275  const TMatrixDBase &m = GetDecompMatrix();
276  R__ASSERT(m.IsValid() && B.IsValid());
277 
278  const Int_t colLwb = B.GetColLwb();
279  const Int_t colUpb = B.GetColUpb();
280  Bool_t status = kTRUE;
281  for (Int_t icol = colLwb; icol <= colUpb && status; icol++) {
282  TMatrixDColumn b(B,icol);
283  status &= Solve(b);
284  }
285 
286  return status;
287 }
288 
289 ////////////////////////////////////////////////////////////////////////////////
290 /// Matrix determinant det = d1*TMath::Power(2.,d2)
291 
293 {
294  if ( !TestBit(kDetermined) ) {
295  if ( !TestBit(kDecomposed) )
296  Decompose();
297  if (TestBit(kSingular) ) {
298  fDet1 = 0.0;
299  fDet2 = 0.0;
300  } else {
301  const TMatrixDBase &m = GetDecompMatrix();
302  R__ASSERT(m.IsValid());
303  TVectorD diagv(m.GetNrows());
304  for (Int_t i = 0; i < diagv.GetNrows(); i++)
305  diagv(i) = m(i,i);
306  DiagProd(diagv,fTol,fDet1,fDet2);
307  }
309  }
310  d1 = fDet1;
311  d2 = fDet2;
312 }
313 
314 ////////////////////////////////////////////////////////////////////////////////
315 /// Print class members
316 
317 void TDecompBase::Print(Option_t * /*opt*/) const
318 {
319  printf("fTol = %.4e\n",fTol);
320  printf("fDet1 = %.4e\n",fDet1);
321  printf("fDet2 = %.4e\n",fDet2);
322  printf("fCondition = %.4e\n",fCondition);
323  printf("fRowLwb = %d\n",fRowLwb);
324  printf("fColLwb = %d\n",fColLwb);
325 }
326 
327 ////////////////////////////////////////////////////////////////////////////////
328 /// Assignment operator
329 
331 {
332  if (this != &source) {
333  TObject::operator=(source);
334  fTol = source.fTol;
335  fDet1 = source.fDet1;
336  fDet2 = source.fDet2;
337  fCondition = source.fCondition;
338  fRowLwb = source.fRowLwb;
339  fColLwb = source.fColLwb;
340  }
341  return *this;
342 }
343 
344 ////////////////////////////////////////////////////////////////////////////////
345 /// Define a Householder-transformation through the parameters up and b .
346 
348  Double_t tol)
349 {
350  const Int_t n = vc.GetNrows();
351  const Double_t * const vp = vc.GetMatrixArray();
352 
353  Double_t c = TMath::Abs(vp[lp]);
354  Int_t i;
355  for (i = l; i < n; i++)
356  c = TMath::Max(TMath::Abs(vp[i]),c);
357 
358  up = 0.0;
359  beta = 0.0;
360  if (c <= tol) {
361 // Warning("DefHouseHolder","max vector=%.4e < %.4e",c,tol);
362  return kFALSE;
363  }
364 
365  Double_t sd = vp[lp]/c; sd *= sd;
366  for (i = l; i < n; i++) {
367  const Double_t tmp = vp[i]/c;
368  sd += tmp*tmp;
369  }
370 
371  Double_t vpprim = c*TMath::Sqrt(sd);
372  if (vp[lp] > 0.) vpprim = -vpprim;
373  up = vp[lp]-vpprim;
374  beta = 1./(vpprim*up);
375 
376  return kTRUE;
377 }
378 
379 ////////////////////////////////////////////////////////////////////////////////
380 /// Apply Householder-transformation.
381 
383  Int_t lp,Int_t l,TMatrixDRow &cr)
384 {
385  const Int_t nv = vc.GetNrows();
386  const Int_t nc = (cr.GetMatrix())->GetNcols();
387 
388  if (nv > nc) {
389  Error("ApplyHouseHolder(const TVectorD &,..,TMatrixDRow &)","matrix row too short");
390  return;
391  }
392 
393  const Int_t inc_c = cr.GetInc();
394  const Double_t * const vp = vc.GetMatrixArray();
395  Double_t * cp = cr.GetPtr();
396 
397  Double_t s = cp[lp*inc_c]*up;
398  Int_t i;
399  for (i = l; i < nv; i++)
400  s += cp[i*inc_c]*vp[i];
401 
402  s = s*beta;
403  cp[lp*inc_c] += s*up;
404  for (i = l; i < nv; i++)
405  cp[i*inc_c] += s*vp[i];
406 }
407 
408 ////////////////////////////////////////////////////////////////////////////////
409 /// Apply Householder-transformation.
410 
412  Int_t lp,Int_t l,TMatrixDColumn &cc)
413 {
414  const Int_t nv = vc.GetNrows();
415  const Int_t nc = (cc.GetMatrix())->GetNrows();
416 
417  if (nv > nc) {
418  Error("ApplyHouseHolder(const TVectorD &,..,TMatrixDRow &)","matrix column too short");
419  return;
420  }
421 
422  const Int_t inc_c = cc.GetInc();
423  const Double_t * const vp = vc.GetMatrixArray();
424  Double_t * cp = cc.GetPtr();
425 
426  Double_t s = cp[lp*inc_c]*up;
427  Int_t i;
428  for (i = l; i < nv; i++)
429  s += cp[i*inc_c]*vp[i];
430 
431  s = s*beta;
432  cp[lp*inc_c] += s*up;
433  for (i = l; i < nv; i++)
434  cp[i*inc_c] += s*vp[i];
435 }
436 
437 ////////////////////////////////////////////////////////////////////////////////
438 /// Apply Householder-transformation.
439 
441  Int_t lp,Int_t l,TVectorD &cv)
442 {
443  const Int_t nv = vc.GetNrows();
444  const Int_t nc = cv.GetNrows();
445 
446  if (nv > nc) {
447  Error("ApplyHouseHolder(const TVectorD &,..,TVectorD &)","vector too short");
448  return;
449  }
450 
451  const Double_t * const vp = vc.GetMatrixArray();
452  Double_t * cp = cv.GetMatrixArray();
453 
454  Double_t s = cp[lp]*up;
455  Int_t i;
456  for (i = l; i < nv; i++)
457  s += cp[i]*vp[i];
458 
459  s = s*beta;
460  cp[lp] += s*up;
461  for (i = l; i < nv; i++)
462  cp[i] += s*vp[i];
463 }
464 
465 ////////////////////////////////////////////////////////////////////////////////
466 /// Defines a Givens-rotation by calculating 2 rotation parameters c and s.
467 /// The rotation is defined with the vector components v1 and v2.
468 
470 {
471  const Double_t a1 = TMath::Abs(v1);
472  const Double_t a2 = TMath::Abs(v2);
473  if (a1 > a2) {
474  const Double_t w = v2/v1;
475  const Double_t q = TMath::Hypot(1.,w);
476  c = 1./q;
477  if (v1 < 0.) c = -c;
478  s = c*w;
479  } else {
480  if (v2 != 0) {
481  const Double_t w = v1/v2;
482  const Double_t q = TMath::Hypot(1.,w);
483  s = 1./q;
484  if (v2 < 0.) s = -s;
485  c = s*w;
486  } else {
487  c = 1.;
488  s = 0.;
489  }
490  }
491 }
492 
493 ////////////////////////////////////////////////////////////////////////////////
494 /// Define and apply a Givens-rotation by calculating 2 rotation
495 /// parameters c and s. The rotation is defined with and applied to the vector
496 /// components v1 and v2.
497 
499 {
500  const Double_t a1 = TMath::Abs(v1);
501  const Double_t a2 = TMath::Abs(v2);
502  if (a1 > a2) {
503  const Double_t w = v2/v1;
504  const Double_t q = TMath::Hypot(1.,w);
505  c = 1./q;
506  if (v1 < 0.) c = -c;
507  s = c*w;
508  v1 = a1*q;
509  v2 = 0.;
510  } else {
511  if (v2 != 0) {
512  const Double_t w = v1/v2;
513  const Double_t q = TMath::Hypot(1.,w);
514  s = 1./q;
515  if (v2 < 0.) s = -s;
516  c = s*w;
517  v1 = a2*q;
518  v2 = 0.;
519  } else {
520  c = 1.;
521  s = 0.;
522  }
523  }
524 }
525 
526 ////////////////////////////////////////////////////////////////////////////////
527 /// Apply a Givens transformation as defined by c and s to the vector compenents
528 /// v1 and v2 .
529 
531 {
532  const Double_t w = z1*c+z2*s;
533  z2 = -z1*s+z2*c;
534  z1 = w;
535 }
static double B[]
virtual Bool_t MultiSolve(TMatrixD &B)
Solve set of equations with RHS in columns of B.
Element * GetPtr() const
const Double_t * v1
Definition: TArcBall.cxx:33
virtual Bool_t TransSolve(TVectorD &b)=0
Int_t Hager(Double_t &est, Int_t iter=5)
Element * GetPtr() const
Double_t fCondition
Definition: TDecompBase.h:42
Int_t fRowLwb
Definition: TDecompBase.h:43
return c
const char Option_t
Definition: RtypesCore.h:62
Int_t GetNrows() const
Definition: TVectorT.h:81
void DefAplGivens(Double_t &v1, Double_t &v2, Double_t &c, Double_t &s)
Define and apply a Givens-rotation by calculating 2 rotation parameters c and s.
void ApplyHouseHolder(const TVectorD &vc, Double_t up, Double_t beta, Int_t lp, Int_t l, TMatrixDRow &cr)
Apply Householder-transformation.
virtual Bool_t Decompose()=0
void Print(Option_t *opt="") const
Print class members.
#define R__ASSERT(e)
Definition: TError.h:98
virtual const TMatrixDBase & GetDecompMatrix() const =0
TAlienJobStatus * status
Definition: TAlienJob.cxx:51
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
TLatex * t1
Definition: textangle.C:20
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
const TMatrixTBase< Element > * GetMatrix() const
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:732
double beta(double x, double y)
Calculates the beta function.
Double_t fTol
Definition: TDecompBase.h:39
Bool_t IsValid() const
Definition: TMatrixTBase.h:157
Bool_t DefHouseHolder(const TVectorD &vc, Int_t lp, Int_t l, Double_t &up, Double_t &beta, Double_t tol)
Define a Householder-transformation through the parameters up and b .
TDecompBase & operator=(const TDecompBase &source)
Assignment operator.
std::map< std::string, std::string >::const_iterator iter
Definition: TAlienJob.cxx:54
Int_t fColLwb
Definition: TDecompBase.h:44
Int_t GetInc() const
TObject & operator=(const TObject &rhs)
TObject assignment operator.
Definition: TObject.cxx:102
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:918
Float_t z[5]
Definition: Ifit.C:16
void Error(const char *location, const char *msgfmt,...)
static void DiagProd(const TVectorD &diag, Double_t tol, Double_t &d1, Double_t &d2)
Int_t GetColUpb() const
Definition: TMatrixTBase.h:136
Int_t GetNrows() const
Definition: TMatrixTBase.h:134
Element * GetMatrixArray()
Definition: TVectorT.h:84
const double tol
virtual Bool_t Solve(TVectorD &b)=0
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:173
TMarker * m
Definition: textangle.C:8
tuple w
Definition: qtexample.py:51
TLine * l
Definition: textangle.C:4
void ApplyGivens(Double_t &z1, Double_t &z2, Double_t c, Double_t s)
Apply a Givens transformation as defined by c and s to the vector compenents v1 and v2 ...
void DefGivens(Double_t v1, Double_t v2, Double_t &c, Double_t &s)
Defines a Givens-rotation by calculating 2 rotation parameters c and s.
const TMatrixTBase< Element > * GetMatrix() const
Int_t GetColLwb() const
Definition: TMatrixTBase.h:135
REAL epsilon
Definition: triangle.c:617
Double_t fDet2
Definition: TDecompBase.h:41
TText * t2
Definition: rootenv.C:28
Double_t fDet1
Definition: TDecompBase.h:40
double Double_t
Definition: RtypesCore.h:55
Element Norm1() const
Compute the 1-norm of the vector SUM{ |v[i]| }.
Definition: TVectorT.cxx:564
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
Double_t y[n]
Definition: legend1.C:17
Double_t Hypot(Double_t x, Double_t y)
Definition: TMath.cxx:60
virtual void Det(Double_t &d1, Double_t &d2)
Matrix determinant det = d1*TMath::Power(2.,d2)
Mother of all ROOT objects.
Definition: TObject.h:58
virtual Double_t Condition()
Matrix condition number.
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
ClassImp(TDecompBase) TDecompBase
Default constructor.
Int_t GetInc() const
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
const Bool_t kTRUE
Definition: Rtypes.h:91
float * q
Definition: THbookFile.cxx:87
const Int_t n
Definition: legend1.C:16