Salome HOME
677f467de9fa18915f173c5dd0f974e8b438cf52
[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 #include <limits>
27
28 #include <BRepAdaptor_Surface.hxx>
29 #include <BRepClass_FaceClassifier.hxx>
30 #include <BRep_Tool.hxx>
31
32 #include <TopAbs.hxx>
33 #include <TopoDS.hxx>
34 #include <TopoDS_Edge.hxx>
35 #include <TopoDS_Face.hxx>
36 #include <TopoDS_Shape.hxx>
37 #include <TopoDS_Vertex.hxx>
38 #include <TopoDS_Iterator.hxx>
39
40 #include <Geom_CylindricalSurface.hxx>
41 #include <Geom_Plane.hxx>
42 #include <Geom_Surface.hxx>
43
44 #include <Precision.hxx>
45 #include <TColStd_MapIteratorOfMapOfInteger.hxx>
46 #include <TColStd_MapOfInteger.hxx>
47 #include <TColStd_SequenceOfAsciiString.hxx>
48 #include <TColgp_Array1OfXYZ.hxx>
49
50 #include <gp_Ax3.hxx>
51 #include <gp_Cylinder.hxx>
52 #include <gp_Dir.hxx>
53 #include <gp_Pln.hxx>
54 #include <gp_Pnt.hxx>
55 #include <gp_Vec.hxx>
56 #include <gp_XYZ.hxx>
57
58 #include "SMDS_Mesh.hxx"
59 #include "SMDS_Iterator.hxx"
60 #include "SMDS_MeshElement.hxx"
61 #include "SMDS_MeshNode.hxx"
62 #include "SMDS_VolumeTool.hxx"
63 #include "SMDS_QuadraticFaceOfNodes.hxx"
64 #include "SMDS_QuadraticEdge.hxx"
65
66 #include "SMESHDS_Mesh.hxx"
67 #include "SMESHDS_GroupBase.hxx"
68
69 /*
70                             AUXILIARY METHODS
71 */
72
73 namespace{
74
75   inline gp_XYZ gpXYZ(const SMDS_MeshNode* aNode )
76   {
77     return gp_XYZ(aNode->X(), aNode->Y(), aNode->Z() );
78   }
79
80   inline double getAngle( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
81   {
82     gp_Vec v1( P1 - P2 ), v2( P3 - P2 );
83
84     return v1.Magnitude() < gp::Resolution() ||
85       v2.Magnitude() < gp::Resolution() ? 0 : v1.Angle( v2 );
86   }
87
88   inline double getArea( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
89   {
90     gp_Vec aVec1( P2 - P1 );
91     gp_Vec aVec2( P3 - P1 );
92     return ( aVec1 ^ aVec2 ).Magnitude() * 0.5;
93   }
94
95   inline double getArea( const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3 )
96   {
97     return getArea( P1.XYZ(), P2.XYZ(), P3.XYZ() );
98   }
99
100
101
102   inline double getDistance( const gp_XYZ& P1, const gp_XYZ& P2 )
103   {
104     double aDist = gp_Pnt( P1 ).Distance( gp_Pnt( P2 ) );
105     return aDist;
106   }
107
108   int getNbMultiConnection( const SMDS_Mesh* theMesh, const int theId )
109   {
110     if ( theMesh == 0 )
111       return 0;
112
113     const SMDS_MeshElement* anEdge = theMesh->FindElement( theId );
114     if ( anEdge == 0 || anEdge->GetType() != SMDSAbs_Edge/* || anEdge->NbNodes() != 2 */)
115       return 0;
116
117     // for each pair of nodes in anEdge (there are 2 pairs in a quadratic edge)
118     // count elements containing both nodes of the pair.
119     // Note that there may be such cases for a quadratic edge (a horizontal line):
120     //
121     //  Case 1          Case 2
122     //  |     |      |        |      |
123     //  |     |      |        |      |
124     //  +-----+------+  +-----+------+ 
125     //  |            |  |            |
126     //  |            |  |            |
127     // result sould be 2 in both cases
128     //
129     int aResult0 = 0, aResult1 = 0;
130      // last node, it is a medium one in a quadratic edge
131     const SMDS_MeshNode* aLastNode = anEdge->GetNode( anEdge->NbNodes() - 1 );
132     const SMDS_MeshNode* aNode0 = anEdge->GetNode( 0 );
133     const SMDS_MeshNode* aNode1 = anEdge->GetNode( 1 );
134     if ( aNode1 == aLastNode ) aNode1 = 0;
135
136     SMDS_ElemIteratorPtr anElemIter = aLastNode->GetInverseElementIterator();
137     while( anElemIter->more() ) {
138       const SMDS_MeshElement* anElem = anElemIter->next();
139       if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
140         SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
141         while ( anIter->more() ) {
142           if ( const SMDS_MeshElement* anElemNode = anIter->next() ) {
143             if ( anElemNode == aNode0 ) {
144               aResult0++;
145               if ( !aNode1 ) break; // not a quadratic edge
146             }
147             else if ( anElemNode == aNode1 )
148               aResult1++;
149           }
150         }
151       }
152     }
153     int aResult = std::max ( aResult0, aResult1 );
154
155 //     TColStd_MapOfInteger aMap;
156
157 //     SMDS_ElemIteratorPtr anIter = anEdge->nodesIterator();
158 //     if ( anIter != 0 ) {
159 //       while( anIter->more() ) {
160 //      const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
161 //      if ( aNode == 0 )
162 //        return 0;
163 //      SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
164 //      while( anElemIter->more() ) {
165 //        const SMDS_MeshElement* anElem = anElemIter->next();
166 //        if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
167 //          int anId = anElem->GetID();
168
169 //          if ( anIter->more() )              // i.e. first node
170 //            aMap.Add( anId );
171 //          else if ( aMap.Contains( anId ) )
172 //            aResult++;
173 //        }
174 //      }
175 //       }
176 //     }
177
178     return aResult;
179   }
180
181   gp_XYZ getNormale( const SMDS_MeshFace* theFace, bool* ok=0 )
182   {
183     int aNbNode = theFace->NbNodes();
184
185     gp_XYZ q1 = gpXYZ( theFace->GetNode(1)) - gpXYZ( theFace->GetNode(0));
186     gp_XYZ q2 = gpXYZ( theFace->GetNode(2)) - gpXYZ( theFace->GetNode(0));
187     gp_XYZ n  = q1 ^ q2;
188     if ( aNbNode > 3 ) {
189       gp_XYZ q3 = gpXYZ( theFace->GetNode(3)) - gpXYZ( theFace->GetNode(0));
190       n += q2 ^ q3;
191     }
192     double len = n.Modulus();
193     bool zeroLen = ( len <= numeric_limits<double>::min());
194     if ( !zeroLen )
195       n /= len;
196
197     if (ok) *ok = !zeroLen;
198
199     return n;
200   }
201 }
202
203
204
205 using namespace SMESH::Controls;
206
207 /*
208  *                               FUNCTORS
209  */
210
211 /*
212   Class       : NumericalFunctor
213   Description : Base class for numerical functors
214 */
215 NumericalFunctor::NumericalFunctor():
216   myMesh(NULL)
217 {
218   myPrecision = -1;
219 }
220
221 void NumericalFunctor::SetMesh( const SMDS_Mesh* theMesh )
222 {
223   myMesh = theMesh;
224 }
225
226 bool NumericalFunctor::GetPoints(const int theId,
227                                  TSequenceOfXYZ& theRes ) const
228 {
229   theRes.clear();
230
231   if ( myMesh == 0 )
232     return false;
233
234   return GetPoints( myMesh->FindElement( theId ), theRes );
235 }
236
237 bool NumericalFunctor::GetPoints(const SMDS_MeshElement* anElem,
238                                  TSequenceOfXYZ&         theRes )
239 {
240   theRes.clear();
241
242   if ( anElem == 0)
243     return false;
244
245   theRes.reserve( anElem->NbNodes() );
246
247   // Get nodes of the element
248   SMDS_ElemIteratorPtr anIter;
249
250   if ( anElem->IsQuadratic() ) {
251     switch ( anElem->GetType() ) {
252     case SMDSAbs_Edge:
253       anIter = dynamic_cast<const SMDS_VtkEdge*>
254         (anElem)->interlacedNodesElemIterator();
255       break;
256     case SMDSAbs_Face:
257       anIter = dynamic_cast<const SMDS_VtkFace*>
258         (anElem)->interlacedNodesElemIterator();
259       break;
260     default:
261       anIter = anElem->nodesIterator();
262       //return false;
263     }
264   }
265   else {
266     anIter = anElem->nodesIterator();
267   }
268
269   if ( anIter ) {
270     while( anIter->more() ) {
271       if ( const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>( anIter->next() ))
272         theRes.push_back( gp_XYZ( aNode->X(), aNode->Y(), aNode->Z() ) );
273     }
274   }
275
276   return true;
277 }
278
279 long  NumericalFunctor::GetPrecision() const
280 {
281   return myPrecision;
282 }
283
284 void  NumericalFunctor::SetPrecision( const long thePrecision )
285 {
286   myPrecision = thePrecision;
287 }
288
289 double NumericalFunctor::GetValue( long theId )
290 {
291   myCurrElement = myMesh->FindElement( theId );
292   TSequenceOfXYZ P;
293   if ( GetPoints( theId, P ))
294   {
295     double aVal = GetValue( P );
296     if ( myPrecision >= 0 )
297     {
298       double prec = pow( 10., (double)( myPrecision ) );
299       aVal = floor( aVal * prec + 0.5 ) / prec;
300     }
301     return aVal;
302   }
303
304   return 0.;
305 }
306
307 //================================================================================
308 /*!
309  * \brief Return histogram of functor values
310  *  \param nbIntervals - number of intervals
311  *  \param nbEvents - number of mesh elements having values within i-th interval
312  *  \param funValues - boundaries of intervals
313  */
314 //================================================================================
315
316 void NumericalFunctor::GetHistogram(int                  nbIntervals,
317                                     std::vector<int>&    nbEvents,
318                                     std::vector<double>& funValues)
319 {
320   if ( nbIntervals < 1 ||
321        !myMesh ||
322        !myMesh->GetMeshInfo().NbElements( GetType() ))
323     return;
324   nbEvents.resize( nbIntervals, 0 );
325   funValues.resize( nbIntervals+1 );
326
327   // get all values sorted
328   std::multiset< double > values;
329   SMDS_ElemIteratorPtr elemIt = myMesh->elementsIterator(GetType());
330   while ( elemIt->more() )
331     values.insert( GetValue( elemIt->next()->GetID() ));
332
333   // case nbIntervals == 1
334   funValues[0] = *values.begin();
335   funValues[nbIntervals] = *values.rbegin();
336   if ( nbIntervals == 1 )
337   {
338     nbEvents[0] = values.size();
339     return;
340   }
341   // case of 1 value
342   if (funValues.front() == funValues.back())
343   {
344     nbEvents.resize( 1 );
345     nbEvents[0] = values.size();
346     funValues[1] = funValues.back();
347     funValues.resize( 2 );
348   }
349   // generic case
350   std::multiset< double >::iterator min = values.begin(), max;
351   for ( int i = 0; i < nbIntervals; ++i )
352   {
353     double r = (i+1) / double( nbIntervals );
354     funValues[i+1] = funValues.front() * (1-r) + funValues.back() * r;
355     if ( min != values.end() && *min <= funValues[i+1] )
356     {
357       max = values.upper_bound( funValues[i+1] ); // greater than funValues[i+1], or end()
358       nbEvents[i] = std::distance( min, max );
359       min = max;
360     }
361   }
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_VtkFace* F =
1626         dynamic_cast<const SMDS_VtkFace*>(anElem);
1627       // use special nodes iterator
1628       SMDS_ElemIteratorPtr anIter = F->interlacedNodesElemIterator();
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 = dynamic_cast<const SMDS_VtkFace*>
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
1907 /*
1908   Class       : FreeBorders
1909   Description : Predicate for free borders
1910 */
1911
1912 FreeBorders::FreeBorders()
1913 {
1914   myMesh = 0;
1915 }
1916
1917 void FreeBorders::SetMesh( const SMDS_Mesh* theMesh )
1918 {
1919   myMesh = theMesh;
1920 }
1921
1922 bool FreeBorders::IsSatisfy( long theId )
1923 {
1924   return getNbMultiConnection( myMesh, theId ) == 1;
1925 }
1926
1927 SMDSAbs_ElementType FreeBorders::GetType() const
1928 {
1929   return SMDSAbs_Edge;
1930 }
1931
1932
1933 /*
1934   Class       : FreeEdges
1935   Description : Predicate for free Edges
1936 */
1937 FreeEdges::FreeEdges()
1938 {
1939   myMesh = 0;
1940 }
1941
1942 void FreeEdges::SetMesh( const SMDS_Mesh* theMesh )
1943 {
1944   myMesh = theMesh;
1945 }
1946
1947 bool FreeEdges::IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId  )
1948 {
1949   TColStd_MapOfInteger aMap;
1950   for ( int i = 0; i < 2; i++ )
1951   {
1952     SMDS_ElemIteratorPtr anElemIter = theNodes[ i ]->GetInverseElementIterator();
1953     while( anElemIter->more() )
1954     {
1955       const SMDS_MeshElement* anElem = anElemIter->next();
1956       if ( anElem != 0 && anElem->GetType() == SMDSAbs_Face )
1957       {
1958         int anId = anElem->GetID();
1959
1960         if ( i == 0 )
1961           aMap.Add( anId );
1962         else if ( aMap.Contains( anId ) && anId != theFaceId )
1963           return false;
1964       }
1965     }
1966   }
1967   return true;
1968 }
1969
1970 bool FreeEdges::IsSatisfy( long theId )
1971 {
1972   if ( myMesh == 0 )
1973     return false;
1974
1975   const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
1976   if ( aFace == 0 || aFace->GetType() != SMDSAbs_Face || aFace->NbNodes() < 3 )
1977     return false;
1978
1979   SMDS_ElemIteratorPtr anIter;
1980   if ( aFace->IsQuadratic() ) {
1981     anIter = dynamic_cast<const SMDS_VtkFace*>
1982       (aFace)->interlacedNodesElemIterator();
1983   }
1984   else {
1985     anIter = aFace->nodesIterator();
1986   }
1987   if ( anIter == 0 )
1988     return false;
1989
1990   int i = 0, nbNodes = aFace->NbNodes();
1991   std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
1992   while( anIter->more() )
1993   {
1994     const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
1995     if ( aNode == 0 )
1996       return false;
1997     aNodes[ i++ ] = aNode;
1998   }
1999   aNodes[ nbNodes ] = aNodes[ 0 ];
2000
2001   for ( i = 0; i < nbNodes; i++ )
2002     if ( IsFreeEdge( &aNodes[ i ], theId ) )
2003       return true;
2004
2005   return false;
2006 }
2007
2008 SMDSAbs_ElementType FreeEdges::GetType() const
2009 {
2010   return SMDSAbs_Face;
2011 }
2012
2013 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2014   myElemId(theElemId)
2015 {
2016   myPntId[0] = thePntId1;  myPntId[1] = thePntId2;
2017   if(thePntId1 > thePntId2){
2018     myPntId[1] = thePntId1;  myPntId[0] = thePntId2;
2019   }
2020 }
2021
2022 bool FreeEdges::Border::operator<(const FreeEdges::Border& x) const{
2023   if(myPntId[0] < x.myPntId[0]) return true;
2024   if(myPntId[0] == x.myPntId[0])
2025     if(myPntId[1] < x.myPntId[1]) return true;
2026   return false;
2027 }
2028
2029 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2030                           FreeEdges::TBorders& theRegistry,
2031                           FreeEdges::TBorders& theContainer)
2032 {
2033   if(theRegistry.find(theBorder) == theRegistry.end()){
2034     theRegistry.insert(theBorder);
2035     theContainer.insert(theBorder);
2036   }else{
2037     theContainer.erase(theBorder);
2038   }
2039 }
2040
2041 void FreeEdges::GetBoreders(TBorders& theBorders)
2042 {
2043   TBorders aRegistry;
2044   SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2045   for(; anIter->more(); ){
2046     const SMDS_MeshFace* anElem = anIter->next();
2047     long anElemId = anElem->GetID();
2048     SMDS_ElemIteratorPtr aNodesIter;
2049     if ( anElem->IsQuadratic() )
2050       aNodesIter = static_cast<const SMDS_VtkFace*>(anElem)->
2051         interlacedNodesElemIterator();
2052     else
2053       aNodesIter = anElem->nodesIterator();
2054     long aNodeId[2];
2055     const SMDS_MeshElement* aNode;
2056     if(aNodesIter->more()){
2057       aNode = aNodesIter->next();
2058       aNodeId[0] = aNodeId[1] = aNode->GetID();
2059     }
2060     for(; aNodesIter->more(); ){
2061       aNode = aNodesIter->next();
2062       long anId = aNode->GetID();
2063       Border aBorder(anElemId,aNodeId[1],anId);
2064       aNodeId[1] = anId;
2065       //std::cout<<aBorder.myPntId[0]<<"; "<<aBorder.myPntId[1]<<"; "<<aBorder.myElemId<<endl;
2066       UpdateBorders(aBorder,aRegistry,theBorders);
2067     }
2068     Border aBorder(anElemId,aNodeId[0],aNodeId[1]);
2069     //std::cout<<aBorder.myPntId[0]<<"; "<<aBorder.myPntId[1]<<"; "<<aBorder.myElemId<<endl;
2070     UpdateBorders(aBorder,aRegistry,theBorders);
2071   }
2072   //std::cout<<"theBorders.size() = "<<theBorders.size()<<endl;
2073 }
2074
2075
2076 /*
2077   Class       : FreeNodes
2078   Description : Predicate for free nodes
2079 */
2080
2081 FreeNodes::FreeNodes()
2082 {
2083   myMesh = 0;
2084 }
2085
2086 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2087 {
2088   myMesh = theMesh;
2089 }
2090
2091 bool FreeNodes::IsSatisfy( long theNodeId )
2092 {
2093   const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2094   if (!aNode)
2095     return false;
2096
2097   return (aNode->NbInverseElements() < 1);
2098 }
2099
2100 SMDSAbs_ElementType FreeNodes::GetType() const
2101 {
2102   return SMDSAbs_Node;
2103 }
2104
2105
2106 /*
2107   Class       : FreeFaces
2108   Description : Predicate for free faces
2109 */
2110
2111 FreeFaces::FreeFaces()
2112 {
2113   myMesh = 0;
2114 }
2115
2116 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2117 {
2118   myMesh = theMesh;
2119 }
2120
2121 bool FreeFaces::IsSatisfy( long theId )
2122 {
2123   if (!myMesh) return false;
2124   // check that faces nodes refers to less than two common volumes
2125   const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2126   if ( !aFace || aFace->GetType() != SMDSAbs_Face )
2127     return false;
2128
2129   int nbNode = aFace->NbNodes();
2130
2131   // collect volumes check that number of volumss with count equal nbNode not less than 2
2132   typedef map< SMDS_MeshElement*, int > TMapOfVolume; // map of volume counters
2133   typedef map< SMDS_MeshElement*, int >::iterator TItrMapOfVolume; // iterator
2134   TMapOfVolume mapOfVol;
2135
2136   SMDS_ElemIteratorPtr nodeItr = aFace->nodesIterator();
2137   while ( nodeItr->more() ) {
2138     const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>(nodeItr->next());
2139     if ( !aNode ) continue;
2140     SMDS_ElemIteratorPtr volItr = aNode->GetInverseElementIterator(SMDSAbs_Volume);
2141     while ( volItr->more() ) {
2142       SMDS_MeshElement* aVol = (SMDS_MeshElement*)volItr->next();
2143       TItrMapOfVolume itr = mapOfVol.insert(make_pair(aVol, 0)).first;
2144       (*itr).second++;
2145     } 
2146   }
2147   int nbVol = 0;
2148   TItrMapOfVolume volItr = mapOfVol.begin();
2149   TItrMapOfVolume volEnd = mapOfVol.end();
2150   for ( ; volItr != volEnd; ++volItr )
2151     if ( (*volItr).second >= nbNode )
2152        nbVol++;
2153   // face is not free if number of volumes constructed on thier nodes more than one
2154   return (nbVol < 2);
2155 }
2156
2157 SMDSAbs_ElementType FreeFaces::GetType() const
2158 {
2159   return SMDSAbs_Face;
2160 }
2161
2162 /*
2163   Class       : LinearOrQuadratic
2164   Description : Predicate to verify whether a mesh element is linear
2165 */
2166
2167 LinearOrQuadratic::LinearOrQuadratic()
2168 {
2169   myMesh = 0;
2170 }
2171
2172 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2173 {
2174   myMesh = theMesh;
2175 }
2176
2177 bool LinearOrQuadratic::IsSatisfy( long theId )
2178 {
2179   if (!myMesh) return false;
2180   const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2181   if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2182     return false;
2183   return (!anElem->IsQuadratic());
2184 }
2185
2186 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2187 {
2188   myType = theType;
2189 }
2190
2191 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2192 {
2193   return myType;
2194 }
2195
2196 /*
2197   Class       : GroupColor
2198   Description : Functor for check color of group to whic mesh element belongs to
2199 */
2200
2201 GroupColor::GroupColor()
2202 {
2203 }
2204
2205 bool GroupColor::IsSatisfy( long theId )
2206 {
2207   return (myIDs.find( theId ) != myIDs.end());
2208 }
2209
2210 void GroupColor::SetType( SMDSAbs_ElementType theType )
2211 {
2212   myType = theType;
2213 }
2214
2215 SMDSAbs_ElementType GroupColor::GetType() const
2216 {
2217   return myType;
2218 }
2219
2220 static bool isEqual( const Quantity_Color& theColor1,
2221                      const Quantity_Color& theColor2 )
2222 {
2223   // tolerance to compare colors
2224   const double tol = 5*1e-3;
2225   return ( fabs( theColor1.Red() - theColor2.Red() ) < tol &&
2226            fabs( theColor1.Green() - theColor2.Green() ) < tol &&
2227            fabs( theColor1.Blue() - theColor2.Blue() ) < tol );
2228 }
2229
2230
2231 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2232 {
2233   myIDs.clear();
2234   
2235   const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2236   if ( !aMesh )
2237     return;
2238
2239   int nbGrp = aMesh->GetNbGroups();
2240   if ( !nbGrp )
2241     return;
2242   
2243   // iterates on groups and find necessary elements ids
2244   const std::set<SMESHDS_GroupBase*>& aGroups = aMesh->GetGroups();
2245   set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
2246   for (; GrIt != aGroups.end(); GrIt++) {
2247     SMESHDS_GroupBase* aGrp = (*GrIt);
2248     if ( !aGrp )
2249       continue;
2250     // check type and color of group
2251     if ( !isEqual( myColor, aGrp->GetColor() ) )
2252       continue;
2253     if ( myType != SMDSAbs_All && myType != (SMDSAbs_ElementType)aGrp->GetType() )
2254       continue;
2255
2256     SMDSAbs_ElementType aGrpElType = (SMDSAbs_ElementType)aGrp->GetType();
2257     if ( myType == aGrpElType || (myType == SMDSAbs_All && aGrpElType != SMDSAbs_Node) ) {
2258       // add elements IDS into control
2259       int aSize = aGrp->Extent();
2260       for (int i = 0; i < aSize; i++)
2261         myIDs.insert( aGrp->GetID(i+1) );
2262     }
2263   }
2264 }
2265
2266 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
2267 {
2268   TCollection_AsciiString aStr = theStr;
2269   aStr.RemoveAll( ' ' );
2270   aStr.RemoveAll( '\t' );
2271   for ( int aPos = aStr.Search( ";;" ); aPos != -1; aPos = aStr.Search( ";;" ) )
2272     aStr.Remove( aPos, 2 );
2273   Standard_Real clr[3];
2274   clr[0] = clr[1] = clr[2] = 0.;
2275   for ( int i = 0; i < 3; i++ ) {
2276     TCollection_AsciiString tmpStr = aStr.Token( ";", i+1 );
2277     if ( !tmpStr.IsEmpty() && tmpStr.IsRealValue() )
2278       clr[i] = tmpStr.RealValue();
2279   }
2280   myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
2281 }
2282
2283 //=======================================================================
2284 // name    : GetRangeStr
2285 // Purpose : Get range as a string.
2286 //           Example: "1,2,3,50-60,63,67,70-"
2287 //=======================================================================
2288 void GroupColor::GetColorStr( TCollection_AsciiString& theResStr ) const
2289 {
2290   theResStr.Clear();
2291   theResStr += TCollection_AsciiString( myColor.Red() );
2292   theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2293   theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2294 }
2295
2296 /*
2297   Class       : ElemGeomType
2298   Description : Predicate to check element geometry type
2299 */
2300
2301 ElemGeomType::ElemGeomType()
2302 {
2303   myMesh = 0;
2304   myType = SMDSAbs_All;
2305   myGeomType = SMDSGeom_TRIANGLE;
2306 }
2307
2308 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2309 {
2310   myMesh = theMesh;
2311 }
2312
2313 bool ElemGeomType::IsSatisfy( long theId )
2314 {
2315   if (!myMesh) return false;
2316   const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2317   if ( !anElem )
2318     return false;
2319   const SMDSAbs_ElementType anElemType = anElem->GetType();
2320   if ( myType != SMDSAbs_All && anElemType != myType )
2321     return false;
2322   const int aNbNode = anElem->NbNodes();
2323   bool isOk = false;
2324   switch( anElemType )
2325   {
2326   case SMDSAbs_Node:
2327     isOk = (myGeomType == SMDSGeom_POINT);
2328     break;
2329
2330   case SMDSAbs_Edge:
2331     isOk = (myGeomType == SMDSGeom_EDGE);
2332     break;
2333
2334   case SMDSAbs_Face:
2335     if ( myGeomType == SMDSGeom_TRIANGLE )
2336       isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 6 : aNbNode == 3));
2337     else if ( myGeomType == SMDSGeom_QUADRANGLE )
2338       isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 8 : aNbNode == 4));
2339     else if ( myGeomType == SMDSGeom_POLYGON )
2340       isOk = anElem->IsPoly();
2341     break;
2342
2343   case SMDSAbs_Volume:
2344     if ( myGeomType == SMDSGeom_TETRA )
2345       isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 10 : aNbNode == 4));
2346     else if ( myGeomType == SMDSGeom_PYRAMID )
2347       isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 13 : aNbNode == 5));
2348     else if ( myGeomType == SMDSGeom_PENTA )
2349       isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 15 : aNbNode == 6));
2350     else if ( myGeomType == SMDSGeom_HEXA )
2351       isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 20 : aNbNode == 8));
2352      else if ( myGeomType == SMDSGeom_POLYHEDRA )
2353       isOk = anElem->IsPoly();
2354     break;
2355     default: break;
2356   }
2357   return isOk;
2358 }
2359
2360 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2361 {
2362   myType = theType;
2363 }
2364
2365 SMDSAbs_ElementType ElemGeomType::GetType() const
2366 {
2367   return myType;
2368 }
2369
2370 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2371 {
2372   myGeomType = theType;
2373 }
2374
2375 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2376 {
2377   return myGeomType;
2378 }
2379
2380 //================================================================================
2381 /*!
2382  * \brief Class CoplanarFaces
2383  */
2384 //================================================================================
2385
2386 CoplanarFaces::CoplanarFaces()
2387   : myMesh(0), myFaceID(0), myToler(0)
2388 {
2389 }
2390 bool CoplanarFaces::IsSatisfy( long theElementId )
2391 {
2392   if ( myCoplanarIDs.empty() )
2393   {
2394     // Build a set of coplanar face ids
2395
2396     if ( !myMesh || !myFaceID || !myToler )
2397       return false;
2398
2399     const SMDS_MeshElement* face = myMesh->FindElement( myFaceID );
2400     if ( !face || face->GetType() != SMDSAbs_Face )
2401       return false;
2402
2403     bool normOK;
2404     gp_Vec myNorm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2405     if (!normOK)
2406       return false;
2407
2408     const double radianTol = myToler * PI180;
2409     typedef SMDS_StdIterator< const SMDS_MeshElement*, SMDS_ElemIteratorPtr > TFaceIt;
2410     std::set<const SMDS_MeshElement*> checkedFaces, checkedNodes;
2411     std::list<const SMDS_MeshElement*> faceQueue( 1, face );
2412     while ( !faceQueue.empty() )
2413     {
2414       face = faceQueue.front();
2415       if ( checkedFaces.insert( face ).second )
2416       {
2417         gp_Vec norm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2418         if (!normOK || myNorm.Angle( norm ) <= radianTol)
2419         {
2420           myCoplanarIDs.insert( face->GetID() );
2421           std::set<const SMDS_MeshElement*> neighborFaces;
2422           for ( int i = 0; i < face->NbCornerNodes(); ++i )
2423           {
2424             const SMDS_MeshNode* n = face->GetNode( i );
2425             if ( checkedNodes.insert( n ).second )
2426               neighborFaces.insert( TFaceIt( n->GetInverseElementIterator(SMDSAbs_Face)),
2427                                     TFaceIt());
2428           }
2429           faceQueue.insert( faceQueue.end(), neighborFaces.begin(), neighborFaces.end() );
2430         }
2431       }
2432       faceQueue.pop_front();
2433     }
2434   }
2435   return myCoplanarIDs.count( theElementId );
2436 }
2437
2438 /*
2439   *Class       : RangeOfIds
2440   *Description : Predicate for Range of Ids.
2441   *              Range may be specified with two ways.
2442   *              1. Using AddToRange method
2443   *              2. With SetRangeStr method. Parameter of this method is a string
2444   *                 like as "1,2,3,50-60,63,67,70-"
2445 */
2446
2447 //=======================================================================
2448 // name    : RangeOfIds
2449 // Purpose : Constructor
2450 //=======================================================================
2451 RangeOfIds::RangeOfIds()
2452 {
2453   myMesh = 0;
2454   myType = SMDSAbs_All;
2455 }
2456
2457 //=======================================================================
2458 // name    : SetMesh
2459 // Purpose : Set mesh
2460 //=======================================================================
2461 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
2462 {
2463   myMesh = theMesh;
2464 }
2465
2466 //=======================================================================
2467 // name    : AddToRange
2468 // Purpose : Add ID to the range
2469 //=======================================================================
2470 bool RangeOfIds::AddToRange( long theEntityId )
2471 {
2472   myIds.Add( theEntityId );
2473   return true;
2474 }
2475
2476 //=======================================================================
2477 // name    : GetRangeStr
2478 // Purpose : Get range as a string.
2479 //           Example: "1,2,3,50-60,63,67,70-"
2480 //=======================================================================
2481 void RangeOfIds::GetRangeStr( TCollection_AsciiString& theResStr )
2482 {
2483   theResStr.Clear();
2484
2485   TColStd_SequenceOfInteger     anIntSeq;
2486   TColStd_SequenceOfAsciiString aStrSeq;
2487
2488   TColStd_MapIteratorOfMapOfInteger anIter( myIds );
2489   for ( ; anIter.More(); anIter.Next() )
2490   {
2491     int anId = anIter.Key();
2492     TCollection_AsciiString aStr( anId );
2493     anIntSeq.Append( anId );
2494     aStrSeq.Append( aStr );
2495   }
2496
2497   for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2498   {
2499     int aMinId = myMin( i );
2500     int aMaxId = myMax( i );
2501
2502     TCollection_AsciiString aStr;
2503     if ( aMinId != IntegerFirst() )
2504       aStr += aMinId;
2505
2506     aStr += "-";
2507
2508     if ( aMaxId != IntegerLast() )
2509       aStr += aMaxId;
2510
2511     // find position of the string in result sequence and insert string in it
2512     if ( anIntSeq.Length() == 0 )
2513     {
2514       anIntSeq.Append( aMinId );
2515       aStrSeq.Append( aStr );
2516     }
2517     else
2518     {
2519       if ( aMinId < anIntSeq.First() )
2520       {
2521         anIntSeq.Prepend( aMinId );
2522         aStrSeq.Prepend( aStr );
2523       }
2524       else if ( aMinId > anIntSeq.Last() )
2525       {
2526         anIntSeq.Append( aMinId );
2527         aStrSeq.Append( aStr );
2528       }
2529       else
2530         for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
2531           if ( aMinId < anIntSeq( j ) )
2532           {
2533             anIntSeq.InsertBefore( j, aMinId );
2534             aStrSeq.InsertBefore( j, aStr );
2535             break;
2536           }
2537     }
2538   }
2539
2540   if ( aStrSeq.Length() == 0 )
2541     return;
2542
2543   theResStr = aStrSeq( 1 );
2544   for ( int j = 2, k = aStrSeq.Length(); j <= k; j++  )
2545   {
2546     theResStr += ",";
2547     theResStr += aStrSeq( j );
2548   }
2549 }
2550
2551 //=======================================================================
2552 // name    : SetRangeStr
2553 // Purpose : Define range with string
2554 //           Example of entry string: "1,2,3,50-60,63,67,70-"
2555 //=======================================================================
2556 bool RangeOfIds::SetRangeStr( const TCollection_AsciiString& theStr )
2557 {
2558   myMin.Clear();
2559   myMax.Clear();
2560   myIds.Clear();
2561
2562   TCollection_AsciiString aStr = theStr;
2563   aStr.RemoveAll( ' ' );
2564   aStr.RemoveAll( '\t' );
2565
2566   for ( int aPos = aStr.Search( ",," ); aPos != -1; aPos = aStr.Search( ",," ) )
2567     aStr.Remove( aPos, 2 );
2568
2569   TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
2570   int i = 1;
2571   while ( tmpStr != "" )
2572   {
2573     tmpStr = aStr.Token( ",", i++ );
2574     int aPos = tmpStr.Search( '-' );
2575
2576     if ( aPos == -1 )
2577     {
2578       if ( tmpStr.IsIntegerValue() )
2579         myIds.Add( tmpStr.IntegerValue() );
2580       else
2581         return false;
2582     }
2583     else
2584     {
2585       TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
2586       TCollection_AsciiString aMinStr = tmpStr;
2587
2588       while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
2589       while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
2590
2591       if ( (!aMinStr.IsEmpty() && !aMinStr.IsIntegerValue()) ||
2592            (!aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue()) )
2593         return false;
2594
2595       myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
2596       myMax.Append( aMaxStr.IsEmpty() ? IntegerLast()  : aMaxStr.IntegerValue() );
2597     }
2598   }
2599
2600   return true;
2601 }
2602
2603 //=======================================================================
2604 // name    : GetType
2605 // Purpose : Get type of supported entities
2606 //=======================================================================
2607 SMDSAbs_ElementType RangeOfIds::GetType() const
2608 {
2609   return myType;
2610 }
2611
2612 //=======================================================================
2613 // name    : SetType
2614 // Purpose : Set type of supported entities
2615 //=======================================================================
2616 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
2617 {
2618   myType = theType;
2619 }
2620
2621 //=======================================================================
2622 // name    : IsSatisfy
2623 // Purpose : Verify whether entity satisfies to this rpedicate
2624 //=======================================================================
2625 bool RangeOfIds::IsSatisfy( long theId )
2626 {
2627   if ( !myMesh )
2628     return false;
2629
2630   if ( myType == SMDSAbs_Node )
2631   {
2632     if ( myMesh->FindNode( theId ) == 0 )
2633       return false;
2634   }
2635   else
2636   {
2637     const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2638     if ( anElem == 0 || (myType != anElem->GetType() && myType != SMDSAbs_All ))
2639       return false;
2640   }
2641
2642   if ( myIds.Contains( theId ) )
2643     return true;
2644
2645   for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2646     if ( theId >= myMin( i ) && theId <= myMax( i ) )
2647       return true;
2648
2649   return false;
2650 }
2651
2652 /*
2653   Class       : Comparator
2654   Description : Base class for comparators
2655 */
2656 Comparator::Comparator():
2657   myMargin(0)
2658 {}
2659
2660 Comparator::~Comparator()
2661 {}
2662
2663 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
2664 {
2665   if ( myFunctor )
2666     myFunctor->SetMesh( theMesh );
2667 }
2668
2669 void Comparator::SetMargin( double theValue )
2670 {
2671   myMargin = theValue;
2672 }
2673
2674 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
2675 {
2676   myFunctor = theFunct;
2677 }
2678
2679 SMDSAbs_ElementType Comparator::GetType() const
2680 {
2681   return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
2682 }
2683
2684 double Comparator::GetMargin()
2685 {
2686   return myMargin;
2687 }
2688
2689
2690 /*
2691   Class       : LessThan
2692   Description : Comparator "<"
2693 */
2694 bool LessThan::IsSatisfy( long theId )
2695 {
2696   return myFunctor && myFunctor->GetValue( theId ) < myMargin;
2697 }
2698
2699
2700 /*
2701   Class       : MoreThan
2702   Description : Comparator ">"
2703 */
2704 bool MoreThan::IsSatisfy( long theId )
2705 {
2706   return myFunctor && myFunctor->GetValue( theId ) > myMargin;
2707 }
2708
2709
2710 /*
2711   Class       : EqualTo
2712   Description : Comparator "="
2713 */
2714 EqualTo::EqualTo():
2715   myToler(Precision::Confusion())
2716 {}
2717
2718 bool EqualTo::IsSatisfy( long theId )
2719 {
2720   return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
2721 }
2722
2723 void EqualTo::SetTolerance( double theToler )
2724 {
2725   myToler = theToler;
2726 }
2727
2728 double EqualTo::GetTolerance()
2729 {
2730   return myToler;
2731 }
2732
2733 /*
2734   Class       : LogicalNOT
2735   Description : Logical NOT predicate
2736 */
2737 LogicalNOT::LogicalNOT()
2738 {}
2739
2740 LogicalNOT::~LogicalNOT()
2741 {}
2742
2743 bool LogicalNOT::IsSatisfy( long theId )
2744 {
2745   return myPredicate && !myPredicate->IsSatisfy( theId );
2746 }
2747
2748 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
2749 {
2750   if ( myPredicate )
2751     myPredicate->SetMesh( theMesh );
2752 }
2753
2754 void LogicalNOT::SetPredicate( PredicatePtr thePred )
2755 {
2756   myPredicate = thePred;
2757 }
2758
2759 SMDSAbs_ElementType LogicalNOT::GetType() const
2760 {
2761   return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
2762 }
2763
2764
2765 /*
2766   Class       : LogicalBinary
2767   Description : Base class for binary logical predicate
2768 */
2769 LogicalBinary::LogicalBinary()
2770 {}
2771
2772 LogicalBinary::~LogicalBinary()
2773 {}
2774
2775 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
2776 {
2777   if ( myPredicate1 )
2778     myPredicate1->SetMesh( theMesh );
2779
2780   if ( myPredicate2 )
2781     myPredicate2->SetMesh( theMesh );
2782 }
2783
2784 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
2785 {
2786   myPredicate1 = thePredicate;
2787 }
2788
2789 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
2790 {
2791   myPredicate2 = thePredicate;
2792 }
2793
2794 SMDSAbs_ElementType LogicalBinary::GetType() const
2795 {
2796   if ( !myPredicate1 || !myPredicate2 )
2797     return SMDSAbs_All;
2798
2799   SMDSAbs_ElementType aType1 = myPredicate1->GetType();
2800   SMDSAbs_ElementType aType2 = myPredicate2->GetType();
2801
2802   return aType1 == aType2 ? aType1 : SMDSAbs_All;
2803 }
2804
2805
2806 /*
2807   Class       : LogicalAND
2808   Description : Logical AND
2809 */
2810 bool LogicalAND::IsSatisfy( long theId )
2811 {
2812   return
2813     myPredicate1 &&
2814     myPredicate2 &&
2815     myPredicate1->IsSatisfy( theId ) &&
2816     myPredicate2->IsSatisfy( theId );
2817 }
2818
2819
2820 /*
2821   Class       : LogicalOR
2822   Description : Logical OR
2823 */
2824 bool LogicalOR::IsSatisfy( long theId )
2825 {
2826   return
2827     myPredicate1 &&
2828     myPredicate2 &&
2829     (myPredicate1->IsSatisfy( theId ) ||
2830     myPredicate2->IsSatisfy( theId ));
2831 }
2832
2833
2834 /*
2835                               FILTER
2836 */
2837
2838 Filter::Filter()
2839 {}
2840
2841 Filter::~Filter()
2842 {}
2843
2844 void Filter::SetPredicate( PredicatePtr thePredicate )
2845 {
2846   myPredicate = thePredicate;
2847 }
2848
2849 template<class TElement, class TIterator, class TPredicate>
2850 inline void FillSequence(const TIterator& theIterator,
2851                          TPredicate& thePredicate,
2852                          Filter::TIdSequence& theSequence)
2853 {
2854   if ( theIterator ) {
2855     while( theIterator->more() ) {
2856       TElement anElem = theIterator->next();
2857       long anId = anElem->GetID();
2858       if ( thePredicate->IsSatisfy( anId ) )
2859         theSequence.push_back( anId );
2860     }
2861   }
2862 }
2863
2864 void
2865 Filter::
2866 GetElementsId( const SMDS_Mesh* theMesh,
2867                PredicatePtr thePredicate,
2868                TIdSequence& theSequence )
2869 {
2870   theSequence.clear();
2871
2872   if ( !theMesh || !thePredicate )
2873     return;
2874
2875   thePredicate->SetMesh( theMesh );
2876
2877   SMDSAbs_ElementType aType = thePredicate->GetType();
2878   switch(aType){
2879   case SMDSAbs_Node:
2880     FillSequence<const SMDS_MeshNode*>(theMesh->nodesIterator(),thePredicate,theSequence);
2881     break;
2882   case SMDSAbs_Edge:
2883     FillSequence<const SMDS_MeshElement*>(theMesh->edgesIterator(),thePredicate,theSequence);
2884     break;
2885   case SMDSAbs_Face:
2886     FillSequence<const SMDS_MeshElement*>(theMesh->facesIterator(),thePredicate,theSequence);
2887     break;
2888   case SMDSAbs_Volume:
2889     FillSequence<const SMDS_MeshElement*>(theMesh->volumesIterator(),thePredicate,theSequence);
2890     break;
2891   case SMDSAbs_All:
2892     FillSequence<const SMDS_MeshElement*>(theMesh->edgesIterator(),thePredicate,theSequence);
2893     FillSequence<const SMDS_MeshElement*>(theMesh->facesIterator(),thePredicate,theSequence);
2894     FillSequence<const SMDS_MeshElement*>(theMesh->volumesIterator(),thePredicate,theSequence);
2895     break;
2896   }
2897 }
2898
2899 void
2900 Filter::GetElementsId( const SMDS_Mesh* theMesh,
2901                        Filter::TIdSequence& theSequence )
2902 {
2903   GetElementsId(theMesh,myPredicate,theSequence);
2904 }
2905
2906 /*
2907                               ManifoldPart
2908 */
2909
2910 typedef std::set<SMDS_MeshFace*>                    TMapOfFacePtr;
2911
2912 /*
2913    Internal class Link
2914 */
2915
2916 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
2917                           SMDS_MeshNode* theNode2 )
2918 {
2919   myNode1 = theNode1;
2920   myNode2 = theNode2;
2921 }
2922
2923 ManifoldPart::Link::~Link()
2924 {
2925   myNode1 = 0;
2926   myNode2 = 0;
2927 }
2928
2929 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
2930 {
2931   if ( myNode1 == theLink.myNode1 &&
2932        myNode2 == theLink.myNode2 )
2933     return true;
2934   else if ( myNode1 == theLink.myNode2 &&
2935             myNode2 == theLink.myNode1 )
2936     return true;
2937   else
2938     return false;
2939 }
2940
2941 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
2942 {
2943   if(myNode1 < x.myNode1) return true;
2944   if(myNode1 == x.myNode1)
2945     if(myNode2 < x.myNode2) return true;
2946   return false;
2947 }
2948
2949 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
2950                             const ManifoldPart::Link& theLink2 )
2951 {
2952   return theLink1.IsEqual( theLink2 );
2953 }
2954
2955 ManifoldPart::ManifoldPart()
2956 {
2957   myMesh = 0;
2958   myAngToler = Precision::Angular();
2959   myIsOnlyManifold = true;
2960 }
2961
2962 ManifoldPart::~ManifoldPart()
2963 {
2964   myMesh = 0;
2965 }
2966
2967 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
2968 {
2969   myMesh = theMesh;
2970   process();
2971 }
2972
2973 SMDSAbs_ElementType ManifoldPart::GetType() const
2974 { return SMDSAbs_Face; }
2975
2976 bool ManifoldPart::IsSatisfy( long theElementId )
2977 {
2978   return myMapIds.Contains( theElementId );
2979 }
2980
2981 void ManifoldPart::SetAngleTolerance( const double theAngToler )
2982 { myAngToler = theAngToler; }
2983
2984 double ManifoldPart::GetAngleTolerance() const
2985 { return myAngToler; }
2986
2987 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
2988 { myIsOnlyManifold = theIsOnly; }
2989
2990 void ManifoldPart::SetStartElem( const long  theStartId )
2991 { myStartElemId = theStartId; }
2992
2993 bool ManifoldPart::process()
2994 {
2995   myMapIds.Clear();
2996   myMapBadGeomIds.Clear();
2997
2998   myAllFacePtr.clear();
2999   myAllFacePtrIntDMap.clear();
3000   if ( !myMesh )
3001     return false;
3002
3003   // collect all faces into own map
3004   SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3005   for (; anFaceItr->more(); )
3006   {
3007     SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3008     myAllFacePtr.push_back( aFacePtr );
3009     myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3010   }
3011
3012   SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3013   if ( !aStartFace )
3014     return false;
3015
3016   // the map of non manifold links and bad geometry
3017   TMapOfLink aMapOfNonManifold;
3018   TColStd_MapOfInteger aMapOfTreated;
3019
3020   // begin cycle on faces from start index and run on vector till the end
3021   //  and from begin to start index to cover whole vector
3022   const int aStartIndx = myAllFacePtrIntDMap[aStartFace];
3023   bool isStartTreat = false;
3024   for ( int fi = aStartIndx; !isStartTreat || fi != aStartIndx ; fi++ )
3025   {
3026     if ( fi == aStartIndx )
3027       isStartTreat = true;
3028     // as result next time when fi will be equal to aStartIndx
3029
3030     SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3031     if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3032       continue;
3033
3034     aMapOfTreated.Add( aFacePtr->GetID() );
3035     TColStd_MapOfInteger aResFaces;
3036     if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3037                          aMapOfNonManifold, aResFaces ) )
3038       continue;
3039     TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3040     for ( ; anItr.More(); anItr.Next() )
3041     {
3042       int aFaceId = anItr.Key();
3043       aMapOfTreated.Add( aFaceId );
3044       myMapIds.Add( aFaceId );
3045     }
3046
3047     if ( fi == ( myAllFacePtr.size() - 1 ) )
3048       fi = 0;
3049   } // end run on vector of faces
3050   return !myMapIds.IsEmpty();
3051 }
3052
3053 static void getLinks( const SMDS_MeshFace* theFace,
3054                       ManifoldPart::TVectorOfLink& theLinks )
3055 {
3056   int aNbNode = theFace->NbNodes();
3057   SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3058   int i = 1;
3059   SMDS_MeshNode* aNode = 0;
3060   for ( ; aNodeItr->more() && i <= aNbNode; )
3061   {
3062
3063     SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3064     if ( i == 1 )
3065       aNode = aN1;
3066     i++;
3067     SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3068     i++;
3069     ManifoldPart::Link aLink( aN1, aN2 );
3070     theLinks.push_back( aLink );
3071   }
3072 }
3073
3074 bool ManifoldPart::findConnected
3075                  ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3076                   SMDS_MeshFace*                           theStartFace,
3077                   ManifoldPart::TMapOfLink&                theNonManifold,
3078                   TColStd_MapOfInteger&                    theResFaces )
3079 {
3080   theResFaces.Clear();
3081   if ( !theAllFacePtrInt.size() )
3082     return false;
3083
3084   if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3085   {
3086     myMapBadGeomIds.Add( theStartFace->GetID() );
3087     return false;
3088   }
3089
3090   ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3091   ManifoldPart::TVectorOfLink aSeqOfBoundary;
3092   theResFaces.Add( theStartFace->GetID() );
3093   ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3094
3095   expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3096                  aDMapLinkFace, theNonManifold, theStartFace );
3097
3098   bool isDone = false;
3099   while ( !isDone && aMapOfBoundary.size() != 0 )
3100   {
3101     bool isToReset = false;
3102     ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3103     for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3104     {
3105       ManifoldPart::Link aLink = *pLink;
3106       if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3107         continue;
3108       // each link could be treated only once
3109       aMapToSkip.insert( aLink );
3110
3111       ManifoldPart::TVectorOfFacePtr aFaces;
3112       // find next
3113       if ( myIsOnlyManifold &&
3114            (theNonManifold.find( aLink ) != theNonManifold.end()) )
3115         continue;
3116       else
3117       {
3118         getFacesByLink( aLink, aFaces );
3119         // filter the element to keep only indicated elements
3120         ManifoldPart::TVectorOfFacePtr aFiltered;
3121         ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3122         for ( ; pFace != aFaces.end(); ++pFace )
3123         {
3124           SMDS_MeshFace* aFace = *pFace;
3125           if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3126             aFiltered.push_back( aFace );
3127         }
3128         aFaces = aFiltered;
3129         if ( aFaces.size() < 2 )  // no neihgbour faces
3130           continue;
3131         else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3132         {
3133           theNonManifold.insert( aLink );
3134           continue;
3135         }
3136       }
3137
3138       // compare normal with normals of neighbor element
3139       SMDS_MeshFace* aPrevFace = aDMapLinkFace[ aLink ];
3140       ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3141       for ( ; pFace != aFaces.end(); ++pFace )
3142       {
3143         SMDS_MeshFace* aNextFace = *pFace;
3144         if ( aPrevFace == aNextFace )
3145           continue;
3146         int anNextFaceID = aNextFace->GetID();
3147         if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3148          // should not be with non manifold restriction. probably bad topology
3149           continue;
3150         // check if face was treated and skipped
3151         if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3152              !isInPlane( aPrevFace, aNextFace ) )
3153           continue;
3154         // add new element to connected and extend the boundaries.
3155         theResFaces.Add( anNextFaceID );
3156         expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3157                         aDMapLinkFace, theNonManifold, aNextFace );
3158         isToReset = true;
3159       }
3160     }
3161     isDone = !isToReset;
3162   }
3163
3164   return !theResFaces.IsEmpty();
3165 }
3166
3167 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3168                               const SMDS_MeshFace* theFace2 )
3169 {
3170   gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3171   gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3172   if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3173   {
3174     myMapBadGeomIds.Add( theFace2->GetID() );
3175     return false;
3176   }
3177   if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
3178     return true;
3179
3180   return false;
3181 }
3182
3183 void ManifoldPart::expandBoundary
3184                    ( ManifoldPart::TMapOfLink&            theMapOfBoundary,
3185                      ManifoldPart::TVectorOfLink&         theSeqOfBoundary,
3186                      ManifoldPart::TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
3187                      ManifoldPart::TMapOfLink&            theNonManifold,
3188                      SMDS_MeshFace*                       theNextFace ) const
3189 {
3190   ManifoldPart::TVectorOfLink aLinks;
3191   getLinks( theNextFace, aLinks );
3192   int aNbLink = (int)aLinks.size();
3193   for ( int i = 0; i < aNbLink; i++ )
3194   {
3195     ManifoldPart::Link aLink = aLinks[ i ];
3196     if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3197       continue;
3198     if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3199     {
3200       if ( myIsOnlyManifold )
3201       {
3202         // remove from boundary
3203         theMapOfBoundary.erase( aLink );
3204         ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3205         for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3206         {
3207           ManifoldPart::Link aBoundLink = *pLink;
3208           if ( aBoundLink.IsEqual( aLink ) )
3209           {
3210             theSeqOfBoundary.erase( pLink );
3211             break;
3212           }
3213         }
3214       }
3215     }
3216     else
3217     {
3218       theMapOfBoundary.insert( aLink );
3219       theSeqOfBoundary.push_back( aLink );
3220       theDMapLinkFacePtr[ aLink ] = theNextFace;
3221     }
3222   }
3223 }
3224
3225 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3226                                    ManifoldPart::TVectorOfFacePtr& theFaces ) const
3227 {
3228   std::set<SMDS_MeshCell *> aSetOfFaces;
3229   // take all faces that shared first node
3230   SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3231   for ( ; anItr->more(); )
3232   {
3233     SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3234     if ( !aFace )
3235       continue;
3236     aSetOfFaces.insert( aFace );
3237   }
3238   // take all faces that shared second node
3239   anItr = theLink.myNode2->facesIterator();
3240   // find the common part of two sets
3241   for ( ; anItr->more(); )
3242   {
3243     SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3244     if ( aSetOfFaces.count( aFace ) )
3245       theFaces.push_back( aFace );
3246   }
3247 }
3248
3249
3250 /*
3251    ElementsOnSurface
3252 */
3253
3254 ElementsOnSurface::ElementsOnSurface()
3255 {
3256   myMesh = 0;
3257   myIds.Clear();
3258   myType = SMDSAbs_All;
3259   mySurf.Nullify();
3260   myToler = Precision::Confusion();
3261   myUseBoundaries = false;
3262 }
3263
3264 ElementsOnSurface::~ElementsOnSurface()
3265 {
3266   myMesh = 0;
3267 }
3268
3269 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
3270 {
3271   if ( myMesh == theMesh )
3272     return;
3273   myMesh = theMesh;
3274   process();
3275 }
3276
3277 bool ElementsOnSurface::IsSatisfy( long theElementId )
3278 {
3279   return myIds.Contains( theElementId );
3280 }
3281
3282 SMDSAbs_ElementType ElementsOnSurface::GetType() const
3283 { return myType; }
3284
3285 void ElementsOnSurface::SetTolerance( const double theToler )
3286 {
3287   if ( myToler != theToler )
3288     myIds.Clear();
3289   myToler = theToler;
3290 }
3291
3292 double ElementsOnSurface::GetTolerance() const
3293 { return myToler; }
3294
3295 void ElementsOnSurface::SetUseBoundaries( bool theUse )
3296 {
3297   if ( myUseBoundaries != theUse ) {
3298     myUseBoundaries = theUse;
3299     SetSurface( mySurf, myType );
3300   }
3301 }
3302
3303 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
3304                                     const SMDSAbs_ElementType theType )
3305 {
3306   myIds.Clear();
3307   myType = theType;
3308   mySurf.Nullify();
3309   if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
3310     return;
3311   mySurf = TopoDS::Face( theShape );
3312   BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
3313   Standard_Real
3314     u1 = SA.FirstUParameter(),
3315     u2 = SA.LastUParameter(),
3316     v1 = SA.FirstVParameter(),
3317     v2 = SA.LastVParameter();
3318   Handle(Geom_Surface) surf = BRep_Tool::Surface( mySurf );
3319   myProjector.Init( surf, u1,u2, v1,v2 );
3320   process();
3321 }
3322
3323 void ElementsOnSurface::process()
3324 {
3325   myIds.Clear();
3326   if ( mySurf.IsNull() )
3327     return;
3328
3329   if ( myMesh == 0 )
3330     return;
3331
3332   if ( myType == SMDSAbs_Face || myType == SMDSAbs_All )
3333   {
3334     myIds.ReSize( myMesh->NbFaces() );
3335     SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
3336     for(; anIter->more(); )
3337       process( anIter->next() );
3338   }
3339
3340   if ( myType == SMDSAbs_Edge || myType == SMDSAbs_All )
3341   {
3342     myIds.ReSize( myIds.Extent() + myMesh->NbEdges() );
3343     SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
3344     for(; anIter->more(); )
3345       process( anIter->next() );
3346   }
3347
3348   if ( myType == SMDSAbs_Node )
3349   {
3350     myIds.ReSize( myMesh->NbNodes() );
3351     SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
3352     for(; anIter->more(); )
3353       process( anIter->next() );
3354   }
3355 }
3356
3357 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
3358 {
3359   SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3360   bool isSatisfy = true;
3361   for ( ; aNodeItr->more(); )
3362   {
3363     SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3364     if ( !isOnSurface( aNode ) )
3365     {
3366       isSatisfy = false;
3367       break;
3368     }
3369   }
3370   if ( isSatisfy )
3371     myIds.Add( theElemPtr->GetID() );
3372 }
3373
3374 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
3375 {
3376   if ( mySurf.IsNull() )
3377     return false;
3378
3379   gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
3380   //  double aToler2 = myToler * myToler;
3381 //   if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
3382 //   {
3383 //     gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
3384 //     if ( aPln.SquareDistance( aPnt ) > aToler2 )
3385 //       return false;
3386 //   }
3387 //   else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
3388 //   {
3389 //     gp_Cylinder aCyl = Handle(Geom_CylindricalSurface)::DownCast(mySurf)->Cylinder();
3390 //     double aRad = aCyl.Radius();
3391 //     gp_Ax3 anAxis = aCyl.Position();
3392 //     gp_XYZ aLoc = aCyl.Location().XYZ();
3393 //     double aXDist = anAxis.XDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3394 //     double aYDist = anAxis.YDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3395 //     if ( fabs(aXDist*aXDist + aYDist*aYDist - aRad*aRad) > aToler2 )
3396 //       return false;
3397 //   }
3398 //   else
3399 //     return false;
3400   myProjector.Perform( aPnt );
3401   bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
3402
3403   return isOn;
3404 }
3405
3406
3407 /*
3408   ElementsOnShape
3409 */
3410
3411 ElementsOnShape::ElementsOnShape()
3412   : myMesh(0),
3413     myType(SMDSAbs_All),
3414     myToler(Precision::Confusion()),
3415     myAllNodesFlag(false)
3416 {
3417   myCurShapeType = TopAbs_SHAPE;
3418 }
3419
3420 ElementsOnShape::~ElementsOnShape()
3421 {
3422 }
3423
3424 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
3425 {
3426   if (myMesh != theMesh) {
3427     myMesh = theMesh;
3428     SetShape(myShape, myType);
3429   }
3430 }
3431
3432 bool ElementsOnShape::IsSatisfy (long theElementId)
3433 {
3434   return myIds.Contains(theElementId);
3435 }
3436
3437 SMDSAbs_ElementType ElementsOnShape::GetType() const
3438 {
3439   return myType;
3440 }
3441
3442 void ElementsOnShape::SetTolerance (const double theToler)
3443 {
3444   if (myToler != theToler) {
3445     myToler = theToler;
3446     SetShape(myShape, myType);
3447   }
3448 }
3449
3450 double ElementsOnShape::GetTolerance() const
3451 {
3452   return myToler;
3453 }
3454
3455 void ElementsOnShape::SetAllNodes (bool theAllNodes)
3456 {
3457   if (myAllNodesFlag != theAllNodes) {
3458     myAllNodesFlag = theAllNodes;
3459     SetShape(myShape, myType);
3460   }
3461 }
3462
3463 void ElementsOnShape::SetShape (const TopoDS_Shape&       theShape,
3464                                 const SMDSAbs_ElementType theType)
3465 {
3466   myType = theType;
3467   myShape = theShape;
3468   myIds.Clear();
3469
3470   if (myMesh == 0) return;
3471
3472   switch (myType)
3473   {
3474   case SMDSAbs_All:
3475     myIds.ReSize(myMesh->NbEdges() + myMesh->NbFaces() + myMesh->NbVolumes());
3476     break;
3477   case SMDSAbs_Node:
3478     myIds.ReSize(myMesh->NbNodes());
3479     break;
3480   case SMDSAbs_Edge:
3481     myIds.ReSize(myMesh->NbEdges());
3482     break;
3483   case SMDSAbs_Face:
3484     myIds.ReSize(myMesh->NbFaces());
3485     break;
3486   case SMDSAbs_Volume:
3487     myIds.ReSize(myMesh->NbVolumes());
3488     break;
3489   default:
3490     break;
3491   }
3492
3493   myShapesMap.Clear();
3494   addShape(myShape);
3495 }
3496
3497 void ElementsOnShape::addShape (const TopoDS_Shape& theShape)
3498 {
3499   if (theShape.IsNull() || myMesh == 0)
3500     return;
3501
3502   if (!myShapesMap.Add(theShape)) return;
3503
3504   myCurShapeType = theShape.ShapeType();
3505   switch (myCurShapeType)
3506   {
3507   case TopAbs_COMPOUND:
3508   case TopAbs_COMPSOLID:
3509   case TopAbs_SHELL:
3510   case TopAbs_WIRE:
3511     {
3512       TopoDS_Iterator anIt (theShape, Standard_True, Standard_True);
3513       for (; anIt.More(); anIt.Next()) addShape(anIt.Value());
3514     }
3515     break;
3516   case TopAbs_SOLID:
3517     {
3518       myCurSC.Load(theShape);
3519       process();
3520     }
3521     break;
3522   case TopAbs_FACE:
3523     {
3524       TopoDS_Face aFace = TopoDS::Face(theShape);
3525       BRepAdaptor_Surface SA (aFace, true);
3526       Standard_Real
3527         u1 = SA.FirstUParameter(),
3528         u2 = SA.LastUParameter(),
3529         v1 = SA.FirstVParameter(),
3530         v2 = SA.LastVParameter();
3531       Handle(Geom_Surface) surf = BRep_Tool::Surface(aFace);
3532       myCurProjFace.Init(surf, u1,u2, v1,v2);
3533       myCurFace = aFace;
3534       process();
3535     }
3536     break;
3537   case TopAbs_EDGE:
3538     {
3539       TopoDS_Edge anEdge = TopoDS::Edge(theShape);
3540       Standard_Real u1, u2;
3541       Handle(Geom_Curve) curve = BRep_Tool::Curve(anEdge, u1, u2);
3542       myCurProjEdge.Init(curve, u1, u2);
3543       process();
3544     }
3545     break;
3546   case TopAbs_VERTEX:
3547     {
3548       TopoDS_Vertex aV = TopoDS::Vertex(theShape);
3549       myCurPnt = BRep_Tool::Pnt(aV);
3550       process();
3551     }
3552     break;
3553   default:
3554     break;
3555   }
3556 }
3557
3558 void ElementsOnShape::process()
3559 {
3560   if (myShape.IsNull() || myMesh == 0)
3561     return;
3562
3563   if (myType == SMDSAbs_Node)
3564   {
3565     SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
3566     while (anIter->more())
3567       process(anIter->next());
3568   }
3569   else
3570   {
3571     if (myType == SMDSAbs_Edge || myType == SMDSAbs_All)
3572     {
3573       SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
3574       while (anIter->more())
3575         process(anIter->next());
3576     }
3577
3578     if (myType == SMDSAbs_Face || myType == SMDSAbs_All)
3579     {
3580       SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
3581       while (anIter->more()) {
3582         process(anIter->next());
3583       }
3584     }
3585
3586     if (myType == SMDSAbs_Volume || myType == SMDSAbs_All)
3587     {
3588       SMDS_VolumeIteratorPtr anIter = myMesh->volumesIterator();
3589       while (anIter->more())
3590         process(anIter->next());
3591     }
3592   }
3593 }
3594
3595 void ElementsOnShape::process (const SMDS_MeshElement* theElemPtr)
3596 {
3597   if (myShape.IsNull())
3598     return;
3599
3600   SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3601   bool isSatisfy = myAllNodesFlag;
3602
3603   gp_XYZ centerXYZ (0, 0, 0);
3604
3605   while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
3606   {
3607     SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3608     gp_Pnt aPnt (aNode->X(), aNode->Y(), aNode->Z());
3609     centerXYZ += aPnt.XYZ();
3610
3611     switch (myCurShapeType)
3612     {
3613     case TopAbs_SOLID:
3614       {
3615         myCurSC.Perform(aPnt, myToler);
3616         isSatisfy = (myCurSC.State() == TopAbs_IN || myCurSC.State() == TopAbs_ON);
3617       }
3618       break;
3619     case TopAbs_FACE:
3620       {
3621         myCurProjFace.Perform(aPnt);
3622         isSatisfy = (myCurProjFace.IsDone() && myCurProjFace.LowerDistance() <= myToler);
3623         if (isSatisfy)
3624         {
3625           // check relatively the face
3626           Quantity_Parameter u, v;
3627           myCurProjFace.LowerDistanceParameters(u, v);
3628           gp_Pnt2d aProjPnt (u, v);
3629           BRepClass_FaceClassifier aClsf (myCurFace, aProjPnt, myToler);
3630           isSatisfy = (aClsf.State() == TopAbs_IN || aClsf.State() == TopAbs_ON);
3631         }
3632       }
3633       break;
3634     case TopAbs_EDGE:
3635       {
3636         myCurProjEdge.Perform(aPnt);
3637         isSatisfy = (myCurProjEdge.NbPoints() > 0 && myCurProjEdge.LowerDistance() <= myToler);
3638       }
3639       break;
3640     case TopAbs_VERTEX:
3641       {
3642         isSatisfy = (aPnt.Distance(myCurPnt) <= myToler);
3643       }
3644       break;
3645     default:
3646       {
3647         isSatisfy = false;
3648       }
3649     }
3650   }
3651
3652   if (isSatisfy && myCurShapeType == TopAbs_SOLID) { // Check the center point for volumes MantisBug 0020168
3653     centerXYZ /= theElemPtr->NbNodes();
3654     gp_Pnt aCenterPnt (centerXYZ);
3655     myCurSC.Perform(aCenterPnt, myToler);
3656     if ( !(myCurSC.State() == TopAbs_IN || myCurSC.State() == TopAbs_ON))
3657       isSatisfy = false;
3658   }
3659
3660   if (isSatisfy)
3661     myIds.Add(theElemPtr->GetID());
3662 }
3663
3664 TSequenceOfXYZ::TSequenceOfXYZ()
3665 {}
3666
3667 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n)
3668 {}
3669
3670 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t)
3671 {}
3672
3673 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray)
3674 {}
3675
3676 template <class InputIterator>
3677 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd)
3678 {}
3679
3680 TSequenceOfXYZ::~TSequenceOfXYZ()
3681 {}
3682
3683 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
3684 {
3685   myArray = theSequenceOfXYZ.myArray;
3686   return *this;
3687 }
3688
3689 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
3690 {
3691   return myArray[n-1];
3692 }
3693
3694 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
3695 {
3696   return myArray[n-1];
3697 }
3698
3699 void TSequenceOfXYZ::clear()
3700 {
3701   myArray.clear();
3702 }
3703
3704 void TSequenceOfXYZ::reserve(size_type n)
3705 {
3706   myArray.reserve(n);
3707 }
3708
3709 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
3710 {
3711   myArray.push_back(v);
3712 }
3713
3714 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
3715 {
3716   return myArray.size();
3717 }