ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/math/AlphaHull.cpp
Revision: 1850
Committed: Wed Feb 20 15:39:39 2013 UTC (12 years, 2 months ago) by gezelter
File size: 16466 byte(s)
Log Message:
Fixed a widespread typo in the license 

File Contents

# Content
1 /* Copyright (c) 2008, 2009, 2010 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. Redistributions of source code must retain the above copyright
9 * notice, this list of conditions and the following disclaimer.
10 *
11 * 2. Redistributions in binary form must reproduce the above copyright
12 * notice, this list of conditions and the following disclaimer in the
13 * documentation and/or other materials provided with the
14 * distribution.
15 *
16 * This software is provided "AS IS," without a warranty of any
17 * kind. All express or implied conditions, representations and
18 * warranties, including any implied warranty of merchantability,
19 * fitness for a particular purpose or non-infringement, are hereby
20 * excluded. The University of Notre Dame and its licensors shall not
21 * be liable for any damages suffered by licensee as a result of
22 * using, modifying or distributing the software or its
23 * derivatives. In no event will the University of Notre Dame or its
24 * licensors be liable for any lost revenue, profit or data, or for
25 * direct, indirect, special, consequential, incidental or punitive
26 * damages, however caused and regardless of the theory of liability,
27 * arising out of the use of or inability to use software, even if the
28 * University of Notre Dame has been advised of the possibility of
29 * such damages.
30 *
31 * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
32 * research, please cite the appropriate papers when you publish your
33 * work. Good starting points are:
34 *
35 * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
36 * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
37 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
38 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
39 * [4] , Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011). *
40 *
41 * AlphaHull.cpp
42 *
43 * Purpose: To calculate Alpha hull, hull volume libqhull.
44 *
45 * Created by Charles F. Vardeman II on 11 Dec 2006.
46 * @author Charles F. Vardeman II
47 * @version $Id$
48 *
49 */
50
51 /* Standard includes independent of library */
52
53 #include <iostream>
54 #include <fstream>
55 #include <list>
56 #include <algorithm>
57 #include <iterator>
58 #include <utility>
59 #include "math/AlphaHull.hpp"
60 #include "utils/simError.h"
61 #ifdef IS_MPI
62 #include <mpi.h>
63 #endif
64 #include "math/qhull.hpp"
65
66 using namespace OpenMD;
67
68 #ifdef HAVE_QHULL
69 double calculate_circumradius(pointT* p0,pointT* p1,pointT* p2, int dim);
70
71 AlphaHull::AlphaHull(double alpha) : Hull(), dim_(3), alpha_(alpha), options_("qhull d QJ Tcv Pp") {
72 }
73
74 void AlphaHull::computeHull(std::vector<StuntDouble*> bodydoubles) {
75
76 int numpoints = bodydoubles.size();
77 // bool alphashape=true;
78
79 Triangles_.clear();
80
81 vertexT *vertex;
82 // vertexT **vertexp;
83 facetT *facet, *neighbor;
84 // setT *vertices, *verticestop, *verticesbottom;
85 int curlong, totlong;
86 pointT *interiorPoint;
87
88 std::vector<double> ptArray(numpoints*dim_);
89
90 // Copy the positon vector into a points vector for qhull.
91 std::vector<StuntDouble*>::iterator SD;
92 int i = 0;
93 for (SD =bodydoubles.begin(); SD != bodydoubles.end(); ++SD){
94 Vector3d pos = (*SD)->getPos();
95 ptArray[dim_ * i] = pos.x();
96 ptArray[dim_ * i + 1] = pos.y();
97 ptArray[dim_ * i + 2] = pos.z();
98 i++;
99 }
100
101 /* Clean up memory from previous convex hull calculations*/
102 boolT ismalloc = False;
103
104 int ridgesCount=0;
105 if (qh_new_qhull(dim_, numpoints, &ptArray[0], ismalloc,
106 const_cast<char *>(options_.c_str()), NULL, stderr)) {
107
108 sprintf(painCave.errMsg, "AlphaHull: Qhull failed to compute convex hull");
109 painCave.isFatal = 1;
110 simError();
111
112 } //qh_new_qhull
113
114
115 #ifdef IS_MPI
116 //If we are doing the mpi version, set up some vectors for data communication
117
118 int nproc = MPI::COMM_WORLD.Get_size();
119 int myrank = MPI::COMM_WORLD.Get_rank();
120 int localHullSites = 0;
121
122 std::vector<int> hullSitesOnProc(nproc, 0);
123 std::vector<int> coordsOnProc(nproc, 0);
124 std::vector<int> displacements(nproc, 0);
125 std::vector<int> vectorDisplacements(nproc, 0);
126
127 std::vector<double> coords;
128 std::vector<double> vels;
129 std::vector<int> indexMap;
130 std::vector<double> masses;
131
132 FORALLvertices{
133 localHullSites++;
134
135 int idx = qh_pointid(vertex->point);
136
137 indexMap.push_back(idx);
138
139 coords.push_back(ptArray[dim_ * idx]);
140 coords.push_back(ptArray[dim_ * idx + 1]);
141 coords.push_back(ptArray[dim_ * idx + 2]);
142
143 StuntDouble* sd = bodydoubles[idx];
144
145 Vector3d vel = sd->getVel();
146 vels.push_back(vel.x());
147 vels.push_back(vel.y());
148 vels.push_back(vel.z());
149
150 masses.push_back(sd->getMass());
151 }
152
153 MPI::COMM_WORLD.Allgather(&localHullSites, 1, MPI::INT, &hullSitesOnProc[0],
154 1, MPI::INT);
155
156 int globalHullSites = 0;
157 for (int iproc = 0; iproc < nproc; iproc++){
158 globalHullSites += hullSitesOnProc[iproc];
159 coordsOnProc[iproc] = dim_ * hullSitesOnProc[iproc];
160 }
161
162 displacements[0] = 0;
163 vectorDisplacements[0] = 0;
164
165 for (int iproc = 1; iproc < nproc; iproc++){
166 displacements[iproc] = displacements[iproc-1] + hullSitesOnProc[iproc-1];
167 vectorDisplacements[iproc] = vectorDisplacements[iproc-1] + coordsOnProc[iproc-1];
168 }
169
170 std::vector<double> globalCoords(dim_ * globalHullSites);
171 std::vector<double> globalVels(dim_ * globalHullSites);
172 std::vector<double> globalMasses(globalHullSites);
173
174 int count = coordsOnProc[myrank];
175
176 MPI::COMM_WORLD.Allgatherv(&coords[0], count, MPI::DOUBLE, &globalCoords[0],
177 &coordsOnProc[0], &vectorDisplacements[0],
178 MPI::DOUBLE);
179
180 MPI::COMM_WORLD.Allgatherv(&vels[0], count, MPI::DOUBLE, &globalVels[0],
181 &coordsOnProc[0], &vectorDisplacements[0],
182 MPI::DOUBLE);
183
184 MPI::COMM_WORLD.Allgatherv(&masses[0], localHullSites, MPI::DOUBLE,
185 &globalMasses[0], &hullSitesOnProc[0],
186 &displacements[0], MPI::DOUBLE);
187
188 // Free previous hull
189 qh_freeqhull(!qh_ALL);
190 qh_memfreeshort(&curlong, &totlong);
191 if (curlong || totlong)
192 std::cerr << "qhull internal warning (main): did not free %d bytes of long memory (%d pieces) "
193 << totlong << curlong << std::endl;
194
195 if (qh_new_qhull(dim_, globalHullSites, &globalCoords[0], ismalloc,
196 const_cast<char *>(options_.c_str()), NULL, stderr)){
197
198 sprintf(painCave.errMsg, "AlphaHull: Qhull failed to compute global convex hull");
199 painCave.isFatal = 1;
200 simError();
201
202 } //qh_new_qhull
203
204
205 #endif
206
207 //Set facet->center as the Voronoi center
208 qh_setvoronoi_all();
209
210
211 int convexNumVert = qh_setsize(qh_facetvertices (qh facet_list, NULL, false));
212 //Insert all the sample points, because, even with alpha=0, the alpha shape/alpha complex will
213 //contain them.
214
215 // tri::Allocator<CMeshO>::AddVertices(pm.cm,convexNumVert);
216
217 /*ivp length is 'qh num_vertices' because each vertex is accessed through its ID whose range is
218 0<=qh_pointid(vertex->point)<qh num_vertices*/
219 // vector<tri::Allocator<CMeshO>::VertexPointer> ivp(qh num_vertices);
220 /*i=0;
221 FORALLvertices{
222 if ((*vertex).point){
223 // pm.cm.vert[i].P()[0] = (*vertex).point[0];
224 // pm.cm.vert[i].P()[1] = (*vertex).point[1];
225 //pm.cm.vert[i].P()[2] = (*vertex).point[2];
226 // ivp[qh_pointid(vertex->point)] = &pm.cm.vert[i];
227 i++;
228 }
229 }
230 */
231 //Set of alpha complex triangles for alphashape filtering
232 setT* set= qh_settemp(4* qh num_facets);
233
234 qh visit_id++;
235 int numFacets=0;
236 std::vector<std::vector <int> > facetlist;
237 interiorPoint = qh interior_point;
238 FORALLfacet_(qh facet_list) {
239 numFacets++;
240 if (!facet->upperdelaunay) {
241 //For all facets (that are tetrahedrons)calculate the radius of the empty circumsphere considering
242 //the distance between the circumcenter and a vertex of the facet
243 vertexT* vertex = (vertexT *)(facet->vertices->e[0].p);
244 double* center = facet->center;
245 double radius = qh_pointdist(vertex->point,center,dim_);
246
247 if (radius>alpha_) // if the facet is not good consider the ridges
248 {
249 //if calculating the alphashape, unmark the facet ('good' is used as 'marked').
250 facet->good=false;
251
252 //Compute each ridge (triangle) once and test the cironference radius with alpha
253 facet->visitid= qh visit_id;
254 qh_makeridges(facet);
255 ridgeT *ridge, **ridgep;
256 int goodTriangles=0;
257 FOREACHridge_(facet->ridges) {
258 neighbor= otherfacet_(ridge, facet);
259 if (( neighbor->visitid != qh visit_id)){
260 //Calculate the radius of the circumference
261 pointT* p0 = ((vertexT*) (ridge->vertices->e[0].p))->point;
262 pointT* p1 = ((vertexT*) (ridge->vertices->e[1].p))->point;
263 pointT* p2 = ((vertexT*) (ridge->vertices->e[2].p))->point;
264
265 radius = calculate_circumradius(p0,p1,p2, dim_);
266
267 if(radius <=alpha_){
268 goodTriangles++;
269 //save the triangle (ridge) for subsequent filtering
270 qh_setappend(&set, ridge);
271 }
272 }
273 }
274
275 //If calculating the alphashape, mark the facet('good' is used as 'marked').
276 //This facet will have some triangles hidden by the facet's neighbor.
277 if(goodTriangles==4)
278 facet->good=true;
279
280 }
281 else //the facet is good. Put all the triangles of the tetrahedron in the mesh
282 {
283 //Compute each ridge (triangle) once
284 facet->visitid= qh visit_id;
285 //If calculating the alphashape, mark the facet('good' is used as 'marked').
286 //This facet will have some triangles hidden by the facet's neighbor.
287 facet->good=true;
288 qh_makeridges(facet);
289 ridgeT *ridge, **ridgep;
290 FOREACHridge_(facet->ridges) {
291 neighbor= otherfacet_(ridge, facet);
292 if ((neighbor->visitid != qh visit_id)){
293 qh_setappend(&set, ridge);
294 }
295 }
296 }
297 }
298 }
299 //assert(numFacets== qh num_facets);
300
301 //Filter the triangles (only the ones on the boundary of the alpha complex) and build the mesh
302
303
304
305 ridgeT *ridge, **ridgep;
306 FOREACHridge_(set) {
307 if ((!ridge->top->good || !ridge->bottom->good || ridge->top->upperdelaunay || ridge->bottom->upperdelaunay)){
308 // tri::Allocator<CMeshO>::FaceIterator fi=tri::Allocator<CMeshO>::AddFaces(pm.cm,1);
309 ridgesCount++;
310 int vertex_n, vertex_i;
311 Triangle face;
312
313 // Vector3d V3dNormal(facet->normal[0], facet->normal[1], facet->normal[2]);
314 //face.setNormal(V3dNormal);
315
316
317 //coordT *center = qh_getcenter(ridge->vertices);
318 //std::cout << "Centers are " << center[0] << " " <<center[1] << " " << center[2] << std::endl;
319 //Vector3d V3dCentroid(center[0], center[1], center[2]);
320 //face.setCentroid(V3dCentroid);
321
322
323 Vector3d faceVel = V3Zero;
324 Vector3d p[3];
325 RealType faceMass = 0.0;
326
327 int ver = 0;
328 std::vector<int> virtexlist;
329 FOREACHvertex_i_(ridge->vertices){
330 int id = qh_pointid(vertex->point);
331 p[ver][0] = vertex->point[0];
332 p[ver][1] = vertex->point[1];
333 p[ver][2] = vertex->point[2];
334 Vector3d vel;
335 RealType mass;
336 ver++;
337 virtexlist.push_back(id);
338 // std::cout << "Ridge: " << ridgesCount << " Vertex " << id << std::endl;
339
340 vel = bodydoubles[id]->getVel();
341 mass = bodydoubles[id]->getMass();
342 face.addVertexSD(bodydoubles[id]);
343
344
345 faceVel = faceVel + vel;
346 faceMass = faceMass + mass;
347 } //FOREACH Vertex
348 facetlist.push_back(virtexlist);
349 face.addVertices(p[0],p[1],p[2]);
350 face.setFacetMass(faceMass);
351 face.setFacetVelocity(faceVel / RealType(3.0));
352
353 RealType area = face.getArea();
354 area_ += area;
355 Vector3d normal = face.getUnitNormal();
356 RealType offset = ((0.0-p[0][0])*normal[0] + (0.0-p[0][1])*normal[1] + (0.0-p[0][2])*normal[2]);
357 RealType dist = normal[0] * interiorPoint[0] + normal[1]*interiorPoint[1] + normal[2]*interiorPoint[2];
358 std::cout << "Dist and normal and area are: " << normal << std::endl;
359 volume_ += dist *area/qh hull_dim;
360
361 Triangles_.push_back(face);
362 }
363 }
364
365 std::cout << "Volume is: " << volume_ << std::endl;
366
367 //assert(pm.cm.fn == ridgesCount);
368 /*
369 std::cout <<"OFF"<<std::endl;
370 std::cout << bodydoubles.size() << " " << facetlist.size() << " " << 3*facetlist.size() << std::endl;
371 for (SD =bodydoubles.begin(); SD != bodydoubles.end(); ++SD){
372 Vector3d pos = (*SD)->getPos();
373 std::cout << pos.x() << " " << pos.y() << " " << pos.z() << std::endl;
374 }
375
376
377 std::vector<std::vector<int> >::iterator thisfacet;
378 std::vector<int>::iterator thisvertex;
379
380 for (thisfacet = facetlist.begin(); thisfacet != facetlist.end(); thisfacet++){
381 std::cout << (*thisfacet).size();
382 for (thisvertex = (*thisfacet).begin(); thisvertex != (*thisfacet).end(); thisvertex++){
383 std::cout << " " << *thisvertex;
384 }
385 std::cout << std::endl;
386 }
387 */
388
389
390
391 /*
392 FORALLfacets {
393 Triangle face;
394
395 Vector3d V3dNormal(facet->normal[0], facet->normal[1], facet->normal[2]);
396 face.setNormal(V3dNormal);
397
398 RealType faceArea = qh_facetarea(facet);
399 face.setArea(faceArea);
400
401 vertices = qh_facet3vertex(facet);
402
403 coordT *center = qh_getcenter(vertices);
404 Vector3d V3dCentroid(center[0], center[1], center[2]);
405 face.setCentroid(V3dCentroid);
406
407 Vector3d faceVel = V3Zero;
408 Vector3d p[3];
409 RealType faceMass = 0.0;
410
411 int ver = 0;
412
413 FOREACHvertex_(vertices){
414 int id = qh_pointid(vertex->point);
415 p[ver][0] = vertex->point[0];
416 p[ver][1] = vertex->point[1];
417 p[ver][2] = vertex->point[2];
418
419 Vector3d vel;
420 RealType mass;
421
422 #ifdef IS_MPI
423 vel = Vector3d(globalVels[dim_ * id],
424 globalVels[dim_ * id + 1],
425 globalVels[dim_ * id + 2]);
426 mass = globalMasses[id];
427
428 // localID will be between 0 and hullSitesOnProc[myrank] if we
429 // own this guy.
430
431 int localID = id - displacements[myrank];
432
433 if (localID >= 0 && localID < hullSitesOnProc[myrank])
434 face.addVertexSD(bodydoubles[indexMap[localID]]);
435
436 #else
437 vel = bodydoubles[id]->getVel();
438 mass = bodydoubles[id]->getMass();
439 face.addVertexSD(bodydoubles[id]);
440 #endif
441
442 faceVel = faceVel + vel;
443 faceMass = faceMass + mass;
444 ver++;
445 } //Foreachvertex
446
447 face.addVertices(p[0], p[1], p[2]);
448 face.setFacetMass(faceMass);
449 face.setFacetVelocity(faceVel/3.0);
450 Triangles_.push_back(face);
451 qh_settempfree(&vertices);
452
453 } //FORALLfacets
454 */
455 // qh_getarea(qh facet_list);
456 //volume_ = qh totvol;
457 // area_ = qh totarea;
458
459 qh_freeqhull(!qh_ALL);
460 qh_memfreeshort(&curlong, &totlong);
461 if (curlong || totlong)
462 std::cerr << "qhull internal warning (main): did not free %d bytes of long memory (%d pieces) "
463 << totlong << curlong << std::endl;
464 }
465
466 void AlphaHull::printHull(const std::string& geomFileName) {
467
468 #ifdef IS_MPI
469 if (worldRank == 0) {
470 #endif
471 FILE *newGeomFile;
472
473 //create new .md file based on old .md file
474 newGeomFile = fopen(geomFileName.c_str(), "w");
475 qh_findgood_all(qh facet_list);
476 for (int i = 0; i < qh_PRINTEND; i++)
477 qh_printfacets(newGeomFile, qh PRINTout[i], qh facet_list, NULL, !qh_ALL);
478
479 fclose(newGeomFile);
480 #ifdef IS_MPI
481 }
482 #endif
483 }
484
485 double calculate_circumradius(pointT* p0,pointT* p1,pointT* p2, int dim){
486 coordT a = qh_pointdist(p0,p1,dim);
487 coordT b = qh_pointdist(p1,p2,dim);
488 coordT c = qh_pointdist(p2,p0,dim);
489
490 coordT sum =(a + b + c)*0.5;
491 coordT area = sum*(a+b-sum)*(a+c-sum)*(b+c-sum);
492 return (double) (a*b*c)/(4*sqrt(area));
493 }
494
495 #endif //QHULL

Properties

Name Value
svn:keywords Author Id Revision Date