proxsuite-nlp  0.10.0
A primal-dual augmented Lagrangian-type solver for nonlinear programming on manifolds.
Loading...
Searching...
No Matches
expose-manifold.cpp
Go to the documentation of this file.
1#ifdef PROXSUITE_NLP_WITH_PINOCCHIO
2#include <pinocchio/fwd.hpp>
5#endif
6
12
13#include <eigenpy/std-vector.hpp>
14
15#include <vector>
16
17namespace proxsuite {
18namespace nlp {
19namespace python {
20using context::ConstVectorRef;
22using context::MatrixRef;
23using context::Scalar;
24using context::VectorRef;
27
29 using context::MatrixXs;
30 using context::VectorXs;
32
33 using BinaryFunTypeRet = VectorXs (Manifold::*)(const ConstVectorRef &,
34 const ConstVectorRef &) const;
35 using BinaryFunType = void (Manifold::*)(
36 const ConstVectorRef &, const ConstVectorRef &, VectorRef) const;
37 using JacobianFunType = void (Manifold::*)(
38 const ConstVectorRef &, const ConstVectorRef &, MatrixRef, int) const;
39
40 bp::class_<Manifold, boost::noncopyable>(
41 "ManifoldAbstract", "Manifold abstract class.", bp::no_init)
42 .add_property("nx", &Manifold::nx, "Manifold representation dimension.")
43 .add_property("ndx", &Manifold::ndx, "Tangent space dimension.")
44 .def("neutral", &Manifold::neutral, "self"_a,
45 "Get the neutral point from the manifold (if a Lie group).")
46 .def("rand", &Manifold::rand, "self"_a,
47 "Sample a random point from the manifold.")
48 .def("isNormalized", &Manifold::isNormalized, ("self"_a, "x"),
49 "Check if the input vector :math:`x` is a viable element of the "
50 "manifold.")
51 .def<BinaryFunType>("integrate", &Manifold::integrate,
52 ("self"_a, "x", "v", "out"))
53 .def<BinaryFunType>("difference", &Manifold::difference,
54 ("self"_a, "x0", "x1", "out"))
55 .def<BinaryFunTypeRet>("integrate", &Manifold::integrate,
56 ("self"_a, "x", "v"))
57 .def<BinaryFunTypeRet>("difference", &Manifold::difference,
58 ("self"_a, "x0", "x1"))
59 .def("interpolate",
60 (void(Manifold::*)(const ConstVectorRef &, const ConstVectorRef &,
61 const Scalar &, VectorRef)
62 const)(&Manifold::interpolate),
63 ("self"_a, "x0", "x1", "u", "out"))
64 .def(
65 "interpolate",
66 (VectorXs(Manifold::*)(const ConstVectorRef &, const ConstVectorRef &,
67 const Scalar &) const)(&Manifold::interpolate),
68 ("self"_a, "x0", "x1", "u"),
69 "Interpolate between two points on the manifold. Allocated version.")
70 .def<JacobianFunType>("Jintegrate", &Manifold::Jintegrate,
71 ("self"_a, "x", "v", "Jout", "arg"),
72 "Compute the Jacobian of the exp operator.")
73 .def<JacobianFunType>("Jdifference", &Manifold::Jdifference,
74 ("self"_a, "x0", "x1", "Jout", "arg"),
75 "Compute the Jacobian of the log operator.")
76 .def(
77 "Jintegrate",
78 +[](const Manifold &m, const ConstVectorRef x,
79 const ConstVectorRef &v, int arg) {
80 MatrixXs Jout(m.ndx(), m.ndx());
81 m.Jintegrate(x, v, Jout, arg);
82 return Jout;
83 },
84 ("self"_a, "x", "v", "arg"),
85 "Compute and return the Jacobian of the exp.")
86 .def("JintegrateTransport", &Manifold::JintegrateTransport,
87 ("self"_a, "x", "v", "J", "arg"),
88 "Perform parallel transport of matrix J expressed at point x+v to "
89 "point x.")
90 .def(
91 "Jdifference",
92 +[](const Manifold &m, const ConstVectorRef x0,
93 const ConstVectorRef &x1, int arg) {
94 MatrixXs Jout(m.ndx(), m.ndx());
95 m.Jdifference(x0, x1, Jout, arg);
96 return Jout;
97 },
98 ("self"_a, "x0", "x1", "arg"),
99 "Compute and return the Jacobian of the log.")
100 .def("tangent_space", &Manifold::tangentSpace, bp::args("self"),
101 "Returns an object representing the tangent space to this manifold.")
102 .def(
103 "__mul__",
104 +[](const PolyManifold &a, const PolyManifold &b) { return a * b; })
105 .def(
106 "__mul__", +[](const PolyManifold &a,
107 const CartesianProduct &b) { return a * b; })
108 .def(
109 "__rmul__", +[](const PolyManifold &a, const CartesianProduct &b) {
110 return a * b;
111 });
112
113 eigenpy::StdVectorPythonVisitor<std::vector<PolyManifold>>::expose(
114 "StdVec_Manifold",
115 eigenpy::details::overload_base_get_item_for_std_vector<
116 std::vector<PolyManifold>>());
117}
118
120template <typename M>
121bp::class_<TangentBundleTpl<M>, bp::bases<Manifold>>
122exposeTangentBundle(const char *name, const char *docstring) {
123 using OutType = TangentBundleTpl<M>;
124 return bp::class_<OutType, bp::bases<Manifold>>(
125 name, docstring, bp::init<M>(("self"_a, "base")))
126 .add_property("base",
127 bp::make_function(&OutType::getBaseSpace,
128 bp::return_internal_reference<>()),
129 "Get the base space.")
131}
132
134template <typename M, class Init>
135bp::class_<TangentBundleTpl<M>, bp::bases<Manifold>>
136exposeTangentBundle(const char *name, const char *docstring, Init init) {
137 return exposeTangentBundle<M>(name, docstring).def(init);
138}
139
141
142#ifdef PROXSUITE_NLP_WITH_PINOCCHIO
145template <typename LieGroup>
146void exposeLieGroup(const char *name, const char *docstring) {
147 bp::class_<PinocchioLieGroup<LieGroup>, bp::bases<Manifold>> cl(
148 name, docstring, bp::init<>("self"_a));
149 cl.def(PolymorphicVisitor<PolyManifold>());
150 cl.enable_pickling_(true);
151}
152
153void exposePinocchioSpaces() {
154 namespace pin = pinocchio;
155 using pin::ModelTpl;
156 using pin::SpecialEuclideanOperationTpl;
157 using pin::SpecialOrthogonalOperationTpl;
158 using pin::VectorSpaceOperationTpl;
159
160 using DynSizeEuclideanSpace = VectorSpaceOperationTpl<Eigen::Dynamic, Scalar>;
161 bp::class_<PinocchioLieGroup<DynSizeEuclideanSpace>, bp::bases<Manifold>>(
162 "EuclideanSpace", "Pinocchio's n-dimensional Euclidean vector space.",
163 bp::no_init)
164 .def(bp::init<int>(("self"_a, "dim")))
165 .def(PolymorphicVisitor<PolyManifold>());
166
167 exposeLieGroup<VectorSpaceOperationTpl<1, Scalar>>(
168 "R", "One-dimensional Euclidean space AKA real number line.");
169 exposeLieGroup<VectorSpaceOperationTpl<2, Scalar>>(
170 "R2", "Two-dimensional Euclidean space.");
171 exposeLieGroup<VectorSpaceOperationTpl<3, Scalar>>(
172 "R3", "Three-dimensional Euclidean space.");
173 exposeLieGroup<VectorSpaceOperationTpl<4, Scalar>>(
174 "R4", "Four-dimensional Euclidean space.");
175 exposeLieGroup<SpecialOrthogonalOperationTpl<2, Scalar>>(
176 "SO2", "SO(2) special orthogonal group.");
177 exposeLieGroup<SpecialOrthogonalOperationTpl<3, Scalar>>(
178 "SO3", "SO(3) special orthogonal group.");
179 exposeLieGroup<SpecialEuclideanOperationTpl<2, Scalar>>(
180 "SE2", "SE(2) special Euclidean group.");
181 exposeLieGroup<SpecialEuclideanOperationTpl<3, Scalar>>(
182 "SE3", "SE(3) special Euclidean group.");
183
184 using SO2 = PinocchioLieGroup<SpecialOrthogonalOperationTpl<2, Scalar>>;
185 using SO3 = PinocchioLieGroup<SpecialOrthogonalOperationTpl<3, Scalar>>;
186 using SE2 = PinocchioLieGroup<SpecialEuclideanOperationTpl<2, Scalar>>;
187 using SE3 = PinocchioLieGroup<SpecialEuclideanOperationTpl<3, Scalar>>;
188
189 /* Expose tangent bundles */
190
191 exposeTangentBundle<SO2>("TSO2", "Tangent bundle of the SO(2) group.")
192 .def(bp::init<>(bp::args("self")));
193 exposeTangentBundle<SO3>("TSO3", "Tangent bundle of the SO(3) group.")
194 .def(bp::init<>(bp::args("self")));
195 exposeTangentBundle<SE2>("TSE2", "Tangent bundle of the SE(2) group.")
196 .def(bp::init<>(bp::args("self")));
197 exposeTangentBundle<SE3>("TSE3", "Tangent bundle of the SE(3) group.")
198 .def(bp::init<>(bp::args("self")));
199
200 /* Groups associated w/ Pinocchio models */
201 using Multibody = MultibodyConfiguration<Scalar>;
202 using Model = ModelTpl<Scalar>;
203 bp::class_<Multibody, bp::bases<Manifold>>(
204 "MultibodyConfiguration", "Configuration group of a multibody",
205 bp::init<const Model &>(bp::args("self", "model")))
206 .add_property("model",
207 bp::make_function(&Multibody::getModel,
208 bp::return_internal_reference<>()),
209 "Return the Pinocchio model instance.")
210 .def(PolymorphicVisitor<PolyManifold>());
211
212 using MultiPhase = MultibodyPhaseSpace<Scalar>;
213 bp::class_<MultiPhase, bp::bases<Manifold>>(
214 "MultibodyPhaseSpace",
215 "Tangent space of the multibody configuration group.",
216 bp::init<const Model &>(("self"_a, "model")))
217 .add_property("model",
218 bp::make_function(&MultiPhase::getModel,
219 bp::return_internal_reference<>()),
220 "Return the Pinocchio model instance.")
221 .add_property("base", bp::make_function(
222 +[](const MultiPhase &m) -> const Multibody & {
223 return m.getBaseSpace();
224 },
225 bp::return_internal_reference<>()))
226 .def(PolymorphicVisitor<PolyManifold>());
227}
228#endif
229
231
233
234 /* Basic vector space */
235 bp::class_<VectorSpaceTpl<Scalar>, bp::bases<Manifold>>(
236 "VectorSpace", "Basic Euclidean vector space.", bp::no_init)
237 .def(bp::init<const int>(("self"_a, "dim")))
239
241#ifdef PROXSUITE_NLP_WITH_PINOCCHIO
242 exposePinocchioSpaces();
243#endif
244}
245
246} // namespace python
247} // namespace nlp
248} // namespace proxsuite
void register_polymorphic_to_python()
Expose a polymorphic value type, e.g. xyz::polymorphic<T, A>.
bp::class_< TangentBundleTpl< M >, bp::bases< Manifold > > exposeTangentBundle(const char *name, const char *docstring)
Expose the tangent bundle of a manifold type M.
Main package namespace.
Definition bcl-params.hpp:5
The cartesian product of two or more manifolds.
Tangent bundle of a base manifold M.