autogenu-jupyter
An automatic code generator and the continuation/GMRES (C/GMRES) based numerical solvers for nonlinear MPC
Loading...
Searching...
No Matches
ocp.hpp
Go to the documentation of this file.
1// This file was automatically generated by autogenu-jupyter (https://github.com/mayataka/autogenu-jupyter).
2// The autogenu-jupyter copyright holders make no ownership claim of its contents.
3
4#ifndef CGMRES__OCP_MOBILEROBOT_HPP_
5#define CGMRES__OCP_MOBILEROBOT_HPP_
6
7#define _USE_MATH_DEFINES
8
9#include <cmath>
10#include <array>
11#include <iostream>
12
13#include "cgmres/types.hpp"
15
16namespace cgmres {
17
23public:
24
28 static constexpr int nx = 3;
29
33 static constexpr int nu = 2;
34
38 static constexpr int nc = 2;
39
43 static constexpr int nh = 2;
44
48 static constexpr int nuc = nu + nc;
49
53 static constexpr int nub = 2;
54
55 double vx_ref = 0.4;
56 double v_min = -0.5;
57 double v_max = 0.5;
58 double w_min = -0.75;
59 double w_max = 0.75;
60 double X_1 = 1;
61 double Y_1 = 0.25;
62 double R_1 = 0.5;
63 double X_2 = 2;
64 double Y_2 = -0.25;
65 double R_2 = 0.5;
66
67 std::array<double, 3> q = {10, 1, 0.01};
68 std::array<double, 2> r = {0.1, 0.1};
69 std::array<double, 3> x_ref = {0, 0, 0};
70
71 static constexpr std::array<int, nub> ubound_indices = {0, 1};
72 std::array<double, nub> umin = {-1.0, -1.0};
73 std::array<double, nub> umax = {1.0, 1.0};
74 std::array<double, nub> dummy_weight = {0.1, 0.1};
75
76 std::array<double, nh> fb_eps = {0.01, 0.01};
77
78 void disp(std::ostream& os) const {
79 os << "OCP_mobilerobot:" << std::endl;
80 os << " nx: " << nx << std::endl;
81 os << " nu: " << nu << std::endl;
82 os << " nc: " << nc << std::endl;
83 os << " nh: " << nh << std::endl;
84 os << " nuc: " << nuc << std::endl;
85 os << " nub: " << nub << std::endl;
86 os << std::endl;
87 os << " vx_ref: " << vx_ref << std::endl;
88 os << " v_min: " << v_min << std::endl;
89 os << " v_max: " << v_max << std::endl;
90 os << " w_min: " << w_min << std::endl;
91 os << " w_max: " << w_max << std::endl;
92 os << " X_1: " << X_1 << std::endl;
93 os << " Y_1: " << Y_1 << std::endl;
94 os << " R_1: " << R_1 << std::endl;
95 os << " X_2: " << X_2 << std::endl;
96 os << " Y_2: " << Y_2 << std::endl;
97 os << " R_2: " << R_2 << std::endl;
98 os << std::endl;
99 Eigen::IOFormat fmt(4, 0, ", ", "", "[", "]");
100 Eigen::IOFormat intfmt(1, 0, ", ", "", "[", "]");
101 os << " q: " << Map<const VectorX>(q.data(), q.size()).transpose().format(fmt) << std::endl;
102 os << " r: " << Map<const VectorX>(r.data(), r.size()).transpose().format(fmt) << std::endl;
103 os << " x_ref: " << Map<const VectorX>(x_ref.data(), x_ref.size()).transpose().format(fmt) << std::endl;
104 os << std::endl;
105 os << " ubound_indices: " << Map<const VectorXi>(ubound_indices.data(), ubound_indices.size()).transpose().format(intfmt) << std::endl;
106 os << " umin: " << Map<const VectorX>(umin.data(), umin.size()).transpose().format(fmt) << std::endl;
107 os << " umax: " << Map<const VectorX>(umax.data(), umax.size()).transpose().format(fmt) << std::endl;
108 os << " dummy_weight: " << Map<const VectorX>(dummy_weight.data(), dummy_weight.size()).transpose().format(fmt) << std::endl;
109 os << std::endl;
110 os << " fb_eps: " << Map<const VectorX>(fb_eps.data(), fb_eps.size()).transpose().format(fmt) << std::endl;
111 }
112
113 friend std::ostream& operator<<(std::ostream& os, const OCP_mobilerobot& ocp) {
114 ocp.disp(os);
115 return os;
116 }
117
118
123 void synchronize() {
124 }
125
135 void eval_f(const double t, const double* x, const double* u,
136 double* dx) const {
137 dx[0] = u[0]*cos(x[2]);
138 dx[1] = u[0]*sin(x[2]);
139 dx[2] = u[1];
140
141 }
142
152 void eval_phix(const double t, const double* x, double* phix) const {
153 phix[0] = (1.0/2.0)*q[0]*(-2*t*vx_ref + 2*x[0]);
154 phix[1] = q[1]*x[1];
155 phix[2] = q[2]*x[2];
156
157 }
158
170 void eval_hx(const double t, const double* x, const double* u,
171 const double* lmd, double* hx) const {
172 const double x0 = -2*x[0];
173 const double x1 = -2*x[1];
174 const double x2 = u[0]*sin(x[2]);
175 const double x3 = cos(x[2]);
176 hx[0] = (1.0/2.0)*q[0]*(-2*t*vx_ref - x0) + u[2]*(2*X_1 + x0) + u[3]*(2*X_2 + x0);
177 hx[1] = q[1]*x[1] + u[2]*(2*Y_1 + x1) + u[3]*(2*Y_2 + x1);
178 hx[2] = -lmd[0]*x2 + lmd[1]*u[0]*x3 + q[2]*x[2] - r[0]*x2*(u[0]*x3 - vx_ref);
179
180 }
181
193 void eval_hu(const double t, const double* x, const double* u,
194 const double* lmd, double* hu) const {
195 const double x0 = cos(x[2]);
196 const double x1 = -x[0];
197 const double x2 = -x[1];
198 const double x3 = -pow(R_1, 2) + pow(-X_1 - x1, 2) + pow(-Y_1 - x2, 2);
199 const double x4 = -pow(R_2, 2) + pow(-X_2 - x1, 2) + pow(-Y_2 - x2, 2);
200 hu[0] = lmd[0]*x0 + lmd[1]*sin(x[2]) + r[0]*x0*(u[0]*x0 - vx_ref);
201 hu[1] = lmd[2] + r[1]*u[1];
202 hu[2] = -u[2] - x3 + sqrt(fb_eps[0] + pow(u[2], 2) + pow(x3, 2));
203 hu[3] = -u[3] - x4 + sqrt(fb_eps[1] + pow(u[3], 2) + pow(x4, 2));
204
205 }
206
214 template <typename VectorType1, typename VectorType2, typename VectorType3>
215 void eval_f(const double t, const MatrixBase<VectorType1>& x,
216 const MatrixBase<VectorType2>& u,
217 const MatrixBase<VectorType3>& dx) const {
218 if (x.size() != nx) {
219 throw std::invalid_argument("[OCP]: x.size() must be " + std::to_string(nx));
220 }
221 if (u.size() != nu) {
222 throw std::invalid_argument("[OCP]: u.size() must be " + std::to_string(nu));
223 }
224 if (dx.size() != nx) {
225 throw std::invalid_argument("[OCP]: dx.size() must be " + std::to_string(nx));
226 }
227 eval_f(t, x.derived().data(), u.derived().data(), CGMRES_EIGEN_CONST_CAST(VectorType3, dx).data());
228 }
229
237 template <typename VectorType1, typename VectorType2>
238 void eval_phix(const double t, const MatrixBase<VectorType1>& x,
239 const MatrixBase<VectorType2>& phix) const {
240 if (x.size() != nx) {
241 throw std::invalid_argument("[OCP]: x.size() must be " + std::to_string(nx));
242 }
243 if (phix.size() != nx) {
244 throw std::invalid_argument("[OCP]: phix.size() must be " + std::to_string(nx));
245 }
246 eval_phix(t, x.derived().data(), CGMRES_EIGEN_CONST_CAST(VectorType2, phix).data());
247 }
248
258 template <typename VectorType1, typename VectorType2, typename VectorType3, typename VectorType4>
259 void eval_hx(const double t, const MatrixBase<VectorType1>& x,
260 const MatrixBase<VectorType2>& uc,
261 const MatrixBase<VectorType3>& lmd,
262 const MatrixBase<VectorType4>& hx) const {
263 if (x.size() != nx) {
264 throw std::invalid_argument("[OCP]: x.size() must be " + std::to_string(nx));
265 }
266 if (uc.size() != nuc) {
267 throw std::invalid_argument("[OCP]: uc.size() must be " + std::to_string(nuc));
268 }
269 if (lmd.size() != nx) {
270 throw std::invalid_argument("[OCP]: lmd.size() must be " + std::to_string(nx));
271 }
272 if (hx.size() != nuc) {
273 throw std::invalid_argument("[OCP]: hx.size() must be " + std::to_string(nx));
274 }
275 eval_hx(t, x.derived().data(), uc.derived().data(), lmd.derived().data(), CGMRES_EIGEN_CONST_CAST(VectorType4, hx).data());
276 }
277
287 template <typename VectorType1, typename VectorType2, typename VectorType3, typename VectorType4>
288 void eval_hu(const double t, const MatrixBase<VectorType1>& x,
289 const MatrixBase<VectorType2>& uc,
290 const MatrixBase<VectorType3>& lmd,
291 const MatrixBase<VectorType4>& hu) const {
292 if (x.size() != nx) {
293 throw std::invalid_argument("[OCP]: x.size() must be " + std::to_string(nx));
294 }
295 if (uc.size() != nuc) {
296 throw std::invalid_argument("[OCP]: uc.size() must be " + std::to_string(nuc));
297 }
298 if (lmd.size() != nx) {
299 throw std::invalid_argument("[OCP]: lmd.size() must be " + std::to_string(nx));
300 }
301 if (hu.size() != nuc) {
302 throw std::invalid_argument("[OCP]: hu.size() must be " + std::to_string(nuc));
303 }
304 eval_hu(t, x.derived().data(), uc.derived().data(), lmd.derived().data(), CGMRES_EIGEN_CONST_CAST(VectorType4, hu).data());
305 }
306
307};
308
309} // namespace cgmres
310
311#endif // CGMRES_OCP_HPP_
Definition of the optimal control problem (OCP) of mobilerobot.
Definition: ocp.hpp:22
double w_min
Definition: ocp.hpp:58
double X_1
Definition: ocp.hpp:60
double X_2
Definition: ocp.hpp:63
void disp(std::ostream &os) const
Definition: ocp.hpp:78
std::array< double, nh > fb_eps
Definition: ocp.hpp:76
static constexpr int nx
Dimension of the state.
Definition: ocp.hpp:28
double Y_1
Definition: ocp.hpp:61
void eval_f(const double t, const MatrixBase< VectorType1 > &x, const MatrixBase< VectorType2 > &u, const MatrixBase< VectorType3 > &dx) const
Computes the state equation dx = f(t, x, u).
Definition: ocp.hpp:215
std::array< double, nub > umin
Definition: ocp.hpp:72
static constexpr int nuc
Dimension of the concatenation of the control input and equality constraints.
Definition: ocp.hpp:48
double R_2
Definition: ocp.hpp:65
friend std::ostream & operator<<(std::ostream &os, const OCP_mobilerobot &ocp)
Definition: ocp.hpp:113
void eval_phix(const double t, const double *x, double *phix) const
Computes the partial derivative of terminal cost with respect to state, i.e., phix = dphi/dx(t,...
Definition: ocp.hpp:152
double vx_ref
Definition: ocp.hpp:55
double w_max
Definition: ocp.hpp:59
static constexpr std::array< int, nub > ubound_indices
Definition: ocp.hpp:71
static constexpr int nc
Dimension of the equality constraints.
Definition: ocp.hpp:38
std::array< double, nub > dummy_weight
Definition: ocp.hpp:74
void eval_hu(const double t, const MatrixBase< VectorType1 > &x, const MatrixBase< VectorType2 > &uc, const MatrixBase< VectorType3 > &lmd, const MatrixBase< VectorType4 > &hu) const
Computes the partial derivative of the Hamiltonian with respect to control input and the equality con...
Definition: ocp.hpp:288
double R_1
Definition: ocp.hpp:62
void eval_f(const double t, const double *x, const double *u, double *dx) const
Computes the state equation dx = f(t, x, u).
Definition: ocp.hpp:135
double v_max
Definition: ocp.hpp:57
double v_min
Definition: ocp.hpp:56
std::array< double, 3 > x_ref
Definition: ocp.hpp:69
std::array< double, 2 > r
Definition: ocp.hpp:68
void eval_hu(const double t, const double *x, const double *u, const double *lmd, double *hu) const
Computes the partial derivative of the Hamiltonian with respect to control input and the equality con...
Definition: ocp.hpp:193
double Y_2
Definition: ocp.hpp:64
void synchronize()
Synchrozies the internal parameters of this OCP with the external references. This method is called a...
Definition: ocp.hpp:123
void eval_hx(const double t, const MatrixBase< VectorType1 > &x, const MatrixBase< VectorType2 > &uc, const MatrixBase< VectorType3 > &lmd, const MatrixBase< VectorType4 > &hx) const
Computes the partial derivative of the Hamiltonian with respect to the state, i.e....
Definition: ocp.hpp:259
static constexpr int nu
Dimension of the control input.
Definition: ocp.hpp:33
static constexpr int nh
Dimension of the Fischer-Burmeister function (already counded in nc).
Definition: ocp.hpp:43
static constexpr int nub
Dimension of the bound constraints on the control input.
Definition: ocp.hpp:53
void eval_phix(const double t, const MatrixBase< VectorType1 > &x, const MatrixBase< VectorType2 > &phix) const
Computes the partial derivative of terminal cost with respect to state, i.e., phix = dphi/dx(t,...
Definition: ocp.hpp:238
std::array< double, nub > umax
Definition: ocp.hpp:73
void eval_hx(const double t, const double *x, const double *u, const double *lmd, double *hx) const
Computes the partial derivative of the Hamiltonian with respect to state, i.e., hx = dH/dx(t,...
Definition: ocp.hpp:170
std::array< double, 3 > q
Definition: ocp.hpp:67
#define CGMRES_EIGEN_CONST_CAST(TYPE, OBJ)
Definition: macros.hpp:7
Definition: continuation_gmres.hpp:11
Eigen::Map< MatrixType > Map
Alias of Eigen::Map.
Definition: types.hpp:50
Eigen::MatrixBase< MatrixType > MatrixBase
Alias of Eigen::MatrixBase.
Definition: types.hpp:29