use a matrix with dynamic rows and columns are run time for points, because it is a little faster than fixed columns

This commit is contained in:
Konstantinos Katrioplas 2018-05-18 16:22:01 +02:00
parent 820e9c6098
commit 0c215dfecb
4 changed files with 16 additions and 7 deletions

View File

@ -33,7 +33,7 @@ namespace Optimal_bounding_box {
template <typename Linear_algebra_traits>
class Evolution
{
typedef typename Linear_algebra_traits::MatrixX3d MatrixXd;
typedef typename Linear_algebra_traits::MatrixXd MatrixXd;
typedef typename Linear_algebra_traits::Matrix3d Matrix3d;
typedef typename Linear_algebra_traits::Vector3d Vector3d;

View File

@ -105,7 +105,7 @@ void find_obb(std::vector<Point>& points,
CGAL_assertion(obb_points.size() == 8);
// eigen linear algebra traits
typedef typename LinearAlgebraTraits::MatrixX3d MatrixXd;
typedef typename LinearAlgebraTraits::MatrixXd MatrixXd;
typedef typename LinearAlgebraTraits::Matrix3d Matrix3d;
MatrixXd points_mat;

View File

@ -20,7 +20,7 @@ void test_nelder_mead()
{
typedef CGAL::Eigen_linear_algebra_traits Linear_algebra_traits;
typedef Linear_algebra_traits::Matrix3d Matrix3d;
typedef Linear_algebra_traits::MatrixX3d MatrixXd;
typedef Linear_algebra_traits::MatrixXd MatrixXd;
MatrixXd data_points(4,3);
data_points(0,0) = 0.866802;
@ -143,7 +143,7 @@ void test_nelder_mead()
void test_genetic_algorithm()
{
CGAL::Eigen_dense_matrix<double, -1, 3> data_points(4, 3); // -1 : dynamic size at run time
CGAL::Eigen_dense_matrix<double, -1, -1> data_points(4, 3); // -1 : dynamic size at run time
data_points(0,0) = 0.866802;
data_points(0,1) = 0.740808,
data_points(0,2) = 0.895304,
@ -170,7 +170,7 @@ void test_genetic_algorithm()
void test_random_unit_tetra()
{
// this is dynamic at run times
CGAL::Eigen_dense_matrix<double, -1, 3> data_points(4, 3);
CGAL::Eigen_dense_matrix<double, -1, -1> data_points(4, 3);
// points are on their convex hull
data_points(0,0) = 0.866802;
@ -244,7 +244,7 @@ void test_reference_tetrahedron(const char* fname)
}
typedef CGAL::Eigen_linear_algebra_traits Linear_algebra_traits;
typedef Linear_algebra_traits::MatrixX3d MatrixXd;
typedef Linear_algebra_traits::MatrixXd MatrixXd;
typedef Linear_algebra_traits::Matrix3d Matrix3d;
// points in a matrix
@ -279,7 +279,7 @@ void test_long_tetrahedron(std::string fname)
}
typedef CGAL::Eigen_linear_algebra_traits Linear_algebra_traits;
typedef Linear_algebra_traits::MatrixX3d MatrixXd;
typedef Linear_algebra_traits::MatrixXd MatrixXd;
typedef Linear_algebra_traits::Matrix3d Matrix3d;
// points in a matrix

View File

@ -174,6 +174,15 @@ const CGAL::Eigen_dense_matrix<NT, D1, D3> operator* (const CGAL::Eigen_dense_ma
return CGAL::Eigen_dense_matrix<NT, D1, D3>(A.m_matrix * B.m_matrix);
}
// D2 and D3 may not be equal at compile time, but equal at run time!
// This overload returns a dynamic matrix.
template <class NT, int D1, int D2, int D3, int D4>
const CGAL::Eigen_dense_matrix<NT> operator* (const CGAL::Eigen_dense_matrix<NT, D1, D2 >& A,
const CGAL::Eigen_dense_matrix<NT, D3, D4 >& B)
{
return CGAL::Eigen_dense_matrix<NT>(A.m_matrix * B.m_matrix);
}
// scalar - matrix multiplication
template <class NT, int D1, int D2>
const CGAL::Eigen_dense_matrix<NT, D1, D2> operator* (const NT& scalar,