1 // Copyright (C) 2007-2013 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.
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_Mesh.hxx"
35 #include "SMESH_OctreeNode.hxx"
37 #include <BRepAdaptor_Surface.hxx>
38 #include <BRepClass_FaceClassifier.hxx>
39 #include <BRep_Tool.hxx>
40 #include <Geom_CylindricalSurface.hxx>
41 #include <Geom_Plane.hxx>
42 #include <Geom_Surface.hxx>
43 #include <Precision.hxx>
44 #include <TColStd_MapIteratorOfMapOfInteger.hxx>
45 #include <TColStd_MapOfInteger.hxx>
46 #include <TColStd_SequenceOfAsciiString.hxx>
47 #include <TColgp_Array1OfXYZ.hxx>
50 #include <TopoDS_Edge.hxx>
51 #include <TopoDS_Face.hxx>
52 #include <TopoDS_Iterator.hxx>
53 #include <TopoDS_Shape.hxx>
54 #include <TopoDS_Vertex.hxx>
56 #include <gp_Cylinder.hxx>
63 #include <vtkMeshQuality.h>
75 const double theEps = 1e-100;
76 const double theInf = 1e+100;
78 inline gp_XYZ gpXYZ(const SMDS_MeshNode* aNode )
80 return gp_XYZ(aNode->X(), aNode->Y(), aNode->Z() );
83 inline double getAngle( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
85 gp_Vec v1( P1 - P2 ), v2( P3 - P2 );
87 return v1.Magnitude() < gp::Resolution() ||
88 v2.Magnitude() < gp::Resolution() ? 0 : v1.Angle( v2 );
91 inline double getArea( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
93 gp_Vec aVec1( P2 - P1 );
94 gp_Vec aVec2( P3 - P1 );
95 return ( aVec1 ^ aVec2 ).Magnitude() * 0.5;
98 inline double getArea( const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3 )
100 return getArea( P1.XYZ(), P2.XYZ(), P3.XYZ() );
105 inline double getDistance( const gp_XYZ& P1, const gp_XYZ& P2 )
107 double aDist = gp_Pnt( P1 ).Distance( gp_Pnt( P2 ) );
111 int getNbMultiConnection( const SMDS_Mesh* theMesh, const int theId )
116 const SMDS_MeshElement* anEdge = theMesh->FindElement( theId );
117 if ( anEdge == 0 || anEdge->GetType() != SMDSAbs_Edge/* || anEdge->NbNodes() != 2 */)
120 // for each pair of nodes in anEdge (there are 2 pairs in a quadratic edge)
121 // count elements containing both nodes of the pair.
122 // Note that there may be such cases for a quadratic edge (a horizontal line):
127 // +-----+------+ +-----+------+
130 // result sould be 2 in both cases
132 int aResult0 = 0, aResult1 = 0;
133 // last node, it is a medium one in a quadratic edge
134 const SMDS_MeshNode* aLastNode = anEdge->GetNode( anEdge->NbNodes() - 1 );
135 const SMDS_MeshNode* aNode0 = anEdge->GetNode( 0 );
136 const SMDS_MeshNode* aNode1 = anEdge->GetNode( 1 );
137 if ( aNode1 == aLastNode ) aNode1 = 0;
139 SMDS_ElemIteratorPtr anElemIter = aLastNode->GetInverseElementIterator();
140 while( anElemIter->more() ) {
141 const SMDS_MeshElement* anElem = anElemIter->next();
142 if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
143 SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
144 while ( anIter->more() ) {
145 if ( const SMDS_MeshElement* anElemNode = anIter->next() ) {
146 if ( anElemNode == aNode0 ) {
148 if ( !aNode1 ) break; // not a quadratic edge
150 else if ( anElemNode == aNode1 )
156 int aResult = std::max ( aResult0, aResult1 );
158 // TColStd_MapOfInteger aMap;
160 // SMDS_ElemIteratorPtr anIter = anEdge->nodesIterator();
161 // if ( anIter != 0 ) {
162 // while( anIter->more() ) {
163 // const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
166 // SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
167 // while( anElemIter->more() ) {
168 // const SMDS_MeshElement* anElem = anElemIter->next();
169 // if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
170 // int anId = anElem->GetID();
172 // if ( anIter->more() ) // i.e. first node
174 // else if ( aMap.Contains( anId ) )
184 gp_XYZ getNormale( const SMDS_MeshFace* theFace, bool* ok=0 )
186 int aNbNode = theFace->NbNodes();
188 gp_XYZ q1 = gpXYZ( theFace->GetNode(1)) - gpXYZ( theFace->GetNode(0));
189 gp_XYZ q2 = gpXYZ( theFace->GetNode(2)) - gpXYZ( theFace->GetNode(0));
192 gp_XYZ q3 = gpXYZ( theFace->GetNode(3)) - gpXYZ( theFace->GetNode(0));
195 double len = n.Modulus();
196 bool zeroLen = ( len <= numeric_limits<double>::min());
200 if (ok) *ok = !zeroLen;
208 using namespace SMESH::Controls;
214 //================================================================================
216 Class : NumericalFunctor
217 Description : Base class for numerical functors
219 //================================================================================
221 NumericalFunctor::NumericalFunctor():
227 void NumericalFunctor::SetMesh( const SMDS_Mesh* theMesh )
232 bool NumericalFunctor::GetPoints(const int theId,
233 TSequenceOfXYZ& theRes ) const
240 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
241 if ( !anElem || anElem->GetType() != this->GetType() )
244 return GetPoints( anElem, theRes );
247 bool NumericalFunctor::GetPoints(const SMDS_MeshElement* anElem,
248 TSequenceOfXYZ& theRes )
255 theRes.reserve( anElem->NbNodes() );
257 // Get nodes of the element
258 SMDS_ElemIteratorPtr anIter;
260 if ( anElem->IsQuadratic() ) {
261 switch ( anElem->GetType() ) {
263 anIter = dynamic_cast<const SMDS_VtkEdge*>
264 (anElem)->interlacedNodesElemIterator();
267 anIter = dynamic_cast<const SMDS_VtkFace*>
268 (anElem)->interlacedNodesElemIterator();
271 anIter = anElem->nodesIterator();
276 anIter = anElem->nodesIterator();
280 while( anIter->more() ) {
281 if ( const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>( anIter->next() ))
282 theRes.push_back( gp_XYZ( aNode->X(), aNode->Y(), aNode->Z() ) );
289 long NumericalFunctor::GetPrecision() const
294 void NumericalFunctor::SetPrecision( const long thePrecision )
296 myPrecision = thePrecision;
297 myPrecisionValue = pow( 10., (double)( myPrecision ) );
300 double NumericalFunctor::GetValue( long theId )
304 myCurrElement = myMesh->FindElement( theId );
307 if ( GetPoints( theId, P ))
308 aVal = Round( GetValue( P ));
313 double NumericalFunctor::Round( const double & aVal )
315 return ( myPrecision >= 0 ) ? floor( aVal * myPrecisionValue + 0.5 ) / myPrecisionValue : aVal;
318 //================================================================================
320 * \brief Return histogram of functor values
321 * \param nbIntervals - number of intervals
322 * \param nbEvents - number of mesh elements having values within i-th interval
323 * \param funValues - boundaries of intervals
324 * \param elements - elements to check vulue of; empty list means "of all"
325 * \param minmax - boundaries of diapason of values to divide into intervals
327 //================================================================================
329 void NumericalFunctor::GetHistogram(int nbIntervals,
330 std::vector<int>& nbEvents,
331 std::vector<double>& funValues,
332 const vector<int>& elements,
333 const double* minmax,
334 const bool isLogarithmic)
336 if ( nbIntervals < 1 ||
338 !myMesh->GetMeshInfo().NbElements( GetType() ))
340 nbEvents.resize( nbIntervals, 0 );
341 funValues.resize( nbIntervals+1 );
343 // get all values sorted
344 std::multiset< double > values;
345 if ( elements.empty() )
347 SMDS_ElemIteratorPtr elemIt = myMesh->elementsIterator(GetType());
348 while ( elemIt->more() )
349 values.insert( GetValue( elemIt->next()->GetID() ));
353 vector<int>::const_iterator id = elements.begin();
354 for ( ; id != elements.end(); ++id )
355 values.insert( GetValue( *id ));
360 funValues[0] = minmax[0];
361 funValues[nbIntervals] = minmax[1];
365 funValues[0] = *values.begin();
366 funValues[nbIntervals] = *values.rbegin();
368 // case nbIntervals == 1
369 if ( nbIntervals == 1 )
371 nbEvents[0] = values.size();
375 if (funValues.front() == funValues.back())
377 nbEvents.resize( 1 );
378 nbEvents[0] = values.size();
379 funValues[1] = funValues.back();
380 funValues.resize( 2 );
383 std::multiset< double >::iterator min = values.begin(), max;
384 for ( int i = 0; i < nbIntervals; ++i )
386 // find end value of i-th interval
387 double r = (i+1) / double(nbIntervals);
388 if (isLogarithmic && funValues.front() > 1e-07 && funValues.back() > 1e-07) {
389 double logmin = log10(funValues.front());
390 double lval = logmin + r * (log10(funValues.back()) - logmin);
391 funValues[i+1] = pow(10.0, lval);
394 funValues[i+1] = funValues.front() * (1-r) + funValues.back() * r;
397 // count values in the i-th interval if there are any
398 if ( min != values.end() && *min <= funValues[i+1] )
400 // find the first value out of the interval
401 max = values.upper_bound( funValues[i+1] ); // max is greater than funValues[i+1], or end()
402 nbEvents[i] = std::distance( min, max );
406 // add values larger than minmax[1]
407 nbEvents.back() += std::distance( min, values.end() );
410 //=======================================================================
413 Description : Functor calculating volume of a 3D element
415 //================================================================================
417 double Volume::GetValue( long theElementId )
419 if ( theElementId && myMesh ) {
420 SMDS_VolumeTool aVolumeTool;
421 if ( aVolumeTool.Set( myMesh->FindElement( theElementId )))
422 return aVolumeTool.GetSize();
427 double Volume::GetBadRate( double Value, int /*nbNodes*/ ) const
432 SMDSAbs_ElementType Volume::GetType() const
434 return SMDSAbs_Volume;
437 //=======================================================================
439 Class : MaxElementLength2D
440 Description : Functor calculating maximum length of 2D element
442 //================================================================================
444 double MaxElementLength2D::GetValue( const TSequenceOfXYZ& P )
450 if( len == 3 ) { // triangles
451 double L1 = getDistance(P( 1 ),P( 2 ));
452 double L2 = getDistance(P( 2 ),P( 3 ));
453 double L3 = getDistance(P( 3 ),P( 1 ));
454 aVal = Max(L1,Max(L2,L3));
456 else if( len == 4 ) { // quadrangles
457 double L1 = getDistance(P( 1 ),P( 2 ));
458 double L2 = getDistance(P( 2 ),P( 3 ));
459 double L3 = getDistance(P( 3 ),P( 4 ));
460 double L4 = getDistance(P( 4 ),P( 1 ));
461 double D1 = getDistance(P( 1 ),P( 3 ));
462 double D2 = getDistance(P( 2 ),P( 4 ));
463 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
465 else if( len == 6 ) { // quadratic triangles
466 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
467 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
468 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
469 aVal = Max(L1,Max(L2,L3));
471 else if( len == 8 || len == 9 ) { // quadratic quadrangles
472 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
473 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
474 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
475 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
476 double D1 = getDistance(P( 1 ),P( 5 ));
477 double D2 = getDistance(P( 3 ),P( 7 ));
478 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
481 if( myPrecision >= 0 )
483 double prec = pow( 10., (double)myPrecision );
484 aVal = floor( aVal * prec + 0.5 ) / prec;
489 double MaxElementLength2D::GetValue( long theElementId )
492 return GetPoints( theElementId, P ) ? GetValue(P) : 0.0;
495 double MaxElementLength2D::GetBadRate( double Value, int /*nbNodes*/ ) const
500 SMDSAbs_ElementType MaxElementLength2D::GetType() const
505 //=======================================================================
507 Class : MaxElementLength3D
508 Description : Functor calculating maximum length of 3D element
510 //================================================================================
512 double MaxElementLength3D::GetValue( long theElementId )
515 if( GetPoints( theElementId, P ) ) {
517 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
518 SMDSAbs_ElementType aType = aElem->GetType();
522 if( len == 4 ) { // tetras
523 double L1 = getDistance(P( 1 ),P( 2 ));
524 double L2 = getDistance(P( 2 ),P( 3 ));
525 double L3 = getDistance(P( 3 ),P( 1 ));
526 double L4 = getDistance(P( 1 ),P( 4 ));
527 double L5 = getDistance(P( 2 ),P( 4 ));
528 double L6 = getDistance(P( 3 ),P( 4 ));
529 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
532 else if( len == 5 ) { // pyramids
533 double L1 = getDistance(P( 1 ),P( 2 ));
534 double L2 = getDistance(P( 2 ),P( 3 ));
535 double L3 = getDistance(P( 3 ),P( 4 ));
536 double L4 = getDistance(P( 4 ),P( 1 ));
537 double L5 = getDistance(P( 1 ),P( 5 ));
538 double L6 = getDistance(P( 2 ),P( 5 ));
539 double L7 = getDistance(P( 3 ),P( 5 ));
540 double L8 = getDistance(P( 4 ),P( 5 ));
541 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
542 aVal = Max(aVal,Max(L7,L8));
545 else if( len == 6 ) { // pentas
546 double L1 = getDistance(P( 1 ),P( 2 ));
547 double L2 = getDistance(P( 2 ),P( 3 ));
548 double L3 = getDistance(P( 3 ),P( 1 ));
549 double L4 = getDistance(P( 4 ),P( 5 ));
550 double L5 = getDistance(P( 5 ),P( 6 ));
551 double L6 = getDistance(P( 6 ),P( 4 ));
552 double L7 = getDistance(P( 1 ),P( 4 ));
553 double L8 = getDistance(P( 2 ),P( 5 ));
554 double L9 = getDistance(P( 3 ),P( 6 ));
555 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
556 aVal = Max(aVal,Max(Max(L7,L8),L9));
559 else if( len == 8 ) { // hexas
560 double L1 = getDistance(P( 1 ),P( 2 ));
561 double L2 = getDistance(P( 2 ),P( 3 ));
562 double L3 = getDistance(P( 3 ),P( 4 ));
563 double L4 = getDistance(P( 4 ),P( 1 ));
564 double L5 = getDistance(P( 5 ),P( 6 ));
565 double L6 = getDistance(P( 6 ),P( 7 ));
566 double L7 = getDistance(P( 7 ),P( 8 ));
567 double L8 = getDistance(P( 8 ),P( 5 ));
568 double L9 = getDistance(P( 1 ),P( 5 ));
569 double L10= getDistance(P( 2 ),P( 6 ));
570 double L11= getDistance(P( 3 ),P( 7 ));
571 double L12= getDistance(P( 4 ),P( 8 ));
572 double D1 = getDistance(P( 1 ),P( 7 ));
573 double D2 = getDistance(P( 2 ),P( 8 ));
574 double D3 = getDistance(P( 3 ),P( 5 ));
575 double D4 = getDistance(P( 4 ),P( 6 ));
576 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
577 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
578 aVal = Max(aVal,Max(L11,L12));
579 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
582 else if( len == 12 ) { // hexagonal prism
583 for ( int i1 = 1; i1 < 12; ++i1 )
584 for ( int i2 = i1+1; i1 <= 12; ++i1 )
585 aVal = Max( aVal, getDistance(P( i1 ),P( i2 )));
588 else if( len == 10 ) { // quadratic tetras
589 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
590 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
591 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
592 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
593 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
594 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
595 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
598 else if( len == 13 ) { // quadratic pyramids
599 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
600 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
601 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
602 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
603 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
604 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
605 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
606 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
607 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
608 aVal = Max(aVal,Max(L7,L8));
611 else if( len == 15 ) { // quadratic pentas
612 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
613 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
614 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
615 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
616 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
617 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
618 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
619 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
620 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
621 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
622 aVal = Max(aVal,Max(Max(L7,L8),L9));
625 else if( len == 20 || len == 27 ) { // quadratic hexas
626 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
627 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
628 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
629 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
630 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
631 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
632 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
633 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
634 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
635 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
636 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
637 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
638 double D1 = getDistance(P( 1 ),P( 7 ));
639 double D2 = getDistance(P( 2 ),P( 8 ));
640 double D3 = getDistance(P( 3 ),P( 5 ));
641 double D4 = getDistance(P( 4 ),P( 6 ));
642 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
643 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
644 aVal = Max(aVal,Max(L11,L12));
645 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
648 else if( len > 1 && aElem->IsPoly() ) { // polys
649 // get the maximum distance between all pairs of nodes
650 for( int i = 1; i <= len; i++ ) {
651 for( int j = 1; j <= len; j++ ) {
652 if( j > i ) { // optimization of the loop
653 double D = getDistance( P(i), P(j) );
654 aVal = Max( aVal, D );
661 if( myPrecision >= 0 )
663 double prec = pow( 10., (double)myPrecision );
664 aVal = floor( aVal * prec + 0.5 ) / prec;
671 double MaxElementLength3D::GetBadRate( double Value, int /*nbNodes*/ ) const
676 SMDSAbs_ElementType MaxElementLength3D::GetType() const
678 return SMDSAbs_Volume;
681 //=======================================================================
684 Description : Functor for calculation of minimum angle
686 //================================================================================
688 double MinimumAngle::GetValue( const TSequenceOfXYZ& P )
695 aMin = getAngle(P( P.size() ), P( 1 ), P( 2 ));
696 aMin = Min(aMin,getAngle(P( P.size()-1 ), P( P.size() ), P( 1 )));
698 for (int i=2; i<P.size();i++){
699 double A0 = getAngle( P( i-1 ), P( i ), P( i+1 ) );
703 return aMin * 180.0 / M_PI;
706 double MinimumAngle::GetBadRate( double Value, int nbNodes ) const
708 //const double aBestAngle = PI / nbNodes;
709 const double aBestAngle = 180.0 - ( 360.0 / double(nbNodes) );
710 return ( fabs( aBestAngle - Value ));
713 SMDSAbs_ElementType MinimumAngle::GetType() const
719 //================================================================================
722 Description : Functor for calculating aspect ratio
724 //================================================================================
726 double AspectRatio::GetValue( long theId )
729 myCurrElement = myMesh->FindElement( theId );
730 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_QUAD )
733 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
734 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
735 aVal = Round( vtkMeshQuality::QuadAspectRatio( avtkCell ));
740 if ( GetPoints( myCurrElement, P ))
741 aVal = Round( GetValue( P ));
746 double AspectRatio::GetValue( const TSequenceOfXYZ& P )
748 // According to "Mesh quality control" by Nadir Bouhamau referring to
749 // Pascal Jean Frey and Paul-Louis George. Maillages, applications aux elements finis.
750 // Hermes Science publications, Paris 1999 ISBN 2-7462-0024-4
753 int nbNodes = P.size();
758 // Compute aspect ratio
760 if ( nbNodes == 3 ) {
761 // Compute lengths of the sides
762 std::vector< double > aLen (nbNodes);
763 for ( int i = 0; i < nbNodes - 1; i++ )
764 aLen[ i ] = getDistance( P( i + 1 ), P( i + 2 ) );
765 aLen[ nbNodes - 1 ] = getDistance( P( 1 ), P( nbNodes ) );
766 // Q = alfa * h * p / S, where
768 // alfa = sqrt( 3 ) / 6
769 // h - length of the longest edge
770 // p - half perimeter
771 // S - triangle surface
772 const double alfa = sqrt( 3. ) / 6.;
773 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
774 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
775 double anArea = getArea( P( 1 ), P( 2 ), P( 3 ) );
776 if ( anArea <= theEps )
778 return alfa * maxLen * half_perimeter / anArea;
780 else if ( nbNodes == 6 ) { // quadratic triangles
781 // Compute lengths of the sides
782 std::vector< double > aLen (3);
783 aLen[0] = getDistance( P(1), P(3) );
784 aLen[1] = getDistance( P(3), P(5) );
785 aLen[2] = getDistance( P(5), P(1) );
786 // Q = alfa * h * p / S, where
788 // alfa = sqrt( 3 ) / 6
789 // h - length of the longest edge
790 // p - half perimeter
791 // S - triangle surface
792 const double alfa = sqrt( 3. ) / 6.;
793 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
794 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
795 double anArea = getArea( P(1), P(3), P(5) );
796 if ( anArea <= theEps )
798 return alfa * maxLen * half_perimeter / anArea;
800 else if( nbNodes == 4 ) { // quadrangle
801 // Compute lengths of the sides
802 std::vector< double > aLen (4);
803 aLen[0] = getDistance( P(1), P(2) );
804 aLen[1] = getDistance( P(2), P(3) );
805 aLen[2] = getDistance( P(3), P(4) );
806 aLen[3] = getDistance( P(4), P(1) );
807 // Compute lengths of the diagonals
808 std::vector< double > aDia (2);
809 aDia[0] = getDistance( P(1), P(3) );
810 aDia[1] = getDistance( P(2), P(4) );
811 // Compute areas of all triangles which can be built
812 // taking three nodes of the quadrangle
813 std::vector< double > anArea (4);
814 anArea[0] = getArea( P(1), P(2), P(3) );
815 anArea[1] = getArea( P(1), P(2), P(4) );
816 anArea[2] = getArea( P(1), P(3), P(4) );
817 anArea[3] = getArea( P(2), P(3), P(4) );
818 // Q = alpha * L * C1 / C2, where
820 // alpha = sqrt( 1/32 )
821 // L = max( L1, L2, L3, L4, D1, D2 )
822 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
823 // C2 = min( S1, S2, S3, S4 )
824 // Li - lengths of the edges
825 // Di - lengths of the diagonals
826 // Si - areas of the triangles
827 const double alpha = sqrt( 1 / 32. );
828 double L = Max( aLen[ 0 ],
832 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
833 double C1 = sqrt( ( aLen[0] * aLen[0] +
836 aLen[3] * aLen[3] ) / 4. );
837 double C2 = Min( anArea[ 0 ],
839 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
842 return alpha * L * C1 / C2;
844 else if( nbNodes == 8 || nbNodes == 9 ) { // nbNodes==8 - quadratic quadrangle
845 // Compute lengths of the sides
846 std::vector< double > aLen (4);
847 aLen[0] = getDistance( P(1), P(3) );
848 aLen[1] = getDistance( P(3), P(5) );
849 aLen[2] = getDistance( P(5), P(7) );
850 aLen[3] = getDistance( P(7), P(1) );
851 // Compute lengths of the diagonals
852 std::vector< double > aDia (2);
853 aDia[0] = getDistance( P(1), P(5) );
854 aDia[1] = getDistance( P(3), P(7) );
855 // Compute areas of all triangles which can be built
856 // taking three nodes of the quadrangle
857 std::vector< double > anArea (4);
858 anArea[0] = getArea( P(1), P(3), P(5) );
859 anArea[1] = getArea( P(1), P(3), P(7) );
860 anArea[2] = getArea( P(1), P(5), P(7) );
861 anArea[3] = getArea( P(3), P(5), P(7) );
862 // Q = alpha * L * C1 / C2, where
864 // alpha = sqrt( 1/32 )
865 // L = max( L1, L2, L3, L4, D1, D2 )
866 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
867 // C2 = min( S1, S2, S3, S4 )
868 // Li - lengths of the edges
869 // Di - lengths of the diagonals
870 // Si - areas of the triangles
871 const double alpha = sqrt( 1 / 32. );
872 double L = Max( aLen[ 0 ],
876 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
877 double C1 = sqrt( ( aLen[0] * aLen[0] +
880 aLen[3] * aLen[3] ) / 4. );
881 double C2 = Min( anArea[ 0 ],
883 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
886 return alpha * L * C1 / C2;
891 double AspectRatio::GetBadRate( double Value, int /*nbNodes*/ ) const
893 // the aspect ratio is in the range [1.0,infinity]
894 // < 1.0 = very bad, zero area
897 return ( Value < 0.9 ) ? 1000 : Value / 1000.;
900 SMDSAbs_ElementType AspectRatio::GetType() const
906 //================================================================================
908 Class : AspectRatio3D
909 Description : Functor for calculating aspect ratio
911 //================================================================================
915 inline double getHalfPerimeter(double theTria[3]){
916 return (theTria[0] + theTria[1] + theTria[2])/2.0;
919 inline double getArea(double theHalfPerim, double theTria[3]){
920 return sqrt(theHalfPerim*
921 (theHalfPerim-theTria[0])*
922 (theHalfPerim-theTria[1])*
923 (theHalfPerim-theTria[2]));
926 inline double getVolume(double theLen[6]){
927 double a2 = theLen[0]*theLen[0];
928 double b2 = theLen[1]*theLen[1];
929 double c2 = theLen[2]*theLen[2];
930 double d2 = theLen[3]*theLen[3];
931 double e2 = theLen[4]*theLen[4];
932 double f2 = theLen[5]*theLen[5];
933 double P = 4.0*a2*b2*d2;
934 double Q = a2*(b2+d2-e2)-b2*(a2+d2-f2)-d2*(a2+b2-c2);
935 double R = (b2+d2-e2)*(a2+d2-f2)*(a2+d2-f2);
936 return sqrt(P-Q+R)/12.0;
939 inline double getVolume2(double theLen[6]){
940 double a2 = theLen[0]*theLen[0];
941 double b2 = theLen[1]*theLen[1];
942 double c2 = theLen[2]*theLen[2];
943 double d2 = theLen[3]*theLen[3];
944 double e2 = theLen[4]*theLen[4];
945 double f2 = theLen[5]*theLen[5];
947 double P = a2*e2*(b2+c2+d2+f2-a2-e2);
948 double Q = b2*f2*(a2+c2+d2+e2-b2-f2);
949 double R = c2*d2*(a2+b2+e2+f2-c2-d2);
950 double S = a2*b2*d2+b2*c2*e2+a2*c2*f2+d2*e2*f2;
952 return sqrt(P+Q+R-S)/12.0;
955 inline double getVolume(const TSequenceOfXYZ& P){
956 gp_Vec aVec1( P( 2 ) - P( 1 ) );
957 gp_Vec aVec2( P( 3 ) - P( 1 ) );
958 gp_Vec aVec3( P( 4 ) - P( 1 ) );
959 gp_Vec anAreaVec( aVec1 ^ aVec2 );
960 return fabs(aVec3 * anAreaVec) / 6.0;
963 inline double getMaxHeight(double theLen[6])
965 double aHeight = std::max(theLen[0],theLen[1]);
966 aHeight = std::max(aHeight,theLen[2]);
967 aHeight = std::max(aHeight,theLen[3]);
968 aHeight = std::max(aHeight,theLen[4]);
969 aHeight = std::max(aHeight,theLen[5]);
975 double AspectRatio3D::GetValue( long theId )
978 myCurrElement = myMesh->FindElement( theId );
979 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_TETRA )
981 // Action from CoTech | ACTION 31.3:
982 // EURIWARE BO: Homogenize the formulas used to calculate the Controls in SMESH to fit with
983 // those of ParaView. The library used by ParaView for those calculations can be reused in SMESH.
984 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
985 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
986 aVal = Round( vtkMeshQuality::TetAspectRatio( avtkCell ));
991 if ( GetPoints( myCurrElement, P ))
992 aVal = Round( GetValue( P ));
997 double AspectRatio3D::GetValue( const TSequenceOfXYZ& P )
999 double aQuality = 0.0;
1000 if(myCurrElement->IsPoly()) return aQuality;
1002 int nbNodes = P.size();
1004 if(myCurrElement->IsQuadratic()) {
1005 if(nbNodes==10) nbNodes=4; // quadratic tetrahedron
1006 else if(nbNodes==13) nbNodes=5; // quadratic pyramid
1007 else if(nbNodes==15) nbNodes=6; // quadratic pentahedron
1008 else if(nbNodes==20) nbNodes=8; // quadratic hexahedron
1009 else if(nbNodes==27) nbNodes=8; // quadratic hexahedron
1010 else return aQuality;
1016 getDistance(P( 1 ),P( 2 )), // a
1017 getDistance(P( 2 ),P( 3 )), // b
1018 getDistance(P( 3 ),P( 1 )), // c
1019 getDistance(P( 2 ),P( 4 )), // d
1020 getDistance(P( 3 ),P( 4 )), // e
1021 getDistance(P( 1 ),P( 4 )) // f
1023 double aTria[4][3] = {
1024 {aLen[0],aLen[1],aLen[2]}, // abc
1025 {aLen[0],aLen[3],aLen[5]}, // adf
1026 {aLen[1],aLen[3],aLen[4]}, // bde
1027 {aLen[2],aLen[4],aLen[5]} // cef
1029 double aSumArea = 0.0;
1030 double aHalfPerimeter = getHalfPerimeter(aTria[0]);
1031 double anArea = getArea(aHalfPerimeter,aTria[0]);
1033 aHalfPerimeter = getHalfPerimeter(aTria[1]);
1034 anArea = getArea(aHalfPerimeter,aTria[1]);
1036 aHalfPerimeter = getHalfPerimeter(aTria[2]);
1037 anArea = getArea(aHalfPerimeter,aTria[2]);
1039 aHalfPerimeter = getHalfPerimeter(aTria[3]);
1040 anArea = getArea(aHalfPerimeter,aTria[3]);
1042 double aVolume = getVolume(P);
1043 //double aVolume = getVolume(aLen);
1044 double aHeight = getMaxHeight(aLen);
1045 static double aCoeff = sqrt(2.0)/12.0;
1046 if ( aVolume > DBL_MIN )
1047 aQuality = aCoeff*aHeight*aSumArea/aVolume;
1052 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 3 ),P( 5 )};
1053 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1056 gp_XYZ aXYZ[4] = {P( 1 ),P( 3 ),P( 4 ),P( 5 )};
1057 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1060 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 5 )};
1061 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1064 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 4 ),P( 5 )};
1065 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1071 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 6 )};
1072 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1075 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 3 )};
1076 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1079 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 6 )};
1080 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1083 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1084 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1087 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 6 )};
1088 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1091 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 3 )};
1092 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1098 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1099 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1102 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 4 )};
1103 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1106 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 7 )};
1107 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1110 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 8 )};
1111 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1114 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 3 )};
1115 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1118 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 4 )};
1119 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1122 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 7 )};
1123 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1126 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 8 )};
1127 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1130 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 3 )};
1131 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1134 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 4 )};
1135 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1138 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 7 )};
1139 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1142 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 8 )};
1143 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1146 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 1 )};
1147 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1150 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 2 )};
1151 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1154 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 5 )};
1155 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1158 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 6 )};
1159 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1162 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 1 )};
1163 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1166 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 2 )};
1167 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1170 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 5 )};
1171 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1174 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 6 )};
1175 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1178 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 1 )};
1179 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1182 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1183 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1186 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 5 )};
1187 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1190 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 6 )};
1191 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1194 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1195 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1198 gp_XYZ aXYZ[4] = {P( 4 ),P( 5 ),P( 8 ),P( 2 )};
1199 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1202 gp_XYZ aXYZ[4] = {P( 1 ),P( 4 ),P( 5 ),P( 3 )};
1203 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1206 gp_XYZ aXYZ[4] = {P( 3 ),P( 6 ),P( 7 ),P( 1 )};
1207 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1210 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 6 ),P( 4 )};
1211 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1214 gp_XYZ aXYZ[4] = {P( 5 ),P( 6 ),P( 8 ),P( 3 )};
1215 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1218 gp_XYZ aXYZ[4] = {P( 7 ),P( 8 ),P( 6 ),P( 1 )};
1219 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1222 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 7 )};
1223 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1226 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 2 ),P( 5 )};
1227 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1233 gp_XYZ aXYZ[8] = {P( 1 ),P( 2 ),P( 4 ),P( 5 ),P( 7 ),P( 8 ),P( 10 ),P( 11 )};
1234 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1237 gp_XYZ aXYZ[8] = {P( 2 ),P( 3 ),P( 5 ),P( 6 ),P( 8 ),P( 9 ),P( 11 ),P( 12 )};
1238 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1241 gp_XYZ aXYZ[8] = {P( 3 ),P( 4 ),P( 6 ),P( 1 ),P( 9 ),P( 10 ),P( 12 ),P( 7 )};
1242 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1245 } // switch(nbNodes)
1247 if ( nbNodes > 4 ) {
1248 // avaluate aspect ratio of quadranle faces
1249 AspectRatio aspect2D;
1250 SMDS_VolumeTool::VolumeType type = SMDS_VolumeTool::GetType( nbNodes );
1251 int nbFaces = SMDS_VolumeTool::NbFaces( type );
1252 TSequenceOfXYZ points(4);
1253 for ( int i = 0; i < nbFaces; ++i ) { // loop on faces of a volume
1254 if ( SMDS_VolumeTool::NbFaceNodes( type, i ) != 4 )
1256 const int* pInd = SMDS_VolumeTool::GetFaceNodesIndices( type, i, true );
1257 for ( int p = 0; p < 4; ++p ) // loop on nodes of a quadranle face
1258 points( p + 1 ) = P( pInd[ p ] + 1 );
1259 aQuality = std::max( aQuality, aspect2D.GetValue( points ));
1265 double AspectRatio3D::GetBadRate( double Value, int /*nbNodes*/ ) const
1267 // the aspect ratio is in the range [1.0,infinity]
1270 return Value / 1000.;
1273 SMDSAbs_ElementType AspectRatio3D::GetType() const
1275 return SMDSAbs_Volume;
1279 //================================================================================
1282 Description : Functor for calculating warping
1284 //================================================================================
1286 double Warping::GetValue( const TSequenceOfXYZ& P )
1288 if ( P.size() != 4 )
1291 gp_XYZ G = ( P( 1 ) + P( 2 ) + P( 3 ) + P( 4 ) ) / 4.;
1293 double A1 = ComputeA( P( 1 ), P( 2 ), P( 3 ), G );
1294 double A2 = ComputeA( P( 2 ), P( 3 ), P( 4 ), G );
1295 double A3 = ComputeA( P( 3 ), P( 4 ), P( 1 ), G );
1296 double A4 = ComputeA( P( 4 ), P( 1 ), P( 2 ), G );
1298 double val = Max( Max( A1, A2 ), Max( A3, A4 ) );
1300 const double eps = 0.1; // val is in degrees
1302 return val < eps ? 0. : val;
1305 double Warping::ComputeA( const gp_XYZ& thePnt1,
1306 const gp_XYZ& thePnt2,
1307 const gp_XYZ& thePnt3,
1308 const gp_XYZ& theG ) const
1310 double aLen1 = gp_Pnt( thePnt1 ).Distance( gp_Pnt( thePnt2 ) );
1311 double aLen2 = gp_Pnt( thePnt2 ).Distance( gp_Pnt( thePnt3 ) );
1312 double L = Min( aLen1, aLen2 ) * 0.5;
1316 gp_XYZ GI = ( thePnt2 + thePnt1 ) / 2. - theG;
1317 gp_XYZ GJ = ( thePnt3 + thePnt2 ) / 2. - theG;
1318 gp_XYZ N = GI.Crossed( GJ );
1320 if ( N.Modulus() < gp::Resolution() )
1325 double H = ( thePnt2 - theG ).Dot( N );
1326 return asin( fabs( H / L ) ) * 180. / M_PI;
1329 double Warping::GetBadRate( double Value, int /*nbNodes*/ ) const
1331 // the warp is in the range [0.0,PI/2]
1332 // 0.0 = good (no warp)
1333 // PI/2 = bad (face pliee)
1337 SMDSAbs_ElementType Warping::GetType() const
1339 return SMDSAbs_Face;
1343 //================================================================================
1346 Description : Functor for calculating taper
1348 //================================================================================
1350 double Taper::GetValue( const TSequenceOfXYZ& P )
1352 if ( P.size() != 4 )
1356 double J1 = getArea( P( 4 ), P( 1 ), P( 2 ) ) / 2.;
1357 double J2 = getArea( P( 3 ), P( 1 ), P( 2 ) ) / 2.;
1358 double J3 = getArea( P( 2 ), P( 3 ), P( 4 ) ) / 2.;
1359 double J4 = getArea( P( 3 ), P( 4 ), P( 1 ) ) / 2.;
1361 double JA = 0.25 * ( J1 + J2 + J3 + J4 );
1365 double T1 = fabs( ( J1 - JA ) / JA );
1366 double T2 = fabs( ( J2 - JA ) / JA );
1367 double T3 = fabs( ( J3 - JA ) / JA );
1368 double T4 = fabs( ( J4 - JA ) / JA );
1370 double val = Max( Max( T1, T2 ), Max( T3, T4 ) );
1372 const double eps = 0.01;
1374 return val < eps ? 0. : val;
1377 double Taper::GetBadRate( double Value, int /*nbNodes*/ ) const
1379 // the taper is in the range [0.0,1.0]
1380 // 0.0 = good (no taper)
1381 // 1.0 = bad (les cotes opposes sont allignes)
1385 SMDSAbs_ElementType Taper::GetType() const
1387 return SMDSAbs_Face;
1390 //================================================================================
1393 Description : Functor for calculating skew in degrees
1395 //================================================================================
1397 static inline double skewAngle( const gp_XYZ& p1, const gp_XYZ& p2, const gp_XYZ& p3 )
1399 gp_XYZ p12 = ( p2 + p1 ) / 2.;
1400 gp_XYZ p23 = ( p3 + p2 ) / 2.;
1401 gp_XYZ p31 = ( p3 + p1 ) / 2.;
1403 gp_Vec v1( p31 - p2 ), v2( p12 - p23 );
1405 return v1.Magnitude() < gp::Resolution() || v2.Magnitude() < gp::Resolution() ? 0. : v1.Angle( v2 );
1408 double Skew::GetValue( const TSequenceOfXYZ& P )
1410 if ( P.size() != 3 && P.size() != 4 )
1414 const double PI2 = M_PI / 2.;
1415 if ( P.size() == 3 )
1417 double A0 = fabs( PI2 - skewAngle( P( 3 ), P( 1 ), P( 2 ) ) );
1418 double A1 = fabs( PI2 - skewAngle( P( 1 ), P( 2 ), P( 3 ) ) );
1419 double A2 = fabs( PI2 - skewAngle( P( 2 ), P( 3 ), P( 1 ) ) );
1421 return Max( A0, Max( A1, A2 ) ) * 180. / M_PI;
1425 gp_XYZ p12 = ( P( 1 ) + P( 2 ) ) / 2.;
1426 gp_XYZ p23 = ( P( 2 ) + P( 3 ) ) / 2.;
1427 gp_XYZ p34 = ( P( 3 ) + P( 4 ) ) / 2.;
1428 gp_XYZ p41 = ( P( 4 ) + P( 1 ) ) / 2.;
1430 gp_Vec v1( p34 - p12 ), v2( p23 - p41 );
1431 double A = v1.Magnitude() <= gp::Resolution() || v2.Magnitude() <= gp::Resolution()
1432 ? 0. : fabs( PI2 - v1.Angle( v2 ) );
1434 double val = A * 180. / M_PI;
1436 const double eps = 0.1; // val is in degrees
1438 return val < eps ? 0. : val;
1442 double Skew::GetBadRate( double Value, int /*nbNodes*/ ) const
1444 // the skew is in the range [0.0,PI/2].
1450 SMDSAbs_ElementType Skew::GetType() const
1452 return SMDSAbs_Face;
1456 //================================================================================
1459 Description : Functor for calculating area
1461 //================================================================================
1463 double Area::GetValue( const TSequenceOfXYZ& P )
1466 if ( P.size() > 2 ) {
1467 gp_Vec aVec1( P(2) - P(1) );
1468 gp_Vec aVec2( P(3) - P(1) );
1469 gp_Vec SumVec = aVec1 ^ aVec2;
1470 for (int i=4; i<=P.size(); i++) {
1471 gp_Vec aVec1( P(i-1) - P(1) );
1472 gp_Vec aVec2( P(i) - P(1) );
1473 gp_Vec tmp = aVec1 ^ aVec2;
1476 val = SumVec.Magnitude() * 0.5;
1481 double Area::GetBadRate( double Value, int /*nbNodes*/ ) const
1483 // meaningless as it is not a quality control functor
1487 SMDSAbs_ElementType Area::GetType() const
1489 return SMDSAbs_Face;
1492 //================================================================================
1495 Description : Functor for calculating length of edge
1497 //================================================================================
1499 double Length::GetValue( const TSequenceOfXYZ& P )
1501 switch ( P.size() ) {
1502 case 2: return getDistance( P( 1 ), P( 2 ) );
1503 case 3: return getDistance( P( 1 ), P( 2 ) ) + getDistance( P( 2 ), P( 3 ) );
1508 double Length::GetBadRate( double Value, int /*nbNodes*/ ) const
1510 // meaningless as it is not quality control functor
1514 SMDSAbs_ElementType Length::GetType() const
1516 return SMDSAbs_Edge;
1519 //================================================================================
1522 Description : Functor for calculating length of edge
1524 //================================================================================
1526 double Length2D::GetValue( long theElementId)
1530 //cout<<"Length2D::GetValue"<<endl;
1531 if (GetPoints(theElementId,P)){
1532 //for(int jj=1; jj<=P.size(); jj++)
1533 // cout<<"jj="<<jj<<" P("<<P(jj).X()<<","<<P(jj).Y()<<","<<P(jj).Z()<<")"<<endl;
1535 double aVal;// = GetValue( P );
1536 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
1537 SMDSAbs_ElementType aType = aElem->GetType();
1546 aVal = getDistance( P( 1 ), P( 2 ) );
1549 else if (len == 3){ // quadratic edge
1550 aVal = getDistance(P( 1 ),P( 3 )) + getDistance(P( 3 ),P( 2 ));
1554 if (len == 3){ // triangles
1555 double L1 = getDistance(P( 1 ),P( 2 ));
1556 double L2 = getDistance(P( 2 ),P( 3 ));
1557 double L3 = getDistance(P( 3 ),P( 1 ));
1558 aVal = Max(L1,Max(L2,L3));
1561 else if (len == 4){ // quadrangles
1562 double L1 = getDistance(P( 1 ),P( 2 ));
1563 double L2 = getDistance(P( 2 ),P( 3 ));
1564 double L3 = getDistance(P( 3 ),P( 4 ));
1565 double L4 = getDistance(P( 4 ),P( 1 ));
1566 aVal = Max(Max(L1,L2),Max(L3,L4));
1569 if (len == 6){ // quadratic triangles
1570 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1571 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1572 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
1573 aVal = Max(L1,Max(L2,L3));
1574 //cout<<"L1="<<L1<<" L2="<<L2<<"L3="<<L3<<" aVal="<<aVal<<endl;
1577 else if (len == 8){ // quadratic quadrangles
1578 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1579 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1580 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
1581 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
1582 aVal = Max(Max(L1,L2),Max(L3,L4));
1585 case SMDSAbs_Volume:
1586 if (len == 4){ // tetraidrs
1587 double L1 = getDistance(P( 1 ),P( 2 ));
1588 double L2 = getDistance(P( 2 ),P( 3 ));
1589 double L3 = getDistance(P( 3 ),P( 1 ));
1590 double L4 = getDistance(P( 1 ),P( 4 ));
1591 double L5 = getDistance(P( 2 ),P( 4 ));
1592 double L6 = getDistance(P( 3 ),P( 4 ));
1593 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1596 else if (len == 5){ // piramids
1597 double L1 = getDistance(P( 1 ),P( 2 ));
1598 double L2 = getDistance(P( 2 ),P( 3 ));
1599 double L3 = getDistance(P( 3 ),P( 4 ));
1600 double L4 = getDistance(P( 4 ),P( 1 ));
1601 double L5 = getDistance(P( 1 ),P( 5 ));
1602 double L6 = getDistance(P( 2 ),P( 5 ));
1603 double L7 = getDistance(P( 3 ),P( 5 ));
1604 double L8 = getDistance(P( 4 ),P( 5 ));
1606 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1607 aVal = Max(aVal,Max(L7,L8));
1610 else if (len == 6){ // pentaidres
1611 double L1 = getDistance(P( 1 ),P( 2 ));
1612 double L2 = getDistance(P( 2 ),P( 3 ));
1613 double L3 = getDistance(P( 3 ),P( 1 ));
1614 double L4 = getDistance(P( 4 ),P( 5 ));
1615 double L5 = getDistance(P( 5 ),P( 6 ));
1616 double L6 = getDistance(P( 6 ),P( 4 ));
1617 double L7 = getDistance(P( 1 ),P( 4 ));
1618 double L8 = getDistance(P( 2 ),P( 5 ));
1619 double L9 = getDistance(P( 3 ),P( 6 ));
1621 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1622 aVal = Max(aVal,Max(Max(L7,L8),L9));
1625 else if (len == 8){ // hexaider
1626 double L1 = getDistance(P( 1 ),P( 2 ));
1627 double L2 = getDistance(P( 2 ),P( 3 ));
1628 double L3 = getDistance(P( 3 ),P( 4 ));
1629 double L4 = getDistance(P( 4 ),P( 1 ));
1630 double L5 = getDistance(P( 5 ),P( 6 ));
1631 double L6 = getDistance(P( 6 ),P( 7 ));
1632 double L7 = getDistance(P( 7 ),P( 8 ));
1633 double L8 = getDistance(P( 8 ),P( 5 ));
1634 double L9 = getDistance(P( 1 ),P( 5 ));
1635 double L10= getDistance(P( 2 ),P( 6 ));
1636 double L11= getDistance(P( 3 ),P( 7 ));
1637 double L12= getDistance(P( 4 ),P( 8 ));
1639 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1640 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1641 aVal = Max(aVal,Max(L11,L12));
1646 if (len == 10){ // quadratic tetraidrs
1647 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
1648 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
1649 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
1650 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1651 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
1652 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
1653 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1656 else if (len == 13){ // quadratic piramids
1657 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
1658 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
1659 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1660 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1661 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1662 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
1663 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
1664 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
1665 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1666 aVal = Max(aVal,Max(L7,L8));
1669 else if (len == 15){ // quadratic pentaidres
1670 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
1671 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
1672 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1673 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1674 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
1675 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
1676 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
1677 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
1678 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
1679 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1680 aVal = Max(aVal,Max(Max(L7,L8),L9));
1683 else if (len == 20){ // quadratic hexaider
1684 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
1685 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
1686 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
1687 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
1688 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
1689 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
1690 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
1691 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
1692 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
1693 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
1694 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
1695 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
1696 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1697 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1698 aVal = Max(aVal,Max(L11,L12));
1710 if ( myPrecision >= 0 )
1712 double prec = pow( 10., (double)( myPrecision ) );
1713 aVal = floor( aVal * prec + 0.5 ) / prec;
1722 double Length2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1724 // meaningless as it is not quality control functor
1728 SMDSAbs_ElementType Length2D::GetType() const
1730 return SMDSAbs_Face;
1733 Length2D::Value::Value(double theLength,long thePntId1, long thePntId2):
1736 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1737 if(thePntId1 > thePntId2){
1738 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1742 bool Length2D::Value::operator<(const Length2D::Value& x) const{
1743 if(myPntId[0] < x.myPntId[0]) return true;
1744 if(myPntId[0] == x.myPntId[0])
1745 if(myPntId[1] < x.myPntId[1]) return true;
1749 void Length2D::GetValues(TValues& theValues){
1751 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1752 for(; anIter->more(); ){
1753 const SMDS_MeshFace* anElem = anIter->next();
1755 if(anElem->IsQuadratic()) {
1756 const SMDS_VtkFace* F =
1757 dynamic_cast<const SMDS_VtkFace*>(anElem);
1758 // use special nodes iterator
1759 SMDS_ElemIteratorPtr anIter = F->interlacedNodesElemIterator();
1764 const SMDS_MeshElement* aNode;
1766 aNode = anIter->next();
1767 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1768 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1769 aNodeId[0] = aNodeId[1] = aNode->GetID();
1772 for(; anIter->more(); ){
1773 const SMDS_MeshNode* N1 = static_cast<const SMDS_MeshNode*> (anIter->next());
1774 P[2] = gp_Pnt(N1->X(),N1->Y(),N1->Z());
1775 aNodeId[2] = N1->GetID();
1776 aLength = P[1].Distance(P[2]);
1777 if(!anIter->more()) break;
1778 const SMDS_MeshNode* N2 = static_cast<const SMDS_MeshNode*> (anIter->next());
1779 P[3] = gp_Pnt(N2->X(),N2->Y(),N2->Z());
1780 aNodeId[3] = N2->GetID();
1781 aLength += P[2].Distance(P[3]);
1782 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1783 Value aValue2(aLength,aNodeId[2],aNodeId[3]);
1785 aNodeId[1] = aNodeId[3];
1786 theValues.insert(aValue1);
1787 theValues.insert(aValue2);
1789 aLength += P[2].Distance(P[0]);
1790 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1791 Value aValue2(aLength,aNodeId[2],aNodeId[0]);
1792 theValues.insert(aValue1);
1793 theValues.insert(aValue2);
1796 SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
1801 const SMDS_MeshElement* aNode;
1802 if(aNodesIter->more()){
1803 aNode = aNodesIter->next();
1804 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1805 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1806 aNodeId[0] = aNodeId[1] = aNode->GetID();
1809 for(; aNodesIter->more(); ){
1810 aNode = aNodesIter->next();
1811 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1812 long anId = aNode->GetID();
1814 P[2] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1816 aLength = P[1].Distance(P[2]);
1818 Value aValue(aLength,aNodeId[1],anId);
1821 theValues.insert(aValue);
1824 aLength = P[0].Distance(P[1]);
1826 Value aValue(aLength,aNodeId[0],aNodeId[1]);
1827 theValues.insert(aValue);
1832 //================================================================================
1834 Class : MultiConnection
1835 Description : Functor for calculating number of faces conneted to the edge
1837 //================================================================================
1839 double MultiConnection::GetValue( const TSequenceOfXYZ& P )
1843 double MultiConnection::GetValue( long theId )
1845 return getNbMultiConnection( myMesh, theId );
1848 double MultiConnection::GetBadRate( double Value, int /*nbNodes*/ ) const
1850 // meaningless as it is not quality control functor
1854 SMDSAbs_ElementType MultiConnection::GetType() const
1856 return SMDSAbs_Edge;
1859 //================================================================================
1861 Class : MultiConnection2D
1862 Description : Functor for calculating number of faces conneted to the edge
1864 //================================================================================
1866 double MultiConnection2D::GetValue( const TSequenceOfXYZ& P )
1871 double MultiConnection2D::GetValue( long theElementId )
1875 const SMDS_MeshElement* aFaceElem = myMesh->FindElement(theElementId);
1876 SMDSAbs_ElementType aType = aFaceElem->GetType();
1881 int i = 0, len = aFaceElem->NbNodes();
1882 SMDS_ElemIteratorPtr anIter = aFaceElem->nodesIterator();
1885 const SMDS_MeshNode *aNode, *aNode0;
1886 TColStd_MapOfInteger aMap, aMapPrev;
1888 for (i = 0; i <= len; i++) {
1893 if (anIter->more()) {
1894 aNode = (SMDS_MeshNode*)anIter->next();
1902 if (i == 0) aNode0 = aNode;
1904 SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
1905 while (anElemIter->more()) {
1906 const SMDS_MeshElement* anElem = anElemIter->next();
1907 if (anElem != 0 && anElem->GetType() == SMDSAbs_Face) {
1908 int anId = anElem->GetID();
1911 if (aMapPrev.Contains(anId)) {
1916 aResult = Max(aResult, aNb);
1927 double MultiConnection2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1929 // meaningless as it is not quality control functor
1933 SMDSAbs_ElementType MultiConnection2D::GetType() const
1935 return SMDSAbs_Face;
1938 MultiConnection2D::Value::Value(long thePntId1, long thePntId2)
1940 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1941 if(thePntId1 > thePntId2){
1942 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1946 bool MultiConnection2D::Value::operator<(const MultiConnection2D::Value& x) const{
1947 if(myPntId[0] < x.myPntId[0]) return true;
1948 if(myPntId[0] == x.myPntId[0])
1949 if(myPntId[1] < x.myPntId[1]) return true;
1953 void MultiConnection2D::GetValues(MValues& theValues){
1954 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1955 for(; anIter->more(); ){
1956 const SMDS_MeshFace* anElem = anIter->next();
1957 SMDS_ElemIteratorPtr aNodesIter;
1958 if ( anElem->IsQuadratic() )
1959 aNodesIter = dynamic_cast<const SMDS_VtkFace*>
1960 (anElem)->interlacedNodesElemIterator();
1962 aNodesIter = anElem->nodesIterator();
1965 //int aNbConnects=0;
1966 const SMDS_MeshNode* aNode0;
1967 const SMDS_MeshNode* aNode1;
1968 const SMDS_MeshNode* aNode2;
1969 if(aNodesIter->more()){
1970 aNode0 = (SMDS_MeshNode*) aNodesIter->next();
1972 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode1;
1973 aNodeId[0] = aNodeId[1] = aNodes->GetID();
1975 for(; aNodesIter->more(); ) {
1976 aNode2 = (SMDS_MeshNode*) aNodesIter->next();
1977 long anId = aNode2->GetID();
1980 Value aValue(aNodeId[1],aNodeId[2]);
1981 MValues::iterator aItr = theValues.find(aValue);
1982 if (aItr != theValues.end()){
1987 theValues[aValue] = 1;
1990 //cout << "NodeIds: "<<aNodeId[1]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1991 aNodeId[1] = aNodeId[2];
1994 Value aValue(aNodeId[0],aNodeId[2]);
1995 MValues::iterator aItr = theValues.find(aValue);
1996 if (aItr != theValues.end()) {
2001 theValues[aValue] = 1;
2004 //cout << "NodeIds: "<<aNodeId[0]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
2009 //================================================================================
2011 Class : BallDiameter
2012 Description : Functor returning diameter of a ball element
2014 //================================================================================
2016 double BallDiameter::GetValue( long theId )
2018 double diameter = 0;
2020 if ( const SMDS_BallElement* ball =
2021 dynamic_cast<const SMDS_BallElement*>( myMesh->FindElement( theId )))
2023 diameter = ball->GetDiameter();
2028 double BallDiameter::GetBadRate( double Value, int /*nbNodes*/ ) const
2030 // meaningless as it is not a quality control functor
2034 SMDSAbs_ElementType BallDiameter::GetType() const
2036 return SMDSAbs_Ball;
2044 //================================================================================
2046 Class : BadOrientedVolume
2047 Description : Predicate bad oriented volumes
2049 //================================================================================
2051 BadOrientedVolume::BadOrientedVolume()
2056 void BadOrientedVolume::SetMesh( const SMDS_Mesh* theMesh )
2061 bool BadOrientedVolume::IsSatisfy( long theId )
2066 SMDS_VolumeTool vTool( myMesh->FindElement( theId ));
2067 return !vTool.IsForward();
2070 SMDSAbs_ElementType BadOrientedVolume::GetType() const
2072 return SMDSAbs_Volume;
2076 Class : BareBorderVolume
2079 bool BareBorderVolume::IsSatisfy(long theElementId )
2081 SMDS_VolumeTool myTool;
2082 if ( myTool.Set( myMesh->FindElement(theElementId)))
2084 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2085 if ( myTool.IsFreeFace( iF ))
2087 const SMDS_MeshNode** n = myTool.GetFaceNodes(iF);
2088 vector< const SMDS_MeshNode*> nodes( n, n+myTool.NbFaceNodes(iF));
2089 if ( !myMesh->FindElement( nodes, SMDSAbs_Face, /*Nomedium=*/false))
2096 //================================================================================
2098 Class : BareBorderFace
2100 //================================================================================
2102 bool BareBorderFace::IsSatisfy(long theElementId )
2105 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2107 if ( face->GetType() == SMDSAbs_Face )
2109 int nbN = face->NbCornerNodes();
2110 for ( int i = 0; i < nbN && !ok; ++i )
2112 // check if a link is shared by another face
2113 const SMDS_MeshNode* n1 = face->GetNode( i );
2114 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2115 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2116 bool isShared = false;
2117 while ( !isShared && fIt->more() )
2119 const SMDS_MeshElement* f = fIt->next();
2120 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2124 const int iQuad = face->IsQuadratic();
2125 myLinkNodes.resize( 2 + iQuad);
2126 myLinkNodes[0] = n1;
2127 myLinkNodes[1] = n2;
2129 myLinkNodes[2] = face->GetNode( i+nbN );
2130 ok = !myMesh->FindElement( myLinkNodes, SMDSAbs_Edge, /*noMedium=*/false);
2138 //================================================================================
2140 Class : OverConstrainedVolume
2142 //================================================================================
2144 bool OverConstrainedVolume::IsSatisfy(long theElementId )
2146 // An element is over-constrained if it has N-1 free borders where
2147 // N is the number of edges/faces for a 2D/3D element.
2148 SMDS_VolumeTool myTool;
2149 if ( myTool.Set( myMesh->FindElement(theElementId)))
2151 int nbSharedFaces = 0;
2152 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2153 if ( !myTool.IsFreeFace( iF ) && ++nbSharedFaces > 1 )
2155 return ( nbSharedFaces == 1 );
2160 //================================================================================
2162 Class : OverConstrainedFace
2164 //================================================================================
2166 bool OverConstrainedFace::IsSatisfy(long theElementId )
2168 // An element is over-constrained if it has N-1 free borders where
2169 // N is the number of edges/faces for a 2D/3D element.
2170 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2171 if ( face->GetType() == SMDSAbs_Face )
2173 int nbSharedBorders = 0;
2174 int nbN = face->NbCornerNodes();
2175 for ( int i = 0; i < nbN; ++i )
2177 // check if a link is shared by another face
2178 const SMDS_MeshNode* n1 = face->GetNode( i );
2179 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2180 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2181 bool isShared = false;
2182 while ( !isShared && fIt->more() )
2184 const SMDS_MeshElement* f = fIt->next();
2185 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2187 if ( isShared && ++nbSharedBorders > 1 )
2190 return ( nbSharedBorders == 1 );
2195 //================================================================================
2197 Class : CoincidentNodes
2198 Description : Predicate of Coincident nodes
2200 //================================================================================
2202 CoincidentNodes::CoincidentNodes()
2207 bool CoincidentNodes::IsSatisfy( long theElementId )
2209 return myCoincidentIDs.Contains( theElementId );
2212 SMDSAbs_ElementType CoincidentNodes::GetType() const
2214 return SMDSAbs_Node;
2217 void CoincidentNodes::SetMesh( const SMDS_Mesh* theMesh )
2219 myMeshModifTracer.SetMesh( theMesh );
2220 if ( myMeshModifTracer.IsMeshModified() )
2222 TIDSortedNodeSet nodesToCheck;
2223 SMDS_NodeIteratorPtr nIt = theMesh->nodesIterator(/*idInceasingOrder=*/true);
2224 while ( nIt->more() )
2225 nodesToCheck.insert( nodesToCheck.end(), nIt->next() );
2227 list< list< const SMDS_MeshNode*> > nodeGroups;
2228 SMESH_OctreeNode::FindCoincidentNodes ( nodesToCheck, &nodeGroups, myToler );
2230 myCoincidentIDs.Clear();
2231 list< list< const SMDS_MeshNode*> >::iterator groupIt = nodeGroups.begin();
2232 for ( ; groupIt != nodeGroups.end(); ++groupIt )
2234 list< const SMDS_MeshNode*>& coincNodes = *groupIt;
2235 list< const SMDS_MeshNode*>::iterator n = coincNodes.begin();
2236 for ( ; n != coincNodes.end(); ++n )
2237 myCoincidentIDs.Add( (*n)->GetID() );
2242 //================================================================================
2244 Class : CoincidentElements
2245 Description : Predicate of Coincident Elements
2246 Note : This class is suitable only for visualization of Coincident Elements
2248 //================================================================================
2250 CoincidentElements::CoincidentElements()
2255 void CoincidentElements::SetMesh( const SMDS_Mesh* theMesh )
2260 bool CoincidentElements::IsSatisfy( long theElementId )
2262 if ( !myMesh ) return false;
2264 if ( const SMDS_MeshElement* e = myMesh->FindElement( theElementId ))
2266 if ( e->GetType() != GetType() ) return false;
2267 set< const SMDS_MeshNode* > elemNodes( e->begin_nodes(), e->end_nodes() );
2268 const int nbNodes = e->NbNodes();
2269 SMDS_ElemIteratorPtr invIt = (*elemNodes.begin())->GetInverseElementIterator( GetType() );
2270 while ( invIt->more() )
2272 const SMDS_MeshElement* e2 = invIt->next();
2273 if ( e2 == e || e2->NbNodes() != nbNodes ) continue;
2275 bool sameNodes = true;
2276 for ( size_t i = 0; i < elemNodes.size() && sameNodes; ++i )
2277 sameNodes = ( elemNodes.count( e2->GetNode( i )));
2285 SMDSAbs_ElementType CoincidentElements1D::GetType() const
2287 return SMDSAbs_Edge;
2289 SMDSAbs_ElementType CoincidentElements2D::GetType() const
2291 return SMDSAbs_Face;
2293 SMDSAbs_ElementType CoincidentElements3D::GetType() const
2295 return SMDSAbs_Volume;
2299 //================================================================================
2302 Description : Predicate for free borders
2304 //================================================================================
2306 FreeBorders::FreeBorders()
2311 void FreeBorders::SetMesh( const SMDS_Mesh* theMesh )
2316 bool FreeBorders::IsSatisfy( long theId )
2318 return getNbMultiConnection( myMesh, theId ) == 1;
2321 SMDSAbs_ElementType FreeBorders::GetType() const
2323 return SMDSAbs_Edge;
2327 //================================================================================
2330 Description : Predicate for free Edges
2332 //================================================================================
2334 FreeEdges::FreeEdges()
2339 void FreeEdges::SetMesh( const SMDS_Mesh* theMesh )
2344 bool FreeEdges::IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId )
2346 TColStd_MapOfInteger aMap;
2347 for ( int i = 0; i < 2; i++ )
2349 SMDS_ElemIteratorPtr anElemIter = theNodes[ i ]->GetInverseElementIterator(SMDSAbs_Face);
2350 while( anElemIter->more() )
2352 if ( const SMDS_MeshElement* anElem = anElemIter->next())
2354 const int anId = anElem->GetID();
2355 if ( anId != theFaceId && !aMap.Add( anId ))
2363 bool FreeEdges::IsSatisfy( long theId )
2368 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2369 if ( aFace == 0 || aFace->GetType() != SMDSAbs_Face || aFace->NbNodes() < 3 )
2372 SMDS_ElemIteratorPtr anIter;
2373 if ( aFace->IsQuadratic() ) {
2374 anIter = dynamic_cast<const SMDS_VtkFace*>
2375 (aFace)->interlacedNodesElemIterator();
2378 anIter = aFace->nodesIterator();
2383 int i = 0, nbNodes = aFace->NbNodes();
2384 std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
2385 while( anIter->more() )
2387 const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
2390 aNodes[ i++ ] = aNode;
2392 aNodes[ nbNodes ] = aNodes[ 0 ];
2394 for ( i = 0; i < nbNodes; i++ )
2395 if ( IsFreeEdge( &aNodes[ i ], theId ) )
2401 SMDSAbs_ElementType FreeEdges::GetType() const
2403 return SMDSAbs_Face;
2406 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2409 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
2410 if(thePntId1 > thePntId2){
2411 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
2415 bool FreeEdges::Border::operator<(const FreeEdges::Border& x) const{
2416 if(myPntId[0] < x.myPntId[0]) return true;
2417 if(myPntId[0] == x.myPntId[0])
2418 if(myPntId[1] < x.myPntId[1]) return true;
2422 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2423 FreeEdges::TBorders& theRegistry,
2424 FreeEdges::TBorders& theContainer)
2426 if(theRegistry.find(theBorder) == theRegistry.end()){
2427 theRegistry.insert(theBorder);
2428 theContainer.insert(theBorder);
2430 theContainer.erase(theBorder);
2434 void FreeEdges::GetBoreders(TBorders& theBorders)
2437 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2438 for(; anIter->more(); ){
2439 const SMDS_MeshFace* anElem = anIter->next();
2440 long anElemId = anElem->GetID();
2441 SMDS_ElemIteratorPtr aNodesIter;
2442 if ( anElem->IsQuadratic() )
2443 aNodesIter = static_cast<const SMDS_VtkFace*>(anElem)->
2444 interlacedNodesElemIterator();
2446 aNodesIter = anElem->nodesIterator();
2448 const SMDS_MeshElement* aNode;
2449 if(aNodesIter->more()){
2450 aNode = aNodesIter->next();
2451 aNodeId[0] = aNodeId[1] = aNode->GetID();
2453 for(; aNodesIter->more(); ){
2454 aNode = aNodesIter->next();
2455 long anId = aNode->GetID();
2456 Border aBorder(anElemId,aNodeId[1],anId);
2458 UpdateBorders(aBorder,aRegistry,theBorders);
2460 Border aBorder(anElemId,aNodeId[0],aNodeId[1]);
2461 UpdateBorders(aBorder,aRegistry,theBorders);
2465 //================================================================================
2468 Description : Predicate for free nodes
2470 //================================================================================
2472 FreeNodes::FreeNodes()
2477 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2482 bool FreeNodes::IsSatisfy( long theNodeId )
2484 const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2488 return (aNode->NbInverseElements() < 1);
2491 SMDSAbs_ElementType FreeNodes::GetType() const
2493 return SMDSAbs_Node;
2497 //================================================================================
2500 Description : Predicate for free faces
2502 //================================================================================
2504 FreeFaces::FreeFaces()
2509 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2514 bool FreeFaces::IsSatisfy( long theId )
2516 if (!myMesh) return false;
2517 // check that faces nodes refers to less than two common volumes
2518 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2519 if ( !aFace || aFace->GetType() != SMDSAbs_Face )
2522 int nbNode = aFace->NbNodes();
2524 // collect volumes check that number of volumss with count equal nbNode not less than 2
2525 typedef map< SMDS_MeshElement*, int > TMapOfVolume; // map of volume counters
2526 typedef map< SMDS_MeshElement*, int >::iterator TItrMapOfVolume; // iterator
2527 TMapOfVolume mapOfVol;
2529 SMDS_ElemIteratorPtr nodeItr = aFace->nodesIterator();
2530 while ( nodeItr->more() ) {
2531 const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>(nodeItr->next());
2532 if ( !aNode ) continue;
2533 SMDS_ElemIteratorPtr volItr = aNode->GetInverseElementIterator(SMDSAbs_Volume);
2534 while ( volItr->more() ) {
2535 SMDS_MeshElement* aVol = (SMDS_MeshElement*)volItr->next();
2536 TItrMapOfVolume itr = mapOfVol.insert(make_pair(aVol, 0)).first;
2541 TItrMapOfVolume volItr = mapOfVol.begin();
2542 TItrMapOfVolume volEnd = mapOfVol.end();
2543 for ( ; volItr != volEnd; ++volItr )
2544 if ( (*volItr).second >= nbNode )
2546 // face is not free if number of volumes constructed on thier nodes more than one
2550 SMDSAbs_ElementType FreeFaces::GetType() const
2552 return SMDSAbs_Face;
2555 //================================================================================
2557 Class : LinearOrQuadratic
2558 Description : Predicate to verify whether a mesh element is linear
2560 //================================================================================
2562 LinearOrQuadratic::LinearOrQuadratic()
2567 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2572 bool LinearOrQuadratic::IsSatisfy( long theId )
2574 if (!myMesh) return false;
2575 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2576 if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2578 return (!anElem->IsQuadratic());
2581 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2586 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2591 //================================================================================
2594 Description : Functor for check color of group to whic mesh element belongs to
2596 //================================================================================
2598 GroupColor::GroupColor()
2602 bool GroupColor::IsSatisfy( long theId )
2604 return (myIDs.find( theId ) != myIDs.end());
2607 void GroupColor::SetType( SMDSAbs_ElementType theType )
2612 SMDSAbs_ElementType GroupColor::GetType() const
2617 static bool isEqual( const Quantity_Color& theColor1,
2618 const Quantity_Color& theColor2 )
2620 // tolerance to compare colors
2621 const double tol = 5*1e-3;
2622 return ( fabs( theColor1.Red() - theColor2.Red() ) < tol &&
2623 fabs( theColor1.Green() - theColor2.Green() ) < tol &&
2624 fabs( theColor1.Blue() - theColor2.Blue() ) < tol );
2628 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2632 const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2636 int nbGrp = aMesh->GetNbGroups();
2640 // iterates on groups and find necessary elements ids
2641 const std::set<SMESHDS_GroupBase*>& aGroups = aMesh->GetGroups();
2642 set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
2643 for (; GrIt != aGroups.end(); GrIt++) {
2644 SMESHDS_GroupBase* aGrp = (*GrIt);
2647 // check type and color of group
2648 if ( !isEqual( myColor, aGrp->GetColor() ) )
2650 if ( myType != SMDSAbs_All && myType != (SMDSAbs_ElementType)aGrp->GetType() )
2653 SMDSAbs_ElementType aGrpElType = (SMDSAbs_ElementType)aGrp->GetType();
2654 if ( myType == aGrpElType || (myType == SMDSAbs_All && aGrpElType != SMDSAbs_Node) ) {
2655 // add elements IDS into control
2656 int aSize = aGrp->Extent();
2657 for (int i = 0; i < aSize; i++)
2658 myIDs.insert( aGrp->GetID(i+1) );
2663 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
2665 TCollection_AsciiString aStr = theStr;
2666 aStr.RemoveAll( ' ' );
2667 aStr.RemoveAll( '\t' );
2668 for ( int aPos = aStr.Search( ";;" ); aPos != -1; aPos = aStr.Search( ";;" ) )
2669 aStr.Remove( aPos, 2 );
2670 Standard_Real clr[3];
2671 clr[0] = clr[1] = clr[2] = 0.;
2672 for ( int i = 0; i < 3; i++ ) {
2673 TCollection_AsciiString tmpStr = aStr.Token( ";", i+1 );
2674 if ( !tmpStr.IsEmpty() && tmpStr.IsRealValue() )
2675 clr[i] = tmpStr.RealValue();
2677 myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
2680 //=======================================================================
2681 // name : GetRangeStr
2682 // Purpose : Get range as a string.
2683 // Example: "1,2,3,50-60,63,67,70-"
2684 //=======================================================================
2686 void GroupColor::GetColorStr( TCollection_AsciiString& theResStr ) const
2689 theResStr += TCollection_AsciiString( myColor.Red() );
2690 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2691 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2694 //================================================================================
2696 Class : ElemGeomType
2697 Description : Predicate to check element geometry type
2699 //================================================================================
2701 ElemGeomType::ElemGeomType()
2704 myType = SMDSAbs_All;
2705 myGeomType = SMDSGeom_TRIANGLE;
2708 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2713 bool ElemGeomType::IsSatisfy( long theId )
2715 if (!myMesh) return false;
2716 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2719 const SMDSAbs_ElementType anElemType = anElem->GetType();
2720 if ( myType != SMDSAbs_All && anElemType != myType )
2722 bool isOk = ( anElem->GetGeomType() == myGeomType );
2726 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2731 SMDSAbs_ElementType ElemGeomType::GetType() const
2736 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2738 myGeomType = theType;
2741 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2746 //================================================================================
2748 Class : ElemEntityType
2749 Description : Predicate to check element entity type
2751 //================================================================================
2753 ElemEntityType::ElemEntityType():
2755 myType( SMDSAbs_All ),
2756 myEntityType( SMDSEntity_0D )
2760 void ElemEntityType::SetMesh( const SMDS_Mesh* theMesh )
2765 bool ElemEntityType::IsSatisfy( long theId )
2767 if ( !myMesh ) return false;
2768 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2770 myEntityType == anElem->GetEntityType() &&
2771 ( myType == SMDSAbs_Edge || myType == SMDSAbs_Face || myType == SMDSAbs_Volume ));
2774 void ElemEntityType::SetType( SMDSAbs_ElementType theType )
2779 SMDSAbs_ElementType ElemEntityType::GetType() const
2784 void ElemEntityType::SetElemEntityType( SMDSAbs_EntityType theEntityType )
2786 myEntityType = theEntityType;
2789 SMDSAbs_EntityType ElemEntityType::GetElemEntityType() const
2791 return myEntityType;
2794 //================================================================================
2796 * \brief Class CoplanarFaces
2798 //================================================================================
2800 CoplanarFaces::CoplanarFaces()
2801 : myFaceID(0), myToler(0)
2804 void CoplanarFaces::SetMesh( const SMDS_Mesh* theMesh )
2806 myMeshModifTracer.SetMesh( theMesh );
2807 if ( myMeshModifTracer.IsMeshModified() )
2809 // Build a set of coplanar face ids
2811 myCoplanarIDs.clear();
2813 if ( !myMeshModifTracer.GetMesh() || !myFaceID || !myToler )
2816 const SMDS_MeshElement* face = myMeshModifTracer.GetMesh()->FindElement( myFaceID );
2817 if ( !face || face->GetType() != SMDSAbs_Face )
2821 gp_Vec myNorm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2825 const double radianTol = myToler * M_PI / 180.;
2826 std::set< SMESH_TLink > checkedLinks;
2828 std::list< pair< const SMDS_MeshElement*, gp_Vec > > faceQueue;
2829 faceQueue.push_back( make_pair( face, myNorm ));
2830 while ( !faceQueue.empty() )
2832 face = faceQueue.front().first;
2833 myNorm = faceQueue.front().second;
2834 faceQueue.pop_front();
2836 for ( int i = 0, nbN = face->NbCornerNodes(); i < nbN; ++i )
2838 const SMDS_MeshNode* n1 = face->GetNode( i );
2839 const SMDS_MeshNode* n2 = face->GetNode(( i+1 )%nbN);
2840 if ( !checkedLinks.insert( SMESH_TLink( n1, n2 )).second )
2842 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator(SMDSAbs_Face);
2843 while ( fIt->more() )
2845 const SMDS_MeshElement* f = fIt->next();
2846 if ( f->GetNodeIndex( n2 ) > -1 )
2848 gp_Vec norm = getNormale( static_cast<const SMDS_MeshFace*>(f), &normOK );
2849 if (!normOK || myNorm.Angle( norm ) <= radianTol)
2851 myCoplanarIDs.insert( f->GetID() );
2852 faceQueue.push_back( make_pair( f, norm ));
2860 bool CoplanarFaces::IsSatisfy( long theElementId )
2862 return myCoplanarIDs.count( theElementId );
2867 *Description : Predicate for Range of Ids.
2868 * Range may be specified with two ways.
2869 * 1. Using AddToRange method
2870 * 2. With SetRangeStr method. Parameter of this method is a string
2871 * like as "1,2,3,50-60,63,67,70-"
2874 //=======================================================================
2875 // name : RangeOfIds
2876 // Purpose : Constructor
2877 //=======================================================================
2878 RangeOfIds::RangeOfIds()
2881 myType = SMDSAbs_All;
2884 //=======================================================================
2886 // Purpose : Set mesh
2887 //=======================================================================
2888 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
2893 //=======================================================================
2894 // name : AddToRange
2895 // Purpose : Add ID to the range
2896 //=======================================================================
2897 bool RangeOfIds::AddToRange( long theEntityId )
2899 myIds.Add( theEntityId );
2903 //=======================================================================
2904 // name : GetRangeStr
2905 // Purpose : Get range as a string.
2906 // Example: "1,2,3,50-60,63,67,70-"
2907 //=======================================================================
2908 void RangeOfIds::GetRangeStr( TCollection_AsciiString& theResStr )
2912 TColStd_SequenceOfInteger anIntSeq;
2913 TColStd_SequenceOfAsciiString aStrSeq;
2915 TColStd_MapIteratorOfMapOfInteger anIter( myIds );
2916 for ( ; anIter.More(); anIter.Next() )
2918 int anId = anIter.Key();
2919 TCollection_AsciiString aStr( anId );
2920 anIntSeq.Append( anId );
2921 aStrSeq.Append( aStr );
2924 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2926 int aMinId = myMin( i );
2927 int aMaxId = myMax( i );
2929 TCollection_AsciiString aStr;
2930 if ( aMinId != IntegerFirst() )
2935 if ( aMaxId != IntegerLast() )
2938 // find position of the string in result sequence and insert string in it
2939 if ( anIntSeq.Length() == 0 )
2941 anIntSeq.Append( aMinId );
2942 aStrSeq.Append( aStr );
2946 if ( aMinId < anIntSeq.First() )
2948 anIntSeq.Prepend( aMinId );
2949 aStrSeq.Prepend( aStr );
2951 else if ( aMinId > anIntSeq.Last() )
2953 anIntSeq.Append( aMinId );
2954 aStrSeq.Append( aStr );
2957 for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
2958 if ( aMinId < anIntSeq( j ) )
2960 anIntSeq.InsertBefore( j, aMinId );
2961 aStrSeq.InsertBefore( j, aStr );
2967 if ( aStrSeq.Length() == 0 )
2970 theResStr = aStrSeq( 1 );
2971 for ( int j = 2, k = aStrSeq.Length(); j <= k; j++ )
2974 theResStr += aStrSeq( j );
2978 //=======================================================================
2979 // name : SetRangeStr
2980 // Purpose : Define range with string
2981 // Example of entry string: "1,2,3,50-60,63,67,70-"
2982 //=======================================================================
2983 bool RangeOfIds::SetRangeStr( const TCollection_AsciiString& theStr )
2989 TCollection_AsciiString aStr = theStr;
2990 aStr.RemoveAll( ' ' );
2991 aStr.RemoveAll( '\t' );
2993 for ( int aPos = aStr.Search( ",," ); aPos != -1; aPos = aStr.Search( ",," ) )
2994 aStr.Remove( aPos, 2 );
2996 TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
2998 while ( tmpStr != "" )
3000 tmpStr = aStr.Token( ",", i++ );
3001 int aPos = tmpStr.Search( '-' );
3005 if ( tmpStr.IsIntegerValue() )
3006 myIds.Add( tmpStr.IntegerValue() );
3012 TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
3013 TCollection_AsciiString aMinStr = tmpStr;
3015 while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
3016 while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
3018 if ( (!aMinStr.IsEmpty() && !aMinStr.IsIntegerValue()) ||
3019 (!aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue()) )
3022 myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
3023 myMax.Append( aMaxStr.IsEmpty() ? IntegerLast() : aMaxStr.IntegerValue() );
3030 //=======================================================================
3032 // Purpose : Get type of supported entities
3033 //=======================================================================
3034 SMDSAbs_ElementType RangeOfIds::GetType() const
3039 //=======================================================================
3041 // Purpose : Set type of supported entities
3042 //=======================================================================
3043 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
3048 //=======================================================================
3050 // Purpose : Verify whether entity satisfies to this rpedicate
3051 //=======================================================================
3052 bool RangeOfIds::IsSatisfy( long theId )
3057 if ( myType == SMDSAbs_Node )
3059 if ( myMesh->FindNode( theId ) == 0 )
3064 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
3065 if ( anElem == 0 || (myType != anElem->GetType() && myType != SMDSAbs_All ))
3069 if ( myIds.Contains( theId ) )
3072 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
3073 if ( theId >= myMin( i ) && theId <= myMax( i ) )
3081 Description : Base class for comparators
3083 Comparator::Comparator():
3087 Comparator::~Comparator()
3090 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
3093 myFunctor->SetMesh( theMesh );
3096 void Comparator::SetMargin( double theValue )
3098 myMargin = theValue;
3101 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
3103 myFunctor = theFunct;
3106 SMDSAbs_ElementType Comparator::GetType() const
3108 return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
3111 double Comparator::GetMargin()
3119 Description : Comparator "<"
3121 bool LessThan::IsSatisfy( long theId )
3123 return myFunctor && myFunctor->GetValue( theId ) < myMargin;
3129 Description : Comparator ">"
3131 bool MoreThan::IsSatisfy( long theId )
3133 return myFunctor && myFunctor->GetValue( theId ) > myMargin;
3139 Description : Comparator "="
3142 myToler(Precision::Confusion())
3145 bool EqualTo::IsSatisfy( long theId )
3147 return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
3150 void EqualTo::SetTolerance( double theToler )
3155 double EqualTo::GetTolerance()
3162 Description : Logical NOT predicate
3164 LogicalNOT::LogicalNOT()
3167 LogicalNOT::~LogicalNOT()
3170 bool LogicalNOT::IsSatisfy( long theId )
3172 return myPredicate && !myPredicate->IsSatisfy( theId );
3175 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
3178 myPredicate->SetMesh( theMesh );
3181 void LogicalNOT::SetPredicate( PredicatePtr thePred )
3183 myPredicate = thePred;
3186 SMDSAbs_ElementType LogicalNOT::GetType() const
3188 return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
3193 Class : LogicalBinary
3194 Description : Base class for binary logical predicate
3196 LogicalBinary::LogicalBinary()
3199 LogicalBinary::~LogicalBinary()
3202 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
3205 myPredicate1->SetMesh( theMesh );
3208 myPredicate2->SetMesh( theMesh );
3211 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
3213 myPredicate1 = thePredicate;
3216 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
3218 myPredicate2 = thePredicate;
3221 SMDSAbs_ElementType LogicalBinary::GetType() const
3223 if ( !myPredicate1 || !myPredicate2 )
3226 SMDSAbs_ElementType aType1 = myPredicate1->GetType();
3227 SMDSAbs_ElementType aType2 = myPredicate2->GetType();
3229 return aType1 == aType2 ? aType1 : SMDSAbs_All;
3235 Description : Logical AND
3237 bool LogicalAND::IsSatisfy( long theId )
3242 myPredicate1->IsSatisfy( theId ) &&
3243 myPredicate2->IsSatisfy( theId );
3249 Description : Logical OR
3251 bool LogicalOR::IsSatisfy( long theId )
3256 (myPredicate1->IsSatisfy( theId ) ||
3257 myPredicate2->IsSatisfy( theId ));
3271 void Filter::SetPredicate( PredicatePtr thePredicate )
3273 myPredicate = thePredicate;
3276 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3277 PredicatePtr thePredicate,
3278 TIdSequence& theSequence )
3280 theSequence.clear();
3282 if ( !theMesh || !thePredicate )
3285 thePredicate->SetMesh( theMesh );
3287 SMDS_ElemIteratorPtr elemIt = theMesh->elementsIterator( thePredicate->GetType() );
3289 while ( elemIt->more() ) {
3290 const SMDS_MeshElement* anElem = elemIt->next();
3291 long anId = anElem->GetID();
3292 if ( thePredicate->IsSatisfy( anId ) )
3293 theSequence.push_back( anId );
3298 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3299 Filter::TIdSequence& theSequence )
3301 GetElementsId(theMesh,myPredicate,theSequence);
3308 typedef std::set<SMDS_MeshFace*> TMapOfFacePtr;
3314 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
3315 SMDS_MeshNode* theNode2 )
3321 ManifoldPart::Link::~Link()
3327 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
3329 if ( myNode1 == theLink.myNode1 &&
3330 myNode2 == theLink.myNode2 )
3332 else if ( myNode1 == theLink.myNode2 &&
3333 myNode2 == theLink.myNode1 )
3339 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
3341 if(myNode1 < x.myNode1) return true;
3342 if(myNode1 == x.myNode1)
3343 if(myNode2 < x.myNode2) return true;
3347 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
3348 const ManifoldPart::Link& theLink2 )
3350 return theLink1.IsEqual( theLink2 );
3353 ManifoldPart::ManifoldPart()
3356 myAngToler = Precision::Angular();
3357 myIsOnlyManifold = true;
3360 ManifoldPart::~ManifoldPart()
3365 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
3371 SMDSAbs_ElementType ManifoldPart::GetType() const
3372 { return SMDSAbs_Face; }
3374 bool ManifoldPart::IsSatisfy( long theElementId )
3376 return myMapIds.Contains( theElementId );
3379 void ManifoldPart::SetAngleTolerance( const double theAngToler )
3380 { myAngToler = theAngToler; }
3382 double ManifoldPart::GetAngleTolerance() const
3383 { return myAngToler; }
3385 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
3386 { myIsOnlyManifold = theIsOnly; }
3388 void ManifoldPart::SetStartElem( const long theStartId )
3389 { myStartElemId = theStartId; }
3391 bool ManifoldPart::process()
3394 myMapBadGeomIds.Clear();
3396 myAllFacePtr.clear();
3397 myAllFacePtrIntDMap.clear();
3401 // collect all faces into own map
3402 SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3403 for (; anFaceItr->more(); )
3405 SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3406 myAllFacePtr.push_back( aFacePtr );
3407 myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3410 SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3414 // the map of non manifold links and bad geometry
3415 TMapOfLink aMapOfNonManifold;
3416 TColStd_MapOfInteger aMapOfTreated;
3418 // begin cycle on faces from start index and run on vector till the end
3419 // and from begin to start index to cover whole vector
3420 const int aStartIndx = myAllFacePtrIntDMap[aStartFace];
3421 bool isStartTreat = false;
3422 for ( int fi = aStartIndx; !isStartTreat || fi != aStartIndx ; fi++ )
3424 if ( fi == aStartIndx )
3425 isStartTreat = true;
3426 // as result next time when fi will be equal to aStartIndx
3428 SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3429 if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3432 aMapOfTreated.Add( aFacePtr->GetID() );
3433 TColStd_MapOfInteger aResFaces;
3434 if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3435 aMapOfNonManifold, aResFaces ) )
3437 TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3438 for ( ; anItr.More(); anItr.Next() )
3440 int aFaceId = anItr.Key();
3441 aMapOfTreated.Add( aFaceId );
3442 myMapIds.Add( aFaceId );
3445 if ( fi == ( myAllFacePtr.size() - 1 ) )
3447 } // end run on vector of faces
3448 return !myMapIds.IsEmpty();
3451 static void getLinks( const SMDS_MeshFace* theFace,
3452 ManifoldPart::TVectorOfLink& theLinks )
3454 int aNbNode = theFace->NbNodes();
3455 SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3457 SMDS_MeshNode* aNode = 0;
3458 for ( ; aNodeItr->more() && i <= aNbNode; )
3461 SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3465 SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3467 ManifoldPart::Link aLink( aN1, aN2 );
3468 theLinks.push_back( aLink );
3472 bool ManifoldPart::findConnected
3473 ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3474 SMDS_MeshFace* theStartFace,
3475 ManifoldPart::TMapOfLink& theNonManifold,
3476 TColStd_MapOfInteger& theResFaces )
3478 theResFaces.Clear();
3479 if ( !theAllFacePtrInt.size() )
3482 if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3484 myMapBadGeomIds.Add( theStartFace->GetID() );
3488 ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3489 ManifoldPart::TVectorOfLink aSeqOfBoundary;
3490 theResFaces.Add( theStartFace->GetID() );
3491 ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3493 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3494 aDMapLinkFace, theNonManifold, theStartFace );
3496 bool isDone = false;
3497 while ( !isDone && aMapOfBoundary.size() != 0 )
3499 bool isToReset = false;
3500 ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3501 for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3503 ManifoldPart::Link aLink = *pLink;
3504 if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3506 // each link could be treated only once
3507 aMapToSkip.insert( aLink );
3509 ManifoldPart::TVectorOfFacePtr aFaces;
3511 if ( myIsOnlyManifold &&
3512 (theNonManifold.find( aLink ) != theNonManifold.end()) )
3516 getFacesByLink( aLink, aFaces );
3517 // filter the element to keep only indicated elements
3518 ManifoldPart::TVectorOfFacePtr aFiltered;
3519 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3520 for ( ; pFace != aFaces.end(); ++pFace )
3522 SMDS_MeshFace* aFace = *pFace;
3523 if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3524 aFiltered.push_back( aFace );
3527 if ( aFaces.size() < 2 ) // no neihgbour faces
3529 else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3531 theNonManifold.insert( aLink );
3536 // compare normal with normals of neighbor element
3537 SMDS_MeshFace* aPrevFace = aDMapLinkFace[ aLink ];
3538 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3539 for ( ; pFace != aFaces.end(); ++pFace )
3541 SMDS_MeshFace* aNextFace = *pFace;
3542 if ( aPrevFace == aNextFace )
3544 int anNextFaceID = aNextFace->GetID();
3545 if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3546 // should not be with non manifold restriction. probably bad topology
3548 // check if face was treated and skipped
3549 if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3550 !isInPlane( aPrevFace, aNextFace ) )
3552 // add new element to connected and extend the boundaries.
3553 theResFaces.Add( anNextFaceID );
3554 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3555 aDMapLinkFace, theNonManifold, aNextFace );
3559 isDone = !isToReset;
3562 return !theResFaces.IsEmpty();
3565 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3566 const SMDS_MeshFace* theFace2 )
3568 gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3569 gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3570 if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3572 myMapBadGeomIds.Add( theFace2->GetID() );
3575 if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
3581 void ManifoldPart::expandBoundary
3582 ( ManifoldPart::TMapOfLink& theMapOfBoundary,
3583 ManifoldPart::TVectorOfLink& theSeqOfBoundary,
3584 ManifoldPart::TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
3585 ManifoldPart::TMapOfLink& theNonManifold,
3586 SMDS_MeshFace* theNextFace ) const
3588 ManifoldPart::TVectorOfLink aLinks;
3589 getLinks( theNextFace, aLinks );
3590 int aNbLink = (int)aLinks.size();
3591 for ( int i = 0; i < aNbLink; i++ )
3593 ManifoldPart::Link aLink = aLinks[ i ];
3594 if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3596 if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3598 if ( myIsOnlyManifold )
3600 // remove from boundary
3601 theMapOfBoundary.erase( aLink );
3602 ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3603 for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3605 ManifoldPart::Link aBoundLink = *pLink;
3606 if ( aBoundLink.IsEqual( aLink ) )
3608 theSeqOfBoundary.erase( pLink );
3616 theMapOfBoundary.insert( aLink );
3617 theSeqOfBoundary.push_back( aLink );
3618 theDMapLinkFacePtr[ aLink ] = theNextFace;
3623 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3624 ManifoldPart::TVectorOfFacePtr& theFaces ) const
3626 std::set<SMDS_MeshCell *> aSetOfFaces;
3627 // take all faces that shared first node
3628 SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3629 for ( ; anItr->more(); )
3631 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3634 aSetOfFaces.insert( aFace );
3636 // take all faces that shared second node
3637 anItr = theLink.myNode2->facesIterator();
3638 // find the common part of two sets
3639 for ( ; anItr->more(); )
3641 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3642 if ( aSetOfFaces.count( aFace ) )
3643 theFaces.push_back( aFace );
3652 ElementsOnSurface::ElementsOnSurface()
3655 myType = SMDSAbs_All;
3657 myToler = Precision::Confusion();
3658 myUseBoundaries = false;
3661 ElementsOnSurface::~ElementsOnSurface()
3665 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
3667 myMeshModifTracer.SetMesh( theMesh );
3668 if ( myMeshModifTracer.IsMeshModified())
3672 bool ElementsOnSurface::IsSatisfy( long theElementId )
3674 return myIds.Contains( theElementId );
3677 SMDSAbs_ElementType ElementsOnSurface::GetType() const
3680 void ElementsOnSurface::SetTolerance( const double theToler )
3682 if ( myToler != theToler )
3687 double ElementsOnSurface::GetTolerance() const
3690 void ElementsOnSurface::SetUseBoundaries( bool theUse )
3692 if ( myUseBoundaries != theUse ) {
3693 myUseBoundaries = theUse;
3694 SetSurface( mySurf, myType );
3698 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
3699 const SMDSAbs_ElementType theType )
3704 if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
3706 mySurf = TopoDS::Face( theShape );
3707 BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
3709 u1 = SA.FirstUParameter(),
3710 u2 = SA.LastUParameter(),
3711 v1 = SA.FirstVParameter(),
3712 v2 = SA.LastVParameter();
3713 Handle(Geom_Surface) surf = BRep_Tool::Surface( mySurf );
3714 myProjector.Init( surf, u1,u2, v1,v2 );
3718 void ElementsOnSurface::process()
3721 if ( mySurf.IsNull() )
3724 if ( !myMeshModifTracer.GetMesh() )
3727 myIds.ReSize( myMeshModifTracer.GetMesh()->GetMeshInfo().NbElements( myType ));
3729 SMDS_ElemIteratorPtr anIter = myMeshModifTracer.GetMesh()->elementsIterator( myType );
3730 for(; anIter->more(); )
3731 process( anIter->next() );
3734 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
3736 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3737 bool isSatisfy = true;
3738 for ( ; aNodeItr->more(); )
3740 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3741 if ( !isOnSurface( aNode ) )
3748 myIds.Add( theElemPtr->GetID() );
3751 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
3753 if ( mySurf.IsNull() )
3756 gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
3757 // double aToler2 = myToler * myToler;
3758 // if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
3760 // gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
3761 // if ( aPln.SquareDistance( aPnt ) > aToler2 )
3764 // else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
3766 // gp_Cylinder aCyl = Handle(Geom_CylindricalSurface)::DownCast(mySurf)->Cylinder();
3767 // double aRad = aCyl.Radius();
3768 // gp_Ax3 anAxis = aCyl.Position();
3769 // gp_XYZ aLoc = aCyl.Location().XYZ();
3770 // double aXDist = anAxis.XDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3771 // double aYDist = anAxis.YDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3772 // if ( fabs(aXDist*aXDist + aYDist*aYDist - aRad*aRad) > aToler2 )
3777 myProjector.Perform( aPnt );
3778 bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
3788 ElementsOnShape::ElementsOnShape()
3790 myType(SMDSAbs_All),
3791 myToler(Precision::Confusion()),
3792 myAllNodesFlag(false)
3796 ElementsOnShape::~ElementsOnShape()
3801 SMDSAbs_ElementType ElementsOnShape::GetType() const
3806 void ElementsOnShape::SetTolerance (const double theToler)
3808 if (myToler != theToler) {
3810 SetShape(myShape, myType);
3814 double ElementsOnShape::GetTolerance() const
3819 void ElementsOnShape::SetAllNodes (bool theAllNodes)
3821 myAllNodesFlag = theAllNodes;
3824 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
3829 void ElementsOnShape::SetShape (const TopoDS_Shape& theShape,
3830 const SMDSAbs_ElementType theType)
3834 if ( myShape.IsNull() ) return;
3836 TopTools_IndexedMapOfShape shapesMap;
3837 TopAbs_ShapeEnum shapeTypes[4] = { TopAbs_SOLID, TopAbs_FACE, TopAbs_EDGE, TopAbs_VERTEX };
3838 TopExp_Explorer sub;
3839 for ( int i = 0; i < 4; ++i )
3841 if ( shapesMap.IsEmpty() )
3842 for ( sub.Init( myShape, shapeTypes[i] ); sub.More(); sub.Next() )
3843 shapesMap.Add( sub.Current() );
3845 for ( sub.Init( myShape, shapeTypes[i], shapeTypes[i-1] ); sub.More(); sub.Next() )
3846 shapesMap.Add( sub.Current() );
3850 myClassifiers.resize( shapesMap.Extent() );
3851 for ( int i = 0; i < shapesMap.Extent(); ++i )
3852 myClassifiers[ i ] = new TClassifier( shapesMap( i+1 ), myToler );
3855 void ElementsOnShape::clearClassifiers()
3857 for ( size_t i = 0; i < myClassifiers.size(); ++i )
3858 delete myClassifiers[ i ];
3859 myClassifiers.clear();
3862 bool ElementsOnShape::IsSatisfy (long elemId)
3864 const SMDS_MeshElement* elem =
3865 ( myType == SMDSAbs_Node ? myMesh->FindNode( elemId ) : myMesh->FindElement( elemId ));
3866 if ( !elem || myClassifiers.empty() )
3869 for ( size_t i = 0; i < myClassifiers.size(); ++i )
3871 SMDS_ElemIteratorPtr aNodeItr = elem->nodesIterator();
3872 bool isSatisfy = myAllNodesFlag;
3874 gp_XYZ centerXYZ (0, 0, 0);
3876 while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
3878 SMESH_TNodeXYZ aPnt ( aNodeItr->next() );
3880 isSatisfy = ! myClassifiers[i]->IsOut( aPnt );
3883 // Check the center point for volumes MantisBug 0020168
3886 myClassifiers[i]->ShapeType() == TopAbs_SOLID)
3888 centerXYZ /= elem->NbNodes();
3889 isSatisfy = ! myClassifiers[i]->IsOut( centerXYZ );
3898 TopAbs_ShapeEnum ElementsOnShape::TClassifier::ShapeType() const
3900 return myShape.ShapeType();
3903 bool ElementsOnShape::TClassifier::IsOut(const gp_Pnt& p)
3905 return (this->*myIsOutFun)( p );
3908 void ElementsOnShape::TClassifier::Init (const TopoDS_Shape& theShape, double theTol)
3912 switch ( myShape.ShapeType() )
3914 case TopAbs_SOLID: {
3915 mySolidClfr.Load(theShape);
3916 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfSolid;
3920 Standard_Real u1,u2,v1,v2;
3921 Handle(Geom_Surface) surf = BRep_Tool::Surface( TopoDS::Face( theShape ));
3922 surf->Bounds( u1,u2,v1,v2 );
3923 myProjFace.Init(surf, u1,u2, v1,v2, myTol );
3924 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfFace;
3928 Standard_Real u1, u2;
3929 Handle(Geom_Curve) curve = BRep_Tool::Curve( TopoDS::Edge(theShape), u1, u2);
3930 myProjEdge.Init(curve, u1, u2);
3931 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfEdge;
3934 case TopAbs_VERTEX:{
3935 myVertexXYZ = BRep_Tool::Pnt( TopoDS::Vertex( theShape ) );
3936 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfVertex;
3940 throw SALOME_Exception("Programmer error in usage of ElementsOnShape::TClassifier");
3944 bool ElementsOnShape::TClassifier::isOutOfSolid (const gp_Pnt& p)
3946 mySolidClfr.Perform( p, myTol );
3947 return ( mySolidClfr.State() != TopAbs_IN && mySolidClfr.State() != TopAbs_ON );
3950 bool ElementsOnShape::TClassifier::isOutOfFace (const gp_Pnt& p)
3952 myProjFace.Perform( p );
3953 if ( myProjFace.IsDone() && myProjFace.LowerDistance() <= myTol )
3955 // check relatively to the face
3956 Quantity_Parameter u, v;
3957 myProjFace.LowerDistanceParameters(u, v);
3958 gp_Pnt2d aProjPnt (u, v);
3959 BRepClass_FaceClassifier aClsf ( TopoDS::Face( myShape ), aProjPnt, myTol );
3960 if ( aClsf.State() == TopAbs_IN || aClsf.State() == TopAbs_ON )
3966 bool ElementsOnShape::TClassifier::isOutOfEdge (const gp_Pnt& p)
3968 myProjEdge.Perform( p );
3969 return ! ( myProjEdge.NbPoints() > 0 && myProjEdge.LowerDistance() <= myTol );
3972 bool ElementsOnShape::TClassifier::isOutOfVertex(const gp_Pnt& p)
3974 return ( myVertexXYZ.Distance( p ) > myTol );
3978 TSequenceOfXYZ::TSequenceOfXYZ()
3981 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n)
3984 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t)
3987 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray)
3990 template <class InputIterator>
3991 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd)
3994 TSequenceOfXYZ::~TSequenceOfXYZ()
3997 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
3999 myArray = theSequenceOfXYZ.myArray;
4003 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
4005 return myArray[n-1];
4008 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
4010 return myArray[n-1];
4013 void TSequenceOfXYZ::clear()
4018 void TSequenceOfXYZ::reserve(size_type n)
4023 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
4025 myArray.push_back(v);
4028 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
4030 return myArray.size();
4033 TMeshModifTracer::TMeshModifTracer():
4034 myMeshModifTime(0), myMesh(0)
4037 void TMeshModifTracer::SetMesh( const SMDS_Mesh* theMesh )
4039 if ( theMesh != myMesh )
4040 myMeshModifTime = 0;
4043 bool TMeshModifTracer::IsMeshModified()
4045 bool modified = false;
4048 modified = ( myMeshModifTime != myMesh->GetMTime() );
4049 myMeshModifTime = myMesh->GetMTime();