22830: EDF 9557 SMESH: Quadratic conversion of a mesh fails

In SMESH_MeshEditor copy node positions in Duplicate*()
In SMESH_MesherHelper
- cache face tolerance
- in GetNodeUV() prevent SIGSEGV if a node is on EDGE not belonging to a FACE
- in GetNodeUV() check if a VERTEX is a seam in a given FACE, not in another
- avoid returning myShape from GetMediumPos() when it can be wrong

52592: Invalid values of Length 2D control

- avoid using not initialized values for showing control values of a group
This commit is contained in:
eap 2014-12-29 18:22:30 +03:00
parent a00f90ebab
commit 133cb38139
6 changed files with 185 additions and 92 deletions

View File

@ -398,7 +398,8 @@ SMESH_DeviceActor
anIdList->SetNumberOfIds(2); anIdList->SetNumberOfIds(2);
Length2D::TValues::const_iterator anIter = aValues.begin(); Length2D::TValues::const_iterator anIter = aValues.begin();
for(vtkIdType aVtkId = 0; anIter != aValues.end(); anIter++,aVtkId++){ aNbCells = 0;
for(; anIter != aValues.end(); anIter++){
const Length2D::Value& aValue = *anIter; const Length2D::Value& aValue = *anIter;
int aNode[2] = { int aNode[2] = {
myVisualObj->GetNodeVTKId(aValue.myPntId[0]), myVisualObj->GetNodeVTKId(aValue.myPntId[0]),
@ -409,9 +410,12 @@ SMESH_DeviceActor
anIdList->SetId( 1, aNode[1] ); anIdList->SetId( 1, aNode[1] );
aConnectivity->InsertNextCell( anIdList ); aConnectivity->InsertNextCell( anIdList );
aCellTypesArray->InsertNextValue( VTK_LINE ); aCellTypesArray->InsertNextValue( VTK_LINE );
aScalars->SetValue(aVtkId,aValue.myLength); aScalars->SetValue(aNbCells,aValue.myLength);
aNbCells++;
} }
} }
aCellTypesArray->SetNumberOfTuples( aNbCells );
aScalars->SetNumberOfTuples( aNbCells );
VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New(); VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New();
aCellLocationsArray->SetNumberOfComponents( 1 ); aCellLocationsArray->SetNumberOfComponents( 1 );
@ -421,7 +425,7 @@ SMESH_DeviceActor
for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ ) for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) ); aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
aDataSet->SetCells( aCellTypesArray, aCellLocationsArray,aConnectivity ); aDataSet->SetCells( aCellTypesArray, aCellLocationsArray, aConnectivity );
SetUnstructuredGrid(aDataSet); SetUnstructuredGrid(aDataSet);
aDataSet->GetCellData()->SetScalars(aScalars); aDataSet->GetCellData()->SetScalars(aScalars);
@ -458,7 +462,8 @@ SMESH_DeviceActor
anIdList->SetNumberOfIds(2); anIdList->SetNumberOfIds(2);
MultiConnection2D::MValues::const_iterator anIter = aValues.begin(); MultiConnection2D::MValues::const_iterator anIter = aValues.begin();
for(vtkIdType aVtkId = 0; anIter != aValues.end(); anIter++,aVtkId++){ aNbCells = 0;
for(; anIter != aValues.end(); anIter++){
const MultiConnection2D::Value& aValue = (*anIter).first; const MultiConnection2D::Value& aValue = (*anIter).first;
int aNode[2] = { int aNode[2] = {
myVisualObj->GetNodeVTKId(aValue.myPntId[0]), myVisualObj->GetNodeVTKId(aValue.myPntId[0]),
@ -469,9 +474,12 @@ SMESH_DeviceActor
anIdList->SetId( 1, aNode[1] ); anIdList->SetId( 1, aNode[1] );
aConnectivity->InsertNextCell( anIdList ); aConnectivity->InsertNextCell( anIdList );
aCellTypesArray->InsertNextValue( VTK_LINE ); aCellTypesArray->InsertNextValue( VTK_LINE );
aScalars->SetValue(aVtkId,(*anIter).second); aScalars->SetValue( aNbCells,(*anIter).second);
aNbCells++;
} }
} }
aCellTypesArray->SetNumberOfTuples( aNbCells );
aScalars->SetNumberOfTuples( aNbCells );
VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New(); VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New();
aCellLocationsArray->SetNumberOfComponents( 1 ); aCellLocationsArray->SetNumberOfComponents( 1 );

View File

@ -10228,6 +10228,7 @@ bool SMESH_MeshEditor::doubleNodes( SMESHDS_Mesh* theMeshDS,
{ {
// duplicate node // duplicate node
aNewNode = theMeshDS->AddNode( aCurrNode->X(), aCurrNode->Y(), aCurrNode->Z() ); aNewNode = theMeshDS->AddNode( aCurrNode->X(), aCurrNode->Y(), aCurrNode->Z() );
copyPosition( aCurrNode, aNewNode );
theNodeNodeMap[ aCurrNode ] = aNewNode; theNodeNodeMap[ aCurrNode ] = aNewNode;
myLastCreatedNodes.Append( aNewNode ); myLastCreatedNodes.Append( aNewNode );
} }
@ -10240,10 +10241,8 @@ bool SMESH_MeshEditor::doubleNodes( SMESHDS_Mesh* theMeshDS,
if ( theIsDoubleElem ) if ( theIsDoubleElem )
AddElement(newNodes, anElem->GetType(), anElem->IsPoly()); AddElement(newNodes, anElem->GetType(), anElem->IsPoly());
else else
{
MESSAGE("ChangeElementNodes");
theMeshDS->ChangeElementNodes( anElem, &newNodes[ 0 ], anElem->NbNodes() ); theMeshDS->ChangeElementNodes( anElem, &newNodes[ 0 ], anElem->NbNodes() );
}
res = true; res = true;
} }
return res; return res;
@ -10254,8 +10253,8 @@ bool SMESH_MeshEditor::doubleNodes( SMESHDS_Mesh* theMeshDS,
\brief Creates a hole in a mesh by doubling the nodes of some particular elements \brief Creates a hole in a mesh by doubling the nodes of some particular elements
\param theNodes - identifiers of nodes to be doubled \param theNodes - identifiers of nodes to be doubled
\param theModifiedElems - identifiers of elements to be updated by the new (doubled) \param theModifiedElems - identifiers of elements to be updated by the new (doubled)
nodes. If list of element identifiers is empty then nodes are doubled but nodes. If list of element identifiers is empty then nodes are doubled but
they not assigned to elements they not assigned to elements
\return TRUE if operation has been completed successfully, FALSE otherwise \return TRUE if operation has been completed successfully, FALSE otherwise
*/ */
//================================================================================ //================================================================================
@ -10291,6 +10290,7 @@ bool SMESH_MeshEditor::DoubleNodes( const std::list< int >& theListOfNodes,
const SMDS_MeshNode* aNewNode = aMeshDS->AddNode( aNode->X(), aNode->Y(), aNode->Z() ); const SMDS_MeshNode* aNewNode = aMeshDS->AddNode( aNode->X(), aNode->Y(), aNode->Z() );
if ( aNewNode ) if ( aNewNode )
{ {
copyPosition( aNode, aNewNode );
anOldNodeToNewNode[ aNode ] = aNewNode; anOldNodeToNewNode[ aNode ] = aNewNode;
myLastCreatedNodes.Append( aNewNode ); myLastCreatedNodes.Append( aNewNode );
} }
@ -10909,6 +10909,7 @@ bool SMESH_MeshEditor::DoubleNodesOnGroupBoundaries( const std::vector<TIDSorted
} }
double *coords = grid->GetPoint(oldId); double *coords = grid->GetPoint(oldId);
SMDS_MeshNode *newNode = meshDS->AddNode(coords[0], coords[1], coords[2]); SMDS_MeshNode *newNode = meshDS->AddNode(coords[0], coords[1], coords[2]);
copyPosition( meshDS->FindNodeVtk( oldId ), newNode );
int newId = newNode->getVtkId(); int newId = newNode->getVtkId();
nodeDomains[oldId][idom] = newId; // cloned node for other domains nodeDomains[oldId][idom] = newId; // cloned node for other domains
//MESSAGE("-+-+-c oldNode " << oldId << " domain " << idomain << " newNode " << newId << " domain " << idom << " size=" <<nodeDomains[oldId].size()); //MESSAGE("-+-+-c oldNode " << oldId << " domain " << idomain << " newNode " << newId << " domain " << idom << " size=" <<nodeDomains[oldId].size());
@ -11322,6 +11323,7 @@ bool SMESH_MeshEditor::CreateFlatElementsOnFacesGroups(const std::vector<TIDSort
if (!clonedNodes.count(node)) if (!clonedNodes.count(node))
{ {
clone = meshDS->AddNode(node->X(), node->Y(), node->Z()); clone = meshDS->AddNode(node->X(), node->Y(), node->Z());
copyPosition( node, clone );
clonedNodes[node] = clone; clonedNodes[node] = clone;
} }
else else
@ -11338,6 +11340,7 @@ bool SMESH_MeshEditor::CreateFlatElementsOnFacesGroups(const std::vector<TIDSort
if (!intermediateNodes.count(node)) if (!intermediateNodes.count(node))
{ {
inter = meshDS->AddNode(node->X(), node->Y(), node->Z()); inter = meshDS->AddNode(node->X(), node->Y(), node->Z());
copyPosition( node, inter );
intermediateNodes[node] = inter; intermediateNodes[node] = inter;
} }
else else
@ -12282,3 +12285,45 @@ int SMESH_MeshEditor::MakeBoundaryMesh(const TIDSortedElemSet& elements,
} }
return nbAddedBnd; return nbAddedBnd;
} }
//================================================================================
/*!
* \brief Copy node position and set \a to node on the same geometry
*/
//================================================================================
void SMESH_MeshEditor::copyPosition( const SMDS_MeshNode* from,
const SMDS_MeshNode* to )
{
if ( !from || !to ) return;
SMDS_PositionPtr pos = from->GetPosition();
if ( !pos || from->getshapeId() < 1 ) return;
switch ( pos->GetTypeOfPosition() )
{
case SMDS_TOP_3DSPACE: break;
case SMDS_TOP_FACE:
{
const SMDS_FacePosition* fPos = static_cast< const SMDS_FacePosition* >( pos );
GetMeshDS()->SetNodeOnFace( to, from->getshapeId(),
fPos->GetUParameter(), fPos->GetVParameter() );
break;
}
case SMDS_TOP_EDGE:
{
// WARNING: it is dangerous to set equal nodes on one EDGE!!!!!!!!
const SMDS_EdgePosition* ePos = static_cast< const SMDS_EdgePosition* >( pos );
GetMeshDS()->SetNodeOnEdge( to, from->getshapeId(), ePos->GetUParameter() );
break;
}
case SMDS_TOP_VERTEX:
{
GetMeshDS()->SetNodeOnVertex( to, from->getshapeId() );
break;
}
case SMDS_TOP_UNSPEC:
default:;
}
}

View File

@ -686,6 +686,9 @@ public:
std::map< const SMDS_MeshNode*, const SMDS_MeshNode* >& theNodeNodeMap, std::map< const SMDS_MeshNode*, const SMDS_MeshNode* >& theNodeNodeMap,
const bool theIsDoubleElem ); const bool theIsDoubleElem );
void copyPosition( const SMDS_MeshNode* from,
const SMDS_MeshNode* to );
private: private:
SMESH_Mesh * myMesh; SMESH_Mesh * myMesh;

View File

@ -259,6 +259,7 @@ void SMESH_MesherHelper::SetSubShape(const TopoDS_Shape& aSh)
{ {
// look for a "seam" edge, a real seam or an edge on period boundary // look for a "seam" edge, a real seam or an edge on period boundary
TopoDS_Edge edge = TopoDS::Edge( exp.Current() ); TopoDS_Edge edge = TopoDS::Edge( exp.Current() );
const int edgeID = meshDS->ShapeToIndex( edge );
if ( myParIndex ) if ( myParIndex )
{ {
BRep_Tool::UVPoints( edge, face, uv1, uv2 ); BRep_Tool::UVPoints( edge, face, uv1, uv2 );
@ -305,7 +306,6 @@ void SMESH_MesherHelper::SetSubShape(const TopoDS_Shape& aSh)
if ( isSeam ) if ( isSeam )
{ {
// store seam shape indices, negative if shape encounters twice // store seam shape indices, negative if shape encounters twice
int edgeID = meshDS->ShapeToIndex( edge );
mySeamShapeIds.insert( IsSeamShape( edgeID ) ? -edgeID : edgeID ); mySeamShapeIds.insert( IsSeamShape( edgeID ) ? -edgeID : edgeID );
for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() ) { for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() ) {
int vertexID = meshDS->ShapeToIndex( v.Current() ); int vertexID = meshDS->ShapeToIndex( v.Current() );
@ -315,10 +315,15 @@ void SMESH_MesherHelper::SetSubShape(const TopoDS_Shape& aSh)
} }
// look for a degenerated edge // look for a degenerated edge
if ( SMESH_Algo::isDegenerated( edge )) { if ( SMESH_Algo::isDegenerated( edge )) {
myDegenShapeIds.insert( meshDS->ShapeToIndex( edge )); myDegenShapeIds.insert( edgeID );
for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() ) for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() )
myDegenShapeIds.insert( meshDS->ShapeToIndex( v.Current() )); myDegenShapeIds.insert( meshDS->ShapeToIndex( v.Current() ));
} }
if ( !BRep_Tool::SameParameter( edge ) ||
!BRep_Tool::SameRange( edge ))
{
setPosOnShapeValidity( edgeID, false );
}
} }
} }
} }
@ -527,11 +532,11 @@ void SMESH_MesherHelper::ToFixNodeParameters(bool toFix)
//======================================================================= //=======================================================================
//function : GetUVOnSeam //function : getUVOnSeam
//purpose : Select UV on either of 2 pcurves of a seam edge, closest to the given UV //purpose : Select UV on either of 2 pcurves of a seam edge, closest to the given UV
//======================================================================= //=======================================================================
gp_Pnt2d SMESH_MesherHelper::GetUVOnSeam( const gp_Pnt2d& uv1, const gp_Pnt2d& uv2 ) const gp_Pnt2d SMESH_MesherHelper::getUVOnSeam( const gp_Pnt2d& uv1, const gp_Pnt2d& uv2 ) const
{ {
gp_Pnt2d result = uv1; gp_Pnt2d result = uv1;
for ( int i = U_periodic; i <= V_periodic ; ++i ) for ( int i = U_periodic; i <= V_periodic ; ++i )
@ -568,38 +573,34 @@ gp_XY SMESH_MesherHelper::GetNodeUV(const TopoDS_Face& F,
const SMDS_PositionPtr Pos = n->GetPosition(); const SMDS_PositionPtr Pos = n->GetPosition();
bool uvOK = false; bool uvOK = false;
if(Pos->GetTypeOfPosition()==SMDS_TOP_FACE) if ( Pos->GetTypeOfPosition() == SMDS_TOP_FACE )
{ {
// node has position on face // node has position on face
const SMDS_FacePosition* fpos = const SMDS_FacePosition* fpos = static_cast<const SMDS_FacePosition*>( Pos );
static_cast<const SMDS_FacePosition*>( Pos ); uv.SetCoord( fpos->GetUParameter(), fpos->GetVParameter() );
uv.SetCoord(fpos->GetUParameter(),fpos->GetVParameter());
if ( check ) if ( check )
uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), 10*MaxTolerance( F )); uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), 10*getFaceMaxTol( F ));
} }
else if(Pos->GetTypeOfPosition()==SMDS_TOP_EDGE) else if ( Pos->GetTypeOfPosition() == SMDS_TOP_EDGE )
{ {
// node has position on edge => it is needed to find // node has position on EDGE => it is needed to find
// corresponding edge from face, get pcurve for this // corresponding EDGE from FACE, get pcurve for this
// edge and retrieve value from this pcurve // EDGE and retrieve value from this pcurve
const SMDS_EdgePosition* epos = const SMDS_EdgePosition* epos = static_cast<const SMDS_EdgePosition*>( Pos );
static_cast<const SMDS_EdgePosition*>( Pos ); const int edgeID = n->getshapeId();
int edgeID = n->getshapeId(); const TopoDS_Edge& E = TopoDS::Edge( GetMeshDS()->IndexToShape( edgeID ));
TopoDS_Edge E = TopoDS::Edge(GetMeshDS()->IndexToShape(edgeID));
double f, l, u = epos->GetUParameter(); double f, l, u = epos->GetUParameter();
Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(E, F, f, l); Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface( E, F, f, l );
bool validU = ( f < u && u < l ); bool validU = ( !C2d.IsNull() && ( f < u ) && ( u < l ));
if ( validU ) if ( validU ) uv = C2d->Value( u );
uv = C2d->Value( u ); else uv.SetCoord( Precision::Infinite(),0.);
else
uv.SetCoord( Precision::Infinite(),0.);
if ( check || !validU ) if ( check || !validU )
uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), 10*MaxTolerance( F ),/*force=*/ !validU ); uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), 10*getFaceMaxTol( F ),/*force=*/ !validU );
// for a node on a seam edge select one of UVs on 2 pcurves // for a node on a seam EDGE select one of UVs on 2 pcurves
if ( n2 && IsSeamShape( edgeID ) ) if ( n2 && IsSeamShape( edgeID ))
{ {
uv = GetUVOnSeam( uv, GetNodeUV( F, n2, 0, check )); uv = getUVOnSeam( uv, GetNodeUV( F, n2, 0, check ));
} }
else else
{ // adjust uv to period { // adjust uv to period
@ -611,23 +612,22 @@ gp_XY SMESH_MesherHelper::GetNodeUV(const TopoDS_Face& F,
if ( isUPeriodic || isVPeriodic ) { if ( isUPeriodic || isVPeriodic ) {
Standard_Real UF,UL,VF,VL; Standard_Real UF,UL,VF,VL;
S->Bounds(UF,UL,VF,VL); S->Bounds(UF,UL,VF,VL);
if ( isUPeriodic ) if ( isUPeriodic ) newUV.SetX( uv.X() + ShapeAnalysis::AdjustToPeriod(uv.X(),UF,UL));
newUV.SetX( uv.X() + ShapeAnalysis::AdjustToPeriod(uv.X(),UF,UL)); if ( isVPeriodic ) newUV.SetY( uv.Y() + ShapeAnalysis::AdjustToPeriod(uv.Y(),VF,VL));
if ( isVPeriodic )
newUV.SetY( uv.Y() + ShapeAnalysis::AdjustToPeriod(uv.Y(),VF,VL)); if ( n2 )
} {
if ( n2 ) gp_Pnt2d uv2 = GetNodeUV( F, n2, 0, check );
{ if ( isUPeriodic && Abs( uv.X()-uv2.X() ) < Abs( newUV.X()-uv2.X() ))
gp_Pnt2d uv2 = GetNodeUV( F, n2, 0, check ); newUV.SetX( uv.X() );
if ( isUPeriodic && Abs( uv.X()-uv2.X() ) < Abs( newUV.X()-uv2.X() )) if ( isVPeriodic && Abs( uv.Y()-uv2.Y() ) < Abs( newUV.Y()-uv2.Y() ))
newUV.SetX( uv.X() ); newUV.SetY( uv.Y() );
if ( isVPeriodic && Abs( uv.Y()-uv2.Y() ) < Abs( newUV.Y()-uv2.Y() )) }
newUV.SetY( uv.Y() );
} }
uv = newUV; uv = newUV;
} }
} }
else if(Pos->GetTypeOfPosition()==SMDS_TOP_VERTEX) else if ( Pos->GetTypeOfPosition() == SMDS_TOP_VERTEX )
{ {
if ( int vertexID = n->getshapeId() ) { if ( int vertexID = n->getshapeId() ) {
const TopoDS_Vertex& V = TopoDS::Vertex(GetMeshDS()->IndexToShape(vertexID)); const TopoDS_Vertex& V = TopoDS::Vertex(GetMeshDS()->IndexToShape(vertexID));
@ -646,7 +646,7 @@ gp_XY SMESH_MesherHelper::GetNodeUV(const TopoDS_Face& F,
// get UV of a vertex closest to the node // get UV of a vertex closest to the node
double dist = 1e100; double dist = 1e100;
gp_Pnt pn = XYZ( n ); gp_Pnt pn = XYZ( n );
for ( TopExp_Explorer vert(F,TopAbs_VERTEX); !uvOK && vert.More(); vert.Next() ) { for ( TopExp_Explorer vert( F,TopAbs_VERTEX ); !uvOK && vert.More(); vert.Next() ) {
TopoDS_Vertex curV = TopoDS::Vertex( vert.Current() ); TopoDS_Vertex curV = TopoDS::Vertex( vert.Current() );
gp_Pnt p = BRep_Tool::Pnt( curV ); gp_Pnt p = BRep_Tool::Pnt( curV );
double curDist = p.SquareDistance( pn ); double curDist = p.SquareDistance( pn );
@ -675,13 +675,23 @@ gp_XY SMESH_MesherHelper::GetNodeUV(const TopoDS_Face& F,
} }
} }
} }
if ( n2 && IsSeamShape( vertexID ) ) if ( n2 && IsSeamShape( vertexID ))
uv = GetUVOnSeam( uv, GetNodeUV( F, n2, 0 )); {
bool isSeam = ( myShape.IsSame( F ));
if ( !isSeam ) {
SMESH_MesherHelper h( *myMesh );
h.SetSubShape( F );
isSeam = IsSeamShape( vertexID );
}
if ( isSeam )
uv = getUVOnSeam( uv, GetNodeUV( F, n2, 0 ));
}
} }
} }
else else
{ {
uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), 10*MaxTolerance( F )); uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), 10*getFaceMaxTol( F ));
} }
if ( check ) if ( check )
@ -750,7 +760,7 @@ bool SMESH_MesherHelper::CheckNodeUV(const TopoDS_Face& F,
const_cast<SMDS_MeshNode*>(n)->SetPosition const_cast<SMDS_MeshNode*>(n)->SetPosition
( SMDS_PositionPtr( new SMDS_FacePosition( U, V ))); ( SMDS_PositionPtr( new SMDS_FacePosition( U, V )));
} }
else if ( uv.Modulus() > numeric_limits<double>::min() ) else if ( myShape.IsSame(F) && uv.Modulus() > numeric_limits<double>::min() )
{ {
setPosOnShapeValidity( shapeID, true ); setPosOnShapeValidity( shapeID, true );
} }
@ -1037,9 +1047,12 @@ bool SMESH_MesherHelper::CheckNodeU(const TopoDS_Edge& E,
//======================================================================= //=======================================================================
//function : GetMediumPos //function : GetMediumPos
//purpose : Return index and type of the shape (EDGE or FACE only) to //purpose : Return index and type of the shape (EDGE or FACE only) to
// set a medium node on // set a medium node on
//param : useCurSubShape - if true, returns the shape set via SetSubShape() //param : useCurSubShape - if true, returns the shape set via SetSubShape()
// if any // if any
// calling GetMediumPos() with useCurSubShape=true is OK only for the
// case where the lower dim mesh is already constructed and converted to quadratic,
// else, nodes on EDGEs are assigned to FACE, for example.
//======================================================================= //=======================================================================
std::pair<int, TopAbs_ShapeEnum> std::pair<int, TopAbs_ShapeEnum>
@ -1442,11 +1455,9 @@ const SMDS_MeshNode* SMESH_MesherHelper::GetMediumNode(const SMDS_MeshNode* n1,
TopoDS_Edge E; double u [2]; TopoDS_Edge E; double u [2];
TopoDS_Face F; gp_XY uv[2]; TopoDS_Face F; gp_XY uv[2];
bool uvOK[2] = { false, false }; bool uvOK[2] = { false, false };
const bool useCurSubShape = ( !myShape.IsNull() && myShape.ShapeType() == TopAbs_EDGE );
pair<int, TopAbs_ShapeEnum> pos = GetMediumPos( n1, n2, mySetElemOnShape ); pair<int, TopAbs_ShapeEnum> pos = GetMediumPos( n1, n2, useCurSubShape );
// calling GetMediumPos() with useCurSubShape=mySetElemOnShape is OK only for the
// case where the lower dim mesh is already constructed, else, nodes on EDGEs are
// assigned to FACE, for example.
// get positions of the given nodes on shapes // get positions of the given nodes on shapes
if ( pos.second == TopAbs_FACE ) if ( pos.second == TopAbs_FACE )
@ -2848,6 +2859,24 @@ double SMESH_MesherHelper::MaxTolerance( const TopoDS_Shape& shape )
return tol; return tol;
} }
//================================================================================
/*!
* \brief Return MaxTolerance( face ), probably cached
*/
//================================================================================
double SMESH_MesherHelper::getFaceMaxTol( const TopoDS_Shape& face ) const
{
int faceID = GetMeshDS()->ShapeToIndex( face );
SMESH_MesherHelper* me = const_cast< SMESH_MesherHelper* >( this );
double & tol = me->myFaceMaxTol.insert( make_pair( faceID, -1. )).first->second;
if ( tol < 0 )
tol = MaxTolerance( face );
return tol;
}
//================================================================================ //================================================================================
/*! /*!
* \brief Return an angle between two EDGEs sharing a common VERTEX with reference * \brief Return an angle between two EDGEs sharing a common VERTEX with reference
@ -4606,6 +4635,8 @@ void SMESH_MesherHelper::FixQuadraticElements(SMESH_ComputeErrorPtr& compError,
// 3. Compute displacement of medium nodes // 3. Compute displacement of medium nodes
// --------------------------------------- // ---------------------------------------
SMESH_MesherHelper faceHlp(*myMesh);
// two loops on QFaces: the first is to treat boundary links, the second is for internal ones. // two loops on QFaces: the first is to treat boundary links, the second is for internal ones.
TopLoc_Location loc; TopLoc_Location loc;
bool checkUV; bool checkUV;
@ -4689,22 +4720,23 @@ void SMESH_MesherHelper::FixQuadraticElements(SMESH_ComputeErrorPtr& compError,
TopoDS_Face face; TopoDS_Face face;
if ( !isInside ) if ( !isInside )
{ {
// compute node displacement of end links of chain in parametric space of face // compute node displacement of end links of chain in parametric space of FACE
TChainLink& linkOnFace = *(++chain.begin()); TChainLink& linkOnFace = *(++chain.begin());
const SMDS_MeshNode* nodeOnFace = linkOnFace->_mediumNode; const SMDS_MeshNode* nodeOnFace = linkOnFace->_mediumNode;
TopoDS_Shape f = GetSubShapeByNode( nodeOnFace, GetMeshDS() ); TopoDS_Shape f = GetSubShapeByNode( nodeOnFace, GetMeshDS() );
if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE ) if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE )
{ {
face = TopoDS::Face( f ); face = TopoDS::Face( f );
faceHlp.SetSubShape( face );
Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc); Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
bool isStraight[2]; bool isStraight[2];
for ( int is1 = 0; is1 < 2; ++is1 ) // move0 or move1 for ( int is1 = 0; is1 < 2; ++is1 ) // move0 or move1
{ {
TChainLink& link = is1 ? chain.back() : chain.front(); TChainLink& link = is1 ? chain.back() : chain.front();
gp_XY uvm = GetNodeUV( face, link->_mediumNode, nodeOnFace, &checkUV); gp_XY uvm = faceHlp.GetNodeUV( face, link->_mediumNode, nodeOnFace, &checkUV );
gp_XY uv1 = GetNodeUV( face, link->node1(), nodeOnFace, &checkUV); gp_XY uv1 = faceHlp.GetNodeUV( face, link->node1(), nodeOnFace, &checkUV );
gp_XY uv2 = GetNodeUV( face, link->node2(), nodeOnFace, &checkUV); gp_XY uv2 = faceHlp.GetNodeUV( face, link->node2(), nodeOnFace, &checkUV );
gp_XY uv12 = GetMiddleUV( surf, uv1, uv2); gp_XY uv12 = faceHlp.GetMiddleUV( surf, uv1, uv2 );
// uvMove = uvm - uv12 // uvMove = uvm - uv12
gp_XY uvMove = applyIn2D(surf, uvm, uv12, gp_XY_Subtracted, /*inPeriod=*/false); gp_XY uvMove = applyIn2D(surf, uvm, uv12, gp_XY_Subtracted, /*inPeriod=*/false);
( is1 ? move1 : move0 ).SetCoord( uvMove.X(), uvMove.Y(), 0 ); ( is1 ? move1 : move0 ).SetCoord( uvMove.X(), uvMove.Y(), 0 );
@ -4719,10 +4751,10 @@ void SMESH_MesherHelper::FixQuadraticElements(SMESH_ComputeErrorPtr& compError,
} }
// check if a chain is already fixed // check if a chain is already fixed
gp_XY uvm = GetNodeUV( face, linkOnFace->_mediumNode, 0, &checkUV); gp_XY uvm = faceHlp.GetNodeUV( face, linkOnFace->_mediumNode, 0, &checkUV );
gp_XY uv1 = GetNodeUV( face, linkOnFace->node1(), nodeOnFace, &checkUV); gp_XY uv1 = faceHlp.GetNodeUV( face, linkOnFace->node1(), nodeOnFace, &checkUV );
gp_XY uv2 = GetNodeUV( face, linkOnFace->node2(), nodeOnFace, &checkUV); gp_XY uv2 = faceHlp.GetNodeUV( face, linkOnFace->node2(), nodeOnFace, &checkUV );
gp_XY uv12 = GetMiddleUV( surf, uv1, uv2); gp_XY uv12 = faceHlp.GetMiddleUV( surf, uv1, uv2 );
if (( uvm - uv12 ).SquareModulus() > 1e-10 ) if (( uvm - uv12 ).SquareModulus() > 1e-10 )
{ {
MSG("Already fixed - ignore"); MSG("Already fixed - ignore");
@ -4767,8 +4799,8 @@ void SMESH_MesherHelper::FixQuadraticElements(SMESH_ComputeErrorPtr& compError,
else { else {
// compute 3D displacement by 2D one // compute 3D displacement by 2D one
Handle(Geom_Surface) s = BRep_Tool::Surface(face,loc); Handle(Geom_Surface) s = BRep_Tool::Surface(face,loc);
gp_XY oldUV = GetNodeUV( face, (*link1)->_mediumNode, 0, &checkUV); gp_XY oldUV = faceHlp.GetNodeUV( face, (*link1)->_mediumNode, 0, &checkUV );
gp_XY newUV = applyIn2D( s, oldUV, gp_XY( move.X(),move.Y()), gp_XY_Added); gp_XY newUV = applyIn2D( s, oldUV, gp_XY( move.X(),move.Y()), gp_XY_Added );
gp_Pnt newPnt = s->Value( newUV.X(), newUV.Y()); gp_Pnt newPnt = s->Value( newUV.X(), newUV.Y());
move = gp_Vec( XYZ((*link1)->_mediumNode), newPnt.Transformed(loc) ); move = gp_Vec( XYZ((*link1)->_mediumNode), newPnt.Transformed(loc) );
if ( SMDS_FacePosition* nPos = if ( SMDS_FacePosition* nPos =
@ -4778,8 +4810,8 @@ void SMESH_MesherHelper::FixQuadraticElements(SMESH_ComputeErrorPtr& compError,
if ( (XYZ((*link1)->node1()) - XYZ((*link1)->node2())).SquareModulus() < if ( (XYZ((*link1)->node1()) - XYZ((*link1)->node2())).SquareModulus() <
move.SquareMagnitude()) move.SquareMagnitude())
{ {
gp_XY uv0 = GetNodeUV( face, (*link0)->_mediumNode, 0, &checkUV); gp_XY uv0 = faceHlp.GetNodeUV( face, (*link0)->_mediumNode, 0, &checkUV );
gp_XY uv2 = GetNodeUV( face, (*link2)->_mediumNode, 0, &checkUV); gp_XY uv2 = faceHlp.GetNodeUV( face, (*link2)->_mediumNode, 0, &checkUV );
MSG( "TOO LONG MOVE \t" << MSG( "TOO LONG MOVE \t" <<
"uv0: "<<uv0.X()<<", "<<uv0.Y()<<" \t" << "uv0: "<<uv0.X()<<", "<<uv0.Y()<<" \t" <<
"uv2: "<<uv2.X()<<", "<<uv2.Y()<<" \t" << "uv2: "<<uv2.X()<<", "<<uv2.Y()<<" \t" <<

View File

@ -671,11 +671,14 @@ public:
* \param uv2 - UV within a face * \param uv2 - UV within a face
* \retval gp_Pnt2d - selected UV * \retval gp_Pnt2d - selected UV
*/ */
gp_Pnt2d GetUVOnSeam( const gp_Pnt2d& uv1, const gp_Pnt2d& uv2 ) const; gp_Pnt2d getUVOnSeam( const gp_Pnt2d& uv1, const gp_Pnt2d& uv2 ) const;
const SMDS_MeshNode* getMediumNodeOnComposedWire(const SMDS_MeshNode* n1, const SMDS_MeshNode* getMediumNodeOnComposedWire(const SMDS_MeshNode* n1,
const SMDS_MeshNode* n2, const SMDS_MeshNode* n2,
bool force3d); bool force3d);
double getFaceMaxTol( const TopoDS_Shape& face ) const;
private: private:
// Forbiden copy constructor // Forbiden copy constructor
@ -710,9 +713,11 @@ public:
double myPar1[2], myPar2[2]; // U and V bounds of a closed periodic surface double myPar1[2], myPar2[2]; // U and V bounds of a closed periodic surface
int myParIndex; // bounds' index (1-U, 2-V, 3-both) int myParIndex; // bounds' index (1-U, 2-V, 3-both)
typedef std::map< int, GeomAPI_ProjectPointOnSurf* > TID2ProjectorOnSurf; std::map< int, double > myFaceMaxTol;
TID2ProjectorOnSurf myFace2Projector;
typedef std::map< int, GeomAPI_ProjectPointOnSurf* > TID2ProjectorOnSurf;
typedef std::map< int, GeomAPI_ProjectPointOnCurve* > TID2ProjectorOnCurve; typedef std::map< int, GeomAPI_ProjectPointOnCurve* > TID2ProjectorOnCurve;
TID2ProjectorOnSurf myFace2Projector;
TID2ProjectorOnCurve myEdge2Projector; TID2ProjectorOnCurve myEdge2Projector;
TopoDS_Shape myShape; TopoDS_Shape myShape;

View File

@ -213,7 +213,7 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
existingNodes.insert( n ); existingNodes.insert( n );
} }
// get EDGESs and their ids and get existing nodes on EDGEs // get EDGEs and their ids and get existing nodes on EDGEs
vector< TopoDS_Edge > edges; vector< TopoDS_Edge > edges;
for ( exp.Init( theShape, TopAbs_EDGE ); exp.More(); exp.Next() ) for ( exp.Init( theShape, TopAbs_EDGE ); exp.More(); exp.Next() )
{ {