qpSWIFT
A Sparse Quadratic Programming Solver
Auxilary.c File Reference
#include "Auxilary.h"
Include dependency graph for Auxilary.c:

Functions

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 *kkt)
 Assembles the KKT matrix. 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 updatevariables (qp_real *x, qp_real *delta_x, qp_real alpha, qp_int count)
 Performs scalar vector addition. 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...
 
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 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_int checksign (qp_real *x, qp_real *delta_x, qp_real alpha, qp_int count)
 Checks if x + alpha*delta_x < 0. More...
 
qp_real norm (qp_real *p, qp_int n)
 Calculates the Eucledian two norm of vector. More...
 
qp_real innerproduct (qp_real *x, qp_real *y, qp_int n)
 Calculates the inner product of two vectors. 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...
 
qp_int ldlinitialsolve (kkt *mykkt, qp_real *delta)
 Solves the kkt linear system to find initial conditions. More...
 
void formlambda (qp_real *lambda, qp_real *s, qp_real *z, qp_int n)
 Computes the scaling point lambda. 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 SparseMatrixTransMultiply (smat *A, qp_real *x, qp_real *y, qp_int start)
 Performs Sparse Matrix Transpose Vector Multiplication as. More...
 
void SparseMatrixMultiply (smat *A, qp_real *x, qp_real *y, qp_int start)
 Performs Sparse Matrix Vector Multiplication as. 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...
 
void SparseMatrixTranspose (smat *A, smat *At)
 Computes the sparse matrix transpose. More...
 
void findminmax (qp_real *z, long n, qp_real *min, qp_real *max)
 Computes the minimum and maximum of the vector. More...
 
qp_int kkt_initialize (QP *myQP)
 Computes the initial condition for the QP problem. 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...
 
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...
 
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...
 

Function Documentation

◆ checksign()

qp_int checksign ( qp_real *  x,
qp_real *  delta_x,
qp_real  alpha,
qp_int  count 
)

Checks if x + alpha*delta_x < 0.

Parameters
[out]Flag0 : Failure 1 : Success
[in]xprimal variable
[in]delta_xsearch direction
[in]alphascalar value

◆ computeresiduals()

void computeresiduals ( QP myQP)

Computes the residuals rx, ry and rz.

Parameters
[out]myQPQP Structure
rx = Px + G'z +c
ry = Ax - b
rz = -s - Gx + h
mu = s'z/m

◆ densetosparse()

void densetosparse ( qp_int  m,
qp_int  n,
qp_real *  pr,
smat A 
)

Converts dense matrix in column major format to CCS format.

Parameters
[out]ASparse Matrix A
[in]mnumber of rows of A
[in]nnumber of columns of A
[in]prpointer to matrix in column major format

◆ densetosparse_ROWMAJOR()

void densetosparse_ROWMAJOR ( qp_int  m,
qp_int  n,
qp_real *  pr,
smat A 
)

Converts dense matrix in row major format to CCS format.

Parameters
[out]ASparse Matrix A
[in]mnumber of rows of A
[in]nnumber of columns of A
[in]prpointer to matrix in row major format

◆ findminmax()

void findminmax ( qp_real *  z,
long  n,
qp_real *  min,
qp_real *  max 
)

Computes the minimum and maximum of the vector.

Parameters
[out]minMinimum of the vector z
[in]maxMaximum of the vector z
[in]nVector z
[in]nnzLength of the vector z

◆ findsteplength()

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.

Parameters
[out]alpha_pPrimal step length
[out]alpha_dDual step length
[in]sPrimal slack variable
[in]delta_sSearch Direction
[in]zDual slack variable
[in]delta_zSearch Direction
[in]m# of inequality constraints

◆ form_ds()

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.

Parameters
[out]dsvector
[in]lambdaScaling Point
[in]delta_sSearch Direction
[in]delta_zSearch Direction
[in]sigmacentering parameter
[in]mucomplementary condition
[in]m# of inequality constraints
[in]selectorselector
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

◆ formkktmatrix_full()

void formkktmatrix_full ( smat P,
smat G,
smat A,
smat Gt,
smat At,
smat kkt 
)

Assembles the KKT matrix.

Parameters
[out]kktKKT Matrix
[in]PCost Function Matrix
[in]AEquality Constraint Matrix
[in]A'Transpose of A
[in]GInequality Constraint Matrix
[in]G'Transpose of G
KKT = [P    A'   G']
      [A    0    0]
      [G    0   -I]

◆ formkktmatrix_U()

void formkktmatrix_U ( smat P,
smat G,
smat Gt,
smat kkt 
)

Forms the Upper triangular part of the kkt matrix

Status: Inactive

◆ formlambda()

void formlambda ( qp_real *  lambda,
qp_real *  s,
qp_real *  z,
qp_int  n 
)

Computes the scaling point lambda.

Parameters
[out]lambdaScaling Point
[in]sPrimal slack variable
[in]zDual variable
[in]nlength of the vector
lambda = sqrt(s/z)

◆ formrho()

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.

Parameters
[out]rhoscalar value rho
[in]sPrimal slack variable
[in]delta_sdelta_s
[in]zDual Variable
[in]delta_zdelta_z
[in]alpha_pPrimal Step length
[in]alpha_dDual Step length
[in]nlength 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;

◆ innerproduct()

qp_real innerproduct ( qp_real *  x,
qp_real *  y,
qp_int  n 
)

Calculates the inner product of two vectors.

Parameters
[out]valueInner Product of the vectors x and y i.e, <x,y>
[in]xvector x
[in]yvector y
[in]nlength of vectors x and y

◆ kkt_initialize()

qp_int kkt_initialize ( QP myQP)

Computes the initial condition for the QP problem.

Parameters
[in]myQPQP Matrix Structure
Solves a closely related QP of the form
min. 0.5*x'Px + c'x + s's
s.t     Ax = b
        Gx + s = h

◆ kktsolve_1()

qp_int kktsolve_1 ( QP myQP)

Solves the KKT linear systems and updates delta_z and delta_s.

Parameters
[out]FlagFlag indicating status of the function; 0 : Failure 1 : Successful
[in]myQPQP Structure

◆ kktsolve_2()

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.

Parameters
[in]myQPQP Structure

◆ ldlinitialsolve()

qp_int ldlinitialsolve ( kkt mykkt,
qp_real *  delta 
)

Solves the kkt linear system to find initial conditions.

Parameters
[in]mykktKKT Structure
[in]deltadelta 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

◆ norm()

qp_real norm ( qp_real *  p,
qp_int  n 
)

Calculates the Eucledian two norm of vector.

Parameters
[out]valueNorm of the vector
[in]pVector
[in]nLength of the vector

◆ obj_value()

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.

Parameters
[out]fvalObjective Value of QP
[in]Pcost Function : quadratic part
[in]ccost function : linear term
[in]xprimal solution
[in]temptemporary workspace variable

◆ SparseMatrixMultiply()

void SparseMatrixMultiply ( smat A,
qp_real *  x,
qp_real *  y,
qp_int  start 
)

Performs Sparse Matrix Vector Multiplication as.

Parameters
[out]yOutut vector y
[in]ASparse Matrix structure
[in]xNumber of rows of the matrix
[in]startSelection index
Computes y = y - Ax
start = 0 ; do nothing
start !=0 ; intialize y=0

◆ SparseMatrixSetup()

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.

Parameters
[out]sparseSparse Matrix structure
[in]mNumber of rows of the matrix
[in]nNumber of Columns of the matrix
[in]nnzNumber of Non zeros of the matrix
[in]jcVector to store column count ; Dim [n+1]
[in]irVector to store row indices in column major format ; Dim[nnz]
[in]prVector to store matrix values in column major format ; Dim[nnz]

◆ SparseMatrixTransMultiply()

void SparseMatrixTransMultiply ( smat A,
qp_real *  x,
qp_real *  y,
qp_int  start 
)

Performs Sparse Matrix Transpose Vector Multiplication as.

Parameters
[out]yOutut vector y
[in]ASparse Matrix structure
[in]xNumber of rows of the matrix
[in]startSelection index
Computes y = y - A'x
start = 0 ; do nothing
start !=0 ; intialize y=0

◆ SparseMatrixTranspose()

void SparseMatrixTranspose ( smat A,
smat At 
)

Computes the sparse matrix transpose.

Parameters
[out]ASparse Matrix
[in]AtTranspose of the Sparse Matrix

◆ test_reach()

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.

Parameters
[out]UPatternReach 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]ParentParent tree of the KKT Matrix
[in]PinvInverse of Permutation vector

◆ Transpose_Row_Count()

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.

Parameters
[out]Ltiir vector of the transpose of sparse choelsky matrix L (Lt->ir)
[out]Ltpjc vector of the transpose of sparse choelsky matrix L (Lt->jc)
[in]mNumber of rows the matrix L
[in]nNumber of columns of the matrix L
[in]Liir vector of the sparse choelsky matrix L (L->ir)
[in]Lpjc vector of the sparse choelsky matrix L (L->jc)

◆ updatekktmatrix()

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,.

Parameters
[out]kktKKT Matrix
[in]sPrimal Slack Variables
[in]zDual Variables
[in]delta_sPrimal slack direction
[in]delta_zDual slack direction
[in]alpha_pPrimal step Length
[in]alpha_dDual step length
[in]m# of ineqaulity constraints
[in]n# of decision variables
[in]p# of equality constraints
[in]indicatorselection vector
indicator = 1 : Pure Affine Direction
indicator = 2 : Pure Centering Direction
indicator = 3 : Pure Newton Direction

◆ updatekktmatrix_b()

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.

Parameters
[out]bKKT Matrix right-hand side vector
[in]rxResidual
[in]ryResidual
[in]rzResidual
[in]dsvector ds
[in]zResidual
[in]m# of inequality constraints
[in]n# of decision variables
[in]p# of equality constraints
b = [rx]
    [ry]
    [rz] - [ds]/[z]

◆ updatevariables()

void updatevariables ( qp_real *  x,
qp_real *  delta_x,
qp_real  alpha,
qp_int  count 
)

Performs scalar vector addition.

Parameters
[out]xvector
[in]xvector
[in]alphascalar value
[in]delta_xsearch direction
x = x + alpha*delta_x