mirror of
https://git.salome-platform.org/gitpub/modules/smesh.git
synced 2025-01-03 21:30:34 +05:00
0020809: [CEA 400] ParaMEDSPLITTER failure on blade.med
* Fix SMESH_MesherHelper::GetNodeU() for a node on vertex of a closed geom edge double GetNodeU(const TopoDS_Edge& theEdge, const SMDS_MeshNode* theNode, + const SMDS_MeshNode* inEdgeNode=0, bool* check=0);
This commit is contained in:
parent
34aebed370
commit
00d55cb030
@ -597,20 +597,33 @@ gp_XY SMESH_MesherHelper::GetMiddleUV(const Handle(Geom_Surface)& surface,
|
|||||||
|
|
||||||
double SMESH_MesherHelper::GetNodeU(const TopoDS_Edge& E,
|
double SMESH_MesherHelper::GetNodeU(const TopoDS_Edge& E,
|
||||||
const SMDS_MeshNode* n,
|
const SMDS_MeshNode* n,
|
||||||
|
const SMDS_MeshNode* inEdgeNode,
|
||||||
bool* check)
|
bool* check)
|
||||||
{
|
{
|
||||||
double param = 0;
|
double param = 0;
|
||||||
const SMDS_PositionPtr Pos = n->GetPosition();
|
const SMDS_PositionPtr Pos = n->GetPosition();
|
||||||
if(Pos->GetTypeOfPosition()==SMDS_TOP_EDGE) {
|
if ( Pos->GetTypeOfPosition()==SMDS_TOP_EDGE )
|
||||||
|
{
|
||||||
const SMDS_EdgePosition* epos =
|
const SMDS_EdgePosition* epos =
|
||||||
static_cast<const SMDS_EdgePosition*>(n->GetPosition().get());
|
static_cast<const SMDS_EdgePosition*>(n->GetPosition().get());
|
||||||
param = epos->GetUParameter();
|
param = epos->GetUParameter();
|
||||||
}
|
}
|
||||||
else if(Pos->GetTypeOfPosition()==SMDS_TOP_VERTEX) {
|
else if( Pos->GetTypeOfPosition() == SMDS_TOP_VERTEX )
|
||||||
SMESHDS_Mesh * meshDS = GetMeshDS();
|
{
|
||||||
int vertexID = n->GetPosition()->GetShapeId();
|
if ( inEdgeNode && TopExp::FirstVertex( E ).IsSame( TopExp::LastVertex( E ))) // issue 0020809
|
||||||
const TopoDS_Vertex& V = TopoDS::Vertex(meshDS->IndexToShape(vertexID));
|
{
|
||||||
param = BRep_Tool::Parameter( V, E );
|
Standard_Real f,l;
|
||||||
|
BRep_Tool::Range( E, f,l );
|
||||||
|
double uInEdge = GetNodeU( E, inEdgeNode );
|
||||||
|
param = ( fabs( uInEdge - f ) < fabs( l - uInEdge )) ? f : l;
|
||||||
|
}
|
||||||
|
else
|
||||||
|
{
|
||||||
|
SMESHDS_Mesh * meshDS = GetMeshDS();
|
||||||
|
int vertexID = n->GetPosition()->GetShapeId();
|
||||||
|
const TopoDS_Vertex& V = TopoDS::Vertex(meshDS->IndexToShape(vertexID));
|
||||||
|
param = BRep_Tool::Parameter( V, E );
|
||||||
|
}
|
||||||
}
|
}
|
||||||
if ( check )
|
if ( check )
|
||||||
*check = CheckNodeU( E, n, param, BRep_Tool::Tolerance( E ));
|
*check = CheckNodeU( E, n, param, BRep_Tool::Tolerance( E ));
|
||||||
@ -739,8 +752,8 @@ const SMDS_MeshNode* SMESH_MesherHelper::GetMediumNode(const SMDS_MeshNode* n1,
|
|||||||
E = TopoDS::Edge(myShape);
|
E = TopoDS::Edge(myShape);
|
||||||
edgeID = myShapeID;
|
edgeID = myShapeID;
|
||||||
}
|
}
|
||||||
u[0] = GetNodeU(E,n1, force3d ? 0 : &uvOK[0]);
|
u[0] = GetNodeU(E,n1,n2, force3d ? 0 : &uvOK[0]);
|
||||||
u[1] = GetNodeU(E,n2, force3d ? 0 : &uvOK[1]);
|
u[1] = GetNodeU(E,n2,n1, force3d ? 0 : &uvOK[1]);
|
||||||
}
|
}
|
||||||
if(!force3d)
|
if(!force3d)
|
||||||
{
|
{
|
||||||
|
@ -269,6 +269,7 @@ public:
|
|||||||
*/
|
*/
|
||||||
double GetNodeU(const TopoDS_Edge& theEdge,
|
double GetNodeU(const TopoDS_Edge& theEdge,
|
||||||
const SMDS_MeshNode* theNode,
|
const SMDS_MeshNode* theNode,
|
||||||
|
const SMDS_MeshNode* inEdgeNode=0,
|
||||||
bool* check=0);
|
bool* check=0);
|
||||||
/*!
|
/*!
|
||||||
* \brief Return node UV on face
|
* \brief Return node UV on face
|
||||||
|
@ -270,7 +270,7 @@ const vector<UVPtStruct>& StdMeshers_FaceSide::GetUVPtStruct(bool isXConst,
|
|||||||
const SMDS_MeshNode* node = nItr->next();
|
const SMDS_MeshNode* node = nItr->next();
|
||||||
if ( myIgnoreMediumNodes && SMESH_MeshEditor::IsMedium( node, SMDSAbs_Edge ))
|
if ( myIgnoreMediumNodes && SMESH_MeshEditor::IsMedium( node, SMDSAbs_Edge ))
|
||||||
continue;
|
continue;
|
||||||
double u = helper.GetNodeU( myEdge[i], node, ¶mOK );
|
double u = helper.GetNodeU( myEdge[i], node, 0, ¶mOK );
|
||||||
double aLenU = GCPnts_AbscissaPoint::Length
|
double aLenU = GCPnts_AbscissaPoint::Length
|
||||||
( const_cast<GeomAdaptor_Curve&>( myC3dAdaptor[i]), myFirst[i], u );
|
( const_cast<GeomAdaptor_Curve&>( myC3dAdaptor[i]), myFirst[i], u );
|
||||||
if ( myEdgeLength[i] < aLenU ) // nonregression test "3D_mesh_NETGEN/G6"
|
if ( myEdgeLength[i] < aLenU ) // nonregression test "3D_mesh_NETGEN/G6"
|
||||||
@ -288,7 +288,7 @@ const vector<UVPtStruct>& StdMeshers_FaceSide::GetUVPtStruct(bool isXConst,
|
|||||||
const SMDS_MeshNode* node = nItr->next();
|
const SMDS_MeshNode* node = nItr->next();
|
||||||
if ( myIgnoreMediumNodes && SMESH_MeshEditor::IsMedium( node, SMDSAbs_Edge ))
|
if ( myIgnoreMediumNodes && SMESH_MeshEditor::IsMedium( node, SMDSAbs_Edge ))
|
||||||
continue;
|
continue;
|
||||||
double u = helper.GetNodeU( myEdge[i], node, ¶mOK );
|
double u = helper.GetNodeU( myEdge[i], node, 0, ¶mOK );
|
||||||
|
|
||||||
// paramSize is signed so orientation is taken into account
|
// paramSize is signed so orientation is taken into account
|
||||||
double normPar = prevNormPar + r * ( u - myFirst[i] ) / paramSize;
|
double normPar = prevNormPar + r * ( u - myFirst[i] ) / paramSize;
|
||||||
|
Loading…
Reference in New Issue
Block a user