ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/math/ConvexHull.cpp
Revision: 1378
Committed: Wed Oct 21 15:48:12 2009 UTC (15 years, 6 months ago) by gezelter
File size: 9358 byte(s)
Log Message:
fixed mapping  bug

File Contents

# User Rev Content
1 chuckv 1365 /* Copyright (c) 2008, 2009 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     * 1. Acknowledgement of the program authors must be made in any
9     * publication of scientific results based in part on use of the
10     * program. An acceptable form of acknowledgement is citation of
11     * the article in which the program was described (Matthew
12     * A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
13     * J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
14     * Parallel Simulation Engine for Molecular Dynamics,"
15     * J. Comput. Chem. 26, pp. 252-271 (2005))
16     *
17     * 2. Redistributions of source code must retain the above copyright
18     * notice, this list of conditions and the following disclaimer.
19     *
20     * 3. Redistributions in binary form must reproduce the above copyright
21     * notice, this list of conditions and the following disclaimer in the
22     * documentation and/or other materials provided with the
23     * distribution.
24     *
25     * This software is provided "AS IS," without a warranty of any
26     * kind. All express or implied conditions, representations and
27     * warranties, including any implied warranty of merchantability,
28     * fitness for a particular purpose or non-infringement, are hereby
29     * excluded. The University of Notre Dame and its licensors shall not
30     * be liable for any damages suffered by licensee as a result of
31     * using, modifying or distributing the software or its
32     * derivatives. In no event will the University of Notre Dame or its
33     * licensors be liable for any lost revenue, profit or data, or for
34     * direct, indirect, special, consequential, incidental or punitive
35     * damages, however caused and regardless of the theory of liability,
36     * arising out of the use of or inability to use software, even if the
37     * University of Notre Dame has been advised of the possibility of
38     * such damages.
39     *
40     *
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 1378 * @version $Id: ConvexHull.cpp,v 1.19 2009-10-21 15:48:12 gezelter Exp $
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 chuckv 1374 using namespace oopse;
66 chuckv 1293
67 chuckv 1374 #ifdef HAVE_QHULL
68     extern "C"
69 chuckv 1242 {
70 chuckv 1374 #include <qhull/qhull.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 chuckv 1293 }
79    
80 chuckv 1242 /* Old options Qt Qu Qg QG0 FA */
81 chuckv 1304 /* More old opts Qc Qi Pp*/
82 chuckv 1188
83 chuckv 1374 ConvexHull::ConvexHull() : Hull(), dim_(3), options_("qhull Qt Pp") {
84     }
85 chuckv 1293
86 chuckv 1374 void ConvexHull::computeHull(std::vector<StuntDouble*> bodydoubles) {
87 chuckv 1293
88 chuckv 1374 int numpoints = bodydoubles.size();
89 chuckv 1293
90 chuckv 1374 Triangles_.clear();
91 chuckv 1242
92 chuckv 1293 vertexT *vertex, **vertexp;
93     facetT *facet;
94     setT *vertices;
95 chuckv 1374 int curlong, totlong;
96 chuckv 1242
97 gezelter 1378 std::vector<double> ptArray(numpoints*dim_);
98 chuckv 1293
99     // Copy the positon vector into a points vector for qhull.
100     std::vector<StuntDouble*>::iterator SD;
101     int i = 0;
102 chuckv 1374 for (SD =bodydoubles.begin(); SD != bodydoubles.end(); ++SD){
103     Vector3d pos = (*SD)->getPos();
104     ptArray[dim_ * i] = pos.x();
105     ptArray[dim_ * i + 1] = pos.y();
106     ptArray[dim_ * i + 2] = pos.z();
107     i++;
108     }
109 chuckv 1242
110 chuckv 1293 boolT ismalloc = False;
111 chuckv 1302 /* Clean up memory from previous convex hull calculations*/
112 chuckv 1308
113 chuckv 1302 if (qh_new_qhull(dim_, numpoints, &ptArray[0], ismalloc,
114 chuckv 1374 const_cast<char *>(options_.c_str()), NULL, stderr)) {
115 chuckv 1188
116 chuckv 1374 sprintf(painCave.errMsg, "ConvexHull: Qhull failed to compute convex hull");
117     painCave.isFatal = 1;
118     simError();
119    
120 chuckv 1293 } //qh_new_qhull
121    
122    
123     #ifdef IS_MPI
124 chuckv 1374 //If we are doing the mpi version, set up some vectors for data communication
125 chuckv 1242
126 chuckv 1374 int nproc = MPI::COMM_WORLD.Get_size();
127     int myrank = MPI::COMM_WORLD.Get_rank();
128     int localHullSites = 0;
129     int* hullSitesOnProc = new int[nproc];
130     int* coordsOnProc = new int[nproc];
131     int* displacements = new int[nproc];
132     int* vectorDisplacements = new int[nproc];
133 chuckv 1293
134 chuckv 1374 std::vector<double> coords;
135     std::vector<double> vels;
136 gezelter 1378 std::vector<int> indexMap;
137 chuckv 1374 std::vector<double> masses;
138    
139     FORALLvertices{
140     localHullSites++;
141 chuckv 1242
142 chuckv 1374 int idx = qh_pointid(vertex->point);
143 gezelter 1378
144     indexMap.push_back(idx);
145    
146 chuckv 1374 coords.push_back(ptArray[dim_ * idx]);
147     coords.push_back(ptArray[dim_ * idx + 1]);
148     coords.push_back(ptArray[dim_ * idx + 2]);
149 chuckv 1293
150 chuckv 1374 StuntDouble* sd = bodydoubles[idx];
151 chuckv 1188
152 chuckv 1374 Vector3d vel = sd->getVel();
153     vels.push_back(vel.x());
154     vels.push_back(vel.y());
155     vels.push_back(vel.z());
156 chuckv 1188
157 chuckv 1374 masses.push_back(sd->getMass());
158     }
159 chuckv 1188
160 chuckv 1374 MPI::COMM_WORLD.Allgather(&localHullSites, 1, MPI::INT, &hullSitesOnProc[0],
161     1, MPI::INT);
162 chuckv 1307
163 chuckv 1374 int globalHullSites = 0;
164 gezelter 1376 for (int iproc = 0; iproc < nproc; iproc++){
165 chuckv 1374 globalHullSites += hullSitesOnProc[iproc];
166     coordsOnProc[iproc] = dim_ * hullSitesOnProc[iproc];
167 chuckv 1302 }
168 chuckv 1293
169 chuckv 1374 displacements[0] = 0;
170     vectorDisplacements[0] = 0;
171 chuckv 1316
172 gezelter 1376 for (int iproc = 1; iproc < nproc; iproc++){
173 chuckv 1374 displacements[iproc] = displacements[iproc-1] + hullSitesOnProc[iproc-1];
174     vectorDisplacements[iproc] = vectorDisplacements[iproc-1] + coordsOnProc[iproc-1];
175 chuckv 1293 }
176    
177 gezelter 1377 std::vector<double> globalCoords(dim_ * globalHullSites);
178     std::vector<double> globalVels(dim_ * globalHullSites);
179 chuckv 1374 std::vector<double> globalMasses(globalHullSites);
180 gezelter 1377
181 chuckv 1374 int count = coordsOnProc[myrank];
182 chuckv 1365
183 gezelter 1377 MPI::COMM_WORLD.Allgatherv(&coords[0], count, MPI::DOUBLE, &globalCoords[0],
184     &coordsOnProc[0], &vectorDisplacements[0],
185 chuckv 1374 MPI::DOUBLE);
186 chuckv 1302
187 gezelter 1377 MPI::COMM_WORLD.Allgatherv(&vels[0], count, MPI::DOUBLE, &globalVels[0],
188     &coordsOnProc[0], &vectorDisplacements[0],
189 chuckv 1374 MPI::DOUBLE);
190 chuckv 1302
191 chuckv 1374 MPI::COMM_WORLD.Allgatherv(&masses[0], localHullSites, MPI::DOUBLE,
192 gezelter 1377 &globalMasses[0], &hullSitesOnProc[0],
193     &displacements[0], MPI::DOUBLE);
194 chuckv 1365
195 chuckv 1293 // Free previous hull
196 chuckv 1242 qh_freeqhull(!qh_ALL);
197     qh_memfreeshort(&curlong, &totlong);
198     if (curlong || totlong)
199     std::cerr << "qhull internal warning (main): did not free %d bytes of long memory (%d pieces) "
200     << totlong << curlong << std::endl;
201 chuckv 1374
202     if (qh_new_qhull(dim_, globalHullSites, &globalCoords[0], ismalloc,
203     const_cast<char *>(options_.c_str()), NULL, stderr)){
204    
205     sprintf(painCave.errMsg, "ConvexHull: Qhull failed to compute global convex hull");
206     painCave.isFatal = 1;
207     simError();
208    
209 chuckv 1293 } //qh_new_qhull
210    
211     #endif
212    
213 chuckv 1374 FORALLfacets {
214     Triangle face;
215 chuckv 1293
216 chuckv 1374 Vector3d V3dNormal(facet->normal[0], facet->normal[1], facet->normal[2]);
217     face.setNormal(V3dNormal);
218 chuckv 1371
219 chuckv 1374 RealType faceArea = qh_facetarea(facet);
220     face.setArea(faceArea);
221    
222     vertices = qh_facet3vertex(facet);
223 chuckv 1293
224 chuckv 1374 coordT *center = qh_getcenter(vertices);
225     Vector3d V3dCentroid(center[0], center[1], center[2]);
226     face.setCentroid(V3dCentroid);
227 chuckv 1304
228 chuckv 1374 Vector3d faceVel = V3Zero;
229     Vector3d p[3];
230     RealType faceMass = 0.0;
231 gezelter 1377
232 chuckv 1374 int ver = 0;
233 chuckv 1293
234 chuckv 1374 FOREACHvertex_(vertices){
235     int id = qh_pointid(vertex->point);
236     p[ver][0] = vertex->point[0];
237     p[ver][1] = vertex->point[1];
238     p[ver][2] = vertex->point[2];
239    
240     Vector3d vel;
241     RealType mass;
242 chuckv 1293
243 chuckv 1302 #ifdef IS_MPI
244 chuckv 1374 vel = Vector3d(globalVels[dim_ * id],
245     globalVels[dim_ * id + 1],
246     globalVels[dim_ * id + 2]);
247     mass = globalMasses[id];
248    
249 gezelter 1377 // localID will be between 0 and hullSitesOnProc[myrank] if we
250     // own this guy.
251    
252 chuckv 1374 int localID = id - displacements[myrank];
253 gezelter 1377
254     if (localID >= 0 && localID < hullSitesOnProc[myrank])
255 gezelter 1378 face.addVertexSD(bodydoubles[indexMap[localID]]);
256 gezelter 1377
257 chuckv 1307 #else
258 chuckv 1374 vel = bodydoubles[id]->getVel();
259     mass = bodydoubles[id]->getMass();
260     face.addVertexSD(bodydoubles[id]);
261 chuckv 1302 #endif
262 chuckv 1374
263     faceVel = faceVel + vel;
264     faceMass = faceMass + mass;
265     ver++;
266     } //Foreachvertex
267 chuckv 1302
268 chuckv 1374 face.addVertices(p[0], p[1], p[2]);
269     face.setFacetMass(faceMass);
270     face.setFacetVelocity(faceVel/3.0);
271     Triangles_.push_back(face);
272     qh_settempfree(&vertices);
273    
274     } //FORALLfacets
275    
276     qh_getarea(qh facet_list);
277     volume_ = qh totvol;
278     area_ = qh totarea;
279    
280 chuckv 1302 #ifdef IS_MPI
281 chuckv 1374 delete [] hullSitesOnProc;
282     delete [] coordsOnProc;
283     delete [] displacements;
284     delete [] vectorDisplacements;
285 chuckv 1302 #endif
286 chuckv 1374
287     qh_freeqhull(!qh_ALL);
288     qh_memfreeshort(&curlong, &totlong);
289     if (curlong || totlong)
290     std::cerr << "qhull internal warning (main): did not free %d bytes of long memory (%d pieces) "
291     << totlong << curlong << std::endl;
292 chuckv 1097 }
293    
294 chuckv 1374 void ConvexHull::printHull(const std::string& geomFileName) {
295 chuckv 1242 FILE *newGeomFile;
296    
297     //create new .md file based on old .md file
298     newGeomFile = fopen(geomFileName.c_str(), "w");
299     qh_findgood_all(qh facet_list);
300     for (int i = 0; i < qh_PRINTEND; i++)
301     qh_printfacets(newGeomFile, qh PRINTout[i], qh facet_list, NULL, !qh_ALL);
302    
303     fclose(newGeomFile);
304 chuckv 1097 }
305 chuckv 1242 #endif //QHULL