PIPS-NLP
Public Member Functions | Public Attributes | Protected Attributes | List of all members
NlpGenLinsys Class Referenceabstract

#include <NlpGenLinsys.h>

Inheritance diagram for NlpGenLinsys:
LinearSystem NlpGenSparseLinsys sLinsys sDummyLinsys sLinsysLeaf sLinsysRoot sLinsysRootAggregation sLinsysRootAug

Public Member Functions

 NlpGenLinsys (NlpGen *factory, NlpGenData *data, LinearAlgebraPackage *la)
 
 NlpGenLinsys ()
 
 ~NlpGenLinsys ()
 
virtual void factor (Data *prob, Variables *vars, RegularizationAlg *RegInfo)
 
virtual void factor (Data *prob, Variables *vars)
 
virtual void solve (Data *prob, Variables *vars, Residuals *res, Variables *step)
 
virtual void solve_NTsteps (Data *prob, Variables *vars, Residuals *resids, Variables *Nstep, Variables *Tstep, Variables *NTstep)
 
virtual void solve_IterRefine (Data *prob_in, Variables *vars_in, Residuals *res_in, Variables *step_in, Residuals *KKT_Resid_in, Variables *KKT_sol_in)
 
virtual void joinRHS (OoqpVector &rhs, OoqpVector &rhs1, OoqpVector &rhs2, OoqpVector &rhs3)
 
virtual void separateVars (OoqpVector &vars1, OoqpVector &vars2, OoqpVector &vars3, OoqpVector &vars)
 
virtual void solveXYZS (OoqpVector &stepx, OoqpVector &stepy, OoqpVector &stepz, OoqpVector &steps, OoqpVector &ztemp, NlpGenData *data)
 
virtual void solveCompressed (OoqpVector &rhs)=0
 
virtual void solveCompressedAugXSYZ (OoqpVector &stepx, OoqpVector &steps, OoqpVector &stepy, OoqpVector &stepz, NlpGenData *prob)
 
virtual void solveCompressedAugXSYZ_PETSC (OoqpVector &stepx, OoqpVector &steps, OoqpVector &stepy, OoqpVector &stepz, NlpGenData *prob)
 
virtual void solveBiCGStab (OoqpVector &stepx, OoqpVector &steps, OoqpVector &stepy, OoqpVector &stepz, NlpGenData *data)
 
virtual void putXDiagonal (OoqpVector &xdiag)=0
 
virtual void putSDiagonal (OoqpVector &xdiag)=0
 
virtual void putZDiagonal (OoqpVector &sdiag)=0
 
virtual void computeDiagonals (OoqpVector &dd, OoqpVector &omega, OoqpVector &t, OoqpVector &lambda, OoqpVector &u, OoqpVector &pi, OoqpVector &v, OoqpVector &gamma, OoqpVector &w, OoqpVector &phi)
 
virtual void solveCompressedBiCGStab (OoqpVector &stepx, OoqpVector &stepy, OoqpVector &stepz, NlpGenData *data)
 
virtual void solveCompressedIterRefin (OoqpVector &stepx, OoqpVector &stepy, OoqpVector &stepz, NlpGenData *data)
 
virtual void matXYZMult (double beta, OoqpVector &res, double alpha, OoqpVector &sol, NlpGenData *data, OoqpVector &solx, OoqpVector &soly, OoqpVector &solz)
 
virtual void matXSYZMult (double beta, OoqpVector &res, double alpha, OoqpVector &sol, NlpGenData *data, OoqpVector &solx, OoqpVector &sols, OoqpVector &soly, OoqpVector &solz)
 
virtual void putYDualDiagonal (OoqpVector &ydiag)=0
 
virtual void joinRHSXSYZ (OoqpVector &rhs_in, OoqpVector &rhs1_in, OoqpVector &rhs2_in, OoqpVector &rhs3_in, OoqpVector &rhs4_in)
 
virtual void separateVarsXSYZ (OoqpVector &x_in, OoqpVector &s_in, OoqpVector &y_in, OoqpVector &z_in, OoqpVector &vars_in)
 
virtual void factorNoMatChange (Data *prob_in, Variables *vars_in, RegularizationAlg *RegInfo)
 
virtual void UpdateMatrices (Data *prob_in, int const updateLevel=2)
 
virtual double computeResidual (NlpGenData *data, OoqpVector &res_, OoqpVector &sol_, OoqpVector &solx_, OoqpVector &sols_, OoqpVector &soly_, OoqpVector &solz_)
 
virtual double computeResidual_FullKKT (NlpGenData *data, NlpGenResiduals *res_, NlpGenVars *sol_, NlpGenVars *var_)
 
virtual double eval_xWx (NlpGenData *prob, NlpGenResiduals *resid, NlpGenVars *steps)
 
virtual void computeQuantitiesForDualReg (NlpGenData *prob, NlpGenVars *vars, NlpGenResiduals *resid, NlpGenVars *steps, double *dualRegQuantities)
 
virtual void setXDiagonal (OoqpVector &xdiag)=0
 
virtual void setSDiagonal (OoqpVector &sdiag)=0
 
virtual void setYDiagonal (OoqpVector &ydiag)=0
 
virtual void setZDiagonal (OoqpVector &zdiag)=0
 
virtual void copyXSYZ_fromArray (OoqpVector &vec_xsyz, double *array_in, const int nb_col)
 
virtual void copyXSYZ_toArray (OoqpVector &vec_xsyz, double *array_in, const int nb_col)
 
- Public Member Functions inherited from LinearSystem
virtual ~LinearSystem ()
 

Public Attributes

OoqpVectornomegaInv
 
OoqpVectorrhs
 
OoqpVectorrhs_back
 
NlpGenfactory
 
long long nx
 
long long my
 
long long mz
 
double priReg
 
double dualReg
 
OoqpVectordd
 
OoqpVectordq
 
OoqpVectortemp_diagX
 
OoqpVectortemp_diagS
 
OoqpVectortemp_diagY
 
OoqpVectortemp_diagZ
 
OoqpVectoradditiveDiag
 
OoqpVectorixupp
 
OoqpVectoricupp
 
OoqpVectorixlow
 
OoqpVectoriclow
 
int nxupp
 
int nxlow
 
int mcupp
 
int mclow
 
int num_slacks
 
int useRefs
 
OoqpVectorsol
 
OoqpVectorres
 
OoqpVectorresx
 
OoqpVectorress
 
OoqpVectorresy
 
OoqpVectorresz
 
OoqpVectorsol2
 
OoqpVectorres2
 
OoqpVectorres3
 
OoqpVectorres4
 
OoqpVectorres5
 
OoqpVectorsol2Bicg
 
OoqpVectorres2Bicg
 
OoqpVectorres3Bicg
 
OoqpVectorres4Bicg
 
OoqpVectorres5Bicg
 
int KryIter
 
bool allocateSpace
 

Protected Attributes

bool fullQ
 

Detailed Description

Linear System solvers for the general NLP formulation. This class contains definitions of methods and data common to the sparse and dense special cases of the general formulation. The derived classes NlpGenSparseLinsys and NlpGenDenseLinsys contain the aspects that are specific to the sparse and dense forms.

See also
NlpGenSparseLinsys
NlpGenDenseLinsys

Constructor & Destructor Documentation

NlpGenLinsys::NlpGenLinsys ( NlpGen factory,
NlpGenData data,
LinearAlgebraPackage la 
)
NlpGenLinsys::NlpGenLinsys ( )
inline
NlpGenLinsys::~NlpGenLinsys ( )

Member Function Documentation

void NlpGenLinsys::computeDiagonals ( OoqpVector dd,
OoqpVector omega,
OoqpVector t,
OoqpVector lambda,
OoqpVector u,
OoqpVector pi,
OoqpVector v,
OoqpVector gamma,
OoqpVector w,
OoqpVector phi 
)
virtual

computes the diagonal matrices in the augmented system from the current set of variables

void NlpGenLinsys::computeQuantitiesForDualReg ( NlpGenData prob,
NlpGenVars vars,
NlpGenResiduals resid,
NlpGenVars steps,
double *  dualRegQuantities 
)
virtual

compute Quantities For the usage of Dual Reg

double NlpGenLinsys::computeResidual ( NlpGenData data,
OoqpVector res_,
OoqpVector sol_,
OoqpVector solx_,
OoqpVector sols_,
OoqpVector soly_,
OoqpVector solz_ 
)
virtual
double NlpGenLinsys::computeResidual_FullKKT ( NlpGenData data,
NlpGenResiduals res_,
NlpGenVars sol_,
NlpGenVars var_ 
)
virtual
void NlpGenLinsys::copyXSYZ_fromArray ( OoqpVector vec_xsyz,
double *  array_in,
const int  nb_col 
)
virtual
void NlpGenLinsys::copyXSYZ_toArray ( OoqpVector vec_xsyz,
double *  array_in,
const int  nb_col 
)
virtual
double NlpGenLinsys::eval_xWx ( NlpGenData prob,
NlpGenResiduals resid,
NlpGenVars steps 
)
virtual

check curvature xWx

void NlpGenLinsys::factor ( Data prob,
Variables vars,
RegularizationAlg RegInfo 
)
virtual

sets up the matrix for the main linear system in "augmented system" form. The actual factorization is performed by a routine specific to either the sparse or dense case.

See also
NlpGenSparseLinsys::factor
NlpGenDenseLinsys::factor

Reimplemented in NlpGenSparseLinsys, and sLinsys.

void NlpGenLinsys::factor ( Data prob,
Variables vars 
)
virtual

factorizes the matrix, stores data related to the factorization to prepare for later calls to "solve"

Implements LinearSystem.

Reimplemented in NlpGenSparseLinsys, and sLinsys.

void NlpGenLinsys::factorNoMatChange ( Data prob_in,
Variables vars_in,
RegularizationAlg RegInfo 
)
virtual
void NlpGenLinsys::joinRHS ( OoqpVector rhs,
OoqpVector rhs1,
OoqpVector rhs2,
OoqpVector rhs3 
)
virtual

assembles a single vector object from three given vectors

Parameters
rhs(output) final joined vector
rhs1(input) first part of rhs
rhs2(input) middle part of rhs
rhs3(input) last part of rhs

Reimplemented in sLinsys, and sDummyLinsys.

void NlpGenLinsys::joinRHSXSYZ ( OoqpVector rhs_in,
OoqpVector rhs1_in,
OoqpVector rhs2_in,
OoqpVector rhs3_in,
OoqpVector rhs4_in 
)
virtual

Reimplemented in sLinsys, and sDummyLinsys.

void NlpGenLinsys::matXSYZMult ( double  beta,
OoqpVector res_,
double  alpha,
OoqpVector sol_,
NlpGenData data,
OoqpVector solx_,
OoqpVector sols_,
OoqpVector soly_,
OoqpVector solz_ 
)
virtual

res = beta*res + alpha*mat*sol stepx, stepy, stepz are used as temporary buffers

void NlpGenLinsys::matXYZMult ( double  beta,
OoqpVector res,
double  alpha,
OoqpVector sol,
NlpGenData data,
OoqpVector solx,
OoqpVector soly,
OoqpVector solz 
)
virtual

res = beta*res - alpha*mat*sol stepx, stepy, stepz are used as temporary buffers

virtual void NlpGenLinsys::putSDiagonal ( OoqpVector xdiag)
pure virtual

places the diagonal resulting from the bounds on Cx into the augmented system matrix — corresponding to slack variables S

Implemented in sLinsysLeaf, sLinsysRoot, sLinsys, NlpGenSparseLinsys, and sDummyLinsys.

virtual void NlpGenLinsys::putXDiagonal ( OoqpVector xdiag)
pure virtual

places the diagonal resulting from the bounds on x into the augmented system matrix — corresponding to pri variables X

Implemented in sLinsysLeaf, sLinsysRoot, sLinsys, NlpGenSparseLinsys, and sDummyLinsys.

virtual void NlpGenLinsys::putYDualDiagonal ( OoqpVector ydiag)
pure virtual

places the diagonal resulting from the bounds on x into the augmented system matrix — for Regularization dual Y

Implemented in sLinsysLeaf, sLinsysRoot, sLinsys, NlpGenSparseLinsys, and sDummyLinsys.

virtual void NlpGenLinsys::putZDiagonal ( OoqpVector sdiag)
pure virtual

places the diagonal resulting from the bounds on Cx into the augmented system matrix — corresponding to slack variables S, z is the dual for s

Implemented in sLinsysLeaf, sLinsysRoot, sLinsys, NlpGenSparseLinsys, and sDummyLinsys.

void NlpGenLinsys::separateVars ( OoqpVector vars1,
OoqpVector vars2,
OoqpVector vars3,
OoqpVector vars 
)
virtual

extracts three component vectors from a given aggregated vector.

Parameters
vars(input) aggregated vector
vars1(output) first part of vars
vars2(output) middle part of vars
vars3(output) last part of vars

Reimplemented in sLinsys, and sDummyLinsys.

void NlpGenLinsys::separateVarsXSYZ ( OoqpVector x_in,
OoqpVector s_in,
OoqpVector y_in,
OoqpVector z_in,
OoqpVector vars_in 
)
virtual

Reimplemented in sLinsys, and sDummyLinsys.

virtual void NlpGenLinsys::setSDiagonal ( OoqpVector sdiag)
pure virtual
virtual void NlpGenLinsys::setXDiagonal ( OoqpVector xdiag)
pure virtual
virtual void NlpGenLinsys::setYDiagonal ( OoqpVector ydiag)
pure virtual
virtual void NlpGenLinsys::setZDiagonal ( OoqpVector zdiag)
pure virtual
void NlpGenLinsys::solve ( Data prob,
Variables vars,
Residuals res,
Variables step 
)
virtual

solves the system for a given set of residuals. Assembles the right-hand side appropriate to the matrix factored in factor, solves the system using the factorization produced there, partitions the solution vector into step components, then recovers the step components eliminated during the block elimination that produced the augmented system form

See also
NlpGenSparseLinsys::solveCompressed
NlpGenDenseLinsys::solveCompressed

Implements LinearSystem.

void NlpGenLinsys::solve_IterRefine ( Data prob_in,
Variables vars_in,
Residuals res_in,
Variables step_in,
Residuals KKT_Resid_in,
Variables KKT_sol_in 
)
virtual
void NlpGenLinsys::solve_NTsteps ( Data prob,
Variables vars,
Residuals resids,
Variables Nstep,
Variables Tstep,
Variables NTstep 
)
virtual
void NlpGenLinsys::solveBiCGStab ( OoqpVector stepx,
OoqpVector steps,
OoqpVector stepy,
OoqpVector stepz,
NlpGenData data 
)
virtual
virtual void NlpGenLinsys::solveCompressed ( OoqpVector rhs)
pure virtual

perform the actual solve using the factors produced in factor.

Parameters
rhson input contains the aggregated right-hand side of the augmented system; on output contains the solution in aggregated form
See also
NlpGenSparseLinsys::solveCompressed
NlpGenDenseLinsys::solveCompressed

Implemented in sLinsys, NlpGenSparseLinsys, and sDummyLinsys.

void NlpGenLinsys::solveCompressedAugXSYZ ( OoqpVector stepx,
OoqpVector steps,
OoqpVector stepy,
OoqpVector stepz,
NlpGenData prob 
)
virtual
void NlpGenLinsys::solveCompressedAugXSYZ_PETSC ( OoqpVector stepx,
OoqpVector steps,
OoqpVector stepy,
OoqpVector stepz,
NlpGenData prob 
)
virtual
void NlpGenLinsys::solveCompressedBiCGStab ( OoqpVector stepx,
OoqpVector stepy,
OoqpVector stepz,
NlpGenData data 
)
virtual
void NlpGenLinsys::solveCompressedIterRefin ( OoqpVector stepx,
OoqpVector stepy,
OoqpVector stepz,
NlpGenData data 
)
virtual
void NlpGenLinsys::solveXYZS ( OoqpVector stepx,
OoqpVector stepy,
OoqpVector stepz,
OoqpVector steps,
OoqpVector ztemp,
NlpGenData data 
)
virtual

assemble right-hand side of augmented system and call solveCompressed to solve it

virtual void NlpGenLinsys::UpdateMatrices ( Data prob_in,
int const  updateLevel = 2 
)
inlinevirtual

Member Data Documentation

OoqpVector* NlpGenLinsys::additiveDiag

diag part: X^{-1}Z

bool NlpGenLinsys::allocateSpace
OoqpVector* NlpGenLinsys::dd

temporary storage vectors

OoqpVector * NlpGenLinsys::dq
double NlpGenLinsys::dualReg
NlpGen* NlpGenLinsys::factory
bool NlpGenLinsys::fullQ
protected
OoqpVector * NlpGenLinsys::iclow
OoqpVector * NlpGenLinsys::icupp
OoqpVector * NlpGenLinsys::ixlow
OoqpVector* NlpGenLinsys::ixupp

index matrices for the upper and lower bounds on x and Cx

int NlpGenLinsys::KryIter
int NlpGenLinsys::mclow
int NlpGenLinsys::mcupp
long long NlpGenLinsys::my
long long NlpGenLinsys::mz
OoqpVector* NlpGenLinsys::nomegaInv

stores a critical diagonal matrix as a vector

int NlpGenLinsys::num_slacks
long long NlpGenLinsys::nx

dimensions of the vectors in the general NLP formulation

int NlpGenLinsys::nxlow
int NlpGenLinsys::nxupp

dimensions of the upper and lower bound vectors

double NlpGenLinsys::priReg

primal reg and dual reg

OoqpVector * NlpGenLinsys::res
OoqpVector * NlpGenLinsys::res2
OoqpVector * NlpGenLinsys::res2Bicg
OoqpVector * NlpGenLinsys::res3
OoqpVector * NlpGenLinsys::res3Bicg
OoqpVector * NlpGenLinsys::res4
OoqpVector * NlpGenLinsys::res4Bicg
OoqpVector * NlpGenLinsys::res5
OoqpVector * NlpGenLinsys::res5Bicg
OoqpVector * NlpGenLinsys::ress
OoqpVector * NlpGenLinsys::resx
OoqpVector * NlpGenLinsys::resy
OoqpVector * NlpGenLinsys::resz
OoqpVector* NlpGenLinsys::rhs

right-hand side of the system

OoqpVector* NlpGenLinsys::rhs_back

copy of right-hand side of the system, to compute residual

OoqpVector* NlpGenLinsys::sol

Work vectors for iterative refinement of the XYZ linear system

OoqpVector* NlpGenLinsys::sol2

Work vectors for BiCGStab and iterative refinement

OoqpVector * NlpGenLinsys::sol2Bicg
OoqpVector * NlpGenLinsys::temp_diagS
OoqpVector * NlpGenLinsys::temp_diagX
OoqpVector * NlpGenLinsys::temp_diagY
OoqpVector * NlpGenLinsys::temp_diagZ
int NlpGenLinsys::useRefs

The documentation for this class was generated from the following files: