1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244
|
// Copyright (c) 2005 Rijksuniversiteit Groningen (Netherlands)
// 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: svn+ssh://scm.gforge.inria.fr/svn/cgal/branches/CGAL-3.5-branch/Skin_surface_3/include/CGAL/Triangulation_incremental_builder_3.h $
// $Id: Triangulation_incremental_builder_3.h 49057 2009-04-30 14:03:52Z spion $
//
//
// Author(s) : Nico Kruithof <Nico@cs.rug.nl>
#ifndef CGAL_TDS_INCREMENTAL_BUILDER_3_H
#define CGAL_TDS_INCREMENTAL_BUILDER_3_H 1
#include <CGAL/basic.h>
#include <CGAL/Triangulation_data_structure_3.h>
#include <CGAL/array.h>
#include <set>
#include <list>
#include <vector>
CGAL_BEGIN_NAMESPACE
template < class Triangulation_3 >
class Triangulation_incremental_builder_3 {
public:
typedef Triangulation_3 T;
typedef typename T::Vertex_handle Vertex_handle;
typedef typename T::Cell_handle Cell_handle;
typedef typename T::Facet Facet;
public:
Triangulation_incremental_builder_3( T &t, bool verbose = false)
: t(t), m_verbose(verbose) {
}
~Triangulation_incremental_builder_3() {}
void begin_triangulation(int dim) {
t.clear();
t.tds().delete_cell(t.infinite_vertex()->cell());
// t.infinite = add_vertex();
t.tds().set_dimension(dim);
}
void end_triangulation() {
construct_infinite_cells();
CGAL_assertion(t.infinite_vertex()->cell() != Cell_handle());
}
Vertex_handle add_vertex();
Cell_handle add_cell(Vertex_handle vh0, Vertex_handle vh1,
Vertex_handle vh2, Vertex_handle vh3);
private:
void construct_infinite_cells();
Cell_handle add_infinite_cell(Cell_handle ch, int i);
void glue_cells(Cell_handle ch0, int ind0, Cell_handle ch1, int ind1);
// Interior facets of the simplical cell:
typedef std::pair < Vertex_handle, Vertex_handle > Vpair;
typedef std::map < Vpair, Facet > MapPair;
typedef typename MapPair::iterator MapPairIt;
typedef cpp0x::array < Vertex_handle, 3 > Vtriple;
typedef std::map < Vtriple, Facet > MapTriple;
typedef typename MapTriple::iterator MapTripleIt;
Vtriple facet(Vertex_handle vh1, Vertex_handle vh2, Vertex_handle vh3) {
if (vh1 < vh2) {
if (vh2 < vh3) {
return CGAL::make_array(vh1,vh2,vh3);
} else {
if (vh1 < vh3) {
return CGAL::make_array(vh1,vh3,vh2);
} else {
return CGAL::make_array(vh3,vh1,vh2);
}
}
}
if (vh1 < vh3) {
return CGAL::make_array(vh2,vh1,vh3);
} else {
if (vh2 < vh3) {
return CGAL::make_array(vh2,vh3,vh1);
} else {
return CGAL::make_array(vh3,vh2,vh1);
}
}
}
MapTriple facets;
T &t;
bool m_verbose;
};
template < class TDS_>
typename Triangulation_incremental_builder_3< TDS_ >::Vertex_handle
Triangulation_incremental_builder_3< TDS_ >::add_vertex() {
return t.tds().create_vertex();
}
template < class TDS_>
typename Triangulation_incremental_builder_3< TDS_ >::Cell_handle
Triangulation_incremental_builder_3< TDS_ >::add_cell(
Vertex_handle vh0, Vertex_handle vh1, Vertex_handle vh2, Vertex_handle vh3)
{
CGAL_assertion(vh0 != NULL); CGAL_assertion(vh1 != NULL);
CGAL_assertion(vh2 != NULL); CGAL_assertion(vh3 != NULL);
CGAL_assertion(vh0 != vh1); CGAL_assertion(vh0 != vh2); CGAL_assertion(vh0 != vh3);
CGAL_assertion(vh1 != vh2); CGAL_assertion(vh1 != vh3); CGAL_assertion(vh2 != vh3);
Cell_handle ch = t.tds().create_cell(vh0, vh1, vh2, vh3);
// Neighbors are by default set to NULL
vh0->set_cell(ch); vh1->set_cell(ch);
vh2->set_cell(ch); vh3->set_cell(ch);
MapTripleIt neighbIt;
for (int i=0; i<4; i++) {
Vtriple vtriple=facet(
ch->vertex((i+1)&3),
ch->vertex((i+2)&3),
ch->vertex((i+3)&3));
neighbIt = facets.find(vtriple);
if (neighbIt != facets.end()) {
Facet f = (*neighbIt).second;
glue_cells(f.first, f.second, ch, i);
facets.erase(neighbIt);
CGAL_assertion(f.first->neighbor(f.second) != NULL);
CGAL_assertion(ch->neighbor(i) != NULL);
} else {
facets[vtriple] = Facet(ch, i);
CGAL_assertion(ch->neighbor(i) == NULL);
}
}
return ch;
}
template < class TDS_>
typename Triangulation_incremental_builder_3< TDS_ >::Cell_handle
Triangulation_incremental_builder_3< TDS_ >::add_infinite_cell(
Cell_handle ch0, int i)
{
CGAL_assertion(ch0->neighbor(i) == NULL);
Vertex_handle vh[4];
vh[i] = t.infinite_vertex();
vh[(i+1)&3] = ch0->vertex((i+1)&3);
vh[(i+2)&3] = ch0->vertex((i+3)&3);
vh[(i+3)&3] = ch0->vertex((i+2)&3);
Cell_handle ch1 = t.tds().create_cell(vh[0], vh[1], vh[2], vh[3]);
// Neighbors are set to NULL
// Do not set points to the infinite cell. All finite vertices point to
// finite cells.
vh[i]->set_cell(ch1);
glue_cells(ch0, i, ch1, i);
return ch1;
}
template < class TDS_>
void
Triangulation_incremental_builder_3< TDS_ >::glue_cells(
Cell_handle ch0, int ind0, Cell_handle ch1, int ind1)
{
CGAL_assertion(ch0->has_vertex(ch1->vertex((ind1+1)&3)));
CGAL_assertion(ch0->has_vertex(ch1->vertex((ind1+2)&3)));
CGAL_assertion(ch0->has_vertex(ch1->vertex((ind1+3)&3)));
CGAL_assertion(ch1->has_vertex(ch0->vertex((ind0+1)&3)));
CGAL_assertion(ch1->has_vertex(ch0->vertex((ind0+2)&3)));
CGAL_assertion(ch1->has_vertex(ch0->vertex((ind0+3)&3)));
CGAL_assertion(ch0->index(ch1->vertex((ind1+1)&3)) != ind0);
CGAL_assertion(ch0->index(ch1->vertex((ind1+2)&3)) != ind0);
CGAL_assertion(ch0->index(ch1->vertex((ind1+3)&3)) != ind0);
CGAL_assertion(ch1->index(ch0->vertex((ind0+1)&3)) != ind1);
CGAL_assertion(ch1->index(ch0->vertex((ind0+2)&3)) != ind1);
CGAL_assertion(ch1->index(ch0->vertex((ind0+3)&3)) != ind1);
ch0->set_neighbor(ind0, ch1);
ch1->set_neighbor(ind1, ch0);
}
// Adds infinite cells to the facets on the convex hull
template < class TDS_>
void
Triangulation_incremental_builder_3< TDS_ >::construct_infinite_cells() {
MapTripleIt ch_facet_it;
MapPair ch_edges;
MapPairIt ch_edge_it;
Vertex_handle vh1, vh2;
for (ch_facet_it = facets.begin();
ch_facet_it != facets.end();
ch_facet_it ++) {
Facet ch_facet = (*ch_facet_it).second;
Cell_handle ch0 = ch_facet.first;
int ind0 = ch_facet.second;
Cell_handle ch1 = add_infinite_cell(ch0, ind0);
// Index of ch1 is also ind0
CGAL_assertion(ch0->neighbor(ind0) != NULL);
CGAL_assertion(ch1->neighbor(ind0) != NULL);
for (int i=1; i<4; i++) {
int i1 = (i==1?2:1);
int i2 = (i==3?2:3);
if (ch1->vertex((ind0+i1)&3) < ch1->vertex((ind0+i2)&3)) {
vh1 = ch1->vertex((ind0+i1)&3);
vh2 = ch1->vertex((ind0+i2)&3);
} else {
vh1 = ch1->vertex((ind0+i2)&3);
vh2 = ch1->vertex((ind0+i1)&3);
}
ch_edge_it = ch_edges.find(Vpair(vh1,vh2));
if (ch_edge_it != ch_edges.end()) {
Facet f_opp = (*ch_edge_it).second;
glue_cells(f_opp.first, f_opp.second, ch1, (ind0+i)&3);
ch_edges.erase(ch_edge_it);
CGAL_assertion(f_opp.first->neighbor(f_opp.second) != NULL);
CGAL_assertion(ch1->neighbor((ind0+i)&3) != NULL);
} else {
ch_edges[Vpair(vh1,vh2)] = Facet(ch1, (ind0+i)&3);
CGAL_assertion(ch1->neighbor((ind0+i)&3) == NULL);
}
}
}
CGAL_assertion(ch_edges.empty());
}
CGAL_END_NAMESPACE
#endif // CGAL_TDS_INCREMENTAL_BUILDER_3_H //
|