1 // Copyright (C) 2007-2011 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
23 #include "SMESH_ControlsDef.hxx"
28 #include <BRepAdaptor_Surface.hxx>
29 #include <BRepClass_FaceClassifier.hxx>
30 #include <BRep_Tool.hxx>
34 #include <TopoDS_Edge.hxx>
35 #include <TopoDS_Face.hxx>
36 #include <TopoDS_Shape.hxx>
37 #include <TopoDS_Vertex.hxx>
38 #include <TopoDS_Iterator.hxx>
40 #include <Geom_CylindricalSurface.hxx>
41 #include <Geom_Plane.hxx>
42 #include <Geom_Surface.hxx>
44 #include <Precision.hxx>
45 #include <TColStd_MapIteratorOfMapOfInteger.hxx>
46 #include <TColStd_MapOfInteger.hxx>
47 #include <TColStd_SequenceOfAsciiString.hxx>
48 #include <TColgp_Array1OfXYZ.hxx>
51 #include <gp_Cylinder.hxx>
58 #include "SMDS_Mesh.hxx"
59 #include "SMDS_Iterator.hxx"
60 #include "SMDS_MeshElement.hxx"
61 #include "SMDS_MeshNode.hxx"
62 #include "SMDS_VolumeTool.hxx"
63 #include "SMDS_QuadraticFaceOfNodes.hxx"
64 #include "SMDS_QuadraticEdge.hxx"
66 #include "SMESHDS_Mesh.hxx"
67 #include "SMESHDS_GroupBase.hxx"
69 #include <vtkMeshQuality.h>
77 inline gp_XYZ gpXYZ(const SMDS_MeshNode* aNode )
79 return gp_XYZ(aNode->X(), aNode->Y(), aNode->Z() );
82 inline double getAngle( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
84 gp_Vec v1( P1 - P2 ), v2( P3 - P2 );
86 return v1.Magnitude() < gp::Resolution() ||
87 v2.Magnitude() < gp::Resolution() ? 0 : v1.Angle( v2 );
90 inline double getArea( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
92 gp_Vec aVec1( P2 - P1 );
93 gp_Vec aVec2( P3 - P1 );
94 return ( aVec1 ^ aVec2 ).Magnitude() * 0.5;
97 inline double getArea( const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3 )
99 return getArea( P1.XYZ(), P2.XYZ(), P3.XYZ() );
104 inline double getDistance( const gp_XYZ& P1, const gp_XYZ& P2 )
106 double aDist = gp_Pnt( P1 ).Distance( gp_Pnt( P2 ) );
110 int getNbMultiConnection( const SMDS_Mesh* theMesh, const int theId )
115 const SMDS_MeshElement* anEdge = theMesh->FindElement( theId );
116 if ( anEdge == 0 || anEdge->GetType() != SMDSAbs_Edge/* || anEdge->NbNodes() != 2 */)
119 // for each pair of nodes in anEdge (there are 2 pairs in a quadratic edge)
120 // count elements containing both nodes of the pair.
121 // Note that there may be such cases for a quadratic edge (a horizontal line):
126 // +-----+------+ +-----+------+
129 // result sould be 2 in both cases
131 int aResult0 = 0, aResult1 = 0;
132 // last node, it is a medium one in a quadratic edge
133 const SMDS_MeshNode* aLastNode = anEdge->GetNode( anEdge->NbNodes() - 1 );
134 const SMDS_MeshNode* aNode0 = anEdge->GetNode( 0 );
135 const SMDS_MeshNode* aNode1 = anEdge->GetNode( 1 );
136 if ( aNode1 == aLastNode ) aNode1 = 0;
138 SMDS_ElemIteratorPtr anElemIter = aLastNode->GetInverseElementIterator();
139 while( anElemIter->more() ) {
140 const SMDS_MeshElement* anElem = anElemIter->next();
141 if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
142 SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
143 while ( anIter->more() ) {
144 if ( const SMDS_MeshElement* anElemNode = anIter->next() ) {
145 if ( anElemNode == aNode0 ) {
147 if ( !aNode1 ) break; // not a quadratic edge
149 else if ( anElemNode == aNode1 )
155 int aResult = std::max ( aResult0, aResult1 );
157 // TColStd_MapOfInteger aMap;
159 // SMDS_ElemIteratorPtr anIter = anEdge->nodesIterator();
160 // if ( anIter != 0 ) {
161 // while( anIter->more() ) {
162 // const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
165 // SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
166 // while( anElemIter->more() ) {
167 // const SMDS_MeshElement* anElem = anElemIter->next();
168 // if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
169 // int anId = anElem->GetID();
171 // if ( anIter->more() ) // i.e. first node
173 // else if ( aMap.Contains( anId ) )
183 gp_XYZ getNormale( const SMDS_MeshFace* theFace, bool* ok=0 )
185 int aNbNode = theFace->NbNodes();
187 gp_XYZ q1 = gpXYZ( theFace->GetNode(1)) - gpXYZ( theFace->GetNode(0));
188 gp_XYZ q2 = gpXYZ( theFace->GetNode(2)) - gpXYZ( theFace->GetNode(0));
191 gp_XYZ q3 = gpXYZ( theFace->GetNode(3)) - gpXYZ( theFace->GetNode(0));
194 double len = n.Modulus();
195 bool zeroLen = ( len <= numeric_limits<double>::min());
199 if (ok) *ok = !zeroLen;
207 using namespace SMESH::Controls;
214 Class : NumericalFunctor
215 Description : Base class for numerical functors
217 NumericalFunctor::NumericalFunctor():
223 void NumericalFunctor::SetMesh( const SMDS_Mesh* theMesh )
228 bool NumericalFunctor::GetPoints(const int theId,
229 TSequenceOfXYZ& theRes ) const
236 return GetPoints( myMesh->FindElement( theId ), theRes );
239 bool NumericalFunctor::GetPoints(const SMDS_MeshElement* anElem,
240 TSequenceOfXYZ& theRes )
247 theRes.reserve( anElem->NbNodes() );
249 // Get nodes of the element
250 SMDS_ElemIteratorPtr anIter;
252 if ( anElem->IsQuadratic() ) {
253 switch ( anElem->GetType() ) {
255 anIter = dynamic_cast<const SMDS_VtkEdge*>
256 (anElem)->interlacedNodesElemIterator();
259 anIter = dynamic_cast<const SMDS_VtkFace*>
260 (anElem)->interlacedNodesElemIterator();
263 anIter = anElem->nodesIterator();
268 anIter = anElem->nodesIterator();
272 while( anIter->more() ) {
273 if ( const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>( anIter->next() ))
274 theRes.push_back( gp_XYZ( aNode->X(), aNode->Y(), aNode->Z() ) );
281 long NumericalFunctor::GetPrecision() const
286 void NumericalFunctor::SetPrecision( const long thePrecision )
288 myPrecision = thePrecision;
289 myPrecisionValue = pow( 10., (double)( myPrecision ) );
292 double NumericalFunctor::GetValue( long theId )
296 myCurrElement = myMesh->FindElement( theId );
299 if ( GetPoints( theId, P ))
300 aVal = Round( GetValue( P ));
305 double NumericalFunctor::Round( const double & aVal )
307 return ( myPrecision >= 0 ) ? floor( aVal * myPrecisionValue + 0.5 ) / myPrecisionValue : aVal;
310 //================================================================================
312 * \brief Return histogram of functor values
313 * \param nbIntervals - number of intervals
314 * \param nbEvents - number of mesh elements having values within i-th interval
315 * \param funValues - boundaries of intervals
316 * \param elements - elements to check vulue of; empty list means "of all"
317 * \param minmax - boundaries of diapason of values to divide into intervals
319 //================================================================================
321 void NumericalFunctor::GetHistogram(int nbIntervals,
322 std::vector<int>& nbEvents,
323 std::vector<double>& funValues,
324 const vector<int>& elements,
325 const double* minmax)
327 if ( nbIntervals < 1 ||
329 !myMesh->GetMeshInfo().NbElements( GetType() ))
331 nbEvents.resize( nbIntervals, 0 );
332 funValues.resize( nbIntervals+1 );
334 // get all values sorted
335 std::multiset< double > values;
336 if ( elements.empty() )
338 SMDS_ElemIteratorPtr elemIt = myMesh->elementsIterator(GetType());
339 while ( elemIt->more() )
340 values.insert( GetValue( elemIt->next()->GetID() ));
344 vector<int>::const_iterator id = elements.begin();
345 for ( ; id != elements.end(); ++id )
346 values.insert( GetValue( *id ));
351 funValues[0] = minmax[0];
352 funValues[nbIntervals] = minmax[1];
356 funValues[0] = *values.begin();
357 funValues[nbIntervals] = *values.rbegin();
359 // case nbIntervals == 1
360 if ( nbIntervals == 1 )
362 nbEvents[0] = values.size();
366 if (funValues.front() == funValues.back())
368 nbEvents.resize( 1 );
369 nbEvents[0] = values.size();
370 funValues[1] = funValues.back();
371 funValues.resize( 2 );
374 std::multiset< double >::iterator min = values.begin(), max;
375 for ( int i = 0; i < nbIntervals; ++i )
377 // find end value of i-th interval
378 double r = (i+1) / double( nbIntervals );
379 funValues[i+1] = funValues.front() * (1-r) + funValues.back() * r;
381 // count values in the i-th interval if there are any
382 if ( min != values.end() && *min <= funValues[i+1] )
384 // find the first value out of the interval
385 max = values.upper_bound( funValues[i+1] ); // max is greater than funValues[i+1], or end()
386 nbEvents[i] = std::distance( min, max );
390 // add values larger than minmax[1]
391 nbEvents.back() += std::distance( min, values.end() );
394 //=======================================================================
395 //function : GetValue
397 //=======================================================================
399 double Volume::GetValue( long theElementId )
401 if ( theElementId && myMesh ) {
402 SMDS_VolumeTool aVolumeTool;
403 if ( aVolumeTool.Set( myMesh->FindElement( theElementId )))
404 return aVolumeTool.GetSize();
409 //=======================================================================
410 //function : GetBadRate
411 //purpose : meaningless as it is not quality control functor
412 //=======================================================================
414 double Volume::GetBadRate( double Value, int /*nbNodes*/ ) const
419 //=======================================================================
422 //=======================================================================
424 SMDSAbs_ElementType Volume::GetType() const
426 return SMDSAbs_Volume;
431 Class : MaxElementLength2D
432 Description : Functor calculating maximum length of 2D element
435 double MaxElementLength2D::GetValue( long theElementId )
438 if( GetPoints( theElementId, P ) ) {
440 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
441 SMDSAbs_ElementType aType = aElem->GetType();
445 if( len == 3 ) { // triangles
446 double L1 = getDistance(P( 1 ),P( 2 ));
447 double L2 = getDistance(P( 2 ),P( 3 ));
448 double L3 = getDistance(P( 3 ),P( 1 ));
449 aVal = Max(L1,Max(L2,L3));
452 else if( len == 4 ) { // quadrangles
453 double L1 = getDistance(P( 1 ),P( 2 ));
454 double L2 = getDistance(P( 2 ),P( 3 ));
455 double L3 = getDistance(P( 3 ),P( 4 ));
456 double L4 = getDistance(P( 4 ),P( 1 ));
457 double D1 = getDistance(P( 1 ),P( 3 ));
458 double D2 = getDistance(P( 2 ),P( 4 ));
459 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
462 else if( len == 6 ) { // quadratic triangles
463 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
464 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
465 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
466 aVal = Max(L1,Max(L2,L3));
469 else if( len == 8 ) { // quadratic quadrangles
470 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
471 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
472 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
473 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
474 double D1 = getDistance(P( 1 ),P( 5 ));
475 double D2 = getDistance(P( 3 ),P( 7 ));
476 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
481 if( myPrecision >= 0 )
483 double prec = pow( 10., (double)myPrecision );
484 aVal = floor( aVal * prec + 0.5 ) / prec;
491 double MaxElementLength2D::GetBadRate( double Value, int /*nbNodes*/ ) const
496 SMDSAbs_ElementType MaxElementLength2D::GetType() const
502 Class : MaxElementLength3D
503 Description : Functor calculating maximum length of 3D element
506 double MaxElementLength3D::GetValue( long theElementId )
509 if( GetPoints( theElementId, P ) ) {
511 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
512 SMDSAbs_ElementType aType = aElem->GetType();
516 if( len == 4 ) { // tetras
517 double L1 = getDistance(P( 1 ),P( 2 ));
518 double L2 = getDistance(P( 2 ),P( 3 ));
519 double L3 = getDistance(P( 3 ),P( 1 ));
520 double L4 = getDistance(P( 1 ),P( 4 ));
521 double L5 = getDistance(P( 2 ),P( 4 ));
522 double L6 = getDistance(P( 3 ),P( 4 ));
523 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
526 else if( len == 5 ) { // pyramids
527 double L1 = getDistance(P( 1 ),P( 2 ));
528 double L2 = getDistance(P( 2 ),P( 3 ));
529 double L3 = getDistance(P( 3 ),P( 4 ));
530 double L4 = getDistance(P( 4 ),P( 1 ));
531 double L5 = getDistance(P( 1 ),P( 5 ));
532 double L6 = getDistance(P( 2 ),P( 5 ));
533 double L7 = getDistance(P( 3 ),P( 5 ));
534 double L8 = getDistance(P( 4 ),P( 5 ));
535 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
536 aVal = Max(aVal,Max(L7,L8));
539 else if( len == 6 ) { // pentas
540 double L1 = getDistance(P( 1 ),P( 2 ));
541 double L2 = getDistance(P( 2 ),P( 3 ));
542 double L3 = getDistance(P( 3 ),P( 1 ));
543 double L4 = getDistance(P( 4 ),P( 5 ));
544 double L5 = getDistance(P( 5 ),P( 6 ));
545 double L6 = getDistance(P( 6 ),P( 4 ));
546 double L7 = getDistance(P( 1 ),P( 4 ));
547 double L8 = getDistance(P( 2 ),P( 5 ));
548 double L9 = getDistance(P( 3 ),P( 6 ));
549 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
550 aVal = Max(aVal,Max(Max(L7,L8),L9));
553 else if( len == 8 ) { // hexas
554 double L1 = getDistance(P( 1 ),P( 2 ));
555 double L2 = getDistance(P( 2 ),P( 3 ));
556 double L3 = getDistance(P( 3 ),P( 4 ));
557 double L4 = getDistance(P( 4 ),P( 1 ));
558 double L5 = getDistance(P( 5 ),P( 6 ));
559 double L6 = getDistance(P( 6 ),P( 7 ));
560 double L7 = getDistance(P( 7 ),P( 8 ));
561 double L8 = getDistance(P( 8 ),P( 5 ));
562 double L9 = getDistance(P( 1 ),P( 5 ));
563 double L10= getDistance(P( 2 ),P( 6 ));
564 double L11= getDistance(P( 3 ),P( 7 ));
565 double L12= getDistance(P( 4 ),P( 8 ));
566 double D1 = getDistance(P( 1 ),P( 7 ));
567 double D2 = getDistance(P( 2 ),P( 8 ));
568 double D3 = getDistance(P( 3 ),P( 5 ));
569 double D4 = getDistance(P( 4 ),P( 6 ));
570 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
571 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
572 aVal = Max(aVal,Max(L11,L12));
573 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
576 else if( len == 10 ) { // quadratic tetras
577 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
578 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
579 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
580 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
581 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
582 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
583 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
586 else if( len == 13 ) { // quadratic pyramids
587 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
588 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
589 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
590 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
591 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
592 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
593 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
594 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
595 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
596 aVal = Max(aVal,Max(L7,L8));
599 else if( len == 15 ) { // quadratic pentas
600 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
601 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
602 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
603 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
604 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
605 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
606 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
607 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
608 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
609 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
610 aVal = Max(aVal,Max(Max(L7,L8),L9));
613 else if( len == 20 ) { // quadratic hexas
614 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
615 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
616 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
617 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
618 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
619 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
620 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
621 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
622 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
623 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
624 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
625 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
626 double D1 = getDistance(P( 1 ),P( 7 ));
627 double D2 = getDistance(P( 2 ),P( 8 ));
628 double D3 = getDistance(P( 3 ),P( 5 ));
629 double D4 = getDistance(P( 4 ),P( 6 ));
630 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
631 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
632 aVal = Max(aVal,Max(L11,L12));
633 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
636 else if( len > 1 && aElem->IsPoly() ) { // polys
637 // get the maximum distance between all pairs of nodes
638 for( int i = 1; i <= len; i++ ) {
639 for( int j = 1; j <= len; j++ ) {
640 if( j > i ) { // optimization of the loop
641 double D = getDistance( P(i), P(j) );
642 aVal = Max( aVal, D );
649 if( myPrecision >= 0 )
651 double prec = pow( 10., (double)myPrecision );
652 aVal = floor( aVal * prec + 0.5 ) / prec;
659 double MaxElementLength3D::GetBadRate( double Value, int /*nbNodes*/ ) const
664 SMDSAbs_ElementType MaxElementLength3D::GetType() const
666 return SMDSAbs_Volume;
672 Description : Functor for calculation of minimum angle
675 double MinimumAngle::GetValue( const TSequenceOfXYZ& P )
682 aMin = getAngle(P( P.size() ), P( 1 ), P( 2 ));
683 aMin = Min(aMin,getAngle(P( P.size()-1 ), P( P.size() ), P( 1 )));
685 for (int i=2; i<P.size();i++){
686 double A0 = getAngle( P( i-1 ), P( i ), P( i+1 ) );
690 return aMin * 180.0 / PI;
693 double MinimumAngle::GetBadRate( double Value, int nbNodes ) const
695 //const double aBestAngle = PI / nbNodes;
696 const double aBestAngle = 180.0 - ( 360.0 / double(nbNodes) );
697 return ( fabs( aBestAngle - Value ));
700 SMDSAbs_ElementType MinimumAngle::GetType() const
708 Description : Functor for calculating aspect ratio
710 double AspectRatio::GetValue( const TSequenceOfXYZ& P )
712 // According to "Mesh quality control" by Nadir Bouhamau referring to
713 // Pascal Jean Frey and Paul-Louis George. Maillages, applications aux elements finis.
714 // Hermes Science publications, Paris 1999 ISBN 2-7462-0024-4
717 int nbNodes = P.size();
722 // Compute aspect ratio
724 if ( nbNodes == 3 ) {
725 // Compute lengths of the sides
726 std::vector< double > aLen (nbNodes);
727 for ( int i = 0; i < nbNodes - 1; i++ )
728 aLen[ i ] = getDistance( P( i + 1 ), P( i + 2 ) );
729 aLen[ nbNodes - 1 ] = getDistance( P( 1 ), P( nbNodes ) );
730 // Q = alfa * h * p / S, where
732 // alfa = sqrt( 3 ) / 6
733 // h - length of the longest edge
734 // p - half perimeter
735 // S - triangle surface
736 const double alfa = sqrt( 3. ) / 6.;
737 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
738 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
739 double anArea = getArea( P( 1 ), P( 2 ), P( 3 ) );
740 if ( anArea <= Precision::Confusion() )
742 return alfa * maxLen * half_perimeter / anArea;
744 else if ( nbNodes == 6 ) { // quadratic triangles
745 // Compute lengths of the sides
746 std::vector< double > aLen (3);
747 aLen[0] = getDistance( P(1), P(3) );
748 aLen[1] = getDistance( P(3), P(5) );
749 aLen[2] = getDistance( P(5), P(1) );
750 // Q = alfa * h * p / S, where
752 // alfa = sqrt( 3 ) / 6
753 // h - length of the longest edge
754 // p - half perimeter
755 // S - triangle surface
756 const double alfa = sqrt( 3. ) / 6.;
757 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
758 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
759 double anArea = getArea( P(1), P(3), P(5) );
760 if ( anArea <= Precision::Confusion() )
762 return alfa * maxLen * half_perimeter / anArea;
764 else if( nbNodes == 4 ) { // quadrangle
765 // Compute lengths of the sides
766 std::vector< double > aLen (4);
767 aLen[0] = getDistance( P(1), P(2) );
768 aLen[1] = getDistance( P(2), P(3) );
769 aLen[2] = getDistance( P(3), P(4) );
770 aLen[3] = getDistance( P(4), P(1) );
771 // Compute lengths of the diagonals
772 std::vector< double > aDia (2);
773 aDia[0] = getDistance( P(1), P(3) );
774 aDia[1] = getDistance( P(2), P(4) );
775 // Compute areas of all triangles which can be built
776 // taking three nodes of the quadrangle
777 std::vector< double > anArea (4);
778 anArea[0] = getArea( P(1), P(2), P(3) );
779 anArea[1] = getArea( P(1), P(2), P(4) );
780 anArea[2] = getArea( P(1), P(3), P(4) );
781 anArea[3] = getArea( P(2), P(3), P(4) );
782 // Q = alpha * L * C1 / C2, where
784 // alpha = sqrt( 1/32 )
785 // L = max( L1, L2, L3, L4, D1, D2 )
786 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
787 // C2 = min( S1, S2, S3, S4 )
788 // Li - lengths of the edges
789 // Di - lengths of the diagonals
790 // Si - areas of the triangles
791 const double alpha = sqrt( 1 / 32. );
792 double L = Max( aLen[ 0 ],
796 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
797 double C1 = sqrt( ( aLen[0] * aLen[0] +
800 aLen[3] * aLen[3] ) / 4. );
801 double C2 = Min( anArea[ 0 ],
803 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
804 if ( C2 <= Precision::Confusion() )
806 return alpha * L * C1 / C2;
808 else if( nbNodes == 8 ){ // nbNodes==8 - quadratic quadrangle
809 // Compute lengths of the sides
810 std::vector< double > aLen (4);
811 aLen[0] = getDistance( P(1), P(3) );
812 aLen[1] = getDistance( P(3), P(5) );
813 aLen[2] = getDistance( P(5), P(7) );
814 aLen[3] = getDistance( P(7), P(1) );
815 // Compute lengths of the diagonals
816 std::vector< double > aDia (2);
817 aDia[0] = getDistance( P(1), P(5) );
818 aDia[1] = getDistance( P(3), P(7) );
819 // Compute areas of all triangles which can be built
820 // taking three nodes of the quadrangle
821 std::vector< double > anArea (4);
822 anArea[0] = getArea( P(1), P(3), P(5) );
823 anArea[1] = getArea( P(1), P(3), P(7) );
824 anArea[2] = getArea( P(1), P(5), P(7) );
825 anArea[3] = getArea( P(3), P(5), P(7) );
826 // Q = alpha * L * C1 / C2, where
828 // alpha = sqrt( 1/32 )
829 // L = max( L1, L2, L3, L4, D1, D2 )
830 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
831 // C2 = min( S1, S2, S3, S4 )
832 // Li - lengths of the edges
833 // Di - lengths of the diagonals
834 // Si - areas of the triangles
835 const double alpha = sqrt( 1 / 32. );
836 double L = Max( aLen[ 0 ],
840 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
841 double C1 = sqrt( ( aLen[0] * aLen[0] +
844 aLen[3] * aLen[3] ) / 4. );
845 double C2 = Min( anArea[ 0 ],
847 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
848 if ( C2 <= Precision::Confusion() )
850 return alpha * L * C1 / C2;
855 double AspectRatio::GetBadRate( double Value, int /*nbNodes*/ ) const
857 // the aspect ratio is in the range [1.0,infinity]
860 return Value / 1000.;
863 SMDSAbs_ElementType AspectRatio::GetType() const
870 Class : AspectRatio3D
871 Description : Functor for calculating aspect ratio
875 inline double getHalfPerimeter(double theTria[3]){
876 return (theTria[0] + theTria[1] + theTria[2])/2.0;
879 inline double getArea(double theHalfPerim, double theTria[3]){
880 return sqrt(theHalfPerim*
881 (theHalfPerim-theTria[0])*
882 (theHalfPerim-theTria[1])*
883 (theHalfPerim-theTria[2]));
886 inline double getVolume(double theLen[6]){
887 double a2 = theLen[0]*theLen[0];
888 double b2 = theLen[1]*theLen[1];
889 double c2 = theLen[2]*theLen[2];
890 double d2 = theLen[3]*theLen[3];
891 double e2 = theLen[4]*theLen[4];
892 double f2 = theLen[5]*theLen[5];
893 double P = 4.0*a2*b2*d2;
894 double Q = a2*(b2+d2-e2)-b2*(a2+d2-f2)-d2*(a2+b2-c2);
895 double R = (b2+d2-e2)*(a2+d2-f2)*(a2+d2-f2);
896 return sqrt(P-Q+R)/12.0;
899 inline double getVolume2(double theLen[6]){
900 double a2 = theLen[0]*theLen[0];
901 double b2 = theLen[1]*theLen[1];
902 double c2 = theLen[2]*theLen[2];
903 double d2 = theLen[3]*theLen[3];
904 double e2 = theLen[4]*theLen[4];
905 double f2 = theLen[5]*theLen[5];
907 double P = a2*e2*(b2+c2+d2+f2-a2-e2);
908 double Q = b2*f2*(a2+c2+d2+e2-b2-f2);
909 double R = c2*d2*(a2+b2+e2+f2-c2-d2);
910 double S = a2*b2*d2+b2*c2*e2+a2*c2*f2+d2*e2*f2;
912 return sqrt(P+Q+R-S)/12.0;
915 inline double getVolume(const TSequenceOfXYZ& P){
916 gp_Vec aVec1( P( 2 ) - P( 1 ) );
917 gp_Vec aVec2( P( 3 ) - P( 1 ) );
918 gp_Vec aVec3( P( 4 ) - P( 1 ) );
919 gp_Vec anAreaVec( aVec1 ^ aVec2 );
920 return fabs(aVec3 * anAreaVec) / 6.0;
923 inline double getMaxHeight(double theLen[6])
925 double aHeight = std::max(theLen[0],theLen[1]);
926 aHeight = std::max(aHeight,theLen[2]);
927 aHeight = std::max(aHeight,theLen[3]);
928 aHeight = std::max(aHeight,theLen[4]);
929 aHeight = std::max(aHeight,theLen[5]);
935 double AspectRatio3D::GetValue( long theId )
938 myCurrElement = myMesh->FindElement( theId );
939 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_TETRA )
941 // Action from CoTech | ACTION 31.3:
942 // EURIWARE BO: Homogenize the formulas used to calculate the Controls in SMESH to fit with
943 // those of ParaView. The library used by ParaView for those calculations can be reused in SMESH.
944 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
945 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
946 aVal = Round( vtkMeshQuality::TetAspectRatio( avtkCell ));
951 if ( GetPoints( myCurrElement, P ))
952 aVal = Round( GetValue( P ));
957 double AspectRatio3D::GetValue( const TSequenceOfXYZ& P )
959 double aQuality = 0.0;
960 if(myCurrElement->IsPoly()) return aQuality;
962 int nbNodes = P.size();
964 if(myCurrElement->IsQuadratic()) {
965 if(nbNodes==10) nbNodes=4; // quadratic tetrahedron
966 else if(nbNodes==13) nbNodes=5; // quadratic pyramid
967 else if(nbNodes==15) nbNodes=6; // quadratic pentahedron
968 else if(nbNodes==20) nbNodes=8; // quadratic hexahedron
969 else return aQuality;
975 getDistance(P( 1 ),P( 2 )), // a
976 getDistance(P( 2 ),P( 3 )), // b
977 getDistance(P( 3 ),P( 1 )), // c
978 getDistance(P( 2 ),P( 4 )), // d
979 getDistance(P( 3 ),P( 4 )), // e
980 getDistance(P( 1 ),P( 4 )) // f
982 double aTria[4][3] = {
983 {aLen[0],aLen[1],aLen[2]}, // abc
984 {aLen[0],aLen[3],aLen[5]}, // adf
985 {aLen[1],aLen[3],aLen[4]}, // bde
986 {aLen[2],aLen[4],aLen[5]} // cef
988 double aSumArea = 0.0;
989 double aHalfPerimeter = getHalfPerimeter(aTria[0]);
990 double anArea = getArea(aHalfPerimeter,aTria[0]);
992 aHalfPerimeter = getHalfPerimeter(aTria[1]);
993 anArea = getArea(aHalfPerimeter,aTria[1]);
995 aHalfPerimeter = getHalfPerimeter(aTria[2]);
996 anArea = getArea(aHalfPerimeter,aTria[2]);
998 aHalfPerimeter = getHalfPerimeter(aTria[3]);
999 anArea = getArea(aHalfPerimeter,aTria[3]);
1001 double aVolume = getVolume(P);
1002 //double aVolume = getVolume(aLen);
1003 double aHeight = getMaxHeight(aLen);
1004 static double aCoeff = sqrt(2.0)/12.0;
1005 if ( aVolume > DBL_MIN )
1006 aQuality = aCoeff*aHeight*aSumArea/aVolume;
1011 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 3 ),P( 5 )};
1012 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1015 gp_XYZ aXYZ[4] = {P( 1 ),P( 3 ),P( 4 ),P( 5 )};
1016 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1019 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 5 )};
1020 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1023 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 4 ),P( 5 )};
1024 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1030 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 6 )};
1031 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1034 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 3 )};
1035 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1038 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 6 )};
1039 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1042 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1043 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1046 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 6 )};
1047 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1050 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 3 )};
1051 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1057 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),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( 4 )};
1062 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1065 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 7 )};
1066 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1069 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 8 )};
1070 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1073 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 3 )};
1074 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1077 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 4 )};
1078 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1081 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 7 )};
1082 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1085 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 8 )};
1086 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1089 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 3 )};
1090 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1093 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 4 )};
1094 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1097 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 7 )};
1098 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1101 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 8 )};
1102 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1105 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 1 )};
1106 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1109 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 2 )};
1110 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1113 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 5 )};
1114 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1117 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 6 )};
1118 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1121 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 1 )};
1122 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1125 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 2 )};
1126 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1129 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 5 )};
1130 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1133 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 6 )};
1134 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1137 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 1 )};
1138 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1141 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1142 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1145 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 5 )};
1146 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1149 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 6 )};
1150 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1153 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1154 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1157 gp_XYZ aXYZ[4] = {P( 4 ),P( 5 ),P( 8 ),P( 2 )};
1158 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1161 gp_XYZ aXYZ[4] = {P( 1 ),P( 4 ),P( 5 ),P( 3 )};
1162 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1165 gp_XYZ aXYZ[4] = {P( 3 ),P( 6 ),P( 7 ),P( 1 )};
1166 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1169 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 6 ),P( 4 )};
1170 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1173 gp_XYZ aXYZ[4] = {P( 5 ),P( 6 ),P( 8 ),P( 3 )};
1174 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1177 gp_XYZ aXYZ[4] = {P( 7 ),P( 8 ),P( 6 ),P( 1 )};
1178 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1181 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 7 )};
1182 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1185 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 2 ),P( 5 )};
1186 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1191 if ( nbNodes > 4 ) {
1192 // avaluate aspect ratio of quadranle faces
1193 AspectRatio aspect2D;
1194 SMDS_VolumeTool::VolumeType type = SMDS_VolumeTool::GetType( nbNodes );
1195 int nbFaces = SMDS_VolumeTool::NbFaces( type );
1196 TSequenceOfXYZ points(4);
1197 for ( int i = 0; i < nbFaces; ++i ) { // loop on faces of a volume
1198 if ( SMDS_VolumeTool::NbFaceNodes( type, i ) != 4 )
1200 const int* pInd = SMDS_VolumeTool::GetFaceNodesIndices( type, i, true );
1201 for ( int p = 0; p < 4; ++p ) // loop on nodes of a quadranle face
1202 points( p + 1 ) = P( pInd[ p ] + 1 );
1203 aQuality = std::max( aQuality, aspect2D.GetValue( points ));
1209 double AspectRatio3D::GetBadRate( double Value, int /*nbNodes*/ ) const
1211 // the aspect ratio is in the range [1.0,infinity]
1214 return Value / 1000.;
1217 SMDSAbs_ElementType AspectRatio3D::GetType() const
1219 return SMDSAbs_Volume;
1225 Description : Functor for calculating warping
1227 double Warping::GetValue( const TSequenceOfXYZ& P )
1229 if ( P.size() != 4 )
1232 gp_XYZ G = ( P( 1 ) + P( 2 ) + P( 3 ) + P( 4 ) ) / 4.;
1234 double A1 = ComputeA( P( 1 ), P( 2 ), P( 3 ), G );
1235 double A2 = ComputeA( P( 2 ), P( 3 ), P( 4 ), G );
1236 double A3 = ComputeA( P( 3 ), P( 4 ), P( 1 ), G );
1237 double A4 = ComputeA( P( 4 ), P( 1 ), P( 2 ), G );
1239 return Max( Max( A1, A2 ), Max( A3, A4 ) );
1242 double Warping::ComputeA( const gp_XYZ& thePnt1,
1243 const gp_XYZ& thePnt2,
1244 const gp_XYZ& thePnt3,
1245 const gp_XYZ& theG ) const
1247 double aLen1 = gp_Pnt( thePnt1 ).Distance( gp_Pnt( thePnt2 ) );
1248 double aLen2 = gp_Pnt( thePnt2 ).Distance( gp_Pnt( thePnt3 ) );
1249 double L = Min( aLen1, aLen2 ) * 0.5;
1250 if ( L < Precision::Confusion())
1253 gp_XYZ GI = ( thePnt2 + thePnt1 ) / 2. - theG;
1254 gp_XYZ GJ = ( thePnt3 + thePnt2 ) / 2. - theG;
1255 gp_XYZ N = GI.Crossed( GJ );
1257 if ( N.Modulus() < gp::Resolution() )
1262 double H = ( thePnt2 - theG ).Dot( N );
1263 return asin( fabs( H / L ) ) * 180. / PI;
1266 double Warping::GetBadRate( double Value, int /*nbNodes*/ ) const
1268 // the warp is in the range [0.0,PI/2]
1269 // 0.0 = good (no warp)
1270 // PI/2 = bad (face pliee)
1274 SMDSAbs_ElementType Warping::GetType() const
1276 return SMDSAbs_Face;
1282 Description : Functor for calculating taper
1284 double Taper::GetValue( const TSequenceOfXYZ& P )
1286 if ( P.size() != 4 )
1290 double J1 = getArea( P( 4 ), P( 1 ), P( 2 ) ) / 2.;
1291 double J2 = getArea( P( 3 ), P( 1 ), P( 2 ) ) / 2.;
1292 double J3 = getArea( P( 2 ), P( 3 ), P( 4 ) ) / 2.;
1293 double J4 = getArea( P( 3 ), P( 4 ), P( 1 ) ) / 2.;
1295 double JA = 0.25 * ( J1 + J2 + J3 + J4 );
1296 if ( JA <= Precision::Confusion() )
1299 double T1 = fabs( ( J1 - JA ) / JA );
1300 double T2 = fabs( ( J2 - JA ) / JA );
1301 double T3 = fabs( ( J3 - JA ) / JA );
1302 double T4 = fabs( ( J4 - JA ) / JA );
1304 return Max( Max( T1, T2 ), Max( T3, T4 ) );
1307 double Taper::GetBadRate( double Value, int /*nbNodes*/ ) const
1309 // the taper is in the range [0.0,1.0]
1310 // 0.0 = good (no taper)
1311 // 1.0 = bad (les cotes opposes sont allignes)
1315 SMDSAbs_ElementType Taper::GetType() const
1317 return SMDSAbs_Face;
1323 Description : Functor for calculating skew in degrees
1325 static inline double skewAngle( const gp_XYZ& p1, const gp_XYZ& p2, const gp_XYZ& p3 )
1327 gp_XYZ p12 = ( p2 + p1 ) / 2.;
1328 gp_XYZ p23 = ( p3 + p2 ) / 2.;
1329 gp_XYZ p31 = ( p3 + p1 ) / 2.;
1331 gp_Vec v1( p31 - p2 ), v2( p12 - p23 );
1333 return v1.Magnitude() < gp::Resolution() || v2.Magnitude() < gp::Resolution() ? 0. : v1.Angle( v2 );
1336 double Skew::GetValue( const TSequenceOfXYZ& P )
1338 if ( P.size() != 3 && P.size() != 4 )
1342 static double PI2 = PI / 2.;
1343 if ( P.size() == 3 )
1345 double A0 = fabs( PI2 - skewAngle( P( 3 ), P( 1 ), P( 2 ) ) );
1346 double A1 = fabs( PI2 - skewAngle( P( 1 ), P( 2 ), P( 3 ) ) );
1347 double A2 = fabs( PI2 - skewAngle( P( 2 ), P( 3 ), P( 1 ) ) );
1349 return Max( A0, Max( A1, A2 ) ) * 180. / PI;
1353 gp_XYZ p12 = ( P( 1 ) + P( 2 ) ) / 2.;
1354 gp_XYZ p23 = ( P( 2 ) + P( 3 ) ) / 2.;
1355 gp_XYZ p34 = ( P( 3 ) + P( 4 ) ) / 2.;
1356 gp_XYZ p41 = ( P( 4 ) + P( 1 ) ) / 2.;
1358 gp_Vec v1( p34 - p12 ), v2( p23 - p41 );
1359 double A = v1.Magnitude() <= gp::Resolution() || v2.Magnitude() <= gp::Resolution()
1360 ? 0. : fabs( PI2 - v1.Angle( v2 ) );
1363 if ( A < Precision::Angular() )
1366 return A * 180. / PI;
1370 double Skew::GetBadRate( double Value, int /*nbNodes*/ ) const
1372 // the skew is in the range [0.0,PI/2].
1378 SMDSAbs_ElementType Skew::GetType() const
1380 return SMDSAbs_Face;
1386 Description : Functor for calculating area
1388 double Area::GetValue( const TSequenceOfXYZ& P )
1391 if ( P.size() > 2 ) {
1392 gp_Vec aVec1( P(2) - P(1) );
1393 gp_Vec aVec2( P(3) - P(1) );
1394 gp_Vec SumVec = aVec1 ^ aVec2;
1395 for (int i=4; i<=P.size(); i++) {
1396 gp_Vec aVec1( P(i-1) - P(1) );
1397 gp_Vec aVec2( P(i) - P(1) );
1398 gp_Vec tmp = aVec1 ^ aVec2;
1401 val = SumVec.Magnitude() * 0.5;
1406 double Area::GetBadRate( double Value, int /*nbNodes*/ ) const
1408 // meaningless as it is not a quality control functor
1412 SMDSAbs_ElementType Area::GetType() const
1414 return SMDSAbs_Face;
1420 Description : Functor for calculating length off edge
1422 double Length::GetValue( const TSequenceOfXYZ& P )
1424 switch ( P.size() ) {
1425 case 2: return getDistance( P( 1 ), P( 2 ) );
1426 case 3: return getDistance( P( 1 ), P( 2 ) ) + getDistance( P( 2 ), P( 3 ) );
1431 double Length::GetBadRate( double Value, int /*nbNodes*/ ) const
1433 // meaningless as it is not quality control functor
1437 SMDSAbs_ElementType Length::GetType() const
1439 return SMDSAbs_Edge;
1444 Description : Functor for calculating length of edge
1447 double Length2D::GetValue( long theElementId)
1451 //cout<<"Length2D::GetValue"<<endl;
1452 if (GetPoints(theElementId,P)){
1453 //for(int jj=1; jj<=P.size(); jj++)
1454 // cout<<"jj="<<jj<<" P("<<P(jj).X()<<","<<P(jj).Y()<<","<<P(jj).Z()<<")"<<endl;
1456 double aVal;// = GetValue( P );
1457 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
1458 SMDSAbs_ElementType aType = aElem->GetType();
1467 aVal = getDistance( P( 1 ), P( 2 ) );
1470 else if (len == 3){ // quadratic edge
1471 aVal = getDistance(P( 1 ),P( 3 )) + getDistance(P( 3 ),P( 2 ));
1475 if (len == 3){ // triangles
1476 double L1 = getDistance(P( 1 ),P( 2 ));
1477 double L2 = getDistance(P( 2 ),P( 3 ));
1478 double L3 = getDistance(P( 3 ),P( 1 ));
1479 aVal = Max(L1,Max(L2,L3));
1482 else if (len == 4){ // quadrangles
1483 double L1 = getDistance(P( 1 ),P( 2 ));
1484 double L2 = getDistance(P( 2 ),P( 3 ));
1485 double L3 = getDistance(P( 3 ),P( 4 ));
1486 double L4 = getDistance(P( 4 ),P( 1 ));
1487 aVal = Max(Max(L1,L2),Max(L3,L4));
1490 if (len == 6){ // quadratic triangles
1491 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1492 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1493 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
1494 aVal = Max(L1,Max(L2,L3));
1495 //cout<<"L1="<<L1<<" L2="<<L2<<"L3="<<L3<<" aVal="<<aVal<<endl;
1498 else if (len == 8){ // quadratic quadrangles
1499 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1500 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1501 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
1502 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
1503 aVal = Max(Max(L1,L2),Max(L3,L4));
1506 case SMDSAbs_Volume:
1507 if (len == 4){ // tetraidrs
1508 double L1 = getDistance(P( 1 ),P( 2 ));
1509 double L2 = getDistance(P( 2 ),P( 3 ));
1510 double L3 = getDistance(P( 3 ),P( 1 ));
1511 double L4 = getDistance(P( 1 ),P( 4 ));
1512 double L5 = getDistance(P( 2 ),P( 4 ));
1513 double L6 = getDistance(P( 3 ),P( 4 ));
1514 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1517 else if (len == 5){ // piramids
1518 double L1 = getDistance(P( 1 ),P( 2 ));
1519 double L2 = getDistance(P( 2 ),P( 3 ));
1520 double L3 = getDistance(P( 3 ),P( 4 ));
1521 double L4 = getDistance(P( 4 ),P( 1 ));
1522 double L5 = getDistance(P( 1 ),P( 5 ));
1523 double L6 = getDistance(P( 2 ),P( 5 ));
1524 double L7 = getDistance(P( 3 ),P( 5 ));
1525 double L8 = getDistance(P( 4 ),P( 5 ));
1527 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1528 aVal = Max(aVal,Max(L7,L8));
1531 else if (len == 6){ // pentaidres
1532 double L1 = getDistance(P( 1 ),P( 2 ));
1533 double L2 = getDistance(P( 2 ),P( 3 ));
1534 double L3 = getDistance(P( 3 ),P( 1 ));
1535 double L4 = getDistance(P( 4 ),P( 5 ));
1536 double L5 = getDistance(P( 5 ),P( 6 ));
1537 double L6 = getDistance(P( 6 ),P( 4 ));
1538 double L7 = getDistance(P( 1 ),P( 4 ));
1539 double L8 = getDistance(P( 2 ),P( 5 ));
1540 double L9 = getDistance(P( 3 ),P( 6 ));
1542 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1543 aVal = Max(aVal,Max(Max(L7,L8),L9));
1546 else if (len == 8){ // hexaider
1547 double L1 = getDistance(P( 1 ),P( 2 ));
1548 double L2 = getDistance(P( 2 ),P( 3 ));
1549 double L3 = getDistance(P( 3 ),P( 4 ));
1550 double L4 = getDistance(P( 4 ),P( 1 ));
1551 double L5 = getDistance(P( 5 ),P( 6 ));
1552 double L6 = getDistance(P( 6 ),P( 7 ));
1553 double L7 = getDistance(P( 7 ),P( 8 ));
1554 double L8 = getDistance(P( 8 ),P( 5 ));
1555 double L9 = getDistance(P( 1 ),P( 5 ));
1556 double L10= getDistance(P( 2 ),P( 6 ));
1557 double L11= getDistance(P( 3 ),P( 7 ));
1558 double L12= getDistance(P( 4 ),P( 8 ));
1560 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1561 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1562 aVal = Max(aVal,Max(L11,L12));
1567 if (len == 10){ // quadratic tetraidrs
1568 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
1569 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
1570 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
1571 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1572 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
1573 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
1574 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1577 else if (len == 13){ // quadratic piramids
1578 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
1579 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
1580 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1581 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1582 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1583 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
1584 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
1585 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
1586 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1587 aVal = Max(aVal,Max(L7,L8));
1590 else if (len == 15){ // quadratic pentaidres
1591 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
1592 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
1593 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1594 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1595 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
1596 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
1597 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
1598 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
1599 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
1600 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1601 aVal = Max(aVal,Max(Max(L7,L8),L9));
1604 else if (len == 20){ // quadratic hexaider
1605 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
1606 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
1607 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
1608 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
1609 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
1610 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
1611 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
1612 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
1613 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
1614 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
1615 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
1616 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
1617 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1618 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1619 aVal = Max(aVal,Max(L11,L12));
1631 if ( myPrecision >= 0 )
1633 double prec = pow( 10., (double)( myPrecision ) );
1634 aVal = floor( aVal * prec + 0.5 ) / prec;
1643 double Length2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1645 // meaningless as it is not quality control functor
1649 SMDSAbs_ElementType Length2D::GetType() const
1651 return SMDSAbs_Face;
1654 Length2D::Value::Value(double theLength,long thePntId1, long thePntId2):
1657 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1658 if(thePntId1 > thePntId2){
1659 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1663 bool Length2D::Value::operator<(const Length2D::Value& x) const{
1664 if(myPntId[0] < x.myPntId[0]) return true;
1665 if(myPntId[0] == x.myPntId[0])
1666 if(myPntId[1] < x.myPntId[1]) return true;
1670 void Length2D::GetValues(TValues& theValues){
1672 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1673 for(; anIter->more(); ){
1674 const SMDS_MeshFace* anElem = anIter->next();
1676 if(anElem->IsQuadratic()) {
1677 const SMDS_VtkFace* F =
1678 dynamic_cast<const SMDS_VtkFace*>(anElem);
1679 // use special nodes iterator
1680 SMDS_ElemIteratorPtr anIter = F->interlacedNodesElemIterator();
1685 const SMDS_MeshElement* aNode;
1687 aNode = anIter->next();
1688 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1689 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1690 aNodeId[0] = aNodeId[1] = aNode->GetID();
1693 for(; anIter->more(); ){
1694 const SMDS_MeshNode* N1 = static_cast<const SMDS_MeshNode*> (anIter->next());
1695 P[2] = gp_Pnt(N1->X(),N1->Y(),N1->Z());
1696 aNodeId[2] = N1->GetID();
1697 aLength = P[1].Distance(P[2]);
1698 if(!anIter->more()) break;
1699 const SMDS_MeshNode* N2 = static_cast<const SMDS_MeshNode*> (anIter->next());
1700 P[3] = gp_Pnt(N2->X(),N2->Y(),N2->Z());
1701 aNodeId[3] = N2->GetID();
1702 aLength += P[2].Distance(P[3]);
1703 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1704 Value aValue2(aLength,aNodeId[2],aNodeId[3]);
1706 aNodeId[1] = aNodeId[3];
1707 theValues.insert(aValue1);
1708 theValues.insert(aValue2);
1710 aLength += P[2].Distance(P[0]);
1711 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1712 Value aValue2(aLength,aNodeId[2],aNodeId[0]);
1713 theValues.insert(aValue1);
1714 theValues.insert(aValue2);
1717 SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
1722 const SMDS_MeshElement* aNode;
1723 if(aNodesIter->more()){
1724 aNode = aNodesIter->next();
1725 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1726 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1727 aNodeId[0] = aNodeId[1] = aNode->GetID();
1730 for(; aNodesIter->more(); ){
1731 aNode = aNodesIter->next();
1732 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1733 long anId = aNode->GetID();
1735 P[2] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1737 aLength = P[1].Distance(P[2]);
1739 Value aValue(aLength,aNodeId[1],anId);
1742 theValues.insert(aValue);
1745 aLength = P[0].Distance(P[1]);
1747 Value aValue(aLength,aNodeId[0],aNodeId[1]);
1748 theValues.insert(aValue);
1754 Class : MultiConnection
1755 Description : Functor for calculating number of faces conneted to the edge
1757 double MultiConnection::GetValue( const TSequenceOfXYZ& P )
1761 double MultiConnection::GetValue( long theId )
1763 return getNbMultiConnection( myMesh, theId );
1766 double MultiConnection::GetBadRate( double Value, int /*nbNodes*/ ) const
1768 // meaningless as it is not quality control functor
1772 SMDSAbs_ElementType MultiConnection::GetType() const
1774 return SMDSAbs_Edge;
1778 Class : MultiConnection2D
1779 Description : Functor for calculating number of faces conneted to the edge
1781 double MultiConnection2D::GetValue( const TSequenceOfXYZ& P )
1786 double MultiConnection2D::GetValue( long theElementId )
1790 const SMDS_MeshElement* aFaceElem = myMesh->FindElement(theElementId);
1791 SMDSAbs_ElementType aType = aFaceElem->GetType();
1796 int i = 0, len = aFaceElem->NbNodes();
1797 SMDS_ElemIteratorPtr anIter = aFaceElem->nodesIterator();
1800 const SMDS_MeshNode *aNode, *aNode0;
1801 TColStd_MapOfInteger aMap, aMapPrev;
1803 for (i = 0; i <= len; i++) {
1808 if (anIter->more()) {
1809 aNode = (SMDS_MeshNode*)anIter->next();
1817 if (i == 0) aNode0 = aNode;
1819 SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
1820 while (anElemIter->more()) {
1821 const SMDS_MeshElement* anElem = anElemIter->next();
1822 if (anElem != 0 && anElem->GetType() == SMDSAbs_Face) {
1823 int anId = anElem->GetID();
1826 if (aMapPrev.Contains(anId)) {
1831 aResult = Max(aResult, aNb);
1842 double MultiConnection2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1844 // meaningless as it is not quality control functor
1848 SMDSAbs_ElementType MultiConnection2D::GetType() const
1850 return SMDSAbs_Face;
1853 MultiConnection2D::Value::Value(long thePntId1, long thePntId2)
1855 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1856 if(thePntId1 > thePntId2){
1857 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1861 bool MultiConnection2D::Value::operator<(const MultiConnection2D::Value& x) const{
1862 if(myPntId[0] < x.myPntId[0]) return true;
1863 if(myPntId[0] == x.myPntId[0])
1864 if(myPntId[1] < x.myPntId[1]) return true;
1868 void MultiConnection2D::GetValues(MValues& theValues){
1869 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1870 for(; anIter->more(); ){
1871 const SMDS_MeshFace* anElem = anIter->next();
1872 SMDS_ElemIteratorPtr aNodesIter;
1873 if ( anElem->IsQuadratic() )
1874 aNodesIter = dynamic_cast<const SMDS_VtkFace*>
1875 (anElem)->interlacedNodesElemIterator();
1877 aNodesIter = anElem->nodesIterator();
1880 //int aNbConnects=0;
1881 const SMDS_MeshNode* aNode0;
1882 const SMDS_MeshNode* aNode1;
1883 const SMDS_MeshNode* aNode2;
1884 if(aNodesIter->more()){
1885 aNode0 = (SMDS_MeshNode*) aNodesIter->next();
1887 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode1;
1888 aNodeId[0] = aNodeId[1] = aNodes->GetID();
1890 for(; aNodesIter->more(); ) {
1891 aNode2 = (SMDS_MeshNode*) aNodesIter->next();
1892 long anId = aNode2->GetID();
1895 Value aValue(aNodeId[1],aNodeId[2]);
1896 MValues::iterator aItr = theValues.find(aValue);
1897 if (aItr != theValues.end()){
1902 theValues[aValue] = 1;
1905 //cout << "NodeIds: "<<aNodeId[1]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1906 aNodeId[1] = aNodeId[2];
1909 Value aValue(aNodeId[0],aNodeId[2]);
1910 MValues::iterator aItr = theValues.find(aValue);
1911 if (aItr != theValues.end()) {
1916 theValues[aValue] = 1;
1919 //cout << "NodeIds: "<<aNodeId[0]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1929 Class : BadOrientedVolume
1930 Description : Predicate bad oriented volumes
1933 BadOrientedVolume::BadOrientedVolume()
1938 void BadOrientedVolume::SetMesh( const SMDS_Mesh* theMesh )
1943 bool BadOrientedVolume::IsSatisfy( long theId )
1948 SMDS_VolumeTool vTool( myMesh->FindElement( theId ));
1949 return !vTool.IsForward();
1952 SMDSAbs_ElementType BadOrientedVolume::GetType() const
1954 return SMDSAbs_Volume;
1958 Class : BareBorderVolume
1961 bool BareBorderVolume::IsSatisfy(long theElementId )
1963 SMDS_VolumeTool myTool;
1964 if ( myTool.Set( myMesh->FindElement(theElementId)))
1966 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
1967 if ( myTool.IsFreeFace( iF ))
1969 const SMDS_MeshNode** n = myTool.GetFaceNodes(iF);
1970 vector< const SMDS_MeshNode*> nodes( n, n+myTool.NbFaceNodes(iF));
1971 if ( !myMesh->FindElement( nodes, SMDSAbs_Face, /*Nomedium=*/false))
1979 Class : BareBorderFace
1982 bool BareBorderFace::IsSatisfy(long theElementId )
1985 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
1987 if ( face->GetType() == SMDSAbs_Face )
1989 int nbN = face->NbCornerNodes();
1990 for ( int i = 0; i < nbN && !ok; ++i )
1992 // check if a link is shared by another face
1993 const SMDS_MeshNode* n1 = face->GetNode( i );
1994 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
1995 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
1996 bool isShared = false;
1997 while ( !isShared && fIt->more() )
1999 const SMDS_MeshElement* f = fIt->next();
2000 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2004 myLinkNodes.resize( 2 + face->IsQuadratic());
2005 myLinkNodes[0] = n1;
2006 myLinkNodes[1] = n2;
2007 if ( face->IsQuadratic() )
2008 myLinkNodes[2] = face->GetNode( i+nbN );
2009 ok = !myMesh->FindElement( myLinkNodes, SMDSAbs_Edge, /*noMedium=*/false);
2018 Class : OverConstrainedVolume
2021 bool OverConstrainedVolume::IsSatisfy(long theElementId )
2023 // An element is over-constrained if it has N-1 free borders where
2024 // N is the number of edges/faces for a 2D/3D element.
2025 SMDS_VolumeTool myTool;
2026 if ( myTool.Set( myMesh->FindElement(theElementId)))
2028 int nbSharedFaces = 0;
2029 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2030 if ( !myTool.IsFreeFace( iF ) && ++nbSharedFaces > 1 )
2032 return ( nbSharedFaces == 1 );
2038 Class : OverConstrainedFace
2041 bool OverConstrainedFace::IsSatisfy(long theElementId )
2043 // An element is over-constrained if it has N-1 free borders where
2044 // N is the number of edges/faces for a 2D/3D element.
2045 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2046 if ( face->GetType() == SMDSAbs_Face )
2048 int nbSharedBorders = 0;
2049 int nbN = face->NbCornerNodes();
2050 for ( int i = 0; i < nbN; ++i )
2052 // check if a link is shared by another face
2053 const SMDS_MeshNode* n1 = face->GetNode( i );
2054 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2055 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2056 bool isShared = false;
2057 while ( !isShared && fIt->more() )
2059 const SMDS_MeshElement* f = fIt->next();
2060 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2062 if ( isShared && ++nbSharedBorders > 1 )
2065 return ( nbSharedBorders == 1 );
2072 Description : Predicate for free borders
2075 FreeBorders::FreeBorders()
2080 void FreeBorders::SetMesh( const SMDS_Mesh* theMesh )
2085 bool FreeBorders::IsSatisfy( long theId )
2087 return getNbMultiConnection( myMesh, theId ) == 1;
2090 SMDSAbs_ElementType FreeBorders::GetType() const
2092 return SMDSAbs_Edge;
2098 Description : Predicate for free Edges
2100 FreeEdges::FreeEdges()
2105 void FreeEdges::SetMesh( const SMDS_Mesh* theMesh )
2110 bool FreeEdges::IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId )
2112 TColStd_MapOfInteger aMap;
2113 for ( int i = 0; i < 2; i++ )
2115 SMDS_ElemIteratorPtr anElemIter = theNodes[ i ]->GetInverseElementIterator();
2116 while( anElemIter->more() )
2118 const SMDS_MeshElement* anElem = anElemIter->next();
2119 if ( anElem != 0 && anElem->GetType() == SMDSAbs_Face )
2121 int anId = anElem->GetID();
2125 else if ( aMap.Contains( anId ) && anId != theFaceId )
2133 bool FreeEdges::IsSatisfy( long theId )
2138 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2139 if ( aFace == 0 || aFace->GetType() != SMDSAbs_Face || aFace->NbNodes() < 3 )
2142 SMDS_ElemIteratorPtr anIter;
2143 if ( aFace->IsQuadratic() ) {
2144 anIter = dynamic_cast<const SMDS_VtkFace*>
2145 (aFace)->interlacedNodesElemIterator();
2148 anIter = aFace->nodesIterator();
2153 int i = 0, nbNodes = aFace->NbNodes();
2154 std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
2155 while( anIter->more() )
2157 const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
2160 aNodes[ i++ ] = aNode;
2162 aNodes[ nbNodes ] = aNodes[ 0 ];
2164 for ( i = 0; i < nbNodes; i++ )
2165 if ( IsFreeEdge( &aNodes[ i ], theId ) )
2171 SMDSAbs_ElementType FreeEdges::GetType() const
2173 return SMDSAbs_Face;
2176 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2179 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
2180 if(thePntId1 > thePntId2){
2181 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
2185 bool FreeEdges::Border::operator<(const FreeEdges::Border& x) const{
2186 if(myPntId[0] < x.myPntId[0]) return true;
2187 if(myPntId[0] == x.myPntId[0])
2188 if(myPntId[1] < x.myPntId[1]) return true;
2192 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2193 FreeEdges::TBorders& theRegistry,
2194 FreeEdges::TBorders& theContainer)
2196 if(theRegistry.find(theBorder) == theRegistry.end()){
2197 theRegistry.insert(theBorder);
2198 theContainer.insert(theBorder);
2200 theContainer.erase(theBorder);
2204 void FreeEdges::GetBoreders(TBorders& theBorders)
2207 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2208 for(; anIter->more(); ){
2209 const SMDS_MeshFace* anElem = anIter->next();
2210 long anElemId = anElem->GetID();
2211 SMDS_ElemIteratorPtr aNodesIter;
2212 if ( anElem->IsQuadratic() )
2213 aNodesIter = static_cast<const SMDS_VtkFace*>(anElem)->
2214 interlacedNodesElemIterator();
2216 aNodesIter = anElem->nodesIterator();
2218 const SMDS_MeshElement* aNode;
2219 if(aNodesIter->more()){
2220 aNode = aNodesIter->next();
2221 aNodeId[0] = aNodeId[1] = aNode->GetID();
2223 for(; aNodesIter->more(); ){
2224 aNode = aNodesIter->next();
2225 long anId = aNode->GetID();
2226 Border aBorder(anElemId,aNodeId[1],anId);
2228 //std::cout<<aBorder.myPntId[0]<<"; "<<aBorder.myPntId[1]<<"; "<<aBorder.myElemId<<endl;
2229 UpdateBorders(aBorder,aRegistry,theBorders);
2231 Border aBorder(anElemId,aNodeId[0],aNodeId[1]);
2232 //std::cout<<aBorder.myPntId[0]<<"; "<<aBorder.myPntId[1]<<"; "<<aBorder.myElemId<<endl;
2233 UpdateBorders(aBorder,aRegistry,theBorders);
2235 //std::cout<<"theBorders.size() = "<<theBorders.size()<<endl;
2241 Description : Predicate for free nodes
2244 FreeNodes::FreeNodes()
2249 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2254 bool FreeNodes::IsSatisfy( long theNodeId )
2256 const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2260 return (aNode->NbInverseElements() < 1);
2263 SMDSAbs_ElementType FreeNodes::GetType() const
2265 return SMDSAbs_Node;
2271 Description : Predicate for free faces
2274 FreeFaces::FreeFaces()
2279 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2284 bool FreeFaces::IsSatisfy( long theId )
2286 if (!myMesh) return false;
2287 // check that faces nodes refers to less than two common volumes
2288 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2289 if ( !aFace || aFace->GetType() != SMDSAbs_Face )
2292 int nbNode = aFace->NbNodes();
2294 // collect volumes check that number of volumss with count equal nbNode not less than 2
2295 typedef map< SMDS_MeshElement*, int > TMapOfVolume; // map of volume counters
2296 typedef map< SMDS_MeshElement*, int >::iterator TItrMapOfVolume; // iterator
2297 TMapOfVolume mapOfVol;
2299 SMDS_ElemIteratorPtr nodeItr = aFace->nodesIterator();
2300 while ( nodeItr->more() ) {
2301 const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>(nodeItr->next());
2302 if ( !aNode ) continue;
2303 SMDS_ElemIteratorPtr volItr = aNode->GetInverseElementIterator(SMDSAbs_Volume);
2304 while ( volItr->more() ) {
2305 SMDS_MeshElement* aVol = (SMDS_MeshElement*)volItr->next();
2306 TItrMapOfVolume itr = mapOfVol.insert(make_pair(aVol, 0)).first;
2311 TItrMapOfVolume volItr = mapOfVol.begin();
2312 TItrMapOfVolume volEnd = mapOfVol.end();
2313 for ( ; volItr != volEnd; ++volItr )
2314 if ( (*volItr).second >= nbNode )
2316 // face is not free if number of volumes constructed on thier nodes more than one
2320 SMDSAbs_ElementType FreeFaces::GetType() const
2322 return SMDSAbs_Face;
2326 Class : LinearOrQuadratic
2327 Description : Predicate to verify whether a mesh element is linear
2330 LinearOrQuadratic::LinearOrQuadratic()
2335 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2340 bool LinearOrQuadratic::IsSatisfy( long theId )
2342 if (!myMesh) return false;
2343 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2344 if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2346 return (!anElem->IsQuadratic());
2349 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2354 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2361 Description : Functor for check color of group to whic mesh element belongs to
2364 GroupColor::GroupColor()
2368 bool GroupColor::IsSatisfy( long theId )
2370 return (myIDs.find( theId ) != myIDs.end());
2373 void GroupColor::SetType( SMDSAbs_ElementType theType )
2378 SMDSAbs_ElementType GroupColor::GetType() const
2383 static bool isEqual( const Quantity_Color& theColor1,
2384 const Quantity_Color& theColor2 )
2386 // tolerance to compare colors
2387 const double tol = 5*1e-3;
2388 return ( fabs( theColor1.Red() - theColor2.Red() ) < tol &&
2389 fabs( theColor1.Green() - theColor2.Green() ) < tol &&
2390 fabs( theColor1.Blue() - theColor2.Blue() ) < tol );
2394 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2398 const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2402 int nbGrp = aMesh->GetNbGroups();
2406 // iterates on groups and find necessary elements ids
2407 const std::set<SMESHDS_GroupBase*>& aGroups = aMesh->GetGroups();
2408 set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
2409 for (; GrIt != aGroups.end(); GrIt++) {
2410 SMESHDS_GroupBase* aGrp = (*GrIt);
2413 // check type and color of group
2414 if ( !isEqual( myColor, aGrp->GetColor() ) )
2416 if ( myType != SMDSAbs_All && myType != (SMDSAbs_ElementType)aGrp->GetType() )
2419 SMDSAbs_ElementType aGrpElType = (SMDSAbs_ElementType)aGrp->GetType();
2420 if ( myType == aGrpElType || (myType == SMDSAbs_All && aGrpElType != SMDSAbs_Node) ) {
2421 // add elements IDS into control
2422 int aSize = aGrp->Extent();
2423 for (int i = 0; i < aSize; i++)
2424 myIDs.insert( aGrp->GetID(i+1) );
2429 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
2431 TCollection_AsciiString aStr = theStr;
2432 aStr.RemoveAll( ' ' );
2433 aStr.RemoveAll( '\t' );
2434 for ( int aPos = aStr.Search( ";;" ); aPos != -1; aPos = aStr.Search( ";;" ) )
2435 aStr.Remove( aPos, 2 );
2436 Standard_Real clr[3];
2437 clr[0] = clr[1] = clr[2] = 0.;
2438 for ( int i = 0; i < 3; i++ ) {
2439 TCollection_AsciiString tmpStr = aStr.Token( ";", i+1 );
2440 if ( !tmpStr.IsEmpty() && tmpStr.IsRealValue() )
2441 clr[i] = tmpStr.RealValue();
2443 myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
2446 //=======================================================================
2447 // name : GetRangeStr
2448 // Purpose : Get range as a string.
2449 // Example: "1,2,3,50-60,63,67,70-"
2450 //=======================================================================
2451 void GroupColor::GetColorStr( TCollection_AsciiString& theResStr ) const
2454 theResStr += TCollection_AsciiString( myColor.Red() );
2455 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2456 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2460 Class : ElemGeomType
2461 Description : Predicate to check element geometry type
2464 ElemGeomType::ElemGeomType()
2467 myType = SMDSAbs_All;
2468 myGeomType = SMDSGeom_TRIANGLE;
2471 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2476 bool ElemGeomType::IsSatisfy( long theId )
2478 if (!myMesh) return false;
2479 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2482 const SMDSAbs_ElementType anElemType = anElem->GetType();
2483 if ( myType != SMDSAbs_All && anElemType != myType )
2485 const int aNbNode = anElem->NbNodes();
2487 switch( anElemType )
2490 isOk = (myGeomType == SMDSGeom_POINT);
2494 isOk = (myGeomType == SMDSGeom_EDGE);
2498 if ( myGeomType == SMDSGeom_TRIANGLE )
2499 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 6 : aNbNode == 3));
2500 else if ( myGeomType == SMDSGeom_QUADRANGLE )
2501 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 8 : aNbNode == 4));
2502 else if ( myGeomType == SMDSGeom_POLYGON )
2503 isOk = anElem->IsPoly();
2506 case SMDSAbs_Volume:
2507 if ( myGeomType == SMDSGeom_TETRA )
2508 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 10 : aNbNode == 4));
2509 else if ( myGeomType == SMDSGeom_PYRAMID )
2510 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 13 : aNbNode == 5));
2511 else if ( myGeomType == SMDSGeom_PENTA )
2512 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 15 : aNbNode == 6));
2513 else if ( myGeomType == SMDSGeom_HEXA )
2514 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 20 : aNbNode == 8));
2515 else if ( myGeomType == SMDSGeom_POLYHEDRA )
2516 isOk = anElem->IsPoly();
2523 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2528 SMDSAbs_ElementType ElemGeomType::GetType() const
2533 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2535 myGeomType = theType;
2538 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2543 //================================================================================
2545 * \brief Class CoplanarFaces
2547 //================================================================================
2549 CoplanarFaces::CoplanarFaces()
2550 : myMesh(0), myFaceID(0), myToler(0)
2553 bool CoplanarFaces::IsSatisfy( long theElementId )
2555 if ( myCoplanarIDs.empty() )
2557 // Build a set of coplanar face ids
2559 if ( !myMesh || !myFaceID || !myToler )
2562 const SMDS_MeshElement* face = myMesh->FindElement( myFaceID );
2563 if ( !face || face->GetType() != SMDSAbs_Face )
2567 gp_Vec myNorm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2571 const double radianTol = myToler * PI180;
2572 typedef SMDS_StdIterator< const SMDS_MeshElement*, SMDS_ElemIteratorPtr > TFaceIt;
2573 std::set<const SMDS_MeshElement*> checkedFaces, checkedNodes;
2574 std::list<const SMDS_MeshElement*> faceQueue( 1, face );
2575 while ( !faceQueue.empty() )
2577 face = faceQueue.front();
2578 if ( checkedFaces.insert( face ).second )
2580 gp_Vec norm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2581 if (!normOK || myNorm.Angle( norm ) <= radianTol)
2583 myCoplanarIDs.insert( face->GetID() );
2584 std::set<const SMDS_MeshElement*> neighborFaces;
2585 for ( int i = 0; i < face->NbCornerNodes(); ++i )
2587 const SMDS_MeshNode* n = face->GetNode( i );
2588 if ( checkedNodes.insert( n ).second )
2589 neighborFaces.insert( TFaceIt( n->GetInverseElementIterator(SMDSAbs_Face)),
2592 faceQueue.insert( faceQueue.end(), neighborFaces.begin(), neighborFaces.end() );
2595 faceQueue.pop_front();
2598 return myCoplanarIDs.count( theElementId );
2603 *Description : Predicate for Range of Ids.
2604 * Range may be specified with two ways.
2605 * 1. Using AddToRange method
2606 * 2. With SetRangeStr method. Parameter of this method is a string
2607 * like as "1,2,3,50-60,63,67,70-"
2610 //=======================================================================
2611 // name : RangeOfIds
2612 // Purpose : Constructor
2613 //=======================================================================
2614 RangeOfIds::RangeOfIds()
2617 myType = SMDSAbs_All;
2620 //=======================================================================
2622 // Purpose : Set mesh
2623 //=======================================================================
2624 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
2629 //=======================================================================
2630 // name : AddToRange
2631 // Purpose : Add ID to the range
2632 //=======================================================================
2633 bool RangeOfIds::AddToRange( long theEntityId )
2635 myIds.Add( theEntityId );
2639 //=======================================================================
2640 // name : GetRangeStr
2641 // Purpose : Get range as a string.
2642 // Example: "1,2,3,50-60,63,67,70-"
2643 //=======================================================================
2644 void RangeOfIds::GetRangeStr( TCollection_AsciiString& theResStr )
2648 TColStd_SequenceOfInteger anIntSeq;
2649 TColStd_SequenceOfAsciiString aStrSeq;
2651 TColStd_MapIteratorOfMapOfInteger anIter( myIds );
2652 for ( ; anIter.More(); anIter.Next() )
2654 int anId = anIter.Key();
2655 TCollection_AsciiString aStr( anId );
2656 anIntSeq.Append( anId );
2657 aStrSeq.Append( aStr );
2660 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2662 int aMinId = myMin( i );
2663 int aMaxId = myMax( i );
2665 TCollection_AsciiString aStr;
2666 if ( aMinId != IntegerFirst() )
2671 if ( aMaxId != IntegerLast() )
2674 // find position of the string in result sequence and insert string in it
2675 if ( anIntSeq.Length() == 0 )
2677 anIntSeq.Append( aMinId );
2678 aStrSeq.Append( aStr );
2682 if ( aMinId < anIntSeq.First() )
2684 anIntSeq.Prepend( aMinId );
2685 aStrSeq.Prepend( aStr );
2687 else if ( aMinId > anIntSeq.Last() )
2689 anIntSeq.Append( aMinId );
2690 aStrSeq.Append( aStr );
2693 for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
2694 if ( aMinId < anIntSeq( j ) )
2696 anIntSeq.InsertBefore( j, aMinId );
2697 aStrSeq.InsertBefore( j, aStr );
2703 if ( aStrSeq.Length() == 0 )
2706 theResStr = aStrSeq( 1 );
2707 for ( int j = 2, k = aStrSeq.Length(); j <= k; j++ )
2710 theResStr += aStrSeq( j );
2714 //=======================================================================
2715 // name : SetRangeStr
2716 // Purpose : Define range with string
2717 // Example of entry string: "1,2,3,50-60,63,67,70-"
2718 //=======================================================================
2719 bool RangeOfIds::SetRangeStr( const TCollection_AsciiString& theStr )
2725 TCollection_AsciiString aStr = theStr;
2726 aStr.RemoveAll( ' ' );
2727 aStr.RemoveAll( '\t' );
2729 for ( int aPos = aStr.Search( ",," ); aPos != -1; aPos = aStr.Search( ",," ) )
2730 aStr.Remove( aPos, 2 );
2732 TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
2734 while ( tmpStr != "" )
2736 tmpStr = aStr.Token( ",", i++ );
2737 int aPos = tmpStr.Search( '-' );
2741 if ( tmpStr.IsIntegerValue() )
2742 myIds.Add( tmpStr.IntegerValue() );
2748 TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
2749 TCollection_AsciiString aMinStr = tmpStr;
2751 while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
2752 while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
2754 if ( (!aMinStr.IsEmpty() && !aMinStr.IsIntegerValue()) ||
2755 (!aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue()) )
2758 myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
2759 myMax.Append( aMaxStr.IsEmpty() ? IntegerLast() : aMaxStr.IntegerValue() );
2766 //=======================================================================
2768 // Purpose : Get type of supported entities
2769 //=======================================================================
2770 SMDSAbs_ElementType RangeOfIds::GetType() const
2775 //=======================================================================
2777 // Purpose : Set type of supported entities
2778 //=======================================================================
2779 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
2784 //=======================================================================
2786 // Purpose : Verify whether entity satisfies to this rpedicate
2787 //=======================================================================
2788 bool RangeOfIds::IsSatisfy( long theId )
2793 if ( myType == SMDSAbs_Node )
2795 if ( myMesh->FindNode( theId ) == 0 )
2800 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2801 if ( anElem == 0 || (myType != anElem->GetType() && myType != SMDSAbs_All ))
2805 if ( myIds.Contains( theId ) )
2808 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2809 if ( theId >= myMin( i ) && theId <= myMax( i ) )
2817 Description : Base class for comparators
2819 Comparator::Comparator():
2823 Comparator::~Comparator()
2826 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
2829 myFunctor->SetMesh( theMesh );
2832 void Comparator::SetMargin( double theValue )
2834 myMargin = theValue;
2837 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
2839 myFunctor = theFunct;
2842 SMDSAbs_ElementType Comparator::GetType() const
2844 return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
2847 double Comparator::GetMargin()
2855 Description : Comparator "<"
2857 bool LessThan::IsSatisfy( long theId )
2859 return myFunctor && myFunctor->GetValue( theId ) < myMargin;
2865 Description : Comparator ">"
2867 bool MoreThan::IsSatisfy( long theId )
2869 return myFunctor && myFunctor->GetValue( theId ) > myMargin;
2875 Description : Comparator "="
2878 myToler(Precision::Confusion())
2881 bool EqualTo::IsSatisfy( long theId )
2883 return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
2886 void EqualTo::SetTolerance( double theToler )
2891 double EqualTo::GetTolerance()
2898 Description : Logical NOT predicate
2900 LogicalNOT::LogicalNOT()
2903 LogicalNOT::~LogicalNOT()
2906 bool LogicalNOT::IsSatisfy( long theId )
2908 return myPredicate && !myPredicate->IsSatisfy( theId );
2911 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
2914 myPredicate->SetMesh( theMesh );
2917 void LogicalNOT::SetPredicate( PredicatePtr thePred )
2919 myPredicate = thePred;
2922 SMDSAbs_ElementType LogicalNOT::GetType() const
2924 return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
2929 Class : LogicalBinary
2930 Description : Base class for binary logical predicate
2932 LogicalBinary::LogicalBinary()
2935 LogicalBinary::~LogicalBinary()
2938 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
2941 myPredicate1->SetMesh( theMesh );
2944 myPredicate2->SetMesh( theMesh );
2947 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
2949 myPredicate1 = thePredicate;
2952 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
2954 myPredicate2 = thePredicate;
2957 SMDSAbs_ElementType LogicalBinary::GetType() const
2959 if ( !myPredicate1 || !myPredicate2 )
2962 SMDSAbs_ElementType aType1 = myPredicate1->GetType();
2963 SMDSAbs_ElementType aType2 = myPredicate2->GetType();
2965 return aType1 == aType2 ? aType1 : SMDSAbs_All;
2971 Description : Logical AND
2973 bool LogicalAND::IsSatisfy( long theId )
2978 myPredicate1->IsSatisfy( theId ) &&
2979 myPredicate2->IsSatisfy( theId );
2985 Description : Logical OR
2987 bool LogicalOR::IsSatisfy( long theId )
2992 (myPredicate1->IsSatisfy( theId ) ||
2993 myPredicate2->IsSatisfy( theId ));
3007 void Filter::SetPredicate( PredicatePtr thePredicate )
3009 myPredicate = thePredicate;
3012 template<class TElement, class TIterator, class TPredicate>
3013 inline void FillSequence(const TIterator& theIterator,
3014 TPredicate& thePredicate,
3015 Filter::TIdSequence& theSequence)
3017 if ( theIterator ) {
3018 while( theIterator->more() ) {
3019 TElement anElem = theIterator->next();
3020 long anId = anElem->GetID();
3021 if ( thePredicate->IsSatisfy( anId ) )
3022 theSequence.push_back( anId );
3029 GetElementsId( const SMDS_Mesh* theMesh,
3030 PredicatePtr thePredicate,
3031 TIdSequence& theSequence )
3033 theSequence.clear();
3035 if ( !theMesh || !thePredicate )
3038 thePredicate->SetMesh( theMesh );
3040 SMDSAbs_ElementType aType = thePredicate->GetType();
3043 FillSequence<const SMDS_MeshNode*>(theMesh->nodesIterator(),thePredicate,theSequence);
3046 FillSequence<const SMDS_MeshElement*>(theMesh->edgesIterator(),thePredicate,theSequence);
3049 FillSequence<const SMDS_MeshElement*>(theMesh->facesIterator(),thePredicate,theSequence);
3051 case SMDSAbs_Volume:
3052 FillSequence<const SMDS_MeshElement*>(theMesh->volumesIterator(),thePredicate,theSequence);
3055 FillSequence<const SMDS_MeshElement*>(theMesh->edgesIterator(),thePredicate,theSequence);
3056 FillSequence<const SMDS_MeshElement*>(theMesh->facesIterator(),thePredicate,theSequence);
3057 FillSequence<const SMDS_MeshElement*>(theMesh->volumesIterator(),thePredicate,theSequence);
3063 Filter::GetElementsId( const SMDS_Mesh* theMesh,
3064 Filter::TIdSequence& theSequence )
3066 GetElementsId(theMesh,myPredicate,theSequence);
3073 typedef std::set<SMDS_MeshFace*> TMapOfFacePtr;
3079 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
3080 SMDS_MeshNode* theNode2 )
3086 ManifoldPart::Link::~Link()
3092 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
3094 if ( myNode1 == theLink.myNode1 &&
3095 myNode2 == theLink.myNode2 )
3097 else if ( myNode1 == theLink.myNode2 &&
3098 myNode2 == theLink.myNode1 )
3104 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
3106 if(myNode1 < x.myNode1) return true;
3107 if(myNode1 == x.myNode1)
3108 if(myNode2 < x.myNode2) return true;
3112 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
3113 const ManifoldPart::Link& theLink2 )
3115 return theLink1.IsEqual( theLink2 );
3118 ManifoldPart::ManifoldPart()
3121 myAngToler = Precision::Angular();
3122 myIsOnlyManifold = true;
3125 ManifoldPart::~ManifoldPart()
3130 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
3136 SMDSAbs_ElementType ManifoldPart::GetType() const
3137 { return SMDSAbs_Face; }
3139 bool ManifoldPart::IsSatisfy( long theElementId )
3141 return myMapIds.Contains( theElementId );
3144 void ManifoldPart::SetAngleTolerance( const double theAngToler )
3145 { myAngToler = theAngToler; }
3147 double ManifoldPart::GetAngleTolerance() const
3148 { return myAngToler; }
3150 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
3151 { myIsOnlyManifold = theIsOnly; }
3153 void ManifoldPart::SetStartElem( const long theStartId )
3154 { myStartElemId = theStartId; }
3156 bool ManifoldPart::process()
3159 myMapBadGeomIds.Clear();
3161 myAllFacePtr.clear();
3162 myAllFacePtrIntDMap.clear();
3166 // collect all faces into own map
3167 SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3168 for (; anFaceItr->more(); )
3170 SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3171 myAllFacePtr.push_back( aFacePtr );
3172 myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3175 SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3179 // the map of non manifold links and bad geometry
3180 TMapOfLink aMapOfNonManifold;
3181 TColStd_MapOfInteger aMapOfTreated;
3183 // begin cycle on faces from start index and run on vector till the end
3184 // and from begin to start index to cover whole vector
3185 const int aStartIndx = myAllFacePtrIntDMap[aStartFace];
3186 bool isStartTreat = false;
3187 for ( int fi = aStartIndx; !isStartTreat || fi != aStartIndx ; fi++ )
3189 if ( fi == aStartIndx )
3190 isStartTreat = true;
3191 // as result next time when fi will be equal to aStartIndx
3193 SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3194 if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3197 aMapOfTreated.Add( aFacePtr->GetID() );
3198 TColStd_MapOfInteger aResFaces;
3199 if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3200 aMapOfNonManifold, aResFaces ) )
3202 TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3203 for ( ; anItr.More(); anItr.Next() )
3205 int aFaceId = anItr.Key();
3206 aMapOfTreated.Add( aFaceId );
3207 myMapIds.Add( aFaceId );
3210 if ( fi == ( myAllFacePtr.size() - 1 ) )
3212 } // end run on vector of faces
3213 return !myMapIds.IsEmpty();
3216 static void getLinks( const SMDS_MeshFace* theFace,
3217 ManifoldPart::TVectorOfLink& theLinks )
3219 int aNbNode = theFace->NbNodes();
3220 SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3222 SMDS_MeshNode* aNode = 0;
3223 for ( ; aNodeItr->more() && i <= aNbNode; )
3226 SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3230 SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3232 ManifoldPart::Link aLink( aN1, aN2 );
3233 theLinks.push_back( aLink );
3237 bool ManifoldPart::findConnected
3238 ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3239 SMDS_MeshFace* theStartFace,
3240 ManifoldPart::TMapOfLink& theNonManifold,
3241 TColStd_MapOfInteger& theResFaces )
3243 theResFaces.Clear();
3244 if ( !theAllFacePtrInt.size() )
3247 if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3249 myMapBadGeomIds.Add( theStartFace->GetID() );
3253 ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3254 ManifoldPart::TVectorOfLink aSeqOfBoundary;
3255 theResFaces.Add( theStartFace->GetID() );
3256 ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3258 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3259 aDMapLinkFace, theNonManifold, theStartFace );
3261 bool isDone = false;
3262 while ( !isDone && aMapOfBoundary.size() != 0 )
3264 bool isToReset = false;
3265 ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3266 for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3268 ManifoldPart::Link aLink = *pLink;
3269 if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3271 // each link could be treated only once
3272 aMapToSkip.insert( aLink );
3274 ManifoldPart::TVectorOfFacePtr aFaces;
3276 if ( myIsOnlyManifold &&
3277 (theNonManifold.find( aLink ) != theNonManifold.end()) )
3281 getFacesByLink( aLink, aFaces );
3282 // filter the element to keep only indicated elements
3283 ManifoldPart::TVectorOfFacePtr aFiltered;
3284 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3285 for ( ; pFace != aFaces.end(); ++pFace )
3287 SMDS_MeshFace* aFace = *pFace;
3288 if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3289 aFiltered.push_back( aFace );
3292 if ( aFaces.size() < 2 ) // no neihgbour faces
3294 else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3296 theNonManifold.insert( aLink );
3301 // compare normal with normals of neighbor element
3302 SMDS_MeshFace* aPrevFace = aDMapLinkFace[ aLink ];
3303 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3304 for ( ; pFace != aFaces.end(); ++pFace )
3306 SMDS_MeshFace* aNextFace = *pFace;
3307 if ( aPrevFace == aNextFace )
3309 int anNextFaceID = aNextFace->GetID();
3310 if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3311 // should not be with non manifold restriction. probably bad topology
3313 // check if face was treated and skipped
3314 if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3315 !isInPlane( aPrevFace, aNextFace ) )
3317 // add new element to connected and extend the boundaries.
3318 theResFaces.Add( anNextFaceID );
3319 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3320 aDMapLinkFace, theNonManifold, aNextFace );
3324 isDone = !isToReset;
3327 return !theResFaces.IsEmpty();
3330 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3331 const SMDS_MeshFace* theFace2 )
3333 gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3334 gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3335 if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3337 myMapBadGeomIds.Add( theFace2->GetID() );
3340 if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
3346 void ManifoldPart::expandBoundary
3347 ( ManifoldPart::TMapOfLink& theMapOfBoundary,
3348 ManifoldPart::TVectorOfLink& theSeqOfBoundary,
3349 ManifoldPart::TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
3350 ManifoldPart::TMapOfLink& theNonManifold,
3351 SMDS_MeshFace* theNextFace ) const
3353 ManifoldPart::TVectorOfLink aLinks;
3354 getLinks( theNextFace, aLinks );
3355 int aNbLink = (int)aLinks.size();
3356 for ( int i = 0; i < aNbLink; i++ )
3358 ManifoldPart::Link aLink = aLinks[ i ];
3359 if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3361 if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3363 if ( myIsOnlyManifold )
3365 // remove from boundary
3366 theMapOfBoundary.erase( aLink );
3367 ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3368 for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3370 ManifoldPart::Link aBoundLink = *pLink;
3371 if ( aBoundLink.IsEqual( aLink ) )
3373 theSeqOfBoundary.erase( pLink );
3381 theMapOfBoundary.insert( aLink );
3382 theSeqOfBoundary.push_back( aLink );
3383 theDMapLinkFacePtr[ aLink ] = theNextFace;
3388 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3389 ManifoldPart::TVectorOfFacePtr& theFaces ) const
3391 std::set<SMDS_MeshCell *> aSetOfFaces;
3392 // take all faces that shared first node
3393 SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3394 for ( ; anItr->more(); )
3396 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3399 aSetOfFaces.insert( aFace );
3401 // take all faces that shared second node
3402 anItr = theLink.myNode2->facesIterator();
3403 // find the common part of two sets
3404 for ( ; anItr->more(); )
3406 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3407 if ( aSetOfFaces.count( aFace ) )
3408 theFaces.push_back( aFace );
3417 ElementsOnSurface::ElementsOnSurface()
3421 myType = SMDSAbs_All;
3423 myToler = Precision::Confusion();
3424 myUseBoundaries = false;
3427 ElementsOnSurface::~ElementsOnSurface()
3432 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
3434 if ( myMesh == theMesh )
3440 bool ElementsOnSurface::IsSatisfy( long theElementId )
3442 return myIds.Contains( theElementId );
3445 SMDSAbs_ElementType ElementsOnSurface::GetType() const
3448 void ElementsOnSurface::SetTolerance( const double theToler )
3450 if ( myToler != theToler )
3455 double ElementsOnSurface::GetTolerance() const
3458 void ElementsOnSurface::SetUseBoundaries( bool theUse )
3460 if ( myUseBoundaries != theUse ) {
3461 myUseBoundaries = theUse;
3462 SetSurface( mySurf, myType );
3466 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
3467 const SMDSAbs_ElementType theType )
3472 if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
3474 mySurf = TopoDS::Face( theShape );
3475 BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
3477 u1 = SA.FirstUParameter(),
3478 u2 = SA.LastUParameter(),
3479 v1 = SA.FirstVParameter(),
3480 v2 = SA.LastVParameter();
3481 Handle(Geom_Surface) surf = BRep_Tool::Surface( mySurf );
3482 myProjector.Init( surf, u1,u2, v1,v2 );
3486 void ElementsOnSurface::process()
3489 if ( mySurf.IsNull() )
3495 if ( myType == SMDSAbs_Face || myType == SMDSAbs_All )
3497 myIds.ReSize( myMesh->NbFaces() );
3498 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
3499 for(; anIter->more(); )
3500 process( anIter->next() );
3503 if ( myType == SMDSAbs_Edge || myType == SMDSAbs_All )
3505 myIds.ReSize( myIds.Extent() + myMesh->NbEdges() );
3506 SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
3507 for(; anIter->more(); )
3508 process( anIter->next() );
3511 if ( myType == SMDSAbs_Node )
3513 myIds.ReSize( myMesh->NbNodes() );
3514 SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
3515 for(; anIter->more(); )
3516 process( anIter->next() );
3520 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
3522 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3523 bool isSatisfy = true;
3524 for ( ; aNodeItr->more(); )
3526 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3527 if ( !isOnSurface( aNode ) )
3534 myIds.Add( theElemPtr->GetID() );
3537 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
3539 if ( mySurf.IsNull() )
3542 gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
3543 // double aToler2 = myToler * myToler;
3544 // if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
3546 // gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
3547 // if ( aPln.SquareDistance( aPnt ) > aToler2 )
3550 // else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
3552 // gp_Cylinder aCyl = Handle(Geom_CylindricalSurface)::DownCast(mySurf)->Cylinder();
3553 // double aRad = aCyl.Radius();
3554 // gp_Ax3 anAxis = aCyl.Position();
3555 // gp_XYZ aLoc = aCyl.Location().XYZ();
3556 // double aXDist = anAxis.XDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3557 // double aYDist = anAxis.YDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3558 // if ( fabs(aXDist*aXDist + aYDist*aYDist - aRad*aRad) > aToler2 )
3563 myProjector.Perform( aPnt );
3564 bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
3574 ElementsOnShape::ElementsOnShape()
3576 myType(SMDSAbs_All),
3577 myToler(Precision::Confusion()),
3578 myAllNodesFlag(false)
3580 myCurShapeType = TopAbs_SHAPE;
3583 ElementsOnShape::~ElementsOnShape()
3587 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
3589 if (myMesh != theMesh) {
3591 SetShape(myShape, myType);
3595 bool ElementsOnShape::IsSatisfy (long theElementId)
3597 return myIds.Contains(theElementId);
3600 SMDSAbs_ElementType ElementsOnShape::GetType() const
3605 void ElementsOnShape::SetTolerance (const double theToler)
3607 if (myToler != theToler) {
3609 SetShape(myShape, myType);
3613 double ElementsOnShape::GetTolerance() const
3618 void ElementsOnShape::SetAllNodes (bool theAllNodes)
3620 if (myAllNodesFlag != theAllNodes) {
3621 myAllNodesFlag = theAllNodes;
3622 SetShape(myShape, myType);
3626 void ElementsOnShape::SetShape (const TopoDS_Shape& theShape,
3627 const SMDSAbs_ElementType theType)
3633 if (myMesh == 0) return;
3638 myIds.ReSize(myMesh->NbEdges() + myMesh->NbFaces() + myMesh->NbVolumes());
3641 myIds.ReSize(myMesh->NbNodes());
3644 myIds.ReSize(myMesh->NbEdges());
3647 myIds.ReSize(myMesh->NbFaces());
3649 case SMDSAbs_Volume:
3650 myIds.ReSize(myMesh->NbVolumes());
3656 myShapesMap.Clear();
3660 void ElementsOnShape::addShape (const TopoDS_Shape& theShape)
3662 if (theShape.IsNull() || myMesh == 0)
3665 if (!myShapesMap.Add(theShape)) return;
3667 myCurShapeType = theShape.ShapeType();
3668 switch (myCurShapeType)
3670 case TopAbs_COMPOUND:
3671 case TopAbs_COMPSOLID:
3675 TopoDS_Iterator anIt (theShape, Standard_True, Standard_True);
3676 for (; anIt.More(); anIt.Next()) addShape(anIt.Value());
3681 myCurSC.Load(theShape);
3687 TopoDS_Face aFace = TopoDS::Face(theShape);
3688 BRepAdaptor_Surface SA (aFace, true);
3690 u1 = SA.FirstUParameter(),
3691 u2 = SA.LastUParameter(),
3692 v1 = SA.FirstVParameter(),
3693 v2 = SA.LastVParameter();
3694 Handle(Geom_Surface) surf = BRep_Tool::Surface(aFace);
3695 myCurProjFace.Init(surf, u1,u2, v1,v2);
3702 TopoDS_Edge anEdge = TopoDS::Edge(theShape);
3703 Standard_Real u1, u2;
3704 Handle(Geom_Curve) curve = BRep_Tool::Curve(anEdge, u1, u2);
3705 myCurProjEdge.Init(curve, u1, u2);
3711 TopoDS_Vertex aV = TopoDS::Vertex(theShape);
3712 myCurPnt = BRep_Tool::Pnt(aV);
3721 void ElementsOnShape::process()
3723 if (myShape.IsNull() || myMesh == 0)
3726 if (myType == SMDSAbs_Node)
3728 SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
3729 while (anIter->more())
3730 process(anIter->next());
3734 if (myType == SMDSAbs_Edge || myType == SMDSAbs_All)
3736 SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
3737 while (anIter->more())
3738 process(anIter->next());
3741 if (myType == SMDSAbs_Face || myType == SMDSAbs_All)
3743 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
3744 while (anIter->more()) {
3745 process(anIter->next());
3749 if (myType == SMDSAbs_Volume || myType == SMDSAbs_All)
3751 SMDS_VolumeIteratorPtr anIter = myMesh->volumesIterator();
3752 while (anIter->more())
3753 process(anIter->next());
3758 void ElementsOnShape::process (const SMDS_MeshElement* theElemPtr)
3760 if (myShape.IsNull())
3763 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3764 bool isSatisfy = myAllNodesFlag;
3766 gp_XYZ centerXYZ (0, 0, 0);
3768 while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
3770 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3771 gp_Pnt aPnt (aNode->X(), aNode->Y(), aNode->Z());
3772 centerXYZ += aPnt.XYZ();
3774 switch (myCurShapeType)
3778 myCurSC.Perform(aPnt, myToler);
3779 isSatisfy = (myCurSC.State() == TopAbs_IN || myCurSC.State() == TopAbs_ON);
3784 myCurProjFace.Perform(aPnt);
3785 isSatisfy = (myCurProjFace.IsDone() && myCurProjFace.LowerDistance() <= myToler);
3788 // check relatively the face
3789 Quantity_Parameter u, v;
3790 myCurProjFace.LowerDistanceParameters(u, v);
3791 gp_Pnt2d aProjPnt (u, v);
3792 BRepClass_FaceClassifier aClsf (myCurFace, aProjPnt, myToler);
3793 isSatisfy = (aClsf.State() == TopAbs_IN || aClsf.State() == TopAbs_ON);
3799 myCurProjEdge.Perform(aPnt);
3800 isSatisfy = (myCurProjEdge.NbPoints() > 0 && myCurProjEdge.LowerDistance() <= myToler);
3805 isSatisfy = (aPnt.Distance(myCurPnt) <= myToler);
3815 if (isSatisfy && myCurShapeType == TopAbs_SOLID) { // Check the center point for volumes MantisBug 0020168
3816 centerXYZ /= theElemPtr->NbNodes();
3817 gp_Pnt aCenterPnt (centerXYZ);
3818 myCurSC.Perform(aCenterPnt, myToler);
3819 if ( !(myCurSC.State() == TopAbs_IN || myCurSC.State() == TopAbs_ON))
3824 myIds.Add(theElemPtr->GetID());
3827 TSequenceOfXYZ::TSequenceOfXYZ()
3830 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n)
3833 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t)
3836 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray)
3839 template <class InputIterator>
3840 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd)
3843 TSequenceOfXYZ::~TSequenceOfXYZ()
3846 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
3848 myArray = theSequenceOfXYZ.myArray;
3852 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
3854 return myArray[n-1];
3857 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
3859 return myArray[n-1];
3862 void TSequenceOfXYZ::clear()
3867 void TSequenceOfXYZ::reserve(size_type n)
3872 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
3874 myArray.push_back(v);
3877 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
3879 return myArray.size();