#include "include/utils.h"
#include "include/Saver.h"
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
using FT = typename Kernel::FT;
using Segments = std::vector<Segment_2>;
using Indices = std::vector<std::size_t>;
using Neighbor_query =
using Angle_regularization =
using Offset_regularization =
using Quadratic_program =
using Quadratic_angle_regularizer =
Kernel, Segments, Neighbor_query, Angle_regularization, Quadratic_program>;
using Quadratic_offset_regularizer =
Kernel, Segments, Neighbor_query, Offset_regularization, Quadratic_program>;
int main(int argc, char *argv[]) {
std::string path = "";
if (argc > 1) path = argv[1];
Saver<Kernel> saver;
std::vector<Segment_2> segments;
create_example_15(segments);
std::vector<Indices> groups(3);
groups[0] = {0, 1, 2, 3, 4, 5, 6};
groups[1] = {7, 8, 9, 10};
groups[2] = {11, 12, 13, 14};
if (path != "") {
const std::string full_path = path + "regularize_15_segments_before";
saver.export_eps_segments(segments, full_path, FT(100));
}
const FT max_angle_2 = FT(10);
Quadratic_program qp_angles;
Neighbor_query neighbor_query(segments);
Angle_regularization angle_regularization(
segments, CGAL::parameters::maximum_angle(max_angle_2));
for (const auto& group : groups) {
neighbor_query.add_group(group);
angle_regularization.add_group(group);
}
Quadratic_angle_regularizer qp_angle_regularizer(
segments, neighbor_query, angle_regularization, qp_angles);
qp_angle_regularizer.regularize();
std::cout << "* number of modified segments (angles) = " <<
angle_regularization.number_of_modified_segments() << std::endl;
const FT max_offset_2 = FT(1) / FT(10);
std::vector<Indices> pgroups;
angle_regularization.parallel_groups(
std::back_inserter(pgroups));
Quadratic_program qp_offsets;
Offset_regularization offset_regularization(
segments, CGAL::parameters::maximum_offset(max_offset_2));
neighbor_query.clear();
for (const auto& pgroup : pgroups) {
neighbor_query.add_group(pgroup);
offset_regularization.add_group(pgroup);
}
Quadratic_offset_regularizer qp_offset_regularizer(
segments, neighbor_query, offset_regularization, qp_offsets);
qp_offset_regularizer.regularize();
std::cout << "* number of modified segments (offsets) = " <<
offset_regularization.number_of_modified_segments() << std::endl;
if (path != "") {
const std::string full_path = path + "regularize_15_segments_after";
saver.export_eps_segments(segments, full_path, FT(100));
}
}
A convenience header that includes all free functions and classes for shape regularization.
Shape regularization algorithm based on the quadratic programming global optimization.
Definition: QP_regularization.h:70
An angle-based regularization type for 2D segments that reinforces parallelism and orthogonality rela...
Definition: Angle_regularization_2.h:55
A neighbor query based on a Delaunay triangulation, which enables to find the nearest neighbors in a ...
Definition: Delaunay_neighbor_query_2.h:61
An offset-based regularization type for 2D segments that reinforces collinearity relationships.
Definition: Offset_regularization_2.h:60