PAL13615 (EDF PAL 315/31 GEOM SMESH : meshing of a "5 edges quadrangle")

implementation
This commit is contained in:
eap 2007-02-20 07:14:12 +00:00
parent 1863381f98
commit 54d60d7615
4 changed files with 448 additions and 1248 deletions

View File

@ -27,12 +27,15 @@
// Module : SMESH // Module : SMESH
// $Header$ // $Header$
using namespace std;
#include "StdMeshers_MEFISTO_2D.hxx" #include "StdMeshers_MEFISTO_2D.hxx"
#include "SMESH_Gen.hxx" #include "SMESH_Gen.hxx"
#include "SMESH_Mesh.hxx" #include "SMESH_Mesh.hxx"
#include "SMESH_subMesh.hxx" #include "SMESH_subMesh.hxx"
#include "SMESH_Block.hxx"
#include "SMESH_MesherHelper.hxx"
#include "StdMeshers_FaceSide.hxx"
#include "StdMeshers_MaxElementArea.hxx" #include "StdMeshers_MaxElementArea.hxx"
#include "StdMeshers_LengthFromEdges.hxx" #include "StdMeshers_LengthFromEdges.hxx"
@ -48,19 +51,16 @@ using namespace std;
#include <TopoDS_Face.hxx> #include <TopoDS_Face.hxx>
#include <TopoDS_Edge.hxx> #include <TopoDS_Edge.hxx>
#include <TopoDS_Shape.hxx>
#include <Geom_Surface.hxx> #include <Geom_Surface.hxx>
#include <GeomAdaptor_Curve.hxx>
#include <Geom2d_Curve.hxx> #include <Geom2d_Curve.hxx>
#include <gp_Pnt2d.hxx> #include <gp_Pnt2d.hxx>
#include <BRep_Tool.hxx> #include <BRep_Tool.hxx>
#include <BRepTools.hxx> #include <BRepTools.hxx>
#include <BRepTools_WireExplorer.hxx>
#include <GCPnts_AbscissaPoint.hxx>
#include <GCPnts_UniformAbscissa.hxx>
#include <TColStd_ListIteratorOfListOfInteger.hxx>
#include <TopTools_ListIteratorOfListOfShape.hxx> #include <TopTools_ListIteratorOfListOfShape.hxx>
#include <TopTools_ListOfShape.hxx> #include <TopTools_ListOfShape.hxx>
#include <TopTools_MapOfShape.hxx>
using namespace std;
//============================================================================= //=============================================================================
/*! /*!
@ -68,8 +68,8 @@ using namespace std;
*/ */
//============================================================================= //=============================================================================
StdMeshers_MEFISTO_2D::StdMeshers_MEFISTO_2D(int hypId, int studyId, StdMeshers_MEFISTO_2D::StdMeshers_MEFISTO_2D(int hypId, int studyId, SMESH_Gen * gen):
SMESH_Gen * gen):SMESH_2D_Algo(hypId, studyId, gen) SMESH_2D_Algo(hypId, studyId, gen)
{ {
MESSAGE("StdMeshers_MEFISTO_2D::StdMeshers_MEFISTO_2D"); MESSAGE("StdMeshers_MEFISTO_2D::StdMeshers_MEFISTO_2D");
_name = "MEFISTO_2D"; _name = "MEFISTO_2D";
@ -102,75 +102,69 @@ StdMeshers_MEFISTO_2D::~StdMeshers_MEFISTO_2D()
//============================================================================= //=============================================================================
bool StdMeshers_MEFISTO_2D::CheckHypothesis bool StdMeshers_MEFISTO_2D::CheckHypothesis
(SMESH_Mesh& aMesh, (SMESH_Mesh& aMesh,
const TopoDS_Shape& aShape, const TopoDS_Shape& aShape,
SMESH_Hypothesis::Hypothesis_Status& aStatus) SMESH_Hypothesis::Hypothesis_Status& aStatus)
{ {
//MESSAGE("StdMeshers_MEFISTO_2D::CheckHypothesis"); _hypMaxElementArea = NULL;
_hypLengthFromEdges = NULL;
_hypMaxElementArea = NULL; list <const SMESHDS_Hypothesis * >::const_iterator itl;
_hypLengthFromEdges = NULL; const SMESHDS_Hypothesis *theHyp;
list <const SMESHDS_Hypothesis * >::const_iterator itl; const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
const SMESHDS_Hypothesis *theHyp; int nbHyp = hyps.size();
if (!nbHyp)
{
aStatus = SMESH_Hypothesis::HYP_MISSING;
return false; // can't work with no hypothesis
}
const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape); itl = hyps.begin();
int nbHyp = hyps.size(); theHyp = (*itl); // use only the first hypothesis
if (!nbHyp)
{
aStatus = SMESH_Hypothesis::HYP_MISSING;
return false; // can't work with no hypothesis
}
itl = hyps.begin(); string hypName = theHyp->GetName();
theHyp = (*itl); // use only the first hypothesis
string hypName = theHyp->GetName(); bool isOk = false;
//int hypId = theHyp->GetID();
//SCRUTE(hypName);
bool isOk = false; if (hypName == "MaxElementArea")
{
_hypMaxElementArea = static_cast<const StdMeshers_MaxElementArea *>(theHyp);
ASSERT(_hypMaxElementArea);
_maxElementArea = _hypMaxElementArea->GetMaxArea();
_edgeLength = 0;
isOk = true;
aStatus = SMESH_Hypothesis::HYP_OK;
}
if (hypName == "MaxElementArea") else if (hypName == "LengthFromEdges")
{ {
_hypMaxElementArea = static_cast<const StdMeshers_MaxElementArea *>(theHyp); _hypLengthFromEdges = static_cast<const StdMeshers_LengthFromEdges *>(theHyp);
ASSERT(_hypMaxElementArea); ASSERT(_hypLengthFromEdges);
_maxElementArea = _hypMaxElementArea->GetMaxArea(); _edgeLength = 0;
_edgeLength = 0; _maxElementArea = 0;
isOk = true; isOk = true;
aStatus = SMESH_Hypothesis::HYP_OK; aStatus = SMESH_Hypothesis::HYP_OK;
} }
else
aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
else if (hypName == "LengthFromEdges") if (isOk)
{ {
_hypLengthFromEdges = static_cast<const StdMeshers_LengthFromEdges *>(theHyp); isOk = false;
ASSERT(_hypLengthFromEdges); if (_maxElementArea > 0)
_edgeLength = 0; {
_maxElementArea = 0; //_edgeLength = 2 * sqrt(_maxElementArea); // triangles : minorant
isOk = true; _edgeLength = sqrt(2. * _maxElementArea/sqrt(3.0));
aStatus = SMESH_Hypothesis::HYP_OK; isOk = true;
} }
else else
aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE; isOk = (_hypLengthFromEdges != NULL); // **** check mode
if (!isOk)
aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
}
if (isOk) return isOk;
{
isOk = false;
if (_maxElementArea > 0)
{
// _edgeLength = 2 * sqrt(_maxElementArea); // triangles : minorant
_edgeLength = sqrt(2. * _maxElementArea/sqrt(3.0));
isOk = true;
}
else
isOk = (_hypLengthFromEdges != NULL); // **** check mode
if (!isOk)
aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
}
//SCRUTE(_edgeLength);
//SCRUTE(_maxElementArea);
return isOk;
} }
//============================================================================= //=============================================================================
@ -183,21 +177,58 @@ bool StdMeshers_MEFISTO_2D::Compute(SMESH_Mesh & aMesh, const TopoDS_Shape & aSh
{ {
MESSAGE("StdMeshers_MEFISTO_2D::Compute"); MESSAGE("StdMeshers_MEFISTO_2D::Compute");
if (_hypLengthFromEdges) TopoDS_Face F = TopoDS::Face(aShape.Oriented(TopAbs_FORWARD));
_edgeLength = ComputeEdgeElementLength(aMesh, aShape);
bool isOk = false;
//const SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
//SMESH_subMesh *theSubMesh = aMesh.GetSubMesh(aShape);
const TopoDS_Face & FF = TopoDS::Face(aShape); // get all edges of a face
bool faceIsForward = (FF.Orientation() == TopAbs_FORWARD); TopoDS_Vertex V1;
TopoDS_Face F = TopoDS::Face(FF.Oriented(TopAbs_FORWARD)); list< TopoDS_Edge > edges;
list< int > nbEdgesInWires;
int nbWires = SMESH_Block::GetOrderedEdges (F, V1, edges, nbEdgesInWires);
Z nblf; //nombre de lignes fermees (enveloppe en tete) if (_hypLengthFromEdges) _edgeLength = 0;
Z *nudslf = NULL; //numero du dernier sommet de chaque ligne fermee
R2 *uvslf = NULL; // split list of all edges into separate wires
Z nbpti = 0; //nombre points internes futurs sommets de la triangulation TWireVector wires ( nbWires );
list< int >::iterator nbE = nbEdgesInWires.begin();
list< TopoDS_Edge >::iterator from, to;
from = to = edges.begin();
for ( int iW = 0; iW < nbWires; ++iW )
{
std::advance( to, *nbE++ );
list< TopoDS_Edge > wireEdges( from, to );
// assure that there is a node on the first vertex
// as StdMeshers_FaceSide::GetUVPtStruct() requires
while ( !VertexNode( TopExp::FirstVertex( wireEdges.front(), true),
aMesh.GetMeshDS()))
{
wireEdges.splice(wireEdges.end(), wireEdges,
wireEdges.begin(), ++wireEdges.begin());
if ( from->IsSame( wireEdges.front() )) {
MESSAGE( "No nodes on vertices on wire " << iW+1);
return false;
}
}
StdMeshers_FaceSide* wire = new StdMeshers_FaceSide( F, wireEdges, &aMesh, true );
wires[ iW ] = StdMeshers_FaceSidePtr( wire );
if (_hypLengthFromEdges && wire->NbSegments() )
_edgeLength += wire->Length() / wire->NbSegments();
from = to;
}
if ( wires[0]->NbSegments() < 3 ) // ex: a circle with 2 segments
return false;
if (_hypLengthFromEdges && _edgeLength < DBL_MIN )
_edgeLength = 100;
// helper builds a quadratic mesh if necessary
myTool = new SMESH_MesherHelper(aMesh);
auto_ptr<SMESH_MesherHelper> helperDeleter( myTool );
_quadraticMesh = myTool->IsQuadraticSubMesh(aShape);
Z nblf; //nombre de lignes fermees (enveloppe en tete)
Z *nudslf = NULL; //numero du dernier sommet de chaque ligne fermee
R2 *uvslf = NULL;
Z nbpti = 0; //nombre points internes futurs sommets de la triangulation
R2 *uvpti = NULL; R2 *uvpti = NULL;
Z nbst; Z nbst;
@ -206,97 +237,55 @@ bool StdMeshers_MEFISTO_2D::Compute(SMESH_Mesh & aMesh, const TopoDS_Shape & aSh
Z *nust = NULL; Z *nust = NULL;
Z ierr = 0; Z ierr = 0;
Z nutysu = 1; // 1: il existe un fonction areteideale_() Z nutysu = 1; // 1: il existe un fonction areteideale_()
// Z nutysu=0; // 0: on utilise aretmx // Z nutysu=0; // 0: on utilise aretmx
R aretmx = _edgeLength; // longueur max aretes future triangulation R aretmx = _edgeLength; // longueur max aretes future triangulation
nblf = NumberOfWires(F); nblf = nbWires;
nudslf = new Z[1 + nblf]; nudslf = new Z[1 + nblf];
nudslf[0] = 0; nudslf[0] = 0;
int iw = 1; int iw = 1;
int nbpnt = 0; int nbpnt = 0;
myTool = new SMESH_MesherHelper(aMesh); // count nb of input points
_quadraticMesh = myTool->IsQuadraticSubMesh(aShape); for ( int iW = 0; iW < nbWires; ++iW )
{
if ( _quadraticMesh && _hypLengthFromEdges ) nbpnt += wires[iW]->NbSegments();
aretmx *= 2.; nudslf[iw++] = nbpnt;
myOuterWire = BRepTools::OuterWire(F);
nbpnt += NumberOfPoints(aMesh, myOuterWire);
if ( nbpnt < 3 ) { // ex: a circle with 2 segments
delete myTool; myTool = 0;
return false;
} }
nudslf[iw++] = nbpnt;
for (TopExp_Explorer exp(F, TopAbs_WIRE); exp.More(); exp.Next()) {
const TopoDS_Wire & W = TopoDS::Wire(exp.Current());
if (!myOuterWire.IsSame(W)) {
nbpnt += NumberOfPoints(aMesh, W);
nudslf[iw++] = nbpnt;
}
}
// avoid passing same uv points for a vertex common to 2 wires
TopTools_IndexedDataMapOfShapeListOfShape VWMap;
if ( iw - 1 > 1 ) // nbofWires > 1
TopExp::MapShapesAndAncestors( F , TopAbs_VERTEX, TopAbs_WIRE, VWMap );
uvslf = new R2[nudslf[nblf]]; uvslf = new R2[nudslf[nblf]];
int m = 0;
double scalex, scaley; double scalex, scaley;
ComputeScaleOnFace(aMesh, F, scalex, scaley); ComputeScaleOnFace(aMesh, F, scalex, scaley);
map<int, const SMDS_MeshNode*> mefistoToDS; // correspondence mefisto index--> points IDNodes // correspondence mefisto index --> Nodes
if ( !LoadPoints(aMesh, F, myOuterWire, uvslf, m, vector< const SMDS_MeshNode*> mefistoToDS(nbpnt, (const SMDS_MeshNode*)0);
mefistoToDS, scalex, scaley, VWMap) ) {
delete myTool; myTool = 0;
return false;
}
for (TopExp_Explorer exp(F, TopAbs_WIRE); exp.More(); exp.Next()) // fill input points UV
{ LoadPoints(wires, uvslf, mefistoToDS, scalex, scaley);
const TopoDS_Wire & W = TopoDS::Wire(exp.Current());
if (!myOuterWire.IsSame(W))
{
if (! LoadPoints(aMesh, F, W, uvslf, m,
mefistoToDS, scalex, scaley, VWMap )) {
delete myTool; myTool = 0;
return false;
}
}
}
uvst = NULL; // Compute
nust = NULL;
aptrte(nutysu, aretmx, aptrte(nutysu, aretmx,
nblf, nudslf, uvslf, nbpti, uvpti, nbst, uvst, nbt, nust, ierr); nblf, nudslf, uvslf, nbpti, uvpti, nbst, uvst, nbt, nust, ierr);
bool isOk = false;
if (ierr == 0) if (ierr == 0)
{ {
MESSAGE("... End Triangulation Generated Triangle Number " << nbt); MESSAGE("... End Triangulation Generated Triangle Number " << nbt);
MESSAGE(" Node Number " << nbst); MESSAGE(" Node Number " << nbst);
StoreResult(aMesh, nbst, uvst, nbt, nust, F, StoreResult(nbst, uvst, nbt, nust, mefistoToDS, scalex, scaley);
faceIsForward, mefistoToDS, scalex, scaley);
isOk = true; isOk = true;
} }
else else
{ {
MESSAGE("Error in Triangulation"); MESSAGE("Error in Triangulation");
isOk = false;
} }
if (nudslf != NULL) if (nudslf != NULL) delete[]nudslf;
delete[]nudslf; if (uvslf != NULL) delete[]uvslf;
if (uvslf != NULL) if (uvst != NULL) delete[]uvst;
delete[]uvslf; if (nust != NULL) delete[]nust;
if (uvst != NULL)
delete[]uvst;
if (nust != NULL)
delete[]nust;
delete myTool; myTool = 0;
return isOk; return isOk;
} }
@ -357,25 +346,30 @@ static bool fixOverlappedLinkUV( R2& uv0, const R2& uv1, const R2& uv2 )
//purpose : //purpose :
//======================================================================= //=======================================================================
static bool fixCommonVertexUV (gp_Pnt2d & theUV, static bool fixCommonVertexUV (R2 & theUV,
const TopoDS_Vertex& theV, const TopoDS_Vertex& theV,
const TopoDS_Wire& theW,
const TopoDS_Wire& theOW,
const TopoDS_Face& theF, const TopoDS_Face& theF,
const TopTools_IndexedDataMapOfShapeListOfShape & theVWMap, const TopTools_IndexedDataMapOfShapeListOfShape & theVWMap,
SMESH_Mesh & theMesh, SMESH_Mesh & theMesh,
bool CreateQuadratic) const double theScaleX,
const double theScaleY,
const bool theCreateQuadratic)
{ {
if( theW.IsSame( theOW ) || if( !theVWMap.Contains( theV )) return false;
!theVWMap.Contains( theV )) return false;
// check if there is another wire sharing theV // check if there is another wire sharing theV
const TopTools_ListOfShape& WList = theVWMap.FindFromKey( theV ); const TopTools_ListOfShape& WList = theVWMap.FindFromKey( theV );
TopTools_ListIteratorOfListOfShape aWIt; TopTools_ListIteratorOfListOfShape aWIt;
TopTools_MapOfShape aWires;
for ( aWIt.Initialize( WList ); aWIt.More(); aWIt.Next() ) for ( aWIt.Initialize( WList ); aWIt.More(); aWIt.Next() )
if ( !theW.IsSame( aWIt.Value() )) aWires.Add( aWIt.Value() );
break; if ( aWires.Extent() < 2 ) return false;
if ( !aWIt.More() ) return false;
TopoDS_Shape anOuterWire = BRepTools::OuterWire(theF);
TopoDS_Shape anInnerWire;
for ( aWIt.Initialize( WList ); aWIt.More() && anInnerWire.IsNull(); aWIt.Next() )
if ( !anOuterWire.IsSame( aWIt.Value() ))
anInnerWire = aWIt.Value();
TopTools_ListOfShape EList; TopTools_ListOfShape EList;
list< double > UList; list< double > UList;
@ -383,7 +377,7 @@ static bool fixCommonVertexUV (gp_Pnt2d & theUV,
// find edges of theW sharing theV // find edges of theW sharing theV
// and find 2d normal to them at theV // and find 2d normal to them at theV
gp_Vec2d N(0.,0.); gp_Vec2d N(0.,0.);
TopoDS_Iterator itE( theW ); TopoDS_Iterator itE( anInnerWire );
for ( ; itE.More(); itE.Next() ) for ( ; itE.More(); itE.Next() )
{ {
const TopoDS_Edge& E = TopoDS::Edge( itE.Value() ); const TopoDS_Edge& E = TopoDS::Edge( itE.Value() );
@ -401,7 +395,7 @@ static bool fixCommonVertexUV (gp_Pnt2d & theUV,
gp_Vec2d d1; gp_Vec2d d1;
gp_Pnt2d p; gp_Pnt2d p;
C2d->D1( u, p, d1 ); C2d->D1( u, p, d1 );
gp_Vec2d n( d1.Y(), -d1.X() ); gp_Vec2d n( d1.Y() * theScaleX, -d1.X() * theScaleY);
if ( E.Orientation() == TopAbs_REVERSED ) if ( E.Orientation() == TopAbs_REVERSED )
n.Reverse(); n.Reverse();
N += n.Normalized(); N += n.Normalized();
@ -411,6 +405,7 @@ static bool fixCommonVertexUV (gp_Pnt2d & theUV,
// define step size by which to move theUV // define step size by which to move theUV
gp_Pnt2d nextUV; // uv of next node on edge, most distant of the four gp_Pnt2d nextUV; // uv of next node on edge, most distant of the four
gp_Pnt2d thisUV( theUV.x, theUV.y );
double maxDist = -DBL_MAX; double maxDist = -DBL_MAX;
TopTools_ListIteratorOfListOfShape aEIt (EList); TopTools_ListIteratorOfListOfShape aEIt (EList);
list< double >::iterator aUIt = UList.begin(); list< double >::iterator aUIt = UList.begin();
@ -431,7 +426,7 @@ static bool fixCommonVertexUV (gp_Pnt2d & theUV,
while ( nIt->more() ) { while ( nIt->more() ) {
const SMDS_MeshNode* node = nIt->next(); const SMDS_MeshNode* node = nIt->next();
// check if node is medium // check if node is medium
if ( CreateQuadratic && SMESH_MesherHelper::IsMedium( node, SMDSAbs_Edge )) if ( theCreateQuadratic && SMESH_MesherHelper::IsMedium( node, SMDSAbs_Edge ))
continue; continue;
const SMDS_EdgePosition* epos = const SMDS_EdgePosition* epos =
static_cast<const SMDS_EdgePosition*>(node->GetPosition().get()); static_cast<const SMDS_EdgePosition*>(node->GetPosition().get());
@ -444,29 +439,33 @@ static bool fixCommonVertexUV (gp_Pnt2d & theUV,
} }
bool isFirstCommon = ( *aUIt == f ); bool isFirstCommon = ( *aUIt == f );
gp_Pnt2d uv = C2d->Value( isFirstCommon ? umin : umax ); gp_Pnt2d uv = C2d->Value( isFirstCommon ? umin : umax );
double dist = theUV.SquareDistance( uv ); double dist = thisUV.SquareDistance( uv );
if ( dist > maxDist ) { if ( dist > maxDist ) {
maxDist = dist; maxDist = dist;
nextUV = uv; nextUV = uv;
} }
} }
R2 uv0, uv1, uv2; R2 uv0, uv1, uv2;
uv0.x = theUV.X(); uv0.y = theUV.Y(); uv0.x = thisUV.X(); uv0.y = thisUV.Y();
uv1.x = nextUV.X(); uv1.y = nextUV.Y(); uv1.x = nextUV.X(); uv1.y = nextUV.Y();
uv2.x = uv0.x; uv2.y = uv0.y; uv2.x = thisUV.X(); uv2.y = thisUV.Y();
uv1.x *= theScaleX; uv1.y *= theScaleY;
if ( fixOverlappedLinkUV( uv0, uv1, uv2 )) if ( fixOverlappedLinkUV( uv0, uv1, uv2 ))
{ {
double step = theUV.Distance( gp_Pnt2d( uv0.x, uv0.y )); double step = thisUV.Distance( gp_Pnt2d( uv0.x, uv0.y ));
// move theUV along the normal by the step // move theUV along the normal by the step
N *= step; N *= step;
MESSAGE("--fixCommonVertexUV move(" << theUV.X() << " " << theUV.Y() MESSAGE("--fixCommonVertexUV move(" << theUV.x << " " << theUV.x
<< ") by (" << N.X() << " " << N.Y() << ")" << ") by (" << N.X() << " " << N.Y() << ")"
<< endl << "--- MAX DIST " << maxDist); << endl << "--- MAX DIST " << maxDist);
theUV.SetXY( theUV.XY() + N.XY() ); theUV.x += N.X();
theUV.y += N.Y();
return true; return true;
} }
@ -479,111 +478,71 @@ static bool fixCommonVertexUV (gp_Pnt2d & theUV,
*/ */
//============================================================================= //=============================================================================
bool StdMeshers_MEFISTO_2D::LoadPoints(SMESH_Mesh & aMesh, bool StdMeshers_MEFISTO_2D::LoadPoints(TWireVector & wires,
const TopoDS_Face & FF, R2 * uvslf,
const TopoDS_Wire & WW, vector<const SMDS_MeshNode*>& mefistoToDS,
R2 * uvslf, double scalex,
int & m, double scaley)
map<int, const SMDS_MeshNode*>&mefistoToDS,
double scalex,
double scaley,
const TopTools_IndexedDataMapOfShapeListOfShape& VWMap)
{ {
// MESSAGE("StdMeshers_MEFISTO_2D::LoadPoints"); // to avoid passing same uv points for a vertex common to 2 wires
TopoDS_Face F;
TopTools_IndexedDataMapOfShapeListOfShape VWMap;
if ( wires.size() > 1 )
{
F = TopoDS::Face( myTool->GetSubShape() );
TopExp::MapShapesAndAncestors( F, TopAbs_VERTEX, TopAbs_WIRE, VWMap );
int nbVertices = 0;
for ( int iW = 0; iW < wires.size(); ++iW )
nbVertices += wires[ iW ]->NbSegments();
if ( nbVertices == VWMap.Extent() )
VWMap.Clear(); // wires have no common vertices
}
TopoDS_Face F = TopoDS::Face(FF.Oriented(TopAbs_FORWARD)); const bool isXConst = false; // meaningles here
const double constValue = 0; // meaningles here
int m = 0;
list< int > mOnVertex; list< int > mOnVertex;
gp_XY scale( scalex, scaley ); for ( int iW = 0; iW < wires.size(); ++iW )
TopoDS_Wire W = TopoDS::Wire(WW.Oriented(TopAbs_FORWARD));
BRepTools_WireExplorer wexp(W, F);
for ( wexp.Init(W, F); wexp.More(); wexp.Next() )
{ {
const TopoDS_Edge & E = wexp.Current(); const vector<UVPtStruct>& uvPtVec = wires[ iW ]->GetUVPtStruct(isXConst,constValue);
bool isForward = (E.Orientation() == TopAbs_FORWARD);
// --- IDNodes of first and last Vertex vector<UVPtStruct>::const_iterator uvPt = uvPtVec.begin();
for ( ++uvPt; uvPt != uvPtVec.end(); ++uvPt )
TopoDS_Vertex VFirst, VLast, V;
TopExp::Vertices(E, VFirst, VLast); // corresponds to f and l
V = isForward ? VFirst : VLast;
ASSERT(!V.IsNull());
SMDS_NodeIteratorPtr lid= aMesh.GetSubMesh(V)->GetSubMeshDS()->GetNodes();
if ( !lid->more() ) {
MESSAGE (" NO NODE BUILT ON VERTEX ");
return false;
}
const SMDS_MeshNode* idFirst = lid->next();
// --- edge internal IDNodes (relies on good order storage, not checked)
map<double, const SMDS_MeshNode*> params;
const SMDS_MeshNode * node;
SMDS_NodeIteratorPtr nodeIt= aMesh.GetSubMesh(E)->GetSubMeshDS()->GetNodes();
while ( nodeIt->more() )
{ {
node = nodeIt->next(); // bind mefisto ID to node
if ( _quadraticMesh && SMESH_MesherHelper::IsMedium( node, SMDSAbs_Edge )) mefistoToDS[m] = uvPt->node;
continue; // set UV
const SMDS_EdgePosition* epos = uvslf[m].x = uvPt->u * scalex;
static_cast<const SMDS_EdgePosition*>(node->GetPosition().get()); uvslf[m].y = uvPt->v * scaley;
double param = epos->GetUParameter(); if ( uvPt->node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX )
if ( !isForward ) param = -param; mOnVertex.push_back( m );
if ( !params.insert( make_pair( param, node )).second )
{
MESSAGE( "BAD NODE ON EDGE POSITIONS" );
return false;
}
}
// --- load 2D values into MEFISTO structure,
// add IDNodes in mefistoToDS map
double f, l, uFirst, u;
Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(E, F, f, l);
uFirst = isForward ? f : l;
// vertex node
gp_Pnt2d p = C2d->Value( uFirst ).XY().Multiplied( scale );
if ( fixCommonVertexUV( p, V, W, myOuterWire, F, VWMap, aMesh, _quadraticMesh ))
myNodesOnCommonV.push_back( idFirst );
mOnVertex.push_back( m );
uvslf[m].x = p.X();
uvslf[m].y = p.Y();
mefistoToDS[m + 1] = idFirst;
//MESSAGE(" "<<m<<" "<<mefistoToDS[m+1]);
//MESSAGE("__ f "<<uFirst<<" "<<uvslf[m].x <<" "<<uvslf[m].y);
m++;
// internal nodes
map<double, const SMDS_MeshNode*>::iterator u_n = params.begin();
for ( int i = 0; u_n != params.end(); ++u_n, ++i )
{
u = isForward ? u_n->first : - u_n->first;
gp_Pnt2d p = C2d->Value( u ).XY().Multiplied( scale );
uvslf[m].x = p.X();
uvslf[m].y = p.Y();
mefistoToDS[m + 1] = u_n->second;
//MESSAGE(" "<<m<<" "<<mefistoToDS[m+1]);
//MESSAGE("__ "<<i<<" "<<param<<" "<<uvslf[m].x <<" "<<uvslf[m].y);
m++; m++;
} }
} // for wexp
// prevent failure on overlapped adjacent links, int mFirst = mOnVertex.front(), mLast = m - 1;
// check only links ending in vertex nodes list< int >::iterator mIt = mOnVertex.begin();
int mFirst = mOnVertex.front(), mLast = m - 1; for ( ; mIt != mOnVertex.end(); ++mIt)
list< int >::iterator mIt = mOnVertex.begin(); {
for ( ; mIt != mOnVertex.end(); ++mIt ) { int m = *mIt;
int i = *mIt; if ( iW && !VWMap.IsEmpty()) { // except outer wire
int iB = i - 1, iA = i + 1; // indices Before and After // avoid passing same uv point for a vertex common to 2 wires
if ( iB < mFirst ) iB = mLast; int vID = mefistoToDS[m]->GetPosition()->GetShapeId();
if ( iA > mLast ) iA = mFirst; TopoDS_Vertex V = TopoDS::Vertex( myTool->GetMeshDS()->IndexToShape( vID ));
fixOverlappedLinkUV (uvslf[ iB ], uvslf[ i ], uvslf[ iA ]); if ( fixCommonVertexUV( uvslf[m], V, F, VWMap, *myTool->GetMesh(),
scalex, scaley, _quadraticMesh )) {
myNodesOnCommonV.push_back( mefistoToDS[m] );
continue;
}
}
// prevent failure on overlapped adjacent links,
// check only links ending in vertex nodes
int mB = m - 1, mA = m + 1; // indices Before and After
if ( mB < mFirst ) mB = mLast;
if ( mA > mLast ) mA = mFirst;
fixOverlappedLinkUV (uvslf[ mB ], uvslf[ m ], uvslf[ mA ]);
}
} }
return true; return true;
@ -595,12 +554,12 @@ bool StdMeshers_MEFISTO_2D::LoadPoints(SMESH_Mesh & aMesh,
*/ */
//============================================================================= //=============================================================================
void StdMeshers_MEFISTO_2D::ComputeScaleOnFace(SMESH_Mesh & aMesh, void StdMeshers_MEFISTO_2D::ComputeScaleOnFace(SMESH_Mesh & aMesh,
const TopoDS_Face & aFace, double &scalex, double &scaley) const TopoDS_Face & aFace,
double & scalex,
double & scaley)
{ {
//MESSAGE("StdMeshers_MEFISTO_2D::ComputeScaleOnFace"); TopoDS_Wire W = BRepTools::OuterWire(aFace);
TopoDS_Face F = TopoDS::Face(aFace.Oriented(TopAbs_FORWARD));
TopoDS_Wire W = BRepTools::OuterWire(F);
double xmin = 1.e300; // min & max of face 2D parametric coord. double xmin = 1.e300; // min & max of face 2D parametric coord.
double xmax = -1.e300; double xmax = -1.e300;
@ -615,7 +574,7 @@ void StdMeshers_MEFISTO_2D::ComputeScaleOnFace(SMESH_Mesh & aMesh,
{ {
const TopoDS_Edge & E = TopoDS::Edge( wexp.Current() ); const TopoDS_Edge & E = TopoDS::Edge( wexp.Current() );
double f, l; double f, l;
Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(E, F, f, l); Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(E, aFace, f, l);
if ( C2d.IsNull() ) continue; if ( C2d.IsNull() ) continue;
double du = (l - f) / double (nbp); double du = (l - f) / double (nbp);
for (int i = 0; i <= nbp; i++) for (int i = 0; i <= nbp; i++)
@ -642,7 +601,8 @@ void StdMeshers_MEFISTO_2D::ComputeScaleOnFace(SMESH_Mesh & aMesh,
double xsize = xmax - xmin; double xsize = xmax - xmin;
double ysize = ymax - ymin; double ysize = ymax - ymin;
Handle(Geom_Surface) S = BRep_Tool::Surface(F); // 3D surface TopLoc_Location L;
Handle(Geom_Surface) S = BRep_Tool::Surface(aFace,L); // 3D surface
double length_x = 0; double length_x = 0;
double length_y = 0; double length_y = 0;
@ -688,21 +648,21 @@ void StdMeshers_MEFISTO_2D::ComputeScaleOnFace(SMESH_Mesh & aMesh,
*/ */
//============================================================================= //=============================================================================
void StdMeshers_MEFISTO_2D::StoreResult(SMESH_Mesh & aMesh, void StdMeshers_MEFISTO_2D::StoreResult(Z nbst, R2 * uvst, Z nbt, Z * nust,
Z nbst, R2 * uvst, Z nbt, Z * nust, vector< const SMDS_MeshNode*>&mefistoToDS,
const TopoDS_Face & F, bool faceIsForward,
map<int, const SMDS_MeshNode*>&mefistoToDS,
double scalex, double scaley) double scalex, double scaley)
{ {
SMESHDS_Mesh * meshDS = aMesh.GetMeshDS(); SMESHDS_Mesh * meshDS = myTool->GetMeshDS();
int faceID = meshDS->ShapeToIndex( F ); int faceID = myTool->GetSubShapeID();
Z n, m; TopoDS_Face F = TopoDS::Face( myTool->GetSubShape() );
Handle(Geom_Surface) S = BRep_Tool::Surface(F); Handle(Geom_Surface) S = BRep_Tool::Surface( F );
for (n = 0; n < nbst; n++) Z n = mefistoToDS.size(); // nb input points
mefistoToDS.resize( nbst );
for ( ; n < nbst; n++)
{ {
if (mefistoToDS.find(n + 1) == mefistoToDS.end()) if (!mefistoToDS[n])
{ {
double u = uvst[n][0] / scalex; double u = uvst[n][0] / scalex;
double v = uvst[n][1] / scaley; double v = uvst[n][1] / scaley;
@ -712,30 +672,28 @@ void StdMeshers_MEFISTO_2D::StoreResult(SMESH_Mesh & aMesh,
meshDS->SetNodeOnFace(node, faceID, u, v); meshDS->SetNodeOnFace(node, faceID, u, v);
//MESSAGE(P.X()<<" "<<P.Y()<<" "<<P.Z()); //MESSAGE(P.X()<<" "<<P.Y()<<" "<<P.Z());
mefistoToDS[n + 1] = node; mefistoToDS[n] = node;
//MESSAGE("NEW: "<<n<<" "<<mefistoToDS[n+1]); //MESSAGE("NEW: "<<n<<" "<<mefistoToDS[n+1]);
} }
} }
m = 0; Z m = 0;
// triangle points must be in trigonometric order if face is Forward // triangle points must be in trigonometric order if face is Forward
// else they must be put clockwise // else they must be put clockwise
bool triangleIsWellOriented = faceIsForward; bool triangleIsWellOriented = ( F.Orientation() == TopAbs_FORWARD );
for (n = 1; n <= nbt; n++) for (n = 1; n <= nbt; n++)
{ {
const SMDS_MeshNode * n1 = mefistoToDS[ nust[m++] ]; const SMDS_MeshNode * n1 = mefistoToDS[ nust[m++] - 1 ];
const SMDS_MeshNode * n2 = mefistoToDS[ nust[m++] ]; const SMDS_MeshNode * n2 = mefistoToDS[ nust[m++] - 1 ];
const SMDS_MeshNode * n3 = mefistoToDS[ nust[m++] ]; const SMDS_MeshNode * n3 = mefistoToDS[ nust[m++] - 1 ];
SMDS_MeshElement * elt; SMDS_MeshElement * elt;
if (triangleIsWellOriented) if (triangleIsWellOriented)
//elt = meshDS->AddFace(n1, n2, n3);
elt = myTool->AddFace(n1, n2, n3); elt = myTool->AddFace(n1, n2, n3);
else else
//elt = meshDS->AddFace(n1, n3, n2);
elt = myTool->AddFace(n1, n3, n2); elt = myTool->AddFace(n1, n3, n2);
meshDS->SetMeshElementOnShape(elt, faceID); meshDS->SetMeshElementOnShape(elt, faceID);
@ -763,82 +721,4 @@ void StdMeshers_MEFISTO_2D::StoreResult(SMESH_Mesh & aMesh,
} }
} }
} }
}
//=============================================================================
/*!
*
*/
//=============================================================================
double StdMeshers_MEFISTO_2D::ComputeEdgeElementLength(SMESH_Mesh & aMesh,
const TopoDS_Shape & aShape)
{
//MESSAGE("StdMeshers_MEFISTO_2D::ComputeEdgeElementLength");
// **** a mettre dans SMESH_2D_Algo ?
//const TopoDS_Face & FF = TopoDS::Face(aShape);
//bool faceIsForward = (FF.Orientation() == TopAbs_FORWARD);
//TopoDS_Face F = TopoDS::Face(FF.Oriented(TopAbs_FORWARD));
double meanElementLength = 100;
double wireLength = 0;
int wireElementsNumber = 0;
for (TopExp_Explorer expe(aShape, TopAbs_EDGE); expe.More(); expe.Next())
{
const TopoDS_Edge & E = TopoDS::Edge(expe.Current());
int nb = aMesh.GetSubMesh(E)->GetSubMeshDS()->NbNodes();
double length = EdgeLength(E);
wireLength += length;
wireElementsNumber += nb;
}
if (wireElementsNumber)
meanElementLength = wireLength / wireElementsNumber;
//SCRUTE(meanElementLength);
return meanElementLength;
}
//=============================================================================
/*!
*
*/
//=============================================================================
ostream & StdMeshers_MEFISTO_2D::SaveTo(ostream & save)
{
return save;
}
//=============================================================================
/*!
*
*/
//=============================================================================
istream & StdMeshers_MEFISTO_2D::LoadFrom(istream & load)
{
return load;
}
//=============================================================================
/*!
*
*/
//=============================================================================
ostream & operator <<(ostream & save, StdMeshers_MEFISTO_2D & hyp)
{
return hyp.SaveTo( save );
}
//=============================================================================
/*!
*
*/
//=============================================================================
istream & operator >>(istream & load, StdMeshers_MEFISTO_2D & hyp)
{
return hyp.LoadFrom( load );
} }

View File

@ -31,24 +31,19 @@
#define _StdMeshers_MEFISTO_2D_HXX_ #define _StdMeshers_MEFISTO_2D_HXX_
#include "SMESH_2D_Algo.hxx" #include "SMESH_2D_Algo.hxx"
#include <TopoDS_Wire.hxx>
#include "SMESH_MesherHelper.hxx"
class SMDS_MeshNode;
class TopTools_IndexedDataMapOfShapeListOfShape;
class TopoDS_Face; class TopoDS_Face;
class TopoDS_WIre;
class StdMeshers_MaxElementArea; class StdMeshers_MaxElementArea;
class StdMeshers_LengthFromEdges; class StdMeshers_LengthFromEdges;
class SMDS_MeshNode; class SMDS_MeshNode;
class SMESH_MesherHelper;
class StdMeshers_FaceSide;
#include <vector>
#include <list> #include <list>
#include <map>
#include "Rn.h" #include "Rn.h"
class StdMeshers_MEFISTO_2D: class StdMeshers_MEFISTO_2D: public SMESH_2D_Algo
public SMESH_2D_Algo
{ {
public: public:
StdMeshers_MEFISTO_2D(int hypId, int studyId, SMESH_Gen* gen); StdMeshers_MEFISTO_2D(int hypId, int studyId, SMESH_Gen* gen);
@ -61,41 +56,29 @@ public:
virtual bool Compute(SMESH_Mesh& aMesh, virtual bool Compute(SMESH_Mesh& aMesh,
const TopoDS_Shape& aShape); const TopoDS_Shape& aShape);
double ComputeEdgeElementLength(SMESH_Mesh& aMesh, typedef boost::shared_ptr< StdMeshers_FaceSide> StdMeshers_FaceSidePtr;
const TopoDS_Shape& aShape); typedef std::vector< StdMeshers_FaceSidePtr > TWireVector;
bool LoadPoints(SMESH_Mesh& aMesh, bool LoadPoints(TWireVector & wires,
const TopoDS_Face& F, R2* uvslf,
const TopoDS_Wire& W, std::vector< const SMDS_MeshNode*>& mefistoToDS,
R2* uvslf, double scalex, double scaley);
int& m,
map<int,const SMDS_MeshNode*>& mefistoToDS,
double scalex, double scaley,
const TopTools_IndexedDataMapOfShapeListOfShape& VWMap);
void ComputeScaleOnFace(SMESH_Mesh& aMesh, void ComputeScaleOnFace(SMESH_Mesh& aMesh,
const TopoDS_Face& aFace, const TopoDS_Face& aFace,
double& scalex, double& scalex,
double& scaley); double& scaley);
void StoreResult (SMESH_Mesh& aMesh, void StoreResult (Z nbst, R2* uvst, Z nbt, Z* nust,
Z nbst, R2* uvst, Z nbt, Z* nust, std::vector< const SMDS_MeshNode*>& mefistoToDS,
const TopoDS_Face& F, bool faceIsForward,
map<int,const SMDS_MeshNode*>& mefistoToDS,
double scalex, double scaley); double scalex, double scaley);
ostream & SaveTo(ostream & save);
istream & LoadFrom(istream & load);
friend ostream & operator << (ostream & save, StdMeshers_MEFISTO_2D & hyp);
friend istream & operator >> (istream & load, StdMeshers_MEFISTO_2D & hyp);
protected: protected:
double _edgeLength; double _edgeLength;
double _maxElementArea; double _maxElementArea;
const StdMeshers_MaxElementArea* _hypMaxElementArea; const StdMeshers_MaxElementArea* _hypMaxElementArea;
const StdMeshers_LengthFromEdges* _hypLengthFromEdges; const StdMeshers_LengthFromEdges* _hypLengthFromEdges;
TopoDS_Wire myOuterWire;
std::list<const SMDS_MeshNode*> myNodesOnCommonV; std::list<const SMDS_MeshNode*> myNodesOnCommonV;
SMESH_MesherHelper* myTool; // toll for working with quadratic elements SMESH_MesherHelper* myTool; // toll for working with quadratic elements

File diff suppressed because it is too large Load Diff

View File

@ -31,40 +31,33 @@
#define _SMESH_QUADRANGLE_2D_HXX_ #define _SMESH_QUADRANGLE_2D_HXX_
#include "SMESH_2D_Algo.hxx" #include "SMESH_2D_Algo.hxx"
#include "SMESH_Mesh.hxx"
#include "Utils_SALOME_Exception.hxx" #include "Utils_SALOME_Exception.hxx"
#include "gp_XY.hxx" class SMESH_Mesh;
class SMESH_MesherHelper;
#include "SMESH_MesherHelper.hxx" class StdMeshers_FaceSide;
struct uvPtStruct;
//class SMDS_MeshNode; //class SMDS_MeshNode;
typedef struct uvPtStruct enum TSideID { BOTTOM_SIDE=0, RIGHT_SIDE, TOP_SIDE, LEFT_SIDE, NB_SIDES };
{
double param;
double normParam;
double u; // original 2d parameter
double v;
double x; // 2d parameter, normalized [0,1]
double y;
const SMDS_MeshNode * node;
} UVPtStruct;
typedef uvPtStruct UVPtStruct;
typedef struct faceQuadStruct typedef struct faceQuadStruct
{ {
int nbPts[4]; //int nbPts[4];
TopoDS_Edge edge[4]; //TopoDS_Edge edge[4];
double first[4]; StdMeshers_FaceSide* side[4];
double last[4]; //double first[4];
bool isEdgeForward[4]; //double last[4];
//bool isEdgeForward[4];
bool isEdgeOut[4]; // true, if an edge has more nodes, than the opposite bool isEdgeOut[4]; // true, if an edge has more nodes, than the opposite
UVPtStruct* uv_edges[4]; //UVPtStruct* uv_edges[4];
UVPtStruct* uv_grid; UVPtStruct* uv_grid;
~faceQuadStruct();
} FaceQuadStruct; } FaceQuadStruct;
class StdMeshers_Quadrangle_2D: class StdMeshers_Quadrangle_2D: public SMESH_2D_Algo
public SMESH_2D_Algo
{ {
public: public:
StdMeshers_Quadrangle_2D(int hypId, int studyId, SMESH_Gen* gen); StdMeshers_Quadrangle_2D(int hypId, int studyId, SMESH_Gen* gen);
@ -83,25 +76,13 @@ public:
const bool CreateQuadratic) const bool CreateQuadratic)
throw (SALOME_Exception); throw (SALOME_Exception);
static void QuadDelete(FaceQuadStruct* quad);
/**
* Returns NLinkNodeMap from myTool
*/
const NLinkNodeMap& GetNLinkNodeMap() { return myTool->GetNLinkNodeMap(); }
ostream & SaveTo(ostream & save);
istream & LoadFrom(istream & load);
friend ostream & operator << (ostream & save, StdMeshers_Quadrangle_2D & hyp);
friend istream & operator >> (istream & load, StdMeshers_Quadrangle_2D & hyp);
protected: protected:
FaceQuadStruct* CheckNbEdges(SMESH_Mesh& aMesh, FaceQuadStruct* CheckNbEdges(SMESH_Mesh& aMesh,
const TopoDS_Shape& aShape) const TopoDS_Shape& aShape)
throw (SALOME_Exception); throw (SALOME_Exception);
void SetNormalizedGrid(SMESH_Mesh& aMesh, bool SetNormalizedGrid(SMESH_Mesh& aMesh,
const TopoDS_Shape& aShape, const TopoDS_Shape& aShape,
FaceQuadStruct*& quad) FaceQuadStruct*& quad)
throw (SALOME_Exception); throw (SALOME_Exception);
@ -131,7 +112,7 @@ protected:
// is not the same in the case where the global number of nodes on edges is even // is not the same in the case where the global number of nodes on edges is even
bool myQuadranglePreference; bool myQuadranglePreference;
SMESH_MesherHelper* myTool; // toll for working with quadratic elements SMESH_MesherHelper* myTool; // tool for working with quadratic elements
}; };
#endif #endif