//********************************************************************************** //test of polyhedron intersection callable from python shell bool do_Polyhedras_Intersect(const shared_ptr<Shape>& cm1,const shared_ptr<Shape>& cm2,const State& state1,const State& state2){ const Se3r& se31=state1.se3; const Se3r& se32=state2.se3; Polyhedra* A = static_cast<Polyhedra*>(cm1.get()); Polyhedra* B = static_cast<Polyhedra*>(cm2.get()); //move and rotate 1st the CGAL structure Polyhedron Matrix3r rot_mat = (se31.orientation).toRotationMatrix(); Vector3r trans_vec = se31.position; Transformation t_rot_trans(rot_mat(0,0),rot_mat(0,1),rot_mat(0,2), trans_vec[0],rot_mat(1,0),rot_mat(1,1),rot_mat(1,2),trans_vec[1],rot_mat(2,0),rot_mat(2,1),rot_mat(2,2),trans_vec[2],1.); Polyhedron PA = A->GetPolyhedron(); std::transform( PA.points_begin(), PA.points_end(), PA.points_begin(), t_rot_trans); //move and rotate 2st the CGAL structure Polyhedron rot_mat = (se32.orientation).toRotationMatrix(); trans_vec = se32.position; t_rot_trans = Transformation(rot_mat(0,0),rot_mat(0,1),rot_mat(0,2), trans_vec[0],rot_mat(1,0),rot_mat(1,1),rot_mat(1,2),trans_vec[1],rot_mat(2,0),rot_mat(2,1),rot_mat(2,2),trans_vec[2],1.); Polyhedron PB = B->GetPolyhedron(); std::transform( PB.points_begin(), PB.points_end(), PB.points_begin(), t_rot_trans); //calculate plane equations std::transform( PA.facets_begin(), PA.facets_end(), PA.planes_begin(),Plane_equation()); std::transform( PB.facets_begin(), PB.facets_end(), PB.planes_begin(),Plane_equation()); //call test return do_intersect(PA,PB); }
void mergeCoplanar(Polyhedron& p,bool step2) { int facetsBefore = p.size_of_facets(); p.normalize_border(); if(!p.size_of_border_halfedges()) { // Calculate normals only once in advace! so all tris should be coplanar with the original std::transform( p.facets_begin(), p.facets_end(),p.planes_begin(),Plane_equation()); // Calculate plane equations (only works on tri<- = bs)=true bool coplanarFound = true; std::vector<Polyhedron::Halfedge_handle> skipHEs; while (coplanarFound) { coplanarFound = false; // Set coplanarFound false int percCount = 1; for (Polyhedron::Halfedge_iterator hit = p.halfedges_begin(); hit != p.halfedges_end(); ++hit,++percCount){ // Loop through all halfedges if (is_coplanar(hit,true)){ // If coplanar and equals semantics Polyhedron::Halfedge_handle removeMe = hit; while (CGAL::circulator_size(removeMe->vertex_begin()) < 3) // Mover handle to beginning of linestring removeMe = removeMe->next(); bool jcnh = false; if (!step2) jcnh = joinCreatesNoHole (hit); else jcnh = joinCreatesNoHole2(hit); if (jcnh){ // If no holes will be created std::cout << "\rFacets before/after: "<<facetsBefore<<" -> "<< p.size_of_facets()<<". ("<<100*percCount/p.size_of_halfedges()<<"%)"; while (CGAL::circulator_size(removeMe->opposite()->vertex_begin()) < 3) // Join vertexes until at the other end of linestring if (removeMe->facet_degree()>3 && removeMe->opposite()->facet_degree()>3) removeMe = (p.join_vertex(removeMe))->next()->opposite(); else // One of the faces turned into a triangle ->remove center vertex break; if (CGAL::circulator_size(removeMe->opposite()->vertex_begin()) < 3) // Remove remained of the border p.erase_center_vertex(removeMe->opposite()); // if two segments remain remove center point else p.join_facet(removeMe); // if one segment remains join facets coplanarFound = true; break; } else { // simplify border, but how to do this safely? not optimal solution implemented. Should do: add inward offseted point of intersection etc. if (std::find(skipHEs.begin(), skipHEs.end(),hit)!=skipHEs.end()) { // Skip if hit in skipList while (CGAL::circulator_size(removeMe->opposite()->vertex_begin()) < 3) { // Join vertexes until at the other end of linestring if (removeMe->facet_degree()>3 && removeMe->opposite()->facet_degree()>3) if (triDoesNotIntersectFacet(removeMe)) // if tri reMe,reME->prev does not intersect left or right facet removeMe = (p.join_vertex(removeMe))->next()->opposite(); // remove removeME else { skipHEs.push_back(removeMe); skipHEs.push_back(removeMe->opposite()); removeMe = removeMe->prev(); // move removeME one halfedge back } else break; // stop if only a triangle remains or at other end } skipHEs.push_back(removeMe); skipHEs.push_back(removeMe->opposite()); } } } } } } //if (!step2) mergeCoplanar(p,true); //else std::cout << "\rFacets before/after: "<<facetsBefore<<" -> "<< p.size_of_facets()<<". (100%)"<<std::endl; }
int main() { Point_3 p( 1, 0, 0); Point_3 q( 0, 1, 0); Point_3 r( 0, 0, 1); Point_3 s( 0, 0, 0); Polyhedron P; P.make_tetrahedron( p, q, r, s); std::transform( P.facets_begin(), P.facets_end(), P.planes_begin(), Plane_equation()); CGAL::set_pretty_mode( std::cout); std::copy( P.planes_begin(), P.planes_end(), std::ostream_iterator<Plane_3>( std::cout, "\n")); return 0; }
int main(int argc, char* argv[]) { CGAL_assertion(argc==2); std::ifstream in1(argv[1]); Polyhedron P1; in1 >> P1; std::transform( P1.facets_begin(), P1.facets_end(), P1.planes_begin(), Plane_equation()); CGAL_assertion(is_strongly_convex_3(P1)); Gausian_map G1(P1); G1.visualize(); }
void Polyhedra::Initialize(){ if (init) return; bool isRandom = false; //get vertices int N = (int) v.size(); if (N==0) { //generate randomly while ((int) v.size()<4) GenerateRandomGeometry(); N = (int) v.size(); isRandom = true; } //compute convex hull of vertices std::vector<CGALpoint> points; points.resize(v.size()); for(int i=0;i<N;i++) { points[i] = CGALpoint(v[i][0],v[i][1],v[i][2]); } CGAL::convex_hull_3(points.begin(), points.end(), P); //connect triagular facets if possible std::transform(P.facets_begin(), P.facets_end(), P.planes_begin(),Plane_equation()); P = Simplify(P, 1E-9); //modify order of v according to CGAl polyhedron int i = 0; v.clear(); for (Polyhedron::Vertex_iterator vIter = P.vertices_begin(); vIter != P.vertices_end(); ++vIter, i++){ v.push_back(Vector3r(vIter->point().x(),vIter->point().y(),vIter->point().z())); } //list surface triangles for plotting faceTri.clear(); std::transform(P.facets_begin(), P.facets_end(), P.planes_begin(),Plane_equation()); for (Polyhedron::Facet_iterator fIter = P.facets_begin(); fIter != P.facets_end(); fIter++){ Polyhedron::Halfedge_around_facet_circulator hfc0; int n = fIter->facet_degree(); hfc0 = fIter->facet_begin(); int a = std::distance(P.vertices_begin(), hfc0->vertex()); for (int i=2; i<n; i++){ ++hfc0; faceTri.push_back(a); faceTri.push_back(std::distance(P.vertices_begin(), hfc0->vertex())); faceTri.push_back(std::distance(P.vertices_begin(), hfc0->next()->vertex())); } } //compute centroid and volume P_volume_centroid(P, &volume, ¢roid); //check vierd behavior of CGAL in tessalation if(isRandom && volume*1.75<4./3.*3.14*size[0]/2.*size[1]/2.*size[2]/2.) { v.clear(); seed = rand(); Initialize(); } Vector3r translation((-1)*centroid); //set centroid to be [0,0,0] for(int i=0;i<N;i++) { v[i] = v[i]-centroid; } if(isRandom) centroid = Vector3r::Zero(); Vector3r origin(0,0,0); //move and rotate also the CGAL structure Polyhedron Transformation t_trans(1.,0.,0.,translation[0],0.,1.,0.,translation[1],0.,0.,1.,translation[2],1.); std::transform( P.points_begin(), P.points_end(), P.points_begin(), t_trans); //compute inertia Real vtet; Vector3r ctet; Matrix3r Itet1, Itet2; Matrix3r inertia_tensor(Matrix3r::Zero()); for(int i=0; i<(int) faceTri.size(); i+=3){ vtet = std::abs((origin-v[faceTri[i+2]]).dot((v[faceTri[i]]-v[faceTri[i+2]]).cross(v[faceTri[i+1]]-v[faceTri[i+2]]))/6.); ctet = (origin+v[faceTri[i]]+v[faceTri[i+1]]+v[faceTri[i+2]]) / 4.; Itet1 = TetraInertiaTensor(origin-ctet, v[faceTri[i]]-ctet, v[faceTri[i+1]]-ctet, v[faceTri[i+2]]-ctet); ctet = ctet-origin; Itet2<< ctet[1]*ctet[1]+ctet[2]*ctet[2], -ctet[0]*ctet[1], -ctet[0]*ctet[2], -ctet[0]*ctet[1], ctet[0]*ctet[0]+ctet[2]*ctet[2], -ctet[2]*ctet[1], -ctet[0]*ctet[2], -ctet[2]*ctet[1], ctet[1]*ctet[1]+ctet[0]*ctet[0]; inertia_tensor = inertia_tensor + Itet1 + Itet2*vtet; } if(std::abs(inertia_tensor(0,1))+std::abs(inertia_tensor(0,2))+std::abs(inertia_tensor(1,2)) < 1E-13){ // no need to rotate, inertia already diagonal orientation = Quaternionr::Identity(); inertia = Vector3r(inertia_tensor(0,0),inertia_tensor(1,1),inertia_tensor(2,2)); }else{ // calculate eigenvectors of I Vector3r rot; Matrix3r I_rot(Matrix3r::Zero()), I_new(Matrix3r::Zero()); matrixEigenDecomposition(inertia_tensor,I_rot,I_new); // I_rot = eigenvectors of inertia_tensors in columns // I_new = eigenvalues on diagonal // set positove direction of vectors - otherwise reloading does not work Matrix3r sign(Matrix3r::Zero()); Real max_v_signed = I_rot(0,0); Real max_v = std::abs(I_rot(0,0)); if (max_v < std::abs(I_rot(1,0))) {max_v_signed = I_rot(1,0); max_v = std::abs(I_rot(1,0));} if (max_v < std::abs(I_rot(2,0))) {max_v_signed = I_rot(2,0); max_v = std::abs(I_rot(2,0));} sign(0,0) = max_v_signed/max_v; max_v_signed = I_rot(0,1); max_v = std::abs(I_rot(0,1)); if (max_v < std::abs(I_rot(1,1))) {max_v_signed = I_rot(1,1); max_v = std::abs(I_rot(1,1));} if (max_v < std::abs(I_rot(2,1))) {max_v_signed = I_rot(2,1); max_v = std::abs(I_rot(2,1));} sign(1,1) = max_v_signed/max_v; sign(2,2) = 1.; I_rot = I_rot*sign; // force the eigenvectors to be right-hand oriented Vector3r third = (I_rot.col(0)).cross(I_rot.col(1)); I_rot(0,2) = third[0]; I_rot(1,2) = third[1]; I_rot(2,2) = third[2]; inertia = Vector3r(I_new(0,0),I_new(1,1),I_new(2,2)); orientation = Quaternionr(I_rot); //rotate the voronoi cell so that x - is maximal inertia axis and z - is minimal inertia axis //orientation.normalize(); //not needed for(int i=0; i< (int) v.size();i++) { v[i] = orientation.conjugate()*v[i]; } //rotate also the CGAL structure Polyhedron Matrix3r rot_mat = (orientation.conjugate()).toRotationMatrix(); Transformation t_rot(rot_mat(0,0),rot_mat(0,1),rot_mat(0,2),rot_mat(1,0),rot_mat(1,1),rot_mat(1,2),rot_mat(2,0),rot_mat(2,1),rot_mat(2,2),1.); std::transform( P.points_begin(), P.points_end(), P.points_begin(), t_rot); } //initialization done init = 1; }
void set_semantic_AABB_C2V(Polyhedron& exteriorPolyhe,PolVector& polyVec) { if (exteriorPolyhe.is_pure_triangle()) { std::transform( exteriorPolyhe.facets_begin(), exteriorPolyhe.facets_end(),exteriorPolyhe.planes_begin(),Plane_equation()); std::vector<std::string> semList; std::vector<std::shared_ptr<AAbbTree>> treeList; // Build Trees. One for each semantic for(PolVector::iterator pvIt = polyVec.begin();pvIt!=polyVec.end();++pvIt) {// Get AABB trees of all semantics if (pvIt->is_pure_triangle()) { std::string semP = pvIt->facets_begin()->semanticBLA; std::vector<std::string>::iterator strIt = std::find(semList.begin(), semList.end(),semP); if (strIt==semList.end()) { // If new sematic semList.push_back(semP); // Add sem std::shared_ptr<AAbbTree> tree = std::make_shared<AAbbTree>(pvIt->facets_begin(),pvIt->facets_end()); // Create tree tree->accelerate_distance_queries(); // accelerate treeList.push_back(tree); // Add tree } else // If not new treeList[strIt-semList.begin()]->insert(pvIt->facets_begin(),pvIt->facets_end()); // Append to tree } else std::cerr << "ERROR: Not pure triangle (set_semantic_AABB2C2V)" << std::endl; } // For each facet calculate the least distance to each tree std::string semListStr = boost::algorithm::join((semList), " "); int percCount = 1; Polyhedron::Facet_iterator exfIt; // Iterate over exterior faces for (exfIt = exteriorPolyhe.facets_begin(); exfIt != exteriorPolyhe.facets_end(); ++exfIt,++percCount) { std::cout << "\r"<<semListStr<<". ("<<100*percCount/exteriorPolyhe.size_of_facets()<<"%)"; Vector_3 orthVec = exfIt->plane().orthogonal_vector(); normalizeVector(orthVec); //if (!normalizeVector(ortVec)) continue; std::vector<distSemFace> dsfList(semList.size()); Point_3 centerPoint = comp_facetCentroid(exfIt); // Compute centroid std::vector<Kernel::FT> leastSemDistances; for (int intIt=0;intIt<(int)treeList.size();++intIt) { // Loop over all trees AAbbTree::Point_and_primitive_id pp = treeList[intIt]->closest_point_and_primitive(centerPoint); dsfList[intIt].dist = CGAL::squared_distance(centerPoint,pp.first); // Store distance semantic and facet for each tree dsfList[intIt].sem = semList[intIt]; dsfList[intIt].fh = pp.second; } std::sort(dsfList.begin(),dsfList.end(),by_dist()); exfIt->leastSqDistance = dsfList[0].dist; // least sqrt distance if (exfIt->isMinkFacet = dsfList[0].dist > SEMANTIC_DISTANCE_THRESHOLD) { exfIt->semanticBLA = TO_DIST_SEMANTIC; // Default semantic if too distant continue; } else exfIt->semanticBLA = dsfList[0].sem; // Semantics of closest Vector_3 faceNormal; Kernel::FT faceSqArea; double minAngle = 10; Kernel::FT maxArea= 0; for (std::vector<distSemFace>::iterator slIt = dsfList.begin();slIt != dsfList.end();++slIt)// HANDLE ANYTHING AS LESS IMPORTANT if (slIt->dist < dsfList[0].dist+OVERLAP_DIST_THRESHOLD) { // Check if Equidistant pointVector facetPoints = comp_facetPoints(exfIt); CGAL::normal_vector_newell_3(facetPoints.begin(),facetPoints.end(),faceNormal); // Calculate normal vector, ortVec set to zero in newell double angle = comp_angle(orthVec,faceNormal); if (angle!=-1 && angle < minAngle+OVERLAP_ANGLE_THRESHOLD) { if (minAngle >= angle+OVERLAP_ANGLE_THRESHOLD) exfIt->equidistSems.clear(); if (angle < minAngle) minAngle = angle; faceSqArea = comp_facetSquaredArea(facetPoints); if (faceSqArea>maxArea-OVERLAP_AREA_THRESHOLD) { if (maxArea<=faceSqArea-OVERLAP_AREA_THRESHOLD) exfIt->equidistSems.clear(); if (faceSqArea>maxArea) maxArea = faceSqArea; exfIt->equidistSems.push_back(slIt->sem); // Add equidist semantics } } } } std::cout << "\r"<<semListStr<<". (100%)" << std::endl; }else std::cerr << "ERROR: Not pure triangle (set_semantic_AABB2C2V)" << std::endl; }
//********************************************************************************** //generate "packing" of non-overlapping balls vector<Vector3r> fillBoxByBalls_cpp(Vector3r minCoord, Vector3r maxCoord, Vector3r sizemin, Vector3r sizemax, Vector3r ratio, int seed, shared_ptr<Material> mat, int NumPoints){ vector<Vector3r> v; Polyhedra trialP; Polyhedron trial, trial_moved; srand(seed); int it = 0; vector<Polyhedron> polyhedrons; vector<vector<Vector3r> > vv; Vector3r position; bool intersection; int count = 0; Vector3r radii; bool fixed_ratio = 0; if (ratio[0] > 0 && ratio[1] > 0 && ratio[2]>0){ fixed_ratio = 1; sizemax[0] = min(min(sizemax[0]/ratio[0], sizemax[1]/ratio[1]), sizemax[2]/ratio[2]); sizemin[0] = max(max(sizemin[0]/ratio[0], sizemin[1]/ratio[1]), sizemin[2]/ratio[2]); } fixed_ratio = 1; //force spherical //it - number of trials to make packing possibly more/less dense Vector3r random_size; while (it<1000){ it = it+1; if (it == 1){ if (fixed_ratio) { double rrr = (rand()*(sizemax[0]-sizemin[0])/RAND_MAX + sizemin[0])/2.; radii = Vector3r(rrr,rrr,rrr); }else { radii = Vector3r(rand()*(sizemax[0]-sizemin[0])/2.,rand()*(sizemax[1]-sizemin[1])/2.,rand()*(sizemax[2]-sizemin[2])/2.)/RAND_MAX + sizemin/2.; } trialP.v = BallPoints(radii,NumPoints,rand()); trialP.Initialize(); trial = trialP.GetPolyhedron(); Matrix3r rot_mat = (trialP.GetOri()).toRotationMatrix(); Transformation t_rot(rot_mat(0,0),rot_mat(0,1),rot_mat(0,2),rot_mat(1,0),rot_mat(1,1),rot_mat(1,2),rot_mat(2,0),rot_mat(2,1),rot_mat(2,2),1.); std::transform( trial.points_begin(), trial.points_end(), trial.points_begin(), t_rot); } position = Vector3r(rand()*(maxCoord[0]-minCoord[0]),rand()*(maxCoord[1]-minCoord[1]),rand()*(maxCoord[2]-minCoord[2]))/RAND_MAX + minCoord; //move CGAL structure Polyhedron Transformation transl(CGAL::TRANSLATION, ToCGALVector(position)); trial_moved = trial; std::transform( trial_moved.points_begin(), trial_moved.points_end(), trial_moved.points_begin(), transl); //calculate plane equations std::transform( trial_moved.facets_begin(), trial_moved.facets_end(), trial_moved.planes_begin(),Plane_equation()); intersection = false; //call test with boundary for(Polyhedron::Vertex_iterator vi = trial_moved.vertices_begin(); (vi != trial_moved.vertices_end()) && (!intersection); vi++){ intersection = (vi->point().x()<minCoord[0]) || (vi->point().x()>maxCoord[0]) || (vi->point().y()<minCoord[1]) || (vi->point().y()>maxCoord[1]) || (vi->point().z()<minCoord[2]) || (vi->point().z()>maxCoord[2]); } //call test with other polyhedrons for(vector<Polyhedron>::iterator a = polyhedrons.begin(); (a != polyhedrons.end()) && (!intersection); a++){ intersection = do_intersect(*a,trial_moved); if (intersection) break; } if (!intersection){ polyhedrons.push_back(trial_moved); v.clear(); for(Polyhedron::Vertex_iterator vi = trial_moved.vertices_begin(); vi != trial_moved.vertices_end(); vi++){ v.push_back(FromCGALPoint(vi->point())); } vv.push_back(v); it = 0; count ++; } } cout << "generated " << count << " polyhedrons"<< endl; //can't be used - no information about material Scene* scene=Omega::instance().getScene().get(); for(vector<vector<Vector3r> >::iterator p=vv.begin(); p!=vv.end(); ++p){ shared_ptr<Body> BP = NewPolyhedra(*p, mat); BP->shape->color = Vector3r(double(rand())/RAND_MAX,double(rand())/RAND_MAX,double(rand())/RAND_MAX); scene->bodies->insert(BP); } return v; }