00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027 #include "TDecompChol.h"
00028 #include "TMath.h"
00029
00030 ClassImp(TDecompChol)
00031
00032
00033 TDecompChol::TDecompChol(Int_t nrows)
00034 {
00035
00036
00037 fU.ResizeTo(nrows,nrows);
00038 }
00039
00040
00041 TDecompChol::TDecompChol(Int_t row_lwb,Int_t row_upb)
00042 {
00043
00044
00045 const Int_t nrows = row_upb-row_lwb+1;
00046 fRowLwb = row_lwb;
00047 fColLwb = row_lwb;
00048 fU.ResizeTo(row_lwb,row_lwb+nrows-1,row_lwb,row_lwb+nrows-1);
00049 }
00050
00051
00052 TDecompChol::TDecompChol(const TMatrixDSym &a,Double_t tol)
00053 {
00054
00055
00056 R__ASSERT(a.IsValid());
00057
00058 SetBit(kMatrixSet);
00059 fCondition = a.Norm1();
00060 fTol = a.GetTol();
00061 if (tol > 0)
00062 fTol = tol;
00063
00064 fRowLwb = a.GetRowLwb();
00065 fColLwb = a.GetColLwb();
00066 fU.ResizeTo(a);
00067 fU = a;
00068 }
00069
00070
00071 TDecompChol::TDecompChol(const TMatrixD &a,Double_t tol)
00072 {
00073
00074
00075 R__ASSERT(a.IsValid());
00076
00077 if (a.GetNrows() != a.GetNcols() || a.GetRowLwb() != a.GetColLwb()) {
00078 Error("TDecompChol(const TMatrixD &","matrix should be square");
00079 return;
00080 }
00081
00082 SetBit(kMatrixSet);
00083 fCondition = a.Norm1();
00084 fTol = a.GetTol();
00085 if (tol > 0)
00086 fTol = tol;
00087
00088 fRowLwb = a.GetRowLwb();
00089 fColLwb = a.GetColLwb();
00090 fU.ResizeTo(a);
00091 fU = a;
00092 }
00093
00094
00095 TDecompChol::TDecompChol(const TDecompChol &another) : TDecompBase(another)
00096 {
00097
00098
00099 *this = another;
00100 }
00101
00102
00103 Bool_t TDecompChol::Decompose()
00104 {
00105
00106
00107
00108 if (TestBit(kDecomposed)) return kTRUE;
00109
00110 if ( !TestBit(kMatrixSet) ) {
00111 Error("Decompose()","Matrix has not been set");
00112 return kFALSE;
00113 }
00114
00115 Int_t i,j,icol,irow;
00116 const Int_t n = fU.GetNrows();
00117 Double_t *pU = fU.GetMatrixArray();
00118 for (icol = 0; icol < n; icol++) {
00119 const Int_t rowOff = icol*n;
00120
00121
00122 Double_t ujj = pU[rowOff+icol];
00123 for (irow = 0; irow < icol; irow++) {
00124 const Int_t pos_ij = irow*n+icol;
00125 ujj -= pU[pos_ij]*pU[pos_ij];
00126 }
00127 if (ujj <= 0) {
00128 Error("Decompose()","matrix not positive definite");
00129 return kFALSE;
00130 }
00131 ujj = TMath::Sqrt(ujj);
00132 pU[rowOff+icol] = ujj;
00133
00134 if (icol < n-1) {
00135 for (j = icol+1; j < n; j++) {
00136 for (i = 0; i < icol; i++) {
00137 const Int_t rowOff2 = i*n;
00138 pU[rowOff+j] -= pU[rowOff2+j]*pU[rowOff2+icol];
00139 }
00140 }
00141 for (j = icol+1; j < n; j++)
00142 pU[rowOff+j] /= ujj;
00143 }
00144 }
00145
00146 for (irow = 0; irow < n; irow++) {
00147 const Int_t rowOff = irow*n;
00148 for (icol = 0; icol < irow; icol++)
00149 pU[rowOff+icol] = 0.;
00150 }
00151
00152 SetBit(kDecomposed);
00153
00154 return kTRUE;
00155 }
00156
00157
00158 const TMatrixDSym TDecompChol::GetMatrix()
00159 {
00160
00161
00162 if (TestBit(kSingular)) {
00163 Error("GetMatrix()","Matrix is singular");
00164 return TMatrixDSym();
00165 }
00166 if ( !TestBit(kDecomposed) ) {
00167 if (!Decompose()) {
00168 Error("GetMatrix()","Decomposition failed");
00169 return TMatrixDSym();
00170 }
00171 }
00172
00173 return TMatrixDSym(TMatrixDSym::kAtA,fU);
00174 }
00175
00176
00177 void TDecompChol::SetMatrix(const TMatrixDSym &a)
00178 {
00179
00180
00181 R__ASSERT(a.IsValid());
00182
00183 ResetStatus();
00184 if (a.GetNrows() != a.GetNcols() || a.GetRowLwb() != a.GetColLwb()) {
00185 Error("SetMatrix(const TMatrixDSym &","matrix should be square");
00186 return;
00187 }
00188
00189 SetBit(kMatrixSet);
00190 fCondition = -1.0;
00191
00192 fRowLwb = a.GetRowLwb();
00193 fColLwb = a.GetColLwb();
00194 fU.ResizeTo(a);
00195 fU = a;
00196 }
00197
00198
00199 Bool_t TDecompChol::Solve(TVectorD &b)
00200 {
00201
00202
00203
00204
00205 R__ASSERT(b.IsValid());
00206 if (TestBit(kSingular)) {
00207 Error("Solve()","Matrix is singular");
00208 return kFALSE;
00209 }
00210 if ( !TestBit(kDecomposed) ) {
00211 if (!Decompose()) {
00212 Error("Solve()","Decomposition failed");
00213 return kFALSE;
00214 }
00215 }
00216
00217 if (fU.GetNrows() != b.GetNrows() || fU.GetRowLwb() != b.GetLwb()) {
00218 Error("Solve(TVectorD &","vector and matrix incompatible");
00219 return kFALSE;
00220 }
00221
00222 const Int_t n = fU.GetNrows();
00223
00224 const Double_t *pU = fU.GetMatrixArray();
00225 Double_t *pb = b.GetMatrixArray();
00226
00227 Int_t i;
00228
00229 for (i = 0; i < n; i++) {
00230 const Int_t off_i = i*n;
00231 if (pU[off_i+i] < fTol)
00232 {
00233 Error("Solve(TVectorD &b)","u[%d,%d]=%.4e < %.4e",i,i,pU[off_i+i],fTol);
00234 return kFALSE;
00235 }
00236 Double_t r = pb[i];
00237 for (Int_t j = 0; j < i; j++) {
00238 const Int_t off_j = j*n;
00239 r -= pU[off_j+i]*pb[j];
00240 }
00241 pb[i] = r/pU[off_i+i];
00242 }
00243
00244
00245 for (i = n-1; i >= 0; i--) {
00246 const Int_t off_i = i*n;
00247 Double_t r = pb[i];
00248 for (Int_t j = i+1; j < n; j++)
00249 r -= pU[off_i+j]*pb[j];
00250 pb[i] = r/pU[off_i+i];
00251 }
00252
00253 return kTRUE;
00254 }
00255
00256
00257 Bool_t TDecompChol::Solve(TMatrixDColumn &cb)
00258 {
00259
00260
00261
00262
00263 TMatrixDBase *b = const_cast<TMatrixDBase *>(cb.GetMatrix());
00264 R__ASSERT(b->IsValid());
00265 if (TestBit(kSingular)) {
00266 Error("Solve()","Matrix is singular");
00267 return kFALSE;
00268 }
00269 if ( !TestBit(kDecomposed) ) {
00270 if (!Decompose()) {
00271 Error("Solve()","Decomposition failed");
00272 return kFALSE;
00273 }
00274 }
00275
00276 if (fU.GetNrows() != b->GetNrows() || fU.GetRowLwb() != b->GetRowLwb())
00277 {
00278 Error("Solve(TMatrixDColumn &cb","vector and matrix incompatible");
00279 return kFALSE;
00280 }
00281
00282 const Int_t n = fU.GetNrows();
00283
00284 const Double_t *pU = fU.GetMatrixArray();
00285 Double_t *pcb = cb.GetPtr();
00286 const Int_t inc = cb.GetInc();
00287
00288 Int_t i;
00289
00290 for (i = 0; i < n; i++) {
00291 const Int_t off_i = i*n;
00292 const Int_t off_i2 = i*inc;
00293 if (pU[off_i+i] < fTol)
00294 {
00295 Error("Solve(TMatrixDColumn &cb)","u[%d,%d]=%.4e < %.4e",i,i,pU[off_i+i],fTol);
00296 return kFALSE;
00297 }
00298 Double_t r = pcb[off_i2];
00299 for (Int_t j = 0; j < i; j++) {
00300 const Int_t off_j = j*n;
00301 r -= pU[off_j+i]*pcb[j*inc];
00302 }
00303 pcb[off_i2] = r/pU[off_i+i];
00304 }
00305
00306
00307 for (i = n-1; i >= 0; i--) {
00308 const Int_t off_i = i*n;
00309 const Int_t off_i2 = i*inc;
00310 Double_t r = pcb[off_i2];
00311 for (Int_t j = i+1; j < n; j++)
00312 r -= pU[off_i+j]*pcb[j*inc];
00313 pcb[off_i2] = r/pU[off_i+i];
00314 }
00315
00316 return kTRUE;
00317 }
00318
00319
00320 void TDecompChol::Det(Double_t &d1,Double_t &d2)
00321 {
00322
00323
00324
00325 if ( !TestBit(kDetermined) ) {
00326 if ( !TestBit(kDecomposed) )
00327 Decompose();
00328 TDecompBase::Det(d1,d2);
00329
00330 fDet1 *= fDet1;
00331 fDet2 += fDet2;
00332 SetBit(kDetermined);
00333 }
00334 d1 = fDet1;
00335 d2 = fDet2;
00336 }
00337
00338
00339 Bool_t TDecompChol::Invert(TMatrixDSym &inv)
00340 {
00341
00342
00343 if (inv.GetNrows() != GetNrows() || inv.GetRowLwb() != GetRowLwb()) {
00344 Error("Invert(TMatrixDSym &","Input matrix has wrong shape");
00345 return kFALSE;
00346 }
00347
00348 inv.UnitMatrix();
00349
00350 const Int_t colLwb = inv.GetColLwb();
00351 const Int_t colUpb = inv.GetColUpb();
00352 Bool_t status = kTRUE;
00353 for (Int_t icol = colLwb; icol <= colUpb && status; icol++) {
00354 TMatrixDColumn b(inv,icol);
00355 status &= Solve(b);
00356 }
00357
00358 return status;
00359 }
00360
00361
00362 TMatrixDSym TDecompChol::Invert(Bool_t &status)
00363 {
00364
00365
00366 const Int_t rowLwb = GetRowLwb();
00367 const Int_t rowUpb = rowLwb+GetNrows()-1;
00368
00369 TMatrixDSym inv(rowLwb,rowUpb);
00370 inv.UnitMatrix();
00371 status = Invert(inv);
00372
00373 return inv;
00374 }
00375
00376
00377 void TDecompChol::Print(Option_t *opt) const
00378 {
00379
00380
00381 TDecompBase::Print(opt);
00382 fU.Print("fU");
00383 }
00384
00385
00386 TDecompChol &TDecompChol::operator=(const TDecompChol &source)
00387 {
00388
00389
00390 if (this != &source) {
00391 TDecompBase::operator=(source);
00392 fU.ResizeTo(source.fU);
00393 fU = source.fU;
00394 }
00395 return *this;
00396 }
00397
00398
00399 TVectorD NormalEqn(const TMatrixD &A,const TVectorD &b)
00400 {
00401
00402
00403
00404
00405
00406 TDecompChol ch(TMatrixDSym(TMatrixDSym::kAtA,A));
00407 Bool_t ok;
00408 return ch.Solve(TMatrixD(TMatrixD::kTransposed,A)*b,ok);
00409 }
00410
00411
00412 TVectorD NormalEqn(const TMatrixD &A,const TVectorD &b,const TVectorD &std)
00413 {
00414
00415
00416
00417
00418
00419
00420
00421 if (!AreCompatible(b,std)) {
00422 ::Error("NormalEqn","vectors b and std are not compatible");
00423 return TVectorD();
00424 }
00425
00426 TMatrixD mAw = A;
00427 TVectorD mBw = b;
00428 for (Int_t irow = 0; irow < A.GetNrows(); irow++) {
00429 TMatrixDRow(mAw,irow) *= 1/std(irow);
00430 mBw(irow) /= std(irow);
00431 }
00432 TDecompChol ch(TMatrixDSym(TMatrixDSym::kAtA,mAw));
00433 Bool_t ok;
00434 return ch.Solve(TMatrixD(TMatrixD::kTransposed,mAw)*mBw,ok);
00435 }
00436
00437
00438 TMatrixD NormalEqn(const TMatrixD &A,const TMatrixD &B)
00439 {
00440
00441
00442
00443
00444
00445
00446 TDecompChol ch(TMatrixDSym(TMatrixDSym::kAtA,A));
00447 TMatrixD mX(A,TMatrixD::kTransposeMult,B);
00448 ch.MultiSolve(mX);
00449 return mX;
00450 }
00451
00452
00453 TMatrixD NormalEqn(const TMatrixD &A,const TMatrixD &B,const TVectorD &std)
00454 {
00455
00456
00457
00458
00459
00460
00461
00462
00463 TMatrixD mAw = A;
00464 TMatrixD mBw = B;
00465 for (Int_t irow = 0; irow < A.GetNrows(); irow++) {
00466 TMatrixDRow(mAw,irow) *= 1/std(irow);
00467 TMatrixDRow(mBw,irow) *= 1/std(irow);
00468 }
00469
00470 TDecompChol ch(TMatrixDSym(TMatrixDSym::kAtA,mAw));
00471 TMatrixD mX(mAw,TMatrixD::kTransposeMult,mBw);
00472 ch.MultiSolve(mX);
00473 return mX;
00474 }