1 // Copyright (C) 2007-2011 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"
69 #include <vtkMeshQuality.h>
77 inline gp_XYZ gpXYZ(const SMDS_MeshNode* aNode )
79 return gp_XYZ(aNode->X(), aNode->Y(), aNode->Z() );
82 inline double getAngle( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
84 gp_Vec v1( P1 - P2 ), v2( P3 - P2 );
86 return v1.Magnitude() < gp::Resolution() ||
87 v2.Magnitude() < gp::Resolution() ? 0 : v1.Angle( v2 );
90 inline double getArea( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
92 gp_Vec aVec1( P2 - P1 );
93 gp_Vec aVec2( P3 - P1 );
94 return ( aVec1 ^ aVec2 ).Magnitude() * 0.5;
97 inline double getArea( const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3 )
99 return getArea( P1.XYZ(), P2.XYZ(), P3.XYZ() );
104 inline double getDistance( const gp_XYZ& P1, const gp_XYZ& P2 )
106 double aDist = gp_Pnt( P1 ).Distance( gp_Pnt( P2 ) );
110 int getNbMultiConnection( const SMDS_Mesh* theMesh, const int theId )
115 const SMDS_MeshElement* anEdge = theMesh->FindElement( theId );
116 if ( anEdge == 0 || anEdge->GetType() != SMDSAbs_Edge/* || anEdge->NbNodes() != 2 */)
119 // for each pair of nodes in anEdge (there are 2 pairs in a quadratic edge)
120 // count elements containing both nodes of the pair.
121 // Note that there may be such cases for a quadratic edge (a horizontal line):
126 // +-----+------+ +-----+------+
129 // result sould be 2 in both cases
131 int aResult0 = 0, aResult1 = 0;
132 // last node, it is a medium one in a quadratic edge
133 const SMDS_MeshNode* aLastNode = anEdge->GetNode( anEdge->NbNodes() - 1 );
134 const SMDS_MeshNode* aNode0 = anEdge->GetNode( 0 );
135 const SMDS_MeshNode* aNode1 = anEdge->GetNode( 1 );
136 if ( aNode1 == aLastNode ) aNode1 = 0;
138 SMDS_ElemIteratorPtr anElemIter = aLastNode->GetInverseElementIterator();
139 while( anElemIter->more() ) {
140 const SMDS_MeshElement* anElem = anElemIter->next();
141 if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
142 SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
143 while ( anIter->more() ) {
144 if ( const SMDS_MeshElement* anElemNode = anIter->next() ) {
145 if ( anElemNode == aNode0 ) {
147 if ( !aNode1 ) break; // not a quadratic edge
149 else if ( anElemNode == aNode1 )
155 int aResult = std::max ( aResult0, aResult1 );
157 // TColStd_MapOfInteger aMap;
159 // SMDS_ElemIteratorPtr anIter = anEdge->nodesIterator();
160 // if ( anIter != 0 ) {
161 // while( anIter->more() ) {
162 // const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
165 // SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
166 // while( anElemIter->more() ) {
167 // const SMDS_MeshElement* anElem = anElemIter->next();
168 // if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
169 // int anId = anElem->GetID();
171 // if ( anIter->more() ) // i.e. first node
173 // else if ( aMap.Contains( anId ) )
183 gp_XYZ getNormale( const SMDS_MeshFace* theFace, bool* ok=0 )
185 int aNbNode = theFace->NbNodes();
187 gp_XYZ q1 = gpXYZ( theFace->GetNode(1)) - gpXYZ( theFace->GetNode(0));
188 gp_XYZ q2 = gpXYZ( theFace->GetNode(2)) - gpXYZ( theFace->GetNode(0));
191 gp_XYZ q3 = gpXYZ( theFace->GetNode(3)) - gpXYZ( theFace->GetNode(0));
194 double len = n.Modulus();
195 bool zeroLen = ( len <= numeric_limits<double>::min());
199 if (ok) *ok = !zeroLen;
207 using namespace SMESH::Controls;
214 Class : NumericalFunctor
215 Description : Base class for numerical functors
217 NumericalFunctor::NumericalFunctor():
223 void NumericalFunctor::SetMesh( const SMDS_Mesh* theMesh )
228 bool NumericalFunctor::GetPoints(const int theId,
229 TSequenceOfXYZ& theRes ) const
236 return GetPoints( myMesh->FindElement( theId ), theRes );
239 bool NumericalFunctor::GetPoints(const SMDS_MeshElement* anElem,
240 TSequenceOfXYZ& theRes )
247 theRes.reserve( anElem->NbNodes() );
249 // Get nodes of the element
250 SMDS_ElemIteratorPtr anIter;
252 if ( anElem->IsQuadratic() ) {
253 switch ( anElem->GetType() ) {
255 anIter = dynamic_cast<const SMDS_VtkEdge*>
256 (anElem)->interlacedNodesElemIterator();
259 anIter = dynamic_cast<const SMDS_VtkFace*>
260 (anElem)->interlacedNodesElemIterator();
263 anIter = anElem->nodesIterator();
268 anIter = anElem->nodesIterator();
272 while( anIter->more() ) {
273 if ( const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>( anIter->next() ))
274 theRes.push_back( gp_XYZ( aNode->X(), aNode->Y(), aNode->Z() ) );
281 long NumericalFunctor::GetPrecision() const
286 void NumericalFunctor::SetPrecision( const long thePrecision )
288 myPrecision = thePrecision;
289 myPrecisionValue = pow( 10., (double)( myPrecision ) );
292 double NumericalFunctor::GetValue( long theId )
296 myCurrElement = myMesh->FindElement( theId );
299 if ( GetPoints( theId, P ))
300 aVal = Round( GetValue( P ));
305 double NumericalFunctor::Round( const double & aVal )
307 return ( myPrecision >= 0 ) ? floor( aVal * myPrecisionValue + 0.5 ) / myPrecisionValue : aVal;
310 //================================================================================
312 * \brief Return histogram of functor values
313 * \param nbIntervals - number of intervals
314 * \param nbEvents - number of mesh elements having values within i-th interval
315 * \param funValues - boundaries of intervals
316 * \param elements - elements to check vulue of; empty list means "of all"
317 * \param minmax - boundaries of diapason of values to divide into intervals
319 //================================================================================
321 void NumericalFunctor::GetHistogram(int nbIntervals,
322 std::vector<int>& nbEvents,
323 std::vector<double>& funValues,
324 const vector<int>& elements,
325 const double* minmax)
327 if ( nbIntervals < 1 ||
329 !myMesh->GetMeshInfo().NbElements( GetType() ))
331 nbEvents.resize( nbIntervals, 0 );
332 funValues.resize( nbIntervals+1 );
334 // get all values sorted
335 std::multiset< double > values;
336 if ( elements.empty() )
338 SMDS_ElemIteratorPtr elemIt = myMesh->elementsIterator(GetType());
339 while ( elemIt->more() )
340 values.insert( GetValue( elemIt->next()->GetID() ));
344 vector<int>::const_iterator id = elements.begin();
345 for ( ; id != elements.end(); ++id )
346 values.insert( GetValue( *id ));
351 funValues[0] = minmax[0];
352 funValues[nbIntervals] = minmax[1];
356 funValues[0] = *values.begin();
357 funValues[nbIntervals] = *values.rbegin();
359 // case nbIntervals == 1
360 if ( nbIntervals == 1 )
362 nbEvents[0] = values.size();
366 if (funValues.front() == funValues.back())
368 nbEvents.resize( 1 );
369 nbEvents[0] = values.size();
370 funValues[1] = funValues.back();
371 funValues.resize( 2 );
374 std::multiset< double >::iterator min = values.begin(), max;
375 for ( int i = 0; i < nbIntervals; ++i )
377 // find end value of i-th interval
378 double r = (i+1) / double( nbIntervals );
379 funValues[i+1] = funValues.front() * (1-r) + funValues.back() * r;
381 // count values in the i-th interval if there are any
382 if ( min != values.end() && *min <= funValues[i+1] )
384 // find the first value out of the interval
385 max = values.upper_bound( funValues[i+1] ); // max is greater than funValues[i+1], or end()
386 nbEvents[i] = std::distance( min, max );
390 // add values larger than minmax[1]
391 nbEvents.back() += std::distance( min, values.end() );
394 //=======================================================================
395 //function : GetValue
397 //=======================================================================
399 double Volume::GetValue( long theElementId )
401 if ( theElementId && myMesh ) {
402 SMDS_VolumeTool aVolumeTool;
403 if ( aVolumeTool.Set( myMesh->FindElement( theElementId )))
404 return aVolumeTool.GetSize();
409 //=======================================================================
410 //function : GetBadRate
411 //purpose : meaningless as it is not quality control functor
412 //=======================================================================
414 double Volume::GetBadRate( double Value, int /*nbNodes*/ ) const
419 //=======================================================================
422 //=======================================================================
424 SMDSAbs_ElementType Volume::GetType() const
426 return SMDSAbs_Volume;
431 Class : MaxElementLength2D
432 Description : Functor calculating maximum length of 2D element
435 double MaxElementLength2D::GetValue( long theElementId )
438 if( GetPoints( theElementId, P ) ) {
440 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
441 SMDSAbs_ElementType aType = aElem->GetType();
445 if( len == 3 ) { // triangles
446 double L1 = getDistance(P( 1 ),P( 2 ));
447 double L2 = getDistance(P( 2 ),P( 3 ));
448 double L3 = getDistance(P( 3 ),P( 1 ));
449 aVal = Max(L1,Max(L2,L3));
452 else if( len == 4 ) { // quadrangles
453 double L1 = getDistance(P( 1 ),P( 2 ));
454 double L2 = getDistance(P( 2 ),P( 3 ));
455 double L3 = getDistance(P( 3 ),P( 4 ));
456 double L4 = getDistance(P( 4 ),P( 1 ));
457 double D1 = getDistance(P( 1 ),P( 3 ));
458 double D2 = getDistance(P( 2 ),P( 4 ));
459 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
462 else if( len == 6 ) { // quadratic triangles
463 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
464 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
465 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
466 aVal = Max(L1,Max(L2,L3));
469 else if( len == 8 || len == 9 ) { // quadratic quadrangles
470 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
471 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
472 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
473 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
474 double D1 = getDistance(P( 1 ),P( 5 ));
475 double D2 = getDistance(P( 3 ),P( 7 ));
476 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
481 if( myPrecision >= 0 )
483 double prec = pow( 10., (double)myPrecision );
484 aVal = floor( aVal * prec + 0.5 ) / prec;
491 double MaxElementLength2D::GetBadRate( double Value, int /*nbNodes*/ ) const
496 SMDSAbs_ElementType MaxElementLength2D::GetType() const
502 Class : MaxElementLength3D
503 Description : Functor calculating maximum length of 3D element
506 double MaxElementLength3D::GetValue( long theElementId )
509 if( GetPoints( theElementId, P ) ) {
511 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
512 SMDSAbs_ElementType aType = aElem->GetType();
516 if( len == 4 ) { // tetras
517 double L1 = getDistance(P( 1 ),P( 2 ));
518 double L2 = getDistance(P( 2 ),P( 3 ));
519 double L3 = getDistance(P( 3 ),P( 1 ));
520 double L4 = getDistance(P( 1 ),P( 4 ));
521 double L5 = getDistance(P( 2 ),P( 4 ));
522 double L6 = getDistance(P( 3 ),P( 4 ));
523 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
526 else if( len == 5 ) { // pyramids
527 double L1 = getDistance(P( 1 ),P( 2 ));
528 double L2 = getDistance(P( 2 ),P( 3 ));
529 double L3 = getDistance(P( 3 ),P( 4 ));
530 double L4 = getDistance(P( 4 ),P( 1 ));
531 double L5 = getDistance(P( 1 ),P( 5 ));
532 double L6 = getDistance(P( 2 ),P( 5 ));
533 double L7 = getDistance(P( 3 ),P( 5 ));
534 double L8 = getDistance(P( 4 ),P( 5 ));
535 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
536 aVal = Max(aVal,Max(L7,L8));
539 else if( len == 6 ) { // pentas
540 double L1 = getDistance(P( 1 ),P( 2 ));
541 double L2 = getDistance(P( 2 ),P( 3 ));
542 double L3 = getDistance(P( 3 ),P( 1 ));
543 double L4 = getDistance(P( 4 ),P( 5 ));
544 double L5 = getDistance(P( 5 ),P( 6 ));
545 double L6 = getDistance(P( 6 ),P( 4 ));
546 double L7 = getDistance(P( 1 ),P( 4 ));
547 double L8 = getDistance(P( 2 ),P( 5 ));
548 double L9 = getDistance(P( 3 ),P( 6 ));
549 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
550 aVal = Max(aVal,Max(Max(L7,L8),L9));
553 else if( len == 8 ) { // hexas
554 double L1 = getDistance(P( 1 ),P( 2 ));
555 double L2 = getDistance(P( 2 ),P( 3 ));
556 double L3 = getDistance(P( 3 ),P( 4 ));
557 double L4 = getDistance(P( 4 ),P( 1 ));
558 double L5 = getDistance(P( 5 ),P( 6 ));
559 double L6 = getDistance(P( 6 ),P( 7 ));
560 double L7 = getDistance(P( 7 ),P( 8 ));
561 double L8 = getDistance(P( 8 ),P( 5 ));
562 double L9 = getDistance(P( 1 ),P( 5 ));
563 double L10= getDistance(P( 2 ),P( 6 ));
564 double L11= getDistance(P( 3 ),P( 7 ));
565 double L12= getDistance(P( 4 ),P( 8 ));
566 double D1 = getDistance(P( 1 ),P( 7 ));
567 double D2 = getDistance(P( 2 ),P( 8 ));
568 double D3 = getDistance(P( 3 ),P( 5 ));
569 double D4 = getDistance(P( 4 ),P( 6 ));
570 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
571 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
572 aVal = Max(aVal,Max(L11,L12));
573 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
576 else if( len == 12 ) { // hexagonal prism
577 for ( int i1 = 1; i1 < 12; ++i1 )
578 for ( int i2 = i1+1; i1 <= 12; ++i1 )
579 aVal = Max( aVal, getDistance(P( i1 ),P( i2 )));
582 else if( len == 10 ) { // quadratic tetras
583 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
584 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
585 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
586 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
587 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
588 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
589 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
592 else if( len == 13 ) { // quadratic pyramids
593 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
594 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
595 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
596 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
597 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
598 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
599 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
600 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
601 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
602 aVal = Max(aVal,Max(L7,L8));
605 else if( len == 15 ) { // quadratic pentas
606 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
607 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
608 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
609 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
610 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
611 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
612 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
613 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
614 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
615 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
616 aVal = Max(aVal,Max(Max(L7,L8),L9));
619 else if( len == 20 || len == 27 ) { // quadratic hexas
620 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
621 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
622 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
623 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
624 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
625 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
626 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
627 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
628 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
629 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
630 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
631 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
632 double D1 = getDistance(P( 1 ),P( 7 ));
633 double D2 = getDistance(P( 2 ),P( 8 ));
634 double D3 = getDistance(P( 3 ),P( 5 ));
635 double D4 = getDistance(P( 4 ),P( 6 ));
636 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
637 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
638 aVal = Max(aVal,Max(L11,L12));
639 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
642 else if( len > 1 && aElem->IsPoly() ) { // polys
643 // get the maximum distance between all pairs of nodes
644 for( int i = 1; i <= len; i++ ) {
645 for( int j = 1; j <= len; j++ ) {
646 if( j > i ) { // optimization of the loop
647 double D = getDistance( P(i), P(j) );
648 aVal = Max( aVal, D );
655 if( myPrecision >= 0 )
657 double prec = pow( 10., (double)myPrecision );
658 aVal = floor( aVal * prec + 0.5 ) / prec;
665 double MaxElementLength3D::GetBadRate( double Value, int /*nbNodes*/ ) const
670 SMDSAbs_ElementType MaxElementLength3D::GetType() const
672 return SMDSAbs_Volume;
678 Description : Functor for calculation of minimum angle
681 double MinimumAngle::GetValue( const TSequenceOfXYZ& P )
688 aMin = getAngle(P( P.size() ), P( 1 ), P( 2 ));
689 aMin = Min(aMin,getAngle(P( P.size()-1 ), P( P.size() ), P( 1 )));
691 for (int i=2; i<P.size();i++){
692 double A0 = getAngle( P( i-1 ), P( i ), P( i+1 ) );
696 return aMin * 180.0 / PI;
699 double MinimumAngle::GetBadRate( double Value, int nbNodes ) const
701 //const double aBestAngle = PI / nbNodes;
702 const double aBestAngle = 180.0 - ( 360.0 / double(nbNodes) );
703 return ( fabs( aBestAngle - Value ));
706 SMDSAbs_ElementType MinimumAngle::GetType() const
714 Description : Functor for calculating aspect ratio
716 double AspectRatio::GetValue( const TSequenceOfXYZ& P )
718 // According to "Mesh quality control" by Nadir Bouhamau referring to
719 // Pascal Jean Frey and Paul-Louis George. Maillages, applications aux elements finis.
720 // Hermes Science publications, Paris 1999 ISBN 2-7462-0024-4
723 int nbNodes = P.size();
728 // Compute aspect ratio
730 if ( nbNodes == 3 ) {
731 // Compute lengths of the sides
732 std::vector< double > aLen (nbNodes);
733 for ( int i = 0; i < nbNodes - 1; i++ )
734 aLen[ i ] = getDistance( P( i + 1 ), P( i + 2 ) );
735 aLen[ nbNodes - 1 ] = getDistance( P( 1 ), P( nbNodes ) );
736 // Q = alfa * h * p / S, where
738 // alfa = sqrt( 3 ) / 6
739 // h - length of the longest edge
740 // p - half perimeter
741 // S - triangle surface
742 const double alfa = sqrt( 3. ) / 6.;
743 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
744 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
745 double anArea = getArea( P( 1 ), P( 2 ), P( 3 ) );
746 if ( anArea <= Precision::Confusion() )
748 return alfa * maxLen * half_perimeter / anArea;
750 else if ( nbNodes == 6 ) { // quadratic triangles
751 // Compute lengths of the sides
752 std::vector< double > aLen (3);
753 aLen[0] = getDistance( P(1), P(3) );
754 aLen[1] = getDistance( P(3), P(5) );
755 aLen[2] = getDistance( P(5), P(1) );
756 // Q = alfa * h * p / S, where
758 // alfa = sqrt( 3 ) / 6
759 // h - length of the longest edge
760 // p - half perimeter
761 // S - triangle surface
762 const double alfa = sqrt( 3. ) / 6.;
763 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
764 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
765 double anArea = getArea( P(1), P(3), P(5) );
766 if ( anArea <= Precision::Confusion() )
768 return alfa * maxLen * half_perimeter / anArea;
770 else if( nbNodes == 4 ) { // quadrangle
771 // Compute lengths of the sides
772 std::vector< double > aLen (4);
773 aLen[0] = getDistance( P(1), P(2) );
774 aLen[1] = getDistance( P(2), P(3) );
775 aLen[2] = getDistance( P(3), P(4) );
776 aLen[3] = getDistance( P(4), P(1) );
777 // Compute lengths of the diagonals
778 std::vector< double > aDia (2);
779 aDia[0] = getDistance( P(1), P(3) );
780 aDia[1] = getDistance( P(2), P(4) );
781 // Compute areas of all triangles which can be built
782 // taking three nodes of the quadrangle
783 std::vector< double > anArea (4);
784 anArea[0] = getArea( P(1), P(2), P(3) );
785 anArea[1] = getArea( P(1), P(2), P(4) );
786 anArea[2] = getArea( P(1), P(3), P(4) );
787 anArea[3] = getArea( P(2), P(3), P(4) );
788 // Q = alpha * L * C1 / C2, where
790 // alpha = sqrt( 1/32 )
791 // L = max( L1, L2, L3, L4, D1, D2 )
792 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
793 // C2 = min( S1, S2, S3, S4 )
794 // Li - lengths of the edges
795 // Di - lengths of the diagonals
796 // Si - areas of the triangles
797 const double alpha = sqrt( 1 / 32. );
798 double L = Max( aLen[ 0 ],
802 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
803 double C1 = sqrt( ( aLen[0] * aLen[0] +
806 aLen[3] * aLen[3] ) / 4. );
807 double C2 = Min( anArea[ 0 ],
809 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
810 if ( C2 <= Precision::Confusion() )
812 return alpha * L * C1 / C2;
814 else if( nbNodes == 8 || nbNodes == 9 ) { // nbNodes==8 - quadratic quadrangle
815 // Compute lengths of the sides
816 std::vector< double > aLen (4);
817 aLen[0] = getDistance( P(1), P(3) );
818 aLen[1] = getDistance( P(3), P(5) );
819 aLen[2] = getDistance( P(5), P(7) );
820 aLen[3] = getDistance( P(7), P(1) );
821 // Compute lengths of the diagonals
822 std::vector< double > aDia (2);
823 aDia[0] = getDistance( P(1), P(5) );
824 aDia[1] = getDistance( P(3), P(7) );
825 // Compute areas of all triangles which can be built
826 // taking three nodes of the quadrangle
827 std::vector< double > anArea (4);
828 anArea[0] = getArea( P(1), P(3), P(5) );
829 anArea[1] = getArea( P(1), P(3), P(7) );
830 anArea[2] = getArea( P(1), P(5), P(7) );
831 anArea[3] = getArea( P(3), P(5), P(7) );
832 // Q = alpha * L * C1 / C2, where
834 // alpha = sqrt( 1/32 )
835 // L = max( L1, L2, L3, L4, D1, D2 )
836 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
837 // C2 = min( S1, S2, S3, S4 )
838 // Li - lengths of the edges
839 // Di - lengths of the diagonals
840 // Si - areas of the triangles
841 const double alpha = sqrt( 1 / 32. );
842 double L = Max( aLen[ 0 ],
846 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
847 double C1 = sqrt( ( aLen[0] * aLen[0] +
850 aLen[3] * aLen[3] ) / 4. );
851 double C2 = Min( anArea[ 0 ],
853 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
854 if ( C2 <= Precision::Confusion() )
856 return alpha * L * C1 / C2;
861 double AspectRatio::GetBadRate( double Value, int /*nbNodes*/ ) const
863 // the aspect ratio is in the range [1.0,infinity]
864 // < 1.0 = very bad, zero area
867 return ( Value < 0.9 ) ? 1000 : Value / 1000.;
870 SMDSAbs_ElementType AspectRatio::GetType() const
877 Class : AspectRatio3D
878 Description : Functor for calculating aspect ratio
882 inline double getHalfPerimeter(double theTria[3]){
883 return (theTria[0] + theTria[1] + theTria[2])/2.0;
886 inline double getArea(double theHalfPerim, double theTria[3]){
887 return sqrt(theHalfPerim*
888 (theHalfPerim-theTria[0])*
889 (theHalfPerim-theTria[1])*
890 (theHalfPerim-theTria[2]));
893 inline double getVolume(double theLen[6]){
894 double a2 = theLen[0]*theLen[0];
895 double b2 = theLen[1]*theLen[1];
896 double c2 = theLen[2]*theLen[2];
897 double d2 = theLen[3]*theLen[3];
898 double e2 = theLen[4]*theLen[4];
899 double f2 = theLen[5]*theLen[5];
900 double P = 4.0*a2*b2*d2;
901 double Q = a2*(b2+d2-e2)-b2*(a2+d2-f2)-d2*(a2+b2-c2);
902 double R = (b2+d2-e2)*(a2+d2-f2)*(a2+d2-f2);
903 return sqrt(P-Q+R)/12.0;
906 inline double getVolume2(double theLen[6]){
907 double a2 = theLen[0]*theLen[0];
908 double b2 = theLen[1]*theLen[1];
909 double c2 = theLen[2]*theLen[2];
910 double d2 = theLen[3]*theLen[3];
911 double e2 = theLen[4]*theLen[4];
912 double f2 = theLen[5]*theLen[5];
914 double P = a2*e2*(b2+c2+d2+f2-a2-e2);
915 double Q = b2*f2*(a2+c2+d2+e2-b2-f2);
916 double R = c2*d2*(a2+b2+e2+f2-c2-d2);
917 double S = a2*b2*d2+b2*c2*e2+a2*c2*f2+d2*e2*f2;
919 return sqrt(P+Q+R-S)/12.0;
922 inline double getVolume(const TSequenceOfXYZ& P){
923 gp_Vec aVec1( P( 2 ) - P( 1 ) );
924 gp_Vec aVec2( P( 3 ) - P( 1 ) );
925 gp_Vec aVec3( P( 4 ) - P( 1 ) );
926 gp_Vec anAreaVec( aVec1 ^ aVec2 );
927 return fabs(aVec3 * anAreaVec) / 6.0;
930 inline double getMaxHeight(double theLen[6])
932 double aHeight = std::max(theLen[0],theLen[1]);
933 aHeight = std::max(aHeight,theLen[2]);
934 aHeight = std::max(aHeight,theLen[3]);
935 aHeight = std::max(aHeight,theLen[4]);
936 aHeight = std::max(aHeight,theLen[5]);
942 double AspectRatio3D::GetValue( long theId )
945 myCurrElement = myMesh->FindElement( theId );
946 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_TETRA )
948 // Action from CoTech | ACTION 31.3:
949 // EURIWARE BO: Homogenize the formulas used to calculate the Controls in SMESH to fit with
950 // those of ParaView. The library used by ParaView for those calculations can be reused in SMESH.
951 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
952 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
953 aVal = Round( vtkMeshQuality::TetAspectRatio( avtkCell ));
958 if ( GetPoints( myCurrElement, P ))
959 aVal = Round( GetValue( P ));
964 double AspectRatio3D::GetValue( const TSequenceOfXYZ& P )
966 double aQuality = 0.0;
967 if(myCurrElement->IsPoly()) return aQuality;
969 int nbNodes = P.size();
971 if(myCurrElement->IsQuadratic()) {
972 if(nbNodes==10) nbNodes=4; // quadratic tetrahedron
973 else if(nbNodes==13) nbNodes=5; // quadratic pyramid
974 else if(nbNodes==15) nbNodes=6; // quadratic pentahedron
975 else if(nbNodes==20) nbNodes=8; // quadratic hexahedron
976 else if(nbNodes==27) nbNodes=8; // quadratic hexahedron
977 else return aQuality;
983 getDistance(P( 1 ),P( 2 )), // a
984 getDistance(P( 2 ),P( 3 )), // b
985 getDistance(P( 3 ),P( 1 )), // c
986 getDistance(P( 2 ),P( 4 )), // d
987 getDistance(P( 3 ),P( 4 )), // e
988 getDistance(P( 1 ),P( 4 )) // f
990 double aTria[4][3] = {
991 {aLen[0],aLen[1],aLen[2]}, // abc
992 {aLen[0],aLen[3],aLen[5]}, // adf
993 {aLen[1],aLen[3],aLen[4]}, // bde
994 {aLen[2],aLen[4],aLen[5]} // cef
996 double aSumArea = 0.0;
997 double aHalfPerimeter = getHalfPerimeter(aTria[0]);
998 double anArea = getArea(aHalfPerimeter,aTria[0]);
1000 aHalfPerimeter = getHalfPerimeter(aTria[1]);
1001 anArea = getArea(aHalfPerimeter,aTria[1]);
1003 aHalfPerimeter = getHalfPerimeter(aTria[2]);
1004 anArea = getArea(aHalfPerimeter,aTria[2]);
1006 aHalfPerimeter = getHalfPerimeter(aTria[3]);
1007 anArea = getArea(aHalfPerimeter,aTria[3]);
1009 double aVolume = getVolume(P);
1010 //double aVolume = getVolume(aLen);
1011 double aHeight = getMaxHeight(aLen);
1012 static double aCoeff = sqrt(2.0)/12.0;
1013 if ( aVolume > DBL_MIN )
1014 aQuality = aCoeff*aHeight*aSumArea/aVolume;
1019 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 3 ),P( 5 )};
1020 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1023 gp_XYZ aXYZ[4] = {P( 1 ),P( 3 ),P( 4 ),P( 5 )};
1024 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1027 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 5 )};
1028 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1031 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 4 ),P( 5 )};
1032 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1038 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 6 )};
1039 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1042 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 3 )};
1043 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1046 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 6 )};
1047 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1050 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1051 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1054 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 6 )};
1055 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1058 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 3 )};
1059 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1065 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1066 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1069 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 4 )};
1070 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1073 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 7 )};
1074 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1077 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 8 )};
1078 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1081 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 3 )};
1082 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1085 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 4 )};
1086 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1089 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 7 )};
1090 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1093 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 8 )};
1094 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1097 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 3 )};
1098 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1101 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 4 )};
1102 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1105 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 7 )};
1106 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1109 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 8 )};
1110 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1113 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 1 )};
1114 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1117 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 2 )};
1118 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1121 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 5 )};
1122 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1125 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 6 )};
1126 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1129 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 1 )};
1130 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1133 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 2 )};
1134 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1137 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 5 )};
1138 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1141 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 6 )};
1142 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1145 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 1 )};
1146 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1149 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1150 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1153 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 5 )};
1154 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1157 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 6 )};
1158 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1161 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1162 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1165 gp_XYZ aXYZ[4] = {P( 4 ),P( 5 ),P( 8 ),P( 2 )};
1166 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1169 gp_XYZ aXYZ[4] = {P( 1 ),P( 4 ),P( 5 ),P( 3 )};
1170 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1173 gp_XYZ aXYZ[4] = {P( 3 ),P( 6 ),P( 7 ),P( 1 )};
1174 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1177 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 6 ),P( 4 )};
1178 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1181 gp_XYZ aXYZ[4] = {P( 5 ),P( 6 ),P( 8 ),P( 3 )};
1182 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1185 gp_XYZ aXYZ[4] = {P( 7 ),P( 8 ),P( 6 ),P( 1 )};
1186 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1189 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 7 )};
1190 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1193 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 2 ),P( 5 )};
1194 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1200 gp_XYZ aXYZ[8] = {P( 1 ),P( 2 ),P( 4 ),P( 5 ),P( 7 ),P( 8 ),P( 10 ),P( 11 )};
1201 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1204 gp_XYZ aXYZ[8] = {P( 2 ),P( 3 ),P( 5 ),P( 6 ),P( 8 ),P( 9 ),P( 11 ),P( 12 )};
1205 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1208 gp_XYZ aXYZ[8] = {P( 3 ),P( 4 ),P( 6 ),P( 1 ),P( 9 ),P( 10 ),P( 12 ),P( 7 )};
1209 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1212 } // switch(nbNodes)
1214 if ( nbNodes > 4 ) {
1215 // avaluate aspect ratio of quadranle faces
1216 AspectRatio aspect2D;
1217 SMDS_VolumeTool::VolumeType type = SMDS_VolumeTool::GetType( nbNodes );
1218 int nbFaces = SMDS_VolumeTool::NbFaces( type );
1219 TSequenceOfXYZ points(4);
1220 for ( int i = 0; i < nbFaces; ++i ) { // loop on faces of a volume
1221 if ( SMDS_VolumeTool::NbFaceNodes( type, i ) != 4 )
1223 const int* pInd = SMDS_VolumeTool::GetFaceNodesIndices( type, i, true );
1224 for ( int p = 0; p < 4; ++p ) // loop on nodes of a quadranle face
1225 points( p + 1 ) = P( pInd[ p ] + 1 );
1226 aQuality = std::max( aQuality, aspect2D.GetValue( points ));
1232 double AspectRatio3D::GetBadRate( double Value, int /*nbNodes*/ ) const
1234 // the aspect ratio is in the range [1.0,infinity]
1237 return Value / 1000.;
1240 SMDSAbs_ElementType AspectRatio3D::GetType() const
1242 return SMDSAbs_Volume;
1248 Description : Functor for calculating warping
1250 double Warping::GetValue( const TSequenceOfXYZ& P )
1252 if ( P.size() != 4 )
1255 gp_XYZ G = ( P( 1 ) + P( 2 ) + P( 3 ) + P( 4 ) ) / 4.;
1257 double A1 = ComputeA( P( 1 ), P( 2 ), P( 3 ), G );
1258 double A2 = ComputeA( P( 2 ), P( 3 ), P( 4 ), G );
1259 double A3 = ComputeA( P( 3 ), P( 4 ), P( 1 ), G );
1260 double A4 = ComputeA( P( 4 ), P( 1 ), P( 2 ), G );
1262 return Max( Max( A1, A2 ), Max( A3, A4 ) );
1265 double Warping::ComputeA( const gp_XYZ& thePnt1,
1266 const gp_XYZ& thePnt2,
1267 const gp_XYZ& thePnt3,
1268 const gp_XYZ& theG ) const
1270 double aLen1 = gp_Pnt( thePnt1 ).Distance( gp_Pnt( thePnt2 ) );
1271 double aLen2 = gp_Pnt( thePnt2 ).Distance( gp_Pnt( thePnt3 ) );
1272 double L = Min( aLen1, aLen2 ) * 0.5;
1273 if ( L < Precision::Confusion())
1276 gp_XYZ GI = ( thePnt2 + thePnt1 ) / 2. - theG;
1277 gp_XYZ GJ = ( thePnt3 + thePnt2 ) / 2. - theG;
1278 gp_XYZ N = GI.Crossed( GJ );
1280 if ( N.Modulus() < gp::Resolution() )
1285 double H = ( thePnt2 - theG ).Dot( N );
1286 return asin( fabs( H / L ) ) * 180. / PI;
1289 double Warping::GetBadRate( double Value, int /*nbNodes*/ ) const
1291 // the warp is in the range [0.0,PI/2]
1292 // 0.0 = good (no warp)
1293 // PI/2 = bad (face pliee)
1297 SMDSAbs_ElementType Warping::GetType() const
1299 return SMDSAbs_Face;
1305 Description : Functor for calculating taper
1307 double Taper::GetValue( const TSequenceOfXYZ& P )
1309 if ( P.size() != 4 )
1313 double J1 = getArea( P( 4 ), P( 1 ), P( 2 ) ) / 2.;
1314 double J2 = getArea( P( 3 ), P( 1 ), P( 2 ) ) / 2.;
1315 double J3 = getArea( P( 2 ), P( 3 ), P( 4 ) ) / 2.;
1316 double J4 = getArea( P( 3 ), P( 4 ), P( 1 ) ) / 2.;
1318 double JA = 0.25 * ( J1 + J2 + J3 + J4 );
1319 if ( JA <= Precision::Confusion() )
1322 double T1 = fabs( ( J1 - JA ) / JA );
1323 double T2 = fabs( ( J2 - JA ) / JA );
1324 double T3 = fabs( ( J3 - JA ) / JA );
1325 double T4 = fabs( ( J4 - JA ) / JA );
1327 return Max( Max( T1, T2 ), Max( T3, T4 ) );
1330 double Taper::GetBadRate( double Value, int /*nbNodes*/ ) const
1332 // the taper is in the range [0.0,1.0]
1333 // 0.0 = good (no taper)
1334 // 1.0 = bad (les cotes opposes sont allignes)
1338 SMDSAbs_ElementType Taper::GetType() const
1340 return SMDSAbs_Face;
1346 Description : Functor for calculating skew in degrees
1348 static inline double skewAngle( const gp_XYZ& p1, const gp_XYZ& p2, const gp_XYZ& p3 )
1350 gp_XYZ p12 = ( p2 + p1 ) / 2.;
1351 gp_XYZ p23 = ( p3 + p2 ) / 2.;
1352 gp_XYZ p31 = ( p3 + p1 ) / 2.;
1354 gp_Vec v1( p31 - p2 ), v2( p12 - p23 );
1356 return v1.Magnitude() < gp::Resolution() || v2.Magnitude() < gp::Resolution() ? 0. : v1.Angle( v2 );
1359 double Skew::GetValue( const TSequenceOfXYZ& P )
1361 if ( P.size() != 3 && P.size() != 4 )
1365 static double PI2 = PI / 2.;
1366 if ( P.size() == 3 )
1368 double A0 = fabs( PI2 - skewAngle( P( 3 ), P( 1 ), P( 2 ) ) );
1369 double A1 = fabs( PI2 - skewAngle( P( 1 ), P( 2 ), P( 3 ) ) );
1370 double A2 = fabs( PI2 - skewAngle( P( 2 ), P( 3 ), P( 1 ) ) );
1372 return Max( A0, Max( A1, A2 ) ) * 180. / PI;
1376 gp_XYZ p12 = ( P( 1 ) + P( 2 ) ) / 2.;
1377 gp_XYZ p23 = ( P( 2 ) + P( 3 ) ) / 2.;
1378 gp_XYZ p34 = ( P( 3 ) + P( 4 ) ) / 2.;
1379 gp_XYZ p41 = ( P( 4 ) + P( 1 ) ) / 2.;
1381 gp_Vec v1( p34 - p12 ), v2( p23 - p41 );
1382 double A = v1.Magnitude() <= gp::Resolution() || v2.Magnitude() <= gp::Resolution()
1383 ? 0. : fabs( PI2 - v1.Angle( v2 ) );
1386 if ( A < Precision::Angular() )
1389 return A * 180. / PI;
1393 double Skew::GetBadRate( double Value, int /*nbNodes*/ ) const
1395 // the skew is in the range [0.0,PI/2].
1401 SMDSAbs_ElementType Skew::GetType() const
1403 return SMDSAbs_Face;
1409 Description : Functor for calculating area
1411 double Area::GetValue( const TSequenceOfXYZ& P )
1414 if ( P.size() > 2 ) {
1415 gp_Vec aVec1( P(2) - P(1) );
1416 gp_Vec aVec2( P(3) - P(1) );
1417 gp_Vec SumVec = aVec1 ^ aVec2;
1418 for (int i=4; i<=P.size(); i++) {
1419 gp_Vec aVec1( P(i-1) - P(1) );
1420 gp_Vec aVec2( P(i) - P(1) );
1421 gp_Vec tmp = aVec1 ^ aVec2;
1424 val = SumVec.Magnitude() * 0.5;
1429 double Area::GetBadRate( double Value, int /*nbNodes*/ ) const
1431 // meaningless as it is not a quality control functor
1435 SMDSAbs_ElementType Area::GetType() const
1437 return SMDSAbs_Face;
1443 Description : Functor for calculating length of edge
1445 double Length::GetValue( const TSequenceOfXYZ& P )
1447 switch ( P.size() ) {
1448 case 2: return getDistance( P( 1 ), P( 2 ) );
1449 case 3: return getDistance( P( 1 ), P( 2 ) ) + getDistance( P( 2 ), P( 3 ) );
1454 double Length::GetBadRate( double Value, int /*nbNodes*/ ) const
1456 // meaningless as it is not quality control functor
1460 SMDSAbs_ElementType Length::GetType() const
1462 return SMDSAbs_Edge;
1467 Description : Functor for calculating length of edge
1470 double Length2D::GetValue( long theElementId)
1474 //cout<<"Length2D::GetValue"<<endl;
1475 if (GetPoints(theElementId,P)){
1476 //for(int jj=1; jj<=P.size(); jj++)
1477 // cout<<"jj="<<jj<<" P("<<P(jj).X()<<","<<P(jj).Y()<<","<<P(jj).Z()<<")"<<endl;
1479 double aVal;// = GetValue( P );
1480 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
1481 SMDSAbs_ElementType aType = aElem->GetType();
1490 aVal = getDistance( P( 1 ), P( 2 ) );
1493 else if (len == 3){ // quadratic edge
1494 aVal = getDistance(P( 1 ),P( 3 )) + getDistance(P( 3 ),P( 2 ));
1498 if (len == 3){ // triangles
1499 double L1 = getDistance(P( 1 ),P( 2 ));
1500 double L2 = getDistance(P( 2 ),P( 3 ));
1501 double L3 = getDistance(P( 3 ),P( 1 ));
1502 aVal = Max(L1,Max(L2,L3));
1505 else if (len == 4){ // quadrangles
1506 double L1 = getDistance(P( 1 ),P( 2 ));
1507 double L2 = getDistance(P( 2 ),P( 3 ));
1508 double L3 = getDistance(P( 3 ),P( 4 ));
1509 double L4 = getDistance(P( 4 ),P( 1 ));
1510 aVal = Max(Max(L1,L2),Max(L3,L4));
1513 if (len == 6){ // quadratic triangles
1514 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1515 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1516 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
1517 aVal = Max(L1,Max(L2,L3));
1518 //cout<<"L1="<<L1<<" L2="<<L2<<"L3="<<L3<<" aVal="<<aVal<<endl;
1521 else if (len == 8){ // quadratic quadrangles
1522 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1523 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1524 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
1525 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
1526 aVal = Max(Max(L1,L2),Max(L3,L4));
1529 case SMDSAbs_Volume:
1530 if (len == 4){ // tetraidrs
1531 double L1 = getDistance(P( 1 ),P( 2 ));
1532 double L2 = getDistance(P( 2 ),P( 3 ));
1533 double L3 = getDistance(P( 3 ),P( 1 ));
1534 double L4 = getDistance(P( 1 ),P( 4 ));
1535 double L5 = getDistance(P( 2 ),P( 4 ));
1536 double L6 = getDistance(P( 3 ),P( 4 ));
1537 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1540 else if (len == 5){ // piramids
1541 double L1 = getDistance(P( 1 ),P( 2 ));
1542 double L2 = getDistance(P( 2 ),P( 3 ));
1543 double L3 = getDistance(P( 3 ),P( 4 ));
1544 double L4 = getDistance(P( 4 ),P( 1 ));
1545 double L5 = getDistance(P( 1 ),P( 5 ));
1546 double L6 = getDistance(P( 2 ),P( 5 ));
1547 double L7 = getDistance(P( 3 ),P( 5 ));
1548 double L8 = getDistance(P( 4 ),P( 5 ));
1550 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1551 aVal = Max(aVal,Max(L7,L8));
1554 else if (len == 6){ // pentaidres
1555 double L1 = getDistance(P( 1 ),P( 2 ));
1556 double L2 = getDistance(P( 2 ),P( 3 ));
1557 double L3 = getDistance(P( 3 ),P( 1 ));
1558 double L4 = getDistance(P( 4 ),P( 5 ));
1559 double L5 = getDistance(P( 5 ),P( 6 ));
1560 double L6 = getDistance(P( 6 ),P( 4 ));
1561 double L7 = getDistance(P( 1 ),P( 4 ));
1562 double L8 = getDistance(P( 2 ),P( 5 ));
1563 double L9 = getDistance(P( 3 ),P( 6 ));
1565 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1566 aVal = Max(aVal,Max(Max(L7,L8),L9));
1569 else if (len == 8){ // hexaider
1570 double L1 = getDistance(P( 1 ),P( 2 ));
1571 double L2 = getDistance(P( 2 ),P( 3 ));
1572 double L3 = getDistance(P( 3 ),P( 4 ));
1573 double L4 = getDistance(P( 4 ),P( 1 ));
1574 double L5 = getDistance(P( 5 ),P( 6 ));
1575 double L6 = getDistance(P( 6 ),P( 7 ));
1576 double L7 = getDistance(P( 7 ),P( 8 ));
1577 double L8 = getDistance(P( 8 ),P( 5 ));
1578 double L9 = getDistance(P( 1 ),P( 5 ));
1579 double L10= getDistance(P( 2 ),P( 6 ));
1580 double L11= getDistance(P( 3 ),P( 7 ));
1581 double L12= getDistance(P( 4 ),P( 8 ));
1583 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1584 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1585 aVal = Max(aVal,Max(L11,L12));
1590 if (len == 10){ // quadratic tetraidrs
1591 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
1592 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
1593 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
1594 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1595 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
1596 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
1597 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1600 else if (len == 13){ // quadratic piramids
1601 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
1602 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
1603 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1604 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1605 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1606 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
1607 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
1608 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
1609 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1610 aVal = Max(aVal,Max(L7,L8));
1613 else if (len == 15){ // quadratic pentaidres
1614 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
1615 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
1616 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1617 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1618 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
1619 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
1620 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
1621 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
1622 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
1623 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1624 aVal = Max(aVal,Max(Max(L7,L8),L9));
1627 else if (len == 20){ // quadratic hexaider
1628 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
1629 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
1630 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
1631 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
1632 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
1633 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
1634 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
1635 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
1636 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
1637 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
1638 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
1639 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
1640 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1641 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1642 aVal = Max(aVal,Max(L11,L12));
1654 if ( myPrecision >= 0 )
1656 double prec = pow( 10., (double)( myPrecision ) );
1657 aVal = floor( aVal * prec + 0.5 ) / prec;
1666 double Length2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1668 // meaningless as it is not quality control functor
1672 SMDSAbs_ElementType Length2D::GetType() const
1674 return SMDSAbs_Face;
1677 Length2D::Value::Value(double theLength,long thePntId1, long thePntId2):
1680 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1681 if(thePntId1 > thePntId2){
1682 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1686 bool Length2D::Value::operator<(const Length2D::Value& x) const{
1687 if(myPntId[0] < x.myPntId[0]) return true;
1688 if(myPntId[0] == x.myPntId[0])
1689 if(myPntId[1] < x.myPntId[1]) return true;
1693 void Length2D::GetValues(TValues& theValues){
1695 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1696 for(; anIter->more(); ){
1697 const SMDS_MeshFace* anElem = anIter->next();
1699 if(anElem->IsQuadratic()) {
1700 const SMDS_VtkFace* F =
1701 dynamic_cast<const SMDS_VtkFace*>(anElem);
1702 // use special nodes iterator
1703 SMDS_ElemIteratorPtr anIter = F->interlacedNodesElemIterator();
1708 const SMDS_MeshElement* aNode;
1710 aNode = anIter->next();
1711 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1712 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1713 aNodeId[0] = aNodeId[1] = aNode->GetID();
1716 for(; anIter->more(); ){
1717 const SMDS_MeshNode* N1 = static_cast<const SMDS_MeshNode*> (anIter->next());
1718 P[2] = gp_Pnt(N1->X(),N1->Y(),N1->Z());
1719 aNodeId[2] = N1->GetID();
1720 aLength = P[1].Distance(P[2]);
1721 if(!anIter->more()) break;
1722 const SMDS_MeshNode* N2 = static_cast<const SMDS_MeshNode*> (anIter->next());
1723 P[3] = gp_Pnt(N2->X(),N2->Y(),N2->Z());
1724 aNodeId[3] = N2->GetID();
1725 aLength += P[2].Distance(P[3]);
1726 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1727 Value aValue2(aLength,aNodeId[2],aNodeId[3]);
1729 aNodeId[1] = aNodeId[3];
1730 theValues.insert(aValue1);
1731 theValues.insert(aValue2);
1733 aLength += P[2].Distance(P[0]);
1734 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1735 Value aValue2(aLength,aNodeId[2],aNodeId[0]);
1736 theValues.insert(aValue1);
1737 theValues.insert(aValue2);
1740 SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
1745 const SMDS_MeshElement* aNode;
1746 if(aNodesIter->more()){
1747 aNode = aNodesIter->next();
1748 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1749 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1750 aNodeId[0] = aNodeId[1] = aNode->GetID();
1753 for(; aNodesIter->more(); ){
1754 aNode = aNodesIter->next();
1755 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1756 long anId = aNode->GetID();
1758 P[2] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1760 aLength = P[1].Distance(P[2]);
1762 Value aValue(aLength,aNodeId[1],anId);
1765 theValues.insert(aValue);
1768 aLength = P[0].Distance(P[1]);
1770 Value aValue(aLength,aNodeId[0],aNodeId[1]);
1771 theValues.insert(aValue);
1777 Class : MultiConnection
1778 Description : Functor for calculating number of faces conneted to the edge
1780 double MultiConnection::GetValue( const TSequenceOfXYZ& P )
1784 double MultiConnection::GetValue( long theId )
1786 return getNbMultiConnection( myMesh, theId );
1789 double MultiConnection::GetBadRate( double Value, int /*nbNodes*/ ) const
1791 // meaningless as it is not quality control functor
1795 SMDSAbs_ElementType MultiConnection::GetType() const
1797 return SMDSAbs_Edge;
1801 Class : MultiConnection2D
1802 Description : Functor for calculating number of faces conneted to the edge
1804 double MultiConnection2D::GetValue( const TSequenceOfXYZ& P )
1809 double MultiConnection2D::GetValue( long theElementId )
1813 const SMDS_MeshElement* aFaceElem = myMesh->FindElement(theElementId);
1814 SMDSAbs_ElementType aType = aFaceElem->GetType();
1819 int i = 0, len = aFaceElem->NbNodes();
1820 SMDS_ElemIteratorPtr anIter = aFaceElem->nodesIterator();
1823 const SMDS_MeshNode *aNode, *aNode0;
1824 TColStd_MapOfInteger aMap, aMapPrev;
1826 for (i = 0; i <= len; i++) {
1831 if (anIter->more()) {
1832 aNode = (SMDS_MeshNode*)anIter->next();
1840 if (i == 0) aNode0 = aNode;
1842 SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
1843 while (anElemIter->more()) {
1844 const SMDS_MeshElement* anElem = anElemIter->next();
1845 if (anElem != 0 && anElem->GetType() == SMDSAbs_Face) {
1846 int anId = anElem->GetID();
1849 if (aMapPrev.Contains(anId)) {
1854 aResult = Max(aResult, aNb);
1865 double MultiConnection2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1867 // meaningless as it is not quality control functor
1871 SMDSAbs_ElementType MultiConnection2D::GetType() const
1873 return SMDSAbs_Face;
1876 MultiConnection2D::Value::Value(long thePntId1, long thePntId2)
1878 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1879 if(thePntId1 > thePntId2){
1880 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1884 bool MultiConnection2D::Value::operator<(const MultiConnection2D::Value& x) const{
1885 if(myPntId[0] < x.myPntId[0]) return true;
1886 if(myPntId[0] == x.myPntId[0])
1887 if(myPntId[1] < x.myPntId[1]) return true;
1891 void MultiConnection2D::GetValues(MValues& theValues){
1892 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1893 for(; anIter->more(); ){
1894 const SMDS_MeshFace* anElem = anIter->next();
1895 SMDS_ElemIteratorPtr aNodesIter;
1896 if ( anElem->IsQuadratic() )
1897 aNodesIter = dynamic_cast<const SMDS_VtkFace*>
1898 (anElem)->interlacedNodesElemIterator();
1900 aNodesIter = anElem->nodesIterator();
1903 //int aNbConnects=0;
1904 const SMDS_MeshNode* aNode0;
1905 const SMDS_MeshNode* aNode1;
1906 const SMDS_MeshNode* aNode2;
1907 if(aNodesIter->more()){
1908 aNode0 = (SMDS_MeshNode*) aNodesIter->next();
1910 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode1;
1911 aNodeId[0] = aNodeId[1] = aNodes->GetID();
1913 for(; aNodesIter->more(); ) {
1914 aNode2 = (SMDS_MeshNode*) aNodesIter->next();
1915 long anId = aNode2->GetID();
1918 Value aValue(aNodeId[1],aNodeId[2]);
1919 MValues::iterator aItr = theValues.find(aValue);
1920 if (aItr != theValues.end()){
1925 theValues[aValue] = 1;
1928 //cout << "NodeIds: "<<aNodeId[1]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1929 aNodeId[1] = aNodeId[2];
1932 Value aValue(aNodeId[0],aNodeId[2]);
1933 MValues::iterator aItr = theValues.find(aValue);
1934 if (aItr != theValues.end()) {
1939 theValues[aValue] = 1;
1942 //cout << "NodeIds: "<<aNodeId[0]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1952 Class : BadOrientedVolume
1953 Description : Predicate bad oriented volumes
1956 BadOrientedVolume::BadOrientedVolume()
1961 void BadOrientedVolume::SetMesh( const SMDS_Mesh* theMesh )
1966 bool BadOrientedVolume::IsSatisfy( long theId )
1971 SMDS_VolumeTool vTool( myMesh->FindElement( theId ));
1972 return !vTool.IsForward();
1975 SMDSAbs_ElementType BadOrientedVolume::GetType() const
1977 return SMDSAbs_Volume;
1981 Class : BareBorderVolume
1984 bool BareBorderVolume::IsSatisfy(long theElementId )
1986 SMDS_VolumeTool myTool;
1987 if ( myTool.Set( myMesh->FindElement(theElementId)))
1989 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
1990 if ( myTool.IsFreeFace( iF ))
1992 const SMDS_MeshNode** n = myTool.GetFaceNodes(iF);
1993 vector< const SMDS_MeshNode*> nodes( n, n+myTool.NbFaceNodes(iF));
1994 if ( !myMesh->FindElement( nodes, SMDSAbs_Face, /*Nomedium=*/false))
2002 Class : BareBorderFace
2005 bool BareBorderFace::IsSatisfy(long theElementId )
2008 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2010 if ( face->GetType() == SMDSAbs_Face )
2012 int nbN = face->NbCornerNodes();
2013 for ( int i = 0; i < nbN && !ok; ++i )
2015 // check if a link is shared by another face
2016 const SMDS_MeshNode* n1 = face->GetNode( i );
2017 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2018 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2019 bool isShared = false;
2020 while ( !isShared && fIt->more() )
2022 const SMDS_MeshElement* f = fIt->next();
2023 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2027 myLinkNodes.resize( 2 + face->IsQuadratic());
2028 myLinkNodes[0] = n1;
2029 myLinkNodes[1] = n2;
2030 if ( face->IsQuadratic() )
2031 myLinkNodes[2] = face->GetNode( i+nbN );
2032 ok = !myMesh->FindElement( myLinkNodes, SMDSAbs_Edge, /*noMedium=*/false);
2041 Class : OverConstrainedVolume
2044 bool OverConstrainedVolume::IsSatisfy(long theElementId )
2046 // An element is over-constrained if it has N-1 free borders where
2047 // N is the number of edges/faces for a 2D/3D element.
2048 SMDS_VolumeTool myTool;
2049 if ( myTool.Set( myMesh->FindElement(theElementId)))
2051 int nbSharedFaces = 0;
2052 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2053 if ( !myTool.IsFreeFace( iF ) && ++nbSharedFaces > 1 )
2055 return ( nbSharedFaces == 1 );
2061 Class : OverConstrainedFace
2064 bool OverConstrainedFace::IsSatisfy(long theElementId )
2066 // An element is over-constrained if it has N-1 free borders where
2067 // N is the number of edges/faces for a 2D/3D element.
2068 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2069 if ( face->GetType() == SMDSAbs_Face )
2071 int nbSharedBorders = 0;
2072 int nbN = face->NbCornerNodes();
2073 for ( int i = 0; i < nbN; ++i )
2075 // check if a link is shared by another face
2076 const SMDS_MeshNode* n1 = face->GetNode( i );
2077 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2078 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2079 bool isShared = false;
2080 while ( !isShared && fIt->more() )
2082 const SMDS_MeshElement* f = fIt->next();
2083 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2085 if ( isShared && ++nbSharedBorders > 1 )
2088 return ( nbSharedBorders == 1 );
2095 Description : Predicate for free borders
2098 FreeBorders::FreeBorders()
2103 void FreeBorders::SetMesh( const SMDS_Mesh* theMesh )
2108 bool FreeBorders::IsSatisfy( long theId )
2110 return getNbMultiConnection( myMesh, theId ) == 1;
2113 SMDSAbs_ElementType FreeBorders::GetType() const
2115 return SMDSAbs_Edge;
2121 Description : Predicate for free Edges
2123 FreeEdges::FreeEdges()
2128 void FreeEdges::SetMesh( const SMDS_Mesh* theMesh )
2133 bool FreeEdges::IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId )
2135 TColStd_MapOfInteger aMap;
2136 for ( int i = 0; i < 2; i++ )
2138 SMDS_ElemIteratorPtr anElemIter = theNodes[ i ]->GetInverseElementIterator(SMDSAbs_Face);
2139 while( anElemIter->more() )
2141 if ( const SMDS_MeshElement* anElem = anElemIter->next())
2143 const int anId = anElem->GetID();
2144 if ( anId != theFaceId && !aMap.Add( anId ))
2152 bool FreeEdges::IsSatisfy( long theId )
2157 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2158 if ( aFace == 0 || aFace->GetType() != SMDSAbs_Face || aFace->NbNodes() < 3 )
2161 SMDS_ElemIteratorPtr anIter;
2162 if ( aFace->IsQuadratic() ) {
2163 anIter = dynamic_cast<const SMDS_VtkFace*>
2164 (aFace)->interlacedNodesElemIterator();
2167 anIter = aFace->nodesIterator();
2172 int i = 0, nbNodes = aFace->NbNodes();
2173 std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
2174 while( anIter->more() )
2176 const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
2179 aNodes[ i++ ] = aNode;
2181 aNodes[ nbNodes ] = aNodes[ 0 ];
2183 for ( i = 0; i < nbNodes; i++ )
2184 if ( IsFreeEdge( &aNodes[ i ], theId ) )
2190 SMDSAbs_ElementType FreeEdges::GetType() const
2192 return SMDSAbs_Face;
2195 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2198 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
2199 if(thePntId1 > thePntId2){
2200 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
2204 bool FreeEdges::Border::operator<(const FreeEdges::Border& x) const{
2205 if(myPntId[0] < x.myPntId[0]) return true;
2206 if(myPntId[0] == x.myPntId[0])
2207 if(myPntId[1] < x.myPntId[1]) return true;
2211 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2212 FreeEdges::TBorders& theRegistry,
2213 FreeEdges::TBorders& theContainer)
2215 if(theRegistry.find(theBorder) == theRegistry.end()){
2216 theRegistry.insert(theBorder);
2217 theContainer.insert(theBorder);
2219 theContainer.erase(theBorder);
2223 void FreeEdges::GetBoreders(TBorders& theBorders)
2226 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2227 for(; anIter->more(); ){
2228 const SMDS_MeshFace* anElem = anIter->next();
2229 long anElemId = anElem->GetID();
2230 SMDS_ElemIteratorPtr aNodesIter;
2231 if ( anElem->IsQuadratic() )
2232 aNodesIter = static_cast<const SMDS_VtkFace*>(anElem)->
2233 interlacedNodesElemIterator();
2235 aNodesIter = anElem->nodesIterator();
2237 const SMDS_MeshElement* aNode;
2238 if(aNodesIter->more()){
2239 aNode = aNodesIter->next();
2240 aNodeId[0] = aNodeId[1] = aNode->GetID();
2242 for(; aNodesIter->more(); ){
2243 aNode = aNodesIter->next();
2244 long anId = aNode->GetID();
2245 Border aBorder(anElemId,aNodeId[1],anId);
2247 UpdateBorders(aBorder,aRegistry,theBorders);
2249 Border aBorder(anElemId,aNodeId[0],aNodeId[1]);
2250 UpdateBorders(aBorder,aRegistry,theBorders);
2257 Description : Predicate for free nodes
2260 FreeNodes::FreeNodes()
2265 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2270 bool FreeNodes::IsSatisfy( long theNodeId )
2272 const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2276 return (aNode->NbInverseElements() < 1);
2279 SMDSAbs_ElementType FreeNodes::GetType() const
2281 return SMDSAbs_Node;
2287 Description : Predicate for free faces
2290 FreeFaces::FreeFaces()
2295 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2300 bool FreeFaces::IsSatisfy( long theId )
2302 if (!myMesh) return false;
2303 // check that faces nodes refers to less than two common volumes
2304 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2305 if ( !aFace || aFace->GetType() != SMDSAbs_Face )
2308 int nbNode = aFace->NbNodes();
2310 // collect volumes check that number of volumss with count equal nbNode not less than 2
2311 typedef map< SMDS_MeshElement*, int > TMapOfVolume; // map of volume counters
2312 typedef map< SMDS_MeshElement*, int >::iterator TItrMapOfVolume; // iterator
2313 TMapOfVolume mapOfVol;
2315 SMDS_ElemIteratorPtr nodeItr = aFace->nodesIterator();
2316 while ( nodeItr->more() ) {
2317 const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>(nodeItr->next());
2318 if ( !aNode ) continue;
2319 SMDS_ElemIteratorPtr volItr = aNode->GetInverseElementIterator(SMDSAbs_Volume);
2320 while ( volItr->more() ) {
2321 SMDS_MeshElement* aVol = (SMDS_MeshElement*)volItr->next();
2322 TItrMapOfVolume itr = mapOfVol.insert(make_pair(aVol, 0)).first;
2327 TItrMapOfVolume volItr = mapOfVol.begin();
2328 TItrMapOfVolume volEnd = mapOfVol.end();
2329 for ( ; volItr != volEnd; ++volItr )
2330 if ( (*volItr).second >= nbNode )
2332 // face is not free if number of volumes constructed on thier nodes more than one
2336 SMDSAbs_ElementType FreeFaces::GetType() const
2338 return SMDSAbs_Face;
2342 Class : LinearOrQuadratic
2343 Description : Predicate to verify whether a mesh element is linear
2346 LinearOrQuadratic::LinearOrQuadratic()
2351 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2356 bool LinearOrQuadratic::IsSatisfy( long theId )
2358 if (!myMesh) return false;
2359 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2360 if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2362 return (!anElem->IsQuadratic());
2365 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2370 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2377 Description : Functor for check color of group to whic mesh element belongs to
2380 GroupColor::GroupColor()
2384 bool GroupColor::IsSatisfy( long theId )
2386 return (myIDs.find( theId ) != myIDs.end());
2389 void GroupColor::SetType( SMDSAbs_ElementType theType )
2394 SMDSAbs_ElementType GroupColor::GetType() const
2399 static bool isEqual( const Quantity_Color& theColor1,
2400 const Quantity_Color& theColor2 )
2402 // tolerance to compare colors
2403 const double tol = 5*1e-3;
2404 return ( fabs( theColor1.Red() - theColor2.Red() ) < tol &&
2405 fabs( theColor1.Green() - theColor2.Green() ) < tol &&
2406 fabs( theColor1.Blue() - theColor2.Blue() ) < tol );
2410 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2414 const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2418 int nbGrp = aMesh->GetNbGroups();
2422 // iterates on groups and find necessary elements ids
2423 const std::set<SMESHDS_GroupBase*>& aGroups = aMesh->GetGroups();
2424 set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
2425 for (; GrIt != aGroups.end(); GrIt++) {
2426 SMESHDS_GroupBase* aGrp = (*GrIt);
2429 // check type and color of group
2430 if ( !isEqual( myColor, aGrp->GetColor() ) )
2432 if ( myType != SMDSAbs_All && myType != (SMDSAbs_ElementType)aGrp->GetType() )
2435 SMDSAbs_ElementType aGrpElType = (SMDSAbs_ElementType)aGrp->GetType();
2436 if ( myType == aGrpElType || (myType == SMDSAbs_All && aGrpElType != SMDSAbs_Node) ) {
2437 // add elements IDS into control
2438 int aSize = aGrp->Extent();
2439 for (int i = 0; i < aSize; i++)
2440 myIDs.insert( aGrp->GetID(i+1) );
2445 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
2447 TCollection_AsciiString aStr = theStr;
2448 aStr.RemoveAll( ' ' );
2449 aStr.RemoveAll( '\t' );
2450 for ( int aPos = aStr.Search( ";;" ); aPos != -1; aPos = aStr.Search( ";;" ) )
2451 aStr.Remove( aPos, 2 );
2452 Standard_Real clr[3];
2453 clr[0] = clr[1] = clr[2] = 0.;
2454 for ( int i = 0; i < 3; i++ ) {
2455 TCollection_AsciiString tmpStr = aStr.Token( ";", i+1 );
2456 if ( !tmpStr.IsEmpty() && tmpStr.IsRealValue() )
2457 clr[i] = tmpStr.RealValue();
2459 myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
2462 //=======================================================================
2463 // name : GetRangeStr
2464 // Purpose : Get range as a string.
2465 // Example: "1,2,3,50-60,63,67,70-"
2466 //=======================================================================
2467 void GroupColor::GetColorStr( TCollection_AsciiString& theResStr ) const
2470 theResStr += TCollection_AsciiString( myColor.Red() );
2471 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2472 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2476 Class : ElemGeomType
2477 Description : Predicate to check element geometry type
2480 ElemGeomType::ElemGeomType()
2483 myType = SMDSAbs_All;
2484 myGeomType = SMDSGeom_TRIANGLE;
2487 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2492 bool ElemGeomType::IsSatisfy( long theId )
2494 if (!myMesh) return false;
2495 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2498 const SMDSAbs_ElementType anElemType = anElem->GetType();
2499 if ( myType != SMDSAbs_All && anElemType != myType )
2501 const int aNbNode = anElem->NbNodes();
2503 switch( anElemType )
2506 isOk = (myGeomType == SMDSGeom_POINT);
2510 isOk = (myGeomType == SMDSGeom_EDGE);
2514 if ( myGeomType == SMDSGeom_TRIANGLE )
2515 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 6 : aNbNode == 3));
2516 else if ( myGeomType == SMDSGeom_QUADRANGLE )
2517 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? ( aNbNode == 8 || aNbNode == 9 ) : aNbNode == 4));
2518 else if ( myGeomType == SMDSGeom_POLYGON )
2519 isOk = anElem->IsPoly();
2522 case SMDSAbs_Volume:
2523 if ( myGeomType == SMDSGeom_TETRA )
2524 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 10 : aNbNode == 4));
2525 else if ( myGeomType == SMDSGeom_PYRAMID )
2526 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 13 : aNbNode == 5));
2527 else if ( myGeomType == SMDSGeom_PENTA )
2528 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 15 : aNbNode == 6));
2529 else if ( myGeomType == SMDSGeom_HEXA )
2530 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? ( aNbNode == 20 || aNbNode == 27 ): aNbNode == 8));
2531 else if ( myGeomType == SMDSGeom_HEXAGONAL_PRISM )
2532 isOk = (anElem->GetEntityType() == SMDSEntity_Hexagonal_Prism );
2533 else if ( myGeomType == SMDSGeom_POLYHEDRA )
2534 isOk = anElem->IsPoly();
2541 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2546 SMDSAbs_ElementType ElemGeomType::GetType() const
2551 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2553 myGeomType = theType;
2556 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2561 //================================================================================
2563 * \brief Class CoplanarFaces
2565 //================================================================================
2567 CoplanarFaces::CoplanarFaces()
2568 : myMesh(0), myFaceID(0), myToler(0)
2571 bool CoplanarFaces::IsSatisfy( long theElementId )
2573 if ( myCoplanarIDs.empty() )
2575 // Build a set of coplanar face ids
2577 if ( !myMesh || !myFaceID || !myToler )
2580 const SMDS_MeshElement* face = myMesh->FindElement( myFaceID );
2581 if ( !face || face->GetType() != SMDSAbs_Face )
2585 gp_Vec myNorm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2589 const double radianTol = myToler * PI180;
2590 typedef SMDS_StdIterator< const SMDS_MeshElement*, SMDS_ElemIteratorPtr > TFaceIt;
2591 std::set<const SMDS_MeshElement*> checkedFaces, checkedNodes;
2592 std::list<const SMDS_MeshElement*> faceQueue( 1, face );
2593 while ( !faceQueue.empty() )
2595 face = faceQueue.front();
2596 if ( checkedFaces.insert( face ).second )
2598 gp_Vec norm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2599 if (!normOK || myNorm.Angle( norm ) <= radianTol)
2601 myCoplanarIDs.insert( face->GetID() );
2602 std::set<const SMDS_MeshElement*> neighborFaces;
2603 for ( int i = 0; i < face->NbCornerNodes(); ++i )
2605 const SMDS_MeshNode* n = face->GetNode( i );
2606 if ( checkedNodes.insert( n ).second )
2607 neighborFaces.insert( TFaceIt( n->GetInverseElementIterator(SMDSAbs_Face)),
2610 faceQueue.insert( faceQueue.end(), neighborFaces.begin(), neighborFaces.end() );
2613 faceQueue.pop_front();
2616 return myCoplanarIDs.count( theElementId );
2621 *Description : Predicate for Range of Ids.
2622 * Range may be specified with two ways.
2623 * 1. Using AddToRange method
2624 * 2. With SetRangeStr method. Parameter of this method is a string
2625 * like as "1,2,3,50-60,63,67,70-"
2628 //=======================================================================
2629 // name : RangeOfIds
2630 // Purpose : Constructor
2631 //=======================================================================
2632 RangeOfIds::RangeOfIds()
2635 myType = SMDSAbs_All;
2638 //=======================================================================
2640 // Purpose : Set mesh
2641 //=======================================================================
2642 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
2647 //=======================================================================
2648 // name : AddToRange
2649 // Purpose : Add ID to the range
2650 //=======================================================================
2651 bool RangeOfIds::AddToRange( long theEntityId )
2653 myIds.Add( theEntityId );
2657 //=======================================================================
2658 // name : GetRangeStr
2659 // Purpose : Get range as a string.
2660 // Example: "1,2,3,50-60,63,67,70-"
2661 //=======================================================================
2662 void RangeOfIds::GetRangeStr( TCollection_AsciiString& theResStr )
2666 TColStd_SequenceOfInteger anIntSeq;
2667 TColStd_SequenceOfAsciiString aStrSeq;
2669 TColStd_MapIteratorOfMapOfInteger anIter( myIds );
2670 for ( ; anIter.More(); anIter.Next() )
2672 int anId = anIter.Key();
2673 TCollection_AsciiString aStr( anId );
2674 anIntSeq.Append( anId );
2675 aStrSeq.Append( aStr );
2678 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2680 int aMinId = myMin( i );
2681 int aMaxId = myMax( i );
2683 TCollection_AsciiString aStr;
2684 if ( aMinId != IntegerFirst() )
2689 if ( aMaxId != IntegerLast() )
2692 // find position of the string in result sequence and insert string in it
2693 if ( anIntSeq.Length() == 0 )
2695 anIntSeq.Append( aMinId );
2696 aStrSeq.Append( aStr );
2700 if ( aMinId < anIntSeq.First() )
2702 anIntSeq.Prepend( aMinId );
2703 aStrSeq.Prepend( aStr );
2705 else if ( aMinId > anIntSeq.Last() )
2707 anIntSeq.Append( aMinId );
2708 aStrSeq.Append( aStr );
2711 for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
2712 if ( aMinId < anIntSeq( j ) )
2714 anIntSeq.InsertBefore( j, aMinId );
2715 aStrSeq.InsertBefore( j, aStr );
2721 if ( aStrSeq.Length() == 0 )
2724 theResStr = aStrSeq( 1 );
2725 for ( int j = 2, k = aStrSeq.Length(); j <= k; j++ )
2728 theResStr += aStrSeq( j );
2732 //=======================================================================
2733 // name : SetRangeStr
2734 // Purpose : Define range with string
2735 // Example of entry string: "1,2,3,50-60,63,67,70-"
2736 //=======================================================================
2737 bool RangeOfIds::SetRangeStr( const TCollection_AsciiString& theStr )
2743 TCollection_AsciiString aStr = theStr;
2744 aStr.RemoveAll( ' ' );
2745 aStr.RemoveAll( '\t' );
2747 for ( int aPos = aStr.Search( ",," ); aPos != -1; aPos = aStr.Search( ",," ) )
2748 aStr.Remove( aPos, 2 );
2750 TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
2752 while ( tmpStr != "" )
2754 tmpStr = aStr.Token( ",", i++ );
2755 int aPos = tmpStr.Search( '-' );
2759 if ( tmpStr.IsIntegerValue() )
2760 myIds.Add( tmpStr.IntegerValue() );
2766 TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
2767 TCollection_AsciiString aMinStr = tmpStr;
2769 while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
2770 while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
2772 if ( (!aMinStr.IsEmpty() && !aMinStr.IsIntegerValue()) ||
2773 (!aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue()) )
2776 myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
2777 myMax.Append( aMaxStr.IsEmpty() ? IntegerLast() : aMaxStr.IntegerValue() );
2784 //=======================================================================
2786 // Purpose : Get type of supported entities
2787 //=======================================================================
2788 SMDSAbs_ElementType RangeOfIds::GetType() const
2793 //=======================================================================
2795 // Purpose : Set type of supported entities
2796 //=======================================================================
2797 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
2802 //=======================================================================
2804 // Purpose : Verify whether entity satisfies to this rpedicate
2805 //=======================================================================
2806 bool RangeOfIds::IsSatisfy( long theId )
2811 if ( myType == SMDSAbs_Node )
2813 if ( myMesh->FindNode( theId ) == 0 )
2818 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2819 if ( anElem == 0 || (myType != anElem->GetType() && myType != SMDSAbs_All ))
2823 if ( myIds.Contains( theId ) )
2826 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2827 if ( theId >= myMin( i ) && theId <= myMax( i ) )
2835 Description : Base class for comparators
2837 Comparator::Comparator():
2841 Comparator::~Comparator()
2844 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
2847 myFunctor->SetMesh( theMesh );
2850 void Comparator::SetMargin( double theValue )
2852 myMargin = theValue;
2855 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
2857 myFunctor = theFunct;
2860 SMDSAbs_ElementType Comparator::GetType() const
2862 return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
2865 double Comparator::GetMargin()
2873 Description : Comparator "<"
2875 bool LessThan::IsSatisfy( long theId )
2877 return myFunctor && myFunctor->GetValue( theId ) < myMargin;
2883 Description : Comparator ">"
2885 bool MoreThan::IsSatisfy( long theId )
2887 return myFunctor && myFunctor->GetValue( theId ) > myMargin;
2893 Description : Comparator "="
2896 myToler(Precision::Confusion())
2899 bool EqualTo::IsSatisfy( long theId )
2901 return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
2904 void EqualTo::SetTolerance( double theToler )
2909 double EqualTo::GetTolerance()
2916 Description : Logical NOT predicate
2918 LogicalNOT::LogicalNOT()
2921 LogicalNOT::~LogicalNOT()
2924 bool LogicalNOT::IsSatisfy( long theId )
2926 return myPredicate && !myPredicate->IsSatisfy( theId );
2929 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
2932 myPredicate->SetMesh( theMesh );
2935 void LogicalNOT::SetPredicate( PredicatePtr thePred )
2937 myPredicate = thePred;
2940 SMDSAbs_ElementType LogicalNOT::GetType() const
2942 return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
2947 Class : LogicalBinary
2948 Description : Base class for binary logical predicate
2950 LogicalBinary::LogicalBinary()
2953 LogicalBinary::~LogicalBinary()
2956 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
2959 myPredicate1->SetMesh( theMesh );
2962 myPredicate2->SetMesh( theMesh );
2965 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
2967 myPredicate1 = thePredicate;
2970 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
2972 myPredicate2 = thePredicate;
2975 SMDSAbs_ElementType LogicalBinary::GetType() const
2977 if ( !myPredicate1 || !myPredicate2 )
2980 SMDSAbs_ElementType aType1 = myPredicate1->GetType();
2981 SMDSAbs_ElementType aType2 = myPredicate2->GetType();
2983 return aType1 == aType2 ? aType1 : SMDSAbs_All;
2989 Description : Logical AND
2991 bool LogicalAND::IsSatisfy( long theId )
2996 myPredicate1->IsSatisfy( theId ) &&
2997 myPredicate2->IsSatisfy( theId );
3003 Description : Logical OR
3005 bool LogicalOR::IsSatisfy( long theId )
3010 (myPredicate1->IsSatisfy( theId ) ||
3011 myPredicate2->IsSatisfy( theId ));
3025 void Filter::SetPredicate( PredicatePtr thePredicate )
3027 myPredicate = thePredicate;
3030 template<class TElement, class TIterator, class TPredicate>
3031 inline void FillSequence(const TIterator& theIterator,
3032 TPredicate& thePredicate,
3033 Filter::TIdSequence& theSequence)
3035 if ( theIterator ) {
3036 while( theIterator->more() ) {
3037 TElement anElem = theIterator->next();
3038 long anId = anElem->GetID();
3039 if ( thePredicate->IsSatisfy( anId ) )
3040 theSequence.push_back( anId );
3047 GetElementsId( const SMDS_Mesh* theMesh,
3048 PredicatePtr thePredicate,
3049 TIdSequence& theSequence )
3051 theSequence.clear();
3053 if ( !theMesh || !thePredicate )
3056 thePredicate->SetMesh( theMesh );
3058 SMDSAbs_ElementType aType = thePredicate->GetType();
3061 FillSequence<const SMDS_MeshNode*>(theMesh->nodesIterator(),thePredicate,theSequence);
3064 FillSequence<const SMDS_MeshElement*>(theMesh->edgesIterator(),thePredicate,theSequence);
3067 FillSequence<const SMDS_MeshElement*>(theMesh->facesIterator(),thePredicate,theSequence);
3069 case SMDSAbs_Volume:
3070 FillSequence<const SMDS_MeshElement*>(theMesh->volumesIterator(),thePredicate,theSequence);
3073 FillSequence<const SMDS_MeshElement*>(theMesh->edgesIterator(),thePredicate,theSequence);
3074 FillSequence<const SMDS_MeshElement*>(theMesh->facesIterator(),thePredicate,theSequence);
3075 FillSequence<const SMDS_MeshElement*>(theMesh->volumesIterator(),thePredicate,theSequence);
3081 Filter::GetElementsId( const SMDS_Mesh* theMesh,
3082 Filter::TIdSequence& theSequence )
3084 GetElementsId(theMesh,myPredicate,theSequence);
3091 typedef std::set<SMDS_MeshFace*> TMapOfFacePtr;
3097 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
3098 SMDS_MeshNode* theNode2 )
3104 ManifoldPart::Link::~Link()
3110 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
3112 if ( myNode1 == theLink.myNode1 &&
3113 myNode2 == theLink.myNode2 )
3115 else if ( myNode1 == theLink.myNode2 &&
3116 myNode2 == theLink.myNode1 )
3122 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
3124 if(myNode1 < x.myNode1) return true;
3125 if(myNode1 == x.myNode1)
3126 if(myNode2 < x.myNode2) return true;
3130 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
3131 const ManifoldPart::Link& theLink2 )
3133 return theLink1.IsEqual( theLink2 );
3136 ManifoldPart::ManifoldPart()
3139 myAngToler = Precision::Angular();
3140 myIsOnlyManifold = true;
3143 ManifoldPart::~ManifoldPart()
3148 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
3154 SMDSAbs_ElementType ManifoldPart::GetType() const
3155 { return SMDSAbs_Face; }
3157 bool ManifoldPart::IsSatisfy( long theElementId )
3159 return myMapIds.Contains( theElementId );
3162 void ManifoldPart::SetAngleTolerance( const double theAngToler )
3163 { myAngToler = theAngToler; }
3165 double ManifoldPart::GetAngleTolerance() const
3166 { return myAngToler; }
3168 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
3169 { myIsOnlyManifold = theIsOnly; }
3171 void ManifoldPart::SetStartElem( const long theStartId )
3172 { myStartElemId = theStartId; }
3174 bool ManifoldPart::process()
3177 myMapBadGeomIds.Clear();
3179 myAllFacePtr.clear();
3180 myAllFacePtrIntDMap.clear();
3184 // collect all faces into own map
3185 SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3186 for (; anFaceItr->more(); )
3188 SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3189 myAllFacePtr.push_back( aFacePtr );
3190 myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3193 SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3197 // the map of non manifold links and bad geometry
3198 TMapOfLink aMapOfNonManifold;
3199 TColStd_MapOfInteger aMapOfTreated;
3201 // begin cycle on faces from start index and run on vector till the end
3202 // and from begin to start index to cover whole vector
3203 const int aStartIndx = myAllFacePtrIntDMap[aStartFace];
3204 bool isStartTreat = false;
3205 for ( int fi = aStartIndx; !isStartTreat || fi != aStartIndx ; fi++ )
3207 if ( fi == aStartIndx )
3208 isStartTreat = true;
3209 // as result next time when fi will be equal to aStartIndx
3211 SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3212 if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3215 aMapOfTreated.Add( aFacePtr->GetID() );
3216 TColStd_MapOfInteger aResFaces;
3217 if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3218 aMapOfNonManifold, aResFaces ) )
3220 TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3221 for ( ; anItr.More(); anItr.Next() )
3223 int aFaceId = anItr.Key();
3224 aMapOfTreated.Add( aFaceId );
3225 myMapIds.Add( aFaceId );
3228 if ( fi == ( myAllFacePtr.size() - 1 ) )
3230 } // end run on vector of faces
3231 return !myMapIds.IsEmpty();
3234 static void getLinks( const SMDS_MeshFace* theFace,
3235 ManifoldPart::TVectorOfLink& theLinks )
3237 int aNbNode = theFace->NbNodes();
3238 SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3240 SMDS_MeshNode* aNode = 0;
3241 for ( ; aNodeItr->more() && i <= aNbNode; )
3244 SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3248 SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3250 ManifoldPart::Link aLink( aN1, aN2 );
3251 theLinks.push_back( aLink );
3255 bool ManifoldPart::findConnected
3256 ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3257 SMDS_MeshFace* theStartFace,
3258 ManifoldPart::TMapOfLink& theNonManifold,
3259 TColStd_MapOfInteger& theResFaces )
3261 theResFaces.Clear();
3262 if ( !theAllFacePtrInt.size() )
3265 if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3267 myMapBadGeomIds.Add( theStartFace->GetID() );
3271 ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3272 ManifoldPart::TVectorOfLink aSeqOfBoundary;
3273 theResFaces.Add( theStartFace->GetID() );
3274 ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3276 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3277 aDMapLinkFace, theNonManifold, theStartFace );
3279 bool isDone = false;
3280 while ( !isDone && aMapOfBoundary.size() != 0 )
3282 bool isToReset = false;
3283 ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3284 for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3286 ManifoldPart::Link aLink = *pLink;
3287 if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3289 // each link could be treated only once
3290 aMapToSkip.insert( aLink );
3292 ManifoldPart::TVectorOfFacePtr aFaces;
3294 if ( myIsOnlyManifold &&
3295 (theNonManifold.find( aLink ) != theNonManifold.end()) )
3299 getFacesByLink( aLink, aFaces );
3300 // filter the element to keep only indicated elements
3301 ManifoldPart::TVectorOfFacePtr aFiltered;
3302 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3303 for ( ; pFace != aFaces.end(); ++pFace )
3305 SMDS_MeshFace* aFace = *pFace;
3306 if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3307 aFiltered.push_back( aFace );
3310 if ( aFaces.size() < 2 ) // no neihgbour faces
3312 else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3314 theNonManifold.insert( aLink );
3319 // compare normal with normals of neighbor element
3320 SMDS_MeshFace* aPrevFace = aDMapLinkFace[ aLink ];
3321 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3322 for ( ; pFace != aFaces.end(); ++pFace )
3324 SMDS_MeshFace* aNextFace = *pFace;
3325 if ( aPrevFace == aNextFace )
3327 int anNextFaceID = aNextFace->GetID();
3328 if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3329 // should not be with non manifold restriction. probably bad topology
3331 // check if face was treated and skipped
3332 if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3333 !isInPlane( aPrevFace, aNextFace ) )
3335 // add new element to connected and extend the boundaries.
3336 theResFaces.Add( anNextFaceID );
3337 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3338 aDMapLinkFace, theNonManifold, aNextFace );
3342 isDone = !isToReset;
3345 return !theResFaces.IsEmpty();
3348 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3349 const SMDS_MeshFace* theFace2 )
3351 gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3352 gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3353 if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3355 myMapBadGeomIds.Add( theFace2->GetID() );
3358 if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
3364 void ManifoldPart::expandBoundary
3365 ( ManifoldPart::TMapOfLink& theMapOfBoundary,
3366 ManifoldPart::TVectorOfLink& theSeqOfBoundary,
3367 ManifoldPart::TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
3368 ManifoldPart::TMapOfLink& theNonManifold,
3369 SMDS_MeshFace* theNextFace ) const
3371 ManifoldPart::TVectorOfLink aLinks;
3372 getLinks( theNextFace, aLinks );
3373 int aNbLink = (int)aLinks.size();
3374 for ( int i = 0; i < aNbLink; i++ )
3376 ManifoldPart::Link aLink = aLinks[ i ];
3377 if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3379 if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3381 if ( myIsOnlyManifold )
3383 // remove from boundary
3384 theMapOfBoundary.erase( aLink );
3385 ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3386 for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3388 ManifoldPart::Link aBoundLink = *pLink;
3389 if ( aBoundLink.IsEqual( aLink ) )
3391 theSeqOfBoundary.erase( pLink );
3399 theMapOfBoundary.insert( aLink );
3400 theSeqOfBoundary.push_back( aLink );
3401 theDMapLinkFacePtr[ aLink ] = theNextFace;
3406 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3407 ManifoldPart::TVectorOfFacePtr& theFaces ) const
3409 std::set<SMDS_MeshCell *> aSetOfFaces;
3410 // take all faces that shared first node
3411 SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3412 for ( ; anItr->more(); )
3414 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3417 aSetOfFaces.insert( aFace );
3419 // take all faces that shared second node
3420 anItr = theLink.myNode2->facesIterator();
3421 // find the common part of two sets
3422 for ( ; anItr->more(); )
3424 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3425 if ( aSetOfFaces.count( aFace ) )
3426 theFaces.push_back( aFace );
3435 ElementsOnSurface::ElementsOnSurface()
3439 myType = SMDSAbs_All;
3441 myToler = Precision::Confusion();
3442 myUseBoundaries = false;
3445 ElementsOnSurface::~ElementsOnSurface()
3450 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
3452 if ( myMesh == theMesh )
3458 bool ElementsOnSurface::IsSatisfy( long theElementId )
3460 return myIds.Contains( theElementId );
3463 SMDSAbs_ElementType ElementsOnSurface::GetType() const
3466 void ElementsOnSurface::SetTolerance( const double theToler )
3468 if ( myToler != theToler )
3473 double ElementsOnSurface::GetTolerance() const
3476 void ElementsOnSurface::SetUseBoundaries( bool theUse )
3478 if ( myUseBoundaries != theUse ) {
3479 myUseBoundaries = theUse;
3480 SetSurface( mySurf, myType );
3484 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
3485 const SMDSAbs_ElementType theType )
3490 if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
3492 mySurf = TopoDS::Face( theShape );
3493 BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
3495 u1 = SA.FirstUParameter(),
3496 u2 = SA.LastUParameter(),
3497 v1 = SA.FirstVParameter(),
3498 v2 = SA.LastVParameter();
3499 Handle(Geom_Surface) surf = BRep_Tool::Surface( mySurf );
3500 myProjector.Init( surf, u1,u2, v1,v2 );
3504 void ElementsOnSurface::process()
3507 if ( mySurf.IsNull() )
3513 if ( myType == SMDSAbs_Face || myType == SMDSAbs_All )
3515 myIds.ReSize( myMesh->NbFaces() );
3516 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
3517 for(; anIter->more(); )
3518 process( anIter->next() );
3521 if ( myType == SMDSAbs_Edge || myType == SMDSAbs_All )
3523 myIds.ReSize( myIds.Extent() + myMesh->NbEdges() );
3524 SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
3525 for(; anIter->more(); )
3526 process( anIter->next() );
3529 if ( myType == SMDSAbs_Node )
3531 myIds.ReSize( myMesh->NbNodes() );
3532 SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
3533 for(; anIter->more(); )
3534 process( anIter->next() );
3538 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
3540 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3541 bool isSatisfy = true;
3542 for ( ; aNodeItr->more(); )
3544 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3545 if ( !isOnSurface( aNode ) )
3552 myIds.Add( theElemPtr->GetID() );
3555 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
3557 if ( mySurf.IsNull() )
3560 gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
3561 // double aToler2 = myToler * myToler;
3562 // if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
3564 // gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
3565 // if ( aPln.SquareDistance( aPnt ) > aToler2 )
3568 // else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
3570 // gp_Cylinder aCyl = Handle(Geom_CylindricalSurface)::DownCast(mySurf)->Cylinder();
3571 // double aRad = aCyl.Radius();
3572 // gp_Ax3 anAxis = aCyl.Position();
3573 // gp_XYZ aLoc = aCyl.Location().XYZ();
3574 // double aXDist = anAxis.XDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3575 // double aYDist = anAxis.YDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3576 // if ( fabs(aXDist*aXDist + aYDist*aYDist - aRad*aRad) > aToler2 )
3581 myProjector.Perform( aPnt );
3582 bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
3592 ElementsOnShape::ElementsOnShape()
3594 myType(SMDSAbs_All),
3595 myToler(Precision::Confusion()),
3596 myAllNodesFlag(false)
3598 myCurShapeType = TopAbs_SHAPE;
3601 ElementsOnShape::~ElementsOnShape()
3605 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
3607 if (myMesh != theMesh) {
3609 SetShape(myShape, myType);
3613 bool ElementsOnShape::IsSatisfy (long theElementId)
3615 return myIds.Contains(theElementId);
3618 SMDSAbs_ElementType ElementsOnShape::GetType() const
3623 void ElementsOnShape::SetTolerance (const double theToler)
3625 if (myToler != theToler) {
3627 SetShape(myShape, myType);
3631 double ElementsOnShape::GetTolerance() const
3636 void ElementsOnShape::SetAllNodes (bool theAllNodes)
3638 if (myAllNodesFlag != theAllNodes) {
3639 myAllNodesFlag = theAllNodes;
3640 SetShape(myShape, myType);
3644 void ElementsOnShape::SetShape (const TopoDS_Shape& theShape,
3645 const SMDSAbs_ElementType theType)
3651 if (myMesh == 0) return;
3656 myIds.ReSize(myMesh->NbEdges() + myMesh->NbFaces() + myMesh->NbVolumes());
3659 myIds.ReSize(myMesh->NbNodes());
3662 myIds.ReSize(myMesh->NbEdges());
3665 myIds.ReSize(myMesh->NbFaces());
3667 case SMDSAbs_Volume:
3668 myIds.ReSize(myMesh->NbVolumes());
3674 myShapesMap.Clear();
3678 void ElementsOnShape::addShape (const TopoDS_Shape& theShape)
3680 if (theShape.IsNull() || myMesh == 0)
3683 if (!myShapesMap.Add(theShape)) return;
3685 myCurShapeType = theShape.ShapeType();
3686 switch (myCurShapeType)
3688 case TopAbs_COMPOUND:
3689 case TopAbs_COMPSOLID:
3693 TopoDS_Iterator anIt (theShape, Standard_True, Standard_True);
3694 for (; anIt.More(); anIt.Next()) addShape(anIt.Value());
3699 myCurSC.Load(theShape);
3705 TopoDS_Face aFace = TopoDS::Face(theShape);
3706 BRepAdaptor_Surface SA (aFace, true);
3708 u1 = SA.FirstUParameter(),
3709 u2 = SA.LastUParameter(),
3710 v1 = SA.FirstVParameter(),
3711 v2 = SA.LastVParameter();
3712 Handle(Geom_Surface) surf = BRep_Tool::Surface(aFace);
3713 myCurProjFace.Init(surf, u1,u2, v1,v2);
3720 TopoDS_Edge anEdge = TopoDS::Edge(theShape);
3721 Standard_Real u1, u2;
3722 Handle(Geom_Curve) curve = BRep_Tool::Curve(anEdge, u1, u2);
3723 myCurProjEdge.Init(curve, u1, u2);
3729 TopoDS_Vertex aV = TopoDS::Vertex(theShape);
3730 myCurPnt = BRep_Tool::Pnt(aV);
3739 void ElementsOnShape::process()
3741 if (myShape.IsNull() || myMesh == 0)
3744 if (myType == SMDSAbs_Node)
3746 SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
3747 while (anIter->more())
3748 process(anIter->next());
3752 if (myType == SMDSAbs_Edge || myType == SMDSAbs_All)
3754 SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
3755 while (anIter->more())
3756 process(anIter->next());
3759 if (myType == SMDSAbs_Face || myType == SMDSAbs_All)
3761 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
3762 while (anIter->more()) {
3763 process(anIter->next());
3767 if (myType == SMDSAbs_Volume || myType == SMDSAbs_All)
3769 SMDS_VolumeIteratorPtr anIter = myMesh->volumesIterator();
3770 while (anIter->more())
3771 process(anIter->next());
3776 void ElementsOnShape::process (const SMDS_MeshElement* theElemPtr)
3778 if (myShape.IsNull())
3781 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3782 bool isSatisfy = myAllNodesFlag;
3784 gp_XYZ centerXYZ (0, 0, 0);
3786 while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
3788 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3789 gp_Pnt aPnt (aNode->X(), aNode->Y(), aNode->Z());
3790 centerXYZ += aPnt.XYZ();
3792 switch (myCurShapeType)
3796 myCurSC.Perform(aPnt, myToler);
3797 isSatisfy = (myCurSC.State() == TopAbs_IN || myCurSC.State() == TopAbs_ON);
3802 myCurProjFace.Perform(aPnt);
3803 isSatisfy = (myCurProjFace.IsDone() && myCurProjFace.LowerDistance() <= myToler);
3806 // check relatively the face
3807 Quantity_Parameter u, v;
3808 myCurProjFace.LowerDistanceParameters(u, v);
3809 gp_Pnt2d aProjPnt (u, v);
3810 BRepClass_FaceClassifier aClsf (myCurFace, aProjPnt, myToler);
3811 isSatisfy = (aClsf.State() == TopAbs_IN || aClsf.State() == TopAbs_ON);
3817 myCurProjEdge.Perform(aPnt);
3818 isSatisfy = (myCurProjEdge.NbPoints() > 0 && myCurProjEdge.LowerDistance() <= myToler);
3823 isSatisfy = (aPnt.Distance(myCurPnt) <= myToler);
3833 if (isSatisfy && myCurShapeType == TopAbs_SOLID) { // Check the center point for volumes MantisBug 0020168
3834 centerXYZ /= theElemPtr->NbNodes();
3835 gp_Pnt aCenterPnt (centerXYZ);
3836 myCurSC.Perform(aCenterPnt, myToler);
3837 if ( !(myCurSC.State() == TopAbs_IN || myCurSC.State() == TopAbs_ON))
3842 myIds.Add(theElemPtr->GetID());
3845 TSequenceOfXYZ::TSequenceOfXYZ()
3848 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n)
3851 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t)
3854 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray)
3857 template <class InputIterator>
3858 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd)
3861 TSequenceOfXYZ::~TSequenceOfXYZ()
3864 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
3866 myArray = theSequenceOfXYZ.myArray;
3870 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
3872 return myArray[n-1];
3875 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
3877 return myArray[n-1];
3880 void TSequenceOfXYZ::clear()
3885 void TSequenceOfXYZ::reserve(size_type n)
3890 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
3892 myArray.push_back(v);
3895 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
3897 return myArray.size();