MatrixX< T > Class Template Reference
[Dynamic-size classesMatrices]

Dynamic-size matrix. More...

#include <matrix.h>

Inheritance diagram for MatrixX< T >:

MatrixBase< T, Derived, VectorType, LUDecompositionType >

List of all members.

Public Types

typedef T ScalType
typedef VectorType VecType

Public Member Functions

Derived adjoint () const
T * array ()
const T * array () const
VectorType column (int col) const
T * columnPtr (int col)
const T * columnPtr (int col) const
void computeAdjoint (Derived *result) const
bool computeBasisKer (Derived *result) const
void computeInverse (Derived *result) const
void computeInverseSafely (Derived *result, bool *invertible) const
bool computeSomeAntecedent (const VectorType &v, VectorType *result) const
determinant () const
int dimKer () const
void findBiggestEntry (int *res_row, int *res_col, int skip=0) const
void getColumn (int column, VectorType *ret) const
void getColumn (int column, T *ret) const
void getRow (int row, VectorType *ret) const
void getRow (int row, T *ret) const
Derived inverse () const
bool isApprox (const Derived &other, const T &precision=Util::epsilon< T >()) const
bool isInvertible () const
bool isNegligible (const T &other, const T &precision=Util::epsilon< T >()) const
bool isZero (const T &precision=Util::epsilon< T >()) const
void leftmultiply (const VectorType &vector, VectorType *res) const
Derived & loadDiagonal (const VectorType &coeffs)
Derived & loadDiagonal (const T &coeff)
Derived & loadIdentity ()
Derived & loadOrthoBasis (const VectorType &vector)
Derived & loadRandom ()
Derived & loadRotation2 (const T &angle)
Derived & loadRotation3 (const T &angle, const VectorType &axis)
Derived & loadScaling (const VectorType &coeffs)
Derived & loadScaling (const T &coeff)
Derived & loadZero ()
 MatrixX (int size, const T *array)
 MatrixX (const MatrixX &other)
 MatrixX (int size=1)
void multiply (const VectorType &vector, VectorType *res) const
void multiply (const Derived &other, Derived *res) const
void multiplyEntries (const Derived &other, Derived *res) const
bool operator!= (const Derived &other) const
T & operator() (int row, int col)
const T & operator() (int row, int col) const
VectorType operator* (const VectorType &vector) const
Derived operator* (const T &factor) const
Derived operator* (const Derived &other) const
Derived & operator*= (const T &factor)
Derived & operator*= (const Derived &other)
Derived operator+ (const Derived &other) const
Derived & operator+= (const Derived &other)
Derived operator- () const
Derived operator- (const Derived &other) const
Derived & operator-= (const Derived &other)
Derived operator/ (const T &factor) const
Derived & operator/= (const T &factor)
bool operator== (const Derived &other) const
T & operator[] (int i)
const T & operator[] (int i) const
Derived & prerotate2 (const T &angle)
Derived & prerotate3 (const T &angle, const VectorType &axis)
Derived & prescale (const VectorType &coeffs)
int rank () const
void readArray (const T *src)
void readRows (const T *rows)
void recursiveGramSchmidt (int n)
void replaceWithAdjoint ()
void replaceWithOpposite ()
void resize (int newsize)
Derived & rotate2 (const T &angle)
Derived & rotate3 (const T &angle, const VectorType &axis)
VectorType row (int row) const
Derived & scale (const VectorType &coeffs)
void setColumn (int column, const VectorType &src)
void setColumn (int column, const T *src)
void setRow (int row, const VectorType &src)
void setRow (int row, const T *src)
int size () const
 ~MatrixX ()

Static Public Member Functions

static bool hasDynamicSize ()


Detailed Description

template<typename T>
class Eigen::MatrixX< T >

Dynamic-size matrix.

A class for dynamic-size square matrices.

The template parameter, T, is the type of the entries of the matrix. It can be any type representing either real or complex numbers. The following typedefs are provided to cover the usual cases:

    typedef MatrixX<double>                 MatrixXd;
    typedef MatrixX<float>                  MatrixXf;
    typedef MatrixX< std::complex<double> > MatrixXcd;
    typedef MatrixX< std::complex<float> >  MatrixXcf;

If you prefer fixed-size matrices (they are faster), see the Matrix class template, which provides exactly the same functionality and API in fixed-size version.

The MatrixX class template provides all the usual operators and methods to manipulate square matrices, and can do some more advanced stuff using involving a LU decomposition, like invert a matrix.

Here are some examples of usage of MatrixX:

    using namespace Eigen;
    using namespace std; // we'll use cout for outputting matrices

    MatrixXd mat1(3), mat2(3); // constructs two new uninitialized 3x3 matrices
    double array[] = { 2.4, 3.1, -0.7,
                       -1.1, 2.9, 4.3,
                       6.7, -3.7, 2.7 };
    mat1.readArray(array);
        // now mat1 is the following matrix:
        // (  2.4  -1.1   6.7 )
        // (  3.1   2.9  -3.7 )
        // ( -0.7   4.3   2.7 )
    mat2.readRows(array);
        // now mat2 is the following matrix:
        // (  2.4   3.1  -0.7 )
        // ( -1.1   2.9   4.3 )
        // (  6.7  -3.7   2.7 )

    mat1 *= mat2; // computes the matrix product mat1 * mat2, stores it in mat1
    mat1 = mat2 + mat1 * mat2; // there are also non-assignment operators
    mat1 = 0.9 * mat1 + mat2 / 2.6; // you can also multiply/divide by numbers

    MatrixXd mat3(5); // construct a new uninitialized 5x5 matrix
    mat3.loadIdentity(); // loads the identity matrix of size 5 into mat3

    mat1(2,3) = -1.4; // Stores the value -1.4 at row 2, column 3 of mat1.
    cout << mat1 << endl;

    double vec_coords[] = {1.1, -2.5, 0.8};
    VectorXd vec( 3, vec_coords ); // construct a vector of size 3.
    vec = mat1 * vec; // multiply vec by mat1

    cout << "determinant of mat1: " << mat1.determinant() << endl;
        // as mat1 has size 3, the determinant is computed by brute force.
        // For larger matrices, this would use a LU decomposition.

Member Typedef Documentation

typedef T ScalType [inherited]

typedef VectorType VecType [inherited]


Constructor & Destructor Documentation

MatrixX ( int  size = 1  )  [inline, explicit]

Constructs an uninitialized square matrix with given size (dimension). The default value for size is 1.

MatrixX ( const MatrixX< T > &  other  )  [inline]

Copy constructor.

MatrixX ( int  size,
const T *  array 
) [inline]

Constructs a matrix with given size from an array. The matrix entries must be stored in column-dominant order in the array.

Parameters:
array the array. It must have dimension at least size*size.
size the size of the matrix (number of rows, or equivalently number of columns).

~MatrixX (  )  [inline]

Destructor, frees the memory allocated for the matrix's array


Member Function Documentation

Derived adjoint (  )  const [inline, inherited]

Returns the adjoint (conjugate transpose, equals transpose unless T is complex numbers) of *this.

This method returns an object by value, which is inefficient. For better performance, use computeAdjoint() or replaceWithAdjoint()

See also:
computeAdjoint(), replaceWithAdjoint()

T* array (  )  [inline, inherited]

Returns:
the array of the matrix, as non-constant.
See also:
operator()(int,int), operator[](int)

const T* array (  )  const [inline, inherited]

Returns:
the array of the matrix, as constant.
See also:
operator()(int,int) const, operator[](int) const

VectorType column ( int  col  )  const [inline, inherited]

Returns a column of the matrix.

This method returns an object by value, which is inefficient. For better performance, use another method.

See also:
columnPtr(int) const, getColumn(int,T*) const, getColumn(int, VectorType *) const

T* columnPtr ( int  col  )  [inline, inherited]

Returns the array data of a column of the matrix, as non-constant

See also:
columnPtr(int) const, column(), setColumn(int, const T*), setColumn(int, const VectorType &)

const T* columnPtr ( int  col  )  const [inline, inherited]

Returns:
the array data of a column of the matrix, as constant
See also:
columnPtr(int), column(), getColumn(int,T*) const, getColumn(int, VectorType *) const

void computeAdjoint ( Derived *  result  )  const [inline, inherited]

Computes the adjoint (conjugate transpose, equals transpose unless T is complex numbers) of *this and stores it in *result.

Don't try to pass this as result , this won't work. To replace *this by its adjoint, use replaceWithAdjoint() instead.

See also:
adjoint(), replaceWithAdjoint()

bool computeBasisKer ( Derived *  result  )  const [inline, inherited]

This method computes a basis of the kernel of the matrix *this.

The result argument is a pointer to the matrix in which the result will be stored. After calling this method, the dimKer() first columns of result form a basis of the kernel of *this. If *this is invertible, then dimKer()==0 and this method does nothing.

Returns:
true if the matrix is non-invertible, hence has a nonzero kernel; false if the matrix is invertible.
This function always performs a LU decomposition, even for small matrices. If you perform other operations that use a LU decomposition, it's better for performance to compute the LU decomposition once and then do the computations from it. See class LUDecomposition and class LUDecompositionX.

See also:
dimKer(), LUDecompositionBase::computeBasisKer()

void computeInverse ( Derived *  result  )  const [inline, inherited]

Computes the inverse matrix of *this, and stores it in *result.

If *this is non-invertible, the *result matrix is left with an undefined value. Therefore, only call this methods on matrices that you know are invertible. You can check by calling isInvertible().

If the matrix has fixed size <= 3, the inverse is obtained by direct computation. If the matrix has size >= 4 and/or has dynamic size, a LU decomposition is computed, and the inverse is obtained from it. If you perform other operations that use a LU decomposition, it's better for performance to compute the LU decomposition once and then do the computations from it. See class LUDecomposition and class LUDecompositionX.

See also:
inverse(), computeInverseSafely()

void computeInverseSafely ( Derived *  result,
bool *  invertible 
) const [inline, inherited]

Safely computes the inverse matrix of *this.

If *this is invertible, its inverse is computed and stored in *result. Moreover, *invertible is set to true.

If *this is non-invertible, the *result matrix is left unchanged and *invertible is set to false.

If the matrix has size <= 3, the inverse is obtained by direct computation. If the matrix has size >= 4, a LU decomposition is computed, and the inverse is obtained from it. If you perform other operations that use a LU decomposition, it's better for performance to compute the LU decomposition once and then do the computations from it. See class LUDecomposition and class LUDecompositionX.

See also:
inverse(), computeInverse()

bool computeSomeAntecedent ( const VectorType &  v,
VectorType *  result 
) const [inline, inherited]

Computes an antecedent of v by *this, that is, a vector u such that Au=v, where A is the matrix *this. If such an antecedent doesn't exist, this method does nothing.

Returns:
true if an antecedent exists, false if no antecedent exists.
Notes:

1. The returned vector u is only one solution of the equation Au=v, which can have more than one solution if A is non-invertible. To get a basis of the whole (affine) space of solutions, use computeBasisKer().

2. This function always performs a LU decomposition, even for small matrices. If you perform other operations that use a LU decomposition, it's better for performance to compute the LU decomposition once and then do the computations from it. See class LUDecomposition and class LUDecompositionX.

See also:
computeBasisKer(), LUDecompositionBase::computeSomeAntecedent()

T determinant (  )  const [inline, inherited]

Returns:
the determinant of the matrix (must be a square matrix).
If the matrix has size <= 3, the determinant is obtained by direct computation. If the matrix has size >= 4, a LU decomposition is computed, and the determinant is obtained from it.

If you perform other operations that use a LU decomposition, it's better for performance to compute the LU decomposition once and then do the computations from it. See class LUDecomposition and class LUDecompositionX.

See also:
LUDecompositionBase::determinant()

int dimKer (  )  const [inline, inherited]

Returns:
the dimension of the kernel of the matrix.
This function always performs a LU decomposition, even for small matrices. If you perform other operations that use a LU decomposition, it's better for performance to compute the LU decomposition once and then do the computations from it. See class LUDecomposition and class LUDecompositionX.

See also:
rank(), isInvertible(), LUDecompositionBase::dimKer()

void findBiggestEntry ( int *  res_row,
int *  res_col,
int  skip = 0 
) const [inline, inherited]

Stores in *res_row and *res_col the position of the entry of the matrix that has biggest absolute value. Skips the skip first rows and columns (default value for skip is 0).

This method is not very fast currently, because it has a nested for loop that doesn't get unrolled and for which we don't provide an unrolled version.

void getColumn ( int  column,
VectorType *  ret 
) const [inline, inherited]

Copies a column of the matrix into the vector ret

See also:
getColumn(int, T *) const, column(), columnPtr(int) const

void getColumn ( int  column,
T *  ret 
) const [inline, inherited]

Copies a column of the matrix into the array ret

See also:
getColumn(int, VectorType *) const, column(), columnPtr(int) const

void getRow ( int  row,
VectorType *  ret 
) const [inline, inherited]

Copies a row of the matrix into the vector ret

See also:
getRow(int, T *) const, row()

void getRow ( int  row,
T *  ret 
) const [inline, inherited]

Copies a row of the matrix into the array ret

See also:
getRow(int, VectorType *) const, row()

static bool hasDynamicSize (  )  [inline, static, inherited]

Returns:
true if the matrix has dynamic size (i.e. is an object of class MatrixX), false if the matrix has fixed size (i.e. is an object of class Matrix).
See also:
size(), resize()

Derived inverse (  )  const [inline, inherited]

This methods returns the inverse matrix of *this. If *this is non-invertible, the returned value is undefined.

This method calls computeInverse(), so the same remarks as for computeInverse() apply here. If you perform other operations that use a LU decomposition, it's better for performance to compute the LU decomposition once and then do the computations from it. See class LUDecomposition and class LUDecompositionX.

This method returns an object by value, which is inefficient.

See also:
computeInverse(), computeInverseSafely()

bool isApprox ( const Derived &  other,
const T &  precision = Util::epsilon<T>() 
) const [inline, inherited]

Returns true if *this and other are approximately equal.

The optional parameter precision allows to control the number of significant digits of precision. For instance, setting precision to 1e-5 results in a precision of 5 decimal digits.

This test is for nonzero matrices. If either of the two matrices being compared is zero, then it returns true if, and only if the other one is also zero -- which is not what one typically wants.

To compare a matrix with the zero matrix, i.e. to check whether a matrix is approximately zero, use isZero() instead.

See also:
operator==(), operator!=(), isZero()

bool isInvertible (  )  const [inline, inherited]

Returns:
true if the matrix is invertible, false if it is singular.
If the matrix haS size <= 3. If the matrix has size >= 4, a LU decomposition is computed, and the invertibility is obtained from it.

If you perform other operations that use a LU decomposition, it's better for performance to compute the LU decomposition once and then do the computations from it. See class LUDecomposition and class LUDecompositionX.

See also:
rank(), dimKer(), LUDecompositionBase::isInvertible(), computeInverseSafely()

bool isNegligible ( const T &  other,
const T &  precision = Util::epsilon<T>() 
) const [inline, inherited]

Returns true if all entries of *this are smaller (in absolute value) than other*precision. In other words, returns true if all entries are much smaller than other. For the meaning of precision, see isApprox().

See also:
isApprox(), isZero()

bool isZero ( const T &  precision = Util::epsilon<T>()  )  const [inline, inherited]

Tests whether *this is approximately equal to the zero matrix.

Equivalent to isNegligible(1). In other words, returns true if all entries of *this are approximately zero, in the sense that they have absolute value smaller than epsilon.

See also:
isNegligible(), isApprox()

void leftmultiply ( const VectorType &  vector,
VectorType *  res 
) const [inline, inherited]

Computes the product vector * *this and stores the result into *res. Here, vector is regarded as a row vector.

For dynamic-size classes, this method resizes *res if necessary. For fixed-size classes, it is required that *res already has the right size, that is: res->size() == this->size().

This method is faster than operator*.

Note for dynamic-size classes: For optimal performance, the vector *res should be initialized with the correct size before calling this method, otherwise it will have to be resized, which is costly.

Note that for the product to make sense, it is required that this->size() == vector.size().

See also:
multiply( const VectorType &, VectorType *) const

Derived & loadDiagonal ( const VectorType &  coeffs  )  [inline, inherited]

Sets *this to be the diagonal matrix with diagonal entries given by the coords of the vector coeffs. If *this doesn't have the same size as coeffs, it gets resized (if it has dynamic size).

See also:
loadScaling(const VectorType &), loadDiagonal(const T &)

Derived & loadDiagonal ( const T &  coeff  )  [inline, inherited]

Sets *this to be coeff times the identity matrix. In other words, the diagonal entries are set to coeff and the non-diagonal entries are set to zero.

See also:
loadScaling(const T&), loadDiagonal(const VectorType &)

Derived& loadIdentity (  )  [inline, inherited]

Sets *this to be the identity matrix.

See also:
loadZero(), loadRandom()

Derived & loadOrthoBasis ( const VectorType &  vector  )  [inline, inherited]

This method loads vector into the first column of *this, and fills the other columns with unit vectors such that all the columns of *this together form an orthogonal basis. This is an orthonormal basis if, and only if vector is a unit vector.

If *this has dynamic size, it gets resized to match the size of vector. If it has fixed size, then it is required that its size equals the size of vector.

Note:
this method assumes that vector is nonzero.
See also:
recursiveGramSchmidt(), VectorBase::makeOrthoVector()

Derived& loadRandom (  )  [inline, inherited]

Sets all entries to random values between -1.0 and 1.0. For complex numbers, both the real and imaginary parts can range from -1.0 to 1.0.

See also:
loadIdentity(), loadZero()

Derived & loadRotation2 ( const T &  angle  )  [inline, inherited]

Sets *this to be a 2-dimensional rotation of given angle.

Parameters:
angle the angle expressed in radians, counter-clockwise
See also:
rotate2(), prerotate2(), loadRotation3()

Derived & loadRotation3 ( const T &  angle,
const VectorType &  axis 
) [inline, inherited]

Sets *this to be a 3-dimensional rotation of given angle in radians.

Parameters:
angle the angle expressed in radians, counter-clockwise if the axis vector it oriented towards the observer.
axis the axis vector around which to rotate. Must be a unit vector, i.e. it is required that axis.norm() == 1.
Note:
if axis.norm() is different from 1, the result is undefined, and certainly won't be a rotation matrix.
See also:
rotate3(), prerotate3(), loadRotation2()

Derived& loadScaling ( const VectorType &  coeffs  )  [inline, inherited]

Derived& loadScaling ( const T &  coeff  )  [inline, inherited]

Derived & loadZero (  )  [inline, inherited]

Sets *this to be the zero matrix.

See also:
loadIdentity(), loadRandom()

void multiply ( const VectorType &  vector,
VectorType *  res 
) const [inline, inherited]

Computes the product *this * vector and stores the result into *res.

For dynamic-size classes, this method resizes *res if necessary. For fixed-size classes, it is required that *res already has the right size, that is: res->size() == this->size().

This method is faster than operator*.

Note for dynamic-size classes: For optimal performance, the vector *res should be initialized with the correct size before calling this method, otherwise it will have to be resized, which is costly.

Note that for the product to make sense, it is required that this->size() == vector.size().

See also:
operator*(const VectorType &) const, leftmultiply()

void multiply ( const Derived &  other,
Derived *  res 
) const [inline, inherited]

Computes the matrix product *this * other and stores the result into *res.

For dynamic-size matrices, this method resizes *res if necessary. For fixed-size matrices, it is required that *res already has the right size, that is: res->size() == this->size().

This method is faster than operator* and operator*=.

Note for dynamic-size matrices: For optimal performance, the matrix res should be initialized with the correct size before calling this method, otherwise it will have to be resized, which is costly.

Note that for the product to make sense, it is required that this->size() == other.size().

See also:
operator*(const Derived &) const, operator*=(const Derived &)

void multiplyEntries ( const Derived &  other,
Derived *  res 
) const [inline, inherited]

Entry-wise multiplication of matrices. This is NOT the standard matrix-matrix product.

*res becomes the matrix such that (*res)(i,j) = (*this)(i,j) * other(i,j).

See also:
multiply(const Derived &, Derived *) const

bool operator!= ( const Derived &  other  )  const [inline, inherited]

Equivalent to !isApprox() with the default precision.

Note:
Despite the name, this operator does a fuzzy compare! It is not equivalent to operator!= on each entry.
See also:
isApprox(),operator==(),isZero()

T& operator() ( int  row,
int  col 
) [inline, inherited]

Returns:
a non-constant reference to the entry of the matrix at given row and column.
See also:
array(int), operator[](int,int)

const T& operator() ( int  row,
int  col 
) const [inline, inherited]

Returns:
a constant reference to the entry of the matrix at given row and column.
See also:
array() const, operator[](int) const

VectorType operator* ( const VectorType &  vector  )  const [inline, inherited]

Returns *this * vector (multiplication of vector by matrix). The size of *this must equal the size of vector.

This method returns an object by value, which is inefficient. For better performance, use multiply( const VectorType &, VectorType *) const

See also:
multiply( const VectorType &, VectorType *) const

Derived operator* ( const T &  factor  )  const [inline, inherited]

Returns *this * factor (multiplication of each coord).

This method returns an object by value, which is inefficient.

See also:
operator*=(const T &)

Derived operator* ( const Derived &  other  )  const [inline, inherited]

Returns the matrix product *this * other.

Note that for the product to make sense, it is required that this->size() == other.size().

This method returns an object by value, which is inefficient. For better performance, use multiply(const Derived &, Derived *) const.

See also:
multiply(const Derived &, Derived *) const, operator*=(const Derived &)

Derived& operator*= ( const T &  factor  )  [inline, inherited]

Stores *this * factor into *this (multiplication of each entry).

See also:
operator*(const T&) const

Derived& operator*= ( const Derived &  other  )  [inline, inherited]

Matrix multiplication on the right: does *this = (*this) * other.

This only makes sense if this->size() == other.size().

See also:
multiply(const Derived &, Derived *) const, operator*(const Derived &) const

Derived operator+ ( const Derived &  other  )  const [inline, inherited]

Returns *this + other (entry-wise addition). The matrices *this and other must have the same size.

This method returns an object by value, which is inefficient.

See also:
operator+=()

Derived& operator+= ( const Derived &  other  )  [inline, inherited]

Stores *this + other into *this (entry-wise addition).

*this and other must have the same size.

See also:
operator+()

Derived operator- ( void   )  const [inline, inherited]

Returns (-(*this)).

This method returns an object by value, which is inefficient.

See also:
replaceWithOpposite(), operator-(const Derived &) const

Derived operator- ( const Derived &  other  )  const [inline, inherited]

Returns *this - other (entry-wise substraction). The matrices *this and other must have the same size.

This method returns an object by value, which is inefficient.

See also:
operator-=(), operator-(void) const

Derived& operator-= ( const Derived &  other  )  [inline, inherited]

Stores *this - other into *this (entry-wise substraction).

*this and other must have the same size.

See also:
operator-(const Derived &) const

Derived operator/ ( const T &  factor  )  const [inline, inherited]

Returns *this / factor (division of each coord).

This method returns an object by value, which is inefficient.

See also:
operator*(const T&) const, operator/=(const T&)

Derived& operator/= ( const T &  factor  )  [inline, inherited]

Stores *this / factor into *this (division of each entry).

See also:
operator*=(const T&)

bool operator== ( const Derived &  other  )  const [inline, inherited]

Equivalent to isApprox() with the default precision.

Note:
Despite the name, this operator does a fuzzy compare! It is not equivalent to operator== on each entry.
See also:
isApprox(),operator!=(),isZero()

T& operator[] ( int  i  )  [inline, inherited]

Returns:
a non-constant reference to the i-th entry of the array of the matrix (which stores the matrix entries in column-major order).
See also:
array(), operator()(int,int)

const T& operator[] ( int  i  )  const [inline, inherited]

Returns:
a constant reference to the i-th entry of the array of the matrix (which stores the matrix entries in column-major order).
See also:
array() const, operator()(int,int) const

Derived& prerotate2 ( const T &  angle  )  [inline, inherited]

Multiplies *this on the left by the rotation matrix of given angle in radians. The template parameter Size must equal 2.

See also:
loadRotation2(), rotate2()

Derived& prerotate3 ( const T &  angle,
const VectorType &  axis 
) [inline, inherited]

Multiplies *this on the left by the rotation matrix of given angle in radians around given axis vector. Only applicable to matrices of size 3.

See also:
rotate3(), loadRotation3()

Derived& prescale ( const VectorType &  coeffs  )  [inline, inherited]

Multiplies *this on the left by the scaling matrix with given vector of scaling coefficients.

See also:
loadScaling(const VectorType &), scale(const VectorType &)

int rank (  )  const [inline, inherited]

Returns:
the rank of the matrix.
This function always performs a LU decomposition, even for small matrices. If you perform other operations that use a LU decomposition, it's better for performance to compute the LU decomposition once and then do the computations from it. See class LUDecomposition and class LUDecompositionX.

See also:
dimKer(), isInvertible(), LUDecompositionBase::rank()

void readArray ( const T *  src  )  [inline, inherited]

Reads the entries of *this from an array. The number of entries read from the array is equal to size()*size().

See also:
operator=(), readRows()

void readRows ( const T *  rows  )  [inline, inherited]

Reads the rows of the matrix from a row-dominant array. For instance, C/C++ two-dimensional arrays are stored in row-dominant order.

For better performance, use readArray()

See also:
operator=(), readArray()

void recursiveGramSchmidt ( int  n  )  [inline, inherited]

Gram-Schmidt algorithm.

Assuming that the n first columns of *this are mutually orthogonal unit vectors, this method fills the remaining columns with unit vectors such that all the columns are mutually orthogonal. Thus, after this method has returned, the columns of *this form an orthonormal basis. In other words, after this function has returned, if we let $A$ denote the matrix *this, then

\[ AA^* = \mathrm{Id}, \]

which is to say that $A$ is an orthogonal matrix. Note that there is an inconsistency in the mathematical language: an orthogonal matrix is a matrix whose columns form an orthonormal basis, not just an orthogonal basis.

Note:
this method assumes that the n first column vectors are unit vectors. Otherwise, it doesn't work.
See also:
loadOrthoBasis()

void replaceWithAdjoint (  )  [inline, inherited]

Replaces *this with its adjoint (conjugate transpose, equals transpose unless T is complex numbers).

See also:
adjoint(), computeAdjoint()

void replaceWithOpposite (  )  [inline, inherited]

Replaces *this with (-(*this)).

See also:
operator-(void) const

void resize ( int  newsize  )  [inline, inherited]

Resizes the matrix. That is only possible if the matrix has dynamic size, i.e. is an object of class MatrixX.

Resizing a fixed-size matrix is not possible, and attempting to do so will only generate a debug message (unless the new size equals the old one).

The matrix entries are not kept, they are left with undefined values after resizing.

See also:
hasDynamicSize(), size()

Derived& rotate2 ( const T &  angle  )  [inline, inherited]

Multiplies *this on the right by the rotation matrix of given angle in radians. The template parameter Size must equal 2.

See also:
loadRotation2(), prerotate2()

Derived& rotate3 ( const T &  angle,
const VectorType &  axis 
) [inline, inherited]

Multiplies *this on the right by the rotation matrix of given angle in radians around given axis vector. Only applicable to matrices of size 3.

See also:
prerotate3(), loadRotation3()

VectorType row ( int  row  )  const [inline, inherited]

Returns a row of the matrix.

This method returns an object by value, which is inefficient. For better performance, use another method.

See also:
getRow(int,T*) const, getRow(int, VectorType *) const

Derived& scale ( const VectorType &  coeffs  )  [inline, inherited]

Multiplies *this on the right by the scaling matrix with given vector of scaling coefficients.

See also:
loadScaling(const VectorType &), prescale(const VectorType &)

void setColumn ( int  column,
const VectorType &  src 
) [inline, inherited]

Copies the vector src into a column of the matrix

See also:
setColumn( int, const T *)

void setColumn ( int  column,
const T *  src 
) [inline, inherited]

Copies the array src into a column of the matrix

See also:
setColumn( int, const VectorType &)

void setRow ( int  row,
const VectorType &  src 
) [inline, inherited]

Copies the vector src into a row of the matrix

See also:
setRow( int, const T *)

void setRow ( int  row,
const T *  src 
) [inline, inherited]

Copies the array src into a row of the matrix

See also:
setRow( int, const VectorType &)

int size (  )  const [inline, inherited]

Returns:
the size (number of rows, or equivalently number of columns) of the matrix.
See also:
hasDynamicSize(), resize()


The documentation for this class was generated from the following file:

Generated on Sun May 17 11:14:22 2009 for Eigen by  doxygen 1.5.8