64 TDecompBK::TDecompBK()
73 TDecompBK::TDecompBK(Int_t nrows)
76 fIpiv =
new Int_t[fNIpiv];
77 memset(fIpiv,0,fNIpiv*
sizeof(Int_t));
78 fU.ResizeTo(nrows,nrows);
84 TDecompBK::TDecompBK(Int_t row_lwb,Int_t row_upb)
86 const Int_t nrows = row_upb-row_lwb+1;
88 fIpiv =
new Int_t[fNIpiv];
89 memset(fIpiv,0,fNIpiv*
sizeof(Int_t));
90 fColLwb = fRowLwb = row_lwb;
91 fU.ResizeTo(nrows,nrows);
97 TDecompBK::TDecompBK(
const TMatrixDSym &a,Double_t tol)
99 R__ASSERT(a.IsValid());
102 fCondition = a.Norm1();
107 fNIpiv = a.GetNcols();
108 fIpiv =
new Int_t[fNIpiv];
109 memset(fIpiv,0,fNIpiv*
sizeof(Int_t));
111 const Int_t nRows = a.GetNrows();
112 fColLwb = fRowLwb = a.GetRowLwb();
113 fU.ResizeTo(nRows,nRows);
114 memcpy(fU.GetMatrixArray(),a.GetMatrixArray(),nRows*nRows*
sizeof(Double_t));
120 TDecompBK::TDecompBK(
const TDecompBK &another) : TDecompBase(another)
131 Bool_t TDecompBK::Decompose()
133 if (TestBit(kDecomposed))
return kTRUE;
135 if ( !TestBit(kMatrixSet) ) {
136 Error(
"Decompose()",
"Matrix has not been set");
143 const Double_t alpha = (1.+TMath::Sqrt(17.))/8.;
148 const Int_t n = fU.GetNcols();
149 Double_t *pU = fU.GetMatrixArray();
150 TMatrixDDiag diag(fU);
160 const Double_t absakk = TMath::Abs(diag(k));
167 TVectorD vcol = TMatrixDColumn_const(fU,k);
169 imax = TMath::LocMax(k,vcol.GetMatrixArray());
176 if (TMath::Max(absakk,colmax) <= fTol) {
181 if (absakk >= alpha*colmax) {
187 TVectorD vrow = TMatrixDRow_const(fU,imax);
189 Int_t jmax = imax+1+TMath::LocMax(k-imax,vrow.GetMatrixArray()+imax+1);
190 Double_t rowmax = vrow[jmax];
192 TVectorD vcol = TMatrixDColumn_const(fU,imax);
194 jmax = TMath::LocMax(imax,vcol.GetMatrixArray());
195 rowmax = TMath::Max(rowmax,vcol[jmax]);
198 if (absakk >= alpha*colmax*(colmax/rowmax)) {
201 }
else if( TMath::Abs(diag(imax)) >= alpha*rowmax) {
211 const Int_t kk = k-kstep+1;
214 Double_t *c_kk = pU+kk;
215 Double_t *c_kp = pU+kp;
216 for (Int_t irow = 0; irow < kp; irow++) {
217 const Double_t t = *c_kk;
224 c_kk = pU+(kp+1)*n+kk;
225 Double_t *r_kp = pU+kp*n+kp+1;
226 for (Int_t icol = 0; icol < kk-kp-1; icol++) {
227 const Double_t t = *c_kk;
234 Double_t t = diag(kk);
239 pU[(k-1)*n+k] = pU[kp*n+k];
246 if (kstep == 1 && k > 0) {
253 const Double_t r1 = 1./diag(k);
254 TMatrixDSub sub1(fU,0,k-1,0,k-1);
255 sub1.Rank1Update(TMatrixDColumn_const(fU,k),-r1);
258 TMatrixDSub sub2(fU,0,k-1,k,k);
270 Double_t *pU_k1 = pU+(k-1)*n;
271 Double_t d12 = pU_k1[k];
272 const Double_t d22 = pU_k1[k-1]/d12;
273 const Double_t d11 = diag(k)/d12;
274 const Double_t t = 1./(d11*d22-1.);
277 for (Int_t j = k-2; j >= 0; j--) {
278 Double_t *pU_j = pU+j*n;
279 const Double_t wkm1 = d12*(d11*pU_j[k-1]-pU_j[k]);
280 const Double_t wk = d12*(d22*pU_j[k]-pU_j[k-1]);
281 for (Int_t i = j; i >= 0; i--) {
282 Double_t *pU_i = pU+i*n;
283 pU_i[j] -= (pU_i[k]*wk+pU_i[k-1]*wkm1);
296 fIpiv[k-1] = -(kp+1);
303 if (!ok) SetBit(kSingular);
304 else SetBit(kDecomposed);
306 fU.Shift(fRowLwb,fRowLwb);
314 void TDecompBK::SetMatrix(
const TMatrixDSym &a)
316 R__ASSERT(a.IsValid());
321 fCondition = a.Norm1();
323 if (fNIpiv != a.GetNcols()) {
324 fNIpiv = a.GetNcols();
326 fIpiv =
new Int_t[fNIpiv];
327 memset(fIpiv,0,fNIpiv*
sizeof(Int_t));
330 const Int_t nRows = a.GetNrows();
331 fColLwb = fRowLwb = a.GetRowLwb();
332 fU.ResizeTo(nRows,nRows);
333 memcpy(fU.GetMatrixArray(),a.GetMatrixArray(),nRows*nRows*
sizeof(Double_t));
339 Bool_t TDecompBK::Solve(TVectorD &b)
341 R__ASSERT(b.IsValid());
342 if (TestBit(kSingular)) {
343 Error(
"Solve()",
"Matrix is singular");
346 if ( !TestBit(kDecomposed) ) {
348 Error(
"Solve()",
"Decomposition failed");
353 if (fU.GetNrows() != b.GetNrows() || fU.GetRowLwb() != b.GetLwb()) {
354 Error(
"Solve(TVectorD &",
"vector and matrix incompatible");
358 const Int_t n = fU.GetNrows();
360 TMatrixDDiag_const diag(fU);
361 const Double_t *pU = fU.GetMatrixArray();
362 Double_t *pb = b.GetMatrixArray();
375 const Int_t kp = fIpiv[k]-1;
377 const Double_t tmp = pb[k];
384 for (Int_t i = 0; i < k; i++)
385 pb[i] -= pU[i*n+k]*pb[k];
394 const Int_t kp = -fIpiv[k]-1;
396 const Double_t tmp = pb[k-1];
404 for (i = 0; i < k-1; i++)
405 pb[i] -= pU[i*n+k]*pb[k];
406 for (i = 0; i < k-1; i++)
407 pb[i] -= pU[i*n+k-1]*pb[k-1];
410 const Double_t *pU_k1 = pU+(k-1)*n;
411 const Double_t ukm1k = pU_k1[k];
412 const Double_t ukm1 = pU_k1[k-1]/ukm1k;
413 const Double_t uk = diag(k)/ukm1k;
414 const Double_t denom = ukm1*uk-1.;
415 const Double_t bkm1 = pb[k-1]/ukm1k;
416 const Double_t bk = pb[k]/ukm1k;
417 pb[k-1] = (uk*bkm1-bk)/denom;
418 pb[k] = (ukm1*bk-bkm1)/denom;
435 for (Int_t i = 0; i < k; i++)
436 pb[k] -= pU[i*n+k]*pb[i];
439 const Int_t kp = fIpiv[k]-1;
441 const Double_t tmp = pb[k];
451 for (i = 0; i < k; i++)
452 pb[k] -= pU[i*n+k]*pb[i];
453 for (i = 0; i < k; i++)
454 pb[k+1] -= pU[i*n+k+1]*pb[i];
457 const Int_t kp = -fIpiv[k]-1;
459 const Double_t tmp = pb[k];
473 Bool_t TDecompBK::Solve(TMatrixDColumn &cb)
475 TMatrixDBase *b =
const_cast<TMatrixDBase *
>(cb.GetMatrix());
476 R__ASSERT(b->IsValid());
477 if (TestBit(kSingular)) {
478 Error(
"Solve()",
"Matrix is singular");
481 if ( !TestBit(kDecomposed) ) {
483 Error(
"Solve()",
"Decomposition failed");
488 if (fU.GetNrows() != b->GetNrows() || fU.GetRowLwb() != b->GetRowLwb()) {
489 Error(
"Solve(TMatrixDColumn &",
"vector and matrix incompatible");
493 const Int_t n = fU.GetNrows();
495 TMatrixDDiag_const diag(fU);
496 const Double_t *pU = fU.GetMatrixArray();
497 Double_t *pcb = cb.GetPtr();
498 const Int_t inc = cb.GetInc();
512 const Int_t kp = fIpiv[k]-1;
514 const Double_t tmp = pcb[k*inc];
515 pcb[k*inc] = pcb[kp*inc];
521 for (Int_t i = 0; i < k; i++)
522 pcb[i*inc] -= pU[i*n+k]*pcb[k*inc];
525 pcb[k*inc] /= diag(k);
531 const Int_t kp = -fIpiv[k]-1;
533 const Double_t tmp = pcb[(k-1)*inc];
534 pcb[(k-1)*inc] = pcb[kp*inc];
541 for (i = 0; i < k-1; i++)
542 pcb[i*inc] -= pU[i*n+k]*pcb[k*inc];
543 for (i = 0; i < k-1; i++)
544 pcb[i*inc] -= pU[i*n+k-1]*pcb[(k-1)*inc];
547 const Double_t *pU_k1 = pU+(k-1)*n;
548 const Double_t ukm1k = pU_k1[k];
549 const Double_t ukm1 = pU_k1[k-1]/ukm1k;
550 const Double_t uk = diag(k)/ukm1k;
551 const Double_t denom = ukm1*uk-1.;
552 const Double_t bkm1 = pcb[(k-1)*inc]/ukm1k;
553 const Double_t bk = pcb[k*inc]/ukm1k;
554 pcb[(k-1)*inc] = (uk*bkm1-bk)/denom;
555 pcb[k*inc] = (ukm1*bk-bkm1)/denom;
572 for (Int_t i = 0; i < k; i++)
573 pcb[k*inc] -= pU[i*n+k]*pcb[i*inc];
576 const Int_t kp = fIpiv[k]-1;
578 const Double_t tmp = pcb[k*inc];
579 pcb[k*inc] = pcb[kp*inc];
588 for (i = 0; i < k; i++)
589 pcb[k*inc] -= pU[i*n+k]*pcb[i*inc];
590 for (i = 0; i < k; i++)
591 pcb[(k+1)*inc] -= pU[i*n+k+1]*pcb[i*inc];
594 const Int_t kp = -fIpiv[k]-1;
596 const Double_t tmp = pcb[k*inc];
597 pcb[k*inc] = pcb[kp*inc];
610 Bool_t TDecompBK::Invert(TMatrixDSym &inv)
612 if (inv.GetNrows() != GetNrows() || inv.GetRowLwb() != GetRowLwb()) {
613 Error(
"Invert(TMatrixDSym &",
"Input matrix has wrong shape");
619 const Int_t colLwb = inv.GetColLwb();
620 const Int_t colUpb = inv.GetColUpb();
621 Bool_t status = kTRUE;
622 for (Int_t icol = colLwb; icol <= colUpb && status; icol++) {
623 TMatrixDColumn b(inv,icol);
633 TMatrixDSym TDecompBK::Invert(Bool_t &status)
635 const Int_t rowLwb = GetRowLwb();
636 const Int_t rowUpb = rowLwb+GetNrows()-1;
638 TMatrixDSym inv(rowLwb,rowUpb);
640 status = Invert(inv);
648 void TDecompBK::Print(Option_t *opt)
const
650 TDecompBase::Print(opt);
652 for (Int_t i = 0; i < fNIpiv; i++)
653 printf(
"[%d] = %d\n",i,fIpiv[i]);
660 TDecompBK &TDecompBK::operator=(
const TDecompBK &source)
662 if (
this != &source) {
663 TDecompBase::operator=(source);
664 fU.ResizeTo(source.fU);
666 if (fNIpiv != source.fNIpiv) {
669 fNIpiv = source.fNIpiv;
670 fIpiv =
new Int_t[fNIpiv];
672 if (fIpiv) memcpy(fIpiv,source.fIpiv,fNIpiv*
sizeof(Int_t));