1 // Copyright (C) 2007-2010 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"
75 inline gp_XYZ gpXYZ(const SMDS_MeshNode* aNode )
77 return gp_XYZ(aNode->X(), aNode->Y(), aNode->Z() );
80 inline double getAngle( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
82 gp_Vec v1( P1 - P2 ), v2( P3 - P2 );
84 return v1.Magnitude() < gp::Resolution() ||
85 v2.Magnitude() < gp::Resolution() ? 0 : v1.Angle( v2 );
88 inline double getArea( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
90 gp_Vec aVec1( P2 - P1 );
91 gp_Vec aVec2( P3 - P1 );
92 return ( aVec1 ^ aVec2 ).Magnitude() * 0.5;
95 inline double getArea( const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3 )
97 return getArea( P1.XYZ(), P2.XYZ(), P3.XYZ() );
102 inline double getDistance( const gp_XYZ& P1, const gp_XYZ& P2 )
104 double aDist = gp_Pnt( P1 ).Distance( gp_Pnt( P2 ) );
108 int getNbMultiConnection( const SMDS_Mesh* theMesh, const int theId )
113 const SMDS_MeshElement* anEdge = theMesh->FindElement( theId );
114 if ( anEdge == 0 || anEdge->GetType() != SMDSAbs_Edge/* || anEdge->NbNodes() != 2 */)
117 // for each pair of nodes in anEdge (there are 2 pairs in a quadratic edge)
118 // count elements containing both nodes of the pair.
119 // Note that there may be such cases for a quadratic edge (a horizontal line):
124 // +-----+------+ +-----+------+
127 // result sould be 2 in both cases
129 int aResult0 = 0, aResult1 = 0;
130 // last node, it is a medium one in a quadratic edge
131 const SMDS_MeshNode* aLastNode = anEdge->GetNode( anEdge->NbNodes() - 1 );
132 const SMDS_MeshNode* aNode0 = anEdge->GetNode( 0 );
133 const SMDS_MeshNode* aNode1 = anEdge->GetNode( 1 );
134 if ( aNode1 == aLastNode ) aNode1 = 0;
136 SMDS_ElemIteratorPtr anElemIter = aLastNode->GetInverseElementIterator();
137 while( anElemIter->more() ) {
138 const SMDS_MeshElement* anElem = anElemIter->next();
139 if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
140 SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
141 while ( anIter->more() ) {
142 if ( const SMDS_MeshElement* anElemNode = anIter->next() ) {
143 if ( anElemNode == aNode0 ) {
145 if ( !aNode1 ) break; // not a quadratic edge
147 else if ( anElemNode == aNode1 )
153 int aResult = std::max ( aResult0, aResult1 );
155 // TColStd_MapOfInteger aMap;
157 // SMDS_ElemIteratorPtr anIter = anEdge->nodesIterator();
158 // if ( anIter != 0 ) {
159 // while( anIter->more() ) {
160 // const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
163 // SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
164 // while( anElemIter->more() ) {
165 // const SMDS_MeshElement* anElem = anElemIter->next();
166 // if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
167 // int anId = anElem->GetID();
169 // if ( anIter->more() ) // i.e. first node
171 // else if ( aMap.Contains( anId ) )
181 gp_XYZ getNormale( const SMDS_MeshFace* theFace, bool* ok=0 )
183 int aNbNode = theFace->NbNodes();
185 gp_XYZ q1 = gpXYZ( theFace->GetNode(1)) - gpXYZ( theFace->GetNode(0));
186 gp_XYZ q2 = gpXYZ( theFace->GetNode(2)) - gpXYZ( theFace->GetNode(0));
189 gp_XYZ q3 = gpXYZ( theFace->GetNode(3)) - gpXYZ( theFace->GetNode(0));
192 double len = n.Modulus();
193 bool zeroLen = ( len <= numeric_limits<double>::min());
197 if (ok) *ok = !zeroLen;
205 using namespace SMESH::Controls;
212 Class : NumericalFunctor
213 Description : Base class for numerical functors
215 NumericalFunctor::NumericalFunctor():
221 void NumericalFunctor::SetMesh( const SMDS_Mesh* theMesh )
226 bool NumericalFunctor::GetPoints(const int theId,
227 TSequenceOfXYZ& theRes ) const
234 return GetPoints( myMesh->FindElement( theId ), theRes );
237 bool NumericalFunctor::GetPoints(const SMDS_MeshElement* anElem,
238 TSequenceOfXYZ& theRes )
245 theRes.reserve( anElem->NbNodes() );
247 // Get nodes of the element
248 SMDS_ElemIteratorPtr anIter;
250 if ( anElem->IsQuadratic() ) {
251 switch ( anElem->GetType() ) {
253 anIter = dynamic_cast<const SMDS_VtkEdge*>
254 (anElem)->interlacedNodesElemIterator();
257 anIter = dynamic_cast<const SMDS_VtkFace*>
258 (anElem)->interlacedNodesElemIterator();
261 anIter = anElem->nodesIterator();
266 anIter = anElem->nodesIterator();
270 while( anIter->more() ) {
271 if ( const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>( anIter->next() ))
272 theRes.push_back( gp_XYZ( aNode->X(), aNode->Y(), aNode->Z() ) );
279 long NumericalFunctor::GetPrecision() const
284 void NumericalFunctor::SetPrecision( const long thePrecision )
286 myPrecision = thePrecision;
289 double NumericalFunctor::GetValue( long theId )
291 myCurrElement = myMesh->FindElement( theId );
293 if ( GetPoints( theId, P ))
295 double aVal = GetValue( P );
296 if ( myPrecision >= 0 )
298 double prec = pow( 10., (double)( myPrecision ) );
299 aVal = floor( aVal * prec + 0.5 ) / prec;
307 //================================================================================
309 * \brief Return histogram of functor values
310 * \param nbIntervals - number of intervals
311 * \param nbEvents - number of mesh elements having values within i-th interval
312 * \param funValues - boundaries of intervals
313 * \param elements - elements to check vulue of; empty list means "of all"
314 * \param minmax - boundaries of diapason of values to divide into intervals
316 //================================================================================
318 void NumericalFunctor::GetHistogram(int nbIntervals,
319 std::vector<int>& nbEvents,
320 std::vector<double>& funValues,
321 const vector<int>& elements,
322 const double* minmax)
324 if ( nbIntervals < 1 ||
326 !myMesh->GetMeshInfo().NbElements( GetType() ))
328 nbEvents.resize( nbIntervals, 0 );
329 funValues.resize( nbIntervals+1 );
331 // get all values sorted
332 std::multiset< double > values;
333 if ( elements.empty() )
335 SMDS_ElemIteratorPtr elemIt = myMesh->elementsIterator(GetType());
336 while ( elemIt->more() )
337 values.insert( GetValue( elemIt->next()->GetID() ));
341 vector<int>::const_iterator id = elements.begin();
342 for ( ; id != elements.end(); ++id )
343 values.insert( GetValue( *id ));
348 funValues[0] = minmax[0];
349 funValues[nbIntervals] = minmax[1];
353 funValues[0] = *values.begin();
354 funValues[nbIntervals] = *values.rbegin();
356 // case nbIntervals == 1
357 if ( nbIntervals == 1 )
359 nbEvents[0] = values.size();
363 if (funValues.front() == funValues.back())
365 nbEvents.resize( 1 );
366 nbEvents[0] = values.size();
367 funValues[1] = funValues.back();
368 funValues.resize( 2 );
371 std::multiset< double >::iterator min = values.begin(), max;
372 for ( int i = 0; i < nbIntervals; ++i )
374 // find end value of i-th interval
375 double r = (i+1) / double( nbIntervals );
376 funValues[i+1] = funValues.front() * (1-r) + funValues.back() * r;
378 // count values in the i-th interval if there are any
379 if ( min != values.end() && *min <= funValues[i+1] )
381 // find the first value out of the interval
382 max = values.upper_bound( funValues[i+1] ); // max is greater than funValues[i+1], or end()
383 nbEvents[i] = std::distance( min, max );
387 // add values larger than minmax[1]
388 nbEvents.back() += std::distance( min, values.end() );
391 //=======================================================================
392 //function : GetValue
394 //=======================================================================
396 double Volume::GetValue( long theElementId )
398 if ( theElementId && myMesh ) {
399 SMDS_VolumeTool aVolumeTool;
400 if ( aVolumeTool.Set( myMesh->FindElement( theElementId )))
401 return aVolumeTool.GetSize();
406 //=======================================================================
407 //function : GetBadRate
408 //purpose : meaningless as it is not quality control functor
409 //=======================================================================
411 double Volume::GetBadRate( double Value, int /*nbNodes*/ ) const
416 //=======================================================================
419 //=======================================================================
421 SMDSAbs_ElementType Volume::GetType() const
423 return SMDSAbs_Volume;
428 Class : MaxElementLength2D
429 Description : Functor calculating maximum length of 2D element
432 double MaxElementLength2D::GetValue( long theElementId )
435 if( GetPoints( theElementId, P ) ) {
437 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
438 SMDSAbs_ElementType aType = aElem->GetType();
442 if( len == 3 ) { // triangles
443 double L1 = getDistance(P( 1 ),P( 2 ));
444 double L2 = getDistance(P( 2 ),P( 3 ));
445 double L3 = getDistance(P( 3 ),P( 1 ));
446 aVal = Max(L1,Max(L2,L3));
449 else if( len == 4 ) { // quadrangles
450 double L1 = getDistance(P( 1 ),P( 2 ));
451 double L2 = getDistance(P( 2 ),P( 3 ));
452 double L3 = getDistance(P( 3 ),P( 4 ));
453 double L4 = getDistance(P( 4 ),P( 1 ));
454 double D1 = getDistance(P( 1 ),P( 3 ));
455 double D2 = getDistance(P( 2 ),P( 4 ));
456 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
459 else if( len == 6 ) { // quadratic triangles
460 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
461 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
462 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
463 aVal = Max(L1,Max(L2,L3));
466 else if( len == 8 ) { // quadratic quadrangles
467 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
468 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
469 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
470 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
471 double D1 = getDistance(P( 1 ),P( 5 ));
472 double D2 = getDistance(P( 3 ),P( 7 ));
473 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
478 if( myPrecision >= 0 )
480 double prec = pow( 10., (double)myPrecision );
481 aVal = floor( aVal * prec + 0.5 ) / prec;
488 double MaxElementLength2D::GetBadRate( double Value, int /*nbNodes*/ ) const
493 SMDSAbs_ElementType MaxElementLength2D::GetType() const
499 Class : MaxElementLength3D
500 Description : Functor calculating maximum length of 3D element
503 double MaxElementLength3D::GetValue( long theElementId )
506 if( GetPoints( theElementId, P ) ) {
508 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
509 SMDSAbs_ElementType aType = aElem->GetType();
513 if( len == 4 ) { // tetras
514 double L1 = getDistance(P( 1 ),P( 2 ));
515 double L2 = getDistance(P( 2 ),P( 3 ));
516 double L3 = getDistance(P( 3 ),P( 1 ));
517 double L4 = getDistance(P( 1 ),P( 4 ));
518 double L5 = getDistance(P( 2 ),P( 4 ));
519 double L6 = getDistance(P( 3 ),P( 4 ));
520 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
523 else if( len == 5 ) { // pyramids
524 double L1 = getDistance(P( 1 ),P( 2 ));
525 double L2 = getDistance(P( 2 ),P( 3 ));
526 double L3 = getDistance(P( 3 ),P( 4 ));
527 double L4 = getDistance(P( 4 ),P( 1 ));
528 double L5 = getDistance(P( 1 ),P( 5 ));
529 double L6 = getDistance(P( 2 ),P( 5 ));
530 double L7 = getDistance(P( 3 ),P( 5 ));
531 double L8 = getDistance(P( 4 ),P( 5 ));
532 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
533 aVal = Max(aVal,Max(L7,L8));
536 else if( len == 6 ) { // pentas
537 double L1 = getDistance(P( 1 ),P( 2 ));
538 double L2 = getDistance(P( 2 ),P( 3 ));
539 double L3 = getDistance(P( 3 ),P( 1 ));
540 double L4 = getDistance(P( 4 ),P( 5 ));
541 double L5 = getDistance(P( 5 ),P( 6 ));
542 double L6 = getDistance(P( 6 ),P( 4 ));
543 double L7 = getDistance(P( 1 ),P( 4 ));
544 double L8 = getDistance(P( 2 ),P( 5 ));
545 double L9 = getDistance(P( 3 ),P( 6 ));
546 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
547 aVal = Max(aVal,Max(Max(L7,L8),L9));
550 else if( len == 8 ) { // hexas
551 double L1 = getDistance(P( 1 ),P( 2 ));
552 double L2 = getDistance(P( 2 ),P( 3 ));
553 double L3 = getDistance(P( 3 ),P( 4 ));
554 double L4 = getDistance(P( 4 ),P( 1 ));
555 double L5 = getDistance(P( 5 ),P( 6 ));
556 double L6 = getDistance(P( 6 ),P( 7 ));
557 double L7 = getDistance(P( 7 ),P( 8 ));
558 double L8 = getDistance(P( 8 ),P( 5 ));
559 double L9 = getDistance(P( 1 ),P( 5 ));
560 double L10= getDistance(P( 2 ),P( 6 ));
561 double L11= getDistance(P( 3 ),P( 7 ));
562 double L12= getDistance(P( 4 ),P( 8 ));
563 double D1 = getDistance(P( 1 ),P( 7 ));
564 double D2 = getDistance(P( 2 ),P( 8 ));
565 double D3 = getDistance(P( 3 ),P( 5 ));
566 double D4 = getDistance(P( 4 ),P( 6 ));
567 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
568 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
569 aVal = Max(aVal,Max(L11,L12));
570 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
573 else if( len == 10 ) { // quadratic tetras
574 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
575 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
576 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
577 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
578 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
579 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
580 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
583 else if( len == 13 ) { // quadratic pyramids
584 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
585 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
586 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
587 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
588 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
589 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
590 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
591 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
592 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
593 aVal = Max(aVal,Max(L7,L8));
596 else if( len == 15 ) { // quadratic pentas
597 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
598 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
599 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
600 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
601 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
602 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
603 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
604 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
605 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
606 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
607 aVal = Max(aVal,Max(Max(L7,L8),L9));
610 else if( len == 20 ) { // quadratic hexas
611 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
612 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
613 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
614 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
615 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
616 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
617 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
618 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
619 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
620 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
621 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
622 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
623 double D1 = getDistance(P( 1 ),P( 7 ));
624 double D2 = getDistance(P( 2 ),P( 8 ));
625 double D3 = getDistance(P( 3 ),P( 5 ));
626 double D4 = getDistance(P( 4 ),P( 6 ));
627 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
628 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
629 aVal = Max(aVal,Max(L11,L12));
630 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
633 else if( len > 1 && aElem->IsPoly() ) { // polys
634 // get the maximum distance between all pairs of nodes
635 for( int i = 1; i <= len; i++ ) {
636 for( int j = 1; j <= len; j++ ) {
637 if( j > i ) { // optimization of the loop
638 double D = getDistance( P(i), P(j) );
639 aVal = Max( aVal, D );
646 if( myPrecision >= 0 )
648 double prec = pow( 10., (double)myPrecision );
649 aVal = floor( aVal * prec + 0.5 ) / prec;
656 double MaxElementLength3D::GetBadRate( double Value, int /*nbNodes*/ ) const
661 SMDSAbs_ElementType MaxElementLength3D::GetType() const
663 return SMDSAbs_Volume;
669 Description : Functor for calculation of minimum angle
672 double MinimumAngle::GetValue( const TSequenceOfXYZ& P )
679 aMin = getAngle(P( P.size() ), P( 1 ), P( 2 ));
680 aMin = Min(aMin,getAngle(P( P.size()-1 ), P( P.size() ), P( 1 )));
682 for (int i=2; i<P.size();i++){
683 double A0 = getAngle( P( i-1 ), P( i ), P( i+1 ) );
687 return aMin * 180.0 / PI;
690 double MinimumAngle::GetBadRate( double Value, int nbNodes ) const
692 //const double aBestAngle = PI / nbNodes;
693 const double aBestAngle = 180.0 - ( 360.0 / double(nbNodes) );
694 return ( fabs( aBestAngle - Value ));
697 SMDSAbs_ElementType MinimumAngle::GetType() const
705 Description : Functor for calculating aspect ratio
707 double AspectRatio::GetValue( const TSequenceOfXYZ& P )
709 // According to "Mesh quality control" by Nadir Bouhamau referring to
710 // Pascal Jean Frey and Paul-Louis George. Maillages, applications aux elements finis.
711 // Hermes Science publications, Paris 1999 ISBN 2-7462-0024-4
714 int nbNodes = P.size();
719 // Compute aspect ratio
721 if ( nbNodes == 3 ) {
722 // Compute lengths of the sides
723 std::vector< double > aLen (nbNodes);
724 for ( int i = 0; i < nbNodes - 1; i++ )
725 aLen[ i ] = getDistance( P( i + 1 ), P( i + 2 ) );
726 aLen[ nbNodes - 1 ] = getDistance( P( 1 ), P( nbNodes ) );
727 // Q = alfa * h * p / S, where
729 // alfa = sqrt( 3 ) / 6
730 // h - length of the longest edge
731 // p - half perimeter
732 // S - triangle surface
733 const double alfa = sqrt( 3. ) / 6.;
734 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
735 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
736 double anArea = getArea( P( 1 ), P( 2 ), P( 3 ) );
737 if ( anArea <= Precision::Confusion() )
739 return alfa * maxLen * half_perimeter / anArea;
741 else if ( nbNodes == 6 ) { // quadratic triangles
742 // Compute lengths of the sides
743 std::vector< double > aLen (3);
744 aLen[0] = getDistance( P(1), P(3) );
745 aLen[1] = getDistance( P(3), P(5) );
746 aLen[2] = getDistance( P(5), P(1) );
747 // Q = alfa * h * p / S, where
749 // alfa = sqrt( 3 ) / 6
750 // h - length of the longest edge
751 // p - half perimeter
752 // S - triangle surface
753 const double alfa = sqrt( 3. ) / 6.;
754 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
755 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
756 double anArea = getArea( P(1), P(3), P(5) );
757 if ( anArea <= Precision::Confusion() )
759 return alfa * maxLen * half_perimeter / anArea;
761 else if( nbNodes == 4 ) { // quadrangle
762 // Compute lengths of the sides
763 std::vector< double > aLen (4);
764 aLen[0] = getDistance( P(1), P(2) );
765 aLen[1] = getDistance( P(2), P(3) );
766 aLen[2] = getDistance( P(3), P(4) );
767 aLen[3] = getDistance( P(4), P(1) );
768 // Compute lengths of the diagonals
769 std::vector< double > aDia (2);
770 aDia[0] = getDistance( P(1), P(3) );
771 aDia[1] = getDistance( P(2), P(4) );
772 // Compute areas of all triangles which can be built
773 // taking three nodes of the quadrangle
774 std::vector< double > anArea (4);
775 anArea[0] = getArea( P(1), P(2), P(3) );
776 anArea[1] = getArea( P(1), P(2), P(4) );
777 anArea[2] = getArea( P(1), P(3), P(4) );
778 anArea[3] = getArea( P(2), P(3), P(4) );
779 // Q = alpha * L * C1 / C2, where
781 // alpha = sqrt( 1/32 )
782 // L = max( L1, L2, L3, L4, D1, D2 )
783 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
784 // C2 = min( S1, S2, S3, S4 )
785 // Li - lengths of the edges
786 // Di - lengths of the diagonals
787 // Si - areas of the triangles
788 const double alpha = sqrt( 1 / 32. );
789 double L = Max( aLen[ 0 ],
793 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
794 double C1 = sqrt( ( aLen[0] * aLen[0] +
797 aLen[3] * aLen[3] ) / 4. );
798 double C2 = Min( anArea[ 0 ],
800 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
801 if ( C2 <= Precision::Confusion() )
803 return alpha * L * C1 / C2;
805 else if( nbNodes == 8 ){ // nbNodes==8 - quadratic quadrangle
806 // Compute lengths of the sides
807 std::vector< double > aLen (4);
808 aLen[0] = getDistance( P(1), P(3) );
809 aLen[1] = getDistance( P(3), P(5) );
810 aLen[2] = getDistance( P(5), P(7) );
811 aLen[3] = getDistance( P(7), P(1) );
812 // Compute lengths of the diagonals
813 std::vector< double > aDia (2);
814 aDia[0] = getDistance( P(1), P(5) );
815 aDia[1] = getDistance( P(3), P(7) );
816 // Compute areas of all triangles which can be built
817 // taking three nodes of the quadrangle
818 std::vector< double > anArea (4);
819 anArea[0] = getArea( P(1), P(3), P(5) );
820 anArea[1] = getArea( P(1), P(3), P(7) );
821 anArea[2] = getArea( P(1), P(5), P(7) );
822 anArea[3] = getArea( P(3), P(5), P(7) );
823 // Q = alpha * L * C1 / C2, where
825 // alpha = sqrt( 1/32 )
826 // L = max( L1, L2, L3, L4, D1, D2 )
827 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
828 // C2 = min( S1, S2, S3, S4 )
829 // Li - lengths of the edges
830 // Di - lengths of the diagonals
831 // Si - areas of the triangles
832 const double alpha = sqrt( 1 / 32. );
833 double L = Max( aLen[ 0 ],
837 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
838 double C1 = sqrt( ( aLen[0] * aLen[0] +
841 aLen[3] * aLen[3] ) / 4. );
842 double C2 = Min( anArea[ 0 ],
844 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
845 if ( C2 <= Precision::Confusion() )
847 return alpha * L * C1 / C2;
852 double AspectRatio::GetBadRate( double Value, int /*nbNodes*/ ) const
854 // the aspect ratio is in the range [1.0,infinity]
857 return Value / 1000.;
860 SMDSAbs_ElementType AspectRatio::GetType() const
867 Class : AspectRatio3D
868 Description : Functor for calculating aspect ratio
872 inline double getHalfPerimeter(double theTria[3]){
873 return (theTria[0] + theTria[1] + theTria[2])/2.0;
876 inline double getArea(double theHalfPerim, double theTria[3]){
877 return sqrt(theHalfPerim*
878 (theHalfPerim-theTria[0])*
879 (theHalfPerim-theTria[1])*
880 (theHalfPerim-theTria[2]));
883 inline double getVolume(double theLen[6]){
884 double a2 = theLen[0]*theLen[0];
885 double b2 = theLen[1]*theLen[1];
886 double c2 = theLen[2]*theLen[2];
887 double d2 = theLen[3]*theLen[3];
888 double e2 = theLen[4]*theLen[4];
889 double f2 = theLen[5]*theLen[5];
890 double P = 4.0*a2*b2*d2;
891 double Q = a2*(b2+d2-e2)-b2*(a2+d2-f2)-d2*(a2+b2-c2);
892 double R = (b2+d2-e2)*(a2+d2-f2)*(a2+d2-f2);
893 return sqrt(P-Q+R)/12.0;
896 inline double getVolume2(double theLen[6]){
897 double a2 = theLen[0]*theLen[0];
898 double b2 = theLen[1]*theLen[1];
899 double c2 = theLen[2]*theLen[2];
900 double d2 = theLen[3]*theLen[3];
901 double e2 = theLen[4]*theLen[4];
902 double f2 = theLen[5]*theLen[5];
904 double P = a2*e2*(b2+c2+d2+f2-a2-e2);
905 double Q = b2*f2*(a2+c2+d2+e2-b2-f2);
906 double R = c2*d2*(a2+b2+e2+f2-c2-d2);
907 double S = a2*b2*d2+b2*c2*e2+a2*c2*f2+d2*e2*f2;
909 return sqrt(P+Q+R-S)/12.0;
912 inline double getVolume(const TSequenceOfXYZ& P){
913 gp_Vec aVec1( P( 2 ) - P( 1 ) );
914 gp_Vec aVec2( P( 3 ) - P( 1 ) );
915 gp_Vec aVec3( P( 4 ) - P( 1 ) );
916 gp_Vec anAreaVec( aVec1 ^ aVec2 );
917 return fabs(aVec3 * anAreaVec) / 6.0;
920 inline double getMaxHeight(double theLen[6])
922 double aHeight = std::max(theLen[0],theLen[1]);
923 aHeight = std::max(aHeight,theLen[2]);
924 aHeight = std::max(aHeight,theLen[3]);
925 aHeight = std::max(aHeight,theLen[4]);
926 aHeight = std::max(aHeight,theLen[5]);
932 double AspectRatio3D::GetValue( const TSequenceOfXYZ& P )
934 double aQuality = 0.0;
935 if(myCurrElement->IsPoly()) return aQuality;
937 int nbNodes = P.size();
939 if(myCurrElement->IsQuadratic()) {
940 if(nbNodes==10) nbNodes=4; // quadratic tetrahedron
941 else if(nbNodes==13) nbNodes=5; // quadratic pyramid
942 else if(nbNodes==15) nbNodes=6; // quadratic pentahedron
943 else if(nbNodes==20) nbNodes=8; // quadratic hexahedron
944 else return aQuality;
950 getDistance(P( 1 ),P( 2 )), // a
951 getDistance(P( 2 ),P( 3 )), // b
952 getDistance(P( 3 ),P( 1 )), // c
953 getDistance(P( 2 ),P( 4 )), // d
954 getDistance(P( 3 ),P( 4 )), // e
955 getDistance(P( 1 ),P( 4 )) // f
957 double aTria[4][3] = {
958 {aLen[0],aLen[1],aLen[2]}, // abc
959 {aLen[0],aLen[3],aLen[5]}, // adf
960 {aLen[1],aLen[3],aLen[4]}, // bde
961 {aLen[2],aLen[4],aLen[5]} // cef
963 double aSumArea = 0.0;
964 double aHalfPerimeter = getHalfPerimeter(aTria[0]);
965 double anArea = getArea(aHalfPerimeter,aTria[0]);
967 aHalfPerimeter = getHalfPerimeter(aTria[1]);
968 anArea = getArea(aHalfPerimeter,aTria[1]);
970 aHalfPerimeter = getHalfPerimeter(aTria[2]);
971 anArea = getArea(aHalfPerimeter,aTria[2]);
973 aHalfPerimeter = getHalfPerimeter(aTria[3]);
974 anArea = getArea(aHalfPerimeter,aTria[3]);
976 double aVolume = getVolume(P);
977 //double aVolume = getVolume(aLen);
978 double aHeight = getMaxHeight(aLen);
979 static double aCoeff = sqrt(2.0)/12.0;
980 if ( aVolume > DBL_MIN )
981 aQuality = aCoeff*aHeight*aSumArea/aVolume;
986 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 3 ),P( 5 )};
987 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
990 gp_XYZ aXYZ[4] = {P( 1 ),P( 3 ),P( 4 ),P( 5 )};
991 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
994 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 5 )};
995 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
998 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 4 ),P( 5 )};
999 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1005 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 6 )};
1006 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1009 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 3 )};
1010 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1013 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 6 )};
1014 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1017 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1018 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1021 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 6 )};
1022 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1025 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 3 )};
1026 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1032 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1033 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1036 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 4 )};
1037 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1040 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 7 )};
1041 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1044 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 8 )};
1045 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1048 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 3 )};
1049 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1052 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 4 )};
1053 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1056 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 7 )};
1057 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1060 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 8 )};
1061 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1064 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 3 )};
1065 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1068 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 4 )};
1069 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1072 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 7 )};
1073 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1076 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 8 )};
1077 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1080 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 1 )};
1081 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1084 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 2 )};
1085 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1088 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 5 )};
1089 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1092 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 6 )};
1093 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1096 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 1 )};
1097 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1100 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 2 )};
1101 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1104 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 5 )};
1105 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1108 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 6 )};
1109 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1112 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 1 )};
1113 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1116 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1117 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1120 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 5 )};
1121 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1124 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 6 )};
1125 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1128 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1129 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1132 gp_XYZ aXYZ[4] = {P( 4 ),P( 5 ),P( 8 ),P( 2 )};
1133 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1136 gp_XYZ aXYZ[4] = {P( 1 ),P( 4 ),P( 5 ),P( 3 )};
1137 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1140 gp_XYZ aXYZ[4] = {P( 3 ),P( 6 ),P( 7 ),P( 1 )};
1141 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1144 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 6 ),P( 4 )};
1145 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1148 gp_XYZ aXYZ[4] = {P( 5 ),P( 6 ),P( 8 ),P( 3 )};
1149 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1152 gp_XYZ aXYZ[4] = {P( 7 ),P( 8 ),P( 6 ),P( 1 )};
1153 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1156 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 7 )};
1157 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1160 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 2 ),P( 5 )};
1161 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1166 if ( nbNodes > 4 ) {
1167 // avaluate aspect ratio of quadranle faces
1168 AspectRatio aspect2D;
1169 SMDS_VolumeTool::VolumeType type = SMDS_VolumeTool::GetType( nbNodes );
1170 int nbFaces = SMDS_VolumeTool::NbFaces( type );
1171 TSequenceOfXYZ points(4);
1172 for ( int i = 0; i < nbFaces; ++i ) { // loop on faces of a volume
1173 if ( SMDS_VolumeTool::NbFaceNodes( type, i ) != 4 )
1175 const int* pInd = SMDS_VolumeTool::GetFaceNodesIndices( type, i, true );
1176 for ( int p = 0; p < 4; ++p ) // loop on nodes of a quadranle face
1177 points( p + 1 ) = P( pInd[ p ] + 1 );
1178 aQuality = std::max( aQuality, aspect2D.GetValue( points ));
1184 double AspectRatio3D::GetBadRate( double Value, int /*nbNodes*/ ) const
1186 // the aspect ratio is in the range [1.0,infinity]
1189 return Value / 1000.;
1192 SMDSAbs_ElementType AspectRatio3D::GetType() const
1194 return SMDSAbs_Volume;
1200 Description : Functor for calculating warping
1202 double Warping::GetValue( const TSequenceOfXYZ& P )
1204 if ( P.size() != 4 )
1207 gp_XYZ G = ( P( 1 ) + P( 2 ) + P( 3 ) + P( 4 ) ) / 4.;
1209 double A1 = ComputeA( P( 1 ), P( 2 ), P( 3 ), G );
1210 double A2 = ComputeA( P( 2 ), P( 3 ), P( 4 ), G );
1211 double A3 = ComputeA( P( 3 ), P( 4 ), P( 1 ), G );
1212 double A4 = ComputeA( P( 4 ), P( 1 ), P( 2 ), G );
1214 return Max( Max( A1, A2 ), Max( A3, A4 ) );
1217 double Warping::ComputeA( const gp_XYZ& thePnt1,
1218 const gp_XYZ& thePnt2,
1219 const gp_XYZ& thePnt3,
1220 const gp_XYZ& theG ) const
1222 double aLen1 = gp_Pnt( thePnt1 ).Distance( gp_Pnt( thePnt2 ) );
1223 double aLen2 = gp_Pnt( thePnt2 ).Distance( gp_Pnt( thePnt3 ) );
1224 double L = Min( aLen1, aLen2 ) * 0.5;
1225 if ( L < Precision::Confusion())
1228 gp_XYZ GI = ( thePnt2 + thePnt1 ) / 2. - theG;
1229 gp_XYZ GJ = ( thePnt3 + thePnt2 ) / 2. - theG;
1230 gp_XYZ N = GI.Crossed( GJ );
1232 if ( N.Modulus() < gp::Resolution() )
1237 double H = ( thePnt2 - theG ).Dot( N );
1238 return asin( fabs( H / L ) ) * 180. / PI;
1241 double Warping::GetBadRate( double Value, int /*nbNodes*/ ) const
1243 // the warp is in the range [0.0,PI/2]
1244 // 0.0 = good (no warp)
1245 // PI/2 = bad (face pliee)
1249 SMDSAbs_ElementType Warping::GetType() const
1251 return SMDSAbs_Face;
1257 Description : Functor for calculating taper
1259 double Taper::GetValue( const TSequenceOfXYZ& P )
1261 if ( P.size() != 4 )
1265 double J1 = getArea( P( 4 ), P( 1 ), P( 2 ) ) / 2.;
1266 double J2 = getArea( P( 3 ), P( 1 ), P( 2 ) ) / 2.;
1267 double J3 = getArea( P( 2 ), P( 3 ), P( 4 ) ) / 2.;
1268 double J4 = getArea( P( 3 ), P( 4 ), P( 1 ) ) / 2.;
1270 double JA = 0.25 * ( J1 + J2 + J3 + J4 );
1271 if ( JA <= Precision::Confusion() )
1274 double T1 = fabs( ( J1 - JA ) / JA );
1275 double T2 = fabs( ( J2 - JA ) / JA );
1276 double T3 = fabs( ( J3 - JA ) / JA );
1277 double T4 = fabs( ( J4 - JA ) / JA );
1279 return Max( Max( T1, T2 ), Max( T3, T4 ) );
1282 double Taper::GetBadRate( double Value, int /*nbNodes*/ ) const
1284 // the taper is in the range [0.0,1.0]
1285 // 0.0 = good (no taper)
1286 // 1.0 = bad (les cotes opposes sont allignes)
1290 SMDSAbs_ElementType Taper::GetType() const
1292 return SMDSAbs_Face;
1298 Description : Functor for calculating skew in degrees
1300 static inline double skewAngle( const gp_XYZ& p1, const gp_XYZ& p2, const gp_XYZ& p3 )
1302 gp_XYZ p12 = ( p2 + p1 ) / 2.;
1303 gp_XYZ p23 = ( p3 + p2 ) / 2.;
1304 gp_XYZ p31 = ( p3 + p1 ) / 2.;
1306 gp_Vec v1( p31 - p2 ), v2( p12 - p23 );
1308 return v1.Magnitude() < gp::Resolution() || v2.Magnitude() < gp::Resolution() ? 0. : v1.Angle( v2 );
1311 double Skew::GetValue( const TSequenceOfXYZ& P )
1313 if ( P.size() != 3 && P.size() != 4 )
1317 static double PI2 = PI / 2.;
1318 if ( P.size() == 3 )
1320 double A0 = fabs( PI2 - skewAngle( P( 3 ), P( 1 ), P( 2 ) ) );
1321 double A1 = fabs( PI2 - skewAngle( P( 1 ), P( 2 ), P( 3 ) ) );
1322 double A2 = fabs( PI2 - skewAngle( P( 2 ), P( 3 ), P( 1 ) ) );
1324 return Max( A0, Max( A1, A2 ) ) * 180. / PI;
1328 gp_XYZ p12 = ( P( 1 ) + P( 2 ) ) / 2.;
1329 gp_XYZ p23 = ( P( 2 ) + P( 3 ) ) / 2.;
1330 gp_XYZ p34 = ( P( 3 ) + P( 4 ) ) / 2.;
1331 gp_XYZ p41 = ( P( 4 ) + P( 1 ) ) / 2.;
1333 gp_Vec v1( p34 - p12 ), v2( p23 - p41 );
1334 double A = v1.Magnitude() <= gp::Resolution() || v2.Magnitude() <= gp::Resolution()
1335 ? 0. : fabs( PI2 - v1.Angle( v2 ) );
1338 if ( A < Precision::Angular() )
1341 return A * 180. / PI;
1345 double Skew::GetBadRate( double Value, int /*nbNodes*/ ) const
1347 // the skew is in the range [0.0,PI/2].
1353 SMDSAbs_ElementType Skew::GetType() const
1355 return SMDSAbs_Face;
1361 Description : Functor for calculating area
1363 double Area::GetValue( const TSequenceOfXYZ& P )
1366 if ( P.size() > 2 ) {
1367 gp_Vec aVec1( P(2) - P(1) );
1368 gp_Vec aVec2( P(3) - P(1) );
1369 gp_Vec SumVec = aVec1 ^ aVec2;
1370 for (int i=4; i<=P.size(); i++) {
1371 gp_Vec aVec1( P(i-1) - P(1) );
1372 gp_Vec aVec2( P(i) - P(1) );
1373 gp_Vec tmp = aVec1 ^ aVec2;
1376 val = SumVec.Magnitude() * 0.5;
1381 double Area::GetBadRate( double Value, int /*nbNodes*/ ) const
1383 // meaningless as it is not a quality control functor
1387 SMDSAbs_ElementType Area::GetType() const
1389 return SMDSAbs_Face;
1395 Description : Functor for calculating length off edge
1397 double Length::GetValue( const TSequenceOfXYZ& P )
1399 switch ( P.size() ) {
1400 case 2: return getDistance( P( 1 ), P( 2 ) );
1401 case 3: return getDistance( P( 1 ), P( 2 ) ) + getDistance( P( 2 ), P( 3 ) );
1406 double Length::GetBadRate( double Value, int /*nbNodes*/ ) const
1408 // meaningless as it is not quality control functor
1412 SMDSAbs_ElementType Length::GetType() const
1414 return SMDSAbs_Edge;
1419 Description : Functor for calculating length of edge
1422 double Length2D::GetValue( long theElementId)
1426 //cout<<"Length2D::GetValue"<<endl;
1427 if (GetPoints(theElementId,P)){
1428 //for(int jj=1; jj<=P.size(); jj++)
1429 // cout<<"jj="<<jj<<" P("<<P(jj).X()<<","<<P(jj).Y()<<","<<P(jj).Z()<<")"<<endl;
1431 double aVal;// = GetValue( P );
1432 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
1433 SMDSAbs_ElementType aType = aElem->GetType();
1442 aVal = getDistance( P( 1 ), P( 2 ) );
1445 else if (len == 3){ // quadratic edge
1446 aVal = getDistance(P( 1 ),P( 3 )) + getDistance(P( 3 ),P( 2 ));
1450 if (len == 3){ // triangles
1451 double L1 = getDistance(P( 1 ),P( 2 ));
1452 double L2 = getDistance(P( 2 ),P( 3 ));
1453 double L3 = getDistance(P( 3 ),P( 1 ));
1454 aVal = Max(L1,Max(L2,L3));
1457 else if (len == 4){ // quadrangles
1458 double L1 = getDistance(P( 1 ),P( 2 ));
1459 double L2 = getDistance(P( 2 ),P( 3 ));
1460 double L3 = getDistance(P( 3 ),P( 4 ));
1461 double L4 = getDistance(P( 4 ),P( 1 ));
1462 aVal = Max(Max(L1,L2),Max(L3,L4));
1465 if (len == 6){ // quadratic triangles
1466 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1467 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1468 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
1469 aVal = Max(L1,Max(L2,L3));
1470 //cout<<"L1="<<L1<<" L2="<<L2<<"L3="<<L3<<" aVal="<<aVal<<endl;
1473 else if (len == 8){ // quadratic quadrangles
1474 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1475 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1476 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
1477 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
1478 aVal = Max(Max(L1,L2),Max(L3,L4));
1481 case SMDSAbs_Volume:
1482 if (len == 4){ // tetraidrs
1483 double L1 = getDistance(P( 1 ),P( 2 ));
1484 double L2 = getDistance(P( 2 ),P( 3 ));
1485 double L3 = getDistance(P( 3 ),P( 1 ));
1486 double L4 = getDistance(P( 1 ),P( 4 ));
1487 double L5 = getDistance(P( 2 ),P( 4 ));
1488 double L6 = getDistance(P( 3 ),P( 4 ));
1489 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1492 else if (len == 5){ // piramids
1493 double L1 = getDistance(P( 1 ),P( 2 ));
1494 double L2 = getDistance(P( 2 ),P( 3 ));
1495 double L3 = getDistance(P( 3 ),P( 4 ));
1496 double L4 = getDistance(P( 4 ),P( 1 ));
1497 double L5 = getDistance(P( 1 ),P( 5 ));
1498 double L6 = getDistance(P( 2 ),P( 5 ));
1499 double L7 = getDistance(P( 3 ),P( 5 ));
1500 double L8 = getDistance(P( 4 ),P( 5 ));
1502 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1503 aVal = Max(aVal,Max(L7,L8));
1506 else if (len == 6){ // pentaidres
1507 double L1 = getDistance(P( 1 ),P( 2 ));
1508 double L2 = getDistance(P( 2 ),P( 3 ));
1509 double L3 = getDistance(P( 3 ),P( 1 ));
1510 double L4 = getDistance(P( 4 ),P( 5 ));
1511 double L5 = getDistance(P( 5 ),P( 6 ));
1512 double L6 = getDistance(P( 6 ),P( 4 ));
1513 double L7 = getDistance(P( 1 ),P( 4 ));
1514 double L8 = getDistance(P( 2 ),P( 5 ));
1515 double L9 = getDistance(P( 3 ),P( 6 ));
1517 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1518 aVal = Max(aVal,Max(Max(L7,L8),L9));
1521 else if (len == 8){ // hexaider
1522 double L1 = getDistance(P( 1 ),P( 2 ));
1523 double L2 = getDistance(P( 2 ),P( 3 ));
1524 double L3 = getDistance(P( 3 ),P( 4 ));
1525 double L4 = getDistance(P( 4 ),P( 1 ));
1526 double L5 = getDistance(P( 5 ),P( 6 ));
1527 double L6 = getDistance(P( 6 ),P( 7 ));
1528 double L7 = getDistance(P( 7 ),P( 8 ));
1529 double L8 = getDistance(P( 8 ),P( 5 ));
1530 double L9 = getDistance(P( 1 ),P( 5 ));
1531 double L10= getDistance(P( 2 ),P( 6 ));
1532 double L11= getDistance(P( 3 ),P( 7 ));
1533 double L12= getDistance(P( 4 ),P( 8 ));
1535 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1536 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1537 aVal = Max(aVal,Max(L11,L12));
1542 if (len == 10){ // quadratic tetraidrs
1543 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
1544 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
1545 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
1546 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1547 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
1548 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
1549 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1552 else if (len == 13){ // quadratic piramids
1553 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
1554 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
1555 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1556 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1557 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1558 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
1559 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
1560 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
1561 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1562 aVal = Max(aVal,Max(L7,L8));
1565 else if (len == 15){ // quadratic pentaidres
1566 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
1567 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
1568 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1569 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1570 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
1571 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
1572 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
1573 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
1574 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
1575 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1576 aVal = Max(aVal,Max(Max(L7,L8),L9));
1579 else if (len == 20){ // quadratic hexaider
1580 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
1581 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
1582 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
1583 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
1584 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
1585 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
1586 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
1587 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
1588 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
1589 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
1590 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
1591 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
1592 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1593 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1594 aVal = Max(aVal,Max(L11,L12));
1606 if ( myPrecision >= 0 )
1608 double prec = pow( 10., (double)( myPrecision ) );
1609 aVal = floor( aVal * prec + 0.5 ) / prec;
1618 double Length2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1620 // meaningless as it is not quality control functor
1624 SMDSAbs_ElementType Length2D::GetType() const
1626 return SMDSAbs_Face;
1629 Length2D::Value::Value(double theLength,long thePntId1, long thePntId2):
1632 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1633 if(thePntId1 > thePntId2){
1634 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1638 bool Length2D::Value::operator<(const Length2D::Value& x) const{
1639 if(myPntId[0] < x.myPntId[0]) return true;
1640 if(myPntId[0] == x.myPntId[0])
1641 if(myPntId[1] < x.myPntId[1]) return true;
1645 void Length2D::GetValues(TValues& theValues){
1647 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1648 for(; anIter->more(); ){
1649 const SMDS_MeshFace* anElem = anIter->next();
1651 if(anElem->IsQuadratic()) {
1652 const SMDS_VtkFace* F =
1653 dynamic_cast<const SMDS_VtkFace*>(anElem);
1654 // use special nodes iterator
1655 SMDS_ElemIteratorPtr anIter = F->interlacedNodesElemIterator();
1660 const SMDS_MeshElement* aNode;
1662 aNode = anIter->next();
1663 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1664 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1665 aNodeId[0] = aNodeId[1] = aNode->GetID();
1668 for(; anIter->more(); ){
1669 const SMDS_MeshNode* N1 = static_cast<const SMDS_MeshNode*> (anIter->next());
1670 P[2] = gp_Pnt(N1->X(),N1->Y(),N1->Z());
1671 aNodeId[2] = N1->GetID();
1672 aLength = P[1].Distance(P[2]);
1673 if(!anIter->more()) break;
1674 const SMDS_MeshNode* N2 = static_cast<const SMDS_MeshNode*> (anIter->next());
1675 P[3] = gp_Pnt(N2->X(),N2->Y(),N2->Z());
1676 aNodeId[3] = N2->GetID();
1677 aLength += P[2].Distance(P[3]);
1678 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1679 Value aValue2(aLength,aNodeId[2],aNodeId[3]);
1681 aNodeId[1] = aNodeId[3];
1682 theValues.insert(aValue1);
1683 theValues.insert(aValue2);
1685 aLength += P[2].Distance(P[0]);
1686 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1687 Value aValue2(aLength,aNodeId[2],aNodeId[0]);
1688 theValues.insert(aValue1);
1689 theValues.insert(aValue2);
1692 SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
1697 const SMDS_MeshElement* aNode;
1698 if(aNodesIter->more()){
1699 aNode = aNodesIter->next();
1700 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1701 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1702 aNodeId[0] = aNodeId[1] = aNode->GetID();
1705 for(; aNodesIter->more(); ){
1706 aNode = aNodesIter->next();
1707 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1708 long anId = aNode->GetID();
1710 P[2] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1712 aLength = P[1].Distance(P[2]);
1714 Value aValue(aLength,aNodeId[1],anId);
1717 theValues.insert(aValue);
1720 aLength = P[0].Distance(P[1]);
1722 Value aValue(aLength,aNodeId[0],aNodeId[1]);
1723 theValues.insert(aValue);
1729 Class : MultiConnection
1730 Description : Functor for calculating number of faces conneted to the edge
1732 double MultiConnection::GetValue( const TSequenceOfXYZ& P )
1736 double MultiConnection::GetValue( long theId )
1738 return getNbMultiConnection( myMesh, theId );
1741 double MultiConnection::GetBadRate( double Value, int /*nbNodes*/ ) const
1743 // meaningless as it is not quality control functor
1747 SMDSAbs_ElementType MultiConnection::GetType() const
1749 return SMDSAbs_Edge;
1753 Class : MultiConnection2D
1754 Description : Functor for calculating number of faces conneted to the edge
1756 double MultiConnection2D::GetValue( const TSequenceOfXYZ& P )
1761 double MultiConnection2D::GetValue( long theElementId )
1765 const SMDS_MeshElement* aFaceElem = myMesh->FindElement(theElementId);
1766 SMDSAbs_ElementType aType = aFaceElem->GetType();
1771 int i = 0, len = aFaceElem->NbNodes();
1772 SMDS_ElemIteratorPtr anIter = aFaceElem->nodesIterator();
1775 const SMDS_MeshNode *aNode, *aNode0;
1776 TColStd_MapOfInteger aMap, aMapPrev;
1778 for (i = 0; i <= len; i++) {
1783 if (anIter->more()) {
1784 aNode = (SMDS_MeshNode*)anIter->next();
1792 if (i == 0) aNode0 = aNode;
1794 SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
1795 while (anElemIter->more()) {
1796 const SMDS_MeshElement* anElem = anElemIter->next();
1797 if (anElem != 0 && anElem->GetType() == SMDSAbs_Face) {
1798 int anId = anElem->GetID();
1801 if (aMapPrev.Contains(anId)) {
1806 aResult = Max(aResult, aNb);
1817 double MultiConnection2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1819 // meaningless as it is not quality control functor
1823 SMDSAbs_ElementType MultiConnection2D::GetType() const
1825 return SMDSAbs_Face;
1828 MultiConnection2D::Value::Value(long thePntId1, long thePntId2)
1830 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1831 if(thePntId1 > thePntId2){
1832 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1836 bool MultiConnection2D::Value::operator<(const MultiConnection2D::Value& x) const{
1837 if(myPntId[0] < x.myPntId[0]) return true;
1838 if(myPntId[0] == x.myPntId[0])
1839 if(myPntId[1] < x.myPntId[1]) return true;
1843 void MultiConnection2D::GetValues(MValues& theValues){
1844 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1845 for(; anIter->more(); ){
1846 const SMDS_MeshFace* anElem = anIter->next();
1847 SMDS_ElemIteratorPtr aNodesIter;
1848 if ( anElem->IsQuadratic() )
1849 aNodesIter = dynamic_cast<const SMDS_VtkFace*>
1850 (anElem)->interlacedNodesElemIterator();
1852 aNodesIter = anElem->nodesIterator();
1855 //int aNbConnects=0;
1856 const SMDS_MeshNode* aNode0;
1857 const SMDS_MeshNode* aNode1;
1858 const SMDS_MeshNode* aNode2;
1859 if(aNodesIter->more()){
1860 aNode0 = (SMDS_MeshNode*) aNodesIter->next();
1862 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode1;
1863 aNodeId[0] = aNodeId[1] = aNodes->GetID();
1865 for(; aNodesIter->more(); ) {
1866 aNode2 = (SMDS_MeshNode*) aNodesIter->next();
1867 long anId = aNode2->GetID();
1870 Value aValue(aNodeId[1],aNodeId[2]);
1871 MValues::iterator aItr = theValues.find(aValue);
1872 if (aItr != theValues.end()){
1877 theValues[aValue] = 1;
1880 //cout << "NodeIds: "<<aNodeId[1]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1881 aNodeId[1] = aNodeId[2];
1884 Value aValue(aNodeId[0],aNodeId[2]);
1885 MValues::iterator aItr = theValues.find(aValue);
1886 if (aItr != theValues.end()) {
1891 theValues[aValue] = 1;
1894 //cout << "NodeIds: "<<aNodeId[0]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1904 Class : BadOrientedVolume
1905 Description : Predicate bad oriented volumes
1908 BadOrientedVolume::BadOrientedVolume()
1913 void BadOrientedVolume::SetMesh( const SMDS_Mesh* theMesh )
1918 bool BadOrientedVolume::IsSatisfy( long theId )
1923 SMDS_VolumeTool vTool( myMesh->FindElement( theId ));
1924 return !vTool.IsForward();
1927 SMDSAbs_ElementType BadOrientedVolume::GetType() const
1929 return SMDSAbs_Volume;
1933 Class : BareBorderVolume
1936 bool BareBorderVolume::IsSatisfy(long theElementId )
1938 SMDS_VolumeTool myTool;
1939 if ( myTool.Set( myMesh->FindElement(theElementId)))
1941 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
1942 if ( myTool.IsFreeFace( iF ))
1944 const SMDS_MeshNode** n = myTool.GetFaceNodes(iF);
1945 vector< const SMDS_MeshNode*> nodes( n, n+myTool.NbFaceNodes(iF));
1946 if ( !myMesh->FindElement( nodes, SMDSAbs_Face, /*Nomedium=*/false))
1954 Class : BareBorderFace
1957 bool BareBorderFace::IsSatisfy(long theElementId )
1959 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
1960 if ( face->GetType() == SMDSAbs_Face )
1962 int nbN = face->NbCornerNodes();
1963 for ( int i = 0; i < nbN; ++i )
1965 // check if a link is shared by another face
1966 const SMDS_MeshNode* n1 = face->GetNode( i );
1967 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
1968 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
1969 bool isShared = false;
1970 while ( !isShared && fIt->more() )
1972 const SMDS_MeshElement* f = fIt->next();
1973 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
1977 myLinkNodes.resize( 2 + face->IsQuadratic());
1978 myLinkNodes[0] = n1;
1979 myLinkNodes[1] = n2;
1980 if ( face->IsQuadratic() )
1981 myLinkNodes[2] = face->GetNode( i+nbN );
1982 return !myMesh->FindElement( myLinkNodes, SMDSAbs_Edge, /*noMedium=*/false);
1990 Class : OverConstrainedVolume
1993 bool OverConstrainedVolume::IsSatisfy(long theElementId )
1995 // An element is over-constrained if it has N-1 free borders where
1996 // N is the number of edges/faces for a 2D/3D element.
1997 SMDS_VolumeTool myTool;
1998 if ( myTool.Set( myMesh->FindElement(theElementId)))
2000 int nbSharedFaces = 0;
2001 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2002 if ( !myTool.IsFreeFace( iF ) && ++nbSharedFaces > 1 )
2004 return ( nbSharedFaces == 1 );
2010 Class : OverConstrainedFace
2013 bool OverConstrainedFace::IsSatisfy(long theElementId )
2015 // An element is over-constrained if it has N-1 free borders where
2016 // N is the number of edges/faces for a 2D/3D element.
2017 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2018 if ( face->GetType() == SMDSAbs_Face )
2020 int nbSharedBorders = 0;
2021 int nbN = face->NbCornerNodes();
2022 for ( int i = 0; i < nbN; ++i )
2024 // check if a link is shared by another face
2025 const SMDS_MeshNode* n1 = face->GetNode( i );
2026 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2027 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2028 bool isShared = false;
2029 while ( !isShared && fIt->more() )
2031 const SMDS_MeshElement* f = fIt->next();
2032 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2034 if ( isShared && ++nbSharedBorders > 1 )
2037 return ( nbSharedBorders == 1 );
2044 Description : Predicate for free borders
2047 FreeBorders::FreeBorders()
2052 void FreeBorders::SetMesh( const SMDS_Mesh* theMesh )
2057 bool FreeBorders::IsSatisfy( long theId )
2059 return getNbMultiConnection( myMesh, theId ) == 1;
2062 SMDSAbs_ElementType FreeBorders::GetType() const
2064 return SMDSAbs_Edge;
2070 Description : Predicate for free Edges
2072 FreeEdges::FreeEdges()
2077 void FreeEdges::SetMesh( const SMDS_Mesh* theMesh )
2082 bool FreeEdges::IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId )
2084 TColStd_MapOfInteger aMap;
2085 for ( int i = 0; i < 2; i++ )
2087 SMDS_ElemIteratorPtr anElemIter = theNodes[ i ]->GetInverseElementIterator();
2088 while( anElemIter->more() )
2090 const SMDS_MeshElement* anElem = anElemIter->next();
2091 if ( anElem != 0 && anElem->GetType() == SMDSAbs_Face )
2093 int anId = anElem->GetID();
2097 else if ( aMap.Contains( anId ) && anId != theFaceId )
2105 bool FreeEdges::IsSatisfy( long theId )
2110 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2111 if ( aFace == 0 || aFace->GetType() != SMDSAbs_Face || aFace->NbNodes() < 3 )
2114 SMDS_ElemIteratorPtr anIter;
2115 if ( aFace->IsQuadratic() ) {
2116 anIter = dynamic_cast<const SMDS_VtkFace*>
2117 (aFace)->interlacedNodesElemIterator();
2120 anIter = aFace->nodesIterator();
2125 int i = 0, nbNodes = aFace->NbNodes();
2126 std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
2127 while( anIter->more() )
2129 const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
2132 aNodes[ i++ ] = aNode;
2134 aNodes[ nbNodes ] = aNodes[ 0 ];
2136 for ( i = 0; i < nbNodes; i++ )
2137 if ( IsFreeEdge( &aNodes[ i ], theId ) )
2143 SMDSAbs_ElementType FreeEdges::GetType() const
2145 return SMDSAbs_Face;
2148 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2151 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
2152 if(thePntId1 > thePntId2){
2153 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
2157 bool FreeEdges::Border::operator<(const FreeEdges::Border& x) const{
2158 if(myPntId[0] < x.myPntId[0]) return true;
2159 if(myPntId[0] == x.myPntId[0])
2160 if(myPntId[1] < x.myPntId[1]) return true;
2164 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2165 FreeEdges::TBorders& theRegistry,
2166 FreeEdges::TBorders& theContainer)
2168 if(theRegistry.find(theBorder) == theRegistry.end()){
2169 theRegistry.insert(theBorder);
2170 theContainer.insert(theBorder);
2172 theContainer.erase(theBorder);
2176 void FreeEdges::GetBoreders(TBorders& theBorders)
2179 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2180 for(; anIter->more(); ){
2181 const SMDS_MeshFace* anElem = anIter->next();
2182 long anElemId = anElem->GetID();
2183 SMDS_ElemIteratorPtr aNodesIter;
2184 if ( anElem->IsQuadratic() )
2185 aNodesIter = static_cast<const SMDS_VtkFace*>(anElem)->
2186 interlacedNodesElemIterator();
2188 aNodesIter = anElem->nodesIterator();
2190 const SMDS_MeshElement* aNode;
2191 if(aNodesIter->more()){
2192 aNode = aNodesIter->next();
2193 aNodeId[0] = aNodeId[1] = aNode->GetID();
2195 for(; aNodesIter->more(); ){
2196 aNode = aNodesIter->next();
2197 long anId = aNode->GetID();
2198 Border aBorder(anElemId,aNodeId[1],anId);
2200 //std::cout<<aBorder.myPntId[0]<<"; "<<aBorder.myPntId[1]<<"; "<<aBorder.myElemId<<endl;
2201 UpdateBorders(aBorder,aRegistry,theBorders);
2203 Border aBorder(anElemId,aNodeId[0],aNodeId[1]);
2204 //std::cout<<aBorder.myPntId[0]<<"; "<<aBorder.myPntId[1]<<"; "<<aBorder.myElemId<<endl;
2205 UpdateBorders(aBorder,aRegistry,theBorders);
2207 //std::cout<<"theBorders.size() = "<<theBorders.size()<<endl;
2213 Description : Predicate for free nodes
2216 FreeNodes::FreeNodes()
2221 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2226 bool FreeNodes::IsSatisfy( long theNodeId )
2228 const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2232 return (aNode->NbInverseElements() < 1);
2235 SMDSAbs_ElementType FreeNodes::GetType() const
2237 return SMDSAbs_Node;
2243 Description : Predicate for free faces
2246 FreeFaces::FreeFaces()
2251 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2256 bool FreeFaces::IsSatisfy( long theId )
2258 if (!myMesh) return false;
2259 // check that faces nodes refers to less than two common volumes
2260 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2261 if ( !aFace || aFace->GetType() != SMDSAbs_Face )
2264 int nbNode = aFace->NbNodes();
2266 // collect volumes check that number of volumss with count equal nbNode not less than 2
2267 typedef map< SMDS_MeshElement*, int > TMapOfVolume; // map of volume counters
2268 typedef map< SMDS_MeshElement*, int >::iterator TItrMapOfVolume; // iterator
2269 TMapOfVolume mapOfVol;
2271 SMDS_ElemIteratorPtr nodeItr = aFace->nodesIterator();
2272 while ( nodeItr->more() ) {
2273 const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>(nodeItr->next());
2274 if ( !aNode ) continue;
2275 SMDS_ElemIteratorPtr volItr = aNode->GetInverseElementIterator(SMDSAbs_Volume);
2276 while ( volItr->more() ) {
2277 SMDS_MeshElement* aVol = (SMDS_MeshElement*)volItr->next();
2278 TItrMapOfVolume itr = mapOfVol.insert(make_pair(aVol, 0)).first;
2283 TItrMapOfVolume volItr = mapOfVol.begin();
2284 TItrMapOfVolume volEnd = mapOfVol.end();
2285 for ( ; volItr != volEnd; ++volItr )
2286 if ( (*volItr).second >= nbNode )
2288 // face is not free if number of volumes constructed on thier nodes more than one
2292 SMDSAbs_ElementType FreeFaces::GetType() const
2294 return SMDSAbs_Face;
2298 Class : LinearOrQuadratic
2299 Description : Predicate to verify whether a mesh element is linear
2302 LinearOrQuadratic::LinearOrQuadratic()
2307 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2312 bool LinearOrQuadratic::IsSatisfy( long theId )
2314 if (!myMesh) return false;
2315 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2316 if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2318 return (!anElem->IsQuadratic());
2321 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2326 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2333 Description : Functor for check color of group to whic mesh element belongs to
2336 GroupColor::GroupColor()
2340 bool GroupColor::IsSatisfy( long theId )
2342 return (myIDs.find( theId ) != myIDs.end());
2345 void GroupColor::SetType( SMDSAbs_ElementType theType )
2350 SMDSAbs_ElementType GroupColor::GetType() const
2355 static bool isEqual( const Quantity_Color& theColor1,
2356 const Quantity_Color& theColor2 )
2358 // tolerance to compare colors
2359 const double tol = 5*1e-3;
2360 return ( fabs( theColor1.Red() - theColor2.Red() ) < tol &&
2361 fabs( theColor1.Green() - theColor2.Green() ) < tol &&
2362 fabs( theColor1.Blue() - theColor2.Blue() ) < tol );
2366 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2370 const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2374 int nbGrp = aMesh->GetNbGroups();
2378 // iterates on groups and find necessary elements ids
2379 const std::set<SMESHDS_GroupBase*>& aGroups = aMesh->GetGroups();
2380 set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
2381 for (; GrIt != aGroups.end(); GrIt++) {
2382 SMESHDS_GroupBase* aGrp = (*GrIt);
2385 // check type and color of group
2386 if ( !isEqual( myColor, aGrp->GetColor() ) )
2388 if ( myType != SMDSAbs_All && myType != (SMDSAbs_ElementType)aGrp->GetType() )
2391 SMDSAbs_ElementType aGrpElType = (SMDSAbs_ElementType)aGrp->GetType();
2392 if ( myType == aGrpElType || (myType == SMDSAbs_All && aGrpElType != SMDSAbs_Node) ) {
2393 // add elements IDS into control
2394 int aSize = aGrp->Extent();
2395 for (int i = 0; i < aSize; i++)
2396 myIDs.insert( aGrp->GetID(i+1) );
2401 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
2403 TCollection_AsciiString aStr = theStr;
2404 aStr.RemoveAll( ' ' );
2405 aStr.RemoveAll( '\t' );
2406 for ( int aPos = aStr.Search( ";;" ); aPos != -1; aPos = aStr.Search( ";;" ) )
2407 aStr.Remove( aPos, 2 );
2408 Standard_Real clr[3];
2409 clr[0] = clr[1] = clr[2] = 0.;
2410 for ( int i = 0; i < 3; i++ ) {
2411 TCollection_AsciiString tmpStr = aStr.Token( ";", i+1 );
2412 if ( !tmpStr.IsEmpty() && tmpStr.IsRealValue() )
2413 clr[i] = tmpStr.RealValue();
2415 myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
2418 //=======================================================================
2419 // name : GetRangeStr
2420 // Purpose : Get range as a string.
2421 // Example: "1,2,3,50-60,63,67,70-"
2422 //=======================================================================
2423 void GroupColor::GetColorStr( TCollection_AsciiString& theResStr ) const
2426 theResStr += TCollection_AsciiString( myColor.Red() );
2427 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2428 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2432 Class : ElemGeomType
2433 Description : Predicate to check element geometry type
2436 ElemGeomType::ElemGeomType()
2439 myType = SMDSAbs_All;
2440 myGeomType = SMDSGeom_TRIANGLE;
2443 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2448 bool ElemGeomType::IsSatisfy( long theId )
2450 if (!myMesh) return false;
2451 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2454 const SMDSAbs_ElementType anElemType = anElem->GetType();
2455 if ( myType != SMDSAbs_All && anElemType != myType )
2457 const int aNbNode = anElem->NbNodes();
2459 switch( anElemType )
2462 isOk = (myGeomType == SMDSGeom_POINT);
2466 isOk = (myGeomType == SMDSGeom_EDGE);
2470 if ( myGeomType == SMDSGeom_TRIANGLE )
2471 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 6 : aNbNode == 3));
2472 else if ( myGeomType == SMDSGeom_QUADRANGLE )
2473 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 8 : aNbNode == 4));
2474 else if ( myGeomType == SMDSGeom_POLYGON )
2475 isOk = anElem->IsPoly();
2478 case SMDSAbs_Volume:
2479 if ( myGeomType == SMDSGeom_TETRA )
2480 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 10 : aNbNode == 4));
2481 else if ( myGeomType == SMDSGeom_PYRAMID )
2482 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 13 : aNbNode == 5));
2483 else if ( myGeomType == SMDSGeom_PENTA )
2484 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 15 : aNbNode == 6));
2485 else if ( myGeomType == SMDSGeom_HEXA )
2486 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 20 : aNbNode == 8));
2487 else if ( myGeomType == SMDSGeom_POLYHEDRA )
2488 isOk = anElem->IsPoly();
2495 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2500 SMDSAbs_ElementType ElemGeomType::GetType() const
2505 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2507 myGeomType = theType;
2510 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2515 //================================================================================
2517 * \brief Class CoplanarFaces
2519 //================================================================================
2521 CoplanarFaces::CoplanarFaces()
2522 : myMesh(0), myFaceID(0), myToler(0)
2525 bool CoplanarFaces::IsSatisfy( long theElementId )
2527 if ( myCoplanarIDs.empty() )
2529 // Build a set of coplanar face ids
2531 if ( !myMesh || !myFaceID || !myToler )
2534 const SMDS_MeshElement* face = myMesh->FindElement( myFaceID );
2535 if ( !face || face->GetType() != SMDSAbs_Face )
2539 gp_Vec myNorm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2543 const double radianTol = myToler * PI180;
2544 typedef SMDS_StdIterator< const SMDS_MeshElement*, SMDS_ElemIteratorPtr > TFaceIt;
2545 std::set<const SMDS_MeshElement*> checkedFaces, checkedNodes;
2546 std::list<const SMDS_MeshElement*> faceQueue( 1, face );
2547 while ( !faceQueue.empty() )
2549 face = faceQueue.front();
2550 if ( checkedFaces.insert( face ).second )
2552 gp_Vec norm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2553 if (!normOK || myNorm.Angle( norm ) <= radianTol)
2555 myCoplanarIDs.insert( face->GetID() );
2556 std::set<const SMDS_MeshElement*> neighborFaces;
2557 for ( int i = 0; i < face->NbCornerNodes(); ++i )
2559 const SMDS_MeshNode* n = face->GetNode( i );
2560 if ( checkedNodes.insert( n ).second )
2561 neighborFaces.insert( TFaceIt( n->GetInverseElementIterator(SMDSAbs_Face)),
2564 faceQueue.insert( faceQueue.end(), neighborFaces.begin(), neighborFaces.end() );
2567 faceQueue.pop_front();
2570 return myCoplanarIDs.count( theElementId );
2575 *Description : Predicate for Range of Ids.
2576 * Range may be specified with two ways.
2577 * 1. Using AddToRange method
2578 * 2. With SetRangeStr method. Parameter of this method is a string
2579 * like as "1,2,3,50-60,63,67,70-"
2582 //=======================================================================
2583 // name : RangeOfIds
2584 // Purpose : Constructor
2585 //=======================================================================
2586 RangeOfIds::RangeOfIds()
2589 myType = SMDSAbs_All;
2592 //=======================================================================
2594 // Purpose : Set mesh
2595 //=======================================================================
2596 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
2601 //=======================================================================
2602 // name : AddToRange
2603 // Purpose : Add ID to the range
2604 //=======================================================================
2605 bool RangeOfIds::AddToRange( long theEntityId )
2607 myIds.Add( theEntityId );
2611 //=======================================================================
2612 // name : GetRangeStr
2613 // Purpose : Get range as a string.
2614 // Example: "1,2,3,50-60,63,67,70-"
2615 //=======================================================================
2616 void RangeOfIds::GetRangeStr( TCollection_AsciiString& theResStr )
2620 TColStd_SequenceOfInteger anIntSeq;
2621 TColStd_SequenceOfAsciiString aStrSeq;
2623 TColStd_MapIteratorOfMapOfInteger anIter( myIds );
2624 for ( ; anIter.More(); anIter.Next() )
2626 int anId = anIter.Key();
2627 TCollection_AsciiString aStr( anId );
2628 anIntSeq.Append( anId );
2629 aStrSeq.Append( aStr );
2632 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2634 int aMinId = myMin( i );
2635 int aMaxId = myMax( i );
2637 TCollection_AsciiString aStr;
2638 if ( aMinId != IntegerFirst() )
2643 if ( aMaxId != IntegerLast() )
2646 // find position of the string in result sequence and insert string in it
2647 if ( anIntSeq.Length() == 0 )
2649 anIntSeq.Append( aMinId );
2650 aStrSeq.Append( aStr );
2654 if ( aMinId < anIntSeq.First() )
2656 anIntSeq.Prepend( aMinId );
2657 aStrSeq.Prepend( aStr );
2659 else if ( aMinId > anIntSeq.Last() )
2661 anIntSeq.Append( aMinId );
2662 aStrSeq.Append( aStr );
2665 for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
2666 if ( aMinId < anIntSeq( j ) )
2668 anIntSeq.InsertBefore( j, aMinId );
2669 aStrSeq.InsertBefore( j, aStr );
2675 if ( aStrSeq.Length() == 0 )
2678 theResStr = aStrSeq( 1 );
2679 for ( int j = 2, k = aStrSeq.Length(); j <= k; j++ )
2682 theResStr += aStrSeq( j );
2686 //=======================================================================
2687 // name : SetRangeStr
2688 // Purpose : Define range with string
2689 // Example of entry string: "1,2,3,50-60,63,67,70-"
2690 //=======================================================================
2691 bool RangeOfIds::SetRangeStr( const TCollection_AsciiString& theStr )
2697 TCollection_AsciiString aStr = theStr;
2698 aStr.RemoveAll( ' ' );
2699 aStr.RemoveAll( '\t' );
2701 for ( int aPos = aStr.Search( ",," ); aPos != -1; aPos = aStr.Search( ",," ) )
2702 aStr.Remove( aPos, 2 );
2704 TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
2706 while ( tmpStr != "" )
2708 tmpStr = aStr.Token( ",", i++ );
2709 int aPos = tmpStr.Search( '-' );
2713 if ( tmpStr.IsIntegerValue() )
2714 myIds.Add( tmpStr.IntegerValue() );
2720 TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
2721 TCollection_AsciiString aMinStr = tmpStr;
2723 while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
2724 while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
2726 if ( (!aMinStr.IsEmpty() && !aMinStr.IsIntegerValue()) ||
2727 (!aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue()) )
2730 myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
2731 myMax.Append( aMaxStr.IsEmpty() ? IntegerLast() : aMaxStr.IntegerValue() );
2738 //=======================================================================
2740 // Purpose : Get type of supported entities
2741 //=======================================================================
2742 SMDSAbs_ElementType RangeOfIds::GetType() const
2747 //=======================================================================
2749 // Purpose : Set type of supported entities
2750 //=======================================================================
2751 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
2756 //=======================================================================
2758 // Purpose : Verify whether entity satisfies to this rpedicate
2759 //=======================================================================
2760 bool RangeOfIds::IsSatisfy( long theId )
2765 if ( myType == SMDSAbs_Node )
2767 if ( myMesh->FindNode( theId ) == 0 )
2772 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2773 if ( anElem == 0 || (myType != anElem->GetType() && myType != SMDSAbs_All ))
2777 if ( myIds.Contains( theId ) )
2780 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2781 if ( theId >= myMin( i ) && theId <= myMax( i ) )
2789 Description : Base class for comparators
2791 Comparator::Comparator():
2795 Comparator::~Comparator()
2798 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
2801 myFunctor->SetMesh( theMesh );
2804 void Comparator::SetMargin( double theValue )
2806 myMargin = theValue;
2809 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
2811 myFunctor = theFunct;
2814 SMDSAbs_ElementType Comparator::GetType() const
2816 return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
2819 double Comparator::GetMargin()
2827 Description : Comparator "<"
2829 bool LessThan::IsSatisfy( long theId )
2831 return myFunctor && myFunctor->GetValue( theId ) < myMargin;
2837 Description : Comparator ">"
2839 bool MoreThan::IsSatisfy( long theId )
2841 return myFunctor && myFunctor->GetValue( theId ) > myMargin;
2847 Description : Comparator "="
2850 myToler(Precision::Confusion())
2853 bool EqualTo::IsSatisfy( long theId )
2855 return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
2858 void EqualTo::SetTolerance( double theToler )
2863 double EqualTo::GetTolerance()
2870 Description : Logical NOT predicate
2872 LogicalNOT::LogicalNOT()
2875 LogicalNOT::~LogicalNOT()
2878 bool LogicalNOT::IsSatisfy( long theId )
2880 return myPredicate && !myPredicate->IsSatisfy( theId );
2883 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
2886 myPredicate->SetMesh( theMesh );
2889 void LogicalNOT::SetPredicate( PredicatePtr thePred )
2891 myPredicate = thePred;
2894 SMDSAbs_ElementType LogicalNOT::GetType() const
2896 return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
2901 Class : LogicalBinary
2902 Description : Base class for binary logical predicate
2904 LogicalBinary::LogicalBinary()
2907 LogicalBinary::~LogicalBinary()
2910 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
2913 myPredicate1->SetMesh( theMesh );
2916 myPredicate2->SetMesh( theMesh );
2919 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
2921 myPredicate1 = thePredicate;
2924 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
2926 myPredicate2 = thePredicate;
2929 SMDSAbs_ElementType LogicalBinary::GetType() const
2931 if ( !myPredicate1 || !myPredicate2 )
2934 SMDSAbs_ElementType aType1 = myPredicate1->GetType();
2935 SMDSAbs_ElementType aType2 = myPredicate2->GetType();
2937 return aType1 == aType2 ? aType1 : SMDSAbs_All;
2943 Description : Logical AND
2945 bool LogicalAND::IsSatisfy( long theId )
2950 myPredicate1->IsSatisfy( theId ) &&
2951 myPredicate2->IsSatisfy( theId );
2957 Description : Logical OR
2959 bool LogicalOR::IsSatisfy( long theId )
2964 (myPredicate1->IsSatisfy( theId ) ||
2965 myPredicate2->IsSatisfy( theId ));
2979 void Filter::SetPredicate( PredicatePtr thePredicate )
2981 myPredicate = thePredicate;
2984 template<class TElement, class TIterator, class TPredicate>
2985 inline void FillSequence(const TIterator& theIterator,
2986 TPredicate& thePredicate,
2987 Filter::TIdSequence& theSequence)
2989 if ( theIterator ) {
2990 while( theIterator->more() ) {
2991 TElement anElem = theIterator->next();
2992 long anId = anElem->GetID();
2993 if ( thePredicate->IsSatisfy( anId ) )
2994 theSequence.push_back( anId );
3001 GetElementsId( const SMDS_Mesh* theMesh,
3002 PredicatePtr thePredicate,
3003 TIdSequence& theSequence )
3005 theSequence.clear();
3007 if ( !theMesh || !thePredicate )
3010 thePredicate->SetMesh( theMesh );
3012 SMDSAbs_ElementType aType = thePredicate->GetType();
3015 FillSequence<const SMDS_MeshNode*>(theMesh->nodesIterator(),thePredicate,theSequence);
3018 FillSequence<const SMDS_MeshElement*>(theMesh->edgesIterator(),thePredicate,theSequence);
3021 FillSequence<const SMDS_MeshElement*>(theMesh->facesIterator(),thePredicate,theSequence);
3023 case SMDSAbs_Volume:
3024 FillSequence<const SMDS_MeshElement*>(theMesh->volumesIterator(),thePredicate,theSequence);
3027 FillSequence<const SMDS_MeshElement*>(theMesh->edgesIterator(),thePredicate,theSequence);
3028 FillSequence<const SMDS_MeshElement*>(theMesh->facesIterator(),thePredicate,theSequence);
3029 FillSequence<const SMDS_MeshElement*>(theMesh->volumesIterator(),thePredicate,theSequence);
3035 Filter::GetElementsId( const SMDS_Mesh* theMesh,
3036 Filter::TIdSequence& theSequence )
3038 GetElementsId(theMesh,myPredicate,theSequence);
3045 typedef std::set<SMDS_MeshFace*> TMapOfFacePtr;
3051 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
3052 SMDS_MeshNode* theNode2 )
3058 ManifoldPart::Link::~Link()
3064 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
3066 if ( myNode1 == theLink.myNode1 &&
3067 myNode2 == theLink.myNode2 )
3069 else if ( myNode1 == theLink.myNode2 &&
3070 myNode2 == theLink.myNode1 )
3076 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
3078 if(myNode1 < x.myNode1) return true;
3079 if(myNode1 == x.myNode1)
3080 if(myNode2 < x.myNode2) return true;
3084 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
3085 const ManifoldPart::Link& theLink2 )
3087 return theLink1.IsEqual( theLink2 );
3090 ManifoldPart::ManifoldPart()
3093 myAngToler = Precision::Angular();
3094 myIsOnlyManifold = true;
3097 ManifoldPart::~ManifoldPart()
3102 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
3108 SMDSAbs_ElementType ManifoldPart::GetType() const
3109 { return SMDSAbs_Face; }
3111 bool ManifoldPart::IsSatisfy( long theElementId )
3113 return myMapIds.Contains( theElementId );
3116 void ManifoldPart::SetAngleTolerance( const double theAngToler )
3117 { myAngToler = theAngToler; }
3119 double ManifoldPart::GetAngleTolerance() const
3120 { return myAngToler; }
3122 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
3123 { myIsOnlyManifold = theIsOnly; }
3125 void ManifoldPart::SetStartElem( const long theStartId )
3126 { myStartElemId = theStartId; }
3128 bool ManifoldPart::process()
3131 myMapBadGeomIds.Clear();
3133 myAllFacePtr.clear();
3134 myAllFacePtrIntDMap.clear();
3138 // collect all faces into own map
3139 SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3140 for (; anFaceItr->more(); )
3142 SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3143 myAllFacePtr.push_back( aFacePtr );
3144 myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3147 SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3151 // the map of non manifold links and bad geometry
3152 TMapOfLink aMapOfNonManifold;
3153 TColStd_MapOfInteger aMapOfTreated;
3155 // begin cycle on faces from start index and run on vector till the end
3156 // and from begin to start index to cover whole vector
3157 const int aStartIndx = myAllFacePtrIntDMap[aStartFace];
3158 bool isStartTreat = false;
3159 for ( int fi = aStartIndx; !isStartTreat || fi != aStartIndx ; fi++ )
3161 if ( fi == aStartIndx )
3162 isStartTreat = true;
3163 // as result next time when fi will be equal to aStartIndx
3165 SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3166 if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3169 aMapOfTreated.Add( aFacePtr->GetID() );
3170 TColStd_MapOfInteger aResFaces;
3171 if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3172 aMapOfNonManifold, aResFaces ) )
3174 TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3175 for ( ; anItr.More(); anItr.Next() )
3177 int aFaceId = anItr.Key();
3178 aMapOfTreated.Add( aFaceId );
3179 myMapIds.Add( aFaceId );
3182 if ( fi == ( myAllFacePtr.size() - 1 ) )
3184 } // end run on vector of faces
3185 return !myMapIds.IsEmpty();
3188 static void getLinks( const SMDS_MeshFace* theFace,
3189 ManifoldPart::TVectorOfLink& theLinks )
3191 int aNbNode = theFace->NbNodes();
3192 SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3194 SMDS_MeshNode* aNode = 0;
3195 for ( ; aNodeItr->more() && i <= aNbNode; )
3198 SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3202 SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3204 ManifoldPart::Link aLink( aN1, aN2 );
3205 theLinks.push_back( aLink );
3209 bool ManifoldPart::findConnected
3210 ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3211 SMDS_MeshFace* theStartFace,
3212 ManifoldPart::TMapOfLink& theNonManifold,
3213 TColStd_MapOfInteger& theResFaces )
3215 theResFaces.Clear();
3216 if ( !theAllFacePtrInt.size() )
3219 if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3221 myMapBadGeomIds.Add( theStartFace->GetID() );
3225 ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3226 ManifoldPart::TVectorOfLink aSeqOfBoundary;
3227 theResFaces.Add( theStartFace->GetID() );
3228 ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3230 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3231 aDMapLinkFace, theNonManifold, theStartFace );
3233 bool isDone = false;
3234 while ( !isDone && aMapOfBoundary.size() != 0 )
3236 bool isToReset = false;
3237 ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3238 for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3240 ManifoldPart::Link aLink = *pLink;
3241 if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3243 // each link could be treated only once
3244 aMapToSkip.insert( aLink );
3246 ManifoldPart::TVectorOfFacePtr aFaces;
3248 if ( myIsOnlyManifold &&
3249 (theNonManifold.find( aLink ) != theNonManifold.end()) )
3253 getFacesByLink( aLink, aFaces );
3254 // filter the element to keep only indicated elements
3255 ManifoldPart::TVectorOfFacePtr aFiltered;
3256 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3257 for ( ; pFace != aFaces.end(); ++pFace )
3259 SMDS_MeshFace* aFace = *pFace;
3260 if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3261 aFiltered.push_back( aFace );
3264 if ( aFaces.size() < 2 ) // no neihgbour faces
3266 else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3268 theNonManifold.insert( aLink );
3273 // compare normal with normals of neighbor element
3274 SMDS_MeshFace* aPrevFace = aDMapLinkFace[ aLink ];
3275 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3276 for ( ; pFace != aFaces.end(); ++pFace )
3278 SMDS_MeshFace* aNextFace = *pFace;
3279 if ( aPrevFace == aNextFace )
3281 int anNextFaceID = aNextFace->GetID();
3282 if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3283 // should not be with non manifold restriction. probably bad topology
3285 // check if face was treated and skipped
3286 if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3287 !isInPlane( aPrevFace, aNextFace ) )
3289 // add new element to connected and extend the boundaries.
3290 theResFaces.Add( anNextFaceID );
3291 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3292 aDMapLinkFace, theNonManifold, aNextFace );
3296 isDone = !isToReset;
3299 return !theResFaces.IsEmpty();
3302 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3303 const SMDS_MeshFace* theFace2 )
3305 gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3306 gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3307 if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3309 myMapBadGeomIds.Add( theFace2->GetID() );
3312 if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
3318 void ManifoldPart::expandBoundary
3319 ( ManifoldPart::TMapOfLink& theMapOfBoundary,
3320 ManifoldPart::TVectorOfLink& theSeqOfBoundary,
3321 ManifoldPart::TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
3322 ManifoldPart::TMapOfLink& theNonManifold,
3323 SMDS_MeshFace* theNextFace ) const
3325 ManifoldPart::TVectorOfLink aLinks;
3326 getLinks( theNextFace, aLinks );
3327 int aNbLink = (int)aLinks.size();
3328 for ( int i = 0; i < aNbLink; i++ )
3330 ManifoldPart::Link aLink = aLinks[ i ];
3331 if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3333 if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3335 if ( myIsOnlyManifold )
3337 // remove from boundary
3338 theMapOfBoundary.erase( aLink );
3339 ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3340 for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3342 ManifoldPart::Link aBoundLink = *pLink;
3343 if ( aBoundLink.IsEqual( aLink ) )
3345 theSeqOfBoundary.erase( pLink );
3353 theMapOfBoundary.insert( aLink );
3354 theSeqOfBoundary.push_back( aLink );
3355 theDMapLinkFacePtr[ aLink ] = theNextFace;
3360 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3361 ManifoldPart::TVectorOfFacePtr& theFaces ) const
3363 std::set<SMDS_MeshCell *> aSetOfFaces;
3364 // take all faces that shared first node
3365 SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3366 for ( ; anItr->more(); )
3368 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3371 aSetOfFaces.insert( aFace );
3373 // take all faces that shared second node
3374 anItr = theLink.myNode2->facesIterator();
3375 // find the common part of two sets
3376 for ( ; anItr->more(); )
3378 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3379 if ( aSetOfFaces.count( aFace ) )
3380 theFaces.push_back( aFace );
3389 ElementsOnSurface::ElementsOnSurface()
3393 myType = SMDSAbs_All;
3395 myToler = Precision::Confusion();
3396 myUseBoundaries = false;
3399 ElementsOnSurface::~ElementsOnSurface()
3404 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
3406 if ( myMesh == theMesh )
3412 bool ElementsOnSurface::IsSatisfy( long theElementId )
3414 return myIds.Contains( theElementId );
3417 SMDSAbs_ElementType ElementsOnSurface::GetType() const
3420 void ElementsOnSurface::SetTolerance( const double theToler )
3422 if ( myToler != theToler )
3427 double ElementsOnSurface::GetTolerance() const
3430 void ElementsOnSurface::SetUseBoundaries( bool theUse )
3432 if ( myUseBoundaries != theUse ) {
3433 myUseBoundaries = theUse;
3434 SetSurface( mySurf, myType );
3438 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
3439 const SMDSAbs_ElementType theType )
3444 if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
3446 mySurf = TopoDS::Face( theShape );
3447 BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
3449 u1 = SA.FirstUParameter(),
3450 u2 = SA.LastUParameter(),
3451 v1 = SA.FirstVParameter(),
3452 v2 = SA.LastVParameter();
3453 Handle(Geom_Surface) surf = BRep_Tool::Surface( mySurf );
3454 myProjector.Init( surf, u1,u2, v1,v2 );
3458 void ElementsOnSurface::process()
3461 if ( mySurf.IsNull() )
3467 if ( myType == SMDSAbs_Face || myType == SMDSAbs_All )
3469 myIds.ReSize( myMesh->NbFaces() );
3470 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
3471 for(; anIter->more(); )
3472 process( anIter->next() );
3475 if ( myType == SMDSAbs_Edge || myType == SMDSAbs_All )
3477 myIds.ReSize( myIds.Extent() + myMesh->NbEdges() );
3478 SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
3479 for(; anIter->more(); )
3480 process( anIter->next() );
3483 if ( myType == SMDSAbs_Node )
3485 myIds.ReSize( myMesh->NbNodes() );
3486 SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
3487 for(; anIter->more(); )
3488 process( anIter->next() );
3492 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
3494 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3495 bool isSatisfy = true;
3496 for ( ; aNodeItr->more(); )
3498 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3499 if ( !isOnSurface( aNode ) )
3506 myIds.Add( theElemPtr->GetID() );
3509 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
3511 if ( mySurf.IsNull() )
3514 gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
3515 // double aToler2 = myToler * myToler;
3516 // if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
3518 // gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
3519 // if ( aPln.SquareDistance( aPnt ) > aToler2 )
3522 // else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
3524 // gp_Cylinder aCyl = Handle(Geom_CylindricalSurface)::DownCast(mySurf)->Cylinder();
3525 // double aRad = aCyl.Radius();
3526 // gp_Ax3 anAxis = aCyl.Position();
3527 // gp_XYZ aLoc = aCyl.Location().XYZ();
3528 // double aXDist = anAxis.XDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3529 // double aYDist = anAxis.YDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3530 // if ( fabs(aXDist*aXDist + aYDist*aYDist - aRad*aRad) > aToler2 )
3535 myProjector.Perform( aPnt );
3536 bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
3546 ElementsOnShape::ElementsOnShape()
3548 myType(SMDSAbs_All),
3549 myToler(Precision::Confusion()),
3550 myAllNodesFlag(false)
3552 myCurShapeType = TopAbs_SHAPE;
3555 ElementsOnShape::~ElementsOnShape()
3559 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
3561 if (myMesh != theMesh) {
3563 SetShape(myShape, myType);
3567 bool ElementsOnShape::IsSatisfy (long theElementId)
3569 return myIds.Contains(theElementId);
3572 SMDSAbs_ElementType ElementsOnShape::GetType() const
3577 void ElementsOnShape::SetTolerance (const double theToler)
3579 if (myToler != theToler) {
3581 SetShape(myShape, myType);
3585 double ElementsOnShape::GetTolerance() const
3590 void ElementsOnShape::SetAllNodes (bool theAllNodes)
3592 if (myAllNodesFlag != theAllNodes) {
3593 myAllNodesFlag = theAllNodes;
3594 SetShape(myShape, myType);
3598 void ElementsOnShape::SetShape (const TopoDS_Shape& theShape,
3599 const SMDSAbs_ElementType theType)
3605 if (myMesh == 0) return;
3610 myIds.ReSize(myMesh->NbEdges() + myMesh->NbFaces() + myMesh->NbVolumes());
3613 myIds.ReSize(myMesh->NbNodes());
3616 myIds.ReSize(myMesh->NbEdges());
3619 myIds.ReSize(myMesh->NbFaces());
3621 case SMDSAbs_Volume:
3622 myIds.ReSize(myMesh->NbVolumes());
3628 myShapesMap.Clear();
3632 void ElementsOnShape::addShape (const TopoDS_Shape& theShape)
3634 if (theShape.IsNull() || myMesh == 0)
3637 if (!myShapesMap.Add(theShape)) return;
3639 myCurShapeType = theShape.ShapeType();
3640 switch (myCurShapeType)
3642 case TopAbs_COMPOUND:
3643 case TopAbs_COMPSOLID:
3647 TopoDS_Iterator anIt (theShape, Standard_True, Standard_True);
3648 for (; anIt.More(); anIt.Next()) addShape(anIt.Value());
3653 myCurSC.Load(theShape);
3659 TopoDS_Face aFace = TopoDS::Face(theShape);
3660 BRepAdaptor_Surface SA (aFace, true);
3662 u1 = SA.FirstUParameter(),
3663 u2 = SA.LastUParameter(),
3664 v1 = SA.FirstVParameter(),
3665 v2 = SA.LastVParameter();
3666 Handle(Geom_Surface) surf = BRep_Tool::Surface(aFace);
3667 myCurProjFace.Init(surf, u1,u2, v1,v2);
3674 TopoDS_Edge anEdge = TopoDS::Edge(theShape);
3675 Standard_Real u1, u2;
3676 Handle(Geom_Curve) curve = BRep_Tool::Curve(anEdge, u1, u2);
3677 myCurProjEdge.Init(curve, u1, u2);
3683 TopoDS_Vertex aV = TopoDS::Vertex(theShape);
3684 myCurPnt = BRep_Tool::Pnt(aV);
3693 void ElementsOnShape::process()
3695 if (myShape.IsNull() || myMesh == 0)
3698 if (myType == SMDSAbs_Node)
3700 SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
3701 while (anIter->more())
3702 process(anIter->next());
3706 if (myType == SMDSAbs_Edge || myType == SMDSAbs_All)
3708 SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
3709 while (anIter->more())
3710 process(anIter->next());
3713 if (myType == SMDSAbs_Face || myType == SMDSAbs_All)
3715 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
3716 while (anIter->more()) {
3717 process(anIter->next());
3721 if (myType == SMDSAbs_Volume || myType == SMDSAbs_All)
3723 SMDS_VolumeIteratorPtr anIter = myMesh->volumesIterator();
3724 while (anIter->more())
3725 process(anIter->next());
3730 void ElementsOnShape::process (const SMDS_MeshElement* theElemPtr)
3732 if (myShape.IsNull())
3735 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3736 bool isSatisfy = myAllNodesFlag;
3738 gp_XYZ centerXYZ (0, 0, 0);
3740 while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
3742 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3743 gp_Pnt aPnt (aNode->X(), aNode->Y(), aNode->Z());
3744 centerXYZ += aPnt.XYZ();
3746 switch (myCurShapeType)
3750 myCurSC.Perform(aPnt, myToler);
3751 isSatisfy = (myCurSC.State() == TopAbs_IN || myCurSC.State() == TopAbs_ON);
3756 myCurProjFace.Perform(aPnt);
3757 isSatisfy = (myCurProjFace.IsDone() && myCurProjFace.LowerDistance() <= myToler);
3760 // check relatively the face
3761 Quantity_Parameter u, v;
3762 myCurProjFace.LowerDistanceParameters(u, v);
3763 gp_Pnt2d aProjPnt (u, v);
3764 BRepClass_FaceClassifier aClsf (myCurFace, aProjPnt, myToler);
3765 isSatisfy = (aClsf.State() == TopAbs_IN || aClsf.State() == TopAbs_ON);
3771 myCurProjEdge.Perform(aPnt);
3772 isSatisfy = (myCurProjEdge.NbPoints() > 0 && myCurProjEdge.LowerDistance() <= myToler);
3777 isSatisfy = (aPnt.Distance(myCurPnt) <= myToler);
3787 if (isSatisfy && myCurShapeType == TopAbs_SOLID) { // Check the center point for volumes MantisBug 0020168
3788 centerXYZ /= theElemPtr->NbNodes();
3789 gp_Pnt aCenterPnt (centerXYZ);
3790 myCurSC.Perform(aCenterPnt, myToler);
3791 if ( !(myCurSC.State() == TopAbs_IN || myCurSC.State() == TopAbs_ON))
3796 myIds.Add(theElemPtr->GetID());
3799 TSequenceOfXYZ::TSequenceOfXYZ()
3802 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n)
3805 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t)
3808 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray)
3811 template <class InputIterator>
3812 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd)
3815 TSequenceOfXYZ::~TSequenceOfXYZ()
3818 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
3820 myArray = theSequenceOfXYZ.myArray;
3824 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
3826 return myArray[n-1];
3829 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
3831 return myArray[n-1];
3834 void TSequenceOfXYZ::clear()
3839 void TSequenceOfXYZ::reserve(size_type n)
3844 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
3846 myArray.push_back(v);
3849 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
3851 return myArray.size();