ACADO Toolkit  1.2.0beta
Toolkit for Automatic Control and Dynamic Optimization
Public Member Functions | Protected Member Functions | Protected Attributes | Friends

Implements a rudimentary dense matrix class. More...

#include <matrix.hpp>

Inheritance diagram for Matrix:
Inheritance graph
[legend]

List of all members.

Public Member Functions

 Matrix ()
 Default constructor.
 Matrix (uint _nRows, uint _nCols)
 Constructor which takes the dimensions of the matrix.
 Matrix (uint _nRows, uint _nCols, const double *const _values)
 Constructor which takes dimension of the matrix and a double array
of appropriate size (_nRows*_nCols) containing the (initial) values
of the matrix.
 Matrix (FILE *file)
 Constructor which takes a file.
 Matrix (int value)
 Constructor which takes an integer to construct a 1x1-matrix.
 Matrix (double value)
 Constructor which takes a double to construct a 1x1-matrix.
 Matrix (const Vector &value, BooleanType doTranspose=BT_FALSE)
 Constructor to convert a Vector into a matrix.
 Matrix (const VariablesGrid &value)
 Constructor to convert a VariabledGrid into a matrix.
 Matrix (const Matrix &rhs)
 Copy constructor (deep copy).
virtual ~Matrix ()
 Destructor.
returnValue init (uint _nRows, uint _nCols)
 Initializer, which takes the dimensions of the matrix.
returnValue init (uint _nRows, uint _nCols, double *_values)
 Initializer,which takes dimension of the matrix and a double array
of appropriate size (_nRows*_nCols) containing the (initial)
values of the matrix.
returnValue appendRows (const Matrix &arg)
 Appends rows at the end of the matrix.
returnValue appendCols (const Matrix &arg)
 Appends columns at the end of the matrix.
Matrixoperator= (const Matrix &rhs)
 Assignment operator (deep copy).
Matrixoperator= (FILE *rhs)
 Assignment operator, which loads a matrix from a file.
Matrixoperator^= (const double *rhs)
 Assignment operator, which assigns the elements from a double*.
double & operator() (uint rowIdx, uint colIdx)
 Access operator that return the value of a certain component.
double operator() (uint rowIdx, uint colIdx) const
 Access operator that return the value of a certain component (const variant).
Matrix operator+ (const Matrix &arg) const
 Adds (element-wise) two matrices to a temporary object.
Matrixoperator+= (const Matrix &arg)
 Adds (element-wise) a matrix to object.
Vector sumCol (const Matrix &arg)
 Computes the column-wise sum the Matrix.
Vector sumRow (const Matrix &arg)
 Computes the row-wise sum the Matrix.
Matrix operator- (const Matrix &arg) const
 Subtracts (element-wise) a matrix from the object and and stores the result to a temporary object.
Matrixoperator-= (const Matrix &arg)
 Subtracts (element-wise) a matrix from the object.
Matrixoperator*= (double scalar)
 Multiplies each component of the object with a given scalar.
Matrixoperator/= (double scalar)
 Divides each component of the object by a given non-zero scalar.
BooleanType operator== (const Matrix &arg) const
 Tests for equality.
Matrix operator* (const Matrix &arg) const
 Multiplies a matrix from the right to the matrix object and stores the result to a temporary object.
Matrix operator^ (const Matrix &arg) const
 Multiplies a matrix from the right to the transposed matrix object and stores the result to a temporary object.
Vector operator* (const Vector &arg) const
 Multiplies a vector from the right to the matrix object and stores the result to a temporary object.
Vector operator^ (const Vector &arg) const
 Multiplies a vector from the right to the transposed matrix object and stores the result to a temporary object.
Matrix transpose () const
Matrix negativeTranspose () const
MatrixmakeVector ()
uint getNumRows () const
 Returns number of rows of the matrix object.
uint getNumCols () const
 Returns number of columns of the matrix object.
Vector getRow (uint idx) const
 Returns a given row of the matrix object.
Vector getCol (uint idx) const
 Returns a given column of the matrix object.
returnValue setRow (uint idx, const Vector &arg)
 Assigns new values to a given row of the matrix object.
returnValue setCol (uint idx, const Vector &arg)
 Assigns new values to a given column of the matrix object.
Matrix getRows (uint idx1, uint idx2) const
 Returns given rows of the matrix object.
Matrix getCols (uint idx1, uint idx2) const
 Returns given columns of the matrix object.
Vector getDiag () const
 Returns a vector containing the diagonal elements of a square matrix.
returnValue setIdentity ()
 Sets object to the identity matrix.
BooleanType isSquare () const
 Tests if object is a square matrix.
BooleanType isSymmetric () const
 Tests if object is a symmetric matrix.
BooleanType isDiagonal () const
 Tests if object is a diagonal matrix.
returnValue symmetrize ()
BooleanType isPositiveSemiDefinite () const
 Tests if object is a positive semi-definite matrix.
BooleanType isPositiveDefinite () const
 Tests if object is a (strictly) positive definite matrix.
Matrix absolute ()
 Returns the a matrix whose components are the absolute values of the components of this object.
Matrix positive ()
 Returns the a matrix whose components are equal to the components of this object, if they are positive or zero, but zero otherwise.
Matrix negative ()
 Returns the a matrix whose components are equal to the components of this object, if they are negative or zero, but zero otherwise.
Matrix minus ()
 Switches the sign of all components.
double getNorm (MatrixNorm norm) const
 Returns specified norm of the matrix.
double getTrace () const
 Returns trace of the matrix.
Vector getEigenvalues () const
 Returns a vector containing the eigenvalues of the matrix.
Vector getEigenvalues (Matrix &Q) const
 Returns a vector containing the eigenvalues of the matrix as
well as the corresponding eigenvectors.
returnValue printEigenvalues () const
 Prints the eigenvalues of the standard output stream.
returnValue getSingularValueDecomposition (Matrix &U, Vector &D, Matrix &V) const
 Computes the Singular Value Decomposition of the matrix and
stores the result in the form

A = U D V^T ,

where the matrices U and V are orthogonal.
Matrix getInverse () const
 Computes the inverse matrix.
Matrix getCholeskyDecomposition () const
 Computes the Cholesky Decomposition of the matrix

A = L L^T

.
Matrix getCholeskyDecomposition (Vector &D) const
 Computes the Cholesky Decomposition of the matrix

A = L D L^T

.
Matrix getCholeskyInverse () const
 Computes the inverse of a positive definite matrix based on
the Cholesky decomposition.
returnValue computeQRdecomposition ()
 Computes the QR decomposition of the matrix.
Vector solveQR (const Vector &b) const
 Solves the system A x = b provided that the routine
computeQRdecomposition() has been used before.
Vector solveTransposeQR (const Vector &b) const
 Solves the system A^T x = b provided that the routine
computeQRdecomposition() has been used before.
returnValue computeSparseLUdecomposition ()
 Computes the sparse LU decomposition of the matrix.
Vector solveSparseLU (const Vector &b) const
 Solves the system A x = b provided that the routine
computeQRdecomposition() has been used before.
Vector solveTransposeSparseLU (const Vector &b) const
 Solves the system A^T x = b provided that the routine
computeQRdecomposition() has been used before.
virtual returnValue printToFile (const char *const filename, const char *const name=DEFAULT_LABEL, const char *const startString=DEFAULT_START_STRING, const char *const endString=DEFAULT_END_STRING, uint width=DEFAULT_WIDTH, uint precision=DEFAULT_PRECISION, const char *const colSeparator=DEFAULT_COL_SEPARATOR, const char *const rowSeparator=DEFAULT_ROW_SEPARATOR) const
 Prints object to file with given name.
virtual returnValue printToFile (FILE *file, const char *const name=DEFAULT_LABEL, const char *const startString=DEFAULT_START_STRING, const char *const endString=DEFAULT_END_STRING, uint width=DEFAULT_WIDTH, uint precision=DEFAULT_PRECISION, const char *const colSeparator=DEFAULT_COL_SEPARATOR, const char *const rowSeparator=DEFAULT_ROW_SEPARATOR) const
 Prints object to given file.
virtual returnValue printToFile (const char *const filename, const char *const name, PrintScheme printScheme) const
 Prints object to file with given name.
virtual returnValue printToFile (FILE *file, const char *const name, PrintScheme printScheme) const
 Prints object to given file.
virtual returnValue printToString (char **string, const char *const name=DEFAULT_LABEL, const char *const startString=DEFAULT_START_STRING, const char *const endString=DEFAULT_END_STRING, uint width=DEFAULT_WIDTH, uint precision=DEFAULT_PRECISION, const char *const colSeparator=DEFAULT_COL_SEPARATOR, const char *const rowSeparator=DEFAULT_ROW_SEPARATOR, BooleanType allocateMemory=BT_TRUE) const
 Prints object to given string.
virtual returnValue printToString (char **string, const char *const name, PrintScheme printScheme, BooleanType allocateMemory=BT_TRUE) const
 Prints object to given string.
virtual uint determineStringLength (const char *const name, const char *const startString, const char *const endString, uint width, uint precision, const char *const colSeparator, const char *const rowSeparator) const
 Determines length of string required for printing object with given settings defining its output format.
VectorspaceElementoperator<< (double *rhs)
 Assignment operator (deep copy).
returnValue convert (double *lhs) const
returnValue init (uint _dim=0, double *_values=0)
 Initializes vector space element with values taken from a double array of appropriate size.
double & operator() (uint idx)
 Access operator that return the value of a certain component.
double operator() (uint idx) const
 Access operator that return the value of a certain component (const variant).
BooleanType operator== (const VectorspaceElement &arg) const
 Tests for element-wise equality.
BooleanType operator!= (const VectorspaceElement &arg) const
 Tests for non-equality.
BooleanType operator< (const VectorspaceElement &arg) const
 Tests for element-wise smaller-than.
BooleanType operator<= (const VectorspaceElement &arg) const
 Tests for element-wise smaller/equal-than.
BooleanType operator> (const VectorspaceElement &arg) const
 Tests for element-wise greater-than.
BooleanType operator>= (const VectorspaceElement &arg) const
 Tests for element-wise greater/equal-than.
returnValue append (const VectorspaceElement &arg)
uint getDim () const
 Returns dimension of vector space.
BooleanType isEmpty () const
 Returns whether the vectorspace element is empty.
BooleanType isEqualTo (double _value) const
 Returns whether all components are equal to _value.
BooleanType isGreaterThan (double _value) const
BooleanType isSmallerThan (double _value) const
BooleanType isZero () const
BooleanType isPositive () const
BooleanType isNegative () const
BooleanType isFinite () const
BooleanType hasNaN () const
BooleanType hasEqualComponents () const
 Returns whether all components have same value.
returnValue setZero ()
 Sets all component to zero.
returnValue setAll (double _value)
 Sets all component to given value.
double getMax () const
 Returns maximum element.
double getMin () const
 Returns minimum element.
double getMean () const
 Returns mean value of all elements.
double getNorm (VectorNorm norm) const
 Returns specified norm of the vectorspace element
interpreted as a vector.
double getNorm (VectorNorm norm, const VectorspaceElement &scale) const
 Returns specified norm of the vectorspace element
interpreted as a vector (with scaling).
returnValue print (const char *const name=DEFAULT_LABEL, const char *const startString=DEFAULT_START_STRING, const char *const endString=DEFAULT_END_STRING, uint width=DEFAULT_WIDTH, uint precision=DEFAULT_PRECISION, const char *const colSeparator=DEFAULT_COL_SEPARATOR, const char *const rowSeparator=DEFAULT_ROW_SEPARATOR) const
 Prints object to standard ouput stream.
returnValue print (const char *const name, PrintScheme printScheme) const
 Prints object to standard ouput stream.
double * getDoublePointer ()

Protected Member Functions

void householdersReductionToBidiagonalForm (double *A, int nrows, int ncols, double *U, double *V, double *diagonal, double *superdiagonal) const
int givensReductionToDiagonalForm (int nrows, int ncols, double *U, double *V, double *diagonal, double *superdiagonal) const
void sortByDecreasingSingularValues (int nrows, int ncols, double *singular_value, double *U, double *V) const
void lowerTriangularInverse (double *L, int n) const

Protected Attributes

uint nRows
 Number of rows.
uint nCols
 Number of columns.
SparseSolversolver
 Optional sparse solver.
double * element
 Element of vector space.
unsigned int dim
 Vector space dimension.

Friends

double * operator^= (double *lhs, Matrix &rhs)
 Assignment operator, which assigns the elements from a double*.
returnValue operator<< (FILE *file, Matrix &arg)
 Prints the matrix into a file.
double * operator<< (double *lhs, VectorspaceElement &rhs)
 Assignment operator (deep copy).
returnValue operator<< (FILE *file, VectorspaceElement &arg)
 Prints the object into a file.

Detailed Description

It is intended to provide a convenient way to deal with linear algebra
objects. The operators that are defined for the Matrix and Vector class
allow a convenient usage of these objects.

Example Code:

See also:
Vector
Author:
Hans Joachim Ferreau, Boris Houska

Constructor & Destructor Documentation

References nCols, nRows, and solver.

Referenced by MatrixVariable::getMatrix().

Matrix::Matrix ( uint  _nRows,
uint  _nCols 
)



Parameters:
_nRowsNumber of rows.
_nColsNumber of columns.

References nCols, nRows, and solver.

Matrix::Matrix ( uint  _nRows,
uint  _nCols,
const double *const  _values 
)


Parameters:
_nRowsNumber of rows.
_nColsNumber of columns.
_nColsDouble array.

References nCols, nRows, and solver.

Matrix::Matrix ( FILE *  file)



Example:

            Matrix A = fopen("matrix.dat", "r");  // loading a matrix from the file "matrix.dat"
            A.print("A");                         // printing the matrix.
            



Parameters:
fileA file from which the matrix should be loaded.

References operator=(), and solver.

Matrix::Matrix ( int  value)
Parameters:
valueThe value of the matrix element

References nCols, nRows, operator()(), and solver.

Matrix::Matrix ( double  value)
Parameters:
valueThe value of the matrix element

References nCols, nRows, operator()(), and solver.

Matrix::Matrix ( const Vector value,
BooleanType  doTranspose = BT_FALSE 
)
Parameters:
valueThe vector containing the elements for the dim-by-1 matrix.

References BT_FALSE, VectorspaceElement::getDim(), nCols, nRows, and solver.

Matrix::Matrix ( const VariablesGrid value)
Parameters:
valueThe grid values for the matrix.

References Grid::getNumPoints(), MatrixVariablesGrid::getNumValues(), Grid::getTime(), nCols, nRows, operator()(), solver, and uint.

Matrix::Matrix ( const Matrix rhs)
Parameters:
rhsRight-hand side object.

References SparseSolver::clone(), nCols, nRows, and solver.

Matrix::~Matrix ( ) [virtual]

References solver.


Member Function Documentation

Matrix Matrix::absolute ( ) [inline]


The result of the decomposition will be strored in this
matrix which will be increased by one column. The matrix
can then be used via the routine "solveQR".

Returns:
The matrix decomposition A = QR in an efficient
storage format.

References ACADOERROR, appendRows(), ASSERT, EPS, getNumCols(), getNumRows(), operator()(), RET_DIV_BY_ZERO, solver, sqrt(), and SUCCESSFUL_RETURN.

Referenced by IntegratorBDF::decomposeJacobian().

returnValue VectorspaceElement::convert ( double *  lhs) const [inline, inherited]
uint Matrix::determineStringLength ( const char *const  name,
const char *const  startString,
const char *const  endString,
uint  width,
uint  precision,
const char *const  colSeparator,
const char *const  rowSeparator 
) const [virtual]
Parameters:
[in]nameName label to be printed before the numerical values.
[in]startStringPrefix before printing the numerical values.
[in]endStringSuffix after printing the numerical values.
[in]widthTotal number of digits per single numerical value.
[in]precisionNumber of decimals per single numerical value.
[in]colSeparatorSeparator between the columns of the numerical values.
[in]rowSeparatorSeparator between the rows of the numerical values.
Returns:
Length of string required for printing object

Reimplemented from VectorspaceElement.

References getNumCols(), getNumRows(), getStringLength(), and uint.

Referenced by LogRecordItem::determineStringLength(), and printToString().



Returns:
The matrix A^{-1} = L^{-T}L^{-1}

References ASSERT, getCholeskyDecomposition(), getNumCols(), getNumRows(), and lowerTriangularInverse().

Vector Matrix::getCol ( uint  idx) const [inline]
Returns:
Temporary object containing the column.
Parameters:
idxIndex of the column to be returned.

References ASSERT, getNumCols(), getNumRows(), operator()(), and uint.

Referenced by ShootingMethod::differentiateForward(), ShootingMethod::differentiateForwardBackward(), IntegrationAlgorithm::getForwardSensitivities(), and VariablesGrid::getVector().

Matrix Matrix::getCols ( uint  idx1,
uint  idx2 
) const [inline]
Returns:
Temporary object containing the columns.
Parameters:
idx1Start index of the columns to be returned.
idx2End index of the columns to be returned.

Reimplemented in MatrixVariable.

References getNumCols(), getNumRows(), init(), operator()(), and uint.

Vector Matrix::getDiag ( ) const [inline]
Returns:
Vector containing the diagonal elements.

References ASSERT, BT_TRUE, getNumRows(), isSquare(), operator()(), and uint.

Referenced by ExportGaussNewtonForces::setupObjectiveEvaluation().

unsigned int VectorspaceElement::getDim ( ) const [inline, inherited]
Returns:
Dimension of vector space.

References VectorspaceElement::dim.

Referenced by Vector::absolute(), absolute(), Function::AD_backward(), Curve::add(), Constraint::add(), DifferentialEquation::addAlgebraic(), DifferentialEquation::addDifferential(), VectorspaceElement::append(), Process::calculateOutput(), ClippingFunctionality::clipSignals(), ImplicitRungeKuttaExport::computeCombinations(), CondensingBasedCPsolver::computeCondensingOperator(), CondensingBasedCPsolver::condense(), VectorspaceElement::convert(), Expression::convert(), EvaluationPoint::copy(), IntegratorBDF::copyBackward(), VectorspaceElement::determineStringLength(), diag(), Integrator::diffTransitionBackward(), Integrator::diffTransitionForward(), IntegratorBDF::evaluate(), IntegratorRK::evaluate(), IntegratorLYAPUNOV::evaluate(), SCPevaluation::evaluateLagrangeGradient(), LSQTerm::evaluateSensitivities(), LSQTerm::evaluateSensitivitiesGN(), LSQEndTerm::evaluateSensitivitiesGN(), Integrator::evaluateTransition(), CondensingBasedCPsolver::expand(), SIMexport::exportAcadoHeader(), Expression::Expression(), RealTimeAlgorithm::feedbackStep(), GaussianNoise::GaussianNoise(), WeightGeneration::generateWeights(), Vector::getAbsolute(), Constraint::getBackwardSensitivities(), Integrator::getBackwardSensitivities(), Actuator::getControlDeadTime(), CondensingBasedCPsolver::getFirstControl(), Constraint::getForwardSensitivities(), VectorspaceElement::getMax(), VectorspaceElement::getMean(), VectorspaceElement::getMin(), GnuplotWindow::getMouseEvent(), VectorspaceElement::getNorm(), Estimator::getNP(), Estimator::getNU(), ClippingFunctionality::getNumControlLimits(), ClippingFunctionality::getNumParameterLimits(), MatrixVariablesGrid::getNumValues(), DenseCP::getNV(), Estimator::getNW(), Estimator::getNX(), Estimator::getNXA(), Sensor::getOutputDeadTime(), Actuator::getParameterDeadTime(), CondensingBasedCPsolver::getParameters(), MultiObjectiveAlgorithm::getParetoFront(), MultiObjectiveAlgorithm::getParetoFrontWithFilter(), IntegratorLYAPUNOV::getProtectedBackwardSensitivities(), IntegratorRK::getProtectedBackwardSensitivities(), BFGSupdate::getSubBlockLine(), MultiObjectiveAlgorithm::getWeightsWithFilter(), TransferDevice::hasDeadTime(), VectorspaceElement::hasEqualComponents(), VectorspaceElement::hasNaN(), Grid::init(), TransferDevice::init(), PIDcontroller::init(), UniformNoise::init(), GaussianNoise::init(), Process::init(), RungeKuttaExport::initializeButcherTableau(), ImplicitRungeKuttaExport::initializeCoefficients(), OptimizationAlgorithmBase::initializeOCPiterate(), Process::initializeStartValues(), Integrator::integrate(), Integrator::integrateSensitivities(), IntermediateState::IntermediateState(), VectorspaceElement::isEmpty(), VectorspaceElement::isEqualTo(), VectorspaceElement::isFinite(), VectorspaceElement::isGreaterThan(), VectorspaceElement::isSmallerThan(), DenseQPsolver::makeBoundsConsistent(), makeVector(), Matrix(), negative(), OCP::OCP(), VectorspaceElement::operator!=(), Vector::operator%(), VectorspaceElement::operator()(), Vector::operator*(), operator*(), Vector::operator*=(), operator*=(), Vector::operator+(), operator+(), Vector::operator+=(), operator+=(), Vector::operator-(), operator-(), operator-(), Vector::operator-=(), operator-=(), Vector::operator/=(), operator/=(), VectorspaceElement::operator<(), operator<<(), VectorspaceElement::operator<=(), MatrixVariable::operator=(), Transition::operator==(), VectorspaceElement::operator==(), operator==(), VectorspaceElement::operator>(), VectorspaceElement::operator>=(), Vector::operator^(), operator^(), positive(), VectorspaceElement::printToString(), CondensingBasedCPsolver::projectHessian(), Process::projectToComponents(), Process::run(), VectorspaceElement::setAll(), Integrator::setBackwardSeed(), IntegratorBDF::setBackwardSeed2(), IntegratorRK::setBackwardSeed2(), IntegratorLYAPUNOV::setBackwardSeed2(), setCol(), ExportNLPSolver::setConstraints(), Actuator::setControlDeadTime(), Actuator::setControlDeadTimes(), ClippingFunctionality::setControlLowerLimits(), ClippingFunctionality::setControlUpperLimits(), PIDcontroller::setDerivativeWeights(), Integrator::setForwardSeed(), IntegratorBDF::setForwardSeed2(), IntegratorLYAPUNOV::setForwardSeed2(), IntegratorRK::setForwardSeed2(), LSQTerm::setGrid(), PIDcontroller::setIntegralWeights(), ModelData::setIntegrationGrid(), ConstraintComponent::setLB(), UniformNoise::setLimits(), VariableSettings::setLowerBounds(), GaussianNoise::setMeans(), Sensor::setOutputDeadTime(), Sensor::setOutputDeadTimes(), Actuator::setParameterDeadTime(), Actuator::setParameterDeadTimes(), ClippingFunctionality::setParameterLowerLimits(), ClippingFunctionality::setParameterUpperLimits(), PIDcontroller::setProportionalWeights(), IntegratorBDF::setProtectedBackwardSeed(), IntegratorLYAPUNOV::setProtectedBackwardSeed(), IntegratorRK::setProtectedBackwardSeed(), IntegratorBDF::setProtectedForwardSeed(), IntegratorLYAPUNOV::setProtectedForwardSeed(), IntegratorRK::setProtectedForwardSeed(), DenseCP::setQPsolution(), setRow(), VariableSettings::setScaling(), BFGSupdate::setSubBlockLine(), ConstraintComponent::setUB(), Vector::setUnitVector(), ExplicitRungeKuttaExport::setup(), ExportGaussNewtonQpDunes::setupConstraintsEvaluation(), ExportGaussNewtonForces::setupConstraintsEvaluation(), VariableSettings::setUpperBounds(), GaussianNoise::setVariances(), VariablesGrid::setVector(), VariablesGrid::shiftBackwards(), Process::simulate(), solveQR(), MultiObjectiveAlgorithm::solveSingleObjective(), solveTransposeQR(), FeedforwardLaw::step(), LinearStateFeedback::step(), PIDcontroller::step(), UniformNoise::step(), GaussianNoise::step(), Process::step(), UniformNoise::UniformNoise(), and MatFile::write().

double * VectorspaceElement::getDoublePointer ( ) [inherited]


Note that this routine requires that the matrix is symmetric.

Returns:
The eigenvalues of the matrix.

Referenced by isPositiveDefinite(), isPositiveSemiDefinite(), printEigenvalues(), and CondensingBasedCPsolver::projectHessian().

The eigenvectors are
stored in the matrix Q (orthogonal matrix).
Note that this routine requires that the matrix is symmetric.

Returns:
The eigenvalues of the matrix.

References ASSERT, VectorspaceElement::dim, VectorspaceElement::element, EPS, fabs(), getNumCols(), getNumRows(), init(), max(), and sqrt().

double VectorspaceElement::getMax ( ) const [inline, inherited]
double VectorspaceElement::getMean ( ) const [inline, inherited]
Returns:
Mean value of all elements.

References VectorspaceElement::element, VectorspaceElement::getDim(), and uint.

double VectorspaceElement::getMin ( ) const [inline, inherited]
double VectorspaceElement::getNorm ( VectorNorm  norm) const [inherited]



Parameters:
normthe type of norm to be computed.

Returns:
Norm of the vector.

References VectorspaceElement::dim, and VectorspaceElement::setAll().

Referenced by IntegratorBDF::applyNewtonStep().

double VectorspaceElement::getNorm ( VectorNorm  norm,
const VectorspaceElement scale 
) const [inherited]



Parameters:
normthe type of norm to be computed.
scalethe elementwise scale.

Returns:
Norm of the vector.

References fabs(), VectorspaceElement::getDim(), pow(), sqrt(), uint, VN_L1, VN_L2, and VN_LINF.

double Matrix::getNorm ( MatrixNorm  norm) const
Returns:
Norm of the matrix.

References fabs(), getNumCols(), getNumRows(), MN_COLUMN_SUM, MN_FROBENIUS, MN_ROW_SUM, sqrt(), and uint.

uint Matrix::getNumCols ( ) const [inline]
Returns:
Number of columns.

References nCols.

Referenced by DifferentialEquation::addAlgebraic(), DifferentialEquation::addDifferential(), appendCols(), appendRows(), PointConstraint::computeForwardSensitivityBlock(), AlgebraicConsistencyConstraint::computeForwardSensitivityBlock(), computeQRdecomposition(), computeSparseLUdecomposition(), CondensingBasedCPsolver::condense(), Expression::convert(), determineStringLength(), ShootingMethod::differentiateForward(), ShootingMethod::differentiateForwardBackward(), Expression::Expression(), DiagonallyImplicitRKExport::formMatrix(), ImplicitRungeKuttaExport::formMatrix(), WeightGeneration::generateWeights(), getCholeskyDecomposition(), getCholeskyInverse(), getCol(), MatrixVariable::getCols(), getCols(), getEigenvalues(), IntegrationAlgorithm::getForwardSensitivities(), getInverse(), getNorm(), BlockMatrix::getNumCols(), MatrixVariablesGrid::getNumCols(), LinearStateFeedback::getNX(), getRow(), MatrixVariable::getRows(), getRows(), getSingularValueDecomposition(), getTrace(), IntermediateState::IntermediateState(), isDiagonal(), isSquare(), minus(), negativeTranspose(), operator()(), Vector::operator*(), operator*(), operator+(), operator+=(), BlockMatrix::operator-(), operator-(), operator-(), operator-=(), MatrixVariablesGrid::operator=(), Transition::operator==(), operator==(), operator^(), operator^=(), MultiObjectiveAlgorithm::printAuxiliaryRoutine(), printToString(), setCol(), IntegratorExport::setLinearInput(), ModelData::setLinearInput(), IntegratorExport::setLinearOutput(), ModelData::setLinearOutput(), NARXExport::setLinearOutput(), ModelData::setNARXmodel(), setRow(), MultiObjectiveAlgorithm::solve(), solveQR(), MultiObjectiveAlgorithm::solveSingleObjective(), solveTransposeQR(), sumCol(), transpose(), ShootingMethod::update(), VariablesGrid::VariablesGrid(), and MatFile::write().

uint Matrix::getNumRows ( ) const [inline]
Returns:
Number of rows.

References nRows.

Referenced by DifferentialEquation::addAlgebraic(), DifferentialEquation::addDifferential(), ShootingMethod::addStage(), appendCols(), appendRows(), PointConstraint::computeForwardSensitivityBlock(), AlgebraicConsistencyConstraint::computeForwardSensitivityBlock(), computeQRdecomposition(), computeSparseLUdecomposition(), CondensingBasedCPsolver::condense(), Expression::convert(), determineStringLength(), diag(), ShootingMethod::differentiateBackward(), PointConstraint::evaluateSensitivities(), BoundaryConstraint::evaluateSensitivities(), CoupledPathConstraint::evaluateSensitivities(), CondensingBasedCPsolver::expand(), Expression::Expression(), DiagonallyImplicitRKExport::formMatrix(), ImplicitRungeKuttaExport::formMatrix(), WeightGeneration::generateWeights(), IntegrationAlgorithm::getBackwardSensitivities(), getCholeskyDecomposition(), getCholeskyInverse(), getCol(), MatrixVariable::getCols(), getCols(), getDiag(), getEigenvalues(), getInverse(), DenseCP::getNC(), getNorm(), LinearStateFeedback::getNU(), BlockMatrix::getNumRows(), MatrixVariablesGrid::getNumRows(), DenseCP::getNV(), getRow(), MatrixVariable::getRows(), getRows(), getSingularValueDecomposition(), MultiObjectiveAlgorithm::getWeightsWithFilter(), RungeKuttaExport::initializeButcherTableau(), IntermediateState::IntermediateState(), isDiagonal(), isSquare(), isSymmetric(), minus(), negativeTranspose(), operator()(), Vector::operator*(), operator*(), operator+(), operator+=(), BlockMatrix::operator-(), operator-(), operator-(), operator-=(), MatrixVariablesGrid::operator=(), Transition::operator==(), operator==(), operator^(), operator^=(), printToString(), setCol(), setIdentity(), IntegratorExport::setLinearInput(), ModelContainer::setLinearInput(), IntegratorExport::setLinearOutput(), ModelContainer::setLinearOutput(), NARXExport::setLinearOutput(), IntegratorExport::setModelData(), ModelData::setNARXmodel(), NARXExport::setNARXmodel(), setRow(), MultiObjectiveAlgorithm::solve(), CondensingBasedCPsolver::solveCPsubproblem(), solveSparseLU(), solveTransposeSparseLU(), sumRow(), symmetrize(), transpose(), VariablesGrid::VariablesGrid(), and MatFile::write().

Vector Matrix::getRow ( uint  idx) const [inline]
Matrix Matrix::getRows ( uint  idx1,
uint  idx2 
) const [inline]
Returns:
Temporary object containing the rows.
Parameters:
idx1Start index of the rows to be returned.
idx2End index of the rows to be returned.

Reimplemented in MatrixVariable.

References getNumCols(), getNumRows(), init(), operator()(), and uint.

Referenced by DiagonallyImplicitRKExport::evaluateRhsImplicitSystem(), ImplicitRungeKuttaExport::evaluateRhsImplicitSystem(), and ExportGaussNewtonCondensed::setupConstraintsEvaluation().

Note that
diagonal matrix D will be stored as a vector containing the
diagonal elements.

Returns:
SUCCESSFUL_RETURN

References VectorspaceElement::dim, VectorspaceElement::element, getNumCols(), getNumRows(), givensReductionToDiagonalForm(), householdersReductionToBidiagonalForm(), VectorspaceElement::init(), init(), sortByDecreasingSingularValues(), and SUCCESSFUL_RETURN.

Referenced by getInverse().

double Matrix::getTrace ( ) const
Returns:
Trace of the matrix.

References BT_FALSE, getNumCols(), INFTY, isSquare(), and uint.

int Matrix::givensReductionToDiagonalForm ( int  nrows,
int  ncols,
double *  U,
double *  V,
double *  diagonal,
double *  superdiagonal 
) const [protected]

References EPS, fabs(), and sqrt().

Referenced by getSingularValueDecomposition().

BooleanType VectorspaceElement::hasEqualComponents ( ) const [inline, inherited]
Returns:
BT_TRUE iff all components have same value, BT_FALSE otherwise.

References acadoIsEqual(), BT_FALSE, BT_TRUE, VectorspaceElement::getDim(), VectorspaceElement::operator()(), and uint.

BooleanType VectorspaceElement::hasNaN ( ) const [inline, inherited]
void Matrix::householdersReductionToBidiagonalForm ( double *  A,
int  nrows,
int  ncols,
double *  U,
double *  V,
double *  diagonal,
double *  superdiagonal 
) const [protected]

References EPS, fabs(), and sqrt().

Referenced by getSingularValueDecomposition().

returnValue VectorspaceElement::init ( uint  _dim = 0,
double *  _values = 0 
) [inherited]

Previously allocated internal memory is freed.

Returns:
SUCCESSFUL_RETURN
Parameters:
_dimVector space dimension.

References VectorspaceElement::dim, VectorspaceElement::element, SUCCESSFUL_RETURN, and uint.

Referenced by IntegratorLYAPUNOV::allocateMemory(), IntegratorRK::allocateMemory(), IntegratorBDF::allocateMemory(), VectorspaceElement::append(), VariableSettings::appendSettings(), ClippingFunctionality::ClippingFunctionality(), CondensingBasedCPsolver::condense(), ConstraintComponent::ConstraintComponent(), IntegratorBDF::constructAll(), PIDcontroller::determineControlAction(), Curve::evaluate(), CondensingBasedCPsolver::expand(), GaussianNoise::GaussianNoise(), IntegrationAlgorithm::getBackwardSensitivities(), getCholeskyDecomposition(), QPsolver_qpOASES::getDualSolution(), QPsolver_qpOASES::getPrimalSolution(), getSingularValueDecomposition(), Integrator::getX(), Integrator::getXA(), FeedforwardLaw::init(), PIDcontroller::init(), RealTimeAlgorithm::init(), Process::init(), CondensingBasedCPsolver::initializeCPsolver(), MatrixVariablesGrid::initMatrixVariables(), Integrator::integrate(), Integrator::integrateSensitivities(), Vector::operator=(), VariableSettings::operator=(), PIDcontroller::PIDcontroller(), Integrator::setBackwardSeed(), IntegratorBDF::setBackwardSeed2(), IntegratorRK::setBackwardSeed2(), IntegratorLYAPUNOV::setBackwardSeed2(), Integrator::setForwardSeed(), IntegratorBDF::setForwardSeed2(), IntegratorRK::setForwardSeed2(), IntegratorLYAPUNOV::setForwardSeed2(), ConstraintComponent::setLB(), VariableSettings::setLowerBound(), IntegratorBDF::setProtectedBackwardSeed(), IntegratorLYAPUNOV::setProtectedBackwardSeed(), IntegratorRK::setProtectedBackwardSeed(), IntegratorBDF::setProtectedForwardSeed(), IntegratorLYAPUNOV::setProtectedForwardSeed(), IntegratorRK::setProtectedForwardSeed(), VariableSettings::setScaling(), ConstraintComponent::setUB(), VariableSettings::setUpperBound(), CondensingBasedCPsolver::setupRelaxedQPdataL1(), CondensingBasedCPsolver::setupRelaxedQPdataL2(), PIDcontroller::step(), TransferDevice::TransferDevice(), UniformNoise::UniformNoise(), and VariableSettings::VariableSettings().

returnValue Matrix::init ( uint  _nRows,
uint  _nCols 
)


Parameters:
_nRowsNumber of rows.
_nColsNumber of columns.

References nCols, nRows, solver, and SUCCESSFUL_RETURN.

Referenced by IntegratorBDF::allocateMemory(), appendCols(), appendRows(), ShootingMethod::clear(), CondensingBasedCPsolver::condense(), IntegratorBDF::constructAll(), IntegratorBDF::determineCorrector(), ShootingMethod::differentiateBackward(), ShootingMethod::differentiateForward(), ShootingMethod::differentiateForwardBackward(), BoxConstraint::evaluateBounds(), MayerTerm::evaluateSensitivities(), BoundaryConstraint::evaluateSensitivities(), PointConstraint::evaluateSensitivities(), LSQTerm::evaluateSensitivities(), PathConstraint::evaluateSensitivities(), CoupledPathConstraint::evaluateSensitivities(), AlgebraicConsistencyConstraint::evaluateSensitivities(), LSQTerm::evaluateSensitivitiesGN(), LSQEndTerm::evaluateSensitivitiesGN(), ShootingMethod::evaluateSensitivitiesLifted(), CondensingBasedCPsolver::expand(), ConstraintElement::get(), getCols(), getEigenvalues(), ParameterEstimationAlgorithm::getParameterVarianceCovariance(), Constraint::getPointConstraint(), getRows(), getSingularValueDecomposition(), BlockMatrix::getSubBlock(), QPsolver_qpOASES::getVarianceCovariance(), DenseCP::init(), MatrixVariable::init(), init(), operator=(), IntegratorBDF::rk_start_solve(), IntegratorBDF::setBackwardSeed2(), ExportNLPSolver::setConstraints(), IntegratorBDF::setForwardSeed2(), LSQTerm::setGrid(), BlockMatrix::setIdentity(), IntegratorBDF::setProtectedBackwardSeed(), IntegratorBDF::setProtectedForwardSeed(), CondensingBasedCPsolver::setupRelaxedQPdataL1(), CondensingBasedCPsolver::setupRelaxedQPdataL2(), MultiObjectiveAlgorithm::solve(), and MultiObjectiveAlgorithm::solveSingleObjective().

returnValue Matrix::init ( uint  _nRows,
uint  _nCols,
double *  _values 
)



Parameters:
_nRowsNumber of rows.
_nColsNumber of columns.
_nColsDouble array.

References init(), nCols, nRows, solver, and SUCCESSFUL_RETURN.

Returns:
BT_TRUE iff matrix object is diagonal.

References acadoIsZero(), BT_FALSE, BT_TRUE, getNumCols(), getNumRows(), and isSquare().

Referenced by ExportGaussNewtonForces::setupObjectiveEvaluation().

BooleanType VectorspaceElement::isEmpty ( ) const [inline, inherited]
Returns:
BT_TRUE iff vectorspace element has dimension zero, BT_FALSE otherwise.

References BT_FALSE, BT_TRUE, and VectorspaceElement::getDim().

Referenced by appendCols(), appendRows(), VariableSettings::appendSettings(), CondensingBasedCPsolver::areRealTimeParametersDefined(), SCPmethod::checkForRealTimeMode(), ShootingMethod::differentiateForward(), ShootingMethod::differentiateForwardBackward(), IntegratorExport::equidistantControlGrid(), IntegratorBDF::evaluate(), IntegratorRK::evaluate(), IntegratorLYAPUNOV::evaluate(), RealTimeAlgorithm::feedbackStep(), VariableSettings::getLowerBound(), VariableSettings::getLowerBounds(), ModelData::getNX(), VariableSettings::getScaling(), VariableSettings::getUpperBound(), VariableSettings::getUpperBounds(), ModelData::hasEquidistantControlGrid(), OCP::hasEquidistantGrid(), VariableSettings::hasLowerBounds(), VariableSettings::hasScaling(), VariableSettings::hasUpperBounds(), RealTimeAlgorithm::init(), Actuator::init(), RungeKuttaExport::initializeButcherTableau(), DenseCP::isLP(), IntegratorBDF::logCurrentIntegratorStep(), VariableSettings::operator=(), IntegratorExport::setLinearInput(), IntegratorExport::setLinearOutput(), NARXExport::setLinearOutput(), VariableSettings::setLowerBound(), IntegratorExport::setModelData(), VariableSettings::setScaling(), VariableSettings::setUpperBound(), SCPmethod::setupRealTimeParameters(), VariablesGrid::shiftBackwards(), MatrixVariablesGrid::shiftBackwards(), CondensingBasedCPsolver::solve(), and VariableSettings::VariableSettings().

BooleanType VectorspaceElement::isEqualTo ( double  _value) const [inline, inherited]
Returns:
BT_TRUE iff all components are equal to _value, BT_FALSE otherwise.

References acadoIsEqual(), BT_FALSE, BT_TRUE, VectorspaceElement::getDim(), and uint.

Referenced by VectorspaceElement::isZero().

BooleanType VectorspaceElement::isFinite ( ) const [inline, inherited]
BooleanType VectorspaceElement::isGreaterThan ( double  _value) const [inline, inherited]
BooleanType VectorspaceElement::isNegative ( ) const [inline, inherited]
BooleanType VectorspaceElement::isPositive ( ) const [inline, inherited]

Note that this test involves a Cholesky decomposition and thus is computationally expensive for larger matrix dimensions.

Returns:
BT_TRUE iff matrix object is (strictly) positive definite.

References acadoIsPositive(), BT_FALSE, getEigenvalues(), VectorspaceElement::getMin(), and isSquare().

Referenced by ExportNLPSolver::setObjective().

Note that this test involves a Cholesky decomposition and thus is computationally expensive for larger matrix dimensions.

Returns:
BT_TRUE iff matrix object is positive definite.

References BT_FALSE, BT_TRUE, getEigenvalues(), VectorspaceElement::getMin(), and isSquare().

Referenced by OCP::minimizeLSQ(), and OCP::minimizeLSQEndTerm().

BooleanType VectorspaceElement::isSmallerThan ( double  _value) const [inline, inherited]
BooleanType Matrix::isSquare ( ) const [inline]
Returns:
BT_TRUE iff matrix object is symmetric.

References acadoIsEqual(), BT_FALSE, BT_TRUE, getNumRows(), isSquare(), and uint.

Referenced by CondensingBasedCPsolver::solveCPsubproblem().

BooleanType VectorspaceElement::isZero ( ) const [inline, inherited]
void Matrix::lowerTriangularInverse ( double *  L,
int  n 
) const [protected]

References ASSERT, and EPS.

Referenced by getCholeskyInverse().

Matrix & Matrix::makeVector ( ) [inline]
Matrix Matrix::minus ( ) [inline]
Matrix Matrix::negative ( ) [inline]
Matrix Matrix::negativeTranspose ( ) const [inline]
BooleanType VectorspaceElement::operator!= ( const VectorspaceElement arg) const [inline, inherited]
Returns:
BT_TRUE iff both objects differ in at least one component.
Parameters:
argObject of comparison.

References ASSERT, BT_FALSE, BT_TRUE, and VectorspaceElement::getDim().

BEGIN_NAMESPACE_ACADO double & VectorspaceElement::operator() ( uint  idx) [inline, inherited]
Returns:
Value of component <idx>
Parameters:
idxIndex of the component to be returned.

References ASSERT, VectorspaceElement::element, and VectorspaceElement::getDim().

Referenced by VectorspaceElement::hasEqualComponents(), Vector::operator%(), Vector::operator*(), and Vector::setUnitVector().

double VectorspaceElement::operator() ( uint  idx) const [inline, inherited]
Returns:
Value of component <idx>
Parameters:
idxIndex of the component to be returned.

References ASSERT, VectorspaceElement::element, and VectorspaceElement::getDim().

BEGIN_NAMESPACE_ACADO double & Matrix::operator() ( uint  rowIdx,
uint  colIdx 
) [inline]
Returns:
Value of component (<rowIdx>,<colIdx>)
Parameters:
rowIdxRow index of the component to be returned.
colIdxColumn index of the component to be returned.

References ASSERT, VectorspaceElement::element, getNumCols(), and getNumRows().

Referenced by computeQRdecomposition(), computeSparseLUdecomposition(), getCol(), getCols(), getDiag(), getRow(), getRows(), Matrix(), minus(), negativeTranspose(), operator*(), operator^(), solveQR(), solveTransposeQR(), symmetrize(), and transpose().

double Matrix::operator() ( uint  rowIdx,
uint  colIdx 
) const [inline]
Returns:
Value of component (<rowIdx>,<colIdx>)
Parameters:
rowIdxRow index of the component to be returned.
colIdxColumn index of the component to be returned.

References ASSERT, VectorspaceElement::element, getNumCols(), and getNumRows().

Matrix Matrix::operator* ( const Matrix arg) const [inline]
Returns:
Temporary object containing result of multiplication.
Parameters:
argMatrix factor.

References ASSERT, getNumCols(), getNumRows(), operator()(), VectorspaceElement::setZero(), and uint.

Vector Matrix::operator* ( const Vector arg) const [inline]
Returns:
Temporary object containing result of multiplication.
Parameters:
argVector factor.

References ASSERT, VectorspaceElement::getDim(), getNumCols(), getNumRows(), operator()(), VectorspaceElement::setZero(), and uint.

Matrix & Matrix::operator*= ( double  scalar) [inline]
Returns:
Reference to object after multiplication.
Parameters:
scalarScalar factor.

References VectorspaceElement::element, fabs(), VectorspaceElement::getDim(), VectorspaceElement::setZero(), uint, and ZERO_EPS.

Matrix Matrix::operator+ ( const Matrix arg) const [inline]
Returns:
Temporary object containing the sum of the matrices.
Parameters:
argSecond summand.

References ASSERT, VectorspaceElement::element, VectorspaceElement::getDim(), getNumCols(), getNumRows(), and uint.

Matrix & Matrix::operator+= ( const Matrix arg) [inline]
Returns:
Reference to object after addition.
Parameters:
argSecond summand.

References ASSERT, VectorspaceElement::element, VectorspaceElement::getDim(), getNumCols(), getNumRows(), and uint.

Matrix Matrix::operator- ( const Matrix arg) const [inline]
Returns:
Temporary object containing the difference of the matrices.
Parameters:
argSubtrahend.

References ASSERT, VectorspaceElement::element, VectorspaceElement::getDim(), getNumCols(), getNumRows(), and uint.

Matrix & Matrix::operator-= ( const Matrix arg) [inline]
Returns:
Reference to object after subtraction.
Parameters:
argSubtrahend.

References ASSERT, VectorspaceElement::element, VectorspaceElement::getDim(), getNumCols(), getNumRows(), and uint.

Matrix & Matrix::operator/= ( double  scalar) [inline]
Returns:
Reference to object after division.
Parameters:
scalarScalar divisor. If it is zero, nothing is done.

References VectorspaceElement::element, fabs(), VectorspaceElement::getDim(), uint, and ZERO_EPS.

BooleanType VectorspaceElement::operator< ( const VectorspaceElement arg) const [inline, inherited]
Returns:
BT_TRUE iff left object is element-wise smaller than the right one.
Parameters:
argObject of comparison.

References ASSERT, BT_FALSE, BT_TRUE, VectorspaceElement::element, EPS, VectorspaceElement::getDim(), and uint.

VectorspaceElement & VectorspaceElement::operator<< ( double *  rhs) [inherited]

References VectorspaceElement::dim, and uint.

BooleanType VectorspaceElement::operator<= ( const VectorspaceElement arg) const [inline, inherited]
Returns:
BT_TRUE iff left object is element-wise smaller/equal than the right one.
Parameters:
argObject of comparison.

References ASSERT, BT_FALSE, BT_TRUE, VectorspaceElement::element, EPS, VectorspaceElement::getDim(), and uint.

Matrix & Matrix::operator= ( const Matrix rhs)
Parameters:
rhsThe right-hand side object.

Reimplemented in MatrixVariable.

References SparseSolver::clone(), nCols, nRows, and solver.

Referenced by appendCols(), appendRows(), and Matrix().

Matrix & Matrix::operator= ( FILE *  rhs)
Parameters:
rhsA file containing the matrix data.

References ACADOINFO, allocateDoublePointerFromFile(), init(), and SUCCESSFUL_RETURN.

BooleanType VectorspaceElement::operator== ( const VectorspaceElement arg) const [inline, inherited]
Returns:
BT_TRUE iff both objects are element-wise equal.
Parameters:
argObject of comparison.

References BT_FALSE, BT_TRUE, VectorspaceElement::element, EPS, fabs(), VectorspaceElement::getDim(), and uint.

BooleanType Matrix::operator== ( const Matrix arg) const [inline]
Parameters:
[in]rhsObject of comparison.
Returns:
BT_TRUE iff both objects are equal,
BT_FALSE otherwise

References acadoIsEqual(), BT_FALSE, BT_TRUE, VectorspaceElement::element, VectorspaceElement::getDim(), getNumCols(), getNumRows(), and uint.

BooleanType VectorspaceElement::operator> ( const VectorspaceElement arg) const [inline, inherited]
Returns:
BT_TRUE iff left object is element-wise greater than the right one.
Parameters:
argObject of comparison.

References ASSERT, BT_FALSE, BT_TRUE, VectorspaceElement::element, EPS, VectorspaceElement::getDim(), and uint.

BooleanType VectorspaceElement::operator>= ( const VectorspaceElement arg) const [inline, inherited]
Returns:
BT_TRUE iff left object is element-wise greater/equal than the right one.
Parameters:
argObject of comparison.

References ASSERT, BT_FALSE, BT_TRUE, VectorspaceElement::element, EPS, VectorspaceElement::getDim(), and uint.

Matrix Matrix::operator^ ( const Matrix arg) const [inline]
Returns:
Temporary object containing result of multiplication.
Parameters:
argMatrix factor.

References ASSERT, getNumCols(), getNumRows(), operator()(), VectorspaceElement::setZero(), and uint.

Vector Matrix::operator^ ( const Vector arg) const [inline]
Returns:
Temporary object containing result of multiplication.
Parameters:
argVector factor.

References ASSERT, VectorspaceElement::getDim(), getNumCols(), getNumRows(), operator()(), VectorspaceElement::setZero(), and uint.

Matrix & Matrix::operator^= ( const double *  rhs)

References nCols, nRows, and uint.

Matrix Matrix::positive ( ) [inline]
returnValue VectorspaceElement::print ( const char *const  name = DEFAULT_LABEL,
const char *const  startString = DEFAULT_START_STRING,
const char *const  endString = DEFAULT_END_STRING,
uint  width = DEFAULT_WIDTH,
uint  precision = DEFAULT_PRECISION,
const char *const  colSeparator = DEFAULT_COL_SEPARATOR,
const char *const  rowSeparator = DEFAULT_ROW_SEPARATOR 
) const [inherited]

Various settings can be specified defining its output format.

Parameters:
[in]nameName label to be printed before the numerical values.
[in]startStringPrefix before printing the numerical values.
[in]endStringSuffix after printing the numerical values.
[in]widthTotal number of digits per single numerical value.
[in]precisionNumber of decimals per single numerical value.
[in]colSeparatorSeparator between the columns of the numerical values.
[in]rowSeparatorSeparator between the rows of the numerical values.
Returns:
SUCCESSFUL_RETURN,
RET_UNKNOWN_BUG

References acadoPrintf(), VectorspaceElement::printToString(), and SUCCESSFUL_RETURN.

Referenced by LSQTerm::evaluate(), SCPmethod::feedbackStep(), Controller::feedbackStep(), Process::init(), Controller::preparationStep(), DenseCP::print(), BlockMatrix::print(), printEigenvalues(), DenseCP::printSolution(), and SimulationEnvironment::step().

returnValue VectorspaceElement::print ( const char *const  name,
PrintScheme  printScheme 
) const [inherited]

Various settings can be specified defining its output format.

Parameters:
[in]nameName label to be printed before the numerical values.
[in]printSchemePrint scheme defining the output format of the information.
Returns:
SUCCESSFUL_RETURN,
RET_UNKNOWN_BUG

References acadoPrintf(), VectorspaceElement::printToString(), and SUCCESSFUL_RETURN.



Returns:
SUCCESSFUL_RETURN.

References getEigenvalues(), and VectorspaceElement::print().

returnValue Matrix::printToFile ( const char *const  filename,
const char *const  name = DEFAULT_LABEL,
const char *const  startString = DEFAULT_START_STRING,
const char *const  endString = DEFAULT_END_STRING,
uint  width = DEFAULT_WIDTH,
uint  precision = DEFAULT_PRECISION,
const char *const  colSeparator = DEFAULT_COL_SEPARATOR,
const char *const  rowSeparator = DEFAULT_ROW_SEPARATOR 
) const [virtual]

Various settings can be specified defining its output format.

Parameters:
[in]filenameFilename for printing.
[in]nameName label to be printed before the numerical values.
[in]startStringPrefix before printing the numerical values.
[in]endStringSuffix after printing the numerical values.
[in]widthTotal number of digits per single numerical value.
[in]precisionNumber of decimals per single numerical value.
[in]colSeparatorSeparator between the columns of the numerical values.
[in]rowSeparatorSeparator between the rows of the numerical values.
Returns:
SUCCESSFUL_RETURN,
RET_FILE_CAN_NOT_BE_OPENED,
RET_UNKNOWN_BUG

Reimplemented from VectorspaceElement.

Referenced by operator<<(), DenseCP::printToFile(), and printToFile().

returnValue Matrix::printToFile ( FILE *  file,
const char *const  name = DEFAULT_LABEL,
const char *const  startString = DEFAULT_START_STRING,
const char *const  endString = DEFAULT_END_STRING,
uint  width = DEFAULT_WIDTH,
uint  precision = DEFAULT_PRECISION,
const char *const  colSeparator = DEFAULT_COL_SEPARATOR,
const char *const  rowSeparator = DEFAULT_ROW_SEPARATOR 
) const [virtual]

Various settings can be specified defining its output format.

Parameters:
[in]fileFile for printing.
[in]nameName label to be printed before the numerical values.
[in]startStringPrefix before printing the numerical values.
[in]endStringSuffix after printing the numerical values.
[in]widthTotal number of digits per single numerical value.
[in]precisionNumber of decimals per single numerical value.
[in]colSeparatorSeparator between the columns of the numerical values.
[in]rowSeparatorSeparator between the rows of the numerical values.
Returns:
SUCCESSFUL_RETURN,
RET_FILE_CAN_NOT_BE_OPENED,
RET_UNKNOWN_BUG

Reimplemented from VectorspaceElement.

References printToFile().

returnValue Matrix::printToFile ( const char *const  filename,
const char *const  name,
PrintScheme  printScheme 
) const [virtual]

Various settings can be specified defining its output format.

Parameters:
[in]filenameFilename for printing.
[in]nameName label to be printed before the numerical values.
[in]printSchemePrint scheme defining the output format of the information.
Returns:
SUCCESSFUL_RETURN,
RET_FILE_CAN_NOT_BE_OPENED,
RET_UNKNOWN_BUG
Parameters:
filenameFilename for printing

Reimplemented from VectorspaceElement.

References ACADOERROR, MatFile::close(), MatFile::open(), printToFile(), PS_MATLAB_BINARY, RET_FILE_CAN_NOT_BE_OPENED, SUCCESSFUL_RETURN, and MatFile::write().

returnValue Matrix::printToFile ( FILE *  file,
const char *const  name,
PrintScheme  printScheme 
) const [virtual]

Various settings can be specified defining its output format.

Parameters:
[in]filenFile for printing.
[in]nameName label to be printed before the numerical values.
[in]printSchemePrint scheme defining the output format of the information.
Returns:
SUCCESSFUL_RETURN,
RET_FILE_CAN_NOT_BE_OPENED,
RET_UNKNOWN_BUG

Reimplemented from VectorspaceElement.

References printToFile().

returnValue Matrix::printToString ( char **  string,
const char *const  name = DEFAULT_LABEL,
const char *const  startString = DEFAULT_START_STRING,
const char *const  endString = DEFAULT_END_STRING,
uint  width = DEFAULT_WIDTH,
uint  precision = DEFAULT_PRECISION,
const char *const  colSeparator = DEFAULT_COL_SEPARATOR,
const char *const  rowSeparator = DEFAULT_ROW_SEPARATOR,
BooleanType  allocateMemory = BT_TRUE 
) const [virtual]

Various settings can be specified defining its output format.

Parameters:
[in,out]stringFile for printing.
[in]nameName label to be printed before the numerical values.
[in]startStringPrefix before printing the numerical values.
[in]endStringSuffix after printing the numerical values.
[in]widthTotal number of digits per single numerical value.
[in]precisionNumber of decimals per single numerical value.
[in]colSeparatorSeparator between the columns of the numerical values.
[in]rowSeparatorSeparator between the rows of the numerical values.
[in]allocateMemoryFlag indicating whether memory for string shall be allocated.
Note:
If 'allocateMemory' flag is set to BT_TRUE, it is assumed that no memory is allocated to the 'string' pointer. Otherwise, it is assumed that sufficient memory has been allocated, e.g. by using the determineStringLength() member function.
Returns:
SUCCESSFUL_RETURN,
RET_FILE_CAN_NOT_BE_OPENED,
RET_UNKNOWN_BUG

Reimplemented from VectorspaceElement.

References ACADOERROR, BT_TRUE, determineStringLength(), getNumCols(), getNumRows(), getStringLength(), RET_UNKNOWN_BUG, SUCCESSFUL_RETURN, and uint.

Referenced by LogRecordItem::getValueString(), printToString(), and MatrixVariablesGrid::printToString().

returnValue Matrix::printToString ( char **  string,
const char *const  name,
PrintScheme  printScheme,
BooleanType  allocateMemory = BT_TRUE 
) const [virtual]

Various settings can be specified defining its output format.

Parameters:
[in,out]stringFile for printing.
[in]nameName label to be printed before the numerical values.
[in]printSchemePrint scheme defining the output format of the information.
[in]allocateMemoryFlag indicating whether memory for string shall be allocated.
Note:
If 'allocateMemory' flag is set to BT_TRUE, it is assumed that no memory is allocated to the 'string' pointer. Otherwise, it is assumed that sufficient memory has been allocated, e.g. by using the determineStringLength() member function.
Returns:
SUCCESSFUL_RETURN,
RET_FILE_CAN_NOT_BE_OPENED,
RET_UNKNOWN_BUG

Reimplemented from VectorspaceElement.

References printToString().

returnValue VectorspaceElement::setAll ( double  _value) [inline, inherited]
returnValue Matrix::setCol ( uint  idx,
const Vector arg 
) [inline]
Returns:
SUCCESSFUL_RETURN
Parameters:
idxColumn index.
argNew values of the column.

References ASSERT, VectorspaceElement::getDim(), getNumCols(), getNumRows(), SUCCESSFUL_RETURN, and uint.

Referenced by ShootingMethod::differentiateForward(), ShootingMethod::differentiateForwardBackward(), and VariablesGrid::shiftBackwards().

returnValue Matrix::setRow ( uint  idx,
const Vector arg 
) [inline]
returnValue VectorspaceElement::setZero ( ) [inline, inherited]
Returns:
SUCCESSFUL_RETURN

References VectorspaceElement::setAll().

Referenced by IntegratorBDF::allocateMemory(), IntegratorBDF::applyMTranspose(), IntegratorBDF::applyNewtonStep(), CondensingBasedCPsolver::condense(), IntegratorBDF::constructAll(), PIDcontroller::determineControlAction(), LSQTerm::evaluateSensitivities(), AlgebraicConsistencyConstraint::evaluateSensitivities(), IntegratorBDF::evaluateSensitivities(), LSQTerm::evaluateSensitivitiesGN(), LSQEndTerm::evaluateSensitivitiesGN(), ShootingMethod::evaluateSensitivitiesLifted(), CondensingBasedCPsolver::expand(), MultiObjectiveAlgorithm::formulateOCP(), Integrator::getBackwardSensitivities(), getCholeskyDecomposition(), VariablesGrid::getIntegral(), BlockMatrix::getSubBlock(), VariablesGrid::getSum(), WeightGeneration::getWeights(), MultiObjectiveAlgorithm::getWeights(), FeedforwardLaw::init(), PIDcontroller::init(), RealTimeAlgorithm::init(), Process::init(), Integrator::integrate(), Integrator::integrateSensitivities(), OCP::minimizeLSQ(), Vector::operator*(), operator*(), Vector::operator*=(), operator*=(), BlockMatrix::operator-(), operator^(), PIDcontroller::PIDcontroller(), IntegratorBDF::setBackwardSeed2(), IntegratorLYAPUNOV::setBackwardSeed2(), IntegratorRK::setBackwardSeed2(), IntegratorBDF::setForwardSeed2(), IntegratorRK::setForwardSeed2(), IntegratorLYAPUNOV::setForwardSeed2(), setIdentity(), IntegratorBDF::setProtectedBackwardSeed(), IntegratorRK::setProtectedBackwardSeed(), IntegratorLYAPUNOV::setProtectedBackwardSeed(), IntegratorBDF::setProtectedForwardSeed(), IntegratorRK::setProtectedForwardSeed(), IntegratorLYAPUNOV::setProtectedForwardSeed(), Vector::setUnitVector(), CondensingBasedCPsolver::setupRelaxedQPdataL1(), CondensingBasedCPsolver::setupRelaxedQPdataL2(), BlockMatrix::setZero(), MultiObjectiveAlgorithm::solve(), LinearStateFeedback::step(), and PIDcontroller::step().

Vector Matrix::solveQR ( const Vector b) const



Returns:
The solution x.

References ASSERT, VectorspaceElement::getDim(), getNumCols(), and operator()().

Referenced by IntegratorBDF::applyNewtonStep().

Vector Matrix::solveSparseLU ( const Vector b) const



Returns:
The solution x.

References ASSERT, getNumRows(), SparseSolver::getX(), SparseSolver::solve(), and solver.

Referenced by IntegratorBDF::applyNewtonStep().

Vector Matrix::solveTransposeQR ( const Vector b) const



Returns:
The solution x.

References ASSERT, VectorspaceElement::getDim(), getNumCols(), and operator()().

Referenced by IntegratorBDF::applyMTranspose().

void Matrix::sortByDecreasingSingularValues ( int  nrows,
int  ncols,
double *  singular_value,
double *  U,
double *  V 
) const [protected]
Vector Matrix::sumCol ( const Matrix arg) [inline]
Returns:
sum

Example:

 a   |  b
 c   |  d

returns [a+b;c+d]

Vector Matrix::sumRow ( const Matrix arg) [inline]
Returns:
sum Example:
 a   |  b
 c   |  d
returns [a+c|b+d]
Matrix Matrix::transpose ( ) const [inline]

Friends And Related Function Documentation

double* operator<< ( double *  lhs,
VectorspaceElement rhs 
) [friend, inherited]
returnValue operator<< ( FILE *  file,
VectorspaceElement arg 
) [friend, inherited]


Returns:
SUCCESSFUL_RETURN
RET_CAN_NOT_WRITE_INTO_FILE
Parameters:
filethe file to print to
argthe element to print
returnValue operator<< ( FILE *  file,
Matrix arg 
) [friend]


Returns:
SUCCESSFUL_RETURN
RET_CAN_NOT_WRITE_INTO_FILE
Parameters:
filethe file to print to
argthe matrix to print
double* operator^= ( double *  lhs,
Matrix rhs 
) [friend]


(friend version)


Member Data Documentation

unsigned int VectorspaceElement::dim [protected, inherited]
double* VectorspaceElement::element [protected, inherited]
uint Matrix::nCols [protected]
uint Matrix::nRows [protected]

The documentation for this class was generated from the following files:
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Friends Defines