ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TMatrixTSym.cxx
Go to the documentation of this file.
1 // @(#)root/matrix:$Id$
2 // Authors: Fons Rademakers, Eddy Offermann Nov 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 //////////////////////////////////////////////////////////////////////////
13 // //
14 // TMatrixTSym //
15 // //
16 // Template class of a symmetric matrix in the linear algebra package //
17 // //
18 // Note that in this implementation both matrix element m[i][j] and //
19 // m[j][i] are updated and stored in memory . However, when making the //
20 // object persistent only the upper right triangle is stored . //
21 // //
22 //////////////////////////////////////////////////////////////////////////
23 
24 #include "TMatrixTSym.h"
25 #include "TBuffer.h"
26 #include "TMatrixTLazy.h"
27 #include "TMatrixTSymCramerInv.h"
28 #include "TDecompLU.h"
29 #include "TMatrixDSymEigen.h"
30 #include "TClass.h"
31 #include "TMath.h"
32 
34 
35 
36 ////////////////////////////////////////////////////////////////////////////////
37 
38 template<class Element>
39 TMatrixTSym<Element>::TMatrixTSym(Int_t no_rows)
40 {
41  Allocate(no_rows,no_rows,0,0,1);
42 }
43 
44 ////////////////////////////////////////////////////////////////////////////////
45 
46 template<class Element>
48 {
49  const Int_t no_rows = row_upb-row_lwb+1;
50  Allocate(no_rows,no_rows,row_lwb,row_lwb,1);
51 }
52 
53 ////////////////////////////////////////////////////////////////////////////////
54 /// option="F": array elements contains the matrix stored column-wise
55 /// like in Fortran, so a[i,j] = elements[i+no_rows*j],
56 /// else it is supposed that array elements are stored row-wise
57 /// a[i,j] = elements[i*no_cols+j]
58 ///
59 /// array elements are copied
60 
61 template<class Element>
62 TMatrixTSym<Element>::TMatrixTSym(Int_t no_rows,const Element *elements,Option_t *option)
63 {
64  Allocate(no_rows,no_rows);
65  SetMatrixArray(elements,option);
66  if (!this->IsSymmetric()) {
67  Error("TMatrixTSym(Int_t,Element*,Option_t*)","matrix not symmetric");
68  }
69 }
70 
71 ////////////////////////////////////////////////////////////////////////////////
72 /// array elements are copied
73 
74 template<class Element>
75 TMatrixTSym<Element>::TMatrixTSym(Int_t row_lwb,Int_t row_upb,const Element *elements,Option_t *option)
76 {
77  const Int_t no_rows = row_upb-row_lwb+1;
78  Allocate(no_rows,no_rows,row_lwb,row_lwb);
79  SetMatrixArray(elements,option);
80  if (!this->IsSymmetric()) {
81  Error("TMatrixTSym(Int_t,Int_t,Element*,Option_t*)","matrix not symmetric");
82  }
83 }
84 
85 ////////////////////////////////////////////////////////////////////////////////
86 
87 template<class Element>
89 {
90  R__ASSERT(another.IsValid());
91  Allocate(another.GetNrows(),another.GetNcols(),another.GetRowLwb(),another.GetColLwb());
92  *this = another;
93 }
94 
95 ////////////////////////////////////////////////////////////////////////////////
96 /// Create a matrix applying a specific operation to the prototype.
97 /// Example: TMatrixTSym<Element> a(10,12); ...; TMatrixTSym<Element> b(TMatrixT::kTransposed, a);
98 /// Supported operations are: kZero, kUnit, kTransposed, kInverted and kAtA.
99 
100 template<class Element>
102 {
103  R__ASSERT(prototype.IsValid());
104 
105  switch(op) {
106  case kZero:
107  Allocate(prototype.GetNrows(),prototype.GetNcols(),
108  prototype.GetRowLwb(),prototype.GetColLwb(),1);
109  break;
110 
111  case kUnit:
112  Allocate(prototype.GetNrows(),prototype.GetNcols(),
113  prototype.GetRowLwb(),prototype.GetColLwb(),1);
114  this->UnitMatrix();
115  break;
116 
117  case kTransposed:
118  Allocate(prototype.GetNcols(), prototype.GetNrows(),
119  prototype.GetColLwb(),prototype.GetRowLwb());
120  Transpose(prototype);
121  break;
122 
123  case kInverted:
124  {
125  Allocate(prototype.GetNrows(),prototype.GetNcols(),
126  prototype.GetRowLwb(),prototype.GetColLwb(),1);
127  *this = prototype;
128  // Since the user can not control the tolerance of this newly created matrix
129  // we put it to the smallest possible number
130  const Element oldTol = this->SetTol(std::numeric_limits<Element>::min());
131  this->Invert();
132  this->SetTol(oldTol);
133  break;
134  }
135 
136  case kAtA:
137  Allocate(prototype.GetNcols(),prototype.GetNcols(),prototype.GetColLwb(),prototype.GetColLwb(),1);
138  TMult(prototype);
139  break;
140 
141  default:
142  Error("TMatrixTSym(EMatrixCreatorOp1,const TMatrixTSym)",
143  "operation %d not yet implemented", op);
144  }
145 }
146 
147 ////////////////////////////////////////////////////////////////////////////////
148 
149 template<class Element>
151 {
152  R__ASSERT(prototype.IsValid());
153 
154  switch(op) {
155  case kAtA:
156  Allocate(prototype.GetNcols(),prototype.GetNcols(),prototype.GetColLwb(),prototype.GetColLwb(),1);
157  TMult(prototype);
158  break;
159 
160  default:
161  Error("TMatrixTSym(EMatrixCreatorOp1,const TMatrixT)",
162  "operation %d not yet implemented", op);
163  }
164 }
165 
166 ////////////////////////////////////////////////////////////////////////////////
167 
168 template<class Element>
170 {
171  R__ASSERT(a.IsValid());
172  R__ASSERT(b.IsValid());
173 
174  switch(op) {
175  case kPlus:
176  {
177  Allocate(a.GetNcols(),a.GetNcols(),a.GetColLwb(),a.GetColLwb(),1);
178  Plus(a,b);
179  break;
180  }
181 
182  case kMinus:
183  {
184  Allocate(a.GetNcols(),a.GetNcols(),a.GetColLwb(),a.GetColLwb(),1);
185  Minus(a,b);
186  break;
187  }
188 
189  default:
190  Error("TMatrixTSym(EMatrixCreatorOp2)", "operation %d not yet implemented", op);
191  }
192 }
193 
194 ////////////////////////////////////////////////////////////////////////////////
195 
196 template<class Element>
198 {
199  Allocate(lazy_constructor.GetRowUpb()-lazy_constructor.GetRowLwb()+1,
200  lazy_constructor.GetRowUpb()-lazy_constructor.GetRowLwb()+1,
201  lazy_constructor.GetRowLwb(),lazy_constructor.GetRowLwb(),1);
202  lazy_constructor.FillIn(*this);
203  if (!this->IsSymmetric()) {
204  Error("TMatrixTSym(TMatrixTSymLazy)","matrix not symmetric");
205  }
206 }
207 
208 ////////////////////////////////////////////////////////////////////////////////
209 /// delete data pointer m, if it was assigned on the heap
210 
211 template<class Element>
213 {
214  if (m) {
215  if (size > this->kSizeMax)
216  delete [] m;
217  m = 0;
218  }
219 }
220 
221 ////////////////////////////////////////////////////////////////////////////////
222 /// return data pointer . if requested size <= kSizeMax, assign pointer
223 /// to the stack space
224 
225 template<class Element>
227 {
228  if (size == 0) return 0;
229  else {
230  if ( size <= this->kSizeMax )
231  return fDataStack;
232  else {
233  Element *heap = new Element[size];
234  return heap;
235  }
236  }
237 }
238 
239 ////////////////////////////////////////////////////////////////////////////////
240 /// copy copySize doubles from *oldp to *newp . However take care of the
241 /// situation where both pointers are assigned to the same stack space
242 
243 template<class Element>
244 Int_t TMatrixTSym<Element>::Memcpy_m(Element *newp,const Element *oldp,Int_t copySize,
245  Int_t newSize,Int_t oldSize)
246 {
247  if (copySize == 0 || oldp == newp)
248  return 0;
249  else {
250  if ( newSize <= this->kSizeMax && oldSize <= this->kSizeMax ) {
251  // both pointers are inside fDataStack, be careful with copy direction !
252  if (newp > oldp) {
253  for (Int_t i = copySize-1; i >= 0; i--)
254  newp[i] = oldp[i];
255  } else {
256  for (Int_t i = 0; i < copySize; i++)
257  newp[i] = oldp[i];
258  }
259  }
260  else
261  memcpy(newp,oldp,copySize*sizeof(Element));
262  }
263  return 0;
264 }
265 
266 ////////////////////////////////////////////////////////////////////////////////
267 /// Allocate new matrix. Arguments are number of rows, columns, row
268 /// lowerbound (0 default) and column lowerbound (0 default).
269 
270 template<class Element>
271 void TMatrixTSym<Element>::Allocate(Int_t no_rows,Int_t no_cols,Int_t row_lwb,Int_t col_lwb,
272  Int_t init,Int_t /*nr_nonzeros*/)
273 {
274  this->fIsOwner = kTRUE;
276  this->fNrows = 0;
277  this->fNcols = 0;
278  this->fRowLwb = 0;
279  this->fColLwb = 0;
280  this->fNelems = 0;
281  fElements = 0;
282 
283  if (no_rows < 0 || no_cols < 0)
284  {
285  Error("Allocate","no_rows=%d no_cols=%d",no_rows,no_cols);
286  this->Invalidate();
287  return;
288  }
289 
290  this->MakeValid();
291  this->fNrows = no_rows;
292  this->fNcols = no_cols;
293  this->fRowLwb = row_lwb;
294  this->fColLwb = col_lwb;
295  this->fNelems = this->fNrows*this->fNcols;
296 
297  if (this->fNelems > 0) {
298  fElements = New_m(this->fNelems);
299  if (init)
300  memset(fElements,0,this->fNelems*sizeof(Element));
301  } else
302  fElements = 0;
303 }
304 
305 ////////////////////////////////////////////////////////////////////////////////
306 /// Symmetric matrix summation. Create a matrix C such that C = A + B.
307 
308 template<class Element>
310 {
311  if (gMatrixCheck) {
312  if (!AreCompatible(a,b)) {
313  Error("Plus","matrices not compatible");
314  return;
315  }
316 
317  if (this->GetMatrixArray() == a.GetMatrixArray()) {
318  Error("Plus","this->GetMatrixArray() == a.GetMatrixArray()");
319  return;
320  }
321 
322  if (this->GetMatrixArray() == b.GetMatrixArray()) {
323  Error("Plus","this->GetMatrixArray() == b.GetMatrixArray()");
324  return;
325  }
326  }
327 
328  const Element * ap = a.GetMatrixArray();
329  const Element * bp = b.GetMatrixArray();
330  Element * cp = this->GetMatrixArray();
331  const Element * const cp_last = cp+this->fNelems;
332 
333  while (cp < cp_last) {
334  *cp = *ap++ + *bp++;
335  cp++;
336  }
337 }
338 
339 ////////////////////////////////////////////////////////////////////////////////
340 /// Symmetric matrix summation. Create a matrix C such that C = A + B.
341 
342 template<class Element>
344 {
345  if (gMatrixCheck) {
346  if (!AreCompatible(a,b)) {
347  Error("Minus","matrices not compatible");
348  return;
349  }
350 
351  if (this->GetMatrixArray() == a.GetMatrixArray()) {
352  Error("Minus","this->GetMatrixArray() == a.GetMatrixArray()");
353  return;
354  }
355 
356  if (this->GetMatrixArray() == b.GetMatrixArray()) {
357  Error("Minus","this->GetMatrixArray() == b.GetMatrixArray()");
358  return;
359  }
360  }
361 
362  const Element * ap = a.GetMatrixArray();
363  const Element * bp = b.GetMatrixArray();
364  Element * cp = this->GetMatrixArray();
365  const Element * const cp_last = cp+this->fNelems;
366 
367  while (cp < cp_last) {
368  *cp = *ap++ - *bp++;
369  cp++;
370  }
371 }
372 
373 ////////////////////////////////////////////////////////////////////////////////
374 /// Create a matrix C such that C = A' * A. In other words,
375 /// c[i,j] = SUM{ a[k,i] * a[k,j] }.
376 
377 template<class Element>
379 {
380  R__ASSERT(a.IsValid());
381 
382 #ifdef CBLAS
383  const Element *ap = a.GetMatrixArray();
384  Element *cp = this->GetMatrixArray();
385  if (typeid(Element) == typeid(Double_t))
386  cblas_dgemm (CblasRowMajor,CblasTrans,CblasNoTrans,this->fNrows,this->fNcols,a.GetNrows(),
387  1.0,ap,a.GetNcols(),ap,a.GetNcols(),1.0,cp,this->fNcols);
388  else if (typeid(Element) != typeid(Float_t))
389  cblas_sgemm (CblasRowMajor,CblasTrans,CblasNoTrans,fNrows,fNcols,a.GetNrows(),
390  1.0,ap,a.GetNcols(),ap,a.GetNcols(),1.0,cp,fNcols);
391  else
392  Error("TMult","type %s not implemented in BLAS library",typeid(Element));
393 #else
394  const Int_t nb = a.GetNoElements();
395  const Int_t ncolsa = a.GetNcols();
396  const Int_t ncolsb = ncolsa;
397  const Element * const ap = a.GetMatrixArray();
398  const Element * const bp = ap;
399  Element * cp = this->GetMatrixArray();
400 
401  const Element *acp0 = ap; // Pointer to A[i,0];
402  while (acp0 < ap+a.GetNcols()) {
403  for (const Element *bcp = bp; bcp < bp+ncolsb; ) { // Pointer to the j-th column of A, Start bcp = A[0,0]
404  const Element *acp = acp0; // Pointer to the i-th column of A, reset to A[0,i]
405  Element cij = 0;
406  while (bcp < bp+nb) { // Scan the i-th column of A and
407  cij += *acp * *bcp; // the j-th col of A
408  acp += ncolsa;
409  bcp += ncolsb;
410  }
411  *cp++ = cij;
412  bcp -= nb-1; // Set bcp to the (j+1)-th col
413  }
414  acp0++; // Set acp0 to the (i+1)-th col
415  }
416 
417  R__ASSERT(cp == this->GetMatrixArray()+this->fNelems && acp0 == ap+ncolsa);
418 #endif
419 }
420 
421 ////////////////////////////////////////////////////////////////////////////////
422 /// Matrix multiplication, with A symmetric
423 /// Create a matrix C such that C = A' * A = A * A = A * A'
424 
425 template<class Element>
427 {
428  R__ASSERT(a.IsValid());
429 
430 #ifdef CBLAS
431  const Element *ap = a.GetMatrixArray();
432  Element *cp = this->GetMatrixArray();
433  if (typeid(Element) == typeid(Double_t))
434  cblas_dsymm (CblasRowMajor,CblasLeft,CblasUpper,this->fNrows,this->fNcols,1.0,
435  ap,a.GetNcols(),ap,a.GetNcols(),0.0,cp,this->fNcols);
436  else if (typeid(Element) != typeid(Float_t))
437  cblas_ssymm (CblasRowMajor,CblasLeft,CblasUpper,fNrows,fNcols,1.0,
438  ap1,a.GetNcols(),ap,a.GetNcols(),0.0,cp,fNcols);
439  else
440  Error("TMult","type %s not implemented in BLAS library",typeid(Element));
441 #else
442  const Int_t nb = a.GetNoElements();
443  const Int_t ncolsa = a.GetNcols();
444  const Int_t ncolsb = ncolsa;
445  const Element * const ap = a.GetMatrixArray();
446  const Element * const bp = ap;
447  Element * cp = this->GetMatrixArray();
448 
449  const Element *acp0 = ap; // Pointer to A[i,0];
450  while (acp0 < ap+a.GetNcols()) {
451  for (const Element *bcp = bp; bcp < bp+ncolsb; ) { // Pointer to the j-th column of A, Start bcp = A[0,0]
452  const Element *acp = acp0; // Pointer to the i-th column of A, reset to A[0,i]
453  Element cij = 0;
454  while (bcp < bp+nb) { // Scan the i-th column of A and
455  cij += *acp * *bcp; // the j-th col of A
456  acp += ncolsa;
457  bcp += ncolsb;
458  }
459  *cp++ = cij;
460  bcp -= nb-1; // Set bcp to the (j+1)-th col
461  }
462  acp0++; // Set acp0 to the (i+1)-th col
463  }
464 
465  R__ASSERT(cp == this->GetMatrixArray()+this->fNelems && acp0 == ap+ncolsa);
466 #endif
467 }
468 
469 ////////////////////////////////////////////////////////////////////////////////
470 
471 template<class Element>
473 {
474  if (gMatrixCheck && row_upb < row_lwb)
475  {
476  Error("Use","row_upb=%d < row_lwb=%d",row_upb,row_lwb);
477  return *this;
478  }
479 
480  this->Clear();
481  this->fNrows = row_upb-row_lwb+1;
482  this->fNcols = this->fNrows;
483  this->fRowLwb = row_lwb;
484  this->fColLwb = row_lwb;
485  this->fNelems = this->fNrows*this->fNcols;
486  fElements = data;
487  this->fIsOwner = kFALSE;
488 
489  return *this;
490 }
491 
492 ////////////////////////////////////////////////////////////////////////////////
493 /// Get submatrix [row_lwb..row_upb][row_lwb..row_upb]; The indexing range of the
494 /// returned matrix depends on the argument option:
495 ///
496 /// option == "S" : return [0..row_upb-row_lwb+1][0..row_upb-row_lwb+1] (default)
497 /// else : return [row_lwb..row_upb][row_lwb..row_upb]
498 
499 template<class Element>
501  TMatrixTSym<Element> &target,Option_t *option) const
502 {
503  if (gMatrixCheck) {
504  R__ASSERT(this->IsValid());
505  if (row_lwb < this->fRowLwb || row_lwb > this->fRowLwb+this->fNrows-1) {
506  Error("GetSub","row_lwb out of bounds");
507  return target;
508  }
509  if (row_upb < this->fRowLwb || row_upb > this->fRowLwb+this->fNrows-1) {
510  Error("GetSub","row_upb out of bounds");
511  return target;
512  }
513  if (row_upb < row_lwb) {
514  Error("GetSub","row_upb < row_lwb");
515  return target;
516  }
517  }
518 
519  TString opt(option);
520  opt.ToUpper();
521  const Int_t shift = (opt.Contains("S")) ? 1 : 0;
522 
523  Int_t row_lwb_sub;
524  Int_t row_upb_sub;
525  if (shift) {
526  row_lwb_sub = 0;
527  row_upb_sub = row_upb-row_lwb;
528  } else {
529  row_lwb_sub = row_lwb;
530  row_upb_sub = row_upb;
531  }
532 
533  target.ResizeTo(row_lwb_sub,row_upb_sub,row_lwb_sub,row_upb_sub);
534  const Int_t nrows_sub = row_upb_sub-row_lwb_sub+1;
535 
536  if (target.GetRowIndexArray() && target.GetColIndexArray()) {
537  for (Int_t irow = 0; irow < nrows_sub; irow++) {
538  for (Int_t icol = 0; icol < nrows_sub; icol++) {
539  target(irow+row_lwb_sub,icol+row_lwb_sub) = (*this)(row_lwb+irow,row_lwb+icol);
540  }
541  }
542  } else {
543  const Element *ap = this->GetMatrixArray()+(row_lwb-this->fRowLwb)*this->fNrows+(row_lwb-this->fRowLwb);
544  Element *bp = target.GetMatrixArray();
545 
546  for (Int_t irow = 0; irow < nrows_sub; irow++) {
547  const Element *ap_sub = ap;
548  for (Int_t icol = 0; icol < nrows_sub; icol++) {
549  *bp++ = *ap_sub++;
550  }
551  ap += this->fNrows;
552  }
553  }
554 
555  return target;
556 }
557 
558 ////////////////////////////////////////////////////////////////////////////////
559 /// Get submatrix [row_lwb..row_upb][col_lwb..col_upb]; The indexing range of the
560 /// returned matrix depends on the argument option:
561 ///
562 /// option == "S" : return [0..row_upb-row_lwb+1][0..col_upb-col_lwb+1] (default)
563 /// else : return [row_lwb..row_upb][col_lwb..col_upb]
564 
565 template<class Element>
567  TMatrixTBase<Element> &target,Option_t *option) const
568 {
569  if (gMatrixCheck) {
570  R__ASSERT(this->IsValid());
571  if (row_lwb < this->fRowLwb || row_lwb > this->fRowLwb+this->fNrows-1) {
572  Error("GetSub","row_lwb out of bounds");
573  return target;
574  }
575  if (col_lwb < this->fColLwb || col_lwb > this->fColLwb+this->fNcols-1) {
576  Error("GetSub","col_lwb out of bounds");
577  return target;
578  }
579  if (row_upb < this->fRowLwb || row_upb > this->fRowLwb+this->fNrows-1) {
580  Error("GetSub","row_upb out of bounds");
581  return target;
582  }
583  if (col_upb < this->fColLwb || col_upb > this->fColLwb+this->fNcols-1) {
584  Error("GetSub","col_upb out of bounds");
585  return target;
586  }
587  if (row_upb < row_lwb || col_upb < col_lwb) {
588  Error("GetSub","row_upb < row_lwb || col_upb < col_lwb");
589  return target;
590  }
591  }
592 
593  TString opt(option);
594  opt.ToUpper();
595  const Int_t shift = (opt.Contains("S")) ? 1 : 0;
596 
597  const Int_t row_lwb_sub = (shift) ? 0 : row_lwb;
598  const Int_t row_upb_sub = (shift) ? row_upb-row_lwb : row_upb;
599  const Int_t col_lwb_sub = (shift) ? 0 : col_lwb;
600  const Int_t col_upb_sub = (shift) ? col_upb-col_lwb : col_upb;
601 
602  target.ResizeTo(row_lwb_sub,row_upb_sub,col_lwb_sub,col_upb_sub);
603  const Int_t nrows_sub = row_upb_sub-row_lwb_sub+1;
604  const Int_t ncols_sub = col_upb_sub-col_lwb_sub+1;
605 
606  if (target.GetRowIndexArray() && target.GetColIndexArray()) {
607  for (Int_t irow = 0; irow < nrows_sub; irow++) {
608  for (Int_t icol = 0; icol < ncols_sub; icol++) {
609  target(irow+row_lwb_sub,icol+col_lwb_sub) = (*this)(row_lwb+irow,col_lwb+icol);
610  }
611  }
612  } else {
613  const Element *ap = this->GetMatrixArray()+(row_lwb-this->fRowLwb)*this->fNcols+(col_lwb-this->fColLwb);
614  Element *bp = target.GetMatrixArray();
615 
616  for (Int_t irow = 0; irow < nrows_sub; irow++) {
617  const Element *ap_sub = ap;
618  for (Int_t icol = 0; icol < ncols_sub; icol++) {
619  *bp++ = *ap_sub++;
620  }
621  ap += this->fNcols;
622  }
623  }
624 
625  return target;
626 }
627 
628 ////////////////////////////////////////////////////////////////////////////////
629 /// Insert matrix source starting at [row_lwb][row_lwb], thereby overwriting the part
630 /// [row_lwb..row_lwb+nrows_source][row_lwb..row_lwb+nrows_source];
631 
632 template<class Element>
634 {
635  if (gMatrixCheck) {
636  R__ASSERT(this->IsValid());
637  R__ASSERT(source.IsValid());
638 
639  if (!source.IsSymmetric()) {
640  Error("SetSub","source matrix is not symmetric");
641  return *this;
642  }
643  if (row_lwb < this->fRowLwb || row_lwb > this->fRowLwb+this->fNrows-1) {
644  Error("SetSub","row_lwb outof bounds");
645  return *this;
646  }
647  if (row_lwb+source.GetNrows() > this->fRowLwb+this->fNrows) {
648  Error("SetSub","source matrix too large");
649  return *this;
650  }
651  }
652 
653  const Int_t nRows_source = source.GetNrows();
654 
655  if (source.GetRowIndexArray() && source.GetColIndexArray()) {
656  const Int_t rowlwb_s = source.GetRowLwb();
657  for (Int_t irow = 0; irow < nRows_source; irow++) {
658  for (Int_t icol = 0; icol < nRows_source; icol++) {
659  (*this)(row_lwb+irow,row_lwb+icol) = source(rowlwb_s+irow,rowlwb_s+icol);
660  }
661  }
662  } else {
663  const Element *bp = source.GetMatrixArray();
664  Element *ap = this->GetMatrixArray()+(row_lwb-this->fRowLwb)*this->fNrows+(row_lwb-this->fRowLwb);
665 
666  for (Int_t irow = 0; irow < nRows_source; irow++) {
667  Element *ap_sub = ap;
668  for (Int_t icol = 0; icol < nRows_source; icol++) {
669  *ap_sub++ = *bp++;
670  }
671  ap += this->fNrows;
672  }
673  }
674 
675  return *this;
676 }
677 
678 ////////////////////////////////////////////////////////////////////////////////
679 /// Insert matrix source starting at [row_lwb][col_lwb] in a symmetric fashion, thereby overwriting the part
680 /// [row_lwb..row_lwb+nrows_source][row_lwb..row_lwb+nrows_source];
681 
682 template<class Element>
684 {
685  if (gMatrixCheck) {
686  R__ASSERT(this->IsValid());
687  R__ASSERT(source.IsValid());
688 
689  if (row_lwb < this->fRowLwb || row_lwb > this->fRowLwb+this->fNrows-1) {
690  Error("SetSub","row_lwb out of bounds");
691  return *this;
692  }
693  if (col_lwb < this->fColLwb || col_lwb > this->fColLwb+this->fNcols-1) {
694  Error("SetSub","col_lwb out of bounds");
695  return *this;
696  }
697 
698  if (row_lwb+source.GetNrows() > this->fRowLwb+this->fNrows || col_lwb+source.GetNcols() > this->fRowLwb+this->fNrows) {
699  Error("SetSub","source matrix too large");
700  return *this;
701  }
702  if (col_lwb+source.GetNcols() > this->fRowLwb+this->fNrows || row_lwb+source.GetNrows() > this->fRowLwb+this->fNrows) {
703  Error("SetSub","source matrix too large");
704  return *this;
705  }
706  }
707 
708  const Int_t nRows_source = source.GetNrows();
709  const Int_t nCols_source = source.GetNcols();
710 
711  const Int_t rowlwb_s = source.GetRowLwb();
712  const Int_t collwb_s = source.GetColLwb();
713  if (row_lwb >= col_lwb) {
714  // lower triangle
715  Int_t irow;
716  for (irow = 0; irow < nRows_source; irow++) {
717  for (Int_t icol = 0; col_lwb+icol <= row_lwb+irow &&
718  icol < nCols_source; icol++) {
719  (*this)(row_lwb+irow,col_lwb+icol) = source(irow+rowlwb_s,icol+collwb_s);
720  }
721  }
722 
723  // upper triangle
724  for (irow = 0; irow < nCols_source; irow++) {
725  for (Int_t icol = nRows_source-1; row_lwb+icol > irow+col_lwb &&
726  icol >= 0; icol--) {
727  (*this)(col_lwb+irow,row_lwb+icol) = source(icol+rowlwb_s,irow+collwb_s);
728  }
729  }
730  } else {
731 
732  }
733 
734  return *this;
735 }
736 
737 ////////////////////////////////////////////////////////////////////////////////
738 
739 template<class Element>
741 {
743  if (!this->IsSymmetric()) {
744  Error("SetMatrixArray","Matrix is not symmetric after Set");
745  }
746 
747  return *this;
748 }
749 
750 ////////////////////////////////////////////////////////////////////////////////
751 
752 template<class Element>
754 {
755  if (row_shift != col_shift) {
756  Error("Shift","row_shift != col_shift");
757  return *this;
758  }
759  return TMatrixTBase<Element>::Shift(row_shift,col_shift);
760 }
761 
762 ////////////////////////////////////////////////////////////////////////////////
763 /// Set size of the matrix to nrows x ncols
764 /// New dynamic elements are created, the overlapping part of the old ones are
765 /// copied to the new structures, then the old elements are deleted.
766 
767 template<class Element>
769 {
770  R__ASSERT(this->IsValid());
771  if (!this->fIsOwner) {
772  Error("ResizeTo(Int_t,Int_t)","Not owner of data array,cannot resize");
773  return *this;
774  }
775 
776  if (nrows != ncols) {
777  Error("ResizeTo(Int_t,Int_t)","nrows != ncols");
778  return *this;
779  }
780 
781  if (this->fNelems > 0) {
782  if (this->fNrows == nrows && this->fNcols == ncols)
783  return *this;
784  else if (nrows == 0 || ncols == 0) {
785  this->fNrows = nrows; this->fNcols = ncols;
786  this->Clear();
787  return *this;
788  }
789 
790  Element *elements_old = GetMatrixArray();
791  const Int_t nelems_old = this->fNelems;
792  const Int_t nrows_old = this->fNrows;
793  const Int_t ncols_old = this->fNcols;
794 
795  Allocate(nrows,ncols);
796  R__ASSERT(this->IsValid());
797 
798  Element *elements_new = GetMatrixArray();
799  // new memory should be initialized but be careful not to wipe out the stack
800  // storage. Initialize all when old or new storage was on the heap
801  if (this->fNelems > this->kSizeMax || nelems_old > this->kSizeMax)
802  memset(elements_new,0,this->fNelems*sizeof(Element));
803  else if (this->fNelems > nelems_old)
804  memset(elements_new+nelems_old,0,(this->fNelems-nelems_old)*sizeof(Element));
805 
806  // Copy overlap
807  const Int_t ncols_copy = TMath::Min(this->fNcols,ncols_old);
808  const Int_t nrows_copy = TMath::Min(this->fNrows,nrows_old);
809 
810  const Int_t nelems_new = this->fNelems;
811  if (ncols_old < this->fNcols) {
812  for (Int_t i = nrows_copy-1; i >= 0; i--) {
813  Memcpy_m(elements_new+i*this->fNcols,elements_old+i*ncols_old,ncols_copy,
814  nelems_new,nelems_old);
815  if (this->fNelems <= this->kSizeMax && nelems_old <= this->kSizeMax)
816  memset(elements_new+i*this->fNcols+ncols_copy,0,(this->fNcols-ncols_copy)*sizeof(Element));
817  }
818  } else {
819  for (Int_t i = 0; i < nrows_copy; i++)
820  Memcpy_m(elements_new+i*this->fNcols,elements_old+i*ncols_old,ncols_copy,
821  nelems_new,nelems_old);
822  }
823 
824  Delete_m(nelems_old,elements_old);
825  } else {
826  Allocate(nrows,ncols,0,0,1);
827  }
828 
829  return *this;
830 }
831 
832 ////////////////////////////////////////////////////////////////////////////////
833 /// Set size of the matrix to [row_lwb:row_upb] x [col_lwb:col_upb]
834 /// New dynamic elemenst are created, the overlapping part of the old ones are
835 /// copied to the new structures, then the old elements are deleted.
836 
837 template<class Element>
839  Int_t /*nr_nonzeros*/)
840 {
841  R__ASSERT(this->IsValid());
842  if (!this->fIsOwner) {
843  Error("ResizeTo(Int_t,Int_t,Int_t,Int_t)","Not owner of data array,cannot resize");
844  return *this;
845  }
846 
847  if (row_lwb != col_lwb) {
848  Error("ResizeTo(Int_t,Int_t,Int_t,Int_t)","row_lwb != col_lwb");
849  return *this;
850  }
851  if (row_upb != col_upb) {
852  Error("ResizeTo(Int_t,Int_t,Int_t,Int_t)","row_upb != col_upb");
853  return *this;
854  }
855 
856  const Int_t new_nrows = row_upb-row_lwb+1;
857  const Int_t new_ncols = col_upb-col_lwb+1;
858 
859  if (this->fNelems > 0) {
860 
861  if (this->fNrows == new_nrows && this->fNcols == new_ncols &&
862  this->fRowLwb == row_lwb && this->fColLwb == col_lwb)
863  return *this;
864  else if (new_nrows == 0 || new_ncols == 0) {
865  this->fNrows = new_nrows; this->fNcols = new_ncols;
866  this->fRowLwb = row_lwb; this->fColLwb = col_lwb;
867  this->Clear();
868  return *this;
869  }
870 
871  Element *elements_old = GetMatrixArray();
872  const Int_t nelems_old = this->fNelems;
873  const Int_t nrows_old = this->fNrows;
874  const Int_t ncols_old = this->fNcols;
875  const Int_t rowLwb_old = this->fRowLwb;
876  const Int_t colLwb_old = this->fColLwb;
877 
878  Allocate(new_nrows,new_ncols,row_lwb,col_lwb);
879  R__ASSERT(this->IsValid());
880 
881  Element *elements_new = GetMatrixArray();
882  // new memory should be initialized but be careful ot to wipe out the stack
883  // storage. Initialize all when old or new storage was on the heap
884  if (this->fNelems > this->kSizeMax || nelems_old > this->kSizeMax)
885  memset(elements_new,0,this->fNelems*sizeof(Element));
886  else if (this->fNelems > nelems_old)
887  memset(elements_new+nelems_old,0,(this->fNelems-nelems_old)*sizeof(Element));
888 
889  // Copy overlap
890  const Int_t rowLwb_copy = TMath::Max(this->fRowLwb,rowLwb_old);
891  const Int_t colLwb_copy = TMath::Max(this->fColLwb,colLwb_old);
892  const Int_t rowUpb_copy = TMath::Min(this->fRowLwb+this->fNrows-1,rowLwb_old+nrows_old-1);
893  const Int_t colUpb_copy = TMath::Min(this->fColLwb+this->fNcols-1,colLwb_old+ncols_old-1);
894 
895  const Int_t nrows_copy = rowUpb_copy-rowLwb_copy+1;
896  const Int_t ncols_copy = colUpb_copy-colLwb_copy+1;
897 
898  if (nrows_copy > 0 && ncols_copy > 0) {
899  const Int_t colOldOff = colLwb_copy-colLwb_old;
900  const Int_t colNewOff = colLwb_copy-this->fColLwb;
901  if (ncols_old < this->fNcols) {
902  for (Int_t i = nrows_copy-1; i >= 0; i--) {
903  const Int_t iRowOld = rowLwb_copy+i-rowLwb_old;
904  const Int_t iRowNew = rowLwb_copy+i-this->fRowLwb;
905  Memcpy_m(elements_new+iRowNew*this->fNcols+colNewOff,
906  elements_old+iRowOld*ncols_old+colOldOff,ncols_copy,this->fNelems,nelems_old);
907  if (this->fNelems <= this->kSizeMax && nelems_old <= this->kSizeMax)
908  memset(elements_new+iRowNew*this->fNcols+colNewOff+ncols_copy,0,
909  (this->fNcols-ncols_copy)*sizeof(Element));
910  }
911  } else {
912  for (Int_t i = 0; i < nrows_copy; i++) {
913  const Int_t iRowOld = rowLwb_copy+i-rowLwb_old;
914  const Int_t iRowNew = rowLwb_copy+i-this->fRowLwb;
915  Memcpy_m(elements_new+iRowNew*this->fNcols+colNewOff,
916  elements_old+iRowOld*ncols_old+colOldOff,ncols_copy,this->fNelems,nelems_old);
917  }
918  }
919  }
920 
921  Delete_m(nelems_old,elements_old);
922  } else {
923  Allocate(new_nrows,new_ncols,row_lwb,col_lwb,1);
924  }
925 
926  return *this;
927 }
928 
929 ////////////////////////////////////////////////////////////////////////////////
930 
931 template<class Element>
933 {
934  const TMatrixT<Element> &tmp = *this;
935  TDecompLU lu(tmp,this->fTol);
936  Double_t d1,d2;
937  lu.Det(d1,d2);
938  return d1*TMath::Power(2.0,d2);
939 }
940 
941 ////////////////////////////////////////////////////////////////////////////////
942 
943 template<class Element>
945 {
946  const TMatrixT<Element> &tmp = *this;
947  TDecompLU lu(tmp,this->fTol);
948  lu.Det(d1,d2);
949 }
950 
951 ////////////////////////////////////////////////////////////////////////////////
952 /// Invert the matrix and calculate its determinant
953 /// Notice that the LU decomposition is used instead of Bunch-Kaufman
954 /// Bunch-Kaufman guarantees a symmetric inverted matrix but is slower than LU .
955 /// The user can access Bunch-Kaufman through the TDecompBK class .
956 
957 template<class Element>
959 {
960  R__ASSERT(this->IsValid());
961  TMatrixD tmp(*this);
962  if (TDecompLU::InvertLU(tmp,Double_t(this->fTol),det)) {
963  const Double_t *p1 = tmp.GetMatrixArray();
964  Element *p2 = this->GetMatrixArray();
965  for (Int_t i = 0; i < this->GetNoElements(); i++)
966  p2[i] = p1[i];
967  }
968 
969  return *this;
970 }
971 
972 ////////////////////////////////////////////////////////////////////////////////
973 /// Invert the matrix and calculate its determinant
974 
975 template<class Element>
977 {
978  R__ASSERT(this->IsValid());
979 
980  const Char_t nRows = Char_t(this->GetNrows());
981  switch (nRows) {
982  case 1:
983  {
984  Element *pM = this->GetMatrixArray();
985  if (*pM == 0.) {
986  Error("InvertFast","matrix is singular");
987  *det = 0;
988  } else {
989  *det = *pM;
990  *pM = 1.0/(*pM);
991  }
992  return *this;
993  }
994  case 2:
995  {
996  TMatrixTSymCramerInv::Inv2x2<Element>(*this,det);
997  return *this;
998  }
999  case 3:
1000  {
1001  TMatrixTSymCramerInv::Inv3x3<Element>(*this,det);
1002  return *this;
1003  }
1004  case 4:
1005  {
1006  TMatrixTSymCramerInv::Inv4x4<Element>(*this,det);
1007  return *this;
1008  }
1009  case 5:
1010  {
1011  TMatrixTSymCramerInv::Inv5x5<Element>(*this,det);
1012  return *this;
1013  }
1014  case 6:
1015  {
1016  TMatrixTSymCramerInv::Inv6x6<Element>(*this,det);
1017  return *this;
1018  }
1019 
1020  default:
1021  {
1022  TMatrixD tmp(*this);
1023  if (TDecompLU::InvertLU(tmp,Double_t(this->fTol),det)) {
1024  const Double_t *p1 = tmp.GetMatrixArray();
1025  Element *p2 = this->GetMatrixArray();
1026  for (Int_t i = 0; i < this->GetNoElements(); i++)
1027  p2[i] = p1[i];
1028  }
1029  return *this;
1030  }
1031  }
1032 }
1033 
1034 ////////////////////////////////////////////////////////////////////////////////
1035 /// Transpose a matrix.
1036 
1037 template<class Element>
1039 {
1040  if (gMatrixCheck) {
1041  R__ASSERT(this->IsValid());
1042  R__ASSERT(source.IsValid());
1043 
1044  if (this->fNrows != source.GetNcols() || this->fRowLwb != source.GetColLwb())
1045  {
1046  Error("Transpose","matrix has wrong shape");
1047  return *this;
1048  }
1049  }
1050 
1051  *this = source;
1052  return *this;
1053 }
1054 
1055 ////////////////////////////////////////////////////////////////////////////////
1056 /// Perform a rank 1 operation on the matrix:
1057 /// A += alpha * v * v^T
1058 
1059 template<class Element>
1061 {
1062  if (gMatrixCheck) {
1063  R__ASSERT(this->IsValid());
1064  R__ASSERT(v.IsValid());
1065  if (v.GetNoElements() < this->fNrows) {
1066  Error("Rank1Update","vector too short");
1067  return *this;
1068  }
1069  }
1070 
1071  const Element * const pv = v.GetMatrixArray();
1072  Element *trp = this->GetMatrixArray(); // pointer to UR part and diagonal, traverse row-wise
1073  Element *tcp = trp; // pointer to LL part, traverse col-wise
1074  for (Int_t i = 0; i < this->fNrows; i++) {
1075  trp += i; // point to [i,i]
1076  tcp += i*this->fNcols; // point to [i,i]
1077  const Element tmp = alpha*pv[i];
1078  for (Int_t j = i; j < this->fNcols; j++) {
1079  if (j > i) *tcp += tmp*pv[j];
1080  *trp++ += tmp*pv[j];
1081  tcp += this->fNcols;
1082  }
1083  tcp -= this->fNelems-1; // point to [0,i]
1084  }
1085 
1086  return *this;
1087 }
1088 
1089 ////////////////////////////////////////////////////////////////////////////////
1090 /// Calculate B * (*this) * B^T , final matrix will be (nrowsb x nrowsb)
1091 /// This is a similarity transform when B is orthogonal . It is more
1092 /// efficient than applying the actual multiplication because this
1093 /// routine realizes that the final matrix is symmetric .
1094 
1095 template<class Element>
1097 {
1098  if (gMatrixCheck) {
1099  R__ASSERT(this->IsValid());
1100  R__ASSERT(b.IsValid());
1101  if (this->fNcols != b.GetNcols() || this->fColLwb != b.GetColLwb()) {
1102  Error("Similarity(const TMatrixT &)","matrices incompatible");
1103  return *this;
1104  }
1105  }
1106 
1107  const Int_t ncolsa = this->fNcols;
1108  const Int_t nb = b.GetNoElements();
1109  const Int_t nrowsb = b.GetNrows();
1110  const Int_t ncolsb = b.GetNcols();
1111 
1112  const Element * const bp = b.GetMatrixArray();
1113 
1114  Element work[kWorkMax];
1115  Bool_t isAllocated = kFALSE;
1116  Element *bap = work;
1117  if (nrowsb*ncolsa > kWorkMax) {
1118  isAllocated = kTRUE;
1119  bap = new Element[nrowsb*ncolsa];
1120  }
1121 
1122  AMultB(bp,nb,ncolsb,this->fElements,this->fNelems,this->fNcols,bap);
1123 
1124  if (nrowsb != this->fNrows)
1125  this->ResizeTo(nrowsb,nrowsb);
1126 
1127 #ifdef CBLAS
1128  Element *cp = this->GetMatrixArray();
1129  if (typeid(Element) == typeid(Double_t))
1130  cblas_dgemm (CblasRowMajor,CblasNoTrans,CblasTrans,this->fNrows,this->fNcols,ba.GetNcols(),
1131  1.0,bap,ba.GetNcols(),bp,b.GetNcols(),1.0,cp,this->fNcols);
1132  else if (typeid(Element) != typeid(Float_t))
1133  cblas_sgemm (CblasRowMajor,CblasNoTrans,CblasTrans,this->fNrows,this->fNcols,ba.GetNcols(),
1134  1.0,bap,ba.GetNcols(),bp,b.GetNcols(),1.0,cp,this->fNcols);
1135  else
1136  Error("Similarity","type %s not implemented in BLAS library",typeid(Element));
1137 #else
1138  const Int_t nba = nrowsb*ncolsa;
1139  const Int_t ncolsba = ncolsa;
1140  const Element * bi1p = bp;
1141  Element * cp = this->GetMatrixArray();
1142  Element * const cp0 = cp;
1143 
1144  Int_t ishift = 0;
1145  const Element *barp0 = bap;
1146  while (barp0 < bap+nba) {
1147  const Element *brp0 = bi1p;
1148  while (brp0 < bp+nb) {
1149  const Element *barp = barp0;
1150  const Element *brp = brp0;
1151  Element cij = 0;
1152  while (brp < brp0+ncolsb)
1153  cij += *barp++ * *brp++;
1154  *cp++ = cij;
1155  brp0 += ncolsb;
1156  }
1157  barp0 += ncolsba;
1158  bi1p += ncolsb;
1159  cp += ++ishift;
1160  }
1161 
1162  R__ASSERT(cp == cp0+this->fNelems+ishift && barp0 == bap+nba);
1163 
1164  cp = cp0;
1165  for (Int_t irow = 0; irow < this->fNrows; irow++) {
1166  const Int_t rowOff1 = irow*this->fNrows;
1167  for (Int_t icol = 0; icol < irow; icol++) {
1168  const Int_t rowOff2 = icol*this->fNrows;
1169  cp[rowOff1+icol] = cp[rowOff2+irow];
1170  }
1171  }
1172 #endif
1173 
1174  if (isAllocated)
1175  delete [] bap;
1176 
1177  return *this;
1178 }
1179 
1180 ////////////////////////////////////////////////////////////////////////////////
1181 /// Calculate B * (*this) * B^T , final matrix will be (nrowsb x nrowsb)
1182 /// This is a similarity transform when B is orthogonal . It is more
1183 /// efficient than applying the actual multiplication because this
1184 /// routine realizes that the final matrix is symmetric .
1185 
1186 template<class Element>
1188 {
1189  if (gMatrixCheck) {
1190  R__ASSERT(this->IsValid());
1191  R__ASSERT(b.IsValid());
1192  if (this->fNcols != b.GetNcols() || this->fColLwb != b.GetColLwb()) {
1193  Error("Similarity(const TMatrixTSym &)","matrices incompatible");
1194  return *this;
1195  }
1196  }
1197 
1198 #ifdef CBLAS
1199  const Int_t nrowsb = b.GetNrows();
1200  const Int_t ncolsa = this->GetNcols();
1201 
1202  Element work[kWorkMax];
1203  Bool_t isAllocated = kFALSE;
1204  Element *abtp = work;
1205  if (this->fNcols > kWorkMax) {
1206  isAllocated = kTRUE;
1207  abtp = new Element[this->fNcols];
1208  }
1209 
1211 
1212  const Element *bp = b.GetMatrixArray();
1213  Element *cp = this->GetMatrixArray();
1214  if (typeid(Element) == typeid(Double_t))
1215  cblas_dsymm (CblasRowMajor,CblasLeft,CblasUpper,this->fNrows,this->fNcols,1.0,
1216  bp,b.GetNcols(),abtp,abt.GetNcols(),0.0,cp,this->fNcols);
1217  else if (typeid(Element) != typeid(Float_t))
1218  cblas_ssymm (CblasRowMajor,CblasLeft,CblasUpper,this->fNrows,this->fNcols,1.0,
1219  bp,b.GetNcols(),abtp,abt.GetNcols(),0.0,cp,this->fNcols);
1220  else
1221  Error("Similarity","type %s not implemented in BLAS library",typeid(Element));
1222 
1223  if (isAllocated)
1224  delete [] abtp;
1225 #else
1226  const Int_t ncolsa = this->GetNcols();
1227  const Int_t nb = b.GetNoElements();
1228  const Int_t nrowsb = b.GetNrows();
1229  const Int_t ncolsb = b.GetNcols();
1230 
1231  const Element * const bp = b.GetMatrixArray();
1232 
1233  Element work[kWorkMax];
1234  Bool_t isAllocated = kFALSE;
1235  Element *bap = work;
1236  if (nrowsb*ncolsa > kWorkMax) {
1237  isAllocated = kTRUE;
1238  bap = new Element[nrowsb*ncolsa];
1239  }
1240 
1241  AMultB(bp,nb,ncolsb,this->fElements,this->fNelems,this->fNcols,bap);
1242 
1243  const Int_t nba = nrowsb*ncolsa;
1244  const Int_t ncolsba = ncolsa;
1245  const Element * bi1p = bp;
1246  Element * cp = this->GetMatrixArray();
1247  Element * const cp0 = cp;
1248 
1249  Int_t ishift = 0;
1250  const Element *barp0 = bap;
1251  while (barp0 < bap+nba) {
1252  const Element *brp0 = bi1p;
1253  while (brp0 < bp+nb) {
1254  const Element *barp = barp0;
1255  const Element *brp = brp0;
1256  Element cij = 0;
1257  while (brp < brp0+ncolsb)
1258  cij += *barp++ * *brp++;
1259  *cp++ = cij;
1260  brp0 += ncolsb;
1261  }
1262  barp0 += ncolsba;
1263  bi1p += ncolsb;
1264  cp += ++ishift;
1265  }
1266 
1267  R__ASSERT(cp == cp0+this->fNelems+ishift && barp0 == bap+nba);
1268 
1269  cp = cp0;
1270  for (Int_t irow = 0; irow < this->fNrows; irow++) {
1271  const Int_t rowOff1 = irow*this->fNrows;
1272  for (Int_t icol = 0; icol < irow; icol++) {
1273  const Int_t rowOff2 = icol*this->fNrows;
1274  cp[rowOff1+icol] = cp[rowOff2+irow];
1275  }
1276  }
1277 
1278  if (isAllocated)
1279  delete [] bap;
1280 #endif
1281 
1282  return *this;
1283 }
1284 
1285 ////////////////////////////////////////////////////////////////////////////////
1286 /// Calculate scalar v * (*this) * v^T
1287 
1288 template<class Element>
1290 {
1291  if (gMatrixCheck) {
1292  R__ASSERT(this->IsValid());
1293  R__ASSERT(v.IsValid());
1294  if (this->fNcols != v.GetNrows() || this->fColLwb != v.GetLwb()) {
1295  Error("Similarity(const TVectorT &)","vector and matrix incompatible");
1296  return -1.;
1297  }
1298  }
1299 
1300  const Element *mp = this->GetMatrixArray(); // Matrix row ptr
1301  const Element *vp = v.GetMatrixArray(); // vector ptr
1302 
1303  Element sum1 = 0;
1304  const Element * const vp_first = vp;
1305  const Element * const vp_last = vp+v.GetNrows();
1306  while (vp < vp_last) {
1307  Element sum2 = 0;
1308  for (const Element *sp = vp_first; sp < vp_last; )
1309  sum2 += *mp++ * *sp++;
1310  sum1 += sum2 * *vp++;
1311  }
1312 
1313  R__ASSERT(mp == this->GetMatrixArray()+this->GetNoElements());
1314 
1315  return sum1;
1316 }
1317 
1318 ////////////////////////////////////////////////////////////////////////////////
1319 /// Calculate B^T * (*this) * B , final matrix will be (ncolsb x ncolsb)
1320 /// It is more efficient than applying the actual multiplication because this
1321 /// routine realizes that the final matrix is symmetric .
1322 
1323 template<class Element>
1325 {
1326  if (gMatrixCheck) {
1327  R__ASSERT(this->IsValid());
1328  R__ASSERT(b.IsValid());
1329  if (this->fNrows != b.GetNrows() || this->fRowLwb != b.GetRowLwb()) {
1330  Error("SimilarityT(const TMatrixT &)","matrices incompatible");
1331  return *this;
1332  }
1333  }
1334 
1335  const Int_t ncolsb = b.GetNcols();
1336  const Int_t ncolsa = this->GetNcols();
1337 
1338  Element work[kWorkMax];
1339  Bool_t isAllocated = kFALSE;
1340  Element *btap = work;
1341  if (ncolsb*ncolsa > kWorkMax) {
1342  isAllocated = kTRUE;
1343  btap = new Element[ncolsb*ncolsa];
1344  }
1345 
1346  TMatrixT<Element> bta; bta.Use(ncolsb,ncolsa,btap);
1347  bta.TMult(b,*this);
1348 
1349  if (ncolsb != this->fNcols)
1350  this->ResizeTo(ncolsb,ncolsb);
1351 
1352 #ifdef CBLAS
1353  const Element *bp = b.GetMatrixArray();
1354  Element *cp = this->GetMatrixArray();
1355  if (typeid(Element) == typeid(Double_t))
1356  cblas_dgemm (CblasRowMajor,CblasNoTrans,CblasNoTrans,this->fNrows,this->fNcols,bta.GetNcols(),
1357  1.0,btap,bta.GetNcols(),bp,b.GetNcols(),1.0,cp,this->fNcols);
1358  else if (typeid(Element) != typeid(Float_t))
1359  cblas_sgemm (CblasRowMajor,CblasNoTrans,CblasNoTrans,this->fNrows,this->fNcols,bta.GetNcols(),
1360  1.0,btap,bta.GetNcols(),bp,b.GetNcols(),1.0,cp,this->fNcols);
1361  else
1362  Error("similarityT","type %s not implemented in BLAS library",typeid(Element));
1363 #else
1364  const Int_t nbta = bta.GetNoElements();
1365  const Int_t nb = b.GetNoElements();
1366  const Int_t ncolsbta = bta.GetNcols();
1367  const Element * const bp = b.GetMatrixArray();
1368  Element * cp = this->GetMatrixArray();
1369  Element * const cp0 = cp;
1370 
1371  Int_t ishift = 0;
1372  const Element *btarp0 = btap; // Pointer to A[i,0];
1373  const Element *bcp0 = bp;
1374  while (btarp0 < btap+nbta) {
1375  for (const Element *bcp = bcp0; bcp < bp+ncolsb; ) { // Pointer to the j-th column of B, Start bcp = B[0,0]
1376  const Element *btarp = btarp0; // Pointer to the i-th row of A, reset to A[i,0]
1377  Element cij = 0;
1378  while (bcp < bp+nb) { // Scan the i-th row of A and
1379  cij += *btarp++ * *bcp; // the j-th col of B
1380  bcp += ncolsb;
1381  }
1382  *cp++ = cij;
1383  bcp -= nb-1; // Set bcp to the (j+1)-th col
1384  }
1385  btarp0 += ncolsbta; // Set ap to the (i+1)-th row
1386  bcp0++;
1387  cp += ++ishift;
1388  }
1389 
1390  R__ASSERT(cp == cp0+this->fNelems+ishift && btarp0 == btap+nbta);
1391 
1392  cp = cp0;
1393  for (Int_t irow = 0; irow < this->fNrows; irow++) {
1394  const Int_t rowOff1 = irow*this->fNrows;
1395  for (Int_t icol = 0; icol < irow; icol++) {
1396  const Int_t rowOff2 = icol*this->fNrows;
1397  cp[rowOff1+icol] = cp[rowOff2+irow];
1398  }
1399  }
1400 #endif
1401 
1402  if (isAllocated)
1403  delete [] btap;
1404 
1405  return *this;
1406 }
1407 
1408 ////////////////////////////////////////////////////////////////////////////////
1409 
1410 template<class Element>
1412 {
1413  if (gMatrixCheck && !AreCompatible(*this,source)) {
1414  Error("operator=","matrices not compatible");
1415  return *this;
1416  }
1417 
1418  if (this->GetMatrixArray() != source.GetMatrixArray()) {
1419  TObject::operator=(source);
1420  memcpy(this->GetMatrixArray(),source.fElements,this->fNelems*sizeof(Element));
1421  }
1422  return *this;
1423 }
1424 
1425 ////////////////////////////////////////////////////////////////////////////////
1426 
1427 template<class Element>
1429 {
1430  R__ASSERT(this->IsValid());
1431 
1432  if (lazy_constructor.fRowUpb != this->GetRowUpb() ||
1433  lazy_constructor.fRowLwb != this->GetRowLwb()) {
1434  Error("operator=(const TMatrixTSymLazy&)", "matrix is incompatible with "
1435  "the assigned Lazy matrix");
1436  return *this;
1437  }
1438 
1439  lazy_constructor.FillIn(*this);
1440  return *this;
1441 }
1442 
1443 ////////////////////////////////////////////////////////////////////////////////
1444 /// Assign val to every element of the matrix.
1445 
1446 template<class Element>
1448 {
1449  R__ASSERT(this->IsValid());
1450 
1451  Element *ep = fElements;
1452  const Element * const ep_last = ep+this->fNelems;
1453  while (ep < ep_last)
1454  *ep++ = val;
1455 
1456  return *this;
1457 }
1458 
1459 ////////////////////////////////////////////////////////////////////////////////
1460 /// Add val to every element of the matrix.
1461 
1462 template<class Element>
1464 {
1465  R__ASSERT(this->IsValid());
1466 
1467  Element *ep = fElements;
1468  const Element * const ep_last = ep+this->fNelems;
1469  while (ep < ep_last)
1470  *ep++ += val;
1471 
1472  return *this;
1473 }
1474 
1475 ////////////////////////////////////////////////////////////////////////////////
1476 /// Subtract val from every element of the matrix.
1477 
1478 template<class Element>
1480 {
1481  R__ASSERT(this->IsValid());
1482 
1483  Element *ep = fElements;
1484  const Element * const ep_last = ep+this->fNelems;
1485  while (ep < ep_last)
1486  *ep++ -= val;
1487 
1488  return *this;
1489 }
1490 
1491 ////////////////////////////////////////////////////////////////////////////////
1492 /// Multiply every element of the matrix with val.
1493 
1494 template<class Element>
1496 {
1497  R__ASSERT(this->IsValid());
1498 
1499  Element *ep = fElements;
1500  const Element * const ep_last = ep+this->fNelems;
1501  while (ep < ep_last)
1502  *ep++ *= val;
1503 
1504  return *this;
1505 }
1506 
1507 ////////////////////////////////////////////////////////////////////////////////
1508 /// Add the source matrix.
1509 
1510 template<class Element>
1512 {
1513  if (gMatrixCheck && !AreCompatible(*this,source)) {
1514  Error("operator+=","matrices not compatible");
1515  return *this;
1516  }
1517 
1518  const Element *sp = source.GetMatrixArray();
1519  Element *tp = this->GetMatrixArray();
1520  const Element * const tp_last = tp+this->fNelems;
1521  while (tp < tp_last)
1522  *tp++ += *sp++;
1523 
1524  return *this;
1525 }
1526 
1527 ////////////////////////////////////////////////////////////////////////////////
1528 /// Subtract the source matrix.
1529 
1530 template<class Element>
1532 {
1533  if (gMatrixCheck && !AreCompatible(*this,source)) {
1534  Error("operator-=","matrices not compatible");
1535  return *this;
1536  }
1537 
1538  const Element *sp = source.GetMatrixArray();
1539  Element *tp = this->GetMatrixArray();
1540  const Element * const tp_last = tp+this->fNelems;
1541  while (tp < tp_last)
1542  *tp++ -= *sp++;
1543 
1544  return *this;
1545 }
1546 
1547 ////////////////////////////////////////////////////////////////////////////////
1548 
1549 template<class Element>
1551 {
1552  R__ASSERT(this->IsValid());
1553 
1554  Element val = 0;
1555  Element *trp = this->GetMatrixArray(); // pointer to UR part and diagonal, traverse row-wise
1556  Element *tcp = trp; // pointer to LL part, traverse col-wise
1557  for (Int_t i = 0; i < this->fNrows; i++) {
1558  trp += i; // point to [i,i]
1559  tcp += i*this->fNcols; // point to [i,i]
1560  for (Int_t j = i; j < this->fNcols; j++) {
1561  action.Operation(val);
1562  if (j > i) *tcp = val;
1563  *trp++ = val;
1564  tcp += this->fNcols;
1565  }
1566  tcp -= this->fNelems-1; // point to [0,i]
1567  }
1568 
1569  return *this;
1570 }
1571 
1572 ////////////////////////////////////////////////////////////////////////////////
1573 /// Apply action to each element of the matrix. To action the location
1574 /// of the current element is passed.
1575 
1576 template<class Element>
1578 {
1579  R__ASSERT(this->IsValid());
1580 
1581  Element val = 0;
1582  Element *trp = this->GetMatrixArray(); // pointer to UR part and diagonal, traverse row-wise
1583  Element *tcp = trp; // pointer to LL part, traverse col-wise
1584  for (Int_t i = 0; i < this->fNrows; i++) {
1585  action.fI = i+this->fRowLwb;
1586  trp += i; // point to [i,i]
1587  tcp += i*this->fNcols; // point to [i,i]
1588  for (Int_t j = i; j < this->fNcols; j++) {
1589  action.fJ = j+this->fColLwb;
1590  action.Operation(val);
1591  if (j > i) *tcp = val;
1592  *trp++ = val;
1593  tcp += this->fNcols;
1594  }
1595  tcp -= this->fNelems-1; // point to [0,i]
1596  }
1597 
1598  return *this;
1599 }
1600 
1601 ////////////////////////////////////////////////////////////////////////////////
1602 /// randomize matrix element values but keep matrix symmetric
1603 
1604 template<class Element>
1606 {
1607  if (gMatrixCheck) {
1608  R__ASSERT(this->IsValid());
1609  if (this->fNrows != this->fNcols || this->fRowLwb != this->fColLwb) {
1610  Error("Randomize(Element,Element,Element &","matrix should be square");
1611  return *this;
1612  }
1613  }
1614 
1615  const Element scale = beta-alpha;
1616  const Element shift = alpha/scale;
1617 
1618  Element *ep = GetMatrixArray();
1619  for (Int_t i = 0; i < this->fNrows; i++) {
1620  const Int_t off = i*this->fNcols;
1621  for (Int_t j = 0; j <= i; j++) {
1622  ep[off+j] = scale*(Drand(seed)+shift);
1623  if (i != j) {
1624  ep[j*this->fNcols+i] = ep[off+j];
1625  }
1626  }
1627  }
1628 
1629  return *this;
1630 }
1631 
1632 ////////////////////////////////////////////////////////////////////////////////
1633 /// randomize matrix element values but keep matrix symmetric positive definite
1634 
1635 template<class Element>
1637 {
1638  if (gMatrixCheck) {
1639  R__ASSERT(this->IsValid());
1640  if (this->fNrows != this->fNcols || this->fRowLwb != this->fColLwb) {
1641  Error("RandomizeSym(Element,Element,Element &","matrix should be square");
1642  return *this;
1643  }
1644  }
1645 
1646  const Element scale = beta-alpha;
1647  const Element shift = alpha/scale;
1648 
1649  Element *ep = GetMatrixArray();
1650  Int_t i;
1651  for (i = 0; i < this->fNrows; i++) {
1652  const Int_t off = i*this->fNcols;
1653  for (Int_t j = 0; j <= i; j++)
1654  ep[off+j] = scale*(Drand(seed)+shift);
1655  }
1656 
1657  for (i = this->fNrows-1; i >= 0; i--) {
1658  const Int_t off1 = i*this->fNcols;
1659  for (Int_t j = i; j >= 0; j--) {
1660  const Int_t off2 = j*this->fNcols;
1661  ep[off1+j] *= ep[off2+j];
1662  for (Int_t k = j-1; k >= 0; k--) {
1663  ep[off1+j] += ep[off1+k]*ep[off2+k];
1664  }
1665  if (i != j)
1666  ep[off2+i] = ep[off1+j];
1667  }
1668  }
1669 
1670  return *this;
1671 }
1672 
1673 ////////////////////////////////////////////////////////////////////////////////
1674 /// Return a matrix containing the eigen-vectors ordered by descending eigen-values.
1675 /// For full functionality use TMatrixDSymEigen .
1676 
1677 template<class Element>
1679 {
1680  TMatrixDSym tmp = *this;
1681  TMatrixDSymEigen eigen(tmp);
1682  eigenValues.ResizeTo(this->fNrows);
1683  eigenValues = eigen.GetEigenValues();
1684  return eigen.GetEigenVectors();
1685 }
1686 
1687 ////////////////////////////////////////////////////////////////////////////////
1688 /// Check to see if two matrices are identical.
1689 
1690 template<class Element>
1692 {
1693  if (!AreCompatible(m1,m2)) return kFALSE;
1694  return (memcmp(m1.GetMatrixArray(),m2.GetMatrixArray(),
1695  m1.GetNoElements()*sizeof(Element)) == 0);
1696 }
1697 
1698 ////////////////////////////////////////////////////////////////////////////////
1699 
1700 template<class Element>
1702 {
1703  TMatrixTSym<Element> target(source1);
1704  target += source2;
1705  return target;
1706 }
1707 
1708 ////////////////////////////////////////////////////////////////////////////////
1709 
1710 template<class Element>
1712 {
1713  TMatrixTSym<Element> target(source1);
1714  target += val;
1715  return target;
1716 }
1717 
1718 ////////////////////////////////////////////////////////////////////////////////
1719 
1720 template<class Element>
1722 {
1723  return operator+(source1,val);
1724 }
1725 
1726 ////////////////////////////////////////////////////////////////////////////////
1727 
1728 template<class Element>
1730 {
1731  TMatrixTSym<Element> target(source1);
1732  target -= source2;
1733  return target;
1734 }
1735 
1736 ////////////////////////////////////////////////////////////////////////////////
1737 
1738 template<class Element>
1740 {
1741  TMatrixTSym<Element> target(source1);
1742  target -= val;
1743  return target;
1744 }
1745 
1746 ////////////////////////////////////////////////////////////////////////////////
1747 
1748 template<class Element>
1750 {
1751  return Element(-1.0)*operator-(source1,val);
1752 }
1753 
1754 ////////////////////////////////////////////////////////////////////////////////
1755 
1756 template<class Element>
1758 {
1759  TMatrixTSym<Element> target(source1);
1760  target *= val;
1761  return target;
1762 }
1763 
1764 ////////////////////////////////////////////////////////////////////////////////
1765 
1766 template<class Element>
1768 {
1769  return operator*(source1,val);
1770 }
1771 
1772 ////////////////////////////////////////////////////////////////////////////////
1773 /// Logical AND
1774 
1775 template<class Element>
1777 {
1778  TMatrixTSym<Element> target;
1779 
1780  if (gMatrixCheck && !AreCompatible(source1,source2)) {
1781  Error("operator&&(const TMatrixTSym&,const TMatrixTSym&)","matrices not compatible");
1782  return target;
1783  }
1784 
1785  target.ResizeTo(source1);
1786 
1787  const Element *sp1 = source1.GetMatrixArray();
1788  const Element *sp2 = source2.GetMatrixArray();
1789  Element *tp = target.GetMatrixArray();
1790  const Element * const tp_last = tp+target.GetNoElements();
1791  while (tp < tp_last)
1792  *tp++ = (*sp1++ != 0.0 && *sp2++ != 0.0);
1793 
1794  return target;
1795 }
1796 
1797 ////////////////////////////////////////////////////////////////////////////////
1798 /// Logical Or
1799 
1800 template<class Element>
1802 {
1803  TMatrixTSym<Element> target;
1804 
1805  if (gMatrixCheck && !AreCompatible(source1,source2)) {
1806  Error("operator||(const TMatrixTSym&,const TMatrixTSym&)","matrices not compatible");
1807  return target;
1808  }
1809 
1810  target.ResizeTo(source1);
1811 
1812  const Element *sp1 = source1.GetMatrixArray();
1813  const Element *sp2 = source2.GetMatrixArray();
1814  Element *tp = target.GetMatrixArray();
1815  const Element * const tp_last = tp+target.GetNoElements();
1816  while (tp < tp_last)
1817  *tp++ = (*sp1++ != 0.0 || *sp2++ != 0.0);
1818 
1819  return target;
1820 }
1821 
1822 ////////////////////////////////////////////////////////////////////////////////
1823 /// source1 > source2
1824 
1825 template<class Element>
1827 {
1828  TMatrixTSym<Element> target;
1829 
1830  if (gMatrixCheck && !AreCompatible(source1,source2)) {
1831  Error("operator>(const TMatrixTSym&,const TMatrixTSym&)","matrices not compatible");
1832  return target;
1833  }
1834 
1835  target.ResizeTo(source1);
1836 
1837  const Element *sp1 = source1.GetMatrixArray();
1838  const Element *sp2 = source2.GetMatrixArray();
1839  Element *tp = target.GetMatrixArray();
1840  const Element * const tp_last = tp+target.GetNoElements();
1841  while (tp < tp_last) {
1842  *tp++ = (*sp1) > (*sp2); sp1++; sp2++;
1843  }
1844 
1845  return target;
1846 }
1847 
1848 ////////////////////////////////////////////////////////////////////////////////
1849 /// source1 >= source2
1850 
1851 template<class Element>
1853 {
1854  TMatrixTSym<Element> target;
1855 
1856  if (gMatrixCheck && !AreCompatible(source1,source2)) {
1857  Error("operator>=(const TMatrixTSym&,const TMatrixTSym&)","matrices not compatible");
1858  return target;
1859  }
1860 
1861  target.ResizeTo(source1);
1862 
1863  const Element *sp1 = source1.GetMatrixArray();
1864  const Element *sp2 = source2.GetMatrixArray();
1865  Element *tp = target.GetMatrixArray();
1866  const Element * const tp_last = tp+target.GetNoElements();
1867  while (tp < tp_last) {
1868  *tp++ = (*sp1) >= (*sp2); sp1++; sp2++;
1869  }
1870 
1871  return target;
1872 }
1873 
1874 ////////////////////////////////////////////////////////////////////////////////
1875 /// source1 <= source2
1876 
1877 template<class Element>
1878 TMatrixTSym<Element> operator<=(const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2)
1879 {
1880  TMatrixTSym<Element> target;
1881 
1882  if (gMatrixCheck && !AreCompatible(source1,source2)) {
1883  Error("operator<=(const TMatrixTSym&,const TMatrixTSym&)","matrices not compatible");
1884  return target;
1885  }
1886 
1887  target.ResizeTo(source1);
1888 
1889  const Element *sp1 = source1.GetMatrixArray();
1890  const Element *sp2 = source2.GetMatrixArray();
1891  Element *tp = target.GetMatrixArray();
1892  const Element * const tp_last = tp+target.GetNoElements();
1893  while (tp < tp_last) {
1894  *tp++ = (*sp1) <= (*sp2); sp1++; sp2++;
1895  }
1896 
1897  return target;
1898 }
1899 
1900 ////////////////////////////////////////////////////////////////////////////////
1901 /// source1 < source2
1902 
1903 template<class Element>
1904 TMatrixTSym<Element> operator<(const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2)
1905 {
1906  TMatrixTSym<Element> target;
1907 
1908  if (gMatrixCheck && !AreCompatible(source1,source2)) {
1909  Error("operator<(const TMatrixTSym&,const TMatrixTSym&)","matrices not compatible");
1910  return target;
1911  }
1912 
1913  target.ResizeTo(source1);
1914 
1915  const Element *sp1 = source1.GetMatrixArray();
1916  const Element *sp2 = source2.GetMatrixArray();
1917  Element *tp = target.GetMatrixArray();
1918  const Element * const tp_last = tp+target.GetNoElements();
1919  while (tp < tp_last) {
1920  *tp++ = (*sp1) < (*sp2); sp1++; sp2++;
1921  }
1922 
1923  return target;
1924 }
1925 
1926 ////////////////////////////////////////////////////////////////////////////////
1927 /// Modify addition: target += scalar * source.
1928 
1929 template<class Element>
1931 {
1932  if (gMatrixCheck && !AreCompatible(target,source)) {
1933  ::Error("Add","matrices not compatible");
1934  return target;
1935  }
1936 
1937  const Int_t nrows = target.GetNrows();
1938  const Int_t ncols = target.GetNcols();
1939  const Int_t nelems = target.GetNoElements();
1940  const Element *sp = source.GetMatrixArray();
1941  Element *trp = target.GetMatrixArray(); // pointer to UR part and diagonal, traverse row-wise
1942  Element *tcp = target.GetMatrixArray(); // pointer to LL part, traverse col-wise
1943  for (Int_t i = 0; i < nrows; i++) {
1944  sp += i;
1945  trp += i; // point to [i,i]
1946  tcp += i*ncols; // point to [i,i]
1947  for (Int_t j = i; j < ncols; j++) {
1948  const Element tmp = scalar * *sp++;
1949  if (j > i) *tcp += tmp;
1950  *trp++ += tmp;
1951  tcp += ncols;
1952  }
1953  tcp -= nelems-1; // point to [0,i]
1954  }
1955 
1956  return target;
1957 }
1958 
1959 ////////////////////////////////////////////////////////////////////////////////
1960 /// Multiply target by the source, element-by-element.
1961 
1962 template<class Element>
1964 {
1965  if (gMatrixCheck && !AreCompatible(target,source)) {
1966  ::Error("ElementMult","matrices not compatible");
1967  return target;
1968  }
1969 
1970  const Int_t nrows = target.GetNrows();
1971  const Int_t ncols = target.GetNcols();
1972  const Int_t nelems = target.GetNoElements();
1973  const Element *sp = source.GetMatrixArray();
1974  Element *trp = target.GetMatrixArray(); // pointer to UR part and diagonal, traverse row-wise
1975  Element *tcp = target.GetMatrixArray(); // pointer to LL part, traverse col-wise
1976  for (Int_t i = 0; i < nrows; i++) {
1977  sp += i;
1978  trp += i; // point to [i,i]
1979  tcp += i*ncols; // point to [i,i]
1980  for (Int_t j = i; j < ncols; j++) {
1981  if (j > i) *tcp *= *sp;
1982  *trp++ *= *sp++;
1983  tcp += ncols;
1984  }
1985  tcp -= nelems-1; // point to [0,i]
1986  }
1987 
1988  return target;
1989 }
1990 
1991 ////////////////////////////////////////////////////////////////////////////////
1992 /// Multiply target by the source, element-by-element.
1993 
1994 template<class Element>
1996 {
1997  if (gMatrixCheck && !AreCompatible(target,source)) {
1998  ::Error("ElementDiv","matrices not compatible");
1999  return target;
2000  }
2001 
2002  const Int_t nrows = target.GetNrows();
2003  const Int_t ncols = target.GetNcols();
2004  const Int_t nelems = target.GetNoElements();
2005  const Element *sp = source.GetMatrixArray();
2006  Element *trp = target.GetMatrixArray(); // pointer to UR part and diagonal, traverse row-wise
2007  Element *tcp = target.GetMatrixArray(); // pointer to LL part, traverse col-wise
2008  for (Int_t i = 0; i < nrows; i++) {
2009  sp += i;
2010  trp += i; // point to [i,i]
2011  tcp += i*ncols; // point to [i,i]
2012  for (Int_t j = i; j < ncols; j++) {
2013  if (*sp != 0.0) {
2014  if (j > i) *tcp /= *sp;
2015  *trp++ /= *sp++;
2016  } else {
2017  const Int_t irow = (sp-source.GetMatrixArray())/source.GetNcols();
2018  const Int_t icol = (sp-source.GetMatrixArray())%source.GetNcols();
2019  Error("ElementDiv","source (%d,%d) is zero",irow,icol);
2020  trp++;
2021  }
2022  tcp += ncols;
2023  }
2024  tcp -= nelems-1; // point to [0,i]
2025  }
2026 
2027  return target;
2028 }
2029 
2030 ////////////////////////////////////////////////////////////////////////////////
2031 /// Stream an object of class TMatrixTSym.
2032 
2033 template<class Element>
2035 {
2036  if (R__b.IsReading()) {
2037  UInt_t R__s, R__c;
2038  Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
2039  Clear();
2040  R__b.ReadClassBuffer(TMatrixTBase<Element>::Class(),this,R__v,R__s,R__c);
2041  fElements = new Element[this->fNelems];
2042  Int_t i;
2043  for (i = 0; i < this->fNrows; i++) {
2044  R__b.ReadFastArray(fElements+i*this->fNcols+i,this->fNcols-i);
2045  }
2046  // copy to Lower left triangle
2047  for (i = 0; i < this->fNrows; i++) {
2048  for (Int_t j = 0; j < i; j++) {
2049  fElements[i*this->fNcols+j] = fElements[j*this->fNrows+i];
2050  }
2051  }
2052  if (this->fNelems <= this->kSizeMax) {
2053  memcpy(fDataStack,fElements,this->fNelems*sizeof(Element));
2054  delete [] fElements;
2055  fElements = fDataStack;
2056  }
2057  } else {
2059  // Only write the Upper right triangle
2060  for (Int_t i = 0; i < this->fNrows; i++) {
2061  R__b.WriteFastArray(fElements+i*this->fNcols+i,this->fNcols-i);
2062  }
2063  }
2064 }
2065 
2066 #ifndef ROOT_TMatrixFSymfwd
2067 #include "TMatrixFSymfwd.h"
2068 #endif
2069 
2070 template class TMatrixTSym<Float_t>;
2071 
2072 template Bool_t operator== <Float_t>(const TMatrixFSym &source1,const TMatrixFSym &source2);
2073 template TMatrixFSym operator+ <Float_t>(const TMatrixFSym &source1,const TMatrixFSym &source2);
2074 template TMatrixFSym operator+ <Float_t>(const TMatrixFSym &source1, Float_t val);
2075 template TMatrixFSym operator+ <Float_t>( Float_t val ,const TMatrixFSym &source2);
2076 template TMatrixFSym operator- <Float_t>(const TMatrixFSym &source1,const TMatrixFSym &source2);
2077 template TMatrixFSym operator- <Float_t>(const TMatrixFSym &source1, Float_t val);
2078 template TMatrixFSym operator- <Float_t>( Float_t val ,const TMatrixFSym &source2);
2079 template TMatrixFSym operator* <Float_t>(const TMatrixFSym &source, Float_t val );
2080 template TMatrixFSym operator* <Float_t>( Float_t val, const TMatrixFSym &source );
2081 template TMatrixFSym operator&& <Float_t>(const TMatrixFSym &source1,const TMatrixFSym &source2);
2082 template TMatrixFSym operator|| <Float_t>(const TMatrixFSym &source1,const TMatrixFSym &source2);
2083 template TMatrixFSym operator> <Float_t>(const TMatrixFSym &source1,const TMatrixFSym &source2);
2084 template TMatrixFSym operator>= <Float_t>(const TMatrixFSym &source1,const TMatrixFSym &source2);
2085 template TMatrixFSym operator<= <Float_t>(const TMatrixFSym &source1,const TMatrixFSym &source2);
2086 template TMatrixFSym operator< <Float_t>(const TMatrixFSym &source1,const TMatrixFSym &source2);
2087 
2088 template TMatrixFSym &Add <Float_t>(TMatrixFSym &target, Float_t scalar,const TMatrixFSym &source);
2089 template TMatrixFSym &ElementMult<Float_t>(TMatrixFSym &target,const TMatrixFSym &source);
2090 template TMatrixFSym &ElementDiv <Float_t>(TMatrixFSym &target,const TMatrixFSym &source);
2091 
2092 #ifndef ROOT_TMatrixDSymfwd
2093 #include "TMatrixDSymfwd.h"
2094 #endif
2095 
2096 template class TMatrixTSym<Double_t>;
2097 
2098 template Bool_t operator== <Double_t>(const TMatrixDSym &source1,const TMatrixDSym &source2);
2099 template TMatrixDSym operator+ <Double_t>(const TMatrixDSym &source1,const TMatrixDSym &source2);
2100 template TMatrixDSym operator+ <Double_t>(const TMatrixDSym &source1, Double_t val);
2101 template TMatrixDSym operator+ <Double_t>( Double_t val ,const TMatrixDSym &source2);
2102 template TMatrixDSym operator- <Double_t>(const TMatrixDSym &source1,const TMatrixDSym &source2);
2103 template TMatrixDSym operator- <Double_t>(const TMatrixDSym &source1, Double_t val);
2104 template TMatrixDSym operator- <Double_t>( Double_t val ,const TMatrixDSym &source2);
2105 template TMatrixDSym operator* <Double_t>(const TMatrixDSym &source, Double_t val );
2106 template TMatrixDSym operator* <Double_t>( Double_t val, const TMatrixDSym &source );
2107 template TMatrixDSym operator&& <Double_t>(const TMatrixDSym &source1,const TMatrixDSym &source2);
2108 template TMatrixDSym operator|| <Double_t>(const TMatrixDSym &source1,const TMatrixDSym &source2);
2109 template TMatrixDSym operator> <Double_t>(const TMatrixDSym &source1,const TMatrixDSym &source2);
2110 template TMatrixDSym operator>= <Double_t>(const TMatrixDSym &source1,const TMatrixDSym &source2);
2111 template TMatrixDSym operator<= <Double_t>(const TMatrixDSym &source1,const TMatrixDSym &source2);
2112 template TMatrixDSym operator< <Double_t>(const TMatrixDSym &source1,const TMatrixDSym &source2);
2113 
2114 template TMatrixDSym &Add <Double_t>(TMatrixDSym &target, Double_t scalar,const TMatrixDSym &source);
2115 template TMatrixDSym &ElementMult<Double_t>(TMatrixDSym &target,const TMatrixDSym &source);
2116 template TMatrixDSym &ElementDiv <Double_t>(TMatrixDSym &target,const TMatrixDSym &source);
Int_t GetRowLwb() const
Definition: TMatrixTBase.h:132
TMatrixTSym< Element > & SimilarityT(const TMatrixT< Element > &n)
Calculate B^T * (*this) * B , final matrix will be (ncolsb x ncolsb) It is more efficient than applyi...
const TVectorD & GetEigenValues() const
virtual const Element * GetMatrixArray() const =0
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...
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...
template TMatrixFSym operator<=< Float_t >(const TMatrixFSym &source1, const TMatrixFSym &source2)
int Invert(LASymMatrix &)
Definition: LaInverse.cxx:22
TVectorT< Element > & ResizeTo(Int_t lwb, Int_t upb)
Resize the vector to [lwb:upb] .
Definition: TVectorT.cxx:291
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
static Vc_ALWAYS_INLINE int_v min(const int_v &x, const int_v &y)
Definition: vector.h:433
template TMatrixDSym operator<< Double_t >(const TMatrixDSym &source1, const TMatrixDSym &source2)
TMatrixTSym< Element > & ElementMult(TMatrixTSym< Element > &target, const TMatrixTSym< Element > &source)
Multiply target by the source, element-by-element.
virtual TMatrixTBase< Element > & ResizeTo(Int_t nrows, Int_t ncols, Int_t=-1)
Set size of the matrix to nrows x ncols New dynamic elements are created, the overlapping part of the...
TMatrixTSym< Element > & Use(Int_t row_lwb, Int_t row_upb, Element *data)
Small helper to encapsulate whether to return the value pointed to by the iterator or its address...
virtual TMatrixTBase< Element > & SetMatrixArray(const Element *data, Option_t *option="")
Copy array data to matrix .
Int_t GetRowLwb() const
Definition: TMatrixTLazy.h:112
Bool_t IsReading() const
Definition: TBuffer.h:83
short Version_t
Definition: RtypesCore.h:61
TMatrixTSym< Element > & Add(TMatrixTSym< Element > &target, Element scalar, const TMatrixTSym< Element > &source)
Modify addition: target += scalar * source.
Bool_t IsValid() const
Definition: TVectorT.h:89
void Allocate(Int_t nrows, Int_t ncols, Int_t row_lwb=0, Int_t col_lwb=0, Int_t init=0, Int_t=-1)
Allocate new matrix.
float Float_t
Definition: RtypesCore.h:53
Int_t GetRowUpb() const
Definition: TMatrixTLazy.h:113
TMatrixTSym< Element > & operator=(const TMatrixTSym< Element > &source)
const char Option_t
Definition: RtypesCore.h:62
templateClassImp(TMatrixTSym) template< class Element > TMatrixTSym< Element >
Definition: TMatrixTSym.cxx:33
Int_t GetNrows() const
Definition: TVectorT.h:81
const TMatrixT< Element > EigenVectors(TVectorT< Element > &eigenValues) const
Return a matrix containing the eigen-vectors ordered by descending eigen-values.
TMatrixTSym< Element > & ElementDiv(TMatrixTSym< Element > &target, const TMatrixTSym< Element > &source)
Multiply target by the source, element-by-element.
TMatrixTSym< Element > operator>=(const TMatrixTSym< Element > &source1, const TMatrixTSym< Element > &source2)
source1 >= source2
void ToUpper()
Change string to upper case.
Definition: TString.cxx:1088
Buffer base class used for serializing objects.
Definition: TBuffer.h:42
void TMult(const TMatrixT< Element > &a)
Create a matrix C such that C = A' * A.
#define R__ASSERT(e)
Definition: TError.h:98
Int_t GetNoElements() const
Definition: TMatrixTBase.h:138
Expr< TransposeOp< SMatrix< T, D1, D2, R >, T, D1, D2 >, T, D2, D1, typename TranspPolicy< T, D1, D2, R >::RepType > Transpose(const SMatrix< T, D1, D2, R > &rhs)
Matrix Transpose B(i,j) = A(j,i) returning a matrix expression.
static Bool_t InvertLU(TMatrixD &a, Double_t tol, Double_t *det=0)
Calculate matrix inversion through in place forward/backward substitution.
Definition: TDecompLU.cxx:775
Basic string class.
Definition: TString.h:137
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:170
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
virtual void Det(Double_t &d1, Double_t &d2)
Calculate determinant det = d1*TMath::Power(2.,d2)
Definition: TDecompLU.cxx:509
const Bool_t kFALSE
Definition: Rtypes.h:92
TMatrixTSym< Element > & operator*=(Element val)
Multiply every element of the matrix with val.
virtual TMatrixTBase< Element > & ResizeTo(Int_t nrows, Int_t ncols, Int_t nr_nonzeros=-1)=0
template TMatrixFSym operator<< Float_t >(const TMatrixFSym &source1, const TMatrixFSym &source2)
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
TMatrixTSym< Element > & Rank1Update(const TVectorT< Element > &v, Element alpha=1.0)
Perform a rank 1 operation on the matrix: A += alpha * v * v^T.
double beta(double x, double y)
Calculates the beta function.
virtual Bool_t IsSymmetric() const
Check whether matrix is symmetric.
Bool_t IsValid() const
Definition: TMatrixTBase.h:157
template TMatrixDSym & ElementDiv< Double_t >(TMatrixDSym &target, const TMatrixDSym &source)
TMatrixTSym< Element > & operator+=(Element val)
Add val to every element of the matrix.
virtual void Operation(Element &element) const =0
virtual Double_t Determinant() const
void AMultB(int n, int m, int k, const double *A, const double *B, double *C)
Definition: cblas.cxx:11
template TMatrixDSym & ElementMult< Double_t >(TMatrixDSym &target, const TMatrixDSym &source)
TMatrixTSym< Element > operator&&(const TMatrixTSym< Element > &source1, const TMatrixTSym< Element > &source2)
Logical AND.
static double p2(double t, double a, double b, double c)
TObject & operator=(const TObject &rhs)
TObject assignment operator.
Definition: TObject.cxx:102
Bool_t AreCompatible(const TMatrixTBase< Element1 > &m1, const TMatrixTBase< Element2 > &m2, Int_t verbose=0)
Check that matrice sm1 and m2 areboth valid and have identical shapes .
void TMult(const TMatrixT< Element > &a, const TMatrixT< Element > &b)
Create a matrix C such that C = A' * B.
Definition: TMatrixT.cxx:844
virtual TMatrixTBase< Element > & SetMatrixArray(const Element *data, Option_t *option="")
Copy array data to matrix .
virtual void FillIn(TMatrixTSym< Element > &m) const =0
void Error(const char *location, const char *msgfmt,...)
void Minus(const TMatrixTSym< Element > &a, const TMatrixTSym< Element > &b)
Symmetric matrix summation. Create a matrix C such that C = A + B.
virtual const Element * GetMatrixArray() const
Definition: TMatrixTSym.h:192
TMatrixTSym< Element > operator>(const TMatrixTSym< Element > &source1, const TMatrixTSym< Element > &source2)
source1 > source2
Int_t GetNrows() const
Definition: TMatrixTBase.h:134
Int_t Memcpy_m(Element *newp, const Element *oldp, Int_t copySize, Int_t newSize, Int_t oldSize)
copy copySize doubles from *oldp to *newp .
Element * GetMatrixArray()
Definition: TVectorT.h:84
TMatrixTBase< Element > & Apply(const TElementActionT< Element > &action)
Apply action to each matrix element.
SVector< double, 2 > v
Definition: Dict.h:5
template TMatrixFSym & Add< Float_t >(TMatrixFSym &target, Float_t scalar, const TMatrixFSym &source)
template TMatrixDSym operator<=< Double_t >(const TMatrixDSym &source1, const TMatrixDSym &source2)
unsigned int UInt_t
Definition: RtypesCore.h:42
TMarker * m
Definition: textangle.C:8
template TMatrixFSym & ElementDiv< Float_t >(TMatrixFSym &target, const TMatrixFSym &source)
virtual const Element * GetMatrixArray() const
Definition: TMatrixT.h:222
static double p1(double t, double a, double b)
Bool_t operator==(const TMatrixTSym< Element > &m1, const TMatrixTSym< Element > &m2)
Check to see if two matrices are identical.
R__EXTERN Int_t gMatrixCheck
Definition: TMatrixTBase.h:91
virtual TMatrixTBase< Element > & Randomize(Element alpha, Element beta, Double_t &seed)
randomize matrix element values but keep matrix symmetric
virtual void ReadFastArray(Bool_t *b, Int_t n)=0
Int_t GetColLwb() const
Definition: TMatrixTBase.h:135
virtual void WriteFastArray(const Bool_t *b, Int_t n)=0
REAL epsilon
Definition: triangle.c:617
TMatrixTSym< Element > & Transpose(const TMatrixTSym< Element > &source)
Transpose a matrix.
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
void Plus(const TMatrixTSym< Element > &a, const TMatrixTSym< Element > &b)
Symmetric matrix summation. Create a matrix C such that C = A + B.
Int_t GetLwb() const
Definition: TVectorT.h:79
static Int_t init()
double Double_t
Definition: RtypesCore.h:55
TMatrixTSym< Element > & operator-=(Element val)
Subtract val from every element of the matrix.
const TMatrixD & GetEigenVectors() const
Int_t GetNcols() const
Definition: TMatrixTBase.h:137
virtual TMatrixTSym< Element > & RandomizePD(Element alpha, Element beta, Double_t &seed)
randomize matrix element values but keep matrix symmetric positive definite
TMatrixTSym< Element > & InvertFast(Double_t *det=0)
Invert the matrix and calculate its determinant.
TMatrixTSym< Element > & SetSub(Int_t row_lwb, const TMatrixTBase< Element > &source)
Insert matrix source starting at [row_lwb][row_lwb], thereby overwriting the part [row_lwb...
TMatrixTSym< Element > & Similarity(const TMatrixT< Element > &n)
Calculate B * (*this) * B^T , final matrix will be (nrowsb x nrowsb) This is a similarity transform w...
virtual const Int_t * GetRowIndexArray() const
Definition: TMatrixTSym.h:88
void Delete_m(Int_t size, Element *&)
delete data pointer m, if it was assigned on the heap
virtual void Operation(Element &element) const =0
TCppObject_t Allocate(TCppType_t type)
Definition: Cppyy.cxx:228
TMatrixT< Element > & Use(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb, Element *data)
Use the array data to fill the matrix ([row_lwb..row_upb] x [col_lwb..col_upb])
Definition: TMatrixT.cxx:1045
TMatrixTSym< Element > & Invert(Double_t *det=0)
Invert the matrix and calculate its determinant Notice that the LU decomposition is used instead of B...
TMatrixTSym< Element > operator||(const TMatrixTSym< Element > &source1, const TMatrixTSym< Element > &source2)
Logical Or.
char Char_t
Definition: RtypesCore.h:29
Int_t GetNoElements() const
Definition: TVectorT.h:82
virtual const Int_t * GetColIndexArray() const
Definition: TMatrixTSym.h:90
TMatrixTSym< Element > operator+(const TMatrixTSym< Element > &source1, const TMatrixTSym< Element > &source2)
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
template TMatrixDSym & Add< Double_t >(TMatrixDSym &target, Double_t scalar, const TMatrixDSym &source)
Element * New_m(Int_t size)
return data pointer .
TMatrixTSym< Element > operator*(const TMatrixTSym< Element > &source1, Element val)
const Bool_t kTRUE
Definition: Rtypes.h:91
template TMatrixFSym & ElementMult< Float_t >(TMatrixFSym &target, const TMatrixFSym &source)
Element * fElements
data container
Definition: TMatrixTSym.h:43
Double_t Drand(Double_t &ix)
Random number generator [0....1] with seed ix.
TMatrixTSym< Element > operator-(const TMatrixTSym< Element > &source1, const TMatrixTSym< Element > &source2)
TMatrixTSym< Element > & GetSub(Int_t row_lwb, Int_t row_upb, TMatrixTSym< Element > &target, Option_t *option="S") const
Get submatrix [row_lwb..row_upb][row_lwb..row_upb]; The indexing range of the returned matrix depends...
virtual const Int_t * GetColIndexArray() const =0
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0
virtual const Int_t * GetRowIndexArray() const =0