]> SALOME platform Git repositories - modules/smesh.git/commitdiff
Salome HOME
23395: EDF 13855 - Crash SALOME when creating a mesh
authoreap <eap@opencascade.com>
Wed, 23 Nov 2016 18:04:46 +0000 (21:04 +0300)
committereap <eap@opencascade.com>
Wed, 23 Nov 2016 18:04:46 +0000 (21:04 +0300)
Pb with a prism whose base face includes 3 wires and only one wire in internal

src/SMESH/SMESH_Mesh.cxx
src/SMESH/SMESH_MesherHelper.cxx
src/StdMeshers/StdMeshers_FaceSide.cxx
src/StdMeshers/StdMeshers_Prism_3D.cxx
src/StdMeshers/StdMeshers_Projection_2D.cxx
src/StdMeshers_I/StdMeshers_ProjectionSource1D_i.cxx
src/StdMeshers_I/StdMeshers_ProjectionSource2D_i.cxx
src/StdMeshers_I/StdMeshers_ProjectionSource3D_i.cxx
src/StdMeshers_I/StdMeshers_QuadrangleParams_i.cxx

index 76d2dbc05d7938820cfa847df0bb9a808f489c49..ff127dece72588a29e297dde56a32e072bd7efb2 100644 (file)
@@ -2368,7 +2368,7 @@ bool SMESH_Mesh::SortByMeshOrder(std::vector<SMESH_subMesh*>& theListToSort) con
   vector<SMESH_subMesh*>::iterator onlyBIt = onlyOrderedList.begin();
   vector<SMESH_subMesh*>::iterator onlyEIt = onlyOrderedList.end();
 
-  // iterate on ordered submeshes and insert them in detected positions
+  // iterate on ordered sub-meshes and insert them in detected positions
   map< int, TPosInList >::iterator i_pos = sortedPos.begin();
   for ( ; onlyBIt != onlyEIt; ++onlyBIt, ++i_pos )
     *(i_pos->second) = *onlyBIt;
@@ -2386,18 +2386,27 @@ bool SMESH_Mesh::IsOrderOK( const SMESH_subMesh* smBefore,
                             const SMESH_subMesh* smAfter ) const
 {
   TListOfListOfInt::const_iterator listIdsIt = _mySubMeshOrder.begin();
-  TListOfInt::const_iterator idBef, idAft;
   for( ; listIdsIt != _mySubMeshOrder.end(); listIdsIt++)
   {
     const TListOfInt& listOfId = *listIdsIt;
-    idBef = std::find( listOfId.begin(), listOfId.end(), smBefore->GetId() );
-    if ( idBef != listOfId.end() )
-      idAft = std::find( listOfId.begin(), listOfId.end(), smAfter->GetId() );
-    if ( idAft != listOfId.end () )
-      return ( std::distance( listOfId.begin(), idBef ) <
-               std::distance( listOfId.begin(), idAft )   );
+    int iB = -1, iA = -1, i = 0;
+    for ( TListOfInt::const_iterator id = listOfId.begin(); id != listOfId.end(); ++id, ++i )
+    {
+      if ( *id == smBefore->GetId() )
+      {
+        iB = i;
+        if ( iA > -1 )
+          return iB < iA;
+      }
+      else if ( *id == smAfter->GetId() )
+      {
+        iA = i;
+        if ( iB > -1 )
+          return iB < iA;
+      }
+    }
   }
-  return true; // no order imposed to given submeshes
+  return true; // no order imposed to given sub-meshes
 } 
 
 //=============================================================================
index e9e9f58b62b6a572fbe351c7693e6f1828c0e235..67df8daebb388f8127ad0e94026a6bcb89fea7ec 100644 (file)
@@ -2930,7 +2930,7 @@ bool SMESH_MesherHelper::IsReversedSubMesh (const TopoDS_Face& theFace)
         TopoDS_Shape s0 = GetSubShapeByNode( nn[0], GetMeshDS() );
         TopoDS_Shape s1 = GetSubShapeByNode( nn[1], GetMeshDS() );
         TopoDS_Shape  E = GetCommonAncestor( s0, s1, *myMesh, TopAbs_EDGE );
-        if ( !E.IsNull() && !s0.IsSame( s1 ))
+        if ( !E.IsNull() && !s0.IsSame( s1 ) && E.Orientation() != TopAbs_INTERNAL )
         {
           // is E seam edge?
           int nb = 0;
index 66bae01e026929239ce2c94a07d007c5e86c18a5..b7a492654a7cae3da753594a2ff3f7d5fbe55ed9 100644 (file)
@@ -448,8 +448,11 @@ const std::vector<UVPtStruct>& StdMeshers_FaceSide::GetUVPtStruct(bool   isXCons
     if ((int) u2node.size() + nbProxyNodes != myNbPonits &&
         (int) u2node.size() + nbProxyNodes != NbPoints( /*update=*/true ))
     {
-      MESSAGE("Wrong node parameters on edges, u2node.size():"
-              <<u2node.size()<<" !=  myNbPonits:"<<myNbPonits);
+      return myPoints;
+    }
+    if (( myNbPonits > 0 ) &&
+        ( u2node.begin()->first < 0 || u2node.rbegin()->first > 1 ))
+    {
       return myPoints;
     }
 
index 3d6249a79570ac2574e98f7b9b75965bef54b8cc..75219db95b4374617c7a473e00266b303a097d4b 100644 (file)
@@ -2544,8 +2544,8 @@ namespace // utils used by StdMeshers_Prism_3D::IsApplicable()
     bool        _isBase;  /* is used in a base FACE */
     EdgeWithNeighbors(const TopoDS_Edge& E, int iE, int nbE, int shift, bool isBase ):
       _edge( E ), _iBase( iE + shift ),
-      _iL( SMESH_MesherHelper::WrapIndex( iE-1, nbE ) + shift ),
-      _iR( SMESH_MesherHelper::WrapIndex( iE+1, nbE ) + shift ),
+      _iL( SMESH_MesherHelper::WrapIndex( iE-1, Max( 1, nbE )) + shift ),
+      _iR( SMESH_MesherHelper::WrapIndex( iE+1, Max( 1, nbE )) + shift ),
       _isBase( isBase )
     {
     }
@@ -2674,16 +2674,16 @@ namespace // utils used by StdMeshers_Prism_3D::IsApplicable()
           edges[ iFirst ]._iL = edges[ iFirst ]._iBase; // connect to self
           edges[ iLast  ]._iR = edges[ iLast ]._iBase;
 
-          // look for an EDGE of the outer WIRE connected to vv
+          // look for an EDGE of the outer WIREs connected to vv
           TopoDS_Vertex v0, v1;
-          for ( iE = 0; iE < nbEdgesInWires.front(); ++iE )
+          for ( iE = 0; iE < iFirst; ++iE )
           {
             v0 = SMESH_MesherHelper::IthVertex( 0, edges[ iE ]._edge );
             v1 = SMESH_MesherHelper::IthVertex( 1, edges[ iE ]._edge );
             if ( vv[0].IsSame( v0 ) || vv[0].IsSame( v1 ))
               edges[ iFirst ]._iL = edges[ iE ]._iBase;
             if ( vv[1].IsSame( v0 ) || vv[1].IsSame( v1 ))
-              edges[ iLast ]._iR = edges[ iE ]._iBase;
+              edges[ iLast  ]._iR = edges[ iE ]._iBase;
           }
         }
         iFirst += *nbE;
index 6ef4c57320b90fe5de63ab40d0d8d96109feea3c..31db57e967fd3f2702c6bd6fbf8d17131acccf48 100644 (file)
 
 #include <utilities.h>
 
+#include <BRepAdaptor_Curve.hxx>
 #include <BRepAdaptor_Surface.hxx>
 #include <BRepMesh_Delaun.hxx>
 #include <BRep_Tool.hxx>
 #include <Bnd_B2d.hxx>
 #include <GeomAPI_ProjectPointOnSurf.hxx>
 #include <GeomLib_IsPlanarSurface.hxx>
+#include <Precision.hxx>
 #include <TopExp.hxx>
 #include <TopExp_Explorer.hxx>
 #include <TopTools_DataMapIteratorOfDataMapOfShapeShape.hxx>
@@ -656,6 +658,8 @@ namespace {
 
     SMESH_MesherHelper srcHelper( *srcMesh );
     srcHelper.SetSubShape( srcFace );
+    SMESH_MesherHelper edgeHelper( *tgtMesh );
+    edgeHelper.ToFixNodeParameters( true );
 
     const SMDS_MeshNode* nullNode = 0;
     TAssocTool::TNodeNodeMap::iterator srcN_tgtN;
@@ -694,10 +698,30 @@ namespace {
           }
           case SMDS_TOP_EDGE:
           {
-            const TopoDS_Shape & srcE = srcMeshDS->IndexToShape( srcNode->getshapeId() );
-            const TopoDS_Shape & tgtE = shape2ShapeMap( srcE, /*isSrc=*/true );
-            double srcU = srcHelper.GetNodeU( TopoDS::Edge( srcE ), srcNode );
-            tgtMeshDS->SetNodeOnEdge( n, TopoDS::Edge( tgtE ), srcU );
+            const TopoDS_Edge& srcE = TopoDS::Edge( srcMeshDS->IndexToShape( srcNode->getshapeId()));
+            const TopoDS_Edge& tgtE = TopoDS::Edge( shape2ShapeMap( srcE, /*isSrc=*/true ));
+            double srcU = srcHelper.GetNodeU( srcE, srcNode );
+            tgtMeshDS->SetNodeOnEdge( n, tgtE, srcU );
+            if ( !tgtFace.IsPartner( srcFace ))
+            {
+              bool isOk = true;
+              edgeHelper.SetSubShape( tgtE );
+              edgeHelper.GetNodeU( tgtE, n, 0, &isOk );
+              if ( !isOk ) // projection of n to tgtE failed (23395)
+              {
+                double sF, sL, tF, tL;
+                BRep_Tool::Range( srcE, sF, sL );
+                BRep_Tool::Range( tgtE, tF, tL );
+                double srcR = ( srcU - sF ) / ( sL - sF );
+                double tgtU  = tF + srcR * ( tL - tF );
+                tgtMeshDS->SetNodeOnEdge( n, tgtE, tgtU );
+                gp_Pnt newP = BRepAdaptor_Curve( tgtE ).Value( tgtU );
+                double dist = newP.Distance( tgtP );
+                double tol = BRep_Tool::Tolerance( tgtE );
+                if ( tol < dist && dist < 1000*tol )
+                  tgtMeshDS->MoveNode( n, newP.X(), newP.Y(), newP.Z() );
+              }
+            }
             break;
           }
           case SMDS_TOP_VERTEX:
@@ -728,12 +752,9 @@ namespace {
 
     if ( !tgtFace.IsPartner( srcFace ) )
     {
-      SMESH_MesherHelper edgeHelper( *tgtMesh );
-      edgeHelper.ToFixNodeParameters( true );
       helper.ToFixNodeParameters( true );
 
       int nbOkPos = 0;
-      bool toCheck = true;
       const double tol2d = 1e-12;
       srcN_tgtN = src2tgtNodes.begin();
       for ( ; srcN_tgtN != src2tgtNodes.end(); ++srcN_tgtN )
@@ -754,9 +775,9 @@ namespace {
         }
         case SMDS_TOP_EDGE:
         {
-          const TopoDS_Edge & tgtE = TopoDS::Edge( tgtMeshDS->IndexToShape( n->getshapeId() ));
-          edgeHelper.SetSubShape( tgtE );
-          edgeHelper.GetNodeU( tgtE, n, 0, &toCheck );
+          // const TopoDS_Edge & tgtE = TopoDS::Edge( tgtMeshDS->IndexToShape( n->getshapeId() ));
+          // edgeHelper.SetSubShape( tgtE );
+          // edgeHelper.GetNodeU( tgtE, n, 0, &toCheck );
           break;
         }
         default:;
index 279e556f3ff22536730db97e7694bb73f6771938..8c0d9778b3803e30017311df07924dc37e5c28ee 100644 (file)
@@ -287,11 +287,14 @@ void StdMeshers_ProjectionSource1D_i::LoadFrom( const char* theStream )
   }
 
   myCorbaMesh = SMESH::SMESH_Mesh::_duplicate( mesh );
-  GetImpl()->SetSourceMesh       ( meshImpl );
-  GetImpl()->SetSourceEdge       ( shapes[ SRC_EDGE ] );
-  GetImpl()->SetVertexAssociation( shapes[ SRC_VERTEX ],
-                                   shapes[ TGT_VERTEX ]);
-
+  try {
+    GetImpl()->SetSourceMesh       ( meshImpl );
+    GetImpl()->SetSourceEdge       ( shapes[ SRC_EDGE ] );
+    GetImpl()->SetVertexAssociation( shapes[ SRC_VERTEX ],
+                                     shapes[ TGT_VERTEX ]);
+  }
+  catch (...) {
+  }
   myBaseImpl->LoadFrom( is );
 
   std::istringstream str( theStream );
index 93cd392c0f14609ede6d0eebaba4651f9b684eed..0c030a040e9d69c63118b0efa09b1cd27eb0c190 100644 (file)
@@ -298,12 +298,16 @@ void StdMeshers_ProjectionSource2D_i::LoadFrom( const char* theStream )
 
   myCorbaMesh = SMESH::SMESH_Mesh::_duplicate( mesh );
 
-  GetImpl()->SetSourceMesh       ( meshImpl );
-  GetImpl()->SetSourceFace       ( shapes[ SRC_FACE ] );
-  GetImpl()->SetVertexAssociation( shapes[ SRC_VERTEX1 ],
-                                   shapes[ SRC_VERTEX2 ],
-                                   shapes[ TGT_VERTEX1 ],
-                                   shapes[ TGT_VERTEX2 ]);
+  try {
+    GetImpl()->SetSourceMesh       ( meshImpl );
+    GetImpl()->SetSourceFace       ( shapes[ SRC_FACE ] );
+    GetImpl()->SetVertexAssociation( shapes[ SRC_VERTEX1 ],
+                                     shapes[ SRC_VERTEX2 ],
+                                     shapes[ TGT_VERTEX1 ],
+                                     shapes[ TGT_VERTEX2 ]);
+  }
+  catch( ... ) {
+  }
   myBaseImpl->LoadFrom( is );
 
   std::istringstream str( theStream );
index a97d95104366e28687341289edfafe334daa61b1..84ad6ca62d67acd3b8e16db247ba11a026b60fc8 100644 (file)
@@ -252,7 +252,7 @@ CORBA::Boolean StdMeshers_ProjectionSource3D_i::IsDimSupported( SMESH::Dimension
 //================================================================================
 /*!
  * \brief Write parameters in a string
 * \retval char* - resulting string
+ * \retval char* - resulting string
  */
 //================================================================================
 
@@ -273,7 +273,7 @@ char* StdMeshers_ProjectionSource3D_i::SaveTo()
 //================================================================================
 /*!
  * \brief Retrieve parameters from the string
 * \param theStream - the input string
+ * \param theStream - the input string
  */
 //================================================================================
 
@@ -299,13 +299,16 @@ void StdMeshers_ProjectionSource3D_i::LoadFrom( const char* theStream )
 
   myCorbaMesh = SMESH::SMESH_Mesh::_duplicate( mesh );
 
-  GetImpl()->SetSourceMesh       ( meshImpl );
-  GetImpl()->SetSource3DShape    ( shapes[ SRC_SHAPE3D ] );
-  GetImpl()->SetVertexAssociation( shapes[ SRC_VERTEX1 ],
-                                   shapes[ SRC_VERTEX2 ],
-                                   shapes[ TGT_VERTEX1 ],
-                                   shapes[ TGT_VERTEX2 ]);
-
+  try {
+    GetImpl()->SetSourceMesh       ( meshImpl );
+    GetImpl()->SetSource3DShape    ( shapes[ SRC_SHAPE3D ] );
+    GetImpl()->SetVertexAssociation( shapes[ SRC_VERTEX1 ],
+                                     shapes[ SRC_VERTEX2 ],
+                                     shapes[ TGT_VERTEX1 ],
+                                     shapes[ TGT_VERTEX2 ]);
+  }
+  catch (...) {
+  }
   myBaseImpl->LoadFrom( is );
 
   std::istringstream str( theStream );
index ebff90789f44fad934783cb8c1ef1af5a23d347f..7f17b00f6ba44f15daa0366a03315350e9eefa54 100644 (file)
@@ -379,7 +379,11 @@ void StdMeshers_QuadrangleParams_i::LoadFrom( const char* theStream )
       for ( int i = 0; i < nb; ++i )
         shapes.push_back( StdMeshers_ObjRefUlils::LoadFromStream( is, & myShapeEntries[i] ));
 
-      GetImpl()->SetEnforcedNodes( shapes, points );
+      try {
+        GetImpl()->SetEnforcedNodes( shapes, points );
+      }
+      catch (...) {
+      }
     }
   }