smesh/src/StdMeshers/StdMeshers_Import_1D2D.cxx

645 lines
22 KiB
C++
Raw Normal View History

2010-11-25 17:44:43 +05:00
// Copyright (C) 2007-2010 CEA/DEN, EDF R&D, OPEN CASCADE
//
// Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
// CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
//
// This library is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 2.1 of the License.
//
// This library is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
// Lesser General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public
// License along with this library; if not, write to the Free Software
// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
//
// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
//
// SMESH SMESH : implementaion of SMESH idl descriptions
// File : StdMeshers_Import_1D2D.cxx
// Module : SMESH
//
#include "StdMeshers_Import_1D2D.hxx"
#include "StdMeshers_Import_1D.hxx"
#include "StdMeshers_ImportSource.hxx"
#include "SMDS_MeshElement.hxx"
#include "SMDS_MeshNode.hxx"
#include "SMESHDS_Group.hxx"
#include "SMESHDS_Mesh.hxx"
#include "SMESH_Comment.hxx"
#include "SMESH_Gen.hxx"
#include "SMESH_Group.hxx"
#include "SMESH_Mesh.hxx"
#include "SMESH_MesherHelper.hxx"
#include "SMESH_subMesh.hxx"
#include "Utils_SALOME_Exception.hxx"
#include "utilities.h"
#include <BRep_Builder.hxx>
#include <BRep_Tool.hxx>
#include <TopExp.hxx>
#include <TopExp_Explorer.hxx>
#include <TopoDS.hxx>
#include <TopoDS_Compound.hxx>
#include <TopoDS_Edge.hxx>
#include <TopoDS_Vertex.hxx>
#include <numeric>
using namespace std;
//=============================================================================
/*!
* Creates StdMeshers_Import_1D2D
*/
//=============================================================================
StdMeshers_Import_1D2D::StdMeshers_Import_1D2D(int hypId, int studyId, SMESH_Gen * gen)
:SMESH_2D_Algo(hypId, studyId, gen), _sourceHyp(0)
{
MESSAGE("StdMeshers_Import_1D2D::StdMeshers_Import_1D2D");
_name = "Import_1D2D";
_shapeType = (1 << TopAbs_FACE);
_compatibleHypothesis.push_back("ImportSource2D");
_requireDescretBoundary = false;
}
//=============================================================================
/*!
* Check presence of a hypothesis
*/
//=============================================================================
bool StdMeshers_Import_1D2D::CheckHypothesis
(SMESH_Mesh& aMesh,
const TopoDS_Shape& aShape,
SMESH_Hypothesis::Hypothesis_Status& aStatus)
{
_sourceHyp = 0;
const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
if ( hyps.size() == 0 )
{
aStatus = SMESH_Hypothesis::HYP_MISSING;
return false; // can't work with no hypothesis
}
if ( hyps.size() > 1 )
{
aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
return false;
}
const SMESHDS_Hypothesis *theHyp = hyps.front();
string hypName = theHyp->GetName();
if (hypName == _compatibleHypothesis.front())
{
_sourceHyp = (StdMeshers_ImportSource1D *)theHyp;
aStatus = SMESH_Hypothesis::HYP_OK;
return true;
}
aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
return true;
}
namespace
{
/*!
* \brief OrientedLink additionally storing a medium node
*/
struct TLink : public SMESH_OrientedLink
{
const SMDS_MeshNode* _medium;
TLink( const SMDS_MeshNode* n1,
const SMDS_MeshNode* n2,
const SMDS_MeshNode* medium=0)
: SMESH_OrientedLink( n1,n2 ), _medium( medium ) {}
};
}
//=============================================================================
/*!
* Import elements from the other mesh
*/
//=============================================================================
bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape & theShape)
{
if ( !_sourceHyp ) return false;
const vector<SMESH_Group*>& srcGroups = _sourceHyp->GetGroups();
if ( srcGroups.empty() )
return error("Invalid source groups");
SMESH_MesherHelper helper(theMesh);
helper.SetSubShape(theShape);
SMESHDS_Mesh* tgtMesh = theMesh.GetMeshDS();
const TopoDS_Face& geomFace = TopoDS::Face( theShape );
const double faceTol = helper.MaxTolerance( geomFace );
const int shapeID = tgtMesh->ShapeToIndex( geomFace );
const bool toCheckOri = (helper.NbAncestors( geomFace, theMesh, TopAbs_SOLID ) == 1 );
Handle(Geom_Surface) surface = BRep_Tool::Surface( geomFace );
2010-12-24 13:18:34 +05:00
const bool reverse =
( helper.GetSubShapeOri( tgtMesh->ShapeToMesh(), geomFace) == TopAbs_REVERSED );
2010-11-25 17:44:43 +05:00
gp_Pnt p; gp_Vec du, dv;
set<int> subShapeIDs;
subShapeIDs.insert( shapeID );
// get nodes on vertices
list < SMESH_TNodeXYZ > vertexNodes;
list < SMESH_TNodeXYZ >::iterator vNIt;
2010-11-25 17:44:43 +05:00
TopExp_Explorer exp( theShape, TopAbs_VERTEX );
for ( ; exp.More(); exp.Next() )
{
const TopoDS_Vertex& v = TopoDS::Vertex( exp.Current() );
if ( !subShapeIDs.insert( tgtMesh->ShapeToIndex( v )).second )
continue;
const SMDS_MeshNode* n = SMESH_Algo::VertexNode( v, tgtMesh );
if ( !n )
{
_gen->Compute(theMesh,v,/*anUpward=*/true);
n = SMESH_Algo::VertexNode( v, tgtMesh );
if ( !n ) return false; // very strange
}
vertexNodes.push_back( SMESH_TNodeXYZ( n ));
2010-11-25 17:44:43 +05:00
}
// to count now many times a link between nodes encounters
map<TLink, int> linkCount;
map<TLink, int>::iterator link2Nb;
// =========================
// Import faces from groups
// =========================
StdMeshers_Import_1D::TNodeNodeMap* n2n;
StdMeshers_Import_1D::TElemElemMap* e2e;
vector<const SMDS_MeshNode*> newNodes;
for ( int iG = 0; iG < srcGroups.size(); ++iG )
{
const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
const int meshID = srcGroup->GetMesh()->GetPersistentId();
const SMESH_Mesh* srcMesh = GetMeshByPersistentID( meshID );
if ( !srcMesh ) continue;
StdMeshers_Import_1D::getMaps( srcMesh, &theMesh, n2n, e2e );
SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
SMDS_MeshNode *tmpNode = helper.AddNode(0,0,0);
2010-11-25 17:44:43 +05:00
gp_XY uv;
while ( srcElems->more() ) // loop on group contents
{
const SMDS_MeshElement* face = srcElems->next();
// find or create nodes of a new face
newNodes.resize( face->NbNodes() );
newNodes.back() = 0;
int nbCreatedNodes = 0;
SMDS_MeshElement::iterator node = face->begin_nodes();
for ( unsigned i = 0; i < newNodes.size(); ++i, ++node )
{
TNodeNodeMap::iterator n2nIt = n2n->insert( make_pair( *node, (SMDS_MeshNode*)0 )).first;
if ( n2nIt->second )
{
if ( !subShapeIDs.count( n2nIt->second->getshapeId() ))
2010-11-25 17:44:43 +05:00
break;
}
else
{
// find an existing vertex node
for ( vNIt = vertexNodes.begin(); vNIt != vertexNodes.end(); ++vNIt)
if ( vNIt->SquareDistance( *node ) < 10 * faceTol * faceTol)
{
(*n2nIt).second = vNIt->_node;
vertexNodes.erase( vNIt );
break;
}
}
if ( !n2nIt->second )
{
// find out if node lies on theShape
tmpNode->setXYZ( (*node)->X(), (*node)->Y(), (*node)->Z());
if ( helper.CheckNodeUV( geomFace, tmpNode, uv, 10 * faceTol, /*force=*/true ))
2010-11-25 17:44:43 +05:00
{
SMDS_MeshNode* newNode = tgtMesh->AddNode( (*node)->X(), (*node)->Y(), (*node)->Z());
n2nIt->second = newNode;
tgtMesh->SetNodeOnFace( newNode, shapeID, uv.X(), uv.Y() );
nbCreatedNodes++;
}
}
if ( !(newNodes[i] = n2nIt->second ))
break;
}
if ( !newNodes.back() )
continue; // not all nodes of the face lie on theShape
// try to find already created face
SMDS_MeshElement * newFace = 0;
if ( nbCreatedNodes == 0 &&
tgtMesh->FindElement(newNodes, SMDSAbs_Face, /*noMedium=*/false))
continue; // repeated face in source groups already created
// check future face orientation
if ( toCheckOri )
{
int iNode = -1;
gp_Vec geomNorm;
do
{
uv = helper.GetNodeUV( geomFace, newNodes[++iNode] );
surface->D1( uv.X(),uv.Y(), p, du,dv );
2010-12-24 13:18:34 +05:00
geomNorm = reverse ? dv^du : du^dv;
2010-11-25 17:44:43 +05:00
}
while ( geomNorm.SquareMagnitude() < 1e-6 && iNode+1 < face->NbCornerNodes());
int iNext = helper.WrapIndex( iNode+1, face->NbCornerNodes() );
int iPrev = helper.WrapIndex( iNode-1, face->NbCornerNodes() );
SMESH_TNodeXYZ prevNode( newNodes[iPrev] );
SMESH_TNodeXYZ curNode ( newNodes[iNode] );
SMESH_TNodeXYZ nextNode( newNodes[iNext] );
2010-11-25 17:44:43 +05:00
gp_Vec n1n0( prevNode - curNode);
gp_Vec n1n2( nextNode - curNode );
gp_Vec meshNorm = n1n2 ^ n1n0;
if ( geomNorm * meshNorm < 0 )
std::reverse( newNodes.begin(), newNodes.end() );
}
// make a new face
switch ( newNodes.size() )
{
case 3:
newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2] );
break;
case 4:
newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2], newNodes[3] );
break;
case 6:
newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2],
newNodes[3], newNodes[4], newNodes[5]);
break;
case 8:
newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2], newNodes[3],
newNodes[4], newNodes[5], newNodes[6], newNodes[7]);
break;
default: continue;
}
tgtMesh->SetMeshElementOnShape( newFace, shapeID );
e2e->insert( make_pair( face, newFace ));
// collect links
int nbNodes = face->NbCornerNodes();
const SMDS_MeshNode* medium = 0;
for ( int i = 0; i < nbNodes; ++i )
{
const SMDS_MeshNode* n1 = newNodes[i];
const SMDS_MeshNode* n2 = newNodes[ (i+1)%nbNodes ];
if ( newFace->IsQuadratic() )
medium = newNodes[i+nbNodes];
link2Nb = linkCount.insert( make_pair( TLink( n1, n2, medium ), 0)).first;
++link2Nb->second;
}
}
helper.GetMeshDS()->RemoveNode(tmpNode);
2010-11-25 17:44:43 +05:00
}
// ==========================================================
// Put nodes on geom edges and create edges on them;
// check if the whole geom face is covered by imported faces
// ==========================================================
vector< TopoDS_Edge > edges;
for ( exp.Init( theShape, TopAbs_EDGE ); exp.More(); exp.Next() )
if ( subShapeIDs.insert( tgtMesh->ShapeToIndex( exp.Current() )).second )
edges.push_back( TopoDS::Edge( exp.Current() ));
bool isFaceMeshed = false;
if ( SMESHDS_SubMesh* tgtSM = tgtMesh->MeshElements( theShape ))
{
// the imported mesh is valid if all external links (encountered once)
// lie on geom edges
subShapeIDs.erase( shapeID ); // to contain edges and vertices only
double u, f, l;
for ( link2Nb = linkCount.begin(); link2Nb != linkCount.end(); ++link2Nb)
{
const TLink& link = (*link2Nb).first;
int nbFaces = link2Nb->second;
if ( nbFaces == 1 )
{
2010-12-24 13:18:34 +05:00
// check if a not shared link lie on face boundary
2010-11-25 17:44:43 +05:00
bool nodesOnBoundary = true;
list< TopoDS_Shape > bndShapes;
for ( int is1stN = 0; is1stN < 2 && nodesOnBoundary; ++is1stN )
{
const SMDS_MeshNode* n = is1stN ? link.node1() : link.node2();
if ( !subShapeIDs.count( n->getshapeId() ))
2010-11-25 17:44:43 +05:00
{
for ( unsigned iE = 0; iE < edges.size(); ++iE )
if ( helper.CheckNodeU( edges[iE], n, u, 10 * faceTol, /*force=*/true ))
{
BRep_Tool::Range(edges[iE],f,l);
if ( Abs(u-f) < 2 * faceTol || Abs(u-l) < 2 * faceTol )
// duplicated node on vertex
return error("Source elements overlap one another");
2010-12-24 13:18:34 +05:00
tgtSM->RemoveNode( n, /*isNodeDeleted=*/false );
2010-11-25 17:44:43 +05:00
tgtMesh->SetNodeOnEdge( (SMDS_MeshNode*)n, edges[iE], u );
break;
}
nodesOnBoundary = subShapeIDs.count( n->getshapeId());
2010-11-25 17:44:43 +05:00
}
if ( nodesOnBoundary )
{
TopoDS_Shape s = helper.GetSubShapeByNode( n, tgtMesh );
if ( s.ShapeType() == TopAbs_VERTEX )
bndShapes.push_front( s ); // vertex first
else
bndShapes.push_back( s ); // edges last
}
}
if ( !nodesOnBoundary )
2010-12-24 13:18:34 +05:00
break; // error: free internal link
2010-11-25 17:44:43 +05:00
if ( bndShapes.front().ShapeType() == TopAbs_EDGE &&
bndShapes.front() != bndShapes.back() )
2010-12-24 13:18:34 +05:00
break; // error: link nodes on different geom edges
2010-11-25 17:44:43 +05:00
// find geom edge the link is on
if ( bndShapes.back().ShapeType() != TopAbs_EDGE )
{
// find geom edge by two vertices
TopoDS_Shape geomEdge;
PShapeIteratorPtr edgeIt = helper.GetAncestors( bndShapes.back(), theMesh, TopAbs_EDGE );
while ( edgeIt->more() )
{
geomEdge = *(edgeIt->next());
if ( !helper.IsSubShape( bndShapes.front(), geomEdge ))
geomEdge.Nullify();
}
if ( geomEdge.IsNull() )
2010-12-24 13:18:34 +05:00
break; // vertices belong to different edges -> error: free internal link
2010-11-25 17:44:43 +05:00
bndShapes.push_back( geomEdge );
}
// create an edge if not yet exists
newNodes.resize(2);
newNodes[0] = link.node1(), newNodes[1] = link.node2();
const SMDS_MeshElement* edge = tgtMesh->FindElement( newNodes, SMDSAbs_Edge );
if ( edge ) continue;
if ( link._reversed ) std::swap( newNodes[0], newNodes[1] );
if ( link._medium )
{
newNodes.push_back( link._medium );
edge = tgtMesh->AddEdge( newNodes[0], newNodes[1], newNodes[2] );
TopoDS_Edge geomEdge = TopoDS::Edge(bndShapes.back());
helper.CheckNodeU( geomEdge, link._medium, u, 10*faceTol, /*force=*/true );
2010-12-24 13:18:34 +05:00
tgtSM->RemoveNode( link._medium, /*isNodeDeleted=*/false );
2010-11-25 17:44:43 +05:00
tgtMesh->SetNodeOnEdge( (SMDS_MeshNode*)link._medium, geomEdge, u );
}
else
{
edge = tgtMesh->AddEdge( newNodes[0], newNodes[1]);
}
if ( !edge )
return false;
tgtMesh->SetMeshElementOnShape( edge, bndShapes.back() );
}
else if ( nbFaces > 2 )
{
return error( "Non-manifold source mesh");
}
}
isFaceMeshed = ( link2Nb == linkCount.end() && !linkCount.empty());
if ( isFaceMeshed )
{
// check that source faces do not overlap:
// there must be only two edges sharing each vertex and bound to sub-edges of theShape
SMESH_MeshEditor editor( &theMesh );
set<int>::iterator subID = subShapeIDs.begin();
for ( ; subID != subShapeIDs.end(); ++subID )
{
const TopoDS_Shape& s = tgtMesh->IndexToShape( *subID );
if ( s.ShapeType() != TopAbs_VERTEX ) continue;
const SMDS_MeshNode* n = SMESH_Algo::VertexNode( TopoDS::Vertex(s), tgtMesh );
SMDS_ElemIteratorPtr eIt = n->GetInverseElementIterator(SMDSAbs_Edge);
int nbEdges = 0;
while ( eIt->more() )
{
const SMDS_MeshElement* edge = eIt->next();
int sId = editor.FindShape( edge );
nbEdges += subShapeIDs.count( sId );
}
if ( nbEdges < 2 )
return false; // weird
if ( nbEdges > 2 )
return error( "Source elements overlap one another");
}
}
}
if ( !isFaceMeshed )
return error( "Source elements don't cover totally the geometrical face" );
// notify sub-meshes of edges on computation
for ( unsigned iE = 0; iE < edges.size(); ++iE )
theMesh.GetSubMesh( edges[iE] )->ComputeStateEngine(SMESH_subMesh::CHECK_COMPUTE_STATE);
// ============
// Copy meshes
// ============
vector<SMESH_Mesh*> srcMeshes = _sourceHyp->GetSourceMeshes();
for ( unsigned i = 0; i < srcMeshes.size(); ++i )
StdMeshers_Import_1D::importMesh( srcMeshes[i], theMesh, _sourceHyp, theShape );
return true;
}
//=============================================================================
/*!
* \brief Set needed event listeners and create a submesh for a copied mesh
*
* This method is called only if a submesh has HYP_OK algo_state.
*/
//=============================================================================
void StdMeshers_Import_1D2D::SetEventListener(SMESH_subMesh* subMesh)
{
if ( !_sourceHyp )
{
const TopoDS_Shape& tgtShape = subMesh->GetSubShape();
SMESH_Mesh* tgtMesh = subMesh->GetFather();
Hypothesis_Status aStatus;
CheckHypothesis( *tgtMesh, tgtShape, aStatus );
}
StdMeshers_Import_1D::setEventListener( subMesh, _sourceHyp );
}
void StdMeshers_Import_1D2D::SubmeshRestored(SMESH_subMesh* subMesh)
{
SetEventListener(subMesh);
}
//=============================================================================
/*!
* Predict nb of mesh entities created by Compute()
*/
//=============================================================================
bool StdMeshers_Import_1D2D::Evaluate(SMESH_Mesh & theMesh,
const TopoDS_Shape & theShape,
MapShapeNbElems& aResMap)
{
if ( !_sourceHyp ) return false;
const vector<SMESH_Group*>& srcGroups = _sourceHyp->GetGroups();
if ( srcGroups.empty() )
return error("Invalid source groups");
vector<int> aVec(SMDSEntity_Last,0);
bool toCopyMesh, toCopyGroups;
_sourceHyp->GetCopySourceMesh(toCopyMesh, toCopyGroups);
if ( toCopyMesh ) // the whole mesh is copied
{
vector<SMESH_Mesh*> srcMeshes = _sourceHyp->GetSourceMeshes();
for ( unsigned i = 0; i < srcMeshes.size(); ++i )
{
SMESH_subMesh* sm = StdMeshers_Import_1D::getSubMeshOfCopiedMesh( theMesh, *srcMeshes[i]);
if ( !sm || aResMap.count( sm )) continue; // already counted
const SMDS_MeshInfo& aMeshInfo = srcMeshes[i]->GetMeshDS()->GetMeshInfo();
for (int i = 0; i < SMDSEntity_Last; i++)
aVec[i] = aMeshInfo.NbEntities((SMDSAbs_EntityType)i);
}
}
else
{
// std-like iterator used to get coordinates of nodes of mesh element
typedef SMDS_StdIterator< SMESH_TNodeXYZ, SMDS_ElemIteratorPtr > TXyzIterator;
2010-11-25 17:44:43 +05:00
SMESH_MesherHelper helper(theMesh);
helper.SetSubShape(theShape);
const TopoDS_Face& geomFace = TopoDS::Face( theShape );
const double faceTol = helper.MaxTolerance( geomFace );
// take into account nodes on vertices
TopExp_Explorer exp( theShape, TopAbs_VERTEX );
for ( ; exp.More(); exp.Next() )
theMesh.GetSubMesh( exp.Current())->Evaluate( aResMap );
// to count now many times a link between nodes encounters,
// negative nb additionally means that a link is quadratic
map<SMESH_TLink, int> linkCount;
map<SMESH_TLink, int>::iterator link2Nb;
// count faces and nodes imported from groups
set<const SMDS_MeshNode* > allNodes;
gp_XY uv;
for ( int iG = 0; iG < srcGroups.size(); ++iG )
{
const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
SMDS_MeshNode *tmpNode =helper.AddNode(0,0,0);
2010-11-25 17:44:43 +05:00
while ( srcElems->more() ) // loop on group contents
{
const SMDS_MeshElement* face = srcElems->next();
// find out if face is located on geomEdge by projecting
// a gravity center of face to geomFace
gp_XYZ gc(0,0,0);
gc = accumulate( TXyzIterator(face->nodesIterator()), TXyzIterator(), gc)/face->NbNodes();
tmpNode->setXYZ( gc.X(), gc.Y(), gc.Z());
if ( helper.CheckNodeUV( geomFace, tmpNode, uv, 10 * faceTol, /*force=*/true ))
2010-11-25 17:44:43 +05:00
{
++aVec[ face->GetEntityType() ];
// collect links
int nbConers = face->NbCornerNodes();
for ( int i = 0; i < face->NbNodes(); ++i )
{
const SMDS_MeshNode* n1 = face->GetNode(i);
allNodes.insert( n1 );
if ( i < nbConers )
{
const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbConers );
link2Nb = linkCount.insert( make_pair( SMESH_TLink( n1, n2 ), 0)).first;
if ( (*link2Nb).second )
link2Nb->second += (link2Nb->second < 0 ) ? -1 : 1;
else
link2Nb->second += ( face->IsQuadratic() ) ? -1 : 1;
}
}
}
}
helper.GetMeshDS()->RemoveNode(tmpNode);
2010-11-25 17:44:43 +05:00
}
int nbNodes = allNodes.size();
allNodes.clear();
// count nodes and edges on geom edges
double u;
for ( exp.Init(theShape, TopAbs_EDGE); exp.More(); exp.Next() )
{
TopoDS_Edge geomEdge = TopoDS::Edge( exp.Current() );
SMESH_subMesh* sm = theMesh.GetSubMesh( geomEdge );
vector<int>& edgeVec = aResMap[sm];
if ( edgeVec.empty() )
{
edgeVec.resize(SMDSEntity_Last,0);
for ( link2Nb = linkCount.begin(); link2Nb != linkCount.end(); )
{
const SMESH_TLink& link = (*link2Nb).first;
int nbFacesOfLink = Abs( link2Nb->second );
bool eraseLink = ( nbFacesOfLink != 1 );
if ( nbFacesOfLink == 1 )
{
if ( helper.CheckNodeU( geomEdge, link.node1(), u, 10*faceTol, /*force=*/true )&&
helper.CheckNodeU( geomEdge, link.node2(), u, 10*faceTol, /*force=*/true ))
{
bool isQuadratic = ( link2Nb->second < 0 );
++edgeVec[ isQuadratic ? SMDSEntity_Quad_Edge : SMDSEntity_Edge ];
++edgeVec[ SMDSEntity_Node ];
--nbNodes;
eraseLink = true;
}
}
if ( eraseLink )
linkCount.erase(link2Nb++);
else
link2Nb++;
}
if ( edgeVec[ SMDSEntity_Node] > 0 )
--edgeVec[ SMDSEntity_Node ]; // for one node on vertex
}
else if ( !helper.IsSeamShape( geomEdge ) ||
geomEdge.Orientation() == TopAbs_FORWARD )
{
nbNodes -= 1+edgeVec[ SMDSEntity_Node ];
}
}
aVec[SMDSEntity_Node] = nbNodes;
}
SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
aResMap.insert(make_pair(sm,aVec));
return true;
}