Simbody  3.7
LinearAlgebra.h
Go to the documentation of this file.
1 #ifndef SimTK_LINEAR_ALGEBRA_H_
2 #define SimTK_LINEAR_ALGEBRA_H_
3 
4 /* -------------------------------------------------------------------------- *
5  * Simbody(tm): SimTKmath *
6  * -------------------------------------------------------------------------- *
7  * This is part of the SimTK biosimulation toolkit originating from *
8  * Simbios, the NIH National Center for Physics-Based Simulation of *
9  * Biological Structures at Stanford, funded under the NIH Roadmap for *
10  * Medical Research, grant U54 GM072970. See https://simtk.org/home/simbody. *
11  * *
12  * Portions copyright (c) 2006-12 Stanford University and the Authors. *
13  * Authors: Jack Middleton *
14  * Contributors: Michael Sherman *
15  * *
16  * Licensed under the Apache License, Version 2.0 (the "License"); you may *
17  * not use this file except in compliance with the License. You may obtain a *
18  * copy of the License at http://www.apache.org/licenses/LICENSE-2.0. *
19  * *
20  * Unless required by applicable law or agreed to in writing, software *
21  * distributed under the License is distributed on an "AS IS" BASIS, *
22  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. *
23  * See the License for the specific language governing permissions and *
24  * limitations under the License. *
25  * -------------------------------------------------------------------------- */
26 
33 #include "SimTKcommon.h"
35 
36 
37 namespace SimTK {
38 
39 // default for reciprocal of the condition number
40 // TODO: sherm 080128 I changed this from 0.01 to a more reasonable
41 // value but it is still wrong because the default should depend
42 // on the matrix size, something like max(m,n)*eps^(7/8) where
43 // eps is machine precision for float or double as appropriate.
44 static const double DefaultRecpCondition = 1e-12;
45 
50 public:
51 
52  Factor() {}
54  template <class ELT> Factor( Matrix_<ELT> m );
56  template <class ELT> void solve( const Vector_<ELT>& b, Vector_<ELT>& x ) const;
58  template <class ELT> void solve( const Matrix_<ELT>& b, Matrix_<ELT>& x ) const;
59 
60 }; // class Factor
61 
62 class FactorLURepBase;
63 
68  public:
69 
71 
73  FactorLU( const FactorLU& c );
74  FactorLU& operator=(const FactorLU& rhs);
75 
76  template <class ELT> FactorLU( const Matrix_<ELT>& m );
78  template <class ELT> void factor( const Matrix_<ELT>& m );
80  template <class ELT> void solve( const Vector_<ELT>& b, Vector_<ELT>& x ) const;
82  template <class ELT> void solve( const Matrix_<ELT>& b, Matrix_<ELT>& x ) const;
83 
85  template <class ELT> void getL( Matrix_<ELT>& l ) const;
87  template <class ELT> void getU( Matrix_<ELT>& u ) const;
89  template < class ELT > void inverse( Matrix_<ELT>& m ) const;
90 
92  bool isSingular() const;
94  int getSingularIndex() const;
95 
96 
97  protected:
98  class FactorLURepBase *rep;
99 
100 }; // class FactorLU
101 
102 
103 class FactorQTZRepBase;
108  public:
109 
111 
113  FactorQTZ( const FactorQTZ& c );
116  template <typename ELT> FactorQTZ( const Matrix_<ELT>& m);
118  template <typename ELT> FactorQTZ( const Matrix_<ELT>& m, double rcond );
120  template <typename ELT> FactorQTZ( const Matrix_<ELT>& m, float rcond );
122  template <typename ELT> void factor( const Matrix_<ELT>& m);
124  template <typename ELT> void factor( const Matrix_<ELT>& m, float rcond );
126  template <typename ELT> void factor( const Matrix_<ELT>& m, double rcond );
128  template <typename ELT> void solve( const Vector_<ELT>& b, Vector_<ELT>& x ) const;
130  template <typename ELT> void solve( const Matrix_<ELT>& b, Matrix_<ELT>& x ) const;
131 
132  template < class ELT > void inverse( Matrix_<ELT>& m ) const;
133 
135  int getRank() const;
137  double getRCondEstimate() const;
138 // void setRank(int rank); TBD
139 
140  protected:
141  class FactorQTZRepBase *rep;
142 }; // class FactorQTZ
147  public:
148 
150 
151  Eigen();
152  Eigen( const Eigen& c );
153  Eigen& operator=(const Eigen& rhs);
154 
156  template <class ELT> Eigen( const Matrix_<ELT>& m );
158  template <class ELT> void factor( const Matrix_<ELT>& m );
160  template <class VAL, class VEC> void getAllEigenValuesAndVectors( Vector_<VAL>& values, Matrix_<VEC>& vectors);
162  template <class T> void getAllEigenValues( Vector_<T>& values);
163 
165  template <class VAL, class VEC> void getFewEigenValuesAndVectors( Vector_<VAL>& values, Matrix_<VEC>& vectors, int ilow, int ihi);
167  template <class T> void getFewEigenVectors( Matrix_<T>& vectors, int ilow, int ihi );
169  template <class T> void getFewEigenValues( Vector_<T>& values, int ilow, int ihi );
170 
172  template <class VAL, class VEC> void getFewEigenValuesAndVectors( Vector_<VAL>& values, Matrix_<VEC>& vectors, typename CNT<VAL>::TReal rlow, typename CNT<VAL>::TReal rhi);
174  template <class T> void getFewEigenVectors( Matrix_<T>& vectors, typename CNT<T>::TReal rlow, typename CNT<T>::TReal rhi );
176  template <class T> void getFewEigenValues( Vector_<T>& values, typename CNT<T>::TReal rlow, typename CNT<T>::TReal rhi );
177 
178 
179  protected:
180  class EigenRepBase *rep;
181 
182 }; // class Eigen
187  public:
188 
193  FactorSVD( const FactorSVD& c );
196 
198  template < class ELT > FactorSVD( const Matrix_<ELT>& m );
201  template < class ELT > FactorSVD( const Matrix_<ELT>& m, float rcond );
204  template < class ELT > FactorSVD( const Matrix_<ELT>& m, double rcond );
206  template < class ELT > void factor( const Matrix_<ELT>& m );
209  template < class ELT > void factor( const Matrix_<ELT>& m, float rcond );
212  template < class ELT > void factor( const Matrix_<ELT>& m, double rcond );
213 
215  template < class T > void getSingularValuesAndVectors( Vector_<typename CNT<T>::TReal>& values,
216  Matrix_<T>& leftVectors, Matrix_<T>& rightVectors );
218  template < class T > void getSingularValues( Vector_<T>& values);
219 
221  int getRank();
223  template < class ELT > void inverse( Matrix_<ELT>& m );
225  template <class ELT> void solve( const Vector_<ELT>& b, Vector_<ELT>& x );
228  template <class ELT> void solve( const Matrix_<ELT>& b, Matrix_<ELT>& x );
229 
230  protected:
231  class FactorSVDRepBase *rep;
232 
233 }; // class FactorSVD
234 
235 } // namespace SimTK
236 
237 #endif //SimTK_LINEAR_ALGEBRA_H_
SimTK::FactorSVD::solve
void solve(const Matrix_< ELT > &b, Matrix_< ELT > &x)
solve for a set of x vectors given multiple right hand side vectors using the singular value decompos...
SimTK::FactorSVD::~FactorSVD
~FactorSVD()
SimTK::Factor::solve
void solve(const Vector_< ELT > &b, Vector_< ELT > &x) const
solves a single right hand side using a factorization
SimTK::FactorLU::FactorLU
FactorLU(const Matrix_< ELT > &m)
SimTK::Eigen::factor
void factor(const Matrix_< ELT > &m)
supply matrix which eigen values will be computed for
SimTK::Eigen::Eigen
Eigen(const Matrix_< ELT > &m)
create a default eigen class
SimTK::FactorLU::solve
void solve(const Vector_< ELT > &b, Vector_< ELT > &x) const
solves a single right hand side
SimTK::FactorQTZ::factor
void factor(const Matrix_< ELT > &m, double rcond)
do QTZ factorization of a matrix for a given reciprocal condition number
SimTK::Eigen::getFewEigenVectors
void getFewEigenVectors(Matrix_< T > &vectors, typename CNT< T >::TReal rlow, typename CNT< T >::TReal rhi)
get a few eigen vectors of a symmetric matrix which are within a range of eigen values
SimTK::FactorQTZ::getRCondEstimate
double getRCondEstimate() const
returns the actual reciprocal condition number at this rank
SimTK::Eigen::getFewEigenValuesAndVectors
void getFewEigenValuesAndVectors(Vector_< VAL > &values, Matrix_< VEC > &vectors, typename CNT< VAL >::TReal rlow, typename CNT< VAL >::TReal rhi)
get a few eigen values and eigen vectors of a symmetric matrix which are within a range of eigen valu...
SimTK::FactorQTZ::~FactorQTZ
~FactorQTZ()
SimTK::FactorQTZ::solve
void solve(const Matrix_< ELT > &b, Matrix_< ELT > &x) const
solve for an array of vectors given multiple right hand sides
SimTK::FactorSVD::factor
void factor(const Matrix_< ELT > &m, float rcond)
supply the matrix to do a singular value decomposition using the specified reciprocal of the conditio...
SimTK_SIMMATH_EXPORT
#define SimTK_SIMMATH_EXPORT
Definition: SimTKmath/include/simmath/internal/common.h:64
SimTK
This is a System that represents the dynamics of a particle moving along a smooth surface.
Definition: Assembler.h:37
SimTK::FactorSVD::solve
void solve(const Vector_< ELT > &b, Vector_< ELT > &x)
solve for x given a right hand side vector using the singular value decomposition
SimTK::FactorSVD::getSingularValuesAndVectors
void getSingularValuesAndVectors(Vector_< typename CNT< T >::TReal > &values, Matrix_< T > &leftVectors, Matrix_< T > &rightVectors)
get the singular values and singular vectors of the matrix
SimTK::Eigen::getFewEigenVectors
void getFewEigenVectors(Matrix_< T > &vectors, int ilow, int ihi)
get a few eigen vectors of a symmetric matrix which are within a range of indices
SimTK::FactorLU::~FactorLU
~FactorLU()
SimTK::FactorSVD::factor
void factor(const Matrix_< ELT > &m)
supply the matrix to do a singular value decomposition
SimTK::FactorQTZ::FactorQTZ
FactorQTZ(const FactorQTZ &c)
SimTK::FactorLU::FactorLU
FactorLU(const FactorLU &c)
SimTK::Factor::Factor
Factor(Matrix_< ELT > m)
creates an factorization of a matrix
SimTK::FactorSVD::factor
void factor(const Matrix_< ELT > &m, double rcond)
supply the matrix to do a singular value decomposition using the specified reciprocal of the conditio...
SimTK::FactorSVD::FactorSVD
FactorSVD(const Matrix_< ELT > &m, double rcond)
singular value decomposition of a matrix using the specified reciprocal of the condition number rcond
SimTK::Eigen
Class to compute Eigen values and Eigen vectors of a matrix.
Definition: LinearAlgebra.h:146
SimTKcommon.h
common.h
SimTK::FactorQTZ::factor
void factor(const Matrix_< ELT > &m, float rcond)
do QTZ factorization of a matrix for a given reciprocal condition number
SimTK::Factor
Base class for the various matrix factorizations.
Definition: LinearAlgebra.h:49
SimTK::Eigen::Eigen
Eigen()
SimTK::FactorLU::operator=
FactorLU & operator=(const FactorLU &rhs)
SimTK::Eigen::~Eigen
~Eigen()
SimTK::FactorSVD::inverse
void inverse(Matrix_< ELT > &m)
get inverse of the matrix using singular value decomposition (sometimes called the pseudo inverse)
SimTK::FactorLU
Class for performing LU matrix factorizations.
Definition: LinearAlgebra.h:67
SimTK::FactorSVD::FactorSVD
FactorSVD()
default constructor
SimTK::FactorQTZ::FactorQTZ
FactorQTZ(const Matrix_< ELT > &m)
do QTZ factorization of a matrix
SimTK::FactorQTZ::factor
void factor(const Matrix_< ELT > &m)
do QTZ factorization of a matrix
SimTK::Eigen::getAllEigenValuesAndVectors
void getAllEigenValuesAndVectors(Vector_< VAL > &values, Matrix_< VEC > &vectors)
get all the eigen values and eigen vectors of a matrix
SimTK::Eigen::getAllEigenValues
void getAllEigenValues(Vector_< T > &values)
get all the eigen values of a matrix
SimTK::FactorQTZ::operator=
FactorQTZ & operator=(const FactorQTZ &rhs)
SimTK::FactorLU::getSingularIndex
int getSingularIndex() const
returns the first diagonal which was found to be singular
SimTK::FactorSVD::FactorSVD
FactorSVD(const Matrix_< ELT > &m)
constructor
SimTK::FactorLU::getU
void getU(Matrix_< ELT > &u) const
returns the upper triangle of an LU factorization
SimTK::FactorQTZ::FactorQTZ
FactorQTZ()
SimTK::FactorQTZ::FactorQTZ
FactorQTZ(const Matrix_< ELT > &m, double rcond)
do QTZ factorization of a matrix for a given reciprocal condition number
SimTK::FactorLU::getL
void getL(Matrix_< ELT > &l) const
returns the lower triangle of an LU factorization
SimTK::FactorLU::isSingular
bool isSingular() const
returns true if matrix was singular
SimTK::FactorLU::solve
void solve(const Matrix_< ELT > &b, Matrix_< ELT > &x) const
solves multiple right hand sides
SimTK::FactorQTZ::rep
class FactorQTZRepBase * rep
Definition: LinearAlgebra.h:141
SimTK::FactorQTZ
Class to perform a QTZ (linear least squares) factorization.
Definition: LinearAlgebra.h:107
SimTK::FactorSVD::getSingularValues
void getSingularValues(Vector_< T > &values)
get just the singular values of the matrix
SimTK::FactorSVD::FactorSVD
FactorSVD(const FactorSVD &c)
copy constructor
SimTK::Factor::Factor
Factor()
Definition: LinearAlgebra.h:52
SimTK::DefaultRecpCondition
static const double DefaultRecpCondition
Definition: LinearAlgebra.h:44
SimTK::FactorQTZ::getRank
int getRank() const
returns the rank of the matrix
SimTK::FactorSVD::FactorSVD
FactorSVD(const Matrix_< ELT > &m, float rcond)
singular value decomposition of a matrix using the specified reciprocal of the condition number rcond
SimTK::Vector_
This is the vector class intended to appear in user code for large, variable size column vectors.
Definition: BigMatrix.h:171
SimTK::Eigen::getFewEigenValuesAndVectors
void getFewEigenValuesAndVectors(Vector_< VAL > &values, Matrix_< VEC > &vectors, int ilow, int ihi)
get a few eigen values and eigen vectors of a symmetric matrix which are within a range of indices
SimTK::Factor::solve
void solve(const Matrix_< ELT > &b, Matrix_< ELT > &x) const
solves multiple right hand sides using a factorization
SimTK::Eigen::getFewEigenValues
void getFewEigenValues(Vector_< T > &values, typename CNT< T >::TReal rlow, typename CNT< T >::TReal rhi)
get a few eigen values of a symmetric matrix which are within a range of eigen values
SimTK::Eigen::getFewEigenValues
void getFewEigenValues(Vector_< T > &values, int ilow, int ihi)
get a few eigen values of a symmetric matrix which are within a range of indices
SimTK::FactorLU::factor
void factor(const Matrix_< ELT > &m)
factors a matrix
SimTK::FactorQTZ::FactorQTZ
FactorQTZ(const Matrix_< ELT > &m, float rcond)
do QTZ factorization of a matrix for a given reciprocal condition number
SimTK::FactorLU::FactorLU
FactorLU()
SimTK::FactorSVD::getRank
int getRank()
get rank of the matrix
SimTK::CNT::TReal
K::TReal TReal
Definition: CompositeNumericalTypes.h:141
SimTK::FactorQTZ::inverse
void inverse(Matrix_< ELT > &m) const
SimTK::Eigen::operator=
Eigen & operator=(const Eigen &rhs)
SimTK::FactorLU::inverse
void inverse(Matrix_< ELT > &m) const
returns the inverse of a matrix using an LU factorization
SimTK::Eigen::Eigen
Eigen(const Eigen &c)
SimTK::Eigen::rep
class EigenRepBase * rep
Definition: LinearAlgebra.h:180
SimTK::FactorSVD::operator=
FactorSVD & operator=(const FactorSVD &rhs)
copy assign
SimTK::Matrix_
This is the matrix class intended to appear in user code for large, variable size matrices.
Definition: BigMatrix.h:168
SimTK::FactorSVD
Class to compute a singular value decomposition of a matrix.
Definition: LinearAlgebra.h:186
SimTK::FactorQTZ::solve
void solve(const Vector_< ELT > &b, Vector_< ELT > &x) const
solve for a vector x given a right hand side vector b
SimTK::FactorSVD::rep
class FactorSVDRepBase * rep
Definition: LinearAlgebra.h:231
SimTK::FactorLU::rep
class FactorLURepBase * rep
Definition: LinearAlgebra.h:98