OpenMD 3.0
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
ConvexHull.cpp
1/*
2 * Copyright (c) 2004-present, The University of Notre Dame. All rights
3 * reserved.
4 *
5 * Redistribution and use in source and binary forms, with or without
6 * modification, are permitted provided that the following conditions are met:
7 *
8 * 1. Redistributions of source code must retain the above copyright notice,
9 * this list of conditions and the following disclaimer.
10 *
11 * 2. Redistributions in binary form must reproduce the above copyright notice,
12 * this list of conditions and the following disclaimer in the documentation
13 * and/or other materials provided with the distribution.
14 *
15 * 3. Neither the name of the copyright holder nor the names of its
16 * contributors may be used to endorse or promote products derived from
17 * this software without specific prior written permission.
18 *
19 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
20 * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
21 * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
22 * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
23 * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
24 * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
25 * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
26 * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
27 * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
28 * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
29 * POSSIBILITY OF SUCH DAMAGE.
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, 234107 (2008).
38 * [4] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
39 * [5] Kuang & Gezelter, Mol. Phys., 110, 691-701 (2012).
40 * [6] Lamichhane, Gezelter & Newman, J. Chem. Phys. 141, 134109 (2014).
41 * [7] Lamichhane, Newman & Gezelter, J. Chem. Phys. 141, 134110 (2014).
42 * [8] Bhattarai, Newman & Gezelter, Phys. Rev. B 99, 094106 (2019).
43 */
44
45#include "math/ConvexHull.hpp"
46
47#include <algorithm>
48#include <fstream>
49#include <iostream>
50#include <iterator>
51#include <list>
52
53#ifdef IS_MPI
54#include <mpi.h>
55#endif
56
57#include "math/qhull.hpp"
58#include "utils/simError.h"
59
60#ifdef HAVE_QHULL
61using namespace OpenMD;
62using namespace std;
63
64ConvexHull::ConvexHull() : Hull(), options_("qhull FA Qt Pp QJ"), dim_(3) {}
65
66void ConvexHull::computeHull(vector<StuntDouble*> bodydoubles) {
67#ifdef HAVE_QHULL_REENTRANT
68 qhT qh_qh;
69 qhT* qh = &qh_qh;
70 QHULL_LIB_CHECK
71#endif
72
73 int numpoints = bodydoubles.size();
74
75 Triangles_.clear();
76
77 vertexT *vertex, **vertexp;
78 facetT* facet;
79 setT* vertices;
80 int curlong, totlong;
81
82 vector<double> ptArray(numpoints * dim_);
83
84 // Copy the positon vector into a points vector for qhull.
85 vector<StuntDouble*>::iterator SD;
86 int i = 0;
87
88 for (SD = bodydoubles.begin(); SD != bodydoubles.end(); ++SD) {
89 Vector3d pos = (*SD)->getPos();
90 ptArray[dim_ * i] = pos.x();
91 ptArray[dim_ * i + 1] = pos.y();
92 ptArray[dim_ * i + 2] = pos.z();
93 i++;
94 }
95
96 /* Clean up memory from previous convex hull calculations */
97 boolT ismalloc = False;
98
99 /* compute the hull for our local points (or all the points for single
100 processor versions) */
101#ifdef HAVE_QHULL_REENTRANT
102 qh_init_A(qh, NULL, NULL, stderr, 0, NULL);
103 int exitcode = setjmp(qh->errexit);
104 if (!exitcode) {
105 qh->NOerrexit = False;
106 qh_initflags(qh, const_cast<char*>(options_.c_str()));
107 qh_init_B(qh, &ptArray[0], numpoints, dim_, ismalloc);
108 qh_qhull(qh);
109 qh_check_output(qh);
110 exitcode = qh_ERRnone;
111 qh->NOerrexit = True;
112 } else {
113 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
114 "ConvexHull: Qhull failed to compute convex hull");
115 painCave.isFatal = 1;
116 simError();
117 }
118#else
119 qh_init_A(NULL, NULL, stderr, 0, NULL);
120 int exitcode = setjmp(qh errexit);
121 if (!exitcode) {
122 qh_initflags(const_cast<char*>(options_.c_str()));
123 qh_init_B(&ptArray[0], numpoints, dim_, ismalloc);
124 qh_qhull();
125 qh_check_output();
126 exitcode = qh_ERRnone;
127 qh NOerrexit = True;
128 } else {
129 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
130 "ConvexHull: Qhull failed to compute convex hull");
131 painCave.isFatal = 1;
132 simError();
133 }
134#endif
135
136#ifdef IS_MPI
137 // If we are doing the mpi version, set up some vectors for data communication
138
139 int nproc;
140 int myrank;
141
142 MPI_Comm_size(MPI_COMM_WORLD, &nproc);
143 MPI_Comm_rank(MPI_COMM_WORLD, &myrank);
144
145 int localHullSites = 0;
146
147 vector<int> hullSitesOnProc(nproc, 0);
148 vector<int> coordsOnProc(nproc, 0);
149 vector<int> displacements(nproc, 0);
150 vector<int> vectorDisplacements(nproc, 0);
151
152 vector<double> coords;
153 vector<double> vels;
154 vector<int> indexMap;
155 vector<double> masses;
156
157 FORALLvertices {
158 localHullSites++;
159
160#ifdef HAVE_QHULL_REENTRANT
161 int idx = qh_pointid(qh, vertex->point);
162#else
163 int idx = qh_pointid(vertex->point);
164#endif
165
166 indexMap.push_back(idx);
167
168 coords.push_back(ptArray[dim_ * idx]);
169 coords.push_back(ptArray[dim_ * idx + 1]);
170 coords.push_back(ptArray[dim_ * idx + 2]);
171
172 StuntDouble* sd = bodydoubles[idx];
173
174 Vector3d vel = sd->getVel();
175 vels.push_back(vel.x());
176 vels.push_back(vel.y());
177 vels.push_back(vel.z());
178
179 masses.push_back(sd->getMass());
180 }
181
182 MPI_Allgather(&localHullSites, 1, MPI_INT, &hullSitesOnProc[0], 1, MPI_INT,
183 MPI_COMM_WORLD);
184
185 int globalHullSites = 0;
186 for (int iproc = 0; iproc < nproc; iproc++) {
187 globalHullSites += hullSitesOnProc[iproc];
188 coordsOnProc[iproc] = dim_ * hullSitesOnProc[iproc];
189 }
190
191 displacements[0] = 0;
192 vectorDisplacements[0] = 0;
193
194 for (int iproc = 1; iproc < nproc; iproc++) {
195 displacements[iproc] =
196 displacements[iproc - 1] + hullSitesOnProc[iproc - 1];
197 vectorDisplacements[iproc] =
198 vectorDisplacements[iproc - 1] + coordsOnProc[iproc - 1];
199 }
200
201 vector<double> globalCoords(dim_ * globalHullSites);
202 vector<double> globalVels(dim_ * globalHullSites);
203 vector<double> globalMasses(globalHullSites);
204
205 int count = coordsOnProc[myrank];
206
207 MPI_Allgatherv(&coords[0], count, MPI_DOUBLE, &globalCoords[0],
208 &coordsOnProc[0], &vectorDisplacements[0], MPI_DOUBLE,
209 MPI_COMM_WORLD);
210
211 MPI_Allgatherv(&vels[0], count, MPI_DOUBLE, &globalVels[0], &coordsOnProc[0],
212 &vectorDisplacements[0], MPI_DOUBLE, MPI_COMM_WORLD);
213
214 MPI_Allgatherv(&masses[0], localHullSites, MPI_DOUBLE, &globalMasses[0],
215 &hullSitesOnProc[0], &displacements[0], MPI_DOUBLE,
216 MPI_COMM_WORLD);
217
218 // Free previous hull
219#ifdef HAVE_QHULL_REENTRANT
220 qh_freeqhull(qh, !qh_ALL);
221 qh_memfreeshort(qh, &curlong, &totlong);
222#else
223 qh_freeqhull(!qh_ALL);
224 qh_memfreeshort(&curlong, &totlong);
225#endif
226 if (curlong || totlong) {
227 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
228 "ConvexHull: qhull internal warning:\n"
229 "\tdid not free %d bytes of long memory (%d pieces)",
230 totlong, curlong);
231 painCave.isFatal = 1;
232 simError();
233 }
234
235#ifdef HAVE_QHULL_REENTRANT
236 qh_init_A(qh, NULL, NULL, stderr, 0, NULL);
237 exitcode = setjmp(qh->errexit);
238 if (!exitcode) {
239 qh->NOerrexit = False;
240 qh_initflags(qh, const_cast<char*>(options_.c_str()));
241 qh_init_B(qh, &globalCoords[0], globalHullSites, dim_, ismalloc);
242 qh_qhull(qh);
243 qh_check_output(qh);
244 exitcode = qh_ERRnone;
245 qh->NOerrexit = True;
246 } else {
247 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
248 "ConvexHull: Qhull failed to compute convex hull");
249 painCave.isFatal = 1;
250 simError();
251 }
252#else
253 qh_init_A(NULL, NULL, stderr, 0, NULL);
254 exitcode = setjmp(qh errexit);
255 if (!exitcode) {
256 qh NOerrexit = False;
257 qh_initflags(const_cast<char*>(options_.c_str()));
258 qh_init_B(&globalCoords[0], globalHullSites, dim_, ismalloc);
259 qh_qhull();
260 qh_check_output();
261 exitcode = qh_ERRnone;
262 qh NOerrexit = True; /* no more setjmp */
263 } else {
264 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
265 "ConvexHull: Qhull failed to compute global convex hull");
266 painCave.isFatal = 1;
267 simError();
268
269 } // qh_new_qhull
270#endif
271
272#endif
273 // commented out below, so comment out here also.
274 // intPoint = qh interior_point;
275 // RealType calcvol = 0.0;
276
277#ifdef HAVE_QHULL_REENTRANT
278 qh_triangulate(qh);
279#else
280 qh_triangulate();
281#endif
282
283 FORALLfacets {
284 Triangle face;
285 // Qhull sets the unit normal in facet->normal
286 Vector3d V3dNormal(facet->normal[0], facet->normal[1], facet->normal[2]);
287 face.setUnitNormal(V3dNormal);
288
289#ifdef HAVE_QHULL_REENTRANT
290 RealType faceArea = qh_facetarea(qh, facet);
291 face.setArea(faceArea);
292 vertices = qh_facet3vertex(qh, facet);
293 coordT* center = qh_getcenter(qh, vertices);
294#else
295 RealType faceArea = qh_facetarea(facet);
296 face.setArea(faceArea);
297 vertices = qh_facet3vertex(facet);
298 coordT* center = qh_getcenter(vertices);
299#endif
300 Vector3d V3dCentroid(center[0], center[1], center[2]);
301 face.setCentroid(V3dCentroid);
302
303 Vector3d faceVel = V3Zero;
304 Vector3d p[3];
305 RealType faceMass = 0.0;
306
307 int ver = 0;
308
309 FOREACHvertex_(vertices) {
310#ifdef HAVE_QHULL_REENTRANT
311 int id = qh_pointid(qh, vertex->point);
312#else
313 int id = qh_pointid(vertex->point);
314#endif
315 p[ver][0] = vertex->point[0];
316 p[ver][1] = vertex->point[1];
317 p[ver][2] = vertex->point[2];
318 Vector3d vel;
319 RealType mass;
320
321#ifdef IS_MPI
322 vel = Vector3d(globalVels[dim_ * id], globalVels[dim_ * id + 1],
323 globalVels[dim_ * id + 2]);
324 mass = globalMasses[id];
325
326 // localID will be between 0 and hullSitesOnProc[myrank] if we
327 // own this guy.
328
329 int localID = id - displacements[myrank];
330
331 if (localID >= 0 && localID < hullSitesOnProc[myrank]) {
332 face.addVertexSD(bodydoubles[indexMap[localID]]);
333 } else {
334 face.addVertexSD(NULL);
335 }
336#else
337 vel = bodydoubles[id]->getVel();
338 mass = bodydoubles[id]->getMass();
339 face.addVertexSD(bodydoubles[id]);
340#endif
341 faceVel = faceVel + vel;
342 faceMass = faceMass + mass;
343 ver++;
344 } // Foreachvertex
345
346 face.addVertices(p[0], p[1], p[2]);
347 face.setFacetMass(faceMass);
348 face.setFacetVelocity(faceVel / RealType(3.0));
349 /*
350 RealType comparea = face.computeArea();
351 realT calcarea = qh_facetarea (facet);
352 Vector3d V3dCompNorm = -face.computeUnitNormal();
353 RealType thisOffset = ((0.0-p[0][0])*V3dCompNorm[0] +
354 (0.0-p[0][1])*V3dCompNorm[1] + (0.0-p[0][2])*V3dCompNorm[2]); RealType dist
355 = facet->offset + intPoint[0]*V3dNormal[0] + intPoint[1]*V3dNormal[1] +
356 intPoint[2]*V3dNormal[2]; cout
357 << "facet offset and computed offset: " << facet->offset << " " <<
358 thisOffset << endl; calcvol += -dist*comparea/qh hull_dim;
359 */
360 Triangles_.push_back(face);
361#ifdef HAVE_QHULL_REENTRANT
362 qh_settempfree(qh, &vertices);
363#else
364 qh_settempfree(&vertices);
365#endif
366 } // FORALLfacets
367
368#ifdef HAVE_QHULL_REENTRANT
369 qh_getarea(qh, qh->facet_list);
370 volume_ = qh->totvol;
371 area_ = qh->totarea;
372 qh_freeqhull(qh, !qh_ALL);
373 qh_memfreeshort(qh, &curlong, &totlong);
374
375#else
376 qh_getarea(qh facet_list);
377 volume_ = qh totvol;
378 area_ = qh totarea;
379
380 qh_freeqhull(!qh_ALL);
381 qh_memfreeshort(&curlong, &totlong);
382#endif
383 if (curlong || totlong) {
384 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
385 "ConvexHull: qhull internal warning:\n"
386 "\tdid not free %d bytes of long memory (%d pieces)",
387 totlong, curlong);
388 painCave.isFatal = 1;
389 simError();
390 }
391}
392
393// void ConvexHull::printHull(const string& geomFileName) {
394
395// #ifdef IS_MPI
396// if (worldRank == 0) {
397// #endif
398// FILE *newGeomFile;
399
400// //create new .omd file based on old .omd file
401// newGeomFile = fopen(geomFileName.c_str(), "w");
402// #ifdef HAVE_QHULL_REENTRANT
403// qh_findgood_all(qh, qh->facet_list);
404// #else
405// qh_findgood_all(qh facet_list);
406// #endif
407// for (int i = 0; i < qh_PRINTEND; i++)
408// #ifdef HAVE_QHULL_REENTRANT
409// qh_printfacets(qh, newGeomFile, qh->PRINTout[i], qh->facet_list, NULL,
410// !qh_ALL);
411// #else
412// qh_printfacets(newGeomFile, qh PRINTout[i], qh facet_list, NULL,
413// !qh_ALL);
414// #endif
415// fclose(newGeomFile);
416// #ifdef IS_MPI
417// }
418// #endif
419// }
420#endif // QHULL
"Don't move, or you're dead! Stand up! Captain, we've got them!"
Vector3d getVel()
Returns the current velocity of this stuntDouble.
RealType getMass()
Returns the mass of this stuntDouble.
Triangle provides geometric data to OpenMD.
Definition Triangle.hpp:65
Real & z()
Returns reference of the third element of Vector3.
Definition Vector3.hpp:120
Real & x()
Returns reference of the first element of Vector3.
Definition Vector3.hpp:96
Real & y()
Returns reference of the second element of Vector3.
Definition Vector3.hpp:108
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.