Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TMatrixTSparse.h
Go to the documentation of this file.
1// @(#)root/matrix:$Id$
2// Authors: Fons Rademakers, Eddy Offermann Feb 2004
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_TMatrixTSparse
13#define ROOT_TMatrixTSparse
14
15#include "TMatrixTBase.h"
16#include "TMatrixTUtils.h"
17
18#include <cstring>
19
20#ifdef CBLAS
21#include <vecLib/vBLAS.h>
22//#include <cblas.h>
23#endif
24
25//////////////////////////////////////////////////////////////////////////
26// //
27// TMatrixTSparse //
28// //
29// Template class of a general sparse matrix in the Harwell-Boeing //
30// format //
31// //
32//////////////////////////////////////////////////////////////////////////
33
34template<class Element> class TMatrixT;
35
36template<class Element> class TMatrixTSparse : public TMatrixTBase<Element> {
37
38protected:
39
40 Int_t *fRowIndex; //[fNrowIndex] row index
41 Int_t *fColIndex; //[fNelems] column index
42 Element *fElements; //[fNelems]
43
44 void Allocate(Int_t nrows,Int_t ncols,Int_t row_lwb = 0,Int_t col_lwb = 0,
45 Int_t init = 0,Int_t nr_nonzeros = 0);
46
47 // Elementary constructors
50 }
51 void AMultB (const TMatrixTSparse<Element> &a,const TMatrixT<Element> &b,Int_t constr=0) {
52 const TMatrixTSparse<Element> bsp = b;
54 void AMultB (const TMatrixT<Element> &a,const TMatrixTSparse<Element> &b,Int_t constr=0) {
56
58 {
61 }
62 void AMultBt(const TMatrixTSparse<Element> &a,const TMatrixT<Element> &b,Int_t constr=0);
63 void AMultBt(const TMatrixT<Element> &a,const TMatrixTSparse<Element> &b,Int_t constr=0);
64
65 void APlusB (const TMatrixTSparse<Element> &a,const TMatrixTSparse<Element> &b,Int_t constr=0);
66 void APlusB (const TMatrixTSparse<Element> &a,const TMatrixT<Element> &b,Int_t constr=0);
67 void APlusB (const TMatrixT<Element> &a,const TMatrixTSparse<Element> &b,Int_t constr=0) { APlusB(b,a,constr); }
68
70 void AMinusB(const TMatrixTSparse<Element> &a,const TMatrixT<Element> &b,Int_t constr=0);
71 void AMinusB(const TMatrixT<Element> &a,const TMatrixTSparse<Element> &b,Int_t constr=0);
72
74 Int_t constr = 0);
75
76 Int_t ReduceSparseMatrix(Int_t nr, Int_t *row, Int_t *col, Element *data);
77
78public:
79
82
83 TMatrixTSparse() { fElements = nullptr; fRowIndex = nullptr; fColIndex = nullptr; }
84 TMatrixTSparse(Int_t nrows,Int_t ncols);
85 TMatrixTSparse(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb, Int_t nr_nonzeros = 0);
86 TMatrixTSparse(Int_t row_lwb,Int_t row_upb,Int_t col_lwb,Int_t col_upb,Int_t nr_nonzeros,
87 Int_t *row, Int_t *col,Element *data);
88 TMatrixTSparse(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb, Int_t *rowptr, Int_t *col, Element *data);
90 TMatrixTSparse(const TMatrixT<Element> &another);
91
94 TMatrixTSparse(const TMatrixTSparse<Element> &a,EMatrixCreatorsOp2 op,const TMatrixT <Element> &b);
95 TMatrixTSparse(const TMatrixT <Element> &a,EMatrixCreatorsOp2 op,const TMatrixTSparse<Element> &b);
96
98
99 const Element *GetMatrixArray () const override;
100 Element *GetMatrixArray () override;
101 const Int_t *GetRowIndexArray() const override;
103 const Int_t *GetColIndexArray() const override;
105
106 TMatrixTBase<Element> &SetRowIndexArray(Int_t *data) override { memmove(fRowIndex,data,(this->fNrows+1)*sizeof(Int_t)); return *this; }
107 TMatrixTBase<Element> &SetColIndexArray(Int_t *data) override { memmove(fColIndex,data,this->fNelems*sizeof(Int_t)); return *this; }
108
112 TMatrixTSparse<Element> &SetSparseIndexAB(const TMatrixT <Element> &a,const TMatrixTSparse<Element> &b);
114 { return SetSparseIndexAB(b,a); }
115
116 void GetMatrix2Array (Element *data,Option_t * /*option*/ ="") const override;
117 TMatrixTBase<Element> &SetMatrixArray (const Element *data,Option_t * /*option*/="") override
118 { memcpy(fElements,data,this->fNelems*sizeof(Element)); return *this; }
119 virtual TMatrixTBase<Element> &SetMatrixArray (Int_t nr_nonzeros,Int_t *irow,Int_t *icol,Element *data);
120 virtual TMatrixTBase<Element> &
121 SetMatrixArray(Int_t nr_nonzeros, Int_t nrows, Int_t ncols, Int_t *irow, Int_t *icol, Element *data);
122 TMatrixTBase<Element> &InsertRow (Int_t row,Int_t col,const Element *v,Int_t n=-1) override;
123 void ExtractRow (Int_t row,Int_t col, Element *v,Int_t n=-1) const override;
124
125 TMatrixTBase<Element> &ResizeTo(Int_t nrows,Int_t ncols,Int_t nr_nonzeros=-1) override;
126 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;
127 inline TMatrixTBase<Element> &ResizeTo(const TMatrixTSparse<Element> &m) {return ResizeTo(m.GetRowLwb(),m.GetRowUpb(),m.GetColLwb(),
128 m.GetColUpb(),m.GetNoElements()); }
129
130 void Clear(Option_t * /*option*/ ="") override { if (this->fIsOwner) {
131 if (fElements) { delete [] fElements; fElements = nullptr; }
132 if (fRowIndex) { delete [] fRowIndex; fRowIndex = nullptr; }
133 if (fColIndex) { delete [] fColIndex; fColIndex = nullptr; }
134 }
135 this->fNelems = 0;
136 this->fNrowIndex = 0;
137 }
138
139 TMatrixTSparse<Element> &Use (Int_t row_lwb,Int_t row_upb,Int_t col_lwb,Int_t col_upb,Int_t nr_nonzeros,
140 Int_t *pRowIndex,Int_t *pColIndex,Element *pData);
141 const TMatrixTSparse<Element> &Use (Int_t row_lwb,Int_t row_upb,Int_t col_lwb,Int_t col_upb,Int_t nr_nonzeros,
142 const Int_t *pRowIndex,const Int_t *pColIndex,const Element *pData) const
143 { return (const TMatrixTSparse<Element>&)
144 ((const_cast<TMatrixTSparse<Element> *>(this))->Use(row_lwb,row_upb,col_lwb,col_upb,nr_nonzeros,
145 const_cast<Int_t *>(pRowIndex),
146 const_cast<Int_t *>(pColIndex),
147 const_cast<Element *>(pData))); }
148 TMatrixTSparse<Element> &Use (Int_t nrows,Int_t ncols,Int_t nr_nonzeros,
149 Int_t *pRowIndex,Int_t *pColIndex,Element *pData);
150 const TMatrixTSparse<Element> &Use (Int_t nrows,Int_t ncols,Int_t nr_nonzeros,
151 const Int_t *pRowIndex,const Int_t *pColIndex,const Element *pData) const;
154
155 TMatrixTBase<Element> &GetSub(Int_t row_lwb,Int_t row_upb,Int_t col_lwb,Int_t col_upb,
156 TMatrixTBase<Element> &target,Option_t *option="S") const override;
157 TMatrixTSparse<Element> GetSub(Int_t row_lwb,Int_t row_upb,Int_t col_lwb,Int_t col_upb,Option_t *option="S") const;
158 TMatrixTBase<Element> &SetSub(Int_t row_lwb,Int_t col_lwb,const TMatrixTBase<Element> &source) override;
159
160 Bool_t IsSymmetric() const override { return (*this == TMatrixTSparse<Element>(kTransposed,*this)); }
162 inline TMatrixTSparse<Element> &T () { return this->Transpose(*this); }
163
164 inline void Mult(const TMatrixTSparse<Element> &a,const TMatrixTSparse<Element> &b) { AMultB(a,b,0); }
165
166 TMatrixTBase<Element> &Zero () override;
168
169 Element RowNorm () const override;
170 Element ColNorm () const override;
171 Int_t NonZeros() const override { return this->fNelems; }
172
173 TMatrixTBase<Element> &NormByDiag(const TVectorT<Element> &/*v*/,Option_t * /*option*/) override
174 { MayNotUse("NormByDiag"); return *this; }
175
176 // Either access a_ij as a(i,j)
177 Element operator()(Int_t rown,Int_t coln) const override;
178 Element &operator()(Int_t rown,Int_t coln) override;
179
180 // or as a[i][j]
182 inline TMatrixTSparseRow <Element> operator[](Int_t rown) { return TMatrixTSparseRow <Element>(*this,rown); }
183
186
191
193 if (this == &source) APlusB (tmp,tmp,1);
194 else APlusB (tmp,source,1);
195 return *this; }
197 APlusB(tmp,source,1); return *this; }
199 if (this == &source) AMinusB (tmp,tmp,1);
200 else AMinusB(tmp,source,1);
201 return *this; }
203 AMinusB(tmp,source,1); return *this; }
205 if (this == &source) AMultB (tmp,tmp,1);
206 else AMultB (tmp,source,1);
207 return *this; }
209 AMultB(tmp,source,1);
210 return *this; }
211
212 TMatrixTBase <Element> &Randomize (Element alpha,Element beta,Double_t &seed) override;
213 virtual TMatrixTSparse<Element> &RandomizePD(Element alpha,Element beta,Double_t &seed);
214
215 ClassDefOverride(TMatrixTSparse,3) // Template of Sparse Matrix class
216};
217
218// When building with -fmodules, it instantiates all pending instantiations,
219// instead of delaying them until the end of the translation unit.
220// We 'got away with' probably because the use and the definition of the
221// explicit specialization do not occur in the same TU.
222//
223// In case we are building with -fmodules, we need to forward declare the
224// specialization in order to compile the dictionary G__Matrix.cxx.
226
227template <class Element> inline const Element *TMatrixTSparse<Element>::GetMatrixArray () const { return fElements; }
228template <class Element> inline Element *TMatrixTSparse<Element>::GetMatrixArray () { return fElements; }
229template <class Element> inline const Int_t *TMatrixTSparse<Element>::GetRowIndexArray() const { return fRowIndex; }
230template <class Element> inline Int_t *TMatrixTSparse<Element>::GetRowIndexArray() { return fRowIndex; }
231template <class Element> inline const Int_t *TMatrixTSparse<Element>::GetColIndexArray() const { return fColIndex; }
232template <class Element> inline Int_t *TMatrixTSparse<Element>::GetColIndexArray() { return fColIndex; }
233
234template <class Element>
236 Int_t *pRowIndex,Int_t *pColIndex,Element *pData)
237 { return Use(0,nrows-1,0,ncols-1,nr_nonzeros,pRowIndex,pColIndex,pData); }
238template <class Element>
240 const Int_t *pRowIndex,const Int_t *pColIndex,const Element *pData) const
241 { return Use(0,nrows-1,0,ncols-1,nr_nonzeros,pRowIndex,pColIndex,pData); }
242template <class Element>
244 { R__ASSERT(a.IsValid());
245 return Use(a.GetRowLwb(),a.GetRowUpb(),a.GetColLwb(),a.GetColUpb(),
246 a.GetNoElements(),a.GetRowIndexArray(),
247 a.GetColIndexArray(),a.GetMatrixArray()); }
248template <class Element>
250 { R__ASSERT(a.IsValid());
251 return Use(a.GetRowLwb(),a.GetRowUpb(),a.GetColLwb(),a.GetColUpb(),
252 a.GetNoElements(),a.GetRowIndexArray(),
253 a.GetColIndexArray(),a.GetMatrixArray()); }
254
255template <class Element>
257 Option_t *option) const
258 {
260 this->GetSub(row_lwb,row_upb,col_lwb,col_upb,tmp,option);
261 return tmp;
262 }
263
264template <class Element> TMatrixTSparse<Element> operator+ (const TMatrixTSparse<Element> &source1,const TMatrixTSparse<Element> &source2);
265template <class Element> TMatrixTSparse<Element> operator+ (const TMatrixTSparse<Element> &source1,const TMatrixT<Element> &source2);
266template <class Element> TMatrixTSparse<Element> operator+ (const TMatrixT<Element> &source1,const TMatrixTSparse<Element> &source2);
267template <class Element> TMatrixTSparse<Element> operator+ (const TMatrixTSparse<Element> &source , Element val );
268template <class Element> TMatrixTSparse<Element> operator+ ( Element val ,const TMatrixTSparse<Element> &source );
269template <class Element> TMatrixTSparse<Element> operator- (const TMatrixTSparse<Element> &source1,const TMatrixTSparse<Element> &source2);
270template <class Element> TMatrixTSparse<Element> operator- (const TMatrixTSparse<Element> &source1,const TMatrixT<Element> &source2);
271template <class Element> TMatrixTSparse<Element> operator- (const TMatrixT<Element> &source1,const TMatrixTSparse<Element> &source2);
272template <class Element> TMatrixTSparse<Element> operator- (const TMatrixTSparse<Element> &source , Element val );
273template <class Element> TMatrixTSparse<Element> operator- ( Element val ,const TMatrixTSparse<Element> &source );
274template <class Element> TMatrixTSparse<Element> operator* (const TMatrixTSparse<Element> &source1,const TMatrixTSparse<Element> &source2);
275template <class Element> TMatrixTSparse<Element> operator* (const TMatrixTSparse<Element> &source1,const TMatrixT<Element> &source2);
276template <class Element> TMatrixTSparse<Element> operator* (const TMatrixT<Element> &source1,const TMatrixTSparse<Element> &source2);
277template <class Element> TMatrixTSparse<Element> operator* ( Element val ,const TMatrixTSparse<Element> &source );
278template <class Element> TMatrixTSparse<Element> operator* (const TMatrixTSparse<Element> &source, Element val );
279
280template <class Element> TMatrixTSparse<Element> &Add (TMatrixTSparse<Element> &target, Element scalar,
281 const TMatrixTSparse<Element> &source);
284
285template <class Element> Bool_t AreCompatible(const TMatrixTSparse<Element> &m1,const TMatrixTSparse<Element> &m2,Int_t verbose=0);
286
287#endif
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
bool Bool_t
Definition RtypesCore.h:63
int Int_t
Definition RtypesCore.h:45
double Double_t
Definition RtypesCore.h:59
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
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
TMatrixTSparse< Element > & ElementDiv(TMatrixTSparse< Element > &target, const TMatrixTSparse< Element > &source)
Divide target by the source, element-by-element.
TMatrixTSparse< Element > operator*(const TMatrixTSparse< Element > &source1, const TMatrixTSparse< Element > &source2)
TMatrixTSparse< Element > operator-(const TMatrixTSparse< Element > &source1, const TMatrixTSparse< Element > &source2)
Bool_t AreCompatible(const TMatrixTSparse< Element > &m1, const TMatrixTSparse< Element > &m2, Int_t verbose=0)
TMatrixTSparse< Element > & Add(TMatrixTSparse< Element > &target, Element scalar, const TMatrixTSparse< Element > &source)
Modify addition: target += scalar * source.
TMatrixTSparse< Element > operator+(const TMatrixTSparse< Element > &source1, const TMatrixTSparse< Element > &source2)
TMatrixTSparse< Element > & ElementMult(TMatrixTSparse< Element > &target, const TMatrixTSparse< Element > &source)
Multiply target by the source, element-by-element.
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:81
TMatrixTBase.
TMatrixTSparse.
TMatrixTBase< Element > & Randomize(Element alpha, Element beta, Double_t &seed) override
randomize matrix element values
TMatrixTBase< Element > & SetMatrixArray(const Element *data, Option_t *="") override
Copy array data to matrix .
Int_t * GetRowIndexArray() override
Int_t ReduceSparseMatrix(Int_t nr, Int_t *row, Int_t *col, Element *data)
Sum matrix entries corresponding to the same matrix element (i,j).
Int_t * GetColIndexArray() override
TMatrixTBase< Element > & GetSub(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb, TMatrixTBase< Element > &target, Option_t *option="S") const override
Get submatrix [row_lwb..row_upb][col_lwb..col_upb]; The indexing range of the returned matrix depends...
TMatrixTSparse< Element > & operator+=(Element val)
Add val to every element of the matrix.
TMatrixTBase< Element > & ResizeTo(const TMatrixTSparse< Element > &m)
static TClass * Class()
TMatrixTSparse< Element > & SetSparseIndex(Int_t nelem_new)
Increase/decrease the number of non-zero elements to nelems_new.
~TMatrixTSparse() override
TMatrixTSparse< Element > & Use(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb, Int_t nr_nonzeros, Int_t *pRowIndex, Int_t *pColIndex, Element *pData)
TMatrixTBase< Element > & InsertRow(Int_t row, Int_t col, const Element *v, Int_t n=-1) override
Insert in row rown, n elements of array v at column coln.
TMatrixTSparse< Element > & operator+=(const TMatrixTSparse< Element > &source)
Element ColNorm() const override
Column matrix norm, MAX{ SUM{ |M(i,j)|, over i}, over j}.
const Int_t * GetRowIndexArray() const override
TMatrixTSparse< Element > & operator-=(const TMatrixTSparse< Element > &source)
void AMinusB(const TMatrixTSparse< Element > &a, const TMatrixTSparse< Element > &b, Int_t constr=0)
General matrix subtraction.
TMatrixTSparse< Element > & operator*=(Element val)
Multiply every element of the matrix with val.
void AMultB(const TMatrixTSparse< Element > &a, const TMatrixTSparse< Element > &b, Int_t constr=0)
void AMultBt(const TMatrixTSparse< Element > &a, const TMatrixTSparse< Element > &b, Int_t constr=0)
void GetMatrix2Array(Element *data, Option_t *="") const override
Copy matrix data to array . It is assumed that array is of size >= fNelems.
TMatrixTSparse< Element > & Transpose(const TMatrixTSparse< Element > &source)
Transpose a matrix. Set the matrix to ncols x nrows if nrows != ncols.
TMatrixTSparse< Element > & operator*=(const TMatrixTSparse< Element > &source)
TMatrixTSparse< Element > & operator-=(Element val)
Subtract val from every element of the matrix.
TMatrixTSparse< Element > & operator+=(const TMatrixT< Element > &source)
TMatrixTBase< Element > & SetColIndexArray(Int_t *data) override
void APlusB(const TMatrixT< Element > &a, const TMatrixTSparse< Element > &b, Int_t constr=0)
Int_t NonZeros() const override
Compute the number of elements != 0.0.
TMatrixTSparseRow< Element > operator[](Int_t rown)
TMatrixTSparse< Element > & operator*=(const TMatrixT< Element > &source)
const TMatrixTSparse< Element > & Use(const TMatrixTSparse< Element > &a) const
Element * fElements
TMatrixTBase< Element > & SetSub(Int_t row_lwb, Int_t col_lwb, const TMatrixTBase< Element > &source) override
Insert matrix source starting at [row_lwb][col_lwb], thereby overwriting the part [row_lwb....
void Allocate(Int_t nrows, Int_t ncols, Int_t row_lwb=0, Int_t col_lwb=0, Int_t init=0, Int_t nr_nonzeros=0)
Allocate new matrix.
void Clear(Option_t *="") override
void AMultB(const TMatrixT< Element > &a, const TMatrixTSparse< Element > &b, Int_t constr=0)
const TMatrixTSparseRow_const< Element > operator[](Int_t rown) const
TMatrixTBase< Element > & SetRowIndexArray(Int_t *data) override
TMatrixTBase< Element > & ResizeTo(Int_t nrows, Int_t ncols, Int_t nr_nonzeros=-1) override
Set size of the matrix to nrows x ncols with nr_nonzeros non-zero entries if nr_nonzeros > 0 .
void APlusB(const TMatrixTSparse< Element > &a, const TMatrixTSparse< Element > &b, Int_t constr=0)
General matrix addition.
const TMatrixTSparse< Element > & Use(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb, Int_t nr_nonzeros, const Int_t *pRowIndex, const Int_t *pColIndex, const Element *pData) const
virtual TMatrixTSparse< Element > & RandomizePD(Element alpha, Element beta, Double_t &seed)
randomize matrix element values but keep matrix symmetric positive definite
TMatrixTBase< Element > & Zero() override
Set matrix elements to zero.
TMatrixTSparse< Element > & Use(TMatrixTSparse< Element > &a)
TMatrixTSparse< Element > & operator=(const TMatrixT< Element > &source)
Notice that the sparsity of the matrix is NOT changed : its fRowIndex/fColIndex are used !
TMatrixTBase< Element > & NormByDiag(const TVectorT< Element > &, Option_t *) override
option:
TMatrixTSparse< Element > & T()
const Int_t * GetColIndexArray() const override
void AMultB(const TMatrixTSparse< Element > &a, const TMatrixT< Element > &b, Int_t constr=0)
Element operator()(Int_t rown, Int_t coln) const override
const TMatrixTSparse< Element > & Use(Int_t nrows, Int_t ncols, Int_t nr_nonzeros, const Int_t *pRowIndex, const Int_t *pColIndex, const Element *pData) const
TMatrixTSparse< Element > & SetSparseIndexAB(const TMatrixTSparse< Element > &a, const TMatrixT< Element > &b)
TMatrixTSparse< Element > GetSub(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb, Option_t *option="S") const
void ExtractRow(Int_t row, Int_t col, Element *v, Int_t n=-1) const override
Store in array v, n matrix elements of row rown starting at column coln.
TMatrixTSparse< Element > & operator-=(const TMatrixT< Element > &source)
Bool_t IsSymmetric() const override
Check whether matrix is symmetric.
TMatrixTBase< Element > & UnitMatrix() override
Make a unit matrix (matrix need not be a square one).
void conservative_sparse_sparse_product_impl(const TMatrixTSparse< Element > &lhs, const TMatrixTSparse< Element > &rhs, Int_t constr=0)
General Sparse Matrix Multiplication (SpMM).
Element RowNorm() const override
Row matrix norm, MAX{ SUM{ |M(i,j)|, over j}, over i}.
TMatrixTSparse< Element > & SetSparseIndexAB(const TMatrixTSparse< Element > &a, const TMatrixTSparse< Element > &b)
Set the row/column indices to the "sum" of matrices a and b It is checked that enough space has been ...
void Mult(const TMatrixTSparse< Element > &a, const TMatrixTSparse< Element > &b)
const Element * GetMatrixArray() const override
Element * GetMatrixArray() override
TMatrixTSparse< Element > & Use(Int_t nrows, Int_t ncols, Int_t nr_nonzeros, Int_t *pRowIndex, Int_t *pColIndex, Element *pData)
TMatrixT.
Definition TMatrixT.h:40
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:1046
TVectorT.
Definition TVectorT.h:27
const Int_t n
Definition legend1.C:16
TMarker m
Definition textangle.C:8