ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/math/AlphaHull.cpp
Revision: 2056
Committed: Fri Feb 20 15:12:07 2015 UTC (10 years, 2 months ago) by gezelter
File size: 16346 byte(s)
Log Message:
Fixes to HullFinder (and by extension to RNEMD) for getSurfaceArea() call.
Adding Multipass Correlation Function (unused right now).

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, 234107 (2008).
38 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
39 * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
40 *
41 * AlphaHull.cpp
42 *
43 * Purpose: To calculate an alpha-shape hull.
44 */
45
46 #ifdef IS_MPI
47 #include <mpi.h>
48 #endif
49
50 /* Standard includes independent of library */
51
52 #include <iostream>
53 #include <fstream>
54 #include <list>
55 #include <algorithm>
56 #include <iterator>
57 #include "math/AlphaHull.hpp"
58 #include "utils/simError.h"
59
60 #include "math/qhull.hpp"
61
62 #ifdef HAVE_QHULL
63 using namespace std;
64 using namespace OpenMD;
65
66 double calculate_circumradius(pointT* p0, pointT* p1, pointT* p2, int dim);
67
68 AlphaHull::AlphaHull(double alpha) : Hull(), dim_(3), alpha_(alpha),
69 options_("qhull d QJ Tcv Pp") {
70 }
71
72 void AlphaHull::computeHull(vector<StuntDouble*> bodydoubles) {
73
74 int numpoints = bodydoubles.size();
75
76 Triangles_.clear();
77
78 vertexT *vertex;
79 facetT *facet, *neighbor;
80 pointT *interiorPoint;
81 int curlong, totlong;
82
83
84 vector<double> ptArray(numpoints*dim_);
85
86 // Copy the positon vector into a points vector for qhull.
87 vector<StuntDouble*>::iterator SD;
88 int i = 0;
89 for (SD =bodydoubles.begin(); SD != bodydoubles.end(); ++SD){
90 Vector3d pos = (*SD)->getPos();
91 ptArray[dim_ * i] = pos.x();
92 ptArray[dim_ * i + 1] = pos.y();
93 ptArray[dim_ * i + 2] = pos.z();
94 i++;
95 }
96
97 /* Clean up memory from previous convex hull calculations*/
98 boolT ismalloc = False;
99
100 /* compute the hull for our local points (or all the points for single
101 processor versions) */
102 if (qh_new_qhull(dim_, numpoints, &ptArray[0], ismalloc,
103 const_cast<char *>(options_.c_str()), NULL, stderr)) {
104
105 sprintf(painCave.errMsg, "AlphaHull: Qhull failed to compute convex hull");
106 painCave.isFatal = 1;
107 simError();
108
109 } //qh_new_qhull
110
111
112 #ifdef IS_MPI
113 //If we are doing the mpi version, set up some vectors for data communication
114
115 int nproc;
116 int myrank;
117 MPI_Comm_size( MPI_COMM_WORLD, &nproc);
118 MPI_Comm_rank( MPI_COMM_WORLD, &myrank);
119
120 int localHullSites = 0;
121
122 vector<int> hullSitesOnProc(nproc, 0);
123 vector<int> coordsOnProc(nproc, 0);
124 vector<int> displacements(nproc, 0);
125 vector<int> vectorDisplacements(nproc, 0);
126
127 vector<double> coords;
128 vector<double> vels;
129 vector<int> indexMap;
130 vector<double> masses;
131
132 FORALLvertices{
133 localHullSites++;
134
135 int idx = qh_pointid(vertex->point);
136
137 indexMap.push_back(idx);
138
139 coords.push_back(ptArray[dim_ * idx]);
140 coords.push_back(ptArray[dim_ * idx + 1]);
141 coords.push_back(ptArray[dim_ * idx + 2]);
142
143 StuntDouble* sd = bodydoubles[idx];
144
145 Vector3d vel = sd->getVel();
146 vels.push_back(vel.x());
147 vels.push_back(vel.y());
148 vels.push_back(vel.z());
149
150 masses.push_back(sd->getMass());
151 }
152
153 MPI_Allgather(&localHullSites, 1, MPI_INT, &hullSitesOnProc[0],
154 1, MPI_INT, MPI_COMM_WORLD);
155
156 int globalHullSites = 0;
157 for (int iproc = 0; iproc < nproc; iproc++){
158 globalHullSites += hullSitesOnProc[iproc];
159 coordsOnProc[iproc] = dim_ * hullSitesOnProc[iproc];
160 }
161
162 displacements[0] = 0;
163 vectorDisplacements[0] = 0;
164
165 for (int iproc = 1; iproc < nproc; iproc++){
166 displacements[iproc] = displacements[iproc-1] + hullSitesOnProc[iproc-1];
167 vectorDisplacements[iproc] = vectorDisplacements[iproc-1] + coordsOnProc[iproc-1];
168 }
169
170 vector<double> globalCoords(dim_ * globalHullSites);
171 vector<double> globalVels(dim_ * globalHullSites);
172 vector<double> globalMasses(globalHullSites);
173
174 int count = coordsOnProc[myrank];
175
176 MPI_Allgatherv(&coords[0], count, MPI_DOUBLE, &globalCoords[0],
177 &coordsOnProc[0], &vectorDisplacements[0],
178 MPI_DOUBLE, MPI_COMM_WORLD);
179
180 MPI_Allgatherv(&vels[0], count, MPI_DOUBLE, &globalVels[0],
181 &coordsOnProc[0], &vectorDisplacements[0],
182 MPI_DOUBLE, MPI_COMM_WORLD);
183
184 MPI_Allgatherv(&masses[0], localHullSites, MPI_DOUBLE,
185 &globalMasses[0], &hullSitesOnProc[0],
186 &displacements[0], MPI_DOUBLE, MPI_COMM_WORLD);
187
188 // Free previous hull
189 qh_freeqhull(!qh_ALL);
190 qh_memfreeshort(&curlong, &totlong);
191 if (curlong || totlong) {
192 sprintf(painCave.errMsg, "AlphaHull: qhull internal warning:\n"
193 "\tdid not free %d bytes of long memory (%d pieces)",
194 totlong, curlong);
195 painCave.isFatal = 1;
196 simError();
197 }
198
199 if (qh_new_qhull(dim_, globalHullSites, &globalCoords[0], ismalloc,
200 const_cast<char *>(options_.c_str()), NULL, stderr)){
201
202 sprintf(painCave.errMsg,
203 "AlphaHull: Qhull failed to compute global convex hull");
204 painCave.isFatal = 1;
205 simError();
206
207 } //qh_new_qhull
208
209 #endif
210
211 //Set facet->center as the Voronoi center
212 qh_setvoronoi_all();
213
214
215 int convexNumVert = qh_setsize(qh_facetvertices (qh facet_list, NULL, false));
216 //Insert all the sample points, because, even with alpha=0, the
217 //alpha shape/alpha complex will contain them.
218
219 // tri::Allocator<CMeshO>::AddVertices(pm.cm,convexNumVert);
220
221 /*ivp length is 'qh num_vertices' because each vertex is accessed through its ID whose range is
222 0<=qh_pointid(vertex->point)<qh num_vertices*/
223 // vector<tri::Allocator<CMeshO>::VertexPointer> ivp(qh num_vertices);
224 /*i=0;
225 FORALLvertices{
226 if ((*vertex).point){
227 // pm.cm.vert[i].P()[0] = (*vertex).point[0];
228 // pm.cm.vert[i].P()[1] = (*vertex).point[1];
229 //pm.cm.vert[i].P()[2] = (*vertex).point[2];
230 // ivp[qh_pointid(vertex->point)] = &pm.cm.vert[i];
231 i++;
232 }
233 }
234 */
235 //Set of alpha complex triangles for alphashape filtering
236 setT* set= qh_settemp(4* qh num_facets);
237
238 qh visit_id++;
239 int numFacets=0;
240 vector<vector <int> > facetlist;
241 interiorPoint = qh interior_point;
242 FORALLfacet_(qh facet_list) {
243 numFacets++;
244 if (!facet->upperdelaunay) {
245 //For all facets (that are tetrahedrons)calculate the radius of
246 //the empty circumsphere considering the distance between the
247 //circumcenter and a vertex of the facet
248 vertexT* vertex = (vertexT *)(facet->vertices->e[0].p);
249 double* center = facet->center;
250 double radius = qh_pointdist(vertex->point,center,dim_);
251
252 if (radius>alpha_) // if the facet is not good consider the ridges
253 {
254 //if calculating the alphashape, unmark the facet ('good' is
255 //used as 'marked').
256 facet->good=false;
257
258 //Compute each ridge (triangle) once and test the
259 //cironference radius with alpha
260 facet->visitid= qh visit_id;
261 qh_makeridges(facet);
262 ridgeT *ridge, **ridgep;
263 int goodTriangles=0;
264 FOREACHridge_(facet->ridges) {
265 neighbor= otherfacet_(ridge, facet);
266 if (( neighbor->visitid != qh visit_id)){
267 //Calculate the radius of the circumference
268 pointT* p0 = ((vertexT*) (ridge->vertices->e[0].p))->point;
269 pointT* p1 = ((vertexT*) (ridge->vertices->e[1].p))->point;
270 pointT* p2 = ((vertexT*) (ridge->vertices->e[2].p))->point;
271
272 radius = calculate_circumradius(p0,p1,p2, dim_);
273
274 if(radius <=alpha_){
275 goodTriangles++;
276 //save the triangle (ridge) for subsequent filtering
277 qh_setappend(&set, ridge);
278 }
279 }
280 }
281
282 //If calculating the alphashape, mark the facet('good' is
283 //used as 'marked'). This facet will have some triangles
284 //hidden by the facet's neighbor.
285 if(goodTriangles==4)
286 facet->good=true;
287
288 }
289 else //the facet is good. Put all the triangles of the
290 //tetrahedron in the mesh
291 {
292 //Compute each ridge (triangle) once
293 facet->visitid= qh visit_id;
294 //If calculating the alphashape, mark the facet('good' is
295 //used as 'marked'). This facet will have some triangles
296 //hidden by the facet's neighbor.
297 facet->good=true;
298 qh_makeridges(facet);
299 ridgeT *ridge, **ridgep;
300 FOREACHridge_(facet->ridges) {
301 neighbor= otherfacet_(ridge, facet);
302 if ((neighbor->visitid != qh visit_id)){
303 qh_setappend(&set, ridge);
304 }
305 }
306 }
307 }
308 }
309 //assert(numFacets== qh num_facets);
310
311 //Filter the triangles (only the ones on the boundary of the alpha
312 //complex) and build the mesh
313
314 int ridgesCount=0;
315
316 ridgeT *ridge, **ridgep;
317 FOREACHridge_(set) {
318 if ((!ridge->top->good || !ridge->bottom->good || ridge->top->upperdelaunay || ridge->bottom->upperdelaunay)){
319 // tri::Allocator<CMeshO>::FaceIterator fi=tri::Allocator<CMeshO>::AddFaces(pm.cm,1);
320 ridgesCount++;
321 int vertex_n, vertex_i;
322 Triangle face;
323
324 // Vector3d V3dNormal(facet->normal[0], facet->normal[1], facet->normal[2]);
325 //face.setNormal(V3dNormal);
326
327
328 //coordT *center = qh_getcenter(ridge->vertices);
329 //cout << "Centers are " << center[0] << " " <<center[1] << " " << center[2] << endl;
330 //Vector3d V3dCentroid(center[0], center[1], center[2]);
331 //face.setCentroid(V3dCentroid);
332
333
334 Vector3d faceVel = V3Zero;
335 Vector3d p[3];
336 RealType faceMass = 0.0;
337
338 int ver = 0;
339 vector<int> virtexlist;
340 FOREACHvertex_i_(ridge->vertices){
341 int id = qh_pointid(vertex->point);
342 p[ver][0] = vertex->point[0];
343 p[ver][1] = vertex->point[1];
344 p[ver][2] = vertex->point[2];
345 Vector3d vel;
346 RealType mass;
347 ver++;
348 virtexlist.push_back(id);
349 // cout << "Ridge: " << ridgesCount << " Vertex " << id << endl;
350
351 vel = bodydoubles[id]->getVel();
352 mass = bodydoubles[id]->getMass();
353 face.addVertexSD(bodydoubles[id]);
354
355
356 faceVel = faceVel + vel;
357 faceMass = faceMass + mass;
358 } //FOREACH Vertex
359 facetlist.push_back(virtexlist);
360 face.addVertices(p[0],p[1],p[2]);
361 face.setFacetMass(faceMass);
362 face.setFacetVelocity(faceVel / RealType(3.0));
363
364 RealType area = face.getArea();
365 area_ += area;
366 Vector3d normal = face.getUnitNormal();
367 RealType offset = ((0.0-p[0][0])*normal[0] + (0.0-p[0][1])*normal[1] + (0.0-p[0][2])*normal[2]);
368 RealType dist = normal[0] * interiorPoint[0] + normal[1]*interiorPoint[1] + normal[2]*interiorPoint[2];
369 volume_ += dist *area/qh hull_dim;
370
371 Triangles_.push_back(face);
372 }
373 }
374
375
376 //assert(pm.cm.fn == ridgesCount);
377 /*
378 std::cout <<"OFF"<<std::endl;
379 std::cout << bodydoubles.size() << " " << facetlist.size() << " " << 3*facetlist.size() << std::endl;
380 for (SD =bodydoubles.begin(); SD != bodydoubles.end(); ++SD){
381 Vector3d pos = (*SD)->getPos();
382 std::cout << pos.x() << " " << pos.y() << " " << pos.z() << std::endl;
383 }
384
385
386 std::vector<std::vector<int> >::iterator thisfacet;
387 std::vector<int>::iterator thisvertex;
388
389 for (thisfacet = facetlist.begin(); thisfacet != facetlist.end(); thisfacet++){
390 std::cout << (*thisfacet).size();
391 for (thisvertex = (*thisfacet).begin(); thisvertex != (*thisfacet).end(); thisvertex++){
392 std::cout << " " << *thisvertex;
393 }
394 std::cout << std::endl;
395 }
396 */
397
398
399
400 /*
401 FORALLfacets {
402 Triangle face;
403
404 Vector3d V3dNormal(facet->normal[0], facet->normal[1], facet->normal[2]);
405 face.setNormal(V3dNormal);
406
407 RealType faceArea = qh_facetarea(facet);
408 face.setArea(faceArea);
409
410 vertices = qh_facet3vertex(facet);
411
412 coordT *center = qh_getcenter(vertices);
413 Vector3d V3dCentroid(center[0], center[1], center[2]);
414 face.setCentroid(V3dCentroid);
415
416 Vector3d faceVel = V3Zero;
417 Vector3d p[3];
418 RealType faceMass = 0.0;
419
420 int ver = 0;
421
422 FOREACHvertex_(vertices){
423 int id = qh_pointid(vertex->point);
424 p[ver][0] = vertex->point[0];
425 p[ver][1] = vertex->point[1];
426 p[ver][2] = vertex->point[2];
427
428 Vector3d vel;
429 RealType mass;
430
431 #ifdef IS_MPI
432 vel = Vector3d(globalVels[dim_ * id],
433 globalVels[dim_ * id + 1],
434 globalVels[dim_ * id + 2]);
435 mass = globalMasses[id];
436
437 // localID will be between 0 and hullSitesOnProc[myrank] if we
438 // own this guy.
439
440 int localID = id - displacements[myrank];
441
442 if (localID >= 0 && localID < hullSitesOnProc[myrank])
443 face.addVertexSD(bodydoubles[indexMap[localID]]);
444
445 #else
446 vel = bodydoubles[id]->getVel();
447 mass = bodydoubles[id]->getMass();
448 face.addVertexSD(bodydoubles[id]);
449 #endif
450
451 faceVel = faceVel + vel;
452 faceMass = faceMass + mass;
453 ver++;
454 } //Foreachvertex
455
456 face.addVertices(p[0], p[1], p[2]);
457 face.setFacetMass(faceMass);
458 face.setFacetVelocity(faceVel/3.0);
459 Triangles_.push_back(face);
460 qh_settempfree(&vertices);
461
462 } //FORALLfacets
463 */
464 // qh_getarea(qh facet_list);
465 //volume_ = qh totvol;
466 // area_ = qh totarea;
467
468 qh_freeqhull(!qh_ALL);
469 qh_memfreeshort(&curlong, &totlong);
470 if (curlong || totlong) {
471 sprintf(painCave.errMsg, "AlphaHull: qhull internal warning:\n"
472 "\tdid not free %d bytes of long memory (%d pieces)",
473 totlong, curlong);
474 painCave.isFatal = 1;
475 simError();
476 }
477 }
478
479 void AlphaHull::printHull(const string& geomFileName) {
480
481 #ifdef IS_MPI
482 if (worldRank == 0) {
483 #endif
484 FILE *newGeomFile;
485
486 //create new .md file based on old .md file
487 newGeomFile = fopen(geomFileName.c_str(), "w");
488 qh_findgood_all(qh facet_list);
489 for (int i = 0; i < qh_PRINTEND; i++)
490 qh_printfacets(newGeomFile, qh PRINTout[i], qh facet_list, NULL, !qh_ALL);
491
492 fclose(newGeomFile);
493 #ifdef IS_MPI
494 }
495 #endif
496 }
497
498 double calculate_circumradius(pointT* p0,pointT* p1,pointT* p2, int dim){
499 coordT a = qh_pointdist(p0,p1,dim);
500 coordT b = qh_pointdist(p1,p2,dim);
501 coordT c = qh_pointdist(p2,p0,dim);
502
503 coordT sum =(a + b + c)*0.5;
504 coordT area = sum*(a+b-sum)*(a+c-sum)*(b+c-sum);
505 return (double) (a*b*c)/(4*sqrt(area));
506 }
507
508 #endif //QHULL

Properties

Name Value
svn:keywords Author Id Revision Date