mirror of https://github.com/CGAL/cgal
Change after review 1
https://cgal.geometryfactory.com/CGAL/Members/wiki/Features/Optimal_bounding_box/1st_round/Answer_to_review_1
This commit is contained in:
parent
a9b1701980
commit
ae19de584a
|
|
@ -2,24 +2,24 @@
|
|||
\ingroup PkgOptimalBoundingBoxConcepts
|
||||
\cgalConcept
|
||||
|
||||
The concept `OrientedBoundingBoxTraits` describes the requirements of the traits class
|
||||
The concept `OrientedBoundingBoxTraits_3` describes the requirements of the traits class
|
||||
used in the function `CGAL::oriented_bounding_box()`, and in particular the need for
|
||||
a 3x3 matrix type.
|
||||
|
||||
\cgalHasModel `CGAL::Oriented_bounding_box_traits`
|
||||
\cgalHasModel `CGAL::Oriented_bounding_box_traits_3`
|
||||
|
||||
*/
|
||||
class OrientedBoundingBoxTraits
|
||||
class OrientedBoundingBoxTraits_3
|
||||
{
|
||||
public:
|
||||
/// The field number type; must be a model of the concept `FieldNumberType`.
|
||||
/// The field number type; must be a model of the concept `FieldNumberType`
|
||||
typedef unspecified_type FT;
|
||||
|
||||
/// The 3D point type; must be model of `Point_3`
|
||||
typedef unspecified_type Point_3;
|
||||
|
||||
/// The 3D affine transformation type; the template parameter `K` must be a model of `Kernel`
|
||||
/// and be compatible with the type `Point_3`.
|
||||
/// and be compatible with the type `Point_3`
|
||||
typedef CGAL::Aff_transformation_3<K> Aff_transformation_3;
|
||||
|
||||
/// A construction object that must provide the function operator:
|
||||
|
|
@ -28,9 +28,12 @@ public:
|
|||
typedef unspecified_type Construct_bbox_3;
|
||||
|
||||
/// A 3x3 matrix type; model of the concept `SvdTraits::Matrix` and which supports
|
||||
/// matrix-matrix and scalar-matrix multiplication, as well as matrix-matrix addition.
|
||||
/// matrix-matrix and scalar-matrix multiplication, as well as matrix-matrix addition
|
||||
typedef unspecified_type Matrix;
|
||||
|
||||
/// Returns the unitary matrix `Q` obtained in the QR-decomposition of the matrix `m`.
|
||||
/// Returns the unitary matrix `Q` obtained in the QR-decomposition of the matrix `m`
|
||||
Matrix get_Q(const Matrix& m) const;
|
||||
|
||||
/// Gives access to the `Construct_bbox_3` construction object
|
||||
Construct_bbox_3 construct_bbox_3_object() const { return Construct_bbox_3(); }
|
||||
};
|
||||
|
|
|
|||
|
|
@ -40,11 +40,11 @@ that are used in this package.
|
|||
|
||||
\cgalCRPSection{Concepts}
|
||||
|
||||
- `OrientedBoundingBoxTraits`
|
||||
- `OrientedBoundingBoxTraits_3`
|
||||
|
||||
\cgalCRPSection{Classes}
|
||||
|
||||
- `CGAL::Oriented_bounding_box_traits`
|
||||
- `CGAL::Oriented_bounding_box_traits_3`
|
||||
|
||||
\cgalCRPSection{Methods}
|
||||
|
||||
|
|
|
|||
|
|
@ -28,17 +28,17 @@ namespace CGAL {
|
|||
|
||||
/// \ingroup PkgOptimalBoundingBoxClasses
|
||||
///
|
||||
/// The class `CGAL::Oriented_bounding_box_traits` is a traits type
|
||||
/// to be used with the functions `CGAL::oriented_bounding_box()`.
|
||||
/// The class `CGAL::Oriented_bounding_box_traits_3` is a traits type
|
||||
/// to be used with the overloads of the function `CGAL::oriented_bounding_box()`.
|
||||
/// It uses the third party library \ref thirdpartyEigen "Eigen", which must therefore
|
||||
/// be available on the system for this class to be used.
|
||||
///
|
||||
/// \tparam K must be a model of `Kernel`
|
||||
///
|
||||
/// \cgalModels `OrientedBoundingBoxTraits`
|
||||
/// \cgalModels `OrientedBoundingBoxTraits_3`
|
||||
///
|
||||
template <typename K>
|
||||
class Oriented_bounding_box_traits
|
||||
class Oriented_bounding_box_traits_3
|
||||
{
|
||||
public:
|
||||
/// The field number type
|
||||
|
|
@ -60,7 +60,7 @@ private:
|
|||
typedef typename Matrix::EigenType EigenType;
|
||||
|
||||
public:
|
||||
/// Offers `construct_bbox_3_object()(const Point_3&)`
|
||||
/// Returns a default-constructed construction object
|
||||
Construct_bbox_3 construct_bbox_3_object() const { return Construct_bbox_3(); }
|
||||
|
||||
/// Performs QR-decomposition of matrix `m` to a unitary matrix and an upper triagonal
|
||||
|
|
@ -75,7 +75,7 @@ public:
|
|||
std::generate(group4.begin(), group4.end(), [&]{ return m_rng.get_int(0, im); });
|
||||
|
||||
// crossover I
|
||||
double bias = 0.1;
|
||||
FT bias = FT(0.1);
|
||||
|
||||
std::vector<Simplex> new_simplices(m);
|
||||
|
||||
|
|
@ -84,10 +84,10 @@ public:
|
|||
std::array<Matrix, 4> offspring;
|
||||
for(int j=0; j<4; ++j)
|
||||
{
|
||||
const double r = m_rng.get_double();
|
||||
const double fitnessA = compute_fitness<Traits>(m_population[group1[i]][j], m_points);
|
||||
const double fitnessB = compute_fitness<Traits>(m_population[group2[i]][j], m_points);
|
||||
const double threshold = (fitnessA < fitnessB) ? (0.5 + bias) : (0.5 - bias);
|
||||
const FT r = FT(m_rng.get_double());
|
||||
const FT fitnessA = compute_fitness<Traits>(m_population[group1[i]][j], m_points);
|
||||
const FT fitnessB = compute_fitness<Traits>(m_population[group2[i]][j], m_points);
|
||||
const FT threshold = (fitnessA < fitnessB) ? (0.5 + bias) : (0.5 - bias);
|
||||
|
||||
if(r < threshold)
|
||||
offspring[j] = m_population[group1[i]][j];
|
||||
|
|
@ -106,10 +106,10 @@ public:
|
|||
std::array<Matrix, 4> offspring;
|
||||
for(int j=0; j<4; ++j)
|
||||
{
|
||||
const double fitnessA = compute_fitness<Traits>(m_population[group3[i]][j], m_points);
|
||||
const double fitnessB = compute_fitness<Traits>(m_population[group4[i]][j], m_points);
|
||||
const double lambda = (fitnessA < fitnessB) ? (0.5 + bias) : (0.5 - bias);
|
||||
const double rambda = 1 - lambda; // the 'l' in 'lambda' stands for left
|
||||
const FT fitnessA = compute_fitness<Traits>(m_population[group3[i]][j], m_points);
|
||||
const FT fitnessB = compute_fitness<Traits>(m_population[group4[i]][j], m_points);
|
||||
const FT lambda = (fitnessA < fitnessB) ? (0.5 + bias) : (0.5 - bias);
|
||||
const FT rambda = 1 - lambda; // the 'l' in 'lambda' stands for left
|
||||
|
||||
// combine information from A and B
|
||||
Matrix new_vertex(3, 3);
|
||||
|
|
@ -148,9 +148,9 @@ public:
|
|||
const std::size_t nelder_mead_iterations = 150;
|
||||
|
||||
// stopping criteria prameters
|
||||
double prev_fit_value = 0.;
|
||||
double new_fit_value = 0.;
|
||||
const double tolerance = 1e-9;
|
||||
FT prev_fit_value = 0.;
|
||||
FT new_fit_value = 0.;
|
||||
const FT tolerance = 1e-9;
|
||||
int stale = 0;
|
||||
|
||||
for(std::size_t t=0; t<generations; ++t)
|
||||
|
|
@ -180,7 +180,7 @@ public:
|
|||
#endif
|
||||
|
||||
new_fit_value = fitness_map.get_best_fitness_value();
|
||||
const double difference = new_fit_value - prev_fit_value;
|
||||
const FT difference = new_fit_value - prev_fit_value;
|
||||
|
||||
#ifdef CGAL_OPTIMAL_BOUNDING_BOX_DEBUG
|
||||
const Matrix& R_now = fitness_map.get_best();
|
||||
|
|
|
|||
|
|
@ -50,7 +50,7 @@ Matrix mean(const Matrix& m1,
|
|||
{
|
||||
typedef typename Traits::FT FT;
|
||||
|
||||
const Matrix reduction = 0.5 * (m1 + m2);
|
||||
const Matrix reduction = FT(0.5) * (m1 + m2);
|
||||
|
||||
return traits.get_Q(reduction);
|
||||
}
|
||||
|
|
@ -63,7 +63,7 @@ const Matrix nm_centroid(const Matrix& S1,
|
|||
{
|
||||
typedef typename Traits::FT FT;
|
||||
|
||||
const Matrix mean = (1./3.) * (S1 + S2 + S3);
|
||||
const Matrix mean = (FT(1) / FT(3)) * (S1 + S2 + S3);
|
||||
|
||||
return traits.get_Q(mean);
|
||||
}
|
||||
|
|
|
|||
|
|
@ -17,7 +17,7 @@
|
|||
|
||||
#include <CGAL/Optimal_bounding_box/internal/evolution.h>
|
||||
#include <CGAL/Optimal_bounding_box/internal/population.h>
|
||||
#include <CGAL/Optimal_bounding_box/Oriented_bounding_box_traits.h>
|
||||
#include <CGAL/Optimal_bounding_box/Oriented_bounding_box_traits_3.h>
|
||||
|
||||
#include <CGAL/boost/graph/Named_function_parameters.h>
|
||||
#include <CGAL/boost/graph/named_params_helper.h>
|
||||
|
|
@ -58,7 +58,6 @@ void construct_oriented_bounding_box(std::array<typename Traits::Point_3, 8>& ob
|
|||
const PointRange& points,
|
||||
const Traits& traits)
|
||||
{
|
||||
typedef typename Traits::FT FT;
|
||||
typedef typename Traits::Point_3 Point;
|
||||
|
||||
// Construct the bbox of the transformed point set
|
||||
|
|
@ -209,12 +208,16 @@ void construct_oriented_bounding_box(Output& output,
|
|||
|
||||
/// \ingroup PkgOptimalBoundingBox_Oriented_bounding_box
|
||||
///
|
||||
/// See above.
|
||||
/// The function `oriented_bounding_box` computes an approximation of the <i>optimal bounding box</i>,
|
||||
/// which is defined as the rectangular box with smallest volume of all the rectangular boxes containing
|
||||
/// the input points.
|
||||
///
|
||||
/// \tparam PointRange a model of `Range`. The value type may not be equal to the traits' `Point_3` type
|
||||
/// See \ref PkgOptimalBoundingBox_Oriented_bounding_box for more information.
|
||||
///
|
||||
/// \tparam PointRange a model of `Range`. The value type may not be equal to the type `%Point_3` of the traits class
|
||||
/// if a point map is provided via named parameters (see below) to access points.
|
||||
/// \tparam Output either the traits' `Aff_transformation_3` type,
|
||||
/// or `std::array<Point, 8>` with `Point` being equivalent to the traits' `Point_3` type
|
||||
/// \tparam Output either the type `Aff_transformation_3` of the traits class,
|
||||
/// or `std::array<Point, 8>` with `Point` being equivalent to the type `%Point_3` of the traits class
|
||||
/// \tparam NamedParameters a sequence of \ref obb_namedparameters "Named Parameters"
|
||||
///
|
||||
/// \param points the input range
|
||||
|
|
@ -223,12 +226,12 @@ void construct_oriented_bounding_box(Output& output,
|
|||
///
|
||||
/// \cgalNamedParamsBegin
|
||||
/// \cgalParamBegin{point_map}
|
||||
/// a model of `ReadablePropertyMap` with value type `geom_traits::Point_3`. If this parameter is omitted,
|
||||
/// `CGAL::Identity_property_map<geom_traits::Point_3>` is used.
|
||||
/// a model of `ReadablePropertyMap` with value type the type`%Point_3` of the traits class.
|
||||
/// If this parameter is omitted, `CGAL::Identity_property_map<%Point_3>` is used.
|
||||
/// \cgalParamEnd
|
||||
/// \cgalParamBegin{geom_traits}
|
||||
/// a geometric traits class instance, model of the concept `OrientedBoundingBoxTraits`.
|
||||
/// %Default is `CGAL::Oriented_bounding_box_traits<K>` where `K` is deduced from the point type.
|
||||
/// a geometric traits class instance, model of the concept `OrientedBoundingBoxTraits_3`.
|
||||
/// %Default is `CGAL::Oriented_bounding_box_traits_3<K>` where `K` is deduced from the point type.
|
||||
/// \cgalParamEnd
|
||||
/// \cgalParamBegin{use_convex_hull}
|
||||
/// a Boolean value to indicate whether the algorithm should first extract the so-called extreme
|
||||
|
|
@ -256,7 +259,7 @@ void oriented_bounding_box(const PointRange& points,
|
|||
#if defined(CGAL_EIGEN3_ENABLED)
|
||||
typedef typename boost::range_value<PointRange>::type Point;
|
||||
typedef typename CGAL::Kernel_traits<Point>::type K;
|
||||
typedef Oriented_bounding_box_traits<K> Default_traits;
|
||||
typedef Oriented_bounding_box_traits_3<K> Default_traits;
|
||||
#else
|
||||
typedef void Default_traits;
|
||||
#endif
|
||||
|
|
@ -288,11 +291,11 @@ void oriented_bounding_box(const PointRange& points,
|
|||
|
||||
/// \ingroup PkgOptimalBoundingBox_Oriented_bounding_box
|
||||
///
|
||||
/// Extracts the vertices of the mesh as a point range and calls the overload above.
|
||||
/// Extracts the vertices of the mesh as a point range and calls the overload using points as input.
|
||||
///
|
||||
/// \tparam PolygonMesh a model of `VertexListGraph`
|
||||
/// \tparam Output either the traits' `Aff_transformation_3` type,
|
||||
/// or `std::array<Point, 8>` with `Point` being equivalent to the traits' `Point_3` type
|
||||
/// \tparam Output either the type `Aff_transformation_3` of the traits class,
|
||||
/// or `std::array<Point, 8>` with `Point` being equivalent to the type `%Point_3` of the traits class
|
||||
/// \tparam NamedParameters a sequence of \ref obb_namedparameters "Named Parameters"
|
||||
///
|
||||
/// \param pmesh the input mesh
|
||||
|
|
@ -306,8 +309,8 @@ void oriented_bounding_box(const PointRange& points,
|
|||
/// `CGAL::vertex_point_t` must be available in `PolygonMesh`
|
||||
/// \cgalParamEnd
|
||||
/// \cgalParamBegin{geom_traits}
|
||||
/// a geometric traits class instance, model of the concept `OrientedBoundingBoxTraits`.
|
||||
/// %Default is `CGAL::Oriented_bounding_box_traits<K>` where `K` is deduced from the point type.
|
||||
/// a geometric traits class instance, model of the concept `OrientedBoundingBoxTraits_3`.
|
||||
/// %Default is `CGAL::Oriented_bounding_box_traits_3<K>` where `K` is deduced from the point type.
|
||||
/// \cgalParamEnd
|
||||
/// \cgalParamBegin{use_convex_hull}
|
||||
/// a Boolean value to indicate whether the algorithm should first extract the so-called extreme
|
||||
|
|
|
|||
|
|
@ -24,7 +24,7 @@
|
|||
* Convenience header file including the headers for all the classes and functions of this package.
|
||||
*/
|
||||
|
||||
#include <CGAL/Optimal_bounding_box/Oriented_bounding_box_traits.h>
|
||||
#include <CGAL/Optimal_bounding_box/Oriented_bounding_box_traits_3.h>
|
||||
#include <CGAL/Optimal_bounding_box/oriented_bounding_box.h>
|
||||
|
||||
#endif // CGAL_OPTIMAL_BOUNDING_BOX_OBB_H
|
||||
|
|
|
|||
|
|
@ -7,7 +7,7 @@
|
|||
typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
|
||||
typedef K::Point_3 Point_3;
|
||||
|
||||
typedef CGAL::Oriented_bounding_box_traits<K> Traits;
|
||||
typedef CGAL::Oriented_bounding_box_traits_3<K> Traits;
|
||||
typedef Traits::Matrix Matrix;
|
||||
|
||||
bool are_equal(double d1, double d2, double epsilon)
|
||||
|
|
|
|||
Loading…
Reference in New Issue