11#include <eigen3/Eigen/Core>
12#include <eigen3/Eigen/Dense>
18inline Eigen::IOFormat
__CleanFmt(4, 0,
", ",
"\n",
"[",
"]");
20template <
typename T,
typename Scalar>
23 for (std::size_t i = 0; i < lb.size(); i++) {
24 Scalar dWidth = ub[i] - lb[i];
32template <
typename T,
typename T2,
typename Scalar>
34 int &ii, Scalar &width) {
35 for (std::size_t i = 0; i < lb.size(); i++) {
36 Scalar dWidth = (ub[i] - lb[i]) * weight(i);
45template <
typename Scalar,
int Dimensions = -1>
struct RnL1 {
47 using cref_t =
const Eigen::Ref<const Eigen::Matrix<Scalar, Dimensions, 1>> &;
48 using ref_t = Eigen::Ref<Eigen::Matrix<Scalar, Dimensions, 1>>;
50 Eigen::Matrix<Scalar, Dimensions, 1>
lb;
51 Eigen::Matrix<Scalar, Dimensions, 1>
ub;
52 Eigen::Matrix<Scalar, Dimensions, 1>
weights;
61 assert(weights_.size() ==
weights.size());
67 out <<
"State Space: RnL1" <<
" RuntimeDIM: " <<
lb.size()
68 <<
" CompileTimeDIM: " << Dimensions << std::endl
70 <<
"ub: " <<
ub.transpose().format(
__CleanFmt) << std::endl;
78 for (
size_t i = 0; i < x.size(); i++) {
93 x =
lb + (
ub -
lb).cwiseProduct(x);
99 out = from + t * (to - from);
103 Scalar &width)
const {
115 if constexpr (Dimensions == Eigen::Dynamic) {
120 assert(x.size() ==
ub.size());
121 assert(x.size() ==
lb.size());
123 for (
size_t i = 0; i < x.size(); i++) {
124 Scalar xx = std::max(
lb(i), std::min(
ub(i), x(i)));
125 Scalar dif = xx - x(i);
129 for (
size_t i = 0; i < Dimensions; i++) {
130 Scalar xx = std::max(
lb(i), std::min(
ub(i), x(i)));
131 Scalar dif = xx - x(i);
144 return (x - y).cwiseAbs().cwiseProduct(
weights).sum();
146 return (x - y).cwiseAbs().sum();
150template <
typename Scalar>
struct Time {
152 using cref_t =
const Eigen::Ref<const Eigen::Matrix<Scalar, 1, 1>> &;
153 using ref_t = Eigen::Ref<Eigen::Matrix<Scalar, 1, 1>>;
154 Eigen::Matrix<Scalar, 1, 1>
lb;
155 Eigen::Matrix<Scalar, 1, 1>
ub;
162 for (
size_t i = 0; i < x.size(); i++) {
174 out <<
"Time: " <<
lb(0) <<
" " <<
ub(0) << std::endl;
179 x(0) = double(rand()) / RAND_MAX * (
ub(0) -
lb(0)) +
lb(0);
183 assert(lb_.size() == 1);
184 assert(ub_.size() == 1);
185 assert(
ub(0) >=
lb(0));
194 assert(to(0) >= from(0));
196 Eigen::Matrix<Scalar, 1, 1> d = to - from;
208 return std::numeric_limits<Scalar>::max();
209 }
else if (x(0) >
lb(0)) {
221 return std::numeric_limits<Scalar>::max();
228template <
typename Scalar>
struct SO2 {
230 using cref_t =
const Eigen::Ref<const Eigen::Matrix<Scalar, 1, 1>> &;
231 using ref_t = Eigen::Ref<Eigen::Matrix<Scalar, 1, 1>>;
237 for (
size_t i = 0; i < x.size(); i++) {
249 out <<
"SO2: " << std::endl
250 <<
"weight: " <<
weight << std::endl
256 x(0) = (double(rand()) / (RAND_MAX + 1.)) * 2. * M_PI - M_PI;
265 assert(weights_.size() == 1);
275 Eigen::Matrix<Scalar, 1, 1> d = to - from;
279 }
else if (d(0) < -M_PI) {
287 }
else if (out(0) < -M_PI) {
298 assert(x(0) >= -M_PI);
299 assert(x(0) <= M_PI);
301 assert(lb(0) >= -M_PI);
302 assert(lb(0) <= M_PI);
304 assert(ub(0) >= -M_PI);
305 assert(ub(0) <= M_PI);
307 if (x(0) >= lb(0) && x(0) <= ub(0)) {
309 }
else if (x(0) > ub(0)) {
310 Scalar d1 = x(0) - ub(0);
311 Scalar d2 = lb(0) - (x(0) - 2 * M_PI);
315 }
else if (x(0) < lb(0)) {
316 Scalar d1 = lb(0) - x(0);
317 Scalar d2 = (x(0) + 2 * M_PI) - ub(0);
329 assert(x(0) >= -M_PI);
330 assert(y(0) >= -M_PI);
332 assert(x(0) <= M_PI);
333 assert(y(0) <= M_PI);
335 Scalar dif = x(0) - y(0);
338 }
else if (dif < -M_PI) {
341 Scalar out = std::abs(dif);
348 using cref_t =
const Eigen::Ref<const Eigen::Matrix<Scalar, 1, 1>> &;
349 using ref_t = Eigen::Ref<Eigen::Matrix<Scalar, 1, 1>>;
351 Eigen::Matrix<Scalar, 1, 1>
lb;
352 Eigen::Matrix<Scalar, 1, 1>
ub;
356 void print(std::ostream &out) { out <<
"SO2Squared: " << std::endl; }
359 for (
size_t i = 0; i < x.size(); i++) {
371 so2.interpolate(from, to, t, out);
392 Scalar d =
so2.distance_to_rectangle(x,
lb,
ub);
398 Scalar d =
so2.distance(x, y);
403template <
typename Scalar,
int Dimensions = -1>
struct RnSquared {
405 using cref_t =
const Eigen::Ref<const Eigen::Matrix<Scalar, Dimensions, 1>> &;
406 using ref_t = Eigen::Ref<Eigen::Matrix<Scalar, Dimensions, 1>>;
407 using vec_t = Eigen::Matrix<Scalar, Dimensions, 1>;
409 Eigen::Matrix<Scalar, Dimensions, 1>
lb;
410 Eigen::Matrix<Scalar, Dimensions, 1>
ub;
415 out <<
"State Space: RnSquared" <<
" RuntimeDIM: " <<
lb.size()
416 <<
" CompileTimeDIM: " << Dimensions << std::endl
418 <<
"ub: " <<
ub.transpose().format(
__CleanFmt) << std::endl;
426 for (
size_t i = 0; i < x.size(); i++) {
446 out = from + t * (to - from);
455 Scalar &width)
const {
466 x =
lb + (
ub -
lb).cwiseProduct(x);
474 if constexpr (Dimensions == Eigen::Dynamic) {
479 assert(x.size() ==
ub.size());
480 assert(x.size() ==
lb.size());
482 for (
size_t i = 0; i < x.size(); i++) {
483 Scalar xx = std::max(
lb(i), std::min(
ub(i), x(i)));
484 Scalar dif = xx - x(i);
488 for (
size_t i = 0; i < Dimensions; i++) {
489 Scalar xx = std::max(
lb(i), std::min(
ub(i), x(i)));
490 Scalar dif = xx - x(i);
503 return (x - y).cwiseProduct(
weights).squaredNorm();
505 return (x - y).squaredNorm();
509template <
typename Scalar,
int Dimensions = -1>
struct Rn {
510 using cref_t =
const Eigen::Ref<const Eigen::Matrix<Scalar, Dimensions, 1>> &;
511 using ref_t = Eigen::Ref<Eigen::Matrix<Scalar, Dimensions, 1>>;
512 using vec_t = Eigen::Matrix<Scalar, Dimensions, 1>;
513 using DIM = std::integral_constant<int, Dimensions>;
516 Eigen::Matrix<Scalar, Dimensions, 1>
lb;
517 Eigen::Matrix<Scalar, Dimensions, 1>
ub;
522 out <<
"State Space: Rn" <<
" RuntimeDIM: " <<
lb.size()
523 <<
" CompileTimeDIM: " << Dimensions << std::endl
525 <<
"ub: " <<
ub.transpose().format(
__CleanFmt) << std::endl;
544 for (
size_t i = 0; i < x.size(); i++) {
558 out = from + t * (to - from);
565 x =
lb + (
ub -
lb).cwiseProduct(x);
569 Scalar &width)
const {
592 using cref_t =
const Eigen::Ref<
const Eigen::Matrix<double, -1, 1>> &;
593 using ref_t = Eigen::Ref<Eigen::Matrix<double, -1, 1>>;
600 ref_t out)
const = 0;
617 using cref_t =
const Eigen::Ref<
const Eigen::Matrix<double, -1, 1>> &;
618 using ref_t = Eigen::Ref<Eigen::Matrix<double, -1, 1>>;
625 ref_t out)
const override {
639 cref_t &ub)
const override {
652 using cref_t =
const Eigen::Ref<
const Eigen::Matrix<double, -1, 1>> &;
653 using ref_t = Eigen::Ref<Eigen::Matrix<double, -1, 1>>;
655 std::shared_ptr<Vpure>
s4;
660 s4->interpolate(from, to, t, out);
667 s4->choose_split_dimension(lb, ub, ii, width);
671 return s4->distance_to_rectangle(x, lb, ub);
675 return s4->distance(x, y);
681 static constexpr int value =
id + 1;
691template <
typename Scalar,
int Dimensions>
struct RnTime {
697 const Eigen::Ref<const Eigen::Matrix<Scalar, effective_dim, 1>> &;
698 using ref_t = Eigen::Ref<Eigen::Matrix<Scalar, effective_dim, 1>>;
707 assert(lambda_t_ >= 0);
708 assert(lambda_r_ >= 0);
714 out <<
"RnTime: " << std::endl;
723 time.interpolate(from.template tail<1>(), to.template tail<1>(), t,
724 out.template tail<1>());
726 if constexpr (Dimensions == Eigen::Dynamic) {
727 size_t n = from.size() - 1;
734 to.template head<Dimensions>(), t,
735 out.template head<Dimensions>());
741 assert(lb_.size() == ub_.size());
743 time.set_bounds(lb_.template tail<1>(), ub_.template tail<1>());
745 if constexpr (Dimensions == Eigen::Dynamic) {
746 size_t n = lb_.size() - 1;
750 ub_.template head<Dimensions>());
755 Scalar &width)
const {
761 time.sample_uniform(x.template tail<1>());
763 if constexpr (Dimensions == Eigen::Dynamic) {
764 size_t n = x.size() - 1;
773 double dt =
time.distance_to_rectangle(
774 x.template tail<1>(), lb.template tail<1>(), ub.template tail<1>());
777 if (dt == std::numeric_limits<Scalar>::max()) {
778 return std::numeric_limits<Scalar>::max();
781 if constexpr (Dimensions == Eigen::Dynamic) {
782 size_t n = x.size() - 1;
786 lb.template head<Dimensions>(),
787 ub.template head<Dimensions>());
794 double dt =
time.distance(x.template tail<1>(), y.template tail<1>());
797 if (dt == std::numeric_limits<Scalar>::max()) {
798 return std::numeric_limits<Scalar>::max();
801 if constexpr (Dimensions == Eigen::Dynamic) {
802 size_t n = x.size() - 1;
805 dr =
rn.
distance(x.template head<Dimensions>(),
806 y.template head<Dimensions>());
812template <
typename Scalar>
struct R2SO2 {
814 using cref_t =
const Eigen::Ref<const Eigen::Matrix<Scalar, 3, 1>> &;
815 using ref_t = Eigen::Ref<Eigen::Matrix<Scalar, 3, 1>>;
816 using vec_t = Eigen::Matrix<Scalar, 3, 1>;
818 using cref2_t =
const Eigen::Ref<const Eigen::Matrix<Scalar, 2, 1>> &;
819 using ref2_t = Eigen::Ref<Eigen::Matrix<Scalar, 2, 1>>;
834 out <<
"R2SO2: " << std::endl;
843 weights.template head<2>() = wr2;
846 so2.set_weights(wso2);
853 so2.check_bounds(x.template tail<1>());
858 so2.sample_uniform(x.template tail<1>());
864 l2.
interpolate(from.template head<2>(), to.template head<2>(), t,
865 out.template head<2>());
866 so2.interpolate(from.template tail<1>(), to.template tail<1>(), t,
867 out.template tail<1>());
873 x.template head<2>(), lb.template head<2>(), ub.template head<2>());
874 Scalar d2 =
so2.distance_to_rectangle(
875 x.template tail<1>(), lb.template tail<1>(), ub.template tail<1>());
881 Scalar d1 =
l2.
distance(x.template head<2>(), y.template head<2>());
882 Scalar d2 =
so2.distance(x.template tail<1>(), y.template tail<1>());
889 using cref_t =
const Eigen::Ref<const Eigen::Matrix<Scalar, 3, 1>> &;
890 using ref_t = Eigen::Ref<Eigen::Matrix<Scalar, 3, 1>>;
892 using cref2_t =
const Eigen::Ref<const Eigen::Matrix<Scalar, 2, 1>> &;
893 using ref2_t = Eigen::Ref<Eigen::Matrix<Scalar, 2, 1>>;
905 out <<
"R2SO2Squared: " << std::endl;
914 so2squared.sample_uniform(x.template tail<1>());
920 x.template head<2>(), lb.template head<2>(), ub.template head<2>());
922 x.template tail<1>(), lb.template tail<1>(), ub.template tail<1>());
929 Scalar d2 =
so2squared.distance(x.template tail<1>(), y.template tail<1>());
936 using cref_t =
const Eigen::Ref<const Eigen::Matrix<Scalar, 4, 1>> &;
937 using ref_t = Eigen::Ref<Eigen::Matrix<Scalar, 4, 1>>;
942 x = Eigen::Quaterniond().UnitRandom().coeffs();
948 out <<
"SO3Squared: " << std::endl;
970 assert(std::abs(x.norm() - 1) < 1e-6);
974 return std::min(d1, d2);
979 assert(x.size() == 4);
980 assert(y.size() == 4);
982 assert(std::abs(x.norm() - 1) < 1e-6);
983 assert(std::abs(y.norm() - 1) < 1e-6);
987 return std::min(d1, d2);
991template <
typename Scalar>
struct SO3 {
993 using cref_t =
const Eigen::Ref<const Eigen::Matrix<Scalar, 4, 1>> &;
994 using ref_t = Eigen::Ref<Eigen::Matrix<Scalar, 4, 1>>;
999 out <<
"SO3: " << std::endl;
1020 out = Eigen::Quaterniond(from).slerp(t, Eigen::Quaterniond(to)).coeffs();
1025 return std::sqrt(
so3squared.distance_to_rectangle(x, lb, ub));
1040 using cref_t =
const Eigen::Ref<const Eigen::Matrix<Scalar, 7, 1>> &;
1041 using ref_t = Eigen::Ref<Eigen::Matrix<Scalar, 7, 1>>;
1042 using cref3_t =
const Eigen::Ref<const Eigen::Matrix<Scalar, 3, 1>> &;
1052 out <<
"R3SO3Squared: " << std::endl;
1060 so3.sample_uniform(x.template tail<4>());
1066 x.template head<3>(), lb.template head<3>(), ub.template head<3>());
1068 Scalar d2 =
so3.distance_to_rectangle(
1069 x.template tail<4>(), lb.template tail<4>(), ub.template tail<4>());
1076 Scalar d1 =
l2.
distance(x.template head<3>(), y.template head<3>());
1077 Scalar d2 =
so3.distance(x.template tail<4>(), y.template tail<4>());
1082template <
typename Scalar>
struct R3SO3 {
1084 using cref_t =
const Eigen::Ref<const Eigen::Matrix<Scalar, 7, 1>> &;
1085 using ref_t = Eigen::Ref<Eigen::Matrix<Scalar, 7, 1>>;
1087 using cref3_t =
const Eigen::Ref<const Eigen::Matrix<Scalar, 3, 1>> &;
1088 using ref3_t = Eigen::Ref<Eigen::Matrix<Scalar, 3, 1>>;
1095 l2.
interpolate(from.template head<3>(), to.template head<3>(), t,
1096 out.template head<3>());
1097 so3.interpolate(from.template tail<4>(), to.template tail<4>(), t,
1098 out.template tail<4>());
1105 out <<
"R3SO3: " << std::endl;
1114 so3.sample_uniform(x.template tail<4>());
1120 x.template head<3>(), lb.template head<3>(), ub.template head<3>());
1122 Scalar d2 =
so3.distance_to_rectangle(
1123 x.template tail<4>(), lb.template tail<4>(), ub.template tail<4>());
1130 Scalar d1 =
l2.
distance(x.template head<3>(), y.template head<3>());
1131 Scalar d2 =
so3.distance(x.template tail<4>(), y.template tail<4>());
1146inline bool starts_with(
const std::string &str,
const std::string &prefix) {
1147 return str.size() >= prefix.size() &&
1148 str.compare(0, prefix.size(), prefix) == 0;
1153 std::string delimiter =
":";
1157 pos = str.find(delimiter);
1158 if (pos == std::string::npos) {
1162 token = str.substr(pos + delimiter.length(), str.size());
1163 int out = std::stoi(token);
1164 std::cout <<
"out " << out << std::endl;
1169 using cref_t =
const Eigen::Ref<
const Eigen::Matrix<Scalar, -1, 1>> &;
1170 using ref_t = Eigen::Ref<Eigen::Matrix<Scalar, -1, 1>>;
1179 Eigen::Matrix<Scalar, -1, 1>
1192 for (
size_t i = 0; i <
spaces.size(); i++) {
1195 obj.set_weights(weights_.segment(counter,
dims[i]));
1204 for (
size_t i = 0; i <
spaces.size(); i++) {
1224 out <<
"Combined: " << std::endl;
1226 out << a << std::endl;
1228 std::visit([&](
auto &obj) { obj.print(out); }, s);
1229 for (
auto &d :
dims)
1230 out << d << std::endl;
1231 out <<
"lb " <<
lb.transpose().format(
__CleanFmt) << std::endl;
1232 out <<
"ub " <<
ub.transpose().format(
__CleanFmt) << std::endl;
1238 for (
size_t i = 0; i < spaces_str.size(); i++) {
1239 if (spaces_str.at(i) ==
"SO2") {
1243 }
else if (spaces_str.at(i) ==
"SO2Squared") {
1247 }
else if (spaces_str.at(i) ==
"SO3") {
1251 }
else if (spaces_str.at(i) ==
"SO3Squared") {
1255 }
else if (
starts_with(spaces_str.at(i),
"RnL1")) {
1259 dims.push_back(dim);
1265 dims.push_back(dim);
1266 }
else if (
starts_with(spaces_str.at(i),
"RnSquared")) {
1270 dims.push_back(dim);
1287 for (
size_t i = 0; i <
spaces.size(); i++) {
1289 [&](
const auto &obj) {
1290 return obj.check_bounds(x.segment(counter,
dims[i]));
1302 assert(lbs.size() == ubs.size());
1310 if (space_name ==
"SO2")
1312 if (space_name ==
"SO2Squared")
1314 if (space_name ==
"SO3")
1316 if (space_name ==
"SO3Squared")
1322 obj.set_bounds(lbs.segment(counter,
dims[i]),
1323 ubs.segment(counter,
dims[i]));
1331 const std::vector<Eigen::VectorXd> &ubs) {
1333 assert(lbs.size() == ubs.size());
1335 for (
size_t i = 0; i < lbs.size(); i++) {
1339 obj.set_bounds(lbs[i], ubs[i]);
1350 for (
size_t i = 0; i <
spaces.size(); i++) {
1352 [&](
const auto &obj) {
1353 obj.sample_uniform(x.segment(counter,
dims[i]));
1366 for (
size_t i = 0; i <
spaces.size(); i++) {
1368 [&](
const auto &obj) {
1369 obj.interpolate(from.segment(counter,
dims[i]),
1370 to.segment(counter,
dims[i]), t,
1371 out.segment(counter,
dims[i]));
1385 auto caller = [&](
const auto &obj) {
1386 return obj.distance(x.segment(counter,
dims[dim_index]),
1387 y.segment(counter,
dims[dim_index]));
1390 for (
size_t i = 0; i <
spaces.size(); i++) {
1392 d += std::visit(caller,
spaces[i]);
1407 auto caller = [&](
const auto &obj) {
1408 return obj.distance_to_rectangle(x.segment(counter,
dims[dim_index]),
1409 lb.segment(counter,
dims[dim_index]),
1410 ub.segment(counter,
dims[dim_index]));
1413 for (
size_t i = 0; i <
spaces.size(); i++) {
1415 d += std::visit(caller,
spaces[i]);
#define THROW_PRETTY_DYNOTREE(arg)
Definition dynotree_macros.h:33
#define CHECK_PRETTY_DYNOTREE(condition, arg)
Definition dynotree_macros.h:36
#define CHECK_PRETTY_DYNOTREE__(condition)
Definition dynotree_macros.h:41
Definition dynotree_macros.h:9
Eigen::IOFormat __CleanFmt(4, 0, ", ", "\n", "[", "]")
DistanceType
Definition StateSpace.h:1136
void choose_split_dimension_default(const T &lb, const T &ub, int &ii, Scalar &width)
Definition StateSpace.h:21
int get_number(const std::string &str)
Definition StateSpace.h:1152
void choose_split_dimension_weights(const T &lb, const T &ub, const T2 &weight, int &ii, Scalar &width)
Definition StateSpace.h:33
bool starts_with(const std::string &str, const std::string &prefix)
Definition StateSpace.h:1146
Definition StateSpace.h:679
static constexpr int value
Definition StateSpace.h:681
Definition StateSpace.h:1168
Combined(const std::vector< std::string > &spaces_str)
Definition StateSpace.h:1236
std::vector< std::string > spaces_names
Definition StateSpace.h:1182
std::variant< RnL1< Scalar >, Rn< Scalar >, RnSquared< Scalar >, SO2< Scalar >, SO2Squared< Scalar >, SO3< Scalar >, SO3Squared< Scalar > > Space
Definition StateSpace.h:1172
Eigen::VectorXd ub
Definition StateSpace.h:1184
Scalar distance(cref_t x, cref_t y) const
Definition StateSpace.h:1378
void set_weights(cref_t weights_)
Definition StateSpace.h:1186
void set_bounds(const std::vector< Eigen::VectorXd > &lbs, const std::vector< Eigen::VectorXd > &ubs)
Definition StateSpace.h:1330
const Eigen::Ref< const Eigen::Matrix< Scalar, -1, 1 > > & cref_t
Definition StateSpace.h:1169
Eigen::VectorXd lb
Definition StateSpace.h:1183
void interpolate(cref_t from, cref_t to, Scalar t, ref_t out) const
Definition StateSpace.h:1360
std::vector< Space > spaces
Definition StateSpace.h:1175
std::vector< int > dims
Definition StateSpace.h:1176
Combined(const std::vector< Space > &spaces, const std::vector< int > &dims)
Definition StateSpace.h:1217
Eigen::Matrix< Scalar, -1, 1 > weights
Definition StateSpace.h:1180
Scalar distance_to_rectangle(cref_t x, cref_t lb, cref_t ub) const
Definition StateSpace.h:1398
void set_bounds(cref_t &lbs, cref_t &ubs)
Definition StateSpace.h:1300
bool use_weights
Definition StateSpace.h:1181
bool check_bounds(cref_t x) const
Definition StateSpace.h:1285
int get_runtime_dim()
Definition StateSpace.h:1202
void print(std::ostream &out)
Definition StateSpace.h:1222
void choose_split_dimension(cref_t lb, cref_t ub, int &ii, Scalar &width)
Definition StateSpace.h:1278
void sample_uniform(ref_t x) const
Definition StateSpace.h:1345
Eigen::Ref< Eigen::Matrix< Scalar, -1, 1 > > ref_t
Definition StateSpace.h:1170
Definition StateSpace.h:887
Scalar distance(cref_t x, cref_t y) const
Definition StateSpace.h:926
const Eigen::Ref< const Eigen::Matrix< Scalar, 2, 1 > > & cref2_t
Definition StateSpace.h:892
const Eigen::Ref< const Eigen::Matrix< Scalar, 3, 1 > > & cref_t
Definition StateSpace.h:889
Eigen::Ref< Eigen::Matrix< Scalar, 3, 1 > > ref_t
Definition StateSpace.h:890
SO2Squared< Scalar > so2squared
Definition StateSpace.h:902
Scalar distance_to_rectangle(cref_t x, cref_t lb, cref_t ub) const
Definition StateSpace.h:917
void print(std::ostream &out)
Definition StateSpace.h:904
void sample_uniform(ref_t x) const
Definition StateSpace.h:912
void set_bounds(cref2_t lb_, cref2_t ub_)
Definition StateSpace.h:910
Eigen::Ref< Eigen::Matrix< Scalar, 2, 1 > > ref2_t
Definition StateSpace.h:893
Scalar angular_weight
Definition StateSpace.h:899
void choose_split_dimension(cref_t lb, cref_t ub, int &ii, Scalar &width)
Definition StateSpace.h:895
RnSquared< Scalar, 2 > rn_squared
Definition StateSpace.h:901
Definition StateSpace.h:812
Eigen::Ref< Eigen::Matrix< Scalar, 2, 1 > > ref2_t
Definition StateSpace.h:819
Scalar distance_to_rectangle(cref_t x, cref_t lb, cref_t ub) const
Definition StateSpace.h:870
Eigen::Ref< Eigen::Matrix< Scalar, 3, 1 > > ref_t
Definition StateSpace.h:815
Scalar angular_weight
Definition StateSpace.h:828
bool check_bounds(cref_t x) const
Definition StateSpace.h:851
void print(std::ostream &out)
Definition StateSpace.h:833
bool use_weights
Definition StateSpace.h:840
void set_bounds(cref2_t lb_, cref2_t ub_)
Definition StateSpace.h:849
Scalar distance(cref_t x, cref_t y) const
Definition StateSpace.h:879
Rn< Scalar, 2 > l2
Definition StateSpace.h:830
vec_t weights
Definition StateSpace.h:839
void set_weights(cref_t wr2, double wso2)
Definition StateSpace.h:842
void choose_split_dimension(cref_t lb, cref_t ub, int &ii, Scalar &width)
Definition StateSpace.h:821
void sample_uniform(ref_t x) const
Definition StateSpace.h:856
void interpolate(cref_t from, cref_t to, Scalar t, ref_t out) const
Definition StateSpace.h:861
Eigen::Matrix< Scalar, 3, 1 > vec_t
Definition StateSpace.h:816
const Eigen::Ref< const Eigen::Matrix< Scalar, 2, 1 > > & cref2_t
Definition StateSpace.h:818
SO2< Scalar > so2
Definition StateSpace.h:831
const Eigen::Ref< const Eigen::Matrix< Scalar, 3, 1 > > & cref_t
Definition StateSpace.h:814
Definition StateSpace.h:1038
Scalar distance_to_rectangle(cref_t &x, cref_t &lb, cref_t &ub) const
Definition StateSpace.h:1063
void sample_uniform(cref3_t lb, cref3_t ub, ref_t x) const
Definition StateSpace.h:1058
const Eigen::Ref< const Eigen::Matrix< Scalar, 7, 1 > > & cref_t
Definition StateSpace.h:1040
void choose_split_dimension(cref_t lb, cref_t ub, int &ii, Scalar &width)
Definition StateSpace.h:1044
RnSquared< Scalar, 3 > l2
Definition StateSpace.h:1048
const Eigen::Ref< const Eigen::Matrix< Scalar, 3, 1 > > & cref3_t
Definition StateSpace.h:1042
SO3Squared< Scalar > so3
Definition StateSpace.h:1049
void set_bounds(cref3_t lb_, cref3_t ub_)
Definition StateSpace.h:1056
Eigen::Ref< Eigen::Matrix< Scalar, 7, 1 > > ref_t
Definition StateSpace.h:1041
void print(std::ostream &out)
Definition StateSpace.h:1051
Scalar distance(cref_t x, cref_t y) const
Definition StateSpace.h:1074
Definition StateSpace.h:1082
void sample_uniform(ref_t x) const
Definition StateSpace.h:1112
SO3< Scalar > so3
Definition StateSpace.h:1102
const Eigen::Ref< const Eigen::Matrix< Scalar, 7, 1 > > & cref_t
Definition StateSpace.h:1084
const Eigen::Ref< const Eigen::Matrix< Scalar, 3, 1 > > & cref3_t
Definition StateSpace.h:1087
void choose_split_dimension(cref_t lb, cref_t ub, int &ii, Scalar &width)
Definition StateSpace.h:1090
Rn< Scalar, 3 > l2
Definition StateSpace.h:1101
Scalar distance(cref_t x, cref_t y) const
Definition StateSpace.h:1128
void set_bounds(cref3_t lb_, cref3_t ub_)
Definition StateSpace.h:1110
void interpolate(cref_t from, cref_t to, Scalar t, ref_t out) const
Definition StateSpace.h:1094
void print(std::ostream &out)
Definition StateSpace.h:1104
Eigen::Ref< Eigen::Matrix< Scalar, 7, 1 > > ref_t
Definition StateSpace.h:1085
Eigen::Ref< Eigen::Matrix< Scalar, 3, 1 > > ref3_t
Definition StateSpace.h:1088
Scalar distance_to_rectangle(cref_t &x, cref_t &lb, cref_t &ub) const
Definition StateSpace.h:1117
Definition StateSpace.h:1035
Definition StateSpace.h:45
void set_bounds(cref_t lb_, cref_t ub_)
Definition StateSpace.h:55
Eigen::Ref< Eigen::Matrix< Scalar, Dimensions, 1 > > ref_t
Definition StateSpace.h:48
Eigen::Matrix< Scalar, Dimensions, 1 > ub
Definition StateSpace.h:51
void choose_split_dimension(cref_t lb, cref_t ub, int &ii, Scalar &width) const
Definition StateSpace.h:102
bool use_weights
Definition StateSpace.h:53
void sample_uniform(ref_t x) const
Definition StateSpace.h:89
void interpolate(cref_t from, cref_t to, Scalar t, ref_t out) const
Definition StateSpace.h:96
void set_weights(cref_t weights_)
Definition StateSpace.h:60
void print(std::ostream &out)
Definition StateSpace.h:66
bool check_bounds(cref_t x) const
Definition StateSpace.h:73
Scalar distance_to_rectangle(cref_t &x, cref_t &lb, cref_t &ub) const
Definition StateSpace.h:110
Eigen::Matrix< Scalar, Dimensions, 1 > weights
Definition StateSpace.h:52
Scalar distance(cref_t x, cref_t y) const
Definition StateSpace.h:138
const Eigen::Ref< const Eigen::Matrix< Scalar, Dimensions, 1 > > & cref_t
Definition StateSpace.h:47
Eigen::Matrix< Scalar, Dimensions, 1 > lb
Definition StateSpace.h:50
Definition StateSpace.h:403
void set_bounds(cref_t lb_, cref_t ub_)
Definition StateSpace.h:449
void interpolate(cref_t from, cref_t to, Scalar t, ref_t out) const
Definition StateSpace.h:443
void choose_split_dimension(cref_t lb, cref_t ub, int &ii, Scalar &width) const
Definition StateSpace.h:454
bool use_weights
Definition StateSpace.h:412
void set_weights(cref_t weights_)
Definition StateSpace.h:437
Eigen::Matrix< Scalar, Dimensions, 1 > vec_t
Definition StateSpace.h:407
bool check_bounds(cref_t x) const
Definition StateSpace.h:421
Eigen::Matrix< Scalar, Dimensions, 1 > lb
Definition StateSpace.h:409
const Eigen::Ref< const Eigen::Matrix< Scalar, Dimensions, 1 > > & cref_t
Definition StateSpace.h:405
void sample_uniform(ref_t x) const
Definition StateSpace.h:462
Eigen::Ref< Eigen::Matrix< Scalar, Dimensions, 1 > > ref_t
Definition StateSpace.h:406
Scalar distance(cref_t &x, cref_t &y) const
Definition StateSpace.h:497
Eigen::Matrix< Scalar, Dimensions, 1 > ub
Definition StateSpace.h:410
vec_t weights
Definition StateSpace.h:411
Scalar distance_to_rectangle(cref_t x, cref_t lb, cref_t ub) const
Definition StateSpace.h:469
void print(std::ostream &out)
Definition StateSpace.h:414
Definition StateSpace.h:691
void interpolate(cref_t from, cref_t to, Scalar t, ref_t out) const
Definition StateSpace.h:719
Eigen::Ref< Eigen::Matrix< Scalar, effective_dim, 1 > > ref_t
Definition StateSpace.h:698
void set_bounds(cref_t lb_, cref_t ub_)
Definition StateSpace.h:739
double lambda_t
Definition StateSpace.h:703
static constexpr int effective_dim
Definition StateSpace.h:695
void print(std::ostream &out)
Definition StateSpace.h:713
const Eigen::Ref< const Eigen::Matrix< Scalar, effective_dim, 1 > > & cref_t
Definition StateSpace.h:696
Time< Scalar > time
Definition StateSpace.h:700
void set_lambda(double lambda_t_, double lambda_r_)
Definition StateSpace.h:706
Rn< Scalar, Dimensions > rn
Definition StateSpace.h:701
Scalar distance(cref_t x, cref_t y) const
Definition StateSpace.h:792
void choose_split_dimension(cref_t lb, cref_t ub, int &ii, Scalar &width) const
Definition StateSpace.h:754
void sample_uniform(ref_t x) const
Definition StateSpace.h:759
Scalar distance_to_rectangle(cref_t &x, cref_t &lb, cref_t &ub) const
Definition StateSpace.h:771
double lambda_r
Definition StateSpace.h:704
Definition StateSpace.h:509
void sample_uniform(ref_t x) const
Definition StateSpace.h:561
bool check_bounds(cref_t x) const
Definition StateSpace.h:539
Scalar distance(cref_t &x, cref_t &y) const
Definition StateSpace.h:582
bool use_weights
Definition StateSpace.h:519
Eigen::Matrix< Scalar, Dimensions, 1 > lb
Definition StateSpace.h:516
std::integral_constant< int, Dimensions > DIM
Definition StateSpace.h:513
const Eigen::Ref< const Eigen::Matrix< Scalar, Dimensions, 1 > > & cref_t
Definition StateSpace.h:510
void choose_split_dimension(cref_t lb, cref_t ub, int &ii, Scalar &width) const
Definition StateSpace.h:568
RnSquared< Scalar, Dimensions > rn_squared
Definition StateSpace.h:515
void set_bounds(cref_t lb_, cref_t ub_)
Definition StateSpace.h:534
Eigen::Matrix< Scalar, Dimensions, 1 > ub
Definition StateSpace.h:517
void print(std::ostream &out)
Definition StateSpace.h:521
Eigen::Matrix< Scalar, Dimensions, 1 > vec_t
Definition StateSpace.h:512
void set_weights(cref_t weights_)
Definition StateSpace.h:528
vec_t weights
Definition StateSpace.h:518
Scalar distance_to_rectangle(cref_t &x, cref_t &lb, cref_t &ub) const
Definition StateSpace.h:576
void interpolate(cref_t from, cref_t to, Scalar t, ref_t out) const
Definition StateSpace.h:555
Eigen::Ref< Eigen::Matrix< Scalar, Dimensions, 1 > > ref_t
Definition StateSpace.h:511
Definition StateSpace.h:613
Rn< double, 4 > rn
Definition StateSpace.h:620
void set_bounds(cref_t lb_, cref_t ub_) override
Definition StateSpace.h:622
double Scalar
Definition StateSpace.h:615
Eigen::Ref< Eigen::Matrix< double, -1, 1 > > ref_t
Definition StateSpace.h:618
virtual Scalar distance(cref_t &x, cref_t &y) const override
Definition StateSpace.h:643
virtual void sample_uniform(ref_t x) const override
Definition StateSpace.h:629
const Eigen::Ref< const Eigen::Matrix< double, -1, 1 > > & cref_t
Definition StateSpace.h:617
virtual void interpolate(cref_t from, cref_t to, Scalar t, ref_t out) const override
Definition StateSpace.h:624
virtual void choose_split_dimension(cref_t lb, cref_t ub, int &ii, Scalar &width) override
Definition StateSpace.h:633
virtual Scalar distance_to_rectangle(cref_t &x, cref_t &lb, cref_t &ub) const override
Definition StateSpace.h:638
Definition StateSpace.h:346
Eigen::Matrix< Scalar, 1, 1 > ub
Definition StateSpace.h:352
void print(std::ostream &out)
Definition StateSpace.h:356
bool check_bounds(cref_t x) const
Definition StateSpace.h:358
SO2< Scalar > so2
Definition StateSpace.h:354
void set_weights(cref_t weights_)
Definition StateSpace.h:374
void choose_split_dimension(cref_t lb, cref_t ub, int &ii, Scalar &width)
Definition StateSpace.h:385
void set_bounds(cref_t lb_, cref_t ub_)
Definition StateSpace.h:380
Scalar distance(cref_t x, cref_t y) const
Definition StateSpace.h:396
void sample_uniform(ref_t x) const
Definition StateSpace.h:378
Scalar distance_to_rectangle(cref_t x, cref_t lb, cref_t ub) const
Definition StateSpace.h:390
void interpolate(cref_t from, cref_t to, Scalar t, ref_t out) const
Definition StateSpace.h:370
const Eigen::Ref< const Eigen::Matrix< Scalar, 1, 1 > > & cref_t
Definition StateSpace.h:348
Eigen::Ref< Eigen::Matrix< Scalar, 1, 1 > > ref_t
Definition StateSpace.h:349
Eigen::Matrix< Scalar, 1, 1 > lb
Definition StateSpace.h:351
Definition StateSpace.h:228
Eigen::Ref< Eigen::Matrix< Scalar, 1, 1 > > ref_t
Definition StateSpace.h:231
void print(std::ostream &out)
Definition StateSpace.h:248
bool use_weights
Definition StateSpace.h:234
void interpolate(cref_t from, cref_t to, Scalar t, ref_t out) const
Definition StateSpace.h:270
Scalar distance(cref_t x, cref_t y) const
Definition StateSpace.h:327
void set_bounds(cref_t lb_, cref_t ub_)
Definition StateSpace.h:259
double weight
Definition StateSpace.h:233
void sample_uniform(ref_t x) const
Definition StateSpace.h:254
bool check_bounds(cref_t x) const
Definition StateSpace.h:236
void choose_split_dimension(cref_t lb, cref_t ub, int &ii, Scalar &width)
Definition StateSpace.h:292
Scalar distance_to_rectangle(cref_t x, cref_t lb, cref_t ub) const
Definition StateSpace.h:296
const Eigen::Ref< const Eigen::Matrix< Scalar, 1, 1 > > & cref_t
Definition StateSpace.h:230
void set_weights(cref_t weights_)
Definition StateSpace.h:264
Definition StateSpace.h:934
void print(std::ostream &out)
Definition StateSpace.h:947
const Eigen::Ref< const Eigen::Matrix< Scalar, 4, 1 > > & cref_t
Definition StateSpace.h:936
Scalar distance(cref_t x, cref_t y) const
Definition StateSpace.h:977
void sample_uniform(ref_t x) const
Definition StateSpace.h:941
bool check_bounds(cref_t x) const
Definition StateSpace.h:945
void interpolate(cref_t from, cref_t to, Scalar t, ref_t out) const
Definition StateSpace.h:960
RnSquared< Scalar, 4 > rn_squared
Definition StateSpace.h:939
Scalar distance_to_rectangle(cref_t &x, cref_t &lb, cref_t &ub) const
Definition StateSpace.h:968
void choose_split_dimension(cref_t lb, cref_t ub, int &ii, Scalar &width)
Definition StateSpace.h:964
void set_weights(cref_t weights_)
Definition StateSpace.h:951
Eigen::Ref< Eigen::Matrix< Scalar, 4, 1 > > ref_t
Definition StateSpace.h:937
void set_bounds(cref_t lb_, cref_t ub_)
Definition StateSpace.h:955
Definition StateSpace.h:991
void sample_uniform(ref_t x) const
Definition StateSpace.h:1005
void interpolate(cref_t from, cref_t to, Scalar t, ref_t out) const
Definition StateSpace.h:1019
void print(std::ostream &out)
Definition StateSpace.h:998
void set_weights(cref_t weights_)
Definition StateSpace.h:1011
void set_bounds(cref_t lb_, cref_t ub_)
Definition StateSpace.h:1007
SO3Squared< Scalar > so3squared
Definition StateSpace.h:996
const Eigen::Ref< const Eigen::Matrix< Scalar, 4, 1 > > & cref_t
Definition StateSpace.h:993
void choose_split_dimension(cref_t lb, cref_t ub, int &ii, Scalar &width)
Definition StateSpace.h:1015
bool check_bounds(cref_t x) const
Definition StateSpace.h:1003
Scalar distance(cref_t x, cref_t y) const
Definition StateSpace.h:1028
Eigen::Ref< Eigen::Matrix< Scalar, 4, 1 > > ref_t
Definition StateSpace.h:994
Scalar distance_to_rectangle(cref_t &x, cref_t &lb, cref_t &ub) const
Definition StateSpace.h:1023
Definition StateSpace.h:150
Eigen::Matrix< Scalar, 1, 1 > lb
Definition StateSpace.h:154
void sample_uniform(ref_t x) const
Definition StateSpace.h:177
const Eigen::Ref< const Eigen::Matrix< Scalar, 1, 1 > > & cref_t
Definition StateSpace.h:152
Scalar distance_to_rectangle(cref_t x, cref_t lb, cref_t ub) const
Definition StateSpace.h:205
void choose_split_dimension(cref_t lb, cref_t ub, int &ii, Scalar &width)
Definition StateSpace.h:200
void set_bounds(cref_t lb_, cref_t ub_)
Definition StateSpace.h:182
Scalar distance(cref_t x, cref_t y) const
Definition StateSpace.h:218
Eigen::Matrix< Scalar, 1, 1 > ub
Definition StateSpace.h:155
void print(std::ostream &out)
Definition StateSpace.h:173
Eigen::Ref< Eigen::Matrix< Scalar, 1, 1 > > ref_t
Definition StateSpace.h:153
void interpolate(cref_t from, cref_t to, Scalar t, ref_t out) const
Definition StateSpace.h:190
bool check_bounds(cref_t x) const
Definition StateSpace.h:157
Definition StateSpace.h:588
virtual void sample_uniform(ref_t x) const =0
virtual Scalar distance_to_rectangle(cref_t &x, cref_t &lb, cref_t &ub) const =0
virtual void interpolate(cref_t from, cref_t to, Scalar t, ref_t out) const =0
virtual Scalar distance(cref_t &x, cref_t &y) const =0
double Scalar
Definition StateSpace.h:590
virtual void set_bounds(cref_t lb_, cref_t ub_)=0
const Eigen::Ref< const Eigen::Matrix< double, -1, 1 > > & cref_t
Definition StateSpace.h:592
Eigen::Ref< Eigen::Matrix< double, -1, 1 > > ref_t
Definition StateSpace.h:593
virtual void choose_split_dimension(cref_t lb, cref_t ub, int &ii, Scalar &width)=0
Rn< double, 4 > rn
Definition StateSpace.h:595
Definition StateSpace.h:648
std::shared_ptr< Vpure > s4
Definition StateSpace.h:655
Eigen::Ref< Eigen::Matrix< double, -1, 1 > > ref_t
Definition StateSpace.h:653
double Scalar
Definition StateSpace.h:650
Scalar distance(cref_t &x, cref_t &y) const
Definition StateSpace.h:674
void sample_uniform(ref_t x) const
Definition StateSpace.h:663
const Eigen::Ref< const Eigen::Matrix< double, -1, 1 > > & cref_t
Definition StateSpace.h:652
void set_bounds(cref_t lb_, cref_t ub_)
Definition StateSpace.h:657
void choose_split_dimension(cref_t lb, cref_t ub, int &ii, Scalar &width)
Definition StateSpace.h:665
void interpolate(cref_t from, cref_t to, Scalar t, ref_t out) const
Definition StateSpace.h:659
Scalar distance_to_rectangle(cref_t &x, cref_t &lb, cref_t &ub) const
Definition StateSpace.h:670