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 #include "SMESH_ControlsDef.hxx"
25 #include "SMDS_BallElement.hxx"
26 #include "SMDS_Iterator.hxx"
27 #include "SMDS_Mesh.hxx"
28 #include "SMDS_MeshElement.hxx"
29 #include "SMDS_MeshNode.hxx"
30 #include "SMDS_QuadraticEdge.hxx"
31 #include "SMDS_QuadraticFaceOfNodes.hxx"
32 #include "SMDS_VolumeTool.hxx"
33 #include "SMESHDS_GroupBase.hxx"
34 #include "SMESHDS_GroupOnFilter.hxx"
35 #include "SMESHDS_Mesh.hxx"
36 #include "SMESH_MeshAlgos.hxx"
37 #include "SMESH_OctreeNode.hxx"
39 #include <Basics_Utils.hxx>
41 #include <BRepAdaptor_Surface.hxx>
42 #include <BRepBndLib.hxx>
43 #include <BRepBuilderAPI_Copy.hxx>
44 #include <BRepClass_FaceClassifier.hxx>
45 #include <BRep_Tool.hxx>
46 #include <Geom_CylindricalSurface.hxx>
47 #include <Geom_Plane.hxx>
48 #include <Geom_Surface.hxx>
49 #include <NCollection_Map.hxx>
50 #include <Precision.hxx>
51 #include <TColStd_MapIteratorOfMapOfInteger.hxx>
52 #include <TColStd_MapOfInteger.hxx>
53 #include <TColStd_SequenceOfAsciiString.hxx>
54 #include <TColgp_Array1OfXYZ.hxx>
58 #include <TopoDS_Edge.hxx>
59 #include <TopoDS_Face.hxx>
60 #include <TopoDS_Iterator.hxx>
61 #include <TopoDS_Shape.hxx>
62 #include <TopoDS_Vertex.hxx>
64 #include <gp_Cylinder.hxx>
71 #include <vtkMeshQuality.h>
82 const double theEps = 1e-100;
83 const double theInf = 1e+100;
85 inline gp_XYZ gpXYZ(const SMDS_MeshNode* aNode )
87 return gp_XYZ(aNode->X(), aNode->Y(), aNode->Z() );
90 inline double getAngle( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
92 gp_Vec v1( P1 - P2 ), v2( P3 - P2 );
94 return v1.Magnitude() < gp::Resolution() ||
95 v2.Magnitude() < gp::Resolution() ? 0 : v1.Angle( v2 );
98 inline double getArea( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
100 gp_Vec aVec1( P2 - P1 );
101 gp_Vec aVec2( P3 - P1 );
102 return ( aVec1 ^ aVec2 ).Magnitude() * 0.5;
105 inline double getArea( const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3 )
107 return getArea( P1.XYZ(), P2.XYZ(), P3.XYZ() );
112 inline double getDistance( const gp_XYZ& P1, const gp_XYZ& P2 )
114 double aDist = gp_Pnt( P1 ).Distance( gp_Pnt( P2 ) );
118 int getNbMultiConnection( const SMDS_Mesh* theMesh, const int theId )
123 const SMDS_MeshElement* anEdge = theMesh->FindElement( theId );
124 if ( anEdge == 0 || anEdge->GetType() != SMDSAbs_Edge/* || anEdge->NbNodes() != 2 */)
127 // for each pair of nodes in anEdge (there are 2 pairs in a quadratic edge)
128 // count elements containing both nodes of the pair.
129 // Note that there may be such cases for a quadratic edge (a horizontal line):
134 // +-----+------+ +-----+------+
137 // result sould be 2 in both cases
139 int aResult0 = 0, aResult1 = 0;
140 // last node, it is a medium one in a quadratic edge
141 const SMDS_MeshNode* aLastNode = anEdge->GetNode( anEdge->NbNodes() - 1 );
142 const SMDS_MeshNode* aNode0 = anEdge->GetNode( 0 );
143 const SMDS_MeshNode* aNode1 = anEdge->GetNode( 1 );
144 if ( aNode1 == aLastNode ) aNode1 = 0;
146 SMDS_ElemIteratorPtr anElemIter = aLastNode->GetInverseElementIterator();
147 while( anElemIter->more() ) {
148 const SMDS_MeshElement* anElem = anElemIter->next();
149 if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
150 SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
151 while ( anIter->more() ) {
152 if ( const SMDS_MeshElement* anElemNode = anIter->next() ) {
153 if ( anElemNode == aNode0 ) {
155 if ( !aNode1 ) break; // not a quadratic edge
157 else if ( anElemNode == aNode1 )
163 int aResult = std::max ( aResult0, aResult1 );
168 gp_XYZ getNormale( const SMDS_MeshFace* theFace, bool* ok=0 )
170 int aNbNode = theFace->NbNodes();
172 gp_XYZ q1 = gpXYZ( theFace->GetNode(1)) - gpXYZ( theFace->GetNode(0));
173 gp_XYZ q2 = gpXYZ( theFace->GetNode(2)) - gpXYZ( theFace->GetNode(0));
176 gp_XYZ q3 = gpXYZ( theFace->GetNode(3)) - gpXYZ( theFace->GetNode(0));
179 double len = n.Modulus();
180 bool zeroLen = ( len <= std::numeric_limits<double>::min());
184 if (ok) *ok = !zeroLen;
192 using namespace SMESH::Controls;
198 //================================================================================
200 Class : NumericalFunctor
201 Description : Base class for numerical functors
203 //================================================================================
205 NumericalFunctor::NumericalFunctor():
211 void NumericalFunctor::SetMesh( const SMDS_Mesh* theMesh )
216 bool NumericalFunctor::GetPoints(const int theId,
217 TSequenceOfXYZ& theRes ) const
224 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
225 if ( !anElem || anElem->GetType() != this->GetType() )
228 return GetPoints( anElem, theRes );
231 bool NumericalFunctor::GetPoints(const SMDS_MeshElement* anElem,
232 TSequenceOfXYZ& theRes )
239 theRes.reserve( anElem->NbNodes() );
240 theRes.setElement( anElem );
242 // Get nodes of the element
243 SMDS_ElemIteratorPtr anIter;
245 if ( anElem->IsQuadratic() ) {
246 switch ( anElem->GetType() ) {
248 anIter = dynamic_cast<const SMDS_VtkEdge*>
249 (anElem)->interlacedNodesElemIterator();
252 anIter = dynamic_cast<const SMDS_VtkFace*>
253 (anElem)->interlacedNodesElemIterator();
256 anIter = anElem->nodesIterator();
260 anIter = anElem->nodesIterator();
265 while( anIter->more() ) {
266 if ( const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>( anIter->next() ))
268 aNode->GetXYZ( xyz );
269 theRes.push_back( gp_XYZ( xyz[0], xyz[1], xyz[2] ));
277 long NumericalFunctor::GetPrecision() const
282 void NumericalFunctor::SetPrecision( const long thePrecision )
284 myPrecision = thePrecision;
285 myPrecisionValue = pow( 10., (double)( myPrecision ) );
288 double NumericalFunctor::GetValue( long theId )
292 myCurrElement = myMesh->FindElement( theId );
295 if ( GetPoints( theId, P )) // elem type is checked here
296 aVal = Round( GetValue( P ));
301 double NumericalFunctor::Round( const double & aVal )
303 return ( myPrecision >= 0 ) ? floor( aVal * myPrecisionValue + 0.5 ) / myPrecisionValue : aVal;
306 //================================================================================
308 * \brief Return histogram of functor values
309 * \param nbIntervals - number of intervals
310 * \param nbEvents - number of mesh elements having values within i-th interval
311 * \param funValues - boundaries of intervals
312 * \param elements - elements to check vulue of; empty list means "of all"
313 * \param minmax - boundaries of diapason of values to divide into intervals
315 //================================================================================
317 void NumericalFunctor::GetHistogram(int nbIntervals,
318 std::vector<int>& nbEvents,
319 std::vector<double>& funValues,
320 const std::vector<int>& elements,
321 const double* minmax,
322 const bool isLogarithmic)
324 if ( nbIntervals < 1 ||
326 !myMesh->GetMeshInfo().NbElements( GetType() ))
328 nbEvents.resize( nbIntervals, 0 );
329 funValues.resize( nbIntervals+1 );
331 // get all values sorted
332 std::multiset< double > values;
333 if ( elements.empty() )
335 SMDS_ElemIteratorPtr elemIt = myMesh->elementsIterator( GetType() );
336 while ( elemIt->more() )
337 values.insert( GetValue( elemIt->next()->GetID() ));
341 std::vector<int>::const_iterator id = elements.begin();
342 for ( ; id != elements.end(); ++id )
343 values.insert( GetValue( *id ));
348 funValues[0] = minmax[0];
349 funValues[nbIntervals] = minmax[1];
353 funValues[0] = *values.begin();
354 funValues[nbIntervals] = *values.rbegin();
356 // case nbIntervals == 1
357 if ( nbIntervals == 1 )
359 nbEvents[0] = values.size();
363 if (funValues.front() == funValues.back())
365 nbEvents.resize( 1 );
366 nbEvents[0] = values.size();
367 funValues[1] = funValues.back();
368 funValues.resize( 2 );
371 std::multiset< double >::iterator min = values.begin(), max;
372 for ( int i = 0; i < nbIntervals; ++i )
374 // find end value of i-th interval
375 double r = (i+1) / double(nbIntervals);
376 if (isLogarithmic && funValues.front() > 1e-07 && funValues.back() > 1e-07) {
377 double logmin = log10(funValues.front());
378 double lval = logmin + r * (log10(funValues.back()) - logmin);
379 funValues[i+1] = pow(10.0, lval);
382 funValues[i+1] = funValues.front() * (1-r) + funValues.back() * r;
385 // count values in the i-th interval if there are any
386 if ( min != values.end() && *min <= funValues[i+1] )
388 // find the first value out of the interval
389 max = values.upper_bound( funValues[i+1] ); // max is greater than funValues[i+1], or end()
390 nbEvents[i] = std::distance( min, max );
394 // add values larger than minmax[1]
395 nbEvents.back() += std::distance( min, values.end() );
398 //=======================================================================
401 Description : Functor calculating volume of a 3D element
403 //================================================================================
405 double Volume::GetValue( long theElementId )
407 if ( theElementId && myMesh ) {
408 SMDS_VolumeTool aVolumeTool;
409 if ( aVolumeTool.Set( myMesh->FindElement( theElementId )))
410 return aVolumeTool.GetSize();
415 double Volume::GetBadRate( double Value, int /*nbNodes*/ ) const
420 SMDSAbs_ElementType Volume::GetType() const
422 return SMDSAbs_Volume;
425 //=======================================================================
427 Class : MaxElementLength2D
428 Description : Functor calculating maximum length of 2D element
430 //================================================================================
432 double MaxElementLength2D::GetValue( const TSequenceOfXYZ& P )
438 if( len == 3 ) { // triangles
439 double L1 = getDistance(P( 1 ),P( 2 ));
440 double L2 = getDistance(P( 2 ),P( 3 ));
441 double L3 = getDistance(P( 3 ),P( 1 ));
442 aVal = Max(L1,Max(L2,L3));
444 else if( len == 4 ) { // quadrangles
445 double L1 = getDistance(P( 1 ),P( 2 ));
446 double L2 = getDistance(P( 2 ),P( 3 ));
447 double L3 = getDistance(P( 3 ),P( 4 ));
448 double L4 = getDistance(P( 4 ),P( 1 ));
449 double D1 = getDistance(P( 1 ),P( 3 ));
450 double D2 = getDistance(P( 2 ),P( 4 ));
451 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
453 else if( len == 6 ) { // quadratic triangles
454 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
455 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
456 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
457 aVal = Max(L1,Max(L2,L3));
459 else if( len == 8 || len == 9 ) { // quadratic quadrangles
460 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
461 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
462 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
463 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
464 double D1 = getDistance(P( 1 ),P( 5 ));
465 double D2 = getDistance(P( 3 ),P( 7 ));
466 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
468 // Diagonals are undefined for concave polygons
469 // else if ( P.getElementEntity() == SMDSEntity_Quad_Polygon && P.size() > 2 ) // quad polygon
472 // aVal = getDistance( P( 1 ), P( P.size() )) + getDistance( P( P.size() ), P( P.size()-1 ));
473 // for ( size_t i = 1; i < P.size()-1; i += 2 )
475 // double L = getDistance( P( i ), P( i+1 )) + getDistance( P( i+1 ), P( i+2 ));
476 // aVal = Max( aVal, L );
479 // for ( int i = P.size()-5; i > 0; i -= 2 )
480 // for ( int j = i + 4; j < P.size() + i - 2; i += 2 )
482 // double D = getDistance( P( i ), P( j ));
483 // aVal = Max( aVal, D );
490 if( myPrecision >= 0 )
492 double prec = pow( 10., (double)myPrecision );
493 aVal = floor( aVal * prec + 0.5 ) / prec;
498 double MaxElementLength2D::GetValue( long theElementId )
501 return GetPoints( theElementId, P ) ? GetValue(P) : 0.0;
504 double MaxElementLength2D::GetBadRate( double Value, int /*nbNodes*/ ) const
509 SMDSAbs_ElementType MaxElementLength2D::GetType() const
514 //=======================================================================
516 Class : MaxElementLength3D
517 Description : Functor calculating maximum length of 3D element
519 //================================================================================
521 double MaxElementLength3D::GetValue( long theElementId )
524 if( GetPoints( theElementId, P ) ) {
526 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
527 SMDSAbs_EntityType aType = aElem->GetEntityType();
530 case SMDSEntity_Tetra: { // tetras
531 double L1 = getDistance(P( 1 ),P( 2 ));
532 double L2 = getDistance(P( 2 ),P( 3 ));
533 double L3 = getDistance(P( 3 ),P( 1 ));
534 double L4 = getDistance(P( 1 ),P( 4 ));
535 double L5 = getDistance(P( 2 ),P( 4 ));
536 double L6 = getDistance(P( 3 ),P( 4 ));
537 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
540 case SMDSEntity_Pyramid: { // pyramids
541 double L1 = getDistance(P( 1 ),P( 2 ));
542 double L2 = getDistance(P( 2 ),P( 3 ));
543 double L3 = getDistance(P( 3 ),P( 4 ));
544 double L4 = getDistance(P( 4 ),P( 1 ));
545 double L5 = getDistance(P( 1 ),P( 5 ));
546 double L6 = getDistance(P( 2 ),P( 5 ));
547 double L7 = getDistance(P( 3 ),P( 5 ));
548 double L8 = getDistance(P( 4 ),P( 5 ));
549 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
550 aVal = Max(aVal,Max(L7,L8));
553 case SMDSEntity_Penta: { // pentas
554 double L1 = getDistance(P( 1 ),P( 2 ));
555 double L2 = getDistance(P( 2 ),P( 3 ));
556 double L3 = getDistance(P( 3 ),P( 1 ));
557 double L4 = getDistance(P( 4 ),P( 5 ));
558 double L5 = getDistance(P( 5 ),P( 6 ));
559 double L6 = getDistance(P( 6 ),P( 4 ));
560 double L7 = getDistance(P( 1 ),P( 4 ));
561 double L8 = getDistance(P( 2 ),P( 5 ));
562 double L9 = getDistance(P( 3 ),P( 6 ));
563 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
564 aVal = Max(aVal,Max(Max(L7,L8),L9));
567 case SMDSEntity_Hexa: { // hexas
568 double L1 = getDistance(P( 1 ),P( 2 ));
569 double L2 = getDistance(P( 2 ),P( 3 ));
570 double L3 = getDistance(P( 3 ),P( 4 ));
571 double L4 = getDistance(P( 4 ),P( 1 ));
572 double L5 = getDistance(P( 5 ),P( 6 ));
573 double L6 = getDistance(P( 6 ),P( 7 ));
574 double L7 = getDistance(P( 7 ),P( 8 ));
575 double L8 = getDistance(P( 8 ),P( 5 ));
576 double L9 = getDistance(P( 1 ),P( 5 ));
577 double L10= getDistance(P( 2 ),P( 6 ));
578 double L11= getDistance(P( 3 ),P( 7 ));
579 double L12= getDistance(P( 4 ),P( 8 ));
580 double D1 = getDistance(P( 1 ),P( 7 ));
581 double D2 = getDistance(P( 2 ),P( 8 ));
582 double D3 = getDistance(P( 3 ),P( 5 ));
583 double D4 = getDistance(P( 4 ),P( 6 ));
584 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
585 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
586 aVal = Max(aVal,Max(L11,L12));
587 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
590 case SMDSEntity_Hexagonal_Prism: { // hexagonal prism
591 for ( int i1 = 1; i1 < 12; ++i1 )
592 for ( int i2 = i1+1; i1 <= 12; ++i1 )
593 aVal = Max( aVal, getDistance(P( i1 ),P( i2 )));
596 case SMDSEntity_Quad_Tetra: { // quadratic tetras
597 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
598 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
599 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
600 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
601 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
602 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
603 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
606 case SMDSEntity_Quad_Pyramid: { // quadratic pyramids
607 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
608 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
609 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
610 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
611 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
612 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
613 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
614 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
615 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
616 aVal = Max(aVal,Max(L7,L8));
619 case SMDSEntity_Quad_Penta: { // quadratic pentas
620 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
621 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
622 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
623 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
624 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
625 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
626 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
627 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
628 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
629 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
630 aVal = Max(aVal,Max(Max(L7,L8),L9));
633 case SMDSEntity_Quad_Hexa:
634 case SMDSEntity_TriQuad_Hexa: { // quadratic hexas
635 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
636 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
637 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
638 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
639 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
640 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
641 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
642 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
643 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
644 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
645 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
646 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
647 double D1 = getDistance(P( 1 ),P( 7 ));
648 double D2 = getDistance(P( 2 ),P( 8 ));
649 double D3 = getDistance(P( 3 ),P( 5 ));
650 double D4 = getDistance(P( 4 ),P( 6 ));
651 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
652 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
653 aVal = Max(aVal,Max(L11,L12));
654 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
657 case SMDSEntity_Quad_Polyhedra:
658 case SMDSEntity_Polyhedra: { // polys
659 // get the maximum distance between all pairs of nodes
660 for( int i = 1; i <= len; i++ ) {
661 for( int j = 1; j <= len; j++ ) {
662 if( j > i ) { // optimization of the loop
663 double D = getDistance( P(i), P(j) );
664 aVal = Max( aVal, D );
670 case SMDSEntity_Node:
672 case SMDSEntity_Edge:
673 case SMDSEntity_Quad_Edge:
674 case SMDSEntity_Triangle:
675 case SMDSEntity_Quad_Triangle:
676 case SMDSEntity_BiQuad_Triangle:
677 case SMDSEntity_Quadrangle:
678 case SMDSEntity_Quad_Quadrangle:
679 case SMDSEntity_BiQuad_Quadrangle:
680 case SMDSEntity_Polygon:
681 case SMDSEntity_Quad_Polygon:
682 case SMDSEntity_Ball:
683 case SMDSEntity_Last: return 0;
684 } // switch ( aType )
686 if( myPrecision >= 0 )
688 double prec = pow( 10., (double)myPrecision );
689 aVal = floor( aVal * prec + 0.5 ) / prec;
696 double MaxElementLength3D::GetBadRate( double Value, int /*nbNodes*/ ) const
701 SMDSAbs_ElementType MaxElementLength3D::GetType() const
703 return SMDSAbs_Volume;
706 //=======================================================================
709 Description : Functor for calculation of minimum angle
711 //================================================================================
713 double MinimumAngle::GetValue( const TSequenceOfXYZ& P )
720 aMin = getAngle(P( P.size() ), P( 1 ), P( 2 ));
721 aMin = Min(aMin,getAngle(P( P.size()-1 ), P( P.size() ), P( 1 )));
723 for ( size_t i = 2; i < P.size(); i++ )
725 double A0 = getAngle( P( i-1 ), P( i ), P( i+1 ) );
729 return aMin * 180.0 / M_PI;
732 double MinimumAngle::GetBadRate( double Value, int nbNodes ) const
734 //const double aBestAngle = PI / nbNodes;
735 const double aBestAngle = 180.0 - ( 360.0 / double(nbNodes) );
736 return ( fabs( aBestAngle - Value ));
739 SMDSAbs_ElementType MinimumAngle::GetType() const
745 //================================================================================
748 Description : Functor for calculating aspect ratio
750 //================================================================================
752 double AspectRatio::GetValue( long theId )
755 myCurrElement = myMesh->FindElement( theId );
756 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_QUAD )
759 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
760 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
761 aVal = Round( vtkMeshQuality::QuadAspectRatio( avtkCell ));
766 if ( GetPoints( myCurrElement, P ))
767 aVal = Round( GetValue( P ));
772 double AspectRatio::GetValue( const TSequenceOfXYZ& P )
774 // According to "Mesh quality control" by Nadir Bouhamau referring to
775 // Pascal Jean Frey and Paul-Louis George. Maillages, applications aux elements finis.
776 // Hermes Science publications, Paris 1999 ISBN 2-7462-0024-4
779 int nbNodes = P.size();
784 // Compute aspect ratio
786 if ( nbNodes == 3 ) {
787 // Compute lengths of the sides
788 std::vector< double > aLen (nbNodes);
789 for ( int i = 0; i < nbNodes - 1; i++ )
790 aLen[ i ] = getDistance( P( i + 1 ), P( i + 2 ) );
791 aLen[ nbNodes - 1 ] = getDistance( P( 1 ), P( nbNodes ) );
792 // Q = alfa * h * p / S, where
794 // alfa = sqrt( 3 ) / 6
795 // h - length of the longest edge
796 // p - half perimeter
797 // S - triangle surface
798 const double alfa = sqrt( 3. ) / 6.;
799 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
800 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
801 double anArea = getArea( P( 1 ), P( 2 ), P( 3 ) );
802 if ( anArea <= theEps )
804 return alfa * maxLen * half_perimeter / anArea;
806 else if ( nbNodes == 6 ) { // quadratic triangles
807 // Compute lengths of the sides
808 std::vector< double > aLen (3);
809 aLen[0] = getDistance( P(1), P(3) );
810 aLen[1] = getDistance( P(3), P(5) );
811 aLen[2] = getDistance( P(5), P(1) );
812 // Q = alfa * h * p / S, where
814 // alfa = sqrt( 3 ) / 6
815 // h - length of the longest edge
816 // p - half perimeter
817 // S - triangle surface
818 const double alfa = sqrt( 3. ) / 6.;
819 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
820 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
821 double anArea = getArea( P(1), P(3), P(5) );
822 if ( anArea <= theEps )
824 return alfa * maxLen * half_perimeter / anArea;
826 else if( nbNodes == 4 ) { // quadrangle
827 // Compute lengths of the sides
828 std::vector< double > aLen (4);
829 aLen[0] = getDistance( P(1), P(2) );
830 aLen[1] = getDistance( P(2), P(3) );
831 aLen[2] = getDistance( P(3), P(4) );
832 aLen[3] = getDistance( P(4), P(1) );
833 // Compute lengths of the diagonals
834 std::vector< double > aDia (2);
835 aDia[0] = getDistance( P(1), P(3) );
836 aDia[1] = getDistance( P(2), P(4) );
837 // Compute areas of all triangles which can be built
838 // taking three nodes of the quadrangle
839 std::vector< double > anArea (4);
840 anArea[0] = getArea( P(1), P(2), P(3) );
841 anArea[1] = getArea( P(1), P(2), P(4) );
842 anArea[2] = getArea( P(1), P(3), P(4) );
843 anArea[3] = getArea( P(2), P(3), P(4) );
844 // Q = alpha * L * C1 / C2, where
846 // alpha = sqrt( 1/32 )
847 // L = max( L1, L2, L3, L4, D1, D2 )
848 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
849 // C2 = min( S1, S2, S3, S4 )
850 // Li - lengths of the edges
851 // Di - lengths of the diagonals
852 // Si - areas of the triangles
853 const double alpha = sqrt( 1 / 32. );
854 double L = Max( aLen[ 0 ],
858 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
859 double C1 = sqrt( ( aLen[0] * aLen[0] +
862 aLen[3] * aLen[3] ) / 4. );
863 double C2 = Min( anArea[ 0 ],
865 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
868 return alpha * L * C1 / C2;
870 else if( nbNodes == 8 || nbNodes == 9 ) { // nbNodes==8 - quadratic quadrangle
871 // Compute lengths of the sides
872 std::vector< double > aLen (4);
873 aLen[0] = getDistance( P(1), P(3) );
874 aLen[1] = getDistance( P(3), P(5) );
875 aLen[2] = getDistance( P(5), P(7) );
876 aLen[3] = getDistance( P(7), P(1) );
877 // Compute lengths of the diagonals
878 std::vector< double > aDia (2);
879 aDia[0] = getDistance( P(1), P(5) );
880 aDia[1] = getDistance( P(3), P(7) );
881 // Compute areas of all triangles which can be built
882 // taking three nodes of the quadrangle
883 std::vector< double > anArea (4);
884 anArea[0] = getArea( P(1), P(3), P(5) );
885 anArea[1] = getArea( P(1), P(3), P(7) );
886 anArea[2] = getArea( P(1), P(5), P(7) );
887 anArea[3] = getArea( P(3), P(5), P(7) );
888 // Q = alpha * L * C1 / C2, where
890 // alpha = sqrt( 1/32 )
891 // L = max( L1, L2, L3, L4, D1, D2 )
892 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
893 // C2 = min( S1, S2, S3, S4 )
894 // Li - lengths of the edges
895 // Di - lengths of the diagonals
896 // Si - areas of the triangles
897 const double alpha = sqrt( 1 / 32. );
898 double L = Max( aLen[ 0 ],
902 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
903 double C1 = sqrt( ( aLen[0] * aLen[0] +
906 aLen[3] * aLen[3] ) / 4. );
907 double C2 = Min( anArea[ 0 ],
909 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
912 return alpha * L * C1 / C2;
917 double AspectRatio::GetBadRate( double Value, int /*nbNodes*/ ) const
919 // the aspect ratio is in the range [1.0,infinity]
920 // < 1.0 = very bad, zero area
923 return ( Value < 0.9 ) ? 1000 : Value / 1000.;
926 SMDSAbs_ElementType AspectRatio::GetType() const
932 //================================================================================
934 Class : AspectRatio3D
935 Description : Functor for calculating aspect ratio
937 //================================================================================
941 inline double getHalfPerimeter(double theTria[3]){
942 return (theTria[0] + theTria[1] + theTria[2])/2.0;
945 inline double getArea(double theHalfPerim, double theTria[3]){
946 return sqrt(theHalfPerim*
947 (theHalfPerim-theTria[0])*
948 (theHalfPerim-theTria[1])*
949 (theHalfPerim-theTria[2]));
952 inline double getVolume(double theLen[6]){
953 double a2 = theLen[0]*theLen[0];
954 double b2 = theLen[1]*theLen[1];
955 double c2 = theLen[2]*theLen[2];
956 double d2 = theLen[3]*theLen[3];
957 double e2 = theLen[4]*theLen[4];
958 double f2 = theLen[5]*theLen[5];
959 double P = 4.0*a2*b2*d2;
960 double Q = a2*(b2+d2-e2)-b2*(a2+d2-f2)-d2*(a2+b2-c2);
961 double R = (b2+d2-e2)*(a2+d2-f2)*(a2+d2-f2);
962 return sqrt(P-Q+R)/12.0;
965 inline double getVolume2(double theLen[6]){
966 double a2 = theLen[0]*theLen[0];
967 double b2 = theLen[1]*theLen[1];
968 double c2 = theLen[2]*theLen[2];
969 double d2 = theLen[3]*theLen[3];
970 double e2 = theLen[4]*theLen[4];
971 double f2 = theLen[5]*theLen[5];
973 double P = a2*e2*(b2+c2+d2+f2-a2-e2);
974 double Q = b2*f2*(a2+c2+d2+e2-b2-f2);
975 double R = c2*d2*(a2+b2+e2+f2-c2-d2);
976 double S = a2*b2*d2+b2*c2*e2+a2*c2*f2+d2*e2*f2;
978 return sqrt(P+Q+R-S)/12.0;
981 inline double getVolume(const TSequenceOfXYZ& P){
982 gp_Vec aVec1( P( 2 ) - P( 1 ) );
983 gp_Vec aVec2( P( 3 ) - P( 1 ) );
984 gp_Vec aVec3( P( 4 ) - P( 1 ) );
985 gp_Vec anAreaVec( aVec1 ^ aVec2 );
986 return fabs(aVec3 * anAreaVec) / 6.0;
989 inline double getMaxHeight(double theLen[6])
991 double aHeight = std::max(theLen[0],theLen[1]);
992 aHeight = std::max(aHeight,theLen[2]);
993 aHeight = std::max(aHeight,theLen[3]);
994 aHeight = std::max(aHeight,theLen[4]);
995 aHeight = std::max(aHeight,theLen[5]);
1001 double AspectRatio3D::GetValue( long theId )
1004 myCurrElement = myMesh->FindElement( theId );
1005 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_TETRA )
1007 // Action from CoTech | ACTION 31.3:
1008 // EURIWARE BO: Homogenize the formulas used to calculate the Controls in SMESH to fit with
1009 // those of ParaView. The library used by ParaView for those calculations can be reused in SMESH.
1010 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
1011 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
1012 aVal = Round( vtkMeshQuality::TetAspectRatio( avtkCell ));
1017 if ( GetPoints( myCurrElement, P ))
1018 aVal = Round( GetValue( P ));
1023 double AspectRatio3D::GetValue( const TSequenceOfXYZ& P )
1025 double aQuality = 0.0;
1026 if(myCurrElement->IsPoly()) return aQuality;
1028 int nbNodes = P.size();
1030 if(myCurrElement->IsQuadratic()) {
1031 if(nbNodes==10) nbNodes=4; // quadratic tetrahedron
1032 else if(nbNodes==13) nbNodes=5; // quadratic pyramid
1033 else if(nbNodes==15) nbNodes=6; // quadratic pentahedron
1034 else if(nbNodes==20) nbNodes=8; // quadratic hexahedron
1035 else if(nbNodes==27) nbNodes=8; // quadratic hexahedron
1036 else return aQuality;
1042 getDistance(P( 1 ),P( 2 )), // a
1043 getDistance(P( 2 ),P( 3 )), // b
1044 getDistance(P( 3 ),P( 1 )), // c
1045 getDistance(P( 2 ),P( 4 )), // d
1046 getDistance(P( 3 ),P( 4 )), // e
1047 getDistance(P( 1 ),P( 4 )) // f
1049 double aTria[4][3] = {
1050 {aLen[0],aLen[1],aLen[2]}, // abc
1051 {aLen[0],aLen[3],aLen[5]}, // adf
1052 {aLen[1],aLen[3],aLen[4]}, // bde
1053 {aLen[2],aLen[4],aLen[5]} // cef
1055 double aSumArea = 0.0;
1056 double aHalfPerimeter = getHalfPerimeter(aTria[0]);
1057 double anArea = getArea(aHalfPerimeter,aTria[0]);
1059 aHalfPerimeter = getHalfPerimeter(aTria[1]);
1060 anArea = getArea(aHalfPerimeter,aTria[1]);
1062 aHalfPerimeter = getHalfPerimeter(aTria[2]);
1063 anArea = getArea(aHalfPerimeter,aTria[2]);
1065 aHalfPerimeter = getHalfPerimeter(aTria[3]);
1066 anArea = getArea(aHalfPerimeter,aTria[3]);
1068 double aVolume = getVolume(P);
1069 //double aVolume = getVolume(aLen);
1070 double aHeight = getMaxHeight(aLen);
1071 static double aCoeff = sqrt(2.0)/12.0;
1072 if ( aVolume > DBL_MIN )
1073 aQuality = aCoeff*aHeight*aSumArea/aVolume;
1078 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 3 ),P( 5 )};
1079 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1082 gp_XYZ aXYZ[4] = {P( 1 ),P( 3 ),P( 4 ),P( 5 )};
1083 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1086 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 5 )};
1087 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1090 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 4 ),P( 5 )};
1091 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1097 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 6 )};
1098 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1101 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 3 )};
1102 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1105 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 6 )};
1106 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1109 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1110 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1113 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 6 )};
1114 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1117 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 3 )};
1118 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1124 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1125 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1128 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 4 )};
1129 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1132 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 7 )};
1133 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1136 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 8 )};
1137 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1140 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 3 )};
1141 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1144 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 4 )};
1145 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1148 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 7 )};
1149 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1152 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 8 )};
1153 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1156 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 3 )};
1157 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1160 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 4 )};
1161 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1164 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 7 )};
1165 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1168 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 8 )};
1169 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1172 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 1 )};
1173 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1176 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 2 )};
1177 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1180 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 5 )};
1181 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1184 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 6 )};
1185 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1188 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 1 )};
1189 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1192 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 2 )};
1193 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1196 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 5 )};
1197 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1200 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 6 )};
1201 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1204 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 1 )};
1205 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1208 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1209 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1212 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 5 )};
1213 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1216 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 6 )};
1217 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1220 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1221 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1224 gp_XYZ aXYZ[4] = {P( 4 ),P( 5 ),P( 8 ),P( 2 )};
1225 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1228 gp_XYZ aXYZ[4] = {P( 1 ),P( 4 ),P( 5 ),P( 3 )};
1229 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1232 gp_XYZ aXYZ[4] = {P( 3 ),P( 6 ),P( 7 ),P( 1 )};
1233 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1236 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 6 ),P( 4 )};
1237 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1240 gp_XYZ aXYZ[4] = {P( 5 ),P( 6 ),P( 8 ),P( 3 )};
1241 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1244 gp_XYZ aXYZ[4] = {P( 7 ),P( 8 ),P( 6 ),P( 1 )};
1245 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1248 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 7 )};
1249 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1252 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 2 ),P( 5 )};
1253 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1259 gp_XYZ aXYZ[8] = {P( 1 ),P( 2 ),P( 4 ),P( 5 ),P( 7 ),P( 8 ),P( 10 ),P( 11 )};
1260 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1263 gp_XYZ aXYZ[8] = {P( 2 ),P( 3 ),P( 5 ),P( 6 ),P( 8 ),P( 9 ),P( 11 ),P( 12 )};
1264 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1267 gp_XYZ aXYZ[8] = {P( 3 ),P( 4 ),P( 6 ),P( 1 ),P( 9 ),P( 10 ),P( 12 ),P( 7 )};
1268 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1271 } // switch(nbNodes)
1273 if ( nbNodes > 4 ) {
1274 // avaluate aspect ratio of quadranle faces
1275 AspectRatio aspect2D;
1276 SMDS_VolumeTool::VolumeType type = SMDS_VolumeTool::GetType( nbNodes );
1277 int nbFaces = SMDS_VolumeTool::NbFaces( type );
1278 TSequenceOfXYZ points(4);
1279 for ( int i = 0; i < nbFaces; ++i ) { // loop on faces of a volume
1280 if ( SMDS_VolumeTool::NbFaceNodes( type, i ) != 4 )
1282 const int* pInd = SMDS_VolumeTool::GetFaceNodesIndices( type, i, true );
1283 for ( int p = 0; p < 4; ++p ) // loop on nodes of a quadranle face
1284 points( p + 1 ) = P( pInd[ p ] + 1 );
1285 aQuality = std::max( aQuality, aspect2D.GetValue( points ));
1291 double AspectRatio3D::GetBadRate( double Value, int /*nbNodes*/ ) const
1293 // the aspect ratio is in the range [1.0,infinity]
1296 return Value / 1000.;
1299 SMDSAbs_ElementType AspectRatio3D::GetType() const
1301 return SMDSAbs_Volume;
1305 //================================================================================
1308 Description : Functor for calculating warping
1310 //================================================================================
1312 double Warping::GetValue( const TSequenceOfXYZ& P )
1314 if ( P.size() != 4 )
1317 gp_XYZ G = ( P( 1 ) + P( 2 ) + P( 3 ) + P( 4 ) ) / 4.;
1319 double A1 = ComputeA( P( 1 ), P( 2 ), P( 3 ), G );
1320 double A2 = ComputeA( P( 2 ), P( 3 ), P( 4 ), G );
1321 double A3 = ComputeA( P( 3 ), P( 4 ), P( 1 ), G );
1322 double A4 = ComputeA( P( 4 ), P( 1 ), P( 2 ), G );
1324 double val = Max( Max( A1, A2 ), Max( A3, A4 ) );
1326 const double eps = 0.1; // val is in degrees
1328 return val < eps ? 0. : val;
1331 double Warping::ComputeA( const gp_XYZ& thePnt1,
1332 const gp_XYZ& thePnt2,
1333 const gp_XYZ& thePnt3,
1334 const gp_XYZ& theG ) const
1336 double aLen1 = gp_Pnt( thePnt1 ).Distance( gp_Pnt( thePnt2 ) );
1337 double aLen2 = gp_Pnt( thePnt2 ).Distance( gp_Pnt( thePnt3 ) );
1338 double L = Min( aLen1, aLen2 ) * 0.5;
1342 gp_XYZ GI = ( thePnt2 + thePnt1 ) / 2. - theG;
1343 gp_XYZ GJ = ( thePnt3 + thePnt2 ) / 2. - theG;
1344 gp_XYZ N = GI.Crossed( GJ );
1346 if ( N.Modulus() < gp::Resolution() )
1351 double H = ( thePnt2 - theG ).Dot( N );
1352 return asin( fabs( H / L ) ) * 180. / M_PI;
1355 double Warping::GetBadRate( double Value, int /*nbNodes*/ ) const
1357 // the warp is in the range [0.0,PI/2]
1358 // 0.0 = good (no warp)
1359 // PI/2 = bad (face pliee)
1363 SMDSAbs_ElementType Warping::GetType() const
1365 return SMDSAbs_Face;
1369 //================================================================================
1372 Description : Functor for calculating taper
1374 //================================================================================
1376 double Taper::GetValue( const TSequenceOfXYZ& P )
1378 if ( P.size() != 4 )
1382 double J1 = getArea( P( 4 ), P( 1 ), P( 2 ) );
1383 double J2 = getArea( P( 3 ), P( 1 ), P( 2 ) );
1384 double J3 = getArea( P( 2 ), P( 3 ), P( 4 ) );
1385 double J4 = getArea( P( 3 ), P( 4 ), P( 1 ) );
1387 double JA = 0.25 * ( J1 + J2 + J3 + J4 );
1391 double T1 = fabs( ( J1 - JA ) / JA );
1392 double T2 = fabs( ( J2 - JA ) / JA );
1393 double T3 = fabs( ( J3 - JA ) / JA );
1394 double T4 = fabs( ( J4 - JA ) / JA );
1396 double val = Max( Max( T1, T2 ), Max( T3, T4 ) );
1398 const double eps = 0.01;
1400 return val < eps ? 0. : val;
1403 double Taper::GetBadRate( double Value, int /*nbNodes*/ ) const
1405 // the taper is in the range [0.0,1.0]
1406 // 0.0 = good (no taper)
1407 // 1.0 = bad (les cotes opposes sont allignes)
1411 SMDSAbs_ElementType Taper::GetType() const
1413 return SMDSAbs_Face;
1416 //================================================================================
1419 Description : Functor for calculating skew in degrees
1421 //================================================================================
1423 static inline double skewAngle( const gp_XYZ& p1, const gp_XYZ& p2, const gp_XYZ& p3 )
1425 gp_XYZ p12 = ( p2 + p1 ) / 2.;
1426 gp_XYZ p23 = ( p3 + p2 ) / 2.;
1427 gp_XYZ p31 = ( p3 + p1 ) / 2.;
1429 gp_Vec v1( p31 - p2 ), v2( p12 - p23 );
1431 return v1.Magnitude() < gp::Resolution() || v2.Magnitude() < gp::Resolution() ? 0. : v1.Angle( v2 );
1434 double Skew::GetValue( const TSequenceOfXYZ& P )
1436 if ( P.size() != 3 && P.size() != 4 )
1440 const double PI2 = M_PI / 2.;
1441 if ( P.size() == 3 )
1443 double A0 = fabs( PI2 - skewAngle( P( 3 ), P( 1 ), P( 2 ) ) );
1444 double A1 = fabs( PI2 - skewAngle( P( 1 ), P( 2 ), P( 3 ) ) );
1445 double A2 = fabs( PI2 - skewAngle( P( 2 ), P( 3 ), P( 1 ) ) );
1447 return Max( A0, Max( A1, A2 ) ) * 180. / M_PI;
1451 gp_XYZ p12 = ( P( 1 ) + P( 2 ) ) / 2.;
1452 gp_XYZ p23 = ( P( 2 ) + P( 3 ) ) / 2.;
1453 gp_XYZ p34 = ( P( 3 ) + P( 4 ) ) / 2.;
1454 gp_XYZ p41 = ( P( 4 ) + P( 1 ) ) / 2.;
1456 gp_Vec v1( p34 - p12 ), v2( p23 - p41 );
1457 double A = v1.Magnitude() <= gp::Resolution() || v2.Magnitude() <= gp::Resolution()
1458 ? 0. : fabs( PI2 - v1.Angle( v2 ) );
1460 double val = A * 180. / M_PI;
1462 const double eps = 0.1; // val is in degrees
1464 return val < eps ? 0. : val;
1468 double Skew::GetBadRate( double Value, int /*nbNodes*/ ) const
1470 // the skew is in the range [0.0,PI/2].
1476 SMDSAbs_ElementType Skew::GetType() const
1478 return SMDSAbs_Face;
1482 //================================================================================
1485 Description : Functor for calculating area
1487 //================================================================================
1489 double Area::GetValue( const TSequenceOfXYZ& P )
1494 gp_Vec aVec1( P(2) - P(1) );
1495 gp_Vec aVec2( P(3) - P(1) );
1496 gp_Vec SumVec = aVec1 ^ aVec2;
1498 for (size_t i=4; i<=P.size(); i++)
1500 gp_Vec aVec1( P(i-1) - P(1) );
1501 gp_Vec aVec2( P(i ) - P(1) );
1502 gp_Vec tmp = aVec1 ^ aVec2;
1505 val = SumVec.Magnitude() * 0.5;
1510 double Area::GetBadRate( double Value, int /*nbNodes*/ ) const
1512 // meaningless as it is not a quality control functor
1516 SMDSAbs_ElementType Area::GetType() const
1518 return SMDSAbs_Face;
1521 //================================================================================
1524 Description : Functor for calculating length of edge
1526 //================================================================================
1528 double Length::GetValue( const TSequenceOfXYZ& P )
1530 switch ( P.size() ) {
1531 case 2: return getDistance( P( 1 ), P( 2 ) );
1532 case 3: return getDistance( P( 1 ), P( 2 ) ) + getDistance( P( 2 ), P( 3 ) );
1537 double Length::GetBadRate( double Value, int /*nbNodes*/ ) const
1539 // meaningless as it is not quality control functor
1543 SMDSAbs_ElementType Length::GetType() const
1545 return SMDSAbs_Edge;
1548 //================================================================================
1551 Description : Functor for calculating minimal length of edge
1553 //================================================================================
1555 double Length2D::GetValue( long theElementId )
1559 if ( GetPoints( theElementId, P ))
1563 SMDSAbs_EntityType aType = P.getElementEntity();
1566 case SMDSEntity_Edge:
1568 aVal = getDistance( P( 1 ), P( 2 ) );
1570 case SMDSEntity_Quad_Edge:
1571 if (len == 3) // quadratic edge
1572 aVal = getDistance(P( 1 ),P( 3 )) + getDistance(P( 3 ),P( 2 ));
1574 case SMDSEntity_Triangle:
1575 if (len == 3){ // triangles
1576 double L1 = getDistance(P( 1 ),P( 2 ));
1577 double L2 = getDistance(P( 2 ),P( 3 ));
1578 double L3 = getDistance(P( 3 ),P( 1 ));
1579 aVal = Min(L1,Min(L2,L3));
1582 case SMDSEntity_Quadrangle:
1583 if (len == 4){ // quadrangles
1584 double L1 = getDistance(P( 1 ),P( 2 ));
1585 double L2 = getDistance(P( 2 ),P( 3 ));
1586 double L3 = getDistance(P( 3 ),P( 4 ));
1587 double L4 = getDistance(P( 4 ),P( 1 ));
1588 aVal = Min(Min(L1,L2),Min(L3,L4));
1591 case SMDSEntity_Quad_Triangle:
1592 case SMDSEntity_BiQuad_Triangle:
1593 if (len >= 6){ // quadratic triangles
1594 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1595 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1596 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
1597 aVal = Min(L1,Min(L2,L3));
1600 case SMDSEntity_Quad_Quadrangle:
1601 case SMDSEntity_BiQuad_Quadrangle:
1602 if (len >= 8){ // quadratic quadrangles
1603 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1604 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1605 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
1606 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
1607 aVal = Min(Min(L1,L2),Min(L3,L4));
1610 case SMDSEntity_Tetra:
1611 if (len == 4){ // tetrahedra
1612 double L1 = getDistance(P( 1 ),P( 2 ));
1613 double L2 = getDistance(P( 2 ),P( 3 ));
1614 double L3 = getDistance(P( 3 ),P( 1 ));
1615 double L4 = getDistance(P( 1 ),P( 4 ));
1616 double L5 = getDistance(P( 2 ),P( 4 ));
1617 double L6 = getDistance(P( 3 ),P( 4 ));
1618 aVal = Min(Min(Min(L1,L2),Min(L3,L4)),Min(L5,L6));
1621 case SMDSEntity_Pyramid:
1622 if (len == 5){ // piramids
1623 double L1 = getDistance(P( 1 ),P( 2 ));
1624 double L2 = getDistance(P( 2 ),P( 3 ));
1625 double L3 = getDistance(P( 3 ),P( 4 ));
1626 double L4 = getDistance(P( 4 ),P( 1 ));
1627 double L5 = getDistance(P( 1 ),P( 5 ));
1628 double L6 = getDistance(P( 2 ),P( 5 ));
1629 double L7 = getDistance(P( 3 ),P( 5 ));
1630 double L8 = getDistance(P( 4 ),P( 5 ));
1632 aVal = Min(Min(Min(L1,L2),Min(L3,L4)),Min(L5,L6));
1633 aVal = Min(aVal,Min(L7,L8));
1636 case SMDSEntity_Penta:
1637 if (len == 6) { // pentaidres
1638 double L1 = getDistance(P( 1 ),P( 2 ));
1639 double L2 = getDistance(P( 2 ),P( 3 ));
1640 double L3 = getDistance(P( 3 ),P( 1 ));
1641 double L4 = getDistance(P( 4 ),P( 5 ));
1642 double L5 = getDistance(P( 5 ),P( 6 ));
1643 double L6 = getDistance(P( 6 ),P( 4 ));
1644 double L7 = getDistance(P( 1 ),P( 4 ));
1645 double L8 = getDistance(P( 2 ),P( 5 ));
1646 double L9 = getDistance(P( 3 ),P( 6 ));
1648 aVal = Min(Min(Min(L1,L2),Min(L3,L4)),Min(L5,L6));
1649 aVal = Min(aVal,Min(Min(L7,L8),L9));
1652 case SMDSEntity_Hexa:
1653 if (len == 8){ // hexahedron
1654 double L1 = getDistance(P( 1 ),P( 2 ));
1655 double L2 = getDistance(P( 2 ),P( 3 ));
1656 double L3 = getDistance(P( 3 ),P( 4 ));
1657 double L4 = getDistance(P( 4 ),P( 1 ));
1658 double L5 = getDistance(P( 5 ),P( 6 ));
1659 double L6 = getDistance(P( 6 ),P( 7 ));
1660 double L7 = getDistance(P( 7 ),P( 8 ));
1661 double L8 = getDistance(P( 8 ),P( 5 ));
1662 double L9 = getDistance(P( 1 ),P( 5 ));
1663 double L10= getDistance(P( 2 ),P( 6 ));
1664 double L11= getDistance(P( 3 ),P( 7 ));
1665 double L12= getDistance(P( 4 ),P( 8 ));
1667 aVal = Min(Min(Min(L1,L2),Min(L3,L4)),Min(L5,L6));
1668 aVal = Min(aVal,Min(Min(L7,L8),Min(L9,L10)));
1669 aVal = Min(aVal,Min(L11,L12));
1672 case SMDSEntity_Quad_Tetra:
1673 if (len == 10){ // quadratic tetraidrs
1674 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
1675 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
1676 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
1677 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1678 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
1679 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
1680 aVal = Min(Min(Min(L1,L2),Min(L3,L4)),Min(L5,L6));
1683 case SMDSEntity_Quad_Pyramid:
1684 if (len == 13){ // quadratic piramids
1685 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
1686 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
1687 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1688 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1689 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1690 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
1691 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
1692 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
1693 aVal = Min(Min(Min(L1,L2),Min(L3,L4)),Min(L5,L6));
1694 aVal = Min(aVal,Min(L7,L8));
1697 case SMDSEntity_Quad_Penta:
1698 if (len == 15){ // quadratic pentaidres
1699 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
1700 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
1701 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1702 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1703 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
1704 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
1705 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
1706 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
1707 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
1708 aVal = Min(Min(Min(L1,L2),Min(L3,L4)),Min(L5,L6));
1709 aVal = Min(aVal,Min(Min(L7,L8),L9));
1712 case SMDSEntity_Quad_Hexa:
1713 case SMDSEntity_TriQuad_Hexa:
1714 if (len >= 20) { // quadratic hexaider
1715 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
1716 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
1717 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
1718 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
1719 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
1720 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
1721 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
1722 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
1723 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
1724 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
1725 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
1726 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
1727 aVal = Min(Min(Min(L1,L2),Min(L3,L4)),Min(L5,L6));
1728 aVal = Min(aVal,Min(Min(L7,L8),Min(L9,L10)));
1729 aVal = Min(aVal,Min(L11,L12));
1732 case SMDSEntity_Polygon:
1734 aVal = getDistance( P(1), P( P.size() ));
1735 for ( size_t i = 1; i < P.size(); ++i )
1736 aVal = Min( aVal, getDistance( P( i ), P( i+1 )));
1739 case SMDSEntity_Quad_Polygon:
1741 aVal = getDistance( P(1), P( P.size() )) + getDistance( P(P.size()), P( P.size()-1 ));
1742 for ( size_t i = 1; i < P.size()-1; i += 2 )
1743 aVal = Min( aVal, getDistance( P( i ), P( i+1 )) + getDistance( P( i+1 ), P( i+2 )));
1746 case SMDSEntity_Hexagonal_Prism:
1747 if (len == 12) { // hexagonal prism
1748 double L1 = getDistance(P( 1 ),P( 2 ));
1749 double L2 = getDistance(P( 2 ),P( 3 ));
1750 double L3 = getDistance(P( 3 ),P( 4 ));
1751 double L4 = getDistance(P( 4 ),P( 5 ));
1752 double L5 = getDistance(P( 5 ),P( 6 ));
1753 double L6 = getDistance(P( 6 ),P( 1 ));
1755 double L7 = getDistance(P( 7 ), P( 8 ));
1756 double L8 = getDistance(P( 8 ), P( 9 ));
1757 double L9 = getDistance(P( 9 ), P( 10 ));
1758 double L10= getDistance(P( 10 ),P( 11 ));
1759 double L11= getDistance(P( 11 ),P( 12 ));
1760 double L12= getDistance(P( 12 ),P( 7 ));
1762 double L13 = getDistance(P( 1 ),P( 7 ));
1763 double L14 = getDistance(P( 2 ),P( 8 ));
1764 double L15 = getDistance(P( 3 ),P( 9 ));
1765 double L16 = getDistance(P( 4 ),P( 10 ));
1766 double L17 = getDistance(P( 5 ),P( 11 ));
1767 double L18 = getDistance(P( 6 ),P( 12 ));
1768 aVal = Min(Min(Min(L1,L2),Min(L3,L4)),Min(L5,L6));
1769 aVal = Min(aVal, Min(Min(Min(L7,L8),Min(L9,L10)),Min(L11,L12)));
1770 aVal = Min(aVal, Min(Min(Min(L13,L14),Min(L15,L16)),Min(L17,L18)));
1773 case SMDSEntity_Polyhedra:
1785 if ( myPrecision >= 0 )
1787 double prec = pow( 10., (double)( myPrecision ) );
1788 aVal = floor( aVal * prec + 0.5 ) / prec;
1797 double Length2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1799 // meaningless as it is not a quality control functor
1803 SMDSAbs_ElementType Length2D::GetType() const
1805 return SMDSAbs_Face;
1808 Length2D::Value::Value(double theLength,long thePntId1, long thePntId2):
1811 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1812 if(thePntId1 > thePntId2){
1813 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1817 bool Length2D::Value::operator<(const Length2D::Value& x) const
1819 if(myPntId[0] < x.myPntId[0]) return true;
1820 if(myPntId[0] == x.myPntId[0])
1821 if(myPntId[1] < x.myPntId[1]) return true;
1825 void Length2D::GetValues(TValues& theValues)
1828 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1829 for(; anIter->more(); ){
1830 const SMDS_MeshFace* anElem = anIter->next();
1832 if(anElem->IsQuadratic()) {
1833 const SMDS_VtkFace* F =
1834 dynamic_cast<const SMDS_VtkFace*>(anElem);
1835 // use special nodes iterator
1836 SMDS_ElemIteratorPtr anIter = F->interlacedNodesElemIterator();
1837 long aNodeId[4] = { 0,0,0,0 };
1841 const SMDS_MeshElement* aNode;
1843 aNode = anIter->next();
1844 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1845 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1846 aNodeId[0] = aNodeId[1] = aNode->GetID();
1849 for(; anIter->more(); ){
1850 const SMDS_MeshNode* N1 = static_cast<const SMDS_MeshNode*> (anIter->next());
1851 P[2] = gp_Pnt(N1->X(),N1->Y(),N1->Z());
1852 aNodeId[2] = N1->GetID();
1853 aLength = P[1].Distance(P[2]);
1854 if(!anIter->more()) break;
1855 const SMDS_MeshNode* N2 = static_cast<const SMDS_MeshNode*> (anIter->next());
1856 P[3] = gp_Pnt(N2->X(),N2->Y(),N2->Z());
1857 aNodeId[3] = N2->GetID();
1858 aLength += P[2].Distance(P[3]);
1859 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1860 Value aValue2(aLength,aNodeId[2],aNodeId[3]);
1862 aNodeId[1] = aNodeId[3];
1863 theValues.insert(aValue1);
1864 theValues.insert(aValue2);
1866 aLength += P[2].Distance(P[0]);
1867 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1868 Value aValue2(aLength,aNodeId[2],aNodeId[0]);
1869 theValues.insert(aValue1);
1870 theValues.insert(aValue2);
1873 SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
1874 long aNodeId[2] = {0,0};
1878 const SMDS_MeshElement* aNode;
1879 if(aNodesIter->more()){
1880 aNode = aNodesIter->next();
1881 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1882 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1883 aNodeId[0] = aNodeId[1] = aNode->GetID();
1886 for(; aNodesIter->more(); ){
1887 aNode = aNodesIter->next();
1888 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1889 long anId = aNode->GetID();
1891 P[2] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1893 aLength = P[1].Distance(P[2]);
1895 Value aValue(aLength,aNodeId[1],anId);
1898 theValues.insert(aValue);
1901 aLength = P[0].Distance(P[1]);
1903 Value aValue(aLength,aNodeId[0],aNodeId[1]);
1904 theValues.insert(aValue);
1909 //================================================================================
1911 Class : MultiConnection
1912 Description : Functor for calculating number of faces conneted to the edge
1914 //================================================================================
1916 double MultiConnection::GetValue( const TSequenceOfXYZ& P )
1920 double MultiConnection::GetValue( long theId )
1922 return getNbMultiConnection( myMesh, theId );
1925 double MultiConnection::GetBadRate( double Value, int /*nbNodes*/ ) const
1927 // meaningless as it is not quality control functor
1931 SMDSAbs_ElementType MultiConnection::GetType() const
1933 return SMDSAbs_Edge;
1936 //================================================================================
1938 Class : MultiConnection2D
1939 Description : Functor for calculating number of faces conneted to the edge
1941 //================================================================================
1943 double MultiConnection2D::GetValue( const TSequenceOfXYZ& P )
1948 double MultiConnection2D::GetValue( long theElementId )
1952 const SMDS_MeshElement* aFaceElem = myMesh->FindElement(theElementId);
1953 SMDSAbs_ElementType aType = aFaceElem->GetType();
1958 int i = 0, len = aFaceElem->NbNodes();
1959 SMDS_ElemIteratorPtr anIter = aFaceElem->nodesIterator();
1962 const SMDS_MeshNode *aNode, *aNode0 = 0;
1963 TColStd_MapOfInteger aMap, aMapPrev;
1965 for (i = 0; i <= len; i++) {
1970 if (anIter->more()) {
1971 aNode = (SMDS_MeshNode*)anIter->next();
1979 if (i == 0) aNode0 = aNode;
1981 SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
1982 while (anElemIter->more()) {
1983 const SMDS_MeshElement* anElem = anElemIter->next();
1984 if (anElem != 0 && anElem->GetType() == SMDSAbs_Face) {
1985 int anId = anElem->GetID();
1988 if (aMapPrev.Contains(anId)) {
1993 aResult = Max(aResult, aNb);
2004 double MultiConnection2D::GetBadRate( double Value, int /*nbNodes*/ ) const
2006 // meaningless as it is not quality control functor
2010 SMDSAbs_ElementType MultiConnection2D::GetType() const
2012 return SMDSAbs_Face;
2015 MultiConnection2D::Value::Value(long thePntId1, long thePntId2)
2017 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
2018 if(thePntId1 > thePntId2){
2019 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
2023 bool MultiConnection2D::Value::operator<(const MultiConnection2D::Value& x) const
2025 if(myPntId[0] < x.myPntId[0]) return true;
2026 if(myPntId[0] == x.myPntId[0])
2027 if(myPntId[1] < x.myPntId[1]) return true;
2031 void MultiConnection2D::GetValues(MValues& theValues)
2033 if ( !myMesh ) return;
2034 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2035 for(; anIter->more(); ){
2036 const SMDS_MeshFace* anElem = anIter->next();
2037 SMDS_ElemIteratorPtr aNodesIter;
2038 if ( anElem->IsQuadratic() )
2039 aNodesIter = dynamic_cast<const SMDS_VtkFace*>
2040 (anElem)->interlacedNodesElemIterator();
2042 aNodesIter = anElem->nodesIterator();
2043 long aNodeId[3] = {0,0,0};
2045 //int aNbConnects=0;
2046 const SMDS_MeshNode* aNode0;
2047 const SMDS_MeshNode* aNode1;
2048 const SMDS_MeshNode* aNode2;
2049 if(aNodesIter->more()){
2050 aNode0 = (SMDS_MeshNode*) aNodesIter->next();
2052 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode1;
2053 aNodeId[0] = aNodeId[1] = aNodes->GetID();
2055 for(; aNodesIter->more(); ) {
2056 aNode2 = (SMDS_MeshNode*) aNodesIter->next();
2057 long anId = aNode2->GetID();
2060 Value aValue(aNodeId[1],aNodeId[2]);
2061 MValues::iterator aItr = theValues.find(aValue);
2062 if (aItr != theValues.end()){
2067 theValues[aValue] = 1;
2070 //cout << "NodeIds: "<<aNodeId[1]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
2071 aNodeId[1] = aNodeId[2];
2074 Value aValue(aNodeId[0],aNodeId[2]);
2075 MValues::iterator aItr = theValues.find(aValue);
2076 if (aItr != theValues.end()) {
2081 theValues[aValue] = 1;
2084 //cout << "NodeIds: "<<aNodeId[0]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
2089 //================================================================================
2091 Class : BallDiameter
2092 Description : Functor returning diameter of a ball element
2094 //================================================================================
2096 double BallDiameter::GetValue( long theId )
2098 double diameter = 0;
2100 if ( const SMDS_BallElement* ball =
2101 dynamic_cast<const SMDS_BallElement*>( myMesh->FindElement( theId )))
2103 diameter = ball->GetDiameter();
2108 double BallDiameter::GetBadRate( double Value, int /*nbNodes*/ ) const
2110 // meaningless as it is not a quality control functor
2114 SMDSAbs_ElementType BallDiameter::GetType() const
2116 return SMDSAbs_Ball;
2119 //================================================================================
2121 Class : NodeConnectivityNumber
2122 Description : Functor returning number of elements connected to a node
2124 //================================================================================
2126 double NodeConnectivityNumber::GetValue( long theId )
2130 if ( const SMDS_MeshNode* node = myMesh->FindNode( theId ))
2132 SMDSAbs_ElementType type;
2133 if ( myMesh->NbVolumes() > 0 )
2134 type = SMDSAbs_Volume;
2135 else if ( myMesh->NbFaces() > 0 )
2136 type = SMDSAbs_Face;
2137 else if ( myMesh->NbEdges() > 0 )
2138 type = SMDSAbs_Edge;
2141 nb = node->NbInverseElements( type );
2146 double NodeConnectivityNumber::GetBadRate( double Value, int /*nbNodes*/ ) const
2151 SMDSAbs_ElementType NodeConnectivityNumber::GetType() const
2153 return SMDSAbs_Node;
2160 //================================================================================
2162 Class : BadOrientedVolume
2163 Description : Predicate bad oriented volumes
2165 //================================================================================
2167 BadOrientedVolume::BadOrientedVolume()
2172 void BadOrientedVolume::SetMesh( const SMDS_Mesh* theMesh )
2177 bool BadOrientedVolume::IsSatisfy( long theId )
2182 SMDS_VolumeTool vTool( myMesh->FindElement( theId ));
2183 return !vTool.IsForward();
2186 SMDSAbs_ElementType BadOrientedVolume::GetType() const
2188 return SMDSAbs_Volume;
2192 Class : BareBorderVolume
2195 bool BareBorderVolume::IsSatisfy(long theElementId )
2197 SMDS_VolumeTool myTool;
2198 if ( myTool.Set( myMesh->FindElement(theElementId)))
2200 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2201 if ( myTool.IsFreeFace( iF ))
2203 const SMDS_MeshNode** n = myTool.GetFaceNodes(iF);
2204 std::vector< const SMDS_MeshNode*> nodes( n, n+myTool.NbFaceNodes(iF));
2205 if ( !myMesh->FindElement( nodes, SMDSAbs_Face, /*Nomedium=*/false))
2212 //================================================================================
2214 Class : BareBorderFace
2216 //================================================================================
2218 bool BareBorderFace::IsSatisfy(long theElementId )
2221 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2223 if ( face->GetType() == SMDSAbs_Face )
2225 int nbN = face->NbCornerNodes();
2226 for ( int i = 0; i < nbN && !ok; ++i )
2228 // check if a link is shared by another face
2229 const SMDS_MeshNode* n1 = face->GetNode( i );
2230 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2231 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2232 bool isShared = false;
2233 while ( !isShared && fIt->more() )
2235 const SMDS_MeshElement* f = fIt->next();
2236 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2240 const int iQuad = face->IsQuadratic();
2241 myLinkNodes.resize( 2 + iQuad);
2242 myLinkNodes[0] = n1;
2243 myLinkNodes[1] = n2;
2245 myLinkNodes[2] = face->GetNode( i+nbN );
2246 ok = !myMesh->FindElement( myLinkNodes, SMDSAbs_Edge, /*noMedium=*/false);
2254 //================================================================================
2256 Class : OverConstrainedVolume
2258 //================================================================================
2260 bool OverConstrainedVolume::IsSatisfy(long theElementId )
2262 // An element is over-constrained if it has N-1 free borders where
2263 // N is the number of edges/faces for a 2D/3D element.
2264 SMDS_VolumeTool myTool;
2265 if ( myTool.Set( myMesh->FindElement(theElementId)))
2267 int nbSharedFaces = 0;
2268 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2269 if ( !myTool.IsFreeFace( iF ) && ++nbSharedFaces > 1 )
2271 return ( nbSharedFaces == 1 );
2276 //================================================================================
2278 Class : OverConstrainedFace
2280 //================================================================================
2282 bool OverConstrainedFace::IsSatisfy(long theElementId )
2284 // An element is over-constrained if it has N-1 free borders where
2285 // N is the number of edges/faces for a 2D/3D element.
2286 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2287 if ( face->GetType() == SMDSAbs_Face )
2289 int nbSharedBorders = 0;
2290 int nbN = face->NbCornerNodes();
2291 for ( int i = 0; i < nbN; ++i )
2293 // check if a link is shared by another face
2294 const SMDS_MeshNode* n1 = face->GetNode( i );
2295 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2296 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2297 bool isShared = false;
2298 while ( !isShared && fIt->more() )
2300 const SMDS_MeshElement* f = fIt->next();
2301 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2303 if ( isShared && ++nbSharedBorders > 1 )
2306 return ( nbSharedBorders == 1 );
2311 //================================================================================
2313 Class : CoincidentNodes
2314 Description : Predicate of Coincident nodes
2316 //================================================================================
2318 CoincidentNodes::CoincidentNodes()
2323 bool CoincidentNodes::IsSatisfy( long theElementId )
2325 return myCoincidentIDs.Contains( theElementId );
2328 SMDSAbs_ElementType CoincidentNodes::GetType() const
2330 return SMDSAbs_Node;
2333 void CoincidentNodes::SetMesh( const SMDS_Mesh* theMesh )
2335 myMeshModifTracer.SetMesh( theMesh );
2336 if ( myMeshModifTracer.IsMeshModified() )
2338 TIDSortedNodeSet nodesToCheck;
2339 SMDS_NodeIteratorPtr nIt = theMesh->nodesIterator(/*idInceasingOrder=*/true);
2340 while ( nIt->more() )
2341 nodesToCheck.insert( nodesToCheck.end(), nIt->next() );
2343 std::list< std::list< const SMDS_MeshNode*> > nodeGroups;
2344 SMESH_OctreeNode::FindCoincidentNodes ( nodesToCheck, &nodeGroups, myToler );
2346 myCoincidentIDs.Clear();
2347 std::list< std::list< const SMDS_MeshNode*> >::iterator groupIt = nodeGroups.begin();
2348 for ( ; groupIt != nodeGroups.end(); ++groupIt )
2350 std::list< const SMDS_MeshNode*>& coincNodes = *groupIt;
2351 std::list< const SMDS_MeshNode*>::iterator n = coincNodes.begin();
2352 for ( ; n != coincNodes.end(); ++n )
2353 myCoincidentIDs.Add( (*n)->GetID() );
2358 //================================================================================
2360 Class : CoincidentElements
2361 Description : Predicate of Coincident Elements
2362 Note : This class is suitable only for visualization of Coincident Elements
2364 //================================================================================
2366 CoincidentElements::CoincidentElements()
2371 void CoincidentElements::SetMesh( const SMDS_Mesh* theMesh )
2376 bool CoincidentElements::IsSatisfy( long theElementId )
2378 if ( !myMesh ) return false;
2380 if ( const SMDS_MeshElement* e = myMesh->FindElement( theElementId ))
2382 if ( e->GetType() != GetType() ) return false;
2383 std::set< const SMDS_MeshNode* > elemNodes( e->begin_nodes(), e->end_nodes() );
2384 const int nbNodes = e->NbNodes();
2385 SMDS_ElemIteratorPtr invIt = (*elemNodes.begin())->GetInverseElementIterator( GetType() );
2386 while ( invIt->more() )
2388 const SMDS_MeshElement* e2 = invIt->next();
2389 if ( e2 == e || e2->NbNodes() != nbNodes ) continue;
2391 bool sameNodes = true;
2392 for ( size_t i = 0; i < elemNodes.size() && sameNodes; ++i )
2393 sameNodes = ( elemNodes.count( e2->GetNode( i )));
2401 SMDSAbs_ElementType CoincidentElements1D::GetType() const
2403 return SMDSAbs_Edge;
2405 SMDSAbs_ElementType CoincidentElements2D::GetType() const
2407 return SMDSAbs_Face;
2409 SMDSAbs_ElementType CoincidentElements3D::GetType() const
2411 return SMDSAbs_Volume;
2415 //================================================================================
2418 Description : Predicate for free borders
2420 //================================================================================
2422 FreeBorders::FreeBorders()
2427 void FreeBorders::SetMesh( const SMDS_Mesh* theMesh )
2432 bool FreeBorders::IsSatisfy( long theId )
2434 return getNbMultiConnection( myMesh, theId ) == 1;
2437 SMDSAbs_ElementType FreeBorders::GetType() const
2439 return SMDSAbs_Edge;
2443 //================================================================================
2446 Description : Predicate for free Edges
2448 //================================================================================
2450 FreeEdges::FreeEdges()
2455 void FreeEdges::SetMesh( const SMDS_Mesh* theMesh )
2460 bool FreeEdges::IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId )
2462 TColStd_MapOfInteger aMap;
2463 for ( int i = 0; i < 2; i++ )
2465 SMDS_ElemIteratorPtr anElemIter = theNodes[ i ]->GetInverseElementIterator(SMDSAbs_Face);
2466 while( anElemIter->more() )
2468 if ( const SMDS_MeshElement* anElem = anElemIter->next())
2470 const int anId = anElem->GetID();
2471 if ( anId != theFaceId && !aMap.Add( anId ))
2479 bool FreeEdges::IsSatisfy( long theId )
2484 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2485 if ( aFace == 0 || aFace->GetType() != SMDSAbs_Face || aFace->NbNodes() < 3 )
2488 SMDS_NodeIteratorPtr anIter = aFace->interlacedNodesIterator();
2492 int i = 0, nbNodes = aFace->NbNodes();
2493 std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
2494 while( anIter->more() )
2495 if ( ! ( aNodes[ i++ ] = anIter->next() ))
2497 aNodes[ nbNodes ] = aNodes[ 0 ];
2499 for ( i = 0; i < nbNodes; i++ )
2500 if ( IsFreeEdge( &aNodes[ i ], theId ) )
2506 SMDSAbs_ElementType FreeEdges::GetType() const
2508 return SMDSAbs_Face;
2511 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2514 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
2515 if(thePntId1 > thePntId2){
2516 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
2520 bool FreeEdges::Border::operator<(const FreeEdges::Border& x) const{
2521 if(myPntId[0] < x.myPntId[0]) return true;
2522 if(myPntId[0] == x.myPntId[0])
2523 if(myPntId[1] < x.myPntId[1]) return true;
2527 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2528 FreeEdges::TBorders& theRegistry,
2529 FreeEdges::TBorders& theContainer)
2531 if(theRegistry.find(theBorder) == theRegistry.end()){
2532 theRegistry.insert(theBorder);
2533 theContainer.insert(theBorder);
2535 theContainer.erase(theBorder);
2539 void FreeEdges::GetBoreders(TBorders& theBorders)
2542 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2543 for(; anIter->more(); ){
2544 const SMDS_MeshFace* anElem = anIter->next();
2545 long anElemId = anElem->GetID();
2546 SMDS_ElemIteratorPtr aNodesIter;
2547 if ( anElem->IsQuadratic() )
2548 aNodesIter = static_cast<const SMDS_VtkFace*>(anElem)->
2549 interlacedNodesElemIterator();
2551 aNodesIter = anElem->nodesIterator();
2552 long aNodeId[2] = {0,0};
2553 const SMDS_MeshElement* aNode;
2554 if(aNodesIter->more()){
2555 aNode = aNodesIter->next();
2556 aNodeId[0] = aNodeId[1] = aNode->GetID();
2558 for(; aNodesIter->more(); ){
2559 aNode = aNodesIter->next();
2560 long anId = aNode->GetID();
2561 Border aBorder(anElemId,aNodeId[1],anId);
2563 UpdateBorders(aBorder,aRegistry,theBorders);
2565 Border aBorder(anElemId,aNodeId[0],aNodeId[1]);
2566 UpdateBorders(aBorder,aRegistry,theBorders);
2570 //================================================================================
2573 Description : Predicate for free nodes
2575 //================================================================================
2577 FreeNodes::FreeNodes()
2582 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2587 bool FreeNodes::IsSatisfy( long theNodeId )
2589 const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2593 return (aNode->NbInverseElements() < 1);
2596 SMDSAbs_ElementType FreeNodes::GetType() const
2598 return SMDSAbs_Node;
2602 //================================================================================
2605 Description : Predicate for free faces
2607 //================================================================================
2609 FreeFaces::FreeFaces()
2614 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2619 bool FreeFaces::IsSatisfy( long theId )
2621 if (!myMesh) return false;
2622 // check that faces nodes refers to less than two common volumes
2623 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2624 if ( !aFace || aFace->GetType() != SMDSAbs_Face )
2627 int nbNode = aFace->NbNodes();
2629 // collect volumes to check that number of volumes with count equal nbNode not less than 2
2630 typedef std::map< SMDS_MeshElement*, int > TMapOfVolume; // map of volume counters
2631 typedef std::map< SMDS_MeshElement*, int >::iterator TItrMapOfVolume; // iterator
2632 TMapOfVolume mapOfVol;
2634 SMDS_ElemIteratorPtr nodeItr = aFace->nodesIterator();
2635 while ( nodeItr->more() )
2637 const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>(nodeItr->next());
2638 if ( !aNode ) continue;
2639 SMDS_ElemIteratorPtr volItr = aNode->GetInverseElementIterator(SMDSAbs_Volume);
2640 while ( volItr->more() )
2642 SMDS_MeshElement* aVol = (SMDS_MeshElement*)volItr->next();
2643 TItrMapOfVolume itr = mapOfVol.insert( std::make_pair( aVol, 0 )).first;
2648 TItrMapOfVolume volItr = mapOfVol.begin();
2649 TItrMapOfVolume volEnd = mapOfVol.end();
2650 for ( ; volItr != volEnd; ++volItr )
2651 if ( (*volItr).second >= nbNode )
2653 // face is not free if number of volumes constructed on thier nodes more than one
2657 SMDSAbs_ElementType FreeFaces::GetType() const
2659 return SMDSAbs_Face;
2662 //================================================================================
2664 Class : LinearOrQuadratic
2665 Description : Predicate to verify whether a mesh element is linear
2667 //================================================================================
2669 LinearOrQuadratic::LinearOrQuadratic()
2674 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2679 bool LinearOrQuadratic::IsSatisfy( long theId )
2681 if (!myMesh) return false;
2682 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2683 if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2685 return (!anElem->IsQuadratic());
2688 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2693 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2698 //================================================================================
2701 Description : Functor for check color of group to whic mesh element belongs to
2703 //================================================================================
2705 GroupColor::GroupColor()
2709 bool GroupColor::IsSatisfy( long theId )
2711 return myIDs.count( theId );
2714 void GroupColor::SetType( SMDSAbs_ElementType theType )
2719 SMDSAbs_ElementType GroupColor::GetType() const
2724 static bool isEqual( const Quantity_Color& theColor1,
2725 const Quantity_Color& theColor2 )
2727 // tolerance to compare colors
2728 const double tol = 5*1e-3;
2729 return ( fabs( theColor1.Red() - theColor2.Red() ) < tol &&
2730 fabs( theColor1.Green() - theColor2.Green() ) < tol &&
2731 fabs( theColor1.Blue() - theColor2.Blue() ) < tol );
2734 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2738 const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2742 int nbGrp = aMesh->GetNbGroups();
2746 // iterates on groups and find necessary elements ids
2747 const std::set<SMESHDS_GroupBase*>& aGroups = aMesh->GetGroups();
2748 std::set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
2749 for (; GrIt != aGroups.end(); GrIt++)
2751 SMESHDS_GroupBase* aGrp = (*GrIt);
2754 // check type and color of group
2755 if ( !isEqual( myColor, aGrp->GetColor() ))
2758 // IPAL52867 (prevent infinite recursion via GroupOnFilter)
2759 if ( SMESHDS_GroupOnFilter * gof = dynamic_cast< SMESHDS_GroupOnFilter* >( aGrp ))
2760 if ( gof->GetPredicate().get() == this )
2763 SMDSAbs_ElementType aGrpElType = (SMDSAbs_ElementType)aGrp->GetType();
2764 if ( myType == aGrpElType || (myType == SMDSAbs_All && aGrpElType != SMDSAbs_Node) ) {
2765 // add elements IDS into control
2766 int aSize = aGrp->Extent();
2767 for (int i = 0; i < aSize; i++)
2768 myIDs.insert( aGrp->GetID(i+1) );
2773 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
2775 Kernel_Utils::Localizer loc;
2776 TCollection_AsciiString aStr = theStr;
2777 aStr.RemoveAll( ' ' );
2778 aStr.RemoveAll( '\t' );
2779 for ( int aPos = aStr.Search( ";;" ); aPos != -1; aPos = aStr.Search( ";;" ) )
2780 aStr.Remove( aPos, 2 );
2781 Standard_Real clr[3];
2782 clr[0] = clr[1] = clr[2] = 0.;
2783 for ( int i = 0; i < 3; i++ ) {
2784 TCollection_AsciiString tmpStr = aStr.Token( ";", i+1 );
2785 if ( !tmpStr.IsEmpty() && tmpStr.IsRealValue() )
2786 clr[i] = tmpStr.RealValue();
2788 myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
2791 //=======================================================================
2792 // name : GetRangeStr
2793 // Purpose : Get range as a string.
2794 // Example: "1,2,3,50-60,63,67,70-"
2795 //=======================================================================
2797 void GroupColor::GetColorStr( TCollection_AsciiString& theResStr ) const
2800 theResStr += TCollection_AsciiString( myColor.Red() );
2801 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2802 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2805 //================================================================================
2807 Class : ElemGeomType
2808 Description : Predicate to check element geometry type
2810 //================================================================================
2812 ElemGeomType::ElemGeomType()
2815 myType = SMDSAbs_All;
2816 myGeomType = SMDSGeom_TRIANGLE;
2819 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2824 bool ElemGeomType::IsSatisfy( long theId )
2826 if (!myMesh) return false;
2827 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2830 const SMDSAbs_ElementType anElemType = anElem->GetType();
2831 if ( myType != SMDSAbs_All && anElemType != myType )
2833 bool isOk = ( anElem->GetGeomType() == myGeomType );
2837 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2842 SMDSAbs_ElementType ElemGeomType::GetType() const
2847 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2849 myGeomType = theType;
2852 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2857 //================================================================================
2859 Class : ElemEntityType
2860 Description : Predicate to check element entity type
2862 //================================================================================
2864 ElemEntityType::ElemEntityType():
2866 myType( SMDSAbs_All ),
2867 myEntityType( SMDSEntity_0D )
2871 void ElemEntityType::SetMesh( const SMDS_Mesh* theMesh )
2876 bool ElemEntityType::IsSatisfy( long theId )
2878 if ( !myMesh ) return false;
2879 if ( myType == SMDSAbs_Node )
2880 return myMesh->FindNode( theId );
2881 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2883 myEntityType == anElem->GetEntityType() );
2886 void ElemEntityType::SetType( SMDSAbs_ElementType theType )
2891 SMDSAbs_ElementType ElemEntityType::GetType() const
2896 void ElemEntityType::SetElemEntityType( SMDSAbs_EntityType theEntityType )
2898 myEntityType = theEntityType;
2901 SMDSAbs_EntityType ElemEntityType::GetElemEntityType() const
2903 return myEntityType;
2906 //================================================================================
2908 * \brief Class ConnectedElements
2910 //================================================================================
2912 ConnectedElements::ConnectedElements():
2913 myNodeID(0), myType( SMDSAbs_All ), myOkIDsReady( false ) {}
2915 SMDSAbs_ElementType ConnectedElements::GetType() const
2918 int ConnectedElements::GetNode() const
2919 { return myXYZ.empty() ? myNodeID : 0; } // myNodeID can be found by myXYZ
2921 std::vector<double> ConnectedElements::GetPoint() const
2924 void ConnectedElements::clearOkIDs()
2925 { myOkIDsReady = false; myOkIDs.clear(); }
2927 void ConnectedElements::SetType( SMDSAbs_ElementType theType )
2929 if ( myType != theType || myMeshModifTracer.IsMeshModified() )
2934 void ConnectedElements::SetMesh( const SMDS_Mesh* theMesh )
2936 myMeshModifTracer.SetMesh( theMesh );
2937 if ( myMeshModifTracer.IsMeshModified() )
2940 if ( !myXYZ.empty() )
2941 SetPoint( myXYZ[0], myXYZ[1], myXYZ[2] ); // find a node near myXYZ it in a new mesh
2945 void ConnectedElements::SetNode( int nodeID )
2950 bool isSameDomain = false;
2951 if ( myOkIDsReady && myMeshModifTracer.GetMesh() && !myMeshModifTracer.IsMeshModified() )
2952 if ( const SMDS_MeshNode* n = myMeshModifTracer.GetMesh()->FindNode( myNodeID ))
2954 SMDS_ElemIteratorPtr eIt = n->GetInverseElementIterator( myType );
2955 while ( !isSameDomain && eIt->more() )
2956 isSameDomain = IsSatisfy( eIt->next()->GetID() );
2958 if ( !isSameDomain )
2962 void ConnectedElements::SetPoint( double x, double y, double z )
2970 bool isSameDomain = false;
2972 // find myNodeID by myXYZ if possible
2973 if ( myMeshModifTracer.GetMesh() )
2975 SMESHUtils::Deleter<SMESH_ElementSearcher> searcher
2976 ( SMESH_MeshAlgos::GetElementSearcher( (SMDS_Mesh&) *myMeshModifTracer.GetMesh() ));
2978 std::vector< const SMDS_MeshElement* > foundElems;
2979 searcher->FindElementsByPoint( gp_Pnt(x,y,z), SMDSAbs_All, foundElems );
2981 if ( !foundElems.empty() )
2983 myNodeID = foundElems[0]->GetNode(0)->GetID();
2984 if ( myOkIDsReady && !myMeshModifTracer.IsMeshModified() )
2985 isSameDomain = IsSatisfy( foundElems[0]->GetID() );
2988 if ( !isSameDomain )
2992 bool ConnectedElements::IsSatisfy( long theElementId )
2994 // Here we do NOT check if the mesh has changed, we do it in Set...() only!!!
2996 if ( !myOkIDsReady )
2998 if ( !myMeshModifTracer.GetMesh() )
3000 const SMDS_MeshNode* node0 = myMeshModifTracer.GetMesh()->FindNode( myNodeID );
3004 std::list< const SMDS_MeshNode* > nodeQueue( 1, node0 );
3005 std::set< int > checkedNodeIDs;
3007 // foreach node in nodeQueue:
3008 // foreach element sharing a node:
3009 // add ID of an element of myType to myOkIDs;
3010 // push all element nodes absent from checkedNodeIDs to nodeQueue;
3011 while ( !nodeQueue.empty() )
3013 const SMDS_MeshNode* node = nodeQueue.front();
3014 nodeQueue.pop_front();
3016 // loop on elements sharing the node
3017 SMDS_ElemIteratorPtr eIt = node->GetInverseElementIterator();
3018 while ( eIt->more() )
3020 // keep elements of myType
3021 const SMDS_MeshElement* element = eIt->next();
3022 if ( element->GetType() == myType )
3023 myOkIDs.insert( myOkIDs.end(), element->GetID() );
3025 // enqueue nodes of the element
3026 SMDS_ElemIteratorPtr nIt = element->nodesIterator();
3027 while ( nIt->more() )
3029 const SMDS_MeshNode* n = static_cast< const SMDS_MeshNode* >( nIt->next() );
3030 if ( checkedNodeIDs.insert( n->GetID() ).second )
3031 nodeQueue.push_back( n );
3035 if ( myType == SMDSAbs_Node )
3036 std::swap( myOkIDs, checkedNodeIDs );
3038 size_t totalNbElems = myMeshModifTracer.GetMesh()->GetMeshInfo().NbElements( myType );
3039 if ( myOkIDs.size() == totalNbElems )
3042 myOkIDsReady = true;
3045 return myOkIDs.empty() ? true : myOkIDs.count( theElementId );
3048 //================================================================================
3050 * \brief Class CoplanarFaces
3052 //================================================================================
3056 inline bool isLessAngle( const gp_Vec& v1, const gp_Vec& v2, const double cos )
3058 double dot = v1 * v2; // cos * |v1| * |v2|
3059 double l1 = v1.SquareMagnitude();
3060 double l2 = v2.SquareMagnitude();
3061 return (( dot * cos >= 0 ) &&
3062 ( dot * dot ) / l1 / l2 >= ( cos * cos ));
3065 CoplanarFaces::CoplanarFaces()
3066 : myFaceID(0), myToler(0)
3069 void CoplanarFaces::SetMesh( const SMDS_Mesh* theMesh )
3071 myMeshModifTracer.SetMesh( theMesh );
3072 if ( myMeshModifTracer.IsMeshModified() )
3074 // Build a set of coplanar face ids
3076 myCoplanarIDs.Clear();
3078 if ( !myMeshModifTracer.GetMesh() || !myFaceID || !myToler )
3081 const SMDS_MeshElement* face = myMeshModifTracer.GetMesh()->FindElement( myFaceID );
3082 if ( !face || face->GetType() != SMDSAbs_Face )
3086 gp_Vec myNorm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
3090 const double cosTol = Cos( myToler * M_PI / 180. );
3091 NCollection_Map< SMESH_TLink, SMESH_TLink > checkedLinks;
3093 std::list< std::pair< const SMDS_MeshElement*, gp_Vec > > faceQueue;
3094 faceQueue.push_back( std::make_pair( face, myNorm ));
3095 while ( !faceQueue.empty() )
3097 face = faceQueue.front().first;
3098 myNorm = faceQueue.front().second;
3099 faceQueue.pop_front();
3101 for ( int i = 0, nbN = face->NbCornerNodes(); i < nbN; ++i )
3103 const SMDS_MeshNode* n1 = face->GetNode( i );
3104 const SMDS_MeshNode* n2 = face->GetNode(( i+1 )%nbN);
3105 if ( !checkedLinks.Add( SMESH_TLink( n1, n2 )))
3107 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator(SMDSAbs_Face);
3108 while ( fIt->more() )
3110 const SMDS_MeshElement* f = fIt->next();
3111 if ( f->GetNodeIndex( n2 ) > -1 )
3113 gp_Vec norm = getNormale( static_cast<const SMDS_MeshFace*>(f), &normOK );
3114 if (!normOK || isLessAngle( myNorm, norm, cosTol))
3116 myCoplanarIDs.Add( f->GetID() );
3117 faceQueue.push_back( std::make_pair( f, norm ));
3125 bool CoplanarFaces::IsSatisfy( long theElementId )
3127 return myCoplanarIDs.Contains( theElementId );
3132 *Description : Predicate for Range of Ids.
3133 * Range may be specified with two ways.
3134 * 1. Using AddToRange method
3135 * 2. With SetRangeStr method. Parameter of this method is a string
3136 * like as "1,2,3,50-60,63,67,70-"
3139 //=======================================================================
3140 // name : RangeOfIds
3141 // Purpose : Constructor
3142 //=======================================================================
3143 RangeOfIds::RangeOfIds()
3146 myType = SMDSAbs_All;
3149 //=======================================================================
3151 // Purpose : Set mesh
3152 //=======================================================================
3153 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
3158 //=======================================================================
3159 // name : AddToRange
3160 // Purpose : Add ID to the range
3161 //=======================================================================
3162 bool RangeOfIds::AddToRange( long theEntityId )
3164 myIds.Add( theEntityId );
3168 //=======================================================================
3169 // name : GetRangeStr
3170 // Purpose : Get range as a string.
3171 // Example: "1,2,3,50-60,63,67,70-"
3172 //=======================================================================
3173 void RangeOfIds::GetRangeStr( TCollection_AsciiString& theResStr )
3177 TColStd_SequenceOfInteger anIntSeq;
3178 TColStd_SequenceOfAsciiString aStrSeq;
3180 TColStd_MapIteratorOfMapOfInteger anIter( myIds );
3181 for ( ; anIter.More(); anIter.Next() )
3183 int anId = anIter.Key();
3184 TCollection_AsciiString aStr( anId );
3185 anIntSeq.Append( anId );
3186 aStrSeq.Append( aStr );
3189 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
3191 int aMinId = myMin( i );
3192 int aMaxId = myMax( i );
3194 TCollection_AsciiString aStr;
3195 if ( aMinId != IntegerFirst() )
3200 if ( aMaxId != IntegerLast() )
3203 // find position of the string in result sequence and insert string in it
3204 if ( anIntSeq.Length() == 0 )
3206 anIntSeq.Append( aMinId );
3207 aStrSeq.Append( aStr );
3211 if ( aMinId < anIntSeq.First() )
3213 anIntSeq.Prepend( aMinId );
3214 aStrSeq.Prepend( aStr );
3216 else if ( aMinId > anIntSeq.Last() )
3218 anIntSeq.Append( aMinId );
3219 aStrSeq.Append( aStr );
3222 for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
3223 if ( aMinId < anIntSeq( j ) )
3225 anIntSeq.InsertBefore( j, aMinId );
3226 aStrSeq.InsertBefore( j, aStr );
3232 if ( aStrSeq.Length() == 0 )
3235 theResStr = aStrSeq( 1 );
3236 for ( int j = 2, k = aStrSeq.Length(); j <= k; j++ )
3239 theResStr += aStrSeq( j );
3243 //=======================================================================
3244 // name : SetRangeStr
3245 // Purpose : Define range with string
3246 // Example of entry string: "1,2,3,50-60,63,67,70-"
3247 //=======================================================================
3248 bool RangeOfIds::SetRangeStr( const TCollection_AsciiString& theStr )
3254 TCollection_AsciiString aStr = theStr;
3255 for ( int i = 1; i <= aStr.Length(); ++i )
3257 char c = aStr.Value( i );
3258 if ( !isdigit( c ) && c != ',' && c != '-' )
3259 aStr.SetValue( i, ',');
3261 aStr.RemoveAll( ' ' );
3263 TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
3265 while ( tmpStr != "" )
3267 tmpStr = aStr.Token( ",", i++ );
3268 int aPos = tmpStr.Search( '-' );
3272 if ( tmpStr.IsIntegerValue() )
3273 myIds.Add( tmpStr.IntegerValue() );
3279 TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
3280 TCollection_AsciiString aMinStr = tmpStr;
3282 while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
3283 while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
3285 if ( (!aMinStr.IsEmpty() && !aMinStr.IsIntegerValue()) ||
3286 (!aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue()) )
3289 myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
3290 myMax.Append( aMaxStr.IsEmpty() ? IntegerLast() : aMaxStr.IntegerValue() );
3297 //=======================================================================
3299 // Purpose : Get type of supported entities
3300 //=======================================================================
3301 SMDSAbs_ElementType RangeOfIds::GetType() const
3306 //=======================================================================
3308 // Purpose : Set type of supported entities
3309 //=======================================================================
3310 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
3315 //=======================================================================
3317 // Purpose : Verify whether entity satisfies to this rpedicate
3318 //=======================================================================
3319 bool RangeOfIds::IsSatisfy( long theId )
3324 if ( myType == SMDSAbs_Node )
3326 if ( myMesh->FindNode( theId ) == 0 )
3331 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
3332 if ( anElem == 0 || (myType != anElem->GetType() && myType != SMDSAbs_All ))
3336 if ( myIds.Contains( theId ) )
3339 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
3340 if ( theId >= myMin( i ) && theId <= myMax( i ) )
3348 Description : Base class for comparators
3350 Comparator::Comparator():
3354 Comparator::~Comparator()
3357 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
3360 myFunctor->SetMesh( theMesh );
3363 void Comparator::SetMargin( double theValue )
3365 myMargin = theValue;
3368 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
3370 myFunctor = theFunct;
3373 SMDSAbs_ElementType Comparator::GetType() const
3375 return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
3378 double Comparator::GetMargin()
3386 Description : Comparator "<"
3388 bool LessThan::IsSatisfy( long theId )
3390 return myFunctor && myFunctor->GetValue( theId ) < myMargin;
3396 Description : Comparator ">"
3398 bool MoreThan::IsSatisfy( long theId )
3400 return myFunctor && myFunctor->GetValue( theId ) > myMargin;
3406 Description : Comparator "="
3409 myToler(Precision::Confusion())
3412 bool EqualTo::IsSatisfy( long theId )
3414 return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
3417 void EqualTo::SetTolerance( double theToler )
3422 double EqualTo::GetTolerance()
3429 Description : Logical NOT predicate
3431 LogicalNOT::LogicalNOT()
3434 LogicalNOT::~LogicalNOT()
3437 bool LogicalNOT::IsSatisfy( long theId )
3439 return myPredicate && !myPredicate->IsSatisfy( theId );
3442 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
3445 myPredicate->SetMesh( theMesh );
3448 void LogicalNOT::SetPredicate( PredicatePtr thePred )
3450 myPredicate = thePred;
3453 SMDSAbs_ElementType LogicalNOT::GetType() const
3455 return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
3460 Class : LogicalBinary
3461 Description : Base class for binary logical predicate
3463 LogicalBinary::LogicalBinary()
3466 LogicalBinary::~LogicalBinary()
3469 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
3472 myPredicate1->SetMesh( theMesh );
3475 myPredicate2->SetMesh( theMesh );
3478 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
3480 myPredicate1 = thePredicate;
3483 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
3485 myPredicate2 = thePredicate;
3488 SMDSAbs_ElementType LogicalBinary::GetType() const
3490 if ( !myPredicate1 || !myPredicate2 )
3493 SMDSAbs_ElementType aType1 = myPredicate1->GetType();
3494 SMDSAbs_ElementType aType2 = myPredicate2->GetType();
3496 return aType1 == aType2 ? aType1 : SMDSAbs_All;
3502 Description : Logical AND
3504 bool LogicalAND::IsSatisfy( long theId )
3509 myPredicate1->IsSatisfy( theId ) &&
3510 myPredicate2->IsSatisfy( theId );
3516 Description : Logical OR
3518 bool LogicalOR::IsSatisfy( long theId )
3523 (myPredicate1->IsSatisfy( theId ) ||
3524 myPredicate2->IsSatisfy( theId ));
3533 // #include <tbb/parallel_for.h>
3534 // #include <tbb/enumerable_thread_specific.h>
3536 // namespace Parallel
3538 // typedef tbb::enumerable_thread_specific< TIdSequence > TIdSeq;
3542 // const SMDS_Mesh* myMesh;
3543 // PredicatePtr myPredicate;
3544 // TIdSeq & myOKIds;
3545 // Predicate( const SMDS_Mesh* m, PredicatePtr p, TIdSeq & ids ):
3546 // myMesh(m), myPredicate(p->Duplicate()), myOKIds(ids) {}
3547 // void operator() ( const tbb::blocked_range<size_t>& r ) const
3549 // for ( size_t i = r.begin(); i != r.end(); ++i )
3550 // if ( myPredicate->IsSatisfy( i ))
3551 // myOKIds.local().push_back();
3563 void Filter::SetPredicate( PredicatePtr thePredicate )
3565 myPredicate = thePredicate;
3568 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3569 PredicatePtr thePredicate,
3570 TIdSequence& theSequence )
3572 theSequence.clear();
3574 if ( !theMesh || !thePredicate )
3577 thePredicate->SetMesh( theMesh );
3579 SMDS_ElemIteratorPtr elemIt = theMesh->elementsIterator( thePredicate->GetType() );
3581 while ( elemIt->more() ) {
3582 const SMDS_MeshElement* anElem = elemIt->next();
3583 long anId = anElem->GetID();
3584 if ( thePredicate->IsSatisfy( anId ) )
3585 theSequence.push_back( anId );
3590 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3591 Filter::TIdSequence& theSequence )
3593 GetElementsId(theMesh,myPredicate,theSequence);
3600 typedef std::set<SMDS_MeshFace*> TMapOfFacePtr;
3606 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
3607 SMDS_MeshNode* theNode2 )
3613 ManifoldPart::Link::~Link()
3619 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
3621 if ( myNode1 == theLink.myNode1 &&
3622 myNode2 == theLink.myNode2 )
3624 else if ( myNode1 == theLink.myNode2 &&
3625 myNode2 == theLink.myNode1 )
3631 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
3633 if(myNode1 < x.myNode1) return true;
3634 if(myNode1 == x.myNode1)
3635 if(myNode2 < x.myNode2) return true;
3639 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
3640 const ManifoldPart::Link& theLink2 )
3642 return theLink1.IsEqual( theLink2 );
3645 ManifoldPart::ManifoldPart()
3648 myAngToler = Precision::Angular();
3649 myIsOnlyManifold = true;
3652 ManifoldPart::~ManifoldPart()
3657 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
3663 SMDSAbs_ElementType ManifoldPart::GetType() const
3664 { return SMDSAbs_Face; }
3666 bool ManifoldPart::IsSatisfy( long theElementId )
3668 return myMapIds.Contains( theElementId );
3671 void ManifoldPart::SetAngleTolerance( const double theAngToler )
3672 { myAngToler = theAngToler; }
3674 double ManifoldPart::GetAngleTolerance() const
3675 { return myAngToler; }
3677 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
3678 { myIsOnlyManifold = theIsOnly; }
3680 void ManifoldPart::SetStartElem( const long theStartId )
3681 { myStartElemId = theStartId; }
3683 bool ManifoldPart::process()
3686 myMapBadGeomIds.Clear();
3688 myAllFacePtr.clear();
3689 myAllFacePtrIntDMap.clear();
3693 // collect all faces into own map
3694 SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3695 for (; anFaceItr->more(); )
3697 SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3698 myAllFacePtr.push_back( aFacePtr );
3699 myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3702 SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3706 // the map of non manifold links and bad geometry
3707 TMapOfLink aMapOfNonManifold;
3708 TColStd_MapOfInteger aMapOfTreated;
3710 // begin cycle on faces from start index and run on vector till the end
3711 // and from begin to start index to cover whole vector
3712 const int aStartIndx = myAllFacePtrIntDMap[aStartFace];
3713 bool isStartTreat = false;
3714 for ( int fi = aStartIndx; !isStartTreat || fi != aStartIndx ; fi++ )
3716 if ( fi == aStartIndx )
3717 isStartTreat = true;
3718 // as result next time when fi will be equal to aStartIndx
3720 SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3721 if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3724 aMapOfTreated.Add( aFacePtr->GetID() );
3725 TColStd_MapOfInteger aResFaces;
3726 if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3727 aMapOfNonManifold, aResFaces ) )
3729 TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3730 for ( ; anItr.More(); anItr.Next() )
3732 int aFaceId = anItr.Key();
3733 aMapOfTreated.Add( aFaceId );
3734 myMapIds.Add( aFaceId );
3737 if ( fi == int( myAllFacePtr.size() - 1 ))
3739 } // end run on vector of faces
3740 return !myMapIds.IsEmpty();
3743 static void getLinks( const SMDS_MeshFace* theFace,
3744 ManifoldPart::TVectorOfLink& theLinks )
3746 int aNbNode = theFace->NbNodes();
3747 SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3749 SMDS_MeshNode* aNode = 0;
3750 for ( ; aNodeItr->more() && i <= aNbNode; )
3753 SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3757 SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3759 ManifoldPart::Link aLink( aN1, aN2 );
3760 theLinks.push_back( aLink );
3764 bool ManifoldPart::findConnected
3765 ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3766 SMDS_MeshFace* theStartFace,
3767 ManifoldPart::TMapOfLink& theNonManifold,
3768 TColStd_MapOfInteger& theResFaces )
3770 theResFaces.Clear();
3771 if ( !theAllFacePtrInt.size() )
3774 if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3776 myMapBadGeomIds.Add( theStartFace->GetID() );
3780 ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3781 ManifoldPart::TVectorOfLink aSeqOfBoundary;
3782 theResFaces.Add( theStartFace->GetID() );
3783 ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3785 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3786 aDMapLinkFace, theNonManifold, theStartFace );
3788 bool isDone = false;
3789 while ( !isDone && aMapOfBoundary.size() != 0 )
3791 bool isToReset = false;
3792 ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3793 for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3795 ManifoldPart::Link aLink = *pLink;
3796 if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3798 // each link could be treated only once
3799 aMapToSkip.insert( aLink );
3801 ManifoldPart::TVectorOfFacePtr aFaces;
3803 if ( myIsOnlyManifold &&
3804 (theNonManifold.find( aLink ) != theNonManifold.end()) )
3808 getFacesByLink( aLink, aFaces );
3809 // filter the element to keep only indicated elements
3810 ManifoldPart::TVectorOfFacePtr aFiltered;
3811 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3812 for ( ; pFace != aFaces.end(); ++pFace )
3814 SMDS_MeshFace* aFace = *pFace;
3815 if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3816 aFiltered.push_back( aFace );
3819 if ( aFaces.size() < 2 ) // no neihgbour faces
3821 else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3823 theNonManifold.insert( aLink );
3828 // compare normal with normals of neighbor element
3829 SMDS_MeshFace* aPrevFace = aDMapLinkFace[ aLink ];
3830 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3831 for ( ; pFace != aFaces.end(); ++pFace )
3833 SMDS_MeshFace* aNextFace = *pFace;
3834 if ( aPrevFace == aNextFace )
3836 int anNextFaceID = aNextFace->GetID();
3837 if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3838 // should not be with non manifold restriction. probably bad topology
3840 // check if face was treated and skipped
3841 if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3842 !isInPlane( aPrevFace, aNextFace ) )
3844 // add new element to connected and extend the boundaries.
3845 theResFaces.Add( anNextFaceID );
3846 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3847 aDMapLinkFace, theNonManifold, aNextFace );
3851 isDone = !isToReset;
3854 return !theResFaces.IsEmpty();
3857 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3858 const SMDS_MeshFace* theFace2 )
3860 gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3861 gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3862 if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3864 myMapBadGeomIds.Add( theFace2->GetID() );
3867 if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
3873 void ManifoldPart::expandBoundary
3874 ( ManifoldPart::TMapOfLink& theMapOfBoundary,
3875 ManifoldPart::TVectorOfLink& theSeqOfBoundary,
3876 ManifoldPart::TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
3877 ManifoldPart::TMapOfLink& theNonManifold,
3878 SMDS_MeshFace* theNextFace ) const
3880 ManifoldPart::TVectorOfLink aLinks;
3881 getLinks( theNextFace, aLinks );
3882 int aNbLink = (int)aLinks.size();
3883 for ( int i = 0; i < aNbLink; i++ )
3885 ManifoldPart::Link aLink = aLinks[ i ];
3886 if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3888 if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3890 if ( myIsOnlyManifold )
3892 // remove from boundary
3893 theMapOfBoundary.erase( aLink );
3894 ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3895 for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3897 ManifoldPart::Link aBoundLink = *pLink;
3898 if ( aBoundLink.IsEqual( aLink ) )
3900 theSeqOfBoundary.erase( pLink );
3908 theMapOfBoundary.insert( aLink );
3909 theSeqOfBoundary.push_back( aLink );
3910 theDMapLinkFacePtr[ aLink ] = theNextFace;
3915 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3916 ManifoldPart::TVectorOfFacePtr& theFaces ) const
3918 std::set<SMDS_MeshCell *> aSetOfFaces;
3919 // take all faces that shared first node
3920 SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3921 for ( ; anItr->more(); )
3923 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3926 aSetOfFaces.insert( aFace );
3928 // take all faces that shared second node
3929 anItr = theLink.myNode2->facesIterator();
3930 // find the common part of two sets
3931 for ( ; anItr->more(); )
3933 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3934 if ( aSetOfFaces.count( aFace ) )
3935 theFaces.push_back( aFace );
3940 Class : BelongToMeshGroup
3941 Description : Verify whether a mesh element is included into a mesh group
3943 BelongToMeshGroup::BelongToMeshGroup(): myGroup( 0 )
3947 void BelongToMeshGroup::SetGroup( SMESHDS_GroupBase* g )
3952 void BelongToMeshGroup::SetStoreName( const std::string& sn )
3957 void BelongToMeshGroup::SetMesh( const SMDS_Mesh* theMesh )
3959 if ( myGroup && myGroup->GetMesh() != theMesh )
3963 if ( !myGroup && !myStoreName.empty() )
3965 if ( const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh))
3967 const std::set<SMESHDS_GroupBase*>& grps = aMesh->GetGroups();
3968 std::set<SMESHDS_GroupBase*>::const_iterator g = grps.begin();
3969 for ( ; g != grps.end() && !myGroup; ++g )
3970 if ( *g && myStoreName == (*g)->GetStoreName() )
3976 myGroup->IsEmpty(); // make GroupOnFilter update its predicate
3980 bool BelongToMeshGroup::IsSatisfy( long theElementId )
3982 return myGroup ? myGroup->Contains( theElementId ) : false;
3985 SMDSAbs_ElementType BelongToMeshGroup::GetType() const
3987 return myGroup ? myGroup->GetType() : SMDSAbs_All;
3990 //================================================================================
3991 // ElementsOnSurface
3992 //================================================================================
3994 ElementsOnSurface::ElementsOnSurface()
3997 myType = SMDSAbs_All;
3999 myToler = Precision::Confusion();
4000 myUseBoundaries = false;
4003 ElementsOnSurface::~ElementsOnSurface()
4007 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
4009 myMeshModifTracer.SetMesh( theMesh );
4010 if ( myMeshModifTracer.IsMeshModified())
4014 bool ElementsOnSurface::IsSatisfy( long theElementId )
4016 return myIds.Contains( theElementId );
4019 SMDSAbs_ElementType ElementsOnSurface::GetType() const
4022 void ElementsOnSurface::SetTolerance( const double theToler )
4024 if ( myToler != theToler )
4029 double ElementsOnSurface::GetTolerance() const
4032 void ElementsOnSurface::SetUseBoundaries( bool theUse )
4034 if ( myUseBoundaries != theUse ) {
4035 myUseBoundaries = theUse;
4036 SetSurface( mySurf, myType );
4040 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
4041 const SMDSAbs_ElementType theType )
4046 if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
4048 mySurf = TopoDS::Face( theShape );
4049 BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
4051 u1 = SA.FirstUParameter(),
4052 u2 = SA.LastUParameter(),
4053 v1 = SA.FirstVParameter(),
4054 v2 = SA.LastVParameter();
4055 Handle(Geom_Surface) surf = BRep_Tool::Surface( mySurf );
4056 myProjector.Init( surf, u1,u2, v1,v2 );
4060 void ElementsOnSurface::process()
4063 if ( mySurf.IsNull() )
4066 if ( !myMeshModifTracer.GetMesh() )
4069 myIds.ReSize( myMeshModifTracer.GetMesh()->GetMeshInfo().NbElements( myType ));
4071 SMDS_ElemIteratorPtr anIter = myMeshModifTracer.GetMesh()->elementsIterator( myType );
4072 for(; anIter->more(); )
4073 process( anIter->next() );
4076 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
4078 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
4079 bool isSatisfy = true;
4080 for ( ; aNodeItr->more(); )
4082 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
4083 if ( !isOnSurface( aNode ) )
4090 myIds.Add( theElemPtr->GetID() );
4093 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
4095 if ( mySurf.IsNull() )
4098 gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
4099 // double aToler2 = myToler * myToler;
4100 // if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
4102 // gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
4103 // if ( aPln.SquareDistance( aPnt ) > aToler2 )
4106 // else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
4108 // gp_Cylinder aCyl = Handle(Geom_CylindricalSurface)::DownCast(mySurf)->Cylinder();
4109 // double aRad = aCyl.Radius();
4110 // gp_Ax3 anAxis = aCyl.Position();
4111 // gp_XYZ aLoc = aCyl.Location().XYZ();
4112 // double aXDist = anAxis.XDirection().XYZ() * ( aPnt.XYZ() - aLoc );
4113 // double aYDist = anAxis.YDirection().XYZ() * ( aPnt.XYZ() - aLoc );
4114 // if ( fabs(aXDist*aXDist + aYDist*aYDist - aRad*aRad) > aToler2 )
4119 myProjector.Perform( aPnt );
4120 bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
4126 //================================================================================
4128 //================================================================================
4131 const int theIsCheckedFlag = 0x0000100;
4134 struct ElementsOnShape::Classifier
4136 Classifier() { mySolidClfr = 0; myFlags = 0; }
4138 void Init(const TopoDS_Shape& s, double tol, const Bnd_B3d* box = 0 );
4139 bool IsOut(const gp_Pnt& p) { return SetChecked( true ), (this->*myIsOutFun)( p ); }
4140 TopAbs_ShapeEnum ShapeType() const { return myShape.ShapeType(); }
4141 const TopoDS_Shape& Shape() const { return myShape; }
4142 const Bnd_B3d* GetBndBox() const { return & myBox; }
4143 bool IsChecked() { return myFlags & theIsCheckedFlag; }
4144 bool IsSetFlag( int flag ) const { return myFlags & flag; }
4145 void SetChecked( bool is ) { is ? SetFlag( theIsCheckedFlag ) : UnsetFlag( theIsCheckedFlag ); }
4146 void SetFlag ( int flag ) { myFlags |= flag; }
4147 void UnsetFlag( int flag ) { myFlags &= ~flag; }
4150 bool isOutOfSolid (const gp_Pnt& p);
4151 bool isOutOfBox (const gp_Pnt& p);
4152 bool isOutOfFace (const gp_Pnt& p);
4153 bool isOutOfEdge (const gp_Pnt& p);
4154 bool isOutOfVertex(const gp_Pnt& p);
4155 bool isBox (const TopoDS_Shape& s);
4157 bool (Classifier::* myIsOutFun)(const gp_Pnt& p);
4158 BRepClass3d_SolidClassifier* mySolidClfr; // ptr because of a run-time forbidden copy-constructor
4160 GeomAPI_ProjectPointOnSurf myProjFace;
4161 GeomAPI_ProjectPointOnCurve myProjEdge;
4163 TopoDS_Shape myShape;
4168 struct ElementsOnShape::OctreeClassifier : public SMESH_Octree
4170 OctreeClassifier( const std::vector< ElementsOnShape::Classifier* >& classifiers );
4171 OctreeClassifier( const OctreeClassifier* otherTree,
4172 const std::vector< ElementsOnShape::Classifier >& clsOther,
4173 std::vector< ElementsOnShape::Classifier >& cls );
4174 void GetClassifiersAtPoint( const gp_XYZ& p,
4175 std::vector< ElementsOnShape::Classifier* >& classifiers );
4177 OctreeClassifier() {}
4178 SMESH_Octree* newChild() const { return new OctreeClassifier; }
4179 void buildChildrenData();
4180 Bnd_B3d* buildRootBox();
4182 std::vector< ElementsOnShape::Classifier* > myClassifiers;
4186 ElementsOnShape::ElementsOnShape():
4188 myType(SMDSAbs_All),
4189 myToler(Precision::Confusion()),
4190 myAllNodesFlag(false)
4194 ElementsOnShape::~ElementsOnShape()
4199 Predicate* ElementsOnShape::clone() const
4201 ElementsOnShape* cln = new ElementsOnShape();
4202 cln->SetAllNodes ( myAllNodesFlag );
4203 cln->SetTolerance( myToler );
4204 cln->SetMesh ( myMeshModifTracer.GetMesh() );
4205 cln->myShape = myShape; // avoid creation of myClassifiers
4206 cln->SetShape ( myShape, myType );
4207 cln->myClassifiers.resize( myClassifiers.size() );
4208 for ( size_t i = 0; i < myClassifiers.size(); ++i )
4209 cln->myClassifiers[ i ].Init( BRepBuilderAPI_Copy( myClassifiers[ i ].Shape()),
4210 myToler, myClassifiers[ i ].GetBndBox() );
4211 if ( myOctree ) // copy myOctree
4213 cln->myOctree = new OctreeClassifier( myOctree, myClassifiers, cln->myClassifiers );
4218 SMDSAbs_ElementType ElementsOnShape::GetType() const
4223 void ElementsOnShape::SetTolerance (const double theToler)
4225 if (myToler != theToler) {
4227 SetShape(myShape, myType);
4231 double ElementsOnShape::GetTolerance() const
4236 void ElementsOnShape::SetAllNodes (bool theAllNodes)
4238 myAllNodesFlag = theAllNodes;
4241 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
4243 myMeshModifTracer.SetMesh( theMesh );
4244 if ( myMeshModifTracer.IsMeshModified())
4246 size_t nbNodes = theMesh ? theMesh->NbNodes() : 0;
4247 if ( myNodeIsChecked.size() == nbNodes )
4249 std::fill( myNodeIsChecked.begin(), myNodeIsChecked.end(), false );
4253 SMESHUtils::FreeVector( myNodeIsChecked );
4254 SMESHUtils::FreeVector( myNodeIsOut );
4255 myNodeIsChecked.resize( nbNodes, false );
4256 myNodeIsOut.resize( nbNodes );
4261 bool ElementsOnShape::getNodeIsOut( const SMDS_MeshNode* n, bool& isOut )
4263 if ( n->GetID() >= (int) myNodeIsChecked.size() ||
4264 !myNodeIsChecked[ n->GetID() ])
4267 isOut = myNodeIsOut[ n->GetID() ];
4271 void ElementsOnShape::setNodeIsOut( const SMDS_MeshNode* n, bool isOut )
4273 if ( n->GetID() < (int) myNodeIsChecked.size() )
4275 myNodeIsChecked[ n->GetID() ] = true;
4276 myNodeIsOut [ n->GetID() ] = isOut;
4280 void ElementsOnShape::SetShape (const TopoDS_Shape& theShape,
4281 const SMDSAbs_ElementType theType)
4283 bool shapeChanges = ( myShape != theShape );
4286 if ( myShape.IsNull() ) return;
4290 TopTools_IndexedMapOfShape shapesMap;
4291 TopAbs_ShapeEnum shapeTypes[4] = { TopAbs_SOLID, TopAbs_FACE, TopAbs_EDGE, TopAbs_VERTEX };
4292 TopExp_Explorer sub;
4293 for ( int i = 0; i < 4; ++i )
4295 if ( shapesMap.IsEmpty() )
4296 for ( sub.Init( myShape, shapeTypes[i] ); sub.More(); sub.Next() )
4297 shapesMap.Add( sub.Current() );
4299 for ( sub.Init( myShape, shapeTypes[i], shapeTypes[i-1] ); sub.More(); sub.Next() )
4300 shapesMap.Add( sub.Current() );
4304 myClassifiers.resize( shapesMap.Extent() );
4305 for ( int i = 0; i < shapesMap.Extent(); ++i )
4306 myClassifiers[ i ].Init( shapesMap( i+1 ), myToler );
4309 if ( theType == SMDSAbs_Node )
4311 SMESHUtils::FreeVector( myNodeIsChecked );
4312 SMESHUtils::FreeVector( myNodeIsOut );
4316 std::fill( myNodeIsChecked.begin(), myNodeIsChecked.end(), false );
4320 void ElementsOnShape::clearClassifiers()
4322 // for ( size_t i = 0; i < myClassifiers.size(); ++i )
4323 // delete myClassifiers[ i ];
4324 myClassifiers.clear();
4330 bool ElementsOnShape::IsSatisfy( long elemId )
4332 const SMDS_Mesh* mesh = myMeshModifTracer.GetMesh();
4333 const SMDS_MeshElement* elem =
4334 ( myType == SMDSAbs_Node ? mesh->FindNode( elemId ) : mesh->FindElement( elemId ));
4335 if ( !elem || myClassifiers.empty() )
4338 bool isSatisfy = myAllNodesFlag, isNodeOut;
4340 gp_XYZ centerXYZ (0, 0, 0);
4342 if ( !myOctree && myClassifiers.size() > 5 )
4344 myWorkClassifiers.resize( myClassifiers.size() );
4345 for ( size_t i = 0; i < myClassifiers.size(); ++i )
4346 myWorkClassifiers[ i ] = & myClassifiers[ i ];
4347 myOctree = new OctreeClassifier( myWorkClassifiers );
4350 SMDS_ElemIteratorPtr aNodeItr = elem->nodesIterator();
4351 while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
4353 SMESH_TNodeXYZ aPnt( aNodeItr->next() );
4357 if ( !getNodeIsOut( aPnt._node, isNodeOut ))
4361 myWorkClassifiers.clear();
4362 myOctree->GetClassifiersAtPoint( aPnt, myWorkClassifiers );
4364 for ( size_t i = 0; i < myWorkClassifiers.size(); ++i )
4365 myWorkClassifiers[i]->SetChecked( false );
4367 for ( size_t i = 0; i < myWorkClassifiers.size() && isNodeOut; ++i )
4368 if ( !myWorkClassifiers[i]->IsChecked() )
4369 isNodeOut = myWorkClassifiers[i]->IsOut( aPnt );
4373 for ( size_t i = 0; i < myClassifiers.size() && isNodeOut; ++i )
4374 isNodeOut = myClassifiers[i].IsOut( aPnt );
4376 setNodeIsOut( aPnt._node, isNodeOut );
4378 isSatisfy = !isNodeOut;
4381 // Check the center point for volumes MantisBug 0020168
4384 myClassifiers[0].ShapeType() == TopAbs_SOLID )
4386 centerXYZ /= elem->NbNodes();
4389 for ( size_t i = 0; i < myWorkClassifiers.size() && !isSatisfy; ++i )
4390 isSatisfy = ! myWorkClassifiers[i]->IsOut( centerXYZ );
4392 for ( size_t i = 0; i < myClassifiers.size() && !isSatisfy; ++i )
4393 isSatisfy = ! myClassifiers[i].IsOut( centerXYZ );
4399 void ElementsOnShape::Classifier::Init( const TopoDS_Shape& theShape,
4401 const Bnd_B3d* theBox )
4407 bool isShapeBox = false;
4408 switch ( myShape.ShapeType() )
4412 if (( isShapeBox = isBox( theShape )))
4414 myIsOutFun = & ElementsOnShape::Classifier::isOutOfBox;
4418 mySolidClfr = new BRepClass3d_SolidClassifier(theShape);
4419 myIsOutFun = & ElementsOnShape::Classifier::isOutOfSolid;
4425 Standard_Real u1,u2,v1,v2;
4426 Handle(Geom_Surface) surf = BRep_Tool::Surface( TopoDS::Face( theShape ));
4427 surf->Bounds( u1,u2,v1,v2 );
4428 myProjFace.Init(surf, u1,u2, v1,v2, myTol );
4429 myIsOutFun = & ElementsOnShape::Classifier::isOutOfFace;
4434 Standard_Real u1, u2;
4435 Handle(Geom_Curve) curve = BRep_Tool::Curve( TopoDS::Edge( theShape ), u1, u2);
4436 myProjEdge.Init(curve, u1, u2);
4437 myIsOutFun = & ElementsOnShape::Classifier::isOutOfEdge;
4442 myVertexXYZ = BRep_Tool::Pnt( TopoDS::Vertex( theShape ) );
4443 myIsOutFun = & ElementsOnShape::Classifier::isOutOfVertex;
4447 throw SALOME_Exception("Programmer error in usage of ElementsOnShape::Classifier");
4459 BRepBndLib::Add( myShape, box );
4461 myBox.Add( box.CornerMin() );
4462 myBox.Add( box.CornerMax() );
4463 gp_XYZ halfSize = 0.5 * ( box.CornerMax().XYZ() - box.CornerMin().XYZ() );
4464 for ( int iDim = 1; iDim <= 3; ++iDim )
4466 double x = halfSize.Coord( iDim );
4467 halfSize.SetCoord( iDim, x + Max( myTol, 1e-2 * x ));
4469 myBox.SetHSize( halfSize );
4474 ElementsOnShape::Classifier::~Classifier()
4476 delete mySolidClfr; mySolidClfr = 0;
4479 bool ElementsOnShape::Classifier::isOutOfSolid (const gp_Pnt& p)
4481 mySolidClfr->Perform( p, myTol );
4482 return ( mySolidClfr->State() != TopAbs_IN && mySolidClfr->State() != TopAbs_ON );
4485 bool ElementsOnShape::Classifier::isOutOfBox (const gp_Pnt& p)
4487 return myBox.IsOut( p.XYZ() );
4490 bool ElementsOnShape::Classifier::isOutOfFace (const gp_Pnt& p)
4492 myProjFace.Perform( p );
4493 if ( myProjFace.IsDone() && myProjFace.LowerDistance() <= myTol )
4495 // check relatively to the face
4496 Quantity_Parameter u, v;
4497 myProjFace.LowerDistanceParameters(u, v);
4498 gp_Pnt2d aProjPnt (u, v);
4499 BRepClass_FaceClassifier aClsf ( TopoDS::Face( myShape ), aProjPnt, myTol );
4500 if ( aClsf.State() == TopAbs_IN || aClsf.State() == TopAbs_ON )
4506 bool ElementsOnShape::Classifier::isOutOfEdge (const gp_Pnt& p)
4508 myProjEdge.Perform( p );
4509 return ! ( myProjEdge.NbPoints() > 0 && myProjEdge.LowerDistance() <= myTol );
4512 bool ElementsOnShape::Classifier::isOutOfVertex(const gp_Pnt& p)
4514 return ( myVertexXYZ.Distance( p ) > myTol );
4517 bool ElementsOnShape::Classifier::isBox (const TopoDS_Shape& theShape)
4519 TopTools_IndexedMapOfShape vMap;
4520 TopExp::MapShapes( theShape, TopAbs_VERTEX, vMap );
4521 if ( vMap.Extent() != 8 )
4525 for ( int i = 1; i <= 8; ++i )
4526 myBox.Add( BRep_Tool::Pnt( TopoDS::Vertex( vMap( i ))).XYZ() );
4528 gp_XYZ pMin = myBox.CornerMin(), pMax = myBox.CornerMax();
4529 for ( int i = 1; i <= 8; ++i )
4531 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( vMap( i )));
4532 for ( int iC = 1; iC <= 3; ++ iC )
4534 double d1 = Abs( pMin.Coord( iC ) - p.Coord( iC ));
4535 double d2 = Abs( pMax.Coord( iC ) - p.Coord( iC ));
4536 if ( Min( d1, d2 ) > myTol )
4540 myBox.Enlarge( myTol );
4545 OctreeClassifier::OctreeClassifier( const std::vector< ElementsOnShape::Classifier* >& classifiers )
4546 :SMESH_Octree( new SMESH_TreeLimit )
4548 myClassifiers = classifiers;
4553 OctreeClassifier::OctreeClassifier( const OctreeClassifier* otherTree,
4554 const std::vector< ElementsOnShape::Classifier >& clsOther,
4555 std::vector< ElementsOnShape::Classifier >& cls )
4556 :SMESH_Octree( new SMESH_TreeLimit )
4558 myBox = new Bnd_B3d( *otherTree->getBox() );
4560 if (( myIsLeaf = otherTree->isLeaf() ))
4562 myClassifiers.resize( otherTree->myClassifiers.size() );
4563 for ( size_t i = 0; i < otherTree->myClassifiers.size(); ++i )
4565 int ind = otherTree->myClassifiers[i] - & clsOther[0];
4566 myClassifiers[ i ] = & cls[ ind ];
4569 else if ( otherTree->myChildren )
4571 myChildren = new SMESH_Tree< Bnd_B3d, 8 > * [ 8 ];
4572 for ( int i = 0; i < nbChildren(); i++ )
4574 new OctreeClassifier( static_cast<const OctreeClassifier*>( otherTree->myChildren[i]),
4579 void ElementsOnShape::
4580 OctreeClassifier::GetClassifiersAtPoint( const gp_XYZ& point,
4581 std::vector< ElementsOnShape::Classifier* >& result )
4583 if ( getBox()->IsOut( point ))
4588 for ( size_t i = 0; i < myClassifiers.size(); ++i )
4589 if ( !myClassifiers[i]->GetBndBox()->IsOut( point ))
4590 result.push_back( myClassifiers[i] );
4594 for (int i = 0; i < nbChildren(); i++)
4595 ((OctreeClassifier*) myChildren[i])->GetClassifiersAtPoint( point, result );
4599 void ElementsOnShape::OctreeClassifier::buildChildrenData()
4601 // distribute myClassifiers among myChildren
4603 const int childFlag[8] = { 0x0000001,
4611 int nbInChild[8] = { 0, 0, 0, 0, 0, 0, 0, 0 };
4613 for ( size_t i = 0; i < myClassifiers.size(); ++i )
4615 for ( int j = 0; j < nbChildren(); j++ )
4617 if ( !myClassifiers[i]->GetBndBox()->IsOut( *myChildren[j]->getBox() ))
4619 myClassifiers[i]->SetFlag( childFlag[ j ]);
4625 for ( int j = 0; j < nbChildren(); j++ )
4627 OctreeClassifier* child = static_cast<OctreeClassifier*>( myChildren[ j ]);
4628 child->myClassifiers.resize( nbInChild[ j ]);
4629 for ( size_t i = 0; nbInChild[ j ] && i < myClassifiers.size(); ++i )
4631 if ( myClassifiers[ i ]->IsSetFlag( childFlag[ j ]))
4634 child->myClassifiers[ nbInChild[ j ]] = myClassifiers[ i ];
4635 myClassifiers[ i ]->UnsetFlag( childFlag[ j ]);
4639 SMESHUtils::FreeVector( myClassifiers );
4641 // define if a child isLeaf()
4642 for ( int i = 0; i < nbChildren(); i++ )
4644 OctreeClassifier* child = static_cast<OctreeClassifier*>( myChildren[ i ]);
4645 child->myIsLeaf = ( child->myClassifiers.size() <= 5 );
4649 Bnd_B3d* ElementsOnShape::OctreeClassifier::buildRootBox()
4651 Bnd_B3d* box = new Bnd_B3d;
4652 for ( size_t i = 0; i < myClassifiers.size(); ++i )
4653 box->Add( *myClassifiers[i]->GetBndBox() );
4658 Class : BelongToGeom
4659 Description : Predicate for verifying whether entity belongs to
4660 specified geometrical support
4663 BelongToGeom::BelongToGeom()
4665 myType(SMDSAbs_NbElementTypes),
4666 myIsSubshape(false),
4667 myTolerance(Precision::Confusion())
4670 Predicate* BelongToGeom::clone() const
4672 BelongToGeom* cln = new BelongToGeom( *this );
4673 cln->myElementsOnShapePtr.reset( static_cast<ElementsOnShape*>( myElementsOnShapePtr->clone() ));
4677 void BelongToGeom::SetMesh( const SMDS_Mesh* theMesh )
4679 if ( myMeshDS != theMesh )
4681 myMeshDS = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
4686 void BelongToGeom::SetGeom( const TopoDS_Shape& theShape )
4688 if ( myShape != theShape )
4695 static bool IsSubShape (const TopTools_IndexedMapOfShape& theMap,
4696 const TopoDS_Shape& theShape)
4698 if (theMap.Contains(theShape)) return true;
4700 if (theShape.ShapeType() == TopAbs_COMPOUND ||
4701 theShape.ShapeType() == TopAbs_COMPSOLID)
4703 TopoDS_Iterator anIt (theShape, Standard_True, Standard_True);
4704 for (; anIt.More(); anIt.Next())
4706 if (!IsSubShape(theMap, anIt.Value())) {
4716 void BelongToGeom::init()
4718 if ( !myMeshDS || myShape.IsNull() ) return;
4720 // is sub-shape of main shape?
4721 TopoDS_Shape aMainShape = myMeshDS->ShapeToMesh();
4722 if (aMainShape.IsNull()) {
4723 myIsSubshape = false;
4726 TopTools_IndexedMapOfShape aMap;
4727 TopExp::MapShapes( aMainShape, aMap );
4728 myIsSubshape = IsSubShape( aMap, myShape );
4732 TopExp::MapShapes( myShape, aMap );
4733 mySubShapesIDs.Clear();
4734 for ( int i = 1; i <= aMap.Extent(); ++i )
4736 int subID = myMeshDS->ShapeToIndex( aMap( i ));
4738 mySubShapesIDs.Add( subID );
4743 //if (!myIsSubshape) // to be always ready to check an element not bound to geometry
4745 if ( !myElementsOnShapePtr )
4746 myElementsOnShapePtr.reset( new ElementsOnShape() );
4747 myElementsOnShapePtr->SetTolerance( myTolerance );
4748 myElementsOnShapePtr->SetAllNodes( true ); // "belong", while false means "lays on"
4749 myElementsOnShapePtr->SetMesh( myMeshDS );
4750 myElementsOnShapePtr->SetShape( myShape, myType );
4754 bool BelongToGeom::IsSatisfy (long theId)
4756 if (myMeshDS == 0 || myShape.IsNull())
4761 return myElementsOnShapePtr->IsSatisfy(theId);
4766 if (myType == SMDSAbs_Node)
4768 if ( const SMDS_MeshNode* aNode = myMeshDS->FindNode( theId ))
4770 if ( aNode->getshapeId() < 1 )
4771 return myElementsOnShapePtr->IsSatisfy(theId);
4773 return mySubShapesIDs.Contains( aNode->getshapeId() );
4778 if ( const SMDS_MeshElement* anElem = myMeshDS->FindElement( theId ))
4780 if ( anElem->GetType() == myType )
4782 if ( anElem->getshapeId() < 1 )
4783 return myElementsOnShapePtr->IsSatisfy(theId);
4785 return mySubShapesIDs.Contains( anElem->getshapeId() );
4793 void BelongToGeom::SetType (SMDSAbs_ElementType theType)
4795 if ( myType != theType )
4802 SMDSAbs_ElementType BelongToGeom::GetType() const
4807 TopoDS_Shape BelongToGeom::GetShape()
4812 const SMESHDS_Mesh* BelongToGeom::GetMeshDS() const
4817 void BelongToGeom::SetTolerance (double theTolerance)
4819 myTolerance = theTolerance;
4823 double BelongToGeom::GetTolerance()
4830 Description : Predicate for verifying whether entiy lying or partially lying on
4831 specified geometrical support
4834 LyingOnGeom::LyingOnGeom()
4836 myType(SMDSAbs_NbElementTypes),
4837 myIsSubshape(false),
4838 myTolerance(Precision::Confusion())
4841 Predicate* LyingOnGeom::clone() const
4843 LyingOnGeom* cln = new LyingOnGeom( *this );
4844 cln->myElementsOnShapePtr.reset( static_cast<ElementsOnShape*>( myElementsOnShapePtr->clone() ));
4848 void LyingOnGeom::SetMesh( const SMDS_Mesh* theMesh )
4850 if ( myMeshDS != theMesh )
4852 myMeshDS = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
4857 void LyingOnGeom::SetGeom( const TopoDS_Shape& theShape )
4859 if ( myShape != theShape )
4866 void LyingOnGeom::init()
4868 if (!myMeshDS || myShape.IsNull()) return;
4870 // is sub-shape of main shape?
4871 TopoDS_Shape aMainShape = myMeshDS->ShapeToMesh();
4872 if (aMainShape.IsNull()) {
4873 myIsSubshape = false;
4876 myIsSubshape = myMeshDS->IsGroupOfSubShapes( myShape );
4881 TopTools_IndexedMapOfShape shapes;
4882 TopExp::MapShapes( myShape, shapes );
4883 mySubShapesIDs.Clear();
4884 for ( int i = 1; i <= shapes.Extent(); ++i )
4886 int subID = myMeshDS->ShapeToIndex( shapes( i ));
4888 mySubShapesIDs.Add( subID );
4891 // else // to be always ready to check an element not bound to geometry
4893 if ( !myElementsOnShapePtr )
4894 myElementsOnShapePtr.reset( new ElementsOnShape() );
4895 myElementsOnShapePtr->SetTolerance( myTolerance );
4896 myElementsOnShapePtr->SetAllNodes( false ); // lays on, while true means "belong"
4897 myElementsOnShapePtr->SetMesh( myMeshDS );
4898 myElementsOnShapePtr->SetShape( myShape, myType );
4902 bool LyingOnGeom::IsSatisfy( long theId )
4904 if ( myMeshDS == 0 || myShape.IsNull() )
4909 return myElementsOnShapePtr->IsSatisfy(theId);
4914 const SMDS_MeshElement* elem =
4915 ( myType == SMDSAbs_Node ) ? myMeshDS->FindNode( theId ) : myMeshDS->FindElement( theId );
4917 if ( mySubShapesIDs.Contains( elem->getshapeId() ))
4920 if ( elem->GetType() != SMDSAbs_Node && elem->GetType() == myType )
4922 SMDS_ElemIteratorPtr nodeItr = elem->nodesIterator();
4923 while ( nodeItr->more() )
4925 const SMDS_MeshElement* aNode = nodeItr->next();
4926 if ( mySubShapesIDs.Contains( aNode->getshapeId() ))
4934 void LyingOnGeom::SetType( SMDSAbs_ElementType theType )
4936 if ( myType != theType )
4943 SMDSAbs_ElementType LyingOnGeom::GetType() const
4948 TopoDS_Shape LyingOnGeom::GetShape()
4953 const SMESHDS_Mesh* LyingOnGeom::GetMeshDS() const
4958 void LyingOnGeom::SetTolerance (double theTolerance)
4960 myTolerance = theTolerance;
4964 double LyingOnGeom::GetTolerance()
4969 TSequenceOfXYZ::TSequenceOfXYZ(): myElem(0)
4972 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n), myElem(0)
4975 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t), myElem(0)
4978 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray), myElem(theSequenceOfXYZ.myElem)
4981 template <class InputIterator>
4982 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd), myElem(0)
4985 TSequenceOfXYZ::~TSequenceOfXYZ()
4988 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
4990 myArray = theSequenceOfXYZ.myArray;
4991 myElem = theSequenceOfXYZ.myElem;
4995 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
4997 return myArray[n-1];
5000 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
5002 return myArray[n-1];
5005 void TSequenceOfXYZ::clear()
5010 void TSequenceOfXYZ::reserve(size_type n)
5015 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
5017 myArray.push_back(v);
5020 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
5022 return myArray.size();
5025 SMDSAbs_EntityType TSequenceOfXYZ::getElementEntity() const
5027 return myElem ? myElem->GetEntityType() : SMDSEntity_Last;
5030 TMeshModifTracer::TMeshModifTracer():
5031 myMeshModifTime(0), myMesh(0)
5034 void TMeshModifTracer::SetMesh( const SMDS_Mesh* theMesh )
5036 if ( theMesh != myMesh )
5037 myMeshModifTime = 0;
5040 bool TMeshModifTracer::IsMeshModified()
5042 bool modified = false;
5045 modified = ( myMeshModifTime != myMesh->GetMTime() );
5046 myMeshModifTime = myMesh->GetMTime();