ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/math/ConvexHull.cpp
Revision: 1704
Committed: Tue Apr 24 20:40:04 2012 UTC (13 years ago) by gezelter
File size: 10162 byte(s)
Log Message:
qhull modifications for compilation 


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

Properties

Name Value
svn:keywords Author Id Revision Date