ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/math/ConvexHull.cpp
Revision: 1665
Committed: Tue Nov 22 20:38:56 2011 UTC (13 years, 5 months ago) by gezelter
File size: 10333 byte(s)
Log Message:
updated copyright notices

File Contents

# Content
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. Redistributions of source code must retain the above copyright
9 * notice, this list of conditions and the following disclaimer.
10 *
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.
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 * 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] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
39 * [4] , 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$
48 *
49 */
50
51 /* Standard includes independent of library */
52
53 #include <iostream>
54 #include <fstream>
55 #include <list>
56 #include <algorithm>
57 #include <iterator>
58 #include "math/ConvexHull.hpp"
59 #include "utils/simError.h"
60
61 #ifdef IS_MPI
62 #include <mpi.h>
63 #endif
64
65 using namespace OpenMD;
66
67 #ifdef HAVE_QHULL
68 extern "C"
69 {
70 #include <qhull/libqhull.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 }
79
80 ConvexHull::ConvexHull() : Hull(), dim_(3), options_("qhull Qt Pp") {
81 }
82
83 void ConvexHull::computeHull(std::vector<StuntDouble*> bodydoubles) {
84
85 int numpoints = bodydoubles.size();
86
87 Triangles_.clear();
88
89 vertexT *vertex, **vertexp;
90 facetT *facet;
91 setT *vertices;
92 int curlong, totlong;
93 pointT *intPoint;
94
95 std::vector<double> ptArray(numpoints*dim_);
96
97 // Copy the positon vector into a points vector for qhull.
98 std::vector<StuntDouble*>::iterator SD;
99 int i = 0;
100 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
108 /* Clean up memory from previous convex hull calculations */
109 boolT ismalloc = False;
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
116 sprintf(painCave.errMsg, "ConvexHull: Qhull failed to compute convex hull");
117 painCave.isFatal = 1;
118 simError();
119
120 } //qh_new_qhull
121
122
123 #ifdef IS_MPI
124 //If we are doing the mpi version, set up some vectors for data communication
125
126 int nproc = MPI::COMM_WORLD.Get_size();
127 int myrank = MPI::COMM_WORLD.Get_rank();
128 int localHullSites = 0;
129
130 std::vector<int> hullSitesOnProc(nproc, 0);
131 std::vector<int> coordsOnProc(nproc, 0);
132 std::vector<int> displacements(nproc, 0);
133 std::vector<int> vectorDisplacements(nproc, 0);
134
135 std::vector<double> coords;
136 std::vector<double> vels;
137 std::vector<int> indexMap;
138 std::vector<double> masses;
139
140 FORALLvertices{
141 localHullSites++;
142
143 int idx = qh_pointid(vertex->point);
144
145 indexMap.push_back(idx);
146
147 coords.push_back(ptArray[dim_ * idx]);
148 coords.push_back(ptArray[dim_ * idx + 1]);
149 coords.push_back(ptArray[dim_ * idx + 2]);
150
151 StuntDouble* sd = bodydoubles[idx];
152
153 Vector3d vel = sd->getVel();
154 vels.push_back(vel.x());
155 vels.push_back(vel.y());
156 vels.push_back(vel.z());
157
158 masses.push_back(sd->getMass());
159 }
160
161 MPI::COMM_WORLD.Allgather(&localHullSites, 1, MPI::INT, &hullSitesOnProc[0],
162 1, MPI::INT);
163
164 int globalHullSites = 0;
165 for (int iproc = 0; iproc < nproc; iproc++){
166 globalHullSites += hullSitesOnProc[iproc];
167 coordsOnProc[iproc] = dim_ * hullSitesOnProc[iproc];
168 }
169
170 displacements[0] = 0;
171 vectorDisplacements[0] = 0;
172
173 for (int iproc = 1; iproc < nproc; iproc++){
174 displacements[iproc] = displacements[iproc-1] + hullSitesOnProc[iproc-1];
175 vectorDisplacements[iproc] = vectorDisplacements[iproc-1] + coordsOnProc[iproc-1];
176 }
177
178 std::vector<double> globalCoords(dim_ * globalHullSites);
179 std::vector<double> globalVels(dim_ * globalHullSites);
180 std::vector<double> globalMasses(globalHullSites);
181
182 int count = coordsOnProc[myrank];
183
184 MPI::COMM_WORLD.Allgatherv(&coords[0], count, MPI::DOUBLE, &globalCoords[0],
185 &coordsOnProc[0], &vectorDisplacements[0],
186 MPI::DOUBLE);
187
188 MPI::COMM_WORLD.Allgatherv(&vels[0], count, MPI::DOUBLE, &globalVels[0],
189 &coordsOnProc[0], &vectorDisplacements[0],
190 MPI::DOUBLE);
191
192 MPI::COMM_WORLD.Allgatherv(&masses[0], localHullSites, MPI::DOUBLE,
193 &globalMasses[0], &hullSitesOnProc[0],
194 &displacements[0], MPI::DOUBLE);
195
196 // Free previous hull
197 qh_freeqhull(!qh_ALL);
198 qh_memfreeshort(&curlong, &totlong);
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)){
209
210 sprintf(painCave.errMsg, "ConvexHull: Qhull failed to compute global convex hull");
211 painCave.isFatal = 1;
212 simError();
213
214 } //qh_new_qhull
215
216 #endif
217 intPoint = qh interior_point;
218 RealType calcvol = 0.0;
219 FORALLfacets {
220 Triangle face;
221 //Qhull sets the unit normal in facet->normal
222 Vector3d V3dNormal(facet->normal[0], facet->normal[1], facet->normal[2]);
223 face.setUnitNormal(V3dNormal);
224
225 RealType faceArea = qh_facetarea(facet);
226 face.setArea(faceArea);
227
228 vertices = qh_facet3vertex(facet);
229
230 coordT *center = qh_getcenter(vertices);
231 Vector3d V3dCentroid(center[0], center[1], center[2]);
232 face.setCentroid(V3dCentroid);
233
234 Vector3d faceVel = V3Zero;
235 Vector3d p[3];
236 RealType faceMass = 0.0;
237
238 int ver = 0;
239
240 FOREACHvertex_(vertices){
241 int id = qh_pointid(vertex->point);
242 p[ver][0] = vertex->point[0];
243 p[ver][1] = vertex->point[1];
244 p[ver][2] = vertex->point[2];
245 Vector3d vel;
246 RealType mass;
247
248 #ifdef IS_MPI
249 vel = Vector3d(globalVels[dim_ * id],
250 globalVels[dim_ * id + 1],
251 globalVels[dim_ * id + 2]);
252 mass = globalMasses[id];
253
254 // localID will be between 0 and hullSitesOnProc[myrank] if we
255 // own this guy.
256
257 int localID = id - displacements[myrank];
258
259
260 if (localID >= 0 && localID < hullSitesOnProc[myrank]){
261 face.addVertexSD(bodydoubles[indexMap[localID]]);
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
270 faceVel = faceVel + vel;
271 faceMass = faceMass + mass;
272 ver++;
273 } //Foreachvertex
274
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
290 } //FORALLfacets
291
292 qh_getarea(qh facet_list);
293 volume_ = qh totvol;
294 area_ = qh totarea;
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 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) {
308
309 #ifdef IS_MPI
310 if (worldRank == 0) {
311 #endif
312 FILE *newGeomFile;
313
314 //create new .md file based on old .md file
315 newGeomFile = fopen(geomFileName.c_str(), "w");
316 qh_findgood_all(qh facet_list);
317 for (int i = 0; i < qh_PRINTEND; i++)
318 qh_printfacets(newGeomFile, qh PRINTout[i], qh facet_list, NULL, !qh_ALL);
319
320 fclose(newGeomFile);
321 #ifdef IS_MPI
322 }
323 #endif
324 }
325 #endif //QHULL

Properties

Name Value
svn:keywords Author Id Revision Date