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 1404 by chuckv, Fri Jan 15 20:46:49 2010 UTC vs.
branches/development/src/math/ConvexHull.cpp (file contents), Revision 1665 by gezelter, Tue Nov 22 20:38:56 2011 UTC

# Line 1 | Line 1
1 < /* Copyright (c) 2008, 2009, 2010 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 67 | Line 67 | extern "C"
67   #ifdef HAVE_QHULL
68   extern "C"
69   {
70 < #include <qhull/qhull.h>
70 > #include <qhull/libqhull.h>
71   #include <qhull/mem.h>
72   #include <qhull/qset.h>
73   #include <qhull/geom.h>
# Line 105 | Line 105 | void ConvexHull::computeHull(std::vector<StuntDouble*>
105      i++;
106    }
107    
108 +  /* Clean up memory from previous convex hull calculations */
109    boolT ismalloc = False;
109  /* Clean up memory from previous convex hull calculations*/
110    
111 +  /* compute the hull for our local points (or all the points for single
112 +     processor versions) */
113    if (qh_new_qhull(dim_, numpoints, &ptArray[0], ismalloc,
114                     const_cast<char *>(options_.c_str()), NULL, stderr)) {
115 <
115 >    
116      sprintf(painCave.errMsg, "ConvexHull: Qhull failed to compute convex hull");
117      painCave.isFatal = 1;
118      simError();
# Line 194 | Line 196 | void ConvexHull::computeHull(std::vector<StuntDouble*>
196    // Free previous hull
197    qh_freeqhull(!qh_ALL);
198    qh_memfreeshort(&curlong, &totlong);
199 <  if (curlong || totlong)
200 <    std::cerr << "qhull internal warning (main): did not free %d bytes of long memory (%d pieces) "
201 <              << totlong << curlong << std::endl;
199 >  if (curlong || totlong) {
200 >    sprintf(painCave.errMsg, "ConvexHull: qhull internal warning:\n"
201 >            "\tdid not free %d bytes of long memory (%d pieces)",
202 >            totlong, curlong);
203 >    painCave.isFatal = 1;
204 >    simError();
205 >  }
206    
207    if (qh_new_qhull(dim_, globalHullSites, &globalCoords[0], ismalloc,
208                     const_cast<char *>(options_.c_str()), NULL, stderr)){
# Line 236 | Line 242 | void ConvexHull::computeHull(std::vector<StuntDouble*>
242        p[ver][0] = vertex->point[0];
243        p[ver][1] = vertex->point[1];
244        p[ver][2] = vertex->point[2];
239      
245        Vector3d vel;
246        RealType mass;
247  
# Line 251 | Line 256 | void ConvexHull::computeHull(std::vector<StuntDouble*>
256  
257        int localID = id - displacements[myrank];
258  
259 <      if (localID >= 0 && localID < hullSitesOnProc[myrank])
259 >
260 >      if (localID >= 0 && localID < hullSitesOnProc[myrank]){
261          face.addVertexSD(bodydoubles[indexMap[localID]]);
262 <      
262 >      }else{
263 >        face.addVertexSD(NULL);
264 >      }
265   #else
266        vel = bodydoubles[id]->getVel();
267        mass = bodydoubles[id]->getMass();
268        face.addVertexSD(bodydoubles[id]);      
269 < #endif
262 <        
269 > #endif  
270        faceVel = faceVel + vel;
271        faceMass = faceMass + mass;
272        ver++;      
# Line 288 | Line 295 | void ConvexHull::computeHull(std::vector<StuntDouble*>
295    //  std::cout << "My volume is: " << calcvol << " qhull volume is:" << volume_ << std::endl;
296    qh_freeqhull(!qh_ALL);
297    qh_memfreeshort(&curlong, &totlong);
298 <  if (curlong || totlong)
299 <    std::cerr << "qhull internal warning (main): did not free %d bytes of long memory (%d pieces) "
300 <              << totlong << curlong << std::endl;    
298 >  if (curlong || totlong) {
299 >    sprintf(painCave.errMsg, "ConvexHull: qhull internal warning:\n"
300 >            "\tdid not free %d bytes of long memory (%d pieces)",
301 >            totlong, curlong);
302 >    painCave.isFatal = 1;
303 >    simError();
304 >  }
305   }
306  
307   void ConvexHull::printHull(const std::string& geomFileName) {

Comparing:
trunk/src/math/ConvexHull.cpp (property svn:keywords), Revision 1404 by chuckv, Fri Jan 15 20:46:49 2010 UTC vs.
branches/development/src/math/ConvexHull.cpp (property svn:keywords), Revision 1665 by gezelter, Tue Nov 22 20:38:56 2011 UTC

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

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines