G4ErrorMatrix.hh

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id: G4ErrorMatrix.hh 69766 2013-05-14 14:33:55Z gcosmo $
00027 //
00028 // Class Description:
00029 //
00030 // Simplified version of CLHEP HepMatrix class
00031 
00032 // History:
00033 // - Imported from CLHEP and modified:   P. Arce    May 2007
00034 // --------------------------------------------------------------------
00035 
00036 #ifndef G4ErrorMatrix_hh
00037 #define G4ErrorMatrix_hh
00038 
00039 #include <vector>
00040 
00041 class G4ErrorSymMatrix;
00042 
00043 typedef std::vector<G4double >::iterator G4ErrorMatrixIter;
00044 typedef std::vector<G4double >::const_iterator G4ErrorMatrixConstIter;
00045 
00046 class G4ErrorMatrix
00047 {
00048 
00049  public:  // with description
00050 
00051    G4ErrorMatrix();
00052    // Default constructor. Gives 0 x 0 matrix.
00053    // Another G4ErrorMatrix can be assigned to it.
00054 
00055    G4ErrorMatrix(G4int p, G4int q);
00056    // Constructor. Gives an unitialized p x q matrix.
00057 
00058    G4ErrorMatrix(G4int p, G4int q, G4int i);
00059    // Constructor. Gives an initialized p x q matrix. 
00060    // If i=0, it is initialized to all 0. If i=1, the diagonal elements
00061    // are set to 1.0.
00062 
00063    G4ErrorMatrix(const G4ErrorMatrix &m1);
00064    // Copy constructor.
00065 
00066    G4ErrorMatrix(const G4ErrorSymMatrix &m1);
00067    // Constructors from G4ErrorSymG4ErrorMatrix, DiagG4ErrorMatrix and Vector.
00068 
00069    virtual ~G4ErrorMatrix();
00070    // Destructor.
00071 
00072    inline virtual G4int num_row() const;
00073    // Returns the number of rows.
00074 
00075    inline virtual G4int num_col() const;
00076    // Returns the number of columns.
00077 
00078    inline virtual const G4double & operator()(G4int row, G4int col) const;
00079    inline virtual G4double & operator()(G4int row, G4int col);
00080    // Read or write a matrix element. 
00081    // ** Note that the indexing starts from (1,1). **
00082 
00083    G4ErrorMatrix & operator *= (G4double t);
00084    // Multiply a G4ErrorMatrix by a floating number.
00085 
00086    G4ErrorMatrix & operator /= (G4double t); 
00087    // Divide a G4ErrorMatrix by a floating number.
00088 
00089    G4ErrorMatrix & operator += ( const G4ErrorMatrix &m2);
00090    G4ErrorMatrix & operator += ( const G4ErrorSymMatrix &m2);
00091    G4ErrorMatrix & operator -= ( const G4ErrorMatrix &m2);
00092    G4ErrorMatrix & operator -= ( const G4ErrorSymMatrix &m2);
00093    // Add or subtract a G4ErrorMatrix. 
00094    // When adding/subtracting Vector, G4ErrorMatrix must have num_col of one.
00095 
00096    G4ErrorMatrix & operator = ( const G4ErrorMatrix &m2);
00097    G4ErrorMatrix & operator = ( const G4ErrorSymMatrix &m2);
00098    // Assignment operators.
00099 
00100    G4ErrorMatrix operator- () const;
00101    // unary minus, ie. flip the sign of each element.
00102 
00103    G4ErrorMatrix apply(G4double (*f)(G4double, G4int, G4int)) const;
00104    // Apply a function to all elements of the matrix.
00105 
00106    G4ErrorMatrix T() const;
00107    // Returns the transpose of a G4ErrorMatrix.
00108 
00109    G4ErrorMatrix sub(G4int min_row, G4int max_row,
00110                      G4int min_col, G4int max_col) const;
00111    // Returns a sub matrix of a G4ErrorMatrix.
00112    // WARNING: rows and columns are numbered from 1
00113 
00114    void sub(G4int row, G4int col, const G4ErrorMatrix &m1);
00115    // Sub matrix of this G4ErrorMatrix is replaced with m1.
00116    // WARNING: rows and columns are numbered from 1
00117 
00118    inline G4ErrorMatrix inverse(G4int& ierr) const;
00119    // Invert a G4ErrorMatrix. G4ErrorMatrix must be square and is not changed.
00120    // Returns ierr = 0 (zero) when successful, otherwise non-zero.
00121 
00122    virtual void invert(G4int& ierr);
00123    // Invert a G4ErrorMatrix. G4ErrorMatrix must be square.
00124    // N.B. the contents of the matrix are replaced by the inverse.
00125    // Returns ierr = 0 (zero) when successful, otherwise non-zero. 
00126    // This method has less overhead then inverse().
00127 
00128    G4double determinant() const;
00129    // calculate the determinant of the matrix.
00130 
00131    G4double trace() const;
00132    // calculate the trace of the matrix (sum of diagonal elements).
00133 
00134    class G4ErrorMatrix_row
00135    {
00136      typedef std::vector<G4double >::const_iterator G4ErrorMatrixConstIter;
00137      public:
00138        inline G4ErrorMatrix_row(G4ErrorMatrix&,G4int);
00139        G4double & operator[](G4int);
00140      private:
00141        G4ErrorMatrix& _a;
00142        G4int _r;
00143    };
00144    class G4ErrorMatrix_row_const
00145    {
00146      public:
00147        inline G4ErrorMatrix_row_const (const G4ErrorMatrix&,G4int);
00148        const G4double & operator[](G4int) const;
00149      private:
00150        const G4ErrorMatrix& _a;
00151        G4int _r;
00152    };
00153    // helper classes for implementing m[i][j]
00154 
00155    inline G4ErrorMatrix_row operator[] (G4int);
00156    inline const G4ErrorMatrix_row_const operator[] (G4int) const;
00157    // Read or write a matrix element.
00158    // While it may not look like it, you simply do m[i][j] to get an element. 
00159    // ** Note that the indexing starts from [0][0]. **
00160 
00161  protected:
00162 
00163    virtual inline G4int num_size() const;
00164    virtual void invertHaywood4(G4int& ierr);
00165    virtual void invertHaywood5(G4int& ierr);
00166    virtual void invertHaywood6(G4int& ierr);
00167 
00168  public:
00169 
00170    static void error(const char *s);
00171 
00172  private:
00173 
00174    friend class G4ErrorMatrix_row;
00175    friend class G4ErrorMatrix_row_const;
00176    friend class G4ErrorSymMatrix;
00177    // Friend classes.
00178 
00179    friend G4ErrorMatrix operator+(const G4ErrorMatrix &m1,
00180                                   const G4ErrorMatrix &m2);
00181    friend G4ErrorMatrix operator-(const G4ErrorMatrix &m1,
00182                                   const G4ErrorMatrix &m2);
00183    friend G4ErrorMatrix operator*(const G4ErrorMatrix &m1,
00184                                   const G4ErrorMatrix &m2);
00185    friend G4ErrorMatrix operator*(const G4ErrorMatrix &m1,
00186                                   const G4ErrorSymMatrix &m2);
00187    friend G4ErrorMatrix operator*(const G4ErrorSymMatrix &m1,
00188                                   const G4ErrorMatrix &m2);
00189    friend G4ErrorMatrix operator*(const G4ErrorSymMatrix &m1,
00190                                   const G4ErrorSymMatrix &m2);
00191    // Multiply a G4ErrorMatrix by a G4ErrorMatrix or Vector.
00192 
00193    // solve the system of linear eq
00194 
00195    friend G4ErrorMatrix qr_solve(G4ErrorMatrix *, const G4ErrorMatrix &b);
00196    friend void tridiagonal(G4ErrorSymMatrix *a,G4ErrorMatrix *hsm);
00197    friend void row_house(G4ErrorMatrix *,const G4ErrorMatrix &, G4double,
00198                          G4int, G4int, G4int, G4int);
00199    friend void back_solve(const G4ErrorMatrix &R, G4ErrorMatrix *b);
00200    friend void col_givens(G4ErrorMatrix *A, G4double c,
00201                           G4double s, G4int k1, G4int k2, 
00202                           G4int rowmin, G4int rowmax);
00203    // Does a column Givens update.
00204 
00205    friend void row_givens(G4ErrorMatrix *A, G4double c,
00206                           G4double s, G4int k1, G4int k2, 
00207                           G4int colmin, G4int colmax);
00208    friend void col_house(G4ErrorMatrix *,const G4ErrorMatrix &, G4double,
00209                          G4int, G4int, G4int, G4int);
00210    friend void house_with_update(G4ErrorMatrix *a,G4int row,G4int col);
00211    friend void house_with_update(G4ErrorMatrix *a,G4ErrorMatrix *v,
00212                                  G4int row, G4int col);
00213    friend void house_with_update2(G4ErrorSymMatrix *a,G4ErrorMatrix *v,
00214                                   G4int row, G4int col); 
00215 
00216    G4int dfact_matrix(G4double &det, G4int *ir);
00217    // factorize the matrix. If successful, the return code is 0. On
00218    // return, det is the determinant and ir[] is row-interchange
00219    // matrix. See CERNLIB's DFACT routine.
00220 
00221    G4int dfinv_matrix(G4int *ir);
00222    // invert the matrix. See CERNLIB DFINV.
00223 
00224    std::vector<G4double > m;
00225 
00226    G4int nrow, ncol;
00227    G4int size;
00228 };
00229 
00230 // Operations other than member functions for G4ErrorMatrix
00231 
00232 G4ErrorMatrix operator*(const G4ErrorMatrix &m1, const G4ErrorMatrix &m2);
00233 G4ErrorMatrix operator*(G4double t, const G4ErrorMatrix &m1);
00234 G4ErrorMatrix operator*(const G4ErrorMatrix &m1, G4double t);
00235 // Multiplication operators
00236 // Note that m *= m1 is always faster than m = m * m1.
00237 
00238 G4ErrorMatrix operator/(const G4ErrorMatrix &m1, G4double t);
00239 // m = m1 / t. (m /= t is faster if you can use it.)
00240 
00241 G4ErrorMatrix operator+(const G4ErrorMatrix &m1, const G4ErrorMatrix &m2);
00242 // m = m1 + m2;
00243 // Note that m += m1 is always faster than m = m + m1.
00244 
00245 G4ErrorMatrix operator-(const G4ErrorMatrix &m1, const G4ErrorMatrix &m2);
00246 // m = m1 - m2;
00247 // Note that m -= m1 is always faster than m = m - m1.
00248 
00249 G4ErrorMatrix dsum(const G4ErrorMatrix&, const G4ErrorMatrix&);
00250 // Direct sum of two matrices. The direct sum of A and B is the matrix 
00251 //        A 0
00252 //        0 B
00253 
00254 std::ostream& operator<<(std::ostream &s, const G4ErrorMatrix &q);
00255 // Read in, write out G4ErrorMatrix into a stream.
00256  
00257 //
00258 // Specialized linear algebra functions
00259 //
00260 
00261 G4ErrorMatrix qr_solve(const G4ErrorMatrix &A, const G4ErrorMatrix &b);
00262 G4ErrorMatrix qr_solve(G4ErrorMatrix *A, const G4ErrorMatrix &b);
00263 // Works like backsolve, except matrix does not need to be upper
00264 // triangular. For nonsquare matrix, it solves in the least square sense.
00265 
00266 G4ErrorMatrix qr_inverse(const G4ErrorMatrix &A);
00267 G4ErrorMatrix qr_inverse(G4ErrorMatrix *A);
00268 // Finds the inverse of a matrix using QR decomposition.  Note, often what
00269 // you really want is solve or backsolve, they can be much quicker than
00270 // inverse in many calculations.
00271 
00272 
00273 void qr_decomp(G4ErrorMatrix *A, G4ErrorMatrix *hsm);
00274 G4ErrorMatrix qr_decomp(G4ErrorMatrix *A);
00275 // Does a QR decomposition of a matrix.
00276 
00277 void back_solve(const G4ErrorMatrix &R, G4ErrorMatrix *b);
00278 // Solves R*x = b where R is upper triangular.  Also has a variation that
00279 // solves a number of equations of this form in one step, where b is a matrix
00280 // with each column a different vector. See also solve.
00281 
00282 void col_house(G4ErrorMatrix *a, const G4ErrorMatrix &v, G4double vnormsq,
00283                G4int row, G4int col, G4int row_start, G4int col_start);
00284 void col_house(G4ErrorMatrix *a, const G4ErrorMatrix &v, G4int row, G4int col,
00285                G4int row_start, G4int col_start);
00286 // Does a column Householder update.
00287 
00288 void col_givens(G4ErrorMatrix *A, G4double c, G4double s,
00289                 G4int k1, G4int k2, G4int row_min=1, G4int row_max=0);
00290 // do a column Givens update
00291 
00292 void row_givens(G4ErrorMatrix *A, G4double c, G4double s,
00293                 G4int k1, G4int k2, G4int col_min=1, G4int col_max=0);
00294 // do a row Givens update
00295 
00296 void givens(G4double a, G4double b, G4double *c, G4double *s);
00297 // algorithm 5.1.5 in Golub and Van Loan
00298 
00299 // Returns a Householder vector to zero elements.
00300 
00301 void house_with_update(G4ErrorMatrix *a, G4int row=1, G4int col=1);
00302 void house_with_update(G4ErrorMatrix *a, G4ErrorMatrix *v,
00303                        G4int row=1, G4int col=1);
00304 // Finds and does Householder reflection on matrix.
00305 
00306 void row_house(G4ErrorMatrix *a, const G4ErrorMatrix &v, G4double vnormsq,
00307                G4int row, G4int col, G4int row_start, G4int col_start);
00308 void row_house(G4ErrorMatrix *a, const G4ErrorMatrix &v, G4int row, G4int col,
00309                G4int row_start, G4int col_start);
00310 // Does a row Householder update.
00311 
00312 #include "G4ErrorMatrix.icc"
00313 
00314 #endif

Generated on Mon May 27 17:48:11 2013 for Geant4 by  doxygen 1.4.7