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;
2124 //================================================================================
2126 Class : BadOrientedVolume
2127 Description : Predicate bad oriented volumes
2129 //================================================================================
2131 BadOrientedVolume::BadOrientedVolume()
2136 void BadOrientedVolume::SetMesh( const SMDS_Mesh* theMesh )
2141 bool BadOrientedVolume::IsSatisfy( long theId )
2146 SMDS_VolumeTool vTool( myMesh->FindElement( theId ));
2147 return !vTool.IsForward();
2150 SMDSAbs_ElementType BadOrientedVolume::GetType() const
2152 return SMDSAbs_Volume;
2156 Class : BareBorderVolume
2159 bool BareBorderVolume::IsSatisfy(long theElementId )
2161 SMDS_VolumeTool myTool;
2162 if ( myTool.Set( myMesh->FindElement(theElementId)))
2164 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2165 if ( myTool.IsFreeFace( iF ))
2167 const SMDS_MeshNode** n = myTool.GetFaceNodes(iF);
2168 std::vector< const SMDS_MeshNode*> nodes( n, n+myTool.NbFaceNodes(iF));
2169 if ( !myMesh->FindElement( nodes, SMDSAbs_Face, /*Nomedium=*/false))
2176 //================================================================================
2178 Class : BareBorderFace
2180 //================================================================================
2182 bool BareBorderFace::IsSatisfy(long theElementId )
2185 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2187 if ( face->GetType() == SMDSAbs_Face )
2189 int nbN = face->NbCornerNodes();
2190 for ( int i = 0; i < nbN && !ok; ++i )
2192 // check if a link is shared by another face
2193 const SMDS_MeshNode* n1 = face->GetNode( i );
2194 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2195 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2196 bool isShared = false;
2197 while ( !isShared && fIt->more() )
2199 const SMDS_MeshElement* f = fIt->next();
2200 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2204 const int iQuad = face->IsQuadratic();
2205 myLinkNodes.resize( 2 + iQuad);
2206 myLinkNodes[0] = n1;
2207 myLinkNodes[1] = n2;
2209 myLinkNodes[2] = face->GetNode( i+nbN );
2210 ok = !myMesh->FindElement( myLinkNodes, SMDSAbs_Edge, /*noMedium=*/false);
2218 //================================================================================
2220 Class : OverConstrainedVolume
2222 //================================================================================
2224 bool OverConstrainedVolume::IsSatisfy(long theElementId )
2226 // An element is over-constrained if it has N-1 free borders where
2227 // N is the number of edges/faces for a 2D/3D element.
2228 SMDS_VolumeTool myTool;
2229 if ( myTool.Set( myMesh->FindElement(theElementId)))
2231 int nbSharedFaces = 0;
2232 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2233 if ( !myTool.IsFreeFace( iF ) && ++nbSharedFaces > 1 )
2235 return ( nbSharedFaces == 1 );
2240 //================================================================================
2242 Class : OverConstrainedFace
2244 //================================================================================
2246 bool OverConstrainedFace::IsSatisfy(long theElementId )
2248 // An element is over-constrained if it has N-1 free borders where
2249 // N is the number of edges/faces for a 2D/3D element.
2250 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2251 if ( face->GetType() == SMDSAbs_Face )
2253 int nbSharedBorders = 0;
2254 int nbN = face->NbCornerNodes();
2255 for ( int i = 0; i < nbN; ++i )
2257 // check if a link is shared by another face
2258 const SMDS_MeshNode* n1 = face->GetNode( i );
2259 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2260 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2261 bool isShared = false;
2262 while ( !isShared && fIt->more() )
2264 const SMDS_MeshElement* f = fIt->next();
2265 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2267 if ( isShared && ++nbSharedBorders > 1 )
2270 return ( nbSharedBorders == 1 );
2275 //================================================================================
2277 Class : CoincidentNodes
2278 Description : Predicate of Coincident nodes
2280 //================================================================================
2282 CoincidentNodes::CoincidentNodes()
2287 bool CoincidentNodes::IsSatisfy( long theElementId )
2289 return myCoincidentIDs.Contains( theElementId );
2292 SMDSAbs_ElementType CoincidentNodes::GetType() const
2294 return SMDSAbs_Node;
2297 void CoincidentNodes::SetMesh( const SMDS_Mesh* theMesh )
2299 myMeshModifTracer.SetMesh( theMesh );
2300 if ( myMeshModifTracer.IsMeshModified() )
2302 TIDSortedNodeSet nodesToCheck;
2303 SMDS_NodeIteratorPtr nIt = theMesh->nodesIterator(/*idInceasingOrder=*/true);
2304 while ( nIt->more() )
2305 nodesToCheck.insert( nodesToCheck.end(), nIt->next() );
2307 std::list< std::list< const SMDS_MeshNode*> > nodeGroups;
2308 SMESH_OctreeNode::FindCoincidentNodes ( nodesToCheck, &nodeGroups, myToler );
2310 myCoincidentIDs.Clear();
2311 std::list< std::list< const SMDS_MeshNode*> >::iterator groupIt = nodeGroups.begin();
2312 for ( ; groupIt != nodeGroups.end(); ++groupIt )
2314 std::list< const SMDS_MeshNode*>& coincNodes = *groupIt;
2315 std::list< const SMDS_MeshNode*>::iterator n = coincNodes.begin();
2316 for ( ; n != coincNodes.end(); ++n )
2317 myCoincidentIDs.Add( (*n)->GetID() );
2322 //================================================================================
2324 Class : CoincidentElements
2325 Description : Predicate of Coincident Elements
2326 Note : This class is suitable only for visualization of Coincident Elements
2328 //================================================================================
2330 CoincidentElements::CoincidentElements()
2335 void CoincidentElements::SetMesh( const SMDS_Mesh* theMesh )
2340 bool CoincidentElements::IsSatisfy( long theElementId )
2342 if ( !myMesh ) return false;
2344 if ( const SMDS_MeshElement* e = myMesh->FindElement( theElementId ))
2346 if ( e->GetType() != GetType() ) return false;
2347 std::set< const SMDS_MeshNode* > elemNodes( e->begin_nodes(), e->end_nodes() );
2348 const int nbNodes = e->NbNodes();
2349 SMDS_ElemIteratorPtr invIt = (*elemNodes.begin())->GetInverseElementIterator( GetType() );
2350 while ( invIt->more() )
2352 const SMDS_MeshElement* e2 = invIt->next();
2353 if ( e2 == e || e2->NbNodes() != nbNodes ) continue;
2355 bool sameNodes = true;
2356 for ( size_t i = 0; i < elemNodes.size() && sameNodes; ++i )
2357 sameNodes = ( elemNodes.count( e2->GetNode( i )));
2365 SMDSAbs_ElementType CoincidentElements1D::GetType() const
2367 return SMDSAbs_Edge;
2369 SMDSAbs_ElementType CoincidentElements2D::GetType() const
2371 return SMDSAbs_Face;
2373 SMDSAbs_ElementType CoincidentElements3D::GetType() const
2375 return SMDSAbs_Volume;
2379 //================================================================================
2382 Description : Predicate for free borders
2384 //================================================================================
2386 FreeBorders::FreeBorders()
2391 void FreeBorders::SetMesh( const SMDS_Mesh* theMesh )
2396 bool FreeBorders::IsSatisfy( long theId )
2398 return getNbMultiConnection( myMesh, theId ) == 1;
2401 SMDSAbs_ElementType FreeBorders::GetType() const
2403 return SMDSAbs_Edge;
2407 //================================================================================
2410 Description : Predicate for free Edges
2412 //================================================================================
2414 FreeEdges::FreeEdges()
2419 void FreeEdges::SetMesh( const SMDS_Mesh* theMesh )
2424 bool FreeEdges::IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId )
2426 TColStd_MapOfInteger aMap;
2427 for ( int i = 0; i < 2; i++ )
2429 SMDS_ElemIteratorPtr anElemIter = theNodes[ i ]->GetInverseElementIterator(SMDSAbs_Face);
2430 while( anElemIter->more() )
2432 if ( const SMDS_MeshElement* anElem = anElemIter->next())
2434 const int anId = anElem->GetID();
2435 if ( anId != theFaceId && !aMap.Add( anId ))
2443 bool FreeEdges::IsSatisfy( long theId )
2448 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2449 if ( aFace == 0 || aFace->GetType() != SMDSAbs_Face || aFace->NbNodes() < 3 )
2452 SMDS_NodeIteratorPtr anIter = aFace->interlacedNodesIterator();
2456 int i = 0, nbNodes = aFace->NbNodes();
2457 std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
2458 while( anIter->more() )
2459 if ( ! ( aNodes[ i++ ] = anIter->next() ))
2461 aNodes[ nbNodes ] = aNodes[ 0 ];
2463 for ( i = 0; i < nbNodes; i++ )
2464 if ( IsFreeEdge( &aNodes[ i ], theId ) )
2470 SMDSAbs_ElementType FreeEdges::GetType() const
2472 return SMDSAbs_Face;
2475 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2478 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
2479 if(thePntId1 > thePntId2){
2480 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
2484 bool FreeEdges::Border::operator<(const FreeEdges::Border& x) const{
2485 if(myPntId[0] < x.myPntId[0]) return true;
2486 if(myPntId[0] == x.myPntId[0])
2487 if(myPntId[1] < x.myPntId[1]) return true;
2491 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2492 FreeEdges::TBorders& theRegistry,
2493 FreeEdges::TBorders& theContainer)
2495 if(theRegistry.find(theBorder) == theRegistry.end()){
2496 theRegistry.insert(theBorder);
2497 theContainer.insert(theBorder);
2499 theContainer.erase(theBorder);
2503 void FreeEdges::GetBoreders(TBorders& theBorders)
2506 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2507 for(; anIter->more(); ){
2508 const SMDS_MeshFace* anElem = anIter->next();
2509 long anElemId = anElem->GetID();
2510 SMDS_ElemIteratorPtr aNodesIter;
2511 if ( anElem->IsQuadratic() )
2512 aNodesIter = static_cast<const SMDS_VtkFace*>(anElem)->
2513 interlacedNodesElemIterator();
2515 aNodesIter = anElem->nodesIterator();
2516 long aNodeId[2] = {0,0};
2517 const SMDS_MeshElement* aNode;
2518 if(aNodesIter->more()){
2519 aNode = aNodesIter->next();
2520 aNodeId[0] = aNodeId[1] = aNode->GetID();
2522 for(; aNodesIter->more(); ){
2523 aNode = aNodesIter->next();
2524 long anId = aNode->GetID();
2525 Border aBorder(anElemId,aNodeId[1],anId);
2527 UpdateBorders(aBorder,aRegistry,theBorders);
2529 Border aBorder(anElemId,aNodeId[0],aNodeId[1]);
2530 UpdateBorders(aBorder,aRegistry,theBorders);
2534 //================================================================================
2537 Description : Predicate for free nodes
2539 //================================================================================
2541 FreeNodes::FreeNodes()
2546 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2551 bool FreeNodes::IsSatisfy( long theNodeId )
2553 const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2557 return (aNode->NbInverseElements() < 1);
2560 SMDSAbs_ElementType FreeNodes::GetType() const
2562 return SMDSAbs_Node;
2566 //================================================================================
2569 Description : Predicate for free faces
2571 //================================================================================
2573 FreeFaces::FreeFaces()
2578 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2583 bool FreeFaces::IsSatisfy( long theId )
2585 if (!myMesh) return false;
2586 // check that faces nodes refers to less than two common volumes
2587 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2588 if ( !aFace || aFace->GetType() != SMDSAbs_Face )
2591 int nbNode = aFace->NbNodes();
2593 // collect volumes to check that number of volumes with count equal nbNode not less than 2
2594 typedef std::map< SMDS_MeshElement*, int > TMapOfVolume; // map of volume counters
2595 typedef std::map< SMDS_MeshElement*, int >::iterator TItrMapOfVolume; // iterator
2596 TMapOfVolume mapOfVol;
2598 SMDS_ElemIteratorPtr nodeItr = aFace->nodesIterator();
2599 while ( nodeItr->more() )
2601 const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>(nodeItr->next());
2602 if ( !aNode ) continue;
2603 SMDS_ElemIteratorPtr volItr = aNode->GetInverseElementIterator(SMDSAbs_Volume);
2604 while ( volItr->more() )
2606 SMDS_MeshElement* aVol = (SMDS_MeshElement*)volItr->next();
2607 TItrMapOfVolume itr = mapOfVol.insert( std::make_pair( aVol, 0 )).first;
2612 TItrMapOfVolume volItr = mapOfVol.begin();
2613 TItrMapOfVolume volEnd = mapOfVol.end();
2614 for ( ; volItr != volEnd; ++volItr )
2615 if ( (*volItr).second >= nbNode )
2617 // face is not free if number of volumes constructed on thier nodes more than one
2621 SMDSAbs_ElementType FreeFaces::GetType() const
2623 return SMDSAbs_Face;
2626 //================================================================================
2628 Class : LinearOrQuadratic
2629 Description : Predicate to verify whether a mesh element is linear
2631 //================================================================================
2633 LinearOrQuadratic::LinearOrQuadratic()
2638 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2643 bool LinearOrQuadratic::IsSatisfy( long theId )
2645 if (!myMesh) return false;
2646 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2647 if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2649 return (!anElem->IsQuadratic());
2652 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2657 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2662 //================================================================================
2665 Description : Functor for check color of group to whic mesh element belongs to
2667 //================================================================================
2669 GroupColor::GroupColor()
2673 bool GroupColor::IsSatisfy( long theId )
2675 return myIDs.count( theId );
2678 void GroupColor::SetType( SMDSAbs_ElementType theType )
2683 SMDSAbs_ElementType GroupColor::GetType() const
2688 static bool isEqual( const Quantity_Color& theColor1,
2689 const Quantity_Color& theColor2 )
2691 // tolerance to compare colors
2692 const double tol = 5*1e-3;
2693 return ( fabs( theColor1.Red() - theColor2.Red() ) < tol &&
2694 fabs( theColor1.Green() - theColor2.Green() ) < tol &&
2695 fabs( theColor1.Blue() - theColor2.Blue() ) < tol );
2698 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2702 const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2706 int nbGrp = aMesh->GetNbGroups();
2710 // iterates on groups and find necessary elements ids
2711 const std::set<SMESHDS_GroupBase*>& aGroups = aMesh->GetGroups();
2712 std::set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
2713 for (; GrIt != aGroups.end(); GrIt++)
2715 SMESHDS_GroupBase* aGrp = (*GrIt);
2718 // check type and color of group
2719 if ( !isEqual( myColor, aGrp->GetColor() ))
2722 // IPAL52867 (prevent infinite recursion via GroupOnFilter)
2723 if ( SMESHDS_GroupOnFilter * gof = dynamic_cast< SMESHDS_GroupOnFilter* >( aGrp ))
2724 if ( gof->GetPredicate().get() == this )
2727 SMDSAbs_ElementType aGrpElType = (SMDSAbs_ElementType)aGrp->GetType();
2728 if ( myType == aGrpElType || (myType == SMDSAbs_All && aGrpElType != SMDSAbs_Node) ) {
2729 // add elements IDS into control
2730 int aSize = aGrp->Extent();
2731 for (int i = 0; i < aSize; i++)
2732 myIDs.insert( aGrp->GetID(i+1) );
2737 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
2739 Kernel_Utils::Localizer loc;
2740 TCollection_AsciiString aStr = theStr;
2741 aStr.RemoveAll( ' ' );
2742 aStr.RemoveAll( '\t' );
2743 for ( int aPos = aStr.Search( ";;" ); aPos != -1; aPos = aStr.Search( ";;" ) )
2744 aStr.Remove( aPos, 2 );
2745 Standard_Real clr[3];
2746 clr[0] = clr[1] = clr[2] = 0.;
2747 for ( int i = 0; i < 3; i++ ) {
2748 TCollection_AsciiString tmpStr = aStr.Token( ";", i+1 );
2749 if ( !tmpStr.IsEmpty() && tmpStr.IsRealValue() )
2750 clr[i] = tmpStr.RealValue();
2752 myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
2755 //=======================================================================
2756 // name : GetRangeStr
2757 // Purpose : Get range as a string.
2758 // Example: "1,2,3,50-60,63,67,70-"
2759 //=======================================================================
2761 void GroupColor::GetColorStr( TCollection_AsciiString& theResStr ) const
2764 theResStr += TCollection_AsciiString( myColor.Red() );
2765 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2766 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2769 //================================================================================
2771 Class : ElemGeomType
2772 Description : Predicate to check element geometry type
2774 //================================================================================
2776 ElemGeomType::ElemGeomType()
2779 myType = SMDSAbs_All;
2780 myGeomType = SMDSGeom_TRIANGLE;
2783 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2788 bool ElemGeomType::IsSatisfy( long theId )
2790 if (!myMesh) return false;
2791 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2794 const SMDSAbs_ElementType anElemType = anElem->GetType();
2795 if ( myType != SMDSAbs_All && anElemType != myType )
2797 bool isOk = ( anElem->GetGeomType() == myGeomType );
2801 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2806 SMDSAbs_ElementType ElemGeomType::GetType() const
2811 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2813 myGeomType = theType;
2816 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2821 //================================================================================
2823 Class : ElemEntityType
2824 Description : Predicate to check element entity type
2826 //================================================================================
2828 ElemEntityType::ElemEntityType():
2830 myType( SMDSAbs_All ),
2831 myEntityType( SMDSEntity_0D )
2835 void ElemEntityType::SetMesh( const SMDS_Mesh* theMesh )
2840 bool ElemEntityType::IsSatisfy( long theId )
2842 if ( !myMesh ) return false;
2843 if ( myType == SMDSAbs_Node )
2844 return myMesh->FindNode( theId );
2845 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2847 myEntityType == anElem->GetEntityType() );
2850 void ElemEntityType::SetType( SMDSAbs_ElementType theType )
2855 SMDSAbs_ElementType ElemEntityType::GetType() const
2860 void ElemEntityType::SetElemEntityType( SMDSAbs_EntityType theEntityType )
2862 myEntityType = theEntityType;
2865 SMDSAbs_EntityType ElemEntityType::GetElemEntityType() const
2867 return myEntityType;
2870 //================================================================================
2872 * \brief Class ConnectedElements
2874 //================================================================================
2876 ConnectedElements::ConnectedElements():
2877 myNodeID(0), myType( SMDSAbs_All ), myOkIDsReady( false ) {}
2879 SMDSAbs_ElementType ConnectedElements::GetType() const
2882 int ConnectedElements::GetNode() const
2883 { return myXYZ.empty() ? myNodeID : 0; } // myNodeID can be found by myXYZ
2885 std::vector<double> ConnectedElements::GetPoint() const
2888 void ConnectedElements::clearOkIDs()
2889 { myOkIDsReady = false; myOkIDs.clear(); }
2891 void ConnectedElements::SetType( SMDSAbs_ElementType theType )
2893 if ( myType != theType || myMeshModifTracer.IsMeshModified() )
2898 void ConnectedElements::SetMesh( const SMDS_Mesh* theMesh )
2900 myMeshModifTracer.SetMesh( theMesh );
2901 if ( myMeshModifTracer.IsMeshModified() )
2904 if ( !myXYZ.empty() )
2905 SetPoint( myXYZ[0], myXYZ[1], myXYZ[2] ); // find a node near myXYZ it in a new mesh
2909 void ConnectedElements::SetNode( int nodeID )
2914 bool isSameDomain = false;
2915 if ( myOkIDsReady && myMeshModifTracer.GetMesh() && !myMeshModifTracer.IsMeshModified() )
2916 if ( const SMDS_MeshNode* n = myMeshModifTracer.GetMesh()->FindNode( myNodeID ))
2918 SMDS_ElemIteratorPtr eIt = n->GetInverseElementIterator( myType );
2919 while ( !isSameDomain && eIt->more() )
2920 isSameDomain = IsSatisfy( eIt->next()->GetID() );
2922 if ( !isSameDomain )
2926 void ConnectedElements::SetPoint( double x, double y, double z )
2934 bool isSameDomain = false;
2936 // find myNodeID by myXYZ if possible
2937 if ( myMeshModifTracer.GetMesh() )
2939 SMESHUtils::Deleter<SMESH_ElementSearcher> searcher
2940 ( SMESH_MeshAlgos::GetElementSearcher( (SMDS_Mesh&) *myMeshModifTracer.GetMesh() ));
2942 std::vector< const SMDS_MeshElement* > foundElems;
2943 searcher->FindElementsByPoint( gp_Pnt(x,y,z), SMDSAbs_All, foundElems );
2945 if ( !foundElems.empty() )
2947 myNodeID = foundElems[0]->GetNode(0)->GetID();
2948 if ( myOkIDsReady && !myMeshModifTracer.IsMeshModified() )
2949 isSameDomain = IsSatisfy( foundElems[0]->GetID() );
2952 if ( !isSameDomain )
2956 bool ConnectedElements::IsSatisfy( long theElementId )
2958 // Here we do NOT check if the mesh has changed, we do it in Set...() only!!!
2960 if ( !myOkIDsReady )
2962 if ( !myMeshModifTracer.GetMesh() )
2964 const SMDS_MeshNode* node0 = myMeshModifTracer.GetMesh()->FindNode( myNodeID );
2968 std::list< const SMDS_MeshNode* > nodeQueue( 1, node0 );
2969 std::set< int > checkedNodeIDs;
2971 // foreach node in nodeQueue:
2972 // foreach element sharing a node:
2973 // add ID of an element of myType to myOkIDs;
2974 // push all element nodes absent from checkedNodeIDs to nodeQueue;
2975 while ( !nodeQueue.empty() )
2977 const SMDS_MeshNode* node = nodeQueue.front();
2978 nodeQueue.pop_front();
2980 // loop on elements sharing the node
2981 SMDS_ElemIteratorPtr eIt = node->GetInverseElementIterator();
2982 while ( eIt->more() )
2984 // keep elements of myType
2985 const SMDS_MeshElement* element = eIt->next();
2986 if ( element->GetType() == myType )
2987 myOkIDs.insert( myOkIDs.end(), element->GetID() );
2989 // enqueue nodes of the element
2990 SMDS_ElemIteratorPtr nIt = element->nodesIterator();
2991 while ( nIt->more() )
2993 const SMDS_MeshNode* n = static_cast< const SMDS_MeshNode* >( nIt->next() );
2994 if ( checkedNodeIDs.insert( n->GetID() ).second )
2995 nodeQueue.push_back( n );
2999 if ( myType == SMDSAbs_Node )
3000 std::swap( myOkIDs, checkedNodeIDs );
3002 size_t totalNbElems = myMeshModifTracer.GetMesh()->GetMeshInfo().NbElements( myType );
3003 if ( myOkIDs.size() == totalNbElems )
3006 myOkIDsReady = true;
3009 return myOkIDs.empty() ? true : myOkIDs.count( theElementId );
3012 //================================================================================
3014 * \brief Class CoplanarFaces
3016 //================================================================================
3020 inline bool isLessAngle( const gp_Vec& v1, const gp_Vec& v2, const double cos )
3022 double dot = v1 * v2; // cos * |v1| * |v2|
3023 double l1 = v1.SquareMagnitude();
3024 double l2 = v2.SquareMagnitude();
3025 return (( dot * cos >= 0 ) &&
3026 ( dot * dot ) / l1 / l2 >= ( cos * cos ));
3029 CoplanarFaces::CoplanarFaces()
3030 : myFaceID(0), myToler(0)
3033 void CoplanarFaces::SetMesh( const SMDS_Mesh* theMesh )
3035 myMeshModifTracer.SetMesh( theMesh );
3036 if ( myMeshModifTracer.IsMeshModified() )
3038 // Build a set of coplanar face ids
3040 myCoplanarIDs.Clear();
3042 if ( !myMeshModifTracer.GetMesh() || !myFaceID || !myToler )
3045 const SMDS_MeshElement* face = myMeshModifTracer.GetMesh()->FindElement( myFaceID );
3046 if ( !face || face->GetType() != SMDSAbs_Face )
3050 gp_Vec myNorm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
3054 const double cosTol = Cos( myToler * M_PI / 180. );
3055 NCollection_Map< SMESH_TLink, SMESH_TLink > checkedLinks;
3057 std::list< std::pair< const SMDS_MeshElement*, gp_Vec > > faceQueue;
3058 faceQueue.push_back( std::make_pair( face, myNorm ));
3059 while ( !faceQueue.empty() )
3061 face = faceQueue.front().first;
3062 myNorm = faceQueue.front().second;
3063 faceQueue.pop_front();
3065 for ( int i = 0, nbN = face->NbCornerNodes(); i < nbN; ++i )
3067 const SMDS_MeshNode* n1 = face->GetNode( i );
3068 const SMDS_MeshNode* n2 = face->GetNode(( i+1 )%nbN);
3069 if ( !checkedLinks.Add( SMESH_TLink( n1, n2 )))
3071 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator(SMDSAbs_Face);
3072 while ( fIt->more() )
3074 const SMDS_MeshElement* f = fIt->next();
3075 if ( f->GetNodeIndex( n2 ) > -1 )
3077 gp_Vec norm = getNormale( static_cast<const SMDS_MeshFace*>(f), &normOK );
3078 if (!normOK || isLessAngle( myNorm, norm, cosTol))
3080 myCoplanarIDs.Add( f->GetID() );
3081 faceQueue.push_back( std::make_pair( f, norm ));
3089 bool CoplanarFaces::IsSatisfy( long theElementId )
3091 return myCoplanarIDs.Contains( theElementId );
3096 *Description : Predicate for Range of Ids.
3097 * Range may be specified with two ways.
3098 * 1. Using AddToRange method
3099 * 2. With SetRangeStr method. Parameter of this method is a string
3100 * like as "1,2,3,50-60,63,67,70-"
3103 //=======================================================================
3104 // name : RangeOfIds
3105 // Purpose : Constructor
3106 //=======================================================================
3107 RangeOfIds::RangeOfIds()
3110 myType = SMDSAbs_All;
3113 //=======================================================================
3115 // Purpose : Set mesh
3116 //=======================================================================
3117 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
3122 //=======================================================================
3123 // name : AddToRange
3124 // Purpose : Add ID to the range
3125 //=======================================================================
3126 bool RangeOfIds::AddToRange( long theEntityId )
3128 myIds.Add( theEntityId );
3132 //=======================================================================
3133 // name : GetRangeStr
3134 // Purpose : Get range as a string.
3135 // Example: "1,2,3,50-60,63,67,70-"
3136 //=======================================================================
3137 void RangeOfIds::GetRangeStr( TCollection_AsciiString& theResStr )
3141 TColStd_SequenceOfInteger anIntSeq;
3142 TColStd_SequenceOfAsciiString aStrSeq;
3144 TColStd_MapIteratorOfMapOfInteger anIter( myIds );
3145 for ( ; anIter.More(); anIter.Next() )
3147 int anId = anIter.Key();
3148 TCollection_AsciiString aStr( anId );
3149 anIntSeq.Append( anId );
3150 aStrSeq.Append( aStr );
3153 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
3155 int aMinId = myMin( i );
3156 int aMaxId = myMax( i );
3158 TCollection_AsciiString aStr;
3159 if ( aMinId != IntegerFirst() )
3164 if ( aMaxId != IntegerLast() )
3167 // find position of the string in result sequence and insert string in it
3168 if ( anIntSeq.Length() == 0 )
3170 anIntSeq.Append( aMinId );
3171 aStrSeq.Append( aStr );
3175 if ( aMinId < anIntSeq.First() )
3177 anIntSeq.Prepend( aMinId );
3178 aStrSeq.Prepend( aStr );
3180 else if ( aMinId > anIntSeq.Last() )
3182 anIntSeq.Append( aMinId );
3183 aStrSeq.Append( aStr );
3186 for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
3187 if ( aMinId < anIntSeq( j ) )
3189 anIntSeq.InsertBefore( j, aMinId );
3190 aStrSeq.InsertBefore( j, aStr );
3196 if ( aStrSeq.Length() == 0 )
3199 theResStr = aStrSeq( 1 );
3200 for ( int j = 2, k = aStrSeq.Length(); j <= k; j++ )
3203 theResStr += aStrSeq( j );
3207 //=======================================================================
3208 // name : SetRangeStr
3209 // Purpose : Define range with string
3210 // Example of entry string: "1,2,3,50-60,63,67,70-"
3211 //=======================================================================
3212 bool RangeOfIds::SetRangeStr( const TCollection_AsciiString& theStr )
3218 TCollection_AsciiString aStr = theStr;
3219 for ( int i = 1; i <= aStr.Length(); ++i )
3221 char c = aStr.Value( i );
3222 if ( !isdigit( c ) && c != ',' && c != '-' )
3223 aStr.SetValue( i, ' ');
3225 aStr.RemoveAll( ' ' );
3227 TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
3229 while ( tmpStr != "" )
3231 tmpStr = aStr.Token( ",", i++ );
3232 int aPos = tmpStr.Search( '-' );
3236 if ( tmpStr.IsIntegerValue() )
3237 myIds.Add( tmpStr.IntegerValue() );
3243 TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
3244 TCollection_AsciiString aMinStr = tmpStr;
3246 while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
3247 while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
3249 if ( (!aMinStr.IsEmpty() && !aMinStr.IsIntegerValue()) ||
3250 (!aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue()) )
3253 myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
3254 myMax.Append( aMaxStr.IsEmpty() ? IntegerLast() : aMaxStr.IntegerValue() );
3261 //=======================================================================
3263 // Purpose : Get type of supported entities
3264 //=======================================================================
3265 SMDSAbs_ElementType RangeOfIds::GetType() const
3270 //=======================================================================
3272 // Purpose : Set type of supported entities
3273 //=======================================================================
3274 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
3279 //=======================================================================
3281 // Purpose : Verify whether entity satisfies to this rpedicate
3282 //=======================================================================
3283 bool RangeOfIds::IsSatisfy( long theId )
3288 if ( myType == SMDSAbs_Node )
3290 if ( myMesh->FindNode( theId ) == 0 )
3295 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
3296 if ( anElem == 0 || (myType != anElem->GetType() && myType != SMDSAbs_All ))
3300 if ( myIds.Contains( theId ) )
3303 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
3304 if ( theId >= myMin( i ) && theId <= myMax( i ) )
3312 Description : Base class for comparators
3314 Comparator::Comparator():
3318 Comparator::~Comparator()
3321 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
3324 myFunctor->SetMesh( theMesh );
3327 void Comparator::SetMargin( double theValue )
3329 myMargin = theValue;
3332 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
3334 myFunctor = theFunct;
3337 SMDSAbs_ElementType Comparator::GetType() const
3339 return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
3342 double Comparator::GetMargin()
3350 Description : Comparator "<"
3352 bool LessThan::IsSatisfy( long theId )
3354 return myFunctor && myFunctor->GetValue( theId ) < myMargin;
3360 Description : Comparator ">"
3362 bool MoreThan::IsSatisfy( long theId )
3364 return myFunctor && myFunctor->GetValue( theId ) > myMargin;
3370 Description : Comparator "="
3373 myToler(Precision::Confusion())
3376 bool EqualTo::IsSatisfy( long theId )
3378 return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
3381 void EqualTo::SetTolerance( double theToler )
3386 double EqualTo::GetTolerance()
3393 Description : Logical NOT predicate
3395 LogicalNOT::LogicalNOT()
3398 LogicalNOT::~LogicalNOT()
3401 bool LogicalNOT::IsSatisfy( long theId )
3403 return myPredicate && !myPredicate->IsSatisfy( theId );
3406 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
3409 myPredicate->SetMesh( theMesh );
3412 void LogicalNOT::SetPredicate( PredicatePtr thePred )
3414 myPredicate = thePred;
3417 SMDSAbs_ElementType LogicalNOT::GetType() const
3419 return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
3424 Class : LogicalBinary
3425 Description : Base class for binary logical predicate
3427 LogicalBinary::LogicalBinary()
3430 LogicalBinary::~LogicalBinary()
3433 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
3436 myPredicate1->SetMesh( theMesh );
3439 myPredicate2->SetMesh( theMesh );
3442 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
3444 myPredicate1 = thePredicate;
3447 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
3449 myPredicate2 = thePredicate;
3452 SMDSAbs_ElementType LogicalBinary::GetType() const
3454 if ( !myPredicate1 || !myPredicate2 )
3457 SMDSAbs_ElementType aType1 = myPredicate1->GetType();
3458 SMDSAbs_ElementType aType2 = myPredicate2->GetType();
3460 return aType1 == aType2 ? aType1 : SMDSAbs_All;
3466 Description : Logical AND
3468 bool LogicalAND::IsSatisfy( long theId )
3473 myPredicate1->IsSatisfy( theId ) &&
3474 myPredicate2->IsSatisfy( theId );
3480 Description : Logical OR
3482 bool LogicalOR::IsSatisfy( long theId )
3487 (myPredicate1->IsSatisfy( theId ) ||
3488 myPredicate2->IsSatisfy( theId ));
3497 // #include <tbb/parallel_for.h>
3498 // #include <tbb/enumerable_thread_specific.h>
3500 // namespace Parallel
3502 // typedef tbb::enumerable_thread_specific< TIdSequence > TIdSeq;
3506 // const SMDS_Mesh* myMesh;
3507 // PredicatePtr myPredicate;
3508 // TIdSeq & myOKIds;
3509 // Predicate( const SMDS_Mesh* m, PredicatePtr p, TIdSeq & ids ):
3510 // myMesh(m), myPredicate(p->Duplicate()), myOKIds(ids) {}
3511 // void operator() ( const tbb::blocked_range<size_t>& r ) const
3513 // for ( size_t i = r.begin(); i != r.end(); ++i )
3514 // if ( myPredicate->IsSatisfy( i ))
3515 // myOKIds.local().push_back();
3527 void Filter::SetPredicate( PredicatePtr thePredicate )
3529 myPredicate = thePredicate;
3532 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3533 PredicatePtr thePredicate,
3534 TIdSequence& theSequence )
3536 theSequence.clear();
3538 if ( !theMesh || !thePredicate )
3541 thePredicate->SetMesh( theMesh );
3543 SMDS_ElemIteratorPtr elemIt = theMesh->elementsIterator( thePredicate->GetType() );
3545 while ( elemIt->more() ) {
3546 const SMDS_MeshElement* anElem = elemIt->next();
3547 long anId = anElem->GetID();
3548 if ( thePredicate->IsSatisfy( anId ) )
3549 theSequence.push_back( anId );
3554 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3555 Filter::TIdSequence& theSequence )
3557 GetElementsId(theMesh,myPredicate,theSequence);
3564 typedef std::set<SMDS_MeshFace*> TMapOfFacePtr;
3570 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
3571 SMDS_MeshNode* theNode2 )
3577 ManifoldPart::Link::~Link()
3583 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
3585 if ( myNode1 == theLink.myNode1 &&
3586 myNode2 == theLink.myNode2 )
3588 else if ( myNode1 == theLink.myNode2 &&
3589 myNode2 == theLink.myNode1 )
3595 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
3597 if(myNode1 < x.myNode1) return true;
3598 if(myNode1 == x.myNode1)
3599 if(myNode2 < x.myNode2) return true;
3603 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
3604 const ManifoldPart::Link& theLink2 )
3606 return theLink1.IsEqual( theLink2 );
3609 ManifoldPart::ManifoldPart()
3612 myAngToler = Precision::Angular();
3613 myIsOnlyManifold = true;
3616 ManifoldPart::~ManifoldPart()
3621 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
3627 SMDSAbs_ElementType ManifoldPart::GetType() const
3628 { return SMDSAbs_Face; }
3630 bool ManifoldPart::IsSatisfy( long theElementId )
3632 return myMapIds.Contains( theElementId );
3635 void ManifoldPart::SetAngleTolerance( const double theAngToler )
3636 { myAngToler = theAngToler; }
3638 double ManifoldPart::GetAngleTolerance() const
3639 { return myAngToler; }
3641 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
3642 { myIsOnlyManifold = theIsOnly; }
3644 void ManifoldPart::SetStartElem( const long theStartId )
3645 { myStartElemId = theStartId; }
3647 bool ManifoldPart::process()
3650 myMapBadGeomIds.Clear();
3652 myAllFacePtr.clear();
3653 myAllFacePtrIntDMap.clear();
3657 // collect all faces into own map
3658 SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3659 for (; anFaceItr->more(); )
3661 SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3662 myAllFacePtr.push_back( aFacePtr );
3663 myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3666 SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3670 // the map of non manifold links and bad geometry
3671 TMapOfLink aMapOfNonManifold;
3672 TColStd_MapOfInteger aMapOfTreated;
3674 // begin cycle on faces from start index and run on vector till the end
3675 // and from begin to start index to cover whole vector
3676 const int aStartIndx = myAllFacePtrIntDMap[aStartFace];
3677 bool isStartTreat = false;
3678 for ( int fi = aStartIndx; !isStartTreat || fi != aStartIndx ; fi++ )
3680 if ( fi == aStartIndx )
3681 isStartTreat = true;
3682 // as result next time when fi will be equal to aStartIndx
3684 SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3685 if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3688 aMapOfTreated.Add( aFacePtr->GetID() );
3689 TColStd_MapOfInteger aResFaces;
3690 if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3691 aMapOfNonManifold, aResFaces ) )
3693 TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3694 for ( ; anItr.More(); anItr.Next() )
3696 int aFaceId = anItr.Key();
3697 aMapOfTreated.Add( aFaceId );
3698 myMapIds.Add( aFaceId );
3701 if ( fi == int( myAllFacePtr.size() - 1 ))
3703 } // end run on vector of faces
3704 return !myMapIds.IsEmpty();
3707 static void getLinks( const SMDS_MeshFace* theFace,
3708 ManifoldPart::TVectorOfLink& theLinks )
3710 int aNbNode = theFace->NbNodes();
3711 SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3713 SMDS_MeshNode* aNode = 0;
3714 for ( ; aNodeItr->more() && i <= aNbNode; )
3717 SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3721 SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3723 ManifoldPart::Link aLink( aN1, aN2 );
3724 theLinks.push_back( aLink );
3728 bool ManifoldPart::findConnected
3729 ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3730 SMDS_MeshFace* theStartFace,
3731 ManifoldPart::TMapOfLink& theNonManifold,
3732 TColStd_MapOfInteger& theResFaces )
3734 theResFaces.Clear();
3735 if ( !theAllFacePtrInt.size() )
3738 if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3740 myMapBadGeomIds.Add( theStartFace->GetID() );
3744 ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3745 ManifoldPart::TVectorOfLink aSeqOfBoundary;
3746 theResFaces.Add( theStartFace->GetID() );
3747 ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3749 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3750 aDMapLinkFace, theNonManifold, theStartFace );
3752 bool isDone = false;
3753 while ( !isDone && aMapOfBoundary.size() != 0 )
3755 bool isToReset = false;
3756 ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3757 for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3759 ManifoldPart::Link aLink = *pLink;
3760 if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3762 // each link could be treated only once
3763 aMapToSkip.insert( aLink );
3765 ManifoldPart::TVectorOfFacePtr aFaces;
3767 if ( myIsOnlyManifold &&
3768 (theNonManifold.find( aLink ) != theNonManifold.end()) )
3772 getFacesByLink( aLink, aFaces );
3773 // filter the element to keep only indicated elements
3774 ManifoldPart::TVectorOfFacePtr aFiltered;
3775 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3776 for ( ; pFace != aFaces.end(); ++pFace )
3778 SMDS_MeshFace* aFace = *pFace;
3779 if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3780 aFiltered.push_back( aFace );
3783 if ( aFaces.size() < 2 ) // no neihgbour faces
3785 else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3787 theNonManifold.insert( aLink );
3792 // compare normal with normals of neighbor element
3793 SMDS_MeshFace* aPrevFace = aDMapLinkFace[ aLink ];
3794 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3795 for ( ; pFace != aFaces.end(); ++pFace )
3797 SMDS_MeshFace* aNextFace = *pFace;
3798 if ( aPrevFace == aNextFace )
3800 int anNextFaceID = aNextFace->GetID();
3801 if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3802 // should not be with non manifold restriction. probably bad topology
3804 // check if face was treated and skipped
3805 if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3806 !isInPlane( aPrevFace, aNextFace ) )
3808 // add new element to connected and extend the boundaries.
3809 theResFaces.Add( anNextFaceID );
3810 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3811 aDMapLinkFace, theNonManifold, aNextFace );
3815 isDone = !isToReset;
3818 return !theResFaces.IsEmpty();
3821 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3822 const SMDS_MeshFace* theFace2 )
3824 gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3825 gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3826 if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3828 myMapBadGeomIds.Add( theFace2->GetID() );
3831 if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
3837 void ManifoldPart::expandBoundary
3838 ( ManifoldPart::TMapOfLink& theMapOfBoundary,
3839 ManifoldPart::TVectorOfLink& theSeqOfBoundary,
3840 ManifoldPart::TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
3841 ManifoldPart::TMapOfLink& theNonManifold,
3842 SMDS_MeshFace* theNextFace ) const
3844 ManifoldPart::TVectorOfLink aLinks;
3845 getLinks( theNextFace, aLinks );
3846 int aNbLink = (int)aLinks.size();
3847 for ( int i = 0; i < aNbLink; i++ )
3849 ManifoldPart::Link aLink = aLinks[ i ];
3850 if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3852 if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3854 if ( myIsOnlyManifold )
3856 // remove from boundary
3857 theMapOfBoundary.erase( aLink );
3858 ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3859 for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3861 ManifoldPart::Link aBoundLink = *pLink;
3862 if ( aBoundLink.IsEqual( aLink ) )
3864 theSeqOfBoundary.erase( pLink );
3872 theMapOfBoundary.insert( aLink );
3873 theSeqOfBoundary.push_back( aLink );
3874 theDMapLinkFacePtr[ aLink ] = theNextFace;
3879 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3880 ManifoldPart::TVectorOfFacePtr& theFaces ) const
3882 std::set<SMDS_MeshCell *> aSetOfFaces;
3883 // take all faces that shared first node
3884 SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3885 for ( ; anItr->more(); )
3887 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3890 aSetOfFaces.insert( aFace );
3892 // take all faces that shared second node
3893 anItr = theLink.myNode2->facesIterator();
3894 // find the common part of two sets
3895 for ( ; anItr->more(); )
3897 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3898 if ( aSetOfFaces.count( aFace ) )
3899 theFaces.push_back( aFace );
3904 Class : BelongToMeshGroup
3905 Description : Verify whether a mesh element is included into a mesh group
3907 BelongToMeshGroup::BelongToMeshGroup(): myGroup( 0 )
3911 void BelongToMeshGroup::SetGroup( SMESHDS_GroupBase* g )
3916 void BelongToMeshGroup::SetStoreName( const std::string& sn )
3921 void BelongToMeshGroup::SetMesh( const SMDS_Mesh* theMesh )
3923 if ( myGroup && myGroup->GetMesh() != theMesh )
3927 if ( !myGroup && !myStoreName.empty() )
3929 if ( const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh))
3931 const std::set<SMESHDS_GroupBase*>& grps = aMesh->GetGroups();
3932 std::set<SMESHDS_GroupBase*>::const_iterator g = grps.begin();
3933 for ( ; g != grps.end() && !myGroup; ++g )
3934 if ( *g && myStoreName == (*g)->GetStoreName() )
3940 myGroup->IsEmpty(); // make GroupOnFilter update its predicate
3944 bool BelongToMeshGroup::IsSatisfy( long theElementId )
3946 return myGroup ? myGroup->Contains( theElementId ) : false;
3949 SMDSAbs_ElementType BelongToMeshGroup::GetType() const
3951 return myGroup ? myGroup->GetType() : SMDSAbs_All;
3954 //================================================================================
3955 // ElementsOnSurface
3956 //================================================================================
3958 ElementsOnSurface::ElementsOnSurface()
3961 myType = SMDSAbs_All;
3963 myToler = Precision::Confusion();
3964 myUseBoundaries = false;
3967 ElementsOnSurface::~ElementsOnSurface()
3971 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
3973 myMeshModifTracer.SetMesh( theMesh );
3974 if ( myMeshModifTracer.IsMeshModified())
3978 bool ElementsOnSurface::IsSatisfy( long theElementId )
3980 return myIds.Contains( theElementId );
3983 SMDSAbs_ElementType ElementsOnSurface::GetType() const
3986 void ElementsOnSurface::SetTolerance( const double theToler )
3988 if ( myToler != theToler )
3993 double ElementsOnSurface::GetTolerance() const
3996 void ElementsOnSurface::SetUseBoundaries( bool theUse )
3998 if ( myUseBoundaries != theUse ) {
3999 myUseBoundaries = theUse;
4000 SetSurface( mySurf, myType );
4004 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
4005 const SMDSAbs_ElementType theType )
4010 if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
4012 mySurf = TopoDS::Face( theShape );
4013 BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
4015 u1 = SA.FirstUParameter(),
4016 u2 = SA.LastUParameter(),
4017 v1 = SA.FirstVParameter(),
4018 v2 = SA.LastVParameter();
4019 Handle(Geom_Surface) surf = BRep_Tool::Surface( mySurf );
4020 myProjector.Init( surf, u1,u2, v1,v2 );
4024 void ElementsOnSurface::process()
4027 if ( mySurf.IsNull() )
4030 if ( !myMeshModifTracer.GetMesh() )
4033 myIds.ReSize( myMeshModifTracer.GetMesh()->GetMeshInfo().NbElements( myType ));
4035 SMDS_ElemIteratorPtr anIter = myMeshModifTracer.GetMesh()->elementsIterator( myType );
4036 for(; anIter->more(); )
4037 process( anIter->next() );
4040 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
4042 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
4043 bool isSatisfy = true;
4044 for ( ; aNodeItr->more(); )
4046 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
4047 if ( !isOnSurface( aNode ) )
4054 myIds.Add( theElemPtr->GetID() );
4057 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
4059 if ( mySurf.IsNull() )
4062 gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
4063 // double aToler2 = myToler * myToler;
4064 // if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
4066 // gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
4067 // if ( aPln.SquareDistance( aPnt ) > aToler2 )
4070 // else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
4072 // gp_Cylinder aCyl = Handle(Geom_CylindricalSurface)::DownCast(mySurf)->Cylinder();
4073 // double aRad = aCyl.Radius();
4074 // gp_Ax3 anAxis = aCyl.Position();
4075 // gp_XYZ aLoc = aCyl.Location().XYZ();
4076 // double aXDist = anAxis.XDirection().XYZ() * ( aPnt.XYZ() - aLoc );
4077 // double aYDist = anAxis.YDirection().XYZ() * ( aPnt.XYZ() - aLoc );
4078 // if ( fabs(aXDist*aXDist + aYDist*aYDist - aRad*aRad) > aToler2 )
4083 myProjector.Perform( aPnt );
4084 bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
4090 //================================================================================
4092 //================================================================================
4095 const int theIsCheckedFlag = 0x0000100;
4098 struct ElementsOnShape::Classifier
4100 //Classifier(const TopoDS_Shape& s, double tol) { Init(s,tol); }
4101 void Init(const TopoDS_Shape& s, double tol, const Bnd_B3d* box = 0 );
4102 bool IsOut(const gp_Pnt& p) { return SetChecked( true ), (this->*myIsOutFun)( p ); }
4103 TopAbs_ShapeEnum ShapeType() const { return myShape.ShapeType(); }
4104 const TopoDS_Shape& Shape() const { return myShape; }
4105 const Bnd_B3d* GetBndBox() const { return & myBox; }
4106 bool IsChecked() { return myFlags & theIsCheckedFlag; }
4107 bool IsSetFlag( int flag ) const { return myFlags & flag; }
4108 void SetChecked( bool is ) { is ? SetFlag( theIsCheckedFlag ) : UnsetFlag( theIsCheckedFlag ); }
4109 void SetFlag ( int flag ) { myFlags |= flag; }
4110 void UnsetFlag( int flag ) { myFlags &= ~flag; }
4112 bool isOutOfSolid (const gp_Pnt& p);
4113 bool isOutOfBox (const gp_Pnt& p);
4114 bool isOutOfFace (const gp_Pnt& p);
4115 bool isOutOfEdge (const gp_Pnt& p);
4116 bool isOutOfVertex(const gp_Pnt& p);
4117 bool isBox (const TopoDS_Shape& s);
4119 bool (Classifier::* myIsOutFun)(const gp_Pnt& p);
4120 BRepClass3d_SolidClassifier mySolidClfr;
4122 GeomAPI_ProjectPointOnSurf myProjFace;
4123 GeomAPI_ProjectPointOnCurve myProjEdge;
4125 TopoDS_Shape myShape;
4130 struct ElementsOnShape::OctreeClassifier : public SMESH_Octree
4132 OctreeClassifier( const std::vector< ElementsOnShape::Classifier* >& classifiers );
4133 OctreeClassifier( const OctreeClassifier* otherTree,
4134 const std::vector< ElementsOnShape::Classifier >& clsOther,
4135 std::vector< ElementsOnShape::Classifier >& cls );
4136 void GetClassifiersAtPoint( const gp_XYZ& p,
4137 std::vector< ElementsOnShape::Classifier* >& classifiers );
4139 OctreeClassifier() {}
4140 SMESH_Octree* newChild() const { return new OctreeClassifier; }
4141 void buildChildrenData();
4142 Bnd_B3d* buildRootBox();
4144 std::vector< ElementsOnShape::Classifier* > myClassifiers;
4148 ElementsOnShape::ElementsOnShape():
4150 myType(SMDSAbs_All),
4151 myToler(Precision::Confusion()),
4152 myAllNodesFlag(false)
4156 ElementsOnShape::~ElementsOnShape()
4161 Predicate* ElementsOnShape::clone() const
4163 ElementsOnShape* cln = new ElementsOnShape();
4164 cln->SetAllNodes ( myAllNodesFlag );
4165 cln->SetTolerance( myToler );
4166 cln->SetMesh ( myMeshModifTracer.GetMesh() );
4167 cln->myShape = myShape; // avoid creation of myClassifiers
4168 cln->SetShape ( myShape, myType );
4169 cln->myClassifiers.resize( myClassifiers.size() );
4170 for ( size_t i = 0; i < myClassifiers.size(); ++i )
4171 cln->myClassifiers[ i ].Init( BRepBuilderAPI_Copy( myClassifiers[ i ].Shape()),
4172 myToler, myClassifiers[ i ].GetBndBox() );
4173 if ( myOctree ) // copy myOctree
4175 cln->myOctree = new OctreeClassifier( myOctree, myClassifiers, cln->myClassifiers );
4180 SMDSAbs_ElementType ElementsOnShape::GetType() const
4185 void ElementsOnShape::SetTolerance (const double theToler)
4187 if (myToler != theToler) {
4189 SetShape(myShape, myType);
4193 double ElementsOnShape::GetTolerance() const
4198 void ElementsOnShape::SetAllNodes (bool theAllNodes)
4200 myAllNodesFlag = theAllNodes;
4203 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
4205 myMeshModifTracer.SetMesh( theMesh );
4206 if ( myMeshModifTracer.IsMeshModified())
4208 size_t nbNodes = theMesh ? theMesh->NbNodes() : 0;
4209 if ( myNodeIsChecked.size() == nbNodes )
4211 std::fill( myNodeIsChecked.begin(), myNodeIsChecked.end(), false );
4215 SMESHUtils::FreeVector( myNodeIsChecked );
4216 SMESHUtils::FreeVector( myNodeIsOut );
4217 myNodeIsChecked.resize( nbNodes, false );
4218 myNodeIsOut.resize( nbNodes );
4223 bool ElementsOnShape::getNodeIsOut( const SMDS_MeshNode* n, bool& isOut )
4225 if ( n->GetID() >= (int) myNodeIsChecked.size() ||
4226 !myNodeIsChecked[ n->GetID() ])
4229 isOut = myNodeIsOut[ n->GetID() ];
4233 void ElementsOnShape::setNodeIsOut( const SMDS_MeshNode* n, bool isOut )
4235 if ( n->GetID() < (int) myNodeIsChecked.size() )
4237 myNodeIsChecked[ n->GetID() ] = true;
4238 myNodeIsOut [ n->GetID() ] = isOut;
4242 void ElementsOnShape::SetShape (const TopoDS_Shape& theShape,
4243 const SMDSAbs_ElementType theType)
4245 bool shapeChanges = ( myShape != theShape );
4248 if ( myShape.IsNull() ) return;
4252 TopTools_IndexedMapOfShape shapesMap;
4253 TopAbs_ShapeEnum shapeTypes[4] = { TopAbs_SOLID, TopAbs_FACE, TopAbs_EDGE, TopAbs_VERTEX };
4254 TopExp_Explorer sub;
4255 for ( int i = 0; i < 4; ++i )
4257 if ( shapesMap.IsEmpty() )
4258 for ( sub.Init( myShape, shapeTypes[i] ); sub.More(); sub.Next() )
4259 shapesMap.Add( sub.Current() );
4261 for ( sub.Init( myShape, shapeTypes[i], shapeTypes[i-1] ); sub.More(); sub.Next() )
4262 shapesMap.Add( sub.Current() );
4266 myClassifiers.resize( shapesMap.Extent() );
4267 for ( int i = 0; i < shapesMap.Extent(); ++i )
4268 myClassifiers[ i ].Init( shapesMap( i+1 ), myToler );
4271 if ( theType == SMDSAbs_Node )
4273 SMESHUtils::FreeVector( myNodeIsChecked );
4274 SMESHUtils::FreeVector( myNodeIsOut );
4278 std::fill( myNodeIsChecked.begin(), myNodeIsChecked.end(), false );
4282 void ElementsOnShape::clearClassifiers()
4284 // for ( size_t i = 0; i < myClassifiers.size(); ++i )
4285 // delete myClassifiers[ i ];
4286 myClassifiers.clear();
4292 bool ElementsOnShape::IsSatisfy( long elemId )
4294 const SMDS_Mesh* mesh = myMeshModifTracer.GetMesh();
4295 const SMDS_MeshElement* elem =
4296 ( myType == SMDSAbs_Node ? mesh->FindNode( elemId ) : mesh->FindElement( elemId ));
4297 if ( !elem || myClassifiers.empty() )
4300 bool isSatisfy = myAllNodesFlag, isNodeOut;
4302 gp_XYZ centerXYZ (0, 0, 0);
4304 if ( !myOctree && myClassifiers.size() > 5 )
4306 myWorkClassifiers.resize( myClassifiers.size() );
4307 for ( size_t i = 0; i < myClassifiers.size(); ++i )
4308 myWorkClassifiers[ i ] = & myClassifiers[ i ];
4309 myOctree = new OctreeClassifier( myWorkClassifiers );
4312 SMDS_ElemIteratorPtr aNodeItr = elem->nodesIterator();
4313 while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
4315 SMESH_TNodeXYZ aPnt( aNodeItr->next() );
4319 if ( !getNodeIsOut( aPnt._node, isNodeOut ))
4323 myWorkClassifiers.clear();
4324 myOctree->GetClassifiersAtPoint( aPnt, myWorkClassifiers );
4326 for ( size_t i = 0; i < myWorkClassifiers.size(); ++i )
4327 myWorkClassifiers[i]->SetChecked( false );
4329 for ( size_t i = 0; i < myWorkClassifiers.size() && isNodeOut; ++i )
4330 if ( !myWorkClassifiers[i]->IsChecked() )
4331 isNodeOut = myWorkClassifiers[i]->IsOut( aPnt );
4335 for ( size_t i = 0; i < myClassifiers.size() && isNodeOut; ++i )
4336 isNodeOut = myClassifiers[i].IsOut( aPnt );
4338 setNodeIsOut( aPnt._node, isNodeOut );
4340 isSatisfy = !isNodeOut;
4343 // Check the center point for volumes MantisBug 0020168
4346 myClassifiers[0].ShapeType() == TopAbs_SOLID )
4348 centerXYZ /= elem->NbNodes();
4351 for ( size_t i = 0; i < myWorkClassifiers.size() && !isSatisfy; ++i )
4352 isSatisfy = ! myWorkClassifiers[i]->IsOut( centerXYZ );
4354 for ( size_t i = 0; i < myClassifiers.size() && !isSatisfy; ++i )
4355 isSatisfy = ! myClassifiers[i].IsOut( centerXYZ );
4361 void ElementsOnShape::Classifier::Init( const TopoDS_Shape& theShape,
4363 const Bnd_B3d* theBox )
4368 bool isShapeBox = false;
4369 switch ( myShape.ShapeType() )
4373 if (( isShapeBox = isBox( theShape )))
4375 myIsOutFun = & ElementsOnShape::Classifier::isOutOfBox;
4379 mySolidClfr.Load(theShape);
4380 myIsOutFun = & ElementsOnShape::Classifier::isOutOfSolid;
4386 Standard_Real u1,u2,v1,v2;
4387 Handle(Geom_Surface) surf = BRep_Tool::Surface( TopoDS::Face( theShape ));
4388 surf->Bounds( u1,u2,v1,v2 );
4389 myProjFace.Init(surf, u1,u2, v1,v2, myTol );
4390 myIsOutFun = & ElementsOnShape::Classifier::isOutOfFace;
4395 Standard_Real u1, u2;
4396 Handle(Geom_Curve) curve = BRep_Tool::Curve( TopoDS::Edge( theShape ), u1, u2);
4397 myProjEdge.Init(curve, u1, u2);
4398 myIsOutFun = & ElementsOnShape::Classifier::isOutOfEdge;
4403 myVertexXYZ = BRep_Tool::Pnt( TopoDS::Vertex( theShape ) );
4404 myIsOutFun = & ElementsOnShape::Classifier::isOutOfVertex;
4408 throw SALOME_Exception("Programmer error in usage of ElementsOnShape::Classifier");
4420 BRepBndLib::Add( myShape, box );
4422 myBox.Add( box.CornerMin() );
4423 myBox.Add( box.CornerMax() );
4424 gp_XYZ halfSize = 0.5 * ( box.CornerMax().XYZ() - box.CornerMin().XYZ() );
4425 for ( int iDim = 1; iDim <= 3; ++iDim )
4427 double x = halfSize.Coord( iDim );
4428 halfSize.SetCoord( iDim, x + Max( myTol, 1e-2 * x ));
4430 myBox.SetHSize( halfSize );
4435 bool ElementsOnShape::Classifier::isOutOfSolid (const gp_Pnt& p)
4437 mySolidClfr.Perform( p, myTol );
4438 return ( mySolidClfr.State() != TopAbs_IN && mySolidClfr.State() != TopAbs_ON );
4441 bool ElementsOnShape::Classifier::isOutOfBox (const gp_Pnt& p)
4443 return myBox.IsOut( p.XYZ() );
4446 bool ElementsOnShape::Classifier::isOutOfFace (const gp_Pnt& p)
4448 myProjFace.Perform( p );
4449 if ( myProjFace.IsDone() && myProjFace.LowerDistance() <= myTol )
4451 // check relatively to the face
4452 Quantity_Parameter u, v;
4453 myProjFace.LowerDistanceParameters(u, v);
4454 gp_Pnt2d aProjPnt (u, v);
4455 BRepClass_FaceClassifier aClsf ( TopoDS::Face( myShape ), aProjPnt, myTol );
4456 if ( aClsf.State() == TopAbs_IN || aClsf.State() == TopAbs_ON )
4462 bool ElementsOnShape::Classifier::isOutOfEdge (const gp_Pnt& p)
4464 myProjEdge.Perform( p );
4465 return ! ( myProjEdge.NbPoints() > 0 && myProjEdge.LowerDistance() <= myTol );
4468 bool ElementsOnShape::Classifier::isOutOfVertex(const gp_Pnt& p)
4470 return ( myVertexXYZ.Distance( p ) > myTol );
4473 bool ElementsOnShape::Classifier::isBox (const TopoDS_Shape& theShape)
4475 TopTools_IndexedMapOfShape vMap;
4476 TopExp::MapShapes( theShape, TopAbs_VERTEX, vMap );
4477 if ( vMap.Extent() != 8 )
4481 for ( int i = 1; i <= 8; ++i )
4482 myBox.Add( BRep_Tool::Pnt( TopoDS::Vertex( vMap( i ))).XYZ() );
4484 gp_XYZ pMin = myBox.CornerMin(), pMax = myBox.CornerMax();
4485 for ( int i = 1; i <= 8; ++i )
4487 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( vMap( i )));
4488 for ( int iC = 1; iC <= 3; ++ iC )
4490 double d1 = Abs( pMin.Coord( iC ) - p.Coord( iC ));
4491 double d2 = Abs( pMax.Coord( iC ) - p.Coord( iC ));
4492 if ( Min( d1, d2 ) > myTol )
4496 myBox.Enlarge( myTol );
4501 OctreeClassifier::OctreeClassifier( const std::vector< ElementsOnShape::Classifier* >& classifiers )
4502 :SMESH_Octree( new SMESH_TreeLimit )
4504 myClassifiers = classifiers;
4509 OctreeClassifier::OctreeClassifier( const OctreeClassifier* otherTree,
4510 const std::vector< ElementsOnShape::Classifier >& clsOther,
4511 std::vector< ElementsOnShape::Classifier >& cls )
4512 :SMESH_Octree( new SMESH_TreeLimit )
4514 myBox = new Bnd_B3d( *otherTree->getBox() );
4516 if (( myIsLeaf = otherTree->isLeaf() ))
4518 myClassifiers.resize( otherTree->myClassifiers.size() );
4519 for ( size_t i = 0; i < otherTree->myClassifiers.size(); ++i )
4521 int ind = otherTree->myClassifiers[i] - & clsOther[0];
4522 myClassifiers[ i ] = & cls[ ind ];
4525 else if ( otherTree->myChildren )
4527 myChildren = new SMESH_Tree < Bnd_B3d, 8 >*[ 8 ];
4528 for ( int i = 0; i < nbChildren(); i++ )
4530 new OctreeClassifier( static_cast<const OctreeClassifier*>( otherTree->myChildren[i]),
4535 void ElementsOnShape::
4536 OctreeClassifier::GetClassifiersAtPoint( const gp_XYZ& point,
4537 std::vector< ElementsOnShape::Classifier* >& result )
4539 if ( getBox()->IsOut( point ))
4544 for ( size_t i = 0; i < myClassifiers.size(); ++i )
4545 if ( !myClassifiers[i]->GetBndBox()->IsOut( point ))
4546 result.push_back( myClassifiers[i] );
4550 for (int i = 0; i < nbChildren(); i++)
4551 ((OctreeClassifier*) myChildren[i])->GetClassifiersAtPoint( point, result );
4555 void ElementsOnShape::OctreeClassifier::buildChildrenData()
4557 // distribute myClassifiers among myChildren
4559 const int childFlag[8] = { 0x0000001,
4567 int nbInChild[8] = { 0, 0, 0, 0, 0, 0, 0, 0 };
4569 for ( size_t i = 0; i < myClassifiers.size(); ++i )
4571 for ( int j = 0; j < nbChildren(); j++ )
4573 if ( !myClassifiers[i]->GetBndBox()->IsOut( *myChildren[j]->getBox() ))
4575 myClassifiers[i]->SetFlag( childFlag[ j ]);
4581 for ( int j = 0; j < nbChildren(); j++ )
4583 OctreeClassifier* child = static_cast<OctreeClassifier*>( myChildren[ j ]);
4584 child->myClassifiers.resize( nbInChild[ j ]);
4585 for ( size_t i = 0; nbInChild[ j ] && i < myClassifiers.size(); ++i )
4587 if ( myClassifiers[ i ]->IsSetFlag( childFlag[ j ]))
4590 child->myClassifiers[ nbInChild[ j ]] = myClassifiers[ i ];
4591 myClassifiers[ i ]->UnsetFlag( childFlag[ j ]);
4595 SMESHUtils::FreeVector( myClassifiers );
4597 // define if a child isLeaf()
4598 for ( int i = 0; i < nbChildren(); i++ )
4600 OctreeClassifier* child = static_cast<OctreeClassifier*>( myChildren[ i ]);
4601 child->myIsLeaf = ( child->myClassifiers.size() <= 5 );
4605 Bnd_B3d* ElementsOnShape::OctreeClassifier::buildRootBox()
4607 Bnd_B3d* box = new Bnd_B3d;
4608 for ( size_t i = 0; i < myClassifiers.size(); ++i )
4609 box->Add( *myClassifiers[i]->GetBndBox() );
4614 Class : BelongToGeom
4615 Description : Predicate for verifying whether entity belongs to
4616 specified geometrical support
4619 BelongToGeom::BelongToGeom()
4621 myType(SMDSAbs_NbElementTypes),
4622 myIsSubshape(false),
4623 myTolerance(Precision::Confusion())
4626 Predicate* BelongToGeom::clone() const
4628 BelongToGeom* cln = new BelongToGeom( *this );
4629 cln->myElementsOnShapePtr.reset( static_cast<ElementsOnShape*>( myElementsOnShapePtr->clone() ));
4633 void BelongToGeom::SetMesh( const SMDS_Mesh* theMesh )
4635 if ( myMeshDS != theMesh )
4637 myMeshDS = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
4642 void BelongToGeom::SetGeom( const TopoDS_Shape& theShape )
4644 if ( myShape != theShape )
4651 static bool IsSubShape (const TopTools_IndexedMapOfShape& theMap,
4652 const TopoDS_Shape& theShape)
4654 if (theMap.Contains(theShape)) return true;
4656 if (theShape.ShapeType() == TopAbs_COMPOUND ||
4657 theShape.ShapeType() == TopAbs_COMPSOLID)
4659 TopoDS_Iterator anIt (theShape, Standard_True, Standard_True);
4660 for (; anIt.More(); anIt.Next())
4662 if (!IsSubShape(theMap, anIt.Value())) {
4672 void BelongToGeom::init()
4674 if ( !myMeshDS || myShape.IsNull() ) return;
4676 // is sub-shape of main shape?
4677 TopoDS_Shape aMainShape = myMeshDS->ShapeToMesh();
4678 if (aMainShape.IsNull()) {
4679 myIsSubshape = false;
4682 TopTools_IndexedMapOfShape aMap;
4683 TopExp::MapShapes( aMainShape, aMap );
4684 myIsSubshape = IsSubShape( aMap, myShape );
4688 TopExp::MapShapes( myShape, aMap );
4689 mySubShapesIDs.Clear();
4690 for ( int i = 1; i <= aMap.Extent(); ++i )
4692 int subID = myMeshDS->ShapeToIndex( aMap( i ));
4694 mySubShapesIDs.Add( subID );
4699 //if (!myIsSubshape) // to be always ready to check an element not bound to geometry
4701 if ( !myElementsOnShapePtr )
4702 myElementsOnShapePtr.reset( new ElementsOnShape() );
4703 myElementsOnShapePtr->SetTolerance( myTolerance );
4704 myElementsOnShapePtr->SetAllNodes( true ); // "belong", while false means "lays on"
4705 myElementsOnShapePtr->SetMesh( myMeshDS );
4706 myElementsOnShapePtr->SetShape( myShape, myType );
4710 bool BelongToGeom::IsSatisfy (long theId)
4712 if (myMeshDS == 0 || myShape.IsNull())
4717 return myElementsOnShapePtr->IsSatisfy(theId);
4722 if (myType == SMDSAbs_Node)
4724 if ( const SMDS_MeshNode* aNode = myMeshDS->FindNode( theId ))
4726 if ( aNode->getshapeId() < 1 )
4727 return myElementsOnShapePtr->IsSatisfy(theId);
4729 return mySubShapesIDs.Contains( aNode->getshapeId() );
4734 if ( const SMDS_MeshElement* anElem = myMeshDS->FindElement( theId ))
4736 if ( anElem->GetType() == myType )
4738 if ( anElem->getshapeId() < 1 )
4739 return myElementsOnShapePtr->IsSatisfy(theId);
4741 return mySubShapesIDs.Contains( anElem->getshapeId() );
4749 void BelongToGeom::SetType (SMDSAbs_ElementType theType)
4751 if ( myType != theType )
4758 SMDSAbs_ElementType BelongToGeom::GetType() const
4763 TopoDS_Shape BelongToGeom::GetShape()
4768 const SMESHDS_Mesh* BelongToGeom::GetMeshDS() const
4773 void BelongToGeom::SetTolerance (double theTolerance)
4775 myTolerance = theTolerance;
4779 double BelongToGeom::GetTolerance()
4786 Description : Predicate for verifying whether entiy lying or partially lying on
4787 specified geometrical support
4790 LyingOnGeom::LyingOnGeom()
4792 myType(SMDSAbs_NbElementTypes),
4793 myIsSubshape(false),
4794 myTolerance(Precision::Confusion())
4797 Predicate* LyingOnGeom::clone() const
4799 LyingOnGeom* cln = new LyingOnGeom( *this );
4800 cln->myElementsOnShapePtr.reset( static_cast<ElementsOnShape*>( myElementsOnShapePtr->clone() ));
4804 void LyingOnGeom::SetMesh( const SMDS_Mesh* theMesh )
4806 if ( myMeshDS != theMesh )
4808 myMeshDS = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
4813 void LyingOnGeom::SetGeom( const TopoDS_Shape& theShape )
4815 if ( myShape != theShape )
4822 void LyingOnGeom::init()
4824 if (!myMeshDS || myShape.IsNull()) return;
4826 // is sub-shape of main shape?
4827 TopoDS_Shape aMainShape = myMeshDS->ShapeToMesh();
4828 if (aMainShape.IsNull()) {
4829 myIsSubshape = false;
4832 myIsSubshape = myMeshDS->IsGroupOfSubShapes( myShape );
4837 TopTools_IndexedMapOfShape shapes;
4838 TopExp::MapShapes( myShape, shapes );
4839 mySubShapesIDs.Clear();
4840 for ( int i = 1; i <= shapes.Extent(); ++i )
4842 int subID = myMeshDS->ShapeToIndex( shapes( i ));
4844 mySubShapesIDs.Add( subID );
4847 // else // to be always ready to check an element not bound to geometry
4849 if ( !myElementsOnShapePtr )
4850 myElementsOnShapePtr.reset( new ElementsOnShape() );
4851 myElementsOnShapePtr->SetTolerance( myTolerance );
4852 myElementsOnShapePtr->SetAllNodes( false ); // lays on, while true means "belong"
4853 myElementsOnShapePtr->SetMesh( myMeshDS );
4854 myElementsOnShapePtr->SetShape( myShape, myType );
4858 bool LyingOnGeom::IsSatisfy( long theId )
4860 if ( myMeshDS == 0 || myShape.IsNull() )
4865 return myElementsOnShapePtr->IsSatisfy(theId);
4870 const SMDS_MeshElement* elem =
4871 ( myType == SMDSAbs_Node ) ? myMeshDS->FindNode( theId ) : myMeshDS->FindElement( theId );
4873 if ( mySubShapesIDs.Contains( elem->getshapeId() ))
4876 if ( elem->GetType() != SMDSAbs_Node && elem->GetType() == myType )
4878 SMDS_ElemIteratorPtr nodeItr = elem->nodesIterator();
4879 while ( nodeItr->more() )
4881 const SMDS_MeshElement* aNode = nodeItr->next();
4882 if ( mySubShapesIDs.Contains( aNode->getshapeId() ))
4890 void LyingOnGeom::SetType( SMDSAbs_ElementType theType )
4892 if ( myType != theType )
4899 SMDSAbs_ElementType LyingOnGeom::GetType() const
4904 TopoDS_Shape LyingOnGeom::GetShape()
4909 const SMESHDS_Mesh* LyingOnGeom::GetMeshDS() const
4914 void LyingOnGeom::SetTolerance (double theTolerance)
4916 myTolerance = theTolerance;
4920 double LyingOnGeom::GetTolerance()
4925 TSequenceOfXYZ::TSequenceOfXYZ(): myElem(0)
4928 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n), myElem(0)
4931 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t), myElem(0)
4934 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray), myElem(theSequenceOfXYZ.myElem)
4937 template <class InputIterator>
4938 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd), myElem(0)
4941 TSequenceOfXYZ::~TSequenceOfXYZ()
4944 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
4946 myArray = theSequenceOfXYZ.myArray;
4947 myElem = theSequenceOfXYZ.myElem;
4951 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
4953 return myArray[n-1];
4956 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
4958 return myArray[n-1];
4961 void TSequenceOfXYZ::clear()
4966 void TSequenceOfXYZ::reserve(size_type n)
4971 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
4973 myArray.push_back(v);
4976 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
4978 return myArray.size();
4981 SMDSAbs_EntityType TSequenceOfXYZ::getElementEntity() const
4983 return myElem ? myElem->GetEntityType() : SMDSEntity_Last;
4986 TMeshModifTracer::TMeshModifTracer():
4987 myMeshModifTime(0), myMesh(0)
4990 void TMeshModifTracer::SetMesh( const SMDS_Mesh* theMesh )
4992 if ( theMesh != myMesh )
4993 myMeshModifTime = 0;
4996 bool TMeshModifTracer::IsMeshModified()
4998 bool modified = false;
5001 modified = ( myMeshModifTime != myMesh->GetMTime() );
5002 myMeshModifTime = myMesh->GetMTime();