Google OR-Tools v9.11
a fast and portable software suite for combinatorial optimization
Loading...
Searching...
No Matches
operations_research::pdlp::ShardedQuadraticProgram Class Reference

#include <sharded_quadratic_program.h>

Public Member Functions

 ShardedQuadraticProgram (QuadraticProgram qp, int num_threads, int num_shards, operations_research::SolverLogger *logger=nullptr)
 
 ShardedQuadraticProgram (const ShardedQuadraticProgram &)=delete
 Movable but not copyable.
 
ShardedQuadraticProgramoperator= (const ShardedQuadraticProgram &)=delete
 
 ShardedQuadraticProgram (ShardedQuadraticProgram &&)=default
 
ShardedQuadraticProgramoperator= (ShardedQuadraticProgram &&)=default
 
const QuadraticProgramQp () const
 
const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > & TransposedConstraintMatrix () const
 Returns a reference to the transpose of the QP's constraint matrix.
 
const SharderConstraintMatrixSharder () const
 Returns a Sharder intended for the columns of the QP's constraint matrix.
 
const SharderTransposedConstraintMatrixSharder () const
 Returns a Sharder intended for the rows of the QP's constraint matrix.
 
const SharderPrimalSharder () const
 Returns a Sharder intended for primal vectors.
 
const SharderDualSharder () const
 Returns a Sharder intended for dual vectors.
 
int64_t PrimalSize () const
 
int64_t DualSize () const
 
void RescaleQuadraticProgram (const Eigen::VectorXd &col_scaling_vec, const Eigen::VectorXd &row_scaling_vec)
 
void SwapVariableBounds (Eigen::VectorXd &variable_lower_bounds, Eigen::VectorXd &variable_upper_bounds)
 
void SwapConstraintBounds (Eigen::VectorXd &constraint_lower_bounds, Eigen::VectorXd &constraint_upper_bounds)
 
void SetConstraintBounds (int64_t constraint_index, std::optional< double > lower_bound, std::optional< double > upper_bound)
 
void SwapObjectiveVector (Eigen::VectorXd &objective)
 
void ReplaceLargeConstraintBoundsWithInfinity (double threshold)
 

Detailed Description

This class stores:

  • A QuadraticProgram (QP)
  • A transposed version of the QP's constraint matrix
  • A thread pool
  • Various Sharder objects for doing sharded matrix and vector computations.

Definition at line 37 of file sharded_quadratic_program.h.

Constructor & Destructor Documentation

◆ ShardedQuadraticProgram() [1/3]

operations_research::pdlp::ShardedQuadraticProgram::ShardedQuadraticProgram ( QuadraticProgram qp,
int num_threads,
int num_shards,
operations_research::SolverLogger * logger = nullptr )

Requires num_shards >= num_threads >= 1.

Note
the qp is intentionally passed by value. If logger is not nullptr, warns about unbalanced matrices using it; otherwise warns via Google standard logging.

Definition at line 77 of file sharded_quadratic_program.cc.

◆ ShardedQuadraticProgram() [2/3]

operations_research::pdlp::ShardedQuadraticProgram::ShardedQuadraticProgram ( const ShardedQuadraticProgram & )
delete

Movable but not copyable.

◆ ShardedQuadraticProgram() [3/3]

operations_research::pdlp::ShardedQuadraticProgram::ShardedQuadraticProgram ( ShardedQuadraticProgram && )
default

Member Function Documentation

◆ ConstraintMatrixSharder()

const Sharder & operations_research::pdlp::ShardedQuadraticProgram::ConstraintMatrixSharder ( ) const
inline

Returns a Sharder intended for the columns of the QP's constraint matrix.

Definition at line 61 of file sharded_quadratic_program.h.

◆ DualSharder()

const Sharder & operations_research::pdlp::ShardedQuadraticProgram::DualSharder ( ) const
inline

Returns a Sharder intended for dual vectors.

Definition at line 71 of file sharded_quadratic_program.h.

◆ DualSize()

int64_t operations_research::pdlp::ShardedQuadraticProgram::DualSize ( ) const
inline

Definition at line 74 of file sharded_quadratic_program.h.

◆ operator=() [1/2]

ShardedQuadraticProgram & operations_research::pdlp::ShardedQuadraticProgram::operator= ( const ShardedQuadraticProgram & )
delete

◆ operator=() [2/2]

ShardedQuadraticProgram & operations_research::pdlp::ShardedQuadraticProgram::operator= ( ShardedQuadraticProgram && )
default

◆ PrimalSharder()

const Sharder & operations_research::pdlp::ShardedQuadraticProgram::PrimalSharder ( ) const
inline

Returns a Sharder intended for primal vectors.

Definition at line 69 of file sharded_quadratic_program.h.

◆ PrimalSize()

int64_t operations_research::pdlp::ShardedQuadraticProgram::PrimalSize ( ) const
inline

Definition at line 73 of file sharded_quadratic_program.h.

◆ Qp()

const QuadraticProgram & operations_research::pdlp::ShardedQuadraticProgram::Qp ( ) const
inline

Definition at line 52 of file sharded_quadratic_program.h.

◆ ReplaceLargeConstraintBoundsWithInfinity()

void operations_research::pdlp::ShardedQuadraticProgram::ReplaceLargeConstraintBoundsWithInfinity ( double threshold)

Replaces constraint bounds whose absolute value is >= threshold with the corresponding infinity.

Definition at line 183 of file sharded_quadratic_program.cc.

◆ RescaleQuadraticProgram()

void operations_research::pdlp::ShardedQuadraticProgram::RescaleQuadraticProgram ( const Eigen::VectorXd & col_scaling_vec,
const Eigen::VectorXd & row_scaling_vec )

Rescales the QP (including objective, variable bounds, constraint bounds, constraint matrix, and transposed constraint matrix) based on col_scaling_vec and row_scaling_vec. That is, rescale the problem so that each variable is rescaled as variable[i] <- variable[i] / col_scaling_vec[i], and the j-th constraint is multiplied by row_scaling_vec[j]. col_scaling_vec and row_scaling_vec must be positive.

Definition at line 148 of file sharded_quadratic_program.cc.

◆ SetConstraintBounds()

void operations_research::pdlp::ShardedQuadraticProgram::SetConstraintBounds ( int64_t constraint_index,
std::optional< double > lower_bound,
std::optional< double > upper_bound )

Definition at line 191 of file sharded_quadratic_program.cc.

◆ SwapConstraintBounds()

void operations_research::pdlp::ShardedQuadraticProgram::SwapConstraintBounds ( Eigen::VectorXd & constraint_lower_bounds,
Eigen::VectorXd & constraint_upper_bounds )
inline

Definition at line 92 of file sharded_quadratic_program.h.

◆ SwapObjectiveVector()

void operations_research::pdlp::ShardedQuadraticProgram::SwapObjectiveVector ( Eigen::VectorXd & objective)
inline

Swaps objective with the objective_vector in the quadratic program. Swapping objective_matrix is not yet supported because it hasn't been needed.

Definition at line 105 of file sharded_quadratic_program.h.

◆ SwapVariableBounds()

void operations_research::pdlp::ShardedQuadraticProgram::SwapVariableBounds ( Eigen::VectorXd & variable_lower_bounds,
Eigen::VectorXd & variable_upper_bounds )
inline

Definition at line 86 of file sharded_quadratic_program.h.

◆ TransposedConstraintMatrix()

const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > & operations_research::pdlp::ShardedQuadraticProgram::TransposedConstraintMatrix ( ) const
inline

Returns a reference to the transpose of the QP's constraint matrix.

Definition at line 56 of file sharded_quadratic_program.h.

◆ TransposedConstraintMatrixSharder()

const Sharder & operations_research::pdlp::ShardedQuadraticProgram::TransposedConstraintMatrixSharder ( ) const
inline

Returns a Sharder intended for the rows of the QP's constraint matrix.

Definition at line 65 of file sharded_quadratic_program.h.


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