31 #ifndef SHOOTINGCONSTRAINT_HPP_
32 #define SHOOTINGCONSTRAINT_HPP_
34 #include <Eigen/Dense>
37 #include <boost/numeric/odeint.hpp>
38 #include <dynamical_systems/tools/eigenIntegration.h>
40 namespace DirectTrajectoryOptimization {
41 template <
typename DIMENSIONS>
42 class DTOMultipleShooting;
48 template<
class DIMENSIONS>
53 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
55 typedef typename DIMENSIONS::state_vector_t state_vector_t;
56 typedef typename DIMENSIONS::control_vector_t control_vector_t;
57 typedef typename DIMENSIONS::state_vector_array_t state_vector_array_t;
58 typedef typename DIMENSIONS::control_vector_array_t control_vector_array_t;
69 virtual Eigen::VectorXd
Evalfx(
const Eigen::VectorXd & in_vect)
override;
71 void forwardIntegrateForTimeRange(
const state_vector_t & curr_state,
const control_vector_t & curr_control,
72 const control_vector_t & next_control,
const double & total_t, state_vector_t & integrated_state);
74 void evalDefectGradient(
const Eigen::VectorXd & inputspace, Eigen::MatrixXd &jacobianMatrix);
76 Eigen::VectorXd evalDefectVector(
const state_vector_array_t & y_trajectory,
77 const control_vector_array_t & u_trajectory,
const Eigen::VectorXd & h_trajectory);
83 _ms_defect_vector.resize(_dms_->num_defects);
88 Eigen::VectorXd _ms_defect_vector;
90 std::shared_ptr<DTOMultipleShooting<DIMENSIONS> > _dms_;
91 std::shared_ptr<DynamicsBase<DIMENSIONS> > _local_dynamics_;
94 boost::numeric::odeint::runge_kutta4<state_vector_t, double,
95 state_vector_t, double,
96 boost::numeric::odeint::vector_space_algebra> ms_stepper_;
99 template<
class DIMENSIONS>
102 state_vector_array_t y_trajectory;
103 control_vector_array_t u_trajectory;
104 Eigen::VectorXd h_trajectory;
106 _dms_->getStateControlIncrementsTrajectoriesFromDecisionVariables(in_vector,y_trajectory,u_trajectory,h_trajectory);
108 Eigen::VectorXd defect = evalDefectVector(y_trajectory,u_trajectory,h_trajectory);
113 template<
class DIMENSIONS>
115 const control_vector_array_t & u_trajectory ,
const Eigen::VectorXd & h_trajectory) {
117 state_vector_t integrated_state;
118 _ms_defect_vector.setZero();
120 int current_defect_index = 0;
122 for (
int k = 0; k < _dms_->n_inc ; k++) {
124 _local_dynamics_ = _dms_->_dynamics[_dms_->node_to_phase_map[k]];
126 forwardIntegrateForTimeRange(y_trajectory[k], u_trajectory[k],
127 u_trajectory[k+1], h_trajectory[k],integrated_state);
129 state_vector_t defect = (y_trajectory[k+1] - integrated_state);
131 _ms_defect_vector.segment<DIMENSIONS::kTotalStates>(current_defect_index) = defect;
133 current_defect_index += DIMENSIONS::kTotalStates;
136 return(_ms_defect_vector);
139 template<
class DIMENSIONS>
140 void ShootingConstraint<DIMENSIONS>::evalDefectGradient(
const Eigen::VectorXd& input, Eigen::MatrixXd &defect_jacobian) {
143 defect_jacobian = this->EvalFxJacobian(input);
146 template<
class DIMENSIONS>
147 void ShootingConstraint<DIMENSIONS>::forwardIntegrateForTimeRange(
148 const state_vector_t & curr_state,
const control_vector_t & curr_control,
149 const control_vector_t & next_control,
const double & total_t, state_vector_t & integrated_state) {
152 integrated_state = curr_state;
154 auto f = [=](
const state_vector_t &x, state_vector_t &dxdt,
const double t) {
156 double t_percentage = t/total_t;
157 control_vector_t control =
158 next_control * t_percentage +
159 curr_control * (1-t_percentage);
161 dxdt = _local_dynamics_->systemDynamics(x, control);
165 boost::numeric::odeint::integrate_const(ms_stepper_, f, integrated_state,
void initialize_num_diff()
Overload this method to define the local pointers.
Definition: GenericConstraintsBase.hpp:151
double ms_dt_
Definition: shooting_constraint.hpp:86
virtual Eigen::VectorXd Evalfx(const Eigen::VectorXd &in_vect) override
To be overloaded by the derived class.
Definition: shooting_constraint.hpp:100
This class implements Direct Multiple Shooting.
Definition: multiple_shooting.hpp:44
Base Class to define vector of constraints and its derivatives.
Definition: GenericConstraintsBase.hpp:44
EIGEN_MAKE_ALIGNED_OPERATOR_NEW GenericConstraintsBase(int num_eq, int num_var)
Generic Constraint Constructor.
Definition: GenericConstraintsBase.hpp:55
Class handling the defects of a multiple shooting method.
Definition: shooting_constraint.hpp:49
Base class for vector constraints and numdiff derivatives.
ShootingConstraint(const int &n_func, const int &n_vars, std::shared_ptr< DTOMultipleShooting< DIMENSIONS >> my_owner)
Constructor of ShootingConstraint.
Definition: shooting_constraint.hpp:66
Class implementing multiple shooting for direct trajectory optimizaion.