Fix test SALOME_TESTS/Grids/smesh/2D_mesh_QuadranglePreference_01/B6

Case of a ring with sub-meshes on both wires
This commit is contained in:
eap 2016-06-16 16:23:35 +03:00
parent aae848792f
commit 466da2436e
8 changed files with 352 additions and 209 deletions

Binary file not shown.

After

Width:  |  Height:  |  Size: 8.5 KiB

View File

@ -66,7 +66,7 @@ objects.
There is also a number of more specific algorithms:
<ul>
<li>\subpage prism_3d_algo_page "for meshing prismatic 3D shapes with hexahedra and prisms"</li>
<li>\subpage quad_from_ma_algo_page "for quadrangle meshing of faces with sinuous borders"</li>
<li>\subpage quad_from_ma_algo_page "for quadrangle meshing of faces with sinuous borders and rings"</li>
<li> <b>Polygon per Face</b> meshing algorithm - generates one mesh
face (either a triangle, a quadrangle or a polygon) per a geometrical
face using all nodes from the face boundary.</li>

View File

@ -32,4 +32,5 @@ The Medial Axis is used in two ways:
borders to find positions of nodes.</li>
</ol>
\image html quad_from_ma_ring_mesh.png "Mesh depends on defined sub-meshes: to the left - sub-meshes on both wires, to the right - a sub-mesh on internal wire only"
*/

View File

@ -1504,6 +1504,12 @@ class Mesh:
pass
return ok
## Return a list of error messages (SMESH.ComputeError) of the last Compute()
def GetComputeErrors(self, shape=0 ):
if shape == 0:
shape = self.mesh.GetShapeToMesh()
return self.smeshpyD.GetComputeErrors( self.mesh, shape )
## Return a name of a sub-shape by its ID
# @param subShapeID a unique ID of a sub-shape
# @return a string describing the sub-shape; possible variants:

View File

@ -548,8 +548,8 @@ const vector<UVPtStruct>& StdMeshers_FaceSide::SimulateUVPtStruct(int nbSeg,
bool isXConst,
double constValue) const
{
if ( myFalsePoints.empty() ) {
if ( myFalsePoints.empty() )
{
if ( NbEdges() == 0 ) return myFalsePoints;
vector<uvPtStruct>* points = const_cast<vector<uvPtStruct>*>( &myFalsePoints );
@ -557,29 +557,30 @@ const vector<UVPtStruct>& StdMeshers_FaceSide::SimulateUVPtStruct(int nbSeg,
int EdgeIndex = 0;
double prevNormPar = 0, paramSize = myNormPar[ EdgeIndex ];
for ( size_t i = 0 ; i < myFalsePoints.size(); ++i ) {
gp_Pnt2d p;
for ( size_t i = 0 ; i < myFalsePoints.size(); ++i )
{
double normPar = double(i) / double(nbSeg);
UVPtStruct & uvPt = (*points)[i];
uvPt.node = 0;
uvPt.x = uvPt.y = uvPt.param = uvPt.normParam = normPar;
if ( isXConst ) uvPt.x = constValue;
else uvPt.y = constValue;
if ( myNormPar[ EdgeIndex ] < normPar ) {
if ( myNormPar[ EdgeIndex ] < normPar )
{
prevNormPar = myNormPar[ EdgeIndex ];
++EdgeIndex;
paramSize = myNormPar[ EdgeIndex ] - prevNormPar;
}
double r = ( normPar - prevNormPar )/ paramSize;
uvPt.param = myFirst[EdgeIndex] * ( 1 - r ) + myLast[EdgeIndex] * r;
if ( !myC2d[ EdgeIndex ].IsNull() ) {
gp_Pnt2d p = myC2d[ EdgeIndex ]->Value( uvPt.param );
if ( !myC2d[ EdgeIndex ].IsNull() )
p = myC2d[ EdgeIndex ]->Value( uvPt.param );
else
p = Value2d( normPar );
uvPt.u = p.X();
uvPt.v = p.Y();
}
else {
uvPt.u = uvPt.v = 1e+100;
}
}
}
return myFalsePoints;
}

View File

@ -1363,6 +1363,116 @@ namespace
return;
} // separateNodes()
//================================================================================
/*!
* \brief Find association of nodes existing on the sinuous sides of a ring
*
* TMAPar2NPoints filled here is used in setQuadSides() only if theSinuFace.IsRing()
* to find most distant nodes of the inner and the outer wires
*/
//================================================================================
void assocNodes( SMESH_MesherHelper& theHelper,
SinuousFace& theSinuFace,
const SMESH_MAT2d::MedialAxis& theMA,
TMAPar2NPoints & thePointsOnE )
{
list< TopoDS_Edge > ee1( theSinuFace._sinuSide [0].begin(), theSinuFace._sinuSide [0].end() );
list< TopoDS_Edge > ee2( theSinuFace._sinuSide [1].begin(), theSinuFace._sinuSide [1].end() );
StdMeshers_FaceSide sideOut( theSinuFace.Face(), ee1, theHelper.GetMesh(), true, true );
StdMeshers_FaceSide sideIn ( theSinuFace.Face(), ee2, theHelper.GetMesh(), true, true );
const UVPtStructVec& uvsOut = sideOut.GetUVPtStruct();
const UVPtStructVec& uvsIn = sideIn.GetUVPtStruct();
// if ( uvs1.size() != uvs2.size() )
// return;
const SMESH_MAT2d::Branch& branch = *theMA.getBranch(0);
SMESH_MAT2d::BoundaryPoint bp[2];
SMESH_MAT2d::BranchPoint brp;
SMESHDS_Mesh* meshDS = theHelper.GetMeshDS();
map< double, const SMDS_MeshNode* > nodeParams; // params of existing nodes
map< double, const SMDS_MeshNode* >::iterator u2n;
// find a node of sideOut most distant from sideIn
vector< BRepAdaptor_Curve > curvesIn( theSinuFace._sinuSide[1].size() );
for ( size_t iE = 0; iE < theSinuFace._sinuSide[1].size(); ++iE )
curvesIn[ iE ].Initialize( theSinuFace._sinuSide[1][iE] );
double maxDist = 0;
SMESH_MAT2d::BoundaryPoint bpIn; // closest IN point
const SMDS_MeshNode* nOut = 0;
const size_t nbEOut = theSinuFace._sinuSide[0].size();
for ( size_t iE = 0; iE < nbEOut; ++iE )
{
const TopoDS_Edge& E = theSinuFace._sinuSide[0][iE];
if ( !SMESH_Algo::GetSortedNodesOnEdge( meshDS, E, /*skipMedium=*/true, nodeParams ))
return;
for ( u2n = nodeParams.begin(); u2n != nodeParams.end(); ++u2n )
{
// point on EDGE (u2n) --> MA point (brp)
if ( !theMA.getBoundary().getBranchPoint( iE, u2n->first, brp ) ||
!branch.getBoundaryPoints( brp, bp[0], bp[1] ))
return;
gp_Pnt pOut = SMESH_TNodeXYZ( u2n->second );
gp_Pnt pIn = curvesIn[ bp[1]._edgeIndex - nbEOut ].Value( bp[1]._param );
double dist = pOut.SquareDistance( pIn );
if ( dist > maxDist )
{
maxDist = dist;
nOut = u2n->second;
bpIn = bp[1];
}
}
}
const SMDS_MeshNode* nIn = 0;
if ( !SMESH_Algo::GetSortedNodesOnEdge( meshDS,
theSinuFace._sinuEdges[ bpIn._edgeIndex ],
/*skipMedium=*/true,
nodeParams ))
return;
u2n = nodeParams.lower_bound( bpIn._param );
if ( u2n == nodeParams.end() )
nIn = nodeParams.rbegin()->second;
else
nIn = u2n->second;
// find position of distant nodes in uvsOut and uvsIn
size_t iDistOut, iDistIn;
for ( iDistOut = 0; iDistOut < uvsOut.size(); ++iDistOut )
{
if ( uvsOut[iDistOut].node == nOut )
break;
}
for ( iDistIn = 0; iDistIn < uvsIn.size(); ++iDistIn )
{
if ( uvsIn[iDistIn].node == nIn )
break;
}
if ( iDistOut == uvsOut.size() || iDistIn == uvsIn.size() )
return;
// store opposite nodes in thePointsOnE (param and EDGE have no sense)
pair< NodePoint, NodePoint > oppNodes( NodePoint( nOut, 0, 0 ), NodePoint( nIn, 0, 0));
thePointsOnE.insert( make_pair( uvsOut[ iDistOut ].normParam, oppNodes ));
int iOut = iDistOut, iIn = iDistIn;
int i, nbNodes = std::min( uvsOut.size(), uvsIn.size() );
if ( nbNodes > 5 ) nbNodes = 5;
for ( i = 0, ++iOut, --iIn; i < nbNodes; ++iOut, --iIn, ++i )
{
iOut = theHelper.WrapIndex( iOut, uvsOut.size() );
iIn = theHelper.WrapIndex( iIn, uvsIn.size() );
oppNodes.first._node = uvsOut[ iOut ].node;
oppNodes.second._node = uvsIn[ iIn ].node;
thePointsOnE.insert( make_pair( uvsOut[ iOut ].normParam, oppNodes ));
}
return;
} // assocNodes()
//================================================================================
/*!
* \brief Setup sides of SinuousFace::_quad
@ -1406,6 +1516,11 @@ namespace
if ( thePointsOnEdges.size() < 4 )
return false;
int nbOut = theFace._quad->side[ 1 ].GetUVPtStruct().size();
int nbIn = theFace._quad->side[ 3 ].GetUVPtStruct().size();
if ( nbOut == 0 || nbIn == 0 )
return false;
// find most distant opposite nodes
double maxDist = 0, dist;
TMAPar2NPoints::const_iterator u2NPdist, u2NP = thePointsOnEdges.begin();
@ -1428,6 +1543,8 @@ namespace
params );
// add a radial quad side
theHelper.SetElementsOnShape( true );
u2NP = thePointsOnEdges.begin();
const SMDS_MeshNode* nOut = u2NP->second.first._node;
const SMDS_MeshNode* nIn = u2NP->second.second._node;
@ -1457,13 +1574,10 @@ namespace
theFace._quad->side[ 0 ] = StdMeshers_FaceSide::New( uvsNew );
theFace._quad->side[ 2 ] = theFace._quad->side[ 0 ];
if ( theFace._quad->side[ 1 ].GetUVPtStruct().empty() ||
theFace._quad->side[ 3 ].GetUVPtStruct().empty() )
return false;
if ( nbIn != nbOut )
theFace._quad->side[ 2 ] = StdMeshers_FaceSide::New( uvsNew );
// assure that the outer sinuous side starts at nOut
if ( theFace._sinuSide[0].size() > 1 )
{
const UVPtStructVec& uvsOut = theFace._quad->side[ 3 ].GetUVPtStruct(); // _sinuSide[0]
size_t i; // find UVPtStruct holding nOut
@ -1485,6 +1599,7 @@ namespace
}
// rotate the IN side if opposite nodes of IN and OUT sides don't match
const SMDS_MeshNode * nIn0 = theFace._quad->side[ 1 ].First().node;
if ( nIn0 != nIn )
{
@ -1571,9 +1686,24 @@ namespace
const SMESH_MAT2d::Branch& branch = *theMA.getBranch(0);
SMESH_MAT2d::BoundaryPoint bp[2];
TMAPar2NPoints pointsOnE;
// check that computed EDGEs are opposite and equally meshed
if ( allComputed )
{
// int nbNodes[2] = { 0, 0 };
// for ( int iSide = 0; iSide < 2; ++iSide ) // loop on two sinuous sides
// nbNodes[ iSide ] += meshDS->MeshElements( theSinuFace._sinuSide[ iSide ])->NbNodes() - 1;
// if ( nbNodes[0] != nbNodes[1] )
// return false;
if ( theSinuFace.IsRing() )
assocNodes( theHelper, theSinuFace, theMA, pointsOnE );
}
else
{
vector< std::size_t > edgeIDs1, edgeIDs2; // indices in theSinuEdges
vector< SMESH_MAT2d::BranchPoint > divPoints;
if ( !allComputed )
branch.getOppositeGeomEdges( edgeIDs1, edgeIDs2, divPoints );
for ( size_t i = 0; i < edgeIDs1.size(); ++i )
@ -1584,7 +1714,7 @@ namespace
int nbNodes2 = meshDS->MeshElements(edgeIDs[ edgeIDs2[i]] )->NbNodes();
if ( nbNodes1 != nbNodes2 )
return false;
if (( i-1 >= 0 ) &&
if (( int(i)-1 >= 0 ) &&
( edgeIDs1[i-1] == edgeIDs1[i] ||
edgeIDs2[i-1] == edgeIDs2[i] ))
return false;
@ -1595,7 +1725,6 @@ namespace
}
// map (param on MA) to (parameters of nodes on a pair of theSinuEdges)
TMAPar2NPoints pointsOnE;
vector<double> maParams;
set<int> projectedEdges; // treated EDGEs which 'isComputed'
@ -1774,6 +1903,7 @@ namespace
mesh->GetSubMesh( theHelper.IthVertex( 1, theSinuEdges[i] ))
->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
}
}
// Setup sides of a quadrangle
if ( !setQuadSides( theHelper, pointsOnE, theSinuFace, the1dAlgo ))
@ -1969,21 +2099,25 @@ bool StdMeshers_QuadFromMedialAxis_1D2D::computeQuads( SMESH_MesherHelper& theHe
int nbNodesShort0 = theQuad->side[0].NbPoints();
int nbNodesShort1 = theQuad->side[2].NbPoints();
int nbNodesSinu0 = theQuad->side[1].NbPoints();
int nbNodesSinu1 = theQuad->side[3].NbPoints();
// compute UV of internal points
myQuadList.push_back( theQuad );
if ( !StdMeshers_Quadrangle_2D::setNormalizedGrid( theQuad ))
return false;
// if ( !StdMeshers_Quadrangle_2D::setNormalizedGrid( theQuad ))
// return false;
// elliptic smooth of internal points to get boundary cell normal to the boundary
bool isRing = theQuad->side[0].grid->Edge(0).IsNull();
if ( !isRing )
if ( !isRing ) {
if ( !StdMeshers_Quadrangle_2D::setNormalizedGrid( theQuad ))
return false;
ellipticSmooth( theQuad, 1 );
}
// create quadrangles
bool ok;
theHelper.SetElementsOnShape( true );
if ( nbNodesShort0 == nbNodesShort1 )
if ( nbNodesShort0 == nbNodesShort1 && nbNodesSinu0 == nbNodesSinu1 )
ok = StdMeshers_Quadrangle_2D::computeQuadDominant( *theHelper.GetMesh(),
theQuad->face, theQuad );
else

View File

@ -609,7 +609,7 @@ bool StdMeshers_Quadrangle_2D::computeQuadDominant(SMESH_Mesh& aMesh,
int stop = 0;
if ( quad->side[3].grid->Edge(0).IsNull() ) // left side is simulated one
{
// quad divided at I but not at J, as nbvertic==nbright==2
if ( nbright == 2 ) // quad divided at I but not at J (2D_mesh_QuadranglePreference_01/B1)
stop++; // we stop at a second node
}
else
@ -657,7 +657,8 @@ bool StdMeshers_Quadrangle_2D::computeQuadDominant(SMESH_Mesh& aMesh,
}
// for each node of the up edge find nearest node
// in the first row of the regular grid and link them
for ( ; i > stop; i--) {
for ( ; i > stop; i--)
{
a = uv_e2[i].node;
b = uv_e2[i - 1].node;
gp_Pnt pb = SMESH_TNodeXYZ( b );
@ -791,8 +792,8 @@ bool StdMeshers_Quadrangle_2D::computeQuadDominant(SMESH_Mesh& aMesh,
}
}
} else {
if (quad->nbNodeOut(3) && nbhoriz == 2) {
// MESSAGE("left edge is out");
if (quad->nbNodeOut(3) && nbhoriz == 2)
{
int g = nbvertic - 1; // last processed node in the grid
int stop = 0;
i = quad->side[ QUAD_LEFT_SIDE ].to-1; // nbleft - 1;