mirror of https://github.com/CGAL/cgal
Adding the simplex class to the Triangulation_3 package.
This commit is contained in:
parent
1fb9c68da2
commit
7ce6d89b94
|
|
@ -533,6 +533,15 @@ as the following example shows.
|
|||
|
||||
\ccIncludeExampleCode{Triangulation_3/example_adding_handles.cpp}
|
||||
|
||||
\subsection{The Simplex class}
|
||||
\label{Triangulation3-sec-simplex}
|
||||
The triangulation defines a \ccc{Simplex} class that represents a
|
||||
simplex (vertex, edge, facet or cell). This example demonstrates how
|
||||
simplices can be stored in a set.
|
||||
|
||||
\ccIncludeExampleCode{Triangulation_3/example_simplex.cpp}
|
||||
|
||||
|
||||
\subsection{Use of the Delaunay hierarchy}
|
||||
|
||||
\ccIncludeExampleCode{Triangulation_3/example_hierarchy.cpp}
|
||||
|
|
@ -592,5 +601,8 @@ In 2005, Christophe Delage implemented the vertex removal function for
|
|||
regular triangulations, which allowed to release this functionality in
|
||||
CGAL 3.2.
|
||||
|
||||
In 2006, Nico Kruithof wrote the \ccc{Triangulations_simplex_3} class
|
||||
that can store simplices of any dimension.
|
||||
|
||||
The authors also wish to thank Jean-Daniel Boissonnat for helpful
|
||||
discussions \cite{bdty-tcgal-00}.
|
||||
|
|
|
|||
|
|
@ -82,6 +82,9 @@ the user can pass them directly as arguments to the functions.
|
|||
\ccTypedef{typedef TriangulationDataStructure_3::Cell_handle Cell_handle;}
|
||||
{handle to a cell}
|
||||
|
||||
\ccTypedef{typedef Triangulation_simplex_3<Self> Simplex;}
|
||||
{Reference to a simplex (vertex, edge, facet or cell) of the triangulation}
|
||||
|
||||
\ccTypedef{typedef TriangulationDataStructure_3::size_type size_type;}
|
||||
{Size type (an unsigned integral type)}
|
||||
\ccGlue
|
||||
|
|
|
|||
|
|
@ -0,0 +1,126 @@
|
|||
% +------------------------------------------------------------------------+
|
||||
% | Reference manual page: Triangulation_simplex_3.tex
|
||||
% +------------------------------------------------------------------------+
|
||||
% | 27.09.2005 Nico Kruithof
|
||||
% | Package: Triangulation_3
|
||||
% |
|
||||
\RCSdef{\RCSTriangulationsimplexRev}{$Id$}
|
||||
\RCSdefDate{\RCSTriangulationsimplexDate}{$Date$}
|
||||
% +------------------------------------------------------------------------+
|
||||
|
||||
\ccRefPageBegin
|
||||
|
||||
%%RefPage: end of header, begin of main body
|
||||
% +------------------------------------------------------------------------+
|
||||
|
||||
|
||||
\begin{ccRefClass}{Triangulation_simplex_3<Triangulation_3>}
|
||||
\label{refTriangulationSimplex}
|
||||
|
||||
\ccDefinition The class \ccRefName\ stores a simplex of any dimension
|
||||
defined by the \ccc{Triangulation_3} class. It also defines the
|
||||
operator less such that simplices can be stored in a \ccc{map} or a
|
||||
\ccc{set} of simplices. The simplex is invalidated by any change in
|
||||
the triangulation.
|
||||
|
||||
\ccInclude{CGAL/Triangulation_simplex_3.h}
|
||||
|
||||
\ccParameters
|
||||
%
|
||||
It is parameterized by the triangulation it derives the simplices
|
||||
from.
|
||||
|
||||
\ccTypes
|
||||
\ccThree{typedef Triangulation_simplex_3<Triangulation_3> }{Triangulated}{}
|
||||
\ccThreeToTwo
|
||||
|
||||
\ccTypedef{typedef Triangulation_simplex_3<Triangulation_3> Simplex;}{
|
||||
The simplex class itself.}
|
||||
|
||||
\ccTypedef{typedef Triangulation_3::Vertex_handle Vertex_handle;}{}
|
||||
\ccGlue
|
||||
\ccTypedef{typedef Triangulation_3::Edge Edge;}{}
|
||||
\ccGlue
|
||||
\ccTypedef{typedef Triangulation_3::Facet Facet;}{}
|
||||
\ccGlue
|
||||
\ccTypedef{typedef Triangulation_3::Cell_handle Cell_handle;}{}
|
||||
|
||||
\ccTypedef{typedef Triangulation_3::Cell_circulator Cell_circulator;}{}
|
||||
\ccGlue
|
||||
\ccTypedef{typedef Triangulation_3::Facet_circulator Facet_circulator;}{}
|
||||
|
||||
\ccTypedef{typedef Triangulation_3::Edge_iterator Edge_iterator;}{}
|
||||
\ccGlue
|
||||
\ccTypedef{typedef Triangulation_3::Facet_iterator Facet_iterator;}{}
|
||||
|
||||
\ccTypedef{typedef Triangulation_3::Finite_vertices_iterator Finite_vertices_iterator;}{}
|
||||
\ccGlue
|
||||
\ccTypedef{typedef Triangulation_3::Finite_edges_iterator Finite_edges_iterator;}{}
|
||||
\ccGlue
|
||||
\ccTypedef{typedef Triangulation_3::Finite_facets_iterator Finite_facets_iterator;}{}
|
||||
\ccGlue
|
||||
\ccTypedef{typedef Triangulation_3::Finite_cells_iterator Finite_cells_iterator;}{}
|
||||
|
||||
|
||||
\ccCreation
|
||||
\ccCreationVariable{simplex} %% choose variable name
|
||||
\ccThree{Triangulation_simplex_3<Triangulation_3>}{simplex}{}
|
||||
\ccThreeToTwo
|
||||
|
||||
\ccConstructor{Triangulation_simplex_3();}{Initialises the simplex to
|
||||
an invalid simplex.}
|
||||
|
||||
\ccConstructor{Triangulation_simplex_3(Vertex_handle vh);}{}
|
||||
\ccGlue
|
||||
\ccConstructor{Triangulation_simplex_3(Edge e);}{}
|
||||
\ccGlue
|
||||
\ccConstructor{Triangulation_simplex_3(Facet f);}{}
|
||||
\ccGlue
|
||||
\ccConstructor{Triangulation_simplex_3(Cell_handle ch);}{}
|
||||
|
||||
\ccConstructor{Triangulation_simplex_3(Cell_circulator ccir);}{}
|
||||
\ccGlue
|
||||
\ccConstructor{Triangulation_simplex_3(Facet_circulator fcir);}{}
|
||||
|
||||
\ccConstructor{Triangulation_simplex_3(Edge_iterator eit);}{}
|
||||
\ccGlue
|
||||
\ccConstructor{Triangulation_simplex_3(Facet_iterator fit);}{}
|
||||
|
||||
\ccOperations
|
||||
\ccThree{Vertex_handle}{Vertex_handle(simplex)xx}{}
|
||||
\ccThreeToTwo
|
||||
|
||||
\ccMethod{int dimension () const;}{returns the dimension of the
|
||||
simplex.}
|
||||
%
|
||||
\ccMethod{operator Vertex_handle () const;}{Returns the \ccc{Vertex_handle}
|
||||
stored in the simplex. \ccPrecond{dimension() == 0}}
|
||||
\ccGlue
|
||||
\ccMethod{operator Edge () const;}{Returns the \ccc{Edge}
|
||||
stored in the simplex. \ccPrecond{dimension() == 1}}
|
||||
\ccGlue
|
||||
\ccMethod{operator Facet () const;}{Returns the \ccc{Facet}
|
||||
stored in the simplex. \ccPrecond{dimension() == 2}}
|
||||
\ccGlue
|
||||
\ccMethod{operator Cell_handle () const;}{Returns the \ccc{Cell_handle}
|
||||
stored in the simplex. \ccPrecond{dimension() == 3}}
|
||||
|
||||
\ccMethod{Cell_handle incident_cell () const;}{Returns a cell incident
|
||||
to the simplex.}
|
||||
|
||||
|
||||
\ccMethod{bool operator==(const
|
||||
Triangulation_simplex_3<Triangulation_3> &s1);}{Test whether two
|
||||
simplices are equal.}
|
||||
%
|
||||
\ccGlue
|
||||
%
|
||||
\ccMethod{bool operator< (const
|
||||
Triangulation_simplex_3<Triangulation_3> &s1);}{Defines a ordering
|
||||
on the simplices. This ordering depends on the memory layout and is
|
||||
independent of the geometry. Therefore, the ordering is not intrinsic}
|
||||
|
||||
\ccSeeAlso
|
||||
\ccRefIdfierPage{CGAL::Triangulation_3<TriangulationTraits_3,TriangulationDataStructure_3>}
|
||||
|
||||
\end{ccRefClass}
|
||||
|
|
@ -85,6 +85,8 @@ Figure~\ref{Triangulation3-fig-orient}.
|
|||
|
||||
\ccRefIdfierPage{CGAL::Regular_triangulation_cell_base_3<Traits,Cb>}
|
||||
|
||||
\ccRefIdfierPage{Triangulation_simplex_3<Triangulation_3>}
|
||||
|
||||
\subsubsection*{Traits Classes}
|
||||
|
||||
\ccRefIdfierPage{CGAL::Regular_triangulation_euclidean_traits_3<K,Weight>}\\
|
||||
|
|
|
|||
|
|
@ -21,6 +21,7 @@
|
|||
|
||||
\input{Triangulation_3_ref/TriangulationCellBase_3.tex}
|
||||
\input{Triangulation_3_ref/TriangulationVertexBase_3.tex}
|
||||
\input{Triangulation_3_ref/Triangulation_simplex_3.tex}
|
||||
\input{Triangulation_3_ref/TriangulationHierarchyVertexBase_3.tex}
|
||||
|
||||
\input{Triangulation_3_ref/RegularCellBase_3.tex}
|
||||
|
|
|
|||
|
|
@ -0,0 +1,56 @@
|
|||
// examples/Triangulation_3/example_simplex.C
|
||||
|
||||
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
|
||||
#include <CGAL/Triangulation_3.h>
|
||||
|
||||
#include <iostream>
|
||||
#include <fstream>
|
||||
#include <cassert>
|
||||
#include <list>
|
||||
#include <vector>
|
||||
|
||||
typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
|
||||
|
||||
typedef CGAL::Triangulation_3<K> Triangulation;
|
||||
|
||||
typedef Triangulation::Finite_vertices_iterator Finite_vertices_iterator;
|
||||
typedef Triangulation::Finite_edges_iterator Finite_edges_iterator;
|
||||
typedef Triangulation::Finite_facets_iterator Finite_facets_iterator;
|
||||
typedef Triangulation::Finite_cells_iterator Finite_cells_iterator;
|
||||
typedef Triangulation::Simplex Simplex;
|
||||
typedef Triangulation::Locate_type Locate_type;
|
||||
typedef Triangulation::Point Point;
|
||||
|
||||
int main()
|
||||
{
|
||||
// construction from a list of points :
|
||||
std::list<Point> L;
|
||||
L.push_front(Point(0,0,0));
|
||||
L.push_front(Point(1,0,0));
|
||||
L.push_front(Point(0,1,0));
|
||||
L.push_front(Point(0,1,1));
|
||||
|
||||
Triangulation T(L.begin(), L.end());
|
||||
|
||||
std::set<Simplex> simplices;
|
||||
|
||||
Finite_vertices_iterator vit = T.finite_vertices_begin();
|
||||
simplices.insert(Simplex(vit));
|
||||
|
||||
Finite_cells_iterator cit = T.finite_cells_begin();
|
||||
simplices.insert(Simplex(cit));
|
||||
|
||||
Finite_edges_iterator eit = T.finite_edges_begin();
|
||||
simplices.insert(Simplex(*eit));
|
||||
|
||||
Finite_facets_iterator fit = T.finite_facets_begin();
|
||||
simplices.insert(Simplex(*fit));
|
||||
|
||||
|
||||
for (std::set<Simplex>::iterator it = simplices.begin();
|
||||
it != simplices.end(); it++) {
|
||||
std::cout << it->dimension() << std::endl;
|
||||
}
|
||||
|
||||
return 0;
|
||||
}
|
||||
|
|
@ -101,6 +101,7 @@ public:
|
|||
typedef Edge_iterator All_edges_iterator;
|
||||
typedef Vertex_iterator All_vertices_iterator;
|
||||
|
||||
typedef typename Tds::Simplex Simplex;
|
||||
private:
|
||||
// This class is used to generate the Finite_*_iterators.
|
||||
class Infinite_tester
|
||||
|
|
|
|||
|
|
@ -48,6 +48,7 @@
|
|||
|
||||
#include <CGAL/Triangulation_ds_cell_base_3.h>
|
||||
#include <CGAL/Triangulation_ds_vertex_base_3.h>
|
||||
#include <CGAL/Triangulation_simplex_3.h>
|
||||
|
||||
#include <CGAL/Triangulation_ds_iterators_3.h>
|
||||
#include <CGAL/Triangulation_ds_circulators_3.h>
|
||||
|
|
@ -179,6 +180,8 @@ public:
|
|||
typedef std::pair<Cell_handle, int> Facet;
|
||||
typedef Triple<Cell_handle, int, int> Edge;
|
||||
|
||||
typedef Triangulation_simplex_3<Tds> Simplex;
|
||||
|
||||
public:
|
||||
Triangulation_data_structure_3()
|
||||
: _dimension(-2)
|
||||
|
|
|
|||
|
|
@ -0,0 +1,284 @@
|
|||
// Copyright (c) 2005 INRIA Sophia-Antipolis (France).
|
||||
// All rights reserved.
|
||||
//
|
||||
// This file is part of CGAL (www.cgal.org); you may redistribute it under
|
||||
// the terms of the Q Public License version 1.0.
|
||||
// See the file LICENSE.QPL distributed with CGAL.
|
||||
//
|
||||
// Licensees holding a valid commercial license may use this file in
|
||||
// accordance with the commercial license agreement provided with the software.
|
||||
//
|
||||
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
|
||||
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
|
||||
//
|
||||
// $URL$
|
||||
// $Id$
|
||||
//
|
||||
//
|
||||
// Author(s) : Nico Kruithof <Nico@cs.rug.nl>
|
||||
// Developed at Rijksuniversiteit Groningen (Netherlands)
|
||||
|
||||
#ifndef CGAL_TRIANGULATION_SIMPLEX_3_H
|
||||
#define CGAL_TRIANGULATION_SIMPLEX_3_H
|
||||
|
||||
CGAL_BEGIN_NAMESPACE
|
||||
|
||||
template < class TriangulationDataStructure_3 >
|
||||
class Triangulation_simplex_3 {
|
||||
typedef TriangulationDataStructure_3 TDS;
|
||||
typedef Triangulation_simplex_3<TDS> Self;
|
||||
public:
|
||||
typedef Self Simplex;
|
||||
|
||||
typedef typename TDS::Vertex_handle Vertex_handle;
|
||||
typedef typename TDS::Edge Edge;
|
||||
typedef typename TDS::Facet Facet;
|
||||
typedef typename TDS::Cell_handle Cell_handle;
|
||||
|
||||
typedef typename TDS::Cell_circulator Cell_circulator;
|
||||
typedef typename TDS::Facet_circulator Facet_circulator;
|
||||
|
||||
typedef typename TDS::Edge_iterator Edge_iterator;
|
||||
typedef typename TDS::Facet_iterator Facet_iterator;
|
||||
|
||||
// Constructors
|
||||
|
||||
// Default constructor initialises to undefined simplex:
|
||||
Triangulation_simplex_3() : ref(-1), ch() { }
|
||||
|
||||
Triangulation_simplex_3(Vertex_handle vh) {
|
||||
set_vertex(vh);
|
||||
}
|
||||
Triangulation_simplex_3(const Edge &e) {
|
||||
set_edge(e);
|
||||
}
|
||||
Triangulation_simplex_3(const Facet &f) {
|
||||
set_facet(f);
|
||||
}
|
||||
Triangulation_simplex_3(Cell_handle ch_) {
|
||||
set_cell(ch_);
|
||||
}
|
||||
|
||||
Triangulation_simplex_3(Cell_circulator ccir) {
|
||||
set_cell(ccir);
|
||||
}
|
||||
Triangulation_simplex_3(Facet_circulator fcir) {
|
||||
set_facet(*fcir);
|
||||
}
|
||||
|
||||
Triangulation_simplex_3(Edge_iterator eit) {
|
||||
set_edge(*eit);
|
||||
}
|
||||
Triangulation_simplex_3(Facet_iterator fit) {
|
||||
set_facet(*fit);
|
||||
}
|
||||
|
||||
// Conversions:
|
||||
operator Vertex_handle () const
|
||||
{
|
||||
CGAL_assertion(dimension() == 0);
|
||||
return ch->vertex(index(0));
|
||||
}
|
||||
operator Edge () const
|
||||
{
|
||||
CGAL_assertion(dimension() == 1);
|
||||
return Edge(ch,index(0),index(1));
|
||||
}
|
||||
operator Facet () const
|
||||
{
|
||||
CGAL_assertion(dimension() == 2);
|
||||
return Facet(ch,index(0));
|
||||
}
|
||||
operator Cell_handle () const
|
||||
{
|
||||
CGAL_assertion(dimension() == 3);
|
||||
return ch;
|
||||
}
|
||||
|
||||
// returns the dimension of the simplex
|
||||
int dimension () const {
|
||||
return (ref & 3);
|
||||
}
|
||||
// returns an incident cell:
|
||||
Cell_handle incident_cell() {
|
||||
return ch;
|
||||
}
|
||||
|
||||
template < class TDS2 >
|
||||
friend bool operator==(Triangulation_simplex_3<TDS2> s0,
|
||||
Triangulation_simplex_3<TDS2> s1);
|
||||
template < class TDS2 >
|
||||
friend bool operator< (Triangulation_simplex_3<TDS2> s0,
|
||||
Triangulation_simplex_3<TDS2> s1);
|
||||
|
||||
private:
|
||||
void set_vertex(const Vertex_handle vh) {
|
||||
ch = vh->cell();
|
||||
ref = (ch->index(vh) << 2); /* dim == 0 */
|
||||
CGAL_assertion (ch != Cell_handle());
|
||||
}
|
||||
void set_edge(const Edge &e) {
|
||||
ch = e.first;
|
||||
ref = (((e.third<< 2) + e.second) << 2) + 1; /* dim */
|
||||
CGAL_assertion (ch != Cell_handle());
|
||||
}
|
||||
void set_facet(const Facet &f) {
|
||||
ch = f.first;
|
||||
ref = (f.second << 2) + 2; /* dim */
|
||||
CGAL_assertion (ch != Cell_handle());
|
||||
}
|
||||
void set_cell(Cell_handle ch_) {
|
||||
ch = ch_;
|
||||
ref = 3; /* dim */
|
||||
CGAL_assertion (ch != Cell_handle());
|
||||
}
|
||||
|
||||
inline int index(int i) const {
|
||||
CGAL_assertion (i==0 || ((i==1) && (dimension()==1)));
|
||||
return (ref >> (2*(i+1))) & 3;
|
||||
}
|
||||
|
||||
int ref; // storage iijjdd (index i, index j, dimension of simplex)
|
||||
Cell_handle ch; // Corresponding cell handle
|
||||
};
|
||||
|
||||
///////////////////////////////
|
||||
// Simplex functions
|
||||
///////////////////////////////
|
||||
template < class TriangulationDataStructure_3 >
|
||||
bool
|
||||
operator!=(Triangulation_simplex_3<TriangulationDataStructure_3> s0,
|
||||
Triangulation_simplex_3<TriangulationDataStructure_3> s1) {
|
||||
return !(s0==s1);
|
||||
}
|
||||
|
||||
template < class TriangulationDataStructure_3 >
|
||||
bool
|
||||
operator==(Triangulation_simplex_3<TriangulationDataStructure_3> s0,
|
||||
Triangulation_simplex_3<TriangulationDataStructure_3> s1) {
|
||||
typedef Triangulation_simplex_3<TriangulationDataStructure_3> Sim;
|
||||
if (s0.dimension() != s1.dimension()) return false;
|
||||
|
||||
typename Sim::Cell_handle neighbor;
|
||||
|
||||
switch (s0.dimension()) {
|
||||
case (0): // Vertex
|
||||
return (s0.ch->vertex(s0.index(0)) == s1.ch->vertex(s1.index(0)));
|
||||
case (1): // Edge
|
||||
return ((s0.ch->vertex(s0.index(0)) == s1.ch->vertex(s1.index(0)) &&
|
||||
s0.ch->vertex(s0.index(1)) == s1.ch->vertex(s1.index(1))) ||
|
||||
(s0.ch->vertex(s0.index(1)) == s1.ch->vertex(s1.index(0)) &&
|
||||
s0.ch->vertex(s0.index(0)) == s1.ch->vertex(s1.index(1))));
|
||||
case (2):
|
||||
if (s0.ch == s1.ch && s0.index(0) == s1.index(0)) {
|
||||
return true;
|
||||
}
|
||||
|
||||
neighbor = s0.ch->neighbor(s0.index(0));
|
||||
if (neighbor == s1.ch &&
|
||||
neighbor->index(s0.ch) == s1.index(0)) {
|
||||
return true;
|
||||
}
|
||||
return false;
|
||||
case (3):
|
||||
return (&(*s0.ch) == &(*s1.ch));
|
||||
}
|
||||
CGAL_assertion(false);
|
||||
return false;
|
||||
}
|
||||
|
||||
template < class TriangulationDataStructure_3 >
|
||||
bool
|
||||
operator<(Triangulation_simplex_3<TriangulationDataStructure_3> s0,
|
||||
Triangulation_simplex_3<TriangulationDataStructure_3> s1) {
|
||||
typedef Triangulation_simplex_3<TriangulationDataStructure_3> Sim;
|
||||
|
||||
if (s0 == s1) return false;
|
||||
if (s0.dimension() < s1.dimension()) return true;
|
||||
if (s0.dimension() > s1.dimension()) return false;
|
||||
|
||||
// Dimensions are equal, compare the memory addresses of the simplices
|
||||
typename Sim::Cell_handle ch1, ch2;
|
||||
typename Sim::Vertex_handle vh1, vh2, vh3, vh4;
|
||||
switch (s0.dimension()) {
|
||||
case (0): // Vertex
|
||||
// Vertextices are not equal
|
||||
return (&(*s0.ch->vertex(s0.index(0))) <
|
||||
&(*s1.ch->vertex(s1.index(0))));
|
||||
case (1): // Edge
|
||||
vh1 = s0.ch->vertex(s0.index(0));
|
||||
vh2 = s0.ch->vertex(s0.index(1));
|
||||
vh3 = s1.ch->vertex(s1.index(0));
|
||||
vh4 = s1.ch->vertex(s1.index(1));
|
||||
|
||||
if ((std::min)(&(*vh1), &(*vh2)) < (std::min)(&(*vh3), &(*vh4)))
|
||||
return true;
|
||||
|
||||
if ((std::min)(&(*vh1), &(*vh2)) > (std::min)(&(*vh3), &(*vh4)))
|
||||
return false;
|
||||
|
||||
if ((std::max)(&(*vh1), &(*vh2)) < (std::max)(&(*vh3), &(*vh4)))
|
||||
return true;
|
||||
|
||||
return false;
|
||||
case (2): // Facet
|
||||
ch1 = s0.ch->neighbor(s0.index(0));
|
||||
ch2 = s1.ch->neighbor(s1.index(0));
|
||||
|
||||
if ((std::min)(&(*s0.ch), &(*ch1)) < (std::min)(&(*s1.ch), &(*ch2)))
|
||||
return true;
|
||||
|
||||
if ((std::min)(&(*s0.ch), &(*ch1)) > (std::min)(&(*s1.ch), &(*ch2)))
|
||||
return false;
|
||||
|
||||
if ((std::max)(&(*s0.ch), &(*ch1)) < (std::max)(&(*s1.ch), &(*ch2)))
|
||||
return true;
|
||||
|
||||
return false;
|
||||
case (3): // Cell
|
||||
return (&(*s0.ch) < &(*s1.ch));
|
||||
}
|
||||
CGAL_assertion(0);
|
||||
return false;
|
||||
}
|
||||
|
||||
template < class TriangulationDataStructure_3 >
|
||||
std::ostream &
|
||||
operator<< (std::ostream& os,
|
||||
const Triangulation_simplex_3<TriangulationDataStructure_3> &s)
|
||||
{
|
||||
typename TriangulationDataStructure_3::Vertex_handle vh;
|
||||
typename TriangulationDataStructure_3::Edge e;
|
||||
typename TriangulationDataStructure_3::Facet f;
|
||||
typename TriangulationDataStructure_3::Cell_handle ch;
|
||||
switch (s.dimension()) {
|
||||
case 0:
|
||||
vh = s;
|
||||
os << &*vh;
|
||||
break;
|
||||
case 1:
|
||||
e = s;
|
||||
os << &*(e.first->vertex(e.second)) << " "
|
||||
<< &*(e.first->vertex(e.third));
|
||||
break;
|
||||
case 2:
|
||||
f = s;
|
||||
os << &*(f.first->vertex((f.second+1)&3)) << " "
|
||||
<< &*(f.first->vertex((f.second+2)&3)) << " "
|
||||
<< &*(f.first->vertex((f.second+3)&3));
|
||||
break;
|
||||
case 3:
|
||||
ch = s;
|
||||
os << &*(ch->vertex(0)) << " "
|
||||
<< &*(ch->vertex(1)) << " "
|
||||
<< &*(ch->vertex(2)) << " "
|
||||
<< &*(ch->vertex(3));
|
||||
break;
|
||||
}
|
||||
return os;
|
||||
}
|
||||
|
||||
|
||||
CGAL_END_NAMESPACE
|
||||
|
||||
#endif // CGAL_TRIANGULATION_SIMPLEX_3_H
|
||||
|
|
@ -0,0 +1,147 @@
|
|||
// Copyright (c) 1998 INRIA Sophia-Antipolis (France).
|
||||
// All rights reserved.
|
||||
//
|
||||
// This file is part of CGAL (www.cgal.org); you may redistribute it under
|
||||
// the terms of the Q Public License version 1.0.
|
||||
// See the file LICENSE.QPL distributed with CGAL.
|
||||
//
|
||||
// Licensees holding a valid commercial license may use this file in
|
||||
// accordance with the commercial license agreement provided with the software.
|
||||
//
|
||||
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
|
||||
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
|
||||
//
|
||||
// $URL$
|
||||
// $Id$
|
||||
//
|
||||
//
|
||||
// Author(s) : Nico Kruithof
|
||||
|
||||
#include <cassert>
|
||||
#include <iostream>
|
||||
|
||||
template <class Triangulation>
|
||||
void
|
||||
_test_cls_triangulation_simplex_3(const Triangulation &)
|
||||
{
|
||||
typedef Triangulation Cls;
|
||||
|
||||
typedef typename Cls::Simplex Simplex;
|
||||
|
||||
typedef typename Cls::Point Point;
|
||||
|
||||
typedef typename Cls::Vertex Vertex;
|
||||
typedef typename Cls::Cell Cell;
|
||||
typedef typename Cls::Facet Facet;
|
||||
typedef typename Cls::Edge Edge;
|
||||
|
||||
typedef typename Cls::size_type size_type;
|
||||
typedef typename Cls::difference_type difference_type;
|
||||
|
||||
typedef typename Cls::Vertex_handle Vertex_handle;
|
||||
typedef typename Cls::Cell_handle Cell_handle;
|
||||
|
||||
typedef typename Cls::Cell_circulator Cell_circulator;
|
||||
typedef typename Cls::Facet_circulator Facet_circulator;
|
||||
|
||||
typedef typename Cls::Cell_iterator Cell_iterator;
|
||||
typedef typename Cls::Facet_iterator Facet_iterator;
|
||||
typedef typename Cls::Edge_iterator Edge_iterator;
|
||||
typedef typename Cls::Vertex_iterator Vertex_iterator;
|
||||
|
||||
typedef typename Cls::Finite_vertices_iterator Finite_vertices_iterator;
|
||||
typedef typename Cls::Finite_edges_iterator Finite_edges_iterator;
|
||||
typedef typename Cls::Finite_facets_iterator Finite_facets_iterator;
|
||||
typedef typename Cls::Finite_cells_iterator Finite_cells_iterator;
|
||||
|
||||
//########################################################################
|
||||
Cls t;
|
||||
|
||||
// Initialise to a 3D triangulation:
|
||||
t.insert(Point(0,0,0));
|
||||
t.insert(Point(1,0,0));
|
||||
t.insert(Point(0,1,0));
|
||||
t.insert(Point(0,0,1));
|
||||
|
||||
{ // Check vertices:
|
||||
Finite_vertices_iterator vit = t.finite_vertices_begin();
|
||||
Vertex_handle vh = vit;
|
||||
|
||||
Simplex s1 = vh;
|
||||
|
||||
Simplex s2(vit);
|
||||
Simplex s3(vh);
|
||||
|
||||
Simplex s4(s1);
|
||||
|
||||
CGAL_assertion(s1.dimension() == 0);
|
||||
CGAL_assertion(s1 == s2);
|
||||
CGAL_assertion(s1 == s3);
|
||||
CGAL_assertion(s1 == s4);
|
||||
Vertex_handle vh2 = s1;
|
||||
CGAL_assertion(vh == vh2);
|
||||
}
|
||||
|
||||
{ // Check edges
|
||||
Finite_edges_iterator eit = t.finite_edges_begin();
|
||||
Edge e = *eit;
|
||||
|
||||
Simplex s1 = *eit;
|
||||
Simplex s2 = e;
|
||||
|
||||
Simplex s3(*eit);
|
||||
Simplex s4(e);
|
||||
|
||||
Simplex s5(s1);
|
||||
|
||||
CGAL_assertion(s1.dimension() == 1);
|
||||
CGAL_assertion(s1 == s2);
|
||||
CGAL_assertion(s1 == s3);
|
||||
CGAL_assertion(s1 == s4);
|
||||
CGAL_assertion(s1 == s5);
|
||||
Edge e2 = s1;
|
||||
CGAL_assertion(e == e2);
|
||||
}
|
||||
|
||||
{ // Check facets
|
||||
Finite_facets_iterator fit = t.finite_facets_begin();
|
||||
Facet f = *fit;
|
||||
|
||||
Simplex s1 = *fit;
|
||||
Simplex s2 = f;
|
||||
|
||||
Simplex s3(*fit);
|
||||
Simplex s4(f);
|
||||
|
||||
Simplex s5(s1);
|
||||
|
||||
CGAL_assertion(s1.dimension() == 2);
|
||||
CGAL_assertion(s1 == s2);
|
||||
CGAL_assertion(s1 == s3);
|
||||
CGAL_assertion(s1 == s4);
|
||||
CGAL_assertion(s1 == s5);
|
||||
Facet f2 = s1;
|
||||
CGAL_assertion(f == f2);
|
||||
}
|
||||
|
||||
{ // Check cells
|
||||
Finite_cells_iterator cit = t.finite_cells_begin();
|
||||
Cell_handle ch = cit;
|
||||
|
||||
Simplex s1 = Simplex(cit);
|
||||
Simplex s2 = ch;
|
||||
|
||||
Simplex s3(cit);
|
||||
Simplex s4(ch);
|
||||
|
||||
Simplex s5(s1);
|
||||
|
||||
CGAL_assertion(s1.dimension() == 3);
|
||||
CGAL_assertion(s1 == s2);
|
||||
CGAL_assertion(s1 == s3);
|
||||
CGAL_assertion(s1 == s4);
|
||||
CGAL_assertion(s1 == s5);
|
||||
Cell_handle ch2 = s1;
|
||||
CGAL_assertion(ch == ch2);
|
||||
}
|
||||
}
|
||||
|
|
@ -0,0 +1,39 @@
|
|||
// Copyright (c) 1998 INRIA Sophia-Antipolis (France).
|
||||
// All rights reserved.
|
||||
//
|
||||
// This file is part of CGAL (www.cgal.org); you may redistribute it under
|
||||
// the terms of the Q Public License version 1.0.
|
||||
// See the file LICENSE.QPL distributed with CGAL.
|
||||
//
|
||||
// Licensees holding a valid commercial license may use this file in
|
||||
// accordance with the commercial license agreement provided with the software.
|
||||
//
|
||||
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
|
||||
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
|
||||
//
|
||||
// $URL$
|
||||
// $Id$
|
||||
//
|
||||
//
|
||||
// Author(s) : Nico Kruithof
|
||||
|
||||
//#define CGAL_TRIANGULATION_DONT_USE_SHORT_NAMES
|
||||
|
||||
#include <CGAL/Triangulation_3.h>
|
||||
|
||||
bool del = false;
|
||||
|
||||
#include <CGAL/_test_types.h>
|
||||
#include <CGAL/_test_cls_triangulation_simplex_3.C>
|
||||
|
||||
// Explicit instantiation of the whole class :
|
||||
template class CGAL::Triangulation_3<K>;
|
||||
|
||||
int main()
|
||||
{
|
||||
typedef CGAL::Triangulation_3<K> Cls3;
|
||||
|
||||
_test_cls_triangulation_simplex_3( Cls3() );
|
||||
|
||||
return 0;
|
||||
}
|
||||
Loading…
Reference in New Issue