mirror of
https://git.salome-platform.org/gitpub/modules/smesh.git
synced 2025-01-14 10:40:34 +05:00
0022216: EDF 2613 SMESH: Projection 1D with multi-dimensional algo (Netgen 1D-2D or BLSurf...)
Improve MakeComputed to compute a source sub-mesh using an all-dimensional algo of dim greater than dim of the given sub-mesh
This commit is contained in:
parent
6b2c537fcd
commit
a541c7e70e
@ -2031,8 +2031,8 @@ TopoDS_Shape StdMeshers_ProjectionUtils::OuterShape( const TopoDS_Face& face,
|
|||||||
|
|
||||||
//================================================================================
|
//================================================================================
|
||||||
/*
|
/*
|
||||||
* Check that submesh is computed and try to compute it if is not
|
* Check that sub-mesh is computed and try to compute it if is not
|
||||||
* \param sm - submesh to compute
|
* \param sm - sub-mesh to compute
|
||||||
* \param iterationNb - int used to stop infinite recursive call
|
* \param iterationNb - int used to stop infinite recursive call
|
||||||
* \retval bool - true if computed
|
* \retval bool - true if computed
|
||||||
*/
|
*/
|
||||||
@ -2047,30 +2047,65 @@ bool StdMeshers_ProjectionUtils::MakeComputed(SMESH_subMesh * sm, const int iter
|
|||||||
if ( sm->IsMeshComputed() )
|
if ( sm->IsMeshComputed() )
|
||||||
return true;
|
return true;
|
||||||
|
|
||||||
SMESH_Mesh* mesh = sm->GetFather();
|
SMESH_Mesh* mesh = sm->GetFather();
|
||||||
SMESH_Gen* gen = mesh->GetGen();
|
SMESH_Gen* gen = mesh->GetGen();
|
||||||
SMESH_Algo* algo = sm->GetAlgo();
|
SMESH_Algo* algo = sm->GetAlgo();
|
||||||
|
TopoDS_Shape shape = sm->GetSubShape();
|
||||||
if ( !algo )
|
if ( !algo )
|
||||||
{
|
{
|
||||||
if ( sm->GetSubShape().ShapeType() != TopAbs_COMPOUND )
|
if ( shape.ShapeType() != TopAbs_COMPOUND )
|
||||||
RETURN_BAD_RESULT("No algo assigned to submesh " << sm->GetId());
|
{
|
||||||
// group
|
// No algo assigned to a non-compound sub-mesh.
|
||||||
bool computed = true;
|
// Try to find an all-dimensional algo of an upper dimension
|
||||||
for ( TopoDS_Iterator grMember( sm->GetSubShape() ); grMember.More(); grMember.Next())
|
int dim = gen->GetShapeDim( shape );
|
||||||
if ( SMESH_subMesh* grSub = mesh->GetSubMesh( grMember.Value() ))
|
for ( ++dim; ( dim <= 3 && !algo ); ++dim )
|
||||||
if ( !MakeComputed( grSub, iterationNb + 1 ))
|
{
|
||||||
computed = false;
|
SMESH_HypoFilter hypoFilter( SMESH_HypoFilter::IsAlgo() );
|
||||||
return computed;
|
hypoFilter.And( SMESH_HypoFilter::HasDim( dim ));
|
||||||
|
list <const SMESHDS_Hypothesis * > hyps;
|
||||||
|
list< TopoDS_Shape > assignedTo;
|
||||||
|
int nbAlgos =
|
||||||
|
mesh->GetHypotheses( shape, hypoFilter, hyps, true, &assignedTo );
|
||||||
|
if ( nbAlgos > 1 ) // concurrent algos
|
||||||
|
{
|
||||||
|
list<SMESH_subMesh*> smList; // where an algo is assigned
|
||||||
|
list< TopoDS_Shape >::iterator shapeIt = assignedTo.begin();
|
||||||
|
for ( ; shapeIt != assignedTo.end(); ++shapeIt )
|
||||||
|
smList.push_back( mesh->GetSubMesh( *shapeIt ));
|
||||||
|
|
||||||
|
mesh->SortByMeshOrder( smList );
|
||||||
|
algo = smList.front()->GetAlgo();
|
||||||
|
shape = smList.front()->GetSubShape();
|
||||||
|
}
|
||||||
|
else if ( nbAlgos == 1 )
|
||||||
|
{
|
||||||
|
algo = (SMESH_Algo*) hyps.front();
|
||||||
|
shape = assignedTo.front();
|
||||||
|
}
|
||||||
|
}
|
||||||
|
if ( !algo )
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
else
|
||||||
|
{
|
||||||
|
// group
|
||||||
|
bool computed = true;
|
||||||
|
for ( TopoDS_Iterator grMember( shape ); grMember.More(); grMember.Next())
|
||||||
|
if ( SMESH_subMesh* grSub = mesh->GetSubMesh( grMember.Value() ))
|
||||||
|
if ( !MakeComputed( grSub, iterationNb + 1 ))
|
||||||
|
computed = false;
|
||||||
|
return computed;
|
||||||
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
string algoType = algo->GetName();
|
string algoType = algo->GetName();
|
||||||
if ( algoType.substr(0, 11) != "Projection_")
|
if ( algoType.substr(0, 11) != "Projection_")
|
||||||
return gen->Compute( *mesh, sm->GetSubShape() );
|
return gen->Compute( *mesh, shape );
|
||||||
|
|
||||||
// try to compute source mesh
|
// try to compute source mesh
|
||||||
|
|
||||||
const list <const SMESHDS_Hypothesis *> & hyps =
|
const list <const SMESHDS_Hypothesis *> & hyps =
|
||||||
algo->GetUsedHypothesis( *mesh, sm->GetSubShape() );
|
algo->GetUsedHypothesis( *mesh, shape );
|
||||||
|
|
||||||
TopoDS_Shape srcShape;
|
TopoDS_Shape srcShape;
|
||||||
SMESH_Mesh* srcMesh = 0;
|
SMESH_Mesh* srcMesh = 0;
|
||||||
@ -2097,16 +2132,16 @@ bool StdMeshers_ProjectionUtils::MakeComputed(SMESH_subMesh * sm, const int iter
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
if ( srcShape.IsNull() ) // no projection source defined
|
if ( srcShape.IsNull() ) // no projection source defined
|
||||||
return gen->Compute( *mesh, sm->GetSubShape() );
|
return gen->Compute( *mesh, shape );
|
||||||
|
|
||||||
if ( srcShape.IsSame( sm->GetSubShape() ))
|
if ( srcShape.IsSame( shape ))
|
||||||
RETURN_BAD_RESULT("Projection from self");
|
RETURN_BAD_RESULT("Projection from self");
|
||||||
|
|
||||||
if ( !srcMesh )
|
if ( !srcMesh )
|
||||||
srcMesh = mesh;
|
srcMesh = mesh;
|
||||||
|
|
||||||
if ( MakeComputed( srcMesh->GetSubMesh( srcShape ), iterationNb + 1 ) &&
|
if ( MakeComputed( srcMesh->GetSubMesh( srcShape ), iterationNb + 1 ) &&
|
||||||
gen->Compute( *mesh, sm->GetSubShape() ))
|
gen->Compute( *mesh, shape ))
|
||||||
return sm->IsMeshComputed();
|
return sm->IsMeshComputed();
|
||||||
|
|
||||||
return false;
|
return false;
|
||||||
|
Loading…
Reference in New Issue
Block a user