mirror of
https://git.salome-platform.org/gitpub/modules/smesh.git
synced 2024-12-28 02:10:35 +05:00
Regression of SALOME_TESTS/Grids/smesh/imps_09/K2
This commit is contained in:
parent
70e4ab70fb
commit
b1f58b701e
@ -157,6 +157,46 @@ namespace {
|
|||||||
return algo;
|
return algo;
|
||||||
}
|
}
|
||||||
};
|
};
|
||||||
|
//=======================================================================
|
||||||
|
/*!
|
||||||
|
* \brief Returns already computed EDGEs
|
||||||
|
*/
|
||||||
|
void getPrecomputedEdges( SMESH_MesherHelper& theHelper,
|
||||||
|
const TopoDS_Shape& theShape,
|
||||||
|
vector< TopoDS_Edge >& theEdges)
|
||||||
|
{
|
||||||
|
theEdges.clear();
|
||||||
|
|
||||||
|
SMESHDS_Mesh* meshDS = theHelper.GetMeshDS();
|
||||||
|
SMESHDS_SubMesh* sm;
|
||||||
|
|
||||||
|
TopTools_IndexedMapOfShape edges;
|
||||||
|
TopExp::MapShapes( theShape, TopAbs_EDGE, edges );
|
||||||
|
for ( int iE = 1; iE <= edges.Extent(); ++iE )
|
||||||
|
{
|
||||||
|
const TopoDS_Shape edge = edges( iE );
|
||||||
|
if (( ! ( sm = meshDS->MeshElements( edge ))) ||
|
||||||
|
( sm->NbElements() == 0 ))
|
||||||
|
continue;
|
||||||
|
|
||||||
|
// there must not be FACEs meshed with triangles and sharing a computed EDGE
|
||||||
|
// as the precomputed EDGEs are used for propagation other to 'vertical' EDGEs
|
||||||
|
bool faceFound = false;
|
||||||
|
PShapeIteratorPtr faceIt =
|
||||||
|
theHelper.GetAncestors( edge, *theHelper.GetMesh(), TopAbs_FACE );
|
||||||
|
while ( const TopoDS_Shape* face = faceIt->next() )
|
||||||
|
|
||||||
|
if (( sm = meshDS->MeshElements( *face )) &&
|
||||||
|
( sm->NbElements() > 0 ) &&
|
||||||
|
( !theHelper.IsSameElemGeometry( sm, SMDSGeom_QUADRANGLE ) ))
|
||||||
|
{
|
||||||
|
faceFound;
|
||||||
|
break;
|
||||||
|
}
|
||||||
|
if ( !faceFound )
|
||||||
|
theEdges.push_back( TopoDS::Edge( edge ));
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
//================================================================================
|
//================================================================================
|
||||||
/*!
|
/*!
|
||||||
@ -623,6 +663,21 @@ bool StdMeshers_Prism_3D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theSh
|
|||||||
compute( prism ));
|
compute( prism ));
|
||||||
}
|
}
|
||||||
|
|
||||||
|
// find propagation chains from already computed EDGEs
|
||||||
|
vector< TopoDS_Edge > computedEdges;
|
||||||
|
getPrecomputedEdges( helper, theShape, computedEdges );
|
||||||
|
myPropagChains = new TopTools_IndexedMapOfShape[ computedEdges.size() + 1 ];
|
||||||
|
SMESHUtils::ArrayDeleter< TopTools_IndexedMapOfShape > pcDel( myPropagChains );
|
||||||
|
for ( size_t i = 0, nb = 0; i < computedEdges.size(); ++i )
|
||||||
|
{
|
||||||
|
StdMeshers_ProjectionUtils::GetPropagationEdge( &theMesh, TopoDS_Edge(),
|
||||||
|
computedEdges[i], myPropagChains + nb );
|
||||||
|
if ( myPropagChains[ nb ].Extent() < 2 ) // an empty map is a termination sign
|
||||||
|
myPropagChains[ nb ].Clear();
|
||||||
|
else
|
||||||
|
nb++;
|
||||||
|
}
|
||||||
|
|
||||||
TopTools_MapOfShape meshedSolids;
|
TopTools_MapOfShape meshedSolids;
|
||||||
list< Prism_3D::TPrismTopo > meshedPrism;
|
list< Prism_3D::TPrismTopo > meshedPrism;
|
||||||
TopTools_ListIteratorOfListOfShape solidIt;
|
TopTools_ListIteratorOfListOfShape solidIt;
|
||||||
@ -650,7 +705,11 @@ bool StdMeshers_Prism_3D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theSh
|
|||||||
!compute( prism ))
|
!compute( prism ))
|
||||||
return false;
|
return false;
|
||||||
|
|
||||||
meshedFaces.push_front( prism.myTop );
|
SMESHDS_SubMesh* smDS = theMesh.GetMeshDS()->MeshElements( prism.myTop );
|
||||||
|
if ( !myHelper->IsSameElemGeometry( smDS, SMDSGeom_QUADRANGLE ))
|
||||||
|
{
|
||||||
|
meshedFaces.push_front( prism.myTop );
|
||||||
|
}
|
||||||
meshedPrism.push_back( prism );
|
meshedPrism.push_back( prism );
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
@ -691,6 +750,7 @@ bool StdMeshers_Prism_3D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theSh
|
|||||||
prism.myBottom = candidateF;
|
prism.myBottom = candidateF;
|
||||||
mySetErrorToSM = false;
|
mySetErrorToSM = false;
|
||||||
if ( !myHelper->IsSubShape( candidateF, prismIt->myShape3D ) &&
|
if ( !myHelper->IsSubShape( candidateF, prismIt->myShape3D ) &&
|
||||||
|
myHelper->IsSubShape( candidateF, solid ) &&
|
||||||
!myHelper->GetMesh()->GetSubMesh( candidateF )->IsMeshComputed() &&
|
!myHelper->GetMesh()->GetSubMesh( candidateF )->IsMeshComputed() &&
|
||||||
initPrism( prism, solid ) &&
|
initPrism( prism, solid ) &&
|
||||||
project2dMesh( prismIt->myBottom, candidateF))
|
project2dMesh( prismIt->myBottom, candidateF))
|
||||||
@ -698,8 +758,12 @@ bool StdMeshers_Prism_3D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theSh
|
|||||||
mySetErrorToSM = true;
|
mySetErrorToSM = true;
|
||||||
if ( !compute( prism ))
|
if ( !compute( prism ))
|
||||||
return false;
|
return false;
|
||||||
meshedFaces.push_front( prism.myTop );
|
SMESHDS_SubMesh* smDS = theMesh.GetMeshDS()->MeshElements( prism.myTop );
|
||||||
meshedFaces.push_front( prism.myBottom );
|
if ( !myHelper->IsSameElemGeometry( smDS, SMDSGeom_QUADRANGLE ))
|
||||||
|
{
|
||||||
|
meshedFaces.push_front( prism.myTop );
|
||||||
|
meshedFaces.push_front( prism.myBottom );
|
||||||
|
}
|
||||||
meshedPrism.push_back( prism );
|
meshedPrism.push_back( prism );
|
||||||
meshedSolids.Add( solid );
|
meshedSolids.Add( solid );
|
||||||
}
|
}
|
||||||
@ -1187,7 +1251,7 @@ bool StdMeshers_Prism_3D::computeWalls(const Prism_3D::TPrismTopo& thePrism)
|
|||||||
SMESHDS_Mesh* meshDS = myHelper->GetMeshDS();
|
SMESHDS_Mesh* meshDS = myHelper->GetMeshDS();
|
||||||
DBGOUT( endl << "COMPUTE Prism " << meshDS->ShapeToIndex( thePrism.myShape3D ));
|
DBGOUT( endl << "COMPUTE Prism " << meshDS->ShapeToIndex( thePrism.myShape3D ));
|
||||||
|
|
||||||
TProjction1dAlgo* projector1D = TProjction1dAlgo::instance( this );
|
TProjction1dAlgo* projector1D = TProjction1dAlgo::instance( this );
|
||||||
StdMeshers_Quadrangle_2D* quadAlgo = TQuadrangleAlgo::instance( this, myHelper );
|
StdMeshers_Quadrangle_2D* quadAlgo = TQuadrangleAlgo::instance( this, myHelper );
|
||||||
|
|
||||||
SMESH_HypoFilter hyp1dFilter( SMESH_HypoFilter::IsAlgo(),/*not=*/true);
|
SMESH_HypoFilter hyp1dFilter( SMESH_HypoFilter::IsAlgo(),/*not=*/true);
|
||||||
@ -1252,8 +1316,17 @@ bool StdMeshers_Prism_3D::computeWalls(const Prism_3D::TPrismTopo& thePrism)
|
|||||||
SMESH_subMesh* srcSM = mesh->GetSubMesh( srcE );
|
SMESH_subMesh* srcSM = mesh->GetSubMesh( srcE );
|
||||||
if ( !srcSM->IsMeshComputed() ) {
|
if ( !srcSM->IsMeshComputed() ) {
|
||||||
DBGOUT( "COMPUTE V edge " << srcSM->GetId() );
|
DBGOUT( "COMPUTE V edge " << srcSM->GetId() );
|
||||||
srcSM->ComputeSubMeshStateEngine( SMESH_subMesh::COMPUTE );
|
TopoDS_Edge prpgSrcE = findPropagationSource( srcE );
|
||||||
srcSM->ComputeStateEngine ( SMESH_subMesh::COMPUTE );
|
if ( !prpgSrcE.IsNull() ) {
|
||||||
|
srcSM->ComputeSubMeshStateEngine( SMESH_subMesh::COMPUTE );
|
||||||
|
projector1D->myHyp.SetSourceEdge( prpgSrcE );
|
||||||
|
projector1D->Compute( *mesh, srcE );
|
||||||
|
srcSM->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
|
||||||
|
}
|
||||||
|
else {
|
||||||
|
srcSM->ComputeSubMeshStateEngine( SMESH_subMesh::COMPUTE );
|
||||||
|
srcSM->ComputeStateEngine ( SMESH_subMesh::COMPUTE );
|
||||||
|
}
|
||||||
if ( !srcSM->IsMeshComputed() )
|
if ( !srcSM->IsMeshComputed() )
|
||||||
return toSM( error( "Can't compute 1D mesh" ));
|
return toSM( error( "Can't compute 1D mesh" ));
|
||||||
}
|
}
|
||||||
@ -1464,6 +1537,21 @@ bool StdMeshers_Prism_3D::computeWalls(const Prism_3D::TPrismTopo& thePrism)
|
|||||||
return true;
|
return true;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
//=======================================================================
|
||||||
|
/*!
|
||||||
|
* \brief Returns a source EDGE of propagation to a given EDGE
|
||||||
|
*/
|
||||||
|
//=======================================================================
|
||||||
|
|
||||||
|
TopoDS_Edge StdMeshers_Prism_3D::findPropagationSource( const TopoDS_Edge& E )
|
||||||
|
{
|
||||||
|
for ( size_t i = 0; !myPropagChains[i].IsEmpty(); ++i )
|
||||||
|
if ( myPropagChains[i].Contains( E ))
|
||||||
|
return TopoDS::Edge( myPropagChains[i].FindKey( 1 ));
|
||||||
|
|
||||||
|
return TopoDS_Edge();
|
||||||
|
}
|
||||||
|
|
||||||
//=======================================================================
|
//=======================================================================
|
||||||
//function : Evaluate
|
//function : Evaluate
|
||||||
//purpose :
|
//purpose :
|
||||||
|
@ -460,6 +460,11 @@ public:
|
|||||||
*/
|
*/
|
||||||
bool computeWalls(const Prism_3D::TPrismTopo& thePrism);
|
bool computeWalls(const Prism_3D::TPrismTopo& thePrism);
|
||||||
|
|
||||||
|
/*!
|
||||||
|
* \brief Returns a source EDGE of propagation to a given EDGE
|
||||||
|
*/
|
||||||
|
TopoDS_Edge findPropagationSource( const TopoDS_Edge& E );
|
||||||
|
|
||||||
/*!
|
/*!
|
||||||
* \brief Find correspondence between bottom and top nodes.
|
* \brief Find correspondence between bottom and top nodes.
|
||||||
* If elements on the bottom and top faces are topologically different,
|
* If elements on the bottom and top faces are topologically different,
|
||||||
@ -514,6 +519,8 @@ private:
|
|||||||
// (the column includes the bottom node)
|
// (the column includes the bottom node)
|
||||||
TNode2ColumnMap myBotToColumnMap;
|
TNode2ColumnMap myBotToColumnMap;
|
||||||
|
|
||||||
|
TopTools_IndexedMapOfShape* myPropagChains;
|
||||||
|
|
||||||
}; // class StdMeshers_Prism_3D
|
}; // class StdMeshers_Prism_3D
|
||||||
|
|
||||||
#endif
|
#endif
|
||||||
|
@ -1613,79 +1613,87 @@ TopoDS_Vertex StdMeshers_ProjectionUtils::GetNextVertex(const TopoDS_Edge& edg
|
|||||||
/*
|
/*
|
||||||
* Return a propagation edge
|
* Return a propagation edge
|
||||||
* \param aMesh - mesh
|
* \param aMesh - mesh
|
||||||
* \param theEdge - edge to find by propagation
|
* \param anEdge - edge to find by propagation
|
||||||
* \param fromEdge - start edge for propagation
|
* \param fromEdge - start edge for propagation
|
||||||
|
* \param chain - return, if !NULL, a propagation chain passed till
|
||||||
|
* anEdge; if anEdge.IsNull() then a full propagation chain is returned;
|
||||||
|
* fromEdge is the 1st in the chain
|
||||||
* \retval pair<int,TopoDS_Edge> - propagation step and found edge
|
* \retval pair<int,TopoDS_Edge> - propagation step and found edge
|
||||||
*/
|
*/
|
||||||
//================================================================================
|
//================================================================================
|
||||||
|
|
||||||
pair<int,TopoDS_Edge>
|
pair<int,TopoDS_Edge>
|
||||||
StdMeshers_ProjectionUtils::GetPropagationEdge( SMESH_Mesh* aMesh,
|
StdMeshers_ProjectionUtils::GetPropagationEdge( SMESH_Mesh* aMesh,
|
||||||
const TopoDS_Edge& theEdge,
|
const TopoDS_Edge& anEdge,
|
||||||
const TopoDS_Edge& fromEdge)
|
const TopoDS_Edge& fromEdge,
|
||||||
|
TopTools_IndexedMapOfShape* chain)
|
||||||
{
|
{
|
||||||
TopTools_IndexedMapOfShape aChain;
|
TopTools_IndexedMapOfShape locChain;
|
||||||
|
TopTools_IndexedMapOfShape& aChain = chain ? *chain : locChain;
|
||||||
int step = 0;
|
int step = 0;
|
||||||
|
|
||||||
|
//TopTools_IndexedMapOfShape checkedWires;
|
||||||
|
BRepTools_WireExplorer aWE;
|
||||||
|
TopoDS_Shape fourEdges[4];
|
||||||
|
|
||||||
// List of edges, added to chain on the previous cycle pass
|
// List of edges, added to chain on the previous cycle pass
|
||||||
TopTools_ListOfShape listPrevEdges;
|
TopTools_ListOfShape listPrevEdges;
|
||||||
listPrevEdges.Append(fromEdge);
|
listPrevEdges.Append( fromEdge );
|
||||||
|
aChain.Add( fromEdge );
|
||||||
|
|
||||||
// Collect all edges pass by pass
|
// Collect all edges pass by pass
|
||||||
while (listPrevEdges.Extent() > 0) {
|
while (listPrevEdges.Extent() > 0)
|
||||||
|
{
|
||||||
step++;
|
step++;
|
||||||
// List of edges, added to chain on this cycle pass
|
// List of edges, added to chain on this cycle pass
|
||||||
TopTools_ListOfShape listCurEdges;
|
TopTools_ListOfShape listCurEdges;
|
||||||
|
|
||||||
// Find the next portion of edges
|
// Find the next portion of edges
|
||||||
TopTools_ListIteratorOfListOfShape itE (listPrevEdges);
|
TopTools_ListIteratorOfListOfShape itE (listPrevEdges);
|
||||||
for (; itE.More(); itE.Next()) {
|
for (; itE.More(); itE.Next())
|
||||||
TopoDS_Shape anE = itE.Value();
|
{
|
||||||
|
const TopoDS_Shape& anE = itE.Value();
|
||||||
|
|
||||||
// Iterate on faces, having edge <anE>
|
// Iterate on faces, having edge <anE>
|
||||||
TopTools_ListIteratorOfListOfShape itA (aMesh->GetAncestors(anE));
|
TopTools_ListIteratorOfListOfShape itA (aMesh->GetAncestors(anE));
|
||||||
for (; itA.More(); itA.Next()) {
|
for (; itA.More(); itA.Next())
|
||||||
TopoDS_Shape aW = itA.Value();
|
{
|
||||||
|
const TopoDS_Shape& aW = itA.Value();
|
||||||
|
|
||||||
// There are objects of different type among the ancestors of edge
|
// There are objects of different type among the ancestors of edge
|
||||||
if (aW.ShapeType() == TopAbs_WIRE) {
|
if ( aW.ShapeType() == TopAbs_WIRE /*&& checkedWires.Add( aW )*/)
|
||||||
TopoDS_Shape anOppE;
|
{
|
||||||
|
Standard_Integer nb = 0, found = -1;
|
||||||
BRepTools_WireExplorer aWE (TopoDS::Wire(aW));
|
for ( aWE.Init( TopoDS::Wire( aW )); aWE.More(); aWE.Next() ) {
|
||||||
Standard_Integer nb = 1, found = 0;
|
if (nb+1 > 4) {
|
||||||
TopTools_Array1OfShape anEdges (1,4);
|
found = -1;
|
||||||
for (; aWE.More(); aWE.Next(), nb++) {
|
|
||||||
if (nb > 4) {
|
|
||||||
found = 0;
|
|
||||||
break;
|
break;
|
||||||
}
|
}
|
||||||
anEdges(nb) = aWE.Current();
|
fourEdges[ nb ] = aWE.Current();
|
||||||
if (anEdges(nb).IsSame(anE)) found = nb;
|
if ( aWE.Current().IsSame( anE )) found = nb;
|
||||||
|
nb++;
|
||||||
}
|
}
|
||||||
|
if (nb == 4 && found >= 0) {
|
||||||
if (nb == 5 && found > 0) {
|
|
||||||
// Quadrangle face found, get an opposite edge
|
// Quadrangle face found, get an opposite edge
|
||||||
Standard_Integer opp = found + 2;
|
TopoDS_Shape& anOppE = fourEdges[( found + 2 ) % 4 ];
|
||||||
if (opp > 4) opp -= 4;
|
|
||||||
anOppE = anEdges(opp);
|
|
||||||
|
|
||||||
// add anOppE to aChain if ...
|
// add anOppE to aChain if ...
|
||||||
if (!aChain.Contains(anOppE)) { // ... anOppE is not in aChain
|
int prevChainSize = aChain.Extent();
|
||||||
|
if ( aChain.Add(anOppE) > prevChainSize ) { // ... anOppE is not in aChain
|
||||||
// Add found edge to the chain oriented so that to
|
// Add found edge to the chain oriented so that to
|
||||||
// have it co-directed with a forward MainEdge
|
// have it co-directed with a forward MainEdge
|
||||||
TopAbs_Orientation ori = anE.Orientation();
|
TopAbs_Orientation ori = anE.Orientation();
|
||||||
if ( anEdges(opp).Orientation() == anEdges(found).Orientation() )
|
if ( anOppE.Orientation() == fourEdges[found].Orientation() )
|
||||||
ori = TopAbs::Reverse( ori );
|
ori = TopAbs::Reverse( ori );
|
||||||
anOppE.Orientation( ori );
|
anOppE.Orientation( ori );
|
||||||
if ( anOppE.IsSame( theEdge ))
|
if ( anOppE.IsSame( anEdge ))
|
||||||
return make_pair( step, TopoDS::Edge( anOppE ));
|
return make_pair( step, TopoDS::Edge( anOppE ));
|
||||||
aChain.Add(anOppE);
|
|
||||||
listCurEdges.Append(anOppE);
|
listCurEdges.Append(anOppE);
|
||||||
}
|
}
|
||||||
} // if (nb == 5 && found > 0)
|
} // if (nb == 4 && found >= 0)
|
||||||
} // if (aF.ShapeType() == TopAbs_WIRE)
|
} // if (aF.ShapeType() == TopAbs_WIRE)
|
||||||
} // for (; itF.More(); itF.Next())
|
} // loop on ancestors of anE
|
||||||
} // for (; itE.More(); itE.Next())
|
} // loop on listPrevEdges
|
||||||
|
|
||||||
listPrevEdges = listCurEdges;
|
listPrevEdges = listCurEdges;
|
||||||
} // while (listPrevEdges.Extent() > 0)
|
} // while (listPrevEdges.Extent() > 0)
|
||||||
|
@ -44,6 +44,7 @@ class SMESH_Hypothesis;
|
|||||||
class SMESH_Mesh;
|
class SMESH_Mesh;
|
||||||
class SMESH_subMesh;
|
class SMESH_subMesh;
|
||||||
class TopTools_IndexedDataMapOfShapeListOfShape;
|
class TopTools_IndexedDataMapOfShapeListOfShape;
|
||||||
|
class TopTools_IndexedMapOfShape;
|
||||||
class TopoDS_Shape;
|
class TopoDS_Shape;
|
||||||
|
|
||||||
/*!
|
/*!
|
||||||
@ -156,11 +157,14 @@ namespace StdMeshers_ProjectionUtils
|
|||||||
* \brief Return an oriented propagation edge
|
* \brief Return an oriented propagation edge
|
||||||
* \param aMesh - mesh
|
* \param aMesh - mesh
|
||||||
* \param fromEdge - start edge for propagation
|
* \param fromEdge - start edge for propagation
|
||||||
|
* \param chain - return, if provided, a propagation chain passed till
|
||||||
|
* anEdge; if anEdge.IsNull() then a full propagation chain is returned
|
||||||
* \retval pair<int,TopoDS_Edge> - propagation step and found edge
|
* \retval pair<int,TopoDS_Edge> - propagation step and found edge
|
||||||
*/
|
*/
|
||||||
std::pair<int,TopoDS_Edge> GetPropagationEdge( SMESH_Mesh* aMesh,
|
std::pair<int,TopoDS_Edge> GetPropagationEdge( SMESH_Mesh* aMesh,
|
||||||
const TopoDS_Edge& anEdge,
|
const TopoDS_Edge& anEdge,
|
||||||
const TopoDS_Edge& fromEdge);
|
const TopoDS_Edge& fromEdge,
|
||||||
|
TopTools_IndexedMapOfShape* chain=0);
|
||||||
|
|
||||||
/*!
|
/*!
|
||||||
* \brief Find corresponding nodes on two faces
|
* \brief Find corresponding nodes on two faces
|
||||||
|
Loading…
Reference in New Issue
Block a user