1 // Copyright (C) 2007-2012 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
22 #include "SMESH_ControlsDef.hxx"
24 #include "SMDS_BallElement.hxx"
25 #include "SMDS_Iterator.hxx"
26 #include "SMDS_Mesh.hxx"
27 #include "SMDS_MeshElement.hxx"
28 #include "SMDS_MeshNode.hxx"
29 #include "SMDS_QuadraticEdge.hxx"
30 #include "SMDS_QuadraticFaceOfNodes.hxx"
31 #include "SMDS_VolumeTool.hxx"
32 #include "SMESHDS_GroupBase.hxx"
33 #include "SMESHDS_Mesh.hxx"
34 #include "SMESH_OctreeNode.hxx"
36 #include <BRepAdaptor_Surface.hxx>
37 #include <BRepClass_FaceClassifier.hxx>
38 #include <BRep_Tool.hxx>
39 #include <Geom_CylindricalSurface.hxx>
40 #include <Geom_Plane.hxx>
41 #include <Geom_Surface.hxx>
42 #include <Precision.hxx>
43 #include <TColStd_MapIteratorOfMapOfInteger.hxx>
44 #include <TColStd_MapOfInteger.hxx>
45 #include <TColStd_SequenceOfAsciiString.hxx>
46 #include <TColgp_Array1OfXYZ.hxx>
49 #include <TopoDS_Edge.hxx>
50 #include <TopoDS_Face.hxx>
51 #include <TopoDS_Iterator.hxx>
52 #include <TopoDS_Shape.hxx>
53 #include <TopoDS_Vertex.hxx>
55 #include <gp_Cylinder.hxx>
62 #include <vtkMeshQuality.h>
74 const double theEps = 1e-100;
75 const double theInf = 1e+100;
77 inline gp_XYZ gpXYZ(const SMDS_MeshNode* aNode )
79 return gp_XYZ(aNode->X(), aNode->Y(), aNode->Z() );
82 inline double getAngle( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
84 gp_Vec v1( P1 - P2 ), v2( P3 - P2 );
86 return v1.Magnitude() < gp::Resolution() ||
87 v2.Magnitude() < gp::Resolution() ? 0 : v1.Angle( v2 );
90 inline double getArea( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
92 gp_Vec aVec1( P2 - P1 );
93 gp_Vec aVec2( P3 - P1 );
94 return ( aVec1 ^ aVec2 ).Magnitude() * 0.5;
97 inline double getArea( const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3 )
99 return getArea( P1.XYZ(), P2.XYZ(), P3.XYZ() );
104 inline double getDistance( const gp_XYZ& P1, const gp_XYZ& P2 )
106 double aDist = gp_Pnt( P1 ).Distance( gp_Pnt( P2 ) );
110 int getNbMultiConnection( const SMDS_Mesh* theMesh, const int theId )
115 const SMDS_MeshElement* anEdge = theMesh->FindElement( theId );
116 if ( anEdge == 0 || anEdge->GetType() != SMDSAbs_Edge/* || anEdge->NbNodes() != 2 */)
119 // for each pair of nodes in anEdge (there are 2 pairs in a quadratic edge)
120 // count elements containing both nodes of the pair.
121 // Note that there may be such cases for a quadratic edge (a horizontal line):
126 // +-----+------+ +-----+------+
129 // result sould be 2 in both cases
131 int aResult0 = 0, aResult1 = 0;
132 // last node, it is a medium one in a quadratic edge
133 const SMDS_MeshNode* aLastNode = anEdge->GetNode( anEdge->NbNodes() - 1 );
134 const SMDS_MeshNode* aNode0 = anEdge->GetNode( 0 );
135 const SMDS_MeshNode* aNode1 = anEdge->GetNode( 1 );
136 if ( aNode1 == aLastNode ) aNode1 = 0;
138 SMDS_ElemIteratorPtr anElemIter = aLastNode->GetInverseElementIterator();
139 while( anElemIter->more() ) {
140 const SMDS_MeshElement* anElem = anElemIter->next();
141 if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
142 SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
143 while ( anIter->more() ) {
144 if ( const SMDS_MeshElement* anElemNode = anIter->next() ) {
145 if ( anElemNode == aNode0 ) {
147 if ( !aNode1 ) break; // not a quadratic edge
149 else if ( anElemNode == aNode1 )
155 int aResult = std::max ( aResult0, aResult1 );
157 // TColStd_MapOfInteger aMap;
159 // SMDS_ElemIteratorPtr anIter = anEdge->nodesIterator();
160 // if ( anIter != 0 ) {
161 // while( anIter->more() ) {
162 // const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
165 // SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
166 // while( anElemIter->more() ) {
167 // const SMDS_MeshElement* anElem = anElemIter->next();
168 // if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
169 // int anId = anElem->GetID();
171 // if ( anIter->more() ) // i.e. first node
173 // else if ( aMap.Contains( anId ) )
183 gp_XYZ getNormale( const SMDS_MeshFace* theFace, bool* ok=0 )
185 int aNbNode = theFace->NbNodes();
187 gp_XYZ q1 = gpXYZ( theFace->GetNode(1)) - gpXYZ( theFace->GetNode(0));
188 gp_XYZ q2 = gpXYZ( theFace->GetNode(2)) - gpXYZ( theFace->GetNode(0));
191 gp_XYZ q3 = gpXYZ( theFace->GetNode(3)) - gpXYZ( theFace->GetNode(0));
194 double len = n.Modulus();
195 bool zeroLen = ( len <= numeric_limits<double>::min());
199 if (ok) *ok = !zeroLen;
207 using namespace SMESH::Controls;
213 //================================================================================
215 Class : NumericalFunctor
216 Description : Base class for numerical functors
218 //================================================================================
220 NumericalFunctor::NumericalFunctor():
226 void NumericalFunctor::SetMesh( const SMDS_Mesh* theMesh )
231 bool NumericalFunctor::GetPoints(const int theId,
232 TSequenceOfXYZ& theRes ) const
239 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
240 if ( !anElem || anElem->GetType() != this->GetType() )
243 return GetPoints( anElem, theRes );
246 bool NumericalFunctor::GetPoints(const SMDS_MeshElement* anElem,
247 TSequenceOfXYZ& theRes )
254 theRes.reserve( anElem->NbNodes() );
256 // Get nodes of the element
257 SMDS_ElemIteratorPtr anIter;
259 if ( anElem->IsQuadratic() ) {
260 switch ( anElem->GetType() ) {
262 anIter = dynamic_cast<const SMDS_VtkEdge*>
263 (anElem)->interlacedNodesElemIterator();
266 anIter = dynamic_cast<const SMDS_VtkFace*>
267 (anElem)->interlacedNodesElemIterator();
270 anIter = anElem->nodesIterator();
275 anIter = anElem->nodesIterator();
279 while( anIter->more() ) {
280 if ( const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>( anIter->next() ))
281 theRes.push_back( gp_XYZ( aNode->X(), aNode->Y(), aNode->Z() ) );
288 long NumericalFunctor::GetPrecision() const
293 void NumericalFunctor::SetPrecision( const long thePrecision )
295 myPrecision = thePrecision;
296 myPrecisionValue = pow( 10., (double)( myPrecision ) );
299 double NumericalFunctor::GetValue( long theId )
303 myCurrElement = myMesh->FindElement( theId );
306 if ( GetPoints( theId, P ))
307 aVal = Round( GetValue( P ));
312 double NumericalFunctor::Round( const double & aVal )
314 return ( myPrecision >= 0 ) ? floor( aVal * myPrecisionValue + 0.5 ) / myPrecisionValue : aVal;
317 //================================================================================
319 * \brief Return histogram of functor values
320 * \param nbIntervals - number of intervals
321 * \param nbEvents - number of mesh elements having values within i-th interval
322 * \param funValues - boundaries of intervals
323 * \param elements - elements to check vulue of; empty list means "of all"
324 * \param minmax - boundaries of diapason of values to divide into intervals
326 //================================================================================
328 void NumericalFunctor::GetHistogram(int nbIntervals,
329 std::vector<int>& nbEvents,
330 std::vector<double>& funValues,
331 const vector<int>& elements,
332 const double* minmax,
333 const bool isLogarithmic)
335 if ( nbIntervals < 1 ||
337 !myMesh->GetMeshInfo().NbElements( GetType() ))
339 nbEvents.resize( nbIntervals, 0 );
340 funValues.resize( nbIntervals+1 );
342 // get all values sorted
343 std::multiset< double > values;
344 if ( elements.empty() )
346 SMDS_ElemIteratorPtr elemIt = myMesh->elementsIterator(GetType());
347 while ( elemIt->more() )
348 values.insert( GetValue( elemIt->next()->GetID() ));
352 vector<int>::const_iterator id = elements.begin();
353 for ( ; id != elements.end(); ++id )
354 values.insert( GetValue( *id ));
359 funValues[0] = minmax[0];
360 funValues[nbIntervals] = minmax[1];
364 funValues[0] = *values.begin();
365 funValues[nbIntervals] = *values.rbegin();
367 // case nbIntervals == 1
368 if ( nbIntervals == 1 )
370 nbEvents[0] = values.size();
374 if (funValues.front() == funValues.back())
376 nbEvents.resize( 1 );
377 nbEvents[0] = values.size();
378 funValues[1] = funValues.back();
379 funValues.resize( 2 );
382 std::multiset< double >::iterator min = values.begin(), max;
383 for ( int i = 0; i < nbIntervals; ++i )
385 // find end value of i-th interval
386 double r = (i+1) / double(nbIntervals);
387 if (isLogarithmic && funValues.front() > 1e-07 && funValues.back() > 1e-07) {
388 double logmin = log10(funValues.front());
389 double lval = logmin + r * (log10(funValues.back()) - logmin);
390 funValues[i+1] = pow(10.0, lval);
393 funValues[i+1] = funValues.front() * (1-r) + funValues.back() * r;
396 // count values in the i-th interval if there are any
397 if ( min != values.end() && *min <= funValues[i+1] )
399 // find the first value out of the interval
400 max = values.upper_bound( funValues[i+1] ); // max is greater than funValues[i+1], or end()
401 nbEvents[i] = std::distance( min, max );
405 // add values larger than minmax[1]
406 nbEvents.back() += std::distance( min, values.end() );
409 //=======================================================================
412 Description : Functor calculating volume of a 3D element
414 //================================================================================
416 double Volume::GetValue( long theElementId )
418 if ( theElementId && myMesh ) {
419 SMDS_VolumeTool aVolumeTool;
420 if ( aVolumeTool.Set( myMesh->FindElement( theElementId )))
421 return aVolumeTool.GetSize();
426 double Volume::GetBadRate( double Value, int /*nbNodes*/ ) const
431 SMDSAbs_ElementType Volume::GetType() const
433 return SMDSAbs_Volume;
436 //=======================================================================
438 Class : MaxElementLength2D
439 Description : Functor calculating maximum length of 2D element
441 //================================================================================
443 double MaxElementLength2D::GetValue( const TSequenceOfXYZ& P )
449 if( len == 3 ) { // triangles
450 double L1 = getDistance(P( 1 ),P( 2 ));
451 double L2 = getDistance(P( 2 ),P( 3 ));
452 double L3 = getDistance(P( 3 ),P( 1 ));
453 aVal = Max(L1,Max(L2,L3));
455 else if( len == 4 ) { // quadrangles
456 double L1 = getDistance(P( 1 ),P( 2 ));
457 double L2 = getDistance(P( 2 ),P( 3 ));
458 double L3 = getDistance(P( 3 ),P( 4 ));
459 double L4 = getDistance(P( 4 ),P( 1 ));
460 double D1 = getDistance(P( 1 ),P( 3 ));
461 double D2 = getDistance(P( 2 ),P( 4 ));
462 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
464 else if( len == 6 ) { // quadratic triangles
465 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
466 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
467 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
468 aVal = Max(L1,Max(L2,L3));
470 else if( len == 8 || len == 9 ) { // quadratic quadrangles
471 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
472 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
473 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
474 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
475 double D1 = getDistance(P( 1 ),P( 5 ));
476 double D2 = getDistance(P( 3 ),P( 7 ));
477 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
480 if( myPrecision >= 0 )
482 double prec = pow( 10., (double)myPrecision );
483 aVal = floor( aVal * prec + 0.5 ) / prec;
488 double MaxElementLength2D::GetValue( long theElementId )
491 return GetPoints( theElementId, P ) ? GetValue(P) : 0.0;
494 double MaxElementLength2D::GetBadRate( double Value, int /*nbNodes*/ ) const
499 SMDSAbs_ElementType MaxElementLength2D::GetType() const
504 //=======================================================================
506 Class : MaxElementLength3D
507 Description : Functor calculating maximum length of 3D element
509 //================================================================================
511 double MaxElementLength3D::GetValue( long theElementId )
514 if( GetPoints( theElementId, P ) ) {
516 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
517 SMDSAbs_ElementType aType = aElem->GetType();
521 if( len == 4 ) { // tetras
522 double L1 = getDistance(P( 1 ),P( 2 ));
523 double L2 = getDistance(P( 2 ),P( 3 ));
524 double L3 = getDistance(P( 3 ),P( 1 ));
525 double L4 = getDistance(P( 1 ),P( 4 ));
526 double L5 = getDistance(P( 2 ),P( 4 ));
527 double L6 = getDistance(P( 3 ),P( 4 ));
528 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
531 else if( len == 5 ) { // pyramids
532 double L1 = getDistance(P( 1 ),P( 2 ));
533 double L2 = getDistance(P( 2 ),P( 3 ));
534 double L3 = getDistance(P( 3 ),P( 4 ));
535 double L4 = getDistance(P( 4 ),P( 1 ));
536 double L5 = getDistance(P( 1 ),P( 5 ));
537 double L6 = getDistance(P( 2 ),P( 5 ));
538 double L7 = getDistance(P( 3 ),P( 5 ));
539 double L8 = getDistance(P( 4 ),P( 5 ));
540 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
541 aVal = Max(aVal,Max(L7,L8));
544 else if( len == 6 ) { // pentas
545 double L1 = getDistance(P( 1 ),P( 2 ));
546 double L2 = getDistance(P( 2 ),P( 3 ));
547 double L3 = getDistance(P( 3 ),P( 1 ));
548 double L4 = getDistance(P( 4 ),P( 5 ));
549 double L5 = getDistance(P( 5 ),P( 6 ));
550 double L6 = getDistance(P( 6 ),P( 4 ));
551 double L7 = getDistance(P( 1 ),P( 4 ));
552 double L8 = getDistance(P( 2 ),P( 5 ));
553 double L9 = getDistance(P( 3 ),P( 6 ));
554 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
555 aVal = Max(aVal,Max(Max(L7,L8),L9));
558 else if( len == 8 ) { // hexas
559 double L1 = getDistance(P( 1 ),P( 2 ));
560 double L2 = getDistance(P( 2 ),P( 3 ));
561 double L3 = getDistance(P( 3 ),P( 4 ));
562 double L4 = getDistance(P( 4 ),P( 1 ));
563 double L5 = getDistance(P( 5 ),P( 6 ));
564 double L6 = getDistance(P( 6 ),P( 7 ));
565 double L7 = getDistance(P( 7 ),P( 8 ));
566 double L8 = getDistance(P( 8 ),P( 5 ));
567 double L9 = getDistance(P( 1 ),P( 5 ));
568 double L10= getDistance(P( 2 ),P( 6 ));
569 double L11= getDistance(P( 3 ),P( 7 ));
570 double L12= getDistance(P( 4 ),P( 8 ));
571 double D1 = getDistance(P( 1 ),P( 7 ));
572 double D2 = getDistance(P( 2 ),P( 8 ));
573 double D3 = getDistance(P( 3 ),P( 5 ));
574 double D4 = getDistance(P( 4 ),P( 6 ));
575 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
576 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
577 aVal = Max(aVal,Max(L11,L12));
578 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
581 else if( len == 12 ) { // hexagonal prism
582 for ( int i1 = 1; i1 < 12; ++i1 )
583 for ( int i2 = i1+1; i1 <= 12; ++i1 )
584 aVal = Max( aVal, getDistance(P( i1 ),P( i2 )));
587 else if( len == 10 ) { // quadratic tetras
588 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
589 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
590 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
591 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
592 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
593 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
594 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
597 else if( len == 13 ) { // quadratic pyramids
598 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
599 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
600 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
601 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
602 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
603 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
604 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
605 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
606 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
607 aVal = Max(aVal,Max(L7,L8));
610 else if( len == 15 ) { // quadratic pentas
611 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
612 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
613 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
614 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
615 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
616 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
617 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
618 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
619 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
620 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
621 aVal = Max(aVal,Max(Max(L7,L8),L9));
624 else if( len == 20 || len == 27 ) { // quadratic hexas
625 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
626 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
627 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
628 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
629 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
630 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
631 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
632 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
633 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
634 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
635 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
636 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
637 double D1 = getDistance(P( 1 ),P( 7 ));
638 double D2 = getDistance(P( 2 ),P( 8 ));
639 double D3 = getDistance(P( 3 ),P( 5 ));
640 double D4 = getDistance(P( 4 ),P( 6 ));
641 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
642 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
643 aVal = Max(aVal,Max(L11,L12));
644 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
647 else if( len > 1 && aElem->IsPoly() ) { // polys
648 // get the maximum distance between all pairs of nodes
649 for( int i = 1; i <= len; i++ ) {
650 for( int j = 1; j <= len; j++ ) {
651 if( j > i ) { // optimization of the loop
652 double D = getDistance( P(i), P(j) );
653 aVal = Max( aVal, D );
660 if( myPrecision >= 0 )
662 double prec = pow( 10., (double)myPrecision );
663 aVal = floor( aVal * prec + 0.5 ) / prec;
670 double MaxElementLength3D::GetBadRate( double Value, int /*nbNodes*/ ) const
675 SMDSAbs_ElementType MaxElementLength3D::GetType() const
677 return SMDSAbs_Volume;
680 //=======================================================================
683 Description : Functor for calculation of minimum angle
685 //================================================================================
687 double MinimumAngle::GetValue( const TSequenceOfXYZ& P )
694 aMin = getAngle(P( P.size() ), P( 1 ), P( 2 ));
695 aMin = Min(aMin,getAngle(P( P.size()-1 ), P( P.size() ), P( 1 )));
697 for (int i=2; i<P.size();i++){
698 double A0 = getAngle( P( i-1 ), P( i ), P( i+1 ) );
702 return aMin * 180.0 / M_PI;
705 double MinimumAngle::GetBadRate( double Value, int nbNodes ) const
707 //const double aBestAngle = PI / nbNodes;
708 const double aBestAngle = 180.0 - ( 360.0 / double(nbNodes) );
709 return ( fabs( aBestAngle - Value ));
712 SMDSAbs_ElementType MinimumAngle::GetType() const
718 //================================================================================
721 Description : Functor for calculating aspect ratio
723 //================================================================================
725 double AspectRatio::GetValue( long theId )
728 myCurrElement = myMesh->FindElement( theId );
729 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_QUAD )
732 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
733 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
734 aVal = Round( vtkMeshQuality::QuadAspectRatio( avtkCell ));
739 if ( GetPoints( myCurrElement, P ))
740 aVal = Round( GetValue( P ));
745 double AspectRatio::GetValue( const TSequenceOfXYZ& P )
747 // According to "Mesh quality control" by Nadir Bouhamau referring to
748 // Pascal Jean Frey and Paul-Louis George. Maillages, applications aux elements finis.
749 // Hermes Science publications, Paris 1999 ISBN 2-7462-0024-4
752 int nbNodes = P.size();
757 // Compute aspect ratio
759 if ( nbNodes == 3 ) {
760 // Compute lengths of the sides
761 std::vector< double > aLen (nbNodes);
762 for ( int i = 0; i < nbNodes - 1; i++ )
763 aLen[ i ] = getDistance( P( i + 1 ), P( i + 2 ) );
764 aLen[ nbNodes - 1 ] = getDistance( P( 1 ), P( nbNodes ) );
765 // Q = alfa * h * p / S, where
767 // alfa = sqrt( 3 ) / 6
768 // h - length of the longest edge
769 // p - half perimeter
770 // S - triangle surface
771 const double alfa = sqrt( 3. ) / 6.;
772 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
773 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
774 double anArea = getArea( P( 1 ), P( 2 ), P( 3 ) );
775 if ( anArea <= theEps )
777 return alfa * maxLen * half_perimeter / anArea;
779 else if ( nbNodes == 6 ) { // quadratic triangles
780 // Compute lengths of the sides
781 std::vector< double > aLen (3);
782 aLen[0] = getDistance( P(1), P(3) );
783 aLen[1] = getDistance( P(3), P(5) );
784 aLen[2] = getDistance( P(5), P(1) );
785 // Q = alfa * h * p / S, where
787 // alfa = sqrt( 3 ) / 6
788 // h - length of the longest edge
789 // p - half perimeter
790 // S - triangle surface
791 const double alfa = sqrt( 3. ) / 6.;
792 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
793 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
794 double anArea = getArea( P(1), P(3), P(5) );
795 if ( anArea <= theEps )
797 return alfa * maxLen * half_perimeter / anArea;
799 else if( nbNodes == 4 ) { // quadrangle
800 // Compute lengths of the sides
801 std::vector< double > aLen (4);
802 aLen[0] = getDistance( P(1), P(2) );
803 aLen[1] = getDistance( P(2), P(3) );
804 aLen[2] = getDistance( P(3), P(4) );
805 aLen[3] = getDistance( P(4), P(1) );
806 // Compute lengths of the diagonals
807 std::vector< double > aDia (2);
808 aDia[0] = getDistance( P(1), P(3) );
809 aDia[1] = getDistance( P(2), P(4) );
810 // Compute areas of all triangles which can be built
811 // taking three nodes of the quadrangle
812 std::vector< double > anArea (4);
813 anArea[0] = getArea( P(1), P(2), P(3) );
814 anArea[1] = getArea( P(1), P(2), P(4) );
815 anArea[2] = getArea( P(1), P(3), P(4) );
816 anArea[3] = getArea( P(2), P(3), P(4) );
817 // Q = alpha * L * C1 / C2, where
819 // alpha = sqrt( 1/32 )
820 // L = max( L1, L2, L3, L4, D1, D2 )
821 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
822 // C2 = min( S1, S2, S3, S4 )
823 // Li - lengths of the edges
824 // Di - lengths of the diagonals
825 // Si - areas of the triangles
826 const double alpha = sqrt( 1 / 32. );
827 double L = Max( aLen[ 0 ],
831 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
832 double C1 = sqrt( ( aLen[0] * aLen[0] +
835 aLen[3] * aLen[3] ) / 4. );
836 double C2 = Min( anArea[ 0 ],
838 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
841 return alpha * L * C1 / C2;
843 else if( nbNodes == 8 || nbNodes == 9 ) { // nbNodes==8 - quadratic quadrangle
844 // Compute lengths of the sides
845 std::vector< double > aLen (4);
846 aLen[0] = getDistance( P(1), P(3) );
847 aLen[1] = getDistance( P(3), P(5) );
848 aLen[2] = getDistance( P(5), P(7) );
849 aLen[3] = getDistance( P(7), P(1) );
850 // Compute lengths of the diagonals
851 std::vector< double > aDia (2);
852 aDia[0] = getDistance( P(1), P(5) );
853 aDia[1] = getDistance( P(3), P(7) );
854 // Compute areas of all triangles which can be built
855 // taking three nodes of the quadrangle
856 std::vector< double > anArea (4);
857 anArea[0] = getArea( P(1), P(3), P(5) );
858 anArea[1] = getArea( P(1), P(3), P(7) );
859 anArea[2] = getArea( P(1), P(5), P(7) );
860 anArea[3] = getArea( P(3), P(5), P(7) );
861 // Q = alpha * L * C1 / C2, where
863 // alpha = sqrt( 1/32 )
864 // L = max( L1, L2, L3, L4, D1, D2 )
865 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
866 // C2 = min( S1, S2, S3, S4 )
867 // Li - lengths of the edges
868 // Di - lengths of the diagonals
869 // Si - areas of the triangles
870 const double alpha = sqrt( 1 / 32. );
871 double L = Max( aLen[ 0 ],
875 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
876 double C1 = sqrt( ( aLen[0] * aLen[0] +
879 aLen[3] * aLen[3] ) / 4. );
880 double C2 = Min( anArea[ 0 ],
882 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
885 return alpha * L * C1 / C2;
890 double AspectRatio::GetBadRate( double Value, int /*nbNodes*/ ) const
892 // the aspect ratio is in the range [1.0,infinity]
893 // < 1.0 = very bad, zero area
896 return ( Value < 0.9 ) ? 1000 : Value / 1000.;
899 SMDSAbs_ElementType AspectRatio::GetType() const
905 //================================================================================
907 Class : AspectRatio3D
908 Description : Functor for calculating aspect ratio
910 //================================================================================
914 inline double getHalfPerimeter(double theTria[3]){
915 return (theTria[0] + theTria[1] + theTria[2])/2.0;
918 inline double getArea(double theHalfPerim, double theTria[3]){
919 return sqrt(theHalfPerim*
920 (theHalfPerim-theTria[0])*
921 (theHalfPerim-theTria[1])*
922 (theHalfPerim-theTria[2]));
925 inline double getVolume(double theLen[6]){
926 double a2 = theLen[0]*theLen[0];
927 double b2 = theLen[1]*theLen[1];
928 double c2 = theLen[2]*theLen[2];
929 double d2 = theLen[3]*theLen[3];
930 double e2 = theLen[4]*theLen[4];
931 double f2 = theLen[5]*theLen[5];
932 double P = 4.0*a2*b2*d2;
933 double Q = a2*(b2+d2-e2)-b2*(a2+d2-f2)-d2*(a2+b2-c2);
934 double R = (b2+d2-e2)*(a2+d2-f2)*(a2+d2-f2);
935 return sqrt(P-Q+R)/12.0;
938 inline double getVolume2(double theLen[6]){
939 double a2 = theLen[0]*theLen[0];
940 double b2 = theLen[1]*theLen[1];
941 double c2 = theLen[2]*theLen[2];
942 double d2 = theLen[3]*theLen[3];
943 double e2 = theLen[4]*theLen[4];
944 double f2 = theLen[5]*theLen[5];
946 double P = a2*e2*(b2+c2+d2+f2-a2-e2);
947 double Q = b2*f2*(a2+c2+d2+e2-b2-f2);
948 double R = c2*d2*(a2+b2+e2+f2-c2-d2);
949 double S = a2*b2*d2+b2*c2*e2+a2*c2*f2+d2*e2*f2;
951 return sqrt(P+Q+R-S)/12.0;
954 inline double getVolume(const TSequenceOfXYZ& P){
955 gp_Vec aVec1( P( 2 ) - P( 1 ) );
956 gp_Vec aVec2( P( 3 ) - P( 1 ) );
957 gp_Vec aVec3( P( 4 ) - P( 1 ) );
958 gp_Vec anAreaVec( aVec1 ^ aVec2 );
959 return fabs(aVec3 * anAreaVec) / 6.0;
962 inline double getMaxHeight(double theLen[6])
964 double aHeight = std::max(theLen[0],theLen[1]);
965 aHeight = std::max(aHeight,theLen[2]);
966 aHeight = std::max(aHeight,theLen[3]);
967 aHeight = std::max(aHeight,theLen[4]);
968 aHeight = std::max(aHeight,theLen[5]);
974 double AspectRatio3D::GetValue( long theId )
977 myCurrElement = myMesh->FindElement( theId );
978 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_TETRA )
980 // Action from CoTech | ACTION 31.3:
981 // EURIWARE BO: Homogenize the formulas used to calculate the Controls in SMESH to fit with
982 // those of ParaView. The library used by ParaView for those calculations can be reused in SMESH.
983 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
984 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
985 aVal = Round( vtkMeshQuality::TetAspectRatio( avtkCell ));
990 if ( GetPoints( myCurrElement, P ))
991 aVal = Round( GetValue( P ));
996 double AspectRatio3D::GetValue( const TSequenceOfXYZ& P )
998 double aQuality = 0.0;
999 if(myCurrElement->IsPoly()) return aQuality;
1001 int nbNodes = P.size();
1003 if(myCurrElement->IsQuadratic()) {
1004 if(nbNodes==10) nbNodes=4; // quadratic tetrahedron
1005 else if(nbNodes==13) nbNodes=5; // quadratic pyramid
1006 else if(nbNodes==15) nbNodes=6; // quadratic pentahedron
1007 else if(nbNodes==20) nbNodes=8; // quadratic hexahedron
1008 else if(nbNodes==27) nbNodes=8; // quadratic hexahedron
1009 else return aQuality;
1015 getDistance(P( 1 ),P( 2 )), // a
1016 getDistance(P( 2 ),P( 3 )), // b
1017 getDistance(P( 3 ),P( 1 )), // c
1018 getDistance(P( 2 ),P( 4 )), // d
1019 getDistance(P( 3 ),P( 4 )), // e
1020 getDistance(P( 1 ),P( 4 )) // f
1022 double aTria[4][3] = {
1023 {aLen[0],aLen[1],aLen[2]}, // abc
1024 {aLen[0],aLen[3],aLen[5]}, // adf
1025 {aLen[1],aLen[3],aLen[4]}, // bde
1026 {aLen[2],aLen[4],aLen[5]} // cef
1028 double aSumArea = 0.0;
1029 double aHalfPerimeter = getHalfPerimeter(aTria[0]);
1030 double anArea = getArea(aHalfPerimeter,aTria[0]);
1032 aHalfPerimeter = getHalfPerimeter(aTria[1]);
1033 anArea = getArea(aHalfPerimeter,aTria[1]);
1035 aHalfPerimeter = getHalfPerimeter(aTria[2]);
1036 anArea = getArea(aHalfPerimeter,aTria[2]);
1038 aHalfPerimeter = getHalfPerimeter(aTria[3]);
1039 anArea = getArea(aHalfPerimeter,aTria[3]);
1041 double aVolume = getVolume(P);
1042 //double aVolume = getVolume(aLen);
1043 double aHeight = getMaxHeight(aLen);
1044 static double aCoeff = sqrt(2.0)/12.0;
1045 if ( aVolume > DBL_MIN )
1046 aQuality = aCoeff*aHeight*aSumArea/aVolume;
1051 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 3 ),P( 5 )};
1052 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1055 gp_XYZ aXYZ[4] = {P( 1 ),P( 3 ),P( 4 ),P( 5 )};
1056 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1059 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 5 )};
1060 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1063 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 4 ),P( 5 )};
1064 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1070 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 6 )};
1071 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1074 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 3 )};
1075 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1078 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 6 )};
1079 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1082 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1083 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1086 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 6 )};
1087 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1090 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 3 )};
1091 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1097 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1098 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1101 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 4 )};
1102 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1105 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 7 )};
1106 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1109 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 8 )};
1110 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1113 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 3 )};
1114 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1117 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 4 )};
1118 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1121 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 7 )};
1122 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1125 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 8 )};
1126 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1129 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 3 )};
1130 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1133 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 4 )};
1134 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1137 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 7 )};
1138 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1141 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 8 )};
1142 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1145 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 1 )};
1146 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1149 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 2 )};
1150 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1153 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 5 )};
1154 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1157 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 6 )};
1158 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1161 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 1 )};
1162 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1165 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 2 )};
1166 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1169 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 5 )};
1170 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1173 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 6 )};
1174 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1177 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 1 )};
1178 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1181 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1182 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1185 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 5 )};
1186 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1189 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 6 )};
1190 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1193 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1194 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1197 gp_XYZ aXYZ[4] = {P( 4 ),P( 5 ),P( 8 ),P( 2 )};
1198 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1201 gp_XYZ aXYZ[4] = {P( 1 ),P( 4 ),P( 5 ),P( 3 )};
1202 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1205 gp_XYZ aXYZ[4] = {P( 3 ),P( 6 ),P( 7 ),P( 1 )};
1206 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1209 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 6 ),P( 4 )};
1210 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1213 gp_XYZ aXYZ[4] = {P( 5 ),P( 6 ),P( 8 ),P( 3 )};
1214 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1217 gp_XYZ aXYZ[4] = {P( 7 ),P( 8 ),P( 6 ),P( 1 )};
1218 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1221 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 7 )};
1222 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1225 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 2 ),P( 5 )};
1226 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1232 gp_XYZ aXYZ[8] = {P( 1 ),P( 2 ),P( 4 ),P( 5 ),P( 7 ),P( 8 ),P( 10 ),P( 11 )};
1233 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1236 gp_XYZ aXYZ[8] = {P( 2 ),P( 3 ),P( 5 ),P( 6 ),P( 8 ),P( 9 ),P( 11 ),P( 12 )};
1237 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1240 gp_XYZ aXYZ[8] = {P( 3 ),P( 4 ),P( 6 ),P( 1 ),P( 9 ),P( 10 ),P( 12 ),P( 7 )};
1241 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1244 } // switch(nbNodes)
1246 if ( nbNodes > 4 ) {
1247 // avaluate aspect ratio of quadranle faces
1248 AspectRatio aspect2D;
1249 SMDS_VolumeTool::VolumeType type = SMDS_VolumeTool::GetType( nbNodes );
1250 int nbFaces = SMDS_VolumeTool::NbFaces( type );
1251 TSequenceOfXYZ points(4);
1252 for ( int i = 0; i < nbFaces; ++i ) { // loop on faces of a volume
1253 if ( SMDS_VolumeTool::NbFaceNodes( type, i ) != 4 )
1255 const int* pInd = SMDS_VolumeTool::GetFaceNodesIndices( type, i, true );
1256 for ( int p = 0; p < 4; ++p ) // loop on nodes of a quadranle face
1257 points( p + 1 ) = P( pInd[ p ] + 1 );
1258 aQuality = std::max( aQuality, aspect2D.GetValue( points ));
1264 double AspectRatio3D::GetBadRate( double Value, int /*nbNodes*/ ) const
1266 // the aspect ratio is in the range [1.0,infinity]
1269 return Value / 1000.;
1272 SMDSAbs_ElementType AspectRatio3D::GetType() const
1274 return SMDSAbs_Volume;
1278 //================================================================================
1281 Description : Functor for calculating warping
1283 //================================================================================
1285 double Warping::GetValue( const TSequenceOfXYZ& P )
1287 if ( P.size() != 4 )
1290 gp_XYZ G = ( P( 1 ) + P( 2 ) + P( 3 ) + P( 4 ) ) / 4.;
1292 double A1 = ComputeA( P( 1 ), P( 2 ), P( 3 ), G );
1293 double A2 = ComputeA( P( 2 ), P( 3 ), P( 4 ), G );
1294 double A3 = ComputeA( P( 3 ), P( 4 ), P( 1 ), G );
1295 double A4 = ComputeA( P( 4 ), P( 1 ), P( 2 ), G );
1297 double val = Max( Max( A1, A2 ), Max( A3, A4 ) );
1299 const double eps = 0.1; // val is in degrees
1301 return val < eps ? 0. : val;
1304 double Warping::ComputeA( const gp_XYZ& thePnt1,
1305 const gp_XYZ& thePnt2,
1306 const gp_XYZ& thePnt3,
1307 const gp_XYZ& theG ) const
1309 double aLen1 = gp_Pnt( thePnt1 ).Distance( gp_Pnt( thePnt2 ) );
1310 double aLen2 = gp_Pnt( thePnt2 ).Distance( gp_Pnt( thePnt3 ) );
1311 double L = Min( aLen1, aLen2 ) * 0.5;
1315 gp_XYZ GI = ( thePnt2 + thePnt1 ) / 2. - theG;
1316 gp_XYZ GJ = ( thePnt3 + thePnt2 ) / 2. - theG;
1317 gp_XYZ N = GI.Crossed( GJ );
1319 if ( N.Modulus() < gp::Resolution() )
1324 double H = ( thePnt2 - theG ).Dot( N );
1325 return asin( fabs( H / L ) ) * 180. / M_PI;
1328 double Warping::GetBadRate( double Value, int /*nbNodes*/ ) const
1330 // the warp is in the range [0.0,PI/2]
1331 // 0.0 = good (no warp)
1332 // PI/2 = bad (face pliee)
1336 SMDSAbs_ElementType Warping::GetType() const
1338 return SMDSAbs_Face;
1342 //================================================================================
1345 Description : Functor for calculating taper
1347 //================================================================================
1349 double Taper::GetValue( const TSequenceOfXYZ& P )
1351 if ( P.size() != 4 )
1355 double J1 = getArea( P( 4 ), P( 1 ), P( 2 ) ) / 2.;
1356 double J2 = getArea( P( 3 ), P( 1 ), P( 2 ) ) / 2.;
1357 double J3 = getArea( P( 2 ), P( 3 ), P( 4 ) ) / 2.;
1358 double J4 = getArea( P( 3 ), P( 4 ), P( 1 ) ) / 2.;
1360 double JA = 0.25 * ( J1 + J2 + J3 + J4 );
1364 double T1 = fabs( ( J1 - JA ) / JA );
1365 double T2 = fabs( ( J2 - JA ) / JA );
1366 double T3 = fabs( ( J3 - JA ) / JA );
1367 double T4 = fabs( ( J4 - JA ) / JA );
1369 double val = Max( Max( T1, T2 ), Max( T3, T4 ) );
1371 const double eps = 0.01;
1373 return val < eps ? 0. : val;
1376 double Taper::GetBadRate( double Value, int /*nbNodes*/ ) const
1378 // the taper is in the range [0.0,1.0]
1379 // 0.0 = good (no taper)
1380 // 1.0 = bad (les cotes opposes sont allignes)
1384 SMDSAbs_ElementType Taper::GetType() const
1386 return SMDSAbs_Face;
1389 //================================================================================
1392 Description : Functor for calculating skew in degrees
1394 //================================================================================
1396 static inline double skewAngle( const gp_XYZ& p1, const gp_XYZ& p2, const gp_XYZ& p3 )
1398 gp_XYZ p12 = ( p2 + p1 ) / 2.;
1399 gp_XYZ p23 = ( p3 + p2 ) / 2.;
1400 gp_XYZ p31 = ( p3 + p1 ) / 2.;
1402 gp_Vec v1( p31 - p2 ), v2( p12 - p23 );
1404 return v1.Magnitude() < gp::Resolution() || v2.Magnitude() < gp::Resolution() ? 0. : v1.Angle( v2 );
1407 double Skew::GetValue( const TSequenceOfXYZ& P )
1409 if ( P.size() != 3 && P.size() != 4 )
1413 const double PI2 = M_PI / 2.;
1414 if ( P.size() == 3 )
1416 double A0 = fabs( PI2 - skewAngle( P( 3 ), P( 1 ), P( 2 ) ) );
1417 double A1 = fabs( PI2 - skewAngle( P( 1 ), P( 2 ), P( 3 ) ) );
1418 double A2 = fabs( PI2 - skewAngle( P( 2 ), P( 3 ), P( 1 ) ) );
1420 return Max( A0, Max( A1, A2 ) ) * 180. / M_PI;
1424 gp_XYZ p12 = ( P( 1 ) + P( 2 ) ) / 2.;
1425 gp_XYZ p23 = ( P( 2 ) + P( 3 ) ) / 2.;
1426 gp_XYZ p34 = ( P( 3 ) + P( 4 ) ) / 2.;
1427 gp_XYZ p41 = ( P( 4 ) + P( 1 ) ) / 2.;
1429 gp_Vec v1( p34 - p12 ), v2( p23 - p41 );
1430 double A = v1.Magnitude() <= gp::Resolution() || v2.Magnitude() <= gp::Resolution()
1431 ? 0. : fabs( PI2 - v1.Angle( v2 ) );
1433 double val = A * 180. / M_PI;
1435 const double eps = 0.1; // val is in degrees
1437 return val < eps ? 0. : val;
1441 double Skew::GetBadRate( double Value, int /*nbNodes*/ ) const
1443 // the skew is in the range [0.0,PI/2].
1449 SMDSAbs_ElementType Skew::GetType() const
1451 return SMDSAbs_Face;
1455 //================================================================================
1458 Description : Functor for calculating area
1460 //================================================================================
1462 double Area::GetValue( const TSequenceOfXYZ& P )
1465 if ( P.size() > 2 ) {
1466 gp_Vec aVec1( P(2) - P(1) );
1467 gp_Vec aVec2( P(3) - P(1) );
1468 gp_Vec SumVec = aVec1 ^ aVec2;
1469 for (int i=4; i<=P.size(); i++) {
1470 gp_Vec aVec1( P(i-1) - P(1) );
1471 gp_Vec aVec2( P(i) - P(1) );
1472 gp_Vec tmp = aVec1 ^ aVec2;
1475 val = SumVec.Magnitude() * 0.5;
1480 double Area::GetBadRate( double Value, int /*nbNodes*/ ) const
1482 // meaningless as it is not a quality control functor
1486 SMDSAbs_ElementType Area::GetType() const
1488 return SMDSAbs_Face;
1491 //================================================================================
1494 Description : Functor for calculating length of edge
1496 //================================================================================
1498 double Length::GetValue( const TSequenceOfXYZ& P )
1500 switch ( P.size() ) {
1501 case 2: return getDistance( P( 1 ), P( 2 ) );
1502 case 3: return getDistance( P( 1 ), P( 2 ) ) + getDistance( P( 2 ), P( 3 ) );
1507 double Length::GetBadRate( double Value, int /*nbNodes*/ ) const
1509 // meaningless as it is not quality control functor
1513 SMDSAbs_ElementType Length::GetType() const
1515 return SMDSAbs_Edge;
1518 //================================================================================
1521 Description : Functor for calculating length of edge
1523 //================================================================================
1525 double Length2D::GetValue( long theElementId)
1529 //cout<<"Length2D::GetValue"<<endl;
1530 if (GetPoints(theElementId,P)){
1531 //for(int jj=1; jj<=P.size(); jj++)
1532 // cout<<"jj="<<jj<<" P("<<P(jj).X()<<","<<P(jj).Y()<<","<<P(jj).Z()<<")"<<endl;
1534 double aVal;// = GetValue( P );
1535 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
1536 SMDSAbs_ElementType aType = aElem->GetType();
1545 aVal = getDistance( P( 1 ), P( 2 ) );
1548 else if (len == 3){ // quadratic edge
1549 aVal = getDistance(P( 1 ),P( 3 )) + getDistance(P( 3 ),P( 2 ));
1553 if (len == 3){ // triangles
1554 double L1 = getDistance(P( 1 ),P( 2 ));
1555 double L2 = getDistance(P( 2 ),P( 3 ));
1556 double L3 = getDistance(P( 3 ),P( 1 ));
1557 aVal = Max(L1,Max(L2,L3));
1560 else if (len == 4){ // quadrangles
1561 double L1 = getDistance(P( 1 ),P( 2 ));
1562 double L2 = getDistance(P( 2 ),P( 3 ));
1563 double L3 = getDistance(P( 3 ),P( 4 ));
1564 double L4 = getDistance(P( 4 ),P( 1 ));
1565 aVal = Max(Max(L1,L2),Max(L3,L4));
1568 if (len == 6){ // quadratic triangles
1569 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1570 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1571 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
1572 aVal = Max(L1,Max(L2,L3));
1573 //cout<<"L1="<<L1<<" L2="<<L2<<"L3="<<L3<<" aVal="<<aVal<<endl;
1576 else if (len == 8){ // quadratic quadrangles
1577 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1578 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1579 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
1580 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
1581 aVal = Max(Max(L1,L2),Max(L3,L4));
1584 case SMDSAbs_Volume:
1585 if (len == 4){ // tetraidrs
1586 double L1 = getDistance(P( 1 ),P( 2 ));
1587 double L2 = getDistance(P( 2 ),P( 3 ));
1588 double L3 = getDistance(P( 3 ),P( 1 ));
1589 double L4 = getDistance(P( 1 ),P( 4 ));
1590 double L5 = getDistance(P( 2 ),P( 4 ));
1591 double L6 = getDistance(P( 3 ),P( 4 ));
1592 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1595 else if (len == 5){ // piramids
1596 double L1 = getDistance(P( 1 ),P( 2 ));
1597 double L2 = getDistance(P( 2 ),P( 3 ));
1598 double L3 = getDistance(P( 3 ),P( 4 ));
1599 double L4 = getDistance(P( 4 ),P( 1 ));
1600 double L5 = getDistance(P( 1 ),P( 5 ));
1601 double L6 = getDistance(P( 2 ),P( 5 ));
1602 double L7 = getDistance(P( 3 ),P( 5 ));
1603 double L8 = getDistance(P( 4 ),P( 5 ));
1605 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1606 aVal = Max(aVal,Max(L7,L8));
1609 else if (len == 6){ // pentaidres
1610 double L1 = getDistance(P( 1 ),P( 2 ));
1611 double L2 = getDistance(P( 2 ),P( 3 ));
1612 double L3 = getDistance(P( 3 ),P( 1 ));
1613 double L4 = getDistance(P( 4 ),P( 5 ));
1614 double L5 = getDistance(P( 5 ),P( 6 ));
1615 double L6 = getDistance(P( 6 ),P( 4 ));
1616 double L7 = getDistance(P( 1 ),P( 4 ));
1617 double L8 = getDistance(P( 2 ),P( 5 ));
1618 double L9 = getDistance(P( 3 ),P( 6 ));
1620 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1621 aVal = Max(aVal,Max(Max(L7,L8),L9));
1624 else if (len == 8){ // hexaider
1625 double L1 = getDistance(P( 1 ),P( 2 ));
1626 double L2 = getDistance(P( 2 ),P( 3 ));
1627 double L3 = getDistance(P( 3 ),P( 4 ));
1628 double L4 = getDistance(P( 4 ),P( 1 ));
1629 double L5 = getDistance(P( 5 ),P( 6 ));
1630 double L6 = getDistance(P( 6 ),P( 7 ));
1631 double L7 = getDistance(P( 7 ),P( 8 ));
1632 double L8 = getDistance(P( 8 ),P( 5 ));
1633 double L9 = getDistance(P( 1 ),P( 5 ));
1634 double L10= getDistance(P( 2 ),P( 6 ));
1635 double L11= getDistance(P( 3 ),P( 7 ));
1636 double L12= getDistance(P( 4 ),P( 8 ));
1638 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1639 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1640 aVal = Max(aVal,Max(L11,L12));
1645 if (len == 10){ // quadratic tetraidrs
1646 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
1647 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
1648 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
1649 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1650 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
1651 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
1652 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1655 else if (len == 13){ // quadratic piramids
1656 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
1657 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
1658 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1659 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1660 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1661 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
1662 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
1663 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
1664 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1665 aVal = Max(aVal,Max(L7,L8));
1668 else if (len == 15){ // quadratic pentaidres
1669 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
1670 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
1671 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1672 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1673 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
1674 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
1675 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
1676 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
1677 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
1678 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1679 aVal = Max(aVal,Max(Max(L7,L8),L9));
1682 else if (len == 20){ // quadratic hexaider
1683 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
1684 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
1685 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
1686 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
1687 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
1688 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
1689 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
1690 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
1691 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
1692 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
1693 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
1694 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
1695 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1696 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1697 aVal = Max(aVal,Max(L11,L12));
1709 if ( myPrecision >= 0 )
1711 double prec = pow( 10., (double)( myPrecision ) );
1712 aVal = floor( aVal * prec + 0.5 ) / prec;
1721 double Length2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1723 // meaningless as it is not quality control functor
1727 SMDSAbs_ElementType Length2D::GetType() const
1729 return SMDSAbs_Face;
1732 Length2D::Value::Value(double theLength,long thePntId1, long thePntId2):
1735 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1736 if(thePntId1 > thePntId2){
1737 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1741 bool Length2D::Value::operator<(const Length2D::Value& x) const{
1742 if(myPntId[0] < x.myPntId[0]) return true;
1743 if(myPntId[0] == x.myPntId[0])
1744 if(myPntId[1] < x.myPntId[1]) return true;
1748 void Length2D::GetValues(TValues& theValues){
1750 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1751 for(; anIter->more(); ){
1752 const SMDS_MeshFace* anElem = anIter->next();
1754 if(anElem->IsQuadratic()) {
1755 const SMDS_VtkFace* F =
1756 dynamic_cast<const SMDS_VtkFace*>(anElem);
1757 // use special nodes iterator
1758 SMDS_ElemIteratorPtr anIter = F->interlacedNodesElemIterator();
1763 const SMDS_MeshElement* aNode;
1765 aNode = anIter->next();
1766 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1767 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1768 aNodeId[0] = aNodeId[1] = aNode->GetID();
1771 for(; anIter->more(); ){
1772 const SMDS_MeshNode* N1 = static_cast<const SMDS_MeshNode*> (anIter->next());
1773 P[2] = gp_Pnt(N1->X(),N1->Y(),N1->Z());
1774 aNodeId[2] = N1->GetID();
1775 aLength = P[1].Distance(P[2]);
1776 if(!anIter->more()) break;
1777 const SMDS_MeshNode* N2 = static_cast<const SMDS_MeshNode*> (anIter->next());
1778 P[3] = gp_Pnt(N2->X(),N2->Y(),N2->Z());
1779 aNodeId[3] = N2->GetID();
1780 aLength += P[2].Distance(P[3]);
1781 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1782 Value aValue2(aLength,aNodeId[2],aNodeId[3]);
1784 aNodeId[1] = aNodeId[3];
1785 theValues.insert(aValue1);
1786 theValues.insert(aValue2);
1788 aLength += P[2].Distance(P[0]);
1789 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1790 Value aValue2(aLength,aNodeId[2],aNodeId[0]);
1791 theValues.insert(aValue1);
1792 theValues.insert(aValue2);
1795 SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
1800 const SMDS_MeshElement* aNode;
1801 if(aNodesIter->more()){
1802 aNode = aNodesIter->next();
1803 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1804 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1805 aNodeId[0] = aNodeId[1] = aNode->GetID();
1808 for(; aNodesIter->more(); ){
1809 aNode = aNodesIter->next();
1810 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1811 long anId = aNode->GetID();
1813 P[2] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1815 aLength = P[1].Distance(P[2]);
1817 Value aValue(aLength,aNodeId[1],anId);
1820 theValues.insert(aValue);
1823 aLength = P[0].Distance(P[1]);
1825 Value aValue(aLength,aNodeId[0],aNodeId[1]);
1826 theValues.insert(aValue);
1831 //================================================================================
1833 Class : MultiConnection
1834 Description : Functor for calculating number of faces conneted to the edge
1836 //================================================================================
1838 double MultiConnection::GetValue( const TSequenceOfXYZ& P )
1842 double MultiConnection::GetValue( long theId )
1844 return getNbMultiConnection( myMesh, theId );
1847 double MultiConnection::GetBadRate( double Value, int /*nbNodes*/ ) const
1849 // meaningless as it is not quality control functor
1853 SMDSAbs_ElementType MultiConnection::GetType() const
1855 return SMDSAbs_Edge;
1858 //================================================================================
1860 Class : MultiConnection2D
1861 Description : Functor for calculating number of faces conneted to the edge
1863 //================================================================================
1865 double MultiConnection2D::GetValue( const TSequenceOfXYZ& P )
1870 double MultiConnection2D::GetValue( long theElementId )
1874 const SMDS_MeshElement* aFaceElem = myMesh->FindElement(theElementId);
1875 SMDSAbs_ElementType aType = aFaceElem->GetType();
1880 int i = 0, len = aFaceElem->NbNodes();
1881 SMDS_ElemIteratorPtr anIter = aFaceElem->nodesIterator();
1884 const SMDS_MeshNode *aNode, *aNode0;
1885 TColStd_MapOfInteger aMap, aMapPrev;
1887 for (i = 0; i <= len; i++) {
1892 if (anIter->more()) {
1893 aNode = (SMDS_MeshNode*)anIter->next();
1901 if (i == 0) aNode0 = aNode;
1903 SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
1904 while (anElemIter->more()) {
1905 const SMDS_MeshElement* anElem = anElemIter->next();
1906 if (anElem != 0 && anElem->GetType() == SMDSAbs_Face) {
1907 int anId = anElem->GetID();
1910 if (aMapPrev.Contains(anId)) {
1915 aResult = Max(aResult, aNb);
1926 double MultiConnection2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1928 // meaningless as it is not quality control functor
1932 SMDSAbs_ElementType MultiConnection2D::GetType() const
1934 return SMDSAbs_Face;
1937 MultiConnection2D::Value::Value(long thePntId1, long thePntId2)
1939 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1940 if(thePntId1 > thePntId2){
1941 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1945 bool MultiConnection2D::Value::operator<(const MultiConnection2D::Value& x) const{
1946 if(myPntId[0] < x.myPntId[0]) return true;
1947 if(myPntId[0] == x.myPntId[0])
1948 if(myPntId[1] < x.myPntId[1]) return true;
1952 void MultiConnection2D::GetValues(MValues& theValues){
1953 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1954 for(; anIter->more(); ){
1955 const SMDS_MeshFace* anElem = anIter->next();
1956 SMDS_ElemIteratorPtr aNodesIter;
1957 if ( anElem->IsQuadratic() )
1958 aNodesIter = dynamic_cast<const SMDS_VtkFace*>
1959 (anElem)->interlacedNodesElemIterator();
1961 aNodesIter = anElem->nodesIterator();
1964 //int aNbConnects=0;
1965 const SMDS_MeshNode* aNode0;
1966 const SMDS_MeshNode* aNode1;
1967 const SMDS_MeshNode* aNode2;
1968 if(aNodesIter->more()){
1969 aNode0 = (SMDS_MeshNode*) aNodesIter->next();
1971 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode1;
1972 aNodeId[0] = aNodeId[1] = aNodes->GetID();
1974 for(; aNodesIter->more(); ) {
1975 aNode2 = (SMDS_MeshNode*) aNodesIter->next();
1976 long anId = aNode2->GetID();
1979 Value aValue(aNodeId[1],aNodeId[2]);
1980 MValues::iterator aItr = theValues.find(aValue);
1981 if (aItr != theValues.end()){
1986 theValues[aValue] = 1;
1989 //cout << "NodeIds: "<<aNodeId[1]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1990 aNodeId[1] = aNodeId[2];
1993 Value aValue(aNodeId[0],aNodeId[2]);
1994 MValues::iterator aItr = theValues.find(aValue);
1995 if (aItr != theValues.end()) {
2000 theValues[aValue] = 1;
2003 //cout << "NodeIds: "<<aNodeId[0]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
2008 //================================================================================
2010 Class : BallDiameter
2011 Description : Functor returning diameter of a ball element
2013 //================================================================================
2015 double BallDiameter::GetValue( long theId )
2017 double diameter = 0;
2019 if ( const SMDS_BallElement* ball =
2020 dynamic_cast<const SMDS_BallElement*>( myMesh->FindElement( theId )))
2022 diameter = ball->GetDiameter();
2027 double BallDiameter::GetBadRate( double Value, int /*nbNodes*/ ) const
2029 // meaningless as it is not a quality control functor
2033 SMDSAbs_ElementType BallDiameter::GetType() const
2035 return SMDSAbs_Ball;
2043 //================================================================================
2045 Class : BadOrientedVolume
2046 Description : Predicate bad oriented volumes
2048 //================================================================================
2050 BadOrientedVolume::BadOrientedVolume()
2055 void BadOrientedVolume::SetMesh( const SMDS_Mesh* theMesh )
2060 bool BadOrientedVolume::IsSatisfy( long theId )
2065 SMDS_VolumeTool vTool( myMesh->FindElement( theId ));
2066 return !vTool.IsForward();
2069 SMDSAbs_ElementType BadOrientedVolume::GetType() const
2071 return SMDSAbs_Volume;
2075 Class : BareBorderVolume
2078 bool BareBorderVolume::IsSatisfy(long theElementId )
2080 SMDS_VolumeTool myTool;
2081 if ( myTool.Set( myMesh->FindElement(theElementId)))
2083 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2084 if ( myTool.IsFreeFace( iF ))
2086 const SMDS_MeshNode** n = myTool.GetFaceNodes(iF);
2087 vector< const SMDS_MeshNode*> nodes( n, n+myTool.NbFaceNodes(iF));
2088 if ( !myMesh->FindElement( nodes, SMDSAbs_Face, /*Nomedium=*/false))
2095 //================================================================================
2097 Class : BareBorderFace
2099 //================================================================================
2101 bool BareBorderFace::IsSatisfy(long theElementId )
2104 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2106 if ( face->GetType() == SMDSAbs_Face )
2108 int nbN = face->NbCornerNodes();
2109 for ( int i = 0; i < nbN && !ok; ++i )
2111 // check if a link is shared by another face
2112 const SMDS_MeshNode* n1 = face->GetNode( i );
2113 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2114 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2115 bool isShared = false;
2116 while ( !isShared && fIt->more() )
2118 const SMDS_MeshElement* f = fIt->next();
2119 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2123 const int iQuad = face->IsQuadratic();
2124 myLinkNodes.resize( 2 + iQuad);
2125 myLinkNodes[0] = n1;
2126 myLinkNodes[1] = n2;
2128 myLinkNodes[2] = face->GetNode( i+nbN );
2129 ok = !myMesh->FindElement( myLinkNodes, SMDSAbs_Edge, /*noMedium=*/false);
2137 //================================================================================
2139 Class : OverConstrainedVolume
2141 //================================================================================
2143 bool OverConstrainedVolume::IsSatisfy(long theElementId )
2145 // An element is over-constrained if it has N-1 free borders where
2146 // N is the number of edges/faces for a 2D/3D element.
2147 SMDS_VolumeTool myTool;
2148 if ( myTool.Set( myMesh->FindElement(theElementId)))
2150 int nbSharedFaces = 0;
2151 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2152 if ( !myTool.IsFreeFace( iF ) && ++nbSharedFaces > 1 )
2154 return ( nbSharedFaces == 1 );
2159 //================================================================================
2161 Class : OverConstrainedFace
2163 //================================================================================
2165 bool OverConstrainedFace::IsSatisfy(long theElementId )
2167 // An element is over-constrained if it has N-1 free borders where
2168 // N is the number of edges/faces for a 2D/3D element.
2169 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2170 if ( face->GetType() == SMDSAbs_Face )
2172 int nbSharedBorders = 0;
2173 int nbN = face->NbCornerNodes();
2174 for ( int i = 0; i < nbN; ++i )
2176 // check if a link is shared by another face
2177 const SMDS_MeshNode* n1 = face->GetNode( i );
2178 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2179 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2180 bool isShared = false;
2181 while ( !isShared && fIt->more() )
2183 const SMDS_MeshElement* f = fIt->next();
2184 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2186 if ( isShared && ++nbSharedBorders > 1 )
2189 return ( nbSharedBorders == 1 );
2194 //================================================================================
2196 Class : CoincidentNodes
2197 Description : Predicate of Coincident nodes
2199 //================================================================================
2201 CoincidentNodes::CoincidentNodes()
2206 bool CoincidentNodes::IsSatisfy( long theElementId )
2208 return myCoincidentIDs.Contains( theElementId );
2211 SMDSAbs_ElementType CoincidentNodes::GetType() const
2213 return SMDSAbs_Node;
2216 void CoincidentNodes::SetMesh( const SMDS_Mesh* theMesh )
2218 myMeshModifTracer.SetMesh( theMesh );
2219 if ( myMeshModifTracer.IsMeshModified() )
2221 TIDSortedNodeSet nodesToCheck;
2222 SMDS_NodeIteratorPtr nIt = theMesh->nodesIterator(/*idInceasingOrder=*/true);
2223 while ( nIt->more() )
2224 nodesToCheck.insert( nodesToCheck.end(), nIt->next() );
2226 list< list< const SMDS_MeshNode*> > nodeGroups;
2227 SMESH_OctreeNode::FindCoincidentNodes ( nodesToCheck, &nodeGroups, myToler );
2229 myCoincidentIDs.Clear();
2230 list< list< const SMDS_MeshNode*> >::iterator groupIt = nodeGroups.begin();
2231 for ( ; groupIt != nodeGroups.end(); ++groupIt )
2233 list< const SMDS_MeshNode*>& coincNodes = *groupIt;
2234 list< const SMDS_MeshNode*>::iterator n = coincNodes.begin();
2235 for ( ; n != coincNodes.end(); ++n )
2236 myCoincidentIDs.Add( (*n)->GetID() );
2241 //================================================================================
2243 Class : CoincidentElements
2244 Description : Predicate of Coincident Elements
2245 Note : This class is suitable only for visualization of Coincident Elements
2247 //================================================================================
2249 CoincidentElements::CoincidentElements()
2254 void CoincidentElements::SetMesh( const SMDS_Mesh* theMesh )
2259 bool CoincidentElements::IsSatisfy( long theElementId )
2261 if ( !myMesh ) return false;
2263 if ( const SMDS_MeshElement* e = myMesh->FindElement( theElementId ))
2265 if ( e->GetType() != GetType() ) return false;
2266 set< const SMDS_MeshNode* > elemNodes( e->begin_nodes(), e->end_nodes() );
2267 const int nbNodes = e->NbNodes();
2268 SMDS_ElemIteratorPtr invIt = (*elemNodes.begin())->GetInverseElementIterator( GetType() );
2269 while ( invIt->more() )
2271 const SMDS_MeshElement* e2 = invIt->next();
2272 if ( e2 == e || e2->NbNodes() != nbNodes ) continue;
2274 bool sameNodes = true;
2275 for ( size_t i = 0; i < elemNodes.size() && sameNodes; ++i )
2276 sameNodes = ( elemNodes.count( e2->GetNode( i )));
2284 SMDSAbs_ElementType CoincidentElements1D::GetType() const
2286 return SMDSAbs_Edge;
2288 SMDSAbs_ElementType CoincidentElements2D::GetType() const
2290 return SMDSAbs_Face;
2292 SMDSAbs_ElementType CoincidentElements3D::GetType() const
2294 return SMDSAbs_Volume;
2298 //================================================================================
2301 Description : Predicate for free borders
2303 //================================================================================
2305 FreeBorders::FreeBorders()
2310 void FreeBorders::SetMesh( const SMDS_Mesh* theMesh )
2315 bool FreeBorders::IsSatisfy( long theId )
2317 return getNbMultiConnection( myMesh, theId ) == 1;
2320 SMDSAbs_ElementType FreeBorders::GetType() const
2322 return SMDSAbs_Edge;
2326 //================================================================================
2329 Description : Predicate for free Edges
2331 //================================================================================
2333 FreeEdges::FreeEdges()
2338 void FreeEdges::SetMesh( const SMDS_Mesh* theMesh )
2343 bool FreeEdges::IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId )
2345 TColStd_MapOfInteger aMap;
2346 for ( int i = 0; i < 2; i++ )
2348 SMDS_ElemIteratorPtr anElemIter = theNodes[ i ]->GetInverseElementIterator(SMDSAbs_Face);
2349 while( anElemIter->more() )
2351 if ( const SMDS_MeshElement* anElem = anElemIter->next())
2353 const int anId = anElem->GetID();
2354 if ( anId != theFaceId && !aMap.Add( anId ))
2362 bool FreeEdges::IsSatisfy( long theId )
2367 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2368 if ( aFace == 0 || aFace->GetType() != SMDSAbs_Face || aFace->NbNodes() < 3 )
2371 SMDS_ElemIteratorPtr anIter;
2372 if ( aFace->IsQuadratic() ) {
2373 anIter = dynamic_cast<const SMDS_VtkFace*>
2374 (aFace)->interlacedNodesElemIterator();
2377 anIter = aFace->nodesIterator();
2382 int i = 0, nbNodes = aFace->NbNodes();
2383 std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
2384 while( anIter->more() )
2386 const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
2389 aNodes[ i++ ] = aNode;
2391 aNodes[ nbNodes ] = aNodes[ 0 ];
2393 for ( i = 0; i < nbNodes; i++ )
2394 if ( IsFreeEdge( &aNodes[ i ], theId ) )
2400 SMDSAbs_ElementType FreeEdges::GetType() const
2402 return SMDSAbs_Face;
2405 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2408 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
2409 if(thePntId1 > thePntId2){
2410 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
2414 bool FreeEdges::Border::operator<(const FreeEdges::Border& x) const{
2415 if(myPntId[0] < x.myPntId[0]) return true;
2416 if(myPntId[0] == x.myPntId[0])
2417 if(myPntId[1] < x.myPntId[1]) return true;
2421 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2422 FreeEdges::TBorders& theRegistry,
2423 FreeEdges::TBorders& theContainer)
2425 if(theRegistry.find(theBorder) == theRegistry.end()){
2426 theRegistry.insert(theBorder);
2427 theContainer.insert(theBorder);
2429 theContainer.erase(theBorder);
2433 void FreeEdges::GetBoreders(TBorders& theBorders)
2436 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2437 for(; anIter->more(); ){
2438 const SMDS_MeshFace* anElem = anIter->next();
2439 long anElemId = anElem->GetID();
2440 SMDS_ElemIteratorPtr aNodesIter;
2441 if ( anElem->IsQuadratic() )
2442 aNodesIter = static_cast<const SMDS_VtkFace*>(anElem)->
2443 interlacedNodesElemIterator();
2445 aNodesIter = anElem->nodesIterator();
2447 const SMDS_MeshElement* aNode;
2448 if(aNodesIter->more()){
2449 aNode = aNodesIter->next();
2450 aNodeId[0] = aNodeId[1] = aNode->GetID();
2452 for(; aNodesIter->more(); ){
2453 aNode = aNodesIter->next();
2454 long anId = aNode->GetID();
2455 Border aBorder(anElemId,aNodeId[1],anId);
2457 UpdateBorders(aBorder,aRegistry,theBorders);
2459 Border aBorder(anElemId,aNodeId[0],aNodeId[1]);
2460 UpdateBorders(aBorder,aRegistry,theBorders);
2464 //================================================================================
2467 Description : Predicate for free nodes
2469 //================================================================================
2471 FreeNodes::FreeNodes()
2476 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2481 bool FreeNodes::IsSatisfy( long theNodeId )
2483 const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2487 return (aNode->NbInverseElements() < 1);
2490 SMDSAbs_ElementType FreeNodes::GetType() const
2492 return SMDSAbs_Node;
2496 //================================================================================
2499 Description : Predicate for free faces
2501 //================================================================================
2503 FreeFaces::FreeFaces()
2508 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2513 bool FreeFaces::IsSatisfy( long theId )
2515 if (!myMesh) return false;
2516 // check that faces nodes refers to less than two common volumes
2517 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2518 if ( !aFace || aFace->GetType() != SMDSAbs_Face )
2521 int nbNode = aFace->NbNodes();
2523 // collect volumes check that number of volumss with count equal nbNode not less than 2
2524 typedef map< SMDS_MeshElement*, int > TMapOfVolume; // map of volume counters
2525 typedef map< SMDS_MeshElement*, int >::iterator TItrMapOfVolume; // iterator
2526 TMapOfVolume mapOfVol;
2528 SMDS_ElemIteratorPtr nodeItr = aFace->nodesIterator();
2529 while ( nodeItr->more() ) {
2530 const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>(nodeItr->next());
2531 if ( !aNode ) continue;
2532 SMDS_ElemIteratorPtr volItr = aNode->GetInverseElementIterator(SMDSAbs_Volume);
2533 while ( volItr->more() ) {
2534 SMDS_MeshElement* aVol = (SMDS_MeshElement*)volItr->next();
2535 TItrMapOfVolume itr = mapOfVol.insert(make_pair(aVol, 0)).first;
2540 TItrMapOfVolume volItr = mapOfVol.begin();
2541 TItrMapOfVolume volEnd = mapOfVol.end();
2542 for ( ; volItr != volEnd; ++volItr )
2543 if ( (*volItr).second >= nbNode )
2545 // face is not free if number of volumes constructed on thier nodes more than one
2549 SMDSAbs_ElementType FreeFaces::GetType() const
2551 return SMDSAbs_Face;
2554 //================================================================================
2556 Class : LinearOrQuadratic
2557 Description : Predicate to verify whether a mesh element is linear
2559 //================================================================================
2561 LinearOrQuadratic::LinearOrQuadratic()
2566 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2571 bool LinearOrQuadratic::IsSatisfy( long theId )
2573 if (!myMesh) return false;
2574 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2575 if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2577 return (!anElem->IsQuadratic());
2580 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2585 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2590 //================================================================================
2593 Description : Functor for check color of group to whic mesh element belongs to
2595 //================================================================================
2597 GroupColor::GroupColor()
2601 bool GroupColor::IsSatisfy( long theId )
2603 return (myIDs.find( theId ) != myIDs.end());
2606 void GroupColor::SetType( SMDSAbs_ElementType theType )
2611 SMDSAbs_ElementType GroupColor::GetType() const
2616 static bool isEqual( const Quantity_Color& theColor1,
2617 const Quantity_Color& theColor2 )
2619 // tolerance to compare colors
2620 const double tol = 5*1e-3;
2621 return ( fabs( theColor1.Red() - theColor2.Red() ) < tol &&
2622 fabs( theColor1.Green() - theColor2.Green() ) < tol &&
2623 fabs( theColor1.Blue() - theColor2.Blue() ) < tol );
2627 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2631 const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2635 int nbGrp = aMesh->GetNbGroups();
2639 // iterates on groups and find necessary elements ids
2640 const std::set<SMESHDS_GroupBase*>& aGroups = aMesh->GetGroups();
2641 set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
2642 for (; GrIt != aGroups.end(); GrIt++) {
2643 SMESHDS_GroupBase* aGrp = (*GrIt);
2646 // check type and color of group
2647 if ( !isEqual( myColor, aGrp->GetColor() ) )
2649 if ( myType != SMDSAbs_All && myType != (SMDSAbs_ElementType)aGrp->GetType() )
2652 SMDSAbs_ElementType aGrpElType = (SMDSAbs_ElementType)aGrp->GetType();
2653 if ( myType == aGrpElType || (myType == SMDSAbs_All && aGrpElType != SMDSAbs_Node) ) {
2654 // add elements IDS into control
2655 int aSize = aGrp->Extent();
2656 for (int i = 0; i < aSize; i++)
2657 myIDs.insert( aGrp->GetID(i+1) );
2662 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
2664 TCollection_AsciiString aStr = theStr;
2665 aStr.RemoveAll( ' ' );
2666 aStr.RemoveAll( '\t' );
2667 for ( int aPos = aStr.Search( ";;" ); aPos != -1; aPos = aStr.Search( ";;" ) )
2668 aStr.Remove( aPos, 2 );
2669 Standard_Real clr[3];
2670 clr[0] = clr[1] = clr[2] = 0.;
2671 for ( int i = 0; i < 3; i++ ) {
2672 TCollection_AsciiString tmpStr = aStr.Token( ";", i+1 );
2673 if ( !tmpStr.IsEmpty() && tmpStr.IsRealValue() )
2674 clr[i] = tmpStr.RealValue();
2676 myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
2679 //=======================================================================
2680 // name : GetRangeStr
2681 // Purpose : Get range as a string.
2682 // Example: "1,2,3,50-60,63,67,70-"
2683 //=======================================================================
2685 void GroupColor::GetColorStr( TCollection_AsciiString& theResStr ) const
2688 theResStr += TCollection_AsciiString( myColor.Red() );
2689 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2690 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2693 //================================================================================
2695 Class : ElemGeomType
2696 Description : Predicate to check element geometry type
2698 //================================================================================
2700 ElemGeomType::ElemGeomType()
2703 myType = SMDSAbs_All;
2704 myGeomType = SMDSGeom_TRIANGLE;
2707 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2712 bool ElemGeomType::IsSatisfy( long theId )
2714 if (!myMesh) return false;
2715 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2718 const SMDSAbs_ElementType anElemType = anElem->GetType();
2719 if ( myType != SMDSAbs_All && anElemType != myType )
2721 bool isOk = ( anElem->GetGeomType() == myGeomType );
2725 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2730 SMDSAbs_ElementType ElemGeomType::GetType() const
2735 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2737 myGeomType = theType;
2740 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2745 //================================================================================
2747 Class : ElemEntityType
2748 Description : Predicate to check element entity type
2750 //================================================================================
2752 ElemEntityType::ElemEntityType():
2754 myType( SMDSAbs_All ),
2755 myEntityType( SMDSEntity_0D )
2759 void ElemEntityType::SetMesh( const SMDS_Mesh* theMesh )
2764 bool ElemEntityType::IsSatisfy( long theId )
2766 if ( !myMesh ) return false;
2767 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2769 myEntityType == anElem->GetEntityType() &&
2770 ( myType == SMDSAbs_Edge || myType == SMDSAbs_Face || myType == SMDSAbs_Volume ));
2773 void ElemEntityType::SetType( SMDSAbs_ElementType theType )
2778 SMDSAbs_ElementType ElemEntityType::GetType() const
2783 void ElemEntityType::SetElemEntityType( SMDSAbs_EntityType theEntityType )
2785 myEntityType = theEntityType;
2788 SMDSAbs_EntityType ElemEntityType::GetElemEntityType() const
2790 return myEntityType;
2793 //================================================================================
2795 * \brief Class CoplanarFaces
2797 //================================================================================
2799 CoplanarFaces::CoplanarFaces()
2800 : myFaceID(0), myToler(0)
2803 void CoplanarFaces::SetMesh( const SMDS_Mesh* theMesh )
2805 myMeshModifTracer.SetMesh( theMesh );
2806 if ( myMeshModifTracer.IsMeshModified() )
2808 // Build a set of coplanar face ids
2810 myCoplanarIDs.clear();
2812 if ( !myMeshModifTracer.GetMesh() || !myFaceID || !myToler )
2815 const SMDS_MeshElement* face = myMeshModifTracer.GetMesh()->FindElement( myFaceID );
2816 if ( !face || face->GetType() != SMDSAbs_Face )
2820 gp_Vec myNorm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2824 const double radianTol = myToler * M_PI / 180.;
2825 std::set< SMESH_TLink > checkedLinks;
2827 std::list< pair< const SMDS_MeshElement*, gp_Vec > > faceQueue;
2828 faceQueue.push_back( make_pair( face, myNorm ));
2829 while ( !faceQueue.empty() )
2831 face = faceQueue.front().first;
2832 myNorm = faceQueue.front().second;
2833 faceQueue.pop_front();
2835 for ( int i = 0, nbN = face->NbCornerNodes(); i < nbN; ++i )
2837 const SMDS_MeshNode* n1 = face->GetNode( i );
2838 const SMDS_MeshNode* n2 = face->GetNode(( i+1 )%nbN);
2839 if ( !checkedLinks.insert( SMESH_TLink( n1, n2 )).second )
2841 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator(SMDSAbs_Face);
2842 while ( fIt->more() )
2844 const SMDS_MeshElement* f = fIt->next();
2845 if ( f->GetNodeIndex( n2 ) > -1 )
2847 gp_Vec norm = getNormale( static_cast<const SMDS_MeshFace*>(f), &normOK );
2848 if (!normOK || myNorm.Angle( norm ) <= radianTol)
2850 myCoplanarIDs.insert( f->GetID() );
2851 faceQueue.push_back( make_pair( f, norm ));
2859 bool CoplanarFaces::IsSatisfy( long theElementId )
2861 return myCoplanarIDs.count( theElementId );
2866 *Description : Predicate for Range of Ids.
2867 * Range may be specified with two ways.
2868 * 1. Using AddToRange method
2869 * 2. With SetRangeStr method. Parameter of this method is a string
2870 * like as "1,2,3,50-60,63,67,70-"
2873 //=======================================================================
2874 // name : RangeOfIds
2875 // Purpose : Constructor
2876 //=======================================================================
2877 RangeOfIds::RangeOfIds()
2880 myType = SMDSAbs_All;
2883 //=======================================================================
2885 // Purpose : Set mesh
2886 //=======================================================================
2887 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
2892 //=======================================================================
2893 // name : AddToRange
2894 // Purpose : Add ID to the range
2895 //=======================================================================
2896 bool RangeOfIds::AddToRange( long theEntityId )
2898 myIds.Add( theEntityId );
2902 //=======================================================================
2903 // name : GetRangeStr
2904 // Purpose : Get range as a string.
2905 // Example: "1,2,3,50-60,63,67,70-"
2906 //=======================================================================
2907 void RangeOfIds::GetRangeStr( TCollection_AsciiString& theResStr )
2911 TColStd_SequenceOfInteger anIntSeq;
2912 TColStd_SequenceOfAsciiString aStrSeq;
2914 TColStd_MapIteratorOfMapOfInteger anIter( myIds );
2915 for ( ; anIter.More(); anIter.Next() )
2917 int anId = anIter.Key();
2918 TCollection_AsciiString aStr( anId );
2919 anIntSeq.Append( anId );
2920 aStrSeq.Append( aStr );
2923 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2925 int aMinId = myMin( i );
2926 int aMaxId = myMax( i );
2928 TCollection_AsciiString aStr;
2929 if ( aMinId != IntegerFirst() )
2934 if ( aMaxId != IntegerLast() )
2937 // find position of the string in result sequence and insert string in it
2938 if ( anIntSeq.Length() == 0 )
2940 anIntSeq.Append( aMinId );
2941 aStrSeq.Append( aStr );
2945 if ( aMinId < anIntSeq.First() )
2947 anIntSeq.Prepend( aMinId );
2948 aStrSeq.Prepend( aStr );
2950 else if ( aMinId > anIntSeq.Last() )
2952 anIntSeq.Append( aMinId );
2953 aStrSeq.Append( aStr );
2956 for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
2957 if ( aMinId < anIntSeq( j ) )
2959 anIntSeq.InsertBefore( j, aMinId );
2960 aStrSeq.InsertBefore( j, aStr );
2966 if ( aStrSeq.Length() == 0 )
2969 theResStr = aStrSeq( 1 );
2970 for ( int j = 2, k = aStrSeq.Length(); j <= k; j++ )
2973 theResStr += aStrSeq( j );
2977 //=======================================================================
2978 // name : SetRangeStr
2979 // Purpose : Define range with string
2980 // Example of entry string: "1,2,3,50-60,63,67,70-"
2981 //=======================================================================
2982 bool RangeOfIds::SetRangeStr( const TCollection_AsciiString& theStr )
2988 TCollection_AsciiString aStr = theStr;
2989 aStr.RemoveAll( ' ' );
2990 aStr.RemoveAll( '\t' );
2992 for ( int aPos = aStr.Search( ",," ); aPos != -1; aPos = aStr.Search( ",," ) )
2993 aStr.Remove( aPos, 2 );
2995 TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
2997 while ( tmpStr != "" )
2999 tmpStr = aStr.Token( ",", i++ );
3000 int aPos = tmpStr.Search( '-' );
3004 if ( tmpStr.IsIntegerValue() )
3005 myIds.Add( tmpStr.IntegerValue() );
3011 TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
3012 TCollection_AsciiString aMinStr = tmpStr;
3014 while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
3015 while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
3017 if ( (!aMinStr.IsEmpty() && !aMinStr.IsIntegerValue()) ||
3018 (!aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue()) )
3021 myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
3022 myMax.Append( aMaxStr.IsEmpty() ? IntegerLast() : aMaxStr.IntegerValue() );
3029 //=======================================================================
3031 // Purpose : Get type of supported entities
3032 //=======================================================================
3033 SMDSAbs_ElementType RangeOfIds::GetType() const
3038 //=======================================================================
3040 // Purpose : Set type of supported entities
3041 //=======================================================================
3042 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
3047 //=======================================================================
3049 // Purpose : Verify whether entity satisfies to this rpedicate
3050 //=======================================================================
3051 bool RangeOfIds::IsSatisfy( long theId )
3056 if ( myType == SMDSAbs_Node )
3058 if ( myMesh->FindNode( theId ) == 0 )
3063 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
3064 if ( anElem == 0 || (myType != anElem->GetType() && myType != SMDSAbs_All ))
3068 if ( myIds.Contains( theId ) )
3071 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
3072 if ( theId >= myMin( i ) && theId <= myMax( i ) )
3080 Description : Base class for comparators
3082 Comparator::Comparator():
3086 Comparator::~Comparator()
3089 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
3092 myFunctor->SetMesh( theMesh );
3095 void Comparator::SetMargin( double theValue )
3097 myMargin = theValue;
3100 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
3102 myFunctor = theFunct;
3105 SMDSAbs_ElementType Comparator::GetType() const
3107 return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
3110 double Comparator::GetMargin()
3118 Description : Comparator "<"
3120 bool LessThan::IsSatisfy( long theId )
3122 return myFunctor && myFunctor->GetValue( theId ) < myMargin;
3128 Description : Comparator ">"
3130 bool MoreThan::IsSatisfy( long theId )
3132 return myFunctor && myFunctor->GetValue( theId ) > myMargin;
3138 Description : Comparator "="
3141 myToler(Precision::Confusion())
3144 bool EqualTo::IsSatisfy( long theId )
3146 return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
3149 void EqualTo::SetTolerance( double theToler )
3154 double EqualTo::GetTolerance()
3161 Description : Logical NOT predicate
3163 LogicalNOT::LogicalNOT()
3166 LogicalNOT::~LogicalNOT()
3169 bool LogicalNOT::IsSatisfy( long theId )
3171 return myPredicate && !myPredicate->IsSatisfy( theId );
3174 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
3177 myPredicate->SetMesh( theMesh );
3180 void LogicalNOT::SetPredicate( PredicatePtr thePred )
3182 myPredicate = thePred;
3185 SMDSAbs_ElementType LogicalNOT::GetType() const
3187 return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
3192 Class : LogicalBinary
3193 Description : Base class for binary logical predicate
3195 LogicalBinary::LogicalBinary()
3198 LogicalBinary::~LogicalBinary()
3201 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
3204 myPredicate1->SetMesh( theMesh );
3207 myPredicate2->SetMesh( theMesh );
3210 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
3212 myPredicate1 = thePredicate;
3215 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
3217 myPredicate2 = thePredicate;
3220 SMDSAbs_ElementType LogicalBinary::GetType() const
3222 if ( !myPredicate1 || !myPredicate2 )
3225 SMDSAbs_ElementType aType1 = myPredicate1->GetType();
3226 SMDSAbs_ElementType aType2 = myPredicate2->GetType();
3228 return aType1 == aType2 ? aType1 : SMDSAbs_All;
3234 Description : Logical AND
3236 bool LogicalAND::IsSatisfy( long theId )
3241 myPredicate1->IsSatisfy( theId ) &&
3242 myPredicate2->IsSatisfy( theId );
3248 Description : Logical OR
3250 bool LogicalOR::IsSatisfy( long theId )
3255 (myPredicate1->IsSatisfy( theId ) ||
3256 myPredicate2->IsSatisfy( theId ));
3270 void Filter::SetPredicate( PredicatePtr thePredicate )
3272 myPredicate = thePredicate;
3275 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3276 PredicatePtr thePredicate,
3277 TIdSequence& theSequence )
3279 theSequence.clear();
3281 if ( !theMesh || !thePredicate )
3284 thePredicate->SetMesh( theMesh );
3286 SMDS_ElemIteratorPtr elemIt = theMesh->elementsIterator( thePredicate->GetType() );
3288 while ( elemIt->more() ) {
3289 const SMDS_MeshElement* anElem = elemIt->next();
3290 long anId = anElem->GetID();
3291 if ( thePredicate->IsSatisfy( anId ) )
3292 theSequence.push_back( anId );
3297 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3298 Filter::TIdSequence& theSequence )
3300 GetElementsId(theMesh,myPredicate,theSequence);
3307 typedef std::set<SMDS_MeshFace*> TMapOfFacePtr;
3313 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
3314 SMDS_MeshNode* theNode2 )
3320 ManifoldPart::Link::~Link()
3326 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
3328 if ( myNode1 == theLink.myNode1 &&
3329 myNode2 == theLink.myNode2 )
3331 else if ( myNode1 == theLink.myNode2 &&
3332 myNode2 == theLink.myNode1 )
3338 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
3340 if(myNode1 < x.myNode1) return true;
3341 if(myNode1 == x.myNode1)
3342 if(myNode2 < x.myNode2) return true;
3346 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
3347 const ManifoldPart::Link& theLink2 )
3349 return theLink1.IsEqual( theLink2 );
3352 ManifoldPart::ManifoldPart()
3355 myAngToler = Precision::Angular();
3356 myIsOnlyManifold = true;
3359 ManifoldPart::~ManifoldPart()
3364 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
3370 SMDSAbs_ElementType ManifoldPart::GetType() const
3371 { return SMDSAbs_Face; }
3373 bool ManifoldPart::IsSatisfy( long theElementId )
3375 return myMapIds.Contains( theElementId );
3378 void ManifoldPart::SetAngleTolerance( const double theAngToler )
3379 { myAngToler = theAngToler; }
3381 double ManifoldPart::GetAngleTolerance() const
3382 { return myAngToler; }
3384 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
3385 { myIsOnlyManifold = theIsOnly; }
3387 void ManifoldPart::SetStartElem( const long theStartId )
3388 { myStartElemId = theStartId; }
3390 bool ManifoldPart::process()
3393 myMapBadGeomIds.Clear();
3395 myAllFacePtr.clear();
3396 myAllFacePtrIntDMap.clear();
3400 // collect all faces into own map
3401 SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3402 for (; anFaceItr->more(); )
3404 SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3405 myAllFacePtr.push_back( aFacePtr );
3406 myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3409 SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3413 // the map of non manifold links and bad geometry
3414 TMapOfLink aMapOfNonManifold;
3415 TColStd_MapOfInteger aMapOfTreated;
3417 // begin cycle on faces from start index and run on vector till the end
3418 // and from begin to start index to cover whole vector
3419 const int aStartIndx = myAllFacePtrIntDMap[aStartFace];
3420 bool isStartTreat = false;
3421 for ( int fi = aStartIndx; !isStartTreat || fi != aStartIndx ; fi++ )
3423 if ( fi == aStartIndx )
3424 isStartTreat = true;
3425 // as result next time when fi will be equal to aStartIndx
3427 SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3428 if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3431 aMapOfTreated.Add( aFacePtr->GetID() );
3432 TColStd_MapOfInteger aResFaces;
3433 if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3434 aMapOfNonManifold, aResFaces ) )
3436 TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3437 for ( ; anItr.More(); anItr.Next() )
3439 int aFaceId = anItr.Key();
3440 aMapOfTreated.Add( aFaceId );
3441 myMapIds.Add( aFaceId );
3444 if ( fi == ( myAllFacePtr.size() - 1 ) )
3446 } // end run on vector of faces
3447 return !myMapIds.IsEmpty();
3450 static void getLinks( const SMDS_MeshFace* theFace,
3451 ManifoldPart::TVectorOfLink& theLinks )
3453 int aNbNode = theFace->NbNodes();
3454 SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3456 SMDS_MeshNode* aNode = 0;
3457 for ( ; aNodeItr->more() && i <= aNbNode; )
3460 SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3464 SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3466 ManifoldPart::Link aLink( aN1, aN2 );
3467 theLinks.push_back( aLink );
3471 bool ManifoldPart::findConnected
3472 ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3473 SMDS_MeshFace* theStartFace,
3474 ManifoldPart::TMapOfLink& theNonManifold,
3475 TColStd_MapOfInteger& theResFaces )
3477 theResFaces.Clear();
3478 if ( !theAllFacePtrInt.size() )
3481 if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3483 myMapBadGeomIds.Add( theStartFace->GetID() );
3487 ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3488 ManifoldPart::TVectorOfLink aSeqOfBoundary;
3489 theResFaces.Add( theStartFace->GetID() );
3490 ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3492 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3493 aDMapLinkFace, theNonManifold, theStartFace );
3495 bool isDone = false;
3496 while ( !isDone && aMapOfBoundary.size() != 0 )
3498 bool isToReset = false;
3499 ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3500 for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3502 ManifoldPart::Link aLink = *pLink;
3503 if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3505 // each link could be treated only once
3506 aMapToSkip.insert( aLink );
3508 ManifoldPart::TVectorOfFacePtr aFaces;
3510 if ( myIsOnlyManifold &&
3511 (theNonManifold.find( aLink ) != theNonManifold.end()) )
3515 getFacesByLink( aLink, aFaces );
3516 // filter the element to keep only indicated elements
3517 ManifoldPart::TVectorOfFacePtr aFiltered;
3518 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3519 for ( ; pFace != aFaces.end(); ++pFace )
3521 SMDS_MeshFace* aFace = *pFace;
3522 if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3523 aFiltered.push_back( aFace );
3526 if ( aFaces.size() < 2 ) // no neihgbour faces
3528 else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3530 theNonManifold.insert( aLink );
3535 // compare normal with normals of neighbor element
3536 SMDS_MeshFace* aPrevFace = aDMapLinkFace[ aLink ];
3537 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3538 for ( ; pFace != aFaces.end(); ++pFace )
3540 SMDS_MeshFace* aNextFace = *pFace;
3541 if ( aPrevFace == aNextFace )
3543 int anNextFaceID = aNextFace->GetID();
3544 if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3545 // should not be with non manifold restriction. probably bad topology
3547 // check if face was treated and skipped
3548 if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3549 !isInPlane( aPrevFace, aNextFace ) )
3551 // add new element to connected and extend the boundaries.
3552 theResFaces.Add( anNextFaceID );
3553 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3554 aDMapLinkFace, theNonManifold, aNextFace );
3558 isDone = !isToReset;
3561 return !theResFaces.IsEmpty();
3564 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3565 const SMDS_MeshFace* theFace2 )
3567 gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3568 gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3569 if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3571 myMapBadGeomIds.Add( theFace2->GetID() );
3574 if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
3580 void ManifoldPart::expandBoundary
3581 ( ManifoldPart::TMapOfLink& theMapOfBoundary,
3582 ManifoldPart::TVectorOfLink& theSeqOfBoundary,
3583 ManifoldPart::TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
3584 ManifoldPart::TMapOfLink& theNonManifold,
3585 SMDS_MeshFace* theNextFace ) const
3587 ManifoldPart::TVectorOfLink aLinks;
3588 getLinks( theNextFace, aLinks );
3589 int aNbLink = (int)aLinks.size();
3590 for ( int i = 0; i < aNbLink; i++ )
3592 ManifoldPart::Link aLink = aLinks[ i ];
3593 if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3595 if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3597 if ( myIsOnlyManifold )
3599 // remove from boundary
3600 theMapOfBoundary.erase( aLink );
3601 ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3602 for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3604 ManifoldPart::Link aBoundLink = *pLink;
3605 if ( aBoundLink.IsEqual( aLink ) )
3607 theSeqOfBoundary.erase( pLink );
3615 theMapOfBoundary.insert( aLink );
3616 theSeqOfBoundary.push_back( aLink );
3617 theDMapLinkFacePtr[ aLink ] = theNextFace;
3622 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3623 ManifoldPart::TVectorOfFacePtr& theFaces ) const
3625 std::set<SMDS_MeshCell *> aSetOfFaces;
3626 // take all faces that shared first node
3627 SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3628 for ( ; anItr->more(); )
3630 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3633 aSetOfFaces.insert( aFace );
3635 // take all faces that shared second node
3636 anItr = theLink.myNode2->facesIterator();
3637 // find the common part of two sets
3638 for ( ; anItr->more(); )
3640 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3641 if ( aSetOfFaces.count( aFace ) )
3642 theFaces.push_back( aFace );
3651 ElementsOnSurface::ElementsOnSurface()
3654 myType = SMDSAbs_All;
3656 myToler = Precision::Confusion();
3657 myUseBoundaries = false;
3660 ElementsOnSurface::~ElementsOnSurface()
3664 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
3666 myMeshModifTracer.SetMesh( theMesh );
3667 if ( myMeshModifTracer.IsMeshModified())
3671 bool ElementsOnSurface::IsSatisfy( long theElementId )
3673 return myIds.Contains( theElementId );
3676 SMDSAbs_ElementType ElementsOnSurface::GetType() const
3679 void ElementsOnSurface::SetTolerance( const double theToler )
3681 if ( myToler != theToler )
3686 double ElementsOnSurface::GetTolerance() const
3689 void ElementsOnSurface::SetUseBoundaries( bool theUse )
3691 if ( myUseBoundaries != theUse ) {
3692 myUseBoundaries = theUse;
3693 SetSurface( mySurf, myType );
3697 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
3698 const SMDSAbs_ElementType theType )
3703 if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
3705 mySurf = TopoDS::Face( theShape );
3706 BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
3708 u1 = SA.FirstUParameter(),
3709 u2 = SA.LastUParameter(),
3710 v1 = SA.FirstVParameter(),
3711 v2 = SA.LastVParameter();
3712 Handle(Geom_Surface) surf = BRep_Tool::Surface( mySurf );
3713 myProjector.Init( surf, u1,u2, v1,v2 );
3717 void ElementsOnSurface::process()
3720 if ( mySurf.IsNull() )
3723 if ( !myMeshModifTracer.GetMesh() )
3726 myIds.ReSize( myMeshModifTracer.GetMesh()->GetMeshInfo().NbElements( myType ));
3728 SMDS_ElemIteratorPtr anIter = myMeshModifTracer.GetMesh()->elementsIterator( myType );
3729 for(; anIter->more(); )
3730 process( anIter->next() );
3733 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
3735 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3736 bool isSatisfy = true;
3737 for ( ; aNodeItr->more(); )
3739 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3740 if ( !isOnSurface( aNode ) )
3747 myIds.Add( theElemPtr->GetID() );
3750 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
3752 if ( mySurf.IsNull() )
3755 gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
3756 // double aToler2 = myToler * myToler;
3757 // if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
3759 // gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
3760 // if ( aPln.SquareDistance( aPnt ) > aToler2 )
3763 // else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
3765 // gp_Cylinder aCyl = Handle(Geom_CylindricalSurface)::DownCast(mySurf)->Cylinder();
3766 // double aRad = aCyl.Radius();
3767 // gp_Ax3 anAxis = aCyl.Position();
3768 // gp_XYZ aLoc = aCyl.Location().XYZ();
3769 // double aXDist = anAxis.XDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3770 // double aYDist = anAxis.YDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3771 // if ( fabs(aXDist*aXDist + aYDist*aYDist - aRad*aRad) > aToler2 )
3776 myProjector.Perform( aPnt );
3777 bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
3787 ElementsOnShape::ElementsOnShape()
3789 myType(SMDSAbs_All),
3790 myToler(Precision::Confusion()),
3791 myAllNodesFlag(false)
3795 ElementsOnShape::~ElementsOnShape()
3800 SMDSAbs_ElementType ElementsOnShape::GetType() const
3805 void ElementsOnShape::SetTolerance (const double theToler)
3807 if (myToler != theToler) {
3809 SetShape(myShape, myType);
3813 double ElementsOnShape::GetTolerance() const
3818 void ElementsOnShape::SetAllNodes (bool theAllNodes)
3820 myAllNodesFlag = theAllNodes;
3823 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
3828 void ElementsOnShape::SetShape (const TopoDS_Shape& theShape,
3829 const SMDSAbs_ElementType theType)
3833 if ( myShape.IsNull() ) return;
3835 TopTools_IndexedMapOfShape shapesMap;
3836 TopAbs_ShapeEnum shapeTypes[4] = { TopAbs_SOLID, TopAbs_FACE, TopAbs_EDGE, TopAbs_VERTEX };
3837 TopExp_Explorer sub;
3838 for ( int i = 0; i < 4; ++i )
3840 if ( shapesMap.IsEmpty() )
3841 for ( sub.Init( myShape, shapeTypes[i] ); sub.More(); sub.Next() )
3842 shapesMap.Add( sub.Current() );
3844 for ( sub.Init( myShape, shapeTypes[i], shapeTypes[i-1] ); sub.More(); sub.Next() )
3845 shapesMap.Add( sub.Current() );
3849 myClassifiers.resize( shapesMap.Extent() );
3850 for ( int i = 0; i < shapesMap.Extent(); ++i )
3851 myClassifiers[ i ] = new TClassifier( shapesMap( i+1 ), myToler );
3854 void ElementsOnShape::clearClassifiers()
3856 for ( size_t i = 0; i < myClassifiers.size(); ++i )
3857 delete myClassifiers[ i ];
3858 myClassifiers.clear();
3861 bool ElementsOnShape::IsSatisfy (long elemId)
3863 const SMDS_MeshElement* elem =
3864 ( myType == SMDSAbs_Node ? myMesh->FindNode( elemId ) : myMesh->FindElement( elemId ));
3865 if ( !elem || myClassifiers.empty() )
3868 for ( size_t i = 0; i < myClassifiers.size(); ++i )
3870 SMDS_ElemIteratorPtr aNodeItr = elem->nodesIterator();
3871 bool isSatisfy = myAllNodesFlag;
3873 gp_XYZ centerXYZ (0, 0, 0);
3875 while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
3877 SMESH_TNodeXYZ aPnt ( aNodeItr->next() );
3879 isSatisfy = ! myClassifiers[i]->IsOut( aPnt );
3882 // Check the center point for volumes MantisBug 0020168
3885 myClassifiers[i]->ShapeType() == TopAbs_SOLID)
3887 centerXYZ /= elem->NbNodes();
3888 isSatisfy = ! myClassifiers[i]->IsOut( centerXYZ );
3897 TopAbs_ShapeEnum ElementsOnShape::TClassifier::ShapeType() const
3899 return myShape.ShapeType();
3902 bool ElementsOnShape::TClassifier::IsOut(const gp_Pnt& p)
3904 return (this->*myIsOutFun)( p );
3907 void ElementsOnShape::TClassifier::Init (const TopoDS_Shape& theShape, double theTol)
3911 switch ( myShape.ShapeType() )
3913 case TopAbs_SOLID: {
3914 mySolidClfr.Load(theShape);
3915 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfSolid;
3919 Standard_Real u1,u2,v1,v2;
3920 Handle(Geom_Surface) surf = BRep_Tool::Surface( TopoDS::Face( theShape ));
3921 surf->Bounds( u1,u2,v1,v2 );
3922 myProjFace.Init(surf, u1,u2, v1,v2, myTol );
3923 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfFace;
3927 Standard_Real u1, u2;
3928 Handle(Geom_Curve) curve = BRep_Tool::Curve( TopoDS::Edge(theShape), u1, u2);
3929 myProjEdge.Init(curve, u1, u2);
3930 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfEdge;
3933 case TopAbs_VERTEX:{
3934 myVertexXYZ = BRep_Tool::Pnt( TopoDS::Vertex( theShape ) );
3935 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfVertex;
3939 throw SALOME_Exception("Programmer error in usage of ElementsOnShape::TClassifier");
3943 bool ElementsOnShape::TClassifier::isOutOfSolid (const gp_Pnt& p)
3945 mySolidClfr.Perform( p, myTol );
3946 return ( mySolidClfr.State() != TopAbs_IN && mySolidClfr.State() != TopAbs_ON );
3949 bool ElementsOnShape::TClassifier::isOutOfFace (const gp_Pnt& p)
3951 myProjFace.Perform( p );
3952 if ( myProjFace.IsDone() && myProjFace.LowerDistance() <= myTol )
3954 // check relatively to the face
3955 Quantity_Parameter u, v;
3956 myProjFace.LowerDistanceParameters(u, v);
3957 gp_Pnt2d aProjPnt (u, v);
3958 BRepClass_FaceClassifier aClsf ( TopoDS::Face( myShape ), aProjPnt, myTol );
3959 if ( aClsf.State() == TopAbs_IN || aClsf.State() == TopAbs_ON )
3965 bool ElementsOnShape::TClassifier::isOutOfEdge (const gp_Pnt& p)
3967 myProjEdge.Perform( p );
3968 return ! ( myProjEdge.NbPoints() > 0 && myProjEdge.LowerDistance() <= myTol );
3971 bool ElementsOnShape::TClassifier::isOutOfVertex(const gp_Pnt& p)
3973 return ( myVertexXYZ.Distance( p ) > myTol );
3977 TSequenceOfXYZ::TSequenceOfXYZ()
3980 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n)
3983 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t)
3986 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray)
3989 template <class InputIterator>
3990 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd)
3993 TSequenceOfXYZ::~TSequenceOfXYZ()
3996 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
3998 myArray = theSequenceOfXYZ.myArray;
4002 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
4004 return myArray[n-1];
4007 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
4009 return myArray[n-1];
4012 void TSequenceOfXYZ::clear()
4017 void TSequenceOfXYZ::reserve(size_type n)
4022 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
4024 myArray.push_back(v);
4027 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
4029 return myArray.size();
4032 TMeshModifTracer::TMeshModifTracer():
4033 myMeshModifTime(0), myMesh(0)
4036 void TMeshModifTracer::SetMesh( const SMDS_Mesh* theMesh )
4038 if ( theMesh != myMesh )
4039 myMeshModifTime = 0;
4042 bool TMeshModifTracer::IsMeshModified()
4044 bool modified = false;
4047 modified = ( myMeshModifTime != myMesh->GetMTime() );
4048 myMeshModifTime = myMesh->GetMTime();