ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/math/ConvexHull.cpp
(Generate patch)

Comparing:
trunk/src/math/ConvexHull.cpp (file contents), Revision 1390 by gezelter, Wed Nov 25 20:02:06 2009 UTC vs.
branches/development/src/math/ConvexHull.cpp (file contents), Revision 1704 by gezelter, Tue Apr 24 20:40:04 2012 UTC

# Line 1 | Line 1
1 < /* Copyright (c) 2008, 2009 The University of Notre Dame. All Rights Reserved.
1 > /* Copyright (c) 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
# Line 35 | Line 35
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]  Vardeman & Gezelter, in progress (2009).                        
38 > * [4] Kuang & Gezelter,  J. Chem. Phys. 133, 164101 (2010).
39 > * [4] , Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011). *
40   *
40 *
41   *  ConvexHull.cpp
42   *
43   *  Purpose: To calculate convexhull, hull volume libqhull.
44   *
45   *  Created by Charles F. Vardeman II on 11 Dec 2006.
46   *  @author  Charles F. Vardeman II
47 < *  @version $Id: ConvexHull.cpp,v 1.21 2009-11-25 20:02:01 gezelter Exp $
47 > *  @version $Id$
48   *
49   */
50  
# Line 62 | Line 62
62   #include <mpi.h>
63   #endif
64  
65 < using namespace OpenMD;
65 > #include "math/qhull.hpp"
66  
67   #ifdef HAVE_QHULL
68 < extern "C"
69 < {
70 < #include <qhull/qhull.h>
71 < #include <qhull/mem.h>
72 < #include <qhull/qset.h>
73 < #include <qhull/geom.h>
74 < #include <qhull/merge.h>
75 < #include <qhull/poly.h>
76 < #include <qhull/io.h>
77 < #include <qhull/stat.h>
78 < }
68 > using namespace OpenMD;
69  
70   ConvexHull::ConvexHull() : Hull(), dim_(3), options_("qhull Qt Pp") {
71   }
# Line 90 | Line 80 | void ConvexHull::computeHull(std::vector<StuntDouble*>
80    facetT *facet;
81    setT *vertices;
82    int curlong, totlong;
83 +  pointT *intPoint;
84    
85    std::vector<double> ptArray(numpoints*dim_);
86  
# Line 104 | Line 95 | void ConvexHull::computeHull(std::vector<StuntDouble*>
95      i++;
96    }
97    
98 +  /* Clean up memory from previous convex hull calculations */
99    boolT ismalloc = False;
108  /* Clean up memory from previous convex hull calculations*/
100    
101 +  /* compute the hull for our local points (or all the points for single
102 +     processor versions) */
103    if (qh_new_qhull(dim_, numpoints, &ptArray[0], ismalloc,
104                     const_cast<char *>(options_.c_str()), NULL, stderr)) {
105 <
105 >    
106      sprintf(painCave.errMsg, "ConvexHull: Qhull failed to compute convex hull");
107      painCave.isFatal = 1;
108      simError();
# Line 193 | Line 186 | void ConvexHull::computeHull(std::vector<StuntDouble*>
186    // Free previous hull
187    qh_freeqhull(!qh_ALL);
188    qh_memfreeshort(&curlong, &totlong);
189 <  if (curlong || totlong)
190 <    std::cerr << "qhull internal warning (main): did not free %d bytes of long memory (%d pieces) "
191 <              << totlong << curlong << std::endl;
189 >  if (curlong || totlong) {
190 >    sprintf(painCave.errMsg, "ConvexHull: qhull internal warning:\n"
191 >            "\tdid not free %d bytes of long memory (%d pieces)",
192 >            totlong, curlong);
193 >    painCave.isFatal = 1;
194 >    simError();
195 >  }
196    
197    if (qh_new_qhull(dim_, globalHullSites, &globalCoords[0], ismalloc,
198                     const_cast<char *>(options_.c_str()), NULL, stderr)){
# Line 207 | Line 204 | void ConvexHull::computeHull(std::vector<StuntDouble*>
204    } //qh_new_qhull
205  
206   #endif
207 <
207 >  intPoint = qh interior_point;
208 >  RealType calcvol = 0.0;
209    FORALLfacets {  
210      Triangle face;
211 <
211 >    //Qhull sets the unit normal in facet->normal
212      Vector3d V3dNormal(facet->normal[0], facet->normal[1], facet->normal[2]);
213 <    face.setNormal(V3dNormal);
213 >    face.setUnitNormal(V3dNormal);
214      
215      RealType faceArea = qh_facetarea(facet);
216      face.setArea(faceArea);
# Line 234 | Line 232 | void ConvexHull::computeHull(std::vector<StuntDouble*>
232        p[ver][0] = vertex->point[0];
233        p[ver][1] = vertex->point[1];
234        p[ver][2] = vertex->point[2];
237      
235        Vector3d vel;
236        RealType mass;
237  
# Line 249 | Line 246 | void ConvexHull::computeHull(std::vector<StuntDouble*>
246  
247        int localID = id - displacements[myrank];
248  
249 <      if (localID >= 0 && localID < hullSitesOnProc[myrank])
249 >
250 >      if (localID >= 0 && localID < hullSitesOnProc[myrank]){
251          face.addVertexSD(bodydoubles[indexMap[localID]]);
252 <      
252 >      }else{
253 >        face.addVertexSD(NULL);
254 >      }
255   #else
256        vel = bodydoubles[id]->getVel();
257        mass = bodydoubles[id]->getMass();
258        face.addVertexSD(bodydoubles[id]);      
259 < #endif
260 <        
259 > #endif  
260        faceVel = faceVel + vel;
261        faceMass = faceMass + mass;
262        ver++;      
# Line 265 | Line 264 | void ConvexHull::computeHull(std::vector<StuntDouble*>
264  
265      face.addVertices(p[0], p[1], p[2]);
266      face.setFacetMass(faceMass);
267 <    face.setFacetVelocity(faceVel/3.0);
267 >    face.setFacetVelocity(faceVel / RealType(3.0));
268 >    /*
269 >    RealType comparea = face.computeArea();
270 >    realT calcarea = qh_facetarea (facet);
271 >    Vector3d V3dCompNorm = -face.computeUnitNormal();
272 >    RealType thisOffset = ((0.0-p[0][0])*V3dCompNorm[0] + (0.0-p[0][1])*V3dCompNorm[1] + (0.0-p[0][2])*V3dCompNorm[2]);
273 >    RealType dist = facet->offset + intPoint[0]*V3dNormal[0] + intPoint[1]*V3dNormal[1] + intPoint[2]*V3dNormal[2];
274 >    std::cout << "facet offset and computed offset: " << facet->offset << "  " << thisOffset <<  std::endl;
275 >    calcvol +=  -dist*comparea/qh hull_dim;
276 >    */
277      Triangles_.push_back(face);
278      qh_settempfree(&vertices);      
279  
# Line 274 | Line 282 | void ConvexHull::computeHull(std::vector<StuntDouble*>
282    qh_getarea(qh facet_list);
283    volume_ = qh totvol;
284    area_ = qh totarea;
285 <
285 >  //  std::cout << "My volume is: " << calcvol << " qhull volume is:" << volume_ << std::endl;
286    qh_freeqhull(!qh_ALL);
287    qh_memfreeshort(&curlong, &totlong);
288 <  if (curlong || totlong)
289 <    std::cerr << "qhull internal warning (main): did not free %d bytes of long memory (%d pieces) "
290 <              << totlong << curlong << std::endl;    
288 >  if (curlong || totlong) {
289 >    sprintf(painCave.errMsg, "ConvexHull: qhull internal warning:\n"
290 >            "\tdid not free %d bytes of long memory (%d pieces)",
291 >            totlong, curlong);
292 >    painCave.isFatal = 1;
293 >    simError();
294 >  }
295   }
296  
297   void ConvexHull::printHull(const std::string& geomFileName) {

Comparing:
trunk/src/math/ConvexHull.cpp (property svn:keywords), Revision 1390 by gezelter, Wed Nov 25 20:02:06 2009 UTC vs.
branches/development/src/math/ConvexHull.cpp (property svn:keywords), Revision 1704 by gezelter, Tue Apr 24 20:40:04 2012 UTC

# Line 0 | Line 1
1 + Author Id Revision Date

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines