0021893: EDF 2133 SMESH : Improvement of 3D extrusion algorithm

-  int NbPoints() const { return myNbPonits; }
+  int NbPoints(const bool update = false) const;

-  int NbSegments() const { return myNbSegments; }
+  int NbSegments(const bool update = false) const;

+  gp_Pnt   Value3d(double U) const;

arg theFirstVertex of SMESH_Block::GetOrderedEdges() became optional
This commit is contained in:
eap 2013-01-28 08:29:06 +00:00
parent 4b17835d08
commit c8d7fe2beb
2 changed files with 116 additions and 13 deletions

View File

@ -160,11 +160,10 @@ StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face& theFace,
double d4 = GCPnts_AbscissaPoint::Length( A2dC, myFirst[i], p4 );
//cout<<"len = "<<len<<" d2 = "<<d2<<" fabs(2*d2/len-1.0) = "<<fabs(2*d2/len-1.0)<<endl;
myIsUniform[i] = !( fabs(2*d2/myEdgeLength[i]-1.0) > 0.01 || fabs(2*d4/d2-1.0) > 0.01 );
if ( !myIsUniform[i] )
//if ( !myIsUniform[i] ) to implement Value3d(u)
{
double fp,lp;
TopLoc_Location L;
Handle(Geom_Curve) C3d = BRep_Tool::Curve(myEdge[i],L,fp,lp);
Handle(Geom_Curve) C3d = BRep_Tool::Curve(myEdge[i],fp,lp);
myC3dAdaptor[i].Load( C3d, fp,lp );
}
}
@ -208,8 +207,8 @@ StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face& theFace,
*/
//================================================================================
StdMeshers_FaceSide::StdMeshers_FaceSide(const SMDS_MeshNode* theNode,
const gp_Pnt2d thePnt2d,
StdMeshers_FaceSide::StdMeshers_FaceSide(const SMDS_MeshNode* theNode,
const gp_Pnt2d thePnt2d,
const StdMeshers_FaceSide* theSide)
{
myC2d.resize(1);
@ -218,7 +217,8 @@ StdMeshers_FaceSide::StdMeshers_FaceSide(const SMDS_MeshNode* theNode,
myDefaultPnt2d = thePnt2d;
myPoints = theSide->GetUVPtStruct();
myNbPonits = myNbSegments = myPoints.size();
myNbPonits = myPoints.size();
myNbSegments = theSide->myNbSegments;
std::vector<uvPtStruct>::iterator it = myPoints.begin();
for(; it!=myPoints.end(); it++) {
(*it).u = thePnt2d.X();
@ -336,7 +336,8 @@ const vector<UVPtStruct>& StdMeshers_FaceSide::GetUVPtStruct(bool isXConst,
}
} // loop on myEdge's
if ( u2node.size() + nbProxyNodes != myNbPonits )
if ( u2node.size() + nbProxyNodes != myNbPonits &&
u2node.size() + nbProxyNodes != NbPoints( /*update=*/true ))
{
MESSAGE("Wrong node parameters on edges, u2node.size():"
<<u2node.size()<<" != myNbPonits:"<<myNbPonits);
@ -601,6 +602,73 @@ void StdMeshers_FaceSide::Reverse()
for ( size_t i = 0; i < myEdge.size(); ++i )
reverseProxySubmesh( myEdge[i] );
}
for ( size_t i = 0; i < myEdge.size(); ++i )
{
double fp,lp;
Handle(Geom_Curve) C3d = BRep_Tool::Curve(myEdge[i],fp,lp);
myC3dAdaptor[i].Load( C3d, fp,lp );
}
}
//=======================================================================
//function : NbPoints
//purpose : Return nb nodes on edges and vertices (+1 to be == GetUVPtStruct().size() )
// Call it with update == true if mesh of this side can be recomputed
// since creation of this side
//=======================================================================
int StdMeshers_FaceSide::NbPoints(const bool update) const
{
if ( !myPoints.empty() )
return myPoints.size();
// if ( !myFalsePoints.empty() )
// return myFalsePoints.size();
if ( update && myEdge.size() > 0 )
{
StdMeshers_FaceSide* me = (StdMeshers_FaceSide*) this;
me->myNbPonits = 0;
me->myNbSegments = 0;
me->myMissingVertexNodes = false;
for ( int i = 0; i < NbEdges(); ++i )
{
TopoDS_Vertex v1 = SMESH_MesherHelper::IthVertex( 0, myEdge[i] );
if ( SMESH_Algo::VertexNode( v1, myProxyMesh->GetMeshDS() ))
me->myNbPonits += 1; // for the first end
else
me->myMissingVertexNodes = true;
if ( const SMESHDS_SubMesh* sm = myProxyMesh->GetSubMesh( Edge(i) )) {
int nbN = sm->NbNodes();
if ( myIgnoreMediumNodes ) {
SMDS_ElemIteratorPtr elemIt = sm->GetElements();
if ( elemIt->more() && elemIt->next()->IsQuadratic() )
nbN -= sm->NbElements();
}
me->myNbPonits += nbN;
}
}
TopoDS_Vertex v1 = SMESH_MesherHelper::IthVertex( 1, Edge( NbEdges()-1 ));
if ( SMESH_Algo::VertexNode( v1, myProxyMesh->GetMeshDS() ))
me->myNbPonits++; // for the last end
else
me->myMissingVertexNodes = true;
}
return myNbPonits;
}
//=======================================================================
//function : NbSegments
//purpose : Return nb edges
// Call it with update == true if mesh of this side can be recomputed
// since creation of this side
//=======================================================================
int StdMeshers_FaceSide::NbSegments(const bool update) const
{
return Max( 0, NbPoints( update ) - 1 );
}
//================================================================================
@ -709,7 +777,6 @@ BRepAdaptor_CompCurve* StdMeshers_FaceSide::GetCurve3d() const
return new BRepAdaptor_CompCurve( aWire );
}
//================================================================================
/*!
* \brief Return 2D point by normalized parameter
@ -743,6 +810,35 @@ gp_Pnt2d StdMeshers_FaceSide::Value2d(double U) const
return myDefaultPnt2d;
}
//================================================================================
/*!
* \brief Return XYZ by normalized parameter
* \param U - normalized parameter value
* \retval gp_Pnt - point
*/
//================================================================================
gp_Pnt StdMeshers_FaceSide::Value3d(double U) const
{
int i = EdgeIndex( U );
double prevU = i ? myNormPar[ i-1 ] : 0;
double r = ( U - prevU )/ ( myNormPar[ i ] - prevU );
double par = myFirst[i] * ( 1 - r ) + myLast[i] * r;
// check parametrization of curve
if( !myIsUniform[i] )
{
double aLen3dU = r * myEdgeLength[i] * ( myFirst[i]>myLast[i] ? -1. : 1.);
GCPnts_AbscissaPoint AbPnt
( const_cast<GeomAdaptor_Curve&>( myC3dAdaptor[i]), aLen3dU, myFirst[i] );
if( AbPnt.IsDone() ) {
par = AbPnt.Parameter();
}
}
return myC3dAdaptor[ i ].Value(par);
}
//================================================================================
/*!
* \brief Return wires of a face as StdMeshers_FaceSide's
@ -755,10 +851,9 @@ TSideVector StdMeshers_FaceSide::GetFaceWires(const TopoDS_Face& theFace,
TError & theError,
SMESH_ProxyMesh::Ptr theProxyMesh)
{
TopoDS_Vertex V1;
list< TopoDS_Edge > edges, internalEdges;
list< int > nbEdgesInWires;
int nbWires = SMESH_Block::GetOrderedEdges (theFace, V1, edges, nbEdgesInWires);
int nbWires = SMESH_Block::GetOrderedEdges (theFace, edges, nbEdgesInWires);
// split list of all edges into separate wires
TSideVector wires( nbWires );

View File

@ -102,13 +102,17 @@ public:
*/
void Reverse();
/*!
* \brief Return nb nodes on edges and vertices (+1 to be == GetUVPtStruct().size() )
* \brief Return nb nodes on edges and vertices (+1 to be == GetUVPtStruct().size() ).
* Call it with update == true if mesh of this side can be recomputed
* since creation of this side
*/
int NbPoints() const { return myNbPonits; }
int NbPoints(const bool update = false) const;
/*!
* \brief Return nb edges
* Call it with update == true if mesh of this side can be recomputed
* since creation of this side
*/
int NbSegments() const { return myNbSegments; }
int NbSegments(const bool update = false) const;
/*!
* \brief Return mesh
*/
@ -147,6 +151,10 @@ public:
* \brief Return UV by normalized parameter
*/
gp_Pnt2d Value2d(double U) const;
/*!
* \brief Return XYZ by normalized parameter
*/
gp_Pnt Value3d(double U) const;
/*!
* \brief Creates a Adaptor2d_Curve2d to be used in SMESH_Block
*/