1 |
< |
/* Copyright (c) 2008, 2009 The University of Notre Dame. All Rights Reserved. |
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. 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 |
8 |
> |
* 1. Redistributions of source code must retain the above copyright |
9 |
|
* notice, this list of conditions and the following disclaimer. |
10 |
|
* |
11 |
< |
* 3. Redistributions in binary form must reproduce the above copyright |
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. |
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 |
|
* ConvexHull.cpp |
42 |
|
* |
43 |
< |
* Purpose: To calculate convexhull, hull volume libqhull. |
44 |
< |
* |
45 |
< |
* Created by Charles F. Vardeman II on 11 Dec 2006. |
46 |
< |
* @author Charles F. Vardeman II |
47 |
< |
* @version $Id: ConvexHull.cpp,v 1.14 2009-10-12 20:11:29 chuckv Exp $ |
48 |
< |
* |
43 |
> |
* Purpose: To calculate a convex hull. |
44 |
|
*/ |
45 |
|
|
46 |
|
/* Standard includes independent of library */ |
47 |
+ |
|
48 |
|
#include <iostream> |
49 |
|
#include <fstream> |
50 |
|
#include <list> |
53 |
|
#include "math/ConvexHull.hpp" |
54 |
|
#include "utils/simError.h" |
55 |
|
|
56 |
+ |
#ifdef IS_MPI |
57 |
+ |
#include <mpi.h> |
58 |
+ |
#endif |
59 |
|
|
60 |
< |
using namespace oopse; |
60 |
> |
#include "math/qhull.hpp" |
61 |
|
|
62 |
< |
/* CGAL version of convex hull first then QHULL */ |
63 |
< |
#ifdef HAVE_CGAL |
64 |
< |
//#include <CGAL/Homogeneous.h> |
66 |
< |
#include <CGAL/basic.h> |
67 |
< |
//#include <CGAL/Simple_cartesian.h> |
68 |
< |
#include <CGAL/Cartesian.h> |
69 |
< |
#include <CGAL/Origin.h> |
70 |
< |
#include <CGAL/Exact_predicates_exact_constructions_kernel.h> |
71 |
< |
#include <CGAL/Convex_hull_traits_3.h> |
72 |
< |
#include <CGAL/convex_hull_3.h> |
73 |
< |
#include <CGAL/Polyhedron_traits_with_normals_3.h> |
74 |
< |
#include <CGAL/Polyhedron_3.h> |
75 |
< |
#include <CGAL/double.h> |
76 |
< |
#include <CGAL/number_utils.h> |
62 |
> |
#ifdef HAVE_QHULL |
63 |
> |
using namespace OpenMD; |
64 |
> |
using namespace std; |
65 |
|
|
66 |
+ |
ConvexHull::ConvexHull() : Hull(), dim_(3), options_("qhull Qt Pp") { |
67 |
+ |
} |
68 |
|
|
69 |
< |
//#include <CGAL/Quotient.h> |
70 |
< |
#include <CGAL/MP_Float.h> |
71 |
< |
//#include <CGAL/Lazy_exact_nt.h> |
69 |
> |
void ConvexHull::computeHull(vector<StuntDouble*> bodydoubles) { |
70 |
> |
|
71 |
> |
int numpoints = bodydoubles.size(); |
72 |
> |
|
73 |
> |
Triangles_.clear(); |
74 |
> |
|
75 |
> |
vertexT *vertex, **vertexp; |
76 |
> |
facetT *facet; |
77 |
> |
setT *vertices; |
78 |
> |
int curlong, totlong; |
79 |
|
|
80 |
+ |
Vector3d boxMax; |
81 |
+ |
Vector3d boxMin; |
82 |
+ |
|
83 |
+ |
vector<double> ptArray(numpoints*dim_); |
84 |
|
|
85 |
+ |
// Copy the positon vector into a points vector for qhull. |
86 |
+ |
vector<StuntDouble*>::iterator SD; |
87 |
+ |
int i = 0; |
88 |
+ |
|
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, "ConvexHull: Qhull failed to compute convex hull"); |
106 |
+ |
painCave.isFatal = 1; |
107 |
+ |
simError(); |
108 |
+ |
|
109 |
+ |
} //qh_new_qhull |
110 |
|
|
85 |
– |
typedef CGAL::MP_Float RT; |
86 |
– |
//typedef double RT; |
87 |
– |
//typedef CGAL::Homogeneous<RT> K; |
88 |
– |
typedef CGAL::Exact_predicates_exact_constructions_kernel K; |
89 |
– |
typedef K::Vector_3 Vector_3; |
90 |
– |
//typedef CGAL::Convex_hull_traits_3<K> Traits; |
91 |
– |
typedef CGAL::Polyhedron_traits_with_normals_3<K> Traits; |
92 |
– |
//typedef Traits::Polyhedron_3 Polyhedron_3; |
93 |
– |
typedef CGAL::Polyhedron_3<Traits> Polyhedron_3; |
94 |
– |
typedef K::Point_3 Point_3; |
111 |
|
|
112 |
+ |
#ifdef IS_MPI |
113 |
+ |
//If we are doing the mpi version, set up some vectors for data communication |
114 |
+ |
|
115 |
+ |
int nproc = MPI::COMM_WORLD.Get_size(); |
116 |
+ |
int myrank = MPI::COMM_WORLD.Get_rank(); |
117 |
+ |
int localHullSites = 0; |
118 |
|
|
119 |
< |
typedef Polyhedron_3::HalfedgeDS HalfedgeDS; |
120 |
< |
typedef Polyhedron_3::Facet_iterator Facet_iterator; |
121 |
< |
typedef Polyhedron_3::Halfedge_around_facet_circulator Halfedge_facet_circulator; |
122 |
< |
typedef Polyhedron_3::Halfedge_handle Halfedge_handle; |
101 |
< |
typedef Polyhedron_3::Facet_iterator Facet_iterator; |
102 |
< |
typedef Polyhedron_3::Plane_iterator Plane_iterator; |
103 |
< |
typedef Polyhedron_3::Vertex_iterator Vertex_iterator; |
104 |
< |
typedef Polyhedron_3::Vertex_handle Vertex_handle; |
105 |
< |
typedef Polyhedron_3::Point_iterator Point_iterator; |
106 |
< |
|
119 |
> |
vector<int> hullSitesOnProc(nproc, 0); |
120 |
> |
vector<int> coordsOnProc(nproc, 0); |
121 |
> |
vector<int> displacements(nproc, 0); |
122 |
> |
vector<int> vectorDisplacements(nproc, 0); |
123 |
|
|
124 |
+ |
vector<double> coords; |
125 |
+ |
vector<double> vels; |
126 |
+ |
vector<int> indexMap; |
127 |
+ |
vector<double> masses; |
128 |
|
|
129 |
< |
class Enriched_Point_3 : public K::Point_3{ |
130 |
< |
public: |
131 |
< |
Enriched_Point_3(double x,double y,double z) : K::Point_3(x,y,z), yupMyPoint(false), mySD(NULL) {} |
129 |
> |
FORALLvertices{ |
130 |
> |
localHullSites++; |
131 |
> |
|
132 |
> |
int idx = qh_pointid(vertex->point); |
133 |
|
|
134 |
< |
bool isMyPoint() const{ return yupMyPoint; } |
114 |
< |
void myPoint(){ yupMyPoint = true; } |
115 |
< |
void setSD(StuntDouble* SD){mySD = SD;} |
116 |
< |
StuntDouble* getStuntDouble(){return mySD;} |
117 |
< |
private: |
118 |
< |
bool yupMyPoint; |
119 |
< |
StuntDouble* mySD; |
134 |
> |
indexMap.push_back(idx); |
135 |
|
|
136 |
< |
}; |
136 |
> |
coords.push_back(ptArray[dim_ * idx]); |
137 |
> |
coords.push_back(ptArray[dim_ * idx + 1]); |
138 |
> |
coords.push_back(ptArray[dim_ * idx + 2]); |
139 |
|
|
140 |
+ |
StuntDouble* sd = bodydoubles[idx]; |
141 |
|
|
142 |
+ |
Vector3d vel = sd->getVel(); |
143 |
+ |
vels.push_back(vel.x()); |
144 |
+ |
vels.push_back(vel.y()); |
145 |
+ |
vels.push_back(vel.z()); |
146 |
|
|
147 |
+ |
masses.push_back(sd->getMass()); |
148 |
+ |
} |
149 |
|
|
150 |
+ |
MPI::COMM_WORLD.Allgather(&localHullSites, 1, MPI::INT, &hullSitesOnProc[0], |
151 |
+ |
1, MPI::INT); |
152 |
|
|
153 |
< |
// compare Point_3's... used in setting up the STL map from points to indices |
154 |
< |
template <typename Pt3> |
155 |
< |
struct Point_3_comp { |
156 |
< |
bool operator() (const Pt3 & p, const Pt3 & q) const { |
131 |
< |
return CGAL::lexicographically_xyz_smaller(p,q); // this is defined inline & hence we had to create fn object & not ptrfun |
153 |
> |
int globalHullSites = 0; |
154 |
> |
for (int iproc = 0; iproc < nproc; iproc++){ |
155 |
> |
globalHullSites += hullSitesOnProc[iproc]; |
156 |
> |
coordsOnProc[iproc] = dim_ * hullSitesOnProc[iproc]; |
157 |
|
} |
133 |
– |
}; |
158 |
|
|
159 |
< |
// coordinate-based hashing inefficient but can we do better if pts are copied? |
160 |
< |
typedef std::map<Point_3, StuntDouble* ,Point_3_comp<Point_3> > ptMapType; |
159 |
> |
displacements[0] = 0; |
160 |
> |
vectorDisplacements[0] = 0; |
161 |
> |
|
162 |
> |
for (int iproc = 1; iproc < nproc; iproc++){ |
163 |
> |
displacements[iproc] = displacements[iproc-1] + hullSitesOnProc[iproc-1]; |
164 |
> |
vectorDisplacements[iproc] = vectorDisplacements[iproc-1] + coordsOnProc[iproc-1]; |
165 |
> |
} |
166 |
|
|
167 |
< |
#ifdef IS_MPI |
168 |
< |
struct { |
169 |
< |
double x,y,z; |
141 |
< |
} surfacePt; |
142 |
< |
#endif |
167 |
> |
vector<double> globalCoords(dim_ * globalHullSites); |
168 |
> |
vector<double> globalVels(dim_ * globalHullSites); |
169 |
> |
vector<double> globalMasses(globalHullSites); |
170 |
|
|
171 |
< |
ConvexHull::ConvexHull() : Hull(){ |
172 |
< |
//If we are doing the mpi version, set up some vectors for data communication |
173 |
< |
#ifdef IS_MPI |
171 |
> |
int count = coordsOnProc[myrank]; |
172 |
> |
|
173 |
> |
MPI::COMM_WORLD.Allgatherv(&coords[0], count, MPI::DOUBLE, &globalCoords[0], |
174 |
> |
&coordsOnProc[0], &vectorDisplacements[0], |
175 |
> |
MPI::DOUBLE); |
176 |
|
|
177 |
+ |
MPI::COMM_WORLD.Allgatherv(&vels[0], count, MPI::DOUBLE, &globalVels[0], |
178 |
+ |
&coordsOnProc[0], &vectorDisplacements[0], |
179 |
+ |
MPI::DOUBLE); |
180 |
|
|
181 |
< |
nproc_ = MPI::COMM_WORLD.Get_size(); |
182 |
< |
myrank_ = MPI::COMM_WORLD.Get_rank(); |
183 |
< |
NstoProc_ = new int[nproc_]; |
152 |
< |
vecdispls_ = new int[nproc_]; |
153 |
< |
displs_ = new int[nproc_]; |
154 |
< |
// Create a surface point type in MPI to send |
155 |
< |
surfacePtType = MPI::DOUBLE.Create_contiguous(3); |
156 |
< |
surfacePtType.Commit(); |
157 |
< |
|
181 |
> |
MPI::COMM_WORLD.Allgatherv(&masses[0], localHullSites, MPI::DOUBLE, |
182 |
> |
&globalMasses[0], &hullSitesOnProc[0], |
183 |
> |
&displacements[0], MPI::DOUBLE); |
184 |
|
|
185 |
< |
#endif |
186 |
< |
} |
187 |
< |
|
188 |
< |
void ConvexHull::computeHull(std::vector<StuntDouble*> bodydoubles) |
189 |
< |
{ |
190 |
< |
|
191 |
< |
std::vector<Enriched_Point_3> points; |
192 |
< |
ptMapType myMap; |
193 |
< |
Point_iterator hc; |
168 |
< |
|
169 |
< |
// Copy the positon vector into a points vector for cgal. |
170 |
< |
std::vector<StuntDouble*>::iterator SD; |
171 |
< |
|
172 |
< |
for (SD =bodydoubles.begin(); SD != bodydoubles.end(); ++SD) |
173 |
< |
{ |
174 |
< |
Vector3d pos = (*SD)->getPos(); |
175 |
< |
Enriched_Point_3* pt = new Enriched_Point_3(pos.x(),pos.y(),pos.z()); |
176 |
< |
pt->setSD(*SD); |
177 |
< |
points.push_back(*pt); |
178 |
< |
// myMap[pt]=(*SD); |
179 |
< |
} |
180 |
< |
|
181 |
< |
// define object to hold convex hull |
182 |
< |
CGAL::Object ch_object_; |
183 |
< |
Polyhedron_3 polyhedron; |
184 |
< |
|
185 |
< |
// compute convex hull |
186 |
< |
|
187 |
< |
std::vector<Enriched_Point_3>::iterator testpt; |
188 |
< |
|
189 |
< |
|
190 |
< |
|
191 |
< |
CGAL::convex_hull_3(points.begin(), points.end(), polyhedron); |
192 |
< |
|
193 |
< |
|
194 |
< |
|
195 |
< |
Ns_ = polyhedron.size_of_vertices(); |
196 |
< |
|
197 |
< |
#ifdef IS_MPI |
198 |
< |
/* Gather an array of the number of verticies on each processor */ |
199 |
< |
|
200 |
< |
|
201 |
< |
surfacePtsGlobal_.clear(); |
202 |
< |
surfacePtsLocal_.clear(); |
203 |
< |
|
204 |
< |
MPI::COMM_WORLD.Allgather(&Ns_,1,MPI::INT,&NstoProc_[0],1,MPI::INT); |
205 |
< |
|
206 |
< |
for (int i = 0; i < nproc_; i++){ |
207 |
< |
Nsglobal_ += NstoProc_[i]; |
185 |
> |
// Free previous hull |
186 |
> |
qh_freeqhull(!qh_ALL); |
187 |
> |
qh_memfreeshort(&curlong, &totlong); |
188 |
> |
if (curlong || totlong) { |
189 |
> |
sprintf(painCave.errMsg, "ConvexHull: qhull internal warning:\n" |
190 |
> |
"\tdid not free %d bytes of long memory (%d pieces)", |
191 |
> |
totlong, curlong); |
192 |
> |
painCave.isFatal = 1; |
193 |
> |
simError(); |
194 |
|
} |
209 |
– |
/*Reminder ideally, we would like to reserve size for the vectors here*/ |
210 |
– |
surfacePtsLocal_.reserve(Ns_); |
211 |
– |
surfacePtsGlobal_.resize(Nsglobal_); |
212 |
– |
// std::fill(surfacePtsGlobal_.begin(),surfacePtsGlobal_.end(),0); |
213 |
– |
|
214 |
– |
/* Build a displacements array */ |
215 |
– |
for (int i = 1; i < nproc_; i++){ |
216 |
– |
vecdispls_[i] = vecdispls_[i-1] + NstoProc_[i-1]; |
217 |
– |
} |
195 |
|
|
196 |
< |
int noffset = vecdispls_[myrank_]; |
197 |
< |
/* gather the potential hull */ |
198 |
< |
|
199 |
< |
|
200 |
< |
for (hc =polyhedron.points_begin();hc != polyhedron.points_end(); ++hc){ |
201 |
< |
Point_3 mypoint = *hc; |
202 |
< |
surfacePt_ mpiSurfacePt; |
203 |
< |
mpiSurfacePt.x = CGAL::to_double(mypoint.x()); |
204 |
< |
mpiSurfacePt.y = CGAL::to_double(mypoint.y()); |
228 |
< |
mpiSurfacePt.z = CGAL::to_double(mypoint.z()); |
229 |
< |
surfacePtsLocal_.push_back(mpiSurfacePt); |
230 |
< |
} |
196 |
> |
if (qh_new_qhull(dim_, globalHullSites, &globalCoords[0], ismalloc, |
197 |
> |
const_cast<char *>(options_.c_str()), NULL, stderr)){ |
198 |
> |
|
199 |
> |
sprintf(painCave.errMsg, |
200 |
> |
"ConvexHull: Qhull failed to compute global convex hull"); |
201 |
> |
painCave.isFatal = 1; |
202 |
> |
simError(); |
203 |
> |
|
204 |
> |
} //qh_new_qhull |
205 |
|
|
232 |
– |
MPI::COMM_WORLD.Allgatherv(&surfacePtsLocal_[0],Ns_,surfacePtType,&surfacePtsGlobal_[0],NstoProc_,vecdispls_,surfacePtType); |
233 |
– |
std::vector<surfacePt_>::iterator spt; |
234 |
– |
std::vector<Enriched_Point_3> gblpoints; |
235 |
– |
|
236 |
– |
int mine = 0; |
237 |
– |
int pointidx = 0; |
238 |
– |
for (spt = surfacePtsGlobal_.begin(); spt != surfacePtsGlobal_.end(); ++spt) |
239 |
– |
{ |
240 |
– |
surfacePt_ thispos = *spt; |
241 |
– |
Enriched_Point_3 ept(thispos.x,thispos.y,thispos.z); |
242 |
– |
if (mine >= noffset && mine < noffset + Ns_){ |
243 |
– |
ept.myPoint(); |
244 |
– |
ept.setSD(points[pointidx].getStuntDouble()); |
245 |
– |
pointidx++; |
246 |
– |
} |
247 |
– |
gblpoints.push_back(ept); |
248 |
– |
|
249 |
– |
mine++; |
250 |
– |
} |
251 |
– |
|
252 |
– |
/* Compute the global hull */ |
253 |
– |
polyhedron.clear(); |
254 |
– |
CGAL::convex_hull_3(gblpoints.begin(), gblpoints.end(), polyhedron); |
255 |
– |
|
256 |
– |
|
206 |
|
#endif |
207 |
+ |
// commented out below, so comment out here also. |
208 |
+ |
// intPoint = qh interior_point; |
209 |
+ |
// RealType calcvol = 0.0; |
210 |
+ |
FORALLfacets { |
211 |
+ |
Triangle face; |
212 |
+ |
//Qhull sets the unit normal in facet->normal |
213 |
+ |
Vector3d V3dNormal(facet->normal[0], facet->normal[1], facet->normal[2]); |
214 |
+ |
face.setUnitNormal(V3dNormal); |
215 |
+ |
|
216 |
+ |
RealType faceArea = qh_facetarea(facet); |
217 |
+ |
face.setArea(faceArea); |
218 |
+ |
|
219 |
+ |
vertices = qh_facet3vertex(facet); |
220 |
+ |
|
221 |
+ |
coordT *center = qh_getcenter(vertices); |
222 |
+ |
Vector3d V3dCentroid(center[0], center[1], center[2]); |
223 |
+ |
face.setCentroid(V3dCentroid); |
224 |
|
|
225 |
+ |
Vector3d faceVel = V3Zero; |
226 |
+ |
Vector3d p[3]; |
227 |
+ |
RealType faceMass = 0.0; |
228 |
|
|
229 |
< |
|
261 |
< |
/* Loop over all of the surface triangles and build data structures for atoms and normals*/ |
262 |
< |
Facet_iterator j; |
263 |
< |
area_ = 0; |
264 |
< |
for ( j = polyhedron.facets_begin(); j !=polyhedron.facets_end(); ++j) { |
265 |
< |
Halfedge_handle h = j->halfedge(); |
229 |
> |
int ver = 0; |
230 |
|
|
231 |
< |
Point_3 r0=h->vertex()->point(); |
232 |
< |
Point_3 r1=h->next()->vertex()->point(); |
233 |
< |
Point_3 r2=h->next()->next()->vertex()->point(); |
231 |
> |
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 |
|
|
239 |
< |
Point_3* pr0 = &r0; |
240 |
< |
Point_3* pr1 = &r1; |
241 |
< |
Point_3* pr2 = &r2; |
242 |
< |
|
243 |
< |
Enriched_Point_3* er0 = static_cast<Enriched_Point_3*>(pr0); |
276 |
< |
Enriched_Point_3* er1 = static_cast<Enriched_Point_3*>(pr1); |
277 |
< |
Enriched_Point_3* er2 = static_cast<Enriched_Point_3*>(pr2); |
278 |
< |
|
279 |
< |
// StuntDouble* sd = er0->getStuntDouble(); |
280 |
< |
std::cerr << "sd globalIndex = " << to_double(er0->x()) << "\n"; |
281 |
< |
|
282 |
< |
Point_3 thisCentroid = CGAL::centroid(r0,r1,r2); |
283 |
< |
|
284 |
< |
Vector_3 normal = CGAL::cross_product(r1-r0,r2-r0); |
285 |
< |
|
286 |
< |
Triangle* face = new Triangle(); |
287 |
< |
Vector3d V3dNormal(CGAL::to_double(normal.x()),CGAL::to_double(normal.y()),CGAL::to_double(normal.z())); |
288 |
< |
Vector3d V3dCentroid(CGAL::to_double(thisCentroid.x()),CGAL::to_double(thisCentroid.y()),CGAL::to_double(thisCentroid.z())); |
289 |
< |
face->setNormal(V3dNormal); |
290 |
< |
face->setCentroid(V3dCentroid); |
291 |
< |
RealType faceArea = 0.5*V3dNormal.length(); |
292 |
< |
face->setArea(faceArea); |
293 |
< |
area_ += faceArea; |
294 |
< |
Triangles_.push_back(face); |
295 |
< |
// ptMapType::const_iterator locn=myMap.find(mypoint); |
296 |
< |
// int myIndex = locn->second; |
297 |
< |
|
298 |
< |
} |
299 |
< |
|
300 |
< |
|
301 |
< |
|
302 |
< |
|
303 |
< |
} |
304 |
< |
void ConvexHull::printHull(const std::string& geomFileName) |
305 |
< |
{ |
306 |
< |
/* |
307 |
< |
std::ofstream newGeomFile; |
308 |
< |
|
309 |
< |
//create new .md file based on old .md file |
310 |
< |
newGeomFile.open("testhull.off"); |
311 |
< |
|
312 |
< |
// Write polyhedron in Object File Format (OFF). |
313 |
< |
CGAL::set_ascii_mode( std::cout); |
314 |
< |
newGeomFile << "OFF" << std::endl << polyhedron.size_of_vertices() << ' ' |
315 |
< |
<< polyhedron.size_of_facets() << " 0" << std::endl; |
316 |
< |
std::copy( polyhedron.points_begin(), polyhedron.points_end(), |
317 |
< |
std::ostream_iterator<Point_3>( newGeomFile, "\n")); |
318 |
< |
for ( Facet_iterator i = polyhedron.facets_begin(); i != polyhedron.facets_end(); ++i) { |
319 |
< |
Halfedge_facet_circulator j = i->facet_begin(); |
320 |
< |
// Facets in polyhedral surfaces are at least triangles. |
321 |
< |
CGAL_assertion( CGAL::circulator_size(j) >= 3); |
322 |
< |
newGeomFile << CGAL::circulator_size(j) << ' '; |
323 |
< |
do { |
324 |
< |
newGeomFile << ' ' << std::distance(polyhedron.vertices_begin(), j->vertex()); |
325 |
< |
} while ( ++j != i->facet_begin()); |
326 |
< |
newGeomFile << std::endl; |
327 |
< |
} |
328 |
< |
|
329 |
< |
newGeomFile.close(); |
330 |
< |
*/ |
331 |
< |
/* |
332 |
< |
std::ofstream newGeomFile; |
333 |
< |
|
334 |
< |
//create new .md file based on old .md file |
335 |
< |
newGeomFile.open(geomFileName.c_str()); |
336 |
< |
|
337 |
< |
// Write polyhedron in Object File Format (OFF). |
338 |
< |
CGAL::set_ascii_mode( std::cout); |
339 |
< |
newGeomFile << "OFF" << std::endl << ch_polyhedron.size_of_vertices() << ' ' |
340 |
< |
<< ch_polyhedron.size_of_facets() << " 0" << std::endl; |
341 |
< |
std::copy( ch_polyhedron.points_begin(), ch_polyhedron.points_end(), |
342 |
< |
std::ostream_iterator<Point_3>( newGeomFile, "\n")); |
343 |
< |
for ( Facet_iterator i = ch_polyhedron.facets_begin(); i != ch_polyhedron.facets_end(); ++i) |
344 |
< |
{ |
345 |
< |
Halfedge_facet_circulator j = i->facet_begin(); |
346 |
< |
// Facets in polyhedral surfaces are at least triangles. |
347 |
< |
CGAL_assertion( CGAL::circulator_size(j) >= 3); |
348 |
< |
newGeomFile << CGAL::circulator_size(j) << ' '; |
349 |
< |
do |
350 |
< |
{ |
351 |
< |
newGeomFile << ' ' << std::distance(ch_polyhedron.vertices_begin(), j->vertex()); |
352 |
< |
} |
353 |
< |
while ( ++j != i->facet_begin()); |
354 |
< |
newGeomFile << std::endl; |
355 |
< |
} |
356 |
< |
|
357 |
< |
newGeomFile.close(); |
358 |
< |
*/ |
359 |
< |
|
360 |
< |
} |
239 |
> |
#ifdef IS_MPI |
240 |
> |
vel = Vector3d(globalVels[dim_ * id], |
241 |
> |
globalVels[dim_ * id + 1], |
242 |
> |
globalVels[dim_ * id + 2]); |
243 |
> |
mass = globalMasses[id]; |
244 |
|
|
245 |
+ |
// localID will be between 0 and hullSitesOnProc[myrank] if we |
246 |
+ |
// own this guy. |
247 |
|
|
248 |
+ |
int localID = id - displacements[myrank]; |
249 |
|
|
250 |
|
|
251 |
< |
|
252 |
< |
|
253 |
< |
|
251 |
> |
if (localID >= 0 && localID < hullSitesOnProc[myrank]){ |
252 |
> |
face.addVertexSD(bodydoubles[indexMap[localID]]); |
253 |
> |
}else{ |
254 |
> |
face.addVertexSD(NULL); |
255 |
> |
} |
256 |
|
#else |
257 |
< |
#ifdef HAVE_QHULL |
258 |
< |
/* Old options Qt Qu Qg QG0 FA */ |
259 |
< |
/* More old opts Qc Qi Pp*/ |
260 |
< |
ConvexHull::ConvexHull() : Hull(), dim_(3), options_("qhull Qt Pp"), Ns_(200), nTriangles_(0) { |
261 |
< |
//If we are doing the mpi version, set up some vectors for data communication |
262 |
< |
#ifdef IS_MPI |
257 |
> |
vel = bodydoubles[id]->getVel(); |
258 |
> |
mass = bodydoubles[id]->getMass(); |
259 |
> |
face.addVertexSD(bodydoubles[id]); |
260 |
> |
#endif |
261 |
> |
faceVel = faceVel + vel; |
262 |
> |
faceMass = faceMass + mass; |
263 |
> |
ver++; |
264 |
> |
} //Foreachvertex |
265 |
|
|
266 |
< |
|
267 |
< |
nproc_ = MPI::COMM_WORLD.Get_size(); |
268 |
< |
myrank_ = MPI::COMM_WORLD.Get_rank(); |
269 |
< |
NstoProc_ = new int[nproc_]; |
270 |
< |
vecdispls_ = new int[nproc_]; |
271 |
< |
vecNstoProc_ = new int[nproc_]; |
272 |
< |
displs_ = new int[nproc_]; |
273 |
< |
|
274 |
< |
// Create a surface point type in MPI to send |
275 |
< |
//surfacePtType = MPI::DOUBLE.Create_contiguous(3); |
276 |
< |
// surfacePtType.Commit(); |
277 |
< |
|
278 |
< |
|
389 |
< |
#endif |
390 |
< |
} |
391 |
< |
|
392 |
< |
|
393 |
< |
|
394 |
< |
void ConvexHull::computeHull(std::vector<StuntDouble*> bodydoubles) |
395 |
< |
{ |
396 |
< |
|
397 |
< |
std::vector<int> surfaceIDs; |
398 |
< |
std::vector<int> surfaceIDsGlobal; |
399 |
< |
std::vector<int> localPtsMap; |
400 |
< |
int numpoints = bodydoubles.size(); |
401 |
< |
|
402 |
< |
//coordT* pt_array; |
403 |
< |
coordT* surfpt_array; |
404 |
< |
vertexT *vertex, **vertexp; |
405 |
< |
facetT *facet; |
406 |
< |
setT *vertices; |
407 |
< |
int curlong,totlong; |
408 |
< |
int id; |
409 |
< |
|
410 |
< |
coordT *point,**pointp; |
411 |
< |
|
412 |
< |
|
413 |
< |
FILE *outdummy = NULL; |
414 |
< |
FILE *errdummy = NULL; |
415 |
< |
|
416 |
< |
//pt_array = (coordT*) malloc(sizeof(coordT) * (numpoints * dim_)); |
417 |
< |
|
418 |
< |
// double* ptArray = new double[numpoints * 3]; |
419 |
< |
std::vector<double> ptArray(numpoints*3); |
420 |
< |
std::vector<bool> isSurfaceID(numpoints); |
421 |
< |
|
422 |
< |
// Copy the positon vector into a points vector for qhull. |
423 |
< |
std::vector<StuntDouble*>::iterator SD; |
424 |
< |
int i = 0; |
425 |
< |
for (SD =bodydoubles.begin(); SD != bodydoubles.end(); ++SD) |
426 |
< |
{ |
427 |
< |
Vector3d pos = (*SD)->getPos(); |
428 |
< |
|
429 |
< |
ptArray[dim_ * i] = pos.x(); |
430 |
< |
ptArray[dim_ * i + 1] = pos.y(); |
431 |
< |
ptArray[dim_ * i + 2] = pos.z(); |
432 |
< |
i++; |
433 |
< |
} |
434 |
< |
|
435 |
< |
|
436 |
< |
|
437 |
< |
|
438 |
< |
|
439 |
< |
|
440 |
< |
boolT ismalloc = False; |
441 |
< |
/* Clean up memory from previous convex hull calculations*/ |
442 |
< |
|
443 |
< |
Triangles_.clear(); |
444 |
< |
surfaceSDs_.clear(); |
445 |
< |
surfaceSDs_.reserve(Ns_); |
446 |
< |
|
447 |
< |
if (qh_new_qhull(dim_, numpoints, &ptArray[0], ismalloc, |
448 |
< |
const_cast<char *>(options_.c_str()), NULL, stderr)) { |
449 |
< |
|
450 |
< |
sprintf(painCave.errMsg, "ConvexHull: Qhull failed to compute convex hull"); |
451 |
< |
painCave.isFatal = 1; |
452 |
< |
simError(); |
453 |
< |
|
454 |
< |
} //qh_new_qhull |
455 |
< |
|
456 |
< |
|
457 |
< |
#ifdef IS_MPI |
458 |
< |
std::vector<double> localPts; |
459 |
< |
std::vector<double> localVel; |
460 |
< |
std::vector<double> localMass; |
461 |
< |
int localPtArraySize; |
462 |
< |
|
463 |
< |
|
464 |
< |
std::fill(isSurfaceID.begin(),isSurfaceID.end(),false); |
465 |
< |
|
466 |
< |
|
467 |
< |
FORALLfacets { |
468 |
< |
|
469 |
< |
if (!facet->simplicial){ |
470 |
< |
// should never happen with Qt |
471 |
< |
sprintf(painCave.errMsg, "ConvexHull: non-simplicaial facet detected"); |
472 |
< |
painCave.isFatal = 1; |
473 |
< |
simError(); |
474 |
< |
} |
475 |
< |
|
476 |
< |
|
477 |
< |
vertices = qh_facet3vertex(facet); |
478 |
< |
FOREACHvertex_(vertices){ |
479 |
< |
id = qh_pointid(vertex->point); |
480 |
< |
|
481 |
< |
if( !isSurfaceID[id] ){ |
482 |
< |
isSurfaceID[id] = true; |
483 |
< |
} |
484 |
< |
} |
266 |
> |
face.addVertices(p[0], p[1], p[2]); |
267 |
> |
face.setFacetMass(faceMass); |
268 |
> |
face.setFacetVelocity(faceVel / RealType(3.0)); |
269 |
> |
/* |
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 |
> |
cout << "facet offset and computed offset: " << facet->offset << " " << thisOffset << endl; |
276 |
> |
calcvol += -dist*comparea/qh hull_dim; |
277 |
> |
*/ |
278 |
> |
Triangles_.push_back(face); |
279 |
|
qh_settempfree(&vertices); |
486 |
– |
|
487 |
– |
} //FORALLfacets |
280 |
|
|
281 |
< |
|
490 |
< |
|
491 |
< |
|
492 |
< |
int idx = 0; |
493 |
< |
int nIsIts = 0; |
494 |
< |
FORALLvertices { |
495 |
< |
idx = qh_pointid(vertex->point); |
496 |
< |
localPts.push_back(ptArray[dim_ * idx]); |
497 |
< |
localPts.push_back(ptArray[dim_ * idx + 1]); |
498 |
< |
localPts.push_back(ptArray[dim_ * idx + 2]); |
499 |
< |
|
500 |
< |
Vector3d vel = bodydoubles[idx]->getVel(); |
501 |
< |
localVel.push_back(vel.x()); |
502 |
< |
localVel.push_back(vel.y()); |
503 |
< |
localVel.push_back(vel.z()); |
504 |
< |
|
505 |
< |
|
506 |
< |
RealType bdmass = bodydoubles[idx]->getMass(); |
507 |
< |
localMass.push_back(bdmass); |
508 |
< |
|
509 |
< |
localPtsMap.push_back(idx); |
510 |
< |
|
511 |
< |
} |
512 |
< |
|
513 |
< |
|
514 |
< |
localPtArraySize = int(localPts.size()/3.0); |
515 |
< |
|
516 |
< |
MPI::COMM_WORLD.Allgather(&localPtArraySize,1,MPI::INT,&NstoProc_[0],1,MPI::INT); |
281 |
> |
} //FORALLfacets |
282 |
|
|
283 |
< |
Nsglobal_=0; |
284 |
< |
for (int i = 0; i < nproc_; i++){ |
285 |
< |
Nsglobal_ += NstoProc_[i]; |
521 |
< |
vecNstoProc_[i] = NstoProc_[i]*3; |
522 |
< |
} |
523 |
< |
|
524 |
< |
|
525 |
< |
int nglobalPts = Nsglobal_*3; |
526 |
< |
|
283 |
> |
qh_getarea(qh facet_list); |
284 |
> |
volume_ = qh totvol; |
285 |
> |
area_ = qh totarea; |
286 |
|
|
287 |
< |
std::vector<double> globalPts(nglobalPts); |
288 |
< |
std::vector<double> globalVel(nglobalPts); |
289 |
< |
std::vector<double> globalMass(Nsglobal_); |
290 |
< |
|
291 |
< |
|
292 |
< |
|
293 |
< |
isSurfaceID.resize(nglobalPts); |
294 |
< |
|
295 |
< |
|
296 |
< |
std::fill(globalPts.begin(),globalPts.end(),0.0); |
297 |
< |
|
539 |
< |
vecdispls_[0] = 0; |
540 |
< |
/* Build a displacements array */ |
541 |
< |
for (int i = 1; i < nproc_; i++){ |
542 |
< |
vecdispls_[i] = vecdispls_[i-1] + vecNstoProc_[i-1]; |
543 |
< |
} |
544 |
< |
|
545 |
< |
displs_[0] = 0; |
546 |
< |
for (int i = 1; i < nproc_; i++){ |
547 |
< |
displs_[i] = displs_[i-1] + NstoProc_[i-1]; |
548 |
< |
} |
549 |
< |
|
550 |
< |
int noffset = vecdispls_[myrank_]; |
551 |
< |
/* gather the potential hull */ |
552 |
< |
|
553 |
< |
MPI::COMM_WORLD.Allgatherv(&localPts[0],localPtArraySize*3,MPI::DOUBLE,&globalPts[0],&vecNstoProc_[0],&vecdispls_[0],MPI::DOUBLE); |
554 |
< |
MPI::COMM_WORLD.Allgatherv(&localVel[0],localPtArraySize*3,MPI::DOUBLE,&globalVel[0],&vecNstoProc_[0],&vecdispls_[0],MPI::DOUBLE); |
555 |
< |
MPI::COMM_WORLD.Allgatherv(&localMass[0],localPtArraySize,MPI::DOUBLE,&globalMass[0],&NstoProc_[0],&displs_[0],MPI::DOUBLE); |
556 |
< |
|
557 |
< |
/* |
558 |
< |
int tmpidx = 0; |
559 |
< |
|
560 |
< |
if (myrank_ == 0){ |
561 |
< |
for (i = 0; i < nglobalPts-3; i++){ |
562 |
< |
std::cout << "Au " << globalPts[tmpidx] << " " << globalPts[tmpidx+1] << " " << globalPts[tmpidx +2] << std::endl; |
563 |
< |
tmpidx = tmpidx + 3; |
287 |
> |
int index = 0; |
288 |
> |
FORALLvertices { |
289 |
> |
Vector3d point(vertex->point[0], vertex->point[1], vertex->point[2]); |
290 |
> |
if (index == 0) { |
291 |
> |
boxMax = point; |
292 |
> |
boxMin = point; |
293 |
> |
} else { |
294 |
> |
for (int i = 0; i < 3; i++) { |
295 |
> |
boxMax[i] = max(boxMax[i], point[i]); |
296 |
> |
boxMin[i] = min(boxMin[i], point[i]); |
297 |
> |
} |
298 |
|
} |
299 |
+ |
index++; |
300 |
|
} |
301 |
< |
*/ |
302 |
< |
|
303 |
< |
// Free previous hull |
301 |
> |
boundingBox_ = Mat3x3d(0.0); |
302 |
> |
boundingBox_(0,0) = boxMax[0] - boxMin[0]; |
303 |
> |
boundingBox_(1,1) = boxMax[1] - boxMin[1]; |
304 |
> |
boundingBox_(2,2) = boxMax[2] - boxMin[2]; |
305 |
> |
|
306 |
|
qh_freeqhull(!qh_ALL); |
307 |
|
qh_memfreeshort(&curlong, &totlong); |
308 |
< |
if (curlong || totlong) |
309 |
< |
std::cerr << "qhull internal warning (main): did not free %d bytes of long memory (%d pieces) " |
310 |
< |
<< totlong << curlong << std::endl; |
308 |
> |
if (curlong || totlong) { |
309 |
> |
sprintf(painCave.errMsg, "ConvexHull: qhull internal warning:\n" |
310 |
> |
"\tdid not free %d bytes of long memory (%d pieces)", |
311 |
> |
totlong, curlong); |
312 |
> |
painCave.isFatal = 1; |
313 |
> |
simError(); |
314 |
> |
} |
315 |
> |
} |
316 |
|
|
317 |
< |
if (qh_new_qhull(dim_, Nsglobal_, &globalPts[0], ismalloc, |
318 |
< |
const_cast<char *>(options_.c_str()), NULL, stderr)){ |
577 |
< |
|
578 |
< |
sprintf(painCave.errMsg, "ConvexHull: Qhull failed to compute global convex hull"); |
579 |
< |
painCave.isFatal = 1; |
580 |
< |
simError(); |
581 |
< |
|
582 |
< |
} //qh_new_qhull |
583 |
< |
|
584 |
< |
#endif |
585 |
< |
|
586 |
< |
|
587 |
< |
|
588 |
< |
|
589 |
< |
|
590 |
< |
|
591 |
< |
unsigned int nf = qh num_facets; |
592 |
< |
|
593 |
< |
/* Build Surface SD list first */ |
594 |
< |
|
595 |
< |
std::fill(isSurfaceID.begin(),isSurfaceID.end(),false); |
596 |
< |
|
597 |
< |
FORALLfacets { |
598 |
< |
|
599 |
< |
if (!facet->simplicial){ |
600 |
< |
// should never happen with Qt |
601 |
< |
sprintf(painCave.errMsg, "ConvexHull: non-simplicaial facet detected"); |
602 |
< |
painCave.isFatal = 1; |
603 |
< |
simError(); |
604 |
< |
} //simplicical |
605 |
< |
|
606 |
< |
Triangle face; |
607 |
< |
Vector3d V3dNormal(facet->normal[0],facet->normal[1],facet->normal[2]); |
608 |
< |
face.setNormal(V3dNormal); |
609 |
< |
|
610 |
< |
|
611 |
< |
|
612 |
< |
//RealType faceArea = 0.5*V3dNormal.length(); |
613 |
< |
RealType faceArea = qh_facetarea(facet); |
614 |
< |
face.setArea(faceArea); |
615 |
< |
|
616 |
< |
|
617 |
< |
vertices = qh_facet3vertex(facet); |
618 |
< |
|
619 |
< |
coordT *center = qh_getcenter(vertices); |
620 |
< |
Vector3d V3dCentroid(center[0], center[1], center[2]); |
621 |
< |
face.setCentroid(V3dCentroid); |
622 |
< |
Vector3d faceVel = V3Zero; |
623 |
< |
Vector3d p[3]; |
624 |
< |
RealType faceMass = 0.0; |
625 |
< |
int ver = 0; |
626 |
< |
FOREACHvertex_(vertices){ |
627 |
< |
id = qh_pointid(vertex->point); |
628 |
< |
p[ver][0] = vertex->point[0]; |
629 |
< |
p[ver][1] = vertex->point[1]; |
630 |
< |
p[ver][2] = vertex->point[2]; |
631 |
< |
int localindex = id; |
317 |
> |
void ConvexHull::printHull(const string& geomFileName) { |
318 |
> |
|
319 |
|
#ifdef IS_MPI |
320 |
< |
Vector3d velVector(globalVel[dim_ * id],globalVel[dim_ * id + 1], globalVel[dim_ * id + 1]); |
634 |
< |
|
635 |
< |
faceVel = faceVel + velVector; |
636 |
< |
faceMass = faceMass + globalMass[id]; |
637 |
< |
if (id >= noffset/3 && id < (noffset + localPtArraySize)/3 ){ |
638 |
< |
localindex = localPtsMap[id-noffset/3]; |
639 |
< |
#else |
640 |
< |
faceVel = faceVel + bodydoubles[localindex]->getVel(); |
641 |
< |
faceMass = faceMass + bodydoubles[localindex]->getMass(); |
320 |
> |
if (worldRank == 0) { |
321 |
|
#endif |
322 |
< |
face.addVertexSD(bodydoubles[localindex]); |
644 |
< |
if( !isSurfaceID[id] ){ |
645 |
< |
isSurfaceID[id] = true; |
646 |
< |
#ifdef IS_MPI |
647 |
< |
|
648 |
< |
#endif |
649 |
< |
|
650 |
< |
surfaceSDs_.push_back(bodydoubles[localindex]); |
651 |
< |
|
652 |
< |
} //IF isSurfaceID |
653 |
< |
|
654 |
< |
#ifdef IS_MPI |
655 |
< |
|
656 |
< |
}else{ |
657 |
< |
face.addVertexSD(NULL); |
658 |
< |
} |
659 |
< |
#endif |
660 |
< |
ver++; |
661 |
< |
} //Foreachvertex |
662 |
< |
/* |
663 |
< |
if (!SETempty_(facet->coplanarset)){ |
664 |
< |
FOREACHpoint_(facet->coplanarset){ |
665 |
< |
id = qh_pointid(point); |
666 |
< |
surfaceSDs_.push_back(bodydoubles[id]); |
667 |
< |
} |
668 |
< |
} |
669 |
< |
*/ |
670 |
< |
face.addVertices(p[0],p[1],p[2]); |
671 |
< |
face.setFacetMass(faceMass); |
672 |
< |
face.setFacetVelocity(faceVel/3.0); |
673 |
< |
Triangles_.push_back(face); |
674 |
< |
qh_settempfree(&vertices); |
675 |
< |
|
676 |
< |
} //FORALLfacets |
677 |
< |
|
678 |
< |
/* |
679 |
< |
std::cout << surfaceSDs_.size() << std::endl; |
680 |
< |
for (SD = surfaceSDs_.begin(); SD != surfaceSDs_.end(); ++SD){ |
681 |
< |
Vector3d thisatom = (*SD)->getPos(); |
682 |
< |
std::cout << "Au " << thisatom.x() << " " << thisatom.y() << " " << thisatom.z() << std::endl; |
683 |
< |
} |
684 |
< |
*/ |
685 |
< |
|
686 |
< |
|
687 |
< |
|
688 |
< |
Ns_ = surfaceSDs_.size(); |
689 |
< |
nTriangles_ = Triangles_.size(); |
322 |
> |
FILE *newGeomFile; |
323 |
|
|
324 |
< |
qh_getarea(qh facet_list); |
325 |
< |
volume_ = qh totvol; |
326 |
< |
area_ = qh totarea; |
324 |
> |
//create new .md file based on old .md file |
325 |
> |
newGeomFile = fopen(geomFileName.c_str(), "w"); |
326 |
> |
qh_findgood_all(qh facet_list); |
327 |
> |
for (int i = 0; i < qh_PRINTEND; i++) |
328 |
> |
qh_printfacets(newGeomFile, qh PRINTout[i], qh facet_list, NULL, !qh_ALL); |
329 |
|
|
330 |
< |
|
331 |
< |
|
332 |
< |
qh_freeqhull(!qh_ALL); |
333 |
< |
qh_memfreeshort(&curlong, &totlong); |
699 |
< |
if (curlong || totlong) |
700 |
< |
std::cerr << "qhull internal warning (main): did not free %d bytes of long memory (%d pieces) " |
701 |
< |
<< totlong << curlong << std::endl; |
702 |
< |
|
703 |
< |
|
704 |
< |
|
330 |
> |
fclose(newGeomFile); |
331 |
> |
#ifdef IS_MPI |
332 |
> |
} |
333 |
> |
#endif |
334 |
|
} |
706 |
– |
|
707 |
– |
|
708 |
– |
|
709 |
– |
void ConvexHull::printHull(const std::string& geomFileName) |
710 |
– |
{ |
711 |
– |
|
712 |
– |
FILE *newGeomFile; |
713 |
– |
|
714 |
– |
//create new .md file based on old .md file |
715 |
– |
newGeomFile = fopen(geomFileName.c_str(), "w"); |
716 |
– |
qh_findgood_all(qh facet_list); |
717 |
– |
for (int i = 0; i < qh_PRINTEND; i++) |
718 |
– |
qh_printfacets(newGeomFile, qh PRINTout[i], qh facet_list, NULL, !qh_ALL); |
719 |
– |
|
720 |
– |
fclose(newGeomFile); |
721 |
– |
} |
335 |
|
#endif //QHULL |
723 |
– |
#endif //CGAL |
724 |
– |
|
725 |
– |
|
726 |
– |
|