|
qpSWIFT
A Sparse Quadratic Programming Solver
|
Go to the source code of this file.
Data Structures | |
| struct | smat |
| struct | kkt |
| struct | stats |
| struct | settings |
| struct | QP |
Typedefs | |
| typedef struct smat | smat |
| typedef struct kkt | kkt |
| typedef struct stats | stats |
| typedef struct settings | settings |
| typedef struct QP | QP |
Functions | |
| qp_int | kkt_initialize (QP *myQP) |
| Computes the initial condition for the QP problem. More... | |
| qp_int | kktsolve_1 (QP *myQP) |
| Solves the KKT linear systems and updates delta_z and delta_s. More... | |
| void | kktsolve_2 (QP *myQP) |
| Solves the kktlinear system from results of kktsolve_1 and updates delta_x, delta_y, delta_z and delta_s. More... | |
| void | SparseMatrixSetup (qp_int m, qp_int n, qp_int nnz, qp_int *jc, qp_int *ir, qp_real *pr, smat *sparse) |
| Sets up the Sparse Matrix in Column Compressed Storage Format based on inputs. More... | |
| void | Transpose_Row_Count (qp_int m, qp_int n, qp_int *Li, qp_int *Lp, qp_int *Lti, qp_int *Ltp) |
| Computes the ir and jc of transpose of a Matrix. More... | |
| void | computeresiduals (QP *myQP) |
| Computes the residuals rx, ry and rz. More... | |
| void | SparseMatrixMultiply (smat *A, qp_real *x, qp_real *y, qp_int start) |
| Performs Sparse Matrix Vector Multiplication as. More... | |
| void | SparseMatrixTransMultiply (smat *A, qp_real *x, qp_real *y, qp_int start) |
| Performs Sparse Matrix Transpose Vector Multiplication as. More... | |
| void | form_ds (qp_real *ds, qp_real *lambda, qp_real *delta_s, qp_real *delta_z, qp_real sigma, qp_real mu, qp_int m, qp_int selector) |
| Updates the ds vector based on the selector. More... | |
| void | formkktmatrix_U (smat *P, smat *G, smat *Gt, smat *kkt) |
| void | formkktmatrix_full (smat *P, smat *G, smat *A, smat *Gt, smat *At, smat *kktmatrix) |
| Assembles the KKT matrix. More... | |
| void | SparseMatrixTranspose (smat *A, smat *At) |
| Computes the sparse matrix transpose. More... | |
| void | updatekktmatrix (smat *kkt, qp_real *s, qp_real *z, qp_real *delta_s, qp_real *delta_z, qp_real alpha_p, qp_real alpha_d, qp_int m, qp_int n, qp_int p, qp_int indicator) |
| Updates the lower diagonal part of the kkt Matrix,. More... | |
| void | updatekktmatrix_b (qp_real *b, qp_real *rx, qp_real *ry, qp_real *rz, qp_real *ds, qp_real *z, qp_int n, qp_int m, qp_int p) |
| Updates the right hand side of the KKT linear system of equations. More... | |
| qp_int | checksign (qp_real *s, qp_real *delta_s, qp_real alpha, qp_int count) |
| Checks if x + alpha*delta_x < 0. More... | |
| void | updatevariables (qp_real *x, qp_real *delta_x, qp_real alpha, qp_int count) |
| Performs scalar vector addition. More... | |
| void | formlambda (qp_real *lambda, qp_real *s, qp_real *z, qp_int n) |
| Computes the scaling point lambda. More... | |
| qp_real | formrho (qp_real *s, qp_real *delta_s, qp_real *z, qp_real *delta_z, qp_real alpha_p, qp_real alpha_d, qp_int n) |
| Computes the scalar rho as. More... | |
| qp_int | ldlinitialsolve (kkt *mykkt, qp_real *delta) |
| Solves the kkt linear system to find initial conditions. More... | |
| void | test_reach (qp_int *Parent, qp_int *Pinv, qp_int *UPattern, qp_int n, qp_int m, qp_int p) |
| Computes the nodes to be updated at each iteration. More... | |
| void | findsteplength (qp_real *s, qp_real *delta_s, qp_real *z, qp_real *delta_z, qp_int m, qp_real *alpha_p, qp_real *alpha_d) |
| Calculates the step length. More... | |
| qp_real | obj_value (smat *P, qp_real *c, qp_real *x, qp_real *temp) |
| Computes the objective function value f = x'Px + c'x. More... | |
| qp_real | innerproduct (qp_real *x, qp_real *y, qp_int n) |
| Calculates the inner product of two vectors. More... | |
| void | densetosparse (qp_int m, qp_int n, qp_real *pr, smat *A) |
| Converts dense matrix in column major format to CCS format. More... | |
| void | densetosparse_ROWMAJOR (qp_int m, qp_int n, qp_real *pr, smat *A) |
| Converts dense matrix in row major format to CCS format. More... | |
KKT Structure
Used to store everything related to KKT matrix and its factorization
QP Structure
The main qpSWIFT structure; contains references to all other structures
Settings Structure
Used to store everything related to qpSWIFT settings
Sparse Matrix Storage Format
Use Sparse Matrix Setup to Initialise a Sparse Matrix
Statistics Structure
Used to store everything related to qpSWIFT statistics
| qp_int checksign | ( | qp_real * | x, |
| qp_real * | delta_x, | ||
| qp_real | alpha, | ||
| qp_int | count | ||
| ) |
Checks if x + alpha*delta_x < 0.
| [out] | Flag | 0 : Failure 1 : Success |
| [in] | x | primal variable |
| [in] | delta_x | search direction |
| [in] | alpha | scalar value |
| void computeresiduals | ( | QP * | myQP | ) |
Computes the residuals rx, ry and rz.
| [out] | myQP | QP Structure |
rx = Px + G'z +c ry = Ax - b rz = -s - Gx + h mu = s'z/m
| void densetosparse | ( | qp_int | m, |
| qp_int | n, | ||
| qp_real * | pr, | ||
| smat * | A | ||
| ) |
Converts dense matrix in column major format to CCS format.
| [out] | A | Sparse Matrix A |
| [in] | m | number of rows of A |
| [in] | n | number of columns of A |
| [in] | pr | pointer to matrix in column major format |
| void densetosparse_ROWMAJOR | ( | qp_int | m, |
| qp_int | n, | ||
| qp_real * | pr, | ||
| smat * | A | ||
| ) |
Converts dense matrix in row major format to CCS format.
| [out] | A | Sparse Matrix A |
| [in] | m | number of rows of A |
| [in] | n | number of columns of A |
| [in] | pr | pointer to matrix in row major format |
| void findsteplength | ( | qp_real * | s, |
| qp_real * | delta_s, | ||
| qp_real * | z, | ||
| qp_real * | delta_z, | ||
| qp_int | m, | ||
| qp_real * | alpha_p, | ||
| qp_real * | alpha_d | ||
| ) |
Calculates the step length.
| [out] | alpha_p | Primal step length |
| [out] | alpha_d | Dual step length |
| [in] | s | Primal slack variable |
| [in] | delta_s | Search Direction |
| [in] | z | Dual slack variable |
| [in] | delta_z | Search Direction |
| [in] | m | # of inequality constraints |
| void form_ds | ( | qp_real * | ds, |
| qp_real * | lambda, | ||
| qp_real * | delta_s, | ||
| qp_real * | delta_z, | ||
| qp_real | sigma, | ||
| qp_real | mu, | ||
| qp_int | m, | ||
| qp_int | selector | ||
| ) |
Updates the ds vector based on the selector.
| [out] | ds | vector |
| [in] | lambda | Scaling Point |
| [in] | delta_s | Search Direction |
| [in] | delta_z | Search Direction |
| [in] | sigma | centering parameter |
| [in] | mu | complementary condition |
| [in] | m | # of inequality constraints |
| [in] | selector | selector |
selector = 1 => ds = -lambda*lambda; \n selector = 2 => ds = -lambda*lambda - (delta_s*delta_z) + (sigma*mu); \n selector = 3 => ds = -lambda*lambda + (sigma*mu); \n
Assembles the KKT matrix.
| [out] | kkt | KKT Matrix |
| [in] | P | Cost Function Matrix |
| [in] | A | Equality Constraint Matrix |
| [in] | A' | Transpose of A |
| [in] | G | Inequality Constraint Matrix |
| [in] | G' | Transpose of G |
KKT = [P A' G']
[A 0 0]
[G 0 -I]
Forms the Upper triangular part of the kkt matrix
Status: Inactive
| void formlambda | ( | qp_real * | lambda, |
| qp_real * | s, | ||
| qp_real * | z, | ||
| qp_int | n | ||
| ) |
Computes the scaling point lambda.
| [out] | lambda | Scaling Point |
| [in] | s | Primal slack variable |
| [in] | z | Dual variable |
| [in] | n | length of the vector |
lambda = sqrt(s/z)
| qp_real formrho | ( | qp_real * | s, |
| qp_real * | delta_s, | ||
| qp_real * | z, | ||
| qp_real * | delta_z, | ||
| qp_real | alpha_p, | ||
| qp_real | alpha_d, | ||
| qp_int | n | ||
| ) |
Computes the scalar rho as.
| [out] | rho | scalar value rho |
| [in] | s | Primal slack variable |
| [in] | delta_s | delta_s |
| [in] | z | Dual Variable |
| [in] | delta_z | delta_z |
| [in] | alpha_p | Primal Step length |
| [in] | alpha_d | Dual Step length |
| [in] | n | length of vectors s,z,delta_s and delta_z (# of ineq. constraints) |
Computes the scalar rho as rho = (s+alpha_p*delta_s)*(z+alpha_d*delta_z)/s'z;
| qp_real innerproduct | ( | qp_real * | x, |
| qp_real * | y, | ||
| qp_int | n | ||
| ) |
Calculates the inner product of two vectors.
| [out] | value | Inner Product of the vectors x and y i.e, <x,y> |
| [in] | x | vector x |
| [in] | y | vector y |
| [in] | n | length of vectors x and y |
| qp_int kkt_initialize | ( | QP * | myQP | ) |
| qp_int kktsolve_1 | ( | QP * | myQP | ) |
Solves the KKT linear systems and updates delta_z and delta_s.
| [out] | Flag | Flag indicating status of the function; 0 : Failure 1 : Successful |
| [in] | myQP | QP Structure |
| void kktsolve_2 | ( | QP * | myQP | ) |
Solves the kktlinear system from results of kktsolve_1 and updates delta_x, delta_y, delta_z and delta_s.
| [in] | myQP | QP Structure |
| qp_int ldlinitialsolve | ( | kkt * | mykkt, |
| qp_real * | delta | ||
| ) |
Solves the kkt linear system to find initial conditions.
| [in] | mykkt | KKT Structure |
| [in] | delta | delta varible |
Invoked by kkt_initialize Creates and updates the LDL workspace variables Performs ldl_symbolic and stores the results Also Performs LDL_numeric ; LDL_perm; LDL_lsolve; LDL_dsolve; LDL_ltsolve; LDL_permt in the same order
| qp_real obj_value | ( | smat * | P, |
| qp_real * | c, | ||
| qp_real * | x, | ||
| qp_real * | temp | ||
| ) |
Computes the objective function value f = x'Px + c'x.
| [out] | fval | Objective Value of QP |
| [in] | P | cost Function : quadratic part |
| [in] | c | cost function : linear term |
| [in] | x | primal solution |
| [in] | temp | temporary workspace variable |
| void SparseMatrixMultiply | ( | smat * | A, |
| qp_real * | x, | ||
| qp_real * | y, | ||
| qp_int | start | ||
| ) |
Performs Sparse Matrix Vector Multiplication as.
| [out] | y | Outut vector y |
| [in] | A | Sparse Matrix structure |
| [in] | x | Number of rows of the matrix |
| [in] | start | Selection index |
Computes y = y - Ax start = 0 ; do nothing start !=0 ; intialize y=0
| void SparseMatrixSetup | ( | qp_int | m, |
| qp_int | n, | ||
| qp_int | nnz, | ||
| qp_int * | jc, | ||
| qp_int * | ir, | ||
| qp_real * | pr, | ||
| smat * | sparse | ||
| ) |
Sets up the Sparse Matrix in Column Compressed Storage Format based on inputs.
| [out] | sparse | Sparse Matrix structure |
| [in] | m | Number of rows of the matrix |
| [in] | n | Number of Columns of the matrix |
| [in] | nnz | Number of Non zeros of the matrix |
| [in] | jc | Vector to store column count ; Dim [n+1] |
| [in] | ir | Vector to store row indices in column major format ; Dim[nnz] |
| [in] | pr | Vector to store matrix values in column major format ; Dim[nnz] |
| void SparseMatrixTransMultiply | ( | smat * | A, |
| qp_real * | x, | ||
| qp_real * | y, | ||
| qp_int | start | ||
| ) |
Performs Sparse Matrix Transpose Vector Multiplication as.
| [out] | y | Outut vector y |
| [in] | A | Sparse Matrix structure |
| [in] | x | Number of rows of the matrix |
| [in] | start | Selection index |
Computes y = y - A'x start = 0 ; do nothing start !=0 ; intialize y=0
Computes the sparse matrix transpose.
| [out] | A | Sparse Matrix |
| [in] | At | Transpose of the Sparse Matrix |
| void test_reach | ( | qp_int * | Parent, |
| qp_int * | Pinv, | ||
| qp_int * | UPattern, | ||
| qp_int | n, | ||
| qp_int | m, | ||
| qp_int | p | ||
| ) |
Computes the nodes to be updated at each iteration.
| [out] | UPattern | Reach of the nodes of the lower diagonal part of KKT Matrix |
| [in] | n | # of decision variables |
| [in] | m | # of inequality constraints |
| [in] | p | # of equality constraints |
| [in] | Parent | Parent tree of the KKT Matrix |
| [in] | Pinv | Inverse of Permutation vector |
| void Transpose_Row_Count | ( | qp_int | m, |
| qp_int | n, | ||
| qp_int * | Li, | ||
| qp_int * | Lp, | ||
| qp_int * | Lti, | ||
| qp_int * | Ltp | ||
| ) |
Computes the ir and jc of transpose of a Matrix.
| [out] | Lti | ir vector of the transpose of sparse choelsky matrix L (Lt->ir) |
| [out] | Ltp | jc vector of the transpose of sparse choelsky matrix L (Lt->jc) |
| [in] | m | Number of rows the matrix L |
| [in] | n | Number of columns of the matrix L |
| [in] | Li | ir vector of the sparse choelsky matrix L (L->ir) |
| [in] | Lp | jc vector of the sparse choelsky matrix L (L->jc) |
| void updatekktmatrix | ( | smat * | kkt, |
| qp_real * | s, | ||
| qp_real * | z, | ||
| qp_real * | delta_s, | ||
| qp_real * | delta_z, | ||
| qp_real | alpha_p, | ||
| qp_real | alpha_d, | ||
| qp_int | m, | ||
| qp_int | n, | ||
| qp_int | p, | ||
| qp_int | indicator | ||
| ) |
Updates the lower diagonal part of the kkt Matrix,.
| [out] | kkt | KKT Matrix |
| [in] | s | Primal Slack Variables |
| [in] | z | Dual Variables |
| [in] | delta_s | Primal slack direction |
| [in] | delta_z | Dual slack direction |
| [in] | alpha_p | Primal step Length |
| [in] | alpha_d | Dual step length |
| [in] | m | # of ineqaulity constraints |
| [in] | n | # of decision variables |
| [in] | p | # of equality constraints |
| [in] | indicator | selection vector |
indicator = 1 : Pure Affine Direction indicator = 2 : Pure Centering Direction indicator = 3 : Pure Newton Direction
| void updatekktmatrix_b | ( | qp_real * | b, |
| qp_real * | rx, | ||
| qp_real * | ry, | ||
| qp_real * | rz, | ||
| qp_real * | ds, | ||
| qp_real * | z, | ||
| qp_int | n, | ||
| qp_int | m, | ||
| qp_int | p | ||
| ) |
Updates the right hand side of the KKT linear system of equations.
| [out] | b | KKT Matrix right-hand side vector |
| [in] | rx | Residual |
| [in] | ry | Residual |
| [in] | rz | Residual |
| [in] | ds | vector ds |
| [in] | z | Residual |
| [in] | m | # of inequality constraints |
| [in] | n | # of decision variables |
| [in] | p | # of equality constraints |
b = [rx]
[ry]
[rz] - [ds]/[z]
| void updatevariables | ( | qp_real * | x, |
| qp_real * | delta_x, | ||
| qp_real | alpha, | ||
| qp_int | count | ||
| ) |
Performs scalar vector addition.
| [out] | x | vector |
| [in] | x | vector |
| [in] | alpha | scalar value |
| [in] | delta_x | search direction |
x = x + alpha*delta_x