DirectTrajectoryOptimization  v0.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
multiple_shooting.hpp
Go to the documentation of this file.
1 /***********************************************************************************
2 Copyright (c) 2017, Diego Pardo. All rights reserved.
3 
4 Redistribution and use in source and binary forms, with or without modification,
5 are permitted provided that the following conditions are met:
6  * Redistributions of source code must retain the above copyright notice,
7  this list of conditions and the following disclaimer.
8  * Redistributions in binary form must reproduce the above copyright notice,
9  this list of conditions and the following disclaimer in the documentation
10  and/or other materials provided with the distribution.
11  * Neither the name of ETH ZURICH nor the names of its contributors may be used
12  to endorse or promote products derived from this software without specific
13  prior written permission.
14 
15 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY
16 EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES
17 OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT
18 SHALL ETH ZURICH BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY,
19 OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
20 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
21 CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
22 (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE,
23 EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
24 ***************************************************************************************/
25 
31 #ifndef INCLUDE_METHOD_INTERFACES_MULTIPLE_SHOOTING_HPP_
32 #define INCLUDE_METHOD_INTERFACES_MULTIPLE_SHOOTING_HPP_
33 
34 #define infy 1e20
35 
36 #include <Eigen/Dense>
37 #include <boost/numeric/odeint.hpp>
40 
41 namespace DirectTrajectoryOptimization {
42 
43 template <typename DIMENSIONS>
45 
50 template<typename DIMENSIONS>
51 class DTOMultipleShooting: public DTOMethodInterface<DIMENSIONS>, public std::enable_shared_from_this<DTOMultipleShooting<DIMENSIONS>>{
52 
53 public:
54  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
55 
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;
60 
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,
65  std::vector<std::shared_ptr<BaseClass::ConstraintsBase<DIMENSIONS> > > constraints,
66  Eigen::Vector2d duration, int number_of_nodes, bool methodVerbose) :
67  DTOMethodInterface<DIMENSIONS>(dynamics, derivatives, costFunction,
68  constraints, duration, number_of_nodes, methodVerbose) { }
69 
70  virtual ~DTOMultipleShooting() {}
71 
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()));
75  _shootingconstraint_->initialize();
76 
77  }
78 
79  // Sparsity structure is the same as the one in direct transcription
80  virtual void evalSparseJacobianConstraint(double* val) override;
81  virtual void evaluateDefects(Eigen::VectorXd & defect_vector) override;
82 
83  std::shared_ptr<ShootingConstraint<DIMENSIONS> > _shootingconstraint_;
84 };
85 
86 template<class DIMENSIONS>
88 
89  Eigen::Map<Eigen::VectorXd> val_vec(val,this->myNLP.lenG);
90  val_vec.setZero();
91  Eigen::MatrixXd Jacobian;
92 
93  Eigen::VectorXd x_local;
94  this->SetDecisionVariablesFromStateControlIncrementsTrajectories(x_local,this->_y_trajectory,
95  this->_u_trajectory,
96  this->_h_trajectory);
97 
98  _shootingconstraint_->evalDefectGradient(x_local,Jacobian);
99 
100 
101  int constraint_row = 0; // row of the constraint
102  int G_index = 0; // current index of the sparsity structure
103  int nnz_per_row;
104 
105  /* defect constraints */
106 
107  for (int k = 0 ; k < this->n_inc ; k++) {
108 
109  int set_of_defects = k * this->_numStates;
110  // for each defect
111  for (int jj = 0; jj < this->_numStates; jj++) {
112  int current_state = set_of_defects + jj;
113  // w.r.t. h_k
114  val_vec(G_index) = Jacobian(set_of_defects + jj , this->h_inds(k));//Jac_h(jj);
115  G_index += 1;
116 
117  // w.r.t. y_k
118  val_vec.segment(G_index, this->_numStates) = Jacobian.block(current_state , this->y_inds(0,k) , 1 , this->_numStates).transpose();//Jac_yk.row(jj);
119  G_index += this->_numStates;
120 
121  // w.r.t. y_k+1
122  val_vec.segment(G_index, this->_numStates) = Jacobian.block(current_state , this->y_inds(0,k+1) , 1 , this->_numStates).transpose();//Jac_ykp1.row(jj);
123  G_index += this->_numStates;
124 
125  // w.r.t. u_k
126  val_vec.segment(G_index, this->_numControls) = Jacobian.block(current_state , this->u_inds(0,k) , 1 , this->_numControls).transpose();//Jac_uk.row(jj);
127  G_index += this->_numControls;
128 
129  // w.r.t. u_k+1
130  val_vec.segment(G_index, this->_numControls) = Jacobian.block(current_state , this->u_inds(0,k+1) , 1 , this->_numControls).transpose();//row(jj);
131  G_index += this->_numControls;
132 
133  constraint_row++;
134  }
135  }
136 
137 
138  if (constraint_row != this->num_defects) {
139  std::cout << "Error evaluating the Jacobian of constraints (DEFECTS)"
140  << std::endl;
141  exit(EXIT_FAILURE);
142  }
143 
144  // ToDo: The following constraint jacobian is common to direct_transcription, but here is replicated.
145  // This should be centralized in base_method_interface
146 
147  /* parametric constraints */
148  // Time constraint
149  val_vec.segment(G_index, this->n_inc) = Eigen::VectorXd::Ones(this->n_inc);
150  G_index += this->n_inc;
151  constraint_row++;
152 
153  // increments constraint
154  nnz_per_row = 2;
155  if (this->_evenTimeIncr) {
156 
157  for (int k = 0 ; k < this->n_inc -1 ; k++) {
158  //iGfun.segment(G_index,nnz_per_row) = Eigen::VectorXi::Constant(nnz_per_row,constraint_row);
159 
160  Eigen::VectorXd aux_value(nnz_per_row);
161  aux_value << 1 , -1;
162  val_vec.segment(G_index,nnz_per_row) = aux_value;
163 
164  G_index+=nnz_per_row;
165  constraint_row++;
166  }
167  }
168 
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;
173  exit(EXIT_FAILURE);
174  }
175 
176  /* User Constraints */
177 
178  for (int k=0; k< this->_numberOfNodes; k++) {
179 
180  //val_vec.segment(G_index,this->num_nonzero_jacobian_user_per_node) = this->_constraints[this->node_to_phase_map[k]]->evalConstraintsFctDerivatives(y_trajectory[k], u_trajectory[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]];
183 
184  constraint_row += this->num_user_constraints_per_node_in_phase[this->node_to_phase_map[k]];
185  }
186 
187  // ToDo Guard Constraints
188  if (this->_usingGuardFunction) {
189  //we have this number already!
190  //ToDo: Move this!
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;
197  }
198 
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;
204  exit(EXIT_FAILURE);
205  }
206 }
207 
208 template<class DIMENSIONS>
209 void DTOMultipleShooting<DIMENSIONS>::evaluateDefects(
210  Eigen::VectorXd & defect_vector) {
211 
212  Eigen::VectorXd x_local;
213  this->SetDecisionVariablesFromStateControlIncrementsTrajectories(x_local,this->_y_trajectory,this->_u_trajectory,this->_h_trajectory);
214 
215  defect_vector = _shootingconstraint_->Evalfx(x_local);
216 }
217 
218 } // namespace DirectTrajectoryOptimization
219 
220 #endif /* INCLUDE_METHOD_INTERFACES_MULTIPLE_SHOOTING_HPP_ */
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.