robotoc
robotoc - efficient ROBOT Optimal Control solvers
Loading...
Searching...
No Matches
cost_function.hpp
Go to the documentation of this file.
1#ifndef ROBOTOC_COST_FUNCTION_HPP_
2#define ROBOTOC_COST_FUNCTION_HPP_
3
4#include <memory>
5#include <cmath>
6#include <vector>
7#include <unordered_map>
8#include <iostream>
9
10#include "Eigen/Core"
11
21
22
23namespace robotoc {
24
31public:
33 = std::shared_ptr<CostFunctionComponentBase>;
34
42 CostFunction(const double discount_factor, const double discount_time_step);
43
48
52 ~CostFunction() = default;
53
57 CostFunction(const CostFunction&) = default;
58
63
67 CostFunction(CostFunction&&) noexcept = default;
68
72 CostFunction& operator=(CostFunction&&) noexcept = default;
73
82 void setDiscountFactor(const double discount_factor,
83 const double discount_time_step);
84
89 double discountFactor() const;
90
95 double discountTimeStep() const;
96
102 bool exist(const std::string& name) const;
103
110 void add(const std::string& name, const CostFunctionComponentBasePtr& cost);
111
117 void erase(const std::string& name);
118
125 CostFunctionComponentBasePtr get(const std::string& name) const;
126
130 void clear();
131
139
149 double evalStageCost(Robot& robot, const ContactStatus& contact_status,
150 CostFunctionData& data, const GridInfo& grid_info,
151 const SplitSolution& s) const;
152
164 double linearizeStageCost(Robot& robot, const ContactStatus& contact_status,
165 CostFunctionData& data, const GridInfo& grid_info,
166 const SplitSolution& s,
167 SplitKKTResidual& kkt_residual) const;
168
183 double quadratizeStageCost(Robot& robot, const ContactStatus& contact_status,
184 CostFunctionData& data, const GridInfo& grid_info,
185 const SplitSolution& s,
186 SplitKKTResidual& kkt_residual,
187 SplitKKTMatrix& kkt_matrix) const;
188
198 const GridInfo& grid_info,
199 const SplitSolution& s) const;
200
212 const GridInfo& grid_info, const SplitSolution& s,
213 SplitKKTResidual& kkt_residual) const;
214
229 const GridInfo& grid_info, const SplitSolution& s,
230 SplitKKTResidual& kkt_residual,
231 SplitKKTMatrix& kkt_matrix) const;
232
242 double evalImpactCost(Robot& robot, const ImpactStatus& impact_status,
243 CostFunctionData& data, const GridInfo& grid_info,
244 const SplitSolution& s) const;
245
257 double linearizeImpactCost(Robot& robot, const ImpactStatus& impact_status,
258 CostFunctionData& data, const GridInfo& grid_info,
259 const SplitSolution& s,
260 SplitKKTResidual& kkt_residual) const;
261
276 double quadratizeImpactCost(Robot& robot, const ImpactStatus& impact_status,
277 CostFunctionData& data, const GridInfo& grid_info,
278 const SplitSolution& s,
279 SplitKKTResidual& kkt_residual,
280 SplitKKTMatrix& kkt_matrix) const;
281
286 std::vector<std::string> getCostComponentList() const;
287
291 void disp(std::ostream& os) const;
292
293 friend std::ostream& operator<<(std::ostream& os,
294 const CostFunction& cost_function);
295
296 friend std::ostream& operator<<(std::ostream& os,
297 const std::shared_ptr<CostFunction>& cost_function);
298
299private:
300 std::vector<CostFunctionComponentBasePtr> costs_;
301 std::unordered_map<std::string, size_t> cost_names_;
302
303 double discount(const double t0, const double t) const {
304 assert(t >= t0);
305 return std::pow(discount_factor_, ((t-t0)/discount_time_step_));
306 }
307
308 double discount_factor_, discount_time_step_;
309 bool discounted_cost_;
310};
311
312} // namespace robotoc
313
314#endif // ROBOTOC_COST_FUNCTION_HPP_
Contact status of robot model.
Definition: contact_status.hpp:32
Stack of the cost function. Composed by cost function components that inherits CostFunctionComponentB...
Definition: cost_function.hpp:30
CostFunctionData createCostFunctionData(const Robot &robot) const
Creates CostFunctionData according to robot model and cost function components.
CostFunction(const double discount_factor, const double discount_time_step)
Constructor with discount factor.
CostFunction & operator=(const CostFunction &)=default
Default copy operator.
void clear()
Clear cost function by removing all components.
std::vector< std::string > getCostComponentList() const
Gets a list of the cost components.
double linearizeImpactCost(Robot &robot, const ImpactStatus &impact_status, CostFunctionData &data, const GridInfo &grid_info, const SplitSolution &s, SplitKKTResidual &kkt_residual) const
Computes the impact cost and its first-order partial derivatives.
void disp(std::ostream &os) const
Displays the cost function onto a ostream.
void add(const std::string &name, const CostFunctionComponentBasePtr &cost)
Adds a cost function component. If a component of the same name exists, throws an exeption.
void setDiscountFactor(const double discount_factor, const double discount_time_step)
Sets the discount facor.
~CostFunction()=default
Destructor.
CostFunctionComponentBasePtr get(const std::string &name) const
Gets a cost function component. If a component of the specified name does not exist,...
double linearizeTerminalCost(Robot &robot, CostFunctionData &data, const GridInfo &grid_info, const SplitSolution &s, SplitKKTResidual &kkt_residual) const
Computes the terminal cost and its first-order partial derivatives.
CostFunction(const CostFunction &)=default
Default copy constructor.
void erase(const std::string &name)
Erases a cost function component. If a component of the specified name does not exist,...
bool exist(const std::string &name) const
Checks if thsi has a cost function component of the specified name.
double evalTerminalCost(Robot &robot, CostFunctionData &data, const GridInfo &grid_info, const SplitSolution &s) const
Computes the terminal cost.
double quadratizeStageCost(Robot &robot, const ContactStatus &contact_status, CostFunctionData &data, const GridInfo &grid_info, const SplitSolution &s, SplitKKTResidual &kkt_residual, SplitKKTMatrix &kkt_matrix) const
Computes the stage cost, its first-order partial derivatives, and its Hessian, i.e....
std::shared_ptr< CostFunctionComponentBase > CostFunctionComponentBasePtr
Definition: cost_function.hpp:33
double quadratizeImpactCost(Robot &robot, const ImpactStatus &impact_status, CostFunctionData &data, const GridInfo &grid_info, const SplitSolution &s, SplitKKTResidual &kkt_residual, SplitKKTMatrix &kkt_matrix) const
Computes the impact cost, its first-order partial derivatives, and its Hessian, i....
double discountFactor() const
Gets the discount facor.
CostFunction()
Default constructor.
double discountTimeStep() const
Gets the discount time step.
double quadratizeTerminalCost(Robot &robot, CostFunctionData &data, const GridInfo &grid_info, const SplitSolution &s, SplitKKTResidual &kkt_residual, SplitKKTMatrix &kkt_matrix) const
Computes the terminal cost, its first-order partial derivatives, and its Hessian, i....
double evalStageCost(Robot &robot, const ContactStatus &contact_status, CostFunctionData &data, const GridInfo &grid_info, const SplitSolution &s) const
Computes the stage cost.
double evalImpactCost(Robot &robot, const ImpactStatus &impact_status, CostFunctionData &data, const GridInfo &grid_info, const SplitSolution &s) const
Computes the impact cost.
CostFunction(CostFunction &&) noexcept=default
Default move constructor.
double linearizeStageCost(Robot &robot, const ContactStatus &contact_status, CostFunctionData &data, const GridInfo &grid_info, const SplitSolution &s, SplitKKTResidual &kkt_residual) const
Computes the stage cost and its first-order partial derivatives.
Impact status of robot model. Wrapper of ContactStatus to treat impacts.
Definition: impact_status.hpp:21
Dynamics and kinematics model of robots. Wraps pinocchio::Model and pinocchio::Data....
Definition: robot.hpp:32
The KKT matrix split into a time stage.
Definition: split_kkt_matrix.hpp:18
KKT residual split into each time stage.
Definition: split_kkt_residual.hpp:18
Solution to the optimal control problem split into a time stage.
Definition: split_solution.hpp:20
Definition: constraint_component_base.hpp:17
Composed of data used to compute the cost function and its derivatives.
Definition: cost_function_data.hpp:17
Grid information.
Definition: grid_info.hpp:24