ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/math/AlphaShape.cpp
Revision: 1242
Committed: Wed May 14 14:31:48 2008 UTC (16 years, 11 months ago) by chuckv
File size: 7967 byte(s)
Log Message:
Changes....

File Contents

# Content
1 /* Copyright (c) 2007 The University of Notre Dame. All Rights Reserved.
2 *
3 * The University of Notre Dame grants you ("Licensee") a
4 * non-exclusive, royalty free, license to use, modify and
5 * redistribute this software in source and binary code form, provided
6 * that the following conditions are met:
7 *
8 * 1. Acknowledgement of the program authors must be made in any
9 * publication of scientific results based in part on use of the
10 * program. An acceptable form of acknowledgement is citation of
11 * the article in which the program was described (Matthew
12 * A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
13 * J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
14 * Parallel Simulation Engine for Molecular Dynamics,"
15 * J. Comput. Chem. 26, pp. 252-271 (2005))
16 *
17 * 2. Redistributions of source code must retain the above copyright
18 * notice, this list of conditions and the following disclaimer.
19 *
20 * 3. Redistributions in binary form must reproduce the above copyright
21 * notice, this list of conditions and the following disclaimer in the
22 * documentation and/or other materials provided with the
23 * distribution.
24 *
25 * This software is provided "AS IS," without a warranty of any
26 * kind. All express or implied conditions, representations and
27 * warranties, including any implied warranty of merchantability,
28 * fitness for a particular purpose or non-infringement, are hereby
29 * excluded. The University of Notre Dame and its licensors shall not
30 * be liable for any damages suffered by licensee as a result of
31 * using, modifying or distributing the software or its
32 * derivatives. In no event will the University of Notre Dame or its
33 * licensors be liable for any lost revenue, profit or data, or for
34 * direct, indirect, special, consequential, incidental or punitive
35 * damages, however caused and regardless of the theory of liability,
36 * arising out of the use of or inability to use software, even if the
37 * University of Notre Dame has been advised of the possibility of
38 * such damages.
39 *
40 *
41 * AlphaShape.cpp
42 *
43 * Purpose: To calculate convexhull, hull volume and radius
44 * using the CGAL library.
45 *
46 * Created by Charles F. Vardeman II on 11 Dec 2006.
47 * @author Charles F. Vardeman II
48 * @version $Id: AlphaShape.cpp,v 1.2 2008-05-14 14:31:48 chuckv Exp $
49 *
50 */
51
52 #include "math/AlphaShape.hpp"
53
54 #include <iostream>
55 #include <list>
56 #include <fstream>
57 #include <CGAL/IO/Geomview_stream.h>
58 #include <CGAL/IO/alpha_shape_geomview_ostream_3.h>
59
60 #include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
61 #include <CGAL/Delaunay_triangulation_3.h>
62 #include <CGAL/Triangulation_hierarchy_3.h>
63 #include <CGAL/Alpha_shape_3.h>
64 #include <CGAL/Tetrahedron_3.h>
65
66
67 struct K : CGAL::Exact_predicates_inexact_constructions_kernel {};
68
69 typedef CGAL::Alpha_shape_vertex_base_3<K> Vb;
70 typedef CGAL::Triangulation_hierarchy_vertex_base_3<Vb> Vbh;
71 typedef CGAL::Alpha_shape_cell_base_3<K> Fb;
72 typedef CGAL::Triangulation_data_structure_3<Vbh,Fb> Tds;
73 typedef CGAL::Delaunay_triangulation_3<K,Tds> Delaunay;
74 typedef CGAL::Triangulation_hierarchy_3<Delaunay> Delaunay_hierarchy;
75 typedef CGAL::Alpha_shape_3<Delaunay_hierarchy> Alpha_shape_3;
76
77 typedef K::Point_3 Point;
78 typedef Alpha_shape_3::Alpha_iterator Alpha_iterator;
79 typedef Alpha_shape_3::NT NT;
80 typedef Alpha_shape_3::Cell_handle Cell_handle;
81 typedef Alpha_shape_3::Vertex_handle Vertex_handle;
82 typedef Alpha_shape_3::Facet Facet;
83 typedef Alpha_shape_3::Edge Edge;
84 typedef CGAL::Tetrahedron_3<K> Tetrahedron;
85
86
87
88
89 using namespace oopse;
90
91 AlphaShape::AlphaShape(){}
92
93
94 bool AlphaShape::genHull(std::vector<Vector3d> pos)
95 {
96 Delaunay_hierarchy dt;
97 //points.reserve(pos.size());
98 // Copy the positon vector into a points vector for cgal.
99 for (unsigned int i = 0; i < pos.size(); ++i)
100 {
101 Point pt(pos[i][0],pos[i][1],pos[i][2]);
102 dt.insert(pt);
103 }
104
105 /* Generate Alpha Shape */
106 std::cout << "Generating alpha shape" << std::endl;
107 Alpha_shape_3 ashape(dt);
108
109 /*
110 CGAL::Geomview_stream gv(CGAL::Bbox_3(0,0,0, 2, 2, 2));
111 gv.set_line_width(4);
112 gv.set_trace(false);
113 gv.set_bg_color(CGAL::Color(0, 200, 200));
114 gv.clear();
115 */
116 Alpha_shape_3::NT alpha_solid = ashape.find_alpha_solid();
117 Alpha_iterator opt = ashape.find_optimal_alpha(1);
118 std::cout << "Smallest alpha value to get a solid through data points is "
119 << alpha_solid << std::endl;
120 std::cout << "Optimal alpha value to get one connected component is "
121 << *opt << std::endl;
122 ashape.set_alpha(*opt);
123 assert(ashape.number_of_solid_components() == 1);
124 std::list<Cell_handle> cells;
125 std::list<Facet> facets;
126 std::list<Edge> edges;
127 std::list<Vertex_handle> vertices;
128
129 ashape.get_alpha_shape_cells(std::back_inserter(cells),
130 Alpha_shape_3::INTERIOR);
131 ashape.get_alpha_shape_facets(std::back_inserter(facets),
132 Alpha_shape_3::REGULAR);
133 ashape.get_alpha_shape_facets(std::back_inserter(facets),
134 Alpha_shape_3::SINGULAR);
135 ashape.get_alpha_shape_edges(std::back_inserter(edges),
136 Alpha_shape_3::SINGULAR);
137 ashape.get_alpha_shape_vertices(std::back_inserter(vertices),
138 Alpha_shape_3::REGULAR);
139 std::cout << " The 0-shape has : " << std::endl;
140 std::cout << cells.size() << " interior tetrahedra" << std::endl;
141 std::cout << facets.size() << " boundary facets" << std::endl;
142 std::cout << edges.size() << " singular edges" << std::endl;
143 std::cout << vertices.size() << " singular vertices" << std::endl;
144
145 RealType volume_;
146 std::list<Cell_handle>::iterator thiscell;
147
148 for(Alpha_shape_3::Cell_iterator c_it = ashape.cells_begin(); c_it != ashape.cells_end(); ++c_it)
149 {
150 volume_ += ashape.tetrahedron(c_it).volume();
151 }
152
153 //gv << (Delaunay) ashape;
154 //std::cout << ashape;
155
156 }
157
158 RealType AlphaShape::getVolume()
159 {
160 /*
161 std::list<Point> L;
162 L.push_front(Point(0,0,0));
163 L.push_front(Point(1,0,0));
164 L.push_front(Point(0,1,0));
165
166 Triangulation T(L.begin(), L.end());
167
168 int n = T.number_of_vertices();
169
170 // insertion from a vector :
171 std::vector<Point> V(3);
172 V[0] = Point(0,0,1);
173 V[1] = Point(1,1,1);
174 V[2] = Point(2,2,2);
175
176 n = n + T.insert(V.begin(), V.end());
177
178 assert( n == 6 ); // 6 points have been inserted
179 assert( T.is_valid() ); // checking validity of T
180
181 double sum_v = 0;
182 for(Triangulation::Cell_iterator c_it = T.cells_begin(); c_it != T.cells_end(); ++c_it)
183 {
184 sum_v += T.tetrahedron(c_it).volume();
185 }
186 std::cout << "sum_v " << sum_v << std::endl;
187 */
188 return 0.0;
189 }
190
191 void AlphaShape::geomviewHull(const std::string& geomFileName)
192 {
193 /*
194 std::ofstream newGeomFile;
195
196 //create new .md file based on old .md file
197 newGeomFile.open(geomFileName.c_str());
198
199 // Write polyhedron in Object File Format (OFF).
200 CGAL::set_ascii_mode( std::cout);
201 newGeomFile << "OFF" << std::endl << ch_polyhedron.size_of_vertices() << ' '
202 << ch_polyhedron.size_of_facets() << " 0" << std::endl;
203 std::copy( ch_polyhedron.points_begin(), ch_polyhedron.points_end(),
204 std::ostream_iterator<Point_3>( newGeomFile, "\n"));
205 for ( Facet_iterator i = ch_polyhedron.facets_begin(); i != ch_polyhedron.facets_end(); ++i)
206 {
207 Halfedge_facet_circulator j = i->facet_begin();
208 // Facets in polyhedral surfaces are at least triangles.
209 CGAL_assertion( CGAL::circulator_size(j) >= 3);
210 newGeomFile << CGAL::circulator_size(j) << ' ';
211 do
212 {
213 newGeomFile << ' ' << std::distance(ch_polyhedron.vertices_begin(), j->vertex());
214 }
215 while ( ++j != i->facet_begin());
216 newGeomFile << std::endl;
217 }
218
219 newGeomFile.close();
220 */
221
222 }