LinearAlgebra.h

Go to the documentation of this file.
00001 #ifndef SimTK_LINEAR_ALGEBRA_H_
00002 #define SimTK_LINEAR_ALGEBRA_H_
00003 /* -------------------------------------------------------------------------- *
00004  *                      SimTK Core: SimTK Simmath(tm)                         *
00005  * -------------------------------------------------------------------------- *
00006  * This is part of the SimTK Core biosimulation toolkit originating from      *
00007  * Simbios, the NIH National Center for Physics-Based Simulation of           *
00008  * Biological Structures at Stanford, funded under the NIH Roadmap for        *
00009  * Medical Research, grant U54 GM072970. See https://simtk.org.               *
00010 
00011  * Portions copyright (c) 2006-2007 Stanford University and Jack Middleton.   *
00012  * Contributors:                                                              *
00013  *
00014  * Permission is hereby granted, free of charge, to any person obtaining      *
00015  * a copy of this software and associated documentation files (the            *
00016  * "Software"), to deal in the Software without restriction, including        *
00017  * without limitation the rights to use, copy, modify, merge, publish,        *
00018  * distribute, sublicense, and/or sell copies of the Software, and to         *
00019  * permit persons to whom the Software is furnished to do so, subject         *
00020  * to the following conditions:                                               *
00021  *
00022  * The above copyright notice and this permission notice shall be included    *
00023  * in all copies or substantial portions of the Software.                     *
00024  *
00025  * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS    *
00026  * OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF                 *
00027  * MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT.     *
00028  * IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY       *
00029  * CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT,       *
00030  * TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE          *
00031  * SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.                     *
00032  */
00033 
00034 
00043 #include <limits.h>
00044 #include "SimTKcommon.h"
00045 #include "SimTKmath.h"
00046 #include "SimTKcommon/internal/BigMatrix.h"
00047 #include "internal/common.h"
00048 
00049 namespace SimTK {
00050 
00051 //  default for reciprocal of the condition number
00052 // TODO: sherm 080128 I changed this from 0.01 to a more reasonable
00053 // value but it is still wrong because the default should depend
00054 // on the matrix size, something like max(m,n)*eps^(7/8) where
00055 // eps is machine precision for float or double as appropriate.
00056 static const double DefaultRecpCondition = 1e-12;
00060 class SimTK_SIMMATH_EXPORT Factor {
00061 public:
00062 
00063   Factor() {}
00065   template <class ELT> Factor( Matrix_<ELT> m );
00067   template <class ELT> void solve( const Vector_<ELT>& b, Vector_<ELT>& x ) const;
00069   template <class ELT> void solve( const Matrix_<ELT>& b, Matrix_<ELT>& x ) const;
00070   
00071 }; // class Factor
00072 
00073 class FactorLURepBase;
00074 
00078 class SimTK_SIMMATH_EXPORT FactorLU: public Factor {
00079     public:
00080 
00081     ~FactorLU();
00082 
00083     FactorLU();
00084     FactorLU( const FactorLU& c );
00085     FactorLU& operator=(const FactorLU& rhs);
00086 
00087     template <class ELT> FactorLU( const Matrix_<ELT>& m );
00089     template <class ELT> void factor( const Matrix_<ELT>& m );
00091     template <class ELT> void solve( const Vector_<ELT>& b, Vector_<ELT>& x ) const;
00093     template <class ELT> void solve( const Matrix_<ELT>& b, Matrix_<ELT>& x ) const;
00094 
00096     template <class ELT> void getL( Matrix_<ELT>& l ) const;
00098     template <class ELT> void getU( Matrix_<ELT>& u ) const;
00100     template < class ELT > void inverse(  Matrix_<ELT>& m ) const;
00101 
00103     bool isSingular() const;
00105     int getSingularIndex() const;
00106 
00107 
00108     protected:
00109     class FactorLURepBase *rep;
00110 
00111 }; // class FactorLU
00112 
00113 
00114 class FactorQTZRepBase;
00118 class SimTK_SIMMATH_EXPORT FactorQTZ: public Factor {
00119     public:
00120 
00121     ~FactorQTZ();
00122 
00123     FactorQTZ();
00124     FactorQTZ( const FactorQTZ& c );
00125     FactorQTZ& operator=(const FactorQTZ& rhs);
00127     template <typename ELT> FactorQTZ( const Matrix_<ELT>& m);
00129     template <typename ELT> FactorQTZ( const Matrix_<ELT>& m, double rcond );
00131     template <typename ELT> FactorQTZ( const Matrix_<ELT>& m, float rcond );
00133     template <typename ELT> void factor( const Matrix_<ELT>& m);
00135     template <typename ELT> void factor( const Matrix_<ELT>& m, float rcond );
00137     template <typename ELT> void factor( const Matrix_<ELT>& m, double rcond );
00139     template <typename ELT> void solve( const Vector_<ELT>& b, Vector_<ELT>& x ) const;
00141     template <typename ELT> void solve( const Matrix_<ELT>& b, Matrix_<ELT>& x ) const;
00142 
00143     template < class ELT > void inverse(  Matrix_<ELT>& m ) const;
00144 
00146     int getRank() const;
00147 //    void setRank(int rank); TBD
00149     void setReciprocalConditionNumber( Real rcond ); 
00150 
00151     protected:
00152     class FactorQTZRepBase *rep;
00153 
00154 }; // class FactorQTZ
00158 class SimTK_SIMMATH_EXPORT Eigen {
00159     public:
00160 
00161     ~Eigen();
00162 
00163     Eigen();
00164     Eigen( const Eigen& c );
00165     Eigen& operator=(const Eigen& rhs);
00166 
00168     template <class ELT> Eigen( const Matrix_<ELT>& m );
00170     template <class ELT> void factor( const Matrix_<ELT>& m );
00172     template <class VAL, class VEC> void getAllEigenValuesAndVectors( Vector_<VAL>& values, Matrix_<VEC>& vectors);
00174     template <class T> void getAllEigenValues( Vector_<T>& values);
00175 
00177     template <class VAL, class VEC> void getFewEigenValuesAndVectors( Vector_<VAL>& values, Matrix_<VEC>& vectors, int ilow, int ihi);
00179     template <class T> void getFewEigenVectors( Matrix_<T>& vectors, int ilow, int ihi );
00181     template <class T> void getFewEigenValues( Vector_<T>& values, int ilow, int ihi );
00182 
00184     template <class VAL, class VEC> void getFewEigenValuesAndVectors( Vector_<VAL>& values, Matrix_<VEC>& vectors, typename CNT<VAL>::TReal rlow, typename CNT<VAL>::TReal rhi);
00186     template <class T> void getFewEigenVectors( Matrix_<T>& vectors, typename CNT<T>::TReal rlow, typename CNT<T>::TReal rhi );
00188     template <class T> void getFewEigenValues( Vector_<T>& values, typename CNT<T>::TReal rlow, typename CNT<T>::TReal rhi );
00189 
00190      
00191     protected:
00192     class EigenRepBase *rep;
00193 
00194 }; // class Eigen
00198 class SimTK_SIMMATH_EXPORT FactorSVD: public Factor {
00199     public:
00200 
00201     ~FactorSVD();
00203     FactorSVD();
00205     FactorSVD( const FactorSVD& c );
00207     FactorSVD& operator=(const FactorSVD& rhs);
00208 
00210     template < class ELT > FactorSVD( const Matrix_<ELT>& m );
00213     template < class ELT > FactorSVD( const Matrix_<ELT>& m, float rcond );
00216     template < class ELT > FactorSVD( const Matrix_<ELT>& m, double rcond );
00218     template < class ELT > void factor( const Matrix_<ELT>& m );
00221     template < class ELT > void factor( const Matrix_<ELT>& m, float rcond );
00224     template < class ELT > void factor( const Matrix_<ELT>& m, double rcond );
00225 
00227     template < class T > void getSingularValuesAndVectors( Vector_<typename CNT<T>::TReal>& values, 
00228                               Matrix_<T>& leftVectors,  Matrix_<T>& rightVectors );
00230     template < class T > void getSingularValues( Vector_<T>& values);
00231 
00233     int getRank();
00235     template < class ELT > void inverse(  Matrix_<ELT>& m );
00237     template <class ELT> void solve( const Vector_<ELT>& b, Vector_<ELT>& x );
00240     template <class ELT> void solve( const Matrix_<ELT>& b, Matrix_<ELT>& x );
00241 
00242     protected:
00243     class FactorSVDRepBase *rep;
00244 
00245 }; // class FactorSVD
00246 
00247 } // namespace SimTK 
00248 
00249 #endif //SimTK_LINEAR_ALGEBRA_H_

Generated on Fri Sep 26 07:44:13 2008 for SimTKcore by  doxygen 1.5.6