Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4ErrorMatrix.hh
이 파일의 문서화 페이지로 가기
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: G4ErrorMatrix.hh 108074 2017-12-19 15:35:08Z gcosmo $
27 //
28 // Class Description:
29 //
30 // Simplified version of CLHEP HepMatrix class
31 
32 // History:
33 // - Imported from CLHEP and modified: P. Arce May 2007
34 // --------------------------------------------------------------------
35 
36 #ifndef G4ErrorMatrix_hh
37 #define G4ErrorMatrix_hh
38 
39 #include <vector>
40 #include "G4Types.hh"
41 
43 
44 typedef std::vector<G4double >::iterator G4ErrorMatrixIter;
45 typedef std::vector<G4double >::const_iterator G4ErrorMatrixConstIter;
46 
48 {
49 
50  public: // with description
51 
52  G4ErrorMatrix();
53  // Default constructor. Gives 0 x 0 matrix.
54  // Another G4ErrorMatrix can be assigned to it.
55 
57  // Constructor. Gives an unitialized p x q matrix.
58 
59  G4ErrorMatrix(G4int p, G4int q, G4int i);
60  // Constructor. Gives an initialized p x q matrix.
61  // If i=0, it is initialized to all 0. If i=1, the diagonal elements
62  // are set to 1.0.
63 
64  G4ErrorMatrix(const G4ErrorMatrix &m1);
65  // Copy 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 
99  // Assignment operators.
100 
101  G4ErrorMatrix operator- () const;
102  // unary minus, ie. flip the sign of each element.
103 
105  // Apply a function to all elements of the matrix.
106 
107  G4ErrorMatrix T() const;
108  // Returns the transpose of a G4ErrorMatrix.
109 
110  G4ErrorMatrix sub(G4int min_row, G4int max_row,
111  G4int min_col, G4int max_col) const;
112  // Returns a sub matrix of a G4ErrorMatrix.
113  // WARNING: rows and columns are numbered from 1
114 
115  void sub(G4int row, G4int col, const G4ErrorMatrix &m1);
116  // Sub matrix of this G4ErrorMatrix is replaced with m1.
117  // WARNING: rows and columns are numbered from 1
118 
119  inline G4ErrorMatrix inverse(G4int& ierr) const;
120  // Invert a G4ErrorMatrix. G4ErrorMatrix must be square and is not changed.
121  // Returns ierr = 0 (zero) when successful, otherwise non-zero.
122 
123  virtual void invert(G4int& ierr);
124  // Invert a G4ErrorMatrix. G4ErrorMatrix must be square.
125  // N.B. the contents of the matrix are replaced by the inverse.
126  // Returns ierr = 0 (zero) when successful, otherwise non-zero.
127  // This method has less overhead then inverse().
128 
129  G4double determinant() const;
130  // calculate the determinant of the matrix.
131 
132  G4double trace() const;
133  // calculate the trace of the matrix (sum of diagonal elements).
134 
136  {
137  typedef std::vector<G4double >::const_iterator G4ErrorMatrixConstIter;
138  public:
141  private:
144  };
146  {
147  public:
149  const G4double & operator[](G4int) const;
150  private:
153  };
154  // helper classes for implementing m[i][j]
155 
157  inline const G4ErrorMatrix_row_const operator[] (G4int) const;
158  // Read or write a matrix element.
159  // While it may not look like it, you simply do m[i][j] to get an element.
160  // ** Note that the indexing starts from [0][0]. **
161 
162  protected:
163 
164  virtual inline G4int num_size() const;
165  virtual void invertHaywood4(G4int& ierr);
166  virtual void invertHaywood5(G4int& ierr);
167  virtual void invertHaywood6(G4int& ierr);
168 
169  public:
170 
171  static void error(const char *s);
172 
173  private:
174 
175  friend class G4ErrorMatrix_row;
177  friend class G4ErrorSymMatrix;
178  // Friend classes.
179 
180  friend G4ErrorMatrix operator+(const G4ErrorMatrix &m1,
181  const G4ErrorMatrix &m2);
182  friend G4ErrorMatrix operator-(const G4ErrorMatrix &m1,
183  const G4ErrorMatrix &m2);
184  friend G4ErrorMatrix operator*(const G4ErrorMatrix &m1,
185  const G4ErrorMatrix &m2);
186  friend G4ErrorMatrix operator*(const G4ErrorMatrix &m1,
187  const G4ErrorSymMatrix &m2);
188  friend G4ErrorMatrix operator*(const G4ErrorSymMatrix &m1,
189  const G4ErrorMatrix &m2);
190  friend G4ErrorMatrix operator*(const G4ErrorSymMatrix &m1,
191  const G4ErrorSymMatrix &m2);
192  // Multiply a G4ErrorMatrix by a G4ErrorMatrix or Vector.
193 
194  // solve the system of linear eq
195 
197  friend void tridiagonal(G4ErrorSymMatrix *a,G4ErrorMatrix *hsm);
198  friend void row_house(G4ErrorMatrix *,const G4ErrorMatrix &, G4double,
199  G4int, G4int, G4int, G4int);
200  friend void back_solve(const G4ErrorMatrix &R, G4ErrorMatrix *b);
201  friend void col_givens(G4ErrorMatrix *A, G4double c,
202  G4double s, G4int k1, G4int k2,
203  G4int rowmin, G4int rowmax);
204  // Does a column Givens update.
205 
206  friend void row_givens(G4ErrorMatrix *A, G4double c,
207  G4double s, G4int k1, G4int k2,
208  G4int colmin, G4int colmax);
209  friend void col_house(G4ErrorMatrix *,const G4ErrorMatrix &, G4double,
210  G4int, G4int, G4int, G4int);
211  friend void house_with_update(G4ErrorMatrix *a,G4int row,G4int col);
213  G4int row, G4int col);
215  G4int row, G4int col);
216 
217  G4int dfact_matrix(G4double &det, G4int *ir);
218  // factorize the matrix. If successful, the return code is 0. On
219  // return, det is the determinant and ir[] is row-interchange
220  // matrix. See CERNLIB's DFACT routine.
221 
222  G4int dfinv_matrix(G4int *ir);
223  // invert the matrix. See CERNLIB DFINV.
224 
225  std::vector<G4double > m;
226 
229 };
230 
231 // Operations other than member functions for G4ErrorMatrix
232 
236 // Multiplication operators
237 // Note that m *= m1 is always faster than m = m * m1.
238 
240 // m = m1 / t. (m /= t is faster if you can use it.)
241 
243 // m = m1 + m2;
244 // Note that m += m1 is always faster than m = m + m1.
245 
247 // m = m1 - m2;
248 // Note that m -= m1 is always faster than m = m - m1.
249 
251 // Direct sum of two matrices. The direct sum of A and B is the matrix
252 // A 0
253 // 0 B
254 
255 std::ostream& operator<<(std::ostream &s, const G4ErrorMatrix &q);
256 // Read in, write out G4ErrorMatrix into a stream.
257 
258 //
259 // Specialized linear algebra functions
260 //
261 
264 // Works like backsolve, except matrix does not need to be upper
265 // triangular. For nonsquare matrix, it solves in the least square sense.
266 
269 // Finds the inverse of a matrix using QR decomposition. Note, often what
270 // you really want is solve or backsolve, they can be much quicker than
271 // inverse in many calculations.
272 
273 
276 // Does a QR decomposition of a matrix.
277 
278 void back_solve(const G4ErrorMatrix &R, G4ErrorMatrix *b);
279 // Solves R*x = b where R is upper triangular. Also has a variation that
280 // solves a number of equations of this form in one step, where b is a matrix
281 // with each column a different vector. See also solve.
282 
283 void col_house(G4ErrorMatrix *a, const G4ErrorMatrix &v, G4double vnormsq,
284  G4int row, G4int col, G4int row_start, G4int col_start);
285 void col_house(G4ErrorMatrix *a, const G4ErrorMatrix &v, G4int row, G4int col,
286  G4int row_start, G4int col_start);
287 // Does a column Householder update.
288 
290  G4int k1, G4int k2, G4int row_min=1, G4int row_max=0);
291 // do a column Givens update
292 
294  G4int k1, G4int k2, G4int col_min=1, G4int col_max=0);
295 // do a row Givens update
296 
297 // void givens(G4double a, G4double b, G4double *c, G4double *s);
298 // algorithm 5.1.5 in Golub and Van Loan
299 
300 // Returns a Householder vector to zero elements.
301 
304  G4int row=1, G4int col=1);
305 // Finds and does Householder reflection on matrix.
306 
307 void row_house(G4ErrorMatrix *a, const G4ErrorMatrix &v, G4double vnormsq,
308  G4int row, G4int col, G4int row_start, G4int col_start);
309 void row_house(G4ErrorMatrix *a, const G4ErrorMatrix &v, G4int row, G4int col,
310  G4int row_start, G4int col_start);
311 // Does a row Householder update.
312 
313 #include "G4ErrorMatrix.icc"
314 
315 #endif
friend G4ErrorMatrix operator+(const G4ErrorMatrix &m1, const G4ErrorMatrix &m2)
G4double trace() const
friend G4ErrorMatrix operator*(const G4ErrorMatrix &m1, const G4ErrorMatrix &m2)
G4ErrorMatrix_row operator[](G4int)
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
friend void col_givens(G4ErrorMatrix *A, G4double c, G4double s, G4int k1, G4int k2, G4int rowmin, G4int rowmax)
G4ErrorMatrix inverse(G4int &ierr) const
virtual G4int num_row() const
std::ostream & operator<<(std::ostream &, const BasicVector3D< float > &)
G4ErrorMatrix & operator=(const G4ErrorMatrix &m2)
G4ErrorMatrix & operator*=(G4double t)
G4ErrorMatrix qr_solve(const G4ErrorMatrix &A, const G4ErrorMatrix &b)
static constexpr double m2
Definition: G4SIunits.hh:130
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)
friend void house_with_update(G4ErrorMatrix *a, G4int row, G4int col)
const char * p
Definition: xmltok.h:285
G4ErrorMatrix T() const
void col_house(G4ErrorMatrix *a, const G4ErrorMatrix &v, G4double vnormsq, G4int row, G4int col, G4int row_start, G4int col_start)
virtual void invertHaywood6(G4int &ierr)
G4ErrorMatrix_row_const(const G4ErrorMatrix &, G4int)
G4ErrorMatrix sub(G4int min_row, G4int max_row, G4int min_col, G4int max_col) const
friend G4ErrorMatrix qr_solve(G4ErrorMatrix *, const G4ErrorMatrix &b)
virtual void invertHaywood5(G4int &ierr)
BasicVector3D< float > operator*(const BasicVector3D< float > &v, double a)
static void error(const char *s)
const XML_Char * s
Definition: expat.h:262
double G4double
Definition: G4Types.hh:76
void house_with_update(G4ErrorMatrix *a, G4int row=1, G4int col=1)
std::vector< G4double >::const_iterator G4ErrorMatrixConstIter
G4int dfinv_matrix(G4int *ir)
virtual void invertHaywood4(G4int &ierr)
G4ErrorMatrix dsum(const G4ErrorMatrix &, const G4ErrorMatrix &)
friend void back_solve(const G4ErrorMatrix &R, G4ErrorMatrix *b)
G4ErrorMatrix apply(G4double(*f)(G4double, G4int, G4int)) const
Int_t col[ntarg]
Definition: Style.C:29
BasicVector3D< float > operator+(const BasicVector3D< float > &v)
void col_givens(G4ErrorMatrix *A, G4double c, G4double s, G4int k1, G4int k2, G4int row_min=1, G4int row_max=0)
friend void tridiagonal(G4ErrorSymMatrix *a, G4ErrorMatrix *hsm)
Double_t R
friend void row_house(G4ErrorMatrix *, const G4ErrorMatrix &, G4double, G4int, G4int, G4int, G4int)
double A(double temperature)
virtual G4int num_size() const
std::vector< G4double >::const_iterator G4ErrorMatrixConstIter
virtual ~G4ErrorMatrix()
std::vector< G4double >::iterator G4ErrorMatrixIter
int G4int
Definition: G4Types.hh:78
G4double determinant() const
std::vector< G4double > m
void back_solve(const G4ErrorMatrix &R, G4ErrorMatrix *b)
G4ErrorMatrix & operator/=(G4double t)
G4ErrorMatrix operator-() const
friend void col_house(G4ErrorMatrix *, const G4ErrorMatrix &, G4double, G4int, G4int, G4int, G4int)
G4ErrorMatrix & operator+=(const G4ErrorMatrix &m2)
G4ErrorMatrix_row(G4ErrorMatrix &, G4int)
friend void house_with_update2(G4ErrorSymMatrix *a, G4ErrorMatrix *v, G4int row, G4int col)
BasicVector3D< float > operator/(const BasicVector3D< float > &v, double a)
void row_house(G4ErrorMatrix *a, const G4ErrorMatrix &v, G4double vnormsq, G4int row, G4int col, G4int row_start, G4int col_start)
virtual const G4double & operator()(G4int row, G4int col) const
void qr_decomp(G4ErrorMatrix *A, G4ErrorMatrix *hsm)
virtual G4int num_col() const
G4ErrorMatrix & operator-=(const G4ErrorMatrix &m2)
friend void row_givens(G4ErrorMatrix *A, G4double c, G4double s, G4int k1, G4int k2, G4int colmin, G4int colmax)
G4int dfact_matrix(G4double &det, G4int *ir)
virtual void invert(G4int &ierr)
BasicVector3D< float > operator-(const BasicVector3D< float > &v)
const G4double & operator[](G4int) const