ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/math/AlphaHull.cpp
Revision: 1704
Committed: Tue Apr 24 20:40:04 2012 UTC (13 years ago) by gezelter
File size: 16446 byte(s)
Log Message:
qhull modifications for compilation 


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, 24107 (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, **vertexp;
82 facetT *facet, *neighbor;
83 setT *vertices, *verticestop, *verticesbottom;
84 int curlong, totlong;
85 pointT *interiorPoint;
86
87 std::vector<double> ptArray(numpoints*dim_);
88
89 // Copy the positon vector into a points vector for qhull.
90 std::vector<StuntDouble*>::iterator SD;
91 int i = 0;
92 for (SD =bodydoubles.begin(); SD != bodydoubles.end(); ++SD){
93 Vector3d pos = (*SD)->getPos();
94 ptArray[dim_ * i] = pos.x();
95 ptArray[dim_ * i + 1] = pos.y();
96 ptArray[dim_ * i + 2] = pos.z();
97 i++;
98 }
99
100 /* Clean up memory from previous convex hull calculations*/
101 boolT ismalloc = False;
102
103 int ridgesCount=0;
104 if (qh_new_qhull(dim_, numpoints, &ptArray[0], ismalloc,
105 const_cast<char *>(options_.c_str()), NULL, stderr)) {
106
107 sprintf(painCave.errMsg, "AlphaHull: Qhull failed to compute convex hull");
108 painCave.isFatal = 1;
109 simError();
110
111 } //qh_new_qhull
112
113
114 #ifdef IS_MPI
115 //If we are doing the mpi version, set up some vectors for data communication
116
117 int nproc = MPI::COMM_WORLD.Get_size();
118 int myrank = MPI::COMM_WORLD.Get_rank();
119 int localHullSites = 0;
120
121 std::vector<int> hullSitesOnProc(nproc, 0);
122 std::vector<int> coordsOnProc(nproc, 0);
123 std::vector<int> displacements(nproc, 0);
124 std::vector<int> vectorDisplacements(nproc, 0);
125
126 std::vector<double> coords;
127 std::vector<double> vels;
128 std::vector<int> indexMap;
129 std::vector<double> masses;
130
131 FORALLvertices{
132 localHullSites++;
133
134 int idx = qh_pointid(vertex->point);
135
136 indexMap.push_back(idx);
137
138 coords.push_back(ptArray[dim_ * idx]);
139 coords.push_back(ptArray[dim_ * idx + 1]);
140 coords.push_back(ptArray[dim_ * idx + 2]);
141
142 StuntDouble* sd = bodydoubles[idx];
143
144 Vector3d vel = sd->getVel();
145 vels.push_back(vel.x());
146 vels.push_back(vel.y());
147 vels.push_back(vel.z());
148
149 masses.push_back(sd->getMass());
150 }
151
152 MPI::COMM_WORLD.Allgather(&localHullSites, 1, MPI::INT, &hullSitesOnProc[0],
153 1, MPI::INT);
154
155 int globalHullSites = 0;
156 for (int iproc = 0; iproc < nproc; iproc++){
157 globalHullSites += hullSitesOnProc[iproc];
158 coordsOnProc[iproc] = dim_ * hullSitesOnProc[iproc];
159 }
160
161 displacements[0] = 0;
162 vectorDisplacements[0] = 0;
163
164 for (int iproc = 1; iproc < nproc; iproc++){
165 displacements[iproc] = displacements[iproc-1] + hullSitesOnProc[iproc-1];
166 vectorDisplacements[iproc] = vectorDisplacements[iproc-1] + coordsOnProc[iproc-1];
167 }
168
169 std::vector<double> globalCoords(dim_ * globalHullSites);
170 std::vector<double> globalVels(dim_ * globalHullSites);
171 std::vector<double> globalMasses(globalHullSites);
172
173 int count = coordsOnProc[myrank];
174
175 MPI::COMM_WORLD.Allgatherv(&coords[0], count, MPI::DOUBLE, &globalCoords[0],
176 &coordsOnProc[0], &vectorDisplacements[0],
177 MPI::DOUBLE);
178
179 MPI::COMM_WORLD.Allgatherv(&vels[0], count, MPI::DOUBLE, &globalVels[0],
180 &coordsOnProc[0], &vectorDisplacements[0],
181 MPI::DOUBLE);
182
183 MPI::COMM_WORLD.Allgatherv(&masses[0], localHullSites, MPI::DOUBLE,
184 &globalMasses[0], &hullSitesOnProc[0],
185 &displacements[0], MPI::DOUBLE);
186
187 // Free previous hull
188 qh_freeqhull(!qh_ALL);
189 qh_memfreeshort(&curlong, &totlong);
190 if (curlong || totlong)
191 std::cerr << "qhull internal warning (main): did not free %d bytes of long memory (%d pieces) "
192 << totlong << curlong << std::endl;
193
194 if (qh_new_qhull(dim_, globalHullSites, &globalCoords[0], ismalloc,
195 const_cast<char *>(options_.c_str()), NULL, stderr)){
196
197 sprintf(painCave.errMsg, "AlphaHull: Qhull failed to compute global convex hull");
198 painCave.isFatal = 1;
199 simError();
200
201 } //qh_new_qhull
202
203
204 #endif
205
206 //Set facet->center as the Voronoi center
207 qh_setvoronoi_all();
208
209
210 int convexNumVert = qh_setsize(qh_facetvertices (qh facet_list, NULL, false));
211 //Insert all the sample points, because, even with alpha=0, the alpha shape/alpha complex will
212 //contain them.
213
214 // tri::Allocator<CMeshO>::AddVertices(pm.cm,convexNumVert);
215
216 /*ivp length is 'qh num_vertices' because each vertex is accessed through its ID whose range is
217 0<=qh_pointid(vertex->point)<qh num_vertices*/
218 // vector<tri::Allocator<CMeshO>::VertexPointer> ivp(qh num_vertices);
219 /*i=0;
220 FORALLvertices{
221 if ((*vertex).point){
222 // pm.cm.vert[i].P()[0] = (*vertex).point[0];
223 // pm.cm.vert[i].P()[1] = (*vertex).point[1];
224 //pm.cm.vert[i].P()[2] = (*vertex).point[2];
225 // ivp[qh_pointid(vertex->point)] = &pm.cm.vert[i];
226 i++;
227 }
228 }
229 */
230 //Set of alpha complex triangles for alphashape filtering
231 setT* set= qh_settemp(4* qh num_facets);
232
233 qh visit_id++;
234 int numFacets=0;
235 std::vector<std::vector <int> > facetlist;
236 interiorPoint = qh interior_point;
237 FORALLfacet_(qh facet_list) {
238 numFacets++;
239 if (!facet->upperdelaunay) {
240 //For all facets (that are tetrahedrons)calculate the radius of the empty circumsphere considering
241 //the distance between the circumcenter and a vertex of the facet
242 vertexT* vertex = (vertexT *)(facet->vertices->e[0].p);
243 double* center = facet->center;
244 double radius = qh_pointdist(vertex->point,center,dim_);
245
246 if (radius>alpha_) // if the facet is not good consider the ridges
247 {
248 //if calculating the alphashape, unmark the facet ('good' is used as 'marked').
249 facet->good=false;
250
251 //Compute each ridge (triangle) once and test the cironference radius with alpha
252 facet->visitid= qh visit_id;
253 qh_makeridges(facet);
254 ridgeT *ridge, **ridgep;
255 int goodTriangles=0;
256 FOREACHridge_(facet->ridges) {
257 neighbor= otherfacet_(ridge, facet);
258 if (( neighbor->visitid != qh visit_id)){
259 //Calculate the radius of the circumference
260 pointT* p0 = ((vertexT*) (ridge->vertices->e[0].p))->point;
261 pointT* p1 = ((vertexT*) (ridge->vertices->e[1].p))->point;
262 pointT* p2 = ((vertexT*) (ridge->vertices->e[2].p))->point;
263
264 radius = calculate_circumradius(p0,p1,p2, dim_);
265
266 if(radius <=alpha_){
267 goodTriangles++;
268 //save the triangle (ridge) for subsequent filtering
269 qh_setappend(&set, ridge);
270 }
271 }
272 }
273
274 //If calculating the alphashape, mark the facet('good' is used as 'marked').
275 //This facet will have some triangles hidden by the facet's neighbor.
276 if(goodTriangles==4)
277 facet->good=true;
278
279 }
280 else //the facet is good. Put all the triangles of the tetrahedron in the mesh
281 {
282 //Compute each ridge (triangle) once
283 facet->visitid= qh visit_id;
284 //If calculating the alphashape, mark the facet('good' is used as 'marked').
285 //This facet will have some triangles hidden by the facet's neighbor.
286 facet->good=true;
287 qh_makeridges(facet);
288 ridgeT *ridge, **ridgep;
289 FOREACHridge_(facet->ridges) {
290 neighbor= otherfacet_(ridge, facet);
291 if ((neighbor->visitid != qh visit_id)){
292 qh_setappend(&set, ridge);
293 }
294 }
295 }
296 }
297 }
298 //assert(numFacets== qh num_facets);
299
300 //Filter the triangles (only the ones on the boundary of the alpha complex) and build the mesh
301
302
303
304 ridgeT *ridge, **ridgep;
305 FOREACHridge_(set) {
306 if ((!ridge->top->good || !ridge->bottom->good || ridge->top->upperdelaunay || ridge->bottom->upperdelaunay)){
307 // tri::Allocator<CMeshO>::FaceIterator fi=tri::Allocator<CMeshO>::AddFaces(pm.cm,1);
308 ridgesCount++;
309 int vertex_n, vertex_i;
310 Triangle face;
311
312 // Vector3d V3dNormal(facet->normal[0], facet->normal[1], facet->normal[2]);
313 //face.setNormal(V3dNormal);
314
315
316 //coordT *center = qh_getcenter(ridge->vertices);
317 //std::cout << "Centers are " << center[0] << " " <<center[1] << " " << center[2] << std::endl;
318 //Vector3d V3dCentroid(center[0], center[1], center[2]);
319 //face.setCentroid(V3dCentroid);
320
321
322 Vector3d faceVel = V3Zero;
323 Vector3d p[3];
324 RealType faceMass = 0.0;
325
326 int ver = 0;
327 std::vector<int> virtexlist;
328 FOREACHvertex_i_(ridge->vertices){
329 int id = qh_pointid(vertex->point);
330 p[ver][0] = vertex->point[0];
331 p[ver][1] = vertex->point[1];
332 p[ver][2] = vertex->point[2];
333 Vector3d vel;
334 RealType mass;
335 ver++;
336 virtexlist.push_back(id);
337 // std::cout << "Ridge: " << ridgesCount << " Vertex " << id << std::endl;
338
339 vel = bodydoubles[id]->getVel();
340 mass = bodydoubles[id]->getMass();
341 face.addVertexSD(bodydoubles[id]);
342
343
344 faceVel = faceVel + vel;
345 faceMass = faceMass + mass;
346 } //FOREACH Vertex
347 facetlist.push_back(virtexlist);
348 face.addVertices(p[0],p[1],p[2]);
349 face.setFacetMass(faceMass);
350 face.setFacetVelocity(faceVel / RealType(3.0));
351
352 RealType area = face.getArea();
353 area_ += area;
354 Vector3d normal = face.getUnitNormal();
355 RealType offset = ((0.0-p[0][0])*normal[0] + (0.0-p[0][1])*normal[1] + (0.0-p[0][2])*normal[2]);
356 RealType dist = normal[0] * interiorPoint[0] + normal[1]*interiorPoint[1] + normal[2]*interiorPoint[2];
357 std::cout << "Dist and normal and area are: " << normal << std::endl;
358 volume_ += dist *area/qh hull_dim;
359
360 Triangles_.push_back(face);
361 }
362 }
363
364 std::cout << "Volume is: " << volume_ << std::endl;
365
366 //assert(pm.cm.fn == ridgesCount);
367 /*
368 std::cout <<"OFF"<<std::endl;
369 std::cout << bodydoubles.size() << " " << facetlist.size() << " " << 3*facetlist.size() << std::endl;
370 for (SD =bodydoubles.begin(); SD != bodydoubles.end(); ++SD){
371 Vector3d pos = (*SD)->getPos();
372 std::cout << pos.x() << " " << pos.y() << " " << pos.z() << std::endl;
373 }
374
375
376 std::vector<std::vector<int> >::iterator thisfacet;
377 std::vector<int>::iterator thisvertex;
378
379 for (thisfacet = facetlist.begin(); thisfacet != facetlist.end(); thisfacet++){
380 std::cout << (*thisfacet).size();
381 for (thisvertex = (*thisfacet).begin(); thisvertex != (*thisfacet).end(); thisvertex++){
382 std::cout << " " << *thisvertex;
383 }
384 std::cout << std::endl;
385 }
386 */
387
388
389
390 /*
391 FORALLfacets {
392 Triangle face;
393
394 Vector3d V3dNormal(facet->normal[0], facet->normal[1], facet->normal[2]);
395 face.setNormal(V3dNormal);
396
397 RealType faceArea = qh_facetarea(facet);
398 face.setArea(faceArea);
399
400 vertices = qh_facet3vertex(facet);
401
402 coordT *center = qh_getcenter(vertices);
403 Vector3d V3dCentroid(center[0], center[1], center[2]);
404 face.setCentroid(V3dCentroid);
405
406 Vector3d faceVel = V3Zero;
407 Vector3d p[3];
408 RealType faceMass = 0.0;
409
410 int ver = 0;
411
412 FOREACHvertex_(vertices){
413 int id = qh_pointid(vertex->point);
414 p[ver][0] = vertex->point[0];
415 p[ver][1] = vertex->point[1];
416 p[ver][2] = vertex->point[2];
417
418 Vector3d vel;
419 RealType mass;
420
421 #ifdef IS_MPI
422 vel = Vector3d(globalVels[dim_ * id],
423 globalVels[dim_ * id + 1],
424 globalVels[dim_ * id + 2]);
425 mass = globalMasses[id];
426
427 // localID will be between 0 and hullSitesOnProc[myrank] if we
428 // own this guy.
429
430 int localID = id - displacements[myrank];
431
432 if (localID >= 0 && localID < hullSitesOnProc[myrank])
433 face.addVertexSD(bodydoubles[indexMap[localID]]);
434
435 #else
436 vel = bodydoubles[id]->getVel();
437 mass = bodydoubles[id]->getMass();
438 face.addVertexSD(bodydoubles[id]);
439 #endif
440
441 faceVel = faceVel + vel;
442 faceMass = faceMass + mass;
443 ver++;
444 } //Foreachvertex
445
446 face.addVertices(p[0], p[1], p[2]);
447 face.setFacetMass(faceMass);
448 face.setFacetVelocity(faceVel/3.0);
449 Triangles_.push_back(face);
450 qh_settempfree(&vertices);
451
452 } //FORALLfacets
453 */
454 // qh_getarea(qh facet_list);
455 //volume_ = qh totvol;
456 // area_ = qh totarea;
457
458 qh_freeqhull(!qh_ALL);
459 qh_memfreeshort(&curlong, &totlong);
460 if (curlong || totlong)
461 std::cerr << "qhull internal warning (main): did not free %d bytes of long memory (%d pieces) "
462 << totlong << curlong << std::endl;
463 }
464
465 void AlphaHull::printHull(const std::string& geomFileName) {
466
467 #ifdef IS_MPI
468 if (worldRank == 0) {
469 #endif
470 FILE *newGeomFile;
471
472 //create new .md file based on old .md file
473 newGeomFile = fopen(geomFileName.c_str(), "w");
474 qh_findgood_all(qh facet_list);
475 for (int i = 0; i < qh_PRINTEND; i++)
476 qh_printfacets(newGeomFile, qh PRINTout[i], qh facet_list, NULL, !qh_ALL);
477
478 fclose(newGeomFile);
479 #ifdef IS_MPI
480 }
481 #endif
482 }
483
484 double calculate_circumradius(pointT* p0,pointT* p1,pointT* p2, int dim){
485 coordT a = qh_pointdist(p0,p1,dim);
486 coordT b = qh_pointdist(p1,p2,dim);
487 coordT c = qh_pointdist(p2,p0,dim);
488
489 coordT sum =(a + b + c)*0.5;
490 coordT area = sum*(a+b-sum)*(a+c-sum)*(b+c-sum);
491 return (double) (a*b*c)/(4*sqrt(area));
492 }
493
494 #endif //QHULL

Properties

Name Value
svn:keywords Author Id Revision Date