Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TDecompQRH.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/** \class TDecompQRH
13 \ingroup Matrix
14
15 QR Decomposition class
16
17 Decompose a general (m x n) matrix A into A = fQ' fR H where
18
19~~~
20 fQ : (m x n) - internal Q' matrix (not orthoginal)
21 fR : (n x n) - upper triangular matrix
22 H : HouseHolder matrix which is stored through
23 fUp: (n) - vector with Householder up's
24 fW : (n) - vector with Householder beta's
25~~~
26
27 If row/column index of A starts at (rowLwb,colLwb) then
28 the decomposed matrices start from :
29~~~
30 fQ' : (rowLwb,0)
31 fR : (0,colLwb)
32 and the decomposed vectors start from :
33 fUp : (0)
34 fW : (0)
35~~~
36
37 In order to get the QR dcomposition of A (i.e. A = QR )
38 The orthoginal matrix Q needs to be computed from the internal Q' and
39 the up's and beta's vector defining the Householder transformation
40
41 The orthogonal Q matrix is returned to the user by calling the
42 function TDecompQRH::GetOrthogonalMatrix()
43
44 Errors arise from formation of reflectors i.e. singularity .
45 Note it attempts to handle the cases where the nRow <= nCol .
46*/
47
48#include "TDecompQRH.h"
49#include "TError.h" // For R__ASSERT
50
51////////////////////////////////////////////////////////////////////////////////
52/// Constructor for (nrows x ncols) matrix
53
55{
56 if (nrows < ncols) {
57 Error("TDecompQRH(Int_t,Int_t","matrix rows should be >= columns");
58 return;
59 }
60
63 if (nrows <= ncols) {
66 } else {
69 }
70}
71
72////////////////////////////////////////////////////////////////////////////////
73/// Constructor for ([row_lwb..row_upb] x [col_lwb..col_upb]) matrix
74
76{
77 const Int_t nrows = row_upb-row_lwb+1;
78 const Int_t ncols = col_upb-col_lwb+1;
79
80 if (nrows < ncols) {
81 Error("TDecompQRH(Int_t,Int_t,Int_t,Int_t","matrix rows should be >= columns");
82 return;
83 }
84
87
90 if (nrows <= ncols) {
93 } else {
96 }
97}
98
99////////////////////////////////////////////////////////////////////////////////
100/// Constructor for general matrix A .
101
103{
104 R__ASSERT(a.IsValid());
105 if (a.GetNrows() < a.GetNcols()) {
106 Error("TDecompQRH(const TMatrixD &","matrix rows should be >= columns");
107 return;
108 }
109
111 fCondition = a.Norm1();
112 fTol = a.GetTol();
113 if (tol > 0.0)
114 fTol = tol;
115
116 fRowLwb = a.GetRowLwb();
117 fColLwb = a.GetColLwb();
118 const Int_t nRow = a.GetNrows();
119 const Int_t nCol = a.GetNcols();
120
122 memcpy(fQ.GetMatrixArray(),a.GetMatrixArray(),nRow*nCol*sizeof(Double_t));
124 if (nRow <= nCol) {
127 } else {
130 }
131}
132
133////////////////////////////////////////////////////////////////////////////////
134/// Copy constructor
135
140
141////////////////////////////////////////////////////////////////////////////////
142/// QR decomposition of matrix a by Householder transformations,
143/// see Golub & Loan first edition p41 & Sec 6.2.
144/// First fR is returned in upper triang of fQ and diagR. fQ returned in
145/// 'u-form' in lower triang of fQ and fW, the latter containing the
146/// "Householder betas".
147/// If the decomposition succeeds, bit kDecomposed is set , otherwise kSingular
148
150{
151 if (TestBit(kDecomposed)) return kTRUE;
152
153 if ( !TestBit(kMatrixSet) ) {
154 Error("Decompose()","Matrix has not been set");
155 return kFALSE;
156 }
157
158 const Int_t nRow = this->GetNrows();
159 const Int_t nCol = this->GetNcols();
160 const Int_t rowLwb = this->GetRowLwb();
161 const Int_t colLwb = this->GetColLwb();
162
165 if (nCol > kWorkMax) diagR.ResizeTo(nCol);
166 else diagR.Use(nCol,work);
167
168 if (QRH(fQ,diagR,fUp,fW,fTol)) {
169 for (Int_t i = 0; i < nRow; i++) {
170 const Int_t ic = (i < nCol) ? i : nCol;
171 for (Int_t j = ic ; j < nCol; j++)
172 fR(i,j) = fQ(i,j);
173 }
175
176 fQ.Shift(rowLwb,0);
177 fR.Shift(0,colLwb);
178
180 }
181
182 return kTRUE;
183}
184
185////////////////////////////////////////////////////////////////////////////////
186/// Decomposition function .
187
189{
190 const Int_t nRow = q.GetNrows();
191 const Int_t nCol = q.GetNcols();
192
193 const Int_t n = (nRow <= nCol) ? nRow-1 : nCol;
194
195 for (Int_t k = 0 ; k < n ; k++) {
197 if (!DefHouseHolder(qc_k,k,k+1,up(k),w(k),tol))
198 return kFALSE;
199 diagR(k) = qc_k(k)-up(k);
200 if (k < nCol-1) {
201 // Apply HouseHolder to sub-matrix
202 for (Int_t j = k+1; j < nCol; j++) {
204 ApplyHouseHolder(qc_k,up(k),w(k),k,k+1,qc_j);
205 }
206 }
207 }
208
209 if (nRow <= nCol) {
210 diagR(nRow-1) = q(nRow-1,nRow-1);
211 up(nRow-1) = 0;
212 w(nRow-1) = 0;
213 }
214
215 return kTRUE;
216}
217
218////////////////////////////////////////////////////////////////////////////////
219/// Set matrix to be decomposed
220
222{
223 R__ASSERT(a.IsValid());
224
225 ResetStatus();
226 if (a.GetNrows() < a.GetNcols()) {
227 Error("TDecompQRH(const TMatrixD &","matrix rows should be >= columns");
228 return;
229 }
230
232 fCondition = a.Norm1();
233
234 fRowLwb = a.GetRowLwb();
235 fColLwb = a.GetColLwb();
236 const Int_t nRow = a.GetNrows();
237 const Int_t nCol = a.GetNcols();
238
240 memcpy(fQ.GetMatrixArray(),a.GetMatrixArray(),nRow*nCol*sizeof(Double_t));
242 if (nRow <= nCol) {
245 } else {
248 }
249}
250
251////////////////////////////////////////////////////////////////////////////////
252/// Solve Ax=b assuming the QR form of A is stored in fR,fQ and fW, but assume b
253/// has *not* been transformed. Solution returned in b.
254
256{
257 R__ASSERT(b.IsValid());
258 if (TestBit(kSingular)) {
259 Error("Solve()","Matrix is singular");
260 return kFALSE;
261 }
262 if ( !TestBit(kDecomposed) ) {
263 if (!Decompose()) {
264 Error("Solve()","Decomposition failed");
265 return kFALSE;
266 }
267 }
268
269 if (fQ.GetNrows() != b.GetNrows() || fQ.GetRowLwb() != b.GetLwb()) {
270 Error("Solve(TVectorD &","vector and matrix incompatible");
271 return kFALSE;
272 }
273
274 const Int_t nQRow = fQ.GetNrows();
275 const Int_t nQCol = fQ.GetNcols();
276
277 // Calculate Q^T.b
278 const Int_t nQ = (nQRow <= nQCol) ? nQRow-1 : nQCol;
279 for (Int_t k = 0; k < nQ; k++) {
281 ApplyHouseHolder(qc_k,fUp(k),fW(k),k,k+1,b);
282 }
283
284 const Int_t nRCol = fR.GetNcols();
285
286 const Double_t *pR = fR.GetMatrixArray();
287 Double_t *pb = b.GetMatrixArray();
288
289 // Backward substitution
290 for (Int_t i = nRCol-1; i >= 0; i--) {
291 const Int_t off_i = i*nRCol;
292 Double_t r = pb[i];
293 for (Int_t j = i+1; j < nRCol; j++)
294 r -= pR[off_i+j]*pb[j];
295 if (TMath::Abs(pR[off_i+i]) < fTol)
296 {
297 Error("Solve(TVectorD &)","R[%d,%d]=%.4e < %.4e",i,i,pR[off_i+i],fTol);
298 return kFALSE;
299 }
300 pb[i] = r/pR[off_i+i];
301 }
302
303 return kTRUE;
304}
305
306////////////////////////////////////////////////////////////////////////////////
307/// Solve Ax=b assuming the QR form of A is stored in fR,fQ and fW, but assume b
308/// has *not* been transformed. Solution returned in b.
309
311{
312 TMatrixDBase *b = const_cast<TMatrixDBase *>(cb.GetMatrix());
313 R__ASSERT(b->IsValid());
314 if (TestBit(kSingular)) {
315 Error("Solve()","Matrix is singular");
316 return kFALSE;
317 }
318 if ( !TestBit(kDecomposed) ) {
319 if (!Decompose()) {
320 Error("Solve()","Decomposition failed");
321 return kFALSE;
322 }
323 }
324
325 if (fQ.GetNrows() != b->GetNrows() || fQ.GetRowLwb() != b->GetRowLwb())
326 {
327 Error("Solve(TMatrixDColumn &","vector and matrix incompatible");
328 return kFALSE;
329 }
330
331 const Int_t nQRow = fQ.GetNrows();
332 const Int_t nQCol = fQ.GetNcols();
333
334 // Calculate Q^T.b
335 const Int_t nQ = (nQRow <= nQCol) ? nQRow-1 : nQCol;
336 for (Int_t k = 0; k < nQ; k++) {
338 ApplyHouseHolder(qc_k,fUp(k),fW(k),k,k+1,cb);
339 }
340
341 const Int_t nRCol = fR.GetNcols();
342
343 const Double_t *pR = fR.GetMatrixArray();
344 Double_t *pcb = cb.GetPtr();
345 const Int_t inc = cb.GetInc();
346
347 // Backward substitution
348 for (Int_t i = nRCol-1; i >= 0; i--) {
349 const Int_t off_i = i*nRCol;
350 const Int_t off_i2 = i*inc;
351 Double_t r = pcb[off_i2];
352 for (Int_t j = i+1; j < nRCol; j++)
353 r -= pR[off_i+j]*pcb[j*inc];
354 if (TMath::Abs(pR[off_i+i]) < fTol)
355 {
356 Error("Solve(TMatrixDColumn &)","R[%d,%d]=%.4e < %.4e",i,i,pR[off_i+i],fTol);
357 return kFALSE;
358 }
359 pcb[off_i2] = r/pR[off_i+i];
360 }
361
362 return kTRUE;
363}
364
365////////////////////////////////////////////////////////////////////////////////
366/// Solve A^T x=b assuming the QR form of A is stored in fR,fQ and fW, but assume b
367/// has *not* been transformed. Solution returned in b.
368
370{
371 R__ASSERT(b.IsValid());
372 if (TestBit(kSingular)) {
373 Error("TransSolve()","Matrix is singular");
374 return kFALSE;
375 }
376 if ( !TestBit(kDecomposed) ) {
377 if (!Decompose()) {
378 Error("TransSolve()","Decomposition failed");
379 return kFALSE;
380 }
381 }
382
383 if (fQ.GetNrows() != fQ.GetNcols() || fQ.GetRowLwb() != fQ.GetColLwb()) {
384 Error("TransSolve(TVectorD &","matrix should be square");
385 return kFALSE;
386 }
387
388 if (fR.GetNrows() != b.GetNrows() || fR.GetRowLwb() != b.GetLwb()) {
389 Error("TransSolve(TVectorD &","vector and matrix incompatible");
390 return kFALSE;
391 }
392
393 const Double_t *pR = fR.GetMatrixArray();
394 Double_t *pb = b.GetMatrixArray();
395
396 const Int_t nRCol = fR.GetNcols();
397
398 // Backward substitution
399 for (Int_t i = 0; i < nRCol; i++) {
400 const Int_t off_i = i*nRCol;
401 Double_t r = pb[i];
402 for (Int_t j = 0; j < i; j++) {
403 const Int_t off_j = j*nRCol;
404 r -= pR[off_j+i]*pb[j];
405 }
406 if (TMath::Abs(pR[off_i+i]) < fTol)
407 {
408 Error("TransSolve(TVectorD &)","R[%d,%d]=%.4e < %.4e",i,i,pR[off_i+i],fTol);
409 return kFALSE;
410 }
411 pb[i] = r/pR[off_i+i];
412 }
413
414 const Int_t nQRow = fQ.GetNrows();
415
416 // Calculate Q.b; it was checked nQRow == nQCol
417 for (Int_t k = nQRow-1; k >= 0; k--) {
419 ApplyHouseHolder(qc_k,fUp(k),fW(k),k,k+1,b);
420 }
421
422 return kTRUE;
423}
424
425////////////////////////////////////////////////////////////////////////////////
426/// Solve A^T x=b assuming the QR form of A is stored in fR,fQ and fW, but assume b
427/// has *not* been transformed. Solution returned in b.
428
430{
431 TMatrixDBase *b = const_cast<TMatrixDBase *>(cb.GetMatrix());
432 R__ASSERT(b->IsValid());
433 if (TestBit(kSingular)) {
434 Error("TransSolve()","Matrix is singular");
435 return kFALSE;
436 }
437 if ( !TestBit(kDecomposed) ) {
438 if (!Decompose()) {
439 Error("TransSolve()","Decomposition failed");
440 return kFALSE;
441 }
442 }
443
444 if (fQ.GetNrows() != fQ.GetNcols() || fQ.GetRowLwb() != fQ.GetColLwb()) {
445 Error("TransSolve(TMatrixDColumn &","matrix should be square");
446 return kFALSE;
447 }
448
449 if (fR.GetNrows() != b->GetNrows() || fR.GetRowLwb() != b->GetRowLwb()) {
450 Error("TransSolve(TMatrixDColumn &","vector and matrix incompatible");
451 return kFALSE;
452 }
453
454 const Double_t *pR = fR.GetMatrixArray();
455 Double_t *pcb = cb.GetPtr();
456 const Int_t inc = cb.GetInc();
457
458 const Int_t nRCol = fR.GetNcols();
459
460 // Backward substitution
461 for (Int_t i = 0; i < nRCol; i++) {
462 const Int_t off_i = i*nRCol;
463 const Int_t off_i2 = i*inc;
464 Double_t r = pcb[off_i2];
465 for (Int_t j = 0; j < i; j++) {
466 const Int_t off_j = j*nRCol;
467 r -= pR[off_j+i]*pcb[j*inc];
468 }
469 if (TMath::Abs(pR[off_i+i]) < fTol)
470 {
471 Error("TransSolve(TMatrixDColumn &)","R[%d,%d]=%.4e < %.4e",i,i,pR[off_i+i],fTol);
472 return kFALSE;
473 }
474 pcb[off_i2] = r/pR[off_i+i];
475 }
476
477 const Int_t nQRow = fQ.GetNrows();
478
479 // Calculate Q.b; it was checked nQRow == nQCol
480 for (Int_t k = nQRow-1; k >= 0; k--) {
482 ApplyHouseHolder(qc_k,fUp(k),fW(k),k,k+1,cb);
483 }
484
485 return kTRUE;
486}
487
488////////////////////////////////////////////////////////////////////////////////
489/// This routine calculates the absolute (!) value of the determinant
490/// |det| = d1*TMath::Power(2.,d2)
491
493{
494 if ( !TestBit(kDetermined) ) {
495 if ( !TestBit(kDecomposed) )
496 Decompose();
497 if (TestBit(kSingular)) {
498 fDet1 = 0.0;
499 fDet2 = 0.0;
500 } else
503 }
504 d1 = fDet1;
505 d2 = fDet2;
506}
507
508////////////////////////////////////////////////////////////////////////////////
509/// For a matrix A(m,n), return the OtrhogonalMatrix Q such as
510/// A = Q * R
511///
512/// Note that this Q is not th einternal fQ matrix obtained in the QRH decomposition, but can be computed
513/// from the fQ and the up and beta vector's defining the Householder transformation
514
516{
517 // apply HouseHolder transformation starting from the identity
518 // Calculate Q.b; it was checked nQRow == nQCol
519
520 const Int_t nRow = this->GetNrows();
521 const Int_t nCol = this->GetNcols();
522 // remmber nCol <= nRow
524 // start from identity matrix
525 for (int i = 0; i < nCol; ++i)
526 orthogQ(i, i) = 1;
527
528
529 // apply the HouseHolder transformations for each column of Q
530 for (int j = 0; j < nCol; ++j) {
532 int nQRow = fQ.GetNrows();
533 for (Int_t k = nQRow - 1; k >= 0; k--) {
535 ApplyHouseHolder(qc_k, fUp(k), fW(k), k, k + 1, b);
536 }
537 }
538 return orthogQ;
539}
540////////////////////////////////////////////////////////////////////////////////
541/// For a matrix A(m,n), its inverse A_inv is defined as A * A_inv = A_inv * A = unit
542/// The user should always supply a matrix of size (m x m) !
543/// If m > n , only the (n x m) part of the returned (pseudo inverse) matrix
544/// should be used .
545
547{
548 if (inv.GetNrows() != GetNrows() || inv.GetNcols() != GetNrows() ||
549 inv.GetRowLwb() != GetRowLwb() || inv.GetColLwb() != GetColLwb()) {
550 Error("Invert(TMatrixD &","Input matrix has wrong shape");
551 return kFALSE;
552 }
553
554 inv.UnitMatrix();
555 const Bool_t status = MultiSolve(inv);
556
557 return status;
558}
559
560////////////////////////////////////////////////////////////////////////////////
561/// For a matrix A(m,n), its inverse A_inv is defined as A * A_inv = A_inv * A = unit
562/// (n x m) Ainv is returned .
563
565{
566 const Int_t rowLwb = GetRowLwb();
567 const Int_t colLwb = GetColLwb();
568 const Int_t rowUpb = rowLwb+GetNrows()-1;
570 inv.UnitMatrix();
571 status = MultiSolve(inv);
572 inv.ResizeTo(rowLwb,rowLwb+GetNcols()-1,colLwb,colLwb+GetNrows()-1);
573
574 return inv;
575}
576
577////////////////////////////////////////////////////////////////////////////////
578/// Print the class members
579
581{
583 fQ.Print("fQ");
584 fR.Print("fR");
585 fUp.Print("fUp");
586 fW.Print("fW");
587}
588
589////////////////////////////////////////////////////////////////////////////////
590/// Assignment operator
591
593{
594 if (this != &source) {
596 fQ.ResizeTo(source.fQ);
597 fR.ResizeTo(source.fR);
598 fUp.ResizeTo(source.fUp);
599 fW.ResizeTo(source.fW);
600 fQ = source.fQ;
601 fR = source.fR;
602 fUp = source.fUp;
603 fW = source.fW;
604 }
605 return *this;
606}
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Bool_t DefHouseHolder(const TVectorD &vc, Int_t lp, Int_t l, Double_t &up, Double_t &b, Double_t tol=0.0)
Define a Householder-transformation through the parameters up and b .
void ApplyHouseHolder(const TVectorD &vc, Double_t up, Double_t b, Int_t lp, Int_t l, TMatrixDRow &cr)
Apply Householder-transformation.
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
float * q
TMatrixTColumn_const< Double_t > TMatrixDColumn_const
TMatrixTColumn< Double_t > TMatrixDColumn
Decomposition Base class.
Definition TDecompBase.h:34
Int_t GetRowLwb() const
Definition TDecompBase.h:73
Double_t fDet1
Definition TDecompBase.h:37
Double_t fDet2
Definition TDecompBase.h:38
void ResetStatus()
Definition TDecompBase.h:43
Int_t GetColLwb() const
Definition TDecompBase.h:74
virtual Bool_t MultiSolve(TMatrixD &B)
Solve set of equations with RHS in columns of B.
Int_t fRowLwb
Definition TDecompBase.h:40
Double_t fTol
Definition TDecompBase.h:36
TDecompBase & operator=(const TDecompBase &source)
Assignment operator.
Double_t fCondition
Definition TDecompBase.h:39
Int_t fColLwb
Definition TDecompBase.h:41
void Print(Option_t *opt="") const override
Print class members.
virtual void Det(Double_t &d1, Double_t &d2)
Matrix determinant det = d1*TMath::Power(2.,d2)
QR Decomposition class.
Definition TDecompQRH.h:26
virtual void SetMatrix(const TMatrixD &a)
Set matrix to be decomposed.
TMatrixD Invert()
Definition TDecompQRH.h:77
Bool_t Solve(TVectorD &b) override
Solve Ax=b assuming the QR form of A is stored in fR,fQ and fW, but assume b has not been transformed...
Int_t GetNrows() const override
Definition TDecompQRH.h:50
void Print(Option_t *opt="") const override
Print the class members.
void Det(Double_t &d1, Double_t &d2) override
This routine calculates the absolute (!) value of the determinant |det| = d1*TMath::Power(2....
TDecompQRH & operator=(const TDecompQRH &source)
Assignment operator.
TMatrixD GetOrthogonalMatrix() const
For a matrix A(m,n), return the OtrhogonalMatrix Q such as A = Q * R.
TMatrixD fQ
Definition TDecompQRH.h:30
TMatrixD fR
Definition TDecompQRH.h:31
static Bool_t QRH(TMatrixD &q, TVectorD &diagR, TVectorD &up, TVectorD &w, Double_t tol)
Decomposition function .
Bool_t Decompose() override
QR decomposition of matrix a by Householder transformations, see Golub & Loan first edition p41 & Sec...
TVectorD fW
Definition TDecompQRH.h:33
Int_t GetNcols() const override
Definition TDecompQRH.h:51
TVectorD fUp
Definition TDecompQRH.h:32
Bool_t TransSolve(TVectorD &b) override
Solve A^T x=b assuming the QR form of A is stored in fR,fQ and fW, but assume b has not been transfor...
Int_t GetNrows() const
void Print(Option_t *name="") const override
Print the matrix as a table of elements.
Int_t GetRowLwb() const
Int_t GetColLwb() const
Int_t GetNcols() const
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.
const TMatrixTBase< Element > * GetMatrix() const
Element * GetPtr() const
const Element * GetMatrixArray() const override
Definition TMatrixT.h:228
TMatrixTBase< Element > & ResizeTo(Int_t nrows, Int_t ncols, Int_t=-1) override
Set size of the matrix to nrows x ncols New dynamic elements are created, the overlapping part of the...
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:202
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:864
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1071
TVectorT< Element > & ResizeTo(Int_t lwb, Int_t upb)
Resize the vector to [lwb:upb] .
Definition TVectorT.cxx:293
void Print(Option_t *option="") const override
Print the vector as a list of elements.
const Int_t n
Definition legend1.C:16
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:124
void inv(rsa_NUMBER *, rsa_NUMBER *, rsa_NUMBER *)
Definition rsaaux.cxx:949