ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/math/AlphaHull.cpp
Revision: 1879
Committed: Sun Jun 16 15:15:42 2013 UTC (11 years, 10 months ago) by gezelter
File size: 16422 byte(s)
Log Message:
MERGE OpenMD development 1783:1878 into trunk

File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date