1 // Copyright (C) 2007-2013 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"
36 #include "SMESH_MeshAlgos.hxx"
38 #include <BRepAdaptor_Surface.hxx>
39 #include <BRepClass_FaceClassifier.hxx>
40 #include <BRep_Tool.hxx>
41 #include <Geom_CylindricalSurface.hxx>
42 #include <Geom_Plane.hxx>
43 #include <Geom_Surface.hxx>
44 #include <Precision.hxx>
45 #include <TColStd_MapIteratorOfMapOfInteger.hxx>
46 #include <TColStd_MapOfInteger.hxx>
47 #include <TColStd_SequenceOfAsciiString.hxx>
48 #include <TColgp_Array1OfXYZ.hxx>
51 #include <TopoDS_Edge.hxx>
52 #include <TopoDS_Face.hxx>
53 #include <TopoDS_Iterator.hxx>
54 #include <TopoDS_Shape.hxx>
55 #include <TopoDS_Vertex.hxx>
57 #include <gp_Cylinder.hxx>
64 #include <vtkMeshQuality.h>
69 #include <Basics_Utils.hxx>
77 const double theEps = 1e-100;
78 const double theInf = 1e+100;
80 inline gp_XYZ gpXYZ(const SMDS_MeshNode* aNode )
82 return gp_XYZ(aNode->X(), aNode->Y(), aNode->Z() );
85 inline double getAngle( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
87 gp_Vec v1( P1 - P2 ), v2( P3 - P2 );
89 return v1.Magnitude() < gp::Resolution() ||
90 v2.Magnitude() < gp::Resolution() ? 0 : v1.Angle( v2 );
93 inline double getArea( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
95 gp_Vec aVec1( P2 - P1 );
96 gp_Vec aVec2( P3 - P1 );
97 return ( aVec1 ^ aVec2 ).Magnitude() * 0.5;
100 inline double getArea( const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3 )
102 return getArea( P1.XYZ(), P2.XYZ(), P3.XYZ() );
107 inline double getDistance( const gp_XYZ& P1, const gp_XYZ& P2 )
109 double aDist = gp_Pnt( P1 ).Distance( gp_Pnt( P2 ) );
113 int getNbMultiConnection( const SMDS_Mesh* theMesh, const int theId )
118 const SMDS_MeshElement* anEdge = theMesh->FindElement( theId );
119 if ( anEdge == 0 || anEdge->GetType() != SMDSAbs_Edge/* || anEdge->NbNodes() != 2 */)
122 // for each pair of nodes in anEdge (there are 2 pairs in a quadratic edge)
123 // count elements containing both nodes of the pair.
124 // Note that there may be such cases for a quadratic edge (a horizontal line):
129 // +-----+------+ +-----+------+
132 // result sould be 2 in both cases
134 int aResult0 = 0, aResult1 = 0;
135 // last node, it is a medium one in a quadratic edge
136 const SMDS_MeshNode* aLastNode = anEdge->GetNode( anEdge->NbNodes() - 1 );
137 const SMDS_MeshNode* aNode0 = anEdge->GetNode( 0 );
138 const SMDS_MeshNode* aNode1 = anEdge->GetNode( 1 );
139 if ( aNode1 == aLastNode ) aNode1 = 0;
141 SMDS_ElemIteratorPtr anElemIter = aLastNode->GetInverseElementIterator();
142 while( anElemIter->more() ) {
143 const SMDS_MeshElement* anElem = anElemIter->next();
144 if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
145 SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
146 while ( anIter->more() ) {
147 if ( const SMDS_MeshElement* anElemNode = anIter->next() ) {
148 if ( anElemNode == aNode0 ) {
150 if ( !aNode1 ) break; // not a quadratic edge
152 else if ( anElemNode == aNode1 )
158 int aResult = std::max ( aResult0, aResult1 );
160 // TColStd_MapOfInteger aMap;
162 // SMDS_ElemIteratorPtr anIter = anEdge->nodesIterator();
163 // if ( anIter != 0 ) {
164 // while( anIter->more() ) {
165 // const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
168 // SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
169 // while( anElemIter->more() ) {
170 // const SMDS_MeshElement* anElem = anElemIter->next();
171 // if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
172 // int anId = anElem->GetID();
174 // if ( anIter->more() ) // i.e. first node
176 // else if ( aMap.Contains( anId ) )
186 gp_XYZ getNormale( const SMDS_MeshFace* theFace, bool* ok=0 )
188 int aNbNode = theFace->NbNodes();
190 gp_XYZ q1 = gpXYZ( theFace->GetNode(1)) - gpXYZ( theFace->GetNode(0));
191 gp_XYZ q2 = gpXYZ( theFace->GetNode(2)) - gpXYZ( theFace->GetNode(0));
194 gp_XYZ q3 = gpXYZ( theFace->GetNode(3)) - gpXYZ( theFace->GetNode(0));
197 double len = n.Modulus();
198 bool zeroLen = ( len <= numeric_limits<double>::min());
202 if (ok) *ok = !zeroLen;
210 using namespace SMESH::Controls;
216 //================================================================================
218 Class : NumericalFunctor
219 Description : Base class for numerical functors
221 //================================================================================
223 NumericalFunctor::NumericalFunctor():
229 void NumericalFunctor::SetMesh( const SMDS_Mesh* theMesh )
234 bool NumericalFunctor::GetPoints(const int theId,
235 TSequenceOfXYZ& theRes ) const
242 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
243 if ( !anElem || anElem->GetType() != this->GetType() )
246 return GetPoints( anElem, theRes );
249 bool NumericalFunctor::GetPoints(const SMDS_MeshElement* anElem,
250 TSequenceOfXYZ& theRes )
257 theRes.reserve( anElem->NbNodes() );
259 // Get nodes of the element
260 SMDS_ElemIteratorPtr anIter;
262 if ( anElem->IsQuadratic() ) {
263 switch ( anElem->GetType() ) {
265 anIter = dynamic_cast<const SMDS_VtkEdge*>
266 (anElem)->interlacedNodesElemIterator();
269 anIter = dynamic_cast<const SMDS_VtkFace*>
270 (anElem)->interlacedNodesElemIterator();
273 anIter = anElem->nodesIterator();
278 anIter = anElem->nodesIterator();
282 while( anIter->more() ) {
283 if ( const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>( anIter->next() ))
284 theRes.push_back( gp_XYZ( aNode->X(), aNode->Y(), aNode->Z() ) );
291 long NumericalFunctor::GetPrecision() const
296 void NumericalFunctor::SetPrecision( const long thePrecision )
298 myPrecision = thePrecision;
299 myPrecisionValue = pow( 10., (double)( myPrecision ) );
302 double NumericalFunctor::GetValue( long theId )
306 myCurrElement = myMesh->FindElement( theId );
309 if ( GetPoints( theId, P ))
310 aVal = Round( GetValue( P ));
315 double NumericalFunctor::Round( const double & aVal )
317 return ( myPrecision >= 0 ) ? floor( aVal * myPrecisionValue + 0.5 ) / myPrecisionValue : aVal;
320 //================================================================================
322 * \brief Return histogram of functor values
323 * \param nbIntervals - number of intervals
324 * \param nbEvents - number of mesh elements having values within i-th interval
325 * \param funValues - boundaries of intervals
326 * \param elements - elements to check vulue of; empty list means "of all"
327 * \param minmax - boundaries of diapason of values to divide into intervals
329 //================================================================================
331 void NumericalFunctor::GetHistogram(int nbIntervals,
332 std::vector<int>& nbEvents,
333 std::vector<double>& funValues,
334 const vector<int>& elements,
335 const double* minmax,
336 const bool isLogarithmic)
338 if ( nbIntervals < 1 ||
340 !myMesh->GetMeshInfo().NbElements( GetType() ))
342 nbEvents.resize( nbIntervals, 0 );
343 funValues.resize( nbIntervals+1 );
345 // get all values sorted
346 std::multiset< double > values;
347 if ( elements.empty() )
349 SMDS_ElemIteratorPtr elemIt = myMesh->elementsIterator(GetType());
350 while ( elemIt->more() )
351 values.insert( GetValue( elemIt->next()->GetID() ));
355 vector<int>::const_iterator id = elements.begin();
356 for ( ; id != elements.end(); ++id )
357 values.insert( GetValue( *id ));
362 funValues[0] = minmax[0];
363 funValues[nbIntervals] = minmax[1];
367 funValues[0] = *values.begin();
368 funValues[nbIntervals] = *values.rbegin();
370 // case nbIntervals == 1
371 if ( nbIntervals == 1 )
373 nbEvents[0] = values.size();
377 if (funValues.front() == funValues.back())
379 nbEvents.resize( 1 );
380 nbEvents[0] = values.size();
381 funValues[1] = funValues.back();
382 funValues.resize( 2 );
385 std::multiset< double >::iterator min = values.begin(), max;
386 for ( int i = 0; i < nbIntervals; ++i )
388 // find end value of i-th interval
389 double r = (i+1) / double(nbIntervals);
390 if (isLogarithmic && funValues.front() > 1e-07 && funValues.back() > 1e-07) {
391 double logmin = log10(funValues.front());
392 double lval = logmin + r * (log10(funValues.back()) - logmin);
393 funValues[i+1] = pow(10.0, lval);
396 funValues[i+1] = funValues.front() * (1-r) + funValues.back() * r;
399 // count values in the i-th interval if there are any
400 if ( min != values.end() && *min <= funValues[i+1] )
402 // find the first value out of the interval
403 max = values.upper_bound( funValues[i+1] ); // max is greater than funValues[i+1], or end()
404 nbEvents[i] = std::distance( min, max );
408 // add values larger than minmax[1]
409 nbEvents.back() += std::distance( min, values.end() );
412 //=======================================================================
415 Description : Functor calculating volume of a 3D element
417 //================================================================================
419 double Volume::GetValue( long theElementId )
421 if ( theElementId && myMesh ) {
422 SMDS_VolumeTool aVolumeTool;
423 if ( aVolumeTool.Set( myMesh->FindElement( theElementId )))
424 return aVolumeTool.GetSize();
429 double Volume::GetBadRate( double Value, int /*nbNodes*/ ) const
434 SMDSAbs_ElementType Volume::GetType() const
436 return SMDSAbs_Volume;
439 //=======================================================================
441 Class : MaxElementLength2D
442 Description : Functor calculating maximum length of 2D element
444 //================================================================================
446 double MaxElementLength2D::GetValue( const TSequenceOfXYZ& P )
452 if( len == 3 ) { // triangles
453 double L1 = getDistance(P( 1 ),P( 2 ));
454 double L2 = getDistance(P( 2 ),P( 3 ));
455 double L3 = getDistance(P( 3 ),P( 1 ));
456 aVal = Max(L1,Max(L2,L3));
458 else if( len == 4 ) { // quadrangles
459 double L1 = getDistance(P( 1 ),P( 2 ));
460 double L2 = getDistance(P( 2 ),P( 3 ));
461 double L3 = getDistance(P( 3 ),P( 4 ));
462 double L4 = getDistance(P( 4 ),P( 1 ));
463 double D1 = getDistance(P( 1 ),P( 3 ));
464 double D2 = getDistance(P( 2 ),P( 4 ));
465 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
467 else if( len == 6 ) { // quadratic triangles
468 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
469 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
470 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
471 aVal = Max(L1,Max(L2,L3));
473 else if( len == 8 || len == 9 ) { // quadratic quadrangles
474 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
475 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
476 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
477 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
478 double D1 = getDistance(P( 1 ),P( 5 ));
479 double D2 = getDistance(P( 3 ),P( 7 ));
480 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
483 if( myPrecision >= 0 )
485 double prec = pow( 10., (double)myPrecision );
486 aVal = floor( aVal * prec + 0.5 ) / prec;
491 double MaxElementLength2D::GetValue( long theElementId )
494 return GetPoints( theElementId, P ) ? GetValue(P) : 0.0;
497 double MaxElementLength2D::GetBadRate( double Value, int /*nbNodes*/ ) const
502 SMDSAbs_ElementType MaxElementLength2D::GetType() const
507 //=======================================================================
509 Class : MaxElementLength3D
510 Description : Functor calculating maximum length of 3D element
512 //================================================================================
514 double MaxElementLength3D::GetValue( long theElementId )
517 if( GetPoints( theElementId, P ) ) {
519 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
520 SMDSAbs_ElementType aType = aElem->GetType();
524 if( len == 4 ) { // tetras
525 double L1 = getDistance(P( 1 ),P( 2 ));
526 double L2 = getDistance(P( 2 ),P( 3 ));
527 double L3 = getDistance(P( 3 ),P( 1 ));
528 double L4 = getDistance(P( 1 ),P( 4 ));
529 double L5 = getDistance(P( 2 ),P( 4 ));
530 double L6 = getDistance(P( 3 ),P( 4 ));
531 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
534 else if( len == 5 ) { // pyramids
535 double L1 = getDistance(P( 1 ),P( 2 ));
536 double L2 = getDistance(P( 2 ),P( 3 ));
537 double L3 = getDistance(P( 3 ),P( 4 ));
538 double L4 = getDistance(P( 4 ),P( 1 ));
539 double L5 = getDistance(P( 1 ),P( 5 ));
540 double L6 = getDistance(P( 2 ),P( 5 ));
541 double L7 = getDistance(P( 3 ),P( 5 ));
542 double L8 = getDistance(P( 4 ),P( 5 ));
543 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
544 aVal = Max(aVal,Max(L7,L8));
547 else if( len == 6 ) { // pentas
548 double L1 = getDistance(P( 1 ),P( 2 ));
549 double L2 = getDistance(P( 2 ),P( 3 ));
550 double L3 = getDistance(P( 3 ),P( 1 ));
551 double L4 = getDistance(P( 4 ),P( 5 ));
552 double L5 = getDistance(P( 5 ),P( 6 ));
553 double L6 = getDistance(P( 6 ),P( 4 ));
554 double L7 = getDistance(P( 1 ),P( 4 ));
555 double L8 = getDistance(P( 2 ),P( 5 ));
556 double L9 = getDistance(P( 3 ),P( 6 ));
557 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
558 aVal = Max(aVal,Max(Max(L7,L8),L9));
561 else if( len == 8 ) { // hexas
562 double L1 = getDistance(P( 1 ),P( 2 ));
563 double L2 = getDistance(P( 2 ),P( 3 ));
564 double L3 = getDistance(P( 3 ),P( 4 ));
565 double L4 = getDistance(P( 4 ),P( 1 ));
566 double L5 = getDistance(P( 5 ),P( 6 ));
567 double L6 = getDistance(P( 6 ),P( 7 ));
568 double L7 = getDistance(P( 7 ),P( 8 ));
569 double L8 = getDistance(P( 8 ),P( 5 ));
570 double L9 = getDistance(P( 1 ),P( 5 ));
571 double L10= getDistance(P( 2 ),P( 6 ));
572 double L11= getDistance(P( 3 ),P( 7 ));
573 double L12= getDistance(P( 4 ),P( 8 ));
574 double D1 = getDistance(P( 1 ),P( 7 ));
575 double D2 = getDistance(P( 2 ),P( 8 ));
576 double D3 = getDistance(P( 3 ),P( 5 ));
577 double D4 = getDistance(P( 4 ),P( 6 ));
578 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
579 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
580 aVal = Max(aVal,Max(L11,L12));
581 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
584 else if( len == 12 ) { // hexagonal prism
585 for ( int i1 = 1; i1 < 12; ++i1 )
586 for ( int i2 = i1+1; i1 <= 12; ++i1 )
587 aVal = Max( aVal, getDistance(P( i1 ),P( i2 )));
590 else if( len == 10 ) { // quadratic tetras
591 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
592 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
593 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
594 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
595 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
596 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
597 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
600 else if( len == 13 ) { // quadratic pyramids
601 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
602 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
603 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
604 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
605 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
606 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
607 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
608 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
609 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
610 aVal = Max(aVal,Max(L7,L8));
613 else if( len == 15 ) { // quadratic pentas
614 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
615 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
616 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
617 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
618 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
619 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
620 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
621 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
622 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
623 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
624 aVal = Max(aVal,Max(Max(L7,L8),L9));
627 else if( len == 20 || len == 27 ) { // quadratic hexas
628 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
629 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
630 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
631 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
632 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
633 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
634 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
635 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
636 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
637 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
638 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
639 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
640 double D1 = getDistance(P( 1 ),P( 7 ));
641 double D2 = getDistance(P( 2 ),P( 8 ));
642 double D3 = getDistance(P( 3 ),P( 5 ));
643 double D4 = getDistance(P( 4 ),P( 6 ));
644 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
645 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
646 aVal = Max(aVal,Max(L11,L12));
647 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
650 else if( len > 1 && aElem->IsPoly() ) { // polys
651 // get the maximum distance between all pairs of nodes
652 for( int i = 1; i <= len; i++ ) {
653 for( int j = 1; j <= len; j++ ) {
654 if( j > i ) { // optimization of the loop
655 double D = getDistance( P(i), P(j) );
656 aVal = Max( aVal, D );
663 if( myPrecision >= 0 )
665 double prec = pow( 10., (double)myPrecision );
666 aVal = floor( aVal * prec + 0.5 ) / prec;
673 double MaxElementLength3D::GetBadRate( double Value, int /*nbNodes*/ ) const
678 SMDSAbs_ElementType MaxElementLength3D::GetType() const
680 return SMDSAbs_Volume;
683 //=======================================================================
686 Description : Functor for calculation of minimum angle
688 //================================================================================
690 double MinimumAngle::GetValue( const TSequenceOfXYZ& P )
697 aMin = getAngle(P( P.size() ), P( 1 ), P( 2 ));
698 aMin = Min(aMin,getAngle(P( P.size()-1 ), P( P.size() ), P( 1 )));
700 for (int i=2; i<P.size();i++){
701 double A0 = getAngle( P( i-1 ), P( i ), P( i+1 ) );
705 return aMin * 180.0 / M_PI;
708 double MinimumAngle::GetBadRate( double Value, int nbNodes ) const
710 //const double aBestAngle = PI / nbNodes;
711 const double aBestAngle = 180.0 - ( 360.0 / double(nbNodes) );
712 return ( fabs( aBestAngle - Value ));
715 SMDSAbs_ElementType MinimumAngle::GetType() const
721 //================================================================================
724 Description : Functor for calculating aspect ratio
726 //================================================================================
728 double AspectRatio::GetValue( long theId )
731 myCurrElement = myMesh->FindElement( theId );
732 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_QUAD )
735 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
736 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
737 aVal = Round( vtkMeshQuality::QuadAspectRatio( avtkCell ));
742 if ( GetPoints( myCurrElement, P ))
743 aVal = Round( GetValue( P ));
748 double AspectRatio::GetValue( const TSequenceOfXYZ& P )
750 // According to "Mesh quality control" by Nadir Bouhamau referring to
751 // Pascal Jean Frey and Paul-Louis George. Maillages, applications aux elements finis.
752 // Hermes Science publications, Paris 1999 ISBN 2-7462-0024-4
755 int nbNodes = P.size();
760 // Compute aspect ratio
762 if ( nbNodes == 3 ) {
763 // Compute lengths of the sides
764 std::vector< double > aLen (nbNodes);
765 for ( int i = 0; i < nbNodes - 1; i++ )
766 aLen[ i ] = getDistance( P( i + 1 ), P( i + 2 ) );
767 aLen[ nbNodes - 1 ] = getDistance( P( 1 ), P( nbNodes ) );
768 // Q = alfa * h * p / S, where
770 // alfa = sqrt( 3 ) / 6
771 // h - length of the longest edge
772 // p - half perimeter
773 // S - triangle surface
774 const double alfa = sqrt( 3. ) / 6.;
775 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
776 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
777 double anArea = getArea( P( 1 ), P( 2 ), P( 3 ) );
778 if ( anArea <= theEps )
780 return alfa * maxLen * half_perimeter / anArea;
782 else if ( nbNodes == 6 ) { // quadratic triangles
783 // Compute lengths of the sides
784 std::vector< double > aLen (3);
785 aLen[0] = getDistance( P(1), P(3) );
786 aLen[1] = getDistance( P(3), P(5) );
787 aLen[2] = getDistance( P(5), P(1) );
788 // Q = alfa * h * p / S, where
790 // alfa = sqrt( 3 ) / 6
791 // h - length of the longest edge
792 // p - half perimeter
793 // S - triangle surface
794 const double alfa = sqrt( 3. ) / 6.;
795 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
796 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
797 double anArea = getArea( P(1), P(3), P(5) );
798 if ( anArea <= theEps )
800 return alfa * maxLen * half_perimeter / anArea;
802 else if( nbNodes == 4 ) { // quadrangle
803 // Compute lengths of the sides
804 std::vector< double > aLen (4);
805 aLen[0] = getDistance( P(1), P(2) );
806 aLen[1] = getDistance( P(2), P(3) );
807 aLen[2] = getDistance( P(3), P(4) );
808 aLen[3] = getDistance( P(4), P(1) );
809 // Compute lengths of the diagonals
810 std::vector< double > aDia (2);
811 aDia[0] = getDistance( P(1), P(3) );
812 aDia[1] = getDistance( P(2), P(4) );
813 // Compute areas of all triangles which can be built
814 // taking three nodes of the quadrangle
815 std::vector< double > anArea (4);
816 anArea[0] = getArea( P(1), P(2), P(3) );
817 anArea[1] = getArea( P(1), P(2), P(4) );
818 anArea[2] = getArea( P(1), P(3), P(4) );
819 anArea[3] = getArea( P(2), P(3), P(4) );
820 // Q = alpha * L * C1 / C2, where
822 // alpha = sqrt( 1/32 )
823 // L = max( L1, L2, L3, L4, D1, D2 )
824 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
825 // C2 = min( S1, S2, S3, S4 )
826 // Li - lengths of the edges
827 // Di - lengths of the diagonals
828 // Si - areas of the triangles
829 const double alpha = sqrt( 1 / 32. );
830 double L = Max( aLen[ 0 ],
834 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
835 double C1 = sqrt( ( aLen[0] * aLen[0] +
838 aLen[3] * aLen[3] ) / 4. );
839 double C2 = Min( anArea[ 0 ],
841 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
844 return alpha * L * C1 / C2;
846 else if( nbNodes == 8 || nbNodes == 9 ) { // nbNodes==8 - quadratic quadrangle
847 // Compute lengths of the sides
848 std::vector< double > aLen (4);
849 aLen[0] = getDistance( P(1), P(3) );
850 aLen[1] = getDistance( P(3), P(5) );
851 aLen[2] = getDistance( P(5), P(7) );
852 aLen[3] = getDistance( P(7), P(1) );
853 // Compute lengths of the diagonals
854 std::vector< double > aDia (2);
855 aDia[0] = getDistance( P(1), P(5) );
856 aDia[1] = getDistance( P(3), P(7) );
857 // Compute areas of all triangles which can be built
858 // taking three nodes of the quadrangle
859 std::vector< double > anArea (4);
860 anArea[0] = getArea( P(1), P(3), P(5) );
861 anArea[1] = getArea( P(1), P(3), P(7) );
862 anArea[2] = getArea( P(1), P(5), P(7) );
863 anArea[3] = getArea( P(3), P(5), P(7) );
864 // Q = alpha * L * C1 / C2, where
866 // alpha = sqrt( 1/32 )
867 // L = max( L1, L2, L3, L4, D1, D2 )
868 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
869 // C2 = min( S1, S2, S3, S4 )
870 // Li - lengths of the edges
871 // Di - lengths of the diagonals
872 // Si - areas of the triangles
873 const double alpha = sqrt( 1 / 32. );
874 double L = Max( aLen[ 0 ],
878 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
879 double C1 = sqrt( ( aLen[0] * aLen[0] +
882 aLen[3] * aLen[3] ) / 4. );
883 double C2 = Min( anArea[ 0 ],
885 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
888 return alpha * L * C1 / C2;
893 double AspectRatio::GetBadRate( double Value, int /*nbNodes*/ ) const
895 // the aspect ratio is in the range [1.0,infinity]
896 // < 1.0 = very bad, zero area
899 return ( Value < 0.9 ) ? 1000 : Value / 1000.;
902 SMDSAbs_ElementType AspectRatio::GetType() const
908 //================================================================================
910 Class : AspectRatio3D
911 Description : Functor for calculating aspect ratio
913 //================================================================================
917 inline double getHalfPerimeter(double theTria[3]){
918 return (theTria[0] + theTria[1] + theTria[2])/2.0;
921 inline double getArea(double theHalfPerim, double theTria[3]){
922 return sqrt(theHalfPerim*
923 (theHalfPerim-theTria[0])*
924 (theHalfPerim-theTria[1])*
925 (theHalfPerim-theTria[2]));
928 inline double getVolume(double theLen[6]){
929 double a2 = theLen[0]*theLen[0];
930 double b2 = theLen[1]*theLen[1];
931 double c2 = theLen[2]*theLen[2];
932 double d2 = theLen[3]*theLen[3];
933 double e2 = theLen[4]*theLen[4];
934 double f2 = theLen[5]*theLen[5];
935 double P = 4.0*a2*b2*d2;
936 double Q = a2*(b2+d2-e2)-b2*(a2+d2-f2)-d2*(a2+b2-c2);
937 double R = (b2+d2-e2)*(a2+d2-f2)*(a2+d2-f2);
938 return sqrt(P-Q+R)/12.0;
941 inline double getVolume2(double theLen[6]){
942 double a2 = theLen[0]*theLen[0];
943 double b2 = theLen[1]*theLen[1];
944 double c2 = theLen[2]*theLen[2];
945 double d2 = theLen[3]*theLen[3];
946 double e2 = theLen[4]*theLen[4];
947 double f2 = theLen[5]*theLen[5];
949 double P = a2*e2*(b2+c2+d2+f2-a2-e2);
950 double Q = b2*f2*(a2+c2+d2+e2-b2-f2);
951 double R = c2*d2*(a2+b2+e2+f2-c2-d2);
952 double S = a2*b2*d2+b2*c2*e2+a2*c2*f2+d2*e2*f2;
954 return sqrt(P+Q+R-S)/12.0;
957 inline double getVolume(const TSequenceOfXYZ& P){
958 gp_Vec aVec1( P( 2 ) - P( 1 ) );
959 gp_Vec aVec2( P( 3 ) - P( 1 ) );
960 gp_Vec aVec3( P( 4 ) - P( 1 ) );
961 gp_Vec anAreaVec( aVec1 ^ aVec2 );
962 return fabs(aVec3 * anAreaVec) / 6.0;
965 inline double getMaxHeight(double theLen[6])
967 double aHeight = std::max(theLen[0],theLen[1]);
968 aHeight = std::max(aHeight,theLen[2]);
969 aHeight = std::max(aHeight,theLen[3]);
970 aHeight = std::max(aHeight,theLen[4]);
971 aHeight = std::max(aHeight,theLen[5]);
977 double AspectRatio3D::GetValue( long theId )
980 myCurrElement = myMesh->FindElement( theId );
981 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_TETRA )
983 // Action from CoTech | ACTION 31.3:
984 // EURIWARE BO: Homogenize the formulas used to calculate the Controls in SMESH to fit with
985 // those of ParaView. The library used by ParaView for those calculations can be reused in SMESH.
986 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
987 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
988 aVal = Round( vtkMeshQuality::TetAspectRatio( avtkCell ));
993 if ( GetPoints( myCurrElement, P ))
994 aVal = Round( GetValue( P ));
999 double AspectRatio3D::GetValue( const TSequenceOfXYZ& P )
1001 double aQuality = 0.0;
1002 if(myCurrElement->IsPoly()) return aQuality;
1004 int nbNodes = P.size();
1006 if(myCurrElement->IsQuadratic()) {
1007 if(nbNodes==10) nbNodes=4; // quadratic tetrahedron
1008 else if(nbNodes==13) nbNodes=5; // quadratic pyramid
1009 else if(nbNodes==15) nbNodes=6; // quadratic pentahedron
1010 else if(nbNodes==20) nbNodes=8; // quadratic hexahedron
1011 else if(nbNodes==27) nbNodes=8; // quadratic hexahedron
1012 else return aQuality;
1018 getDistance(P( 1 ),P( 2 )), // a
1019 getDistance(P( 2 ),P( 3 )), // b
1020 getDistance(P( 3 ),P( 1 )), // c
1021 getDistance(P( 2 ),P( 4 )), // d
1022 getDistance(P( 3 ),P( 4 )), // e
1023 getDistance(P( 1 ),P( 4 )) // f
1025 double aTria[4][3] = {
1026 {aLen[0],aLen[1],aLen[2]}, // abc
1027 {aLen[0],aLen[3],aLen[5]}, // adf
1028 {aLen[1],aLen[3],aLen[4]}, // bde
1029 {aLen[2],aLen[4],aLen[5]} // cef
1031 double aSumArea = 0.0;
1032 double aHalfPerimeter = getHalfPerimeter(aTria[0]);
1033 double anArea = getArea(aHalfPerimeter,aTria[0]);
1035 aHalfPerimeter = getHalfPerimeter(aTria[1]);
1036 anArea = getArea(aHalfPerimeter,aTria[1]);
1038 aHalfPerimeter = getHalfPerimeter(aTria[2]);
1039 anArea = getArea(aHalfPerimeter,aTria[2]);
1041 aHalfPerimeter = getHalfPerimeter(aTria[3]);
1042 anArea = getArea(aHalfPerimeter,aTria[3]);
1044 double aVolume = getVolume(P);
1045 //double aVolume = getVolume(aLen);
1046 double aHeight = getMaxHeight(aLen);
1047 static double aCoeff = sqrt(2.0)/12.0;
1048 if ( aVolume > DBL_MIN )
1049 aQuality = aCoeff*aHeight*aSumArea/aVolume;
1054 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 3 ),P( 5 )};
1055 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1058 gp_XYZ aXYZ[4] = {P( 1 ),P( 3 ),P( 4 ),P( 5 )};
1059 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1062 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 5 )};
1063 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1066 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 4 ),P( 5 )};
1067 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1073 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 6 )};
1074 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1077 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 3 )};
1078 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1081 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 6 )};
1082 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( 2 ),P( 5 ),P( 4 ),P( 6 )};
1090 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1093 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 3 )};
1094 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1100 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1101 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1104 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 4 )};
1105 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1108 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 7 )};
1109 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1112 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 8 )};
1113 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1116 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 3 )};
1117 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1120 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 4 )};
1121 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1124 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 7 )};
1125 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1128 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 8 )};
1129 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1132 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 3 )};
1133 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1136 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 4 )};
1137 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1140 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 7 )};
1141 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1144 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 8 )};
1145 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1148 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 1 )};
1149 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1152 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 2 )};
1153 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1156 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 5 )};
1157 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1160 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 6 )};
1161 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1164 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 1 )};
1165 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1168 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 2 )};
1169 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1172 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 5 )};
1173 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1176 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 6 )};
1177 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1180 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 1 )};
1181 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1184 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1185 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1188 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 5 )};
1189 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1192 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 6 )};
1193 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1196 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1197 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1200 gp_XYZ aXYZ[4] = {P( 4 ),P( 5 ),P( 8 ),P( 2 )};
1201 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1204 gp_XYZ aXYZ[4] = {P( 1 ),P( 4 ),P( 5 ),P( 3 )};
1205 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1208 gp_XYZ aXYZ[4] = {P( 3 ),P( 6 ),P( 7 ),P( 1 )};
1209 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1212 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 6 ),P( 4 )};
1213 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1216 gp_XYZ aXYZ[4] = {P( 5 ),P( 6 ),P( 8 ),P( 3 )};
1217 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1220 gp_XYZ aXYZ[4] = {P( 7 ),P( 8 ),P( 6 ),P( 1 )};
1221 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1224 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 7 )};
1225 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1228 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 2 ),P( 5 )};
1229 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1235 gp_XYZ aXYZ[8] = {P( 1 ),P( 2 ),P( 4 ),P( 5 ),P( 7 ),P( 8 ),P( 10 ),P( 11 )};
1236 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1239 gp_XYZ aXYZ[8] = {P( 2 ),P( 3 ),P( 5 ),P( 6 ),P( 8 ),P( 9 ),P( 11 ),P( 12 )};
1240 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1243 gp_XYZ aXYZ[8] = {P( 3 ),P( 4 ),P( 6 ),P( 1 ),P( 9 ),P( 10 ),P( 12 ),P( 7 )};
1244 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1247 } // switch(nbNodes)
1249 if ( nbNodes > 4 ) {
1250 // avaluate aspect ratio of quadranle faces
1251 AspectRatio aspect2D;
1252 SMDS_VolumeTool::VolumeType type = SMDS_VolumeTool::GetType( nbNodes );
1253 int nbFaces = SMDS_VolumeTool::NbFaces( type );
1254 TSequenceOfXYZ points(4);
1255 for ( int i = 0; i < nbFaces; ++i ) { // loop on faces of a volume
1256 if ( SMDS_VolumeTool::NbFaceNodes( type, i ) != 4 )
1258 const int* pInd = SMDS_VolumeTool::GetFaceNodesIndices( type, i, true );
1259 for ( int p = 0; p < 4; ++p ) // loop on nodes of a quadranle face
1260 points( p + 1 ) = P( pInd[ p ] + 1 );
1261 aQuality = std::max( aQuality, aspect2D.GetValue( points ));
1267 double AspectRatio3D::GetBadRate( double Value, int /*nbNodes*/ ) const
1269 // the aspect ratio is in the range [1.0,infinity]
1272 return Value / 1000.;
1275 SMDSAbs_ElementType AspectRatio3D::GetType() const
1277 return SMDSAbs_Volume;
1281 //================================================================================
1284 Description : Functor for calculating warping
1286 //================================================================================
1288 double Warping::GetValue( const TSequenceOfXYZ& P )
1290 if ( P.size() != 4 )
1293 gp_XYZ G = ( P( 1 ) + P( 2 ) + P( 3 ) + P( 4 ) ) / 4.;
1295 double A1 = ComputeA( P( 1 ), P( 2 ), P( 3 ), G );
1296 double A2 = ComputeA( P( 2 ), P( 3 ), P( 4 ), G );
1297 double A3 = ComputeA( P( 3 ), P( 4 ), P( 1 ), G );
1298 double A4 = ComputeA( P( 4 ), P( 1 ), P( 2 ), G );
1300 double val = Max( Max( A1, A2 ), Max( A3, A4 ) );
1302 const double eps = 0.1; // val is in degrees
1304 return val < eps ? 0. : val;
1307 double Warping::ComputeA( const gp_XYZ& thePnt1,
1308 const gp_XYZ& thePnt2,
1309 const gp_XYZ& thePnt3,
1310 const gp_XYZ& theG ) const
1312 double aLen1 = gp_Pnt( thePnt1 ).Distance( gp_Pnt( thePnt2 ) );
1313 double aLen2 = gp_Pnt( thePnt2 ).Distance( gp_Pnt( thePnt3 ) );
1314 double L = Min( aLen1, aLen2 ) * 0.5;
1318 gp_XYZ GI = ( thePnt2 + thePnt1 ) / 2. - theG;
1319 gp_XYZ GJ = ( thePnt3 + thePnt2 ) / 2. - theG;
1320 gp_XYZ N = GI.Crossed( GJ );
1322 if ( N.Modulus() < gp::Resolution() )
1327 double H = ( thePnt2 - theG ).Dot( N );
1328 return asin( fabs( H / L ) ) * 180. / M_PI;
1331 double Warping::GetBadRate( double Value, int /*nbNodes*/ ) const
1333 // the warp is in the range [0.0,PI/2]
1334 // 0.0 = good (no warp)
1335 // PI/2 = bad (face pliee)
1339 SMDSAbs_ElementType Warping::GetType() const
1341 return SMDSAbs_Face;
1345 //================================================================================
1348 Description : Functor for calculating taper
1350 //================================================================================
1352 double Taper::GetValue( const TSequenceOfXYZ& P )
1354 if ( P.size() != 4 )
1358 double J1 = getArea( P( 4 ), P( 1 ), P( 2 ) ) / 2.;
1359 double J2 = getArea( P( 3 ), P( 1 ), P( 2 ) ) / 2.;
1360 double J3 = getArea( P( 2 ), P( 3 ), P( 4 ) ) / 2.;
1361 double J4 = getArea( P( 3 ), P( 4 ), P( 1 ) ) / 2.;
1363 double JA = 0.25 * ( J1 + J2 + J3 + J4 );
1367 double T1 = fabs( ( J1 - JA ) / JA );
1368 double T2 = fabs( ( J2 - JA ) / JA );
1369 double T3 = fabs( ( J3 - JA ) / JA );
1370 double T4 = fabs( ( J4 - JA ) / JA );
1372 double val = Max( Max( T1, T2 ), Max( T3, T4 ) );
1374 const double eps = 0.01;
1376 return val < eps ? 0. : val;
1379 double Taper::GetBadRate( double Value, int /*nbNodes*/ ) const
1381 // the taper is in the range [0.0,1.0]
1382 // 0.0 = good (no taper)
1383 // 1.0 = bad (les cotes opposes sont allignes)
1387 SMDSAbs_ElementType Taper::GetType() const
1389 return SMDSAbs_Face;
1392 //================================================================================
1395 Description : Functor for calculating skew in degrees
1397 //================================================================================
1399 static inline double skewAngle( const gp_XYZ& p1, const gp_XYZ& p2, const gp_XYZ& p3 )
1401 gp_XYZ p12 = ( p2 + p1 ) / 2.;
1402 gp_XYZ p23 = ( p3 + p2 ) / 2.;
1403 gp_XYZ p31 = ( p3 + p1 ) / 2.;
1405 gp_Vec v1( p31 - p2 ), v2( p12 - p23 );
1407 return v1.Magnitude() < gp::Resolution() || v2.Magnitude() < gp::Resolution() ? 0. : v1.Angle( v2 );
1410 double Skew::GetValue( const TSequenceOfXYZ& P )
1412 if ( P.size() != 3 && P.size() != 4 )
1416 const double PI2 = M_PI / 2.;
1417 if ( P.size() == 3 )
1419 double A0 = fabs( PI2 - skewAngle( P( 3 ), P( 1 ), P( 2 ) ) );
1420 double A1 = fabs( PI2 - skewAngle( P( 1 ), P( 2 ), P( 3 ) ) );
1421 double A2 = fabs( PI2 - skewAngle( P( 2 ), P( 3 ), P( 1 ) ) );
1423 return Max( A0, Max( A1, A2 ) ) * 180. / M_PI;
1427 gp_XYZ p12 = ( P( 1 ) + P( 2 ) ) / 2.;
1428 gp_XYZ p23 = ( P( 2 ) + P( 3 ) ) / 2.;
1429 gp_XYZ p34 = ( P( 3 ) + P( 4 ) ) / 2.;
1430 gp_XYZ p41 = ( P( 4 ) + P( 1 ) ) / 2.;
1432 gp_Vec v1( p34 - p12 ), v2( p23 - p41 );
1433 double A = v1.Magnitude() <= gp::Resolution() || v2.Magnitude() <= gp::Resolution()
1434 ? 0. : fabs( PI2 - v1.Angle( v2 ) );
1436 double val = A * 180. / M_PI;
1438 const double eps = 0.1; // val is in degrees
1440 return val < eps ? 0. : val;
1444 double Skew::GetBadRate( double Value, int /*nbNodes*/ ) const
1446 // the skew is in the range [0.0,PI/2].
1452 SMDSAbs_ElementType Skew::GetType() const
1454 return SMDSAbs_Face;
1458 //================================================================================
1461 Description : Functor for calculating area
1463 //================================================================================
1465 double Area::GetValue( const TSequenceOfXYZ& P )
1468 if ( P.size() > 2 ) {
1469 gp_Vec aVec1( P(2) - P(1) );
1470 gp_Vec aVec2( P(3) - P(1) );
1471 gp_Vec SumVec = aVec1 ^ aVec2;
1472 for (int i=4; i<=P.size(); i++) {
1473 gp_Vec aVec1( P(i-1) - P(1) );
1474 gp_Vec aVec2( P(i) - P(1) );
1475 gp_Vec tmp = aVec1 ^ aVec2;
1478 val = SumVec.Magnitude() * 0.5;
1483 double Area::GetBadRate( double Value, int /*nbNodes*/ ) const
1485 // meaningless as it is not a quality control functor
1489 SMDSAbs_ElementType Area::GetType() const
1491 return SMDSAbs_Face;
1494 //================================================================================
1497 Description : Functor for calculating length of edge
1499 //================================================================================
1501 double Length::GetValue( const TSequenceOfXYZ& P )
1503 switch ( P.size() ) {
1504 case 2: return getDistance( P( 1 ), P( 2 ) );
1505 case 3: return getDistance( P( 1 ), P( 2 ) ) + getDistance( P( 2 ), P( 3 ) );
1510 double Length::GetBadRate( double Value, int /*nbNodes*/ ) const
1512 // meaningless as it is not quality control functor
1516 SMDSAbs_ElementType Length::GetType() const
1518 return SMDSAbs_Edge;
1521 //================================================================================
1524 Description : Functor for calculating length of edge
1526 //================================================================================
1528 double Length2D::GetValue( long theElementId)
1532 //cout<<"Length2D::GetValue"<<endl;
1533 if (GetPoints(theElementId,P)){
1534 //for(int jj=1; jj<=P.size(); jj++)
1535 // cout<<"jj="<<jj<<" P("<<P(jj).X()<<","<<P(jj).Y()<<","<<P(jj).Z()<<")"<<endl;
1537 double aVal;// = GetValue( P );
1538 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
1539 SMDSAbs_ElementType aType = aElem->GetType();
1548 aVal = getDistance( P( 1 ), P( 2 ) );
1551 else if (len == 3){ // quadratic edge
1552 aVal = getDistance(P( 1 ),P( 3 )) + getDistance(P( 3 ),P( 2 ));
1556 if (len == 3){ // triangles
1557 double L1 = getDistance(P( 1 ),P( 2 ));
1558 double L2 = getDistance(P( 2 ),P( 3 ));
1559 double L3 = getDistance(P( 3 ),P( 1 ));
1560 aVal = Max(L1,Max(L2,L3));
1563 else if (len == 4){ // quadrangles
1564 double L1 = getDistance(P( 1 ),P( 2 ));
1565 double L2 = getDistance(P( 2 ),P( 3 ));
1566 double L3 = getDistance(P( 3 ),P( 4 ));
1567 double L4 = getDistance(P( 4 ),P( 1 ));
1568 aVal = Max(Max(L1,L2),Max(L3,L4));
1571 if (len == 6){ // quadratic triangles
1572 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1573 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1574 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
1575 aVal = Max(L1,Max(L2,L3));
1576 //cout<<"L1="<<L1<<" L2="<<L2<<"L3="<<L3<<" aVal="<<aVal<<endl;
1579 else if (len == 8){ // quadratic quadrangles
1580 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1581 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1582 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
1583 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
1584 aVal = Max(Max(L1,L2),Max(L3,L4));
1587 case SMDSAbs_Volume:
1588 if (len == 4){ // tetraidrs
1589 double L1 = getDistance(P( 1 ),P( 2 ));
1590 double L2 = getDistance(P( 2 ),P( 3 ));
1591 double L3 = getDistance(P( 3 ),P( 1 ));
1592 double L4 = getDistance(P( 1 ),P( 4 ));
1593 double L5 = getDistance(P( 2 ),P( 4 ));
1594 double L6 = getDistance(P( 3 ),P( 4 ));
1595 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1598 else if (len == 5){ // piramids
1599 double L1 = getDistance(P( 1 ),P( 2 ));
1600 double L2 = getDistance(P( 2 ),P( 3 ));
1601 double L3 = getDistance(P( 3 ),P( 4 ));
1602 double L4 = getDistance(P( 4 ),P( 1 ));
1603 double L5 = getDistance(P( 1 ),P( 5 ));
1604 double L6 = getDistance(P( 2 ),P( 5 ));
1605 double L7 = getDistance(P( 3 ),P( 5 ));
1606 double L8 = getDistance(P( 4 ),P( 5 ));
1608 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1609 aVal = Max(aVal,Max(L7,L8));
1612 else if (len == 6){ // pentaidres
1613 double L1 = getDistance(P( 1 ),P( 2 ));
1614 double L2 = getDistance(P( 2 ),P( 3 ));
1615 double L3 = getDistance(P( 3 ),P( 1 ));
1616 double L4 = getDistance(P( 4 ),P( 5 ));
1617 double L5 = getDistance(P( 5 ),P( 6 ));
1618 double L6 = getDistance(P( 6 ),P( 4 ));
1619 double L7 = getDistance(P( 1 ),P( 4 ));
1620 double L8 = getDistance(P( 2 ),P( 5 ));
1621 double L9 = getDistance(P( 3 ),P( 6 ));
1623 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1624 aVal = Max(aVal,Max(Max(L7,L8),L9));
1627 else if (len == 8){ // hexaider
1628 double L1 = getDistance(P( 1 ),P( 2 ));
1629 double L2 = getDistance(P( 2 ),P( 3 ));
1630 double L3 = getDistance(P( 3 ),P( 4 ));
1631 double L4 = getDistance(P( 4 ),P( 1 ));
1632 double L5 = getDistance(P( 5 ),P( 6 ));
1633 double L6 = getDistance(P( 6 ),P( 7 ));
1634 double L7 = getDistance(P( 7 ),P( 8 ));
1635 double L8 = getDistance(P( 8 ),P( 5 ));
1636 double L9 = getDistance(P( 1 ),P( 5 ));
1637 double L10= getDistance(P( 2 ),P( 6 ));
1638 double L11= getDistance(P( 3 ),P( 7 ));
1639 double L12= getDistance(P( 4 ),P( 8 ));
1641 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1642 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1643 aVal = Max(aVal,Max(L11,L12));
1648 if (len == 10){ // quadratic tetraidrs
1649 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
1650 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
1651 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
1652 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1653 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
1654 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
1655 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1658 else if (len == 13){ // quadratic piramids
1659 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
1660 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
1661 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1662 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1663 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1664 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
1665 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
1666 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
1667 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1668 aVal = Max(aVal,Max(L7,L8));
1671 else if (len == 15){ // quadratic pentaidres
1672 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
1673 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
1674 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1675 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1676 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
1677 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
1678 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
1679 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
1680 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
1681 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1682 aVal = Max(aVal,Max(Max(L7,L8),L9));
1685 else if (len == 20){ // quadratic hexaider
1686 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
1687 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
1688 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
1689 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
1690 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
1691 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
1692 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
1693 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
1694 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
1695 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
1696 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
1697 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
1698 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1699 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1700 aVal = Max(aVal,Max(L11,L12));
1712 if ( myPrecision >= 0 )
1714 double prec = pow( 10., (double)( myPrecision ) );
1715 aVal = floor( aVal * prec + 0.5 ) / prec;
1724 double Length2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1726 // meaningless as it is not quality control functor
1730 SMDSAbs_ElementType Length2D::GetType() const
1732 return SMDSAbs_Face;
1735 Length2D::Value::Value(double theLength,long thePntId1, long thePntId2):
1738 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1739 if(thePntId1 > thePntId2){
1740 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1744 bool Length2D::Value::operator<(const Length2D::Value& x) const{
1745 if(myPntId[0] < x.myPntId[0]) return true;
1746 if(myPntId[0] == x.myPntId[0])
1747 if(myPntId[1] < x.myPntId[1]) return true;
1751 void Length2D::GetValues(TValues& theValues){
1753 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1754 for(; anIter->more(); ){
1755 const SMDS_MeshFace* anElem = anIter->next();
1757 if(anElem->IsQuadratic()) {
1758 const SMDS_VtkFace* F =
1759 dynamic_cast<const SMDS_VtkFace*>(anElem);
1760 // use special nodes iterator
1761 SMDS_ElemIteratorPtr anIter = F->interlacedNodesElemIterator();
1766 const SMDS_MeshElement* aNode;
1768 aNode = anIter->next();
1769 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1770 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1771 aNodeId[0] = aNodeId[1] = aNode->GetID();
1774 for(; anIter->more(); ){
1775 const SMDS_MeshNode* N1 = static_cast<const SMDS_MeshNode*> (anIter->next());
1776 P[2] = gp_Pnt(N1->X(),N1->Y(),N1->Z());
1777 aNodeId[2] = N1->GetID();
1778 aLength = P[1].Distance(P[2]);
1779 if(!anIter->more()) break;
1780 const SMDS_MeshNode* N2 = static_cast<const SMDS_MeshNode*> (anIter->next());
1781 P[3] = gp_Pnt(N2->X(),N2->Y(),N2->Z());
1782 aNodeId[3] = N2->GetID();
1783 aLength += P[2].Distance(P[3]);
1784 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1785 Value aValue2(aLength,aNodeId[2],aNodeId[3]);
1787 aNodeId[1] = aNodeId[3];
1788 theValues.insert(aValue1);
1789 theValues.insert(aValue2);
1791 aLength += P[2].Distance(P[0]);
1792 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1793 Value aValue2(aLength,aNodeId[2],aNodeId[0]);
1794 theValues.insert(aValue1);
1795 theValues.insert(aValue2);
1798 SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
1803 const SMDS_MeshElement* aNode;
1804 if(aNodesIter->more()){
1805 aNode = aNodesIter->next();
1806 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1807 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1808 aNodeId[0] = aNodeId[1] = aNode->GetID();
1811 for(; aNodesIter->more(); ){
1812 aNode = aNodesIter->next();
1813 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1814 long anId = aNode->GetID();
1816 P[2] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1818 aLength = P[1].Distance(P[2]);
1820 Value aValue(aLength,aNodeId[1],anId);
1823 theValues.insert(aValue);
1826 aLength = P[0].Distance(P[1]);
1828 Value aValue(aLength,aNodeId[0],aNodeId[1]);
1829 theValues.insert(aValue);
1834 //================================================================================
1836 Class : MultiConnection
1837 Description : Functor for calculating number of faces conneted to the edge
1839 //================================================================================
1841 double MultiConnection::GetValue( const TSequenceOfXYZ& P )
1845 double MultiConnection::GetValue( long theId )
1847 return getNbMultiConnection( myMesh, theId );
1850 double MultiConnection::GetBadRate( double Value, int /*nbNodes*/ ) const
1852 // meaningless as it is not quality control functor
1856 SMDSAbs_ElementType MultiConnection::GetType() const
1858 return SMDSAbs_Edge;
1861 //================================================================================
1863 Class : MultiConnection2D
1864 Description : Functor for calculating number of faces conneted to the edge
1866 //================================================================================
1868 double MultiConnection2D::GetValue( const TSequenceOfXYZ& P )
1873 double MultiConnection2D::GetValue( long theElementId )
1877 const SMDS_MeshElement* aFaceElem = myMesh->FindElement(theElementId);
1878 SMDSAbs_ElementType aType = aFaceElem->GetType();
1883 int i = 0, len = aFaceElem->NbNodes();
1884 SMDS_ElemIteratorPtr anIter = aFaceElem->nodesIterator();
1887 const SMDS_MeshNode *aNode, *aNode0;
1888 TColStd_MapOfInteger aMap, aMapPrev;
1890 for (i = 0; i <= len; i++) {
1895 if (anIter->more()) {
1896 aNode = (SMDS_MeshNode*)anIter->next();
1904 if (i == 0) aNode0 = aNode;
1906 SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
1907 while (anElemIter->more()) {
1908 const SMDS_MeshElement* anElem = anElemIter->next();
1909 if (anElem != 0 && anElem->GetType() == SMDSAbs_Face) {
1910 int anId = anElem->GetID();
1913 if (aMapPrev.Contains(anId)) {
1918 aResult = Max(aResult, aNb);
1929 double MultiConnection2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1931 // meaningless as it is not quality control functor
1935 SMDSAbs_ElementType MultiConnection2D::GetType() const
1937 return SMDSAbs_Face;
1940 MultiConnection2D::Value::Value(long thePntId1, long thePntId2)
1942 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1943 if(thePntId1 > thePntId2){
1944 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1948 bool MultiConnection2D::Value::operator<(const MultiConnection2D::Value& x) const{
1949 if(myPntId[0] < x.myPntId[0]) return true;
1950 if(myPntId[0] == x.myPntId[0])
1951 if(myPntId[1] < x.myPntId[1]) return true;
1955 void MultiConnection2D::GetValues(MValues& theValues){
1956 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1957 for(; anIter->more(); ){
1958 const SMDS_MeshFace* anElem = anIter->next();
1959 SMDS_ElemIteratorPtr aNodesIter;
1960 if ( anElem->IsQuadratic() )
1961 aNodesIter = dynamic_cast<const SMDS_VtkFace*>
1962 (anElem)->interlacedNodesElemIterator();
1964 aNodesIter = anElem->nodesIterator();
1967 //int aNbConnects=0;
1968 const SMDS_MeshNode* aNode0;
1969 const SMDS_MeshNode* aNode1;
1970 const SMDS_MeshNode* aNode2;
1971 if(aNodesIter->more()){
1972 aNode0 = (SMDS_MeshNode*) aNodesIter->next();
1974 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode1;
1975 aNodeId[0] = aNodeId[1] = aNodes->GetID();
1977 for(; aNodesIter->more(); ) {
1978 aNode2 = (SMDS_MeshNode*) aNodesIter->next();
1979 long anId = aNode2->GetID();
1982 Value aValue(aNodeId[1],aNodeId[2]);
1983 MValues::iterator aItr = theValues.find(aValue);
1984 if (aItr != theValues.end()){
1989 theValues[aValue] = 1;
1992 //cout << "NodeIds: "<<aNodeId[1]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1993 aNodeId[1] = aNodeId[2];
1996 Value aValue(aNodeId[0],aNodeId[2]);
1997 MValues::iterator aItr = theValues.find(aValue);
1998 if (aItr != theValues.end()) {
2003 theValues[aValue] = 1;
2006 //cout << "NodeIds: "<<aNodeId[0]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
2011 //================================================================================
2013 Class : BallDiameter
2014 Description : Functor returning diameter of a ball element
2016 //================================================================================
2018 double BallDiameter::GetValue( long theId )
2020 double diameter = 0;
2022 if ( const SMDS_BallElement* ball =
2023 dynamic_cast<const SMDS_BallElement*>( myMesh->FindElement( theId )))
2025 diameter = ball->GetDiameter();
2030 double BallDiameter::GetBadRate( double Value, int /*nbNodes*/ ) const
2032 // meaningless as it is not a quality control functor
2036 SMDSAbs_ElementType BallDiameter::GetType() const
2038 return SMDSAbs_Ball;
2046 //================================================================================
2048 Class : BadOrientedVolume
2049 Description : Predicate bad oriented volumes
2051 //================================================================================
2053 BadOrientedVolume::BadOrientedVolume()
2058 void BadOrientedVolume::SetMesh( const SMDS_Mesh* theMesh )
2063 bool BadOrientedVolume::IsSatisfy( long theId )
2068 SMDS_VolumeTool vTool( myMesh->FindElement( theId ));
2069 return !vTool.IsForward();
2072 SMDSAbs_ElementType BadOrientedVolume::GetType() const
2074 return SMDSAbs_Volume;
2078 Class : BareBorderVolume
2081 bool BareBorderVolume::IsSatisfy(long theElementId )
2083 SMDS_VolumeTool myTool;
2084 if ( myTool.Set( myMesh->FindElement(theElementId)))
2086 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2087 if ( myTool.IsFreeFace( iF ))
2089 const SMDS_MeshNode** n = myTool.GetFaceNodes(iF);
2090 vector< const SMDS_MeshNode*> nodes( n, n+myTool.NbFaceNodes(iF));
2091 if ( !myMesh->FindElement( nodes, SMDSAbs_Face, /*Nomedium=*/false))
2098 //================================================================================
2100 Class : BareBorderFace
2102 //================================================================================
2104 bool BareBorderFace::IsSatisfy(long theElementId )
2107 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2109 if ( face->GetType() == SMDSAbs_Face )
2111 int nbN = face->NbCornerNodes();
2112 for ( int i = 0; i < nbN && !ok; ++i )
2114 // check if a link is shared by another face
2115 const SMDS_MeshNode* n1 = face->GetNode( i );
2116 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2117 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2118 bool isShared = false;
2119 while ( !isShared && fIt->more() )
2121 const SMDS_MeshElement* f = fIt->next();
2122 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2126 const int iQuad = face->IsQuadratic();
2127 myLinkNodes.resize( 2 + iQuad);
2128 myLinkNodes[0] = n1;
2129 myLinkNodes[1] = n2;
2131 myLinkNodes[2] = face->GetNode( i+nbN );
2132 ok = !myMesh->FindElement( myLinkNodes, SMDSAbs_Edge, /*noMedium=*/false);
2140 //================================================================================
2142 Class : OverConstrainedVolume
2144 //================================================================================
2146 bool OverConstrainedVolume::IsSatisfy(long theElementId )
2148 // An element is over-constrained if it has N-1 free borders where
2149 // N is the number of edges/faces for a 2D/3D element.
2150 SMDS_VolumeTool myTool;
2151 if ( myTool.Set( myMesh->FindElement(theElementId)))
2153 int nbSharedFaces = 0;
2154 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2155 if ( !myTool.IsFreeFace( iF ) && ++nbSharedFaces > 1 )
2157 return ( nbSharedFaces == 1 );
2162 //================================================================================
2164 Class : OverConstrainedFace
2166 //================================================================================
2168 bool OverConstrainedFace::IsSatisfy(long theElementId )
2170 // An element is over-constrained if it has N-1 free borders where
2171 // N is the number of edges/faces for a 2D/3D element.
2172 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2173 if ( face->GetType() == SMDSAbs_Face )
2175 int nbSharedBorders = 0;
2176 int nbN = face->NbCornerNodes();
2177 for ( int i = 0; i < nbN; ++i )
2179 // check if a link is shared by another face
2180 const SMDS_MeshNode* n1 = face->GetNode( i );
2181 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2182 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2183 bool isShared = false;
2184 while ( !isShared && fIt->more() )
2186 const SMDS_MeshElement* f = fIt->next();
2187 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2189 if ( isShared && ++nbSharedBorders > 1 )
2192 return ( nbSharedBorders == 1 );
2197 //================================================================================
2199 Class : CoincidentNodes
2200 Description : Predicate of Coincident nodes
2202 //================================================================================
2204 CoincidentNodes::CoincidentNodes()
2209 bool CoincidentNodes::IsSatisfy( long theElementId )
2211 return myCoincidentIDs.Contains( theElementId );
2214 SMDSAbs_ElementType CoincidentNodes::GetType() const
2216 return SMDSAbs_Node;
2219 void CoincidentNodes::SetMesh( const SMDS_Mesh* theMesh )
2221 myMeshModifTracer.SetMesh( theMesh );
2222 if ( myMeshModifTracer.IsMeshModified() )
2224 TIDSortedNodeSet nodesToCheck;
2225 SMDS_NodeIteratorPtr nIt = theMesh->nodesIterator(/*idInceasingOrder=*/true);
2226 while ( nIt->more() )
2227 nodesToCheck.insert( nodesToCheck.end(), nIt->next() );
2229 list< list< const SMDS_MeshNode*> > nodeGroups;
2230 SMESH_OctreeNode::FindCoincidentNodes ( nodesToCheck, &nodeGroups, myToler );
2232 myCoincidentIDs.Clear();
2233 list< list< const SMDS_MeshNode*> >::iterator groupIt = nodeGroups.begin();
2234 for ( ; groupIt != nodeGroups.end(); ++groupIt )
2236 list< const SMDS_MeshNode*>& coincNodes = *groupIt;
2237 list< const SMDS_MeshNode*>::iterator n = coincNodes.begin();
2238 for ( ; n != coincNodes.end(); ++n )
2239 myCoincidentIDs.Add( (*n)->GetID() );
2244 //================================================================================
2246 Class : CoincidentElements
2247 Description : Predicate of Coincident Elements
2248 Note : This class is suitable only for visualization of Coincident Elements
2250 //================================================================================
2252 CoincidentElements::CoincidentElements()
2257 void CoincidentElements::SetMesh( const SMDS_Mesh* theMesh )
2262 bool CoincidentElements::IsSatisfy( long theElementId )
2264 if ( !myMesh ) return false;
2266 if ( const SMDS_MeshElement* e = myMesh->FindElement( theElementId ))
2268 if ( e->GetType() != GetType() ) return false;
2269 set< const SMDS_MeshNode* > elemNodes( e->begin_nodes(), e->end_nodes() );
2270 const int nbNodes = e->NbNodes();
2271 SMDS_ElemIteratorPtr invIt = (*elemNodes.begin())->GetInverseElementIterator( GetType() );
2272 while ( invIt->more() )
2274 const SMDS_MeshElement* e2 = invIt->next();
2275 if ( e2 == e || e2->NbNodes() != nbNodes ) continue;
2277 bool sameNodes = true;
2278 for ( size_t i = 0; i < elemNodes.size() && sameNodes; ++i )
2279 sameNodes = ( elemNodes.count( e2->GetNode( i )));
2287 SMDSAbs_ElementType CoincidentElements1D::GetType() const
2289 return SMDSAbs_Edge;
2291 SMDSAbs_ElementType CoincidentElements2D::GetType() const
2293 return SMDSAbs_Face;
2295 SMDSAbs_ElementType CoincidentElements3D::GetType() const
2297 return SMDSAbs_Volume;
2301 //================================================================================
2304 Description : Predicate for free borders
2306 //================================================================================
2308 FreeBorders::FreeBorders()
2313 void FreeBorders::SetMesh( const SMDS_Mesh* theMesh )
2318 bool FreeBorders::IsSatisfy( long theId )
2320 return getNbMultiConnection( myMesh, theId ) == 1;
2323 SMDSAbs_ElementType FreeBorders::GetType() const
2325 return SMDSAbs_Edge;
2329 //================================================================================
2332 Description : Predicate for free Edges
2334 //================================================================================
2336 FreeEdges::FreeEdges()
2341 void FreeEdges::SetMesh( const SMDS_Mesh* theMesh )
2346 bool FreeEdges::IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId )
2348 TColStd_MapOfInteger aMap;
2349 for ( int i = 0; i < 2; i++ )
2351 SMDS_ElemIteratorPtr anElemIter = theNodes[ i ]->GetInverseElementIterator(SMDSAbs_Face);
2352 while( anElemIter->more() )
2354 if ( const SMDS_MeshElement* anElem = anElemIter->next())
2356 const int anId = anElem->GetID();
2357 if ( anId != theFaceId && !aMap.Add( anId ))
2365 bool FreeEdges::IsSatisfy( long theId )
2370 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2371 if ( aFace == 0 || aFace->GetType() != SMDSAbs_Face || aFace->NbNodes() < 3 )
2374 SMDS_ElemIteratorPtr anIter;
2375 if ( aFace->IsQuadratic() ) {
2376 anIter = dynamic_cast<const SMDS_VtkFace*>
2377 (aFace)->interlacedNodesElemIterator();
2380 anIter = aFace->nodesIterator();
2385 int i = 0, nbNodes = aFace->NbNodes();
2386 std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
2387 while( anIter->more() )
2389 const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
2392 aNodes[ i++ ] = aNode;
2394 aNodes[ nbNodes ] = aNodes[ 0 ];
2396 for ( i = 0; i < nbNodes; i++ )
2397 if ( IsFreeEdge( &aNodes[ i ], theId ) )
2403 SMDSAbs_ElementType FreeEdges::GetType() const
2405 return SMDSAbs_Face;
2408 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2411 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
2412 if(thePntId1 > thePntId2){
2413 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
2417 bool FreeEdges::Border::operator<(const FreeEdges::Border& x) const{
2418 if(myPntId[0] < x.myPntId[0]) return true;
2419 if(myPntId[0] == x.myPntId[0])
2420 if(myPntId[1] < x.myPntId[1]) return true;
2424 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2425 FreeEdges::TBorders& theRegistry,
2426 FreeEdges::TBorders& theContainer)
2428 if(theRegistry.find(theBorder) == theRegistry.end()){
2429 theRegistry.insert(theBorder);
2430 theContainer.insert(theBorder);
2432 theContainer.erase(theBorder);
2436 void FreeEdges::GetBoreders(TBorders& theBorders)
2439 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2440 for(; anIter->more(); ){
2441 const SMDS_MeshFace* anElem = anIter->next();
2442 long anElemId = anElem->GetID();
2443 SMDS_ElemIteratorPtr aNodesIter;
2444 if ( anElem->IsQuadratic() )
2445 aNodesIter = static_cast<const SMDS_VtkFace*>(anElem)->
2446 interlacedNodesElemIterator();
2448 aNodesIter = anElem->nodesIterator();
2450 const SMDS_MeshElement* aNode;
2451 if(aNodesIter->more()){
2452 aNode = aNodesIter->next();
2453 aNodeId[0] = aNodeId[1] = aNode->GetID();
2455 for(; aNodesIter->more(); ){
2456 aNode = aNodesIter->next();
2457 long anId = aNode->GetID();
2458 Border aBorder(anElemId,aNodeId[1],anId);
2460 UpdateBorders(aBorder,aRegistry,theBorders);
2462 Border aBorder(anElemId,aNodeId[0],aNodeId[1]);
2463 UpdateBorders(aBorder,aRegistry,theBorders);
2467 //================================================================================
2470 Description : Predicate for free nodes
2472 //================================================================================
2474 FreeNodes::FreeNodes()
2479 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2484 bool FreeNodes::IsSatisfy( long theNodeId )
2486 const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2490 return (aNode->NbInverseElements() < 1);
2493 SMDSAbs_ElementType FreeNodes::GetType() const
2495 return SMDSAbs_Node;
2499 //================================================================================
2502 Description : Predicate for free faces
2504 //================================================================================
2506 FreeFaces::FreeFaces()
2511 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2516 bool FreeFaces::IsSatisfy( long theId )
2518 if (!myMesh) return false;
2519 // check that faces nodes refers to less than two common volumes
2520 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2521 if ( !aFace || aFace->GetType() != SMDSAbs_Face )
2524 int nbNode = aFace->NbNodes();
2526 // collect volumes check that number of volumss with count equal nbNode not less than 2
2527 typedef map< SMDS_MeshElement*, int > TMapOfVolume; // map of volume counters
2528 typedef map< SMDS_MeshElement*, int >::iterator TItrMapOfVolume; // iterator
2529 TMapOfVolume mapOfVol;
2531 SMDS_ElemIteratorPtr nodeItr = aFace->nodesIterator();
2532 while ( nodeItr->more() ) {
2533 const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>(nodeItr->next());
2534 if ( !aNode ) continue;
2535 SMDS_ElemIteratorPtr volItr = aNode->GetInverseElementIterator(SMDSAbs_Volume);
2536 while ( volItr->more() ) {
2537 SMDS_MeshElement* aVol = (SMDS_MeshElement*)volItr->next();
2538 TItrMapOfVolume itr = mapOfVol.insert(make_pair(aVol, 0)).first;
2543 TItrMapOfVolume volItr = mapOfVol.begin();
2544 TItrMapOfVolume volEnd = mapOfVol.end();
2545 for ( ; volItr != volEnd; ++volItr )
2546 if ( (*volItr).second >= nbNode )
2548 // face is not free if number of volumes constructed on thier nodes more than one
2552 SMDSAbs_ElementType FreeFaces::GetType() const
2554 return SMDSAbs_Face;
2557 //================================================================================
2559 Class : LinearOrQuadratic
2560 Description : Predicate to verify whether a mesh element is linear
2562 //================================================================================
2564 LinearOrQuadratic::LinearOrQuadratic()
2569 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2574 bool LinearOrQuadratic::IsSatisfy( long theId )
2576 if (!myMesh) return false;
2577 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2578 if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2580 return (!anElem->IsQuadratic());
2583 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2588 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2593 //================================================================================
2596 Description : Functor for check color of group to whic mesh element belongs to
2598 //================================================================================
2600 GroupColor::GroupColor()
2604 bool GroupColor::IsSatisfy( long theId )
2606 return (myIDs.find( theId ) != myIDs.end());
2609 void GroupColor::SetType( SMDSAbs_ElementType theType )
2614 SMDSAbs_ElementType GroupColor::GetType() const
2619 static bool isEqual( const Quantity_Color& theColor1,
2620 const Quantity_Color& theColor2 )
2622 // tolerance to compare colors
2623 const double tol = 5*1e-3;
2624 return ( fabs( theColor1.Red() - theColor2.Red() ) < tol &&
2625 fabs( theColor1.Green() - theColor2.Green() ) < tol &&
2626 fabs( theColor1.Blue() - theColor2.Blue() ) < tol );
2630 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2634 const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2638 int nbGrp = aMesh->GetNbGroups();
2642 // iterates on groups and find necessary elements ids
2643 const std::set<SMESHDS_GroupBase*>& aGroups = aMesh->GetGroups();
2644 set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
2645 for (; GrIt != aGroups.end(); GrIt++) {
2646 SMESHDS_GroupBase* aGrp = (*GrIt);
2649 // check type and color of group
2650 if ( !isEqual( myColor, aGrp->GetColor() ) )
2652 if ( myType != SMDSAbs_All && myType != (SMDSAbs_ElementType)aGrp->GetType() )
2655 SMDSAbs_ElementType aGrpElType = (SMDSAbs_ElementType)aGrp->GetType();
2656 if ( myType == aGrpElType || (myType == SMDSAbs_All && aGrpElType != SMDSAbs_Node) ) {
2657 // add elements IDS into control
2658 int aSize = aGrp->Extent();
2659 for (int i = 0; i < aSize; i++)
2660 myIDs.insert( aGrp->GetID(i+1) );
2665 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
2667 Kernel_Utils::Localizer loc;
2668 TCollection_AsciiString aStr = theStr;
2669 aStr.RemoveAll( ' ' );
2670 aStr.RemoveAll( '\t' );
2671 for ( int aPos = aStr.Search( ";;" ); aPos != -1; aPos = aStr.Search( ";;" ) )
2672 aStr.Remove( aPos, 2 );
2673 Standard_Real clr[3];
2674 clr[0] = clr[1] = clr[2] = 0.;
2675 for ( int i = 0; i < 3; i++ ) {
2676 TCollection_AsciiString tmpStr = aStr.Token( ";", i+1 );
2677 if ( !tmpStr.IsEmpty() && tmpStr.IsRealValue() )
2678 clr[i] = tmpStr.RealValue();
2680 myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
2683 //=======================================================================
2684 // name : GetRangeStr
2685 // Purpose : Get range as a string.
2686 // Example: "1,2,3,50-60,63,67,70-"
2687 //=======================================================================
2689 void GroupColor::GetColorStr( TCollection_AsciiString& theResStr ) const
2692 theResStr += TCollection_AsciiString( myColor.Red() );
2693 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2694 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2697 //================================================================================
2699 Class : ElemGeomType
2700 Description : Predicate to check element geometry type
2702 //================================================================================
2704 ElemGeomType::ElemGeomType()
2707 myType = SMDSAbs_All;
2708 myGeomType = SMDSGeom_TRIANGLE;
2711 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2716 bool ElemGeomType::IsSatisfy( long theId )
2718 if (!myMesh) return false;
2719 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2722 const SMDSAbs_ElementType anElemType = anElem->GetType();
2723 if ( myType != SMDSAbs_All && anElemType != myType )
2725 bool isOk = ( anElem->GetGeomType() == myGeomType );
2729 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2734 SMDSAbs_ElementType ElemGeomType::GetType() const
2739 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2741 myGeomType = theType;
2744 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2749 //================================================================================
2751 Class : ElemEntityType
2752 Description : Predicate to check element entity type
2754 //================================================================================
2756 ElemEntityType::ElemEntityType():
2758 myType( SMDSAbs_All ),
2759 myEntityType( SMDSEntity_0D )
2763 void ElemEntityType::SetMesh( const SMDS_Mesh* theMesh )
2768 bool ElemEntityType::IsSatisfy( long theId )
2770 if ( !myMesh ) return false;
2771 if ( myType == SMDSAbs_Node )
2772 return myMesh->FindNode( theId );
2773 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2775 myEntityType == anElem->GetEntityType() );
2778 void ElemEntityType::SetType( SMDSAbs_ElementType theType )
2783 SMDSAbs_ElementType ElemEntityType::GetType() const
2788 void ElemEntityType::SetElemEntityType( SMDSAbs_EntityType theEntityType )
2790 myEntityType = theEntityType;
2793 SMDSAbs_EntityType ElemEntityType::GetElemEntityType() const
2795 return myEntityType;
2798 //================================================================================
2800 * \brief Class ConnectedElements
2802 //================================================================================
2804 ConnectedElements::ConnectedElements():
2805 myNodeID(0), myType( SMDSAbs_All ), myOkIDsReady( false ) {}
2807 SMDSAbs_ElementType ConnectedElements::GetType() const
2810 int ConnectedElements::GetNode() const
2811 { return myXYZ.empty() ? myNodeID : 0; } // myNodeID can be found by myXYZ
2813 std::vector<double> ConnectedElements::GetPoint() const
2816 void ConnectedElements::clearOkIDs()
2817 { myOkIDsReady = false; myOkIDs.clear(); }
2819 void ConnectedElements::SetType( SMDSAbs_ElementType theType )
2821 if ( myType != theType || myMeshModifTracer.IsMeshModified() )
2826 void ConnectedElements::SetMesh( const SMDS_Mesh* theMesh )
2828 myMeshModifTracer.SetMesh( theMesh );
2829 if ( myMeshModifTracer.IsMeshModified() )
2832 if ( !myXYZ.empty() )
2833 SetPoint( myXYZ[0], myXYZ[1], myXYZ[2] ); // find a node near myXYZ it in a new mesh
2837 void ConnectedElements::SetNode( int nodeID )
2842 bool isSameDomain = false;
2843 if ( myOkIDsReady && myMeshModifTracer.GetMesh() && !myMeshModifTracer.IsMeshModified() )
2844 if ( const SMDS_MeshNode* n = myMeshModifTracer.GetMesh()->FindNode( myNodeID ))
2846 SMDS_ElemIteratorPtr eIt = n->GetInverseElementIterator( myType );
2847 while ( !isSameDomain && eIt->more() )
2848 isSameDomain = IsSatisfy( eIt->next()->GetID() );
2850 if ( !isSameDomain )
2854 void ConnectedElements::SetPoint( double x, double y, double z )
2862 bool isSameDomain = false;
2864 // find myNodeID by myXYZ if possible
2865 if ( myMeshModifTracer.GetMesh() )
2867 auto_ptr<SMESH_ElementSearcher> searcher
2868 ( SMESH_MeshAlgos::GetElementSearcher( (SMDS_Mesh&) *myMeshModifTracer.GetMesh() ));
2870 vector< const SMDS_MeshElement* > foundElems;
2871 searcher->FindElementsByPoint( gp_Pnt(x,y,z), SMDSAbs_All, foundElems );
2873 if ( !foundElems.empty() )
2875 myNodeID = foundElems[0]->GetNode(0)->GetID();
2876 if ( myOkIDsReady && !myMeshModifTracer.IsMeshModified() )
2877 isSameDomain = IsSatisfy( foundElems[0]->GetID() );
2880 if ( !isSameDomain )
2884 bool ConnectedElements::IsSatisfy( long theElementId )
2886 // Here we do NOT check if the mesh has changed, we do it in Set...() only!!!
2888 if ( !myOkIDsReady )
2890 if ( !myMeshModifTracer.GetMesh() )
2892 const SMDS_MeshNode* node0 = myMeshModifTracer.GetMesh()->FindNode( myNodeID );
2896 list< const SMDS_MeshNode* > nodeQueue( 1, node0 );
2897 std::set< int > checkedNodeIDs;
2899 // foreach node in nodeQueue:
2900 // foreach element sharing a node:
2901 // add ID of an element of myType to myOkIDs;
2902 // push all element nodes absent from checkedNodeIDs to nodeQueue;
2903 while ( !nodeQueue.empty() )
2905 const SMDS_MeshNode* node = nodeQueue.front();
2906 nodeQueue.pop_front();
2908 // loop on elements sharing the node
2909 SMDS_ElemIteratorPtr eIt = node->GetInverseElementIterator();
2910 while ( eIt->more() )
2912 // keep elements of myType
2913 const SMDS_MeshElement* element = eIt->next();
2914 if ( element->GetType() == myType )
2915 myOkIDs.insert( myOkIDs.end(), element->GetID() );
2917 // enqueue nodes of the element
2918 SMDS_ElemIteratorPtr nIt = element->nodesIterator();
2919 while ( nIt->more() )
2921 const SMDS_MeshNode* n = static_cast< const SMDS_MeshNode* >( nIt->next() );
2922 if ( checkedNodeIDs.insert( n->GetID() ).second )
2923 nodeQueue.push_back( n );
2927 if ( myType == SMDSAbs_Node )
2928 std::swap( myOkIDs, checkedNodeIDs );
2930 size_t totalNbElems = myMeshModifTracer.GetMesh()->GetMeshInfo().NbElements( myType );
2931 if ( myOkIDs.size() == totalNbElems )
2934 myOkIDsReady = true;
2937 return myOkIDs.empty() ? true : myOkIDs.count( theElementId );
2940 //================================================================================
2942 * \brief Class CoplanarFaces
2944 //================================================================================
2946 CoplanarFaces::CoplanarFaces()
2947 : myFaceID(0), myToler(0)
2950 void CoplanarFaces::SetMesh( const SMDS_Mesh* theMesh )
2952 myMeshModifTracer.SetMesh( theMesh );
2953 if ( myMeshModifTracer.IsMeshModified() )
2955 // Build a set of coplanar face ids
2957 myCoplanarIDs.clear();
2959 if ( !myMeshModifTracer.GetMesh() || !myFaceID || !myToler )
2962 const SMDS_MeshElement* face = myMeshModifTracer.GetMesh()->FindElement( myFaceID );
2963 if ( !face || face->GetType() != SMDSAbs_Face )
2967 gp_Vec myNorm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2971 const double radianTol = myToler * M_PI / 180.;
2972 std::set< SMESH_TLink > checkedLinks;
2974 std::list< pair< const SMDS_MeshElement*, gp_Vec > > faceQueue;
2975 faceQueue.push_back( make_pair( face, myNorm ));
2976 while ( !faceQueue.empty() )
2978 face = faceQueue.front().first;
2979 myNorm = faceQueue.front().second;
2980 faceQueue.pop_front();
2982 for ( int i = 0, nbN = face->NbCornerNodes(); i < nbN; ++i )
2984 const SMDS_MeshNode* n1 = face->GetNode( i );
2985 const SMDS_MeshNode* n2 = face->GetNode(( i+1 )%nbN);
2986 if ( !checkedLinks.insert( SMESH_TLink( n1, n2 )).second )
2988 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator(SMDSAbs_Face);
2989 while ( fIt->more() )
2991 const SMDS_MeshElement* f = fIt->next();
2992 if ( f->GetNodeIndex( n2 ) > -1 )
2994 gp_Vec norm = getNormale( static_cast<const SMDS_MeshFace*>(f), &normOK );
2995 if (!normOK || myNorm.Angle( norm ) <= radianTol)
2997 myCoplanarIDs.insert( f->GetID() );
2998 faceQueue.push_back( make_pair( f, norm ));
3006 bool CoplanarFaces::IsSatisfy( long theElementId )
3008 return myCoplanarIDs.count( theElementId );
3013 *Description : Predicate for Range of Ids.
3014 * Range may be specified with two ways.
3015 * 1. Using AddToRange method
3016 * 2. With SetRangeStr method. Parameter of this method is a string
3017 * like as "1,2,3,50-60,63,67,70-"
3020 //=======================================================================
3021 // name : RangeOfIds
3022 // Purpose : Constructor
3023 //=======================================================================
3024 RangeOfIds::RangeOfIds()
3027 myType = SMDSAbs_All;
3030 //=======================================================================
3032 // Purpose : Set mesh
3033 //=======================================================================
3034 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
3039 //=======================================================================
3040 // name : AddToRange
3041 // Purpose : Add ID to the range
3042 //=======================================================================
3043 bool RangeOfIds::AddToRange( long theEntityId )
3045 myIds.Add( theEntityId );
3049 //=======================================================================
3050 // name : GetRangeStr
3051 // Purpose : Get range as a string.
3052 // Example: "1,2,3,50-60,63,67,70-"
3053 //=======================================================================
3054 void RangeOfIds::GetRangeStr( TCollection_AsciiString& theResStr )
3058 TColStd_SequenceOfInteger anIntSeq;
3059 TColStd_SequenceOfAsciiString aStrSeq;
3061 TColStd_MapIteratorOfMapOfInteger anIter( myIds );
3062 for ( ; anIter.More(); anIter.Next() )
3064 int anId = anIter.Key();
3065 TCollection_AsciiString aStr( anId );
3066 anIntSeq.Append( anId );
3067 aStrSeq.Append( aStr );
3070 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
3072 int aMinId = myMin( i );
3073 int aMaxId = myMax( i );
3075 TCollection_AsciiString aStr;
3076 if ( aMinId != IntegerFirst() )
3081 if ( aMaxId != IntegerLast() )
3084 // find position of the string in result sequence and insert string in it
3085 if ( anIntSeq.Length() == 0 )
3087 anIntSeq.Append( aMinId );
3088 aStrSeq.Append( aStr );
3092 if ( aMinId < anIntSeq.First() )
3094 anIntSeq.Prepend( aMinId );
3095 aStrSeq.Prepend( aStr );
3097 else if ( aMinId > anIntSeq.Last() )
3099 anIntSeq.Append( aMinId );
3100 aStrSeq.Append( aStr );
3103 for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
3104 if ( aMinId < anIntSeq( j ) )
3106 anIntSeq.InsertBefore( j, aMinId );
3107 aStrSeq.InsertBefore( j, aStr );
3113 if ( aStrSeq.Length() == 0 )
3116 theResStr = aStrSeq( 1 );
3117 for ( int j = 2, k = aStrSeq.Length(); j <= k; j++ )
3120 theResStr += aStrSeq( j );
3124 //=======================================================================
3125 // name : SetRangeStr
3126 // Purpose : Define range with string
3127 // Example of entry string: "1,2,3,50-60,63,67,70-"
3128 //=======================================================================
3129 bool RangeOfIds::SetRangeStr( const TCollection_AsciiString& theStr )
3135 TCollection_AsciiString aStr = theStr;
3136 aStr.RemoveAll( ' ' );
3137 aStr.RemoveAll( '\t' );
3139 for ( int aPos = aStr.Search( ",," ); aPos != -1; aPos = aStr.Search( ",," ) )
3140 aStr.Remove( aPos, 2 );
3142 TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
3144 while ( tmpStr != "" )
3146 tmpStr = aStr.Token( ",", i++ );
3147 int aPos = tmpStr.Search( '-' );
3151 if ( tmpStr.IsIntegerValue() )
3152 myIds.Add( tmpStr.IntegerValue() );
3158 TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
3159 TCollection_AsciiString aMinStr = tmpStr;
3161 while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
3162 while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
3164 if ( (!aMinStr.IsEmpty() && !aMinStr.IsIntegerValue()) ||
3165 (!aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue()) )
3168 myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
3169 myMax.Append( aMaxStr.IsEmpty() ? IntegerLast() : aMaxStr.IntegerValue() );
3176 //=======================================================================
3178 // Purpose : Get type of supported entities
3179 //=======================================================================
3180 SMDSAbs_ElementType RangeOfIds::GetType() const
3185 //=======================================================================
3187 // Purpose : Set type of supported entities
3188 //=======================================================================
3189 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
3194 //=======================================================================
3196 // Purpose : Verify whether entity satisfies to this rpedicate
3197 //=======================================================================
3198 bool RangeOfIds::IsSatisfy( long theId )
3203 if ( myType == SMDSAbs_Node )
3205 if ( myMesh->FindNode( theId ) == 0 )
3210 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
3211 if ( anElem == 0 || (myType != anElem->GetType() && myType != SMDSAbs_All ))
3215 if ( myIds.Contains( theId ) )
3218 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
3219 if ( theId >= myMin( i ) && theId <= myMax( i ) )
3227 Description : Base class for comparators
3229 Comparator::Comparator():
3233 Comparator::~Comparator()
3236 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
3239 myFunctor->SetMesh( theMesh );
3242 void Comparator::SetMargin( double theValue )
3244 myMargin = theValue;
3247 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
3249 myFunctor = theFunct;
3252 SMDSAbs_ElementType Comparator::GetType() const
3254 return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
3257 double Comparator::GetMargin()
3265 Description : Comparator "<"
3267 bool LessThan::IsSatisfy( long theId )
3269 return myFunctor && myFunctor->GetValue( theId ) < myMargin;
3275 Description : Comparator ">"
3277 bool MoreThan::IsSatisfy( long theId )
3279 return myFunctor && myFunctor->GetValue( theId ) > myMargin;
3285 Description : Comparator "="
3288 myToler(Precision::Confusion())
3291 bool EqualTo::IsSatisfy( long theId )
3293 return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
3296 void EqualTo::SetTolerance( double theToler )
3301 double EqualTo::GetTolerance()
3308 Description : Logical NOT predicate
3310 LogicalNOT::LogicalNOT()
3313 LogicalNOT::~LogicalNOT()
3316 bool LogicalNOT::IsSatisfy( long theId )
3318 return myPredicate && !myPredicate->IsSatisfy( theId );
3321 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
3324 myPredicate->SetMesh( theMesh );
3327 void LogicalNOT::SetPredicate( PredicatePtr thePred )
3329 myPredicate = thePred;
3332 SMDSAbs_ElementType LogicalNOT::GetType() const
3334 return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
3339 Class : LogicalBinary
3340 Description : Base class for binary logical predicate
3342 LogicalBinary::LogicalBinary()
3345 LogicalBinary::~LogicalBinary()
3348 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
3351 myPredicate1->SetMesh( theMesh );
3354 myPredicate2->SetMesh( theMesh );
3357 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
3359 myPredicate1 = thePredicate;
3362 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
3364 myPredicate2 = thePredicate;
3367 SMDSAbs_ElementType LogicalBinary::GetType() const
3369 if ( !myPredicate1 || !myPredicate2 )
3372 SMDSAbs_ElementType aType1 = myPredicate1->GetType();
3373 SMDSAbs_ElementType aType2 = myPredicate2->GetType();
3375 return aType1 == aType2 ? aType1 : SMDSAbs_All;
3381 Description : Logical AND
3383 bool LogicalAND::IsSatisfy( long theId )
3388 myPredicate1->IsSatisfy( theId ) &&
3389 myPredicate2->IsSatisfy( theId );
3395 Description : Logical OR
3397 bool LogicalOR::IsSatisfy( long theId )
3402 (myPredicate1->IsSatisfy( theId ) ||
3403 myPredicate2->IsSatisfy( theId ));
3417 void Filter::SetPredicate( PredicatePtr thePredicate )
3419 myPredicate = thePredicate;
3422 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3423 PredicatePtr thePredicate,
3424 TIdSequence& theSequence )
3426 theSequence.clear();
3428 if ( !theMesh || !thePredicate )
3431 thePredicate->SetMesh( theMesh );
3433 SMDS_ElemIteratorPtr elemIt = theMesh->elementsIterator( thePredicate->GetType() );
3435 while ( elemIt->more() ) {
3436 const SMDS_MeshElement* anElem = elemIt->next();
3437 long anId = anElem->GetID();
3438 if ( thePredicate->IsSatisfy( anId ) )
3439 theSequence.push_back( anId );
3444 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3445 Filter::TIdSequence& theSequence )
3447 GetElementsId(theMesh,myPredicate,theSequence);
3454 typedef std::set<SMDS_MeshFace*> TMapOfFacePtr;
3460 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
3461 SMDS_MeshNode* theNode2 )
3467 ManifoldPart::Link::~Link()
3473 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
3475 if ( myNode1 == theLink.myNode1 &&
3476 myNode2 == theLink.myNode2 )
3478 else if ( myNode1 == theLink.myNode2 &&
3479 myNode2 == theLink.myNode1 )
3485 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
3487 if(myNode1 < x.myNode1) return true;
3488 if(myNode1 == x.myNode1)
3489 if(myNode2 < x.myNode2) return true;
3493 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
3494 const ManifoldPart::Link& theLink2 )
3496 return theLink1.IsEqual( theLink2 );
3499 ManifoldPart::ManifoldPart()
3502 myAngToler = Precision::Angular();
3503 myIsOnlyManifold = true;
3506 ManifoldPart::~ManifoldPart()
3511 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
3517 SMDSAbs_ElementType ManifoldPart::GetType() const
3518 { return SMDSAbs_Face; }
3520 bool ManifoldPart::IsSatisfy( long theElementId )
3522 return myMapIds.Contains( theElementId );
3525 void ManifoldPart::SetAngleTolerance( const double theAngToler )
3526 { myAngToler = theAngToler; }
3528 double ManifoldPart::GetAngleTolerance() const
3529 { return myAngToler; }
3531 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
3532 { myIsOnlyManifold = theIsOnly; }
3534 void ManifoldPart::SetStartElem( const long theStartId )
3535 { myStartElemId = theStartId; }
3537 bool ManifoldPart::process()
3540 myMapBadGeomIds.Clear();
3542 myAllFacePtr.clear();
3543 myAllFacePtrIntDMap.clear();
3547 // collect all faces into own map
3548 SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3549 for (; anFaceItr->more(); )
3551 SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3552 myAllFacePtr.push_back( aFacePtr );
3553 myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3556 SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3560 // the map of non manifold links and bad geometry
3561 TMapOfLink aMapOfNonManifold;
3562 TColStd_MapOfInteger aMapOfTreated;
3564 // begin cycle on faces from start index and run on vector till the end
3565 // and from begin to start index to cover whole vector
3566 const int aStartIndx = myAllFacePtrIntDMap[aStartFace];
3567 bool isStartTreat = false;
3568 for ( int fi = aStartIndx; !isStartTreat || fi != aStartIndx ; fi++ )
3570 if ( fi == aStartIndx )
3571 isStartTreat = true;
3572 // as result next time when fi will be equal to aStartIndx
3574 SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3575 if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3578 aMapOfTreated.Add( aFacePtr->GetID() );
3579 TColStd_MapOfInteger aResFaces;
3580 if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3581 aMapOfNonManifold, aResFaces ) )
3583 TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3584 for ( ; anItr.More(); anItr.Next() )
3586 int aFaceId = anItr.Key();
3587 aMapOfTreated.Add( aFaceId );
3588 myMapIds.Add( aFaceId );
3591 if ( fi == ( myAllFacePtr.size() - 1 ) )
3593 } // end run on vector of faces
3594 return !myMapIds.IsEmpty();
3597 static void getLinks( const SMDS_MeshFace* theFace,
3598 ManifoldPart::TVectorOfLink& theLinks )
3600 int aNbNode = theFace->NbNodes();
3601 SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3603 SMDS_MeshNode* aNode = 0;
3604 for ( ; aNodeItr->more() && i <= aNbNode; )
3607 SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3611 SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3613 ManifoldPart::Link aLink( aN1, aN2 );
3614 theLinks.push_back( aLink );
3618 bool ManifoldPart::findConnected
3619 ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3620 SMDS_MeshFace* theStartFace,
3621 ManifoldPart::TMapOfLink& theNonManifold,
3622 TColStd_MapOfInteger& theResFaces )
3624 theResFaces.Clear();
3625 if ( !theAllFacePtrInt.size() )
3628 if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3630 myMapBadGeomIds.Add( theStartFace->GetID() );
3634 ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3635 ManifoldPart::TVectorOfLink aSeqOfBoundary;
3636 theResFaces.Add( theStartFace->GetID() );
3637 ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3639 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3640 aDMapLinkFace, theNonManifold, theStartFace );
3642 bool isDone = false;
3643 while ( !isDone && aMapOfBoundary.size() != 0 )
3645 bool isToReset = false;
3646 ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3647 for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3649 ManifoldPart::Link aLink = *pLink;
3650 if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3652 // each link could be treated only once
3653 aMapToSkip.insert( aLink );
3655 ManifoldPart::TVectorOfFacePtr aFaces;
3657 if ( myIsOnlyManifold &&
3658 (theNonManifold.find( aLink ) != theNonManifold.end()) )
3662 getFacesByLink( aLink, aFaces );
3663 // filter the element to keep only indicated elements
3664 ManifoldPart::TVectorOfFacePtr aFiltered;
3665 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3666 for ( ; pFace != aFaces.end(); ++pFace )
3668 SMDS_MeshFace* aFace = *pFace;
3669 if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3670 aFiltered.push_back( aFace );
3673 if ( aFaces.size() < 2 ) // no neihgbour faces
3675 else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3677 theNonManifold.insert( aLink );
3682 // compare normal with normals of neighbor element
3683 SMDS_MeshFace* aPrevFace = aDMapLinkFace[ aLink ];
3684 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3685 for ( ; pFace != aFaces.end(); ++pFace )
3687 SMDS_MeshFace* aNextFace = *pFace;
3688 if ( aPrevFace == aNextFace )
3690 int anNextFaceID = aNextFace->GetID();
3691 if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3692 // should not be with non manifold restriction. probably bad topology
3694 // check if face was treated and skipped
3695 if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3696 !isInPlane( aPrevFace, aNextFace ) )
3698 // add new element to connected and extend the boundaries.
3699 theResFaces.Add( anNextFaceID );
3700 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3701 aDMapLinkFace, theNonManifold, aNextFace );
3705 isDone = !isToReset;
3708 return !theResFaces.IsEmpty();
3711 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3712 const SMDS_MeshFace* theFace2 )
3714 gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3715 gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3716 if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3718 myMapBadGeomIds.Add( theFace2->GetID() );
3721 if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
3727 void ManifoldPart::expandBoundary
3728 ( ManifoldPart::TMapOfLink& theMapOfBoundary,
3729 ManifoldPart::TVectorOfLink& theSeqOfBoundary,
3730 ManifoldPart::TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
3731 ManifoldPart::TMapOfLink& theNonManifold,
3732 SMDS_MeshFace* theNextFace ) const
3734 ManifoldPart::TVectorOfLink aLinks;
3735 getLinks( theNextFace, aLinks );
3736 int aNbLink = (int)aLinks.size();
3737 for ( int i = 0; i < aNbLink; i++ )
3739 ManifoldPart::Link aLink = aLinks[ i ];
3740 if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3742 if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3744 if ( myIsOnlyManifold )
3746 // remove from boundary
3747 theMapOfBoundary.erase( aLink );
3748 ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3749 for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3751 ManifoldPart::Link aBoundLink = *pLink;
3752 if ( aBoundLink.IsEqual( aLink ) )
3754 theSeqOfBoundary.erase( pLink );
3762 theMapOfBoundary.insert( aLink );
3763 theSeqOfBoundary.push_back( aLink );
3764 theDMapLinkFacePtr[ aLink ] = theNextFace;
3769 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3770 ManifoldPart::TVectorOfFacePtr& theFaces ) const
3772 std::set<SMDS_MeshCell *> aSetOfFaces;
3773 // take all faces that shared first node
3774 SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3775 for ( ; anItr->more(); )
3777 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3780 aSetOfFaces.insert( aFace );
3782 // take all faces that shared second node
3783 anItr = theLink.myNode2->facesIterator();
3784 // find the common part of two sets
3785 for ( ; anItr->more(); )
3787 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3788 if ( aSetOfFaces.count( aFace ) )
3789 theFaces.push_back( aFace );
3798 ElementsOnSurface::ElementsOnSurface()
3801 myType = SMDSAbs_All;
3803 myToler = Precision::Confusion();
3804 myUseBoundaries = false;
3807 ElementsOnSurface::~ElementsOnSurface()
3811 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
3813 myMeshModifTracer.SetMesh( theMesh );
3814 if ( myMeshModifTracer.IsMeshModified())
3818 bool ElementsOnSurface::IsSatisfy( long theElementId )
3820 return myIds.Contains( theElementId );
3823 SMDSAbs_ElementType ElementsOnSurface::GetType() const
3826 void ElementsOnSurface::SetTolerance( const double theToler )
3828 if ( myToler != theToler )
3833 double ElementsOnSurface::GetTolerance() const
3836 void ElementsOnSurface::SetUseBoundaries( bool theUse )
3838 if ( myUseBoundaries != theUse ) {
3839 myUseBoundaries = theUse;
3840 SetSurface( mySurf, myType );
3844 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
3845 const SMDSAbs_ElementType theType )
3850 if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
3852 mySurf = TopoDS::Face( theShape );
3853 BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
3855 u1 = SA.FirstUParameter(),
3856 u2 = SA.LastUParameter(),
3857 v1 = SA.FirstVParameter(),
3858 v2 = SA.LastVParameter();
3859 Handle(Geom_Surface) surf = BRep_Tool::Surface( mySurf );
3860 myProjector.Init( surf, u1,u2, v1,v2 );
3864 void ElementsOnSurface::process()
3867 if ( mySurf.IsNull() )
3870 if ( !myMeshModifTracer.GetMesh() )
3873 myIds.ReSize( myMeshModifTracer.GetMesh()->GetMeshInfo().NbElements( myType ));
3875 SMDS_ElemIteratorPtr anIter = myMeshModifTracer.GetMesh()->elementsIterator( myType );
3876 for(; anIter->more(); )
3877 process( anIter->next() );
3880 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
3882 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3883 bool isSatisfy = true;
3884 for ( ; aNodeItr->more(); )
3886 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3887 if ( !isOnSurface( aNode ) )
3894 myIds.Add( theElemPtr->GetID() );
3897 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
3899 if ( mySurf.IsNull() )
3902 gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
3903 // double aToler2 = myToler * myToler;
3904 // if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
3906 // gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
3907 // if ( aPln.SquareDistance( aPnt ) > aToler2 )
3910 // else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
3912 // gp_Cylinder aCyl = Handle(Geom_CylindricalSurface)::DownCast(mySurf)->Cylinder();
3913 // double aRad = aCyl.Radius();
3914 // gp_Ax3 anAxis = aCyl.Position();
3915 // gp_XYZ aLoc = aCyl.Location().XYZ();
3916 // double aXDist = anAxis.XDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3917 // double aYDist = anAxis.YDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3918 // if ( fabs(aXDist*aXDist + aYDist*aYDist - aRad*aRad) > aToler2 )
3923 myProjector.Perform( aPnt );
3924 bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
3934 ElementsOnShape::ElementsOnShape()
3936 myType(SMDSAbs_All),
3937 myToler(Precision::Confusion()),
3938 myAllNodesFlag(false)
3942 ElementsOnShape::~ElementsOnShape()
3947 SMDSAbs_ElementType ElementsOnShape::GetType() const
3952 void ElementsOnShape::SetTolerance (const double theToler)
3954 if (myToler != theToler) {
3956 SetShape(myShape, myType);
3960 double ElementsOnShape::GetTolerance() const
3965 void ElementsOnShape::SetAllNodes (bool theAllNodes)
3967 myAllNodesFlag = theAllNodes;
3970 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
3975 void ElementsOnShape::SetShape (const TopoDS_Shape& theShape,
3976 const SMDSAbs_ElementType theType)
3980 if ( myShape.IsNull() ) return;
3982 TopTools_IndexedMapOfShape shapesMap;
3983 TopAbs_ShapeEnum shapeTypes[4] = { TopAbs_SOLID, TopAbs_FACE, TopAbs_EDGE, TopAbs_VERTEX };
3984 TopExp_Explorer sub;
3985 for ( int i = 0; i < 4; ++i )
3987 if ( shapesMap.IsEmpty() )
3988 for ( sub.Init( myShape, shapeTypes[i] ); sub.More(); sub.Next() )
3989 shapesMap.Add( sub.Current() );
3991 for ( sub.Init( myShape, shapeTypes[i], shapeTypes[i-1] ); sub.More(); sub.Next() )
3992 shapesMap.Add( sub.Current() );
3996 myClassifiers.resize( shapesMap.Extent() );
3997 for ( int i = 0; i < shapesMap.Extent(); ++i )
3998 myClassifiers[ i ] = new TClassifier( shapesMap( i+1 ), myToler );
4001 void ElementsOnShape::clearClassifiers()
4003 for ( size_t i = 0; i < myClassifiers.size(); ++i )
4004 delete myClassifiers[ i ];
4005 myClassifiers.clear();
4008 bool ElementsOnShape::IsSatisfy (long elemId)
4010 const SMDS_MeshElement* elem =
4011 ( myType == SMDSAbs_Node ? myMesh->FindNode( elemId ) : myMesh->FindElement( elemId ));
4012 if ( !elem || myClassifiers.empty() )
4015 for ( size_t i = 0; i < myClassifiers.size(); ++i )
4017 SMDS_ElemIteratorPtr aNodeItr = elem->nodesIterator();
4018 bool isSatisfy = myAllNodesFlag;
4020 gp_XYZ centerXYZ (0, 0, 0);
4022 while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
4024 SMESH_TNodeXYZ aPnt ( aNodeItr->next() );
4026 isSatisfy = ! myClassifiers[i]->IsOut( aPnt );
4029 // Check the center point for volumes MantisBug 0020168
4032 myClassifiers[i]->ShapeType() == TopAbs_SOLID)
4034 centerXYZ /= elem->NbNodes();
4035 isSatisfy = ! myClassifiers[i]->IsOut( centerXYZ );
4044 TopAbs_ShapeEnum ElementsOnShape::TClassifier::ShapeType() const
4046 return myShape.ShapeType();
4049 bool ElementsOnShape::TClassifier::IsOut(const gp_Pnt& p)
4051 return (this->*myIsOutFun)( p );
4054 void ElementsOnShape::TClassifier::Init (const TopoDS_Shape& theShape, double theTol)
4058 switch ( myShape.ShapeType() )
4060 case TopAbs_SOLID: {
4061 mySolidClfr.Load(theShape);
4062 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfSolid;
4066 Standard_Real u1,u2,v1,v2;
4067 Handle(Geom_Surface) surf = BRep_Tool::Surface( TopoDS::Face( theShape ));
4068 surf->Bounds( u1,u2,v1,v2 );
4069 myProjFace.Init(surf, u1,u2, v1,v2, myTol );
4070 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfFace;
4074 Standard_Real u1, u2;
4075 Handle(Geom_Curve) curve = BRep_Tool::Curve( TopoDS::Edge(theShape), u1, u2);
4076 myProjEdge.Init(curve, u1, u2);
4077 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfEdge;
4080 case TopAbs_VERTEX:{
4081 myVertexXYZ = BRep_Tool::Pnt( TopoDS::Vertex( theShape ) );
4082 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfVertex;
4086 throw SALOME_Exception("Programmer error in usage of ElementsOnShape::TClassifier");
4090 bool ElementsOnShape::TClassifier::isOutOfSolid (const gp_Pnt& p)
4092 mySolidClfr.Perform( p, myTol );
4093 return ( mySolidClfr.State() != TopAbs_IN && mySolidClfr.State() != TopAbs_ON );
4096 bool ElementsOnShape::TClassifier::isOutOfFace (const gp_Pnt& p)
4098 myProjFace.Perform( p );
4099 if ( myProjFace.IsDone() && myProjFace.LowerDistance() <= myTol )
4101 // check relatively to the face
4102 Quantity_Parameter u, v;
4103 myProjFace.LowerDistanceParameters(u, v);
4104 gp_Pnt2d aProjPnt (u, v);
4105 BRepClass_FaceClassifier aClsf ( TopoDS::Face( myShape ), aProjPnt, myTol );
4106 if ( aClsf.State() == TopAbs_IN || aClsf.State() == TopAbs_ON )
4112 bool ElementsOnShape::TClassifier::isOutOfEdge (const gp_Pnt& p)
4114 myProjEdge.Perform( p );
4115 return ! ( myProjEdge.NbPoints() > 0 && myProjEdge.LowerDistance() <= myTol );
4118 bool ElementsOnShape::TClassifier::isOutOfVertex(const gp_Pnt& p)
4120 return ( myVertexXYZ.Distance( p ) > myTol );
4124 TSequenceOfXYZ::TSequenceOfXYZ()
4127 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n)
4130 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t)
4133 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray)
4136 template <class InputIterator>
4137 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd)
4140 TSequenceOfXYZ::~TSequenceOfXYZ()
4143 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
4145 myArray = theSequenceOfXYZ.myArray;
4149 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
4151 return myArray[n-1];
4154 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
4156 return myArray[n-1];
4159 void TSequenceOfXYZ::clear()
4164 void TSequenceOfXYZ::reserve(size_type n)
4169 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
4171 myArray.push_back(v);
4174 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
4176 return myArray.size();
4179 TMeshModifTracer::TMeshModifTracer():
4180 myMeshModifTime(0), myMesh(0)
4183 void TMeshModifTracer::SetMesh( const SMDS_Mesh* theMesh )
4185 if ( theMesh != myMesh )
4186 myMeshModifTime = 0;
4189 bool TMeshModifTracer::IsMeshModified()
4191 bool modified = false;
4194 modified = ( myMeshModifTime != myMesh->GetMTime() );
4195 myMeshModifTime = myMesh->GetMTime();