ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/math/ConvexHull.cpp
Revision: 1825
Committed: Wed Jan 9 19:27:52 2013 UTC (12 years, 3 months ago) by gezelter
File size: 10224 byte(s)
Log Message:
Deleting unused variables

File Contents

# User Rev Content
1 gezelter 1411 /* Copyright (c) 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 gezelter 1665 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
39     * [4] , Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011). *
40 chuckv 1097 *
41     * 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 1442 * @version $Id$
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 1704 #include "math/qhull.hpp"
66 chuckv 1293
67 chuckv 1374 #ifdef HAVE_QHULL
68 gezelter 1704 using namespace OpenMD;
69 chuckv 1293
70 chuckv 1374 ConvexHull::ConvexHull() : Hull(), dim_(3), options_("qhull Qt Pp") {
71     }
72 chuckv 1293
73 chuckv 1374 void ConvexHull::computeHull(std::vector<StuntDouble*> bodydoubles) {
74 chuckv 1293
75 chuckv 1374 int numpoints = bodydoubles.size();
76 chuckv 1293
77 chuckv 1374 Triangles_.clear();
78 chuckv 1242
79 chuckv 1293 vertexT *vertex, **vertexp;
80     facetT *facet;
81     setT *vertices;
82 chuckv 1374 int curlong, totlong;
83 gezelter 1825 // pointT *intPoint;
84 chuckv 1242
85 gezelter 1378 std::vector<double> ptArray(numpoints*dim_);
86 chuckv 1293
87     // Copy the positon vector into a points vector for qhull.
88     std::vector<StuntDouble*>::iterator SD;
89     int i = 0;
90 gezelter 1812
91 chuckv 1374 for (SD =bodydoubles.begin(); SD != bodydoubles.end(); ++SD){
92     Vector3d pos = (*SD)->getPos();
93     ptArray[dim_ * i] = pos.x();
94     ptArray[dim_ * i + 1] = pos.y();
95     ptArray[dim_ * i + 2] = pos.z();
96     i++;
97     }
98 chuckv 1242
99 gezelter 1411 /* Clean up memory from previous convex hull calculations */
100 chuckv 1293 boolT ismalloc = False;
101 chuckv 1308
102 gezelter 1411 /* compute the hull for our local points (or all the points for single
103     processor versions) */
104 chuckv 1302 if (qh_new_qhull(dim_, numpoints, &ptArray[0], ismalloc,
105 chuckv 1374 const_cast<char *>(options_.c_str()), NULL, stderr)) {
106 gezelter 1411
107 chuckv 1374 sprintf(painCave.errMsg, "ConvexHull: Qhull failed to compute convex hull");
108     painCave.isFatal = 1;
109     simError();
110    
111 chuckv 1293 } //qh_new_qhull
112    
113    
114     #ifdef IS_MPI
115 chuckv 1374 //If we are doing the mpi version, set up some vectors for data communication
116 chuckv 1242
117 chuckv 1374 int nproc = MPI::COMM_WORLD.Get_size();
118     int myrank = MPI::COMM_WORLD.Get_rank();
119     int localHullSites = 0;
120 chuckv 1293
121 gezelter 1384 std::vector<int> hullSitesOnProc(nproc, 0);
122     std::vector<int> coordsOnProc(nproc, 0);
123     std::vector<int> displacements(nproc, 0);
124     std::vector<int> vectorDisplacements(nproc, 0);
125    
126 chuckv 1374 std::vector<double> coords;
127     std::vector<double> vels;
128 gezelter 1378 std::vector<int> indexMap;
129 chuckv 1374 std::vector<double> masses;
130    
131     FORALLvertices{
132     localHullSites++;
133 chuckv 1242
134 chuckv 1374 int idx = qh_pointid(vertex->point);
135 gezelter 1378
136     indexMap.push_back(idx);
137    
138 chuckv 1374 coords.push_back(ptArray[dim_ * idx]);
139     coords.push_back(ptArray[dim_ * idx + 1]);
140     coords.push_back(ptArray[dim_ * idx + 2]);
141 chuckv 1293
142 chuckv 1374 StuntDouble* sd = bodydoubles[idx];
143 chuckv 1188
144 chuckv 1374 Vector3d vel = sd->getVel();
145     vels.push_back(vel.x());
146     vels.push_back(vel.y());
147     vels.push_back(vel.z());
148 chuckv 1188
149 chuckv 1374 masses.push_back(sd->getMass());
150     }
151 chuckv 1188
152 chuckv 1374 MPI::COMM_WORLD.Allgather(&localHullSites, 1, MPI::INT, &hullSitesOnProc[0],
153     1, MPI::INT);
154 chuckv 1307
155 chuckv 1374 int globalHullSites = 0;
156 gezelter 1376 for (int iproc = 0; iproc < nproc; iproc++){
157 chuckv 1374 globalHullSites += hullSitesOnProc[iproc];
158     coordsOnProc[iproc] = dim_ * hullSitesOnProc[iproc];
159 chuckv 1302 }
160 chuckv 1293
161 chuckv 1374 displacements[0] = 0;
162     vectorDisplacements[0] = 0;
163 chuckv 1316
164 gezelter 1376 for (int iproc = 1; iproc < nproc; iproc++){
165 chuckv 1374 displacements[iproc] = displacements[iproc-1] + hullSitesOnProc[iproc-1];
166     vectorDisplacements[iproc] = vectorDisplacements[iproc-1] + coordsOnProc[iproc-1];
167 chuckv 1293 }
168    
169 gezelter 1377 std::vector<double> globalCoords(dim_ * globalHullSites);
170     std::vector<double> globalVels(dim_ * globalHullSites);
171 chuckv 1374 std::vector<double> globalMasses(globalHullSites);
172 gezelter 1377
173 chuckv 1374 int count = coordsOnProc[myrank];
174 chuckv 1365
175 gezelter 1377 MPI::COMM_WORLD.Allgatherv(&coords[0], count, MPI::DOUBLE, &globalCoords[0],
176     &coordsOnProc[0], &vectorDisplacements[0],
177 chuckv 1374 MPI::DOUBLE);
178 chuckv 1302
179 gezelter 1377 MPI::COMM_WORLD.Allgatherv(&vels[0], count, MPI::DOUBLE, &globalVels[0],
180     &coordsOnProc[0], &vectorDisplacements[0],
181 chuckv 1374 MPI::DOUBLE);
182 chuckv 1302
183 chuckv 1374 MPI::COMM_WORLD.Allgatherv(&masses[0], localHullSites, MPI::DOUBLE,
184 gezelter 1377 &globalMasses[0], &hullSitesOnProc[0],
185     &displacements[0], MPI::DOUBLE);
186 chuckv 1365
187 chuckv 1293 // Free previous hull
188 chuckv 1242 qh_freeqhull(!qh_ALL);
189     qh_memfreeshort(&curlong, &totlong);
190 gezelter 1411 if (curlong || totlong) {
191     sprintf(painCave.errMsg, "ConvexHull: qhull internal warning:\n"
192     "\tdid not free %d bytes of long memory (%d pieces)",
193     totlong, curlong);
194     painCave.isFatal = 1;
195     simError();
196     }
197 chuckv 1374
198     if (qh_new_qhull(dim_, globalHullSites, &globalCoords[0], ismalloc,
199     const_cast<char *>(options_.c_str()), NULL, stderr)){
200    
201     sprintf(painCave.errMsg, "ConvexHull: Qhull failed to compute global convex hull");
202     painCave.isFatal = 1;
203     simError();
204    
205 chuckv 1293 } //qh_new_qhull
206    
207     #endif
208 gezelter 1825 // commented out below, so comment out here also.
209     // intPoint = qh interior_point;
210     // RealType calcvol = 0.0;
211 chuckv 1374 FORALLfacets {
212     Triangle face;
213 chuckv 1404 //Qhull sets the unit normal in facet->normal
214 chuckv 1374 Vector3d V3dNormal(facet->normal[0], facet->normal[1], facet->normal[2]);
215 chuckv 1404 face.setUnitNormal(V3dNormal);
216 chuckv 1371
217 chuckv 1374 RealType faceArea = qh_facetarea(facet);
218     face.setArea(faceArea);
219    
220     vertices = qh_facet3vertex(facet);
221 chuckv 1293
222 chuckv 1374 coordT *center = qh_getcenter(vertices);
223     Vector3d V3dCentroid(center[0], center[1], center[2]);
224     face.setCentroid(V3dCentroid);
225 chuckv 1304
226 chuckv 1374 Vector3d faceVel = V3Zero;
227     Vector3d p[3];
228     RealType faceMass = 0.0;
229 gezelter 1377
230 chuckv 1374 int ver = 0;
231 chuckv 1293
232 chuckv 1374 FOREACHvertex_(vertices){
233     int id = qh_pointid(vertex->point);
234     p[ver][0] = vertex->point[0];
235     p[ver][1] = vertex->point[1];
236     p[ver][2] = vertex->point[2];
237     Vector3d vel;
238     RealType mass;
239 chuckv 1293
240 chuckv 1302 #ifdef IS_MPI
241 chuckv 1374 vel = Vector3d(globalVels[dim_ * id],
242     globalVels[dim_ * id + 1],
243     globalVels[dim_ * id + 2]);
244     mass = globalMasses[id];
245    
246 gezelter 1377 // localID will be between 0 and hullSitesOnProc[myrank] if we
247     // own this guy.
248    
249 chuckv 1374 int localID = id - displacements[myrank];
250 gezelter 1377
251 chuckv 1410
252     if (localID >= 0 && localID < hullSitesOnProc[myrank]){
253 gezelter 1378 face.addVertexSD(bodydoubles[indexMap[localID]]);
254 chuckv 1410 }else{
255     face.addVertexSD(NULL);
256     }
257 chuckv 1307 #else
258 chuckv 1374 vel = bodydoubles[id]->getVel();
259     mass = bodydoubles[id]->getMass();
260     face.addVertexSD(bodydoubles[id]);
261 gezelter 1411 #endif
262 chuckv 1374 faceVel = faceVel + vel;
263     faceMass = faceMass + mass;
264     ver++;
265     } //Foreachvertex
266 chuckv 1302
267 chuckv 1374 face.addVertices(p[0], p[1], p[2]);
268     face.setFacetMass(faceMass);
269 gezelter 1668 face.setFacetVelocity(faceVel / RealType(3.0));
270 chuckv 1404 /*
271     RealType comparea = face.computeArea();
272     realT calcarea = qh_facetarea (facet);
273     Vector3d V3dCompNorm = -face.computeUnitNormal();
274     RealType thisOffset = ((0.0-p[0][0])*V3dCompNorm[0] + (0.0-p[0][1])*V3dCompNorm[1] + (0.0-p[0][2])*V3dCompNorm[2]);
275     RealType dist = facet->offset + intPoint[0]*V3dNormal[0] + intPoint[1]*V3dNormal[1] + intPoint[2]*V3dNormal[2];
276     std::cout << "facet offset and computed offset: " << facet->offset << " " << thisOffset << std::endl;
277     calcvol += -dist*comparea/qh hull_dim;
278     */
279 chuckv 1374 Triangles_.push_back(face);
280     qh_settempfree(&vertices);
281    
282     } //FORALLfacets
283    
284     qh_getarea(qh facet_list);
285     volume_ = qh totvol;
286     area_ = qh totarea;
287 chuckv 1404 // std::cout << "My volume is: " << calcvol << " qhull volume is:" << volume_ << std::endl;
288 chuckv 1374 qh_freeqhull(!qh_ALL);
289     qh_memfreeshort(&curlong, &totlong);
290 gezelter 1411 if (curlong || totlong) {
291     sprintf(painCave.errMsg, "ConvexHull: qhull internal warning:\n"
292     "\tdid not free %d bytes of long memory (%d pieces)",
293     totlong, curlong);
294     painCave.isFatal = 1;
295     simError();
296     }
297 chuckv 1097 }
298    
299 chuckv 1374 void ConvexHull::printHull(const std::string& geomFileName) {
300 gezelter 1384
301     #ifdef IS_MPI
302     if (worldRank == 0) {
303     #endif
304 chuckv 1242 FILE *newGeomFile;
305    
306     //create new .md file based on old .md file
307     newGeomFile = fopen(geomFileName.c_str(), "w");
308     qh_findgood_all(qh facet_list);
309     for (int i = 0; i < qh_PRINTEND; i++)
310     qh_printfacets(newGeomFile, qh PRINTout[i], qh facet_list, NULL, !qh_ALL);
311    
312     fclose(newGeomFile);
313 gezelter 1384 #ifdef IS_MPI
314     }
315     #endif
316 chuckv 1097 }
317 chuckv 1242 #endif //QHULL

Properties

Name Value
svn:keywords Author Id Revision Date