algebra.h File Reference

Introduction

Header of the lineal algebra routines used in the CuikSuite.

This file includes high level routines which require LU or QR decompositions, and similar procedures.

See also
algebra.c

Definition in file algebra.h.

Data Structures

struct  TBroyden
 Broyden method. More...
 

Macros

#define DAMPED_NEWTON   0
 Damped Newton process. More...
 
#define DN_F   0.1
 Damping Newton initial factor. More...
 
#define DN_MF   0.95
 Damping Newton maximum factor. More...
 
#define DN_W   0.25
 Damping Newton update ratio. More...
 

Functions

int MatrixDeterminantSgn (double epsilon, unsigned int n, double *A)
 Sign of the determinant of a matrix. More...
 
double MatrixDeterminant (unsigned int n, double *A)
 Determinant of a matrix. More...
 
unsigned int MatrixExponential (unsigned int m, double *A, double t, double *eA)
 Exponential of a matrix. More...
 
unsigned int MatrixExponentialVector (unsigned int m, double *A, double t, double *v)
 Exponential of a matrix and product by a vector. More...
 
unsigned int FindRank (double epsilon, unsigned int nr, unsigned int nc, double *mT)
 Determines the row-rank of a matrix. More...
 
unsigned int FindKernel (unsigned int nr, unsigned int nc, double *mT, unsigned int dof, boolean check, double epsilon, double **T)
 Computes the kernel of a matrix. More...
 
unsigned int FindKernelAndRank (unsigned int nr, unsigned int nc, double *mT, unsigned int *rank, double epsilon, double **T)
 Computes the kernel and the rank of a matrix. More...
 
unsigned int FindKernelAndIndependentRows (unsigned int nr, unsigned int nc, double *mT, unsigned int dof, double epsilon, boolean *singular, boolean **IR, double **T)
 Computes the kernel of a matrix and determines the independent rows of this matrix. More...
 
unsigned int FindIndependentRows (unsigned int nr, unsigned int nc, double *mT, unsigned int dof, double epsilon, boolean *singular, boolean **IR)
 Computes the independent rows of this matrix. More...
 
boolean InvertMatrix (unsigned int nr, double *A)
 Inverse of a matrix. More...
 
void InitNewton (unsigned int nr, unsigned int nc, TNewton *n)
 Defines a Newton structure. More...
 
double * GetNewtonMatrixBuffer (TNewton *n)
 Buffer to store the Newton matrix. More...
 
double * GetNewtonRHBuffer (TNewton *n)
 Buffer to store the Newton right hand. More...
 
void NewtonSetMatrix (unsigned int i, unsigned int j, double v, TNewton *n)
 Defines the matrix being used in a Newton step. More...
 
void NewtonSetRH (unsigned int i, double v, TNewton *n)
 Defines the vector being used in a Newton step. More...
 
int NewtonStep (double nullSingularValue, double *x, double *dif, TNewton *n)
 One step in a Newton iteration. More...
 
void DeleteNewton (TNewton *n)
 Releases a Newton structure. More...
 
void InitBroyden (unsigned int nr, unsigned int nc, TBroyden *b)
 Defines a Broyden structure. More...
 
void ResetBroyden (TBroyden *b)
 Resets a Broyden structure. More...
 
double * GetBroydenMatrixBuffer (TBroyden *b)
 Buffer to store the Broyden matrix. More...
 
double * GetBroydenRHBuffer (TBroyden *b)
 Buffer to store the Broyden right hand. More...
 
int BroydenStep (double epsilon, unsigned int *tp, double *x, TBroyden *b)
 One step in a Broyden iteration. More...
 
void BroydenUpdateJacobian (double epsilon, TBroyden *b)
 Updates the approximated Jacobian. More...
 
void DeleteBroyden (TBroyden *b)
 Releases a Broyden structure. More...
 
void InitLS (unsigned int nr, unsigned int nc, unsigned int nrh, TLinearSystem *ls)
 Defines a linear system structure. More...
 
double * GetLSMatrixBuffer (TLinearSystem *ls)
 Buffer to store the A matrix. More...
 
double * GetLSRHBuffer (TLinearSystem *ls)
 Buffer to store the linear system right hand (RH). More...
 
double * GetLSSolutionBuffer (TLinearSystem *ls)
 Buffer to store the linear system solution. More...
 
void LSSetMatrix (unsigned int i, unsigned int j, double v, TLinearSystem *ls)
 Defines the matrix being used in a linear system. More...
 
void LSSetRH (unsigned int i, unsigned int j, double v, TLinearSystem *ls)
 Defines the vector being used in a linear system. More...
 
int LSSolve (TLinearSystem *ls)
 Solves the linear sytem. More...
 
int LSSolveCond (double *rcond, TLinearSystem *ls)
 Solves the linear sytem. More...
 
void DeleteLS (TLinearSystem *ls)
 Releases a linear system structure. More...
 

Macro Definition Documentation

◆ DAMPED_NEWTON

#define DAMPED_NEWTON   0

If 1, the Newton correction step is damped: it is multiplied by a factor alpha in [0,1] to reduce overshooting and to improve converged.

The alpha factor is initialized to DN_F and increased towards DN_MF using the DN_W weight at each step.

Only experimental. Not actually used.

Definition at line 37 of file algebra.h.

◆ DN_F

#define DN_F   0.1

Initial value for the constant used to damp the Newton step.

Definition at line 44 of file algebra.h.

◆ DN_MF

#define DN_MF   0.95

Maximum value for the constant used to damp the Newton step.

Definition at line 51 of file algebra.h.

◆ DN_W

#define DN_W   0.25

Factor used to damp the Newton step.

The update rule is alpha=alpha+DM_W*(DM_MF-alpha) with alpha initialized to DN_F.

Definition at line 61 of file algebra.h.

Function Documentation

◆ MatrixDeterminantSgn()

int MatrixDeterminantSgn ( double  epsilon,
unsigned int  n,
double *  A 
)

Sign of the determinant of a squared matrix.

IMPORTANT: The intput matrix may be modified inside this function. The caller CAN NOT assume that it remains unchanged.

Parameters
epsilonNumerical accuracy.
nNumber of rows/columsn of the matrix.
AThe matrix.
Returns
The matrix determinant sign (-1,0,1).

Referenced by GetChartDegree().

◆ MatrixDeterminant()

double MatrixDeterminant ( unsigned int  n,
double *  A 
)

Determinant of a squared matrix.

This function must be used with caution (only for small matrices) since the determinant can easily overflow.

IMPORTANT: The intput matrix may be modified inside this function. The caller CAN NOT assume that it remains unchanged.

Parameters
nNumber of rows/columsn of the matrix.
AThe matrix.
Returns
The matrix determinant.

Referenced by CompareTangentSpaces().

◆ MatrixExponential()

unsigned int MatrixExponential ( unsigned int  m,
double *  A,
double  t,
double *  eA 
)

Computes the exponential of a matrix (at a given time step).

For a constant matrix A, the error of this function grows with the time 't'. Use with caution for large t's (i.e., 5, 10, is already large).

Parameters
mNumber of rows/columsn of the matrix.
AThe matrix.
tThe time.
eAThe exponential of the matrix: exp(A*t).
Returns
1 if the process was succesful.

Definition at line 1583 of file algebra.c.

References dgpadm_(), IdentityMatrix(), and NEW.

Referenced by dGt(), and LQRPolicy().

◆ MatrixExponentialVector()

unsigned int MatrixExponentialVector ( unsigned int  m,
double *  A,
double  t,
double *  v 
)

Computes the product of the exponential of a matrix with a vector at a given time step.

For a constant matrix A, the error of this function grows with the time 't'. Use with caution for large t's because the error in this case scales much faster than when using MatrixExponential. In this case the error can be significative just for t=1. If accuracy is required (over speed) just use MatrixExponential and multiply the result by vector v.

Parameters
mNumber of rows/columsn of the matrix.
AThe matrix.
tThe time.
vThe input vector. Re-used to store the output: exp(A*t)*v.
Returns
1 if the process was succesful.

Definition at line 1660 of file algebra.c.

References dgchbv_(), and NEW.

◆ FindRank()

unsigned int FindRank ( double  epsilon,
unsigned int  nr,
unsigned int  nc,
double *  mT 
)

Determines the rank of a matrix, i.e. the dimension of the space spanned by the rows/column of the matrix.

For a given problem, the number of variables minus the rank of the Jacobian gives the dimensionality of the solution space, assuming that the Jacobian is evaluated in a regular point. The dimensionality of the solution space is the same as that of its tangent space.

IMPORTANT: The intput matrix may be modified inside this function. The caller CAN NOT assume that it remains unchanged.

Parameters
epsilonNumerical accuracy.
nrNumber of rows of the matrix (no transposed).
ncNumber of columns of the matrix (no transposed).
mTThe TRANSPOSED matrix. This matrix is overwriten inside this function!!
Returns
The rank of the matrix. NO_UINT if the rank can not be computed.

Definition at line 1724 of file algebra.c.

References AnalyzeKernel(), FALSE, and TRUE.

Referenced by ManifoldDimension().

◆ FindKernel()

unsigned int FindKernel ( unsigned int  nr,
unsigned int  nc,
double *  mT,
unsigned int  dof,
boolean  check,
double  epsilon,
double **  T 
)

Defines a basis of the null space of a matrix.

IMPORTANT: The intput matrix may be modified inside this function. The caller CAN NOT assume that it remains unchanged.

Parameters
nrNumber of rows of the matrix (no transposed).
ncNumbe of columns of the matrix (no transposed).
mTThe TRANSPOSED matrix stored as a vector. This matrix is overwriten inside this function!!
dofExpected dimension of the kernel.
checkIf TRUE the function introduces some consistancy checks (whether the kernel dimensionality is larger or smaller than the expected one).
epsilonValues below epsilon are taken as zero.
TThe output kernel. This is a (nc x dof) matrix (stored as a vector). The space for this matrix is allocated in this function but must be de-allocated externally.
Returns
0 if all the operations are correct, 1 if there the kernel is too large (i.e., the given point is singular), 2 if the chart could not be defined since the kernel is too small at the given point, and 3 if there is an error while performing QR decomposition. These outputs come directly from AnalyzeKernel.

Definition at line 1738 of file algebra.c.

References AnalyzeKernel(), FALSE, and TRUE.

Referenced by FindRightNullVector(), and RefineSingularPoint().

◆ FindKernelAndRank()

unsigned int FindKernelAndRank ( unsigned int  nr,
unsigned int  nc,
double *  mT,
unsigned int *  rank,
double  epsilon,
double **  T 
)

Defines a basis of the null space of a matrix, without any clue about the dimension of this null space. The returne rank is actually the number of variables minus the dimension of such space.

IMPORTANT: The intput matrix may be modified inside this function. The caller CAN NOT assume that it remains unchanged.

Parameters
nrNumber of rows of the matrix (no transposed).
ncNumbe of columns of the matrix (no transposed).
mTThe TRANSPOSED matrix stored as a vector. This matrix is overwriten inside this function!!
rankThe rank of the input matrix: number of variables (columns) minus dimension of the output null space (num columns of the matrix representing the matrix).
epsilonValues below epsilon are taken as zero.
TThe output kernel. This is a (nc x dof) matrix (stored as a vector). The space for this matrix is allocated in this function but must be de-allocated externally.
Returns
0 if all the operations are correct and 3 if there is an error while performing QR decomposition. These outputs come directly from AnalyzeKernel.

Definition at line 1752 of file algebra.c.

References AnalyzeKernel(), FALSE, and TRUE.

◆ FindKernelAndIndependentRows()

unsigned int FindKernelAndIndependentRows ( unsigned int  nr,
unsigned int  nc,
double *  mT,
unsigned int  dof,
double  epsilon,
boolean singular,
boolean **  IR,
double **  T 
)

Defines a basis of the null space of a matrix and determines a subset of the rows of the matrix that are independent.

This is useful because in our case most (all?) the matrices have redundancy (i.e., rows that are linearly dependent on other rows). However, for some purposes we need to determine a subset of the rows that are linearly independent.

IMPORTANT: The intput matrix may be modified inside this function. The caller CAN NOT assume that it remains unchanged.

Parameters
nrNumber of rows of the matrix (no transposed).
ncNumbe of columns of the matrix (no transposed).
mTThe TRANSPOSED matrix stored as a vector. This matrix is overwriten inside this function!!
dofExpected dimension of the kernel. If zero, the function tries to determine the rank automatically.
epsilonValues below epsilon are taken as zero.
singularTRUE if the matrix is singular (has more null eigen values than the expected ones). Output.
IRThe set of independent rows as a boolean vector with as many entriees as rows in the input matrix and TRUE for the independent rows. The space for this vector is allocated here but must be deallocated externally. If the matrix is singular this contains the most likely basis of the matrix (up to the numerical accuracy). Caution must be taken to use this output in this case.
TThe output kernel. This is a (nc x dof) matrix (stored as a vector). The space for this matrix is allocated in this function but must be de-allocated externally.
Returns
0 if all the operations are correct, 1 if there the kernel is too large (i.e., the given point is singular), 2 if the chart could not be defined since the kernel is too small at the given point, and 3 if there is an error while performing QR decomposition. These outputs come directly from AnalyzeKernel.

Definition at line 1765 of file algebra.c.

References AnalyzeKernel(), FALSE, and TRUE.

Referenced by ComputeJacobianKernelBasis().

◆ FindIndependentRows()

unsigned int FindIndependentRows ( unsigned int  nr,
unsigned int  nc,
double *  mT,
unsigned int  dof,
double  epsilon,
boolean singular,
boolean **  IR 
)

Determines a subset of the rows of a matrix that are independent.

IMPORTANT: The intput matrix may be modified inside this function. The caller CAN NOT assume that it remains unchanged.

Parameters
nrNumber of rows of the matrix (no transposed).
ncNumbe of columns of the matrix (no transposed).
mTThe TRANSPOSED matrix stored as a vector. This matrix is overwriten inside this function!!
dofExpected dimension of the kernel. If zero, the function tries to determine the rank automatically.
epsilonValues below epsilon are taken as zero.
singularTRUE if the matrix is singular (has more null eigen values than the expected ones). Output.
IRThe set of independent rows as a boolean vector with as many entriees as rows in the input matrix and TRUE for the independent rows. The space for this vector is allocated here but must be deallocated externally. If the matrix is singular this contains the most likely basis of the matrix (up to the numerical accuracy). Caution must be taken to use this output in this case.
Returns
0 if all the operations are correct, 1 if there the kernel is too large (i.e., the given point is singular), 2 if the chart could not be defined since the kernel is too small at the given point, and 3 if there is an error while performing QR decomposition. These outputs come directly from AnalyzeKernel.

Definition at line 1779 of file algebra.c.

References AnalyzeKernel(), FALSE, and TRUE.

Referenced by GetPositionJacobian().

◆ InvertMatrix()

boolean InvertMatrix ( unsigned int  nr,
double *  A 
)

Computes the inverse of a matrix.

Parameters
nrNumber of rows (and columns) of the matrix.
AA square matrix.
Returns
TRUE if the inverse is computed correctly.

◆ InitNewton()

void InitNewton ( unsigned int  nr,
unsigned int  nc,
TNewton *  n 
)

Initializes the structure to be used in a iterative Newton process.

Note that the dimensionality of the solution set can be deduced at each Newton step (using the QR with pivoting) but it is more efficient if we know it beforehand.

Parameters
nrNumber of rows.
ncNumber of columns.
nThe Newton structure to initialize.

Referenced by CuikNewtonInBox(), CuikNewtonSimp(), Newton2ManifoldPlane(), and RefineSingularPoint().

◆ GetNewtonMatrixBuffer()

double* GetNewtonMatrixBuffer ( TNewton *  n)
inline

Buffer to store the Newton matrix.

This buffer must be accessed using the RC2INDEX macro since the matrix can be stored row major or column major depending on the underlying lineal algebra library being used.

Parameters
nThe Newton structure.
Returns
A pointer to the buffer where to store the matrix.

Definition at line 1794 of file algebra.c.

Referenced by CuikNewtonInBox(), CuikNewtonSimp(), Newton2ManifoldPlane(), and RefineSingularPoint().

◆ GetNewtonRHBuffer()

double* GetNewtonRHBuffer ( TNewton *  n)
inline

Buffer to store the Newton RH.

Parameters
nThe Newton structure.
Returns
A pointer to the buffer where to store the RH.

Definition at line 1799 of file algebra.c.

Referenced by CuikNewtonInBox(), CuikNewtonSimp(), Newton2ManifoldPlane(), and RefineSingularPoint().

◆ NewtonSetMatrix()

void NewtonSetMatrix ( unsigned int  i,
unsigned int  j,
double  v,
TNewton *  n 
)
inline

Sets one element of the matrix to be used in one Newton step. This matrix is typically initilized externally, but here we provide a mehtod to set it.

Parameters
iThe row.
jThe column.
vThe new value.
nThe Newton structure to set.

Definition at line 1804 of file algebra.c.

References RC2INDEX.

Referenced by CuikNewtonInBox().

◆ NewtonSetRH()

void NewtonSetRH ( unsigned int  i,
double  v,
TNewton *  n 
)
inline

Sets one element of the vector to be used in one Newton step. This vector is typically initilized externally, but here we provide a mehtod to set it.

Parameters
iThe index in the vector.
vThe new value.
nThe Newton structure to set.

Definition at line 1809 of file algebra.c.

Referenced by CuikNewtonInBox().

◆ NewtonStep()

int NewtonStep ( double  nullSingularValue,
double *  x,
double *  dif,
TNewton *  n 
)

Computes and applies the correction in one step of a Newton iteration.

This function allows to easily apply Newton processes just defining the input data (A, b), setting the initial value for x, and allocating space for the process. The space for the process is not allocated inside this function to avoid continous malloc/free calls since this function is typically used many times in a row.

IMPORTANT: The intput matrix may be modified inside this function. The caller CAN NOT assume that it remains unchanged.

Parameters
nullSingularValueSingular values below this are set to zero.
xThe current approximation to the system solution. This is updated internally.
difNorm of the change applied to x.
nThe Newton information.
Returns
0 if no error is encountered in the Newton step.

Referenced by CuikNewtonInBox(), CuikNewtonSimp(), Newton2ManifoldPlane(), and RefineSingularPoint().

◆ DeleteNewton()

void DeleteNewton ( TNewton *  n)

Releases the memory allocated for a Newton step.

Parameters
nA Newton structure.

Referenced by CuikNewtonInBox(), CuikNewtonSimp(), Newton2ManifoldPlane(), and RefineSingularPoint().

◆ InitBroyden()

void InitBroyden ( unsigned int  nr,
unsigned int  nc,
TBroyden b 
)

Initializes the structure to be used in a iterative Broyden process.

Note that the dimensionality of the solution set can be deduced at each Broyden step (using the QR with pivoting) but it is more efficient if we know it beforehand.

Parameters
nrNumber of rows.
ncNumber of columns.
bThe Broyden structure to initialize.

Definition at line 1839 of file algebra.c.

References Error(), GetLSMatrixBuffer(), GetLSRHBuffer(), GetLSSolutionBuffer(), InitLS(), and NEW.

Referenced by InitDynamicSpace().

◆ ResetBroyden()

void ResetBroyden ( TBroyden b)
inline

Resets the Broyden structure (the internal iteration) without releasing memory.

Parameters
bThe Broyden structure to initialize.

Definition at line 1862 of file algebra.c.

Referenced by NextDynamicState().

◆ GetBroydenMatrixBuffer()

double* GetBroydenMatrixBuffer ( TBroyden b)
inline

Buffer to store the Broyden matrix.

This buffer must be accessed using the RC2INDEX macro since the matrix can be stored row major or column major depending on the underlying lineal algebra library being used.

Parameters
bThe Broyden structure.
Returns
A pointer to the buffer where to store the matrix.

Definition at line 1867 of file algebra.c.

References Error().

Referenced by InitDynamicSpace().

◆ GetBroydenRHBuffer()

double* GetBroydenRHBuffer ( TBroyden b)
inline

Buffer to store the Broyden RH.

Parameters
bThe Broyden structure.
Returns
A pointer to the buffer where to store the RH.

Definition at line 1874 of file algebra.c.

Referenced by InitDynamicSpace().

◆ BroydenStep()

int BroydenStep ( double  epsilon,
unsigned int *  tp,
double *  x,
TBroyden b 
)

Computes and applies the correction in one step of a Broyden iteration. Calls the Newton iteration with the current Jacobian matrix and updates this matrix.

The result is stored in the right-hand side buffer.

Parameters
epsilonNumerical accuracy.
tpTopology of the variables.
xThe current approximation to the system solution. This is updated internally.
bThe Broyden information.
Returns
0 if no error is encountered in the Broyden step.

Definition at line 1879 of file algebra.c.

References LSSolve(), NEW, Norm(), PI2PI, PrintMatrix(), PrintVector(), and TOPOLOGY_S.

Referenced by NextDynamicState().

◆ BroydenUpdateJacobian()

void BroydenUpdateJacobian ( double  epsilon,
TBroyden b 
)

Applies the Broyden's Jacobian update procedure. This can only be used after using BroydenStep at least once.

Parameters
epsilonThe numerical accuracy.
bThe Broyden structure.

Definition at line 1952 of file algebra.c.

References AccumulateVector(), DifferenceVector(), GeneralDotProduct(), MatrixVectorProduct(), and RankOneUpdate().

Referenced by NextDynamicState().

◆ DeleteBroyden()

void DeleteBroyden ( TBroyden b)

Releases the memory allocated for a Broyden step.

Parameters
bA Broyden structure.

Definition at line 1987 of file algebra.c.

References DeleteLS().

Referenced by DeleteDynamicSpace().

◆ InitLS()

void InitLS ( unsigned int  nr,
unsigned int  nc,
unsigned int  nrh,
TLinearSystem *  ls 
)

Initializes the structure to be used when solving a linear system.

Parameters
nrNumber of rows.
ncNumber of columns.
nrhNumber of right hand terms.
lsThe linear system structure to initialize.

Referenced by Chart2Manifold(), InitBroyden(), InitDynamicSpace(), and SetLinearizedDynamics().

◆ GetLSMatrixBuffer()

double* GetLSMatrixBuffer ( TLinearSystem *  ls)
inline

Buffer to store the A matrix.

This buffer must be accessed using the RC2INDEX macro since the matrix can be stored row major or column major depending on the underlying lineal algebra library being used.

Parameters
lsThe linear system structure.
Returns
A pointer to the buffer where to store the matrix.

Definition at line 1814 of file algebra.c.

Referenced by Chart2Manifold(), InitBroyden(), InitDynamicSpace(), and SetLinearizedDynamics().

◆ GetLSRHBuffer()

double* GetLSRHBuffer ( TLinearSystem *  ls)
inline

Buffer to store the linear system RH.

Parameters
lsThe linear system structure.
Returns
A pointer to the buffer where to store the RH.

Definition at line 1819 of file algebra.c.

Referenced by Chart2Manifold(), InitBroyden(), InitDynamicSpace(), and SetLinearizedDynamics().

◆ GetLSSolutionBuffer()

double* GetLSSolutionBuffer ( TLinearSystem *  ls)
inline

Buffer to store the linear sytem solution. Note that this maybe the same as the buffer for the rhs.

Parameters
lsThe linear system structure.
Returns
A pointer to the buffer where to store the solution.

Definition at line 1824 of file algebra.c.

Referenced by Chart2Manifold(), ComputeAcceleration(), InitBroyden(), InitDynamicSpace(), and SetLinearizedDynamics().

◆ LSSetMatrix()

void LSSetMatrix ( unsigned int  i,
unsigned int  j,
double  v,
TLinearSystem *  ls 
)
inline

Sets one element of the matrix to be used in a linear system.

Parameters
iThe row.
jThe column.
vThe new value.
lsThe linear system structure.

Definition at line 1829 of file algebra.c.

References RC2INDEX.

◆ LSSetRH()

void LSSetRH ( unsigned int  i,
unsigned int  j,
double  v,
TLinearSystem *  ls 
)
inline

Sets one element of the right hand vector of a linear system.

Parameters
iThe row index.
jThe column index.
vThe new value.
lsThe linear system structure.

Definition at line 1834 of file algebra.c.

References RC2INDEX.

◆ LSSolve()

int LSSolve ( TLinearSystem *  ls)

Computes the solution of a linear system.

IMPORTANT: The intput matrix may be modified inside this function. The caller CAN NOT assume that it remains unchanged.

Parameters
lsThe linear system to solve.
Returns
0 if no error is encountered when solving the system.

Referenced by BroydenStep(), Chart2Manifold(), ComputeAcceleration(), ComputeInverseDynamics(), and SetLinearizedDynamics().

◆ LSSolveCond()

int LSSolveCond ( double *  rcond,
TLinearSystem *  ls 
)

Computes the solution of a linear system only and returns the the reciprocal of the condition number of the system matrix. This provides information about the confidence in the solution (solutions with low condition number are not reliable).

Since we compute 1 / ( norm(A) * norm(inv(A)) ), we only implement this option for squared systems (otherwise inv(A) would not be defined).

IMPORTANT: The intput matrix may be modified inside this function. The caller CAN NOT assume that it remains unchanged.

Parameters
rcondThe condition number. Returned mainly for debugging purposes.
lsThe linear system to solve.
Returns
0 if no error is encountered when solving the system.

Referenced by LQRComputePolicy(), and LQRComputePolicy_t().

◆ DeleteLS()

void DeleteLS ( TLinearSystem *  ls)

Releases the memory allocated for a linear system.

Parameters
lsA linear system structure.

Referenced by Chart2Manifold(), DeleteBroyden(), DeleteDynamicSpace(), and SetLinearizedDynamics().