aligator  0.14.0
A primal-dual augmented Lagrangian-type solver for nonlinear trajectory optimization.
 
Loading...
Searching...
No Matches
expose-pinocchio-manifolds.cpp
Go to the documentation of this file.
3
6
7namespace aligator::python {
10using PolyManifold = xyz::polymorphic<Manifold>;
11
14template <typename LieGroup>
15void exposeLieGroup(const char *name, const char *docstring) {
16 bp::class_<PinocchioLieGroup<LieGroup>, bp::bases<Manifold>>(
17 name, docstring, bp::init<>("self"_a))
19 .enable_pickling_(true);
20}
21
23template <typename M>
24bp::class_<TangentBundleTpl<M>, bp::bases<Manifold>>
25exposeTangentBundle(const char *name, const char *docstring) {
26 using OutType = TangentBundleTpl<M>;
27 return bp::class_<OutType, bp::bases<Manifold>>(
28 name, docstring, bp::init<M>(("self"_a, "base")))
29 .add_property("base",
30 bp::make_function(&OutType::getBaseSpace,
31 bp::return_internal_reference<>()),
32 "Get the base space.")
34}
35
37template <typename M, class Init>
38bp::class_<TangentBundleTpl<M>, bp::bases<Manifold>>
39exposeTangentBundle(const char *name, const char *docstring, Init init) {
40 return exposeTangentBundle<M>(name, docstring).def(init);
41}
42
44 bp::scope scope = get_namespace("manifolds");
45
46 namespace pin = pinocchio;
47 using pin::ModelTpl;
48 using pin::SpecialEuclideanOperationTpl;
49 using pin::SpecialOrthogonalOperationTpl;
50 using pin::VectorSpaceOperationTpl;
51
52 using DynSizeEuclideanSpace = VectorSpaceOperationTpl<Eigen::Dynamic, Scalar>;
53 bp::class_<PinocchioLieGroup<DynSizeEuclideanSpace>, bp::bases<Manifold>>(
54 "EuclideanSpace", "Pinocchio's n-dimensional Euclidean vector space.",
55 bp::no_init)
56 .def(bp::init<int>(("self"_a, "dim")))
58
60 "R", "One-dimensional Euclidean space AKA real number line.");
62 "R2", "Two-dimensional Euclidean space.");
64 "R3", "Three-dimensional Euclidean space.");
66 "R4", "Four-dimensional Euclidean space.");
68 "SO2", "SO(2) special orthogonal group.");
70 "SO3", "SO(3) special orthogonal group.");
72 "SE2", "SE(2) special Euclidean group.");
74 "SE3", "SE(3) special Euclidean group.");
75
80
81 /* Expose tangent bundles */
82
83 exposeTangentBundle<SO2>("TSO2", "Tangent bundle of the SO(2) group.")
84 .def(bp::init<>(bp::args("self")));
85 exposeTangentBundle<SO3>("TSO3", "Tangent bundle of the SO(3) group.")
86 .def(bp::init<>(bp::args("self")));
87 exposeTangentBundle<SE2>("TSE2", "Tangent bundle of the SE(2) group.")
88 .def(bp::init<>(bp::args("self")));
89 exposeTangentBundle<SE3>("TSE3", "Tangent bundle of the SE(3) group.")
90 .def(bp::init<>(bp::args("self")));
91
92 /* Groups associated w/ Pinocchio models */
93 using Multibody = MultibodyConfiguration<Scalar>;
94 using Model = ModelTpl<Scalar>;
95 bp::class_<Multibody, bp::bases<Manifold>>(
96 "MultibodyConfiguration", "Configuration group of a multibody",
97 bp::init<const Model &>(bp::args("self", "model")))
98 .add_property("model",
99 bp::make_function(&Multibody::getModel,
100 bp::return_internal_reference<>()),
101 "Return the Pinocchio model instance.")
103 .enable_pickling_(true);
104
105 using MultiPhase = MultibodyPhaseSpace<Scalar>;
106 bp::class_<MultiPhase, bp::bases<Manifold>>(
107 "MultibodyPhaseSpace",
108 "Tangent space of the multibody configuration group.",
109 bp::init<const Model &>(("self"_a, "model")))
110 .add_property("model",
111 bp::make_function(&MultiPhase::getModel,
112 bp::return_internal_reference<>()),
113 "Return the Pinocchio model instance.")
114 .add_property("base", bp::make_function(
115 +[](const MultiPhase &m) -> const Multibody & {
116 return m.getBaseSpace();
117 },
118 bp::return_internal_reference<>()))
120 .enable_pickling_(true);
121}
122
123} // namespace aligator::python
ManifoldAbstractTpl< Scalar > Manifold
Definition context.hpp:14
The Python bindings.
Definition blk-matrix.hpp:5
xyz::polymorphic< Manifold > PolyManifold
void exposeLieGroup(const char *name, const char *docstring)
bp::class_< TangentBundleTpl< M >, bp::bases< Manifold > > exposeTangentBundle(const char *name, const char *docstring)
Expose the tangent bundle of a manifold type M.
bp::object get_namespace(const std::string &name)
Create or retrieve a Python scope (that is, a class or module namespace).
Definition utils.hpp:22
Multibody configuration group , defined using the Pinocchio library.
Definition multibody.hpp:17
The tangent bundle of a multibody configuration group.
Wrap a Pinocchio Lie group into a ManifoldAbstractTpl object.
Tangent bundle of a base manifold M.