mirror of https://github.com/CGAL/cgal
bugfixes (min/max, typos, proper boolean return for failure case)
This commit is contained in:
parent
ce78896f96
commit
93f6177906
|
|
@ -107,7 +107,7 @@ struct Cartesian_grid_location<GeomTraits, Cache_vertex_locations>
|
|||
for(std::size_t k=0; k<dims[2]; ++k)
|
||||
for(std::size_t j=0; j<dims[1]; ++j)
|
||||
for(std::size_t i=0; i<dims[0]; ++i)
|
||||
m_points.emplace_back(span.min() + Vector_3(i * spacing[0], j * spacing[1], k * spacing[2]));
|
||||
m_points.emplace_back((span.min)() + Vector_3(i * spacing[0], j * spacing[1], k * spacing[2]));
|
||||
}
|
||||
|
||||
template <typename Grid>
|
||||
|
|
|
|||
|
|
@ -153,28 +153,28 @@ template <typename K, typename S1, typename S2>
|
|||
typename K::FT
|
||||
shape_union(const S1& s1, const S2& s2, const typename K::Point_3& q)
|
||||
{
|
||||
return std::min(s1(q), s2(q));
|
||||
return (std::min)(s1(q), s2(q));
|
||||
}
|
||||
|
||||
template <typename K, typename S1, typename S2>
|
||||
typename K::FT
|
||||
shape_difference(const S1& s1, const S2& s2, const typename K::Point_3& q)
|
||||
{
|
||||
return std::max(s1(q), -s2(q));
|
||||
return (std::max)(s1(q), -s2(q));
|
||||
}
|
||||
|
||||
template <typename K, typename S1, typename S2>
|
||||
typename K::FT
|
||||
shape_intersection(const S1& s1, const S2& s2, const typename K::Point_3& q)
|
||||
{
|
||||
return std::max(s1(q), s2(q));
|
||||
return (std::max)(s1(q), s2(q));
|
||||
}
|
||||
|
||||
template <typename K, typename S1, typename S2>
|
||||
typename K::FT
|
||||
shape_symmetric_difference(const S1& s1, const S2& s2, const typename K::Point_3& q)
|
||||
{
|
||||
return std::max(-std::min(s1(q), s2(q)), std::max(s1(q), s2(q)));
|
||||
return (std::max)(-(std::min)(s1(q), s2(q)), (std::max)(s1(q), s2(q)));
|
||||
}
|
||||
|
||||
|
||||
|
|
|
|||
|
|
@ -313,17 +313,14 @@ private:
|
|||
d = b * b - FT(4) * a * c;
|
||||
}
|
||||
|
||||
void calc_coordinates(const std::array<FT, 8>& values, const std::size_t idx, const FT i0, const FT a, const FT b, const FT c, const FT d, const std::vector<bool> &f_flag, unsigned char &q_sol, FT ui[2], FT vi[2], FT wi[2]) {
|
||||
bool calc_coordinates(const std::array<FT, 8>& values, const std::size_t idx, const FT i0, const FT a, const FT b, const FT c, const FT d, const std::vector<bool> &f_flag, unsigned char &q_sol, FT ui[2], FT vi[2], FT wi[2]) {
|
||||
const int* remap = internal::Cube_table::asymptotic_remap[idx];
|
||||
|
||||
if (values[remap[0]] == values[remap[2]] && values[remap[1]] == values[remap[3]]) {
|
||||
ui[0] = FT(1) / FT(0);
|
||||
return;
|
||||
}
|
||||
if (values[remap[0]] == values[remap[4]] && values[remap[1]] == values[remap[5]]) {
|
||||
ui[0] = FT(1) / FT(0);
|
||||
return;
|
||||
}
|
||||
if (values[remap[0]] == values[remap[2]] && values[remap[1]] == values[remap[3]])
|
||||
return false;
|
||||
|
||||
if (values[remap[0]] == values[remap[4]] && values[remap[1]] == values[remap[5]])
|
||||
return false;
|
||||
|
||||
FT d2 = sqrt(CGAL::to_double(d));
|
||||
|
||||
|
|
@ -460,6 +457,8 @@ private:
|
|||
|
||||
if (0 < wi[1] && wi[1] < 1)
|
||||
q_sol |= 32;
|
||||
|
||||
return true;
|
||||
}
|
||||
|
||||
bool p_slice(const cell_descriptor& cell,
|
||||
|
|
@ -872,7 +871,8 @@ private:
|
|||
if (a == 0 || d < 0)
|
||||
continue;
|
||||
|
||||
calc_coordinates(values, idx, i0, a, b, c, d, f_flag, q_sol, ui, vi, wi);
|
||||
if (!calc_coordinates(values, idx, i0, a, b, c, d, f_flag, q_sol, ui, vi, wi))
|
||||
continue;
|
||||
|
||||
if (!std::isfinite(CGAL::to_double(ui[0])) || !std::isfinite(CGAL::to_double(ui[1])))
|
||||
continue;
|
||||
|
|
|
|||
|
|
@ -24,7 +24,7 @@
|
|||
#include <unordered_map>
|
||||
#include <unordered_set>
|
||||
|
||||
#define CGAL_TESTUISTE_ISOSURFACING_OUTPUT
|
||||
#define CGAL_TESTSUITE_ISOSURFACING_OUTPUT
|
||||
|
||||
namespace IS = CGAL::Isosurfacing;
|
||||
|
||||
|
|
@ -68,9 +68,9 @@ bool check_closed_not_empty(const Grid& grid, const IS::Interpolated_discrete_va
|
|||
for(std::size_t y=0; y<grid.ydim(); ++y) {
|
||||
for(std::size_t z=0; z<grid.zdim(); ++z)
|
||||
{
|
||||
total_min = std::min(total_min, values(x, y, z));
|
||||
total_min = (std::min)(total_min, values(x, y, z));
|
||||
if (x == 0 || y == 0 || z == 0 || x == grid.xdim() - 1 || y == grid.ydim() - 1 || z == grid.zdim() - 1)
|
||||
boundary_min = std::min(boundary_min, values(x, y, z));
|
||||
boundary_min = (std::min)(boundary_min, values(x, y, z));
|
||||
}
|
||||
}
|
||||
}
|
||||
|
|
@ -155,12 +155,12 @@ void embed_in_positive_grid(const Grid& grid, IS::Interpolated_discrete_values_3
|
|||
|
||||
CGAL::Random rand;
|
||||
|
||||
FT min = std::numeric_limits<FT>::max();
|
||||
FT min = (std::numeric_limits<FT>::max)();
|
||||
FT max = std::numeric_limits<FT>::lowest();
|
||||
for (const FT v : center_values)
|
||||
{
|
||||
min = std::min(min, v);
|
||||
max = std::max(max, v);
|
||||
min = (std::min)(min, v);
|
||||
max = (std::max)(max, v);
|
||||
}
|
||||
max += 1e-3;
|
||||
|
||||
|
|
@ -495,7 +495,7 @@ void compare_tmc_mc_trilinear(const std::array<typename KERNEL::FT, 8>& case_val
|
|||
Triangle_range triangles_high_res;
|
||||
IS::marching_cubes<CGAL::Parallel_if_available_tag>(domain_high_res, iso, points_high_res, triangles_high_res, CGAL::parameters::use_topologically_correct_marching_cubes(true));
|
||||
|
||||
#ifdef CGAL_TESTUISTE_ISOSURFACING_OUTPUT
|
||||
#ifdef CGAL_TESTSUITE_ISOSURFACING_OUTPUT
|
||||
CGAL::IO::write_polygon_soup("trilinear.off", points_high_res, triangles_high_res);
|
||||
#endif
|
||||
|
||||
|
|
@ -513,7 +513,7 @@ void compare_tmc_mc_trilinear(const std::array<typename KERNEL::FT, 8>& case_val
|
|||
Triangle_range triangles_mc;
|
||||
IS::marching_cubes<CGAL::Parallel_if_available_tag>(domain, iso, points_mc, triangles_mc, CGAL::parameters::use_topologically_correct_marching_cubes(false));
|
||||
|
||||
#ifdef CGAL_TESTUISTE_ISOSURFACING_OUTPUT
|
||||
#ifdef CGAL_TESTSUITE_ISOSURFACING_OUTPUT
|
||||
CGAL::IO::write_polygon_soup("mc.off", points_mc, triangles_mc);
|
||||
#endif
|
||||
}
|
||||
|
|
@ -522,7 +522,7 @@ void compare_tmc_mc_trilinear(const std::array<typename KERNEL::FT, 8>& case_val
|
|||
Triangle_range triangles;
|
||||
IS::marching_cubes<CGAL::Sequential_tag>(domain, iso, points, triangles, CGAL::parameters::use_topologically_correct_marching_cubes(true));
|
||||
|
||||
#ifdef CGAL_TESTUISTE_ISOSURFACING_OUTPUT
|
||||
#ifdef CGAL_TESTSUITE_ISOSURFACING_OUTPUT
|
||||
CGAL::IO::write_polygon_soup("tmc.off", points, triangles);
|
||||
#endif
|
||||
}
|
||||
|
|
|
|||
Loading…
Reference in New Issue