Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ErrorMatrix.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 HepMatrix class
30
31// History:
32// - Imported from CLHEP and modified: P. Arce May 2007
33// --------------------------------------------------------------------
34
35#ifndef G4ErrorMatrix_hh
36#define G4ErrorMatrix_hh
37
38#include <vector>
39#include "G4Types.hh"
40
42
43typedef std::vector<G4double >::iterator G4ErrorMatrixIter;
44typedef std::vector<G4double >::const_iterator G4ErrorMatrixConstIter;
45
47{
48
49 public: // with description
50
52 // Default constructor. Gives 0 x 0 matrix.
53 // Another G4ErrorMatrix can be assigned to it.
54
56 // Constructor. Gives an unitialized p x q matrix.
57
59 // Constructor. Gives an initialized p x q matrix.
60 // If i=0, it is initialized to all 0. If i=1, the diagonal elements
61 // are set to 1.0.
62
65 // Copy and move constructor.
66
68 // Constructors from G4ErrorSymG4ErrorMatrix, DiagG4ErrorMatrix and Vector.
69
70 virtual ~G4ErrorMatrix();
71 // Destructor.
72
73 inline virtual G4int num_row() const;
74 // Returns the number of rows.
75
76 inline virtual G4int num_col() const;
77 // Returns the number of columns.
78
79 inline virtual const G4double & operator()(G4int row, G4int col) const;
80 inline virtual G4double & operator()(G4int row, G4int col);
81 // Read or write a matrix element.
82 // ** Note that the indexing starts from (1,1). **
83
85 // Multiply a G4ErrorMatrix by a floating number.
86
88 // Divide a G4ErrorMatrix by a floating number.
89
94 // Add or subtract a G4ErrorMatrix.
95 // When adding/subtracting Vector, G4ErrorMatrix must have num_col of one.
96
100 // Assignment and move operators.
101
103 // unary minus, ie. flip the sign of each element.
104
106 // Apply a function to all elements of the matrix.
107
108 G4ErrorMatrix T() const;
109 // Returns the transpose of a G4ErrorMatrix.
110
111 G4ErrorMatrix sub(G4int min_row, G4int max_row,
112 G4int min_col, G4int max_col) const;
113 // Returns a sub matrix of a G4ErrorMatrix.
114 // WARNING: rows and columns are numbered from 1
115
116 void sub(G4int row, G4int col, const G4ErrorMatrix &m1);
117 // Sub matrix of this G4ErrorMatrix is replaced with m1.
118 // WARNING: rows and columns are numbered from 1
119
120 inline G4ErrorMatrix inverse(G4int& ierr) const;
121 // Invert a G4ErrorMatrix. G4ErrorMatrix must be square and is not changed.
122 // Returns ierr = 0 (zero) when successful, otherwise non-zero.
123
124 virtual void invert(G4int& ierr);
125 // Invert a G4ErrorMatrix. G4ErrorMatrix must be square.
126 // N.B. the contents of the matrix are replaced by the inverse.
127 // Returns ierr = 0 (zero) when successful, otherwise non-zero.
128 // This method has less overhead then inverse().
129
130 G4double determinant() const;
131 // calculate the determinant of the matrix.
132
133 G4double trace() const;
134 // calculate the trace of the matrix (sum of diagonal elements).
135
137 {
138 typedef std::vector<G4double >::const_iterator G4ErrorMatrixConstIter;
139 public:
142 private:
143 G4ErrorMatrix& _a;
144 G4int _r;
145 };
147 {
148 public:
150 const G4double & operator[](G4int) const;
151 private:
152 const G4ErrorMatrix& _a;
153 G4int _r;
154 };
155 // helper classes for implementing m[i][j]
156
159 // Read or write a matrix element.
160 // While it may not look like it, you simply do m[i][j] to get an element.
161 // ** Note that the indexing starts from [0][0]. **
162
163 protected:
164
165 virtual inline G4int num_size() const;
166 virtual void invertHaywood4(G4int& ierr);
167 virtual void invertHaywood5(G4int& ierr);
168 virtual void invertHaywood6(G4int& ierr);
169
170 public:
171
172 static void error(const char *s);
173
174 private:
175
176 friend class G4ErrorMatrix_row;
178 friend class G4ErrorSymMatrix;
179 // Friend classes.
180
181 friend G4ErrorMatrix operator+(const G4ErrorMatrix &m1,
182 const G4ErrorMatrix &m2);
183 friend G4ErrorMatrix operator-(const G4ErrorMatrix &m1,
184 const G4ErrorMatrix &m2);
185 friend G4ErrorMatrix operator*(const G4ErrorMatrix &m1,
186 const G4ErrorMatrix &m2);
187 friend G4ErrorMatrix operator*(const G4ErrorMatrix &m1,
188 const G4ErrorSymMatrix &m2);
190 const G4ErrorMatrix &m2);
192 const G4ErrorSymMatrix &m2);
193 // Multiply a G4ErrorMatrix by a G4ErrorMatrix or Vector.
194
195 // solve the system of linear eq
196
201 friend void back_solve(const G4ErrorMatrix &R, G4ErrorMatrix *b);
203 G4double s, G4int k1, G4int k2,
204 G4int rowmin, G4int rowmax);
205 // Does a column Givens update.
206
208 G4double s, G4int k1, G4int k2,
209 G4int colmin, G4int colmax);
214 G4int row, G4int col);
216 G4int row, G4int col);
217
218 G4int dfact_matrix(G4double &det, G4int *ir);
219 // factorize the matrix. If successful, the return code is 0. On
220 // return, det is the determinant and ir[] is row-interchange
221 // matrix. See CERNLIB's DFACT routine.
222
223 G4int dfinv_matrix(G4int *ir);
224 // invert the matrix. See CERNLIB DFINV.
225
226 std::vector<G4double > m;
227
228 G4int nrow, ncol;
229 G4int size;
230};
231
232// Operations other than member functions for G4ErrorMatrix
233
237// Multiplication operators
238// Note that m *= m1 is always faster than m = m * m1.
239
241// m = m1 / t. (m /= t is faster if you can use it.)
242
244// m = m1 + m2;
245// Note that m += m1 is always faster than m = m + m1.
246
248// m = m1 - m2;
249// Note that m -= m1 is always faster than m = m - m1.
250
252// Direct sum of two matrices. The direct sum of A and B is the matrix
253// A 0
254// 0 B
255
256std::ostream& operator<<(std::ostream &s, const G4ErrorMatrix &q);
257// Read in, write out G4ErrorMatrix into a stream.
258
259//
260// Specialized linear algebra functions
261//
262
265// Works like backsolve, except matrix does not need to be upper
266// triangular. For nonsquare matrix, it solves in the least square sense.
267
270// Finds the inverse of a matrix using QR decomposition. Note, often what
271// you really want is solve or backsolve, they can be much quicker than
272// inverse in many calculations.
273
274
277// Does a QR decomposition of a matrix.
278
280// Solves R*x = b where R is upper triangular. Also has a variation that
281// solves a number of equations of this form in one step, where b is a matrix
282// with each column a different vector. See also solve.
283
285 G4int row, G4int col, G4int row_start, G4int col_start);
286void col_house(G4ErrorMatrix *a, const G4ErrorMatrix &v, G4int row, G4int col,
287 G4int row_start, G4int col_start);
288// Does a column Householder update.
289
291 G4int k1, G4int k2, G4int row_min=1, G4int row_max=0);
292// do a column Givens update
293
295 G4int k1, G4int k2, G4int col_min=1, G4int col_max=0);
296// do a row Givens update
297
298// void givens(G4double a, G4double b, G4double *c, G4double *s);
299// algorithm 5.1.5 in Golub and Van Loan
300
301// Returns a Householder vector to zero elements.
302
305 G4int row=1, G4int col=1);
306// Finds and does Householder reflection on matrix.
307
309 G4int row, G4int col, G4int row_start, G4int col_start);
310void row_house(G4ErrorMatrix *a, const G4ErrorMatrix &v, G4int row, G4int col,
311 G4int row_start, G4int col_start);
312// Does a row Householder update.
313
314#include "G4ErrorMatrix.icc"
315
316#endif
void qr_decomp(G4ErrorMatrix *A, G4ErrorMatrix *hsm)
void col_house(G4ErrorMatrix *a, const G4ErrorMatrix &v, G4double vnormsq, G4int row, G4int col, G4int row_start, G4int col_start)
G4ErrorMatrix qr_solve(const G4ErrorMatrix &A, const G4ErrorMatrix &b)
std::vector< G4double >::iterator G4ErrorMatrixIter
G4ErrorMatrix qr_inverse(const G4ErrorMatrix &A)
void row_givens(G4ErrorMatrix *A, G4double c, G4double s, G4int k1, G4int k2, G4int col_min=1, G4int col_max=0)
G4ErrorMatrix operator-(const G4ErrorMatrix &m1, const G4ErrorMatrix &m2)
void col_givens(G4ErrorMatrix *A, G4double c, G4double s, G4int k1, G4int k2, G4int row_min=1, G4int row_max=0)
G4ErrorMatrix dsum(const G4ErrorMatrix &, const G4ErrorMatrix &)
G4ErrorMatrix operator*(const G4ErrorMatrix &m1, const G4ErrorMatrix &m2)
G4ErrorMatrix operator+(const G4ErrorMatrix &m1, const G4ErrorMatrix &m2)
void row_house(G4ErrorMatrix *a, const G4ErrorMatrix &v, G4double vnormsq, G4int row, G4int col, G4int row_start, G4int col_start)
void house_with_update(G4ErrorMatrix *a, G4int row=1, G4int col=1)
std::ostream & operator<<(std::ostream &s, const G4ErrorMatrix &q)
std::vector< G4double >::const_iterator G4ErrorMatrixConstIter
void back_solve(const G4ErrorMatrix &R, G4ErrorMatrix *b)
G4ErrorMatrix operator/(const G4ErrorMatrix &m1, G4double t)
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4double & operator[](G4int) const
G4ErrorMatrix_row_const(const G4ErrorMatrix &, G4int)
G4ErrorMatrix_row(G4ErrorMatrix &, G4int)
friend void row_house(G4ErrorMatrix *, const G4ErrorMatrix &, G4double, G4int, G4int, G4int, G4int)
G4ErrorMatrix apply(G4double(*f)(G4double, G4int, G4int)) const
virtual void invertHaywood4(G4int &ierr)
G4ErrorMatrix operator-() const
friend void col_house(G4ErrorMatrix *, const G4ErrorMatrix &, G4double, G4int, G4int, G4int, G4int)
friend void row_givens(G4ErrorMatrix *A, G4double c, G4double s, G4int k1, G4int k2, G4int colmin, G4int colmax)
virtual void invert(G4int &ierr)
G4ErrorMatrix & operator/=(G4double t)
friend void house_with_update(G4ErrorMatrix *a, G4ErrorMatrix *v, G4int row, G4int col)
friend void tridiagonal(G4ErrorSymMatrix *a, G4ErrorMatrix *hsm)
G4ErrorMatrix T() const
friend void house_with_update(G4ErrorMatrix *a, G4int row, G4int col)
virtual ~G4ErrorMatrix()
G4double determinant() const
friend G4ErrorMatrix operator*(const G4ErrorMatrix &m1, const G4ErrorMatrix &m2)
friend G4ErrorMatrix operator+(const G4ErrorMatrix &m1, const G4ErrorMatrix &m2)
virtual void invertHaywood5(G4int &ierr)
G4ErrorMatrix & operator=(const G4ErrorMatrix &m2)
virtual G4int num_col() const
virtual void invertHaywood6(G4int &ierr)
G4ErrorMatrix & operator*=(G4double t)
G4ErrorMatrix & operator-=(const G4ErrorMatrix &m2)
virtual G4int num_size() const
G4ErrorMatrix inverse(G4int &ierr) const
virtual G4int num_row() const
static void error(const char *s)
G4ErrorMatrix_row operator[](G4int)
G4ErrorMatrix(G4ErrorMatrix &&)=default
G4ErrorMatrix & operator+=(const G4ErrorMatrix &m2)
friend void house_with_update2(G4ErrorSymMatrix *a, G4ErrorMatrix *v, G4int row, G4int col)
G4double trace() const
friend void col_givens(G4ErrorMatrix *A, G4double c, G4double s, G4int k1, G4int k2, G4int rowmin, G4int rowmax)
virtual const G4double & operator()(G4int row, G4int col) const
virtual G4double & operator()(G4int row, G4int col)
G4ErrorMatrix sub(G4int min_row, G4int max_row, G4int min_col, G4int max_col) const
friend void back_solve(const G4ErrorMatrix &R, G4ErrorMatrix *b)
friend G4ErrorMatrix qr_solve(G4ErrorMatrix *, const G4ErrorMatrix &b)