OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
AlphaHull.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/AlphaHull.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 std;
62using namespace OpenMD;
63
64double calculate_circumradius(pointT* p0, pointT* p1, pointT* p2, int dim);
65
66AlphaHull::AlphaHull(double alpha) :
67 Hull(), dim_(4), alpha_(alpha), options_("qhull d QJ Tcv Pp") {}
68
69void AlphaHull::computeHull(vector<StuntDouble*> bodydoubles) {
70#ifdef HAVE_QHULL_REENTRANT
71 qhT qh_qh;
72 qhT* qh = &qh_qh;
73 QHULL_LIB_CHECK
74#endif
75
76 int numpoints = bodydoubles.size();
77
78 Triangles_.clear();
79
80 vertexT* vertex;
81 facetT *facet, *neighbor;
82 pointT* interiorPoint;
83 int curlong, totlong;
84
85 vector<double> ptArray(numpoints * dim_);
86
87 // Copy the positon vector into a points vector for qhull.
88 vector<StuntDouble*>::iterator SD;
89 int i = 0;
90
91 for (SD = bodydoubles.begin(); SD != bodydoubles.end(); ++SD) {
92 Vector3d pos = (*SD)->getPos();
93 ptArray[dim_ * i] = pos.x();
94 ptArray[dim_ * i + 1] = pos.y();
95 ptArray[dim_ * i + 2] = pos.z();
96 ptArray[dim_ * i + 3] = pos.lengthSquare();
97 i++;
98 }
99 /* Clean up memory from previous convex hull calculations*/
100 boolT ismalloc = False;
101
102 /* compute the hull for our local points (or all the points for single
103 processor versions) */
104#ifdef HAVE_QHULL_REENTRANT
105 qh_init_A(qh, NULL, NULL, stderr, 0, NULL);
106 int exitcode = setjmp(qh->errexit);
107 if (!exitcode) {
108 qh->NOerrexit = False;
109 qh_initflags(qh, const_cast<char*>(options_.c_str()));
110 qh_init_B(qh, &ptArray[0], numpoints, dim_, ismalloc);
111 qh_qhull(qh);
112 qh_check_output(qh);
113 exitcode = qh_ERRnone;
114 qh->NOerrexit = True;
115 } else {
116 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
117 "AlphaHull: Qhull failed to compute convex hull");
118 painCave.isFatal = 1;
119 simError();
120 }
121#else
122 qh_init_A(NULL, NULL, stderr, 0, NULL);
123 int exitcode = setjmp(qh errexit);
124 if (!exitcode) {
125 qh_initflags(const_cast<char*>(options_.c_str()));
126 qh_init_B(&ptArray[0], numpoints, dim_, ismalloc);
127 qh_qhull();
128 qh_check_output();
129 exitcode = qh_ERRnone;
130 qh NOerrexit = True;
131 } else {
132 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
133 "AlphaHull: Qhull failed to compute convex hull");
134 painCave.isFatal = 1;
135 simError();
136 }
137#endif
138
139#ifdef IS_MPI
140 // If we are doing the mpi version, set up some vectors for data communication
141
142 int nproc;
143 int myrank;
144 MPI_Comm_size(MPI_COMM_WORLD, &nproc);
145 MPI_Comm_rank(MPI_COMM_WORLD, &myrank);
146
147 int localHullSites = 0;
148
149 vector<int> hullSitesOnProc(nproc, 0);
150 vector<int> coordsOnProc(nproc, 0);
151 vector<int> displacements(nproc, 0);
152 vector<int> vectorDisplacements(nproc, 0);
153
154 vector<double> coords;
155 vector<double> vels;
156 vector<int> indexMap;
157 vector<double> masses;
158
159 FORALLvertices {
160 localHullSites++;
161
162#ifdef HAVE_QHULL_REENTRANT
163 int idx = qh_pointid(qh, vertex->point);
164#else
165 int idx = qh_pointid(vertex->point);
166#endif
167
168 indexMap.push_back(idx);
169
170 coords.push_back(ptArray[dim_ * idx]);
171 coords.push_back(ptArray[dim_ * idx + 1]);
172 coords.push_back(ptArray[dim_ * idx + 2]);
173 coords.push_back(ptArray[dim_ * idx + 3]);
174
175 StuntDouble* sd = bodydoubles[idx];
176
177 Vector3d vel = sd->getVel();
178 vels.push_back(vel.x());
179 vels.push_back(vel.y());
180 vels.push_back(vel.z());
181 vels.push_back(0.0);
182
183 masses.push_back(sd->getMass());
184 }
185
186 MPI_Allgather(&localHullSites, 1, MPI_INT, &hullSitesOnProc[0], 1, MPI_INT,
187 MPI_COMM_WORLD);
188
189 int globalHullSites = 0;
190 for (int iproc = 0; iproc < nproc; iproc++) {
191 globalHullSites += hullSitesOnProc[iproc];
192 coordsOnProc[iproc] = dim_ * hullSitesOnProc[iproc];
193 }
194
195 displacements[0] = 0;
196 vectorDisplacements[0] = 0;
197
198 for (int iproc = 1; iproc < nproc; iproc++) {
199 displacements[iproc] =
200 displacements[iproc - 1] + hullSitesOnProc[iproc - 1];
201 vectorDisplacements[iproc] =
202 vectorDisplacements[iproc - 1] + coordsOnProc[iproc - 1];
203 }
204
205 vector<double> globalCoords(dim_ * globalHullSites);
206 vector<double> globalVels(dim_ * globalHullSites);
207 vector<double> globalMasses(globalHullSites);
208
209 int count = coordsOnProc[myrank];
210
211 MPI_Allgatherv(&coords[0], count, MPI_DOUBLE, &globalCoords[0],
212 &coordsOnProc[0], &vectorDisplacements[0], MPI_DOUBLE,
213 MPI_COMM_WORLD);
214
215 MPI_Allgatherv(&vels[0], count, MPI_DOUBLE, &globalVels[0], &coordsOnProc[0],
216 &vectorDisplacements[0], MPI_DOUBLE, MPI_COMM_WORLD);
217
218 MPI_Allgatherv(&masses[0], localHullSites, MPI_DOUBLE, &globalMasses[0],
219 &hullSitesOnProc[0], &displacements[0], MPI_DOUBLE,
220 MPI_COMM_WORLD);
221
222 // Free previous hull
223#ifdef HAVE_QHULL_REENTRANT
224 qh_freeqhull(qh, !qh_ALL);
225 qh_memfreeshort(qh, &curlong, &totlong);
226#else
227 qh_freeqhull(!qh_ALL);
228 qh_memfreeshort(&curlong, &totlong);
229#endif
230 if (curlong || totlong) {
231 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
232 "AlphaHull: qhull internal warning:\n"
233 "\tdid not free %d bytes of long memory (%d pieces)",
234 totlong, curlong);
235 painCave.isFatal = 1;
236 simError();
237 }
238
239#ifdef HAVE_QHULL_REENTRANT
240 qh_init_A(qh, NULL, NULL, stderr, 0, NULL);
241 exitcode = setjmp(qh->errexit);
242 if (!exitcode) {
243 qh->NOerrexit = False;
244 qh_initflags(qh, const_cast<char*>(options_.c_str()));
245 qh_init_B(qh, &globalCoords[0], globalHullSites, dim_, ismalloc);
246 qh_qhull(qh);
247 qh_check_output(qh);
248 exitcode = qh_ERRnone;
249 qh->NOerrexit = True;
250 } else {
251 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
252 "AlphaHull: Qhull failed to compute convex hull");
253 painCave.isFatal = 1;
254 simError();
255 }
256#else
257 qh_init_A(NULL, NULL, stderr, 0, NULL);
258 exitcode = setjmp(qh errexit);
259 if (!exitcode) {
260 qh NOerrexit = False;
261 qh_initflags(const_cast<char*>(options_.c_str()));
262 qh_init_B(&globalCoords[0], globalHullSites, dim_, ismalloc);
263 qh_qhull();
264 qh_check_output();
265 exitcode = qh_ERRnone;
266 qh NOerrexit = True;
267 } else {
268 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
269 "AlphaHull: Qhull failed to compute convex hull");
270 painCave.isFatal = 1;
271 simError();
272 }
273#endif
274
275#endif
276
277 // Set facet->center as the Voronoi center
278#ifdef HAVE_QHULL_REENTRANT
279 qh_setvoronoi_all(qh);
280#else
281 qh_setvoronoi_all();
282#endif
283
284 // Set of alpha complex triangles for alphashape filtering
285 vector<vector<int>> facetlist;
286 int numFacets = 0;
287
288#ifdef HAVE_QHULL_REENTRANT
289 setT* set = qh_settemp(qh, 4 * qh->num_facets);
290 qh->visit_id++;
291 interiorPoint = qh->interior_point;
292#else
293 setT* set = qh_settemp(4 * qh num_facets);
294 qh visit_id++;
295 interiorPoint = qh interior_point;
296#endif
297
298#ifdef HAVE_QHULL_REENTRANT
299 FORALLfacet_(qh->facet_list) {
300#else
301 FORALLfacet_(qh facet_list) {
302#endif
303 numFacets++;
304 if (!facet->upperdelaunay) {
305 // For all facets (that are tetrahedrons)calculate the radius of
306 // the empty circumsphere considering the distance between the
307 // circumcenter and a vertex of the facet
308 vertexT* vertex = (vertexT*)(facet->vertices->e[0].p);
309
310#ifdef HAVE_QHULL_REENTRANT
311 double* center = qh_facetcenter(qh, facet->vertices);
312#else
313 double* center = qh_facetcenter(facet->vertices);
314#endif
315 double radius = qh_pointdist(center, vertex->point, dim_ - 1);
316 // If radius is bigger than alpha, remove the tetrahedron
317
318 if (radius > alpha_) {
319 // if calculating the alphashape, unmark the facet ('good' is
320 // used as 'marked').
321 facet->good = false;
322
323 // Compute each ridge (triangle) once and test the
324 // cironference radius with alpha
325#ifdef HAVE_QHULL_REENTRANT
326 facet->visitid = qh->visit_id;
327 qh_makeridges(qh, facet);
328#else
329 facet->visitid = qh visit_id;
330 qh_makeridges(facet);
331#endif
332 ridgeT *ridge, **ridgep;
333 int goodTriangles = 0;
334 FOREACHridge_(facet->ridges) {
335 neighbor = otherfacet_(ridge, facet);
336#ifdef HAVE_QHULL_REENTRANT
337 if ((neighbor->visitid != qh->visit_id)) {
338#else
339 if ((neighbor->visitid != qh visit_id)) {
340#endif
341 // Calculate the radius of the circumference
342 pointT* p0 = ((vertexT*)(ridge->vertices->e[0].p))->point;
343 pointT* p1 = ((vertexT*)(ridge->vertices->e[1].p))->point;
344 pointT* p2 = ((vertexT*)(ridge->vertices->e[2].p))->point;
345
346 radius = calculate_circumradius(p0, p1, p2, dim_ - 1);
347
348 if (radius <= alpha_) {
349 goodTriangles++;
350 // save the triangle (ridge) for subsequent filtering
351#ifdef HAVE_QHULL_REENTRANT
352 qh_setappend(qh, &set, ridge);
353#else
354 qh_setappend(&set, ridge);
355#endif
356 }
357 }
358 }
359
360 // If calculating the alphashape, mark the facet('good' is
361 // used as 'marked'). This facet will have some triangles
362 // hidden by the facet's neighbor.
363 if (goodTriangles == 4) facet->good = true;
364
365 } else // the facet is good. Put all the triangles of the
366 // tetrahedron in the mesh
367 {
368 // Compute each ridge (triangle) once
369#ifdef HAVE_QHULL_REENTRANT
370 facet->visitid = qh->visit_id;
371#else
372 facet->visitid = qh visit_id;
373#endif
374 // Mark the facet('good' is used as 'marked'). This facet
375 // will have some triangles hidden by the facet's neighbor.
376 facet->good = true;
377#ifdef HAVE_QHULL_REENTRANT
378 qh_makeridges(qh, facet);
379#else
380 qh_makeridges(facet);
381#endif
382 ridgeT *ridge, **ridgep;
383 FOREACHridge_(facet->ridges) {
384 neighbor = otherfacet_(ridge, facet);
385#ifdef HAVE_QHULL_REENTRANT
386 if ((neighbor->visitid != qh->visit_id)) {
387 qh_setappend(qh, &set, ridge);
388 }
389#else
390 if ((neighbor->visitid != qh visit_id)) { qh_setappend(&set, ridge); }
391#endif
392 }
393 }
394 }
395 }
396
397 int ridgesCount = 0;
398
399 ridgeT *ridge, **ridgep;
400 FOREACHridge_(set) {
401 if ((!ridge->top->good || !ridge->bottom->good ||
402 ridge->top->upperdelaunay || ridge->bottom->upperdelaunay)) {
403 ridgesCount++;
404 int vertex_n, vertex_i;
405 Triangle face;
406
407 Vector3d faceVel = V3Zero;
408 Vector3d p[3];
409 RealType faceMass = 0.0;
410
411 int ver = 0;
412 vector<int> vertexlist;
413
414#ifdef HAVE_QHULL_REENTRANT
415 FOREACHvertex_i_(qh, ridge->vertices) {
416#else
417 FOREACHvertex_i_(ridge->vertices) {
418#endif
419#ifdef HAVE_QHULL_REENTRANT
420 int id = qh_pointid(qh, vertex->point);
421#else
422 int id = qh_pointid(vertex->point);
423#endif
424 p[ver][0] = vertex->point[0];
425 p[ver][1] = vertex->point[1];
426 p[ver][2] = vertex->point[2];
427 Vector3d vel;
428 RealType mass;
429 ver++;
430 vertexlist.push_back(id);
431
432 vel = bodydoubles[id]->getVel();
433 mass = bodydoubles[id]->getMass();
434 face.addVertexSD(bodydoubles[id]);
435
436 faceVel = faceVel + vel;
437 faceMass = faceMass + mass;
438 } // FOREACH Vertex
439 facetlist.push_back(vertexlist);
440 face.addVertices(p[0], p[1], p[2]);
441 face.setFacetMass(faceMass);
442 face.setFacetVelocity(faceVel / RealType(3.0));
443
444 RealType area = face.getArea();
445 area_ += area;
446 Vector3d normal = face.getUnitNormal();
447 RealType dist = normal[0] * interiorPoint[0] +
448 normal[1] * interiorPoint[1] +
449 normal[2] * interiorPoint[2];
450#ifdef HAVE_QHULL_REENTRANT
451 volume_ += dist * area / qh->hull_dim;
452#else
453 volume_ += dist * area / qh hull_dim;
454#endif
455
456 Triangles_.push_back(face);
457 }
458 }
459
460#ifdef HAVE_QHULL_REENTRANT
461 qh_freeqhull(qh, !qh_ALL);
462 qh_memfreeshort(qh, &curlong, &totlong);
463#else
464 qh_freeqhull(!qh_ALL);
465 qh_memfreeshort(&curlong, &totlong);
466#endif
467 if (curlong || totlong) {
468 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
469 "AlphaHull: qhull internal warning:\n"
470 "\tdid not free %d bytes of long memory (%d pieces)",
471 totlong, curlong);
472 painCave.isFatal = 1;
473 simError();
474 }
475}
476
477double calculate_circumradius(pointT* p0, pointT* p1, pointT* p2, int dim) {
478 coordT a = qh_pointdist(p0, p1, dim);
479 coordT b = qh_pointdist(p1, p2, dim);
480 coordT c = qh_pointdist(p2, p0, dim);
481
482 coordT sum = (a + b + c) * 0.5;
483 coordT area = sum * (a + b - sum) * (a + c - sum) * (b + c - sum);
484 return (double)(a * b * c) / (4 * sqrt(area));
485}
486
487#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
Real lengthSquare()
Returns the squared length of this vector.
Definition Vector.hpp:399
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.