rsLQR  0.1
rsLQR Solver

Data Structures

struct  UnitRange
 Represents a range of consecutive integers. More...
 
struct  BinaryNode_s
 
struct  OrderedBinaryTree
 The binary tree for the rsLQR solver. More...
 
struct  NdLqrCholeskyFactors
 Stores a list of CholeskyInfo structs for the rsLQR solver. More...
 
struct  NdFactor
 A chunk of memory for a single time step. More...
 
struct  NdData
 Core storage container for the rsLQR solver. More...
 
struct  NdLqrProfile
 A struct describing how long each part of the solve took, in milliseconds. More...
 
struct  NdLqrSolver
 Main solver for rsLQR. More...
 

Typedefs

typedef struct BinaryNode_s BinaryNode
 One of the nodes in the binary tree for the rsLQR solver. More...
 

Functions

OrderedBinaryTree ndlqr_BuildTree (int N)
 Construct a new binary tree for a horizon of length N. More...
 
int ndlqr_FreeTree (OrderedBinaryTree *tree)
 Frees the data in tree. More...
 
int ndlqr_GetIndexFromLeaf (const OrderedBinaryTree *tree, int leaf, int level)
 Get the knot point index given the leaf index at a given level. More...
 
int ndlqr_GetIndexLevel (const OrderedBinaryTree *tree, int index)
 Get the level for a given knot point index. More...
 
int ndlqr_GetIndexAtLevel (const OrderedBinaryTree *tree, int index, int level)
 Get the index in level that corresponds to index. More...
 
NdLqrCholeskyFactorsndlqr_NewCholeskyFactors (int depth, int nhorizon)
 Initialize a new NdLqrCholeskyFactors object. More...
 
int ndlqr_FreeCholeskyFactors (NdLqrCholeskyFactors *cholfacts)
 Free the memory of a CholeskyFactors object. More...
 
int ndlqr_GetQFactorizon (NdLqrCholeskyFactors *cholfacts, int index, CholeskyInfo **cholfact)
 Get the CholeskyInfo for the matrix Q at index index. More...
 
int ndlqr_GetRFactorizon (NdLqrCholeskyFactors *cholfacts, int index, CholeskyInfo **cholfact)
 Get the CholeskyInfo for the matrix R at index index. More...
 
int ndlqr_GetSFactorization (NdLqrCholeskyFactors *cholfacts, int leaf, int level, CholeskyInfo **cholfact)
 Get the CholeskyInfo for the output of ndlqr_SolveCholeskyFactor(). More...
 
Matrix ndlqr_GetLambdaFactor (NdFactor *factor)
 
Matrix ndlqr_GetStateFactor (NdFactor *factor)
 
Matrix ndlqr_GetInputFactor (NdFactor *factor)
 
NdDatandlqr_NewNdData (int nstates, int ninputs, int nhorizon, int width)
 Initialize the NdData structure. More...
 
int ndlqr_FreeNdData (NdData *nddata)
 Frees the memory allocated in an NdData structure. More...
 
int ndlqr_GetNdFactor (NdData *nddata, int index, int level, NdFactor **factor)
 Retrieve an individual NdFactor out of the NdData. More...
 
void ndlqr_ResetNdData (NdData *nddata)
 Resets all of the memory for an NdData to zero. More...
 
int ndlqr_SolveLeaf (NdLqrSolver *solver, int index)
 Solve all the equations for the lowest-level diagonal blocks, by timestep. More...
 
int ndlqr_SolveLeaves (NdLqrSolver *solver)
 
int ndlqr_FactorInnerProduct (NdData *data, NdData *fact, int index, int data_level, int fact_level)
 Calculates one of the inner products needed at level data_level. More...
 
int ndlqr_SolveCholeskyFactor (NdData *fact, CholeskyInfo *cholinfo, int index, int level, int upper_level)
 Use the precomputed Cholesky factorization to solve for y at each parent level. More...
 
bool ndlqr_ShouldCalcLambda (OrderedBinaryTree *tree, int index, int i)
 Determines if the \( \Lambda \) should be updated during ndlqr_UpdateSchurFactor() More...
 
int ndlqr_UpdateShurFactor (NdData *fact, NdData *soln, int index, int i, int level, int upper_level, bool calc_lambda)
 Calculates \( x \) and \( z \) to complete the factorization at the current level. More...
 
int ndlqr_ComputeShurCompliment (NdLqrSolver *solver, int index, int level, int upper_level)
 
int ndlqr_Solve (NdLqrSolver *solver)
 Solve an LQR problem using rsLQR. More...
 
Matrix ndlqr_GetSolution (NdLqrSolver *solver)
 Return the solution vector. More...
 
int ndlqr_CopySolution (NdLqrSolver *solver, double *soln)
 Copies the solution vector to a user-supplied array. More...
 
NdLqrProfile ndlqr_NewNdLqrProfile ()
 Create a profile initialized with zeros.
 
void ndlqr_ResetProfile (NdLqrProfile *prof)
 Reset the profile to its initialized state. More...
 
void ndlqr_CopyProfile (NdLqrProfile *dest, NdLqrProfile *src)
 Copy the profile information to a new profile. More...
 
void ndlqr_PrintProfile (NdLqrProfile *profile)
 Print a summary fo the profile. More...
 
void ndlqr_CompareProfile (NdLqrProfile *base, NdLqrProfile *prof)
 Compare two profiles, printing the comparison to stdout. More...
 
NdLqrSolverndlqr_NewNdLqrSolver (int nstates, int ninputs, int nhorizon)
 Create a new solver, allocating all the required memory. More...
 
int ndlqr_FreeNdLqrSolver (NdLqrSolver *solver)
 Deallocates the memory for the solver. More...
 
int ndlqr_InitializeWithLQRProblem (const LQRProblem *lqrprob, NdLqrSolver *solver)
 Initialize the solver with data from an LQR Problem. More...
 
void ndlqr_ResetSolver (NdLqrSolver *solver)
 Resets the rsLQR solver. More...
 
void ndlqr_PrintSolveSummary (NdLqrSolver *solver)
 Prints a summary of the solve. More...
 
int ndlqr_GetNumVars (NdLqrSolver *solver)
 Gets the total number of decision variables for the problem. More...
 
int ndlqr_SetNumThreads (NdLqrSolver *solver, int num_threads)
 Set the number of threads to be used during the solve. More...
 
int ndlqr_GetNumThreads (NdLqrSolver *solver)
 Get the number of threads used during the rsLQR solve. More...
 
int ndlqr_PrintSolveProfile (NdLqrSolver *solver)
 Prints a summary of how long individual components took. More...
 
NdLqrProfile ndlqr_GetProfile (NdLqrSolver *solver)
 Ge the internal profile data from a solve. More...
 

Detailed Description

Typedef Documentation

◆ BinaryNode

typedef struct BinaryNode_s BinaryNode

One of the nodes in the binary tree for the rsLQR solver.

Function Documentation

◆ ndlqr_BuildTree()

OrderedBinaryTree ndlqr_BuildTree ( int  N)

Construct a new binary tree for a horizon of length N.

Must be paired with a corresponding call to ndlqr_FreeTree().

Parameters
Nhorizon length. Must be a power of 2.
Returns
A new binary tree

◆ ndlqr_CompareProfile()

void ndlqr_CompareProfile ( NdLqrProfile base,
NdLqrProfile prof 
)

Compare two profiles, printing the comparison to stdout.

Parameters
baseThe baseline profile
profThe "new" profile

◆ ndlqr_CopyProfile()

void ndlqr_CopyProfile ( NdLqrProfile dest,
NdLqrProfile src 
)

Copy the profile information to a new profile.

Parameters
destNew location for data. Existing data will be overwritten.
srcData to be copied.

◆ ndlqr_CopySolution()

int ndlqr_CopySolution ( NdLqrSolver solver,
double *  soln 
)

Copies the solution vector to a user-supplied array.

See ndlqr_GetSolution() for variable order.

Precondition
ndlqr_Solve() has already been called
Parameters
solverThe solver with the solution.
solnOutput location for the solution vector. Should have length at least equal to the output of ndlqr_GetNumVars().
Returns

◆ ndlqr_FactorInnerProduct()

int ndlqr_FactorInnerProduct ( NdData data,
NdData fact,
int  index,
int  data_level,
int  fact_level 
)

Calculates one of the inner products needed at level data_level.

Calculates the following:

\[ \bar{\Lambda_{k+1}^{(p)}} = \bar{\Lambda_k^{(p)}}S + \langle X_k^{(j)}, \bar{X}_k^{(p)} \rangle + \langle U_k^{(j)}, \bar{U}_k^{(p)} \rangle + \langle X_{k+1}^{(j)}, \bar{X}_{k+1}^{(p)} \rangle + \langle U_{k+1}^{(j)}, \bar{U}_{k+1}^{(p)} \rangle \]

where \( j \) is data_level, \( p \) is fact_level, and \( k \) is index.

Parameters
dataThe data for the original KKT matrix
factThe current data for the factorization
indexKnot point index
data_levelLevel index for the original data, or equivalently the current level being processed by the outer solve
fact_levelLevel index for the factorization data, or equivalently the parent or upper level. fact_level >= data_level.
Returns
0 if successful

◆ ndlqr_FreeCholeskyFactors()

int ndlqr_FreeCholeskyFactors ( NdLqrCholeskyFactors cholfacts)

Free the memory of a CholeskyFactors object.

Parameters
cholfactsAn initialized NdLqrCholeskyFactors object
Postcondition
cholfacts = NULL
Returns
0 if successful

◆ ndlqr_FreeNdData()

int ndlqr_FreeNdData ( NdData nddata)

Frees the memory allocated in an NdData structure.

Parameters
nddataInitialized NdData structure
Returns
0 if successful

◆ ndlqr_FreeNdLqrSolver()

int ndlqr_FreeNdLqrSolver ( NdLqrSolver solver)

Deallocates the memory for the solver.

Parameters
solverAn initialized solver.
Returns
0 if successful.
Postcondition
solver == NULL

◆ ndlqr_FreeTree()

int ndlqr_FreeTree ( OrderedBinaryTree tree)

Frees the data in tree.

Parameters
treeAn initialized binary tree
Returns
0 if successful

◆ ndlqr_GetIndexAtLevel()

int ndlqr_GetIndexAtLevel ( const OrderedBinaryTree tree,
int  index,
int  level 
)

Get the index in level that corresponds to index.

If the level is higher than the level of the given index, it's simply the parent at that level. If it's lower, then it's the index that's closest to the given one, with ties broken by choosing the left (or smaller) of the two.

Parameters
treePrecomputed binary tree
indexStart index of the search. The result will be the index closest to this index.
levelThe level in which the returned index should belong to.
Returns
int The index closest to the provided one, in the given level. -1 if unsucessful.

◆ ndlqr_GetIndexFromLeaf()

int ndlqr_GetIndexFromLeaf ( const OrderedBinaryTree tree,
int  leaf,
int  level 
)

Get the knot point index given the leaf index at a given level.

Parameters
treeAn initialized binary tree for the problem horizon
leafLeaf index
levelLevel of tree from which to get the index
Returns

◆ ndlqr_GetIndexLevel()

int ndlqr_GetIndexLevel ( const OrderedBinaryTree tree,
int  index 
)

Get the level for a given knot point index.

Parameters
treeAn initialized binary tree for the problem horizon
indexKnot point index
Returns
The level for the given index

◆ ndlqr_GetNdFactor()

int ndlqr_GetNdFactor ( NdData nddata,
int  index,
int  level,
NdFactor **  factor 
)

Retrieve an individual NdFactor out of the NdData.

Retrieves a block of memory out of NdData, stored as an NdFactor. Typical usage will look like:

NdData* ndata = ndlqr_NewNdData(nstates, ninputs, nhorizon, nstates);
NdFactor* factor;
int index = 0; // set to desired index
int level = 1; // set to desired level
ndlqr_GetNdFactor(nddata, index, level, &factor);
Parameters
nddataStorage location for the desired block of memory
indexTime step of the factor to extract
levelLevel (or column in the NdData) of the desired factor
factorStorage location for the factor.
Returns
0 if successful

◆ ndlqr_GetNumThreads()

int ndlqr_GetNumThreads ( NdLqrSolver solver)

Get the number of threads used during the rsLQR solve.

Parameters
solverA solver which has already been initialized and solved
Returns
number of OpenMP threads used the by solver

◆ ndlqr_GetNumVars()

int ndlqr_GetNumVars ( NdLqrSolver solver)

Gets the total number of decision variables for the problem.

Parameters
solver

◆ ndlqr_GetProfile()

NdLqrProfile ndlqr_GetProfile ( NdLqrSolver solver)

Ge the internal profile data from a solve.

Parameters
solverA solver which has already been initialized and solved
Returns
A profile object containing timing information about the solve See NdLqrProfile for more info.

◆ ndlqr_GetQFactorizon()

int ndlqr_GetQFactorizon ( NdLqrCholeskyFactors cholfacts,
int  index,
CholeskyInfo **  cholfact 
)

Get the CholeskyInfo for the matrix Q at index index.

Parameters
[in]cholfactsAll the stored info for the Cholesky solves
[in]indexKnot point index
[out]cholfactLocation for the retrieved CholeskyInfo
Returns
0 if successful

◆ ndlqr_GetRFactorizon()

int ndlqr_GetRFactorizon ( NdLqrCholeskyFactors cholfacts,
int  index,
CholeskyInfo **  cholfact 
)

Get the CholeskyInfo for the matrix R at index index.

Parameters
[in]cholfactsAll the stored info for the Cholesky solves
[in]indexKnot point index
[out]cholfactLocation for the retrieved CholeskyInfo
Returns
0 if successful

◆ ndlqr_GetSFactorization()

int ndlqr_GetSFactorization ( NdLqrCholeskyFactors cholfacts,
int  leaf,
int  level,
CholeskyInfo **  cholfact 
)

Get the CholeskyInfo for the output of ndlqr_SolveCholeskyFactor().

Gets the CholeskyInfo for \( \bar{B}_i^{(j)} \).

Parameters
[in]cholfactsAll the stored info for the Cholesky solves
[in]leafThe leaf index for the desired factor
[in]levelThe level of the binary tree for the factor
[out]cholfactLocation for the retrieved CholeskyInfo
Returns
0 if successful

◆ ndlqr_GetSolution()

Matrix ndlqr_GetSolution ( NdLqrSolver solver)

Return the solution vector.

Returns the solution vector as a Matrix object, which is a simple wrapper around a raw pointer which points to the data actually stored by the solver. The user must not free the data, as it is owned by the solver. To get a solution vector owned by the caller, use ndlqr_CopySolution() instead.

Variable Ordering

The variables are ordered in the solution vector as follows:

\[ \begin{bmatrix} \lambda_1^T & x_1^T & u_1^T & \lambda_2^T & \dots & x_{N-1}^T & u_{N-1}^T & \lambda_N^T & x_N^T \end{bmatrix}^T \]

Precondition
ndlqr_Solve() has already been called
Parameters
solver
Returns

◆ ndlqr_InitializeWithLQRProblem()

int ndlqr_InitializeWithLQRProblem ( const LQRProblem lqrprob,
NdLqrSolver solver 
)

Initialize the solver with data from an LQR Problem.

Precondition
Solver has already been initialized via ndlqr_NewNdLqrSolver()
Parameters
lqrprobAn initialized LQR problem with the data to be be solved.
solverAn initialized solver.
Returns
0 if successful

◆ ndlqr_NewCholeskyFactors()

NdLqrCholeskyFactors* ndlqr_NewCholeskyFactors ( int  depth,
int  nhorizon 
)

Initialize a new NdLqrCholeskyFactors object.

Must be paired with a call to ndlqr_FreeCholeskyFactors().

Parameters
depthDepth of the binary tree
nhorizonLength of the time horizon. log2(@depth) = nhorizon.
Returns
An initialized NdLqrCholeskyFactors object.

◆ ndlqr_NewNdData()

NdData* ndlqr_NewNdData ( int  nstates,
int  ninputs,
int  nhorizon,
int  width 
)

Initialize the NdData structure.

Note this allocates a large block of memory. This should be followed by a single call to ndlqr_FreeNdData().

Parameters
nstatesNumber of variables in the state vector
ninputsNumber of control inputs
nhorizonLength of the time horizon
widthWith of each factor. Should be nstates for KKT matrix data, or 1 for the right-hand side vector.
Returns
The initialized NdData structure

◆ ndlqr_NewNdLqrSolver()

NdLqrSolver* ndlqr_NewNdLqrSolver ( int  nstates,
int  ninputs,
int  nhorizon 
)

Create a new solver, allocating all the required memory.

Must be followed by a later call to ndlqr_FreeNdLqrSolver().

Parameters
nstatesNumber of elements in the state vector
ninputsNumber of control inputs
nhorizonLength of the time horizon. Must be a power of 2.
Returns
A pointer to the new solver

◆ ndlqr_PrintProfile()

void ndlqr_PrintProfile ( NdLqrProfile profile)

Print a summary fo the profile.

Parameters
profile

◆ ndlqr_PrintSolveProfile()

int ndlqr_PrintSolveProfile ( NdLqrSolver solver)

Prints a summary of how long individual components took.

Precondition
ndlqr_Solve() has already been called
Parameters
solverA solver which has already been initialized and solved
Returns
0 if successful

◆ ndlqr_PrintSolveSummary()

void ndlqr_PrintSolveSummary ( NdLqrSolver solver)

Prints a summary of the solve.

Prints solve time, the residual norm, and the number of theads.

Precondition
ndlqr_Solve() has already been called
Parameters
solver

◆ ndlqr_ResetNdData()

void ndlqr_ResetNdData ( NdData nddata)

Resets all of the memory for an NdData to zero.

Parameters
nddataInitialized NdData structure

◆ ndlqr_ResetProfile()

void ndlqr_ResetProfile ( NdLqrProfile prof)

Reset the profile to its initialized state.

Parameters
profA profile

◆ ndlqr_ResetSolver()

void ndlqr_ResetSolver ( NdLqrSolver solver)

Resets the rsLQR solver.

Resets all of the data in the solver to how it was when it was first initialized.

Parameters
solver

◆ ndlqr_SetNumThreads()

int ndlqr_SetNumThreads ( NdLqrSolver solver,
int  num_threads 
)

Set the number of threads to be used during the solve.

Does not guarantee that the specified number of threads will be used. To query the actual number of threads used during the solve, use the ndlqr_GetNumThreads() function after the solve.

Parameters
solverrsLQR solver
num_threadsrequested number of threads
Returns
0 if successful

◆ ndlqr_ShouldCalcLambda()

bool ndlqr_ShouldCalcLambda ( OrderedBinaryTree tree,
int  index,
int  i 
)

Determines if the \( \Lambda \) should be updated during ndlqr_UpdateSchurFactor()

Basically checks to see if the \( \Lambda \) for the given index is a $B_i^{(p)}$ for any level greater than or equal to the current level.

Parameters
treeBinary tree with cached information about the structure of the problem
indexKnot point index calculated using ndlqr_GetIndexAtLevel()
iKnot point index being processed
Returns
0 if successful

◆ ndlqr_Solve()

int ndlqr_Solve ( NdLqrSolver solver)

Solve an LQR problem using rsLQR.

This is the main method of the package. Using the data already allocated in the solver, this first calculates the factorization of the matrix data and then solves for solution vector. Note that calling this method twice will not result in the same data since it replaces the original right-hand-side data with the solution vector.

Extracting the solution

The solution can be retrieved after the solve by calling either ndlqr_GetSolution() or ndlqr_CopySolution().

Calling solve multiple times

If you want to call this method multiple times on the same data (e.g. when benchmarking solve times), use ndlqr_ResetNdData() between solves to reset the factorization data. Then re-initialize the solver with the problem data via ndlqr_InitializeWithLQRProblem().

Parameters
solverAn nsLQR solver that has been initialized with the desired problem data.
Returns
0 if successful.

◆ ndlqr_SolveCholeskyFactor()

int ndlqr_SolveCholeskyFactor ( NdData fact,
CholeskyInfo cholinfo,
int  index,
int  level,
int  upper_level 
)

Use the precomputed Cholesky factorization to solve for y at each parent level.

Solve the following linear system of equations, overwriting the right-hand-side.

\[ \bar{\Lambda}_{k+1}^{(j)} y = \bar{\Lambda}_{k+1}^{(p)} \]

where \( j \) is level, \( p \) is upper_level, and \( k \) is index.

This is the same as solving

\[ -\bar{B}_i^{(j)} y_i^{(j,p)} = \bar{b}_i^{(j,p)} \]

using the notation from the paper.

Parameters
factData for the factorization
cholinfoCached Cholesky factorization for \( \bar{\Lambda}_{k+1}^{(j)} \).
indexKnot point index. Should be calculated using ndlqr_GetIndexFromLeaf().
levelLevel index for the level currently being processed by the upper-level solve.
upper_levelLevel index for the right-hand-side. upper_level > level.
Returns
0 if successful

◆ ndlqr_SolveLeaf()

int ndlqr_SolveLeaf ( NdLqrSolver solver,
int  index 
)

Solve all the equations for the lowest-level diagonal blocks, by timestep.

Called at the beginning of the solve, this method calculates the initial pieces of the factorization needed to start the recursion up through the rest of the levels. For and LQR problem, it calculates

\[ Q_k^{-1} A_k^T Q_k^{-1} q_k R_k^{-1} B_k^T R_k^{-1} r_k \]

for each index \( k \), with special-casing applied to the first and last indices.

Parameters
solverAn initialized rsLQR solver
indexKnotpoint index
Returns
0 if successful

◆ ndlqr_UpdateShurFactor()

int ndlqr_UpdateShurFactor ( NdData fact,
NdData soln,
int  index,
int  i,
int  level,
int  upper_level,
bool  calc_lambda 
)

Calculates \( x \) and \( z \) to complete the factorization at the current level.

Calculates the following

\[ \bar{\Lambda}_k^{(p)} = \bar{\Lambda}_k^{(p)} - \bar{\Lambda}_k^{(j)} \bar{\Lambda}_{k_\text{sep}}^{(p)} \bar{X}_k^{(p)} = \bar{X}_k^{(p)} - \bar{X}_k^{(j)} \bar{\Lambda}_{k_\text{sep}}^{(p)} \bar{U}_k^{(p)} = \bar{U}_k^{(p)} - \bar{U}_k^{(j)} \bar{\Lambda}_{k_\text{sep}}^{(p)} \]

where \( \bar{\Lambda}_{k_\text{sep}}^{(p)} \) is equivalent to $y_i^{(j,p)}$ from the paper and is the result of the ndlqr_SolveCholeskyFactor() for index, level, and upper_level.

Parameters
factNdData for the factorization data
solnNdData for the factorization data or solution vector
indexKnot point index of the "separator" at level level. Should be calculated using ndlqr_GetIndexAtlevel().
iKnot point index to be processed
levelLevel index for the level currently being processed
upper_levelLevel index of the upper level. upper_level > level
calc_lambdaWhether the \( \Lambda \) factor should be updated. This should be calculated using ndlqr_ShouldCalcLambda().
Returns
ndlqr_GetNdFactor
int ndlqr_GetNdFactor(NdData *nddata, int index, int level, NdFactor **factor)
Retrieve an individual NdFactor out of the NdData.
Definition: nddata.c:82
ndlqr_NewNdData
NdData * ndlqr_NewNdData(int nstates, int ninputs, int nhorizon, int width)
Initialize the NdData structure.
Definition: nddata.c:15
NdData
Core storage container for the rsLQR solver.
Definition: nddata.h:83
NdFactor
A chunk of memory for a single time step.
Definition: nddata.h:40