ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/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.
Revision 1411 by gezelter, Wed Feb 24 16:09:21 2010 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 90 | Line 90 | void ConvexHull::computeHull(std::vector<StuntDouble*>
90    facetT *facet;
91    setT *vertices;
92    int curlong, totlong;
93 +  pointT *intPoint;
94    
95    std::vector<double> ptArray(numpoints*dim_);
96  
# Line 104 | 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;
108  /* 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 193 | 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 207 | Line 214 | void ConvexHull::computeHull(std::vector<StuntDouble*>
214    } //qh_new_qhull
215  
216   #endif
217 <
217 >  intPoint = qh interior_point;
218 >  RealType calcvol = 0.0;
219    FORALLfacets {  
220      Triangle face;
221 <
221 >    //Qhull sets the unit normal in facet->normal
222      Vector3d V3dNormal(facet->normal[0], facet->normal[1], facet->normal[2]);
223 <    face.setNormal(V3dNormal);
223 >    face.setUnitNormal(V3dNormal);
224      
225      RealType faceArea = qh_facetarea(facet);
226      face.setArea(faceArea);
# Line 234 | 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];
237      
245        Vector3d vel;
246        RealType mass;
247  
# Line 249 | 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
260 <        
269 > #endif  
270        faceVel = faceVel + vel;
271        faceMass = faceMass + mass;
272        ver++;      
# Line 266 | Line 275 | void ConvexHull::computeHull(std::vector<StuntDouble*>
275      face.addVertices(p[0], p[1], p[2]);
276      face.setFacetMass(faceMass);
277      face.setFacetVelocity(faceVel/3.0);
278 +    /*
279 +    RealType comparea = face.computeArea();
280 +    realT calcarea = qh_facetarea (facet);
281 +    Vector3d V3dCompNorm = -face.computeUnitNormal();
282 +    RealType thisOffset = ((0.0-p[0][0])*V3dCompNorm[0] + (0.0-p[0][1])*V3dCompNorm[1] + (0.0-p[0][2])*V3dCompNorm[2]);
283 +    RealType dist = facet->offset + intPoint[0]*V3dNormal[0] + intPoint[1]*V3dNormal[1] + intPoint[2]*V3dNormal[2];
284 +    std::cout << "facet offset and computed offset: " << facet->offset << "  " << thisOffset <<  std::endl;
285 +    calcvol +=  -dist*comparea/qh hull_dim;
286 +    */
287      Triangles_.push_back(face);
288      qh_settempfree(&vertices);      
289  
# Line 274 | Line 292 | void ConvexHull::computeHull(std::vector<StuntDouble*>
292    qh_getarea(qh facet_list);
293    volume_ = qh totvol;
294    area_ = qh totarea;
295 <
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) {

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines