Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TMatrixTSym.h
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#ifndef ROOT_TMatrixTSym
13#define ROOT_TMatrixTSym
14
15//////////////////////////////////////////////////////////////////////////
16// //
17// TMatrixTSym //
18// //
19// Implementation of a symmetric matrix in the linear algebra package //
20// //
21// Note that in this implementation both matrix element m[i][j] and //
22// m[j][i] are updated and stored in memory. However, when making the //
23// object persistent only the upper right triangle is stored. //
24// //
25//////////////////////////////////////////////////////////////////////////
26
27#include "TMatrixTBase.h"
28#include "TMatrixTUtils.h"
29
30#include <cassert>
31
32template<class Element>class TMatrixT;
33template<class Element>class TMatrixTSymLazy;
34template<class Element>class TVectorT;
35
36template<class Element> class TMatrixTSym : public TMatrixTBase<Element> {
37
38protected:
39
40 Element fDataStack[TMatrixTBase<Element>::kSizeMax]; //! data container
41 Element *fElements; //[fNelems] elements themselves
42
43 Element *New_m (Int_t size);
44 void Delete_m(Int_t size,Element*&);
45 Int_t Memcpy_m(Element *newp,const Element *oldp,Int_t copySize,
46 Int_t newSize,Int_t oldSize);
47 void Allocate(Int_t nrows,Int_t ncols,Int_t row_lwb = 0,Int_t col_lwb = 0,Int_t init = 0,
48 Int_t /*nr_nonzeros*/ = -1);
49
50public:
51
52 enum {kWorkMax = 100}; // size of work array
55
56 TMatrixTSym() { fElements = nullptr; }
57 explicit TMatrixTSym(Int_t nrows);
58 TMatrixTSym(Int_t row_lwb,Int_t row_upb);
59 TMatrixTSym(Int_t nrows,const Element *data,Option_t *option="");
60 TMatrixTSym(Int_t row_lwb,Int_t row_upb,const Element *data,Option_t *option="");
61 TMatrixTSym(const TMatrixTSym<Element> &another);
62 template <class Element2> TMatrixTSym(const TMatrixTSym<Element2> &another)
63 {
64 R__ASSERT(another.IsValid());
65 Allocate(another.GetNrows(),another.GetNcols(),another.GetRowLwb(),another.GetColLwb());
66 *this = another;
67 }
68
70 TMatrixTSym(EMatrixCreatorsOp1 op,const TMatrixT <Element> &prototype);
72 TMatrixTSym(const TMatrixTSymLazy<Element> &lazy_constructor);
73
75
76 // Elementary constructors
77 void TMult(const TMatrixT <Element> &a);
78 void TMult(const TMatrixTSym<Element> &a);
79 void Mult (const TMatrixTSym<Element> &a) { TMult(a); }
80
83
84 inline void SetElement(Int_t rown, Int_t coln, Element val);
85
86 const Element *GetMatrixArray () const override;
87 Element *GetMatrixArray () override;
88 const Int_t *GetRowIndexArray() const override { return nullptr; }
89 Int_t *GetRowIndexArray() override { return nullptr; }
90 const Int_t *GetColIndexArray() const override { return nullptr; }
91 Int_t *GetColIndexArray() override { return nullptr; }
92
93 TMatrixTBase<Element> &SetRowIndexArray(Int_t * /*data*/) override { MayNotUse("SetRowIndexArray(Int_t *)"); return *this; }
94 TMatrixTBase<Element> &SetColIndexArray(Int_t * /*data*/) override { MayNotUse("SetColIndexArray(Int_t *)"); return *this; }
95
96 void Clear (Option_t * /*option*/ ="") override { if (this->fIsOwner) Delete_m(this->fNelems,fElements);
97 else fElements = nullptr;
98 this->fNelems = 0; }
99 Bool_t IsSymmetric() const override { return kTRUE; }
100
101 TMatrixTSym <Element> &Use (Int_t row_lwb,Int_t row_upb,Element *data);
102 const TMatrixTSym <Element> &Use (Int_t row_lwb,Int_t row_upb,const Element *data) const
103 { return (const TMatrixTSym<Element>&)
104 ((const_cast<TMatrixTSym<Element> *>(this))->Use(row_lwb,row_upb,const_cast<Element *>(data))); }
105 TMatrixTSym <Element> &Use (Int_t nrows,Element *data);
106 const TMatrixTSym <Element> &Use (Int_t nrows,const Element *data) const;
107 TMatrixTSym <Element> &Use (TMatrixTSym<Element> &a);
108 const TMatrixTSym <Element> &Use (const TMatrixTSym<Element> &a) const;
109
110 TMatrixTSym <Element> &GetSub (Int_t row_lwb,Int_t row_upb,TMatrixTSym<Element> &target,Option_t *option="S") const;
111 TMatrixTBase<Element> &GetSub (Int_t row_lwb,Int_t row_upb,Int_t col_lwb,Int_t col_upb,
112 TMatrixTBase<Element> &target,Option_t *option="S") const override;
113 TMatrixTSym <Element> GetSub (Int_t row_lwb,Int_t row_upb,Int_t col_lwb,Int_t col_upb,Option_t *option="S") const;
114 TMatrixTSym <Element> &SetSub (Int_t row_lwb,const TMatrixTBase<Element> &source);
115 TMatrixTBase<Element> &SetSub (Int_t row_lwb,Int_t col_lwb,const TMatrixTBase<Element> &source) override;
116
117 TMatrixTBase<Element> &SetMatrixArray(const Element *data, Option_t *option="") override;
118
119 TMatrixTBase<Element> &Shift (Int_t row_shift,Int_t col_shift) override;
120 TMatrixTBase<Element> &ResizeTo (Int_t nrows,Int_t ncols,Int_t /*nr_nonzeros*/ =-1) override;
121 TMatrixTBase<Element> &ResizeTo (Int_t row_lwb,Int_t row_upb,Int_t col_lwb,Int_t col_upb,Int_t /*nr_nonzeros*/ =-1) override;
123 return ResizeTo(m.GetRowLwb(),m.GetRowUpb(),m.GetColLwb(),m.GetColUpb()); }
124
125 Double_t Determinant () const override;
126 void Determinant (Double_t &d1,Double_t &d2) const override;
127
128 TMatrixTSym<Element> &Invert (Double_t *det=nullptr);
131 inline TMatrixTSym<Element> &T () { return this->Transpose(*this); }
132 TMatrixTSym<Element> &Rank1Update (const TVectorT <Element> &v,Element alpha=1.0);
133 TMatrixTSym<Element> &Similarity (const TMatrixT <Element> &n);
135 Element Similarity (const TVectorT <Element> &v) const;
136 TMatrixTSym<Element> &SimilarityT (const TMatrixT <Element> &n);
137
138 // Either access a_ij as a(i,j)
139 inline Element operator()(Int_t rown,Int_t coln) const override;
140 inline Element &operator()(Int_t rown,Int_t coln) override; ///< Access element a_ij where i=rown and j=coln. \warning Modifying this element by the caller breaks the symmetry of the matrix if a_ji is not modified accordingly.
141
142 // or as a[i][j]
143 inline const TMatrixTRow_const<Element> operator[](Int_t rown) const { return TMatrixTRow_const<Element>(*this,rown); }
144 inline TMatrixTRow <Element> operator[](Int_t rown) { return TMatrixTRow <Element>(*this,rown); } ///< Access row a_i where i=rown. \note A concatenated call to [coln] allows to access element a_ij where i=rown and j=coln \warning Modifying this row by the caller breaks the symmetry of the matrix if a_j is not modified accordingly.
145
146 TMatrixTSym<Element> &operator= (const TMatrixTSym <Element> &source);
148 template <class Element2> TMatrixTSym<Element> &operator= (const TMatrixTSym<Element2> &source)
149 {
150 if (!AreCompatible(*this,source)) {
151 Error("operator=(const TMatrixTSym2 &)","matrices not compatible");
152 return *this;
153 }
154
155 TObject::operator=(source);
156 const Element2 * const ps = source.GetMatrixArray();
157 Element * const pt = this->GetMatrixArray();
158 for (Int_t i = 0; i < this->fNelems; i++)
159 pt[i] = ps[i];
160 this->fTol = source.GetTol();
161 return *this;
162 }
163
164 TMatrixTSym<Element> &operator= (Element val);
165 TMatrixTSym<Element> &operator-=(Element val);
166 TMatrixTSym<Element> &operator+=(Element val);
167 TMatrixTSym<Element> &operator*=(Element val);
168
169 TMatrixTSym &operator+=(const TMatrixTSym &source);
170 TMatrixTSym &operator-=(const TMatrixTSym &source);
171
172 TMatrixTBase<Element> &Apply(const TElementActionT <Element> &action) override;
174
175 TMatrixTBase<Element> &Randomize (Element alpha,Element beta,Double_t &seed) override;
176 virtual TMatrixTSym <Element> &RandomizePD(Element alpha,Element beta,Double_t &seed);
177
178 const TMatrixT<Element> EigenVectors(TVectorT<Element> &eigenValues) const;
179
180 ClassDefOverride(TMatrixTSym,2) // Template of Symmetric Matrix class
181};
182// When building with -fmodules, it instantiates all pending instantiations,
183// instead of delaying them until the end of the translation unit.
184// We 'got away with' probably because the use and the definition of the
185// explicit specialization do not occur in the same TU.
186//
187// In case we are building with -fmodules, we need to forward declare the
188// specialization in order to compile the dictionary G__Matrix.cxx.
190
191template <class Element> inline const Element *TMatrixTSym<Element>::GetMatrixArray() const { return fElements; }
192template <class Element> inline Element *TMatrixTSym<Element>::GetMatrixArray() { return fElements; }
193
194template <class Element> inline TMatrixTSym<Element> &TMatrixTSym<Element>::Use (Int_t nrows,Element *data) { return Use(0,nrows-1,data); }
195template <class Element> inline const TMatrixTSym<Element> &TMatrixTSym<Element>::Use (Int_t nrows,const Element *data) const
196 { return Use(0,nrows-1,data); }
198 { return Use(a.GetRowLwb(),a.GetRowUpb(),a.GetMatrixArray()); }
199template <class Element> inline const TMatrixTSym<Element> &TMatrixTSym<Element>::Use (const TMatrixTSym<Element> &a) const
200 { return Use(a.GetRowLwb(),a.GetRowUpb(),a.GetMatrixArray()); }
201
202template <class Element> inline TMatrixTSym<Element> TMatrixTSym<Element>::GetSub (Int_t row_lwb,Int_t row_upb,Int_t col_lwb,Int_t col_upb,
203 Option_t *option) const
204 {
206 this->GetSub(row_lwb,row_upb,col_lwb,col_upb,tmp,option);
207 return tmp;
208 }
209
210template <class Element> inline Element TMatrixTSym<Element>::operator()(Int_t rown,Int_t coln) const
211{
212 R__ASSERT(this->IsValid());
213 const Int_t arown = rown-this->fRowLwb;
214 const Int_t acoln = coln-this->fColLwb;
215 if (arown >= this->fNrows || arown < 0) {
216 Error("operator()","Request row(%d) outside matrix range of %d - %d",rown,this->fRowLwb,this->fRowLwb+this->fNrows);
218 }
219 if (acoln >= this->fNcols || acoln < 0) {
220 Error("operator()","Request column(%d) outside matrix range of %d - %d",coln,this->fColLwb,this->fColLwb+this->fNcols);
222 }
223 return (fElements[arown*this->fNcols+acoln]);
224}
225
226template <class Element> inline Element &TMatrixTSym<Element>::operator()(Int_t rown,Int_t coln)
227{
228 R__ASSERT(this->IsValid());
229 const Int_t arown = rown-this->fRowLwb;
230 const Int_t acoln = coln-this->fColLwb;
231 if (arown >= this->fNrows || arown < 0) {
232 Error("operator()","Request row(%d) outside matrix range of %d - %d",rown,this->fRowLwb,this->fRowLwb+this->fNrows);
234 }
235 if (acoln >= this->fNcols || acoln < 0) {
236 Error("operator()","Request column(%d) outside matrix range of %d - %d",coln,this->fColLwb,this->fColLwb+this->fNcols);
238 }
239 return (fElements[arown*this->fNcols+acoln]);
240}
241
242////////////////////////////////////////////////////////////////////////////////
243/// Efficiently sets element (rown,coln) equal to val
244/// Index bound checks can be deactivated by defining NDEBUG
245
246template <class Element>
247inline void TMatrixTSym<Element>::SetElement(Int_t rown, Int_t coln, Element val)
248{
249 assert(this->IsValid());
250 rown = rown - this->fRowLwb;
251 coln = coln - this->fColLwb;
252 assert((rown < this->fNrows && rown >= 0) && "SetElement() error: row index outside matrix range");
253 assert((coln < this->fNcols && coln >= 0) && "SetElement() error: column index outside matrix range");
254 fElements[rown * this->fNcols + coln] = val;
255}
256
257template <class Element> Bool_t operator== (const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2);
258template <class Element> TMatrixTSym<Element> operator+ (const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2);
259template <class Element> TMatrixTSym<Element> operator+ (const TMatrixTSym<Element> &source1, Element val);
260template <class Element> TMatrixTSym<Element> operator+ ( Element val ,const TMatrixTSym<Element> &source2);
261template <class Element> TMatrixTSym<Element> operator- (const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2);
262template <class Element> TMatrixTSym<Element> operator- (const TMatrixTSym<Element> &source1, Element val);
263template <class Element> TMatrixTSym<Element> operator- ( Element val ,const TMatrixTSym<Element> &source2);
264template <class Element> TMatrixTSym<Element> operator* (const TMatrixTSym<Element> &source, Element val );
265template <class Element> TMatrixTSym<Element> operator* ( Element val, const TMatrixTSym<Element> &source );
266// Preventing warnings with -Weffc++ in GCC since overloading the || and && operators was a design choice.
267#if (__GNUC__ * 10000 + __GNUC_MINOR__ * 100 + __GNUC_PATCHLEVEL__) >= 40600
268#pragma GCC diagnostic push
269#pragma GCC diagnostic ignored "-Weffc++"
270#endif
271template <class Element> TMatrixTSym<Element> operator&& (const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2);
272template <class Element> TMatrixTSym<Element> operator|| (const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2);
273#if (__GNUC__ * 10000 + __GNUC_MINOR__ * 100 + __GNUC_PATCHLEVEL__) >= 40600
274#pragma GCC diagnostic pop
275#endif
276template <class Element> TMatrixTSym<Element> operator> (const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2);
277template <class Element> TMatrixTSym<Element> operator>= (const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2);
278template <class Element> TMatrixTSym<Element> operator<= (const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2);
279template <class Element> TMatrixTSym<Element> operator< (const TMatrixTSym<Element> &source1,const TMatrixTSym<Element> &source2);
280
281template <class Element> TMatrixTSym<Element> &Add (TMatrixTSym<Element> &target, Element scalar,const TMatrixTSym<Element> &source);
283template <class Element> TMatrixTSym<Element> &ElementDiv (TMatrixTSym<Element> &target,const TMatrixTSym<Element> &source);
284
285#endif
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
bool Bool_t
Definition RtypesCore.h:63
int Int_t
Definition RtypesCore.h:45
double Double_t
Definition RtypesCore.h:59
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
const char Option_t
Definition RtypesCore.h:66
#define ClassDefOverride(name, id)
Definition Rtypes.h:346
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:185
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t target
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
TMatrixTSym< Element > operator+(const TMatrixTSym< Element > &source1, const TMatrixTSym< Element > &source2)
TMatrixTSym< Element > operator-(const TMatrixTSym< Element > &source1, const TMatrixTSym< Element > &source2)
TMatrixTSym< Element > & Add(TMatrixTSym< Element > &target, Element scalar, const TMatrixTSym< Element > &source)
Modify addition: target += scalar * source.
TMatrixTSym< Element > & ElementMult(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
TMatrixTSym< Element > operator>=(const TMatrixTSym< Element > &source1, const TMatrixTSym< Element > &source2)
source1 >= source2
TMatrixTSym< Element > operator<=(const TMatrixTSym< Element > &source1, const TMatrixTSym< Element > &source2)
source1 <= source2
TMatrixTSym< Element > operator||(const TMatrixTSym< Element > &source1, const TMatrixTSym< Element > &source2)
Logical Or.
TMatrixTSym< Element > operator*(const TMatrixTSym< Element > &source, Element val)
TMatrixTSym< Element > operator&&(const TMatrixTSym< Element > &source1, const TMatrixTSym< Element > &source2)
Logical AND.
Bool_t operator==(const TMatrixTSym< Element > &source1, const TMatrixTSym< Element > &source2)
Check to see if two matrices are identical.
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:81
TMatrixTBase.
Int_t GetNrows() const
Int_t GetRowLwb() const
Int_t GetColLwb() const
static Element & NaNValue()
Bool_t IsValid() const
Int_t GetNcols() const
Element GetTol() const
TMatrixTSym.
Definition TMatrixTSym.h:36
TMatrixTSym< Element > GetSub(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb, Option_t *option="S") const
Element & operator()(Int_t rown, Int_t coln) override
Access element a_ij where i=rown and j=coln.
TMatrixTBase< Element > & ResizeTo(const TMatrixTSym< Element > &m)
TMatrixTSym< Element > & operator*=(Element val)
Multiply every element of the matrix with val.
TMatrixTBase< Element > & ResizeTo(Int_t nrows, Int_t ncols, Int_t=-1) override
Set size of the matrix to nrows x ncols New dynamic elements are created, the overlapping part of the...
TMatrixTSym(const TMatrixTSym< Element2 > &another)
Definition TMatrixTSym.h:62
TMatrixTRow< Element > operator[](Int_t rown)
Access row a_i where i=rown.
const Int_t * GetColIndexArray() const override
Definition TMatrixTSym.h:90
TMatrixTSym< Element > & operator+=(Element val)
Add val to every element of the matrix.
const TMatrixTSym< Element > & Use(Int_t nrows, const Element *data) const
const TMatrixTSym< Element > & Use(Int_t row_lwb, Int_t row_upb, const Element *data) const
virtual TMatrixTSym< Element > & RandomizePD(Element alpha, Element beta, Double_t &seed)
randomize matrix element values but keep matrix symmetric positive definite
TMatrixTBase< Element > & Randomize(Element alpha, Element beta, Double_t &seed) override
randomize matrix element values but keep matrix symmetric
Int_t * GetRowIndexArray() override
Definition TMatrixTSym.h:89
TMatrixTBase< Element > & SetColIndexArray(Int_t *) override
Definition TMatrixTSym.h:94
TMatrixTSym< Element > & Transpose(const TMatrixTSym< Element > &source)
Transpose a matrix.
Element * New_m(Int_t size)
return data pointer .
TMatrixTSym< Element > & Use(Int_t nrows, Element *data)
TMatrixTSym< Element > & Use(Int_t row_lwb, Int_t row_upb, Element *data)
const TMatrixTRow_const< Element > operator[](Int_t rown) const
TMatrixTSym< Element > & Rank1Update(const TVectorT< Element > &v, Element alpha=1.0)
Perform a rank 1 operation on the matrix: A += alpha * v * v^T.
TMatrixTSym< Element > & Use(TMatrixTSym< Element > &a)
TMatrixTSym< Element > & T()
void Minus(const TMatrixTSym< Element > &a, const TMatrixTSym< Element > &b)
Symmetric matrix subtraction. Replace this matrix with C such that C = A - B.
void TMult(const TMatrixT< Element > &a)
Replace this matrix with C such that C = A' * A.
const Int_t * GetRowIndexArray() const override
Definition TMatrixTSym.h:88
void Delete_m(Int_t size, Element *&)
delete data pointer m, if it was assigned on the heap
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.
Element * GetMatrixArray() override
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...
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 .
const TMatrixT< Element > EigenVectors(TVectorT< Element > &eigenValues) const
Return a matrix containing the eigen-vectors ordered by descending eigen-values.
TMatrixTSym< Element > & operator-=(Element val)
Subtract val from every element of the matrix.
Bool_t IsSymmetric() const override
Check whether matrix is symmetric.
Definition TMatrixTSym.h:99
void Plus(const TMatrixTSym< Element > &a, const TMatrixTSym< Element > &b)
Symmetric matrix summation. Replace this matrix with C such that C = A + B.
Element fDataStack[TMatrixTBase< Element >::kSizeMax]
Definition TMatrixTSym.h:40
TMatrixTSym< Element > & InvertFast(Double_t *det=nullptr)
Invert the matrix and calculate its determinant.
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...
void SetElement(Int_t rown, Int_t coln, Element val)
Efficiently sets element (rown,coln) equal to val Index bound checks can be deactivated by defining N...
static TClass * Class()
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...
TMatrixTSym< Element > & operator=(const TMatrixTSym< Element > &source)
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....
TMatrixTBase< Element > & SetMatrixArray(const Element *data, Option_t *option="") override
Copy array data to matrix .
void Clear(Option_t *="") override
Definition TMatrixTSym.h:96
void Mult(const TMatrixTSym< Element > &a)
Definition TMatrixTSym.h:79
Element * fElements
data container
Definition TMatrixTSym.h:41
TMatrixTBase< Element > & Shift(Int_t row_shift, Int_t col_shift) override
Shift the row index by adding row_shift and the column index by adding col_shift, respectively.
TMatrixTBase< Element > & Apply(const TElementActionT< Element > &action) override
Apply action to each matrix element.
~TMatrixTSym() override
Definition TMatrixTSym.h:74
Double_t Determinant() const override
TMatrixTSym< Element > & Invert(Double_t *det=nullptr)
Invert the matrix and calculate its determinant Notice that the LU decomposition is used instead of B...
const TMatrixTSym< Element > & Use(const TMatrixTSym< Element > &a) const
Element operator()(Int_t rown, Int_t coln) const override
const Element * GetMatrixArray() const override
TMatrixTBase< Element > & SetRowIndexArray(Int_t *) override
Definition TMatrixTSym.h:93
Int_t * GetColIndexArray() override
Definition TMatrixTSym.h:91
TMatrixT.
Definition TMatrixT.h:40
TObject & operator=(const TObject &rhs)
TObject assignment operator.
Definition TObject.h:296
void MayNotUse(const char *method) const
Use this method to signal that a method (defined in a base class) may not be called in a derived clas...
Definition TObject.cxx:1058
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1005
TVectorT.
Definition TVectorT.h:29
TPaveText * pt
const Int_t n
Definition legend1.C:16
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 .
TMarker m
Definition textangle.C:8