robotoc
robotoc - efficient ROBOT Optimal Control solvers
Loading...
Searching...
No Matches
unconstr_parnmpc_solver.hpp
Go to the documentation of this file.
1#ifndef ROBOTOC_UNCONSTR_PARNMPC_SOLVER_HPP_
2#define ROBOTOC_UNCONSTR_PARNMPC_SOLVER_HPP_
3
4#include <vector>
5#include <memory>
6
7#include "Eigen/Core"
8
19#include "robotoc/ocp/ocp.hpp"
24
25
26namespace robotoc {
27
35public:
43 const SolverOptions& solver_options=SolverOptions());
44
49
54
59
64
69
73 UnconstrParNMPCSolver& operator=(UnconstrParNMPCSolver&&) noexcept = default;
74
79 void setSolverOptions(const SolverOptions& solver_options);
80
85 void discretize(const double t);
86
92
97
108 void solve(const double t, const Eigen::VectorXd& q, const Eigen::VectorXd& v,
109 const bool init_solver=true);
110
116
124 const SplitSolution& getSolution(const int stage) const;
125
131 std::vector<Eigen::VectorXd> getSolution(const std::string& name) const;
132
138 void setSolution(const std::string& name, const Eigen::VectorXd& value);
139
148 double KKTError(const double t, const Eigen::VectorXd& q,
149 const Eigen::VectorXd& v);
150
156 double KKTError() const;
157
162 const std::vector<GridInfo>& getTimeDiscretization() const;
163
169 void setRobotProperties(const RobotProperties& properties);
170
171private:
172 aligned_vector<Robot> robots_;
173 std::vector<GridInfo> time_discretization_;
174 UnconstrBackwardCorrection backward_correction_;
175 UnconstrLineSearch line_search_;
176 OCP ocp_;
177 KKTMatrix kkt_matrix_;
178 KKTResidual kkt_residual_;
179 Solution s_;
180 Direction d_;
181 SolverOptions solver_options_;
182 SolverStatistics solver_statistics_;
183 Timer timer_;
184
192 void updateSolution(const double t, const Eigen::VectorXd& q,
193 const Eigen::VectorXd& v);
194
195};
196
197} // namespace robotoc
198
199#endif // ROBOTOC_UNCONSTR_PARNMPC_SOLVER_HPP_
Dynamics and kinematics model of robots. Wraps pinocchio::Model and pinocchio::Data....
Definition: robot.hpp:32
Solution to the optimal control problem split into a time stage.
Definition: split_solution.hpp:20
A timer class to take benchmarks.
Definition: timer.hpp:12
Backward correction for optimal control problems of unconstrained rigid-body systems.
Definition: unconstr_backward_correction.hpp:33
Line search for optimal control problems for unconstrained rigid-body systems.
Definition: unconstr_line_search.hpp:29
Optimal control problem solver of unconstrained rigid-body systems by ParNMPC algorithm....
Definition: unconstr_parnmpc_solver.hpp:34
void initBackwardCorrection()
Initializes the backward correction solver.
~UnconstrParNMPCSolver()=default
Destructor.
void solve(const double t, const Eigen::VectorXd &q, const Eigen::VectorXd &v, const bool init_solver=true)
Solves the optimal control problem. Internally calls updateSolution().
UnconstrParNMPCSolver & operator=(const UnconstrParNMPCSolver &)=default
Default copy assign operator.
double KKTError() const
Returns the l2-norm of the KKT residuals using the results of UnconstrParNMPCsolver::updateSolution()...
UnconstrParNMPCSolver(const UnconstrParNMPCSolver &)=default
Default copy constructor.
UnconstrParNMPCSolver(const OCP &ocp, const SolverOptions &solver_options=SolverOptions())
Construct optimal control problem solver.
void setRobotProperties(const RobotProperties &properties)
Sets a collection of the properties for robot model in this solver.
void discretize(const double t)
Discretizes the problem and reiszes the data structures.
UnconstrParNMPCSolver()
Default constructor.
void initConstraints()
Initializes the priaml-dual interior point method for inequality constraints.
UnconstrParNMPCSolver(UnconstrParNMPCSolver &&) noexcept=default
Default move constructor.
void setSolverOptions(const SolverOptions &solver_options)
Sets the solver option.
void setSolution(const std::string &name, const Eigen::VectorXd &value)
Sets the solution over the horizon.
const SolverStatistics & getSolverStatistics() const
Gets the solver statistics.
const std::vector< GridInfo > & getTimeDiscretization() const
Gets the discretization.
const SplitSolution & getSolution(const int stage) const
Get the split solution of a time stage. For example, the control input torques at the initial stage c...
Definition: constraint_component_base.hpp:17
aligned_vector< SplitKKTMatrix > KKTMatrix
The KKT matrix of the optimal control problem.
Definition: kkt_matrix.hpp:16
aligned_vector< SplitSolution > Solution
Solution to the optimal control problem.
Definition: solution.hpp:16
aligned_vector< SplitDirection > Direction
Newton direction of the solution to the optimal control problem.
Definition: direction.hpp:16
std::vector< T, Eigen::aligned_allocator< T > > aligned_vector
std vector with Eigen::aligned_allocator.
Definition: aligned_vector.hpp:14
aligned_vector< SplitKKTResidual > KKTResidual
The KKT residual of the optimal control problem.
Definition: kkt_residual.hpp:16
Grid information.
Definition: grid_info.hpp:24
The optimal control problem.
Definition: ocp.hpp:22
Collection of the robot properties, which can change after constructing robot models.
Definition: robot_properties.hpp:30
Options of optimal control solvers.
Definition: solver_options.hpp:17
Statistics of optimal control solvers.
Definition: solver_statistics.hpp:17