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 1374 by chuckv, Tue Oct 20 20:05:28 2009 UTC vs.
branches/development/src/math/ConvexHull.cpp (file contents), Revision 1867 by gezelter, Mon Apr 29 17:53:48 2013 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
5   * redistribute this software in source and binary code form, provided
6   * that the following conditions are met:
7   *
8 < * 1. Acknowledgement of the program authors must be made in any
9 < *    publication of scientific results based in part on use of the
10 < *    program.  An acceptable form of acknowledgement is citation of
11 < *    the article in which the program was described (Matthew
12 < *    A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
13 < *    J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
14 < *    Parallel Simulation Engine for Molecular Dynamics,"
15 < *    J. Comput. Chem. 26, pp. 252-271 (2005))
16 < *
17 < * 2. Redistributions of source code must retain the above copyright
8 > * 1. Redistributions of source code must retain the above copyright
9   *    notice, this list of conditions and the following disclaimer.
10   *
11 < * 3. Redistributions in binary form must reproduce the above copyright
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.
# Line 37 | Line 28
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 + * [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 + * [4] Kuang & Gezelter,  J. Chem. Phys. 133, 164101 (2010).
39 + * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
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.16 2009-10-20 20:05:28 chuckv Exp $
48 < *
43 > *  Purpose: To calculate a convex hull.
44   */
45  
46   /* Standard includes independent of library */
# Line 62 | Line 57
57   #include <mpi.h>
58   #endif
59  
60 < using namespace oopse;
60 > #include "math/qhull.hpp"
61  
62   #ifdef HAVE_QHULL
63 < extern "C"
64 < {
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 < }
63 > using namespace OpenMD;
64 > using namespace std;
65  
66 < /* Old options Qt Qu Qg QG0 FA */
81 < /* More old opts Qc Qi Pp*/
82 <
83 < ConvexHull::ConvexHull() : Hull(), dim_(3), options_("qhull Qt Pp") {
66 > ConvexHull::ConvexHull() : Hull(), dim_(3), options_("qhull FA Qt Pp") {
67   }
68  
69 < void ConvexHull::computeHull(std::vector<StuntDouble*> bodydoubles) {
70 <
69 > void ConvexHull::computeHull(vector<StuntDouble*> bodydoubles) {
70 >  
71    int numpoints = bodydoubles.size();
72 <
72 >  
73    Triangles_.clear();
74    
75    vertexT *vertex, **vertexp;
76    facetT *facet;
77    setT *vertices;
78    int curlong, totlong;
96  
97  std::vector<double> ptArray(numpoints*3);
98  std::vector<bool> isSurfaceID(numpoints);
79  
80 +  vector<double> ptArray(numpoints*dim_);
81 +
82    // Copy the positon vector into a points vector for qhull.
83 <  std::vector<StuntDouble*>::iterator SD;
83 >  vector<StuntDouble*>::iterator SD;
84    int i = 0;
85 +
86    for (SD =bodydoubles.begin(); SD != bodydoubles.end(); ++SD){
87      Vector3d pos = (*SD)->getPos();      
88      ptArray[dim_ * i] = pos.x();
# Line 108 | Line 91 | void ConvexHull::computeHull(std::vector<StuntDouble*>
91      i++;
92    }
93    
94 +  /* Clean up memory from previous convex hull calculations */
95    boolT ismalloc = False;
112  /* Clean up memory from previous convex hull calculations*/
96    
97 +  /* compute the hull for our local points (or all the points for single
98 +     processor versions) */
99    if (qh_new_qhull(dim_, numpoints, &ptArray[0], ismalloc,
100                     const_cast<char *>(options_.c_str()), NULL, stderr)) {
101 <
101 >    
102      sprintf(painCave.errMsg, "ConvexHull: Qhull failed to compute convex hull");
103      painCave.isFatal = 1;
104      simError();
# Line 127 | Line 112 | void ConvexHull::computeHull(std::vector<StuntDouble*>
112    int nproc = MPI::COMM_WORLD.Get_size();
113    int myrank = MPI::COMM_WORLD.Get_rank();
114    int localHullSites = 0;
130  int* hullSitesOnProc = new int[nproc];
131  int* coordsOnProc = new int[nproc];
132  int* displacements = new int[nproc];
133  int* vectorDisplacements = new int[nproc];
115  
116 <  std::vector<double> coords;
117 <  std::vector<double> vels;
118 <  std::vector<int> objectIDs;
119 <  std::vector<double> masses;
116 >  vector<int> hullSitesOnProc(nproc, 0);
117 >  vector<int> coordsOnProc(nproc, 0);
118 >  vector<int> displacements(nproc, 0);
119 >  vector<int> vectorDisplacements(nproc, 0);
120  
121 +  vector<double> coords;
122 +  vector<double> vels;
123 +  vector<int> indexMap;
124 +  vector<double> masses;
125 +
126    FORALLvertices{
127      localHullSites++;
128      
129      int idx = qh_pointid(vertex->point);
130 +
131 +    indexMap.push_back(idx);
132 +
133      coords.push_back(ptArray[dim_  * idx]);
134      coords.push_back(ptArray[dim_  * idx + 1]);
135      coords.push_back(ptArray[dim_  * idx + 2]);
# Line 159 | Line 148 | void ConvexHull::computeHull(std::vector<StuntDouble*>
148                              1, MPI::INT);
149  
150    int globalHullSites = 0;
151 <  for (int iproc = 0; i < nproc; iproc++){
151 >  for (int iproc = 0; iproc < nproc; iproc++){
152      globalHullSites += hullSitesOnProc[iproc];
153      coordsOnProc[iproc] = dim_ * hullSitesOnProc[iproc];
154    }
# Line 167 | Line 156 | void ConvexHull::computeHull(std::vector<StuntDouble*>
156    displacements[0] = 0;
157    vectorDisplacements[0] = 0;
158    
159 <  for (int iproc = 1; i < nproc; iproc++){
159 >  for (int iproc = 1; iproc < nproc; iproc++){
160      displacements[iproc] = displacements[iproc-1] + hullSitesOnProc[iproc-1];
161      vectorDisplacements[iproc] = vectorDisplacements[iproc-1] + coordsOnProc[iproc-1];
162    }
163  
164 <  std::vector<double> globalCoords(dim_*globalHullSites);
165 <  std::vector<double> globalVels(dim_*globalHullSites);
166 <  std::vector<double> globalMasses(globalHullSites);
164 >  vector<double> globalCoords(dim_ * globalHullSites);
165 >  vector<double> globalVels(dim_ * globalHullSites);
166 >  vector<double> globalMasses(globalHullSites);
167 >
168    int count = coordsOnProc[myrank];
169    
170 <  MPI::COMM_WORLD.Allgatherv(&coords[0], count, MPI::DOUBLE,
171 <                             &globalCoords[0], &coordsOnProc[0], &vectorDisplacements[0],
170 >  MPI::COMM_WORLD.Allgatherv(&coords[0], count, MPI::DOUBLE, &globalCoords[0],
171 >                             &coordsOnProc[0], &vectorDisplacements[0],
172                               MPI::DOUBLE);
173  
174 <  MPI::COMM_WORLD.Allgatherv(&vels[0], count, MPI::DOUBLE,
175 <                             &globalVels[0], &coordsOnProc[0], &vectorDisplacements[0],
174 >  MPI::COMM_WORLD.Allgatherv(&vels[0], count, MPI::DOUBLE, &globalVels[0],
175 >                             &coordsOnProc[0], &vectorDisplacements[0],
176                               MPI::DOUBLE);
177  
178    MPI::COMM_WORLD.Allgatherv(&masses[0], localHullSites, MPI::DOUBLE,
179 <                             &globalMasses[0], &hullSitesOnProc[0], &displacements[0],
180 <                             MPI::DOUBLE);
179 >                             &globalMasses[0], &hullSitesOnProc[0],
180 >                             &displacements[0], MPI::DOUBLE);
181  
182    // Free previous hull
183    qh_freeqhull(!qh_ALL);
184    qh_memfreeshort(&curlong, &totlong);
185 <  if (curlong || totlong)
186 <    std::cerr << "qhull internal warning (main): did not free %d bytes of long memory (%d pieces) "
187 <              << totlong << curlong << std::endl;
185 >  if (curlong || totlong) {
186 >    sprintf(painCave.errMsg, "ConvexHull: qhull internal warning:\n"
187 >            "\tdid not free %d bytes of long memory (%d pieces)",
188 >            totlong, curlong);
189 >    painCave.isFatal = 1;
190 >    simError();
191 >  }
192    
193    if (qh_new_qhull(dim_, globalHullSites, &globalCoords[0], ismalloc,
194                     const_cast<char *>(options_.c_str()), NULL, stderr)){
195      
196 <    sprintf(painCave.errMsg, "ConvexHull: Qhull failed to compute global convex hull");
196 >    sprintf(painCave.errMsg,
197 >            "ConvexHull: Qhull failed to compute global convex hull");
198      painCave.isFatal = 1;
199      simError();
200      
201    } //qh_new_qhull
202  
203   #endif
204 +  // commented out below, so comment out here also.
205 +  // intPoint = qh interior_point;
206 +  // RealType calcvol = 0.0;
207 +  
208 +  qh_triangulate ();
209 +  int num_facets = qh num_facets;
210 +  int num_vertices = qh num_vertices;
211  
212    FORALLfacets {  
213      Triangle face;
214 <
214 >    //Qhull sets the unit normal in facet->normal
215      Vector3d V3dNormal(facet->normal[0], facet->normal[1], facet->normal[2]);
216 <    face.setNormal(V3dNormal);
216 >    face.setUnitNormal(V3dNormal);
217      
218      RealType faceArea = qh_facetarea(facet);
219      face.setArea(faceArea);
# Line 225 | Line 227 | void ConvexHull::computeHull(std::vector<StuntDouble*>
227      Vector3d faceVel = V3Zero;
228      Vector3d p[3];
229      RealType faceMass = 0.0;
230 +
231      int ver = 0;
232  
233      FOREACHvertex_(vertices){
# Line 232 | Line 235 | void ConvexHull::computeHull(std::vector<StuntDouble*>
235        p[ver][0] = vertex->point[0];
236        p[ver][1] = vertex->point[1];
237        p[ver][2] = vertex->point[2];
235      
238        Vector3d vel;
239        RealType mass;
240  
# Line 242 | Line 244 | void ConvexHull::computeHull(std::vector<StuntDouble*>
244                       globalVels[dim_ * id + 2]);
245        mass = globalMasses[id];
246  
247 <      // localID will be between 0 and hullSitesOnProc[myrank] if we own this guy.
247 >      // localID will be between 0 and hullSitesOnProc[myrank] if we
248 >      // own this guy.
249 >
250        int localID = id - displacements[myrank];
251 <      if (id >= 0 && id < hullSitesOnProc[myrank])
252 <        face.addVertexSD(bodydoubles[localID]);
253 <      else
251 >
252 >
253 >      if (localID >= 0 && localID < hullSitesOnProc[myrank]){
254 >        face.addVertexSD(bodydoubles[indexMap[localID]]);
255 >      }else{
256          face.addVertexSD(NULL);
257 +      }
258   #else
259        vel = bodydoubles[id]->getVel();
260        mass = bodydoubles[id]->getMass();
261        face.addVertexSD(bodydoubles[id]);      
262 < #endif
256 <        
262 > #endif  
263        faceVel = faceVel + vel;
264        faceMass = faceMass + mass;
265        ver++;      
# Line 261 | Line 267 | void ConvexHull::computeHull(std::vector<StuntDouble*>
267  
268      face.addVertices(p[0], p[1], p[2]);
269      face.setFacetMass(faceMass);
270 <    face.setFacetVelocity(faceVel/3.0);
270 >    face.setFacetVelocity(faceVel / RealType(3.0));
271 >    /*
272 >    RealType comparea = face.computeArea();
273 >    realT calcarea = qh_facetarea (facet);
274 >    Vector3d V3dCompNorm = -face.computeUnitNormal();
275 >    RealType thisOffset = ((0.0-p[0][0])*V3dCompNorm[0] + (0.0-p[0][1])*V3dCompNorm[1] + (0.0-p[0][2])*V3dCompNorm[2]);
276 >    RealType dist = facet->offset + intPoint[0]*V3dNormal[0] + intPoint[1]*V3dNormal[1] + intPoint[2]*V3dNormal[2];
277 >    cout << "facet offset and computed offset: " << facet->offset << "  " << thisOffset <<  endl;
278 >    calcvol +=  -dist*comparea/qh hull_dim;
279 >    */
280      Triangles_.push_back(face);
281      qh_settempfree(&vertices);      
282  
# Line 270 | Line 285 | void ConvexHull::computeHull(std::vector<StuntDouble*>
285    qh_getarea(qh facet_list);
286    volume_ = qh totvol;
287    area_ = qh totarea;
288 <
274 < #ifdef IS_MPI
275 <  delete [] hullSitesOnProc;
276 <  delete [] coordsOnProc;
277 <  delete [] displacements;
278 <  delete [] vectorDisplacements;
279 < #endif
280 <  
288 >    
289    qh_freeqhull(!qh_ALL);
290    qh_memfreeshort(&curlong, &totlong);
291 <  if (curlong || totlong)
292 <    std::cerr << "qhull internal warning (main): did not free %d bytes of long memory (%d pieces) "
293 <              << totlong << curlong << std::endl;    
291 >  if (curlong || totlong) {
292 >    sprintf(painCave.errMsg, "ConvexHull: qhull internal warning:\n"
293 >            "\tdid not free %d bytes of long memory (%d pieces)",
294 >            totlong, curlong);
295 >    painCave.isFatal = 1;
296 >    simError();
297 >  }
298   }
299  
300 <
289 <
290 < void ConvexHull::printHull(const std::string& geomFileName) {
291 <  FILE *newGeomFile;
300 > void ConvexHull::printHull(const string& geomFileName) {
301    
302 <  //create new .md file based on old .md file
303 <  newGeomFile = fopen(geomFileName.c_str(), "w");
304 <  qh_findgood_all(qh facet_list);
305 <  for (int i = 0; i < qh_PRINTEND; i++)
306 <    qh_printfacets(newGeomFile, qh PRINTout[i], qh facet_list, NULL, !qh_ALL);
307 <  
308 <  fclose(newGeomFile);
302 > #ifdef IS_MPI
303 >  if (worldRank == 0)  {
304 > #endif
305 >    FILE *newGeomFile;
306 >    
307 >    //create new .md file based on old .md file
308 >    newGeomFile = fopen(geomFileName.c_str(), "w");
309 >    qh_findgood_all(qh facet_list);
310 >    for (int i = 0; i < qh_PRINTEND; i++)
311 >      qh_printfacets(newGeomFile, qh PRINTout[i], qh facet_list, NULL, !qh_ALL);
312 >    
313 >    fclose(newGeomFile);
314 > #ifdef IS_MPI
315 >  }
316 > #endif
317   }
318   #endif //QHULL

Comparing:
trunk/src/math/ConvexHull.cpp (property svn:keywords), Revision 1374 by chuckv, Tue Oct 20 20:05:28 2009 UTC vs.
branches/development/src/math/ConvexHull.cpp (property svn:keywords), Revision 1867 by gezelter, Mon Apr 29 17:53:48 2013 UTC

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

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines