Geant4 9.6.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// $Id$
27//
28// Class Description:
29//
30// Simplified version of CLHEP HepSymMatrix class.
31
32// History:
33// - Imported from CLHEP and modified: P. Arce May 2007
34// --------------------------------------------------------------------
35
36#ifndef G4ErrorSymMatrix_hh
37#define G4ErrorSymMatrix_hh
38
39#include <vector>
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
58 // Copy constructor.
59
60 // Constructor from DiagMatrix
61
62 virtual ~G4ErrorSymMatrix();
63 // Destructor.
64
65 inline G4int num_row() const;
66 inline G4int num_col() const;
67 // Returns number of rows/columns.
68
69 const G4double & operator()(G4int row, G4int col) const;
71 // Read and write a G4ErrorSymMatrix element.
72 // ** Note that indexing starts from (1,1). **
73
74 const G4double & fast(G4int row, G4int col) const;
75 G4double & fast(G4int row, G4int col);
76 // fast element access.
77 // Must be row>=col;
78 // ** Note that indexing starts from (1,1). **
79
80 void assign(const G4ErrorMatrix &m2);
81 // Assigns m2 to s, assuming m2 is a symmetric matrix.
82
83 void assign(const G4ErrorSymMatrix &m2);
84 // Another form of assignment. For consistency.
85
87 // Multiply a G4ErrorSymMatrix by a floating number.
88
90 // Divide a G4ErrorSymMatrix by a floating number.
91
94 // Add or subtract a G4ErrorSymMatrix.
95
97 // Assignment operators. Notice that there is no G4ErrorSymMatrix = Matrix.
98
100 // unary minus, ie. flip the sign of each element.
101
103 // Returns the transpose of a G4ErrorSymMatrix (which is itself).
104
106 // Apply a function to all elements of the matrix.
107
110 // Returns m1*s*m1.T().
111
113 // temporary. test of new similarity.
114 // Returns m1.T()*s*m1.
115
116 G4ErrorSymMatrix sub(G4int min_row, G4int max_row) const;
117 // Returns a sub matrix of a G4ErrorSymMatrix.
118
119 void sub(G4int row, const G4ErrorSymMatrix &m1);
120 // Sub matrix of this G4ErrorSymMatrix is replaced with m1.
121
122 G4ErrorSymMatrix sub(G4int min_row, G4int max_row);
123 // SGI CC bug. I have to have both with/without const. I should not need
124 // one without const.
125
126 inline G4ErrorSymMatrix inverse(G4int &ifail) const;
127 // Invert a Matrix. The matrix is not changed
128 // Returns 0 when successful, otherwise non-zero.
129
130 void invert(G4int &ifail);
131 // Invert a Matrix.
132 // N.B. the contents of the matrix are replaced by the inverse.
133 // Returns ierr = 0 when successful, otherwise non-zero.
134 // This method has less overhead then inverse().
135
136 G4double determinant() const;
137 // calculate the determinant of the matrix.
138
139 G4double trace() const;
140 // calculate the trace of the matrix (sum of diagonal elements).
141
143 {
144 public:
147 private:
149 G4int _r;
150 };
152 {
153 public:
155 inline const G4double & operator[](G4int) const;
156 private:
157 const G4ErrorSymMatrix& _a;
158 G4int _r;
159 };
160 // helper class to implement m[i][j]
161
164 // Read or write a matrix element.
165 // While it may not look like it, you simply do m[i][j] to get an
166 // element.
167 // ** Note that the indexing starts from [0][0]. **
168
169 // Special-case inversions for 5x5 and 6x6 symmetric positive definite:
170 // These set ifail=0 and invert if the matrix was positive definite;
171 // otherwise ifail=1 and the matrix is left unaltered.
172
173 void invertCholesky5 (G4int &ifail);
174 void invertCholesky6 (G4int &ifail);
175
176 // Inversions for 5x5 and 6x6 forcing use of specific methods: The
177 // behavior (though not the speed) will be identical to invert(ifail).
178
179 void invertHaywood4 (G4int & ifail);
180 void invertHaywood5 (G4int &ifail);
181 void invertHaywood6 (G4int &ifail);
182 void invertBunchKaufman (G4int &ifail);
183
184 protected:
185
186 inline G4int num_size() const;
187
188 private:
189
192 friend class G4ErrorMatrix;
193
196 friend void diag_step(G4ErrorSymMatrix *t, G4int begin, G4int end);
198 G4int begin, G4int end);
201 G4int row, G4int col);
202
204 const G4ErrorSymMatrix &m2);
206 const G4ErrorSymMatrix &m2);
208 const G4ErrorSymMatrix &m2);
210 const G4ErrorMatrix &m2);
211 friend G4ErrorMatrix operator*(const G4ErrorMatrix &m1,
212 const G4ErrorSymMatrix &m2);
213 // Multiply a Matrix by a Matrix or Vector.
214
215 // Returns v * v.T();
216
217 std::vector<G4double > m;
218
219 G4int nrow;
220 G4int size; // total number of elements
221
222 static G4double posDefFraction5x5;
223 static G4double adjustment5x5;
224 static const G4double CHOLESKY_THRESHOLD_5x5;
225 static const G4double CHOLESKY_CREEP_5x5;
226
227 static G4double posDefFraction6x6;
228 static G4double adjustment6x6;
229 static const G4double CHOLESKY_THRESHOLD_6x6;
230 static const G4double CHOLESKY_CREEP_6x6;
231
232 void invert4 (G4int & ifail);
233 void invert5 (G4int & ifail);
234 void invert6 (G4int & ifail);
235};
236
237//
238// Operations other than member functions for Matrix, G4ErrorSymMatrix,
239// DiagMatrix and Vectors
240//
241
242std::ostream& operator<<(std::ostream &s, const G4ErrorSymMatrix &q);
243// Write out Matrix, G4ErrorSymMatrix, DiagMatrix and Vector into ostream.
244
246 const G4ErrorSymMatrix &m2);
248 const G4ErrorMatrix &m2);
250 const G4ErrorSymMatrix &m2);
253// Multiplication operators.
254// Note that m *= m1 is always faster than m = m * m1
255
257// s = s1 / t. (s /= t is faster if you can use it.)
258
260 const G4ErrorSymMatrix &s2);
262 const G4ErrorMatrix &m2);
264 const G4ErrorSymMatrix &s2);
265// Addition operators
266
268 const G4ErrorSymMatrix &s2);
270 const G4ErrorMatrix &m2);
272 const G4ErrorSymMatrix &s2);
273// subtraction operators
274
276 const G4ErrorSymMatrix &s2);
277// Direct sum of two symmetric matrices;
278
280// Find the conditon number of a symmetric matrix.
281
284// Implicit symmetric QR step with Wilkinson Shift
285
287// Diagonalize a symmetric matrix.
288// It returns the matrix U so that s_old = U * s_diag * U.T()
289
291 G4int row=1, G4int col=1);
292// Finds and does Householder reflection on matrix.
293
296// Does a Householder tridiagonalization of a symmetric matrix.
297
298#include "G4ErrorSymMatrix.icc"
299
300#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:64
int G4int
Definition: G4Types.hh:66
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)
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)
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)