mirror of
https://git.salome-platform.org/gitpub/modules/smesh.git
synced 2025-02-05 05:04:16 +05:00
23092: EDF 10836 SMESH: UseExisting2DElements fails when geometry contains more than one face
Decrease 2D tolerance user for UV classification
This commit is contained in:
parent
062d52d5db
commit
22463abfc9
@ -237,6 +237,26 @@ bool SMESH_Mesh::MeshExists( int meshId ) const
|
|||||||
return _myDocument ? bool( _myDocument->GetMesh( meshId )) : false;
|
return _myDocument ? bool( _myDocument->GetMesh( meshId )) : false;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
//================================================================================
|
||||||
|
/*!
|
||||||
|
* \brief Return a mesh by id
|
||||||
|
*/
|
||||||
|
//================================================================================
|
||||||
|
|
||||||
|
SMESH_Mesh* SMESH_Mesh::FindMesh( int meshId ) const
|
||||||
|
{
|
||||||
|
if ( _id == meshId )
|
||||||
|
return (SMESH_Mesh*) this;
|
||||||
|
|
||||||
|
if ( StudyContextStruct *aStudyContext = _gen->GetStudyContext( _studyId ))
|
||||||
|
{
|
||||||
|
std::map < int, SMESH_Mesh * >::iterator i_m = aStudyContext->mapMesh.find( meshId );
|
||||||
|
if ( i_m != aStudyContext->mapMesh.end() )
|
||||||
|
return i_m->second;
|
||||||
|
}
|
||||||
|
return NULL;
|
||||||
|
}
|
||||||
|
|
||||||
//=============================================================================
|
//=============================================================================
|
||||||
/*!
|
/*!
|
||||||
* \brief Set geometry to be meshed
|
* \brief Set geometry to be meshed
|
||||||
|
@ -169,6 +169,8 @@ class SMESH_EXPORT SMESH_Mesh
|
|||||||
|
|
||||||
bool MeshExists( int meshId ) const;
|
bool MeshExists( int meshId ) const;
|
||||||
|
|
||||||
|
SMESH_Mesh* FindMesh( int meshId ) const;
|
||||||
|
|
||||||
SMESHDS_Mesh * GetMeshDS() { return _myMeshDS; }
|
SMESHDS_Mesh * GetMeshDS() { return _myMeshDS; }
|
||||||
|
|
||||||
const SMESHDS_Mesh * GetMeshDS() const { return _myMeshDS; }
|
const SMESHDS_Mesh * GetMeshDS() const { return _myMeshDS; }
|
||||||
|
@ -2470,11 +2470,8 @@ void SMESH_subMesh::deleteOwnListeners()
|
|||||||
list< OwnListenerData >::iterator d;
|
list< OwnListenerData >::iterator d;
|
||||||
for ( d = _ownListeners.begin(); d != _ownListeners.end(); ++d )
|
for ( d = _ownListeners.begin(); d != _ownListeners.end(); ++d )
|
||||||
{
|
{
|
||||||
if ( !_father->MeshExists( d->myMeshID ))
|
SMESH_Mesh* mesh = _father->FindMesh( d->myMeshID );
|
||||||
continue;
|
if ( !mesh || !mesh->GetSubMeshContaining( d->mySubMeshID ))
|
||||||
if ( _father->GetId() == d->myMeshID &&
|
|
||||||
this->GetId() != d->mySubMeshID &&
|
|
||||||
!_father->GetSubMeshContaining( d->mySubMeshID ))
|
|
||||||
continue;
|
continue;
|
||||||
d->mySubMesh->DeleteEventListener( d->myListener );
|
d->mySubMesh->DeleteEventListener( d->myListener );
|
||||||
}
|
}
|
||||||
|
@ -44,11 +44,13 @@
|
|||||||
#include "Utils_SALOME_Exception.hxx"
|
#include "Utils_SALOME_Exception.hxx"
|
||||||
#include "utilities.h"
|
#include "utilities.h"
|
||||||
|
|
||||||
|
#include <BRepBndLib.hxx>
|
||||||
#include <BRepClass_FaceClassifier.hxx>
|
#include <BRepClass_FaceClassifier.hxx>
|
||||||
#include <BRepTools.hxx>
|
#include <BRepTools.hxx>
|
||||||
#include <BRep_Builder.hxx>
|
#include <BRep_Builder.hxx>
|
||||||
#include <BRep_Tool.hxx>
|
#include <BRep_Tool.hxx>
|
||||||
#include <Bnd_B2d.hxx>
|
#include <Bnd_B2d.hxx>
|
||||||
|
#include <Bnd_Box.hxx>
|
||||||
#include <GeomAPI_ProjectPointOnSurf.hxx>
|
#include <GeomAPI_ProjectPointOnSurf.hxx>
|
||||||
#include <GeomAdaptor_Surface.hxx>
|
#include <GeomAdaptor_Surface.hxx>
|
||||||
#include <Precision.hxx>
|
#include <Precision.hxx>
|
||||||
@ -197,6 +199,7 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
|
|||||||
|
|
||||||
BRepClass_FaceClassifier classifier;
|
BRepClass_FaceClassifier classifier;
|
||||||
Bnd_B2d bndBox2d;
|
Bnd_B2d bndBox2d;
|
||||||
|
Bnd_Box bndBox3d;
|
||||||
{
|
{
|
||||||
Standard_Real umin,umax,vmin,vmax;
|
Standard_Real umin,umax,vmin,vmax;
|
||||||
BRepTools::UVBounds(geomFace,umin,umax,vmin,vmax);
|
BRepTools::UVBounds(geomFace,umin,umax,vmin,vmax);
|
||||||
@ -212,6 +215,9 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
|
|||||||
bndBox2d.Add( pmax );
|
bndBox2d.Add( pmax );
|
||||||
}
|
}
|
||||||
bndBox2d.Enlarge( 1e-2 * Sqrt( bndBox2d.SquareExtent() ));
|
bndBox2d.Enlarge( 1e-2 * Sqrt( bndBox2d.SquareExtent() ));
|
||||||
|
|
||||||
|
BRepBndLib::Add( geomFace, bndBox3d );
|
||||||
|
bndBox3d.Enlarge( 1e-4 * sqrt( bndBox3d.SquareExtent() ));
|
||||||
}
|
}
|
||||||
|
|
||||||
set<int> subShapeIDs;
|
set<int> subShapeIDs;
|
||||||
@ -284,9 +290,10 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
|
|||||||
const double groupTol = 0.5 * sqrt( getMinElemSize2( srcGroup ));
|
const double groupTol = 0.5 * sqrt( getMinElemSize2( srcGroup ));
|
||||||
minGroupTol = std::min( groupTol, minGroupTol );
|
minGroupTol = std::min( groupTol, minGroupTol );
|
||||||
|
|
||||||
GeomAdaptor_Surface S( surface );
|
//GeomAdaptor_Surface S( surface );
|
||||||
const double clsfTol = Min( S.UResolution( 0.1 * groupTol ),
|
// const double clsfTol = Min( S.UResolution( 0.1 * groupTol ), -- issue 0023092
|
||||||
S.VResolution( 0.1 * groupTol ));
|
// S.VResolution( 0.1 * groupTol ));
|
||||||
|
const double clsfTol = BRep_Tool::Tolerance( geomFace );
|
||||||
|
|
||||||
StdMeshers_Import_1D::TNodeNodeMap::iterator n2nIt;
|
StdMeshers_Import_1D::TNodeNodeMap::iterator n2nIt;
|
||||||
pair< StdMeshers_Import_1D::TNodeNodeMap::iterator, bool > it_isnew;
|
pair< StdMeshers_Import_1D::TNodeNodeMap::iterator, bool > it_isnew;
|
||||||
@ -296,13 +303,16 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
|
|||||||
{
|
{
|
||||||
const SMDS_MeshElement* face = srcElems->next();
|
const SMDS_MeshElement* face = srcElems->next();
|
||||||
|
|
||||||
|
SMDS_MeshElement::iterator node = face->begin_nodes();
|
||||||
|
if ( bndBox3d.IsOut( SMESH_TNodeXYZ( *node )))
|
||||||
|
continue;
|
||||||
|
|
||||||
// find or create nodes of a new face
|
// find or create nodes of a new face
|
||||||
nodeState.resize( face->NbNodes() );
|
nodeState.resize( face->NbNodes() );
|
||||||
newNodes.resize( nodeState.size() );
|
newNodes.resize( nodeState.size() );
|
||||||
newNodes.back() = 0;
|
newNodes.back() = 0;
|
||||||
int nbCreatedNodes = 0;
|
int nbCreatedNodes = 0;
|
||||||
bool isOut = false, isIn = false; // if at least one node isIn - do not classify other nodes
|
bool isOut = false, isIn = false; // if at least one node isIn - do not classify other nodes
|
||||||
SMDS_MeshElement::iterator node = face->begin_nodes();
|
|
||||||
for ( size_t i = 0; i < newNodes.size(); ++i, ++node )
|
for ( size_t i = 0; i < newNodes.size(); ++i, ++node )
|
||||||
{
|
{
|
||||||
SMESH_TNodeXYZ nXYZ = *node;
|
SMESH_TNodeXYZ nXYZ = *node;
|
||||||
|
Loading…
Reference in New Issue
Block a user