30 ClassImp(TDecompChol);
35 TDecompChol::TDecompChol(Int_t nrows)
37 fU.ResizeTo(nrows,nrows);
43 TDecompChol::TDecompChol(Int_t row_lwb,Int_t row_upb)
45 const Int_t nrows = row_upb-row_lwb+1;
48 fU.ResizeTo(row_lwb,row_lwb+nrows-1,row_lwb,row_lwb+nrows-1);
54 TDecompChol::TDecompChol(
const TMatrixDSym &a,Double_t tol)
56 R__ASSERT(a.IsValid());
59 fCondition = a.Norm1();
64 fRowLwb = a.GetRowLwb();
65 fColLwb = a.GetColLwb();
73 TDecompChol::TDecompChol(
const TMatrixD &a,Double_t tol)
75 R__ASSERT(a.IsValid());
77 if (a.GetNrows() != a.GetNcols() || a.GetRowLwb() != a.GetColLwb()) {
78 Error(
"TDecompChol(const TMatrixD &",
"matrix should be square");
83 fCondition = a.Norm1();
88 fRowLwb = a.GetRowLwb();
89 fColLwb = a.GetColLwb();
97 TDecompChol::TDecompChol(
const TDecompChol &another) : TDecompBase(another)
106 Bool_t TDecompChol::Decompose()
108 if (TestBit(kDecomposed))
return kTRUE;
110 if ( !TestBit(kMatrixSet) ) {
111 Error(
"Decompose()",
"Matrix has not been set");
116 const Int_t n = fU.GetNrows();
117 Double_t *pU = fU.GetMatrixArray();
118 for (icol = 0; icol < n; icol++) {
119 const Int_t rowOff = icol*n;
122 Double_t ujj = pU[rowOff+icol];
123 for (irow = 0; irow < icol; irow++) {
124 const Int_t pos_ij = irow*n+icol;
125 ujj -= pU[pos_ij]*pU[pos_ij];
128 Error(
"Decompose()",
"matrix not positive definite");
131 ujj = TMath::Sqrt(ujj);
132 pU[rowOff+icol] = ujj;
135 for (j = icol+1; j < n; j++) {
136 for (i = 0; i < icol; i++) {
137 const Int_t rowOff2 = i*n;
138 pU[rowOff+j] -= pU[rowOff2+j]*pU[rowOff2+icol];
141 for (j = icol+1; j < n; j++)
146 for (irow = 0; irow < n; irow++) {
147 const Int_t rowOff = irow*n;
148 for (icol = 0; icol < irow; icol++)
149 pU[rowOff+icol] = 0.;
160 const TMatrixDSym TDecompChol::GetMatrix()
162 if (TestBit(kSingular)) {
163 Error(
"GetMatrix()",
"Matrix is singular");
164 return TMatrixDSym();
166 if ( !TestBit(kDecomposed) ) {
168 Error(
"GetMatrix()",
"Decomposition failed");
169 return TMatrixDSym();
173 return TMatrixDSym(TMatrixDSym::kAtA,fU);
179 void TDecompChol::SetMatrix(
const TMatrixDSym &a)
181 R__ASSERT(a.IsValid());
184 if (a.GetNrows() != a.GetNcols() || a.GetRowLwb() != a.GetColLwb()) {
185 Error(
"SetMatrix(const TMatrixDSym &",
"matrix should be square");
192 fRowLwb = a.GetRowLwb();
193 fColLwb = a.GetColLwb();
203 Bool_t TDecompChol::Solve(TVectorD &b)
205 R__ASSERT(b.IsValid());
206 if (TestBit(kSingular)) {
207 Error(
"Solve()",
"Matrix is singular");
210 if ( !TestBit(kDecomposed) ) {
212 Error(
"Solve()",
"Decomposition failed");
217 if (fU.GetNrows() != b.GetNrows() || fU.GetRowLwb() != b.GetLwb()) {
218 Error(
"Solve(TVectorD &",
"vector and matrix incompatible");
222 const Int_t n = fU.GetNrows();
224 const Double_t *pU = fU.GetMatrixArray();
225 Double_t *pb = b.GetMatrixArray();
229 for (i = 0; i < n; i++) {
230 const Int_t off_i = i*n;
231 if (pU[off_i+i] < fTol)
233 Error(
"Solve(TVectorD &b)",
"u[%d,%d]=%.4e < %.4e",i,i,pU[off_i+i],fTol);
237 for (Int_t j = 0; j < i; j++) {
238 const Int_t off_j = j*n;
239 r -= pU[off_j+i]*pb[j];
241 pb[i] = r/pU[off_i+i];
245 for (i = n-1; i >= 0; i--) {
246 const Int_t off_i = i*n;
248 for (Int_t j = i+1; j < n; j++)
249 r -= pU[off_i+j]*pb[j];
250 pb[i] = r/pU[off_i+i];
261 Bool_t TDecompChol::Solve(TMatrixDColumn &cb)
263 TMatrixDBase *b =
const_cast<TMatrixDBase *
>(cb.GetMatrix());
264 R__ASSERT(b->IsValid());
265 if (TestBit(kSingular)) {
266 Error(
"Solve()",
"Matrix is singular");
269 if ( !TestBit(kDecomposed) ) {
271 Error(
"Solve()",
"Decomposition failed");
276 if (fU.GetNrows() != b->GetNrows() || fU.GetRowLwb() != b->GetRowLwb())
278 Error(
"Solve(TMatrixDColumn &cb",
"vector and matrix incompatible");
282 const Int_t n = fU.GetNrows();
284 const Double_t *pU = fU.GetMatrixArray();
285 Double_t *pcb = cb.GetPtr();
286 const Int_t inc = cb.GetInc();
290 for (i = 0; i < n; i++) {
291 const Int_t off_i = i*n;
292 const Int_t off_i2 = i*inc;
293 if (pU[off_i+i] < fTol)
295 Error(
"Solve(TMatrixDColumn &cb)",
"u[%d,%d]=%.4e < %.4e",i,i,pU[off_i+i],fTol);
298 Double_t r = pcb[off_i2];
299 for (Int_t j = 0; j < i; j++) {
300 const Int_t off_j = j*n;
301 r -= pU[off_j+i]*pcb[j*inc];
303 pcb[off_i2] = r/pU[off_i+i];
307 for (i = n-1; i >= 0; i--) {
308 const Int_t off_i = i*n;
309 const Int_t off_i2 = i*inc;
310 Double_t r = pcb[off_i2];
311 for (Int_t j = i+1; j < n; j++)
312 r -= pU[off_i+j]*pcb[j*inc];
313 pcb[off_i2] = r/pU[off_i+i];
323 void TDecompChol::Det(Double_t &d1,Double_t &d2)
325 if ( !TestBit(kDetermined) ) {
326 if ( !TestBit(kDecomposed) )
328 TDecompBase::Det(d1,d2);
341 Bool_t TDecompChol::Invert(TMatrixDSym &inv)
343 if (inv.GetNrows() != GetNrows() || inv.GetRowLwb() != GetRowLwb()) {
344 Error(
"Invert(TMatrixDSym &",
"Input matrix has wrong shape");
350 const Int_t colLwb = inv.GetColLwb();
351 const Int_t colUpb = inv.GetColUpb();
352 Bool_t status = kTRUE;
353 for (Int_t icol = colLwb; icol <= colUpb && status; icol++) {
354 TMatrixDColumn b(inv,icol);
364 TMatrixDSym TDecompChol::Invert(Bool_t &status)
366 const Int_t rowLwb = GetRowLwb();
367 const Int_t rowUpb = rowLwb+GetNrows()-1;
369 TMatrixDSym inv(rowLwb,rowUpb);
371 status = Invert(inv);
379 void TDecompChol::Print(Option_t *opt)
const
381 TDecompBase::Print(opt);
388 TDecompChol &TDecompChol::operator=(
const TDecompChol &source)
390 if (
this != &source) {
391 TDecompBase::operator=(source);
392 fU.ResizeTo(source.fU);
404 TVectorD NormalEqn(
const TMatrixD &A,
const TVectorD &b)
406 TDecompChol ch(TMatrixDSym(TMatrixDSym::kAtA,A));
408 return ch.Solve(TMatrixD(TMatrixD::kTransposed,A)*b,ok);
419 TVectorD NormalEqn(
const TMatrixD &A,
const TVectorD &b,
const TVectorD &std)
421 if (!AreCompatible(b,std)) {
422 ::Error(
"NormalEqn",
"vectors b and std are not compatible");
428 for (Int_t irow = 0; irow < A.GetNrows(); irow++) {
429 TMatrixDRow(mAw,irow) *= 1/std(irow);
430 mBw(irow) /= std(irow);
432 TDecompChol ch(TMatrixDSym(TMatrixDSym::kAtA,mAw));
434 return ch.Solve(TMatrixD(TMatrixD::kTransposed,mAw)*mBw,ok);
444 TMatrixD NormalEqn(
const TMatrixD &A,
const TMatrixD &B)
446 TDecompChol ch(TMatrixDSym(TMatrixDSym::kAtA,A));
447 TMatrixD mX(A,TMatrixD::kTransposeMult,B);
463 TMatrixD NormalEqn(
const TMatrixD &A,
const TMatrixD &B,
const TVectorD &std)
467 for (Int_t irow = 0; irow < A.GetNrows(); irow++) {
468 TMatrixDRow(mAw,irow) *= 1/std(irow);
469 TMatrixDRow(mBw,irow) *= 1/std(irow);
472 TDecompChol ch(TMatrixDSym(TMatrixDSym::kAtA,mAw));
473 TMatrixD mX(mAw,TMatrixD::kTransposeMult,mBw);