mirror of https://github.com/CGAL/cgal
add another version of tetrahedron_soup_to_triangulation_3
based on indices, similarly to polygon_soup_to_polygon_mesh and using named parameters for optional parameters use it in an example
This commit is contained in:
parent
2bbdca32fb
commit
595c969757
|
|
@ -198,6 +198,9 @@ CGAL_add_named_parameter(cell_selector_t, cell_selector, cell_selector)
|
|||
CGAL_add_named_parameter(facet_is_constrained_t, facet_is_constrained, facet_is_constrained_map)
|
||||
CGAL_add_named_parameter(smooth_constrained_edges_t, smooth_constrained_edges, smooth_constrained_edges)
|
||||
|
||||
// MDS_3 parameters
|
||||
CGAL_add_named_parameter(surface_facets_t, surface_facets, surface_facets)
|
||||
|
||||
// output parameters
|
||||
CGAL_add_named_parameter(face_proxy_map_t, face_proxy_map, face_proxy_map)
|
||||
CGAL_add_named_parameter(proxies_t, proxies, proxies)
|
||||
|
|
|
|||
|
|
@ -8,35 +8,69 @@
|
|||
#include <CGAL/tetrahedron_soup_to_triangulation_3.h>
|
||||
|
||||
#include <vector>
|
||||
#include <map>
|
||||
|
||||
typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
|
||||
|
||||
typedef CGAL::Delaunay_triangulation_3<K> DT3;
|
||||
typedef CGAL::Tetrahedral_remeshing::Remeshing_triangulation_3<K> Remeshing_triangulation;
|
||||
typedef CGAL::Mesh_complex_3_in_triangulation_3<Remeshing_triangulation> C3T3;
|
||||
|
||||
typedef K::Point_3 Point_3;
|
||||
typedef K::Tetrahedron_3 Tetrahedron_3;
|
||||
|
||||
typedef DT3::Vertex_handle Vertex_handle;
|
||||
|
||||
int main(int argc, char* argv[])
|
||||
{
|
||||
//build a triangulation
|
||||
DT3 delaunay;
|
||||
CGAL::Random_points_in_cube_3<Point_3> randp(2.);
|
||||
while (delaunay.number_of_vertices() < 100)
|
||||
delaunay.insert(*randp++);
|
||||
const int nbv = 100;
|
||||
|
||||
//collect tetrahedra
|
||||
std::vector<Tetrahedron_3> tetrahedra(delaunay.number_of_finite_cells());
|
||||
//a triangulation
|
||||
DT3 delaunay;
|
||||
std::map<Vertex_handle, int> v2i;
|
||||
std::vector<DT3::Point> points(nbv);
|
||||
std::vector<Tetrahedron_3> tetrahedra;
|
||||
std::vector<std::array<int, 5> > tets_by_indices;
|
||||
|
||||
//insert random points
|
||||
CGAL::Random_points_in_cube_3<Point_3> randp(2.);
|
||||
int i = 0;
|
||||
while (i < nbv)
|
||||
{
|
||||
points[i] = *randp++;
|
||||
Vertex_handle v = delaunay.insert(points[i]);
|
||||
v2i[v] = i++;
|
||||
}
|
||||
|
||||
tetrahedra.reserve(delaunay.number_of_finite_cells());
|
||||
tets_by_indices.reserve(delaunay.number_of_finite_cells());
|
||||
for (DT3::Cell_handle c : delaunay.finite_cell_handles())
|
||||
{
|
||||
tetrahedra.push_back(delaunay.tetrahedron(c));
|
||||
|
||||
//build triangulation
|
||||
std::array<int, 5> tet;
|
||||
tet[0] = v2i.at(c->vertex(0));
|
||||
tet[1] = v2i.at(c->vertex(1));
|
||||
tet[2] = v2i.at(c->vertex(2));
|
||||
tet[3] = v2i.at(c->vertex(3));
|
||||
tet[4] = Remeshing_triangulation::Cell::Subdomain_index(1);
|
||||
|
||||
tets_by_indices.push_back(tet);
|
||||
}
|
||||
|
||||
//build triangulation from tetrahedra
|
||||
Remeshing_triangulation tr;
|
||||
CGAL::tetrahedron_soup_to_triangulation_3(tetrahedra, tr);
|
||||
|
||||
std::ofstream os("dt_rebuilt.mesh");
|
||||
CGAL::write_MEDIT(os, tr);
|
||||
//buid triangulation from indices
|
||||
Remeshing_triangulation tr2;
|
||||
CGAL::tetrahedron_soup_to_triangulation_3(points, tets_by_indices, tr2);
|
||||
|
||||
//build a C3T3
|
||||
C3T3 c3t3;
|
||||
c3t3.triangulation() = tr;
|
||||
|
||||
std::ofstream ofs("c3t3_output.mesh");
|
||||
c3t3.output_to_medit(ofs);
|
||||
|
||||
return EXIT_SUCCESS;
|
||||
}
|
||||
|
|
|
|||
|
|
@ -341,7 +341,6 @@ bool build_triangulation(Tr& tr,
|
|||
const std::vector<typename Tr::Point>& points,
|
||||
const std::vector<std::array<int,5> >& finite_cells,
|
||||
const std::map<std::array<int,3>, typename Tr::Cell::Surface_patch_index>& border_facets,
|
||||
std::vector<typename Tr::Vertex_handle>& vertex_handle_vector,
|
||||
const bool verbose = false,
|
||||
bool replace_domain_0 = false)
|
||||
{
|
||||
|
|
@ -354,6 +353,7 @@ bool build_triangulation(Tr& tr,
|
|||
typedef boost::unordered_map<Facet_vvv, std::vector<Incident_cell> > Incident_cells_map;
|
||||
|
||||
Incident_cells_map incident_cells_map;
|
||||
std::vector<Vertex_handle> vertex_handle_vector;
|
||||
vertex_handle_vector.resize(points.size() + 1); // id to vertex_handle
|
||||
//index 0 is for infinite vertex
|
||||
// 1 to n for points in `points`
|
||||
|
|
@ -484,9 +484,8 @@ bool build_triangulation_from_file(std::istream& is,
|
|||
if(finite_cells.empty())
|
||||
return false;
|
||||
|
||||
std::vector<typename Tr::Vertex_handle> vertices(points.size() + 1);
|
||||
bool is_well_built = build_triangulation(tr,
|
||||
points, finite_cells, border_facets, vertices, false, replace_domain_0);
|
||||
points, finite_cells, border_facets, false, replace_domain_0);
|
||||
return is_well_built;
|
||||
}
|
||||
|
||||
|
|
|
|||
|
|
@ -21,6 +21,7 @@
|
|||
#include <CGAL/license/MDS_3.h>
|
||||
|
||||
#include <CGAL/MDS_3/tet_soup_to_c3t3.h>
|
||||
#include <CGAL/boost/graph/named_params_helper.h>
|
||||
|
||||
#include <vector>
|
||||
#include <array>
|
||||
|
|
@ -28,10 +29,8 @@
|
|||
|
||||
namespace CGAL {
|
||||
|
||||
/*!
|
||||
\ingroup PkgMDS3Functions
|
||||
|
||||
* @brief
|
||||
/** \ingroup PkgMDS3Functions
|
||||
* builds a 3D triangulation from a soup of tetrahedra.
|
||||
*
|
||||
* @tparam TetrahedronRange
|
||||
* @tparam Triangulation model of
|
||||
|
|
@ -45,47 +44,99 @@ namespace CGAL {
|
|||
Triangulation& tr)
|
||||
{
|
||||
typedef Triangulation Tr;
|
||||
typedef typename Tr::Cell_handle Cell_handle;
|
||||
typedef typename Tr::Vertex_handle Vertex_handle;
|
||||
typedef typename Tr::Point Point;
|
||||
|
||||
std::vector<Point> points;
|
||||
std::vector<std::array<int, 5> > finite_cells;
|
||||
std::map<std::array<int, 3>, typename Tr::Cell::Surface_patch_index> border_facets;
|
||||
std::vector<Vertex_handle> vertex_handle_vector;
|
||||
std::map<Vertex_handle, int> v2i;
|
||||
std::map<Point, int> p2i;
|
||||
|
||||
for (typename TetrahedronRange::value_type tet : tets)
|
||||
{
|
||||
CGAL_assertion(tet.orientation() != CGAL::NEGATIVE);
|
||||
CGAL_assertion(tet.orientation() == CGAL::POSITIVE);
|
||||
std::array<int, 5> cell;
|
||||
|
||||
Cell_handle hint = Cell_handle();
|
||||
for (int i = 0; i < 4; ++i)
|
||||
{
|
||||
const Point& pi = tet[i];
|
||||
typename Tr::Locate_type lt;
|
||||
int li, lj;
|
||||
hint = tr.locate(pi, lt, li, lj, hint);
|
||||
if (lt != Tr::Locate_type::VERTEX)
|
||||
if (p2i.find(pi) == p2i.end())
|
||||
{
|
||||
points.push_back(pi);
|
||||
Vertex_handle newv = tr.insert(pi, lt, hint, li, lj);
|
||||
vertex_handle_vector.push_back(newv);
|
||||
|
||||
CGAL_assertion(points.size() == vertex_handle_vector.size());
|
||||
v2i.insert(std::make_pair(newv, points.size() - 1));
|
||||
p2i.insert(std::make_pair(pi, points.size() - 1));
|
||||
cell[i] = static_cast<int>(points.size() - 1);
|
||||
}
|
||||
else
|
||||
cell[i] = v2i.at(hint->vertex(li));
|
||||
cell[i] = p2i.at(pi);
|
||||
}
|
||||
cell[4] = 1;
|
||||
|
||||
CGAL_assertion(CGAL::orientation(points[cell[0]],
|
||||
points[cell[1]], points[cell[2]], points[cell[3]]) == CGAL::POSITIVE);
|
||||
|
||||
finite_cells.push_back(cell);
|
||||
}
|
||||
|
||||
CGAL::MDS_3::build_triangulation(tr, points, finite_cells,
|
||||
border_facets, vertex_handle_vector);
|
||||
CGAL::MDS_3::build_triangulation(tr, points, finite_cells, border_facets);
|
||||
}
|
||||
|
||||
/** \ingroup PkgMDS3Functions
|
||||
* builds a 3D triangulation from a soup of tetrahedra.
|
||||
*
|
||||
* @tparam PointRange a model of the concept `RandomAccessContainer`
|
||||
* whose value type is the point type
|
||||
* @tparam TetrahedronRange a model of the concept `RandomAccessContainer` whose
|
||||
* value type is a model of the concept `RandomAccessContainer` whose value type is `std::size_t`
|
||||
*
|
||||
* @param points points of the soup of tetrahedra
|
||||
* @param tets each element in the range describes a tetrahedron using the indices of the points
|
||||
* in `points` (indices 0 to 3), and the associated `Subdomain_index` (index 4)
|
||||
* @param tr the 3D triangulation to be built
|
||||
* @param np an optional sequence of \ref mds3_namedparameters "Named Parameters" among the ones listed below
|
||||
*
|
||||
* \cgalNamedParamsBegin
|
||||
* \cgalParamNBegin{surface_facets}
|
||||
* \cgalParamDescription{each element in the range describes a surface facet using the indices of points
|
||||
* in `points` (indices 0 to 2), and the associated `Surface_patch_index` (index 3)}
|
||||
* \cgalParamType{a class model of `AssociativeContainer`
|
||||
* whose key type is model of `RandomAccessContainer`
|
||||
* and mapped type is `Tr::Cell::Surface_patch_index`}
|
||||
* \cgalParamDefault{An empty `std::map<std::array<int, 3>, typename Tr::Cell::Surface_patch_index>`}
|
||||
* \cgalParamExtra{If this parameter is omitted, an internal property map for `CGAL::vertex_point_t`
|
||||
* must be available in `PolygonMesh`.}
|
||||
* \cgalParamNEnd
|
||||
* \cgalNamedParamsEnd
|
||||
*
|
||||
* @pre `points` contains each point only once
|
||||
*
|
||||
* @sa `CGAL::Polygon_mesh_processing::polygon_soup_to_polygon_mesh()`
|
||||
*/
|
||||
template<typename PointRange,
|
||||
typename TetrahedronRange,
|
||||
typename Triangulation,
|
||||
typename NamedParameters>
|
||||
void tetrahedron_soup_to_triangulation_3(const PointRange& points,
|
||||
const TetrahedronRange& tets,
|
||||
Triangulation& tr,
|
||||
const NamedParameters& np)
|
||||
{
|
||||
using parameters::choose_parameter;
|
||||
using parameters::get_parameter;
|
||||
|
||||
std::map<std::array<int, 3>,
|
||||
typename Triangulation::Cell::Surface_patch_index> empty_map;
|
||||
auto facets = choose_parameter(get_parameter(np, internal_np::surface_facets), empty_map);
|
||||
|
||||
CGAL::MDS_3::build_triangulation(tr, points, tets, facets);
|
||||
}
|
||||
|
||||
template<typename PointRange,
|
||||
typename TetrahedronRange,
|
||||
typename Triangulation>
|
||||
void tetrahedron_soup_to_triangulation_3(const PointRange& points,
|
||||
const TetrahedronRange& tets,
|
||||
Triangulation& tr)
|
||||
{
|
||||
tetrahedron_soup_to_triangulation_3(points, tets, tr, parameters::all_default());
|
||||
}
|
||||
|
||||
} //namespace CGAL
|
||||
|
|
|
|||
Loading…
Reference in New Issue