1 // Copyright (C) 2007-2012 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
23 #include "SMESH_ControlsDef.hxx"
25 #include "SMDS_BallElement.hxx"
26 #include "SMDS_Iterator.hxx"
27 #include "SMDS_Mesh.hxx"
28 #include "SMDS_MeshElement.hxx"
29 #include "SMDS_MeshNode.hxx"
30 #include "SMDS_QuadraticEdge.hxx"
31 #include "SMDS_QuadraticFaceOfNodes.hxx"
32 #include "SMDS_VolumeTool.hxx"
33 #include "SMESHDS_GroupBase.hxx"
34 #include "SMESHDS_Mesh.hxx"
35 #include "SMESH_OctreeNode.hxx"
37 #include <BRepAdaptor_Surface.hxx>
38 #include <BRepClass_FaceClassifier.hxx>
39 #include <BRep_Tool.hxx>
40 #include <Geom_CylindricalSurface.hxx>
41 #include <Geom_Plane.hxx>
42 #include <Geom_Surface.hxx>
43 #include <Precision.hxx>
44 #include <TColStd_MapIteratorOfMapOfInteger.hxx>
45 #include <TColStd_MapOfInteger.hxx>
46 #include <TColStd_SequenceOfAsciiString.hxx>
47 #include <TColgp_Array1OfXYZ.hxx>
50 #include <TopoDS_Edge.hxx>
51 #include <TopoDS_Face.hxx>
52 #include <TopoDS_Iterator.hxx>
53 #include <TopoDS_Shape.hxx>
54 #include <TopoDS_Vertex.hxx>
56 #include <gp_Cylinder.hxx>
63 #include <vtkMeshQuality.h>
75 inline gp_XYZ gpXYZ(const SMDS_MeshNode* aNode )
77 return gp_XYZ(aNode->X(), aNode->Y(), aNode->Z() );
80 inline double getAngle( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
82 gp_Vec v1( P1 - P2 ), v2( P3 - P2 );
84 return v1.Magnitude() < gp::Resolution() ||
85 v2.Magnitude() < gp::Resolution() ? 0 : v1.Angle( v2 );
88 inline double getArea( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
90 gp_Vec aVec1( P2 - P1 );
91 gp_Vec aVec2( P3 - P1 );
92 return ( aVec1 ^ aVec2 ).Magnitude() * 0.5;
95 inline double getArea( const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3 )
97 return getArea( P1.XYZ(), P2.XYZ(), P3.XYZ() );
102 inline double getDistance( const gp_XYZ& P1, const gp_XYZ& P2 )
104 double aDist = gp_Pnt( P1 ).Distance( gp_Pnt( P2 ) );
108 int getNbMultiConnection( const SMDS_Mesh* theMesh, const int theId )
113 const SMDS_MeshElement* anEdge = theMesh->FindElement( theId );
114 if ( anEdge == 0 || anEdge->GetType() != SMDSAbs_Edge/* || anEdge->NbNodes() != 2 */)
117 // for each pair of nodes in anEdge (there are 2 pairs in a quadratic edge)
118 // count elements containing both nodes of the pair.
119 // Note that there may be such cases for a quadratic edge (a horizontal line):
124 // +-----+------+ +-----+------+
127 // result sould be 2 in both cases
129 int aResult0 = 0, aResult1 = 0;
130 // last node, it is a medium one in a quadratic edge
131 const SMDS_MeshNode* aLastNode = anEdge->GetNode( anEdge->NbNodes() - 1 );
132 const SMDS_MeshNode* aNode0 = anEdge->GetNode( 0 );
133 const SMDS_MeshNode* aNode1 = anEdge->GetNode( 1 );
134 if ( aNode1 == aLastNode ) aNode1 = 0;
136 SMDS_ElemIteratorPtr anElemIter = aLastNode->GetInverseElementIterator();
137 while( anElemIter->more() ) {
138 const SMDS_MeshElement* anElem = anElemIter->next();
139 if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
140 SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
141 while ( anIter->more() ) {
142 if ( const SMDS_MeshElement* anElemNode = anIter->next() ) {
143 if ( anElemNode == aNode0 ) {
145 if ( !aNode1 ) break; // not a quadratic edge
147 else if ( anElemNode == aNode1 )
153 int aResult = std::max ( aResult0, aResult1 );
155 // TColStd_MapOfInteger aMap;
157 // SMDS_ElemIteratorPtr anIter = anEdge->nodesIterator();
158 // if ( anIter != 0 ) {
159 // while( anIter->more() ) {
160 // const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
163 // SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
164 // while( anElemIter->more() ) {
165 // const SMDS_MeshElement* anElem = anElemIter->next();
166 // if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
167 // int anId = anElem->GetID();
169 // if ( anIter->more() ) // i.e. first node
171 // else if ( aMap.Contains( anId ) )
181 gp_XYZ getNormale( const SMDS_MeshFace* theFace, bool* ok=0 )
183 int aNbNode = theFace->NbNodes();
185 gp_XYZ q1 = gpXYZ( theFace->GetNode(1)) - gpXYZ( theFace->GetNode(0));
186 gp_XYZ q2 = gpXYZ( theFace->GetNode(2)) - gpXYZ( theFace->GetNode(0));
189 gp_XYZ q3 = gpXYZ( theFace->GetNode(3)) - gpXYZ( theFace->GetNode(0));
192 double len = n.Modulus();
193 bool zeroLen = ( len <= numeric_limits<double>::min());
197 if (ok) *ok = !zeroLen;
205 using namespace SMESH::Controls;
212 Class : NumericalFunctor
213 Description : Base class for numerical functors
215 NumericalFunctor::NumericalFunctor():
221 void NumericalFunctor::SetMesh( const SMDS_Mesh* theMesh )
226 bool NumericalFunctor::GetPoints(const int theId,
227 TSequenceOfXYZ& theRes ) const
234 return GetPoints( myMesh->FindElement( theId ), theRes );
237 bool NumericalFunctor::GetPoints(const SMDS_MeshElement* anElem,
238 TSequenceOfXYZ& theRes )
245 theRes.reserve( anElem->NbNodes() );
247 // Get nodes of the element
248 SMDS_ElemIteratorPtr anIter;
250 if ( anElem->IsQuadratic() ) {
251 switch ( anElem->GetType() ) {
253 anIter = dynamic_cast<const SMDS_VtkEdge*>
254 (anElem)->interlacedNodesElemIterator();
257 anIter = dynamic_cast<const SMDS_VtkFace*>
258 (anElem)->interlacedNodesElemIterator();
261 anIter = anElem->nodesIterator();
266 anIter = anElem->nodesIterator();
270 while( anIter->more() ) {
271 if ( const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>( anIter->next() ))
272 theRes.push_back( gp_XYZ( aNode->X(), aNode->Y(), aNode->Z() ) );
279 long NumericalFunctor::GetPrecision() const
284 void NumericalFunctor::SetPrecision( const long thePrecision )
286 myPrecision = thePrecision;
287 myPrecisionValue = pow( 10., (double)( myPrecision ) );
290 double NumericalFunctor::GetValue( long theId )
294 myCurrElement = myMesh->FindElement( theId );
297 if ( GetPoints( theId, P ))
298 aVal = Round( GetValue( P ));
303 double NumericalFunctor::Round( const double & aVal )
305 return ( myPrecision >= 0 ) ? floor( aVal * myPrecisionValue + 0.5 ) / myPrecisionValue : aVal;
308 //================================================================================
310 * \brief Return histogram of functor values
311 * \param nbIntervals - number of intervals
312 * \param nbEvents - number of mesh elements having values within i-th interval
313 * \param funValues - boundaries of intervals
314 * \param elements - elements to check vulue of; empty list means "of all"
315 * \param minmax - boundaries of diapason of values to divide into intervals
317 //================================================================================
319 void NumericalFunctor::GetHistogram(int nbIntervals,
320 std::vector<int>& nbEvents,
321 std::vector<double>& funValues,
322 const vector<int>& elements,
323 const double* minmax)
325 if ( nbIntervals < 1 ||
327 !myMesh->GetMeshInfo().NbElements( GetType() ))
329 nbEvents.resize( nbIntervals, 0 );
330 funValues.resize( nbIntervals+1 );
332 // get all values sorted
333 std::multiset< double > values;
334 if ( elements.empty() )
336 SMDS_ElemIteratorPtr elemIt = myMesh->elementsIterator(GetType());
337 while ( elemIt->more() )
338 values.insert( GetValue( elemIt->next()->GetID() ));
342 vector<int>::const_iterator id = elements.begin();
343 for ( ; id != elements.end(); ++id )
344 values.insert( GetValue( *id ));
349 funValues[0] = minmax[0];
350 funValues[nbIntervals] = minmax[1];
354 funValues[0] = *values.begin();
355 funValues[nbIntervals] = *values.rbegin();
357 // case nbIntervals == 1
358 if ( nbIntervals == 1 )
360 nbEvents[0] = values.size();
364 if (funValues.front() == funValues.back())
366 nbEvents.resize( 1 );
367 nbEvents[0] = values.size();
368 funValues[1] = funValues.back();
369 funValues.resize( 2 );
372 std::multiset< double >::iterator min = values.begin(), max;
373 for ( int i = 0; i < nbIntervals; ++i )
375 // find end value of i-th interval
376 double r = (i+1) / double( nbIntervals );
377 funValues[i+1] = funValues.front() * (1-r) + funValues.back() * r;
379 // count values in the i-th interval if there are any
380 if ( min != values.end() && *min <= funValues[i+1] )
382 // find the first value out of the interval
383 max = values.upper_bound( funValues[i+1] ); // max is greater than funValues[i+1], or end()
384 nbEvents[i] = std::distance( min, max );
388 // add values larger than minmax[1]
389 nbEvents.back() += std::distance( min, values.end() );
392 //=======================================================================
393 //function : GetValue
395 //=======================================================================
397 double Volume::GetValue( long theElementId )
399 if ( theElementId && myMesh ) {
400 SMDS_VolumeTool aVolumeTool;
401 if ( aVolumeTool.Set( myMesh->FindElement( theElementId )))
402 return aVolumeTool.GetSize();
407 //=======================================================================
408 //function : GetBadRate
409 //purpose : meaningless as it is not quality control functor
410 //=======================================================================
412 double Volume::GetBadRate( double Value, int /*nbNodes*/ ) const
417 //=======================================================================
420 //=======================================================================
422 SMDSAbs_ElementType Volume::GetType() const
424 return SMDSAbs_Volume;
429 Class : MaxElementLength2D
430 Description : Functor calculating maximum length of 2D element
432 double MaxElementLength2D::GetValue( const TSequenceOfXYZ& P )
438 if( len == 3 ) { // triangles
439 double L1 = getDistance(P( 1 ),P( 2 ));
440 double L2 = getDistance(P( 2 ),P( 3 ));
441 double L3 = getDistance(P( 3 ),P( 1 ));
442 aVal = Max(L1,Max(L2,L3));
444 else if( len == 4 ) { // quadrangles
445 double L1 = getDistance(P( 1 ),P( 2 ));
446 double L2 = getDistance(P( 2 ),P( 3 ));
447 double L3 = getDistance(P( 3 ),P( 4 ));
448 double L4 = getDistance(P( 4 ),P( 1 ));
449 double D1 = getDistance(P( 1 ),P( 3 ));
450 double D2 = getDistance(P( 2 ),P( 4 ));
451 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
453 else if( len == 6 ) { // quadratic triangles
454 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
455 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
456 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
457 aVal = Max(L1,Max(L2,L3));
459 else if( len == 8 || len == 9 ) { // quadratic quadrangles
460 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
461 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
462 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
463 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
464 double D1 = getDistance(P( 1 ),P( 5 ));
465 double D2 = getDistance(P( 3 ),P( 7 ));
466 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
469 if( myPrecision >= 0 )
471 double prec = pow( 10., (double)myPrecision );
472 aVal = floor( aVal * prec + 0.5 ) / prec;
477 double MaxElementLength2D::GetValue( long theElementId )
480 if( GetPoints( theElementId, P ) && myMesh->FindElement( theElementId )->GetType() == SMDSAbs_Face) {
486 double MaxElementLength2D::GetBadRate( double Value, int /*nbNodes*/ ) const
491 SMDSAbs_ElementType MaxElementLength2D::GetType() const
497 Class : MaxElementLength3D
498 Description : Functor calculating maximum length of 3D element
501 double MaxElementLength3D::GetValue( long theElementId )
504 if( GetPoints( theElementId, P ) ) {
506 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
507 SMDSAbs_ElementType aType = aElem->GetType();
511 if( len == 4 ) { // tetras
512 double L1 = getDistance(P( 1 ),P( 2 ));
513 double L2 = getDistance(P( 2 ),P( 3 ));
514 double L3 = getDistance(P( 3 ),P( 1 ));
515 double L4 = getDistance(P( 1 ),P( 4 ));
516 double L5 = getDistance(P( 2 ),P( 4 ));
517 double L6 = getDistance(P( 3 ),P( 4 ));
518 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
521 else if( len == 5 ) { // pyramids
522 double L1 = getDistance(P( 1 ),P( 2 ));
523 double L2 = getDistance(P( 2 ),P( 3 ));
524 double L3 = getDistance(P( 3 ),P( 4 ));
525 double L4 = getDistance(P( 4 ),P( 1 ));
526 double L5 = getDistance(P( 1 ),P( 5 ));
527 double L6 = getDistance(P( 2 ),P( 5 ));
528 double L7 = getDistance(P( 3 ),P( 5 ));
529 double L8 = getDistance(P( 4 ),P( 5 ));
530 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
531 aVal = Max(aVal,Max(L7,L8));
534 else if( len == 6 ) { // pentas
535 double L1 = getDistance(P( 1 ),P( 2 ));
536 double L2 = getDistance(P( 2 ),P( 3 ));
537 double L3 = getDistance(P( 3 ),P( 1 ));
538 double L4 = getDistance(P( 4 ),P( 5 ));
539 double L5 = getDistance(P( 5 ),P( 6 ));
540 double L6 = getDistance(P( 6 ),P( 4 ));
541 double L7 = getDistance(P( 1 ),P( 4 ));
542 double L8 = getDistance(P( 2 ),P( 5 ));
543 double L9 = getDistance(P( 3 ),P( 6 ));
544 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
545 aVal = Max(aVal,Max(Max(L7,L8),L9));
548 else if( len == 8 ) { // hexas
549 double L1 = getDistance(P( 1 ),P( 2 ));
550 double L2 = getDistance(P( 2 ),P( 3 ));
551 double L3 = getDistance(P( 3 ),P( 4 ));
552 double L4 = getDistance(P( 4 ),P( 1 ));
553 double L5 = getDistance(P( 5 ),P( 6 ));
554 double L6 = getDistance(P( 6 ),P( 7 ));
555 double L7 = getDistance(P( 7 ),P( 8 ));
556 double L8 = getDistance(P( 8 ),P( 5 ));
557 double L9 = getDistance(P( 1 ),P( 5 ));
558 double L10= getDistance(P( 2 ),P( 6 ));
559 double L11= getDistance(P( 3 ),P( 7 ));
560 double L12= getDistance(P( 4 ),P( 8 ));
561 double D1 = getDistance(P( 1 ),P( 7 ));
562 double D2 = getDistance(P( 2 ),P( 8 ));
563 double D3 = getDistance(P( 3 ),P( 5 ));
564 double D4 = getDistance(P( 4 ),P( 6 ));
565 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
566 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
567 aVal = Max(aVal,Max(L11,L12));
568 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
571 else if( len == 12 ) { // hexagonal prism
572 for ( int i1 = 1; i1 < 12; ++i1 )
573 for ( int i2 = i1+1; i1 <= 12; ++i1 )
574 aVal = Max( aVal, getDistance(P( i1 ),P( i2 )));
577 else if( len == 10 ) { // quadratic tetras
578 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
579 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
580 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
581 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
582 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
583 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
584 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
587 else if( len == 13 ) { // quadratic pyramids
588 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
589 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
590 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
591 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
592 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
593 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
594 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
595 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
596 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
597 aVal = Max(aVal,Max(L7,L8));
600 else if( len == 15 ) { // quadratic pentas
601 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
602 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
603 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
604 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
605 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
606 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
607 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
608 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
609 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
610 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
611 aVal = Max(aVal,Max(Max(L7,L8),L9));
614 else if( len == 20 || len == 27 ) { // quadratic hexas
615 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
616 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
617 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
618 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
619 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
620 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
621 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
622 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
623 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
624 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
625 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
626 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
627 double D1 = getDistance(P( 1 ),P( 7 ));
628 double D2 = getDistance(P( 2 ),P( 8 ));
629 double D3 = getDistance(P( 3 ),P( 5 ));
630 double D4 = getDistance(P( 4 ),P( 6 ));
631 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
632 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
633 aVal = Max(aVal,Max(L11,L12));
634 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
637 else if( len > 1 && aElem->IsPoly() ) { // polys
638 // get the maximum distance between all pairs of nodes
639 for( int i = 1; i <= len; i++ ) {
640 for( int j = 1; j <= len; j++ ) {
641 if( j > i ) { // optimization of the loop
642 double D = getDistance( P(i), P(j) );
643 aVal = Max( aVal, D );
650 if( myPrecision >= 0 )
652 double prec = pow( 10., (double)myPrecision );
653 aVal = floor( aVal * prec + 0.5 ) / prec;
660 double MaxElementLength3D::GetBadRate( double Value, int /*nbNodes*/ ) const
665 SMDSAbs_ElementType MaxElementLength3D::GetType() const
667 return SMDSAbs_Volume;
673 Description : Functor for calculation of minimum angle
676 double MinimumAngle::GetValue( const TSequenceOfXYZ& P )
683 aMin = getAngle(P( P.size() ), P( 1 ), P( 2 ));
684 aMin = Min(aMin,getAngle(P( P.size()-1 ), P( P.size() ), P( 1 )));
686 for (int i=2; i<P.size();i++){
687 double A0 = getAngle( P( i-1 ), P( i ), P( i+1 ) );
691 return aMin * 180.0 / M_PI;
694 double MinimumAngle::GetBadRate( double Value, int nbNodes ) const
696 //const double aBestAngle = PI / nbNodes;
697 const double aBestAngle = 180.0 - ( 360.0 / double(nbNodes) );
698 return ( fabs( aBestAngle - Value ));
701 SMDSAbs_ElementType MinimumAngle::GetType() const
709 Description : Functor for calculating aspect ratio
711 double AspectRatio::GetValue( long theId )
714 myCurrElement = myMesh->FindElement( theId );
715 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_QUAD )
718 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
719 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
720 aVal = Round( vtkMeshQuality::QuadAspectRatio( avtkCell ));
725 if ( GetPoints( myCurrElement, P ))
726 aVal = Round( GetValue( P ));
731 double AspectRatio::GetValue( const TSequenceOfXYZ& P )
733 // According to "Mesh quality control" by Nadir Bouhamau referring to
734 // Pascal Jean Frey and Paul-Louis George. Maillages, applications aux elements finis.
735 // Hermes Science publications, Paris 1999 ISBN 2-7462-0024-4
738 int nbNodes = P.size();
743 // Compute aspect ratio
745 if ( nbNodes == 3 ) {
746 // Compute lengths of the sides
747 std::vector< double > aLen (nbNodes);
748 for ( int i = 0; i < nbNodes - 1; i++ )
749 aLen[ i ] = getDistance( P( i + 1 ), P( i + 2 ) );
750 aLen[ nbNodes - 1 ] = getDistance( P( 1 ), P( nbNodes ) );
751 // Q = alfa * h * p / S, where
753 // alfa = sqrt( 3 ) / 6
754 // h - length of the longest edge
755 // p - half perimeter
756 // S - triangle surface
757 const double alfa = sqrt( 3. ) / 6.;
758 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
759 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
760 double anArea = getArea( P( 1 ), P( 2 ), P( 3 ) );
761 if ( anArea <= Precision::Confusion() )
763 return alfa * maxLen * half_perimeter / anArea;
765 else if ( nbNodes == 6 ) { // quadratic triangles
766 // Compute lengths of the sides
767 std::vector< double > aLen (3);
768 aLen[0] = getDistance( P(1), P(3) );
769 aLen[1] = getDistance( P(3), P(5) );
770 aLen[2] = getDistance( P(5), P(1) );
771 // Q = alfa * h * p / S, where
773 // alfa = sqrt( 3 ) / 6
774 // h - length of the longest edge
775 // p - half perimeter
776 // S - triangle surface
777 const double alfa = sqrt( 3. ) / 6.;
778 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
779 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
780 double anArea = getArea( P(1), P(3), P(5) );
781 if ( anArea <= Precision::Confusion() )
783 return alfa * maxLen * half_perimeter / anArea;
785 else if( nbNodes == 4 ) { // quadrangle
786 // Compute lengths of the sides
787 std::vector< double > aLen (4);
788 aLen[0] = getDistance( P(1), P(2) );
789 aLen[1] = getDistance( P(2), P(3) );
790 aLen[2] = getDistance( P(3), P(4) );
791 aLen[3] = getDistance( P(4), P(1) );
792 // Compute lengths of the diagonals
793 std::vector< double > aDia (2);
794 aDia[0] = getDistance( P(1), P(3) );
795 aDia[1] = getDistance( P(2), P(4) );
796 // Compute areas of all triangles which can be built
797 // taking three nodes of the quadrangle
798 std::vector< double > anArea (4);
799 anArea[0] = getArea( P(1), P(2), P(3) );
800 anArea[1] = getArea( P(1), P(2), P(4) );
801 anArea[2] = getArea( P(1), P(3), P(4) );
802 anArea[3] = getArea( P(2), P(3), P(4) );
803 // Q = alpha * L * C1 / C2, where
805 // alpha = sqrt( 1/32 )
806 // L = max( L1, L2, L3, L4, D1, D2 )
807 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
808 // C2 = min( S1, S2, S3, S4 )
809 // Li - lengths of the edges
810 // Di - lengths of the diagonals
811 // Si - areas of the triangles
812 const double alpha = sqrt( 1 / 32. );
813 double L = Max( aLen[ 0 ],
817 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
818 double C1 = sqrt( ( aLen[0] * aLen[0] +
821 aLen[3] * aLen[3] ) / 4. );
822 double C2 = Min( anArea[ 0 ],
824 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
825 if ( C2 <= Precision::Confusion() )
827 return alpha * L * C1 / C2;
829 else if( nbNodes == 8 || nbNodes == 9 ) { // nbNodes==8 - quadratic quadrangle
830 // Compute lengths of the sides
831 std::vector< double > aLen (4);
832 aLen[0] = getDistance( P(1), P(3) );
833 aLen[1] = getDistance( P(3), P(5) );
834 aLen[2] = getDistance( P(5), P(7) );
835 aLen[3] = getDistance( P(7), P(1) );
836 // Compute lengths of the diagonals
837 std::vector< double > aDia (2);
838 aDia[0] = getDistance( P(1), P(5) );
839 aDia[1] = getDistance( P(3), P(7) );
840 // Compute areas of all triangles which can be built
841 // taking three nodes of the quadrangle
842 std::vector< double > anArea (4);
843 anArea[0] = getArea( P(1), P(3), P(5) );
844 anArea[1] = getArea( P(1), P(3), P(7) );
845 anArea[2] = getArea( P(1), P(5), P(7) );
846 anArea[3] = getArea( P(3), P(5), P(7) );
847 // Q = alpha * L * C1 / C2, where
849 // alpha = sqrt( 1/32 )
850 // L = max( L1, L2, L3, L4, D1, D2 )
851 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
852 // C2 = min( S1, S2, S3, S4 )
853 // Li - lengths of the edges
854 // Di - lengths of the diagonals
855 // Si - areas of the triangles
856 const double alpha = sqrt( 1 / 32. );
857 double L = Max( aLen[ 0 ],
861 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
862 double C1 = sqrt( ( aLen[0] * aLen[0] +
865 aLen[3] * aLen[3] ) / 4. );
866 double C2 = Min( anArea[ 0 ],
868 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
869 if ( C2 <= Precision::Confusion() )
871 return alpha * L * C1 / C2;
876 double AspectRatio::GetBadRate( double Value, int /*nbNodes*/ ) const
878 // the aspect ratio is in the range [1.0,infinity]
879 // < 1.0 = very bad, zero area
882 return ( Value < 0.9 ) ? 1000 : Value / 1000.;
885 SMDSAbs_ElementType AspectRatio::GetType() const
892 Class : AspectRatio3D
893 Description : Functor for calculating aspect ratio
897 inline double getHalfPerimeter(double theTria[3]){
898 return (theTria[0] + theTria[1] + theTria[2])/2.0;
901 inline double getArea(double theHalfPerim, double theTria[3]){
902 return sqrt(theHalfPerim*
903 (theHalfPerim-theTria[0])*
904 (theHalfPerim-theTria[1])*
905 (theHalfPerim-theTria[2]));
908 inline double getVolume(double theLen[6]){
909 double a2 = theLen[0]*theLen[0];
910 double b2 = theLen[1]*theLen[1];
911 double c2 = theLen[2]*theLen[2];
912 double d2 = theLen[3]*theLen[3];
913 double e2 = theLen[4]*theLen[4];
914 double f2 = theLen[5]*theLen[5];
915 double P = 4.0*a2*b2*d2;
916 double Q = a2*(b2+d2-e2)-b2*(a2+d2-f2)-d2*(a2+b2-c2);
917 double R = (b2+d2-e2)*(a2+d2-f2)*(a2+d2-f2);
918 return sqrt(P-Q+R)/12.0;
921 inline double getVolume2(double theLen[6]){
922 double a2 = theLen[0]*theLen[0];
923 double b2 = theLen[1]*theLen[1];
924 double c2 = theLen[2]*theLen[2];
925 double d2 = theLen[3]*theLen[3];
926 double e2 = theLen[4]*theLen[4];
927 double f2 = theLen[5]*theLen[5];
929 double P = a2*e2*(b2+c2+d2+f2-a2-e2);
930 double Q = b2*f2*(a2+c2+d2+e2-b2-f2);
931 double R = c2*d2*(a2+b2+e2+f2-c2-d2);
932 double S = a2*b2*d2+b2*c2*e2+a2*c2*f2+d2*e2*f2;
934 return sqrt(P+Q+R-S)/12.0;
937 inline double getVolume(const TSequenceOfXYZ& P){
938 gp_Vec aVec1( P( 2 ) - P( 1 ) );
939 gp_Vec aVec2( P( 3 ) - P( 1 ) );
940 gp_Vec aVec3( P( 4 ) - P( 1 ) );
941 gp_Vec anAreaVec( aVec1 ^ aVec2 );
942 return fabs(aVec3 * anAreaVec) / 6.0;
945 inline double getMaxHeight(double theLen[6])
947 double aHeight = std::max(theLen[0],theLen[1]);
948 aHeight = std::max(aHeight,theLen[2]);
949 aHeight = std::max(aHeight,theLen[3]);
950 aHeight = std::max(aHeight,theLen[4]);
951 aHeight = std::max(aHeight,theLen[5]);
957 double AspectRatio3D::GetValue( long theId )
960 myCurrElement = myMesh->FindElement( theId );
961 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_TETRA )
963 // Action from CoTech | ACTION 31.3:
964 // EURIWARE BO: Homogenize the formulas used to calculate the Controls in SMESH to fit with
965 // those of ParaView. The library used by ParaView for those calculations can be reused in SMESH.
966 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
967 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
968 aVal = Round( vtkMeshQuality::TetAspectRatio( avtkCell ));
973 if ( GetPoints( myCurrElement, P ))
974 aVal = Round( GetValue( P ));
979 double AspectRatio3D::GetValue( const TSequenceOfXYZ& P )
981 double aQuality = 0.0;
982 if(myCurrElement->IsPoly()) return aQuality;
984 int nbNodes = P.size();
986 if(myCurrElement->IsQuadratic()) {
987 if(nbNodes==10) nbNodes=4; // quadratic tetrahedron
988 else if(nbNodes==13) nbNodes=5; // quadratic pyramid
989 else if(nbNodes==15) nbNodes=6; // quadratic pentahedron
990 else if(nbNodes==20) nbNodes=8; // quadratic hexahedron
991 else if(nbNodes==27) nbNodes=8; // quadratic hexahedron
992 else return aQuality;
998 getDistance(P( 1 ),P( 2 )), // a
999 getDistance(P( 2 ),P( 3 )), // b
1000 getDistance(P( 3 ),P( 1 )), // c
1001 getDistance(P( 2 ),P( 4 )), // d
1002 getDistance(P( 3 ),P( 4 )), // e
1003 getDistance(P( 1 ),P( 4 )) // f
1005 double aTria[4][3] = {
1006 {aLen[0],aLen[1],aLen[2]}, // abc
1007 {aLen[0],aLen[3],aLen[5]}, // adf
1008 {aLen[1],aLen[3],aLen[4]}, // bde
1009 {aLen[2],aLen[4],aLen[5]} // cef
1011 double aSumArea = 0.0;
1012 double aHalfPerimeter = getHalfPerimeter(aTria[0]);
1013 double anArea = getArea(aHalfPerimeter,aTria[0]);
1015 aHalfPerimeter = getHalfPerimeter(aTria[1]);
1016 anArea = getArea(aHalfPerimeter,aTria[1]);
1018 aHalfPerimeter = getHalfPerimeter(aTria[2]);
1019 anArea = getArea(aHalfPerimeter,aTria[2]);
1021 aHalfPerimeter = getHalfPerimeter(aTria[3]);
1022 anArea = getArea(aHalfPerimeter,aTria[3]);
1024 double aVolume = getVolume(P);
1025 //double aVolume = getVolume(aLen);
1026 double aHeight = getMaxHeight(aLen);
1027 static double aCoeff = sqrt(2.0)/12.0;
1028 if ( aVolume > DBL_MIN )
1029 aQuality = aCoeff*aHeight*aSumArea/aVolume;
1034 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 3 ),P( 5 )};
1035 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1038 gp_XYZ aXYZ[4] = {P( 1 ),P( 3 ),P( 4 ),P( 5 )};
1039 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1042 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 5 )};
1043 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1046 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 4 ),P( 5 )};
1047 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1053 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 6 )};
1054 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1057 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 3 )};
1058 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1061 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 6 )};
1062 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1065 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1066 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1069 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 6 )};
1070 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1073 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 3 )};
1074 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1080 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1081 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1084 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 4 )};
1085 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1088 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 7 )};
1089 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1092 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 8 )};
1093 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1096 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 3 )};
1097 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1100 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 4 )};
1101 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1104 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 7 )};
1105 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1108 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 8 )};
1109 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1112 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 3 )};
1113 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1116 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 4 )};
1117 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1120 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 7 )};
1121 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1124 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 8 )};
1125 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1128 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 1 )};
1129 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1132 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 2 )};
1133 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1136 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 5 )};
1137 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1140 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 6 )};
1141 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1144 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 1 )};
1145 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1148 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 2 )};
1149 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1152 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 5 )};
1153 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1156 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 6 )};
1157 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1160 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 1 )};
1161 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1164 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1165 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1168 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 5 )};
1169 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1172 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 6 )};
1173 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1176 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1177 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1180 gp_XYZ aXYZ[4] = {P( 4 ),P( 5 ),P( 8 ),P( 2 )};
1181 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1184 gp_XYZ aXYZ[4] = {P( 1 ),P( 4 ),P( 5 ),P( 3 )};
1185 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1188 gp_XYZ aXYZ[4] = {P( 3 ),P( 6 ),P( 7 ),P( 1 )};
1189 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1192 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 6 ),P( 4 )};
1193 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1196 gp_XYZ aXYZ[4] = {P( 5 ),P( 6 ),P( 8 ),P( 3 )};
1197 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1200 gp_XYZ aXYZ[4] = {P( 7 ),P( 8 ),P( 6 ),P( 1 )};
1201 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1204 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 7 )};
1205 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1208 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 2 ),P( 5 )};
1209 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1215 gp_XYZ aXYZ[8] = {P( 1 ),P( 2 ),P( 4 ),P( 5 ),P( 7 ),P( 8 ),P( 10 ),P( 11 )};
1216 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1219 gp_XYZ aXYZ[8] = {P( 2 ),P( 3 ),P( 5 ),P( 6 ),P( 8 ),P( 9 ),P( 11 ),P( 12 )};
1220 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1223 gp_XYZ aXYZ[8] = {P( 3 ),P( 4 ),P( 6 ),P( 1 ),P( 9 ),P( 10 ),P( 12 ),P( 7 )};
1224 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1227 } // switch(nbNodes)
1229 if ( nbNodes > 4 ) {
1230 // avaluate aspect ratio of quadranle faces
1231 AspectRatio aspect2D;
1232 SMDS_VolumeTool::VolumeType type = SMDS_VolumeTool::GetType( nbNodes );
1233 int nbFaces = SMDS_VolumeTool::NbFaces( type );
1234 TSequenceOfXYZ points(4);
1235 for ( int i = 0; i < nbFaces; ++i ) { // loop on faces of a volume
1236 if ( SMDS_VolumeTool::NbFaceNodes( type, i ) != 4 )
1238 const int* pInd = SMDS_VolumeTool::GetFaceNodesIndices( type, i, true );
1239 for ( int p = 0; p < 4; ++p ) // loop on nodes of a quadranle face
1240 points( p + 1 ) = P( pInd[ p ] + 1 );
1241 aQuality = std::max( aQuality, aspect2D.GetValue( points ));
1247 double AspectRatio3D::GetBadRate( double Value, int /*nbNodes*/ ) const
1249 // the aspect ratio is in the range [1.0,infinity]
1252 return Value / 1000.;
1255 SMDSAbs_ElementType AspectRatio3D::GetType() const
1257 return SMDSAbs_Volume;
1263 Description : Functor for calculating warping
1265 double Warping::GetValue( const TSequenceOfXYZ& P )
1267 if ( P.size() != 4 )
1270 gp_XYZ G = ( P( 1 ) + P( 2 ) + P( 3 ) + P( 4 ) ) / 4.;
1272 double A1 = ComputeA( P( 1 ), P( 2 ), P( 3 ), G );
1273 double A2 = ComputeA( P( 2 ), P( 3 ), P( 4 ), G );
1274 double A3 = ComputeA( P( 3 ), P( 4 ), P( 1 ), G );
1275 double A4 = ComputeA( P( 4 ), P( 1 ), P( 2 ), G );
1277 return Max( Max( A1, A2 ), Max( A3, A4 ) );
1280 double Warping::ComputeA( const gp_XYZ& thePnt1,
1281 const gp_XYZ& thePnt2,
1282 const gp_XYZ& thePnt3,
1283 const gp_XYZ& theG ) const
1285 double aLen1 = gp_Pnt( thePnt1 ).Distance( gp_Pnt( thePnt2 ) );
1286 double aLen2 = gp_Pnt( thePnt2 ).Distance( gp_Pnt( thePnt3 ) );
1287 double L = Min( aLen1, aLen2 ) * 0.5;
1288 if ( L < Precision::Confusion())
1291 gp_XYZ GI = ( thePnt2 + thePnt1 ) / 2. - theG;
1292 gp_XYZ GJ = ( thePnt3 + thePnt2 ) / 2. - theG;
1293 gp_XYZ N = GI.Crossed( GJ );
1295 if ( N.Modulus() < gp::Resolution() )
1300 double H = ( thePnt2 - theG ).Dot( N );
1301 return asin( fabs( H / L ) ) * 180. / M_PI;
1304 double Warping::GetBadRate( double Value, int /*nbNodes*/ ) const
1306 // the warp is in the range [0.0,PI/2]
1307 // 0.0 = good (no warp)
1308 // PI/2 = bad (face pliee)
1312 SMDSAbs_ElementType Warping::GetType() const
1314 return SMDSAbs_Face;
1320 Description : Functor for calculating taper
1322 double Taper::GetValue( const TSequenceOfXYZ& P )
1324 if ( P.size() != 4 )
1328 double J1 = getArea( P( 4 ), P( 1 ), P( 2 ) ) / 2.;
1329 double J2 = getArea( P( 3 ), P( 1 ), P( 2 ) ) / 2.;
1330 double J3 = getArea( P( 2 ), P( 3 ), P( 4 ) ) / 2.;
1331 double J4 = getArea( P( 3 ), P( 4 ), P( 1 ) ) / 2.;
1333 double JA = 0.25 * ( J1 + J2 + J3 + J4 );
1334 if ( JA <= Precision::Confusion() )
1337 double T1 = fabs( ( J1 - JA ) / JA );
1338 double T2 = fabs( ( J2 - JA ) / JA );
1339 double T3 = fabs( ( J3 - JA ) / JA );
1340 double T4 = fabs( ( J4 - JA ) / JA );
1342 return Max( Max( T1, T2 ), Max( T3, T4 ) );
1345 double Taper::GetBadRate( double Value, int /*nbNodes*/ ) const
1347 // the taper is in the range [0.0,1.0]
1348 // 0.0 = good (no taper)
1349 // 1.0 = bad (les cotes opposes sont allignes)
1353 SMDSAbs_ElementType Taper::GetType() const
1355 return SMDSAbs_Face;
1361 Description : Functor for calculating skew in degrees
1363 static inline double skewAngle( const gp_XYZ& p1, const gp_XYZ& p2, const gp_XYZ& p3 )
1365 gp_XYZ p12 = ( p2 + p1 ) / 2.;
1366 gp_XYZ p23 = ( p3 + p2 ) / 2.;
1367 gp_XYZ p31 = ( p3 + p1 ) / 2.;
1369 gp_Vec v1( p31 - p2 ), v2( p12 - p23 );
1371 return v1.Magnitude() < gp::Resolution() || v2.Magnitude() < gp::Resolution() ? 0. : v1.Angle( v2 );
1374 double Skew::GetValue( const TSequenceOfXYZ& P )
1376 if ( P.size() != 3 && P.size() != 4 )
1380 static double PI2 = M_PI / 2.;
1381 if ( P.size() == 3 )
1383 double A0 = fabs( PI2 - skewAngle( P( 3 ), P( 1 ), P( 2 ) ) );
1384 double A1 = fabs( PI2 - skewAngle( P( 1 ), P( 2 ), P( 3 ) ) );
1385 double A2 = fabs( PI2 - skewAngle( P( 2 ), P( 3 ), P( 1 ) ) );
1387 return Max( A0, Max( A1, A2 ) ) * 180. / M_PI;
1391 gp_XYZ p12 = ( P( 1 ) + P( 2 ) ) / 2.;
1392 gp_XYZ p23 = ( P( 2 ) + P( 3 ) ) / 2.;
1393 gp_XYZ p34 = ( P( 3 ) + P( 4 ) ) / 2.;
1394 gp_XYZ p41 = ( P( 4 ) + P( 1 ) ) / 2.;
1396 gp_Vec v1( p34 - p12 ), v2( p23 - p41 );
1397 double A = v1.Magnitude() <= gp::Resolution() || v2.Magnitude() <= gp::Resolution()
1398 ? 0. : fabs( PI2 - v1.Angle( v2 ) );
1401 if ( A < Precision::Angular() )
1404 return A * 180. / M_PI;
1408 double Skew::GetBadRate( double Value, int /*nbNodes*/ ) const
1410 // the skew is in the range [0.0,PI/2].
1416 SMDSAbs_ElementType Skew::GetType() const
1418 return SMDSAbs_Face;
1424 Description : Functor for calculating area
1426 double Area::GetValue( const TSequenceOfXYZ& P )
1429 if ( P.size() > 2 ) {
1430 gp_Vec aVec1( P(2) - P(1) );
1431 gp_Vec aVec2( P(3) - P(1) );
1432 gp_Vec SumVec = aVec1 ^ aVec2;
1433 for (int i=4; i<=P.size(); i++) {
1434 gp_Vec aVec1( P(i-1) - P(1) );
1435 gp_Vec aVec2( P(i) - P(1) );
1436 gp_Vec tmp = aVec1 ^ aVec2;
1439 val = SumVec.Magnitude() * 0.5;
1444 double Area::GetBadRate( double Value, int /*nbNodes*/ ) const
1446 // meaningless as it is not a quality control functor
1450 SMDSAbs_ElementType Area::GetType() const
1452 return SMDSAbs_Face;
1458 Description : Functor for calculating length of edge
1460 double Length::GetValue( const TSequenceOfXYZ& P )
1462 switch ( P.size() ) {
1463 case 2: return getDistance( P( 1 ), P( 2 ) );
1464 case 3: return getDistance( P( 1 ), P( 2 ) ) + getDistance( P( 2 ), P( 3 ) );
1469 double Length::GetBadRate( double Value, int /*nbNodes*/ ) const
1471 // meaningless as it is not quality control functor
1475 SMDSAbs_ElementType Length::GetType() const
1477 return SMDSAbs_Edge;
1482 Description : Functor for calculating length of edge
1485 double Length2D::GetValue( long theElementId)
1489 //cout<<"Length2D::GetValue"<<endl;
1490 if (GetPoints(theElementId,P)){
1491 //for(int jj=1; jj<=P.size(); jj++)
1492 // cout<<"jj="<<jj<<" P("<<P(jj).X()<<","<<P(jj).Y()<<","<<P(jj).Z()<<")"<<endl;
1494 double aVal;// = GetValue( P );
1495 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
1496 SMDSAbs_ElementType aType = aElem->GetType();
1505 aVal = getDistance( P( 1 ), P( 2 ) );
1508 else if (len == 3){ // quadratic edge
1509 aVal = getDistance(P( 1 ),P( 3 )) + getDistance(P( 3 ),P( 2 ));
1513 if (len == 3){ // triangles
1514 double L1 = getDistance(P( 1 ),P( 2 ));
1515 double L2 = getDistance(P( 2 ),P( 3 ));
1516 double L3 = getDistance(P( 3 ),P( 1 ));
1517 aVal = Max(L1,Max(L2,L3));
1520 else if (len == 4){ // quadrangles
1521 double L1 = getDistance(P( 1 ),P( 2 ));
1522 double L2 = getDistance(P( 2 ),P( 3 ));
1523 double L3 = getDistance(P( 3 ),P( 4 ));
1524 double L4 = getDistance(P( 4 ),P( 1 ));
1525 aVal = Max(Max(L1,L2),Max(L3,L4));
1528 if (len == 6){ // quadratic triangles
1529 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1530 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1531 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
1532 aVal = Max(L1,Max(L2,L3));
1533 //cout<<"L1="<<L1<<" L2="<<L2<<"L3="<<L3<<" aVal="<<aVal<<endl;
1536 else if (len == 8){ // quadratic quadrangles
1537 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1538 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1539 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
1540 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
1541 aVal = Max(Max(L1,L2),Max(L3,L4));
1544 case SMDSAbs_Volume:
1545 if (len == 4){ // tetraidrs
1546 double L1 = getDistance(P( 1 ),P( 2 ));
1547 double L2 = getDistance(P( 2 ),P( 3 ));
1548 double L3 = getDistance(P( 3 ),P( 1 ));
1549 double L4 = getDistance(P( 1 ),P( 4 ));
1550 double L5 = getDistance(P( 2 ),P( 4 ));
1551 double L6 = getDistance(P( 3 ),P( 4 ));
1552 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1555 else if (len == 5){ // piramids
1556 double L1 = getDistance(P( 1 ),P( 2 ));
1557 double L2 = getDistance(P( 2 ),P( 3 ));
1558 double L3 = getDistance(P( 3 ),P( 4 ));
1559 double L4 = getDistance(P( 4 ),P( 1 ));
1560 double L5 = getDistance(P( 1 ),P( 5 ));
1561 double L6 = getDistance(P( 2 ),P( 5 ));
1562 double L7 = getDistance(P( 3 ),P( 5 ));
1563 double L8 = getDistance(P( 4 ),P( 5 ));
1565 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1566 aVal = Max(aVal,Max(L7,L8));
1569 else if (len == 6){ // pentaidres
1570 double L1 = getDistance(P( 1 ),P( 2 ));
1571 double L2 = getDistance(P( 2 ),P( 3 ));
1572 double L3 = getDistance(P( 3 ),P( 1 ));
1573 double L4 = getDistance(P( 4 ),P( 5 ));
1574 double L5 = getDistance(P( 5 ),P( 6 ));
1575 double L6 = getDistance(P( 6 ),P( 4 ));
1576 double L7 = getDistance(P( 1 ),P( 4 ));
1577 double L8 = getDistance(P( 2 ),P( 5 ));
1578 double L9 = getDistance(P( 3 ),P( 6 ));
1580 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1581 aVal = Max(aVal,Max(Max(L7,L8),L9));
1584 else if (len == 8){ // hexaider
1585 double L1 = getDistance(P( 1 ),P( 2 ));
1586 double L2 = getDistance(P( 2 ),P( 3 ));
1587 double L3 = getDistance(P( 3 ),P( 4 ));
1588 double L4 = getDistance(P( 4 ),P( 1 ));
1589 double L5 = getDistance(P( 5 ),P( 6 ));
1590 double L6 = getDistance(P( 6 ),P( 7 ));
1591 double L7 = getDistance(P( 7 ),P( 8 ));
1592 double L8 = getDistance(P( 8 ),P( 5 ));
1593 double L9 = getDistance(P( 1 ),P( 5 ));
1594 double L10= getDistance(P( 2 ),P( 6 ));
1595 double L11= getDistance(P( 3 ),P( 7 ));
1596 double L12= getDistance(P( 4 ),P( 8 ));
1598 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1599 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1600 aVal = Max(aVal,Max(L11,L12));
1605 if (len == 10){ // quadratic tetraidrs
1606 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
1607 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
1608 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
1609 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1610 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
1611 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
1612 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1615 else if (len == 13){ // quadratic piramids
1616 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
1617 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
1618 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1619 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1620 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1621 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
1622 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
1623 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
1624 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1625 aVal = Max(aVal,Max(L7,L8));
1628 else if (len == 15){ // quadratic pentaidres
1629 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
1630 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
1631 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1632 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1633 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
1634 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
1635 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
1636 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
1637 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
1638 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1639 aVal = Max(aVal,Max(Max(L7,L8),L9));
1642 else if (len == 20){ // quadratic hexaider
1643 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
1644 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
1645 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
1646 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
1647 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
1648 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
1649 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
1650 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
1651 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
1652 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
1653 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
1654 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
1655 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1656 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1657 aVal = Max(aVal,Max(L11,L12));
1669 if ( myPrecision >= 0 )
1671 double prec = pow( 10., (double)( myPrecision ) );
1672 aVal = floor( aVal * prec + 0.5 ) / prec;
1681 double Length2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1683 // meaningless as it is not quality control functor
1687 SMDSAbs_ElementType Length2D::GetType() const
1689 return SMDSAbs_Face;
1692 Length2D::Value::Value(double theLength,long thePntId1, long thePntId2):
1695 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1696 if(thePntId1 > thePntId2){
1697 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1701 bool Length2D::Value::operator<(const Length2D::Value& x) const{
1702 if(myPntId[0] < x.myPntId[0]) return true;
1703 if(myPntId[0] == x.myPntId[0])
1704 if(myPntId[1] < x.myPntId[1]) return true;
1708 void Length2D::GetValues(TValues& theValues){
1710 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1711 for(; anIter->more(); ){
1712 const SMDS_MeshFace* anElem = anIter->next();
1714 if(anElem->IsQuadratic()) {
1715 const SMDS_VtkFace* F =
1716 dynamic_cast<const SMDS_VtkFace*>(anElem);
1717 // use special nodes iterator
1718 SMDS_ElemIteratorPtr anIter = F->interlacedNodesElemIterator();
1723 const SMDS_MeshElement* aNode;
1725 aNode = anIter->next();
1726 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1727 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1728 aNodeId[0] = aNodeId[1] = aNode->GetID();
1731 for(; anIter->more(); ){
1732 const SMDS_MeshNode* N1 = static_cast<const SMDS_MeshNode*> (anIter->next());
1733 P[2] = gp_Pnt(N1->X(),N1->Y(),N1->Z());
1734 aNodeId[2] = N1->GetID();
1735 aLength = P[1].Distance(P[2]);
1736 if(!anIter->more()) break;
1737 const SMDS_MeshNode* N2 = static_cast<const SMDS_MeshNode*> (anIter->next());
1738 P[3] = gp_Pnt(N2->X(),N2->Y(),N2->Z());
1739 aNodeId[3] = N2->GetID();
1740 aLength += P[2].Distance(P[3]);
1741 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1742 Value aValue2(aLength,aNodeId[2],aNodeId[3]);
1744 aNodeId[1] = aNodeId[3];
1745 theValues.insert(aValue1);
1746 theValues.insert(aValue2);
1748 aLength += P[2].Distance(P[0]);
1749 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1750 Value aValue2(aLength,aNodeId[2],aNodeId[0]);
1751 theValues.insert(aValue1);
1752 theValues.insert(aValue2);
1755 SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
1760 const SMDS_MeshElement* aNode;
1761 if(aNodesIter->more()){
1762 aNode = aNodesIter->next();
1763 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1764 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1765 aNodeId[0] = aNodeId[1] = aNode->GetID();
1768 for(; aNodesIter->more(); ){
1769 aNode = aNodesIter->next();
1770 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1771 long anId = aNode->GetID();
1773 P[2] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1775 aLength = P[1].Distance(P[2]);
1777 Value aValue(aLength,aNodeId[1],anId);
1780 theValues.insert(aValue);
1783 aLength = P[0].Distance(P[1]);
1785 Value aValue(aLength,aNodeId[0],aNodeId[1]);
1786 theValues.insert(aValue);
1792 Class : MultiConnection
1793 Description : Functor for calculating number of faces conneted to the edge
1795 double MultiConnection::GetValue( const TSequenceOfXYZ& P )
1799 double MultiConnection::GetValue( long theId )
1801 return getNbMultiConnection( myMesh, theId );
1804 double MultiConnection::GetBadRate( double Value, int /*nbNodes*/ ) const
1806 // meaningless as it is not quality control functor
1810 SMDSAbs_ElementType MultiConnection::GetType() const
1812 return SMDSAbs_Edge;
1816 Class : MultiConnection2D
1817 Description : Functor for calculating number of faces conneted to the edge
1819 double MultiConnection2D::GetValue( const TSequenceOfXYZ& P )
1824 double MultiConnection2D::GetValue( long theElementId )
1828 const SMDS_MeshElement* aFaceElem = myMesh->FindElement(theElementId);
1829 SMDSAbs_ElementType aType = aFaceElem->GetType();
1834 int i = 0, len = aFaceElem->NbNodes();
1835 SMDS_ElemIteratorPtr anIter = aFaceElem->nodesIterator();
1838 const SMDS_MeshNode *aNode, *aNode0;
1839 TColStd_MapOfInteger aMap, aMapPrev;
1841 for (i = 0; i <= len; i++) {
1846 if (anIter->more()) {
1847 aNode = (SMDS_MeshNode*)anIter->next();
1855 if (i == 0) aNode0 = aNode;
1857 SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
1858 while (anElemIter->more()) {
1859 const SMDS_MeshElement* anElem = anElemIter->next();
1860 if (anElem != 0 && anElem->GetType() == SMDSAbs_Face) {
1861 int anId = anElem->GetID();
1864 if (aMapPrev.Contains(anId)) {
1869 aResult = Max(aResult, aNb);
1880 double MultiConnection2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1882 // meaningless as it is not quality control functor
1886 SMDSAbs_ElementType MultiConnection2D::GetType() const
1888 return SMDSAbs_Face;
1891 MultiConnection2D::Value::Value(long thePntId1, long thePntId2)
1893 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1894 if(thePntId1 > thePntId2){
1895 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1899 bool MultiConnection2D::Value::operator<(const MultiConnection2D::Value& x) const{
1900 if(myPntId[0] < x.myPntId[0]) return true;
1901 if(myPntId[0] == x.myPntId[0])
1902 if(myPntId[1] < x.myPntId[1]) return true;
1906 void MultiConnection2D::GetValues(MValues& theValues){
1907 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1908 for(; anIter->more(); ){
1909 const SMDS_MeshFace* anElem = anIter->next();
1910 SMDS_ElemIteratorPtr aNodesIter;
1911 if ( anElem->IsQuadratic() )
1912 aNodesIter = dynamic_cast<const SMDS_VtkFace*>
1913 (anElem)->interlacedNodesElemIterator();
1915 aNodesIter = anElem->nodesIterator();
1918 //int aNbConnects=0;
1919 const SMDS_MeshNode* aNode0;
1920 const SMDS_MeshNode* aNode1;
1921 const SMDS_MeshNode* aNode2;
1922 if(aNodesIter->more()){
1923 aNode0 = (SMDS_MeshNode*) aNodesIter->next();
1925 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode1;
1926 aNodeId[0] = aNodeId[1] = aNodes->GetID();
1928 for(; aNodesIter->more(); ) {
1929 aNode2 = (SMDS_MeshNode*) aNodesIter->next();
1930 long anId = aNode2->GetID();
1933 Value aValue(aNodeId[1],aNodeId[2]);
1934 MValues::iterator aItr = theValues.find(aValue);
1935 if (aItr != theValues.end()){
1940 theValues[aValue] = 1;
1943 //cout << "NodeIds: "<<aNodeId[1]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1944 aNodeId[1] = aNodeId[2];
1947 Value aValue(aNodeId[0],aNodeId[2]);
1948 MValues::iterator aItr = theValues.find(aValue);
1949 if (aItr != theValues.end()) {
1954 theValues[aValue] = 1;
1957 //cout << "NodeIds: "<<aNodeId[0]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1963 Class : BallDiameter
1964 Description : Functor returning diameter of a ball element
1966 double BallDiameter::GetValue( long theId )
1968 double diameter = 0;
1970 if ( const SMDS_BallElement* ball =
1971 dynamic_cast<const SMDS_BallElement*>( myMesh->FindElement( theId )))
1973 diameter = ball->GetDiameter();
1978 double BallDiameter::GetBadRate( double Value, int /*nbNodes*/ ) const
1980 // meaningless as it is not a quality control functor
1984 SMDSAbs_ElementType BallDiameter::GetType() const
1986 return SMDSAbs_Ball;
1995 Class : BadOrientedVolume
1996 Description : Predicate bad oriented volumes
1999 BadOrientedVolume::BadOrientedVolume()
2004 void BadOrientedVolume::SetMesh( const SMDS_Mesh* theMesh )
2009 bool BadOrientedVolume::IsSatisfy( long theId )
2014 SMDS_VolumeTool vTool( myMesh->FindElement( theId ));
2015 return !vTool.IsForward();
2018 SMDSAbs_ElementType BadOrientedVolume::GetType() const
2020 return SMDSAbs_Volume;
2024 Class : BareBorderVolume
2027 bool BareBorderVolume::IsSatisfy(long theElementId )
2029 SMDS_VolumeTool myTool;
2030 if ( myTool.Set( myMesh->FindElement(theElementId)))
2032 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2033 if ( myTool.IsFreeFace( iF ))
2035 const SMDS_MeshNode** n = myTool.GetFaceNodes(iF);
2036 vector< const SMDS_MeshNode*> nodes( n, n+myTool.NbFaceNodes(iF));
2037 if ( !myMesh->FindElement( nodes, SMDSAbs_Face, /*Nomedium=*/false))
2045 Class : BareBorderFace
2048 bool BareBorderFace::IsSatisfy(long theElementId )
2051 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2053 if ( face->GetType() == SMDSAbs_Face )
2055 int nbN = face->NbCornerNodes();
2056 for ( int i = 0; i < nbN && !ok; ++i )
2058 // check if a link is shared by another face
2059 const SMDS_MeshNode* n1 = face->GetNode( i );
2060 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2061 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2062 bool isShared = false;
2063 while ( !isShared && fIt->more() )
2065 const SMDS_MeshElement* f = fIt->next();
2066 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2070 const int iQuad = face->IsQuadratic();
2071 myLinkNodes.resize( 2 + iQuad);
2072 myLinkNodes[0] = n1;
2073 myLinkNodes[1] = n2;
2075 myLinkNodes[2] = face->GetNode( i+nbN );
2076 ok = !myMesh->FindElement( myLinkNodes, SMDSAbs_Edge, /*noMedium=*/false);
2085 Class : OverConstrainedVolume
2088 bool OverConstrainedVolume::IsSatisfy(long theElementId )
2090 // An element is over-constrained if it has N-1 free borders where
2091 // N is the number of edges/faces for a 2D/3D element.
2092 SMDS_VolumeTool myTool;
2093 if ( myTool.Set( myMesh->FindElement(theElementId)))
2095 int nbSharedFaces = 0;
2096 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2097 if ( !myTool.IsFreeFace( iF ) && ++nbSharedFaces > 1 )
2099 return ( nbSharedFaces == 1 );
2105 Class : OverConstrainedFace
2108 bool OverConstrainedFace::IsSatisfy(long theElementId )
2110 // An element is over-constrained if it has N-1 free borders where
2111 // N is the number of edges/faces for a 2D/3D element.
2112 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2113 if ( face->GetType() == SMDSAbs_Face )
2115 int nbSharedBorders = 0;
2116 int nbN = face->NbCornerNodes();
2117 for ( int i = 0; i < nbN; ++i )
2119 // check if a link is shared by another face
2120 const SMDS_MeshNode* n1 = face->GetNode( i );
2121 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2122 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2123 bool isShared = false;
2124 while ( !isShared && fIt->more() )
2126 const SMDS_MeshElement* f = fIt->next();
2127 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2129 if ( isShared && ++nbSharedBorders > 1 )
2132 return ( nbSharedBorders == 1 );
2138 Class : CoincidentNodes
2139 Description : Predicate of Coincident nodes
2142 CoincidentNodes::CoincidentNodes()
2147 bool CoincidentNodes::IsSatisfy( long theElementId )
2149 return myCoincidentIDs.Contains( theElementId );
2152 SMDSAbs_ElementType CoincidentNodes::GetType() const
2154 return SMDSAbs_Node;
2157 void CoincidentNodes::SetMesh( const SMDS_Mesh* theMesh )
2159 myMeshModifTracer.SetMesh( theMesh );
2160 if ( myMeshModifTracer.IsMeshModified() )
2162 TIDSortedNodeSet nodesToCheck;
2163 SMDS_NodeIteratorPtr nIt = theMesh->nodesIterator(/*idInceasingOrder=*/true);
2164 while ( nIt->more() )
2165 nodesToCheck.insert( nodesToCheck.end(), nIt->next() );
2167 list< list< const SMDS_MeshNode*> > nodeGroups;
2168 SMESH_OctreeNode::FindCoincidentNodes ( nodesToCheck, &nodeGroups, myToler );
2170 myCoincidentIDs.Clear();
2171 list< list< const SMDS_MeshNode*> >::iterator groupIt = nodeGroups.begin();
2172 for ( ; groupIt != nodeGroups.end(); ++groupIt )
2174 list< const SMDS_MeshNode*>& coincNodes = *groupIt;
2175 list< const SMDS_MeshNode*>::iterator n = coincNodes.begin();
2176 for ( ; n != coincNodes.end(); ++n )
2177 myCoincidentIDs.Add( (*n)->GetID() );
2183 Class : CoincidentElements
2184 Description : Predicate of Coincident Elements
2185 Note : This class is suitable only for visualization of Coincident Elements
2188 CoincidentElements::CoincidentElements()
2193 void CoincidentElements::SetMesh( const SMDS_Mesh* theMesh )
2198 bool CoincidentElements::IsSatisfy( long theElementId )
2200 if ( !myMesh ) return false;
2202 if ( const SMDS_MeshElement* e = myMesh->FindElement( theElementId ))
2204 if ( e->GetType() != GetType() ) return false;
2205 set< const SMDS_MeshNode* > elemNodes( e->begin_nodes(), e->end_nodes() );
2206 const int nbNodes = e->NbNodes();
2207 SMDS_ElemIteratorPtr invIt = (*elemNodes.begin())->GetInverseElementIterator( GetType() );
2208 while ( invIt->more() )
2210 const SMDS_MeshElement* e2 = invIt->next();
2211 if ( e2 == e || e2->NbNodes() != nbNodes ) continue;
2213 bool sameNodes = true;
2214 for ( size_t i = 0; i < elemNodes.size() && sameNodes; ++i )
2215 sameNodes = ( elemNodes.count( e2->GetNode( i )));
2223 SMDSAbs_ElementType CoincidentElements1D::GetType() const
2225 return SMDSAbs_Edge;
2227 SMDSAbs_ElementType CoincidentElements2D::GetType() const
2229 return SMDSAbs_Face;
2231 SMDSAbs_ElementType CoincidentElements3D::GetType() const
2233 return SMDSAbs_Volume;
2239 Description : Predicate for free borders
2242 FreeBorders::FreeBorders()
2247 void FreeBorders::SetMesh( const SMDS_Mesh* theMesh )
2252 bool FreeBorders::IsSatisfy( long theId )
2254 return getNbMultiConnection( myMesh, theId ) == 1;
2257 SMDSAbs_ElementType FreeBorders::GetType() const
2259 return SMDSAbs_Edge;
2265 Description : Predicate for free Edges
2267 FreeEdges::FreeEdges()
2272 void FreeEdges::SetMesh( const SMDS_Mesh* theMesh )
2277 bool FreeEdges::IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId )
2279 TColStd_MapOfInteger aMap;
2280 for ( int i = 0; i < 2; i++ )
2282 SMDS_ElemIteratorPtr anElemIter = theNodes[ i ]->GetInverseElementIterator(SMDSAbs_Face);
2283 while( anElemIter->more() )
2285 if ( const SMDS_MeshElement* anElem = anElemIter->next())
2287 const int anId = anElem->GetID();
2288 if ( anId != theFaceId && !aMap.Add( anId ))
2296 bool FreeEdges::IsSatisfy( long theId )
2301 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2302 if ( aFace == 0 || aFace->GetType() != SMDSAbs_Face || aFace->NbNodes() < 3 )
2305 SMDS_ElemIteratorPtr anIter;
2306 if ( aFace->IsQuadratic() ) {
2307 anIter = dynamic_cast<const SMDS_VtkFace*>
2308 (aFace)->interlacedNodesElemIterator();
2311 anIter = aFace->nodesIterator();
2316 int i = 0, nbNodes = aFace->NbNodes();
2317 std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
2318 while( anIter->more() )
2320 const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
2323 aNodes[ i++ ] = aNode;
2325 aNodes[ nbNodes ] = aNodes[ 0 ];
2327 for ( i = 0; i < nbNodes; i++ )
2328 if ( IsFreeEdge( &aNodes[ i ], theId ) )
2334 SMDSAbs_ElementType FreeEdges::GetType() const
2336 return SMDSAbs_Face;
2339 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2342 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
2343 if(thePntId1 > thePntId2){
2344 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
2348 bool FreeEdges::Border::operator<(const FreeEdges::Border& x) const{
2349 if(myPntId[0] < x.myPntId[0]) return true;
2350 if(myPntId[0] == x.myPntId[0])
2351 if(myPntId[1] < x.myPntId[1]) return true;
2355 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2356 FreeEdges::TBorders& theRegistry,
2357 FreeEdges::TBorders& theContainer)
2359 if(theRegistry.find(theBorder) == theRegistry.end()){
2360 theRegistry.insert(theBorder);
2361 theContainer.insert(theBorder);
2363 theContainer.erase(theBorder);
2367 void FreeEdges::GetBoreders(TBorders& theBorders)
2370 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2371 for(; anIter->more(); ){
2372 const SMDS_MeshFace* anElem = anIter->next();
2373 long anElemId = anElem->GetID();
2374 SMDS_ElemIteratorPtr aNodesIter;
2375 if ( anElem->IsQuadratic() )
2376 aNodesIter = static_cast<const SMDS_VtkFace*>(anElem)->
2377 interlacedNodesElemIterator();
2379 aNodesIter = anElem->nodesIterator();
2381 const SMDS_MeshElement* aNode;
2382 if(aNodesIter->more()){
2383 aNode = aNodesIter->next();
2384 aNodeId[0] = aNodeId[1] = aNode->GetID();
2386 for(; aNodesIter->more(); ){
2387 aNode = aNodesIter->next();
2388 long anId = aNode->GetID();
2389 Border aBorder(anElemId,aNodeId[1],anId);
2391 UpdateBorders(aBorder,aRegistry,theBorders);
2393 Border aBorder(anElemId,aNodeId[0],aNodeId[1]);
2394 UpdateBorders(aBorder,aRegistry,theBorders);
2401 Description : Predicate for free nodes
2404 FreeNodes::FreeNodes()
2409 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2414 bool FreeNodes::IsSatisfy( long theNodeId )
2416 const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2420 return (aNode->NbInverseElements() < 1);
2423 SMDSAbs_ElementType FreeNodes::GetType() const
2425 return SMDSAbs_Node;
2431 Description : Predicate for free faces
2434 FreeFaces::FreeFaces()
2439 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2444 bool FreeFaces::IsSatisfy( long theId )
2446 if (!myMesh) return false;
2447 // check that faces nodes refers to less than two common volumes
2448 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2449 if ( !aFace || aFace->GetType() != SMDSAbs_Face )
2452 int nbNode = aFace->NbNodes();
2454 // collect volumes check that number of volumss with count equal nbNode not less than 2
2455 typedef map< SMDS_MeshElement*, int > TMapOfVolume; // map of volume counters
2456 typedef map< SMDS_MeshElement*, int >::iterator TItrMapOfVolume; // iterator
2457 TMapOfVolume mapOfVol;
2459 SMDS_ElemIteratorPtr nodeItr = aFace->nodesIterator();
2460 while ( nodeItr->more() ) {
2461 const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>(nodeItr->next());
2462 if ( !aNode ) continue;
2463 SMDS_ElemIteratorPtr volItr = aNode->GetInverseElementIterator(SMDSAbs_Volume);
2464 while ( volItr->more() ) {
2465 SMDS_MeshElement* aVol = (SMDS_MeshElement*)volItr->next();
2466 TItrMapOfVolume itr = mapOfVol.insert(make_pair(aVol, 0)).first;
2471 TItrMapOfVolume volItr = mapOfVol.begin();
2472 TItrMapOfVolume volEnd = mapOfVol.end();
2473 for ( ; volItr != volEnd; ++volItr )
2474 if ( (*volItr).second >= nbNode )
2476 // face is not free if number of volumes constructed on thier nodes more than one
2480 SMDSAbs_ElementType FreeFaces::GetType() const
2482 return SMDSAbs_Face;
2486 Class : LinearOrQuadratic
2487 Description : Predicate to verify whether a mesh element is linear
2490 LinearOrQuadratic::LinearOrQuadratic()
2495 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2500 bool LinearOrQuadratic::IsSatisfy( long theId )
2502 if (!myMesh) return false;
2503 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2504 if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2506 return (!anElem->IsQuadratic());
2509 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2514 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2521 Description : Functor for check color of group to whic mesh element belongs to
2524 GroupColor::GroupColor()
2528 bool GroupColor::IsSatisfy( long theId )
2530 return (myIDs.find( theId ) != myIDs.end());
2533 void GroupColor::SetType( SMDSAbs_ElementType theType )
2538 SMDSAbs_ElementType GroupColor::GetType() const
2543 static bool isEqual( const Quantity_Color& theColor1,
2544 const Quantity_Color& theColor2 )
2546 // tolerance to compare colors
2547 const double tol = 5*1e-3;
2548 return ( fabs( theColor1.Red() - theColor2.Red() ) < tol &&
2549 fabs( theColor1.Green() - theColor2.Green() ) < tol &&
2550 fabs( theColor1.Blue() - theColor2.Blue() ) < tol );
2554 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2558 const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2562 int nbGrp = aMesh->GetNbGroups();
2566 // iterates on groups and find necessary elements ids
2567 const std::set<SMESHDS_GroupBase*>& aGroups = aMesh->GetGroups();
2568 set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
2569 for (; GrIt != aGroups.end(); GrIt++) {
2570 SMESHDS_GroupBase* aGrp = (*GrIt);
2573 // check type and color of group
2574 if ( !isEqual( myColor, aGrp->GetColor() ) )
2576 if ( myType != SMDSAbs_All && myType != (SMDSAbs_ElementType)aGrp->GetType() )
2579 SMDSAbs_ElementType aGrpElType = (SMDSAbs_ElementType)aGrp->GetType();
2580 if ( myType == aGrpElType || (myType == SMDSAbs_All && aGrpElType != SMDSAbs_Node) ) {
2581 // add elements IDS into control
2582 int aSize = aGrp->Extent();
2583 for (int i = 0; i < aSize; i++)
2584 myIDs.insert( aGrp->GetID(i+1) );
2589 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
2591 TCollection_AsciiString aStr = theStr;
2592 aStr.RemoveAll( ' ' );
2593 aStr.RemoveAll( '\t' );
2594 for ( int aPos = aStr.Search( ";;" ); aPos != -1; aPos = aStr.Search( ";;" ) )
2595 aStr.Remove( aPos, 2 );
2596 Standard_Real clr[3];
2597 clr[0] = clr[1] = clr[2] = 0.;
2598 for ( int i = 0; i < 3; i++ ) {
2599 TCollection_AsciiString tmpStr = aStr.Token( ";", i+1 );
2600 if ( !tmpStr.IsEmpty() && tmpStr.IsRealValue() )
2601 clr[i] = tmpStr.RealValue();
2603 myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
2606 //=======================================================================
2607 // name : GetRangeStr
2608 // Purpose : Get range as a string.
2609 // Example: "1,2,3,50-60,63,67,70-"
2610 //=======================================================================
2611 void GroupColor::GetColorStr( TCollection_AsciiString& theResStr ) const
2614 theResStr += TCollection_AsciiString( myColor.Red() );
2615 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2616 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2620 Class : ElemGeomType
2621 Description : Predicate to check element geometry type
2624 ElemGeomType::ElemGeomType()
2627 myType = SMDSAbs_All;
2628 myGeomType = SMDSGeom_TRIANGLE;
2631 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2636 bool ElemGeomType::IsSatisfy( long theId )
2638 if (!myMesh) return false;
2639 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2642 const SMDSAbs_ElementType anElemType = anElem->GetType();
2643 if ( myType != SMDSAbs_All && anElemType != myType )
2645 bool isOk = ( anElem->GetGeomType() == myGeomType );
2649 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2654 SMDSAbs_ElementType ElemGeomType::GetType() const
2659 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2661 myGeomType = theType;
2664 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2669 //================================================================================
2671 * \brief Class CoplanarFaces
2673 //================================================================================
2675 CoplanarFaces::CoplanarFaces()
2676 : myFaceID(0), myToler(0)
2679 void CoplanarFaces::SetMesh( const SMDS_Mesh* theMesh )
2681 myMeshModifTracer.SetMesh( theMesh );
2682 if ( myMeshModifTracer.IsMeshModified() )
2684 // Build a set of coplanar face ids
2686 myCoplanarIDs.clear();
2688 if ( !myMeshModifTracer.GetMesh() || !myFaceID || !myToler )
2691 const SMDS_MeshElement* face = myMeshModifTracer.GetMesh()->FindElement( myFaceID );
2692 if ( !face || face->GetType() != SMDSAbs_Face )
2696 gp_Vec myNorm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2700 const double radianTol = myToler * M_PI / 180.;
2701 std::set< SMESH_TLink > checkedLinks;
2703 std::list< pair< const SMDS_MeshElement*, gp_Vec > > faceQueue;
2704 faceQueue.push_back( make_pair( face, myNorm ));
2705 while ( !faceQueue.empty() )
2707 face = faceQueue.front().first;
2708 myNorm = faceQueue.front().second;
2709 faceQueue.pop_front();
2711 for ( int i = 0, nbN = face->NbCornerNodes(); i < nbN; ++i )
2713 const SMDS_MeshNode* n1 = face->GetNode( i );
2714 const SMDS_MeshNode* n2 = face->GetNode(( i+1 )%nbN);
2715 if ( !checkedLinks.insert( SMESH_TLink( n1, n2 )).second )
2717 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator(SMDSAbs_Face);
2718 while ( fIt->more() )
2720 const SMDS_MeshElement* f = fIt->next();
2721 if ( f->GetNodeIndex( n2 ) > -1 )
2723 gp_Vec norm = getNormale( static_cast<const SMDS_MeshFace*>(f), &normOK );
2724 if (!normOK || myNorm.Angle( norm ) <= radianTol)
2726 myCoplanarIDs.insert( f->GetID() );
2727 faceQueue.push_back( make_pair( f, norm ));
2735 bool CoplanarFaces::IsSatisfy( long theElementId )
2737 return myCoplanarIDs.count( theElementId );
2742 *Description : Predicate for Range of Ids.
2743 * Range may be specified with two ways.
2744 * 1. Using AddToRange method
2745 * 2. With SetRangeStr method. Parameter of this method is a string
2746 * like as "1,2,3,50-60,63,67,70-"
2749 //=======================================================================
2750 // name : RangeOfIds
2751 // Purpose : Constructor
2752 //=======================================================================
2753 RangeOfIds::RangeOfIds()
2756 myType = SMDSAbs_All;
2759 //=======================================================================
2761 // Purpose : Set mesh
2762 //=======================================================================
2763 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
2768 //=======================================================================
2769 // name : AddToRange
2770 // Purpose : Add ID to the range
2771 //=======================================================================
2772 bool RangeOfIds::AddToRange( long theEntityId )
2774 myIds.Add( theEntityId );
2778 //=======================================================================
2779 // name : GetRangeStr
2780 // Purpose : Get range as a string.
2781 // Example: "1,2,3,50-60,63,67,70-"
2782 //=======================================================================
2783 void RangeOfIds::GetRangeStr( TCollection_AsciiString& theResStr )
2787 TColStd_SequenceOfInteger anIntSeq;
2788 TColStd_SequenceOfAsciiString aStrSeq;
2790 TColStd_MapIteratorOfMapOfInteger anIter( myIds );
2791 for ( ; anIter.More(); anIter.Next() )
2793 int anId = anIter.Key();
2794 TCollection_AsciiString aStr( anId );
2795 anIntSeq.Append( anId );
2796 aStrSeq.Append( aStr );
2799 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2801 int aMinId = myMin( i );
2802 int aMaxId = myMax( i );
2804 TCollection_AsciiString aStr;
2805 if ( aMinId != IntegerFirst() )
2810 if ( aMaxId != IntegerLast() )
2813 // find position of the string in result sequence and insert string in it
2814 if ( anIntSeq.Length() == 0 )
2816 anIntSeq.Append( aMinId );
2817 aStrSeq.Append( aStr );
2821 if ( aMinId < anIntSeq.First() )
2823 anIntSeq.Prepend( aMinId );
2824 aStrSeq.Prepend( aStr );
2826 else if ( aMinId > anIntSeq.Last() )
2828 anIntSeq.Append( aMinId );
2829 aStrSeq.Append( aStr );
2832 for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
2833 if ( aMinId < anIntSeq( j ) )
2835 anIntSeq.InsertBefore( j, aMinId );
2836 aStrSeq.InsertBefore( j, aStr );
2842 if ( aStrSeq.Length() == 0 )
2845 theResStr = aStrSeq( 1 );
2846 for ( int j = 2, k = aStrSeq.Length(); j <= k; j++ )
2849 theResStr += aStrSeq( j );
2853 //=======================================================================
2854 // name : SetRangeStr
2855 // Purpose : Define range with string
2856 // Example of entry string: "1,2,3,50-60,63,67,70-"
2857 //=======================================================================
2858 bool RangeOfIds::SetRangeStr( const TCollection_AsciiString& theStr )
2864 TCollection_AsciiString aStr = theStr;
2865 aStr.RemoveAll( ' ' );
2866 aStr.RemoveAll( '\t' );
2868 for ( int aPos = aStr.Search( ",," ); aPos != -1; aPos = aStr.Search( ",," ) )
2869 aStr.Remove( aPos, 2 );
2871 TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
2873 while ( tmpStr != "" )
2875 tmpStr = aStr.Token( ",", i++ );
2876 int aPos = tmpStr.Search( '-' );
2880 if ( tmpStr.IsIntegerValue() )
2881 myIds.Add( tmpStr.IntegerValue() );
2887 TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
2888 TCollection_AsciiString aMinStr = tmpStr;
2890 while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
2891 while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
2893 if ( (!aMinStr.IsEmpty() && !aMinStr.IsIntegerValue()) ||
2894 (!aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue()) )
2897 myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
2898 myMax.Append( aMaxStr.IsEmpty() ? IntegerLast() : aMaxStr.IntegerValue() );
2905 //=======================================================================
2907 // Purpose : Get type of supported entities
2908 //=======================================================================
2909 SMDSAbs_ElementType RangeOfIds::GetType() const
2914 //=======================================================================
2916 // Purpose : Set type of supported entities
2917 //=======================================================================
2918 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
2923 //=======================================================================
2925 // Purpose : Verify whether entity satisfies to this rpedicate
2926 //=======================================================================
2927 bool RangeOfIds::IsSatisfy( long theId )
2932 if ( myType == SMDSAbs_Node )
2934 if ( myMesh->FindNode( theId ) == 0 )
2939 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2940 if ( anElem == 0 || (myType != anElem->GetType() && myType != SMDSAbs_All ))
2944 if ( myIds.Contains( theId ) )
2947 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2948 if ( theId >= myMin( i ) && theId <= myMax( i ) )
2956 Description : Base class for comparators
2958 Comparator::Comparator():
2962 Comparator::~Comparator()
2965 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
2968 myFunctor->SetMesh( theMesh );
2971 void Comparator::SetMargin( double theValue )
2973 myMargin = theValue;
2976 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
2978 myFunctor = theFunct;
2981 SMDSAbs_ElementType Comparator::GetType() const
2983 return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
2986 double Comparator::GetMargin()
2994 Description : Comparator "<"
2996 bool LessThan::IsSatisfy( long theId )
2998 return myFunctor && myFunctor->GetValue( theId ) < myMargin;
3004 Description : Comparator ">"
3006 bool MoreThan::IsSatisfy( long theId )
3008 return myFunctor && myFunctor->GetValue( theId ) > myMargin;
3014 Description : Comparator "="
3017 myToler(Precision::Confusion())
3020 bool EqualTo::IsSatisfy( long theId )
3022 return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
3025 void EqualTo::SetTolerance( double theToler )
3030 double EqualTo::GetTolerance()
3037 Description : Logical NOT predicate
3039 LogicalNOT::LogicalNOT()
3042 LogicalNOT::~LogicalNOT()
3045 bool LogicalNOT::IsSatisfy( long theId )
3047 return myPredicate && !myPredicate->IsSatisfy( theId );
3050 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
3053 myPredicate->SetMesh( theMesh );
3056 void LogicalNOT::SetPredicate( PredicatePtr thePred )
3058 myPredicate = thePred;
3061 SMDSAbs_ElementType LogicalNOT::GetType() const
3063 return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
3068 Class : LogicalBinary
3069 Description : Base class for binary logical predicate
3071 LogicalBinary::LogicalBinary()
3074 LogicalBinary::~LogicalBinary()
3077 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
3080 myPredicate1->SetMesh( theMesh );
3083 myPredicate2->SetMesh( theMesh );
3086 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
3088 myPredicate1 = thePredicate;
3091 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
3093 myPredicate2 = thePredicate;
3096 SMDSAbs_ElementType LogicalBinary::GetType() const
3098 if ( !myPredicate1 || !myPredicate2 )
3101 SMDSAbs_ElementType aType1 = myPredicate1->GetType();
3102 SMDSAbs_ElementType aType2 = myPredicate2->GetType();
3104 return aType1 == aType2 ? aType1 : SMDSAbs_All;
3110 Description : Logical AND
3112 bool LogicalAND::IsSatisfy( long theId )
3117 myPredicate1->IsSatisfy( theId ) &&
3118 myPredicate2->IsSatisfy( theId );
3124 Description : Logical OR
3126 bool LogicalOR::IsSatisfy( long theId )
3131 (myPredicate1->IsSatisfy( theId ) ||
3132 myPredicate2->IsSatisfy( theId ));
3146 void Filter::SetPredicate( PredicatePtr thePredicate )
3148 myPredicate = thePredicate;
3151 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3152 PredicatePtr thePredicate,
3153 TIdSequence& theSequence )
3155 theSequence.clear();
3157 if ( !theMesh || !thePredicate )
3160 thePredicate->SetMesh( theMesh );
3162 SMDS_ElemIteratorPtr elemIt = theMesh->elementsIterator( thePredicate->GetType() );
3164 while ( elemIt->more() ) {
3165 const SMDS_MeshElement* anElem = elemIt->next();
3166 long anId = anElem->GetID();
3167 if ( thePredicate->IsSatisfy( anId ) )
3168 theSequence.push_back( anId );
3173 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3174 Filter::TIdSequence& theSequence )
3176 GetElementsId(theMesh,myPredicate,theSequence);
3183 typedef std::set<SMDS_MeshFace*> TMapOfFacePtr;
3189 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
3190 SMDS_MeshNode* theNode2 )
3196 ManifoldPart::Link::~Link()
3202 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
3204 if ( myNode1 == theLink.myNode1 &&
3205 myNode2 == theLink.myNode2 )
3207 else if ( myNode1 == theLink.myNode2 &&
3208 myNode2 == theLink.myNode1 )
3214 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
3216 if(myNode1 < x.myNode1) return true;
3217 if(myNode1 == x.myNode1)
3218 if(myNode2 < x.myNode2) return true;
3222 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
3223 const ManifoldPart::Link& theLink2 )
3225 return theLink1.IsEqual( theLink2 );
3228 ManifoldPart::ManifoldPart()
3231 myAngToler = Precision::Angular();
3232 myIsOnlyManifold = true;
3235 ManifoldPart::~ManifoldPart()
3240 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
3246 SMDSAbs_ElementType ManifoldPart::GetType() const
3247 { return SMDSAbs_Face; }
3249 bool ManifoldPart::IsSatisfy( long theElementId )
3251 return myMapIds.Contains( theElementId );
3254 void ManifoldPart::SetAngleTolerance( const double theAngToler )
3255 { myAngToler = theAngToler; }
3257 double ManifoldPart::GetAngleTolerance() const
3258 { return myAngToler; }
3260 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
3261 { myIsOnlyManifold = theIsOnly; }
3263 void ManifoldPart::SetStartElem( const long theStartId )
3264 { myStartElemId = theStartId; }
3266 bool ManifoldPart::process()
3269 myMapBadGeomIds.Clear();
3271 myAllFacePtr.clear();
3272 myAllFacePtrIntDMap.clear();
3276 // collect all faces into own map
3277 SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3278 for (; anFaceItr->more(); )
3280 SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3281 myAllFacePtr.push_back( aFacePtr );
3282 myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3285 SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3289 // the map of non manifold links and bad geometry
3290 TMapOfLink aMapOfNonManifold;
3291 TColStd_MapOfInteger aMapOfTreated;
3293 // begin cycle on faces from start index and run on vector till the end
3294 // and from begin to start index to cover whole vector
3295 const int aStartIndx = myAllFacePtrIntDMap[aStartFace];
3296 bool isStartTreat = false;
3297 for ( int fi = aStartIndx; !isStartTreat || fi != aStartIndx ; fi++ )
3299 if ( fi == aStartIndx )
3300 isStartTreat = true;
3301 // as result next time when fi will be equal to aStartIndx
3303 SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3304 if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3307 aMapOfTreated.Add( aFacePtr->GetID() );
3308 TColStd_MapOfInteger aResFaces;
3309 if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3310 aMapOfNonManifold, aResFaces ) )
3312 TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3313 for ( ; anItr.More(); anItr.Next() )
3315 int aFaceId = anItr.Key();
3316 aMapOfTreated.Add( aFaceId );
3317 myMapIds.Add( aFaceId );
3320 if ( fi == ( myAllFacePtr.size() - 1 ) )
3322 } // end run on vector of faces
3323 return !myMapIds.IsEmpty();
3326 static void getLinks( const SMDS_MeshFace* theFace,
3327 ManifoldPart::TVectorOfLink& theLinks )
3329 int aNbNode = theFace->NbNodes();
3330 SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3332 SMDS_MeshNode* aNode = 0;
3333 for ( ; aNodeItr->more() && i <= aNbNode; )
3336 SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3340 SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3342 ManifoldPart::Link aLink( aN1, aN2 );
3343 theLinks.push_back( aLink );
3347 bool ManifoldPart::findConnected
3348 ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3349 SMDS_MeshFace* theStartFace,
3350 ManifoldPart::TMapOfLink& theNonManifold,
3351 TColStd_MapOfInteger& theResFaces )
3353 theResFaces.Clear();
3354 if ( !theAllFacePtrInt.size() )
3357 if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3359 myMapBadGeomIds.Add( theStartFace->GetID() );
3363 ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3364 ManifoldPart::TVectorOfLink aSeqOfBoundary;
3365 theResFaces.Add( theStartFace->GetID() );
3366 ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3368 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3369 aDMapLinkFace, theNonManifold, theStartFace );
3371 bool isDone = false;
3372 while ( !isDone && aMapOfBoundary.size() != 0 )
3374 bool isToReset = false;
3375 ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3376 for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3378 ManifoldPart::Link aLink = *pLink;
3379 if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3381 // each link could be treated only once
3382 aMapToSkip.insert( aLink );
3384 ManifoldPart::TVectorOfFacePtr aFaces;
3386 if ( myIsOnlyManifold &&
3387 (theNonManifold.find( aLink ) != theNonManifold.end()) )
3391 getFacesByLink( aLink, aFaces );
3392 // filter the element to keep only indicated elements
3393 ManifoldPart::TVectorOfFacePtr aFiltered;
3394 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3395 for ( ; pFace != aFaces.end(); ++pFace )
3397 SMDS_MeshFace* aFace = *pFace;
3398 if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3399 aFiltered.push_back( aFace );
3402 if ( aFaces.size() < 2 ) // no neihgbour faces
3404 else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3406 theNonManifold.insert( aLink );
3411 // compare normal with normals of neighbor element
3412 SMDS_MeshFace* aPrevFace = aDMapLinkFace[ aLink ];
3413 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3414 for ( ; pFace != aFaces.end(); ++pFace )
3416 SMDS_MeshFace* aNextFace = *pFace;
3417 if ( aPrevFace == aNextFace )
3419 int anNextFaceID = aNextFace->GetID();
3420 if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3421 // should not be with non manifold restriction. probably bad topology
3423 // check if face was treated and skipped
3424 if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3425 !isInPlane( aPrevFace, aNextFace ) )
3427 // add new element to connected and extend the boundaries.
3428 theResFaces.Add( anNextFaceID );
3429 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3430 aDMapLinkFace, theNonManifold, aNextFace );
3434 isDone = !isToReset;
3437 return !theResFaces.IsEmpty();
3440 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3441 const SMDS_MeshFace* theFace2 )
3443 gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3444 gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3445 if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3447 myMapBadGeomIds.Add( theFace2->GetID() );
3450 if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
3456 void ManifoldPart::expandBoundary
3457 ( ManifoldPart::TMapOfLink& theMapOfBoundary,
3458 ManifoldPart::TVectorOfLink& theSeqOfBoundary,
3459 ManifoldPart::TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
3460 ManifoldPart::TMapOfLink& theNonManifold,
3461 SMDS_MeshFace* theNextFace ) const
3463 ManifoldPart::TVectorOfLink aLinks;
3464 getLinks( theNextFace, aLinks );
3465 int aNbLink = (int)aLinks.size();
3466 for ( int i = 0; i < aNbLink; i++ )
3468 ManifoldPart::Link aLink = aLinks[ i ];
3469 if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3471 if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3473 if ( myIsOnlyManifold )
3475 // remove from boundary
3476 theMapOfBoundary.erase( aLink );
3477 ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3478 for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3480 ManifoldPart::Link aBoundLink = *pLink;
3481 if ( aBoundLink.IsEqual( aLink ) )
3483 theSeqOfBoundary.erase( pLink );
3491 theMapOfBoundary.insert( aLink );
3492 theSeqOfBoundary.push_back( aLink );
3493 theDMapLinkFacePtr[ aLink ] = theNextFace;
3498 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3499 ManifoldPart::TVectorOfFacePtr& theFaces ) const
3501 std::set<SMDS_MeshCell *> aSetOfFaces;
3502 // take all faces that shared first node
3503 SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3504 for ( ; anItr->more(); )
3506 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3509 aSetOfFaces.insert( aFace );
3511 // take all faces that shared second node
3512 anItr = theLink.myNode2->facesIterator();
3513 // find the common part of two sets
3514 for ( ; anItr->more(); )
3516 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3517 if ( aSetOfFaces.count( aFace ) )
3518 theFaces.push_back( aFace );
3527 ElementsOnSurface::ElementsOnSurface()
3530 myType = SMDSAbs_All;
3532 myToler = Precision::Confusion();
3533 myUseBoundaries = false;
3536 ElementsOnSurface::~ElementsOnSurface()
3540 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
3542 myMeshModifTracer.SetMesh( theMesh );
3543 if ( myMeshModifTracer.IsMeshModified())
3547 bool ElementsOnSurface::IsSatisfy( long theElementId )
3549 return myIds.Contains( theElementId );
3552 SMDSAbs_ElementType ElementsOnSurface::GetType() const
3555 void ElementsOnSurface::SetTolerance( const double theToler )
3557 if ( myToler != theToler )
3562 double ElementsOnSurface::GetTolerance() const
3565 void ElementsOnSurface::SetUseBoundaries( bool theUse )
3567 if ( myUseBoundaries != theUse ) {
3568 myUseBoundaries = theUse;
3569 SetSurface( mySurf, myType );
3573 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
3574 const SMDSAbs_ElementType theType )
3579 if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
3581 mySurf = TopoDS::Face( theShape );
3582 BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
3584 u1 = SA.FirstUParameter(),
3585 u2 = SA.LastUParameter(),
3586 v1 = SA.FirstVParameter(),
3587 v2 = SA.LastVParameter();
3588 Handle(Geom_Surface) surf = BRep_Tool::Surface( mySurf );
3589 myProjector.Init( surf, u1,u2, v1,v2 );
3593 void ElementsOnSurface::process()
3596 if ( mySurf.IsNull() )
3599 if ( !myMeshModifTracer.GetMesh() )
3602 myIds.ReSize( myMeshModifTracer.GetMesh()->GetMeshInfo().NbElements( myType ));
3604 SMDS_ElemIteratorPtr anIter = myMeshModifTracer.GetMesh()->elementsIterator( myType );
3605 for(; anIter->more(); )
3606 process( anIter->next() );
3609 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
3611 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3612 bool isSatisfy = true;
3613 for ( ; aNodeItr->more(); )
3615 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3616 if ( !isOnSurface( aNode ) )
3623 myIds.Add( theElemPtr->GetID() );
3626 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
3628 if ( mySurf.IsNull() )
3631 gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
3632 // double aToler2 = myToler * myToler;
3633 // if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
3635 // gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
3636 // if ( aPln.SquareDistance( aPnt ) > aToler2 )
3639 // else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
3641 // gp_Cylinder aCyl = Handle(Geom_CylindricalSurface)::DownCast(mySurf)->Cylinder();
3642 // double aRad = aCyl.Radius();
3643 // gp_Ax3 anAxis = aCyl.Position();
3644 // gp_XYZ aLoc = aCyl.Location().XYZ();
3645 // double aXDist = anAxis.XDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3646 // double aYDist = anAxis.YDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3647 // if ( fabs(aXDist*aXDist + aYDist*aYDist - aRad*aRad) > aToler2 )
3652 myProjector.Perform( aPnt );
3653 bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
3663 ElementsOnShape::ElementsOnShape()
3665 myType(SMDSAbs_All),
3666 myToler(Precision::Confusion()),
3667 myAllNodesFlag(false)
3669 myCurShapeType = TopAbs_SHAPE;
3672 ElementsOnShape::~ElementsOnShape()
3676 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
3678 myMeshModifTracer.SetMesh( theMesh );
3679 if ( myMeshModifTracer.IsMeshModified())
3680 SetShape(myShape, myType);
3683 bool ElementsOnShape::IsSatisfy (long theElementId)
3685 return myIds.Contains(theElementId);
3688 SMDSAbs_ElementType ElementsOnShape::GetType() const
3693 void ElementsOnShape::SetTolerance (const double theToler)
3695 if (myToler != theToler) {
3697 SetShape(myShape, myType);
3701 double ElementsOnShape::GetTolerance() const
3706 void ElementsOnShape::SetAllNodes (bool theAllNodes)
3708 if (myAllNodesFlag != theAllNodes) {
3709 myAllNodesFlag = theAllNodes;
3710 SetShape(myShape, myType);
3714 void ElementsOnShape::SetShape (const TopoDS_Shape& theShape,
3715 const SMDSAbs_ElementType theType)
3721 const SMDS_Mesh* myMesh = myMeshModifTracer.GetMesh();
3723 if ( !myMesh ) return;
3725 myIds.ReSize( myMeshModifTracer.GetMesh()->GetMeshInfo().NbElements( myType ));
3727 myShapesMap.Clear();
3731 void ElementsOnShape::addShape (const TopoDS_Shape& theShape)
3733 if (theShape.IsNull() || myMeshModifTracer.GetMesh() == 0)
3736 if (!myShapesMap.Add(theShape)) return;
3738 myCurShapeType = theShape.ShapeType();
3739 switch (myCurShapeType)
3741 case TopAbs_COMPOUND:
3742 case TopAbs_COMPSOLID:
3746 TopoDS_Iterator anIt (theShape, Standard_True, Standard_True);
3747 for (; anIt.More(); anIt.Next()) addShape(anIt.Value());
3752 myCurSC.Load(theShape);
3758 TopoDS_Face aFace = TopoDS::Face(theShape);
3759 BRepAdaptor_Surface SA (aFace, true);
3761 u1 = SA.FirstUParameter(),
3762 u2 = SA.LastUParameter(),
3763 v1 = SA.FirstVParameter(),
3764 v2 = SA.LastVParameter();
3765 Handle(Geom_Surface) surf = BRep_Tool::Surface(aFace);
3766 myCurProjFace.Init(surf, u1,u2, v1,v2);
3773 TopoDS_Edge anEdge = TopoDS::Edge(theShape);
3774 Standard_Real u1, u2;
3775 Handle(Geom_Curve) curve = BRep_Tool::Curve(anEdge, u1, u2);
3776 myCurProjEdge.Init(curve, u1, u2);
3782 TopoDS_Vertex aV = TopoDS::Vertex(theShape);
3783 myCurPnt = BRep_Tool::Pnt(aV);
3792 void ElementsOnShape::process()
3794 const SMDS_Mesh* myMesh = myMeshModifTracer.GetMesh();
3795 if (myShape.IsNull() || myMesh == 0)
3798 SMDS_ElemIteratorPtr anIter = myMesh->elementsIterator(myType);
3799 while (anIter->more())
3800 process(anIter->next());
3803 void ElementsOnShape::process (const SMDS_MeshElement* theElemPtr)
3805 if (myShape.IsNull())
3808 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3809 bool isSatisfy = myAllNodesFlag;
3811 gp_XYZ centerXYZ (0, 0, 0);
3813 while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
3815 SMESH_TNodeXYZ aPnt ( aNodeItr->next() );
3818 switch (myCurShapeType)
3822 myCurSC.Perform(aPnt, myToler);
3823 isSatisfy = (myCurSC.State() == TopAbs_IN || myCurSC.State() == TopAbs_ON);
3828 myCurProjFace.Perform(aPnt);
3829 isSatisfy = (myCurProjFace.IsDone() && myCurProjFace.LowerDistance() <= myToler);
3832 // check relatively the face
3833 Quantity_Parameter u, v;
3834 myCurProjFace.LowerDistanceParameters(u, v);
3835 gp_Pnt2d aProjPnt (u, v);
3836 BRepClass_FaceClassifier aClsf (myCurFace, aProjPnt, myToler);
3837 isSatisfy = (aClsf.State() == TopAbs_IN || aClsf.State() == TopAbs_ON);
3843 myCurProjEdge.Perform(aPnt);
3844 isSatisfy = (myCurProjEdge.NbPoints() > 0 && myCurProjEdge.LowerDistance() <= myToler);
3849 isSatisfy = (myCurPnt.Distance(aPnt) <= myToler);
3859 if (isSatisfy && myCurShapeType == TopAbs_SOLID) { // Check the center point for volumes MantisBug 0020168
3860 centerXYZ /= theElemPtr->NbNodes();
3861 gp_Pnt aCenterPnt (centerXYZ);
3862 myCurSC.Perform(aCenterPnt, myToler);
3863 if ( !(myCurSC.State() == TopAbs_IN || myCurSC.State() == TopAbs_ON))
3868 myIds.Add(theElemPtr->GetID());
3871 TSequenceOfXYZ::TSequenceOfXYZ()
3874 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n)
3877 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t)
3880 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray)
3883 template <class InputIterator>
3884 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd)
3887 TSequenceOfXYZ::~TSequenceOfXYZ()
3890 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
3892 myArray = theSequenceOfXYZ.myArray;
3896 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
3898 return myArray[n-1];
3901 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
3903 return myArray[n-1];
3906 void TSequenceOfXYZ::clear()
3911 void TSequenceOfXYZ::reserve(size_type n)
3916 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
3918 myArray.push_back(v);
3921 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
3923 return myArray.size();
3926 TMeshModifTracer::TMeshModifTracer():
3927 myMeshModifTime(0), myMesh(0)
3930 void TMeshModifTracer::SetMesh( const SMDS_Mesh* theMesh )
3932 if ( theMesh != myMesh )
3933 myMeshModifTime = 0;
3936 bool TMeshModifTracer::IsMeshModified()
3938 bool modified = false;
3941 modified = ( myMeshModifTime != myMesh->GetMTime() );
3942 myMeshModifTime = myMesh->GetMTime();