ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/math/ConvexHull.cpp
Revision: 1879
Committed: Sun Jun 16 15:15:42 2013 UTC (11 years, 10 months ago) by gezelter
File size: 9982 byte(s)
Log Message:
MERGE OpenMD development 1783:1878 into trunk

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 gezelter 1879 * [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 1782 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
39 gezelter 1879 * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
40 chuckv 1097 *
41     * ConvexHull.cpp
42     *
43 gezelter 1879 * 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 1782 #include "math/qhull.hpp"
61 chuckv 1293
62 chuckv 1374 #ifdef HAVE_QHULL
63 gezelter 1782 using namespace OpenMD;
64 gezelter 1879 using namespace std;
65 chuckv 1293
66 gezelter 1879 ConvexHull::ConvexHull() : Hull(), dim_(3), options_("qhull FA Qt Pp") {
67 chuckv 1374 }
68 chuckv 1293
69 gezelter 1879 void ConvexHull::computeHull(vector<StuntDouble*> bodydoubles) {
70    
71 chuckv 1374 int numpoints = bodydoubles.size();
72 gezelter 1879
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 chuckv 1293
80 gezelter 1879 vector<double> ptArray(numpoints*dim_);
81    
82 chuckv 1293 // Copy the positon vector into a points vector for qhull.
83 gezelter 1879 vector<StuntDouble*>::iterator SD;
84 chuckv 1293 int i = 0;
85 gezelter 1879
86 chuckv 1374 for (SD =bodydoubles.begin(); SD != bodydoubles.end(); ++SD){
87     Vector3d pos = (*SD)->getPos();
88     ptArray[dim_ * i] = pos.x();
89     ptArray[dim_ * i + 1] = pos.y();
90     ptArray[dim_ * i + 2] = pos.z();
91     i++;
92     }
93 chuckv 1242
94 gezelter 1411 /* Clean up memory from previous convex hull calculations */
95 chuckv 1293 boolT ismalloc = False;
96 chuckv 1308
97 gezelter 1411 /* compute the hull for our local points (or all the points for single
98     processor versions) */
99 chuckv 1302 if (qh_new_qhull(dim_, numpoints, &ptArray[0], ismalloc,
100 chuckv 1374 const_cast<char *>(options_.c_str()), NULL, stderr)) {
101 gezelter 1411
102 chuckv 1374 sprintf(painCave.errMsg, "ConvexHull: Qhull failed to compute convex hull");
103     painCave.isFatal = 1;
104     simError();
105    
106 chuckv 1293 } //qh_new_qhull
107    
108    
109     #ifdef IS_MPI
110 chuckv 1374 //If we are doing the mpi version, set up some vectors for data communication
111 chuckv 1242
112 chuckv 1374 int nproc = MPI::COMM_WORLD.Get_size();
113     int myrank = MPI::COMM_WORLD.Get_rank();
114     int localHullSites = 0;
115 chuckv 1293
116 gezelter 1879 vector<int> hullSitesOnProc(nproc, 0);
117     vector<int> coordsOnProc(nproc, 0);
118     vector<int> displacements(nproc, 0);
119     vector<int> vectorDisplacements(nproc, 0);
120 gezelter 1384
121 gezelter 1879 vector<double> coords;
122     vector<double> vels;
123     vector<int> indexMap;
124     vector<double> masses;
125 chuckv 1374
126     FORALLvertices{
127     localHullSites++;
128 chuckv 1242
129 chuckv 1374 int idx = qh_pointid(vertex->point);
130 gezelter 1378
131     indexMap.push_back(idx);
132    
133 chuckv 1374 coords.push_back(ptArray[dim_ * idx]);
134     coords.push_back(ptArray[dim_ * idx + 1]);
135     coords.push_back(ptArray[dim_ * idx + 2]);
136 chuckv 1293
137 chuckv 1374 StuntDouble* sd = bodydoubles[idx];
138 chuckv 1188
139 chuckv 1374 Vector3d vel = sd->getVel();
140     vels.push_back(vel.x());
141     vels.push_back(vel.y());
142     vels.push_back(vel.z());
143 chuckv 1188
144 chuckv 1374 masses.push_back(sd->getMass());
145     }
146 chuckv 1188
147 chuckv 1374 MPI::COMM_WORLD.Allgather(&localHullSites, 1, MPI::INT, &hullSitesOnProc[0],
148     1, MPI::INT);
149 chuckv 1307
150 chuckv 1374 int globalHullSites = 0;
151 gezelter 1376 for (int iproc = 0; iproc < nproc; iproc++){
152 chuckv 1374 globalHullSites += hullSitesOnProc[iproc];
153     coordsOnProc[iproc] = dim_ * hullSitesOnProc[iproc];
154 chuckv 1302 }
155 chuckv 1293
156 chuckv 1374 displacements[0] = 0;
157     vectorDisplacements[0] = 0;
158 chuckv 1316
159 gezelter 1376 for (int iproc = 1; iproc < nproc; iproc++){
160 chuckv 1374 displacements[iproc] = displacements[iproc-1] + hullSitesOnProc[iproc-1];
161     vectorDisplacements[iproc] = vectorDisplacements[iproc-1] + coordsOnProc[iproc-1];
162 chuckv 1293 }
163    
164 gezelter 1879 vector<double> globalCoords(dim_ * globalHullSites);
165     vector<double> globalVels(dim_ * globalHullSites);
166     vector<double> globalMasses(globalHullSites);
167 gezelter 1377
168 chuckv 1374 int count = coordsOnProc[myrank];
169 chuckv 1365
170 gezelter 1377 MPI::COMM_WORLD.Allgatherv(&coords[0], count, MPI::DOUBLE, &globalCoords[0],
171     &coordsOnProc[0], &vectorDisplacements[0],
172 chuckv 1374 MPI::DOUBLE);
173 chuckv 1302
174 gezelter 1377 MPI::COMM_WORLD.Allgatherv(&vels[0], count, MPI::DOUBLE, &globalVels[0],
175     &coordsOnProc[0], &vectorDisplacements[0],
176 chuckv 1374 MPI::DOUBLE);
177 chuckv 1302
178 chuckv 1374 MPI::COMM_WORLD.Allgatherv(&masses[0], localHullSites, MPI::DOUBLE,
179 gezelter 1377 &globalMasses[0], &hullSitesOnProc[0],
180     &displacements[0], MPI::DOUBLE);
181 chuckv 1365
182 chuckv 1293 // Free previous hull
183 chuckv 1242 qh_freeqhull(!qh_ALL);
184     qh_memfreeshort(&curlong, &totlong);
185 gezelter 1411 if (curlong || totlong) {
186     sprintf(painCave.errMsg, "ConvexHull: qhull internal warning:\n"
187     "\tdid not free %d bytes of long memory (%d pieces)",
188     totlong, curlong);
189     painCave.isFatal = 1;
190     simError();
191     }
192 chuckv 1374
193     if (qh_new_qhull(dim_, globalHullSites, &globalCoords[0], ismalloc,
194     const_cast<char *>(options_.c_str()), NULL, stderr)){
195    
196 gezelter 1879 sprintf(painCave.errMsg,
197     "ConvexHull: Qhull failed to compute global convex hull");
198 chuckv 1374 painCave.isFatal = 1;
199     simError();
200    
201 chuckv 1293 } //qh_new_qhull
202    
203     #endif
204 gezelter 1879 // commented out below, so comment out here also.
205     // intPoint = qh interior_point;
206     // RealType calcvol = 0.0;
207    
208     qh_triangulate ();
209 gezelter 1782
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 1782 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 1879 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 1879
287 chuckv 1374 qh_freeqhull(!qh_ALL);
288     qh_memfreeshort(&curlong, &totlong);
289 gezelter 1411 if (curlong || totlong) {
290     sprintf(painCave.errMsg, "ConvexHull: qhull internal warning:\n"
291     "\tdid not free %d bytes of long memory (%d pieces)",
292     totlong, curlong);
293     painCave.isFatal = 1;
294     simError();
295     }
296 chuckv 1097 }
297    
298 gezelter 1879 void ConvexHull::printHull(const string& geomFileName) {
299    
300 gezelter 1384 #ifdef IS_MPI
301     if (worldRank == 0) {
302     #endif
303 gezelter 1879 FILE *newGeomFile;
304    
305     //create new .md file based on old .md file
306     newGeomFile = fopen(geomFileName.c_str(), "w");
307     qh_findgood_all(qh facet_list);
308     for (int i = 0; i < qh_PRINTEND; i++)
309     qh_printfacets(newGeomFile, qh PRINTout[i], qh facet_list, NULL, !qh_ALL);
310    
311     fclose(newGeomFile);
312 gezelter 1384 #ifdef IS_MPI
313     }
314     #endif
315 chuckv 1097 }
316 chuckv 1242 #endif //QHULL

Properties

Name Value
svn:keywords Author Id Revision Date