ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/math/ConvexHull.cpp
Revision: 1404
Committed: Fri Jan 15 20:46:49 2010 UTC (15 years, 3 months ago) by chuckv
Original Path: trunk/src/math/ConvexHull.cpp
File size: 10065 byte(s)
Log Message:
Cleanup of unit normal code in Triangle.cpp
SMIPD now call getUnitNormal instead of getNormal.


File Contents

# User Rev Content
1 chuckv 1404 /* Copyright (c) 2008, 2009, 2010 The University of Notre Dame. All Rights Reserved.
2 chuckv 1097 *
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 gezelter 1390 * 1. Redistributions of source code must retain the above copyright
9 chuckv 1097 * notice, this list of conditions and the following disclaimer.
10     *
11 gezelter 1390 * 2. Redistributions in binary form must reproduce the above copyright
12 chuckv 1097 * 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 gezelter 1390 * 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, 24107 (2008).
38     * [4] Vardeman & Gezelter, in progress (2009).
39 chuckv 1097 *
40 gezelter 1390 *
41 chuckv 1097 * ConvexHull.cpp
42     *
43 chuckv 1137 * Purpose: To calculate convexhull, hull volume libqhull.
44 chuckv 1097 *
45     * Created by Charles F. Vardeman II on 11 Dec 2006.
46     * @author Charles F. Vardeman II
47 gezelter 1390 * @version $Id: ConvexHull.cpp,v 1.21 2009-11-25 20:02:01 gezelter Exp $
48 chuckv 1097 *
49     */
50    
51 chuckv 1242 /* Standard includes independent of library */
52 chuckv 1374
53 chuckv 1097 #include <iostream>
54     #include <fstream>
55 chuckv 1188 #include <list>
56     #include <algorithm>
57     #include <iterator>
58 chuckv 1141 #include "math/ConvexHull.hpp"
59 chuckv 1188 #include "utils/simError.h"
60 chuckv 1293
61     #ifdef IS_MPI
62 chuckv 1374 #include <mpi.h>
63 chuckv 1293 #endif
64    
65 gezelter 1390 using namespace OpenMD;
66 chuckv 1293
67 chuckv 1374 #ifdef HAVE_QHULL
68     extern "C"
69 chuckv 1242 {
70 chuckv 1374 #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 chuckv 1293 }
79    
80 chuckv 1374 ConvexHull::ConvexHull() : Hull(), dim_(3), options_("qhull Qt Pp") {
81     }
82 chuckv 1293
83 chuckv 1374 void ConvexHull::computeHull(std::vector<StuntDouble*> bodydoubles) {
84 chuckv 1293
85 chuckv 1374 int numpoints = bodydoubles.size();
86 chuckv 1293
87 chuckv 1374 Triangles_.clear();
88 chuckv 1242
89 chuckv 1293 vertexT *vertex, **vertexp;
90     facetT *facet;
91     setT *vertices;
92 chuckv 1374 int curlong, totlong;
93 chuckv 1404 pointT *intPoint;
94 chuckv 1242
95 gezelter 1378 std::vector<double> ptArray(numpoints*dim_);
96 chuckv 1293
97     // Copy the positon vector into a points vector for qhull.
98     std::vector<StuntDouble*>::iterator SD;
99     int i = 0;
100 chuckv 1374 for (SD =bodydoubles.begin(); SD != bodydoubles.end(); ++SD){
101     Vector3d pos = (*SD)->getPos();
102     ptArray[dim_ * i] = pos.x();
103     ptArray[dim_ * i + 1] = pos.y();
104     ptArray[dim_ * i + 2] = pos.z();
105     i++;
106     }
107 chuckv 1242
108 chuckv 1293 boolT ismalloc = False;
109 chuckv 1302 /* Clean up memory from previous convex hull calculations*/
110 chuckv 1308
111 chuckv 1302 if (qh_new_qhull(dim_, numpoints, &ptArray[0], ismalloc,
112 chuckv 1374 const_cast<char *>(options_.c_str()), NULL, stderr)) {
113 chuckv 1188
114 chuckv 1374 sprintf(painCave.errMsg, "ConvexHull: Qhull failed to compute convex hull");
115     painCave.isFatal = 1;
116     simError();
117    
118 chuckv 1293 } //qh_new_qhull
119    
120    
121     #ifdef IS_MPI
122 chuckv 1374 //If we are doing the mpi version, set up some vectors for data communication
123 chuckv 1242
124 chuckv 1374 int nproc = MPI::COMM_WORLD.Get_size();
125     int myrank = MPI::COMM_WORLD.Get_rank();
126     int localHullSites = 0;
127 chuckv 1293
128 gezelter 1384 std::vector<int> hullSitesOnProc(nproc, 0);
129     std::vector<int> coordsOnProc(nproc, 0);
130     std::vector<int> displacements(nproc, 0);
131     std::vector<int> vectorDisplacements(nproc, 0);
132    
133 chuckv 1374 std::vector<double> coords;
134     std::vector<double> vels;
135 gezelter 1378 std::vector<int> indexMap;
136 chuckv 1374 std::vector<double> masses;
137    
138     FORALLvertices{
139     localHullSites++;
140 chuckv 1242
141 chuckv 1374 int idx = qh_pointid(vertex->point);
142 gezelter 1378
143     indexMap.push_back(idx);
144    
145 chuckv 1374 coords.push_back(ptArray[dim_ * idx]);
146     coords.push_back(ptArray[dim_ * idx + 1]);
147     coords.push_back(ptArray[dim_ * idx + 2]);
148 chuckv 1293
149 chuckv 1374 StuntDouble* sd = bodydoubles[idx];
150 chuckv 1188
151 chuckv 1374 Vector3d vel = sd->getVel();
152     vels.push_back(vel.x());
153     vels.push_back(vel.y());
154     vels.push_back(vel.z());
155 chuckv 1188
156 chuckv 1374 masses.push_back(sd->getMass());
157     }
158 chuckv 1188
159 chuckv 1374 MPI::COMM_WORLD.Allgather(&localHullSites, 1, MPI::INT, &hullSitesOnProc[0],
160     1, MPI::INT);
161 chuckv 1307
162 chuckv 1374 int globalHullSites = 0;
163 gezelter 1376 for (int iproc = 0; iproc < nproc; iproc++){
164 chuckv 1374 globalHullSites += hullSitesOnProc[iproc];
165     coordsOnProc[iproc] = dim_ * hullSitesOnProc[iproc];
166 chuckv 1302 }
167 chuckv 1293
168 chuckv 1374 displacements[0] = 0;
169     vectorDisplacements[0] = 0;
170 chuckv 1316
171 gezelter 1376 for (int iproc = 1; iproc < nproc; iproc++){
172 chuckv 1374 displacements[iproc] = displacements[iproc-1] + hullSitesOnProc[iproc-1];
173     vectorDisplacements[iproc] = vectorDisplacements[iproc-1] + coordsOnProc[iproc-1];
174 chuckv 1293 }
175    
176 gezelter 1377 std::vector<double> globalCoords(dim_ * globalHullSites);
177     std::vector<double> globalVels(dim_ * globalHullSites);
178 chuckv 1374 std::vector<double> globalMasses(globalHullSites);
179 gezelter 1377
180 chuckv 1374 int count = coordsOnProc[myrank];
181 chuckv 1365
182 gezelter 1377 MPI::COMM_WORLD.Allgatherv(&coords[0], count, MPI::DOUBLE, &globalCoords[0],
183     &coordsOnProc[0], &vectorDisplacements[0],
184 chuckv 1374 MPI::DOUBLE);
185 chuckv 1302
186 gezelter 1377 MPI::COMM_WORLD.Allgatherv(&vels[0], count, MPI::DOUBLE, &globalVels[0],
187     &coordsOnProc[0], &vectorDisplacements[0],
188 chuckv 1374 MPI::DOUBLE);
189 chuckv 1302
190 chuckv 1374 MPI::COMM_WORLD.Allgatherv(&masses[0], localHullSites, MPI::DOUBLE,
191 gezelter 1377 &globalMasses[0], &hullSitesOnProc[0],
192     &displacements[0], MPI::DOUBLE);
193 chuckv 1365
194 chuckv 1293 // Free previous hull
195 chuckv 1242 qh_freeqhull(!qh_ALL);
196     qh_memfreeshort(&curlong, &totlong);
197     if (curlong || totlong)
198     std::cerr << "qhull internal warning (main): did not free %d bytes of long memory (%d pieces) "
199     << totlong << curlong << std::endl;
200 chuckv 1374
201     if (qh_new_qhull(dim_, globalHullSites, &globalCoords[0], ismalloc,
202     const_cast<char *>(options_.c_str()), NULL, stderr)){
203    
204     sprintf(painCave.errMsg, "ConvexHull: Qhull failed to compute global convex hull");
205     painCave.isFatal = 1;
206     simError();
207    
208 chuckv 1293 } //qh_new_qhull
209    
210     #endif
211 chuckv 1404 intPoint = qh interior_point;
212     RealType calcvol = 0.0;
213 chuckv 1374 FORALLfacets {
214     Triangle face;
215 chuckv 1404 //Qhull sets the unit normal in facet->normal
216 chuckv 1374 Vector3d V3dNormal(facet->normal[0], facet->normal[1], facet->normal[2]);
217 chuckv 1404 face.setUnitNormal(V3dNormal);
218 chuckv 1371
219 chuckv 1374 RealType faceArea = qh_facetarea(facet);
220     face.setArea(faceArea);
221    
222     vertices = qh_facet3vertex(facet);
223 chuckv 1293
224 chuckv 1374 coordT *center = qh_getcenter(vertices);
225     Vector3d V3dCentroid(center[0], center[1], center[2]);
226     face.setCentroid(V3dCentroid);
227 chuckv 1304
228 chuckv 1374 Vector3d faceVel = V3Zero;
229     Vector3d p[3];
230     RealType faceMass = 0.0;
231 gezelter 1377
232 chuckv 1374 int ver = 0;
233 chuckv 1293
234 chuckv 1374 FOREACHvertex_(vertices){
235     int id = qh_pointid(vertex->point);
236     p[ver][0] = vertex->point[0];
237     p[ver][1] = vertex->point[1];
238     p[ver][2] = vertex->point[2];
239    
240     Vector3d vel;
241     RealType mass;
242 chuckv 1293
243 chuckv 1302 #ifdef IS_MPI
244 chuckv 1374 vel = Vector3d(globalVels[dim_ * id],
245     globalVels[dim_ * id + 1],
246     globalVels[dim_ * id + 2]);
247     mass = globalMasses[id];
248    
249 gezelter 1377 // localID will be between 0 and hullSitesOnProc[myrank] if we
250     // own this guy.
251    
252 chuckv 1374 int localID = id - displacements[myrank];
253 gezelter 1377
254     if (localID >= 0 && localID < hullSitesOnProc[myrank])
255 gezelter 1378 face.addVertexSD(bodydoubles[indexMap[localID]]);
256 gezelter 1377
257 chuckv 1307 #else
258 chuckv 1374 vel = bodydoubles[id]->getVel();
259     mass = bodydoubles[id]->getMass();
260     face.addVertexSD(bodydoubles[id]);
261 chuckv 1302 #endif
262 chuckv 1374
263     faceVel = faceVel + vel;
264     faceMass = faceMass + mass;
265     ver++;
266     } //Foreachvertex
267 chuckv 1302
268 chuckv 1374 face.addVertices(p[0], p[1], p[2]);
269     face.setFacetMass(faceMass);
270     face.setFacetVelocity(faceVel/3.0);
271 chuckv 1404 /*
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     std::cout << "facet offset and computed offset: " << facet->offset << " " << thisOffset << std::endl;
278     calcvol += -dist*comparea/qh hull_dim;
279     */
280 chuckv 1374 Triangles_.push_back(face);
281     qh_settempfree(&vertices);
282    
283     } //FORALLfacets
284    
285     qh_getarea(qh facet_list);
286     volume_ = qh totvol;
287     area_ = qh totarea;
288 chuckv 1404 // std::cout << "My volume is: " << calcvol << " qhull volume is:" << volume_ << std::endl;
289 chuckv 1374 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;
294 chuckv 1097 }
295    
296 chuckv 1374 void ConvexHull::printHull(const std::string& geomFileName) {
297 gezelter 1384
298     #ifdef IS_MPI
299     if (worldRank == 0) {
300     #endif
301 chuckv 1242 FILE *newGeomFile;
302    
303     //create new .md file based on old .md file
304     newGeomFile = fopen(geomFileName.c_str(), "w");
305     qh_findgood_all(qh facet_list);
306     for (int i = 0; i < qh_PRINTEND; i++)
307     qh_printfacets(newGeomFile, qh PRINTout[i], qh facet_list, NULL, !qh_ALL);
308    
309     fclose(newGeomFile);
310 gezelter 1384 #ifdef IS_MPI
311     }
312     #endif
313 chuckv 1097 }
314 chuckv 1242 #endif //QHULL