robotoc
robotoc - efficient ROBOT Optimal Control solvers
Loading...
Searching...
No Matches
unconstr_backward_correction.hpp
Go to the documentation of this file.
1#ifndef ROBOTOC_UNCONSTR_BACKWARD_CORRECTION_HPP_
2#define ROBOTOC_UNCONSTR_BACKWARD_CORRECTION_HPP_
3
4#include <vector>
5
6#include "Eigen/Core"
7
18#include "robotoc/ocp/ocp.hpp"
24
25
26namespace robotoc {
27
34public:
41 UnconstrBackwardCorrection(const OCP& ocp, const int nthreads);
42
47
52
57
62
67
72
78 void setNumThreads(const int nthreads);
79
90 const std::vector<GridInfo>& time_discretization,
91 const Solution& s, KKTMatrix& kkt_matrix,
92 KKTResidual& kkt_residual);
93
101 const std::vector<GridInfo>& time_discretization,
102 const Solution& s);
103
114 const std::vector<GridInfo>& time_discretization,
115 const Eigen::VectorXd& q, const Eigen::VectorXd& v,
116 const Solution& s, KKTResidual& kkt_residual);
117
129 const std::vector<GridInfo>& time_discretization,
130 const Eigen::VectorXd& q, const Eigen::VectorXd& v,
131 const Solution& s, KKTMatrix& kkt_matrix,
132 KKTResidual& kkt_residual);
133
146 const std::vector<GridInfo>& time_discretization,
147 const Eigen::VectorXd& q, const Eigen::VectorXd& v,
148 const Solution& s, KKTMatrix& kkt_matrix,
149 KKTResidual& kkt_residual);
150
160 void backwardCorrection(const std::vector<GridInfo>& time_discretization,
161 const Solution& s, const KKTMatrix& kkt_matrix,
162 const KKTResidual& kkt_residual, Direction& d);
163
168 double primalStepSize() const;
169
174 double dualStepSize() const;
175
180 const PerformanceIndex& getEval() const {
181 return performance_index_;
182 }
183
194 const std::vector<GridInfo>& time_discretization,
195 const double primal_step_size,
196 const double dual_step_size, Direction& d, Solution& s);
197
207 const std::vector<GridInfo>& time_discretization,
208 const double primal_step_size,
209 const Direction& d, Solution& s);
210
211private:
212 int nthreads_;
213 OCP ocp_;
214 ParNMPCIntermediateStage intermediate_stage_;
215 ParNMPCTerminalStage terminal_stage_;
217 PerformanceIndex performance_index_;
219 Solution s_new_;
220 std::vector<Eigen::MatrixXd> aux_mat_;
221 Eigen::VectorXd primal_step_sizes_, dual_step_sizes_;
222
223};
224
225} // namespace robotoc
226
227#endif // ROBOTOC_UNCONSTR_BACKWARD_CORRECTION_HPP_
The intermediate stage of ParNMPC for unconstrained rigid-body systems.
Definition: parnmpc_intermediate_stage.hpp:28
The terminal stage of ParNMPC for unconstrained rigid-body systems.
Definition: parnmpc_terminal_stage.hpp:28
Dynamics and kinematics model of robots. Wraps pinocchio::Model and pinocchio::Data....
Definition: robot.hpp:32
Backward correction for optimal control problems of unconstrained rigid-body systems.
Definition: unconstr_backward_correction.hpp:33
double dualStepSize() const
Gets the maximum dual step size of the fraction-to-boundary-rule.
void initConstraints(aligned_vector< Robot > &robots, const std::vector< GridInfo > &time_discretization, const Solution &s)
Computes the cost and constraint violations.
UnconstrBackwardCorrection(const OCP &ocp, const int nthreads)
Construct a backward correction.
UnconstrBackwardCorrection()
Default constructor.
void setNumThreads(const int nthreads)
Sets the number of threads of the parallel computations.
UnconstrBackwardCorrection & operator=(const UnconstrBackwardCorrection &)=default
Default copy operator.
double primalStepSize() const
Gets the maximum primal step size of the fraction-to-boundary-rule.
void coarseUpdate(aligned_vector< Robot > &robots, const std::vector< GridInfo > &time_discretization, const Eigen::VectorXd &q, const Eigen::VectorXd &v, const Solution &s, KKTMatrix &kkt_matrix, KKTResidual &kkt_residual)
Eval KKT and coarse updates the solution leveraging the parallel computation.
UnconstrBackwardCorrection(UnconstrBackwardCorrection &&) noexcept=default
Default move constructor.
UnconstrBackwardCorrection(const UnconstrBackwardCorrection &)=default
Default copy constructor.
~UnconstrBackwardCorrection()=default
Destructor.
void integratePrimalSolution(const aligned_vector< Robot > &robots, const std::vector< GridInfo > &time_discretization, const double primal_step_size, const Direction &d, Solution &s)
Integrates the solution.
const PerformanceIndex & getEval() const
Gets the performance index of the evaluation.
Definition: unconstr_backward_correction.hpp:180
void integrateSolution(const aligned_vector< Robot > &robots, const std::vector< GridInfo > &time_discretization, const double primal_step_size, const double dual_step_size, Direction &d, Solution &s)
Integrates the solution.
void initAuxMat(aligned_vector< Robot > &robots, const std::vector< GridInfo > &time_discretization, const Solution &s, KKTMatrix &kkt_matrix, KKTResidual &kkt_residual)
Initializes the auxiliary matrices by the terminal cost Hessian computed by the current solution.
void evalKKT(aligned_vector< Robot > &robots, const std::vector< GridInfo > &time_discretization, const Eigen::VectorXd &q, const Eigen::VectorXd &v, const Solution &s, KKTMatrix &kkt_matrix, KKTResidual &kkt_residual)
Computes the KKT residual and matrix.
void backwardCorrection(const std::vector< GridInfo > &time_discretization, const Solution &s, const KKTMatrix &kkt_matrix, const KKTResidual &kkt_residual, Direction &d)
Performs the backward correction for coarse updated solution and computes the Newton direction.
void evalOCP(aligned_vector< Robot > &robots, const std::vector< GridInfo > &time_discretization, const Eigen::VectorXd &q, const Eigen::VectorXd &v, const Solution &s, KKTResidual &kkt_residual)
Computes the cost and constraint violations.
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
Performance index of optimal control problems.
Definition: performance_index.hpp:13