Readability changes (no real changes)

This commit is contained in:
Mael Rouxel-Labbé 2018-01-11 11:42:45 +01:00
parent 4f02f892cb
commit 42c25d1c84
2 changed files with 133 additions and 126 deletions

View File

@ -45,7 +45,8 @@ struct Data_access
Data_access<Map>(const Map& m): map(m){}
std::pair< Data_type, bool>
operator()(const Key_type& p) const {
operator()(const Key_type& p) const
{
typename Map::const_iterator mit = map.find(p);
if(mit!= map.end())
return std::make_pair(mit->second, true);
@ -68,7 +69,8 @@ linear_interpolation(ForwardIterator first, ForwardIterator beyond,
typedef typename Functor::result_type::first_type Value_type;
Value_type result(0);
typename Functor::result_type val;
for(; first !=beyond; ++first){
for(; first !=beyond; ++first)
{
val = function_value(first->first);
CGAL_assertion(val.second);
result += (first->second/norm) * val.first;
@ -77,7 +79,8 @@ linear_interpolation(ForwardIterator first, ForwardIterator beyond,
}
template < class ForwardIterator, class Functor, class GradFunctor, class Traits, class Point>
template < class ForwardIterator, class Functor, class GradFunctor,
class Traits, class Point >
std::pair< typename Functor::result_type::first_type, bool>
quadratic_interpolation(ForwardIterator first, ForwardIterator beyond,
const typename
@ -89,20 +92,21 @@ quadratic_interpolation(ForwardIterator first, ForwardIterator beyond,
const Traits& traits)
{
CGAL_precondition(norm > 0);
Interpolation::internal::V2P<Traits> v2p(traits);
typedef typename Functor::result_type::first_type Value_type;
Interpolation::internal::V2P<Traits> v2p(traits);
Value_type result(0);
typename Functor::result_type f;
typename GradFunctor::result_type grad;
for(; first !=beyond; ++first){
for(; first !=beyond; ++first)
{
f = function_value(first->first);
grad = function_gradient(first->first);
//test if value and gradient are correctly retrieved:
CGAL_assertion(f.second);
if(!grad.second)
return std::make_pair(Value_type(0), false);
result += (first->second/norm)
*( f.first + grad.first*
result += (first->second/norm) * (f.first + grad.first *
traits.construct_scaled_vector_d_object()
(traits.construct_vector_d_object()(v2p(first->first), p),0.5));
}
@ -110,7 +114,8 @@ quadratic_interpolation(ForwardIterator first, ForwardIterator beyond,
}
template < class ForwardIterator, class Functor, class GradFunctor, class Traits, class Point>
template < class ForwardIterator, class Functor, class GradFunctor,
class Traits, class Point >
std::pair< typename Functor::result_type::first_type, bool>
sibson_c1_interpolation(ForwardIterator first, ForwardIterator beyond,
const typename
@ -122,7 +127,6 @@ sibson_c1_interpolation(ForwardIterator first, ForwardIterator beyond,
const Traits& traits)
{
CGAL_precondition(norm >0);
Interpolation::internal::V2P<Traits> v2p(traits);
typedef typename Functor::result_type::first_type Value_type;
typedef typename Traits::FT Coord_type;
@ -130,26 +134,28 @@ sibson_c1_interpolation(ForwardIterator first, ForwardIterator beyond,
Value_type linear_int(0), gradient_int(0);
typename Functor::result_type f;
typename GradFunctor::result_type grad;
Interpolation::internal::V2P<Traits> v2p(traits);
for(; first !=beyond; ++first){
for(; first !=beyond; ++first)
{
f = function_value(first->first);
grad = function_gradient(first->first);
CGAL_assertion(f.second);
if(!grad.second)
//the values are not correct:
return std::make_pair(Value_type(0), false);
return std::make_pair(Value_type(0), false); //the values are not correct
Coord_type coeff = first->second/norm;
Coord_type squared_dist = traits.
compute_squared_distance_d_object()(v2p(first->first), p);
Coord_type squared_dist = traits.compute_squared_distance_d_object()(v2p(first->first), p);
Coord_type dist = CGAL_NTS sqrt(squared_dist);
if(squared_dist ==0){
if(squared_dist == 0)
{
ForwardIterator it = first;
CGAL_USE(it);
CGAL_assertion(++it == beyond);
return std::make_pair(f.first, true);
}
//three different terms to mix linear and gradient
//interpolation
term1 += coeff / dist;
@ -158,8 +164,7 @@ sibson_c1_interpolation(ForwardIterator first, ForwardIterator beyond,
linear_int += coeff * f.first;
gradient_int += (coeff/dist)
* (f.first + grad.first *
gradient_int += (coeff/dist) * (f.first + grad.first *
traits.construct_vector_d_object()(v2p(first->first), p));
}
@ -170,20 +175,19 @@ sibson_c1_interpolation(ForwardIterator first, ForwardIterator beyond,
(term4 + term2), true);
}
//this method works with rational number types:
// This method works with rational number types:
// modification of Sibson's interpolant without sqrt
// following a proposition by Gunther Rote:
//
// the general scheme:
// The general scheme:
// Coord_type inv_weight = f(dist); //i.e. dist^2
// term1 += coeff/inv_weight;
// term2 += coeff * squared_dist;
// term3 += coeff*(squared_dist/inv_weight);
// gradient_int += (coeff/inv_weight)*
// (vh->get_value()+ vh->get_gradient()
// *(p - vh->point()));
// gradient_int += (coeff/inv_weight) * (vh->get_value()+ vh->get_gradient() * (p - vh->point()));
template < class ForwardIterator, class Functor, class GradFunctor, class Traits, class Point>
template < class ForwardIterator, class Functor, class GradFunctor,
class Traits, class Point >
std::pair< typename Functor::result_type::first_type, bool >
sibson_c1_interpolation_square(ForwardIterator first, ForwardIterator beyond,
const typename
@ -204,23 +208,26 @@ sibson_c1_interpolation_square(ForwardIterator first, ForwardIterator beyond,
typename Functor::result_type f;
typename GradFunctor::result_type grad;
for(; first!=beyond; ++first){
for(; first!=beyond; ++first)
{
f = function_value(first->first);
grad = function_gradient(first->first);
CGAL_assertion(f.second);
if(!grad.second)
//the gradient is not known
return std::make_pair(Value_type(0), false);
return std::make_pair(Value_type(0), false); // the gradient is not known
Coord_type coeff = first->second/norm;
Coord_type squared_dist = traits.compute_squared_distance_d_object()(v2p(first->first), v2p(p));
if(squared_dist ==0){
if(squared_dist ==0)
{
ForwardIterator it = first;
CGAL_USE(it);
CGAL_assertion(++it == beyond);
return std::make_pair(f.first, true);
}
//three different terms to mix linear and gradient
//interpolation
term1 += coeff / squared_dist;
@ -265,7 +272,8 @@ farin_c1_interpolation(RandomAccessIterator first,
typename GradFunctor::result_type grad;
int n = static_cast<int>(beyond - first);
if( n==1){
if(n == 1)
{
f = function_value(first->first);
CGAL_assertion(f.second);
return std::make_pair(f.first, true);
@ -279,10 +287,10 @@ farin_c1_interpolation(RandomAccessIterator first,
Value_type result(0);
const Coord_type fac3(3);
std::vector< std::vector<Value_type> >
ordinates(n,std::vector<Value_type>(n, Value_type(0)));
std::vector< std::vector<Value_type> > ordinates(n, std::vector<Value_type>(n, Value_type(0)));
for(int i =0; i<n; ++i){
for(int i=0; i<n; ++i)
{
it = first + i;
Coord_type coord_i_square = CGAL_NTS square(it->second);
@ -296,14 +304,15 @@ farin_c1_interpolation(RandomAccessIterator first,
// compute tangent plane control point (one 2, one 1 entry)
Value_type res_i(0);
for(int j =0; j<n; ++j){
if(i!=j){
for(int j=0; j<n; ++j)
{
if(i!=j)
{
it2 = first + j;
grad = function_gradient(it->first);
if(!grad.second)
//the gradient is not known
return std::make_pair(Value_type(0), false);
return std::make_pair(Value_type(0), false); // the gradient is not known
//ordinates[i][j] = (p_j - p_i) * g_i
ordinates[i][j] = grad.first *
@ -321,8 +330,11 @@ farin_c1_interpolation(RandomAccessIterator first,
// the third type of control points: three 1 entries i,j,k
for(int i=0; i< n; ++i)
{
for(int j=i+1; j< n; ++j)
for(int k=j+1; k<n; ++k){
{
for(int k=j+1; k<n; ++k)
{
// add 6* (u_i*u_j*u_k) * b_ijk
// b_ijk = 1.5 * a - 0.5*c
// where
@ -338,6 +350,8 @@ farin_c1_interpolation(RandomAccessIterator first,
ordinates[k][j]))
* (first+i)->second *(first+j)->second *(first+k)->second ;
}
}
}
return std::make_pair(result/(CGAL_NTS square(norm)*norm), true);
}

View File

@ -24,17 +24,19 @@
#include <CGAL/license/Interpolation.h>
#include <CGAL/Origin.h>
#include <CGAL/Interpolation/internal/helpers.h>
#include <CGAL/natural_neighbor_coordinates_2.h>
#include <CGAL/regular_neighbor_coordinates_2.h>
#include <iterator>
#include <utility>
#include <CGAL/Origin.h>
#include <boost/mpl/if.hpp>
#include <boost/type_traits/is_same.hpp>
#include <iterator>
#include <utility>
#include <vector>
namespace CGAL {
template < class ForwardIterator, class Functor, class Traits >
@ -42,22 +44,20 @@ typename Traits::Vector_d
sibson_gradient_fitting(ForwardIterator first,
ForwardIterator beyond,
const typename
std::iterator_traits<ForwardIterator>::
value_type::second_type& norm,
std::iterator_traits<ForwardIterator>::value_type::second_type& norm,
const typename Traits::Point_d& bare_p,
const typename Functor::result_type::first_type fn,
Functor function_value,
const Traits& traits)
{
Interpolation::internal::V2P<Traits> v2p(traits);
CGAL_precondition( first != beyond && norm != 0);
Interpolation::internal::V2P<Traits> v2p(traits);
typedef typename Traits::Aff_transformation_d Aff_transformation;
typedef typename Traits::FT Coord_type;
typename Traits::Vector_d pn =
traits.construct_vector_d_object()(NULL_VECTOR);
Aff_transformation scaling, m,
Hn(traits.construct_null_matrix_d_object()());
typename Traits::Vector_d pn = traits.construct_vector_d_object()(NULL_VECTOR);
Aff_transformation scaling, m, Hn(traits.construct_null_matrix_d_object()());
for(; first!=beyond; ++first)
{
@ -70,13 +70,11 @@ sibson_gradient_fitting(ForwardIterator first,
// compute the vector pn:
typename Functor::result_type f = function_value(first->first);
CGAL_assertion(f.second); //function value of first->first is valid
pn = pn + traits.construct_scaled_vector_d_object()
(d,scale * (f.first - fn));
pn = pn + traits.construct_scaled_vector_d_object()(d, scale * (f.first - fn));
//compute the matrix Hn:
m = traits.construct_outer_product_d_object()(d);
scaling = traits.construct_scaling_matrix_d_object()(scale);
Hn = traits.construct_sum_matrix_d_object()(Hn, scaling * m);
}
return Hn.inverse().transform(pn);
@ -88,9 +86,7 @@ typename Traits::Vector_d
sibson_gradient_fitting(const Triangul& tr,
ForwardIterator first,
ForwardIterator beyond,
const typename
std::iterator_traits<ForwardIterator>::
value_type::second_type& norm,
const typename std::iterator_traits<ForwardIterator>::value_type::second_type& norm,
VH vh,
Functor function_value,
const Traits& traits,
@ -116,8 +112,7 @@ sibson_gradient_fitting(const Triangul& tr,
ForwardIterator first,
ForwardIterator beyond,
const typename
std::iterator_traits<ForwardIterator>::
value_type::second_type& norm,
std::iterator_traits<ForwardIterator>::value_type::second_type& norm,
VH vh,
Functor function_value,
const Traits& traits,
@ -142,8 +137,7 @@ sibson_gradient_fitting(const Triangul& tr,
ForwardIterator first,
ForwardIterator beyond,
const typename
std::iterator_traits<ForwardIterator>::
value_type::second_type& norm,
std::iterator_traits<ForwardIterator>::value_type::second_type& norm,
VH vh,
Functor function_value,
const Traits& traits,
@ -205,7 +199,7 @@ sibson_gradient_fitting_internal(const Triangul& tr,
}
//the following functions allow to fit the gradients for all points in
// The following functions allow to fit the gradients for all points in
// a triangulation except the convex hull points.
// -> _nn2: natural_neighbor_coordinates_2
// -> _rn2: regular_neighbor_coordinates_2
@ -263,7 +257,6 @@ sibson_gradient_fitting_rn_2(const Rt& rt,
std::pair< typename Functor::argument_type,
typename Traits::FT > > > CoordInserter;
typedef typename boost::mpl::if_<
boost::is_same<typename Functor::argument_type, typename Rt::Weighted_point>,
Interpolation::internal::Vertex2WPoint<Rt, typename Traits::FT>,