23239: [CEA 1739] Regression : crash trying to create mesh

(StdMesher_Prism_3D.cxx)

+ minor changes
This commit is contained in:
eap 2016-02-17 20:36:53 +03:00 committed by vsr
parent 920fe932b1
commit 04b68235aa
4 changed files with 24 additions and 17 deletions

View File

@ -344,7 +344,7 @@ class SMESH_EXPORT SMESH_Mesh
bool IsOrderOK( const SMESH_subMesh* smBefore,
const SMESH_subMesh* smAfter ) const;
ostream& Dump(ostream & save);
std::ostream& Dump(ostream & save);
private:

View File

@ -47,7 +47,7 @@
#include <TopoDS_Wire.hxx>
#ifdef _DEBUG_
#define _MYDEBUG_
//#define _MYDEBUG_
#include "SMESH_File.hxx"
#include "SMESH_Comment.hxx"
#endif
@ -162,6 +162,8 @@ namespace
if ( inSegments.size() > 1000 )
return;
const char* fileName = "/misc/dn25/salome/eap/salome/misc/Code/C++/MAdebug.txt";
const char* user = getenv("USER");
if ( !user || strcmp( user, "eap" )) return;
SMESH_File file(fileName, false );
file.remove();
file.openForWriting();
@ -1518,9 +1520,9 @@ Adaptor3d_Curve* SMESH_MAT2d::MedialAxis::make3DCurve(const Branch& branch) cons
*/
//================================================================================
void SMESH_MAT2d::Branch::init( vector<const TVDEdge*>& maEdges,
const Boundary* boundary,
map< const TVDVertex*, BranchEndType > endType )
void SMESH_MAT2d::Branch::init( vector<const TVDEdge*>& maEdges,
const Boundary* boundary,
map< const TVDVertex*, BranchEndType >& endType )
{
if ( maEdges.empty() ) return;

View File

@ -120,9 +120,9 @@ namespace SMESH_MAT2d
public: // internal: construction
void init( std::vector<const TVDEdge*>& maEdges,
const Boundary* boundary,
std::map< const TVDVertex*, BranchEndType > endType);
void init( std::vector<const TVDEdge*>& maEdges,
const Boundary* boundary,
std::map< const TVDVertex*, BranchEndType >& endType);
void setBranchesToEnds( const std::vector< Branch >& branches);
BranchPoint getPoint( const TVDVertex* vertex ) const;
void setRemoved( const BranchPoint& proxyPoint );

View File

@ -2512,16 +2512,17 @@ namespace // utils used by StdMeshers_Prism_3D::IsApplicable()
}
EdgeWithNeighbors() {}
};
struct PrismSide
// PrismSide contains all FACEs linking a bottom EDGE with a top one.
struct PrismSide
{
TopoDS_Face _face;
TopTools_IndexedMapOfShape *_faces; // pointer because its copy constructor is private
TopoDS_Edge _topEdge;
vector< EdgeWithNeighbors >*_edges;
int _iBotEdge;
vector< bool > _isCheckedEdge;
TopoDS_Face _face; // a currently treated upper FACE
TopTools_IndexedMapOfShape *_faces; // all FACEs (pointer because of a private copy constructor)
TopoDS_Edge _topEdge; // a current top EDGE
vector< EdgeWithNeighbors >*_edges; // all EDGEs of _face
int _iBotEdge; // index of _topEdge within _edges
vector< bool > _isCheckedEdge; // mark EDGEs whose two owner FACEs found
int _nbCheckedEdges; // nb of EDGEs whose location is defined
PrismSide *_leftSide;
PrismSide *_leftSide; // neighbor sides
PrismSide *_rightSide;
void SetExcluded() { _leftSide = _rightSide = NULL; }
bool IsExcluded() const { return !_leftSide; }
@ -2679,7 +2680,7 @@ bool StdMeshers_Prism_3D::IsApplicable(const TopoDS_Shape & shape, bool toCheckA
typedef vector< EdgeWithNeighbors > TEdgeWithNeighborsVec;
vector< TEdgeWithNeighborsVec > faceEdgesVec( allFaces.Extent() + 1 );
const size_t nbEdgesMax = facesOfEdge.Extent() * 2; // there can be seam EDGES
const size_t nbEdgesMax = facesOfEdge.Extent() * 2; // there can be seam EDGEs
TopTools_IndexedMapOfShape* facesOfSide = new TopTools_IndexedMapOfShape[ nbEdgesMax ];
SMESHUtils::ArrayDeleter<TopTools_IndexedMapOfShape> delFacesOfSide( facesOfSide );
@ -2815,6 +2816,10 @@ bool StdMeshers_Prism_3D::IsApplicable(const TopoDS_Shape & shape, bool toCheckA
side._isCheckedEdge[ side._iBotEdge ] = true;
side._nbCheckedEdges = 1; // bottom EDGE is known
}
else // probably a triangular top face found
{
side._face.Nullify();
}
side._topEdge.Nullify();
isOK = ( !side._edges->empty() || side._faces->Extent() > 1 );