Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ErrorSymMatrix.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// Class Description:
28//
29// Simplified version of CLHEP HepSymMatrix class.
30
31// History:
32// - Imported from CLHEP and modified: P. Arce May 2007
33// --------------------------------------------------------------------
34
35#ifndef G4ErrorSymMatrix_hh
36#define G4ErrorSymMatrix_hh
37
38#include <vector>
39#include "globals.hh"
40
41class G4ErrorMatrix;
42
44{
45 public: //with description
46
48 // Default constructor. Gives 0x0 symmetric matrix.
49 // Another G4ErrorSymMatrix can be assigned to it.
50
51 explicit G4ErrorSymMatrix(G4int p);
53 // Constructor. Gives p x p symmetric matrix.
54 // With a second argument, the matrix is initialized. 0 means a zero
55 // matrix, 1 means the identity matrix.
56
59 // Copy and move constructor.
60
61 // Constructor from DiagMatrix
62
63 virtual ~G4ErrorSymMatrix();
64 // Destructor.
65
66 inline G4int num_row() const;
67 inline G4int num_col() const;
68 // Returns number of rows/columns.
69
70 const G4double & operator()(G4int row, G4int col) const;
72 // Read and write a G4ErrorSymMatrix element.
73 // ** Note that indexing starts from (1,1). **
74
75 const G4double & fast(G4int row, G4int col) const;
76 G4double & fast(G4int row, G4int col);
77 // fast element access.
78 // Must be row>=col;
79 // ** Note that indexing starts from (1,1). **
80
81 void assign(const G4ErrorMatrix &m2);
82 // Assigns m2 to s, assuming m2 is a symmetric matrix.
83
84 void assign(const G4ErrorSymMatrix &m2);
85 // Another form of assignment. For consistency.
86
88 // Multiply a G4ErrorSymMatrix by a floating number.
89
91 // Divide a G4ErrorSymMatrix by a floating number.
92
95 // Add or subtract a G4ErrorSymMatrix.
96
99 // Assignment and move operators.
100 // Notice that there is no G4ErrorSymMatrix = Matrix.
101
103 // unary minus, ie. flip the sign of each element.
104
106 // Returns the transpose of a G4ErrorSymMatrix (which is itself).
107
109 // Apply a function to all elements of the matrix.
110
113 // Returns m1*s*m1.T().
114
116 // temporary. test of new similarity.
117 // Returns m1.T()*s*m1.
118
119 G4ErrorSymMatrix sub(G4int min_row, G4int max_row) const;
120 // Returns a sub matrix of a G4ErrorSymMatrix.
121
122 void sub(G4int row, const G4ErrorSymMatrix &m1);
123 // Sub matrix of this G4ErrorSymMatrix is replaced with m1.
124
125 G4ErrorSymMatrix sub(G4int min_row, G4int max_row);
126 // SGI CC bug. I have to have both with/without const. I should not need
127 // one without const.
128
129 inline G4ErrorSymMatrix inverse(G4int &ifail) const;
130 // Invert a Matrix. The matrix is not changed
131 // Returns 0 when successful, otherwise non-zero.
132
133 void invert(G4int &ifail);
134 // Invert a Matrix.
135 // N.B. the contents of the matrix are replaced by the inverse.
136 // Returns ierr = 0 when successful, otherwise non-zero.
137 // This method has less overhead then inverse().
138
139 G4double determinant() const;
140 // calculate the determinant of the matrix.
141
142 G4double trace() const;
143 // calculate the trace of the matrix (sum of diagonal elements).
144
146 {
147 public:
150 private:
152 G4int _r;
153 };
155 {
156 public:
158 inline const G4double & operator[](G4int) const;
159 private:
160 const G4ErrorSymMatrix& _a;
161 G4int _r;
162 };
163 // helper class to implement m[i][j]
164
167 // Read or write a matrix element.
168 // While it may not look like it, you simply do m[i][j] to get an
169 // element.
170 // ** Note that the indexing starts from [0][0]. **
171
172 // Special-case inversions for 5x5 and 6x6 symmetric positive definite:
173 // These set ifail=0 and invert if the matrix was positive definite;
174 // otherwise ifail=1 and the matrix is left unaltered.
175
176 void invertCholesky5 (G4int &ifail);
177 void invertCholesky6 (G4int &ifail);
178
179 // Inversions for 5x5 and 6x6 forcing use of specific methods: The
180 // behavior (though not the speed) will be identical to invert(ifail).
181
182 void invertHaywood4 (G4int & ifail);
183 void invertHaywood5 (G4int &ifail);
184 void invertHaywood6 (G4int &ifail);
185 void invertBunchKaufman (G4int &ifail);
186
187 protected:
188
189 inline G4int num_size() const;
190
191 private:
192
195 friend class G4ErrorMatrix;
196
199 friend void diag_step(G4ErrorSymMatrix *t, G4int begin, G4int end);
201 G4int begin, G4int end);
204 G4int row, G4int col);
205
207 const G4ErrorSymMatrix &m2);
209 const G4ErrorSymMatrix &m2);
211 const G4ErrorSymMatrix &m2);
213 const G4ErrorMatrix &m2);
214 friend G4ErrorMatrix operator*(const G4ErrorMatrix &m1,
215 const G4ErrorSymMatrix &m2);
216 // Multiply a Matrix by a Matrix or Vector.
217
218 // Returns v * v.T();
219
220 std::vector<G4double > m;
221
222 G4int nrow;
223 G4int size; // total number of elements
224
225 static G4ThreadLocal G4double posDefFraction5x5;
226 static G4ThreadLocal G4double adjustment5x5;
227 static const G4double CHOLESKY_THRESHOLD_5x5;
228 static const G4double CHOLESKY_CREEP_5x5;
229
230 static G4ThreadLocal G4double posDefFraction6x6;
231 static G4ThreadLocal G4double adjustment6x6;
232 static const G4double CHOLESKY_THRESHOLD_6x6;
233 static const G4double CHOLESKY_CREEP_6x6;
234
235 void invert4 (G4int & ifail);
236 void invert5 (G4int & ifail);
237 void invert6 (G4int & ifail);
238};
239
240//
241// Operations other than member functions for Matrix, G4ErrorSymMatrix,
242// DiagMatrix and Vectors
243//
244
245std::ostream& operator<<(std::ostream &s, const G4ErrorSymMatrix &q);
246// Write out Matrix, G4ErrorSymMatrix, DiagMatrix and Vector into ostream.
247
249 const G4ErrorSymMatrix &m2);
251 const G4ErrorMatrix &m2);
253 const G4ErrorSymMatrix &m2);
256// Multiplication operators.
257// Note that m *= m1 is always faster than m = m * m1
258
260// s = s1 / t. (s /= t is faster if you can use it.)
261
263 const G4ErrorSymMatrix &s2);
265 const G4ErrorMatrix &m2);
267 const G4ErrorSymMatrix &s2);
268// Addition operators
269
271 const G4ErrorSymMatrix &s2);
273 const G4ErrorMatrix &m2);
275 const G4ErrorSymMatrix &s2);
276// subtraction operators
277
279 const G4ErrorSymMatrix &s2);
280// Direct sum of two symmetric matrices;
281
283// Find the conditon number of a symmetric matrix.
284
287// Implicit symmetric QR step with Wilkinson Shift
288
290// Diagonalize a symmetric matrix.
291// It returns the matrix U so that s_old = U * s_diag * U.T()
292
294 G4int row=1, G4int col=1);
295// Finds and does Householder reflection on matrix.
296
299// Does a Householder tridiagonalization of a symmetric matrix.
300
301#include "G4ErrorSymMatrix.icc"
302
303#endif
void diag_step(G4ErrorSymMatrix *t, G4int begin, G4int end)
G4ErrorMatrix operator-(const G4ErrorMatrix &m1, const G4ErrorSymMatrix &s2)
G4ErrorSymMatrix operator/(const G4ErrorSymMatrix &m1, G4double t)
void tridiagonal(G4ErrorSymMatrix *a, G4ErrorMatrix *hsm)
G4ErrorMatrix diagonalize(G4ErrorSymMatrix *s)
void house_with_update2(G4ErrorSymMatrix *a, G4ErrorMatrix *v, G4int row=1, G4int col=1)
G4double condition(const G4ErrorSymMatrix &m)
G4ErrorSymMatrix dsum(const G4ErrorSymMatrix &s1, const G4ErrorSymMatrix &s2)
std::ostream & operator<<(std::ostream &s, const G4ErrorSymMatrix &q)
G4ErrorMatrix operator*(const G4ErrorMatrix &m1, const G4ErrorSymMatrix &m2)
G4ErrorMatrix operator+(const G4ErrorMatrix &m1, const G4ErrorSymMatrix &s2)
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4ErrorSymMatrix_row_const(const G4ErrorSymMatrix &, G4int)
const G4double & operator[](G4int) const
G4ErrorSymMatrix_row(G4ErrorSymMatrix &, G4int)
G4int num_row() const
G4ErrorSymMatrix similarityT(const G4ErrorMatrix &m1) const
friend void diag_step(G4ErrorSymMatrix *t, G4int begin, G4int end)
G4ErrorSymMatrix(G4ErrorSymMatrix &&)=default
void invertCholesky6(G4int &ifail)
void assign(const G4ErrorSymMatrix &m2)
friend G4ErrorSymMatrix operator+(const G4ErrorSymMatrix &m1, const G4ErrorSymMatrix &m2)
void invert(G4int &ifail)
void assign(const G4ErrorMatrix &m2)
G4double trace() const
const G4double & fast(G4int row, G4int col) const
void invertHaywood6(G4int &ifail)
friend void tridiagonal(G4ErrorSymMatrix *a, G4ErrorMatrix *hsm)
void invertHaywood5(G4int &ifail)
friend G4ErrorMatrix diagonalize(G4ErrorSymMatrix *s)
G4ErrorSymMatrix sub(G4int min_row, G4int max_row) const
G4ErrorSymMatrix apply(G4double(*f)(G4double, G4int, G4int)) const
G4ErrorSymMatrix operator-() const
const G4double & operator()(G4int row, G4int col) const
G4double & fast(G4int row, G4int col)
G4double determinant() const
G4ErrorSymMatrix & operator/=(G4double t)
G4ErrorSymMatrix inverse(G4int &ifail) const
G4ErrorSymMatrix & operator*=(G4double t)
friend G4double condition(const G4ErrorSymMatrix &m)
G4ErrorSymMatrix & operator=(G4ErrorSymMatrix &&)=default
void invertHaywood4(G4int &ifail)
G4ErrorSymMatrix & operator+=(const G4ErrorSymMatrix &m2)
virtual ~G4ErrorSymMatrix()
G4int num_size() const
friend void house_with_update2(G4ErrorSymMatrix *a, G4ErrorMatrix *v, G4int row, G4int col)
G4ErrorSymMatrix_row operator[](G4int)
G4ErrorSymMatrix & operator-=(const G4ErrorSymMatrix &m2)
G4ErrorSymMatrix similarity(const G4ErrorMatrix &m1) const
void invertCholesky5(G4int &ifail)
G4int num_col() const
friend void diag_step(G4ErrorSymMatrix *t, G4ErrorMatrix *u, G4int begin, G4int end)
G4ErrorSymMatrix T() const
G4ErrorSymMatrix & operator=(const G4ErrorSymMatrix &m2)
friend G4ErrorMatrix operator*(const G4ErrorSymMatrix &m1, const G4ErrorSymMatrix &m2)
G4double & operator()(G4int row, G4int col)
void invertBunchKaufman(G4int &ifail)
#define G4ThreadLocal
Definition: tls.hh:77