#include <BRepTools_WireExplorer.hxx>
#include <BRep_Tool.hxx>
#include <Geom2d_Curve.hxx>
+#include <GeomAPI_ProjectPointOnCurve.hxx>
#include <GeomAPI_ProjectPointOnSurf.hxx>
#include <Geom_Curve.hxx>
#include <Geom_Surface.hxx>
gp_XYZ XYZ(const SMDS_MeshNode* n) { return gp_XYZ(n->X(), n->Y(), n->Z()); }
+ enum { U_periodic = 1, V_periodic = 2 };
}
//================================================================================
//================================================================================
SMESH_MesherHelper::SMESH_MesherHelper(SMESH_Mesh& theMesh)
- : myMesh(&theMesh), myShapeID(0), myCreateQuadratic(false)
+ : myParIndex(0), myMesh(&theMesh), myShapeID(0), myCreateQuadratic(false)
{
+ myPar1[0] = myPar2[0] = myPar1[1] = myPar2[1] = 0;
mySetElemOnShape = ( ! myMesh->HasShapeToMesh() );
}
//=======================================================================
-//function : CheckShape
-//purpose :
+//function : IsQuadraticSubMesh
+//purpose : Check submesh for given shape: if all elements on this shape
+// are quadratic, quadratic elements will be created.
+// Also fill myTLinkNodeMap
//=======================================================================
bool SMESH_MesherHelper::IsQuadraticSubMesh(const TopoDS_Shape& aSh)
return myCreateQuadratic;
}
-//================================================================================
-/*!
- * \brief Set geomerty to make elements on
- * \param aSh - geomertic shape
- */
-//================================================================================
+//=======================================================================
+//function : SetSubShape
+//purpose : Set geomerty to make elements on
+//=======================================================================
void SMESH_MesherHelper::SetSubShape(const int aShID)
{
SetSubShape( TopoDS_Shape() );
}
-//================================================================================
-/*!
- * \brief Set geomerty to make elements on
- * \param aSh - geomertic shape
- */
-//================================================================================
+//=======================================================================
+//function : SetSubShape
+//purpose : Set geomerty to create elements on
+//=======================================================================
void SMESH_MesherHelper::SetSubShape(const TopoDS_Shape& aSh)
{
}
SMESHDS_Mesh* meshDS = GetMeshDS();
myShapeID = meshDS->ShapeToIndex(aSh);
+ myParIndex = 0;
// treatment of periodic faces
for ( TopExp_Explorer eF( aSh, TopAbs_FACE ); eF.More(); eF.Next() )
const TopoDS_Edge& edge = TopoDS::Edge( exp.Current() );
if ( BRep_Tool::IsClosed( edge, face )) {
// initialize myPar1, myPar2 and myParIndex
- if ( mySeamShapeIds.empty() ) {
- gp_Pnt2d uv1, uv2;
- BRep_Tool::UVPoints( edge, face, uv1, uv2 );
- if ( Abs( uv1.Coord(1) - uv2.Coord(1) ) < Abs( uv1.Coord(2) - uv2.Coord(2) ))
- {
- myParIndex = 1; // U periodic
- myPar1 = surface.FirstUParameter();
- myPar2 = surface.LastUParameter();
- }
- else {
- myParIndex = 2; // V periodic
- myPar1 = surface.FirstVParameter();
- myPar2 = surface.LastVParameter();
- }
+ gp_Pnt2d uv1, uv2;
+ BRep_Tool::UVPoints( edge, face, uv1, uv2 );
+ if ( Abs( uv1.Coord(1) - uv2.Coord(1) ) < Abs( uv1.Coord(2) - uv2.Coord(2) ))
+ {
+ myParIndex |= U_periodic;
+ myPar1[0] = surface.FirstUParameter();
+ myPar2[0] = surface.LastUParameter();
+ }
+ else {
+ myParIndex |= V_periodic;
+ myPar1[1] = surface.FirstVParameter();
+ myPar2[1] = surface.LastVParameter();
}
// store seam shape indices, negative if shape encounters twice
int edgeID = meshDS->ShapeToIndex( edge );
}
}
-//================================================================================
- /*!
- * \brief Check if inFaceNode argument is necessary for call GetNodeUV(F,..)
- * \param F - the face
- * \retval bool - return true if the face is periodic
- */
-//================================================================================
+//=======================================================================
+//function : GetNodeUVneedInFaceNode
+//purpose : Check if inFaceNode argument is necessary for call GetNodeUV(F,..)
+// Return true if the face is periodic.
+// If F is Null, answer about subshape set through IsQuadraticSubMesh() or
+// * SetSubShape()
+//=======================================================================
bool SMESH_MesherHelper::GetNodeUVneedInFaceNode(const TopoDS_Face& F) const
{
}
//=======================================================================
-/*!
- * \brief Return support shape of a node
- * \param node - the node
- * \param meshDS - mesh DS
- * \retval TopoDS_Shape - found support shape
- */
+//function : GetSubShapeByNode
+//purpose : Return support shape of a node
//=======================================================================
TopoDS_Shape SMESH_MesherHelper::GetSubShapeByNode(const SMDS_MeshNode* node,
//=======================================================================
//function : AddTLinkNode
-//purpose :
+//purpose : add a link in my data structure
//=======================================================================
-/*!
- * Auxilary function for filling myTLinkNodeMap
- */
+
void SMESH_MesherHelper::AddTLinkNode(const SMDS_MeshNode* n1,
const SMDS_MeshNode* n2,
const SMDS_MeshNode* n12)
}
//=======================================================================
-/*!
- * \brief Select UV on either of 2 pcurves of a seam edge, closest to the given UV
- * \param uv1 - UV on the seam
- * \param uv2 - UV within a face
- * \retval gp_Pnt2d - selected UV
- */
+//function : GetUVOnSeam
+//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
{
- double p1 = uv1.Coord( myParIndex );
- double p2 = uv2.Coord( myParIndex );
- double p3 = ( Abs( p1 - myPar1 ) < Abs( p1 - myPar2 )) ? myPar2 : myPar1;
- if ( Abs( p2 - p1 ) > Abs( p2 - p3 ))
- p1 = p3;
gp_Pnt2d result = uv1;
- result.SetCoord( myParIndex, p1 );
+ for ( int i = U_periodic; i <= V_periodic ; ++i )
+ {
+ if ( myParIndex & i )
+ {
+ double p1 = uv1.Coord( i );
+ double dp1 = Abs( p1-myPar1[i-1]), dp2 = Abs( p1-myPar2[i-1]);
+ if ( myParIndex == i ||
+ dp1 < ( myPar2[i-1] - myPar2[i-1] ) / 100. ||
+ dp2 < ( myPar2[i-1] - myPar2[i-1] ) / 100. )
+ {
+ double p2 = uv2.Coord( i );
+ double p1Alt = ( dp1 < dp2 ) ? myPar2[i-1] : myPar1[i-1];
+ if ( Abs( p2 - p1 ) > Abs( p2 - p1Alt ))
+ result.SetCoord( i, p1Alt );
+ }
+ }
+ }
return result;
}
//=======================================================================
-/*!
- * \brief Return node UV on face
- * \param F - the face
- * \param n - the node
- * \param n2 - a node of element being created located inside a face
- * \param check - optional flag returing false if found UV are invalid
- * \retval gp_XY - resulting UV
- */
+//function : GetNodeUV
+//purpose : Return node UV on face
//=======================================================================
gp_XY SMESH_MesherHelper::GetNodeUV(const TopoDS_Face& F,
const SMDS_MeshNode* n2,
bool* check) const
{
- gp_Pnt2d uv( 1e100, 1e100 );
+ gp_Pnt2d uv( Precision::Infinite(), Precision::Infinite() );
const SMDS_PositionPtr Pos = n->GetPosition();
bool uvOK = false;
if(Pos->GetTypeOfPosition()==SMDS_TOP_FACE)
const SMDS_FacePosition* fpos =
static_cast<const SMDS_FacePosition*>(n->GetPosition().get());
uv.SetCoord(fpos->GetUParameter(),fpos->GetVParameter());
- uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), BRep_Tool::Tolerance( F ));
+ if ( check )
+ uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), BRep_Tool::Tolerance( F ));
}
else if(Pos->GetTypeOfPosition()==SMDS_TOP_EDGE)
{
TopoDS_Edge E = TopoDS::Edge(GetMeshDS()->IndexToShape(edgeID));
double f, l, u = epos->GetUParameter();
Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(E, F, f, l);
- if ( f < u && u < l )
+ bool validU = ( f < u && u < l );
+ if ( validU )
uv = C2d->Value( u );
else
uv.SetCoord(0.,0.);
- uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), BRep_Tool::Tolerance( E ));
+ if ( check || !validU )
+ uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), BRep_Tool::Tolerance( E ),/*force=*/ !validU );
// for a node on a seam edge select one of UVs on 2 pcurves
if ( n2 && IsSeamShape( edgeID ) )
}
}
else {
+ uvOK = false;
TopTools_ListIteratorOfListOfShape it( myMesh->GetAncestors( V ));
for ( ; it.More(); it.Next() ) {
if ( it.Value().ShapeType() == TopAbs_EDGE ) {
if ( !C2d.IsNull() ) {
double u = ( V == TopExp::FirstVertex( edge ) ) ? f : l;
uv = C2d->Value( u );
+ uvOK = true;
break;
}
}
}
//=======================================================================
-/*!
- * \brief Check and fix node UV on a face
- * \retval bool - false if UV is bad and could not be fixed
- */
+//function : CheckNodeUV
+//purpose : Check and fix node UV on a face
//=======================================================================
bool SMESH_MesherHelper::CheckNodeUV(const TopoDS_Face& F,
const SMDS_MeshNode* n,
gp_XY& uv,
- const double tol) const
+ const double tol,
+ const bool force) const
{
- if ( !myOkNodePosShapes.count( n->GetPosition()->GetShapeId() ))
+ if ( force || !myOkNodePosShapes.count( n->GetPosition()->GetShapeId() ))
{
// check that uv is correct
TopLoc_Location loc;
Handle(Geom_Surface) surface = BRep_Tool::Surface( F,loc );
gp_Pnt nodePnt = XYZ( n );
if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
- if ( nodePnt.Distance( surface->Value( uv.X(), uv.Y() )) > tol )
+ if ( Precision::IsInfinite( uv.X() ) ||
+ Precision::IsInfinite( uv.Y() ) ||
+ nodePnt.Distance( surface->Value( uv.X(), uv.Y() )) > tol )
{
// uv incorrect, project the node to surface
GeomAPI_ProjectPointOnSurf projector( nodePnt, surface, tol );
return true;
}
+namespace
+{
+ struct TMiddle
+ {
+ gp_XY operator()(const gp_XY& uv1, const gp_XY& uv2) const { return ( uv1 + uv2 ) / 2.; }
+ };
+ struct TAdd
+ {
+ gp_XY operator()(const gp_XY& uv1, const gp_XY& uv2) const { return ( uv1 + uv2 ); }
+ };
+ struct TSubtract
+ {
+ gp_XY operator()(const gp_XY& uv1, const gp_XY& uv2) const { return ( uv1 - uv2 ); }
+ };
+
+ //================================================================================
+ /*!
+ * \brief Perform given operation on two points in parametric space of given surface
+ * Example: gp_XY uvSum = applyXYFUN( surf, uv1, uv2, gp_XYFun(Added))
+ */
+ //================================================================================
+
+ template<typename FUNC>
+ gp_XY applyFunc(const Handle(Geom_Surface)& surface,
+ const gp_XY& uv1,
+ gp_XY uv2,
+ const bool resultInPeriod=true)
+ {
+ Standard_Boolean isUPeriodic = surface.IsNull() ? false : surface->IsUPeriodic();
+ Standard_Boolean isVPeriodic = surface.IsNull() ? false : surface->IsVPeriodic();
+ if ( !isUPeriodic && !isVPeriodic )
+ return FUNC()(uv1,uv2);
+
+ // move uv2 not far than half-period from uv1
+ if ( isUPeriodic )
+ uv2.SetX( uv2.X()+ShapeAnalysis::AdjustByPeriod(uv2.X(),uv1.X(),surface->UPeriod()) );
+ if ( isVPeriodic )
+ uv2.SetY( uv2.Y()+ShapeAnalysis::AdjustByPeriod(uv2.Y(),uv1.Y(),surface->VPeriod()) );
+
+ // execute operation
+ gp_XY res = FUNC()(uv1,uv2);
+
+ // move result within period
+ if ( resultInPeriod )
+ {
+ Standard_Real UF,UL,VF,VL;
+ surface->Bounds(UF,UL,VF,VL);
+ if ( isUPeriodic )
+ res.SetX( res.X() + ShapeAnalysis::AdjustToPeriod(res.X(),UF,UL));
+ if ( isVPeriodic )
+ res.SetY( res.Y() + ShapeAnalysis::AdjustToPeriod(res.Y(),VF,VL));
+ }
+
+ return res;
+ }
+}
//=======================================================================
-/*!
- * \brief Return middle UV taking in account surface period
- */
+//function : GetMiddleUV
+//purpose : Return middle UV taking in account surface period
//=======================================================================
gp_XY SMESH_MesherHelper::GetMiddleUV(const Handle(Geom_Surface)& surface,
const gp_XY& p1,
const gp_XY& p2)
{
- if ( surface.IsNull() )
- return 0.5 * ( p1 + p2 );
- //checking if surface is periodic
- Standard_Real UF,UL,VF,VL;
- surface->Bounds(UF,UL,VF,VL);
-
- Standard_Real u,v;
- Standard_Boolean isUPeriodic = surface->IsUPeriodic();
- if(isUPeriodic) {
- Standard_Real UPeriod = surface->UPeriod();
- Standard_Real p2x = p2.X()+ShapeAnalysis::AdjustByPeriod(p2.X(),p1.X(),UPeriod);
- Standard_Real pmid = (p1.X()+p2x)/2.;
- u = pmid+ShapeAnalysis::AdjustToPeriod(pmid,UF,UL);
- }
- else {
- u= (p1.X()+p2.X())/2.;
- }
- Standard_Boolean isVPeriodic = surface->IsVPeriodic();
- if(isVPeriodic) {
- Standard_Real VPeriod = surface->VPeriod();
- Standard_Real p2y = p2.Y()+ShapeAnalysis::AdjustByPeriod(p2.Y(),p1.Y(),VPeriod);
- Standard_Real pmid = (p1.Y()+p2y)/2.;
- v = pmid+ShapeAnalysis::AdjustToPeriod(pmid,VF,VL);
- }
- else {
- v = (p1.Y()+p2.Y())/2.;
- }
- return gp_XY( u,v );
+ return applyFunc<TMiddle>( surface, p1, p2 );
}
//=======================================================================
-/*!
- * \brief Return node U on edge
- * \param E - the Edge
- * \param n - the node
- * \retval double - resulting U
- *
- * Auxilary function called form GetMediumNode()
- */
+//function : GetNodeU
+//purpose : Return node U on edge
//=======================================================================
double SMESH_MesherHelper::GetNodeU(const TopoDS_Edge& E,
const TopoDS_Vertex& V = TopoDS::Vertex(meshDS->IndexToShape(vertexID));
param = BRep_Tool::Parameter( V, E );
}
+ if ( check )
+ *check = CheckNodeU( E, n, param, BRep_Tool::Tolerance( E ));
return param;
}
-//================================================================================
-/*!
- * \brief Return existing or create new medium nodes between given ones
- * \param force3d - if true, new node is the middle of n1 and n2,
- * else is located on geom face or geom edge
- */
-//================================================================================
+//=======================================================================
+//function : CheckNodeU
+//purpose : Check and fix node U on an edge
+// Return false if U is bad and could not be fixed
+//=======================================================================
+
+bool SMESH_MesherHelper::CheckNodeU(const TopoDS_Edge& E,
+ const SMDS_MeshNode* n,
+ double& u,
+ const double tol,
+ const bool force) const
+{
+ if ( force || !myOkNodePosShapes.count( n->GetPosition()->GetShapeId() ))
+ {
+ // check that u is correct
+ TopLoc_Location loc; double f,l;
+ Handle(Geom_Curve) curve = BRep_Tool::Curve( E,loc,f,l );
+ if ( curve.IsNull() ) // degenerated edge
+ {
+ if ( u+tol < f || u-tol > l )
+ {
+ double r = Max( 0.5, 1 - tol*n->GetID()); // to get a unique u on edge
+ u = f*r + l*(1-r);
+ }
+ }
+ else
+ {
+ gp_Pnt nodePnt = SMESH_MeshEditor::TNodeXYZ( n );
+ if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
+ if ( nodePnt.Distance( curve->Value( u )) > tol )
+ {
+ // u incorrect, project the node to the curve
+ GeomAPI_ProjectPointOnCurve projector( nodePnt, curve, f, l );
+ if ( projector.NbPoints() < 1 )
+ {
+ MESSAGE( "SMESH_MesherHelper::CheckNodeU() failed to project" );
+ return false;
+ }
+ Quantity_Parameter U = projector.LowerDistanceParameter();
+ if ( nodePnt.Distance( curve->Value( U )) > tol )
+ {
+ MESSAGE( "SMESH_MesherHelper::CheckNodeU(), invalid projection" );
+ return false;
+ }
+ u = double( U );
+ }
+ else if ( fabs( u ) > numeric_limits<double>::min() )
+ {
+ ((SMESH_MesherHelper*) this)->myOkNodePosShapes.insert( n->GetPosition()->GetShapeId() );
+ }
+ }
+ }
+ return true;
+}
+
+//=======================================================================
+//function : GetMediumNode
+//purpose : Return existing or create new medium nodes between given ones
+//=======================================================================
const SMDS_MeshNode* SMESH_MesherHelper::GetMediumNode(const SMDS_MeshNode* n1,
const SMDS_MeshNode* n2,
- bool force3d)
+ bool force3d)
{
+ // Find existing node
+
SMESH_TLink link(n1,n2);
ItTLinkNode itLN = myTLinkNodeMap.find( link );
if ( itLN != myTLinkNodeMap.end() ) {
return (*itLN).second;
}
- // create medium node
+
+ // Create medium node
+
SMDS_MeshNode* n12;
SMESHDS_Mesh* meshDS = GetMeshDS();
+
+ // get type of shape for the new medium node
int faceID = -1, edgeID = -1;
const SMDS_PositionPtr Pos1 = n1->GetPosition();
const SMDS_PositionPtr Pos2 = n2->GetPosition();
edgeID = Pos2->GetShapeId();
}
}
+ // get positions of the given nodes on shapes
+ TopoDS_Edge E; double u [2];
+ TopoDS_Face F; gp_XY uv[2];
+ bool uvOK[2] = { false, false };
+ TopAbs_ShapeEnum shapeType = myShape.IsNull() ? TopAbs_SHAPE : myShape.ShapeType();
+ if ( faceID>0 || shapeType == TopAbs_FACE)
+ {
+ if( myShape.IsNull() )
+ F = TopoDS::Face(meshDS->IndexToShape(faceID));
+ else {
+ F = TopoDS::Face(myShape);
+ faceID = myShapeID;
+ }
+ uv[0] = GetNodeUV(F,n1,n2, force3d ? 0 : &uvOK[0]);
+ uv[1] = GetNodeUV(F,n2,n1, force3d ? 0 : &uvOK[1]);
+ }
+ else if (edgeID>0 || shapeType == TopAbs_EDGE)
+ {
+ if( myShape.IsNull() )
+ E = TopoDS::Edge(meshDS->IndexToShape(edgeID));
+ else {
+ E = TopoDS::Edge(myShape);
+ edgeID = myShapeID;
+ }
+ u[0] = GetNodeU(E,n1, force3d ? 0 : &uvOK[0]);
+ u[1] = GetNodeU(E,n2, force3d ? 0 : &uvOK[1]);
+ }
if(!force3d)
{
// we try to create medium node using UV parameters of
// nodes, else - medium between corresponding 3d points
-
- TopAbs_ShapeEnum shapeType = myShape.IsNull() ? TopAbs_SHAPE : myShape.ShapeType();
- if(faceID>0 || shapeType == TopAbs_FACE) {
- // obtaining a face and 2d points for nodes
- TopoDS_Face F;
- if( myShape.IsNull() )
- F = TopoDS::Face(meshDS->IndexToShape(faceID));
- else {
- F = TopoDS::Face(myShape);
- faceID = myShapeID;
- }
- bool uvOK1, uvOK2;
- gp_XY p1 = GetNodeUV(F,n1,n2, &uvOK1);
- gp_XY p2 = GetNodeUV(F,n2,n1, &uvOK2);
-
- if ( uvOK1 && uvOK2 )
+ if( ! F.IsNull() )
+ {
+ if ( uvOK[0] && uvOK[1] )
{
if ( IsDegenShape( Pos1->GetShapeId() ))
- p1.SetCoord( myParIndex, p2.Coord( myParIndex ));
+ if ( myParIndex & U_periodic ) uv[0].SetCoord( 1, uv[1].Coord( 1 ));
+ else uv[0].SetCoord( 2, uv[1].Coord( 2 ));
else if ( IsDegenShape( Pos2->GetShapeId() ))
- p2.SetCoord( myParIndex, p1.Coord( myParIndex ));
+ if ( myParIndex & U_periodic ) uv[1].SetCoord( 1, uv[0].Coord( 1 ));
+ else uv[1].SetCoord( 2, uv[0].Coord( 2 ));
TopLoc_Location loc;
Handle(Geom_Surface) S = BRep_Tool::Surface(F,loc);
- gp_XY uv = GetMiddleUV( S, p1, p2 );
- gp_Pnt P = S->Value( uv.X(), uv.Y() ).Transformed(loc);
+ gp_XY UV = GetMiddleUV( S, uv[0], uv[1] );
+ gp_Pnt P = S->Value( UV.X(), UV.Y() ).Transformed(loc);
n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
- meshDS->SetNodeOnFace(n12, faceID, uv.X(), uv.Y());
+ meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y());
myTLinkNodeMap.insert(make_pair(link,n12));
return n12;
}
}
- if (edgeID>0 || shapeType == TopAbs_EDGE) {
-
- TopoDS_Edge E;
- if( myShape.IsNull() )
- E = TopoDS::Edge(meshDS->IndexToShape(edgeID));
- else {
- E = TopoDS::Edge(myShape);
- edgeID = myShapeID;
- }
-
- double p1 = GetNodeU(E,n1);
- double p2 = GetNodeU(E,n2);
-
+ else if ( !E.IsNull() )
+ {
double f,l;
Handle(Geom_Curve) C = BRep_Tool::Curve(E, f, l);
- if(!C.IsNull()) {
-
+ if(!C.IsNull())
+ {
Standard_Boolean isPeriodic = C->IsPeriodic();
- double u;
+ double U;
if(isPeriodic) {
Standard_Real Period = C->Period();
- Standard_Real p = p2+ShapeAnalysis::AdjustByPeriod(p2,p1,Period);
- Standard_Real pmid = (p1+p)/2.;
- u = pmid+ShapeAnalysis::AdjustToPeriod(pmid,C->FirstParameter(),C->LastParameter());
+ Standard_Real p = u[1]+ShapeAnalysis::AdjustByPeriod(u[1],u[0],Period);
+ Standard_Real pmid = (u[0]+p)/2.;
+ U = pmid+ShapeAnalysis::AdjustToPeriod(pmid,C->FirstParameter(),C->LastParameter());
}
else
- u = (p1+p2)/2.;
+ U = (u[0]+u[1])/2.;
- gp_Pnt P = C->Value( u );
+ gp_Pnt P = C->Value( U );
n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
- meshDS->SetNodeOnEdge(n12, edgeID, u);
+ meshDS->SetNodeOnEdge(n12, edgeID, U);
myTLinkNodeMap.insert(make_pair(link,n12));
return n12;
}
double y = ( n1->Y() + n2->Y() )/2.;
double z = ( n1->Z() + n2->Z() )/2.;
n12 = meshDS->AddNode(x,y,z);
- if(edgeID>0)
- meshDS->SetNodeOnEdge(n12, edgeID);
- else if(faceID>0)
- meshDS->SetNodeOnFace(n12, faceID);
+ if ( !F.IsNull() )
+ {
+ gp_XY UV = ( uv[0] + uv[1] ) / 2.;
+ CheckNodeUV( F, n12, UV, BRep_Tool::Tolerance( F ), /*force=*/true);
+ meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y() );
+ }
+ else if ( !E.IsNull() )
+ {
+ double U = ( u[0] + u[1] ) / 2.;
+ CheckNodeU( E, n12, U, BRep_Tool::Tolerance( E ), /*force=*/true);
+ meshDS->SetNodeOnEdge(n12, edgeID, U);
+ }
else
+ {
meshDS->SetNodeInVolume(n12, myShapeID);
+ }
myTLinkNodeMap.insert( make_pair( link, n12 ));
return n12;
}
//=======================================================================
-/*!
- * Creates a node
- */
+//function : AddNode
+//purpose : Creates a node
//=======================================================================
SMDS_MeshNode* SMESH_MesherHelper::AddNode(double x, double y, double z, int ID)
}
//=======================================================================
-/*!
- * Creates quadratic or linear edge
- */
+//function : AddEdge
+//purpose : Creates quadratic or linear edge
//=======================================================================
SMDS_MeshEdge* SMESH_MesherHelper::AddEdge(const SMDS_MeshNode* n1,
- const SMDS_MeshNode* n2,
- const int id,
- const bool force3d)
+ const SMDS_MeshNode* n2,
+ const int id,
+ const bool force3d)
{
SMESHDS_Mesh * meshDS = GetMeshDS();
}
//=======================================================================
-/*!
- * Creates quadratic or linear triangle
- */
+//function : AddFace
+//purpose : Creates quadratic or linear triangle
//=======================================================================
SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
}
//=======================================================================
-/*!
- * Creates quadratic or linear quadrangle
- */
+//function : AddFace
+//purpose : Creates quadratic or linear quadrangle
//=======================================================================
SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
const SMDS_MeshNode* n2,
const SMDS_MeshNode* n3,
const SMDS_MeshNode* n4,
- const int id,
- const bool force3d)
+ const int id,
+ const bool force3d)
{
SMESHDS_Mesh * meshDS = GetMeshDS();
SMDS_MeshFace* elem = 0;
}
//=======================================================================
-/*!
- * Creates quadratic or linear volume
- */
+//function : AddVolume
+//purpose : Creates quadratic or linear prism
//=======================================================================
SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
}
//=======================================================================
-/*!
- * Creates quadratic or linear volume
- */
+//function : AddVolume
+//purpose : Creates quadratic or linear tetrahedron
//=======================================================================
SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
}
//=======================================================================
-/*!
- * Creates quadratic or linear pyramid
- */
+//function : AddVolume
+//purpose : Creates quadratic or linear pyramid
//=======================================================================
SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
}
//=======================================================================
-/*!
- * Creates quadratic or linear hexahedron
- */
+//function : AddVolume
+//purpose : Creates quadratic or linear hexahedron
//=======================================================================
SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
}
//=======================================================================
-/*!
- * \brief Load nodes bound to face into a map of node columns
- * \param theParam2ColumnMap - map of node columns to fill
- * \param theFace - the face on which nodes are searched for
- * \param theBaseEdge - the edge nodes of which are columns' bases
- * \param theMesh - the mesh containing nodes
- * \retval bool - false if something is wrong
- *
- * The key of the map is a normalized parameter of each
- * base node on theBaseEdge.
- * This method works in supposition that nodes on the face
- * forms a rectangular grid and elements can be quardrangles or triangles
- */
+//function : LoadNodeColumns
+//purpose : Load nodes bound to face into a map of node columns
//=======================================================================
bool SMESH_MesherHelper::LoadNodeColumns(TParam2ColumnMap & theParam2ColumnMap,
}
//=======================================================================
-/*!
- * \brief Return number of unique ancestors of the shape
- */
+//function : NbAncestors
+//purpose : Return number of unique ancestors of the shape
//=======================================================================
int SMESH_MesherHelper::NbAncestors(const TopoDS_Shape& shape,
}
//=======================================================================
-/**
- * Check mesh without geometry for: if all elements on this shape are quadratic,
- * quadratic elements will be created.
- * Used then generated 3D mesh without geometry.
- */
+//function : GetSubShapeOri
+//purpose : Return orientation of sub-shape in the main shape
+//=======================================================================
+
+TopAbs_Orientation SMESH_MesherHelper::GetSubShapeOri(const TopoDS_Shape& shape,
+ const TopoDS_Shape& subShape)
+{
+ TopAbs_Orientation ori = TopAbs_Orientation(-1);
+ if ( !shape.IsNull() && !subShape.IsNull() )
+ {
+ TopExp_Explorer e( shape, subShape.ShapeType() );
+ if ( shape.Orientation() >= TopAbs_INTERNAL ) // TopAbs_INTERNAL or TopAbs_EXTERNAL
+ e.Init( shape.Oriented(TopAbs_FORWARD), subShape.ShapeType() );
+ for ( ; e.More(); e.Next())
+ if ( subShape.IsSame( e.Current() ))
+ break;
+ if ( e.More() )
+ ori = e.Current().Orientation();
+ }
+ return ori;
+}
+
+//=======================================================================
+//function : IsQuadraticMesh
+//purpose : Check mesh without geometry for: if all elements on this shape are quadratic,
+// quadratic elements will be created.
+// Used then generated 3D mesh without geometry.
//=======================================================================
SMESH_MesherHelper:: MType SMESH_MesherHelper::IsQuadraticMesh()
}
//=======================================================================
-/*!
- * \brief Return an alternative parameter for a node on seam
- */
+//function : GetOtherParam
+//purpose : Return an alternative parameter for a node on seam
//=======================================================================
double SMESH_MesherHelper::GetOtherParam(const double param) const
{
- return fabs(param-myPar1) < fabs(param-myPar2) ? myPar2 : myPar1;
+ int i = myParIndex & U_periodic ? 0 : 1;
+ return fabs(param-myPar1[i]) < fabs(param-myPar2[i]) ? myPar2[i] : myPar1[i];
}
//=======================================================================
mutable vector< const QLink* > _sides;
mutable bool _sideIsAdded[4]; // added in chain of links
gp_Vec _normal;
+#ifdef _DEBUG_
+ mutable const SMDS_MeshElement* _face;
+#endif
- QFace( const vector< const QLink*>& links );
+ QFace( const vector< const QLink*>& links, const SMDS_MeshElement* face=0 );
void SetVolume(const SMDS_MeshElement* v) const { _volumes[ _volumes[0] ? 1 : 0 ] = v; }
for (int i=0; i<_sides.size(); ++i ) if ( _sides[i] == side ) return i;
return -1;
}
- bool GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& error) const;
+ bool GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& err) const;
- bool GetLinkChain( TChainLink& link, TChain& chain, SMDS_TypeOfPosition pos, int& error) const
+ bool GetLinkChain( TChainLink& link, TChain& chain, SMDS_TypeOfPosition pos, int& err) const
{
int i = LinkIndex( link._qlink );
if ( i < 0 ) return true;
_sideIsAdded[i] = true;
link.SetFace( this );
// continue from opposite link
- return GetLinkChain( (i+2)%_sides.size(), chain, pos, error );
+ return GetLinkChain( (i+2)%_sides.size(), chain, pos, err );
}
bool IsBoundary() const { return !_volumes[1]; }
const TChainLink& avoidLink,
TLinkInSet * notBoundaryLink = 0,
const SMDS_MeshNode* nodeToContain = 0,
- bool * isAdjacentUsed = 0) const;
+ bool * isAdjacentUsed = 0,
+ int nbRecursionsLeft = -1) const;
TLinkInSet GetLinkByNode( const TLinkSet& links,
const TChainLink& avoidLink,
*/
//================================================================================
- QFace::QFace( const vector< const QLink*>& links )
+ QFace::QFace( const vector< const QLink*>& links, const SMDS_MeshElement* face )
{
_volumes[0] = _volumes[1] = 0;
_sides = links;
_normal /= sqrt( normSqSize );
else
_normal.SetCoord(1e-33,0,0);
+
+#ifdef _DEBUG_
+ _face = face;
+#endif
}
//================================================================================
/*!
if ( _sides.size() != 4 ) { // triangle - visit all my continous faces
MSGBEG( *this );
- for ( int i = 0; i < _sides.size(); ++i ) {
- if ( !_sideIsAdded[i] && _sides[i] ) {
- _sideIsAdded[i]=true;
- TChain::iterator chLink = chain.insert( chain.begin(), TChainLink(_sides[i]));
- chLink->SetFace( this );
- if ( _sides[i]->MediumPos() >= pos )
- if ( const QFace* f = _sides[i]->GetContinuesFace( this ))
- f->GetLinkChain( *chLink, chain, pos, error );
+ list< const QFace* > faces( 1, this );
+ for (list< const QFace* >::iterator fIt = faces.begin(); fIt != faces.end(); ++fIt ) {
+ const QFace* face = *fIt;
+ for ( int i = 0; i < face->_sides.size(); ++i ) {
+ if ( !face->_sideIsAdded[i] && face->_sides[i] ) {
+ face->_sideIsAdded[i] = true;
+ TChain::iterator chLink = chain.insert( chain.begin(), TChainLink(face->_sides[i]));
+ chLink->SetFace( face );
+ if ( face->_sides[i]->MediumPos() >= pos )
+ if ( const QFace* contFace = face->_sides[i]->GetContinuesFace( face ))
+ faces.push_back( contFace );
+ }
}
}
if ( error < ERR_TRI )
chLink->SetFace( this );
MSGBEG( *this );
- // propagate from rectangle to neighbour faces
+ // propagate from quadrangle to neighbour faces
if ( link->MediumPos() >= pos ) {
int nbLinkFaces = link->_faces.size();
if ( nbLinkFaces == 4 || nbLinkFaces < 4 && link->OnBoundary()) {
* \param nodeToContain - node the returned link must contain; if provided, search
* also performed on adjacent faces
* \param isAdjacentUsed - returns true if link is found in adjacent faces
+ * \param nbRecursionsLeft - to limit recursion
*/
//================================================================================
const TChainLink& avoidLink,
TLinkInSet * notBoundaryLink,
const SMDS_MeshNode* nodeToContain,
- bool * isAdjacentUsed) const
+ bool * isAdjacentUsed,
+ int nbRecursionsLeft) const
{
TLinkInSet linksEnd = links.end(), boundaryLink = linksEnd;
continue;
TLinkInSet link = links.find( _sides[iL] );
if ( link == linksEnd ) continue;
+ if ( (*link)->MediumPos() > SMDS_TOP_FACE )
+ continue; // We work on faces here, don't go into a volume
// check link
if ( link->IsBoundary() ) {
if ( boundaryLink != linksEnd ) break;
}
- if ( boundaryLink == linksEnd && nodeToContain ) // cellect adjacent faces
+ if ( boundaryLink == linksEnd && nodeToContain ) // collect adjacent faces
if ( const QFace* adj = link->NextFace( this ))
if ( adj->Contains( nodeToContain ))
adjacentFaces.push_back( make_pair( adj, link ));
}
if ( isAdjacentUsed ) *isAdjacentUsed = false;
- if ( boundaryLink == linksEnd && nodeToContain ) // check adjacent faces
+ if ( boundaryLink == linksEnd && nodeToContain && nbRecursionsLeft) // check adjacent faces
{
+ if ( nbRecursionsLeft < 0 )
+ nbRecursionsLeft = nodeToContain->NbInverseElements();
TFaceLinkList::iterator adj = adjacentFaces.begin();
for ( ; boundaryLink == linksEnd && adj != adjacentFaces.end(); ++adj )
- boundaryLink = adj->first->GetBoundaryLink( links, *(adj->second),
- 0, nodeToContain, isAdjacentUsed);
+ boundaryLink = adj->first->GetBoundaryLink( links, *(adj->second), 0, nodeToContain,
+ isAdjacentUsed, nbRecursionsLeft-1);
if ( isAdjacentUsed ) *isAdjacentUsed = true;
}
return boundaryLink;
int iL1 = (iL + 1) % 3, iL2 = (iL + 2) % 3; // indices of the two other links of triangle
TLinkInSet link1 = theLinks.find( _sides[iL1] );
TLinkInSet link2 = theLinks.find( _sides[iL2] );
+ if ( link1 == theLinks.end() || link2 == theLinks.end() )
+ return thePrevLen;
const QFace* f1 = link1->NextFace( this ); // adjacent faces
const QFace* f2 = link2->NextFace( this );
for ( TopExp_Explorer f(myShape,TopAbs_FACE,TopAbs_SOLID); f.More(); f.Next() ) {
faces.Add( f.Current() );
}
- for ( TopExp_Explorer v(myShape,TopAbs_SOLID); v.More(); v.Next() ) {
- if ( myMesh->GetSubMesh( v.Current() )->IsEmpty() ) { // get faces of solid
- for ( TopExp_Explorer f( v.Current(), TopAbs_FACE); f.More(); f.Next() )
+ for ( TopExp_Explorer s(myShape,TopAbs_SOLID); s.More(); s.Next() ) {
+ if ( myMesh->GetSubMesh( s.Current() )->IsEmpty() ) { // get faces of solid
+ for ( TopExp_Explorer f( s.Current(), TopAbs_FACE); f.More(); f.Next() )
faces.Add( f.Current() );
}
else { // fix nodes in the solid and its faces
SMESH_MesherHelper h(*myMesh);
- h.SetSubShape( v.Current() );
+ h.SetSubShape( s.Current() );
h.FixQuadraticElements(false);
}
}
for ( TopTools_MapIteratorOfMapOfShape fIt( faces ); fIt.More(); fIt.Next() ) {
SMESH_MesherHelper h(*myMesh);
h.SetSubShape( fIt.Key() );
- h.FixQuadraticElements();
+ h.FixQuadraticElements(true);
}
return;
}
hasRectFaces = hasRectFaces ||
( volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_HEXA ||
volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_PENTA );
+#ifdef _DEBUG_
+ if ( nbN == 6 )
+ pFace->_face = GetMeshDS()->FindFace(faceNodes[0],faceNodes[2],faceNodes[4]);
+ else
+ pFace->_face = GetMeshDS()->FindFace(faceNodes[0],faceNodes[2],
+ faceNodes[4],faceNodes[6] );
+#endif
}
}
set< QLink >::iterator pLink = links.begin();
// not treat boundary of volumic submesh
int isInside = ( elemType == SMDSAbs_Volume && volumeOnly ) ? 1 : 0;
for ( ; isInside < 2; ++isInside ) {
- MSG( "--------------- LOOP " << isInside << " ------------------");
+ MSG( "--------------- LOOP (inside=" << isInside << ") ------------------");
SMDS_TypeOfPosition pos = isInside ? SMDS_TOP_3DSPACE : SMDS_TOP_FACE;
for ( pFace = faces.begin(); pFace != faces.end(); ++pFace ) {
if ( !pFace->GetLinkChain( dir+2, rawChain, pos, error ) && error ==ERR_UNKNOWN ) continue;
vector< TChain > chains;
- if ( error == ERR_OK ) { // chains contains continues rectangles
+ if ( error == ERR_OK ) { // chain contains continues rectangles
chains.resize(1);
chains[0].splice( chains[0].begin(), rawChain );
}
- else if ( error == ERR_TRI ) { // chains contains continues triangles
+ else if ( error == ERR_TRI ) { // chain contains continues triangles
TSplitTriaResult res = splitTrianglesIntoChains( rawChain, chains, pos );
- if ( res != _OK ) { // not rectangles split into triangles
+ if ( res != _OK ) { // not quadrangles split into triangles
fixTriaNearBoundary( rawChain, *this );
break;
}
}
- else if ( error == ERR_PRISM ) { // side faces of prisms
+ else if ( error == ERR_PRISM ) { // quadrangle side faces of prisms
fixPrism( rawChain );
break;
}
// compute node displacement of end links in parametric space of face
const SMDS_MeshNode* nodeOnFace = (*(++chain.begin()))->_mediumNode;
TopoDS_Shape f = GetSubShapeByNode( nodeOnFace, GetMeshDS() );
- if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE ) {
+ if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE )
+ {
face = TopoDS::Face( f );
- for ( int is1 = 0; is1 < 2; ++is1 ) { // move0 or move1
+ Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
+ for ( int is1 = 0; is1 < 2; ++is1 ) // move0 or move1
+ {
TChainLink& link = is1 ? chain.back() : chain.front();
+ gp_XY uvm = GetNodeUV( face, link->_mediumNode, nodeOnFace, &checkUV);
gp_XY uv1 = GetNodeUV( face, link->node1(), nodeOnFace, &checkUV);
gp_XY uv2 = GetNodeUV( face, link->node2(), nodeOnFace, &checkUV);
- gp_XY uvm = GetNodeUV( face, link->_mediumNode, nodeOnFace, &checkUV);
- gp_XY uvMove = uvm - GetMiddleUV( BRep_Tool::Surface(face,loc), uv1, uv2);
- if ( is1 ) move1.SetCoord( uvMove.X(), uvMove.Y(), 0 );
- else move0.SetCoord( uvMove.X(), uvMove.Y(), 0 );
+ gp_XY uv12 = GetMiddleUV( surf, uv1, uv2);
+ // uvMove = uvm - uv12
+ gp_XY uvMove = applyFunc<TSubtract>(surf, uvm, uv12,/*inPeriod=*/false);
+ ( is1 ? move1 : move0 ).SetCoord( uvMove.X(), uvMove.Y(), 0 );
}
if ( move0.SquareMagnitude() < straightTol2 &&
move1.SquareMagnitude() < straightTol2 ) {
}
else {
// compute 3D displacement by 2D one
+ Handle(Geom_Surface) s = BRep_Tool::Surface(face,loc);
gp_XY oldUV = GetNodeUV( face, (*link1)->_mediumNode, 0, &checkUV);
- gp_XY newUV = oldUV + gp_XY( move.X(), move.Y() );
- gp_Pnt newPnt = BRep_Tool::Surface(face,loc)->Value( newUV.X(), newUV.Y());
+ gp_XY newUV = applyFunc<TAdd>( s, oldUV, gp_XY( move.X(),move.Y() ));
+ gp_Pnt newPnt = s->Value( newUV.X(), newUV.Y());
move = gp_Vec( XYZ((*link1)->_mediumNode), newPnt.Transformed(loc) );
#ifdef _DEBUG_
if ( (XYZ((*link1)->node1()) - XYZ((*link1)->node2())).SquareModulus() <
}
}
}
+
+//=======================================================================
+/*!
+ * \brief Iterator on ancestors of the given type
+ */
+//=======================================================================
+
+struct TAncestorsIterator : public SMDS_Iterator<const TopoDS_Shape*>
+{
+ TopTools_ListIteratorOfListOfShape _ancIter;
+ TopAbs_ShapeEnum _type;
+ TAncestorsIterator( const TopTools_ListOfShape& ancestors, TopAbs_ShapeEnum type)
+ : _ancIter( ancestors ), _type( type )
+ {
+ if ( _ancIter.More() && _ancIter.Value().ShapeType() != _type ) next();
+ }
+ virtual bool more()
+ {
+ return _ancIter.More();
+ }
+ virtual const TopoDS_Shape* next()
+ {
+ const TopoDS_Shape* s = _ancIter.More() ? & _ancIter.Value() : 0;
+ if ( _ancIter.More() )
+ for ( _ancIter.Next(); _ancIter.More(); _ancIter.Next())
+ if ( _ancIter.Value().ShapeType() == _type )
+ break;
+ return s;
+ }
+};
+
+//=======================================================================
+/*!
+ * \brief Return iterator on ancestors of the given type
+ */
+//=======================================================================
+
+PShapeIteratorPtr SMESH_MesherHelper::GetAncestors(const TopoDS_Shape& shape,
+ const SMESH_Mesh& mesh,
+ TopAbs_ShapeEnum ancestorType)
+{
+ return PShapeIteratorPtr( new TAncestorsIterator( mesh.GetAncestors(shape), ancestorType));
+}