]> SALOME platform Git repositories - modules/smesh.git/commitdiff
Salome HOME
Porting to OCCT 7.8.0 OCCT780
authorjfa <jfa@opencascade.com>
Mon, 15 Jan 2024 13:40:33 +0000 (13:40 +0000)
committermbs <martin.bernhard@opencascade.com>
Thu, 16 May 2024 07:51:45 +0000 (08:51 +0100)
15 files changed:
src/Controls/SMESH_Controls.cxx
src/DriverSTL/DriverSTL_R_SMDS_Mesh.cxx
src/SMESH/SMESH_MeshEditor.cxx
src/SMESH/SMESH_MesherHelper.cxx
src/SMESHDS/CMakeLists.txt
src/SMESHDS/SMESHDS_Mesh.hxx
src/SMESHGUI/SMESHGUI.cxx
src/SMESHUtils/SMESH_FreeBorders.cxx
src/SMESHUtils/SMESH_MeshAlgos.cxx
src/SMESHUtils/SMESH_Offset.cxx
src/SMESHUtils/SMESH_TypeDefs.hxx
src/StdMeshers/StdMeshers_Quadrangle_2D.cxx
src/Tools/blocFissure/gmu/fissureCoude.py
test/ex_MakePolyLine.py
test/netgen_runner_2D.py

index 57855777584e688899d0d94b68d9c5e5086ef901..55d5c19a809ff3f173828d7fefea0bf0c1d74667 100644 (file)
@@ -3445,7 +3445,7 @@ void CoplanarFaces::SetMesh( const SMDS_Mesh* theMesh )
       return;
 
     const double cosTol = Cos( myToler * M_PI / 180. );
-    NCollection_Map< SMESH_TLink, SMESH_TLink > checkedLinks;
+    NCollection_Map< SMESH_TLink, SMESH_TLinkHasher > checkedLinks;
 
     std::list< std::pair< const SMDS_MeshElement*, gp_Vec > > faceQueue;
     faceQueue.push_back( std::make_pair( face, myNorm ));
index c838ec8065d6e16f07917ed29016c9862558e32b..c2fd53597f643c596c4dbf1d257ed8c9ff2c2a2a 100644 (file)
@@ -23,6 +23,7 @@
 #include "DriverSTL_R_SMDS_Mesh.h"
 
 #include <Basics_Utils.hxx>
+#include <Basics_OCCTVersion.hxx>
 
 #include <gp_Pnt.hxx>
 #include <NCollection_DataMap.hxx>
@@ -41,8 +42,12 @@ namespace
     //function : HashCode
     //purpose  :
     //=======================================================================
+#if OCC_VERSION_LARGE < 0x07080000
     inline static Standard_Integer HashCode
     (const gp_Pnt& point,  Standard_Integer Upper)
+#else
+    size_t operator()(const gp_Pnt& point) const
+#endif
     {
       union
       {
@@ -52,14 +57,22 @@ namespace
 
       point.Coord( U.R[0], U.R[1], U.R[2] );
 
+#if OCC_VERSION_LARGE < 0x07080000
       return ::HashCode(U.I[0]/23+U.I[1]/19+U.I[2]/17+U.I[3]/13+U.I[4]/11+U.I[5]/7,Upper);
+#else
+      return static_cast<size_t>(U.I[0]/23+U.I[1]/19+U.I[2]/17+U.I[3]/13+U.I[4]/11+U.I[5]/7);
+#endif
     }
     //=======================================================================
     //function : IsEqual
     //purpose  :
     //=======================================================================
+#if OCC_VERSION_LARGE < 0x07080000
     inline static Standard_Boolean IsEqual
     (const gp_Pnt& point1, const gp_Pnt& point2)
+#else
+    bool operator()(const gp_Pnt& point1, const gp_Pnt& point2) const
+#endif
     {
       static Standard_Real tab1[3], tab2[3];
       point1.Coord(tab1[0],tab1[1],tab1[2]);
index 8c940322e0229a86bdf12fb802492d234ae6f260..125bbaa2e04c031916233c38deb0fe52ee32ab03 100644 (file)
 #include "SMESH_TryCatch.hxx" // include after OCCT headers!
 
 #include <smIdType.hxx>
+#include <Basics_OCCTVersion.hxx>
 
 #define cast2Node(elem) static_cast<const SMDS_MeshNode*>( elem )
 
@@ -1507,7 +1508,7 @@ int SMESH_MeshEditor::Reorient2D( TIDSortedElemSet &  theFaces,
   vector< const SMDS_MeshElement* > facesNearLink;
   vector< std::pair< int, int > >   nodeIndsOfFace;
   TIDSortedElemSet                  avoidSet, emptySet;
-  NCollection_Map< SMESH_TLink, SMESH_TLink > checkedLinks;
+  NCollection_Map< SMESH_TLink, SMESH_TLinkHasher > checkedLinks;
 
   while ( !theRefFaces.empty() )
   {
@@ -7906,6 +7907,8 @@ bool SMESH_MeshEditor::applyMerge( const SMDS_MeshElement* elem,
 // purpose : allow comparing elements basing on their nodes
 // ========================================================
 
+struct ComparableElementHasher;
+
 class ComparableElement : public boost::container::flat_set< smIdType >
 {
   typedef boost::container::flat_set< smIdType >  int_set;
@@ -7914,6 +7917,8 @@ class ComparableElement : public boost::container::flat_set< smIdType >
   smIdType                mySumID;
   mutable int             myGroupID;
 
+  friend ComparableElementHasher;
+
 public:
 
   ComparableElement( const SMDS_MeshElement* theElem ):
@@ -7942,7 +7947,11 @@ public:
     mySumID   = src.mySumID;
     myGroupID = src.myGroupID;
   }
+};
 
+struct ComparableElementHasher
+{
+#if OCC_VERSION_LARGE < 0x07080000
   static int HashCode(const ComparableElement& se, int limit )
   {
     return ::HashCode( FromSmIdType<int>(se.mySumID), limit );
@@ -7951,7 +7960,17 @@ public:
   {
     return ( se1 == se2 );
   }
+#else
+  size_t operator()(const ComparableElement& se) const
+  {
+    return static_cast<size_t>(FromSmIdType<int>(se.mySumID));
+  }
 
+  bool operator()(const ComparableElement& se1, const ComparableElement& se2) const
+  {
+    return ( se1 == se2 );
+  }
+#endif
 };
 
 //=======================================================================
@@ -7969,8 +7988,8 @@ void SMESH_MeshEditor::FindEqualElements( TIDSortedElemSet &        theElements,
   if ( theElements.empty() ) elemIt = GetMeshDS()->elementsIterator();
   else                       elemIt = SMESHUtils::elemSetIterator( theElements );
 
-  typedef NCollection_Map< ComparableElement, ComparableElement > TMapOfElements;
-  typedef std::list<smIdType>                                     TGroupOfElems;
+  typedef NCollection_Map< ComparableElement, ComparableElementHasher > TMapOfElements;
+  typedef std::list<smIdType>                                           TGroupOfElems;
   TMapOfElements               mapOfElements;
   std::vector< TGroupOfElems > arrayOfGroups;
   TGroupOfElems                groupOfElems;
index 1be97fb72763ae21a0d01f134bd2799f431425ca..49a1152810d7d6f874aa8fc8fb69ad24583b0c41 100644 (file)
@@ -3279,6 +3279,16 @@ double SMESH_MesherHelper::getFaceMaxTol( const TopoDS_Shape& face ) const
   return tol;
 }
 
+bool CheckAlmostZero(gp_Vec & vec1,gp_Vec & vec2, gp_Vec & vecref)
+{
+  auto v1   = gp_Dir(vec1);
+  auto v2   = gp_Dir(vec2);
+  auto vref = gp_Dir(vecref);
+  auto XYZ = v1.Crossed (v2);
+  double cond  = XYZ.X()*vref.X()+XYZ.Y()*vref.Y()+XYZ.Z()*vref.Z();
+  return (Abs(cond) <= 1e-12);
+}
+
 //================================================================================
 /*!
  * \brief Return an angle between two EDGEs sharing a common VERTEX with reference
@@ -3341,7 +3351,9 @@ double SMESH_MesherHelper::GetAngle( const TopoDS_Edge &   theE1,
     if ( theE2.Orientation() /*GetSubShapeOri( F, theE2 )*/ == TopAbs_REVERSED )
       vec2.Reverse();
     angle = vec1.AngleWithRef( vec2, vecRef );
-
+    if ( angle < 0. && CheckAlmostZero(vec1,vec2,vecRef))
+      angle*=-1;
+      
     if ( Abs ( angle ) >= 0.99 * M_PI )
     {
       BRep_Tool::Range( theE1, f, l );
index 4d058b84da8cb3aacb52f215b6ae89bc24e3892e..56e87ef7be3fec908a6908fcf152b5d4046704f8 100644 (file)
@@ -67,7 +67,6 @@ SET(SMESHDS_HEADERS
   SMESHDS_GroupOnGeom.hxx
   SMESHDS_GroupOnFilter.hxx
   SMESH_SMESHDS.hxx
-  SMESHDS_DataMapOfShape.hxx
   SMESH_Controls.hxx
 )
 
index c90ba3e25bec7b0bef4fb88f0235d2535f6bd645..1b8d2e2fd6d18a31c4d3cebc451969b2952d19df 100644 (file)
@@ -32,6 +32,8 @@
 #include "SMDS_Mesh.hxx"
 #include "SMESHDS_SubMesh.hxx"
 
+#include <Basics_OCCTVersion.hxx>
+
 #include <TopTools_IndexedMapOfShape.hxx>
 #include <TopoDS_Shape.hxx>
 
@@ -57,8 +59,34 @@ class SMDS_BallElement;
  * So this functionality implement on new NCollection_DataMap technology
  */
 #include <NCollection_DataMap.hxx>
-#include "SMESHDS_DataMapOfShape.hxx"
 typedef std::list<const SMESHDS_Hypothesis*>                          THypList;
+
+struct SMESHDS_Hasher
+{
+#if OCC_VERSION_LARGE < 0x07080000
+  static inline Standard_Boolean IsEqual(const TopoDS_Shape& S1,
+                                         const TopoDS_Shape& S2)
+  {
+    return S1.IsSame(S2);
+  }
+  static inline Standard_Integer HashCode(const TopoDS_Shape& S,
+                                          const Standard_Integer Upper)
+  {
+    return ::HashCode( S, Upper);
+  }
+#else
+  bool operator()(const TopoDS_Shape& S1, const TopoDS_Shape& S2) const
+  {
+    // for the purpose of ShapeToHypothesis map we don't consider shapes orientation
+    return S1.IsSame(S2);
+  }
+  size_t operator()(const TopoDS_Shape& S) const
+  {
+    return std::hash<TopoDS_Shape>{}(S);
+  }
+#endif
+};
+
 typedef NCollection_DataMap< TopoDS_Shape, THypList, SMESHDS_Hasher > ShapeToHypothesis;
 
 class SMESHDS_GroupBase;
index 5f1d4549a0704076c63ebd7755a371544b2e5cdc..90f0b230d7466fbb7ef4bc954b8220d61bf5737a 100644 (file)
 //  File   : SMESHGUI.cxx
 //  Author : Nicolas REJNERI, Open CASCADE S.A.S.
 
+#include <Basics_OCCTVersion.hxx>
+
+#if OCC_VERSION_LARGE < 0x07080000
+
 #include <Standard_math.hxx>  // E.A. must be included before Python.h to fix compilation on windows
+
+#else
+
+#ifdef _MSC_VER
+#ifndef _USE_MATH_DEFINES
+#define _USE_MATH_DEFINES
+#endif
+#include <math.h>
+#endif
+
+#endif
+
 #ifdef HAVE_FINITE
 #undef HAVE_FINITE            // VSR: avoid compilation warning on Linux : "HAVE_FINITE" redefined
 #endif
index ce7fb230f0c7c42d6d418dca969b9ae421e27b70..f1575ca25122142e0ee7684b4e3829a56d1d6472 100644 (file)
@@ -444,7 +444,7 @@ void SMESH_MeshAlgos::FindCoincidentFreeBorders(SMDS_Mesh&              mesh,
                                                 CoincidentFreeBorders & foundFreeBordes)
 {
   // find free links
-  typedef NCollection_DataMap<SMESH_TLink, const SMDS_MeshElement*, SMESH_TLink > TLink2FaceMap;
+  typedef NCollection_DataMap<SMESH_TLink, const SMDS_MeshElement*, SMESH_TLinkHasher > TLink2FaceMap;
   TLink2FaceMap linkMap;
   int nbSharedLinks = 0;
   SMDS_FaceIteratorPtr faceIt = mesh.facesIterator();
@@ -840,7 +840,7 @@ void SMESH_MeshAlgos::FindFreeBorders(SMDS_Mesh&       theMesh,
   bool isManifold = true;
 
   // find free links
-  typedef NCollection_DataMap<SMESH_TLink, const SMDS_MeshElement*, SMESH_TLink > TLink2FaceMap;
+  typedef NCollection_DataMap<SMESH_TLink, const SMDS_MeshElement*, SMESH_TLinkHasher > TLink2FaceMap;
   TLink2FaceMap linkMap;
   int nbSharedLinks = 0;
   SMDS_FaceIteratorPtr faceIt = theMesh.facesIterator();
index 24958cc81fd77b243c71cd3142d8caf096d5b2d7..fac231ec8507711c65bef041672a8074deb78cd8 100644 (file)
@@ -1989,8 +1989,8 @@ SMESH_MeshAlgos::FindSharpEdges( SMDS_Mesh* theMesh,
   std::vector< Edge > resultEdges;
   if ( !theMesh ) return resultEdges;
 
-  typedef std::pair< bool, const SMDS_MeshNode* >                            TIsSharpAndMedium;
-  typedef NCollection_DataMap< SMESH_TLink, TIsSharpAndMedium, SMESH_TLink > TLinkSharpMap;
+  typedef std::pair< bool, const SMDS_MeshNode* >                                  TIsSharpAndMedium;
+  typedef NCollection_DataMap< SMESH_TLink, TIsSharpAndMedium, SMESH_TLinkHasher > TLinkSharpMap;
 
   TLinkSharpMap linkIsSharp;
   Standard_Integer nbBuckets = FromSmIdType<Standard_Integer>( theMesh->NbFaces() );
@@ -2101,8 +2101,8 @@ SMESH_MeshAlgos::SeparateFacesByEdges( SMDS_Mesh* theMesh, const std::vector< Ed
 
   // build map of face edges (SMESH_TLink) and their faces
 
-  typedef std::vector< const SMDS_MeshElement* >                    TFaceVec;
-  typedef NCollection_DataMap< SMESH_TLink, TFaceVec, SMESH_TLink > TFacesByLinks;
+  typedef std::vector< const SMDS_MeshElement* >                          TFaceVec;
+  typedef NCollection_DataMap< SMESH_TLink, TFaceVec, SMESH_TLinkHasher > TFacesByLinks;
   TFacesByLinks facesByLink;
   Standard_Integer nbBuckets = FromSmIdType<Standard_Integer>( theMesh->NbFaces() );
   if ( nbBuckets > 0 )
@@ -2270,7 +2270,7 @@ std::vector< const SMDS_MeshNode*> SMESH_MeshAlgos::GetCommonNodes(const SMDS_Me
 bool SMESH_MeshAlgos::IsOn2DBoundary( const SMDS_MeshNode*                 theNode,
                                       std::vector< const SMDS_MeshNode*> * theNeibors )
 {
-  typedef NCollection_DataMap< SMESH_TLink, int, SMESH_TLink > TLinkCountMap;
+  typedef NCollection_DataMap< SMESH_TLink, int, SMESH_TLinkHasher > TLinkCountMap;
   TLinkCountMap linkCountMap( 10 );
 
   int nbFreeLinks = 0;
index 754411be86c868dd7184270f1434eb566c1a7e77..3b0ba7a211ac93198e4a0c6a09d226930690cfee 100644 (file)
@@ -29,6 +29,7 @@
 #include "SMDS_Mesh.hxx"
 
 #include <Utils_SALOME_Exception.hxx>
+#include <Basics_OCCTVersion.hxx>
 
 #include <Bnd_B3d.hxx>
 #include <NCollection_Map.hxx>
@@ -43,7 +44,7 @@ namespace
   const int theMaxNbFaces = 256; // max number of faces sharing a node
 
   typedef NCollection_DataMap< const SMDS_MeshNode*, const SMDS_MeshNode*, SMESH_Hasher > TNNMap;
-  typedef NCollection_Map< SMESH_Link, SMESH_Link >                                       TLinkMap;
+  typedef NCollection_Map< SMESH_Link, SMESH_TLinkHasher >                                TLinkMap;
 
   //--------------------------------------------------------------------------------
   /*!
@@ -75,7 +76,11 @@ namespace
     const SMDS_MeshNode* IntNode() const { return myIntNode.Node(); }
     const SMDS_MeshNode* Node1() const { return myNode[ myReverse ]; }
     const SMDS_MeshNode* Node2() const { return myNode[ !myReverse ]; }
+  };
 
+  struct CutLinkHasher
+  {
+#if OCC_VERSION_LARGE < 0x07080000
     static Standard_Integer HashCode(const CutLink&         link,
                                      const Standard_Integer upper)
     {
@@ -90,9 +95,23 @@ namespace
                link1.myNode[1] == link2.myNode[1] &&
                link1.myIndex == link2.myIndex );
     }
+#else
+    size_t operator()(const CutLink& link) const
+    {
+      return size_t( link.myNode[0]->GetID() +
+                     link.myNode[1]->GetID() +
+                     link.myIndex );
+    }
+    bool operator()(const CutLink& link1, const CutLink& link2 ) const
+    {
+      return ( link1.myNode[0] == link2.myNode[0] &&
+               link1.myNode[1] == link2.myNode[1] &&
+               link1.myIndex == link2.myIndex );
+    }
+#endif
   };
 
-  typedef NCollection_Map< CutLink, CutLink > TCutLinkMap;
+  typedef NCollection_Map< CutLink, CutLinkHasher > TCutLinkMap;
 
   //--------------------------------------------------------------------------------
   /*!
@@ -274,7 +293,16 @@ namespace
                       TLinkMap&                    theCutOffCoplanarLinks) const;
     void InitLinks() const;
     bool IsCoplanar( const EdgePart* edge ) const;
+    void Dump() const;
 
+  private:
+
+    EdgePart* getTwin( const EdgePart* edge ) const;
+  };
+
+  struct CutFaceHasher
+  {
+#if OCC_VERSION_LARGE < 0x07080000
     static Standard_Integer HashCode(const CutFace& f, const Standard_Integer upper)
     {
       return ::HashCode( FromSmIdType<int>(f.myInitFace->GetID()), upper );
@@ -283,14 +311,20 @@ namespace
     {
       return f1.myInitFace == f2.myInitFace;
     }
-    void Dump() const;
-
-  private:
+#else
+    size_t operator()(const CutFace& f) const
+    {
+      return FromSmIdType<int>(f.myInitFace->GetID());
+    }
 
-    EdgePart* getTwin( const EdgePart* edge ) const;
+    bool operator()(const CutFace& f1, const CutFace& f2) const
+    {
+      return f1.myInitFace == f2.myInitFace;
+    }
+#endif
   };
 
-  typedef NCollection_Map< CutFace, CutFace > TCutFaceMap;
+  typedef NCollection_Map< CutFace, CutFaceHasher > TCutFaceMap;
 
   //--------------------------------------------------------------------------------
   /*!
index 2e920132d4aecac0e10e544f6e3340ae1f1f101a..9ae5a4ce25bbbc30fccf3dd29138f5a7c00b06f3 100644 (file)
@@ -27,6 +27,8 @@
 #ifndef __SMESH_TypeDefs_HXX__
 #define __SMESH_TypeDefs_HXX__
 
+#include <Basics_OCCTVersion.hxx>
+
 #include "SMESH_Utils.hxx"
 
 #include "SMDS_SetIterator.hxx"
@@ -36,6 +38,7 @@
 
 #include <gp_XYZ.hxx>
 #include <gp_XY.hxx>
+#include <NCollection_Sequence.hxx>
 
 #include <map>
 #include <list>
@@ -156,6 +159,19 @@ struct SMESH_TLink: public NLink
   const SMDS_MeshNode* node2() const { return second; }
 
   // methods for usage of SMESH_TLink as a hasher in NCollection maps
+  //static int HashCode(const SMESH_TLink& link, int aLimit)
+  //{
+  //  return smIdHasher::HashCode( link.node1()->GetID() + link.node2()->GetID(), aLimit );
+  //}
+  //static Standard_Boolean IsEqual(const SMESH_TLink& l1, const SMESH_TLink& l2)
+  //{
+  //  return ( l1.node1() == l2.node1() && l1.node2() == l2.node2() );
+  //}
+};
+// a hasher in NCollection maps
+struct SMESH_TLinkHasher
+{
+#if OCC_VERSION_LARGE < 0x07080000
   static int HashCode(const SMESH_TLink& link, int aLimit)
   {
     return smIdHasher::HashCode( link.node1()->GetID() + link.node2()->GetID(), aLimit );
@@ -164,6 +180,16 @@ struct SMESH_TLink: public NLink
   {
     return ( l1.node1() == l2.node1() && l1.node2() == l2.node2() );
   }
+#else
+  size_t operator()(const SMESH_TLink& link) const
+  {
+    return smIdHasher()( link.node1()->GetID() + link.node2()->GetID() );
+  }
+  bool operator()(const SMESH_TLink& l1, const SMESH_TLink& l2) const
+  {
+    return ( l1.node1() == l2.node1() && l1.node2() == l2.node2() );
+  }
+#endif
 };
 typedef SMESH_TLink SMESH_Link;
 
@@ -218,6 +244,7 @@ typedef SMESH_TNodeXYZ SMESH_NodeXYZ;
 
 struct SMESH_Hasher
 {
+#if OCC_VERSION_LARGE < 0x07080000
   static Standard_Integer HashCode(const SMDS_MeshElement* e, const Standard_Integer upper)
   {
     return smIdHasher::HashCode( e->GetID(), upper );
@@ -226,6 +253,17 @@ struct SMESH_Hasher
   {
     return ( e1 == e2 );
   }
+#else
+  size_t operator()(const SMDS_MeshElement* e) const
+  {
+    return smIdHasher()( e->GetID() );
+  }
+
+  bool operator()(const SMDS_MeshElement* e1, const SMDS_MeshElement* e2) const
+  {
+    return ( e1 == e2 );
+  }
+#endif
 };
 
 //--------------------------------------------------
@@ -262,10 +300,8 @@ typedef std::vector< const SMDS_MeshElement* > SMESH_SequenceOfElemPtr;
 
 // --------------------------------------------------------------------------------
 // class SMESH_SequenceOfNode
-#include <NCollection_DefineSequence.hxx>
-typedef const SMDS_MeshNode* SMDS_MeshNodePtr;
 
-DEFINE_SEQUENCE(SMESH_SequenceOfNode,
-                SMESH_BaseCollectionNodePtr, SMDS_MeshNodePtr)
+typedef const SMDS_MeshNode*                     SMDS_MeshNodePtr;
+typedef NCollection_Sequence< SMDS_MeshNodePtr > SMESH_SequenceOfNode;
 
 #endif
index e05279c3598008534aabff0fae37bece57051c9f..1abd265956b30f3ebef365877d8ec18edeaddab6 100644 (file)
@@ -49,7 +49,7 @@
 #include <Bnd_Box.hxx>
 #include <GeomAPI_ProjectPointOnSurf.hxx>
 #include <Geom_Surface.hxx>
-#include <NCollection_DefineArray2.hxx>
+#include <NCollection_Array2.hxx>
 #include <Precision.hxx>
 #include <ShapeAnalysis.hxx>
 #include <TColStd_SequenceOfInteger.hxx>
index f8ce927f03ec780398840ae7eadd854fcc5470da..014ab56fe5b0ee8b93a257d59fbe27f880e2ea54 100644 (file)
@@ -106,7 +106,10 @@ class fissureCoude(fissureGenerique):
 
     # --- tube coude sain
 
-    geometrieSaine = geompy.MakePartition([tube_1, coude, tube_2, P1, P2], [plan_y], [], [], geompy.ShapeType["SOLID"], 0, [], 1)
+    #Refactor to perform the partition in two steps to avoid break from occt (reported bug#)
+    #geometrieSaine = geompy.MakePartition([tube_1, coude, tube_2, P1, P2], [plan_y], [], [], geompy.ShapeType["SOLID"], 0, [], 1)
+    geometrieSaine0 = geompy.MakePartition([tube_1, coude], [plan_y], [], [], geompy.ShapeType["SOLID"], 0, [], 1)
+    geometrieSaine = geompy.MakePartition([geometrieSaine0, tube_2, P1, P2], [plan_y], [], [], geompy.ShapeType["SOLID"], 0, [], 1)
     geomPublish(initLog.debug, geometrieSaine, self.nomCas, self.numeroCas )
     [P1, P2] = geompy.RestoreGivenSubShapes(geometrieSaine, [P1, P2], GEOM.FSM_GetInPlaceByHistory, False, True)
 
index a7fb6a7e7dcf9f6c6ed55ee3146da0414154c036..46115fccf6108121240ebd2e10ec87418f798143 100644 (file)
@@ -25,12 +25,21 @@ if not isDone:
 # define arguments for MakePolyLine
 
 segments = []
+
 # between nodes 20 and 1, default plane
 segments.append( SMESH.PolySegment( 20, 0, SMESH.PointStruct(-1, -1, -1), 1, 0, SMESH.PointStruct(-1, -1, -1), smesh.MakeDirStruct(0,0,0) ))
+
 # between nodes 1 and 100, default plane
 segments.append( SMESH.PolySegment( 1, 0, SMESH.PointStruct(-1, -1, -1), 200, 0, SMESH.PointStruct(-1, -1, -1), smesh.MakeDirStruct(0,0,0) ))
-# between nodes 200 and edge (578, 577), plane includes vector (1,1,1)
-segments.append( SMESH.PolySegment( 200, 0, SMESH.PointStruct(-1, -1, -1), 578, 577, SMESH.PointStruct(-1, -1, -1), smesh.MakeDirStruct(1,1,1) ))
+
+# between node 200 and edge (578, 577), plane includes vector (1,1,1)
+#segments.append( SMESH.PolySegment( 200, 0, SMESH.PointStruct(-1, -1, -1), 578, 577, SMESH.PointStruct(-1, -1, -1), smesh.MakeDirStruct(1,1,1) ))
+# nodes 578 and 577 are not always neighbour, so, use another approach
+
+# between node 200 and an edge, close to point (200, 90, 70), plane includes vector (1,1,1)
+elems = Mesh_1.FindElementsByPoint(200, 90, 70, SMESH.FACE)
+nodes = Mesh_1.GetElemNodes(elems[0])
+segments.append( SMESH.PolySegment( 200, 0, SMESH.PointStruct(-1, -1, -1), nodes[0], nodes[1], SMESH.PointStruct(-1, -1, -1), smesh.MakeDirStruct(1,1,1) ))
 
 Mesh_1.MakePolyLine( segments, "1D group")
 
index 36be1b7170b22e7bfc39aeb05a4eea01ee087000..255a4a61a71208e4abbc1ed9ac4653d887696edc 100644 (file)
@@ -169,8 +169,8 @@ def test_netgen2dLenghtFromEdge():
     print("Nb Segments:", nb_segments)
     print("Nb Points:", nb_points)
 
-    assert nb_triangles == 12
-    assert nb_points == 8
+    assert nb_triangles == 12*2
+    assert nb_points == 14
     assert nb_segments == 12
 
 if __name__ == "__main__":