1 // Copyright (C) 2007-2012 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"
25 #include "SMDS_BallElement.hxx"
26 #include "SMDS_Iterator.hxx"
27 #include "SMDS_Mesh.hxx"
28 #include "SMDS_MeshElement.hxx"
29 #include "SMDS_MeshNode.hxx"
30 #include "SMDS_QuadraticEdge.hxx"
31 #include "SMDS_QuadraticFaceOfNodes.hxx"
32 #include "SMDS_VolumeTool.hxx"
33 #include "SMESHDS_GroupBase.hxx"
34 #include "SMESHDS_Mesh.hxx"
35 #include "SMESH_OctreeNode.hxx"
37 #include <BRepAdaptor_Surface.hxx>
38 #include <BRepClass_FaceClassifier.hxx>
39 #include <BRep_Tool.hxx>
40 #include <Geom_CylindricalSurface.hxx>
41 #include <Geom_Plane.hxx>
42 #include <Geom_Surface.hxx>
43 #include <Precision.hxx>
44 #include <TColStd_MapIteratorOfMapOfInteger.hxx>
45 #include <TColStd_MapOfInteger.hxx>
46 #include <TColStd_SequenceOfAsciiString.hxx>
47 #include <TColgp_Array1OfXYZ.hxx>
50 #include <TopoDS_Edge.hxx>
51 #include <TopoDS_Face.hxx>
52 #include <TopoDS_Iterator.hxx>
53 #include <TopoDS_Shape.hxx>
54 #include <TopoDS_Vertex.hxx>
56 #include <gp_Cylinder.hxx>
63 #include <vtkMeshQuality.h>
75 const double theEps = 1e-100;
76 const double theInf = 1e+100;
78 inline gp_XYZ gpXYZ(const SMDS_MeshNode* aNode )
80 return gp_XYZ(aNode->X(), aNode->Y(), aNode->Z() );
83 inline double getAngle( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
85 gp_Vec v1( P1 - P2 ), v2( P3 - P2 );
87 return v1.Magnitude() < gp::Resolution() ||
88 v2.Magnitude() < gp::Resolution() ? 0 : v1.Angle( v2 );
91 inline double getArea( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
93 gp_Vec aVec1( P2 - P1 );
94 gp_Vec aVec2( P3 - P1 );
95 return ( aVec1 ^ aVec2 ).Magnitude() * 0.5;
98 inline double getArea( const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3 )
100 return getArea( P1.XYZ(), P2.XYZ(), P3.XYZ() );
105 inline double getDistance( const gp_XYZ& P1, const gp_XYZ& P2 )
107 double aDist = gp_Pnt( P1 ).Distance( gp_Pnt( P2 ) );
111 int getNbMultiConnection( const SMDS_Mesh* theMesh, const int theId )
116 const SMDS_MeshElement* anEdge = theMesh->FindElement( theId );
117 if ( anEdge == 0 || anEdge->GetType() != SMDSAbs_Edge/* || anEdge->NbNodes() != 2 */)
120 // for each pair of nodes in anEdge (there are 2 pairs in a quadratic edge)
121 // count elements containing both nodes of the pair.
122 // Note that there may be such cases for a quadratic edge (a horizontal line):
127 // +-----+------+ +-----+------+
130 // result sould be 2 in both cases
132 int aResult0 = 0, aResult1 = 0;
133 // last node, it is a medium one in a quadratic edge
134 const SMDS_MeshNode* aLastNode = anEdge->GetNode( anEdge->NbNodes() - 1 );
135 const SMDS_MeshNode* aNode0 = anEdge->GetNode( 0 );
136 const SMDS_MeshNode* aNode1 = anEdge->GetNode( 1 );
137 if ( aNode1 == aLastNode ) aNode1 = 0;
139 SMDS_ElemIteratorPtr anElemIter = aLastNode->GetInverseElementIterator();
140 while( anElemIter->more() ) {
141 const SMDS_MeshElement* anElem = anElemIter->next();
142 if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
143 SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
144 while ( anIter->more() ) {
145 if ( const SMDS_MeshElement* anElemNode = anIter->next() ) {
146 if ( anElemNode == aNode0 ) {
148 if ( !aNode1 ) break; // not a quadratic edge
150 else if ( anElemNode == aNode1 )
156 int aResult = std::max ( aResult0, aResult1 );
158 // TColStd_MapOfInteger aMap;
160 // SMDS_ElemIteratorPtr anIter = anEdge->nodesIterator();
161 // if ( anIter != 0 ) {
162 // while( anIter->more() ) {
163 // const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
166 // SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
167 // while( anElemIter->more() ) {
168 // const SMDS_MeshElement* anElem = anElemIter->next();
169 // if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
170 // int anId = anElem->GetID();
172 // if ( anIter->more() ) // i.e. first node
174 // else if ( aMap.Contains( anId ) )
184 gp_XYZ getNormale( const SMDS_MeshFace* theFace, bool* ok=0 )
186 int aNbNode = theFace->NbNodes();
188 gp_XYZ q1 = gpXYZ( theFace->GetNode(1)) - gpXYZ( theFace->GetNode(0));
189 gp_XYZ q2 = gpXYZ( theFace->GetNode(2)) - gpXYZ( theFace->GetNode(0));
192 gp_XYZ q3 = gpXYZ( theFace->GetNode(3)) - gpXYZ( theFace->GetNode(0));
195 double len = n.Modulus();
196 bool zeroLen = ( len <= numeric_limits<double>::min());
200 if (ok) *ok = !zeroLen;
208 using namespace SMESH::Controls;
215 Class : NumericalFunctor
216 Description : Base class for numerical functors
218 NumericalFunctor::NumericalFunctor():
224 void NumericalFunctor::SetMesh( const SMDS_Mesh* theMesh )
229 bool NumericalFunctor::GetPoints(const int theId,
230 TSequenceOfXYZ& theRes ) const
237 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
238 if ( !anElem || anElem->GetType() != this->GetType() )
241 return GetPoints( anElem, theRes );
244 bool NumericalFunctor::GetPoints(const SMDS_MeshElement* anElem,
245 TSequenceOfXYZ& theRes )
252 theRes.reserve( anElem->NbNodes() );
254 // Get nodes of the element
255 SMDS_ElemIteratorPtr anIter;
257 if ( anElem->IsQuadratic() ) {
258 switch ( anElem->GetType() ) {
260 anIter = dynamic_cast<const SMDS_VtkEdge*>
261 (anElem)->interlacedNodesElemIterator();
264 anIter = dynamic_cast<const SMDS_VtkFace*>
265 (anElem)->interlacedNodesElemIterator();
268 anIter = anElem->nodesIterator();
273 anIter = anElem->nodesIterator();
277 while( anIter->more() ) {
278 if ( const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>( anIter->next() ))
279 theRes.push_back( gp_XYZ( aNode->X(), aNode->Y(), aNode->Z() ) );
286 long NumericalFunctor::GetPrecision() const
291 void NumericalFunctor::SetPrecision( const long thePrecision )
293 myPrecision = thePrecision;
294 myPrecisionValue = pow( 10., (double)( myPrecision ) );
297 double NumericalFunctor::GetValue( long theId )
301 myCurrElement = myMesh->FindElement( theId );
304 if ( GetPoints( theId, P ))
305 aVal = Round( GetValue( P ));
310 double NumericalFunctor::Round( const double & aVal )
312 return ( myPrecision >= 0 ) ? floor( aVal * myPrecisionValue + 0.5 ) / myPrecisionValue : aVal;
315 //================================================================================
317 * \brief Return histogram of functor values
318 * \param nbIntervals - number of intervals
319 * \param nbEvents - number of mesh elements having values within i-th interval
320 * \param funValues - boundaries of intervals
321 * \param elements - elements to check vulue of; empty list means "of all"
322 * \param minmax - boundaries of diapason of values to divide into intervals
324 //================================================================================
326 void NumericalFunctor::GetHistogram(int nbIntervals,
327 std::vector<int>& nbEvents,
328 std::vector<double>& funValues,
329 const vector<int>& elements,
330 const double* minmax)
332 if ( nbIntervals < 1 ||
334 !myMesh->GetMeshInfo().NbElements( GetType() ))
336 nbEvents.resize( nbIntervals, 0 );
337 funValues.resize( nbIntervals+1 );
339 // get all values sorted
340 std::multiset< double > values;
341 if ( elements.empty() )
343 SMDS_ElemIteratorPtr elemIt = myMesh->elementsIterator(GetType());
344 while ( elemIt->more() )
345 values.insert( GetValue( elemIt->next()->GetID() ));
349 vector<int>::const_iterator id = elements.begin();
350 for ( ; id != elements.end(); ++id )
351 values.insert( GetValue( *id ));
356 funValues[0] = minmax[0];
357 funValues[nbIntervals] = minmax[1];
361 funValues[0] = *values.begin();
362 funValues[nbIntervals] = *values.rbegin();
364 // case nbIntervals == 1
365 if ( nbIntervals == 1 )
367 nbEvents[0] = values.size();
371 if (funValues.front() == funValues.back())
373 nbEvents.resize( 1 );
374 nbEvents[0] = values.size();
375 funValues[1] = funValues.back();
376 funValues.resize( 2 );
379 std::multiset< double >::iterator min = values.begin(), max;
380 for ( int i = 0; i < nbIntervals; ++i )
382 // find end value of i-th interval
383 double r = (i+1) / double( nbIntervals );
384 funValues[i+1] = funValues.front() * (1-r) + funValues.back() * r;
386 // count values in the i-th interval if there are any
387 if ( min != values.end() && *min <= funValues[i+1] )
389 // find the first value out of the interval
390 max = values.upper_bound( funValues[i+1] ); // max is greater than funValues[i+1], or end()
391 nbEvents[i] = std::distance( min, max );
395 // add values larger than minmax[1]
396 nbEvents.back() += std::distance( min, values.end() );
399 //=======================================================================
400 //function : GetValue
402 //=======================================================================
404 double Volume::GetValue( long theElementId )
406 if ( theElementId && myMesh ) {
407 SMDS_VolumeTool aVolumeTool;
408 if ( aVolumeTool.Set( myMesh->FindElement( theElementId )))
409 return aVolumeTool.GetSize();
414 //=======================================================================
415 //function : GetBadRate
416 //purpose : meaningless as it is not quality control functor
417 //=======================================================================
419 double Volume::GetBadRate( double Value, int /*nbNodes*/ ) const
424 //=======================================================================
427 //=======================================================================
429 SMDSAbs_ElementType Volume::GetType() const
431 return SMDSAbs_Volume;
434 //=======================================================================
436 Class : MaxElementLength2D
437 Description : Functor calculating maximum length of 2D element
439 double MaxElementLength2D::GetValue( const TSequenceOfXYZ& P )
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));
451 else if( len == 4 ) { // quadrangles
452 double L1 = getDistance(P( 1 ),P( 2 ));
453 double L2 = getDistance(P( 2 ),P( 3 ));
454 double L3 = getDistance(P( 3 ),P( 4 ));
455 double L4 = getDistance(P( 4 ),P( 1 ));
456 double D1 = getDistance(P( 1 ),P( 3 ));
457 double D2 = getDistance(P( 2 ),P( 4 ));
458 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
460 else if( len == 6 ) { // quadratic triangles
461 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
462 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
463 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
464 aVal = Max(L1,Max(L2,L3));
466 else if( len == 8 || len == 9 ) { // quadratic quadrangles
467 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
468 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
469 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
470 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
471 double D1 = getDistance(P( 1 ),P( 5 ));
472 double D2 = getDistance(P( 3 ),P( 7 ));
473 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
476 if( myPrecision >= 0 )
478 double prec = pow( 10., (double)myPrecision );
479 aVal = floor( aVal * prec + 0.5 ) / prec;
484 double MaxElementLength2D::GetValue( long theElementId )
487 return GetPoints( theElementId, P ) ? GetValue(P) : 0.0;
490 double MaxElementLength2D::GetBadRate( double Value, int /*nbNodes*/ ) const
495 SMDSAbs_ElementType MaxElementLength2D::GetType() const
500 //=======================================================================
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;
675 //=======================================================================
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 / M_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( long theId )
719 myCurrElement = myMesh->FindElement( theId );
720 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_QUAD )
723 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
724 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
725 aVal = Round( vtkMeshQuality::QuadAspectRatio( avtkCell ));
730 if ( GetPoints( myCurrElement, P ))
731 aVal = Round( GetValue( P ));
736 double AspectRatio::GetValue( const TSequenceOfXYZ& P )
738 // According to "Mesh quality control" by Nadir Bouhamau referring to
739 // Pascal Jean Frey and Paul-Louis George. Maillages, applications aux elements finis.
740 // Hermes Science publications, Paris 1999 ISBN 2-7462-0024-4
743 int nbNodes = P.size();
748 // Compute aspect ratio
750 if ( nbNodes == 3 ) {
751 // Compute lengths of the sides
752 std::vector< double > aLen (nbNodes);
753 for ( int i = 0; i < nbNodes - 1; i++ )
754 aLen[ i ] = getDistance( P( i + 1 ), P( i + 2 ) );
755 aLen[ nbNodes - 1 ] = getDistance( P( 1 ), P( nbNodes ) );
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( 2 ), P( 3 ) );
766 if ( anArea <= theEps )
768 return alfa * maxLen * half_perimeter / anArea;
770 else if ( nbNodes == 6 ) { // quadratic triangles
771 // Compute lengths of the sides
772 std::vector< double > aLen (3);
773 aLen[0] = getDistance( P(1), P(3) );
774 aLen[1] = getDistance( P(3), P(5) );
775 aLen[2] = getDistance( P(5), P(1) );
776 // Q = alfa * h * p / S, where
778 // alfa = sqrt( 3 ) / 6
779 // h - length of the longest edge
780 // p - half perimeter
781 // S - triangle surface
782 const double alfa = sqrt( 3. ) / 6.;
783 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
784 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
785 double anArea = getArea( P(1), P(3), P(5) );
786 if ( anArea <= theEps )
788 return alfa * maxLen * half_perimeter / anArea;
790 else if( nbNodes == 4 ) { // quadrangle
791 // Compute lengths of the sides
792 std::vector< double > aLen (4);
793 aLen[0] = getDistance( P(1), P(2) );
794 aLen[1] = getDistance( P(2), P(3) );
795 aLen[2] = getDistance( P(3), P(4) );
796 aLen[3] = getDistance( P(4), P(1) );
797 // Compute lengths of the diagonals
798 std::vector< double > aDia (2);
799 aDia[0] = getDistance( P(1), P(3) );
800 aDia[1] = getDistance( P(2), P(4) );
801 // Compute areas of all triangles which can be built
802 // taking three nodes of the quadrangle
803 std::vector< double > anArea (4);
804 anArea[0] = getArea( P(1), P(2), P(3) );
805 anArea[1] = getArea( P(1), P(2), P(4) );
806 anArea[2] = getArea( P(1), P(3), P(4) );
807 anArea[3] = getArea( P(2), P(3), P(4) );
808 // Q = alpha * L * C1 / C2, where
810 // alpha = sqrt( 1/32 )
811 // L = max( L1, L2, L3, L4, D1, D2 )
812 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
813 // C2 = min( S1, S2, S3, S4 )
814 // Li - lengths of the edges
815 // Di - lengths of the diagonals
816 // Si - areas of the triangles
817 const double alpha = sqrt( 1 / 32. );
818 double L = Max( aLen[ 0 ],
822 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
823 double C1 = sqrt( ( aLen[0] * aLen[0] +
826 aLen[3] * aLen[3] ) / 4. );
827 double C2 = Min( anArea[ 0 ],
829 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
832 return alpha * L * C1 / C2;
834 else if( nbNodes == 8 || nbNodes == 9 ) { // nbNodes==8 - quadratic quadrangle
835 // Compute lengths of the sides
836 std::vector< double > aLen (4);
837 aLen[0] = getDistance( P(1), P(3) );
838 aLen[1] = getDistance( P(3), P(5) );
839 aLen[2] = getDistance( P(5), P(7) );
840 aLen[3] = getDistance( P(7), P(1) );
841 // Compute lengths of the diagonals
842 std::vector< double > aDia (2);
843 aDia[0] = getDistance( P(1), P(5) );
844 aDia[1] = getDistance( P(3), P(7) );
845 // Compute areas of all triangles which can be built
846 // taking three nodes of the quadrangle
847 std::vector< double > anArea (4);
848 anArea[0] = getArea( P(1), P(3), P(5) );
849 anArea[1] = getArea( P(1), P(3), P(7) );
850 anArea[2] = getArea( P(1), P(5), P(7) );
851 anArea[3] = getArea( P(3), P(5), P(7) );
852 // Q = alpha * L * C1 / C2, where
854 // alpha = sqrt( 1/32 )
855 // L = max( L1, L2, L3, L4, D1, D2 )
856 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
857 // C2 = min( S1, S2, S3, S4 )
858 // Li - lengths of the edges
859 // Di - lengths of the diagonals
860 // Si - areas of the triangles
861 const double alpha = sqrt( 1 / 32. );
862 double L = Max( aLen[ 0 ],
866 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
867 double C1 = sqrt( ( aLen[0] * aLen[0] +
870 aLen[3] * aLen[3] ) / 4. );
871 double C2 = Min( anArea[ 0 ],
873 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
876 return alpha * L * C1 / C2;
881 double AspectRatio::GetBadRate( double Value, int /*nbNodes*/ ) const
883 // the aspect ratio is in the range [1.0,infinity]
884 // < 1.0 = very bad, zero area
887 return ( Value < 0.9 ) ? 1000 : Value / 1000.;
890 SMDSAbs_ElementType AspectRatio::GetType() const
897 Class : AspectRatio3D
898 Description : Functor for calculating aspect ratio
902 inline double getHalfPerimeter(double theTria[3]){
903 return (theTria[0] + theTria[1] + theTria[2])/2.0;
906 inline double getArea(double theHalfPerim, double theTria[3]){
907 return sqrt(theHalfPerim*
908 (theHalfPerim-theTria[0])*
909 (theHalfPerim-theTria[1])*
910 (theHalfPerim-theTria[2]));
913 inline double getVolume(double theLen[6]){
914 double a2 = theLen[0]*theLen[0];
915 double b2 = theLen[1]*theLen[1];
916 double c2 = theLen[2]*theLen[2];
917 double d2 = theLen[3]*theLen[3];
918 double e2 = theLen[4]*theLen[4];
919 double f2 = theLen[5]*theLen[5];
920 double P = 4.0*a2*b2*d2;
921 double Q = a2*(b2+d2-e2)-b2*(a2+d2-f2)-d2*(a2+b2-c2);
922 double R = (b2+d2-e2)*(a2+d2-f2)*(a2+d2-f2);
923 return sqrt(P-Q+R)/12.0;
926 inline double getVolume2(double theLen[6]){
927 double a2 = theLen[0]*theLen[0];
928 double b2 = theLen[1]*theLen[1];
929 double c2 = theLen[2]*theLen[2];
930 double d2 = theLen[3]*theLen[3];
931 double e2 = theLen[4]*theLen[4];
932 double f2 = theLen[5]*theLen[5];
934 double P = a2*e2*(b2+c2+d2+f2-a2-e2);
935 double Q = b2*f2*(a2+c2+d2+e2-b2-f2);
936 double R = c2*d2*(a2+b2+e2+f2-c2-d2);
937 double S = a2*b2*d2+b2*c2*e2+a2*c2*f2+d2*e2*f2;
939 return sqrt(P+Q+R-S)/12.0;
942 inline double getVolume(const TSequenceOfXYZ& P){
943 gp_Vec aVec1( P( 2 ) - P( 1 ) );
944 gp_Vec aVec2( P( 3 ) - P( 1 ) );
945 gp_Vec aVec3( P( 4 ) - P( 1 ) );
946 gp_Vec anAreaVec( aVec1 ^ aVec2 );
947 return fabs(aVec3 * anAreaVec) / 6.0;
950 inline double getMaxHeight(double theLen[6])
952 double aHeight = std::max(theLen[0],theLen[1]);
953 aHeight = std::max(aHeight,theLen[2]);
954 aHeight = std::max(aHeight,theLen[3]);
955 aHeight = std::max(aHeight,theLen[4]);
956 aHeight = std::max(aHeight,theLen[5]);
962 double AspectRatio3D::GetValue( long theId )
965 myCurrElement = myMesh->FindElement( theId );
966 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_TETRA )
968 // Action from CoTech | ACTION 31.3:
969 // EURIWARE BO: Homogenize the formulas used to calculate the Controls in SMESH to fit with
970 // those of ParaView. The library used by ParaView for those calculations can be reused in SMESH.
971 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
972 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
973 aVal = Round( vtkMeshQuality::TetAspectRatio( avtkCell ));
978 if ( GetPoints( myCurrElement, P ))
979 aVal = Round( GetValue( P ));
984 double AspectRatio3D::GetValue( const TSequenceOfXYZ& P )
986 double aQuality = 0.0;
987 if(myCurrElement->IsPoly()) return aQuality;
989 int nbNodes = P.size();
991 if(myCurrElement->IsQuadratic()) {
992 if(nbNodes==10) nbNodes=4; // quadratic tetrahedron
993 else if(nbNodes==13) nbNodes=5; // quadratic pyramid
994 else if(nbNodes==15) nbNodes=6; // quadratic pentahedron
995 else if(nbNodes==20) nbNodes=8; // quadratic hexahedron
996 else if(nbNodes==27) nbNodes=8; // quadratic hexahedron
997 else return aQuality;
1003 getDistance(P( 1 ),P( 2 )), // a
1004 getDistance(P( 2 ),P( 3 )), // b
1005 getDistance(P( 3 ),P( 1 )), // c
1006 getDistance(P( 2 ),P( 4 )), // d
1007 getDistance(P( 3 ),P( 4 )), // e
1008 getDistance(P( 1 ),P( 4 )) // f
1010 double aTria[4][3] = {
1011 {aLen[0],aLen[1],aLen[2]}, // abc
1012 {aLen[0],aLen[3],aLen[5]}, // adf
1013 {aLen[1],aLen[3],aLen[4]}, // bde
1014 {aLen[2],aLen[4],aLen[5]} // cef
1016 double aSumArea = 0.0;
1017 double aHalfPerimeter = getHalfPerimeter(aTria[0]);
1018 double anArea = getArea(aHalfPerimeter,aTria[0]);
1020 aHalfPerimeter = getHalfPerimeter(aTria[1]);
1021 anArea = getArea(aHalfPerimeter,aTria[1]);
1023 aHalfPerimeter = getHalfPerimeter(aTria[2]);
1024 anArea = getArea(aHalfPerimeter,aTria[2]);
1026 aHalfPerimeter = getHalfPerimeter(aTria[3]);
1027 anArea = getArea(aHalfPerimeter,aTria[3]);
1029 double aVolume = getVolume(P);
1030 //double aVolume = getVolume(aLen);
1031 double aHeight = getMaxHeight(aLen);
1032 static double aCoeff = sqrt(2.0)/12.0;
1033 if ( aVolume > DBL_MIN )
1034 aQuality = aCoeff*aHeight*aSumArea/aVolume;
1039 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 3 ),P( 5 )};
1040 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1043 gp_XYZ aXYZ[4] = {P( 1 ),P( 3 ),P( 4 ),P( 5 )};
1044 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1047 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 5 )};
1048 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1051 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 4 ),P( 5 )};
1052 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1058 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 6 )};
1059 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1062 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 3 )};
1063 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1066 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 6 )};
1067 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1070 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1071 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1074 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 6 )};
1075 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1078 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 3 )};
1079 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1085 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1086 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1089 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 4 )};
1090 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1093 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 7 )};
1094 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1097 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 8 )};
1098 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1101 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 3 )};
1102 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1105 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 4 )};
1106 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1109 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 7 )};
1110 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1113 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 8 )};
1114 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1117 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 3 )};
1118 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1121 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 4 )};
1122 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1125 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 7 )};
1126 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1129 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 8 )};
1130 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1133 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 1 )};
1134 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1137 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 2 )};
1138 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1141 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 5 )};
1142 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1145 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 6 )};
1146 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1149 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 1 )};
1150 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1153 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 2 )};
1154 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1157 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 5 )};
1158 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1161 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 6 )};
1162 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1165 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 1 )};
1166 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1169 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1170 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1173 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 5 )};
1174 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1177 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 6 )};
1178 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1181 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1182 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1185 gp_XYZ aXYZ[4] = {P( 4 ),P( 5 ),P( 8 ),P( 2 )};
1186 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1189 gp_XYZ aXYZ[4] = {P( 1 ),P( 4 ),P( 5 ),P( 3 )};
1190 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1193 gp_XYZ aXYZ[4] = {P( 3 ),P( 6 ),P( 7 ),P( 1 )};
1194 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1197 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 6 ),P( 4 )};
1198 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1201 gp_XYZ aXYZ[4] = {P( 5 ),P( 6 ),P( 8 ),P( 3 )};
1202 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1205 gp_XYZ aXYZ[4] = {P( 7 ),P( 8 ),P( 6 ),P( 1 )};
1206 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1209 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 7 )};
1210 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1213 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 2 ),P( 5 )};
1214 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1220 gp_XYZ aXYZ[8] = {P( 1 ),P( 2 ),P( 4 ),P( 5 ),P( 7 ),P( 8 ),P( 10 ),P( 11 )};
1221 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1224 gp_XYZ aXYZ[8] = {P( 2 ),P( 3 ),P( 5 ),P( 6 ),P( 8 ),P( 9 ),P( 11 ),P( 12 )};
1225 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1228 gp_XYZ aXYZ[8] = {P( 3 ),P( 4 ),P( 6 ),P( 1 ),P( 9 ),P( 10 ),P( 12 ),P( 7 )};
1229 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1232 } // switch(nbNodes)
1234 if ( nbNodes > 4 ) {
1235 // avaluate aspect ratio of quadranle faces
1236 AspectRatio aspect2D;
1237 SMDS_VolumeTool::VolumeType type = SMDS_VolumeTool::GetType( nbNodes );
1238 int nbFaces = SMDS_VolumeTool::NbFaces( type );
1239 TSequenceOfXYZ points(4);
1240 for ( int i = 0; i < nbFaces; ++i ) { // loop on faces of a volume
1241 if ( SMDS_VolumeTool::NbFaceNodes( type, i ) != 4 )
1243 const int* pInd = SMDS_VolumeTool::GetFaceNodesIndices( type, i, true );
1244 for ( int p = 0; p < 4; ++p ) // loop on nodes of a quadranle face
1245 points( p + 1 ) = P( pInd[ p ] + 1 );
1246 aQuality = std::max( aQuality, aspect2D.GetValue( points ));
1252 double AspectRatio3D::GetBadRate( double Value, int /*nbNodes*/ ) const
1254 // the aspect ratio is in the range [1.0,infinity]
1257 return Value / 1000.;
1260 SMDSAbs_ElementType AspectRatio3D::GetType() const
1262 return SMDSAbs_Volume;
1268 Description : Functor for calculating warping
1270 double Warping::GetValue( const TSequenceOfXYZ& P )
1272 if ( P.size() != 4 )
1275 gp_XYZ G = ( P( 1 ) + P( 2 ) + P( 3 ) + P( 4 ) ) / 4.;
1277 double A1 = ComputeA( P( 1 ), P( 2 ), P( 3 ), G );
1278 double A2 = ComputeA( P( 2 ), P( 3 ), P( 4 ), G );
1279 double A3 = ComputeA( P( 3 ), P( 4 ), P( 1 ), G );
1280 double A4 = ComputeA( P( 4 ), P( 1 ), P( 2 ), G );
1282 return Max( Max( A1, A2 ), Max( A3, A4 ) );
1285 double Warping::ComputeA( const gp_XYZ& thePnt1,
1286 const gp_XYZ& thePnt2,
1287 const gp_XYZ& thePnt3,
1288 const gp_XYZ& theG ) const
1290 double aLen1 = gp_Pnt( thePnt1 ).Distance( gp_Pnt( thePnt2 ) );
1291 double aLen2 = gp_Pnt( thePnt2 ).Distance( gp_Pnt( thePnt3 ) );
1292 double L = Min( aLen1, aLen2 ) * 0.5;
1296 gp_XYZ GI = ( thePnt2 + thePnt1 ) / 2. - theG;
1297 gp_XYZ GJ = ( thePnt3 + thePnt2 ) / 2. - theG;
1298 gp_XYZ N = GI.Crossed( GJ );
1300 if ( N.Modulus() < gp::Resolution() )
1305 double H = ( thePnt2 - theG ).Dot( N );
1306 return asin( fabs( H / L ) ) * 180. / M_PI;
1309 double Warping::GetBadRate( double Value, int /*nbNodes*/ ) const
1311 // the warp is in the range [0.0,PI/2]
1312 // 0.0 = good (no warp)
1313 // PI/2 = bad (face pliee)
1317 SMDSAbs_ElementType Warping::GetType() const
1319 return SMDSAbs_Face;
1325 Description : Functor for calculating taper
1327 double Taper::GetValue( const TSequenceOfXYZ& P )
1329 if ( P.size() != 4 )
1333 double J1 = getArea( P( 4 ), P( 1 ), P( 2 ) ) / 2.;
1334 double J2 = getArea( P( 3 ), P( 1 ), P( 2 ) ) / 2.;
1335 double J3 = getArea( P( 2 ), P( 3 ), P( 4 ) ) / 2.;
1336 double J4 = getArea( P( 3 ), P( 4 ), P( 1 ) ) / 2.;
1338 double JA = 0.25 * ( J1 + J2 + J3 + J4 );
1342 double T1 = fabs( ( J1 - JA ) / JA );
1343 double T2 = fabs( ( J2 - JA ) / JA );
1344 double T3 = fabs( ( J3 - JA ) / JA );
1345 double T4 = fabs( ( J4 - JA ) / JA );
1347 return Max( Max( T1, T2 ), Max( T3, T4 ) );
1350 double Taper::GetBadRate( double Value, int /*nbNodes*/ ) const
1352 // the taper is in the range [0.0,1.0]
1353 // 0.0 = good (no taper)
1354 // 1.0 = bad (les cotes opposes sont allignes)
1358 SMDSAbs_ElementType Taper::GetType() const
1360 return SMDSAbs_Face;
1366 Description : Functor for calculating skew in degrees
1368 static inline double skewAngle( const gp_XYZ& p1, const gp_XYZ& p2, const gp_XYZ& p3 )
1370 gp_XYZ p12 = ( p2 + p1 ) / 2.;
1371 gp_XYZ p23 = ( p3 + p2 ) / 2.;
1372 gp_XYZ p31 = ( p3 + p1 ) / 2.;
1374 gp_Vec v1( p31 - p2 ), v2( p12 - p23 );
1376 return v1.Magnitude() < gp::Resolution() || v2.Magnitude() < gp::Resolution() ? 0. : v1.Angle( v2 );
1379 double Skew::GetValue( const TSequenceOfXYZ& P )
1381 if ( P.size() != 3 && P.size() != 4 )
1385 static double PI2 = M_PI / 2.;
1386 if ( P.size() == 3 )
1388 double A0 = fabs( PI2 - skewAngle( P( 3 ), P( 1 ), P( 2 ) ) );
1389 double A1 = fabs( PI2 - skewAngle( P( 1 ), P( 2 ), P( 3 ) ) );
1390 double A2 = fabs( PI2 - skewAngle( P( 2 ), P( 3 ), P( 1 ) ) );
1392 return Max( A0, Max( A1, A2 ) ) * 180. / M_PI;
1396 gp_XYZ p12 = ( P( 1 ) + P( 2 ) ) / 2.;
1397 gp_XYZ p23 = ( P( 2 ) + P( 3 ) ) / 2.;
1398 gp_XYZ p34 = ( P( 3 ) + P( 4 ) ) / 2.;
1399 gp_XYZ p41 = ( P( 4 ) + P( 1 ) ) / 2.;
1401 gp_Vec v1( p34 - p12 ), v2( p23 - p41 );
1402 double A = v1.Magnitude() <= gp::Resolution() || v2.Magnitude() <= gp::Resolution()
1403 ? 0. : fabs( PI2 - v1.Angle( v2 ) );
1409 return A * 180. / M_PI;
1413 double Skew::GetBadRate( double Value, int /*nbNodes*/ ) const
1415 // the skew is in the range [0.0,PI/2].
1421 SMDSAbs_ElementType Skew::GetType() const
1423 return SMDSAbs_Face;
1429 Description : Functor for calculating area
1431 double Area::GetValue( const TSequenceOfXYZ& P )
1434 if ( P.size() > 2 ) {
1435 gp_Vec aVec1( P(2) - P(1) );
1436 gp_Vec aVec2( P(3) - P(1) );
1437 gp_Vec SumVec = aVec1 ^ aVec2;
1438 for (int i=4; i<=P.size(); i++) {
1439 gp_Vec aVec1( P(i-1) - P(1) );
1440 gp_Vec aVec2( P(i) - P(1) );
1441 gp_Vec tmp = aVec1 ^ aVec2;
1444 val = SumVec.Magnitude() * 0.5;
1449 double Area::GetBadRate( double Value, int /*nbNodes*/ ) const
1451 // meaningless as it is not a quality control functor
1455 SMDSAbs_ElementType Area::GetType() const
1457 return SMDSAbs_Face;
1463 Description : Functor for calculating length of edge
1465 double Length::GetValue( const TSequenceOfXYZ& P )
1467 switch ( P.size() ) {
1468 case 2: return getDistance( P( 1 ), P( 2 ) );
1469 case 3: return getDistance( P( 1 ), P( 2 ) ) + getDistance( P( 2 ), P( 3 ) );
1474 double Length::GetBadRate( double Value, int /*nbNodes*/ ) const
1476 // meaningless as it is not quality control functor
1480 SMDSAbs_ElementType Length::GetType() const
1482 return SMDSAbs_Edge;
1487 Description : Functor for calculating length of edge
1490 double Length2D::GetValue( long theElementId)
1494 //cout<<"Length2D::GetValue"<<endl;
1495 if (GetPoints(theElementId,P)){
1496 //for(int jj=1; jj<=P.size(); jj++)
1497 // cout<<"jj="<<jj<<" P("<<P(jj).X()<<","<<P(jj).Y()<<","<<P(jj).Z()<<")"<<endl;
1499 double aVal;// = GetValue( P );
1500 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
1501 SMDSAbs_ElementType aType = aElem->GetType();
1510 aVal = getDistance( P( 1 ), P( 2 ) );
1513 else if (len == 3){ // quadratic edge
1514 aVal = getDistance(P( 1 ),P( 3 )) + getDistance(P( 3 ),P( 2 ));
1518 if (len == 3){ // triangles
1519 double L1 = getDistance(P( 1 ),P( 2 ));
1520 double L2 = getDistance(P( 2 ),P( 3 ));
1521 double L3 = getDistance(P( 3 ),P( 1 ));
1522 aVal = Max(L1,Max(L2,L3));
1525 else if (len == 4){ // quadrangles
1526 double L1 = getDistance(P( 1 ),P( 2 ));
1527 double L2 = getDistance(P( 2 ),P( 3 ));
1528 double L3 = getDistance(P( 3 ),P( 4 ));
1529 double L4 = getDistance(P( 4 ),P( 1 ));
1530 aVal = Max(Max(L1,L2),Max(L3,L4));
1533 if (len == 6){ // quadratic triangles
1534 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1535 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1536 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
1537 aVal = Max(L1,Max(L2,L3));
1538 //cout<<"L1="<<L1<<" L2="<<L2<<"L3="<<L3<<" aVal="<<aVal<<endl;
1541 else if (len == 8){ // quadratic quadrangles
1542 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1543 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1544 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
1545 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
1546 aVal = Max(Max(L1,L2),Max(L3,L4));
1549 case SMDSAbs_Volume:
1550 if (len == 4){ // tetraidrs
1551 double L1 = getDistance(P( 1 ),P( 2 ));
1552 double L2 = getDistance(P( 2 ),P( 3 ));
1553 double L3 = getDistance(P( 3 ),P( 1 ));
1554 double L4 = getDistance(P( 1 ),P( 4 ));
1555 double L5 = getDistance(P( 2 ),P( 4 ));
1556 double L6 = getDistance(P( 3 ),P( 4 ));
1557 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1560 else if (len == 5){ // piramids
1561 double L1 = getDistance(P( 1 ),P( 2 ));
1562 double L2 = getDistance(P( 2 ),P( 3 ));
1563 double L3 = getDistance(P( 3 ),P( 4 ));
1564 double L4 = getDistance(P( 4 ),P( 1 ));
1565 double L5 = getDistance(P( 1 ),P( 5 ));
1566 double L6 = getDistance(P( 2 ),P( 5 ));
1567 double L7 = getDistance(P( 3 ),P( 5 ));
1568 double L8 = getDistance(P( 4 ),P( 5 ));
1570 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1571 aVal = Max(aVal,Max(L7,L8));
1574 else if (len == 6){ // pentaidres
1575 double L1 = getDistance(P( 1 ),P( 2 ));
1576 double L2 = getDistance(P( 2 ),P( 3 ));
1577 double L3 = getDistance(P( 3 ),P( 1 ));
1578 double L4 = getDistance(P( 4 ),P( 5 ));
1579 double L5 = getDistance(P( 5 ),P( 6 ));
1580 double L6 = getDistance(P( 6 ),P( 4 ));
1581 double L7 = getDistance(P( 1 ),P( 4 ));
1582 double L8 = getDistance(P( 2 ),P( 5 ));
1583 double L9 = getDistance(P( 3 ),P( 6 ));
1585 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1586 aVal = Max(aVal,Max(Max(L7,L8),L9));
1589 else if (len == 8){ // hexaider
1590 double L1 = getDistance(P( 1 ),P( 2 ));
1591 double L2 = getDistance(P( 2 ),P( 3 ));
1592 double L3 = getDistance(P( 3 ),P( 4 ));
1593 double L4 = getDistance(P( 4 ),P( 1 ));
1594 double L5 = getDistance(P( 5 ),P( 6 ));
1595 double L6 = getDistance(P( 6 ),P( 7 ));
1596 double L7 = getDistance(P( 7 ),P( 8 ));
1597 double L8 = getDistance(P( 8 ),P( 5 ));
1598 double L9 = getDistance(P( 1 ),P( 5 ));
1599 double L10= getDistance(P( 2 ),P( 6 ));
1600 double L11= getDistance(P( 3 ),P( 7 ));
1601 double L12= getDistance(P( 4 ),P( 8 ));
1603 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1604 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1605 aVal = Max(aVal,Max(L11,L12));
1610 if (len == 10){ // quadratic tetraidrs
1611 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
1612 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
1613 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
1614 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1615 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
1616 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
1617 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1620 else if (len == 13){ // quadratic piramids
1621 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
1622 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
1623 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1624 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1625 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1626 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
1627 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
1628 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
1629 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1630 aVal = Max(aVal,Max(L7,L8));
1633 else if (len == 15){ // quadratic pentaidres
1634 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
1635 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
1636 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1637 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1638 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
1639 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
1640 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
1641 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
1642 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
1643 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1644 aVal = Max(aVal,Max(Max(L7,L8),L9));
1647 else if (len == 20){ // quadratic hexaider
1648 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
1649 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
1650 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
1651 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
1652 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
1653 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
1654 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
1655 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
1656 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
1657 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
1658 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
1659 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
1660 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1661 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1662 aVal = Max(aVal,Max(L11,L12));
1674 if ( myPrecision >= 0 )
1676 double prec = pow( 10., (double)( myPrecision ) );
1677 aVal = floor( aVal * prec + 0.5 ) / prec;
1686 double Length2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1688 // meaningless as it is not quality control functor
1692 SMDSAbs_ElementType Length2D::GetType() const
1694 return SMDSAbs_Face;
1697 Length2D::Value::Value(double theLength,long thePntId1, long thePntId2):
1700 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1701 if(thePntId1 > thePntId2){
1702 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1706 bool Length2D::Value::operator<(const Length2D::Value& x) const{
1707 if(myPntId[0] < x.myPntId[0]) return true;
1708 if(myPntId[0] == x.myPntId[0])
1709 if(myPntId[1] < x.myPntId[1]) return true;
1713 void Length2D::GetValues(TValues& theValues){
1715 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1716 for(; anIter->more(); ){
1717 const SMDS_MeshFace* anElem = anIter->next();
1719 if(anElem->IsQuadratic()) {
1720 const SMDS_VtkFace* F =
1721 dynamic_cast<const SMDS_VtkFace*>(anElem);
1722 // use special nodes iterator
1723 SMDS_ElemIteratorPtr anIter = F->interlacedNodesElemIterator();
1728 const SMDS_MeshElement* aNode;
1730 aNode = anIter->next();
1731 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1732 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1733 aNodeId[0] = aNodeId[1] = aNode->GetID();
1736 for(; anIter->more(); ){
1737 const SMDS_MeshNode* N1 = static_cast<const SMDS_MeshNode*> (anIter->next());
1738 P[2] = gp_Pnt(N1->X(),N1->Y(),N1->Z());
1739 aNodeId[2] = N1->GetID();
1740 aLength = P[1].Distance(P[2]);
1741 if(!anIter->more()) break;
1742 const SMDS_MeshNode* N2 = static_cast<const SMDS_MeshNode*> (anIter->next());
1743 P[3] = gp_Pnt(N2->X(),N2->Y(),N2->Z());
1744 aNodeId[3] = N2->GetID();
1745 aLength += P[2].Distance(P[3]);
1746 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1747 Value aValue2(aLength,aNodeId[2],aNodeId[3]);
1749 aNodeId[1] = aNodeId[3];
1750 theValues.insert(aValue1);
1751 theValues.insert(aValue2);
1753 aLength += P[2].Distance(P[0]);
1754 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1755 Value aValue2(aLength,aNodeId[2],aNodeId[0]);
1756 theValues.insert(aValue1);
1757 theValues.insert(aValue2);
1760 SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
1765 const SMDS_MeshElement* aNode;
1766 if(aNodesIter->more()){
1767 aNode = aNodesIter->next();
1768 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1769 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1770 aNodeId[0] = aNodeId[1] = aNode->GetID();
1773 for(; aNodesIter->more(); ){
1774 aNode = aNodesIter->next();
1775 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1776 long anId = aNode->GetID();
1778 P[2] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1780 aLength = P[1].Distance(P[2]);
1782 Value aValue(aLength,aNodeId[1],anId);
1785 theValues.insert(aValue);
1788 aLength = P[0].Distance(P[1]);
1790 Value aValue(aLength,aNodeId[0],aNodeId[1]);
1791 theValues.insert(aValue);
1797 Class : MultiConnection
1798 Description : Functor for calculating number of faces conneted to the edge
1800 double MultiConnection::GetValue( const TSequenceOfXYZ& P )
1804 double MultiConnection::GetValue( long theId )
1806 return getNbMultiConnection( myMesh, theId );
1809 double MultiConnection::GetBadRate( double Value, int /*nbNodes*/ ) const
1811 // meaningless as it is not quality control functor
1815 SMDSAbs_ElementType MultiConnection::GetType() const
1817 return SMDSAbs_Edge;
1821 Class : MultiConnection2D
1822 Description : Functor for calculating number of faces conneted to the edge
1824 double MultiConnection2D::GetValue( const TSequenceOfXYZ& P )
1829 double MultiConnection2D::GetValue( long theElementId )
1833 const SMDS_MeshElement* aFaceElem = myMesh->FindElement(theElementId);
1834 SMDSAbs_ElementType aType = aFaceElem->GetType();
1839 int i = 0, len = aFaceElem->NbNodes();
1840 SMDS_ElemIteratorPtr anIter = aFaceElem->nodesIterator();
1843 const SMDS_MeshNode *aNode, *aNode0;
1844 TColStd_MapOfInteger aMap, aMapPrev;
1846 for (i = 0; i <= len; i++) {
1851 if (anIter->more()) {
1852 aNode = (SMDS_MeshNode*)anIter->next();
1860 if (i == 0) aNode0 = aNode;
1862 SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
1863 while (anElemIter->more()) {
1864 const SMDS_MeshElement* anElem = anElemIter->next();
1865 if (anElem != 0 && anElem->GetType() == SMDSAbs_Face) {
1866 int anId = anElem->GetID();
1869 if (aMapPrev.Contains(anId)) {
1874 aResult = Max(aResult, aNb);
1885 double MultiConnection2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1887 // meaningless as it is not quality control functor
1891 SMDSAbs_ElementType MultiConnection2D::GetType() const
1893 return SMDSAbs_Face;
1896 MultiConnection2D::Value::Value(long thePntId1, long thePntId2)
1898 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1899 if(thePntId1 > thePntId2){
1900 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1904 bool MultiConnection2D::Value::operator<(const MultiConnection2D::Value& x) const{
1905 if(myPntId[0] < x.myPntId[0]) return true;
1906 if(myPntId[0] == x.myPntId[0])
1907 if(myPntId[1] < x.myPntId[1]) return true;
1911 void MultiConnection2D::GetValues(MValues& theValues){
1912 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1913 for(; anIter->more(); ){
1914 const SMDS_MeshFace* anElem = anIter->next();
1915 SMDS_ElemIteratorPtr aNodesIter;
1916 if ( anElem->IsQuadratic() )
1917 aNodesIter = dynamic_cast<const SMDS_VtkFace*>
1918 (anElem)->interlacedNodesElemIterator();
1920 aNodesIter = anElem->nodesIterator();
1923 //int aNbConnects=0;
1924 const SMDS_MeshNode* aNode0;
1925 const SMDS_MeshNode* aNode1;
1926 const SMDS_MeshNode* aNode2;
1927 if(aNodesIter->more()){
1928 aNode0 = (SMDS_MeshNode*) aNodesIter->next();
1930 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode1;
1931 aNodeId[0] = aNodeId[1] = aNodes->GetID();
1933 for(; aNodesIter->more(); ) {
1934 aNode2 = (SMDS_MeshNode*) aNodesIter->next();
1935 long anId = aNode2->GetID();
1938 Value aValue(aNodeId[1],aNodeId[2]);
1939 MValues::iterator aItr = theValues.find(aValue);
1940 if (aItr != theValues.end()){
1945 theValues[aValue] = 1;
1948 //cout << "NodeIds: "<<aNodeId[1]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1949 aNodeId[1] = aNodeId[2];
1952 Value aValue(aNodeId[0],aNodeId[2]);
1953 MValues::iterator aItr = theValues.find(aValue);
1954 if (aItr != theValues.end()) {
1959 theValues[aValue] = 1;
1962 //cout << "NodeIds: "<<aNodeId[0]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1968 Class : BallDiameter
1969 Description : Functor returning diameter of a ball element
1971 double BallDiameter::GetValue( long theId )
1973 double diameter = 0;
1975 if ( const SMDS_BallElement* ball =
1976 dynamic_cast<const SMDS_BallElement*>( myMesh->FindElement( theId )))
1978 diameter = ball->GetDiameter();
1983 double BallDiameter::GetBadRate( double Value, int /*nbNodes*/ ) const
1985 // meaningless as it is not a quality control functor
1989 SMDSAbs_ElementType BallDiameter::GetType() const
1991 return SMDSAbs_Ball;
2000 Class : BadOrientedVolume
2001 Description : Predicate bad oriented volumes
2004 BadOrientedVolume::BadOrientedVolume()
2009 void BadOrientedVolume::SetMesh( const SMDS_Mesh* theMesh )
2014 bool BadOrientedVolume::IsSatisfy( long theId )
2019 SMDS_VolumeTool vTool( myMesh->FindElement( theId ));
2020 return !vTool.IsForward();
2023 SMDSAbs_ElementType BadOrientedVolume::GetType() const
2025 return SMDSAbs_Volume;
2029 Class : BareBorderVolume
2032 bool BareBorderVolume::IsSatisfy(long theElementId )
2034 SMDS_VolumeTool myTool;
2035 if ( myTool.Set( myMesh->FindElement(theElementId)))
2037 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2038 if ( myTool.IsFreeFace( iF ))
2040 const SMDS_MeshNode** n = myTool.GetFaceNodes(iF);
2041 vector< const SMDS_MeshNode*> nodes( n, n+myTool.NbFaceNodes(iF));
2042 if ( !myMesh->FindElement( nodes, SMDSAbs_Face, /*Nomedium=*/false))
2050 Class : BareBorderFace
2053 bool BareBorderFace::IsSatisfy(long theElementId )
2056 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2058 if ( face->GetType() == SMDSAbs_Face )
2060 int nbN = face->NbCornerNodes();
2061 for ( int i = 0; i < nbN && !ok; ++i )
2063 // check if a link is shared by another face
2064 const SMDS_MeshNode* n1 = face->GetNode( i );
2065 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2066 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2067 bool isShared = false;
2068 while ( !isShared && fIt->more() )
2070 const SMDS_MeshElement* f = fIt->next();
2071 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2075 const int iQuad = face->IsQuadratic();
2076 myLinkNodes.resize( 2 + iQuad);
2077 myLinkNodes[0] = n1;
2078 myLinkNodes[1] = n2;
2080 myLinkNodes[2] = face->GetNode( i+nbN );
2081 ok = !myMesh->FindElement( myLinkNodes, SMDSAbs_Edge, /*noMedium=*/false);
2090 Class : OverConstrainedVolume
2093 bool OverConstrainedVolume::IsSatisfy(long theElementId )
2095 // An element is over-constrained if it has N-1 free borders where
2096 // N is the number of edges/faces for a 2D/3D element.
2097 SMDS_VolumeTool myTool;
2098 if ( myTool.Set( myMesh->FindElement(theElementId)))
2100 int nbSharedFaces = 0;
2101 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2102 if ( !myTool.IsFreeFace( iF ) && ++nbSharedFaces > 1 )
2104 return ( nbSharedFaces == 1 );
2110 Class : OverConstrainedFace
2113 bool OverConstrainedFace::IsSatisfy(long theElementId )
2115 // An element is over-constrained if it has N-1 free borders where
2116 // N is the number of edges/faces for a 2D/3D element.
2117 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2118 if ( face->GetType() == SMDSAbs_Face )
2120 int nbSharedBorders = 0;
2121 int nbN = face->NbCornerNodes();
2122 for ( int i = 0; i < nbN; ++i )
2124 // check if a link is shared by another face
2125 const SMDS_MeshNode* n1 = face->GetNode( i );
2126 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2127 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2128 bool isShared = false;
2129 while ( !isShared && fIt->more() )
2131 const SMDS_MeshElement* f = fIt->next();
2132 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2134 if ( isShared && ++nbSharedBorders > 1 )
2137 return ( nbSharedBorders == 1 );
2143 Class : CoincidentNodes
2144 Description : Predicate of Coincident nodes
2147 CoincidentNodes::CoincidentNodes()
2152 bool CoincidentNodes::IsSatisfy( long theElementId )
2154 return myCoincidentIDs.Contains( theElementId );
2157 SMDSAbs_ElementType CoincidentNodes::GetType() const
2159 return SMDSAbs_Node;
2162 void CoincidentNodes::SetMesh( const SMDS_Mesh* theMesh )
2164 myMeshModifTracer.SetMesh( theMesh );
2165 if ( myMeshModifTracer.IsMeshModified() )
2167 TIDSortedNodeSet nodesToCheck;
2168 SMDS_NodeIteratorPtr nIt = theMesh->nodesIterator(/*idInceasingOrder=*/true);
2169 while ( nIt->more() )
2170 nodesToCheck.insert( nodesToCheck.end(), nIt->next() );
2172 list< list< const SMDS_MeshNode*> > nodeGroups;
2173 SMESH_OctreeNode::FindCoincidentNodes ( nodesToCheck, &nodeGroups, myToler );
2175 myCoincidentIDs.Clear();
2176 list< list< const SMDS_MeshNode*> >::iterator groupIt = nodeGroups.begin();
2177 for ( ; groupIt != nodeGroups.end(); ++groupIt )
2179 list< const SMDS_MeshNode*>& coincNodes = *groupIt;
2180 list< const SMDS_MeshNode*>::iterator n = coincNodes.begin();
2181 for ( ; n != coincNodes.end(); ++n )
2182 myCoincidentIDs.Add( (*n)->GetID() );
2188 Class : CoincidentElements
2189 Description : Predicate of Coincident Elements
2190 Note : This class is suitable only for visualization of Coincident Elements
2193 CoincidentElements::CoincidentElements()
2198 void CoincidentElements::SetMesh( const SMDS_Mesh* theMesh )
2203 bool CoincidentElements::IsSatisfy( long theElementId )
2205 if ( !myMesh ) return false;
2207 if ( const SMDS_MeshElement* e = myMesh->FindElement( theElementId ))
2209 if ( e->GetType() != GetType() ) return false;
2210 set< const SMDS_MeshNode* > elemNodes( e->begin_nodes(), e->end_nodes() );
2211 const int nbNodes = e->NbNodes();
2212 SMDS_ElemIteratorPtr invIt = (*elemNodes.begin())->GetInverseElementIterator( GetType() );
2213 while ( invIt->more() )
2215 const SMDS_MeshElement* e2 = invIt->next();
2216 if ( e2 == e || e2->NbNodes() != nbNodes ) continue;
2218 bool sameNodes = true;
2219 for ( size_t i = 0; i < elemNodes.size() && sameNodes; ++i )
2220 sameNodes = ( elemNodes.count( e2->GetNode( i )));
2228 SMDSAbs_ElementType CoincidentElements1D::GetType() const
2230 return SMDSAbs_Edge;
2232 SMDSAbs_ElementType CoincidentElements2D::GetType() const
2234 return SMDSAbs_Face;
2236 SMDSAbs_ElementType CoincidentElements3D::GetType() const
2238 return SMDSAbs_Volume;
2244 Description : Predicate for free borders
2247 FreeBorders::FreeBorders()
2252 void FreeBorders::SetMesh( const SMDS_Mesh* theMesh )
2257 bool FreeBorders::IsSatisfy( long theId )
2259 return getNbMultiConnection( myMesh, theId ) == 1;
2262 SMDSAbs_ElementType FreeBorders::GetType() const
2264 return SMDSAbs_Edge;
2270 Description : Predicate for free Edges
2272 FreeEdges::FreeEdges()
2277 void FreeEdges::SetMesh( const SMDS_Mesh* theMesh )
2282 bool FreeEdges::IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId )
2284 TColStd_MapOfInteger aMap;
2285 for ( int i = 0; i < 2; i++ )
2287 SMDS_ElemIteratorPtr anElemIter = theNodes[ i ]->GetInverseElementIterator(SMDSAbs_Face);
2288 while( anElemIter->more() )
2290 if ( const SMDS_MeshElement* anElem = anElemIter->next())
2292 const int anId = anElem->GetID();
2293 if ( anId != theFaceId && !aMap.Add( anId ))
2301 bool FreeEdges::IsSatisfy( long theId )
2306 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2307 if ( aFace == 0 || aFace->GetType() != SMDSAbs_Face || aFace->NbNodes() < 3 )
2310 SMDS_ElemIteratorPtr anIter;
2311 if ( aFace->IsQuadratic() ) {
2312 anIter = dynamic_cast<const SMDS_VtkFace*>
2313 (aFace)->interlacedNodesElemIterator();
2316 anIter = aFace->nodesIterator();
2321 int i = 0, nbNodes = aFace->NbNodes();
2322 std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
2323 while( anIter->more() )
2325 const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
2328 aNodes[ i++ ] = aNode;
2330 aNodes[ nbNodes ] = aNodes[ 0 ];
2332 for ( i = 0; i < nbNodes; i++ )
2333 if ( IsFreeEdge( &aNodes[ i ], theId ) )
2339 SMDSAbs_ElementType FreeEdges::GetType() const
2341 return SMDSAbs_Face;
2344 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2347 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
2348 if(thePntId1 > thePntId2){
2349 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
2353 bool FreeEdges::Border::operator<(const FreeEdges::Border& x) const{
2354 if(myPntId[0] < x.myPntId[0]) return true;
2355 if(myPntId[0] == x.myPntId[0])
2356 if(myPntId[1] < x.myPntId[1]) return true;
2360 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2361 FreeEdges::TBorders& theRegistry,
2362 FreeEdges::TBorders& theContainer)
2364 if(theRegistry.find(theBorder) == theRegistry.end()){
2365 theRegistry.insert(theBorder);
2366 theContainer.insert(theBorder);
2368 theContainer.erase(theBorder);
2372 void FreeEdges::GetBoreders(TBorders& theBorders)
2375 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2376 for(; anIter->more(); ){
2377 const SMDS_MeshFace* anElem = anIter->next();
2378 long anElemId = anElem->GetID();
2379 SMDS_ElemIteratorPtr aNodesIter;
2380 if ( anElem->IsQuadratic() )
2381 aNodesIter = static_cast<const SMDS_VtkFace*>(anElem)->
2382 interlacedNodesElemIterator();
2384 aNodesIter = anElem->nodesIterator();
2386 const SMDS_MeshElement* aNode;
2387 if(aNodesIter->more()){
2388 aNode = aNodesIter->next();
2389 aNodeId[0] = aNodeId[1] = aNode->GetID();
2391 for(; aNodesIter->more(); ){
2392 aNode = aNodesIter->next();
2393 long anId = aNode->GetID();
2394 Border aBorder(anElemId,aNodeId[1],anId);
2396 UpdateBorders(aBorder,aRegistry,theBorders);
2398 Border aBorder(anElemId,aNodeId[0],aNodeId[1]);
2399 UpdateBorders(aBorder,aRegistry,theBorders);
2406 Description : Predicate for free nodes
2409 FreeNodes::FreeNodes()
2414 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2419 bool FreeNodes::IsSatisfy( long theNodeId )
2421 const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2425 return (aNode->NbInverseElements() < 1);
2428 SMDSAbs_ElementType FreeNodes::GetType() const
2430 return SMDSAbs_Node;
2436 Description : Predicate for free faces
2439 FreeFaces::FreeFaces()
2444 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2449 bool FreeFaces::IsSatisfy( long theId )
2451 if (!myMesh) return false;
2452 // check that faces nodes refers to less than two common volumes
2453 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2454 if ( !aFace || aFace->GetType() != SMDSAbs_Face )
2457 int nbNode = aFace->NbNodes();
2459 // collect volumes check that number of volumss with count equal nbNode not less than 2
2460 typedef map< SMDS_MeshElement*, int > TMapOfVolume; // map of volume counters
2461 typedef map< SMDS_MeshElement*, int >::iterator TItrMapOfVolume; // iterator
2462 TMapOfVolume mapOfVol;
2464 SMDS_ElemIteratorPtr nodeItr = aFace->nodesIterator();
2465 while ( nodeItr->more() ) {
2466 const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>(nodeItr->next());
2467 if ( !aNode ) continue;
2468 SMDS_ElemIteratorPtr volItr = aNode->GetInverseElementIterator(SMDSAbs_Volume);
2469 while ( volItr->more() ) {
2470 SMDS_MeshElement* aVol = (SMDS_MeshElement*)volItr->next();
2471 TItrMapOfVolume itr = mapOfVol.insert(make_pair(aVol, 0)).first;
2476 TItrMapOfVolume volItr = mapOfVol.begin();
2477 TItrMapOfVolume volEnd = mapOfVol.end();
2478 for ( ; volItr != volEnd; ++volItr )
2479 if ( (*volItr).second >= nbNode )
2481 // face is not free if number of volumes constructed on thier nodes more than one
2485 SMDSAbs_ElementType FreeFaces::GetType() const
2487 return SMDSAbs_Face;
2491 Class : LinearOrQuadratic
2492 Description : Predicate to verify whether a mesh element is linear
2495 LinearOrQuadratic::LinearOrQuadratic()
2500 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2505 bool LinearOrQuadratic::IsSatisfy( long theId )
2507 if (!myMesh) return false;
2508 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2509 if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2511 return (!anElem->IsQuadratic());
2514 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2519 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2526 Description : Functor for check color of group to whic mesh element belongs to
2529 GroupColor::GroupColor()
2533 bool GroupColor::IsSatisfy( long theId )
2535 return (myIDs.find( theId ) != myIDs.end());
2538 void GroupColor::SetType( SMDSAbs_ElementType theType )
2543 SMDSAbs_ElementType GroupColor::GetType() const
2548 static bool isEqual( const Quantity_Color& theColor1,
2549 const Quantity_Color& theColor2 )
2551 // tolerance to compare colors
2552 const double tol = 5*1e-3;
2553 return ( fabs( theColor1.Red() - theColor2.Red() ) < tol &&
2554 fabs( theColor1.Green() - theColor2.Green() ) < tol &&
2555 fabs( theColor1.Blue() - theColor2.Blue() ) < tol );
2559 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2563 const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2567 int nbGrp = aMesh->GetNbGroups();
2571 // iterates on groups and find necessary elements ids
2572 const std::set<SMESHDS_GroupBase*>& aGroups = aMesh->GetGroups();
2573 set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
2574 for (; GrIt != aGroups.end(); GrIt++) {
2575 SMESHDS_GroupBase* aGrp = (*GrIt);
2578 // check type and color of group
2579 if ( !isEqual( myColor, aGrp->GetColor() ) )
2581 if ( myType != SMDSAbs_All && myType != (SMDSAbs_ElementType)aGrp->GetType() )
2584 SMDSAbs_ElementType aGrpElType = (SMDSAbs_ElementType)aGrp->GetType();
2585 if ( myType == aGrpElType || (myType == SMDSAbs_All && aGrpElType != SMDSAbs_Node) ) {
2586 // add elements IDS into control
2587 int aSize = aGrp->Extent();
2588 for (int i = 0; i < aSize; i++)
2589 myIDs.insert( aGrp->GetID(i+1) );
2594 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
2596 TCollection_AsciiString aStr = theStr;
2597 aStr.RemoveAll( ' ' );
2598 aStr.RemoveAll( '\t' );
2599 for ( int aPos = aStr.Search( ";;" ); aPos != -1; aPos = aStr.Search( ";;" ) )
2600 aStr.Remove( aPos, 2 );
2601 Standard_Real clr[3];
2602 clr[0] = clr[1] = clr[2] = 0.;
2603 for ( int i = 0; i < 3; i++ ) {
2604 TCollection_AsciiString tmpStr = aStr.Token( ";", i+1 );
2605 if ( !tmpStr.IsEmpty() && tmpStr.IsRealValue() )
2606 clr[i] = tmpStr.RealValue();
2608 myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
2611 //=======================================================================
2612 // name : GetRangeStr
2613 // Purpose : Get range as a string.
2614 // Example: "1,2,3,50-60,63,67,70-"
2615 //=======================================================================
2616 void GroupColor::GetColorStr( TCollection_AsciiString& theResStr ) const
2619 theResStr += TCollection_AsciiString( myColor.Red() );
2620 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2621 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2625 Class : ElemGeomType
2626 Description : Predicate to check element geometry type
2629 ElemGeomType::ElemGeomType()
2632 myType = SMDSAbs_All;
2633 myGeomType = SMDSGeom_TRIANGLE;
2636 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2641 bool ElemGeomType::IsSatisfy( long theId )
2643 if (!myMesh) return false;
2644 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2647 const SMDSAbs_ElementType anElemType = anElem->GetType();
2648 if ( myType != SMDSAbs_All && anElemType != myType )
2650 bool isOk = ( anElem->GetGeomType() == myGeomType );
2654 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2659 SMDSAbs_ElementType ElemGeomType::GetType() const
2664 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2666 myGeomType = theType;
2669 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2674 //================================================================================
2676 * \brief Class CoplanarFaces
2678 //================================================================================
2680 CoplanarFaces::CoplanarFaces()
2681 : myFaceID(0), myToler(0)
2684 void CoplanarFaces::SetMesh( const SMDS_Mesh* theMesh )
2686 myMeshModifTracer.SetMesh( theMesh );
2687 if ( myMeshModifTracer.IsMeshModified() )
2689 // Build a set of coplanar face ids
2691 myCoplanarIDs.clear();
2693 if ( !myMeshModifTracer.GetMesh() || !myFaceID || !myToler )
2696 const SMDS_MeshElement* face = myMeshModifTracer.GetMesh()->FindElement( myFaceID );
2697 if ( !face || face->GetType() != SMDSAbs_Face )
2701 gp_Vec myNorm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2705 const double radianTol = myToler * M_PI / 180.;
2706 std::set< SMESH_TLink > checkedLinks;
2708 std::list< pair< const SMDS_MeshElement*, gp_Vec > > faceQueue;
2709 faceQueue.push_back( make_pair( face, myNorm ));
2710 while ( !faceQueue.empty() )
2712 face = faceQueue.front().first;
2713 myNorm = faceQueue.front().second;
2714 faceQueue.pop_front();
2716 for ( int i = 0, nbN = face->NbCornerNodes(); i < nbN; ++i )
2718 const SMDS_MeshNode* n1 = face->GetNode( i );
2719 const SMDS_MeshNode* n2 = face->GetNode(( i+1 )%nbN);
2720 if ( !checkedLinks.insert( SMESH_TLink( n1, n2 )).second )
2722 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator(SMDSAbs_Face);
2723 while ( fIt->more() )
2725 const SMDS_MeshElement* f = fIt->next();
2726 if ( f->GetNodeIndex( n2 ) > -1 )
2728 gp_Vec norm = getNormale( static_cast<const SMDS_MeshFace*>(f), &normOK );
2729 if (!normOK || myNorm.Angle( norm ) <= radianTol)
2731 myCoplanarIDs.insert( f->GetID() );
2732 faceQueue.push_back( make_pair( f, norm ));
2740 bool CoplanarFaces::IsSatisfy( long theElementId )
2742 return myCoplanarIDs.count( theElementId );
2747 *Description : Predicate for Range of Ids.
2748 * Range may be specified with two ways.
2749 * 1. Using AddToRange method
2750 * 2. With SetRangeStr method. Parameter of this method is a string
2751 * like as "1,2,3,50-60,63,67,70-"
2754 //=======================================================================
2755 // name : RangeOfIds
2756 // Purpose : Constructor
2757 //=======================================================================
2758 RangeOfIds::RangeOfIds()
2761 myType = SMDSAbs_All;
2764 //=======================================================================
2766 // Purpose : Set mesh
2767 //=======================================================================
2768 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
2773 //=======================================================================
2774 // name : AddToRange
2775 // Purpose : Add ID to the range
2776 //=======================================================================
2777 bool RangeOfIds::AddToRange( long theEntityId )
2779 myIds.Add( theEntityId );
2783 //=======================================================================
2784 // name : GetRangeStr
2785 // Purpose : Get range as a string.
2786 // Example: "1,2,3,50-60,63,67,70-"
2787 //=======================================================================
2788 void RangeOfIds::GetRangeStr( TCollection_AsciiString& theResStr )
2792 TColStd_SequenceOfInteger anIntSeq;
2793 TColStd_SequenceOfAsciiString aStrSeq;
2795 TColStd_MapIteratorOfMapOfInteger anIter( myIds );
2796 for ( ; anIter.More(); anIter.Next() )
2798 int anId = anIter.Key();
2799 TCollection_AsciiString aStr( anId );
2800 anIntSeq.Append( anId );
2801 aStrSeq.Append( aStr );
2804 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2806 int aMinId = myMin( i );
2807 int aMaxId = myMax( i );
2809 TCollection_AsciiString aStr;
2810 if ( aMinId != IntegerFirst() )
2815 if ( aMaxId != IntegerLast() )
2818 // find position of the string in result sequence and insert string in it
2819 if ( anIntSeq.Length() == 0 )
2821 anIntSeq.Append( aMinId );
2822 aStrSeq.Append( aStr );
2826 if ( aMinId < anIntSeq.First() )
2828 anIntSeq.Prepend( aMinId );
2829 aStrSeq.Prepend( aStr );
2831 else if ( aMinId > anIntSeq.Last() )
2833 anIntSeq.Append( aMinId );
2834 aStrSeq.Append( aStr );
2837 for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
2838 if ( aMinId < anIntSeq( j ) )
2840 anIntSeq.InsertBefore( j, aMinId );
2841 aStrSeq.InsertBefore( j, aStr );
2847 if ( aStrSeq.Length() == 0 )
2850 theResStr = aStrSeq( 1 );
2851 for ( int j = 2, k = aStrSeq.Length(); j <= k; j++ )
2854 theResStr += aStrSeq( j );
2858 //=======================================================================
2859 // name : SetRangeStr
2860 // Purpose : Define range with string
2861 // Example of entry string: "1,2,3,50-60,63,67,70-"
2862 //=======================================================================
2863 bool RangeOfIds::SetRangeStr( const TCollection_AsciiString& theStr )
2869 TCollection_AsciiString aStr = theStr;
2870 aStr.RemoveAll( ' ' );
2871 aStr.RemoveAll( '\t' );
2873 for ( int aPos = aStr.Search( ",," ); aPos != -1; aPos = aStr.Search( ",," ) )
2874 aStr.Remove( aPos, 2 );
2876 TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
2878 while ( tmpStr != "" )
2880 tmpStr = aStr.Token( ",", i++ );
2881 int aPos = tmpStr.Search( '-' );
2885 if ( tmpStr.IsIntegerValue() )
2886 myIds.Add( tmpStr.IntegerValue() );
2892 TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
2893 TCollection_AsciiString aMinStr = tmpStr;
2895 while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
2896 while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
2898 if ( (!aMinStr.IsEmpty() && !aMinStr.IsIntegerValue()) ||
2899 (!aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue()) )
2902 myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
2903 myMax.Append( aMaxStr.IsEmpty() ? IntegerLast() : aMaxStr.IntegerValue() );
2910 //=======================================================================
2912 // Purpose : Get type of supported entities
2913 //=======================================================================
2914 SMDSAbs_ElementType RangeOfIds::GetType() const
2919 //=======================================================================
2921 // Purpose : Set type of supported entities
2922 //=======================================================================
2923 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
2928 //=======================================================================
2930 // Purpose : Verify whether entity satisfies to this rpedicate
2931 //=======================================================================
2932 bool RangeOfIds::IsSatisfy( long theId )
2937 if ( myType == SMDSAbs_Node )
2939 if ( myMesh->FindNode( theId ) == 0 )
2944 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2945 if ( anElem == 0 || (myType != anElem->GetType() && myType != SMDSAbs_All ))
2949 if ( myIds.Contains( theId ) )
2952 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2953 if ( theId >= myMin( i ) && theId <= myMax( i ) )
2961 Description : Base class for comparators
2963 Comparator::Comparator():
2967 Comparator::~Comparator()
2970 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
2973 myFunctor->SetMesh( theMesh );
2976 void Comparator::SetMargin( double theValue )
2978 myMargin = theValue;
2981 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
2983 myFunctor = theFunct;
2986 SMDSAbs_ElementType Comparator::GetType() const
2988 return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
2991 double Comparator::GetMargin()
2999 Description : Comparator "<"
3001 bool LessThan::IsSatisfy( long theId )
3003 return myFunctor && myFunctor->GetValue( theId ) < myMargin;
3009 Description : Comparator ">"
3011 bool MoreThan::IsSatisfy( long theId )
3013 return myFunctor && myFunctor->GetValue( theId ) > myMargin;
3019 Description : Comparator "="
3022 myToler(Precision::Confusion())
3025 bool EqualTo::IsSatisfy( long theId )
3027 return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
3030 void EqualTo::SetTolerance( double theToler )
3035 double EqualTo::GetTolerance()
3042 Description : Logical NOT predicate
3044 LogicalNOT::LogicalNOT()
3047 LogicalNOT::~LogicalNOT()
3050 bool LogicalNOT::IsSatisfy( long theId )
3052 return myPredicate && !myPredicate->IsSatisfy( theId );
3055 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
3058 myPredicate->SetMesh( theMesh );
3061 void LogicalNOT::SetPredicate( PredicatePtr thePred )
3063 myPredicate = thePred;
3066 SMDSAbs_ElementType LogicalNOT::GetType() const
3068 return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
3073 Class : LogicalBinary
3074 Description : Base class for binary logical predicate
3076 LogicalBinary::LogicalBinary()
3079 LogicalBinary::~LogicalBinary()
3082 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
3085 myPredicate1->SetMesh( theMesh );
3088 myPredicate2->SetMesh( theMesh );
3091 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
3093 myPredicate1 = thePredicate;
3096 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
3098 myPredicate2 = thePredicate;
3101 SMDSAbs_ElementType LogicalBinary::GetType() const
3103 if ( !myPredicate1 || !myPredicate2 )
3106 SMDSAbs_ElementType aType1 = myPredicate1->GetType();
3107 SMDSAbs_ElementType aType2 = myPredicate2->GetType();
3109 return aType1 == aType2 ? aType1 : SMDSAbs_All;
3115 Description : Logical AND
3117 bool LogicalAND::IsSatisfy( long theId )
3122 myPredicate1->IsSatisfy( theId ) &&
3123 myPredicate2->IsSatisfy( theId );
3129 Description : Logical OR
3131 bool LogicalOR::IsSatisfy( long theId )
3136 (myPredicate1->IsSatisfy( theId ) ||
3137 myPredicate2->IsSatisfy( theId ));
3151 void Filter::SetPredicate( PredicatePtr thePredicate )
3153 myPredicate = thePredicate;
3156 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3157 PredicatePtr thePredicate,
3158 TIdSequence& theSequence )
3160 theSequence.clear();
3162 if ( !theMesh || !thePredicate )
3165 thePredicate->SetMesh( theMesh );
3167 SMDS_ElemIteratorPtr elemIt = theMesh->elementsIterator( thePredicate->GetType() );
3169 while ( elemIt->more() ) {
3170 const SMDS_MeshElement* anElem = elemIt->next();
3171 long anId = anElem->GetID();
3172 if ( thePredicate->IsSatisfy( anId ) )
3173 theSequence.push_back( anId );
3178 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3179 Filter::TIdSequence& theSequence )
3181 GetElementsId(theMesh,myPredicate,theSequence);
3188 typedef std::set<SMDS_MeshFace*> TMapOfFacePtr;
3194 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
3195 SMDS_MeshNode* theNode2 )
3201 ManifoldPart::Link::~Link()
3207 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
3209 if ( myNode1 == theLink.myNode1 &&
3210 myNode2 == theLink.myNode2 )
3212 else if ( myNode1 == theLink.myNode2 &&
3213 myNode2 == theLink.myNode1 )
3219 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
3221 if(myNode1 < x.myNode1) return true;
3222 if(myNode1 == x.myNode1)
3223 if(myNode2 < x.myNode2) return true;
3227 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
3228 const ManifoldPart::Link& theLink2 )
3230 return theLink1.IsEqual( theLink2 );
3233 ManifoldPart::ManifoldPart()
3236 myAngToler = Precision::Angular();
3237 myIsOnlyManifold = true;
3240 ManifoldPart::~ManifoldPart()
3245 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
3251 SMDSAbs_ElementType ManifoldPart::GetType() const
3252 { return SMDSAbs_Face; }
3254 bool ManifoldPart::IsSatisfy( long theElementId )
3256 return myMapIds.Contains( theElementId );
3259 void ManifoldPart::SetAngleTolerance( const double theAngToler )
3260 { myAngToler = theAngToler; }
3262 double ManifoldPart::GetAngleTolerance() const
3263 { return myAngToler; }
3265 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
3266 { myIsOnlyManifold = theIsOnly; }
3268 void ManifoldPart::SetStartElem( const long theStartId )
3269 { myStartElemId = theStartId; }
3271 bool ManifoldPart::process()
3274 myMapBadGeomIds.Clear();
3276 myAllFacePtr.clear();
3277 myAllFacePtrIntDMap.clear();
3281 // collect all faces into own map
3282 SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3283 for (; anFaceItr->more(); )
3285 SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3286 myAllFacePtr.push_back( aFacePtr );
3287 myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3290 SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3294 // the map of non manifold links and bad geometry
3295 TMapOfLink aMapOfNonManifold;
3296 TColStd_MapOfInteger aMapOfTreated;
3298 // begin cycle on faces from start index and run on vector till the end
3299 // and from begin to start index to cover whole vector
3300 const int aStartIndx = myAllFacePtrIntDMap[aStartFace];
3301 bool isStartTreat = false;
3302 for ( int fi = aStartIndx; !isStartTreat || fi != aStartIndx ; fi++ )
3304 if ( fi == aStartIndx )
3305 isStartTreat = true;
3306 // as result next time when fi will be equal to aStartIndx
3308 SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3309 if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3312 aMapOfTreated.Add( aFacePtr->GetID() );
3313 TColStd_MapOfInteger aResFaces;
3314 if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3315 aMapOfNonManifold, aResFaces ) )
3317 TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3318 for ( ; anItr.More(); anItr.Next() )
3320 int aFaceId = anItr.Key();
3321 aMapOfTreated.Add( aFaceId );
3322 myMapIds.Add( aFaceId );
3325 if ( fi == ( myAllFacePtr.size() - 1 ) )
3327 } // end run on vector of faces
3328 return !myMapIds.IsEmpty();
3331 static void getLinks( const SMDS_MeshFace* theFace,
3332 ManifoldPart::TVectorOfLink& theLinks )
3334 int aNbNode = theFace->NbNodes();
3335 SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3337 SMDS_MeshNode* aNode = 0;
3338 for ( ; aNodeItr->more() && i <= aNbNode; )
3341 SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3345 SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3347 ManifoldPart::Link aLink( aN1, aN2 );
3348 theLinks.push_back( aLink );
3352 bool ManifoldPart::findConnected
3353 ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3354 SMDS_MeshFace* theStartFace,
3355 ManifoldPart::TMapOfLink& theNonManifold,
3356 TColStd_MapOfInteger& theResFaces )
3358 theResFaces.Clear();
3359 if ( !theAllFacePtrInt.size() )
3362 if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3364 myMapBadGeomIds.Add( theStartFace->GetID() );
3368 ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3369 ManifoldPart::TVectorOfLink aSeqOfBoundary;
3370 theResFaces.Add( theStartFace->GetID() );
3371 ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3373 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3374 aDMapLinkFace, theNonManifold, theStartFace );
3376 bool isDone = false;
3377 while ( !isDone && aMapOfBoundary.size() != 0 )
3379 bool isToReset = false;
3380 ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3381 for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3383 ManifoldPart::Link aLink = *pLink;
3384 if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3386 // each link could be treated only once
3387 aMapToSkip.insert( aLink );
3389 ManifoldPart::TVectorOfFacePtr aFaces;
3391 if ( myIsOnlyManifold &&
3392 (theNonManifold.find( aLink ) != theNonManifold.end()) )
3396 getFacesByLink( aLink, aFaces );
3397 // filter the element to keep only indicated elements
3398 ManifoldPart::TVectorOfFacePtr aFiltered;
3399 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3400 for ( ; pFace != aFaces.end(); ++pFace )
3402 SMDS_MeshFace* aFace = *pFace;
3403 if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3404 aFiltered.push_back( aFace );
3407 if ( aFaces.size() < 2 ) // no neihgbour faces
3409 else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3411 theNonManifold.insert( aLink );
3416 // compare normal with normals of neighbor element
3417 SMDS_MeshFace* aPrevFace = aDMapLinkFace[ aLink ];
3418 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3419 for ( ; pFace != aFaces.end(); ++pFace )
3421 SMDS_MeshFace* aNextFace = *pFace;
3422 if ( aPrevFace == aNextFace )
3424 int anNextFaceID = aNextFace->GetID();
3425 if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3426 // should not be with non manifold restriction. probably bad topology
3428 // check if face was treated and skipped
3429 if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3430 !isInPlane( aPrevFace, aNextFace ) )
3432 // add new element to connected and extend the boundaries.
3433 theResFaces.Add( anNextFaceID );
3434 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3435 aDMapLinkFace, theNonManifold, aNextFace );
3439 isDone = !isToReset;
3442 return !theResFaces.IsEmpty();
3445 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3446 const SMDS_MeshFace* theFace2 )
3448 gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3449 gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3450 if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3452 myMapBadGeomIds.Add( theFace2->GetID() );
3455 if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
3461 void ManifoldPart::expandBoundary
3462 ( ManifoldPart::TMapOfLink& theMapOfBoundary,
3463 ManifoldPart::TVectorOfLink& theSeqOfBoundary,
3464 ManifoldPart::TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
3465 ManifoldPart::TMapOfLink& theNonManifold,
3466 SMDS_MeshFace* theNextFace ) const
3468 ManifoldPart::TVectorOfLink aLinks;
3469 getLinks( theNextFace, aLinks );
3470 int aNbLink = (int)aLinks.size();
3471 for ( int i = 0; i < aNbLink; i++ )
3473 ManifoldPart::Link aLink = aLinks[ i ];
3474 if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3476 if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3478 if ( myIsOnlyManifold )
3480 // remove from boundary
3481 theMapOfBoundary.erase( aLink );
3482 ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3483 for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3485 ManifoldPart::Link aBoundLink = *pLink;
3486 if ( aBoundLink.IsEqual( aLink ) )
3488 theSeqOfBoundary.erase( pLink );
3496 theMapOfBoundary.insert( aLink );
3497 theSeqOfBoundary.push_back( aLink );
3498 theDMapLinkFacePtr[ aLink ] = theNextFace;
3503 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3504 ManifoldPart::TVectorOfFacePtr& theFaces ) const
3506 std::set<SMDS_MeshCell *> aSetOfFaces;
3507 // take all faces that shared first node
3508 SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3509 for ( ; anItr->more(); )
3511 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3514 aSetOfFaces.insert( aFace );
3516 // take all faces that shared second node
3517 anItr = theLink.myNode2->facesIterator();
3518 // find the common part of two sets
3519 for ( ; anItr->more(); )
3521 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3522 if ( aSetOfFaces.count( aFace ) )
3523 theFaces.push_back( aFace );
3532 ElementsOnSurface::ElementsOnSurface()
3535 myType = SMDSAbs_All;
3537 myToler = Precision::Confusion();
3538 myUseBoundaries = false;
3541 ElementsOnSurface::~ElementsOnSurface()
3545 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
3547 myMeshModifTracer.SetMesh( theMesh );
3548 if ( myMeshModifTracer.IsMeshModified())
3552 bool ElementsOnSurface::IsSatisfy( long theElementId )
3554 return myIds.Contains( theElementId );
3557 SMDSAbs_ElementType ElementsOnSurface::GetType() const
3560 void ElementsOnSurface::SetTolerance( const double theToler )
3562 if ( myToler != theToler )
3567 double ElementsOnSurface::GetTolerance() const
3570 void ElementsOnSurface::SetUseBoundaries( bool theUse )
3572 if ( myUseBoundaries != theUse ) {
3573 myUseBoundaries = theUse;
3574 SetSurface( mySurf, myType );
3578 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
3579 const SMDSAbs_ElementType theType )
3584 if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
3586 mySurf = TopoDS::Face( theShape );
3587 BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
3589 u1 = SA.FirstUParameter(),
3590 u2 = SA.LastUParameter(),
3591 v1 = SA.FirstVParameter(),
3592 v2 = SA.LastVParameter();
3593 Handle(Geom_Surface) surf = BRep_Tool::Surface( mySurf );
3594 myProjector.Init( surf, u1,u2, v1,v2 );
3598 void ElementsOnSurface::process()
3601 if ( mySurf.IsNull() )
3604 if ( !myMeshModifTracer.GetMesh() )
3607 myIds.ReSize( myMeshModifTracer.GetMesh()->GetMeshInfo().NbElements( myType ));
3609 SMDS_ElemIteratorPtr anIter = myMeshModifTracer.GetMesh()->elementsIterator( myType );
3610 for(; anIter->more(); )
3611 process( anIter->next() );
3614 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
3616 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3617 bool isSatisfy = true;
3618 for ( ; aNodeItr->more(); )
3620 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3621 if ( !isOnSurface( aNode ) )
3628 myIds.Add( theElemPtr->GetID() );
3631 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
3633 if ( mySurf.IsNull() )
3636 gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
3637 // double aToler2 = myToler * myToler;
3638 // if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
3640 // gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
3641 // if ( aPln.SquareDistance( aPnt ) > aToler2 )
3644 // else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
3646 // gp_Cylinder aCyl = Handle(Geom_CylindricalSurface)::DownCast(mySurf)->Cylinder();
3647 // double aRad = aCyl.Radius();
3648 // gp_Ax3 anAxis = aCyl.Position();
3649 // gp_XYZ aLoc = aCyl.Location().XYZ();
3650 // double aXDist = anAxis.XDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3651 // double aYDist = anAxis.YDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3652 // if ( fabs(aXDist*aXDist + aYDist*aYDist - aRad*aRad) > aToler2 )
3657 myProjector.Perform( aPnt );
3658 bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
3668 ElementsOnShape::ElementsOnShape()
3670 myType(SMDSAbs_All),
3671 myToler(Precision::Confusion()),
3672 myAllNodesFlag(false)
3676 ElementsOnShape::~ElementsOnShape()
3681 SMDSAbs_ElementType ElementsOnShape::GetType() const
3686 void ElementsOnShape::SetTolerance (const double theToler)
3688 if (myToler != theToler) {
3690 SetShape(myShape, myType);
3694 double ElementsOnShape::GetTolerance() const
3699 void ElementsOnShape::SetAllNodes (bool theAllNodes)
3701 myAllNodesFlag = theAllNodes;
3704 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
3709 void ElementsOnShape::SetShape (const TopoDS_Shape& theShape,
3710 const SMDSAbs_ElementType theType)
3714 if ( myShape.IsNull() ) return;
3716 TopTools_IndexedMapOfShape shapesMap;
3717 TopAbs_ShapeEnum shapeTypes[4] = { TopAbs_SOLID, TopAbs_FACE, TopAbs_EDGE, TopAbs_VERTEX };
3718 TopExp_Explorer sub;
3719 for ( int i = 0; i < 4; ++i )
3721 if ( shapesMap.IsEmpty() )
3722 for ( sub.Init( myShape, shapeTypes[i] ); sub.More(); sub.Next() )
3723 shapesMap.Add( sub.Current() );
3725 for ( sub.Init( myShape, shapeTypes[i], shapeTypes[i-1] ); sub.More(); sub.Next() )
3726 shapesMap.Add( sub.Current() );
3730 myClassifiers.resize( shapesMap.Extent() );
3731 for ( int i = 0; i < shapesMap.Extent(); ++i )
3732 myClassifiers[ i ] = new TClassifier( shapesMap( i+1 ), myToler );
3735 void ElementsOnShape::clearClassifiers()
3737 for ( size_t i = 0; i < myClassifiers.size(); ++i )
3738 delete myClassifiers[ i ];
3739 myClassifiers.clear();
3742 bool ElementsOnShape::IsSatisfy (long elemId)
3744 const SMDS_MeshElement* elem =
3745 ( myType == SMDSAbs_Node ? myMesh->FindNode( elemId ) : myMesh->FindElement( elemId ));
3746 if ( !elem || myClassifiers.empty() )
3749 for ( size_t i = 0; i < myClassifiers.size(); ++i )
3751 SMDS_ElemIteratorPtr aNodeItr = elem->nodesIterator();
3752 bool isSatisfy = myAllNodesFlag;
3754 gp_XYZ centerXYZ (0, 0, 0);
3756 while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
3758 SMESH_TNodeXYZ aPnt ( aNodeItr->next() );
3760 isSatisfy = ! myClassifiers[i]->IsOut( aPnt );
3763 // Check the center point for volumes MantisBug 0020168
3764 if (isSatisfy && myClassifiers[i]->ShapeType() == TopAbs_SOLID)
3766 centerXYZ /= elem->NbNodes();
3767 isSatisfy = ! myClassifiers[i]->IsOut( centerXYZ );
3776 TopAbs_ShapeEnum ElementsOnShape::TClassifier::ShapeType() const
3778 return myShape.ShapeType();
3781 bool ElementsOnShape::TClassifier::IsOut(const gp_Pnt& p)
3783 return (this->*myIsOutFun)( p );
3786 void ElementsOnShape::TClassifier::Init (const TopoDS_Shape& theShape, double theTol)
3790 switch ( myShape.ShapeType() )
3792 case TopAbs_SOLID: {
3793 mySolidClfr.Load(theShape);
3794 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfSolid;
3798 Standard_Real u1,u2,v1,v2;
3799 Handle(Geom_Surface) surf = BRep_Tool::Surface( TopoDS::Face( theShape ));
3800 surf->Bounds( u1,u2,v1,v2 );
3801 myProjFace.Init(surf, u1,u2, v1,v2, myTol );
3802 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfFace;
3806 Standard_Real u1, u2;
3807 Handle(Geom_Curve) curve = BRep_Tool::Curve( TopoDS::Edge(theShape), u1, u2);
3808 myProjEdge.Init(curve, u1, u2);
3809 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfEdge;
3812 case TopAbs_VERTEX:{
3813 myVertexXYZ = BRep_Tool::Pnt( TopoDS::Vertex( theShape ) );
3814 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfVertex;
3818 throw SALOME_Exception("Programmer error in usage of ElementsOnShape::TClassifier");
3822 bool ElementsOnShape::TClassifier::isOutOfSolid (const gp_Pnt& p)
3824 mySolidClfr.Perform( p, myTol );
3825 return ( mySolidClfr.State() != TopAbs_IN && mySolidClfr.State() != TopAbs_ON );
3828 bool ElementsOnShape::TClassifier::isOutOfFace (const gp_Pnt& p)
3830 myProjFace.Perform( p );
3831 if ( myProjFace.IsDone() && myProjFace.LowerDistance() <= myTol );
3833 // check relatively to the face
3834 Quantity_Parameter u, v;
3835 myProjFace.LowerDistanceParameters(u, v);
3836 gp_Pnt2d aProjPnt (u, v);
3837 BRepClass_FaceClassifier aClsf ( TopoDS::Face( myShape ), aProjPnt, myTol );
3838 if ( aClsf.State() == TopAbs_IN || aClsf.State() == TopAbs_ON )
3844 bool ElementsOnShape::TClassifier::isOutOfEdge (const gp_Pnt& p)
3846 myProjEdge.Perform( p );
3847 return ! ( myProjEdge.NbPoints() > 0 && myProjEdge.LowerDistance() <= myTol );
3850 bool ElementsOnShape::TClassifier::isOutOfVertex(const gp_Pnt& p)
3852 return ( myVertexXYZ.Distance( p ) > myTol );
3856 TSequenceOfXYZ::TSequenceOfXYZ()
3859 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n)
3862 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t)
3865 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray)
3868 template <class InputIterator>
3869 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd)
3872 TSequenceOfXYZ::~TSequenceOfXYZ()
3875 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
3877 myArray = theSequenceOfXYZ.myArray;
3881 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
3883 return myArray[n-1];
3886 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
3888 return myArray[n-1];
3891 void TSequenceOfXYZ::clear()
3896 void TSequenceOfXYZ::reserve(size_type n)
3901 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
3903 myArray.push_back(v);
3906 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
3908 return myArray.size();
3911 TMeshModifTracer::TMeshModifTracer():
3912 myMeshModifTime(0), myMesh(0)
3915 void TMeshModifTracer::SetMesh( const SMDS_Mesh* theMesh )
3917 if ( theMesh != myMesh )
3918 myMeshModifTime = 0;
3921 bool TMeshModifTracer::IsMeshModified()
3923 bool modified = false;
3926 modified = ( myMeshModifTime != myMesh->GetMTime() );
3927 myMeshModifTime = myMesh->GetMTime();