1 // Copyright (C) 2007-2011 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"
27 #include <BRepAdaptor_Surface.hxx>
28 #include <BRepClass_FaceClassifier.hxx>
29 #include <BRep_Tool.hxx>
33 #include <TopoDS_Edge.hxx>
34 #include <TopoDS_Face.hxx>
35 #include <TopoDS_Shape.hxx>
36 #include <TopoDS_Vertex.hxx>
37 #include <TopoDS_Iterator.hxx>
39 #include <Geom_CylindricalSurface.hxx>
40 #include <Geom_Plane.hxx>
41 #include <Geom_Surface.hxx>
43 #include <Precision.hxx>
44 #include <TColStd_MapIteratorOfMapOfInteger.hxx>
45 #include <TColStd_MapOfInteger.hxx>
46 #include <TColStd_SequenceOfAsciiString.hxx>
47 #include <TColgp_Array1OfXYZ.hxx>
50 #include <gp_Cylinder.hxx>
57 #include "SMDS_Mesh.hxx"
58 #include "SMDS_Iterator.hxx"
59 #include "SMDS_MeshElement.hxx"
60 #include "SMDS_MeshNode.hxx"
61 #include "SMDS_VolumeTool.hxx"
62 #include "SMDS_QuadraticFaceOfNodes.hxx"
63 #include "SMDS_QuadraticEdge.hxx"
65 #include "SMESHDS_Mesh.hxx"
66 #include "SMESHDS_GroupBase.hxx"
68 #include "SMESH_OctreeNode.hxx"
70 #include <vtkMeshQuality.h>
78 inline gp_XYZ gpXYZ(const SMDS_MeshNode* aNode )
80 return gp_XYZ(aNode->X(), aNode->Y(), aNode->Z() );
83 inline double getAngle( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
85 gp_Vec v1( P1 - P2 ), v2( P3 - P2 );
87 return v1.Magnitude() < gp::Resolution() ||
88 v2.Magnitude() < gp::Resolution() ? 0 : v1.Angle( v2 );
91 inline double getArea( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
93 gp_Vec aVec1( P2 - P1 );
94 gp_Vec aVec2( P3 - P1 );
95 return ( aVec1 ^ aVec2 ).Magnitude() * 0.5;
98 inline double getArea( const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3 )
100 return getArea( P1.XYZ(), P2.XYZ(), P3.XYZ() );
105 inline double getDistance( const gp_XYZ& P1, const gp_XYZ& P2 )
107 double aDist = gp_Pnt( P1 ).Distance( gp_Pnt( P2 ) );
111 int getNbMultiConnection( const SMDS_Mesh* theMesh, const int theId )
116 const SMDS_MeshElement* anEdge = theMesh->FindElement( theId );
117 if ( anEdge == 0 || anEdge->GetType() != SMDSAbs_Edge/* || anEdge->NbNodes() != 2 */)
120 // for each pair of nodes in anEdge (there are 2 pairs in a quadratic edge)
121 // count elements containing both nodes of the pair.
122 // Note that there may be such cases for a quadratic edge (a horizontal line):
127 // +-----+------+ +-----+------+
130 // result sould be 2 in both cases
132 int aResult0 = 0, aResult1 = 0;
133 // last node, it is a medium one in a quadratic edge
134 const SMDS_MeshNode* aLastNode = anEdge->GetNode( anEdge->NbNodes() - 1 );
135 const SMDS_MeshNode* aNode0 = anEdge->GetNode( 0 );
136 const SMDS_MeshNode* aNode1 = anEdge->GetNode( 1 );
137 if ( aNode1 == aLastNode ) aNode1 = 0;
139 SMDS_ElemIteratorPtr anElemIter = aLastNode->GetInverseElementIterator();
140 while( anElemIter->more() ) {
141 const SMDS_MeshElement* anElem = anElemIter->next();
142 if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
143 SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
144 while ( anIter->more() ) {
145 if ( const SMDS_MeshElement* anElemNode = anIter->next() ) {
146 if ( anElemNode == aNode0 ) {
148 if ( !aNode1 ) break; // not a quadratic edge
150 else if ( anElemNode == aNode1 )
156 int aResult = std::max ( aResult0, aResult1 );
158 // TColStd_MapOfInteger aMap;
160 // SMDS_ElemIteratorPtr anIter = anEdge->nodesIterator();
161 // if ( anIter != 0 ) {
162 // while( anIter->more() ) {
163 // const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
166 // SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
167 // while( anElemIter->more() ) {
168 // const SMDS_MeshElement* anElem = anElemIter->next();
169 // if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
170 // int anId = anElem->GetID();
172 // if ( anIter->more() ) // i.e. first node
174 // else if ( aMap.Contains( anId ) )
184 gp_XYZ getNormale( const SMDS_MeshFace* theFace, bool* ok=0 )
186 int aNbNode = theFace->NbNodes();
188 gp_XYZ q1 = gpXYZ( theFace->GetNode(1)) - gpXYZ( theFace->GetNode(0));
189 gp_XYZ q2 = gpXYZ( theFace->GetNode(2)) - gpXYZ( theFace->GetNode(0));
192 gp_XYZ q3 = gpXYZ( theFace->GetNode(3)) - gpXYZ( theFace->GetNode(0));
195 double len = n.Modulus();
196 bool zeroLen = ( len <= numeric_limits<double>::min());
200 if (ok) *ok = !zeroLen;
208 using namespace SMESH::Controls;
215 Class : NumericalFunctor
216 Description : Base class for numerical functors
218 NumericalFunctor::NumericalFunctor():
224 void NumericalFunctor::SetMesh( const SMDS_Mesh* theMesh )
229 bool NumericalFunctor::GetPoints(const int theId,
230 TSequenceOfXYZ& theRes ) const
237 return GetPoints( myMesh->FindElement( theId ), theRes );
240 bool NumericalFunctor::GetPoints(const SMDS_MeshElement* anElem,
241 TSequenceOfXYZ& theRes )
248 theRes.reserve( anElem->NbNodes() );
250 // Get nodes of the element
251 SMDS_ElemIteratorPtr anIter;
253 if ( anElem->IsQuadratic() ) {
254 switch ( anElem->GetType() ) {
256 anIter = dynamic_cast<const SMDS_VtkEdge*>
257 (anElem)->interlacedNodesElemIterator();
260 anIter = dynamic_cast<const SMDS_VtkFace*>
261 (anElem)->interlacedNodesElemIterator();
264 anIter = anElem->nodesIterator();
269 anIter = anElem->nodesIterator();
273 while( anIter->more() ) {
274 if ( const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>( anIter->next() ))
275 theRes.push_back( gp_XYZ( aNode->X(), aNode->Y(), aNode->Z() ) );
282 long NumericalFunctor::GetPrecision() const
287 void NumericalFunctor::SetPrecision( const long thePrecision )
289 myPrecision = thePrecision;
290 myPrecisionValue = pow( 10., (double)( myPrecision ) );
293 double NumericalFunctor::GetValue( long theId )
297 myCurrElement = myMesh->FindElement( theId );
300 if ( GetPoints( theId, P ))
301 aVal = Round( GetValue( P ));
306 double NumericalFunctor::Round( const double & aVal )
308 return ( myPrecision >= 0 ) ? floor( aVal * myPrecisionValue + 0.5 ) / myPrecisionValue : aVal;
311 //================================================================================
313 * \brief Return histogram of functor values
314 * \param nbIntervals - number of intervals
315 * \param nbEvents - number of mesh elements having values within i-th interval
316 * \param funValues - boundaries of intervals
317 * \param elements - elements to check vulue of; empty list means "of all"
318 * \param minmax - boundaries of diapason of values to divide into intervals
320 //================================================================================
322 void NumericalFunctor::GetHistogram(int nbIntervals,
323 std::vector<int>& nbEvents,
324 std::vector<double>& funValues,
325 const vector<int>& elements,
326 const double* minmax)
328 if ( nbIntervals < 1 ||
330 !myMesh->GetMeshInfo().NbElements( GetType() ))
332 nbEvents.resize( nbIntervals, 0 );
333 funValues.resize( nbIntervals+1 );
335 // get all values sorted
336 std::multiset< double > values;
337 if ( elements.empty() )
339 SMDS_ElemIteratorPtr elemIt = myMesh->elementsIterator(GetType());
340 while ( elemIt->more() )
341 values.insert( GetValue( elemIt->next()->GetID() ));
345 vector<int>::const_iterator id = elements.begin();
346 for ( ; id != elements.end(); ++id )
347 values.insert( GetValue( *id ));
352 funValues[0] = minmax[0];
353 funValues[nbIntervals] = minmax[1];
357 funValues[0] = *values.begin();
358 funValues[nbIntervals] = *values.rbegin();
360 // case nbIntervals == 1
361 if ( nbIntervals == 1 )
363 nbEvents[0] = values.size();
367 if (funValues.front() == funValues.back())
369 nbEvents.resize( 1 );
370 nbEvents[0] = values.size();
371 funValues[1] = funValues.back();
372 funValues.resize( 2 );
375 std::multiset< double >::iterator min = values.begin(), max;
376 for ( int i = 0; i < nbIntervals; ++i )
378 // find end value of i-th interval
379 double r = (i+1) / double( nbIntervals );
380 funValues[i+1] = funValues.front() * (1-r) + funValues.back() * r;
382 // count values in the i-th interval if there are any
383 if ( min != values.end() && *min <= funValues[i+1] )
385 // find the first value out of the interval
386 max = values.upper_bound( funValues[i+1] ); // max is greater than funValues[i+1], or end()
387 nbEvents[i] = std::distance( min, max );
391 // add values larger than minmax[1]
392 nbEvents.back() += std::distance( min, values.end() );
395 //=======================================================================
396 //function : GetValue
398 //=======================================================================
400 double Volume::GetValue( long theElementId )
402 if ( theElementId && myMesh ) {
403 SMDS_VolumeTool aVolumeTool;
404 if ( aVolumeTool.Set( myMesh->FindElement( theElementId )))
405 return aVolumeTool.GetSize();
410 //=======================================================================
411 //function : GetBadRate
412 //purpose : meaningless as it is not quality control functor
413 //=======================================================================
415 double Volume::GetBadRate( double Value, int /*nbNodes*/ ) const
420 //=======================================================================
423 //=======================================================================
425 SMDSAbs_ElementType Volume::GetType() const
427 return SMDSAbs_Volume;
432 Class : MaxElementLength2D
433 Description : Functor calculating maximum length of 2D element
436 double MaxElementLength2D::GetValue( long theElementId )
439 if( GetPoints( theElementId, P ) ) {
441 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
442 SMDSAbs_ElementType aType = aElem->GetType();
446 if( len == 3 ) { // triangles
447 double L1 = getDistance(P( 1 ),P( 2 ));
448 double L2 = getDistance(P( 2 ),P( 3 ));
449 double L3 = getDistance(P( 3 ),P( 1 ));
450 aVal = Max(L1,Max(L2,L3));
453 else if( len == 4 ) { // quadrangles
454 double L1 = getDistance(P( 1 ),P( 2 ));
455 double L2 = getDistance(P( 2 ),P( 3 ));
456 double L3 = getDistance(P( 3 ),P( 4 ));
457 double L4 = getDistance(P( 4 ),P( 1 ));
458 double D1 = getDistance(P( 1 ),P( 3 ));
459 double D2 = getDistance(P( 2 ),P( 4 ));
460 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
463 else if( len == 6 ) { // quadratic triangles
464 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
465 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
466 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
467 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));
482 if( myPrecision >= 0 )
484 double prec = pow( 10., (double)myPrecision );
485 aVal = floor( aVal * prec + 0.5 ) / prec;
492 double MaxElementLength2D::GetBadRate( double Value, int /*nbNodes*/ ) const
497 SMDSAbs_ElementType MaxElementLength2D::GetType() const
503 Class : MaxElementLength3D
504 Description : Functor calculating maximum length of 3D element
507 double MaxElementLength3D::GetValue( long theElementId )
510 if( GetPoints( theElementId, P ) ) {
512 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
513 SMDSAbs_ElementType aType = aElem->GetType();
517 if( len == 4 ) { // tetras
518 double L1 = getDistance(P( 1 ),P( 2 ));
519 double L2 = getDistance(P( 2 ),P( 3 ));
520 double L3 = getDistance(P( 3 ),P( 1 ));
521 double L4 = getDistance(P( 1 ),P( 4 ));
522 double L5 = getDistance(P( 2 ),P( 4 ));
523 double L6 = getDistance(P( 3 ),P( 4 ));
524 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
527 else if( len == 5 ) { // pyramids
528 double L1 = getDistance(P( 1 ),P( 2 ));
529 double L2 = getDistance(P( 2 ),P( 3 ));
530 double L3 = getDistance(P( 3 ),P( 4 ));
531 double L4 = getDistance(P( 4 ),P( 1 ));
532 double L5 = getDistance(P( 1 ),P( 5 ));
533 double L6 = getDistance(P( 2 ),P( 5 ));
534 double L7 = getDistance(P( 3 ),P( 5 ));
535 double L8 = getDistance(P( 4 ),P( 5 ));
536 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
537 aVal = Max(aVal,Max(L7,L8));
540 else if( len == 6 ) { // pentas
541 double L1 = getDistance(P( 1 ),P( 2 ));
542 double L2 = getDistance(P( 2 ),P( 3 ));
543 double L3 = getDistance(P( 3 ),P( 1 ));
544 double L4 = getDistance(P( 4 ),P( 5 ));
545 double L5 = getDistance(P( 5 ),P( 6 ));
546 double L6 = getDistance(P( 6 ),P( 4 ));
547 double L7 = getDistance(P( 1 ),P( 4 ));
548 double L8 = getDistance(P( 2 ),P( 5 ));
549 double L9 = getDistance(P( 3 ),P( 6 ));
550 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
551 aVal = Max(aVal,Max(Max(L7,L8),L9));
554 else if( len == 8 ) { // hexas
555 double L1 = getDistance(P( 1 ),P( 2 ));
556 double L2 = getDistance(P( 2 ),P( 3 ));
557 double L3 = getDistance(P( 3 ),P( 4 ));
558 double L4 = getDistance(P( 4 ),P( 1 ));
559 double L5 = getDistance(P( 5 ),P( 6 ));
560 double L6 = getDistance(P( 6 ),P( 7 ));
561 double L7 = getDistance(P( 7 ),P( 8 ));
562 double L8 = getDistance(P( 8 ),P( 5 ));
563 double L9 = getDistance(P( 1 ),P( 5 ));
564 double L10= getDistance(P( 2 ),P( 6 ));
565 double L11= getDistance(P( 3 ),P( 7 ));
566 double L12= getDistance(P( 4 ),P( 8 ));
567 double D1 = getDistance(P( 1 ),P( 7 ));
568 double D2 = getDistance(P( 2 ),P( 8 ));
569 double D3 = getDistance(P( 3 ),P( 5 ));
570 double D4 = getDistance(P( 4 ),P( 6 ));
571 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
572 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
573 aVal = Max(aVal,Max(L11,L12));
574 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
577 else if( len == 12 ) { // hexagonal prism
578 for ( int i1 = 1; i1 < 12; ++i1 )
579 for ( int i2 = i1+1; i1 <= 12; ++i1 )
580 aVal = Max( aVal, getDistance(P( i1 ),P( i2 )));
583 else if( len == 10 ) { // quadratic tetras
584 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
585 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
586 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
587 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
588 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
589 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
590 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
593 else if( len == 13 ) { // quadratic pyramids
594 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
595 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
596 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
597 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
598 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
599 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
600 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
601 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
602 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
603 aVal = Max(aVal,Max(L7,L8));
606 else if( len == 15 ) { // quadratic pentas
607 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
608 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
609 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
610 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
611 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
612 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
613 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
614 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
615 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
616 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
617 aVal = Max(aVal,Max(Max(L7,L8),L9));
620 else if( len == 20 || len == 27 ) { // quadratic hexas
621 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
622 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
623 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
624 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
625 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
626 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
627 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
628 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
629 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
630 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
631 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
632 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
633 double D1 = getDistance(P( 1 ),P( 7 ));
634 double D2 = getDistance(P( 2 ),P( 8 ));
635 double D3 = getDistance(P( 3 ),P( 5 ));
636 double D4 = getDistance(P( 4 ),P( 6 ));
637 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
638 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
639 aVal = Max(aVal,Max(L11,L12));
640 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
643 else if( len > 1 && aElem->IsPoly() ) { // polys
644 // get the maximum distance between all pairs of nodes
645 for( int i = 1; i <= len; i++ ) {
646 for( int j = 1; j <= len; j++ ) {
647 if( j > i ) { // optimization of the loop
648 double D = getDistance( P(i), P(j) );
649 aVal = Max( aVal, D );
656 if( myPrecision >= 0 )
658 double prec = pow( 10., (double)myPrecision );
659 aVal = floor( aVal * prec + 0.5 ) / prec;
666 double MaxElementLength3D::GetBadRate( double Value, int /*nbNodes*/ ) const
671 SMDSAbs_ElementType MaxElementLength3D::GetType() const
673 return SMDSAbs_Volume;
679 Description : Functor for calculation of minimum angle
682 double MinimumAngle::GetValue( const TSequenceOfXYZ& P )
689 aMin = getAngle(P( P.size() ), P( 1 ), P( 2 ));
690 aMin = Min(aMin,getAngle(P( P.size()-1 ), P( P.size() ), P( 1 )));
692 for (int i=2; i<P.size();i++){
693 double A0 = getAngle( P( i-1 ), P( i ), P( i+1 ) );
697 return aMin * 180.0 / M_PI;
700 double MinimumAngle::GetBadRate( double Value, int nbNodes ) const
702 //const double aBestAngle = PI / nbNodes;
703 const double aBestAngle = 180.0 - ( 360.0 / double(nbNodes) );
704 return ( fabs( aBestAngle - Value ));
707 SMDSAbs_ElementType MinimumAngle::GetType() const
715 Description : Functor for calculating aspect ratio
717 double AspectRatio::GetValue( const TSequenceOfXYZ& P )
719 // According to "Mesh quality control" by Nadir Bouhamau referring to
720 // Pascal Jean Frey and Paul-Louis George. Maillages, applications aux elements finis.
721 // Hermes Science publications, Paris 1999 ISBN 2-7462-0024-4
724 int nbNodes = P.size();
729 // Compute aspect ratio
731 if ( nbNodes == 3 ) {
732 // Compute lengths of the sides
733 std::vector< double > aLen (nbNodes);
734 for ( int i = 0; i < nbNodes - 1; i++ )
735 aLen[ i ] = getDistance( P( i + 1 ), P( i + 2 ) );
736 aLen[ nbNodes - 1 ] = getDistance( P( 1 ), P( nbNodes ) );
737 // Q = alfa * h * p / S, where
739 // alfa = sqrt( 3 ) / 6
740 // h - length of the longest edge
741 // p - half perimeter
742 // S - triangle surface
743 const double alfa = sqrt( 3. ) / 6.;
744 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
745 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
746 double anArea = getArea( P( 1 ), P( 2 ), P( 3 ) );
747 if ( anArea <= Precision::Confusion() )
749 return alfa * maxLen * half_perimeter / anArea;
751 else if ( nbNodes == 6 ) { // quadratic triangles
752 // Compute lengths of the sides
753 std::vector< double > aLen (3);
754 aLen[0] = getDistance( P(1), P(3) );
755 aLen[1] = getDistance( P(3), P(5) );
756 aLen[2] = getDistance( P(5), P(1) );
757 // Q = alfa * h * p / S, where
759 // alfa = sqrt( 3 ) / 6
760 // h - length of the longest edge
761 // p - half perimeter
762 // S - triangle surface
763 const double alfa = sqrt( 3. ) / 6.;
764 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
765 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
766 double anArea = getArea( P(1), P(3), P(5) );
767 if ( anArea <= Precision::Confusion() )
769 return alfa * maxLen * half_perimeter / anArea;
771 else if( nbNodes == 4 ) { // quadrangle
772 // Compute lengths of the sides
773 std::vector< double > aLen (4);
774 aLen[0] = getDistance( P(1), P(2) );
775 aLen[1] = getDistance( P(2), P(3) );
776 aLen[2] = getDistance( P(3), P(4) );
777 aLen[3] = getDistance( P(4), P(1) );
778 // Compute lengths of the diagonals
779 std::vector< double > aDia (2);
780 aDia[0] = getDistance( P(1), P(3) );
781 aDia[1] = getDistance( P(2), P(4) );
782 // Compute areas of all triangles which can be built
783 // taking three nodes of the quadrangle
784 std::vector< double > anArea (4);
785 anArea[0] = getArea( P(1), P(2), P(3) );
786 anArea[1] = getArea( P(1), P(2), P(4) );
787 anArea[2] = getArea( P(1), P(3), P(4) );
788 anArea[3] = getArea( P(2), P(3), P(4) );
789 // Q = alpha * L * C1 / C2, where
791 // alpha = sqrt( 1/32 )
792 // L = max( L1, L2, L3, L4, D1, D2 )
793 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
794 // C2 = min( S1, S2, S3, S4 )
795 // Li - lengths of the edges
796 // Di - lengths of the diagonals
797 // Si - areas of the triangles
798 const double alpha = sqrt( 1 / 32. );
799 double L = Max( aLen[ 0 ],
803 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
804 double C1 = sqrt( ( aLen[0] * aLen[0] +
807 aLen[3] * aLen[3] ) / 4. );
808 double C2 = Min( anArea[ 0 ],
810 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
811 if ( C2 <= Precision::Confusion() )
813 return alpha * L * C1 / C2;
815 else if( nbNodes == 8 || nbNodes == 9 ) { // nbNodes==8 - quadratic quadrangle
816 // Compute lengths of the sides
817 std::vector< double > aLen (4);
818 aLen[0] = getDistance( P(1), P(3) );
819 aLen[1] = getDistance( P(3), P(5) );
820 aLen[2] = getDistance( P(5), P(7) );
821 aLen[3] = getDistance( P(7), P(1) );
822 // Compute lengths of the diagonals
823 std::vector< double > aDia (2);
824 aDia[0] = getDistance( P(1), P(5) );
825 aDia[1] = getDistance( P(3), P(7) );
826 // Compute areas of all triangles which can be built
827 // taking three nodes of the quadrangle
828 std::vector< double > anArea (4);
829 anArea[0] = getArea( P(1), P(3), P(5) );
830 anArea[1] = getArea( P(1), P(3), P(7) );
831 anArea[2] = getArea( P(1), P(5), P(7) );
832 anArea[3] = getArea( P(3), P(5), P(7) );
833 // Q = alpha * L * C1 / C2, where
835 // alpha = sqrt( 1/32 )
836 // L = max( L1, L2, L3, L4, D1, D2 )
837 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
838 // C2 = min( S1, S2, S3, S4 )
839 // Li - lengths of the edges
840 // Di - lengths of the diagonals
841 // Si - areas of the triangles
842 const double alpha = sqrt( 1 / 32. );
843 double L = Max( aLen[ 0 ],
847 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
848 double C1 = sqrt( ( aLen[0] * aLen[0] +
851 aLen[3] * aLen[3] ) / 4. );
852 double C2 = Min( anArea[ 0 ],
854 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
855 if ( C2 <= Precision::Confusion() )
857 return alpha * L * C1 / C2;
862 double AspectRatio::GetBadRate( double Value, int /*nbNodes*/ ) const
864 // the aspect ratio is in the range [1.0,infinity]
865 // < 1.0 = very bad, zero area
868 return ( Value < 0.9 ) ? 1000 : Value / 1000.;
871 SMDSAbs_ElementType AspectRatio::GetType() const
878 Class : AspectRatio3D
879 Description : Functor for calculating aspect ratio
883 inline double getHalfPerimeter(double theTria[3]){
884 return (theTria[0] + theTria[1] + theTria[2])/2.0;
887 inline double getArea(double theHalfPerim, double theTria[3]){
888 return sqrt(theHalfPerim*
889 (theHalfPerim-theTria[0])*
890 (theHalfPerim-theTria[1])*
891 (theHalfPerim-theTria[2]));
894 inline double getVolume(double theLen[6]){
895 double a2 = theLen[0]*theLen[0];
896 double b2 = theLen[1]*theLen[1];
897 double c2 = theLen[2]*theLen[2];
898 double d2 = theLen[3]*theLen[3];
899 double e2 = theLen[4]*theLen[4];
900 double f2 = theLen[5]*theLen[5];
901 double P = 4.0*a2*b2*d2;
902 double Q = a2*(b2+d2-e2)-b2*(a2+d2-f2)-d2*(a2+b2-c2);
903 double R = (b2+d2-e2)*(a2+d2-f2)*(a2+d2-f2);
904 return sqrt(P-Q+R)/12.0;
907 inline double getVolume2(double theLen[6]){
908 double a2 = theLen[0]*theLen[0];
909 double b2 = theLen[1]*theLen[1];
910 double c2 = theLen[2]*theLen[2];
911 double d2 = theLen[3]*theLen[3];
912 double e2 = theLen[4]*theLen[4];
913 double f2 = theLen[5]*theLen[5];
915 double P = a2*e2*(b2+c2+d2+f2-a2-e2);
916 double Q = b2*f2*(a2+c2+d2+e2-b2-f2);
917 double R = c2*d2*(a2+b2+e2+f2-c2-d2);
918 double S = a2*b2*d2+b2*c2*e2+a2*c2*f2+d2*e2*f2;
920 return sqrt(P+Q+R-S)/12.0;
923 inline double getVolume(const TSequenceOfXYZ& P){
924 gp_Vec aVec1( P( 2 ) - P( 1 ) );
925 gp_Vec aVec2( P( 3 ) - P( 1 ) );
926 gp_Vec aVec3( P( 4 ) - P( 1 ) );
927 gp_Vec anAreaVec( aVec1 ^ aVec2 );
928 return fabs(aVec3 * anAreaVec) / 6.0;
931 inline double getMaxHeight(double theLen[6])
933 double aHeight = std::max(theLen[0],theLen[1]);
934 aHeight = std::max(aHeight,theLen[2]);
935 aHeight = std::max(aHeight,theLen[3]);
936 aHeight = std::max(aHeight,theLen[4]);
937 aHeight = std::max(aHeight,theLen[5]);
943 double AspectRatio3D::GetValue( long theId )
946 myCurrElement = myMesh->FindElement( theId );
947 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_TETRA )
949 // Action from CoTech | ACTION 31.3:
950 // EURIWARE BO: Homogenize the formulas used to calculate the Controls in SMESH to fit with
951 // those of ParaView. The library used by ParaView for those calculations can be reused in SMESH.
952 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
953 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
954 aVal = Round( vtkMeshQuality::TetAspectRatio( avtkCell ));
959 if ( GetPoints( myCurrElement, P ))
960 aVal = Round( GetValue( P ));
965 double AspectRatio3D::GetValue( const TSequenceOfXYZ& P )
967 double aQuality = 0.0;
968 if(myCurrElement->IsPoly()) return aQuality;
970 int nbNodes = P.size();
972 if(myCurrElement->IsQuadratic()) {
973 if(nbNodes==10) nbNodes=4; // quadratic tetrahedron
974 else if(nbNodes==13) nbNodes=5; // quadratic pyramid
975 else if(nbNodes==15) nbNodes=6; // quadratic pentahedron
976 else if(nbNodes==20) nbNodes=8; // quadratic hexahedron
977 else if(nbNodes==27) nbNodes=8; // quadratic hexahedron
978 else return aQuality;
984 getDistance(P( 1 ),P( 2 )), // a
985 getDistance(P( 2 ),P( 3 )), // b
986 getDistance(P( 3 ),P( 1 )), // c
987 getDistance(P( 2 ),P( 4 )), // d
988 getDistance(P( 3 ),P( 4 )), // e
989 getDistance(P( 1 ),P( 4 )) // f
991 double aTria[4][3] = {
992 {aLen[0],aLen[1],aLen[2]}, // abc
993 {aLen[0],aLen[3],aLen[5]}, // adf
994 {aLen[1],aLen[3],aLen[4]}, // bde
995 {aLen[2],aLen[4],aLen[5]} // cef
997 double aSumArea = 0.0;
998 double aHalfPerimeter = getHalfPerimeter(aTria[0]);
999 double anArea = getArea(aHalfPerimeter,aTria[0]);
1001 aHalfPerimeter = getHalfPerimeter(aTria[1]);
1002 anArea = getArea(aHalfPerimeter,aTria[1]);
1004 aHalfPerimeter = getHalfPerimeter(aTria[2]);
1005 anArea = getArea(aHalfPerimeter,aTria[2]);
1007 aHalfPerimeter = getHalfPerimeter(aTria[3]);
1008 anArea = getArea(aHalfPerimeter,aTria[3]);
1010 double aVolume = getVolume(P);
1011 //double aVolume = getVolume(aLen);
1012 double aHeight = getMaxHeight(aLen);
1013 static double aCoeff = sqrt(2.0)/12.0;
1014 if ( aVolume > DBL_MIN )
1015 aQuality = aCoeff*aHeight*aSumArea/aVolume;
1020 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 3 ),P( 5 )};
1021 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1024 gp_XYZ aXYZ[4] = {P( 1 ),P( 3 ),P( 4 ),P( 5 )};
1025 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1028 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 5 )};
1029 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1032 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 4 ),P( 5 )};
1033 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1039 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 6 )};
1040 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1043 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 3 )};
1044 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1047 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 6 )};
1048 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1051 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1052 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1055 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 6 )};
1056 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1059 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 3 )};
1060 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1066 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1067 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1070 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 4 )};
1071 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1074 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 7 )};
1075 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1078 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 8 )};
1079 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1082 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 3 )};
1083 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1086 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 4 )};
1087 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1090 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 7 )};
1091 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1094 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 8 )};
1095 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1098 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 3 )};
1099 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1102 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 4 )};
1103 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1106 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 7 )};
1107 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1110 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 8 )};
1111 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1114 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 1 )};
1115 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1118 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 2 )};
1119 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1122 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 5 )};
1123 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1126 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 6 )};
1127 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1130 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 1 )};
1131 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1134 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 2 )};
1135 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1138 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 5 )};
1139 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1142 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 6 )};
1143 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1146 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 1 )};
1147 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1150 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1151 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1154 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 5 )};
1155 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1158 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 6 )};
1159 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1162 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1163 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1166 gp_XYZ aXYZ[4] = {P( 4 ),P( 5 ),P( 8 ),P( 2 )};
1167 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1170 gp_XYZ aXYZ[4] = {P( 1 ),P( 4 ),P( 5 ),P( 3 )};
1171 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1174 gp_XYZ aXYZ[4] = {P( 3 ),P( 6 ),P( 7 ),P( 1 )};
1175 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1178 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 6 ),P( 4 )};
1179 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1182 gp_XYZ aXYZ[4] = {P( 5 ),P( 6 ),P( 8 ),P( 3 )};
1183 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1186 gp_XYZ aXYZ[4] = {P( 7 ),P( 8 ),P( 6 ),P( 1 )};
1187 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1190 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 7 )};
1191 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1194 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 2 ),P( 5 )};
1195 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1201 gp_XYZ aXYZ[8] = {P( 1 ),P( 2 ),P( 4 ),P( 5 ),P( 7 ),P( 8 ),P( 10 ),P( 11 )};
1202 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1205 gp_XYZ aXYZ[8] = {P( 2 ),P( 3 ),P( 5 ),P( 6 ),P( 8 ),P( 9 ),P( 11 ),P( 12 )};
1206 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1209 gp_XYZ aXYZ[8] = {P( 3 ),P( 4 ),P( 6 ),P( 1 ),P( 9 ),P( 10 ),P( 12 ),P( 7 )};
1210 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1213 } // switch(nbNodes)
1215 if ( nbNodes > 4 ) {
1216 // avaluate aspect ratio of quadranle faces
1217 AspectRatio aspect2D;
1218 SMDS_VolumeTool::VolumeType type = SMDS_VolumeTool::GetType( nbNodes );
1219 int nbFaces = SMDS_VolumeTool::NbFaces( type );
1220 TSequenceOfXYZ points(4);
1221 for ( int i = 0; i < nbFaces; ++i ) { // loop on faces of a volume
1222 if ( SMDS_VolumeTool::NbFaceNodes( type, i ) != 4 )
1224 const int* pInd = SMDS_VolumeTool::GetFaceNodesIndices( type, i, true );
1225 for ( int p = 0; p < 4; ++p ) // loop on nodes of a quadranle face
1226 points( p + 1 ) = P( pInd[ p ] + 1 );
1227 aQuality = std::max( aQuality, aspect2D.GetValue( points ));
1233 double AspectRatio3D::GetBadRate( double Value, int /*nbNodes*/ ) const
1235 // the aspect ratio is in the range [1.0,infinity]
1238 return Value / 1000.;
1241 SMDSAbs_ElementType AspectRatio3D::GetType() const
1243 return SMDSAbs_Volume;
1249 Description : Functor for calculating warping
1251 double Warping::GetValue( const TSequenceOfXYZ& P )
1253 if ( P.size() != 4 )
1256 gp_XYZ G = ( P( 1 ) + P( 2 ) + P( 3 ) + P( 4 ) ) / 4.;
1258 double A1 = ComputeA( P( 1 ), P( 2 ), P( 3 ), G );
1259 double A2 = ComputeA( P( 2 ), P( 3 ), P( 4 ), G );
1260 double A3 = ComputeA( P( 3 ), P( 4 ), P( 1 ), G );
1261 double A4 = ComputeA( P( 4 ), P( 1 ), P( 2 ), G );
1263 return Max( Max( A1, A2 ), Max( A3, A4 ) );
1266 double Warping::ComputeA( const gp_XYZ& thePnt1,
1267 const gp_XYZ& thePnt2,
1268 const gp_XYZ& thePnt3,
1269 const gp_XYZ& theG ) const
1271 double aLen1 = gp_Pnt( thePnt1 ).Distance( gp_Pnt( thePnt2 ) );
1272 double aLen2 = gp_Pnt( thePnt2 ).Distance( gp_Pnt( thePnt3 ) );
1273 double L = Min( aLen1, aLen2 ) * 0.5;
1274 if ( L < Precision::Confusion())
1277 gp_XYZ GI = ( thePnt2 + thePnt1 ) / 2. - theG;
1278 gp_XYZ GJ = ( thePnt3 + thePnt2 ) / 2. - theG;
1279 gp_XYZ N = GI.Crossed( GJ );
1281 if ( N.Modulus() < gp::Resolution() )
1286 double H = ( thePnt2 - theG ).Dot( N );
1287 return asin( fabs( H / L ) ) * 180. / M_PI;
1290 double Warping::GetBadRate( double Value, int /*nbNodes*/ ) const
1292 // the warp is in the range [0.0,PI/2]
1293 // 0.0 = good (no warp)
1294 // PI/2 = bad (face pliee)
1298 SMDSAbs_ElementType Warping::GetType() const
1300 return SMDSAbs_Face;
1306 Description : Functor for calculating taper
1308 double Taper::GetValue( const TSequenceOfXYZ& P )
1310 if ( P.size() != 4 )
1314 double J1 = getArea( P( 4 ), P( 1 ), P( 2 ) ) / 2.;
1315 double J2 = getArea( P( 3 ), P( 1 ), P( 2 ) ) / 2.;
1316 double J3 = getArea( P( 2 ), P( 3 ), P( 4 ) ) / 2.;
1317 double J4 = getArea( P( 3 ), P( 4 ), P( 1 ) ) / 2.;
1319 double JA = 0.25 * ( J1 + J2 + J3 + J4 );
1320 if ( JA <= Precision::Confusion() )
1323 double T1 = fabs( ( J1 - JA ) / JA );
1324 double T2 = fabs( ( J2 - JA ) / JA );
1325 double T3 = fabs( ( J3 - JA ) / JA );
1326 double T4 = fabs( ( J4 - JA ) / JA );
1328 return Max( Max( T1, T2 ), Max( T3, T4 ) );
1331 double Taper::GetBadRate( double Value, int /*nbNodes*/ ) const
1333 // the taper is in the range [0.0,1.0]
1334 // 0.0 = good (no taper)
1335 // 1.0 = bad (les cotes opposes sont allignes)
1339 SMDSAbs_ElementType Taper::GetType() const
1341 return SMDSAbs_Face;
1347 Description : Functor for calculating skew in degrees
1349 static inline double skewAngle( const gp_XYZ& p1, const gp_XYZ& p2, const gp_XYZ& p3 )
1351 gp_XYZ p12 = ( p2 + p1 ) / 2.;
1352 gp_XYZ p23 = ( p3 + p2 ) / 2.;
1353 gp_XYZ p31 = ( p3 + p1 ) / 2.;
1355 gp_Vec v1( p31 - p2 ), v2( p12 - p23 );
1357 return v1.Magnitude() < gp::Resolution() || v2.Magnitude() < gp::Resolution() ? 0. : v1.Angle( v2 );
1360 double Skew::GetValue( const TSequenceOfXYZ& P )
1362 if ( P.size() != 3 && P.size() != 4 )
1366 static double PI2 = M_PI / 2.;
1367 if ( P.size() == 3 )
1369 double A0 = fabs( PI2 - skewAngle( P( 3 ), P( 1 ), P( 2 ) ) );
1370 double A1 = fabs( PI2 - skewAngle( P( 1 ), P( 2 ), P( 3 ) ) );
1371 double A2 = fabs( PI2 - skewAngle( P( 2 ), P( 3 ), P( 1 ) ) );
1373 return Max( A0, Max( A1, A2 ) ) * 180. / M_PI;
1377 gp_XYZ p12 = ( P( 1 ) + P( 2 ) ) / 2.;
1378 gp_XYZ p23 = ( P( 2 ) + P( 3 ) ) / 2.;
1379 gp_XYZ p34 = ( P( 3 ) + P( 4 ) ) / 2.;
1380 gp_XYZ p41 = ( P( 4 ) + P( 1 ) ) / 2.;
1382 gp_Vec v1( p34 - p12 ), v2( p23 - p41 );
1383 double A = v1.Magnitude() <= gp::Resolution() || v2.Magnitude() <= gp::Resolution()
1384 ? 0. : fabs( PI2 - v1.Angle( v2 ) );
1387 if ( A < Precision::Angular() )
1390 return A * 180. / M_PI;
1394 double Skew::GetBadRate( double Value, int /*nbNodes*/ ) const
1396 // the skew is in the range [0.0,PI/2].
1402 SMDSAbs_ElementType Skew::GetType() const
1404 return SMDSAbs_Face;
1410 Description : Functor for calculating area
1412 double Area::GetValue( const TSequenceOfXYZ& P )
1415 if ( P.size() > 2 ) {
1416 gp_Vec aVec1( P(2) - P(1) );
1417 gp_Vec aVec2( P(3) - P(1) );
1418 gp_Vec SumVec = aVec1 ^ aVec2;
1419 for (int i=4; i<=P.size(); i++) {
1420 gp_Vec aVec1( P(i-1) - P(1) );
1421 gp_Vec aVec2( P(i) - P(1) );
1422 gp_Vec tmp = aVec1 ^ aVec2;
1425 val = SumVec.Magnitude() * 0.5;
1430 double Area::GetBadRate( double Value, int /*nbNodes*/ ) const
1432 // meaningless as it is not a quality control functor
1436 SMDSAbs_ElementType Area::GetType() const
1438 return SMDSAbs_Face;
1444 Description : Functor for calculating length of edge
1446 double Length::GetValue( const TSequenceOfXYZ& P )
1448 switch ( P.size() ) {
1449 case 2: return getDistance( P( 1 ), P( 2 ) );
1450 case 3: return getDistance( P( 1 ), P( 2 ) ) + getDistance( P( 2 ), P( 3 ) );
1455 double Length::GetBadRate( double Value, int /*nbNodes*/ ) const
1457 // meaningless as it is not quality control functor
1461 SMDSAbs_ElementType Length::GetType() const
1463 return SMDSAbs_Edge;
1468 Description : Functor for calculating length of edge
1471 double Length2D::GetValue( long theElementId)
1475 //cout<<"Length2D::GetValue"<<endl;
1476 if (GetPoints(theElementId,P)){
1477 //for(int jj=1; jj<=P.size(); jj++)
1478 // cout<<"jj="<<jj<<" P("<<P(jj).X()<<","<<P(jj).Y()<<","<<P(jj).Z()<<")"<<endl;
1480 double aVal;// = GetValue( P );
1481 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
1482 SMDSAbs_ElementType aType = aElem->GetType();
1491 aVal = getDistance( P( 1 ), P( 2 ) );
1494 else if (len == 3){ // quadratic edge
1495 aVal = getDistance(P( 1 ),P( 3 )) + getDistance(P( 3 ),P( 2 ));
1499 if (len == 3){ // triangles
1500 double L1 = getDistance(P( 1 ),P( 2 ));
1501 double L2 = getDistance(P( 2 ),P( 3 ));
1502 double L3 = getDistance(P( 3 ),P( 1 ));
1503 aVal = Max(L1,Max(L2,L3));
1506 else if (len == 4){ // quadrangles
1507 double L1 = getDistance(P( 1 ),P( 2 ));
1508 double L2 = getDistance(P( 2 ),P( 3 ));
1509 double L3 = getDistance(P( 3 ),P( 4 ));
1510 double L4 = getDistance(P( 4 ),P( 1 ));
1511 aVal = Max(Max(L1,L2),Max(L3,L4));
1514 if (len == 6){ // quadratic triangles
1515 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1516 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1517 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
1518 aVal = Max(L1,Max(L2,L3));
1519 //cout<<"L1="<<L1<<" L2="<<L2<<"L3="<<L3<<" aVal="<<aVal<<endl;
1522 else if (len == 8){ // quadratic quadrangles
1523 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1524 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1525 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
1526 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
1527 aVal = Max(Max(L1,L2),Max(L3,L4));
1530 case SMDSAbs_Volume:
1531 if (len == 4){ // tetraidrs
1532 double L1 = getDistance(P( 1 ),P( 2 ));
1533 double L2 = getDistance(P( 2 ),P( 3 ));
1534 double L3 = getDistance(P( 3 ),P( 1 ));
1535 double L4 = getDistance(P( 1 ),P( 4 ));
1536 double L5 = getDistance(P( 2 ),P( 4 ));
1537 double L6 = getDistance(P( 3 ),P( 4 ));
1538 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1541 else if (len == 5){ // piramids
1542 double L1 = getDistance(P( 1 ),P( 2 ));
1543 double L2 = getDistance(P( 2 ),P( 3 ));
1544 double L3 = getDistance(P( 3 ),P( 4 ));
1545 double L4 = getDistance(P( 4 ),P( 1 ));
1546 double L5 = getDistance(P( 1 ),P( 5 ));
1547 double L6 = getDistance(P( 2 ),P( 5 ));
1548 double L7 = getDistance(P( 3 ),P( 5 ));
1549 double L8 = getDistance(P( 4 ),P( 5 ));
1551 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1552 aVal = Max(aVal,Max(L7,L8));
1555 else if (len == 6){ // pentaidres
1556 double L1 = getDistance(P( 1 ),P( 2 ));
1557 double L2 = getDistance(P( 2 ),P( 3 ));
1558 double L3 = getDistance(P( 3 ),P( 1 ));
1559 double L4 = getDistance(P( 4 ),P( 5 ));
1560 double L5 = getDistance(P( 5 ),P( 6 ));
1561 double L6 = getDistance(P( 6 ),P( 4 ));
1562 double L7 = getDistance(P( 1 ),P( 4 ));
1563 double L8 = getDistance(P( 2 ),P( 5 ));
1564 double L9 = getDistance(P( 3 ),P( 6 ));
1566 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1567 aVal = Max(aVal,Max(Max(L7,L8),L9));
1570 else if (len == 8){ // hexaider
1571 double L1 = getDistance(P( 1 ),P( 2 ));
1572 double L2 = getDistance(P( 2 ),P( 3 ));
1573 double L3 = getDistance(P( 3 ),P( 4 ));
1574 double L4 = getDistance(P( 4 ),P( 1 ));
1575 double L5 = getDistance(P( 5 ),P( 6 ));
1576 double L6 = getDistance(P( 6 ),P( 7 ));
1577 double L7 = getDistance(P( 7 ),P( 8 ));
1578 double L8 = getDistance(P( 8 ),P( 5 ));
1579 double L9 = getDistance(P( 1 ),P( 5 ));
1580 double L10= getDistance(P( 2 ),P( 6 ));
1581 double L11= getDistance(P( 3 ),P( 7 ));
1582 double L12= getDistance(P( 4 ),P( 8 ));
1584 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1585 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1586 aVal = Max(aVal,Max(L11,L12));
1591 if (len == 10){ // quadratic tetraidrs
1592 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
1593 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
1594 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
1595 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1596 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
1597 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
1598 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1601 else if (len == 13){ // quadratic piramids
1602 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
1603 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
1604 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1605 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1606 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1607 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
1608 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
1609 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
1610 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1611 aVal = Max(aVal,Max(L7,L8));
1614 else if (len == 15){ // quadratic pentaidres
1615 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
1616 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
1617 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1618 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1619 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
1620 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
1621 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
1622 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
1623 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
1624 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1625 aVal = Max(aVal,Max(Max(L7,L8),L9));
1628 else if (len == 20){ // quadratic hexaider
1629 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
1630 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
1631 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
1632 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
1633 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
1634 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
1635 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
1636 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
1637 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
1638 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
1639 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
1640 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
1641 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1642 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1643 aVal = Max(aVal,Max(L11,L12));
1655 if ( myPrecision >= 0 )
1657 double prec = pow( 10., (double)( myPrecision ) );
1658 aVal = floor( aVal * prec + 0.5 ) / prec;
1667 double Length2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1669 // meaningless as it is not quality control functor
1673 SMDSAbs_ElementType Length2D::GetType() const
1675 return SMDSAbs_Face;
1678 Length2D::Value::Value(double theLength,long thePntId1, long thePntId2):
1681 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1682 if(thePntId1 > thePntId2){
1683 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1687 bool Length2D::Value::operator<(const Length2D::Value& x) const{
1688 if(myPntId[0] < x.myPntId[0]) return true;
1689 if(myPntId[0] == x.myPntId[0])
1690 if(myPntId[1] < x.myPntId[1]) return true;
1694 void Length2D::GetValues(TValues& theValues){
1696 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1697 for(; anIter->more(); ){
1698 const SMDS_MeshFace* anElem = anIter->next();
1700 if(anElem->IsQuadratic()) {
1701 const SMDS_VtkFace* F =
1702 dynamic_cast<const SMDS_VtkFace*>(anElem);
1703 // use special nodes iterator
1704 SMDS_ElemIteratorPtr anIter = F->interlacedNodesElemIterator();
1709 const SMDS_MeshElement* aNode;
1711 aNode = anIter->next();
1712 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1713 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1714 aNodeId[0] = aNodeId[1] = aNode->GetID();
1717 for(; anIter->more(); ){
1718 const SMDS_MeshNode* N1 = static_cast<const SMDS_MeshNode*> (anIter->next());
1719 P[2] = gp_Pnt(N1->X(),N1->Y(),N1->Z());
1720 aNodeId[2] = N1->GetID();
1721 aLength = P[1].Distance(P[2]);
1722 if(!anIter->more()) break;
1723 const SMDS_MeshNode* N2 = static_cast<const SMDS_MeshNode*> (anIter->next());
1724 P[3] = gp_Pnt(N2->X(),N2->Y(),N2->Z());
1725 aNodeId[3] = N2->GetID();
1726 aLength += P[2].Distance(P[3]);
1727 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1728 Value aValue2(aLength,aNodeId[2],aNodeId[3]);
1730 aNodeId[1] = aNodeId[3];
1731 theValues.insert(aValue1);
1732 theValues.insert(aValue2);
1734 aLength += P[2].Distance(P[0]);
1735 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1736 Value aValue2(aLength,aNodeId[2],aNodeId[0]);
1737 theValues.insert(aValue1);
1738 theValues.insert(aValue2);
1741 SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
1746 const SMDS_MeshElement* aNode;
1747 if(aNodesIter->more()){
1748 aNode = aNodesIter->next();
1749 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1750 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1751 aNodeId[0] = aNodeId[1] = aNode->GetID();
1754 for(; aNodesIter->more(); ){
1755 aNode = aNodesIter->next();
1756 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1757 long anId = aNode->GetID();
1759 P[2] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1761 aLength = P[1].Distance(P[2]);
1763 Value aValue(aLength,aNodeId[1],anId);
1766 theValues.insert(aValue);
1769 aLength = P[0].Distance(P[1]);
1771 Value aValue(aLength,aNodeId[0],aNodeId[1]);
1772 theValues.insert(aValue);
1778 Class : MultiConnection
1779 Description : Functor for calculating number of faces conneted to the edge
1781 double MultiConnection::GetValue( const TSequenceOfXYZ& P )
1785 double MultiConnection::GetValue( long theId )
1787 return getNbMultiConnection( myMesh, theId );
1790 double MultiConnection::GetBadRate( double Value, int /*nbNodes*/ ) const
1792 // meaningless as it is not quality control functor
1796 SMDSAbs_ElementType MultiConnection::GetType() const
1798 return SMDSAbs_Edge;
1802 Class : MultiConnection2D
1803 Description : Functor for calculating number of faces conneted to the edge
1805 double MultiConnection2D::GetValue( const TSequenceOfXYZ& P )
1810 double MultiConnection2D::GetValue( long theElementId )
1814 const SMDS_MeshElement* aFaceElem = myMesh->FindElement(theElementId);
1815 SMDSAbs_ElementType aType = aFaceElem->GetType();
1820 int i = 0, len = aFaceElem->NbNodes();
1821 SMDS_ElemIteratorPtr anIter = aFaceElem->nodesIterator();
1824 const SMDS_MeshNode *aNode, *aNode0;
1825 TColStd_MapOfInteger aMap, aMapPrev;
1827 for (i = 0; i <= len; i++) {
1832 if (anIter->more()) {
1833 aNode = (SMDS_MeshNode*)anIter->next();
1841 if (i == 0) aNode0 = aNode;
1843 SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
1844 while (anElemIter->more()) {
1845 const SMDS_MeshElement* anElem = anElemIter->next();
1846 if (anElem != 0 && anElem->GetType() == SMDSAbs_Face) {
1847 int anId = anElem->GetID();
1850 if (aMapPrev.Contains(anId)) {
1855 aResult = Max(aResult, aNb);
1866 double MultiConnection2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1868 // meaningless as it is not quality control functor
1872 SMDSAbs_ElementType MultiConnection2D::GetType() const
1874 return SMDSAbs_Face;
1877 MultiConnection2D::Value::Value(long thePntId1, long thePntId2)
1879 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1880 if(thePntId1 > thePntId2){
1881 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1885 bool MultiConnection2D::Value::operator<(const MultiConnection2D::Value& x) const{
1886 if(myPntId[0] < x.myPntId[0]) return true;
1887 if(myPntId[0] == x.myPntId[0])
1888 if(myPntId[1] < x.myPntId[1]) return true;
1892 void MultiConnection2D::GetValues(MValues& theValues){
1893 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1894 for(; anIter->more(); ){
1895 const SMDS_MeshFace* anElem = anIter->next();
1896 SMDS_ElemIteratorPtr aNodesIter;
1897 if ( anElem->IsQuadratic() )
1898 aNodesIter = dynamic_cast<const SMDS_VtkFace*>
1899 (anElem)->interlacedNodesElemIterator();
1901 aNodesIter = anElem->nodesIterator();
1904 //int aNbConnects=0;
1905 const SMDS_MeshNode* aNode0;
1906 const SMDS_MeshNode* aNode1;
1907 const SMDS_MeshNode* aNode2;
1908 if(aNodesIter->more()){
1909 aNode0 = (SMDS_MeshNode*) aNodesIter->next();
1911 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode1;
1912 aNodeId[0] = aNodeId[1] = aNodes->GetID();
1914 for(; aNodesIter->more(); ) {
1915 aNode2 = (SMDS_MeshNode*) aNodesIter->next();
1916 long anId = aNode2->GetID();
1919 Value aValue(aNodeId[1],aNodeId[2]);
1920 MValues::iterator aItr = theValues.find(aValue);
1921 if (aItr != theValues.end()){
1926 theValues[aValue] = 1;
1929 //cout << "NodeIds: "<<aNodeId[1]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1930 aNodeId[1] = aNodeId[2];
1933 Value aValue(aNodeId[0],aNodeId[2]);
1934 MValues::iterator aItr = theValues.find(aValue);
1935 if (aItr != theValues.end()) {
1940 theValues[aValue] = 1;
1943 //cout << "NodeIds: "<<aNodeId[0]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1953 Class : BadOrientedVolume
1954 Description : Predicate bad oriented volumes
1957 BadOrientedVolume::BadOrientedVolume()
1962 void BadOrientedVolume::SetMesh( const SMDS_Mesh* theMesh )
1967 bool BadOrientedVolume::IsSatisfy( long theId )
1972 SMDS_VolumeTool vTool( myMesh->FindElement( theId ));
1973 return !vTool.IsForward();
1976 SMDSAbs_ElementType BadOrientedVolume::GetType() const
1978 return SMDSAbs_Volume;
1982 Class : BareBorderVolume
1985 bool BareBorderVolume::IsSatisfy(long theElementId )
1987 SMDS_VolumeTool myTool;
1988 if ( myTool.Set( myMesh->FindElement(theElementId)))
1990 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
1991 if ( myTool.IsFreeFace( iF ))
1993 const SMDS_MeshNode** n = myTool.GetFaceNodes(iF);
1994 vector< const SMDS_MeshNode*> nodes( n, n+myTool.NbFaceNodes(iF));
1995 if ( !myMesh->FindElement( nodes, SMDSAbs_Face, /*Nomedium=*/false))
2003 Class : BareBorderFace
2006 bool BareBorderFace::IsSatisfy(long theElementId )
2009 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2011 if ( face->GetType() == SMDSAbs_Face )
2013 int nbN = face->NbCornerNodes();
2014 for ( int i = 0; i < nbN && !ok; ++i )
2016 // check if a link is shared by another face
2017 const SMDS_MeshNode* n1 = face->GetNode( i );
2018 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2019 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2020 bool isShared = false;
2021 while ( !isShared && fIt->more() )
2023 const SMDS_MeshElement* f = fIt->next();
2024 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2028 const int iQuad = face->IsQuadratic();
2029 myLinkNodes.resize( 2 + iQuad);
2030 myLinkNodes[0] = n1;
2031 myLinkNodes[1] = n2;
2033 myLinkNodes[2] = face->GetNode( i+nbN );
2034 ok = !myMesh->FindElement( myLinkNodes, SMDSAbs_Edge, /*noMedium=*/false);
2043 Class : OverConstrainedVolume
2046 bool OverConstrainedVolume::IsSatisfy(long theElementId )
2048 // An element is over-constrained if it has N-1 free borders where
2049 // N is the number of edges/faces for a 2D/3D element.
2050 SMDS_VolumeTool myTool;
2051 if ( myTool.Set( myMesh->FindElement(theElementId)))
2053 int nbSharedFaces = 0;
2054 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2055 if ( !myTool.IsFreeFace( iF ) && ++nbSharedFaces > 1 )
2057 return ( nbSharedFaces == 1 );
2063 Class : OverConstrainedFace
2066 bool OverConstrainedFace::IsSatisfy(long theElementId )
2068 // An element is over-constrained if it has N-1 free borders where
2069 // N is the number of edges/faces for a 2D/3D element.
2070 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2071 if ( face->GetType() == SMDSAbs_Face )
2073 int nbSharedBorders = 0;
2074 int nbN = face->NbCornerNodes();
2075 for ( int i = 0; i < nbN; ++i )
2077 // check if a link is shared by another face
2078 const SMDS_MeshNode* n1 = face->GetNode( i );
2079 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2080 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2081 bool isShared = false;
2082 while ( !isShared && fIt->more() )
2084 const SMDS_MeshElement* f = fIt->next();
2085 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2087 if ( isShared && ++nbSharedBorders > 1 )
2090 return ( nbSharedBorders == 1 );
2096 Class : CoincidentNodes
2097 Description : Predicate of Coincident nodes
2100 CoincidentNodes::CoincidentNodes()
2105 bool CoincidentNodes::IsSatisfy( long theElementId )
2107 return myCoincidentIDs.Contains( theElementId );
2110 SMDSAbs_ElementType CoincidentNodes::GetType() const
2112 return SMDSAbs_Node;
2115 void CoincidentNodes::SetMesh( const SMDS_Mesh* theMesh )
2117 myMeshModifTracer.SetMesh( theMesh );
2118 if ( myMeshModifTracer.IsMeshModified() )
2120 TIDSortedNodeSet nodesToCheck;
2121 SMDS_NodeIteratorPtr nIt = theMesh->nodesIterator(/*idInceasingOrder=*/true);
2122 while ( nIt->more() )
2123 nodesToCheck.insert( nodesToCheck.end(), nIt->next() );
2125 list< list< const SMDS_MeshNode*> > nodeGroups;
2126 SMESH_OctreeNode::FindCoincidentNodes ( nodesToCheck, &nodeGroups, myToler );
2128 myCoincidentIDs.Clear();
2129 list< list< const SMDS_MeshNode*> >::iterator groupIt = nodeGroups.begin();
2130 for ( ; groupIt != nodeGroups.end(); ++groupIt )
2132 list< const SMDS_MeshNode*>& coincNodes = *groupIt;
2133 list< const SMDS_MeshNode*>::iterator n = coincNodes.begin();
2134 for ( ; n != coincNodes.end(); ++n )
2135 myCoincidentIDs.Add( (*n)->GetID() );
2141 Class : CoincidentElements
2142 Description : Predicate of Coincident Elements
2143 Note : This class is suitable only for visualization of Coincident Elements
2146 CoincidentElements::CoincidentElements()
2151 void CoincidentElements::SetMesh( const SMDS_Mesh* theMesh )
2156 bool CoincidentElements::IsSatisfy( long theElementId )
2158 if ( !myMesh ) return false;
2160 if ( const SMDS_MeshElement* e = myMesh->FindElement( theElementId ))
2162 if ( e->GetType() != GetType() ) return false;
2163 set< const SMDS_MeshNode* > elemNodes( e->begin_nodes(), e->end_nodes() );
2164 SMDS_ElemIteratorPtr invIt = (*elemNodes.begin())->GetInverseElementIterator( GetType() );
2165 while ( invIt->more() )
2167 const SMDS_MeshElement* e2 = invIt->next();
2168 if ( e2 == e || e2->NbNodes() != (int)elemNodes.size() ) continue;
2170 bool sameNodes = true;
2171 for ( size_t i = 0; i < elemNodes.size() && sameNodes; ++i )
2172 sameNodes = ( elemNodes.count( e2->GetNode( i )));
2180 SMDSAbs_ElementType CoincidentElements1D::GetType() const
2182 return SMDSAbs_Edge;
2184 SMDSAbs_ElementType CoincidentElements2D::GetType() const
2186 return SMDSAbs_Face;
2188 SMDSAbs_ElementType CoincidentElements3D::GetType() const
2190 return SMDSAbs_Volume;
2196 Description : Predicate for free borders
2199 FreeBorders::FreeBorders()
2204 void FreeBorders::SetMesh( const SMDS_Mesh* theMesh )
2209 bool FreeBorders::IsSatisfy( long theId )
2211 return getNbMultiConnection( myMesh, theId ) == 1;
2214 SMDSAbs_ElementType FreeBorders::GetType() const
2216 return SMDSAbs_Edge;
2222 Description : Predicate for free Edges
2224 FreeEdges::FreeEdges()
2229 void FreeEdges::SetMesh( const SMDS_Mesh* theMesh )
2234 bool FreeEdges::IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId )
2236 TColStd_MapOfInteger aMap;
2237 for ( int i = 0; i < 2; i++ )
2239 SMDS_ElemIteratorPtr anElemIter = theNodes[ i ]->GetInverseElementIterator(SMDSAbs_Face);
2240 while( anElemIter->more() )
2242 if ( const SMDS_MeshElement* anElem = anElemIter->next())
2244 const int anId = anElem->GetID();
2245 if ( anId != theFaceId && !aMap.Add( anId ))
2253 bool FreeEdges::IsSatisfy( long theId )
2258 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2259 if ( aFace == 0 || aFace->GetType() != SMDSAbs_Face || aFace->NbNodes() < 3 )
2262 SMDS_ElemIteratorPtr anIter;
2263 if ( aFace->IsQuadratic() ) {
2264 anIter = dynamic_cast<const SMDS_VtkFace*>
2265 (aFace)->interlacedNodesElemIterator();
2268 anIter = aFace->nodesIterator();
2273 int i = 0, nbNodes = aFace->NbNodes();
2274 std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
2275 while( anIter->more() )
2277 const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
2280 aNodes[ i++ ] = aNode;
2282 aNodes[ nbNodes ] = aNodes[ 0 ];
2284 for ( i = 0; i < nbNodes; i++ )
2285 if ( IsFreeEdge( &aNodes[ i ], theId ) )
2291 SMDSAbs_ElementType FreeEdges::GetType() const
2293 return SMDSAbs_Face;
2296 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2299 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
2300 if(thePntId1 > thePntId2){
2301 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
2305 bool FreeEdges::Border::operator<(const FreeEdges::Border& x) const{
2306 if(myPntId[0] < x.myPntId[0]) return true;
2307 if(myPntId[0] == x.myPntId[0])
2308 if(myPntId[1] < x.myPntId[1]) return true;
2312 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2313 FreeEdges::TBorders& theRegistry,
2314 FreeEdges::TBorders& theContainer)
2316 if(theRegistry.find(theBorder) == theRegistry.end()){
2317 theRegistry.insert(theBorder);
2318 theContainer.insert(theBorder);
2320 theContainer.erase(theBorder);
2324 void FreeEdges::GetBoreders(TBorders& theBorders)
2327 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2328 for(; anIter->more(); ){
2329 const SMDS_MeshFace* anElem = anIter->next();
2330 long anElemId = anElem->GetID();
2331 SMDS_ElemIteratorPtr aNodesIter;
2332 if ( anElem->IsQuadratic() )
2333 aNodesIter = static_cast<const SMDS_VtkFace*>(anElem)->
2334 interlacedNodesElemIterator();
2336 aNodesIter = anElem->nodesIterator();
2338 const SMDS_MeshElement* aNode;
2339 if(aNodesIter->more()){
2340 aNode = aNodesIter->next();
2341 aNodeId[0] = aNodeId[1] = aNode->GetID();
2343 for(; aNodesIter->more(); ){
2344 aNode = aNodesIter->next();
2345 long anId = aNode->GetID();
2346 Border aBorder(anElemId,aNodeId[1],anId);
2348 UpdateBorders(aBorder,aRegistry,theBorders);
2350 Border aBorder(anElemId,aNodeId[0],aNodeId[1]);
2351 UpdateBorders(aBorder,aRegistry,theBorders);
2358 Description : Predicate for free nodes
2361 FreeNodes::FreeNodes()
2366 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2371 bool FreeNodes::IsSatisfy( long theNodeId )
2373 const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2377 return (aNode->NbInverseElements() < 1);
2380 SMDSAbs_ElementType FreeNodes::GetType() const
2382 return SMDSAbs_Node;
2388 Description : Predicate for free faces
2391 FreeFaces::FreeFaces()
2396 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2401 bool FreeFaces::IsSatisfy( long theId )
2403 if (!myMesh) return false;
2404 // check that faces nodes refers to less than two common volumes
2405 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2406 if ( !aFace || aFace->GetType() != SMDSAbs_Face )
2409 int nbNode = aFace->NbNodes();
2411 // collect volumes check that number of volumss with count equal nbNode not less than 2
2412 typedef map< SMDS_MeshElement*, int > TMapOfVolume; // map of volume counters
2413 typedef map< SMDS_MeshElement*, int >::iterator TItrMapOfVolume; // iterator
2414 TMapOfVolume mapOfVol;
2416 SMDS_ElemIteratorPtr nodeItr = aFace->nodesIterator();
2417 while ( nodeItr->more() ) {
2418 const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>(nodeItr->next());
2419 if ( !aNode ) continue;
2420 SMDS_ElemIteratorPtr volItr = aNode->GetInverseElementIterator(SMDSAbs_Volume);
2421 while ( volItr->more() ) {
2422 SMDS_MeshElement* aVol = (SMDS_MeshElement*)volItr->next();
2423 TItrMapOfVolume itr = mapOfVol.insert(make_pair(aVol, 0)).first;
2428 TItrMapOfVolume volItr = mapOfVol.begin();
2429 TItrMapOfVolume volEnd = mapOfVol.end();
2430 for ( ; volItr != volEnd; ++volItr )
2431 if ( (*volItr).second >= nbNode )
2433 // face is not free if number of volumes constructed on thier nodes more than one
2437 SMDSAbs_ElementType FreeFaces::GetType() const
2439 return SMDSAbs_Face;
2443 Class : LinearOrQuadratic
2444 Description : Predicate to verify whether a mesh element is linear
2447 LinearOrQuadratic::LinearOrQuadratic()
2452 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2457 bool LinearOrQuadratic::IsSatisfy( long theId )
2459 if (!myMesh) return false;
2460 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2461 if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2463 return (!anElem->IsQuadratic());
2466 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2471 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2478 Description : Functor for check color of group to whic mesh element belongs to
2481 GroupColor::GroupColor()
2485 bool GroupColor::IsSatisfy( long theId )
2487 return (myIDs.find( theId ) != myIDs.end());
2490 void GroupColor::SetType( SMDSAbs_ElementType theType )
2495 SMDSAbs_ElementType GroupColor::GetType() const
2500 static bool isEqual( const Quantity_Color& theColor1,
2501 const Quantity_Color& theColor2 )
2503 // tolerance to compare colors
2504 const double tol = 5*1e-3;
2505 return ( fabs( theColor1.Red() - theColor2.Red() ) < tol &&
2506 fabs( theColor1.Green() - theColor2.Green() ) < tol &&
2507 fabs( theColor1.Blue() - theColor2.Blue() ) < tol );
2511 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2515 const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2519 int nbGrp = aMesh->GetNbGroups();
2523 // iterates on groups and find necessary elements ids
2524 const std::set<SMESHDS_GroupBase*>& aGroups = aMesh->GetGroups();
2525 set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
2526 for (; GrIt != aGroups.end(); GrIt++) {
2527 SMESHDS_GroupBase* aGrp = (*GrIt);
2530 // check type and color of group
2531 if ( !isEqual( myColor, aGrp->GetColor() ) )
2533 if ( myType != SMDSAbs_All && myType != (SMDSAbs_ElementType)aGrp->GetType() )
2536 SMDSAbs_ElementType aGrpElType = (SMDSAbs_ElementType)aGrp->GetType();
2537 if ( myType == aGrpElType || (myType == SMDSAbs_All && aGrpElType != SMDSAbs_Node) ) {
2538 // add elements IDS into control
2539 int aSize = aGrp->Extent();
2540 for (int i = 0; i < aSize; i++)
2541 myIDs.insert( aGrp->GetID(i+1) );
2546 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
2548 TCollection_AsciiString aStr = theStr;
2549 aStr.RemoveAll( ' ' );
2550 aStr.RemoveAll( '\t' );
2551 for ( int aPos = aStr.Search( ";;" ); aPos != -1; aPos = aStr.Search( ";;" ) )
2552 aStr.Remove( aPos, 2 );
2553 Standard_Real clr[3];
2554 clr[0] = clr[1] = clr[2] = 0.;
2555 for ( int i = 0; i < 3; i++ ) {
2556 TCollection_AsciiString tmpStr = aStr.Token( ";", i+1 );
2557 if ( !tmpStr.IsEmpty() && tmpStr.IsRealValue() )
2558 clr[i] = tmpStr.RealValue();
2560 myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
2563 //=======================================================================
2564 // name : GetRangeStr
2565 // Purpose : Get range as a string.
2566 // Example: "1,2,3,50-60,63,67,70-"
2567 //=======================================================================
2568 void GroupColor::GetColorStr( TCollection_AsciiString& theResStr ) const
2571 theResStr += TCollection_AsciiString( myColor.Red() );
2572 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2573 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2577 Class : ElemGeomType
2578 Description : Predicate to check element geometry type
2581 ElemGeomType::ElemGeomType()
2584 myType = SMDSAbs_All;
2585 myGeomType = SMDSGeom_TRIANGLE;
2588 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2593 bool ElemGeomType::IsSatisfy( long theId )
2595 if (!myMesh) return false;
2596 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2599 const SMDSAbs_ElementType anElemType = anElem->GetType();
2600 if ( myType != SMDSAbs_All && anElemType != myType )
2602 const int aNbNode = anElem->NbNodes();
2604 switch( anElemType )
2607 isOk = (myGeomType == SMDSGeom_POINT);
2611 isOk = (myGeomType == SMDSGeom_EDGE);
2615 if ( myGeomType == SMDSGeom_TRIANGLE )
2616 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 6 : aNbNode == 3));
2617 else if ( myGeomType == SMDSGeom_QUADRANGLE )
2618 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? ( aNbNode == 8 || aNbNode == 9 ) : aNbNode == 4));
2619 else if ( myGeomType == SMDSGeom_POLYGON )
2620 isOk = anElem->IsPoly();
2623 case SMDSAbs_Volume:
2624 if ( myGeomType == SMDSGeom_TETRA )
2625 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 10 : aNbNode == 4));
2626 else if ( myGeomType == SMDSGeom_PYRAMID )
2627 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 13 : aNbNode == 5));
2628 else if ( myGeomType == SMDSGeom_PENTA )
2629 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 15 : aNbNode == 6));
2630 else if ( myGeomType == SMDSGeom_HEXA )
2631 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? ( aNbNode == 20 || aNbNode == 27 ): aNbNode == 8));
2632 else if ( myGeomType == SMDSGeom_HEXAGONAL_PRISM )
2633 isOk = (anElem->GetEntityType() == SMDSEntity_Hexagonal_Prism );
2634 else if ( myGeomType == SMDSGeom_POLYHEDRA )
2635 isOk = anElem->IsPoly();
2642 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2647 SMDSAbs_ElementType ElemGeomType::GetType() const
2652 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2654 myGeomType = theType;
2657 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2662 //================================================================================
2664 * \brief Class CoplanarFaces
2666 //================================================================================
2668 CoplanarFaces::CoplanarFaces()
2669 : myFaceID(0), myToler(0)
2672 void CoplanarFaces::SetMesh( const SMDS_Mesh* theMesh )
2674 myMeshModifTracer.SetMesh( theMesh );
2675 if ( myMeshModifTracer.IsMeshModified() )
2677 // Build a set of coplanar face ids
2679 myCoplanarIDs.clear();
2681 if ( !myMeshModifTracer.GetMesh() || !myFaceID || !myToler )
2684 const SMDS_MeshElement* face = myMeshModifTracer.GetMesh()->FindElement( myFaceID );
2685 if ( !face || face->GetType() != SMDSAbs_Face )
2689 gp_Vec myNorm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2693 const double radianTol = myToler * M_PI / 180.;
2694 typedef SMDS_StdIterator< const SMDS_MeshElement*, SMDS_ElemIteratorPtr > TFaceIt;
2695 std::set<const SMDS_MeshElement*> checkedFaces, checkedNodes;
2696 std::list<const SMDS_MeshElement*> faceQueue( 1, face );
2697 while ( !faceQueue.empty() )
2699 face = faceQueue.front();
2700 if ( checkedFaces.insert( face ).second )
2702 gp_Vec norm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2703 if (!normOK || myNorm.Angle( norm ) <= radianTol)
2705 myCoplanarIDs.insert( face->GetID() );
2706 std::set<const SMDS_MeshElement*> neighborFaces;
2707 for ( int i = 0; i < face->NbCornerNodes(); ++i )
2709 const SMDS_MeshNode* n = face->GetNode( i );
2710 if ( checkedNodes.insert( n ).second )
2711 neighborFaces.insert( TFaceIt( n->GetInverseElementIterator(SMDSAbs_Face)),
2714 faceQueue.insert( faceQueue.end(), neighborFaces.begin(), neighborFaces.end() );
2717 faceQueue.pop_front();
2721 bool CoplanarFaces::IsSatisfy( long theElementId )
2723 return myCoplanarIDs.count( theElementId );
2728 *Description : Predicate for Range of Ids.
2729 * Range may be specified with two ways.
2730 * 1. Using AddToRange method
2731 * 2. With SetRangeStr method. Parameter of this method is a string
2732 * like as "1,2,3,50-60,63,67,70-"
2735 //=======================================================================
2736 // name : RangeOfIds
2737 // Purpose : Constructor
2738 //=======================================================================
2739 RangeOfIds::RangeOfIds()
2742 myType = SMDSAbs_All;
2745 //=======================================================================
2747 // Purpose : Set mesh
2748 //=======================================================================
2749 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
2754 //=======================================================================
2755 // name : AddToRange
2756 // Purpose : Add ID to the range
2757 //=======================================================================
2758 bool RangeOfIds::AddToRange( long theEntityId )
2760 myIds.Add( theEntityId );
2764 //=======================================================================
2765 // name : GetRangeStr
2766 // Purpose : Get range as a string.
2767 // Example: "1,2,3,50-60,63,67,70-"
2768 //=======================================================================
2769 void RangeOfIds::GetRangeStr( TCollection_AsciiString& theResStr )
2773 TColStd_SequenceOfInteger anIntSeq;
2774 TColStd_SequenceOfAsciiString aStrSeq;
2776 TColStd_MapIteratorOfMapOfInteger anIter( myIds );
2777 for ( ; anIter.More(); anIter.Next() )
2779 int anId = anIter.Key();
2780 TCollection_AsciiString aStr( anId );
2781 anIntSeq.Append( anId );
2782 aStrSeq.Append( aStr );
2785 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2787 int aMinId = myMin( i );
2788 int aMaxId = myMax( i );
2790 TCollection_AsciiString aStr;
2791 if ( aMinId != IntegerFirst() )
2796 if ( aMaxId != IntegerLast() )
2799 // find position of the string in result sequence and insert string in it
2800 if ( anIntSeq.Length() == 0 )
2802 anIntSeq.Append( aMinId );
2803 aStrSeq.Append( aStr );
2807 if ( aMinId < anIntSeq.First() )
2809 anIntSeq.Prepend( aMinId );
2810 aStrSeq.Prepend( aStr );
2812 else if ( aMinId > anIntSeq.Last() )
2814 anIntSeq.Append( aMinId );
2815 aStrSeq.Append( aStr );
2818 for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
2819 if ( aMinId < anIntSeq( j ) )
2821 anIntSeq.InsertBefore( j, aMinId );
2822 aStrSeq.InsertBefore( j, aStr );
2828 if ( aStrSeq.Length() == 0 )
2831 theResStr = aStrSeq( 1 );
2832 for ( int j = 2, k = aStrSeq.Length(); j <= k; j++ )
2835 theResStr += aStrSeq( j );
2839 //=======================================================================
2840 // name : SetRangeStr
2841 // Purpose : Define range with string
2842 // Example of entry string: "1,2,3,50-60,63,67,70-"
2843 //=======================================================================
2844 bool RangeOfIds::SetRangeStr( const TCollection_AsciiString& theStr )
2850 TCollection_AsciiString aStr = theStr;
2851 aStr.RemoveAll( ' ' );
2852 aStr.RemoveAll( '\t' );
2854 for ( int aPos = aStr.Search( ",," ); aPos != -1; aPos = aStr.Search( ",," ) )
2855 aStr.Remove( aPos, 2 );
2857 TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
2859 while ( tmpStr != "" )
2861 tmpStr = aStr.Token( ",", i++ );
2862 int aPos = tmpStr.Search( '-' );
2866 if ( tmpStr.IsIntegerValue() )
2867 myIds.Add( tmpStr.IntegerValue() );
2873 TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
2874 TCollection_AsciiString aMinStr = tmpStr;
2876 while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
2877 while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
2879 if ( (!aMinStr.IsEmpty() && !aMinStr.IsIntegerValue()) ||
2880 (!aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue()) )
2883 myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
2884 myMax.Append( aMaxStr.IsEmpty() ? IntegerLast() : aMaxStr.IntegerValue() );
2891 //=======================================================================
2893 // Purpose : Get type of supported entities
2894 //=======================================================================
2895 SMDSAbs_ElementType RangeOfIds::GetType() const
2900 //=======================================================================
2902 // Purpose : Set type of supported entities
2903 //=======================================================================
2904 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
2909 //=======================================================================
2911 // Purpose : Verify whether entity satisfies to this rpedicate
2912 //=======================================================================
2913 bool RangeOfIds::IsSatisfy( long theId )
2918 if ( myType == SMDSAbs_Node )
2920 if ( myMesh->FindNode( theId ) == 0 )
2925 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2926 if ( anElem == 0 || (myType != anElem->GetType() && myType != SMDSAbs_All ))
2930 if ( myIds.Contains( theId ) )
2933 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2934 if ( theId >= myMin( i ) && theId <= myMax( i ) )
2942 Description : Base class for comparators
2944 Comparator::Comparator():
2948 Comparator::~Comparator()
2951 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
2954 myFunctor->SetMesh( theMesh );
2957 void Comparator::SetMargin( double theValue )
2959 myMargin = theValue;
2962 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
2964 myFunctor = theFunct;
2967 SMDSAbs_ElementType Comparator::GetType() const
2969 return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
2972 double Comparator::GetMargin()
2980 Description : Comparator "<"
2982 bool LessThan::IsSatisfy( long theId )
2984 return myFunctor && myFunctor->GetValue( theId ) < myMargin;
2990 Description : Comparator ">"
2992 bool MoreThan::IsSatisfy( long theId )
2994 return myFunctor && myFunctor->GetValue( theId ) > myMargin;
3000 Description : Comparator "="
3003 myToler(Precision::Confusion())
3006 bool EqualTo::IsSatisfy( long theId )
3008 return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
3011 void EqualTo::SetTolerance( double theToler )
3016 double EqualTo::GetTolerance()
3023 Description : Logical NOT predicate
3025 LogicalNOT::LogicalNOT()
3028 LogicalNOT::~LogicalNOT()
3031 bool LogicalNOT::IsSatisfy( long theId )
3033 return myPredicate && !myPredicate->IsSatisfy( theId );
3036 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
3039 myPredicate->SetMesh( theMesh );
3042 void LogicalNOT::SetPredicate( PredicatePtr thePred )
3044 myPredicate = thePred;
3047 SMDSAbs_ElementType LogicalNOT::GetType() const
3049 return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
3054 Class : LogicalBinary
3055 Description : Base class for binary logical predicate
3057 LogicalBinary::LogicalBinary()
3060 LogicalBinary::~LogicalBinary()
3063 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
3066 myPredicate1->SetMesh( theMesh );
3069 myPredicate2->SetMesh( theMesh );
3072 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
3074 myPredicate1 = thePredicate;
3077 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
3079 myPredicate2 = thePredicate;
3082 SMDSAbs_ElementType LogicalBinary::GetType() const
3084 if ( !myPredicate1 || !myPredicate2 )
3087 SMDSAbs_ElementType aType1 = myPredicate1->GetType();
3088 SMDSAbs_ElementType aType2 = myPredicate2->GetType();
3090 return aType1 == aType2 ? aType1 : SMDSAbs_All;
3096 Description : Logical AND
3098 bool LogicalAND::IsSatisfy( long theId )
3103 myPredicate1->IsSatisfy( theId ) &&
3104 myPredicate2->IsSatisfy( theId );
3110 Description : Logical OR
3112 bool LogicalOR::IsSatisfy( long theId )
3117 (myPredicate1->IsSatisfy( theId ) ||
3118 myPredicate2->IsSatisfy( theId ));
3132 void Filter::SetPredicate( PredicatePtr thePredicate )
3134 myPredicate = thePredicate;
3137 template<class TElement, class TIterator, class TPredicate>
3138 inline void FillSequence(const TIterator& theIterator,
3139 TPredicate& thePredicate,
3140 Filter::TIdSequence& theSequence)
3142 if ( theIterator ) {
3143 while( theIterator->more() ) {
3144 TElement anElem = theIterator->next();
3145 long anId = anElem->GetID();
3146 if ( thePredicate->IsSatisfy( anId ) )
3147 theSequence.push_back( anId );
3154 GetElementsId( const SMDS_Mesh* theMesh,
3155 PredicatePtr thePredicate,
3156 TIdSequence& theSequence )
3158 theSequence.clear();
3160 if ( !theMesh || !thePredicate )
3163 thePredicate->SetMesh( theMesh );
3165 SMDSAbs_ElementType aType = thePredicate->GetType();
3168 FillSequence<const SMDS_MeshNode*>(theMesh->nodesIterator(),thePredicate,theSequence);
3171 FillSequence<const SMDS_MeshElement*>(theMesh->edgesIterator(),thePredicate,theSequence);
3174 FillSequence<const SMDS_MeshElement*>(theMesh->facesIterator(),thePredicate,theSequence);
3176 case SMDSAbs_Volume:
3177 FillSequence<const SMDS_MeshElement*>(theMesh->volumesIterator(),thePredicate,theSequence);
3180 FillSequence<const SMDS_MeshElement*>(theMesh->edgesIterator(),thePredicate,theSequence);
3181 FillSequence<const SMDS_MeshElement*>(theMesh->facesIterator(),thePredicate,theSequence);
3182 FillSequence<const SMDS_MeshElement*>(theMesh->volumesIterator(),thePredicate,theSequence);
3188 Filter::GetElementsId( const SMDS_Mesh* theMesh,
3189 Filter::TIdSequence& theSequence )
3191 GetElementsId(theMesh,myPredicate,theSequence);
3198 typedef std::set<SMDS_MeshFace*> TMapOfFacePtr;
3204 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
3205 SMDS_MeshNode* theNode2 )
3211 ManifoldPart::Link::~Link()
3217 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
3219 if ( myNode1 == theLink.myNode1 &&
3220 myNode2 == theLink.myNode2 )
3222 else if ( myNode1 == theLink.myNode2 &&
3223 myNode2 == theLink.myNode1 )
3229 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
3231 if(myNode1 < x.myNode1) return true;
3232 if(myNode1 == x.myNode1)
3233 if(myNode2 < x.myNode2) return true;
3237 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
3238 const ManifoldPart::Link& theLink2 )
3240 return theLink1.IsEqual( theLink2 );
3243 ManifoldPart::ManifoldPart()
3246 myAngToler = Precision::Angular();
3247 myIsOnlyManifold = true;
3250 ManifoldPart::~ManifoldPart()
3255 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
3261 SMDSAbs_ElementType ManifoldPart::GetType() const
3262 { return SMDSAbs_Face; }
3264 bool ManifoldPart::IsSatisfy( long theElementId )
3266 return myMapIds.Contains( theElementId );
3269 void ManifoldPart::SetAngleTolerance( const double theAngToler )
3270 { myAngToler = theAngToler; }
3272 double ManifoldPart::GetAngleTolerance() const
3273 { return myAngToler; }
3275 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
3276 { myIsOnlyManifold = theIsOnly; }
3278 void ManifoldPart::SetStartElem( const long theStartId )
3279 { myStartElemId = theStartId; }
3281 bool ManifoldPart::process()
3284 myMapBadGeomIds.Clear();
3286 myAllFacePtr.clear();
3287 myAllFacePtrIntDMap.clear();
3291 // collect all faces into own map
3292 SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3293 for (; anFaceItr->more(); )
3295 SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3296 myAllFacePtr.push_back( aFacePtr );
3297 myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3300 SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3304 // the map of non manifold links and bad geometry
3305 TMapOfLink aMapOfNonManifold;
3306 TColStd_MapOfInteger aMapOfTreated;
3308 // begin cycle on faces from start index and run on vector till the end
3309 // and from begin to start index to cover whole vector
3310 const int aStartIndx = myAllFacePtrIntDMap[aStartFace];
3311 bool isStartTreat = false;
3312 for ( int fi = aStartIndx; !isStartTreat || fi != aStartIndx ; fi++ )
3314 if ( fi == aStartIndx )
3315 isStartTreat = true;
3316 // as result next time when fi will be equal to aStartIndx
3318 SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3319 if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3322 aMapOfTreated.Add( aFacePtr->GetID() );
3323 TColStd_MapOfInteger aResFaces;
3324 if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3325 aMapOfNonManifold, aResFaces ) )
3327 TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3328 for ( ; anItr.More(); anItr.Next() )
3330 int aFaceId = anItr.Key();
3331 aMapOfTreated.Add( aFaceId );
3332 myMapIds.Add( aFaceId );
3335 if ( fi == ( myAllFacePtr.size() - 1 ) )
3337 } // end run on vector of faces
3338 return !myMapIds.IsEmpty();
3341 static void getLinks( const SMDS_MeshFace* theFace,
3342 ManifoldPart::TVectorOfLink& theLinks )
3344 int aNbNode = theFace->NbNodes();
3345 SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3347 SMDS_MeshNode* aNode = 0;
3348 for ( ; aNodeItr->more() && i <= aNbNode; )
3351 SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3355 SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3357 ManifoldPart::Link aLink( aN1, aN2 );
3358 theLinks.push_back( aLink );
3362 bool ManifoldPart::findConnected
3363 ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3364 SMDS_MeshFace* theStartFace,
3365 ManifoldPart::TMapOfLink& theNonManifold,
3366 TColStd_MapOfInteger& theResFaces )
3368 theResFaces.Clear();
3369 if ( !theAllFacePtrInt.size() )
3372 if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3374 myMapBadGeomIds.Add( theStartFace->GetID() );
3378 ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3379 ManifoldPart::TVectorOfLink aSeqOfBoundary;
3380 theResFaces.Add( theStartFace->GetID() );
3381 ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3383 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3384 aDMapLinkFace, theNonManifold, theStartFace );
3386 bool isDone = false;
3387 while ( !isDone && aMapOfBoundary.size() != 0 )
3389 bool isToReset = false;
3390 ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3391 for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3393 ManifoldPart::Link aLink = *pLink;
3394 if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3396 // each link could be treated only once
3397 aMapToSkip.insert( aLink );
3399 ManifoldPart::TVectorOfFacePtr aFaces;
3401 if ( myIsOnlyManifold &&
3402 (theNonManifold.find( aLink ) != theNonManifold.end()) )
3406 getFacesByLink( aLink, aFaces );
3407 // filter the element to keep only indicated elements
3408 ManifoldPart::TVectorOfFacePtr aFiltered;
3409 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3410 for ( ; pFace != aFaces.end(); ++pFace )
3412 SMDS_MeshFace* aFace = *pFace;
3413 if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3414 aFiltered.push_back( aFace );
3417 if ( aFaces.size() < 2 ) // no neihgbour faces
3419 else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3421 theNonManifold.insert( aLink );
3426 // compare normal with normals of neighbor element
3427 SMDS_MeshFace* aPrevFace = aDMapLinkFace[ aLink ];
3428 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3429 for ( ; pFace != aFaces.end(); ++pFace )
3431 SMDS_MeshFace* aNextFace = *pFace;
3432 if ( aPrevFace == aNextFace )
3434 int anNextFaceID = aNextFace->GetID();
3435 if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3436 // should not be with non manifold restriction. probably bad topology
3438 // check if face was treated and skipped
3439 if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3440 !isInPlane( aPrevFace, aNextFace ) )
3442 // add new element to connected and extend the boundaries.
3443 theResFaces.Add( anNextFaceID );
3444 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3445 aDMapLinkFace, theNonManifold, aNextFace );
3449 isDone = !isToReset;
3452 return !theResFaces.IsEmpty();
3455 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3456 const SMDS_MeshFace* theFace2 )
3458 gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3459 gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3460 if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3462 myMapBadGeomIds.Add( theFace2->GetID() );
3465 if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
3471 void ManifoldPart::expandBoundary
3472 ( ManifoldPart::TMapOfLink& theMapOfBoundary,
3473 ManifoldPart::TVectorOfLink& theSeqOfBoundary,
3474 ManifoldPart::TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
3475 ManifoldPart::TMapOfLink& theNonManifold,
3476 SMDS_MeshFace* theNextFace ) const
3478 ManifoldPart::TVectorOfLink aLinks;
3479 getLinks( theNextFace, aLinks );
3480 int aNbLink = (int)aLinks.size();
3481 for ( int i = 0; i < aNbLink; i++ )
3483 ManifoldPart::Link aLink = aLinks[ i ];
3484 if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3486 if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3488 if ( myIsOnlyManifold )
3490 // remove from boundary
3491 theMapOfBoundary.erase( aLink );
3492 ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3493 for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3495 ManifoldPart::Link aBoundLink = *pLink;
3496 if ( aBoundLink.IsEqual( aLink ) )
3498 theSeqOfBoundary.erase( pLink );
3506 theMapOfBoundary.insert( aLink );
3507 theSeqOfBoundary.push_back( aLink );
3508 theDMapLinkFacePtr[ aLink ] = theNextFace;
3513 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3514 ManifoldPart::TVectorOfFacePtr& theFaces ) const
3516 std::set<SMDS_MeshCell *> aSetOfFaces;
3517 // take all faces that shared first node
3518 SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3519 for ( ; anItr->more(); )
3521 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3524 aSetOfFaces.insert( aFace );
3526 // take all faces that shared second node
3527 anItr = theLink.myNode2->facesIterator();
3528 // find the common part of two sets
3529 for ( ; anItr->more(); )
3531 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3532 if ( aSetOfFaces.count( aFace ) )
3533 theFaces.push_back( aFace );
3542 ElementsOnSurface::ElementsOnSurface()
3546 myType = SMDSAbs_All;
3548 myToler = Precision::Confusion();
3549 myUseBoundaries = false;
3552 ElementsOnSurface::~ElementsOnSurface()
3557 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
3559 if ( myMesh == theMesh )
3565 bool ElementsOnSurface::IsSatisfy( long theElementId )
3567 return myIds.Contains( theElementId );
3570 SMDSAbs_ElementType ElementsOnSurface::GetType() const
3573 void ElementsOnSurface::SetTolerance( const double theToler )
3575 if ( myToler != theToler )
3580 double ElementsOnSurface::GetTolerance() const
3583 void ElementsOnSurface::SetUseBoundaries( bool theUse )
3585 if ( myUseBoundaries != theUse ) {
3586 myUseBoundaries = theUse;
3587 SetSurface( mySurf, myType );
3591 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
3592 const SMDSAbs_ElementType theType )
3597 if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
3599 mySurf = TopoDS::Face( theShape );
3600 BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
3602 u1 = SA.FirstUParameter(),
3603 u2 = SA.LastUParameter(),
3604 v1 = SA.FirstVParameter(),
3605 v2 = SA.LastVParameter();
3606 Handle(Geom_Surface) surf = BRep_Tool::Surface( mySurf );
3607 myProjector.Init( surf, u1,u2, v1,v2 );
3611 void ElementsOnSurface::process()
3614 if ( mySurf.IsNull() )
3620 if ( myType == SMDSAbs_Face || myType == SMDSAbs_All )
3622 myIds.ReSize( myMesh->NbFaces() );
3623 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
3624 for(; anIter->more(); )
3625 process( anIter->next() );
3628 if ( myType == SMDSAbs_Edge || myType == SMDSAbs_All )
3630 myIds.ReSize( myIds.Extent() + myMesh->NbEdges() );
3631 SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
3632 for(; anIter->more(); )
3633 process( anIter->next() );
3636 if ( myType == SMDSAbs_Node )
3638 myIds.ReSize( myMesh->NbNodes() );
3639 SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
3640 for(; anIter->more(); )
3641 process( anIter->next() );
3645 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
3647 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3648 bool isSatisfy = true;
3649 for ( ; aNodeItr->more(); )
3651 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3652 if ( !isOnSurface( aNode ) )
3659 myIds.Add( theElemPtr->GetID() );
3662 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
3664 if ( mySurf.IsNull() )
3667 gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
3668 // double aToler2 = myToler * myToler;
3669 // if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
3671 // gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
3672 // if ( aPln.SquareDistance( aPnt ) > aToler2 )
3675 // else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
3677 // gp_Cylinder aCyl = Handle(Geom_CylindricalSurface)::DownCast(mySurf)->Cylinder();
3678 // double aRad = aCyl.Radius();
3679 // gp_Ax3 anAxis = aCyl.Position();
3680 // gp_XYZ aLoc = aCyl.Location().XYZ();
3681 // double aXDist = anAxis.XDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3682 // double aYDist = anAxis.YDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3683 // if ( fabs(aXDist*aXDist + aYDist*aYDist - aRad*aRad) > aToler2 )
3688 myProjector.Perform( aPnt );
3689 bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
3699 ElementsOnShape::ElementsOnShape()
3701 myType(SMDSAbs_All),
3702 myToler(Precision::Confusion()),
3703 myAllNodesFlag(false)
3705 myCurShapeType = TopAbs_SHAPE;
3708 ElementsOnShape::~ElementsOnShape()
3712 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
3714 myMeshModifTracer.SetMesh( theMesh );
3715 if ( myMeshModifTracer.IsMeshModified())
3716 SetShape(myShape, myType);
3719 bool ElementsOnShape::IsSatisfy (long theElementId)
3721 return myIds.Contains(theElementId);
3724 SMDSAbs_ElementType ElementsOnShape::GetType() const
3729 void ElementsOnShape::SetTolerance (const double theToler)
3731 if (myToler != theToler) {
3733 SetShape(myShape, myType);
3737 double ElementsOnShape::GetTolerance() const
3742 void ElementsOnShape::SetAllNodes (bool theAllNodes)
3744 if (myAllNodesFlag != theAllNodes) {
3745 myAllNodesFlag = theAllNodes;
3746 SetShape(myShape, myType);
3750 void ElementsOnShape::SetShape (const TopoDS_Shape& theShape,
3751 const SMDSAbs_ElementType theType)
3757 const SMDS_Mesh* myMesh = myMeshModifTracer.GetMesh();
3759 if ( !myMesh ) return;
3764 myIds.ReSize(myMesh->NbEdges() + myMesh->NbFaces() + myMesh->NbVolumes());
3767 myIds.ReSize(myMesh->NbNodes());
3770 myIds.ReSize(myMesh->NbEdges());
3773 myIds.ReSize(myMesh->NbFaces());
3775 case SMDSAbs_Volume:
3776 myIds.ReSize(myMesh->NbVolumes());
3782 myShapesMap.Clear();
3786 void ElementsOnShape::addShape (const TopoDS_Shape& theShape)
3788 if (theShape.IsNull() || myMeshModifTracer.GetMesh() == 0)
3791 if (!myShapesMap.Add(theShape)) return;
3793 myCurShapeType = theShape.ShapeType();
3794 switch (myCurShapeType)
3796 case TopAbs_COMPOUND:
3797 case TopAbs_COMPSOLID:
3801 TopoDS_Iterator anIt (theShape, Standard_True, Standard_True);
3802 for (; anIt.More(); anIt.Next()) addShape(anIt.Value());
3807 myCurSC.Load(theShape);
3813 TopoDS_Face aFace = TopoDS::Face(theShape);
3814 BRepAdaptor_Surface SA (aFace, true);
3816 u1 = SA.FirstUParameter(),
3817 u2 = SA.LastUParameter(),
3818 v1 = SA.FirstVParameter(),
3819 v2 = SA.LastVParameter();
3820 Handle(Geom_Surface) surf = BRep_Tool::Surface(aFace);
3821 myCurProjFace.Init(surf, u1,u2, v1,v2);
3828 TopoDS_Edge anEdge = TopoDS::Edge(theShape);
3829 Standard_Real u1, u2;
3830 Handle(Geom_Curve) curve = BRep_Tool::Curve(anEdge, u1, u2);
3831 myCurProjEdge.Init(curve, u1, u2);
3837 TopoDS_Vertex aV = TopoDS::Vertex(theShape);
3838 myCurPnt = BRep_Tool::Pnt(aV);
3847 void ElementsOnShape::process()
3849 const SMDS_Mesh* myMesh = myMeshModifTracer.GetMesh();
3850 if (myShape.IsNull() || myMesh == 0)
3853 if (myType == SMDSAbs_Node)
3855 SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
3856 while (anIter->more())
3857 process(anIter->next());
3861 if (myType == SMDSAbs_Edge || myType == SMDSAbs_All)
3863 SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
3864 while (anIter->more())
3865 process(anIter->next());
3868 if (myType == SMDSAbs_Face || myType == SMDSAbs_All)
3870 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
3871 while (anIter->more()) {
3872 process(anIter->next());
3876 if (myType == SMDSAbs_Volume || myType == SMDSAbs_All)
3878 SMDS_VolumeIteratorPtr anIter = myMesh->volumesIterator();
3879 while (anIter->more())
3880 process(anIter->next());
3885 void ElementsOnShape::process (const SMDS_MeshElement* theElemPtr)
3887 if (myShape.IsNull())
3890 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3891 bool isSatisfy = myAllNodesFlag;
3893 gp_XYZ centerXYZ (0, 0, 0);
3895 while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
3897 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3898 gp_Pnt aPnt (aNode->X(), aNode->Y(), aNode->Z());
3899 centerXYZ += aPnt.XYZ();
3901 switch (myCurShapeType)
3905 myCurSC.Perform(aPnt, myToler);
3906 isSatisfy = (myCurSC.State() == TopAbs_IN || myCurSC.State() == TopAbs_ON);
3911 myCurProjFace.Perform(aPnt);
3912 isSatisfy = (myCurProjFace.IsDone() && myCurProjFace.LowerDistance() <= myToler);
3915 // check relatively the face
3916 Quantity_Parameter u, v;
3917 myCurProjFace.LowerDistanceParameters(u, v);
3918 gp_Pnt2d aProjPnt (u, v);
3919 BRepClass_FaceClassifier aClsf (myCurFace, aProjPnt, myToler);
3920 isSatisfy = (aClsf.State() == TopAbs_IN || aClsf.State() == TopAbs_ON);
3926 myCurProjEdge.Perform(aPnt);
3927 isSatisfy = (myCurProjEdge.NbPoints() > 0 && myCurProjEdge.LowerDistance() <= myToler);
3932 isSatisfy = (aPnt.Distance(myCurPnt) <= myToler);
3942 if (isSatisfy && myCurShapeType == TopAbs_SOLID) { // Check the center point for volumes MantisBug 0020168
3943 centerXYZ /= theElemPtr->NbNodes();
3944 gp_Pnt aCenterPnt (centerXYZ);
3945 myCurSC.Perform(aCenterPnt, myToler);
3946 if ( !(myCurSC.State() == TopAbs_IN || myCurSC.State() == TopAbs_ON))
3951 myIds.Add(theElemPtr->GetID());
3954 TSequenceOfXYZ::TSequenceOfXYZ()
3957 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n)
3960 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t)
3963 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray)
3966 template <class InputIterator>
3967 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd)
3970 TSequenceOfXYZ::~TSequenceOfXYZ()
3973 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
3975 myArray = theSequenceOfXYZ.myArray;
3979 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
3981 return myArray[n-1];
3984 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
3986 return myArray[n-1];
3989 void TSequenceOfXYZ::clear()
3994 void TSequenceOfXYZ::reserve(size_type n)
3999 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
4001 myArray.push_back(v);
4004 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
4006 return myArray.size();
4009 TMeshModifTracer::TMeshModifTracer():
4010 myMeshModifTime(0), myMesh(0)
4013 void TMeshModifTracer::SetMesh( const SMDS_Mesh* theMesh )
4015 if ( theMesh != myMesh )
4016 myMeshModifTime = 0;
4019 bool TMeshModifTracer::IsMeshModified()
4021 bool modified = false;
4024 modified = ( myMeshModifTime != myMesh->GetMTime() );
4025 myMeshModifTime = myMesh->GetMTime();