Project
Loading...
Searching...
No Matches
SymMatrix.h
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
15
16#ifndef ALICEO2_FWDALIGN_SYMMATRIX_H
17#define ALICEO2_FWDALIGN_SYMMATRIX_H
18
19#include <TVectorD.h>
20#include <TString.h>
21
23
24namespace o2
25{
26namespace fwdalign
27{
28
37class SymMatrix : public MatrixSq
38{
39
40 public:
42 SymMatrix();
43
45 SymMatrix(Int_t size);
46
48 SymMatrix(const SymMatrix& mat);
49
51 ~SymMatrix() override;
52
54 void Clear(Option_t* option = "") override;
55 void Reset() override;
56
57 Int_t GetSize() const override { return fNrowIndex; }
58 Int_t GetSizeUsed() const { return fRowLwb; }
59 Int_t GetSizeBooked() const { return fNcols; }
60 Int_t GetSizeAdded() const { return fNrows; }
61
63 Float_t GetDensity() const override;
64
67
70
73
74 Double_t operator()(Int_t rown, Int_t coln) const override;
75 Double_t& operator()(Int_t rown, Int_t coln) override;
76
77 Double_t DiagElem(Int_t r) const override { return (*(const SymMatrix*)this)(r, r); }
78 Double_t& DiagElem(Int_t r) override { return (*this)(r, r); }
79
81 Double_t* GetRow(Int_t r);
82
84 void Print(const Option_t* option = "") const override;
85
87 void AddRows(int nrows = 1);
88
89 void SetSizeUsed(Int_t sz) { fRowLwb = sz; }
90
91 void Scale(Double_t coeff);
92
94 Bool_t Multiply(const SymMatrix& right);
95
99 void MultiplyByVec(const Double_t* vecIn, Double_t* vecOut) const override;
100
101 void MultiplyByVec(const TVectorD& vecIn, TVectorD& vecOut) const override;
102 void AddToRow(Int_t r, Double_t* valc, Int_t* indc, Int_t n) override;
103
104 // ---------------------------------- Dummy methods of MatrixBase
105 const Double_t* GetMatrixArray() const override { return fElems; };
106 Double_t* GetMatrixArray() override { return (Double_t*)fElems; }
107 const Int_t* GetRowIndexArray() const override
108 {
109 Error("GetRowIndexArray", "Dummy");
110 return nullptr;
111 };
112 Int_t* GetRowIndexArray() override
113 {
114 Error("GetRowIndexArray", "Dummy");
115 return nullptr;
116 };
117 const Int_t* GetColIndexArray() const override
118 {
119 Error("GetColIndexArray", "Dummy");
120 return nullptr;
121 };
122 Int_t* GetColIndexArray() override
123 {
124 Error("GetColIndexArray", "Dummy");
125 return nullptr;
126 };
128 {
129 Error("SetRowIndexArray", "Dummy");
130 return *this;
131 }
133 {
134 Error("SetColIndexArray", "Dummy");
135 return *this;
136 }
137 TMatrixDBase& GetSub(Int_t, Int_t, Int_t, Int_t, TMatrixDBase&, Option_t*) const override
138 {
139 Error("GetSub", "Dummy");
140 return *((TMatrixDBase*)this);
141 }
142 TMatrixDBase& SetSub(Int_t, Int_t, const TMatrixDBase&) override
143 {
144 Error("GetSub", "Dummy");
145 return *this;
146 }
147 TMatrixDBase& ResizeTo(Int_t, Int_t, Int_t) override
148 {
149 Error("ResizeTo", "Dummy");
150 return *this;
151 }
152 TMatrixDBase& ResizeTo(Int_t, Int_t, Int_t, Int_t, Int_t) override
153 {
154 Error("ResizeTo", "Dummy");
155 return *this;
156 }
157
158 // ----------------------------- Choleski methods ----------------------------------------
168
170 void InvertChol(SymMatrix* mchol);
171
173 Bool_t InvertChol();
174
183 Bool_t SolveChol(Double_t* brhs, Bool_t invert = kFALSE);
184
185 Bool_t SolveChol(Double_t* brhs, Double_t* bsol, Bool_t invert = kFALSE);
186 Bool_t SolveChol(TVectorD& brhs, Bool_t invert = kFALSE);
187 Bool_t SolveChol(const TVectorD& brhs, TVectorD& bsol, Bool_t invert = kFALSE);
188
198 Bool_t SolveCholN(Double_t* bn, int nRHS, Bool_t invert = kFALSE);
199
204 int SolveSpmInv(double* vecB, Bool_t stabilize = kTRUE);
205
206 protected:
207 virtual Int_t GetIndex(Int_t row, Int_t col) const;
208 Double_t GetEl(Int_t row, Int_t col) const { return operator()(row, col); }
209 void SetEl(Int_t row, Int_t col, Double_t val) { operator()(row, col) = val; }
210
211 protected:
212 Double_t* fElems;
213 Double_t** fElemsAdd;
214
216 static Int_t fgCopyCnt;
217
219};
220
221//___________________________________________________________
222inline Int_t SymMatrix::GetIndex(Int_t row, Int_t col) const
223{
224 // lower triangle is actually filled
225 return ((row * (row + 1)) >> 1) + col;
226}
227
228//___________________________________________________________
229inline Double_t SymMatrix::operator()(Int_t row, Int_t col) const
230{
231 //
232 if (row < col) {
233 Swap(row, col);
234 }
235 if (row >= fNrowIndex) {
236 return 0;
237 }
238 return (const Double_t&)(row < fNcols ? fElems[GetIndex(row, col)] : (fElemsAdd[row - fNcols])[col]);
239}
240
241//___________________________________________________________
242inline Double_t& SymMatrix::operator()(Int_t row, Int_t col)
243{
244 if (row < col) {
245 Swap(row, col);
246 }
247 if (row >= fNrowIndex) {
248 AddRows(row - fNrowIndex + 1);
249 }
250 return (row < fNcols ? fElems[GetIndex(row, col)] : (fElemsAdd[row - fNcols])[col]);
251}
252
253//___________________________________________________________
254inline void SymMatrix::MultiplyByVec(const TVectorD& vecIn, TVectorD& vecOut) const
255{
256 MultiplyByVec(vecIn.GetMatrixArray(), vecOut.GetMatrixArray());
257}
258
259//___________________________________________________________
260inline void SymMatrix::Scale(Double_t coeff)
261{
262 for (int i = fNrowIndex; i--;) {
263 for (int j = i; j--;) {
264 double& el = operator()(i, j);
265 if (el) {
266 el *= coeff;
267 }
268 }
269 }
270}
271
272//___________________________________________________________
273inline void SymMatrix::AddToRow(Int_t r, Double_t* valc, Int_t* indc, Int_t n)
274{
275 for (int i = n; i--;) {
276 (*this)(indc[i], r) += valc[i];
277 }
278}
279
280} // namespace fwdalign
281} // namespace o2
282
283#endif
int32_t i
Abstract class (from AliROOT) for square matrix used for millepede2 operation.
uint32_t j
Definition RawData.h:0
uint32_t col
Definition RawData.h:4
void Swap(int &r, int &c) const
Definition MatrixSq.h:137
Fast symmetric matrix with dynamically expandable size.
Definition SymMatrix.h:38
SymMatrix()
default constructor
Definition SymMatrix.cxx:29
Double_t operator()(Int_t rown, Int_t coln) const override
Definition SymMatrix.h:229
void Scale(Double_t coeff)
Definition SymMatrix.h:260
SymMatrix & operator-=(const SymMatrix &src)
minus operator
Int_t GetSizeUsed() const
Definition SymMatrix.h:58
static Int_t fgCopyCnt
matrix copy counter
Definition SymMatrix.h:216
void Print(const Option_t *option="") const override
print itself
TMatrixDBase & GetSub(Int_t, Int_t, Int_t, Int_t, TMatrixDBase &, Option_t *) const override
Definition SymMatrix.h:137
int SolveSpmInv(double *vecB, Bool_t stabilize=kTRUE)
Obtain solution of a system of linear equations with symmetric matrix and the inverse (using 'singula...
~SymMatrix() override
destructor
Definition SymMatrix.cxx:82
const Double_t * GetMatrixArray() const override
Definition SymMatrix.h:105
Double_t ** fElemsAdd
Elements (rows) added dynamicaly.
Definition SymMatrix.h:213
Double_t & DiagElem(Int_t r) override
Definition SymMatrix.h:78
ClassDefOverride(SymMatrix, 0)
void SetEl(Int_t row, Int_t col, Double_t val)
Definition SymMatrix.h:209
Bool_t SolveChol(Double_t *brhs, Bool_t invert=kFALSE)
Solves the set of n linear equations A x = b.
void AddToRow(Int_t r, Double_t *valc, Int_t *indc, Int_t n) override
Definition SymMatrix.h:273
void MultiplyByVec(const Double_t *vecIn, Double_t *vecOut) const override
fill vecOut by matrix*vecIn
Int_t GetSize() const override
Definition SymMatrix.h:57
Double_t * GetMatrixArray() override
Definition SymMatrix.h:106
Float_t GetDensity() const override
get fraction of non-zero elements
Double_t DiagElem(Int_t r) const override
Definition SymMatrix.h:77
TMatrixDBase & SetColIndexArray(Int_t *) override
Definition SymMatrix.h:132
static SymMatrix * fgBuffer
buffer for fast solution
Definition SymMatrix.h:215
Double_t * fElems
Elements booked by constructor.
Definition SymMatrix.h:212
TMatrixDBase & ResizeTo(Int_t, Int_t, Int_t, Int_t, Int_t) override
Definition SymMatrix.h:152
Int_t GetSizeAdded() const
Definition SymMatrix.h:60
SymMatrix & operator=(const SymMatrix &src)
assignment operator
Definition SymMatrix.cxx:92
const Int_t * GetRowIndexArray() const override
Definition SymMatrix.h:107
Bool_t InvertChol()
Invert matrix using Choleski decomposition.
SymMatrix * DecomposeChol()
Return a matrix with Choleski decomposition.
Double_t GetEl(Int_t row, Int_t col) const
Definition SymMatrix.h:208
Bool_t SolveCholN(Double_t *bn, int nRHS, Bool_t invert=kFALSE)
Solves the set of n linear equations A x = b; this version solve multiple RHSs at once.
virtual Int_t GetIndex(Int_t row, Int_t col) const
Definition SymMatrix.h:222
void Clear(Option_t *option="") override
clear dynamic part
TMatrixDBase & SetRowIndexArray(Int_t *) override
Definition SymMatrix.h:127
void Reset() override
Double_t * GetRow(Int_t r)
get pointer on the row
Int_t * GetColIndexArray() override
Definition SymMatrix.h:122
Int_t GetSizeBooked() const
Definition SymMatrix.h:59
void AddRows(int nrows=1)
add empty rows
const Int_t * GetColIndexArray() const override
Definition SymMatrix.h:117
Int_t * GetRowIndexArray() override
Definition SymMatrix.h:112
SymMatrix & operator+=(const SymMatrix &src)
add operator
TMatrixDBase & ResizeTo(Int_t, Int_t, Int_t) override
Definition SymMatrix.h:147
void SetSizeUsed(Int_t sz)
Definition SymMatrix.h:89
Bool_t Multiply(const SymMatrix &right)
multiply from the right
TMatrixDBase & SetSub(Int_t, Int_t, const TMatrixDBase &) override
Definition SymMatrix.h:142
GLdouble n
Definition glcorearb.h:1982
GLenum src
Definition glcorearb.h:1767
GLsizeiptr size
Definition glcorearb.h:659
GLdouble GLdouble right
Definition glcorearb.h:4077
GLuint GLfloat * val
Definition glcorearb.h:1582
GLboolean r
Definition glcorearb.h:1233
GLboolean invert
Definition glcorearb.h:543
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
std::vector< int > row