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"
27 #include <BRepAdaptor_Surface.hxx>
28 #include <BRepClass_FaceClassifier.hxx>
29 #include <BRep_Tool.hxx>
33 #include <TopoDS_Edge.hxx>
34 #include <TopoDS_Face.hxx>
35 #include <TopoDS_Shape.hxx>
36 #include <TopoDS_Vertex.hxx>
37 #include <TopoDS_Iterator.hxx>
39 #include <Geom_CylindricalSurface.hxx>
40 #include <Geom_Plane.hxx>
41 #include <Geom_Surface.hxx>
43 #include <Precision.hxx>
44 #include <TColStd_MapIteratorOfMapOfInteger.hxx>
45 #include <TColStd_MapOfInteger.hxx>
46 #include <TColStd_SequenceOfAsciiString.hxx>
47 #include <TColgp_Array1OfXYZ.hxx>
50 #include <gp_Cylinder.hxx>
57 #include "SMDS_Mesh.hxx"
58 #include "SMDS_Iterator.hxx"
59 #include "SMDS_MeshElement.hxx"
60 #include "SMDS_MeshNode.hxx"
61 #include "SMDS_VolumeTool.hxx"
62 #include "SMDS_QuadraticFaceOfNodes.hxx"
63 #include "SMDS_QuadraticEdge.hxx"
65 #include "SMESHDS_Mesh.hxx"
66 #include "SMESHDS_GroupBase.hxx"
73 inline double getAngle( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
75 gp_Vec v1( P1 - P2 ), v2( P3 - P2 );
77 return v1.Magnitude() < gp::Resolution() ||
78 v2.Magnitude() < gp::Resolution() ? 0 : v1.Angle( v2 );
81 inline double getArea( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
83 gp_Vec aVec1( P2 - P1 );
84 gp_Vec aVec2( P3 - P1 );
85 return ( aVec1 ^ aVec2 ).Magnitude() * 0.5;
88 inline double getArea( const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3 )
90 return getArea( P1.XYZ(), P2.XYZ(), P3.XYZ() );
95 inline double getDistance( const gp_XYZ& P1, const gp_XYZ& P2 )
97 double aDist = gp_Pnt( P1 ).Distance( gp_Pnt( P2 ) );
101 int getNbMultiConnection( const SMDS_Mesh* theMesh, const int theId )
106 const SMDS_MeshElement* anEdge = theMesh->FindElement( theId );
107 if ( anEdge == 0 || anEdge->GetType() != SMDSAbs_Edge/* || anEdge->NbNodes() != 2 */)
110 // for each pair of nodes in anEdge (there are 2 pairs in a quadratic edge)
111 // count elements containing both nodes of the pair.
112 // Note that there may be such cases for a quadratic edge (a horizontal line):
117 // +-----+------+ +-----+------+
120 // result sould be 2 in both cases
122 int aResult0 = 0, aResult1 = 0;
123 // last node, it is a medium one in a quadratic edge
124 const SMDS_MeshNode* aLastNode = anEdge->GetNode( anEdge->NbNodes() - 1 );
125 const SMDS_MeshNode* aNode0 = anEdge->GetNode( 0 );
126 const SMDS_MeshNode* aNode1 = anEdge->GetNode( 1 );
127 if ( aNode1 == aLastNode ) aNode1 = 0;
129 SMDS_ElemIteratorPtr anElemIter = aLastNode->GetInverseElementIterator();
130 while( anElemIter->more() ) {
131 const SMDS_MeshElement* anElem = anElemIter->next();
132 if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
133 SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
134 while ( anIter->more() ) {
135 if ( const SMDS_MeshElement* anElemNode = anIter->next() ) {
136 if ( anElemNode == aNode0 ) {
138 if ( !aNode1 ) break; // not a quadratic edge
140 else if ( anElemNode == aNode1 )
146 int aResult = std::max ( aResult0, aResult1 );
148 // TColStd_MapOfInteger aMap;
150 // SMDS_ElemIteratorPtr anIter = anEdge->nodesIterator();
151 // if ( anIter != 0 ) {
152 // while( anIter->more() ) {
153 // const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
156 // SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
157 // while( anElemIter->more() ) {
158 // const SMDS_MeshElement* anElem = anElemIter->next();
159 // if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
160 // int anId = anElem->GetID();
162 // if ( anIter->more() ) // i.e. first node
164 // else if ( aMap.Contains( anId ) )
178 using namespace SMESH::Controls;
185 Class : NumericalFunctor
186 Description : Base class for numerical functors
188 NumericalFunctor::NumericalFunctor():
194 void NumericalFunctor::SetMesh( const SMDS_Mesh* theMesh )
199 bool NumericalFunctor::GetPoints(const int theId,
200 TSequenceOfXYZ& theRes ) const
207 return GetPoints( myMesh->FindElement( theId ), theRes );
210 bool NumericalFunctor::GetPoints(const SMDS_MeshElement* anElem,
211 TSequenceOfXYZ& theRes )
218 theRes.reserve( anElem->NbNodes() );
220 // Get nodes of the element
221 SMDS_ElemIteratorPtr anIter;
223 if ( anElem->IsQuadratic() ) {
224 switch ( anElem->GetType() ) {
226 anIter = static_cast<const SMDS_QuadraticEdge*>
227 (anElem)->interlacedNodesElemIterator();
230 anIter = static_cast<const SMDS_QuadraticFaceOfNodes*>
231 (anElem)->interlacedNodesElemIterator();
234 anIter = anElem->nodesIterator();
239 anIter = anElem->nodesIterator();
243 while( anIter->more() ) {
244 if ( const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>( anIter->next() ))
245 theRes.push_back( gp_XYZ( aNode->X(), aNode->Y(), aNode->Z() ) );
252 long NumericalFunctor::GetPrecision() const
257 void NumericalFunctor::SetPrecision( const long thePrecision )
259 myPrecision = thePrecision;
262 double NumericalFunctor::GetValue( long theId )
264 myCurrElement = myMesh->FindElement( theId );
266 if ( GetPoints( theId, P ))
268 double aVal = GetValue( P );
269 if ( myPrecision >= 0 )
271 double prec = pow( 10., (double)( myPrecision ) );
272 aVal = floor( aVal * prec + 0.5 ) / prec;
280 //================================================================================
282 * \brief Return histogram of functor values
283 * \param nbIntervals - number of intervals
284 * \param nbEvents - number of mesh elements having values within i-th interval
285 * \param funValues - boundaries of intervals
286 * \param elements - elements to check vulue of; empty list means "of all"
287 * \param minmax - boundaries of diapason of values to divide into intervals
289 //================================================================================
291 void NumericalFunctor::GetHistogram(int nbIntervals,
292 std::vector<int>& nbEvents,
293 std::vector<double>& funValues,
294 const vector<int>& elements,
295 const double* minmax)
297 if ( nbIntervals < 1 ||
299 !myMesh->GetMeshInfo().NbElements( GetType() ))
301 nbEvents.resize( nbIntervals, 0 );
302 funValues.resize( nbIntervals+1 );
304 // get all values sorted
305 std::multiset< double > values;
306 if ( elements.empty() )
308 SMDS_ElemIteratorPtr elemIt = myMesh->elementsIterator(GetType());
309 while ( elemIt->more() )
310 values.insert( GetValue( elemIt->next()->GetID() ));
314 vector<int>::const_iterator id = elements.begin();
315 for ( ; id != elements.end(); ++id )
316 values.insert( GetValue( *id ));
321 funValues[0] = minmax[0];
322 funValues[nbIntervals] = minmax[1];
326 funValues[0] = *values.begin();
327 funValues[nbIntervals] = *values.rbegin();
329 // case nbIntervals == 1
330 if ( nbIntervals == 1 )
332 nbEvents[0] = values.size();
336 if (funValues.front() == funValues.back())
338 nbEvents.resize( 1 );
339 nbEvents[0] = values.size();
340 funValues[1] = funValues.back();
341 funValues.resize( 2 );
344 std::multiset< double >::iterator min = values.begin(), max;
345 for ( int i = 0; i < nbIntervals; ++i )
347 // find end value of i-th interval
348 double r = (i+1) / double( nbIntervals );
349 funValues[i+1] = funValues.front() * (1-r) + funValues.back() * r;
351 // count values in the i-th interval if there are any
352 if ( min != values.end() && *min <= funValues[i+1] )
354 // find the first value out of the interval
355 max = values.upper_bound( funValues[i+1] ); // max is greater than funValues[i+1], or end()
356 nbEvents[i] = std::distance( min, max );
360 // add values larger than minmax[1]
361 nbEvents.back() += std::distance( min, values.end() );
364 //=======================================================================
365 //function : GetValue
367 //=======================================================================
369 double Volume::GetValue( long theElementId )
371 if ( theElementId && myMesh ) {
372 SMDS_VolumeTool aVolumeTool;
373 if ( aVolumeTool.Set( myMesh->FindElement( theElementId )))
374 return aVolumeTool.GetSize();
379 //=======================================================================
380 //function : GetBadRate
381 //purpose : meaningless as it is not quality control functor
382 //=======================================================================
384 double Volume::GetBadRate( double Value, int /*nbNodes*/ ) const
389 //=======================================================================
392 //=======================================================================
394 SMDSAbs_ElementType Volume::GetType() const
396 return SMDSAbs_Volume;
401 Class : MaxElementLength2D
402 Description : Functor calculating maximum length of 2D element
405 double MaxElementLength2D::GetValue( long theElementId )
408 if( GetPoints( theElementId, P ) ) {
410 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
411 SMDSAbs_ElementType aType = aElem->GetType();
415 if( len == 3 ) { // triangles
416 double L1 = getDistance(P( 1 ),P( 2 ));
417 double L2 = getDistance(P( 2 ),P( 3 ));
418 double L3 = getDistance(P( 3 ),P( 1 ));
419 aVal = Max(L1,Max(L2,L3));
422 else if( len == 4 ) { // quadrangles
423 double L1 = getDistance(P( 1 ),P( 2 ));
424 double L2 = getDistance(P( 2 ),P( 3 ));
425 double L3 = getDistance(P( 3 ),P( 4 ));
426 double L4 = getDistance(P( 4 ),P( 1 ));
427 double D1 = getDistance(P( 1 ),P( 3 ));
428 double D2 = getDistance(P( 2 ),P( 4 ));
429 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
432 else if( len == 6 ) { // quadratic triangles
433 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
434 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
435 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
436 aVal = Max(L1,Max(L2,L3));
439 else if( len == 8 ) { // quadratic quadrangles
440 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
441 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
442 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
443 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
444 double D1 = getDistance(P( 1 ),P( 5 ));
445 double D2 = getDistance(P( 3 ),P( 7 ));
446 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
451 if( myPrecision >= 0 )
453 double prec = pow( 10., (double)myPrecision );
454 aVal = floor( aVal * prec + 0.5 ) / prec;
461 double MaxElementLength2D::GetBadRate( double Value, int /*nbNodes*/ ) const
466 SMDSAbs_ElementType MaxElementLength2D::GetType() const
472 Class : MaxElementLength3D
473 Description : Functor calculating maximum length of 3D element
476 double MaxElementLength3D::GetValue( long theElementId )
479 if( GetPoints( theElementId, P ) ) {
481 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
482 SMDSAbs_ElementType aType = aElem->GetType();
486 if( len == 4 ) { // tetras
487 double L1 = getDistance(P( 1 ),P( 2 ));
488 double L2 = getDistance(P( 2 ),P( 3 ));
489 double L3 = getDistance(P( 3 ),P( 1 ));
490 double L4 = getDistance(P( 1 ),P( 4 ));
491 double L5 = getDistance(P( 2 ),P( 4 ));
492 double L6 = getDistance(P( 3 ),P( 4 ));
493 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
496 else if( len == 5 ) { // pyramids
497 double L1 = getDistance(P( 1 ),P( 2 ));
498 double L2 = getDistance(P( 2 ),P( 3 ));
499 double L3 = getDistance(P( 3 ),P( 4 ));
500 double L4 = getDistance(P( 4 ),P( 1 ));
501 double L5 = getDistance(P( 1 ),P( 5 ));
502 double L6 = getDistance(P( 2 ),P( 5 ));
503 double L7 = getDistance(P( 3 ),P( 5 ));
504 double L8 = getDistance(P( 4 ),P( 5 ));
505 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
506 aVal = Max(aVal,Max(L7,L8));
509 else if( len == 6 ) { // pentas
510 double L1 = getDistance(P( 1 ),P( 2 ));
511 double L2 = getDistance(P( 2 ),P( 3 ));
512 double L3 = getDistance(P( 3 ),P( 1 ));
513 double L4 = getDistance(P( 4 ),P( 5 ));
514 double L5 = getDistance(P( 5 ),P( 6 ));
515 double L6 = getDistance(P( 6 ),P( 4 ));
516 double L7 = getDistance(P( 1 ),P( 4 ));
517 double L8 = getDistance(P( 2 ),P( 5 ));
518 double L9 = getDistance(P( 3 ),P( 6 ));
519 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
520 aVal = Max(aVal,Max(Max(L7,L8),L9));
523 else if( len == 8 ) { // hexas
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( 5 ),P( 6 ));
529 double L6 = getDistance(P( 6 ),P( 7 ));
530 double L7 = getDistance(P( 7 ),P( 8 ));
531 double L8 = getDistance(P( 8 ),P( 5 ));
532 double L9 = getDistance(P( 1 ),P( 5 ));
533 double L10= getDistance(P( 2 ),P( 6 ));
534 double L11= getDistance(P( 3 ),P( 7 ));
535 double L12= getDistance(P( 4 ),P( 8 ));
536 double D1 = getDistance(P( 1 ),P( 7 ));
537 double D2 = getDistance(P( 2 ),P( 8 ));
538 double D3 = getDistance(P( 3 ),P( 5 ));
539 double D4 = getDistance(P( 4 ),P( 6 ));
540 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
541 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
542 aVal = Max(aVal,Max(L11,L12));
543 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
546 else if( len == 10 ) { // quadratic tetras
547 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
548 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
549 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
550 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
551 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
552 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
553 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
556 else if( len == 13 ) { // quadratic pyramids
557 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
558 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
559 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
560 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
561 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
562 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
563 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
564 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
565 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
566 aVal = Max(aVal,Max(L7,L8));
569 else if( len == 15 ) { // quadratic pentas
570 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
571 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
572 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
573 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
574 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
575 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
576 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
577 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
578 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
579 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
580 aVal = Max(aVal,Max(Max(L7,L8),L9));
583 else if( len == 20 ) { // quadratic hexas
584 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
585 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
586 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
587 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
588 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
589 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
590 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
591 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
592 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
593 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
594 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
595 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
596 double D1 = getDistance(P( 1 ),P( 7 ));
597 double D2 = getDistance(P( 2 ),P( 8 ));
598 double D3 = getDistance(P( 3 ),P( 5 ));
599 double D4 = getDistance(P( 4 ),P( 6 ));
600 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
601 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
602 aVal = Max(aVal,Max(L11,L12));
603 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
606 else if( len > 1 && aElem->IsPoly() ) { // polys
607 // get the maximum distance between all pairs of nodes
608 for( int i = 1; i <= len; i++ ) {
609 for( int j = 1; j <= len; j++ ) {
610 if( j > i ) { // optimization of the loop
611 double D = getDistance( P(i), P(j) );
612 aVal = Max( aVal, D );
619 if( myPrecision >= 0 )
621 double prec = pow( 10., (double)myPrecision );
622 aVal = floor( aVal * prec + 0.5 ) / prec;
629 double MaxElementLength3D::GetBadRate( double Value, int /*nbNodes*/ ) const
634 SMDSAbs_ElementType MaxElementLength3D::GetType() const
636 return SMDSAbs_Volume;
642 Description : Functor for calculation of minimum angle
645 double MinimumAngle::GetValue( const TSequenceOfXYZ& P )
652 aMin = getAngle(P( P.size() ), P( 1 ), P( 2 ));
653 aMin = Min(aMin,getAngle(P( P.size()-1 ), P( P.size() ), P( 1 )));
655 for (int i=2; i<P.size();i++){
656 double A0 = getAngle( P( i-1 ), P( i ), P( i+1 ) );
660 return aMin * 180.0 / PI;
663 double MinimumAngle::GetBadRate( double Value, int nbNodes ) const
665 //const double aBestAngle = PI / nbNodes;
666 const double aBestAngle = 180.0 - ( 360.0 / double(nbNodes) );
667 return ( fabs( aBestAngle - Value ));
670 SMDSAbs_ElementType MinimumAngle::GetType() const
678 Description : Functor for calculating aspect ratio
680 double AspectRatio::GetValue( const TSequenceOfXYZ& P )
682 // According to "Mesh quality control" by Nadir Bouhamau referring to
683 // Pascal Jean Frey and Paul-Louis George. Maillages, applications aux elements finis.
684 // Hermes Science publications, Paris 1999 ISBN 2-7462-0024-4
687 int nbNodes = P.size();
692 // Compute aspect ratio
694 if ( nbNodes == 3 ) {
695 // Compute lengths of the sides
696 std::vector< double > aLen (nbNodes);
697 for ( int i = 0; i < nbNodes - 1; i++ )
698 aLen[ i ] = getDistance( P( i + 1 ), P( i + 2 ) );
699 aLen[ nbNodes - 1 ] = getDistance( P( 1 ), P( nbNodes ) );
700 // Q = alfa * h * p / S, where
702 // alfa = sqrt( 3 ) / 6
703 // h - length of the longest edge
704 // p - half perimeter
705 // S - triangle surface
706 const double alfa = sqrt( 3. ) / 6.;
707 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
708 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
709 double anArea = getArea( P( 1 ), P( 2 ), P( 3 ) );
710 if ( anArea <= Precision::Confusion() )
712 return alfa * maxLen * half_perimeter / anArea;
714 else if ( nbNodes == 6 ) { // quadratic triangles
715 // Compute lengths of the sides
716 std::vector< double > aLen (3);
717 aLen[0] = getDistance( P(1), P(3) );
718 aLen[1] = getDistance( P(3), P(5) );
719 aLen[2] = getDistance( P(5), P(1) );
720 // Q = alfa * h * p / S, where
722 // alfa = sqrt( 3 ) / 6
723 // h - length of the longest edge
724 // p - half perimeter
725 // S - triangle surface
726 const double alfa = sqrt( 3. ) / 6.;
727 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
728 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
729 double anArea = getArea( P(1), P(3), P(5) );
730 if ( anArea <= Precision::Confusion() )
732 return alfa * maxLen * half_perimeter / anArea;
734 else if( nbNodes == 4 ) { // quadrangle
735 // Compute lengths of the sides
736 std::vector< double > aLen (4);
737 aLen[0] = getDistance( P(1), P(2) );
738 aLen[1] = getDistance( P(2), P(3) );
739 aLen[2] = getDistance( P(3), P(4) );
740 aLen[3] = getDistance( P(4), P(1) );
741 // Compute lengths of the diagonals
742 std::vector< double > aDia (2);
743 aDia[0] = getDistance( P(1), P(3) );
744 aDia[1] = getDistance( P(2), P(4) );
745 // Compute areas of all triangles which can be built
746 // taking three nodes of the quadrangle
747 std::vector< double > anArea (4);
748 anArea[0] = getArea( P(1), P(2), P(3) );
749 anArea[1] = getArea( P(1), P(2), P(4) );
750 anArea[2] = getArea( P(1), P(3), P(4) );
751 anArea[3] = getArea( P(2), P(3), P(4) );
752 // Q = alpha * L * C1 / C2, where
754 // alpha = sqrt( 1/32 )
755 // L = max( L1, L2, L3, L4, D1, D2 )
756 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
757 // C2 = min( S1, S2, S3, S4 )
758 // Li - lengths of the edges
759 // Di - lengths of the diagonals
760 // Si - areas of the triangles
761 const double alpha = sqrt( 1 / 32. );
762 double L = Max( aLen[ 0 ],
766 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
767 double C1 = sqrt( ( aLen[0] * aLen[0] +
770 aLen[3] * aLen[3] ) / 4. );
771 double C2 = Min( anArea[ 0 ],
773 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
774 if ( C2 <= Precision::Confusion() )
776 return alpha * L * C1 / C2;
778 else if( nbNodes == 8 ){ // nbNodes==8 - quadratic quadrangle
779 // Compute lengths of the sides
780 std::vector< double > aLen (4);
781 aLen[0] = getDistance( P(1), P(3) );
782 aLen[1] = getDistance( P(3), P(5) );
783 aLen[2] = getDistance( P(5), P(7) );
784 aLen[3] = getDistance( P(7), P(1) );
785 // Compute lengths of the diagonals
786 std::vector< double > aDia (2);
787 aDia[0] = getDistance( P(1), P(5) );
788 aDia[1] = getDistance( P(3), P(7) );
789 // Compute areas of all triangles which can be built
790 // taking three nodes of the quadrangle
791 std::vector< double > anArea (4);
792 anArea[0] = getArea( P(1), P(3), P(5) );
793 anArea[1] = getArea( P(1), P(3), P(7) );
794 anArea[2] = getArea( P(1), P(5), P(7) );
795 anArea[3] = getArea( P(3), P(5), P(7) );
796 // Q = alpha * L * C1 / C2, where
798 // alpha = sqrt( 1/32 )
799 // L = max( L1, L2, L3, L4, D1, D2 )
800 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
801 // C2 = min( S1, S2, S3, S4 )
802 // Li - lengths of the edges
803 // Di - lengths of the diagonals
804 // Si - areas of the triangles
805 const double alpha = sqrt( 1 / 32. );
806 double L = Max( aLen[ 0 ],
810 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
811 double C1 = sqrt( ( aLen[0] * aLen[0] +
814 aLen[3] * aLen[3] ) / 4. );
815 double C2 = Min( anArea[ 0 ],
817 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
818 if ( C2 <= Precision::Confusion() )
820 return alpha * L * C1 / C2;
825 double AspectRatio::GetBadRate( double Value, int /*nbNodes*/ ) const
827 // the aspect ratio is in the range [1.0,infinity]
830 return Value / 1000.;
833 SMDSAbs_ElementType AspectRatio::GetType() const
840 Class : AspectRatio3D
841 Description : Functor for calculating aspect ratio
845 inline double getHalfPerimeter(double theTria[3]){
846 return (theTria[0] + theTria[1] + theTria[2])/2.0;
849 inline double getArea(double theHalfPerim, double theTria[3]){
850 return sqrt(theHalfPerim*
851 (theHalfPerim-theTria[0])*
852 (theHalfPerim-theTria[1])*
853 (theHalfPerim-theTria[2]));
856 inline double getVolume(double theLen[6]){
857 double a2 = theLen[0]*theLen[0];
858 double b2 = theLen[1]*theLen[1];
859 double c2 = theLen[2]*theLen[2];
860 double d2 = theLen[3]*theLen[3];
861 double e2 = theLen[4]*theLen[4];
862 double f2 = theLen[5]*theLen[5];
863 double P = 4.0*a2*b2*d2;
864 double Q = a2*(b2+d2-e2)-b2*(a2+d2-f2)-d2*(a2+b2-c2);
865 double R = (b2+d2-e2)*(a2+d2-f2)*(a2+d2-f2);
866 return sqrt(P-Q+R)/12.0;
869 inline double getVolume2(double theLen[6]){
870 double a2 = theLen[0]*theLen[0];
871 double b2 = theLen[1]*theLen[1];
872 double c2 = theLen[2]*theLen[2];
873 double d2 = theLen[3]*theLen[3];
874 double e2 = theLen[4]*theLen[4];
875 double f2 = theLen[5]*theLen[5];
877 double P = a2*e2*(b2+c2+d2+f2-a2-e2);
878 double Q = b2*f2*(a2+c2+d2+e2-b2-f2);
879 double R = c2*d2*(a2+b2+e2+f2-c2-d2);
880 double S = a2*b2*d2+b2*c2*e2+a2*c2*f2+d2*e2*f2;
882 return sqrt(P+Q+R-S)/12.0;
885 inline double getVolume(const TSequenceOfXYZ& P){
886 gp_Vec aVec1( P( 2 ) - P( 1 ) );
887 gp_Vec aVec2( P( 3 ) - P( 1 ) );
888 gp_Vec aVec3( P( 4 ) - P( 1 ) );
889 gp_Vec anAreaVec( aVec1 ^ aVec2 );
890 return fabs(aVec3 * anAreaVec) / 6.0;
893 inline double getMaxHeight(double theLen[6])
895 double aHeight = std::max(theLen[0],theLen[1]);
896 aHeight = std::max(aHeight,theLen[2]);
897 aHeight = std::max(aHeight,theLen[3]);
898 aHeight = std::max(aHeight,theLen[4]);
899 aHeight = std::max(aHeight,theLen[5]);
905 double AspectRatio3D::GetValue( const TSequenceOfXYZ& P )
907 double aQuality = 0.0;
908 if(myCurrElement->IsPoly()) return aQuality;
910 int nbNodes = P.size();
912 if(myCurrElement->IsQuadratic()) {
913 if(nbNodes==10) nbNodes=4; // quadratic tetrahedron
914 else if(nbNodes==13) nbNodes=5; // quadratic pyramid
915 else if(nbNodes==15) nbNodes=6; // quadratic pentahedron
916 else if(nbNodes==20) nbNodes=8; // quadratic hexahedron
917 else return aQuality;
923 getDistance(P( 1 ),P( 2 )), // a
924 getDistance(P( 2 ),P( 3 )), // b
925 getDistance(P( 3 ),P( 1 )), // c
926 getDistance(P( 2 ),P( 4 )), // d
927 getDistance(P( 3 ),P( 4 )), // e
928 getDistance(P( 1 ),P( 4 )) // f
930 double aTria[4][3] = {
931 {aLen[0],aLen[1],aLen[2]}, // abc
932 {aLen[0],aLen[3],aLen[5]}, // adf
933 {aLen[1],aLen[3],aLen[4]}, // bde
934 {aLen[2],aLen[4],aLen[5]} // cef
936 double aSumArea = 0.0;
937 double aHalfPerimeter = getHalfPerimeter(aTria[0]);
938 double anArea = getArea(aHalfPerimeter,aTria[0]);
940 aHalfPerimeter = getHalfPerimeter(aTria[1]);
941 anArea = getArea(aHalfPerimeter,aTria[1]);
943 aHalfPerimeter = getHalfPerimeter(aTria[2]);
944 anArea = getArea(aHalfPerimeter,aTria[2]);
946 aHalfPerimeter = getHalfPerimeter(aTria[3]);
947 anArea = getArea(aHalfPerimeter,aTria[3]);
949 double aVolume = getVolume(P);
950 //double aVolume = getVolume(aLen);
951 double aHeight = getMaxHeight(aLen);
952 static double aCoeff = sqrt(2.0)/12.0;
953 if ( aVolume > DBL_MIN )
954 aQuality = aCoeff*aHeight*aSumArea/aVolume;
959 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 3 ),P( 5 )};
960 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
963 gp_XYZ aXYZ[4] = {P( 1 ),P( 3 ),P( 4 ),P( 5 )};
964 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
967 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 5 )};
968 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
971 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 4 ),P( 5 )};
972 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
978 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 6 )};
979 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
982 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 3 )};
983 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
986 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 6 )};
987 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
990 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
991 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
994 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 6 )};
995 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
998 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 3 )};
999 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1005 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1006 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1009 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 4 )};
1010 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1013 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 7 )};
1014 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1017 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 8 )};
1018 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1021 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 3 )};
1022 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1025 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 4 )};
1026 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1029 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 7 )};
1030 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1033 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 8 )};
1034 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1037 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 3 )};
1038 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1041 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 4 )};
1042 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1045 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 7 )};
1046 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1049 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 8 )};
1050 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1053 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 1 )};
1054 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1057 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 2 )};
1058 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1061 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 5 )};
1062 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1065 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 6 )};
1066 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1069 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 1 )};
1070 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1073 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 2 )};
1074 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1077 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 5 )};
1078 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1081 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 6 )};
1082 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1085 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 1 )};
1086 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1089 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1090 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1093 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 5 )};
1094 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1097 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 6 )};
1098 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1101 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1102 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1105 gp_XYZ aXYZ[4] = {P( 4 ),P( 5 ),P( 8 ),P( 2 )};
1106 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1109 gp_XYZ aXYZ[4] = {P( 1 ),P( 4 ),P( 5 ),P( 3 )};
1110 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1113 gp_XYZ aXYZ[4] = {P( 3 ),P( 6 ),P( 7 ),P( 1 )};
1114 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1117 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 6 ),P( 4 )};
1118 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1121 gp_XYZ aXYZ[4] = {P( 5 ),P( 6 ),P( 8 ),P( 3 )};
1122 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1125 gp_XYZ aXYZ[4] = {P( 7 ),P( 8 ),P( 6 ),P( 1 )};
1126 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1129 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 7 )};
1130 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1133 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 2 ),P( 5 )};
1134 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1139 if ( nbNodes > 4 ) {
1140 // avaluate aspect ratio of quadranle faces
1141 AspectRatio aspect2D;
1142 SMDS_VolumeTool::VolumeType type = SMDS_VolumeTool::GetType( nbNodes );
1143 int nbFaces = SMDS_VolumeTool::NbFaces( type );
1144 TSequenceOfXYZ points(4);
1145 for ( int i = 0; i < nbFaces; ++i ) { // loop on faces of a volume
1146 if ( SMDS_VolumeTool::NbFaceNodes( type, i ) != 4 )
1148 const int* pInd = SMDS_VolumeTool::GetFaceNodesIndices( type, i, true );
1149 for ( int p = 0; p < 4; ++p ) // loop on nodes of a quadranle face
1150 points( p + 1 ) = P( pInd[ p ] + 1 );
1151 aQuality = std::max( aQuality, aspect2D.GetValue( points ));
1157 double AspectRatio3D::GetBadRate( double Value, int /*nbNodes*/ ) const
1159 // the aspect ratio is in the range [1.0,infinity]
1162 return Value / 1000.;
1165 SMDSAbs_ElementType AspectRatio3D::GetType() const
1167 return SMDSAbs_Volume;
1173 Description : Functor for calculating warping
1175 double Warping::GetValue( const TSequenceOfXYZ& P )
1177 if ( P.size() != 4 )
1180 gp_XYZ G = ( P( 1 ) + P( 2 ) + P( 3 ) + P( 4 ) ) / 4.;
1182 double A1 = ComputeA( P( 1 ), P( 2 ), P( 3 ), G );
1183 double A2 = ComputeA( P( 2 ), P( 3 ), P( 4 ), G );
1184 double A3 = ComputeA( P( 3 ), P( 4 ), P( 1 ), G );
1185 double A4 = ComputeA( P( 4 ), P( 1 ), P( 2 ), G );
1187 return Max( Max( A1, A2 ), Max( A3, A4 ) );
1190 double Warping::ComputeA( const gp_XYZ& thePnt1,
1191 const gp_XYZ& thePnt2,
1192 const gp_XYZ& thePnt3,
1193 const gp_XYZ& theG ) const
1195 double aLen1 = gp_Pnt( thePnt1 ).Distance( gp_Pnt( thePnt2 ) );
1196 double aLen2 = gp_Pnt( thePnt2 ).Distance( gp_Pnt( thePnt3 ) );
1197 double L = Min( aLen1, aLen2 ) * 0.5;
1198 if ( L < Precision::Confusion())
1201 gp_XYZ GI = ( thePnt2 + thePnt1 ) / 2. - theG;
1202 gp_XYZ GJ = ( thePnt3 + thePnt2 ) / 2. - theG;
1203 gp_XYZ N = GI.Crossed( GJ );
1205 if ( N.Modulus() < gp::Resolution() )
1210 double H = ( thePnt2 - theG ).Dot( N );
1211 return asin( fabs( H / L ) ) * 180. / PI;
1214 double Warping::GetBadRate( double Value, int /*nbNodes*/ ) const
1216 // the warp is in the range [0.0,PI/2]
1217 // 0.0 = good (no warp)
1218 // PI/2 = bad (face pliee)
1222 SMDSAbs_ElementType Warping::GetType() const
1224 return SMDSAbs_Face;
1230 Description : Functor for calculating taper
1232 double Taper::GetValue( const TSequenceOfXYZ& P )
1234 if ( P.size() != 4 )
1238 double J1 = getArea( P( 4 ), P( 1 ), P( 2 ) ) / 2.;
1239 double J2 = getArea( P( 3 ), P( 1 ), P( 2 ) ) / 2.;
1240 double J3 = getArea( P( 2 ), P( 3 ), P( 4 ) ) / 2.;
1241 double J4 = getArea( P( 3 ), P( 4 ), P( 1 ) ) / 2.;
1243 double JA = 0.25 * ( J1 + J2 + J3 + J4 );
1244 if ( JA <= Precision::Confusion() )
1247 double T1 = fabs( ( J1 - JA ) / JA );
1248 double T2 = fabs( ( J2 - JA ) / JA );
1249 double T3 = fabs( ( J3 - JA ) / JA );
1250 double T4 = fabs( ( J4 - JA ) / JA );
1252 return Max( Max( T1, T2 ), Max( T3, T4 ) );
1255 double Taper::GetBadRate( double Value, int /*nbNodes*/ ) const
1257 // the taper is in the range [0.0,1.0]
1258 // 0.0 = good (no taper)
1259 // 1.0 = bad (les cotes opposes sont allignes)
1263 SMDSAbs_ElementType Taper::GetType() const
1265 return SMDSAbs_Face;
1271 Description : Functor for calculating skew in degrees
1273 static inline double skewAngle( const gp_XYZ& p1, const gp_XYZ& p2, const gp_XYZ& p3 )
1275 gp_XYZ p12 = ( p2 + p1 ) / 2.;
1276 gp_XYZ p23 = ( p3 + p2 ) / 2.;
1277 gp_XYZ p31 = ( p3 + p1 ) / 2.;
1279 gp_Vec v1( p31 - p2 ), v2( p12 - p23 );
1281 return v1.Magnitude() < gp::Resolution() || v2.Magnitude() < gp::Resolution() ? 0. : v1.Angle( v2 );
1284 double Skew::GetValue( const TSequenceOfXYZ& P )
1286 if ( P.size() != 3 && P.size() != 4 )
1290 static double PI2 = PI / 2.;
1291 if ( P.size() == 3 )
1293 double A0 = fabs( PI2 - skewAngle( P( 3 ), P( 1 ), P( 2 ) ) );
1294 double A1 = fabs( PI2 - skewAngle( P( 1 ), P( 2 ), P( 3 ) ) );
1295 double A2 = fabs( PI2 - skewAngle( P( 2 ), P( 3 ), P( 1 ) ) );
1297 return Max( A0, Max( A1, A2 ) ) * 180. / PI;
1301 gp_XYZ p12 = ( P( 1 ) + P( 2 ) ) / 2.;
1302 gp_XYZ p23 = ( P( 2 ) + P( 3 ) ) / 2.;
1303 gp_XYZ p34 = ( P( 3 ) + P( 4 ) ) / 2.;
1304 gp_XYZ p41 = ( P( 4 ) + P( 1 ) ) / 2.;
1306 gp_Vec v1( p34 - p12 ), v2( p23 - p41 );
1307 double A = v1.Magnitude() <= gp::Resolution() || v2.Magnitude() <= gp::Resolution()
1308 ? 0. : fabs( PI2 - v1.Angle( v2 ) );
1311 if ( A < Precision::Angular() )
1314 return A * 180. / PI;
1318 double Skew::GetBadRate( double Value, int /*nbNodes*/ ) const
1320 // the skew is in the range [0.0,PI/2].
1326 SMDSAbs_ElementType Skew::GetType() const
1328 return SMDSAbs_Face;
1334 Description : Functor for calculating area
1336 double Area::GetValue( const TSequenceOfXYZ& P )
1339 if ( P.size() > 2 ) {
1340 gp_Vec aVec1( P(2) - P(1) );
1341 gp_Vec aVec2( P(3) - P(1) );
1342 gp_Vec SumVec = aVec1 ^ aVec2;
1343 for (int i=4; i<=P.size(); i++) {
1344 gp_Vec aVec1( P(i-1) - P(1) );
1345 gp_Vec aVec2( P(i) - P(1) );
1346 gp_Vec tmp = aVec1 ^ aVec2;
1349 val = SumVec.Magnitude() * 0.5;
1354 double Area::GetBadRate( double Value, int /*nbNodes*/ ) const
1356 // meaningless as it is not a quality control functor
1360 SMDSAbs_ElementType Area::GetType() const
1362 return SMDSAbs_Face;
1368 Description : Functor for calculating length off edge
1370 double Length::GetValue( const TSequenceOfXYZ& P )
1372 switch ( P.size() ) {
1373 case 2: return getDistance( P( 1 ), P( 2 ) );
1374 case 3: return getDistance( P( 1 ), P( 2 ) ) + getDistance( P( 2 ), P( 3 ) );
1379 double Length::GetBadRate( double Value, int /*nbNodes*/ ) const
1381 // meaningless as it is not quality control functor
1385 SMDSAbs_ElementType Length::GetType() const
1387 return SMDSAbs_Edge;
1392 Description : Functor for calculating length of edge
1395 double Length2D::GetValue( long theElementId)
1399 //cout<<"Length2D::GetValue"<<endl;
1400 if (GetPoints(theElementId,P)){
1401 //for(int jj=1; jj<=P.size(); jj++)
1402 // cout<<"jj="<<jj<<" P("<<P(jj).X()<<","<<P(jj).Y()<<","<<P(jj).Z()<<")"<<endl;
1404 double aVal;// = GetValue( P );
1405 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
1406 SMDSAbs_ElementType aType = aElem->GetType();
1415 aVal = getDistance( P( 1 ), P( 2 ) );
1418 else if (len == 3){ // quadratic edge
1419 aVal = getDistance(P( 1 ),P( 3 )) + getDistance(P( 3 ),P( 2 ));
1423 if (len == 3){ // triangles
1424 double L1 = getDistance(P( 1 ),P( 2 ));
1425 double L2 = getDistance(P( 2 ),P( 3 ));
1426 double L3 = getDistance(P( 3 ),P( 1 ));
1427 aVal = Max(L1,Max(L2,L3));
1430 else if (len == 4){ // quadrangles
1431 double L1 = getDistance(P( 1 ),P( 2 ));
1432 double L2 = getDistance(P( 2 ),P( 3 ));
1433 double L3 = getDistance(P( 3 ),P( 4 ));
1434 double L4 = getDistance(P( 4 ),P( 1 ));
1435 aVal = Max(Max(L1,L2),Max(L3,L4));
1438 if (len == 6){ // quadratic triangles
1439 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1440 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1441 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
1442 aVal = Max(L1,Max(L2,L3));
1443 //cout<<"L1="<<L1<<" L2="<<L2<<"L3="<<L3<<" aVal="<<aVal<<endl;
1446 else if (len == 8){ // quadratic quadrangles
1447 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1448 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1449 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
1450 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
1451 aVal = Max(Max(L1,L2),Max(L3,L4));
1454 case SMDSAbs_Volume:
1455 if (len == 4){ // tetraidrs
1456 double L1 = getDistance(P( 1 ),P( 2 ));
1457 double L2 = getDistance(P( 2 ),P( 3 ));
1458 double L3 = getDistance(P( 3 ),P( 1 ));
1459 double L4 = getDistance(P( 1 ),P( 4 ));
1460 double L5 = getDistance(P( 2 ),P( 4 ));
1461 double L6 = getDistance(P( 3 ),P( 4 ));
1462 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1465 else if (len == 5){ // piramids
1466 double L1 = getDistance(P( 1 ),P( 2 ));
1467 double L2 = getDistance(P( 2 ),P( 3 ));
1468 double L3 = getDistance(P( 3 ),P( 4 ));
1469 double L4 = getDistance(P( 4 ),P( 1 ));
1470 double L5 = getDistance(P( 1 ),P( 5 ));
1471 double L6 = getDistance(P( 2 ),P( 5 ));
1472 double L7 = getDistance(P( 3 ),P( 5 ));
1473 double L8 = getDistance(P( 4 ),P( 5 ));
1475 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1476 aVal = Max(aVal,Max(L7,L8));
1479 else if (len == 6){ // pentaidres
1480 double L1 = getDistance(P( 1 ),P( 2 ));
1481 double L2 = getDistance(P( 2 ),P( 3 ));
1482 double L3 = getDistance(P( 3 ),P( 1 ));
1483 double L4 = getDistance(P( 4 ),P( 5 ));
1484 double L5 = getDistance(P( 5 ),P( 6 ));
1485 double L6 = getDistance(P( 6 ),P( 4 ));
1486 double L7 = getDistance(P( 1 ),P( 4 ));
1487 double L8 = getDistance(P( 2 ),P( 5 ));
1488 double L9 = getDistance(P( 3 ),P( 6 ));
1490 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1491 aVal = Max(aVal,Max(Max(L7,L8),L9));
1494 else if (len == 8){ // hexaider
1495 double L1 = getDistance(P( 1 ),P( 2 ));
1496 double L2 = getDistance(P( 2 ),P( 3 ));
1497 double L3 = getDistance(P( 3 ),P( 4 ));
1498 double L4 = getDistance(P( 4 ),P( 1 ));
1499 double L5 = getDistance(P( 5 ),P( 6 ));
1500 double L6 = getDistance(P( 6 ),P( 7 ));
1501 double L7 = getDistance(P( 7 ),P( 8 ));
1502 double L8 = getDistance(P( 8 ),P( 5 ));
1503 double L9 = getDistance(P( 1 ),P( 5 ));
1504 double L10= getDistance(P( 2 ),P( 6 ));
1505 double L11= getDistance(P( 3 ),P( 7 ));
1506 double L12= getDistance(P( 4 ),P( 8 ));
1508 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1509 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1510 aVal = Max(aVal,Max(L11,L12));
1515 if (len == 10){ // quadratic tetraidrs
1516 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
1517 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
1518 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
1519 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1520 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
1521 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
1522 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1525 else if (len == 13){ // quadratic piramids
1526 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
1527 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
1528 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1529 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1530 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1531 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
1532 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
1533 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
1534 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1535 aVal = Max(aVal,Max(L7,L8));
1538 else if (len == 15){ // quadratic pentaidres
1539 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
1540 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
1541 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1542 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1543 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
1544 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
1545 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
1546 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
1547 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
1548 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1549 aVal = Max(aVal,Max(Max(L7,L8),L9));
1552 else if (len == 20){ // quadratic hexaider
1553 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
1554 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
1555 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
1556 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
1557 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
1558 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
1559 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
1560 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
1561 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
1562 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
1563 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
1564 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
1565 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1566 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1567 aVal = Max(aVal,Max(L11,L12));
1579 if ( myPrecision >= 0 )
1581 double prec = pow( 10., (double)( myPrecision ) );
1582 aVal = floor( aVal * prec + 0.5 ) / prec;
1591 double Length2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1593 // meaningless as it is not quality control functor
1597 SMDSAbs_ElementType Length2D::GetType() const
1599 return SMDSAbs_Face;
1602 Length2D::Value::Value(double theLength,long thePntId1, long thePntId2):
1605 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1606 if(thePntId1 > thePntId2){
1607 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1611 bool Length2D::Value::operator<(const Length2D::Value& x) const{
1612 if(myPntId[0] < x.myPntId[0]) return true;
1613 if(myPntId[0] == x.myPntId[0])
1614 if(myPntId[1] < x.myPntId[1]) return true;
1618 void Length2D::GetValues(TValues& theValues){
1620 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1621 for(; anIter->more(); ){
1622 const SMDS_MeshFace* anElem = anIter->next();
1624 if(anElem->IsQuadratic()) {
1625 const SMDS_QuadraticFaceOfNodes* F =
1626 static_cast<const SMDS_QuadraticFaceOfNodes*>(anElem);
1627 // use special nodes iterator
1628 SMDS_NodeIteratorPtr anIter = F->interlacedNodesIterator();
1633 const SMDS_MeshElement* aNode;
1635 aNode = anIter->next();
1636 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1637 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1638 aNodeId[0] = aNodeId[1] = aNode->GetID();
1641 for(; anIter->more(); ){
1642 const SMDS_MeshNode* N1 = static_cast<const SMDS_MeshNode*> (anIter->next());
1643 P[2] = gp_Pnt(N1->X(),N1->Y(),N1->Z());
1644 aNodeId[2] = N1->GetID();
1645 aLength = P[1].Distance(P[2]);
1646 if(!anIter->more()) break;
1647 const SMDS_MeshNode* N2 = static_cast<const SMDS_MeshNode*> (anIter->next());
1648 P[3] = gp_Pnt(N2->X(),N2->Y(),N2->Z());
1649 aNodeId[3] = N2->GetID();
1650 aLength += P[2].Distance(P[3]);
1651 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1652 Value aValue2(aLength,aNodeId[2],aNodeId[3]);
1654 aNodeId[1] = aNodeId[3];
1655 theValues.insert(aValue1);
1656 theValues.insert(aValue2);
1658 aLength += P[2].Distance(P[0]);
1659 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1660 Value aValue2(aLength,aNodeId[2],aNodeId[0]);
1661 theValues.insert(aValue1);
1662 theValues.insert(aValue2);
1665 SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
1670 const SMDS_MeshElement* aNode;
1671 if(aNodesIter->more()){
1672 aNode = aNodesIter->next();
1673 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1674 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1675 aNodeId[0] = aNodeId[1] = aNode->GetID();
1678 for(; aNodesIter->more(); ){
1679 aNode = aNodesIter->next();
1680 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1681 long anId = aNode->GetID();
1683 P[2] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1685 aLength = P[1].Distance(P[2]);
1687 Value aValue(aLength,aNodeId[1],anId);
1690 theValues.insert(aValue);
1693 aLength = P[0].Distance(P[1]);
1695 Value aValue(aLength,aNodeId[0],aNodeId[1]);
1696 theValues.insert(aValue);
1702 Class : MultiConnection
1703 Description : Functor for calculating number of faces conneted to the edge
1705 double MultiConnection::GetValue( const TSequenceOfXYZ& P )
1709 double MultiConnection::GetValue( long theId )
1711 return getNbMultiConnection( myMesh, theId );
1714 double MultiConnection::GetBadRate( double Value, int /*nbNodes*/ ) const
1716 // meaningless as it is not quality control functor
1720 SMDSAbs_ElementType MultiConnection::GetType() const
1722 return SMDSAbs_Edge;
1726 Class : MultiConnection2D
1727 Description : Functor for calculating number of faces conneted to the edge
1729 double MultiConnection2D::GetValue( const TSequenceOfXYZ& P )
1734 double MultiConnection2D::GetValue( long theElementId )
1738 const SMDS_MeshElement* aFaceElem = myMesh->FindElement(theElementId);
1739 SMDSAbs_ElementType aType = aFaceElem->GetType();
1744 int i = 0, len = aFaceElem->NbNodes();
1745 SMDS_ElemIteratorPtr anIter = aFaceElem->nodesIterator();
1748 const SMDS_MeshNode *aNode, *aNode0;
1749 TColStd_MapOfInteger aMap, aMapPrev;
1751 for (i = 0; i <= len; i++) {
1756 if (anIter->more()) {
1757 aNode = (SMDS_MeshNode*)anIter->next();
1765 if (i == 0) aNode0 = aNode;
1767 SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
1768 while (anElemIter->more()) {
1769 const SMDS_MeshElement* anElem = anElemIter->next();
1770 if (anElem != 0 && anElem->GetType() == SMDSAbs_Face) {
1771 int anId = anElem->GetID();
1774 if (aMapPrev.Contains(anId)) {
1779 aResult = Max(aResult, aNb);
1790 double MultiConnection2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1792 // meaningless as it is not quality control functor
1796 SMDSAbs_ElementType MultiConnection2D::GetType() const
1798 return SMDSAbs_Face;
1801 MultiConnection2D::Value::Value(long thePntId1, long thePntId2)
1803 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1804 if(thePntId1 > thePntId2){
1805 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1809 bool MultiConnection2D::Value::operator<(const MultiConnection2D::Value& x) const{
1810 if(myPntId[0] < x.myPntId[0]) return true;
1811 if(myPntId[0] == x.myPntId[0])
1812 if(myPntId[1] < x.myPntId[1]) return true;
1816 void MultiConnection2D::GetValues(MValues& theValues){
1817 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1818 for(; anIter->more(); ){
1819 const SMDS_MeshFace* anElem = anIter->next();
1820 SMDS_ElemIteratorPtr aNodesIter;
1821 if ( anElem->IsQuadratic() )
1822 aNodesIter = static_cast<const SMDS_QuadraticFaceOfNodes*>
1823 (anElem)->interlacedNodesElemIterator();
1825 aNodesIter = anElem->nodesIterator();
1828 //int aNbConnects=0;
1829 const SMDS_MeshNode* aNode0;
1830 const SMDS_MeshNode* aNode1;
1831 const SMDS_MeshNode* aNode2;
1832 if(aNodesIter->more()){
1833 aNode0 = (SMDS_MeshNode*) aNodesIter->next();
1835 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode1;
1836 aNodeId[0] = aNodeId[1] = aNodes->GetID();
1838 for(; aNodesIter->more(); ) {
1839 aNode2 = (SMDS_MeshNode*) aNodesIter->next();
1840 long anId = aNode2->GetID();
1843 Value aValue(aNodeId[1],aNodeId[2]);
1844 MValues::iterator aItr = theValues.find(aValue);
1845 if (aItr != theValues.end()){
1850 theValues[aValue] = 1;
1853 //cout << "NodeIds: "<<aNodeId[1]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1854 aNodeId[1] = aNodeId[2];
1857 Value aValue(aNodeId[0],aNodeId[2]);
1858 MValues::iterator aItr = theValues.find(aValue);
1859 if (aItr != theValues.end()) {
1864 theValues[aValue] = 1;
1867 //cout << "NodeIds: "<<aNodeId[0]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1877 Class : BadOrientedVolume
1878 Description : Predicate bad oriented volumes
1881 BadOrientedVolume::BadOrientedVolume()
1886 void BadOrientedVolume::SetMesh( const SMDS_Mesh* theMesh )
1891 bool BadOrientedVolume::IsSatisfy( long theId )
1896 SMDS_VolumeTool vTool( myMesh->FindElement( theId ));
1897 return !vTool.IsForward();
1900 SMDSAbs_ElementType BadOrientedVolume::GetType() const
1902 return SMDSAbs_Volume;
1906 Class : BareBorderVolume
1909 bool BareBorderVolume::IsSatisfy(long theElementId )
1911 SMDS_VolumeTool myTool;
1912 if ( myTool.Set( myMesh->FindElement(theElementId)))
1914 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
1915 if ( myTool.IsFreeFace( iF ))
1917 const SMDS_MeshNode** n = myTool.GetFaceNodes(iF);
1918 vector< const SMDS_MeshNode*> nodes( n, n+myTool.NbFaceNodes(iF));
1919 if ( !myMesh->FindElement( nodes, SMDSAbs_Face, /*Nomedium=*/false))
1927 Class : BareBorderFace
1930 bool BareBorderFace::IsSatisfy(long theElementId )
1932 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
1933 if ( face->GetType() == SMDSAbs_Face )
1935 int nbN = face->NbCornerNodes();
1936 for ( int i = 0; i < nbN; ++i )
1938 // check if a link is shared by another face
1939 const SMDS_MeshNode* n1 = face->GetNode( i );
1940 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
1941 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
1942 bool isShared = false;
1943 while ( !isShared && fIt->more() )
1945 const SMDS_MeshElement* f = fIt->next();
1946 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
1950 myLinkNodes.resize( 2 + face->IsQuadratic());
1951 myLinkNodes[0] = n1;
1952 myLinkNodes[1] = n2;
1953 if ( face->IsQuadratic() )
1954 myLinkNodes[2] = face->GetNode( i+nbN );
1955 return !myMesh->FindElement( myLinkNodes, SMDSAbs_Edge, /*noMedium=*/false);
1963 Class : OverConstrainedVolume
1966 bool OverConstrainedVolume::IsSatisfy(long theElementId )
1968 // An element is over-constrained if it has N-1 free borders where
1969 // N is the number of edges/faces for a 2D/3D element.
1970 SMDS_VolumeTool myTool;
1971 if ( myTool.Set( myMesh->FindElement(theElementId)))
1973 int nbSharedFaces = 0;
1974 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
1975 if ( !myTool.IsFreeFace( iF ) && ++nbSharedFaces > 1 )
1977 return ( nbSharedFaces == 1 );
1983 Class : OverConstrainedFace
1986 bool OverConstrainedFace::IsSatisfy(long theElementId )
1988 // An element is over-constrained if it has N-1 free borders where
1989 // N is the number of edges/faces for a 2D/3D element.
1990 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
1991 if ( face->GetType() == SMDSAbs_Face )
1993 int nbSharedBorders = 0;
1994 int nbN = face->NbCornerNodes();
1995 for ( int i = 0; i < nbN; ++i )
1997 // check if a link is shared by another face
1998 const SMDS_MeshNode* n1 = face->GetNode( i );
1999 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2000 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2001 bool isShared = false;
2002 while ( !isShared && fIt->more() )
2004 const SMDS_MeshElement* f = fIt->next();
2005 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2007 if ( isShared && ++nbSharedBorders > 1 )
2010 return ( nbSharedBorders == 1 );
2017 Description : Predicate for free borders
2020 FreeBorders::FreeBorders()
2025 void FreeBorders::SetMesh( const SMDS_Mesh* theMesh )
2030 bool FreeBorders::IsSatisfy( long theId )
2032 return getNbMultiConnection( myMesh, theId ) == 1;
2035 SMDSAbs_ElementType FreeBorders::GetType() const
2037 return SMDSAbs_Edge;
2043 Description : Predicate for free Edges
2045 FreeEdges::FreeEdges()
2050 void FreeEdges::SetMesh( const SMDS_Mesh* theMesh )
2055 bool FreeEdges::IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId )
2057 TColStd_MapOfInteger aMap;
2058 for ( int i = 0; i < 2; i++ )
2060 SMDS_ElemIteratorPtr anElemIter = theNodes[ i ]->GetInverseElementIterator();
2061 while( anElemIter->more() )
2063 const SMDS_MeshElement* anElem = anElemIter->next();
2064 if ( anElem != 0 && anElem->GetType() == SMDSAbs_Face )
2066 int anId = anElem->GetID();
2070 else if ( aMap.Contains( anId ) && anId != theFaceId )
2078 bool FreeEdges::IsSatisfy( long theId )
2083 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2084 if ( aFace == 0 || aFace->GetType() != SMDSAbs_Face || aFace->NbNodes() < 3 )
2087 SMDS_ElemIteratorPtr anIter;
2088 if ( aFace->IsQuadratic() ) {
2089 anIter = static_cast<const SMDS_QuadraticFaceOfNodes*>
2090 (aFace)->interlacedNodesElemIterator();
2093 anIter = aFace->nodesIterator();
2098 int i = 0, nbNodes = aFace->NbNodes();
2099 std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
2100 while( anIter->more() )
2102 const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
2105 aNodes[ i++ ] = aNode;
2107 aNodes[ nbNodes ] = aNodes[ 0 ];
2109 for ( i = 0; i < nbNodes; i++ )
2110 if ( IsFreeEdge( &aNodes[ i ], theId ) )
2116 SMDSAbs_ElementType FreeEdges::GetType() const
2118 return SMDSAbs_Face;
2121 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2124 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
2125 if(thePntId1 > thePntId2){
2126 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
2130 bool FreeEdges::Border::operator<(const FreeEdges::Border& x) const{
2131 if(myPntId[0] < x.myPntId[0]) return true;
2132 if(myPntId[0] == x.myPntId[0])
2133 if(myPntId[1] < x.myPntId[1]) return true;
2137 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2138 FreeEdges::TBorders& theRegistry,
2139 FreeEdges::TBorders& theContainer)
2141 if(theRegistry.find(theBorder) == theRegistry.end()){
2142 theRegistry.insert(theBorder);
2143 theContainer.insert(theBorder);
2145 theContainer.erase(theBorder);
2149 void FreeEdges::GetBoreders(TBorders& theBorders)
2152 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2153 for(; anIter->more(); ){
2154 const SMDS_MeshFace* anElem = anIter->next();
2155 long anElemId = anElem->GetID();
2156 SMDS_ElemIteratorPtr aNodesIter;
2157 if ( anElem->IsQuadratic() )
2158 aNodesIter = static_cast<const SMDS_QuadraticFaceOfNodes*>(anElem)->
2159 interlacedNodesElemIterator();
2161 aNodesIter = anElem->nodesIterator();
2163 const SMDS_MeshElement* aNode;
2164 if(aNodesIter->more()){
2165 aNode = aNodesIter->next();
2166 aNodeId[0] = aNodeId[1] = aNode->GetID();
2168 for(; aNodesIter->more(); ){
2169 aNode = aNodesIter->next();
2170 long anId = aNode->GetID();
2171 Border aBorder(anElemId,aNodeId[1],anId);
2173 //std::cout<<aBorder.myPntId[0]<<"; "<<aBorder.myPntId[1]<<"; "<<aBorder.myElemId<<endl;
2174 UpdateBorders(aBorder,aRegistry,theBorders);
2176 Border aBorder(anElemId,aNodeId[0],aNodeId[1]);
2177 //std::cout<<aBorder.myPntId[0]<<"; "<<aBorder.myPntId[1]<<"; "<<aBorder.myElemId<<endl;
2178 UpdateBorders(aBorder,aRegistry,theBorders);
2180 //std::cout<<"theBorders.size() = "<<theBorders.size()<<endl;
2186 Description : Predicate for free nodes
2189 FreeNodes::FreeNodes()
2194 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2199 bool FreeNodes::IsSatisfy( long theNodeId )
2201 const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2205 return (aNode->NbInverseElements() < 1);
2208 SMDSAbs_ElementType FreeNodes::GetType() const
2210 return SMDSAbs_Node;
2216 Description : Predicate for free faces
2219 FreeFaces::FreeFaces()
2224 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2229 bool FreeFaces::IsSatisfy( long theId )
2231 if (!myMesh) return false;
2232 // check that faces nodes refers to less than two common volumes
2233 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2234 if ( !aFace || aFace->GetType() != SMDSAbs_Face )
2237 int nbNode = aFace->NbNodes();
2239 // collect volumes check that number of volumss with count equal nbNode not less than 2
2240 typedef map< SMDS_MeshElement*, int > TMapOfVolume; // map of volume counters
2241 typedef map< SMDS_MeshElement*, int >::iterator TItrMapOfVolume; // iterator
2242 TMapOfVolume mapOfVol;
2244 SMDS_ElemIteratorPtr nodeItr = aFace->nodesIterator();
2245 while ( nodeItr->more() ) {
2246 const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>(nodeItr->next());
2247 if ( !aNode ) continue;
2248 SMDS_ElemIteratorPtr volItr = aNode->GetInverseElementIterator(SMDSAbs_Volume);
2249 while ( volItr->more() ) {
2250 SMDS_MeshElement* aVol = (SMDS_MeshElement*)volItr->next();
2251 TItrMapOfVolume itr = mapOfVol.insert(make_pair(aVol, 0)).first;
2256 TItrMapOfVolume volItr = mapOfVol.begin();
2257 TItrMapOfVolume volEnd = mapOfVol.end();
2258 for ( ; volItr != volEnd; ++volItr )
2259 if ( (*volItr).second >= nbNode )
2261 // face is not free if number of volumes constructed on thier nodes more than one
2265 SMDSAbs_ElementType FreeFaces::GetType() const
2267 return SMDSAbs_Face;
2271 Class : LinearOrQuadratic
2272 Description : Predicate to verify whether a mesh element is linear
2275 LinearOrQuadratic::LinearOrQuadratic()
2280 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2285 bool LinearOrQuadratic::IsSatisfy( long theId )
2287 if (!myMesh) return false;
2288 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2289 if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2291 return (!anElem->IsQuadratic());
2294 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2299 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2306 Description : Functor for check color of group to whic mesh element belongs to
2309 GroupColor::GroupColor()
2313 bool GroupColor::IsSatisfy( long theId )
2315 return (myIDs.find( theId ) != myIDs.end());
2318 void GroupColor::SetType( SMDSAbs_ElementType theType )
2323 SMDSAbs_ElementType GroupColor::GetType() const
2328 static bool isEqual( const Quantity_Color& theColor1,
2329 const Quantity_Color& theColor2 )
2331 // tolerance to compare colors
2332 const double tol = 5*1e-3;
2333 return ( fabs( theColor1.Red() - theColor2.Red() ) < tol &&
2334 fabs( theColor1.Green() - theColor2.Green() ) < tol &&
2335 fabs( theColor1.Blue() - theColor2.Blue() ) < tol );
2339 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2343 const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2347 int nbGrp = aMesh->GetNbGroups();
2351 // iterates on groups and find necessary elements ids
2352 const std::set<SMESHDS_GroupBase*>& aGroups = aMesh->GetGroups();
2353 set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
2354 for (; GrIt != aGroups.end(); GrIt++) {
2355 SMESHDS_GroupBase* aGrp = (*GrIt);
2358 // check type and color of group
2359 if ( !isEqual( myColor, aGrp->GetColor() ) )
2361 if ( myType != SMDSAbs_All && myType != (SMDSAbs_ElementType)aGrp->GetType() )
2364 SMDSAbs_ElementType aGrpElType = (SMDSAbs_ElementType)aGrp->GetType();
2365 if ( myType == aGrpElType || (myType == SMDSAbs_All && aGrpElType != SMDSAbs_Node) ) {
2366 // add elements IDS into control
2367 int aSize = aGrp->Extent();
2368 for (int i = 0; i < aSize; i++)
2369 myIDs.insert( aGrp->GetID(i+1) );
2374 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
2376 TCollection_AsciiString aStr = theStr;
2377 aStr.RemoveAll( ' ' );
2378 aStr.RemoveAll( '\t' );
2379 for ( int aPos = aStr.Search( ";;" ); aPos != -1; aPos = aStr.Search( ";;" ) )
2380 aStr.Remove( aPos, 2 );
2381 Standard_Real clr[3];
2382 clr[0] = clr[1] = clr[2] = 0.;
2383 for ( int i = 0; i < 3; i++ ) {
2384 TCollection_AsciiString tmpStr = aStr.Token( ";", i+1 );
2385 if ( !tmpStr.IsEmpty() && tmpStr.IsRealValue() )
2386 clr[i] = tmpStr.RealValue();
2388 myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
2391 //=======================================================================
2392 // name : GetRangeStr
2393 // Purpose : Get range as a string.
2394 // Example: "1,2,3,50-60,63,67,70-"
2395 //=======================================================================
2396 void GroupColor::GetColorStr( TCollection_AsciiString& theResStr ) const
2399 theResStr += TCollection_AsciiString( myColor.Red() );
2400 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2401 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2405 Class : ElemGeomType
2406 Description : Predicate to check element geometry type
2409 ElemGeomType::ElemGeomType()
2412 myType = SMDSAbs_All;
2413 myGeomType = SMDSGeom_TRIANGLE;
2416 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2421 bool ElemGeomType::IsSatisfy( long theId )
2423 if (!myMesh) return false;
2424 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2427 const SMDSAbs_ElementType anElemType = anElem->GetType();
2428 if ( myType != SMDSAbs_All && anElemType != myType )
2430 const int aNbNode = anElem->NbNodes();
2432 switch( anElemType )
2435 isOk = (myGeomType == SMDSGeom_POINT);
2439 isOk = (myGeomType == SMDSGeom_EDGE);
2443 if ( myGeomType == SMDSGeom_TRIANGLE )
2444 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 6 : aNbNode == 3));
2445 else if ( myGeomType == SMDSGeom_QUADRANGLE )
2446 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 8 : aNbNode == 4));
2447 else if ( myGeomType == SMDSGeom_POLYGON )
2448 isOk = anElem->IsPoly();
2451 case SMDSAbs_Volume:
2452 if ( myGeomType == SMDSGeom_TETRA )
2453 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 10 : aNbNode == 4));
2454 else if ( myGeomType == SMDSGeom_PYRAMID )
2455 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 13 : aNbNode == 5));
2456 else if ( myGeomType == SMDSGeom_PENTA )
2457 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 15 : aNbNode == 6));
2458 else if ( myGeomType == SMDSGeom_HEXA )
2459 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 20 : aNbNode == 8));
2460 else if ( myGeomType == SMDSGeom_POLYHEDRA )
2461 isOk = anElem->IsPoly();
2468 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2473 SMDSAbs_ElementType ElemGeomType::GetType() const
2478 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2480 myGeomType = theType;
2483 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2490 Description : Predicate for Range of Ids.
2491 Range may be specified with two ways.
2492 1. Using AddToRange method
2493 2. With SetRangeStr method. Parameter of this method is a string
2494 like as "1,2,3,50-60,63,67,70-"
2497 //=======================================================================
2498 // name : RangeOfIds
2499 // Purpose : Constructor
2500 //=======================================================================
2501 RangeOfIds::RangeOfIds()
2504 myType = SMDSAbs_All;
2507 //=======================================================================
2509 // Purpose : Set mesh
2510 //=======================================================================
2511 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
2516 //=======================================================================
2517 // name : AddToRange
2518 // Purpose : Add ID to the range
2519 //=======================================================================
2520 bool RangeOfIds::AddToRange( long theEntityId )
2522 myIds.Add( theEntityId );
2526 //=======================================================================
2527 // name : GetRangeStr
2528 // Purpose : Get range as a string.
2529 // Example: "1,2,3,50-60,63,67,70-"
2530 //=======================================================================
2531 void RangeOfIds::GetRangeStr( TCollection_AsciiString& theResStr )
2535 TColStd_SequenceOfInteger anIntSeq;
2536 TColStd_SequenceOfAsciiString aStrSeq;
2538 TColStd_MapIteratorOfMapOfInteger anIter( myIds );
2539 for ( ; anIter.More(); anIter.Next() )
2541 int anId = anIter.Key();
2542 TCollection_AsciiString aStr( anId );
2543 anIntSeq.Append( anId );
2544 aStrSeq.Append( aStr );
2547 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2549 int aMinId = myMin( i );
2550 int aMaxId = myMax( i );
2552 TCollection_AsciiString aStr;
2553 if ( aMinId != IntegerFirst() )
2558 if ( aMaxId != IntegerLast() )
2561 // find position of the string in result sequence and insert string in it
2562 if ( anIntSeq.Length() == 0 )
2564 anIntSeq.Append( aMinId );
2565 aStrSeq.Append( aStr );
2569 if ( aMinId < anIntSeq.First() )
2571 anIntSeq.Prepend( aMinId );
2572 aStrSeq.Prepend( aStr );
2574 else if ( aMinId > anIntSeq.Last() )
2576 anIntSeq.Append( aMinId );
2577 aStrSeq.Append( aStr );
2580 for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
2581 if ( aMinId < anIntSeq( j ) )
2583 anIntSeq.InsertBefore( j, aMinId );
2584 aStrSeq.InsertBefore( j, aStr );
2590 if ( aStrSeq.Length() == 0 )
2593 theResStr = aStrSeq( 1 );
2594 for ( int j = 2, k = aStrSeq.Length(); j <= k; j++ )
2597 theResStr += aStrSeq( j );
2601 //=======================================================================
2602 // name : SetRangeStr
2603 // Purpose : Define range with string
2604 // Example of entry string: "1,2,3,50-60,63,67,70-"
2605 //=======================================================================
2606 bool RangeOfIds::SetRangeStr( const TCollection_AsciiString& theStr )
2612 TCollection_AsciiString aStr = theStr;
2613 aStr.RemoveAll( ' ' );
2614 aStr.RemoveAll( '\t' );
2616 for ( int aPos = aStr.Search( ",," ); aPos != -1; aPos = aStr.Search( ",," ) )
2617 aStr.Remove( aPos, 2 );
2619 TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
2621 while ( tmpStr != "" )
2623 tmpStr = aStr.Token( ",", i++ );
2624 int aPos = tmpStr.Search( '-' );
2628 if ( tmpStr.IsIntegerValue() )
2629 myIds.Add( tmpStr.IntegerValue() );
2635 TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
2636 TCollection_AsciiString aMinStr = tmpStr;
2638 while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
2639 while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
2641 if ( !aMinStr.IsEmpty() && !aMinStr.IsIntegerValue() ||
2642 !aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue() )
2645 myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
2646 myMax.Append( aMaxStr.IsEmpty() ? IntegerLast() : aMaxStr.IntegerValue() );
2653 //=======================================================================
2655 // Purpose : Get type of supported entities
2656 //=======================================================================
2657 SMDSAbs_ElementType RangeOfIds::GetType() const
2662 //=======================================================================
2664 // Purpose : Set type of supported entities
2665 //=======================================================================
2666 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
2671 //=======================================================================
2673 // Purpose : Verify whether entity satisfies to this rpedicate
2674 //=======================================================================
2675 bool RangeOfIds::IsSatisfy( long theId )
2680 if ( myType == SMDSAbs_Node )
2682 if ( myMesh->FindNode( theId ) == 0 )
2687 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2688 if ( anElem == 0 || myType != anElem->GetType() && myType != SMDSAbs_All )
2692 if ( myIds.Contains( theId ) )
2695 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2696 if ( theId >= myMin( i ) && theId <= myMax( i ) )
2704 Description : Base class for comparators
2706 Comparator::Comparator():
2710 Comparator::~Comparator()
2713 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
2716 myFunctor->SetMesh( theMesh );
2719 void Comparator::SetMargin( double theValue )
2721 myMargin = theValue;
2724 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
2726 myFunctor = theFunct;
2729 SMDSAbs_ElementType Comparator::GetType() const
2731 return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
2734 double Comparator::GetMargin()
2742 Description : Comparator "<"
2744 bool LessThan::IsSatisfy( long theId )
2746 return myFunctor && myFunctor->GetValue( theId ) < myMargin;
2752 Description : Comparator ">"
2754 bool MoreThan::IsSatisfy( long theId )
2756 return myFunctor && myFunctor->GetValue( theId ) > myMargin;
2762 Description : Comparator "="
2765 myToler(Precision::Confusion())
2768 bool EqualTo::IsSatisfy( long theId )
2770 return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
2773 void EqualTo::SetTolerance( double theToler )
2778 double EqualTo::GetTolerance()
2785 Description : Logical NOT predicate
2787 LogicalNOT::LogicalNOT()
2790 LogicalNOT::~LogicalNOT()
2793 bool LogicalNOT::IsSatisfy( long theId )
2795 return myPredicate && !myPredicate->IsSatisfy( theId );
2798 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
2801 myPredicate->SetMesh( theMesh );
2804 void LogicalNOT::SetPredicate( PredicatePtr thePred )
2806 myPredicate = thePred;
2809 SMDSAbs_ElementType LogicalNOT::GetType() const
2811 return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
2816 Class : LogicalBinary
2817 Description : Base class for binary logical predicate
2819 LogicalBinary::LogicalBinary()
2822 LogicalBinary::~LogicalBinary()
2825 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
2828 myPredicate1->SetMesh( theMesh );
2831 myPredicate2->SetMesh( theMesh );
2834 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
2836 myPredicate1 = thePredicate;
2839 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
2841 myPredicate2 = thePredicate;
2844 SMDSAbs_ElementType LogicalBinary::GetType() const
2846 if ( !myPredicate1 || !myPredicate2 )
2849 SMDSAbs_ElementType aType1 = myPredicate1->GetType();
2850 SMDSAbs_ElementType aType2 = myPredicate2->GetType();
2852 return aType1 == aType2 ? aType1 : SMDSAbs_All;
2858 Description : Logical AND
2860 bool LogicalAND::IsSatisfy( long theId )
2865 myPredicate1->IsSatisfy( theId ) &&
2866 myPredicate2->IsSatisfy( theId );
2872 Description : Logical OR
2874 bool LogicalOR::IsSatisfy( long theId )
2879 myPredicate1->IsSatisfy( theId ) ||
2880 myPredicate2->IsSatisfy( theId );
2894 void Filter::SetPredicate( PredicatePtr thePredicate )
2896 myPredicate = thePredicate;
2899 template<class TElement, class TIterator, class TPredicate>
2900 inline void FillSequence(const TIterator& theIterator,
2901 TPredicate& thePredicate,
2902 Filter::TIdSequence& theSequence)
2904 if ( theIterator ) {
2905 while( theIterator->more() ) {
2906 TElement anElem = theIterator->next();
2907 long anId = anElem->GetID();
2908 if ( thePredicate->IsSatisfy( anId ) )
2909 theSequence.push_back( anId );
2916 GetElementsId( const SMDS_Mesh* theMesh,
2917 PredicatePtr thePredicate,
2918 TIdSequence& theSequence )
2920 theSequence.clear();
2922 if ( !theMesh || !thePredicate )
2925 thePredicate->SetMesh( theMesh );
2927 SMDSAbs_ElementType aType = thePredicate->GetType();
2930 FillSequence<const SMDS_MeshNode*>(theMesh->nodesIterator(),thePredicate,theSequence);
2933 FillSequence<const SMDS_MeshElement*>(theMesh->edgesIterator(),thePredicate,theSequence);
2936 FillSequence<const SMDS_MeshElement*>(theMesh->facesIterator(),thePredicate,theSequence);
2938 case SMDSAbs_Volume:
2939 FillSequence<const SMDS_MeshElement*>(theMesh->volumesIterator(),thePredicate,theSequence);
2942 FillSequence<const SMDS_MeshElement*>(theMesh->edgesIterator(),thePredicate,theSequence);
2943 FillSequence<const SMDS_MeshElement*>(theMesh->facesIterator(),thePredicate,theSequence);
2944 FillSequence<const SMDS_MeshElement*>(theMesh->volumesIterator(),thePredicate,theSequence);
2950 Filter::GetElementsId( const SMDS_Mesh* theMesh,
2951 Filter::TIdSequence& theSequence )
2953 GetElementsId(theMesh,myPredicate,theSequence);
2960 typedef std::set<SMDS_MeshFace*> TMapOfFacePtr;
2966 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
2967 SMDS_MeshNode* theNode2 )
2973 ManifoldPart::Link::~Link()
2979 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
2981 if ( myNode1 == theLink.myNode1 &&
2982 myNode2 == theLink.myNode2 )
2984 else if ( myNode1 == theLink.myNode2 &&
2985 myNode2 == theLink.myNode1 )
2991 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
2993 if(myNode1 < x.myNode1) return true;
2994 if(myNode1 == x.myNode1)
2995 if(myNode2 < x.myNode2) return true;
2999 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
3000 const ManifoldPart::Link& theLink2 )
3002 return theLink1.IsEqual( theLink2 );
3005 ManifoldPart::ManifoldPart()
3008 myAngToler = Precision::Angular();
3009 myIsOnlyManifold = true;
3012 ManifoldPart::~ManifoldPart()
3017 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
3023 SMDSAbs_ElementType ManifoldPart::GetType() const
3024 { return SMDSAbs_Face; }
3026 bool ManifoldPart::IsSatisfy( long theElementId )
3028 return myMapIds.Contains( theElementId );
3031 void ManifoldPart::SetAngleTolerance( const double theAngToler )
3032 { myAngToler = theAngToler; }
3034 double ManifoldPart::GetAngleTolerance() const
3035 { return myAngToler; }
3037 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
3038 { myIsOnlyManifold = theIsOnly; }
3040 void ManifoldPart::SetStartElem( const long theStartId )
3041 { myStartElemId = theStartId; }
3043 bool ManifoldPart::process()
3046 myMapBadGeomIds.Clear();
3048 myAllFacePtr.clear();
3049 myAllFacePtrIntDMap.clear();
3053 // collect all faces into own map
3054 SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3055 for (; anFaceItr->more(); )
3057 SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3058 myAllFacePtr.push_back( aFacePtr );
3059 myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3062 SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3066 // the map of non manifold links and bad geometry
3067 TMapOfLink aMapOfNonManifold;
3068 TColStd_MapOfInteger aMapOfTreated;
3070 // begin cycle on faces from start index and run on vector till the end
3071 // and from begin to start index to cover whole vector
3072 const int aStartIndx = myAllFacePtrIntDMap[aStartFace];
3073 bool isStartTreat = false;
3074 for ( int fi = aStartIndx; !isStartTreat || fi != aStartIndx ; fi++ )
3076 if ( fi == aStartIndx )
3077 isStartTreat = true;
3078 // as result next time when fi will be equal to aStartIndx
3080 SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3081 if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3084 aMapOfTreated.Add( aFacePtr->GetID() );
3085 TColStd_MapOfInteger aResFaces;
3086 if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3087 aMapOfNonManifold, aResFaces ) )
3089 TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3090 for ( ; anItr.More(); anItr.Next() )
3092 int aFaceId = anItr.Key();
3093 aMapOfTreated.Add( aFaceId );
3094 myMapIds.Add( aFaceId );
3097 if ( fi == ( myAllFacePtr.size() - 1 ) )
3099 } // end run on vector of faces
3100 return !myMapIds.IsEmpty();
3103 static void getLinks( const SMDS_MeshFace* theFace,
3104 ManifoldPart::TVectorOfLink& theLinks )
3106 int aNbNode = theFace->NbNodes();
3107 SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3109 SMDS_MeshNode* aNode = 0;
3110 for ( ; aNodeItr->more() && i <= aNbNode; )
3113 SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3117 SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3119 ManifoldPart::Link aLink( aN1, aN2 );
3120 theLinks.push_back( aLink );
3124 static gp_XYZ getNormale( const SMDS_MeshFace* theFace )
3127 int aNbNode = theFace->NbNodes();
3128 TColgp_Array1OfXYZ anArrOfXYZ(1,4);
3129 SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3131 for ( ; aNodeItr->more() && i <= 4; i++ ) {
3132 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3133 anArrOfXYZ.SetValue(i, gp_XYZ( aNode->X(), aNode->Y(), aNode->Z() ) );
3136 gp_XYZ q1 = anArrOfXYZ.Value(2) - anArrOfXYZ.Value(1);
3137 gp_XYZ q2 = anArrOfXYZ.Value(3) - anArrOfXYZ.Value(1);
3139 if ( aNbNode > 3 ) {
3140 gp_XYZ q3 = anArrOfXYZ.Value(4) - anArrOfXYZ.Value(1);
3143 double len = n.Modulus();
3150 bool ManifoldPart::findConnected
3151 ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3152 SMDS_MeshFace* theStartFace,
3153 ManifoldPart::TMapOfLink& theNonManifold,
3154 TColStd_MapOfInteger& theResFaces )
3156 theResFaces.Clear();
3157 if ( !theAllFacePtrInt.size() )
3160 if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3162 myMapBadGeomIds.Add( theStartFace->GetID() );
3166 ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3167 ManifoldPart::TVectorOfLink aSeqOfBoundary;
3168 theResFaces.Add( theStartFace->GetID() );
3169 ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3171 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3172 aDMapLinkFace, theNonManifold, theStartFace );
3174 bool isDone = false;
3175 while ( !isDone && aMapOfBoundary.size() != 0 )
3177 bool isToReset = false;
3178 ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3179 for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3181 ManifoldPart::Link aLink = *pLink;
3182 if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3184 // each link could be treated only once
3185 aMapToSkip.insert( aLink );
3187 ManifoldPart::TVectorOfFacePtr aFaces;
3189 if ( myIsOnlyManifold &&
3190 (theNonManifold.find( aLink ) != theNonManifold.end()) )
3194 getFacesByLink( aLink, aFaces );
3195 // filter the element to keep only indicated elements
3196 ManifoldPart::TVectorOfFacePtr aFiltered;
3197 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3198 for ( ; pFace != aFaces.end(); ++pFace )
3200 SMDS_MeshFace* aFace = *pFace;
3201 if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3202 aFiltered.push_back( aFace );
3205 if ( aFaces.size() < 2 ) // no neihgbour faces
3207 else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3209 theNonManifold.insert( aLink );
3214 // compare normal with normals of neighbor element
3215 SMDS_MeshFace* aPrevFace = aDMapLinkFace[ aLink ];
3216 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3217 for ( ; pFace != aFaces.end(); ++pFace )
3219 SMDS_MeshFace* aNextFace = *pFace;
3220 if ( aPrevFace == aNextFace )
3222 int anNextFaceID = aNextFace->GetID();
3223 if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3224 // should not be with non manifold restriction. probably bad topology
3226 // check if face was treated and skipped
3227 if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3228 !isInPlane( aPrevFace, aNextFace ) )
3230 // add new element to connected and extend the boundaries.
3231 theResFaces.Add( anNextFaceID );
3232 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3233 aDMapLinkFace, theNonManifold, aNextFace );
3237 isDone = !isToReset;
3240 return !theResFaces.IsEmpty();
3243 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3244 const SMDS_MeshFace* theFace2 )
3246 gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3247 gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3248 if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3250 myMapBadGeomIds.Add( theFace2->GetID() );
3253 if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
3259 void ManifoldPart::expandBoundary
3260 ( ManifoldPart::TMapOfLink& theMapOfBoundary,
3261 ManifoldPart::TVectorOfLink& theSeqOfBoundary,
3262 ManifoldPart::TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
3263 ManifoldPart::TMapOfLink& theNonManifold,
3264 SMDS_MeshFace* theNextFace ) const
3266 ManifoldPart::TVectorOfLink aLinks;
3267 getLinks( theNextFace, aLinks );
3268 int aNbLink = (int)aLinks.size();
3269 for ( int i = 0; i < aNbLink; i++ )
3271 ManifoldPart::Link aLink = aLinks[ i ];
3272 if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3274 if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3276 if ( myIsOnlyManifold )
3278 // remove from boundary
3279 theMapOfBoundary.erase( aLink );
3280 ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3281 for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3283 ManifoldPart::Link aBoundLink = *pLink;
3284 if ( aBoundLink.IsEqual( aLink ) )
3286 theSeqOfBoundary.erase( pLink );
3294 theMapOfBoundary.insert( aLink );
3295 theSeqOfBoundary.push_back( aLink );
3296 theDMapLinkFacePtr[ aLink ] = theNextFace;
3301 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3302 ManifoldPart::TVectorOfFacePtr& theFaces ) const
3304 SMDS_Mesh::SetOfFaces aSetOfFaces;
3305 // take all faces that shared first node
3306 SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3307 for ( ; anItr->more(); )
3309 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3312 aSetOfFaces.Add( aFace );
3314 // take all faces that shared second node
3315 anItr = theLink.myNode2->facesIterator();
3316 // find the common part of two sets
3317 for ( ; anItr->more(); )
3319 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3320 if ( aSetOfFaces.Contains( aFace ) )
3321 theFaces.push_back( aFace );
3330 ElementsOnSurface::ElementsOnSurface()
3334 myType = SMDSAbs_All;
3336 myToler = Precision::Confusion();
3337 myUseBoundaries = false;
3340 ElementsOnSurface::~ElementsOnSurface()
3345 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
3347 if ( myMesh == theMesh )
3353 bool ElementsOnSurface::IsSatisfy( long theElementId )
3355 return myIds.Contains( theElementId );
3358 SMDSAbs_ElementType ElementsOnSurface::GetType() const
3361 void ElementsOnSurface::SetTolerance( const double theToler )
3363 if ( myToler != theToler )
3368 double ElementsOnSurface::GetTolerance() const
3371 void ElementsOnSurface::SetUseBoundaries( bool theUse )
3373 if ( myUseBoundaries != theUse ) {
3374 myUseBoundaries = theUse;
3375 SetSurface( mySurf, myType );
3379 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
3380 const SMDSAbs_ElementType theType )
3385 if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
3387 mySurf = TopoDS::Face( theShape );
3388 BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
3390 u1 = SA.FirstUParameter(),
3391 u2 = SA.LastUParameter(),
3392 v1 = SA.FirstVParameter(),
3393 v2 = SA.LastVParameter();
3394 Handle(Geom_Surface) surf = BRep_Tool::Surface( mySurf );
3395 myProjector.Init( surf, u1,u2, v1,v2 );
3399 void ElementsOnSurface::process()
3402 if ( mySurf.IsNull() )
3408 if ( myType == SMDSAbs_Face || myType == SMDSAbs_All )
3410 myIds.ReSize( myMesh->NbFaces() );
3411 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
3412 for(; anIter->more(); )
3413 process( anIter->next() );
3416 if ( myType == SMDSAbs_Edge || myType == SMDSAbs_All )
3418 myIds.ReSize( myIds.Extent() + myMesh->NbEdges() );
3419 SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
3420 for(; anIter->more(); )
3421 process( anIter->next() );
3424 if ( myType == SMDSAbs_Node )
3426 myIds.ReSize( myMesh->NbNodes() );
3427 SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
3428 for(; anIter->more(); )
3429 process( anIter->next() );
3433 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
3435 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3436 bool isSatisfy = true;
3437 for ( ; aNodeItr->more(); )
3439 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3440 if ( !isOnSurface( aNode ) )
3447 myIds.Add( theElemPtr->GetID() );
3450 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
3452 if ( mySurf.IsNull() )
3455 gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
3456 // double aToler2 = myToler * myToler;
3457 // if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
3459 // gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
3460 // if ( aPln.SquareDistance( aPnt ) > aToler2 )
3463 // else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
3465 // gp_Cylinder aCyl = Handle(Geom_CylindricalSurface)::DownCast(mySurf)->Cylinder();
3466 // double aRad = aCyl.Radius();
3467 // gp_Ax3 anAxis = aCyl.Position();
3468 // gp_XYZ aLoc = aCyl.Location().XYZ();
3469 // double aXDist = anAxis.XDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3470 // double aYDist = anAxis.YDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3471 // if ( fabs(aXDist*aXDist + aYDist*aYDist - aRad*aRad) > aToler2 )
3476 myProjector.Perform( aPnt );
3477 bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
3487 ElementsOnShape::ElementsOnShape()
3489 myType(SMDSAbs_All),
3490 myToler(Precision::Confusion()),
3491 myAllNodesFlag(false)
3493 myCurShapeType = TopAbs_SHAPE;
3496 ElementsOnShape::~ElementsOnShape()
3500 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
3502 if (myMesh != theMesh) {
3504 SetShape(myShape, myType);
3508 bool ElementsOnShape::IsSatisfy (long theElementId)
3510 return myIds.Contains(theElementId);
3513 SMDSAbs_ElementType ElementsOnShape::GetType() const
3518 void ElementsOnShape::SetTolerance (const double theToler)
3520 if (myToler != theToler) {
3522 SetShape(myShape, myType);
3526 double ElementsOnShape::GetTolerance() const
3531 void ElementsOnShape::SetAllNodes (bool theAllNodes)
3533 if (myAllNodesFlag != theAllNodes) {
3534 myAllNodesFlag = theAllNodes;
3535 SetShape(myShape, myType);
3539 void ElementsOnShape::SetShape (const TopoDS_Shape& theShape,
3540 const SMDSAbs_ElementType theType)
3546 if (myMesh == 0) return;
3551 myIds.ReSize(myMesh->NbEdges() + myMesh->NbFaces() + myMesh->NbVolumes());
3554 myIds.ReSize(myMesh->NbNodes());
3557 myIds.ReSize(myMesh->NbEdges());
3560 myIds.ReSize(myMesh->NbFaces());
3562 case SMDSAbs_Volume:
3563 myIds.ReSize(myMesh->NbVolumes());
3569 myShapesMap.Clear();
3573 void ElementsOnShape::addShape (const TopoDS_Shape& theShape)
3575 if (theShape.IsNull() || myMesh == 0)
3578 if (!myShapesMap.Add(theShape)) return;
3580 myCurShapeType = theShape.ShapeType();
3581 switch (myCurShapeType)
3583 case TopAbs_COMPOUND:
3584 case TopAbs_COMPSOLID:
3588 TopoDS_Iterator anIt (theShape, Standard_True, Standard_True);
3589 for (; anIt.More(); anIt.Next()) addShape(anIt.Value());
3594 myCurSC.Load(theShape);
3600 TopoDS_Face aFace = TopoDS::Face(theShape);
3601 BRepAdaptor_Surface SA (aFace, true);
3603 u1 = SA.FirstUParameter(),
3604 u2 = SA.LastUParameter(),
3605 v1 = SA.FirstVParameter(),
3606 v2 = SA.LastVParameter();
3607 Handle(Geom_Surface) surf = BRep_Tool::Surface(aFace);
3608 myCurProjFace.Init(surf, u1,u2, v1,v2);
3615 TopoDS_Edge anEdge = TopoDS::Edge(theShape);
3616 Standard_Real u1, u2;
3617 Handle(Geom_Curve) curve = BRep_Tool::Curve(anEdge, u1, u2);
3618 myCurProjEdge.Init(curve, u1, u2);
3624 TopoDS_Vertex aV = TopoDS::Vertex(theShape);
3625 myCurPnt = BRep_Tool::Pnt(aV);
3634 void ElementsOnShape::process()
3636 if (myShape.IsNull() || myMesh == 0)
3639 if (myType == SMDSAbs_Node)
3641 SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
3642 while (anIter->more())
3643 process(anIter->next());
3647 if (myType == SMDSAbs_Edge || myType == SMDSAbs_All)
3649 SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
3650 while (anIter->more())
3651 process(anIter->next());
3654 if (myType == SMDSAbs_Face || myType == SMDSAbs_All)
3656 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
3657 while (anIter->more()) {
3658 process(anIter->next());
3662 if (myType == SMDSAbs_Volume || myType == SMDSAbs_All)
3664 SMDS_VolumeIteratorPtr anIter = myMesh->volumesIterator();
3665 while (anIter->more())
3666 process(anIter->next());
3671 void ElementsOnShape::process (const SMDS_MeshElement* theElemPtr)
3673 if (myShape.IsNull())
3676 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3677 bool isSatisfy = myAllNodesFlag;
3679 gp_XYZ centerXYZ (0, 0, 0);
3681 while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
3683 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3684 gp_Pnt aPnt (aNode->X(), aNode->Y(), aNode->Z());
3685 centerXYZ += aPnt.XYZ();
3687 switch (myCurShapeType)
3691 myCurSC.Perform(aPnt, myToler);
3692 isSatisfy = (myCurSC.State() == TopAbs_IN || myCurSC.State() == TopAbs_ON);
3697 myCurProjFace.Perform(aPnt);
3698 isSatisfy = (myCurProjFace.IsDone() && myCurProjFace.LowerDistance() <= myToler);
3701 // check relatively the face
3702 Quantity_Parameter u, v;
3703 myCurProjFace.LowerDistanceParameters(u, v);
3704 gp_Pnt2d aProjPnt (u, v);
3705 BRepClass_FaceClassifier aClsf (myCurFace, aProjPnt, myToler);
3706 isSatisfy = (aClsf.State() == TopAbs_IN || aClsf.State() == TopAbs_ON);
3712 myCurProjEdge.Perform(aPnt);
3713 isSatisfy = (myCurProjEdge.NbPoints() > 0 && myCurProjEdge.LowerDistance() <= myToler);
3718 isSatisfy = (aPnt.Distance(myCurPnt) <= myToler);
3728 if (isSatisfy && myCurShapeType == TopAbs_SOLID) { // Check the center point for volumes MantisBug 0020168
3729 centerXYZ /= theElemPtr->NbNodes();
3730 gp_Pnt aCenterPnt (centerXYZ);
3731 myCurSC.Perform(aCenterPnt, myToler);
3732 if ( !(myCurSC.State() == TopAbs_IN || myCurSC.State() == TopAbs_ON))
3737 myIds.Add(theElemPtr->GetID());
3740 TSequenceOfXYZ::TSequenceOfXYZ()
3743 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n)
3746 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t)
3749 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray)
3752 template <class InputIterator>
3753 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd)
3756 TSequenceOfXYZ::~TSequenceOfXYZ()
3759 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
3761 myArray = theSequenceOfXYZ.myArray;
3765 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
3767 return myArray[n-1];
3770 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
3772 return myArray[n-1];
3775 void TSequenceOfXYZ::clear()
3780 void TSequenceOfXYZ::reserve(size_type n)
3785 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
3787 myArray.push_back(v);
3790 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
3792 return myArray.size();