1 // Copyright (C) 2007-2016 CEA/DEN, EDF R&D, OPEN CASCADE
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
11 // Lesser General Public License for more details.
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
20 // File : StdMeshers_RadialQuadrangle_1D2D.cxx
23 #include "StdMeshers_RadialQuadrangle_1D2D.hxx"
25 #include "StdMeshers_NumberOfLayers.hxx"
26 #include "StdMeshers_LayerDistribution.hxx"
27 #include "StdMeshers_Regular_1D.hxx"
28 #include "StdMeshers_NumberOfSegments.hxx"
29 #include "StdMeshers_FaceSide.hxx"
31 #include "SMDS_MeshNode.hxx"
32 #include "SMESHDS_Mesh.hxx"
33 #include "SMESHDS_SubMesh.hxx"
34 #include "SMESH_Gen.hxx"
35 #include "SMESH_HypoFilter.hxx"
36 #include "SMESH_Mesh.hxx"
37 #include "SMESH_MesherHelper.hxx"
38 #include "SMESH_subMesh.hxx"
39 #include "SMESH_subMeshEventListener.hxx"
40 #include "SMESH_Block.hxx"
42 #include "utilities.h"
44 #include <BRepAdaptor_CompCurve.hxx>
45 #include <BRepAdaptor_Curve.hxx>
46 #include <BRepBuilderAPI_MakeEdge.hxx>
47 #include <BRep_Builder.hxx>
48 #include <BRep_Tool.hxx>
49 #include <GCPnts_AbscissaPoint.hxx>
50 #include <Geom2d_Line.hxx>
51 #include <Geom2d_TrimmedCurve.hxx>
52 #include <GeomAPI_ProjectPointOnSurf.hxx>
53 #include <Geom_Circle.hxx>
54 #include <Geom_Line.hxx>
55 #include <Geom_TrimmedCurve.hxx>
56 #include <ShapeFix_Edge.hxx>
57 #include <TColgp_SequenceOfPnt.hxx>
58 #include <TColgp_SequenceOfPnt2d.hxx>
60 #include <TopExp_Explorer.hxx>
61 #include <TopTools_ListIteratorOfListOfShape.hxx>
67 #define RETURN_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); return false; }
68 #define gpXYZ(n) gp_XYZ(n->X(),n->Y(),n->Z())
71 //=======================================================================
72 //function : StdMeshers_RadialQuadrangle_1D2D
74 //=======================================================================
76 StdMeshers_RadialQuadrangle_1D2D::StdMeshers_RadialQuadrangle_1D2D(int hypId,
79 :StdMeshers_Quadrangle_2D( hypId, studyId, gen )
81 _name = "RadialQuadrangle_1D2D";
82 _shapeType = (1 << TopAbs_FACE); // 1 bit per shape type
84 _compatibleHypothesis.push_back("LayerDistribution2D");
85 _compatibleHypothesis.push_back("NumberOfLayers2D");
86 _requireDiscreteBoundary = false;
87 _supportSubmeshes = true;
88 _neededLowerHyps[ 1 ] = true; // suppress warning on hiding a global 1D algo
91 myDistributionHypo = 0;
95 //================================================================================
99 //================================================================================
101 StdMeshers_RadialQuadrangle_1D2D::~StdMeshers_RadialQuadrangle_1D2D()
105 //=======================================================================
106 //function : CheckHypothesis
108 //=======================================================================
110 bool StdMeshers_RadialQuadrangle_1D2D::CheckHypothesis
112 const TopoDS_Shape& aShape,
113 SMESH_Hypothesis::Hypothesis_Status& aStatus)
117 myDistributionHypo = 0;
119 list <const SMESHDS_Hypothesis * >::const_iterator itl;
121 const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
122 if ( hyps.size() == 0 ) {
123 aStatus = SMESH_Hypothesis::HYP_OK;
124 return true; // can work with no hypothesis
127 if ( hyps.size() > 1 ) {
128 aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
132 const SMESHDS_Hypothesis *theHyp = hyps.front();
134 string hypName = theHyp->GetName();
136 if (hypName == "NumberOfLayers2D") {
137 myNbLayerHypo = static_cast<const StdMeshers_NumberOfLayers *>(theHyp);
138 aStatus = SMESH_Hypothesis::HYP_OK;
141 if (hypName == "LayerDistribution2D") {
142 myDistributionHypo = static_cast<const StdMeshers_LayerDistribution *>(theHyp);
143 aStatus = SMESH_Hypothesis::HYP_OK;
146 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
152 // ------------------------------------------------------------------------------
154 * \brief Listener used to mark edges meshed by StdMeshers_RadialQuadrangle_1D2D
156 class TEdgeMarker : public SMESH_subMeshEventListener
158 TEdgeMarker(): SMESH_subMeshEventListener(/*isDeletable=*/false,
159 "StdMeshers_RadialQuadrangle_1D2D::TEdgeMarker") {}
161 //!< Return static listener
162 static SMESH_subMeshEventListener* getListener()
164 static TEdgeMarker theEdgeMarker;
165 return &theEdgeMarker;
167 //! Clear face sumbesh if something happens on edges
168 void ProcessEvent(const int event,
170 SMESH_subMesh* edgeSubMesh,
171 EventListenerData* data,
172 const SMESH_Hypothesis* /*hyp*/)
174 if ( data && !data->mySubMeshes.empty() && eventType == SMESH_subMesh::ALGO_EVENT)
176 ASSERT( data->mySubMeshes.front() != edgeSubMesh );
177 SMESH_subMesh* faceSubMesh = data->mySubMeshes.front();
178 faceSubMesh->ComputeStateEngine( SMESH_subMesh::CLEAN );
183 //================================================================================
185 * \brief Mark an edge as computed by StdMeshers_RadialQuadrangle_1D2D
187 //================================================================================
189 void markEdgeAsComputedByMe(const TopoDS_Edge& edge, SMESH_subMesh* faceSubMesh)
191 if ( SMESH_subMesh* edgeSM = faceSubMesh->GetFather()->GetSubMeshContaining( edge ))
193 if ( !edgeSM->GetEventListenerData( TEdgeMarker::getListener() ))
194 faceSubMesh->SetEventListener( TEdgeMarker::getListener(),
195 SMESH_subMeshEventListenerData::MakeData(faceSubMesh),
200 //================================================================================
202 * \brief Return sides of the face connected in the order: aCircEdge, aLinEdge1, aLinEdge2
203 * \retval int - nb of sides
205 //================================================================================
207 int analyseFace(const TopoDS_Shape& aShape,
209 StdMeshers_FaceSidePtr& aCircSide,
210 StdMeshers_FaceSidePtr& aLinSide1,
211 StdMeshers_FaceSidePtr& aLinSide2)
213 const TopoDS_Face& face = TopoDS::Face( aShape );
214 aCircSide.reset(); aLinSide1.reset(); aLinSide2.reset();
216 list< TopoDS_Edge > edges;
217 list< int > nbEdgesInWire;
218 int nbWire = SMESH_Block::GetOrderedEdges ( face, edges, nbEdgesInWire );
219 if ( nbWire > 2 || nbEdgesInWire.front() < 1 ) return 0;
221 // remove degenerated EDGEs
222 list<TopoDS_Edge>::iterator edge = edges.begin();
223 while ( edge != edges.end() )
224 if ( SMESH_Algo::isDegenerated( *edge ))
225 edge = edges.erase( edge );
228 int nbEdges = edges.size();
230 // find VERTEXes between continues EDGEs
231 TopTools_MapOfShape contVV;
234 TopoDS_Edge ePrev = edges.back();
235 for ( edge = edges.begin(); edge != edges.end(); ++edge )
237 if ( SMESH_Algo::IsContinuous( ePrev, *edge ))
238 contVV.Add( SMESH_MesherHelper::IthVertex( 0, *edge ));
242 // make edges start from a non-continues VERTEX
243 if ( 1 < contVV.Extent() && contVV.Extent() < nbEdges )
245 while ( contVV.Contains( SMESH_MesherHelper::IthVertex( 0, edges.front() )))
246 edges.splice( edges.end(), edges, edges.begin() );
251 while ( !edges.empty() )
253 list< TopoDS_Edge > sideEdges;
254 sideEdges.splice( sideEdges.end(), edges, edges.begin() );
255 while ( !edges.empty() &&
256 contVV.Contains( SMESH_MesherHelper::IthVertex( 0, edges.front() )))
257 sideEdges.splice( sideEdges.end(), edges, edges.begin() );
259 StdMeshers_FaceSidePtr side;
261 side = StdMeshers_FaceSide::New( face, sideEdges, aMesh,
262 /*isFwd=*/true, /*skipMedium=*/ true );
263 sides.push_back( side );
266 if ( !aMesh ) // call from IsApplicable()
269 if ( sides.size() > 3 )
272 if ( nbWire == 2 && (( sides.size() != 2 ) ||
273 ( sides[0]->IsClosed() && sides[1]->IsClosed() ) ||
274 ( !sides[0]->IsClosed() && !sides[1]->IsClosed() )))
277 // detect an elliptic side
279 if ( sides.size() == 1 )
281 aCircSide = sides[0];
285 // sort sides by deviation from a straight line
286 multimap< double, int > deviation2sideInd;
287 const double nbSamples = 7;
288 for ( size_t iS = 0; iS < sides.size(); ++iS )
290 gp_Pnt pf = BRep_Tool::Pnt( sides[iS]->FirstVertex() );
291 gp_Pnt pl = BRep_Tool::Pnt( sides[iS]->LastVertex() );
293 double v1Len = v1.Magnitude();
294 if ( v1Len < std::numeric_limits< double >::min() )
296 deviation2sideInd.insert( make_pair( sides[iS]->Length(), iS )); // the side seems closed
300 for ( int i = 0; i < nbSamples; ++i )
302 gp_Pnt pi( sides[iS]->Value3d(( i + 1 ) / nbSamples ));
304 double h = 0.5 * v1.Crossed( vi ).Magnitude() / v1Len;
305 devia = Max( devia, h );
307 deviation2sideInd.insert( make_pair( devia, iS ));
309 double maxDevi = deviation2sideInd.rbegin()->first;
310 if ( maxDevi < 1e-7 && sides.size() == 3 )
312 // a triangle FACE; use a side with the most outstanding length as an elliptic one
313 deviation2sideInd.clear();
314 multimap< double, int > len2sideInd;
315 for ( size_t iS = 0; iS < sides.size(); ++iS )
316 len2sideInd.insert( make_pair( sides[iS]->Length(), iS ));
318 multimap< double, int >::iterator l2i = len2sideInd.begin();
319 double len0 = l2i->first;
320 double len1 = (++l2i)->first;
321 double len2 = (++l2i)->first;
322 if ( len1 - len0 > len2 - len1 )
323 deviation2sideInd.insert( make_pair( 0., len2sideInd.begin()->second ));
325 deviation2sideInd.insert( make_pair( 0., len2sideInd.rbegin()->second ));
328 int iCirc = deviation2sideInd.rbegin()->second;
329 aCircSide = sides[ iCirc ];
330 aLinSide1 = sides[( iCirc + 1 ) % sides.size() ];
331 if ( sides.size() > 2 )
333 aLinSide2 = sides[( iCirc + 2 ) % sides.size() ];
334 aLinSide2->Reverse(); // to be "parallel" to aLinSide1
337 if (( nbWire == 2 && aLinSide1 ) &&
338 ( aLinSide1->Edge(0).Orientation() == TopAbs_INTERNAL ) &&
339 ( aCircSide->IsClosed() ))
341 // assure that aCircSide starts at aLinSide1
342 TopoDS_Vertex v0 = aLinSide1->FirstVertex();
343 TopoDS_Vertex v1 = aLinSide1->LastVertex();
344 if ( ! aCircSide->FirstVertex().IsSame( v0 ) &&
345 ! aCircSide->FirstVertex().IsSame( v1 ))
348 for ( ; iE < aCircSide->NbEdges(); ++iE )
349 if ( aCircSide->FirstVertex(iE).IsSame( v0 ) ||
350 aCircSide->FirstVertex(iE).IsSame( v1 ))
352 if ( iE == aCircSide->NbEdges() )
356 for ( int i = 0; i < aCircSide->NbEdges(); ++i, ++iE )
357 edges.push_back( aCircSide->Edge( iE % aCircSide->NbEdges() ));
359 aCircSide = StdMeshers_FaceSide::New( face, edges, aMesh,
360 /*isFwd=*/true, /*skipMedium=*/ true );
367 //================================================================================
369 * \brief Checks if the common vertex between LinSide's lies inside the circle
371 * \return bool - false if there are 3 EDGEs and the corner is outside
373 //================================================================================
375 bool isCornerInsideCircle(const StdMeshers_FaceSidePtr& CircSide,
376 const StdMeshers_FaceSidePtr& LinSide1,
377 const StdMeshers_FaceSidePtr& LinSide2)
379 // if ( CircSide && LinSide1 && LinSide2 )
381 // Handle(Geom_Circle) aCirc = Handle(Geom_Circle)::DownCast( getCurve( CircSide ));
382 // TopoDS_Vertex aCommonV;
383 // if ( !aCirc.IsNull() &&
384 // TopExp::CommonVertex( LinSide1, LinSide2, aCommonV ))
386 // gp_Pnt aCommonP = BRep_Tool::Pnt( aCommonV );
387 // gp_Pnt aCenter = aCirc->Location();
388 // double dist = aCenter.Distance( aCommonP );
389 // return dist < 0.1 * aCirc->Radius();
395 //================================================================================
397 * \brief Create an EDGE connecting the ellipse center with the most distant point
399 * \param [in] circSide - the elliptic side
400 * \param [in] face - the FACE
401 * \param [out] circNode - a node on circSide most distant from the center
402 * \return TopoDS_Edge - the create EDGE
404 //================================================================================
406 TopoDS_Edge makeEdgeToCenter( StdMeshers_FaceSidePtr& circSide,
407 const TopoDS_Face& face,
408 const SMDS_MeshNode*& circNode)
410 // find the center and a point most distant from it
412 double maxDist = 0, normPar;
414 for ( int i = 0; i < 32; ++i )
416 double u = 0.5 * i / 32.;
417 gp_Pnt2d p1 = circSide->Value2d( u );
418 gp_Pnt2d p2 = circSide->Value2d( u + 0.5 );
419 double dist = p1.SquareDistance( p2 );
420 if ( dist > maxDist )
428 gp_XY center = 0.5 * ( uv1 + uv2 );
429 double len = 0.5 * Sqrt( maxDist );
430 bool isCirc = ( Abs( len - circSide->Value2d( 0 ).Distance( center )) < 1e-3 * len );
432 // find a node closest to the most distant point
435 const UVPtStructVec& circNodes = circSide->GetUVPtStruct();
438 double minDist = 1e100;
439 for ( size_t i = 0; i <= circNodes.size(); ++i )
441 double dist = Abs( circNodes[i].normParam - normPar );
442 if ( dist < minDist )
449 circNode = circNodes[iDist].node;
450 uv1 = circNodes[iDist].UV();
451 len = ( uv1 - center ).Modulus();
453 // make the EDGE between the most distant point and the center
455 Handle(Geom2d_Line) line = new Geom2d_Line( uv1, gp_Dir2d( center - uv1 ) );
456 Handle(Geom2d_Curve) pcu = new Geom2d_TrimmedCurve( line, 0, len );
458 Handle(Geom_Surface) surface = BRep_Tool::Surface( face );
459 TopoDS_Edge edge = BRepBuilderAPI_MakeEdge( pcu, surface, 0, len );
461 BRep_Builder().UpdateEdge( edge, pcu, face, 1e-7 );
462 ShapeFix_Edge().FixAddCurve3d( edge );
465 // assure that circSide starts at circNode
466 if ( iDist != 0 && iDist != circNodes.size()-1 )
468 // create a new circSide
469 UVPtStructVec nodesNew;
470 nodesNew.reserve( circNodes.size() );
471 nodesNew.insert( nodesNew.end(), circNodes.begin() + iDist, circNodes.end() );
472 nodesNew.insert( nodesNew.end(), circNodes.begin() + 1, circNodes.begin() + iDist + 1 );
473 circSide = StdMeshers_FaceSide::New( nodesNew );
479 //================================================================================
481 * \brief Set nodes existing on a linSide to UVPtStructVec and create missing nodes
482 * corresponding to layerPositions
484 //================================================================================
486 void makeMissingMesh( StdMeshers_FaceSidePtr& linSide,
487 UVPtStructVec& nodes,
488 const vector< double >& layerPositions,
489 SMESH_MesherHelper* helper )
491 // tolerance to compare normParam
493 for ( size_t i = 1; i < layerPositions.size(); ++i )
494 tol = Min( tol, layerPositions[i] - layerPositions[i-1] );
497 // merge existing nodes with layerPositions into UVPtStructVec
498 // ------------------------------------------------------------
500 const UVPtStructVec& exiNodes = linSide->GetUVPtStruct();
502 nodes.reserve( layerPositions.size() + exiNodes.size() );
503 vector< double >::const_iterator pos = layerPositions.begin(), posEnd = layerPositions.end();
504 for ( size_t i = 0; i < exiNodes.size(); ++i )
506 switch ( exiNodes[i].node->GetPosition()->GetTypeOfPosition() )
508 case SMDS_TOP_VERTEX:
510 // allocate UVPtStruct's for non-existing nodes
511 while ( pos != posEnd && *pos < exiNodes[i].normParam - tol )
514 uvPS.normParam = *pos++;
515 nodes.push_back( uvPS );
517 // save existing node on a VERTEX
518 nodes.push_back( exiNodes[i] );
523 // save existing nodes on an EDGE
524 while ( i < exiNodes.size() && exiNodes[i].node->GetPosition()->GetDim() == 1 )
526 nodes.push_back( exiNodes[i++] );
528 // save existing node on a VERTEX
529 if ( i < exiNodes.size() && exiNodes[i].node->GetPosition()->GetDim() == 0 )
531 nodes.push_back( exiNodes[i] );
538 // skip layer positions covered by saved nodes
539 while ( pos != posEnd && *pos < nodes.back().normParam + tol )
544 // allocate UVPtStruct's for the rest non-existing nodes
545 while ( pos != posEnd )
548 uvPS.normParam = *pos++;
549 nodes.push_back( uvPS );
552 // create missing nodes
553 // ---------------------
555 SMESHDS_Mesh * meshDS = helper->GetMeshDS();
556 const TopoDS_Face& F = TopoDS::Face( helper->GetSubShape() );
557 Handle(ShapeAnalysis_Surface) surface = helper->GetSurface( F );
559 helper->SetElementsOnShape( false ); // we create nodes on EDGEs, not on the FACE
561 for ( size_t i = 0; i < nodes.size(); ++i )
563 if ( nodes[ i ].node ) continue;
565 gp_Pnt2d uv = linSide->Value2d( nodes[i].normParam );
566 gp_Pnt xyz = surface->Value( uv.X(), uv.Y() );
568 nodes[ i ].SetUV( uv.XY() );
569 nodes[ i ].node = helper->AddNode( xyz.X(), xyz.Y(), xyz.Z() );
572 // set nodes on VERTEXes
573 for ( int iE = 0; iE < linSide->NbEdges(); ++iE )
575 TopoDS_Vertex v = linSide->LastVertex( iE );
576 if ( SMESH_Algo::VertexNode( v, meshDS ))
579 double normPar = linSide->LastParameter( iE );
581 while ( nodes[ i ].normParam < normPar )
583 if (( nodes[ i ].normParam - normPar ) > ( normPar - nodes[ i-1 ].normParam ))
585 meshDS->SetNodeOnVertex( nodes[ i ].node, v );
588 // set nodes on EDGEs
590 for ( size_t i = 0; i < nodes.size(); ++i )
592 if ( nodes[ i ].node->getshapeId() > 0 ) continue;
594 double u = linSide->Parameter( nodes[i].normParam, edgeID );
595 meshDS->SetNodeOnEdge( nodes[ i ].node, edgeID, u );
599 for ( size_t i = 1; i < nodes.size(); ++i )
601 if ( meshDS->FindEdge( nodes[i].node, nodes[i-1].node )) continue;
603 const SMDS_MeshElement* seg = helper->AddEdge( nodes[i].node, nodes[i-1].node );
605 double normParam = 0.5 * ( nodes[i].normParam + nodes[i-1].normParam );
606 edgeID = linSide->EdgeID( linSide->EdgeIndex( normParam ));
607 meshDS->SetMeshElementOnShape( seg, edgeID );
610 helper->SetElementsOnShape( true );
612 } // end makeMissingMesh()
614 //================================================================================
615 //================================================================================
617 * \brief Class computing layers distribution using data of
618 * StdMeshers_LayerDistribution hypothesis
620 //================================================================================
621 //================================================================================
623 class TNodeDistributor: public StdMeshers_Regular_1D
625 list <const SMESHDS_Hypothesis *> myUsedHyps;
627 // -----------------------------------------------------------------------------
628 static TNodeDistributor* GetDistributor(SMESH_Mesh& aMesh)
630 const int myID = -1001;
631 TNodeDistributor* myHyp = dynamic_cast<TNodeDistributor*>( aMesh.GetHypothesis( myID ));
633 myHyp = new TNodeDistributor( myID, 0, aMesh.GetGen() );
636 // -----------------------------------------------------------------------------
637 //! Computes distribution of nodes on a straight line ending at pIn and pOut
638 bool Compute( vector< double > & positions,
639 const TopoDS_Edge& edge,
640 Adaptor3d_Curve& curve,
644 const SMESH_Hypothesis* hyp1d)
646 if ( !hyp1d ) return error( "Invalid LayerDistribution hypothesis");
649 myUsedHyps.push_back( hyp1d );
651 SMESH_Hypothesis::Hypothesis_Status aStatus;
652 if ( !StdMeshers_Regular_1D::CheckHypothesis( mesh, edge, aStatus ))
653 return error( "StdMeshers_Regular_1D::CheckHypothesis() failed "
654 "with LayerDistribution hypothesis");
656 double len = GCPnts_AbscissaPoint::Length( curve, f, l );
658 list< double > params;
659 if ( !StdMeshers_Regular_1D::computeInternalParameters( mesh, curve, len, f, l, params, false ))
660 return error("StdMeshers_Regular_1D failed to compute layers distribution");
662 params.push_front( f );
663 params.push_back ( l );
665 positions.reserve( params.size() );
666 for (list<double>::iterator itU = params.begin(); itU != params.end(); itU++)
667 positions.push_back(( *itU - f ) / ( l - f ));
670 // -----------------------------------------------------------------------------
671 //! Make mesh on an adge using assigned 1d hyp or defaut nb of segments
672 bool ComputeCircularEdge( SMESH_Mesh& aMesh,
673 const StdMeshers_FaceSidePtr& aSide )
676 for ( int i = 0; i < aSide->NbEdges(); ++i )
678 const TopoDS_Edge& edge = aSide->Edge( i );
679 _gen->Compute( aMesh, edge );
680 SMESH_subMesh *sm = aMesh.GetSubMesh( edge );
681 if ( sm->GetComputeState() != SMESH_subMesh::COMPUTE_OK)
683 // find any 1d hyp assigned (there can be a hyp w/o algo)
684 myUsedHyps = SMESH_Algo::GetUsedHypothesis( aMesh, edge, /*ignoreAux=*/true );
685 Hypothesis_Status aStatus;
686 if ( !StdMeshers_Regular_1D::CheckHypothesis( aMesh, edge, aStatus ))
688 // no valid 1d hyp assigned, use default nb of segments
689 _hypType = NB_SEGMENTS;
690 _ivalue[ DISTR_TYPE_IND ] = StdMeshers_NumberOfSegments::DT_Regular;
691 _ivalue[ NB_SEGMENTS_IND ] = _gen->GetDefaultNbSegments();
693 ok &= StdMeshers_Regular_1D::Compute( aMesh, edge );
698 // -----------------------------------------------------------------------------
699 //! Make mesh on an adge using assigned 1d hyp or defaut nb of segments
700 bool EvaluateCircularEdge(SMESH_Mesh& aMesh,
701 const StdMeshers_FaceSidePtr aSide,
702 MapShapeNbElems& aResMap)
705 for ( int i = 0; i < aSide->NbEdges(); ++i )
707 const TopoDS_Edge& anEdge = aSide->Edge( i );
708 _gen->Evaluate( aMesh, anEdge, aResMap );
709 if ( aResMap.count( aMesh.GetSubMesh( anEdge )))
712 // find any 1d hyp assigned
713 myUsedHyps = SMESH_Algo::GetUsedHypothesis(aMesh, anEdge, /*ignoreAux=*/true);
714 Hypothesis_Status aStatus;
715 if ( !StdMeshers_Regular_1D::CheckHypothesis( aMesh, anEdge, aStatus ))
717 // no valid 1d hyp assigned, use default nb of segments
718 _hypType = NB_SEGMENTS;
719 _ivalue[ DISTR_TYPE_IND ] = StdMeshers_NumberOfSegments::DT_Regular;
720 _ivalue[ NB_SEGMENTS_IND ] = _gen->GetDefaultNbSegments();
722 ok &= StdMeshers_Regular_1D::Evaluate( aMesh, anEdge, aResMap );
727 // -----------------------------------------------------------------------------
728 TNodeDistributor( int hypId, int studyId, SMESH_Gen* gen)
729 : StdMeshers_Regular_1D( hypId, studyId, gen)
732 // -----------------------------------------------------------------------------
733 virtual const list <const SMESHDS_Hypothesis *> &
734 GetUsedHypothesis(SMESH_Mesh &, const TopoDS_Shape &, const bool)
738 // -----------------------------------------------------------------------------
742 //=======================================================================
744 * \brief Allow algo to do something after persistent restoration
745 * \param subMesh - restored submesh
747 * call markEdgeAsComputedByMe()
749 //=======================================================================
751 void StdMeshers_RadialQuadrangle_1D2D::SubmeshRestored(SMESH_subMesh* faceSubMesh)
753 if ( !faceSubMesh->IsEmpty() )
755 for ( TopExp_Explorer e( faceSubMesh->GetSubShape(), TopAbs_EDGE ); e.More(); e.Next() )
757 markEdgeAsComputedByMe( TopoDS::Edge( e.Current() ), faceSubMesh );
762 //=======================================================================
765 //=======================================================================
767 bool StdMeshers_RadialQuadrangle_1D2D::Compute(SMESH_Mesh& aMesh,
768 const TopoDS_Shape& aShape)
770 SMESH_MesherHelper helper( aMesh );
771 StdMeshers_Quadrangle_2D::myHelper = & helper;
772 StdMeshers_Quadrangle_2D::myNeedSmooth = false;
773 StdMeshers_Quadrangle_2D::myCheckOri = false;
774 StdMeshers_Quadrangle_2D::myQuadList.clear();
776 StdMeshers_FaceSidePtr circSide, linSide1, linSide2;
777 int nbSides = analyseFace( aShape, &aMesh, circSide, linSide1, linSide2 );
778 if( nbSides > 3 || nbSides < 1 )
779 return error("The face must be a full ellipse or a part of ellipse (i.e. the number "
780 "of edges is less or equal to 3 and one of them is an ellipse curve)");
782 // get not yet computed EDGEs
783 list< TopoDS_Edge > emptyEdges;
784 for ( TopExp_Explorer e( aShape, TopAbs_EDGE ); e.More(); e.Next() )
786 if ( aMesh.GetSubMesh( e.Current() )->IsEmpty() )
787 emptyEdges.push_back( TopoDS::Edge( e.Current() ));
790 TNodeDistributor* algo1d = TNodeDistributor::GetDistributor(aMesh);
792 if ( !algo1d->ComputeCircularEdge( aMesh, circSide ))
793 return error( algo1d->GetComputeError() );
796 TopoDS_Face F = TopoDS::Face(aShape);
797 Handle(Geom_Surface) S = BRep_Tool::Surface(F);
799 myHelper->IsQuadraticSubMesh( aShape );
800 myHelper->SetElementsOnShape( true );
802 vector< double > layerPositions; // [0,1]
804 const SMDS_MeshNode* centerNode = 0;
805 gp_Pnt2d centerUV(0,0);
807 // ------------------------------------------------------------------------------------------
808 if ( nbSides == 1 ) // full ellipse
810 const SMDS_MeshNode* circNode;
811 TopoDS_Edge linEdge = makeEdgeToCenter( circSide, F, circNode );
813 StdMeshers_FaceSidePtr tmpSide =
814 StdMeshers_FaceSide::New( F, linEdge, &aMesh, /*isFrw=*/true, /*skipMedium=*/true );
816 if ( !computeLayerPositions( tmpSide, layerPositions ))
819 UVPtStructVec nodes( layerPositions.size() );
820 nodes[0].node = circNode;
821 for ( size_t i = 0; i < layerPositions.size(); ++i )
823 gp_Pnt2d uv = tmpSide->Value2d( layerPositions[i] );
824 gp_Pnt xyz = S->Value( uv.X(), uv.Y() );
826 nodes[ i ].SetUV( uv.XY() );
827 nodes[ i ].normParam = layerPositions[i];
829 nodes[ i ].node = myHelper->AddNode( xyz.X(), xyz.Y(), xyz.Z(), 0, uv.X(), uv.Y() );
832 linSide1 = StdMeshers_FaceSide::New( nodes );
833 linSide2 = StdMeshers_FaceSide::New( nodes );
835 centerNode = nodes.back().node;
836 centerUV = nodes.back().UV();
838 // ------------------------------------------------------------------------------------------
839 else if ( nbSides == 2 && linSide1->Edge(0).Orientation() == TopAbs_INTERNAL )
841 // full ellipse with an internal radial side
843 // eliminate INTERNAL orientation
844 list< TopoDS_Edge > edges;
845 for ( int iE = 0; iE < linSide1->NbEdges(); ++iE )
847 edges.push_back( linSide1->Edge( iE ));
848 edges.back().Orientation( TopAbs_FORWARD );
851 // orient the internal side
852 bool isVIn0Shared = false;
853 TopoDS_Vertex vIn0 = myHelper->IthVertex( 0, edges.front() );
854 for ( int iE = 0; iE < circSide->NbEdges() && !isVIn0Shared; ++iE )
855 isVIn0Shared = vIn0.IsSame( circSide->FirstVertex( iE ));
857 linSide1 = StdMeshers_FaceSide::New( F, edges, &aMesh,
858 /*isFrw=*/isVIn0Shared, /*skipMedium=*/true );
861 if ( !computeLayerPositions( linSide1, layerPositions, &nbMeshedEdges ))
864 // merge existing nodes with new nodes at layerPositions into a UVPtStructVec
866 if ( nbMeshedEdges != linSide1->NbEdges() )
867 makeMissingMesh( linSide1, nodes, layerPositions, myHelper );
869 nodes = linSide1->GetUVPtStruct();
871 linSide1 = StdMeshers_FaceSide::New( nodes );
872 linSide2 = StdMeshers_FaceSide::New( nodes );
874 centerNode = nodes.back().node;
875 centerUV = nodes.back().UV();
877 // ------------------------------------------------------------------------------------------
878 else if ( nbSides == 2 )
880 // find positions of layers for the first half of linSide1
882 if ( !computeLayerPositions( linSide1, layerPositions, &nbMeshedEdges, /*useHalf=*/true ))
885 // make positions for the whole linSide1
886 for ( size_t i = 0; i < layerPositions.size(); ++i )
888 layerPositions[i] *= 0.5;
890 layerPositions.reserve( layerPositions.size() * 2 );
891 for ( int nb = layerPositions.size()-1; nb > 0; --nb )
892 layerPositions.push_back( layerPositions.back() + layerPositions[nb] - layerPositions[nb-1] );
894 // merge existing nodes with new nodes at layerPositions into a UVPtStructVec
896 if ( nbMeshedEdges != linSide1->NbEdges() )
897 makeMissingMesh( linSide1, nodes, layerPositions, myHelper );
899 nodes = linSide1->GetUVPtStruct();
901 // find a central node
903 while ( nodes[ i ].normParam < 0.5 ) ++i;
904 if (( nodes[ i ].normParam - 0.5 ) > ( 0.5 - nodes[ i-1 ].normParam )) --i;
906 // distribute nodes between two linear sides
907 UVPtStructVec nodes2( nodes.rbegin(), nodes.rbegin() + nodes.size() - i );
908 nodes.resize( i + 1 );
910 linSide1 = StdMeshers_FaceSide::New( nodes );
911 linSide2 = StdMeshers_FaceSide::New( nodes2 );
913 centerNode = nodes.back().node;
914 centerUV = nodes.back().UV();
916 // ------------------------------------------------------------------------------------------
919 // one curve must be a part of ellipse and 2 other curves must be segments of line
921 int nbMeshedEdges1, nbMeshedEdges2;
922 vector< double > layerPositions2;
923 bool ok1 = computeLayerPositions( linSide1, layerPositions, &nbMeshedEdges1 );
924 bool ok2 = computeLayerPositions( linSide2, layerPositions2, &nbMeshedEdges2 );
928 bool linSide1Computed = ( nbMeshedEdges1 == linSide1->NbEdges() );
929 bool linSide2Computed = ( nbMeshedEdges2 == linSide2->NbEdges() );
933 if ( linSide1Computed && !linSide2Computed )
935 // use layer positions of linSide1 to mesh linSide2
936 makeMissingMesh( linSide2, nodes, layerPositions, myHelper );
937 linSide2 = StdMeshers_FaceSide::New( nodes );
939 else if ( linSide2Computed && !linSide1Computed )
941 // use layer positions of linSide2 to mesh linSide1
942 makeMissingMesh( linSide1, nodes, layerPositions2, myHelper );
943 linSide1 = StdMeshers_FaceSide::New( nodes );
945 else if ( !linSide2Computed && !linSide1Computed )
947 // use layer positions of a longer side to mesh the shorter side
948 vector< double >& lp =
949 ( linSide1->Length() > linSide2->Length() ) ? layerPositions : layerPositions2;
951 makeMissingMesh( linSide1, nodes, lp, myHelper );
952 linSide1 = StdMeshers_FaceSide::New( nodes );
953 makeMissingMesh( linSide2, nodes, lp, myHelper );
954 linSide2 = StdMeshers_FaceSide::New( nodes );
957 const UVPtStructVec& nodes2 = linSide2->GetUVPtStruct();
958 centerNode = nodes2.back().node;
959 centerUV = nodes2.back().UV();
962 list< TopoDS_Edge >::iterator ee = emptyEdges.begin();
963 for ( ; ee != emptyEdges.end(); ++ee )
964 markEdgeAsComputedByMe( *ee, aMesh.GetSubMesh( F ));
966 circSide->GetUVPtStruct(); // let sides take into account just computed nodes
967 linSide1->GetUVPtStruct();
968 linSide2->GetUVPtStruct();
970 FaceQuadStruct::Ptr quad( new FaceQuadStruct );
971 quad->side.resize( 4 );
972 quad->side[0] = circSide;
973 quad->side[1] = linSide1;
974 quad->side[2] = StdMeshers_FaceSide::New( circSide.get(), centerNode, ¢erUV );
975 quad->side[3] = linSide2;
977 myQuadList.push_back( quad );
979 // create quadrangles
981 if ( linSide1->NbPoints() == linSide2->NbPoints() )
982 ok = StdMeshers_Quadrangle_2D::computeQuadDominant( aMesh, F, quad );
984 ok = StdMeshers_Quadrangle_2D::computeTriangles( aMesh, F, quad );
986 StdMeshers_Quadrangle_2D::myHelper = 0;
991 //================================================================================
993 * \brief Compute nodes on the radial edge
994 * \retval int - nb of segments on the linSide
996 //================================================================================
998 int StdMeshers_RadialQuadrangle_1D2D::computeLayerPositions(StdMeshers_FaceSidePtr linSide,
999 vector< double > & positions,
1003 // First, try to compute positions of layers
1007 SMESH_Mesh * mesh = myHelper->GetMesh();
1009 const SMESH_Hypothesis* hyp1D = myDistributionHypo ? myDistributionHypo->GetLayerDistribution() : 0;
1010 int nbLayers = myNbLayerHypo ? myNbLayerHypo->GetNumberOfLayers() : 0;
1012 if ( !hyp1D && !nbLayers )
1014 // No own algo hypotheses assigned, so first try to find any 1D hypothesis.
1015 // find a hyp usable by TNodeDistributor
1016 TopoDS_Shape edge = linSide->Edge(0);
1017 const SMESH_HypoFilter* hypKind =
1018 TNodeDistributor::GetDistributor(*mesh)->GetCompatibleHypoFilter(/*ignoreAux=*/true);
1019 hyp1D = mesh->GetHypothesis( edge, *hypKind, /*fromAncestors=*/true);
1021 if ( hyp1D ) // try to compute with hyp1D
1023 BRepAdaptor_CompCurve* curve = linSide->GetCurve3d();
1024 SMESHUtils::Deleter< BRepAdaptor_CompCurve > delCurve( curve );
1025 double f = curve->FirstParameter();
1026 double l = curve->LastParameter();
1029 l = 0.5 * ( f + l ); // first part of linSide is used
1031 if ( !TNodeDistributor::GetDistributor(*mesh)->Compute( positions, linSide->Edge(0),
1032 *curve, f, l, *mesh, hyp1D ))
1034 if ( myDistributionHypo ) { // bad hyp assigned
1035 return error( TNodeDistributor::GetDistributor(*mesh)->GetComputeError() );
1038 // bad hyp found, its Ok, lets try with default nb of segments
1043 if ( positions.empty() ) // try to use nb of layers
1046 nbLayers = _gen->GetDefaultNbSegments();
1050 positions.resize( nbLayers + 1 );
1051 for ( int z = 0; z < nbLayers; ++z )
1052 positions[ z ] = double( z )/ double( nbLayers );
1053 positions.back() = 1;
1057 // Second, check presence of a mesh built by other algo on linSide
1059 int nbEdgesComputed = 0;
1060 for ( int i = 0; i < linSide->NbEdges(); ++i )
1062 nbEdgesComputed += ( !mesh->GetSubMesh( linSide->Edge(i))->IsEmpty() );
1065 if ( nbEdgesComputed == linSide->NbEdges() )
1067 const UVPtStructVec& points = linSide->GetUVPtStruct();
1068 if ( points.size() >= 2 )
1070 positions.resize( points.size() );
1071 for ( size_t i = 0; i < points.size(); ++i )
1072 positions[ i ] = points[i].normParam;
1076 if ( nbMeshedEdges ) *nbMeshedEdges = nbEdgesComputed;
1078 return positions.size();
1082 //=======================================================================
1083 //function : Evaluate
1085 //=======================================================================
1087 bool StdMeshers_RadialQuadrangle_1D2D::Evaluate(SMESH_Mesh& aMesh,
1088 const TopoDS_Shape& aShape,
1089 MapShapeNbElems& aResMap)
1091 if( aShape.ShapeType() != TopAbs_FACE ) {
1094 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
1095 if( aResMap.count(sm) )
1098 vector<int>& aResVec =
1099 aResMap.insert( make_pair(sm, vector<int>(SMDSEntity_Last,0))).first->second;
1101 myHelper = new SMESH_MesherHelper( aMesh );
1102 myHelper->SetSubShape( aShape );
1103 SMESHUtils::Deleter<SMESH_MesherHelper> helperDeleter( myHelper );
1105 TNodeDistributor* algo1d = TNodeDistributor::GetDistributor(aMesh);
1107 StdMeshers_FaceSidePtr circSide, linSide1, linSide2;
1108 int nbSides = analyseFace( aShape, &aMesh, circSide, linSide1, linSide2 );
1109 if( nbSides > 3 || nbSides < 1 )
1112 if ( algo1d->EvaluateCircularEdge( aMesh, circSide, aResMap ))
1115 vector< double > layerPositions; // [0,1]
1117 // ------------------------------------------------------------------------------------------
1120 const TopoDS_Face& F = TopoDS::Face( aShape );
1122 const SMDS_MeshNode* circNode;
1123 TopoDS_Edge linEdge = makeEdgeToCenter( circSide, F, circNode );
1125 StdMeshers_FaceSidePtr tmpSide =
1126 StdMeshers_FaceSide::New( F, linEdge, &aMesh, /*isFrw=*/true, /*skipMedium=*/true );
1128 if ( !computeLayerPositions( tmpSide, layerPositions ))
1131 // ------------------------------------------------------------------------------------------
1132 else if ( nbSides == 2 && linSide1->Edge(0).Orientation() == TopAbs_INTERNAL )
1134 if ( !computeLayerPositions( linSide1, layerPositions ))
1137 // ------------------------------------------------------------------------------------------
1138 else if ( nbSides == 2 )
1140 // find positions of layers for the first half of linSide1
1141 if ( !computeLayerPositions( linSide1, layerPositions, 0, /*useHalf=*/true ))
1144 // ------------------------------------------------------------------------------------------
1145 else // nbSides == 3
1147 if ( !computeLayerPositions(( linSide1->Length() > linSide2->Length() ) ? linSide1 : linSide2,
1152 bool isQuadratic = false;
1153 for ( TopExp_Explorer edge( aShape, TopAbs_EDGE ); edge.More() && !isQuadratic ; edge.Next() )
1155 sm = aMesh.GetSubMesh( edge.Current() );
1156 vector<int>& nbElems = aResMap[ sm ];
1157 if ( SMDSEntity_Quad_Edge < (int) nbElems.size() )
1158 isQuadratic = nbElems[ SMDSEntity_Quad_Edge ];
1161 int nbCircSegments = 0;
1162 for ( int iE = 0; iE < circSide->NbEdges(); ++iE )
1164 sm = aMesh.GetSubMesh( circSide->Edge( iE ));
1165 vector<int>& nbElems = aResMap[ sm ];
1166 if ( SMDSEntity_Quad_Edge < (int) nbElems.size() )
1167 nbCircSegments += ( nbElems[ SMDSEntity_Edge ] + nbElems[ SMDSEntity_Quad_Edge ]);
1170 int nbQuads = nbCircSegments * ( layerPositions.size() - 1 );
1171 int nbTria = nbCircSegments;
1172 int nbNodes = ( nbCircSegments - 1 ) * ( layerPositions.size() - 2 );
1175 nbNodes += (( nbCircSegments - 1 ) * ( layerPositions.size() - 1 ) + // radial
1176 ( nbCircSegments ) * ( layerPositions.size() - 2 )); // circular
1177 aResVec[SMDSEntity_Quad_Triangle ] = nbTria;
1178 aResVec[SMDSEntity_Quad_Quadrangle] = nbQuads;
1182 aResVec[SMDSEntity_Triangle ] = nbTria;
1183 aResVec[SMDSEntity_Quadrangle] = nbQuads;
1185 aResVec[SMDSEntity_Node] = nbNodes;
1189 // evaluation for linSides
1190 vector<int> aResVec(SMDSEntity_Last, 0);
1191 if ( isQuadratic ) {
1192 aResVec[SMDSEntity_Node ] = 2 * ( layerPositions.size() - 1 ) + 1;
1193 aResVec[SMDSEntity_Quad_Edge] = layerPositions.size() - 1;
1196 aResVec[SMDSEntity_Node] = layerPositions.size() - 2;
1197 aResVec[SMDSEntity_Edge] = layerPositions.size() - 1;
1199 sm = aMesh.GetSubMesh( linSide1->Edge(0) );
1200 aResMap[sm] = aResVec;
1203 sm = aMesh.GetSubMesh( linSide2->Edge(0) );
1204 aResMap[sm] = aResVec;
1211 //================================================================================
1213 * \brief Return true if the algorithm can compute mesh on this shape
1215 //================================================================================
1217 bool StdMeshers_RadialQuadrangle_1D2D::IsApplicable( const TopoDS_Shape & aShape, bool toCheckAll )
1219 int nbFoundFaces = 0;
1220 for (TopExp_Explorer exp( aShape, TopAbs_FACE ); exp.More(); exp.Next(), ++nbFoundFaces )
1222 StdMeshers_FaceSidePtr circSide, linSide1, linSide2;
1223 int nbSides = analyseFace( exp.Current(), NULL, circSide, linSide1, linSide2 );
1224 bool ok = ( 0 < nbSides && nbSides <= 3 &&
1225 isCornerInsideCircle( circSide, linSide1, linSide2 ));
1226 if( toCheckAll && !ok ) return false;
1227 if( !toCheckAll && ok ) return true;
1229 if( toCheckAll && nbFoundFaces != 0 ) return true;