Salome HOME
compilation warning
[modules/smesh.git] / src / Controls / SMESH_Controls.cxx
1 //  Copyright (C) 2007-2010  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 //  Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 //  CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
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.
10 //
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.
15 //
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
19 //
20 //  See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22
23 #include "SMESH_ControlsDef.hxx"
24
25 #include <set>
26
27 #include <BRepAdaptor_Surface.hxx>
28 #include <BRepClass_FaceClassifier.hxx>
29 #include <BRep_Tool.hxx>
30
31 #include <TopAbs.hxx>
32 #include <TopoDS.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>
38
39 #include <Geom_CylindricalSurface.hxx>
40 #include <Geom_Plane.hxx>
41 #include <Geom_Surface.hxx>
42
43 #include <Precision.hxx>
44 #include <TColStd_MapIteratorOfMapOfInteger.hxx>
45 #include <TColStd_MapOfInteger.hxx>
46 #include <TColStd_SequenceOfAsciiString.hxx>
47 #include <TColgp_Array1OfXYZ.hxx>
48
49 #include <gp_Ax3.hxx>
50 #include <gp_Cylinder.hxx>
51 #include <gp_Dir.hxx>
52 #include <gp_Pln.hxx>
53 #include <gp_Pnt.hxx>
54 #include <gp_Vec.hxx>
55 #include <gp_XYZ.hxx>
56
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"
64
65 #include "SMESHDS_Mesh.hxx"
66 #include "SMESHDS_GroupBase.hxx"
67
68 /*
69                             AUXILIARY METHODS
70 */
71
72 namespace{
73   inline double getAngle( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
74   {
75     gp_Vec v1( P1 - P2 ), v2( P3 - P2 );
76
77     return v1.Magnitude() < gp::Resolution() ||
78       v2.Magnitude() < gp::Resolution() ? 0 : v1.Angle( v2 );
79   }
80
81   inline double getArea( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
82   {
83     gp_Vec aVec1( P2 - P1 );
84     gp_Vec aVec2( P3 - P1 );
85     return ( aVec1 ^ aVec2 ).Magnitude() * 0.5;
86   }
87
88   inline double getArea( const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3 )
89   {
90     return getArea( P1.XYZ(), P2.XYZ(), P3.XYZ() );
91   }
92
93
94
95   inline double getDistance( const gp_XYZ& P1, const gp_XYZ& P2 )
96   {
97     double aDist = gp_Pnt( P1 ).Distance( gp_Pnt( P2 ) );
98     return aDist;
99   }
100
101   int getNbMultiConnection( const SMDS_Mesh* theMesh, const int theId )
102   {
103     if ( theMesh == 0 )
104       return 0;
105
106     const SMDS_MeshElement* anEdge = theMesh->FindElement( theId );
107     if ( anEdge == 0 || anEdge->GetType() != SMDSAbs_Edge/* || anEdge->NbNodes() != 2 */)
108       return 0;
109
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):
113     //
114     //  Case 1          Case 2
115     //  |     |      |        |      |
116     //  |     |      |        |      |
117     //  +-----+------+  +-----+------+ 
118     //  |            |  |            |
119     //  |            |  |            |
120     // result sould be 2 in both cases
121     //
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;
128
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 ) {
137               aResult0++;
138               if ( !aNode1 ) break; // not a quadratic edge
139             }
140             else if ( anElemNode == aNode1 )
141               aResult1++;
142           }
143         }
144       }
145     }
146     int aResult = std::max ( aResult0, aResult1 );
147
148 //     TColStd_MapOfInteger aMap;
149
150 //     SMDS_ElemIteratorPtr anIter = anEdge->nodesIterator();
151 //     if ( anIter != 0 ) {
152 //       while( anIter->more() ) {
153 //      const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
154 //      if ( aNode == 0 )
155 //        return 0;
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();
161
162 //          if ( anIter->more() )              // i.e. first node
163 //            aMap.Add( anId );
164 //          else if ( aMap.Contains( anId ) )
165 //            aResult++;
166 //        }
167 //      }
168 //       }
169 //     }
170
171     return aResult;
172   }
173
174 }
175
176
177
178 using namespace SMESH::Controls;
179
180 /*
181                                 FUNCTORS
182 */
183
184 /*
185   Class       : NumericalFunctor
186   Description : Base class for numerical functors
187 */
188 NumericalFunctor::NumericalFunctor():
189   myMesh(NULL)
190 {
191   myPrecision = -1;
192 }
193
194 void NumericalFunctor::SetMesh( const SMDS_Mesh* theMesh )
195 {
196   myMesh = theMesh;
197 }
198
199 bool NumericalFunctor::GetPoints(const int theId,
200                                  TSequenceOfXYZ& theRes ) const
201 {
202   theRes.clear();
203
204   if ( myMesh == 0 )
205     return false;
206
207   return GetPoints( myMesh->FindElement( theId ), theRes );
208 }
209
210 bool NumericalFunctor::GetPoints(const SMDS_MeshElement* anElem,
211                                  TSequenceOfXYZ&         theRes )
212 {
213   theRes.clear();
214
215   if ( anElem == 0)
216     return false;
217
218   theRes.reserve( anElem->NbNodes() );
219
220   // Get nodes of the element
221   SMDS_ElemIteratorPtr anIter;
222
223   if ( anElem->IsQuadratic() ) {
224     switch ( anElem->GetType() ) {
225     case SMDSAbs_Edge:
226       anIter = static_cast<const SMDS_QuadraticEdge*>
227         (anElem)->interlacedNodesElemIterator();
228       break;
229     case SMDSAbs_Face:
230       anIter = static_cast<const SMDS_QuadraticFaceOfNodes*>
231         (anElem)->interlacedNodesElemIterator();
232       break;
233     default:
234       anIter = anElem->nodesIterator();
235       //return false;
236     }
237   }
238   else {
239     anIter = anElem->nodesIterator();
240   }
241
242   if ( anIter ) {
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() ) );
246     }
247   }
248
249   return true;
250 }
251
252 long  NumericalFunctor::GetPrecision() const
253 {
254   return myPrecision;
255 }
256
257 void  NumericalFunctor::SetPrecision( const long thePrecision )
258 {
259   myPrecision = thePrecision;
260 }
261
262 double NumericalFunctor::GetValue( long theId )
263 {
264   myCurrElement = myMesh->FindElement( theId );
265   TSequenceOfXYZ P;
266   if ( GetPoints( theId, P ))
267   {
268     double aVal = GetValue( P );
269     if ( myPrecision >= 0 )
270     {
271       double prec = pow( 10., (double)( myPrecision ) );
272       aVal = floor( aVal * prec + 0.5 ) / prec;
273     }
274     return aVal;
275   }
276
277   return 0.;
278 }
279
280 //================================================================================
281 /*!
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
288  */
289 //================================================================================
290
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)
296 {
297   if ( nbIntervals < 1 ||
298        !myMesh ||
299        !myMesh->GetMeshInfo().NbElements( GetType() ))
300     return;
301   nbEvents.resize( nbIntervals, 0 );
302   funValues.resize( nbIntervals+1 );
303
304   // get all values sorted
305   std::multiset< double > values;
306   if ( elements.empty() )
307   {
308     SMDS_ElemIteratorPtr elemIt = myMesh->elementsIterator(GetType());
309     while ( elemIt->more() )
310       values.insert( GetValue( elemIt->next()->GetID() ));
311   }
312   else
313   {
314     vector<int>::const_iterator id = elements.begin();
315     for ( ; id != elements.end(); ++id )
316       values.insert( GetValue( *id ));
317   }
318
319   if ( minmax )
320   {
321     funValues[0] = minmax[0];
322     funValues[nbIntervals] = minmax[1];
323   }
324   else
325   {
326     funValues[0] = *values.begin();
327     funValues[nbIntervals] = *values.rbegin();
328   }
329   // case nbIntervals == 1
330   if ( nbIntervals == 1 )
331   {
332     nbEvents[0] = values.size();
333     return;
334   }
335   // case of 1 value
336   if (funValues.front() == funValues.back())
337   {
338     nbEvents.resize( 1 );
339     nbEvents[0] = values.size();
340     funValues[1] = funValues.back();
341     funValues.resize( 2 );
342   }
343   // generic case
344   std::multiset< double >::iterator min = values.begin(), max;
345   for ( int i = 0; i < nbIntervals; ++i )
346   {
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;
350
351     // count values in the i-th interval if there are any
352     if ( min != values.end() && *min <= funValues[i+1] )
353     {
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 );
357       min = max;
358     }
359   }
360   // add values larger than minmax[1]
361   nbEvents.back() += std::distance( min, values.end() );
362 }
363
364 //=======================================================================
365 //function : GetValue
366 //purpose  : 
367 //=======================================================================
368
369 double Volume::GetValue( long theElementId )
370 {
371   if ( theElementId && myMesh ) {
372     SMDS_VolumeTool aVolumeTool;
373     if ( aVolumeTool.Set( myMesh->FindElement( theElementId )))
374       return aVolumeTool.GetSize();
375   }
376   return 0;
377 }
378
379 //=======================================================================
380 //function : GetBadRate
381 //purpose  : meaningless as it is not quality control functor
382 //=======================================================================
383
384 double Volume::GetBadRate( double Value, int /*nbNodes*/ ) const
385 {
386   return Value;
387 }
388
389 //=======================================================================
390 //function : GetType
391 //purpose  : 
392 //=======================================================================
393
394 SMDSAbs_ElementType Volume::GetType() const
395 {
396   return SMDSAbs_Volume;
397 }
398
399
400 /*
401   Class       : MaxElementLength2D
402   Description : Functor calculating maximum length of 2D element
403 */
404
405 double MaxElementLength2D::GetValue( long theElementId )
406 {
407   TSequenceOfXYZ P;
408   if( GetPoints( theElementId, P ) ) {
409     double aVal = 0;
410     const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
411     SMDSAbs_ElementType aType = aElem->GetType();
412     int len = P.size();
413     switch( aType ) {
414     case SMDSAbs_Face:
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));
420         break;
421       }
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));
430         break;
431       }
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));
437         break;
438       }
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));
447         break;
448       }
449     }
450
451     if( myPrecision >= 0 )
452     {
453       double prec = pow( 10., (double)myPrecision );
454       aVal = floor( aVal * prec + 0.5 ) / prec;
455     }
456     return aVal;
457   }
458   return 0.;
459 }
460
461 double MaxElementLength2D::GetBadRate( double Value, int /*nbNodes*/ ) const
462 {
463   return Value;
464 }
465
466 SMDSAbs_ElementType MaxElementLength2D::GetType() const
467 {
468   return SMDSAbs_Face;
469 }
470
471 /*
472   Class       : MaxElementLength3D
473   Description : Functor calculating maximum length of 3D element
474 */
475
476 double MaxElementLength3D::GetValue( long theElementId )
477 {
478   TSequenceOfXYZ P;
479   if( GetPoints( theElementId, P ) ) {
480     double aVal = 0;
481     const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
482     SMDSAbs_ElementType aType = aElem->GetType();
483     int len = P.size();
484     switch( aType ) {
485     case SMDSAbs_Volume:
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));
494         break;
495       }
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));
507         break;
508       }
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));
521         break;
522       }
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)));
544         break;
545       }
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));
554         break;
555       }
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));
567         break;
568       }
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));
581         break;
582       }
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)));
604         break;
605       }
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 );
613             }
614           }
615         }
616       }
617     }
618
619     if( myPrecision >= 0 )
620     {
621       double prec = pow( 10., (double)myPrecision );
622       aVal = floor( aVal * prec + 0.5 ) / prec;
623     }
624     return aVal;
625   }
626   return 0.;
627 }
628
629 double MaxElementLength3D::GetBadRate( double Value, int /*nbNodes*/ ) const
630 {
631   return Value;
632 }
633
634 SMDSAbs_ElementType MaxElementLength3D::GetType() const
635 {
636   return SMDSAbs_Volume;
637 }
638
639
640 /*
641   Class       : MinimumAngle
642   Description : Functor for calculation of minimum angle
643 */
644
645 double MinimumAngle::GetValue( const TSequenceOfXYZ& P )
646 {
647   double aMin;
648
649   if (P.size() <3)
650     return 0.;
651
652   aMin = getAngle(P( P.size() ), P( 1 ), P( 2 ));
653   aMin = Min(aMin,getAngle(P( P.size()-1 ), P( P.size() ), P( 1 )));
654
655   for (int i=2; i<P.size();i++){
656       double A0 = getAngle( P( i-1 ), P( i ), P( i+1 ) );
657     aMin = Min(aMin,A0);
658   }
659
660   return aMin * 180.0 / PI;
661 }
662
663 double MinimumAngle::GetBadRate( double Value, int nbNodes ) const
664 {
665   //const double aBestAngle = PI / nbNodes;
666   const double aBestAngle = 180.0 - ( 360.0 / double(nbNodes) );
667   return ( fabs( aBestAngle - Value ));
668 }
669
670 SMDSAbs_ElementType MinimumAngle::GetType() const
671 {
672   return SMDSAbs_Face;
673 }
674
675
676 /*
677   Class       : AspectRatio
678   Description : Functor for calculating aspect ratio
679 */
680 double AspectRatio::GetValue( const TSequenceOfXYZ& P )
681 {
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
685   // PAL10872
686
687   int nbNodes = P.size();
688
689   if ( nbNodes < 3 )
690     return 0;
691
692   // Compute aspect ratio
693
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
701     //
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() )
711       return 0.;
712     return alfa * maxLen * half_perimeter / anArea;
713   }
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
721     //
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() )
731       return 0.;
732     return alfa * maxLen * half_perimeter / anArea;
733   }
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
753     //
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 ],
763                  Max( aLen[ 1 ],
764                    Max( aLen[ 2 ],
765                      Max( aLen[ 3 ],
766                        Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
767     double C1 = sqrt( ( aLen[0] * aLen[0] +
768                         aLen[1] * aLen[1] +
769                         aLen[2] * aLen[2] +
770                         aLen[3] * aLen[3] ) / 4. );
771     double C2 = Min( anArea[ 0 ],
772                   Min( anArea[ 1 ],
773                     Min( anArea[ 2 ], anArea[ 3 ] ) ) );
774     if ( C2 <= Precision::Confusion() )
775       return 0.;
776     return alpha * L * C1 / C2;
777   }
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
797     //
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 ],
807                  Max( aLen[ 1 ],
808                    Max( aLen[ 2 ],
809                      Max( aLen[ 3 ],
810                        Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
811     double C1 = sqrt( ( aLen[0] * aLen[0] +
812                         aLen[1] * aLen[1] +
813                         aLen[2] * aLen[2] +
814                         aLen[3] * aLen[3] ) / 4. );
815     double C2 = Min( anArea[ 0 ],
816                   Min( anArea[ 1 ],
817                     Min( anArea[ 2 ], anArea[ 3 ] ) ) );
818     if ( C2 <= Precision::Confusion() )
819       return 0.;
820     return alpha * L * C1 / C2;
821   }
822   return 0;
823 }
824
825 double AspectRatio::GetBadRate( double Value, int /*nbNodes*/ ) const
826 {
827   // the aspect ratio is in the range [1.0,infinity]
828   // 1.0 = good
829   // infinity = bad
830   return Value / 1000.;
831 }
832
833 SMDSAbs_ElementType AspectRatio::GetType() const
834 {
835   return SMDSAbs_Face;
836 }
837
838
839 /*
840   Class       : AspectRatio3D
841   Description : Functor for calculating aspect ratio
842 */
843 namespace{
844
845   inline double getHalfPerimeter(double theTria[3]){
846     return (theTria[0] + theTria[1] + theTria[2])/2.0;
847   }
848
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]));
854   }
855
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;
867   }
868
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];
876
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;
881
882     return sqrt(P+Q+R-S)/12.0;
883   }
884
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;
891   }
892
893   inline double getMaxHeight(double theLen[6])
894   {
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]);
900     return aHeight;
901   }
902
903 }
904
905 double AspectRatio3D::GetValue( const TSequenceOfXYZ& P )
906 {
907   double aQuality = 0.0;
908   if(myCurrElement->IsPoly()) return aQuality;
909
910   int nbNodes = P.size();
911
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;
918   }
919
920   switch(nbNodes){
921   case 4:{
922     double aLen[6] = {
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
929     };
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
935     };
936     double aSumArea = 0.0;
937     double aHalfPerimeter = getHalfPerimeter(aTria[0]);
938     double anArea = getArea(aHalfPerimeter,aTria[0]);
939     aSumArea += anArea;
940     aHalfPerimeter = getHalfPerimeter(aTria[1]);
941     anArea = getArea(aHalfPerimeter,aTria[1]);
942     aSumArea += anArea;
943     aHalfPerimeter = getHalfPerimeter(aTria[2]);
944     anArea = getArea(aHalfPerimeter,aTria[2]);
945     aSumArea += anArea;
946     aHalfPerimeter = getHalfPerimeter(aTria[3]);
947     anArea = getArea(aHalfPerimeter,aTria[3]);
948     aSumArea += anArea;
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;
955     break;
956   }
957   case 5:{
958     {
959       gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 3 ),P( 5 )};
960       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
961     }
962     {
963       gp_XYZ aXYZ[4] = {P( 1 ),P( 3 ),P( 4 ),P( 5 )};
964       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
965     }
966     {
967       gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 5 )};
968       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
969     }
970     {
971       gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 4 ),P( 5 )};
972       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
973     }
974     break;
975   }
976   case 6:{
977     {
978       gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 6 )};
979       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
980     }
981     {
982       gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 3 )};
983       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
984     }
985     {
986       gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 6 )};
987       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
988     }
989     {
990       gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
991       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
992     }
993     {
994       gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 6 )};
995       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
996     }
997     {
998       gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 3 )};
999       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1000     }
1001     break;
1002   }
1003   case 8:{
1004     {
1005       gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1006       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1007     }
1008     {
1009       gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 4 )};
1010       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1011     }
1012     {
1013       gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 7 )};
1014       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1015     }
1016     {
1017       gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 8 )};
1018       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1019     }
1020     {
1021       gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 3 )};
1022       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1023     }
1024     {
1025       gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 4 )};
1026       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1027     }
1028     {
1029       gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 7 )};
1030       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1031     }
1032     {
1033       gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 8 )};
1034       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1035     }
1036     {
1037       gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 3 )};
1038       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1039     }
1040     {
1041       gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 4 )};
1042       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1043     }
1044     {
1045       gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 7 )};
1046       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1047     }
1048     {
1049       gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 8 )};
1050       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1051     }
1052     {
1053       gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 1 )};
1054       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1055     }
1056     {
1057       gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 2 )};
1058       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1059     }
1060     {
1061       gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 5 )};
1062       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1063     }
1064     {
1065       gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 6 )};
1066       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1067     }
1068     {
1069       gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 1 )};
1070       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1071     }
1072     {
1073       gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 2 )};
1074       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1075     }
1076     {
1077       gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 5 )};
1078       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1079     }
1080     {
1081       gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 6 )};
1082       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1083     }
1084     {
1085       gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 1 )};
1086       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1087     }
1088     {
1089       gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1090       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1091     }
1092     {
1093       gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 5 )};
1094       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1095     }
1096     {
1097       gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 6 )};
1098       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1099     }
1100     {
1101       gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1102       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1103     }
1104     {
1105       gp_XYZ aXYZ[4] = {P( 4 ),P( 5 ),P( 8 ),P( 2 )};
1106       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1107     }
1108     {
1109       gp_XYZ aXYZ[4] = {P( 1 ),P( 4 ),P( 5 ),P( 3 )};
1110       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1111     }
1112     {
1113       gp_XYZ aXYZ[4] = {P( 3 ),P( 6 ),P( 7 ),P( 1 )};
1114       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1115     }
1116     {
1117       gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 6 ),P( 4 )};
1118       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1119     }
1120     {
1121       gp_XYZ aXYZ[4] = {P( 5 ),P( 6 ),P( 8 ),P( 3 )};
1122       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1123     }
1124     {
1125       gp_XYZ aXYZ[4] = {P( 7 ),P( 8 ),P( 6 ),P( 1 )};
1126       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1127     }
1128     {
1129       gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 7 )};
1130       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1131     }
1132     {
1133       gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 2 ),P( 5 )};
1134       aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1135     }
1136     break;
1137   }
1138   }
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 )
1147         continue;
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 ));
1152     }
1153   }
1154   return aQuality;
1155 }
1156
1157 double AspectRatio3D::GetBadRate( double Value, int /*nbNodes*/ ) const
1158 {
1159   // the aspect ratio is in the range [1.0,infinity]
1160   // 1.0 = good
1161   // infinity = bad
1162   return Value / 1000.;
1163 }
1164
1165 SMDSAbs_ElementType AspectRatio3D::GetType() const
1166 {
1167   return SMDSAbs_Volume;
1168 }
1169
1170
1171 /*
1172   Class       : Warping
1173   Description : Functor for calculating warping
1174 */
1175 double Warping::GetValue( const TSequenceOfXYZ& P )
1176 {
1177   if ( P.size() != 4 )
1178     return 0;
1179
1180   gp_XYZ G = ( P( 1 ) + P( 2 ) + P( 3 ) + P( 4 ) ) / 4.;
1181
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 );
1186
1187   return Max( Max( A1, A2 ), Max( A3, A4 ) );
1188 }
1189
1190 double Warping::ComputeA( const gp_XYZ& thePnt1,
1191                           const gp_XYZ& thePnt2,
1192                           const gp_XYZ& thePnt3,
1193                           const gp_XYZ& theG ) const
1194 {
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())
1199     return 0.;
1200
1201   gp_XYZ GI = ( thePnt2 + thePnt1 ) / 2. - theG;
1202   gp_XYZ GJ = ( thePnt3 + thePnt2 ) / 2. - theG;
1203   gp_XYZ N  = GI.Crossed( GJ );
1204
1205   if ( N.Modulus() < gp::Resolution() )
1206     return PI / 2;
1207
1208   N.Normalize();
1209
1210   double H = ( thePnt2 - theG ).Dot( N );
1211   return asin( fabs( H / L ) ) * 180. / PI;
1212 }
1213
1214 double Warping::GetBadRate( double Value, int /*nbNodes*/ ) const
1215 {
1216   // the warp is in the range [0.0,PI/2]
1217   // 0.0 = good (no warp)
1218   // PI/2 = bad  (face pliee)
1219   return Value;
1220 }
1221
1222 SMDSAbs_ElementType Warping::GetType() const
1223 {
1224   return SMDSAbs_Face;
1225 }
1226
1227
1228 /*
1229   Class       : Taper
1230   Description : Functor for calculating taper
1231 */
1232 double Taper::GetValue( const TSequenceOfXYZ& P )
1233 {
1234   if ( P.size() != 4 )
1235     return 0.;
1236
1237   // Compute taper
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.;
1242
1243   double JA = 0.25 * ( J1 + J2 + J3 + J4 );
1244   if ( JA <= Precision::Confusion() )
1245     return 0.;
1246
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 );
1251
1252   return Max( Max( T1, T2 ), Max( T3, T4 ) );
1253 }
1254
1255 double Taper::GetBadRate( double Value, int /*nbNodes*/ ) const
1256 {
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)
1260   return Value;
1261 }
1262
1263 SMDSAbs_ElementType Taper::GetType() const
1264 {
1265   return SMDSAbs_Face;
1266 }
1267
1268
1269 /*
1270   Class       : Skew
1271   Description : Functor for calculating skew in degrees
1272 */
1273 static inline double skewAngle( const gp_XYZ& p1, const gp_XYZ& p2, const gp_XYZ& p3 )
1274 {
1275   gp_XYZ p12 = ( p2 + p1 ) / 2.;
1276   gp_XYZ p23 = ( p3 + p2 ) / 2.;
1277   gp_XYZ p31 = ( p3 + p1 ) / 2.;
1278
1279   gp_Vec v1( p31 - p2 ), v2( p12 - p23 );
1280
1281   return v1.Magnitude() < gp::Resolution() || v2.Magnitude() < gp::Resolution() ? 0. : v1.Angle( v2 );
1282 }
1283
1284 double Skew::GetValue( const TSequenceOfXYZ& P )
1285 {
1286   if ( P.size() != 3 && P.size() != 4 )
1287     return 0.;
1288
1289   // Compute skew
1290   static double PI2 = PI / 2.;
1291   if ( P.size() == 3 )
1292   {
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 ) ) );
1296
1297     return Max( A0, Max( A1, A2 ) ) * 180. / PI;
1298   }
1299   else
1300   {
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.;
1305
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 ) );
1309
1310     //BUG SWP12743
1311     if ( A < Precision::Angular() )
1312       return 0.;
1313
1314     return A * 180. / PI;
1315   }
1316 }
1317
1318 double Skew::GetBadRate( double Value, int /*nbNodes*/ ) const
1319 {
1320   // the skew is in the range [0.0,PI/2].
1321   // 0.0 = good
1322   // PI/2 = bad
1323   return Value;
1324 }
1325
1326 SMDSAbs_ElementType Skew::GetType() const
1327 {
1328   return SMDSAbs_Face;
1329 }
1330
1331
1332 /*
1333   Class       : Area
1334   Description : Functor for calculating area
1335 */
1336 double Area::GetValue( const TSequenceOfXYZ& P )
1337 {
1338   double val = 0.0;
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;
1347       SumVec.Add(tmp);
1348     }
1349     val = SumVec.Magnitude() * 0.5;
1350   }
1351   return val;
1352 }
1353
1354 double Area::GetBadRate( double Value, int /*nbNodes*/ ) const
1355 {
1356   // meaningless as it is not a quality control functor
1357   return Value;
1358 }
1359
1360 SMDSAbs_ElementType Area::GetType() const
1361 {
1362   return SMDSAbs_Face;
1363 }
1364
1365
1366 /*
1367   Class       : Length
1368   Description : Functor for calculating length off edge
1369 */
1370 double Length::GetValue( const TSequenceOfXYZ& P )
1371 {
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 ) );
1375   default: return 0.;
1376   }
1377 }
1378
1379 double Length::GetBadRate( double Value, int /*nbNodes*/ ) const
1380 {
1381   // meaningless as it is not quality control functor
1382   return Value;
1383 }
1384
1385 SMDSAbs_ElementType Length::GetType() const
1386 {
1387   return SMDSAbs_Edge;
1388 }
1389
1390 /*
1391   Class       : Length2D
1392   Description : Functor for calculating length of edge
1393 */
1394
1395 double Length2D::GetValue( long theElementId)
1396 {
1397   TSequenceOfXYZ P;
1398
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;
1403
1404     double  aVal;// = GetValue( P );
1405     const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
1406     SMDSAbs_ElementType aType = aElem->GetType();
1407
1408     int len = P.size();
1409
1410     switch (aType){
1411     case SMDSAbs_All:
1412     case SMDSAbs_Node:
1413     case SMDSAbs_Edge:
1414       if (len == 2){
1415         aVal = getDistance( P( 1 ), P( 2 ) );
1416         break;
1417       }
1418       else if (len == 3){ // quadratic edge
1419         aVal = getDistance(P( 1 ),P( 3 )) + getDistance(P( 3 ),P( 2 ));
1420         break;
1421       }
1422     case SMDSAbs_Face:
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));
1428         break;
1429       }
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));
1436         break;
1437       }
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;
1444         break;
1445       }
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));
1452         break;
1453       }
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));
1463         break;
1464       }
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 ));
1474
1475         aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1476         aVal = Max(aVal,Max(L7,L8));
1477         break;
1478       }
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 ));
1489
1490         aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1491         aVal = Max(aVal,Max(Max(L7,L8),L9));
1492         break;
1493       }
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 ));
1507
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));
1511         break;
1512
1513       }
1514
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));
1523         break;
1524       }
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));
1536         break;
1537       }
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));
1550         break;
1551       }
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));
1568         break;
1569
1570       }
1571
1572     default: aVal=-1;
1573     }
1574
1575     if (aVal <0){
1576       return 0.;
1577     }
1578
1579     if ( myPrecision >= 0 )
1580     {
1581       double prec = pow( 10., (double)( myPrecision ) );
1582       aVal = floor( aVal * prec + 0.5 ) / prec;
1583     }
1584
1585     return aVal;
1586
1587   }
1588   return 0.;
1589 }
1590
1591 double Length2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1592 {
1593   // meaningless as it is not quality control functor
1594   return Value;
1595 }
1596
1597 SMDSAbs_ElementType Length2D::GetType() const
1598 {
1599   return SMDSAbs_Face;
1600 }
1601
1602 Length2D::Value::Value(double theLength,long thePntId1, long thePntId2):
1603   myLength(theLength)
1604 {
1605   myPntId[0] = thePntId1;  myPntId[1] = thePntId2;
1606   if(thePntId1 > thePntId2){
1607     myPntId[1] = thePntId1;  myPntId[0] = thePntId2;
1608   }
1609 }
1610
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;
1615   return false;
1616 }
1617
1618 void Length2D::GetValues(TValues& theValues){
1619   TValues aValues;
1620   SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1621   for(; anIter->more(); ){
1622     const SMDS_MeshFace* anElem = anIter->next();
1623
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();
1629       long aNodeId[4];
1630       gp_Pnt P[4];
1631
1632       double aLength;
1633       const SMDS_MeshElement* aNode;
1634       if(anIter->more()){
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();
1639         aLength = 0;
1640       }
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]);
1653         P[1] = P[3];
1654         aNodeId[1] = aNodeId[3];
1655         theValues.insert(aValue1);
1656         theValues.insert(aValue2);
1657       }
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);
1663     }
1664     else {
1665       SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
1666       long aNodeId[2];
1667       gp_Pnt P[3];
1668
1669       double aLength;
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();
1676         aLength = 0;
1677       }
1678       for(; aNodesIter->more(); ){
1679         aNode = aNodesIter->next();
1680         const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1681         long anId = aNode->GetID();
1682         
1683         P[2] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1684         
1685         aLength = P[1].Distance(P[2]);
1686         
1687         Value aValue(aLength,aNodeId[1],anId);
1688         aNodeId[1] = anId;
1689         P[1] = P[2];
1690         theValues.insert(aValue);
1691       }
1692
1693       aLength = P[0].Distance(P[1]);
1694
1695       Value aValue(aLength,aNodeId[0],aNodeId[1]);
1696       theValues.insert(aValue);
1697     }
1698   }
1699 }
1700
1701 /*
1702   Class       : MultiConnection
1703   Description : Functor for calculating number of faces conneted to the edge
1704 */
1705 double MultiConnection::GetValue( const TSequenceOfXYZ& P )
1706 {
1707   return 0;
1708 }
1709 double MultiConnection::GetValue( long theId )
1710 {
1711   return getNbMultiConnection( myMesh, theId );
1712 }
1713
1714 double MultiConnection::GetBadRate( double Value, int /*nbNodes*/ ) const
1715 {
1716   // meaningless as it is not quality control functor
1717   return Value;
1718 }
1719
1720 SMDSAbs_ElementType MultiConnection::GetType() const
1721 {
1722   return SMDSAbs_Edge;
1723 }
1724
1725 /*
1726   Class       : MultiConnection2D
1727   Description : Functor for calculating number of faces conneted to the edge
1728 */
1729 double MultiConnection2D::GetValue( const TSequenceOfXYZ& P )
1730 {
1731   return 0;
1732 }
1733
1734 double MultiConnection2D::GetValue( long theElementId )
1735 {
1736   int aResult = 0;
1737
1738   const SMDS_MeshElement* aFaceElem = myMesh->FindElement(theElementId);
1739   SMDSAbs_ElementType aType = aFaceElem->GetType();
1740
1741   switch (aType) {
1742   case SMDSAbs_Face:
1743     {
1744       int i = 0, len = aFaceElem->NbNodes();
1745       SMDS_ElemIteratorPtr anIter = aFaceElem->nodesIterator();
1746       if (!anIter) break;
1747
1748       const SMDS_MeshNode *aNode, *aNode0;
1749       TColStd_MapOfInteger aMap, aMapPrev;
1750
1751       for (i = 0; i <= len; i++) {
1752         aMapPrev = aMap;
1753         aMap.Clear();
1754
1755         int aNb = 0;
1756         if (anIter->more()) {
1757           aNode = (SMDS_MeshNode*)anIter->next();
1758         } else {
1759           if (i == len)
1760             aNode = aNode0;
1761           else
1762             break;
1763         }
1764         if (!aNode) break;
1765         if (i == 0) aNode0 = aNode;
1766
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();
1772
1773             aMap.Add(anId);
1774             if (aMapPrev.Contains(anId)) {
1775               aNb++;
1776             }
1777           }
1778         }
1779         aResult = Max(aResult, aNb);
1780       }
1781     }
1782     break;
1783   default:
1784     aResult = 0;
1785   }
1786
1787   return aResult;
1788 }
1789
1790 double MultiConnection2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1791 {
1792   // meaningless as it is not quality control functor
1793   return Value;
1794 }
1795
1796 SMDSAbs_ElementType MultiConnection2D::GetType() const
1797 {
1798   return SMDSAbs_Face;
1799 }
1800
1801 MultiConnection2D::Value::Value(long thePntId1, long thePntId2)
1802 {
1803   myPntId[0] = thePntId1;  myPntId[1] = thePntId2;
1804   if(thePntId1 > thePntId2){
1805     myPntId[1] = thePntId1;  myPntId[0] = thePntId2;
1806   }
1807 }
1808
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;
1813   return false;
1814 }
1815
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();
1824     else
1825       aNodesIter = anElem->nodesIterator();
1826     long aNodeId[3];
1827
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();
1834       aNode1 = aNode0;
1835       const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode1;
1836       aNodeId[0] = aNodeId[1] = aNodes->GetID();
1837     }
1838     for(; aNodesIter->more(); ) {
1839       aNode2 = (SMDS_MeshNode*) aNodesIter->next();
1840       long anId = aNode2->GetID();
1841       aNodeId[2] = anId;
1842
1843       Value aValue(aNodeId[1],aNodeId[2]);
1844       MValues::iterator aItr = theValues.find(aValue);
1845       if (aItr != theValues.end()){
1846         aItr->second += 1;
1847         //aNbConnects = nb;
1848       }
1849       else {
1850         theValues[aValue] = 1;
1851         //aNbConnects = 1;
1852       }
1853       //cout << "NodeIds: "<<aNodeId[1]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1854       aNodeId[1] = aNodeId[2];
1855       aNode1 = aNode2;
1856     }
1857     Value aValue(aNodeId[0],aNodeId[2]);
1858     MValues::iterator aItr = theValues.find(aValue);
1859     if (aItr != theValues.end()) {
1860       aItr->second += 1;
1861       //aNbConnects = nb;
1862     }
1863     else {
1864       theValues[aValue] = 1;
1865       //aNbConnects = 1;
1866     }
1867     //cout << "NodeIds: "<<aNodeId[0]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1868   }
1869
1870 }
1871
1872 /*
1873                             PREDICATES
1874 */
1875
1876 /*
1877   Class       : BadOrientedVolume
1878   Description : Predicate bad oriented volumes
1879 */
1880
1881 BadOrientedVolume::BadOrientedVolume()
1882 {
1883   myMesh = 0;
1884 }
1885
1886 void BadOrientedVolume::SetMesh( const SMDS_Mesh* theMesh )
1887 {
1888   myMesh = theMesh;
1889 }
1890
1891 bool BadOrientedVolume::IsSatisfy( long theId )
1892 {
1893   if ( myMesh == 0 )
1894     return false;
1895
1896   SMDS_VolumeTool vTool( myMesh->FindElement( theId ));
1897   return !vTool.IsForward();
1898 }
1899
1900 SMDSAbs_ElementType BadOrientedVolume::GetType() const
1901 {
1902   return SMDSAbs_Volume;
1903 }
1904
1905 /*
1906   Class       : BareBorderVolume
1907 */
1908
1909 bool BareBorderVolume::IsSatisfy(long theElementId )
1910 {
1911   SMDS_VolumeTool  myTool;
1912   if ( myTool.Set( myMesh->FindElement(theElementId)))
1913   {
1914     for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
1915       if ( myTool.IsFreeFace( iF ))
1916       {
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))
1920           return true;
1921       }
1922   }
1923   return false;
1924 }
1925
1926 /*
1927   Class       : BareBorderFace
1928 */
1929
1930 bool BareBorderFace::IsSatisfy(long theElementId )
1931 {
1932   if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
1933     if ( face->GetType() == SMDSAbs_Face )
1934     {
1935       int nbN = face->NbCornerNodes();
1936       for ( int i = 0; i < nbN; ++i )
1937       {
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() )
1944         {
1945           const SMDS_MeshElement* f = fIt->next();
1946           isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
1947         }
1948         if ( !isShared )
1949         {
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);
1956         }
1957       }
1958     }
1959   return false;
1960 }
1961
1962 /*
1963   Class       : OverConstrainedVolume
1964 */
1965
1966 bool OverConstrainedVolume::IsSatisfy(long theElementId )
1967 {
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)))
1972   {
1973     int nbSharedFaces = 0;
1974     for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
1975       if ( !myTool.IsFreeFace( iF ) && ++nbSharedFaces > 1 )
1976         break;
1977     return ( nbSharedFaces == 1 );
1978   }
1979   return false;
1980 }
1981
1982 /*
1983   Class       : OverConstrainedFace
1984 */
1985
1986 bool OverConstrainedFace::IsSatisfy(long theElementId )
1987 {
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 )
1992     {
1993       int nbSharedBorders = 0;
1994       int nbN = face->NbCornerNodes();
1995       for ( int i = 0; i < nbN; ++i )
1996       {
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() )
2003         {
2004           const SMDS_MeshElement* f = fIt->next();
2005           isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2006         }
2007         if ( isShared && ++nbSharedBorders > 1 )
2008           break;
2009       }
2010       return ( nbSharedBorders == 1 );
2011     }
2012   return false;
2013 }
2014
2015 /*
2016   Class       : FreeBorders
2017   Description : Predicate for free borders
2018 */
2019
2020 FreeBorders::FreeBorders()
2021 {
2022   myMesh = 0;
2023 }
2024
2025 void FreeBorders::SetMesh( const SMDS_Mesh* theMesh )
2026 {
2027   myMesh = theMesh;
2028 }
2029
2030 bool FreeBorders::IsSatisfy( long theId )
2031 {
2032   return getNbMultiConnection( myMesh, theId ) == 1;
2033 }
2034
2035 SMDSAbs_ElementType FreeBorders::GetType() const
2036 {
2037   return SMDSAbs_Edge;
2038 }
2039
2040
2041 /*
2042   Class       : FreeEdges
2043   Description : Predicate for free Edges
2044 */
2045 FreeEdges::FreeEdges()
2046 {
2047   myMesh = 0;
2048 }
2049
2050 void FreeEdges::SetMesh( const SMDS_Mesh* theMesh )
2051 {
2052   myMesh = theMesh;
2053 }
2054
2055 bool FreeEdges::IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId  )
2056 {
2057   TColStd_MapOfInteger aMap;
2058   for ( int i = 0; i < 2; i++ )
2059   {
2060     SMDS_ElemIteratorPtr anElemIter = theNodes[ i ]->GetInverseElementIterator();
2061     while( anElemIter->more() )
2062     {
2063       const SMDS_MeshElement* anElem = anElemIter->next();
2064       if ( anElem != 0 && anElem->GetType() == SMDSAbs_Face )
2065       {
2066         int anId = anElem->GetID();
2067
2068         if ( i == 0 )
2069           aMap.Add( anId );
2070         else if ( aMap.Contains( anId ) && anId != theFaceId )
2071           return false;
2072       }
2073     }
2074   }
2075   return true;
2076 }
2077
2078 bool FreeEdges::IsSatisfy( long theId )
2079 {
2080   if ( myMesh == 0 )
2081     return false;
2082
2083   const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2084   if ( aFace == 0 || aFace->GetType() != SMDSAbs_Face || aFace->NbNodes() < 3 )
2085     return false;
2086
2087   SMDS_ElemIteratorPtr anIter;
2088   if ( aFace->IsQuadratic() ) {
2089     anIter = static_cast<const SMDS_QuadraticFaceOfNodes*>
2090       (aFace)->interlacedNodesElemIterator();
2091   }
2092   else {
2093     anIter = aFace->nodesIterator();
2094   }
2095   if ( anIter == 0 )
2096     return false;
2097
2098   int i = 0, nbNodes = aFace->NbNodes();
2099   std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
2100   while( anIter->more() )
2101   {
2102     const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
2103     if ( aNode == 0 )
2104       return false;
2105     aNodes[ i++ ] = aNode;
2106   }
2107   aNodes[ nbNodes ] = aNodes[ 0 ];
2108
2109   for ( i = 0; i < nbNodes; i++ )
2110     if ( IsFreeEdge( &aNodes[ i ], theId ) )
2111       return true;
2112
2113   return false;
2114 }
2115
2116 SMDSAbs_ElementType FreeEdges::GetType() const
2117 {
2118   return SMDSAbs_Face;
2119 }
2120
2121 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2122   myElemId(theElemId)
2123 {
2124   myPntId[0] = thePntId1;  myPntId[1] = thePntId2;
2125   if(thePntId1 > thePntId2){
2126     myPntId[1] = thePntId1;  myPntId[0] = thePntId2;
2127   }
2128 }
2129
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;
2134   return false;
2135 }
2136
2137 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2138                           FreeEdges::TBorders& theRegistry,
2139                           FreeEdges::TBorders& theContainer)
2140 {
2141   if(theRegistry.find(theBorder) == theRegistry.end()){
2142     theRegistry.insert(theBorder);
2143     theContainer.insert(theBorder);
2144   }else{
2145     theContainer.erase(theBorder);
2146   }
2147 }
2148
2149 void FreeEdges::GetBoreders(TBorders& theBorders)
2150 {
2151   TBorders aRegistry;
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();
2160     else
2161       aNodesIter = anElem->nodesIterator();
2162     long aNodeId[2];
2163     const SMDS_MeshElement* aNode;
2164     if(aNodesIter->more()){
2165       aNode = aNodesIter->next();
2166       aNodeId[0] = aNodeId[1] = aNode->GetID();
2167     }
2168     for(; aNodesIter->more(); ){
2169       aNode = aNodesIter->next();
2170       long anId = aNode->GetID();
2171       Border aBorder(anElemId,aNodeId[1],anId);
2172       aNodeId[1] = anId;
2173       //std::cout<<aBorder.myPntId[0]<<"; "<<aBorder.myPntId[1]<<"; "<<aBorder.myElemId<<endl;
2174       UpdateBorders(aBorder,aRegistry,theBorders);
2175     }
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);
2179   }
2180   //std::cout<<"theBorders.size() = "<<theBorders.size()<<endl;
2181 }
2182
2183
2184 /*
2185   Class       : FreeNodes
2186   Description : Predicate for free nodes
2187 */
2188
2189 FreeNodes::FreeNodes()
2190 {
2191   myMesh = 0;
2192 }
2193
2194 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2195 {
2196   myMesh = theMesh;
2197 }
2198
2199 bool FreeNodes::IsSatisfy( long theNodeId )
2200 {
2201   const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2202   if (!aNode)
2203     return false;
2204
2205   return (aNode->NbInverseElements() < 1);
2206 }
2207
2208 SMDSAbs_ElementType FreeNodes::GetType() const
2209 {
2210   return SMDSAbs_Node;
2211 }
2212
2213
2214 /*
2215   Class       : FreeFaces
2216   Description : Predicate for free faces
2217 */
2218
2219 FreeFaces::FreeFaces()
2220 {
2221   myMesh = 0;
2222 }
2223
2224 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2225 {
2226   myMesh = theMesh;
2227 }
2228
2229 bool FreeFaces::IsSatisfy( long theId )
2230 {
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 )
2235     return false;
2236
2237   int nbNode = aFace->NbNodes();
2238
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;
2243
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;
2252       (*itr).second++;
2253     } 
2254   }
2255   int nbVol = 0;
2256   TItrMapOfVolume volItr = mapOfVol.begin();
2257   TItrMapOfVolume volEnd = mapOfVol.end();
2258   for ( ; volItr != volEnd; ++volItr )
2259     if ( (*volItr).second >= nbNode )
2260        nbVol++;
2261   // face is not free if number of volumes constructed on thier nodes more than one
2262   return (nbVol < 2);
2263 }
2264
2265 SMDSAbs_ElementType FreeFaces::GetType() const
2266 {
2267   return SMDSAbs_Face;
2268 }
2269
2270 /*
2271   Class       : LinearOrQuadratic
2272   Description : Predicate to verify whether a mesh element is linear
2273 */
2274
2275 LinearOrQuadratic::LinearOrQuadratic()
2276 {
2277   myMesh = 0;
2278 }
2279
2280 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2281 {
2282   myMesh = theMesh;
2283 }
2284
2285 bool LinearOrQuadratic::IsSatisfy( long theId )
2286 {
2287   if (!myMesh) return false;
2288   const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2289   if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2290     return false;
2291   return (!anElem->IsQuadratic());
2292 }
2293
2294 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2295 {
2296   myType = theType;
2297 }
2298
2299 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2300 {
2301   return myType;
2302 }
2303
2304 /*
2305   Class       : GroupColor
2306   Description : Functor for check color of group to whic mesh element belongs to
2307 */
2308
2309 GroupColor::GroupColor()
2310 {
2311 }
2312
2313 bool GroupColor::IsSatisfy( long theId )
2314 {
2315   return (myIDs.find( theId ) != myIDs.end());
2316 }
2317
2318 void GroupColor::SetType( SMDSAbs_ElementType theType )
2319 {
2320   myType = theType;
2321 }
2322
2323 SMDSAbs_ElementType GroupColor::GetType() const
2324 {
2325   return myType;
2326 }
2327
2328 static bool isEqual( const Quantity_Color& theColor1,
2329                      const Quantity_Color& theColor2 )
2330 {
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 );
2336 }
2337
2338
2339 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2340 {
2341   myIDs.clear();
2342   
2343   const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2344   if ( !aMesh )
2345     return;
2346
2347   int nbGrp = aMesh->GetNbGroups();
2348   if ( !nbGrp )
2349     return;
2350   
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);
2356     if ( !aGrp )
2357       continue;
2358     // check type and color of group
2359     if ( !isEqual( myColor, aGrp->GetColor() ) )
2360       continue;
2361     if ( myType != SMDSAbs_All && myType != (SMDSAbs_ElementType)aGrp->GetType() )
2362       continue;
2363
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) );
2370     }
2371   }
2372 }
2373
2374 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
2375 {
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();
2387   }
2388   myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
2389 }
2390
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
2397 {
2398   theResStr.Clear();
2399   theResStr += TCollection_AsciiString( myColor.Red() );
2400   theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2401   theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2402 }
2403
2404 /*
2405   Class       : ElemGeomType
2406   Description : Predicate to check element geometry type
2407 */
2408
2409 ElemGeomType::ElemGeomType()
2410 {
2411   myMesh = 0;
2412   myType = SMDSAbs_All;
2413   myGeomType = SMDSGeom_TRIANGLE;
2414 }
2415
2416 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2417 {
2418   myMesh = theMesh;
2419 }
2420
2421 bool ElemGeomType::IsSatisfy( long theId )
2422 {
2423   if (!myMesh) return false;
2424   const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2425   if ( !anElem )
2426     return false;
2427   const SMDSAbs_ElementType anElemType = anElem->GetType();
2428   if ( myType != SMDSAbs_All && anElemType != myType )
2429     return false;
2430   const int aNbNode = anElem->NbNodes();
2431   bool isOk = false;
2432   switch( anElemType )
2433   {
2434   case SMDSAbs_Node:
2435     isOk = (myGeomType == SMDSGeom_POINT);
2436     break;
2437
2438   case SMDSAbs_Edge:
2439     isOk = (myGeomType == SMDSGeom_EDGE);
2440     break;
2441
2442   case SMDSAbs_Face:
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();
2449     break;
2450
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();
2462     break;
2463     default: break;
2464   }
2465   return isOk;
2466 }
2467
2468 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2469 {
2470   myType = theType;
2471 }
2472
2473 SMDSAbs_ElementType ElemGeomType::GetType() const
2474 {
2475   return myType;
2476 }
2477
2478 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2479 {
2480   myGeomType = theType;
2481 }
2482
2483 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2484 {
2485   return myGeomType;
2486 }
2487
2488 /*
2489   Class       : RangeOfIds
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-"
2495 */
2496
2497 //=======================================================================
2498 // name    : RangeOfIds
2499 // Purpose : Constructor
2500 //=======================================================================
2501 RangeOfIds::RangeOfIds()
2502 {
2503   myMesh = 0;
2504   myType = SMDSAbs_All;
2505 }
2506
2507 //=======================================================================
2508 // name    : SetMesh
2509 // Purpose : Set mesh
2510 //=======================================================================
2511 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
2512 {
2513   myMesh = theMesh;
2514 }
2515
2516 //=======================================================================
2517 // name    : AddToRange
2518 // Purpose : Add ID to the range
2519 //=======================================================================
2520 bool RangeOfIds::AddToRange( long theEntityId )
2521 {
2522   myIds.Add( theEntityId );
2523   return true;
2524 }
2525
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 )
2532 {
2533   theResStr.Clear();
2534
2535   TColStd_SequenceOfInteger     anIntSeq;
2536   TColStd_SequenceOfAsciiString aStrSeq;
2537
2538   TColStd_MapIteratorOfMapOfInteger anIter( myIds );
2539   for ( ; anIter.More(); anIter.Next() )
2540   {
2541     int anId = anIter.Key();
2542     TCollection_AsciiString aStr( anId );
2543     anIntSeq.Append( anId );
2544     aStrSeq.Append( aStr );
2545   }
2546
2547   for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2548   {
2549     int aMinId = myMin( i );
2550     int aMaxId = myMax( i );
2551
2552     TCollection_AsciiString aStr;
2553     if ( aMinId != IntegerFirst() )
2554       aStr += aMinId;
2555
2556     aStr += "-";
2557
2558     if ( aMaxId != IntegerLast() )
2559       aStr += aMaxId;
2560
2561     // find position of the string in result sequence and insert string in it
2562     if ( anIntSeq.Length() == 0 )
2563     {
2564       anIntSeq.Append( aMinId );
2565       aStrSeq.Append( aStr );
2566     }
2567     else
2568     {
2569       if ( aMinId < anIntSeq.First() )
2570       {
2571         anIntSeq.Prepend( aMinId );
2572         aStrSeq.Prepend( aStr );
2573       }
2574       else if ( aMinId > anIntSeq.Last() )
2575       {
2576         anIntSeq.Append( aMinId );
2577         aStrSeq.Append( aStr );
2578       }
2579       else
2580         for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
2581           if ( aMinId < anIntSeq( j ) )
2582           {
2583             anIntSeq.InsertBefore( j, aMinId );
2584             aStrSeq.InsertBefore( j, aStr );
2585             break;
2586           }
2587     }
2588   }
2589
2590   if ( aStrSeq.Length() == 0 )
2591     return;
2592
2593   theResStr = aStrSeq( 1 );
2594   for ( int j = 2, k = aStrSeq.Length(); j <= k; j++  )
2595   {
2596     theResStr += ",";
2597     theResStr += aStrSeq( j );
2598   }
2599 }
2600
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 )
2607 {
2608   myMin.Clear();
2609   myMax.Clear();
2610   myIds.Clear();
2611
2612   TCollection_AsciiString aStr = theStr;
2613   aStr.RemoveAll( ' ' );
2614   aStr.RemoveAll( '\t' );
2615
2616   for ( int aPos = aStr.Search( ",," ); aPos != -1; aPos = aStr.Search( ",," ) )
2617     aStr.Remove( aPos, 2 );
2618
2619   TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
2620   int i = 1;
2621   while ( tmpStr != "" )
2622   {
2623     tmpStr = aStr.Token( ",", i++ );
2624     int aPos = tmpStr.Search( '-' );
2625
2626     if ( aPos == -1 )
2627     {
2628       if ( tmpStr.IsIntegerValue() )
2629         myIds.Add( tmpStr.IntegerValue() );
2630       else
2631         return false;
2632     }
2633     else
2634     {
2635       TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
2636       TCollection_AsciiString aMinStr = tmpStr;
2637
2638       while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
2639       while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
2640
2641       if ( !aMinStr.IsEmpty() && !aMinStr.IsIntegerValue() ||
2642            !aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue() )
2643         return false;
2644
2645       myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
2646       myMax.Append( aMaxStr.IsEmpty() ? IntegerLast()  : aMaxStr.IntegerValue() );
2647     }
2648   }
2649
2650   return true;
2651 }
2652
2653 //=======================================================================
2654 // name    : GetType
2655 // Purpose : Get type of supported entities
2656 //=======================================================================
2657 SMDSAbs_ElementType RangeOfIds::GetType() const
2658 {
2659   return myType;
2660 }
2661
2662 //=======================================================================
2663 // name    : SetType
2664 // Purpose : Set type of supported entities
2665 //=======================================================================
2666 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
2667 {
2668   myType = theType;
2669 }
2670
2671 //=======================================================================
2672 // name    : IsSatisfy
2673 // Purpose : Verify whether entity satisfies to this rpedicate
2674 //=======================================================================
2675 bool RangeOfIds::IsSatisfy( long theId )
2676 {
2677   if ( !myMesh )
2678     return false;
2679
2680   if ( myType == SMDSAbs_Node )
2681   {
2682     if ( myMesh->FindNode( theId ) == 0 )
2683       return false;
2684   }
2685   else
2686   {
2687     const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2688     if ( anElem == 0 || myType != anElem->GetType() && myType != SMDSAbs_All )
2689       return false;
2690   }
2691
2692   if ( myIds.Contains( theId ) )
2693     return true;
2694
2695   for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2696     if ( theId >= myMin( i ) && theId <= myMax( i ) )
2697       return true;
2698
2699   return false;
2700 }
2701
2702 /*
2703   Class       : Comparator
2704   Description : Base class for comparators
2705 */
2706 Comparator::Comparator():
2707   myMargin(0)
2708 {}
2709
2710 Comparator::~Comparator()
2711 {}
2712
2713 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
2714 {
2715   if ( myFunctor )
2716     myFunctor->SetMesh( theMesh );
2717 }
2718
2719 void Comparator::SetMargin( double theValue )
2720 {
2721   myMargin = theValue;
2722 }
2723
2724 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
2725 {
2726   myFunctor = theFunct;
2727 }
2728
2729 SMDSAbs_ElementType Comparator::GetType() const
2730 {
2731   return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
2732 }
2733
2734 double Comparator::GetMargin()
2735 {
2736   return myMargin;
2737 }
2738
2739
2740 /*
2741   Class       : LessThan
2742   Description : Comparator "<"
2743 */
2744 bool LessThan::IsSatisfy( long theId )
2745 {
2746   return myFunctor && myFunctor->GetValue( theId ) < myMargin;
2747 }
2748
2749
2750 /*
2751   Class       : MoreThan
2752   Description : Comparator ">"
2753 */
2754 bool MoreThan::IsSatisfy( long theId )
2755 {
2756   return myFunctor && myFunctor->GetValue( theId ) > myMargin;
2757 }
2758
2759
2760 /*
2761   Class       : EqualTo
2762   Description : Comparator "="
2763 */
2764 EqualTo::EqualTo():
2765   myToler(Precision::Confusion())
2766 {}
2767
2768 bool EqualTo::IsSatisfy( long theId )
2769 {
2770   return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
2771 }
2772
2773 void EqualTo::SetTolerance( double theToler )
2774 {
2775   myToler = theToler;
2776 }
2777
2778 double EqualTo::GetTolerance()
2779 {
2780   return myToler;
2781 }
2782
2783 /*
2784   Class       : LogicalNOT
2785   Description : Logical NOT predicate
2786 */
2787 LogicalNOT::LogicalNOT()
2788 {}
2789
2790 LogicalNOT::~LogicalNOT()
2791 {}
2792
2793 bool LogicalNOT::IsSatisfy( long theId )
2794 {
2795   return myPredicate && !myPredicate->IsSatisfy( theId );
2796 }
2797
2798 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
2799 {
2800   if ( myPredicate )
2801     myPredicate->SetMesh( theMesh );
2802 }
2803
2804 void LogicalNOT::SetPredicate( PredicatePtr thePred )
2805 {
2806   myPredicate = thePred;
2807 }
2808
2809 SMDSAbs_ElementType LogicalNOT::GetType() const
2810 {
2811   return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
2812 }
2813
2814
2815 /*
2816   Class       : LogicalBinary
2817   Description : Base class for binary logical predicate
2818 */
2819 LogicalBinary::LogicalBinary()
2820 {}
2821
2822 LogicalBinary::~LogicalBinary()
2823 {}
2824
2825 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
2826 {
2827   if ( myPredicate1 )
2828     myPredicate1->SetMesh( theMesh );
2829
2830   if ( myPredicate2 )
2831     myPredicate2->SetMesh( theMesh );
2832 }
2833
2834 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
2835 {
2836   myPredicate1 = thePredicate;
2837 }
2838
2839 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
2840 {
2841   myPredicate2 = thePredicate;
2842 }
2843
2844 SMDSAbs_ElementType LogicalBinary::GetType() const
2845 {
2846   if ( !myPredicate1 || !myPredicate2 )
2847     return SMDSAbs_All;
2848
2849   SMDSAbs_ElementType aType1 = myPredicate1->GetType();
2850   SMDSAbs_ElementType aType2 = myPredicate2->GetType();
2851
2852   return aType1 == aType2 ? aType1 : SMDSAbs_All;
2853 }
2854
2855
2856 /*
2857   Class       : LogicalAND
2858   Description : Logical AND
2859 */
2860 bool LogicalAND::IsSatisfy( long theId )
2861 {
2862   return
2863     myPredicate1 &&
2864     myPredicate2 &&
2865     myPredicate1->IsSatisfy( theId ) &&
2866     myPredicate2->IsSatisfy( theId );
2867 }
2868
2869
2870 /*
2871   Class       : LogicalOR
2872   Description : Logical OR
2873 */
2874 bool LogicalOR::IsSatisfy( long theId )
2875 {
2876   return
2877     myPredicate1 &&
2878     myPredicate2 &&
2879     myPredicate1->IsSatisfy( theId ) ||
2880     myPredicate2->IsSatisfy( theId );
2881 }
2882
2883
2884 /*
2885                               FILTER
2886 */
2887
2888 Filter::Filter()
2889 {}
2890
2891 Filter::~Filter()
2892 {}
2893
2894 void Filter::SetPredicate( PredicatePtr thePredicate )
2895 {
2896   myPredicate = thePredicate;
2897 }
2898
2899 template<class TElement, class TIterator, class TPredicate>
2900 inline void FillSequence(const TIterator& theIterator,
2901                          TPredicate& thePredicate,
2902                          Filter::TIdSequence& theSequence)
2903 {
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 );
2910     }
2911   }
2912 }
2913
2914 void
2915 Filter::
2916 GetElementsId( const SMDS_Mesh* theMesh,
2917                PredicatePtr thePredicate,
2918                TIdSequence& theSequence )
2919 {
2920   theSequence.clear();
2921
2922   if ( !theMesh || !thePredicate )
2923     return;
2924
2925   thePredicate->SetMesh( theMesh );
2926
2927   SMDSAbs_ElementType aType = thePredicate->GetType();
2928   switch(aType){
2929   case SMDSAbs_Node:
2930     FillSequence<const SMDS_MeshNode*>(theMesh->nodesIterator(),thePredicate,theSequence);
2931     break;
2932   case SMDSAbs_Edge:
2933     FillSequence<const SMDS_MeshElement*>(theMesh->edgesIterator(),thePredicate,theSequence);
2934     break;
2935   case SMDSAbs_Face:
2936     FillSequence<const SMDS_MeshElement*>(theMesh->facesIterator(),thePredicate,theSequence);
2937     break;
2938   case SMDSAbs_Volume:
2939     FillSequence<const SMDS_MeshElement*>(theMesh->volumesIterator(),thePredicate,theSequence);
2940     break;
2941   case SMDSAbs_All:
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);
2945     break;
2946   }
2947 }
2948
2949 void
2950 Filter::GetElementsId( const SMDS_Mesh* theMesh,
2951                        Filter::TIdSequence& theSequence )
2952 {
2953   GetElementsId(theMesh,myPredicate,theSequence);
2954 }
2955
2956 /*
2957                               ManifoldPart
2958 */
2959
2960 typedef std::set<SMDS_MeshFace*>                    TMapOfFacePtr;
2961
2962 /*
2963    Internal class Link
2964 */
2965
2966 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
2967                           SMDS_MeshNode* theNode2 )
2968 {
2969   myNode1 = theNode1;
2970   myNode2 = theNode2;
2971 }
2972
2973 ManifoldPart::Link::~Link()
2974 {
2975   myNode1 = 0;
2976   myNode2 = 0;
2977 }
2978
2979 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
2980 {
2981   if ( myNode1 == theLink.myNode1 &&
2982        myNode2 == theLink.myNode2 )
2983     return true;
2984   else if ( myNode1 == theLink.myNode2 &&
2985             myNode2 == theLink.myNode1 )
2986     return true;
2987   else
2988     return false;
2989 }
2990
2991 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
2992 {
2993   if(myNode1 < x.myNode1) return true;
2994   if(myNode1 == x.myNode1)
2995     if(myNode2 < x.myNode2) return true;
2996   return false;
2997 }
2998
2999 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
3000                             const ManifoldPart::Link& theLink2 )
3001 {
3002   return theLink1.IsEqual( theLink2 );
3003 }
3004
3005 ManifoldPart::ManifoldPart()
3006 {
3007   myMesh = 0;
3008   myAngToler = Precision::Angular();
3009   myIsOnlyManifold = true;
3010 }
3011
3012 ManifoldPart::~ManifoldPart()
3013 {
3014   myMesh = 0;
3015 }
3016
3017 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
3018 {
3019   myMesh = theMesh;
3020   process();
3021 }
3022
3023 SMDSAbs_ElementType ManifoldPart::GetType() const
3024 { return SMDSAbs_Face; }
3025
3026 bool ManifoldPart::IsSatisfy( long theElementId )
3027 {
3028   return myMapIds.Contains( theElementId );
3029 }
3030
3031 void ManifoldPart::SetAngleTolerance( const double theAngToler )
3032 { myAngToler = theAngToler; }
3033
3034 double ManifoldPart::GetAngleTolerance() const
3035 { return myAngToler; }
3036
3037 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
3038 { myIsOnlyManifold = theIsOnly; }
3039
3040 void ManifoldPart::SetStartElem( const long  theStartId )
3041 { myStartElemId = theStartId; }
3042
3043 bool ManifoldPart::process()
3044 {
3045   myMapIds.Clear();
3046   myMapBadGeomIds.Clear();
3047
3048   myAllFacePtr.clear();
3049   myAllFacePtrIntDMap.clear();
3050   if ( !myMesh )
3051     return false;
3052
3053   // collect all faces into own map
3054   SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3055   for (; anFaceItr->more(); )
3056   {
3057     SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3058     myAllFacePtr.push_back( aFacePtr );
3059     myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3060   }
3061
3062   SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3063   if ( !aStartFace )
3064     return false;
3065
3066   // the map of non manifold links and bad geometry
3067   TMapOfLink aMapOfNonManifold;
3068   TColStd_MapOfInteger aMapOfTreated;
3069
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++ )
3075   {
3076     if ( fi == aStartIndx )
3077       isStartTreat = true;
3078     // as result next time when fi will be equal to aStartIndx
3079
3080     SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3081     if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3082       continue;
3083
3084     aMapOfTreated.Add( aFacePtr->GetID() );
3085     TColStd_MapOfInteger aResFaces;
3086     if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3087                          aMapOfNonManifold, aResFaces ) )
3088       continue;
3089     TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3090     for ( ; anItr.More(); anItr.Next() )
3091     {
3092       int aFaceId = anItr.Key();
3093       aMapOfTreated.Add( aFaceId );
3094       myMapIds.Add( aFaceId );
3095     }
3096
3097     if ( fi == ( myAllFacePtr.size() - 1 ) )
3098       fi = 0;
3099   } // end run on vector of faces
3100   return !myMapIds.IsEmpty();
3101 }
3102
3103 static void getLinks( const SMDS_MeshFace* theFace,
3104                       ManifoldPart::TVectorOfLink& theLinks )
3105 {
3106   int aNbNode = theFace->NbNodes();
3107   SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3108   int i = 1;
3109   SMDS_MeshNode* aNode = 0;
3110   for ( ; aNodeItr->more() && i <= aNbNode; )
3111   {
3112
3113     SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3114     if ( i == 1 )
3115       aNode = aN1;
3116     i++;
3117     SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3118     i++;
3119     ManifoldPart::Link aLink( aN1, aN2 );
3120     theLinks.push_back( aLink );
3121   }
3122 }
3123
3124 static gp_XYZ getNormale( const SMDS_MeshFace* theFace )
3125 {
3126   gp_XYZ n;
3127   int aNbNode = theFace->NbNodes();
3128   TColgp_Array1OfXYZ anArrOfXYZ(1,4);
3129   SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3130   int i = 1;
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() ) );
3134   }
3135
3136   gp_XYZ q1 = anArrOfXYZ.Value(2) - anArrOfXYZ.Value(1);
3137   gp_XYZ q2 = anArrOfXYZ.Value(3) - anArrOfXYZ.Value(1);
3138   n  = q1 ^ q2;
3139   if ( aNbNode > 3 ) {
3140     gp_XYZ q3 = anArrOfXYZ.Value(4) - anArrOfXYZ.Value(1);
3141     n += q2 ^ q3;
3142   }
3143   double len = n.Modulus();
3144   if ( len > 0 )
3145     n /= len;
3146
3147   return n;
3148 }
3149
3150 bool ManifoldPart::findConnected
3151                  ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3152                   SMDS_MeshFace*                           theStartFace,
3153                   ManifoldPart::TMapOfLink&                theNonManifold,
3154                   TColStd_MapOfInteger&                    theResFaces )
3155 {
3156   theResFaces.Clear();
3157   if ( !theAllFacePtrInt.size() )
3158     return false;
3159
3160   if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3161   {
3162     myMapBadGeomIds.Add( theStartFace->GetID() );
3163     return false;
3164   }
3165
3166   ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3167   ManifoldPart::TVectorOfLink aSeqOfBoundary;
3168   theResFaces.Add( theStartFace->GetID() );
3169   ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3170
3171   expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3172                  aDMapLinkFace, theNonManifold, theStartFace );
3173
3174   bool isDone = false;
3175   while ( !isDone && aMapOfBoundary.size() != 0 )
3176   {
3177     bool isToReset = false;
3178     ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3179     for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3180     {
3181       ManifoldPart::Link aLink = *pLink;
3182       if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3183         continue;
3184       // each link could be treated only once
3185       aMapToSkip.insert( aLink );
3186
3187       ManifoldPart::TVectorOfFacePtr aFaces;
3188       // find next
3189       if ( myIsOnlyManifold &&
3190            (theNonManifold.find( aLink ) != theNonManifold.end()) )
3191         continue;
3192       else
3193       {
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 )
3199         {
3200           SMDS_MeshFace* aFace = *pFace;
3201           if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3202             aFiltered.push_back( aFace );
3203         }
3204         aFaces = aFiltered;
3205         if ( aFaces.size() < 2 )  // no neihgbour faces
3206           continue;
3207         else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3208         {
3209           theNonManifold.insert( aLink );
3210           continue;
3211         }
3212       }
3213
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 )
3218       {
3219         SMDS_MeshFace* aNextFace = *pFace;
3220         if ( aPrevFace == aNextFace )
3221           continue;
3222         int anNextFaceID = aNextFace->GetID();
3223         if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3224          // should not be with non manifold restriction. probably bad topology
3225           continue;
3226         // check if face was treated and skipped
3227         if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3228              !isInPlane( aPrevFace, aNextFace ) )
3229           continue;
3230         // add new element to connected and extend the boundaries.
3231         theResFaces.Add( anNextFaceID );
3232         expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3233                         aDMapLinkFace, theNonManifold, aNextFace );
3234         isToReset = true;
3235       }
3236     }
3237     isDone = !isToReset;
3238   }
3239
3240   return !theResFaces.IsEmpty();
3241 }
3242
3243 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3244                               const SMDS_MeshFace* theFace2 )
3245 {
3246   gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3247   gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3248   if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3249   {
3250     myMapBadGeomIds.Add( theFace2->GetID() );
3251     return false;
3252   }
3253   if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
3254     return true;
3255
3256   return false;
3257 }
3258
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
3265 {
3266   ManifoldPart::TVectorOfLink aLinks;
3267   getLinks( theNextFace, aLinks );
3268   int aNbLink = (int)aLinks.size();
3269   for ( int i = 0; i < aNbLink; i++ )
3270   {
3271     ManifoldPart::Link aLink = aLinks[ i ];
3272     if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3273       continue;
3274     if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3275     {
3276       if ( myIsOnlyManifold )
3277       {
3278         // remove from boundary
3279         theMapOfBoundary.erase( aLink );
3280         ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3281         for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3282         {
3283           ManifoldPart::Link aBoundLink = *pLink;
3284           if ( aBoundLink.IsEqual( aLink ) )
3285           {
3286             theSeqOfBoundary.erase( pLink );
3287             break;
3288           }
3289         }
3290       }
3291     }
3292     else
3293     {
3294       theMapOfBoundary.insert( aLink );
3295       theSeqOfBoundary.push_back( aLink );
3296       theDMapLinkFacePtr[ aLink ] = theNextFace;
3297     }
3298   }
3299 }
3300
3301 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3302                                    ManifoldPart::TVectorOfFacePtr& theFaces ) const
3303 {
3304   SMDS_Mesh::SetOfFaces aSetOfFaces;
3305   // take all faces that shared first node
3306   SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3307   for ( ; anItr->more(); )
3308   {
3309     SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3310     if ( !aFace )
3311       continue;
3312     aSetOfFaces.Add( aFace );
3313   }
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(); )
3318   {
3319     SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3320     if ( aSetOfFaces.Contains( aFace ) )
3321       theFaces.push_back( aFace );
3322   }
3323 }
3324
3325
3326 /*
3327    ElementsOnSurface
3328 */
3329
3330 ElementsOnSurface::ElementsOnSurface()
3331 {
3332   myMesh = 0;
3333   myIds.Clear();
3334   myType = SMDSAbs_All;
3335   mySurf.Nullify();
3336   myToler = Precision::Confusion();
3337   myUseBoundaries = false;
3338 }
3339
3340 ElementsOnSurface::~ElementsOnSurface()
3341 {
3342   myMesh = 0;
3343 }
3344
3345 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
3346 {
3347   if ( myMesh == theMesh )
3348     return;
3349   myMesh = theMesh;
3350   process();
3351 }
3352
3353 bool ElementsOnSurface::IsSatisfy( long theElementId )
3354 {
3355   return myIds.Contains( theElementId );
3356 }
3357
3358 SMDSAbs_ElementType ElementsOnSurface::GetType() const
3359 { return myType; }
3360
3361 void ElementsOnSurface::SetTolerance( const double theToler )
3362 {
3363   if ( myToler != theToler )
3364     myIds.Clear();
3365   myToler = theToler;
3366 }
3367
3368 double ElementsOnSurface::GetTolerance() const
3369 { return myToler; }
3370
3371 void ElementsOnSurface::SetUseBoundaries( bool theUse )
3372 {
3373   if ( myUseBoundaries != theUse ) {
3374     myUseBoundaries = theUse;
3375     SetSurface( mySurf, myType );
3376   }
3377 }
3378
3379 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
3380                                     const SMDSAbs_ElementType theType )
3381 {
3382   myIds.Clear();
3383   myType = theType;
3384   mySurf.Nullify();
3385   if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
3386     return;
3387   mySurf = TopoDS::Face( theShape );
3388   BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
3389   Standard_Real
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 );
3396   process();
3397 }
3398
3399 void ElementsOnSurface::process()
3400 {
3401   myIds.Clear();
3402   if ( mySurf.IsNull() )
3403     return;
3404
3405   if ( myMesh == 0 )
3406     return;
3407
3408   if ( myType == SMDSAbs_Face || myType == SMDSAbs_All )
3409   {
3410     myIds.ReSize( myMesh->NbFaces() );
3411     SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
3412     for(; anIter->more(); )
3413       process( anIter->next() );
3414   }
3415
3416   if ( myType == SMDSAbs_Edge || myType == SMDSAbs_All )
3417   {
3418     myIds.ReSize( myIds.Extent() + myMesh->NbEdges() );
3419     SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
3420     for(; anIter->more(); )
3421       process( anIter->next() );
3422   }
3423
3424   if ( myType == SMDSAbs_Node )
3425   {
3426     myIds.ReSize( myMesh->NbNodes() );
3427     SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
3428     for(; anIter->more(); )
3429       process( anIter->next() );
3430   }
3431 }
3432
3433 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
3434 {
3435   SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3436   bool isSatisfy = true;
3437   for ( ; aNodeItr->more(); )
3438   {
3439     SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3440     if ( !isOnSurface( aNode ) )
3441     {
3442       isSatisfy = false;
3443       break;
3444     }
3445   }
3446   if ( isSatisfy )
3447     myIds.Add( theElemPtr->GetID() );
3448 }
3449
3450 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
3451 {
3452   if ( mySurf.IsNull() )
3453     return false;
3454
3455   gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
3456   //  double aToler2 = myToler * myToler;
3457 //   if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
3458 //   {
3459 //     gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
3460 //     if ( aPln.SquareDistance( aPnt ) > aToler2 )
3461 //       return false;
3462 //   }
3463 //   else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
3464 //   {
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 )
3472 //       return false;
3473 //   }
3474 //   else
3475 //     return false;
3476   myProjector.Perform( aPnt );
3477   bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
3478
3479   return isOn;
3480 }
3481
3482
3483 /*
3484   ElementsOnShape
3485 */
3486
3487 ElementsOnShape::ElementsOnShape()
3488   : myMesh(0),
3489     myType(SMDSAbs_All),
3490     myToler(Precision::Confusion()),
3491     myAllNodesFlag(false)
3492 {
3493   myCurShapeType = TopAbs_SHAPE;
3494 }
3495
3496 ElementsOnShape::~ElementsOnShape()
3497 {
3498 }
3499
3500 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
3501 {
3502   if (myMesh != theMesh) {
3503     myMesh = theMesh;
3504     SetShape(myShape, myType);
3505   }
3506 }
3507
3508 bool ElementsOnShape::IsSatisfy (long theElementId)
3509 {
3510   return myIds.Contains(theElementId);
3511 }
3512
3513 SMDSAbs_ElementType ElementsOnShape::GetType() const
3514 {
3515   return myType;
3516 }
3517
3518 void ElementsOnShape::SetTolerance (const double theToler)
3519 {
3520   if (myToler != theToler) {
3521     myToler = theToler;
3522     SetShape(myShape, myType);
3523   }
3524 }
3525
3526 double ElementsOnShape::GetTolerance() const
3527 {
3528   return myToler;
3529 }
3530
3531 void ElementsOnShape::SetAllNodes (bool theAllNodes)
3532 {
3533   if (myAllNodesFlag != theAllNodes) {
3534     myAllNodesFlag = theAllNodes;
3535     SetShape(myShape, myType);
3536   }
3537 }
3538
3539 void ElementsOnShape::SetShape (const TopoDS_Shape&       theShape,
3540                                 const SMDSAbs_ElementType theType)
3541 {
3542   myType = theType;
3543   myShape = theShape;
3544   myIds.Clear();
3545
3546   if (myMesh == 0) return;
3547
3548   switch (myType)
3549   {
3550   case SMDSAbs_All:
3551     myIds.ReSize(myMesh->NbEdges() + myMesh->NbFaces() + myMesh->NbVolumes());
3552     break;
3553   case SMDSAbs_Node:
3554     myIds.ReSize(myMesh->NbNodes());
3555     break;
3556   case SMDSAbs_Edge:
3557     myIds.ReSize(myMesh->NbEdges());
3558     break;
3559   case SMDSAbs_Face:
3560     myIds.ReSize(myMesh->NbFaces());
3561     break;
3562   case SMDSAbs_Volume:
3563     myIds.ReSize(myMesh->NbVolumes());
3564     break;
3565   default:
3566     break;
3567   }
3568
3569   myShapesMap.Clear();
3570   addShape(myShape);
3571 }
3572
3573 void ElementsOnShape::addShape (const TopoDS_Shape& theShape)
3574 {
3575   if (theShape.IsNull() || myMesh == 0)
3576     return;
3577
3578   if (!myShapesMap.Add(theShape)) return;
3579
3580   myCurShapeType = theShape.ShapeType();
3581   switch (myCurShapeType)
3582   {
3583   case TopAbs_COMPOUND:
3584   case TopAbs_COMPSOLID:
3585   case TopAbs_SHELL:
3586   case TopAbs_WIRE:
3587     {
3588       TopoDS_Iterator anIt (theShape, Standard_True, Standard_True);
3589       for (; anIt.More(); anIt.Next()) addShape(anIt.Value());
3590     }
3591     break;
3592   case TopAbs_SOLID:
3593     {
3594       myCurSC.Load(theShape);
3595       process();
3596     }
3597     break;
3598   case TopAbs_FACE:
3599     {
3600       TopoDS_Face aFace = TopoDS::Face(theShape);
3601       BRepAdaptor_Surface SA (aFace, true);
3602       Standard_Real
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);
3609       myCurFace = aFace;
3610       process();
3611     }
3612     break;
3613   case TopAbs_EDGE:
3614     {
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);
3619       process();
3620     }
3621     break;
3622   case TopAbs_VERTEX:
3623     {
3624       TopoDS_Vertex aV = TopoDS::Vertex(theShape);
3625       myCurPnt = BRep_Tool::Pnt(aV);
3626       process();
3627     }
3628     break;
3629   default:
3630     break;
3631   }
3632 }
3633
3634 void ElementsOnShape::process()
3635 {
3636   if (myShape.IsNull() || myMesh == 0)
3637     return;
3638
3639   if (myType == SMDSAbs_Node)
3640   {
3641     SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
3642     while (anIter->more())
3643       process(anIter->next());
3644   }
3645   else
3646   {
3647     if (myType == SMDSAbs_Edge || myType == SMDSAbs_All)
3648     {
3649       SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
3650       while (anIter->more())
3651         process(anIter->next());
3652     }
3653
3654     if (myType == SMDSAbs_Face || myType == SMDSAbs_All)
3655     {
3656       SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
3657       while (anIter->more()) {
3658         process(anIter->next());
3659       }
3660     }
3661
3662     if (myType == SMDSAbs_Volume || myType == SMDSAbs_All)
3663     {
3664       SMDS_VolumeIteratorPtr anIter = myMesh->volumesIterator();
3665       while (anIter->more())
3666         process(anIter->next());
3667     }
3668   }
3669 }
3670
3671 void ElementsOnShape::process (const SMDS_MeshElement* theElemPtr)
3672 {
3673   if (myShape.IsNull())
3674     return;
3675
3676   SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3677   bool isSatisfy = myAllNodesFlag;
3678
3679   gp_XYZ centerXYZ (0, 0, 0);
3680
3681   while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
3682   {
3683     SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3684     gp_Pnt aPnt (aNode->X(), aNode->Y(), aNode->Z());
3685     centerXYZ += aPnt.XYZ();
3686
3687     switch (myCurShapeType)
3688     {
3689     case TopAbs_SOLID:
3690       {
3691         myCurSC.Perform(aPnt, myToler);
3692         isSatisfy = (myCurSC.State() == TopAbs_IN || myCurSC.State() == TopAbs_ON);
3693       }
3694       break;
3695     case TopAbs_FACE:
3696       {
3697         myCurProjFace.Perform(aPnt);
3698         isSatisfy = (myCurProjFace.IsDone() && myCurProjFace.LowerDistance() <= myToler);
3699         if (isSatisfy)
3700         {
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);
3707         }
3708       }
3709       break;
3710     case TopAbs_EDGE:
3711       {
3712         myCurProjEdge.Perform(aPnt);
3713         isSatisfy = (myCurProjEdge.NbPoints() > 0 && myCurProjEdge.LowerDistance() <= myToler);
3714       }
3715       break;
3716     case TopAbs_VERTEX:
3717       {
3718         isSatisfy = (aPnt.Distance(myCurPnt) <= myToler);
3719       }
3720       break;
3721     default:
3722       {
3723         isSatisfy = false;
3724       }
3725     }
3726   }
3727
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))
3733       isSatisfy = false;
3734   }
3735
3736   if (isSatisfy)
3737     myIds.Add(theElemPtr->GetID());
3738 }
3739
3740 TSequenceOfXYZ::TSequenceOfXYZ()
3741 {}
3742
3743 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n)
3744 {}
3745
3746 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t)
3747 {}
3748
3749 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray)
3750 {}
3751
3752 template <class InputIterator>
3753 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd)
3754 {}
3755
3756 TSequenceOfXYZ::~TSequenceOfXYZ()
3757 {}
3758
3759 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
3760 {
3761   myArray = theSequenceOfXYZ.myArray;
3762   return *this;
3763 }
3764
3765 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
3766 {
3767   return myArray[n-1];
3768 }
3769
3770 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
3771 {
3772   return myArray[n-1];
3773 }
3774
3775 void TSequenceOfXYZ::clear()
3776 {
3777   myArray.clear();
3778 }
3779
3780 void TSequenceOfXYZ::reserve(size_type n)
3781 {
3782   myArray.reserve(n);
3783 }
3784
3785 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
3786 {
3787   myArray.push_back(v);
3788 }
3789
3790 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
3791 {
3792   return myArray.size();
3793 }