mirror of
https://git.salome-platform.org/gitpub/modules/smesh.git
synced 2024-12-27 09:50:34 +05:00
Re-fix regression of SALOME_TESTS/Grids/smesh/3D_mesh_Extrusion_01/B2
This commit is contained in:
parent
8a1ff9ba77
commit
7d8e88b13c
@ -82,13 +82,13 @@ public:
|
||||
inline int GetID() const { return myID; }
|
||||
|
||||
///Return the type of the current element
|
||||
virtual SMDSAbs_ElementType GetType() const = 0;
|
||||
virtual SMDSAbs_EntityType GetEntityType() const = 0;
|
||||
virtual SMDSAbs_ElementType GetType() const = 0;
|
||||
virtual SMDSAbs_EntityType GetEntityType() const = 0;
|
||||
virtual SMDSAbs_GeometryType GetGeomType() const = 0;
|
||||
virtual vtkIdType GetVtkType() const = 0;
|
||||
virtual vtkIdType GetVtkType() const = 0;
|
||||
|
||||
virtual bool IsPoly() const { return false; }
|
||||
virtual bool IsQuadratic() const;
|
||||
|
||||
virtual bool IsMediumNode(const SMDS_MeshNode* node) const;
|
||||
virtual int NbCornerNodes() const;
|
||||
|
||||
@ -143,10 +143,14 @@ public:
|
||||
*/
|
||||
virtual int GetNodeIndex( const SMDS_MeshNode* node ) const;
|
||||
|
||||
inline ShortType getMeshId() const { return myMeshId; }
|
||||
inline LongType getshapeId() const { return myShapeId; }
|
||||
inline int getIdInShape() const { return myIdInShape; }
|
||||
inline int getVtkId() const { return myVtkID; }
|
||||
inline ShortType getMeshId() const { return myMeshId; }
|
||||
inline LongType getshapeId() const { return myShapeId >> BITS_SHIFT; }
|
||||
inline int getIdInShape() const { return myIdInShape; }
|
||||
inline int getVtkId() const { return myVtkID; }
|
||||
|
||||
// mark this element; to be used in algos
|
||||
inline void setIsMarked( bool is ) const;
|
||||
inline bool isMarked() const;
|
||||
|
||||
/*!
|
||||
* \brief Filters of elements, to be used with SMDS_SetIterator
|
||||
@ -180,10 +184,10 @@ public:
|
||||
};
|
||||
|
||||
protected:
|
||||
inline void setId(int id) {myID = id; }
|
||||
inline void setShapeId(LongType shapeId) {myShapeId = shapeId; }
|
||||
inline void setIdInShape(int id) { myIdInShape = id; }
|
||||
inline void setVtkId(int vtkId) { myVtkID = vtkId; }
|
||||
inline void setId(int id) { myID = id; }
|
||||
inline void setVtkId(int vtkId) { myVtkID = vtkId; }
|
||||
inline void setIdInShape(int id) { myIdInShape = id; }
|
||||
inline void setShapeId(LongType shapeId) { myShapeId = ( shapeId << BITS_SHIFT ) | ( myShapeId & BIT_IS_MARKED ); }
|
||||
SMDS_MeshElement(int ID=-1);
|
||||
SMDS_MeshElement(int id, ShortType meshId, LongType shapeId = 0);
|
||||
virtual void init(int id = -1, ShortType meshId = -1, LongType shapeId = 0);
|
||||
@ -195,12 +199,25 @@ protected:
|
||||
int myVtkID;
|
||||
//! SMDS_Mesh identification in SMESH
|
||||
ShortType myMeshId;
|
||||
//! SubShape and SubMesh identification in SMESHDS
|
||||
//! SubShape and SubMesh identification in SMESHDS; one bit is used to mark the element
|
||||
LongType myShapeId;
|
||||
//! Element index in SMESHDS_SubMesh vector
|
||||
int myIdInShape;
|
||||
|
||||
enum Bits { // use the 1st right bit of myShapeId to set/unset a mark
|
||||
BIT_IS_MARKED = 1,
|
||||
BITS_SHIFT = 1
|
||||
};
|
||||
};
|
||||
|
||||
inline void SMDS_MeshElement::setIsMarked( bool is ) const
|
||||
{
|
||||
const_cast< SMDS_MeshElement* >( this )->myShapeId = ( myShapeId & ~BIT_IS_MARKED ) | is;
|
||||
}
|
||||
inline bool SMDS_MeshElement::isMarked() const
|
||||
{
|
||||
return myShapeId & BIT_IS_MARKED;
|
||||
}
|
||||
|
||||
// ============================================================
|
||||
/*!
|
||||
|
@ -256,7 +256,7 @@ int loadVE( const list< TopoDS_Edge > & eList,
|
||||
//purpose :
|
||||
//=======================================================================
|
||||
|
||||
SMESH_Pattern::SMESH_Pattern ()
|
||||
SMESH_Pattern::SMESH_Pattern (): myToKeepNodes(false)
|
||||
{
|
||||
}
|
||||
|
||||
@ -565,11 +565,13 @@ static bool isMeshBoundToShape(SMESHDS_Mesh * aMeshDS,
|
||||
bool SMESH_Pattern::Load (SMESH_Mesh* theMesh,
|
||||
const TopoDS_Face& theFace,
|
||||
bool theProject,
|
||||
TopoDS_Vertex the1stVertex)
|
||||
TopoDS_Vertex the1stVertex,
|
||||
bool theKeepNodes)
|
||||
{
|
||||
MESSAGE(" ::Load(face) " );
|
||||
Clear();
|
||||
myIs2D = true;
|
||||
myToKeepNodes = theKeepNodes;
|
||||
|
||||
SMESHDS_Mesh * aMeshDS = theMesh->GetMeshDS();
|
||||
SMESHDS_SubMesh * fSubMesh = aMeshDS->MeshElements( theFace );
|
||||
@ -966,6 +968,19 @@ bool SMESH_Pattern::Load (SMESH_Mesh* theMesh,
|
||||
myIsBoundaryPointsFound = true;
|
||||
}
|
||||
|
||||
if ( myToKeepNodes )
|
||||
{
|
||||
myInNodes.resize( nodePointIDMap.size() + closeNodePointIDMap.size() );
|
||||
|
||||
TNodePointIDMap::iterator nIdIt = nodePointIDMap.begin();
|
||||
for ( ; nIdIt != nodePointIDMap.end(); nIdIt++ )
|
||||
myInNodes[ nIdIt->second ] = smdsNode( nIdIt->first );
|
||||
|
||||
nIdIt = closeNodePointIDMap.begin();
|
||||
for ( ; nIdIt != closeNodePointIDMap.end(); nIdIt++ )
|
||||
myInNodes[ nIdIt->second ] = smdsNode( nIdIt->first );
|
||||
}
|
||||
|
||||
// Assure that U range is proportional to V range
|
||||
|
||||
Bnd_Box2d bndBox;
|
||||
@ -3181,11 +3196,13 @@ bool SMESH_Pattern::Apply (std::set<const SMDS_MeshVolume*> & theVolumes,
|
||||
//=======================================================================
|
||||
|
||||
bool SMESH_Pattern::Load (SMESH_Mesh* theMesh,
|
||||
const TopoDS_Shell& theBlock)
|
||||
const TopoDS_Shell& theBlock,
|
||||
bool theKeepNodes)
|
||||
{
|
||||
MESSAGE(" ::Load(volume) " );
|
||||
Clear();
|
||||
myIs2D = false;
|
||||
myToKeepNodes = theKeepNodes;
|
||||
SMESHDS_SubMesh * aSubMesh;
|
||||
|
||||
const bool isQuadMesh = theMesh->NbVolumes( ORDER_QUADRATIC );
|
||||
@ -3219,7 +3236,7 @@ bool SMESH_Pattern::Load (SMESH_Mesh* theMesh,
|
||||
SMDS_NodeIteratorPtr nIt = aSubMesh->GetNodes();
|
||||
if ( !nIt->more() ) continue;
|
||||
|
||||
// store a node and a point
|
||||
// store a node and a point
|
||||
while ( nIt->more() ) {
|
||||
const SMDS_MeshNode* node = smdsNode( nIt->next() );
|
||||
if ( isQuadMesh && SMESH_MeshEditor::IsMedium( node, SMDSAbs_Volume ))
|
||||
@ -3297,6 +3314,14 @@ bool SMESH_Pattern::Load (SMESH_Mesh* theMesh,
|
||||
|
||||
myIsBoundaryPointsFound = true;
|
||||
|
||||
if ( myToKeepNodes )
|
||||
{
|
||||
myInNodes.resize( nodePointIDMap.size() );
|
||||
TNodePointIDMap::iterator nIdIt = nodePointIDMap.begin();
|
||||
for ( ; nIdIt != nodePointIDMap.end(); nIdIt++ )
|
||||
myInNodes[ nIdIt->second ] = smdsNode( nIdIt->first );
|
||||
}
|
||||
|
||||
return setErrorCode( ERR_OK );
|
||||
}
|
||||
|
||||
@ -4144,6 +4169,9 @@ bool SMESH_Pattern::MakeMesh(SMESH_Mesh* theMesh,
|
||||
|
||||
aMeshDS->compactMesh();
|
||||
|
||||
if ( myToKeepNodes )
|
||||
myOutNodes.swap( nodesVector );
|
||||
|
||||
// const map<int,SMESHDS_SubMesh*>& sm = aMeshDS->SubMeshes();
|
||||
// map<int,SMESHDS_SubMesh*>::const_iterator i_sm = sm.begin();
|
||||
// for ( ; i_sm != sm.end(); i_sm++ )
|
||||
|
@ -71,13 +71,15 @@ class SMESH_EXPORT SMESH_Pattern {
|
||||
bool Load (SMESH_Mesh* theMesh,
|
||||
const TopoDS_Face& theFace,
|
||||
bool theProject = false,
|
||||
TopoDS_Vertex the1stVertex=TopoDS_Vertex());
|
||||
TopoDS_Vertex the1stVertex=TopoDS_Vertex(),
|
||||
bool theKeepNodes = false );
|
||||
// Create a pattern from the mesh built on <theFace>.
|
||||
// <theProject>==true makes override nodes positions
|
||||
// on <theFace> computed by mesher
|
||||
|
||||
bool Load (SMESH_Mesh* theMesh,
|
||||
const TopoDS_Shell& theBlock);
|
||||
const TopoDS_Shell& theBlock,
|
||||
bool theKeepNodes = false);
|
||||
// Create a pattern from the mesh built on <theBlock>
|
||||
|
||||
bool Save (std::ostream& theFile);
|
||||
@ -216,6 +218,11 @@ class SMESH_EXPORT SMESH_Pattern {
|
||||
{ return myElemXYZIDs.empty() || !applied ? myElemPointIDs : myElemXYZIDs; }
|
||||
// Return nodal connectivity of the elements of the pattern
|
||||
|
||||
void GetInOutNodes( std::vector< const SMDS_MeshNode* > *& inNodes,
|
||||
std::vector< const SMDS_MeshNode* > *& outNodes )
|
||||
{ inNodes = & myInNodes; outNodes = & myOutNodes; }
|
||||
// Return loaded and just created nodes
|
||||
|
||||
void DumpPoints() const;
|
||||
// Debug
|
||||
|
||||
@ -368,6 +375,10 @@ private:
|
||||
// nb of key-points in each of pattern boundaries
|
||||
std::list< int > myNbKeyPntInBoundary;
|
||||
|
||||
// nodes corresponding to myPoints
|
||||
bool myToKeepNodes; // to keep these data
|
||||
std::vector< const SMDS_MeshNode* > myInNodes; // loaded nodes
|
||||
std::vector< const SMDS_MeshNode* > myOutNodes; // created nodes
|
||||
|
||||
// to compute while applying to mesh elements, not to shapes
|
||||
|
||||
@ -377,7 +388,7 @@ private:
|
||||
std::vector<const SMDS_MeshElement*> myElements; // refined elements
|
||||
std::vector<const SMDS_MeshNode*> myOrderedNodes;
|
||||
|
||||
// elements to replace with polygon or polyhedron
|
||||
// elements to replace with polygon or polyhedron
|
||||
std::vector<const SMDS_MeshElement*> myPolyElems;
|
||||
// definitions of new poly elements
|
||||
std::list< TElemDef > myPolyElemXYZIDs;
|
||||
|
@ -41,14 +41,16 @@
|
||||
#include "SMESH_Comment.hxx"
|
||||
#include "SMESH_Gen.hxx"
|
||||
#include "SMESH_Mesh.hxx"
|
||||
#include "SMESH_MeshAlgos.hxx"
|
||||
#include "SMESH_MesherHelper.hxx"
|
||||
#include "SMESH_Pattern.hxx"
|
||||
#include "SMESH_subMesh.hxx"
|
||||
#include "SMESH_subMeshEventListener.hxx"
|
||||
|
||||
#include "utilities.h"
|
||||
#include <utilities.h>
|
||||
|
||||
#include <BRepAdaptor_Surface.hxx>
|
||||
#include <BRepMesh_Delaun.hxx>
|
||||
#include <BRep_Tool.hxx>
|
||||
#include <Bnd_B2d.hxx>
|
||||
#include <GeomAPI_ProjectPointOnSurf.hxx>
|
||||
@ -260,7 +262,7 @@ namespace {
|
||||
//================================================================================
|
||||
/*!
|
||||
* \brief find new nodes belonging to one free border of mesh on face
|
||||
* \param sm - submesh on edge or vertex containg nodes to choose from
|
||||
* \param sm - submesh on edge or vertex containing nodes to choose from
|
||||
* \param face - the face bound by the submesh
|
||||
* \param u2nodes - map to fill with nodes
|
||||
* \param seamNodes - set of found nodes
|
||||
@ -1102,7 +1104,7 @@ namespace {
|
||||
{
|
||||
SMESH_subMesh* faceSM = helper.GetMesh()->GetSubMesh( helper.GetSubShape() );
|
||||
|
||||
if ( helper.IsDistorted2D( faceSM, /*checkUV=*/false ))
|
||||
if ( helper.IsDistorted2D( faceSM, /*checkUV=*/true ))
|
||||
{
|
||||
SMESH_MeshEditor editor( helper.GetMesh() );
|
||||
SMESHDS_SubMesh* smDS = faceSM->GetSubMeshDS();
|
||||
@ -1148,6 +1150,239 @@ namespace {
|
||||
return true;
|
||||
}
|
||||
|
||||
typedef list< pair< const SMDS_MeshNode*, const BRepMesh_Triangle* > > TNodeTriaList;
|
||||
|
||||
//================================================================================
|
||||
/*!
|
||||
* \brief Add in-FACE nodes surrounding a given node to a queue
|
||||
*/
|
||||
//================================================================================
|
||||
|
||||
void addCloseNodes( const SMDS_MeshNode* srcNode,
|
||||
const BRepMesh_Triangle* bmTria,
|
||||
const int srcFaceID,
|
||||
TNodeTriaList & noTriQueue )
|
||||
{
|
||||
// find in-FACE nodes
|
||||
SMDS_ElemIteratorPtr elems = srcNode->GetInverseElementIterator(SMDSAbs_Face);
|
||||
while ( elems->more() )
|
||||
{
|
||||
const SMDS_MeshElement* elem = elems->next();
|
||||
if ( elem->getshapeId() == srcFaceID )
|
||||
{
|
||||
for ( int i = 0, nb = elem->NbNodes(); i < nb; ++i )
|
||||
{
|
||||
const SMDS_MeshNode* n = elem->GetNode( i );
|
||||
if ( !n->isMarked() )
|
||||
noTriQueue.push_back( make_pair( n, bmTria ));
|
||||
}
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
//================================================================================
|
||||
/*!
|
||||
* \brief Find a delauney triangle containing a given 2D point and return
|
||||
* barycentric coordinates within the found triangle
|
||||
*/
|
||||
//================================================================================
|
||||
|
||||
const BRepMesh_Triangle* findTriangle( const gp_XY& uv,
|
||||
const BRepMesh_Triangle* bmTria,
|
||||
Handle(BRepMesh_DataStructureOfDelaun)& triaDS,
|
||||
double bc[3] )
|
||||
{
|
||||
int nodeIDs[3];
|
||||
gp_XY nodeUVs[3];
|
||||
int linkIDs[3];
|
||||
Standard_Boolean ori[3];
|
||||
|
||||
while ( bmTria )
|
||||
{
|
||||
// check bmTria
|
||||
|
||||
triaDS->ElementNodes( *bmTria, nodeIDs );
|
||||
nodeUVs[0] = triaDS->GetNode( nodeIDs[0] ).Coord();
|
||||
nodeUVs[1] = triaDS->GetNode( nodeIDs[1] ).Coord();
|
||||
nodeUVs[2] = triaDS->GetNode( nodeIDs[2] ).Coord();
|
||||
|
||||
SMESH_MeshAlgos::GetBarycentricCoords( uv,
|
||||
nodeUVs[0], nodeUVs[1], nodeUVs[2],
|
||||
bc[0], bc[1] );
|
||||
if ( bc[0] >= 0 && bc[1] >= 0 && bc[0] + bc[1] <= 1 )
|
||||
{
|
||||
bc[2] = 1 - bc[0] - bc[1];
|
||||
return bmTria;
|
||||
}
|
||||
|
||||
// look for a neighbor triangle, which is adjacent to a link intersected
|
||||
// by a segment( triangle center -> uv )
|
||||
|
||||
gp_XY gc = ( nodeUVs[0] + nodeUVs[1] + nodeUVs[2] ) / 3.;
|
||||
gp_XY seg = uv - gc;
|
||||
|
||||
bmTria->Edges( linkIDs, ori );
|
||||
int triaID = triaDS->IndexOf( *bmTria );
|
||||
bmTria = 0;
|
||||
|
||||
for ( int i = 0; i < 3; ++i )
|
||||
{
|
||||
const BRepMesh_PairOfIndex & triIDs = triaDS->ElementsConnectedTo( linkIDs[i] );
|
||||
if ( triIDs.Extent() < 2 )
|
||||
continue; // no neighbor triangle
|
||||
|
||||
// check if a link intersects gc2uv
|
||||
const BRepMesh_Edge & link = triaDS->GetLink( linkIDs[i] );
|
||||
const BRepMesh_Vertex & n1 = triaDS->GetNode( link.FirstNode() );
|
||||
const BRepMesh_Vertex & n2 = triaDS->GetNode( link.LastNode() );
|
||||
gp_XY uv1 = n1.Coord();
|
||||
gp_XY lin = n2.Coord() - uv1; // link direction
|
||||
|
||||
double crossSegLin = seg ^ lin;
|
||||
if ( Abs( crossSegLin ) < std::numeric_limits<double>::min() )
|
||||
continue; // parallel
|
||||
|
||||
double uSeg = ( uv1 - gc ) ^ lin / crossSegLin;
|
||||
if ( 0. <= uSeg && uSeg <= 1. )
|
||||
{
|
||||
bmTria = & triaDS->GetElement( triIDs.Index( 1 + ( triIDs.Index(1) == triaID )));
|
||||
break;
|
||||
}
|
||||
}
|
||||
}
|
||||
return bmTria;
|
||||
}
|
||||
|
||||
//================================================================================
|
||||
/*!
|
||||
* \brief Morph mesh on the target face to lie within FACE boundary w/o distortion
|
||||
*
|
||||
* algo:
|
||||
* - make a CDT on the src FACE
|
||||
* - find a triangle containing a src node and get its barycentric coordinates
|
||||
* - move the node to a point with the same barycentric coordinates in a corresponding
|
||||
* tgt triangle
|
||||
*/
|
||||
//================================================================================
|
||||
|
||||
bool morph( SMESH_MesherHelper& tgtHelper,
|
||||
const TopoDS_Face& tgtFace,
|
||||
const TopoDS_Face& srcFace,
|
||||
const TSideVector& tgtWires,
|
||||
const TSideVector& srcWires,
|
||||
const TAssocTool::TNodeNodeMap& src2tgtNodes )
|
||||
{
|
||||
if ( srcWires.size() != tgtWires.size() ) return false;
|
||||
|
||||
// count boundary points
|
||||
int iP = 1, nbP = 0;
|
||||
for ( size_t iW = 0; iW < srcWires.size(); ++iW )
|
||||
nbP += srcWires[iW]->NbPoints() - 1; // 1st and last points coincide
|
||||
|
||||
// fill boundary points
|
||||
BRepMesh::Array1OfVertexOfDelaun srcVert( 1, 1 + nbP ), tgtVert( 1, 1 + nbP );
|
||||
vector< const SMDS_MeshNode* > bndSrcNodes( nbP + 1 ); bndSrcNodes[0] = 0;
|
||||
BRepMesh_Vertex v( 0, 0, BRepMesh_Frontier );
|
||||
for ( size_t iW = 0; iW < srcWires.size(); ++iW )
|
||||
{
|
||||
const UVPtStructVec& srcPnt = srcWires[iW]->GetUVPtStruct();
|
||||
const UVPtStructVec& tgtPnt = tgtWires[iW]->GetUVPtStruct();
|
||||
if ( srcPnt.size() != tgtPnt.size() ) return false;
|
||||
|
||||
for ( int i = 0, nb = srcPnt.size() - 1; i < nb; ++i, ++iP )
|
||||
{
|
||||
bndSrcNodes[ iP ] = srcPnt[i].node;
|
||||
srcPnt[i].node->setIsMarked( true );
|
||||
|
||||
v.ChangeCoord() = srcPnt[i].UV();
|
||||
srcVert( iP ) = v;
|
||||
v.ChangeCoord() = tgtPnt[i].UV();
|
||||
tgtVert( iP ) = v;
|
||||
}
|
||||
}
|
||||
// triangulate the srcFace in 2D
|
||||
BRepMesh_Delaun delauney( srcVert );
|
||||
Handle(BRepMesh_DataStructureOfDelaun) triaDS = delauney.Result();
|
||||
|
||||
Handle(ShapeAnalysis_Surface) tgtSurface = tgtHelper.GetSurface( tgtFace );
|
||||
SMESHDS_Mesh* srcMesh = srcWires[0]->GetMesh()->GetMeshDS();
|
||||
SMESHDS_Mesh* tgtMesh = tgtHelper.GetMeshDS();
|
||||
const SMDS_MeshNode *srcNode, *tgtNode;
|
||||
const BRepMesh_Triangle *bmTria;
|
||||
|
||||
// un-mark internal src nodes; later we will mark moved nodes
|
||||
SMDS_NodeIteratorPtr nIt = srcMesh->MeshElements( srcFace )->GetNodes();
|
||||
if ( !nIt || !nIt->more() ) return true;
|
||||
while ( nIt->more() )
|
||||
( srcNode = nIt->next() )->setIsMarked( false );
|
||||
|
||||
// initialize a queue of nodes with starting triangles
|
||||
const int srcFaceID = srcNode->getshapeId();
|
||||
TNodeTriaList noTriQueue;
|
||||
size_t iBndSrcN = 1;
|
||||
for ( ; iBndSrcN < bndSrcNodes.size() && noTriQueue.empty(); ++iBndSrcN )
|
||||
{
|
||||
// get a triangle
|
||||
const BRepMesh::ListOfInteger & linkIds = triaDS->LinksConnectedTo( iBndSrcN );
|
||||
const BRepMesh_PairOfIndex & triaIds = triaDS->ElementsConnectedTo( linkIds.First() );
|
||||
const BRepMesh_Triangle& tria = triaDS->GetElement( triaIds.Index(1) );
|
||||
|
||||
addCloseNodes( bndSrcNodes[ iBndSrcN ], &tria, srcFaceID, noTriQueue );
|
||||
}
|
||||
|
||||
// Move tgt nodes
|
||||
|
||||
double bc[3]; // barycentric coordinates
|
||||
int nodeIDs[3];
|
||||
bool checkUV = true;
|
||||
const SMDS_FacePosition* pos;
|
||||
|
||||
while ( !noTriQueue.empty() )
|
||||
{
|
||||
srcNode = noTriQueue.front().first;
|
||||
bmTria = noTriQueue.front().second;
|
||||
noTriQueue.pop_front();
|
||||
if ( srcNode->isMarked() )
|
||||
continue;
|
||||
srcNode->setIsMarked( true );
|
||||
|
||||
// find a delauney triangle containing the src node
|
||||
gp_XY uv = tgtHelper.GetNodeUV( srcFace, srcNode, NULL, &checkUV );
|
||||
bmTria = findTriangle( uv, bmTria, triaDS, bc );
|
||||
if ( !bmTria )
|
||||
continue;
|
||||
|
||||
// compute new coordinates for a corresponding tgt node
|
||||
gp_XY uvNew( 0., 0. ), nodeUV;
|
||||
triaDS->ElementNodes( *bmTria, nodeIDs );
|
||||
for ( int i = 0; i < 3; ++i )
|
||||
uvNew += bc[i] * tgtVert( nodeIDs[i]).Coord();
|
||||
gp_Pnt xyz = tgtSurface->Value( uvNew );
|
||||
|
||||
// find and move tgt node
|
||||
TAssocTool::TNodeNodeMap::const_iterator n2n = src2tgtNodes.find( srcNode );
|
||||
if ( n2n == src2tgtNodes.end() ) continue;
|
||||
tgtNode = n2n->second;
|
||||
tgtMesh->MoveNode( tgtNode, xyz.X(), xyz.Y(), xyz.Z() );
|
||||
|
||||
if (( pos = dynamic_cast< const SMDS_FacePosition* >( tgtNode->GetPosition() )))
|
||||
const_cast<SMDS_FacePosition*>( pos )->SetParameters( uvNew.X(), uvNew.Y() );
|
||||
|
||||
addCloseNodes( srcNode, bmTria, srcFaceID, noTriQueue );
|
||||
|
||||
// assure that all src nodes are visited
|
||||
for ( ; iBndSrcN < bndSrcNodes.size() && noTriQueue.empty(); ++iBndSrcN )
|
||||
{
|
||||
const BRepMesh::ListOfInteger & linkIds = triaDS->LinksConnectedTo( iBndSrcN );
|
||||
const BRepMesh_PairOfIndex & triaIds = triaDS->ElementsConnectedTo( linkIds.First() );
|
||||
const BRepMesh_Triangle& tria = triaDS->GetElement( triaIds.Index(1) );
|
||||
addCloseNodes( bndSrcNodes[ iBndSrcN ], &tria, srcFaceID, noTriQueue );
|
||||
}
|
||||
}
|
||||
|
||||
return true;
|
||||
}
|
||||
|
||||
//=======================================================================
|
||||
/*
|
||||
* Set initial association of VERTEXes for the case of projection
|
||||
@ -1246,23 +1481,6 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
|
||||
}
|
||||
TopoDS_Face srcFace = TopoDS::Face( shape2ShapeMap( tgtFace ).Oriented(TopAbs_FORWARD));
|
||||
|
||||
// orient faces
|
||||
// if ( srcMesh == tgtMesh )
|
||||
// {
|
||||
// TopoDS_Shape solid =
|
||||
// helper.GetCommonAncestor( srcFace, tgtFace, *tgtMesh, TopAbs_SOLID );
|
||||
// if ( !solid.IsNull() )
|
||||
// {
|
||||
// srcFace.Orientation( helper.GetSubShapeOri( solid, srcFace ));
|
||||
// tgtFace.Orientation( helper.GetSubShapeOri( solid, tgtFace ));
|
||||
// }
|
||||
// else if ( helper.NbAncestors( srcFace, *tgtMesh, TopAbs_SOLID ) == 1 &&
|
||||
// helper.NbAncestors( tgtFace, *tgtMesh, TopAbs_SOLID ) == 1 )
|
||||
// {
|
||||
// srcFace.Orientation( helper.GetSubShapeOri( tgtMesh->GetShapeToMesh(), srcFace ));
|
||||
// tgtFace.Orientation( helper.GetSubShapeOri( tgtMesh->GetShapeToMesh(), tgtFace ));
|
||||
// }
|
||||
// }
|
||||
// ----------------------------------------------
|
||||
// Assure that mesh on a source Face is computed
|
||||
// ----------------------------------------------
|
||||
@ -1461,7 +1679,7 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
|
||||
|
||||
// Load pattern from the source face
|
||||
SMESH_Pattern mapper;
|
||||
mapper.Load( srcMesh, srcFace, toProjectNodes, srcV1 );
|
||||
mapper.Load( srcMesh, srcFace, toProjectNodes, srcV1, /*keepNodes=*/true );
|
||||
if ( mapper.GetErrorCode() != SMESH_Pattern::ERR_OK )
|
||||
return error(COMPERR_BAD_INPUT_MESH,"Can't load mesh pattern from the source face");
|
||||
|
||||
@ -1485,6 +1703,15 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
|
||||
if ( mapper.GetErrorCode() != SMESH_Pattern::ERR_OK )
|
||||
return error("Can't make mesh by source mesh pattern");
|
||||
|
||||
// fill _src2tgtNodes
|
||||
std::vector< const SMDS_MeshNode* > *srcNodes, *tgtNodes;
|
||||
mapper.GetInOutNodes( srcNodes, tgtNodes );
|
||||
size_t nbN = std::min( srcNodes->size(), tgtNodes->size() );
|
||||
for ( size_t i = 0; i < nbN; ++i )
|
||||
if ( (*srcNodes)[i] && (*tgtNodes)[i] )
|
||||
_src2tgtNodes.insert( make_pair( (*srcNodes)[i], (*tgtNodes)[i] ));
|
||||
|
||||
|
||||
} // end of projection using Pattern mapping
|
||||
|
||||
{
|
||||
@ -1674,9 +1901,13 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
|
||||
// boundary, also bad face can be created if EDGEs already discretized
|
||||
// --> fix bad faces by smoothing
|
||||
// ----------------------------------------------------------------
|
||||
if ( !fixDistortedFaces( helper, tgtWires ))
|
||||
return error("Invalid mesh generated");
|
||||
if ( helper.IsDistorted2D( tgtSubMesh, /*checkUV=*/false ))
|
||||
{
|
||||
morph( helper, tgtFace, srcFace, tgtWires, srcWires, _src2tgtNodes );
|
||||
|
||||
if ( !fixDistortedFaces( helper, tgtWires ))
|
||||
return error("Invalid mesh generated");
|
||||
}
|
||||
// ---------------------------
|
||||
// Check elements orientation
|
||||
// ---------------------------
|
||||
|
Loading…
Reference in New Issue
Block a user