Drake
Drake C++ Documentation

Detailed Description

Implements a cost of the form

\[ .5 x'Qx + b'x + c \]

.

#include <drake/solvers/cost.h>

Public Member Functions

template<typename DerivedQ , typename Derivedb >
 QuadraticCost (const Eigen::MatrixBase< DerivedQ > &Q, const Eigen::MatrixBase< Derivedb > &b, double c=0., std::optional< bool > is_hessian_psd=std::nullopt)
 Constructs a cost of the form

\[ .5 x'Qx + b'x + c \]

. More...

 
 ~QuadraticCost () override
 
const Eigen::MatrixXd & Q () const
 Returns the symmetric matrix Q, as the Hessian of the cost. More...
 
const Eigen::VectorXd & b () const
 
double c () const
 
bool is_convex () const
 Returns true if this cost is convex. More...
 
template<typename DerivedQ , typename DerivedB >
void UpdateCoefficients (const Eigen::MatrixBase< DerivedQ > &new_Q, const Eigen::MatrixBase< DerivedB > &new_b, double new_c=0., std::optional< bool > is_hessian_psd=std::nullopt)
 Updates the quadratic and linear term of the constraint. More...
 
Does not allow copy, move, or assignment
 QuadraticCost (const QuadraticCost &)=delete
 
QuadraticCostoperator= (const QuadraticCost &)=delete
 
 QuadraticCost (QuadraticCost &&)=delete
 
QuadraticCostoperator= (QuadraticCost &&)=delete
 
- Public Member Functions inherited from Cost
 Cost (const Cost &)=delete
 
Costoperator= (const Cost &)=delete
 
 Cost (Cost &&)=delete
 
Costoperator= (Cost &&)=delete
 
- Public Member Functions inherited from EvaluatorBase
virtual ~EvaluatorBase ()
 
void Eval (const Eigen::Ref< const Eigen::VectorXd > &x, Eigen::VectorXd *y) const
 Evaluates the expression. More...
 
void Eval (const Eigen::Ref< const AutoDiffVecXd > &x, AutoDiffVecXd *y) const
 Evaluates the expression. More...
 
void Eval (const Eigen::Ref< const VectorX< symbolic::Variable >> &x, VectorX< symbolic::Expression > *y) const
 Evaluates the expression. More...
 
void set_description (const std::string &description)
 Set a human-friendly description for the evaluator. More...
 
const std::string & get_description () const
 Getter for a human-friendly description for the evaluator. More...
 
std::ostream & Display (std::ostream &os, const VectorX< symbolic::Variable > &vars) const
 Formats this evaluator into the given stream using vars for the bound decision variable names. More...
 
std::ostream & Display (std::ostream &os) const
 Formats this evaluator into the given stream, without displaying the decision variables it is bound to. More...
 
std::string ToLatex (const VectorX< symbolic::Variable > &vars, int precision=3) const
 Returns a LaTeX string describing this evaluator. More...
 
int num_vars () const
 Getter for the number of variables, namely the number of rows in x, as used in Eval(x, y). More...
 
int num_outputs () const
 Getter for the number of outputs, namely the number of rows in y, as used in Eval(x, y). More...
 
void SetGradientSparsityPattern (const std::vector< std::pair< int, int >> &gradient_sparsity_pattern)
 Set the sparsity pattern of the gradient matrix ∂y/∂x (the gradient of y value in Eval, w.r.t x in Eval) . More...
 
const std::optional< std::vector< std::pair< int, int > > > & gradient_sparsity_pattern () const
 Returns the vector of (row_index, col_index) that contains all the entries in the gradient of Eval function (∂y/∂x) whose value could be non-zero, namely if ∂yᵢ/∂xⱼ could be non-zero, then the pair (i, j) is in gradient_sparsity_pattern. More...
 
bool is_thread_safe () const
 Returns whether it is safe to call Eval in parallel. More...
 
 EvaluatorBase (const EvaluatorBase &)=delete
 
EvaluatorBaseoperator= (const EvaluatorBase &)=delete
 
 EvaluatorBase (EvaluatorBase &&)=delete
 
EvaluatorBaseoperator= (EvaluatorBase &&)=delete
 

Additional Inherited Members

- Protected Member Functions inherited from Cost
 Cost (int num_vars, const std::string &description="")
 Constructs a cost evaluator. More...
 
- Protected Member Functions inherited from EvaluatorBase
 EvaluatorBase (int num_outputs, int num_vars, const std::string &description="")
 Constructs a evaluator. More...
 
void set_num_outputs (int num_outputs)
 
void set_is_thread_safe (bool is_thread_safe)
 

Constructor & Destructor Documentation

◆ QuadraticCost() [1/3]

QuadraticCost ( const QuadraticCost )
delete

◆ QuadraticCost() [2/3]

QuadraticCost ( QuadraticCost &&  )
delete

◆ QuadraticCost() [3/3]

QuadraticCost ( const Eigen::MatrixBase< DerivedQ > &  Q,
const Eigen::MatrixBase< Derivedb > &  b,
double  c = 0.,
std::optional< bool >  is_hessian_psd = std::nullopt 
)

Constructs a cost of the form

\[ .5 x'Qx + b'x + c \]

.

Parameters
QQuadratic term.
bLinear term.
c(optional) Constant term.
is_hessian_psd(optional) Indicates if the Hessian matrix Q is positive semidefinite (psd) or not. If set to true, then the user guarantees that Q is psd; if set to false, then the user guarantees that Q is not psd. If set to std::nullopt, then the constructor will check if Q is psd or not. The default is std::nullopt. To speed up the constructor, set is_hessian_psd to either true or false.

◆ ~QuadraticCost()

~QuadraticCost ( )
override

Member Function Documentation

◆ b()

const Eigen::VectorXd& b ( ) const

◆ c()

double c ( ) const

◆ is_convex()

bool is_convex ( ) const

Returns true if this cost is convex.

A quadratic cost if convex if and only if its Hessian matrix Q is positive semidefinite.

◆ operator=() [1/2]

QuadraticCost& operator= ( const QuadraticCost )
delete

◆ operator=() [2/2]

QuadraticCost& operator= ( QuadraticCost &&  )
delete

◆ Q()

const Eigen::MatrixXd& Q ( ) const

Returns the symmetric matrix Q, as the Hessian of the cost.

◆ UpdateCoefficients()

void UpdateCoefficients ( const Eigen::MatrixBase< DerivedQ > &  new_Q,
const Eigen::MatrixBase< DerivedB > &  new_b,
double  new_c = 0.,
std::optional< bool >  is_hessian_psd = std::nullopt 
)

Updates the quadratic and linear term of the constraint.

The new matrices need to have the same dimension as before.

Parameters
new_QNew quadratic term.
new_bNew linear term.
new_c(optional) New constant term.
is_hessian_psd(optional) Indicates if the Hessian matrix Q is positive semidefinite (psd) or not. If set to true, then the user guarantees that Q is psd; if set to false, then the user guarantees that Q is not psd. If set to std::nullopt, then this function will check if Q is psd or not. The default is std::nullopt. To speed up the computation, set is_hessian_psd to either true or false.

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