30 #include <Eigen/Dense>
36 namespace DirectTrajectoryOptimization {
38 template <
class DIMENSIONS>
39 void DTODirectTranscription<DIMENSIONS>::SetDircolMethod(
const int & dcm) {
44 template <
class DIMENSIONS>
45 void DTODirectTranscription<DIMENSIONS>::evaluateDefects(
46 Eigen::VectorXd & defect_vector) {
48 defect_vector.setZero();
51 state_vector_t y_kmid;
52 state_vector_t F_kmid, ymid_dot;
53 control_vector_t u_kmid;
56 int current_defect_index = 0;
57 int phase_first_node = 0;
60 for(
int p = 0; p < this->num_phases ; p++ ){
62 for(
int n = 0 ; n < this->number_nodes_phase[p] - 1 ; n++) {
64 node = phase_first_node + n;
67 if (dircol_method == 0) {
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] );
72 }
else if(dircol_method == 1) {
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]);
77 u_kmid = 0.5*(this->_u_trajectory[node] + this->_u_trajectory[node+1]);
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]);
83 F_kmid = this->_dynamics[this->node_to_phase_map[node]]->systemDynamics(y_kmid, u_kmid);
86 defect_vector.segment(current_defect_index, this->_numStates) = F_kmid - ymid_dot;
88 current_defect_index += this->_numStates;
92 phase_first_node += this->number_nodes_phase[p];
98 template <
class DIMENSIONS>
99 void DTODirectTranscription<DIMENSIONS>::evalSparseJacobianConstraint(
double* val) {
101 Eigen::Map<Eigen::VectorXd> val_vec(val,this->myNLP.lenG);
103 int constraint_row = 0;
109 if (dircol_method == 0) {
111 state_vector_t Jac_h;
112 state_matrix_t Jac_yk, Jac_ykp1;
113 control_gain_matrix_t Jac_uk, Jac_ukp1;
117 int phase_first_node = 0;
120 for(
int p = 0; p < this->num_phases ; p++ ){
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();
128 for(
int n = 0 ; n < this->number_nodes_phase[p] - 1 ; n++) {
130 node = phase_first_node + n;
132 alpha = 0.5*this->_h_trajectory(inc_index);
137 Jac_h = 0.5*(this->_ydot_trajectory[node]+ this->_ydot_trajectory[node+1]);
139 Jac_yk = I_y + alpha * aux_s_matrix;
141 Jac_uk = alpha * aux_c_matrix;
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();
149 Jac_ykp1 = -I_y + alpha * aux_s_matrix;
151 Jac_ukp1 = alpha * aux_c_matrix;
153 for(
int jj=0; jj<this->_numStates; jj++) {
155 val_vec(G_index) = Jac_h(jj);
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;
170 phase_first_node += this->number_nodes_phase[p];
173 }
else if (dircol_method == 1) {
175 state_vector_t y_kmid, F_kmid;
176 control_vector_t u_kmid;
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;
187 int phase_first_node = 0;
190 for(
int p = 0; p < this->num_phases ; p++ ){
192 for(
int n = 0 ; n < this->number_nodes_phase[p] - 1 ; n++) {
194 node = phase_first_node + n;
196 fac = 3/(2*this->_h_trajectory(inc_index));
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]);
203 u_kmid = 0.5*(this->_u_trajectory[node] + this->_u_trajectory[node+1]);
205 F_kmid = this->_dynamics[this->node_to_phase_map[node]]->systemDynamics(y_kmid, u_kmid);
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();
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();
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();
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]);
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;
227 Jac_uk = this->_h_trajectory(inc_index)*0.125*dFkmid_dykmid*dFk_duk + 0.5*dFkmid_dukmid + 0.25*dFk_duk;
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;
231 Jac_ukp1 = -this->_h_trajectory(inc_index)*0.125*dFkmid_dykmid*dFkp1_dukp1 + 0.5*dFkmid_dukmid + 0.25*dFkp1_dukp1;
234 for (
int jj = 0 ; jj < this->_numStates; jj++) {
235 val_vec(G_index) = Jac_h(jj);
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;
250 phase_first_node += this->number_nodes_phase[p];
254 if (constraint_row!= this->num_defects) {
255 std::cout <<
"Error evaluating the Jacobian of constraints (DEFECTS)" << std::endl;
261 val_vec.segment(G_index, this->n_inc) = Eigen::VectorXd::Ones(this->n_inc);
262 G_index += this->n_inc;
267 if (this->_evenTimeIncr) {
269 for (
int k = 0 ; k < this->n_inc -1 ; k++) {
271 Eigen::VectorXd aux_value(nnz_per_row);
273 val_vec.segment(G_index,nnz_per_row) = aux_value;
275 G_index+=nnz_per_row;
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;
288 for (
int k=0; k< this->_numberOfNodes; k++) {
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]];
296 if (this->_multiphaseproblem) {
300 for(
int i = 0; i < this->num_phases -1 ; i++) {
303 int nnz_per_guard = this->_numStates*this->size_guard[i];
306 state_vector_t y_l = this->_y_trajectory[this->leftphase_nodes(i)];
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];
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];
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)];
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];
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;
Direct transcription class for trajectory optimization.