// Copyright (C) 2007-2014 CEA/DEN, EDF R&D, OPEN CASCADE // // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN, // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS // // This library is free software; you can redistribute it and/or // modify it under the terms of the GNU Lesser General Public // License as published by the Free Software Foundation; either // version 2.1 of the License, or (at your option) any later version. // // This library is distributed in the hope that it will be useful, // but WITHOUT ANY WARRANTY; without even the implied warranty of // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU // Lesser General Public License for more details. // // You should have received a copy of the GNU Lesser General Public // License along with this library; if not, write to the Free Software // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA // // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com // #include #include #include #include #include #include // OCCT Includes #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include // CAREFUL ! position of this file is critic : see Lucien PIGNOLONI / OCC #define STD_SORT_ALGO 1 namespace { /** * This function constructs and returns modified shape from the original one * for singular cases. It is used for the method GetMinDistanceSingular. * * \param theShape the original shape * \param theModifiedShape output parameter. The modified shape. * \param theAddDist output parameter. The added distance for modified shape. * \retval true if the shape is modified; false otherwise. * * \internal */ Standard_Boolean ModifyShape(const TopoDS_Shape &theShape, TopoDS_Shape &theModifiedShape, Standard_Real &theAddDist) { Standard_Boolean isModified = Standard_False; TopExp_Explorer anExp; int nbf = 0; theAddDist = 0.; theModifiedShape.Nullify(); for ( anExp.Init( theShape, TopAbs_FACE ); anExp.More(); anExp.Next() ) { nbf++; theModifiedShape = anExp.Current(); } if(nbf==1) { TopoDS_Shape sh = theShape; while(sh.ShapeType()==TopAbs_COMPOUND) { TopoDS_Iterator it(sh); sh = it.Value(); } Handle(Geom_Surface) S = BRep_Tool::Surface(TopoDS::Face(theModifiedShape)); if( S->IsKind(STANDARD_TYPE(Geom_SphericalSurface)) || S->IsKind(STANDARD_TYPE(Geom_ToroidalSurface)) || S->IsUPeriodic()) { const Standard_Boolean isShell = (sh.ShapeType()==TopAbs_SHELL || sh.ShapeType()==TopAbs_FACE); if( isShell || S->IsUPeriodic() ) { // non solid case or any periodic surface (Mantis 22454). double U1,U2,V1,V2; // changes for 0020677: EDF 1219 GEOM: MinDistance gives 0 instead of 20.88 //S->Bounds(U1,U2,V1,V2); changed by ShapeAnalysis::GetFaceUVBounds(TopoDS::Face(theModifiedShape),U1,U2,V1,V2); // end of changes for 020677 (dmv) Handle(Geom_RectangularTrimmedSurface) TrS1 = new Geom_RectangularTrimmedSurface(S,U1,(U1+U2)/2.,V1,V2); Handle(Geom_RectangularTrimmedSurface) TrS2 = new Geom_RectangularTrimmedSurface(S,(U1+U2)/2.,U2,V1,V2); BRep_Builder B; TopoDS_Face F1,F2; TopoDS_Shape aMShape; if (isShell) { B.MakeCompound(TopoDS::Compound(aMShape)); } else { B.MakeShell(TopoDS::Shell(aMShape)); } B.MakeFace(F1,TrS1,1.e-7); B.Add(aMShape,F1); B.MakeFace(F2,TrS2,1.e-7); B.Add(aMShape,F2); Handle(ShapeFix_Shape) sfs = new ShapeFix_Shape; if (!isShell) { // The original shape is a solid. TopoDS_Solid aSolid; B.MakeSolid(aSolid); B.Add(aSolid, aMShape); aMShape = aSolid; } sfs->Init(aMShape); sfs->SetPrecision(1.e-6); sfs->SetMaxTolerance(1.0); sfs->Perform(); theModifiedShape = sfs->Shape(); isModified = Standard_True; } else { if( S->IsKind(STANDARD_TYPE(Geom_SphericalSurface)) ) { Handle(Geom_SphericalSurface) SS = Handle(Geom_SphericalSurface)::DownCast(S); gp_Pnt PC = SS->Location(); BRep_Builder B; TopoDS_Vertex V; B.MakeVertex(V,PC,1.e-7); theModifiedShape = V; theAddDist = SS->Radius(); isModified = Standard_True; } else { Handle(Geom_ToroidalSurface) TS = Handle(Geom_ToroidalSurface)::DownCast(S); gp_Ax3 ax3 = TS->Position(); Handle(Geom_Circle) C = new Geom_Circle(ax3.Ax2(),TS->MajorRadius()); BRep_Builder B; TopoDS_Edge E; B.MakeEdge(E,C,1.e-7); theModifiedShape = E; theAddDist = TS->MinorRadius(); isModified = Standard_True; } } } else { theModifiedShape = theShape; } } else theModifiedShape = theShape; return isModified; } //======================================================================= //function : ShapeToDouble //purpose : used by CompareShapes::operator() //======================================================================= std::pair ShapeToDouble (const TopoDS_Shape& S, bool isOldSorting) { // Computing of CentreOfMass gp_Pnt GPoint; double Len; if (S.ShapeType() == TopAbs_VERTEX) { GPoint = BRep_Tool::Pnt(TopoDS::Vertex(S)); Len = (double)S.Orientation(); } else { GProp_GProps GPr; // BEGIN: fix for Mantis issue 0020842 if (isOldSorting) { BRepGProp::LinearProperties(S, GPr); } else { if (S.ShapeType() == TopAbs_EDGE || S.ShapeType() == TopAbs_WIRE) { BRepGProp::LinearProperties(S, GPr); } else if (S.ShapeType() == TopAbs_FACE || S.ShapeType() == TopAbs_SHELL) { BRepGProp::SurfaceProperties(S, GPr); } else { BRepGProp::VolumeProperties(S, GPr); } } // END: fix for Mantis issue 0020842 GPoint = GPr.CentreOfMass(); Len = GPr.Mass(); } double dMidXYZ = GPoint.X() * 999.0 + GPoint.Y() * 99.0 + GPoint.Z() * 0.9; return std::make_pair(dMidXYZ, Len); } void parseWard( const GEOMUtils::LevelsList &theLevelList, std::string &treeStr ) { treeStr.append( "{" ); for( GEOMUtils::LevelsList::const_iterator j = theLevelList.begin(); j != theLevelList.end(); ++j ) { if ( j != theLevelList.begin() ) { treeStr.append( ";" ); } GEOMUtils::LevelInfo level = (*j); GEOMUtils::LevelInfo::iterator upIter; for ( upIter = level.begin(); upIter != level.end(); ++upIter ) { if ( upIter != level.begin() ) { treeStr.append( "," ); } treeStr.append( upIter->first ); for ( std::vector::iterator k = upIter->second.begin(); k != upIter->second.end(); ++k ) { treeStr.append( "_" ); treeStr.append( *k ); } } } treeStr.append( "}" ); } GEOMUtils::LevelsList parseWard( const std::string& theData, std::size_t& theCursor ) { std::size_t indexStart = theData.find( "{", theCursor ) + 1; std::size_t indexEnd = theData.find( "}", indexStart ); std::string ward = theData.substr( indexStart, indexEnd - indexStart ); std::stringstream ss(ward); std::string substr; std::vector levelsListStr; while ( std::getline( ss, substr, ';' ) ) { if ( !substr.empty() ) levelsListStr.push_back( substr ); } GEOMUtils::LevelsList levelsListData; for( int level = 0; level < levelsListStr.size(); level++ ) { std::vector namesListStr; std::stringstream ss1( levelsListStr[level] ); while ( std::getline( ss1, substr, ',' ) ) { if ( !substr.empty() ) namesListStr.push_back( substr ); } GEOMUtils::LevelInfo levelInfoData; for( int node = 0; node < namesListStr.size(); node++ ) { std::vector linksListStr; std::stringstream ss2( namesListStr[node] ); while ( std::getline( ss2, substr, '_' ) ) { if ( !substr.empty() ) linksListStr.push_back( substr ); } std::string nodeItem = linksListStr[0]; if( !nodeItem.empty() ) { GEOMUtils::NodeLinks linksListData; for( int link = 1; link < linksListStr.size(); link++ ) { std::string linkItem = linksListStr[link]; linksListData.push_back( linkItem ); }// Links levelInfoData[nodeItem] = linksListData; } }// Level's objects levelsListData.push_back(levelInfoData); }// Levels theCursor = indexEnd + 1; return levelsListData; } } //======================================================================= //function : GetPosition //purpose : //======================================================================= gp_Ax3 GEOMUtils::GetPosition (const TopoDS_Shape& theShape) { gp_Ax3 aResult; if (theShape.IsNull()) return aResult; // Axes aResult.Transform(theShape.Location().Transformation()); if (theShape.ShapeType() == TopAbs_FACE) { Handle(Geom_Surface) aGS = BRep_Tool::Surface(TopoDS::Face(theShape)); if (!aGS.IsNull() && aGS->IsKind(STANDARD_TYPE(Geom_Plane))) { Handle(Geom_Plane) aGPlane = Handle(Geom_Plane)::DownCast(aGS); gp_Pln aPln = aGPlane->Pln(); aResult = aPln.Position(); // In case of reverse orinetation of the face invert the plane normal // (the face's normal does not mathc the plane's normal in this case) if(theShape.Orientation() == TopAbs_REVERSED) { gp_Dir Vx = aResult.XDirection(); gp_Dir N = aResult.Direction().Mirrored(Vx); gp_Pnt P = aResult.Location(); aResult = gp_Ax3(P, N, Vx); } } } // Origin gp_Pnt aPnt; TopAbs_ShapeEnum aShType = theShape.ShapeType(); if (aShType == TopAbs_VERTEX) { aPnt = BRep_Tool::Pnt(TopoDS::Vertex(theShape)); } else { if (aShType == TopAbs_COMPOUND) { aShType = GetTypeOfSimplePart(theShape); } GProp_GProps aSystem; if (aShType == TopAbs_EDGE || aShType == TopAbs_WIRE) BRepGProp::LinearProperties(theShape, aSystem); else if (aShType == TopAbs_FACE || aShType == TopAbs_SHELL) BRepGProp::SurfaceProperties(theShape, aSystem); else BRepGProp::VolumeProperties(theShape, aSystem); aPnt = aSystem.CentreOfMass(); } aResult.SetLocation(aPnt); return aResult; } //======================================================================= //function : GetVector //purpose : //======================================================================= gp_Vec GEOMUtils::GetVector (const TopoDS_Shape& theShape, Standard_Boolean doConsiderOrientation) { if (theShape.IsNull()) Standard_NullObject::Raise("Null shape is given for a vector"); if (theShape.ShapeType() != TopAbs_EDGE) Standard_TypeMismatch::Raise("Invalid shape is given, must be a vector or an edge"); TopoDS_Edge anE = TopoDS::Edge(theShape); TopoDS_Vertex V1, V2; TopExp::Vertices(anE, V1, V2, doConsiderOrientation); if (V1.IsNull() || V2.IsNull()) Standard_NullObject::Raise("Invalid edge is given, it must have two points"); gp_Vec aV (BRep_Tool::Pnt(V1), BRep_Tool::Pnt(V2)); if (aV.Magnitude() < gp::Resolution()) { Standard_ConstructionError::Raise("Vector of zero length is given"); } return aV; } //======================================================================= //function : CompareShapes::operator() //purpose : used by std::sort(), called from SortShapes() //======================================================================= bool GEOMUtils::CompareShapes::operator() (const TopoDS_Shape& theShape1, const TopoDS_Shape& theShape2) { if (!myMap.IsBound(theShape1)) { myMap.Bind(theShape1, ShapeToDouble(theShape1, myIsOldSorting)); } if (!myMap.IsBound(theShape2)) { myMap.Bind(theShape2, ShapeToDouble(theShape2, myIsOldSorting)); } std::pair val1 = myMap.Find(theShape1); std::pair val2 = myMap.Find(theShape2); double tol = Precision::Confusion(); bool exchange = Standard_False; double dMidXYZ = val1.first - val2.first; if (dMidXYZ >= tol) { exchange = Standard_True; } else if (Abs(dMidXYZ) < tol) { double dLength = val1.second - val2.second; if (dLength >= tol) { exchange = Standard_True; } else if (Abs(dLength) < tol && theShape1.ShapeType() <= TopAbs_FACE) { // PAL17233 // equal values possible on shapes such as two halves of a sphere and // a membrane inside the sphere Bnd_Box box1,box2; BRepBndLib::Add(theShape1, box1); if (!box1.IsVoid()) { BRepBndLib::Add(theShape2, box2); Standard_Real dSquareExtent = box1.SquareExtent() - box2.SquareExtent(); if (dSquareExtent >= tol) { exchange = Standard_True; } else if (Abs(dSquareExtent) < tol) { Standard_Real aXmin, aYmin, aZmin, aXmax, aYmax, aZmax, val1, val2; box1.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax); val1 = (aXmin+aXmax)*999.0 + (aYmin+aYmax)*99.0 + (aZmin+aZmax)*0.9; box2.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax); val2 = (aXmin+aXmax)*999.0 + (aYmin+aYmax)*99.0 + (aZmin+aZmax)*0.9; if ((val1 - val2) >= tol) { exchange = Standard_True; } } } } } //return val1 < val2; return !exchange; } //======================================================================= //function : SortShapes //purpose : //======================================================================= void GEOMUtils::SortShapes (TopTools_ListOfShape& SL, const Standard_Boolean isOldSorting) { #ifdef STD_SORT_ALGO std::vector aShapesVec; aShapesVec.reserve(SL.Extent()); TopTools_ListIteratorOfListOfShape it (SL); for (; it.More(); it.Next()) { aShapesVec.push_back(it.Value()); } SL.Clear(); CompareShapes shComp (isOldSorting); std::stable_sort(aShapesVec.begin(), aShapesVec.end(), shComp); //std::sort(aShapesVec.begin(), aShapesVec.end(), shComp); std::vector::const_iterator anIter = aShapesVec.begin(); for (; anIter != aShapesVec.end(); ++anIter) { SL.Append(*anIter); } #else // old implementation Standard_Integer MaxShapes = SL.Extent(); TopTools_Array1OfShape aShapes (1,MaxShapes); TColStd_Array1OfInteger OrderInd(1,MaxShapes); TColStd_Array1OfReal MidXYZ (1,MaxShapes); //X,Y,Z; TColStd_Array1OfReal Length (1,MaxShapes); //X,Y,Z; // Computing of CentreOfMass Standard_Integer Index; GProp_GProps GPr; gp_Pnt GPoint; TopTools_ListIteratorOfListOfShape it(SL); for (Index=1; it.More(); Index++) { TopoDS_Shape S = it.Value(); SL.Remove( it ); // == it.Next() aShapes(Index) = S; OrderInd.SetValue (Index, Index); if (S.ShapeType() == TopAbs_VERTEX) { GPoint = BRep_Tool::Pnt( TopoDS::Vertex( S )); Length.SetValue( Index, (Standard_Real) S.Orientation()); } else { // BEGIN: fix for Mantis issue 0020842 if (isOldSorting) { BRepGProp::LinearProperties (S, GPr); } else { if (S.ShapeType() == TopAbs_EDGE || S.ShapeType() == TopAbs_WIRE) { BRepGProp::LinearProperties (S, GPr); } else if (S.ShapeType() == TopAbs_FACE || S.ShapeType() == TopAbs_SHELL) { BRepGProp::SurfaceProperties(S, GPr); } else { BRepGProp::VolumeProperties(S, GPr); } } // END: fix for Mantis issue 0020842 GPoint = GPr.CentreOfMass(); Length.SetValue(Index, GPr.Mass()); } MidXYZ.SetValue(Index, GPoint.X()*999.0 + GPoint.Y()*99.0 + GPoint.Z()*0.9); //cout << Index << " L: " << Length(Index) << "CG: " << MidXYZ(Index) << endl; } // Sorting Standard_Integer aTemp; Standard_Boolean exchange, Sort = Standard_True; Standard_Real tol = Precision::Confusion(); while (Sort) { Sort = Standard_False; for (Index=1; Index < MaxShapes; Index++) { exchange = Standard_False; Standard_Real dMidXYZ = MidXYZ(OrderInd(Index)) - MidXYZ(OrderInd(Index+1)); Standard_Real dLength = Length(OrderInd(Index)) - Length(OrderInd(Index+1)); if ( dMidXYZ >= tol ) { // cout << "MidXYZ: " << MidXYZ(OrderInd(Index))<< " > " <= tol ) { // cout << "Length: " << Length(OrderInd(Index))<< " > " <= tol ) { // cout << "SquareExtent: " << box1.SquareExtent()<<" > "< val2; if ((val1 - val2) >= tol) { exchange = Standard_True; } //cout << "box: " << val1<<" > "< 1) { Standard_ConstructionError::Raise("Multiple edges near the given point are found"); } else if (nbFound == 0) { Standard_ConstructionError::Raise("There are no edges near the given point"); } else { } return aResult; } //======================================================================= //function : PreciseBoundingBox //purpose : //======================================================================= Standard_Boolean GEOMUtils::PreciseBoundingBox (const TopoDS_Shape &theShape, Bnd_Box &theBox) { if ( theBox.IsVoid() ) BRepBndLib::Add( theShape, theBox ); Standard_Real aBound[6]; theBox.Get(aBound[0], aBound[2], aBound[4], aBound[1], aBound[3], aBound[5]); Standard_Integer i; const gp_Pnt aMid(0.5*(aBound[1] + aBound[0]), // XMid 0.5*(aBound[3] + aBound[2]), // YMid 0.5*(aBound[5] + aBound[4])); // ZMid const gp_XYZ aSize(aBound[1] - aBound[0], // DX aBound[3] - aBound[2], // DY aBound[5] - aBound[4]); // DZ const gp_Pnt aPnt[6] = { gp_Pnt(aBound[0] - (aBound[1] - aBound[0]), aMid.Y(), aMid.Z()), // XMin gp_Pnt(aBound[1] + (aBound[1] - aBound[0]), aMid.Y(), aMid.Z()), // XMax gp_Pnt(aMid.X(), aBound[2] - (aBound[3] - aBound[2]), aMid.Z()), // YMin gp_Pnt(aMid.X(), aBound[3] + (aBound[3] - aBound[2]), aMid.Z()), // YMax gp_Pnt(aMid.X(), aMid.Y(), aBound[4] - (aBound[5] - aBound[4])), // ZMin gp_Pnt(aMid.X(), aMid.Y(), aBound[5] + (aBound[5] - aBound[4])) // ZMax }; const gp_Dir aDir[3] = { gp::DX(), gp::DY(), gp::DZ() }; const Standard_Real aPlnSize[3] = { 0.5*Max(aSize.Y(), aSize.Z()), // XMin, XMax planes 0.5*Max(aSize.X(), aSize.Z()), // YMin, YMax planes 0.5*Max(aSize.X(), aSize.Y()) // ZMin, ZMax planes }; gp_Pnt aPMin[2]; for (i = 0; i < 6; i++) { const Standard_Integer iHalf = i/2; const gp_Pln aPln(aPnt[i], aDir[iHalf]); BRepBuilderAPI_MakeFace aMkFace(aPln, -aPlnSize[iHalf], aPlnSize[iHalf], -aPlnSize[iHalf], aPlnSize[iHalf]); if (!aMkFace.IsDone()) { return Standard_False; } TopoDS_Shape aFace = aMkFace.Shape(); // Get minimal distance between planar face and shape. Standard_Real aMinDist = GEOMUtils::GetMinDistance(aFace, theShape, aPMin[0], aPMin[1]); if (aMinDist < 0.) { return Standard_False; } aBound[i] = aPMin[1].Coord(iHalf + 1); } // Update Bounding box with the new values. theBox.SetVoid(); theBox.Update(aBound[0], aBound[2], aBound[4], aBound[1], aBound[3], aBound[5]); return Standard_True; } //======================================================================= //function : GetMinDistanceSingular //purpose : //======================================================================= double GEOMUtils::GetMinDistanceSingular(const TopoDS_Shape& aSh1, const TopoDS_Shape& aSh2, gp_Pnt& Ptmp1, gp_Pnt& Ptmp2) { TopoDS_Shape tmpSh1; TopoDS_Shape tmpSh2; Standard_Real AddDist1 = 0.; Standard_Real AddDist2 = 0.; Standard_Boolean IsChange1 = ModifyShape(aSh1, tmpSh1, AddDist1); Standard_Boolean IsChange2 = ModifyShape(aSh2, tmpSh2, AddDist2); if( !IsChange1 && !IsChange2 ) return -2.0; BRepExtrema_DistShapeShape dst(tmpSh1,tmpSh2); if (dst.IsDone()) { double MinDist = 1.e9; gp_Pnt PMin1, PMin2, P1, P2; for (int i = 1; i <= dst.NbSolution(); i++) { P1 = dst.PointOnShape1(i); P2 = dst.PointOnShape2(i); Standard_Real Dist = P1.Distance(P2); if (MinDist > Dist) { MinDist = Dist; PMin1 = P1; PMin2 = P2; } } if(MinDist<1.e-7) { Ptmp1 = PMin1; Ptmp2 = PMin2; } else { gp_Dir aDir(gp_Vec(PMin1,PMin2)); if( MinDist > (AddDist1+AddDist2) ) { Ptmp1 = gp_Pnt( PMin1.X() + aDir.X()*AddDist1, PMin1.Y() + aDir.Y()*AddDist1, PMin1.Z() + aDir.Z()*AddDist1 ); Ptmp2 = gp_Pnt( PMin2.X() - aDir.X()*AddDist2, PMin2.Y() - aDir.Y()*AddDist2, PMin2.Z() - aDir.Z()*AddDist2 ); return (MinDist - AddDist1 - AddDist2); } else { if( AddDist1 > 0 ) { Ptmp1 = gp_Pnt( PMin1.X() + aDir.X()*AddDist1, PMin1.Y() + aDir.Y()*AddDist1, PMin1.Z() + aDir.Z()*AddDist1 ); Ptmp2 = Ptmp1; } else { Ptmp2 = gp_Pnt( PMin2.X() - aDir.X()*AddDist2, PMin2.Y() - aDir.Y()*AddDist2, PMin2.Z() - aDir.Z()*AddDist2 ); Ptmp1 = Ptmp2; } } } double res = MinDist - AddDist1 - AddDist2; if(res<0.) res = 0.0; return res; } return -2.0; } //======================================================================= //function : GetMinDistance //purpose : //======================================================================= Standard_Real GEOMUtils::GetMinDistance (const TopoDS_Shape& theShape1, const TopoDS_Shape& theShape2, gp_Pnt& thePnt1, gp_Pnt& thePnt2) { Standard_Real aResult = 1.e9; // Issue 0020231: A min distance bug with torus and vertex. // Make GetMinDistance() return zero if a sole VERTEX is inside any of SOLIDs // which of shapes consists of only one vertex? TopExp_Explorer exp1(theShape1,TopAbs_VERTEX), exp2(theShape2,TopAbs_VERTEX); TopoDS_Shape V1 = exp1.More() ? exp1.Current() : TopoDS_Shape(); TopoDS_Shape V2 = exp2.More() ? exp2.Current() : TopoDS_Shape(); exp1.Next(); exp2.Next(); if ( exp1.More() ) V1.Nullify(); if ( exp2.More() ) V2.Nullify(); // vertex and container of solids TopoDS_Shape V = V1.IsNull() ? V2 : V1; TopoDS_Shape S = V1.IsNull() ? theShape1 : theShape2; if ( !V.IsNull() ) { // classify vertex against solids gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( V ) ); for ( exp1.Init( S, TopAbs_SOLID ); exp1.More(); exp1.Next() ) { BRepClass3d_SolidClassifier classifier( exp1.Current(), p, 1e-6); if ( classifier.State() == TopAbs_IN ) { thePnt1 = p; thePnt2 = p; return 0.0; } } } // End Issue 0020231 // skl 30.06.2008 // additional workaround for bugs 19899, 19908 and 19910 from Mantis double dist = GEOMUtils::GetMinDistanceSingular (theShape1, theShape2, thePnt1, thePnt2); if (dist > -1.0) { return dist; } BRepExtrema_DistShapeShape dst (theShape1, theShape2); if (dst.IsDone()) { gp_Pnt P1, P2; for (int i = 1; i <= dst.NbSolution(); i++) { P1 = dst.PointOnShape1(i); P2 = dst.PointOnShape2(i); Standard_Real Dist = P1.Distance(P2); if (aResult > Dist) { aResult = Dist; thePnt1 = P1; thePnt2 = P2; } } } return aResult; } //======================================================================= // function : ConvertClickToPoint() // purpose : Returns the point clicked in 3D view //======================================================================= gp_Pnt GEOMUtils::ConvertClickToPoint( int x, int y, Handle(V3d_View) aView ) { V3d_Coordinate XEye, YEye, ZEye, XAt, YAt, ZAt; aView->Eye( XEye, YEye, ZEye ); aView->At( XAt, YAt, ZAt ); gp_Pnt EyePoint( XEye, YEye, ZEye ); gp_Pnt AtPoint( XAt, YAt, ZAt ); gp_Vec EyeVector( EyePoint, AtPoint ); gp_Dir EyeDir( EyeVector ); gp_Pln PlaneOfTheView = gp_Pln( AtPoint, EyeDir ); Standard_Real X, Y, Z; aView->Convert( x, y, X, Y, Z ); gp_Pnt ConvertedPoint( X, Y, Z ); gp_Pnt2d ConvertedPointOnPlane = ProjLib::Project( PlaneOfTheView, ConvertedPoint ); gp_Pnt ResultPoint = ElSLib::Value( ConvertedPointOnPlane.X(), ConvertedPointOnPlane.Y(), PlaneOfTheView ); return ResultPoint; } //======================================================================= // function : ConvertTreeToString() // purpose : Returns the string representation of dependency tree //======================================================================= void GEOMUtils::ConvertTreeToString( const TreeModel &tree, std::string &treeStr ) { TreeModel::const_iterator i; for ( i = tree.begin(); i != tree.end(); ++i ) { treeStr.append( i->first ); treeStr.append( "-" ); std::vector upLevelList = i->second.first; treeStr.append( "upward" ); parseWard( upLevelList, treeStr ); std::vector downLevelList = i->second.second; treeStr.append( "downward" ); parseWard( downLevelList, treeStr ); } } //======================================================================= // function : ConvertStringToTree() // purpose : Returns the dependency tree //======================================================================= void GEOMUtils::ConvertStringToTree( const std::string &theData, TreeModel &tree ) { std::size_t cursor = 0; while( theData.find('-',cursor) != std::string::npos ) //find next selected object { std::size_t objectIndex = theData.find( '-', cursor ); std::string objectEntry = theData.substr( cursor, objectIndex - cursor ); cursor = objectIndex; std::size_t upwardIndexBegin = theData.find("{",cursor) + 1; std::size_t upwardIndexFinish = theData.find("}",upwardIndexBegin); LevelsList upwardList = parseWard( theData, cursor ); LevelsList downwardList = parseWard( theData, cursor ); tree[objectEntry] = std::pair( upwardList, downwardList ); } }