DirectTrajectoryOptimization  v0.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
direct_transcription_implementation.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 
30 #include <Eigen/Dense>
31 #include <cstdlib>
32 #include <iostream>
33 
34 #include "direct_transcription.hpp"
35 
36 namespace DirectTrajectoryOptimization {
37 
38 template <class DIMENSIONS>
39 void DTODirectTranscription<DIMENSIONS>::SetDircolMethod(const int & dcm) {
40 
41  dircol_method = dcm;
42 }
43 
44 template <class DIMENSIONS>
45 void DTODirectTranscription<DIMENSIONS>::evaluateDefects(
46  Eigen::VectorXd & defect_vector) {
47 
48  defect_vector.setZero();
49 
50  // ToDo : Move this to the update variables method and validate with : dircol==1
51  state_vector_t y_kmid;
52  state_vector_t F_kmid, ymid_dot;
53  control_vector_t u_kmid;
54 
55  int node = 0;
56  int current_defect_index = 0;
57  int phase_first_node = 0;
58  int inc_index=0;
59 
60  for(int p = 0; p < this->num_phases ; p++ ){
61 
62  for(int n = 0 ; n < this->number_nodes_phase[p] - 1 ; n++) {
63 
64  node = phase_first_node + n;
65 
66  // Trapezodial
67  if (dircol_method == 0) {
68 
69  defect_vector.segment(current_defect_index, this->_numStates) = -this->_y_trajectory[node+1] + this->_y_trajectory[node] +
70  this->_h_trajectory(inc_index)*0.5*( this->_ydot_trajectory[node] + this->_ydot_trajectory[node+1] );
71 
72  } else if(dircol_method == 1) {
73 
74  // Hermitian interpolation of state at intermediate points
75  y_kmid = 0.5*(this->_y_trajectory[node] + this->_y_trajectory[node+1]) + this->_h_trajectory(inc_index)*0.125*(this->_ydot_trajectory[node] - this->_ydot_trajectory[node+1]);
76  // Midpoint rule for control at intermediate points
77  u_kmid = 0.5*(this->_u_trajectory[node] + this->_u_trajectory[node+1]);
78 
79  // approximated dynamics at intermediate point assuming
80  ymid_dot = (1.5/this->_h_trajectory(inc_index))*(this->_y_trajectory[node+1] - this->_y_trajectory[node]) - 0.25*(this->_ydot_trajectory[node] + this->_ydot_trajectory[node+1]);
81 
82  // dynamics at the intermediate points
83  F_kmid = this->_dynamics[this->node_to_phase_map[node]]->systemDynamics(y_kmid, u_kmid);
84 
85  //finally the defect constraints
86  defect_vector.segment(current_defect_index, this->_numStates) = F_kmid - ymid_dot;
87  }
88  current_defect_index += this->_numStates;
89  inc_index++;
90 
91  }
92  phase_first_node += this->number_nodes_phase[p];
93 
94 
95  }
96 }
97 
98 template <class DIMENSIONS>
99 void DTODirectTranscription<DIMENSIONS>::evalSparseJacobianConstraint(double* val) {
100 
101  Eigen::Map<Eigen::VectorXd> val_vec(val,this->myNLP.lenG);
102 
103  int constraint_row = 0; // row of the constraint
104  int G_index = 0; // current index of the sparsity structure
105  int nnz_per_row;
106 
107  // defect constraints
108  // Trapezodial
109  if (dircol_method == 0) {
110  //state_vector_t F_k, F_kp1;
111  state_vector_t Jac_h;
112  state_matrix_t Jac_yk, Jac_ykp1;
113  control_gain_matrix_t Jac_uk, Jac_ukp1;
114 
115  double alpha = 0;
116  int node = 0;
117  int phase_first_node = 0;
118  int inc_index = 0;
119 
120  for(int p = 0; p < this->num_phases ; p++ ){
121 
122  // initialize iteration variables
123  this->_derivatives[this->node_to_phase_map[phase_first_node]]->setCurrentStateAndControl(this->_y_trajectory[phase_first_node],this->_u_trajectory[phase_first_node]);
124  state_matrix_t aux_s_matrix = this->_derivatives[this->node_to_phase_map[phase_first_node]]->getDerivativeState();
125  control_gain_matrix_t aux_c_matrix = this->_derivatives[this->node_to_phase_map[phase_first_node]]->getDerivativeControl();
126 
127 
128  for(int n = 0 ; n < this->number_nodes_phase[p] - 1 ; n++) {
129 
130  node = phase_first_node + n;
131 
132  alpha = 0.5*this->_h_trajectory(inc_index);
133 
134  //F_kp1 = this->_dynamics[this->node_to_phase_map[k+1]]->systemDynamics(y_trajectory[k+1], u_trajectory[k+1]);
135 
136  //w.r.t h (room for efficiency)
137  Jac_h = 0.5*(this->_ydot_trajectory[node]+ this->_ydot_trajectory[node+1]);
138  //w.r.t y_k
139  Jac_yk = I_y + alpha * aux_s_matrix;
140  //w.r.t u_k
141  Jac_uk = alpha * aux_c_matrix;
142 
143  //Very important : Do not move this above previous Jac's.
144  this->_derivatives[this->node_to_phase_map[node+1]]->setCurrentStateAndControl(this->_y_trajectory[node+1],this->_u_trajectory[node+1]);
145  aux_s_matrix = this->_derivatives[this->node_to_phase_map[node+1]]->getDerivativeState();
146  aux_c_matrix = this->_derivatives[this->node_to_phase_map[node+1]]->getDerivativeControl();
147 
148  //w.r.t y_kp1
149  Jac_ykp1 = -I_y + alpha * aux_s_matrix;
150  //w.r.t u_kp1
151  Jac_ukp1 = alpha * aux_c_matrix;
152 
153  for(int jj=0; jj<this->_numStates; jj++) {
154 
155  val_vec(G_index) = Jac_h(jj);
156  G_index += 1;
157  val_vec.segment(G_index, this->_numStates) = Jac_yk.row(jj);
158  G_index += this->_numStates;
159  val_vec.segment(G_index, this->_numStates) = Jac_ykp1.row(jj);
160  G_index += this->_numStates;
161  val_vec.segment(G_index, this->_numControls) = Jac_uk.row(jj);
162  G_index += this->_numControls;
163  val_vec.segment(G_index, this->_numControls) = Jac_ukp1.row(jj);
164  G_index += this->_numControls;
165 
166  constraint_row++;
167  }
168  inc_index++;
169  }
170  phase_first_node += this->number_nodes_phase[p];
171  }
172 
173  } else if (dircol_method == 1) { // Simpson-Hermite
174 
175  state_vector_t y_kmid, F_kmid;
176  control_vector_t u_kmid;
177 
178  state_vector_t Jac_h;
179  state_matrix_t Jac_yk, Jac_ykmid, Jac_ykp1;
180  state_matrix_t dFk_dyk, dFkmid_dykmid, dFkp1_dykp1;
181  control_gain_matrix_t Jac_uk, Jac_ukmid, Jac_ukp1;
182  control_gain_matrix_t dFk_duk, dFkmid_dukmid, dFkp1_dukp1;
183 
184  double fac = 0.0;
185 
186  int node = 0;
187  int phase_first_node = 0;
188  int inc_index = 0;
189 
190  for(int p = 0; p < this->num_phases ; p++ ){
191 
192  for(int n = 0 ; n < this->number_nodes_phase[p] - 1 ; n++) {
193 
194  node = phase_first_node + n;
195 
196  fac = 3/(2*this->_h_trajectory(inc_index));
197 
198  //F_k = this->_dynamics[this->node_to_phase_map[k]]->systemDynamics(y_trajectory[k], u_trajectory[k]);
199  //F_kp1 = this->_dynamics[this->node_to_phase_map[node+1]]->systemDynamics(y_trajectory[node+1], u_trajectory[node+1]);
200 
201  y_kmid = 0.5*(this->_y_trajectory[node] + this->_y_trajectory[node+1]) + this->_h_trajectory(inc_index)*0.125*(this->_ydot_trajectory[node] - this->_ydot_trajectory[node+1]);
202 
203  u_kmid = 0.5*(this->_u_trajectory[node] + this->_u_trajectory[node+1]);
204 
205  F_kmid = this->_dynamics[this->node_to_phase_map[node]]->systemDynamics(y_kmid, u_kmid);
206 
207  // Dynamic derivatives
208  this->_derivatives[this->node_to_phase_map[node]]->setCurrentStateAndControl(this->_y_trajectory[node],this->_u_trajectory[node]);
209  dFk_dyk = this->_derivatives[this->node_to_phase_map[node]]->getDerivativeState();
210  dFk_duk = this->_derivatives[this->node_to_phase_map[node]]->getDerivativeControl();
211 
212  this->_derivatives[this->node_to_phase_map[node]]->setCurrentStateAndControl(y_kmid,u_kmid);
213  dFkmid_dykmid = this->_derivatives[this->node_to_phase_map[node]]->getDerivativeState();
214  dFkmid_dukmid = this->_derivatives[this->node_to_phase_map[node]]->getDerivativeControl();
215 
216  this->_derivatives[this->node_to_phase_map[node+1]]->setCurrentStateAndControl(this->_y_trajectory[node+1],this->_u_trajectory[node+1]);
217  dFkp1_dykp1 = this->_derivatives[this->node_to_phase_map[node+1]]->getDerivativeState();
218  dFkp1_dukp1 = this->_derivatives[this->node_to_phase_map[node+1]]->getDerivativeControl();
219 
220  //Elements of Jacobian
221  //w.r.t h
222  Jac_h = 0.125 * dFkmid_dykmid * (this->_ydot_trajectory[node] - this->_ydot_trajectory[node+1])
223  - (fac/this->_h_trajectory(inc_index))*(this->_y_trajectory[node] - this->_y_trajectory[node+1]); //
224  //w.r.t y_k
225  Jac_yk = this->_h_trajectory(inc_index)*0.125*dFkmid_dykmid*dFk_dyk + 0.5*dFkmid_dykmid + I_y*fac + 0.25*dFk_dyk; //
226  //w.r.t u_k
227  Jac_uk = this->_h_trajectory(inc_index)*0.125*dFkmid_dykmid*dFk_duk + 0.5*dFkmid_dukmid + 0.25*dFk_duk; //
228  //w.r.t y_kp1
229  Jac_ykp1 = -this->_h_trajectory(inc_index)*0.125*dFkmid_dykmid*dFkp1_dykp1 + 0.5*dFkmid_dykmid - fac*I_y + 0.25*dFkp1_dykp1; //
230  //w.r.t u_kp1
231  Jac_ukp1 = -this->_h_trajectory(inc_index)*0.125*dFkmid_dykmid*dFkp1_dukp1 + 0.5*dFkmid_dukmid + 0.25*dFkp1_dukp1; //
232 
233 
234  for (int jj = 0 ; jj < this->_numStates; jj++) {
235  val_vec(G_index) = Jac_h(jj);
236  G_index += 1;
237  val_vec.segment(G_index, this->_numStates) = Jac_yk.row(jj);
238  G_index += this->_numStates;
239  val_vec.segment(G_index, this->_numStates) = Jac_ykp1.row(jj);
240  G_index += this->_numStates;
241  val_vec.segment(G_index, this->_numControls) = Jac_uk.row(jj);
242  G_index += this->_numControls;
243  val_vec.segment(G_index, this->_numControls) = Jac_ukp1.row(jj);
244  G_index += this->_numControls;
245 
246  constraint_row++;
247  }
248  inc_index++;
249  }
250  phase_first_node += this->number_nodes_phase[p];
251  }
252  }
253 
254  if (constraint_row!= this->num_defects) {
255  std::cout << "Error evaluating the Jacobian of constraints (DEFECTS)" << std::endl;
256  exit(EXIT_FAILURE);
257  }
258 
259  // parametric constraints
260  // Time constraint
261  val_vec.segment(G_index, this->n_inc) = Eigen::VectorXd::Ones(this->n_inc);
262  G_index += this->n_inc;
263  constraint_row++;
264 
265  // increments constraint
266  nnz_per_row = 2;
267  if (this->_evenTimeIncr) {
268 
269  for (int k = 0 ; k < this->n_inc -1 ; k++) {
270 
271  Eigen::VectorXd aux_value(nnz_per_row);
272  aux_value << 1 , -1;
273  val_vec.segment(G_index,nnz_per_row) = aux_value;
274 
275  G_index+=nnz_per_row;
276  constraint_row++;
277  }
278  }
279 
280  if (constraint_row != this->num_parametric_constraints + this->num_defects) {
281  std::cout << "Error evaluating the Jacobian of constraints (PARAMETRIC)" << std::endl;
282  std::cout << "contratint_row : " << constraint_row << std::endl;
283  std::cout << "num_parametric_constraints : " << this->num_parametric_constraints << std::endl;
284  exit(EXIT_FAILURE);
285  }
286 
287  // User Constraints
288  for (int k=0; k< this->_numberOfNodes; k++) {
289 
290  val_vec.segment(G_index,this->num_nonzero_jacobian_user_per_node_in_phase[this->node_to_phase_map[k]]) =
291  this->_constraints[this->node_to_phase_map[k]]->evalConstraintsFctDerivatives(this->_y_trajectory[k], this->_u_trajectory[k]);
292  G_index += this->num_nonzero_jacobian_user_per_node_in_phase[this->node_to_phase_map[k]];
293  constraint_row += this->num_user_constraints_per_node_in_phase[this->node_to_phase_map[k]];
294  }
295 
296  if (this->_multiphaseproblem) {
297 
298  state_vector_t gJac;
299  /* Guard Constraints */
300  for(int i = 0; i < this->num_phases -1 ; i++) {
301 
302  // nnz per guard
303  int nnz_per_guard = this->_numStates*this->size_guard[i];
304 
305  // state at the left
306  state_vector_t y_l = this->_y_trajectory[this->leftphase_nodes(i)];
307 
308  // guard->evalJacobian should return a single-row Jacobian
309  val_vec.segment(G_index,nnz_per_guard) = this->_guards[i]->evalConstraintsFctDerivatives(y_l);
310  G_index += nnz_per_guard;
311  constraint_row += this->size_guard[i];
312 
313  }
314 
315  /* Interphase constraints */
316 
317  for(int i = 0 ; i < this->num_phases - 1 ; i++) {
318  int nnz_interphase_constraint = 2 * (this->_numStates + this->_numControls) * this->size_interphase_constraint[i];
319 
320  //ToDo: Set this outside this method!
321  state_vector_t y_l = this->_y_trajectory[this->leftphase_nodes(i)];
322  control_vector_t u_l = this->_u_trajectory[this->leftphase_nodes(i)];
323  state_vector_t y_r = this->_y_trajectory[this->rightphase_nodes(i)];
324  control_vector_t u_r = this->_u_trajectory[this->rightphase_nodes(i)];
325 
326  // _interphase_constraint->getJacobian() should return Jacobian matrix row-wise
327  val_vec.segment(G_index,nnz_interphase_constraint) = this->_interphase_constraint[i]->evalConstraintsFctDerivatives(y_l,y_r,u_l,u_r);
328  G_index += nnz_interphase_constraint;
329  constraint_row += this->size_interphase_constraint[i];
330  }
331 
332  }
333 
334  if (constraint_row!=this->myNLP.num_F) {
335  std::cout << "Error evaluating the Jacobian of constraints (USER)" << std::endl;
336  std::cout << "constraint_row : " << constraint_row << std::endl;
337  std::cout << "num_F : " << this->myNLP.num_F << std::endl;
338  std::cout << "lenG : " << this->myNLP.lenG << std::endl;
339  exit(EXIT_FAILURE);
340  }
341 }
342 
343 } // namespace DirectTrajectoryOptimization
Direct transcription class for trajectory optimization.