1 // Copyright (C) 2007-2016 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License, or (at your option) any later version.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
23 // SMESH SMESH : implementaion of SMESH idl descriptions
24 // File : SMESH_Algo.cxx
25 // Author : Paul RASCLE, EDF
28 #include "SMESH_Algo.hxx"
30 #include "SMDS_EdgePosition.hxx"
31 #include "SMDS_FacePosition.hxx"
32 #include "SMDS_MeshElement.hxx"
33 #include "SMDS_MeshNode.hxx"
34 #include "SMDS_VolumeTool.hxx"
35 #include "SMESHDS_Mesh.hxx"
36 #include "SMESHDS_SubMesh.hxx"
37 #include "SMESH_Comment.hxx"
38 #include "SMESH_Gen.hxx"
39 #include "SMESH_HypoFilter.hxx"
40 #include "SMESH_Mesh.hxx"
41 #include "SMESH_MeshAlgos.hxx"
42 #include "SMESH_TypeDefs.hxx"
43 #include "SMESH_subMesh.hxx"
45 #include <Basics_OCCTVersion.hxx>
47 #include <BRepAdaptor_Curve.hxx>
48 #include <BRepLProp.hxx>
49 #include <BRep_Tool.hxx>
50 #include <GCPnts_AbscissaPoint.hxx>
51 #include <GeomAdaptor_Curve.hxx>
52 #include <Geom_Surface.hxx>
53 #include <LDOMParser.hxx>
55 #include <TopExp_Explorer.hxx>
56 #include <TopLoc_Location.hxx>
57 #include <TopTools_ListIteratorOfListOfShape.hxx>
58 #include <TopTools_ListOfShape.hxx>
60 #include <TopoDS_Edge.hxx>
61 #include <TopoDS_Face.hxx>
62 #include <TopoDS_Vertex.hxx>
63 #include <TopoDS_Wire.hxx>
65 #include <gp_Pnt2d.hxx>
68 #include <Standard_ErrorHandler.hxx>
69 #include <Standard_Failure.hxx>
71 #include "utilities.h"
75 #include "SMESH_ProxyMesh.hxx"
76 #include "SMESH_MesherHelper.hxx"
80 //================================================================================
82 * \brief Returns \a true if two algorithms (described by \a this and the given
83 * algo data) are compatible by their output and input types of elements.
85 //================================================================================
87 bool SMESH_Algo::Features::IsCompatible( const SMESH_Algo::Features& algo2 ) const
89 if ( _dim > algo2._dim ) return algo2.IsCompatible( *this );
90 // algo2 is of highter dimension
91 if ( _outElemTypes.empty() || algo2._inElemTypes.empty() )
93 bool compatible = true;
94 set<SMDSAbs_GeometryType>::const_iterator myOutType = _outElemTypes.begin();
95 for ( ; myOutType != _outElemTypes.end() && compatible; ++myOutType )
96 compatible = algo2._inElemTypes.count( *myOutType );
100 //================================================================================
102 * \brief Return Data of the algorithm
104 //================================================================================
106 const SMESH_Algo::Features& SMESH_Algo::GetFeatures( const std::string& algoType )
108 static map< string, SMESH_Algo::Features > theFeaturesByName;
109 if ( theFeaturesByName.empty() )
111 // Read Plugin.xml files
112 vector< string > xmlPaths = SMESH_Gen::GetPluginXMLPaths();
113 LDOMParser xmlParser;
114 for ( size_t iXML = 0; iXML < xmlPaths.size(); ++iXML )
116 bool error = xmlParser.parse( xmlPaths[iXML].c_str() );
119 TCollection_AsciiString data;
120 INFOS( xmlParser.GetError(data) );
123 // <algorithm type="Regular_1D"
126 // output="QUAD,TRIA">
128 LDOM_Document xmlDoc = xmlParser.getDocument();
129 LDOM_NodeList algoNodeList = xmlDoc.getElementsByTagName( "algorithm" );
130 for ( int i = 0; i < algoNodeList.getLength(); ++i )
132 LDOM_Node algoNode = algoNodeList.item( i );
133 LDOM_Element& algoElem = (LDOM_Element&) algoNode;
134 TCollection_AsciiString algoType = algoElem.getAttribute("type");
135 TCollection_AsciiString input = algoElem.getAttribute("input");
136 TCollection_AsciiString output = algoElem.getAttribute("output");
137 TCollection_AsciiString dim = algoElem.getAttribute("dim");
138 TCollection_AsciiString label = algoElem.getAttribute("label-id");
139 if ( algoType.IsEmpty() ) continue;
141 Features & data = theFeaturesByName[ algoType.ToCString() ];
142 data._dim = dim.IntegerValue();
143 data._label = label.ToCString();
144 for ( int isInput = 0; isInput < 2; ++isInput )
146 TCollection_AsciiString& typeStr = isInput ? input : output;
147 set<SMDSAbs_GeometryType>& typeSet = isInput ? data._inElemTypes : data._outElemTypes;
149 while ( beg <= typeStr.Length() )
151 while ( beg < typeStr.Length() && !isalpha( typeStr.Value( beg ) ))
154 while ( end < typeStr.Length() && isalpha( typeStr.Value( end + 1 ) ))
158 TCollection_AsciiString typeName = typeStr.SubString( beg, end );
159 if ( typeName == "EDGE" ) typeSet.insert( SMDSGeom_EDGE );
160 else if ( typeName == "TRIA" ) typeSet.insert( SMDSGeom_TRIANGLE );
161 else if ( typeName == "QUAD" ) typeSet.insert( SMDSGeom_QUADRANGLE );
169 return theFeaturesByName[ algoType ];
172 //=============================================================================
176 //=============================================================================
178 SMESH_Algo::SMESH_Algo (int hypId, int studyId, SMESH_Gen * gen)
179 : SMESH_Hypothesis(hypId, studyId, gen)
181 _compatibleAllHypFilter = _compatibleNoAuxHypFilter = NULL;
182 _onlyUnaryInput = _requireDiscreteBoundary = _requireShape = true;
183 _quadraticMesh = _supportSubmeshes = false;
185 for ( int i = 0; i < 4; ++i )
186 _neededLowerHyps[ i ] = false;
189 //=============================================================================
193 //=============================================================================
195 SMESH_Algo::~SMESH_Algo()
197 delete _compatibleNoAuxHypFilter;
198 // delete _compatibleAllHypFilter; -- _compatibleNoAuxHypFilter does it!!!
201 //=============================================================================
205 //=============================================================================
207 SMESH_0D_Algo::SMESH_0D_Algo(int hypId, int studyId, SMESH_Gen* gen)
208 : SMESH_Algo(hypId, studyId, gen)
210 _shapeType = (1 << TopAbs_VERTEX);
213 SMESH_1D_Algo::SMESH_1D_Algo(int hypId, int studyId, SMESH_Gen* gen)
214 : SMESH_Algo(hypId, studyId, gen)
216 _shapeType = (1 << TopAbs_EDGE);
219 SMESH_2D_Algo::SMESH_2D_Algo(int hypId, int studyId, SMESH_Gen* gen)
220 : SMESH_Algo(hypId, studyId, gen)
222 _shapeType = (1 << TopAbs_FACE);
225 SMESH_3D_Algo::SMESH_3D_Algo(int hypId, int studyId, SMESH_Gen* gen)
226 : SMESH_Algo(hypId, studyId, gen)
228 _shapeType = (1 << TopAbs_SOLID);
232 //=============================================================================
234 * Usually an algoritm has nothing to save
236 //=============================================================================
238 ostream & SMESH_Algo::SaveTo(ostream & save) { return save; }
239 istream & SMESH_Algo::LoadFrom(istream & load) { return load; }
241 //=============================================================================
245 //=============================================================================
247 const vector < string > &SMESH_Algo::GetCompatibleHypothesis()
249 return _compatibleHypothesis;
252 //=============================================================================
254 * List the hypothesis used by the algorithm associated to the shape.
255 * Hypothesis associated to father shape -are- taken into account (see
256 * GetAppliedHypothesis). Relevant hypothesis have a name (type) listed in
257 * the algorithm. This method could be surcharged by specific algorithms, in
258 * case of several hypothesis simultaneously applicable.
260 //=============================================================================
262 const list <const SMESHDS_Hypothesis *> &
263 SMESH_Algo::GetUsedHypothesis(SMESH_Mesh & aMesh,
264 const TopoDS_Shape & aShape,
265 const bool ignoreAuxiliary) const
267 SMESH_Algo* me = const_cast< SMESH_Algo* >( this );
268 me->_usedHypList.clear();
269 if ( const SMESH_HypoFilter* filter = GetCompatibleHypoFilter( ignoreAuxiliary ))
271 aMesh.GetHypotheses( aShape, *filter, me->_usedHypList, true );
272 if ( ignoreAuxiliary && _usedHypList.size() > 1 )
273 me->_usedHypList.clear(); //only one compatible hypothesis allowed
278 //=============================================================================
280 * List the relevant hypothesis associated to the shape. Relevant hypothesis
281 * have a name (type) listed in the algorithm. Hypothesis associated to
282 * father shape -are not- taken into account (see GetUsedHypothesis)
284 //=============================================================================
286 const list<const SMESHDS_Hypothesis *> &
287 SMESH_Algo::GetAppliedHypothesis(SMESH_Mesh & aMesh,
288 const TopoDS_Shape & aShape,
289 const bool ignoreAuxiliary) const
291 SMESH_Algo* me = const_cast< SMESH_Algo* >( this );
292 me->_appliedHypList.clear();
293 if ( const SMESH_HypoFilter* filter = GetCompatibleHypoFilter( ignoreAuxiliary ))
294 aMesh.GetHypotheses( aShape, *filter, me->_appliedHypList, false );
296 return _appliedHypList;
299 //=============================================================================
301 * Compute length of an edge
303 //=============================================================================
305 double SMESH_Algo::EdgeLength(const TopoDS_Edge & E)
307 double UMin = 0, UMax = 0;
309 Handle(Geom_Curve) C = BRep_Tool::Curve(E, L, UMin, UMax);
312 GeomAdaptor_Curve AdaptCurve(C, UMin, UMax); //range is important for periodic curves
313 double length = GCPnts_AbscissaPoint::Length(AdaptCurve, UMin, UMax);
317 //================================================================================
319 * \brief Just return false as the algorithm does not hold parameters values
321 //================================================================================
323 bool SMESH_Algo::SetParametersByMesh(const SMESH_Mesh* /*theMesh*/,
324 const TopoDS_Shape& /*theShape*/)
328 bool SMESH_Algo::SetParametersByDefaults(const TDefaults& , const SMESH_Mesh*)
332 //================================================================================
334 * \brief Fill vector of node parameters on geometrical edge, including vertex nodes
335 * \param theMesh - The mesh containing nodes
336 * \param theEdge - The geometrical edge of interest
337 * \param theParams - The resulting vector of sorted node parameters
338 * \retval bool - false if not all parameters are OK
340 //================================================================================
342 bool SMESH_Algo::GetNodeParamOnEdge(const SMESHDS_Mesh* theMesh,
343 const TopoDS_Edge& theEdge,
344 vector< double > & theParams)
348 if ( !theMesh || theEdge.IsNull() )
351 SMESHDS_SubMesh * eSubMesh = theMesh->MeshElements( theEdge );
352 if ( !eSubMesh || !eSubMesh->GetElements()->more() )
353 return false; // edge is not meshed
355 //int nbEdgeNodes = 0;
356 set < double > paramSet;
359 // loop on nodes of an edge: sort them by param on edge
360 SMDS_NodeIteratorPtr nIt = eSubMesh->GetNodes();
361 while ( nIt->more() )
363 const SMDS_MeshNode* node = nIt->next();
364 const SMDS_PositionPtr& pos = node->GetPosition();
365 if ( pos->GetTypeOfPosition() != SMDS_TOP_EDGE )
367 const SMDS_EdgePosition* epos =
368 static_cast<const SMDS_EdgePosition*>(node->GetPosition());
369 if ( !paramSet.insert( epos->GetUParameter() ).second )
370 return false; // equal parameters
373 // add vertex nodes params
375 TopExp::Vertices( theEdge, V1, V2);
376 if ( VertexNode( V1, theMesh ) &&
377 !paramSet.insert( BRep_Tool::Parameter(V1,theEdge) ).second )
378 return false; // there are equal parameters
379 if ( VertexNode( V2, theMesh ) &&
380 !paramSet.insert( BRep_Tool::Parameter(V2,theEdge) ).second )
381 return false; // there are equal parameters
384 theParams.resize( paramSet.size() );
385 set < double >::iterator par = paramSet.begin();
386 vector< double >::iterator vecPar = theParams.begin();
387 for ( ; par != paramSet.end(); ++par, ++vecPar )
390 return theParams.size() > 1;
393 //================================================================================
395 * \brief Fill vector of node parameters on geometrical edge, including vertex nodes
396 * \param theMesh - The mesh containing nodes
397 * \param theEdge - The geometrical edge of interest
398 * \param theParams - The resulting vector of sorted node parameters
399 * \retval bool - false if not all parameters are OK
401 //================================================================================
403 bool SMESH_Algo::GetSortedNodesOnEdge(const SMESHDS_Mesh* theMesh,
404 const TopoDS_Edge& theEdge,
405 const bool ignoreMediumNodes,
406 map< double, const SMDS_MeshNode* > & theNodes,
407 const SMDSAbs_ElementType typeToCheck)
411 if ( !theMesh || theEdge.IsNull() )
414 SMESHDS_SubMesh * eSubMesh = theMesh->MeshElements( theEdge );
415 if ( !eSubMesh || ( eSubMesh->NbElements() == 0 && eSubMesh->NbNodes() == 0))
416 return false; // edge is not meshed
419 set < double > paramSet;
422 // loop on nodes of an edge: sort them by param on edge
423 SMDS_NodeIteratorPtr nIt = eSubMesh->GetNodes();
424 while ( nIt->more() )
426 const SMDS_MeshNode* node = nIt->next();
427 if ( ignoreMediumNodes && SMESH_MesherHelper::IsMedium( node, typeToCheck ))
429 const SMDS_PositionPtr& pos = node->GetPosition();
430 if ( pos->GetTypeOfPosition() != SMDS_TOP_EDGE )
432 const SMDS_EdgePosition* epos =
433 static_cast<const SMDS_EdgePosition*>(node->GetPosition());
434 theNodes.insert( theNodes.end(), make_pair( epos->GetUParameter(), node ));
439 TopoDS_Vertex v1, v2;
440 TopExp::Vertices(theEdge, v1, v2);
441 const SMDS_MeshNode* n1 = VertexNode( v1, eSubMesh, 0 );
442 const SMDS_MeshNode* n2 = VertexNode( v2, eSubMesh, 0 );
444 BRep_Tool::Range(theEdge, f, l);
445 if ( v1.Orientation() != TopAbs_FORWARD )
447 if ( n1 && ++nbNodes )
448 theNodes.insert( make_pair( f, n1 ));
449 if ( n2 && ++nbNodes )
450 theNodes.insert( make_pair( l, n2 ));
452 return (int)theNodes.size() == nbNodes;
455 //================================================================================
457 * \brief Returns the filter recognizing only compatible hypotheses
458 * \param ignoreAuxiliary - make filter ignore auxiliary hypotheses
459 * \retval SMESH_HypoFilter* - the filter that can be NULL
461 //================================================================================
463 const SMESH_HypoFilter*
464 SMESH_Algo::GetCompatibleHypoFilter(const bool ignoreAuxiliary) const
466 if ( !_compatibleHypothesis.empty() )
468 if ( !_compatibleAllHypFilter )
470 SMESH_HypoFilter* filter = new SMESH_HypoFilter();
471 filter->Init( filter->HasName( _compatibleHypothesis[0] ));
472 for ( size_t i = 1; i < _compatibleHypothesis.size(); ++i )
473 filter->Or( filter->HasName( _compatibleHypothesis[ i ] ));
475 SMESH_HypoFilter* filterNoAux = new SMESH_HypoFilter( filter );
476 filterNoAux->AndNot( filterNoAux->IsAuxiliary() );
478 // _compatibleNoAuxHypFilter will detele _compatibleAllHypFilter!!!
479 SMESH_Algo* me = const_cast< SMESH_Algo* >( this );
480 me->_compatibleAllHypFilter = filter;
481 me->_compatibleNoAuxHypFilter = filterNoAux;
483 return ignoreAuxiliary ? _compatibleNoAuxHypFilter : _compatibleAllHypFilter;
488 //================================================================================
490 * \brief Return continuity of two edges
491 * \param E1 - the 1st edge
492 * \param E2 - the 2nd edge
493 * \retval GeomAbs_Shape - regularity at the junction between E1 and E2
495 //================================================================================
497 GeomAbs_Shape SMESH_Algo::Continuity(const TopoDS_Edge& theE1,
498 const TopoDS_Edge& theE2)
500 // avoid pb with internal edges
501 TopoDS_Edge E1 = theE1, E2 = theE2;
502 if (E1.Orientation() > TopAbs_REVERSED) // INTERNAL
503 E1.Orientation( TopAbs_FORWARD );
504 if (E2.Orientation() > TopAbs_REVERSED) // INTERNAL
505 E2.Orientation( TopAbs_FORWARD );
507 TopoDS_Vertex V, VV1[2], VV2[2];
508 TopExp::Vertices( E1, VV1[0], VV1[1], true );
509 TopExp::Vertices( E2, VV2[0], VV2[1], true );
510 if ( VV1[1].IsSame( VV2[0] )) { V = VV1[1]; }
511 else if ( VV1[0].IsSame( VV2[1] )) { V = VV1[0]; }
512 else if ( VV1[1].IsSame( VV2[1] )) { V = VV1[1]; E1.Reverse(); }
513 else if ( VV1[0].IsSame( VV2[0] )) { V = VV1[0]; E1.Reverse(); }
514 else { return GeomAbs_C0; }
516 Standard_Real u1 = BRep_Tool::Parameter( V, E1 );
517 Standard_Real u2 = BRep_Tool::Parameter( V, E2 );
518 BRepAdaptor_Curve C1( E1 ), C2( E2 );
519 Standard_Real tol = BRep_Tool::Tolerance( V );
520 Standard_Real angTol = 2e-3;
523 return BRepLProp::Continuity(C1, C2, u1, u2, tol, angTol);
525 catch (Standard_Failure) {
530 //================================================================================
532 * \brief Return true if an edge can be considered straight
534 //================================================================================
536 bool SMESH_Algo::IsStraight( const TopoDS_Edge & E,
537 const bool degenResult)
541 if ( BRep_Tool::Curve( E, f, l ).IsNull())
544 BRepAdaptor_Curve curve( E );
545 switch( curve.GetType() )
550 case GeomAbs_Ellipse:
551 case GeomAbs_Hyperbola:
552 case GeomAbs_Parabola:
554 // case GeomAbs_BezierCurve:
555 // case GeomAbs_BSplineCurve:
556 // case GeomAbs_OtherCurve:
559 const double f = curve.FirstParameter();
560 const double l = curve.LastParameter();
561 const gp_Pnt pf = curve.Value( f );
562 const gp_Pnt pl = curve.Value( l );
563 const gp_Vec v1( pf, pl );
564 const double v1Len = v1.Magnitude();
565 if ( v1Len < std::numeric_limits< double >::min() )
566 return false; // E seems closed
567 const double tol = Min( 10 * curve.Tolerance(), v1Len * 1e-2 );
568 const double nbSamples = 7;
569 for ( int i = 0; i < nbSamples; ++i )
571 const double r = ( i + 1 ) / nbSamples;
572 const gp_Pnt pi = curve.Value( f * r + l * ( 1 - r ));
573 const gp_Vec vi( pf, pi );
574 const double h = 0.5 * v1.Crossed( vi ).Magnitude() / v1Len;
581 //================================================================================
583 * \brief Return true if an edge has no 3D curve
585 //================================================================================
587 bool SMESH_Algo::isDegenerated( const TopoDS_Edge & E, const bool checkLength )
590 return EdgeLength( E ) == 0;
593 Handle(Geom_Curve) C = BRep_Tool::Curve( E, loc, f,l );
597 //================================================================================
599 * \brief Return the node built on a vertex
600 * \param V - the vertex
601 * \param meshDS - mesh
602 * \retval const SMDS_MeshNode* - found node or NULL
603 * \sa SMESH_MesherHelper::GetSubShapeByNode( const SMDS_MeshNode*, SMESHDS_Mesh* )
605 //================================================================================
607 const SMDS_MeshNode* SMESH_Algo::VertexNode(const TopoDS_Vertex& V,
608 const SMESHDS_Mesh* meshDS)
610 if ( SMESHDS_SubMesh* sm = meshDS->MeshElements(V) ) {
611 SMDS_NodeIteratorPtr nIt= sm->GetNodes();
618 //=======================================================================
620 * \brief Return the node built on a vertex.
621 * A node moved to other geometry by MergeNodes() is also returned.
622 * \param V - the vertex
624 * \retval const SMDS_MeshNode* - found node or NULL
626 //=======================================================================
628 const SMDS_MeshNode* SMESH_Algo::VertexNode(const TopoDS_Vertex& V,
629 const SMESH_Mesh* mesh)
631 const SMDS_MeshNode* node = VertexNode( V, mesh->GetMeshDS() );
633 if ( !node && mesh->HasModificationsToDiscard() )
635 PShapeIteratorPtr edgeIt = SMESH_MesherHelper::GetAncestors( V, *mesh, TopAbs_EDGE );
636 while ( const TopoDS_Shape* edge = edgeIt->next() )
637 if ( SMESHDS_SubMesh* edgeSM = mesh->GetMeshDS()->MeshElements( *edge ))
638 if ( edgeSM->NbElements() > 0 )
639 return VertexNode( V, edgeSM, mesh, /*checkV=*/false );
644 //=======================================================================
646 * \brief Return the node built on a vertex.
647 * A node moved to other geometry by MergeNodes() is also returned.
648 * \param V - the vertex
649 * \param edgeSM - sub-mesh of a meshed EDGE sharing the vertex
650 * \param checkV - if \c true, presence of a node on the vertex is checked
651 * \retval const SMDS_MeshNode* - found node or NULL
653 //=======================================================================
655 const SMDS_MeshNode* SMESH_Algo::VertexNode(const TopoDS_Vertex& V,
656 const SMESHDS_SubMesh* edgeSM,
657 const SMESH_Mesh* mesh,
660 const SMDS_MeshNode* node = checkV ? VertexNode( V, edgeSM->GetParent() ) : 0;
662 if ( !node && edgeSM )
664 // find nodes not shared by mesh segments
665 typedef set< const SMDS_MeshNode* > TNodeSet;
666 typedef map< const SMDS_MeshNode*, const SMDS_MeshNode* > TNodeMap;
667 TNodeMap notSharedNodes;
668 TNodeSet otherShapeNodes;
669 vector< const SMDS_MeshNode* > segNodes(3);
670 SMDS_ElemIteratorPtr segIt = edgeSM->GetElements();
671 while ( segIt->more() )
673 const SMDS_MeshElement* seg = segIt->next();
674 if ( seg->GetType() != SMDSAbs_Edge )
676 segNodes.assign( seg->begin_nodes(), seg->end_nodes() );
677 for ( int i = 0; i < 2; ++i )
679 const SMDS_MeshNode* n1 = segNodes[i];
680 const SMDS_MeshNode* n2 = segNodes[1-i];
681 pair<TNodeMap::iterator, bool> it2new = notSharedNodes.insert( make_pair( n1, n2 ));
682 if ( !it2new.second ) // n encounters twice
683 notSharedNodes.erase( it2new.first );
684 if ( n1->getshapeId() != edgeSM->GetID() )
685 otherShapeNodes.insert( n1 );
688 if ( otherShapeNodes.size() == 1 && notSharedNodes.empty() ) // a closed EDGE
689 return *otherShapeNodes.begin();
691 if ( notSharedNodes.size() == 2 ) // two end nodes found
693 SMESHDS_Mesh* meshDS = edgeSM->GetParent();
694 const TopoDS_Shape& E = meshDS->IndexToShape( edgeSM->GetID() );
695 if ( E.IsNull() || E.ShapeType() != TopAbs_EDGE )
697 const SMDS_MeshNode* n1 = notSharedNodes.begin ()->first;
698 const SMDS_MeshNode* n2 = notSharedNodes.rbegin()->first;
699 TopoDS_Shape S1 = SMESH_MesherHelper::GetSubShapeByNode( n1, meshDS );
700 if ( S1.ShapeType() == TopAbs_VERTEX && SMESH_MesherHelper::IsSubShape( S1, E ))
702 TopoDS_Shape S2 = SMESH_MesherHelper::GetSubShapeByNode( n2, meshDS );
703 if ( S2.ShapeType() == TopAbs_VERTEX && SMESH_MesherHelper::IsSubShape( S2, E ))
705 if ( edgeSM->NbElements() <= 2 || !mesh ) // one-two segments
707 gp_Pnt pV = BRep_Tool::Pnt( V );
708 double dist1 = pV.SquareDistance( SMESH_TNodeXYZ( n1 ));
709 double dist2 = pV.SquareDistance( SMESH_TNodeXYZ( n2 ));
710 return dist1 < dist2 ? n1 : n2;
714 SMESH_MesherHelper helper( const_cast<SMESH_Mesh&>( *mesh ));
715 const SMDS_MeshNode* n1i = notSharedNodes.begin ()->second;
716 const SMDS_MeshNode* n2i = notSharedNodes.rbegin()->second;
717 const TopoDS_Edge& edge = TopoDS::Edge( E );
719 double pos1 = helper.GetNodeU( edge, n1i, n2i, &posOK );
720 double pos2 = helper.GetNodeU( edge, n2i, n1i, &posOK );
721 double posV = BRep_Tool::Parameter( V, edge );
722 if ( Abs( pos1 - posV ) < Abs( pos2 - posV )) return n1;
730 //=======================================================================
731 //function : GetMeshError
732 //purpose : Finds topological errors of a sub-mesh
733 //WARNING : 1D check is NOT implemented so far
734 //=======================================================================
736 SMESH_Algo::EMeshError SMESH_Algo::GetMeshError(SMESH_subMesh* subMesh)
738 EMeshError err = MEr_OK;
740 SMESHDS_SubMesh* smDS = subMesh->GetSubMeshDS();
744 switch ( subMesh->GetSubShape().ShapeType() )
746 case TopAbs_FACE: { // ====================== 2D =====================
748 SMDS_ElemIteratorPtr fIt = smDS->GetElements();
752 // We check that olny links on EDGEs encouter once, the rest links, twice
753 set< SMESH_TLink > links;
754 while ( fIt->more() )
756 const SMDS_MeshElement* f = fIt->next();
757 int nbNodes = f->NbCornerNodes(); // ignore medium nodes
758 for ( int i = 0; i < nbNodes; ++i )
760 const SMDS_MeshNode* n1 = f->GetNode( i );
761 const SMDS_MeshNode* n2 = f->GetNode(( i+1 ) % nbNodes);
762 std::pair< set< SMESH_TLink >::iterator, bool > it_added =
763 links.insert( SMESH_TLink( n1, n2 ));
764 if ( !it_added.second )
765 // As we do NOT(!) check if mesh is manifold, we believe that a link can
766 // encounter once or twice only (not three times), we erase a link as soon
767 // as it encounters twice to speed up search in the <links> map.
768 links.erase( it_added.first );
771 // the links remaining in the <links> should all be on EDGE
772 set< SMESH_TLink >::iterator linkIt = links.begin();
773 for ( ; linkIt != links.end(); ++linkIt )
775 const SMESH_TLink& link = *linkIt;
776 if ( link.node1()->GetPosition()->GetTypeOfPosition() > SMDS_TOP_EDGE ||
777 link.node2()->GetPosition()->GetTypeOfPosition() > SMDS_TOP_EDGE )
780 // TODO: to check orientation
783 case TopAbs_SOLID: { // ====================== 3D =====================
785 SMDS_ElemIteratorPtr vIt = smDS->GetElements();
789 SMDS_VolumeTool vTool;
790 while ( !vIt->more() )
792 if (!vTool.Set( vIt->next() ))
795 for ( int iF = 0; iF < vTool.NbFaces(); ++iF )
796 if ( vTool.IsFreeFace( iF ))
798 int nbN = vTool.NbFaceNodes( iF );
799 const SMDS_MeshNode** nodes = vTool.GetFaceNodes( iF );
800 for ( int i = 0; i < nbN; ++i )
801 if ( nodes[i]->GetPosition()->GetTypeOfPosition() > SMDS_TOP_FACE )
812 //================================================================================
814 * \brief Sets event listener to submeshes if necessary
815 * \param subMesh - submesh where algo is set
817 * After being set, event listener is notified on each event of a submesh.
818 * By default non listener is set
820 //================================================================================
822 void SMESH_Algo::SetEventListener(SMESH_subMesh* /*subMesh*/)
826 //================================================================================
828 * \brief Allow algo to do something after persistent restoration
829 * \param subMesh - restored submesh
831 * This method is called only if a submesh has HYP_OK algo_state.
833 //================================================================================
835 void SMESH_Algo::SubmeshRestored(SMESH_subMesh* /*subMesh*/)
839 //================================================================================
841 * \brief Computes mesh without geometry
842 * \param aMesh - the mesh
843 * \param aHelper - helper that must be used for adding elements to \aaMesh
844 * \retval bool - is a success
846 //================================================================================
848 bool SMESH_Algo::Compute(SMESH_Mesh & /*aMesh*/, SMESH_MesherHelper* /*aHelper*/)
850 return error( COMPERR_BAD_INPUT_MESH, "Mesh built on shape expected");
853 //=======================================================================
854 //function : CancelCompute
855 //purpose : Sets _computeCanceled to true. It's usage depends on
856 // * implementation of a particular mesher.
857 //=======================================================================
859 void SMESH_Algo::CancelCompute()
861 _computeCanceled = true;
862 _error = COMPERR_CANCELED;
865 //================================================================================
867 * If possible, returns progress of computation [0.,1.]
869 //================================================================================
871 double SMESH_Algo::GetProgress() const
876 //================================================================================
878 * \brief store error and comment and then return ( error == COMPERR_OK )
880 //================================================================================
882 bool SMESH_Algo::error(int error, const SMESH_Comment& comment)
886 return ( error == COMPERR_OK );
889 //================================================================================
891 * \brief store error and return ( error == COMPERR_OK )
893 //================================================================================
895 bool SMESH_Algo::error(SMESH_ComputeErrorPtr error)
898 _error = error->myName;
899 _comment = error->myComment;
900 _badInputElements = error->myBadElements;
901 return error->IsOK();
906 //================================================================================
908 * \brief return compute error
910 //================================================================================
912 SMESH_ComputeErrorPtr SMESH_Algo::GetComputeError() const
914 SMESH_ComputeErrorPtr err = SMESH_ComputeError::New( _error, _comment, this );
915 // hope this method is called by only SMESH_subMesh after this->Compute()
916 err->myBadElements.splice( err->myBadElements.end(),
917 (list<const SMDS_MeshElement*>&) _badInputElements );
921 //================================================================================
923 * \brief initialize compute error before call of Compute()
925 //================================================================================
927 void SMESH_Algo::InitComputeError()
931 list<const SMDS_MeshElement*>::iterator elem = _badInputElements.begin();
932 for ( ; elem != _badInputElements.end(); ++elem )
933 if ( (*elem)->GetID() < 1 )
935 _badInputElements.clear();
937 _computeCanceled = false;
942 //================================================================================
944 * \brief Return compute progress by nb of calls of this method
946 //================================================================================
948 double SMESH_Algo::GetProgressByTic() const
951 for ( size_t i = 0; i < _smToCompute.size(); ++i )
952 computeCost += _smToCompute[i]->GetComputeCost();
954 const_cast<SMESH_Algo*>( this )->_progressTic++;
956 double x = 5 * _progressTic;
957 x = ( x < computeCost ) ? ( x / computeCost ) : 1.;
958 return 0.9 * sin( x * M_PI / 2 );
961 //================================================================================
963 * \brief store a bad input element preventing computation,
964 * which may be a temporary one i.e. not residing the mesh,
965 * then it will be deleted by InitComputeError()
967 //================================================================================
969 void SMESH_Algo::addBadInputElement(const SMDS_MeshElement* elem)
972 _badInputElements.push_back( elem );
975 //=======================================================================
976 //function : addBadInputElements
977 //purpose : store a bad input elements or nodes preventing computation
978 //=======================================================================
980 void SMESH_Algo::addBadInputElements(const SMESHDS_SubMesh* sm,
987 SMDS_NodeIteratorPtr nIt = sm->GetNodes();
988 while ( nIt->more() ) addBadInputElement( nIt->next() );
992 SMDS_ElemIteratorPtr eIt = sm->GetElements();
993 while ( eIt->more() ) addBadInputElement( eIt->next() );
998 //=============================================================================
1002 //=============================================================================
1004 // int SMESH_Algo::NumberOfWires(const TopoDS_Shape& S)
1007 // for (TopExp_Explorer exp(S,TopAbs_WIRE); exp.More(); exp.Next())
1012 //=============================================================================
1016 //=============================================================================
1018 int SMESH_Algo::NumberOfPoints(SMESH_Mesh& aMesh, const TopoDS_Wire& W)
1021 for (TopExp_Explorer exp(W,TopAbs_EDGE); exp.More(); exp.Next()) {
1022 const TopoDS_Edge& E = TopoDS::Edge(exp.Current());
1023 int nb = aMesh.GetSubMesh(E)->GetSubMeshDS()->NbNodes();
1026 nbPoints += nb + 1; // internal points plus 1 vertex of 2 (last point ?)
1032 //================================================================================
1034 * Method in which an algorithm generating a structured mesh
1035 * fixes positions of in-face nodes after there movement
1036 * due to insertion of viscous layers.
1038 //================================================================================
1040 bool SMESH_2D_Algo::FixInternalNodes(const SMESH_ProxyMesh& mesh,
1041 const TopoDS_Face& face)
1043 const SMESHDS_SubMesh* smDS = mesh.GetSubMesh(face);
1044 if ( !smDS || smDS->NbElements() < 1 )
1047 SMESH_MesherHelper helper( *mesh.GetMesh() );
1049 // get all faces from a proxy sub-mesh
1050 typedef SMDS_StdIterator< const SMDS_MeshElement*, SMDS_ElemIteratorPtr > TIterator;
1051 TIDSortedElemSet allFaces( TIterator( smDS->GetElements() ), TIterator() );
1052 TIDSortedElemSet avoidSet, firstRowQuads;
1054 // indices of nodes to pass to a neighbour quad using SMESH_MeshAlgos::FindFaceInSet()
1057 // get two first rows of nodes by passing through the first row of faces
1058 vector< vector< const SMDS_MeshNode* > > nodeRows;
1059 int iRow1 = 0, iRow2 = 1;
1060 const SMDS_MeshElement* quad;
1062 // look for a corner quadrangle and it's corner node
1063 const SMDS_MeshElement* cornerQuad = 0;
1064 int cornerNodeInd = -1;
1065 SMDS_ElemIteratorPtr fIt = smDS->GetElements();
1066 while ( !cornerQuad && fIt->more() )
1068 cornerQuad = fIt->next();
1069 if ( cornerQuad->NbCornerNodes() != 4 )
1071 SMDS_NodeIteratorPtr nIt = cornerQuad->nodeIterator();
1072 for ( int i = 0; i < 4; ++i )
1074 int nbInverseQuads = 0;
1075 SMDS_ElemIteratorPtr fIt = nIt->next()->GetInverseElementIterator(SMDSAbs_Face);
1076 while ( fIt->more() )
1077 nbInverseQuads += allFaces.count( fIt->next() );
1078 if ( nbInverseQuads == 1 )
1079 cornerNodeInd = i, i = 4;
1081 if ( cornerNodeInd < 0 )
1084 if ( !cornerQuad || cornerNodeInd < 0 )
1087 iN1 = helper.WrapIndex( cornerNodeInd + 1, 4 );
1088 iN2 = helper.WrapIndex( cornerNodeInd + 2, 4 );
1089 int iN3 = helper.WrapIndex( cornerNodeInd + 3, 4 );
1091 nodeRows[iRow1].push_back( cornerQuad->GetNode( cornerNodeInd ));
1092 nodeRows[iRow1].push_back( cornerQuad->GetNode( iN1 ));
1093 nodeRows[iRow2].push_back( cornerQuad->GetNode( iN3 ));
1094 nodeRows[iRow2].push_back( cornerQuad->GetNode( iN2 ));
1095 firstRowQuads.insert( cornerQuad );
1097 // pass through the rest quads in a face row
1102 avoidSet.insert( quad );
1103 if (( quad = SMESH_MeshAlgos::FindFaceInSet( nodeRows[iRow1].back(),
1104 nodeRows[iRow2].back(),
1105 allFaces, avoidSet, &iN1, &iN2)))
1107 nodeRows[iRow1].push_back( quad->GetNode( helper.WrapIndex( iN2 + 2, 4 )));
1108 nodeRows[iRow2].push_back( quad->GetNode( helper.WrapIndex( iN1 + 2, 4 )));
1109 if ( quad->NbCornerNodes() != 4 )
1113 if ( nodeRows[iRow1].size() < 3 )
1114 return true; // there is nothing to fix
1117 nodeRows.reserve( smDS->NbElements() / nodeRows[iRow1].size() );
1119 // get the rest node rows
1124 // get the first quad in the next face row
1125 if (( quad = SMESH_MeshAlgos::FindFaceInSet( nodeRows[iRow1][0],
1127 allFaces, /*avoid=*/firstRowQuads,
1130 if ( quad->NbCornerNodes() != 4 )
1132 nodeRows.resize( iRow2+1 );
1133 nodeRows[iRow2].push_back( quad->GetNode( helper.WrapIndex( iN2 + 2, 4 )));
1134 nodeRows[iRow2].push_back( quad->GetNode( helper.WrapIndex( iN1 + 2, 4 )));
1135 firstRowQuads.insert( quad );
1139 break; // no more rows
1142 // pass through the rest quads in a face row
1146 avoidSet.insert( quad );
1147 if (( quad = SMESH_MeshAlgos::FindFaceInSet( nodeRows[iRow1][ nodeRows[iRow2].size()-1 ],
1148 nodeRows[iRow2].back(),
1149 allFaces, avoidSet, &iN1, &iN2)))
1151 if ( quad->NbCornerNodes() != 4 )
1153 nodeRows[iRow2].push_back( quad->GetNode( helper.WrapIndex( iN1 + 2, 4 )));
1156 if ( nodeRows[iRow1].size() != nodeRows[iRow2].size() )
1159 if ( nodeRows.size() < 3 )
1160 return true; // there is nothing to fix
1162 // get params of the first (bottom) and last (top) node rows
1163 UVPtStructVec uvB( nodeRows[0].size() ), uvT( nodeRows[0].size() );
1164 for ( int isBot = 0; isBot < 2; ++isBot )
1166 UVPtStructVec & uvps = isBot ? uvB : uvT;
1167 vector< const SMDS_MeshNode* >& nodes = nodeRows[ isBot ? 0 : nodeRows.size()-1 ];
1168 for ( size_t i = 0; i < nodes.size(); ++i )
1170 uvps[i].node = nodes[i];
1171 gp_XY uv = helper.GetNodeUV( face, uvps[i].node );
1172 uvps[i].u = uv.Coord(1);
1173 uvps[i].v = uv.Coord(2);
1176 // calculate x (normalized param)
1177 for ( size_t i = 1; i < nodes.size(); ++i )
1178 uvps[i].x = uvps[i-1].x + SMESH_TNodeXYZ( uvps[i-1].node ).Distance( uvps[i].node );
1179 for ( size_t i = 1; i < nodes.size(); ++i )
1180 uvps[i].x /= uvps.back().x;
1183 // get params of the left and right node rows
1184 UVPtStructVec uvL( nodeRows.size() ), uvR( nodeRows.size() );
1185 for ( int isLeft = 0; isLeft < 2; ++isLeft )
1187 UVPtStructVec & uvps = isLeft ? uvL : uvR;
1188 const int iCol = isLeft ? 0 : nodeRows[0].size() - 1;
1189 for ( size_t i = 0; i < nodeRows.size(); ++i )
1191 uvps[i].node = nodeRows[i][iCol];
1192 gp_XY uv = helper.GetNodeUV( face, uvps[i].node );
1193 uvps[i].u = uv.Coord(1);
1194 uvps[i].v = uv.Coord(2);
1197 // calculate y (normalized param)
1198 for ( size_t i = 1; i < nodeRows.size(); ++i )
1199 uvps[i].y = uvps[i-1].y + SMESH_TNodeXYZ( uvps[i-1].node ).Distance( uvps[i].node );
1200 for ( size_t i = 1; i < nodeRows.size(); ++i )
1201 uvps[i].y /= uvps.back().y;
1204 // update node coordinates
1205 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
1206 Handle(Geom_Surface) S = BRep_Tool::Surface( face );
1207 gp_XY a0 ( uvB.front().u, uvB.front().v );
1208 gp_XY a1 ( uvB.back().u, uvB.back().v );
1209 gp_XY a2 ( uvT.back().u, uvT.back().v );
1210 gp_XY a3 ( uvT.front().u, uvT.front().v );
1211 for ( size_t iRow = 1; iRow < nodeRows.size()-1; ++iRow )
1213 gp_XY p1 ( uvR[ iRow ].u, uvR[ iRow ].v );
1214 gp_XY p3 ( uvL[ iRow ].u, uvL[ iRow ].v );
1215 const double y0 = uvL[ iRow ].y;
1216 const double y1 = uvR[ iRow ].y;
1217 for ( size_t iCol = 1; iCol < nodeRows[0].size()-1; ++iCol )
1219 gp_XY p0 ( uvB[ iCol ].u, uvB[ iCol ].v );
1220 gp_XY p2 ( uvT[ iCol ].u, uvT[ iCol ].v );
1221 const double x0 = uvB[ iCol ].x;
1222 const double x1 = uvT[ iCol ].x;
1223 double x = (x0 + y0 * (x1 - x0)) / (1 - (y1 - y0) * (x1 - x0));
1224 double y = y0 + x * (y1 - y0);
1225 gp_XY uv = helper.calcTFI( x, y, a0,a1,a2,a3, p0,p1,p2,p3 );
1226 gp_Pnt p = S->Value( uv.Coord(1), uv.Coord(2));
1227 const SMDS_MeshNode* n = nodeRows[iRow][iCol];
1228 meshDS->MoveNode( n, p.X(), p.Y(), p.Z() );
1229 if ( SMDS_FacePosition* pos = dynamic_cast< SMDS_FacePosition*>( n->GetPosition() ))
1230 pos->SetParameters( uv.Coord(1), uv.Coord(2) );