31 #ifndef INCLUDE_METHOD_INTERFACES_MULTIPLE_SHOOTING_HPP_
32 #define INCLUDE_METHOD_INTERFACES_MULTIPLE_SHOOTING_HPP_
36 #include <Eigen/Dense>
37 #include <boost/numeric/odeint.hpp>
41 namespace DirectTrajectoryOptimization {
43 template <
typename DIMENSIONS>
50 template<
typename DIMENSIONS>
54 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
56 typedef typename DIMENSIONS::state_vector_t state_vector_t;
57 typedef typename DIMENSIONS::control_vector_t control_vector_t;
58 typedef typename DIMENSIONS::state_vector_array_t state_vector_array_t;
59 typedef typename DIMENSIONS::control_vector_array_t control_vector_array_t;
62 std::vector<std::shared_ptr<DynamicsBase<DIMENSIONS> > > dynamics,
63 std::vector<std::shared_ptr<DerivativesBaseDS<DIMENSIONS> > > derivatives,
64 std::shared_ptr<BaseClass::CostFunctionBase<DIMENSIONS> > costFunction,
66 Eigen::Vector2d duration,
int number_of_nodes,
bool methodVerbose) :
68 constraints, duration, number_of_nodes, methodVerbose) { }
72 virtual void Method_specific_initialization()
override {
73 _shootingconstraint_ = std::shared_ptr<ShootingConstraint<DIMENSIONS> >(
new ShootingConstraint<DIMENSIONS>(this->num_defects,this->myNLP.num_vars,
74 std::enable_shared_from_this<DTOMultipleShooting<DIMENSIONS>>::shared_from_this()));
80 virtual void evalSparseJacobianConstraint(
double* val)
override;
81 virtual void evaluateDefects(Eigen::VectorXd & defect_vector)
override;
86 template<
class DIMENSIONS>
89 Eigen::Map<Eigen::VectorXd> val_vec(val,this->myNLP.lenG);
91 Eigen::MatrixXd Jacobian;
93 Eigen::VectorXd x_local;
94 this->SetDecisionVariablesFromStateControlIncrementsTrajectories(x_local,this->_y_trajectory,
98 _shootingconstraint_->evalDefectGradient(x_local,Jacobian);
101 int constraint_row = 0;
107 for (
int k = 0 ; k < this->n_inc ; k++) {
109 int set_of_defects = k * this->_numStates;
111 for (
int jj = 0; jj < this->_numStates; jj++) {
112 int current_state = set_of_defects + jj;
114 val_vec(G_index) = Jacobian(set_of_defects + jj , this->h_inds(k));
118 val_vec.segment(G_index, this->_numStates) = Jacobian.block(current_state , this->y_inds(0,k) , 1 , this->_numStates).transpose();
119 G_index += this->_numStates;
122 val_vec.segment(G_index, this->_numStates) = Jacobian.block(current_state , this->y_inds(0,k+1) , 1 , this->_numStates).transpose();
123 G_index += this->_numStates;
126 val_vec.segment(G_index, this->_numControls) = Jacobian.block(current_state , this->u_inds(0,k) , 1 , this->_numControls).transpose();
127 G_index += this->_numControls;
130 val_vec.segment(G_index, this->_numControls) = Jacobian.block(current_state , this->u_inds(0,k+1) , 1 , this->_numControls).transpose();
131 G_index += this->_numControls;
138 if (constraint_row != this->num_defects) {
139 std::cout <<
"Error evaluating the Jacobian of constraints (DEFECTS)"
149 val_vec.segment(G_index, this->n_inc) = Eigen::VectorXd::Ones(this->n_inc);
150 G_index += this->n_inc;
155 if (this->_evenTimeIncr) {
157 for (
int k = 0 ; k < this->n_inc -1 ; k++) {
160 Eigen::VectorXd aux_value(nnz_per_row);
162 val_vec.segment(G_index,nnz_per_row) = aux_value;
164 G_index+=nnz_per_row;
169 if (constraint_row != this->num_parametric_constraints + this->num_defects) {
170 std::cout <<
"Error evaluating the Jacobian of constraints (PARAMETRIC)" << std::endl;
171 std::cout <<
"contratint_row : " << constraint_row << std::endl;
172 std::cout <<
"num_parametric_constraints : " << this->num_parametric_constraints << std::endl;
178 for (
int k=0; k< this->_numberOfNodes; k++) {
181 val_vec.segment(G_index,this->num_nonzero_jacobian_user_per_node_in_phase[this->node_to_phase_map[k]]) = this->_constraints[this->node_to_phase_map[k]]->evalConstraintsFctDerivatives(this->_y_trajectory[k],this->_u_trajectory[k]);
182 G_index += this->num_nonzero_jacobian_user_per_node_in_phase[this->node_to_phase_map[k]];
184 constraint_row += this->num_user_constraints_per_node_in_phase[this->node_to_phase_map[k]];
188 if (this->_usingGuardFunction) {
191 int num_nnz_jac_guards = this->num_guards * this->_numStates;
192 Eigen::VectorXi guard_derivatives;
193 guard_derivatives.resize(this->_numStates);
194 val_vec.segment(G_index,num_nnz_jac_guards).setZero();
195 G_index +=num_nnz_jac_guards;
196 constraint_row += this->num_guards;
199 if (constraint_row!=this->myNLP.num_F) {
200 std::cout <<
"Error evaluating the Jacobian of constraints (USER)" << std::endl;
201 std::cout <<
"constraint_row : " << constraint_row << std::endl;
202 std::cout <<
"num_F : " << this->myNLP.num_F << std::endl;
203 std::cout <<
"lenG : " << this->myNLP.lenG << std::endl;
208 template<
class DIMENSIONS>
209 void DTOMultipleShooting<DIMENSIONS>::evaluateDefects(
210 Eigen::VectorXd & defect_vector) {
212 Eigen::VectorXd x_local;
213 this->SetDecisionVariablesFromStateControlIncrementsTrajectories(x_local,this->_y_trajectory,this->_u_trajectory,this->_h_trajectory);
215 defect_vector = _shootingconstraint_->Evalfx(x_local);
This class serves as a base interface for inheriting classes.
Definition: base_method_interface.hpp:60
std::shared_ptr< ShootingConstraint< DIMENSIONS > > _shootingconstraint_
Definition: multiple_shooting.hpp:83
Base class for TO constraints.
Definition: ConstraintsBase.hpp:46
This class implements Direct Multiple Shooting.
Definition: multiple_shooting.hpp:44
Multiple shooting defect function including numerical differentiation.
Base class for direct methods.