mirror of
https://git.salome-platform.org/gitpub/modules/smesh.git
synced 2025-01-27 22:20:33 +05:00
Fix regressions
3D_mesh_Extrusion_00/A6 3D_mesh_Polyhedrons_00/A2,A3 groups_07/H1 imps_06/G0 imps_07/H2,H4,H5,H7 imps_08/I5
This commit is contained in:
parent
ae32dcd34f
commit
2246612bb5
@ -1003,19 +1003,20 @@ SMDS_MeshVolume* SMDS_Mesh::AddVolumeWithID(const SMDS_MeshNode * n1,
|
|||||||
}
|
}
|
||||||
else {
|
else {
|
||||||
// --- retrieve nodes ID
|
// --- retrieve nodes ID
|
||||||
|
myNodeIds.resize(12);
|
||||||
myNodeIds[0] = n1->getVtkId();
|
myNodeIds[0] = n1->getVtkId();
|
||||||
myNodeIds[0] = n6->getVtkId();
|
myNodeIds[1] = n6->getVtkId();
|
||||||
myNodeIds[0] = n5->getVtkId();
|
myNodeIds[2] = n5->getVtkId();
|
||||||
myNodeIds[0] = n4->getVtkId();
|
myNodeIds[3] = n4->getVtkId();
|
||||||
myNodeIds[0] = n3->getVtkId();
|
myNodeIds[4] = n3->getVtkId();
|
||||||
myNodeIds[0] = n2->getVtkId();
|
myNodeIds[5] = n2->getVtkId();
|
||||||
|
|
||||||
myNodeIds[0] = n7->getVtkId();
|
myNodeIds[6] = n7->getVtkId();
|
||||||
myNodeIds[0] = n12->getVtkId();
|
myNodeIds[7] = n12->getVtkId();
|
||||||
myNodeIds[0] = n11->getVtkId();
|
myNodeIds[8] = n11->getVtkId();
|
||||||
myNodeIds[0] = n10->getVtkId();
|
myNodeIds[9] = n10->getVtkId();
|
||||||
myNodeIds[0] = n9->getVtkId();
|
myNodeIds[10] = n9->getVtkId();
|
||||||
myNodeIds[0] = n8->getVtkId();
|
myNodeIds[11] = n8->getVtkId();
|
||||||
|
|
||||||
SMDS_VtkVolume *volvtk = myVolumePool->getNew();
|
SMDS_VtkVolume *volvtk = myVolumePool->getNew();
|
||||||
volvtk->init(myNodeIds, this);
|
volvtk->init(myNodeIds, this);
|
||||||
|
@ -3020,6 +3020,7 @@ CORBA::Boolean Filter_i::SetCriteria( const SMESH::Filter::Criteria& theCriteria
|
|||||||
|
|
||||||
SMESH::Predicate_ptr aPrevPredicate = SMESH::Predicate::_nil();
|
SMESH::Predicate_ptr aPrevPredicate = SMESH::Predicate::_nil();
|
||||||
int aPrevBinary = SMESH::FT_Undefined;
|
int aPrevBinary = SMESH::FT_Undefined;
|
||||||
|
aBinaries.back() = SMESH::FT_Undefined;
|
||||||
|
|
||||||
for ( aPredIter = aPredicates.begin(), aBinaryIter = aBinaries.begin();
|
for ( aPredIter = aPredicates.begin(), aBinaryIter = aBinaries.begin();
|
||||||
aPredIter != aPredicates.end() && aBinaryIter != aBinaries.end();
|
aPredIter != aPredicates.end() && aBinaryIter != aBinaries.end();
|
||||||
|
@ -307,33 +307,56 @@ namespace MeshEditor_I {
|
|||||||
*/
|
*/
|
||||||
//================================================================================
|
//================================================================================
|
||||||
|
|
||||||
|
enum IDSource_Error { IDSource_OK, IDSource_INVALID, IDSource_EMPTY };
|
||||||
|
|
||||||
bool idSourceToSet(SMESH::SMESH_IDSource_ptr theIDSource,
|
bool idSourceToSet(SMESH::SMESH_IDSource_ptr theIDSource,
|
||||||
const SMESHDS_Mesh* theMeshDS,
|
const SMESHDS_Mesh* theMeshDS,
|
||||||
TIDSortedElemSet& theElemSet,
|
TIDSortedElemSet& theElemSet,
|
||||||
const SMDSAbs_ElementType theType,
|
const SMDSAbs_ElementType theType,
|
||||||
const bool emptyIfIsMesh=false)
|
const bool emptyIfIsMesh = false,
|
||||||
|
IDSource_Error* error = 0)
|
||||||
|
|
||||||
{
|
{
|
||||||
if ( CORBA::is_nil( theIDSource ) )
|
if ( error ) *error = IDSource_OK;
|
||||||
return false;
|
|
||||||
if ( emptyIfIsMesh && SMESH::DownCast<SMESH_Mesh_i*>( theIDSource ))
|
|
||||||
return true;
|
|
||||||
|
|
||||||
|
if ( CORBA::is_nil( theIDSource ) )
|
||||||
|
{
|
||||||
|
if ( error ) *error = IDSource_INVALID;
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
if ( emptyIfIsMesh && SMESH::DownCast<SMESH_Mesh_i*>( theIDSource ))
|
||||||
|
{
|
||||||
|
if ( error && theMeshDS->GetMeshInfo().NbElements( theType ) == 0 )
|
||||||
|
*error = IDSource_EMPTY;
|
||||||
|
return true;
|
||||||
|
}
|
||||||
SMESH::long_array_var anIDs = theIDSource->GetIDs();
|
SMESH::long_array_var anIDs = theIDSource->GetIDs();
|
||||||
if ( anIDs->length() == 0 )
|
if ( anIDs->length() == 0 )
|
||||||
|
{
|
||||||
|
if ( error ) *error = IDSource_EMPTY;
|
||||||
return false;
|
return false;
|
||||||
|
}
|
||||||
SMESH::array_of_ElementType_var types = theIDSource->GetTypes();
|
SMESH::array_of_ElementType_var types = theIDSource->GetTypes();
|
||||||
if ( types->length() == 1 && types[0] == SMESH::NODE ) // group of nodes
|
if ( types->length() == 1 && types[0] == SMESH::NODE ) // group of nodes
|
||||||
{
|
{
|
||||||
if ( theType == SMDSAbs_All || theType == SMDSAbs_Node )
|
if ( theType == SMDSAbs_All || theType == SMDSAbs_Node )
|
||||||
|
{
|
||||||
arrayToSet( anIDs, theMeshDS, theElemSet, SMDSAbs_Node );
|
arrayToSet( anIDs, theMeshDS, theElemSet, SMDSAbs_Node );
|
||||||
|
}
|
||||||
else
|
else
|
||||||
|
{
|
||||||
|
if ( error ) *error = IDSource_INVALID;
|
||||||
return false;
|
return false;
|
||||||
|
}
|
||||||
}
|
}
|
||||||
else
|
else
|
||||||
{
|
{
|
||||||
arrayToSet( anIDs, theMeshDS, theElemSet, theType);
|
arrayToSet( anIDs, theMeshDS, theElemSet, theType);
|
||||||
return bool(anIDs->length()) == bool(theElemSet.size());
|
if ( bool(anIDs->length()) != bool(theElemSet.size()))
|
||||||
|
{
|
||||||
|
if ( error ) *error = IDSource_INVALID;
|
||||||
|
return false;
|
||||||
|
}
|
||||||
}
|
}
|
||||||
return true;
|
return true;
|
||||||
}
|
}
|
||||||
@ -1644,8 +1667,12 @@ CORBA::Long SMESH_MeshEditor_i::Reorient2D(SMESH::SMESH_IDSource_ptr the2Dgroup,
|
|||||||
|
|
||||||
TIDSortedElemSet elements;
|
TIDSortedElemSet elements;
|
||||||
prepareIdSource( the2Dgroup );
|
prepareIdSource( the2Dgroup );
|
||||||
if ( !idSourceToSet( the2Dgroup, getMeshDS(), elements, SMDSAbs_Face, /*emptyIfIsMesh=*/1))
|
IDSource_Error error;
|
||||||
return 0;//THROW_SALOME_CORBA_EXCEPTION("No faces in given group", SALOME::BAD_PARAM);
|
idSourceToSet( the2Dgroup, getMeshDS(), elements, SMDSAbs_Face, /*emptyIfIsMesh=*/1, &error );
|
||||||
|
if ( error == IDSource_EMPTY )
|
||||||
|
return 0;
|
||||||
|
if ( error == IDSource_INVALID )
|
||||||
|
THROW_SALOME_CORBA_EXCEPTION("No faces in given group", SALOME::BAD_PARAM);
|
||||||
|
|
||||||
|
|
||||||
const SMDS_MeshElement* face = 0;
|
const SMDS_MeshElement* face = 0;
|
||||||
@ -1725,8 +1752,8 @@ CORBA::Long SMESH_MeshEditor_i::Reorient2DBy3D(const SMESH::ListOfIDSources& fac
|
|||||||
|
|
||||||
TIDSortedElemSet volumes;
|
TIDSortedElemSet volumes;
|
||||||
prepareIdSource( volumeGroup );
|
prepareIdSource( volumeGroup );
|
||||||
if ( !idSourceToSet( volumeGroup, getMeshDS(), volumes, SMDSAbs_Volume, /*emptyIfIsMesh=*/1))
|
IDSource_Error volsError;
|
||||||
THROW_SALOME_CORBA_EXCEPTION("No volumes in a given object", SALOME::BAD_PARAM);
|
idSourceToSet( volumeGroup, getMeshDS(), volumes, SMDSAbs_Volume, /*emptyIfMesh=*/1, &volsError);
|
||||||
|
|
||||||
int nbReori = 0;
|
int nbReori = 0;
|
||||||
for ( size_t i = 0; i < faceGroups.length(); ++i )
|
for ( size_t i = 0; i < faceGroups.length(); ++i )
|
||||||
@ -1735,13 +1762,16 @@ CORBA::Long SMESH_MeshEditor_i::Reorient2DBy3D(const SMESH::ListOfIDSources& fac
|
|||||||
prepareIdSource( faceGrp );
|
prepareIdSource( faceGrp );
|
||||||
|
|
||||||
TIDSortedElemSet faces;
|
TIDSortedElemSet faces;
|
||||||
if ( !idSourceToSet( faceGrp, getMeshDS(), faces, SMDSAbs_Face, /*emptyIfIsMesh=*/1) &&
|
IDSource_Error error;
|
||||||
faceGroups.length() == 1 )
|
idSourceToSet( faceGrp, getMeshDS(), faces, SMDSAbs_Face, /*emptyIfIsMesh=*/1, &error );
|
||||||
; //THROW_SALOME_CORBA_EXCEPTION("No faces in a given object", SALOME::BAD_PARAM);
|
if ( error == IDSource_INVALID && faceGroups.length() == 1 )
|
||||||
|
THROW_SALOME_CORBA_EXCEPTION("No faces in a given object", SALOME::BAD_PARAM);
|
||||||
|
if ( error == IDSource_OK && volsError != IDSource_OK )
|
||||||
|
THROW_SALOME_CORBA_EXCEPTION("No volumes in a given object", SALOME::BAD_PARAM);
|
||||||
|
|
||||||
nbReori += getEditor().Reorient2DBy3D( faces, volumes, outsideNormal );
|
nbReori += getEditor().Reorient2DBy3D( faces, volumes, outsideNormal );
|
||||||
|
|
||||||
if ( faces.empty() ) // all faces in the mesh treated
|
if ( error != IDSource_EMPTY && faces.empty() ) // all faces in the mesh treated
|
||||||
break;
|
break;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
@ -1908,26 +1908,17 @@ class Mesh:
|
|||||||
# @return SMESH_GroupOnFilter
|
# @return SMESH_GroupOnFilter
|
||||||
# @ingroup l2_grps_create
|
# @ingroup l2_grps_create
|
||||||
def MakeGroupByCriterion(self, groupName, Criterion):
|
def MakeGroupByCriterion(self, groupName, Criterion):
|
||||||
aFilterMgr = self.smeshpyD.CreateFilterManager()
|
return self.MakeGroupByCriteria( groupName, [Criterion] )
|
||||||
aFilter = aFilterMgr.CreateFilter()
|
|
||||||
aCriteria = []
|
|
||||||
aCriteria.append(Criterion)
|
|
||||||
aFilter.SetCriteria(aCriteria)
|
|
||||||
group = self.MakeGroupByFilter(groupName, aFilter)
|
|
||||||
aFilterMgr.UnRegister()
|
|
||||||
return group
|
|
||||||
|
|
||||||
## Creates a mesh group by the given criteria (list of criteria)
|
## Creates a mesh group by the given criteria (list of criteria)
|
||||||
# @param groupName the name of the mesh group
|
# @param groupName the name of the mesh group
|
||||||
# @param theCriteria the list of criteria
|
# @param theCriteria the list of criteria
|
||||||
|
# @param binOp binary operator used when binary operator of criteria is undefined
|
||||||
# @return SMESH_GroupOnFilter
|
# @return SMESH_GroupOnFilter
|
||||||
# @ingroup l2_grps_create
|
# @ingroup l2_grps_create
|
||||||
def MakeGroupByCriteria(self, groupName, theCriteria):
|
def MakeGroupByCriteria(self, groupName, theCriteria, binOp=SMESH.FT_LogicalAND):
|
||||||
aFilterMgr = self.smeshpyD.CreateFilterManager()
|
aFilter = self.smeshpyD.GetFilterFromCriteria( theCriteria, binOp )
|
||||||
aFilter = aFilterMgr.CreateFilter()
|
|
||||||
aFilter.SetCriteria(theCriteria)
|
|
||||||
group = self.MakeGroupByFilter(groupName, aFilter)
|
group = self.MakeGroupByFilter(groupName, aFilter)
|
||||||
aFilterMgr.UnRegister()
|
|
||||||
return group
|
return group
|
||||||
|
|
||||||
## Creates a mesh group by the given filter
|
## Creates a mesh group by the given filter
|
||||||
@ -2589,7 +2580,7 @@ class Mesh:
|
|||||||
|
|
||||||
## Get measure structure specifying bounding box data of the specified object(s)
|
## Get measure structure specifying bounding box data of the specified object(s)
|
||||||
# @param IDs single source object or list of source objects or list of nodes/elements IDs
|
# @param IDs single source object or list of source objects or list of nodes/elements IDs
|
||||||
# @param isElem if @a objects is a list of IDs, @c True value in this parameters specifies that @a objects are elements,
|
# @param isElem if @a IDs is a list of IDs, @c True value in this parameters specifies that @a objects are elements,
|
||||||
# @c False specifies that @a objects are nodes
|
# @c False specifies that @a objects are nodes
|
||||||
# @return Measure structure
|
# @return Measure structure
|
||||||
# @sa BoundingBox()
|
# @sa BoundingBox()
|
||||||
|
@ -70,47 +70,6 @@ StdMeshers_Import_1D::StdMeshers_Import_1D(int hypId, int studyId, SMESH_Gen * g
|
|||||||
_compatibleHypothesis.push_back("ImportSource1D");
|
_compatibleHypothesis.push_back("ImportSource1D");
|
||||||
}
|
}
|
||||||
|
|
||||||
//=============================================================================
|
|
||||||
/*!
|
|
||||||
* Check presence of a hypothesis
|
|
||||||
*/
|
|
||||||
//=============================================================================
|
|
||||||
|
|
||||||
bool StdMeshers_Import_1D::CheckHypothesis
|
|
||||||
(SMESH_Mesh& aMesh,
|
|
||||||
const TopoDS_Shape& aShape,
|
|
||||||
SMESH_Hypothesis::Hypothesis_Status& aStatus)
|
|
||||||
{
|
|
||||||
_sourceHyp = 0;
|
|
||||||
|
|
||||||
const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
|
|
||||||
if ( hyps.size() == 0 )
|
|
||||||
{
|
|
||||||
aStatus = SMESH_Hypothesis::HYP_MISSING;
|
|
||||||
return false; // can't work with no hypothesis
|
|
||||||
}
|
|
||||||
|
|
||||||
if ( hyps.size() > 1 )
|
|
||||||
{
|
|
||||||
aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
|
|
||||||
return false;
|
|
||||||
}
|
|
||||||
|
|
||||||
const SMESHDS_Hypothesis *theHyp = hyps.front();
|
|
||||||
|
|
||||||
string hypName = theHyp->GetName();
|
|
||||||
|
|
||||||
if (hypName == _compatibleHypothesis.front())
|
|
||||||
{
|
|
||||||
_sourceHyp = (StdMeshers_ImportSource1D *)theHyp;
|
|
||||||
aStatus = SMESH_Hypothesis::HYP_OK;
|
|
||||||
return true;
|
|
||||||
}
|
|
||||||
|
|
||||||
aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
|
|
||||||
return true;
|
|
||||||
}
|
|
||||||
|
|
||||||
//================================================================================
|
//================================================================================
|
||||||
namespace // INTERNAL STUFF
|
namespace // INTERNAL STUFF
|
||||||
//================================================================================
|
//================================================================================
|
||||||
@ -346,7 +305,7 @@ namespace // INTERNAL STUFF
|
|||||||
bool rmGroups = (d->_copyGroupSubM.erase( sm ) && d->_copyGroupSubM.empty()) || rmMesh;
|
bool rmGroups = (d->_copyGroupSubM.erase( sm ) && d->_copyGroupSubM.empty()) || rmMesh;
|
||||||
if ( rmMesh )
|
if ( rmMesh )
|
||||||
d->removeImportedMesh( sm->GetFather()->GetMeshDS() );
|
d->removeImportedMesh( sm->GetFather()->GetMeshDS() );
|
||||||
if ( rmGroups && data )
|
if ( rmGroups && data && data->myType == SRC_HYP )
|
||||||
d->removeGroups( sm, data->_srcHyp );
|
d->removeGroups( sm, data->_srcHyp );
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
@ -383,7 +342,7 @@ namespace // INTERNAL STUFF
|
|||||||
// remove imported mesh and groups
|
// remove imported mesh and groups
|
||||||
d->removeImportedMesh( sm->GetFather()->GetMeshDS() );
|
d->removeImportedMesh( sm->GetFather()->GetMeshDS() );
|
||||||
|
|
||||||
if ( data )
|
if ( data && data->myType == SRC_HYP )
|
||||||
d->removeGroups( sm, data->_srcHyp );
|
d->removeGroups( sm, data->_srcHyp );
|
||||||
|
|
||||||
// clear the rest submeshes
|
// clear the rest submeshes
|
||||||
@ -395,7 +354,7 @@ namespace // INTERNAL STUFF
|
|||||||
{
|
{
|
||||||
SMESH_subMesh* subM = *sub;
|
SMESH_subMesh* subM = *sub;
|
||||||
_ListenerData* hypData = (_ListenerData*) subM->GetEventListenerData( get() );
|
_ListenerData* hypData = (_ListenerData*) subM->GetEventListenerData( get() );
|
||||||
if ( hypData )
|
if ( hypData && hypData->myType == SRC_HYP )
|
||||||
d->removeGroups( sm, hypData->_srcHyp );
|
d->removeGroups( sm, hypData->_srcHyp );
|
||||||
|
|
||||||
subM->ComputeStateEngine( SMESH_subMesh::CLEAN );
|
subM->ComputeStateEngine( SMESH_subMesh::CLEAN );
|
||||||
@ -408,7 +367,7 @@ namespace // INTERNAL STUFF
|
|||||||
if ( sm->GetSubShape().ShapeType() == TopAbs_FACE )
|
if ( sm->GetSubShape().ShapeType() == TopAbs_FACE )
|
||||||
sm->ComputeSubMeshStateEngine( SMESH_subMesh::CLEAN );
|
sm->ComputeSubMeshStateEngine( SMESH_subMesh::CLEAN );
|
||||||
}
|
}
|
||||||
if ( data )
|
if ( data && data->myType == SRC_HYP )
|
||||||
d->trackHypParams( sm, data->_srcHyp );
|
d->trackHypParams( sm, data->_srcHyp );
|
||||||
d->_n2n.clear();
|
d->_n2n.clear();
|
||||||
d->_e2e.clear();
|
d->_e2e.clear();
|
||||||
@ -633,6 +592,48 @@ namespace // INTERNAL STUFF
|
|||||||
|
|
||||||
} // namespace
|
} // namespace
|
||||||
|
|
||||||
|
//=============================================================================
|
||||||
|
/*!
|
||||||
|
* Check presence of a hypothesis
|
||||||
|
*/
|
||||||
|
//=============================================================================
|
||||||
|
|
||||||
|
bool StdMeshers_Import_1D::CheckHypothesis
|
||||||
|
(SMESH_Mesh& aMesh,
|
||||||
|
const TopoDS_Shape& aShape,
|
||||||
|
SMESH_Hypothesis::Hypothesis_Status& aStatus)
|
||||||
|
{
|
||||||
|
_sourceHyp = 0;
|
||||||
|
|
||||||
|
const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
|
||||||
|
if ( hyps.size() == 0 )
|
||||||
|
{
|
||||||
|
aStatus = SMESH_Hypothesis::HYP_MISSING;
|
||||||
|
return false; // can't work with no hypothesis
|
||||||
|
}
|
||||||
|
|
||||||
|
if ( hyps.size() > 1 )
|
||||||
|
{
|
||||||
|
aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
|
||||||
|
const SMESHDS_Hypothesis *theHyp = hyps.front();
|
||||||
|
|
||||||
|
string hypName = theHyp->GetName();
|
||||||
|
|
||||||
|
if (hypName == _compatibleHypothesis.front())
|
||||||
|
{
|
||||||
|
_sourceHyp = (StdMeshers_ImportSource1D *)theHyp;
|
||||||
|
aStatus = _sourceHyp->GetGroups().empty() ? HYP_BAD_PARAMETER : HYP_OK;
|
||||||
|
if ( aStatus == HYP_BAD_PARAMETER )
|
||||||
|
_Listener::waitHypModification( aMesh.GetSubMesh( aShape ));
|
||||||
|
return aStatus == HYP_OK;
|
||||||
|
}
|
||||||
|
|
||||||
|
aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
|
||||||
//=============================================================================
|
//=============================================================================
|
||||||
/*!
|
/*!
|
||||||
|
@ -488,10 +488,13 @@ namespace {
|
|||||||
map<const SMDS_MeshNode* , const SMDS_MeshNode*> src2tgtNodes;
|
map<const SMDS_MeshNode* , const SMDS_MeshNode*> src2tgtNodes;
|
||||||
map<const SMDS_MeshNode* , const SMDS_MeshNode*>::iterator srcN_tgtN;
|
map<const SMDS_MeshNode* , const SMDS_MeshNode*>::iterator srcN_tgtN;
|
||||||
|
|
||||||
|
bool tgtEdgesMeshed = false;
|
||||||
for ( TopExp_Explorer srcExp( srcFace, TopAbs_EDGE); srcExp.More(); srcExp.Next() )
|
for ( TopExp_Explorer srcExp( srcFace, TopAbs_EDGE); srcExp.More(); srcExp.Next() )
|
||||||
{
|
{
|
||||||
const TopoDS_Shape& srcEdge = srcExp.Current();
|
const TopoDS_Shape& srcEdge = srcExp.Current();
|
||||||
const TopoDS_Shape& tgtEdge = shape2ShapeMap( srcEdge, /*isSrc=*/true );
|
const TopoDS_Shape& tgtEdge = shape2ShapeMap( srcEdge, /*isSrc=*/true );
|
||||||
|
tgtEdgesMeshed != tgtMesh->GetSubMesh( tgtEdge )->IsEmpty();
|
||||||
|
|
||||||
if ( srcMesh->GetSubMesh( srcEdge )->IsEmpty() ||
|
if ( srcMesh->GetSubMesh( srcEdge )->IsEmpty() ||
|
||||||
tgtMesh->GetSubMesh( tgtEdge )->IsEmpty() )
|
tgtMesh->GetSubMesh( tgtEdge )->IsEmpty() )
|
||||||
continue;
|
continue;
|
||||||
@ -544,6 +547,31 @@ namespace {
|
|||||||
src2tgtNodes.insert( make_pair( u_sn->second, u_tn->second ));
|
src2tgtNodes.insert( make_pair( u_sn->second, u_tn->second ));
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
// check nodes on VERTEXes for a case of not meshes EDGEs
|
||||||
|
for ( TopExp_Explorer srcExp( srcFace, TopAbs_VERTEX); srcExp.More(); srcExp.Next() )
|
||||||
|
{
|
||||||
|
const TopoDS_Shape& srcV = srcExp.Current();
|
||||||
|
const TopoDS_Shape& tgtV = shape2ShapeMap( srcV, /*isSrc=*/true );
|
||||||
|
const SMDS_MeshNode* srcN = SMESH_Algo::VertexNode( TopoDS::Vertex( srcV ), srcMeshDS );
|
||||||
|
const SMDS_MeshNode* tgtN = SMESH_Algo::VertexNode( TopoDS::Vertex( tgtV ), srcMeshDS );
|
||||||
|
if ( !srcN )
|
||||||
|
continue;
|
||||||
|
if ( !tgtN || tgtV.ShapeType() != TopAbs_VERTEX )
|
||||||
|
return false;
|
||||||
|
|
||||||
|
if ( !tgtV.IsPartner( srcV ))
|
||||||
|
{
|
||||||
|
// check that transformation is OK by three nodes
|
||||||
|
gp_Pnt p0S = SMESH_TNodeXYZ( srcN );
|
||||||
|
gp_Pnt p0T = SMESH_TNodeXYZ( tgtN );
|
||||||
|
if ( p0T.SquareDistance( p0S.Transformed( trsf )) > tol )
|
||||||
|
{
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
src2tgtNodes.insert( make_pair( srcN, tgtN ));
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
// Make new faces
|
// Make new faces
|
||||||
|
|
||||||
@ -552,6 +580,10 @@ namespace {
|
|||||||
helper.SetSubShape( tgtFace );
|
helper.SetSubShape( tgtFace );
|
||||||
helper.IsQuadraticSubMesh( tgtFace );
|
helper.IsQuadraticSubMesh( tgtFace );
|
||||||
|
|
||||||
|
SMESHDS_SubMesh* srcSubDS = srcMeshDS->MeshElements( srcFace );
|
||||||
|
if ( !tgtEdgesMeshed && srcSubDS->NbElements() )
|
||||||
|
helper.SetIsQuadratic( srcSubDS->GetElements()->next()->IsQuadratic() );
|
||||||
|
|
||||||
SMESH_MesherHelper srcHelper( *srcMesh );
|
SMESH_MesherHelper srcHelper( *srcMesh );
|
||||||
srcHelper.SetSubShape( srcFace );
|
srcHelper.SetSubShape( srcFace );
|
||||||
|
|
||||||
@ -563,7 +595,6 @@ namespace {
|
|||||||
if ( isReverse )
|
if ( isReverse )
|
||||||
std::swap( tri1, tri2 ), std::swap( quad1, quad3 );
|
std::swap( tri1, tri2 ), std::swap( quad1, quad3 );
|
||||||
|
|
||||||
SMESHDS_SubMesh* srcSubDS = srcMeshDS->MeshElements( srcFace );
|
|
||||||
SMDS_ElemIteratorPtr elemIt = srcSubDS->GetElements();
|
SMDS_ElemIteratorPtr elemIt = srcSubDS->GetElements();
|
||||||
vector< const SMDS_MeshNode* > tgtNodes;
|
vector< const SMDS_MeshNode* > tgtNodes;
|
||||||
while ( elemIt->more() ) // loop on all mesh faces on srcFace
|
while ( elemIt->more() ) // loop on all mesh faces on srcFace
|
||||||
|
@ -226,6 +226,8 @@ namespace {
|
|||||||
static SMESH_HypoFilter hypo;
|
static SMESH_HypoFilter hypo;
|
||||||
hypo.Init( hypo.HasDim( 1 )).
|
hypo.Init( hypo.HasDim( 1 )).
|
||||||
AndNot ( hypo.IsAlgo() ).
|
AndNot ( hypo.IsAlgo() ).
|
||||||
|
AndNot ( hypo.HasName( StdMeshers_Propagation::GetName() )).
|
||||||
|
AndNot ( hypo.HasName( StdMeshers_PropagOfDistribution::GetName() )).
|
||||||
AndNot ( hypo.IsAssignedTo( theSubMesh->GetFather()->GetShapeToMesh() ));
|
AndNot ( hypo.IsAssignedTo( theSubMesh->GetFather()->GetShapeToMesh() ));
|
||||||
|
|
||||||
return theSubMesh->GetFather()->GetHypothesis( theSubMesh, hypo, true, theSssignedTo );
|
return theSubMesh->GetFather()->GetHypothesis( theSubMesh, hypo, true, theSssignedTo );
|
||||||
|
@ -52,7 +52,7 @@ class STDMESHERS_EXPORT StdMeshers_Propagation : public SMESH_Hypothesis
|
|||||||
virtual std::ostream & SaveTo(std::ostream & save);
|
virtual std::ostream & SaveTo(std::ostream & save);
|
||||||
virtual std::istream & LoadFrom(std::istream & load);
|
virtual std::istream & LoadFrom(std::istream & load);
|
||||||
|
|
||||||
static std::string GetName ();
|
static std::string GetName();
|
||||||
|
|
||||||
/*!
|
/*!
|
||||||
* \brief Returns a filter selecting both StdMeshers_Propagation and
|
* \brief Returns a filter selecting both StdMeshers_Propagation and
|
||||||
|
@ -83,7 +83,7 @@
|
|||||||
#include <cmath>
|
#include <cmath>
|
||||||
#include <limits>
|
#include <limits>
|
||||||
|
|
||||||
//#define __myDEBUG
|
#define __myDEBUG
|
||||||
|
|
||||||
using namespace std;
|
using namespace std;
|
||||||
|
|
||||||
@ -3218,7 +3218,7 @@ bool _ViscousBuilder::inflate(_SolidData& data)
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
avgThick /= data._edges.size();
|
avgThick /= data._edges.size();
|
||||||
debugMsg( "-- Thickness " << avgThick*100 << "% reached" );
|
debugMsg( "-- Thickness " << curThick << " ("<< avgThick*100 << "%) reached" );
|
||||||
|
|
||||||
if ( distToIntersection < tgtThick*avgThick*1.5 )
|
if ( distToIntersection < tgtThick*avgThick*1.5 )
|
||||||
{
|
{
|
||||||
@ -3247,7 +3247,7 @@ bool _ViscousBuilder::inflate(_SolidData& data)
|
|||||||
( new SMESH_ComputeError (COMPERR_WARNING,
|
( new SMESH_ComputeError (COMPERR_WARNING,
|
||||||
SMESH_Comment("Thickness ") << tgtThick <<
|
SMESH_Comment("Thickness ") << tgtThick <<
|
||||||
" of viscous layers not reached,"
|
" of viscous layers not reached,"
|
||||||
" average reached thickness is " << avgThick*100 << "%."));
|
" average reached thickness is " << avgThick*tgtThick));
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
Loading…
Reference in New Issue
Block a user