ROOT  6.06/09
Reference Guide
TDecompBK.cxx
Go to the documentation of this file.
1 // @(#)root/matrix:$Id$
2 // Authors: Fons Rademakers, Eddy Offermann Sep 2004
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 #include "TDecompBK.h"
13 #include "TMath.h"
14 #include "TError.h"
15 
17 
18 ///////////////////////////////////////////////////////////////////////////
19 // //
20 // The Bunch-Kaufman diagonal pivoting method decomposes a real //
21 // symmetric matrix A using //
22 // //
23 // A = U*D*U^T //
24 // //
25 // where U is a product of permutation and unit upper triangular //
26 // matrices, U^T is the transpose of U, and D is symmetric and block //
27 // diagonal with 1-by-1 and 2-by-2 diagonal blocks. //
28 // //
29 // U = P(n-1)*U(n-1)* ... *P(k)U(k)* ..., //
30 // i.e., U is a product of terms P(k)*U(k), where k decreases from n-1 //
31 // to 0 in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1//
32 // and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as //
33 // defined by IPIV(k), and U(k) is a unit upper triangular matrix, such //
34 // that if the diagonal block D(k) is of order s (s = 1 or 2), then //
35 // //
36 // ( I v 0 ) k-s //
37 // U(k) = ( 0 I 0 ) s //
38 // ( 0 0 I ) n-k //
39 // k-s s n-k //
40 // //
41 // If s = 1, D(k) overwrites A(k,k), and v overwrites A(0:k-1,k). //
42 // If s = 2, the upper triangle of D(k) overwrites A(k-1,k-1), A(k-1,k),//
43 // and A(k,k), and v overwrites A(0:k-2,k-1:k). //
44 // //
45 // fU contains on entry the symmetric matrix A of which only the upper //
46 // triangular part is referenced . On exit fU contains the block diagonal//
47 // matrix D and the multipliers used to obtain the factor U, see above . //
48 // //
49 // fIpiv if dimension n contains details of the interchanges and the //
50 // the block structure of D . If (fIPiv(k) > 0, then rows and columns k //
51 // and fIPiv(k) were interchanged and D(k,k) is a 1-by-1 diagonal block. //
52 // If IPiv(k) = fIPiv(k-1) < 0, rows and columns k-1 and -IPiv(k) were //
53 // interchanged and D(k-1:k,k-1:k) is a 2-by-2 diagonal block. //
54 // //
55 ///////////////////////////////////////////////////////////////////////////
56 
57 ////////////////////////////////////////////////////////////////////////////////
58 /// Default constructor
59 
61 {
62  fNIpiv = 0;
63  fIpiv = 0;
64 }
65 
66 ////////////////////////////////////////////////////////////////////////////////
67 /// Constructor for (nrows x nrows) symmetric matrix
68 
70 {
71  fNIpiv = nrows;
72  fIpiv = new Int_t[fNIpiv];
73  memset(fIpiv,0,fNIpiv*sizeof(Int_t));
74  fU.ResizeTo(nrows,nrows);
75 }
76 
77 ////////////////////////////////////////////////////////////////////////////////
78 /// Constructor for ([row_lwb..row_upb] x [row_lwb..row_upb]) symmetric matrix
79 
81 {
82  const Int_t nrows = row_upb-row_lwb+1;
83  fNIpiv = nrows;
84  fIpiv = new Int_t[fNIpiv];
85  memset(fIpiv,0,fNIpiv*sizeof(Int_t));
86  fColLwb = fRowLwb = row_lwb;
87  fU.ResizeTo(nrows,nrows);
88 }
89 
90 ////////////////////////////////////////////////////////////////////////////////
91 /// Constructor for symmetric matrix A
92 
94 {
95  R__ASSERT(a.IsValid());
96 
98  fCondition = a.Norm1();
99  fTol = a.GetTol();
100  if (tol > 0)
101  fTol = tol;
102 
103  fNIpiv = a.GetNcols();
104  fIpiv = new Int_t[fNIpiv];
105  memset(fIpiv,0,fNIpiv*sizeof(Int_t));
106 
107  const Int_t nRows = a.GetNrows();
108  fColLwb = fRowLwb = a.GetRowLwb();
109  fU.ResizeTo(nRows,nRows);
110  memcpy(fU.GetMatrixArray(),a.GetMatrixArray(),nRows*nRows*sizeof(Double_t));
111 }
112 
113 ////////////////////////////////////////////////////////////////////////////////
114 /// Copy constructor
115 
116 TDecompBK::TDecompBK(const TDecompBK &another) : TDecompBase(another)
117 {
118  fNIpiv = 0;
119  fIpiv = 0;
120  *this = another;
121 }
122 
123 ////////////////////////////////////////////////////////////////////////////////
124 /// Matrix A is decomposed in components U and D so that A = U*D*U^T
125 /// If the decomposition succeeds, bit kDecomposed is set , otherwise kSingular
126 
128 {
129  if (TestBit(kDecomposed)) return kTRUE;
130 
131  if ( !TestBit(kMatrixSet) ) {
132  Error("Decompose()","Matrix has not been set");
133  return kFALSE;
134  }
135 
136  Bool_t ok = kTRUE;
137 
138 // Initialize alpha for use in choosing pivot block size.
139  const Double_t alpha = (1.+TMath::Sqrt(17.))/8.;
140 
141 // Factorize a as u*d*u' using the upper triangle of a .
142 // k is the main loop index, decreasing from n-1 to 0 in steps of 1 or 2
143 
144  const Int_t n = fU.GetNcols();
145  Double_t *pU = fU.GetMatrixArray();
146  TMatrixDDiag diag(fU);
147 
148  Int_t imax = 0;
149  Int_t k = n-1;
150  while (k >= 0) {
151  Int_t kstep = 1;
152 
153  // determine rows and columns to be interchanged and whether
154  // a 1-by-1 or 2-by-2 pivot block will be used
155 
156  const Double_t absakk = TMath::Abs(diag(k));
157 
158  // imax is the row-index of the largest off-diagonal element in
159  // column k, and colmax is its absolute value
160 
161  Double_t colmax;
162  if ( k > 0 ) {
163  TVectorD vcol = TMatrixDColumn_const(fU,k);
164  vcol.Abs();
165  imax = TMath::LocMax(k,vcol.GetMatrixArray());
166  colmax = vcol[imax];
167  } else {
168  colmax = 0.;
169  }
170 
171  Int_t kp;
172  if (TMath::Max(absakk,colmax) <= fTol) {
173  // the block diagonal matrix will be singular
174  kp = k;
175  ok = kFALSE;
176  } else {
177  if (absakk >= alpha*colmax) {
178  // no interchange, use 1-by-1 pivot block
179  kp = k;
180  } else {
181  // jmax is the column-index of the largest off-diagonal
182  // element in row imax, and rowmax is its absolute value
183  TVectorD vrow = TMatrixDRow_const(fU,imax);
184  vrow.Abs();
185  Int_t jmax = imax+1+TMath::LocMax(k-imax,vrow.GetMatrixArray()+imax+1);
186  Double_t rowmax = vrow[jmax];
187  if (imax > 0) {
188  TVectorD vcol = TMatrixDColumn_const(fU,imax);
189  vcol.Abs();
190  jmax = TMath::LocMax(imax,vcol.GetMatrixArray());
191  rowmax = TMath::Max(rowmax,vcol[jmax]);
192  }
193 
194  if (absakk >= alpha*colmax*(colmax/rowmax)) {
195  // No interchange, use 1-by-1 pivot block
196  kp = k;
197  } else if( TMath::Abs(diag(imax)) >= alpha*rowmax) {
198  // Interchange rows and columns k and imax, use 1-by-1 pivot block
199  kp = imax;
200  } else {
201  // Interchange rows and columns k-1 and imax, use 2-by-2 pivot block
202  kp = imax;
203  kstep = 2;
204  }
205  }
206 
207  const Int_t kk = k-kstep+1;
208  if (kp != kk) {
209  // Interchange rows and columns kk and kp in the leading submatrix a(0:k,0:k)
210  Double_t *c_kk = pU+kk;
211  Double_t *c_kp = pU+kp;
212  for (Int_t irow = 0; irow < kp; irow++) {
213  const Double_t t = *c_kk;
214  *c_kk = *c_kp;
215  *c_kp = t;
216  c_kk += n;
217  c_kp += n;
218  }
219 
220  c_kk = pU+(kp+1)*n+kk;
221  Double_t *r_kp = pU+kp*n+kp+1;
222  for (Int_t icol = 0; icol < kk-kp-1; icol++) {
223  const Double_t t = *c_kk;
224  *c_kk = *r_kp;
225  *r_kp = t;
226  c_kk += n;
227  r_kp += 1;
228  }
229 
230  Double_t t = diag(kk);
231  diag(kk) = diag(kp);
232  diag(kp) = t;
233  if (kstep == 2) {
234  t = pU[(k-1)*n+k];
235  pU[(k-1)*n+k] = pU[kp*n+k];
236  pU[kp*n+k] = t;
237  }
238  }
239 
240  // Update the leading submatrix
241 
242  if (kstep == 1 && k > 0) {
243  // 1-by-1 pivot block d(k): column k now holds w(k) = u(k)*d(k)
244  // where u(k) is the k-th column of u
245 
246  // perform a rank-1 update of a(0:k-1,0:k-1) as
247  // a := a - u(k)*d(k)*u(k)' = a - w(k)*1/d(k)*w(k)'
248 
249  const Double_t r1 = 1./diag(k);
250  TMatrixDSub sub1(fU,0,k-1,0,k-1);
251  sub1.Rank1Update(TMatrixDColumn_const(fU,k),-r1);
252 
253  // store u(k) in column k
254  TMatrixDSub sub2(fU,0,k-1,k,k);
255  sub2 *= r1;
256  } else {
257  // 2-by-2 pivot block d(k): columns k and k-1 now hold
258  // ( w(k-1) w(k) ) = ( u(k-1) u(k) )*d(k)
259  // where u(k) and u(k-1) are the k-th and (k-1)-th columns of u
260 
261  // perform a rank-2 update of a(0:k-2,0:k-2) as
262  // a := a - ( u(k-1) u(k) )*d(k)*( u(k-1) u(k) )'
263  // = a - ( w(k-1) w(k) )*inv(d(k))*( w(k-1) w(k) )'
264 
265  if ( k > 1 ) {
266  Double_t *pU_k1 = pU+(k-1)*n;
267  Double_t d12 = pU_k1[k];
268  const Double_t d22 = pU_k1[k-1]/d12;
269  const Double_t d11 = diag(k)/d12;
270  const Double_t t = 1./(d11*d22-1.);
271  d12 = t/d12;
272 
273  for (Int_t j = k-2; j >= 0; j--) {
274  Double_t *pU_j = pU+j*n;
275  const Double_t wkm1 = d12*(d11*pU_j[k-1]-pU_j[k]);
276  const Double_t wk = d12*(d22*pU_j[k]-pU_j[k-1]);
277  for (Int_t i = j; i >= 0; i--) {
278  Double_t *pU_i = pU+i*n;
279  pU_i[j] -= (pU_i[k]*wk+pU_i[k-1]*wkm1);
280  }
281  pU_j[k] = wk;
282  pU_j[k-1] = wkm1;
283  }
284  }
285  }
286 
287  // Store details of the interchanges in fIpiv
288  if (kstep == 1) {
289  fIpiv[k] = (kp+1);
290  } else {
291  fIpiv[k] = -(kp+1);
292  fIpiv[k-1] = -(kp+1);
293  }
294  }
295 
296  k -= kstep;
297  }
298 
299  if (!ok) SetBit(kSingular);
300  else SetBit(kDecomposed);
301 
303 
304  return ok;
305 }
306 
307 ////////////////////////////////////////////////////////////////////////////////
308 /// Set the matrix to be decomposed, decomposition status is reset.
309 
311 {
312  R__ASSERT(a.IsValid());
313 
314  ResetStatus();
315 
317  fCondition = a.Norm1();
318 
319  if (fNIpiv != a.GetNcols()) {
320  fNIpiv = a.GetNcols();
321  delete [] fIpiv;
322  fIpiv = new Int_t[fNIpiv];
323  memset(fIpiv,0,fNIpiv*sizeof(Int_t));
324  }
325 
326  const Int_t nRows = a.GetNrows();
327  fColLwb = fRowLwb = a.GetRowLwb();
328  fU.ResizeTo(nRows,nRows);
329  memcpy(fU.GetMatrixArray(),a.GetMatrixArray(),nRows*nRows*sizeof(Double_t));
330 }
331 
332 ////////////////////////////////////////////////////////////////////////////////
333 /// Solve Ax=b assuming the BK form of A is stored in fU . Solution returned in b.
334 
336 {
337  R__ASSERT(b.IsValid());
338  if (TestBit(kSingular)) {
339  Error("Solve()","Matrix is singular");
340  return kFALSE;
341  }
342  if ( !TestBit(kDecomposed) ) {
343  if (!Decompose()) {
344  Error("Solve()","Decomposition failed");
345  return kFALSE;
346  }
347  }
348 
349  if (fU.GetNrows() != b.GetNrows() || fU.GetRowLwb() != b.GetLwb()) {
350  Error("Solve(TVectorD &","vector and matrix incompatible");
351  return kFALSE;
352  }
353 
354  const Int_t n = fU.GetNrows();
355 
356  TMatrixDDiag_const diag(fU);
357  const Double_t *pU = fU.GetMatrixArray();
358  Double_t *pb = b.GetMatrixArray();
359 
360  // solve a*x = b, where a = u*d*u'. First solve u*d*x = b, overwriting b with x.
361  // k is the main loop index, decreasing from n-1 to 0 in steps of 1 or 2,
362  // depending on the size of the diagonal blocks.
363 
364  Int_t k = n-1;
365  while (k >= 0) {
366 
367  if (fIpiv[k] > 0) {
368 
369  // 1 x 1 diagonal block
370  // interchange rows k and ipiv(k).
371  const Int_t kp = fIpiv[k]-1;
372  if (kp != k) {
373  const Double_t tmp = pb[k];
374  pb[k] = pb[kp];
375  pb[kp] = tmp;
376  }
377 
378  // multiply by inv(u(k)), where u(k) is the transformation
379  // stored in column k of a.
380  for (Int_t i = 0; i < k; i++)
381  pb[i] -= pU[i*n+k]*pb[k];
382 
383  // multiply by the inverse of the diagonal block.
384  pb[k] /= diag(k);
385  k--;
386  } else {
387 
388  // 2 x 2 diagonal block
389  // interchange rows k-1 and -ipiv(k).
390  const Int_t kp = -fIpiv[k]-1;
391  if (kp != k-1) {
392  const Double_t tmp = pb[k-1];
393  pb[k-1] = pb[kp];
394  pb[kp] = tmp;
395  }
396 
397  // multiply by inv(u(k)), where u(k) is the transformation
398  // stored in columns k-1 and k of a.
399  Int_t i;
400  for (i = 0; i < k-1; i++)
401  pb[i] -= pU[i*n+k]*pb[k];
402  for (i = 0; i < k-1; i++)
403  pb[i] -= pU[i*n+k-1]*pb[k-1];
404 
405  // multiply by the inverse of the diagonal block.
406  const Double_t *pU_k1 = pU+(k-1)*n;
407  const Double_t ukm1k = pU_k1[k];
408  const Double_t ukm1 = pU_k1[k-1]/ukm1k;
409  const Double_t uk = diag(k)/ukm1k;
410  const Double_t denom = ukm1*uk-1.;
411  const Double_t bkm1 = pb[k-1]/ukm1k;
412  const Double_t bk = pb[k]/ukm1k;
413  pb[k-1] = (uk*bkm1-bk)/denom;
414  pb[k] = (ukm1*bk-bkm1)/denom;
415  k -= 2;
416  }
417  }
418 
419  // Next solve u'*x = b, overwriting b with x.
420  //
421  // k is the main loop index, increasing from 0 to n-1 in steps of
422  // 1 or 2, depending on the size of the diagonal blocks.
423 
424  k = 0;
425  while (k < n) {
426 
427  if (fIpiv[k] > 0) {
428  // 1 x 1 diagonal block
429  // multiply by inv(u'(k)), where u(k) is the transformation
430  // stored in column k of a.
431  for (Int_t i = 0; i < k; i++)
432  pb[k] -= pU[i*n+k]*pb[i];
433 
434  // interchange elements k and ipiv(k).
435  const Int_t kp = fIpiv[k]-1;
436  if (kp != k) {
437  const Double_t tmp = pb[k];
438  pb[k] = pb[kp];
439  pb[kp] = tmp;
440  }
441  k++;
442  } else {
443  // 2 x 2 diagonal block
444  // multiply by inv(u'(k+1)), where u(k+1) is the transformation
445  // stored in columns k and k+1 of a.
446  Int_t i ;
447  for (i = 0; i < k; i++)
448  pb[k] -= pU[i*n+k]*pb[i];
449  for (i = 0; i < k; i++)
450  pb[k+1] -= pU[i*n+k+1]*pb[i];
451 
452  // interchange elements k and -ipiv(k).
453  const Int_t kp = -fIpiv[k]-1;
454  if (kp != k) {
455  const Double_t tmp = pb[k];
456  pb[k] = pb[kp];
457  pb[kp] = tmp;
458  }
459  k += 2;
460  }
461  }
462 
463  return kTRUE;
464 }
465 
466 ////////////////////////////////////////////////////////////////////////////////
467 /// Solve Ax=b assuming the BK form of A is stored in fU . Solution returned in b.
468 
470 {
471  TMatrixDBase *b = const_cast<TMatrixDBase *>(cb.GetMatrix());
472  R__ASSERT(b->IsValid());
473  if (TestBit(kSingular)) {
474  Error("Solve()","Matrix is singular");
475  return kFALSE;
476  }
477  if ( !TestBit(kDecomposed) ) {
478  if (!Decompose()) {
479  Error("Solve()","Decomposition failed");
480  return kFALSE;
481  }
482  }
483 
484  if (fU.GetNrows() != b->GetNrows() || fU.GetRowLwb() != b->GetRowLwb()) {
485  Error("Solve(TMatrixDColumn &","vector and matrix incompatible");
486  return kFALSE;
487  }
488 
489  const Int_t n = fU.GetNrows();
490 
491  TMatrixDDiag_const diag(fU);
492  const Double_t *pU = fU.GetMatrixArray();
493  Double_t *pcb = cb.GetPtr();
494  const Int_t inc = cb.GetInc();
495 
496 
497  // solve a*x = b, where a = u*d*u'. First solve u*d*x = b, overwriting b with x.
498  // k is the main loop index, decreasing from n-1 to 0 in steps of 1 or 2,
499  // depending on the size of the diagonal blocks.
500 
501  Int_t k = n-1;
502  while (k >= 0) {
503 
504  if (fIpiv[k] > 0) {
505 
506  // 1 x 1 diagonal block
507  // interchange rows k and ipiv(k).
508  const Int_t kp = fIpiv[k]-1;
509  if (kp != k) {
510  const Double_t tmp = pcb[k*inc];
511  pcb[k*inc] = pcb[kp*inc];
512  pcb[kp*inc] = tmp;
513  }
514 
515  // multiply by inv(u(k)), where u(k) is the transformation
516  // stored in column k of a.
517  for (Int_t i = 0; i < k; i++)
518  pcb[i*inc] -= pU[i*n+k]*pcb[k*inc];
519 
520  // multiply by the inverse of the diagonal block.
521  pcb[k*inc] /= diag(k);
522  k--;
523  } else {
524 
525  // 2 x 2 diagonal block
526  // interchange rows k-1 and -ipiv(k).
527  const Int_t kp = -fIpiv[k]-1;
528  if (kp != k-1) {
529  const Double_t tmp = pcb[(k-1)*inc];
530  pcb[(k-1)*inc] = pcb[kp*inc];
531  pcb[kp*inc] = tmp;
532  }
533 
534  // multiply by inv(u(k)), where u(k) is the transformation
535  // stored in columns k-1 and k of a.
536  Int_t i;
537  for (i = 0; i < k-1; i++)
538  pcb[i*inc] -= pU[i*n+k]*pcb[k*inc];
539  for (i = 0; i < k-1; i++)
540  pcb[i*inc] -= pU[i*n+k-1]*pcb[(k-1)*inc];
541 
542  // multiply by the inverse of the diagonal block.
543  const Double_t *pU_k1 = pU+(k-1)*n;
544  const Double_t ukm1k = pU_k1[k];
545  const Double_t ukm1 = pU_k1[k-1]/ukm1k;
546  const Double_t uk = diag(k)/ukm1k;
547  const Double_t denom = ukm1*uk-1.;
548  const Double_t bkm1 = pcb[(k-1)*inc]/ukm1k;
549  const Double_t bk = pcb[k*inc]/ukm1k;
550  pcb[(k-1)*inc] = (uk*bkm1-bk)/denom;
551  pcb[k*inc] = (ukm1*bk-bkm1)/denom;
552  k -= 2;
553  }
554  }
555 
556  // Next solve u'*x = b, overwriting b with x.
557  //
558  // k is the main loop index, increasing from 0 to n-1 in steps of
559  // 1 or 2, depending on the size of the diagonal blocks.
560 
561  k = 0;
562  while (k < n) {
563 
564  if (fIpiv[k] > 0) {
565  // 1 x 1 diagonal block
566  // multiply by inv(u'(k)), where u(k) is the transformation
567  // stored in column k of a.
568  for (Int_t i = 0; i < k; i++)
569  pcb[k*inc] -= pU[i*n+k]*pcb[i*inc];
570 
571  // interchange elements k and ipiv(k).
572  const Int_t kp = fIpiv[k]-1;
573  if (kp != k) {
574  const Double_t tmp = pcb[k*inc];
575  pcb[k*inc] = pcb[kp*inc];
576  pcb[kp*inc] = tmp;
577  }
578  k++;
579  } else {
580  // 2 x 2 diagonal block
581  // multiply by inv(u'(k+1)), where u(k+1) is the transformation
582  // stored in columns k and k+1 of a.
583  Int_t i;
584  for (i = 0; i < k; i++)
585  pcb[k*inc] -= pU[i*n+k]*pcb[i*inc];
586  for (i = 0; i < k; i++)
587  pcb[(k+1)*inc] -= pU[i*n+k+1]*pcb[i*inc];
588 
589  // interchange elements k and -ipiv(k).
590  const Int_t kp = -fIpiv[k]-1;
591  if (kp != k) {
592  const Double_t tmp = pcb[k*inc];
593  pcb[k*inc] = pcb[kp*inc];
594  pcb[kp*inc] = tmp;
595  }
596  k += 2;
597  }
598  }
599 
600  return kTRUE;
601 }
602 
603 ////////////////////////////////////////////////////////////////////////////////
604 /// For a symmetric matrix A(m,m), its inverse A_inv(m,m) is returned .
605 
607 {
608  if (inv.GetNrows() != GetNrows() || inv.GetRowLwb() != GetRowLwb()) {
609  Error("Invert(TMatrixDSym &","Input matrix has wrong shape");
610  return kFALSE;
611  }
612 
613  inv.UnitMatrix();
614 
615  const Int_t colLwb = inv.GetColLwb();
616  const Int_t colUpb = inv.GetColUpb();
617  Bool_t status = kTRUE;
618  for (Int_t icol = colLwb; icol <= colUpb && status; icol++) {
619  TMatrixDColumn b(inv,icol);
620  status &= Solve(b);
621  }
622 
623  return status;
624 }
625 
626 ////////////////////////////////////////////////////////////////////////////////
627 /// For a symmetric matrix A(m,m), its inverse A_inv(m,m) is returned .
628 
630 {
631  const Int_t rowLwb = GetRowLwb();
632  const Int_t rowUpb = rowLwb+GetNrows()-1;
633 
634  TMatrixDSym inv(rowLwb,rowUpb);
635  inv.UnitMatrix();
636  status = Invert(inv);
637 
638  return inv;
639 }
640 
641 ////////////////////////////////////////////////////////////////////////////////
642 /// Print the class members
643 
644 void TDecompBK::Print(Option_t *opt) const
645 {
646  TDecompBase::Print(opt);
647  printf("fIpiv:\n");
648  for (Int_t i = 0; i < fNIpiv; i++)
649  printf("[%d] = %d\n",i,fIpiv[i]);
650  fU.Print("fU");
651 }
652 
653 ////////////////////////////////////////////////////////////////////////////////
654 /// Assigment operator
655 
657 {
658  if (this != &source) {
659  TDecompBase::operator=(source);
660  fU.ResizeTo(source.fU);
661  fU = source.fU;
662  if (fNIpiv != source.fNIpiv) {
663  if (fIpiv)
664  delete [] fIpiv;
665  fNIpiv = source.fNIpiv;
666  fIpiv = new Int_t[fNIpiv];
667  }
668  if (fIpiv) memcpy(fIpiv,source.fIpiv,fNIpiv*sizeof(Int_t));
669  }
670  return *this;
671 }
Int_t GetRowLwb() const
Definition: TMatrixTBase.h:132
void Rank1Update(const TVectorT< Element > &vec, Element alpha=1.0)
Perform a rank 1 operation on the matrix: A += alpha * v * v^T.
virtual TMatrixTBase< Element > & UnitMatrix()
Make a unit matrix (matrix need not be a square one).
virtual TMatrixTBase< Element > & Shift(Int_t row_shift, Int_t col_shift)
Shift the row index by adding row_shift and the column index by adding col_shift, respectively...
Element * GetPtr() const
Long64_t LocMax(Long64_t n, const T *a)
Definition: TMath.h:724
Double_t fCondition
Definition: TDecompBase.h:42
Bool_t IsValid() const
Definition: TVectorT.h:89
Int_t fRowLwb
Definition: TDecompBase.h:43
const char Option_t
Definition: RtypesCore.h:62
Int_t fNIpiv
Definition: TDecompBK.h:35
virtual Bool_t Solve(TVectorD &b)
Solve Ax=b assuming the BK form of A is stored in fU . Solution returned in b.
Definition: TDecompBK.cxx:335
Int_t GetNrows() const
Definition: TVectorT.h:81
void Print(Option_t *opt="") const
Print class members.
#define R__ASSERT(e)
Definition: TError.h:98
double inv(double x)
For comparisons.
Definition: inv.h:58
TMatrixTColumn_const< Double_t > TMatrixDColumn_const
TAlienJobStatus * status
Definition: TAlienJob.cxx:51
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
const Bool_t kFALSE
Definition: Rtypes.h:92
ClassImp(TDecompBK) TDecompBK
Default constructor.
Definition: TDecompBK.cxx:16
virtual TMatrixTBase< Element > & ResizeTo(Int_t nrows, Int_t ncols, Int_t=-1)
Set size of the matrix to nrows x ncols New dynamic elements are created, the overlapping part of the...
Definition: TMatrixT.cxx:1201
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:732
Double_t fTol
Definition: TDecompBase.h:39
Element Norm1() const
Definition: TMatrixTBase.h:188
Bool_t IsValid() const
Definition: TMatrixTBase.h:157
Element GetTol() const
Definition: TMatrixTBase.h:139
TDecompBase & operator=(const TDecompBase &source)
Assignment operator.
TMatrixTRow_const< Double_t > TMatrixDRow_const
void Print(Option_t *opt="") const
Print the class members.
Definition: TDecompBK.cxx:644
Int_t fColLwb
Definition: TDecompBase.h:44
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:918
Int_t GetColUpb() const
Definition: TMatrixTBase.h:136
virtual const Element * GetMatrixArray() const
Definition: TMatrixTSym.h:183
Int_t GetNrows() const
Definition: TMatrixTBase.h:134
Element * GetMatrixArray()
Definition: TVectorT.h:84
void Print(Option_t *name="") const
Print the matrix as a table of elements.
const double tol
TDecompBK & operator=(const TDecompBK &source)
Assigment operator.
Definition: TDecompBK.cxx:656
TMatrixDSym Invert()
Definition: TDecompBK.h:69
unsigned int r1[N_CITIES]
Definition: simanTSP.cxx:321
virtual Int_t GetNrows() const
Definition: TDecompBK.h:50
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:173
virtual const Element * GetMatrixArray() const
Definition: TMatrixT.h:211
void ResetStatus()
Definition: TDecompBase.h:46
const TMatrixTBase< Element > * GetMatrix() const
Int_t GetColLwb() const
Definition: TMatrixTBase.h:135
Int_t GetRowLwb() const
Definition: TDecompBase.h:76
TMatrixD fU
Definition: TDecompBK.h:37
Int_t GetLwb() const
Definition: TVectorT.h:79
double Double_t
Definition: RtypesCore.h:55
virtual void SetMatrix(const TMatrixDSym &a)
Set the matrix to be decomposed, decomposition status is reset.
Definition: TDecompBK.cxx:310
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
Int_t * fIpiv
Definition: TDecompBK.h:36
Int_t GetNcols() const
Definition: TMatrixTBase.h:137
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
Int_t GetInc() const
TVectorT< Element > & Abs()
Take an absolute value of a vector, i.e. apply Abs() to each element.
Definition: TVectorT.cxx:461
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
const Bool_t kTRUE
Definition: Rtypes.h:91
virtual Bool_t Decompose()
Matrix A is decomposed in components U and D so that A = U*D*U^T If the decomposition succeeds...
Definition: TDecompBK.cxx:127
const Int_t n
Definition: legend1.C:16