ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/math/ConvexHull.cpp
Revision: 3544
Committed: Mon Oct 19 17:44:18 2009 UTC (15 years, 6 months ago) by chuckv
File size: 20642 byte(s)
Log Message:
More MPI convex hull bug fixes. noffset corresponds to the correct array indexs. Still mismatch in forces between serial and parallel versions.

File Contents

# User Rev Content
1 chuckv 3535 /* Copyright (c) 2008, 2009 The University of Notre Dame. All Rights Reserved.
2 chuckv 3083 *
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 3137 * Purpose: To calculate convexhull, hull volume libqhull.
44 chuckv 3083 *
45     * Created by Charles F. Vardeman II on 11 Dec 2006.
46     * @author Charles F. Vardeman II
47 chuckv 3544 * @version $Id: ConvexHull.cpp,v 1.15 2009-10-19 17:44:18 chuckv Exp $
48 chuckv 3083 *
49     */
50    
51 chuckv 3392 /* Standard includes independent of library */
52 chuckv 3083 #include <iostream>
53     #include <fstream>
54 chuckv 3275 #include <list>
55     #include <algorithm>
56     #include <iterator>
57 chuckv 3141 #include "math/ConvexHull.hpp"
58 chuckv 3275 #include "utils/simError.h"
59 chuckv 3450
60    
61 chuckv 3083 using namespace oopse;
62    
63 chuckv 3392 /* CGAL version of convex hull first then QHULL */
64     #ifdef HAVE_CGAL
65 chuckv 3450 //#include <CGAL/Homogeneous.h>
66 chuckv 3392 #include <CGAL/basic.h>
67 chuckv 3450 //#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 chuckv 3392 #include <CGAL/Convex_hull_traits_3.h>
72     #include <CGAL/convex_hull_3.h>
73 chuckv 3450 #include <CGAL/Polyhedron_traits_with_normals_3.h>
74     #include <CGAL/Polyhedron_3.h>
75 chuckv 3392 #include <CGAL/double.h>
76 chuckv 3450 #include <CGAL/number_utils.h>
77 chuckv 3275
78 chuckv 3450
79     //#include <CGAL/Quotient.h>
80 chuckv 3392 #include <CGAL/MP_Float.h>
81 chuckv 3450 //#include <CGAL/Lazy_exact_nt.h>
82 chuckv 3275
83 chuckv 3450
84    
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 chuckv 3416 typedef K::Point_3 Point_3;
95 chuckv 3275
96    
97 chuckv 3450 typedef Polyhedron_3::HalfedgeDS HalfedgeDS;
98     typedef Polyhedron_3::Facet_iterator Facet_iterator;
99     typedef Polyhedron_3::Halfedge_around_facet_circulator Halfedge_facet_circulator;
100     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    
107 chuckv 3275
108 chuckv 3450
109     class Enriched_Point_3 : public K::Point_3{
110     public:
111     Enriched_Point_3(double x,double y,double z) : K::Point_3(x,y,z), yupMyPoint(false), mySD(NULL) {}
112    
113     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;
120    
121     };
122    
123    
124    
125    
126    
127     // compare Point_3's... used in setting up the STL map from points to indices
128     template <typename Pt3>
129     struct Point_3_comp {
130     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
132     }
133     };
134    
135     // coordinate-based hashing inefficient but can we do better if pts are copied?
136     typedef std::map<Point_3, StuntDouble* ,Point_3_comp<Point_3> > ptMapType;
137    
138     #ifdef IS_MPI
139     struct {
140     double x,y,z;
141     } surfacePt;
142     #endif
143    
144     ConvexHull::ConvexHull() : Hull(){
145     //If we are doing the mpi version, set up some vectors for data communication
146     #ifdef IS_MPI
147    
148    
149     nproc_ = MPI::COMM_WORLD.Get_size();
150     myrank_ = MPI::COMM_WORLD.Get_rank();
151     NstoProc_ = new int[nproc_];
152 chuckv 3473 vecdispls_ = new int[nproc_];
153     displs_ = new int[nproc_];
154 chuckv 3450 // Create a surface point type in MPI to send
155     surfacePtType = MPI::DOUBLE.Create_contiguous(3);
156     surfacePtType.Commit();
157    
158    
159     #endif
160     }
161    
162     void ConvexHull::computeHull(std::vector<StuntDouble*> bodydoubles)
163 chuckv 3392 {
164 chuckv 3450
165     std::vector<Enriched_Point_3> points;
166     ptMapType myMap;
167     Point_iterator hc;
168 chuckv 3392
169     // Copy the positon vector into a points vector for cgal.
170 chuckv 3450 std::vector<StuntDouble*>::iterator SD;
171    
172     for (SD =bodydoubles.begin(); SD != bodydoubles.end(); ++SD)
173 chuckv 3392 {
174 chuckv 3450 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 chuckv 3392 }
180    
181     // define object to hold convex hull
182 chuckv 3416 CGAL::Object ch_object_;
183     Polyhedron_3 polyhedron;
184    
185 chuckv 3392 // compute convex hull
186 chuckv 3450
187     std::vector<Enriched_Point_3>::iterator testpt;
188    
189    
190    
191     CGAL::convex_hull_3(points.begin(), points.end(), polyhedron);
192    
193 chuckv 3416
194 chuckv 3450
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];
208     }
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 chuckv 3473 vecdispls_[i] = vecdispls_[i-1] + NstoProc_[i-1];
217 chuckv 3450 }
218    
219 chuckv 3473 int noffset = vecdispls_[myrank_];
220 chuckv 3450 /* gather the potential hull */
221    
222    
223     for (hc =polyhedron.points_begin();hc != polyhedron.points_end(); ++hc){
224     Point_3 mypoint = *hc;
225     surfacePt_ mpiSurfacePt;
226     mpiSurfacePt.x = CGAL::to_double(mypoint.x());
227     mpiSurfacePt.y = CGAL::to_double(mypoint.y());
228     mpiSurfacePt.z = CGAL::to_double(mypoint.z());
229     surfacePtsLocal_.push_back(mpiSurfacePt);
230     }
231    
232 chuckv 3473 MPI::COMM_WORLD.Allgatherv(&surfacePtsLocal_[0],Ns_,surfacePtType,&surfacePtsGlobal_[0],NstoProc_,vecdispls_,surfacePtType);
233 chuckv 3450 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 chuckv 3416 }
251    
252 chuckv 3450 /* Compute the global hull */
253     polyhedron.clear();
254     CGAL::convex_hull_3(gblpoints.begin(), gblpoints.end(), polyhedron);
255    
256    
257     #endif
258    
259    
260    
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();
266    
267     Point_3 r0=h->vertex()->point();
268     Point_3 r1=h->next()->vertex()->point();
269     Point_3 r2=h->next()->next()->vertex()->point();
270    
271     Point_3* pr0 = &r0;
272     Point_3* pr1 = &r1;
273     Point_3* pr2 = &r2;
274    
275     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 chuckv 3416 /*
307 chuckv 3450 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 chuckv 3392 }
328 chuckv 3450
329     newGeomFile.close();
330 chuckv 3416 */
331 chuckv 3450 /*
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 chuckv 3392 }
361 chuckv 3275
362    
363    
364 chuckv 3450
365    
366    
367    
368 chuckv 3392 #else
369     #ifdef HAVE_QHULL
370     /* Old options Qt Qu Qg QG0 FA */
371 chuckv 3461 /* More old opts Qc Qi Pp*/
372 chuckv 3464 ConvexHull::ConvexHull() : Hull(), dim_(3), options_("qhull Qt Pp"), Ns_(200), nTriangles_(0) {
373 chuckv 3450 //If we are doing the mpi version, set up some vectors for data communication
374     #ifdef IS_MPI
375 chuckv 3275
376 chuckv 3450
377     nproc_ = MPI::COMM_WORLD.Get_size();
378     myrank_ = MPI::COMM_WORLD.Get_rank();
379     NstoProc_ = new int[nproc_];
380 chuckv 3473 vecdispls_ = new int[nproc_];
381     vecNstoProc_ = new int[nproc_];
382     displs_ = new int[nproc_];
383 chuckv 3450
384     // Create a surface point type in MPI to send
385 chuckv 3459 //surfacePtType = MPI::DOUBLE.Create_contiguous(3);
386     // surfacePtType.Commit();
387 chuckv 3450
388    
389     #endif
390     }
391    
392    
393    
394     void ConvexHull::computeHull(std::vector<StuntDouble*> bodydoubles)
395 chuckv 3392 {
396    
397 chuckv 3450 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 chuckv 3392
410 chuckv 3450 coordT *point,**pointp;
411    
412    
413 chuckv 3392 FILE *outdummy = NULL;
414     FILE *errdummy = NULL;
415    
416 chuckv 3450 //pt_array = (coordT*) malloc(sizeof(coordT) * (numpoints * dim_));
417    
418 chuckv 3459 // double* ptArray = new double[numpoints * 3];
419     std::vector<double> ptArray(numpoints*3);
420 chuckv 3450 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 chuckv 3392
435 chuckv 3450
436 chuckv 3392
437    
438    
439    
440 chuckv 3450 boolT ismalloc = False;
441 chuckv 3459 /* Clean up memory from previous convex hull calculations*/
442 chuckv 3465
443 chuckv 3450 Triangles_.clear();
444     surfaceSDs_.clear();
445 chuckv 3459 surfaceSDs_.reserve(Ns_);
446    
447     if (qh_new_qhull(dim_, numpoints, &ptArray[0], ismalloc,
448 chuckv 3450 const_cast<char *>(options_.c_str()), NULL, stderr)) {
449 chuckv 3275
450 chuckv 3450 sprintf(painCave.errMsg, "ConvexHull: Qhull failed to compute convex hull");
451 chuckv 3464 painCave.isFatal = 1;
452 chuckv 3450 simError();
453    
454     } //qh_new_qhull
455    
456    
457     #ifdef IS_MPI
458     std::vector<double> localPts;
459 chuckv 3464 std::vector<double> localVel;
460 chuckv 3473 std::vector<double> localMass;
461 chuckv 3450 int localPtArraySize;
462 chuckv 3392
463 chuckv 3459
464 chuckv 3450 std::fill(isSurfaceID.begin(),isSurfaceID.end(),false);
465 chuckv 3459
466 chuckv 3450
467     FORALLfacets {
468 chuckv 3392
469 chuckv 3450 if (!facet->simplicial){
470 chuckv 3392 // should never happen with Qt
471     sprintf(painCave.errMsg, "ConvexHull: non-simplicaial facet detected");
472 chuckv 3464 painCave.isFatal = 1;
473 chuckv 3392 simError();
474 chuckv 3450 }
475    
476    
477     vertices = qh_facet3vertex(facet);
478     FOREACHvertex_(vertices){
479 chuckv 3392 id = qh_pointid(vertex->point);
480 chuckv 3450
481     if( !isSurfaceID[id] ){
482     isSurfaceID[id] = true;
483 chuckv 3392 }
484 chuckv 3450 }
485     qh_settempfree(&vertices);
486 chuckv 3392
487 chuckv 3450 } //FORALLfacets
488 chuckv 3275
489 chuckv 3450
490 chuckv 3275
491    
492 chuckv 3450 int idx = 0;
493 chuckv 3459 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 chuckv 3464
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 chuckv 3535
506 chuckv 3473 RealType bdmass = bodydoubles[idx]->getMass();
507     localMass.push_back(bdmass);
508    
509 chuckv 3459 localPtsMap.push_back(idx);
510 chuckv 3473
511 chuckv 3459 }
512 chuckv 3450
513 chuckv 3459
514 chuckv 3473 localPtArraySize = int(localPts.size()/3.0);
515 chuckv 3459
516 chuckv 3450 MPI::COMM_WORLD.Allgather(&localPtArraySize,1,MPI::INT,&NstoProc_[0],1,MPI::INT);
517 chuckv 3473
518 chuckv 3459 Nsglobal_=0;
519 chuckv 3450 for (int i = 0; i < nproc_; i++){
520     Nsglobal_ += NstoProc_[i];
521 chuckv 3473 vecNstoProc_[i] = NstoProc_[i]*3;
522 chuckv 3450 }
523 chuckv 3459
524    
525 chuckv 3473 int nglobalPts = Nsglobal_*3;
526 chuckv 3459
527 chuckv 3450
528 chuckv 3473 std::vector<double> globalPts(nglobalPts);
529     std::vector<double> globalVel(nglobalPts);
530     std::vector<double> globalMass(Nsglobal_);
531 chuckv 3459
532 chuckv 3535
533    
534 chuckv 3459 isSurfaceID.resize(nglobalPts);
535    
536    
537 chuckv 3450 std::fill(globalPts.begin(),globalPts.end(),0.0);
538 chuckv 3459
539 chuckv 3473 vecdispls_[0] = 0;
540 chuckv 3450 /* Build a displacements array */
541     for (int i = 1; i < nproc_; i++){
542 chuckv 3473 vecdispls_[i] = vecdispls_[i-1] + vecNstoProc_[i-1];
543     }
544    
545     displs_[0] = 0;
546     for (int i = 1; i < nproc_; i++){
547 chuckv 3450 displs_[i] = displs_[i-1] + NstoProc_[i-1];
548     }
549 chuckv 3459
550 chuckv 3544 int noffset = displs_[myrank_];
551 chuckv 3450 /* gather the potential hull */
552 chuckv 3392
553 chuckv 3535 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 chuckv 3473 MPI::COMM_WORLD.Allgatherv(&localMass[0],localPtArraySize,MPI::DOUBLE,&globalMass[0],&NstoProc_[0],&displs_[0],MPI::DOUBLE);
556 chuckv 3535
557 chuckv 3459 /*
558 chuckv 3535 int tmpidx = 0;
559    
560 chuckv 3459 if (myrank_ == 0){
561 chuckv 3535 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;
564 chuckv 3459 }
565     }
566     */
567 chuckv 3535
568 chuckv 3450 // Free previous hull
569 chuckv 3392 qh_freeqhull(!qh_ALL);
570     qh_memfreeshort(&curlong, &totlong);
571     if (curlong || totlong)
572     std::cerr << "qhull internal warning (main): did not free %d bytes of long memory (%d pieces) "
573     << totlong << curlong << std::endl;
574 chuckv 3450
575 chuckv 3535 if (qh_new_qhull(dim_, Nsglobal_, &globalPts[0], ismalloc,
576 chuckv 3450 const_cast<char *>(options_.c_str()), NULL, stderr)){
577    
578     sprintf(painCave.errMsg, "ConvexHull: Qhull failed to compute global convex hull");
579 chuckv 3459 painCave.isFatal = 1;
580 chuckv 3450 simError();
581    
582     } //qh_new_qhull
583    
584     #endif
585    
586    
587    
588    
589    
590    
591     unsigned int nf = qh num_facets;
592 chuckv 3544
593 chuckv 3450 /* Build Surface SD list first */
594    
595     std::fill(isSurfaceID.begin(),isSurfaceID.end(),false);
596 chuckv 3544 int numvers = 0;
597 chuckv 3450 FORALLfacets {
598    
599     if (!facet->simplicial){
600     // should never happen with Qt
601     sprintf(painCave.errMsg, "ConvexHull: non-simplicaial facet detected");
602 chuckv 3459 painCave.isFatal = 1;
603 chuckv 3450 simError();
604     } //simplicical
605    
606 chuckv 3465 Triangle face;
607 chuckv 3450 Vector3d V3dNormal(facet->normal[0],facet->normal[1],facet->normal[2]);
608 chuckv 3465 face.setNormal(V3dNormal);
609 chuckv 3461
610    
611    
612 chuckv 3473 //RealType faceArea = 0.5*V3dNormal.length();
613     RealType faceArea = qh_facetarea(facet);
614 chuckv 3465 face.setArea(faceArea);
615 chuckv 3450
616    
617     vertices = qh_facet3vertex(facet);
618 chuckv 3461
619     coordT *center = qh_getcenter(vertices);
620     Vector3d V3dCentroid(center[0], center[1], center[2]);
621 chuckv 3465 face.setCentroid(V3dCentroid);
622 chuckv 3464 Vector3d faceVel = V3Zero;
623 chuckv 3473 Vector3d p[3];
624     RealType faceMass = 0.0;
625     int ver = 0;
626 chuckv 3450 FOREACHvertex_(vertices){
627     id = qh_pointid(vertex->point);
628 chuckv 3473 p[ver][0] = vertex->point[0];
629     p[ver][1] = vertex->point[1];
630     p[ver][2] = vertex->point[2];
631 chuckv 3459 int localindex = id;
632     #ifdef IS_MPI
633 chuckv 3544 Vector3d velVector(globalVel[dim_ * id],globalVel[dim_ * id + 1], globalVel[dim_ * id + 2]);
634 chuckv 3473
635 chuckv 3464 faceVel = faceVel + velVector;
636 chuckv 3473 faceMass = faceMass + globalMass[id];
637 chuckv 3544 if (id >= noffset && id < (noffset + localPtArraySize) ){
638    
639     localindex = localPtsMap[id-noffset];
640 chuckv 3464 #else
641     faceVel = faceVel + bodydoubles[localindex]->getVel();
642 chuckv 3473 faceMass = faceMass + bodydoubles[localindex]->getMass();
643 chuckv 3459 #endif
644 chuckv 3473 face.addVertexSD(bodydoubles[localindex]);
645 chuckv 3459 if( !isSurfaceID[id] ){
646     isSurfaceID[id] = true;
647     #ifdef IS_MPI
648    
649     #endif
650    
651     surfaceSDs_.push_back(bodydoubles[localindex]);
652 chuckv 3544 // std::cout <<"This ID is: " << bodydoubles[localindex]->getGlobalIndex() << std::endl;
653 chuckv 3459
654     } //IF isSurfaceID
655    
656     #ifdef IS_MPI
657    
658     }else{
659 chuckv 3473 face.addVertexSD(NULL);
660 chuckv 3459 }
661     #endif
662 chuckv 3544 numvers++;
663 chuckv 3473 ver++;
664 chuckv 3450 } //Foreachvertex
665 chuckv 3461 /*
666     if (!SETempty_(facet->coplanarset)){
667     FOREACHpoint_(facet->coplanarset){
668     id = qh_pointid(point);
669     surfaceSDs_.push_back(bodydoubles[id]);
670     }
671     }
672 chuckv 3464 */
673 chuckv 3473 face.addVertices(p[0],p[1],p[2]);
674     face.setFacetMass(faceMass);
675 chuckv 3465 face.setFacetVelocity(faceVel/3.0);
676 chuckv 3450 Triangles_.push_back(face);
677     qh_settempfree(&vertices);
678 chuckv 3464
679 chuckv 3450 } //FORALLfacets
680    
681 chuckv 3535 /*
682 chuckv 3461 std::cout << surfaceSDs_.size() << std::endl;
683     for (SD = surfaceSDs_.begin(); SD != surfaceSDs_.end(); ++SD){
684     Vector3d thisatom = (*SD)->getPos();
685     std::cout << "Au " << thisatom.x() << " " << thisatom.y() << " " << thisatom.z() << std::endl;
686     }
687     */
688 chuckv 3275
689 chuckv 3450
690 chuckv 3461 Ns_ = surfaceSDs_.size();
691 chuckv 3464 nTriangles_ = Triangles_.size();
692 chuckv 3461
693     qh_getarea(qh facet_list);
694     volume_ = qh totvol;
695     area_ = qh totarea;
696    
697    
698    
699     qh_freeqhull(!qh_ALL);
700     qh_memfreeshort(&curlong, &totlong);
701     if (curlong || totlong)
702     std::cerr << "qhull internal warning (main): did not free %d bytes of long memory (%d pieces) "
703     << totlong << curlong << std::endl;
704    
705    
706    
707 chuckv 3083 }
708    
709 chuckv 3450
710    
711     void ConvexHull::printHull(const std::string& geomFileName)
712 chuckv 3275 {
713 chuckv 3137
714 chuckv 3392 FILE *newGeomFile;
715    
716     //create new .md file based on old .md file
717     newGeomFile = fopen(geomFileName.c_str(), "w");
718     qh_findgood_all(qh facet_list);
719     for (int i = 0; i < qh_PRINTEND; i++)
720     qh_printfacets(newGeomFile, qh PRINTout[i], qh facet_list, NULL, !qh_ALL);
721    
722     fclose(newGeomFile);
723 chuckv 3083 }
724 chuckv 3392 #endif //QHULL
725     #endif //CGAL
726 chuckv 3083
727 chuckv 3275
728