1 // Copyright (C) 2007-2010 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
23 #include "SMESH_ControlsDef.hxx"
28 #include <BRepAdaptor_Surface.hxx>
29 #include <BRepClass_FaceClassifier.hxx>
30 #include <BRep_Tool.hxx>
34 #include <TopoDS_Edge.hxx>
35 #include <TopoDS_Face.hxx>
36 #include <TopoDS_Shape.hxx>
37 #include <TopoDS_Vertex.hxx>
38 #include <TopoDS_Iterator.hxx>
40 #include <Geom_CylindricalSurface.hxx>
41 #include <Geom_Plane.hxx>
42 #include <Geom_Surface.hxx>
44 #include <Precision.hxx>
45 #include <TColStd_MapIteratorOfMapOfInteger.hxx>
46 #include <TColStd_MapOfInteger.hxx>
47 #include <TColStd_SequenceOfAsciiString.hxx>
48 #include <TColgp_Array1OfXYZ.hxx>
51 #include <gp_Cylinder.hxx>
58 #include "SMDS_Mesh.hxx"
59 #include "SMDS_Iterator.hxx"
60 #include "SMDS_MeshElement.hxx"
61 #include "SMDS_MeshNode.hxx"
62 #include "SMDS_VolumeTool.hxx"
63 #include "SMDS_QuadraticFaceOfNodes.hxx"
64 #include "SMDS_QuadraticEdge.hxx"
66 #include "SMESHDS_Mesh.hxx"
67 #include "SMESHDS_GroupBase.hxx"
75 inline gp_XYZ gpXYZ(const SMDS_MeshNode* aNode )
77 return gp_XYZ(aNode->X(), aNode->Y(), aNode->Z() );
80 inline double getAngle( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
82 gp_Vec v1( P1 - P2 ), v2( P3 - P2 );
84 return v1.Magnitude() < gp::Resolution() ||
85 v2.Magnitude() < gp::Resolution() ? 0 : v1.Angle( v2 );
88 inline double getArea( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
90 gp_Vec aVec1( P2 - P1 );
91 gp_Vec aVec2( P3 - P1 );
92 return ( aVec1 ^ aVec2 ).Magnitude() * 0.5;
95 inline double getArea( const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3 )
97 return getArea( P1.XYZ(), P2.XYZ(), P3.XYZ() );
102 inline double getDistance( const gp_XYZ& P1, const gp_XYZ& P2 )
104 double aDist = gp_Pnt( P1 ).Distance( gp_Pnt( P2 ) );
108 int getNbMultiConnection( const SMDS_Mesh* theMesh, const int theId )
113 const SMDS_MeshElement* anEdge = theMesh->FindElement( theId );
114 if ( anEdge == 0 || anEdge->GetType() != SMDSAbs_Edge/* || anEdge->NbNodes() != 2 */)
117 // for each pair of nodes in anEdge (there are 2 pairs in a quadratic edge)
118 // count elements containing both nodes of the pair.
119 // Note that there may be such cases for a quadratic edge (a horizontal line):
124 // +-----+------+ +-----+------+
127 // result sould be 2 in both cases
129 int aResult0 = 0, aResult1 = 0;
130 // last node, it is a medium one in a quadratic edge
131 const SMDS_MeshNode* aLastNode = anEdge->GetNode( anEdge->NbNodes() - 1 );
132 const SMDS_MeshNode* aNode0 = anEdge->GetNode( 0 );
133 const SMDS_MeshNode* aNode1 = anEdge->GetNode( 1 );
134 if ( aNode1 == aLastNode ) aNode1 = 0;
136 SMDS_ElemIteratorPtr anElemIter = aLastNode->GetInverseElementIterator();
137 while( anElemIter->more() ) {
138 const SMDS_MeshElement* anElem = anElemIter->next();
139 if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
140 SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
141 while ( anIter->more() ) {
142 if ( const SMDS_MeshElement* anElemNode = anIter->next() ) {
143 if ( anElemNode == aNode0 ) {
145 if ( !aNode1 ) break; // not a quadratic edge
147 else if ( anElemNode == aNode1 )
153 int aResult = std::max ( aResult0, aResult1 );
155 // TColStd_MapOfInteger aMap;
157 // SMDS_ElemIteratorPtr anIter = anEdge->nodesIterator();
158 // if ( anIter != 0 ) {
159 // while( anIter->more() ) {
160 // const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
163 // SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
164 // while( anElemIter->more() ) {
165 // const SMDS_MeshElement* anElem = anElemIter->next();
166 // if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
167 // int anId = anElem->GetID();
169 // if ( anIter->more() ) // i.e. first node
171 // else if ( aMap.Contains( anId ) )
181 gp_XYZ getNormale( const SMDS_MeshFace* theFace, bool* ok=0 )
183 int aNbNode = theFace->NbNodes();
185 gp_XYZ q1 = gpXYZ( theFace->GetNode(1)) - gpXYZ( theFace->GetNode(0));
186 gp_XYZ q2 = gpXYZ( theFace->GetNode(2)) - gpXYZ( theFace->GetNode(0));
189 gp_XYZ q3 = gpXYZ( theFace->GetNode(3)) - gpXYZ( theFace->GetNode(0));
192 double len = n.Modulus();
193 bool zeroLen = ( len <= numeric_limits<double>::min());
197 if (ok) *ok = !zeroLen;
205 using namespace SMESH::Controls;
212 Class : NumericalFunctor
213 Description : Base class for numerical functors
215 NumericalFunctor::NumericalFunctor():
221 void NumericalFunctor::SetMesh( const SMDS_Mesh* theMesh )
226 bool NumericalFunctor::GetPoints(const int theId,
227 TSequenceOfXYZ& theRes ) const
234 return GetPoints( myMesh->FindElement( theId ), theRes );
237 bool NumericalFunctor::GetPoints(const SMDS_MeshElement* anElem,
238 TSequenceOfXYZ& theRes )
245 theRes.reserve( anElem->NbNodes() );
247 // Get nodes of the element
248 SMDS_ElemIteratorPtr anIter;
250 if ( anElem->IsQuadratic() ) {
251 switch ( anElem->GetType() ) {
253 anIter = static_cast<const SMDS_QuadraticEdge*>
254 (anElem)->interlacedNodesElemIterator();
257 anIter = static_cast<const SMDS_QuadraticFaceOfNodes*>
258 (anElem)->interlacedNodesElemIterator();
261 anIter = anElem->nodesIterator();
266 anIter = anElem->nodesIterator();
270 while( anIter->more() ) {
271 if ( const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>( anIter->next() ))
272 theRes.push_back( gp_XYZ( aNode->X(), aNode->Y(), aNode->Z() ) );
279 long NumericalFunctor::GetPrecision() const
284 void NumericalFunctor::SetPrecision( const long thePrecision )
286 myPrecision = thePrecision;
289 double NumericalFunctor::GetValue( long theId )
291 myCurrElement = myMesh->FindElement( theId );
293 if ( GetPoints( theId, P ))
295 double aVal = GetValue( P );
296 if ( myPrecision >= 0 )
298 double prec = pow( 10., (double)( myPrecision ) );
299 aVal = floor( aVal * prec + 0.5 ) / prec;
307 //================================================================================
309 * \brief Return histogram of functor values
310 * \param nbIntervals - number of intervals
311 * \param nbEvents - number of mesh elements having values within i-th interval
312 * \param funValues - boundaries of intervals
314 //================================================================================
316 void NumericalFunctor::GetHistogram(int nbIntervals,
317 std::vector<int>& nbEvents,
318 std::vector<double>& funValues)
320 if ( nbIntervals < 1 ||
322 !myMesh->GetMeshInfo().NbElements( GetType() ))
324 nbEvents.resize( nbIntervals, 0 );
325 funValues.resize( nbIntervals+1 );
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() ));
333 // case nbIntervals == 1
334 funValues[0] = *values.begin();
335 funValues[nbIntervals] = *values.rbegin();
336 if ( nbIntervals == 1 )
338 nbEvents[0] = values.size();
342 if (funValues.front() == funValues.back())
344 nbEvents.resize( 1 );
345 nbEvents[0] = values.size();
346 funValues[1] = funValues.back();
347 funValues.resize( 2 );
350 std::multiset< double >::iterator min = values.begin(), max;
351 for ( int i = 0; i < nbIntervals; ++i )
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] )
357 max = values.upper_bound( funValues[i+1] ); // greater than funValues[i+1], or end()
358 nbEvents[i] = std::distance( min, max );
364 //=======================================================================
365 //function : GetValue
367 //=======================================================================
369 double Volume::GetValue( long theElementId )
371 if ( theElementId && myMesh ) {
372 SMDS_VolumeTool aVolumeTool;
373 if ( aVolumeTool.Set( myMesh->FindElement( theElementId )))
374 return aVolumeTool.GetSize();
379 //=======================================================================
380 //function : GetBadRate
381 //purpose : meaningless as it is not quality control functor
382 //=======================================================================
384 double Volume::GetBadRate( double Value, int /*nbNodes*/ ) const
389 //=======================================================================
392 //=======================================================================
394 SMDSAbs_ElementType Volume::GetType() const
396 return SMDSAbs_Volume;
401 Class : MaxElementLength2D
402 Description : Functor calculating maximum length of 2D element
405 double MaxElementLength2D::GetValue( long theElementId )
408 if( GetPoints( theElementId, P ) ) {
410 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
411 SMDSAbs_ElementType aType = aElem->GetType();
415 if( len == 3 ) { // triangles
416 double L1 = getDistance(P( 1 ),P( 2 ));
417 double L2 = getDistance(P( 2 ),P( 3 ));
418 double L3 = getDistance(P( 3 ),P( 1 ));
419 aVal = Max(L1,Max(L2,L3));
422 else if( len == 4 ) { // quadrangles
423 double L1 = getDistance(P( 1 ),P( 2 ));
424 double L2 = getDistance(P( 2 ),P( 3 ));
425 double L3 = getDistance(P( 3 ),P( 4 ));
426 double L4 = getDistance(P( 4 ),P( 1 ));
427 double D1 = getDistance(P( 1 ),P( 3 ));
428 double D2 = getDistance(P( 2 ),P( 4 ));
429 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
432 else if( len == 6 ) { // quadratic triangles
433 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
434 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
435 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
436 aVal = Max(L1,Max(L2,L3));
439 else if( len == 8 ) { // quadratic quadrangles
440 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
441 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
442 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
443 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
444 double D1 = getDistance(P( 1 ),P( 5 ));
445 double D2 = getDistance(P( 3 ),P( 7 ));
446 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
451 if( myPrecision >= 0 )
453 double prec = pow( 10., (double)myPrecision );
454 aVal = floor( aVal * prec + 0.5 ) / prec;
461 double MaxElementLength2D::GetBadRate( double Value, int /*nbNodes*/ ) const
466 SMDSAbs_ElementType MaxElementLength2D::GetType() const
472 Class : MaxElementLength3D
473 Description : Functor calculating maximum length of 3D element
476 double MaxElementLength3D::GetValue( long theElementId )
479 if( GetPoints( theElementId, P ) ) {
481 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
482 SMDSAbs_ElementType aType = aElem->GetType();
486 if( len == 4 ) { // tetras
487 double L1 = getDistance(P( 1 ),P( 2 ));
488 double L2 = getDistance(P( 2 ),P( 3 ));
489 double L3 = getDistance(P( 3 ),P( 1 ));
490 double L4 = getDistance(P( 1 ),P( 4 ));
491 double L5 = getDistance(P( 2 ),P( 4 ));
492 double L6 = getDistance(P( 3 ),P( 4 ));
493 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
496 else if( len == 5 ) { // pyramids
497 double L1 = getDistance(P( 1 ),P( 2 ));
498 double L2 = getDistance(P( 2 ),P( 3 ));
499 double L3 = getDistance(P( 3 ),P( 4 ));
500 double L4 = getDistance(P( 4 ),P( 1 ));
501 double L5 = getDistance(P( 1 ),P( 5 ));
502 double L6 = getDistance(P( 2 ),P( 5 ));
503 double L7 = getDistance(P( 3 ),P( 5 ));
504 double L8 = getDistance(P( 4 ),P( 5 ));
505 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
506 aVal = Max(aVal,Max(L7,L8));
509 else if( len == 6 ) { // pentas
510 double L1 = getDistance(P( 1 ),P( 2 ));
511 double L2 = getDistance(P( 2 ),P( 3 ));
512 double L3 = getDistance(P( 3 ),P( 1 ));
513 double L4 = getDistance(P( 4 ),P( 5 ));
514 double L5 = getDistance(P( 5 ),P( 6 ));
515 double L6 = getDistance(P( 6 ),P( 4 ));
516 double L7 = getDistance(P( 1 ),P( 4 ));
517 double L8 = getDistance(P( 2 ),P( 5 ));
518 double L9 = getDistance(P( 3 ),P( 6 ));
519 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
520 aVal = Max(aVal,Max(Max(L7,L8),L9));
523 else if( len == 8 ) { // hexas
524 double L1 = getDistance(P( 1 ),P( 2 ));
525 double L2 = getDistance(P( 2 ),P( 3 ));
526 double L3 = getDistance(P( 3 ),P( 4 ));
527 double L4 = getDistance(P( 4 ),P( 1 ));
528 double L5 = getDistance(P( 5 ),P( 6 ));
529 double L6 = getDistance(P( 6 ),P( 7 ));
530 double L7 = getDistance(P( 7 ),P( 8 ));
531 double L8 = getDistance(P( 8 ),P( 5 ));
532 double L9 = getDistance(P( 1 ),P( 5 ));
533 double L10= getDistance(P( 2 ),P( 6 ));
534 double L11= getDistance(P( 3 ),P( 7 ));
535 double L12= getDistance(P( 4 ),P( 8 ));
536 double D1 = getDistance(P( 1 ),P( 7 ));
537 double D2 = getDistance(P( 2 ),P( 8 ));
538 double D3 = getDistance(P( 3 ),P( 5 ));
539 double D4 = getDistance(P( 4 ),P( 6 ));
540 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
541 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
542 aVal = Max(aVal,Max(L11,L12));
543 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
546 else if( len == 10 ) { // quadratic tetras
547 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
548 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
549 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
550 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
551 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
552 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
553 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
556 else if( len == 13 ) { // quadratic pyramids
557 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
558 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
559 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
560 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
561 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
562 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
563 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
564 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
565 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
566 aVal = Max(aVal,Max(L7,L8));
569 else if( len == 15 ) { // quadratic pentas
570 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
571 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
572 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
573 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
574 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
575 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
576 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
577 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
578 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
579 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
580 aVal = Max(aVal,Max(Max(L7,L8),L9));
583 else if( len == 20 ) { // quadratic hexas
584 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
585 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
586 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
587 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
588 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
589 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
590 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
591 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
592 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
593 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
594 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
595 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
596 double D1 = getDistance(P( 1 ),P( 7 ));
597 double D2 = getDistance(P( 2 ),P( 8 ));
598 double D3 = getDistance(P( 3 ),P( 5 ));
599 double D4 = getDistance(P( 4 ),P( 6 ));
600 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
601 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
602 aVal = Max(aVal,Max(L11,L12));
603 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
606 else if( len > 1 && aElem->IsPoly() ) { // polys
607 // get the maximum distance between all pairs of nodes
608 for( int i = 1; i <= len; i++ ) {
609 for( int j = 1; j <= len; j++ ) {
610 if( j > i ) { // optimization of the loop
611 double D = getDistance( P(i), P(j) );
612 aVal = Max( aVal, D );
619 if( myPrecision >= 0 )
621 double prec = pow( 10., (double)myPrecision );
622 aVal = floor( aVal * prec + 0.5 ) / prec;
629 double MaxElementLength3D::GetBadRate( double Value, int /*nbNodes*/ ) const
634 SMDSAbs_ElementType MaxElementLength3D::GetType() const
636 return SMDSAbs_Volume;
642 Description : Functor for calculation of minimum angle
645 double MinimumAngle::GetValue( const TSequenceOfXYZ& P )
652 aMin = getAngle(P( P.size() ), P( 1 ), P( 2 ));
653 aMin = Min(aMin,getAngle(P( P.size()-1 ), P( P.size() ), P( 1 )));
655 for (int i=2; i<P.size();i++){
656 double A0 = getAngle( P( i-1 ), P( i ), P( i+1 ) );
660 return aMin * 180.0 / PI;
663 double MinimumAngle::GetBadRate( double Value, int nbNodes ) const
665 //const double aBestAngle = PI / nbNodes;
666 const double aBestAngle = 180.0 - ( 360.0 / double(nbNodes) );
667 return ( fabs( aBestAngle - Value ));
670 SMDSAbs_ElementType MinimumAngle::GetType() const
678 Description : Functor for calculating aspect ratio
680 double AspectRatio::GetValue( const TSequenceOfXYZ& P )
682 // According to "Mesh quality control" by Nadir Bouhamau referring to
683 // Pascal Jean Frey and Paul-Louis George. Maillages, applications aux elements finis.
684 // Hermes Science publications, Paris 1999 ISBN 2-7462-0024-4
687 int nbNodes = P.size();
692 // Compute aspect ratio
694 if ( nbNodes == 3 ) {
695 // Compute lengths of the sides
696 std::vector< double > aLen (nbNodes);
697 for ( int i = 0; i < nbNodes - 1; i++ )
698 aLen[ i ] = getDistance( P( i + 1 ), P( i + 2 ) );
699 aLen[ nbNodes - 1 ] = getDistance( P( 1 ), P( nbNodes ) );
700 // Q = alfa * h * p / S, where
702 // alfa = sqrt( 3 ) / 6
703 // h - length of the longest edge
704 // p - half perimeter
705 // S - triangle surface
706 const double alfa = sqrt( 3. ) / 6.;
707 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
708 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
709 double anArea = getArea( P( 1 ), P( 2 ), P( 3 ) );
710 if ( anArea <= Precision::Confusion() )
712 return alfa * maxLen * half_perimeter / anArea;
714 else if ( nbNodes == 6 ) { // quadratic triangles
715 // Compute lengths of the sides
716 std::vector< double > aLen (3);
717 aLen[0] = getDistance( P(1), P(3) );
718 aLen[1] = getDistance( P(3), P(5) );
719 aLen[2] = getDistance( P(5), P(1) );
720 // Q = alfa * h * p / S, where
722 // alfa = sqrt( 3 ) / 6
723 // h - length of the longest edge
724 // p - half perimeter
725 // S - triangle surface
726 const double alfa = sqrt( 3. ) / 6.;
727 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
728 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
729 double anArea = getArea( P(1), P(3), P(5) );
730 if ( anArea <= Precision::Confusion() )
732 return alfa * maxLen * half_perimeter / anArea;
734 else if( nbNodes == 4 ) { // quadrangle
735 // Compute lengths of the sides
736 std::vector< double > aLen (4);
737 aLen[0] = getDistance( P(1), P(2) );
738 aLen[1] = getDistance( P(2), P(3) );
739 aLen[2] = getDistance( P(3), P(4) );
740 aLen[3] = getDistance( P(4), P(1) );
741 // Compute lengths of the diagonals
742 std::vector< double > aDia (2);
743 aDia[0] = getDistance( P(1), P(3) );
744 aDia[1] = getDistance( P(2), P(4) );
745 // Compute areas of all triangles which can be built
746 // taking three nodes of the quadrangle
747 std::vector< double > anArea (4);
748 anArea[0] = getArea( P(1), P(2), P(3) );
749 anArea[1] = getArea( P(1), P(2), P(4) );
750 anArea[2] = getArea( P(1), P(3), P(4) );
751 anArea[3] = getArea( P(2), P(3), P(4) );
752 // Q = alpha * L * C1 / C2, where
754 // alpha = sqrt( 1/32 )
755 // L = max( L1, L2, L3, L4, D1, D2 )
756 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
757 // C2 = min( S1, S2, S3, S4 )
758 // Li - lengths of the edges
759 // Di - lengths of the diagonals
760 // Si - areas of the triangles
761 const double alpha = sqrt( 1 / 32. );
762 double L = Max( aLen[ 0 ],
766 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
767 double C1 = sqrt( ( aLen[0] * aLen[0] +
770 aLen[3] * aLen[3] ) / 4. );
771 double C2 = Min( anArea[ 0 ],
773 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
774 if ( C2 <= Precision::Confusion() )
776 return alpha * L * C1 / C2;
778 else if( nbNodes == 8 ){ // nbNodes==8 - quadratic quadrangle
779 // Compute lengths of the sides
780 std::vector< double > aLen (4);
781 aLen[0] = getDistance( P(1), P(3) );
782 aLen[1] = getDistance( P(3), P(5) );
783 aLen[2] = getDistance( P(5), P(7) );
784 aLen[3] = getDistance( P(7), P(1) );
785 // Compute lengths of the diagonals
786 std::vector< double > aDia (2);
787 aDia[0] = getDistance( P(1), P(5) );
788 aDia[1] = getDistance( P(3), P(7) );
789 // Compute areas of all triangles which can be built
790 // taking three nodes of the quadrangle
791 std::vector< double > anArea (4);
792 anArea[0] = getArea( P(1), P(3), P(5) );
793 anArea[1] = getArea( P(1), P(3), P(7) );
794 anArea[2] = getArea( P(1), P(5), P(7) );
795 anArea[3] = getArea( P(3), P(5), P(7) );
796 // Q = alpha * L * C1 / C2, where
798 // alpha = sqrt( 1/32 )
799 // L = max( L1, L2, L3, L4, D1, D2 )
800 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
801 // C2 = min( S1, S2, S3, S4 )
802 // Li - lengths of the edges
803 // Di - lengths of the diagonals
804 // Si - areas of the triangles
805 const double alpha = sqrt( 1 / 32. );
806 double L = Max( aLen[ 0 ],
810 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
811 double C1 = sqrt( ( aLen[0] * aLen[0] +
814 aLen[3] * aLen[3] ) / 4. );
815 double C2 = Min( anArea[ 0 ],
817 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
818 if ( C2 <= Precision::Confusion() )
820 return alpha * L * C1 / C2;
825 double AspectRatio::GetBadRate( double Value, int /*nbNodes*/ ) const
827 // the aspect ratio is in the range [1.0,infinity]
830 return Value / 1000.;
833 SMDSAbs_ElementType AspectRatio::GetType() const
840 Class : AspectRatio3D
841 Description : Functor for calculating aspect ratio
845 inline double getHalfPerimeter(double theTria[3]){
846 return (theTria[0] + theTria[1] + theTria[2])/2.0;
849 inline double getArea(double theHalfPerim, double theTria[3]){
850 return sqrt(theHalfPerim*
851 (theHalfPerim-theTria[0])*
852 (theHalfPerim-theTria[1])*
853 (theHalfPerim-theTria[2]));
856 inline double getVolume(double theLen[6]){
857 double a2 = theLen[0]*theLen[0];
858 double b2 = theLen[1]*theLen[1];
859 double c2 = theLen[2]*theLen[2];
860 double d2 = theLen[3]*theLen[3];
861 double e2 = theLen[4]*theLen[4];
862 double f2 = theLen[5]*theLen[5];
863 double P = 4.0*a2*b2*d2;
864 double Q = a2*(b2+d2-e2)-b2*(a2+d2-f2)-d2*(a2+b2-c2);
865 double R = (b2+d2-e2)*(a2+d2-f2)*(a2+d2-f2);
866 return sqrt(P-Q+R)/12.0;
869 inline double getVolume2(double theLen[6]){
870 double a2 = theLen[0]*theLen[0];
871 double b2 = theLen[1]*theLen[1];
872 double c2 = theLen[2]*theLen[2];
873 double d2 = theLen[3]*theLen[3];
874 double e2 = theLen[4]*theLen[4];
875 double f2 = theLen[5]*theLen[5];
877 double P = a2*e2*(b2+c2+d2+f2-a2-e2);
878 double Q = b2*f2*(a2+c2+d2+e2-b2-f2);
879 double R = c2*d2*(a2+b2+e2+f2-c2-d2);
880 double S = a2*b2*d2+b2*c2*e2+a2*c2*f2+d2*e2*f2;
882 return sqrt(P+Q+R-S)/12.0;
885 inline double getVolume(const TSequenceOfXYZ& P){
886 gp_Vec aVec1( P( 2 ) - P( 1 ) );
887 gp_Vec aVec2( P( 3 ) - P( 1 ) );
888 gp_Vec aVec3( P( 4 ) - P( 1 ) );
889 gp_Vec anAreaVec( aVec1 ^ aVec2 );
890 return fabs(aVec3 * anAreaVec) / 6.0;
893 inline double getMaxHeight(double theLen[6])
895 double aHeight = std::max(theLen[0],theLen[1]);
896 aHeight = std::max(aHeight,theLen[2]);
897 aHeight = std::max(aHeight,theLen[3]);
898 aHeight = std::max(aHeight,theLen[4]);
899 aHeight = std::max(aHeight,theLen[5]);
905 double AspectRatio3D::GetValue( const TSequenceOfXYZ& P )
907 double aQuality = 0.0;
908 if(myCurrElement->IsPoly()) return aQuality;
910 int nbNodes = P.size();
912 if(myCurrElement->IsQuadratic()) {
913 if(nbNodes==10) nbNodes=4; // quadratic tetrahedron
914 else if(nbNodes==13) nbNodes=5; // quadratic pyramid
915 else if(nbNodes==15) nbNodes=6; // quadratic pentahedron
916 else if(nbNodes==20) nbNodes=8; // quadratic hexahedron
917 else return aQuality;
923 getDistance(P( 1 ),P( 2 )), // a
924 getDistance(P( 2 ),P( 3 )), // b
925 getDistance(P( 3 ),P( 1 )), // c
926 getDistance(P( 2 ),P( 4 )), // d
927 getDistance(P( 3 ),P( 4 )), // e
928 getDistance(P( 1 ),P( 4 )) // f
930 double aTria[4][3] = {
931 {aLen[0],aLen[1],aLen[2]}, // abc
932 {aLen[0],aLen[3],aLen[5]}, // adf
933 {aLen[1],aLen[3],aLen[4]}, // bde
934 {aLen[2],aLen[4],aLen[5]} // cef
936 double aSumArea = 0.0;
937 double aHalfPerimeter = getHalfPerimeter(aTria[0]);
938 double anArea = getArea(aHalfPerimeter,aTria[0]);
940 aHalfPerimeter = getHalfPerimeter(aTria[1]);
941 anArea = getArea(aHalfPerimeter,aTria[1]);
943 aHalfPerimeter = getHalfPerimeter(aTria[2]);
944 anArea = getArea(aHalfPerimeter,aTria[2]);
946 aHalfPerimeter = getHalfPerimeter(aTria[3]);
947 anArea = getArea(aHalfPerimeter,aTria[3]);
949 double aVolume = getVolume(P);
950 //double aVolume = getVolume(aLen);
951 double aHeight = getMaxHeight(aLen);
952 static double aCoeff = sqrt(2.0)/12.0;
953 if ( aVolume > DBL_MIN )
954 aQuality = aCoeff*aHeight*aSumArea/aVolume;
959 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 3 ),P( 5 )};
960 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
963 gp_XYZ aXYZ[4] = {P( 1 ),P( 3 ),P( 4 ),P( 5 )};
964 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
967 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 5 )};
968 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
971 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 4 ),P( 5 )};
972 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
978 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 6 )};
979 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
982 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 3 )};
983 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
986 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 6 )};
987 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
990 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
991 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
994 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 6 )};
995 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
998 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 3 )};
999 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1005 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1006 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1009 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 4 )};
1010 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1013 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 7 )};
1014 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1017 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 8 )};
1018 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1021 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 3 )};
1022 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1025 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 4 )};
1026 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1029 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 7 )};
1030 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1033 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 8 )};
1034 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1037 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 3 )};
1038 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1041 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 4 )};
1042 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1045 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 7 )};
1046 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1049 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 8 )};
1050 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1053 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 1 )};
1054 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1057 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 2 )};
1058 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1061 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 5 )};
1062 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1065 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 6 )};
1066 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1069 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 1 )};
1070 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1073 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 2 )};
1074 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1077 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 5 )};
1078 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1081 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 6 )};
1082 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1085 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 1 )};
1086 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1089 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1090 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1093 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 5 )};
1094 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1097 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 6 )};
1098 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1101 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1102 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1105 gp_XYZ aXYZ[4] = {P( 4 ),P( 5 ),P( 8 ),P( 2 )};
1106 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1109 gp_XYZ aXYZ[4] = {P( 1 ),P( 4 ),P( 5 ),P( 3 )};
1110 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1113 gp_XYZ aXYZ[4] = {P( 3 ),P( 6 ),P( 7 ),P( 1 )};
1114 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1117 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 6 ),P( 4 )};
1118 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1121 gp_XYZ aXYZ[4] = {P( 5 ),P( 6 ),P( 8 ),P( 3 )};
1122 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1125 gp_XYZ aXYZ[4] = {P( 7 ),P( 8 ),P( 6 ),P( 1 )};
1126 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1129 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 7 )};
1130 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1133 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 2 ),P( 5 )};
1134 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1139 if ( nbNodes > 4 ) {
1140 // avaluate aspect ratio of quadranle faces
1141 AspectRatio aspect2D;
1142 SMDS_VolumeTool::VolumeType type = SMDS_VolumeTool::GetType( nbNodes );
1143 int nbFaces = SMDS_VolumeTool::NbFaces( type );
1144 TSequenceOfXYZ points(4);
1145 for ( int i = 0; i < nbFaces; ++i ) { // loop on faces of a volume
1146 if ( SMDS_VolumeTool::NbFaceNodes( type, i ) != 4 )
1148 const int* pInd = SMDS_VolumeTool::GetFaceNodesIndices( type, i, true );
1149 for ( int p = 0; p < 4; ++p ) // loop on nodes of a quadranle face
1150 points( p + 1 ) = P( pInd[ p ] + 1 );
1151 aQuality = std::max( aQuality, aspect2D.GetValue( points ));
1157 double AspectRatio3D::GetBadRate( double Value, int /*nbNodes*/ ) const
1159 // the aspect ratio is in the range [1.0,infinity]
1162 return Value / 1000.;
1165 SMDSAbs_ElementType AspectRatio3D::GetType() const
1167 return SMDSAbs_Volume;
1173 Description : Functor for calculating warping
1175 double Warping::GetValue( const TSequenceOfXYZ& P )
1177 if ( P.size() != 4 )
1180 gp_XYZ G = ( P( 1 ) + P( 2 ) + P( 3 ) + P( 4 ) ) / 4.;
1182 double A1 = ComputeA( P( 1 ), P( 2 ), P( 3 ), G );
1183 double A2 = ComputeA( P( 2 ), P( 3 ), P( 4 ), G );
1184 double A3 = ComputeA( P( 3 ), P( 4 ), P( 1 ), G );
1185 double A4 = ComputeA( P( 4 ), P( 1 ), P( 2 ), G );
1187 return Max( Max( A1, A2 ), Max( A3, A4 ) );
1190 double Warping::ComputeA( const gp_XYZ& thePnt1,
1191 const gp_XYZ& thePnt2,
1192 const gp_XYZ& thePnt3,
1193 const gp_XYZ& theG ) const
1195 double aLen1 = gp_Pnt( thePnt1 ).Distance( gp_Pnt( thePnt2 ) );
1196 double aLen2 = gp_Pnt( thePnt2 ).Distance( gp_Pnt( thePnt3 ) );
1197 double L = Min( aLen1, aLen2 ) * 0.5;
1198 if ( L < Precision::Confusion())
1201 gp_XYZ GI = ( thePnt2 + thePnt1 ) / 2. - theG;
1202 gp_XYZ GJ = ( thePnt3 + thePnt2 ) / 2. - theG;
1203 gp_XYZ N = GI.Crossed( GJ );
1205 if ( N.Modulus() < gp::Resolution() )
1210 double H = ( thePnt2 - theG ).Dot( N );
1211 return asin( fabs( H / L ) ) * 180. / PI;
1214 double Warping::GetBadRate( double Value, int /*nbNodes*/ ) const
1216 // the warp is in the range [0.0,PI/2]
1217 // 0.0 = good (no warp)
1218 // PI/2 = bad (face pliee)
1222 SMDSAbs_ElementType Warping::GetType() const
1224 return SMDSAbs_Face;
1230 Description : Functor for calculating taper
1232 double Taper::GetValue( const TSequenceOfXYZ& P )
1234 if ( P.size() != 4 )
1238 double J1 = getArea( P( 4 ), P( 1 ), P( 2 ) ) / 2.;
1239 double J2 = getArea( P( 3 ), P( 1 ), P( 2 ) ) / 2.;
1240 double J3 = getArea( P( 2 ), P( 3 ), P( 4 ) ) / 2.;
1241 double J4 = getArea( P( 3 ), P( 4 ), P( 1 ) ) / 2.;
1243 double JA = 0.25 * ( J1 + J2 + J3 + J4 );
1244 if ( JA <= Precision::Confusion() )
1247 double T1 = fabs( ( J1 - JA ) / JA );
1248 double T2 = fabs( ( J2 - JA ) / JA );
1249 double T3 = fabs( ( J3 - JA ) / JA );
1250 double T4 = fabs( ( J4 - JA ) / JA );
1252 return Max( Max( T1, T2 ), Max( T3, T4 ) );
1255 double Taper::GetBadRate( double Value, int /*nbNodes*/ ) const
1257 // the taper is in the range [0.0,1.0]
1258 // 0.0 = good (no taper)
1259 // 1.0 = bad (les cotes opposes sont allignes)
1263 SMDSAbs_ElementType Taper::GetType() const
1265 return SMDSAbs_Face;
1271 Description : Functor for calculating skew in degrees
1273 static inline double skewAngle( const gp_XYZ& p1, const gp_XYZ& p2, const gp_XYZ& p3 )
1275 gp_XYZ p12 = ( p2 + p1 ) / 2.;
1276 gp_XYZ p23 = ( p3 + p2 ) / 2.;
1277 gp_XYZ p31 = ( p3 + p1 ) / 2.;
1279 gp_Vec v1( p31 - p2 ), v2( p12 - p23 );
1281 return v1.Magnitude() < gp::Resolution() || v2.Magnitude() < gp::Resolution() ? 0. : v1.Angle( v2 );
1284 double Skew::GetValue( const TSequenceOfXYZ& P )
1286 if ( P.size() != 3 && P.size() != 4 )
1290 static double PI2 = PI / 2.;
1291 if ( P.size() == 3 )
1293 double A0 = fabs( PI2 - skewAngle( P( 3 ), P( 1 ), P( 2 ) ) );
1294 double A1 = fabs( PI2 - skewAngle( P( 1 ), P( 2 ), P( 3 ) ) );
1295 double A2 = fabs( PI2 - skewAngle( P( 2 ), P( 3 ), P( 1 ) ) );
1297 return Max( A0, Max( A1, A2 ) ) * 180. / PI;
1301 gp_XYZ p12 = ( P( 1 ) + P( 2 ) ) / 2.;
1302 gp_XYZ p23 = ( P( 2 ) + P( 3 ) ) / 2.;
1303 gp_XYZ p34 = ( P( 3 ) + P( 4 ) ) / 2.;
1304 gp_XYZ p41 = ( P( 4 ) + P( 1 ) ) / 2.;
1306 gp_Vec v1( p34 - p12 ), v2( p23 - p41 );
1307 double A = v1.Magnitude() <= gp::Resolution() || v2.Magnitude() <= gp::Resolution()
1308 ? 0. : fabs( PI2 - v1.Angle( v2 ) );
1311 if ( A < Precision::Angular() )
1314 return A * 180. / PI;
1318 double Skew::GetBadRate( double Value, int /*nbNodes*/ ) const
1320 // the skew is in the range [0.0,PI/2].
1326 SMDSAbs_ElementType Skew::GetType() const
1328 return SMDSAbs_Face;
1334 Description : Functor for calculating area
1336 double Area::GetValue( const TSequenceOfXYZ& P )
1339 if ( P.size() > 2 ) {
1340 gp_Vec aVec1( P(2) - P(1) );
1341 gp_Vec aVec2( P(3) - P(1) );
1342 gp_Vec SumVec = aVec1 ^ aVec2;
1343 for (int i=4; i<=P.size(); i++) {
1344 gp_Vec aVec1( P(i-1) - P(1) );
1345 gp_Vec aVec2( P(i) - P(1) );
1346 gp_Vec tmp = aVec1 ^ aVec2;
1349 val = SumVec.Magnitude() * 0.5;
1354 double Area::GetBadRate( double Value, int /*nbNodes*/ ) const
1356 // meaningless as it is not a quality control functor
1360 SMDSAbs_ElementType Area::GetType() const
1362 return SMDSAbs_Face;
1368 Description : Functor for calculating length off edge
1370 double Length::GetValue( const TSequenceOfXYZ& P )
1372 switch ( P.size() ) {
1373 case 2: return getDistance( P( 1 ), P( 2 ) );
1374 case 3: return getDistance( P( 1 ), P( 2 ) ) + getDistance( P( 2 ), P( 3 ) );
1379 double Length::GetBadRate( double Value, int /*nbNodes*/ ) const
1381 // meaningless as it is not quality control functor
1385 SMDSAbs_ElementType Length::GetType() const
1387 return SMDSAbs_Edge;
1392 Description : Functor for calculating length of edge
1395 double Length2D::GetValue( long theElementId)
1399 //cout<<"Length2D::GetValue"<<endl;
1400 if (GetPoints(theElementId,P)){
1401 //for(int jj=1; jj<=P.size(); jj++)
1402 // cout<<"jj="<<jj<<" P("<<P(jj).X()<<","<<P(jj).Y()<<","<<P(jj).Z()<<")"<<endl;
1404 double aVal;// = GetValue( P );
1405 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
1406 SMDSAbs_ElementType aType = aElem->GetType();
1415 aVal = getDistance( P( 1 ), P( 2 ) );
1418 else if (len == 3){ // quadratic edge
1419 aVal = getDistance(P( 1 ),P( 3 )) + getDistance(P( 3 ),P( 2 ));
1423 if (len == 3){ // triangles
1424 double L1 = getDistance(P( 1 ),P( 2 ));
1425 double L2 = getDistance(P( 2 ),P( 3 ));
1426 double L3 = getDistance(P( 3 ),P( 1 ));
1427 aVal = Max(L1,Max(L2,L3));
1430 else if (len == 4){ // quadrangles
1431 double L1 = getDistance(P( 1 ),P( 2 ));
1432 double L2 = getDistance(P( 2 ),P( 3 ));
1433 double L3 = getDistance(P( 3 ),P( 4 ));
1434 double L4 = getDistance(P( 4 ),P( 1 ));
1435 aVal = Max(Max(L1,L2),Max(L3,L4));
1438 if (len == 6){ // quadratic triangles
1439 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1440 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1441 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
1442 aVal = Max(L1,Max(L2,L3));
1443 //cout<<"L1="<<L1<<" L2="<<L2<<"L3="<<L3<<" aVal="<<aVal<<endl;
1446 else if (len == 8){ // quadratic quadrangles
1447 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1448 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1449 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
1450 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
1451 aVal = Max(Max(L1,L2),Max(L3,L4));
1454 case SMDSAbs_Volume:
1455 if (len == 4){ // tetraidrs
1456 double L1 = getDistance(P( 1 ),P( 2 ));
1457 double L2 = getDistance(P( 2 ),P( 3 ));
1458 double L3 = getDistance(P( 3 ),P( 1 ));
1459 double L4 = getDistance(P( 1 ),P( 4 ));
1460 double L5 = getDistance(P( 2 ),P( 4 ));
1461 double L6 = getDistance(P( 3 ),P( 4 ));
1462 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1465 else if (len == 5){ // piramids
1466 double L1 = getDistance(P( 1 ),P( 2 ));
1467 double L2 = getDistance(P( 2 ),P( 3 ));
1468 double L3 = getDistance(P( 3 ),P( 4 ));
1469 double L4 = getDistance(P( 4 ),P( 1 ));
1470 double L5 = getDistance(P( 1 ),P( 5 ));
1471 double L6 = getDistance(P( 2 ),P( 5 ));
1472 double L7 = getDistance(P( 3 ),P( 5 ));
1473 double L8 = getDistance(P( 4 ),P( 5 ));
1475 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1476 aVal = Max(aVal,Max(L7,L8));
1479 else if (len == 6){ // pentaidres
1480 double L1 = getDistance(P( 1 ),P( 2 ));
1481 double L2 = getDistance(P( 2 ),P( 3 ));
1482 double L3 = getDistance(P( 3 ),P( 1 ));
1483 double L4 = getDistance(P( 4 ),P( 5 ));
1484 double L5 = getDistance(P( 5 ),P( 6 ));
1485 double L6 = getDistance(P( 6 ),P( 4 ));
1486 double L7 = getDistance(P( 1 ),P( 4 ));
1487 double L8 = getDistance(P( 2 ),P( 5 ));
1488 double L9 = getDistance(P( 3 ),P( 6 ));
1490 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1491 aVal = Max(aVal,Max(Max(L7,L8),L9));
1494 else if (len == 8){ // hexaider
1495 double L1 = getDistance(P( 1 ),P( 2 ));
1496 double L2 = getDistance(P( 2 ),P( 3 ));
1497 double L3 = getDistance(P( 3 ),P( 4 ));
1498 double L4 = getDistance(P( 4 ),P( 1 ));
1499 double L5 = getDistance(P( 5 ),P( 6 ));
1500 double L6 = getDistance(P( 6 ),P( 7 ));
1501 double L7 = getDistance(P( 7 ),P( 8 ));
1502 double L8 = getDistance(P( 8 ),P( 5 ));
1503 double L9 = getDistance(P( 1 ),P( 5 ));
1504 double L10= getDistance(P( 2 ),P( 6 ));
1505 double L11= getDistance(P( 3 ),P( 7 ));
1506 double L12= getDistance(P( 4 ),P( 8 ));
1508 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1509 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1510 aVal = Max(aVal,Max(L11,L12));
1515 if (len == 10){ // quadratic tetraidrs
1516 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
1517 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
1518 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
1519 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1520 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
1521 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
1522 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1525 else if (len == 13){ // quadratic piramids
1526 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
1527 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
1528 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1529 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1530 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1531 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
1532 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
1533 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
1534 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1535 aVal = Max(aVal,Max(L7,L8));
1538 else if (len == 15){ // quadratic pentaidres
1539 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
1540 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
1541 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1542 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1543 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
1544 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
1545 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
1546 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
1547 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
1548 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1549 aVal = Max(aVal,Max(Max(L7,L8),L9));
1552 else if (len == 20){ // quadratic hexaider
1553 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
1554 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
1555 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
1556 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
1557 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
1558 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
1559 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
1560 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
1561 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
1562 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
1563 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
1564 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
1565 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1566 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1567 aVal = Max(aVal,Max(L11,L12));
1579 if ( myPrecision >= 0 )
1581 double prec = pow( 10., (double)( myPrecision ) );
1582 aVal = floor( aVal * prec + 0.5 ) / prec;
1591 double Length2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1593 // meaningless as it is not quality control functor
1597 SMDSAbs_ElementType Length2D::GetType() const
1599 return SMDSAbs_Face;
1602 Length2D::Value::Value(double theLength,long thePntId1, long thePntId2):
1605 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1606 if(thePntId1 > thePntId2){
1607 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1611 bool Length2D::Value::operator<(const Length2D::Value& x) const{
1612 if(myPntId[0] < x.myPntId[0]) return true;
1613 if(myPntId[0] == x.myPntId[0])
1614 if(myPntId[1] < x.myPntId[1]) return true;
1618 void Length2D::GetValues(TValues& theValues){
1620 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1621 for(; anIter->more(); ){
1622 const SMDS_MeshFace* anElem = anIter->next();
1624 if(anElem->IsQuadratic()) {
1625 const SMDS_QuadraticFaceOfNodes* F =
1626 static_cast<const SMDS_QuadraticFaceOfNodes*>(anElem);
1627 // use special nodes iterator
1628 SMDS_NodeIteratorPtr anIter = F->interlacedNodesIterator();
1633 const SMDS_MeshElement* aNode;
1635 aNode = anIter->next();
1636 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1637 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1638 aNodeId[0] = aNodeId[1] = aNode->GetID();
1641 for(; anIter->more(); ){
1642 const SMDS_MeshNode* N1 = static_cast<const SMDS_MeshNode*> (anIter->next());
1643 P[2] = gp_Pnt(N1->X(),N1->Y(),N1->Z());
1644 aNodeId[2] = N1->GetID();
1645 aLength = P[1].Distance(P[2]);
1646 if(!anIter->more()) break;
1647 const SMDS_MeshNode* N2 = static_cast<const SMDS_MeshNode*> (anIter->next());
1648 P[3] = gp_Pnt(N2->X(),N2->Y(),N2->Z());
1649 aNodeId[3] = N2->GetID();
1650 aLength += P[2].Distance(P[3]);
1651 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1652 Value aValue2(aLength,aNodeId[2],aNodeId[3]);
1654 aNodeId[1] = aNodeId[3];
1655 theValues.insert(aValue1);
1656 theValues.insert(aValue2);
1658 aLength += P[2].Distance(P[0]);
1659 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1660 Value aValue2(aLength,aNodeId[2],aNodeId[0]);
1661 theValues.insert(aValue1);
1662 theValues.insert(aValue2);
1665 SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
1670 const SMDS_MeshElement* aNode;
1671 if(aNodesIter->more()){
1672 aNode = aNodesIter->next();
1673 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1674 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1675 aNodeId[0] = aNodeId[1] = aNode->GetID();
1678 for(; aNodesIter->more(); ){
1679 aNode = aNodesIter->next();
1680 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1681 long anId = aNode->GetID();
1683 P[2] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1685 aLength = P[1].Distance(P[2]);
1687 Value aValue(aLength,aNodeId[1],anId);
1690 theValues.insert(aValue);
1693 aLength = P[0].Distance(P[1]);
1695 Value aValue(aLength,aNodeId[0],aNodeId[1]);
1696 theValues.insert(aValue);
1702 Class : MultiConnection
1703 Description : Functor for calculating number of faces conneted to the edge
1705 double MultiConnection::GetValue( const TSequenceOfXYZ& P )
1709 double MultiConnection::GetValue( long theId )
1711 return getNbMultiConnection( myMesh, theId );
1714 double MultiConnection::GetBadRate( double Value, int /*nbNodes*/ ) const
1716 // meaningless as it is not quality control functor
1720 SMDSAbs_ElementType MultiConnection::GetType() const
1722 return SMDSAbs_Edge;
1726 Class : MultiConnection2D
1727 Description : Functor for calculating number of faces conneted to the edge
1729 double MultiConnection2D::GetValue( const TSequenceOfXYZ& P )
1734 double MultiConnection2D::GetValue( long theElementId )
1738 const SMDS_MeshElement* aFaceElem = myMesh->FindElement(theElementId);
1739 SMDSAbs_ElementType aType = aFaceElem->GetType();
1744 int i = 0, len = aFaceElem->NbNodes();
1745 SMDS_ElemIteratorPtr anIter = aFaceElem->nodesIterator();
1748 const SMDS_MeshNode *aNode, *aNode0;
1749 TColStd_MapOfInteger aMap, aMapPrev;
1751 for (i = 0; i <= len; i++) {
1756 if (anIter->more()) {
1757 aNode = (SMDS_MeshNode*)anIter->next();
1765 if (i == 0) aNode0 = aNode;
1767 SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
1768 while (anElemIter->more()) {
1769 const SMDS_MeshElement* anElem = anElemIter->next();
1770 if (anElem != 0 && anElem->GetType() == SMDSAbs_Face) {
1771 int anId = anElem->GetID();
1774 if (aMapPrev.Contains(anId)) {
1779 aResult = Max(aResult, aNb);
1790 double MultiConnection2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1792 // meaningless as it is not quality control functor
1796 SMDSAbs_ElementType MultiConnection2D::GetType() const
1798 return SMDSAbs_Face;
1801 MultiConnection2D::Value::Value(long thePntId1, long thePntId2)
1803 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1804 if(thePntId1 > thePntId2){
1805 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1809 bool MultiConnection2D::Value::operator<(const MultiConnection2D::Value& x) const{
1810 if(myPntId[0] < x.myPntId[0]) return true;
1811 if(myPntId[0] == x.myPntId[0])
1812 if(myPntId[1] < x.myPntId[1]) return true;
1816 void MultiConnection2D::GetValues(MValues& theValues){
1817 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1818 for(; anIter->more(); ){
1819 const SMDS_MeshFace* anElem = anIter->next();
1820 SMDS_ElemIteratorPtr aNodesIter;
1821 if ( anElem->IsQuadratic() )
1822 aNodesIter = static_cast<const SMDS_QuadraticFaceOfNodes*>
1823 (anElem)->interlacedNodesElemIterator();
1825 aNodesIter = anElem->nodesIterator();
1828 //int aNbConnects=0;
1829 const SMDS_MeshNode* aNode0;
1830 const SMDS_MeshNode* aNode1;
1831 const SMDS_MeshNode* aNode2;
1832 if(aNodesIter->more()){
1833 aNode0 = (SMDS_MeshNode*) aNodesIter->next();
1835 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode1;
1836 aNodeId[0] = aNodeId[1] = aNodes->GetID();
1838 for(; aNodesIter->more(); ) {
1839 aNode2 = (SMDS_MeshNode*) aNodesIter->next();
1840 long anId = aNode2->GetID();
1843 Value aValue(aNodeId[1],aNodeId[2]);
1844 MValues::iterator aItr = theValues.find(aValue);
1845 if (aItr != theValues.end()){
1850 theValues[aValue] = 1;
1853 //cout << "NodeIds: "<<aNodeId[1]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1854 aNodeId[1] = aNodeId[2];
1857 Value aValue(aNodeId[0],aNodeId[2]);
1858 MValues::iterator aItr = theValues.find(aValue);
1859 if (aItr != theValues.end()) {
1864 theValues[aValue] = 1;
1867 //cout << "NodeIds: "<<aNodeId[0]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1877 Class : BadOrientedVolume
1878 Description : Predicate bad oriented volumes
1881 BadOrientedVolume::BadOrientedVolume()
1886 void BadOrientedVolume::SetMesh( const SMDS_Mesh* theMesh )
1891 bool BadOrientedVolume::IsSatisfy( long theId )
1896 SMDS_VolumeTool vTool( myMesh->FindElement( theId ));
1897 return !vTool.IsForward();
1900 SMDSAbs_ElementType BadOrientedVolume::GetType() const
1902 return SMDSAbs_Volume;
1909 Description : Predicate for free borders
1912 FreeBorders::FreeBorders()
1917 void FreeBorders::SetMesh( const SMDS_Mesh* theMesh )
1922 bool FreeBorders::IsSatisfy( long theId )
1924 return getNbMultiConnection( myMesh, theId ) == 1;
1927 SMDSAbs_ElementType FreeBorders::GetType() const
1929 return SMDSAbs_Edge;
1935 Description : Predicate for free Edges
1937 FreeEdges::FreeEdges()
1942 void FreeEdges::SetMesh( const SMDS_Mesh* theMesh )
1947 bool FreeEdges::IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId )
1949 TColStd_MapOfInteger aMap;
1950 for ( int i = 0; i < 2; i++ )
1952 SMDS_ElemIteratorPtr anElemIter = theNodes[ i ]->GetInverseElementIterator();
1953 while( anElemIter->more() )
1955 const SMDS_MeshElement* anElem = anElemIter->next();
1956 if ( anElem != 0 && anElem->GetType() == SMDSAbs_Face )
1958 int anId = anElem->GetID();
1962 else if ( aMap.Contains( anId ) && anId != theFaceId )
1970 bool FreeEdges::IsSatisfy( long theId )
1975 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
1976 if ( aFace == 0 || aFace->GetType() != SMDSAbs_Face || aFace->NbNodes() < 3 )
1979 SMDS_ElemIteratorPtr anIter;
1980 if ( aFace->IsQuadratic() ) {
1981 anIter = static_cast<const SMDS_QuadraticFaceOfNodes*>
1982 (aFace)->interlacedNodesElemIterator();
1985 anIter = aFace->nodesIterator();
1990 int i = 0, nbNodes = aFace->NbNodes();
1991 std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
1992 while( anIter->more() )
1994 const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
1997 aNodes[ i++ ] = aNode;
1999 aNodes[ nbNodes ] = aNodes[ 0 ];
2001 for ( i = 0; i < nbNodes; i++ )
2002 if ( IsFreeEdge( &aNodes[ i ], theId ) )
2008 SMDSAbs_ElementType FreeEdges::GetType() const
2010 return SMDSAbs_Face;
2013 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2016 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
2017 if(thePntId1 > thePntId2){
2018 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
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;
2029 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2030 FreeEdges::TBorders& theRegistry,
2031 FreeEdges::TBorders& theContainer)
2033 if(theRegistry.find(theBorder) == theRegistry.end()){
2034 theRegistry.insert(theBorder);
2035 theContainer.insert(theBorder);
2037 theContainer.erase(theBorder);
2041 void FreeEdges::GetBoreders(TBorders& theBorders)
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_QuadraticFaceOfNodes*>(anElem)->
2051 interlacedNodesElemIterator();
2053 aNodesIter = anElem->nodesIterator();
2055 const SMDS_MeshElement* aNode;
2056 if(aNodesIter->more()){
2057 aNode = aNodesIter->next();
2058 aNodeId[0] = aNodeId[1] = aNode->GetID();
2060 for(; aNodesIter->more(); ){
2061 aNode = aNodesIter->next();
2062 long anId = aNode->GetID();
2063 Border aBorder(anElemId,aNodeId[1],anId);
2065 //std::cout<<aBorder.myPntId[0]<<"; "<<aBorder.myPntId[1]<<"; "<<aBorder.myElemId<<endl;
2066 UpdateBorders(aBorder,aRegistry,theBorders);
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);
2072 //std::cout<<"theBorders.size() = "<<theBorders.size()<<endl;
2078 Description : Predicate for free nodes
2081 FreeNodes::FreeNodes()
2086 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2091 bool FreeNodes::IsSatisfy( long theNodeId )
2093 const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2097 return (aNode->NbInverseElements() < 1);
2100 SMDSAbs_ElementType FreeNodes::GetType() const
2102 return SMDSAbs_Node;
2108 Description : Predicate for free faces
2111 FreeFaces::FreeFaces()
2116 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2121 bool FreeFaces::IsSatisfy( long theId )
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 )
2129 int nbNode = aFace->NbNodes();
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;
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;
2148 TItrMapOfVolume volItr = mapOfVol.begin();
2149 TItrMapOfVolume volEnd = mapOfVol.end();
2150 for ( ; volItr != volEnd; ++volItr )
2151 if ( (*volItr).second >= nbNode )
2153 // face is not free if number of volumes constructed on thier nodes more than one
2157 SMDSAbs_ElementType FreeFaces::GetType() const
2159 return SMDSAbs_Face;
2163 Class : LinearOrQuadratic
2164 Description : Predicate to verify whether a mesh element is linear
2167 LinearOrQuadratic::LinearOrQuadratic()
2172 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2177 bool LinearOrQuadratic::IsSatisfy( long theId )
2179 if (!myMesh) return false;
2180 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2181 if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2183 return (!anElem->IsQuadratic());
2186 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2191 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2198 Description : Functor for check color of group to whic mesh element belongs to
2201 GroupColor::GroupColor()
2205 bool GroupColor::IsSatisfy( long theId )
2207 return (myIDs.find( theId ) != myIDs.end());
2210 void GroupColor::SetType( SMDSAbs_ElementType theType )
2215 SMDSAbs_ElementType GroupColor::GetType() const
2220 static bool isEqual( const Quantity_Color& theColor1,
2221 const Quantity_Color& theColor2 )
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 );
2231 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2235 const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2239 int nbGrp = aMesh->GetNbGroups();
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);
2250 // check type and color of group
2251 if ( !isEqual( myColor, aGrp->GetColor() ) )
2253 if ( myType != SMDSAbs_All && myType != (SMDSAbs_ElementType)aGrp->GetType() )
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) );
2266 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
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();
2280 myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
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
2291 theResStr += TCollection_AsciiString( myColor.Red() );
2292 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2293 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2297 Class : ElemGeomType
2298 Description : Predicate to check element geometry type
2301 ElemGeomType::ElemGeomType()
2304 myType = SMDSAbs_All;
2305 myGeomType = SMDSGeom_TRIANGLE;
2308 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2313 bool ElemGeomType::IsSatisfy( long theId )
2315 if (!myMesh) return false;
2316 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2319 const SMDSAbs_ElementType anElemType = anElem->GetType();
2320 if ( myType != SMDSAbs_All && anElemType != myType )
2322 const int aNbNode = anElem->NbNodes();
2324 switch( anElemType )
2327 isOk = (myGeomType == SMDSGeom_POINT);
2331 isOk = (myGeomType == SMDSGeom_EDGE);
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();
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();
2360 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2365 SMDSAbs_ElementType ElemGeomType::GetType() const
2370 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2372 myGeomType = theType;
2375 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2380 //================================================================================
2382 * \brief Class CoplanarFaces
2384 //================================================================================
2386 CoplanarFaces::CoplanarFaces()
2387 : myMesh(0), myFaceID(0), myToler(0)
2390 bool CoplanarFaces::IsSatisfy( long theElementId )
2392 if ( myCoplanarIDs.empty() )
2394 // Build a set of coplanar face ids
2396 if ( !myMesh || !myFaceID || !myToler )
2399 const SMDS_MeshElement* face = myMesh->FindElement( myFaceID );
2400 if ( !face || face->GetType() != SMDSAbs_Face )
2404 gp_Vec myNorm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
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() )
2414 face = faceQueue.front();
2415 if ( checkedFaces.insert( face ).second )
2417 gp_Vec norm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2418 if (!normOK || myNorm.Angle( norm ) <= radianTol)
2420 myCoplanarIDs.insert( face->GetID() );
2421 std::set<const SMDS_MeshElement*> neighborFaces;
2422 for ( int i = 0; i < face->NbCornerNodes(); ++i )
2424 const SMDS_MeshNode* n = face->GetNode( i );
2425 if ( checkedNodes.insert( n ).second )
2426 neighborFaces.insert( TFaceIt( n->GetInverseElementIterator(SMDSAbs_Face)),
2429 faceQueue.insert( faceQueue.end(), neighborFaces.begin(), neighborFaces.end() );
2432 faceQueue.pop_front();
2435 return myCoplanarIDs.count( theElementId );
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-"
2447 //=======================================================================
2448 // name : RangeOfIds
2449 // Purpose : Constructor
2450 //=======================================================================
2451 RangeOfIds::RangeOfIds()
2454 myType = SMDSAbs_All;
2457 //=======================================================================
2459 // Purpose : Set mesh
2460 //=======================================================================
2461 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
2466 //=======================================================================
2467 // name : AddToRange
2468 // Purpose : Add ID to the range
2469 //=======================================================================
2470 bool RangeOfIds::AddToRange( long theEntityId )
2472 myIds.Add( theEntityId );
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 )
2485 TColStd_SequenceOfInteger anIntSeq;
2486 TColStd_SequenceOfAsciiString aStrSeq;
2488 TColStd_MapIteratorOfMapOfInteger anIter( myIds );
2489 for ( ; anIter.More(); anIter.Next() )
2491 int anId = anIter.Key();
2492 TCollection_AsciiString aStr( anId );
2493 anIntSeq.Append( anId );
2494 aStrSeq.Append( aStr );
2497 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2499 int aMinId = myMin( i );
2500 int aMaxId = myMax( i );
2502 TCollection_AsciiString aStr;
2503 if ( aMinId != IntegerFirst() )
2508 if ( aMaxId != IntegerLast() )
2511 // find position of the string in result sequence and insert string in it
2512 if ( anIntSeq.Length() == 0 )
2514 anIntSeq.Append( aMinId );
2515 aStrSeq.Append( aStr );
2519 if ( aMinId < anIntSeq.First() )
2521 anIntSeq.Prepend( aMinId );
2522 aStrSeq.Prepend( aStr );
2524 else if ( aMinId > anIntSeq.Last() )
2526 anIntSeq.Append( aMinId );
2527 aStrSeq.Append( aStr );
2530 for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
2531 if ( aMinId < anIntSeq( j ) )
2533 anIntSeq.InsertBefore( j, aMinId );
2534 aStrSeq.InsertBefore( j, aStr );
2540 if ( aStrSeq.Length() == 0 )
2543 theResStr = aStrSeq( 1 );
2544 for ( int j = 2, k = aStrSeq.Length(); j <= k; j++ )
2547 theResStr += aStrSeq( j );
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 )
2562 TCollection_AsciiString aStr = theStr;
2563 aStr.RemoveAll( ' ' );
2564 aStr.RemoveAll( '\t' );
2566 for ( int aPos = aStr.Search( ",," ); aPos != -1; aPos = aStr.Search( ",," ) )
2567 aStr.Remove( aPos, 2 );
2569 TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
2571 while ( tmpStr != "" )
2573 tmpStr = aStr.Token( ",", i++ );
2574 int aPos = tmpStr.Search( '-' );
2578 if ( tmpStr.IsIntegerValue() )
2579 myIds.Add( tmpStr.IntegerValue() );
2585 TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
2586 TCollection_AsciiString aMinStr = tmpStr;
2588 while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
2589 while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
2591 if ( !aMinStr.IsEmpty() && !aMinStr.IsIntegerValue() ||
2592 !aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue() )
2595 myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
2596 myMax.Append( aMaxStr.IsEmpty() ? IntegerLast() : aMaxStr.IntegerValue() );
2603 //=======================================================================
2605 // Purpose : Get type of supported entities
2606 //=======================================================================
2607 SMDSAbs_ElementType RangeOfIds::GetType() const
2612 //=======================================================================
2614 // Purpose : Set type of supported entities
2615 //=======================================================================
2616 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
2621 //=======================================================================
2623 // Purpose : Verify whether entity satisfies to this rpedicate
2624 //=======================================================================
2625 bool RangeOfIds::IsSatisfy( long theId )
2630 if ( myType == SMDSAbs_Node )
2632 if ( myMesh->FindNode( theId ) == 0 )
2637 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2638 if ( anElem == 0 || myType != anElem->GetType() && myType != SMDSAbs_All )
2642 if ( myIds.Contains( theId ) )
2645 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2646 if ( theId >= myMin( i ) && theId <= myMax( i ) )
2654 Description : Base class for comparators
2656 Comparator::Comparator():
2660 Comparator::~Comparator()
2663 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
2666 myFunctor->SetMesh( theMesh );
2669 void Comparator::SetMargin( double theValue )
2671 myMargin = theValue;
2674 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
2676 myFunctor = theFunct;
2679 SMDSAbs_ElementType Comparator::GetType() const
2681 return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
2684 double Comparator::GetMargin()
2692 Description : Comparator "<"
2694 bool LessThan::IsSatisfy( long theId )
2696 return myFunctor && myFunctor->GetValue( theId ) < myMargin;
2702 Description : Comparator ">"
2704 bool MoreThan::IsSatisfy( long theId )
2706 return myFunctor && myFunctor->GetValue( theId ) > myMargin;
2712 Description : Comparator "="
2715 myToler(Precision::Confusion())
2718 bool EqualTo::IsSatisfy( long theId )
2720 return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
2723 void EqualTo::SetTolerance( double theToler )
2728 double EqualTo::GetTolerance()
2735 Description : Logical NOT predicate
2737 LogicalNOT::LogicalNOT()
2740 LogicalNOT::~LogicalNOT()
2743 bool LogicalNOT::IsSatisfy( long theId )
2745 return myPredicate && !myPredicate->IsSatisfy( theId );
2748 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
2751 myPredicate->SetMesh( theMesh );
2754 void LogicalNOT::SetPredicate( PredicatePtr thePred )
2756 myPredicate = thePred;
2759 SMDSAbs_ElementType LogicalNOT::GetType() const
2761 return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
2766 Class : LogicalBinary
2767 Description : Base class for binary logical predicate
2769 LogicalBinary::LogicalBinary()
2772 LogicalBinary::~LogicalBinary()
2775 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
2778 myPredicate1->SetMesh( theMesh );
2781 myPredicate2->SetMesh( theMesh );
2784 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
2786 myPredicate1 = thePredicate;
2789 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
2791 myPredicate2 = thePredicate;
2794 SMDSAbs_ElementType LogicalBinary::GetType() const
2796 if ( !myPredicate1 || !myPredicate2 )
2799 SMDSAbs_ElementType aType1 = myPredicate1->GetType();
2800 SMDSAbs_ElementType aType2 = myPredicate2->GetType();
2802 return aType1 == aType2 ? aType1 : SMDSAbs_All;
2808 Description : Logical AND
2810 bool LogicalAND::IsSatisfy( long theId )
2815 myPredicate1->IsSatisfy( theId ) &&
2816 myPredicate2->IsSatisfy( theId );
2822 Description : Logical OR
2824 bool LogicalOR::IsSatisfy( long theId )
2829 myPredicate1->IsSatisfy( theId ) ||
2830 myPredicate2->IsSatisfy( theId );
2844 void Filter::SetPredicate( PredicatePtr thePredicate )
2846 myPredicate = thePredicate;
2849 template<class TElement, class TIterator, class TPredicate>
2850 inline void FillSequence(const TIterator& theIterator,
2851 TPredicate& thePredicate,
2852 Filter::TIdSequence& theSequence)
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 );
2866 GetElementsId( const SMDS_Mesh* theMesh,
2867 PredicatePtr thePredicate,
2868 TIdSequence& theSequence )
2870 theSequence.clear();
2872 if ( !theMesh || !thePredicate )
2875 thePredicate->SetMesh( theMesh );
2877 SMDSAbs_ElementType aType = thePredicate->GetType();
2880 FillSequence<const SMDS_MeshNode*>(theMesh->nodesIterator(),thePredicate,theSequence);
2883 FillSequence<const SMDS_MeshElement*>(theMesh->edgesIterator(),thePredicate,theSequence);
2886 FillSequence<const SMDS_MeshElement*>(theMesh->facesIterator(),thePredicate,theSequence);
2888 case SMDSAbs_Volume:
2889 FillSequence<const SMDS_MeshElement*>(theMesh->volumesIterator(),thePredicate,theSequence);
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);
2900 Filter::GetElementsId( const SMDS_Mesh* theMesh,
2901 Filter::TIdSequence& theSequence )
2903 GetElementsId(theMesh,myPredicate,theSequence);
2910 typedef std::set<SMDS_MeshFace*> TMapOfFacePtr;
2916 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
2917 SMDS_MeshNode* theNode2 )
2923 ManifoldPart::Link::~Link()
2929 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
2931 if ( myNode1 == theLink.myNode1 &&
2932 myNode2 == theLink.myNode2 )
2934 else if ( myNode1 == theLink.myNode2 &&
2935 myNode2 == theLink.myNode1 )
2941 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
2943 if(myNode1 < x.myNode1) return true;
2944 if(myNode1 == x.myNode1)
2945 if(myNode2 < x.myNode2) return true;
2949 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
2950 const ManifoldPart::Link& theLink2 )
2952 return theLink1.IsEqual( theLink2 );
2955 ManifoldPart::ManifoldPart()
2958 myAngToler = Precision::Angular();
2959 myIsOnlyManifold = true;
2962 ManifoldPart::~ManifoldPart()
2967 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
2973 SMDSAbs_ElementType ManifoldPart::GetType() const
2974 { return SMDSAbs_Face; }
2976 bool ManifoldPart::IsSatisfy( long theElementId )
2978 return myMapIds.Contains( theElementId );
2981 void ManifoldPart::SetAngleTolerance( const double theAngToler )
2982 { myAngToler = theAngToler; }
2984 double ManifoldPart::GetAngleTolerance() const
2985 { return myAngToler; }
2987 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
2988 { myIsOnlyManifold = theIsOnly; }
2990 void ManifoldPart::SetStartElem( const long theStartId )
2991 { myStartElemId = theStartId; }
2993 bool ManifoldPart::process()
2996 myMapBadGeomIds.Clear();
2998 myAllFacePtr.clear();
2999 myAllFacePtrIntDMap.clear();
3003 // collect all faces into own map
3004 SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3005 for (; anFaceItr->more(); )
3007 SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3008 myAllFacePtr.push_back( aFacePtr );
3009 myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3012 SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3016 // the map of non manifold links and bad geometry
3017 TMapOfLink aMapOfNonManifold;
3018 TColStd_MapOfInteger aMapOfTreated;
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++ )
3026 if ( fi == aStartIndx )
3027 isStartTreat = true;
3028 // as result next time when fi will be equal to aStartIndx
3030 SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3031 if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3034 aMapOfTreated.Add( aFacePtr->GetID() );
3035 TColStd_MapOfInteger aResFaces;
3036 if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3037 aMapOfNonManifold, aResFaces ) )
3039 TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3040 for ( ; anItr.More(); anItr.Next() )
3042 int aFaceId = anItr.Key();
3043 aMapOfTreated.Add( aFaceId );
3044 myMapIds.Add( aFaceId );
3047 if ( fi == ( myAllFacePtr.size() - 1 ) )
3049 } // end run on vector of faces
3050 return !myMapIds.IsEmpty();
3053 static void getLinks( const SMDS_MeshFace* theFace,
3054 ManifoldPart::TVectorOfLink& theLinks )
3056 int aNbNode = theFace->NbNodes();
3057 SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3059 SMDS_MeshNode* aNode = 0;
3060 for ( ; aNodeItr->more() && i <= aNbNode; )
3063 SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3067 SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3069 ManifoldPart::Link aLink( aN1, aN2 );
3070 theLinks.push_back( aLink );
3074 bool ManifoldPart::findConnected
3075 ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3076 SMDS_MeshFace* theStartFace,
3077 ManifoldPart::TMapOfLink& theNonManifold,
3078 TColStd_MapOfInteger& theResFaces )
3080 theResFaces.Clear();
3081 if ( !theAllFacePtrInt.size() )
3084 if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3086 myMapBadGeomIds.Add( theStartFace->GetID() );
3090 ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3091 ManifoldPart::TVectorOfLink aSeqOfBoundary;
3092 theResFaces.Add( theStartFace->GetID() );
3093 ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3095 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3096 aDMapLinkFace, theNonManifold, theStartFace );
3098 bool isDone = false;
3099 while ( !isDone && aMapOfBoundary.size() != 0 )
3101 bool isToReset = false;
3102 ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3103 for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3105 ManifoldPart::Link aLink = *pLink;
3106 if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3108 // each link could be treated only once
3109 aMapToSkip.insert( aLink );
3111 ManifoldPart::TVectorOfFacePtr aFaces;
3113 if ( myIsOnlyManifold &&
3114 (theNonManifold.find( aLink ) != theNonManifold.end()) )
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 )
3124 SMDS_MeshFace* aFace = *pFace;
3125 if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3126 aFiltered.push_back( aFace );
3129 if ( aFaces.size() < 2 ) // no neihgbour faces
3131 else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3133 theNonManifold.insert( aLink );
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 )
3143 SMDS_MeshFace* aNextFace = *pFace;
3144 if ( aPrevFace == aNextFace )
3146 int anNextFaceID = aNextFace->GetID();
3147 if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3148 // should not be with non manifold restriction. probably bad topology
3150 // check if face was treated and skipped
3151 if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3152 !isInPlane( aPrevFace, aNextFace ) )
3154 // add new element to connected and extend the boundaries.
3155 theResFaces.Add( anNextFaceID );
3156 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3157 aDMapLinkFace, theNonManifold, aNextFace );
3161 isDone = !isToReset;
3164 return !theResFaces.IsEmpty();
3167 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3168 const SMDS_MeshFace* theFace2 )
3170 gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3171 gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3172 if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3174 myMapBadGeomIds.Add( theFace2->GetID() );
3177 if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
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
3190 ManifoldPart::TVectorOfLink aLinks;
3191 getLinks( theNextFace, aLinks );
3192 int aNbLink = (int)aLinks.size();
3193 for ( int i = 0; i < aNbLink; i++ )
3195 ManifoldPart::Link aLink = aLinks[ i ];
3196 if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3198 if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3200 if ( myIsOnlyManifold )
3202 // remove from boundary
3203 theMapOfBoundary.erase( aLink );
3204 ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3205 for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3207 ManifoldPart::Link aBoundLink = *pLink;
3208 if ( aBoundLink.IsEqual( aLink ) )
3210 theSeqOfBoundary.erase( pLink );
3218 theMapOfBoundary.insert( aLink );
3219 theSeqOfBoundary.push_back( aLink );
3220 theDMapLinkFacePtr[ aLink ] = theNextFace;
3225 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3226 ManifoldPart::TVectorOfFacePtr& theFaces ) const
3228 SMDS_Mesh::SetOfFaces aSetOfFaces;
3229 // take all faces that shared first node
3230 SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3231 for ( ; anItr->more(); )
3233 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3236 aSetOfFaces.Add( aFace );
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(); )
3243 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3244 if ( aSetOfFaces.Contains( aFace ) )
3245 theFaces.push_back( aFace );
3254 ElementsOnSurface::ElementsOnSurface()
3258 myType = SMDSAbs_All;
3260 myToler = Precision::Confusion();
3261 myUseBoundaries = false;
3264 ElementsOnSurface::~ElementsOnSurface()
3269 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
3271 if ( myMesh == theMesh )
3277 bool ElementsOnSurface::IsSatisfy( long theElementId )
3279 return myIds.Contains( theElementId );
3282 SMDSAbs_ElementType ElementsOnSurface::GetType() const
3285 void ElementsOnSurface::SetTolerance( const double theToler )
3287 if ( myToler != theToler )
3292 double ElementsOnSurface::GetTolerance() const
3295 void ElementsOnSurface::SetUseBoundaries( bool theUse )
3297 if ( myUseBoundaries != theUse ) {
3298 myUseBoundaries = theUse;
3299 SetSurface( mySurf, myType );
3303 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
3304 const SMDSAbs_ElementType theType )
3309 if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
3311 mySurf = TopoDS::Face( theShape );
3312 BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
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 );
3323 void ElementsOnSurface::process()
3326 if ( mySurf.IsNull() )
3332 if ( myType == SMDSAbs_Face || myType == SMDSAbs_All )
3334 myIds.ReSize( myMesh->NbFaces() );
3335 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
3336 for(; anIter->more(); )
3337 process( anIter->next() );
3340 if ( myType == SMDSAbs_Edge || myType == SMDSAbs_All )
3342 myIds.ReSize( myIds.Extent() + myMesh->NbEdges() );
3343 SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
3344 for(; anIter->more(); )
3345 process( anIter->next() );
3348 if ( myType == SMDSAbs_Node )
3350 myIds.ReSize( myMesh->NbNodes() );
3351 SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
3352 for(; anIter->more(); )
3353 process( anIter->next() );
3357 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
3359 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3360 bool isSatisfy = true;
3361 for ( ; aNodeItr->more(); )
3363 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3364 if ( !isOnSurface( aNode ) )
3371 myIds.Add( theElemPtr->GetID() );
3374 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
3376 if ( mySurf.IsNull() )
3379 gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
3380 // double aToler2 = myToler * myToler;
3381 // if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
3383 // gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
3384 // if ( aPln.SquareDistance( aPnt ) > aToler2 )
3387 // else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
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 )
3400 myProjector.Perform( aPnt );
3401 bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
3411 ElementsOnShape::ElementsOnShape()
3413 myType(SMDSAbs_All),
3414 myToler(Precision::Confusion()),
3415 myAllNodesFlag(false)
3417 myCurShapeType = TopAbs_SHAPE;
3420 ElementsOnShape::~ElementsOnShape()
3424 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
3426 if (myMesh != theMesh) {
3428 SetShape(myShape, myType);
3432 bool ElementsOnShape::IsSatisfy (long theElementId)
3434 return myIds.Contains(theElementId);
3437 SMDSAbs_ElementType ElementsOnShape::GetType() const
3442 void ElementsOnShape::SetTolerance (const double theToler)
3444 if (myToler != theToler) {
3446 SetShape(myShape, myType);
3450 double ElementsOnShape::GetTolerance() const
3455 void ElementsOnShape::SetAllNodes (bool theAllNodes)
3457 if (myAllNodesFlag != theAllNodes) {
3458 myAllNodesFlag = theAllNodes;
3459 SetShape(myShape, myType);
3463 void ElementsOnShape::SetShape (const TopoDS_Shape& theShape,
3464 const SMDSAbs_ElementType theType)
3470 if (myMesh == 0) return;
3475 myIds.ReSize(myMesh->NbEdges() + myMesh->NbFaces() + myMesh->NbVolumes());
3478 myIds.ReSize(myMesh->NbNodes());
3481 myIds.ReSize(myMesh->NbEdges());
3484 myIds.ReSize(myMesh->NbFaces());
3486 case SMDSAbs_Volume:
3487 myIds.ReSize(myMesh->NbVolumes());
3493 myShapesMap.Clear();
3497 void ElementsOnShape::addShape (const TopoDS_Shape& theShape)
3499 if (theShape.IsNull() || myMesh == 0)
3502 if (!myShapesMap.Add(theShape)) return;
3504 myCurShapeType = theShape.ShapeType();
3505 switch (myCurShapeType)
3507 case TopAbs_COMPOUND:
3508 case TopAbs_COMPSOLID:
3512 TopoDS_Iterator anIt (theShape, Standard_True, Standard_True);
3513 for (; anIt.More(); anIt.Next()) addShape(anIt.Value());
3518 myCurSC.Load(theShape);
3524 TopoDS_Face aFace = TopoDS::Face(theShape);
3525 BRepAdaptor_Surface SA (aFace, true);
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);
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);
3548 TopoDS_Vertex aV = TopoDS::Vertex(theShape);
3549 myCurPnt = BRep_Tool::Pnt(aV);
3558 void ElementsOnShape::process()
3560 if (myShape.IsNull() || myMesh == 0)
3563 if (myType == SMDSAbs_Node)
3565 SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
3566 while (anIter->more())
3567 process(anIter->next());
3571 if (myType == SMDSAbs_Edge || myType == SMDSAbs_All)
3573 SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
3574 while (anIter->more())
3575 process(anIter->next());
3578 if (myType == SMDSAbs_Face || myType == SMDSAbs_All)
3580 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
3581 while (anIter->more()) {
3582 process(anIter->next());
3586 if (myType == SMDSAbs_Volume || myType == SMDSAbs_All)
3588 SMDS_VolumeIteratorPtr anIter = myMesh->volumesIterator();
3589 while (anIter->more())
3590 process(anIter->next());
3595 void ElementsOnShape::process (const SMDS_MeshElement* theElemPtr)
3597 if (myShape.IsNull())
3600 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3601 bool isSatisfy = myAllNodesFlag;
3603 gp_XYZ centerXYZ (0, 0, 0);
3605 while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
3607 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3608 gp_Pnt aPnt (aNode->X(), aNode->Y(), aNode->Z());
3609 centerXYZ += aPnt.XYZ();
3611 switch (myCurShapeType)
3615 myCurSC.Perform(aPnt, myToler);
3616 isSatisfy = (myCurSC.State() == TopAbs_IN || myCurSC.State() == TopAbs_ON);
3621 myCurProjFace.Perform(aPnt);
3622 isSatisfy = (myCurProjFace.IsDone() && myCurProjFace.LowerDistance() <= myToler);
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);
3636 myCurProjEdge.Perform(aPnt);
3637 isSatisfy = (myCurProjEdge.NbPoints() > 0 && myCurProjEdge.LowerDistance() <= myToler);
3642 isSatisfy = (aPnt.Distance(myCurPnt) <= myToler);
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))
3661 myIds.Add(theElemPtr->GetID());
3664 TSequenceOfXYZ::TSequenceOfXYZ()
3667 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n)
3670 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t)
3673 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray)
3676 template <class InputIterator>
3677 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd)
3680 TSequenceOfXYZ::~TSequenceOfXYZ()
3683 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
3685 myArray = theSequenceOfXYZ.myArray;
3689 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
3691 return myArray[n-1];
3694 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
3696 return myArray[n-1];
3699 void TSequenceOfXYZ::clear()
3704 void TSequenceOfXYZ::reserve(size_type n)
3709 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
3711 myArray.push_back(v);
3714 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
3716 return myArray.size();