ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/math/AlphaHull.cpp
Revision: 1858
Committed: Wed Apr 3 21:32:13 2013 UTC (12 years ago) by gezelter
File size: 16976 byte(s)
Log Message:
Some bugfixes for cell-linked-list-style neighbor lists when the simulation
doesn't use periodic boundary conditions.   Cleaning up some of the Hull
stuff while we're in there.

File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date