00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024 #include "TMatrixTSym.h"
00025 #include "TMatrixTLazy.h"
00026 #include "TMatrixTSymCramerInv.h"
00027 #include "TDecompLU.h"
00028 #include "TMatrixDSymEigen.h"
00029 #include "TClass.h"
00030 #include "TMath.h"
00031
00032 #ifndef R__ALPHA
00033 templateClassImp(TMatrixTSym)
00034 #endif
00035
00036
00037 template<class Element>
00038 TMatrixTSym<Element>::TMatrixTSym(Int_t no_rows)
00039 {
00040 Allocate(no_rows,no_rows,0,0,1);
00041 }
00042
00043
00044 template<class Element>
00045 TMatrixTSym<Element>::TMatrixTSym(Int_t row_lwb,Int_t row_upb)
00046 {
00047 const Int_t no_rows = row_upb-row_lwb+1;
00048 Allocate(no_rows,no_rows,row_lwb,row_lwb,1);
00049 }
00050
00051
00052 template<class Element>
00053 TMatrixTSym<Element>::TMatrixTSym(Int_t no_rows,const Element *elements,Option_t *option)
00054 {
00055
00056
00057
00058
00059
00060
00061
00062 Allocate(no_rows,no_rows);
00063 SetMatrixArray(elements,option);
00064 if (!this->IsSymmetric()) {
00065 Error("TMatrixTSym(Int_t,Element*,Option_t*)","matrix not symmetric");
00066 }
00067 }
00068
00069
00070 template<class Element>
00071 TMatrixTSym<Element>::TMatrixTSym(Int_t row_lwb,Int_t row_upb,const Element *elements,Option_t *option)
00072 {
00073
00074
00075 const Int_t no_rows = row_upb-row_lwb+1;
00076 Allocate(no_rows,no_rows,row_lwb,row_lwb);
00077 SetMatrixArray(elements,option);
00078 if (!this->IsSymmetric()) {
00079 Error("TMatrixTSym(Int_t,Int_t,Element*,Option_t*)","matrix not symmetric");
00080 }
00081 }
00082
00083
00084 template<class Element>
00085 TMatrixTSym<Element>::TMatrixTSym(const TMatrixTSym<Element> &another) : TMatrixTBase<Element>(another)
00086 {
00087 R__ASSERT(another.IsValid());
00088 Allocate(another.GetNrows(),another.GetNcols(),another.GetRowLwb(),another.GetColLwb());
00089 *this = another;
00090 }
00091
00092
00093 template<class Element>
00094 TMatrixTSym<Element>::TMatrixTSym(EMatrixCreatorsOp1 op,const TMatrixTSym<Element> &prototype)
00095 {
00096
00097
00098
00099
00100 R__ASSERT(prototype.IsValid());
00101
00102 switch(op) {
00103 case kZero:
00104 Allocate(prototype.GetNrows(),prototype.GetNcols(),
00105 prototype.GetRowLwb(),prototype.GetColLwb(),1);
00106 break;
00107
00108 case kUnit:
00109 Allocate(prototype.GetNrows(),prototype.GetNcols(),
00110 prototype.GetRowLwb(),prototype.GetColLwb(),1);
00111 this->UnitMatrix();
00112 break;
00113
00114 case kTransposed:
00115 Allocate(prototype.GetNcols(), prototype.GetNrows(),
00116 prototype.GetColLwb(),prototype.GetRowLwb());
00117 Transpose(prototype);
00118 break;
00119
00120 case kInverted:
00121 {
00122 Allocate(prototype.GetNrows(),prototype.GetNcols(),
00123 prototype.GetRowLwb(),prototype.GetColLwb(),1);
00124 *this = prototype;
00125
00126
00127 const Element oldTol = this->SetTol(std::numeric_limits<Element>::min());
00128 this->Invert();
00129 this->SetTol(oldTol);
00130 break;
00131 }
00132
00133 case kAtA:
00134 Allocate(prototype.GetNcols(),prototype.GetNcols(),prototype.GetColLwb(),prototype.GetColLwb(),1);
00135 TMult(prototype);
00136 break;
00137
00138 default:
00139 Error("TMatrixTSym(EMatrixCreatorOp1,const TMatrixTSym)",
00140 "operation %d not yet implemented", op);
00141 }
00142 }
00143
00144
00145 template<class Element>
00146 TMatrixTSym<Element>::TMatrixTSym(EMatrixCreatorsOp1 op,const TMatrixT<Element> &prototype)
00147 {
00148 R__ASSERT(prototype.IsValid());
00149
00150 switch(op) {
00151 case kAtA:
00152 Allocate(prototype.GetNcols(),prototype.GetNcols(),prototype.GetColLwb(),prototype.GetColLwb(),1);
00153 TMult(prototype);
00154 break;
00155
00156 default:
00157 Error("TMatrixTSym(EMatrixCreatorOp1,const TMatrixT)",
00158 "operation %d not yet implemented", op);
00159 }
00160 }
00161
00162
00163 template<class Element>
00164 TMatrixTSym<Element>::TMatrixTSym(const TMatrixTSym<Element> &a,EMatrixCreatorsOp2 op,const TMatrixTSym<Element> &b)
00165 {
00166 R__ASSERT(a.IsValid());
00167 R__ASSERT(b.IsValid());
00168
00169 switch(op) {
00170 case kPlus:
00171 {
00172 Plus(a,b);
00173 break;
00174 }
00175
00176 case kMinus:
00177 {
00178 Minus(a,b);
00179 break;
00180 }
00181
00182 default:
00183 Error("TMatrixTSym(EMatrixCreatorOp2)", "operation %d not yet implemented", op);
00184 }
00185 }
00186
00187
00188 template<class Element>
00189 TMatrixTSym<Element>::TMatrixTSym(const TMatrixTSymLazy<Element> &lazy_constructor)
00190 {
00191 Allocate(lazy_constructor.GetRowUpb()-lazy_constructor.GetRowLwb()+1,
00192 lazy_constructor.GetRowUpb()-lazy_constructor.GetRowLwb()+1,
00193 lazy_constructor.GetRowLwb(),lazy_constructor.GetRowLwb(),1);
00194 lazy_constructor.FillIn(*this);
00195 if (!this->IsSymmetric()) {
00196 Error("TMatrixTSym(TMatrixTSymLazy)","matrix not symmetric");
00197 }
00198 }
00199
00200
00201 template<class Element>
00202 void TMatrixTSym<Element>::Delete_m(Int_t size,Element *&m)
00203 {
00204
00205
00206 if (m) {
00207 if (size > this->kSizeMax)
00208 delete [] m;
00209 m = 0;
00210 }
00211 }
00212
00213
00214 template<class Element>
00215 Element* TMatrixTSym<Element>::New_m(Int_t size)
00216 {
00217
00218
00219
00220 if (size == 0) return 0;
00221 else {
00222 if ( size <= this->kSizeMax )
00223 return fDataStack;
00224 else {
00225 Element *heap = new Element[size];
00226 return heap;
00227 }
00228 }
00229 }
00230
00231
00232 template<class Element>
00233 Int_t TMatrixTSym<Element>::Memcpy_m(Element *newp,const Element *oldp,Int_t copySize,
00234 Int_t newSize,Int_t oldSize)
00235 {
00236
00237
00238
00239 if (copySize == 0 || oldp == newp)
00240 return 0;
00241 else {
00242 if ( newSize <= this->kSizeMax && oldSize <= this->kSizeMax ) {
00243
00244 if (newp > oldp) {
00245 for (Int_t i = copySize-1; i >= 0; i--)
00246 newp[i] = oldp[i];
00247 } else {
00248 for (Int_t i = 0; i < copySize; i++)
00249 newp[i] = oldp[i];
00250 }
00251 }
00252 else
00253 memcpy(newp,oldp,copySize*sizeof(Element));
00254 }
00255 return 0;
00256 }
00257
00258
00259 template<class Element>
00260 void TMatrixTSym<Element>::Allocate(Int_t no_rows,Int_t no_cols,Int_t row_lwb,Int_t col_lwb,
00261 Int_t init,Int_t )
00262 {
00263
00264
00265
00266 this->fIsOwner = kTRUE;
00267 this->fTol = std::numeric_limits<Element>::epsilon();
00268 this->fNrows = 0;
00269 this->fNcols = 0;
00270 this->fRowLwb = 0;
00271 this->fColLwb = 0;
00272 this->fNelems = 0;
00273 fElements = 0;
00274
00275 if (no_rows < 0 || no_cols < 0)
00276 {
00277 Error("Allocate","no_rows=%d no_cols=%d",no_rows,no_cols);
00278 this->Invalidate();
00279 return;
00280 }
00281
00282 this->MakeValid();
00283 this->fNrows = no_rows;
00284 this->fNcols = no_cols;
00285 this->fRowLwb = row_lwb;
00286 this->fColLwb = col_lwb;
00287 this->fNelems = this->fNrows*this->fNcols;
00288
00289 if (this->fNelems > 0) {
00290 fElements = New_m(this->fNelems);
00291 if (init)
00292 memset(fElements,0,this->fNelems*sizeof(Element));
00293 } else
00294 fElements = 0;
00295 }
00296
00297
00298 template<class Element>
00299 void TMatrixTSym<Element>::Plus(const TMatrixTSym<Element> &a,const TMatrixTSym<Element> &b)
00300 {
00301
00302
00303 if (gMatrixCheck) {
00304 if (!AreCompatible(a,b)) {
00305 Error("Plus","matrices not compatible");
00306 return;
00307 }
00308
00309 if (this->GetMatrixArray() == a.GetMatrixArray()) {
00310 Error("Plus","this->GetMatrixArray() == a.GetMatrixArray()");
00311 return;
00312 }
00313
00314 if (this->GetMatrixArray() == b.GetMatrixArray()) {
00315 Error("Plus","this->GetMatrixArray() == b.GetMatrixArray()");
00316 return;
00317 }
00318 }
00319
00320 const Element * ap = a.GetMatrixArray();
00321 const Element * bp = b.GetMatrixArray();
00322 Element * cp = this->GetMatrixArray();
00323 const Element * const cp_last = cp+this->fNelems;
00324
00325 while (cp < cp_last) {
00326 *cp = *ap++ + *bp++;
00327 cp++;
00328 }
00329 }
00330
00331
00332 template<class Element>
00333 void TMatrixTSym<Element>::Minus(const TMatrixTSym<Element> &a,const TMatrixTSym<Element> &b)
00334 {
00335
00336
00337 if (gMatrixCheck) {
00338 if (!AreCompatible(a,b)) {
00339 Error("Minus","matrices not compatible");
00340 return;
00341 }
00342
00343 if (this->GetMatrixArray() == a.GetMatrixArray()) {
00344 Error("Minus","this->GetMatrixArray() == a.GetMatrixArray()");
00345 return;
00346 }
00347
00348 if (this->GetMatrixArray() == b.GetMatrixArray()) {
00349 Error("Minus","this->GetMatrixArray() == b.GetMatrixArray()");
00350 return;
00351 }
00352 }
00353
00354 const Element * ap = a.GetMatrixArray();
00355 const Element * bp = b.GetMatrixArray();
00356 Element * cp = this->GetMatrixArray();
00357 const Element * const cp_last = cp+this->fNelems;
00358
00359 while (cp < cp_last) {
00360 *cp = *ap++ - *bp++;
00361 cp++;
00362 }
00363 }
00364
00365
00366 template<class Element>
00367 void TMatrixTSym<Element>::TMult(const TMatrixT<Element> &a)
00368 {
00369
00370
00371
00372 R__ASSERT(a.IsValid());
00373
00374 #ifdef CBLAS
00375 const Element *ap = a.GetMatrixArray();
00376 Element *cp = this->GetMatrixArray();
00377 if (typeid(Element) == typeid(Double_t))
00378 cblas_dgemm (CblasRowMajor,CblasTrans,CblasNoTrans,this->fNrows,this->fNcols,a.GetNrows(),
00379 1.0,ap,a.GetNcols(),ap,a.GetNcols(),1.0,cp,this->fNcols);
00380 else if (typeid(Element) != typeid(Float_t))
00381 cblas_sgemm (CblasRowMajor,CblasTrans,CblasNoTrans,fNrows,fNcols,a.GetNrows(),
00382 1.0,ap,a.GetNcols(),ap,a.GetNcols(),1.0,cp,fNcols);
00383 else
00384 Error("TMult","type %s not implemented in BLAS library",typeid(Element));
00385 #else
00386 const Int_t nb = a.GetNoElements();
00387 const Int_t ncolsa = a.GetNcols();
00388 const Int_t ncolsb = ncolsa;
00389 const Element * const ap = a.GetMatrixArray();
00390 const Element * const bp = ap;
00391 Element * cp = this->GetMatrixArray();
00392
00393 const Element *acp0 = ap;
00394 while (acp0 < ap+a.GetNcols()) {
00395 for (const Element *bcp = bp; bcp < bp+ncolsb; ) {
00396 const Element *acp = acp0;
00397 Element cij = 0;
00398 while (bcp < bp+nb) {
00399 cij += *acp * *bcp;
00400 acp += ncolsa;
00401 bcp += ncolsb;
00402 }
00403 *cp++ = cij;
00404 bcp -= nb-1;
00405 }
00406 acp0++;
00407 }
00408
00409 R__ASSERT(cp == this->GetMatrixArray()+this->fNelems && acp0 == ap+ncolsa);
00410 #endif
00411 }
00412
00413
00414 template<class Element>
00415 void TMatrixTSym<Element>::TMult(const TMatrixTSym<Element> &a)
00416 {
00417
00418
00419
00420 R__ASSERT(a.IsValid());
00421
00422 #ifdef CBLAS
00423 const Element *ap = a.GetMatrixArray();
00424 Element *cp = this->GetMatrixArray();
00425 if (typeid(Element) == typeid(Double_t))
00426 cblas_dsymm (CblasRowMajor,CblasLeft,CblasUpper,this->fNrows,this->fNcols,1.0,
00427 ap,a.GetNcols(),ap,a.GetNcols(),0.0,cp,this->fNcols);
00428 else if (typeid(Element) != typeid(Float_t))
00429 cblas_ssymm (CblasRowMajor,CblasLeft,CblasUpper,fNrows,fNcols,1.0,
00430 ap1,a.GetNcols(),ap,a.GetNcols(),0.0,cp,fNcols);
00431 else
00432 Error("TMult","type %s not implemented in BLAS library",typeid(Element));
00433 #else
00434 const Int_t nb = a.GetNoElements();
00435 const Int_t ncolsa = a.GetNcols();
00436 const Int_t ncolsb = ncolsa;
00437 const Element * const ap = a.GetMatrixArray();
00438 const Element * const bp = ap;
00439 Element * cp = this->GetMatrixArray();
00440
00441 const Element *acp0 = ap;
00442 while (acp0 < ap+a.GetNcols()) {
00443 for (const Element *bcp = bp; bcp < bp+ncolsb; ) {
00444 const Element *acp = acp0;
00445 Element cij = 0;
00446 while (bcp < bp+nb) {
00447 cij += *acp * *bcp;
00448 acp += ncolsa;
00449 bcp += ncolsb;
00450 }
00451 *cp++ = cij;
00452 bcp -= nb-1;
00453 }
00454 acp0++;
00455 }
00456
00457 R__ASSERT(cp == this->GetMatrixArray()+this->fNelems && acp0 == ap+ncolsa);
00458 #endif
00459 }
00460
00461
00462 template<class Element>
00463 TMatrixTSym<Element> &TMatrixTSym<Element>::Use(Int_t row_lwb,Int_t row_upb,Element *data)
00464 {
00465 if (gMatrixCheck && row_upb < row_lwb)
00466 {
00467 Error("Use","row_upb=%d < row_lwb=%d",row_upb,row_lwb);
00468 return *this;
00469 }
00470
00471 this->Clear();
00472 this->fNrows = row_upb-row_lwb+1;
00473 this->fNcols = this->fNrows;
00474 this->fRowLwb = row_lwb;
00475 this->fColLwb = row_lwb;
00476 this->fNelems = this->fNrows*this->fNcols;
00477 fElements = data;
00478 this->fIsOwner = kFALSE;
00479
00480 return *this;
00481 }
00482
00483
00484 template<class Element>
00485 TMatrixTSym<Element> &TMatrixTSym<Element>::GetSub(Int_t row_lwb,Int_t row_upb,
00486 TMatrixTSym<Element> &target,Option_t *option) const
00487 {
00488
00489
00490
00491
00492
00493
00494 if (gMatrixCheck) {
00495 R__ASSERT(this->IsValid());
00496 if (row_lwb < this->fRowLwb || row_lwb > this->fRowLwb+this->fNrows-1) {
00497 Error("GetSub","row_lwb out of bounds");
00498 return target;
00499 }
00500 if (row_upb < this->fRowLwb || row_upb > this->fRowLwb+this->fNrows-1) {
00501 Error("GetSub","row_upb out of bounds");
00502 return target;
00503 }
00504 if (row_upb < row_lwb) {
00505 Error("GetSub","row_upb < row_lwb");
00506 return target;
00507 }
00508 }
00509
00510 TString opt(option);
00511 opt.ToUpper();
00512 const Int_t shift = (opt.Contains("S")) ? 1 : 0;
00513
00514 Int_t row_lwb_sub;
00515 Int_t row_upb_sub;
00516 if (shift) {
00517 row_lwb_sub = 0;
00518 row_upb_sub = row_upb-row_lwb;
00519 } else {
00520 row_lwb_sub = row_lwb;
00521 row_upb_sub = row_upb;
00522 }
00523
00524 target.ResizeTo(row_lwb_sub,row_upb_sub,row_lwb_sub,row_upb_sub);
00525 const Int_t nrows_sub = row_upb_sub-row_lwb_sub+1;
00526
00527 if (target.GetRowIndexArray() && target.GetColIndexArray()) {
00528 for (Int_t irow = 0; irow < nrows_sub; irow++) {
00529 for (Int_t icol = 0; icol < nrows_sub; icol++) {
00530 target(irow+row_lwb_sub,icol+row_lwb_sub) = (*this)(row_lwb+irow,row_lwb+icol);
00531 }
00532 }
00533 } else {
00534 const Element *ap = this->GetMatrixArray()+(row_lwb-this->fRowLwb)*this->fNrows+(row_lwb-this->fRowLwb);
00535 Element *bp = target.GetMatrixArray();
00536
00537 for (Int_t irow = 0; irow < nrows_sub; irow++) {
00538 const Element *ap_sub = ap;
00539 for (Int_t icol = 0; icol < nrows_sub; icol++) {
00540 *bp++ = *ap_sub++;
00541 }
00542 ap += this->fNrows;
00543 }
00544 }
00545
00546 return target;
00547 }
00548
00549
00550 template<class Element>
00551 TMatrixTBase<Element> &TMatrixTSym<Element>::GetSub(Int_t row_lwb,Int_t row_upb,Int_t col_lwb,Int_t col_upb,
00552 TMatrixTBase<Element> &target,Option_t *option) const
00553 {
00554
00555
00556
00557
00558
00559
00560 if (gMatrixCheck) {
00561 R__ASSERT(this->IsValid());
00562 if (row_lwb < this->fRowLwb || row_lwb > this->fRowLwb+this->fNrows-1) {
00563 Error("GetSub","row_lwb out of bounds");
00564 return target;
00565 }
00566 if (col_lwb < this->fColLwb || col_lwb > this->fColLwb+this->fNcols-1) {
00567 Error("GetSub","col_lwb out of bounds");
00568 return target;
00569 }
00570 if (row_upb < this->fRowLwb || row_upb > this->fRowLwb+this->fNrows-1) {
00571 Error("GetSub","row_upb out of bounds");
00572 return target;
00573 }
00574 if (col_upb < this->fColLwb || col_upb > this->fColLwb+this->fNcols-1) {
00575 Error("GetSub","col_upb out of bounds");
00576 return target;
00577 }
00578 if (row_upb < row_lwb || col_upb < col_lwb) {
00579 Error("GetSub","row_upb < row_lwb || col_upb < col_lwb");
00580 return target;
00581 }
00582 }
00583
00584 TString opt(option);
00585 opt.ToUpper();
00586 const Int_t shift = (opt.Contains("S")) ? 1 : 0;
00587
00588 const Int_t row_lwb_sub = (shift) ? 0 : row_lwb;
00589 const Int_t row_upb_sub = (shift) ? row_upb-row_lwb : row_upb;
00590 const Int_t col_lwb_sub = (shift) ? 0 : col_lwb;
00591 const Int_t col_upb_sub = (shift) ? col_upb-col_lwb : col_upb;
00592
00593 target.ResizeTo(row_lwb_sub,row_upb_sub,col_lwb_sub,col_upb_sub);
00594 const Int_t nrows_sub = row_upb_sub-row_lwb_sub+1;
00595 const Int_t ncols_sub = col_upb_sub-col_lwb_sub+1;
00596
00597 if (target.GetRowIndexArray() && target.GetColIndexArray()) {
00598 for (Int_t irow = 0; irow < nrows_sub; irow++) {
00599 for (Int_t icol = 0; icol < ncols_sub; icol++) {
00600 target(irow+row_lwb_sub,icol+col_lwb_sub) = (*this)(row_lwb+irow,col_lwb+icol);
00601 }
00602 }
00603 } else {
00604 const Element *ap = this->GetMatrixArray()+(row_lwb-this->fRowLwb)*this->fNcols+(col_lwb-this->fColLwb);
00605 Element *bp = target.GetMatrixArray();
00606
00607 for (Int_t irow = 0; irow < nrows_sub; irow++) {
00608 const Element *ap_sub = ap;
00609 for (Int_t icol = 0; icol < ncols_sub; icol++) {
00610 *bp++ = *ap_sub++;
00611 }
00612 ap += this->fNcols;
00613 }
00614 }
00615
00616 return target;
00617 }
00618
00619
00620 template<class Element>
00621 TMatrixTSym<Element> &TMatrixTSym<Element>::SetSub(Int_t row_lwb,const TMatrixTBase<Element> &source)
00622 {
00623
00624
00625
00626 if (gMatrixCheck) {
00627 R__ASSERT(this->IsValid());
00628 R__ASSERT(source.IsValid());
00629
00630 if (!source.IsSymmetric()) {
00631 Error("SetSub","source matrix is not symmetric");
00632 return *this;
00633 }
00634 if (row_lwb < this->fRowLwb || row_lwb > this->fRowLwb+this->fNrows-1) {
00635 Error("SetSub","row_lwb outof bounds");
00636 return *this;
00637 }
00638 if (row_lwb+source.GetNrows() > this->fRowLwb+this->fNrows) {
00639 Error("SetSub","source matrix too large");
00640 return *this;
00641 }
00642 }
00643
00644 const Int_t nRows_source = source.GetNrows();
00645
00646 if (source.GetRowIndexArray() && source.GetColIndexArray()) {
00647 const Int_t rowlwb_s = source.GetRowLwb();
00648 for (Int_t irow = 0; irow < nRows_source; irow++) {
00649 for (Int_t icol = 0; icol < nRows_source; icol++) {
00650 (*this)(row_lwb+irow,row_lwb+icol) = source(rowlwb_s+irow,rowlwb_s+icol);
00651 }
00652 }
00653 } else {
00654 const Element *bp = source.GetMatrixArray();
00655 Element *ap = this->GetMatrixArray()+(row_lwb-this->fRowLwb)*this->fNrows+(row_lwb-this->fRowLwb);
00656
00657 for (Int_t irow = 0; irow < nRows_source; irow++) {
00658 Element *ap_sub = ap;
00659 for (Int_t icol = 0; icol < nRows_source; icol++) {
00660 *ap_sub++ = *bp++;
00661 }
00662 ap += this->fNrows;
00663 }
00664 }
00665
00666 return *this;
00667 }
00668
00669
00670 template<class Element>
00671 TMatrixTBase<Element> &TMatrixTSym<Element>::SetSub(Int_t row_lwb,Int_t col_lwb,const TMatrixTBase<Element> &source)
00672 {
00673
00674
00675
00676 if (gMatrixCheck) {
00677 R__ASSERT(this->IsValid());
00678 R__ASSERT(source.IsValid());
00679
00680 if (row_lwb < this->fRowLwb || row_lwb > this->fRowLwb+this->fNrows-1) {
00681 Error("SetSub","row_lwb out of bounds");
00682 return *this;
00683 }
00684 if (col_lwb < this->fColLwb || col_lwb > this->fColLwb+this->fNcols-1) {
00685 Error("SetSub","col_lwb out of bounds");
00686 return *this;
00687 }
00688
00689 if (row_lwb+source.GetNrows() > this->fRowLwb+this->fNrows || col_lwb+source.GetNcols() > this->fRowLwb+this->fNrows) {
00690 Error("SetSub","source matrix too large");
00691 return *this;
00692 }
00693 if (col_lwb+source.GetNcols() > this->fRowLwb+this->fNrows || row_lwb+source.GetNrows() > this->fRowLwb+this->fNrows) {
00694 Error("SetSub","source matrix too large");
00695 return *this;
00696 }
00697 }
00698
00699 const Int_t nRows_source = source.GetNrows();
00700 const Int_t nCols_source = source.GetNcols();
00701
00702 const Int_t rowlwb_s = source.GetRowLwb();
00703 const Int_t collwb_s = source.GetColLwb();
00704 if (row_lwb >= col_lwb) {
00705
00706 Int_t irow;
00707 for (irow = 0; irow < nRows_source; irow++) {
00708 for (Int_t icol = 0; col_lwb+icol <= row_lwb+irow &&
00709 icol < nCols_source; icol++) {
00710 (*this)(row_lwb+irow,col_lwb+icol) = source(irow+rowlwb_s,icol+collwb_s);
00711 }
00712 }
00713
00714
00715 for (irow = 0; irow < nCols_source; irow++) {
00716 for (Int_t icol = nRows_source-1; row_lwb+icol > irow+col_lwb &&
00717 icol >= 0; icol--) {
00718 (*this)(col_lwb+irow,row_lwb+icol) = source(icol+rowlwb_s,irow+collwb_s);
00719 }
00720 }
00721 } else {
00722
00723 }
00724
00725 return *this;
00726 }
00727
00728
00729 template<class Element>
00730 TMatrixTBase<Element> &TMatrixTSym<Element>::SetMatrixArray(const Element *data,Option_t *option)
00731 {
00732 TMatrixTBase<Element>::SetMatrixArray(data,option);
00733 if (!this->IsSymmetric()) {
00734 Error("SetMatrixArray","Matrix is not symmetric after Set");
00735 }
00736
00737 return *this;
00738 }
00739
00740
00741 template<class Element>
00742 TMatrixTBase<Element> &TMatrixTSym<Element>::Shift(Int_t row_shift,Int_t col_shift)
00743 {
00744 if (row_shift != col_shift) {
00745 Error("Shift","row_shift != col_shift");
00746 return *this;
00747 }
00748 return TMatrixTBase<Element>::Shift(row_shift,col_shift);
00749 }
00750
00751
00752 template<class Element>
00753 TMatrixTBase<Element> &TMatrixTSym<Element>::ResizeTo(Int_t nrows,Int_t ncols,Int_t )
00754 {
00755
00756
00757
00758
00759 R__ASSERT(this->IsValid());
00760 if (!this->fIsOwner) {
00761 Error("ResizeTo(Int_t,Int_t)","Not owner of data array,cannot resize");
00762 return *this;
00763 }
00764
00765 if (nrows != ncols) {
00766 Error("ResizeTo(Int_t,Int_t)","nrows != ncols");
00767 return *this;
00768 }
00769
00770 if (this->fNelems > 0) {
00771 if (this->fNrows == nrows && this->fNcols == ncols)
00772 return *this;
00773 else if (nrows == 0 || ncols == 0) {
00774 this->fNrows = nrows; this->fNcols = ncols;
00775 this->Clear();
00776 return *this;
00777 }
00778
00779 Element *elements_old = GetMatrixArray();
00780 const Int_t nelems_old = this->fNelems;
00781 const Int_t nrows_old = this->fNrows;
00782 const Int_t ncols_old = this->fNcols;
00783
00784 Allocate(nrows,ncols);
00785 R__ASSERT(this->IsValid());
00786
00787 Element *elements_new = GetMatrixArray();
00788
00789
00790 if (this->fNelems > this->kSizeMax || nelems_old > this->kSizeMax)
00791 memset(elements_new,0,this->fNelems*sizeof(Element));
00792 else if (this->fNelems > nelems_old)
00793 memset(elements_new+nelems_old,0,(this->fNelems-nelems_old)*sizeof(Element));
00794
00795
00796 const Int_t ncols_copy = TMath::Min(this->fNcols,ncols_old);
00797 const Int_t nrows_copy = TMath::Min(this->fNrows,nrows_old);
00798
00799 const Int_t nelems_new = this->fNelems;
00800 if (ncols_old < this->fNcols) {
00801 for (Int_t i = nrows_copy-1; i >= 0; i--) {
00802 Memcpy_m(elements_new+i*this->fNcols,elements_old+i*ncols_old,ncols_copy,
00803 nelems_new,nelems_old);
00804 if (this->fNelems <= this->kSizeMax && nelems_old <= this->kSizeMax)
00805 memset(elements_new+i*this->fNcols+ncols_copy,0,(this->fNcols-ncols_copy)*sizeof(Element));
00806 }
00807 } else {
00808 for (Int_t i = 0; i < nrows_copy; i++)
00809 Memcpy_m(elements_new+i*this->fNcols,elements_old+i*ncols_old,ncols_copy,
00810 nelems_new,nelems_old);
00811 }
00812
00813 Delete_m(nelems_old,elements_old);
00814 } else {
00815 Allocate(nrows,ncols,0,0,1);
00816 }
00817
00818 return *this;
00819 }
00820
00821
00822 template<class Element>
00823 TMatrixTBase<Element> &TMatrixTSym<Element>::ResizeTo(Int_t row_lwb,Int_t row_upb,Int_t col_lwb,Int_t col_upb,
00824 Int_t )
00825 {
00826
00827
00828
00829
00830 R__ASSERT(this->IsValid());
00831 if (!this->fIsOwner) {
00832 Error("ResizeTo(Int_t,Int_t,Int_t,Int_t)","Not owner of data array,cannot resize");
00833 return *this;
00834 }
00835
00836 if (row_lwb != col_lwb) {
00837 Error("ResizeTo(Int_t,Int_t,Int_t,Int_t)","row_lwb != col_lwb");
00838 return *this;
00839 }
00840 if (row_upb != col_upb) {
00841 Error("ResizeTo(Int_t,Int_t,Int_t,Int_t)","row_upb != col_upb");
00842 return *this;
00843 }
00844
00845 const Int_t new_nrows = row_upb-row_lwb+1;
00846 const Int_t new_ncols = col_upb-col_lwb+1;
00847
00848 if (this->fNelems > 0) {
00849
00850 if (this->fNrows == new_nrows && this->fNcols == new_ncols &&
00851 this->fRowLwb == row_lwb && this->fColLwb == col_lwb)
00852 return *this;
00853 else if (new_nrows == 0 || new_ncols == 0) {
00854 this->fNrows = new_nrows; this->fNcols = new_ncols;
00855 this->fRowLwb = row_lwb; this->fColLwb = col_lwb;
00856 this->Clear();
00857 return *this;
00858 }
00859
00860 Element *elements_old = GetMatrixArray();
00861 const Int_t nelems_old = this->fNelems;
00862 const Int_t nrows_old = this->fNrows;
00863 const Int_t ncols_old = this->fNcols;
00864 const Int_t rowLwb_old = this->fRowLwb;
00865 const Int_t colLwb_old = this->fColLwb;
00866
00867 Allocate(new_nrows,new_ncols,row_lwb,col_lwb);
00868 R__ASSERT(this->IsValid());
00869
00870 Element *elements_new = GetMatrixArray();
00871
00872
00873 if (this->fNelems > this->kSizeMax || nelems_old > this->kSizeMax)
00874 memset(elements_new,0,this->fNelems*sizeof(Element));
00875 else if (this->fNelems > nelems_old)
00876 memset(elements_new+nelems_old,0,(this->fNelems-nelems_old)*sizeof(Element));
00877
00878
00879 const Int_t rowLwb_copy = TMath::Max(this->fRowLwb,rowLwb_old);
00880 const Int_t colLwb_copy = TMath::Max(this->fColLwb,colLwb_old);
00881 const Int_t rowUpb_copy = TMath::Min(this->fRowLwb+this->fNrows-1,rowLwb_old+nrows_old-1);
00882 const Int_t colUpb_copy = TMath::Min(this->fColLwb+this->fNcols-1,colLwb_old+ncols_old-1);
00883
00884 const Int_t nrows_copy = rowUpb_copy-rowLwb_copy+1;
00885 const Int_t ncols_copy = colUpb_copy-colLwb_copy+1;
00886
00887 if (nrows_copy > 0 && ncols_copy > 0) {
00888 const Int_t colOldOff = colLwb_copy-colLwb_old;
00889 const Int_t colNewOff = colLwb_copy-this->fColLwb;
00890 if (ncols_old < this->fNcols) {
00891 for (Int_t i = nrows_copy-1; i >= 0; i--) {
00892 const Int_t iRowOld = rowLwb_copy+i-rowLwb_old;
00893 const Int_t iRowNew = rowLwb_copy+i-this->fRowLwb;
00894 Memcpy_m(elements_new+iRowNew*this->fNcols+colNewOff,
00895 elements_old+iRowOld*ncols_old+colOldOff,ncols_copy,this->fNelems,nelems_old);
00896 if (this->fNelems <= this->kSizeMax && nelems_old <= this->kSizeMax)
00897 memset(elements_new+iRowNew*this->fNcols+colNewOff+ncols_copy,0,
00898 (this->fNcols-ncols_copy)*sizeof(Element));
00899 }
00900 } else {
00901 for (Int_t i = 0; i < nrows_copy; i++) {
00902 const Int_t iRowOld = rowLwb_copy+i-rowLwb_old;
00903 const Int_t iRowNew = rowLwb_copy+i-this->fRowLwb;
00904 Memcpy_m(elements_new+iRowNew*this->fNcols+colNewOff,
00905 elements_old+iRowOld*ncols_old+colOldOff,ncols_copy,this->fNelems,nelems_old);
00906 }
00907 }
00908 }
00909
00910 Delete_m(nelems_old,elements_old);
00911 } else {
00912 Allocate(new_nrows,new_ncols,row_lwb,col_lwb,1);
00913 }
00914
00915 return *this;
00916 }
00917
00918
00919 template<class Element>
00920 Double_t TMatrixTSym<Element>::Determinant() const
00921 {
00922 const TMatrixT<Element> &tmp = *this;
00923 TDecompLU lu(tmp,this->fTol);
00924 Double_t d1,d2;
00925 lu.Det(d1,d2);
00926 return d1*TMath::Power(2.0,d2);
00927 }
00928
00929
00930 template<class Element>
00931 void TMatrixTSym<Element>::Determinant(Double_t &d1,Double_t &d2) const
00932 {
00933 const TMatrixT<Element> &tmp = *this;
00934 TDecompLU lu(tmp,this->fTol);
00935 lu.Det(d1,d2);
00936 }
00937
00938
00939 template<class Element>
00940 TMatrixTSym<Element> &TMatrixTSym<Element>::Invert(Double_t *det)
00941 {
00942
00943
00944
00945
00946
00947 R__ASSERT(this->IsValid());
00948 TMatrixD tmp(*this);
00949 if (TDecompLU::InvertLU(tmp,Double_t(this->fTol),det)) {
00950 const Double_t *p1 = tmp.GetMatrixArray();
00951 Element *p2 = this->GetMatrixArray();
00952 for (Int_t i = 0; i < this->GetNoElements(); i++)
00953 p2[i] = p1[i];
00954 }
00955
00956 return *this;
00957 }
00958
00959
00960 template<class Element>
00961 TMatrixTSym<Element> &TMatrixTSym<Element>::InvertFast(Double_t *det)
00962 {
00963
00964
00965 R__ASSERT(this->IsValid());
00966
00967 const Char_t nRows = Char_t(this->GetNrows());
00968 switch (nRows) {
00969 case 1:
00970 {
00971 Element *pM = this->GetMatrixArray();
00972 if (*pM == 0.) {
00973 Error("InvertFast","matrix is singular");
00974 *det = 0;
00975 } else {
00976 *det = *pM;
00977 *pM = 1.0/(*pM);
00978 }
00979 return *this;
00980 }
00981 case 2:
00982 {
00983 TMatrixTSymCramerInv::Inv2x2<Element>(*this,det);
00984 return *this;
00985 }
00986 case 3:
00987 {
00988 TMatrixTSymCramerInv::Inv3x3<Element>(*this,det);
00989 return *this;
00990 }
00991 case 4:
00992 {
00993 TMatrixTSymCramerInv::Inv4x4<Element>(*this,det);
00994 return *this;
00995 }
00996 case 5:
00997 {
00998 TMatrixTSymCramerInv::Inv5x5<Element>(*this,det);
00999 return *this;
01000 }
01001 case 6:
01002 {
01003 TMatrixTSymCramerInv::Inv6x6<Element>(*this,det);
01004 return *this;
01005 }
01006
01007 default:
01008 {
01009 TMatrixD tmp(*this);
01010 if (TDecompLU::InvertLU(tmp,Double_t(this->fTol),det)) {
01011 const Double_t *p1 = tmp.GetMatrixArray();
01012 Element *p2 = this->GetMatrixArray();
01013 for (Int_t i = 0; i < this->GetNoElements(); i++)
01014 p2[i] = p1[i];
01015 }
01016 return *this;
01017 }
01018 }
01019 }
01020
01021
01022 template<class Element>
01023 TMatrixTSym<Element> &TMatrixTSym<Element>::Transpose(const TMatrixTSym<Element> &source)
01024 {
01025
01026
01027 if (gMatrixCheck) {
01028 R__ASSERT(this->IsValid());
01029 R__ASSERT(source.IsValid());
01030
01031 if (this->fNrows != source.GetNcols() || this->fRowLwb != source.GetColLwb())
01032 {
01033 Error("Transpose","matrix has wrong shape");
01034 return *this;
01035 }
01036 }
01037
01038 *this = source;
01039 return *this;
01040 }
01041
01042
01043 template<class Element>
01044 TMatrixTSym<Element> &TMatrixTSym<Element>::Rank1Update(const TVectorT<Element> &v,Element alpha)
01045 {
01046
01047
01048
01049 if (gMatrixCheck) {
01050 R__ASSERT(this->IsValid());
01051 R__ASSERT(v.IsValid());
01052 if (v.GetNoElements() < this->fNrows) {
01053 Error("Rank1Update","vector too short");
01054 return *this;
01055 }
01056 }
01057
01058 const Element * const pv = v.GetMatrixArray();
01059 Element *trp = this->GetMatrixArray();
01060 Element *tcp = trp;
01061 for (Int_t i = 0; i < this->fNrows; i++) {
01062 trp += i;
01063 tcp += i*this->fNcols;
01064 const Element tmp = alpha*pv[i];
01065 for (Int_t j = i; j < this->fNcols; j++) {
01066 if (j > i) *tcp += tmp*pv[j];
01067 *trp++ += tmp*pv[j];
01068 tcp += this->fNcols;
01069 }
01070 tcp -= this->fNelems-1;
01071 }
01072
01073 return *this;
01074 }
01075
01076
01077 template<class Element>
01078 TMatrixTSym<Element> &TMatrixTSym<Element>::Similarity(const TMatrixT<Element> &b)
01079 {
01080
01081
01082
01083
01084
01085 if (gMatrixCheck) {
01086 R__ASSERT(this->IsValid());
01087 R__ASSERT(b.IsValid());
01088 if (this->fNcols != b.GetNcols() || this->fColLwb != b.GetColLwb()) {
01089 Error("Similarity(const TMatrixT &)","matrices incompatible");
01090 return *this;
01091 }
01092 }
01093
01094 const Int_t ncolsa = this->fNcols;
01095 const Int_t nb = b.GetNoElements();
01096 const Int_t nrowsb = b.GetNrows();
01097 const Int_t ncolsb = b.GetNcols();
01098
01099 const Element * const bp = b.GetMatrixArray();
01100
01101 Element work[kWorkMax];
01102 Bool_t isAllocated = kFALSE;
01103 Element *bap = work;
01104 if (nrowsb*ncolsa > kWorkMax) {
01105 isAllocated = kTRUE;
01106 bap = new Element[nrowsb*ncolsa];
01107 }
01108
01109 AMultB(bp,nb,ncolsb,this->fElements,this->fNelems,this->fNcols,bap);
01110
01111 if (nrowsb != this->fNrows)
01112 this->ResizeTo(nrowsb,nrowsb);
01113
01114 #ifdef CBLAS
01115 Element *cp = this->GetMatrixArray();
01116 if (typeid(Element) == typeid(Double_t))
01117 cblas_dgemm (CblasRowMajor,CblasNoTrans,CblasTrans,this->fNrows,this->fNcols,ba.GetNcols(),
01118 1.0,bap,ba.GetNcols(),bp,b.GetNcols(),1.0,cp,this->fNcols);
01119 else if (typeid(Element) != typeid(Float_t))
01120 cblas_sgemm (CblasRowMajor,CblasNoTrans,CblasTrans,this->fNrows,this->fNcols,ba.GetNcols(),
01121 1.0,bap,ba.GetNcols(),bp,b.GetNcols(),1.0,cp,this->fNcols);
01122 else
01123 Error("Similarity","type %s not implemented in BLAS library",typeid(Element));
01124 #else
01125 const Int_t nba = nrowsb*ncolsa;
01126 const Int_t ncolsba = ncolsa;
01127 const Element * bi1p = bp;
01128 Element * cp = this->GetMatrixArray();
01129 Element * const cp0 = cp;
01130
01131 Int_t ishift = 0;
01132 const Element *barp0 = bap;
01133 while (barp0 < bap+nba) {
01134 const Element *brp0 = bi1p;
01135 while (brp0 < bp+nb) {
01136 const Element *barp = barp0;
01137 const Element *brp = brp0;
01138 Element cij = 0;
01139 while (brp < brp0+ncolsb)
01140 cij += *barp++ * *brp++;
01141 *cp++ = cij;
01142 brp0 += ncolsb;
01143 }
01144 barp0 += ncolsba;
01145 bi1p += ncolsb;
01146 cp += ++ishift;
01147 }
01148
01149 R__ASSERT(cp == cp0+this->fNelems+ishift && barp0 == bap+nba);
01150
01151 cp = cp0;
01152 for (Int_t irow = 0; irow < this->fNrows; irow++) {
01153 const Int_t rowOff1 = irow*this->fNrows;
01154 for (Int_t icol = 0; icol < irow; icol++) {
01155 const Int_t rowOff2 = icol*this->fNrows;
01156 cp[rowOff1+icol] = cp[rowOff2+irow];
01157 }
01158 }
01159 #endif
01160
01161 if (isAllocated)
01162 delete [] bap;
01163
01164 return *this;
01165 }
01166
01167
01168 template<class Element>
01169 TMatrixTSym<Element> &TMatrixTSym<Element>::Similarity(const TMatrixTSym<Element> &b)
01170 {
01171
01172
01173
01174
01175
01176 if (gMatrixCheck) {
01177 R__ASSERT(this->IsValid());
01178 R__ASSERT(b.IsValid());
01179 if (this->fNcols != b.GetNcols() || this->fColLwb != b.GetColLwb()) {
01180 Error("Similarity(const TMatrixTSym &)","matrices incompatible");
01181 return *this;
01182 }
01183 }
01184
01185 #ifdef CBLAS
01186 const Int_t nrowsb = b.GetNrows();
01187 const Int_t ncolsa = this->GetNcols();
01188
01189 Element work[kWorkMax];
01190 Bool_t isAllocated = kFALSE;
01191 Element *abtp = work;
01192 if (this->fNcols > kWorkMax) {
01193 isAllocated = kTRUE;
01194 abtp = new Element[this->fNcols];
01195 }
01196
01197 const TMatrixT<Element> abt(*this,TMatrixT<Element>::kMultTranspose,b);
01198
01199 const Element *bp = b.GetMatrixArray();
01200 Element *cp = this->GetMatrixArray();
01201 if (typeid(Element) == typeid(Double_t))
01202 cblas_dsymm (CblasRowMajor,CblasLeft,CblasUpper,this->fNrows,this->fNcols,1.0,
01203 bp,b.GetNcols(),abtp,abt.GetNcols(),0.0,cp,this->fNcols);
01204 else if (typeid(Element) != typeid(Float_t))
01205 cblas_ssymm (CblasRowMajor,CblasLeft,CblasUpper,this->fNrows,this->fNcols,1.0,
01206 bp,b.GetNcols(),abtp,abt.GetNcols(),0.0,cp,this->fNcols);
01207 else
01208 Error("Similarity","type %s not implemented in BLAS library",typeid(Element));
01209
01210 if (isAllocated)
01211 delete [] abtp;
01212 #else
01213 const Int_t ncolsa = this->GetNcols();
01214 const Int_t nb = b.GetNoElements();
01215 const Int_t nrowsb = b.GetNrows();
01216 const Int_t ncolsb = b.GetNcols();
01217
01218 const Element * const bp = b.GetMatrixArray();
01219
01220 Element work[kWorkMax];
01221 Bool_t isAllocated = kFALSE;
01222 Element *bap = work;
01223 if (nrowsb*ncolsa > kWorkMax) {
01224 isAllocated = kTRUE;
01225 bap = new Element[nrowsb*ncolsa];
01226 }
01227
01228 AMultB(bp,nb,ncolsb,this->fElements,this->fNelems,this->fNcols,bap);
01229
01230 const Int_t nba = nrowsb*ncolsa;
01231 const Int_t ncolsba = ncolsa;
01232 const Element * bi1p = bp;
01233 Element * cp = this->GetMatrixArray();
01234 Element * const cp0 = cp;
01235
01236 Int_t ishift = 0;
01237 const Element *barp0 = bap;
01238 while (barp0 < bap+nba) {
01239 const Element *brp0 = bi1p;
01240 while (brp0 < bp+nb) {
01241 const Element *barp = barp0;
01242 const Element *brp = brp0;
01243 Element cij = 0;
01244 while (brp < brp0+ncolsb)
01245 cij += *barp++ * *brp++;
01246 *cp++ = cij;
01247 brp0 += ncolsb;
01248 }
01249 barp0 += ncolsba;
01250 bi1p += ncolsb;
01251 cp += ++ishift;
01252 }
01253
01254 R__ASSERT(cp == cp0+this->fNelems+ishift && barp0 == bap+nba);
01255
01256 cp = cp0;
01257 for (Int_t irow = 0; irow < this->fNrows; irow++) {
01258 const Int_t rowOff1 = irow*this->fNrows;
01259 for (Int_t icol = 0; icol < irow; icol++) {
01260 const Int_t rowOff2 = icol*this->fNrows;
01261 cp[rowOff1+icol] = cp[rowOff2+irow];
01262 }
01263 }
01264
01265 if (isAllocated)
01266 delete [] bap;
01267 #endif
01268
01269 return *this;
01270 }
01271
01272
01273 template<class Element>
01274 Element TMatrixTSym<Element>::Similarity(const TVectorT<Element> &v) const
01275 {
01276
01277
01278 if (gMatrixCheck) {
01279 R__ASSERT(this->IsValid());
01280 R__ASSERT(v.IsValid());
01281 if (this->fNcols != v.GetNrows() || this->fColLwb != v.GetLwb()) {
01282 Error("Similarity(const TVectorT &)","vector and matrix incompatible");
01283 return -1.;
01284 }
01285 }
01286
01287 const Element *mp = this->GetMatrixArray();
01288 const Element *vp = v.GetMatrixArray();
01289
01290 Element sum1 = 0;
01291 const Element * const vp_first = vp;
01292 const Element * const vp_last = vp+v.GetNrows();
01293 while (vp < vp_last) {
01294 Element sum2 = 0;
01295 for (const Element *sp = vp_first; sp < vp_last; )
01296 sum2 += *mp++ * *sp++;
01297 sum1 += sum2 * *vp++;
01298 }
01299
01300 R__ASSERT(mp == this->GetMatrixArray()+this->GetNoElements());
01301
01302 return sum1;
01303 }
01304
01305
01306 template<class Element>
01307 TMatrixTSym<Element> &TMatrixTSym<Element>::SimilarityT(const TMatrixT<Element> &b)
01308 {
01309
01310
01311
01312
01313 if (gMatrixCheck) {
01314 R__ASSERT(this->IsValid());
01315 R__ASSERT(b.IsValid());
01316 if (this->fNrows != b.GetNrows() || this->fRowLwb != b.GetRowLwb()) {
01317 Error("SimilarityT(const TMatrixT &)","matrices incompatible");
01318 return *this;
01319 }
01320 }
01321
01322 const Int_t ncolsb = b.GetNcols();
01323 const Int_t ncolsa = this->GetNcols();
01324
01325 Element work[kWorkMax];
01326 Bool_t isAllocated = kFALSE;
01327 Element *btap = work;
01328 if (ncolsb*ncolsa > kWorkMax) {
01329 isAllocated = kTRUE;
01330 btap = new Element[ncolsb*ncolsa];
01331 }
01332
01333 TMatrixT<Element> bta; bta.Use(ncolsb,ncolsa,btap);
01334 bta.TMult(b,*this);
01335
01336 if (ncolsb != this->fNcols)
01337 this->ResizeTo(ncolsb,ncolsb);
01338
01339 #ifdef CBLAS
01340 const Element *bp = b.GetMatrixArray();
01341 Element *cp = this->GetMatrixArray();
01342 if (typeid(Element) == typeid(Double_t))
01343 cblas_dgemm (CblasRowMajor,CblasNoTrans,CblasNoTrans,this->fNrows,this->fNcols,bta.GetNcols(),
01344 1.0,btap,bta.GetNcols(),bp,b.GetNcols(),1.0,cp,this->fNcols);
01345 else if (typeid(Element) != typeid(Float_t))
01346 cblas_sgemm (CblasRowMajor,CblasNoTrans,CblasNoTrans,this->fNrows,this->fNcols,bta.GetNcols(),
01347 1.0,btap,bta.GetNcols(),bp,b.GetNcols(),1.0,cp,this->fNcols);
01348 else
01349 Error("similarityT","type %s not implemented in BLAS library",typeid(Element));
01350 #else
01351 const Int_t nbta = bta.GetNoElements();
01352 const Int_t nb = b.GetNoElements();
01353 const Int_t ncolsbta = bta.GetNcols();
01354 const Element * const bp = b.GetMatrixArray();
01355 Element * cp = this->GetMatrixArray();
01356 Element * const cp0 = cp;
01357
01358 Int_t ishift = 0;
01359 const Element *btarp0 = btap;
01360 const Element *bcp0 = bp;
01361 while (btarp0 < btap+nbta) {
01362 for (const Element *bcp = bcp0; bcp < bp+ncolsb; ) {
01363 const Element *btarp = btarp0;
01364 Element cij = 0;
01365 while (bcp < bp+nb) {
01366 cij += *btarp++ * *bcp;
01367 bcp += ncolsb;
01368 }
01369 *cp++ = cij;
01370 bcp -= nb-1;
01371 }
01372 btarp0 += ncolsbta;
01373 bcp0++;
01374 cp += ++ishift;
01375 }
01376
01377 R__ASSERT(cp == cp0+this->fNelems+ishift && btarp0 == btap+nbta);
01378
01379 cp = cp0;
01380 for (Int_t irow = 0; irow < this->fNrows; irow++) {
01381 const Int_t rowOff1 = irow*this->fNrows;
01382 for (Int_t icol = 0; icol < irow; icol++) {
01383 const Int_t rowOff2 = icol*this->fNrows;
01384 cp[rowOff1+icol] = cp[rowOff2+irow];
01385 }
01386 }
01387 #endif
01388
01389 if (isAllocated)
01390 delete [] btap;
01391
01392 return *this;
01393 }
01394
01395
01396 template<class Element>
01397 TMatrixTSym<Element> &TMatrixTSym<Element>::operator=(const TMatrixTSym<Element> &source)
01398 {
01399 if (gMatrixCheck && !AreCompatible(*this,source)) {
01400 Error("operator=","matrices not compatible");
01401 return *this;
01402 }
01403
01404 if (this->GetMatrixArray() != source.GetMatrixArray()) {
01405 TObject::operator=(source);
01406 memcpy(this->GetMatrixArray(),source.fElements,this->fNelems*sizeof(Element));
01407 }
01408 return *this;
01409 }
01410
01411
01412 template<class Element>
01413 TMatrixTSym<Element> &TMatrixTSym<Element>::operator=(const TMatrixTSymLazy<Element> &lazy_constructor)
01414 {
01415 R__ASSERT(this->IsValid());
01416
01417 if (lazy_constructor.fRowUpb != this->GetRowUpb() ||
01418 lazy_constructor.fRowLwb != this->GetRowLwb()) {
01419 Error("operator=(const TMatrixTSymLazy&)", "matrix is incompatible with "
01420 "the assigned Lazy matrix");
01421 return *this;
01422 }
01423
01424 lazy_constructor.FillIn(*this);
01425 return *this;
01426 }
01427
01428
01429 template<class Element>
01430 TMatrixTSym<Element> &TMatrixTSym<Element>::operator=(Element val)
01431 {
01432
01433
01434 R__ASSERT(this->IsValid());
01435
01436 Element *ep = fElements;
01437 const Element * const ep_last = ep+this->fNelems;
01438 while (ep < ep_last)
01439 *ep++ = val;
01440
01441 return *this;
01442 }
01443
01444
01445 template<class Element>
01446 TMatrixTSym<Element> &TMatrixTSym<Element>::operator+=(Element val)
01447 {
01448
01449
01450 R__ASSERT(this->IsValid());
01451
01452 Element *ep = fElements;
01453 const Element * const ep_last = ep+this->fNelems;
01454 while (ep < ep_last)
01455 *ep++ += val;
01456
01457 return *this;
01458 }
01459
01460
01461 template<class Element>
01462 TMatrixTSym<Element> &TMatrixTSym<Element>::operator-=(Element val)
01463 {
01464
01465
01466 R__ASSERT(this->IsValid());
01467
01468 Element *ep = fElements;
01469 const Element * const ep_last = ep+this->fNelems;
01470 while (ep < ep_last)
01471 *ep++ -= val;
01472
01473 return *this;
01474 }
01475
01476
01477 template<class Element>
01478 TMatrixTSym<Element> &TMatrixTSym<Element>::operator*=(Element val)
01479 {
01480
01481
01482 R__ASSERT(this->IsValid());
01483
01484 Element *ep = fElements;
01485 const Element * const ep_last = ep+this->fNelems;
01486 while (ep < ep_last)
01487 *ep++ *= val;
01488
01489 return *this;
01490 }
01491
01492
01493 template<class Element>
01494 TMatrixTSym<Element> &TMatrixTSym<Element>::operator+=(const TMatrixTSym<Element> &source)
01495 {
01496
01497
01498 if (gMatrixCheck && !AreCompatible(*this,source)) {
01499 Error("operator+=","matrices not compatible");
01500 return *this;
01501 }
01502
01503 const Element *sp = source.GetMatrixArray();
01504 Element *tp = this->GetMatrixArray();
01505 const Element * const tp_last = tp+this->fNelems;
01506 while (tp < tp_last)
01507 *tp++ += *sp++;
01508
01509 return *this;
01510 }
01511
01512
01513 template<class Element>
01514 TMatrixTSym<Element> &TMatrixTSym<Element>::operator-=(const TMatrixTSym<Element> &source)
01515 {
01516
01517
01518 if (gMatrixCheck && !AreCompatible(*this,source)) {
01519 Error("operator-=","matrices not compatible");
01520 return *this;
01521 }
01522
01523 const Element *sp = source.GetMatrixArray();
01524 Element *tp = this->GetMatrixArray();
01525 const Element * const tp_last = tp+this->fNelems;
01526 while (tp < tp_last)
01527 *tp++ -= *sp++;
01528
01529 return *this;
01530 }
01531
01532
01533 template<class Element>
01534 TMatrixTBase<Element> &TMatrixTSym<Element>::Apply(const TElementActionT<Element> &action)
01535 {
01536 R__ASSERT(this->IsValid());
01537
01538 Element val = 0;
01539 Element *trp = this->GetMatrixArray();
01540 Element *tcp = trp;
01541 for (Int_t i = 0; i < this->fNrows; i++) {
01542 trp += i;
01543 tcp += i*this->fNcols;
01544 for (Int_t j = i; j < this->fNcols; j++) {
01545 action.Operation(val);
01546 if (j > i) *tcp = val;
01547 *trp++ = val;
01548 tcp += this->fNcols;
01549 }
01550 tcp -= this->fNelems-1;
01551 }
01552
01553 return *this;
01554 }
01555
01556
01557 template<class Element>
01558 TMatrixTBase<Element> &TMatrixTSym<Element>::Apply(const TElementPosActionT<Element> &action)
01559 {
01560
01561
01562
01563 R__ASSERT(this->IsValid());
01564
01565 Element val = 0;
01566 Element *trp = this->GetMatrixArray();
01567 Element *tcp = trp;
01568 for (Int_t i = 0; i < this->fNrows; i++) {
01569 action.fI = i+this->fRowLwb;
01570 trp += i;
01571 tcp += i*this->fNcols;
01572 for (Int_t j = i; j < this->fNcols; j++) {
01573 action.fJ = j+this->fColLwb;
01574 action.Operation(val);
01575 if (j > i) *tcp = val;
01576 *trp++ = val;
01577 tcp += this->fNcols;
01578 }
01579 tcp -= this->fNelems-1;
01580 }
01581
01582 return *this;
01583 }
01584
01585
01586 template<class Element>
01587 TMatrixTBase<Element> &TMatrixTSym<Element>::Randomize(Element alpha,Element beta,Double_t &seed)
01588 {
01589
01590
01591 if (gMatrixCheck) {
01592 R__ASSERT(this->IsValid());
01593 if (this->fNrows != this->fNcols || this->fRowLwb != this->fColLwb) {
01594 Error("Randomize(Element,Element,Element &","matrix should be square");
01595 return *this;
01596 }
01597 }
01598
01599 const Element scale = beta-alpha;
01600 const Element shift = alpha/scale;
01601
01602 Element *ep = GetMatrixArray();
01603 for (Int_t i = 0; i < this->fNrows; i++) {
01604 const Int_t off = i*this->fNcols;
01605 for (Int_t j = 0; j <= i; j++) {
01606 ep[off+j] = scale*(Drand(seed)+shift);
01607 if (i != j) {
01608 ep[j*this->fNcols+i] = ep[off+j];
01609 }
01610 }
01611 }
01612
01613 return *this;
01614 }
01615
01616
01617 template<class Element>
01618 TMatrixTSym<Element> &TMatrixTSym<Element>::RandomizePD(Element alpha,Element beta,Double_t &seed)
01619 {
01620
01621
01622 if (gMatrixCheck) {
01623 R__ASSERT(this->IsValid());
01624 if (this->fNrows != this->fNcols || this->fRowLwb != this->fColLwb) {
01625 Error("RandomizeSym(Element,Element,Element &","matrix should be square");
01626 return *this;
01627 }
01628 }
01629
01630 const Element scale = beta-alpha;
01631 const Element shift = alpha/scale;
01632
01633 Element *ep = GetMatrixArray();
01634 Int_t i;
01635 for (i = 0; i < this->fNrows; i++) {
01636 const Int_t off = i*this->fNcols;
01637 for (Int_t j = 0; j <= i; j++)
01638 ep[off+j] = scale*(Drand(seed)+shift);
01639 }
01640
01641 for (i = this->fNrows-1; i >= 0; i--) {
01642 const Int_t off1 = i*this->fNcols;
01643 for (Int_t j = i; j >= 0; j--) {
01644 const Int_t off2 = j*this->fNcols;
01645 ep[off1+j] *= ep[off2+j];
01646 for (Int_t k = j-1; k >= 0; k--) {
01647 ep[off1+j] += ep[off1+k]*ep[off2+k];
01648 }
01649 if (i != j)
01650 ep[off2+i] = ep[off1+j];
01651 }
01652 }
01653
01654 return *this;
01655 }
01656
01657
01658 template<class Element>
01659 const TMatrixT<Element> TMatrixTSym<Element>::EigenVectors(TVectorT<Element> &eigenValues) const
01660 {
01661
01662
01663
01664 TMatrixDSym tmp = *this;
01665 TMatrixDSymEigen eigen(tmp);
01666 eigenValues.ResizeTo(this->fNrows);
01667 eigenValues = eigen.GetEigenValues();
01668 return eigen.GetEigenVectors();
01669 }
01670
01671
01672 template<class Element>
01673 Bool_t operator==(const TMatrixTSym<Element> &m1,const TMatrixTSym<Element> &m2)
01674 {
01675
01676
01677 if (!AreCompatible(m1,m2)) return kFALSE;
01678 return (memcmp(m1.GetMatrixArray(),m2.GetMatrixArray(),
01679 m1.GetNoElements()*sizeof(Element)) == 0);
01680 }
01681
01682
01683 template<class Element>
01684 TMatrixTSym<Element> operator+(const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2)
01685 {
01686 TMatrixTSym<Element> target(source1);
01687 target += source2;
01688 return target;
01689 }
01690
01691
01692 template<class Element>
01693 TMatrixTSym<Element> operator+(const TMatrixTSym<Element> &source1,Element val)
01694 {
01695 TMatrixTSym<Element> target(source1);
01696 target += val;
01697 return target;
01698 }
01699
01700
01701 template<class Element>
01702 TMatrixTSym<Element> operator+(Element val,const TMatrixTSym<Element> &source1)
01703 {
01704 return operator+(source1,val);
01705 }
01706
01707
01708 template<class Element>
01709 TMatrixTSym<Element> operator-(const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2)
01710 {
01711 TMatrixTSym<Element> target(source1);
01712 target -= source2;
01713 return target;
01714 }
01715
01716
01717 template<class Element>
01718 TMatrixTSym<Element> operator-(const TMatrixTSym<Element> &source1,Element val)
01719 {
01720 TMatrixTSym<Element> target(source1);
01721 target -= val;
01722 return target;
01723 }
01724
01725
01726 template<class Element>
01727 TMatrixTSym<Element> operator-(Element val,const TMatrixTSym<Element> &source1)
01728 {
01729 return Element(-1.0)*operator-(source1,val);
01730 }
01731
01732
01733 template<class Element>
01734 TMatrixTSym<Element> operator*(const TMatrixTSym<Element> &source1,Element val)
01735 {
01736 TMatrixTSym<Element> target(source1);
01737 target *= val;
01738 return target;
01739 }
01740
01741
01742 template<class Element>
01743 TMatrixTSym<Element> operator*(Element val,const TMatrixTSym<Element> &source1)
01744 {
01745 return operator*(source1,val);
01746 }
01747
01748
01749 template<class Element>
01750 TMatrixTSym<Element> operator&&(const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2)
01751 {
01752
01753
01754 TMatrixTSym<Element> target;
01755
01756 if (gMatrixCheck && !AreCompatible(source1,source2)) {
01757 Error("operator&&(const TMatrixTSym&,const TMatrixTSym&)","matrices not compatible");
01758 return target;
01759 }
01760
01761 target.ResizeTo(source1);
01762
01763 const Element *sp1 = source1.GetMatrixArray();
01764 const Element *sp2 = source2.GetMatrixArray();
01765 Element *tp = target.GetMatrixArray();
01766 const Element * const tp_last = tp+target.GetNoElements();
01767 while (tp < tp_last)
01768 *tp++ = (*sp1++ != 0.0 && *sp2++ != 0.0);
01769
01770 return target;
01771 }
01772
01773
01774 template<class Element>
01775 TMatrixTSym<Element> operator||(const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2)
01776 {
01777
01778
01779 TMatrixTSym<Element> target;
01780
01781 if (gMatrixCheck && !AreCompatible(source1,source2)) {
01782 Error("operator||(const TMatrixTSym&,const TMatrixTSym&)","matrices not compatible");
01783 return target;
01784 }
01785
01786 target.ResizeTo(source1);
01787
01788 const Element *sp1 = source1.GetMatrixArray();
01789 const Element *sp2 = source2.GetMatrixArray();
01790 Element *tp = target.GetMatrixArray();
01791 const Element * const tp_last = tp+target.GetNoElements();
01792 while (tp < tp_last)
01793 *tp++ = (*sp1++ != 0.0 || *sp2++ != 0.0);
01794
01795 return target;
01796 }
01797
01798
01799 template<class Element>
01800 TMatrixTSym<Element> operator>(const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2)
01801 {
01802
01803
01804 TMatrixTSym<Element> target;
01805
01806 if (gMatrixCheck && !AreCompatible(source1,source2)) {
01807 Error("operator>(const TMatrixTSym&,const TMatrixTSym&)","matrices not compatible");
01808 return target;
01809 }
01810
01811 target.ResizeTo(source1);
01812
01813 const Element *sp1 = source1.GetMatrixArray();
01814 const Element *sp2 = source2.GetMatrixArray();
01815 Element *tp = target.GetMatrixArray();
01816 const Element * const tp_last = tp+target.GetNoElements();
01817 while (tp < tp_last) {
01818 *tp++ = (*sp1) > (*sp2); sp1++; sp2++;
01819 }
01820
01821 return target;
01822 }
01823
01824
01825 template<class Element>
01826 TMatrixTSym<Element> operator>=(const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2)
01827 {
01828
01829
01830 TMatrixTSym<Element> target;
01831
01832 if (gMatrixCheck && !AreCompatible(source1,source2)) {
01833 Error("operator>=(const TMatrixTSym&,const TMatrixTSym&)","matrices not compatible");
01834 return target;
01835 }
01836
01837 target.ResizeTo(source1);
01838
01839 const Element *sp1 = source1.GetMatrixArray();
01840 const Element *sp2 = source2.GetMatrixArray();
01841 Element *tp = target.GetMatrixArray();
01842 const Element * const tp_last = tp+target.GetNoElements();
01843 while (tp < tp_last) {
01844 *tp++ = (*sp1) >= (*sp2); sp1++; sp2++;
01845 }
01846
01847 return target;
01848 }
01849
01850
01851 template<class Element>
01852 TMatrixTSym<Element> operator<=(const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2)
01853 {
01854
01855
01856 TMatrixTSym<Element> target;
01857
01858 if (gMatrixCheck && !AreCompatible(source1,source2)) {
01859 Error("operator<=(const TMatrixTSym&,const TMatrixTSym&)","matrices not compatible");
01860 return target;
01861 }
01862
01863 target.ResizeTo(source1);
01864
01865 const Element *sp1 = source1.GetMatrixArray();
01866 const Element *sp2 = source2.GetMatrixArray();
01867 Element *tp = target.GetMatrixArray();
01868 const Element * const tp_last = tp+target.GetNoElements();
01869 while (tp < tp_last) {
01870 *tp++ = (*sp1) <= (*sp2); sp1++; sp2++;
01871 }
01872
01873 return target;
01874 }
01875
01876
01877 template<class Element>
01878 TMatrixTSym<Element> operator<(const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2)
01879 {
01880
01881
01882 TMatrixTSym<Element> target;
01883
01884 if (gMatrixCheck && !AreCompatible(source1,source2)) {
01885 Error("operator<(const TMatrixTSym&,const TMatrixTSym&)","matrices not compatible");
01886 return target;
01887 }
01888
01889 target.ResizeTo(source1);
01890
01891 const Element *sp1 = source1.GetMatrixArray();
01892 const Element *sp2 = source2.GetMatrixArray();
01893 Element *tp = target.GetMatrixArray();
01894 const Element * const tp_last = tp+target.GetNoElements();
01895 while (tp < tp_last) {
01896 *tp++ = (*sp1) < (*sp2); sp1++; sp2++;
01897 }
01898
01899 return target;
01900 }
01901
01902
01903 template<class Element>
01904 TMatrixTSym<Element> &Add(TMatrixTSym<Element> &target,Element scalar,const TMatrixTSym<Element> &source)
01905 {
01906
01907
01908 if (gMatrixCheck && !AreCompatible(target,source)) {
01909 ::Error("Add","matrices not compatible");
01910 return target;
01911 }
01912
01913 const Int_t nrows = target.GetNrows();
01914 const Int_t ncols = target.GetNcols();
01915 const Int_t nelems = target.GetNoElements();
01916 const Element *sp = source.GetMatrixArray();
01917 Element *trp = target.GetMatrixArray();
01918 Element *tcp = target.GetMatrixArray();
01919 for (Int_t i = 0; i < nrows; i++) {
01920 sp += i;
01921 trp += i;
01922 tcp += i*ncols;
01923 for (Int_t j = i; j < ncols; j++) {
01924 const Element tmp = scalar * *sp++;
01925 if (j > i) *tcp += tmp;
01926 *trp++ += tmp;
01927 tcp += ncols;
01928 }
01929 tcp -= nelems-1;
01930 }
01931
01932 return target;
01933 }
01934
01935
01936 template<class Element>
01937 TMatrixTSym<Element> &ElementMult(TMatrixTSym<Element> &target,const TMatrixTSym<Element> &source)
01938 {
01939
01940
01941 if (gMatrixCheck && !AreCompatible(target,source)) {
01942 ::Error("ElementMult","matrices not compatible");
01943 return target;
01944 }
01945
01946 const Int_t nrows = target.GetNrows();
01947 const Int_t ncols = target.GetNcols();
01948 const Int_t nelems = target.GetNoElements();
01949 const Element *sp = source.GetMatrixArray();
01950 Element *trp = target.GetMatrixArray();
01951 Element *tcp = target.GetMatrixArray();
01952 for (Int_t i = 0; i < nrows; i++) {
01953 sp += i;
01954 trp += i;
01955 tcp += i*ncols;
01956 for (Int_t j = i; j < ncols; j++) {
01957 if (j > i) *tcp *= *sp;
01958 *trp++ *= *sp++;
01959 tcp += ncols;
01960 }
01961 tcp -= nelems-1;
01962 }
01963
01964 return target;
01965 }
01966
01967
01968 template<class Element>
01969 TMatrixTSym<Element> &ElementDiv(TMatrixTSym<Element> &target,const TMatrixTSym<Element> &source)
01970 {
01971
01972
01973 if (gMatrixCheck && !AreCompatible(target,source)) {
01974 ::Error("ElementDiv","matrices not compatible");
01975 return target;
01976 }
01977
01978 const Int_t nrows = target.GetNrows();
01979 const Int_t ncols = target.GetNcols();
01980 const Int_t nelems = target.GetNoElements();
01981 const Element *sp = source.GetMatrixArray();
01982 Element *trp = target.GetMatrixArray();
01983 Element *tcp = target.GetMatrixArray();
01984 for (Int_t i = 0; i < nrows; i++) {
01985 sp += i;
01986 trp += i;
01987 tcp += i*ncols;
01988 for (Int_t j = i; j < ncols; j++) {
01989 if (*sp != 0.0) {
01990 if (j > i) *tcp /= *sp;
01991 *trp++ /= *sp++;
01992 } else {
01993 const Int_t irow = (sp-source.GetMatrixArray())/source.GetNcols();
01994 const Int_t icol = (sp-source.GetMatrixArray())%source.GetNcols();
01995 Error("ElementDiv","source (%d,%d) is zero",irow,icol);
01996 trp++;
01997 }
01998 tcp += ncols;
01999 }
02000 tcp -= nelems-1;
02001 }
02002
02003 return target;
02004 }
02005
02006
02007 template<class Element>
02008 void TMatrixTSym<Element>::Streamer(TBuffer &R__b)
02009 {
02010
02011
02012 if (R__b.IsReading()) {
02013 UInt_t R__s, R__c;
02014 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
02015 Clear();
02016 R__b.ReadClassBuffer(TMatrixTBase<Element>::Class(),this,R__v,R__s,R__c);
02017 fElements = new Element[this->fNelems];
02018 Int_t i;
02019 for (i = 0; i < this->fNrows; i++) {
02020 R__b.ReadFastArray(fElements+i*this->fNcols+i,this->fNcols-i);
02021 }
02022
02023 for (i = 0; i < this->fNrows; i++) {
02024 for (Int_t j = 0; j < i; j++) {
02025 fElements[i*this->fNcols+j] = fElements[j*this->fNrows+i];
02026 }
02027 }
02028 if (this->fNelems <= this->kSizeMax) {
02029 memcpy(fDataStack,fElements,this->fNelems*sizeof(Element));
02030 delete [] fElements;
02031 fElements = fDataStack;
02032 }
02033 } else {
02034 R__b.WriteClassBuffer(TMatrixTBase<Element>::Class(),this);
02035
02036 for (Int_t i = 0; i < this->fNrows; i++) {
02037 R__b.WriteFastArray(fElements+i*this->fNcols+i,this->fNcols-i);
02038 }
02039 }
02040 }
02041
02042 #ifndef ROOT_TMatrixFSymfwd
02043 #include "TMatrixFSymfwd.h"
02044 #endif
02045
02046 template class TMatrixTSym<Float_t>;
02047
02048 template Bool_t operator== <Float_t>(const TMatrixFSym &source1,const TMatrixFSym &source2);
02049 template TMatrixFSym operator+ <Float_t>(const TMatrixFSym &source1,const TMatrixFSym &source2);
02050 template TMatrixFSym operator+ <Float_t>(const TMatrixFSym &source1, Float_t val);
02051 template TMatrixFSym operator+ <Float_t>( Float_t val ,const TMatrixFSym &source2);
02052 template TMatrixFSym operator- <Float_t>(const TMatrixFSym &source1,const TMatrixFSym &source2);
02053 template TMatrixFSym operator- <Float_t>(const TMatrixFSym &source1, Float_t val);
02054 template TMatrixFSym operator- <Float_t>( Float_t val ,const TMatrixFSym &source2);
02055 template TMatrixFSym operator* <Float_t>(const TMatrixFSym &source, Float_t val );
02056 template TMatrixFSym operator* <Float_t>( Float_t val, const TMatrixFSym &source );
02057 template TMatrixFSym operator&& <Float_t>(const TMatrixFSym &source1,const TMatrixFSym &source2);
02058 template TMatrixFSym operator|| <Float_t>(const TMatrixFSym &source1,const TMatrixFSym &source2);
02059 template TMatrixFSym operator> <Float_t>(const TMatrixFSym &source1,const TMatrixFSym &source2);
02060 template TMatrixFSym operator>= <Float_t>(const TMatrixFSym &source1,const TMatrixFSym &source2);
02061 template TMatrixFSym operator<= <Float_t>(const TMatrixFSym &source1,const TMatrixFSym &source2);
02062 template TMatrixFSym operator< <Float_t>(const TMatrixFSym &source1,const TMatrixFSym &source2);
02063
02064 template TMatrixFSym &Add <Float_t>(TMatrixFSym &target, Float_t scalar,const TMatrixFSym &source);
02065 template TMatrixFSym &ElementMult<Float_t>(TMatrixFSym &target,const TMatrixFSym &source);
02066 template TMatrixFSym &ElementDiv <Float_t>(TMatrixFSym &target,const TMatrixFSym &source);
02067
02068 #ifndef ROOT_TMatrixDSymfwd
02069 #include "TMatrixDSymfwd.h"
02070 #endif
02071
02072 template class TMatrixTSym<Double_t>;
02073
02074 template Bool_t operator== <Double_t>(const TMatrixDSym &source1,const TMatrixDSym &source2);
02075 template TMatrixDSym operator+ <Double_t>(const TMatrixDSym &source1,const TMatrixDSym &source2);
02076 template TMatrixDSym operator+ <Double_t>(const TMatrixDSym &source1, Double_t val);
02077 template TMatrixDSym operator+ <Double_t>( Double_t val ,const TMatrixDSym &source2);
02078 template TMatrixDSym operator- <Double_t>(const TMatrixDSym &source1,const TMatrixDSym &source2);
02079 template TMatrixDSym operator- <Double_t>(const TMatrixDSym &source1, Double_t val);
02080 template TMatrixDSym operator- <Double_t>( Double_t val ,const TMatrixDSym &source2);
02081 template TMatrixDSym operator* <Double_t>(const TMatrixDSym &source, Double_t val );
02082 template TMatrixDSym operator* <Double_t>( Double_t val, const TMatrixDSym &source );
02083 template TMatrixDSym operator&& <Double_t>(const TMatrixDSym &source1,const TMatrixDSym &source2);
02084 template TMatrixDSym operator|| <Double_t>(const TMatrixDSym &source1,const TMatrixDSym &source2);
02085 template TMatrixDSym operator> <Double_t>(const TMatrixDSym &source1,const TMatrixDSym &source2);
02086 template TMatrixDSym operator>= <Double_t>(const TMatrixDSym &source1,const TMatrixDSym &source2);
02087 template TMatrixDSym operator<= <Double_t>(const TMatrixDSym &source1,const TMatrixDSym &source2);
02088 template TMatrixDSym operator< <Double_t>(const TMatrixDSym &source1,const TMatrixDSym &source2);
02089
02090 template TMatrixDSym &Add <Double_t>(TMatrixDSym &target, Double_t scalar,const TMatrixDSym &source);
02091 template TMatrixDSym &ElementMult<Double_t>(TMatrixDSym &target,const TMatrixDSym &source);
02092 template TMatrixDSym &ElementDiv <Double_t>(TMatrixDSym &target,const TMatrixDSym &source);