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